krotov.propagators module¶
Routines that can be passed as propagator to optimize_pulses()
The numerical effort involved in the optimization is almost entirely within the
simulation of the system dynamics. In every iteration and for every objective,
the system must be “propagated” once forwards in time and once backwards in
time, see also krotov.parallelization
.
The implementation of this time propagation must be inside the user-supplied
routine propagator that is passed to optimize_pulses()
and must
calculate the propagation over a single time step. In particular,
qutip.mesolve.mesolve()
is not automatically used for simulating any
dynamics within the optimization. The signature for any propagator
must be the same as the “reference” expm()
propagator:
>>> str(inspect.signature(krotov.propagators.expm))
'(H, state, dt, c_ops=None, backwards=False, initialize=False)'
The arguments are as follows (cf. Propagator
):
H is the system Hamiltonian or Liouvillian, in a nested-list format similar to that used by
qutip.mesolve.mesolve()
, e.g., for a Hamiltonian \(\Op{H} = \Op{H}_0 + c \Op{H}_1\), where \(c\) is the value of a control field at a particular point in time, propagator would receive a list[H0, [H1, c]]
whereH0
andH1
arequtip.Qobj
operators. The nested-list for H used here, with scalar values for the controls, is obtained internally from the format used bymesolve()
, with time-dependent controls over the entire time grid, viakrotov.conversions.plug_in_pulse_values()
.state is the
qutip.Qobj
state that should be propagated, either a Hilbert space state, or a density matrix.dt is the time step (a float). It is always positive, even for
backwards=True
.c_ops is None, or a list of collapse (Lindblad) operators, where each list element is a
qutip.Qobj
instance (or possibly a nested list, for time-dependent Lindblad operators. Note that is generally preferred for H to be a Liouvillian, for dissipative dynamics.backwards (
bool
): If passed as True, the propagator should propagate backwards in time. In Hilbert space, this means using -dt instead of dt. In Liouville space, there is no difference between forward and backward propagation. In the context of Krotov’s method, the backward propagation uses the conjugate Hamiltonian, respectively Liouvillian. However, the propagator routine does not need to be aware of this fact: it will receive the appropriate H and c_ops.initialize (
bool
): A flag to indicate the beginning of a propagation over a time grid. If False in subsequent calls, the propagator may assume that the input state is the result of the previous call to propagator.
Warning
The routines in this module are provided with no guarantee to be either
general or efficient. The expm()
propagator is exact to machine
precision, but generally extremely slow. For “production use”, it is
recommended to supply a problem-specific propagator that is highly
optimized for speed. You might consider the use of Cython. This is key to
minimize the runtime of the optimization.
The initialize flag enables “stateful” propagators that cache data between
calls. This can significantly improve numerical efficiency.
DensityMatrixODEPropagator
is an example for such a propagator. In
general, any stateful propagator should be an instance of
Propagator
.
Summary¶
Classes:
Propagator for density matrix evolution under a Lindbladian |
|
Abstract base class for stateful propagators |
Functions:
Propagate using matrix exponentiation. |
__all__
: DensityMatrixODEPropagator
, Propagator
, expm
Reference¶
-
krotov.propagators.
expm
(H, state, dt, c_ops=None, backwards=False, initialize=False)[source]¶ Propagate using matrix exponentiation.
This supports H being a Hamiltonian (for a Hilbert space state) or a Liouvillian (for state being a density matrix) in nested-list format. Collapse operators c_ops are not supported. The propagator is not stateful, thus initialize is ignored. The matrix exponentiation is evaluated in single-threaded mode, to prevent accidental nested parallelization.
-
class
krotov.propagators.
Propagator
[source]¶ Bases:
abc.ABC
Abstract base class for stateful propagators
-
abstract
__call__
(H, state, dt, c_ops=None, backwards=False, initialize=False)[source]¶ Evaluation of a single propagation step
- Parameters
H (list) – A Hamiltonian or Liouvillian in qutip’s nested-list format, with a scalar value in the place of a time-dependency. For example,
[H0, [H1, u]]
for a drift HamiltonianH0
, a control HamiltonianH1
, and a scalar valueu
that is a time-dependent control evaluated for a particular point in time.state (qutip.Qobj) – The state to propagate
dt (float) – The time step over which to propagate
c_ops (list or None) – A list of Lindblad operators. Using explicit Lindblad operators should be avoided: it is usually more efficient to convert them into a Lindbladian, passed as H
backwards (bool) – Whether the propagation is forward in time or backward in time
initialize (bool) – Whether the propagator should (re-)initialize for a new propagation, when the propagator is used to advance on a time grid, initialize should be passed as True for the initial time step (0 to dt in a forward propagation, or T to T-dt for a backward propagation), and False otherwise.
Note
A propagator may assume the propagation to be “sequential” when initialize is False. That is, the state to propagate is the result of the previous call to the propagator.
-
abstract
-
class
krotov.propagators.
DensityMatrixODEPropagator
(method='adams', order=12, atol=1e-08, rtol=1e-06, nsteps=1000, first_step=0, min_step=0, max_step=0, reentrant=False)[source]¶ Bases:
krotov.propagators.Propagator
Propagator for density matrix evolution under a Lindbladian
See
qutip.solver.Options
for all arguments except reentrant. Passing True for the reentrant re-initializes the propagator in every time step.Warning
By default, the propagator is not “re-entrant”. That is, you cannot use more than one instance of
DensityMatrixODEPropagator
in the same process at the same time. This limitation is due toscipy.integrate.ode
with the “zvode” integrator not being re-entrant. Passingreentrant=True
side-steps this problem by re-initializatingscipy.integrate.ode
in every time step. This makes it possible to useDensityMatrixODEPropagator
in the optimization of multiple objectives, but creates a significant overhead.-
__call__
(H, state, dt, c_ops=None, backwards=False, initialize=False)[source]¶ Evaluation of a single propagation step
- Parameters
H (list) – A Liouvillian superoperator in qutip’s nested-list format, with a scalar value in the place of a time-dependency. For example,
[L0, [L1, u]]
for a drift LiouvillianL0
, a control LiouvillianH1
, and a scalar valueu
that is a time-dependent control evaluated for a particular point in time. If initialize is False, only the control values are taken into account; any operators are assumed to be identical to the internally cached values of H during initialization.state (qutip.Qobj) – The density matrix to propagate. The passed value is ignored unless initialize is given as True. Otherwise, it is assumed that state matches the (internally stored) state that was the result from the previous propagation step.
dt (float) – The time step over which to propagate
c_ops (list or None) – An empty list, or None. Since this propagator assumes a full Liouvillian, it cannot be combined with Lindblad operators.
backwards (bool) – Whether the propagation is forward in time or backward in time. Since the equation of motion for a Liouvillian and conjugate Liouvillian is the same, this parameter has no effect. Instead, for the backward propagation, the conjugate Liouvillian must be passed for L.
initialize (bool) – Whether to (re-)initialize for a new propagation. This caches H (except for the control values) and state internally.
-