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