The Elk Code
 
Loading...
Searching...
No Matches
wigner3jf.f90
Go to the documentation of this file.
1
2! Copyright (C) 2014 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!BOP
7! !ROUTINE: wigner3jf
8! !INTERFACE:
9real(8) function wigner3jf(j12,j22,j32,m12,m22,m32)
10! !INPUT/OUTPUT PARAMETERS:
11! j12, j22, j32 : angular momentum quantum numbers times 2 (in,integer)
12! m12, m22, m32 : magnetic quantum numbers times 2 (in,integer)
13! !DESCRIPTION:
14! Returns the Wigner $3j$-symbol for the case where the arguments may be
15! fractional, i.e. multiples of $\frac{1}{2}$. The input parameters to this
16! function are taken to be twice their actual values, which allows them to
17! remain integers. The formula used is identical to that in {\tt wigner3j}.
18!
19! !REVISION HISTORY:
20! Created January 2014 (JKD)
21!EOP
22!BOC
23implicit none
24! arguments
25integer, intent(in) :: j12,j22,j32
26integer, intent(in) :: m12,m22,m32
27! local variables
28integer jm1,jm2,jm3,n1,n2
29integer l12,l22,l32,l42
30integer k,k1,k2,l1,l2,l3
31real(8) sgn,sm,t1
32! external functions
33real(8), external :: factn,factr
34! check input variables
35if ((j12 < 0).or.(j22 < 0).or.(j32 < 0).or.(abs(m12) > j12).or. &
36 (abs(m22) > j22).or.(abs(m32) > j32)) then
37 write(*,*)
38 write(*,'("Error(wigner3jf): invalid arguments :")')
39 write(*,'("j12 = ",I8," j22 = ",I8," j32 = ",I8)') j12,j22,j32
40 write(*,'("m12 = ",I8," m22 = ",I8," m32 = ",I8)') m12,m22,m32
41 write(*,*)
42 stop
43end if
44if ((j12 == 0).and.(j22 == 0).and.(j32 == 0)) then
45 wigner3jf=1.d0
46 return
47end if
48if ((j12 > 100).or.(j22 > 100).or.(j32 > 100)) then
49 write(*,*)
50 write(*,'("Error(wigner3jf): angular momenta out of range : ",3I8)') j12, &
51 j22,j32
52 write(*,*)
53 stop
54end if
55jm1=j12+m12
56jm2=j22+m22
57jm3=j32+m32
58if ((mod(jm1,2) /= 0).or.(mod(jm2,2) /= 0).or.(mod(jm3,2) /= 0)) then
59 wigner3jf=0.d0
60 return
61end if
62l12=j22-j12+j32
63l22=j12-j22+j32
64l32=j12+j22-j32
65l42=j12+j22+j32
66if ((mod(l12,2) /= 0).or.(mod(l22,2) /= 0).or.(mod(l32,2) /= 0).or. &
67 (mod(l42,2) /= 0)) then
68 wigner3jf=0.d0
69 return
70end if
71l1=l12/2
72l2=l22/2
73l3=l32/2
74if ((m12+m22+m32 /= 0).or.(l1 < 0).or.(l2 < 0).or.(l3 < 0)) then
75 wigner3jf=0.d0
76 return
77end if
78n1=(j12-m12)/2
79n2=(j22+m22)/2
80k1=max(0,n1-l2,n2-l1)
81k2=min(l3,n1,n2)
82k=k1+(j22-j12+m32)/2
83if (mod(k,2) /= 0) then
84 sgn=-1.d0
85else
86 sgn=1.d0
87end if
88sm=0.d0
89do k=k1,k2
90 t1=sgn*factr(l1,l1-n2+k)*factr(l2,l2-n1+k)*factr(l3,l3-k)
91 sm=sm+t1/(factn(k)*factn(n1-k)*factn(n2-k))
92 sgn=-sgn
93end do
94jm1=jm1/2
95jm2=jm2/2
96jm3=jm3/2
97t1=factr(jm1,l1)*factr(jm2,l2)*factr(jm3,l3)
98jm1=(j12-m12)/2
99jm2=(j22-m22)/2
100jm3=(j32-m32)/2
101t1=t1*factr(jm3,1+l42/2)*factn(jm1)*factn(jm2)
102wigner3jf=sm*sqrt(t1)
103end function
104!EOC
105
elemental real(8) function factn(n)
Definition factn.f90:7
real(8) function factr(n, d)
Definition factr.f90:10
real(8) function wigner3jf(j12, j22, j32, m12, m22, m32)
Definition wigner3jf.f90:10