The Elk Code
hmllolo.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 hmllolo(is,ias,ngp,ld,h)
7 use modmain
8 implicit none
9 ! arguments
10 integer, intent(in) :: is,ias,ngp,ld
11 complex(8), intent(out) :: h(ld,*)
12 ! local variables
13 integer ilo,jlo,i,j
14 integer l0,l1,l2,l3
15 integer lm1,lm3,lma,lmb
16 complex(8) z1
17 do jlo=1,nlorb(is)
18  l3=lorbl(jlo,is)
19  do lm3=l3**2+1,(l3+1)**2
20  j=ngp+idxlo(lm3,jlo,ias)
21  do ilo=1,jlo
22  l1=lorbl(ilo,is)
23  if (mod(l1+l3,2) == 0) then; l0=0; else; l0=1; end if
24  do lm1=l1**2+1,(l1+1)**2
25  i=ngp+idxlo(lm1,ilo,ias)
26  if (i > j) cycle
27  z1=0.d0
28  do l2=l0,lmaxo,2
29  lma=l2**2+1; lmb=(l2+1)**2
30  z1=z1+sum(gntyry(lma:lmb,lm3,lm1)*hlolo(lma:lmb,jlo,ilo,ias))
31  end do
32  h(i,j)=z1
33  end do
34  end do
35  end do
36 end do
37 end subroutine
38 
integer, dimension(maxspecies) nlorb
Definition: modmain.f90:786
integer, dimension(:,:,:), allocatable idxlo
Definition: modmain.f90:855
integer lmaxo
Definition: modmain.f90:201
complex(8), dimension(:,:,:), allocatable gntyry
Definition: modmain.f90:867
integer, dimension(maxlorb, maxspecies) lorbl
Definition: modmain.f90:796
pure subroutine hmllolo(is, ias, ngp, ld, h)
Definition: hmllolo.f90:7
real(8), dimension(:,:,:,:), allocatable hlolo
Definition: modmain.f90:865