The Elk Code
genvchi0.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2010 S. Sharma, J. K. Dewhurst 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 genvchi0(t3hw,ik,lock,vqpl,gclgq,jlgqr,ylmgq,sfacgq,nm,vchi0)
7 use modmain
8 use modomp
9 implicit none
10 ! local variables
11 logical, intent(in) :: t3hw
12 integer, intent(in) :: ik
13 integer(omp_lock_kind), intent(inout) :: lock(nwrf)
14 real(8), intent(in) :: vqpl(3),gclgq(ngrf),jlgqr(njcmax,nspecies,ngrf)
15 complex(8), intent(in) :: ylmgq(lmmaxo,ngrf),sfacgq(ngrf,natmtot)
16 integer, intent(in) :: nm
17 complex(4), intent(inout) :: vchi0(nm,nm,nwrf)
18 ! local variables
19 logical tq0
20 integer isym,jk,jkq,iw,nthd
21 integer nst,nstq,ist,jst,kst,lst
22 integer nm2,ig0,ig1,ig,jg,i,j
23 real(8) vkql(3),ei,ej,eij,t1
24 complex(8) a(3,3),z1
25 complex(4) c1
26 ! automatic arrays
27 integer idx(nstsv),idxq(nstsv)
28 integer ngp(nspnfv),ngpq(nspnfv)
29 ! allocatable arrays
30 integer, allocatable :: igpig(:,:),igpqig(:,:)
31 complex(4), allocatable :: wfmt(:,:,:,:),wfir(:,:,:)
32 complex(4), allocatable :: wfmtq(:,:,:,:),wfirq(:,:,:)
33 complex(4), allocatable :: crhomt(:,:),crhoir(:),cw(:),b(:,:)
34 complex(8), allocatable :: zrhoig(:),pmat(:,:,:)
35 ! check if q = 0
36 tq0=.false.
37 if (sum(abs(vqpl(:))) < epslat) tq0=.true.
38 ! k+q-vector in lattice coordinates
39 vkql(:)=vkl(:,ik)+vqpl(:)
40 ! equivalent reduced k-points for k and k+q
41 jk=ivkik(ivk(1,ik),ivk(2,ik),ivk(3,ik))
42 call findkpt(vkql,isym,jkq)
43 ! count and index states at k and k+q in energy window
44 nst=0
45 do ist=1,nstsv
46  if (abs(evalsv(ist,jk)-efermi) > emaxrf) cycle
47  nst=nst+1
48  idx(nst)=ist
49 end do
50 nstq=0
51 do ist=1,nstsv
52  if (abs(evalsv(ist,jkq)-efermi) > emaxrf) cycle
53  nstq=nstq+1
54  idxq(nstq)=ist
55 end do
56 ! generate the wavefunctions for all states at k and k+q in energy window
57 allocate(igpig(ngkmax,nspnfv))
58 allocate(wfmt(npcmtmax,natmtot,nspinor,nst),wfir(ngtc,nspinor,nst))
59 call genwfsvp_sp(.false.,.false.,nst,idx,ngdgc,igfc,vkl(:,ik),ngp,igpig,wfmt, &
60  ngtc,wfir)
61 deallocate(igpig)
62 allocate(igpqig(ngkmax,nspnfv))
63 allocate(wfmtq(npcmtmax,natmtot,nspinor,nstq),wfirq(ngtc,nspinor,nstq))
64 call genwfsvp_sp(.false.,.false.,nstq,idxq,ngdgc,igfc,vkql,ngpq,igpqig,wfmtq, &
65  ngtc,wfirq)
66 deallocate(igpqig)
67 ! read the momentum matrix elements from file for q = 0
68 if (tq0) then
69  allocate(pmat(nstsv,nstsv,3))
70  call getpmat(vkl(:,ik),pmat)
71 ! divide by unit cell volume
72  t1=1.d0/omega
73  pmat(1:nstsv,1:nstsv,1:3)=t1*pmat(1:nstsv,1:nstsv,1:3)
74 end if
75 nm2=nm**2
76 call holdthd(nst,nthd)
77 !$OMP PARALLEL DEFAULT(SHARED) &
78 !$OMP PRIVATE(crhomt,crhoir,zrhoig,cw,b) &
79 !$OMP PRIVATE(jst,kst,lst,ei,ej,eij,t1,iw) &
80 !$OMP PRIVATE(ig,jg,z1,c1,ig0,ig1,i,j,a) &
81 !$OMP NUM_THREADS(nthd)
82 allocate(crhomt(npcmtmax,natmtot),crhoir(ngtc))
83 allocate(zrhoig(ngrf),cw(nwrf))
84 if (tq0.and.t3hw) then
85  allocate(b(-1:ngrf,-1:ngrf))
86 else
87  allocate(b(ngrf,ngrf))
88 end if
89 !$OMP DO SCHEDULE(DYNAMIC)
90 do ist=1,nst
91  kst=idx(ist)
92  ei=evalsv(kst,jk)
93  do jst=1,nstq
94  lst=idxq(jst)
95  t1=wkptnr*omega*(occsv(kst,jk)-occsv(lst,jkq))
96  if (abs(t1) < 1.d-8) cycle
97  ej=evalsv(lst,jkq)
98  eij=ei-ej
99 ! frequency-dependent part in response function formula for all frequencies
100  do iw=1,nwrf
101  cw(iw)=t1/(eij+wrf(iw))
102  end do
103 ! compute the complex density in G+q-space
104  call gencrho(.true.,.true.,ngtc,wfmt(:,:,:,ist),wfir(:,:,ist), &
105  wfmtq(:,:,:,jst),wfirq(:,:,jst),crhomt,crhoir)
106  call zftcf(ngrf,jlgqr,ylmgq,ngrf,sfacgq,crhomt,crhoir,zrhoig)
107 ! Hermitian part of body
108  do jg=1,ngrf
109  do ig=1,jg-1
110  b(ig,jg)=conjg(b(jg,ig))
111  end do
112  z1=gclgq(jg)*conjg(zrhoig(jg))
113  do ig=jg,ngrf
114  b(ig,jg)=gclgq(ig)*zrhoig(ig)*z1
115  end do
116  end do
117 ! case of q = 0
118  if (tq0) then
119  if (t3hw) then
120  b(-1:1,-1:1)=0.e0
121 ! calculate 3 x ngrf wings of matrix
122  t1=-sqrt(fourpi)/eij
123  do i=-1,1
124  z1=t1*pmat(kst,lst,i+2)
125  b(i,2:ngrf)=z1*conjg(zrhoig(2:ngrf))*gclgq(2:ngrf)
126  do j=2,ngrf
127  b(j,i)=conjg(b(i,j))
128  end do
129  end do
130  else
131 ! use trace of 3 x 3 head of matrix
132  t1=sum(dble(pmat(kst,lst,1:3))**2+aimag(pmat(kst,lst,1:3))**2)/3.d0
133  b(1,1)=(fourpi/eij**2)*t1
134 ! wings of matrix
135  t1=-sqrt(fourpi)/eij
136  z1=(t1/3.d0)*(pmat(kst,lst,1)+pmat(kst,lst,2)+pmat(kst,lst,3))
137  b(1,2:ngrf)=z1*conjg(zrhoig(2:ngrf))*gclgq(2:ngrf)
138  b(2:ngrf,1)=conjg(b(1,2:ngrf))
139  end if
140  end if
141 ! add to body and wings of the response function
142  if (t3hw.or.(mbwgrf < 0)) then
143 ! full matrix
144  do iw=1,nwrf
145  call omp_set_lock(lock(iw))
146  call caxpy(nm2,cw(iw),b,1,vchi0(1,1,iw),1)
147  call omp_unset_lock(lock(iw))
148  end do
149  else
150 ! treat matrix as banded
151  do iw=1,nwrf
152  call omp_set_lock(lock(iw))
153  c1=cw(iw)
154  do jg=1,ngrf
155  ig0=max(jg-mbwgrf,1); ig1=min(jg+mbwgrf,ngrf)
156  vchi0(ig0:ig1,jg,iw)=vchi0(ig0:ig1,jg,iw)+c1*b(ig0:ig1,jg)
157  end do
158  call omp_unset_lock(lock(iw))
159  end do
160  end if
161 ! calculate 3 x 3 head with alternative formula to improve numerical accuracy
162  if (tq0.and.t3hw) then
163  t1=-fourpi/eij
164  cw(1:nwrf)=cw(1:nwrf)/wrf(1:nwrf)
165  do j=1,3
166  do i=1,3
167  a(i,j)=t1*pmat(kst,lst,i)*conjg(pmat(kst,lst,j))
168  end do
169  end do
170 ! add to the head of the response function
171  do iw=1,nwrf
172  call omp_set_lock(lock(iw))
173  vchi0(1:3,1:3,iw)=vchi0(1:3,1:3,iw)+cw(iw)*a(1:3,1:3)
174  call omp_unset_lock(lock(iw))
175  end do
176  end if
177 ! end loop over jst
178  end do
179 ! end loop over ist
180 end do
181 !$OMP END DO
182 deallocate(crhomt,crhoir,zrhoig,cw,b)
183 !$OMP END PARALLEL
184 call freethd(nthd)
185 deallocate(wfmt,wfir,wfmtq,wfirq)
186 if (tq0) deallocate(pmat)
187 end subroutine
188 
real(8) efermi
Definition: modmain.f90:907
integer npcmtmax
Definition: modmain.f90:216
real(8), dimension(:,:), allocatable evalsv
Definition: modmain.f90:921
integer ngtc
Definition: modmain.f90:392
subroutine genwfsvp_sp(tsh, tgp, nst, idx, ngridg_, igfft_, vpl, ngp, igpig, wfmt, ld, wfir)
Definition: genwfsvp_sp.f90:8
integer mbwgrf
Definition: modmain.f90:1166
integer ngkmax
Definition: modmain.f90:499
real(8) omega
Definition: modmain.f90:20
subroutine gencrho(tsh, tspc, ngt, wfmt1, wfir1, wfmt2, wfir2, crhomt, crhoir)
Definition: gencrho.f90:7
subroutine getpmat(vpl, pmat)
Definition: getpmat.f90:7
Definition: modomp.f90:6
complex(8), dimension(:), allocatable wrf
Definition: modmain.f90:1170
integer, dimension(:,:,:), allocatable ivkik
Definition: modmain.f90:467
real(8) emaxrf
Definition: modmain.f90:1162
subroutine genvchi0(t3hw, ik, lock, vqpl, gclgq, jlgqr, ylmgq, sfacgq, nm, vchi0)
Definition: genvchi0.f90:7
subroutine zftcf(ngp, jlgpr, ylmgp, ld, sfacgp, cfmt, cfir, zfgp)
Definition: zftcf.f90:7
integer, dimension(:), allocatable igfc
Definition: modmain.f90:410
real(8), dimension(:,:), allocatable occsv
Definition: modmain.f90:905
integer nspinor
Definition: modmain.f90:267
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
real(8) epslat
Definition: modmain.f90:24
integer, dimension(3) ngdgc
Definition: modmain.f90:388
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
real(8) wkptnr
Definition: modmain.f90:477
real(8), parameter fourpi
Definition: modmain.f90:1234
subroutine findkpt(vpl, isym, ik)
Definition: findkpt.f90:7
integer, dimension(:,:), allocatable ivk
Definition: modmain.f90:465