The Elk Code
 
Loading...
Searching...
No Matches
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
6subroutine hmlaloq(is,ias,ngp,ngpq,apwalm,apwalmq,ld,hq)
7use modmain
8implicit none
9! arguments
10integer, intent(in) :: is,ias
11integer, intent(in) :: ngp,ngpq
12complex(8), intent(in) :: apwalm(ngkmax,apwordmax,lmmaxapw)
13complex(8), intent(in) :: apwalmq(ngkmax,apwordmax,lmmaxapw)
14integer, intent(in) :: ld
15complex(8), intent(inout) :: hq(ld,*)
16! local variables
17integer io,ilo
18integer l0,l1,l2,l3
19integer lm1,lm3,lma,lmb
20integer i0,j0,i,j
21complex(8) z1
22do ilo=1,nlorb(is)
23 l1=lorbl(ilo,is)
24 do lm1=l1**2+1,(l1+1)**2
25 i=idxlo(lm1,ilo,ias)
26 i0=ngpq+i
27 j0=ngp+i
28 do l3=0,lmaxapw
29 if (mod(l1+l3,2) == 0) then; l0=0; else; l0=1; end if
30 do lm3=l3**2+1,(l3+1)**2
31 do io=1,apword(l3,is)
32 z1=0.d0
33 do l2=l0,lmaxo,2
34 lma=l2**2+1; lmb=lma+2*l2
35 z1=z1+sum(gntyry(lma:lmb,lm3,lm1)*hloa(lma:lmb,io,l3,ilo,ias))
36 end do
37 if (abs(z1%re)+abs(z1%im) > 1.d-12) then
38 do i=1,ngpq
39 hq(i,j0)=hq(i,j0)+conjg(z1*apwalmq(i,io,lm3))
40 end do
41 do j=1,ngp
42 hq(i0,j)=hq(i0,j)+z1*apwalm(j,io,lm3)
43 end do
44 end if
45 end do
46 end do
47 end do
48 end do
49end do
50end subroutine
51
subroutine hmlaloq(is, ias, ngp, ngpq, apwalm, apwalmq, ld, hq)
Definition hmlaloq.f90:7
integer, dimension(:,:,:), allocatable idxlo
Definition modmain.f90:850
integer, dimension(0:maxlapw, maxspecies) apword
Definition modmain.f90:758
integer lmaxapw
Definition modmain.f90:197
complex(8), dimension(:,:,:), allocatable gntyry
Definition modmain.f90:862
integer lmaxo
Definition modmain.f90:201
integer, dimension(maxspecies) nlorb
Definition modmain.f90:786
real(8), dimension(:,:,:,:,:), allocatable hloa
Definition modmain.f90:858
integer, dimension(maxlorb, maxspecies) lorbl
Definition modmain.f90:796