krotov.functionals module

Functionals and chi_constructor routines.

Any chi_constructor routine passed to optimize_pulses() must take the following keyword-arguments:

  • fw_states_T (list of Qobj): The list of states resulting from the forward-propagation of each Objective.initial_state under the guess pulses of the current iteration (the optimized pulses of the previous iteration)

  • objectives (list of Objective): A list of the optimization objectives.

  • tau_vals (list of complex or None): The overlaps of the Objective.target and the corresponding fw_states_T, assuming Objective.target contains a quantum state. If the objective defines no target state, a list of Nones

Krotov’s method does not have an explicit dependence on the optimization functional. It only enters through the chi_constructor which calculates the boundary condition for the backward propagation, that is, the states

\[\Ket{\chi_k^{(i)}(T)} = - \left.\frac{\partial J_T} {\partial \Bra{\phi_k}} \right\vert_{\phi^{(i)}(T)}\]

for functionals defined in Hilbert space, or

\[\Op{\chi}_k^{(i)}(T) = - \left.\frac{\partial J_T} {\partial \langle\!\langle\Op{\rho}_k\vert} \right\vert_{\rho^{(i)}(T)}\]

in Liouville space, using the abstract Hilbert-Schmidt notation \(\langle\!\langle a \vert b \rangle\!\rangle \equiv \tr[a^\dagger b]\). Passing a specific chi_constructor results in the minimization of the final time functional from which that chi_constructor was derived.

The functions in this module that evaluate functionals are intended for use inside a function that is passed as an info_hook to optimize_pulses(). Thus, they calculate \(J_T\) from the same keyword arguments as the info_hook. The values for \(J_T\) may be used in a convergence analysis, see krotov.convergence.

Summary

Functions:

F_avg

Average gate fidelity

F_re

Real-part fidelity

F_sm

Square-modulus fidelity

F_ss

State-to-state phase-insensitive fidelity

J_T_hs

Hilbert-Schmidt distance measure functional \(J_{T,\text{hs}}\)

J_T_re

Real-part functional \(J_{T,\text{re}}\)

J_T_sm

Square-modulus functional \(J_{T,\text{sm}}\)

J_T_ss

State-to-state phase-insensitive functional \(J_{T,\text{ss}}\)

chis_hs

States \(\Op{\chi}_k\) for functional \(J_{T,\text{hs}}\)

chis_re

States \(\ket{\chi_k}\) for functional \(J_{T,\text{re}}\)

chis_sm

States \(\ket{\chi_k}\) for functional \(J_{T,\text{sm}}\)

chis_ss

States \(\ket{\chi_k}\) for functional \(J_{T,\text{ss}}\)

f_tau

Average complex overlaps of the target states with the fw_states_T.

gate

Gate that maps basis_states to fw_states_T

mapped_basis

Result of applying the gate O to basis_states

__all__: F_avg, F_re, F_sm, F_ss, J_T_hs, J_T_re, J_T_sm, J_T_ss, chis_hs, chis_re, chis_sm, chis_ss, f_tau, gate, mapped_basis

Reference

krotov.functionals.f_tau(fw_states_T, objectives, tau_vals=None, **kwargs)[source]

Average complex overlaps of the target states with the fw_states_T.

That is,

\[f_{\tau} = \frac{1}{N} \sum_{k=1}^{N} w_k \tau_k\]

where \(\tau_k\) are the elements of tau_vals, assumed to be

\[\tau_k = \Braket{\Psi_k^{\tgt}}{\Psi_k(T)},\]

in Hilbert space, or

\[\tau_k = \tr\left[\Op{\rho}_k^{\tgt\,\dagger}\Op{\rho}_k(T)\right]\]

in Liouville space, where \(\ket{\Psi_k}\) or \(\Op{\rho}_k\) are the elements of fw_states_T, and \(\ket{\Psi_k^{\tgt}}\) or \(\Op{\rho}^{\tgt}\) are the target states from the target attribute of the objectives. If tau_vals are None, they will be calculated internally.

\(N\) is the number of objectives, and \(w_k\) is an optional weight for each objective. For any objective that has a (custom) weight attribute, the \(w_k\) is taken from that attribute; otherwise, \(w_k = 1\). The weights, if present, are not automatically normalized, they are assumed to have values such that the resulting \(f_{\tau}\) lies in the unit circle of the complex plane. Usually, this means that the weights should sum to \(N\). The exception would be for mixed target states, where the weights should compensate for the non-unit purity. The problem may be circumvented by using J_T_hs() for mixed target states.

The kwargs are ignored, allowing the function to be used in an info_hook.

krotov.functionals.F_ss(fw_states_T, objectives, tau_vals=None, **kwargs)[source]

State-to-state phase-insensitive fidelity

\[F_{\text{ss}} = \frac{1}{N} \sum_{k=1}^{N} w_k \Abs{\tau_k}^2 \quad\in [0, 1]\]

with \(N\), \(w_k\) and \(\tau_k\) as in f_tau().

The kwargs are ignored, allowing the function to be used in an info_hook.

krotov.functionals.J_T_ss(fw_states_T, objectives, tau_vals=None, **kwargs)[source]

State-to-state phase-insensitive functional \(J_{T,\text{ss}}\)

\[J_{T,\text{ss}} = 1 - F_{\text{ss}} \quad\in [0, 1].\]

All arguments are passed to F_ss().

krotov.functionals.chis_ss(fw_states_T, objectives, tau_vals)[source]

States \(\ket{\chi_k}\) for functional \(J_{T,\text{ss}}\)

\[\Ket{\chi_k} = -\frac{\partial J_{T,\text{ss}}}{\partial \bra{\Psi_k(T)}} = \frac{1}{N} w_k \tau_k \Ket{\Psi^{\tgt}_k}\]

with \(\tau_k\) and \(w_k\) as defined in f_tau().

krotov.functionals.F_sm(fw_states_T, objectives, tau_vals=None, **kwargs)[source]

Square-modulus fidelity

\[F_{\text{sm}} = \Abs{f_{\tau}}^2 \quad\in [0, 1].\]

All arguments are passed to f_tau() to evaluate \(f_{\tau}\).

krotov.functionals.J_T_sm(fw_states_T, objectives, tau_vals=None, **kwargs)[source]

Square-modulus functional \(J_{T,\text{sm}}\)

\[J_{T,\text{sm}} = 1 - F_{\text{sm}} \quad\in [0, 1]\]

All arguments are passed to f_tau() while evaluating \(F_{\text{sm}}\) in F_sm().

krotov.functionals.chis_sm(fw_states_T, objectives, tau_vals)[source]

States \(\ket{\chi_k}\) for functional \(J_{T,\text{sm}}\)

\[\Ket{\chi_k} = -\frac{\partial J_{T,\text{sm}}}{\partial \bra{\Psi_k(T)}} = \frac{1}{N^2} w_k \sum_{j}^{N} w_j\tau_j\Ket{\Psi^{\tgt}_k}\]

with optional weights \(w_k\), cf. f_tau() (default: \(w_k=1\)). If given, the weights should generally sum to \(N\).

krotov.functionals.F_re(fw_states_T, objectives, tau_vals=None, **kwargs)[source]

Real-part fidelity

\[\begin{split}F_{\text{re}} = \Re[f_{\tau}] \quad\in \begin{cases} [-1, 1] & \text{in Hilbert space} \\ [0, 1] & \text{in Liouville space.} \end{cases}\end{split}\]

All arguments are passed to f_tau() to evaluate \(f_{\tau}\).

krotov.functionals.J_T_re(fw_states_T, objectives, tau_vals=None, **kwargs)[source]

Real-part functional \(J_{T,\text{re}}\)

\[\begin{split}J_{T,\text{re}} = 1 - F_{\text{re}} \quad\in \begin{cases} [0, 2] & \text{in Hilbert space} \\ [0, 1] & \text{in Liouville space.} \end{cases}\end{split}\]

All arguments are passed to f_tau() while evaluating \(F_{\text{re}}\) in F_re().

Note

If the target states are mixed, \(J_{T,\text{re}}\) may take negative values (for fw_states_T that are “in the right direction”, but more pure than the target states). In this case, you may consider using J_T_hs().

krotov.functionals.chis_re(fw_states_T, objectives, tau_vals)[source]

States \(\ket{\chi_k}\) for functional \(J_{T,\text{re}}\)

\[\Ket{\chi_k} = -\frac{\partial J_{T,\text{re}}}{\partial \bra{\Psi_k(T)}} = \frac{1}{2N} w_k \Ket{\Psi^{\tgt}_k}\]

with optional weights \(w_k\), cf. f_tau() (default: \(w_k=1\)). If given, the weights should generally sum to \(N\).

Note: tau_vals are ignored, but are present to satisfy the requirments of the chi_constructor interface.

krotov.functionals.J_T_hs(fw_states_T, objectives, tau_vals=None, **kwargs)[source]

Hilbert-Schmidt distance measure functional \(J_{T,\text{hs}}\)

\[\begin{split}J_{T,\text{hs}} = \frac{1}{2N} \sum_{k=1}^{N} w_k \Norm{\Op{\rho}_k(T) - \Op{\rho}_k^{\tgt}}_{\text{hs}}^2 \quad \in \begin{cases} [0, 2] & \text{in Hilbert space} \\ [0, 1] & \text{in Liouville space} \end{cases}\end{split}\]

in Liouville space (using the Hilbert-Schmidt norm), or equivalently with \(\ket{\Psi_k(T)}\) and \(\ket{\Psi_k^{tgt}}\) in Hilbert space. The functional is evaluated as

\[J_{T,\text{hs}} = \frac{1}{2N} \sum_{k=1}^{N} w_k \left( \Norm{\Op{\rho}_k(T)}_{\text{hs}}^2 + \Norm{\Op{\rho}^{\tgt}}_{\text{hs}}^2 - 2 \Re[\tau_k] \right)\]

where the \(\Op{\rho}_k\) are the elements of fw_states_T, the \(\Op{\rho}_k^{\tgt}\) are the target states from the target attribute of the objectives, and the \(\tau_k\) are the elements of tau_vals (which will be calculated internally if passed as None).

The \(w_k\) are optional weights, cf. f_tau(). If given, the weights should generally sum to \(N\).

The kwargs are ignored, allowing the function to be used in an info_hook.

Note

For pure states (or Hilbert space states), \(J_{T,\text{hs}}\) is equivalent to \(J_{T,\text{re}}\), cf. J_T_re(). However, the backward-propagated states \(\chi_k\) obtained from the two functionals (chis_re() and chis_hs()) are not equivalent. This may result in a vastly different optimization landscape that requires a significantly different value of the \(\lambda_a\) value that regulates the overall magnitude of the pulse updates (given in pulse_options in optimize_pulses()).

krotov.functionals.chis_hs(fw_states_T, objectives, tau_vals)[source]

States \(\Op{\chi}_k\) for functional \(J_{T,\text{hs}}\)

\[\Op{\chi}_k = -\frac{\partial J_{T,\text{sm}}} {\partial \langle\!\langle \Op{\rho}_k(T)\vert} = \frac{1}{2N} w_k \left(\Op{\rho}^{\tgt}_k - \Op{\rho}_k(T)\right)\]

with optional weights \(w_k\), cf. f_tau() (default: \(w_k=1\)).

This is derived from \(J_{T,\text{hs}}\) rewritten in the abstract Hilbert-Schmidt notation \(\langle\!\langle a \vert b \rangle\!\rangle \equiv \tr[a^\dagger b]\):

\[J_{T,\text{hs}} = \frac{-1}{2N} \sum_{k=1}^{N} w_k \big( \underbrace{ \langle\!\langle \Op{\rho}_k(T) \vert \Op{\rho}_k^{\tgt} \rangle\!\rangle + \langle\!\langle \Op{\rho}_k^{\tgt}\vert \Op{\rho}_k(T) \rangle\!\rangle }_{=2\Re[\tau_k]} - \underbrace{ \langle\!\langle \Op{\rho}_k(T) \vert \Op{\rho}_k(T) \rangle\!\rangle }_{=\Norm{\Op{\rho}_k(T)}_{\text{hs}}^2} - \underbrace{ \langle\!\langle \Op{\rho}_k^{\tgt} \vert \Op{\rho}_k^{\tgt} \rangle\!\rangle }_{=\Norm{\Op{\rho}^{\tgt}}_{\text{hs}}^2} \big).\]

Note: tau_vals are ignored, but are present to satisfy the requirments of the chi_constructor interface.

krotov.functionals.F_avg(fw_states_T, basis_states, gate, mapped_basis_states=None, prec=1e-05)[source]

Average gate fidelity

\[F_{\text{avg}} = \int \big\langle \Psi \big\vert \Op{O}^\dagger \DynMap[\ketbra{\Psi}{\Psi}] \Op{O} \big\vert \Psi \big\rangle \dd \Psi\]

where \(\Op{O}\) is the target gate, and \(\DynMap\) represents the dynamical map from time zero to \(T\).

In Liouville space, this is numerically evaluated as

\[F_{\text{avg}} = \frac{1}{N (N+1)} \sum_{i,j=1}^N \left( \big\langle \phi_i \big\vert \Op{O}^\dagger \Op{\rho}_{ij} \Op{O} \big\vert \phi_j \big\rangle + \big\langle \phi_i \big\vert \Op{O}^\dagger \Op{\rho}_{jj} \Op{O} \big\vert \phi_i \big\rangle \right),\]

where \(\ket{\phi_i}\) is the \(i\)’th element of basis_states, and \(\Op{\rho}_{ij}\) is the \((i-1) N + j\)’th element of fw_states_T, that is, \(\Op{\rho}_{ij} = \DynMap[\ketbra{\phi_i}{\phi_j}]\), with \(N\) the dimension of the Hilbert space.

In Hilbert space (unitary dynamics), this simplifies to

\[F_{\text{avg}} = \frac{1}{N (N+1)} \left( \Abs{\tr\left[\Op{O}^\dagger \Op{U}\right]}^2 + \tr\left[\Op{O}^\dagger \Op{U} \Op{U}^\dagger \Op{O}\right] \right),\]

where \(\Op{U}\) the gate that maps basis_states to the result of a forward propagation of those basis states, stored in fw_states_T.

Parameters
  • fw_states_T (list[qutip.Qobj]) – The forward propagated states. For dissipative dynamics, this must be the forward propagation of the full basis of Liouville space, that is, all \(N^2\) dyadic combinations of the Hilbert space logical basis states. For unitary dynamics, the \(N\) forward-propagated basis_states.

  • basis_states (list[qutip.Qobj]) – The \(N\) Hilbert space logical basis states

  • gate (qutip.Qobj) – The \(N \times N\) quantum gate in the logical subspace, e.g. qutip.qip.gates.cnot().

  • mapped_basis_states (None or list[qutip.Qobj]) – If given, the result of applying gate to basis_states. If not given, this will be calculated internally via mapped_basis(). It is recommended to pass pre-calculated mapped_basis_states when evaluating \(F_{\text{avg}}\) repeatedly for the same target.

  • prec (float) – assert that the fidelity is correct at least up to the given precision. Mathematically, \(F_{\text{avg}}\) is a real value. However, errors in the fw_states_T can lead to a small non-zero imaginary part. We assert that this imaginary part is below prec.

krotov.functionals.gate(basis_states, fw_states_T)[source]

Gate that maps basis_states to fw_states_T

Example

>>> from qutip import ket
>>> basis = [ket(nums) for nums in [(0, 0), (0, 1), (1, 0), (1, 1)]]
>>> fw_states_T = mapped_basis(qutip.gates.cnot(), basis)
>>> U = gate(basis, fw_states_T)
>>> assert (U - qutip.gates.cnot()).norm() < 1e-15
krotov.functionals.mapped_basis(O, basis_states)[source]

Result of applying the gate O to basis_states

Example

>>> from qutip import ket
>>> basis = [ket(nums) for nums in [(0, 0), (0, 1), (1, 0), (1, 1)]]
>>> states = mapped_basis(qutip.gates.cnot(), basis)
>>> assert (states[0] - ket((0, 0))).norm() < 1e-15
>>> assert (states[1] - ket((0, 1))).norm() < 1e-15
>>> assert (states[2] - ket((1, 1))).norm() < 1e-15  # swap (1, 1) ...
>>> assert (states[3] - ket((1, 0))).norm() < 1e-15  # ... and (1, 0)