Quantum mechanics

General functions used in quantum mechanics.

nexus.KroneckerDelta(i, j)

Evaluates the Kronecker delta function

\[\begin{split}\delta_{ij} = \begin{cases} 0 &\text{if } i \neq j, \\ 1 &\text{if } i=j. \end{cases}\end{split}\]
Parameters:
  • i (int) – Integer number i.

  • j (int) – Integer number j.

Returns:

Kronecker delta.

Return type:

bool

nexus.Factorial(*args)

Overload 1:

Calculates the factorial \(n!\) for n <= 20.

Parameters:

n (int) – Integer number n.

Returns:

Factorial of n.

Return type:

int


Overload 2:

Calculates the factorial \(n!\) for n <= 20.

Parameters:

n (float) – Integer number. The input is rounded and converted to an integer.

Returns:

Factorial of n.

Return type:

int

Example

fac = nx.Factorial(3)

fac = nx.Factorial(3.0)
nexus.ClebshGordon(j1, j2, m1, m2, J, M)

Calculates the Clebsh-Gordon coefficient to the given angular momentum quantum numbers via Racah’s formula Eq. (3.19), [Rose]. The input follows the convention:

\[<j_1, m_1, j_2, m_2 | J, M> = C(j_1, j_2, J; m_1, m_2, M) = C(j_1, j_2, J; m_1, M-m_1),\]

where C(…) is the notation of [Rose] for the Clebsh-Gordon coefficients.

Parameters:
  • j1 (float) – j1

  • j2 (float) – j2

  • m1 (float) – m1

  • m2 (float) – m2

  • J (float) – J

  • M (float) – M

Returns:

Clebsh-Gordon coefficient.

Return type:

float

Example

clebshgordon = nx.ClebshGordon(1/2, 1, -1/2, 1, 3/2, 1/2)
nexus.ThreeJSymbol(j1, j2, j3, m1, m2, m3)

Calculates the 3j-symbol via

\[\begin{split}\begin{Bmatrix} j_1 & j_2 & j_3 \\ m_1 & m_2 & m_3 \end{Bmatrix} = \frac{(-1)^{j_1 -j_2 - m_3}}{\sqrt{2 j_3 + 1}} <j_1, j_2, m_1, m_2 | j_3, -m_3>.\end{split}\]
Parameters:
  • j1 (float) – j1

  • j2 (float) – j2

  • j3 (float) – j3

  • m1 (float) – m1

  • m2 (float) – m2

  • m3 (float) – m3

Returns:

ThreeJSymbol.

Return type:

float

nexus.VCoefficient(j1, j2, j3, m1, m2, m3)

Calculates the Racah V-coefficient after Eq. (3.20), [Rose].

\[V(j_1, j_2, j_3, m_1, m_2, m_3) = (-1)^{-j_3 + m_3} \sqrt{2j_3+1} <j_1, j_2, m_1, m_2 | j_3, -m_3>.\]
Parameters:
  • j1 (float) – j1

  • j2 (float) – j2

  • j3 (float) – j3

  • m1 (float) – m1

  • m2 (float) – m2

  • m3 (float) – m3

Returns:

V-coefficient.

Return type:

float

nexus.TriangleCoefficient(a, b, c)

Calculates the triangle coefficient after Eq. (6.8), [Rose].

\[\Delta (abc) = \sqrt{\frac{(a+b-c)! (a-b+c)! (-a+b+c)!}{(a+b+c+1)!}}.\]
Parameters:
  • a (float) – a

  • b (float) – b

  • c (float) – c

Returns:

triangle coefficient \(\Delta\).

Return type:

float

nexus.WCoefficient(a, b, c, d, e, f)

Calculates the Racah W coefficient \(W(abcd ; ef)\) after [Rose].

Parameters:
  • a (float) – a

  • b (float) – b

  • c (float) – c

  • d (float) – d

  • e (float) – e

  • f (float) – f

Returns:

Racah W coefficient \(W(abcd ; ef)\).

Return type:

float

nexus.SixJSymbol(a, b, c, d, e, f)

Calculates the 6j-symbol after Eq. (9.11), [Varshalovich].

\[\begin{split}\begin{Bmatrix} a & b & c \\ d & e & f \end{Bmatrix} = (-1)^{a+b+d+e} W(abed;cf).\end{split}\]
Parameters:
  • a (float) – a

  • b (float) – b

  • c (float) – c

  • d (float) – d

  • e (float) – e

  • f (float) – f

Returns:

6j-symbol.

Return type:

float

nexus.SmallWignerDmatrix(j, m1, m2, beta)

Calculates the small Wigner D matrix elements \(d^j_{m',m}(\beta) = d^j_{m1,m2}(\beta)\) after Eq. (4.13), [Rose].

Parameters:
  • j (float) – j

  • m1 (float) – m1

  • m2 (float) – m2

  • beta (float) – \(\beta\) angle (rad).

Returns:

Small Wigner D matrix element.

Return type:

float

Example

matrix_element = nx.SmallWignerDmatrix(1, 1, 1, nx.pi/2.5)
nexus.WignerDmatrix(j, m1, m2, alpha, beta, gamma)

Calculates the Wigner D matrix elements \(D^j_{m',m} (\alpha, \beta, \gamma) = D^j_{m1,m2} (\alpha, \beta, \gamma)\) after Eq. (4.12), [Rose].

Parameters:
  • j (float) – j

  • m1 (float) – m1

  • m2 (float) – m2

  • alpha (float) – \(\alpha\) angle (rad).

  • beta (float) – \(\beta\) angle (rad).

  • gamma (float) – \(\gamma\) angle (rad).

Returns:

Wigner D matrix element.

Return type:

complex

Example

matrix_element = nx.WignerDmatrix(1, 1, 1, nx.pi/2.5, nx.pi/1.5, nx.pi/2)
nexus.Hermite(n, x)

Computes the physicist’s Hermite polynomials of the degree n and argument x.

\[P_l^m(x) = (-1)^n e^{x^2} \frac{d^n}{dx^n} e^{-x^2}.\]
Parameters:
  • n (int) – degree \(n \geq 0\).

  • x (float) – argument x.

Returns:

physicist’s Hermite polynomials.

Return type:

double

nexus.Laguerre(n, x)

Computes the Laguerre polynomials of the degree n and argument x.

\[L_n(x) = \frac{e^x}{n!} \frac{d^n}{dx^n} \left(e^{-x} x^n \right).\]
Parameters:
  • n (int) – degree \(n \geq 0\).

  • x (float) – argument x.

Returns:

Laguerre polynomials.

Return type:

double

nexus.AssociatedLaguerre(n, k, x)

Computes the associated Laguerre polynomials of the degree n, order k and argument x.

\[L_n^k(x) = (-1)^k \frac{d^k}{dx^k} L_{n+k}(x).\]
Parameters:
  • n (int) – degree \(n \geq 0\).

  • k (int) – order \(k \geq 0\).

  • x (float) – argument x.

Returns:

associated Laguerre polynomials.

Return type:

double

nexus.Legendre(n, x)

Computes the unassociated Legendre polynomials of the degree n and argument x.

\[P_n(x) = \frac{1}{2^n n!} \frac{d^n}{dx^n} \left(x^2-1 \right)^n.\]
Parameters:
  • n (int) – degree \(n \geq 0\).

  • x (float) – argument x.

Returns:

unassociated Legendre polynomials.

Return type:

double

nexus.AssociatedLegendre(l, m, x, csp)

Computes the associated Legendre polynomials of the degree l, order m, and argument x. The Condon–Shortley phase can be included.

\[P_\ell^m(x) = \left(1- x^2 \right)^{m/2} \frac{d^m}{dx^m} \left(P_\ell(x) \right).\]

If the Condon–Shortley phase is included a factor \((-1)^m\) is multiplied to \(P_\ell^m(x)\).

Parameters:
  • l (int) – degree \(\ell \geq 0\).

  • m (int) – order \(m \geq 0\).

  • x (float) – argument x.

  • csp (bool) – If True the Condon-Shortley phase is included.

Returns:

associated Legendre polynomials.

Return type:

double

nexus.SphericalHarmonics(l, m, theta, phi, csp, racah_norm)

Computes the spherical harmonic of the degree l, order m, and arguments theta and phi. The Condon–Shortley phase can be included. The normalization factor can be specified.

\[Y_\ell^m(\theta,\varphi) = N P_\ell^m(\cos{\theta}) e^{i m \varphi}.\]

where \(P_\ell^m(x)\) are the associated Legendre polynominals (without the Condon-Shortly phase) and the normalization factor is

\[N = \sqrt{\frac{2 \ell + 1}{4 \pi} \frac{(\ell-m)!}{(\ell+m)!}}\]

If the Condon–Shortley phase is included a factor \((-1)^m\) is multiplied to \(Y_\ell^m(x)\).

If the Racah norm is used, the normalization is \(N = \sqrt{\frac{(\ell-m)!}{(\ell+m)!}}\).

Parameters:
  • l (int) – degree \(\ell \geq 0\).

  • m (int) – order \(m \geq 0\).

  • theta (float) – argument \(\theta\).

  • phi (float) – argument \(\varphi\).

  • csp (bool) – If True the Condon-Shortley phase is included.

  • Racah_norm (bool) – If True the Racah norm is used.

Returns:

spherical harmonics.

Return type:

Complex

nexus.DiagonalizeHermitian(matrix)

Diagonalizes a hermitian matrix. Used to diagonalize the nuclear Hamiltonian. See Eq. (14) [Sturhahn].

Parameters:

matrix (ndarray) – Hermitian matrix to diagonalize.

Returns:

Eigenvalues and eigenvectors of the matrix.

Return type:

Eigensystem

Example

site = nx.BareHyperfine(
  magnetic_field = 33,
  magnetic_theta = nx.DegToRad(90)
)

# returns the ground state Hamiltonian in units of linewidth
hamiltonian = nx.HamiltonianGroundState(site, nx.moessbauer.Fe57)

# returns the real eigenvalues and eigenvectors
# the eigenvalues give the detuning energy in units of linewidth
eigen_system = nx.DiagonalizeHermitian(hamiltonian)

# list of eigenvalues
print(eigen_system.eigenvalues)

for vector in eigen_system.eigenvectors:
  print("Eigenvector")
  print(vector)
class nexus.Eulerangles(alpha=0.0, beta=0.0, gamma=0.0)

Constructor for Eulerangles class.

Parameters:
  • alpha (float) – \(\alpha\) in rad.

  • beta (float) – \(\beta\) in rad.

  • gamma (float) – \(\gamma\) in rad.

alpha

\(\alpha\) in rad.

Type:

float

beta

\(\beta\) in rad.

Type:

float

gamma

\(\gamma\) in rad.

Type:

float

nexus.VecSpherHarmPolarization(*args)

Calculates the vector spherical harmonics in the polarization base of \(\sigma\) and \(\pi\) for given Euler angles. See Eq. (15), [Sturhahn].

Parameters:
  • L (int) – photon angular momentum of the transition, dipole L = 1, quadrupole L = 2, …

  • lambda (int) – transition type, magnetic transition lambda = 0, electronic transition lambda = 1

  • euler (Eulerangles) – Euler angles. Default is 0 for all angles.

Returns:

first dimension has size of 2*L+1, each entry has a .simga and .pi attribute.

Return type:

ndarray

Example

print("E1, L=1, lam= 1")
VSH = nx.VecSpherHarmPolarization(1,1)
print(VSH[0].sigma)
print(VSH[0].pi)
print(VSH[1].sigma)
print(VSH[1].pi)
print(VSH[2].sigma)
print(VSH[2].pi)


print("M1, L=1, lam= 0")
VSH = nx.VecSpherHarmPolarization(1,0)
print(VSH[0].sigma)
print(VSH[0].pi)
print(VSH[1].sigma)
print(VSH[1].pi)
print(VSH[2].sigma)
print(VSH[2].pi)


print("E2, L=2, lam= 1")
VSH = nx.VecSpherHarmPolarization(2,1)
print(VSH[0].sigma)
print(VSH[0].pi)
print(VSH[1].sigma)
print(VSH[1].pi)
print(VSH[2].sigma)
print(VSH[2].pi)
print(VSH[3].sigma)
print(VSH[3].pi)
print(VSH[4].sigma)
print(VSH[4].pi)