Nuclear transitions

Functions to calculate the nuclear transition energies and transition operators.

nexus.HamiltonianGroundState(hyperfine, isotope)

Calculates the matrix Hamiltonian \(H\) of the ground state in the basis of the nuclear spin operator, after Eq. (12) [Sturhahn]. With a ground state spin of \(J\) one obtains a square matrix of dimension \(2J+1\). The matrix elements \(H_{m',m}\) are in the range from \(m = -J,...,J\). In general, the matrix must be diagonalized in order to find the energy detuning (eigenvalues) in units of the linewidth (\(\Gamma\)).

Parameters:
  • hyperfine (BareHyperfine) – Hyperfine parameters acting on the nucleus.

  • isotope (MoessbauerIsotope) – Moessbauer isotope on which the hyperfine parameters act.

Returns:

Complex matrix Hamiltonian. Energy in units of the linewidth (\(\Gamma\)).

Return type:

ndarray

Example

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

# returns the ground state Hamiltonian
hamiltonian = nx.HamiltonianGroundState(site, nx.moessbauer.Fe57)
print(hamiltonian)
nexus.HamiltonianExcitedState(hyperfine, isotope)

Calculates the matrix Hamiltonian \(H\) of the excited state in the basis of the nuclear spin operator, after Eq. (12) [Sturhahn]. With an excited state spin of \(J\) one obtains a square matrix of dimension \(2J+1\). The matrix elements \(H_{m',m}\) are in the range from \(m = -J,...,J\). In general, the matrix must be diagonalized in order to find the energy detuning (eigenvalues) in units of the linewidth (\(\Gamma\)).

Parameters:
  • hyperfine (BareHyperfine) – Hyperfine parameters acting on the nucleus.

  • isotope (MoessbauerIsotope) – Moessbauer isotope on which the hyperfine parameters act.

Returns:

Complex matrix Hamiltonian. Energy in units of the linewidth (\(\Gamma\)).

Return type:

ndarray

Example

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

# returns the excited state Hamiltonian
hamiltonian = nx.HamiltonianExcitedState(site, nx.moessbauer.Fe57)
print(hamiltonian)
nexus.multipol_TLni(eigenvec_ground, eigenvec_excited, isotope, multipolarity2)

Calculates the factor \(\delta_{L,\lambda} T_{Lni}^{(\lambda)}\) of Eq. (7) in [Sturhahn]. \(I_i\) and \(I_n\) are spins of the ground and excited state, respectively. \(L\) is the angular momentum from the multipolarity (E1 or M1 \(L=1\), E2 \(L=2\)). \(\lambda\) indicates electronic (1) or magnetic (0) scattering. The sum over \(M\) (which is \(\Delta m\)) is the multipole selection rule and runs over all possible changes in angular momentum from \(M = -L, ..., L\). The sum over \(m\) runs over all possible ground state spin quantum numbers \(-I_i, ..., I_i`\). The eigenvector components of the ground and excited state, given in bra-ket notation, are summed up for each ground state \(m\) and all possible multipole transition \(M\).

Parameters:
  • eigenvec_ground (ndarray) – Complex eigenvector of the ground state corresponding to a magnetic quantum number \(m_g\).

  • eigenvec_excited (ndarray) – Complex eigenvector of the excited state corresponding to a magnetic quantum number \(m_e\).

  • isotope (MoessbauerIsotope) – isotope for the calculation.

  • multipolarity2 (bool) – Specifies if the first or second multipole component should be calculated, e.g. for a M1E2 transition. If False the first multipole component of the transition is calculated, e.g. M1. If True the second one is calculated, e.g. E2.

Returns:

Returns a 3 dimensional vector. The first two entries are the vector components (sigma and pi) corresponding to the vector spherical harmonics \(Y_{LM}^{\lambda}\). The third component is the corresponding factor for pure isotropic scattering.

Return type:

ndarray

nexus.HyperfineTransitions(hyperfine, isotope)

Calculate all hyperfine-split transitions in a nucleus for a given set of hyperfine parameters. The model is described in [Sturhahn]. The eigenvalues and eigenvectors are calculated by diagonalization of the nuclear Hamiltonian, Eq. (12), and the energy detuning is found by Eq. (14). The calculated polarization dependent transition matrix of the returned Transition objects is the one of the single transition and is the nominator of the last term in Eq. (7)

\[\mathbf{T}_{\mu, \nu} = \left[\sum_{L,\lambda} \delta_{L,\lambda} T_{Lni}^{\lambda}\right]_{\mu}^* \left[\sum_{L,\lambda} \delta_{L',\lambda'}^* T_{L'ni}^{\lambda'}\right]_{\nu}\]

already weighted by the multipole mixing coefficients. The sum over all transitions \(n\) is not calculated here. The angles of the hyperfine parameters are defined from an arbitrary quantization axis of the nucleus. For calculations of photon interactions with the nuclear system, this quantization axis is set along the photon wave vector \(k\) in nexus.

Parameters:
  • hyper (BareHyperfine) – Hyperfine parameters acting on the nucleus.

  • isotope (MoessbauerIsotope) – Moessbauer isotope on which the hyperfine parameters act.

Returns:

List of all possible transitions in the nucleus.

Return type:

Transitions

Example

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

# returns a Transitions object
nuclear_transitions = nx.HyperfineTransitions(site, nx.moessbauer.Fe57)

# a single transition
print(nuclear_transitions[0])

# energy detuning in units of Gamma
print(nuclear_transitions[0].energy_detuning)

# number of all transitions in the nucleus
print(len(nuclear_transitions))

# all transitions
for trans in nuclear_transitions:
  print(trans)