The Elk Code
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 
6 subroutine gendmatk(tspndg,tlmdg,lmin,lmax,ias,nst,idx,ngp,apwalm,evecfv, &
7  evecsv,ld,dmat)
8 use modmain
9 implicit none
10 ! arguments
11 logical, intent(in) :: tspndg,tlmdg
12 integer, intent(in) :: lmin,lmax,ias
13 integer, intent(in) :: nst,idx(*),ngp(nspnfv)
14 complex(8), intent(in) :: apwalm(ngkmax,apwordmax,lmmaxapw,natmtot,nspnfv)
15 complex(8), intent(in) :: evecfv(nmatmax,nstfv,nspnfv),evecsv(nstsv,nstsv)
16 integer, intent(in) :: ld
17 complex(8), intent(out) :: dmat(ld,nspinor,ld,nspinor,nst)
18 ! local variables
19 integer ispn,jspn,ist,is
20 integer nrc,nrci,irco
21 integer l,lma,lmb,lm1,lm2
22 integer npci,ni,no,i1,i2
23 complex(8) zsm
24 ! allocatable arrays
25 complex(8), allocatable :: wfmt(:,:,:)
26 if (lmin < 0) then
27  write(*,*)
28  write(*,'("Error(gendmatk): lmin < 0 : ",I8)') lmin
29  write(*,*)
30  stop
31 end if
32 if (lmax > lmaxo) then
33  write(*,*)
34  write(*,'("Error(gendmatk): lmax > lmaxo : ",2I8)') lmax,lmaxo
35  write(*,*)
36  stop
37 end if
38 is=idxis(ias)
39 nrc=nrcmt(is)
40 nrci=nrcmti(is)
41 irco=nrci+1
42 npci=npcmti(is)
43 ni=npci-1
44 no=npcmt(is)-npci-1
45 ! generate the second-variational wavefunctions
46 allocate(wfmt(npcmtmax,nspinor,nst))
47 call wfmtsv(.true.,lradstp,is,ias,nst,idx,ngp,apwalm,evecfv,evecsv,npcmtmax, &
48  wfmt)
49 ! zero the density matrix
50 dmat(:,:,:,:,:)=0.d0
51 ! loop over second-variational states
52 do 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
77 end do
78 deallocate(wfmt)
79 end subroutine
80 
integer, dimension(maxspecies) npcmt
Definition: modmain.f90:214
integer npcmtmax
Definition: modmain.f90:216
integer lmmaxo
Definition: modmain.f90:203
integer lmaxo
Definition: modmain.f90:201
integer lradstp
Definition: modmain.f90:171
subroutine gendmatk(tspndg, tlmdg, lmin, lmax, ias, nst, idx, ngp, apwalm, evecfv, evecsv, ld, dmat)
Definition: gendmatk.f90:8
real(8), dimension(:,:), allocatable wr2cmt
Definition: modmain.f90:189
integer, dimension(maxspecies) npcmti
Definition: modmain.f90:214
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
integer lmmaxi
Definition: modmain.f90:207
integer, dimension(maxspecies) nrcmt
Definition: modmain.f90:173
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
integer lmaxi
Definition: modmain.f90:205