The Elk Code
genvmatk.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2007 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
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 genvmatk(vmt,vir,ngp,igpig,wfmt,ld,wfgp,vmat)
7 use modmain
8 use moddftu
9 use modomp
10 implicit none
11 ! arguments
12 ! the potential is multiplied by the radial integration weights in the
13 ! muffin-tin and by the characteristic function in the interstitial region
14 real(8), intent(in) :: vmt(npcmtmax,natmtot),vir(ngtc)
15 integer, intent(in) :: ngp(nspnfv),igpig(ngkmax,nspnfv)
16 complex(4), intent(in) :: wfmt(npcmtmax,natmtot,nspinor,nstsv)
17 integer, intent(in) :: ld
18 complex(4), intent(in) :: wfgp(ld,nspinor,nstsv)
19 complex(8), intent(out) :: vmat(nstsv,nstsv)
20 ! local variables
21 integer ist,jst,nj,ispn,jspn
22 integer is,ias,nrc,nrci,nrco
23 integer npc,npc2,ipco
24 integer ld1,ld2,n,igp,nthd
25 ! automatic arrays
26 complex(4) wfmt1(npcmtmax),wfmt2(npcmtmax)
27 complex(4) wfir(ngtc),c(ngkmax),y(nstsv)
28 ! external functions
29 real(4), external :: sdot
30 ld1=npcmtmax*natmtot*nspinor
31 ld2=ld*nspinor
32 ! zero the upper triangular matrix elements
33 do jst=1,nstsv
34  vmat(1:jst,jst)=0.d0
35 end do
36 call holdthd(nstsv,nthd)
37 !$OMP PARALLEL DEFAULT(SHARED) &
38 !$OMP PRIVATE(wfmt1,wfmt2,wfir,c,y) &
39 !$OMP PRIVATE(jst,nj,ispn,jspn) &
40 !$OMP PRIVATE(ias,is,npc,npc2,nrc) &
41 !$OMP PRIVATE(nrci,nrco,ipco,n,igp) &
42 !$OMP NUM_THREADS(nthd)
43 !-------------------------!
44 ! muffin-tin part !
45 !-------------------------!
46 !$OMP DO SCHEDULE(DYNAMIC)
47 do jst=1,nstsv
48  nj=jst-1
49  do ispn=1,nspinor
50  do ias=1,natmtot
51  is=idxis(ias)
52  npc=npcmt(is)
53  npc2=npc*2
54 ! apply local potential to wavefunction
55  wfmt1(1:npc)=vmt(1:npc,ias)*wfmt(1:npc,ias,ispn,jst)
56 ! apply muffin-tin DFT+U potential matrix if required (note that this should be
57 ! used only in the spin-unpolarised case)
58  if (dftu /= 0) then
59  if (any(tvmmt(0:lmaxdm,ias))) then
60  nrc=nrcmt(is)
61  nrci=nrcmti(is)
62  nrco=nrc-nrci
63  ipco=npcmti(is)+1
64  call cgemm('N','N',lmmaxi,nrci,lmmaxi,cone,vmatmti(1,1,1,1,ias), &
65  lmmaxi,wfmt(1,ias,ispn,jst),lmmaxi,czero,wfmt2,lmmaxi)
66  call cgemm('N','N',lmmaxo,nrco,lmmaxo,cone,vmatmto(1,1,1,1,ias), &
67  lmmaxo,wfmt(ipco,ias,ispn,jst),lmmaxo,czero,wfmt2(ipco),lmmaxo)
68  call cfcmtwr(nrc,nrci,wr2cmt(:,is),wfmt2)
69  wfmt1(1:npc)=wfmt1(1:npc)+wfmt2(1:npc)
70  end if
71  end if
72 ! compute the inner products
73  call cgemv('C',npc,nj,cone,wfmt(1,ias,ispn,1),ld1,wfmt1,1,czero,y,1)
74  vmat(1:nj,jst)=vmat(1:nj,jst)+y(1:nj)
75  vmat(jst,jst)=vmat(jst,jst)+sdot(npc2,wfmt(1,ias,ispn,jst),1,wfmt1,1)
76  end do
77  end do
78 end do
79 !$OMP END DO
80 !---------------------------!
81 ! interstitial part !
82 !---------------------------!
83 !$OMP DO SCHEDULE(DYNAMIC)
84 do jst=1,nstsv
85  nj=jst-1
86  do ispn=1,nspinor
87  jspn=jspnfv(ispn)
88  n=ngp(jspn)
89 ! Fourier transform wavefunction to real-space
90  wfir(1:ngtc)=0.e0
91  do igp=1,n
92  wfir(igfc(igpig(igp,jspn)))=wfgp(igp,ispn,jst)
93  end do
94  call cfftifc(3,ngdgc,1,wfir)
95 ! apply potential to wavefunction
96  wfir(1:ngtc)=vir(1:ngtc)*wfir(1:ngtc)
97 ! Fourier transform to G+p-space
98  call cfftifc(3,ngdgc,-1,wfir)
99  do igp=1,n
100  c(igp)=wfir(igfc(igpig(igp,jspn)))
101  end do
102 ! compute the inner products
103  call cgemv('C',n,nj,cone,wfgp(1,ispn,1),ld2,c,1,czero,y,1)
104  vmat(1:nj,jst)=vmat(1:nj,jst)+y(1:nj)
105  vmat(jst,jst)=vmat(jst,jst)+sdot(n*2,wfgp(1,ispn,jst),1,c,1)
106  end do
107 end do
108 !$OMP END DO
109 !$OMP END PARALLEL
110 call freethd(nthd)
111 ! lower triangular part
112 do ist=1,nstsv
113  do jst=1,ist-1
114  vmat(ist,jst)=conjg(vmat(jst,ist))
115  end do
116 end do
117 end subroutine
118 
integer, dimension(maxspecies) npcmt
Definition: modmain.f90:214
integer lmmaxo
Definition: modmain.f90:203
complex(4), dimension(:,:,:,:,:), allocatable vmatmti
Definition: moddftu.f90:22
complex(4), parameter czero
Definition: modmain.f90:1239
Definition: modomp.f90:6
complex(4), parameter cone
Definition: modmain.f90:1239
subroutine cfftifc(nd, n, sgn, c)
Definition: cfftifc_fftw.f90:7
complex(4), dimension(:,:,:,:,:), allocatable vmatmto
Definition: moddftu.f90:22
subroutine genvmatk(vmt, vir, ngp, igpig, wfmt, ld, wfgp, vmat)
Definition: genvmatk.f90:7
integer, dimension(:), allocatable igfc
Definition: modmain.f90:410
logical, dimension(:,:), allocatable tvmmt
Definition: moddftu.f90:26
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 dftu
Definition: moddftu.f90:32
integer, dimension(3) ngdgc
Definition: modmain.f90:388
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer, dimension(maxspecies) nrcmt
Definition: modmain.f90:173
integer, dimension(maxspecies) nrcmti
Definition: modmain.f90:211
integer, parameter lmaxdm
Definition: moddftu.f90:13
pure subroutine cfcmtwr(nr, nri, wr, cfmt)
Definition: cfcmtwr.f90:7
integer, dimension(2) jspnfv
Definition: modmain.f90:291