The following types are exported and available by using Juqbox
.
Juqbox.Working_Arrays
— Typewa = Working_Arrays(params::objparams, nCoeff::Int64)
Constructor for the mutable struct Working_Arrays containing preallocated working arrays.
Arguments
param:: objparams
: Struct with problem definitionnCoeff:: Int64
: Number of parameters in optimization
Juqbox.bcparams
— Typebcpar = bcparams(T, D1, Ncoupled, Nunc, omega, pcof)
General constructor of struct bcparams
for setting up B-splines with carrier waves.
bcpar = bcparams(T, D1, omega, pcof)
Simplified constructor for the case when there are no uncoupled controls and Ncoupled = size(omega,1)
.
Arguments
T:: Float64
: Duration of spline functionD1:: Int64
: Number of basis functions in each segmentNcoupled::Int64
: Number of coupled controls in the simulationNunc::Int64
: Number of uncoupled controls in the simulationomega::Array{Float64,2}
: Carrier wave frequenciespcof:: Array{Float64, 1}
: Coefficient vector. Must have D1*Nseg elements
First dimensions of the omega
array:
- Without uncoupled controls,
Nunc=0
andsize(omega,1) = Ncoupled
. - With uncoupled controls,
Nunc > 0
andsize(omega,1) = Ncoupled + Nunc
.
Second dimension of the omega
array:
size(omega, 2) = Nfreq
Ordering of the pcof
array:
First consider the case without uncoupled control functions, Nunc = 0
: Then the pcof
array then has 2*Ncoupled*Nfreq*D1
elements. Each ctrl ∈ [1,Ncoupled]
and freq ∈ [1,Nfreq]
corresponds to D1
elements in the pcof
vector. For the case Ncoupled = 2
and Nfreq = 2
, the elements are ordered according to
ctrl | freq | α_1 | α_2 |
---|---|---|---|
1 | 1 | 1:D1 | D1+1:2 D1 |
1 | 2 | 2 D1+1: 3 D1 | 3 D1+1:4 D1 |
2 | 1 | 4 D1+1: 5 D1 | 5 D1+1:6 D1 |
2 | 2 | 6 D1+1: 7 D1 | 7 D1+1: 8D1 |
If there are uncoupled controls, Nunc > 0
, the pcof
array should have (2*Ncoupled + Nunc)*Nfreq*D1
elements. The last Nunc*Nfreq*D1
elements correspond to the uncoupled control functions and are ordered in a corresponding way.
External links
- Spline Wavelet on Wikipedia.
Juqbox.lsolver_object
— Typelinear_solver = lsolver_object(; tol = tol,
max_iter = max_iter,
nrhs = nrhs,
solver = NEUMANN_SOLVER)
Constructor for the mutable struct lsolver_object. That allcoates arrays and sets up the function pointers for the different linear solvers supported.
Arguments
tol::Float64 = 1e-10
: Convergence tolerance of the iterative solver (only needed for Jacobi)max_iter::Int64 = 3
: Max number of iterations for the linear solvernrhs::Int64 = 1
: Number of right-hand sides (used for tolerance scaling)solver::Int64 = NEUMANN_SOLVER
: (keyword) ID of the iterative solver. Can take the value of NEUMANNSOLVER (i.e. 1) or JACOBISOLVER (i.e. 2) See examples/cnot2-jacobi-setup.jl
Juqbox.objparams
— Typeparams = objparams(Ne, Ng, T, Nsteps;
Uinit=Uinit,
Utarget=Utarget,
Cfreq=Cfreq,
Rfreq=Rfreq,
Hconst=Hconst [,
Hsym_ops=Hsym_ops,
Hanti_ops=Hanti_ops,
Hunc_ops=Hunc_ops,
wmatScale=wmatScale,
objFuncType=objFuncType,
leak_ubound=leak_ubound,
linear_solver = lsolver_object(),
use_sparse = use_sparse],
dVds = dVds)
Constructor for the mutable struct objparams. The sizes of the arrays in the argument list are based on Ntot = prod(Ne + Ng)
, Ness = prod(Ne)
, Nosc = length(Ne) = length(Ng)
.
Notes: It is assumed that length(Hsym_ops) = length(Hanti_ops) =: Ncoupled
. The matrices Hconst
, Hsym_ops[j]
and Hanti_ops[j]
, for j∈[1,Ncoupled], must all be of size Ntot × Ntot
. The matrices Hsym_ops[j]
must be symmetric and Hanti_ops[j]
must be skew-symmetric. The matrices Hunc_ops[j]
, for j∈[1,Nunc], where Nunc = length(Hunc_ops)
, must also be of size Ntot × Ntot
and either be symmetric or skew-symmetric.
Arguments
Ne::Array{Int64,1}
: Number of essential energy levels for each subsystemNg::Array{Int64,1}
: Number of guard energy levels for each subsystemT::Float64
: Duration of gateNsteps::Int64
: Number of timesteps for integrating Schroedinger's equationUinit::Array{Float64,2}
: (keyword) Matrix holding the initial conditions for the solution matrix of size Uinit[Ntot, Ness]Utarget::Array{Complex{Float64},2}
: (keyword) Matrix holding the target gate matrix of size Uinit[Ntot, Ness]Cfreq::Array{Float64,2}
: (keyword) Carrier wave (angular) frequencies of size Cfreq[Nctrl, Nfreq]Rfreq::Array{Float64,1}
: (keyword) Rotational (regular) frequencies for each control Hamiltonian; size Rfreq[Nctrl]Hconst::Array{Float64,2}
: (keyword) Time-independent part of the Hamiltonian matrix of size Ntot × NtotHsym_ops:: Array{Array{Float64,2},1}
: (keyword) Array of symmetric control Hamiltonians, each of size Ntot × NtotHanti_ops:: Array{Array{Float64,2},1}
: (keyword) Array of anti-symmetric control Hamiltonians, each of size Ntot × NtotHunc_ops:: Array{Array{Float64,2},1}
: (keyword) Array of uncoupled control Hamiltonians, each of size Ntot × NtotwmatScale::Float64 = 1.0
: (keyword) Scaling factor for suppressing guarded energy levelsobjFuncType::Int64 = 1
# 1 = objective function include infidelity and leakage # 2 = objective function only includes infidelity... no leakage in obj function or constraint # 3 = objective function only includes infidelity; leakage treated as inequality constraintleak_ubound::Float64 = 1.0e-3
: The upper bound on the leakage inequality constraint (See examples/cnot2-leakieq-setup.jl )linear_solver::lsolver_object = lsolver_object()
: The linear solver object used to solve the implicit & adjoint systemuse_sparse::Bool = false
: (keyword) Set to true to sparsify all Hamiltonian matricesdVds::Array{Complex{Float64},2}
: (keyword) Matrix holding the complex-valued matrix dV/ds of size Ntot x Ne (for continuation)
Juqbox.splineparams
— Typespar = splineparams(T, D1, Nseg, pcof)
Constructor for struct splineparams, which sets up the parameters for a regular B-spline function (without carrier waves).
Arguments
T:: Float64
: Duration of spline functionD1:: Int64
: Number of basis functions in each splineNseg:: Int64
: Number of splines (real, imaginary, different ctrl func)pcof:: Array{Float64, 1}
: Coefficient vector. Must have D1*Nseg elements
External links
- Spline Wavelet on Wikipedia.