atomic-physics API

core

class atomic_physics.core.Atom(magnetic_field: float, level_data: dict[atomic_physics.core.Level, atomic_physics.core.LevelData], level_states: dict[atomic_physics.core.Level, atomic_physics.core.LevelStates], transitions: dict[str, atomic_physics.core.Transition], nuclear_spin: float)

Represents a specific atom at a given magnetic field.

magnetic_field

the applied magnetic field (Tesla)

Type:

float

level_data

atomic structure dictionary, mapping Levels to LevelData.

Type:

dict[atomic_physics.core.Level, atomic_physics.core.LevelData]

level_states

level states dictionary, mapping Levels to LevelStates.

Type:

dict[atomic_physics.core.Level, atomic_physics.core.LevelStates]

transitions

transitions dictionary, mapping names of transitions to Transition objects.

Type:

dict[str, atomic_physics.core.Transition]

nuclear_spin

Nuclear spin.

Type:

float

num_states

the total number of states contained within the atom.

Type:

int

levels

tuple containing all of the Levels within this atom.

Type:

tuple[atomic_physics.core.Level, …]

state_energies

vector of energies of each state in the atom. States are ordered by decreasing energy, with state_energies[0] corresponding the state with highest energy and state_energies[-1] corresponding to the ground state.

Type:

numpy.ndarray

state_vectors

array of state eigenvectors in the high-field basis. The vector state_vectors[:, state_index] gives the representation of the state with energy state_energies[state_index] in the basis of high-field energy eigenstates (M_I, M_J). NB we use the high-field basis internally because it is the most convenient to represent our operators in. These states are not energy-ordered and are not generally energy eigenstates.

Type:

numpy.ndarray

high_field_M_I

array of values of M_I for each state in the high-field basis used in state_vectors.

Type:

numpy.ndarray

high_field_M_J

array of values of M_J for each state in the high-field basis used in state_vectors.

Type:

numpy.ndarray

M

vector of magnetic quantum numbers for each state.

Type:

numpy.ndarray

F

vector of total (electron + nucleus) angular momentum for each state. This is not generally a good quantum number and should be read as “this is the value of F we would get if we started in this state and adiabatically ramped the magnetic field to 0”.

Type:

numpy.ndarray

M_I

vector of nuclear magnetic quantum numbers for each state. This is not generally a good quantum number and should be read as “this is the value of M_I we would get if we started in this state and adiabatically ramped the magnetic field to infinity”.

Type:

numpy.ndarray

M_J

vector of electron magnetic quantum numbers for each state. This is not generally a good quantum number and should be read as “this is the value of M_J we would get if we started in this state and adiabatically ramped the magnetic field to infinity”.

Type:

numpy.ndarray

_electric_multipoles

electric multiple matrix (scattering amplitudes between states). See get_electric_multipoles(). To keep construction fast, we calculate this lazily when needed.

Type:

numpy.ndarray | None

_magnetic_dipoles

magnetic dipole matrix elements. To keep construction fast, we calculate this lazily when needed.

Type:

numpy.ndarray | None

get_electric_multipoles() ndarray

Returns the electric multi-pole matrix for each transition.

The values we return are scattering amplitudes between states (square root of the spontaneous scattering rates), rather than matrix elements themselves, since they are somewhat more convenient to work with.

Currently we only calculate dipole (E1) and Quadrupole (E2) scattering amplitudes.

To do: define what we mean by a matrix element, and how the amplitudes stored in ePole are related to the matrix elements and to the Rabi frequencies.

To do: double check the signs of the amplitudes.

get_level_for_state(state: int) Level

Returns the level a state lies in.

Parameters:

state – index of the state.

Returns:

the Level the state lies within.

get_magnetic_dipoles() ndarray

Returns an array of magnetic dipole matrix elements.

The matrix elements are given by \(R_{ul} := (-1)^{q+1} \left<u|\mu_q|l\right>\) where:

  • \(\left|u\right>\) (\(\left|l\right>\)) is the state with greater (lower) energy

  • \(q := Mu - Ml\).

  • \(\mu_q\) is the \(q\) th component of the magnetic dipole operator in the spherical basis.

get_rabi_rf(lower: int, upper: int, amplitude: float) float

Returns the Rabi frequency for a magnetic dipole transition.

See also get_magnetic_dipoles().

Parameters:
  • lower – index of the state with lower energy involved in the transition.

  • upper – index of the state with higher energy involved in the transition.

  • amplitude – amplitude of the component of the driving magnetic field ( spherical basis) which couples to this transition (Tesla).

Returns:

the Rabi frequency. We retain phase information, so this can be either positive or negative. We define the Rabi frequency so that \(t_{\pi} = \pi / |\Omega|\)

get_saturation_intensity(transition: str) float

Returns the saturation intensity for a transition.

We adopt a convention whereby, for a resonantly-driven cycling transition, one saturation intensity gives equal stimulated and spontaneous transition rates.

Parameters:

transition – the transition name

Returns:

saturation intensity (W/m^2)

get_slice_for_level(level: Level) slice

Returns a slice object, which can be used to index into a state vector to select only the states which lie within a given level.

Parameters:

level – the level to select.

Returns:

the slice object.

get_state_for_F(level: Level, F: float, M_F: float) int

Returns the index of the state with a given value of \(F\) and \(M_F\) within a given level.

\(F\) is only a good quantum number at zero field. At other fields, this is a “best guess” at the right state (see F).

Parameters:
  • level – the level to look within.

  • F – the value of \(F\) to look for.

  • M_F – the value of \(M_F\) to look for.

Returns:

the index of the corresponding state.

get_state_for_MI_MJ(level: Level, M_I: float, M_J: float) int

Returns the index of the state with a given value of \(M_I\) and \(M_J\) within a given level.

\(M_I\) and \(M_J\) are only good quantum numbers at high field. At other fields, this is a “best guess” at the right state (see M_I and M_J).

Parameters:
  • level – the level to look within.

  • M_I – the value of \(M_I\) to look for.

  • M_J – the value of \(M_J\) to look for.

Returns:

the index of the corresponding state.

get_states_for_M(level: Level, M: float) ndarray

Returns the indicates of the states within a given level with the specified value of \(M\).

Parameters:
  • level – the level to look within.

  • M – the value of \(M\) to look for.

Returns:

array of state indices.

get_states_for_level(level: Level) ndarray

Returns an array containing the indices of states which lie within a given level.

This can be used, for example, to index into a state vector to select only the states which lie within a given level.

Parameters:

level – the level to select.

Returns:

array of indices of states within the level.

get_transition_for_levels(levels: tuple[atomic_physics.core.Level, atomic_physics.core.Level]) Transition

Returns the transition between a pair of Levels.

Parameters:

levels – tuple containing the two Levels involved in the transition.

get_transition_frequency_for_states(states: tuple[int, int], relative: bool = True) float

Returns the frequency (angular units) of the transition between a pair of states.

Parameters:
  • states – tuple of indices of the states involved in the transition.

  • relative – if False the returned frequency is the absolute frequency of the transition between the two states. If True we subtract the centre of gravity frequency of the transition between the two levels the states lie in to calculate the frequency relative to the overall transition’s centre of gravity.

Returns:

the transition frequency (rad/s).

intensity_to_power(transition: str, waist_radius: float, intensity: float) float

Returns the power needed to achieve a given intensity.

Parameters:
  • transition – name of the transition to drive

  • waist_radius – Gaussian beam waist (\(1/e^2\) intensity radius in m)

  • intensity – beam intensity (saturation intensities).

Returns:

beam power (W)

class atomic_physics.core.AtomFactory(level_data: tuple[atomic_physics.core.LevelData, ...], transitions: dict[str, atomic_physics.core.Transition], nuclear_spin: float)

Flexible interface for creating Atoms.

Example usage:

from atomic_physics.core import Atom
from atomic_physics.ions.ca40 import Ca40


atom: Atom = Ca40(magnetic_field=10e-4)
level_data

tuple of atomic structure data for each level in the atom

Type:

tuple[atomic_physics.core.LevelData, …]

transitions

tuple of transitions between levels

Type:

dict[str, atomic_physics.core.Transition]

nuclear_spin

the atom’s nuclear spin

Type:

float

filter_levels(level_filter: tuple[atomic_physics.core.Level, ...]) AtomFactory

Returns a new AtomFactory which contains only a subset of the atomic levels.

This method is used to avoid unnecessary calculations when not all atomic levels are needed.

Parameters:

level_filter – list of Levels to include.

Returns:

a new AtomFactory containing only the selected subset of this AtomFactory's levels.

class atomic_physics.core.Laser(transition: str, polarization: int, intensity: float, detuning: float)

Represents a laser.

transition

string giving the name of the transition driven by the laser.

Type:

str

polarization

the laser’s polarization, defined as the difference in angular momentum between the higher and lower energy states coupled by this laser (q = M_upper - M_lower). q = +1 for σ+ light, q = -1 for σ- light and q = 0 for π light.

Type:

int

intensity

the laser intensity (saturation intensities).

Type:

float

detuning

the laser’s detuning (rad/s) from the transition centre of gravity, defined so that w_laser = w_transition + detuning.

Type:

float

class atomic_physics.core.Level(n: int, L: float, J: float, S: float)

Represents a single level.

n

the level’s principal quantum number

Type:

int

L

the level’s orbital angular momentum

Type:

float

J

the level’s total (spin + orbital) electronic angular momentum

Type:

float

S

the level’s spin angular momentum

Type:

float

class atomic_physics.core.LevelData(level: Level, Ahfs: float, Bhfs: float, g_J: float | None = None, g_I: float = 0.0)

Atomic structure data about a single level.

level

the Level this data is for.

Type:

atomic_physics.core.Level

Ahfs

Nuclear A coefficient.

Type:

float

Bhfs

Nuclear B (quadrupole) coefficient.

Type:

float

g_J

G factor. If None, we use the Lande g factor.

Type:

float | None

g_I

Nuclear g factor.

Type:

float

class atomic_physics.core.LevelStates(frequency: float, start_index: int, stop_index: int, num_states: int)

Stores information about the states within a level.

freq

frequency (rad/s) of the centre-of-gravity transition from the ground-level to this level. These frequencies are calculated automatically by combining data from atom’s transitions. In cases where the the level structure is multiply connected, taking different routes between the same levels may give slightly different values for the level frequencies due to inconsistencies between the measured transition frequencies. We do not make guarantees about which route is taken for these calculations and this data should not be relied on for accurate calculations (see Atom's get_transition_frequency_for_states method for that).

start_index

index into the state vector for the first (highest energy) state within this level.

Type:

int

stop_index

index into the state vector of the laset (lowest energy) state within this level.

Type:

int

num_states

the number of states within the level.

Type:

int

class atomic_physics.core.RFDrive(frequency: float, amplitude: float, polarization: ndarray)

Represents an AC magnetic field, which drives magnetic dipole transitions.

frequency

frequency of the RF drive (rad/s).

Type:

float

amplitude

amplitude of the field (T).

Type:

float

polarization

Jones vector describing the magnetic field’s polarization.

Type:

numpy.ndarray

class atomic_physics.core.Transition(lower: Level, upper: Level, frequency: float, einstein_A: float)

Represents a transition between a pair of Levels.

lower

the Level in the transition with lower energy

Type:

atomic_physics.core.Level

upper

the Level in the transition with greater energy

Type:

atomic_physics.core.Level

freq

the transition frequency (rad/s)

einstein_A

the transition’s Einstein A coefficient

Type:

float

Rate Equations

class atomic_physics.rate_equations.Rates(atom: Atom)

Rate equations calculations.

See the Rate Equations section of the documentation for details.

Example usage:

# Electron shelving simulation in 43Ca+
import matplotlib.pyplot as plt
import numpy as np
from scipy.linalg import expm

from atomic_physics.core import Laser
from atomic_physics.ions import ca43
from atomic_physics.rate_equations import Rates

t_ax = np.linspace(0, 100e-6, 100)  # Scan the duration of the "shelving" pulse
intensity = 0.02  # 393 intensity

ion = ca43.Ca43(magnetic_field=146e-4)

# Ion starts in the F=4, M=+4 "stretched" state within the 4S1/2 ground-level
stretch = ion.get_state_for_F(ca43.ground_level, F=4, M_F=+4)
Vi = np.zeros((ion.num_states, 1))
Vi[stretch] = 1

# Tune the 393nm laser to resonance with the
# 4S1/2(F=4, M_F=+4) <> 4P3/2(F=5, M_F=+5) transition
detuning = ion.get_transition_frequency_for_states(
    (stretch, ion.get_state_for_F(ca43.P32, F=5, M_F=+5))
)
lasers = (
    Laser("393", polarization=+1, intensity=intensity, detuning=detuning),
)  # resonant 393 sigma+

rates = Rates(ion)
transitions = rates.get_transitions_matrix(lasers)

shelved = np.zeros(len(t_ax))  # Population in the 3D5/2 level at the end
for idx, t in np.ndenumerate(t_ax):
    Vf = expm(transitions * t) @ Vi  # NB use of matrix operations here!
    shelved[idx] = np.sum(Vf[ion.get_slice_for_level(ca43.shelf)])

plt.plot(t_ax * 1e6, shelved)
plt.ylabel("Shelved Population")
plt.xlabel("Shelving time (us)")
plt.grid()
plt.show()
get_spont_matrix() ndarray

Returns the spontaneous emission matrix.

get_steady_state_populations(transitions: ndarray) ndarray

Returns the steady-state population vector for a given transitions matrix.

Parameters:

transitions – transitions matrix.

Returns:

steady-state population vector.

get_stim_matrix(lasers: tuple[atomic_physics.core.Laser, ...]) ndarray

Returns the stimulated emission matrix for a set of lasers.

get_transitions_matrix(lasers: tuple[atomic_physics.core.Laser, ...]) ndarray

Returns the complete (spontaneous + stimulated emissions) transitions matrix for a given set of lasers.

Helper Utilities

atomic_physics.utils.ac_zeeman_shift_for_state(atom: Atom, state: int, drive: RFDrive) float

Returns the AC Zeeman shift on a given state resulting from an applied RF field.

The calculated shift includes the counter-rotating term (Bloch-Siegert shift) and is calculated by summing \(\frac{1}{2} \Omega^2 \left(\omega_0 / (w_0^2 - w_{\mathrm{rf}}^2)\right)\) over all magnetic dipole transitions which couple to the state. Where \(\Omega\) is the Rabi frequency for the transition and \(w_0\) is the transition frequency.

Example:

import numpy as np

from atomic_physics.core import RFDrive
from atomic_physics.ions.ca40 import Ca40, ground_level
from atomic_physics.polarization import SIGMA_PLUS_POLARIZATION
from atomic_physics.utils import ac_zeeman_shift_for_state

ion = Ca40(magnetic_field=10e-4)
w_transition = ion.get_transition_frequency_for_states(
    states=(
        ion.get_states_for_M(level=ground_level, M=+1 / 2),
        ion.get_states_for_M(level=ground_level, M=-1 / 2),
    )
)
w_rf = w_transition + 1e6 * 2 * np.pi  # RF is 1 MHz blue of the transition

ac_zeeman_shift_for_state(
    atom=ion,
    state=ion.get_states_for_M(level=ground_level, M=+1 / 2),
    drive=RFDrive(frequency=w_rf, amplitude=1e-6, polarization=SIGMA_PLUS_POLARIZATION),
)
Parameters:
  • atom – the atom.

  • state – index of the state to calculate the shift for.

  • drive – the applied RF drive.

Returns:

the AC Zeeman shift (rad/s).

atomic_physics.utils.ac_zeeman_shift_for_transition(atom: Atom, states: tuple[int, int], drive: RFDrive) float

Returns the AC Zeeman shift on a transition resulting from an applied RF field.

See ac_zeeman_shift_for_state() for details.

Parameters:
  • atom – the atom.

  • states – tuple containing indices of the two states involved in the transition.

  • drive – the applied RF drive.

Returns:

the AC Zeeman shift (rad/s).

atomic_physics.utils.d2f_dB2(atom_factory: AtomFactory, magnetic_field: float, states: tuple[int, int], eps: float = 1e-06) float

Returns the second-order field-sensitivity (\(\frac{\mathrm{d}^2f}{\mathrm{d}B^2}\)) of a transition between two states in the same level at a given magnetic field.

Parameters:
  • atom_factory – factory class for the atom of interest.

  • magnetic_field – the magnetic to calculate the field sensitivity at.

  • states – tuple of indices involved in this transition.

  • eps – field difference (T) to use when calculating derivatives numerically.

Returns:

the transition’s field sensitivity (rad/s/T^2).

atomic_physics.utils.df_dB(atom_factory: AtomFactory, magnetic_field: float, states: tuple[int, int], eps: float = 1e-06) float

Returns the field-sensitivity (\(\frac{\mathrm{d}f}{\mathrm{d}B}\)) of a transition between two states in the same level at a given magnetic field.

Parameters:
  • atom_factory – factory class for the atom of interest.

  • magnetic_field – the magnetic to calculate the field sensitivity at.

  • states – tuple of indices involved in this transition.

  • eps – field difference (T) to use when calculating derivatives numerically.

Returns:

the transition’s field sensitivity (rad/s/T).

atomic_physics.utils.field_insensitive_point(atom_factory: AtomFactory, level_0: Level, F_0: float, M_F_0: float, level_1: Level, F_1: float, M_F_1: float, magnetic_field_guess: float = 0.0001, eps: float = 0.0001) float | None

Returns the magnetic field at which the frequency of a transition between two states in the same level becomes first-order field independent.

Since the energy ordering of states can change with magnetic field, we label states by F and M_F instead of using state indices.

Parameters:
  • atom_factory – factory class for the atom of interest.

  • level_0 – level the first state involved in the transition lies in.

  • F_0 – value of F for the first state involved in the transition.

  • M_F_0 – value of M_F for the first state involved in the transition.

  • level_1 – level the second state involved in the transition lies in.

  • F_1 – value of F for the second state involved in the transition.

  • M_F_1 – value of M_F for the second state involved in the transition.

  • magnetic_field_guess – Initial guess for the magnetic field insensitive point (T). This is used both as a seed for the root finding algorithm and as a scale factor to help numerical accuracy.

  • eps – step size as a fraction of magnetic_field_guess to use when calculating numerical derivatives.

Returns:

the field-independent point (T) or None if none found.

atomic_physics.utils.rayleigh_range(transition: Transition, waist_radius: float) float

Returns the Rayleigh range for a given beam radius.

Parameters:
  • transition – the transition the laser drives.

  • waist_radius – Gaussian beam waist (\(1/e^2\) intensity radius in m).

Returns:

the Rayleigh range (m).

Operators

atomic_physics.operators.AngularMomentumLoweringOp(j: float) ndarray

Angular momentum lowering operator (\(J_-\)) represented in the basis of angular momentum eigenstates.

Basis states are labelled in order of increasing \(M_J\).

The returned operator is defined so that:

\[\begin{split}J_{-\left(n, m\right)} &:= \left<M_J=n|J_-|M_J=m\right>\\ & = \delta\left(n,m-1\right)\sqrt{(J+m)(J-m+1)}\end{split}\]

This operator is related to the spherical basis operator \(J_{-1}\) by:

\[J_{-1} = +\frac{1}{\sqrt{2}} J_-\]
atomic_physics.operators.AngularMomentumProjectionOp(j: float) ndarray

Angular momentum projection operation represented in the basis of angular momentum eigenstates.

Basis states are labelled in order of increasing \(M_J\).

The returned operator is defined so that:

\[\begin{split}J_{z\left(n, m\right)} &:= \left<M_J=n|J_z|M_J=m\right>\\ & = m\delta\left(n,m\right)\end{split}\]
atomic_physics.operators.AngularMomentumRaisingOp(j: float) ndarray

Angular momentum raising operator (\(J_+\)) represented in the basis of angular momentum eigenstates.

Basis states are labelled in order of increasing \(M_J\).

The returned operator is defined so that:

\[\begin{split}J_{+\left(n, m\right)} &:= \left<M_J=n|J_+|M_J=m\right>\\ & = \delta\left(n,m+1\right)\sqrt{(J-m)(J+m+1)}\end{split}\]

This operator is related to the spherical basis operator \(J_{+1}\) by:

\[J_{+1} = -\frac{1}{\sqrt{2}} J_+\]
atomic_physics.operators.expectation_value(state_vectors: ndarray, operator: ndarray) ndarray

Calculates the expectation value of an operator for a given set of states in a given basis.

Parameters:
  • state_vectors – array of shape (num_basis_states, num_states) where num_basis_states is the number of states in the basis that operator is represented in and num_states is the number of states to calculate the expectation value for. The vector state_vector[:, state_index] should give the representation of state state_index in the basis.

  • operator – operator represented in the basis given by state_vectors.

Returns:

a vector, giving the expectation value of operator for each state in state_vectors.

Polarization

Tools for working with polarization vectors.

See Definitions and Conventions for definitions and further discussion about representation of polarization within atomic-physics.

atomic_physics.polarization.PI_POLARIZATION = array([0, 0, 1])

Jones vector for π polarization.

π-polarized radiation drives transitions where the magnetic quantum number is the same in both states (\(M_{\mathrm{upper}} = M_{\mathrm{lower}}\)).

atomic_physics.polarization.SIGMA_MINUS_POLARIZATION = array([-0.70710678+0.j        ,  0.        +0.70710678j,         0.        +0.j        ])

Jones vector for σ- polarization.

σ- polarized radiation drives transitions where the magnetic quantum number in the state with higher energy is 1 lower than in the state with lower energy ( \(M_{\mathrm{upper}} = M_{\mathrm{lower}} - 1\)).

atomic_physics.polarization.SIGMA_PLUS_POLARIZATION = array([0.70710678+0.j        , 0.        +0.70710678j,        0.        +0.j        ])

Jones vector for σ+ polarization.

σ+ polarized radiation drives transitions where the magnetic quantum number in the state with higher energy is 1 greater than in the state with lower energy ( \(M_{\mathrm{upper}} = M_{\mathrm{lower}} + 1\)).

atomic_physics.polarization.X_POLARIZATION: ndarray = array([1, 0, 0])

Jones vector for a field with linear polarization along the x-axis.

atomic_physics.polarization.Y_POLARIZATION: ndarray = array([0, 1, 0])

Jones vector for a field with linear polarization along the y-axis.

atomic_physics.polarization.Z_POLARIZATION: ndarray = array([0, 0, 1])

Jones vector for a field with linear polarization along the z-axis.

atomic_physics.polarization.cartesian_to_spherical(cartesian: ndarray) ndarray

Converts a vector in Cartesian coordinates to the spherical basis.

Parameters:

cartesian – input array in Cartesian coordinates.

Returns:

the input array converted to the spherical basis.

atomic_physics.polarization.dM_for_transition(atom: Atom, states: tuple[int, int]) int

Returns the change in magnetic quantum number for a given transition.

\(\delta M := M_{\mathrm{upper}} - M_{\mathrm{lower}}\) where \(M_{\mathrm{upper}}\) (\(M_{\mathrm{lower}}\)) is the magnetic quantum number for the state with greater (lower) energy.

Parameters:
  • atom – the atom.

  • states – tuple containing indices of the two states involved in the transition.

Returns:

the change in magnetic quantum number.

atomic_physics.polarization.half_wave_plate() ndarray

Returns the Jones matrix for a field propagating along the z-axis passing through a half-wave plate, whose fast axis is aligned along \(\hat{\mathbf{x}}\).

atomic_physics.polarization.quarter_wave_plate() ndarray

Returns the Jones matrix for a field propagating along the z-axis passing through a quarter-wave plate, whose fast axis is aligned along \(\hat{\mathbf{x}}\).

atomic_physics.polarization.retarder_jones(phi: float) ndarray

Returns the Jones matrix for a field propagating along the z-axis passing through a retarder, whose fast axis is aligned along \(\hat{\mathbf{x}}\).

Parameters:

phi – the phase shift added by the retarder.

Returns:

the Jones matrix.

atomic_physics.polarization.rotate_jones_matrix(input_matrix: ndarray, theta: float) ndarray

Rotates a Jones matrix about the z-axis (x-y plane rotation).

the Jones matrix for a component rotated by an angle theta about the z-axis is given by R(theta) @ J @ R(-theta).

Parameters:
  • input_array – Jones matrix to be rotated.

  • theta – rotation angle (radians).

Returns:

the rotated Jones matrix.

atomic_physics.polarization.spherical_dot(vec_a: ndarray, vec_b: ndarray) ndarray

Returns the dot product between two vectors in the spherical basis.

atomic_physics.polarization.spherical_to_cartesian(spherical: ndarray) ndarray

Converts a vector in the spherical basis to Cartesian coordinates.

Parameters:

cartesian – input array in the spherical basis.

Returns:

the input array converted to Cartesian coordinates.

Atoms

Two State Atom

Ideal spin 1/2 atom with two states.

atomic_physics.atoms.two_state.TwoStateAtom = AtomFactory(level_data=(LevelData(level=Level(n=0, L=0, J=0.5, S=0.5), Ahfs=0, Bhfs=0, g_J=2.00231930436092, g_I=0.0),), transitions={}, nuclear_spin=0.0, level_states={Level(n=0, L=0, J=0.5, S=0.5): LevelStates(frequency=0.0, start_index=0, stop_index=2, num_states=2)})

Ideal spin 1/2 atom with two states.

atomic_physics.atoms.two_state.field_for_frequency(frequency: float) float

Returns the B-field needed to produce a given transition frequency.

Parameters:

frequency – the desired transition frequency (rad/s).

Returns:

the required magnetic field (T).

atomic_physics.atoms.two_state.ground_level = Level(n=0, L=0, J=0.5, S=0.5)

The only level within the TwoStateAtom.

atomic_physics.atoms.two_state.lower_state = 1

Index of the M=-1/2 (lower-energy) state.

atomic_physics.atoms.two_state.upper_state = 0

Index of the M=+1/2 (higher-energy) state.

Ions

Barium 133

\(^{133}\mathrm{Ba}^+\)

The transition frequencies are calculated based on the \(^{138}\mathrm{Ba}^+\) frequencies from [1] and the isotope shifts in the second reference listed next to the frequency. Where no references are given, the transition frequencies were calculated based on transition frequencies between other levels. For example the frequency of the 1762 nm transition f(1762) was found from f(1762)=f(455)-f(614).

References:

* [1] A. Kramida, NIST Atomic Spectra Database (ver. 5.9) (2021)
* [2] Zhiqiang Zhang, K. J. Arnold, S. R. Chanu, R. Kaewuam,
  M. S. Safronova, and M. D. Barrett Phys. Rev. A 101, 062515 (2020)
* [3] N. Yu, W. Nagourney, and H. Dehmelt, Phys. Rev. Lett. 78, 4898 (1997)
* [4] N. J. Stone, Table of nuclear magnetic dipole and electric
  quadrupole moments, Atomic Data and Nuclear Data Tables, Volume 90, Issue 1 (2005)
* [5] H. Knab, K. H. Knöll, F. Scheerer and G. Werth, Zeitschrift für
  Physik D Atoms, Molecules and Clusters volume 25, pages205–208 (1993)
* [6] David Hucul, Justin E. Christensen, Eric R. Hudson, and
  Wesley C. Campbell, Phys. Rev. Lett. 119, 100501 (2017)
* [7] - J.E. Christensen, D. Hucul, W.C. Campbell et al., npj Quantum Inf 6, 35
  (2020).
atomic_physics.ions.ba133.Ba133 = AtomFactory(level_data=(LevelData(level=Level(n=6, L=0, J=0.5, S=0.5), Ahfs=-6.5766751523280195e-24, Bhfs=0, g_J=2.0024906, g_I=-1.54334), LevelData(level=Level(n=6, L=1, J=0.5, S=0.5), Ahfs=-1.2191969076e-24, Bhfs=0, g_J=0.6658935652130266, g_I=-1.54334), LevelData(level=Level(n=6, L=1, J=1.5, S=0.5), Ahfs=-2.064020851725e-25, Bhfs=0, g_J=1.3341064347869733, g_I=-1.54334), LevelData(level=Level(n=5, L=2, J=1.5, S=0.5), Ahfs=-3.104313865275e-25, Bhfs=0, g_J=0.799536139127816, g_I=-1.54334), LevelData(level=Level(n=5, L=2, J=2.5, S=0.5), Ahfs=1.8332127415e-26, Bhfs=0, g_J=1.200463860872184, g_I=-1.54334)), transitions={'493': Transition(lower=Level(n=6, L=0, J=0.5, S=0.5), upper=Level(n=6, L=1, J=0.5, S=0.5), frequency=3816572113379737.5, einstein_A=95300000.0), '455': Transition(lower=Level(n=6, L=0, J=0.5, S=0.5), upper=Level(n=6, L=1, J=1.5, S=0.5), frequency=4135068020081962.5, einstein_A=111000000.0), '650': Transition(lower=Level(n=5, L=2, J=1.5, S=0.5), upper=Level(n=6, L=1, J=0.5, S=0.5), frequency=2898508217515499.0, einstein_A=31000000.0), '585': Transition(lower=Level(n=5, L=2, J=1.5, S=0.5), upper=Level(n=6, L=1, J=1.5, S=0.5), frequency=3217004124217724.0, einstein_A=6000000.0), '614': Transition(lower=Level(n=5, L=2, J=2.5, S=0.5), upper=Level(n=6, L=1, J=1.5, S=0.5), frequency=3066132110108545.5, einstein_A=41200000.0), '1762': Transition(lower=Level(n=6, L=0, J=0.5, S=0.5), upper=Level(n=5, L=2, J=2.5, S=0.5), frequency=1068935909973417.1, einstein_A=0.033178500331785), '2051': Transition(lower=Level(n=6, L=0, J=0.5, S=0.5), upper=Level(n=5, L=2, J=1.5, S=0.5), frequency=918063895864238.6, einstein_A=0.0125)}, nuclear_spin=0.5, level_states={Level(n=6, L=1, J=1.5, S=0.5): LevelStates(frequency=4135068020081962.5, start_index=0, stop_index=8, num_states=8), Level(n=6, L=1, J=0.5, S=0.5): LevelStates(frequency=3816572113379737.5, start_index=8, stop_index=12, num_states=4), Level(n=5, L=2, J=2.5, S=0.5): LevelStates(frequency=1068935909973417.1, start_index=12, stop_index=24, num_states=12), Level(n=5, L=2, J=1.5, S=0.5): LevelStates(frequency=918063895864238.6, start_index=24, stop_index=32, num_states=8), Level(n=6, L=0, J=0.5, S=0.5): LevelStates(frequency=0, start_index=32, stop_index=36, num_states=4)})

\(^{133}\mathrm{Ba}^+\) atomic structure.

atomic_physics.ions.ba133.D32 = Level(n=5, L=2, J=1.5, S=0.5)

The \(\left|n=5, S=1/2, L=2, J=3/2\right>\) level.

atomic_physics.ions.ba133.D52 = Level(n=5, L=2, J=2.5, S=0.5)

The \(\left|n=5, S=1/2, L=2, J=5/2\right>\) level.

atomic_physics.ions.ba133.P12 = Level(n=6, L=1, J=0.5, S=0.5)

The \(\left|n=6, S=1/2, L=1, J=1/2\right>\) level.

atomic_physics.ions.ba133.P32 = Level(n=6, L=1, J=1.5, S=0.5)

The \(\left|n=6, S=1/2, L=1, J=3/2\right>\) level.

atomic_physics.ions.ba133.S12 = Level(n=6, L=0, J=0.5, S=0.5)

The \(\left|n=6, S=1/2, L=0, J=1/2\right>\) level.

atomic_physics.ions.ba133.ground_level = Level(n=6, L=0, J=0.5, S=0.5)

Alias for the \(\left|n=6, S=1/2, L=0, J=1/2\right>\) ground level of \(^{133}\mathrm{Ba}^+\).

atomic_physics.ions.ba133.shelf = Level(n=5, L=2, J=2.5, S=0.5)

Alias for the \(\left|n=5, S=1/2, L=2, J=5/2\right>\) “shelf” level of \(^{133}\mathrm{Ba}^+\).

Barium 135

\(^{135}\mathrm{Ba}^+\)

The transition frequencies are calculated based on the \(^{138}\mathrm{Ba}^+\) frequencies from [1] and the isotope shifts in the second reference listed next to the frequency. Where no references are given, the transition frequencies were calculated based on transition frequencies between other levels. For example the frequency of the 1762 nm transition f(1762) was found from f(1762)=f(455)-f(614).

References:

* [1] A. Kramida, NIST Atomic Spectra Database (ver. 5.9) (2021)
* [2] Zhiqiang Zhang, K. J. Arnold, S. R. Chanu, R. Kaewuam,
  M. S. Safronova, and M. D. Barrett Phys. Rev. A 101, 062515 (2020)
* [3] N. Yu, W. Nagourney, and H. Dehmelt, Phys. Rev. Lett. 78, 4898 (1997)
* [4] N. J. Stone, Table of nuclear magnetic dipole and electric
  quadrupole moments, Atomic Data and Nuclear Data Tables, Volume 90, Issue 1 (2005)
* [5] H. Knab, K. H. Knöll, F. Scheerer and G. Werth, Zeitschrift für
  Physik D Atoms, Molecules and Clusters volume 25, pages205–208 (1993)
* [6] W. Becker, G. Werth, Zeitschrift für Physik A Atoms and Nuclei,
  Volume 311, Issue 1-2, pp. 41-47 (1983)
* [7] P Villemoes et al, J. Phys. B: At. Mol. Opt. Phys. 26 4289 (1993)
* [8] Roger E. Silverans, Gustaaf Borghs, Peter De Bisschop, and
  Marleen Van Hove, Phys. Rev. A 33, 2117 (1986)
* [9] K. Wendt, S. A. Ahmad, F. Buchinger, A. C. Mueller, R. Neugart, and
  E. -W. Otten, Zeitschrift für Physik A Atoms and Nuclei volume 318,
  pages 125–129 (1984)
atomic_physics.ions.ba135.Ba135 = AtomFactory(level_data=(LevelData(level=Level(n=6, L=0, J=0.5, S=0.5), Ahfs=2.3798658153882437e-24, Bhfs=0, g_J=2.0024906, g_I=0.5586266666666666), LevelData(level=Level(n=6, L=1, J=0.5, S=0.5), Ahfs=4.40368622169e-25, Bhfs=0, g_J=0.6658935652130266, g_I=0.5586266666666666), LevelData(level=Level(n=6, L=1, J=1.5, S=0.5), Ahfs=7.4874592695e-26, Bhfs=3.9093813885e-26, g_J=1.3341064347869733, g_I=0.5586266666666666), LevelData(level=Level(n=5, L=2, J=1.5, S=0.5), Ahfs=1.12371391152447e-25, Bhfs=1.9184328383892e-32, g_J=0.799536139127816, g_I=0.5586266666666666), LevelData(level=Level(n=5, L=2, J=2.5, S=0.5), Ahfs=-7.113086306024999e-27, Bhfs=2.563759062438e-26, g_J=1.200463860872184, g_I=0.5586266666666666)), transitions={'493': Transition(lower=Level(n=6, L=0, J=0.5, S=0.5), upper=Level(n=6, L=1, J=0.5, S=0.5), frequency=3816572113379584.5, einstein_A=95300000.0), '455': Transition(lower=Level(n=6, L=0, J=0.5, S=0.5), upper=Level(n=6, L=1, J=1.5, S=0.5), frequency=4135068020081993.5, einstein_A=111000000.0), '650': Transition(lower=Level(n=5, L=2, J=1.5, S=0.5), upper=Level(n=6, L=1, J=0.5, S=0.5), frequency=2898508217514774.5, einstein_A=31000000.0), '585': Transition(lower=Level(n=5, L=2, J=1.5, S=0.5), upper=Level(n=6, L=1, J=1.5, S=0.5), frequency=3217004124217184.0, einstein_A=6000000.0), '614': Transition(lower=Level(n=5, L=2, J=2.5, S=0.5), upper=Level(n=6, L=1, J=1.5, S=0.5), frequency=3066132110107825.5, einstein_A=41200000.0), '1762': Transition(lower=Level(n=6, L=0, J=0.5, S=0.5), upper=Level(n=5, L=2, J=2.5, S=0.5), frequency=1068935909974168.0, einstein_A=0.033178500331785), '2051': Transition(lower=Level(n=6, L=0, J=0.5, S=0.5), upper=Level(n=5, L=2, J=1.5, S=0.5), frequency=918063895864810.0, einstein_A=0.0125)}, nuclear_spin=1.5, level_states={Level(n=6, L=1, J=1.5, S=0.5): LevelStates(frequency=4135068020081994.0, start_index=0, stop_index=16, num_states=16), Level(n=6, L=1, J=0.5, S=0.5): LevelStates(frequency=3816572113379584.5, start_index=16, stop_index=24, num_states=8), Level(n=5, L=2, J=2.5, S=0.5): LevelStates(frequency=1068935909974168.0, start_index=24, stop_index=48, num_states=24), Level(n=5, L=2, J=1.5, S=0.5): LevelStates(frequency=918063895864810.0, start_index=48, stop_index=64, num_states=16), Level(n=6, L=0, J=0.5, S=0.5): LevelStates(frequency=0, start_index=64, stop_index=72, num_states=8)})

\(^{135}\mathrm{Ba}^+\) atomic structure.

atomic_physics.ions.ba135.D32 = Level(n=5, L=2, J=1.5, S=0.5)

The \(\left|n=5, S=1/2, L=2, J=3/2\right>\) level.

atomic_physics.ions.ba135.D52 = Level(n=5, L=2, J=2.5, S=0.5)

The \(\left|n=5, S=1/2, L=2, J=5/2\right>\) level.

atomic_physics.ions.ba135.P12 = Level(n=6, L=1, J=0.5, S=0.5)

The \(\left|n=6, S=1/2, L=1, J=1/2\right>\) level.

atomic_physics.ions.ba135.P32 = Level(n=6, L=1, J=1.5, S=0.5)

The \(\left|n=6, S=1/2, L=1, J=3/2\right>\) level.

atomic_physics.ions.ba135.S12 = Level(n=6, L=0, J=0.5, S=0.5)

The \(\left|n=6, S=1/2, L=0, J=1/2\right>\) level.

atomic_physics.ions.ba135.ground_level = Level(n=6, L=0, J=0.5, S=0.5)

Alias for the \(\left|n=6, S=1/2, L=0, J=1/2\right>\) ground level of \(^{135}\mathrm{Ba}^+\).

atomic_physics.ions.ba135.shelf = Level(n=5, L=2, J=2.5, S=0.5)

Alias for the \(\left|n=5, S=1/2, L=2, J=5/2\right>\) “shelf” level of \(^{135}\mathrm{Ba}^+\).

Barium 137

\(^{137}\mathrm{Ba}^+\)

The transition frequencies are calculated based on the \(^{138}\mathrm{Ba}^+\) frequencies from [1] and the isotope shifts in the second reference listed next to the frequency. Where no references are given, the transition frequencies were calculated based on transition frequencies between other levels. For example the frequency of the 1762 nm transition f(1762) was found from f(1762)=f(455)-f(614).

References:

* [1] A. Kramida, NIST Atomic Spectra Database (ver. 5.9) (2021)
* [2] Zhiqiang Zhang, K. J. Arnold, S. R. Chanu, R. Kaewuam,
  M. S. Safronova, and M. D. Barrett Phys. Rev. A 101, 062515 (2020)
* [3] N. Yu, W. Nagourney, and H. Dehmelt, Phys. Rev. Lett. 78, 4898 (1997)
* [4] N. J. Stone, Table of nuclear magnetic dipole and electric
  quadrupole moments, Atomic Data and Nuclear Data Tables, Volume 90, Issue 1 (2005)
* [5] H. Knab, K. H. Knöll, F. Scheerer and G. Werth, Zeitschrift für
  Physik D Atoms, Molecules and Clusters volume 25, pages 205–208 (1993)
* [6] R. Blatt and G. Werth, Phys. Rev. A 25, 1476 (1982)
* [7] P Villemoes et al, J. Phys. B: At. Mol. Opt. Phys. 26 4289 (1993)
* [8] Nicholas C. Lewty, Boon Leng Chuah, Radu Cazan, B. K. Sahoo, and
  M. D. Barrett, Opt. Express 21, 7131-7132 (2013)
* [9] Nicholas C. Lewty, Boon Leng Chuah, Radu Cazan, Murray D. Barrett,
  and B. K. Sahoo, Phys. Rev. A 88, 012518 (2013)
* [10] K. Wendt, S. A. Ahmad, F. Buchinger, A. C. Mueller, R. Neugart, and
  E. -W. Otten, Zeitschrift für Physik A Atoms and Nuclei volume 318,
  pages 125–129 (1984)
atomic_physics.ions.ba137.Ba137 = AtomFactory(level_data=(LevelData(level=Level(n=6, L=0, J=0.5, S=0.5), Ahfs=2.6629320068812833e-24, Bhfs=0, g_J=2.0024906, g_I=0.6249133333333333), LevelData(level=Level(n=6, L=1, J=0.5, S=0.5), Ahfs=4.927808370555e-25, Bhfs=0, g_J=0.6658935652130266, g_I=0.6249133333333333), LevelData(level=Level(n=6, L=1, J=1.5, S=0.5), Ahfs=8.4283612308e-26, Bhfs=6.12911488875e-26, g_J=1.3341064347869733, g_I=0.6249133333333333), LevelData(level=Level(n=5, L=2, J=1.5, S=0.5), Ahfs=1.2571715848627352e-25, Bhfs=2.951027153553318e-26, g_J=0.799536139127816, g_I=0.6249133333333333), LevelData(level=Level(n=5, L=2, J=2.5, S=0.5), Ahfs=-7.97065483347651e-27, Bhfs=3.9442027123522797e-26, g_J=1.200463860872184, g_I=0.6249133333333333)), transitions={'493': Transition(lower=Level(n=6, L=0, J=0.5, S=0.5), upper=Level(n=6, L=1, J=0.5, S=0.5), frequency=3816572113379097.5, einstein_A=95300000.0), '455': Transition(lower=Level(n=6, L=0, J=0.5, S=0.5), upper=Level(n=6, L=1, J=1.5, S=0.5), frequency=4135068020081369.0, einstein_A=111000000.0), '650': Transition(lower=Level(n=5, L=2, J=1.5, S=0.5), upper=Level(n=6, L=1, J=0.5, S=0.5), frequency=2898508217514336.5, einstein_A=31000000.0), '585': Transition(lower=Level(n=5, L=2, J=1.5, S=0.5), upper=Level(n=6, L=1, J=1.5, S=0.5), frequency=3217004124216608.0, einstein_A=6000000.0), '614': Transition(lower=Level(n=5, L=2, J=2.5, S=0.5), upper=Level(n=6, L=1, J=1.5, S=0.5), frequency=3066132110107203.0, einstein_A=41200000.0), '1762': Transition(lower=Level(n=6, L=0, J=0.5, S=0.5), upper=Level(n=5, L=2, J=2.5, S=0.5), frequency=1068935909974166.0, einstein_A=0.033178500331785), '2051': Transition(lower=Level(n=6, L=0, J=0.5, S=0.5), upper=Level(n=5, L=2, J=1.5, S=0.5), frequency=918063895864760.9, einstein_A=0.0125)}, nuclear_spin=1.5, level_states={Level(n=6, L=1, J=1.5, S=0.5): LevelStates(frequency=4135068020081369.0, start_index=0, stop_index=16, num_states=16), Level(n=6, L=1, J=0.5, S=0.5): LevelStates(frequency=3816572113379097.5, start_index=16, stop_index=24, num_states=8), Level(n=5, L=2, J=2.5, S=0.5): LevelStates(frequency=1068935909974166.0, start_index=24, stop_index=48, num_states=24), Level(n=5, L=2, J=1.5, S=0.5): LevelStates(frequency=918063895864760.9, start_index=48, stop_index=64, num_states=16), Level(n=6, L=0, J=0.5, S=0.5): LevelStates(frequency=0, start_index=64, stop_index=72, num_states=8)})

\(^{137}\mathrm{Ba}^+\) atomic structure.

atomic_physics.ions.ba137.D32 = Level(n=5, L=2, J=1.5, S=0.5)

The \(\left|n=5, S=1/2, L=2, J=3/2\right>\) level.

atomic_physics.ions.ba137.D52 = Level(n=5, L=2, J=2.5, S=0.5)

The \(\left|n=5, S=1/2, L=2, J=5/2\right>\) level.

atomic_physics.ions.ba137.P12 = Level(n=6, L=1, J=0.5, S=0.5)

The \(\left|n=6, S=1/2, L=1, J=1/2\right>\) level.

atomic_physics.ions.ba137.P32 = Level(n=6, L=1, J=1.5, S=0.5)

The \(\left|n=6, S=1/2, L=1, J=3/2\right>\) level.

atomic_physics.ions.ba137.S12 = Level(n=6, L=0, J=0.5, S=0.5)

The \(\left|n=6, S=1/2, L=0, J=1/2\right>\) level.

atomic_physics.ions.ba137.ground_level = Level(n=6, L=0, J=0.5, S=0.5)

Alias for the \(\left|n=6, S=1/2, L=0, J=1/2\right>\) ground level of \(^{137}\mathrm{Ba}^+\).

atomic_physics.ions.ba137.shelf = Level(n=5, L=2, J=2.5, S=0.5)

Alias for the \(\left|n=5, S=1/2, L=2, J=5/2\right>\) “shelf” level of \(^{137}\mathrm{Ba}^+\).

Barium 138

\(^{138}\mathrm{Ba}^+\)

Where no references are given next to the transition frequencies, they were calculated based on transition frequencies between other levels. For example the frequency of the 1762 nm transition f(1762) was found from f(1762)=f(455)-f(614).

References:

* [1] A. Kramida, NIST Atomic Spectra Database (ver. 5.9) (2021)
* [2] Zhiqiang Zhang, K. J. Arnold, S. R. Chanu, R. Kaewuam,
  M. S. Safronova, and M. D. Barrett Phys. Rev. A 101, 062515 (2020)
* [3] N. Yu, W. Nagourney, and H. Dehmelt, Phys. Rev. Lett. 78, 4898 (1997)
atomic_physics.ions.ba138.Ba138 = AtomFactory(level_data=(LevelData(level=Level(n=6, L=0, J=0.5, S=0.5), Ahfs=0, Bhfs=0, g_J=2.00231930436092, g_I=0.0), LevelData(level=Level(n=6, L=1, J=0.5, S=0.5), Ahfs=0, Bhfs=0, g_J=0.6658935652130266, g_I=0.0), LevelData(level=Level(n=6, L=1, J=1.5, S=0.5), Ahfs=0, Bhfs=0, g_J=1.3341064347869733, g_I=0.0), LevelData(level=Level(n=5, L=2, J=1.5, S=0.5), Ahfs=0, Bhfs=0, g_J=0.799536139127816, g_I=0.0), LevelData(level=Level(n=5, L=2, J=2.5, S=0.5), Ahfs=0, Bhfs=0, g_J=1.200463860872184, g_I=0.0)), transitions={'493': Transition(lower=Level(n=6, L=0, J=0.5, S=0.5), upper=Level(n=6, L=1, J=0.5, S=0.5), frequency=3816572113377394.0, einstein_A=95300000.0), '455': Transition(lower=Level(n=6, L=0, J=0.5, S=0.5), upper=Level(n=6, L=1, J=1.5, S=0.5), frequency=4135068020079713.5, einstein_A=111000000.0), '650': Transition(lower=Level(n=5, L=2, J=1.5, S=0.5), upper=Level(n=6, L=1, J=0.5, S=0.5), frequency=2898508217514255.0, einstein_A=31000000.0), '585': Transition(lower=Level(n=5, L=2, J=1.5, S=0.5), upper=Level(n=6, L=1, J=1.5, S=0.5), frequency=3217004124216574.5, einstein_A=6000000.0), '614': Transition(lower=Level(n=5, L=2, J=2.5, S=0.5), upper=Level(n=6, L=1, J=1.5, S=0.5), frequency=3066132110107188.5, einstein_A=41200000.0), '1762': Transition(lower=Level(n=6, L=0, J=0.5, S=0.5), upper=Level(n=5, L=2, J=2.5, S=0.5), frequency=1068935909972525.1, einstein_A=0.033178500331785), '2051': Transition(lower=Level(n=6, L=0, J=0.5, S=0.5), upper=Level(n=5, L=2, J=1.5, S=0.5), frequency=918063895863139.0, einstein_A=0.0125)}, nuclear_spin=0.0, level_states={Level(n=6, L=1, J=1.5, S=0.5): LevelStates(frequency=4135068020079713.5, start_index=0, stop_index=4, num_states=4), Level(n=6, L=1, J=0.5, S=0.5): LevelStates(frequency=3816572113377394.0, start_index=4, stop_index=6, num_states=2), Level(n=5, L=2, J=2.5, S=0.5): LevelStates(frequency=1068935909972525.1, start_index=6, stop_index=12, num_states=6), Level(n=5, L=2, J=1.5, S=0.5): LevelStates(frequency=918063895863139.0, start_index=12, stop_index=16, num_states=4), Level(n=6, L=0, J=0.5, S=0.5): LevelStates(frequency=0, start_index=16, stop_index=18, num_states=2)})

\(^{138}\mathrm{Ba}^+\) atomic structure.

atomic_physics.ions.ba138.D32 = Level(n=5, L=2, J=1.5, S=0.5)

The \(\left|n=5, S=1/2, L=2, J=3/2\right>\) level.

atomic_physics.ions.ba138.D52 = Level(n=5, L=2, J=2.5, S=0.5)

The \(\left|n=5, S=1/2, L=2, J=5/2\right>\) level.

atomic_physics.ions.ba138.P12 = Level(n=6, L=1, J=0.5, S=0.5)

The \(\left|n=6, S=1/2, L=1, J=1/2\right>\) level.

atomic_physics.ions.ba138.P32 = Level(n=6, L=1, J=1.5, S=0.5)

The \(\left|n=6, S=1/2, L=1, J=3/2\right>\) level.

atomic_physics.ions.ba138.S12 = Level(n=6, L=0, J=0.5, S=0.5)

The \(\left|n=6, S=1/2, L=0, J=1/2\right>\) level.

atomic_physics.ions.ba138.ground_level = Level(n=6, L=0, J=0.5, S=0.5)

Alias for the \(\left|n=6, S=1/2, L=0, J=1/2\right>\) ground level of \(^{138}\mathrm{Ba}^+\).

atomic_physics.ions.ba138.shelf = Level(n=5, L=2, J=2.5, S=0.5)

Alias for the \(\left|n=5, S=1/2, L=2, J=5/2\right>\) “shelf” level of \(^{138}\mathrm{Ba}^+\).

Calcium 40

\(^{40}\mathrm{Ca}^+\)

References:

* [1] A. Kramida, At. Data Nucl. Data Tables 133-134, 101322 (2020)
* [2] T. P. Harty, DPhil Thesis (2013)
* [3] M. Chwalla et all, PRL 102, 023002 (2009)
atomic_physics.ions.ca40.Ca40 = AtomFactory(level_data=(LevelData(level=Level(n=4, L=0, J=0.5, S=0.5), Ahfs=0, Bhfs=0, g_J=2.00225664, g_I=0.0), LevelData(level=Level(n=4, L=1, J=0.5, S=0.5), Ahfs=0, Bhfs=0, g_J=0.6658935652130266, g_I=0.0), LevelData(level=Level(n=4, L=1, J=1.5, S=0.5), Ahfs=0, Bhfs=0, g_J=1.3341064347869733, g_I=0.0), LevelData(level=Level(n=3, L=2, J=1.5, S=0.5), Ahfs=0, Bhfs=0, g_J=0.799536139127816, g_I=0.0), LevelData(level=Level(n=3, L=2, J=2.5, S=0.5), Ahfs=0, Bhfs=0, g_J=1.200334, g_I=0.0)), transitions={'397': Transition(lower=Level(n=4, L=0, J=0.5, S=0.5), upper=Level(n=4, L=1, J=0.5, S=0.5), frequency=4745204585539877.0, einstein_A=132000000.0), '393': Transition(lower=Level(n=4, L=0, J=0.5, S=0.5), upper=Level(n=4, L=1, J=1.5, S=0.5), frequency=4787190380628514.0, einstein_A=135000000.0), '866': Transition(lower=Level(n=3, L=2, J=1.5, S=0.5), upper=Level(n=4, L=1, J=0.5, S=0.5), frequency=2173983592933215.0, einstein_A=8400000.0), '850': Transition(lower=Level(n=3, L=2, J=1.5, S=0.5), upper=Level(n=4, L=1, J=1.5, S=0.5), frequency=2215969388021852.0, einstein_A=955000.0), '854': Transition(lower=Level(n=3, L=2, J=2.5, S=0.5), upper=Level(n=4, L=1, J=1.5, S=0.5), frequency=2204536510188146.5, einstein_A=8480000.0), '729': Transition(lower=Level(n=4, L=0, J=0.5, S=0.5), upper=Level(n=3, L=2, J=2.5, S=0.5), frequency=2582653870442892.0, einstein_A=0.856), '733': Transition(lower=Level(n=4, L=0, J=0.5, S=0.5), upper=Level(n=3, L=2, J=1.5, S=0.5), frequency=2571220992605833.0, einstein_A=0.85)}, nuclear_spin=0.0, level_states={Level(n=4, L=1, J=1.5, S=0.5): LevelStates(frequency=4787190380627685.0, start_index=0, stop_index=4, num_states=4), Level(n=4, L=1, J=0.5, S=0.5): LevelStates(frequency=4745204585539048.0, start_index=4, stop_index=6, num_states=2), Level(n=3, L=2, J=2.5, S=0.5): LevelStates(frequency=2582653870442892.0, start_index=6, stop_index=12, num_states=6), Level(n=3, L=2, J=1.5, S=0.5): LevelStates(frequency=2571220992605833.0, start_index=12, stop_index=16, num_states=4), Level(n=4, L=0, J=0.5, S=0.5): LevelStates(frequency=0, start_index=16, stop_index=18, num_states=2)})

\(^{40}\mathrm{Ca}^+\) atomic structure.

atomic_physics.ions.ca40.D32 = Level(n=3, L=2, J=1.5, S=0.5)

The \(\left|n=3, S=1/2, L=2, J=3/2\right>\) level.

atomic_physics.ions.ca40.D52 = Level(n=3, L=2, J=2.5, S=0.5)

The \(\left|n=3, S=1/2, L=2, J=5/2\right>\) level.

atomic_physics.ions.ca40.P12 = Level(n=4, L=1, J=0.5, S=0.5)

The \(\left|n=4, S=1/2, L=1, J=1/2\right>\) level.

atomic_physics.ions.ca40.P32 = Level(n=4, L=1, J=1.5, S=0.5)

The \(\left|n=4, S=1/2, L=1, J=3/2\right>\) level.

atomic_physics.ions.ca40.S12 = Level(n=4, L=0, J=0.5, S=0.5)

The \(\left|n=4, S=1/2, L=0, J=1/2\right>\) level.

atomic_physics.ions.ca40.ground_level = Level(n=4, L=0, J=0.5, S=0.5)

Alias for the \(\left|n=4, S=1/2, L=0, J=1/2\right>\) ground level of \(^{40}\mathrm{Ca}^+\).

atomic_physics.ions.ca40.shelf = Level(n=3, L=2, J=2.5, S=0.5)

Alias for the \(\left|n=3, S=1/2, L=2, J=5/2\right>\) “shelf” level of \(^{40}\mathrm{Ca}^+\).

Calcium 43

\(^{43}\mathrm{Ca}^+\)

References:

* [1] F. Arbes, et al., Zeitschrift fur Physik D: Atoms, Molecules and
  Clusters, 31, 27 (1994)
* [2] G. Tommaseo, et al., The European Physical Journal D, 25 (2003)
* [3] T. P. Harty, et al. Phys. Rev. Lett 113, 220501 (2014)
* [4] W.  Nortershauser, et al., The European Physical Journal D, 2 (1998)
* [5] J. Benhelm, et al., PHYSICAL REVIEW A 75, 032506 (2007)
* [6] A. Kramida, At. Data Nucl. Data Tables 133-134, 101322 (2020)
atomic_physics.ions.ca43.Ca43 = AtomFactory(level_data=(LevelData(level=Level(n=4, L=0, J=0.5, S=0.5), Ahfs=-5.343276695526923e-25, Bhfs=0, g_J=2.00225664, g_I=-0.37581371428571425), LevelData(level=Level(n=4, L=1, J=0.5, S=0.5), Ahfs=-9.6343059981e-26, Bhfs=0, g_J=0.6658935652130266, g_I=-0.37581371428571425), LevelData(level=Level(n=4, L=1, J=1.5, S=0.5), Ahfs=-2.0805860271e-26, Bhfs=-4.5719884034999995e-27, g_J=1.3341064347869733, g_I=-0.37581371428571425), LevelData(level=Level(n=3, L=2, J=1.5, S=0.5), Ahfs=-3.13413118095e-26, Bhfs=-2.4516459555e-27, g_J=0.799536139127816, g_I=-0.37581371428571425), LevelData(level=Level(n=3, L=2, J=2.5, S=0.5), Ahfs=-2.5795953700964998e-27, Bhfs=-2.810116350615e-27, g_J=1.2003, g_I=-0.37581371428571425)), transitions={'397': Transition(lower=Level(n=4, L=0, J=0.5, S=0.5), upper=Level(n=4, L=1, J=0.5, S=0.5), frequency=4745208845784560.0, einstein_A=132000000.0), '393': Transition(lower=Level(n=4, L=0, J=0.5, S=0.5), upper=Level(n=4, L=1, J=1.5, S=0.5), frequency=4787194645660984.0, einstein_A=135000000.0), '866': Transition(lower=Level(n=3, L=2, J=1.5, S=0.5), upper=Level(n=4, L=1, J=0.5, S=0.5), frequency=2173961839062849.8, einstein_A=8400000.0), '850': Transition(lower=Level(n=3, L=2, J=1.5, S=0.5), upper=Level(n=4, L=1, J=1.5, S=0.5), frequency=2215947638939274.0, einstein_A=955000.0), '854': Transition(lower=Level(n=3, L=2, J=2.5, S=0.5), upper=Level(n=4, L=1, J=1.5, S=0.5), frequency=2204514796046362.0, einstein_A=8480000.0), '729': Transition(lower=Level(n=4, L=0, J=0.5, S=0.5), upper=Level(n=3, L=2, J=2.5, S=0.5), frequency=2582679849602684.0, einstein_A=0.856), '733': Transition(lower=Level(n=4, L=0, J=0.5, S=0.5), upper=Level(n=3, L=2, J=1.5, S=0.5), frequency=2571247006721710.0, einstein_A=0.85)}, nuclear_spin=3.5, level_states={Level(n=4, L=1, J=1.5, S=0.5): LevelStates(frequency=4787194645660984.0, start_index=0, stop_index=32, num_states=32), Level(n=4, L=1, J=0.5, S=0.5): LevelStates(frequency=4745208845784560.0, start_index=32, stop_index=48, num_states=16), Level(n=3, L=2, J=2.5, S=0.5): LevelStates(frequency=2582679849602684.0, start_index=48, stop_index=96, num_states=48), Level(n=3, L=2, J=1.5, S=0.5): LevelStates(frequency=2571247006721710.0, start_index=96, stop_index=128, num_states=32), Level(n=4, L=0, J=0.5, S=0.5): LevelStates(frequency=0, start_index=128, stop_index=144, num_states=16)})

\(^{43}\mathrm{Ca}^+\) atomic structure.

atomic_physics.ions.ca43.D32 = Level(n=3, L=2, J=1.5, S=0.5)

The \(\left|n=3, S=1/2, L=2, J=3/2\right>\) level.

atomic_physics.ions.ca43.D52 = Level(n=3, L=2, J=2.5, S=0.5)

The \(\left|n=3, S=1/2, L=2, J=5/2\right>\) level.

atomic_physics.ions.ca43.P12 = Level(n=4, L=1, J=0.5, S=0.5)

The \(\left|n=4, S=1/2, L=1, J=1/2\right>\) level.

atomic_physics.ions.ca43.P32 = Level(n=4, L=1, J=1.5, S=0.5)

The \(\left|n=4, S=1/2, L=1, J=3/2\right>\) level.

atomic_physics.ions.ca43.S12 = Level(n=4, L=0, J=0.5, S=0.5)

The \(\left|n=4, S=1/2, L=0, J=1/2\right>\) level.

atomic_physics.ions.ca43.ground_level = Level(n=4, L=0, J=0.5, S=0.5)

Alias for the \(\left|n=4, S=1/2, L=0, J=1/2\right>\) ground level of \(^{43}\mathrm{Ca}^+\).

atomic_physics.ions.ca43.shelf = Level(n=3, L=2, J=2.5, S=0.5)

Alias for the \(\left|n=3, S=1/2, L=2, J=5/2\right>\) “shelf” level of \(^{43}\mathrm{Ca}^+\).

Magnesium 25

\(^{25}\mathrm{Mg}^+\)

References:

* [1] W. M. Itano and D. J. Wineland, Precision measurement of the
  ground-state hyperfine constant of Mg+, PRA, 24, 3 (1981)
* [2] W. H. Yuan et. al., Precision measurement of the light shift of
  25Mg+ions, Phys. Rev. A, 98, 5 (2018)
* [3] G. Clos et. al., Decoherence-Assisted Spectroscopy of a Single
  Mg+ Ion, Phys. Rev. Lett., 112, 11 (2014)
* [4] M. Kaur et. al., Radiative transition properties of singly charged
  magnesium, calcium, strontium and barium ions, Atomic Data and Nuclear
  Data Tables, 137 (2021)
* [5] Z. T. Xu et. al., Precision measurement of the 25Mg+
  ground-state hyperfine constant, Phys. Rev. A, 96, 5, (2017)
* [6] J. Nguyen, The Linewidth and Hyperfine A Constant of the 2P1/2 State
  of a Magnesium Ion Confined in a Linear Paul Trap, Thesis,
  McMaster University (2009) http://hdl.handle.net/11375/17398
* [7]  N. J. Stone, Table of nuclear magnetic dipole and electric
  quadrupole moments, Atomic Data and Nuclear Data Tables, Volume 90,
  Issue 1 (2005)
* [8] L. Toppozini, Trapped-Mg+ Apparatus for Control and Structure Studies,
  Thesis, McMaster University (2006) http://hdl.handle.net/11375/21342
atomic_physics.ions.mg25.Mg25 = AtomFactory(level_data=(LevelData(level=Level(n=3, L=0, J=0.5, S=0.5), Ahfs=-3.9508224791217467e-25, Bhfs=0, g_J=2.002, g_I=-0.34218000000000004), LevelData(level=Level(n=3, L=1, J=0.5, S=0.5), Ahfs=6.769193265239999e-26, Bhfs=0, g_J=0.6658935652130266, g_I=-0.34218000000000004), LevelData(level=Level(n=3, L=1, J=1.5, S=0.5), Ahfs=-1.2653938686858e-26, Bhfs=1.48035021042195e-26, g_J=1.3341064347869733, g_I=-0.34218000000000004)), transitions={'280': Transition(lower=Level(n=3, L=0, J=0.5, S=0.5), upper=Level(n=3, L=1, J=0.5, S=0.5), frequency=6718861106202451.0, einstein_A=558000000.0), '279': Transition(lower=Level(n=3, L=0, J=0.5, S=0.5), upper=Level(n=3, L=1, J=1.5, S=0.5), frequency=6736105873764683.0, einstein_A=260000000.0)}, nuclear_spin=2.5, level_states={Level(n=3, L=1, J=1.5, S=0.5): LevelStates(frequency=6736105873764683.0, start_index=0, stop_index=24, num_states=24), Level(n=3, L=1, J=0.5, S=0.5): LevelStates(frequency=6718861106202451.0, start_index=24, stop_index=36, num_states=12), Level(n=3, L=0, J=0.5, S=0.5): LevelStates(frequency=0, start_index=36, stop_index=48, num_states=12)})

\(^{25}\mathrm{Mg}^+\) atomic structure.

atomic_physics.ions.mg25.P12 = Level(n=3, L=1, J=0.5, S=0.5)

The \(\left|n=3, S=1/2, L=1, J=1/2\right>\) level.

atomic_physics.ions.mg25.P32 = Level(n=3, L=1, J=1.5, S=0.5)

The \(\left|n=3, S=1/2, L=1, J=3/2\right>\) level.

atomic_physics.ions.mg25.S12 = Level(n=3, L=0, J=0.5, S=0.5)

The \(\left|n=3, S=1/2, L=0, J=1/2\right>\) level.

atomic_physics.ions.mg25.ground_level = Level(n=3, L=0, J=0.5, S=0.5)

Alias for the \(\left|n=3, S=1/2, L=0, J=1/2\right>\) ground level of \(^{25}\mathrm{Mg}^+\).

Strontium 88

\(^{88}\mathrm{Sr}^+\)

References:

* [1] A. Kramida, NIST Atomic Spectra Database (ver. 5.9) (2021)
* [2] P. Dubé, Metrologia (2015)
* [3] V. Letchumanan, Phys. Rev. A (2005)
atomic_physics.ions.sr88.D32 = Level(n=5, L=2, J=1.5, S=0.5)

The \(\left|n=4, S=1/2, L=2, J=3/2\right>\) level.

atomic_physics.ions.sr88.D52 = Level(n=4, L=2, J=2.5, S=0.5)

The \(\left|n=4, S=1/2, L=2, J=5/2\right>\) level.

atomic_physics.ions.sr88.P12 = Level(n=5, L=1, J=0.5, S=0.5)

The \(\left|n=5, S=1/2, L=1, J=1/2\right>\) level.

atomic_physics.ions.sr88.P32 = Level(n=5, L=1, J=1.5, S=0.5)

The \(\left|n=5, S=1/2, L=1, J=3/2\right>\) level.

atomic_physics.ions.sr88.S12 = Level(n=5, L=0, J=0.5, S=0.5)

The \(\left|n=5, S=1/2, L=0, J=1/2\right>\) level.

atomic_physics.ions.sr88.Sr88 = AtomFactory(level_data=(LevelData(level=Level(n=5, L=0, J=0.5, S=0.5), Ahfs=0, Bhfs=0, g_J=2.00231930436092, g_I=0.0), LevelData(level=Level(n=5, L=1, J=0.5, S=0.5), Ahfs=0, Bhfs=0, g_J=0.6658935652130266, g_I=0.0), LevelData(level=Level(n=5, L=1, J=1.5, S=0.5), Ahfs=0, Bhfs=0, g_J=1.3341064347869733, g_I=0.0), LevelData(level=Level(n=5, L=2, J=1.5, S=0.5), Ahfs=0, Bhfs=0, g_J=0.799536139127816, g_I=0.0), LevelData(level=Level(n=4, L=2, J=2.5, S=0.5), Ahfs=0, Bhfs=0, g_J=1.200463860872184, g_I=0.0)), transitions={'422': Transition(lower=Level(n=5, L=0, J=0.5, S=0.5), upper=Level(n=5, L=1, J=0.5, S=0.5), frequency=4468368742080117.0, einstein_A=127900000.0), '408': Transition(lower=Level(n=5, L=0, J=0.5, S=0.5), upper=Level(n=5, L=1, J=1.5, S=0.5), frequency=4619381269281887.0, einstein_A=141000000.0), '1092': Transition(lower=Level(n=5, L=2, J=1.5, S=0.5), upper=Level(n=5, L=1, J=0.5, S=0.5), frequency=1725765746181632.2, einstein_A=7460000.0), '1004': Transition(lower=Level(n=5, L=2, J=1.5, S=0.5), upper=Level(n=5, L=1, J=1.5, S=0.5), frequency=1876772445586797.5, einstein_A=1000000.0), '1033': Transition(lower=Level(n=4, L=2, J=2.5, S=0.5), upper=Level(n=5, L=1, J=1.5, S=0.5), frequency=1823951977527592.8, einstein_A=8700000.0), '674': Transition(lower=Level(n=5, L=0, J=0.5, S=0.5), upper=Level(n=4, L=2, J=2.5, S=0.5), frequency=2794629154802134.0, einstein_A=2.55885), '687': Transition(lower=Level(n=5, L=0, J=0.5, S=0.5), upper=Level(n=5, L=2, J=1.5, S=0.5), frequency=2742581055871865.5, einstein_A=2.299)}, nuclear_spin=0, level_states={Level(n=5, L=1, J=1.5, S=0.5): LevelStates(frequency=4619353501458663.0, start_index=0, stop_index=4, num_states=4), Level(n=5, L=1, J=0.5, S=0.5): LevelStates(frequency=4468346802053498.0, start_index=4, stop_index=6, num_states=2), Level(n=4, L=2, J=2.5, S=0.5): LevelStates(frequency=2794629154802134.0, start_index=6, stop_index=12, num_states=6), Level(n=5, L=2, J=1.5, S=0.5): LevelStates(frequency=2742581055871865.5, start_index=12, stop_index=16, num_states=4), Level(n=5, L=0, J=0.5, S=0.5): LevelStates(frequency=0, start_index=16, stop_index=18, num_states=2)})

\(^{88}\mathrm{Sr}^+\) atomic structure.

atomic_physics.ions.sr88.ground_level = Level(n=5, L=0, J=0.5, S=0.5)

Alias for the \(\left|n=5, S=1/2, L=0, J=1/2\right>\) ground level of \(^{88}\mathrm{Sr}^+\).

atomic_physics.ions.sr88.shelf = Level(n=4, L=2, J=2.5, S=0.5)

Alias for the \(\left|n=4, S=1/2, L=2, J=5/2\right>\) “shelf” level of \(^{88}\mathrm{Sr}^+\).