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,ngp,ngpq
12 complex(8), intent(in) :: apwalm(ngkmax,apwordmax,lmmaxapw)
13 complex(8), intent(in) :: apwalmq(ngkmax,apwordmax,lmmaxapw)
14 complex(8), intent(in) :: dapwalm(ngkmax,apwordmax,lmmaxapw)
15 complex(8), intent(in) :: dapwalmq(ngkmax,apwordmax,lmmaxapw)
16 integer, intent(in) :: ld
17 complex(8), intent(inout) :: dh(*)
18 ! local variables
19 integer io,jo,i
20 integer l0,l1,l2,l3
21 integer lm1,lm3,lma,lmb
22 complex(8) z1
23 ! automatic arrays
24 complex(8) y1(ngp),a1(lmoapw(is),ngpq),b1(lmoapw(is),ngp)
25 complex(8) y2(ngp),a2(lmoapw(is),ngpq),b2(lmoapw(is),ngp)
26 i=0
27 do l1=0,lmaxapw
28  do lm1=l1**2+1,(l1+1)**2
29  do io=1,apword(l1,is)
30  i=i+1
31  y1(1:ngp)=0.d0
32  do l3=0,lmaxapw
33  l0=merge(0,1,mod(l1+l3,2) == 0)
34  do lm3=l3**2+1,(l3+1)**2
35  do jo=1,apword(l3,is)
36  z1=0.d0
37  do l2=l0,lmaxo,2
38  lma=l2**2+1; lmb=lma+2*l2
39  z1=z1+sum(gntyyy(lma:lmb,lm3,lm1)*dhaa(lma:lmb,jo,l3,io,l1,ias))
40  end do
41  if (abs(z1%re)+abs(z1%im) > 1.d-12) &
42  call zaxpy(ngp,z1,apwalm(:,jo,lm3),1,y1,1)
43  end do
44  end do
45  end do
46  if (ias == iasph) then
47  y2(1:ngp)=0.d0
48  do l3=0,lmaxapw
49  l0=merge(0,1,mod(l1+l3,2) == 0)
50  do lm3=l3**2+1,(l3+1)**2
51  do jo=1,apword(l3,is)
52  z1=0.d0
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,dapwalm(:,jo,lm3),1,y1,1)
59  call zaxpy(ngp,z1,apwalm(:,jo,lm3),1,y2,1)
60  end if
61  end do
62  end do
63  end do
64  a2(i,1:ngpq)=dapwalmq(1:ngpq,io,lm1)
65  b2(i,1:ngp)=y2(1:ngp)
66  end if
67  a1(i,1:ngpq)=apwalmq(1:ngpq,io,lm1)
68  b1(i,1:ngp)=y1(1:ngp)
69  end do
70  end do
71 end do
72 call zmctm(lmoapw(is),ngpq,ngp,a1,b1,ld,dh)
73 if (ias == iasph) call zmctm(lmoapw(is),ngpq,ngp,a2,b2,ld,dh)
74 end subroutine
75 
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:857
integer, dimension(0:maxlapw, maxspecies) apword
Definition: modmain.f90:756
complex(8), dimension(:,:,:), allocatable gntyry
Definition: modmain.f90:863
subroutine zmctm(l, m, n, a, b, ld, c)
Definition: zmctm.f90:7
complex(8), dimension(:,:,:,:,:,:), allocatable dhaa
Definition: modphonon.f90:120