Single atom calculations

Note

Some of the functions (Wavefunction, AtomSurfaceVdW, OpticalLattice1D, DynamicPolarizability, and optical materials properties) in this module are part of new ARC 3.0.0 version. See more at E. J. Robertson, N. Šibalić, R. M. Potvliege and M. P. A. Jones, arXiv:2007.12016 .

Overview

Wavefunction Methods

Wavefunction.getRtimesPsiSpherical(theta, phi, r) Calculates list of \(r \cdot \psi_{m_s} (\theta, \phi, r)\)
Wavefunction.getRtimesPsi(x, y, z) Calculates list of \(r \cdot \psi_{m_s} (x, y, z)\)
Wavefunction.getPsi(x, y, z) Calculates list of \(\psi_{m_s} (x,y,z)\)
Wavefunction.getRtimesPsiSquaredInPlane([…]) Calculates \(|r \cdot \psi|^2\) on a mesh in a given plane.
Wavefunction.plot2D([plane, pointsPerAxis, …]) 2D colour plot of \(|r \cdot \psi|^2\) wavefunction in a requested plane.
Wavefunction.plot3D([plane, pointsPerAxis, …]) 3D colour surface plot of \(|r \cdot \psi|^2\) wavefunction in a requested plane.

StarkMap Methods

StarkMap.defineBasis(n, l, j, mj, nMin, …) Initializes basis of states around state of interest
StarkMap.diagonalise(eFieldList[, …]) Finds atom eigenstates in a given electric field
StarkMap.plotLevelDiagram([units, …]) Makes a plot of a stark map of energy levels
StarkMap.showPlot([interactive]) Shows plot made by plotLevelDiagram
StarkMap.savePlot([filename]) Saves plot made by plotLevelDiagram
StarkMap.exportData(fileBase[, exportFormat]) Exports StarkMap calculation data.
StarkMap.getPolarizability([maxField, …]) Returns the polarizability of the state (set during the initalization process).
StarkMap.getState(state, electricField, …) Returns basis states and coefficients that make up for a given electric field the eigenstate with largest contribution of the original state.

LevelPlot Methods

LevelPlot is also called Grotrian diagram, or term diagram.

LevelPlot.makeLevels(nFrom, nTo, lFrom, lTo) Constructs energy level diagram in a given range
LevelPlot.drawLevels() Draws a level diagram plot
LevelPlot.showPlot() Shows a level diagram plot

AtomSurfaceVdW Methods

AtomSurfaceVdW.getC3contribution(n1, l1, j1, …) Contribution to \(C_3\) of \(|n_1, \ell_1, j_1\rangle\) state due to dipole coupling to \(|n_2, \ell_2, j_2\rangle\) state.
AtomSurfaceVdW.getStateC3(n, l, j, …[, s, …]) Van der Waals atom-surface interaction coefficient for a given state (\(C_3\) in units of \(\mathrm{J}\cdot\mathrm{m}^3\) )

OpticalLattice1D Methods

OpticalLattice1D.getRecoilEnergy() Recoil energy for atoms in given optical lattice
OpticalLattice1D.getTrappingFrequency(…) Atom’s trapping frequecy for given trapth depth
OpticalLattice1D.defineBasis([lLimit]) Define basis for Bloch band calculations
OpticalLattice1D.diagonalise(…[, …]) Calculates energy levels (Bloch bands) for given quasimomentumList
OpticalLattice1D.plotLevelDiagram() Plots energy level diagram (Bloch bands).
OpticalLattice1D.BlochWavefunction(…) Bloch wavefunction as a function of 1D coordinate.
OpticalLattice1D.getWannierFunction(x[, …]) Gives value at cooridnate x of a Wannier function localized at given lattice index.

DynamicPolarizability Methods

DynamicPolarizability.defineBasis(nMin, nMax) Defines basis for calculation of dynamic polarizability
DynamicPolarizability.getPolarizability(…) Calculates of scalar, vector, tensor, core and pondermotive
DynamicPolarizability.plotPolarizability(…) Plots of polarisability for a range of wavelengths.

Optical material properties

OpticalMaterial() Abstract class implementing calculation of basic properties for optical materials.
Air() Air as an optical material at normal conditions
Sapphire() Sapphire as optical material.

Detailed documentation

This module provides calculations of single-atom properties.

Included calculations are Stark maps, level plot visualisations, lifetimes and radiative decays.

class arc.calculations_atom_single.AtomSurfaceVdW(atom, surfaceMaterial=None)[source]

Calculates atom-surface Van der Waals interaction.

Energy of atom state \(|i\rangle\) at distance \(z\) from the surface of material is offseted in energy by \(V_{\rm VdW}\) at small distances \(z\ll\rm{min}(\lambda_{i,j})\) , where \(\lambda_{i,j}\) are the wavelengths from atom state \(|i \rangle\) to all strongly-coupled states \(j\) , due to (unretarded) atom-surface interaction, also called Van der Waals interaction. The interaction potential can be expressed as

\(V_{\rm VdW} = - \frac{C_3}{z^3}\)

This class calculates \(C_3\) for individual states \(|i\rangle\).

See example atom-surface calculation snippet.

Parameters:
  • atom (AlkaliAtom or DivalentAtom) – specified Alkali or Alkaline Earth atom whose interaction with surface we want to explore
  • material (from arc.materials) – specified surface material

Note

To find frequecy shift of a transition \(|\rm a \rangle\rightarrow |\rm b \rangle\), one needs to calculate difference in \(C_3\) coefficients obtained for the two states \(|\rm a\rangle\) and \(|\rm b\rangle\) respectively. See example TODO (TO-DO)

getC3contribution(n1, l1, j1, n2, l2, j2, s=0.5)[source]

Contribution to \(C_3\) of \(|n_1, \ell_1, j_1\rangle\) state due to dipole coupling to \(|n_2, \ell_2, j_2\rangle\) state.

Calculates \(\frac{1}{4\pi\varepsilon_0}\ \frac{ n(\omega_{\rm ab})^2 - 1}{ n(\omega_{\rm ab})^2 + 1}\ \frac{ \left| \langle a| D_x | b \rangle \right|^2 \ + \left| \langle a | D_y | b \rangle \right|^2 + \ 2 \cdot \left|\langle a |D_z| b \rangle \right|^2}{16}\)

where \(|{\rm a}\rangle \equiv |n_1, \ell_1, j_1\rangle\) , \(|{\rm b}\rangle \equiv |n_2, \ell_2, j_2\rangle\), \(\mathbf{D} \equiv e \cdot \mathbf{r} \ \equiv \hat{x} D_x + \hat{y} D_y\ + \hat{z} D_z\) is atomic dipole operator and \(n(\omega_{\rm ab})\) is refractive index of the considered surface at transition frequency \(\omega_{\rm ab}\) .

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 od state 2
  • l2 (int) – orbital angular momentum of state 2
  • j2 (float) – total angular momentum of state 2
  • s (float) – optional, spin angular momentum of states. Default value of 0.5 is correct for AlkaliAtoms. For DivalentAtom it has to be explicitly stated
Returns:

contribution to VdW coefficient \(C_3\) ,estimated error \(\delta C_3\) (in units of \({\rm J}\cdot{\rm m}^3\)), and refractive index \(n\) of the surface material for the given transition.

Return type:

float, float, float

Warning

This is just contribution of one transition to the level shift of a particular state. To calculate total level shift, check AtomSurfaceVdW.getStateC3

getStateC3(n, l, j, coupledStatesList, s=0.5, debugOutput=False)[source]

Van der Waals atom-surface interaction coefficient for a given state (\(C_3\) in units of \(\mathrm{J}\cdot\mathrm{m}^3\) )

Parameters:
  • n (int) – principal quantum number of the state
  • l (int) – orbital angular momentum of the state
  • j (int) – total angular momentum of state
  • coupledStatesList (array) – array of states that are strongly dipole-coupled to the initial state, whose contribution to \(C_3\) will be take into account. Format [[n1,l1,j1],…]
  • s (float, optional) – total spin angular momentum for the considered state. Default value of 0.5 is correct for AlkaliAtoms, but it has to be explicitly specifiied for DivalentAtom.
  • debugOutput (bool, optional) – prints additional output information, False by default.
Returns:

\(C_3\) (in units of \({\rm J}\cdot {\rm m}^3\) ), estimated error \(\delta C_3\)

Return type:

float, float

class arc.calculations_atom_single.DynamicPolarizability(atom, n, l, j, s=0.5)[source]

Calculations of magic wavelengths and dynamic polarizability (scalar and tensor).

Parameters:
  • atom – alkali or alkaline element of choice
  • n (int) – principal quantum number of the selected stated
  • l (int) – orbital angular momentum of the selected state
  • j (float) – total angular momentum of selected state
  • s (float) – optional, spin state of the atom. Default value of 0.5 is correct for Alkali atoms, but it has to be explicitly specified for DivalentAtom.
defineBasis(nMin, nMax)[source]

Defines basis for calculation of dynamic polarizability

Parameters:
  • nMin (int) – minimal principal quantum number of states to be taken into account for calculation
  • nMax (int) – maxi,al principal quantum number of states to be taken into account for calculation
getPolarizability(driveWavelength, units='SI', accountForStateLifetime=False, mj=None)[source]
Calculates of scalar, vector, tensor, core and pondermotive polarizability, and returns state corresponding to the closest transition resonance.

Note that pondermotive polarisability is calculated as \(\alpha_P = e^2 / (2 m_e \omega^2)\), i.e. assumes that the definition of the energy shift in field \(E\) is \(\frac{1}{2}\alpha_P E^2\). For more datils check the preprint arXiv:2007.12016 that introduced the update.

Args:
driveWavelength (float): wavelength of driving field
(in units of m)
units (string): optional, ‘SI’ or ‘a.u.’ (equivalently ‘au’),
switches between SI units for returned result (\(Hz V^-2 m^2\) ) and atomic units (“\(a_0^3\) “). Defaul ‘SI’
accountForStateLifetime (bool): optional, should we account
for finite transition linewidths caused by finite state lifetimes. By default False.
Returns:
scalar, vector, tensor, pondermotive polarisability of state, core polarisability and atomic state whose resonance is closest in energy. Returned units depend on units parameter (default SI).
plotPolarizability(wavelengthList, mj=None, addToPlotAxis=None, line='b-', units='SI', addCorePolarisability=True, addPondermotivePolarisability=False, accountForStateLifetime=False, debugOutput=False)[source]

Plots of polarisability for a range of wavelengths.

Can be combined for different states to allow finding magic wavelengths for pairs of states. Currently supports only driving with linearly polarised light. See example magic wavelength snippet.

Parameters:
  • wavelengthList (array) – wavelengths for which we want to calculate polarisability (in units of m).
  • mj (float) – optional, mj projection of the total angular momenutum for the states for which we are calculating polarisability. By default it’s +j.
  • line (string) – optional, line style short definition to be passed to matplotlib when plotting calculated polarisabilities
  • units (string) – optional, ‘SI’ or ‘a.u.’ (equivalently ‘au’), switches between SI units for returned result (\(Hz V^-2 m^2\) ) and atomic units (“\(a_0^3\) “). Deafault ‘SI’.
  • addCorePolarisability (bool) – optional, should ionic core polarisability be taken into account. By default True.
  • addPondermotivePolarisability (bool) – optional, should pondermotive polarisability (also called free-electron polarisability) be added to the total polarisability. Default is False. It assumes that there is no significant variation of trapping field intensity over the range of the electric cloud. If this condition is not satisfied, one has to calculate total shift as average over the electron wavefunction.
  • accountForStateLifetime (bool) – optional, should we account for finite transition linewidths caused by finite state lifetimes. By default False.
  • debugOutput (bool) – optonal. Print additional output on resonances Default value False.
class arc.calculations_atom_single.LevelPlot(atomType)[source]

Single atom level plots and decays (a Grotrian diagram, or term diagram)

For an example see Rydberg energy levels example snippet.

Parameters:atom (AlkaliAtom or DivalentAtom) – ={ arc.alkali_atom_data.Lithium6, arc.alkali_atom_data.Lithium7, arc.alkali_atom_data.Sodium, arc.alkali_atom_data.Potassium39, arc.alkali_atom_data.Potassium40, arc.alkali_atom_data.Potassium41, arc.alkali_atom_data.Rubidium85, arc.alkali_atom_data.Rubidium87, arc.alkali_atom_data.Caesium, arc.divalent_atom_data.Strontium88, arc.divalent_atom_data.Calcium40 arc.divalent_atom_data.Ytterbium174 } Alkali atom type whose levels we want to examine
drawLevels()[source]

Draws a level diagram plot

makeLevels(nFrom, nTo, lFrom, lTo, sList=[0.5])[source]

Constructs energy level diagram in a given range

Parameters:
  • nFrom (int) – minimal principal quantum number of the states we are interested in
  • nTo (int) – maximal principal quantum number of the states we are interested in
  • lFrom (int) – minimal orbital angular momentum of the states we are interested in
  • lTo (int) – maximal orbital angular momentum of the states we are interested in
  • sList (float) – optional, spin angular momentum. Default value of [0.5] corresponds to Alkali atoms. For Alkaline Earth it has to be specified. For divalent atoms one can plot either one spin state by setting for example sList=[0]`, or both spin states sList=[0,1]`
showPlot()[source]

Shows a level diagram plot

class arc.calculations_atom_single.OpticalLattice1D(atom, trapWavenegth)[source]

Atom properties in optical lattices in 1D.

See example optical lattice calculations snippet.

Parameters:
  • atom – one of AlkaliAtom or DivalentAtom
  • trapWavenegth (float) – wavelength of trapping laser light (in units of m)
BlochWavefunction(trapPotentialDepth, quasimomentum, blochBandIndex)[source]

Bloch wavefunction as a function of 1D coordinate.

Example

Returns Bloch wavefunction. Use as following:

trapPotentialDepth = 40  # units of recoil energy
quasimomentum = 0
blochBandIndex = 0  # Bloch band lowest in energy is 0
wf = lattice.BlochWavefunction(trapPotentialDepth,
                               quasimomentum,
                               blochBandIndex)
wf(x)  # returns complex number corresponding to value of Bloch
       # wavefunction at point x (cooridnate given in units of
       # 1/k where k = 2 \pi / trapWavenegth )
       # by default k=1, so one full wavelength is 2\pi
Parameters:
  • trapPotentialDepth (float) – (in units of recoil energy OpticalLattice1D.getRecoilEnergy)
  • quasimomentum (float) – (in units of 2 pi / OpticalLattice1D.trapWavenegth; note that reciprocal lattice momentum in this units is 2, and that full range of quasimomentum is from -1 to +1)
Returns:

Bloch wavefunction as a function of coordinate (see call example above)

defineBasis(lLimit=35)[source]

Define basis for Bloch band calculations

Bloch states are calculated suming up all relevant states with momenta in range [-lLimit * 4 * pi /trapWavenegth, +lLimit * 4 * pi /trapWavenegth] Note that factor of 4 occurs since potential lattice period is twice the trapWavelength for standing wave.

Parameters:lLimit (integer) – Optional, defines maximal momentum to be taken for calculation of Bloch States as lLimit * 4 * pi / trapWavenegth . By default set to 35.
diagonalise(trapPotentialDepth, quasimomentumList, saveBandIndex=None)[source]

Calculates energy levels (Bloch bands) for given quasimomentumList

Energy levels and their quasimomentum are saved in internal variables energy and quasimomentum. Energies are saved in units of recoil energy, and quasimomentum in units of The optional parameter saveBandIndex specifies index of the Bloch band for which eigenvectrors should be saved. If provided, eigenvectors for each value quasimomentumList[i] are saved in savedBlochBand[i].

Parameters:
  • latticePotential (float) – lattice depth formed by the standing wave of laser, with wavelength specified during initialisation of the lattice (in units of recoil energy).
  • quasimomentumList (array) – array of quasimomentum values for which energy levels will be calculated (in units of \(\hbar \cdot k\), where \(k\) is trapping laser wavevector; since reciprocal lattice has twice the trapping laser wavevector due to standing wave, full range of quasimomentum is from -1 to +1)
  • saveBandIndex (int) – optional, default None. If provided, specifies for which Bloch band should the eignevectors be also saved. saveBlochBand=0 corresponds to lowest energy band.
energy = []

energy of states obtained by OpticalLattice1D.diagonalise method in format [[energies for quasimomentum1 ], [energies for quasimomentum2 ], …]

getRecoilEnergy()[source]

Recoil energy for atoms in given optical lattice

Returns:recoil energy in units of J
Return type:float
getTrappingFrequency(trapPotentialDepth)[source]

Atom’s trapping frequecy for given trapth depth

Parameters:trapPotentialDepth (float) – lattice depth (in units of J)
Returns:trapping frequency (in Hz)
Return type:float
getWannierFunction(x, latticeIndex=0, k=1)[source]

Gives value at cooridnate x of a Wannier function localized at given lattice index.

Parameters:
  • x (float) – spatial coordinate (in units of \(2\pi/k\) ; for default value of laser drivng wavevecto \(k=1\) , one trappinWavelength is \(2\pi\) ). Coordinate origin is at latticeIndex=0 .
  • latticeIndex (int) – optional, lattice index at which the Wannier function is localised. By defualt 0.
  • k (float) – optional; laser driving wavevector, defines unit of length. Default value is 1, making one trapping laser wavelenth equal to \(2\pi\)
plotLevelDiagram()[source]

Plots energy level diagram (Bloch bands).

Based on diagonalisation of the lattice potential, plots descrete eigen energy spectra obtained for each value of the quasimomentum used in OpticalLattice1D.diagonalise method.

Returns:matploltib figure with a Bloch bands
quasimomentum = []

list of quzimomentum for which the energies of states was calculated by OpticalLattice1D.diagonalise method in format [quasimomentum1, quasimomentum2, …]

savedBlochBand = []

list of saved eigen energy state compositions for each of the Calculated quasimomentums for the selected index of the Bloch band in OpticalLattice1D.diagonalise method in format [[eigen state decomposition for quasimomentum 1], [eigen state decomposition for quasimomentum 2], …]

trapPotentialDepth = 0

save slattice trap potential depth for which calculation OpticalLattice1D.diagonalise was done

class arc.calculations_atom_single.StarkMap(atom)[source]

Calculates Stark maps for single atom in a field

This initializes calculation for the atom of a given type. For details of calculation see Zimmerman [1]. For a quick working example see Stark map example snippet.

Parameters:atom (AlkaliAtom or DivalentAtom) – ={ arc.alkali_atom_data.Lithium6, arc.alkali_atom_data.Lithium7, arc.alkali_atom_data.Sodium, arc.alkali_atom_data.Potassium39, arc.alkali_atom_data.Potassium40, arc.alkali_atom_data.Potassium41, arc.alkali_atom_data.Rubidium85, arc.alkali_atom_data.Rubidium87, arc.alkali_atom_data.Caesium, arc.divalent_atom_data.Strontium88, arc.divalent_atom_data.Calcium40 arc.divalent_atom_data.Ytterbium174 } Select the alkali metal for energy level diagram calculation

Examples

State \(28~S_{1/2}~|m_j|=0.5\) polarizability calculation

>>> from arc import *
>>> calc = StarkMap(Caesium())
>>> calc.defineBasis(28, 0, 0.5, 0.5, 23, 32, 20)
>>> calc.diagonalise(np.linspace(00.,6000,600))
>>> print("%.5f MHz cm^2 / V^2 " % calc.getPolarizability())
0.76705 MHz cm^2 / V^2

Stark map calculation

>>> from arc import *
>>> calc = StarkMap(Caesium())
>>> calc.defineBasis(28, 0, 0.5, 0.5, 23, 32, 20)
>>> calc.diagonalise(np.linspace(00.,60000,600))
>>> calc.plotLevelDiagram()
>>> calc.showPlot()
<< matplotlib plot will open containing a Stark map >>

Examples

Advanced interfacing of Stark map calculations (StarkMap class) Here we show one easy way to obtain the Stark matrix (from diagonal mat1 and off-diagonal part mat2 ) and basis states (stored in basisStates ), if this middle-product of the calculation is needed for some code build on top of the existing ARC package.

>>> from arc import *
>>> calc = StarkMap(Caesium())
>>> calc.defineBasis(28, 0, 0.5, 0.5, 23, 32, 20)
>>> # Now we have matrix and basis states, that we can used in our own code
>>> # Let's say we want Stark map at electric field of 0.2 V/m
>>> eField = 0.2 # V/m
>>> # We can easily extract Stark matrix
>>> # as diagonal matrix (state detunings)
>>> #  + off-diagonal matrix (propotional to electric field)
>>> matrix = calc.mat1+calc.mat2*eField
>>> # and the basis states as array [ [n,l,j,mj] , ...]
>>> basisStates = calc.basisStates
>>> # you can do your own calculation now...

References

[1]M. L. Zimmerman et.al, PRA 20:2251 (1979) https://doi.org/10.1103/PhysRevA.20.2251
ax = None

pointer towards matplotlib figure axis after plotLevelDiagram is called to create figure

basisStates = None

List of basis states for calculation in the form [ [n,l,j,mj], …]. Calculated by defineBasis .

defineBasis(n, l, j, mj, nMin, nMax, maxL, Bz=0, progressOutput=False, debugOutput=False, s=0.5)[source]

Initializes basis of states around state of interest

Defines basis of states for further calculation. \(n,l,j,m_j\) specify state whose neighbourhood and polarizability we want to explore. Other parameters specify basis of calculations. This method stores basis in basisStates, while corresponding interaction matrix is stored in two parts. First part is diagonal electric-field independent part stored in mat1, while the second part mat2 corresponds to off-diagonal elements that are propotional to electric field. Overall interaction matrix for electric field eField can be then obtained as fullStarkMatrix = mat1 + mat2 *eField

Parameters:
  • n (int) – principal quantum number of the state
  • l (int) – angular orbital momentum of the state
  • j (flaot) – total angular momentum of the state
  • mj (float) – projection of total angular momentum of the state
  • nMin (int) – minimal principal quantum number of the states to be included in the basis for calculation
  • nMax (int) – maximal principal quantum number of the states to be included in the basis for calculation
  • maxL (int) – maximal value of orbital angular momentum for the states to be included in the basis for calculation
  • Bz (float) – optional, magnetic field directed along z-axis in units of Tesla. Calculation will be correct only for weak magnetic fields, where paramagnetic term is much stronger then diamagnetic term. Diamagnetic term is neglected.
  • progressOutput (bool, optional) – if True prints the progress of calculation; Set to false by default.
  • debugOutput (bool, optional) – if True prints additional information usefull for debuging. Set to false by default.
  • s (float) – optional. Total spin angular momentum for the state. Default value of 0.5 is correct for Alkaline Atoms, but value has to be specified explicitly for divalent atoms (e.g. s=0 or s=1 for singlet and triplet states, that have total spin angular momenutum equal to 0 or 1 respectively).
diagonalise(eFieldList, drivingFromState=[0, 0, 0, 0, 0], progressOutput=False, debugOutput=False)[source]

Finds atom eigenstates in a given electric field

Eigenstates are calculated for a list of given electric fields. To extract polarizability of the originaly stated state see getPolarizability method. Results are saved in eFieldList, y and highlight.

Parameters:
  • eFieldList (array) – array of electric field strength (in V/m) for which we want to know energy eigenstates
  • progressOutput (bool, optional) – if True prints the progress of calculation; Set to false by default.
  • debugOutput (bool, optional) – if True prints additional information usefull for debuging. Set to false by default.
eFieldList = None

Saves electric field (in units of V/m) for which energy levels are calculated

See also

y, highlight, diagonalise

exportData(fileBase, exportFormat='csv')[source]

Exports StarkMap calculation data.

Only supported format (selected by default) is .csv in a human-readable form with a header that saves details of calculation. Function saves three files: 1) filebase _eField.csv; 2) filebase _energyLevels 3) filebase _highlight

For more details on the format, see header of the saved files.

Parameters:
  • filebase (string) – filebase for the names of the saved files without format extension. Add as a prefix a directory path if necessary (e.g. saving outside the current working directory)
  • exportFormat (string) – optional. Format of the exported file. Currently only .csv is supported but this can be extended in the future.
fig = None

pointer towards matplotlib figure after plotLevelDiagram is called to create figure

getPolarizability(maxField=10000000000.0, showPlot=False, debugOutput=False, minStateContribution=0.0)[source]

Returns the polarizability of the state (set during the initalization process).

Fits offset of the energy level of the state to \(\frac{1}{2} \alpha_0 E^2\), where \(E\) is the applied static electric field, and returns fitted value \(\alpha_0\)

Parameters:
  • maxField (float, optional) – maximum field (in V/m) to be used for fitting the polarizability. By default, max field is very large, so it will use eigenvalues calculated in the whole range.
  • showPlot (bool, optional) – shows plot of calculated eigenValues of the given state (dots), and the fit (solid line) for extracting polarizability
  • debugOutput (bool, optional) – if True prints additional information usefull for debuging. Set to false by default.
Returns:

scalar polarizability in units of MHz cm \(^2\) / V \(^2\)

Return type:

float

getState(state, electricField, minN, maxN, maxL, accountForAmplitude=0.95, debugOutput=False)[source]

Returns basis states and coefficients that make up for a given electric field the eigenstate with largest contribution of the original state.

Parameters:
  • state (array) – target basis state in format \([n,\ell,j,m_j]\) corresponding to the state whose composition we want to track as we apply the electric field
  • electricField (float) – applied DC electric field in units of V/m.
  • minN (int) – minimal principal quantum number to be taken for calculation of the Stark mixing
  • maxN (int) – maximal principal quantum nunber to be take for calculation of the Start mixing
  • maxL (int) – maximal orbital angular momentum of states that should be taken in calculation of the Stark mixing
  • accountForAmplitude (float) – optinal, relative amplitude of state that should be reached with the subset of the eigen states returned. The returned eigen states will be sorted in the declining relative contribution to the final eigen state, and once total accounted amplitude of the state reaches 0.95, further output of additional small contribution of the other basis states to the final states will be supressed. Default value of 0.95 will force output until basis state accounts for 95% of the state amplitude.
  • debugOutput (bool) – optional, prints additional debug information if True. Default False.
Returns:

array of states in format [[n1, l1, j1, mj1], …] and array of complex coefficients in format [c1, c2, …] corresponding the projections of the eigenstate (thas has largest contribution of the original state in the given electric field) on the basis states, and energy of the found state in (eV)

highlight = None

highlight[i] is an array of values measuring highlighted feature in the eigenstates at electric field intensity eFieldList[i]. E.g. highlight[i][j] measures highlighted feature of the state with energy y[i][j] at electric field eFieldList[i]. What will be highlighted feature is defined in the call of diagonalise (see that part of documentation for details).

See also

eFieldList, y, diagonalise

indexOfCoupledState = None

Index of coupled state (initial state passed to defineBasis) in basisStates list of basis states

mat1 = None

diagonal elements of Stark-matrix (detuning of states) calculated by defineBasis in the basis basisStates.

mat2 = None

off-diagonal elements of Stark-matrix divided by electric field value. To get off diagonal elemements multiply this matrix with electric field value. Full Stark matrix is obtained as fullStarkMatrix = mat1 + mat2 *eField. Calculated by defineBasis in the basis basisStates.

plotLevelDiagram(units=1, highlighState=True, progressOutput=False, debugOutput=False, highlightColour='red', addToExistingPlot=False)[source]

Makes a plot of a stark map of energy levels

To save this plot, see savePlot. To print this plot see showPlot. Pointers (handles) towards matplotlib figure and axis used are saved in fig and ax variables respectively.

Parameters:
  • units (int,optional) – possible values {1,2} ; if the value is 1 (default) Stark diagram will be plotted in energy units cm \({}^{-1}\); if value is 2, Stark diagram will be plotted as energy \(/h\) in units of GHz
  • highlightState (bool, optional) – False by default. If True, scatter plot colour map will map in red amount of original state for the given eigenState
  • progressOutput (bool, optional) – if True prints the progress of calculation; Set to False by default.
  • debugOutput (bool, optional) – if True prints additional information usefull for debuging. Set to False by default.
  • addToExistingPlot (bool, optional) – if True adds points to existing old plot. Note that then interactive plotting doesn’t work. False by default.
s = None

spin manifold in which we are working default value of 0.5 is correct for Alkaline Atoms. Otherwise it has to be specified when calling defineBasis as s=0 or s=1 for singlet and triplet states respectively

savePlot(filename='StarkMap.pdf')[source]

Saves plot made by plotLevelDiagram

Parameters:filename (str, optional) – file location where the plot should be saved
showPlot(interactive=True)[source]

Shows plot made by plotLevelDiagram

y = None

y[i] is an array of eigenValues corresponding to the energies of the atom states at the electric field eFieldList[i]. For example y[i][j] is energy of the j eigenvalue (energy of the state) measured in cm \({}^{-1}\) relative to the ionization threshold.

class arc.calculations_atom_single.Wavefunction(atom, basisStates, coefficients)[source]
Calculates and plots electron wavefunctions.

For an example see wavefunction plotting example snippet.

Args:

atom: atom type considered (for example Rubidum87()) basisStates (array): array of states in fine basis that contribute

to the state whose wavefunction is requested. \([[n_1, \ell_1, j_1, m_{j1}], ...]\) For efficient calculation do not pass all the possible basis states, but just the once that have significant contribution to the reqested state.
coefficients (array): array [c1, …] of complex coefficents
\(c_i = \langle \psi_i |\psi\rangle\) corresponding to decomposition of required state \(|\psi\rangle\) on basis states \(|\psi_i \rangle\) .
getPsi(x, y, z)[source]

Calculates list of \(\psi_{m_s} (x,y,z)\)

At point define by Cartesian coordinates returns list of \(\psi_{m_s} (x,y,z)\) wavefunction values corresponding to different electron spin projection values \(m_s\).

Parameters:
  • x (float) – Cartesian coordinates of selected point, relative to the atom core. (in atomic units of Bohr radius \(a_0\))
  • y (float) – Cartesian coordinates of selected point, relative to the atom core. (in atomic units of Bohr radius \(a_0\))
  • z (float) – Cartesian coordinates of selected point, relative to the atom core. (in atomic units of Bohr radius \(a_0\))
Returns:

list of complex values corresponding to \(\psi_{m_s} (\theta, \phi, r)\) for different spin states \(m_s\) contributing to the state in decreasing order of \(m_s\). For example, for arc.AlkaliAtom returns \(\psi_{m_s=+1/2} (\theta, \phi, r)\) and \(\psi_{m_s=-1/2} (\theta, \phi, r)\) . )`.

getRtimesPsi(x, y, z)[source]

Calculates list of \(r \cdot \psi_{m_s} (x, y, z)\)

At a point defined by Cartesian coordinates returns list of \(r \cdot \psi_{m_s} (x, y, z)\) wavefunction values for different electron spin projection values \(m_s\).

Parameters:
  • x (float) – Cartesian coordinates of selected point, relative to the atom core. (in atomic units of Bohr radius \(a_0\))
  • y (float) – Cartesian coordinates of selected point, relative to the atom core. (in atomic units of Bohr radius \(a_0\))
  • z (float) – Cartesian coordinates of selected point, relative to the atom core. (in atomic units of Bohr radius \(a_0\))
Returns:

list of complex values corresponding to \(r \cdot \psi_{m_s} (\theta, \phi, r)\) for different spin states \(m_s\) contributing to the state in decreasing order of \(m_s\). For example, for arc.AlkaliAtom returns \(r \cdot \psi_{m_s=+1/2} (\theta, \phi, r)\) and \(r \cdot \psi_{m_s=-1/2} (\theta, \phi, r)\) . )`, where \(r=\sqrt{x^2+y^2+z^2}\).

getRtimesPsiSpherical(theta, phi, r)[source]

Calculates list of \(r \cdot \psi_{m_s} (\theta, \phi, r)\)

At point defined by spherical coordinates, returns list of \(r \cdot \psi_{m_s} (\theta, \phi, r)\) wavefunction values for different electron spin projection values \(m_s\).

Coordinates are defined relative to atomic core.

Parameters:
  • theta (float) – polar angle (angle between \(z\) axis and vector pointing towards selected point) (in units of radians).
  • phi (float) – azimuthal angle (angle between \(x\) axis and projection at \(x-y\) plane of vector pointing towards selected point) (in units of radians).
  • r (float) – distance between coordinate origin and selected point. (in atomic units of Bohr radius \(a_0\))
Returns:

list of complex values corresponding to \(\psi_{m_s} (\theta, \phi, r)\) for different spin states \(m_s\) contributing to the state in decreasing order of \(m_s\). For example, for arc.AlkaliAtom returns \(r \cdot \psi_{m_s=+1/2} (\theta, \phi, r)\) and \(r \cdot \psi_{m_s=-1/2} (\theta, \phi, r) \)

getRtimesPsiSquaredInPlane(plane='x-z', pointsPerAxis=150, axisLength=None, units='atomic')[source]

Calculates \(|r \cdot \psi|^2\) on a mesh in a given plane.

Parameters:
  • plane (str) – optiona, set’s calculation plane to ‘x-y’ or ‘x-z’. Default value ‘x-y’
  • pointsPerAxis (int) – optional, a number of mesh points per Carthesian axis. Default value of 150, gives a mesh with total size of \(150 \times 150 = 22500\) points.
  • axisLength (float) – optional, length of the square in the selected plane on which wavefunction will be calculated. By default it is largw enough to fit the whole wavefunction (in atomic units of Bohr radius \(a_0\)).
  • units (str) – optional, units of length in which calculated mesh will be returned (note that axisLength is on the other hand always in atomi units.). Supported values are ‘atomic’ or ‘nm’. Default value ‘atomic’ .
Returns:

meshCoordinate1, meshCoordinate2 and \(|r \cdot \psi|^2 = \sum_{m_s} |r \cdot \psi_{m_s}|^2\), where sum is over possible electron spin projection values \(m_s\).

plot2D(plane='x-z', pointsPerAxis=150, axisLength=None, units='atomic', colorbar=True, labels=True)[source]

2D colour plot of \(|r \cdot \psi|^2\) wavefunction in a requested plane.

Parameters:
  • plane (str) – optiona, set’s calculation plane to ‘x-y’ or ‘x-z’. Default value ‘x-y’
  • pointsPerAxis (int) – optional, a number of mesh points per Carthesian axis. Default value of 150, gives a mesh with total size of \(150 \times 150 = 22500\) points.
  • axisLength (float) – optional, length of the square in the selected plane on which wavefunction will be calculated. By default it is large enough to fit the whole wavefunction (in atomic units of Bohr radius \(a_0\)).
  • units (str) – optional, units of length in which calculated mesh will be returned (note that axisLength is on the other hand always in atomi units.). Supported values are ‘atomic’ or ‘nm’. Default value ‘atomic’ .
  • colorbar (bool) – optional, determens if the colour bar scale of should be shown. Default value is True.
  • labels (bool) – optional, determines if the labels on the axis of the plot should be shown. Default value is True.
Returns:

matplotlib.pyplot.figure object with a requested plot. Use show() method to see figure.

plot3D(plane='x-z', pointsPerAxis=150, axisLength=None, units='atomic', labels=True)[source]

3D colour surface plot of \(|r \cdot \psi|^2\) wavefunction in a requested plane.

Parameters:
  • plane (str) – optiona, set’s calculation plane to ‘x-y’ or ‘x-z’. Default value ‘x-y’
  • pointsPerAxis (int) – optional, a number of mesh points per Carthesian axis. Default value of 150, gives a mesh with total size of \(150 \times 150 = 22500\) points.
  • axisLength (float) – optional, length of the square in the selected plane on which wavefunction will be calculated. By default it is large enough to fit the whole wavefunction (in atomic units of Bohr radius \(a_0\)).
  • units (str) – optional, units of length in which calculated mesh will be returned (note that axisLength is on the other hand always in atomi units.). Supported values are ‘atomic’ or ‘nm’. Default value ‘atomic’ .
  • labels (bool) – optional, determines if the labels on the axis of the plot should be shown. Default value is True.
Returns:

matplotlib.figure object with a requested plot. Use show() method to see figure.

class arc.materials.Air[source]

Bases: arc.materials.OpticalMaterial

Air as an optical material at normal conditions

getN(vacuumWavelength=None, *args, **kwargs)[source]

Assumes temperature: 15 °C, pressure: 101325 Pa

class arc.materials.OpticalMaterial[source]

Bases: object

Abstract class implementing calculation of basic properties for optical materials.

getN(*args, **kwargs)[source]

Refractive index of material

name = ''

Human-friendly name of material

sources = []

List of .csv files listing refractive index measurements first column in these files is wavelength (in mu m), the second refractive index

sourcesComment = []

Any notes about measured values

sourcesRange = []

Array of max and minimal wavelegth pairs [lambdaMin, lambdaMax] for each of the sources. Automatically loaded from sources list

class arc.materials.Sapphire[source]

Bases: arc.materials.OpticalMaterial

Sapphire as optical material.

getN(vacuumWavelength=None, airWavelength=None, axis='ordinary', *args, **kwargs)[source]