File: wigner3j.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (89 lines) | stat: -rw-r--r-- 2,624 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89

! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
! This file is distributed under the terms of the GNU Lesser General Public
! License. See the file COPYING for license details.

!BOP
! !ROUTINE: wigner3j
! !INTERFACE:
real(8) function wigner3j(j1,j2,j3,m1,m2,m3)
! !INPUT/OUTPUT PARAMETERS:
!   j1, j2, j3 : angular momentum quantum numbers (in,integer)
!   m1, m2, m3 : magnetic quantum numbers (in,integer)
! !DESCRIPTION:
!   Returns the Wigner $3j$-symbol. There are many equivalent formulae for
!   the $3j$-symbols, the following provides high accuracy for $j\le 50$
!   \begin{align*}
!    &\begin{pmatrix} j_1 & j_2 & j_3 \\ m_1 & m_2 & m_3 \end{pmatrix}= \\
!    &(-1)^{j1+j2+m3}\sqrt{\frac{(j_1+m_1)!\,(j_2+m_2)!\,(j_3+m_3)!\,
!    (j_3-m_3)!\,(j_1-m_1)!\,(j_2-m_2)!}{(j_2-j_1+j_3)!\,(j_1-j_2+j_3)!\,
!    (j_1+j_2-j_3)!\,(1+j_1+j_2+j_3)!}}\,\sum_k(-1)^k \\
!    &\frac{(j_2-j_1+j_3)!\,(j_1-j_2+j_3)!\,(j_1+j_2-j_3)!}{(j_3-j_1-m_2+k)!\,
!    (j_3-j_2+m_1+k)!\,(j_1+j_2-j_3-k)!\,k!\,(j_1-m_1-k)!\,(j_2+m_2-k)!},
!   \end{align*}
!   where the sum is over all integers $k$ for which the factorials in the
!   summand are non-negative.
!
! !REVISION HISTORY:
!   Created November 2002 (JKD)
!EOP
!BOC
implicit none
! arguments
integer, intent(in) :: j1,j2,j3
integer, intent(in) :: m1,m2,m3
! local variables
integer k,k1,k2,l1,l2,l3,n1,n2
real(8) sgn,sum,t1
! external functions
real(8) factnm,factr
external factnm,factr
! check input variables
if ((j1.lt.0).or.(j2.lt.0).or.(j3.lt.0).or.(abs(m1).gt.j1).or.(abs(m2).gt.j2) &
 .or.(abs(m3).gt.j3)) then
  write(*,*)
  write(*,'("Error(wigner3j): invalid arguments :")')
  write(*,'("j1 = ",I8," j2 = ",I8," j3 = ",I8)') j1,j2,j3
  write(*,'("m1 = ",I8," m2 = ",I8," m3 = ",I8)') m1,m2,m3
  write(*,*)
  stop
end if
if ((j1.eq.0).and.(j2.eq.0).and.(j3.eq.0)) then
  wigner3j=1.d0
  return
end if
if ((j1.gt.50).or.(j2.gt.50).or.(j3.gt.50)) then
  write(*,*)
  write(*,'("Error(wigner3j): angular momenta out of range : ",3I8)') j1,j2,j3
  write(*,*)
  stop
end if
l1=j2-j1+j3
l2=j1-j2+j3
l3=j1+j2-j3
if ((m1+m2+m3.ne.0).or.(l1.lt.0).or.(l2.lt.0).or.(l3.lt.0)) then
  wigner3j=0.d0
  return
end if
n1=j1-m1
n2=j2+m2
k1=max(0,n1-l2,n2-l1)
k2=min(l3,n1,n2)
if (mod(k1-j1+j2+m3,2).ne.0) then
  sgn=-1.d0
else
  sgn=1.d0
end if
sum=0.d0
do k=k1,k2
  t1=sgn*factr(l1,l1-n2+k)*factr(l2,l2-n1+k)*factr(l3,l3-k)
  sum=sum+t1/(factnm(k,1)*factnm(n1-k,1)*factnm(n2-k,1))
  sgn=-sgn
end do
t1=factr(j1+m1,l1)*factr(j2+m2,l2)*factr(j3+m3,l3)
t1=t1*factr(j3-m3,1+j1+j2+j3)*factnm(j1-m1,1)*factnm(j2-m2,1)
wigner3j=sum*sqrt(t1)
return
end function
!EOC