The Elk Code
olpaaq.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 
6 subroutine olpaaq(is,ngp,ngpq,apwalm,apwalmq,ld,oq)
7 use modmain
8 implicit none
9 ! arguments
10 integer, intent(in) :: is,ngp,ngpq
11 complex(8), intent(in) :: apwalm(ngkmax,apwordmax,lmmaxapw)
12 complex(8), intent(in) :: apwalmq(ngkmax,apwordmax,lmmaxapw)
13 integer, intent(in) :: ld
14 complex(8), intent(inout) :: oq(*)
15 ! local variables
16 integer io,l,lm,i
17 ! automatic arrays
18 complex(8) a(lmoapw(is),ngpq),b(lmoapw(is),ngp)
19 i=0
20 do l=0,lmaxapw
21  do lm=l**2+1,(l+1)**2
22  do io=1,apword(l,is)
23  i=i+1
24  a(i,1:ngpq)=apwalmq(1:ngpq,io,lm)
25  b(i,1:ngp)=apwalm(1:ngp,io,lm)
26  end do
27  end do
28 end do
29 call zmctm(lmoapw(is),ngpq,ngp,a,b,ld,oq)
30 end subroutine
31 
integer lmaxapw
Definition: modmain.f90:197
subroutine olpaaq(is, ngp, ngpq, apwalm, apwalmq, ld, oq)
Definition: olpaaq.f90:7
integer, dimension(0:maxlapw, maxspecies) apword
Definition: modmain.f90:758
subroutine zmctm(l, m, n, a, b, ld, c)
Definition: zmctm.f90:7