The Elk Code
genjprk.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2018 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 genjprk(ik,jrmt_,jrir_)
7 use modmain
8 implicit none
9 ! arguments
10 integer, intent(in) :: ik
11 real(8), intent(inout) :: jrmt_(npcmtmax,natmtot,3),jrir_(ngtc,3)
12 ! local variables
13 integer ispn,jspn,nst,ist,jst
14 integer is,ia,ias,nrc,nrci,npc
15 integer igk,ifg,i
16 real(8) wo
17 ! automatic arrays
18 integer idx(nstsv)
19 real(8) rfmt(npcmtmax)
20 complex(8) gwfmt(npcmtmax,3),zfmt1(npcmtmax),zfmt2(npcmtmax)
21 ! allocatable arrays
22 complex(8), allocatable :: apwalm(:,:,:,:,:),evecfv(:,:,:),evecsv(:,:)
23 complex(8), allocatable :: wfmt(:,:,:,:),wfgk(:,:,:),zfft1(:),zfft2(:)
24 allocate(apwalm(ngkmax,apwordmax,lmmaxapw,natmtot,nspnfv))
25 allocate(evecfv(nmatmax,nstfv,nspnfv),evecsv(nstsv,nstsv))
26 ! find the matching coefficients
27 do ispn=1,nspnfv
28  call match(ngk(ispn,ik),vgkc(:,:,ispn,ik),gkc(:,ispn,ik), &
29  sfacgk(:,:,ispn,ik),apwalm(:,:,:,:,ispn))
30 end do
31 ! get the eigenvectors from file
32 call getevecfv(filext,ik,vkl(:,ik),vgkl(:,:,:,ik),evecfv)
33 call getevecsv(filext,ik,vkl(:,ik),evecsv)
34 ! count and index the occupied states
35 nst=0
36 do ist=1,nstsv
37  if (abs(occsv(ist,ik)) < epsocc) cycle
38  nst=nst+1
39  idx(nst)=ist
40 end do
41 ! calculate the second-variational wavefunctions for occupied states
42 allocate(wfmt(npcmtmax,natmtot,nspinor,nst),wfgk(ngkmax,nspinor,nst))
43 call genwfsv(.true.,.true.,nst,idx,ngdgc,igfc,ngk(:,ik),igkig(:,:,ik),apwalm, &
44  evecfv,evecsv,wfmt,ngkmax,wfgk)
45 deallocate(apwalm,evecfv,evecsv)
46 !-------------------------------------------------!
47 ! muffin-tin paramagnetic current density !
48 !-------------------------------------------------!
49 do ist=1,nst
50  jst=idx(ist)
51  wo=wkpt(ik)*occsv(jst,ik)
52  do ispn=1,nspinor
53  do is=1,nspecies
54  nrc=nrcmt(is)
55  nrci=nrcmti(is)
56  npc=npcmt(is)
57  do ia=1,natoms(is)
58  ias=idxas(ia,is)
59 ! compute the gradient of the wavefunction
60  call gradzfmt(nrc,nrci,rlcmt(:,-1,is),wcrcmt(:,:,is), &
61  wfmt(:,ias,ispn,ist),npcmtmax,gwfmt)
62 ! convert wavefunction to spherical coordinates and conjugate
63  call zbsht(nrc,nrci,wfmt(:,ias,ispn,ist),zfmt1)
64  zfmt1(1:npc)=conjg(zfmt1(1:npc))
65  do i=1,3
66 ! convert wavefunction gradient to spherical coordinates
67  call zbsht(nrc,nrci,gwfmt(:,i),zfmt2)
68 ! compute the partial current density
69  rfmt(1:npc)=aimag(zfmt1(1:npc)*zfmt2(1:npc))
70  jrmt_(1:npc,ias,i)=jrmt_(1:npc,ias,i)+wo*rfmt(1:npc)
71  end do
72  end do
73  end do
74  end do
75 end do
76 deallocate(wfmt)
77 !---------------------------------------------------!
78 ! interstitial paramagnetic current density !
79 !---------------------------------------------------!
80 allocate(zfft1(ngtc),zfft2(ngtc))
81 do ist=1,nst
82  jst=idx(ist)
83  wo=wkpt(ik)*occsv(jst,ik)/omega
84  do ispn=1,nspinor
85  jspn=jspnfv(ispn)
86 ! Fourier transform to real-space and conjugate
87  zfft1(1:ngtc)=0.d0
88  do igk=1,ngk(jspn,ik)
89  ifg=igfc(igkig(igk,jspn,ik))
90  zfft1(ifg)=wfgk(igk,ispn,ist)
91  end do
92  call zfftifc(3,ngdgc,1,zfft1)
93  zfft1(1:ngtc)=conjg(zfft1(1:ngtc))
94  do i=1,3
95 ! compute the gradient of the wavefunction
96  zfft2(1:ngtc)=0.d0
97  do igk=1,ngk(jspn,ik)
98  ifg=igfc(igkig(igk,jspn,ik))
99  zfft2(ifg)=vgkc(i,igk,jspn,ik)*zi*wfgk(igk,ispn,ist)
100  end do
101  call zfftifc(3,ngdgc,1,zfft2)
102  jrir_(1:ngtc,i)=jrir_(1:ngtc,i)+wo*aimag(zfft1(1:ngtc)*zfft2(1:ngtc))
103  end do
104  end do
105 end do
106 deallocate(wfgk,zfft1,zfft2)
107 end subroutine
108 
integer nmatmax
Definition: modmain.f90:853
integer, dimension(maxspecies) npcmt
Definition: modmain.f90:214
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition: getevecsv.f90:7
character(256) filext
Definition: modmain.f90:1301
subroutine genwfsv(tsh, tgp, nst, idx, ngridg_, igfft_, ngp, igpig, apwalm, evecfv, evecsv, wfmt, ld, wfir)
Definition: genwfsv.f90:11
integer lmmaxapw
Definition: modmain.f90:199
integer, dimension(maxatoms, maxspecies) idxas
Definition: modmain.f90:42
integer ngkmax
Definition: modmain.f90:499
subroutine getevecfv(fext, ikp, vpl, vgpl, evecfv)
Definition: getevecfv.f90:10
real(8) omega
Definition: modmain.f90:20
subroutine match(ngp, vgpc, gpc, sfacgp, apwalm)
Definition: match.f90:10
subroutine gradzfmt(nr, nri, ri, wcr, zfmt, ld, gzfmt)
Definition: gradzfmt.f90:10
real(8) epsocc
Definition: modmain.f90:903
subroutine genjprk(ik, jrmt_, jrir_)
Definition: genjprk.f90:7
complex(8), dimension(:,:,:,:), allocatable sfacgk
Definition: modmain.f90:509
integer, dimension(:,:), allocatable ngk
Definition: modmain.f90:497
subroutine zfftifc(nd, n, sgn, z)
Definition: zfftifc_fftw.f90:7
real(8), dimension(:), allocatable wkpt
Definition: modmain.f90:475
real(8), dimension(:,:,:,:), allocatable vgkl
Definition: modmain.f90:503
integer, dimension(:), allocatable igfc
Definition: modmain.f90:410
real(8), dimension(:,:,:), allocatable rlcmt
Definition: modmain.f90:181
real(8), dimension(:,:), allocatable occsv
Definition: modmain.f90:905
integer nspinor
Definition: modmain.f90:267
real(8), dimension(:,:,:,:), allocatable vgkc
Definition: modmain.f90:505
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
integer, dimension(maxspecies) natoms
Definition: modmain.f90:36
integer apwordmax
Definition: modmain.f90:760
real(8), dimension(:,:,:), allocatable wcrcmt
Definition: modmain.f90:193
subroutine zbsht(nr, nri, zfmt1, zfmt2)
Definition: zbsht.f90:10
real(8), dimension(:,:,:), allocatable gkc
Definition: modmain.f90:507
integer, dimension(3) ngdgc
Definition: modmain.f90:388
integer nspecies
Definition: modmain.f90:34
complex(8), parameter zi
Definition: modmain.f90:1240
integer, dimension(maxspecies) nrcmt
Definition: modmain.f90:173
integer, dimension(maxspecies) nrcmti
Definition: modmain.f90:211
integer, dimension(:,:,:), allocatable igkig
Definition: modmain.f90:501
integer nstfv
Definition: modmain.f90:887
integer, dimension(2) jspnfv
Definition: modmain.f90:291
integer nspnfv
Definition: modmain.f90:289