The Elk Code
 
Loading...
Searching...
No Matches
gendmatk.f90
Go to the documentation of this file.
1
2! Copyright (C) 2007 F. Bultmark, F. Cricchio, L. Nordstrom and J. K. Dewhurst.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine gendmatk(tspndg,tlmdg,lmin,lmax,ias,nst,idx,ngp,apwalm,evecfv, &
7 evecsv,ld,dmat)
8use modmain
9implicit none
10! arguments
11logical, intent(in) :: tspndg,tlmdg
12integer, intent(in) :: lmin,lmax,ias
13integer, intent(in) :: nst,idx(*),ngp(nspnfv)
14complex(8), intent(in) :: apwalm(ngkmax,apwordmax,lmmaxapw,natmtot,nspnfv)
15complex(8), intent(in) :: evecfv(nmatmax,nstfv,nspnfv),evecsv(nstsv,nstsv)
16integer, intent(in) :: ld
17complex(8), intent(out) :: dmat(ld,nspinor,ld,nspinor,nst)
18! local variables
19integer ispn,jspn,ist,is
20integer nrc,nrci,irco
21integer l,lma,lmb,lm1,lm2
22integer npci,ni,no,i1,i2
23complex(8) zsm
24! allocatable arrays
25complex(8), allocatable :: wfmt(:,:,:)
26if (lmin < 0) then
27 write(*,*)
28 write(*,'("Error(gendmatk): lmin < 0 : ",I8)') lmin
29 write(*,*)
30 stop
31end if
32if (lmax > lmaxo) then
33 write(*,*)
34 write(*,'("Error(gendmatk): lmax > lmaxo : ",2I8)') lmax,lmaxo
35 write(*,*)
36 stop
37end if
38is=idxis(ias)
39nrc=nrcmt(is)
40nrci=nrcmti(is)
41irco=nrci+1
42npci=npcmti(is)
43ni=npci-1
44no=npcmt(is)-npci-1
45! generate the second-variational wavefunctions
46allocate(wfmt(npcmtmax,nspinor,nst))
47call wfmtsv(.true.,lradstp,is,ias,nst,idx,ngp,apwalm,evecfv,evecsv,npcmtmax, &
48 wfmt)
49! zero the density matrix
50dmat(:,:,:,:,:)=0.d0
51! loop over second-variational states
52do ist=1,nst
53 do ispn=1,nspinor
54 do jspn=1,nspinor
55 if (tspndg.and.(ispn /= jspn)) cycle
56 do l=lmin,lmax
57 lma=l**2+1; lmb=lma+2*l
58 do lm1=lma,lmb
59 do lm2=lma,lmb
60 if (tlmdg.and.(lm1 /= lm2)) cycle
61 if (l <= lmaxi) then
62 zsm=sum(wfmt(lm1:lm1+ni:lmmaxi,ispn,ist) &
63 *conjg(wfmt(lm2:lm2+ni:lmmaxi,jspn,ist))*wr2cmt(1:nrci,is))
64 else
65 zsm=0.d0
66 end if
67 i1=npci+lm1; i2=npci+lm2
68 zsm=zsm+sum(wfmt(i1:i1+no:lmmaxo,ispn,ist) &
69 *conjg(wfmt(i2:i2+no:lmmaxo,jspn,ist))*wr2cmt(irco:nrc,is))
70 dmat(lm1,ispn,lm2,jspn,ist)=zsm
71 end do
72 end do
73 end do
74 end do
75 end do
76! end loop over second-variational states
77end do
78deallocate(wfmt)
79end subroutine
80
subroutine gendmatk(tspndg, tlmdg, lmin, lmax, ias, nst, idx, ngp, apwalm, evecfv, evecsv, ld, dmat)
Definition gendmatk.f90:8
integer, dimension(maxspecies) npcmti
Definition modmain.f90:214
integer lmmaxi
Definition modmain.f90:207
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer lradstp
Definition modmain.f90:171
real(8), dimension(:,:), allocatable wr2cmt
Definition modmain.f90:189
integer npcmtmax
Definition modmain.f90:216
integer, dimension(maxspecies) npcmt
Definition modmain.f90:214
integer lmaxo
Definition modmain.f90:201
integer lmaxi
Definition modmain.f90:205
integer lmmaxo
Definition modmain.f90:203
integer, dimension(maxspecies) nrcmti
Definition modmain.f90:211
subroutine wfmtsv(tsh, lrstp, is, ias, nst, idx, ngp, apwalm, evecfv, evecsv, ld, wfmt)
Definition wfmtsv.f90:7