The Elk Code
hmlalo.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 pure subroutine hmlalo(is,ias,ngp,apwalm,ld,h)
7 use modmain
8 implicit none
9 ! arguments
10 integer, intent(in) :: is,ias,ngp
11 complex(8), intent(in) :: apwalm(ngkmax,apwordmax,lmmaxapw)
12 integer, intent(in) :: ld
13 complex(8), intent(inout) :: h(ld,*)
14 ! local variables
15 integer io,ilo,j
16 integer l0,l1,l2,l3
17 integer lm1,lm3,lma,lmb
18 complex(8) z1
19 do ilo=1,nlorb(is)
20  l1=lorbl(ilo,is)
21  do lm1=l1**2+1,(l1+1)**2
22  j=ngp+idxlo(lm1,ilo,ias)
23  do l3=0,lmaxapw
24  if (mod(l1+l3,2) == 0) then; l0=0; else; l0=1; end if
25  do lm3=l3**2+1,(l3+1)**2
26  do io=1,apword(l3,is)
27  z1=0.d0
28  do l2=l0,lmaxo,2
29  lma=l2**2+1; lmb=lma+2*l2
30  z1=z1+sum(gntyry(lma:lmb,lm3,lm1)*hloa(lma:lmb,io,l3,ilo,ias))
31  end do
32 ! note that what is actually computed is the Hermitian conjugate of <lo|H|APW>
33  if (abs(z1%re)+abs(z1%im) > 1.d-12) then
34  h(1:ngp,j)=h(1:ngp,j)+conjg(z1*apwalm(1:ngp,io,lm3))
35  end if
36  end do
37  end do
38  end do
39  end do
40 end do
41 end subroutine
42 
integer, dimension(maxspecies) nlorb
Definition: modmain.f90:786
integer, dimension(:,:,:), allocatable idxlo
Definition: modmain.f90:855
integer lmmaxapw
Definition: modmain.f90:199
integer ngkmax
Definition: modmain.f90:499
integer lmaxo
Definition: modmain.f90:201
integer lmaxapw
Definition: modmain.f90:197
integer, dimension(0:maxlapw, maxspecies) apword
Definition: modmain.f90:758
complex(8), dimension(:,:,:), allocatable gntyry
Definition: modmain.f90:867
real(8), dimension(:,:,:,:,:), allocatable hloa
Definition: modmain.f90:863
integer apwordmax
Definition: modmain.f90:760
integer, dimension(maxlorb, maxspecies) lorbl
Definition: modmain.f90:796
pure subroutine hmlalo(is, ias, ngp, apwalm, ld, h)
Definition: hmlalo.f90:7