The Elk Code
dpotks.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2011 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 subroutine dpotks
7 use modmain
8 use modphonon
9 use modomp
10 implicit none
11 ! local variables
12 integer is,ias,np,nthd
13 ! convert density derivative to spherical coordinates
14 call holdthd(natmtot,nthd)
15 !$OMP PARALLEL DO DEFAULT(SHARED) &
16 !$OMP PRIVATE(is) &
17 !$OMP SCHEDULE(DYNAMIC) &
18 !$OMP NUM_THREADS(nthd)
19 do ias=1,natmtot
20  is=idxis(ias)
21  call zbshtip(nrmt(is),nrmti(is),drhomt(:,ias))
22 end do
23 !$OMP END PARALLEL DO
24 call freethd(nthd)
25 ! compute the exchange-correlation potential derivative
26 call dpotxc
27 ! convert density derivative to spherical harmonics
28 call holdthd(natmtot,nthd)
29 !$OMP PARALLEL DO DEFAULT(SHARED) &
30 !$OMP PRIVATE(is) &
31 !$OMP SCHEDULE(DYNAMIC) &
32 !$OMP NUM_THREADS(nthd)
33 do ias=1,natmtot
34  is=idxis(ias)
35  call zfshtip(nrmt(is),nrmti(is),drhomt(:,ias))
36 end do
37 !$OMP END PARALLEL DO
38 call freethd(nthd)
39 ! generate the Coulomb potential derivative
40 call dpotcoul
41 ! add to the Kohn-Sham potential derivative
42 do ias=1,natmtot
43  is=idxis(ias)
44  np=npmt(is)
45  dvsmt(1:np,ias)=dvsmt(1:np,ias)+dvclmt(1:np,ias)
46 end do
48 ! remove the gradient part of the potential derivative for displaced muffin-tin
49 np=npmt(isph)
50 dvsmt(1:np,iasph)=dvsmt(1:np,iasph)+gvsmt(1:np)
51 end subroutine
52 
subroutine dpotks
Definition: dpotks.f90:7
integer ngtot
Definition: modmain.f90:390
complex(8), dimension(:,:), pointer, contiguous drhomt
Definition: modphonon.f90:100
complex(8), dimension(:,:), allocatable dvclmt
Definition: modphonon.f90:104
complex(8), dimension(:), allocatable dvsir
Definition: modphonon.f90:112
integer, dimension(maxspecies) npmt
Definition: modmain.f90:213
subroutine dpotcoul
Definition: dpotcoul.f90:7
Definition: modomp.f90:6
integer iasph
Definition: modphonon.f90:15
complex(8), dimension(:), allocatable gvsmt
Definition: modphonon.f90:106
subroutine zbshtip(nr, nri, zfmt)
Definition: zbshtip.f90:7
subroutine zfshtip(nr, nri, zfmt)
Definition: zfshtip.f90:7
complex(8), dimension(:,:), pointer, contiguous dvsmt
Definition: modphonon.f90:110
integer isph
Definition: modphonon.f90:15
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer natmtot
Definition: modmain.f90:40
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
subroutine dpotxc
Definition: dpotxc.f90:7
complex(8), dimension(:), allocatable dvclir
Definition: modphonon.f90:104
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150