The Elk Code
dhmlaa.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2012 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 
6 subroutine dhmlaa(is,ias,ngp,ngpq,apwalm,apwalmq,dapwalm,dapwalmq,ld,dh)
7 use modmain
8 use modphonon
9 implicit none
10 ! arguments
11 integer, intent(in) :: is,ias
12 integer, intent(in) :: ngp,ngpq
13 complex(8), intent(in) :: apwalm(ngkmax,apwordmax,lmmaxapw)
14 complex(8), intent(in) :: apwalmq(ngkmax,apwordmax,lmmaxapw)
15 complex(8), intent(in) :: dapwalm(ngkmax,apwordmax,lmmaxapw)
16 complex(8), intent(in) :: dapwalmq(ngkmax,apwordmax,lmmaxapw)
17 integer, intent(in) :: ld
18 complex(8), intent(inout) :: dh(*)
19 ! local variables
20 integer io,jo,i
21 integer l0,l1,l2,l3
22 integer lm1,lm3,lma,lmb
23 complex(8) z1
24 ! automatic arrays
25 complex(8) y1(ngp),a1(lmoapw(is),ngpq),b1(lmoapw(is),ngp)
26 complex(8) y2(ngp),a2(lmoapw(is),ngpq),b2(lmoapw(is),ngp)
27 i=0
28 do l1=0,lmaxapw
29  do lm1=l1**2+1,(l1+1)**2
30  do io=1,apword(l1,is)
31  i=i+1
32  y1(1:ngp)=0.d0
33  do l3=0,lmaxapw
34  if (mod(l1+l3,2) == 0) then; l0=0; else; l0=1; end if
35  do lm3=l3**2+1,(l3+1)**2
36  do jo=1,apword(l3,is)
37  z1=0.d0
38  do l2=l0,lmaxo,2
39  lma=l2**2+1; lmb=lma+2*l2
40  z1=z1+sum(gntyyy(lma:lmb,lm3,lm1)*dhaa(lma:lmb,jo,l3,io,l1,ias))
41  end do
42  if (abs(z1%re)+abs(z1%im) > 1.d-12) then
43  call zaxpy(ngp,z1,apwalm(:,jo,lm3),1,y1,1)
44  end if
45  end do
46  end do
47  end do
48  if (ias == iasph) then
49  y2(1:ngp)=0.d0
50  do l3=0,lmaxapw
51  if (mod(l1+l3,2) == 0) then; l0=0; else; l0=1; end if
52  do lm3=l3**2+1,(l3+1)**2
53  do jo=1,apword(l3,is)
54  z1=0.d0
55  do l2=l0,lmaxo,2
56  lma=l2**2+1; lmb=lma+2*l2
57  z1=z1+sum(gntyry(lma:lmb,lm3,lm1)*haa(lma:lmb,jo,l3,io,l1,ias))
58  end do
59  if (abs(z1%re)+abs(z1%im) > 1.d-12) then
60  call zaxpy(ngp,z1,dapwalm(:,jo,lm3),1,y1,1)
61  call zaxpy(ngp,z1,apwalm(:,jo,lm3),1,y2,1)
62  end if
63  end do
64  end do
65  end do
66  a2(i,1:ngpq)=dapwalmq(1:ngpq,io,lm1)
67  b2(i,1:ngp)=y2(1:ngp)
68  end if
69  a1(i,1:ngpq)=apwalmq(1:ngpq,io,lm1)
70  b1(i,1:ngp)=y1(1:ngp)
71  end do
72  end do
73 end do
74 call zmctm(lmoapw(is),ngpq,ngp,a1,b1,ld,dh)
75 if (ias == iasph) then
76  call zmctm(lmoapw(is),ngpq,ngp,a2,b2,ld,dh)
77 end if
78 end subroutine
79 
subroutine dhmlaa(is, ias, ngp, ngpq, apwalm, apwalmq, dapwalm, dapwalmq, ld, dh)
Definition: dhmlaa.f90:7
integer iasph
Definition: modphonon.f90:15
integer lmaxo
Definition: modmain.f90:201
real(8), dimension(:,:,:), allocatable gntyyy
Definition: modphonon.f90:126
integer lmaxapw
Definition: modmain.f90:197
real(8), dimension(:,:,:,:,:,:), allocatable haa
Definition: modmain.f90:861
integer, dimension(0:maxlapw, maxspecies) apword
Definition: modmain.f90:758
complex(8), dimension(:,:,:), allocatable gntyry
Definition: modmain.f90:867
subroutine zmctm(l, m, n, a, b, ld, c)
Definition: zmctm.f90:7
complex(8), dimension(:,:,:,:,:,:), allocatable dhaa
Definition: modphonon.f90:120