12 integer,
intent(in) :: ik
13 complex(8),
intent(inout) :: dynibs(3,natmtot)
15 integer ispn0,ispn1,ispn,jspn
17 integer is,ias,ist,jst,jk
18 integer iv(3),ig,i,j,l
20 complex(8) z1,z2,dt1,dz1,dz2
22 real(8) evalfv(nstfv,nspnfv)
23 complex(8) vh(nmatmax),vo(nmatmax),dvh(nmatmax),dvo(nmatmax)
24 complex(8) ffv(nstfv,nstfv),dffv(nstfv,nstfv),y(nstfv),dy(nstfv)
26 integer,
allocatable :: ijg(:,:),ijgq(:,:)
27 real(8),
allocatable :: dp(:,:),dpq(:,:)
28 complex(8),
allocatable :: apwalm(:,:,:,:),apwalmq(:,:,:,:),dapwalm(:,:,:)
29 complex(8),
allocatable :: evecfv(:,:,:),devecfv(:,:,:)
30 complex(8),
allocatable :: evecsv(:,:),devecsv(:,:)
31 complex(8),
allocatable :: h(:,:),o(:,:),dlh(:,:),dlo(:,:)
32 complex(8),
allocatable :: hq(:,:),oq(:,:),dh(:,:),od(:,:)
33 complex(8),
allocatable :: dlhq(:,:),dloq(:,:),ddlh(:,:),ddlo(:,:)
35 complex(8),
external :: zdotc
37 allocate(ijg(nmatmax,nmatmax),ijgq(nmatmax,nmatmax))
38 allocate(dp(nmatmax,nmatmax),dpq(nmatmax,nmatmax))
42 allocate(evecfv(nmatmax,nstfv,nspnfv))
43 allocate(devecfv(nmatmax,nstfv,nspnfv))
44 allocate(h(nmatmax,nmatmax),o(nmatmax,nmatmax))
45 allocate(dlh(nmatmax,nmatmax),dlo(nmatmax,nmatmax))
46 allocate(hq(nmatmax,nmatmax),oq(nmatmax,nmatmax))
47 allocate(dh(nmatmax,nmatmax),od(nmatmax,nmatmax))
48 allocate(dlhq(nmatmax,nmatmax),dloq(nmatmax,nmatmax))
49 allocate(ddlh(nmatmax,nmatmax),ddlo(nmatmax,nmatmax))
65 ispn0=jspn; ispn1=jspn
77 ijg(i,j)=
ivgig(iv(1),iv(2),iv(3))
78 dp(i,j)=0.5d0*dot_product(
vgkc(1:3,i,jspn,ik),
vgkc(1:3,j,jspn,ik))
85 ijgq(i,j)=
ivgig(iv(1),iv(2),iv(3))
86 dpq(i,j)=0.5d0*dot_product(
vgkqc(1:3,i,jspn,ik),
vgkc(1:3,j,jspn,ik))
90 call match(n,
vgkc(:,:,jspn,ik),
gkc(:,jspn,ik),
sfacgk(:,:,jspn,ik),apwalm)
91 call match(nq,
vgkqc(:,:,jspn,ik),
gkqc(:,jspn,ik),
sfacgkq(:,:,jspn,ik), &
102 call hmlaa(.false.,is,ias,n,apwalm(:,:,:,ias),nmatmax,h)
103 call hmlalo(is,ias,n,apwalm(:,:,:,ias),nmatmax,h)
107 call olpaa(.false.,is,n,apwalm(:,:,:,ias),nmatmax,o)
108 call olpalo(is,ias,n,apwalm(:,:,:,ias),nmatmax,o)
110 call hmlaaq(is,ias,n,nq,apwalm(:,:,:,ias),apwalmq(:,:,:,ias),nmatmax,hq)
111 call hmlaloq(is,ias,n,nq,apwalm(:,:,:,ias),apwalmq(:,:,:,ias),nmatmax,hq)
113 call olpaaq(is,n,nq,apwalm(:,:,:,ias),apwalmq(:,:,:,ias),nmatmax,oq)
114 call olpaloq(is,ias,n,nq,apwalm(:,:,:,ias),apwalmq(:,:,:,ias),nmatmax,oq)
117 call dhmlaa(is,ias,n,n,apwalm(:,:,:,ias),apwalm(:,:,:,ias),dapwalm,dapwalm,&
119 call dhmlalo(is,ias,n,n,apwalm(:,:,:,ias),apwalm(:,:,:,ias),dapwalm, &
122 call dolpaa(is,ias,n,n,apwalm(:,:,:,ias),apwalm(:,:,:,ias),dapwalm,dapwalm,&
124 call dolpalo(is,ias,n,n,dapwalm,dapwalm,nmatmax,od)
133 z2=t1*(dp(i,j)*z1+h(i,j))
134 dlh(i,j)=cmplx(-z2%im,z2%re,8)
136 dlo(i,j)=cmplx(-z2%im,z2%re,8)
144 dlh(i,j)=cmplx(-z1%im,z1%re,8)
146 dlo(i,j)=cmplx(-z1%im,z1%re,8)
161 z2=t1*(dpq(i,j)*z1+hq(i,j))
162 dlhq(i,j)=cmplx(-z2%im,z2%re,8)
164 dloq(i,j)=cmplx(-z2%im,z2%re,8)
167 t1=-
vgkc(l,j,jspn,ik)
170 dlhq(i,j)=cmplx(-z1%im,z1%re,8)
172 dloq(i,j)=cmplx(-z1%im,z1%re,8)
178 t1=
vgkqc(l,i,jspn,ik)
180 dlhq(i,j)=cmplx(-z1%im,z1%re,8)
182 dloq(i,j)=cmplx(-z1%im,z1%re,8)
195 if (ias ==
iasph)
then 197 dz1=
vgc(
ipph,ig)*cmplx(z1%im,-z1%re,8)
201 z2=t1*(dp(i,j)*dz1+dh(i,j))
202 ddlh(i,j)=cmplx(-z2%im,z2%re,8)
204 ddlo(i,j)=cmplx(-z2%im,z2%re,8)
207 t1=-
vgkc(l,j,jspn,ik)
210 ddlh(i,j)=cmplx(-z1%im,z1%re,8)
212 ddlo(i,j)=cmplx(-z1%im,z1%re,8)
220 ddlh(i,j)=cmplx(-z1%im,z1%re,8)
222 ddlo(i,j)=cmplx(-z1%im,z1%re,8)
233 call zhemv(
'U',nm,
zone,dlh,nmatmax,evecfv(:,jst,jspn),1,
zzero,vh,1)
234 call zhemv(
'U',nm,
zone,dlo,nmatmax,evecfv(:,jst,jspn),1,
zzero,vo,1)
237 z1=zdotc(nm,evecfv(:,ist,jspn),1,vh,1)
238 z2=zdotc(nm,evecfv(:,ist,jspn),1,vo,1)
239 ffv(ist,jst)=z1-t1*z2
246 call zhemv(
'U',nm,
zone,dlo,nmatmax,evecfv(:,jst,jspn),1,
zzero,vo,1)
247 call zgemv(
'N',nm,nm,
zone,ddlh,nmatmax,evecfv(:,jst,jspn),1,
zzero,dvh,1)
248 call zgemv(
'N',nm,nm,
zone,ddlo,nmatmax,evecfv(:,jst,jspn),1,
zzero,dvo,1)
252 z2=zdotc(nm,evecfv(:,ist,jspn),1,vo,1)
253 dz1=zdotc(nm,evecfv(:,ist,jspn),1,dvh,1)
254 dz2=zdotc(nm,evecfv(:,ist,jspn),1,dvo,1)
255 dffv(ist,jst)=dffv(ist,jst)+dz1-dt1*z2-t1*dz2
257 call zgemv(
'C',nmq,nm,
zone,dlhq,nmatmax,devecfv(:,jst,jspn),1,
zzero, &
259 call zgemv(
'C',nmq,nm,
zone,dloq,nmatmax,devecfv(:,jst,jspn),1,
zzero, &
262 dz1=2.d0*zdotc(nm,evecfv(:,ist,jspn),1,dvh,1)
263 dz2=2.d0*zdotc(nm,evecfv(:,ist,jspn),1,dvo,1)
264 dffv(ist,jst)=dffv(ist,jst)+dz1-t1*dz2
273 call zgemv(
'N',nstfv,nstfv,
zone,ffv,nstfv,evecsv(i,j),1,
zzero,y,1)
274 call zgemv(
'N',nstfv,nstfv,
zone,dffv,nstfv,evecsv(i,j),1,
zzero,dy,1)
275 call zgemv(
'N',nstfv,nstfv,
zone,ffv,nstfv,devecsv(i,j),1,
zone,dy,1)
276 dz1=zdotc(nstfv,evecsv(i,j),1,dy,1)
277 dz1=dz1+zdotc(nstfv,devecsv(i,j),1,y,1)
278 z1=z1+
occsv(j,jk)*dz1
285 z1=z1+
occsv(j,jk)*dffv(j,j)
287 z1=z1+
doccsv(j,ik)*dble(ffv(j,j))
291 dynibs(l,ias)=dynibs(l,ias)-
wkptnr*z1
298 deallocate(ijg,ijgq,dp,dpq)
299 deallocate(apwalm,apwalmq,dapwalm)
300 deallocate(evecfv,devecfv)
301 deallocate(h,o,dlh,dlo,hq,oq,dh,od)
302 deallocate(dlhq,dloq,ddlh,ddlo)
303 if (
tevecsv)
deallocate(evecsv,devecsv)
complex(8), dimension(:,:), allocatable sfacg
subroutine dhmlaa(is, ias, ngp, ngpq, apwalm, apwalmq, dapwalm, dapwalmq, ld, dh)
subroutine getevecsv(fext, ikp, vpl, evecsv)
integer, dimension(3) ngridg
real(8), dimension(:,:,:), allocatable gkqc
integer, dimension(:,:,:), allocatable ivgig
subroutine getevecfv(fext, ikp, vpl, vgpl, evecfv)
subroutine match(ngp, vgpc, gpc, sfacgp, apwalm)
integer, dimension(:,:,:), allocatable ivkik
subroutine hmlaa(thr, is, ias, ngp, apwalm, ld, h)
complex(8), parameter zone
subroutine getdevecsv(ik, iq, is, ia, ip, devecsv)
complex(8), dimension(:,:,:,:), allocatable sfacgk
subroutine olpaaq(is, ngp, ngpq, apwalm, apwalmq, ld, oq)
subroutine hmlaaq(is, ias, ngp, ngpq, apwalm, apwalmq, ld, hq)
real(8), dimension(:,:), allocatable vgc
subroutine olpaa(tor, is, ngp, apwalm, ld, o)
integer, dimension(:,:), allocatable ngk
complex(8), dimension(:,:,:,:), allocatable sfacgkq
pure subroutine olpalo(is, ias, ngp, apwalm, ld, o)
real(8), dimension(:,:), allocatable doccsv
subroutine olpaloq(is, ias, ngp, ngpq, apwalm, apwalmq, ld, oq)
real(8), dimension(:,:,:,:), allocatable vgkl
real(8), dimension(:,:,:), allocatable devalfv
real(8), dimension(:,:), allocatable ffacg
real(8), dimension(:,:), allocatable occsv
real(8), dimension(:,:), allocatable vgqc
subroutine getevalfv(fext, ikp, vpl, evalfv)
subroutine dolpaa(is, ias, ngp, ngpq, apwalm, apwalmq, dapwalm, dapwalmq, ld, od)
subroutine dforcek(ik, dynibs)
real(8), dimension(:,:,:,:), allocatable vgkc
integer, dimension(:,:), allocatable ivg
complex(8), parameter zzero
real(8), dimension(:,:), allocatable vkl
integer, dimension(2, 3) intgv
integer, dimension(maxatoms *maxspecies) idxis
real(8), dimension(:,:), allocatable ffacgq
real(8), dimension(:,:,:), allocatable gkc
complex(8), dimension(:,:), allocatable sfacgq
subroutine getdevecfv(ik, iq, is, ia, ip, devecfv)
integer, dimension(:,:), allocatable ngkq
real(8), dimension(:,:,:,:), allocatable vgkqc
subroutine hmlaloq(is, ias, ngp, ngpq, apwalm, apwalmq, ld, hq)
pure subroutine dmatch(ias, ip, ngp, vgpc, apwalm, dapwalm)
pure subroutine dolpalo(is, ias, ngp, ngpq, dapwalm, dapwalmq, ld, od)
integer, dimension(:,:,:), allocatable igkig
pure subroutine dhmlalo(is, ias, ngp, ngpq, apwalm, apwalmq, dapwalm, dapwalmq, ld, dh)
integer, dimension(:,:), allocatable ivk
integer, dimension(:,:,:), allocatable igkqig
pure subroutine hmlalo(is, ias, ngp, apwalm, ld, h)