File: rdmdtsdn.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (40 lines) | stat: -rw-r--r-- 923 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40

! Copyright (C) 2008 T. Baldsiefen, S. Sharma, J. K. Dewhurst and E. K. U. Gross.
! This file is distributed under the terms of the GNU Lesser General Public
! License. See the file COPYING for license details.


!BOP
! !ROUTINE: rdmdtsdn
! !INTERFACE:
subroutine rdmdtsdn(dedn)
! !USES:
use modmain
use modrdm
! !INPUT/OUTPUT PARAMETERS:
!   dedn : energy derivative (inout,real(nstsv,nkpt))
! !DESCRIPTION:
!   Calculates the derivative of the entropic contribution to the free energy
!   w.r.t. occupation numbers and adds it to the total.
!
! !REVISION HISTORY:
!   Created 2008 (Baldsiefen)
!EOP
!BOC
implicit none
! arguments
real(8), intent(inout) :: dedn(nstsv,nkpt)
! local variables
integer  ik,ist
real(8)  t1
do ik=1,nkpt
  do ist=1,nstsv
    t1=max(occsv(ist,ik),epsocc)
    t1=min(t1,occmax-epsocc)
    dedn(ist,ik)=dedn(ist,ik)-rdmtemp*kboltz*log(t1/(occmax-t1))
  end do
end do
return
end subroutine
!EOC