The Elk Code
hmlaloq.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2013 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 hmlaloq(is,ias,ngp,ngpq,apwalm,apwalmq,ld,hq)
7 use modmain
8 implicit none
9 ! arguments
10 integer, intent(in) :: is,ias,ngp,ngpq
11 complex(8), intent(in) :: apwalm(ngkmax,apwordmax,lmmaxapw)
12 complex(8), intent(in) :: apwalmq(ngkmax,apwordmax,lmmaxapw)
13 integer, intent(in) :: ld
14 complex(8), intent(inout) :: hq(ld,*)
15 ! local variables
16 integer io,ilo,i,j
17 integer l0,l1,l2,l3
18 integer lm1,lm3,lma,lmb
19 complex(8) z1
20 do ilo=1,nlorb(is)
21  l1=lorbl(ilo,is)
22  do lm1=l1**2+1,(l1+1)**2
23  j=idxlo(lm1,ilo,ias)
24  i=ngpq+j
25  j=ngp+j
26  do l3=0,lmaxapw
27  l0=merge(0,1,mod(l1+l3,2) == 0)
28  do lm3=l3**2+1,(l3+1)**2
29  do io=1,apword(l3,is)
30  z1=0.d0
31  do l2=l0,lmaxo,2
32  lma=l2**2+1; lmb=lma+2*l2
33  z1=z1+sum(gntyry(lma:lmb,lm3,lm1)*hloa(lma:lmb,io,l3,ilo,ias))
34  end do
35  if (abs(z1%re)+abs(z1%im) > 1.d-12) then
36  hq(1:ngpq,j)=hq(1:ngpq,j)+conjg(z1*apwalmq(1:ngpq,io,lm3))
37  hq(i,1:ngp)=hq(i,1:ngp)+z1*apwalm(1:ngp,io,lm3)
38  end if
39  end do
40  end do
41  end do
42  end do
43 end do
44 end subroutine
45 
integer, dimension(maxspecies) nlorb
Definition: modmain.f90:784
integer, dimension(:,:,:), allocatable idxlo
Definition: modmain.f90:851
integer lmaxo
Definition: modmain.f90:201
integer lmaxapw
Definition: modmain.f90:197
integer, dimension(0:maxlapw, maxspecies) apword
Definition: modmain.f90:756
complex(8), dimension(:,:,:), allocatable gntyry
Definition: modmain.f90:863
real(8), dimension(:,:,:,:,:), allocatable hloa
Definition: modmain.f90:859
subroutine hmlaloq(is, ias, ngp, ngpq, apwalm, apwalmq, ld, hq)
Definition: hmlaloq.f90:7
integer, dimension(maxlorb, maxspecies) lorbl
Definition: modmain.f90:794