API

IonSim.IonType
Ion

The physical parameters defining an isolated ion's internal structure.

IonSim.IonConfigurationType
IonConfiguration

Physical configuration of ions. Stores a collection of ions and information about the interactions of their center of mass motion.

IonSim.IonInstanceType

IonInstance(selected_sublevels::Vector{Tuple}[, starkshift::Dict]) Ion instance of some species

selected_sublevels specifies which energy sublevels will be present in the Hilbert space of this Ion instance, as a subset of all possible sublevels.

Each element of selected_sublevels is a 2-element Tuple (level, sublevels), with the first element being the name of a level and the second specifying which sublevels should be included. Allowed sublevels are those whose magnetic quantum number m is in the set {-f, -f+1, -f+2, ... f-1, f}, where f is the total angular momentum quantum number of level. For each level specified there are three allowed options to specify the set of sublevels to include:

  • sublevels::Real: Includes only one m = sublevels
  • sublevels::Vector{Real}: Includes all sublevels whose magnetic quantum number m is in sublevels
  • sublevels = "all": Includes all allowed sublevels

If instead selected_sublevels = "all", then all sublevels of all levels are included.

Omission of a level in selected_sublevels will exclude all sublevels.

Fields

  • species_properties::NamedTuple: Contains constants specifying parameters specific to species
  • sublevels::Vector{Tuple{String,Real}}: List of all sublevels present in the Hilbert space
  • sublevel_aliases::Dict{String,Tuple}: Dict specifying aliases assigned to sublevels, in the format alias => sublevel
  • shape::Vector{Int}: Dimension of the Hilbert space
  • stark_shift::OrderedDict: A dictionary with keys denoting the selected levels and values, a real number for describing a shift of the level's energy. This is just a convenient way to add Stark shifts to the simulation without additional resources
  • ionnumber: When the ion is added to an IonConfiguration, this value keeps track of its order in the chain
  • position: When the ion is added to an IonConfiguration, this value keeps track of its physical position in meters
IonSim.IonPropertiesType

IonProperties type. Required keywords

  • mass: Mass of ion in kg
  • charge: Charge of ion in units of elementary charge
  • full_level_structure: OrderedDict describing properties of each energy level
    • key::String: Name of energy level. Spectroscopic notation encouraged, e.g. "S1/2,f=1"
    • value::NamedTuple(:n, :l, :j, :f, :E): Quantum numbers n, l, j, f, and energy E (in Hz)
  • full_transitions: Dict of all allowed transitions between energy levels
    • key::Tuple{String,String} Pair of levels, ordered (lower, upper) in energy
    • value::NamedTuple(:multipole, :einsteinA): Leading-order multipole of the transition (e.g. "E1", "E2") and Einstein A coefficient (between fine structure levels only; hyperfine factors are calculated when needed)

Optional keywords

  • default_sublevel_selection: Default value of selected_sublevels argument in Ion constructor
  • gfactors: Dict(level::String => g::Real) Custom Landé g-factors, if contributions from higher-than-first-order perturbations are desired
  • nonlinear_zeeman: Dict describing nonlinear contributions to Zeeman shift of certain sublevels
    • key::Tuple{String,Real}: sublevel name
    • value::Function(B::Real): Nonlinear term(s) of Zeeman shift. Full Zeeman shift will be calculated as the sum of the usual linear term and this function
IonSim.IonSimBasisType
IonSimBasis

An abstract type for specialized bases, which are unique to IonSim.

IonSim.LaserType
Laser(;λ=nothing, E=0, Δ=0, ϵ=(x̂+ŷ)/√2, k=ẑ, ϕ=0, pointing::Array{Tuple{Int,Real}})

The physical parameters defining laser light. args

  • λ::Union{Real,Nothing}: the wavelength of the laser in meters
  • E::Union{Function,Real}: magnitude of the E-field in V/m
  • Δ: static detuning from f = c/λ in [Hz]
  • ϵ::NamedTuple: (ϵ.x, ϵ.y, ϵ.z), polarization direction, requires norm of 1
  • k::NamedTuple: (k.x, k.y, k.z), propagation direction, requires norm of 1
  • ϕ::Union{Function,Real}: time-dependent phase. of course, this can also be used to model a time-dependent detuning. Units are in radians. Note: if this is set to a function of time, then when constructing a Hamiltonian with the hamiltonian function, the units of time will be as specified by the timescale keyword argument.
  • pointing: an array of Tuple{Int,Real} for describing ion-laser pointing configuration. (first element of the tuple is the index for an ion and the second element is the scaling factor for the laser's Efield which must be between 0 and 1).
IonSim.LinearChainType
LinearChain(;
        ions::Vector{Ion}, com_frequencies::NamedTuple{(:x,:y,:z)}, 
        vibrational_modes::NamedTuple{(:x,:y,:z),Tuple{Vararg{Vector{VibrationalMode},3}}}
    )

Contains all of the information necessary to describe a collection of ions trapped in a 3D harmonic potential and forming a linear coulomb crystal.

user-defined fields

  • ions::Vector{Ion}: a list of ions that compose the linear Coulomb crystal
  • com_frequencies::NamedTuple{(:x,:y,:z),Tuple{Vararg{Vector{VibrationalMode},3}}}: Describes the COM frequencies (x=ν_x, y=ν_y, z=ν_z). The $z$-axis is taken to be parallel to the crystal's symmetry axis.
  • vibrational_modes::NamedTuple{(:x,:y,:z)}: e.g. vibrational_modes=(x=[1], y=[2], z=[1,2]). Specifies the axis and a list of integers which correspond to the $i^{th}$ farthest mode away from the COM for that axis. For example, vibrational_modes=(z=[2]) would specify the axial stretch mode. These are the modes that will be modeled in the chain. Note: vibrational_modes=(x=[],y=[],z=[1]), vibrational_modes=(y=[],z=[1]) and vibrational_modes=(;z=[1]) are all acceptable and equivalent.

derived fields

  • full_normal_mode_description::NamedTuple{(:x,:y,:z)}: For each axis, this contains an array of tuples where the first element is a vibrational frequency [Hz] and the second element is a vector describing the corresponding normalized normal mode structure.
IonSim.TrapType
Trap(;
        configuration::LinearChain, B::Real=0, Bhat::NamedTuple{(:x,:y,:z)=ẑ, ∇B::Real=0,
         δB::Union{Real,Function}=0, lasers::Vector{Laser}
)

Information necessary to describe the Hamiltonian for a collection of ions in a linear chain interacting with laser light. user-defined fields

  • configuration<:LinearChain
  • B: A real value describing the mean magnitude of the B-field [Tesla].
  • Bhat::NamedTuple{(:x,:y,:z)}: Describes the direction of the B-field (defaults to ẑ).
  • ∇B: Magnitude of the B-field gradient. We assume that the gradient always points along the z-direction. [Tesla / meter]
  • δB::Function: Time-dependence of the B-field [Tesla]
  • lasers::Array{<:Laser}: For each laser in the array, the pointing field should contain an array of Tuple{Int,Real}. The first element specifies the index of an ion in the ions field that the laser interacts with. The second element specifies a scaling factor for the strength of that interaction (to be used, e.g., for modeling cross-talk).

derived fields

  • _cnst_δB::Bool: A Boolean flag signifying whether or not δB is a constant function.

  • basis<:CompositeBasis: The basis for describing the combined system, ions + vibrational modes. If constructing the Hamiltonian explictly (with hamiltonian), then the ordering of the basis is set, by convention, as $ion₁ ⊗ ion₂ ⊗ ... ⊗ ion_N ⊗ mode₁ ⊗ mode₂ ⊗ ... ⊗ mode_N$, where the ion bases are ordered according to the order in T.configuration.ions and the vibrational modes are ordered according to the order in [T.configuration.vibrational_modes.x, T.configuration.vibrational_modes.y, T.configuration.vibrational_modes.z]. E.g. for:

    chain = LinearChain(ions=[C1, C2], com_frequencies=(x=2e6,y=2e6,z=1e6), selected_modes=(x=[1, 2], y=[], z=[1]))

    The ordering of the basis would be

    C1.basis ⊗ C2.basis ⊗ chain.vibrational_modes.x[1].basis ⊗ chain.vibrational_modes.x[2].basis ⊗ chain.vibrational_modes.z[1].basis

    Otherwise, the ordering is according to the form of the initial state used in the solver.

IonSim.VibrationalModeType
VibrationalMode(
        ν::Real, mode_structure::Vector{Real}, δν::Union{Function,Real}=0.; N::Int=10,
        axis::NamedTuple{(:x,:y,:z)}=ẑ
    )

user-defined fields

  • ν::Real: frequency [Hz]
  • mode_structure::Vector{Real}: The normalized eigenvector describing the collective motion of the ions belonging to this mode.
  • δν::Union{Function,Real}: Either a function describing time-dependent fluctuations of ν or a real number which will be converted to the constant function t -> δν.
  • N::Int: Highest level included in the Hilbert space.
  • axis::NamedTuple{(:x,:y,:z)}: The axis of symmetry for the vibration. This must lie along one of the basis vectors , or .

derived fields

  • shape::Vector{Int}: Indicates dimension of used Hilbert space (=[N+1]).
  • _cnst_δν::Bool: A Boolean flag signifying whether or not δν is a constant function.

Note: the iᵗʰ Fock state (|i⟩) can be obtained by indexing as v=VibrationalMode(...); v[i]

IonSim.AnmMethod
Anm(N::Real, com::NamedTuple{(:x,:y,:z)}, axis::NamedTuple{(:x,:y,:z)})

Computes the normal modes and corresponding trap frequencies along a particular axis for a collection of N ions in a linear Coloumb crystal and returns an array of tuples with first element the frequency of the normal mode and 2nd element the corresponding eigenvector.

com should be a NamedTuple of COM frequences for the different axes: (x<:Real, y<:Real, z<:Real), where the $z$-axis is taken to be parallel to the axis of the crystal.

IonSim.Efield_from_pi_time!Method
Efield_from_pi_time!(
        pi_time::Real, Bhat::NamedTuple{(:x,:y,:z)}, laser::Laser, ion::Ion, 
        transition::Union{Tuple{String,String},Vector{<:String}}
    )

Same as Efield_from_pi_time, but updates laser[:E] in-place.

IonSim.Efield_from_pi_timeMethod
Efield_from_pi_time(
    pi_time::Real, Bhat::NamedTuple{(:x,:y,:z)}, laser::Laser, ion::Ion, 
    transition::Union{Tuple{String,String},Vector{<:String}}
)

Compute the E-field needed to get a certain pi_time with a certain laser-ion transition.

Alternatively, one may use

Efield_from_pi_time(
            pi_time::Real, T::Trap, laser_index::Int, ion_index::Int, 
            transition::Union{Tuple{String,String},Vector{<:String}}
        )

which is the same as Efield_from_pi_time(pi_time, T.Bhat, T.lasers[laser_index], T.configuration.ions[ion_index], transition)

IonSim.Efield_from_rabi_frequencyMethod
Efield_from_rabi_frequency(
    Ω::Real, Bhat::NamedTuple{(:x,:y,:z)}, laser::Laser, ion::Ion, 
    transition::Union{Tuple{String,String},Vector{<:String}}
)

Compute the Efield needed to get a certain rabi frequency Ω with a certain laser-ion transition.

Alternatively, one may use

Efield_from_rabi_frequency(
            Ω::Real, T::Trap, laser_index::Int, ion_index::Int, 
            transition::Union{Tuple{String,String},Vector{<:String}}
        )

which is the same as Efield_from_rabi_frequency(pi_time, T.Bhat, T.lasers[laser_index], T.configuration.ions[ion_index], transition)

IonSim.alias2sublevelMethod
alias2sublevel(I::Ion, alias::String)

Returns the sublevel corresponding to the given alias alias of I. Inverse function of sublevelalias.

IonSim.characteristic_length_scaleMethod
characteristic_length_scale(M::Real, ν::Real)

Returns the characteristic length scale for a linear chain of identical ions of mass M and with axial trap frequency 2π × ν.

IonSim.clear_sublevel_alias!Method
clear_sublevel_alias!(I::Ion, sublevel)

Erases the assignment of an alias to sublevel of Ion I. Accepts either the full sublevel Tuple{String,Real} or its alias String. Also accepts a vector of sublevels to clear multiple alias assignments in a single call.

IonSim.coherentstateMethod
coherentstate(v::VibrationalMode, α::Number)

Returns a coherent state on v with complex amplitude $α$.

IonSim.createMethod
create(v::VibrationalMode)

returns the creation operator for v such that: create(v) * v[i] = √(i+1) * v[i+1].

IonSim.destroyMethod
destroy(v::VibrationalMode)

Returns the destruction operator for v such that: destroy(v) * v[i] = √i * v[i-1].

IonSim.einsteinAMethod
einsteinA(I::Ion, Lpair::Tuple)

Returns Einstein A coefficient corresponding to the transition Lpair[1] -> Lpair[2]. The first level must be the lower level and the second must be the upper level.

IonSim.energyMethod
energy(I::Ion, level::String)

Returns the energy of level of I.

IonSim.energyMethod
energy(I::Ion, sublevel; B=0, ignore_starkshift=false)

Returns energy of sublevel of I. A Zeeman shift may be included by setting the value of the magnetic field B. The Stark shift may be omitted by setting ignore_starkshift=true.

IonSim.get_basisMethod
get_basis(T::Trap)

Returns the composite basis describing the Hilbert space for T.

IonSim.get_vibrational_modesMethod
get_vibrational_modes(lc::LinearChain)

Returns an array of all of the selected VibrationalModes in the LinearChain. The order is [lc.x..., lc.y..., lc.z...].

IonSim.get_ηMethod
get_η(V::VibrationalMode, L::Laser, I::Ion)

The Lamb-Dicke parameter: $|k|cos(\theta)\sqrt{\frac{\hbar}{2m\nu}}$ for a given vibrational mode, ion and laser.

IonSim.global_beam!Method
global_beam!(T::Trap, laser::Laser)

Set laser to shine with full intensity on all ions in Trap.

IonSim.hamiltonianMethod
hamiltonian(
        T::Trap; timescale::Real=1e-6, lamb_dicke_order::Union{Vector{Int},Int}=1, 
        rwa_cutoff::Real=Inf, displacement="truncated", time_dependent_eta=false
    )

Constructs the Hamiltonian for T as a function of time. Return type is a function h(t::Real, ψ) that, itself, returns a QuantumOptics.SparseOperator.

args

  • timescale: e.g. a value of 1e-6 will take time to be in $\mu s$

  • lamb_dicke_order: Only consider terms that change the phonon number by up to this value. If this is an Int, then the cutoff is applied to all modes. If this is a Vector{Int}, then lamb_dicke_order[i] is applied to the iᵗʰ mode, according to the order in T.basis. Note: this isn't quite the same thing as the Lamb-Dicke approximation since setting lamb_dicke_order=1 will retain, for example, terms proportional to $a^\dagger a$.

  • rwa_cutoff: drop terms in the Hamiltonian that oscillate faster than this cutoff.

  • displacement: This can be either "truncated"(default) or "analytic".

    When an atom is irradiated, both the atom's energy and its momentum will generally be affected. For an atom in a harmonic potential, the exchange of momentum can be modeled as a displacement operation $D(α=iηe^{-iνt}) = exp[αa^† - α^*a]$, where $η$ is the Lamb-Dicke parameter, which can be described equivalently as either being proportional to the square root of the ratio of the recoil frequency with the ground state energy of the atom's motion or as the ratio of the spread of the ground state wavefunction to the wavelength of the laser.

    When "truncated" is selected, the matrix elements of $D(α)$ are computed by constructing $α^* a, αa^†$ in a truncated basis (according to the dimension specified in your model) and then exponentiating their difference. This has the advantage, amongst other things, of guaranting unitarity.

    If "analytic" is selected, then the matrix elements are computed assuming an infinite- dimensional Hilbert space.

    For small displacements ($η ≪ N$, where $N$ is the dimension of the motion's Hilbert space), both of these methods will be good approximations.

  • time_dependent_eta::Bool: In addition to impacting the vibrational subspace directly, a change in the trap frequency, $δν$, will also change the Lamb-Dicke parameter. Since typically $δν≪ν$, this effect will be small $η ≈ η₀(1 + δν/2ν)$ and doesn't warrant the additional computational resources needed to calculate and update it in time. In this case, we can set time_dependent_eta=false (default), which will set $η(t) = η₀$.

IonSim.ionprojectorMethod
ionprojector(obj, sublevels...; only_ions=false)

If obj<:IonConfiguration this will return $|ψ₁⟩⟨ψ₁|⊗...⊗|ψ\_N⟩⟨ψ\_N|⊗𝟙$ where $|ψᵢ⟩$ = obj.ions[i][sublevels[i]] and the identity operator $𝟙$ is over all of the COM modes considered in obj.

If only_ions=true, then the projector is defined only over the ion subspace.

If instead obj<:Trap, then this is the same as obj = Trap.configuration.

IonSim.ionstateMethod
ionstate(object, sublevel)

For object<:Ion and sublevel<:Tuple{String,Real} (full sublevel name) or sublevel<:String (alias), this returns the ket corresponding to the Ion being in the state $|index⟩$. The object can also be an IonConfiguration or Trap instance, in which case $N$ arguments should be given in place of index, where $N$ equals the number of ions in the IonConfiguration or Trap. This will return the state $|index₁⟩⊗|index₂⟩⊗...⊗|index\_N⟩$.

One may also specify sublevel<:Int. If object<:Ion, this will return the ket given by $|index⟩ = (0 ... 1 ... 0)ᵀ$ where the nonzero element in the column vector is located at index.

IonSim.landegfFunction
landegf(l::Real, j::Real, f::Real, i::Real, s::Real=1//2)

Landé g-factor of hyperfine energy level

args

  • l: orbital angular momentum quantum number
  • j: electron total angular momentum quantum number
  • f: total angular momentum quantum number
  • i: nuclear spin angular momentum quantum number
  • s: electronic spin angular momentum quantum number (defaults to 1/2)
IonSim.landegfMethod
landegf(I::Ion, level::String)

landegf for the quantum numbers of level in I.

IonSim.levelsMethod
levels(I::Ion)

Returns array of all energy levels of I.

IonSim.leveltransitionsMethod
leveltransitions(I::Ion)

Returns all allowed transitions between levels of I as a vector of Tuple{String,String}.

IonSim.lifetimeMethod
lifetime(I::Ion, level::String)

Computes lifetime of level by summing the transition rates out of level.

The sum is taken over all levels that the species may have, rather than the levels present in the instance I.

IonSim.linear_equilibrium_positionsMethod
linear_equilibrium_positions(N::Int)

Returns the scaled equilibrium positions of N ions in a harmonic potential, assuming that all ions have the same mass. ref

IonSim.matrix_elementFunction
matrix_element(I::Ion, transition::Tuple, Efield::Real, khat::NamedTuple, ϵhat::NamedTuple, Bhat::NamedTuple=(;z=1))

Computes the matrix elements (units of Hz) between two energy sublevels args

  • I: Ion undergoing transition
  • transition: Tuple of sublevels (full names or aliases) between which the transition is being calculated. Must be formatted such that energy(transition[2]) > energy(transition[1])
  • Efield: Amplitude of driving electric field
  • khat: Unit vector of light wavevector
  • ϵhat: Unit vector of light polarization
  • Bhat: Unit vector of magnetic field
IonSim.matrix_elementMethod
matrix_element(I::Ion, transition::Tuple, T::Trap, laser::Laser, time::Real)

Calls matrix_element(I::Ion, transition::Tuple, Efield::Real, khat::NamedTuple, ϵhat::NamedTuple, Bhat::NamedTuple=(;z=1)) with Efield, khat, and ϵhat evaluated for laser at time time, and Bhat evaluated for T.

One may alternatively replace ion with ion_index::Int, which instead specifies the index of the intended ion within T.

IonSim.numberMethod
number(v::VibrationalMode)

Returns the number operator for v such that: number(v) * v[i] = i * v[i].

IonSim.quantumnumbersMethod
quantumnumbers(I::Ion, level::String)
quantumnumbers(I::Ion, sublevel)

Returns the quantum numbers of an energy level or sublevel of I as a NamedTuple.

If second argument is a level, returns (:n, :i, :s, :l, :j, :f)

If second argument is a sublevel, returns (:n, :i, :s, :l, :j, :f, :m)

IonSim.set_gradient!Method
set_gradient!(
        T::Trap, ion_indxs::Tuple{Int,Int}, transition::Tuple, f::Real
    )

Sets the Bfield gradient in place to achieve a detuning f between the transition of two ions, which are assumed to be of the same species. ion_indxs refer to the ordering of the ions in the chain.

IonSim.set_stark_shift!Method
set_stark_shift!(I::Ion, stark_shift_dict::Dict)

Applies set_stark_shift(I, sublevel, shift) to all pairs sublevel => shift of the Dict stark_shift_dict.

IonSim.set_stark_shift!Method
set_stark_shift!(I::Ion, sublevel, shift::Real)

Applies a stark shift shift to the chosen sublevel of I (overwriting any previously assigned stark shift).

IonSim.set_sublevel_alias!Method
set_sublevel_alias!(I::Ion, aliasassignments)

aliasassignments specifies which aliases should be assigned to which sublevels. There are two options to do this:

  • aliasassignments is a Vector of Tuples, with the first element of each being the sublevel (sublevel::Tuple{String,Real}) and the second being its assigned alias (alias::String)
  • aliasassignments is a Dict with the format alias::String => sublevel::Tuple{String,Real}

Calls set_sublevel_alias!(I, sublevel, alias) for each pair sublevel, alias.

IonSim.set_sublevel_alias!Method
set_sublevel_alias!(I::Ion, sublevel::Tuple{String,Real}, alias::String)

Assigns an alias alias to sublevel of I. This then allows one to pass alias in place of sublevel (for I only) into any function which accepts a sublevel as an argument.

IonSim.sigmaMethod
sigma(ion::Ion, ψ1::sublevel[, ψ2::sublevel])

Returns $|ψ1\rangle\langle ψ2|$, where $|ψ_i\rangle$ corresponds to the state returned by ion[ψᵢ].

If ψ2 is not given, then $|ψ1\rangle\langle ψ1|$ is returned.

IonSim.stark_shiftMethod
stark_shift(I::Ion, sublevel)

Returns the assigned stark shift of sublevel of Ion I.

IonSim.sublevel2levelMethod
sublevel2level(I::Ion, sublevel)

Retuns the energy level of I corresponding to sublevel.

IonSim.sublevelaliasMethod
sublevelalias(I::Ion, sublevel::Tuple{String,Real})

Returns the alias assined to sublevel of I. If no alias is assigned, returns nothing.

IonSim.subleveltransitionsMethod
subleveltransitions(I::Ion)

Returns all allowed transitions between sublevels of I as a vector of Tuple{S,S} where S=Tuple{String,Real}.

IonSim.transitionfrequencyMethod
transitionfrequency(ion::Ion, transition::Tuple, T::Trap; ignore_starkshift=false)

Retuns The frequency of the transition transition including the Zeeman shift experienced by ion given its position in T.

One may alternatively replace ion with ion_index::Int, which instead specifies the index of the intended ion within T.

IonSim.transitionfrequencyMethod
transitionfrequency(I::Ion, transition::Tuple; B=0, ignore_starkshift=false)

transition is a Tuple of two sublevels or levels.

Computes the absolute values of the difference in energies between transition[1] and transition[2].

If between sublevels, then the Zeeman shift may be included by setting the value of the magnetic field B, and Stark shifts may be omitted by setting ignore_starkshift=true.

IonSim.transitionmultipoleMethod
transitionmultipole(I::Ion, Lpair::Tuple)

Returns the transition multiple ('E1', 'E2', etc.) corresponding to the transition Lpair[1] -> Lpair[2]. The first level must be the lower level and the second must be the upper level.

IonSim.transitionwavelengthMethod
transitionwavelength(ion::Ion, transition::Tuple, T::Trap; ignore_starkshift=false)

Retuns The wavelength of the transition transition including the Zeeman shift experienced by ion given its position in T.

One may alternatively replace ion with ion_index::Int, which instead specifies the index of the intended ion within T.

IonSim.transitionwavelengthMethod
transitionwavelength(I::Ion, transition::Tuple; B=0, ignore_starkshift=false)

Returns the wavelength corresponding to transitionfrequency(I::Ion, transition::Tuple; B=0, ignore_starkshift=false).

IonSim.zeeman_shiftMethod
zeeman_shift(I::Ion, sublevel, B::Real)

Returns the Zeeman shift at a magnetic field of B of sublevel of I.

If sublevel has a custom g-factor defined, then this is used. Otherwise, landegf is used to compute the Landé g-factor.

Zeeman shift calculated as $ΔE = (μ_B/ħ) ⋅ g_f ⋅ B ⋅ m / 2π$

IonSim.zeeman_shiftMethod
zeeman_shift(I::Ion, sublevel, T::Trap)

Calls zeeman_shift(I::Ion, sublevel, B::Real) with B evaluated for ion I in T.

One may alternatively replace ion with ion_index::Int, which instead specifies the index of the intended ion within T.

QuantumOpticsBase.coherentthermalstateMethod
coherentthermalstate(v::VibrationalMode, n̄::Real, α::Number; method="truncated)

Returns a displaced thermal state for v, which is created by applying a displacement operation to a thermal state. The mean occupation of the thermal state is and α is the complex amplitude of the displacement.

method can be either "truncated" or "analytic" and this argument determines how the displacement operator is computed (see: displace) .

QuantumOpticsBase.displaceMethod
displace(v::VibrationalMode, α::Number; method="truncated")

Returns the displacement operator $D(α)$ corresponding to v.

If method="truncated" (default), the matrix elements are computed according to $D(α) = exp[αa^† - α^*a]$ where $a$ and $a^†$ live in a truncated Hilbert space of dimension v.N+1. Otherwise if method="analytic", the matrix elements are computed assuming an infinite-dimension Hilbert space. In general, this option will not return a unitary operator.

QuantumOpticsBase.thermalstateMethod
thermalstate(v::VibrationalMode, n̄::Real; method="truncated")

Returns a thermal density matrix with $⟨a^†a⟩ ≈ n̄$. Note: approximate because we are dealing with a finite dimensional Hilbert space that must be normalized.

method can be set to either "truncated" (default) or "analytic". In the former case, the thermal density matrix is generated according to the formula: $ρ_{th} = exp(-νa^†a/T) / Tr [exp(-νa^†a/T)]$. In the later case, the analytic formula, assuming an infinite-dimensional Hilbert space, is used: $[ρ_{th}]_{ij} = δ_{ij} \frac{nⁱ}{(n+1)^{i+1}}.$

IonSim.analytical.rabi_flopMethod
rabi_flop(tspan, Ω::Real, η::Real, n̄::Real; s::Int=0) <br>

Single ion rabi flop. Returns: $\sum_{n=0}^∞ p_n sin^2(\Omega_n t)$ <br> with $\Omega_n = Ωe^{-η^2/2}η^s\sqrt{\frac{n!}{(n+s)!}}L_{n}^{s}(η^2)$ <br> where $s$ is the order of the (blue) sideband that we are driving and $L_{n}^{s}$ is the associated Laguerre polynomial. ref

IonSim.analytical.two_ion_msMethod
two_ion_ms(tspan, Ω::Real, ν::Real, δ::Real, η::Real, n̄::Real)

ref <br> Assumes vibrational mode starts in a thermal state with: $\langle a^\dagger a\rangle = n̄$ and ions start in doubly ground state. Returns (ρgg, ρee), the population in the doubly ground and doubly excited state, respectively. $[Ω], [ν], [δ] = Hz$

IonSim.PhysicalConstants.c_rank1Constant

c_rank1: Matrix of spherical basis vectors (defined in e.g. Quantum dynamics of cold trapped ions with application to quantum computation, Appl. Phys. B 66, 181-190 (1998).

Useful for converting between coordinates of rank-1 spherical tensors and Cartesian coordinates

IonSim.PhysicalConstants.c_rank2Constant

c_rank2: Matrix of spherical basis rank-2 tensors (defined in e.g. Quantum dynamics of cold trapped ions with application to quantum computation, Appl. Phys. B 66, 181-190 (1998).

Useful for converting between coordinates of rank-2 spherical tensors and Cartesian coordinates