API / Code Reference

Python package for simulating light-pulse atom interferometers.

Submodules

aisim.atoms module

Classes and functions related to the atomic cloud.

class aisim.atoms.AtomicEnsemble(phase_space_vectors, state_kets=[1, 0], time=0, weights=None)

Bases: object

Represents an atomic ensemble consisting of n atoms.

Each atom is is defined by its phase space vector (x0, y0, z0, vx, vy, vz) at time t=0. From this phase space vector the position at later times can be calculated. Optionally, weights can be added for each atom in the ensemble. Slicing along the axis of the n atoms is supported.

Parameters:
  • phase_space_vectors (ndarray) – n × 6 dimensional array representing the phase space vectors (x0, y0, z0, vx, vy, vz) of the atoms in an atomic ensemble
  • state_kets (m × 1 or n × m x 1 array or list, optional) –
    vector(s) representing the m internal degrees of freedom of the
    atoms. If the list or array is one-dimensional, all atoms are initialized in the same internal state. Alternatively, each atom can be initialized with a different state vector by passing an array of state vectors for every atom. E.g. to initialize all atoms in the ground state of a two-level system, pass [1, 0] which is the default.
  • time (float, optional) – the initial time (default 0) when the phase space and state vectors are initialized
  • weights (1darray , optional) – Optional weights for each of the n atoms in the ensemble
phase_space_vectors

n × 6 dimensional array representing the phase space vectors (x0, y0, z0, vx, vy, vz) of the atoms in an atomic ensemble

Type:ndarray
position
velocity
state_kets
state_bras
density_matrices
density_matrix
calc_position(t)

Calculate the positions (x, y, z) of the atoms after propagation.

Parameters:t (float) – time when the positions should be calculated
Returns:pos – n × 3 dimensional array of the positions (x, y, z)
Return type:array
density_matrices

Density matrix of the m level system of the n atoms.

These are pure states.

Type:(n × m x m) array
density_matrix

Density matrix of the ensemble’s m level system.

Type:(m x m) array
fidelity(rho_target)

Calculate fidelity of ensemble’s density matrix and target matrix [1].

Parameters:rho_target (array) – target density matrix as m x m array
Returns:fidelity – fidelity of AtomicEnsemble’s compared to target density matrix
Return type:float

References

[1] https://en.wikipedia.org/wiki/Fidelity_of_quantum_states

position

Positions (x, y, z) of the atoms in the ensemble.

Type:(n × 3) array
state_bras

The bra vectors of the m level system.

Type:(n × 1 x m) array
state_kets

The ket vectors of the m level system.

Type:(n × m x 1) array
state_occupation(state)

Calculate the state population of each atom in the ensemble.

Parameters:state (int or list_like) – Specifies which state population should be calculated. E.g. the excited state of a two-level system can be calculated by passing either 1 or [0, 1].
Returns:occupation – n dimensional array of the state population of each of the n atom
Return type:array
time

Time changes when propagating the atomic ensemble.

Type:float
velocity

Velocities of the atoms in the ensemble.

Type:array
aisim.atoms.combine_grids(pos, vel)

Combine a position and velocity grid into an array of phase space vectors.

The resulting array contains (x, y, z, vx, vy, vz).

Parameters:vel (pos,) – position and velocity grids as generated by make_grid
Returns:phase_space_vectors
Return type:n * m × 1 array
aisim.atoms.combine_weights(pos_weights, vel_weights)

Combine the weights of a position and velocity grid.

Complements _combine_grids.

Parameters:vel_weights (pos_weights,) – weights of a position and velocity grids.
aisim.atoms.create_ensemble_from_grids(pos_params, vel_params, **kwargs)

Create an atomic ensemble from evenly spaced position and velocity grids.

The resulting position and velocity grids are evenly spaced on polar coordinates.

Parameters:
  • vel_params (pos_params,) – Dictionary containing the parameters determining the position and velocity distributions of the atomic ensemble. They each have to contain the arguments described in the docstring of make_grid, i.e. std_rho, std_z (required), n_rho, n_theta, n_z, m_std_rho, m_std_z, weight, optional.
  • **kwargs – Optional keyworded arguments passed to AtomicEnsemble
Returns:

ensemble – Atomic ensemble contains all possible combinations of the position and velocity grid as phase space vectors. They vectors are weighted according to the combined (multiplied) weights of the respective position and velocity distributions ccording to the weight arguments in pos_params and vel_params

Return type:

AtomicEnsemble

aisim.atoms.create_random_ensemble_from_gaussian_distribution(pos_params, vel_params, n_samples, seed=None, **kwargs)

Random atomic ensemble from normal position and velocity distributions.

Parameters:
  • vel_params (pos_params,) – Dictionary containing the parameters determining the position and velocity distributions of the atomic ensemble. Entries for position space are ‘mean_x’,’std_x’ ,’mean_y’, ‘std_y’, ‘mean_z’, ‘std_z’. Entries for velocity space are ‘mean_vx’,’std_vx’, ‘mean_vy’, ‘std_vy’,’mean_vz’, ‘std_vz’.
  • n_samples (float) – number of random samples
  • seed (int or 1-d array_like, optional) – Set the seed of the random number generator to get predictable samples. If set, this number is passed to numpy.random.seed.
  • **kwargs – Optional keyworded arguments passed to AtomicEnsemble
Returns:

ensemble – Atomic ensemble containing the generated phase space vectors.

Return type:

AtomicEnsemble

aisim.atoms.make_grid(std_rho, std_z, n_rho=20, n_theta=36, n_z=1, m_std_rho=3, m_std_z=0)

Evenly spaced grid of positions (or velocities) and weights.

Each of these positions (or velocities) are evenly spaced in polar coordinates and weighted according to a gaussian distribution.

Parameters:
  • std_sigma (std_rho,) – 1/e radius of the position or velocity distribution.
  • n_theta, n_z (n_rho,) – number of grid points per standard deviation along rho and z direction and total number of points along theta, respectively
  • m_std_z (m_std_rho,) – number of standard deviations for the rho and z distribution
Returns:

  • grid (n × 3 array) – Grid of n vectors in carthesian coordinates (x, y, z). In polar coordinates, the grid has this form: [[dRho, 0, -m_std_z*sigma_z/2] [dRho, dTheta, …] [dRho , 2*dTheta, …] [… , …, 0] [dRho , <2*pi, …] [2*dRho , 0, …] [2*dRho , dTheta, …] [… , …, … [m_std_rho*sigma_rho , <2*pi, +m_std_z*sigma_z/2]]
  • weights (1 × n array) – weights for each vector in the grid

aisim.beam module

Classes and functions related to the interferometry lasers.

class aisim.beam.IntensityProfile(r_beam, center_rabi_freq)

Bases: object

Class that defines a intensity profile.

Parameters:
  • r_beam (float) – beam radius in m
  • center_rabi_freq (float) – Rabi frequency at center of intensity profile
r_beam

beam radius in m

Type:float
center_rabi_freq

Rabi frequency at center of intensity profile

Type:float
get_rabi_freq(pos)

Rabi frequency at a position of a Gaussian beam.

Parameters:pos ((n × 2) array or (n × 3) array) – positions of the n atoms in two or three dimensions (x, y, [z]).
Returns:rabi_freqs – the Rabi frequencies for the n positions
Return type:array
class aisim.beam.Wavefront(r_beam, coeff)

Bases: object

Class that defines a wavefront.

Parameters:
  • r_beam (float) – beam radius in m
  • coeff (list) – list of 36 Zernike coefficients in multiples of the wavelength
r_beam

beam radius in m

Type:float
coeff

list of 36 Zernike coefficients in multiples of the wavelength

Type:list
get_value(pos)

Get the wavefront at a position.

Parameters:pos (n × 3 array) – array of position vectors (x, y, z) where the wavefront is probed
Returns:wf – The value of the wavefront at the positions
Return type:nd array
plot(ax=None)

Plot the wavefront.

Parameters:ax (Axis , optional) – If axis is provided, they will be used for the plot. if not provided, a new plot will automatically be created.
plot_coeff(ax=None)

Plot the coefficients as a bar chart.

Parameters:ax (Axis , optional) – If axis is provided, they will be used for the plot. if not provided, a new plot will automatically be created.
classmethod zern_iso(rho, theta, coeff, r_beam)

Calculate the sum of the first 36 Zernike polynomials.

Based on ISO24157:2008.

Parameters:
  • theta (rho,) – Polar coordinates of the position where the sum of Zernike polynomials should be calculated.
  • coeff (array) – first 36 Zernike coefficients
  • r_beam (float) – radius of the wavefront
Returns:

values – value(s) of the wavefront at the probed position

Return type:

float or array of float

class aisim.beam.Wavevectors(k1=8055366, k2=-8055366)

Bases: object

Class that defines the wave vectors of the two Ramen beams.

Parameters:k2 (k1,) – 1D wave vectors (wavenumber) in z-direction of the two Raman beams in rad/m, defaults to 2*pi/780e-9
k1, k2

1D wave vectors (wavenumber) in z-direction of the two Raman beams in rad/m

Type:float
doppler_shift(atoms)

Calculate the Doppler shifts for an atomic ensemble.

Parameters:atoms (AtomicEnsemble) – an atomic enemble with a finite velocity in the z direction
Returns:dopler_shift – Doppler shift of each atom in the ensemble in rad/s
Return type:1d array
aisim.beam.gen_wavefront(r_beam, std=0)

Create an artificial wavefront.

Parameters:
  • r_beam (float) – beam radius
  • std (float) – standard deviation of each Zernike polynomial coefficient in multiples of the wavelength.
Returns:

wf – artificial wavefront

Return type:

Wavefront

aisim.convert module

Functions to convert various quantities.

aisim.convert.cart2pol(cart)

Convert vectors in cartesian coordinates to polar coordinates.

Parameters:cart (n × 3 array) – array of n vectors in cartesian coordinates (x, y, z)
Returns:pol – array of n vectors in polar coordinates (rho, theta, z)
Return type:n × 3 array
aisim.convert.kb = 1.3806488e-23

Boltzmann constant in J/K.

aisim.convert.mass = {'Rb87': 1.443161e-25}

Atomic mass in kg.

aisim.convert.phase_to_grav(phase, T, keff)

Convert a phase shift to gravitational acceleration.

Takes the phase shift measured in a Mach Zehnder atom interferometer and converts it to the corresponding gravtional accleration.

Parameters:
  • phase (float) – Interferometer phase
  • T (float) – interferometer time
  • keff (float) – effective wavenumber
Returns:

gravitational acceleration in m/s²

Return type:

float

aisim.convert.pol2cart(pol)

Convert vectors in polar coordinates to cartesian coordinates.

Parameters:pol (n × 3 array) – array of n vectors in polar coordinates (rho, theta, z)
Returns:cart – array of n vectors in cartesian coordinates (x, y, z)
Return type:n × 3 array
aisim.convert.temp(sigma_v, species='Rb87')

Calculate the temperature of an atomic cloud from its velocity spread.

Parameters:
  • sigma_v (float) – velocity spread (1 sigma) in meters per second
  • species (str) – the atomic species, has to be a key in mass
Returns:

temp – temperature of the cloud in Kelvin

Return type:

float

Raises:

ValueError – If negative velocity spread is passed

aisim.convert.vel_from_temp(temp, species='Rb87')

Calculate the velocity spread (1 sigma) from the temperature of the cloud.

Parameters:
  • temp (float) – temperature of the cloud in Kelvin
  • species (str) – the atomic species, has to be a key in mass
Returns:

vel – velocity spread (1 sigma)

Return type:

float

Raises:

ValueError – If negative temperature is passed

aisim.det module

Classes and functions related to the detection system.

class aisim.det.Detector(t_det, **kwargs)

Bases: object

A generic detection zone.

This is only a template without functionality. Deriving classes have to implement _detected_idx.

Parameters:
  • t_det (float) – time of the detection
  • **kwargs – Additional arguments used by classes that inherit from this class. All keyworded arguments are stored as attribues.
t_det

time of the detection

Type:float
** kwargs

all keyworded arguments that are passed upon creation

detected_atoms(atoms)

Determine wheter a position is within the detection zone.

Returns a new AtomicEnsemble object containing only the detected phase space vectors.

Parameters:atoms (AtomicEnsemble) – the atomic ensemble
Returns:detected_atoms – atomic ensemble containing only phase space vectors that are eventually detected
Return type:AtomicEnsemble
class aisim.det.PolarDetector(t_det, **kwargs)

Bases: aisim.det.Detector

A spherical detection zone.

All atoms within a circle with the specified radius within the x-y plane will be detected.

Parameters:
  • t_det (float) – time of the detection
  • r_det – radius of the spherical detection zone
t_det

time of the detection

Type:float
** kwargs

all keyworded arguments that are passed upon creation

class aisim.det.SphericalDetector(t_det, **kwargs)

Bases: aisim.det.Detector

A spherical detection zone.

All atoms within a sphere of the specified radius will be detected.

Parameters:
  • t_det (float) – time of the detection
  • r_det – radius of the spherical detection zone
t_det

time of the detection

Type:float
** kwargs

all keyworded arguments that are passed upon creation

aisim.prop module

Classes to propagate atomic ensembles.

class aisim.prop.FreePropagator(time_delta, **kwargs)

Bases: aisim.prop.Propagator

Propagator implementing free propagation without light-matter interaction.

Parameters:time_delta (float) – time that should be propagated
class aisim.prop.Propagator(time_delta, **kwargs)

Bases: object

A generic propagator.

This is just a template class without an implemented propagation matrix.

Parameters:
  • time_delta (float) – time that should be propagated
  • **kwargs – Additional arguments used by classes that inherit from this class. All keyworded arguments are stored as attribues.
propagate(atoms)

Propagate an atomic ensemble.

Parameters:atoms (AtomicEnsemble) – atomic ensemble that should be is propagated
class aisim.prop.TwoLevelTransitionPropagator(time_delta, intensity_profile, wave_vectors=None, wf=None, phase_scan=0)

Bases: aisim.prop.Propagator

A time propagator of an effective Raman two-level system.

Parameters:
  • time_delta (float) – length of pulse
  • intensity_profile (IntensityProfile) – Intensity profile of the interferometry lasers
  • wave_vectors (Wavevectors) – wave vectors of the two Raman beams for calculation of Doppler shifts
  • wf (Wavefront , optional) – wavefront aberrations of the interferometry beam
  • phase_scan (float) – effective phase for fringe scans

Notes

The propagator is for example defined in [1].

References

[1] Young, B. C., Kasevich, M., & Chu, S. (1997). Precision atom interferometry with light pulses. In P. R. Berman (Ed.), Atom Interferometry (pp. 363–406). Academic Press. https://doi.org/10.1016/B978-012092460-8/50010-2

aisim.sims module