The Elk Code
rdmdedn.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2007-2008 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: rdmdedn
8 ! !INTERFACE:
9 subroutine rdmdedn(dedn)
10 ! !USES:
11 use modmain
12 use modrdm
13 use modomp
14 ! !INPUT/OUTPUT PARAMETERS:
15 ! dedn : free energy derivative (out,real(nstsv,nkpt))
16 ! !DESCRIPTION:
17 ! Calculates the negative of the derivative of total free energy w.r.t.
18 ! occupation numbers.
19 !
20 ! !REVISION HISTORY:
21 ! Created 2008 (Sharma)
22 !EOP
23 !BOC
24 implicit none
25 ! arguments
26 real(8), intent(out) :: dedn(nstsv,nkpt)
27 ! local variables
28 integer ik,ist,nthd
29 ! allocatable arrays
30 complex(8), allocatable :: evecsv(:,:),c(:,:)
31 call holdthd(nkpt,nthd)
32 !$OMP PARALLEL DEFAULT(SHARED) &
33 !$OMP PRIVATE(evecsv,c,ist) &
34 !$OMP NUM_THREADS(nthd)
35 allocate(evecsv(nstsv,nstsv),c(nstsv,nstsv))
36 !$OMP DO SCHEDULE(DYNAMIC)
37 do ik=1,nkpt
38 ! get eigenvectors from file
39  call getevecsv(filext,ik,vkl(:,ik),evecsv)
40 ! kinetic and Coulomb potential contribution
41  call zgemm('C','N',nstsv,nstsv,nstsv,zone,evecsv,nstsv,dkdc(:,:,ik),nstsv, &
42  zzero,c,nstsv)
43  do ist=1,nstsv
44  dedn(ist,ik)=-(dble(c(ist,ist))+dble(vclmat(ist,ist,ik)))
45  end do
46 end do
47 !$OMP END DO
48 deallocate(evecsv,c)
49 !$OMP END PARALLEL
50 call freethd(nthd)
51 ! add exchange correlation contribution
52 call rdmdexcdn(dedn)
53 ! add entropic contribution if needed
54 if (rdmtemp > 0.d0) call rdmdtsdn(dedn)
55 end subroutine
56 !EOC
57 
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition: getevecsv.f90:7
character(256) filext
Definition: modmain.f90:1301
subroutine rdmdtsdn(dedn)
Definition: rdmdtsdn.f90:10
Definition: modomp.f90:6
complex(8), parameter zone
Definition: modmain.f90:1240
subroutine rdmdexcdn(dedn)
Definition: rdmdexcdn.f90:10
subroutine rdmdedn(dedn)
Definition: rdmdedn.f90:10
complex(8), dimension(:,:,:), allocatable dkdc
Definition: modrdm.f90:15
real(8) rdmtemp
Definition: modrdm.f90:31
Definition: modrdm.f90:6
complex(8), parameter zzero
Definition: modmain.f90:1240
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
complex(8), dimension(:,:,:), allocatable vclmat
Definition: modrdm.f90:13