The Elk Code
genzvbmatk.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2018 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 genzvbmatk(zvmt,zvir,zbmt,zbir,ngp,igpig,wfmt,wfir,wfgp,vbmat)
7 use modmain
8 use modomp
9 implicit none
10 ! arguments
11 ! the potential and field are multiplied by the radial integration weights in
12 ! the muffin-tin and by the characteristic function in the interstitial region
13 complex(8), intent(in) :: zvmt(npcmtmax,natmtot),zvir(ngtc)
14 complex(8), intent(in) :: zbmt(npcmtmax,natmtot,ndmag),zbir(ngtc,ndmag)
15 integer, intent(in) :: ngp,igpig(ngp)
16 complex(4), intent(in) :: wfmt(npcmtmax,natmtot,nspinor,nstsv)
17 ! note that wfir does not have a 1/sqrt(omega) prefactor
18 complex(4), intent(in) :: wfir(ngtc,nspinor,nstsv)
19 complex(4), intent(in) :: wfgp(ngp,nspinor,nstsv)
20 complex(8), intent(out) :: vbmat(nstsv,nstsv)
21 ! local variables
22 integer jst,ispn,nthd
23 integer is,ias,npc
24 integer ld1,ld2,igp
25 ! automatic arrays
26 complex(4) wfmt1(npcmtmax,nspinor),wfir1(ngtc,nspinor),y(nstsv),c(ngp)
27 ld1=npcmtmax*natmtot*nspinor
28 ld2=ngp*nspinor
29 ! zero the matrix elements
30 vbmat(1:nstsv,1:nstsv)=0.d0
31 call holdthd(nstsv,nthd)
32 !$OMP PARALLEL DEFAULT(SHARED) &
33 !$OMP PRIVATE(wfmt1,wfir1,y,c) &
34 !$OMP PRIVATE(jst,ias,is) &
35 !$OMP PRIVATE(npc,ispn,igp) &
36 !$OMP NUM_THREADS(nthd)
37 !-------------------------!
38 ! muffin-tin part !
39 !-------------------------!
40 !$OMP DO SCHEDULE(DYNAMIC)
41 do jst=1,nstsv
42  do ias=1,natmtot
43  is=idxis(ias)
44  npc=npcmt(is)
45 ! apply local potential and magnetic field to spinor wavefunction
46  if (ncmag) then
47 ! non-collinear case
48  call zvbmk1(npc,zvmt(:,ias),zbmt(:,ias,1),zbmt(:,ias,2),zbmt(:,ias,3), &
49  wfmt(:,ias,1,jst),wfmt(:,ias,2,jst),wfmt1,wfmt1(:,2))
50  else
51 ! collinear case
52  call zvbmk2(npc,zvmt(:,ias),zbmt(:,ias,1),wfmt(:,ias,1,jst), &
53  wfmt(:,ias,2,jst),wfmt1,wfmt1(:,2))
54  end if
55 ! compute the inner products
56  call cgemv('C',npc,nstsv,cone,wfmt(1,ias,1,1),ld1,wfmt1(1,1),1,czero,y,1)
57  call cgemv('C',npc,nstsv,cone,wfmt(1,ias,2,1),ld1,wfmt1(1,2),1,cone,y,1)
58  vbmat(1:nstsv,jst)=vbmat(1:nstsv,jst)+y(1:nstsv)
59  end do
60 end do
61 !$OMP END DO
62 !---------------------------!
63 ! interstitial part !
64 !---------------------------!
65 !$OMP DO SCHEDULE(DYNAMIC)
66 do jst=1,nstsv
67 ! apply local potential and magnetic field to spinor wavefunction
68  if (ncmag) then
69 ! non-collinear case
70  call zvbmk1(ngtc,zvir,zbir,zbir(:,2),zbir(:,3),wfir(:,1,jst), &
71  wfir(:,2,jst),wfir1,wfir1(:,2))
72  else
73 ! collinear case
74  call zvbmk2(ngtc,zvir,zbir,wfir(:,1,jst),wfir(:,2,jst),wfir1,wfir1(:,2))
75  end if
76  do ispn=1,nspinor
77 ! Fourier transform to G+p-space
78  call cfftifc(3,ngdgc,-1,wfir1(:,ispn))
79  do igp=1,ngp
80  c(igp)=wfir1(igfc(igpig(igp)),ispn)
81  end do
82  call cgemv('C',ngp,nstsv,cone,wfgp(1,ispn,1),ld2,c,1,czero,y,1)
83  vbmat(1:nstsv,jst)=vbmat(1:nstsv,jst)+y(1:nstsv)
84  end do
85 end do
86 !$OMP END DO
87 !$OMP END PARALLEL
88 call freethd(nthd)
89 
90 contains
91 
92 pure subroutine zvbmk1(n,zv,zb1,zb2,zb3,wf11,wf12,wf21,wf22)
93 implicit none
94 ! arguments
95 integer, intent(in) :: n
96 complex(8), intent(in) :: zv(n),zb1(n),zb2(n),zb3(n)
97 complex(4), intent(in) :: wf11(n),wf12(n)
98 complex(4), intent(out) :: wf21(n),wf22(n)
99 ! local variables
100 integer i
101 complex(8) z1
102 do i=1,n
103  z1=cmplx(-zb2(i)%im,zb2(i)%re,8)
104  wf21(i)=(zv(i)+zb3(i))*wf11(i)+(zb1(i)-z1)*wf12(i)
105  wf22(i)=(zv(i)-zb3(i))*wf12(i)+(zb1(i)+z1)*wf11(i)
106 end do
107 end subroutine
108 
109 pure subroutine zvbmk2(n,zv,zb,wf11,wf12,wf21,wf22)
110 implicit none
111 ! arguments
112 integer, intent(in) :: n
113 complex(8), intent(in) :: zv(n),zb(n)
114 complex(4), intent(in) :: wf11(n),wf12(n)
115 complex(4), intent(out) :: wf21(n),wf22(n)
116 ! local variables
117 integer i
118 do i=1,n
119  wf21(i)=(zv(i)+zb(i))*wf11(i)
120  wf22(i)=(zv(i)-zb(i))*wf12(i)
121 end do
122 end subroutine
123 
124 end subroutine
125 
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
pure subroutine zvbmk2(n, zv, zb, wf11, wf12, wf21, wf22)
Definition: genzvbmatk.f90:110
subroutine genzvbmatk(zvmt, zvir, zbmt, zbir, ngp, igpig, wfmt, wfir, wfgp, vbmat)
Definition: genzvbmatk.f90:7
subroutine cfftifc(nd, n, sgn, c)
Definition: cfftifc_fftw.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
logical ncmag
Definition: modmain.f90:240
pure subroutine zvbmk1(n, zv, zb1, zb2, zb3, wf11, wf12, wf21, wf22)
Definition: genzvbmatk.f90:93