The Elk Code
 
Loading...
Searching...
No Matches
rdmdtsdn.f90
Go to the documentation of this file.
1
2! Copyright (C) 2008 T. Baldsiefen, S. Sharma, J. K. Dewhurst and
3! E. K. U. Gross. This file is distributed under the terms of the GNU Lesser
4! General Public License. See the file COPYING for license details.
5
6!BOP
7! !ROUTINE: rdmdtsdn
8! !INTERFACE:
9subroutine rdmdtsdn(dedn)
10! !USES:
11use modmain
12use modrdm
13! !INPUT/OUTPUT PARAMETERS:
14! dedn : energy derivative (inout,real(nstsv,nkpt))
15! !DESCRIPTION:
16! Calculates the derivative of the entropic contribution to the free energy
17! with respect to the occupation numbers and adds it to the total.
18!
19! !REVISION HISTORY:
20! Created 2008 (Baldsiefen)
21!EOP
22!BOC
23implicit none
24! arguments
25real(8), intent(inout) :: dedn(nstsv,nkpt)
26! local variables
27integer ik,ist
28real(8) t1
29do ik=1,nkpt
30 do ist=1,nstsv
31 t1=max(occsv(ist,ik),epsocc)
32 t1=min(t1,occmax-epsocc)
33 dedn(ist,ik)=dedn(ist,ik)-rdmtemp*kboltz*log(t1/(occmax-t1))
34 end do
35end do
36end subroutine
37!EOC
38
real(8), parameter kboltz
Definition modmain.f90:1262
real(8) epsocc
Definition modmain.f90:900
real(8) occmax
Definition modmain.f90:898
real(8), dimension(:,:), allocatable occsv
Definition modmain.f90:902
real(8) rdmtemp
Definition modrdm.f90:31
subroutine rdmdtsdn(dedn)
Definition rdmdtsdn.f90:10