The Elk Code
rdmdedc.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2007 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: rdmdedc
8 ! !INTERFACE:
9 subroutine rdmdedc(dedc)
10 ! !USES:
11 use modmain
12 use modrdm
13 use modomp
14 ! !INPUT/OUTPUT PARAMETERS:
15 ! dedc : energy derivative (out,complex(nstsv,nstsv,nkpt))
16 ! !DESCRIPTION:
17 ! Calculates the derivative of the total energy w.r.t. the second-variational
18 ! coefficients {\tt evecsv}.
19 !
20 ! !REVISION HISTORY:
21 ! Created 2008 (Sharma)
22 !EOP
23 !BOC
24 implicit none
25 ! arguments
26 complex(8), intent(out) :: dedc(nstsv,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 the eigenvectors from file
39  call getevecsv(filext,ik,vkl(:,ik),evecsv)
40 ! kinetic and Coulomb potential contribution
41  call zgemm('N','N',nstsv,nstsv,nstsv,zone,evecsv,nstsv,vclmat(:,:,ik),nstsv, &
42  zzero,c,nstsv)
43  do ist=1,nstsv
44  dedc(:,ist,ik)=occsv(ist,ik)*(dkdc(:,ist,ik)+c(:,ist))
45  end do
46 end do
47 !$OMP END DO
48 deallocate(evecsv,c)
49 !$OMP END PARALLEL
50 call freethd(nthd)
51 ! exchange-correlation contribution
52 call rdmdexcdc(dedc)
53 end subroutine
54 !EOC
55 
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition: getevecsv.f90:7
character(256) filext
Definition: modmain.f90:1301
Definition: modomp.f90:6
complex(8), parameter zone
Definition: modmain.f90:1240
complex(8), dimension(:,:,:), allocatable dkdc
Definition: modrdm.f90:15
Definition: modrdm.f90:6
subroutine rdmdedc(dedc)
Definition: rdmdedc.f90:10
real(8), dimension(:,:), allocatable occsv
Definition: modmain.f90:905
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
subroutine rdmdexcdc(dedc)
Definition: rdmdexcdc.f90:10