The Elk Code
 
Loading...
Searching...
No Matches
gaunt.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: gaunt
8! !INTERFACE:
9real(8) function gaunt(l1,l2,l3,m1,m2,m3)
10! !INPUT/OUTPUT PARAMETERS:
11! l1, l2, l3 : angular momentum quantum numbers (in,integer)
12! m1, m2, m3 : magnetic quantum numbers (in,integer)
13! !DESCRIPTION:
14! Returns the Gaunt coefficient given by
15! $$ \langle Y^{l_1}_{m_1}|Y^{l_2}_{m_2}|Y^{l_3}_{m_3} \rangle
16! = (-1)^{m_1}\left[\frac{(2l_1+1)(2l_2+1)(2l_3+1)}{4\pi} \right]
17! ^{\frac{1}{2}}
18! \begin{pmatrix} l_1 & l_2 & l_3 \\ 0 & 0 & 0 \end{pmatrix}
19! \begin{pmatrix} l_1 & l_2 & l_3 \\ -m_1 & m_2 & m_3 \end{pmatrix}. $$
20! Suitable for $l_i$ less than 50.
21!
22! !REVISION HISTORY:
23! Created November 2002 (JKD)
24!EOP
25!BOC
26implicit none
27! arguments
28integer, intent(in) :: l1,l2,l3
29integer, intent(in) :: m1,m2,m3
30! local variables
31integer j,j1,j2,j3,jh
32real(8) t1
33! real constant 1/sqrt(4*pi)
34real(8), parameter :: c1=0.28209479177387814347d0
35! external functions
36real(8), external :: wigner3j,factn,factr
37if ((l1 < 0).or.(l2 < 0).or.(l3 < 0).or.(abs(m1) > l1).or.(abs(m2) > l2) &
38 .or.(abs(m3) > l3)) then
39 write(*,*)
40 write(*,'("Error(gaunt): non-physical arguments :")')
41 write(*,'("l1 = ",I8," l2 = ",I8," l3 = ",I8)') l1,l2,l3
42 write(*,'("m1 = ",I8," m2 = ",I8," m3 = ",I8)') m1,m2,m3
43 write(*,*)
44 stop
45end if
46if ((l1 > 50).or.(l2 > 50).or.(l3 > 50)) then
47 write(*,*)
48 write(*,'("Error(gaunt): angular momenta out of range : ",3I8)') l1,l2,l3
49 write(*,*)
50 stop
51end if
52if (m1-m2-m3 /= 0) then
53 gaunt=0.d0
54 return
55end if
56j1=l2-l1+l3
57j2=l1-l2+l3
58j3=l1+l2-l3
59if ((j1 < 0).or.(j2 < 0).or.(j3 < 0)) then
60 gaunt=0.d0
61 return
62end if
63j=l1+l2+l3
64if (mod(j,2) /= 0) then
65 gaunt=0.d0
66 return
67end if
68jh=j/2
69t1=sqrt(dble((2*l1+1)*(2*l2+1)*(2*l3+1))*factr(j1,j+1)*factn(j2)*factn(j3))
70t1=t1*factr(jh,jh-l1)/(factn(jh-l2)*factn(jh-l3))
71gaunt=t1*c1*wigner3j(l1,l2,l3,-m1,m2,m3)
72if (mod(m1+jh,2) /= 0) gaunt=-gaunt
73end function
74!EOC
75
elemental real(8) function factn(n)
Definition factn.f90:7
real(8) function factr(n, d)
Definition factr.f90:10
real(8) function gaunt(l1, l2, l3, m1, m2, m3)
Definition gaunt.f90:10
real(8) function wigner3j(j1, j2, j3, m1, m2, m3)
Definition wigner3j.f90:10