krotov.conversions module¶
Routines for structural conversions.
Conversion between between data structures used by QuTiP’s
mesolve()
and data structures used internally in an
optimization with Krotov’s method. This includes the
time discretization of control fields, and in particular converting between a
discretization defined on the points of the time grid (“controls”) and
piecewise-constant “pulses” defined on the intervals of the time grid.
Summary¶
Functions:
Convert control on time grid to control on time grid intervals. |
|
Discretize the given control onto the tlist time grid. |
|
Extract a list of (unique) controls from the objectives. |
|
Extract a map of where controls are used in objectives. |
|
Plug pulse values into H. |
|
Convert pulse from time-grid intervals to time-grid points. |
|
Convert pulse_options into a list. |
__all__
: control_onto_interval
, discretize
, extract_controls
, extract_controls_mapping
, plug_in_pulse_values
, pulse_onto_tlist
, pulse_options_dict_to_list
Reference¶
-
krotov.conversions.
control_onto_interval
(control)[source]¶ Convert control on time grid to control on time grid intervals.
- Parameters
control (numpy.ndarray) – values of controls on time grid
- Returns
pulse defined on the intervals of the time grid
- Return type
The value for the first and last interval will be identical to the values at
control[0]
andcontrol[-1]
to ensure proper boundary conditions. All other intervals are calculated such that the original values in control are the average of the interval-values before and after that point in time.The
pulse_onto_tlist()
function calculates the inverse to this transformation.Note
For a callable control, call
discretize()
first.
-
krotov.conversions.
discretize
(control, tlist, args=(None), kwargs=None, via_midpoints=False)[source]¶ Discretize the given control onto the tlist time grid.
If control is a callable, return array of values for control evaluated for all points in tlist. If control is already discretized, check that the discretization matches tlist (by size).
- Parameters
control (callable or numpy.ndarray) – control to be discretized. If callable, must take time value t as its first argument.
tlist (numpy.ndarray) – time grid to discretize one
args (tuple or list) – If control is a callable, further positional arguments to pass to control. The default passes a single value None, to match the requirements for a callable control function in QuTiP.
kwargs (None or dict) – If control is callable, further keyword arguments to pass to control. If None, no keyword arguments will be passed.
via_midpoints (bool) – If True, sample control at the midpoints of tlist (except for the initial and final values which are evaluated at
tlist[0]
andtlist[1]
to preseve exact boundary conditions) and then un-average viapulse_onto_tlist()
to fit onto tlist. If False, evaluate directly on tlist.
Note
If
via_midpoints=True
, the discretized values are generally not exactly the result of evaluating control at the values of tlist. Instead, the values are adjusted slightly to guarantee numerical stability when converting between a sampling on the time grid and a sampling on the mid points of the time grid intervals, as required by Krotov’s method, see Time discretization.- Returns
Discretized array of real control values, same length as tlist
- Return type
- Raises
TypeError – If control is not a function that takes two arguments (t, args), or a numpy array
ValueError – If control is numpy array of incorrect size.
-
krotov.conversions.
extract_controls
(objectives)[source]¶ Extract a list of (unique) controls from the objectives.
Controls are unique if they are not the same object, cf. Python’s is keyword.
See
extract_controls_mapping()
for an example.
-
krotov.conversions.
extract_controls_mapping
(objectives, controls)[source]¶ Extract a map of where controls are used in objectives.
The result is a nested list where the first index relates to the objectives, the second index relates to the Hamiltonian (0) or the c_ops (1…), and the third index relates to the controls.
Example
>>> import qutip >>> import krotov >>> X, Y, Z = qutip.Qobj(), qutip.Qobj(), qutip.Qobj() # dummy Hams >>> u1, u2 = np.array([]), np.array([]) # dummy controls >>> psi0, psi_tgt = qutip.Qobj(), qutip.Qobj() # dummy states
>>> H1 = [X, [Y, u1], [Z, u1]] # ham for first objective >>> H2 = [X, [Y, u2]] # ham for second objective >>> c_ops = [[[X, u1]], [[Y, u2]]] >>> objectives = [ ... krotov.Objective( ... initial_state=psi0, ... target=psi_tgt, ... H=H1, ... c_ops=c_ops ... ), ... krotov.Objective( ... initial_state=psi0, ... target=psi_tgt, ... H=H2, ... c_ops=c_ops ... ) ... ] >>> controls = extract_controls(objectives) >>> assert controls == [u1, u2]
>>> controls_mapping = extract_controls_mapping(objectives, controls) >>> controls_mapping [[[[1, 2], []], [[0], []], [[], [0]]], [[[], [1]], [[0], []], [[], [0]]]]
The structure should be read as follows:
For the first objective (0), in the Hamiltonian (0), where is the first pulse (0) used? (answer: in
H1[1]
andH1[2]
)>>> controls_mapping[0][0][0] [1, 2]
For the second objective (1), in the second
c_ops
(2), where is the second pulse (1) used? (answer: inc_ops[1][0]
)>>> controls_mapping[1][2][1] [0]
For the second objective (1), in the Hamiltonian (0), where is the first pulse (0) used? (answer: nowhere)
>>> controls_mapping[1][0][0] []
-
krotov.conversions.
plug_in_pulse_values
(H, pulses, mapping, time_index, conjugate=False)[source]¶ Plug pulse values into H.
- Parameters
H (list) – nested list for a QuTiP-time-dependent operator
pulses (list) – list of pulses in array format
mapping (list) – nested list: for each pulse, a list of indices in H where pulse value should be inserted
time_index (int) – Index of the value of each pulse that should be plugged in
conjugate (bool) – If True, use conjugate complex pulse values
- Returns
a list with the same structure as H that contains the same
Qobj
operators as H, but where every time dependency is replaced by the value of the appropriate pulse at time_index.- Return type
Example
>>> X, Y, Z = 'X', 'Y', 'Z' # dummy Hams, these would normally be Qobjs >>> u1, u2 = np.array([0, 10, 0]), np.array([0, 20, 0]) >>> H = [X, [X, u1], [Y, u1], [Z, u2]] >>> pulses = [u1, u2] >>> mapping = [[1, 2], [3]] # u1 is in H[1] and H[2], u2 is in H[3] >>> plug_in_pulse_values(H, pulses, mapping, time_index=1) ['X', ['X', 10], ['Y', 10], ['Z', 20]]
Note
It is of no consequence whether H contains the pulses, as long as it has the right structure:
>>> H = [X, [X, None], [Y, None], [Z, None]] >>> plug_in_pulse_values(H, pulses, mapping, time_index=1) ['X', ['X', 10], ['Y', 10], ['Z', 20]]
-
krotov.conversions.
pulse_onto_tlist
(pulse)[source]¶ Convert pulse from time-grid intervals to time-grid points.
- Parameters
pulse (numpy.ndarray) – values defined on the interval of a time grid
- Returns
values of the control defined directly on the time grid points. The size of the returned array is one greater than the size of pulse.
- Return type
Inverse of
control_onto_interval()
.The first and last value are also the first and last value of the returned control field. For all other points, the value is the average of the value of the input values before and after the point.
-
krotov.conversions.
pulse_options_dict_to_list
(pulse_options, controls)[source]¶ Convert pulse_options into a list.
Given a dict pulse_options that contains an options-dict for every control in controls (cf.
optimize_pulses()
), return a list of the options-dicts in the same order as controls.- Raises
ValueError – if pulse_options to not contain all of the controls