The Elk Code
 
Loading...
Searching...
No Matches
olpaa.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
6subroutine olpaa(tor,is,ngp,apwalm,ld,o)
7use modmain
8implicit none
9! arguments
10logical, intent(in) :: tor
11integer, intent(in) :: is,ngp
12complex(8), intent(in) :: apwalm(ngkmax,apwordmax,lmmaxapw)
13integer, intent(in) :: ld
14complex(8), intent(inout) :: o(*)
15! local variables
16integer io,l,lm,i
17! automatic arrays
18complex(8) a(lmoapw(is),ngp)
19i=0
20do 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:ngp)=apwalm(1:ngp,io,lm)
25 end do
26 end do
27end do
28if (tor) then
29! matrix O is real
30 call rzmctmu(lmoapw(is),ngp,a,a,ld,o)
31else
32! matrix O is complex
33 call zmctmu(lmoapw(is),ngp,a,a,ld,o)
34end if
35end subroutine
36
integer, dimension(0:maxlapw, maxspecies) apword
Definition modmain.f90:758
integer lmaxapw
Definition modmain.f90:197
subroutine olpaa(tor, is, ngp, apwalm, ld, o)
Definition olpaa.f90:7
subroutine rzmctmu(l, n, a, b, ld, c)
Definition rzmctmu.f90:7
subroutine zmctmu(l, n, a, b, ld, c)
Definition zmctmu.f90:7