The Elk Code
 
Loading...
Searching...
No Matches
genpmatk.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
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: genpmatk
8! !INTERFACE:
9subroutine genpmatk(ngp,igpig,vgpc,wfmt,wfgp,pmat)
10! !USES:
11use modmain
12! !INPUT/OUTPUT PARAMETERS:
13! ngp : number of G+p-vectors (in,integer(nspnfv))
14! igpig : index from G+p-vectors to G-vectors (in,integer(ngkmax,nspnfv))
15! vgpc : G+p-vectors in Cartesian coordinates (in,real(3,ngkmax,nspnfv))
16! wfmt : muffin-tin wavefunction in spherical harmonics
17! (in,complex(npcmtmax,natmtot,nspinor,nstsv))
18! wfgp : interstitial wavefunction in plane wave basis
19! (in,complex(ngkmax,nspinor,nstsv))
20! pmat : momentum matrix elements (out,complex(nstsv,nstsv,3))
21! !DESCRIPTION:
22! Calculates the momentum matrix elements
23! $$ P_{ij}=\int d^3r\,\Psi_{i{\bf k}}^*({\bf r})\left(-i\nabla
24! +\frac{1}{4c^2}\left[\vec{\sigma}\times\nabla V_s({\bf r})\right]\right)
25! \Psi_{j{\bf k}}({\bf r}), $$
26! where $V_s$ is the Kohn-Sham effective potential. The second term in the
27! brackets is only calculated if spin-orbit coupling is enabled. See Rathgen
28! and Katsnelson, {\it Physica Scripta} {\bf T109}, 170 (2004).
29!
30! !REVISION HISTORY:
31! Created November 2003 (Sharma)
32! Fixed bug found by Juergen Spitaler, September 2006 (JKD)
33! Added spin-orbit correction, July 2010 (JKD)
34! Fixed bug found by Koichi Kitahara, January 2014 (JKD)
35!EOP
36!BOC
37implicit none
38! arguments
39integer, intent(in) :: ngp(nspnfv),igpig(ngkmax,nspnfv)
40real(8), intent(in) :: vgpc(3,ngkmax,nspnfv)
41complex(8), intent(in) :: wfmt(npcmtmax,natmtot,nspinor,nstsv)
42complex(8), intent(in) :: wfgp(ngkmax,nspinor,nstsv)
43complex(8), intent(out) :: pmat(nstsv,nstsv,3)
44! local variables
45integer ist,jst,ispn,jspn
46integer is,ia,ias,nrc,nrci,npc
47integer ld1,ld2,igp,ifg,i
48real(8) cso
49complex(8) z1,z2,z11,z12,z21,z22,z31,z32
50! automatic arrays
51real(8) rfmt(npcmtmax)
52complex(8) gwfmt(npcmtmax,3,nspinor),gvmt(npcmtmax,3)
53complex(8) zfmt1(npcmtmax,nspinor),zfmt2(npcmtmax,3,nspinor)
54complex(8) gwfir(ngtc,3),z(ngkmax)
55! coefficient of spin-orbit coupling
56cso=1.d0/(4.d0*solsc**2)
57! leading dimensions
58ld1=npcmtmax*natmtot*nspinor
59ld2=ngkmax*nspinor
60! zero the momentum matrix elements array
61pmat(1:nstsv,1:nstsv,1:3)=0.d0
62!---------------------------------!
63! muffin-tin contribution !
64!---------------------------------!
65do is=1,nspecies
66 nrc=nrcmt(is)
67 nrci=nrcmti(is)
68 npc=npcmt(is)
69 do ia=1,natoms(is)
70 ias=idxas(ia,is)
71! compute gradient of potential for spin-orbit correction if required
72 if (spinorb) then
73 call rfmtftoc(nrc,nrci,vsmt(:,ias),rfmt)
74 call rtozfmt(nrc,nrci,rfmt,zfmt1)
75 call gradzfmt(nrc,nrci,rlcmt(:,-1,is),wcrcmt(:,:,is),zfmt1,npcmtmax,gvmt)
76! convert to spherical coordinates
77 do i=1,3
78 call zbshtip(nrc,nrci,gvmt(:,i))
79 end do
80 end if
81 do jst=1,nstsv
82 do ispn=1,nspinor
83! compute the gradient of the wavefunction
84 call gradzfmt(nrc,nrci,rlcmt(:,-1,is),wcrcmt(:,:,is), &
85 wfmt(:,ias,ispn,jst),npcmtmax,gwfmt(:,:,ispn))
86 end do
87! add spin-orbit correction if required
88 if (spinorb) then
89 do ispn=1,nspinor
90! convert wavefunction to spherical coordinates
91 call zbsht(nrc,nrci,wfmt(:,ias,ispn,jst),zfmt1(:,ispn))
92 end do
93! compute i σ ⨯ (∇ V(r)) φ(r)
94 do i=1,npc
95 z1=zi*zfmt1(i,1)
96 z2=zi*zfmt1(i,2)
97 z11=gvmt(i,1)*z1; z12=gvmt(i,1)*z2
98 z21=gvmt(i,2)*z1; z22=gvmt(i,2)*z2
99 z31=gvmt(i,3)*z1; z32=gvmt(i,3)*z2
100 zfmt2(i,1,1)=zmi*z32-z21
101 zfmt2(i,1,2)=zi*z31+z22
102 zfmt2(i,2,1)=z11-z32
103 zfmt2(i,2,2)=-z12-z31
104 zfmt2(i,3,1)=zi*z12+z22
105 zfmt2(i,3,2)=zmi*z11+z21
106 end do
107! convert to spherical harmonics and add to wavefunction gradient
108 do ispn=1,nspinor
109 do i=1,3
110 call zfsht(nrc,nrci,zfmt2(:,i,ispn),zfmt1)
111 gwfmt(1:npc,i,ispn)=gwfmt(1:npc,i,ispn)+cso*zfmt1(1:npc,1)
112 end do
113 end do
114 end if
115 do i=1,3
116 do ispn=1,nspinor
117! apply the radial integral weights
118 call zfmtwr(nrc,nrci,wr2cmt(:,is),gwfmt(:,i,ispn))
119! compute the overlaps
120 call zgemv('C',npc,jst,zone,wfmt(1,ias,ispn,1),ld1,gwfmt(1,i,ispn),1,&
121 zone,pmat(1,jst,i),1)
122 end do
123 end do
124 end do
125 end do
126end do
127!-----------------------------------!
128! interstitial contribution !
129!-----------------------------------!
130do jst=1,nstsv
131 do ispn=1,nspinor
132 jspn=jspnfv(ispn)
133! compute the gradient
134 gwfir(1:ngtc,1:3)=0.d0
135 do igp=1,ngp(jspn)
136 ifg=igfc(igpig(igp,jspn))
137 z1=wfgp(igp,ispn,jst)
138 gwfir(ifg,1:3)=vgpc(1:3,igp,jspn)*zi*z1
139 end do
140 do i=1,3
141! Fourier transform to real-space
142 call zfftifc(3,ngdgc,1,gwfir(:,i))
143! multiply by coarse characteristic function
144 gwfir(1:ngtc,i)=gwfir(1:ngtc,i)*cfrc(1:ngtc)
145! Fourier transform back to G-space
146 call zfftifc(3,ngdgc,-1,gwfir(:,i))
147 end do
148! find the overlaps
149 do i=1,3
150 do igp=1,ngp(jspn)
151 ifg=igfc(igpig(igp,jspn))
152 z(igp)=gwfir(ifg,i)
153 end do
154 call zgemv('C',ngp(jspn),jst,zone,wfgp(1,ispn,1),ld2,z,1,zone, &
155 pmat(1,jst,i),1)
156 end do
157 end do
158end do
159! multiply by -i and set lower triangular part
160do i=1,3
161 do jst=1,nstsv
162 do ist=1,jst-1
163 z1=zmi*pmat(ist,jst,i)
164 pmat(ist,jst,i)=z1
165 pmat(jst,ist,i)=conjg(z1)
166 end do
167 pmat(jst,jst,i)=aimag(pmat(jst,jst,i))
168 end do
169end do
170end subroutine
171!EOC
172
subroutine genpmatk(ngp, igpig, vgpc, wfmt, wfgp, pmat)
Definition genpmatk.f90:10
subroutine gradzfmt(nr, nri, ri, wcr, zfmt, ld, gzfmt)
Definition gradzfmt.f90:10
complex(8), parameter zmi
Definition modmain.f90:1239
real(8), dimension(:), allocatable cfrc
Definition modmain.f90:438
integer, dimension(3) ngdgc
Definition modmain.f90:388
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
integer, dimension(maxatoms, maxspecies) idxas
Definition modmain.f90:42
complex(8), parameter zi
Definition modmain.f90:1239
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
complex(8), parameter zone
Definition modmain.f90:1238
logical spinorb
Definition modmain.f90:230
integer, dimension(:), allocatable igfc
Definition modmain.f90:410
real(8), dimension(:,:,:), allocatable rlcmt
Definition modmain.f90:181
real(8), dimension(:,:), allocatable wr2cmt
Definition modmain.f90:189
real(8), dimension(:,:,:), allocatable wcrcmt
Definition modmain.f90:193
integer, dimension(maxspecies) npcmt
Definition modmain.f90:214
real(8), dimension(:,:), pointer, contiguous vsmt
Definition modmain.f90:649
real(8) solsc
Definition modmain.f90:1252
integer, dimension(2) jspnfv
Definition modmain.f90:291
integer, dimension(maxspecies) nrcmti
Definition modmain.f90:211
integer nspecies
Definition modmain.f90:34
pure subroutine rfmtftoc(nrc, nrci, rfmt, rfcmt)
Definition rfmtftoc.f90:7
pure subroutine rtozfmt(nr, nri, rfmt, zfmt)
Definition rtozfmt.f90:7
subroutine zbsht(nr, nri, zfmt1, zfmt2)
Definition zbsht.f90:10
subroutine zbshtip(nr, nri, zfmt)
Definition zbshtip.f90:7
subroutine zfftifc(nd, n, sgn, z)
pure subroutine zfmtwr(nr, nri, wr, zfmt)
Definition zfmtwr.f90:7
subroutine zfsht(nr, nri, zfmt1, zfmt2)
Definition zfsht.f90:10