The Elk Code
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:
9 subroutine rdmdtsdn(dedn)
10 ! !USES:
11 use modmain
12 use 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
23 implicit none
24 ! arguments
25 real(8), intent(inout) :: dedn(nstsv,nkpt)
26 ! local variables
27 integer ik,ist
28 real(8) t1
29 do 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
35 end do
36 end subroutine
37 !EOC
38 
subroutine rdmdtsdn(dedn)
Definition: rdmdtsdn.f90:10
real(8), parameter kboltz
Definition: modmain.f90:1263
real(8) epsocc
Definition: modmain.f90:903
real(8) rdmtemp
Definition: modrdm.f90:31
Definition: modrdm.f90:6
real(8) occmax
Definition: modmain.f90:901
real(8), dimension(:,:), allocatable occsv
Definition: modmain.f90:905