The Elk Code
dforce.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2012 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 dforce(dyn)
7 use modmain
8 use modphonon
9 use modmpi
10 use modomp
11 implicit none
12 ! arguments
13 complex(8), intent(out) :: dyn(3,natmtot)
14 ! local variables
15 integer ik,is,ias,jas,ip
16 integer nr,nri,nthd
17 complex(8) z1
18 ! automatic arrays
19 complex(8) dynibs(3,natmtot)
20 ! allocatable arrays
21 complex(8), allocatable :: zrhomt(:,:),zrhoir(:)
22 complex(8), allocatable :: grhomt(:,:,:),grhoir(:,:)
23 complex(8), allocatable :: zvclmt(:,:),zvclir(:)
24 complex(8), allocatable :: gvclmt(:,:,:),gvclir(:,:)
25 complex(8), allocatable :: zfmt(:),gvcln(:,:),gzfmt(:,:)
26 ! external functions
27 complex(8), external :: zfmtinp
28 allocate(zrhomt(npmtmax,natmtot),zrhoir(ngtot))
29 allocate(grhomt(npmtmax,natmtot,3),grhoir(ngtot,3))
30 allocate(zvclmt(npmtmax,natmtot+1),zvclir(ngtot))
31 allocate(gvclmt(npmtmax,natmtot,3),gvclir(ngtot,3))
32 allocate(zfmt(npmtmax),gvcln(npmtmax,3),gzfmt(npmtmax,3))
33 ! make complex copy of the density
34 do ias=1,natmtot
35  is=idxis(ias)
36  call rtozfmt(nrmt(is),nrmti(is),rhomt(:,ias),zrhomt(:,ias))
37 end do
38 zrhoir(1:ngtot)=rhoir(1:ngtot)
39 ! compute the gradient of the density
40 call gradzf(zrhomt,zrhoir,grhomt,grhoir)
41 !--------------------------------------------------------------!
42 ! Hellmann-Feynman force derivative for displaced atom !
43 !--------------------------------------------------------------!
44 ! calculate the gradient of the nuclear potential
45 call gradzvcln(isph,gvcln)
46 do ip=1,3
47 ! compute the q-dependent nuclear Coulomb potential derivative
48  zvclmt(1:npmtmax,1:natmtot)=0.d0
49  zvclmt(1:npmt(isph),iasph)=gvcln(1:npmt(isph),ip)
50 ! zero the interstitial density
51  zrhoir(1:ngtot)=0.d0
52  if (ip == ipph) then; jas=iasph; else; jas=0; endif
54  ngvec,jlgqrmt,ylmgq,sfacgq,zrhoir,npmtmax,zvclmt,zvclir)
55 ! multiply with density derivative and integrate
56  z1=sum(cfunir(1:ngtot)*conjg(zvclir(1:ngtot))*drhoir(1:ngtot))
57  z1=z1*(omega/ngtot)
58  do ias=1,natmtot
59  is=idxis(ias)
60  z1=z1+zfmtinp(nrmt(is),nrmti(is),wr2mt(:,is),zvclmt(:,ias),drhomt(:,ias))
61  end do
62  dyn(ip,iasph)=-z1
63 ! nuclear-nuclear term
64  if (ip == ipph) then; jas=natmtot+1; else; jas=iasph; end if
65  call gradzfmt(nrmt(isph),nrmti(isph),rlmt(:,-1,isph),wcrmt(:,:,isph), &
66  zvclmt(:,jas),npmtmax,gzfmt)
67  z1=spzn(isph)*gzfmt(1,ipph)*y00
68  dyn(ip,iasph)=dyn(ip,iasph)-z1
69 end do
70 ! compute the lattice-periodic nuclear Coulomb potential derivative
71 zvclmt(1:npmtmax,1:natmtot)=0.d0
72 zvclmt(1:npmt(isph),iasph)=gvcln(1:npmt(isph),ipph)
74  ngvec,jlgrmt,ylmg,sfacg,zrhoir,npmtmax,zvclmt,zvclir)
75 do ip=1,3
76 ! multiply with density gradient and integrate
77  z1=sum(cfunir(1:ngtot)*zvclir(1:ngtot)*grhoir(1:ngtot,ip))
78  z1=z1*(omega/ngtot)
79  do ias=1,natmtot
80  is=idxis(ias)
81  z1=z1+zfmtinp(nrmt(is),nrmti(is),wr2mt(:,is),grhomt(:,ias,ip),zvclmt(:,ias))
82  end do
83  dyn(ip,iasph)=dyn(ip,iasph)-z1
84 end do
85 ! nuclear-nuclear term
86 do ip=1,3
87  do ias=1,natmtot
88  is=idxis(ias)
89  if (ias == iasph) then; jas=natmtot+1; else; jas=ias; end if
90  call gradzfmt(nrmt(is),nrmti(is),rlmt(:,-1,is),wcrmt(:,:,is),zvclmt(:,jas),&
91  npmtmax,gzfmt)
92  z1=spzn(is)*gzfmt(1,ip)*y00
93  dyn(ip,iasph)=dyn(ip,iasph)+z1
94  end do
95 end do
96 !-------------------------------------------------------------------!
97 ! Hellmann-Feynman force derivative for non-displaced atoms !
98 !-------------------------------------------------------------------!
99 do ias=1,natmtot
100  if (ias == iasph) cycle
101  is=idxis(ias)
102  nr=nrmt(is)
103  nri=nrmti(is)
104 ! compute the gradient of the Coulomb potential derivative at the nucleus
105  call gradzfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),dvclmt(:,ias),npmtmax,gzfmt)
106  do ip=1,3
107  dyn(ip,ias)=spzn(is)*gzfmt(1,ip)*y00
108  end do
109 end do
110 !--------------------------------------------!
111 ! IBS correction to force derivative !
112 !--------------------------------------------!
113 dynibs(:,:)=0.d0
114 ! k-point dependent part
115 call holdthd(nkptnr/np_mpi,nthd)
116 !$OMP PARALLEL DO DEFAULT(SHARED) &
117 !$OMP REDUCTION(+:dynibs) &
118 !$OMP SCHEDULE(DYNAMIC) &
119 !$OMP NUM_THREADS(nthd)
120 do ik=1,nkptnr
121 ! distribute among MPI processes
122  if (mod(ik-1,np_mpi) /= lp_mpi) cycle
123  call dforcek(ik,dynibs)
124 end do
125 !$OMP END PARALLEL DO
126 call freethd(nthd)
127 ! add IBS force derivatives from each process and redistribute
128 if (np_mpi > 1) then
129  call mpi_allreduce(mpi_in_place,dynibs,3*natmtot,mpi_double_complex,mpi_sum, &
130  mpicom,ierror)
131 end if
132 ! k-point independent part
133 do ias=1,natmtot
134  is=idxis(ias)
135  nr=nrmt(is)
136  nri=nrmti(is)
137  do ip=1,3
138  z1=zfmtinp(nr,nri,wr2mt(:,is),grhomt(:,ias,ip),dvsmt(:,ias))
139  dynibs(ip,ias)=dynibs(ip,ias)-z1
140  end do
141 ! compute the gradient of the density derivative
142  if (ias == iasph) then; jas=natmtot+1; else; jas=ias; end if
143  call gradzfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),drhomt(:,jas),npmtmax,gzfmt)
144 ! convert Kohn-Sham potential to complex spherical harmonics
145  call rtozfmt(nr,nri,vsmt(:,ias),zfmt)
146  do ip=1,3
147  z1=zfmtinp(nr,nri,wr2mt(:,is),zfmt,gzfmt(:,ip))
148  dynibs(ip,ias)=dynibs(ip,ias)-z1
149  end do
150 end do
151 ! add the IBS force derivatives to the total
152 dyn(:,:)=dyn(:,:)+dynibs(:,:)
153 deallocate(zrhomt,zrhoir,grhomt,grhoir)
154 deallocate(zvclmt,zvclir,gvclmt,gvclir,zfmt,gzfmt)
155 end subroutine
156 
complex(8), dimension(:,:), allocatable sfacg
Definition: modmain.f90:430
subroutine dforce(dyn)
Definition: dforce.f90:7
integer, dimension(3) ngridg
Definition: modmain.f90:386
real(8), dimension(:,:,:), allocatable jlgqrmt
Definition: modphonon.f90:68
integer ngtot
Definition: modmain.f90:390
complex(8), dimension(:,:), pointer, contiguous drhomt
Definition: modphonon.f90:100
real(8), dimension(:), pointer, contiguous rhoir
Definition: modmain.f90:614
complex(8), dimension(:,:), allocatable dvclmt
Definition: modphonon.f90:104
subroutine zpotcoul(iash, nrmt_, nrmti_, npmt_, ld1, rl, ngridg_, igfft_, ngp, gpc, gclgp, ld2, jlgprmt, ylmgp, sfacgp, zrhoir, ld3, zvclmt, zvclir)
Definition: zpotcoul.f90:11
real(8), dimension(:,:,:), allocatable rlmt
Definition: modmain.f90:179
integer, dimension(maxspecies) npmt
Definition: modmain.f90:213
real(8) omega
Definition: modmain.f90:20
Definition: modomp.f90:6
subroutine gradzfmt(nr, nri, ri, wcr, zfmt, ld, gzfmt)
Definition: gradzfmt.f90:10
real(8), dimension(:,:), pointer, contiguous rhomt
Definition: modmain.f90:614
integer iasph
Definition: modphonon.f90:15
pure subroutine rtozfmt(nr, nri, rfmt, zfmt)
Definition: rtozfmt.f90:7
integer nkptnr
Definition: modmain.f90:463
integer np_mpi
Definition: modmpi.f90:13
complex(8), dimension(:,:), allocatable ylmg
Definition: modmain.f90:428
integer, dimension(:), allocatable igfft
Definition: modmain.f90:406
subroutine gradzvcln(is, gzfmt)
Definition: gradzvcln.f90:7
complex(8), dimension(:,:), pointer, contiguous dvsmt
Definition: modphonon.f90:110
integer ipph
Definition: modphonon.f90:15
integer ngvec
Definition: modmain.f90:396
real(8), dimension(:), allocatable cfunir
Definition: modmain.f90:436
integer isph
Definition: modphonon.f90:15
subroutine dforcek(ik, dynibs)
Definition: dforcek.f90:7
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
subroutine gradzf(zfmt, zfir, gzfmt, gzfir)
Definition: gradzf.f90:7
Definition: modmpi.f90:6
real(8), dimension(:,:,:), allocatable jlgrmt
Definition: modmain.f90:426
real(8), dimension(:), allocatable gclg
Definition: modmain.f90:424
complex(8), dimension(:), pointer, contiguous drhoir
Definition: modphonon.f90:100
real(8), dimension(:), allocatable gc
Definition: modmain.f90:422
complex(8), dimension(:,:), allocatable sfacgq
Definition: modphonon.f90:72
integer lp_mpi
Definition: modmpi.f90:15
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
real(8), dimension(maxspecies) spzn
Definition: modmain.f90:80
integer npmtmax
Definition: modmain.f90:216
real(8), dimension(:,:,:), allocatable wcrmt
Definition: modmain.f90:187
integer nrmtmax
Definition: modmain.f90:152
real(8), dimension(:), allocatable gqc
Definition: modphonon.f90:64
real(8), dimension(:,:), pointer, contiguous vsmt
Definition: modmain.f90:649
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
real(8), parameter y00
Definition: modmain.f90:1236
integer mpicom
Definition: modmpi.f90:11
real(8), dimension(:), allocatable gclgq
Definition: modphonon.f90:66
complex(8), dimension(:,:), allocatable ylmgq
Definition: modphonon.f90:70
real(8), dimension(:,:), allocatable wr2mt
Definition: modmain.f90:183
integer ierror
Definition: modmpi.f90:19
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150