Alkali atom functions

Overview

Classes and global methods

AlkaliAtom([preferQuantumDefects, cpp_numerov]) Implements general calculations for alkali atoms.
NumerovBack(innerLimit, outerLimit, kfun, ...) Full Python implementation of Numerov integration
saveCalculation(calculation, fileName) Saves calculation for future use.
loadSavedCalculation(fileName) Loads previously saved calculation.
printState(n, l, j) Prints state spectroscopic label for numeric \(n\),
printStateString(n, l, j) Returns state spectroscopic label for numeric \(n\),

AlkaliAtom Methods

AlkaliAtom.getDipoleMatrixElement(n1, l1, ...) Dipole matrix element
AlkaliAtom.getTransitionWavelength(n1, l1, ...) Calculated transition wavelength (in vacuum) in m.
AlkaliAtom.getTransitionFrequency(n1, l1, ...) Calculated transition frequency in Hz
AlkaliAtom.getRabiFrequency(n1, l1, j1, mj1, ...) Returns a Rabi frequency for resonantly driven atom in a
AlkaliAtom.getRabiFrequency2(n1, l1, j1, ...) Returns a Rabi frequency for resonant excitation with a given
AlkaliAtom.getStateLifetime(n, l, j[, ...]) Returns the lifetime of the state (in s)
AlkaliAtom.getTransitionRate(n1, l1, j1, n2, ...) Transition rate due to coupling to vacuum modes (black body included)
AlkaliAtom.getReducedMatrixElementJ_asymmetric(n1, ...) Reduced matrix element in \(J\) basis, defined in asymmetric notation.
AlkaliAtom.getReducedMatrixElementJ(n1, l1, ...) Reduced matrix element in \(J\) basis (symmetric notation)
AlkaliAtom.getReducedMatrixElementL(n1, l1, ...) Reduced matrix element in \(L\) basis (symmetric notation)
AlkaliAtom.getRadialMatrixElement(n1, l1, ...) Radial part of the dipole matrix element
AlkaliAtom.getQuadrupoleMatrixElement(n1, ...) Radial part of the quadrupole matrix element
AlkaliAtom.getPressure(temperature) Vapour pressure (in Pa) at given temperature
AlkaliAtom.getNumberDensity(temperature) Atom number density at given temperature
AlkaliAtom.getAverageInteratomicSpacing(...) Returns average interatomic spacing in atomic vapour
AlkaliAtom.corePotential(l, r) core potential felt by valence electron
AlkaliAtom.effectiveCharge(l, r) effective charge of the core felt by valence electron
AlkaliAtom.potential(l, s, j, r) returns total potential that electron feels
AlkaliAtom.radialWavefunction(l, s, j, ...) Radial part of electron wavefunction
AlkaliAtom.getEnergy(n, l, j) Energy of the level relative to the ionisation level (in eV)
AlkaliAtom.getQuantumDefect(n, l, j) Quantum defect of the level.
AlkaliAtom.getC6term(n, l, j, n1, l1, j1, ...) C6 interaction term for the given two pair-states
AlkaliAtom.getC3term(n, l, j, n1, l1, j1, ...) C3 interaction term for the given two pair-states
AlkaliAtom.getEnergyDefect(n, l, j, n1, l1, ...) Energy defect for the given two pair-states (one of the state has
AlkaliAtom.getEnergyDefect2(n, l, j, nn, ll, ...) Energy defect for the given two pair-states
AlkaliAtom.updateDipoleMatrixElementsFile() Updates the file with pre-calculated dipole matrix elements.
AlkaliAtom.getRadialCoupling(n, l, j, n1, l1, j1) Returns radial part of the coupling between two states (dipole and
AlkaliAtom.getAverageSpeed(temperature) Average (mean) speed at a given temperature
AlkaliAtom.getLiteratureDME(n1, l1, j1, n2, ...) Returns literature information on requested transition.

Detailed documentation

Implements general single-atom calculations

This module calculates single (isolated) atom properties of all alkali metals in general. For example, it calculates dipole matrix elements, quandrupole matrix elements, etc. Also, some helpful general functions are here, e.g. for saving and loading calculations (single-atom and pair-state based), printing state labels etc.

class arc.alkali_atom_functions.AlkaliAtom(preferQuantumDefects=True, cpp_numerov=True)[source]

Implements general calculations for alkali atoms.

This abstract class implements general calculations methods.

Parameters:
  • preferQuantumDefects (bool) – Use quantum defects for energy level calculations. If False, uses NIST ASD values where available. If True, uses quantum defects for energy calculations for principal quantum numbers equal or above minQuantumDefectN which is specified for each element separately. For principal quantum numbers below this value, NIST ASD values are used, since quantum defects don’t reproduce well low-lying states. Default is True.
  • cpp_numerov (bool) – should the wavefunction be calculated with Numerov algorithm implemented in C++; if False, it uses pure Python implementation that is much slower. Default is True.
Z = 0.0

Atomic number

abundance = 1.0

relative isotope abundance

alphaC = 0.0

Core polarizability

corePotential(l, r)[source]

core potential felt by valence electron

For more details about derivation of model potential see Ref. [2].

Parameters:
  • l (int) – orbital angular momentum
  • r (float) – distance from the nucleus (in a.u.)
Returns:

core potential felt by valence electron (in a.u. ???)

Return type:

float

References

[2](1, 2) M. Marinescu, H. R. Sadeghpour, and A. Dalgarno PRA 49, 982 (1994), https://doi.org/10.1103/PhysRevA.49.982
cpp_numerov = True

swich - should the wavefunction be calculated with Numerov algorithm implemented in C++

dipoleMatrixElementFile = ''

location of hard-disk stored dipole matrix elements

effectiveCharge(l, r)[source]

effective charge of the core felt by valence electron

For more details about derivation of model potential see Ref. [2].

Parameters:
  • l (int) – orbital angular momentum
  • r (float) – distance from the nucleus (in a.u.)
Returns:

effective charge (in a.u.)

Return type:

float

elementName = 'elementName'

Human-readable element name

extraLevels = []

levels that are for smaller principal quantum number (n) than ground level, but are above in energy due to angular part

getAverageInteratomicSpacing(temperature)[source]

Returns average interatomic spacing in atomic vapour

See calculation of basic properties example snippet.

Parameters:temperature (float) – temperature of the atomic vapour
Returns:average interatomic spacing in m
Return type:float
getAverageSpeed(temperature)[source]

Average (mean) speed at a given temperature

Parameters:temperature (float) – temperature (K)
Returns:mean speed (m/s)
Return type:float
getC3term(n, l, j, n1, l1, j1, n2, l2, j2)[source]

C3 interaction term for the given two pair-states

Calculates \(C_3\) intaraction term for \(|n,l,j,n,l,j\rangle \leftrightarrow |n_1,l_1,j_1,n_2,l_2,j_2\rangle\)

Parameters:
  • n (int) – principal quantum number
  • l (int) – orbital angular momenutum
  • j (float) – total angular momentum
  • n1 (int) – principal quantum number
  • l1 (int) – orbital angular momentum
  • j1 (float) – total angular momentum
  • n2 (int) – principal quantum number
  • l2 (int) – orbital angular momentum
  • j2 (float) – total angular momentum
Returns:

\(C_3 = \frac{\langle n,l,j |er|n_1,l_1,j_1\rangle \langle n,l,j |er|n_2,l_2,j_2\rangle}{4\pi\varepsilon_0}\) (\(h\) Hz m \({}^3\)).

Return type:

float

getC6term(n, l, j, n1, l1, j1, n2, l2, j2)[source]

C6 interaction term for the given two pair-states

Calculates \(C_6\) intaraction term for \(|n,l,j,n,l,j\rangle \leftrightarrow |n_1,l_1,j_1,n_2,l_2,j_2\rangle\). For details of calculation see Ref. [3].

Parameters:
  • n (int) – principal quantum number
  • l (int) – orbital angular momenutum
  • j (float) – total angular momentum
  • n1 (int) – principal quantum number
  • l1 (int) – orbital angular momentum
  • j1 (float) – total angular momentum
  • n2 (int) – principal quantum number
  • l2 (int) – orbital angular momentum
  • j2 (float) – total angular momentum
Returns:

\(C_6 = \frac{1}{4\pi\varepsilon_0} \frac{|\langle n,l,j |er|n_1,l_1,j_1\rangle|^2| \langle n,l,j |er|n_2,l_2,j_2\rangle|^2} {E(n_1,l_1,j_2,n_2,j_2,j_2)-E(n,l,j,n,l,j)}\) (\(h\) Hz m \({}^6\)).

Return type:

float

Example

We can reproduce values from Ref. [3] for C3 coupling to particular channels. Taking for example channels described by the Eq. (50a-c) we can get the values:

from arc import *

channels = [[70,0,0.5, 70, 1,1.5, 69,1, 1.5],\
            [70,0,0.5, 70, 1,1.5, 69,1, 0.5],\
            [70,0,0.5, 69, 1,1.5, 70,1, 0.5],\
            [70,0,0.5, 70, 1,0.5, 69,1, 0.5]]

print(" = = = Caesium = = = ")
atom = Caesium()
for channel in channels:
    print("%.0f  GHz (mu m)^6" % ( atom.getC6term(*channel)/h*1.e27 ))

print("\n = = = Rubidium  = = =")
atom = Rubidium()
for channel in channels:
    print("%.0f  GHz (mu m)^6" % ( atom.getC6term(*channel)/h*1.e27 ))

Returns:

 = = = Caesium = = =
722  GHz (mu m)^6
316  GHz (mu m)^6
383  GHz (mu m)^6
228  GHz (mu m)^6

 = = = Rubidium  = = =
799  GHz (mu m)^6
543  GHz (mu m)^6
589  GHz (mu m)^6
437  GHz (mu m)^6

which is in good agreement with the values cited in the Ref. [3]. Small discrepancies for Caesium originate from slightly different quantum defects used in calculations.

References

[3](1, 2, 3) T. G. Walker, M. Saffman, PRA 77, 032723 (2008) https://doi.org/10.1103/PhysRevA.77.032723
getDipoleMatrixElement(n1, l1, j1, mj1, n2, l2, j2, mj2, q)[source]

Dipole matrix element \(\langle n_1 l_1 j_1 m_{j_1} |e\mathbf{r}|n_2 l_2 j_2 m_{j_2}\rangle\) in units of \(a_0 e\)

Returns:dipole matrix element( \(a_0 e\))
Return type:float

Example

For example, calculation of \(5 S_{1/2}m_j=-\frac{1}{2} \rightarrow 5 P_{3/2}m_j=-\frac{3}{2}\) transition dipole matrix element for laser driving \(\sigma^-\) transition:

from arc import *
atom = Rubidium()
# transition 5 S_{1/2} m_j=-0.5 -> 5 P_{3/2} m_j=-1.5 for laser
# driving sigma- transition
print(atom.getDipoleMatrixElement(5,0,0.5,-0.5,5,1,1.5,-1.5,-1))
getEnergy(n, l, j)[source]

Energy of the level relative to the ionisation level (in eV)

Returned energies are with respect to the center of gravity of the hyperfine-split states. If preferQuantumDefects =False (set during initialization) program will try use NIST energy value, if such exists, falling back to energy calculation with quantum defects if the measured value doesn’t exist. For preferQuantumDefects =True, program will always calculate energies from quantum defects (useful for comparing quantum defect calculations with measured energy level values).

Parameters:
  • n (int) – principal quantum number
  • l (int) – orbital angular momentum
  • j (float) – total angular momentum
Returns:

state energy (eV)

Return type:

float

getEnergyDefect(n, l, j, n1, l1, j1, n2, l2, j2)[source]

Energy defect for the given two pair-states (one of the state has two atoms in the same state)

Energy difference between the states \(E(n,l,j,n,l,j) - E(n_1,l_1,j_1,n_2,l_2,j_2)\)

Parameters:
  • n (int) – principal quantum number
  • l (int) – orbital angular momenutum
  • j (float) – total angular momentum
  • n1 (int) – principal quantum number
  • l1 (int) – orbital angular momentum
  • j1 (float) – total angular momentum
  • n2 (int) – principal quantum number
  • l2 (int) – orbital angular momentum
  • j2 (float) – total angular momentum
Returns:

energy defect (SI units: J)

Return type:

float

getEnergyDefect2(n, l, j, nn, ll, jj, n1, l1, j1, n2, l2, j2)[source]

Energy defect for the given two pair-states

Energy difference between the states \(E(n,l,j,nn,ll,jj) - E(n_1,l_1,j_1,n_2,l_2,j_2)\)

See pair-state energy defects example snippet.

Parameters:
  • n (int) – principal quantum number
  • l (int) – orbital angular momenutum
  • j (float) – total angular momentum
  • nn (int) – principal quantum number
  • ll (int) – orbital angular momenutum
  • jj (float) – total angular momentum
  • n1 (int) – principal quantum number
  • l1 (int) – orbital angular momentum
  • j1 (float) – total angular momentum
  • n2 (int) – principal quantum number
  • l2 (int) – orbital angular momentum
  • j2 (float) – total angular momentum
Returns:

energy defect (SI units: J)

Return type:

float

getLiteratureDME(n1, l1, j1, n2, l2, j2)[source]

Returns literature information on requested transition.

Parameters:
  • n1,l1,j1 – one of the states we are coupling
  • n2,l2,j2 – the other state to which we are coupling
Returns:

hasLiteratureValue?, dme, referenceInformation

If Boolean value is True, a literature value for dipole matrix element was found and reduced DME in J basis is returned as the number. Third returned argument (array) contains additional information about the literature value in the following order [ typeOfSource, errorEstimate , comment , reference, reference DOI] upon success to find a literature value for dipole matrix element:

  • typeOfSource=1 if the value is theoretical calculation; otherwise, if it is experimentally obtained value typeOfSource=0
  • comment details where within the publication the value can be found
  • errorEstimate is absolute error estimate
  • reference is human-readable formatted reference
  • reference DOI provides link to the publication.

Boolean value is False, followed by zero and an empty array if no literature value for dipole matrix element is found.

Return type:

bool, float, [int,float,string,string,string]

Note

The literature values are stored in /data folder in <element name>_literature_dme.csv files as a ; separated values. Each row in the file consists of one literature entry, that has information in the following order:

  • n1
  • l1
  • j1
  • n2
  • l2
  • j2
  • dipole matrix element reduced l basis (a.u.)
  • comment (e.g. where in the paper value appears?)
  • value origin: 1 for theoretical; 0 for experimental values
  • accuracy
  • source (human readable formatted citation)
  • doi number (e.g. 10.1103/RevModPhys.82.2313 )

If there are several values for a given transition, program will output the value that has smallest error (under column accuracy). The list of values can be expanded - every time program runs this file is read and the list is parsed again for use in calculations.

getNumberDensity(temperature)[source]

Atom number density at given temperature

See calculation of basic properties example snippet.

Parameters:temperature (float) – temperature in K
Returns:atom concentration in \(1/m^3\)
Return type:float
getPressure(temperature)[source]

Vapour pressure (in Pa) at given temperature

Parameters:temperature (float) – temperature in K
Returns:vapour pressure in Pa
Return type:float
getQuadrupoleMatrixElement(n1, l1, j1, n2, l2, j2)[source]

Radial part of the quadrupole matrix element

Calculates \(\int \mathbf{d}r~R_{n_1,l_1,j_1}(r)\cdot R_{n_1,l_1,j_1}(r) \cdot r^4\). See Quadrupole calculation example snippet .

Parameters:
  • n1 (int) – principal quantum number of state 1
  • l1 (int) – orbital angular momentum of state 1
  • j1 (float) – total angular momentum of state 1
  • n2 (int) – principal quantum number of state 2
  • l2 (int) – orbital angular momentum of state 2
  • j2 (float) – total angular momentum of state 2
Returns:

quadrupole matrix element (\(a_0^2 e\)).

Return type:

float

getQuantumDefect(n, l, j)[source]

Quantum defect of the level.

For an example, see Rydberg energy levels example snippet.

Parameters:
  • n (int) – principal quantum number
  • l (int) – orbital angular momentum
  • j (float) – total angular momentum
Returns:

quantum defect

Return type:

float

getRabiFrequency(n1, l1, j1, mj1, n2, l2, j2, q, laserPower, laserWaist)[source]

Returns a Rabi frequency for resonantly driven atom in a center of TEM00 mode of a driving field

Parameters:
  • n1,l1,j1,mj1 – state from which we are driving transition
  • n2,l2,j2 – state to which we are driving transition
  • q – laser polarization (-1,0,1 correspond to \(\sigma^-\), \(\pi\) and \(\sigma^+\) respectively)
  • laserPower – laser power in units of W
  • laserWaist – laser \(1/e^2\) waist (radius) in units of m
Returns:

Frequency in rad \(^{-1}\). If you want frequency in Hz, divide by returned value by \(2\pi\)

Return type:

float

getRabiFrequency2(n1, l1, j1, mj1, n2, l2, j2, q, electricFieldAmplitude)[source]

Returns a Rabi frequency for resonant excitation with a given electric field amplitude

Parameters:
  • n1,l1,j1,mj1 – state from which we are driving transition
  • n2,l2,j2 – state to which we are driving transition
  • q – laser polarization (-1,0,1 correspond to \(\sigma^-\), \(\pi\) and \(\sigma^+\) respectively)
  • electricFieldAmplitude – amplitude of electric field driving (V/m)
Returns:

Frequency in rad \(^{-1}\). If you want frequency in Hz, divide by returned value by \(2\pi\)

Return type:

float

getRadialCoupling(n, l, j, n1, l1, j1)[source]

Returns radial part of the coupling between two states (dipole and quadrupole interactions only)

Parameters:
  • n1 (int) – principal quantum number
  • l1 (int) – orbital angular momentum
  • j1 (float) – total angular momentum
  • n2 (int) – principal quantum number
  • l2 (int) – orbital angular momentum
  • j2 (float) – total angular momentum
Returns:

radial coupling strength (in a.u.), or zero for forbidden transitions in dipole and quadrupole approximation.

Return type:

float

getRadialMatrixElement(n1, l1, j1, n2, l2, j2, useLiterature=True)[source]

Radial part of the dipole matrix element

Calculates \(\int \mathbf{d}r~R_{n_1,l_1,j_1}(r)\cdot R_{n_1,l_1,j_1}(r) \cdot r^3\).

Parameters:
  • n1 (int) – principal quantum number of state 1
  • l1 (int) – orbital angular momentum of state 1
  • j1 (float) – total angular momentum of state 1
  • n2 (int) – principal quantum number of state 2
  • l2 (int) – orbital angular momentum of state 2
  • j2 (float) – total angular momentum of state 2
Returns:

dipole matrix element (\(a_0 e\)).

Return type:

float

getReducedMatrixElementJ(n1, l1, j1, n2, l2, j2)[source]

Reduced matrix element in \(J\) basis (symmetric notation)

Parameters:
  • n1 (int) – principal quantum number of state 1
  • l1 (int) – orbital angular momentum of state 1
  • j1 (float) – total angular momentum of state 1
  • n2 (int) – principal quantum number of state 2
  • l2 (int) – orbital angular momentum of state 2
  • j2 (float) – total angular momentum of state 2
Returns:

reduced dipole matrix element in \(J\) basis \(\langle j || er || j' \rangle\) (\(a_0 e\)).

Return type:

float

getReducedMatrixElementJ_asymmetric(n1, l1, j1, n2, l2, j2)[source]

Reduced matrix element in \(J\) basis, defined in asymmetric notation.

Note that notation for symmetric and asymmetricly defined reduced matrix element is not consistent in the literature. For example, notation is used e.g. in Steck [1] is precisely the oposite.

Note

Note that this notation is asymmetric: \(( j||e r ||j' ) \neq ( j'||e r ||j )\). Relation between the two notation is \(\langle j||er||j'\rangle= \sqrt{2j+1} ( j ||er ||j')\). This function always returns value for transition from lower to higher energy state, independent of the order of states entered in the function call.

Parameters:
  • n1 (int) – principal quantum number of state 1
  • l1 (int) – orbital angular momentum of state 1
  • j1 (float) – total angular momentum of state 1
  • n2 (int) – principal quantum number of state 2
  • l2 (int) – orbital angular momentum of state 2
  • j2 (float) – total angular momentum of state 2
Returns:

reduced dipole matrix element in Steck notation \(( j || er || j' )\) (\(a_0 e\)).

Return type:

float

[1]Daniel A. Steck, “Cesium D Line Data,” (revision 2.0.1, 2 May 2008). http://steck.us/alkalidata
getReducedMatrixElementL(n1, l1, j1, n2, l2, j2)[source]

Reduced matrix element in \(L\) basis (symmetric notation)

Parameters:
  • n1 (int) – principal quantum number of state 1
  • l1 (int) – orbital angular momentum of state 1
  • j1 (float) – total angular momentum of state 1
  • n2 (int) – principal quantum number of state 2
  • l2 (int) – orbital angular momentum of state 2
  • j2 (float) – total angular momentum of state 2
Returns:

reduced dipole matrix element in \(L\) basis \(\langle l || er || l' \rangle\) (\(a_0 e\)).

Return type:

float

getStateLifetime(n, l, j, temperature=0, includeLevelsUpTo=0)[source]

Returns the lifetime of the state (in s)

For non-zero temperatures, user must specify up to which principal quantum number levels, that is above the initial state, should be included in order to account for black-body induced transitions to higher lying states. See Rydberg lifetimes example snippet.

Parameters:
  • l, j (n,) – specifies state whose lifetime we are calculating
  • temperature – optional. Temperature at which the atom environment is, measured in K. If this parameter is non-zero, user has to specify transitions up to which state (due to black-body decay) should be included in calculation.
  • includeLevelsUpTo (int) – optional and not needed for atom lifetimes calculated at zero temperature. At non zero temperatures, this specify maximum principal quantum number of the state to which black-body induced transitions will be included. Minimal value of the parameter in that case is \(n+1\)
Returns:

State lifetime in units of s (seconds)

Return type:

float

See also

getTransitionRate for calculating rates of individual transition rates between the two states

getTransitionFrequency(n1, l1, j1, n2, l2, j2)[source]

Calculated transition frequency in Hz

Returned values is given relative to the centre of gravity of the hyperfine-split states.

Parameters:
  • n1 (int) – principal quantum number of the state from which we are going
  • l1 (int) – orbital angular momentum of the state from which we are going
  • j1 (float) – total angular momentum of the state from which we are going
  • n2 (int) – principal quantum number of the state to which we are going
  • l2 (int) – orbital angular momentum of the state to which we are going
  • j2 (float) – total angular momentum of the state to which we are going
Returns:

transition frequency (in Hz). If the returned value is negative, level from which we are going is above the level to which we are going.

Return type:

float

getTransitionRate(n1, l1, j1, n2, l2, j2, temperature=0.0)[source]

Transition rate due to coupling to vacuum modes (black body included)

Calculates transition rate from the first given state to the second given state \(|n_1,l_1,j_1\rangle \rightarrow |n_2,j_2,j_2\rangle\) at given temperature due to interaction with the vacuum field. For zero temperature this returns Einstein A coefficient. For details of calculation see Ref. [4] and Ref. [5]. See Black-body induced population transfer example snippet.

Parameters:
  • n1 (int) – principal quantum number
  • l1 (int) – orbital angular momentum
  • j1 (float) – total angular momentum
  • n2 (int) – principal quantum number
  • l2 (int) – orbital angular momentum
  • j2 (float) – total angular momentum
  • [temperature] (float) – temperature in K
Returns:

transition rate in s \({}^{-1}\) (SI)

Return type:

float

References

[4]C. E. Theodosiou, PRA 30, 2881 (1984) https://doi.org/10.1103/PhysRevA.30.2881
[5]I. I. Beterov, I. I. Ryabtsev, D. B. Tretyakov, and V. M. Entin, PRA 79, 052504 (2009) https://doi.org/10.1103/PhysRevA.79.052504
getTransitionWavelength(n1, l1, j1, n2, l2, j2)[source]

Calculated transition wavelength (in vacuum) in m.

Returned values is given relative to the centre of gravity of the hyperfine-split states.

Parameters:
  • n1 (int) – principal quantum number of the state from which we are going
  • l1 (int) – orbital angular momentum of the state from which we are going
  • j1 (float) – total angular momentum of the state from which we are going
  • n2 (int) – principal quantum number of the state to which we are going
  • l2 (int) – orbital angular momentum of the state to which we are going
  • j2 (float) – total angular momentum of the state to which we are going
Returns:

transition wavelength (in m). If the returned value is negative, level from which we are going is above the level to which we are going.

Return type:

float

groundStateN = 0

principal quantum number for the ground state

levelDataFromNIST = ''

location of stored NIST values of measured energy levels in eV

literatureDMEfilename = ''

Filename of the additional literature source values of dipole matrix elements.

These additional values should be saved as reduced dipole matrix elements in J basis.

mass = 0.0

atomic mass in kg

minQuantumDefectN = 0

minimal quantum number for which quantum defects can be used; uses measured energy levels otherwise

potential(l, s, j, r)[source]

returns total potential that electron feels

Total potential = core potential + Spin-Orbit interaction

Parameters:
  • l (int) – orbital angular momentum
  • s (float) – spin angular momentum
  • j (float) – total angular momentum
  • r (float) – distance from the nucleus (in a.u.)
Returns:

potential (in a.u.)

Return type:

float

quadrupoleMatrixElementFile = ''

location of hard-disk stored dipole matrix elements

quantumDefect = [[[0.0, 0.0, 0.0, 0.0, 0.0, 0.0], [0.0, 0.0, 0.0, 0.0, 0.0, 0.0], [0.0, 0.0, 0.0, 0.0, 0.0, 0.0], [0.0, 0.0, 0.0, 0.0, 0.0, 0.0], [0.0, 0.0, 0.0, 0.0, 0.0, 0.0]], [[0.0, 0.0, 0.0, 0.0, 0.0, 0.0], [0.0, 0.0, 0.0, 0.0, 0.0, 0.0], [0.0, 0.0, 0.0, 0.0, 0.0, 0.0], [0.0, 0.0, 0.0, 0.0, 0.0, 0.0], [0.0, 0.0, 0.0, 0.0, 0.0, 0.0]]]

Contains list of modified Rydberg-Ritz coefficients for calculating quantum defects for [[ \(S_{1/2},P_{1/2},D_{3/2},F_{5/2}\)], [ \(S_{1/2},P_{3/2},D_{5/2},F_{7/2}\)]].

radialWavefunction(l, s, j, stateEnergy, innerLimit, outerLimit, step)[source]

Radial part of electron wavefunction

Calculates radial function with Numerov (from outside towards the core). Note that wavefunction might not be calculated all the way to the requested innerLimit if the divergence occurs before. In that case third returned argument gives nonzero value, corresponding to the first index in the array for which wavefunction was calculated. For quick example see Rydberg wavefunction calculation snippet.

Parameters:
  • l (int) – orbital angular momentum
  • s (float) – spin angular momentum
  • j (float) – total angular momentum
  • stateEnergy (float) – state energy, relative to ionization threshold, should be given in atomic units (Hatree)
  • innerLimit (float) – inner limit at which wavefunction is requested
  • outerLimit (float) – outer limit at which wavefunction is requested
  • step (flaot) – radial step for integration mesh (a.u.)
Returns:

\(r\)

\(R(r)\cdot r\)

Return type:

List[float], List[flaot], int

Note

Radial wavefunction is not scaled to unity! This normalization condition means that we are using spherical harmonics which are normalized such that \(\int \mathrm{d}\theta~\mathrm{d}\psi~Y(l,m_l)^* \times Y(l',m_{l'}) = \delta (l,l') ~\delta (m_l, m_{l'})\).

Note

Alternative calculation methods can be added here (potenatial package expansion).

scaledRydbergConstant = 0

in eV

updateDipoleMatrixElementsFile()[source]

Updates the file with pre-calculated dipole matrix elements.

This function will add the the file all the elements that have been calculated in the previous run, allowing quick access to them in the future calculations.

arc.alkali_atom_functions.NumerovBack(innerLimit, outerLimit, kfun, step, init1, init2)[source]

Full Python implementation of Numerov integration

Calculates solution function \(rad(r)\) with descrete step in \(r\) size of step, integrating from outerLimit towards the innerLimit (from outside, inwards) equation \(\frac{\mathrm{d}^2 rad(r)}{\mathrm{d} r^2} = kfun(r)\cdot rad(r)\).

Parameters:
  • innerLimit (float) – inner limit of integration
  • outerLimit (flaot) – outer limit of integration
  • kfun (function(double)) – pointer to function used in equation (see longer explanation above)
  • step – descrete step size for integration
  • init1 (float) – initial value, rad`(`outerLimit`+`step)
  • init2 (float) – initial value, rad`(`outerLimit`+:math:`2cdot step)
Returns:

\(r\) (a.u), \(rad(r)\);

Return type:

numpy array of float , numpy array of float, int

Note

Returned function is not normalized!

Note

If AlkaliAtom.cpp_numerov swich is set to True (default option), much faster C implementation of the algorithm will be used instead. That is recommended option. See documentation installation instructions for more details.

arc.alkali_atom_functions.loadSavedCalculation(fileName)[source]

Loads previously saved calculation.

Loads calculations_atom_pairstate.PairStateInteractions and calculations_atom_single.StarkMap calculation instance from file named filename where it was previously saved with saveCalculation .

Example

See example for saveCalculation.

Parameters:fileName – name of the file where calculation will be saved
Returns:saved calculation
arc.alkali_atom_functions.printState(n, l, j)[source]

Prints state spectroscopic label for numeric \(n\), \(l\), \(s\) label of the state

Parameters:
  • n (int) – principal quantum number
  • l (int) – orbital angular momentum
  • j (float) – total angular momentum
arc.alkali_atom_functions.printStateString(n, l, j)[source]

Returns state spectroscopic label for numeric \(n\), \(l\), \(s\) label of the state

Parameters:
  • n (int) – principal quantum number
  • l (int) – orbital angular momentum
  • j (float) – total angular momentum
Returns:

label for the state in standard spectroscopic notation

Return type:

string

arc.alkali_atom_functions.printStateStringLatex(n, l, j)[source]

Returns latex code for spectroscopic label for numeric \(n\), \(l\), \(s\) label of the state

Parameters:
  • n (int) – principal quantum number
  • l (int) – orbital angular momentum
  • j (float) – total angular momentum
Returns:

label for the state in standard spectroscopic notation

Return type:

string

arc.alkali_atom_functions.saveCalculation(calculation, fileName)[source]

Saves calculation for future use.

Saves calculations_atom_pairstate.PairStateInteractions and calculations_atom_single.StarkMap calculations in compact binary format in file named filename. It uses cPickle serialization library in Python, and also zips the final file.

Calculation can be retrieved and used with loadSavedCalculation

Parameters:
  • calculation – class instance of calculations (instance of calculations_atom_pairstate.PairStateInteractions or calculations_atom_single.StarkMap) to be saved.
  • fileName – name of the file where calculation will be saved

Example

Let’s suppose that we did the part of the calculation_atom_pairstate.PairStateInteractions calculation that involves generation of the interaction matrix. After that we can save the full calculation in a single file:

calc = PairStateInteractions(Rubidium(), 60,0,0.5,60,0,0.5, 0.5,0.5)
calc.defineBasis(0,0, 5,5, 25.e9)
calc.diagonalise(np.linspace(0.5,10.0,200),150)
saveCalculation(calc, "mySavedCalculation.pkl")

Then, at a later time, and even on the another machine, we can load that file and continue with calculation. We can for example explore the calculated level diagram:

calc = loadSavedCalculation("mySavedCalculation.pkl")
calc.plotLevelDiagram()
calc.showPlot()
rvdw = calc.getVdwFromLevelDiagram(0.5,14,minStateContribution=0.5,\
                                   showPlot = True)

Or, we can do additional matrix diagonalization, in some new range, then and find C6 by fitting the obtained level diagram:

calc = loadSavedCalculation("mySavedCalculation.pkl")
calc.diagonalise(np.linspace(3,6.0,200),20)
calc.getC6fromLevelDiagram(3,6.0,showPlot=True)

Note that for all loading of saved calculations we’ve been using function loadSavedCalculation .

Note

This doesn’t save results of plotLevelDiagram for the corresponding calculations. Call the plot function before calling showPlot function for the corresponding calculation.