The Elk Code
 
Loading...
Searching...
No Matches
genvcl1221.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: genvcl1221
8! !INTERFACE:
9subroutine genvcl1221(ikp,vcl1221)
10! !USES:
11use modmain
12! !INPUT/OUTPUT PARAMETERS:
13! ikp : k-point from non-reduced set (in,integer)
14! vcl1221 : Coulomb matrix elements (out,real(nstsv,nstsv,nkpt))
15! !DESCRIPTION:
16! Calculates the Coulomb matrix elements of the type $V(1,2,2,1)$. See the
17! routine {\tt genvcl1223} for details.
18!
19! !REVISION HISTORY:
20! Created June 2008 (Sharma)
21!EOP
22!BOC
23implicit none
24! arguments
25integer, intent(in) :: ikp
26real(8), intent(out) :: vcl1221(nstsv,nstsv,nkpt)
27! local variables
28integer ik,ist1,ist2
29integer iv(3),iq,ig
30real(8) vc(3)
31complex(8) z1
32! allocatable arrays
33real(8), allocatable :: vgqc(:,:),gqc(:),gclgq(:),jlgqrmt(:,:,:)
34complex(8), allocatable :: apwalm(:,:,:,:),evecfv(:,:),evecsv(:,:)
35complex(8), allocatable :: ylmgq(:,:),sfacgq(:,:)
36complex(4), allocatable :: wfmt1(:,:,:,:),wfir1(:,:,:)
37complex(4), allocatable :: wfmt2(:,:,:,:),wfir2(:,:,:)
38complex(4), allocatable :: crhomt(:,:),crhoir(:)
39complex(4), allocatable :: cvclmt(:,:),cvclir(:)
40! external functions
41complex(8), external :: zcfinp
42! allocate local arrays
43allocate(vgqc(3,ngvc),gqc(ngvc),gclgq(ngvc))
44allocate(jlgqrmt(0:lnpsd,ngvc,nspecies))
45allocate(apwalm(ngkmax,apwordmax,lmmaxapw,natmtot))
46allocate(evecfv(nmatmax,nstfv),evecsv(nstsv,nstsv))
47allocate(ylmgq(lmmaxo,ngvc),sfacgq(ngvc,natmtot))
48allocate(wfmt1(npcmtmax,natmtot,nspinor,nstsv),wfir1(ngtc,nspinor,nstsv))
49allocate(wfmt2(npcmtmax,natmtot,nspinor,nstsv),wfir2(ngtc,nspinor,nstsv))
50allocate(crhomt(npcmtmax,natmtot),crhoir(ngtc))
51allocate(cvclmt(npcmtmax,natmtot),cvclir(ngtc))
52! get the eigenvectors from file for non-reduced k-point ikp
53call getevecfv(filext,0,vkl(:,ikp),vgkl(:,:,1,ikp),evecfv)
54call getevecsv(filext,0,vkl(:,ikp),evecsv)
55! find the matching coefficients
56call match(ngk(1,ikp),vgkc(:,:,1,ikp),gkc(:,1,ikp),sfacgk(:,:,1,ikp),apwalm)
57! calculate the wavefunctions for all states of passed non-reduced k-point ikp
58call genwfsv_sp(.false.,.false.,nstsv,[0],ngdgc,igfc,ngk(1,ikp),igkig(:,1,ikp),&
59 apwalm,evecfv,evecsv,wfmt2,ngtc,wfir2)
60! start loop over reduced k-point set
61do ik=1,nkpt
62! determine the q-vector
63 iv(1:3)=ivk(1:3,ik)-ivk(1:3,ikp)
64 iv(1:3)=modulo(iv(1:3),ngridk(1:3))
65! check if the q-point is in user-defined set
66 iv(1:3)=iv(1:3)*ngridq(1:3)
67 if (any(mod(iv(1:3),ngridk(1:3)) /= 0)) cycle
68 iv(1:3)=iv(1:3)/ngridk(1:3)
69 iq=ivqiq(iv(1),iv(2),iv(3))
70 vc(1:3)=vkc(1:3,ik)-vkc(1:3,ikp)
71 do ig=1,ngvc
72! determine the G+q-vectors
73 vgqc(1:3,ig)=vgc(1:3,ig)+vc(1:3)
74! G+q-vector length
75 gqc(ig)=sqrt(vgqc(1,ig)**2+vgqc(2,ig)**2+vgqc(3,ig)**2)
76! spherical harmonics for G+q-vectors
77 call genylmv(.true.,lmaxo,vgqc(:,ig),ylmgq(:,ig))
78 end do
79! structure factors for G+q
80 call gensfacgp(ngvc,vgqc,ngvc,sfacgq)
81! generate the regularised Coulomb Green's function in G+q-space
82 call gengclgq(.true.,iq,ngvc,gqc,gclgq)
83! compute the required spherical Bessel functions
84 call genjlgprmt(lnpsd,ngvc,gqc,ngvc,jlgqrmt)
85! find the matching coefficients
86 call match(ngk(1,ik),vgkc(:,:,1,ik),gkc(:,1,ik),sfacgk(:,:,1,ik),apwalm)
87! get the eigenvectors from file
88 call getevecfv(filext,ik,vkl(:,ik),vgkl(:,:,:,ik),evecfv)
89 call getevecsv(filext,ik,vkl(:,ik),evecsv)
90! calculate the wavefunctions for all states of the reduced k-point
91 call genwfsv_sp(.false.,.false.,nstsv,[0],ngdgc,igfc,ngk(1,ik),igkig(:,1,ik),&
92 apwalm,evecfv,evecsv,wfmt1,ngtc,wfir1)
93!----------------------------------------------!
94! valence-valence-valence contribution !
95!----------------------------------------------!
96 do ist1=1,nstsv
97 do ist2=1,nstsv
98! calculate the complex overlap density
99 call gencrho(.true.,.true.,ngtc,wfmt2(:,:,:,ist2),wfir2(:,:,ist2), &
100 wfmt1(:,:,:,ist1),wfir1(:,:,ist1),crhomt,crhoir)
101! compute the potential and G=0 coefficient of the density
104 gclgq,ngvc,jlgqrmt,ylmgq,sfacgq,crhoir,npcmtmax,cvclmt,cvclir)
105 cvclir(1:ngtc)=cvclir(1:ngtc)*cfrc(1:ngtc)
106 z1=zcfinp(crhomt,crhoir,cvclmt,cvclir)
107 vcl1221(ist1,ist2,ik)=wqptnr*z1%re
108! end loop over ist2
109 end do
110! end loop over ist1
111 end do
112! end loop over reduced k-point set
113end do
114deallocate(vgqc,gqc,gclgq,jlgqrmt)
115deallocate(apwalm,evecfv,evecsv,ylmgq,sfacgq)
116deallocate(wfmt1,wfmt2,wfir1,wfir2)
117deallocate(crhomt,crhoir,cvclmt,cvclir)
118end subroutine
119!EOC
120
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 genvcl1221(ikp, vcl1221)
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