The Elk Code
 
Loading...
Searching...
No Matches
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:
9subroutine rdmdedc(dedc)
10! !USES:
11use modmain
12use modrdm
13use 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
24implicit none
25! arguments
26complex(8), intent(out) :: dedc(nstsv,nstsv,nkpt)
27! local variables
28integer ik,ist,nthd
29! allocatable arrays
30complex(8), allocatable :: evecsv(:,:),c(:,:)
31call holdthd(nkpt,nthd)
32!$OMP PARALLEL DEFAULT(SHARED) &
33!$OMP PRIVATE(evecsv,c,ist) &
34!$OMP NUM_THREADS(nthd)
35allocate(evecsv(nstsv,nstsv),c(nstsv,nstsv))
36!$OMP DO SCHEDULE(DYNAMIC)
37do 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
46end do
47!$OMP END DO
48deallocate(evecsv,c)
49!$OMP END PARALLEL
50call freethd(nthd)
51! exchange-correlation contribution
52call rdmdexcdc(dedc)
53end subroutine
54!EOC
55
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition getevecsv.f90:7
complex(8), parameter zzero
Definition modmain.f90:1238
character(256) filext
Definition modmain.f90:1300
complex(8), parameter zone
Definition modmain.f90:1238
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
real(8), dimension(:,:), allocatable occsv
Definition modmain.f90:902
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
complex(8), dimension(:,:,:), allocatable vclmat
Definition modrdm.f90:13
complex(8), dimension(:,:,:), allocatable dkdc
Definition modrdm.f90:15
subroutine rdmdedc(dedc)
Definition rdmdedc.f90:10
subroutine rdmdexcdc(dedc)
Definition rdmdexcdc.f90:10