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::GaussianRandomFunction: 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
Ψ::Vector{<:Complex}: Initial state of the spin system, the length of the vector must be4, default is1/√2(|↑↓⟩-|↓↑⟩)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
Ψ::Vector{<:Complex}: Initial state of the spin system, the length of the vector must be4, default is1/√2(|↑↓⟩-|↓↑⟩)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.dephasingfactor — FunctionCalculate the dephasingfactor according to a special combinator of the noise sequence.
Arguments
model::ShuttlingModel: The spin shuttling modelc::Vector{Int}: The combinator of the noise sequence, which should have the same length as the number of spins.method::Symbol: The method to calculate the characteristic value,:trapezoid,:simpson,:adaptive
Returns
Real: The dephasing factor
Example
model = OneSpinModel(1.0, 1.0, 100, OrnsteinUhlenbeckField([1.0, 1.0, 1.0]), t->t)
c = [1, 1, 1]
dephasingfactor(model, c)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
Quantum Process
SpinShuttling.MixingUnitaryChannel — TypeMixingUnitaryChannel <: QuantumProcess A quantum process that applies a set of unitary operations with given probabilities.
Fields
Us::Vector{Matrix{Complex{Float64}}}: Vector of unitary matrices, each of size2^n × 2^nfornqubits.Ps::Vector{Float64}: Vector of mixing probabilities, which sum to 1
SpinShuttling.processtomography — Functionprocesstomography(KrausOps::Vector{<:Matrix{<:Complex}}, basis::Vector{<:Matrix{<:Complex}}; normalized::Bool=false) Compute the transfer matrix for a quantum channel defined by a set of Kraus operators.
Arguments
KrausOps::Vector{<:Matrix{<:Complex}}: Vector of Kraus operators, each of size2^n × 2^nfornqubits.basis::Vector{<:Matrix{<:Complex}}: Vector of basis matrices, each of size2^n × 2^nfornqubits.normalized::Bool=false: Iftrue, the transfer matrix is normalized by the number of Kraus operators.
Returns
Matrix{Float64}: The process tomography transfer matrix of size `N × N
SpinShuttling.paulitransfermatrix — Functionpaulitransfermatrix(KrausOps::Vector{<:Matrix{<:Complex}}; normalized::Bool=false) Compute the Pauli transfer matrix for a set of Kraus operators.
Arguments
KrausOps::Vector{<:Matrix{<:Complex}}: Vector of Kraus operators, each of size2^n × 2^nfornqubits.normalized::Bool=false: Iftrue, the transfer matrix is normalized by the number of Kraus operators.
Returns
Matrix{Float64}: The Pauli transfer matrix of size4^n × 4^n.
SpinShuttling.KrausOps — FunctionKrausOps(channel::MixingUnitaryChannel)::Vector{Matrix{Complex{Float64}}} Convert a MixingUnitaryChannel to its Kraus operators.
Arguments
channel::MixingUnitaryChannel: The mixing unitary channel.
Returns
Vector{Matrix{Complex{Float64}}}: Vector of Kraus operators,
Quantum Measures
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.concurrence — Functionconcurrence(ρ) -> Float64Compute the Wootters concurrence of a two‑qubit (4×4) density matrix ρ.
Arguments
ρ::AbstractMatrix{<:Complex}: Hermitian, positive‑semidefinite matrix withsize == (4,4)andtr(ρ) ≈ 1.
Returns
Float64in the closed interval [0,1].
Algorithm
- Build the spin‑flip operator
Y ⊗ Y(Kronecker product of Pauli‑Y). - Form spin-flip matrix ρ̃.
- Compute the eigenvalues
λ_iofρ * ρ̃, take their square‑roots, sort descending. - Return
max(0, λ₁ − λ₂ − λ₃ − λ₄).
Reference
W.K. Wootters, Phys. Rev. Lett. 80, 2245 (1998).
SpinShuttling.vonneumannentropy — Functionvonneumannentropy(ρ::AbstractMatrix{<:Number})::Float64Compute the von Neumann entropy of a quantum state represented by a density matrix ρ.
Arguments
ρ::AbstractMatrix{<:Number}: Hermitian, positive-semidefinite matrix representing the quantum state.
Returns
Float64: The von Neumann entropy, a non-negative real number.
SpinShuttling.processfidelity — Functionreturns the process fidelity of a quantum channel defined by a transfer matrix Λ and a target process S.
Arguments
Λ::Matrix{<:Number}: Transfer matrix of the quantum channel.S::Matrix{<:Number}: Target process matrix.d::Int=0: Dimension of the Hilbert space. Ifd
is 0, it is inferred from the size of Λ.
Returns
Float64: The process fidelity, a real number in the closed interval **[0
processfidelity(E::MixingUnitaryChannel, S::Matrix{<:Number}; d::Int=0) Compute the process fidelity of a quantum channel defined by a MixingUnitaryChannel and a target process S.
Arguments
E::MixingUnitaryChannel: The mixing unitary channel.S::Matrix{<:Number}: The target process matrix.d::Int=0: Dimension of the Hilbert space. Ifdis0, it is inferred from the size ofE.
Returns
Float64: The process fidelity, a real number in the closed interval **[0
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,:straightor: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.GaussianRandomFunction — 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.CompositeGaussianRandomFunction — 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::GaussianRandomFunction: 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
GaussianRandomFunction: 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
GaussianRandomFunction: 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.OrnsteinUhlenbeckFieldorPinkLorentzianField