The Elk Code
 
Loading...
Searching...
No Matches
lopzflmn.f90
Go to the documentation of this file.
1
2! Copyright (C) 2020 J. K. Dewhurst and S. Sharma.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6pure subroutine lopzflmn(lmax,n,ld,zflm,zlflm1,zlflm2,zlflm3)
7implicit none
8! arguments
9integer, intent(in) :: lmax,n,ld
10complex(8), intent(in) :: zflm(ld,n)
11complex(8), intent(out) :: zlflm1(ld,n),zlflm2(ld,n),zlflm3(ld,n)
12! local variables
13integer l,m,lm,i
14real(8) t1
15complex(8), parameter :: zi=(0.d0,1.d0),zmi=(0.d0,-1.d0)
16complex(8) z1
17lm=0
18do l=0,lmax
19 do m=-l,l
20 lm=lm+1
21 if (m == -l) then
22 zlflm1(lm,1:n)=0.d0
23 zlflm2(lm,1:n)=0.d0
24 end if
25 if (m < l) then
26 t1=0.5d0*sqrt(dble((l-m)*(l+m+1)))
27 do i=1,n
28 z1=t1*zflm(lm,i)
29 zlflm1(lm+1,i)=z1
30 zlflm2(lm+1,i)=zmi*z1
31 end do
32 end if
33 if (m > -l) then
34 t1=0.5d0*sqrt(dble((l+m)*(l-m+1)))
35 do i=1,n
36 z1=t1*zflm(lm,i)
37 zlflm1(lm-1,i)=zlflm1(lm-1,i)+z1
38 zlflm2(lm-1,i)=zlflm2(lm-1,i)+zi*z1
39 end do
40 end if
41 if (m /= 0) then
42 zlflm3(lm,1:n)=dble(m)*zflm(lm,1:n)
43 else
44 zlflm3(lm,1:n)=0.d0
45 end if
46 end do
47end do
48end subroutine
49
pure subroutine lopzflmn(lmax, n, ld, zflm, zlflm1, zlflm2, zlflm3)
Definition lopzflmn.f90:7