API
QuantumPropagators
The highest-level API of the QuantumPropagators.jl
package (apart from some convenience functions) consists of a single function:
propagate
— Propagate a quantum state over an entire time grid under a given generator
At a slightly lower level, propagation of quantum states in encapsulated by The Propagator interface:
init_prop
— Initialize apropagator
object, which is of some concrete (method-dependent) sub-type ofAbstractPropagator
reinit_prop!
— Re-initialize thepropagator
prop_step!
— Advance thepropagator
by a single time step forward or backward
In some cases, the ability to mutate the propagator after each propagation step can be useful. This can be achieved with the following private methods:
set_state!
— Mutate the current quantumstate
of thepropagator
(not exported)set_t!
— Mutate the current time of thepropagator
(not exported)
The dynamics of a quantum state are determined by a time-dependent dynamical generator (a Hamiltonian or Liouvillian). The QuantumPropagators
package re-exports the two main initialization routines for generators from QuantumPropagators.Generators
:
hamiltonian
— Construct a time-dependent generator for a propagation in Hilbert space under the Schrödinger equationliouvillian
— Construct a time-dependent generator for a propagation in Liouville space under the master equation in Lindblad form
To set up the time-dependent control fields in a Hamiltonian, methods from the submodules QuantumPropagators.Controls
, QuantumPropagators.Shapes
, and QuantumPropagators.Amplitudes
can be used.
The above constitutes the main interface of QuantumPropagators
. At the lowest level, further functionality is provided by sub-modules like QuantumPropagators.Cheby
, which defines a standalone API specifically for the Chebychev propagation method. The full list of sub-modules and their public members is:
QuantumPropagators.Amplitudes
QuantumPropagators.Arnoldi
QuantumPropagators.Cheby
QuantumPropagators.Controls
QuantumPropagators.ExpProp
QuantumPropagators.Generators
QuantumPropagators.Interfaces
QuantumPropagators.Newton
QuantumPropagators.Shapes
QuantumPropagators.SpectralRange
QuantumPropagators.Storage
Convenience functions
There are some "convenience functions" that wrap around propagate
for common tasks:
Reference for top-level QuantumPropagators
module
Public Members:
Re-exported Members:
Private Members:
cheby_get_spectral_envelope
ExpPropagator
PWCPropagator
ChebyPropagator
AbstractPropagator
enable_timings
set_state!
set_t!
timings_enabled
PiecewisePropagator
NewtonPropagator
disable_timings
ode_function
Public members
QuantumPropagators.Propagation
— TypeWrapper around the parameters of a call to propagate
.
Propagation(
generator, tlist;
pre_propagation=nothing, post_propagation=nothing,
kwargs...
)
Propagation(
propagator;
pre_propagation=nothing, post_propagation=nothing,
kwargs...
)
is a wrapper around the arguments for propagate
/ init_prop
, for use within propagate_sequence
.
The positional and keyword arguments are those accepted by the above propagation routines, excluding the initial state. A Propagation
may in addition include the pre_propagation
and post_propagation
keyword arguments recognized by propagate_sequence
.
QuantumPropagators.init_prop
— FunctionInitialize a Propagator
.
propagator = init_prop(
state, generator, tlist;
method, # mandatory keyword argument
backward=false,
inplace=QuantumPropagators.Interfaces.supports_inplace(state),
piecewise=nothing,
pwc=nothing,
kwargs...
)
initializes a propagator for the time propagation of the given state
over a time grid tlist
under the time-dependent generator (Hamiltonian/Liouvillian) generator
.
Arguments
state
: The "initial" state for the propagation. Forbackward=false
, this state is taken to be at initial time (tlist[begin]
); and forbackward=true
, at the final time (tlist[end]
)generator
: The time-dependent generator of the dynamicstlist
: The time grid over which which the propagation is defined. This may or may not be equidistant.
Mandatory keyword arguments
method
: The propagation method to use. May be given as a name (Symbol
), but the recommended usage is to pass a module implementing the propagation method, e.g.,using QuantumPropagators: Cheby; method = Cheby
. Passing a module ensures that the code implementing the method is correctly loaded. This is particularly important for propagators using third-party backends, like withmethod=OrdinaryDiffEq
.
Optional keyword arguments
backward
: Iftrue
, initialize the propagator for a backward propagation. The resultingpropagator.t
will betlist[end]
, and subsequent calls toprop_step!
will move backward ontlist
.inplace
: Iftrue
, thestate
property of the resulting propagator will be changed in-place by any call toprop_step!
. Iffalse
, each call toprop_step!
changes the reference forpropagator.state
, and the propagation will not use any in-place operations. Not all propagation methods may support both in-place and not-in-place propagation. In-place propagation is generally more efficient for larger Hilbert space dimensions, but may not be compatible, e.g., with automatic differentiation.piecewise
: If given as a boolean,true
enforces that the resulting propagator is aPiecewisePropagator
, andfalse
enforces that it not aPiecewisePropagator
. For the defaultpiecewise=nothing
, whatever type of propagation is the default for the givenmethod
will be used. Throw an error if the givenmethod
does not support the required type of propagation.pwc
: Likepiecewise
, but for the strongerPWCPropagator
.
All other kwargs
are method-dependent and are ignored for methods that do not support them.
The type of the returned propagator
is a sub-type of AbstractPropagator
, respectively a sub-type of PiecewisePropagator
if piecewise=true
or a sub-type of PWCPropagator
if pwc=true
.
Internals
Internally, the (mandatory) keyword method
is converted into a fourth positional argument. This allows propagation methods to define their own implementation of init_prop
via multiple dispatch. However, when calling init_prop
in high-level code, method
must always be given as a keyword argument.
See also
reinit_prop!
— Re-initialize a propagatorpropagate
— Higher-level propagation interface- `check_propagator — a function to verify the interface described above.
using QuantumPropagators: Cheby
cheby_propagator = init_prop(
state,
generator,
tlist;
method=Cheby,
inplace=QuantumPropagators.Interfaces.supports_inplace(state),
backward=false,
verbose=false,
parameters=nothing,
control_ranges=nothing,
specrange_method=:auto,
specrange_buffer=0.01,
cheby_coeffs_limit=1e-12,
check_normalization=false,
specrange_kwargs...
)
initializes a ChebyPropagator
.
Method-specific keyword arguments
control_ranges
: a dict the maps the controls ingenerator
(seeget_controls
) to a tuple of min/max values. The Chebychev coefficients will be calculated based on a spectral envelope that assumes that each control can take arbitrary values within the min/max range. If not given, the ranges are determined automatically. Specifying manual control ranges can be useful when the the control amplitudes (parameters
) may change during the propagation, e.g. in a sequential-update control scheme.specrange_method
: Method to pass to thespecrange
functionspecrange_buffer
: An additional factor by which to enlarge the estimated spectral range returned byspecrange
, in order to ensure that Chebychev coefficients are based on an overestimation of the spectral range.cheby_coeffs_limit
: The maximum magnitude of Chebychev coefficients that should be treated as non-zerocheck_normalization
: Check whether the Hamiltonian has been properly normalized, i.e., that the spectral range ofgenerator
has not been underestimated. This slowes down the propagation, but is advisable for novelgenerators
.uniform_dt_tolerance=1e-12
: How much the intervals oftlist
are allowed to vary while still being considered constant.specrange_kwargs
: All further keyword arguments are passed to thespecrange
function. Most notably, with the defaultspecrange_method=:auto
(orspecrange_method=:manual
), passingE_min
andE_max
allows to manually specify the spectral range ofgenerator
.
using QuantumPropagators: Newton
newton_propagator = init_prop(
state,
generator,
tlist;
method=Newton,
inplace=QuantumPropagators.Interfaces.supports_inplace(state),
backward=false,
verbose=false,
parameters=nothing,
m_max=10,
func=(z -> exp(-1im * z)),
norm_min=1e-14,
relerr=1e-12,
max_restarts=50,
_...
)
initializes a NewtonPropagator
.
Method-specific keyword arguments
using QuantumPropagators: ExpProp
exp_propagator = init_prop(
state,
generator,
tlist;
method=ExpProp,
inplace=QuantumPropagators.Interfaces.supports_inplace(state),
backward=false,
verbose=false,
parameters=nothing,
func=(H_dt -> exp(-1im * H_dt))
convert_state=_exp_prop_convert_state(state),
convert_operator=_exp_prop_convert_operator(generator),
_...
)
initializes an ExpPropagator
.
Method-specific keyword arguments
func
: The function to evaluate. The argumentH_dt
is obtained by constructing an operatorH
fromgenerator
via theevaluate
function and the multiplied with the time stepdt
for the current time interval. The propagation then simply multiplies the return value offunc
with the current stateconvert_state
: Type to which to temporarily convert states before multiplying the return value offunc
.convert_operator
: Type to which to convert the operatorH
before multiplying it withdt
and plugging the result intofunc
The convert_state
and convert_operator
parameters are useful for when the generator
and or state
are unusual data structures for which the relevant methods to calculate func
are not defined. Often, it is easier to temporarily convert them to standard complex matrices and vectors than to implement the missing methods.
using OrdinaryDiffEq # or: `using DifferentialEquations`
ode_propagator = init_prop(
state,
generator,
tlist;
method=OrdinaryDiffEq, # or: `method=DifferentialEquations`
inplace=QuantumPropagators.Interfaces.supports_inplace(state),
backward=false,
verbose=false,
parameters=nothing,
piecewise=false,
pwc=false,
alg=OrdinaryDiffEq.Tsit5(),
solver_options...
)
initializes a propagator that uses an ODE solver from the OrdinaryDiffEq.jl package as a backend.
By default, the resulting propagator is for time-continuous controls that can be evaluated with evaluate(control, t)
for any t
in the range of tlist[begin]
to tlist[end]
. The controls may be parametrized, see get_parameters
. Any parameters will be available in the parameters
attribute of the resulting ode_propagator
, as a dictionary mapping controls to a vector of parameter values. Mutating ode_propagator.parameters[control]
will be reflected in any subsequent call to prop_step!
.
If pwc=true
(or, equivalently piecewise=true
), all controls will be discretized with discretize_on_midpoints
and the propagation will be for piecewise constant dynamics. The resulting ode_propagator
will be an instance of PWCPropagator
, with the corresponding semantics. In particular, the ode_propatator.parameters
will be a mapping of controls to discretized pulse values, not the analytical parameters obtained with get_parameters(control)
as in the default case.
Internally, the generator
will be wrapped with QuantumPropagators.ode_function
. The resulting function f
will be called internally as f(ϕ, Ψ, vals_dict, t)
or f(Ψ, vals_dict, t)
depending on the inplace
keyword argument.
Method-specific keyword arguments
pwc
: Whether to propagate for piecewise-constant controls or, with the defaultpwc=false
, for time-continuous controls.piecewise
: Currently equivalent topwc
, but future version may change this to allow for other piecewise (e.g., piecewise-linear) controls.parameters
: If given, a mapping of controls to parameter values (pwc=false
) or pulse values on the intervals of the time grid (pwc=true
). By default, theparameters
are determined automatically usingget_parameters
, respectivelydiscretize_on_midpoints
ifpwc=true
. If they are given manually, they must follow the exact same semantics. In particular, forpwc=false
, any parameters must alias the parameters in the controls, such that mutatingparameters
is automatically reflected inevaluate
. Theparameters
will be available as an attribute of theode_propagator
.alg
: The algorithm to use for the ODE Solver, see the list of solvers in the DifferentialEquations manual. The defaultTsit5()
method is the recommended choice for non-stiff problems.solver_options
: All other keyword arguments are passed to the ODE solver, see the list of Solve Keyword Arguments in the DifferentialEquations manual. Note that the options for "Default Algorithm Hinting" do not apply, sincealg
must be specified manually. Also, the "Output Control" is managed by theode_propagator
, so these options should not be used.
QuantumPropagators.prop_step!
— FunctionAdvance the propagator
by a single time step.
state = prop_step!(propagator)
returns the state obtained from propagating to the next point on the time grid from propagator.t
, respectively the previous point if propagator.backward
is true.
When the propagation would lead out of the time grid, prop_step!
leaves propagator
unchanged and returns nothing
. Thus, a return value of nothing
may be used to signal that a propagation has completed.
QuantumPropagators.propagate
— FunctionPropagate a state over an entire time grid.
state = propagate(
state,
generator,
tlist;
method, # mandatory keyword argument
check=true,
backward=false,
inplace=QuantumPropagators.Interfaces.supports_inplace(state),
verbose=false,
piecewise=nothing,
pwc=nothing,
storage=nothing,
observables=<store state>,
callback=nothing,
show_progress=false,
init_prop_kwargs...)
propagates state
of the entire time grid and returns the propagated states, or a storage array of data collected during the propagation. This high-level routine performs the following three steps:
If
check=true
(default), check thatstate
,generator
, andtlist
are consistent with the required interface.Initialize a
propagator
viainit_prop
:init_prop(state, generator, tlist; method, inplace, init_prop_kwargs...)
Call and return the result of
propagate(propagator; storage, observables, show_progress, callback)
Arguments
state
: The "initial" state for the propagation. Forbackward=false
, this state is taken to be at initial time (tlist[begin]
); and forbackward=true
, at the final time (tlist[end]
)generator
: The time-dependent generator of the dynamicstlist
: The time grid over which which the propagation is defined. This may or may not be equidistant.
Mandatory keyword arguments
method
: The propagation method to use. May be given as a name (Symbol
), but the recommended usage is to pass a module implementing the propagation method, cf.init_prop
.
Optional keyword arguments
check
: iftrue
, check thatstate
,generator
, andtlist
passcheck_state
,check_generator
andcheck_tlist
, respectively.backward
: Iftrue
, propagate backward in timeinplace
: Iftrue
, propagate using in-place operations. Iffalse
, avoid in-place operations. Not all propagation methods support both in-place and not-in-place propagation. Note thatinplace=true
requires thatQuantumPropagators.Interfaces.supports_inplace
forstate
istrue
.piecewise
: If given as a boolean, ensure that the internalpropagator
is an instance ofPiecewisePropagator
, cf.init_prop
.pwc
: If given a a boolean, do a piecewise constant propagation where the generator in each interval is constant (the internalpropagator
is aPWCPropagator
, cf.init_prop
)storage
: Flag whether to store and return the propagated states / observables, or pre-allocated storage array. See Notes below.observables
: Converters for data to be stored instorage
. See Notes below.callback
: Function to call after each propagation step. See Notes below.show_progress
: Whether to show a progress bar. See Notes below.
All remaining keyword arguments are passed to init_prop
to initialize the Propagator
that is used internally to drive the optimization. Unknown keyword arguments will be ignored.
Notes
In general, there is no requirement that tlist
has a constant time step, although some propagation methods (most notably Cheby
) only support a uniform time grid.
If storage
is given as a container pre-allocated via init_storage
, it will be filled with data determined by the observables
. Specifically, after each propagation step,
data = map_observables(observables, tlist, i, state)
write_to_storage!(storage, i, data)
is executed, where state
is defined at time tlist[i]
. See map_observables
and write_to_storage!
for details. The default values for observables
results simply in the propagated states at every point in time being stored.
The storage
parameter may also be given as true
, and a new storage array will be created internally with init_storage
and returned instead of the propagated state:
data = propagate(
state, generator, tlist; method,
backward=false; storage=true, observables=observables,
callback=nothing, show_progress=false, init_prop_kwargs...)
If backward
is true
, the input state is assumed to be at time tlist[end]
, and the propagation progresses backward in time (with a negative time step dt
). If storage
is given, it will be filled back-to-front during the backward propagation.
If callback
is given as a callable, it will be called after each propagation step, as callback(propagator, observables)
where propagator
is Propagator
object driving the propagation. The callback
is called before calculating any observables, or storing the propagated state in storage
. Example usage includes writing data to file, or modifying state
via set_state!
, e.g., removing amplitude from the lowest and highest level to mitigate "truncation error".
If show_progress
is given as true
, a progress bar will be shown for long-running propagation. In order to customize the progress bar, show_progress
may also be a function that receives length(tlist)
and returns a ProgressMeter.Progress
instance.
If in_place=false
is given, the propagation avoids in-place operations. This is slower than inplace=true
, but is often required in the context of automatic differentiation (AD), e.g., with Zygote. That is, use in_place=false
if propagate
is called inside a function to be passed to Zygote.gradient
, Zygote.pullback
, or a similar function. In an AD context, storage
and show_progress
should not be used.
The propagate
routine returns the propagated state at tlist[end]
, respectively tlist[1]
if backward=true
, or a storage array with the stored states / observable data if storage=true
.
state = propagate(
state,
propagator;
storage=nothing,
observables=<store state>,
show_progress=false,
callback=nothing,
reinit_prop_kwargs...
)
re-initializes the given propagator
with state
(see reinit_prop!
) and then calls the lower-level propagate(propagator; ...)
.
state = propagate(
propagator;
storage=nothing,
observables=<store state>,
show_progress=false,
callback=nothing,
)
propagates a freshly initialized propagator
(immediately after init_prop
). Used in the higher-level propagate(state, generator, tlist; kwargs...)
.
QuantumPropagators.propagate_sequence
— FunctionPropagate a state through a sequence of generators.
states = propagate_sequence(
state,
propagations;
storage=nothing,
pre_propagation=nothing,
post_propagation=nothing,
kwargs...
)
takes an initial
state and performs a sequence of propagate
calls using the parameters in propagations
. The initial state for each step in the sequence is the state resulting from the previous step. Optionally, before and after each step, a pre_propagation
and post_propagation
function may modify the state instantaneously, e.g., to perform a frame transformation. Return the vector of states at the end of each step (after any post_propagation
, before any next pre_propagation
of the next step).
Arguments
state
: The initial statepropagations
: A vector ofPropagation
instances, one per step in the sequence, each containing the arguments for the call topropagate
for that step. ThePropagation
contains the generator and time grid for each step as positional parameters, or alternatively a pre-initializedPropagator
, and any keyword arguments forpropagate
that are specific to that step. Note thatpropagate
keyword arguments that are common to all steps can be given directly topropagate_sequence
.storage
: Ifstorage=true
, return a vector of storage objects as returned bypropagate(…, storage=true)
for each propagation step, instead of the state after each step. To use a pre-initializedstorage
, eachPropagation
inpropagations
should have astorage
keyword argument instead.pre_propagation
: If notnothing
, must be a function that receives the same arguments aspropagate
and returns a state. Called immediately before thepropagate
of each step, and the state returned bypre_propagation
will become the initial state for the subsequent call topropagate
. Generally,pre_propagation
would be different in each step of the sequence, and should be given as a keyword argument in a particularPropagation
.post_propagation
: If notnothing
, a function that receives the same arguments aspropagate
and returns a state, seepre_propagation
. The returned state becomes the initial state for the next step in the sequence (and may be further processed by the followingpre_propagation
). Likepre_propagation
, this will generally be set as a keyword argument for a particularPropagation
, not as a global keyword argument topropagate_sequence
.
All other keyword arguments are forwarded to propagate
. Thus, keyword arguments that are common to all steps in the sequence should be given as keyword arguments to propagate_sequence
directly.
QuantumPropagators.reinit_prop!
— FunctionRe-initialize a propagator.
reinit_prop!(propagator, state; kwargs...)
resets the propagator
to state
at the beginning of the time grid, respectively the end of the time grid if propagator.backward
is true.
At a minimum, this is equivalent to a call to set_state!
follow by a call to set_t!
, but some propagators may have additional requirements on re-initialization, such as refreshing expansion coefficients for ChebyPropagator
. In this case, the kwargs
may be additional keyword arguments specific to the concrete type of propagator.
reinit_prop!(
propagator::ChebyPropagator,
state;
transform_control_ranges=((c, ϵ_min, ϵ_max, check) => (ϵ_min, ϵ_max)),
kwargs...
)
re-initializes an existing ChebyPropagator
. This may or may not involve recalculating the Chebychev coefficients based on the current control amplitudes in propagator.parameters
.
Method-specific keyword arguments
transform_control_ranges
: a function(c, ϵ_min, ϵ_max, check) => (ϵ_min′, ϵ_max′)
. For each controlc
, the function is called withcheck=true
andϵ_min
(ϵ_max
) the current minimum (maximum) values for the control frompropagator.parameters
). The Chebychev coefficients will be recalculated if the existing coefficients were obtained assuming a range forc
outside the returnedϵ_min′, ϵ_max′
.If the coefficients do need to be recalculated,
transform_control_ranges
is called a second time withcheck=false
, and the returned(ϵ_min′, ϵ_max′)
are used for estimating the new spectral range.For example,
function transform_control_ranges(c, ϵ_min, ϵ_max, check) if check return (min(ϵ_min, 2 * ϵ_min), max(ϵ_max, 2 * ϵ_max)) else return (min(ϵ_min, 5 * ϵ_min), max(ϵ_max, 5 * ϵ_max)) end end
will re-calculate the Chebychev coefficients only if the current amplitudes differ by more than a factor of two from the ranges that were used when initializing the propagator (
control_ranges
parameter ininit_prop
, which would have had to overestimate the actual amplitudes by at least a factor of two). When re-calculating, thecontrol_ranges
will overestimate the amplitudes by a factor of five. With thistransform_control_ranges
, the propagation will be stable as long as the amplitudes do not change dynamically by more than a factor of 2.5 from their original range, while also not re-calculating coefficients unnecessarily in each pass because of modest changes in the amplitudes.The
transform_control_ranges
argument is only relevant in the context of optimal control, where the samepropagator
will be used for many iterations with changing control field amplitudes.
All other keyword arguments are ignored.
Private members
QuantumPropagators.cheby_get_spectral_envelope
— FunctionDetermine the spectral envelope of a generator
.
E_min, E_max = cheby_get_spectral_envelope(
generator, tlist, control_ranges, method; kwargs...
)
estimates a lower bound E_min
the lowest eigenvalue of the generator for any values of the controls specified by control_ranges
, and an upper bound E_max
for the highest eigenvalue.
This is done by constructing operators from the extremal values for the controls as specified in control_ranges
and taking the smallest/largest return values from specrange
for those operators.
Arguments
generator
: dynamical generator, e.g. a time-dependenttlist
: The time grid for the propagationcontrol_ranges
: a dict that maps controls that occur ingenerator
(cf.get_controls
to a tuple of minimum and maximum amplitude for that controlmethod
: method name to pass tospecrange
kwargs
: Any remaining keyword arguments are passed tospecrange
QuantumPropagators.ExpPropagator
— TypePropagator for propagation via direct exponentiation (method=QuantumPropagators.ExpProp
)
This is a PWCPropagator
.
QuantumPropagators.PWCPropagator
— TypePiecewisePropagator
sub-type for piecewise-constant propagators.
Like the more general PiecewisePropagator
, this is characterized by propagator.parameters
mapping the controls in the generator
to a vector of amplitude value on the midpoints of the time grid intervals.
The propagation will use these values as constant within each interval.
QuantumPropagators.ChebyPropagator
— TypePropagator for Chebychev propagation (method=QuantumPropagators.Cheby
).
This is a PWCPropagator
.
QuantumPropagators.AbstractPropagator
— TypeAbstract base type for all Propagator
objects.
All Propagator
objects must be instantiated via init_prop
and implement the following interface.
Properties
state
(read-only): The current quantum state in the propagationtlist
(read-only): The time grid for the propagationt
(read-only): The time at whichstate
is defined. An element oftlist
.parameters
: parameters that determine the dynamics. The structure of the parameters depends on the concretePropagator
type (i.e., the propagation method). Mutating theparameters
affects subsequent propagation steps.backward
: Boolean flag to indicate whether the propagation moves forward or backward in timeinplace
: Boolean flag to indicate whetherpropagator.state
is modified in-place or is recreated by every call ofprop_step!
orset_state!
. Withinplace=false
, the propagator should generally avoid in-place operations, such as calls toQuantumPropagators.Controls.evaluate!
.
Concrete Propagator
types may have additional properties or fields, but these should be considered private.
Methods
reinit_prop!
— reset the propagator to a new initial state at the beginning of the time grid (or the end, for backward propagation)prop_step!
– advance the propagator by one step forward (or backward) on the time grid.set_state!
— safely mutate the current quantumstate
of the propagation. Note that directly mutating thestate
property is not safe. However,Ψ = propagator.state; foo_mutate!(Ψ), set_state!(propagator, Ψ)
for some mutating functionfoo_mutate!
is guaranteed to be safe and efficient for both in-place and not-in-place propagators.set_t!
— safely mutate the current time (propagator.t
), snapping to the values oftlist
.
See also
PiecewisePropagator
— specialization ofAbstractPropagator
for piecewise propagation methods.PWCPropagator
— specialization ofPiecewisePropagator
for piecewise-constant propagation methods.
QuantumPropagators.enable_timings
— FunctionEnable the collection of TimerOutputs
data.
QuantumPropagators.enable_timings()
enables certain portions of the package to collect TimerOutputs
internally. This aids in profiling and benchmarking propagation methods.
Specifically, after enable_timings()
, for any ChebyPropagator
or NewtonPropagator
, timing data will become available in propagator.wrk.timing_data
(as a TimerOutput
instance). This data is reset when the propagator is re-instantiated with init_prop
or re-initialized with reinit_prop!
. This makes the data local to any call of propagate
.
Note that enable_timings()
triggers recompilation, so propagate
should be called at least twice to avoid compilation overhead in the timing data. There is still a small overhead for collecting the timing data.
The collection of timing data can be disabled again with disable_timings
.
Returns QuantumPropagators.timings_enabled()
, i.e., true
if successful.
QuantumPropagators.set_state!
— FunctionSet the current state
of the propagator
.
set_state!(propagator, state)
sets the propagator.state
property and returns propagator.state
. In order to mutate the current state after a call to prop_step!
, the following pattern is recommended:
Ψ = propagator.state
foo_mutate!(Ψ)
set_state!(propagator, Ψ)
where foo_mutate!
is some function that mutates Ψ
. This is guaranteed to work efficiently both for in-place and not-in-place propagators, without incurring unnecessary copies.
foo_mutate!(propagator.state)
by itself is not a safe operation. Always follow it by
set_state!(propagator, propagator.state)
See also
set_t!
— setpropagator.t
.
QuantumPropagators.set_t!
— FunctionSet the current time for the propagation.
set_t!(propagator, t)
Sets propagator.t
to the given value of t
, where t
must be an element of propagator.tlist
.
See also
set_state!
— setpropagator.state
QuantumPropagators.timings_enabled
— FunctionCheck whether the collection of TimerOutputs
data is active.
QuantumPropagators.timings_enabled()
returns true
if QuantumPropagators.enable_timings()
was called, and false
otherwise or after QuantumPropagators.disable_timings()
.
QuantumPropagators.PiecewisePropagator
— TypeAbstractPropagator
sub-type for piecewise propagators.
A piecewise propagator is determined by a single parameter per control and time grid interval. Consequently, the propagator.parameters
are a dictionary mapping the controls found in the generator
via get_controls
to a vector of values defined on the intervals of the time grid, see discretize_on_midpoints
. This does not necessarily imply that these values are the piecewise-constant amplitudes for the intervals. A general piecewise propagator might use interpolation to obtain actual amplitudes within any given time interval.
When the amplitudes are piecewise-constant, the propagator should be a concrete instantiation of a PWCPropagator
.
QuantumPropagators.NewtonPropagator
— TypePropagator for Newton propagation (method=QuantumPropagators.Newton
).
This is a PWCPropagator
.
QuantumPropagators.disable_timings
— FunctionDisable the collection of TimerOutputs
data.
QuantumPropagators.disable_timings()
disables the collection of timing data previously enabled with enable_timings
. This triggers recompilation to completely remove profiling from the code. That is, there is zero cost when the collection of timing data is disabled.
Returns QuantumPropagators.timings_enabled()
, i.e., false
if successful.
QuantumPropagators.ode_function
— FunctionWrap around a Generator
, for use as an ODE function.
f = ode_function(generator, tlist; c=-1im)
creates a function suitable to be passed to ODEProblem
.
\[\gdef\op#1{\hat{#1}} \gdef\ket#1{\vert{#1}\rangle}\]
With generator
corresponding to $\op{H}(t)$, this implicitly encodes the ODE
\[\frac{\partial}{\partial t} \ket{\Psi(t)} = c \op{H}(t) \ket{\Psi(t)}\]
for the state $\ket{\Psi(t)}$. With the default $c = -i$, this corresponds to the Schrödinger equation, or the Liouville equation with convention=:LvN
.
The resulting f
works both in-place and not-in-place, as
f(ϕ, Ψ, vals_dict, t) # in-place `f(du, u, p, t)`
ϕ = f(Ψ, vals_dict, t) # not-in-place `f(u, p, t)`
Calling f
as above is functionally equivalent to calling evaluate
to obtain an operator H
from the original time-dependent generator
, and then applying H
to the current quantum state Ψ
:
H = evaluate(f.generator, t; vals_dict=vals_dict)
ϕ = c * H * Ψ
where vals_dict
may be a dictionary mapping controls to values (set as the parameters p
of the underlying ODE solver).
If QuantumPropagators.enable_timings()
has been called, profiling data is collected in f.timing_data
.
QuantumPropagators.Generators
The following routines define dynamical generators (Hamiltonians, Liouvillians). This includes the initialization of generators and the methods that define how these generators contain controls and control amplitudes. These methods must be defined for any custom generator or control amplitude.
Reference
Public Members:
Public members
QuantumPropagators.Generators.Generator
— TypeA time-dependent generator.
Generator(ops::Vector{OT}, amplitudes::Vector{AT})
produces an object of type Generator{OT,AT}
that represents
\[Ĥ(t)= Ĥ_0 + \sum_l a_l(\{ϵ_{l'}(t)\}, t) \, Ĥ_l\,,\]
where $Ĥ_l$ are the ops
and $a_l(t)$ are the amplitudes
. $Ĥ(t)$ and $Ĥ_l$ may represent operators in Hilbert space or super-operators in Liouville space. If the number of amplitudes
is less than the number of ops
, the first ops
are considered as drift terms ($Ĥ_0$, respectively subsequent terms with $a_l ≡ 1$). At least one time-dependent amplitude is required. Each amplitude may depend on one or more control functions $ϵ_{l'}(t)$, although most typically $a_l(t) ≡ ϵ_l(t)$, that is, the amplitudes
are simply a vector of the controls. See hamiltonian
for details.
A Generator
object should generally not be instantiated directly, but via hamiltonian
or liouvillian
.
The list of ops
and amplitudes
are properties of the Generator
. They should not be mutated.
See also
QuantumPropagators.Generators.Operator
— TypeA static operator in Hilbert or Liouville space.
Operator(ops::Vector{OT}, coeffs::Vector{CT})
produces an object of type Operator{OT,CT}
that encapsulates the "lazy" sum
\[Ĥ = \sum_l c_l Ĥ_l\,,\]
where $Ĥ_l$ are the ops
and $c_l$ are the coeffs
, which each must be a constant Number
. If the number of coefficients is less than the number of operators, the first ops
are considered to have $c_l = 1$.
An Operator
object would generally not be instantiated directly, but be obtained from a Generator
via evaluate
.
The $Ĥ_l$ in the sum are considered immutable. This implies that an Operator
can be updated in-place with evaluate!
by only changing the coeffs
.
QuantumPropagators.Generators.ScaledOperator
— TypeA static operator with a scalar pre-factor.
op = ScaledOperator(α, Ĥ)
represents the "lazy" product $α Ĥ$ where $Ĥ$ is an operator (typically an Operator
instance) and $α$ is a scalar.
QuantumPropagators.Generators.hamiltonian
— FunctionInitialize a (usually time-dependent) Hamiltonian.
The most common usage is, e.g.,
using QuantumPropagators
H₀ = ComplexF64[0 0; 0 1];
H₁ = ComplexF64[0 1; 1 0];
ϵ₁(t) = 1.0;
hamiltonian(H₀, (H₀, ϵ₁))
# output
Generator with 2 ops and 1 amplitudes
ops::Vector{Matrix{ComplexF64}}:
ComplexF64[0.0 + 0.0im 0.0 + 0.0im; 0.0 + 0.0im 1.0 + 0.0im]
ComplexF64[0.0 + 0.0im 0.0 + 0.0im; 0.0 + 0.0im 1.0 + 0.0im]
amplitudes::Vector{typeof(ϵ₁)}:
ϵ₁
In general,
H = hamiltonian(terms...; check=true)
constructs a Hamiltonian based on the given terms
. Each term must be an operator or a tuple (op, ampl)
of an operator and a control amplitude. Single operators are considered "drift" terms.
In most cases, each control amplitude will simply be a control function or vector of pulse values. In general, ampl
can be an arbitrary object that depends on one or more controls, which must be obtainable via get_controls(ampl)
. See QuantumPropagators.Interfaces.check_amplitude
for the required interface.
The hamiltonian
function will generally return a Generator
instance. However, if none of the given terms are time-dependent, it may return a static operator (e.g., an AbstractMatrix
or Operator
):
hamiltonian(H₀)
# output
2×2 Matrix{ComplexF64}:
0.0+0.0im 0.0+0.0im
0.0+0.0im 1.0+0.0im
hamiltonian(H₀, (H₁, 2.0))
# output
Operator with 2 ops and 1 coeffs
ops::Vector{Matrix{ComplexF64}}:
ComplexF64[0.0 + 0.0im 0.0 + 0.0im; 0.0 + 0.0im 1.0 + 0.0im]
ComplexF64[0.0 + 0.0im 1.0 + 0.0im; 1.0 + 0.0im 0.0 + 0.0im]
coeffs: [2.0]
The hamiltonian
function may generate warnings if the terms
are of an unexpected type or structure. These can be suppressed with check=false
.
QuantumPropagators.Generators.liouvillian
— FunctionConstruct a Liouvillian Generator
.
ℒ = liouvillian(Ĥ, c_ops=(); convention=:LvN, check=true)
calculates the sparse Liouvillian super-operator ℒ
from the Hamiltonian Ĥ
and a list c_ops
of Lindblad operators.
With convention=:LvN
, applying the resulting ℒ
to a vectorized density matrix ρ⃗
calculates $\frac{d}{dt} \vec{\rho}(t) = ℒ \vec{\rho}(t)$ equivalent to the Liouville-von-Neumann equation for the density matrix $ρ̂$,
\[\frac{d}{dt} ρ̂(t) = -i [Ĥ, ρ̂(t)] + \sum_k\left( Â_k ρ̂ Â_k^\dagger - \frac{1}{2} A_k^\dagger Â_k ρ̂ - \frac{1}{2} ρ̂ Â_k^\dagger Â_k \right)\,,\]
where the Lindblad operators $Â_k$ are the elements of c_ops
.
The Hamiltonian $Ĥ$ will generally be time-dependent. For example, it may be a Generator
as returned by hamiltonian
. For example, for a Hamiltonian with the terms (Ĥ₀, (Ĥ₁, ϵ₁), (Ĥ₂, ϵ₂))
, where Ĥ₀
, Ĥ₁
, Ĥ₂
are matrices, and ϵ₁
and ϵ₂
are functions of time, the resulting ℒ
will be a Generator
corresponding to terms (ℒ₀, (ℒ₁, ϵ₁), (ℒ₂, ϵ₂))
, where the initial terms is the superoperator ℒ₀
for the static component of the Liouvillian, i.e., the commutator with the drift Hamiltonian Ĥ₀
, plus the dissipator (sum over $k$), as a sparse matrix. Time-dependent Lindblad operators are not currently supported. The remaining elements are tuples (ℒ₁, ϵ₁)
and (ℒ₂, ϵ₂)
corresponding to the commutators with the two control Hamiltonians, where ℒ₁
and ℒ₂
again are sparse matrices.
If $Ĥ$ is not time-dependent, the resulting ℒ
will likewise be a static operator. Passing H=nothing
with non-empty c_ops
initializes a pure dissipator.
With convention=:TDSE
, the Liouvillian will be constructed for the equation of motion $i \hbar \frac{d}{dt} \vec{\rho}(t) = ℒ \vec{\rho}(t)$ to match exactly the form of the time-dependent Schrödinger equation. While this notation is not standard in the literature of open quantum systems, it has the benefit that the resulting ℒ
can be used in a numerical propagator for a (non-Hermitian) Schrödinger equation without any change. Thus, for numerical applications, convention=:TDSE
is generally preferred. The returned ℒ
between the two conventions differs only by a factor of $i$, since we generally assume $\hbar=1$.
The convention
keyword argument is mandatory, to force a conscious choice.
See Goerz et. al. "Optimal control theory for a unitary operation under dissipative evolution", arXiv 1312.0111v2, Appendix B.2 for the explicit construction of the Liouvillian superoperator as a sparse matrix.
Passing check=false
, suppresses warnings and errors about unexpected types or the structure of the arguments, cf. hamiltonian
.
QuantumPropagators.Shapes
The following routines define useful function $S(t)$ that can be used for control functions or amplitudes.
Reference
Public Members:
Public members
QuantumPropagators.Shapes.blackman
— FunctionBlackman window shape.
blackman(t, t₀, T; a=0.16)
calculates
\[B(t; t_0, T) = \frac{1}{2}\left( 1 - a - \cos\left(2π \frac{t - t_0}{T - t_0}\right) + a \cos\left(4π \frac{t - t_0}{T - t_0}\right) \right)\,,\]
for a scalar t
, with $a$ = 0.16.
See https://en.wikipedia.org/wiki/Window_function#Blackman_windows
A Blackman shape looks nearly identical to a Gaussian with a 6-sigma interval between t₀
and T
. Unlike the Gaussian, however, it will go exactly to zero at the edges. Thus, Blackman pulses are often preferable to Gaussians.
QuantumPropagators.Shapes.box
— FunctionBox shape (Theta-function).
box(t, t₀, T)
evaluates the Heaviside (Theta-) function $\Theta(t) = 1$ for $t_0 \le t \le T$; and $\Theta(t) = 0$ otherwise.
QuantumPropagators.Shapes.flattop
— FunctionFlat shape (amplitude 1.0) with a switch-on/switch-off from zero.
flattop(t; T, t_rise, t₀=0.0, t_fall=t_rise, func=:blackman)
evaluates a shape function that starts at 0 at $t=t₀$, and ramps to to 1 during the t_rise
interval. The function then remains at value 1, before ramping down to 0 again during the interval t_fall
before T
. For $t < t₀$ and $t > T$, the shape is zero.
The default switch-on/-off shape is half of a Blackman window (see blackman
).
For func=:sinsq
, the switch-on/-off shape is a sine-squared curve.
QuantumPropagators.Controls
The following routines define method that must be defined for any control function, or for generators with respect to control functions.
Reference
Public Members:
ParameterizedFunction
discretize
discretize_on_midpoints
evaluate
evaluate!
get_controls
get_parameters
get_tlist_midpoints
substitute
t_mid
Public members
QuantumPropagators.Controls.ParameterizedFunction
— TypeAbstract type for function-like objects with parameters
.
A struct that is an implementation of a ParameterizedFunction
:
- must have a
parameters
field that is anAbstractVector
of floats (e.g., aComponentArrays.ComponentVector
) - must be callable with a single float argument
t
, - may define getters and setters for referencing the values in
parameters
with convenient names.
The parameters
field of any ParameterizedFunction
can be accessed via get_parameters
.
See How to define a parameterized control for an example. You may use the QuantumPropagators.Interfaces.check_parameterized_function
to check the implementation of a ParameterizedFunction
subtype.
QuantumPropagators.Controls.discretize
— FunctionEvaluate control
at every point of tlist
.
values = discretize(control, tlist; via_midpoints=true)
discretizes the given control
to a Vector of values defined on the points of tlist
.
If control
is a function, it is first evaluated at the midpoint of tlist
, see discretize_on_midpoints
, and then the values on the midpoints are converted to values on tlist
. This discretization is more stable than directly evaluating the control function at the values of tlist
, and ensures that repeated round-trips between discretize
and discretize_on_midpoints
can be done safely, see the note in the documentation of discretize_on_midpoints
.
The latter can still be achieved by passing via_midpoints=false
. While such a direct discretization is suitable e.g. for plotting, but it is unsuitable for round-trips between discretize
and discretize_on_midpoints
(constant controls on tlist
may result in a zig-zag on the intervals of tlist
).
If control
is a vector, a copy of control
will be returned if it is of the same length as tlist
. Otherwise, control
must have one less value than tlist
, and is assumed to be defined on the midpoints of tlist
. In that case, discretize
acts as the inverse of discretize_on_midpoints
. See discretize_on_midpoints
for how control values on tlist
and control values on the intervals of tlist
are related.
QuantumPropagators.Controls.discretize_on_midpoints
— FunctionEvaluate control
at the midpoints of tlist
.
values = discretize_on_midpoints(control, tlist)
discretizes the given control
to a Vector of values on the midpoints of tlist
. Hence, the resulting values
will contain one less value than tlist
.
If control
is a vector of values defined on tlist
(i.e., of the same length as tlist
), it will be converted to a vector of values on the intervals of tlist
. The value for the first and last "midpoint" will remain the original values at the beginning and end of tlist
, in order to ensure exact boundary conditions. For all other midpoints, the value for that midpoint will be calculated by "un-averaging".
For example, for a control
and tlist
of length 5, consider the following diagram:
tlist index: 1 2 3 4 5
tlist: ⋅ ⋅ ⋅ ⋅ ⋅ input values cᵢ (i ∈ 1..5)
|̂/ ̄ ̄ ̂\ / ̂\ / ̂ ̄ ̄\|̂
midpoints: x x x x output values pᵢ (i ∈ 1..4)
midpoints index: 1 2 3 4
We will have $p₁=c₁$ for the first value, $p₄=c₅$ for the last value. For all other points, the control values $cᵢ = \frac{p_{i-1} + p_{i}}{2}$ are the average of the values on the midpoints. This implies the "un-averaging" for the midpoint values $pᵢ = 2 c_{i} - p_{i-1}$.
An arbitrary input control
array may not be compatible with the above averaging formula. In this case, the conversion will be "lossy" (discretize
will not recover the original control
array; the difference should be considered a "discretization error"). However, any further round-trip conversions between points and intervals are bijective and preserve the boundary conditions. In this case, the discretize_on_midpoints
and discretize
methods are each other's inverse. This also implies that for an optimal control procedure, it is safe to modify midpoint values. Modifying the the values on the time grid directly on the other hand may accumulate discretization errors.
If control
is a vector of one less length than tlist
, a copy of control
will be returned, under the assumption that the input is already properly discretized.
If control
is a function, the function will be directly evaluated at the midpoints marked as x
in the above diagram..
See also
get_tlist_midpoints
– get all the midpoints on which the control will be discretized.t_mid
– get a particular midpoint.discretize
– discretize directly ontlist
instead of on the midpoints
QuantumPropagators.Controls.evaluate
— FunctionEvaluate all controls.
In general, evaluate(object, args...; vals_dict=IdDict())
evaluates the object
for a specific point in time indicated by the positional args
. Any control in object
is evaluated at the specified point in time. Alternatively, the vals_dict
maps a controls to value ("plug in this value for the given control")
For example,
op = evaluate(generator, t)
evaluates generator
at time t
. This requires that any control in generator
is a callable that takes t
as a single argument.
op = evaluate(generator, tlist, n)
evaluates generator
for the n'th interval of tlist
. This uses the definitions for the midpoints in discretize_on_midpoints
. The controls in generator
may be vectors (see discretize
, discretize_on_midpoints
) or callables of t
.
op = evaluate(generator, t; vals_dict)
op = evaluate(generator, tlist, n; vals_dict)
resolves any explicit time dependencies in generator
at the specified point in time, but uses the value in the given vals_dict
for any control in vals_dict
.
a = evaluate(ampl, tlist, n; vals_dict=IdDict())
a = evaluate(ampl, t; vals_dict=IdDict())
evaluates a control amplitude to a scalar by evaluating any explicit time dependency, and by replacing each control with the corresponding value in vals_dict
.
Calling evaluate
for an object with no implicit or explicit time dependence should return the object unchanged.
For generators without any explicit time dependence,
op = evaluate(generator; vals_dict)
can be used. The vals_dict
in this case must contain values for all controls in generator
.
See also:
evaluate!
— update an existing operator with a re-evaluation of a
generator at a different point in time.
QuantumPropagators.Controls.evaluate!
— FunctionUpdate an existing evaluation of a generator
.
evaluate!(op, generator, args..; vals_dict=IdDict())
performs an in-place update on an op
the was obtained from a previous call to evaluate
with the same generator
, but for a different point in time and/or different values in vals_dict
.
QuantumPropagators.Controls.get_controls
— FunctionExtract a Tuple of controls.
controls = get_controls(generator)
extracts the controls from a single dynamical generator.
For example, if generator = hamiltonian(H0, (H1, ϵ1), (H2, ϵ2))
, extracts (ϵ1, ϵ2)
.
get_controls(operator)
for a static operator (matrix) returns an empty tuple.
QuantumPropagators.Controls.get_parameters
— FunctionObtain analytic parameters of the given control
.
parameters = get_parameters(control)
obtains parameters
as an AbstractVector{Float64}
containing any tunable analytic parameters associated with the control
. The specific type of parameters
depends on how control
is defined, but a ComponentArrays.ComponentVector
should be a common array type.
Mutating the resulting vector must directly affect the control in any subsequent call to evaluate
. That is, the values in parameters
must alias values inside the control
.
Note that the control
must be an object specifically designed to have analytic parameters. Typically, it should be implemented as a subtype of ParameterizedFunction
. For a simple function ϵ(t)
or a vector of pulse values, which are the default types of controls discussed in the documentation of hamiltonian
, the get_parameters
function will return an empty vector.
More generally,
parameters = get_parameters(object)
collects and combines all unique parameter arrays from the controls inside the object
. The object
may be a Generator
, Trajectory
, ControlProblem
, or any other object for which get_controls(object)
is defined. If there are multiple controls with different parameter arrays, these are combined in a RecursiveArrayTools.ArrayPartition
. This requires the RecursiveArrayTools
package to be loaded. Again, mutating parameters
directly affects the underlying controls.
The parameters
may be used as part of the parameters
attribute of a propagator for time-continuous dynamics, like a general ODE solver, or in an optimization that tunes analytic control parameters, e.g., with a Nelder-Mead method. Examples might include the widths, peak amplitudes, and times of a superposition of Gaussians [7], cf. the example of a ParameterizedFunction
, or the amplitudes associated with spectral components in a random truncated basis [11].
The parameters
are not intended for optimization methods such as GRAPE or Krotov that fundamentally use a piecewise-constant control ansatz. In the context of such methods, the "control parameters" are always the amplitudes of the control
at the mid-points of the time grid, as obtained by discretize_on_midpoints
, and get_parameters
is ignored.
QuantumPropagators.Controls.get_tlist_midpoints
— FunctionShift time grid values to the interval midpoints
tlist_midpoints = get_tlist_midpoints(
tlist; preserve_start=true, preserve_end=true
)
takes a vector tlist
of length $n$ and returns a Vector{Float64}
of length $n-1$ containing the midpoint values of each interval. The intervals in tlist
are not required to be uniform.
By default, the first and last point of tlist
is preserved, see discretize_on_midpoints
. This behavior can be disabled by passing preserve_start
and preserve_end
as false
in order to use the midpoints of the first and last interval, respectively.
See also
t_mid
– get a particular midpoint.
QuantumPropagators.Controls.substitute
— FunctionSubstitute inside the given object.
object = substitute(object, replacements)
returns a modified object with the replacements defined in the given replacements
dictionary. Things that can be replaced include operators, controls, and amplitudes. For example,
generator = substitute(generator::Generator, replacements)
operator = substitute(operator::Operator, replacements)
amplitude = substitute(amplitude, controls_replacements)
Note that substitute
cannot be used to replace dynamic quantities, e.g. controls, with static value. Use evaluate
instead for that purpose.
QuantumPropagators.Controls.t_mid
— FunctionMidpoint of n'th interval of tlist.
t = t_mid(tlist, n)
returns the t
that is the midpoint between points tlist[n+1]
and tlist[n]
, but snapping to the beginning/end to follow the convention explained in discretize_on_midpoints
(to preserve exact boundary conditions at the edges of the time grid.)
See also
get_tlist_midpoints
– get all the midpoints in one go.
QuantumPropagators.Amplitudes
The following types define control amplitudes $a(\{ϵ_l(t)\}, t)$ that depend on one or more control functions $ϵ_l(t)$.
Reference
Public Members:
Public members
QuantumPropagators.Amplitudes.LockedAmplitude
— TypeA time-dependent amplitude that is not a control.
ampl = LockedAmplitude(shape)
wraps around shape
, which must be either a vector of values defined on the midpoints of a time grid or a callable shape(t)
.
ampl = LockedAmplitude(shape, tlist)
discretizes shape
to the midpoints of tlist
.
QuantumPropagators.Amplitudes.ShapedAmplitude
— TypeProduct of a fixed shape and a control.
ampl = ShapedAmplitude(control; shape=shape)
produces an amplitude $a(t) = S(t) ϵ(t)$, where $S(t)$ corresponds to shape
and $ϵ(t)$ corresponds to control
. Both control
and shape
should be either a vector of values defined on the midpoints of a time grid or a callable control(t)
, respectively shape(t)
. In the latter case, ampl
will also be callable.
ampl = ShapedAmplitude(control, tlist; shape=shape)
discretizes control
and shape
to the midpoints of tlist
.
QuantumPropagators.Storage
The following routines allow to manage and extend storage arrays storage
parameter in propagate
. See the discussion of Expectation Values for more details.
Reference
Public Members:
Public members
QuantumPropagators.Storage.get_from_storage
— FunctionQuantumPropagators.Storage.get_from_storage!
— FunctionObtain data from storage.
get_from_storage!(data, storage, i)
extracts data from the storage
for the i'th time slot. Inverse of write_to_storage!
. This modifies data
in-place. If get_from_storage!
is implemented for arbitrary observables
, it is the developer's responsibility that init_storage
, write_to_storage!
, and get_from_storage!
are compatible.
To extract immutable data
, the non-in-place version
data = get_from_storage(storage, i)
can be used.
QuantumPropagators.Storage.init_storage
— FunctionCreate a storage
array for propagation.
storage = init_storage(state, tlist)
creates a storage array suitable for storing a state
for each point in tlist
.
storage = init_storage(state, tlist, observables)
creates a storage array suitable for the data generated by the observables
applied to state
, see map_observables
, for each point in tlist
.
storage = init_storage(data, nt)
creates a storage arrays suitable for storing data
nt times, where nt=length(tlist)
. By default, this will be a vector of typeof(data)
and length nt
, or a n × nt
Matrix with the same eltype
as data
if data
is a Vector of length n
.
QuantumPropagators.Storage.map_observable
— FunctionApply a single observable
to state
.
data = map_observable(observable, tlist, i, state)
By default, observable
can be one of the following:
- A function taking the three arguments
state
,tlist
,i
, wherestate
is defined at timetlist[i]
. - A function taking a single argument
state
, under the assumption that the observable is time-independent - A matrix for which to calculate the expectation value with respect to the vector
state
.
The default map_observables
delegates to this function.
QuantumPropagators.Storage.map_observables
— FunctionObtain "observable" data from state
.
data = map_observables(observables, tlist, i, state)
calculates the data for a tuple of observables
applied to state
defined at time tlist[i]
. For a single observable (tuple of length 1), simply return the result of map_observable
.
For multiple observables, return the tuple resulting from applying map_observable
for each observable. If the tuple is "uniform" (all elements are of the same type, e.g. if each observable calculates the expectation value of a Hermitian operator), it is converted to a Vector. This allows for compact storage in a storage array, see init_storage
.
QuantumPropagators.Storage.write_to_storage!
— FunctionPlace data into storage
for time slot i
.
write_to_storage!(storage, i, data)
for a storage
array created by init_storage
stores the data
obtained from map_observables
at time slot i
.
Conceptually, this corresponds roughly to storage[i] = data
, but storage
may have its own idea on how to store data for a specific time slot. For example, with the default init_storage
Vector data will be stored in a matrix, and write_to_storage!
will in this case write data to the i'th column of the matrix.
For a given type of storage
and data
, it is the developer's responsibility that init_storage
and write_to_storage!
are compatible.
QuantumPropagators.Cheby
The following routines implement time propagation via expansion in Chebychev polynomials.
Reference
Public Members:
Public members
QuantumPropagators.Cheby.ChebyWrk
— TypeWorkspace for the Chebychev propagation routine.
ChebyWrk(Ψ, Δ, E_min, dt; limit=1e-12)
initializes the workspace for the propagation of a state similar to Ψ
under a Hamiltonian with eigenvalues between E_min
and E_min + Δ
, and a time step dt
. Chebychev coefficients smaller than the given limit
are discarded.
QuantumPropagators.Cheby.cheby
— FunctionEvaluate Ψ = exp(-𝕚 * H * dt) Ψ
.
Ψ_out = cheby(Ψ, H, dt, wrk; E_min=nothing, check_normalization=false)
acts like cheby!
but does not modify Ψ
in-place.
QuantumPropagators.Cheby.cheby!
— FunctionEvaluate Ψ = exp(-𝕚 * H * dt) Ψ
in-place.
cheby!(Ψ, H, dt, wrk; E_min=nothing, check_normalization=false)
Arguments
Ψ
: on input, initial vector. Will be overwritten with result.H
: Hermitian operatordt
: time stepwrk
: internal workspaceE_min
: minimum eigenvalue of H, to be used instead of theE_min
from the initialization ofwrk
. The samewrk
may be used for different valuesE_min
, as long as the spectra radiusΔ
and the time stepdt
are the same as those used for the initialization ofwrk
.check_normalizataion
: perform checks that the H does not exceed the spectral radius for which the workspace was initialized.
The routine will not allocate any internal storage. This implementation requires copyto!
lmul!
, and axpy!
to be implemented for Ψ
, and the three-argument mul!
for Ψ
and H
.
QuantumPropagators.Cheby.cheby_coeffs
— FunctionCalculate Chebychev coefficients.
a::Vector{Float64} = cheby_coeffs(Δ, dt; limit=1e-12)
return an array of coefficients larger than limit
.
Arguments
Δ
: the spectral radius of the underlying operatordt
: the time step
See also cheby_coeffs!
for an in-place version.
QuantumPropagators.Cheby.cheby_coeffs!
— FunctionCalculate Chebychev coefficients in-place.
n::Int = cheby_coeffs!(coeffs, Δ, dt, limit=1e-12)
overwrites the first n
values in coeffs
with new coefficients larger than limit
for the given new spectral radius Δ
and time step dt
. The coeffs
array will be resized if necessary, and may have a length > n
on exit.
See also cheby_coeffs
for an non-in-place version.
QuantumPropagators.Newton
The following routines implement time propagation via expansion in Newton polynomials using a restarted Arnoldi scheme to determine evaluation points.
Reference
Public Members:
Private Members:
Public members
QuantumPropagators.Newton.NewtonWrk
— TypeNewtonWrk(v0, m_max=10)
Workspace for the Newton-with-restarted-Arnoldi propagation routine.
Initializes the workspace for the propagation of a vector v0
, using a maximum Krylov dimension of m_max
in each restart iteration. Note that m_max
should be smaller than the length of v0
.
QuantumPropagators.Newton.newton!
— Functionnewton!(Ψ, H, dt, wrk; func=(z -> exp(-1im*z)), norm_min=1e-14, relerr=1e-12,
max_restarts=50, _...)
Evaluate Ψ = func(H*dt) Ψ
using a Newton-with-restarted-Arnoldi scheme.
Arguments
Ψ
: The state to propagate, will be overwritten in-place with the propagated stateH
: Operator acting onΨ
. Together withdt
, this is the argument tofunc
dt
: Implicit time step. Together withH
, this is the argument tofunc
wkr
: Work array, initialized withNewtonWrk
func
: The function to apply toH dt
, taking a single (scalar) complex-valued argumentz
in place ofH dt
. The defaultfunc
is to evaluate the time evaluations operator for the Schrödinger equationnorm_min
: the minimum norm at which to consider a state similar toΨ
as zerorelerr
: The relative error defining the convergence condition for the restart iteration. Propagation stops when the norm of the accumulatedΨ
is stable up to the given relative errormax_restarts
: The maximum number of restart iterations. Exceedingmax_restarts
will throw anAssertionError
.
All other keyword arguments are ignored.
Private members
QuantumPropagators.Newton.extend_newton_coeffs!
— Functionextend_newton_coeffs!(a, n_a, leja, func, n_leja, radius)
Extend the array a
of existing Newton coefficients for the expansion of the func
from n_a
coefficients to n_leja
coefficients. Return a new value n_a=n_a+n_leja
with the total number of Newton coefficients in the updated a
.
Arguments
a
: On input, a zero-based array of lengthn_a
or greater, containing Newton coefficients. On output, array containing a totaln_leja
coefficients. The arraya
will be resized if necessary, and may have a length greater thann_leja
on outputn_a
: The number of Newton coefficients ina
, on input. Elements ofa
beyond the firstn_a
elements will be overwritten.leja
: Array of normalized Leja points, containing at leastn_leja
elements.func
: Function for which to calculate Newton coefficientsn_leja
: The number of elements inleja
to use for calculating new coefficients, and the total number of Newton coefficients on outputradius
: Normalization radius for divided differences
QuantumPropagators.Newton.extend_leja!
— Functionextend_leja!(leja, n, newpoints, n_use)
Given an array of n
(ordered) Leja points, extract n_use
points from newpoints
, and append them to the existing Leja points. The array leja
should be sufficiently large to hold the new Leja points, which are appended after index n_old
. It will be re-allocated if necessary and may have a size of up to 2*(n+n_use)
.
Arguments
leja
: Array of leja values. Must contain the "old" leja values to be kept inleja(0:n-1)
. On output,n_use
new leja points will be inleja(n+:n+n_use-1)
, for the original value ofn
. Theleja
array must use zero-based indexing.n
: On input, number of "old" leja points inleja
. On output, total number of leja points (i.e.n=n+n_use
)newpoints
: On input, candidate points for new leja points. Then_use
best values will be chosen and added toleja
. On output, the values ofnew_points
are undefined.n_use
: Number of points that should be added toleja
QuantumPropagators.ExpProp
The following routines implement time propagation via explicit exponentiation of the dynamical generator.
Reference
Public Members:
Private Members:
Public members
QuantumPropagators.ExpProp.ExpPropWrk
— TypeExpPropWrk(v0)
Workspace for propagation via direct matrix exponentiation.
Initializes the workspace for the propagation of a vector v0
QuantumPropagators.ExpProp.expprop!
— Functionexpprop!(Ψ, H, dt, wrk; func=(H_dt -> exp(-1im * H_dt)), _...)
Evaluate Ψ = func(H*dt) Ψ
by directly evaluating U = func(H*dt)
, i.e. by matrix exponentiation for the default func
, and then multiplying U
and Ψ
in-place with mul!
.
The workspace wrk
must be initialized with ExpPropWrk
to provide storage for a temporary state.
Keyword arguments besides func
are ignored.
Private members
QuantumPropagators.ExpProp.expprop
— FunctionΨ_out = expprop(Ψ, H, dt, wrk; func=(H_dt -> exp(-1im * H_dt)), _...)
evaluates Ψ_out = func(H*dt) Ψ
as in expprop!
, but not acting in-place.
QuantumPropagators.SpectralRange
Chebychev propagation relies on estimating the spectral range of the Hamiltonian, which in turn may be done via Arnoldi iteration.
Reference
Public Members:
Private Members:
Public members
QuantumPropagators.SpectralRange.specrange
— FunctionCalculate the spectral range of a Hamiltonian H
on the real axis.
E_min, E_max = specrange(H; method=:auto, kwargs...)
calculates the approximate lowest and highest eigenvalues of H
. Any imaginary part in the eigenvalues is ignored: the routine is intended for (although not strictly limited to) a Hermitian H
.
This delegates to
specrange(H, method; kwargs...)
for the different methods.
The default method=:auto
chooses the best method for the given H
. This is :diag
for small matrices, and :arnoldi
otherwise. If both E_min
and E_max
are given in the kwargs
, those will be returned directly (method=:manual
).
Keyword arguments not relevant to the underlying implementation will be ignored.
E_min, E_max = specrange(
H, :arnoldi;
rng=Random.GLOBAL_RNG,
state=random_state(H; rng),
m_min=20,
m_max=60,
prec=1e-3,
norm_min=1e-15,
enlarge=true
)
uses Arnoldi iteration with state
as the starting vector. It approximates the eigenvalues of H
with between m_min
and m_max
Ritz values, until the lowest and highest eigenvalue are stable to a relative precision of prec
. The norm_min
parameter is passed to the underlying arnoldi!
.
If enlarge=true
(default) the returned E_min
and E_max
will be enlarged via a heuristic to slightly over-estimate the spectral radius instead of under-estimating it.
E_min, E_max = specrange(H, :diag)
uses exact diagonization via the standard eigvals
function to obtain the smallest and largest eigenvalue. This should only be used for relatively small matrices.
E_min, E_max = specrange(H, :manual; E_min, E_max)
directly returns the given E_min
and E_max
without considering H
.
Private members
QuantumPropagators.SpectralRange.ritzvals
— FunctionCalculate a vector for Ritz values converged to a given precision.
R = ritzvals(G, state, m_min, m_max=2*m_min; prec=1e-5, norm_min=1e-15)
calculates a complex vector R
of at least m_min
(assuming a sufficient Krylov dimension) and at most m_max
Ritz values.
QuantumPropagators.SpectralRange.random_state
— FunctionRandom normalized quantum state.
Ψ = random_state(H; rng=Random.GLOBAL_RNG)
returns a random normalized state compatible with the Hamiltonian H
. This is intended to provide a starting vector for estimating the spectral radius of H
via an Arnoldi method.
QuantumPropagators.Arnoldi
Arnoldi iteration is an approximate method to find extremal eigenvalues of a dynamical generator by projecting it into a Krylov subspace. It is used to estimate spectral ranges for Chebychev Propagation and to find appropriate Leja points to support Newton Propagation
Reference
Private Members:
Private members
QuantumPropagators.Arnoldi.extend_arnoldi!
— FunctionExtend dimension of Hessenberg matrix by one.
extend_arnoldi!(Hess, q, m, H, dt; norm_min=1e-15)
extends the entries in Hess
from size (m-1)×(m-1) to size m×m, and the list q
of Arnoldi vectors from m to (m+1). It is assumed that the input Hess
was created by a call to arnoldi!
with extended=false
or a previous call to extend_arnoldi!
. Note that Hess
itself is not resized, so it must be allocated to size m×m or greater on input.
QuantumPropagators.Arnoldi.diagonalize_hessenberg_matrix
— Functiondiagonalize_hessenberg_matrix(Hess, m; accumulate=false)
Diagonalize the m × m top left submatrix of the given Hessenberg matrix.
If accumulate
is true
, return the concatenated eigenvalues for Hess[1:1,1:1]
to Hess[1:m,1:m]
, that is, all sumatrices of size 1 through m
.
QuantumPropagators.Arnoldi.arnoldi!
— Functionm = arnoldi!(Hess, q, m, Ψ, H, dt=1.0; extended=true, norm_min=1e-15)
Calculate the Hessenberg matrix and Arnoldi vectors of H dt
, from Ψ
.
For a given order m
, the m×m
Hessemberg matrix is calculated and stored in in the pre-allocated Hess
. Further an array of m
normalized Arnoldi vectors is stored in in the pre-allocated q
, plus one additional unnormalized Arnoldi vector. The unnormalized m+1
st vector could be used to easily extend a given m×m
Hessenberg matrix to a (m+1)×(m+1)
matrix.
If the extended Hessenberg matrix is requested (extended=true
, default), the m+1
st Arnoldi vector is also normalized, and it's norm will be stored in m+1, m
entry of the (extended) Hessenberg matrix, which is an (m+1)×(m+1)
matrix.
Return the size m
of the calculated Hessenberg matrix. This will usually be the input m
, except when the Krylov dimension of H
starting from Ψ
is less then m
. E.g., if Ψ
is an eigenstate of H
, the returned m
will be 1.
See https://en.wikipedia.org/wiki/Arnoldi_iteration for a description of the algorithm.
Arguments
Hess::Matrix{ComplexF64}
: Pre-allocated storage for the Hessemberg matrix. Can be uninitialized on input. The matrix must be at least of sizem×m
, or(m+1)×(m+1)
ifextended=true
. On output, them×m
sub-matrix ofHess
(with the returned outputm
) will contain the Hessenberg matrix, and all other elements ofHess
be be set to zero.q
: Pre-allocated array of states similar toΨ
, as storage for the calculated Arnoldi vectors. These may be un-initialized on input. Must be at least of lengthm+1
m
: The requested dimensions of the output Hessenberg matrix.Ψ
: The starting vector for the Arnoldi procedure. This can be of any type, as long asΦ = H * Ψ
results in a vector similar toΨ
, there is an inner products ofΦ
andΨ
(Ψ⋅Φ
is defined), andnorm(Ψ)
is defined.H
: The operator (up todt
) for which to calculate the Arnoldi procedure. Can be of any type, as long asH * Ψ
is defined.dt
: The implicit time step; the total operator for which to calculate the Arnoldi procedure isH * dt
extended
: Iftrue
(default), calculate the extended Hessenberg matrix, and normalized the final Arnoldi vectornorm_min
: the minimum value of the norm ofΨ
at whichΨ
should be considered the zero vector
QuantumPropagators.Interfaces
The following routines allow to test whether custom data structures match the interface requirements of QuantumPropagators
.
Reference
Public Members:
check_amplitude
check_control
check_generator
check_operator
check_parameterized
check_parameterized_function
check_propagator
check_state
check_tlist
supports_inplace
Public members
QuantumPropagators.Interfaces.check_amplitude
— FunctionCheck amplitude appearing in Generator
.
@test check_amplitude(ampl; tlist, quiet=false)
verifies that the given ampl
is a valid element in the list of amplitudes
of a Generator
object. Specifically:
get_controls(ampl)
must be defined and return a tuple- all controls in
ampl
must passcheck_control
substitute(ampl, controls_replacements)
must be definedevaluate(ampl, tlist, n)
must be defined and return a Numberevaluate(ampl, tlist, n; vals_dict)
must be defined and return a Number
If for_parameterization
(may require the RecursiveArrayTools
package to be loaded):
get_parameters(ampl)
must be defined and return a vector of floats. Mutating that vector must mutate the controls inside theampl
.
The function returns true
for a valid amplitude and false
for an invalid amplitude. Unless quiet=true
, it will log an error to indicate which of the conditions failed.
QuantumPropagators.Interfaces.check_control
— FunctionCheck that control
can be evaluated on a time grid.
@test check_control(
control;
tlist,
for_parameterization=true,
for_time_continuous=(control isa Function),
quiet=false
)
verifies the given control
(one of the elements of the tuple returned by get_controls
):
evaluate(control, tlist, n)
must be defined and return aFloat64
evaluate(control, tlist, n; vals_dict=IdDict(control => v))
must be defined and returnv
discretize(control, tlist)
must be defined and return a vector of floats of the same size astlist
. Only iflength(tlist) > 2
.- all values in
discretize(control, tlist)
must be finite (isfinite
). discretize_on_midpoints(control, tlist)
must be defined and return a vector of floats with one element less thantlist
. Only iflength(tlist) > 2
.- all values in
discretize_on_midpoints(control, tlist)
must be finite (isfinite
)
If for_time_continuous
:
evaluate(control, t)
must be defined and return aFloat64
evaluate(control, t; vals_dict=IdDict(control => v))
must be defined and returnv
If for_parameterization
:
get_parameters(control)
must be defined and return a vector of floats. Mutating that vector must mutate the control.
The function returns true
for a valid control and false
for an invalid control. Unless quiet=true
, it will log an error to indicate which of the conditions failed.
QuantumPropagators.Interfaces.check_generator
— FunctionCheck the dynamical generator
for propagating state
over tlist
.
@test check_generator(
generator; state, tlist,
for_pwc=true, for_time_continuous=false,
for_expval=true, for_parameterization=false,
atol=1e-14, quiet=false)
verifies the given generator
:
get_controls(generator)
must be defined and return a tuple- all controls returned by
get_controls(generator)
must passcheck_control
substitute(generator, replacements)
must be defined- If
generator
is aGenerator
instance, all elements ofgenerator.amplitudes
must passcheck_amplitude
withfor_parameterization
.
If for_pwc
(default):
op = evaluate(generator, tlist, n)
must return a valid operator (check_operator
), with forwarded keyword arguments (includingfor_expval
)- If
QuantumPropagators.Interfaces.supports_inplace(op)
istrue
,evaluate!(op, generator, tlist, n)
must be defined
If for_time_continuous
:
evaluate(generator, t)
must return a valid operator (check_operator
), with forwarded keyword arguments (includingfor_expval
)If
QuantumPropagators.Interfaces.supports_inplace(op)
istrue
,evaluate!(op, generator, t)
must be defined
If for_parameterization
(may require the RecursiveArrayTools
package to be loaded):
get_parameters(generator)
must be defined and return a vector of floats. Mutating that vector must mutate the controls inside thegenerator
.
The function returns true
for a valid generator and false
for an invalid generator. Unless quiet=true
, it will log an error to indicate which of the conditions failed.
QuantumPropagators.Interfaces.check_operator
— FunctionCheck that op
is a valid operator that can be applied to state
.
@test check_operator(op; state, tlist=[0.0, 1.0],
for_expval=true, atol=1e-14, quiet=false)
verifies the given op
relative to state
. The state
must pass check_state
.
An "operator" is any object that evaluate
returns when evaluating a time-dependent dynamic generator. The specific requirements for op
are:
op
must not be time-dependent:evaluate(op, tlist, 1) ≡ op
op
must not contain any controls:length(get_controls(op)) == 0
op * state
must be defined- The
QuantumPropagators.Interfaces.supports_inplace
method must be defined forop
. If it returnstrue
, it must be possible to evaluate a generator in-place into the existingop
. Seecheck_generator
.
If QuantumPropagators.Interfaces.supports_inplace(state)
:
- The 3-argument
LinearAlgebra.mul!
must applyop
to the givenstate
- The 5-argument
LinearAlgebra.mul!
must applyop
to the givenstate
LinearAlgebra.mul!
must match*
, if applicableLinearAlgebra.mul!
must return the resulting state
If for_expval
(typically required for optimal control):
LinearAlgebra.dot(state, op, state)
must return return a numberdot(state, op, state)
must matchdot(state, op * state)
, if applicable
The function returns true
for a valid operator and false
for an invalid operator. Unless quiet=true
, it will log an error to indicate which of the conditions failed.
QuantumPropagators.Interfaces.check_parameterized
— FunctionCheck that that the object supports the parameterization interface.
@test check_parameterized(object; name="::$typeof(object))", quiet=false)
verifies that the given object
:
- can be passed to
get_parameters
, which must return anAbstractVector
ofFloat64
- is mutated by mutating the
parameters
obtained byget_parameters
See also
check_parameterized_function
isobject
is aParameterizedFunction
QuantumPropagators.Interfaces.check_parameterized_function
— FunctionCheck a ParameterizedFunction
instance.
@test check_parameterized_function(f; tlist; quiet=false)
verifies that the given f
:
- is an instance of
ParameterizedFunction
. - has a field
parameters
that is anAbstractVector{Float64}
. - is a callable as
f(t)
for values oft
intlist
, returning aFloat64
. get_parameters
provides access to theparameters
field.- passes
check_parameterized
See also
check_parameterized
for objects that have parameters (get_parameters
), but are not instances ofParameterizedFunction
QuantumPropagators.Interfaces.check_propagator
— FunctionCheck that the given propagator
implements the required interface.
@test check_propagator(propagator; atol=1e-14, quiet=false)
verifies that the propagator
matches the interface described for an AbstractPropagator
. The propagator
must have been freshly initialized with init_prop
.
propagator
must have the propertiesstate
,tlist
,t
,parameters
,backward
, andinplace
propagator.state
must be a valid state (seecheck_state
)- If
propagator.inplace
is true,supports_inplace
forpropagator.state
must also be true propagator.tlist
must be monotonically increasing.propagator.t
must be the first or last element ofpropagator.tlist
, depending onpropagator.backward
prop_step!(propagator)
must be defined and return a valid state until the time grid is exhausted- For an in-place propagator, the state returned by
prop_step!
must be thepropagator.state
object - For a not-in-place propagator, the state returned by
prop_step!
must be a new object prop_step!
must advancepropagator.t
forward or backward one step on the time gridprop_step!
must returnnothing
when going beyond the time gridset_t!(propagator, t)
must be defined and setpropagator.t
set_state!(propagator, state)
must be defined and setpropagator.state
.set_state!(propagator, state)
for an in-place propagator must overwritepropagator.state
in-place.set_state!
must return the setpropagator.state
- In a
PiecewisePropagator
,propagator.parameters
must be a dict mapping controls to a vector of values, one for each interval onpropagator.tlist
reinit_prop!
must be defined and re-initialize the propagatorreinit_prop!(propagator, state)
must be idempotent. That is, repeated calls toreinit_prop!
leave thepropagator
unchanged.
The function returns true
for a valid propagator and false
for an invalid propagator. Unless quiet=true
, it will log an error to indicate which of the conditions failed.
QuantumPropagators.Interfaces.check_state
— FunctionCheck that state
is a valid element of a Hilbert space.
@test check_state(state; normalized=false, atol=1e-15, quiet=false)
verifies the following requirements:
- The inner product (
LinearAlgebra.dot
) of two states must return a Complex number. - The
LinearAlgebra.norm
ofstate
must be defined via the inner product. This is the definition of a Hilbert space, a.k.a a "complete inner product space" or more precisely a "Banach space (normed vector space) where the norm is induced by an inner product". - The `QuantumPropagators.Interfaces.supports_inplace method must be defined for
state
Any state
must support the following not-in-place operations:
state + state
andstate - state
must be definedcopy(state)
must be defined and return an object of the same type asstate
c * state
for a scalarc
must be definednorm(state + state)
must fulfill the triangle inequalityzero(state)
must be defined and produce a state with norm 00.0 * state
must produce a state with norm 0copy(state) - state
must have norm 0norm(state)
must have absolute homogeneity:norm(s * state) = s * norm(state)
If supports_inplace(state)
is true
, the state
must also support the following:
similar(state)
must be defined and return a valid state of the same type astate
copyto!(other, state)
must be definedfill!(state, c)
must be definedLinearAlgebra.lmul!(c, state)
for a scalarc
must be definedLinearAlgebra.axpy!(c, state, other)
must be definednorm(state)
must fulfill the same general mathematical norm properties as for the non-in-place norm.
If normalized
(not required by default):
LinearAlgebra.norm(state)
must be 1
It is strongly recommended to always support immutable operations (also for mutable states)
The function returns true
for a valid state and false
for an invalid state. Unless quiet=true
, it will log an error to indicate which of the conditions failed.
QuantumPropagators.Interfaces.check_tlist
— FunctionCheck that the given tlist
is valid.
@test check_tlist(tlist; quiet=false)
verifies the given time grid. A valid time grid must
- be a
Vector{Float64}
, - contain at least two points (beginning and end),
- be monotonically increasing
The function returns true
for a valid time grid and false
for an invalid time grid. Unless quiet=true
, it will log an error to indicated which of the conditions failed.
QuantumPropagators.Interfaces.supports_inplace
— FunctionIndicate whether a given state or operator supports in-place operations
supports_inplace(state)
Indicates that propagators can assume that the in-place requirements defined in QuantumPropagators.Interfaces.check_state
hold. States with in-place support must also fulfill specific properties when interacting with operators, see QuantumPropagators.Interfaces.check_operator
.
supports_inplace(op)
Indicates that the operator can be evaluated in-place with evaluate!
, see QuantumPropagators.Interfaces.check_generator
Note that supports_inplace
is not quite the same as Base.ismutable
: When using custom structs for states or operators, even if those structs are not defined as mutable
, they may still define the in-place interface (typically because their components are mutable).
Index
QuantumPropagators.AbstractPropagator
QuantumPropagators.Amplitudes.LockedAmplitude
QuantumPropagators.Amplitudes.ShapedAmplitude
QuantumPropagators.Cheby.ChebyWrk
QuantumPropagators.ChebyPropagator
QuantumPropagators.Controls.ParameterizedFunction
QuantumPropagators.ExpProp.ExpPropWrk
QuantumPropagators.ExpPropagator
QuantumPropagators.Generators.Generator
QuantumPropagators.Generators.Operator
QuantumPropagators.Generators.ScaledOperator
QuantumPropagators.Newton.NewtonWrk
QuantumPropagators.NewtonPropagator
QuantumPropagators.PWCPropagator
QuantumPropagators.PiecewisePropagator
QuantumPropagators.Propagation
QuantumPropagators.Arnoldi.arnoldi!
QuantumPropagators.Arnoldi.diagonalize_hessenberg_matrix
QuantumPropagators.Arnoldi.extend_arnoldi!
QuantumPropagators.Cheby.cheby
QuantumPropagators.Cheby.cheby!
QuantumPropagators.Cheby.cheby_coeffs
QuantumPropagators.Cheby.cheby_coeffs!
QuantumPropagators.Controls.discretize
QuantumPropagators.Controls.discretize_on_midpoints
QuantumPropagators.Controls.evaluate
QuantumPropagators.Controls.evaluate!
QuantumPropagators.Controls.get_controls
QuantumPropagators.Controls.get_parameters
QuantumPropagators.Controls.get_tlist_midpoints
QuantumPropagators.Controls.substitute
QuantumPropagators.Controls.t_mid
QuantumPropagators.ExpProp.expprop
QuantumPropagators.ExpProp.expprop!
QuantumPropagators.Generators.hamiltonian
QuantumPropagators.Generators.liouvillian
QuantumPropagators.Interfaces.check_amplitude
QuantumPropagators.Interfaces.check_control
QuantumPropagators.Interfaces.check_generator
QuantumPropagators.Interfaces.check_operator
QuantumPropagators.Interfaces.check_parameterized
QuantumPropagators.Interfaces.check_parameterized_function
QuantumPropagators.Interfaces.check_propagator
QuantumPropagators.Interfaces.check_state
QuantumPropagators.Interfaces.check_tlist
QuantumPropagators.Interfaces.supports_inplace
QuantumPropagators.Newton.extend_leja!
QuantumPropagators.Newton.extend_newton_coeffs!
QuantumPropagators.Newton.newton!
QuantumPropagators.Shapes.blackman
QuantumPropagators.Shapes.box
QuantumPropagators.Shapes.flattop
QuantumPropagators.SpectralRange.random_state
QuantumPropagators.SpectralRange.ritzvals
QuantumPropagators.SpectralRange.specrange
QuantumPropagators.Storage.get_from_storage
QuantumPropagators.Storage.get_from_storage!
QuantumPropagators.Storage.init_storage
QuantumPropagators.Storage.map_observable
QuantumPropagators.Storage.map_observables
QuantumPropagators.Storage.write_to_storage!
QuantumPropagators.cheby_get_spectral_envelope
QuantumPropagators.disable_timings
QuantumPropagators.enable_timings
QuantumPropagators.init_prop
QuantumPropagators.init_prop
QuantumPropagators.init_prop
QuantumPropagators.init_prop
QuantumPropagators.init_prop
QuantumPropagators.ode_function
QuantumPropagators.prop_step!
QuantumPropagators.propagate
QuantumPropagators.propagate_sequence
QuantumPropagators.reinit_prop!
QuantumPropagators.set_state!
QuantumPropagators.set_t!
QuantumPropagators.timings_enabled