Solver#

class Solver(static_hamiltonian=None, hamiltonian_operators=None, static_dissipators=None, dissipator_operators=None, hamiltonian_channels=None, dissipator_channels=None, channel_carrier_freqs=None, dt=None, rotating_frame=None, in_frame_basis=False, array_library=None, vectorized=None, rwa_cutoff_freq=None, rwa_carrier_freqs=None, validate=True)[source]#

Bases: object

Solver class for simulating both Hamiltonian and Lindblad dynamics, with high level type-handling of input states.

If only Hamiltonian information is provided, this class will internally construct a HamiltonianModel instance, and simulate the model using the Schrodinger equation \(\dot{y}(t) = -iH(t)y(t)\) (see the Solver.solve() method documentation for details on how different initial state types are handled). HamiltonianModel represents a decomposition of the Hamiltonian of the form:

\[H(t) = H_0 + \sum_i s_i(t) H_i,\]

where \(H_0\) is the static component, the \(H_i\) are the time-dependent components of the Hamiltonian, and the \(s_i(t)\) are the time-dependent signals, specifiable as either Signal objects, or constructed from Qiskit Pulse schedules if Solver is configured for Pulse simulation (see below).

If dissipators are specified as part of the model, then a LindbladModel is constructed, and simulations are performed by solving the Lindblad equation:

\[\dot{y}(t) = -i[H(t), y(t)] + \mathcal{D}_0(y(t)) + \mathcal{D}(t)(y(t)),\]

where \(H(t)\) is the Hamiltonian part, specified as above, and \(\mathcal{D}_0\) and \(\mathcal{D}(t)\) are the static and time-dependent portions of the dissipator, given by:

\[\mathcal{D}_0(y(t)) = \sum_j N_j y(t) N_j^\dagger - \frac{1}{2} \{N_j^\dagger N_j, y(t)\},\]

and

\[\mathcal{D}(t)(y(t)) = \sum_j \gamma_j(t) L_j y(t) L_j^\dagger - \frac{1}{2} \{L_j^\dagger L_j, y(t)\},\]

with \(N_j\) the static dissipators, \(L_j\) the time-dependent dissipator operators, and \(\gamma_j(t)\) the time-dependent signals specifiable as either Signal objects, or constructed from Qiskit Pulse schedules if Solver is configured for Pulse simulation (see below).

Transformations on the model can be specified via the optional arguments:

  • rotating_frame: Transforms the model into a rotating frame. Note that the operator specifying the frame will be substracted from the static_hamiltonian. If supplied as a 1d array, rotating_frame is interpreted as the diagonal elements of a diagonal matrix. Given a frame operator \(F = -i H_0\), for the Schrodinger equation entering the rotating frame of \(F\), corresponds to transforming the solution as \(y(t) \mapsto exp(-tF)y(t)\), and for the Lindblad equation it corresponds to transforming the solution as \(y(t) \mapsto exp(-tF)y(t)exp(tF)\). See RotatingFrame for more details.

  • in_frame_basis: Whether to represent the model in the basis in which the frame operator is diagonal, henceforth called the “frame basis”. If rotating_frame is None or was supplied as a 1d array, this kwarg has no effect. If rotating_frame was specified as a 2d array, the frame basis is the diagonalizing basis supplied by np.linalg.eigh. If in_frame_basis==True, this objects behaves as if all operators were supplied in the frame basis: calls to solve will assume the initial state is supplied in the frame basis, and the results will be returned in the frame basis. If in_frame_basis==False, the system will still be solved in the frame basis for efficiency, however the initial state (and final output states) will automatically be transformed into (and, respectively, out of) the frame basis.

  • rwa_cutoff_freq and rwa_carrier_freqs: Performs a rotating wave approximation (RWA) on the model with cutoff frequency rwa_cutoff_freq, assuming the time-dependent coefficients of the model have carrier frequencies specified by rwa_carrier_freqs. If dissipator_operators is None, rwa_carrier_freqs must be a list of floats of length equal to hamiltonian_operators, and if dissipator_operators is not None, rwa_carrier_freqs must be a tuple of lists of floats, with the first entry the list of carrier frequencies for hamiltonian_operators, and the second entry the list of carrier frequencies for dissipator_operators. See rotating_wave_approximation() for details on the mathematical approximation.

Additionally, the array_library argument controls the underlying array representation used to store and evaluate the model. See the model evaluation section of the Models API documentation for a more detailed description of this argument.

Note

When using the rwa_cutoff_freq optional argument, Solver cannot be instantiated within a JAX-transformable function. However, after construction, instances can still be used within JAX-transformable functions regardless of whether an rwa_cutoff_freq is set.

Solver can be configured to simulate Qiskit Pulse schedules by setting all of the following parameters, which determine how Pulse schedules are interpreted:

  • hamiltonian_channels: List of channels in string format corresponding to the time-dependent coefficients of hamiltonian_operators.

  • dissipator_channels: List of channels in string format corresponding to time-dependent coefficients of dissipator_operators.

  • channel_carrier_freqs: Dictionary mapping channel names to frequencies. A frequency must be specified for every channel appearing in hamiltonian_channels and dissipator_channels. When simulating schedules, these frequencies are interpreted as the analog carrier frequencies associated with the channel; deviations from these frequencies due to SetFrequency or ShiftFrequency instructions are implemented by digitally modulating the samples for the channel envelope. If an rwa_cutoff_freq is specified, and no rwa_carrier_freqs is specified, these frequencies will be used for the RWA.

  • dt: The envelope sample width.

If configured to simulate Pulse schedules, and a JAX-based solver method is chosen when calling Solver.solve(), Solver.solve() will automatically attempt to compile a single function to re-use for all schedule simulations.

The evolution given by the model can be simulated by calling Solver.solve(), which calls solve_lmde(), and does various automatic type handling operations for qiskit.quantum_info state and super operator types.

Initialize solver with model information.

Parameters:
  • static_hamiltonian (Union[ndarray, number, int, float, complex, Tracer, Array, Array, spmatrix, BCOO, list, None]) – Constant Hamiltonian term. If a rotating_frame is specified, the frame_operator will be subtracted from the static_hamiltonian.

  • hamiltonian_operators (Union[ndarray, number, int, float, complex, Tracer, Array, Array, spmatrix, BCOO, list, None]) – Hamiltonian operators.

  • static_dissipators (Union[ndarray, number, int, float, complex, Tracer, Array, Array, spmatrix, BCOO, list, None]) – Constant dissipation operators.

  • dissipator_operators (Union[ndarray, number, int, float, complex, Tracer, Array, Array, spmatrix, BCOO, list, None]) – Dissipation operators with time-dependent coefficients.

  • hamiltonian_channels (Optional[List[str]]) – List of channel names in pulse schedules corresponding to Hamiltonian operators.

  • dissipator_channels (Optional[List[str]]) – List of channel names in pulse schedules corresponding to dissipator operators.

  • channel_carrier_freqs (Optional[dict]) – Dictionary mapping channel names to floats which represent the carrier frequency of the pulse channel with the corresponding name.

  • dt (Optional[float]) – Sample rate for simulating pulse schedules.

  • rotating_frame (Union[ndarray, number, int, float, complex, Tracer, Array, Array, spmatrix, BCOO, list, RotatingFrame, None]) – Rotating frame to transform the model into. Rotating frames which are diagonal can be supplied as a 1d array of the diagonal elements, to explicitly indicate that they are diagonal.

  • in_frame_basis (bool) – Whether to represent the model in the basis in which the rotating frame operator is diagonalized. See class documentation for a more detailed explanation on how this argument affects object behaviour.

  • array_library (Optional[str]) – Array library to use for storing operators of underlying model. See the model evaluation section of the Models API documentation for a more detailed description of this argument.

  • vectorized (Optional[bool]) – If including dissipator terms, whether or not to construct the LindbladModel in vectorized form. See the model evaluation section of the Models API documentation for a more detailed description of this argument.

  • rwa_cutoff_freq (Optional[float]) – Rotating wave approximation cutoff frequency. If None, no approximation is made.

  • rwa_carrier_freqs (Union[ndarray, number, int, float, complex, Tracer, Array, Array, spmatrix, BCOO, list, Tuple[Union[ndarray, number, int, float, complex, Tracer, Array, Array, spmatrix, BCOO, list], Union[ndarray, number, int, float, complex, Tracer, Array, Array, spmatrix, BCOO, list]], None]) – Carrier frequencies to use for rotating wave approximation. If no time dependent coefficients in model leave as None, if no time-dependent dissipators specify as a list of frequencies for each Hamiltonian operator, and if time-dependent dissipators present specify as a tuple of lists of frequencies, one for Hamiltonian operators and one for dissipators.

  • validate (bool) – Whether or not to validate Hamiltonian operators as being Hermitian.

Raises:
  • QiskitError – If arguments concerning pulse-schedule interpretation are insufficiently

  • specified.

Methods

solve(t_span, y0, signals=None, convert_results=True, **kwargs)[source]#

Solve a dynamical problem, or a set of dynamical problems.

Calls solve_lmde(), and returns an OdeResult object in the style of scipy.integrate.solve_ivp, with results formatted to be the same types as the input. See Additional Information for special handling of various input types, and for specifying multiple simulations at once.

Parameters:
  • t_span (Union[ndarray, number, int, float, complex, Tracer, Array, Array, spmatrix, BCOO, list]) – Time interval to integrate over.

  • y0 (Union[ndarray, number, int, float, complex, Tracer, Array, Array, spmatrix, BCOO, list, QuantumState, BaseOperator]) – Initial state.

  • signals (Union[List[Union[Schedule, ScheduleBlock]], List[Signal], Tuple[List[Signal], List[Signal]], None]) – Specification of time-dependent coefficients to simulate, either in Signal format or as Qiskit Pulse Pulse schedules. If specifying in Signal format, if dissipator_operators is None, specify as a list of signals for the Hamiltonian component, otherwise specify as a tuple of two lists, one for Hamiltonian components, and one for the dissipator_operators coefficients.

  • convert_results (bool) – If True, convert returned solver state results to the same class as y0. If False, states will be returned in the native array type used by the specified solver method.

  • **kwargs – Keyword args passed to solve_lmde().

Returns:

object with formatted output types.

Return type:

OdeResult

Raises:

QiskitError – Initial state y0 is of invalid shape. If signals specifies Schedule simulation but Solver hasn’t been configured to simulate pulse schedules.

Additional Information:

The behaviour of this method is impacted by the input type of y0 and the internal model, summarized in the following table:

Table 2 Type-based behaviour#

y0 type

Model type

yf type

Description

ArrayLike, np.ndarray, Operator

Any

Same as y0

Solves either the Schrodinger equation or Lindblad equation with initial state y0 as specified.

Statevector

HamiltonianModel

Statevector

Solves the Schrodinger equation with initial state y0.

DensityMatrix

HamiltonianModel

DensityMatrix

Solves the Schrodinger equation with initial state the identity matrix to compute the unitary, then conjugates y0 with the result to solve for the density matrix.

Statevector, DensityMatrix

LindbladModel

DensityMatrix

Solve the Lindblad equation with initial state y0, converting to a DensityMatrix first if y0 is a Statevector.

QuantumChannel

HamiltonianModel

SuperOp

Converts y0 to a SuperOp representation, then solves the Schrodinger equation with initial state the identity matrix to compute the unitary and composes with y0.

QuantumChannel

LindbladModel

SuperOp

Solves the vectorized Lindblad equation with initial state y0. vectorized must be set to True.

In some cases (e.g. if using JAX), wrapping the returned states in the type given in the yf type column above may be undesirable. Setting convert_results=False prevents this wrapping, while still allowing usage of the automatic type-handling for the input state.

In addition to the above, this method can be used to specify multiple simulations simultaneously. This can be done by specifying one or more of the arguments t_span, y0, or signals as a list of valid inputs. For this mode of operation, all of these arguments must be either lists of the same length, or a single valid input, which will be used repeatedly.

For example the following code runs three simulations, returning results in a list:

t_span = [span1, span2, span3]
y0 = [state1, state2, state3]
signals = [signals1, signals2, signals3]

results = solver.solve(t_span=t_span, y0=y0, signals=signals)

The following code block runs three simulations, for different sets of signals, repeatedly using the same t_span and y0:

t_span = [t0, tf]
y0 = state1
signals = [signals1, signals2, signal3]

results = solver.solve(t_span=t_span, y0=y0, signals=signals)

Attributes

model#

The model of the system, either a Hamiltonian or Lindblad model.