The Elk Code
 
Loading...
Searching...
No Matches
hmlaa.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!BOP
7! !ROUTINE: hmlaa
8! !INTERFACE:
9subroutine hmlaa(thr,is,ias,ngp,apwalm,ld,h)
10! !USES:
11use modmain
12! !INPUT/OUTPUT PARAMETERS:
13! thr : .true. if the matrix h is real valued (in,logical)
14! is : species number (in,integer)
15! ias : joint atom and species number (in,integer)
16! ngp : number of G+p-vectors (in,integer)
17! apwalm : APW matching coefficients (in,complex(ngkmax,apwordmax,lmmaxapw))
18! ld : leading dimension of h (in,integer)
19! h : Hamiltonian matrix (inout,complex(*))
20! !DESCRIPTION:
21! Calculates the APW-APW contribution to the Hamiltonian matrix.
22!
23! !REVISION HISTORY:
24! Created October 2002 (JKD)
25!EOP
26!BOC
27implicit none
28! arguments
29logical, intent(in) :: thr
30integer, intent(in) :: is,ias,ngp
31complex(8), intent(in) :: apwalm(ngkmax,apwordmax,lmmaxapw)
32integer, intent(in) :: ld
33complex(8), intent(inout) :: h(*)
34! local variables
35integer io,jo,i
36integer l0,l1,l2,l3
37integer lm1,lm3,lma,lmb
38complex(8) z1
39! automatic arrays
40complex(8) y(ngp),a(lmoapw(is),ngp),b(lmoapw(is),ngp)
41i=0
42do l1=0,lmaxapw
43 do lm1=l1**2+1,(l1+1)**2
44 do io=1,apword(l1,is)
45 i=i+1
46 y(1:ngp)=0.d0
47 do l3=0,lmaxapw
48 if (mod(l1+l3,2) == 0) then; l0=0; else; l0=1; end if
49 do lm3=l3**2+1,(l3+1)**2
50 do jo=1,apword(l3,is)
51 z1=0.d0
52! kinetic and potential contribution
53 do l2=l0,lmaxo,2
54 lma=l2**2+1; lmb=lma+2*l2
55 z1=z1+sum(gntyry(lma:lmb,lm3,lm1)*haa(lma:lmb,jo,l3,io,l1,ias))
56 end do
57 if (abs(z1%re)+abs(z1%im) > 1.d-12) then
58 call zaxpy(ngp,z1,apwalm(:,jo,lm3),1,y,1)
59 end if
60 end do
61 end do
62 end do
63 a(i,1:ngp)=apwalm(1:ngp,io,lm1)
64 b(i,1:ngp)=y(1:ngp)
65 end do
66 end do
67end do
68if (thr) then
69! matrix H is real
70 call rzmctmu(lmoapw(is),ngp,a,b,ld,h)
71else
72! matrix H is complex
73 call zmctmu(lmoapw(is),ngp,a,b,ld,h)
74end if
75end subroutine
76!EOC
77
subroutine hmlaa(thr, is, ias, ngp, apwalm, ld, h)
Definition hmlaa.f90:10
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 rzmctmu(l, n, a, b, ld, c)
Definition rzmctmu.f90:7
subroutine zmctmu(l, n, a, b, ld, c)
Definition zmctmu.f90:7