The Elk Code
 
Loading...
Searching...
No Matches
genvcl1223.f90
Go to the documentation of this file.
1
2! Copyright (C) 2007-2008 J. K. Dewhurst, S. Sharma 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
6!BOP
7! !ROUTINE: genvcl1223
8! !INTERFACE:
9subroutine genvcl1223(ikp,vcl1223)
10! !USES:
11use modmain
12! !INPUT/OUTPUT PARAMETERS:
13! ikp : k-point from non-reduced set (in,integer)
14! vcl1223 : Coulomb matrix elements (out,complex(nstsv,nstsv,nstsv,nkpt))
15! !DESCRIPTION:
16! Calculates Coulomb matrix elements of the type
17! $$ V(1,2,2,3)=\int d^3r\,d^3r'\,\frac{\varphi_{i_1{\bf k}}^*({\bf r})
18! \varphi_{i_2{\bf k}'}({\bf r})\varphi_{i_2{\bf k}'}^*({\bf r}')
19! \varphi_{i_3{\bf k}}({\bf r}')}{|{\bf r}-{\bf r}'|}. $$
20!
21! !REVISION HISTORY:
22! Created 2008 (Sharma)
23!EOP
24!BOC
25implicit none
26! arguments
27integer, intent(in) :: ikp
28complex(8), intent(out) :: vcl1223(nstsv,nstsv,nstsv,nkpt)
29! local variables
30integer ik,ist1,ist2,ist3
31integer iv(3),iq,ig
32real(8) vc(3)
33complex(8) z1
34! allocatable arrays
35real(8), allocatable :: vgqc(:,:),gqc(:),gclgq(:),jlgqrmt(:,:,:)
36complex(8), allocatable :: apwalm(:,:,:,:),evecfv(:,:),evecsv(:,:)
37complex(8), allocatable :: ylmgq(:,:),sfacgq(:,:)
38complex(4), allocatable :: wfmt1(:,:,:,:),wfir1(:,:,:)
39complex(4), allocatable :: wfmt2(:,:,:,:),wfir2(:,:,:)
40complex(4), allocatable :: crhomt(:,:,:),crhoir(:,:)
41complex(4), allocatable :: cvclmt(:,:),cvclir(:)
42! external functions
43complex(8), external :: zcfinp
44! allocate local arrays
45allocate(vgqc(3,ngvc),gqc(ngvc),gclgq(ngvc))
46allocate(jlgqrmt(0:lnpsd,ngvc,nspecies))
47allocate(apwalm(ngkmax,apwordmax,lmmaxapw,natmtot))
48allocate(evecfv(nmatmax,nstfv),evecsv(nstsv,nstsv))
49allocate(ylmgq(lmmaxo,ngvc),sfacgq(ngvc,natmtot))
50allocate(wfmt1(npcmtmax,natmtot,nspinor,nstsv),wfir1(ngtc,nspinor,nstsv))
51allocate(wfmt2(npcmtmax,natmtot,nspinor,nstsv),wfir2(ngtc,nspinor,nstsv))
52allocate(crhomt(npcmtmax,natmtot,nstsv),crhoir(ngtc,nstsv))
53allocate(cvclmt(npcmtmax,natmtot),cvclir(ngtc))
54! get the eigenvectors from file for non-reduced k-point ikp
55call getevecfv(filext,0,vkl(:,ikp),vgkl(:,:,1,ikp),evecfv)
56call getevecsv(filext,0,vkl(:,ikp),evecsv)
57! find the matching coefficients
58call match(ngk(1,ikp),vgkc(:,:,1,ikp),gkc(:,1,ikp),sfacgk(:,:,1,ikp),apwalm)
59! calculate the wavefunctions for all states of passed non-reduced k-point ikp
60call genwfsv_sp(.false.,.false.,nstsv,[0],ngdgc,igfc,ngk(1,ikp),igkig(:,1,ikp),&
61 apwalm,evecfv,evecsv,wfmt2,ngtc,wfir2)
62! start loop over reduced k-point set
63do ik=1,nkpt
64! determine the q-vector
65 iv(:)=ivk(:,ik)-ivk(:,ikp)
66 iv(:)=modulo(iv(:),ngridk(:))
67! check if the q-point is in user-defined set
68 iv(:)=iv(:)*ngridq(:)
69 if (any(mod(iv(:),ngridk(:)) /= 0)) cycle
70 iv(:)=iv(:)/ngridk(:)
71 iq=ivqiq(iv(1),iv(2),iv(3))
72 vc(:)=vkc(:,ik)-vkc(:,ikp)
73 do ig=1,ngvc
74! determine the G+q-vectors
75 vgqc(1:3,ig)=vgc(1:3,ig)+vc(1:3)
76! G+q-vector length
77 gqc(ig)=sqrt(vgqc(1,ig)**2+vgqc(2,ig)**2+vgqc(3,ig)**2)
78! spherical harmonics for G+q-vectors
79 call genylmv(.true.,lmaxo,vgqc(:,ig),ylmgq(:,ig))
80 end do
81! structure factors for G+q
82 call gensfacgp(ngvc,vgqc,ngvc,sfacgq)
83! generate the regularised Coulomb Green's function in G+q-space
84 call gengclgq(.true.,iq,ngvc,gqc,gclgq)
85! compute the required spherical Bessel functions
86 call genjlgprmt(lnpsd,ngvc,gqc,ngvc,jlgqrmt)
87! find the matching coefficients
88 call match(ngk(1,ik),vgkc(:,:,1,ik),gkc(:,1,ik),sfacgk(:,:,1,ik),apwalm)
89! get the eigenvectors from file
90 call getevecfv(filext,ik,vkl(:,ik),vgkl(:,:,:,ik),evecfv)
91 call getevecsv(filext,ik,vkl(:,ik),evecsv)
92! calculate the wavefunctions for all states of the reduced k-point
93 call genwfsv_sp(.false.,.false.,nstsv,[0],ngdgc,igfc,ngk(1,ik),igkig(:,1,ik),&
94 apwalm,evecfv,evecsv,wfmt1,ngtc,wfir1)
95!----------------------------------------------!
96! valence-valence-valence contribution !
97!----------------------------------------------!
98 do ist2=1,nstsv
99 do ist1=1,nstsv
100! calculate the complex overlap density for all states
101 call gencrho(.true.,.true.,ngtc,wfmt2(:,:,:,ist2),wfir2(:,:,ist2), &
102 wfmt1(:,:,:,ist1),wfir1(:,:,ist1),crhomt(:,:,ist1),crhoir(:,ist1))
103 end do
104 do ist3=1,nstsv
105! compute the Coulomb potential
107 crhomt(:,:,ist3),cvclmt)
109 gclgq,ngvc,jlgqrmt,ylmgq,sfacgq,crhoir(:,ist3),npcmtmax,cvclmt,cvclir)
110 cvclir(1:ngtc)=cvclir(1:ngtc)*cfrc(1:ngtc)
111 do ist1=ist3,nstsv
112 z1=zcfinp(crhomt(:,:,ist1),crhoir(:,ist1),cvclmt,cvclir)
113 vcl1223(ist1,ist3,ist2,ik)=wqptnr*z1
114 end do
115 end do
116 end do
117! calculate the lower diagonal
118 do ist1=1,nstsv
119 do ist3=1,ist1-1
120 vcl1223(ist3,ist1,1:nstsv,ik)=conjg(vcl1223(ist1,ist3,1:nstsv,ik))
121 end do
122 end do
123! end loop over reduced k-point set
124end do
125deallocate(vgqc,gqc,gclgq,jlgqrmt)
126deallocate(apwalm,evecfv,evecsv,ylmgq,sfacgq)
127deallocate(wfmt1,wfmt2,wfir1,wfir2)
128deallocate(crhomt,crhoir,cvclmt,cvclir)
129end subroutine
130!EOC
131
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
pure subroutine gengclgq(treg, iq, ngq, gqc, gclgq)
Definition gengclgq.f90:7
subroutine genjlgprmt(lmax, ngp, gpc, ld, jlgprmt)
pure subroutine gensfacgp(ngp, vgpc, ld, sfacgp)
Definition gensfacgp.f90:10
subroutine genvcl1223(ikp, vcl1223)
subroutine genwfsv_sp(tsh, tgp, nst, idx, ngridg_, igfft_, ngp, igpig, apwalm, evecfv, evecsv, wfmt, ld, wfir)
Definition genwfsv_sp.f90:8
pure subroutine genylmv(t4pil, lmax, v, ylm)
Definition genylmv.f90:10
subroutine getevecfv(fext, ikp, vpl, vgpl, evecfv)
Definition getevecfv.f90:10
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition getevecsv.f90:7
subroutine match(ngp, vgpc, gpc, sfacgp, apwalm)
Definition match.f90:10
real(8), dimension(:,:,:,:), allocatable vgkc
Definition modmain.f90:505
real(8), dimension(:,:,:), allocatable gkc
Definition modmain.f90:507
integer nspinor
Definition modmain.f90:267
real(8), dimension(:), allocatable cfrc
Definition modmain.f90:438
integer, dimension(3) ngdgc
Definition modmain.f90:388
character(256) filext
Definition modmain.f90:1300
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 ivk
Definition modmain.f90:465
integer, dimension(:,:,:), allocatable ivqiq
Definition modmain.f90:531
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 natmtot
Definition modmain.f90:40
real(8), dimension(:,:,:,:), allocatable vgkl
Definition modmain.f90:503
integer npcmtmax
Definition modmain.f90:216
integer, dimension(3) ngridk
Definition modmain.f90:448
real(8), dimension(:,:,:), allocatable wprcmt
Definition modmain.f90:191
integer, dimension(maxspecies) npcmt
Definition modmain.f90:214
integer nmatmax
Definition modmain.f90:848
integer lmaxo
Definition modmain.f90:201
real(8) wqptnr
Definition modmain.f90:551
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
integer ngkmax
Definition modmain.f90:499
real(8), dimension(:,:), allocatable vgc
Definition modmain.f90:420
integer nrcmtmax
Definition modmain.f90:175
real(8), dimension(:,:), allocatable vkc
Definition modmain.f90:473
integer lmmaxo
Definition modmain.f90:203
integer ngvc
Definition modmain.f90:398
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
integer, dimension(3) ngridq
Definition modmain.f90:515
integer lnpsd
Definition modmain.f90:628