The Elk Code
gauntyry.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: gauntyry
8 ! !INTERFACE:
9 complex(8) function gauntyry(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 complex Gaunt-like coefficient given by
15 ! $\langle Y^{l_1}_{m_1}|R^{l_2}_{m_2}|Y^{l_3}_{m_3}\rangle$, where $Y_{lm}$
16 ! and $R_{lm}$ are the complex and real spherical harmonics, respectively.
17 ! Suitable for $l_i$ less than 50. See routine {\tt genrlm}.
18 !
19 ! !REVISION HISTORY:
20 ! Created November 2002 (JKD)
21 !EOP
22 !BOC
23 implicit none
24 ! arguments
25 integer, intent(in) :: l1,l2,l3
26 integer, intent(in) :: m1,m2,m3
27 ! local variables
28 ! real constant sqrt(2)/2
29 real(8), parameter :: c1=0.7071067811865475244d0
30 real(8) t1
31 ! external functions
32 real(8), external :: gaunt
33 if (m2 > 0) then
34  if (mod(m2,2) == 0) then
35  t1=c1*(gaunt(l1,l2,l3,m1,m2,m3)+gaunt(l1,l2,l3,m1,-m2,m3))
36  else
37  t1=c1*(gaunt(l1,l2,l3,m1,m2,m3)-gaunt(l1,l2,l3,m1,-m2,m3))
38  end if
39  gauntyry=t1
40 else if (m2 < 0) then
41  if (mod(m2,2) == 0) then
42  t1=c1*(gaunt(l1,l2,l3,m1,m2,m3)-gaunt(l1,l2,l3,m1,-m2,m3))
43  else
44  t1=c1*(gaunt(l1,l2,l3,m1,m2,m3)+gaunt(l1,l2,l3,m1,-m2,m3))
45  end if
46  gauntyry=cmplx(0.d0,-t1,8)
47 else
48  gauntyry=gaunt(l1,l2,l3,m1,m2,m3)
49 end if
50 end function
51 !EOC
52 
real(8) function gaunt(l1, l2, l3, m1, m2, m3)
Definition: gaunt.f90:10
complex(8) function gauntyry(l1, l2, l3, m1, m2, m3)
Definition: gauntyry.f90:10