CollectiveSpins.geometry.boxMethod
geometry.box(a, b, c; Nx=2, Ny=2, Nz=2)

Positions of spins on a orthorhombic lattice.

The lattice starts at the origin and continues into positive x, y and z direction.

CollectiveSpins.geometry.chainMethod
geometry.chain(a, N)

Positions of spins on a chain in x-direction.

The chain starts at the origin and continues into positive x-direction.

Arguments

  • a: Spin-spin distance.
  • N: Number of spins
CollectiveSpins.geometry.cubeMethod
geometry.cube(a; Nx=2, Ny=2, Nz=2)

Positions of spins on a cubic lattice with edge length a.

The lattice starts at the origin and continues into positive x, y and z direction.

CollectiveSpins.geometry.excitation_phasesMethod
geometry.excitation_phases(k, pos)

Calculate the excitation phases created in a geometry by an incident plane wave with wave vector k.

Arguments:
* `k`: Wave vector of the incident plane wave.
* `pos`: List of atomic positions.
CollectiveSpins.geometry.hexagonalMethod
geometry.hexagonal(a; Nr=1)

Positions of spins on a hexagonal lattice in the xy-plane.

The hexagonal lattice consists of Nr rings.

CollectiveSpins.geometry.lhc1Method
geometry.lhc1()

Geometry of the LHC-I light harvesting complex in the xy-plane centered around `z=0` with a mean interatomic distance of `1.0`. Positions `1` to `32` constitute the outer ring, while the inner ring comprises positions `33` to `38`.
CollectiveSpins.geometry.lhc1_pMethod
geometry.lhc1_p()

Normalized polarization vectors in an LHC-I light harvesting complex. Entries `1` to `32` pertain to the outer ring, while entries `33` to `38` give the inner ring's dipole orientations.

Note: This function does not return a geometry, butrather te polarizations of dipoles in the LHC-I complex.
CollectiveSpins.geometry.rectangleMethod
geometry.rectangle(a, b; Nx=2, Ny=2)

Positions of spins on a rectangular lattice in the xy-plane.

The lattice starts at the origin and continues into positive x and y direction.

Arguments

  • a: Spin-spin distance in x-direction.
  • b: Spin-spin distance in y-direction.
  • Nx=2: Number of spins into x direction.
  • Ny=2: Number of spins into y direction.
CollectiveSpins.geometry.ringMethod
geometry.ring(a, N; distance = false)

Ring of N particles with radius a. If distance is set to true, then a gives the distance between particles and the radius is determined accordingly.

CollectiveSpins.geometry.ring_p_tangentialMethod
geometry.ring_p_tangential(N)

Return tangential polarization vectors for a ring of N atoms in the xy-plane created by ring(a, N).

Arguments:
* `N`: Number of emitters.
CollectiveSpins.geometry.squareMethod
geometry.square(a; Nx=2, Ny=2)

Positions of spins on a square lattice in the xy-plane with distance a.

The lattice starts at the origin and continues into positive x and y direction.

CollectiveSpins.geometry.triangleMethod
geometry.triangle(a)

Positions of spins on a equilateral triangle in the xy-plane with edge length a.

The positions are: (0,0,0), (a,0,0), (a/2, h, 0)

CollectiveSpins.meanfield.ProductStateType

Class describing a Meanfield state (Product state).

The data layout is [sx1 sx2 ... sy1 sy2 ... sz1 sz2 ...]

Arguments

  • N: Number of spins.
  • data: Vector of length 3*N.
CollectiveSpins.meanfield.blochstateMethod
meanfield.blochstate(phi, theta[, N=1])

Product state of N single spin Bloch states.

All spins have the same azimuthal angle phi and polar angle theta.

CollectiveSpins.meanfield.timeevolutionMethod
meanfield.timeevolution(T, S::SpinCollection, state0[; fout])

Meanfield time evolution.

Arguments

  • T: Points of time for which output will be generated.
  • S: SpinCollection describing the system.
  • state0: Initial ProductState.
  • fout (optional): Function with signature fout(t, state) that is called whenever output should be generated.
CollectiveSpins.meanfield.timeevolution_symmetricMethod
meanfield.timeevolution_symmetric(T, state0, Ωeff, Γeff[; γ, δ0, fout])

Symmetric meanfield time evolution.

Arguments

  • T: Points of time for which output will be generated.
  • state0: Initial ProductState.
  • Ωeff: Effective dipole-dipole interaction.
  • Γeff: Effective collective decay rate.
  • γ=1: Single spin decay rate.
  • δ0=0: Phase shift for rotated symmetric meanfield time evolution.
  • fout (optional): Function with signature fout(t, state) that is called whenever output should be generated.
CollectiveSpins.CavityModeType

A class representing a single mode in a cavity.

Arguments

  • cutoff: Number of included Fock states.
  • delta=0 Detuning.
  • eta=0: Pump strength.
  • kappa=0: Decay rate.
CollectiveSpins.CavitySpinCollectionType

A class representing a cavity coupled to many spins.

Arguments

  • cavity: A CavityMode.
  • spincollection: A SpinCollection.
  • g: A vector defing the coupling strengths between the i-th spin and the cavity mode. Alternatively a single number can be given for identical coupling for all spins.
CollectiveSpins.SpinType

A class representing a single spin.

A spin is defined by its position and its detuning to a main frequency.

Arguments

  • position: A vector defining a point in R3.
  • delta: Detuning.
CollectiveSpins.SpinCollectionType

A class representing a system consisting of many spins.

Arguments

  • spins: Vector of spins.
  • polarizations: Unit vectors defining the directions of the spins.
  • gammas: Decay rates.
CollectiveSpins.SpinCollectionMethod

Create a SpinCollection without explicitly creating single spins.

Arguments

  • positions: Vector containing the positions of all single spins.
  • polarizations: Unit vectors defining the directions of the spins.
  • deltas=0: Detunings.
  • gammas=1: Decay rates.
CollectiveSpins.SystemType

Abstract base class for all systems defined in this library.

Currently there are the following concrete systems:

  • Spin
  • SpinCollection
  • CavityMode
  • CavitySpinCollection
CollectiveSpins.interaction.FMethod
interaction.F(ri::Vector, rj::Vector, µi::Vector, µj::Vector)

General F function for arbitrary positions and dipole orientations.

Arguments:

  • ri: Position of first spin
  • rj: Position of second spin
  • µi: Dipole orientation of first spin.
  • µj: Dipole orientation of second spin.
CollectiveSpins.interaction.GMethod
interaction.G(ri::Vector, rj::Vector, µi::Vector, µj::Vector)

General G function for arbitrary positions and dipole orientations.

Arguments:

  • ri: Position of first spin
  • rj: Position of second spin
  • µi: Dipole orientation of first spin.
  • µj: Dipole orientation of second spin.
CollectiveSpins.interaction.GammaFunction
interaction.Gamma(ri::Vector, rj::Vector, µi::Vector, µj::Vector, γi::Real=1, γj::Real=1)

Arguments:

  • ri: Position of first spin
  • rj: Position of second spin
  • µi: Dipole orientation of first spin.
  • µj: Dipole orientation of second spin.
  • γi: Decay rate of first spin.
  • γj: Decay rate of second spin.

Note that the dipole moments μi and μj are normalized internally. To account for dipole moments with different lengths you need to scale the decay rates γi and γj, respectively.

CollectiveSpins.interaction.GreenTensorFunction
GreenTensor(r::Vector, k::Number=2π)

Calculate the Green's Tensor at position r for wave number k defined by

\[G = e^{ikr}\Big[\left(\frac{1}{kr} + \frac{i}{(kr)^2} - \frac{1}{(kr)^3}\right)*I - \textbf{r}\textbf{r}^T\left(\frac{1}{kr} + \frac{3i}{(kr)^2} - \frac{3}{(kr)^3}\right)\Big]\]

Choosing k=2π corresponds to the position r being given in units of the wavelength associated with the dipole transition.

Returns a 3×3 complex Matrix.

CollectiveSpins.interaction.OmegaFunction
interaction.Omega(ri::Vector, rj::Vector, µi::Vector, µj::Vector, γi::Real=1, γj::Real=1)

Arguments:

  • ri: Position of first spin
  • rj: Position of second spin
  • µi: Dipole orientation of first spin.
  • µj: Dipole orientation of second spin.
  • γi: Decay rate of first spin.
  • γj: Decay rate of second spin.

Note that the dipole moments μi and μj are normalized internally. To account for dipole moments with different lengths you need to scale the decay rates γi and γj, respectively.

CollectiveSpins.field.intensityMethod
field.intensity(r, S, state)

Arguments:
* `r`: Position (vector in R3).
* `S`: SpinCollection System
* `state`: system state (Independent, Meanfield, MPC or Quantum)
CollectiveSpins.field.intensityMethod
field.intensity(r, S, state, sm)

Arguments:
* `r`: Position (vector in R3)
* `S`: SpinCollection System
* `state`: system state (Independent, Meanfield, MPC or Quantum)
* `sm`: function sm(j::Int) that gives the sigma-minus operator for the j-th particle.
CollectiveSpins.mpc.MPCStateType

Class describing a MPC state (Product state + Correlations).

The data layout is vector that in matrix form looks like

Cxx Cxy Cyy Cxz Czz Cyz

where the Cij are the appropriate correlation matrices. The expectation values sx, sy and sz are the diagonals of the matrices Cxx, Cyy and Czz, respectively.

Arguments

  • N: Number of spins.
  • data: Vector of length (3N)(2*N+1).
CollectiveSpins.mpc.blochstateMethod
mpc.blochstate(phi, theta[, N=1])

Product state of N single spin Bloch states.

All spins have the same azimuthal angle phi and polar angle theta.

CollectiveSpins.mpc.rotateMethod
mpc.rotate(axis, angles, state)

Rotations on the Bloch sphere for the given MPCState.

Arguments

  • axis: Rotation axis.
  • angles: Rotation angle(s).
  • state: MPCState that should be rotated.
CollectiveSpins.mpc.squeezeMethod
mpc.squeeze(axis, χT, state0)

Spin squeezing along an arbitrary axis.

Arguments

  • axis: Squeezing axis.
  • χT: Squeezing strength.
  • state0: MPCState that should be squeezed.
CollectiveSpins.mpc.squeeze_sxMethod
mpc.squeeze_sx(χT, state)

Spin squeezing along sx.

Arguments

  • χT: Squeezing strength.
  • state0: MPCState that should be squeezed.
CollectiveSpins.mpc.timeevolutionMethod
mpc.timeevolution(T, S::SpinCollection, state0[; fout])

MPC time evolution.

Arguments

  • T: Points of time for which output will be generated.
  • S: SpinCollection describing the system.
  • state0: Initial MPCState.
  • fout (optional): Function with signature fout(t, state) that is called whenever output should be generated.
CollectiveSpins.effective_interaction.chainMethod
effective_interaction.chain(a, Θ, N)

Effective Omega and Gamma for an infinite chain.

The calculation is done by adding N spins left and N spins right of a central spin.

Arguments

  • a: Spin-spin distance.
  • θ: Angle between polarization axis and spin chain.
  • N: Number of included spins.
CollectiveSpins.effective_interaction.chain_orthogonalMethod
effective_interaction.chain_orthogonal(a, N)

Effective Omega and Gamma for an infinite chain with orthogonal polarization axis.

The calculation is done by adding N spins left and N spins right of a central spin.

Arguments

  • a: Spin-spin distance.
  • N: Number of included spins.
CollectiveSpins.effective_interaction.cubiclattice_orthogonalMethod
effective_interaction.cubiclattice_orthogonal(a, N)

Effective Omega and Gamma for a infinite cubic lattice.

The polarization axis is orthogonal to the top face of a unit cell and the calculation is done by creating a (2N+1)(2N+1)(2N+1) cubic lattice and calculate the combined interaction for the central spin.

Arguments

  • a: Spin-spin distance.
  • N: Number of included spins.
CollectiveSpins.effective_interaction.hexagonallattice3d_orthogonalMethod
effective_interaction.hexagonallattice3d_orthogonal(a, b, N)

Effective Omega and Gamma for a infinite 3D hexagonal lattice.

The lattice consists of stacked planes of hexagonal lattices where the the polarization axis is orthogonal to the planes. The calculation is done by creating hexagonal lattices with N rings, stacking 2N+1 lattices of this kind above each other and calculating the combined interaction for the central spin.

Arguments

  • a: Spin-spin distance for hexagons.
  • b: Distance between planes of hexagonal lattices
  • N: Number of included spins.
CollectiveSpins.effective_interaction.hexagonallattice_orthogonalMethod
effective_interaction.hexagonallattice_orthogonal(a, N)

Effective Omega and Gamma for a infinite hexagonal lattice.

The polarization axis is orthogonal to the square lattice plane and the calculation is done by creating hexagonal lattice consisting of N rings and calculate the combined interaction for the central spin.

Arguments

  • a: Spin-spin distance.
  • N: Number of included spins.
CollectiveSpins.effective_interaction.squarelattice_orthogonalMethod
effective_interaction.squarelattice_orthogonal(a, N)

Effective Omega and Gamma for a infinite square lattice.

The polarization axis is orthogonal to the square lattice plane and the calculation is done by creating a (2N+1)*(2N+1) square lattice and calculate the combined interaction for the central spin.

Arguments

  • a: Spin-spin distance.
  • N: Number of included spins.
CollectiveSpins.effective_interaction.tetragonallattice_orthogonalMethod
effective_interaction.tetragonallattice_orthogonal(a, b, N)

Effective Omega and Gamma for a infinite tetragonal lattice.

The polarization axis is orthogonal to the top face of a unit cell and the calculation is done by creating a (2N+1)(2N+1)(2N+1) tetragonal lattice and calculate the combined interaction for the central spin.

Arguments

  • a: Spin-spin distance for bottom side square.
  • b: Height of the unit cell.
  • N: Number of included spins.
CollectiveSpins.reducedspin.JumpOperatorsFunction
JumpOperators(S::SpinCollectino, M::Int=1)

Gamma Matix and Jump Operators for dipole-dipole interactino.

Arguments

  • S: Spin collection.
  • M: Number of excitations.
CollectiveSpins.reducedspin.projector_idxMethod
projector_idx(b, j, m)

Find the indices of all basis states in the m excitation manifold involving the atom j in order to build a projector onto the upper state of this atom.

CollectiveSpins.reducedspin.reducedsigmapsigmamMethod
reducedsigmapsigmam(b::ReducedSpinBasis, i, j)

Create the operator product σᵢ⁺σⱼ⁻ directly on a ReducedSpinBasis. Useful especially for a basis where only one excitation is included, since then the single operators are zero (do not conserve excitation number), but the product is not.

CollectiveSpins.reducedspin.reducedspinstateMethod
reducedspinstate(b::ReducedSpinBasis, inds::Vector{Int})

State where the excitations are placed in the atoms given by inds. Note, that b.MS <= length(inds) <= b.M must be satisfied.

Examples

julia> b = CollectiveSpins.ReducedSpinBasis(4,2)
ReducedSpin(N=4, M=2, MS=0)

julia> GS = CollectiveSpins.reducedspinstate(b,[]) # get the ground state
Ket(dim=11)
  basis: ReducedSpin(N=4, M=2, MS=0)
 1.0 + 0.0im
 0.0 + 0.0im
 0.0 + 0.0im
 0.0 + 0.0im
 0.0 + 0.0im
 0.0 + 0.0im
 0.0 + 0.0im
 0.0 + 0.0im
 0.0 + 0.0im
 0.0 + 0.0im
 0.0 + 0.0im

julia> ψ2 = CollectiveSpins.reducedspinstate(b,[1,2]) # First and second atom excited
Ket(dim=11)
  basis: ReducedSpin(N=4, M=2, MS=0)
 0.0 + 0.0im
 0.0 + 0.0im
 0.0 + 0.0im
 0.0 + 0.0im
 0.0 + 0.0im
 1.0 + 0.0im
 0.0 + 0.0im
 0.0 + 0.0im
 0.0 + 0.0im
 0.0 + 0.0im
 0.0 + 0.0im
CollectiveSpins.reducedspin.reducedspintransitionMethod
reducedspintransition(b::ReducedSpinBasis, to::Vector{Int}, from::Vector{Int})

Transition operator $|\mathrm{to}⟩⟨\mathrm{from}|$, where to and from are given as vectors denoting the excitations.

CollectiveSpins.reducedspin.schroedinger_nhMethod
schroedinger_nh(T, S::SpinCollection, psi0; kwargs...)

Time evolution of the non-Hermitian Hamiltonian in the single excitation manifold (M=1) via the Schrödinger equation.
CollectiveSpins.reducedspin.transition_idxMethod
transition_idx(b, j, m)

Find the indices of all basis states where a transition from the m-1 excitation manifold into the m excitation manifold can occur by raising atom j

CollectiveSpins.effective_interaction_rotated.chain_orthogonalMethod
effective_interaction_rotated.chain_orthogonal(a, N, dϕ)

Effective Omega and Gamma for an infinite chain.

The polarization axis is orthogonal to the chain and the calculation is done by adding N spins left and N spins right of a central spin.

Arguments

  • a: Spin-spin distance.
  • N: Number of included spins.
  • : Phase shift between neighboring spins.
CollectiveSpins.quantum.JumpOperators_diagonalMethod
quantum.JumpOperators_diagonal(S)

Jump operators of the given system. (diagonalized)

Diagonalized means that the Gamma matrix is diagonalized and the jump operators are changed accordingly.

CollectiveSpins.quantum.blochstateMethod
quantum.blochstate(phi, theta[, N=1])

Product state of N single spin Bloch states.

All spins have the same azimuthal angle phi and polar angle theta.

CollectiveSpins.quantum.rotateMethod
meanfield.rotate(axis, angles, state)

Rotations on the Bloch sphere for the given density operator.

Arguments

  • axis: Rotation axis.
  • angles: Rotation angle(s).
  • ρ: Density operator that should be rotated.
CollectiveSpins.quantum.squeezeMethod
quantum.squeeze(axis, χT, ρ₀)

Spin squeezing along an arbitrary axis.

Arguments

  • axis: Squeezing axis.
  • χT: Squeezing strength.
  • ρ₀: Operator that should be squeezed.
CollectiveSpins.quantum.squeeze_sxMethod
quantum.squeeze_sx(χT, ρ₀)

Spin squeezing along sx.

Arguments

  • χT: Squeezing strength.
  • ρ₀: Operator that should be squeezed.
CollectiveSpins.quantum.timeevolutionMethod
quantum.timeevolution(T, S, state0[; fout])

Master equation time evolution.

Diagonalized means that the Gamma matrix is diagonalized and the jump operators are changed accordingly.

Arguments

  • T: Points of time for which output will be generated.
  • S: System
  • ρ₀: Initial density operator.
  • fout (optional): Function with signature fout(t, state) that is called whenever output should be generated.
CollectiveSpins.quantum.timeevolution_diagonalMethod
quantum.timeevolution_diagonal(T, S, state0[; fout])

Master equation time evolution. (diagonalized)

Diagonalized means that the Gamma matrix is diagonalized and the jump operators are changed accordingly.

Arguments

  • T: Points of time for which output will be generated.
  • S: System
  • ρ₀: Initial density operator.
  • fout (optional): Function with signature fout(t, state) that is called whenever output should be generated.
CollectiveSpins.collective_modes.Gamma_k_2DMethod
Gamma_k_2D(k_vec, a_vec1, a_vec2, polarization)

Collective decay rate Gammak of in-plane mode kvec for an 2D array of atoms in yz-plane.

WLOG, this calculation scales natural atomic frequency wavelength lambda0=1 and decay rate gamma0=1.

Arguments

  • k_vec: yz-axis quasimomentum of collective mode in first BZ such that |k|<= pi/a.
  • a_vec1, a_vec2: 2D Bravais lattice vectors.
  • polarization: 3D, complex vector of atomic polarization.
CollectiveSpins.collective_modes.Gamma_k_chainMethod
Gamma_k_chain(k, a, polarization)

Collective decay rate Gamma_k of mode k for an infinite chain of atoms along x-axis.

WLOG, this calculation scales natural atomic frequency wavelength lambda0=1 and decay rate gamma0=1.

Arguments

  • k: x-axis quasimomentum of collective mode in first BZ such that |k|<= pi/a.
  • a: Atomic lattice spacing.
  • polarization: 3D, complex vector of atomic polarization.
CollectiveSpins.collective_modes.Green_Tensor_k_2DMethod
Green_Tensor_k_2D(k_vec, a_vec1, a_vec2)

Green's Tensor in reciprocal space for in-plane mode k_vec for a 2D array of atoms in xy-plane. WLOG, this calculation scales natural atomic frequency wavelength lambda0=1 and decay rate gamma0=1.

Arguments

  • k_vec: xy-axis quasimomentum of collective mode in first BZ
  • a_vec1, a_vec2: 2D Bravais lattice vectors (real-space).
CollectiveSpins.collective_modes.N_lightconeMethod
N_lightcone(rec, rec_mag)

INTERNAL FUNCTION. Number of reciprocal lattice vectors that are within the lightcone.

Arguments

  • rec: Magnitude of reciprocal lattice vector overlap along axis.
  • rec_mag: Magnitude of reciprocal lattice vector.
CollectiveSpins.collective_modes.Omega_k_2DMethod
Omega_k_2D(k_vec, a_vec1, a_vec2, polarization)

Collective frequency shift Omegak of in-plane mode kvec for an 2D array of atoms in yz-plane.

WLOG, this calculation scales natural atomic frequency wavelength lambda0=1 and decay rate gamma0=1.

Arguments

  • k_vec: yz-axis quasimomentum of collective mode in first BZ such that |k|<= pi/a.
  • a_vec1, a_vec2: 2D Bravais lattice vectors.
  • polarization: 3D, complex vector of atomic polarization.
  • Lambda: Cutoff frequency of renormalization.
  • N1: Number of terms in a_vec1 reciprocal lattice direction.
  • N2: Number of terms in a_vec2 reciprocal lattice direction.
CollectiveSpins.collective_modes.Omega_k_chainMethod
Omega_k_chain(k, a, polarization)

Collective frequency shift Omega_k of mode k for an infinite chain of atoms along x-axis.

WLOG, this calculation scales natural atomic frequency wavelength lambda0=1 and decay rate gamma0=1.

Arguments

  • k: x-axis quasimomentum of collective mode in first BZ such that |k|<= pi/a.
  • a: Atomic lattice spacing.
  • polarization: 3D, complex vector of atomic polarization.
CollectiveSpins.collective_modes.polarization_renormMethod
polarization_renorm(polarization)

INTERNAL FUNCTION. Normalizes polarization vector and separates in and out-of-lattice-plane components.

Arguments

  • polarization: 3D, complex vector of atomic polarization.
CollectiveSpins.collective_modes.renorm_constsMethod
renorm_consts(V, b1, b2, Lambda, N1, N2)

INTERNAL FUNCTION. Returns default cutoff frequency and number of reciprocal lattice terms to be summed if keyword arguments not specified.

Arguments

  • V: Volume of BZ.
  • b1, b2: Reciprocal lattice vectors from avec1 and avec2.
  • Lambda: Cutoff frequency of renormalization.
  • N1: Number of terms in a_vec1 reciprocal lattice direction.
  • N2: Number of terms in a_vec2 reciprocal lattice direction.
CollectiveSpins.independent.blochstateMethod
independent.blochstate(phi, theta[, N=1])

Product state of N single spin Bloch states.

All spins have the same azimuthal angle phi and polar angle theta.

CollectiveSpins.independent.timeevolutionMethod
independent.timeevolution(T, gamma, state0)

Independent time evolution.

Arguments

  • T: Points of time for which output will be generated.
  • gamma: Decay rate(s).
  • state0: Initial state.
CollectiveSpins.independent.timeevolutionMethod
independent.timeevolution(T, S::SpinCollection, state0)

Independent time evolution.

Arguments

  • T: Points of time for which output will be generated.
  • S: SpinCollection describing the system.
  • state0: Initial state.