CGcoefficient.HalfIntType
HalfInt = Union{Integer, Rational}

Angular momentum quantum numbers may be half integers and integers. With HalfInt type, you can use both 2 and 3//2 as parameters. But the parameter like 5//3, who's denominator is not 2 while gives out error.

CGcoefficient.SqrtRationalType
SqrtRational <: Real

Store exact value of a Rational's square root. It is stored in s√r format, and we do not simplify it during arithmetics. You can use the simplify function to simplify it.

Base.floatMethod
float(x::SqrtRational)::Float64

Although we use BigInt in calculation, we do not convert it to BigFloat. Because this package is designed for numeric calculation, BigFloat is unnecessary.

CGcoefficient.CGMethod
CG(j1::HalfInt, j2::HalfInt, j3::HalfInt, m1::HalfInt, m2::HalfInt, m3::HalfInt)

CG coefficient $\langle j_1m_1 j_2m_2 | j_3m_3 \rangle$

CGcoefficient.CG0Method
CG0(j1::Integer, j2::Integer, j3::Integer)

CG coefficient special case: $\langle j_1 0 j_2 0 | j_3 0 \rangle$.

CGcoefficient.MoshinskyMethod
Moshinsky(N::Integer, L::Integer, n::Integer, l::Integer, n1::Integer, l1::Integer, n2::Integer, l2::Integer, Λ::Integer)

Moshinsky bracket,Ref: Buck et al. Nuc. Phys. A 600 (1996) 387-402. Since this function is designed for demonstration the exact result, It only calculate the case of $\tan(\beta) = 1$.

CGcoefficient.RacahMethod
Racah(j1::HalfInt, j2::HalfInt, j3::HalfInt, j4::HalfInt, j5::HalfInt, j6::HalfInt)

Racah coefficient

\[W(j_1j_2j_3j_4, j_5j_6) = (-1)^{j_1+j_2+j_3+j_4} \begin{Bmatrix} j_1 & j_2 & j_5 \\ j_4 & j_3 & j_6 \end{Bmatrix}\]

CGcoefficient.binomial_data_sizeMethod
binomial_data_size(n::Int)::Int

Store (half) binomial data in the order of

# n - data
  0   1
  1   1
  2   1 2
  3   1 3
  4   1 4 6
  5   1 5 10
  6   1 6 15 20

Return the total number of the stored data.

CGcoefficient.check_coupleMethod
check_couple(dj1::T, dj2::T, dj3::T) where {T <: Integer}

check if three angular monentum number j1, j2, j3 can couple

CGcoefficient.check_jmMethod
check_jm(dj::T, dm::T) where {T <: Integer}

check if the m-quantum number if one of the components of the j-quantum number, in other words, m and j has the same parity, and abs(m) < j

CGcoefficient.d3jMethod
d3j(dj1::Integer, dj2::Integer, dj3::Integer, dm1::Integer, dm2::Integer, dm3::Integer)

3j-symbol function with double angular monentum number parameters, so that the parameters can be integer.

CGcoefficient.d6jMethod
d6j(dj1::Integer, dj2::Integer, dj3::Integer, dj4::Integer, dj5::Integer, dj6::Integer)

6j-symbol function with double angular monentum number parameters, so that the parameters can be integer.

CGcoefficient.d9jMethod
d9j(dj1::Integer, dj2::Integer, dj3::Integer,
    dj4::Integer, dj5::Integer, dj6::Integer,
    dj7::Integer, dj8::Integer, dj9::Integer)

9j-symbol function with double angular monentum number parameters, so that the parameters can be integer.

CGcoefficient.dCGMethod
dCG(dj1::Integer, dj2::Integer, dj3::Integer, dm1::Integer, dm2::Integer, dm3::Integer)

CG coefficient function with double angular monentum number parameters, so that the parameters can be integer. You can calculate dCG(1, 1, 2, 1, 1, 2) to calculate the real CG(1//2, 1//2, 1, 1/2, 1//2, 1)

CGcoefficient.dRacahMethod
dRacah(dj1::Integer, dj2::Integer, dj3::Integer, dj4::Integer, dj5::Integer, dj6::Integer)

Racah coefficient function with double angular momentum parameters, so that the parameters can be integer.

CGcoefficient.dfuncMethod
dfunc(dj::Integer, dm1::Integer, dm2::Integer, β::Float64)

Wigner d-function.

CGcoefficient.f3jMethod
f3j(dj1::Integer, dj2::Integer, dj3::Integer, dm1::Integer, dm2::Integer, dm3::Integer)

float64 and fast Wigner 3j symbol.

CGcoefficient.f6jMethod
f6j(dj1::Integer, dj2::Integer, dj3::Integer, dj4::Integer, dj5::Integer, dj6::Integer)

float64 and fast Wigner 6j symbol.

CGcoefficient.f9jMethod
f9j(dj1::Integer, dj2::Integer, dj3::Integer,
    dj4::Integer, dj5::Integer, dj6::Integer,
    dj7::Integer, dj8::Integer, dj9::Integer)

float64 and fast Wigner 9j symbol.

CGcoefficient.fCGMethod
fCG(dj1::Integer, dj2::Integer, dj3::Integer, dm1::Integer, dm2::Integer, dm3::Integer)

float64 and fast CG coefficient.

CGcoefficient.fCG0Method
fCG0(dj1::Integer, dj2::Integer, dj3::Integer)

float64 and fast CG coefficient for m1 == m2 == m3 == 0.

CGcoefficient.fCGspinMethod
fCGspin(ds1::Integer, ds2::Integer, S::Integer)

float64 and fast CG coefficient for two spin-1/2 coupling.

CGcoefficient.fMoshinskyFunction
fMoshinsky(N::Integer, L::Integer, n::Integer, l::Integer, n1::Integer, l1::Integer, n2::Integer, l2::Integer, tanβ::Float64 = 1.0)

float64 and fast Moshinsky bracket.

CGcoefficient.fRacahMethod
Racah(dj1::Integer, dj2::Integer, dj3::Integer, dj4::Integer, dj5::Integer, dj6::Integer)

float64 and fast Racah coefficient.

CGcoefficient.flsjjMethod
lsjj(l1::Integer, l2::Integer, dj1::Integer, dj2::Integer, L::Integer, S::Integer, J::Integer)

float64 and fast lsjj coefficient.

CGcoefficient.fnorm9jMethod
fnorm9j(dj1::Integer, dj2::Integer, dj3::Integer,
        dj4::Integer, dj5::Integer, dj6::Integer,
        dj7::Integer, dj8::Integer, dj9::Integer)

float64 and fast normalized Wigner 9j symbol.

CGcoefficient.lsjjMethod
lsjj(l1::Integer, l2::Integer, j1::HalfInt, j2::HalfInt, L::Integer, S::Integer, J::Integer)

LS-coupling to jj-coupling transformation coefficient

\[|l_1 l_2 j_1 j_2; J\rangle = \sum_{LS} \langle l_1 l_2 LSJ | l_1 l_2 j_1 j_2; J \rangle |l_1 l_2 LSJ \rangle\]

CGcoefficient.nineJMethod
nineJ(j1::HalfInt, j2::HalfInt, j3::HalfInt,
      j4::HalfInt, j5::HalfInt, j6::HalfInt,
      j7::HalfInt, j8::HalfInt, j9::HalfInt)

Wigner 9j-symbol

\[\begin{Bmatrix} j_1 & j_2 & j_3 \\ j_4 & j_5 & j_6 \\ j_7 & j_8 & j_9 \end{Bmatrix}\]

CGcoefficient.norm9JMethod
norm9J(j1::HalfInt, j2::HalfInt, j3::HalfInt,
       j4::HalfInt, j5::HalfInt, j6::HalfInt,
       j7::HalfInt, j8::HalfInt, j9::HalfInt)

normalized Wigner 9j-symbol

\[\begin{bmatrix} j_1 & j_2 & j_3 \\ j_4 & j_5 & j_6 \\ j_7 & j_8 & j_9 \end{bmatrix} = \sqrt{(2j_3+1)(2j_6+1)(2j_7+1)(2j_8+1)} \begin{Bmatrix} j_1 & j_2 & j_3 \\ j_4 & j_5 & j_6 \\ j_7 & j_8 & j_9 \end{Bmatrix}\]

CGcoefficient.sixJMethod
sixJ(j1::HalfInt, j2::HalfInt, j3::HalfInt, j4::HalfInt, j5::HalfInt, j6::HalfInt)

Wigner 6j-symbol

\[\begin{Bmatrix} j_1 & j_2 & j_3 \\ j_4 & j_5 & j_6 \end{Bmatrix}\]

CGcoefficient.threeJMethod
threeJ(j1::HalfInt, j2::HalfInt, j3::HalfInt, m1::HalfInt, m2::HalfInt, m3::HalfInt)

Wigner 3j-symbol

\[\begin{pmatrix} j_1 & j_2 & j_3 \\ m_1 & m_2 & m_3 \end{pmatrix}\]

CGcoefficient.unsafe_fbinomialMethod
function unsafe_fbinomial(n::Int, k::Int)

Same as fbinomial, but without bounds check. Thus for n < 0 or n > _nmax or k < 0 or k > n, the result is undefined. In Wigner symbol calculation, the mathematics guarantee that unsafe_fbinomial(n, k) is always safe.

CGcoefficient.wigner_init_floatMethod
wigner_init_float(n::Integer, mode::AbstractString, rank::Integer)

This function reserves memory for fbinomial(n, k). In this code, the fbinomial function is only valid in the stored range. If you call a fbinomial function out of the range, it just gives you 0. The __init__() function stores to nmax = 67.

The parameters means

Calculate rangeCG & 3j6j & Racah9j
meaning of typetype\rank369
max angular momentum"Jmax"3*Jmax+14*Jmax+15*Jmax+1
max two-body coupled angular momentum"2bjmax"2*jmax+13*jmax+14*jmax+1
max binomial"nmax"nmaxnamxnmax

The "2bjmax" mode means your calculation only consider two-body coupling, and no three-body coupling. This mode assumes that in all these coefficients, at least one of the angular momentun is just a single particle angular momentum. With this assumption, "2bjmax" mode will use less memory than "Jmax" mode.

"Jmax" means the global maximum angular momentum, for every parameters. It is always safe with out any assumption.

The "nmax" mode directly set nmax, and the rank parameter is ignored.

rank = 6 means you only need to calculate CG and/or 6j symbols, you don't need to calculate 9j symbol.

For example

wigner_init_float(21, "Jmax", 6)

means you calculate CG and 6j symbols, and donot calculate 9j symbol. The maximum angular momentum in your system is 21.

You do not need to rememmber those values in the table. You just need to find the maximum angular momentum in you canculation, then call the function.

The wignerinitfloat function is not thread safe, so you should call it before you start your calculation.