The Elk Code
 
Loading...
Searching...
No Matches
hmlxbsek.f90
Go to the documentation of this file.
1
2! Copyright (C) 2010 S. Sharma, J. K. Dewhurst 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
6subroutine hmlxbsek(ik2)
7use modmain
8implicit none
9! arguments
10integer, intent(in) :: ik2
11! local variables
12integer ik1,ist1,ist2,jst1,jst2
13integer i1,i2,j1,j2,a1,a2,b1,b2
14integer is,ias,l
15real(8) t0
16complex(8) z1
17! automatic arrays
18integer ngp(nspnfv)
19! allocatable arrays
20integer, allocatable :: igpig(:,:)
21complex(4), allocatable :: wfmt1(:,:,:,:),wfir1(:,:,:)
22complex(4), allocatable :: wfmt2(:,:,:,:),wfir2(:,:,:)
23complex(4), allocatable :: crhomt(:,:),crhoir(:)
24complex(4), allocatable :: cvclmt(:,:,:),cvclir(:,:)
25! external functions
26complex(8), external :: zcfinp
27! allocate local arrays
28allocate(igpig(ngkmax,nspnfv))
29allocate(wfmt1(npcmtmax,natmtot,nspinor,nstsv),wfir1(ngtc,nspinor,nstsv))
30allocate(wfmt2(npcmtmax,natmtot,nspinor,nstsv),wfir2(ngtc,nspinor,nstsv))
31allocate(crhomt(npcmtmax,natmtot),crhoir(ngtc))
32allocate(cvclmt(npcmtmax,natmtot,nvcbse),cvclir(ngtc,nvcbse))
33! calculate the wavefunctions for all states of k-point ik2
34call genwfsvp_sp(.false.,.false.,nstsv,[0],ngdgc,igfc,vkl(:,ik2),ngp,igpig, &
35 wfmt2,ngtc,wfir2)
36l=0
37do i2=1,nvbse
38 ist2=istbse(i2,ik2)
39 do j2=1,ncbse
40 jst2=jstbse(j2,ik2)
41 a2=ijkbse(i2,j2,ik2)
42 l=l+1
43! calculate the complex overlap density
44 call gencrho(.true.,.true.,ngtc,wfmt2(:,:,:,ist2),wfir2(:,:,ist2), &
45 wfmt2(:,:,:,jst2),wfir2(:,:,jst2),crhomt,crhoir)
46! compute the Coulomb potential
48 cvclmt(:,:,l))
50 ngvec,jlgrmt,ylmg,sfacg,crhoir,npcmtmax,cvclmt(:,:,l),cvclir(:,l))
51 cvclir(:,l)=cvclir(:,l)*cfrc(:)
52 end do
53end do
55! start loop over ik1
56do ik1=1,nkptnr
57 if (ik1 == ik2) then
58 wfmt1(:,:,:,:)=wfmt2(:,:,:,:)
59 wfir1(:,:,:)=wfir2(:,:,:)
60 else
61 call genwfsvp_sp(.false.,.false.,nstsv,[0],ngdgc,igfc,vkl(:,ik1),ngp,igpig,&
62 wfmt1,ngtc,wfir1)
63 end if
64 do i1=1,nvbse
65 ist1=istbse(i1,ik1)
66 do j1=1,ncbse
67 jst1=jstbse(j1,ik1)
68 a1=ijkbse(i1,j1,ik1)
69! calculate the complex overlap density
70 call gencrho(.true.,.true.,ngtc,wfmt1(:,:,:,ist1),wfir1(:,:,ist1), &
71 wfmt1(:,:,:,jst1),wfir1(:,:,jst1),crhomt,crhoir)
72 l=0
73 do i2=1,nvbse
74 ist2=istbse(i2,ik2)
75 do j2=1,ncbse
76 jst2=jstbse(j2,ik2)
77 a2=ijkbse(i2,j2,ik2)
78 l=l+1
79! compute the matrix element
80 z1=t0*zcfinp(crhomt,crhoir,cvclmt(:,:,l),cvclir(:,l))
81 hmlbse(a1,a2)=hmlbse(a1,a2)+z1
82! compute off-diagonal blocks if required
83 if (bsefull) then
84 b1=a1+nbbse
85 b2=a2+nbbse
86 hmlbse(b1,b2)=hmlbse(b1,b2)-conjg(z1)
87! conjugate the potential
88 do ias=1,natmtot
89 is=idxis(ias)
90 call cfmtconj(nrcmt(is),nrcmti(is),npcmt(is),cvclmt(:,ias,l))
91 end do
92 cvclir(:,l)=conjg(cvclir(:,l))
93 z1=t0*zcfinp(crhomt,crhoir,cvclmt(:,:,l),cvclir(:,l))
94 hmlbse(a1,b2)=hmlbse(a1,b2)+z1
95 hmlbse(b1,a2)=hmlbse(b1,a2)-conjg(z1)
96 end if
97 end do
98 end do
99 end do
100 end do
101end do
102deallocate(igpig,wfmt1,wfmt2,wfir1,wfir2)
103deallocate(crhomt,crhoir,cvclmt,cvclir)
104end subroutine
105
pure subroutine cfmtconj(nr, nri, np, cfmt)
Definition cfmtconj.f90:7
subroutine cpotcoul(nrmt_, nrmti_, npmt_, ld1, rl, ngridg_, igfft_, ngp, gpc, gclgp, ld2, jlgprmt, ylmgp, sfacgp, crhoir, ld3, cvclmt, cvclir)
Definition cpotcoul.f90:8
subroutine gencrho(tsh, tspc, ngt, wfmt1, wfir1, wfmt2, wfir2, crhomt, crhoir)
Definition gencrho.f90:7
subroutine gencvclmt(nrmt_, nrmti_, ld1, rl, wpr, ld2, crhomt, cvclmt)
Definition gencvclmt.f90:7
subroutine genwfsvp_sp(tsh, tgp, nst, idx, ngridg_, igfft_, vpl, ngp, igpig, wfmt, ld, wfir)
subroutine hmlxbsek(ik2)
Definition hmlxbsek.f90:7
logical bsefull
Definition modmain.f90:1203
integer nvcbse
Definition modmain.f90:1186
integer nkptnr
Definition modmain.f90:463
integer nspinor
Definition modmain.f90:267
integer, dimension(:,:), allocatable istbse
Definition modmain.f90:1192
real(8), dimension(:), allocatable cfrc
Definition modmain.f90:438
integer, dimension(3) ngdgc
Definition modmain.f90:388
real(8), dimension(:,:,:), allocatable jlgrmt
Definition modmain.f90:426
integer nbbse
Definition modmain.f90:1188
integer ngvec
Definition modmain.f90:396
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
complex(8), dimension(:,:), allocatable hmlbse
Definition modmain.f90:1198
integer, dimension(:), allocatable igfc
Definition modmain.f90:410
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
complex(8), dimension(:,:), allocatable sfacg
Definition modmain.f90:430
integer ngtc
Definition modmain.f90:392
real(8), dimension(:,:,:), allocatable rlcmt
Definition modmain.f90:181
integer natmtot
Definition modmain.f90:40
integer, dimension(:,:,:), allocatable ijkbse
Definition modmain.f90:1196
integer npcmtmax
Definition modmain.f90:216
integer, dimension(:,:), allocatable jstbse
Definition modmain.f90:1194
real(8), dimension(:,:,:), allocatable wprcmt
Definition modmain.f90:191
integer, dimension(maxspecies) npcmt
Definition modmain.f90:214
integer nstsv
Definition modmain.f90:886
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
integer ngkmax
Definition modmain.f90:499
real(8) occmax
Definition modmain.f90:898
real(8), dimension(:), allocatable gclg
Definition modmain.f90:424
integer nrcmtmax
Definition modmain.f90:175
integer ngvc
Definition modmain.f90:398
complex(8), dimension(:,:), allocatable ylmg
Definition modmain.f90:428
integer, dimension(maxspecies) nrcmti
Definition modmain.f90:211
real(8), dimension(:), allocatable gc
Definition modmain.f90:422
real(8) wkptnr
Definition modmain.f90:477
integer ncbse
Definition modmain.f90:1176
integer nvbse
Definition modmain.f90:1176