The Elk Code
 
Loading...
Searching...
No Matches
clebgor.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 Lesser General Public
4! License. See the file COPYING for license details.
5
6!BOP
7! !ROUTINE: clebgor
8! !INTERFACE:
9real(8) function clebgor(j1,j2,j3,m1,m2,m3)
10! !INPUT/OUTPUT PARAMETERS:
11! j1, j2, j3 : angular momentum quantum numbers (in,integer)
12! m1, m2, m3 : magnetic quantum numbers (in,integer)
13! !DESCRIPTION:
14! Returns the Clebsch-Gordon coefficients using the Wigner $3j$-symbols
15! $$ C(J_1 J_2 J_3 | m_1 m_2 m_3)=(-1)^{J_1-J_2+m_3}\sqrt{2J_3+1}
16! \begin{pmatrix} J_1 & J_2 & J_3 \\ m_1 & m_2 & -m_3 \end{pmatrix}. $$
17! Suitable for $J_i\le 50$. See {\tt wigner3j}.
18!
19! !REVISION HISTORY:
20! Created September 2003 (JKD)
21!EOP
22!BOC
23implicit none
24! arguments
25integer, intent(in) :: j1,j2,j3
26integer, intent(in) :: m1,m2,m3
27! external functions
28real(8), external :: wigner3j
29if ((j1 < 0).or.(j2 < 0).or.(j3 < 0).or.(abs(m1) > j1).or.(abs(m2) > j2) &
30 .or.(abs(m3) > j3)) then
31 write(*,*)
32 write(*,'("Error(clebgor): non-physical arguments :")')
33 write(*,'("j1 = ",I8," j2 = ",I8," j3 = ",I8)') j1,j2,j3
34 write(*,'("m1 = ",I8," m2 = ",I8," m3 = ",I8)') m1,m2,m3
35 write(*,*)
36 stop
37end if
38if ((j1 == 0).and.(j2 == 0).and.(j3 == 0)) then
39 clebgor=1.d0
40 return
41end if
42if ((j1 > 50).or.(j2 > 50).or.(j3 > 50)) then
43 write(*,*)
44 write(*,'("Error(clebgor): angular momenta out of range : ",3I8)') j1,j2,j3
45 write(*,*)
46 stop
47end if
48if ((m1+m2 /= m3).or.(j1+j2 < j3).or.(j2+j3 < j1).or.(j1+j3 < j2)) then
49 clebgor=0.d0
50 return
51end if
52clebgor=sqrt(dble(2*j3+1))*wigner3j(j1,j2,j3,m1,m2,-m3)
53if (mod(j1-j2+m3,2) /= 0) clebgor=-clebgor
54end function
55!EOC
56
real(8) function clebgor(j1, j2, j3, m1, m2, m3)
Definition clebgor.f90:10
real(8) function wigner3j(j1, j2, j3, m1, m2, m3)
Definition wigner3j.f90:10