The Elk Code
 
Loading...
Searching...
No Matches
hmlaaq.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 hmlaaq(is,ias,ngp,ngpq,apwalm,apwalmq,ld,hq)
7use modmain
8implicit none
9integer, intent(in) :: is,ias,ngp,ngpq
10complex(8), intent(in) :: apwalm(ngkmax,apwordmax,lmmaxapw)
11complex(8), intent(in) :: apwalmq(ngkmax,apwordmax,lmmaxapw)
12integer, intent(in) :: ld
13complex(8), intent(inout) :: hq(*)
14! local variables
15integer io,jo,i
16integer l0,l1,l2,l3
17integer lm1,lm3,lma,lmb
18complex(8) z1
19! automatic arrays
20complex(8) y(ngp),a(lmoapw(is),ngpq),b(lmoapw(is),ngp)
21i=0
22do l1=0,lmaxapw
23 do lm1=l1**2+1,(l1+1)**2
24 do io=1,apword(l1,is)
25 i=i+1
26 y(:)=0.d0
27 do l3=0,lmaxapw
28 if (mod(l1+l3,2) == 0) then; l0=0; else; l0=1; end if
29 do lm3=l3**2+1,(l3+1)**2
30 do jo=1,apword(l3,is)
31 z1=0.d0
32! kinetic and potential contribution
33 do l2=l0,lmaxo,2
34 lma=l2**2+1; lmb=lma+2*l2
35 z1=z1+sum(gntyry(lma:lmb,lm3,lm1)*haa(lma:lmb,jo,l3,io,l1,ias))
36 end do
37 if (abs(z1%re)+abs(z1%im) > 1.d-12) then
38 call zaxpy(ngp,z1,apwalm(:,jo,lm3),1,y,1)
39 end if
40 end do
41 end do
42 end do
43 a(i,1:ngpq)=apwalmq(1:ngpq,io,lm1)
44 b(i,1:ngp)=y(1:ngp)
45 end do
46 end do
47end do
48call zmctm(lmoapw(is),ngpq,ngp,a,b,ld,hq)
49end subroutine
50
subroutine hmlaaq(is, ias, ngp, ngpq, apwalm, apwalmq, ld, hq)
Definition hmlaaq.f90:7
integer, dimension(0:maxlapw, maxspecies) apword
Definition modmain.f90:758
real(8), dimension(:,:,:,:,:,:), allocatable haa
Definition modmain.f90:856
integer lmaxapw
Definition modmain.f90:197
complex(8), dimension(:,:,:), allocatable gntyry
Definition modmain.f90:862
integer lmaxo
Definition modmain.f90:201
subroutine zmctm(l, m, n, a, b, ld, c)
Definition zmctm.f90:7