The Elk Code
 
Loading...
Searching...
No Matches
dolpaa.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 dolpaa(is,ias,ngp,ngpq,apwalm,apwalmq,dapwalm,dapwalmq,ld,od)
7use modmain
8use modphonon
9implicit none
10! arguments
11integer, intent(in) :: is,ias
12integer, intent(in) :: ngp,ngpq
13complex(8), intent(in) :: apwalm(ngkmax,apwordmax,lmmaxapw)
14complex(8), intent(in) :: apwalmq(ngkmax,apwordmax,lmmaxapw)
15complex(8), intent(in) :: dapwalm(ngkmax,apwordmax,lmmaxapw)
16complex(8), intent(in) :: dapwalmq(ngkmax,apwordmax,lmmaxapw)
17integer, intent(in) :: ld
18complex(8), intent(inout) :: od(*)
19! local variables
20integer io,l,lm,i
21! automatic arrays
22complex(8) a(lmoapw(is),ngpq),b(lmoapw(is),ngp)
23if (ias /= iasph) return
24i=0
25do l=0,lmaxapw
26 do lm=l**2+1,(l+1)**2
27 do io=1,apword(l,is)
28 i=i+1
29 a(i,1:ngpq)=apwalmq(1:ngpq,io,lm)
30 b(i,1:ngp)=dapwalm(1:ngp,io,lm)
31 end do
32 end do
33end do
34call zmctm(lmoapw(is),ngpq,ngp,a,b,ld,od)
35i=0
36do l=0,lmaxapw
37 do lm=l**2+1,(l+1)**2
38 do io=1,apword(l,is)
39 i=i+1
40 a(i,1:ngpq)=dapwalmq(1:ngpq,io,lm)
41 b(i,1:ngp)=apwalm(1:ngp,io,lm)
42 end do
43 end do
44end do
45call zmctm(lmoapw(is),ngpq,ngp,a,b,ld,od)
46end subroutine
47
subroutine dolpaa(is, ias, ngp, ngpq, apwalm, apwalmq, dapwalm, dapwalmq, ld, od)
Definition dolpaa.f90:7
integer, dimension(0:maxlapw, maxspecies) apword
Definition modmain.f90:758
integer lmaxapw
Definition modmain.f90:197
integer iasph
Definition modphonon.f90:15
subroutine zmctm(l, m, n, a, b, ld, c)
Definition zmctm.f90:7