The Elk Code
 
Loading...
Searching...
No Matches
rdmentropy.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 General
4! Public License. See the file COPYING for license details.
5
6!BOP
7! !ROUTINE: rdmentropy
8! !INTERFACE:
9subroutine rdmentropy
10! !USES:
11use modmain
12use modrdm
13! !DESCRIPTION:
14! Calculates RDMFT entropy $S=-\sum_i n_i\log(n_i/n_{\rm max})
15! +(n_{\rm max}-n_i)\log(1-n_i/n_{\rm max})$, where $n_{\rm max}$ is the
16! maximum allowed occupancy (1 or 2).
17!
18! !REVISION HISTORY:
19! Created 2008 (Baldsiefen)
20!EOP
21!BOC
22implicit none
23! local variables
24integer ik,ist
25real(8) t1
26rdmentrpy=0.d0
27do ik=1,nkpt
28 do ist=1,nstsv
29 t1=max(occsv(ist,ik),epsocc)
30 t1=min(t1,occmax-epsocc)
31 rdmentrpy=rdmentrpy-wkpt(ik)*(t1*log(t1/occmax) &
32 +(occmax-t1)*log(1.d0-t1/occmax))
33 end do
34end do
36end subroutine
37!EOC
38
real(8), dimension(:), allocatable wkpt
Definition modmain.f90:475
integer nkpt
Definition modmain.f90:461
real(8), parameter kboltz
Definition modmain.f90:1262
integer nstsv
Definition modmain.f90:886
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) rdmentrpy
Definition modrdm.f90:33
subroutine rdmentropy