krotov.propagators module

Routines that can be passed as propagator to optimize_pulses()

Summary

Data:

DensityMatrixODEPropagator Propagator for density matrix evolution under a Lindbladian
Propagator Abstract base class for stateful propagators
expm Propagate using matrix exponentiation

__all__: DensityMatrixODEPropagator, Propagator, expm

Reference

krotov.propagators.expm(H, state, dt, c_ops, backwards=False, initialize=False)[source]

Propagate using matrix exponentiation

class krotov.propagators.Propagator[source]

Bases: abc.ABC

Abstract base class for stateful propagators

__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 Hamiltonian H0, a control Hamiltonian H1, and a scalar value u that is a time-dependent control evaluated for a particular point in time.
  • state (qutip.Qobj) – The state to propagation
  • 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.

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)[source]

Bases: krotov.propagators.Propagator

Propagator for density matrix evolution under a Lindbladian

See qutip.solver.Options for arguments.

__call__(L, rho, dt, c_ops=None, backwards=False, initialize=False)[source]

Evaluation of a single propagation step

Parameters:
  • L (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 Liouvillian L0, a control Liouvillian H1, and a scalar value u 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 L during initialization.
  • rho (qutip.Qobj) – The density matrix to propagate. The passed value is ignored unless initialize is given as True. Otherwise, it is assumed that rho 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 L (except for the control values) and rho internally.