The Elk Code
genzvmatk.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2017 T. Mueller, 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 subroutine genzvmatk(zvmt,zvir,ngp,igpig,wfmt,wfir,wfgp,vmat)
7 use modmain
8 use modomp
9 implicit none
10 ! arguments
11 ! the potential is multiplied by the radial integration weights in the
12 ! muffin-tin and by the characteristic function in the interstitial region
13 complex(8), intent(in) :: zvmt(npcmtmax,natmtot),zvir(ngtc)
14 integer, intent(in) :: ngp,igpig(ngp)
15 complex(4), intent(in) :: wfmt(npcmtmax,natmtot,nspinor,nstsv)
16 ! note that wfir does not have a 1/sqrt(omega) prefactor
17 complex(4), intent(in) :: wfir(ngtc,nspinor,nstsv),wfgp(ngp,nspinor,nstsv)
18 complex(8), intent(out) :: vmat(nstsv,nstsv)
19 ! local variables
20 integer jst,ispn,nthd
21 integer is,ias,npc
22 integer ld1,ld2,igp
23 ! automatic arrays
24 complex(4) wfmt1(npcmtmax),wfir1(ngtc),y(nstsv),c(ngp)
25 ld1=npcmtmax*natmtot*nspinor
26 ld2=ngp*nspinor
27 ! zero the matrix elements
28 vmat(1:nstsv,1:nstsv)=0.d0
29 call holdthd(nstsv,nthd)
30 !$OMP PARALLEL DEFAULT(SHARED) &
31 !$OMP PRIVATE(wfmt1,wfir1,y,c) &
32 !$OMP PRIVATE(jst,ispn,ias) &
33 !$OMP PRIVATE(is,npc,igp) &
34 !$OMP NUM_THREADS(nthd)
35 !-------------------------!
36 ! muffin-tin part !
37 !-------------------------!
38 !$OMP DO SCHEDULE(DYNAMIC)
39 do jst=1,nstsv
40  do ispn=1,nspinor
41  do ias=1,natmtot
42  is=idxis(ias)
43  npc=npcmt(is)
44 ! apply complex potential to wavefunction
45  wfmt1(1:npc)=zvmt(1:npc,ias)*wfmt(1:npc,ias,ispn,jst)
46 ! compute the inner products
47  call cgemv('C',npc,nstsv,cone,wfmt(1,ias,ispn,1),ld1,wfmt1,1,czero,y,1)
48  vmat(1:nstsv,jst)=vmat(1:nstsv,jst)+y(1:nstsv)
49  end do
50  end do
51 end do
52 !$OMP END DO
53 !---------------------------!
54 ! interstitial part !
55 !---------------------------!
56 !$OMP DO SCHEDULE(DYNAMIC)
57 do jst=1,nstsv
58  do ispn=1,nspinor
59 ! apply potential to wavefunction
60  wfir1(1:ngtc)=zvir(1:ngtc)*wfir(1:ngtc,ispn,jst)
61 ! Fourier transform to G+p-space
62  call cfftifc(3,ngdgc,-1,wfir1)
63  do igp=1,ngp
64  c(igp)=wfir1(igfc(igpig(igp)))
65  end do
66  call cgemv('C',ngp,nstsv,cone,wfgp(1,ispn,1),ld2,c,1,czero,y,1)
67  vmat(1:nstsv,jst)=vmat(1:nstsv,jst)+y(1:nstsv)
68  end do
69 end do
70 !$OMP END DO
71 !$OMP END PARALLEL
72 call freethd(nthd)
73 end subroutine
74 
integer, dimension(maxspecies) npcmt
Definition: modmain.f90:214
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
subroutine genzvmatk(zvmt, zvir, ngp, igpig, wfmt, wfir, wfgp, vmat)
Definition: genzvmatk.f90:7
integer, dimension(:), allocatable igfc
Definition: modmain.f90:410
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
integer, dimension(3) ngdgc
Definition: modmain.f90:388
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78