The Elk Code
 
Loading...
Searching...
No Matches
tm3rtoz.f90
Go to the documentation of this file.
1
2! Copyright (C) 2022 Leon Kerber, 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
6subroutine tm3rtoz(l,k,p,r,ld,wkpr,zkpr)
7use modmain
8implicit none
9! arguments
10integer, intent(in) :: l,k,p,r,ld
11real(8), intent(in) :: wkpr(-ld:ld)
12complex(8), intent(out) :: zkpr(-ld:ld)
13! local variables
14integer g,t
15real(8) a,b,t0,t1,t2
16! external functions
17real(8), external :: wigner3j,factn,factn2,factr
18! old convention normalisation factors and phase factors
19g=k+p+r
20if (mod(g,2) == 0) then
21 t0=1.d0/wigner3j(k,p,r,0,0,0)
22else
23 t0=sqrt(factr(g+1,g-2*k)/(factn(g-2*p)*factn(g-2*r)))
24 t0=t0*factn2(g-2*k)*factn2(g-2*p)*factn2(g-2*r)/factn2(g)
25end if
26t0=t0/sqrt(dble(2*r+1))
27t0=t0*sqrt(factn(2*l-k)*factn(2*l+k+1))/factn(2*l)
28t0=t0*sqrt(factn(2+p))
29if (mod(k+p,2) /= 0) t0=-t0
30! remove orthonormal convention normalisation factors
31t0=t0/sqrt(dble(2*k+1))
32t0=t0/sqrt(dble(2*p+1))
33t0=t0/2.d0
34do t=-r,r
35 t1=t0*(wkpr(t)+wkpr(-t))
36 t2=t0*(wkpr(t)-wkpr(-t))
37 if (mod(t,2) == 0) then
38 a=t1
39 b=t2
40 else
41 a=-t2
42 b=-t1
43 end if
44 if ((k == r).and.(p == 1)) then
45 if (mod(k,2) == 0) then
46 b=-b
47 else
48 a=-a
49 end if
50 end if
51 zkpr(t)=cmplx(a,b,8)
52end do
53end subroutine
54
elemental real(8) function factn2(n)
Definition factn2.f90:7
elemental real(8) function factn(n)
Definition factn.f90:7
real(8) function factr(n, d)
Definition factr.f90:10
subroutine tm3rtoz(l, k, p, r, ld, wkpr, zkpr)
Definition tm3rtoz.f90:7
real(8) function wigner3j(j1, j2, j3, m1, m2, m3)
Definition wigner3j.f90:10