Problem Setup

The parameters that determine the physics of the problem are held in a SchrodingerProb.

For a state-transfer problem, please give the initial conditions as column matrices. There are currently bugs when giving them as vectors.

QuantumGateDesign.SchrodingerProbType
SchrodingerProb(system_sym, system_asym, sym_operators, asym_operators, u0, v0, tf, nsteps, N_ess_levels)

Set up an object containing the data that defines the physics and numerics of the problem (the Hamiltonians, initial conditions, number of timesteps, etc).

Arguments

  • system_sym::M: the symmetric/real part of the system Hamiltonian.
  • system_asym::M: the antisymmetric/imaginary part of the system Hamiltonian.
  • sym_operators::Vector{M}: a vector whose i-th entry is the symmetric part of the i-th control Hamiltonian.
  • asym_operators::Vector{M}: a vector whose i-th entry is the antisymmetric part of the i-th control Hamiltonian.
  • u0::M: the real part of the initial conditions. The i-th column corresponds to the i-th initial state in the gate basis.
  • v0::M: the imaginary part of the initial conditions. The i-th column corresponds to the i-th initial state in the gate basis.
  • tf::Real: duration of the gate.
  • nsteps::Int64: number of timesteps to take.
  • N_ess_levels::Int64: number of levels in the 'essential' subspace, i.e. the part of the subspace actually used for computation.
  • guard_subspace_projector::Union{M, missing}=missing: matrix projecting a state vector in to the 'guard' subspace.

Keyword Arguments

  • gmres_abstol::Float64: absolute tolerance to use when solving linear systems with GMRES.
  • gmres_reltol::Float64: relative tolerance to use when solving linear systems with GMRES.
  • `preconditioner_type: the type of preconditioner to use in the algorithm.

where M <: AbstractMatrix{Float64}

source
QuantumGateDesign.guard_projectorFunction

Given the size of each subsystem and the number of essential levels in each subsystem, return a matrix which projects the state vector onto the guarded subspace.

Note that the first subsystem corresponds to the leftmost bit of the quantum bitstring.

E.g.

$|n_0 n_1 n_2 \rangle = |n_0\rangle \otimes |n_1\rangle \otimes |n_2\rangle$

Examples ≡≡≡≡≡≡≡≡

``` julia> guard_projector([3], [2]) 6×6 SparseMatrixCSC{Int64, Int64} with 6 stored entries: 0 ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ 0 ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ 1 ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ 0 ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ 0 ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ 1

julia> guard_projector([2,2], [2,1]) 8×8 SparseMatrixCSC{Int64, Int64} with 8 stored entries: 0 ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ 0 ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ 1 ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ 1 ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ 0 ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ 0 ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ 1 ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ 1 ```

source