The Elk Code
 
Loading...
Searching...
No Matches
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
6subroutine genzvbmatk(zvmt,zvir,zbmt,zbir,ngp,igpig,wfmt,wfir,wfgp,vbmat)
7use modmain
8use modomp
9implicit 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
13complex(8), intent(in) :: zvmt(npcmtmax,natmtot),zvir(ngtc)
14complex(8), intent(in) :: zbmt(npcmtmax,natmtot,ndmag),zbir(ngtc,ndmag)
15integer, intent(in) :: ngp,igpig(ngp)
16complex(4), intent(in) :: wfmt(npcmtmax,natmtot,nspinor,nstsv)
17! note that wfir does not have a 1/sqrt(omega) prefactor
18complex(4), intent(in) :: wfir(ngtc,nspinor,nstsv)
19complex(4), intent(in) :: wfgp(ngp,nspinor,nstsv)
20complex(8), intent(out) :: vbmat(nstsv,nstsv)
21! local variables
22integer jst,ispn,nthd
23integer is,ias,npc
24integer ld1,ld2,igp
25! automatic arrays
26complex(4) wfmt1(npcmtmax,nspinor),wfir1(ngtc,nspinor),y(nstsv),c(ngp)
27ld1=npcmtmax*natmtot*nspinor
28ld2=ngp*nspinor
29! zero the matrix elements
30vbmat(1:nstsv,1:nstsv)=0.d0
31call 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)
41do 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
60end do
61!$OMP END DO
62!---------------------------!
63! interstitial part !
64!---------------------------!
65!$OMP DO SCHEDULE(DYNAMIC)
66do 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
85end do
86!$OMP END DO
87!$OMP END PARALLEL
88call freethd(nthd)
89return
90
91contains
92
93pure subroutine zvbmk1(n,zv,zb1,zb2,zb3,wf11,wf12,wf21,wf22)
94implicit none
95! arguments
96integer, intent(in) :: n
97complex(8), intent(in) :: zv(n),zb1(n),zb2(n),zb3(n)
98complex(4), intent(in) :: wf11(n),wf12(n)
99complex(4), intent(out) :: wf21(n),wf22(n)
100! local variables
101integer i
102complex(8) z1
103do i=1,n
104 z1=zi*zb2(i)
105 wf21(i)=(zv(i)+zb3(i))*wf11(i)+(zb1(i)-z1)*wf12(i)
106 wf22(i)=(zv(i)-zb3(i))*wf12(i)+(zb1(i)+z1)*wf11(i)
107end do
108end subroutine
109
110pure subroutine zvbmk2(n,zv,zb,wf11,wf12,wf21,wf22)
111implicit none
112! arguments
113integer, intent(in) :: n
114complex(8), intent(in) :: zv(n),zb(n)
115complex(4), intent(in) :: wf11(n),wf12(n)
116complex(4), intent(out) :: wf21(n),wf22(n)
117! local variables
118integer i
119do i=1,n
120 wf21(i)=(zv(i)+zb(i))*wf11(i)
121 wf22(i)=(zv(i)-zb(i))*wf12(i)
122end do
123end subroutine
124
125end subroutine
126
subroutine cfftifc(nd, n, sgn, c)
pure subroutine zvbmk1(n, zv, zb1, zb2, zb3, wf11, wf12, wf21, wf22)
pure subroutine zvbmk2(n, zv, zb, wf11, wf12, wf21, wf22)
subroutine genzvbmatk(zvmt, zvir, zbmt, zbir, ngp, igpig, wfmt, wfir, wfgp, vbmat)
Definition genzvbmatk.f90:7
logical ncmag
Definition modmain.f90:240
complex(4), parameter czero
Definition modmain.f90:1236
integer, dimension(3) ngdgc
Definition modmain.f90:388
complex(8), parameter zi
Definition modmain.f90:1239
integer, dimension(:), allocatable igfc
Definition modmain.f90:410
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer, dimension(maxspecies) npcmt
Definition modmain.f90:214
complex(4), parameter cone
Definition modmain.f90:1236
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106