APIs
Spin Shuttling Models
SpinShuttling.ShuttlingModel
— TypeSpin shuttling model defined by a stochastic field, the realization of the stochastic field is specified by the paths of the shuttled spins.
Arguments
n::Int
: Number of spinsΨ::Vector{<:Complex}
: Initial state of the spin system, the length of the vector must be `2^nT::Real
: Maximum timeN::Int
: Time discretizationB::GaussianRandomField
: Noise fieldX::Vector{<:Function}
: Shuttling paths, the length of the vector must be `ninitialize::Bool
: Initialize the random function
Fields
R::RandomFunction
: Random function defined by the noise field and the paths
Example
model = ShuttlingModel(1, [1+0im, 0+0im], 1.0, 100, OrnsteinUhlenbeckField([1.0, 1.0, 1.0]), [t->t])
SpinShuttling.OneSpinModel
— FunctionGeneral one spin shuttling model initialized at initial state |Ψ₀⟩, with arbitrary shuttling path x(t).
Arguments
Ψ::Vector{<:Complex}
: Initial state of the spin system, the length of the vector must be `2^nT::Real
: Maximum timeN::Int
: Time discretizationB::GaussianRandomField
: Noise fieldx::Function
: Shuttling pathinitialize::Bool
: Initialize the random function
Example
model = OneSpinModel(1 / √2 * [1+0im, 1+0im], 1.0, 100, OrnsteinUhlenbeckField([1.0, 1.0, 1.0]), t->t)
One spin shuttling model initialzied at |Ψ₀⟩=|+⟩. The qubit is shuttled at constant velocity along the path x(t)=L/T*t
, with total time T
in μs
and length L
in μm
.
Arguments
T::Real
: Maximum timeL::Real
: Length of the pathN::Int
: Time discretizationB::GaussianRandomField
: Noise fieldv::Real
: Velocity of the shuttling
Example
model = OneSpinModel(1.0, 1.0, 100, OrnsteinUhlenbeckField([1.0, 1.0, 1.0]), 1.0)
SpinShuttling.OneSpinForthBackModel
— FunctionOne spin shuttling model initialzied at |Ψ₀⟩=|+⟩. The qubit is shuttled at constant velocity along a forth-back path x(t, T, L) = t<T/2 ? 2L/T*t : 2L/T*(T-t)
, with total time T
in μs
and length L
in μm
.
Arguments
T::Real
: Maximum timeL::Real
: Length of the pathN::Int
: Time discretizationB::GaussianRandomField
: Noise fieldv::Real
: Velocity of the shuttling
Example
model = OneSpinForthBackModel(1.0, 1.0, 100, OrnsteinUhlenbeckField([1.0, 1.0, 1.0]))
SpinShuttling.TwoSpinModel
— FunctionGeneral two spin shuttling model initialized at initial state |Ψ₀⟩, with arbitrary shuttling paths x₁(t), x₂(t).
Arguments
Ψ::Vector{<:Complex}
: Initial state of the spin system, the length of the vector must be `2^nT::Real
: Maximum timeN::Int
: Time discretizationB::GaussianRandomField
: Noise fieldx₁::Function
: Shuttling path for the first spinx₂::Function
: Shuttling path for the second spin
Example
model = TwoSpinModel(1 / √2 * [1+0im, 1+0im, 1+0im, 1+0im], 1.0, 100, OrnsteinUhlenbeckField([1.0, 1.0, 1.0]), t->t, t->t+1)
SpinShuttling.TwoSpinSequentialModel
— FunctionTwo spin shuttling model initialized at the singlet state |Ψ₀⟩=1/√2(|↑↓⟩-|↓↑⟩)
. The qubits are shuttled at constant velocity along the path x₁(t)=L/T₁*t
and x₂(t)=L/T₁*(t-T₀)
. The delay between the them is T₀
and the total shuttling time is T₁+T₀
. It should be noticed that due to the exclusion of fermions, x₁(t)
and x₂(t)
cannot overlap.
Arguments
T₀::Real
: Delay timeT₁::Real
: Shuttling timeL::Real
: Length of the pathN::Int
: Time discretizationB::GaussianRandomField
: Noise field
Example
model = TwoSpinSequentialModel(1.0, 1.0, 1.0, 100, OrnsteinUhlenbeckField([1.0, 1.0, 1.0]))
SpinShuttling.TwoSpinParallelModel
— FunctionTwo spin shuttling model initialized at the singlet state |Ψ₀⟩=1/√2(|↑↓⟩-|↓↑⟩)
. The qubits are shuttled at constant velocity along the 2D path x₁(t)=L/T*t, y₁(t)=0
and x₂(t)=L/T*t, y₂(t)=D
. The total shuttling time is T
and the length of the path is L
in μm
.
Arguments
T::Real
: Total timeD::Real
: Distance between the two qubitsL::Real
: Length of the pathN::Int
: Time discretizationB::GaussianRandomField
: Noise field
Example
model = TwoSpinParallelModel(1.0, 1.0, 1.0, 100, OrnsteinUhlenbeckField([1.0, 1.0, 1.0]))
SpinShuttling.dephasingmatrix
— FunctionCalculate the dephasing matrix of a given spin shuttling model.
Arguments
model::ShuttlingModel
: The spin shuttling modelmethod::Symbol
: The method to calculate the dephasing factor,:trapezoid
,:simpson
,:adaptive
Returns
W::Matrix{<:Real}
: The dephasing matrix
Example
model = OneSpinModel(1.0, 1.0, 100, OrnsteinUhlenbeckField([1.0, 1.0, 1.0]), t->t)
W = dephasingmatrix(model)
Sample the dephasing matrix array for a given normal random vector.
Arguments
model::ShuttlingModel
: The spin shuttling modelrandseq::Vector{<:Real}
: The random sequenceisarray::Bool
: Return the dephasing matrix array for each time step
SpinShuttling.sampling
— FunctionSampling an observable that defines on a specific spin shuttling model
Arguments
model::ShuttlingModel
: The spin shuttling modelobjective::Function
: The objective functionobjective(mode::ShuttlingModel; randseq)
`M::Int
: Monte-Carlo sampling size
SpinShuttling.statefidelity
— FunctionCalculate the state fidelity of a spin shuttling model using numerical integration of the covariance matrix.
Arguments
model::ShuttlingModel
: The spin shuttling modelmethod::Symbol
: The method to calculate the dephasing factor,:trapezoid
,:simpson
,:adaptive
Returns
Real
: The state fidelity
Sample the state fidelity of a spin shuttling model using Monte-Carlo sampling.
Arguments
model::ShuttlingModel
: The spin shuttling modelrandseq::Vector{<:Real}
: The random sequenceisarray::Bool
: Return the dephasing matrix array for each time step
SpinShuttling.W
— FunctionAnalytical dephasing factor of a one-spin shuttling model.
Arguments
T::Real
: Total timeL::Real
: Length of the pathB<:GaussianRandomField
: Noise field, Ornstein-Uhlenbeck or Pink-Lorentzianpath::Symbol
: Path of the shuttling model,:straight
or:forthback
Analytical dephasing factor of a sequenced two-spin EPR pair shuttling model.
Stochastics
SpinShuttling.OrnsteinUhlenbeckField
— TypeOrnstein-Uhlenbeck field, the correlation function of which is σ^2 * exp(-|t₁ - t₂|/θₜ) * exp(-|x₁-x₂|/θₓ)
where t
is time and x
is position.
SpinShuttling.PinkLorentzianField
— TypePink-Lorentzian Field, the correlation function of which is σ^2 * (expinti(-γ[2]abs(t₁ - t₂)) - expinti(-γ[1]abs(t₁ - t₂)))/log(γ[2]/γ[1]) * exp(-κ|x₁-x₂|)
where expinti
is the exponential integral function.
SpinShuttling.PinkPiField
— TypePink-HeavisidePi Field, the correlation function of which is σ^2 * (expinti(-γ[2]abs(t₁ - t₂)) - expinti(-γ[1]abs(t₁ - t₂)))/log(γ[2]/γ[1]) * sinc(-κ|x₁-x₂|)
where expinti
is the exponential integral function.
SpinShuttling.RandomFunction
— TypeGenerate a random time series from a Gaussian random field.
R()
generates a random time series from a Gaussian random field R
R(randseq)
generates a random time series from a Gaussian random field R
with a given random sequence randseq
.
SpinShuttling.CompositeRandomFunction
— FunctionCreate a new random function composed by a linear combination of random processes. The input random function represents the direct sum of these processes. The output random function is a tensor contraction from the input.
Arguments
R::RandomFunction
: a direct sum of random processes R₁⊕ R₂⊕ ... ⊕ Rₙc::Vector{Int}
: a vector of coefficientsinitialize::Bool=true
: whether to initialize the Cholesky decomposition of the covariance matrix
Returns
RandomFunction
: a new random function composed by a linear combination of random processes
SpinShuttling.restriction
— FunctionRestrict the random field along a parameteized curve.
Arguments
X::Vector{<:Function}
: a vector of functions [x₁(t), x₂(t), ...]t::AbstractArray
: time arrayB::GaussianRandomField
: a Gaussian random fieldinitialize::Bool
: initialize the random function
Returns
RandomFunction
: a new random function restricted along the curve
SpinShuttling.initialize!
— FunctionInitialize the Cholesky decomposition of the covariance matrix of a random function.
SpinShuttling.characteristicfunction
— FunctionCompute the characteristic functional of the process from the numerical quadrature of the covariance matrix. Using Simpson's rule by default.
SpinShuttling.characteristicvalue
— FunctionCompute the final phase of the characteristic functional of the process from the numerical quadrature of the covariance matrix. Using Simpson's rule by default.
SpinShuttling.covariancematrix
— FunctionCovariance matrix of a Gaussian random field. When P₁=P₂
, it is the auto-covariance matrix of a Gaussian random process. When P₁!=P₂
, it is the cross-covariance matrix between two Gaussian random processes.
Arguments
P₁::Vector{<:Point}
: time-position arrayP₂::Vector{<:Point}
: time-position arrayprocess::GaussianRandomField
: a Gaussian random field
Auto-Covariance matrix of a Gaussian random process.
Arguments
P::Vector{<:Point}
: time-position arrayprocess::GaussianRandomField
: a Gaussian random field
Returns
Symmetric{Real}
: auto-covariance matrix
SpinShuttling.covariance
— FunctionCovariance function of Gaussian random field.
Arguments
p₁::Point
: time-position arrayp₂::Point
: time-position arrayprocess<:GaussianRandomField
: a Gaussian random field, e.g.OrnsteinUhlenbeckField
orPinkLorentzianField