API

Types

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.

source

Some useful function

CGcoefficient.check_jmFunction
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

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

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

source
CGcoefficient.binomial_data_sizeFunction
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.

source
CGcoefficient.check_CGFunction
check_CG(dj1::T, dj2::T, dj3::T, dm1::T, dm2::T, dm3::T) where {T <: Integer}

Check if the Clebsch-Gordan coefficient is valid.

source
CGcoefficient.check_3jFunction
check_3j(dj1::T, dj2::T, dj3::T, dm1::T, dm2::T, dm3::T) where {T <: Integer}

Check if the Wigner 3-j symbol is valid.

source
CGcoefficient.check_6jFunction
check_6j(dj1::T, dj2::T, dj3::T, dj4::T, dj5::T, dj6::T) where {T <: Integer}

Check if the Wigner 6-j symbol is valid.

source
CGcoefficient.check_9jFunction
check_9j(dj1::T, dj2::T, dj3::T, dj4::T, dj5::T, dj6::T, dj7::T, dj8::T, dj9::T) where {T <: Integer}

Check if the Wigner 9-j symbol is valid.

source
Base.floatMethod
float(x::SqrtRational)::BigFloat

Convert a SqrtRational to BigFloat.

source

Exact functions

The default functions give out exact result in the format of SqrtRational. The results are simplified to give out shotest possible result. Their arguments are Real (but only allow twice of which can be convert into integer).

CGcoefficient.CGFunction
CG(j1::Real, j2::Real, j3::Real, m1::Real, m2::Real, m3::Real)

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

source
CGcoefficient.CG0Function
CG0(j1::Integer, j2::Integer, j3::Integer)

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

source
CGcoefficient.threeJFunction
threeJ(j1::Real, j2::Real, j3::Real, m1::Real, m2::Real, m3::Real)

Wigner 3j-symbol

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

source
CGcoefficient.sixJFunction
sixJ(j1::Real, j2::Real, j3::Real, j4::Real, j5::Real, j6::Real)

Wigner 6j-symbol

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

source
CGcoefficient.nineJFunction
nineJ(j1::Real, j2::Real, j3::Real,
      j4::Real, j5::Real, j6::Real,
      j7::Real, j8::Real, j9::Real)

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}\]

source
CGcoefficient.norm9JFunction
norm9J(j1::Real, j2::Real, j3::Real,
       j4::Real, j5::Real, j6::Real,
       j7::Real, j8::Real, j9::Real)

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}\]

source
CGcoefficient.lsjjFunction
lsjj(l1::Integer, l2::Integer, j1::Real, j2::Real, L::Integer, S::Integer, J::Integer)

LS-coupling to jj-coupling transformation coefficient

\[\langle l_1l_2L,s_1s_2S;J|l_1s_1j_1,l_2s_2j_2;J\rangle = \begin{bmatrix}l_1 & s_1 & j_1 \\ l_2 & s_2 & j_2 \\ L & S & J\end{bmatrix}\]

source
CGcoefficient.RacahFunction
Racah(j1::Real, j2::Real, j3::Real, j4::Real, j5::Real, j6::Real)

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}\]

source
CGcoefficient.MoshinskyFunction
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. This function is designed for demonstration the exact result, It only calculate the case of $\tan(\beta) = \sqrt{m_1\omega_1/m_2\omega_2} = 1$.

source

People often use double of angular momentum quantum number as parameters, so we can use integer as parameters. This package also offers such functions, where the d letter means double. These functions also give out exact SqrtRational results.

CGcoefficient.dCGFunction
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)

source
CGcoefficient.d3jFunction
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.

source
CGcoefficient.d6jFunction
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.

source
CGcoefficient.d9jFunction
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.

source
CGcoefficient.dRacahFunction
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.

source

float version functions

Float version functions is always used for numeric calculation. They are designed for fast calculation. You should call wigner_init_float to reserve the inner binomial table. They only resive Integer arguments, thus fCG, f3j, f6j, fRacha, f9j only resive arguements which are double of the exact quantum number. In some functions, e.g. flsjj, some quantum number can be half-integer, the corresponding arguments start with d (which mean double of the quantum number); while some other quantum number must be integer, the corresponding arguments does not starts with d.

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

We calculate the Wigner Symbols with binomials. We will store some binomials first, then when we need one binomial, we just load it. In this package, the functions starts with f character dependent on the stored binomials. The wigner_init_float function is used to extend the range of the stored binomial table.

The parameters means:

  • mode:
    • "Jmax": set the maximum angular momentum in the system.
    • "2bjmax": set the maximum two-body coupled angular momentum in the system.
    • "nmax": directly set the maximum binomial number.
    • "Moshinsky": set the maximum N = 2n+l in the Moshinsky bracket.
  • n: the maximum angular momentum or the maximum binomial number, dependent on the mode.
  • rank:
    • 3: calculate only CG coefficients and 3j symbols.
    • 6: calculate CG coefficients, 3j symbols and 6j symbols.
    • 9: calculate CG coefficients, 3j symbols, 6j symbols and 9j symbols.

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

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.

The "Moshinsky" mode means you want to calculate the Moshinsky brackets. The n parameter is the maximum N = 2n+l in the Moshinsky bracket. The rank parameter is ignored.

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

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.

The exact values of the maximum 'nmaxfor differentrank` are shown in the table below (details).

Calculate rangeCG & 3j6j & Racah9j
meaning of typetype\rank369
max angular momentum"Jmax"3Jmax+14Jmax+15Jmax+1
max two-body coupled angular momentum"2bjmax"2jmax+13jmax+14jmax+1
max $N = 2n+l$ for Moshinsky bracket"Moshinsky"6Nmax+16Nmax+16Nmax+1
max binomial"nmax"nmaxnamxnmax

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 wigner_init_float function is not thread safe, so you should call it before you start your calculation.

source
CGcoefficient.unsafe_fbinomialFunction
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.

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

float64 and fast CG coefficient.

source
CGcoefficient.fCG0Function
fCG0(dj1::Integer, dj2::Integer, dj3::Integer)

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

source
CGcoefficient.fCGspinFunction
fCGspin(ds1::Integer, ds2::Integer, S::Integer)

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

source
CGcoefficient.fCG3spinFunction
fCG3spin(ds1::Integer, ds2::Integer, ds3::Integer, S12::Integer, dS::Integer)

float64 and fast CG coefficient for three spin-1/2 coupling: <S12,M12|1/2,m1;1/2,m2><S,M|S12,M12;1/2,m3>

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

float64 and fast Wigner 3j symbol.

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

float64 and fast Wigner 6j symbol.

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

float64 and fast Racah coefficient.

source
CGcoefficient.f9jFunction
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.

source
CGcoefficient.fnorm9jFunction
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.

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

float64 and fast lsjj coefficient.

source
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.

source

prime factorization version

I implement a primary version of the prime factorization algorithm, based on the wigxjpf library.

CGcoefficient.PFRationalType
PFRational{T<:Integer} <: Real

A positive rational number in prime factorization form.

\[q = \prod_{i=1}^{n} p_i^{e_i}\]

where $p_i$ is the $i$-th prime number and $e_i$ can be negative.

source
Base.gcdMethod
gcd(x::PFRational, y::PFRational)

gcd of rational means gcd of numerator and lcm of denominator

source
Base.lcmMethod
lcm(x::PFRational, y::PFRational)

lcm of rational means lcm of numerator and gcd of denominator

source
CGcoefficient.wigner_init_pfFunction
wigner_init_pf(n::Integer, mode::AbstractString, rank::Integer)

Initialize the prime factorization tables for Wigner symbols. The parameters are similar to wigner_inif_float. You must call this function before using the e and ef versions of Wigner symbols.

source
CGcoefficient.pf_binomialFunction
pf_binomial(n::Integer, k::Integer)::PFRational{Int16}

Return the stored binomial coefficient C(n, k) in prime factorization form.

source
CGcoefficient.eCGFunction
eCG(dj1::Integer, dj2::Integer, dj3::Integer, dm1::Integer, dm2::Integer, dm3::Integer)::SqrtRational{BigInt}

Exact CG coefficient using prime factorization algorithm.

source
CGcoefficient.eCG0Function
eCG0(j1::Integer, j2::Integer, j3::Integer)::SqrtRational{BigInt}

Exact CG coefficient (m1 = m2 = m3 = 0) using prime factorization algorithm.

source
CGcoefficient.e3jFunction
e3j(dj1::Integer, dj2::Integer, dj3::Integer, dm1::Integer, dm2::Integer, dm3::Integer)::SqrtRational{BigInt}

Exact 3j symbol using prime factorization algorithm.

source
CGcoefficient.e6jFunction
e6j(j1::Integer, j2::Integer, j3::Integer, j4::Integer, j5::Integer, j6::Integer)::SqrtRational{BigInt}

Exact 6j symbol using prime factorization algorithm.

source
CGcoefficient.efCGFunction
efCG(dj1::Integer, dj2::Integer, dj3::Integer, dm1::Integer, dm2::Integer, dm3::Integer)::Float64

Exact CG coefficient (Float64) using prime factorization algorithm.

source
CGcoefficient.efCG0Function
efCG0(j1::Integer, j2::Integer, j3::Integer)::Float64

Exact CG coefficient (m1 = m2 = m3 = 0) (Float64) using prime factorization algorithm.

source
CGcoefficient.ef3jFunction
ef3j(dj1::Integer, dj2::Integer, dj3::Integer, dm1::Integer, dm2::Integer, dm3::Integer)::Float64

Exact 3j symbol (Float64) using prime factorization algorithm.

source
CGcoefficient.ef6jFunction
ef6j(j1::Integer, j2::Integer, j3::Integer, j4::Integer, j5::Integer, j6::Integer)::Float64

Exact 6j symbol (Float64) using prime factorization algorithm.

source