The Elk Code
gradwfcr2.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 subroutine gradwfcr2(gwf2mt)
7 use modmain
8 implicit none
9 ! arguments
10 real(8), intent(inout) :: gwf2mt(npmtmax,natmtot)
11 ! local variables
12 integer ist,is,ias
13 integer nr,nri,iro,ir
14 integer np,l,lm,i
15 ! allocatable arrays
16 complex(8), allocatable :: wfmt(:),gwfmt(:,:),zfmt(:)
17 allocate(wfmt(npmtmax),gwfmt(npmtmax,3),zfmt(npmtmax))
18 do ias=1,natmtot
19  is=idxis(ias)
20  nr=nrmt(is)
21  nri=nrmti(is)
22  iro=nri+1
23  np=npmt(is)
24  do ist=1,nstsp(is)
25  if (spcore(ist,is).and.(ksp(ist,is) == lsp(ist,is)+1)) then
26  l=lsp(ist,is)
27  do lm=l**2+1,(l+1)**2
28  wfmt(1:np)=0.d0
29  i=lm
30  do ir=1,nri
31  wfmt(i)=rwfcr(ir,1,ist,ias)/rsp(ir,is)
32  i=i+lmmaxi
33  end do
34  do ir=iro,nr
35  wfmt(i)=rwfcr(ir,1,ist,ias)/rsp(ir,is)
36  i=i+lmmaxo
37  end do
38  call gradzfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),wfmt,npmtmax,gwfmt)
39  do i=1,3
40  call zbsht(nr,nri,gwfmt(:,i),zfmt)
41 ! factor of 2 from spin
42  gwf2mt(1:np,ias)=gwf2mt(1:np,ias) &
43  +2.d0*(dble(zfmt(1:np))**2+aimag(zfmt(1:np))**2)
44  end do
45  end do
46  end if
47  end do
48 ! end loops over atoms
49 end do
50 deallocate(wfmt,gwfmt,zfmt)
51 end subroutine
52 
integer, dimension(maxstsp, maxspecies) ksp
Definition: modmain.f90:125
integer, dimension(maxstsp, maxspecies) lsp
Definition: modmain.f90:123
integer lmmaxo
Definition: modmain.f90:203
real(8), dimension(:,:,:), allocatable rlmt
Definition: modmain.f90:179
integer, dimension(maxspecies) npmt
Definition: modmain.f90:213
subroutine gradzfmt(nr, nri, ri, wcr, zfmt, ld, gzfmt)
Definition: gradzfmt.f90:10
subroutine gradwfcr2(gwf2mt)
Definition: gradwfcr2.f90:7
logical, dimension(maxstsp, maxspecies) spcore
Definition: modmain.f90:127
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
integer lmmaxi
Definition: modmain.f90:207
subroutine zbsht(nr, nri, zfmt1, zfmt2)
Definition: zbsht.f90:10
integer, dimension(maxspecies) nstsp
Definition: modmain.f90:113
real(8), dimension(:,:,:,:), allocatable rwfcr
Definition: modmain.f90:939
real(8), dimension(:,:), allocatable rsp
Definition: modmain.f90:135
real(8), dimension(:,:,:), allocatable wcrmt
Definition: modmain.f90:187
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150