Dynamics and Time-Evolution

Schrödinger Equation

This module provides solvers for the unitary Schrodinger equation.

sesolve(
H: QobjEvoLike,
psi0: Qobj,
tlist: ArrayLike,
_e_ops=None,
_args=None,
_options=None,
*,
e_ops: EopsLike | list[EopsLike] | dict[Any, EopsLike] = None,
args: dict[str, Any] = None,
options: dict[str, Any] = None,
**kwargs,
) Result[source]

Schrodinger equation evolution of a state vector or unitary matrix for a given Hamiltonian.

Evolve the state vector (psi0) using a given Hamiltonian (H), by integrating the set of ordinary differential equations that define the system. Alternatively evolve a unitary matrix in solving the Schrodinger operator equation.

The output is either the state vector or unitary matrix at arbitrary points in time (tlist), or the expectation values of the supplied operators (e_ops). If e_ops is a callback function, it is invoked for each time in tlist with time and the state as arguments, and the function does not use any return values. e_ops cannot be used in conjunction with solving the Schrodinger operator equation

Time-dependent operators

For time-dependent problems, H and c_ops can be a QobjEvo or object that can be interpreted as QobjEvo such as a list of (Qobj, Coefficient) pairs or a function.

Parameters:
HQobj, QobjEvo, QobjEvo compatible format.

System Hamiltonian as a Qobj or QobjEvo for time-dependent Hamiltonians. List of [Qobj, Coefficient] or callable that can be made into QobjEvo are also accepted.

psi0Qobj

initial state vector (ket) or initial unitary operator psi0 = U

tlistlist / array

list of times for \(t\).

e_opsQobj, callable, list or dict, optional

Single operator, or list or dict of operators, for which to evaluate expectation values. Operator can be Qobj, QobjEvo or callables with the signature f(t: float, state: Qobj) -> Any.

argsdict, optional

dictionary of parameters for time-dependent Hamiltonians

optionsdict, optional

Dictionary of options for the solver.

  • store_final_state : bool
    Whether or not to store the final state of the evolution in the result class.
  • store_states : bool, None
    Whether or not to store the state vectors or density matrices. On None the states will be saved if no expectation operators are given.
  • normalize_output : bool
    Normalize output state to hide ODE numerical errors. Only normalize the state if the initial state is already normalized.
  • progress_bar : str {‘text’, ‘enhanced’, ‘tqdm’, ‘’}
    How to present the solver progress. ‘tqdm’ uses the python module of the same name and raise an error if not installed. Empty string or False will disable the bar.
  • progress_kwargs : dict
    kwargs to pass to the progress_bar. Qutip’s bars use chunk_size.
  • method : str [“adams”, “bdf”, “lsoda”, “dop853”, “vern9”, etc.]
    Which differential equation integration method to use.
  • atol, rtol : float
    Absolute and relative tolerance of the ODE integrator.
  • nsteps : int
    Maximum number of (internally defined) steps allowed in one tlist step.
  • max_step : float
    Maximum lenght of one internal step. When using pulses, it should be less than half the width of the thinnest pulse.

Other options could be supported depending on the integration method, see Integrator.

Returns:
result: Result

An instance of the class Result, which contains a list of array result.expect of expectation values for the times specified by tlist, and/or a list result.states of state vectors or density matrices corresponding to the times in tlist [if e_ops is an empty list of store_states=True in options].

krylovsolve(
H: Qobj,
psi0: Qobj,
tlist: ArrayLike,
krylov_dim: int,
_e_ops=None,
_args=None,
_options=None,
*,
e_ops: dict[Any, Qobj | QobjEvo | Callable[[float, Qobj], Any]] = None,
args: dict[str, Any] = None,
options: dict[str, Any] = None,
) Result[source]

Schrodinger equation evolution of a state vector for time independent Hamiltonians using Krylov method.

Evolve the state vector (“psi0”) finding an approximation for the time evolution operator of Hamiltonian (“H”) by obtaining the projection of the time evolution operator on a set of small dimensional Krylov subspaces (m << dim(H)).

The output is either the state vector or unitary matrix at arbitrary points in time (tlist), or the expectation values of the supplied operators (e_ops). If e_ops is a callback function, it is invoked for each time in tlist with time and the state as arguments, and the function does not use any return values. e_ops cannot be used in conjunction with solving the Schrodinger operator equation

Parameters:
HQobj, QobjEvo, QobjEvo compatible format.

System Hamiltonian as a Qobj or QobjEvo for time-dependent Hamiltonians. List of [Qobj, Coefficient] or callable that can be made into QobjEvo are also accepted.

psi0Qobj

Initial state vector (ket)

tlistlist / array

list of times for \(t\).

krylov_dim: int

Dimension of Krylov approximation subspaces used for the time evolution approximation.

e_opsQobj, callable, or list, optional

Single operator or list of operators for which to evaluate expectation values or callable or list of callable. Callable signature must be, f(t: float, state: Qobj). See expect for more detail of operator expectation.

argsdict, optional

dictionary of parameters for time-dependent Hamiltonians

optionsdict, optional

Dictionary of options for the solver.

  • store_final_state : bool
    Whether or not to store the final state of the evolution in the result class.
  • store_states : bool, None
    Whether or not to store the state vectors or density matrices. On None the states will be saved if no expectation operators are given.
  • normalize_output : bool
    Normalize output state to hide ODE numerical errors. Only normalize the state if the initial state is already normalized.
  • progress_bar : str {‘text’, ‘enhanced’, ‘tqdm’, ‘’}
    How to present the solver progress. ‘tqdm’ uses the python module of the same name and raise an error if not installed. Empty string or False will disable the bar.
  • progress_kwargs : dict
    kwargs to pass to the progress_bar. Qutip’s bars use chunk_size.
  • atol: float
    Absolute tolerance of the ODE integrator.
  • nsteps : int
    Maximum number of (internally defined) steps allowed in one tlist step.
  • min_step, max_step : float
    Miniumum and maximum lenght of one internal step.
  • always_compute_step: bool
    If True, the step lenght is computed each time a new Krylov subspace is computed. Otherwise it is computed only once when creating the integrator.
  • sub_system_tol: float
    Tolerance to detect an happy breakdown. An happy breakdown happens when the initial ket is in a subspace of the Hamiltonian smaller than krylov_dim.
Returns:
result: Result

An instance of the class Result, which contains a list of array result.expect of expectation values for the times specified by tlist, and/or a list result.states of state vectors or density matrices corresponding to the times in tlist [if e_ops is an empty list of store_states=True in options].

class SESolver(
H: Qobj | QobjEvo,
*,
options: dict[str, Any] = None,
)[source]

Bases: Solver

Schrodinger equation evolution of a state vector or unitary matrix for a given Hamiltonian.

Parameters:
HQobj, QobjEvo

System Hamiltonian as a Qobj or QobjEvo for time-dependent Hamiltonians. List of [Qobj, Coefficient] or callable that can be made into QobjEvo are also accepted.

optionsdict, optional

Options for the solver, see SESolver.options and Integrator for a list of all options.

Attributes:
stats: dict

Diverse diagnostic statistics of the evolution.

classmethod ExpectFeedback(
operator: Qobj | QobjEvo,
default: Any = 0.0,
)

Expectation value of the instantaneous state of the evolution to be used by a time-dependent operator.

When used as an args:

QobjEvo([op, func], args={"E0": Solver.ExpectFeedback(oper)})

The func will receive expect(oper, state) as E0 during the evolution.

Parameters:
operatorQobj, QobjEvo

Operator to compute the expectation values of.

defaultfloat, default0.

Initial value to be used at setup.

classmethod StateFeedback(
default: Qobj | Data = None,
raw_data: bool = False,
prop: bool = False,
)[source]

State of the evolution to be used in a time-dependent operator.

When used as an args:

QobjEvo([op, func], args={"state": SESolver.StateFeedback()})

The func will receive the ket as state during the evolution.

Parameters:
defaultQobj or qutip.core.data.Data, defaultNone

Initial value to be used at setup of the system.

propbool, defaultFalse

Set to True when using sesolve for computing propagators.

raw_databool, defaultFalse

If True, the raw matrix will be passed instead of a Qobj. For density matrices, the matrices can be column stacked or square depending on the integration method.

property options: dict

Solver’s options:

store_final_state: bool, default: False

Whether or not to store the final state of the evolution in the result class.

store_states: bool, default: None

Whether or not to store the state vectors or density matrices. On None the states will be saved if no expectation operators are given.

normalize_output: bool, default: True

Normalize output state to hide ODE numerical errors.

progress_bar: str {“text”, “enhanced”, “tqdm”, “”}, default: “”

How to present the solver progress. ‘tqdm’ uses the python module of the same name and raise an error if not installed. Empty string or False will disable the bar.

progress_kwargs: dict, default: {“chunk_size”: 10}

Arguments to pass to the progress_bar. Qutip’s bars use chunk_size.

method: str, default: “adams”

Which ordinary differential equation integration method to use.

run(
state0: Qobj,
tlist: ArrayLike,
*,
e_ops: EopsLike | list[EopsLike] | dict[Any, EopsLike] = None,
args: dict[str, Any] = None,
) Result

Do the evolution of the Quantum system.

For a state0 at time tlist[0] do the evolution as directed by rhs and for each time in tlist store the state and/or expectation values in a Result. The evolution method and stored results are determined by options.

Parameters:
state0Qobj

Initial state of the evolution.

tlistlist of double

Time for which to save the results (state and/or expect) of the evolution. The first element of the list is the initial time of the evolution. Each times of the list must be increasing, but does not need to be uniformy distributed.

argsdict, optional

Change the args of the rhs for the evolution.

e_opsQobj, QobjEvo, callable, list, or dict optional

Single, list or dict of Qobj, QobjEvo or callable to compute the expectation values. Function[s] must have the signature f(t : float, state : Qobj) -> expect.

Returns:
resultsResult

Results of the evolution. States and/or expect will be saved. You can control the saved data in the options.

start(state0: Qobj, t0: Number) None

Set the initial state and time for a step evolution.

Parameters:
state0Qobj

Initial state of the evolution.

t0double

Initial time of the evolution.

step(
t: Number,
*,
args: dict[str, Any] = None,
copy: bool = True,
) Qobj

Evolve the state to t and return the state as a Qobj.

Parameters:
tdouble

Time to evolve to, must be higher than the last call.

argsdict, optional {None}

Update the args of the system. The change is effective from the beginning of the interval. Changing args can slow the evolution.

copybool, optional {True}

Whether to return a copy of the data or the data in the ODE solver.

Notes

The state must be initialized first by calling start or run. If run is called, step will continue from the last time and state obtained.

property sys_dims

Dimensions of the space that the system use:

qutip.basis(sovler.dims) will create a state with proper dimensions for this solver.

Master Equation

This module provides solvers for the Lindblad master equation and von Neumann equation.

mesolve(
H: QobjEvoLike,
rho0: Qobj,
tlist: ArrayLike,
c_ops: Qobj | QobjEvo | list[QobjEvoLike] = None,
_e_ops=None,
_args=None,
_options=None,
*,
e_ops: EopsLike | list[EopsLike] | dict[Any, EopsLike] = None,
args: dict[str, Any] = None,
options: dict[str, Any] = None,
**kwargs,
) Result[source]

Master equation evolution of a density matrix for a given Hamiltonian and set of collapse operators, or a Liouvillian.

Evolve the state vector or density matrix (rho0) using a given Hamiltonian or Liouvillian (H) and an optional set of collapse operators (c_ops), by integrating the set of ordinary differential equations that define the system. In the absence of collapse operators the system is evolved according to the unitary evolution of the Hamiltonian.

The output is either the state vector at arbitrary points in time (tlist), or the expectation values of the supplied operators (e_ops). If e_ops is a callback function, it is invoked for each time in tlist with time and the state as arguments, and the function does not use any return values.

If either H or the Qobj elements in c_ops are superoperators, they will be treated as direct contributions to the total system Liouvillian. This allows the solution of master equations that are not in standard Lindblad form.

Time-dependent operators

For time-dependent problems, H and c_ops can be a QobjEvo or object that can be interpreted as QobjEvo such as a list of (Qobj, Coefficient) pairs or a function.

Additional options

Additional options to mesolve can be set via the options argument. Many ODE integration options can be set this way, and the store_states and store_final_state options can be used to store states even though expectation values are requested via the e_ops argument.

Parameters:
HQobj, QobjEvo, QobjEvo compatible format.

Possibly time-dependent system Liouvillian or Hamiltonian as a Qobj or QobjEvo. List of [Qobj, Coefficient] or callable that can be made into QobjEvo are also accepted.

rho0Qobj

initial density matrix or state vector (ket).

tlistlist / array

list of times for \(t\).

c_opslist of (QobjEvo, QobjEvo compatible format)

Single collapse operator, or list of collapse operators, or a list of Liouvillian superoperators. None is equivalent to an empty list.

e_opsQobj, callable, list or dict, optional

Single operator, or list or dict of operators, for which to evaluate expectation values. Operator can be Qobj, QobjEvo or callables with the signature f(t: float, state: Qobj) -> Any.

argsdict, optional

dictionary of parameters for time-dependent Hamiltonians and collapse operators.

optionsdict, optional

Dictionary of options for the solver.

  • store_final_state : bool
    Whether or not to store the final state of the evolution in the result class.
  • store_states : bool, None
    Whether or not to store the state vectors or density matrices. On None the states will be saved if no expectation operators are given.
  • normalize_output : bool
    Normalize output state to hide ODE numerical errors. Only normalize the state if the initial state is already normalized.
  • progress_bar : str {‘text’, ‘enhanced’, ‘tqdm’, ‘’}
    How to present the solver progress. ‘tqdm’ uses the python module of the same name and raise an error if not installed. Empty string or False will disable the bar.
  • progress_kwargs : dict
    kwargs to pass to the progress_bar. Qutip’s bars use chunk_size.
  • method : str [“adams”, “bdf”, “lsoda”, “dop853”, “vern9”, etc.]
    Which differential equation integration method to use.
  • atol, rtol : float
    Absolute and relative tolerance of the ODE integrator.
  • nsteps : int
    Maximum number of (internally defined) steps allowed in one tlist step.
  • max_step : float
    Maximum lenght of one internal step. When using pulses, it should be less than half the width of the thinnest pulse.

Other options could be supported depending on the integration method, see Integrator.

Returns:
result: Result

An instance of the class Result, which contains a list of array result.expect of expectation values for the times specified by tlist, and/or a list result.states of state vectors or density matrices corresponding to the times in tlist [if e_ops is an empty list of store_states=True in options].

Notes

When no collapse operator are given and the H is not a superoperator, it will defer to sesolve.

class MESolver(
H: Qobj | QobjEvo,
c_ops: Qobj | QobjEvo | list[Qobj | QobjEvo] = None,
*,
options: dict = None,
)[source]

Bases: SESolver

Master equation evolution of a density matrix for a given Hamiltonian and set of collapse operators, or a Liouvillian.

Evolve the density matrix (rho0) using a given Hamiltonian or Liouvillian (H) and an optional set of collapse operators (c_ops), by integrating the set of ordinary differential equations that define the system.

If either H or the Qobj elements in c_ops are superoperators, they will be treated as direct contributions to the total system Liouvillian. This allows the solution of master equations that are not in standard Lindblad form.

Parameters:
HQobj, QobjEvo

Possibly time-dependent system Liouvillian or Hamiltonian as a Qobj or QobjEvo. List of [Qobj, Coefficient] or callable that can be made into QobjEvo are also accepted.

c_opslist of Qobj, QobjEvo

Single collapse operator, or list of collapse operators, or a list of Liouvillian superoperators. None is equivalent to an empty list.

optionsdict, optional

Options for the solver, see MESolver.options and Integrator for a list of all options.

Attributes:
stats: dict

Diverse diagnostic statistics of the evolution.

classmethod ExpectFeedback(
operator: Qobj | QobjEvo,
default: Any = 0.0,
)

Expectation value of the instantaneous state of the evolution to be used by a time-dependent operator.

When used as an args:

QobjEvo([op, func], args={"E0": Solver.ExpectFeedback(oper)})

The func will receive expect(oper, state) as E0 during the evolution.

Parameters:
operatorQobj, QobjEvo

Operator to compute the expectation values of.

defaultfloat, default0.

Initial value to be used at setup.

classmethod StateFeedback(
default: Qobj | Data = None,
raw_data: bool = False,
prop: bool = False,
)[source]

State of the evolution to be used in a time-dependent operator.

When used as an args:

QobjEvo([op, func], args={"state": MESolver.StateFeedback()})

The func will receive the density matrix as state during the evolution.

Parameters:
defaultQobj or qutip.core.data.Data, defaultNone

Initial value to be used at setup of the system.

propbool, defaultFalse

Set to True when computing propagators. The default with take the shape of the propagator instead of a state.

raw_databool, defaultFalse

If True, the raw matrix will be passed instead of a Qobj. For density matrices, the matrices can be column stacked or square depending on the integration method.

property options: dict

Solver’s options:

store_final_state: bool, default: False

Whether or not to store the final state of the evolution in the result class.

store_states: bool, default: None

Whether or not to store the state vectors or density matrices. On None the states will be saved if no expectation operators are given.

normalize_output: bool, default: True

Normalize output state to hide ODE numerical errors.

progress_bar: str {“text”, “enhanced”, “tqdm”, “”}, default: “”

How to present the solver progress. ‘tqdm’ uses the python module of the same name and raise an error if not installed. Empty string or False will disable the bar.

progress_kwargs: dict, default: {“chunk_size”: 10}

Arguments to pass to the progress_bar. Qutip’s bars use chunk_size.

method: str, default: “adams”

Which ordinary differential equation integration method to use.

run(
state0: Qobj,
tlist: ArrayLike,
*,
e_ops: EopsLike | list[EopsLike] | dict[Any, EopsLike] = None,
args: dict[str, Any] = None,
) Result

Do the evolution of the Quantum system.

For a state0 at time tlist[0] do the evolution as directed by rhs and for each time in tlist store the state and/or expectation values in a Result. The evolution method and stored results are determined by options.

Parameters:
state0Qobj

Initial state of the evolution.

tlistlist of double

Time for which to save the results (state and/or expect) of the evolution. The first element of the list is the initial time of the evolution. Each times of the list must be increasing, but does not need to be uniformy distributed.

argsdict, optional

Change the args of the rhs for the evolution.

e_opsQobj, QobjEvo, callable, list, or dict optional

Single, list or dict of Qobj, QobjEvo or callable to compute the expectation values. Function[s] must have the signature f(t : float, state : Qobj) -> expect.

Returns:
resultsResult

Results of the evolution. States and/or expect will be saved. You can control the saved data in the options.

start(state0: Qobj, t0: Number) None

Set the initial state and time for a step evolution.

Parameters:
state0Qobj

Initial state of the evolution.

t0double

Initial time of the evolution.

step(
t: Number,
*,
args: dict[str, Any] = None,
copy: bool = True,
) Qobj

Evolve the state to t and return the state as a Qobj.

Parameters:
tdouble

Time to evolve to, must be higher than the last call.

argsdict, optional {None}

Update the args of the system. The change is effective from the beginning of the interval. Changing args can slow the evolution.

copybool, optional {True}

Whether to return a copy of the data or the data in the ODE solver.

Notes

The state must be initialized first by calling start or run. If run is called, step will continue from the last time and state obtained.

property sys_dims

Dimensions of the space that the system use:

qutip.basis(sovler.dims) will create a state with proper dimensions for this solver.

class Result(
e_ops: dict[Any, Qobj | QobjEvo | Callable[[float, Qobj], Any]],
options: ResultOptions,
*,
solver: str = None,
stats: dict[str, Any] = None,
**kw,
)[source]

Base class for storing solver results.

Parameters:
e_opsQobj, QobjEvo, function or list or dict of these

The e_ops parameter defines the set of values to record at each time step t. If an element is a Qobj or QobjEvo the value recorded is the expectation value of that operator given the state at t. If the element is a function, f, the value recorded is f(t, state).

The values are recorded in the e_data and expect attributes of this result object. e_data is a dictionary and expect is a list, where each item contains the values of the corresponding e_op.

optionsdict

The options for this result class.

solverstr or None

The name of the solver generating these results.

statsdict or None

The stats generated by the solver while producing these results. Note that the solver may update the stats directly while producing results.

kwdict

Additional parameters specific to a result sub-class.

Attributes:
timeslist

A list of the times at which the expectation values and states were recorded.

stateslist of Qobj

The state at each time t (if the recording of the state was requested).

final_stateQobj:

The final state (if the recording of the final state was requested).

expectlist of arrays of expectation values

A list containing the values of each e_op. The list is in the same order in which the e_ops were supplied and empty if no e_ops were given.

Each element is itself a list and contains the values of the corresponding e_op, with one value for each time in .times.

The same lists of values may be accessed via the .e_data dictionary and the original e_ops are available via the .e_ops attribute.

e_datadict

A dictionary containing the values of each e_op. If the e_ops were supplied as a dictionary, the keys are the same as in that dictionary. Otherwise the keys are the index of the e_op in the .expect list.

The lists of expectation values returned are the same lists as those returned by .expect.

e_opsdict

A dictionary containing the supplied e_ops as ExpectOp instances. The keys of the dictionary are the same as for .e_data. Each value is object where .e_ops[k](t, state) calculates the value of e_op k at time t and the given state, and .e_ops[k].op is the original object supplied to create the e_op.

solverstr or None

The name of the solver generating these results.

statsdict or None

The stats generated by the solver while producing these results.

optionsdict

The options for this result class.

Monte Carlo Evolution

mcsolve(
H: QobjEvoLike,
state: Qobj,
tlist: ArrayLike,
c_ops: QobjEvoLike | list[QobjEvoLike] = (),
_e_ops=None,
_ntraj=None,
*,
e_ops: EopsLike | list[EopsLike] | dict[Any, EopsLike] = None,
ntraj: int = 500,
args: dict[str, Any] = None,
options: dict[str, Any] = None,
seeds: int | SeedSequence | list[int | SeedSequence] = None,
target_tol: float | tuple[float, float] | list[tuple[float, float]] = None,
timeout: float = None,
**kwargs,
) McResult[source]

Monte Carlo evolution of a state vector \(|\psi \rangle\) for a given Hamiltonian and sets of collapse operators. Options for the underlying ODE solver are given by the Options class.

Parameters:
HQobj, QobjEvo, list, callable.

System Hamiltonian as a Qobj, QobjEvo. It can also be any input type that QobjEvo accepts (see QobjEvo’s documentation). H can also be a superoperator (liouvillian) if some collapse operators are to be treated deterministically.

stateQobj

Initial state vector or density matrix.

tlistarray_like

Times at which results are recorded.

c_opslist

A list of collapse operators in any input type that QobjEvo accepts (see QobjEvo’s documentation). They must be operators even if H is a superoperator. If none are given, the solver will defer to sesolve or mesolve.

e_opsQobj, callable, list or dict, optional

Single operator, or list or dict of operators, for which to evaluate expectation values. Operator can be Qobj, QobjEvo or callables with the signature f(t: float, state: Qobj) -> Any.

ntrajint, default: 500

Maximum number of trajectories to run. Can be cut short if a time limit is passed with the timeout keyword or if the target tolerance is reached, see target_tol.

argsdict, optional

Arguments for time-dependent Hamiltonian and collapse operator terms.

optionsdict, optional

Dictionary of options for the solver.

  • store_final_state : bool
    Whether or not to store the final state of the evolution in the result class.
  • store_states : bool, None
    Whether or not to store the state vectors or density matrices. On None the states will be saved if no expectation operators are given.
  • progress_bar : str {‘text’, ‘enhanced’, ‘tqdm’, ‘’}
    How to present the solver progress. ‘tqdm’ uses the python module of the same name and raise an error if not installed. Empty string or False will disable the bar.
  • progress_kwargs : dict
    kwargs to pass to the progress_bar. Qutip’s bars use chunk_size.
  • method : str [“adams”, “bdf”, “lsoda”, “dop853”, “vern9”, etc.]
    Which differential equation integration method to use.
  • atol, rtol : float
    Absolute and relative tolerance of the ODE integrator.
  • nsteps : int
    Maximum number of (internally defined) steps allowed in one tlist step.
  • max_step : float
    Maximum length of one internal step. When using pulses, it should be less than half the width of the thinnest pulse.
  • keep_runs_results : bool, [False]
    Whether to store results from all trajectories or just store the averages.
  • map : str {“serial”, “parallel”, “loky”, “mpi”}
    How to run the trajectories. “parallel” uses the multiprocessing module to run in parallel while “loky” and “mpi” use the “loky” and “mpi4py” modules to do so.
  • num_cpus : int
    Number of cpus to use when running in parallel. None detect the number of available cpus.
  • norm_t_tol, norm_tol, norm_steps : float, float, int
    Parameters used to find the collapse location. norm_t_tol and norm_tol are the tolerance in time and norm respectively. An error will be raised if the collapse could not be found within norm_steps tries.
  • mc_corr_eps : float
    Small number used to detect non-physical collapse caused by numerical imprecision.
  • improved_sampling : Bool
    Whether to use the improved sampling algorithm from Abdelhafez et al. PRA (2019)

Additional options are listed under options. More options may be available depending on the selected differential equation integration method, see Integrator.

seedsint, SeedSequence, list, optional

Seed for the random number generator. It can be a single seed used to spawn seeds for each trajectory or a list of seeds, one for each trajectory. Seeds are saved in the result and they can be reused with:

seeds=prev_result.seeds
target_tolfloat, tuple, list, optional

Target tolerance of the evolution. The evolution will compute trajectories until the error on the expectation values is lower than this tolerance. The maximum number of trajectories employed is given by ntraj. The error is computed using jackknife resampling. target_tol can be an absolute tolerance or a pair of absolute and relative tolerance, in that order. Lastly, it can be a list of pairs of (atol, rtol) for each e_ops.

timeoutfloat, optional

Maximum time for the evolution in second. When reached, no more trajectories will be computed.

Returns:
resultsMcResult

Object storing all results from the simulation. Which results is saved depends on the presence of e_ops and the options used. collapse and photocurrent is available to Monte Carlo simulation results. If the initial condition is mixed, the result has additional attributes initial_states and ntraj_per_initial_state.

Notes

The simulation will end when the first end condition is reached between ntraj, timeout and target_tol. If the initial condition is mixed, target_tol is not supported. If the initial condition is mixed, and the end condition is not ntraj, the results returned by this function should be considered invalid.

class MCSolver(
H: Qobj | QobjEvo,
c_ops: Qobj | QobjEvo | list[Qobj | QobjEvo],
*,
options: dict[str, Any] = None,
)[source]

Bases: MultiTrajSolver

Monte Carlo Solver of a state vector \(|\psi \rangle\) for a given Hamiltonian and sets of collapse operators. Options for the underlying ODE solver are given by the Options class.

Parameters:
HQobj, QobjEvo, list, callable.

System Hamiltonian as a Qobj, QobjEvo. It can also be any input type that QobjEvo accepts (see QobjEvo’s documentation). H can also be a superoperator (liouvillian) if some collapse operators are to be treated deterministically.

c_opslist

A list of collapse operators in any input type that QobjEvo accepts (see QobjEvo’s documentation). They must be operators even if H is a superoperator.

optionsdict, [optional]

Options for the evolution.

classmethod CollapseFeedback(default: list = None)[source]

Collapse of the trajectory argument for time dependent systems.

When used as an args:

QobjEvo([op, func], args={"cols": MCSolver.CollapseFeedback()})

The func will receive a list of (time, operator number) for each collapses of the trajectory as cols.

Note

CollapseFeedback can’t be added to a running solver when updating arguments between steps: solver.step(..., args={}).

Parameters:
defaultlist, default[]

Argument value to use outside of solver.

classmethod ExpectFeedback(
operator: Qobj | QobjEvo,
default: Any = 0.0,
)

Expectation value of the instantaneous state of the evolution to be used by a time-dependent operator.

When used as an args:

QobjEvo([op, func], args={"E0": Solver.ExpectFeedback(oper)})

The func will receive expect(oper, state) as E0 during the evolution.

Parameters:
operatorQobj, QobjEvo

Operator to compute the expectation values of.

defaultfloat, default0.

Initial value to be used at setup.

classmethod StateFeedback(
default: Qobj | Data = None,
raw_data: bool = False,
prop: bool = False,
)[source]

State of the evolution to be used in a time-dependent operator.

When used as an args:

QobjEvo([op, func], args={"state": MCSolver.StateFeedback()})

The func will receive the density matrix as state during the evolution.

Parameters:
defaultQobj or qutip.core.data.Data, defaultNone

Initial value to be used at setup of the system.

openbool, default False

Set to True when using the monte carlo solver for open systems.

raw_databool, defaultFalse

If True, the raw matrix will be passed instead of a Qobj. For density matrices, the matrices can be column stacked or square depending on the integration method.

property options: dict

Options for monte carlo solver:

store_final_state: bool, default: False

Whether or not to store the final state of the evolution in the result class.

store_states: bool, default: None

Whether or not to store the state vectors or density matrices. On None the states will be saved if no expectation operators are given.

progress_bar: str {‘text’, ‘enhanced’, ‘tqdm’, ‘’}, default: “text”

How to present the solver progress. ‘tqdm’ uses the python module of the same name and raise an error if not installed. Empty string or False will disable the bar.

progress_kwargs: dict, default: {“chunk_size”:10}

Arguments to pass to the progress_bar. Qutip’s bars use chunk_size.

keep_runs_results: bool, default: False

Whether to store results from all trajectories or just store the averages.

method: str, default: “adams”

Which differential equation integration method to use.

map: str {“serial”, “parallel”, “loky”, “mpi”}, default: “serial”

How to run the trajectories. “parallel” uses the multiprocessing module to run in parallel while “loky” and “mpi” use the “loky” and “mpi4py” modules to do so.

mpi_options: dict, default: {}

Only applies if map is “mpi”. This dictionary will be passed as keyword arguments to the mpi4py.futures.MPIPoolExecutor constructor. Note that the max_workers argument is provided separately through the num_cpus option.

num_cpus: None, int

Number of cpus to use when running in parallel. None detect the number of available cpus.

bitgenerator: {None, “MT19937”, “PCG64”, “PCG64DXSM”, …}

Which of numpy.random’s bitgenerator to use. With None, your numpy version’s default is used.

mc_corr_eps: float, default: 1e-10

Small number used to detect non-physical collapse caused by numerical imprecision.

norm_t_tol: float, default: 1e-6

Tolerance in time used when finding the collapse.

norm_tol: float, default: 1e-4

Tolerance in norm used when finding the collapse.

norm_steps: int, default: 5

Maximum number of tries to find the collapse.

improved_sampling: Bool, default: False

Whether to use the improved sampling algorithm of Abdelhafez et al. PRA (2019)

run(
state: Qobj | list[tuple[Qobj, float]],
tlist: ArrayLike,
ntraj: int | list[int] = None,
*,
args: dict[str, Any] = None,
e_ops: EopsLike | list[EopsLike] | dict[Any, EopsLike] = None,
target_tol: float | tuple[float, float] | list[tuple[float, float]] = None,
timeout: float = None,
seeds: int | SeedSequence | list[int | SeedSequence] = None,
) McResult[source]

Do the evolution of the Quantum system.

For a state at time tlist[0], do up to ntraj simulations of the Monte-Carlo evolution. For each time in tlist store the state and/or expectation values in a MultiTrajResult. The evolution method and stored results are determined by options.

Parameters:
state{Qobj, list of (Qobj, float)}

Initial state of the evolution. May be either a pure state or a statistical ensemble. An ensemble can be provided either as a density matrix, or as a list of tuples. In the latter case, the first element of each tuple is a pure state, and the second element is its weight, i.e., a number between 0 and 1 describing the fraction of the ensemble in that state. The sum of all weights must be one.

tlistlist of double

Time for which to save the results (state and/or expect) of the evolution. The first element of the list is the initial time of the evolution. Time in the list must be in increasing order, but does not need to be uniformly distributed.

ntraj{int, list of int}

Number of trajectories to add. If the initial state is pure, this must be single number. If the initial state is a mixed ensemble, specified as a list of pure states, this parameter may also be a list of numbers with the same number of entries. It then specifies the number of trajectories for each pure state. If the initial state is mixed and this parameter is a single number, it specifies the total number of trajectories, which are distributed over the initial ensemble automatically.

argsdict, optional

Change the args of the rhs for the evolution.

e_opsQobj, callable, list or dict, optional

Single operator, or list or dict of operators, for which to evaluate expectation values. Operator can be Qobj, QobjEvo or callables with the signature f(t: float, state: Qobj) -> Any.

timeoutfloat, optional

Maximum time in seconds for the trajectories to run. Once this time is reached, the simulation will end even if the number of trajectories is less than ntraj. The map function, set in options, can interupt the running trajectory or wait for it to finish. Set to an arbitrary high number to disable.

target_tol{float, tuple, list}, optional

Target tolerance of the evolution. The evolution will compute trajectories until the error on the expectation values is lower than this tolerance. The maximum number of trajectories employed is given by ntraj. The error is computed using jackknife resampling. target_tol can be an absolute tolerance or a pair of absolute and relative tolerance, in that order. Lastly, it can be a list of pairs of (atol, rtol) for each e_ops.

seeds{int, SeedSequence, list}, optional

Seed or list of seeds for each trajectories.

Returns:
resultsMcResult

Results of the evolution. States and/or expect will be saved. You can control the saved data in the options. If the initial condition is mixed, the result has additional attributes initial_states and ntraj_per_initial_state.

start(
state0: Qobj,
t0: float,
seed: int | SeedSequence = None,
)

Set the initial state and time for a step evolution.

Parameters:
stateQobj

Initial state of the evolution.

t0double

Initial time of the evolution.

seedint, SeedSequence, list, optional

Seed for the random number generator. It can be a single seed used to spawn seeds for each trajectory or a list of seed, one for each trajectory.

Notes

When using step evolution, only one trajectory can be computed at once.

step(
t: float,
*,
args: dict[str, Any] = None,
copy: bool = True,
) Qobj

Evolve the state to t and return the state as a Qobj.

Parameters:
tdouble

Time to evolve to, must be higher than the last call.

argsdict, optional

Update the args of the system. The change is effective from the beginning of the interval. Changing args can slow the evolution.

copybool, default: True

Whether to return a copy of the data or the data in the ODE solver.

property sys_dims

Dimensions of the space that the system use:

qutip.basis(sovler.dims) will create a state with proper dimensions for this solver.

nm_mcsolve(
H: QobjEvoLike,
state: Qobj,
tlist: ArrayLike,
ops_and_rates: list[tuple[Qobj, CoefficientLike]] = (),
_e_ops=None,
_ntraj=None,
*,
e_ops: EopsLike | list[EopsLike] | dict[Any, EopsLike] = None,
ntraj: int = 500,
args: dict[str, Any] = None,
options: dict[str, Any] = None,
seeds: int | SeedSequence | list[int | SeedSequence] = None,
target_tol: float | tuple[float, float] | list[tuple[float, float]] = None,
timeout: float = None,
**kwargs,
) NmmcResult[source]

Monte-Carlo evolution corresponding to a Lindblad equation with “rates” that may be negative. Usage of this function is analogous to mcsolve, but the c_ops parameter is replaced by an ops_and_rates parameter to allow for negative rates. Options for the underlying ODE solver are given by the Options class.

Parameters:
HQobj, QobjEvo, list, callable.

System Hamiltonian as a Qobj, QobjEvo. It can also be any input type that QobjEvo accepts (see QobjEvo’s documentation). H can also be a superoperator (liouvillian) if some collapse operators are to be treated deterministically.

stateQobj

Initial state vector or density matrix.

tlistarray_like

Times at which results are recorded.

ops_and_rateslist

A list of tuples (L, Gamma), where the Lindblad operator L is a Qobj and Gamma represents the corresponding rate, which is allowed to be negative. The Lindblad operators must be operators even if H is a superoperator. If none are given, the solver will defer to sesolve or mesolve. Each rate Gamma may be just a number (in the case of a constant rate) or, otherwise, specified using any format accepted by coefficient.

e_opsQobj, callable, list or dict, optional

Single operator, or list or dict of operators, for which to evaluate expectation values. Operator can be Qobj, QobjEvo or callables with the signature f(t: float, state: Qobj) -> Any.

ntrajint, default: 500

Maximum number of trajectories to run. Can be cut short if a time limit is passed with the timeout keyword or if the target tolerance is reached, see target_tol.

argsdict, optional

Arguments for time-dependent Hamiltonian and collapse operator terms.

optionsdict, optional

Dictionary of options for the solver.

  • store_final_state : bool
    Whether or not to store the final state of the evolution in the result class.
  • store_states : bool, None
    Whether or not to store the state vectors or density matrices. On None the states will be saved if no expectation operators are given.
  • progress_bar : str {‘text’, ‘enhanced’, ‘tqdm’, ‘’}
    How to present the solver progress. ‘tqdm’ uses the python module of the same name and raise an error if not installed. Empty string or False will disable the bar.
  • progress_kwargs : dict
    kwargs to pass to the progress_bar. Qutip’s bars use chunk_size.
  • method : str [“adams”, “bdf”, “lsoda”, “dop853”, “vern9”, etc.]
    Which differential equation integration method to use.
  • atol, rtol : float
    Absolute and relative tolerance of the ODE integrator.
  • nsteps : int
    Maximum number of (internally defined) steps allowed in one tlist step.
  • max_step : float
    Maximum length of one internal step. When using pulses, it should be less than half the width of the thinnest pulse.
  • keep_runs_results : bool, [False]
    Whether to store results from all trajectories or just store the averages.
  • map : str {“serial”, “parallel”, “loky”, “mpi”}
    How to run the trajectories. “parallel” uses the multiprocessing module to run in parallel while “loky” and “mpi” use the “loky” and “mpi4py” modules to do so.
  • num_cpus : int
    Number of cpus to use when running in parallel. None detect the number of available cpus.
  • norm_t_tol, norm_tol, norm_steps : float, float, int
    Parameters used to find the collapse location. norm_t_tol and norm_tol are the tolerance in time and norm respectively. An error will be raised if the collapse could not be found within norm_steps tries.
  • improved_sampling : Bool
    Whether to use the improved sampling algorithm from Abdelhafez et al. PRA (2019)
  • mc_corr_eps : float
    Small number used to detect non-physical collapse caused by numerical imprecision.
  • completeness_rtol, completeness_atol : float, float
    Parameters used in determining whether the given Lindblad operators satisfy a certain completeness relation. If they do not, an additional Lindblad operator is added automatically (with zero rate).
  • martingale_quad_limit : float or int
    An upper bound on the number of subintervals used in the adaptive integration of the martingale.

Note that the ‘improved_sampling’ option is not currently supported. Additional options are listed under options. More options may be available depending on the selected differential equation integration method, see Integrator.

seedsint, SeedSequence, list, optional

Seed for the random number generator. It can be a single seed used to spawn seeds for each trajectory or a list of seeds, one for each trajectory. Seeds are saved in the result and they can be reused with:

seeds=prev_result.seeds
target_tolfloat, tuple, list, optional

Target tolerance of the evolution. The evolution will compute trajectories until the error on the expectation values is lower than this tolerance. The maximum number of trajectories employed is given by ntraj. The error is computed using jackknife resampling. target_tol can be an absolute tolerance or a pair of absolute and relative tolerance, in that order. Lastly, it can be a list of pairs of (atol, rtol) for each e_ops.

timeoutfloat, optional

Maximum time for the evolution in seconds. When reached, no more trajectories will be computed.

Returns:
resultsNmmcResult

Object storing all results from the simulation. Compared to a result returned by mcsolve, this result contains the additional field trace (and runs_trace if store_final_state is set). Note that the states on the individual trajectories are not normalized. This field contains the average of their trace, which will converge to one in the limit of sufficiently many trajectories. If the initial condition is mixed, the result has additional attributes initial_states and ntraj_per_initial_state.

class NonMarkovianMCSolver(
H: Qobj | QobjEvo,
ops_and_rates: Sequence[tuple[Qobj, float | Coefficient]],
*,
options: dict[str, Any] = None,
)[source]

Bases: MCSolver

Monte Carlo Solver for Lindblad equations with “rates” that may be negative. The c_ops parameter of MCSolver is replaced by an ops_and_rates parameter to allow for negative rates. Options for the underlying ODE solver are given by the Options class.

Parameters:
HQobj, QobjEvo, list, callable.

System Hamiltonian as a Qobj, QobjEvo. It can also be any input type that QobjEvo accepts (see QobjEvo documentation). H can also be a superoperator (liouvillian) if some collapse operators are to be treated deterministically.

ops_and_rateslist

A list of tuples (L, Gamma), where the Lindblad operator L is a Qobj and Gamma represents the corresponding rate, which is allowed to be negative. The Lindblad operators must be operators even if H is a superoperator. Each rate Gamma may be just a number (in the case of a constant rate) or, otherwise, a Coefficient.

optionsSolverOptions, [optional]

Options for the evolution.

classmethod CollapseFeedback(default: list = None)

Collapse of the trajectory argument for time dependent systems.

When used as an args:

QobjEvo([op, func], args={"cols": MCSolver.CollapseFeedback()})

The func will receive a list of (time, operator number) for each collapses of the trajectory as cols.

Note

CollapseFeedback can’t be added to a running solver when updating arguments between steps: solver.step(..., args={}).

Parameters:
defaultlist, default[]

Argument value to use outside of solver.

classmethod ExpectFeedback(
operator: Qobj | QobjEvo,
default: Any = 0.0,
)

Expectation value of the instantaneous state of the evolution to be used by a time-dependent operator.

When used as an args:

QobjEvo([op, func], args={"E0": Solver.ExpectFeedback(oper)})

The func will receive expect(oper, state) as E0 during the evolution.

Parameters:
operatorQobj, QobjEvo

Operator to compute the expectation values of.

defaultfloat, default0.

Initial value to be used at setup.

classmethod StateFeedback(
default: Qobj | Data = None,
raw_data: bool = False,
prop: bool = False,
)

State of the evolution to be used in a time-dependent operator.

When used as an args:

QobjEvo([op, func], args={"state": MCSolver.StateFeedback()})

The func will receive the density matrix as state during the evolution.

Parameters:
defaultQobj or qutip.core.data.Data, defaultNone

Initial value to be used at setup of the system.

openbool, default False

Set to True when using the monte carlo solver for open systems.

raw_databool, defaultFalse

If True, the raw matrix will be passed instead of a Qobj. For density matrices, the matrices can be column stacked or square depending on the integration method.

current_martingale()[source]

Returns the value of the influence martingale along the current trajectory. The value of the martingale is the product of the continuous and the discrete contribution. The current time and the collapses that have happened are read out from the internal integrator.

property options: dict[str, Any]

Options for non-Markovian Monte Carlo solver:

store_final_state: bool, default: False

Whether or not to store the final state of the evolution in the result class.

store_states: bool, default: None

Whether or not to store the state vectors or density matrices. On None the states will be saved if no expectation operators are given.

progress_bar: str {‘text’, ‘enhanced’, ‘tqdm’, ‘’}, default: “text”

How to present the solver progress. ‘tqdm’ uses the python module of the same name and raise an error if not installed. Empty string or False will disable the bar.

progress_kwargs: dict, default: {“chunk_size”:10}

Arguments to pass to the progress_bar. Qutip’s bars use chunk_size.

keep_runs_results: bool, default: False

Whether to store results from all trajectories or just store the averages.

method: str, default: “adams”

Which differential equation integration method to use.

map: str {“serial”, “parallel”, “loky”, “mpi”}, default: “serial”

How to run the trajectories. “parallel” uses the multiprocessing module to run in parallel while “loky” and “mpi” use the “loky” and “mpi4py” modules to do so.

mpi_options: dict, default: {}

Only applies if map is “mpi”. This dictionary will be passed as keyword arguments to the mpi4py.futures.MPIPoolExecutor constructor. Note that the max_workers argument is provided separately through the num_cpus option.

num_cpus: None, int

Number of cpus to use when running in parallel. None detect the number of available cpus.

bitgenerator: {None, “MT19937”, “PCG64”, “PCG64DXSM”, …}

Which of numpy.random’s bitgenerator to use. With None, your numpy version’s default is used.

mc_corr_eps: float, default: 1e-10

Small number used to detect non-physical collapse caused by numerical imprecision.

norm_t_tol: float, default: 1e-6

Tolerance in time used when finding the collapse.

norm_tol: float, default: 1e-4

Tolerance in norm used when finding the collapse.

norm_steps: int, default: 5

Maximum number of tries to find the collapse.

improved_sampling: Bool, default: False

Whether to use the improved sampling algorithm of Abdelhafez et al. PRA (2019)

completeness_rtol: float, default: 1e-5

Used in determining whether the given Lindblad operators satisfy a certain completeness relation. If they do not, an additional Lindblad operator is added automatically (with zero rate).

completeness_atol: float, default: 1e-8

Used in determining whether the given Lindblad operators satisfy a certain completeness relation. If they do not, an additional Lindblad operator is added automatically (with zero rate).

martingale_quad_limit: float or int, default: 100

An upper bound on the number of subintervals used in the adaptive integration of the martingale.

Note that the ‘improved_sampling’ option is not currently supported.

rate(t, i)[source]

Return the i’th unshifted rate at time t.

Parameters:
tfloat

The time at which to calculate the rate.

iint

Which rate to calculate.

Returns:
ratefloat

The value of rate i at time t.

rate_shift(t)[source]

Return the rate shift at time t.

The rate shift is 2 * abs(min([0, rate_1(t), rate_2(t), ...])).

Parameters:
tfloat

The time at which to calculate the rate shift.

Returns:
rate_shiftfloat

The rate shift amount.

run(
state: Qobj,
tlist: ArrayLike,
ntraj: int = 1,
*,
args: dict[str, Any] = None,
**kwargs,
)[source]

Do the evolution of the Quantum system.

For a state at time tlist[0] do the evolution as directed by rhs and for each time in tlist store the state and/or expectation values in a MultiTrajResult. The evolution method and stored results are determined by options.

Parameters:
stateQobj

Initial state of the evolution.

tlistlist of double

Time for which to save the results (state and/or expect) of the evolution. The first element of the list is the initial time of the evolution. Time in the list must be in increasing order, but does not need to be uniformly distributed.

ntrajint

Number of trajectories to add.

argsdict, optional

Change the args of the rhs for the evolution.

e_opslist

list of Qobj or QobjEvo to compute the expectation values. Alternatively, function[s] with the signature f(t, state) -> expect can be used.

timeoutfloat, optional

Maximum time in seconds for the trajectories to run. Once this time is reached, the simulation will end even if the number of trajectories is less than ntraj. The map function, set in options, can interupt the running trajectory or wait for it to finish. Set to an arbitrary high number to disable.

target_tol{float, tuple, list}, optional

Target tolerance of the evolution. The evolution will compute trajectories until the error on the expectation values is lower than this tolerance. The maximum number of trajectories employed is given by ntraj. The error is computed using jackknife resampling. target_tol can be an absolute tolerance or a pair of absolute and relative tolerance, in that order. Lastly, it can be a list of pairs of (atol, rtol) for each e_ops.

seeds{int, SeedSequence, list}, optional

Seed or list of seeds for each trajectories.

Returns:
resultsMultiTrajResult

Results of the evolution. States and/or expect will be saved. You can control the saved data in the options.

sqrt_shifted_rate(t, i)[source]

Return the square root of the i’th shifted rate at time t.

Parameters:
tfloat

The time at wich to calculate the shifted rate.

iint

Which shifted rate to calculate.

Returns:
ratefloat

The square root of the shifted value of rate i at time t.

start(
state: Qobj,
t0: float,
seed: int | SeedSequence = None,
)[source]

Set the initial state and time for a step evolution.

Parameters:
stateQobj

Initial state of the evolution.

t0double

Initial time of the evolution.

seedint, SeedSequence, list, optional

Seed for the random number generator. It can be a single seed used to spawn seeds for each trajectory or a list of seed, one for each trajectory.

Notes

When using step evolution, only one trajectory can be computed at once.

step(
t: float,
*,
args: dict[str, Any] = None,
copy: bool = True,
) Qobj[source]

Evolve the state to t and return the state as a Qobj.

Parameters:
tdouble

Time to evolve to, must be higher than the last call.

argsdict, optional

Update the args of the system. The change is effective from the beginning of the interval. Changing args can slow the evolution.

copybool, default: True

Whether to return a copy of the data or the data in the ODE solver.

property sys_dims

Dimensions of the space that the system use:

qutip.basis(sovler.dims) will create a state with proper dimensions for this solver.

class McResult(
e_ops,
options: MultiTrajResultOptions,
*,
solver=None,
stats=None,
**kw,
)[source]

Bases: _McBaseResult

Class for storing Monte-Carlo solver results.

Parameters:
e_opsQobj, QobjEvo, function or list or dict of these

The e_ops parameter defines the set of values to record at each time step t. If an element is a Qobj or QobjEvo the value recorded is the expectation value of that operator given the state at t. If the element is a function, f, the value recorded is f(t, state).

The values are recorded in the .expect attribute of this result object. .expect is a list, where each item contains the values of the corresponding e_op.

optionsSolverResultsOptions

The options for this result class.

solverstr or None

The name of the solver generating these results.

statsdict

The stats generated by the solver while producing these results. Note that the solver may update the stats directly while producing results. Must include a value for “num_collapse”.

kwdict

Additional parameters specific to a result sub-class.

Attributes:
collapselist

For each run, a list of every collapse as a tuple of the time it happened and the corresponding c_ops index.

property photocurrent

Average photocurrent or measurement of the evolution.

property runs_photocurrent

Photocurrent or measurement of each runs.

class NmmcResult(
e_ops,
options: MultiTrajResultOptions,
*,
solver=None,
stats=None,
**kw,
)[source]

Bases: _McBaseResult

Class for storing the results of the non-Markovian Monte-Carlo solver.

Parameters:
e_opsQobj, QobjEvo, function or list or dict of these

The e_ops parameter defines the set of values to record at each time step t. If an element is a Qobj or QobjEvo the value recorded is the expectation value of that operator given the state at t. If the element is a function, f, the value recorded is f(t, state).

The values are recorded in the .expect attribute of this result object. .expect is a list, where each item contains the values of the corresponding e_op.

optionsSolverResultsOptions

The options for this result class.

solverstr or None

The name of the solver generating these results.

statsdict

The stats generated by the solver while producing these results. Note that the solver may update the stats directly while producing results. Must include a value for “num_collapse”.

kwdict

Additional parameters specific to a result sub-class.

Attributes:
average_tracelist

Refers to average_trace or runs_trace, depending on whether keep_runs_results is set in the options.

std_tracelist

Refers to average_trace or runs_trace, depending on whether keep_runs_results is set in the options.

runs_tracelist of lists

For each recorded trajectory, the trace at each time. Only present if keep_runs_results is set in the options.

property average_trace

Refers to average_trace or runs_trace, depending on whether keep_runs_results is set in the options.

merge(other, p=None)[source]

Merges two multi-trajectory results.

If this result represent an ensemble \(\rho\), and other represents an ensemble \(\rho'\), then the merged result represents the ensemble

\[\rho_{\mathrm{merge}} = p \rho + (1 - p) \rho'\]

where p is a parameter between 0 and 1. Its default value is \(p_{\textrm{def}} = N / (N + N')\), N and N’ being the number of trajectories in the two result objects.

Parameters:
otherMultiTrajResult

The multi-trajectory result to merge with this one

pfloat [optional]

The relative weight of this result in the combination. By default, will be chosen such that all trajectories contribute equally to the merged result.

property std_trace

Refers to average_trace or runs_trace, depending on whether keep_runs_results is set in the options.

property trace

Refers to average_trace or runs_trace, depending on whether keep_runs_results is set in the options.

Bloch-Redfield Master Equation

This module provides solvers for the Lindblad master equation and von Neumann equation.

brmesolve(
H: QobjEvoLike,
psi0: Qobj,
tlist: ArrayLike,
a_ops: list[tuple[QobjEvoLike, CoefficientLike]] = None,
sec_cutoff: float = 0.1,
*_pos_args,
c_ops: list[QobjEvoLike] = None,
e_ops: EopsLike | list[EopsLike] | dict[Any, EopsLike] = None,
args: dict[str, Any] = None,
options: dict[str, Any] = None,
**kwargs,
)[source]

Solves for the dynamics of a system using the Bloch-Redfield master equation, given an input Hamiltonian, Hermitian bath-coupling terms and their associated spectral functions, as well as possible Lindblad collapse operators.

Parameters:
HQobj, QobjEvo

Possibly time-dependent system Liouvillian or Hamiltonian as a Qobj or QobjEvo. list of [Qobj, Coefficient] or callable that can be made into QobjEvo are also accepted.

psi0: :obj:`.Qobj`

Initial density matrix or state vector (ket).

tlistarray_like

List of times for evaluating evolution

a_opslist of (a_op, spectra)

Nested list of system operators that couple to the environment, and the corresponding bath spectra.

a_opQobj, QobjEvo, QobjEvo compatible format

The operator coupling to the environment. Must be hermitian.

spectraCoefficient, str, func

The corresponding bath spectral responce. Can be a Coefficient using an ‘w’ args, a function of the frequence or a string. Coefficient build from a numpy array are understood as a function of w instead of t. Function are expected to be of the signature f(w) or f(t, w, **args).

The spectra function can depend on t if the corresponding a_op is a QobjEvo.

Example:

a_ops = [
    (a+a.dag(), ('w>0', args={"w": 0})),
    (QobjEvo(a+a.dag()), 'w > exp(-t)'),
    ([[b+b.dag(), lambda t: ...]], lambda w: ...)),
    (c+c.dag(), SpectraCoefficient(coefficient(array, tlist=ws))),
]

Note

Cubic_Spline has been replaced by:

spline = qutip.coefficient(array, tlist=times)

Whether the a_ops is time dependent is decided by the type of the operator: Qobj vs QobjEvo instead of the type of the spectra.

sec_cutofffloat, default: 0.1

Cutoff for secular approximation. Use -1 if secular approximation is not used when evaluating bath-coupling terms.

*_pos_args
Temporary shim to update the signature from
(..., a_ops, e_ops, c_ops, args, sec_cutoff, options)
to
(..., a_ops, sec_cutoff, *, e_ops, c_ops, args, options)
making e_ops, c_ops, args and options keyword only parameter from qutip 5.3.
c_opslist of (QobjEvo, QobjEvo compatible format), optional

List of collapse operators.

argsdict, optional

Dictionary of parameters for time-dependent Hamiltonians and collapse operators. The key w is reserved for the spectra function.

e_opslist, dict, Qobj or callback function, optional

Single operator, or list or dict of operators, for which to evaluate expectation values. Operator can be Qobj, QobjEvo or callables with the signature f(t: float, state: Qobj) -> Any. Callable signature must be, f(t: float, state: Qobj).

optionsdict, optional

Dictionary of options for the solver.

  • store_final_state : bool
    Whether or not to store the final state of the evolution in the result class.
  • store_states : bool, None
    Whether or not to store the state vectors or density matrices. On None the states will be saved if no expectation operators are given.
  • normalize_output : bool
    Normalize output state to hide ODE numerical errors. Only normalize the state if the initial state is already normalized.
  • progress_bar : str {‘text’, ‘enhanced’, ‘tqdm’, ‘’}
    How to present the solver progress. ‘tqdm’ uses the python module of the same name and raise an error if not installed. Empty string or False will disable the bar.
  • progress_kwargs : dict
    kwargs to pass to the progress_bar. Qutip’s bars use chunk_size.
  • tensor_type : str [‘sparse’, ‘dense’, ‘data’]
    Which data type to use when computing the brtensor. With a cutoff ‘sparse’ is usually the most efficient.
  • sparse_eigensolver : bool {False} Whether to use the sparse eigensolver
  • method : str [“adams”, “bdf”, “lsoda”, “dop853”, “vern9”, etc.] Which differential equation integration method to use.
  • atol, rtol : float
    Absolute and relative tolerance of the ODE integrator.
  • nsteps : int
    Maximum number of (internally defined) steps allowed in one tlist step.
  • max_step : float, 0
    Maximum lenght of one internal step. When using pulses, it should be less than half the width of the thinnest pulse.

Other options could be supported depending on the integration method, see Integrator.

Returns:
result: Result

An instance of the class qutip.solver.Result, which contains either an array of expectation values, for operators given in e_ops, or a list of states for the times specified by tlist.

class BRSolver(
H: Qobj | QobjEvo,
a_ops: list[tuple[Qobj | QobjEvo, Coefficient]],
c_ops: Qobj | QobjEvo | list[Qobj | QobjEvo] = None,
sec_cutoff: float = 0.1,
*,
options: dict[str, Any] = None,
)[source]

Bases: Solver

Bloch Redfield equation evolution of a density matrix for a given Hamiltonian and set of bath coupling operators.

Parameters:
HQobj, QobjEvo

Possibly time-dependent system Liouvillian or Hamiltonian as a Qobj or QobjEvo. list of [Qobj, Coefficient] or callable that can be made into QobjEvo are also accepted.

a_opslist of (a_op, spectra)

Nested list of system operators that couple to the environment, and the corresponding bath spectra.

a_opQobj, QobjEvo

The operator coupling to the environment. Must be hermitian.

spectraCoefficient

The corresponding bath spectra. As a Coefficient using an ‘w’ args. Can depend on t only if a_op is a QobjEvo. SpectraCoefficient can be used to conver a coefficient depending on t to one depending on w.

Example:

a_ops = [
    (a+a.dag(), coefficient('w>0', args={'w':0})),
    (QobjEvo([b+b.dag(), lambda t: ...]),
     coefficient(lambda t, w: ...), args={"w": 0}),
    (c+c.dag(), SpectraCoefficient(coefficient(array, tlist=ws))),
]
c_opslist of Qobj, QobjEvo

Single collapse operator, or list of collapse operators, or a list of Lindblad dissipator. None is equivalent to an empty list.

optionsdict, optional

Options for the solver, see BRSolver.options and Integrator for a list of all options.

sec_cutofffloat {0.1}

Cutoff for secular approximation. Use -1 if secular approximation is not used when evaluating bath-coupling terms.

Attributes:
stats: dict

Diverse diagnostic statistics of the evolution.

classmethod ExpectFeedback(
operator: Qobj | QobjEvo,
default: Any = 0.0,
)

Expectation value of the instantaneous state of the evolution to be used by a time-dependent operator.

When used as an args:

QobjEvo([op, func], args={"E0": Solver.ExpectFeedback(oper)})

The func will receive expect(oper, state) as E0 during the evolution.

Parameters:
operatorQobj, QobjEvo

Operator to compute the expectation values of.

defaultfloat, default0.

Initial value to be used at setup.

classmethod StateFeedback(default=None, raw_data=False)[source]

State of the evolution to be used in a time-dependent operator.

When used as an args:

QobjEvo([op, func], args={"state": BRMESolver.StateFeedback()})

The func will receive the density matrix as state during the evolution.

Note

The state will not be in the lab basis, but in the evolution basis.

Parameters:
defaultQobj or qutip.core.data.Data, defaultNone

Initial value to be used at setup of the system.

raw_databool, defaultFalse

If True, the raw matrix will be passed instead of a Qobj. For density matrices, the matrices can be column stacked or square depending on the integration method.

property options

Options for bloch redfield solver:

store_final_state: bool, default: False

Whether or not to store the final state of the evolution in the result class.

store_states: bool, default: None

Whether or not to store the state vectors or density matrices. On None the states will be saved if no expectation operators are given.

normalize_output: bool, default: False

Normalize output state to hide ODE numerical errors.

progress_bar: str {‘text’, ‘enhanced’, ‘tqdm’, ‘’}, default: “”

How to present the solver progress. ‘tqdm’ uses the python module of the same name and raise an error if not installed. Empty string or False will disable the bar.

progress_kwargs: dict, default: {“chunk_size”:10}

Arguments to pass to the progress_bar. Qutip’s bars use chunk_size.

tensor_type: str [‘sparse’, ‘dense’, ‘data’], default: “sparse”

Which data type to use when computing the brtensor. With a cutoff ‘sparse’ is usually the most efficient.

sparse_eigensolver: bool, default: False

Whether to use the sparse eigensolver

method: str, default: “adams”

Which ODE integrator methods are supported.

run(
state0: Qobj,
tlist: ArrayLike,
*,
e_ops: EopsLike | list[EopsLike] | dict[Any, EopsLike] = None,
args: dict[str, Any] = None,
) Result

Do the evolution of the Quantum system.

For a state0 at time tlist[0] do the evolution as directed by rhs and for each time in tlist store the state and/or expectation values in a Result. The evolution method and stored results are determined by options.

Parameters:
state0Qobj

Initial state of the evolution.

tlistlist of double

Time for which to save the results (state and/or expect) of the evolution. The first element of the list is the initial time of the evolution. Each times of the list must be increasing, but does not need to be uniformy distributed.

argsdict, optional

Change the args of the rhs for the evolution.

e_opsQobj, QobjEvo, callable, list, or dict optional

Single, list or dict of Qobj, QobjEvo or callable to compute the expectation values. Function[s] must have the signature f(t : float, state : Qobj) -> expect.

Returns:
resultsResult

Results of the evolution. States and/or expect will be saved. You can control the saved data in the options.

start(state0: Qobj, t0: Number) None

Set the initial state and time for a step evolution.

Parameters:
state0Qobj

Initial state of the evolution.

t0double

Initial time of the evolution.

step(
t: Number,
*,
args: dict[str, Any] = None,
copy: bool = True,
) Qobj

Evolve the state to t and return the state as a Qobj.

Parameters:
tdouble

Time to evolve to, must be higher than the last call.

argsdict, optional {None}

Update the args of the system. The change is effective from the beginning of the interval. Changing args can slow the evolution.

copybool, optional {True}

Whether to return a copy of the data or the data in the ODE solver.

Notes

The state must be initialized first by calling start or run. If run is called, step will continue from the last time and state obtained.

property sys_dims

Dimensions of the space that the system use:

qutip.basis(sovler.dims) will create a state with proper dimensions for this solver.

Floquet States and Floquet-Markov Master Equation

floquet_tensor(
H: QobjEvo | FloquetBasis,
c_ops: list[Qobj],
spectra_cb: list[Callable[[float], complex]],
T: float = 0,
w_th: float = 0.0,
kmax: int = 5,
nT: int = 100,
) Qobj[source]

Construct a tensor that represents the master equation in the floquet basis.

Simplest RWA approximation [Grifoni et al, Phys.Rep. 304 229 (1998)]

Parameters:
HQobjEvo, FloquetBasis

Periodic Hamiltonian a floquet basis system.

c_opslist of Qobj

list of collapse operators.

spectra_cblist callback functions

List of callback functions that compute the noise power spectrum as a function of frequency for the collapse operators in c_ops.

Tfloat, optional

The period of the time-dependence of the hamiltonian. Optional if H is a FloquetBasis object.

w_thfloat, default: 0.0

The temperature in units of frequency.

kmaxint, default: 5

The truncation of the number of sidebands (default 5).

nTint, default: 100

The number of integration steps (for calculating X) within one period.

Returns:
outputarray

The Floquet-Markov master equation tensor R.

fmmesolve(
H: QobjEvoLike | FloquetBasis,
rho0: Qobj,
tlist: ArrayLike,
c_ops: list[Qobj] = None,
spectra_cb: list[Callable[[float], complex]] = None,
T: float = 0.0,
w_th: float = 0.0,
*pos_args,
e_ops: EopsLike | list[EopsLike] | dict[Any, EopsLike] = None,
args: dict[str, Any] = None,
options: dict[str, Any] = None,
) FloquetResult[source]

Solve the dynamics for the system using the Floquet-Markov master equation.

Parameters:
HQobj, QobjEvo, QobjEvo compatible format.

Periodic system Hamiltonian as QobjEvo. List of [Qobj, Coefficient] or callable that can be made into QobjEvo are also accepted.

rho0 / psi0Qobj

Initial density matrix or state vector (ket).

tlistlist / array

List of times for \(t\).

c_opslist of Qobj, optional

List of collapse operators. Time dependent collapse operators are not supported. Fall back on fsesolve if not provided.

e_opslist of Qobj / callback function, optional

Single operator, or list or dict of operators, for which to evaluate expectation values. Operator can be Qobj, QobjEvo or callables with the signature f(t: float, state: Qobj) -> Any. See expect for more detail of operator expectation. The states are reverted to the lab basis before computing the expectation values.

spectra_cblist callback functions, default: lambda w: (w > 0)

List of callback functions that compute the noise power spectrum as a function of frequency for the collapse operators in c_ops.

Tfloat, default=tlist[-1]

The period of the time-dependence of the hamiltonian. The default value 0 indicates that the ‘tlist’ spans a single period of the driving.

w_thfloat, default: 0.0

The temperature of the environment in units of frequency. For example, if the Hamiltonian written in units of 2pi GHz, and the temperature is given in K, use the following conversion:

temperature = 25e-3 # unit K h = 6.626e-34 kB = 1.38e-23 args[‘w_th’] = temperature * (kB / h) * 2 * pi * 1e-9

argsdict, optional

Dictionary of parameters for time-dependent Hamiltonian

optionsdict, optional

Dictionary of options for the solver.

  • store_final_state : bool
    Whether or not to store the final state of the evolution in the result class.
  • store_states : bool, None
    Whether or not to store the state vectors or density matrices. On None the states will be saved if no expectation operators are given.
  • store_floquet_states : bool
    Whether or not to store the density matrices in the floquet basis in result.floquet_states.
  • normalize_output : bool
    Normalize output state to hide ODE numerical errors. Only normalize the state if the initial state is already normalized.
  • progress_bar : str {‘text’, ‘enhanced’, ‘tqdm’, ‘’}
    How to present the solver progress. ‘tqdm’ uses the python module of the same name and raise an error if not installed. Empty string or False will disable the bar.
  • progress_kwargs : dict
    kwargs to pass to the progress_bar. Qutip’s bars use chunk_size.
  • method : str [“adams”, “bdf”, “lsoda”, “dop853”, “vern9”, etc.]
    Which differential equation integration method to use.
  • atol, rtol : float
    Absolute and relative tolerance of the ODE integrator.
  • nsteps : int
    Maximum number of (internally defined) steps allowed in one tlist step.
  • max_step : float
    Maximum lenght of one internal step. When using pulses, it should be less than half the width of the thinnest pulse.

Other options could be supported depending on the integration method, see Integrator.

Returns:
result: Result

An instance of the class Result, which contains the expectation values for the times specified by tlist, and/or the state density matrices corresponding to the times.

fsesolve(
H: QobjEvoLike | FloquetBasis,
psi0: Qobj,
tlist: ArrayLike,
T: float = 0.0,
*pos_args,
e_ops: EopsLike | list[EopsLike] | dict[Any, EopsLike] = None,
args: dict[str, Any] = None,
options: dict[str, Any] = None,
) Result[source]

Solve the Schrodinger equation using the Floquet formalism.

Parameters:
HQobj, QobjEvo, QobjEvo compatible format.

Periodic system Hamiltonian as QobjEvo. List of [Qobj, Coefficient] or callable that can be made into QobjEvo are also accepted.

psi0Qobj

Initial state vector (ket). If an operator is provided,

tlistlist / array

List of times for \(t\).

Tfloat, default=tlist[-1]

The period of the time-dependence of the hamiltonian.

e_opslist or dict of Qobj / callback function, optional

Single operator, or list or dict of operators, for which to evaluate expectation values. Operator can be Qobj, QobjEvo or callables with the signature f(t: float, state: Qobj) -> Any. See expect for more detail of operator expectation.

argsdictionary, optional

Dictionary with variables required to evaluate H.

optionsdict, optional

Options for the results.

  • store_final_state : bool
    Whether or not to store the final state of the evolution in the result class.
  • store_states : bool, None
    Whether or not to store the state vectors or density matrices. On None the states will be saved if no expectation operators are given.
  • normalize_output : bool
    Normalize output state to hide ODE numerical errors. Only normalize the state if the initial state is already normalized.
Returns:
outputResult

An instance of the class Result, which contains either an array of expectation values or an array of state vectors, for the times specified by tlist.

class FMESolver(
floquet_basis: FloquetBasis,
a_ops: list[tuple[Qobj, Callable[[float], float]]],
w_th: float = 0.0,
*,
kmax: int = 5,
nT: int = None,
options: dict[str, Any] = None,
)[source]

Bases: MESolver

Solver for the Floquet-Markov master equation.

Note

Operators (c_ops and e_ops) are in the laboratory basis.

Parameters:
floquet_basisFloquetBasis

The system Hamiltonian wrapped in a FloquetBasis object. Choosing a different integrator for the floquet_basis than for the evolution of the floquet state can improve the performance.

a_opslist of tuple(Qobj, callable)

List of collapse operators and the corresponding function for the noise power spectrum. The collapse operator must be a Qobj and cannot be time dependent. The spectrum function must take and return an numpy array.

w_thfloat

The temperature of the environment in units of Hamiltonian frequency.

kmaxint [5]

The truncation of the number of sidebands..

nTint [20*kmax]

The number of integration steps (for calculating X) within one period.

optionsdict, optional

Options for the solver, see FMESolver.options and Integrator for a list of all options.

classmethod ExpectFeedback()[source]

Expect of the state of the evolution to be used in a time-dependent operator.

Not not implemented for FMESolver

classmethod StateFeedback()[source]

State of the evolution to be used in a time-dependent operator.

Not not implemented for FMESolver

property options: dict

Solver’s options:

store_final_state: bool, default: False

Whether or not to store the final state of the evolution in the result class.

store_states: bool, default: None

Whether or not to store the state vectors or density matrices. On None the states will be saved if no expectation operators are given.

normalize_output: bool, default: True

Normalize output state to hide ODE numerical errors.

progress_bar: str {“text”, “enhanced”, “tqdm”, “”}, default: “”

How to present the solver progress. ‘tqdm’ uses the python module of the same name and raise an error if not installed. Empty string or False will disable the bar.

progress_kwargs: dict, default: {“chunk_size”: 10}

Arguments to pass to the progress_bar. Qutip’s bars use chunk_size.

method: str, default: “adams”

Which ordinary differential equation integration method to use.

run(
state0: Qobj,
tlist: ArrayLike,
*,
floquet: bool = False,
args: dict[str, Any] = None,
e_ops: EopsLike | list[EopsLike] | dict[Any, EopsLike] = None,
) FloquetResult[source]

Calculate the evolution of the quantum system.

For a state0 at time tlist[0] do the evolution as directed by rhs and for each time in tlist store the state and/or expectation values in a Result. The evolution method and stored results are determined by options.

Parameters:
state0Qobj

Initial state of the evolution.

tlistlist of double

Time for which to save the results (state and/or expect) of the evolution. The first element of the list is the initial time of the evolution. Each times of the list must be increasing, but does not need to be uniformy distributed.

floquetbool, optional {False}

Whether the initial state in the floquet basis or laboratory basis.

argsdict, optional

Not supported

e_opslist or dict, optional

List or dict of Qobj, QobjEvo or callable to compute the expectation values. Function[s] must have the signature f(t : float, state : Qobj) -> expect.

Returns:
resultsFloquetResult

Results of the evolution. States and/or expect will be saved. You can control the saved data in the options.

start(
state0: Qobj,
t0: float,
*,
floquet: bool = False,
) None[source]

Set the initial state and time for a step evolution. options for the evolutions are read at this step.

Parameters:
state0Qobj

Initial state of the evolution.

t0double

Initial time of the evolution.

floquetbool, optional {False}

Whether the initial state is in the floquet basis or laboratory basis.

step(
t: float,
*,
args: dict[str, Any] = None,
copy: bool = True,
floquet: bool = False,
) Qobj[source]

Evolve the state to t and return the state as a Qobj.

Parameters:
tdouble

Time to evolve to, must be higher than the last call.

copybool, optional {True}

Whether to return a copy of the data or the data in the ODE solver.

floquetbool, optional {False}

Whether to return the state in the floquet basis or laboratory basis.

argsdict, optional {None}

Not supported

Notes

The state must be initialized first by calling start or run. If run is called, step will continue from the last time and state obtained.

property sys_dims

Dimensions of the space that the system use:

qutip.basis(sovler.dims) will create a state with proper dimensions for this solver.

class FloquetBasis(
H: QobjEvoLike,
T: float,
args: dict[str, Any] = None,
options: dict[str, Any] = None,
sparse: bool = False,
sort: bool = True,
precompute: ArrayLike = None,
times: ArrayLike = None,
)[source]

Utility to compute floquet modes and states.

Attributes:
UPropagator

The propagator of the Hamiltonian over one period.

evecsData

Matrix where each column is an initial Floquet mode.

e_quasinp.ndarray[float]

The quasi energies of the Hamiltonian.

from_floquet_basis(
floquet_basis: QobjOrData,
t: float = 0,
) QobjOrData[source]

Transform a ket or density matrix from the Floquet basis at time t to the lab basis.

Parameters:
floquet_basisQobj, Data

Initial state in the Floquet basis at time t. May be either a ket or density matrix.

tfloat [0]

The time at which to evaluate the Floquet states.

Returns:
outputQobj, Data

The state in the lab basis. The return type is the same as the type of the input state.

mode(
t: float,
data: Literal[False],
) Qobj[source]
mode(
t: float,
data: Literal[True],
) Data

Calculate the Floquet modes at time t.

Parameters:
tfloat

The time for which to evaluate the Floquet mode.

databool [False]

Whether to return the states as a single data matrix or a list of ket states.

Returns:
outputlist[Qobj], Data

A list of Floquet states for the time t or the states as column in a single matrix.

state(
t: float,
data: Literal[False],
) Qobj[source]
state(
t: float,
data: Literal[True],
) Data

Evaluate the floquet states at time t.

Parameters:
tfloat

The time for which to evaluate the Floquet states.

databool [False]

Whether to return the states as a single data matrix or a list of ket states.

Returns:
outputlist[Qobj], Data

A list of Floquet states for the time t or the states as column in a single matrix.

to_floquet_basis(
lab_basis: QobjOrData,
t: float = 0,
) QobjOrData[source]

Transform a ket or density matrix in the lab basis to the Floquet basis at time t.

Parameters:
lab_basisQobj, Data

Initial state in the lab basis.

tfloat [0]

The time at which to evaluate the Floquet states.

Returns:
outputQobj, Data

The state in the Floquet basis. The return type is the same as the type of the input state.

Stochastic Schrödinger Equation and Master Equation

smesolve(
H: QobjEvoLike,
rho0: Qobj,
tlist: ArrayLike,
c_ops: Qobj | QobjEvo | Sequence[QobjEvoLike] = (),
sc_ops: Qobj | QobjEvo | Sequence[QobjEvoLike] = (),
heterodyne: bool = False,
*,
e_ops: EopsLike | list[EopsLike] | dict[Any, EopsLike] = None,
args: dict[str, Any] = None,
ntraj: int = 500,
options: dict[str, Any] = None,
seeds: int | SeedSequence | Sequence[int | SeedSequence] = None,
target_tol: float | tuple[float, float] | list[tuple[float, float]] = None,
timeout: float = None,
**kwargs,
) StochasticResult[source]

Solve stochastic master equation.

Parameters:
HQobj, QobjEvo, QobjEvo compatible format.

System Hamiltonian as a Qobj or QobjEvo for time-dependent Hamiltonians. List of [Qobj, Coefficient] or callable that can be made into QobjEvo are also accepted.

rho0Qobj

Initial density matrix or state vector (ket).

tlistlist / array

List of times for \(t\).

c_opslist of (QobjEvo, QobjEvo compatible format), optional

Deterministic collapse operator which will contribute with a standard Lindblad type of dissipation.

sc_opslist of (QobjEvo, QobjEvo compatible format)

List of stochastic collapse operators.

e_opsQobj, callable, list or dict, optional

Single operator, or list or dict of operators, for which to evaluate expectation values. Operator can be Qobj, QobjEvo or callables with the signature f(t: float, state: Qobj) -> Any.

argsdict, optional

Dictionary of parameters for time-dependent Hamiltonians and collapse operators.

ntrajint, default: 500

Number of trajectories to compute.

heterodynebool, default: False

Whether to use heterodyne or homodyne detection.

seedsint, SeedSequence, list, optional

Seed for the random number generator. It can be a single seed used to spawn seeds for each trajectory or a list of seeds, one for each trajectory. Seeds are saved in the result and they can be reused with:

seeds=prev_result.seeds

When using a parallel map, the trajectories can be re-ordered.

target_tol{float, tuple, list}, optional

Target tolerance of the evolution. The evolution will compute trajectories until the error on the expectation values is lower than this tolerance. The maximum number of trajectories employed is given by ntraj. The error is computed using jackknife resampling. target_tol can be an absolute tolerance or a pair of absolute and relative tolerance, in that order. Lastly, it can be a list of pairs of (atol, rtol) for each e_ops.

timeoutfloat, optional

Maximum time for the evolution in second. When reached, no more trajectories will be computed. Overwrite the option of the same name.

optionsdict, optional

Dictionary of options for the solver.

  • store_final_state : bool
    Whether or not to store the final state of the evolution in the result class.
  • store_states : bool, None
    Whether or not to store the state vectors or density matrices. On None the states will be saved if no expectation operators are given.
  • store_measurement: str, {‘start’, ‘middle’, ‘end’, ‘’}
    Whether and how to store the measurement for each trajectories. ‘start’, ‘middle’, ‘end’ indicate when in the interval the expectation value of the m_ops is taken.
  • keep_runs_results : bool
    Whether to store results from all trajectories or just store the averages.
  • normalize_output : bool
    Normalize output state to hide ODE numerical errors. Only normalize the state if the initial state is already normalized.
  • progress_bar : str {‘text’, ‘enhanced’, ‘tqdm’, ‘’}
    How to present the solver progress. ‘tqdm’ uses the python module of the same name and raise an error if not installed. Empty string or False will disable the bar.
  • progress_kwargs : dict
    kwargs to pass to the progress_bar. Qutip’s bars use chunk_size.
  • method : str
    Which stochastic differential equation integration method to use. Main ones are {“euler”, “rouchon”, “platen”, “taylor1.5_imp”}
  • map : str {“serial”, “parallel”, “loky”, “mpi”}
    How to run the trajectories. “parallel” uses the multiprocessing module to run in parallel while “loky” and “mpi” use the “loky” and “mpi4py” modules to do so.
  • num_cpus : NoneType, int
    Number of cpus to use when running in parallel. None detect the number of available cpus.
  • dt : float
    The finite steps lenght for the Stochastic integration method. Default change depending on the integrator.

Additional options are listed under options. More options may be available depending on the selected differential equation integration method, see SIntegrator.

Returns:
output: Result

An instance of the class Result.

ssesolve(
H: QobjEvoLike,
psi0: Qobj,
tlist: ArrayLike,
sc_ops: QobjEvoLike | Sequence[QobjEvoLike] = (),
heterodyne: bool = False,
*,
e_ops: EopsLike | list[EopsLike] | dict[Any, EopsLike] = None,
args: dict[str, Any] = None,
ntraj: int = 500,
options: dict[str, Any] = None,
seeds: int | SeedSequence | Sequence[int | SeedSequence] = None,
target_tol: float | tuple[float, float] | list[tuple[float, float]] = None,
timeout: float = None,
**kwargs,
) StochasticResult[source]

Solve stochastic Schrodinger equation.

Parameters:
HQobj, QobjEvo, QobjEvo compatible format.

System Hamiltonian as a Qobj or QobjEvo for time-dependent Hamiltonians. List of [Qobj, Coefficient] or callable that can be made into QobjEvo are also accepted.

psi0Qobj

Initial state vector (ket).

tlistlist / array

List of times for \(t\).

sc_opslist of (QobjEvo, QobjEvo compatible format)

List of stochastic collapse operators.

e_opsQobj, callable, list or dict, optional

Single operator, or list or dict of operators, for which to evaluate expectation values. Operator can be Qobj, QobjEvo or callables with the signature f(t: float, state: Qobj) -> Any.

argsdict, optional

Dictionary of parameters for time-dependent Hamiltonians and collapse operators.

ntrajint, default: 500

Number of trajectories to compute.

heterodynebool, default: False

Whether to use heterodyne or homodyne detection.

seedsint, SeedSequence, list, optional

Seed for the random number generator. It can be a single seed used to spawn seeds for each trajectory or a list of seeds, one for each trajectory. Seeds are saved in the result and they can be reused with:

seeds=prev_result.seeds
target_tol{float, tuple, list}, optional

Target tolerance of the evolution. The evolution will compute trajectories until the error on the expectation values is lower than this tolerance. The maximum number of trajectories employed is given by ntraj. The error is computed using jackknife resampling. target_tol can be an absolute tolerance or a pair of absolute and relative tolerance, in that order. Lastly, it can be a list of pairs of (atol, rtol) for each e_ops.

timeoutfloat, optional

Maximum time for the evolution in second. When reached, no more trajectories will be computed. Overwrite the option of the same name.

optionsdict, optional

Dictionary of options for the solver.

  • store_final_state : bool
    Whether or not to store the final state of the evolution in the result class.
  • store_states : bool, None
    Whether or not to store the state vectors or density matrices. On None the states will be saved if no expectation operators are given.
  • store_measurement: str, {‘start’, ‘middle’, ‘end’, ‘’}
    Whether and how to store the measurement for each trajectories. ‘start’, ‘middle’, ‘end’ indicate when in the interval the expectation value of the m_ops is taken.
  • keep_runs_results : bool
    Whether to store results from all trajectories or just store the averages.
  • normalize_output : bool
    Normalize output state to hide ODE numerical errors. Only normalize the state if the initial state is already normalized.
  • progress_bar : str {‘text’, ‘enhanced’, ‘tqdm’, ‘’}
    How to present the solver progress. ‘tqdm’ uses the python module of the same name and raise an error if not installed. Empty string or False will disable the bar.
  • progress_kwargs : dict
    kwargs to pass to the progress_bar. Qutip’s bars use chunk_size.
  • method : str
    Which stochastic differential equation integration method to use. Main ones are {“euler”, “rouchon”, “platen”, “taylor1.5_imp”}
  • map : str {“serial”, “parallel”, “loky”, “mpi”}
    How to run the trajectories. “parallel” uses the multiprocessing module to run in parallel while “loky” and “mpi” use the “loky” and “mpi4py” modules to do so.
  • num_cpus : NoneType, int
    Number of cpus to use when running in parallel. None detect the number of available cpus.
  • dt : float
    The finite steps lenght for the Stochastic integration method. Default change depending on the integrator.

Additional options are listed under options. More options may be available depending on the selected differential equation integration method, see SIntegrator.

Returns:
output: Result

An instance of the class Result.

class SMESolver(
H: Qobj | QobjEvo,
sc_ops: Sequence[Qobj | QobjEvo],
heterodyne: bool,
*,
c_ops: Sequence[Qobj | QobjEvo] = (),
options: dict[str, Any] = None,
)[source]

Stochastic Master Equation Solver.

Parameters:
HQobj, QobjEvo, QobjEvo compatible format.

System Hamiltonian as a Qobj or QobjEvo for time-dependent Hamiltonians. List of [Qobj, Coefficient] or callable that can be made into QobjEvo are also accepted.

sc_opslist of (QobjEvo, QobjEvo compatible format)

List of stochastic collapse operators.

heterodynebool, default: False

Whether to use heterodyne or homodyne detection.

optionsdict, optional

Options for the solver, see SMESolver.options and SIntegrator for a list of all options.

classmethod ExpectFeedback(
operator: Qobj | QobjEvo,
default: Any = 0.0,
)

Expectation value of the instantaneous state of the evolution to be used by a time-dependent operator.

When used as an args:

QobjEvo([op, func], args={"E0": Solver.ExpectFeedback(oper)})

The func will receive expect(oper, state) as E0 during the evolution.

Parameters:
operatorQobj, QobjEvo

Operator to compute the expectation values of.

defaultfloat, default0.

Initial value to be used at setup.

classmethod StateFeedback(
default: Qobj | Data = None,
raw_data: bool = False,
)

State of the evolution to be used in a time-dependent operator.

When used as an args:

QobjEvo([op, func], args={"state": SMESolver.StateFeedback()})

The func will receive the density matrix as state during the evolution.

Note

Not supported by the rouchon mehtod.

Parameters:
defaultQobj or qutip.core.data.Data, defaultNone

Initial value to be used at setup of the system.

raw_databool, defaultFalse

If True, the raw matrix will be passed instead of a Qobj. For density matrices, the matrices can be column stacked or square depending on the integration method.

classmethod WienerFeedback(
default: Callable[[float], ndarray[Any, dtype[float]]] = None,
)

Wiener function of the trajectory argument for time dependent systems.

When used as an args:

QobjEvo([op, func], args={"W": SMESolver.WienerFeedback()})

The func will receive a function as W that return an array of wiener processes values at t. The wiener process for the i-th sc_ops is the i-th element for homodyne detection and the (2i, 2i+1) pairs of process in heterodyne detection. The process is a step function with step of length options["dt"].

Note

WienerFeedback can’t be added to a running solver when updating arguments between steps: solver.step(..., args={}).

Parameters:
defaultcallable, optional

Default function used outside the solver. When not passed, a function returning np.array([0]) is used.

property options: dict[str, Any]

Options for stochastic solver:

store_final_state: bool, default: False

Whether or not to store the final state of the evolution in the result class.

store_states: None, bool, default: None

Whether or not to store the state vectors or density matrices. On None the states will be saved if no expectation operators are given.

store_measurement: str, {‘start’, ‘middle’, ‘end’, ‘’}, default: “”

Whether and how to store the measurement for each trajectories. ‘start’, ‘middle’, ‘end’ indicate when in the interval the expectation value of the m_ops is taken. Storing measurements will also store the wiener process, or brownian noise for each trajectories.

progress_bar: str {‘text’, ‘enhanced’, ‘tqdm’, ‘’}, default: “text”

How to present the solver progress. ‘tqdm’ uses the python module of the same name and raise an error if not installed. Empty string or False will disable the bar.

progress_kwargs: dict, default: {“chunk_size”:10}

Arguments to pass to the progress_bar. Qutip’s bars use chunk_size.

keep_runs_results: bool, default: False

Whether to store results from all trajectories or just store the averages.

normalize_output: bool

Normalize output state to hide ODE numerical errors.

method: str, default: “platen”

Which differential equation integration method to use.

map: str {“serial”, “parallel”, “loky”, “mpi”}, default: “serial”

How to run the trajectories. “parallel” uses the multiprocessing module to run in parallel while “loky” and “mpi” use the “loky” and “mpi4py” modules to do so.

mpi_options: dict, default: {}

Only applies if map is “mpi”. This dictionary will be passed as keyword arguments to the mpi4py.futures.MPIPoolExecutor constructor. Note that the max_workers argument is provided separately through the num_cpus option.

num_cpus: None, int, default: None

Number of cpus to use when running in parallel. None detect the number of available cpus.

bitgenerator: {None, “MT19937”, “PCG64DXSM”, …}, default: None

Which of numpy.random’s bitgenerator to use. With None, your numpy version’s default is used.

run(
state: Qobj,
tlist: ArrayLike,
ntraj: int = 1,
*,
args: dict[str, Any] = None,
e_ops: dict[Any, Qobj | QobjEvo | Callable[[float, Qobj], Any]] = None,
target_tol: float | tuple[float, float] | list[tuple[float, float]] = None,
timeout: float = None,
seeds: int | SeedSequence | list[int | SeedSequence] = None,
) MultiTrajResult

Do the evolution of the Quantum system.

For a state at time tlist[0] do the evolution as directed by rhs and for each time in tlist store the state and/or expectation values in a MultiTrajResult. The evolution method and stored results are determined by options.

Parameters:
stateQobj

Initial state of the evolution.

tlistlist of double

Time for which to save the results (state and/or expect) of the evolution. The first element of the list is the initial time of the evolution. Time in the list must be in increasing order, but does not need to be uniformly distributed.

ntrajint

Number of trajectories to add.

argsdict, optional

Change the args of the rhs for the evolution.

e_opslist

list of Qobj or QobjEvo to compute the expectation values. Alternatively, function[s] with the signature f(t, state) -> expect can be used.

timeoutfloat, optional

Maximum time in seconds for the trajectories to run. Once this time is reached, the simulation will end even if the number of trajectories is less than ntraj. The map function, set in options, can interupt the running trajectory or wait for it to finish. Set to an arbitrary high number to disable.

target_tol{float, tuple, list}, optional

Target tolerance of the evolution. The evolution will compute trajectories until the error on the expectation values is lower than this tolerance. The maximum number of trajectories employed is given by ntraj. The error is computed using jackknife resampling. target_tol can be an absolute tolerance or a pair of absolute and relative tolerance, in that order. Lastly, it can be a list of pairs of (atol, rtol) for each e_ops.

seeds{int, SeedSequence, list}, optional

Seed or list of seeds for each trajectories.

Returns:
resultsMultiTrajResult

Results of the evolution. States and/or expect will be saved. You can control the saved data in the options.

run_from_experiment(
state: Qobj,
tlist: ArrayLike,
noise: Sequence[float],
*,
args: dict[str, Any] = None,
e_ops: EopsLike | list[EopsLike] | dict[Any, EopsLike] = None,
measurement: bool = False,
)

Run a single trajectory from a given state and noise.

Parameters:
stateQobj

Initial state of the system.

tlistarray_like

List of times for which to evaluate the state. The tlist must increase uniformly.

noisearray_like

Noise for each time step and each stochastic collapse operators. For homodyne detection, noise[i, t_idx] is the Wiener increments between tlist[t_idx] and tlist[t_idx+1] for the i-th sc_ops. For heterodyne detection, an extra dimension is added for the pair of measurement: noise[i, j, t_idx]``with ``j in {0,1}.

argsdict, optional

Arguments to pass to the Hamiltonian and collapse operators.

e_opsQobj, callable, list or dict, optional

Single operator, or list or dict of operators, for which to evaluate expectation values. Operator can be Qobj, QobjEvo or callables with the signature f(t: float, state: Qobj) -> Any.

measurementbool, defaultFalse

Whether the passed noise is the Wiener increments dW (gaussian noise with standard derivation of dt**0.5), or the measurement.

Homodyne measurement is:

noise[i][t] = dW/dt + expect(sc_ops[i] + sc_ops[i].dag, state[t])

Heterodyne measurement is:

noise[i][0][t] = dW/dt * 2**0.5
  + expect(sc_ops[i] + sc_ops[i].dag, state[t])

noise[i][1][t] = dW/dt * 2**0.5
  -1j * expect(sc_ops[i] - sc_ops[i].dag, state[t])

Note that this function expects the expectation values to be taken at the start of the time step, corresponding to the “start” setting for the “store_measurements” option.

Only available for limited integration methods.

Returns:
resultStochasticTrajResult

Result of the trajectory.

Notes

Only default values of m_ops and dW_factors are supported.

start(
state0: Qobj,
t0: float,
seed: int | SeedSequence = None,
)

Set the initial state and time for a step evolution.

Parameters:
stateQobj

Initial state of the evolution.

t0double

Initial time of the evolution.

seedint, SeedSequence, list, optional

Seed for the random number generator. It can be a single seed used to spawn seeds for each trajectory or a list of seed, one for each trajectory.

Notes

When using step evolution, only one trajectory can be computed at once.

step(t, *, args=None, copy=True, wiener_increment=False)

Evolve the state to t and return the state as a Qobj.

Parameters:
tdouble

Time to evolve to, must be higher than the last call.

argsdict, optional

Update the args of the system. The change is effective from the beginning of the interval. Changing args can slow the evolution.

copybool, default: True

Whether to return a copy of the data or the data in the ODE solver.

wiener_increment: bool, default: False

Whether to return dW in addition to the state.

property sys_dims

Dimensions of the space that the system use:

qutip.basis(sovler.dims) will create a state with proper dimensions for this solver.

class SSESolver(
H: Qobj | QobjEvo,
sc_ops: Sequence[Qobj | QobjEvo],
heterodyne: bool,
*,
c_ops: Sequence[Qobj | QobjEvo] = (),
options: dict[str, Any] = None,
)[source]

Stochastic Schrodinger Equation Solver.

Parameters:
HQobj, QobjEvo, QobjEvo compatible format.

System Hamiltonian as a Qobj or QobjEvo for time-dependent Hamiltonians. List of [Qobj, Coefficient] or callable that can be made into QobjEvo are also accepted.

c_opslist of (QobjEvo, QobjEvo compatible format)

Deterministic collapse operator which will contribute with a standard Lindblad type of dissipation.

sc_opslist of (QobjEvo, QobjEvo compatible format)

List of stochastic collapse operators.

heterodynebool, default: False

Whether to use heterodyne or homodyne detection.

optionsdict, optional

Options for the solver, see SSESolver.options and SIntegrator for a list of all options.

classmethod ExpectFeedback(
operator: Qobj | QobjEvo,
default: Any = 0.0,
)

Expectation value of the instantaneous state of the evolution to be used by a time-dependent operator.

When used as an args:

QobjEvo([op, func], args={"E0": Solver.ExpectFeedback(oper)})

The func will receive expect(oper, state) as E0 during the evolution.

Parameters:
operatorQobj, QobjEvo

Operator to compute the expectation values of.

defaultfloat, default0.

Initial value to be used at setup.

classmethod StateFeedback(
default: Qobj | Data = None,
raw_data: bool = False,
)

State of the evolution to be used in a time-dependent operator.

When used as an args:

QobjEvo([op, func], args={"state": SMESolver.StateFeedback()})

The func will receive the density matrix as state during the evolution.

Note

Not supported by the rouchon mehtod.

Parameters:
defaultQobj or qutip.core.data.Data, defaultNone

Initial value to be used at setup of the system.

raw_databool, defaultFalse

If True, the raw matrix will be passed instead of a Qobj. For density matrices, the matrices can be column stacked or square depending on the integration method.

classmethod WienerFeedback(
default: Callable[[float], ndarray[Any, dtype[float]]] = None,
)

Wiener function of the trajectory argument for time dependent systems.

When used as an args:

QobjEvo([op, func], args={"W": SMESolver.WienerFeedback()})

The func will receive a function as W that return an array of wiener processes values at t. The wiener process for the i-th sc_ops is the i-th element for homodyne detection and the (2i, 2i+1) pairs of process in heterodyne detection. The process is a step function with step of length options["dt"].

Note

WienerFeedback can’t be added to a running solver when updating arguments between steps: solver.step(..., args={}).

Parameters:
defaultcallable, optional

Default function used outside the solver. When not passed, a function returning np.array([0]) is used.

property options: dict[str, Any]

Options for stochastic solver:

store_final_state: bool, default: False

Whether or not to store the final state of the evolution in the result class.

store_states: None, bool, default: None

Whether or not to store the state vectors or density matrices. On None the states will be saved if no expectation operators are given.

store_measurement: str, {‘start’, ‘middle’, ‘end’, ‘’}, default: “”

Whether and how to store the measurement for each trajectories. ‘start’, ‘middle’, ‘end’ indicate when in the interval the expectation value of the m_ops is taken. Storing measurements will also store the wiener process, or brownian noise for each trajectories.

progress_bar: str {‘text’, ‘enhanced’, ‘tqdm’, ‘’}, default: “text”

How to present the solver progress. ‘tqdm’ uses the python module of the same name and raise an error if not installed. Empty string or False will disable the bar.

progress_kwargs: dict, default: {“chunk_size”:10}

Arguments to pass to the progress_bar. Qutip’s bars use chunk_size.

keep_runs_results: bool, default: False

Whether to store results from all trajectories or just store the averages.

normalize_output: bool

Normalize output state to hide ODE numerical errors.

method: str, default: “platen”

Which differential equation integration method to use.

map: str {“serial”, “parallel”, “loky”, “mpi”}, default: “serial”

How to run the trajectories. “parallel” uses the multiprocessing module to run in parallel while “loky” and “mpi” use the “loky” and “mpi4py” modules to do so.

mpi_options: dict, default: {}

Only applies if map is “mpi”. This dictionary will be passed as keyword arguments to the mpi4py.futures.MPIPoolExecutor constructor. Note that the max_workers argument is provided separately through the num_cpus option.

num_cpus: None, int, default: None

Number of cpus to use when running in parallel. None detect the number of available cpus.

bitgenerator: {None, “MT19937”, “PCG64DXSM”, …}, default: None

Which of numpy.random’s bitgenerator to use. With None, your numpy version’s default is used.

run(
state: Qobj,
tlist: ArrayLike,
ntraj: int = 1,
*,
args: dict[str, Any] = None,
e_ops: dict[Any, Qobj | QobjEvo | Callable[[float, Qobj], Any]] = None,
target_tol: float | tuple[float, float] | list[tuple[float, float]] = None,
timeout: float = None,
seeds: int | SeedSequence | list[int | SeedSequence] = None,
) MultiTrajResult

Do the evolution of the Quantum system.

For a state at time tlist[0] do the evolution as directed by rhs and for each time in tlist store the state and/or expectation values in a MultiTrajResult. The evolution method and stored results are determined by options.

Parameters:
stateQobj

Initial state of the evolution.

tlistlist of double

Time for which to save the results (state and/or expect) of the evolution. The first element of the list is the initial time of the evolution. Time in the list must be in increasing order, but does not need to be uniformly distributed.

ntrajint

Number of trajectories to add.

argsdict, optional

Change the args of the rhs for the evolution.

e_opslist

list of Qobj or QobjEvo to compute the expectation values. Alternatively, function[s] with the signature f(t, state) -> expect can be used.

timeoutfloat, optional

Maximum time in seconds for the trajectories to run. Once this time is reached, the simulation will end even if the number of trajectories is less than ntraj. The map function, set in options, can interupt the running trajectory or wait for it to finish. Set to an arbitrary high number to disable.

target_tol{float, tuple, list}, optional

Target tolerance of the evolution. The evolution will compute trajectories until the error on the expectation values is lower than this tolerance. The maximum number of trajectories employed is given by ntraj. The error is computed using jackknife resampling. target_tol can be an absolute tolerance or a pair of absolute and relative tolerance, in that order. Lastly, it can be a list of pairs of (atol, rtol) for each e_ops.

seeds{int, SeedSequence, list}, optional

Seed or list of seeds for each trajectories.

Returns:
resultsMultiTrajResult

Results of the evolution. States and/or expect will be saved. You can control the saved data in the options.

run_from_experiment(
state: Qobj,
tlist: ArrayLike,
noise: Sequence[float],
*,
args: dict[str, Any] = None,
e_ops: EopsLike | list[EopsLike] | dict[Any, EopsLike] = None,
measurement: bool = False,
)

Run a single trajectory from a given state and noise.

Parameters:
stateQobj

Initial state of the system.

tlistarray_like

List of times for which to evaluate the state. The tlist must increase uniformly.

noisearray_like

Noise for each time step and each stochastic collapse operators. For homodyne detection, noise[i, t_idx] is the Wiener increments between tlist[t_idx] and tlist[t_idx+1] for the i-th sc_ops. For heterodyne detection, an extra dimension is added for the pair of measurement: noise[i, j, t_idx]``with ``j in {0,1}.

argsdict, optional

Arguments to pass to the Hamiltonian and collapse operators.

e_opsQobj, callable, list or dict, optional

Single operator, or list or dict of operators, for which to evaluate expectation values. Operator can be Qobj, QobjEvo or callables with the signature f(t: float, state: Qobj) -> Any.

measurementbool, defaultFalse

Whether the passed noise is the Wiener increments dW (gaussian noise with standard derivation of dt**0.5), or the measurement.

Homodyne measurement is:

noise[i][t] = dW/dt + expect(sc_ops[i] + sc_ops[i].dag, state[t])

Heterodyne measurement is:

noise[i][0][t] = dW/dt * 2**0.5
  + expect(sc_ops[i] + sc_ops[i].dag, state[t])

noise[i][1][t] = dW/dt * 2**0.5
  -1j * expect(sc_ops[i] - sc_ops[i].dag, state[t])

Note that this function expects the expectation values to be taken at the start of the time step, corresponding to the “start” setting for the “store_measurements” option.

Only available for limited integration methods.

Returns:
resultStochasticTrajResult

Result of the trajectory.

Notes

Only default values of m_ops and dW_factors are supported.

start(
state0: Qobj,
t0: float,
seed: int | SeedSequence = None,
)

Set the initial state and time for a step evolution.

Parameters:
stateQobj

Initial state of the evolution.

t0double

Initial time of the evolution.

seedint, SeedSequence, list, optional

Seed for the random number generator. It can be a single seed used to spawn seeds for each trajectory or a list of seed, one for each trajectory.

Notes

When using step evolution, only one trajectory can be computed at once.

step(t, *, args=None, copy=True, wiener_increment=False)

Evolve the state to t and return the state as a Qobj.

Parameters:
tdouble

Time to evolve to, must be higher than the last call.

argsdict, optional

Update the args of the system. The change is effective from the beginning of the interval. Changing args can slow the evolution.

copybool, default: True

Whether to return a copy of the data or the data in the ODE solver.

wiener_increment: bool, default: False

Whether to return dW in addition to the state.

property sys_dims

Dimensions of the space that the system use:

qutip.basis(sovler.dims) will create a state with proper dimensions for this solver.

class MultiTrajResult(
e_ops,
options: MultiTrajResultOptions,
*,
solver=None,
stats=None,
**kw,
)[source]

Base class for storing results for solver using multiple trajectories.

Parameters:
e_opsQobj, QobjEvo, function or list or dict of these

The e_ops parameter defines the set of values to record at each time step t. If an element is a Qobj or QobjEvo the value recorded is the expectation value of that operator given the state at t. If the element is a function, f, the value recorded is f(t, state).

The values are recorded in the .expect attribute of this result object. .expect is a list, where each item contains the values of the corresponding e_op.

Function e_ops must return a number so the average can be computed.

optionsdict

The options for this result class.

solverstr or None

The name of the solver generating these results.

statsdict or None

The stats generated by the solver while producing these results. Note that the solver may update the stats directly while producing results.

kwdict

Additional parameters specific to a result sub-class.

Attributes:
timeslist

A list of the times at which the expectation values and states were recorded.

average_stateslist of Qobj

States averages as density matrices.

runs_stateslist of list of Qobj

States of every runs as states[run][t].

average_final_stateQobj:

Last states of each trajectories averaged into a density matrix.

runs_final_statelist of Qobj

The final state for each trajectory (if the recording of the final state and trajectories was requested).

average_expectlist of array of expectation values

A list containing the values of each e_op averaged over each trajectories. The list is in the same order in which the e_ops were supplied and empty if no e_ops were given.

Each element is itself an array and contains the values of the corresponding e_op, with one value for each time in .times.

std_expectlist of array of expectation values

A list containing the standard derivation of each e_op over each trajectories. The list is in the same order in which the e_ops were supplied and empty if no e_ops were given.

Each element is itself an array and contains the values of the corresponding e_op, with one value for each time in .times.

runs_expectlist of array of expectation values

A list containing the values of each e_op for each trajectories. The list is in the same order in which the e_ops were supplied and empty if no e_ops were given. Only available if the storing of trajectories was requested.

The order of the elements is runs_expect[e_ops][trajectory][time].

Each element is itself an array and contains the values of the corresponding e_op, with one value for each time in .times.

average_e_datadict

A dictionary containing the values of each e_op averaged over each trajectories. If the e_ops were supplied as a dictionary, the keys are the same as in that dictionary. Otherwise the keys are the index of the e_op in the .expect list.

The lists of expectation values returned are the same lists as those returned by .expect.

std_e_datadict

A dictionary containing the standard derivation of each e_op over each trajectories. If the e_ops were supplied as a dictionary, the keys are the same as in that dictionary. Otherwise the keys are the index of the e_op in the .expect list.

The lists of expectation values returned are the same lists as those returned by .expect.

runs_e_datadict

A dictionary containing the values of each e_op for each trajectories. If the e_ops were supplied as a dictionary, the keys are the same as in that dictionary. Otherwise the keys are the index of the e_op in the .expect list. Only available if the storing of trajectories was requested.

The order of the elements is runs_expect[e_ops][trajectory][time].

The lists of expectation values returned are the same lists as those returned by .expect.

runs_weightslist

For each trajectory, the weight with which that trajectory enters averages.

deterministic_weightslist

For each deterministic trajectory, (when using improved_sampling) the weight with which that trajectory enters averages.

solverstr or None

The name of the solver generating these results.

num_trajectories: int

Number of trajectories computed.

seeds: list of SeedSequence

The seeds used to compute each trajectories.

trajectories: list of Result

If the option keep_runs_results is set, a list of all trajectories.

deterministic_trajectories: list of Result

A list of the no-jump trajectories if the option improved_sampling is set.

statsdict or None

The stats generated by the solver while producing these results.

optionsSolverResultsOptions

The options for this result class.

add_deterministic(trajectory, weight)[source]

Add a trajectory that was not randomly generated. The weight provided here is the exact weight that will be used for this trajectory in all averages.

Parameters:
trajectoryResult

The result of the simulation of the deterministic trajectory

weightfloat

Number (usually between 0 and 1), exact weight of this trajectory

property average_final_state

Last states of each trajectories averaged into a density matrix.

property average_states

States averages as density matrices.

property final_state

Runs final states if available, average otherwise.

merge(other, p=None)[source]

Merges two multi-trajectory results.

If this result represent an ensemble \(\rho\), and other represents an ensemble \(\rho'\), then the merged result represents the ensemble

\[\rho_{\mathrm{merge}} = p \rho + (1 - p) \rho'\]

where p is a parameter between 0 and 1. Its default value is \(p_{\textrm{def}} = N / (N + N')\), N and N’ being the number of trajectories in the two result objects.

Parameters:
otherMultiTrajResult

The multi-trajectory result to merge with this one

pfloat [optional]

The relative weight of this result in the combination. By default, will be chosen such that all trajectories contribute equally to the merged result.

property runs_final_states

Last states of each trajectories.

property runs_states

States of every runs as states[run][t].

property states

Runs final states if available, average otherwise.

steady_state(N=0)[source]

Average the states of the last N times of every runs as a density matrix. Should converge to the steady state in the right circumstances.

Parameters:
Nint [optional]

Number of states from the end of tlist to average. Per default all states will be averaged.

Non-Markovian Solvers

This module contains an implementation of the non-Markovian transfer tensor method (TTM), introduced in [1].

[1] Javier Cerrillo and Jianshu Cao, Phys. Rev. Lett 112, 110401 (2014)

ttmsolve(dynmaps, state0, times, e_ops=(), num_learning=0, options=None)[source]

Expand time-evolution using the Transfer Tensor Method [1], based on a set of precomputed dynamical maps.

Parameters:
dynmapslist of Qobj, callable

List of precomputed dynamical maps (superoperators) for the first times of times or a callback function that returns the superoperator at a given time.

state0Qobj

Initial density matrix or state vector (ket).

timesarray_like

List of times \(t_n\) at which to compute results. Must be uniformily spaced.

e_opsQobj, callable, or list, optional

Single operator or list of operators for which to evaluate expectation values or callable or list of callable. Callable signature must be, f(t: float, state: Qobj). See expect for more detail of operator expectation.

num_learningint, default: 0

Number of times used to construct the dynmaps operators when dynmaps is a callable.

optionsdictionary, optional

Dictionary of options for the solver.

  • store_final_state : bool Whether or not to store the final state of the evolution in the result class.

  • store_states : bool, None Whether or not to store the state vectors or density matrices. On None the states will be saved if no expectation operators are given.

  • normalize_output : bool Normalize output state to hide ODE numerical errors.

  • threshold : float Threshold for halting. Halts if \(||T_{n}-T_{n-1}||\) is below treshold.

Returns:
output: Result

An instance of the class Result.

[1]

Javier Cerrillo and Jianshu Cao, Phys. Rev. Lett 112, 110401 (2014) ..

Integrator

Different ODE solver from many sources (scipy, diffrax, home made, etc.) used by qutip solvers. Their options are added to the solver options:

class IntegratorScipyAdams(system, options)[source]

Integrator using Scipy ode with zvode integrator using adams method. Ordinary Differential Equation solver by netlib (https://www.netlib.org/odepack).

Usable with method="adams"

property options

Supported options by zvode integrator:

atolfloat, default: 1e-8

Absolute tolerance.

rtolfloat, default: 1e-6

Relative tolerance.

orderint, default: 12, ‘adams’ or 5, ‘bdf’

Order of integrator <=12 ‘adams’, <=5 ‘bdf’

nstepsint, default: 2500

Max. number of internal steps/call.

first_stepfloat, default: 0

Size of initial step (0 = automatic).

min_stepfloat, default: 0

Minimum step size (0 = automatic).

max_stepfloat, default: 0

Maximum step size (0 = automatic) When using pulses, change to half the thinest pulse otherwise it may be skipped.

class IntegratorScipyBDF(system, options)[source]

Integrator using Scipy ode with zvode integrator using bdf method. Ordinary Differential Equation solver by netlib (https://www.netlib.org/odepack).

Usable with method="bdf"

property options

Supported options by zvode integrator:

atolfloat, default: 1e-8

Absolute tolerance.

rtolfloat, default: 1e-6

Relative tolerance.

orderint, default: 12, ‘adams’ or 5, ‘bdf’

Order of integrator <=12 ‘adams’, <=5 ‘bdf’

nstepsint, default: 2500

Max. number of internal steps/call.

first_stepfloat, default: 0

Size of initial step (0 = automatic).

min_stepfloat, default: 0

Minimum step size (0 = automatic).

max_stepfloat, default: 0

Maximum step size (0 = automatic) When using pulses, change to half the thinest pulse otherwise it may be skipped.

class IntegratorScipylsoda(system, options)[source]

Integrator using Scipy ode with lsoda integrator. ODE solver by netlib (https://www.netlib.org/odepack) Automatically choose between ‘Adams’ and ‘BDF’ methods to solve both stiff and non-stiff systems.

Usable with method="lsoda"

property options

Supported options by lsoda integrator:

atolfloat, default: 1e-8

Absolute tolerance.

rtolfloat, default: 1e-6

Relative tolerance.

nstepsint, default: 2500

Max. number of internal steps/call.

max_order_nsint, default: 12

Maximum order used in the nonstiff case (<= 12).

max_order_sint, default: 5

Maximum order used in the stiff case (<= 5).

first_stepfloat, default: 0

Size of initial step (0 = automatic).

max_stepfloat, default: 0

Maximum step size (0 = automatic) When using pulses, change to half the thinest pulse otherwise it may be skipped.

min_stepfloat, default: 0

Minimum step size (0 = automatic)

class IntegratorScipyDop853(system, options)[source]

Integrator using Scipy ode with dop853 integrator. Eight order runge-kutta method by Dormand & Prince. Use fortran implementation from [E. Hairer, S.P. Norsett and G. Wanner, Solving Ordinary Differential Equations i. Nonstiff Problems. 2nd edition. Springer Series in Computational Mathematics, Springer-Verlag (1993)].

Usable with method="dop853"

property options

Supported options by dop853 integrator:

atolfloat, default: 1e-8

Absolute tolerance.

rtolfloat, default: 1e-6

Relative tolerance.

nstepsint, default: 2500

Max. number of internal steps/call.

first_stepfloat, default: 0

Size of initial step (0 = automatic).

max_stepfloat, default: 0

Maximum step size (0 = automatic)

ifactor, dfactorfloat, default: 6., 0.3

Maximum factor to increase/decrease step size by in one step

betafloat, default: 0

Beta parameter for stabilised step size control.

See scipy.integrate.ode ode for more detail

class IntegratorVern7(system, options)[source]

QuTiP’s implementation of Verner’s “most efficient” Runge-Kutta method of order 7. These are Runge-Kutta methods with variable steps and dense output.

The implementation uses QuTiP’s Data objects for the state, allowing sparse, GPU or other data layer objects to be used efficiently by the solver in their native formats.

See https://www.sfu.ca/~jverner/ for a detailed description of the methods.

Usable with method="vern7"

property options

Supported options by verner method:

atolfloat, default: 1e-8

Absolute tolerance.

rtolfloat, default: 1e-6

Relative tolerance.

nstepsint, default: 1000

Max. number of internal steps/call.

first_stepfloat, default: 0

Size of initial step (0 = automatic).

min_stepfloat, default: 0

Minimum step size (0 = automatic).

max_stepfloat, default: 0

Maximum step size (0 = automatic) When using pulses, change to half the thinest pulse otherwise it may be skipped.

interpolatebool, default: True

Whether to use interpolation step, faster most of the time.

class IntegratorVern9(system, options)[source]

QuTiP’s implementation of Verner’s “most efficient” Runge-Kutta method of order 9. These are Runge-Kutta methods with variable steps and dense output.

The implementation uses QuTiP’s Data objects for the state, allowing sparse, GPU or other data layer objects to be used efficiently by the solver in their native formats.

See https://www.sfu.ca/~jverner/ for a detailed description of the methods.

Usable with method="vern9"

property options

Supported options by verner method:

atolfloat, default: 1e-8

Absolute tolerance.

rtolfloat, default: 1e-6

Relative tolerance.

nstepsint, default: 1000

Max. number of internal steps/call.

first_stepfloat, default: 0

Size of initial step (0 = automatic).

min_stepfloat, default: 0

Minimum step size (0 = automatic).

max_stepfloat, default: 0

Maximum step size (0 = automatic) When using pulses, change to half the thinest pulse otherwise it may be skipped.

interpolatebool, default: True

Whether to use interpolation step, faster most of the time.

class IntegratorDiag(system, options)[source]

Integrator solving the ODE by diagonalizing the system and solving analytically. It can only solve constant system and has a long preparation time, but the integration is fast.

Usable with method="diag"

property options

Supported options by “diag” method:

eigensolver_dtypestr, default: “dense”

Qutip data type {“dense”, “csr”, etc.} to use when computing the eigenstates. The dense eigen solver is usually faster and more stable.

class IntegratorKrylov(system, options)[source]

Evolve the state vector (“psi0”) finding an approximation for the time evolution operator of Hamiltonian (“H”) by obtaining the projection of the time evolution operator on a set of small dimensional Krylov subspaces (m << dim(H)).

property options

Supported options by krylov method:

atolfloat, default: 1e-7

Absolute tolerance.

nstepsint, default: 100

Max. number of internal steps/call.

min_step, max_stepfloat, default: (1e-5, 1e5)

Minimum and maximum step size.

krylov_dim: int, default: 0

Dimension of Krylov approximation subspaces used for the time evolution approximation. If the defaut 0 is given, the dimension is calculated from the system size N, using min(int((N + 100)**0.5), N-1).

sub_system_tol: float, default: 1e-7

Tolerance to detect a happy breakdown. A happy breakdown occurs when the initial ket is in a subspace of the Hamiltonian smaller than krylov_dim.

always_compute_step: bool, default: False

If True, the step length is computed each time a new Krylov subspace is computed. Otherwise it is computed only once when creating the integrator.

Stochastic Integrator

class RouchonSODE(rhs, options)[source]

Stochastic integration method keeping the positivity of the density matrix. See eq. (4) Pierre Rouchon and Jason F. Ralpha, Efficient Quantum Filtering for Quantum Feedback Control, arXiv:1410.5345 [quant-ph], Phys. Rev. A 91, 012118, (2015).

  • Order: strong 1

Notes

This method should be used with very small dt. Unlike other methods that will return unphysical state (negative eigenvalues, Nans) when the time step is too large, this method will return state that seems normal.

property options

Supported options by Rouchon Stochastic Integrators:

dtfloat, default: 0.001

Internal time step.

tolfloat, default: 1e-7

Relative tolerance.

class EulerSODE(rhs, options)[source]

A simple generalization of the Euler method for ordinary differential equations to stochastic differential equations. Only solver which could take non-commuting sc_ops.

  • Order: 0.5

property options

Supported options by Explicit Stochastic Integrators:

dtfloat, default: 0.001

Internal time step.

tolfloat, default: 1e-10

Tolerance for the time steps.

class Milstein_SODE(rhs, options)[source]

An order 1.0 strong Taylor scheme. Better approximate numerical solution to stochastic differential equations. See eq. (3.12) of chapter 10.3 of Peter E. Kloeden and Exkhard Platen, Numerical Solution of Stochastic Differential Equations..

  • Order strong 1.0

property options

Supported options by Explicit Stochastic Integrators:

dtfloat, default: 0.001

Internal time step.

tolfloat, default: 1e-10

Tolerance for the time steps.

class Taylor1_5_SODE(rhs, options)[source]

Order 1.5 strong Taylor scheme. Solver with more terms of the Ito-Taylor expansion. See eq. (4.6) of chapter 10.4 of Peter E. Kloeden and Exkhard Platen, Numerical Solution of Stochastic Differential Equations.

  • Order strong 1.5

property options

Supported options by Order 1.5 strong Taylor Stochastic Integrators:

dtfloat, default: 0.001

Internal time step.

tolfloat, default: 1e-10

Relative tolerance.

derr_dtfloat, default: 1e-6

Finite time difference used to compute the derrivative of the hamiltonian and sc_ops.

class Implicit_Milstein_SODE(rhs, options)[source]

An order 1.0 implicit strong Taylor scheme. Implicit Milstein scheme for the numerical simulation of stiff stochastic differential equations. Eq. (2.11) with alpha=0.5 of chapter 12.2 of Peter E. Kloeden and Exkhard Platen, Numerical Solution of Stochastic Differential Equations.

  • Order strong 1.0

property options

Supported options by Implicit Stochastic Integrators:

dtfloat, default: 0.001

Internal time step.

tolfloat, default: 1e-10

Tolerance for the time steps.

solve_methodstr, default: None

Method used for solver the Ax=b of the implicit step. Accept methods supported by qutip.core.data.solve. When the system is constant, the inverse of the matrix A can be used by entering inv.

solve_optionsdict, default: {}

Options to pass to the call to qutip.core.data.solve.

class Implicit_Taylor1_5_SODE(rhs, options)[source]

Order 1.5 implicit strong Taylor scheme. Solver with more terms of the Ito-Taylor expansion. Eq. (2.18) with alpha=0.5 of chapter 12.2 of Peter E. Kloeden and Exkhard Platen, Numerical Solution of Stochastic Differential Equations.

  • Order strong 1.5

property options

Supported options by Implicit Order 1.5 strong Taylor Stochastic Integrators:

dtfloat, default: 0.001

Internal time step.

tolfloat, default: 1e-10

Tolerance for the time steps.

solve_methodstr, default: None

Method used for solver the Ax=b of the implicit step. Accept methods supported by qutip.core.data.solve. When the system is constant, the inverse of the matrix A can be used by entering inv.

solve_optionsdict, default: {}

Options to pass to the call to qutip.core.data.solve.

derr_dtfloat, default: 1e-6

Finite time difference used to compute the derrivative of the hamiltonian and sc_ops.

class PlatenSODE(rhs, options)[source]

Explicit scheme, creates the Milstein using finite differences instead of analytic derivatives. Also contains some higher order terms, thus converges better than Milstein while staying strong order 1.0. Does not require derivatives. See eq. (7.47) of chapter 7 of H.-P. Breuer and F. Petruccione, The Theory of Open Quantum Systems.

  • Order: strong 1, weak 2

property options

Supported options by Explicit Stochastic Integrators:

dtfloat, default: 0.001

Internal time step.

tolfloat, default: 1e-10

Tolerance for the time steps.

class Explicit1_5_SODE(rhs, options)[source]

Explicit order 1.5 strong schemes. Reproduce the order 1.5 strong Taylor scheme using finite difference instead of derivatives. Slower than taylor15 but usable when derrivatives cannot be analytically obtained. See eq. (2.13) of chapter 11.2 of Peter E. Kloeden and Exkhard Platen, Numerical Solution of Stochastic Differential Equations.

  • Order: strong 1.5

property options

Supported options by Explicit Stochastic Integrators:

dtfloat, default: 0.001

Internal time step.

tolfloat, default: 1e-10

Tolerance for the time steps.

class PredCorr_SODE(rhs, options)[source]

Generalization of the trapezoidal method to stochastic differential equations. More stable than explicit methods. See eq. (5.4) of chapter 15.5 of Peter E. Kloeden and Exkhard Platen, Numerical Solution of Stochastic Differential Equations.

  • Order strong 0.5, weak 1.0

  • Codes to only correct the stochastic part (\(\alpha=0\), \(\eta=1/2\)): 'pred-corr', 'predictor-corrector' or 'pc-euler'

  • Codes to correct both the stochastic and deterministic parts (\(\alpha=1/2\), \(\eta=1/2\)): 'pc-euler-imp', 'pc-euler-2' or 'pred-corr-2'

property options

Supported options by Explicit Stochastic Integrators:

dtfloat, default: 0.001

Internal time step.

tolfloat, default: 1e-10

Tolerance for the time steps.

alphafloat, default: 0.

Implicit factor to the drift. eff_drift ~= drift(t) * (1-alpha) + drift(t+dt) * alpha

etafloat, default: 0.5

Implicit factor to the diffusion. eff_diffusion ~= diffusion(t) * (1-eta) + diffusion(t+dt) * eta

Parallelization

This module provides functions for parallel execution of loops and function mappings, using the builtin Python module multiprocessing or the loky parallel execution library.

loky_pmap(
task,
values,
task_args=None,
task_kwargs=None,
reduce_func=None,
map_kw=None,
progress_bar=None,
progress_bar_kwargs={},
)[source]

Parallel execution of a mapping of values to the function task. This is functionally equivalent to:

result = [task(value, *task_args, **task_kwargs) for value in values]

Use the loky module instead of multiprocessing.

Parameters:
taska Python function

The function that is to be called for each value in task_vec.

valuesarray / list

The list or array of values for which the task function is to be evaluated.

task_argslist, optional

The optional additional arguments to the task function.

task_kwargsdictionary, optional

The optional additional keyword arguments to the task function.

reduce_funcfunc, optional

If provided, it will be called with the output of each task instead of storing them in a list. Note that the order in which results are passed to reduce_func is not defined. It should return None or a number. When returning a number, it represents the estimation of the number of tasks left. On a return <= 0, the map will end early.

progress_barstr, optional

Progress bar options’s string for showing progress.

progress_bar_kwargsdict, optional

Options for the progress bar.

map_kw: dict, optional

Dictionary containing entry for: - timeout: float, Maximum time (sec) for the whole map. - num_cpus: int, Number of jobs to run at once. - fail_fast: bool, Abort at the first error.

Returns:
resultlist

The result list contains the value of task(value, *task_args, **task_kwargs) for each value in values. If a reduce_func is provided, and empty list will be returned.

mpi_pmap(
task,
values,
task_args=None,
task_kwargs=None,
reduce_func=None,
map_kw=None,
progress_bar=None,
progress_bar_kwargs={},
)[source]

Parallel execution of a mapping of values to the function task. This is functionally equivalent to:

result = [task(value, *task_args, **task_kwargs) for value in values]

Uses the mpi4py module to execute the tasks asynchronously with MPI processes. For more information, consult the documentation of mpi4py and the mpi4py.MPIPoolExecutor class.

Note: in keeping consistent with the API of parallel_map, the parameter determining the number of requested worker processes is called num_cpus. The value of map_kw[‘num_cpus’] is passed to the MPIPoolExecutor as its max_workers argument. If this parameter is not provided, the environment variable QUTIP_NUM_PROCESSES is used instead. If this environment variable is not set either, QuTiP will use default values that might be unsuitable for MPI applications.

Parameters:
taska Python function

The function that is to be called for each value in task_vec.

valuesarray / list

The list or array of values for which the task function is to be evaluated.

task_argslist, optional

The optional additional arguments to the task function.

task_kwargsdictionary, optional

The optional additional keyword arguments to the task function.

reduce_funcfunc, optional

If provided, it will be called with the output of each task instead of storing them in a list. Note that the order in which results are passed to reduce_func is not defined. It should return None or a number. When returning a number, it represents the estimation of the number of tasks left. On a return <= 0, the map will end early.

progress_barstr, optional

Progress bar options’s string for showing progress.

progress_bar_kwargsdict, optional

Options for the progress bar.

map_kw: dict, optional

Dictionary containing entry for: - timeout: float, Maximum time (sec) for the whole map. - num_cpus: int, Number of jobs to run at once. - fail_fast: bool, Abort at the first error. All remaining entries of map_kw will be passed to the mpi4py.MPIPoolExecutor constructor.

Returns:
resultlist

The result list contains the value of task(value, *task_args, **task_kwargs) for each value in values. If a reduce_func is provided, and empty list will be returned.

parallel_map(
task,
values,
task_args=None,
task_kwargs=None,
reduce_func=None,
map_kw=None,
progress_bar=None,
progress_bar_kwargs={},
)[source]

Parallel execution of a mapping of values to the function task. This is functionally equivalent to:

result = [task(value, *task_args, **task_kwargs) for value in values]
Parameters:
taska Python function

The function that is to be called for each value in task_vec.

valuesarray / list

The list or array of values for which the task function is to be evaluated.

task_argslist, optional

The optional additional arguments to the task function.

task_kwargsdictionary, optional

The optional additional keyword arguments to the task function.

reduce_funcfunc, optional

If provided, it will be called with the output of each task instead of storing them in a list. Note that the order in which results are passed to reduce_func is not defined. It should return None or a number. When returning a number, it represents the estimation of the number of tasks left. On a return <= 0, the map will end early.

progress_barstr, optional

Progress bar options’s string for showing progress.

progress_bar_kwargsdict, optional

Options for the progress bar.

map_kw: dict, optional

Dictionary containing entry for: - timeout: float, Maximum time (sec) for the whole map. - num_cpus: int, Number of jobs to run at once. - fail_fast: bool, Abort at the first error.

Returns:
resultlist

The result list contains the value of task(value, *task_args, **task_kwargs) for each value in values. If a reduce_func is provided, and empty list will be returned.

serial_map(
task,
values,
task_args=None,
task_kwargs=None,
reduce_func=None,
map_kw=None,
progress_bar=None,
progress_bar_kwargs={},
)[source]

Serial mapping function with the same call signature as parallel_map, for easy switching between serial and parallel execution. This is functionally equivalent to:

result = [task(value, *task_args, **task_kwargs) for value in values]

This function work as a drop-in replacement of parallel_map.

Parameters:
taska Python function

The function that is to be called for each value in task_vec.

valuesarray / list

The list or array of values for which the task function is to be evaluated.

task_argslist, optional

The optional additional argument to the task function.

task_kwargsdictionary, optional

The optional additional keyword argument to the task function.

reduce_funcfunc, optional

If provided, it will be called with the output of each tasks instead of storing a them in a list. It should return None or a number. When returning a number, it represent the estimation of the number of task left. On a return <= 0, the map will end early.

progress_barstr, optional

Progress bar options’s string for showing progress.

progress_bar_kwargsdict, optional

Options for the progress bar.

map_kw: dict, optional

Dictionary containing: - timeout: float, Maximum time (sec) for the whole map. - fail_fast: bool, Raise an error at the first.

Returns:
resultlist

The result list contains the value of task(value, *task_args, **task_kwargs) for each value in values. If a reduce_func is provided, and empty list will be returned.

Propagators

propagator(
H: QobjEvoLike,
t: Number,
c_ops: QobjEvoLike | list[QobjEvoLike] = None,
args: dict[str, Any] = None,
options: dict[str, Any] = None,
**kwargs,
) Qobj | list[Qobj][source]

Calculate the propagator U(t) for the density matrix or wave function such that \(\psi(t) = U(t)\psi(0)\) or \(\rho_{\mathrm vec}(t) = U(t) \rho_{\mathrm vec}(0)\) where \(\rho_{\mathrm vec}\) is the vector representation of the density matrix.

Parameters:
HQobj, QobjEvo, QobjEvo compatible format

Possibly time-dependent system Liouvillian or Hamiltonian as a Qobj or QobjEvo. list of [Qobj, Coefficient] or callable that can be made into QobjEvo are also accepted.

tfloat or array-like

Time or list of times for which to evaluate the propagator. If a single time t is passed, the propagator from 0 to t is computed. When t is a list, the propagators from the first time in the list to each elements in t is returned. In that case, the first output will always be the identity matrix.

c_opslist, optional

List of collapse operators as Qobj, QobjEvo or list that can be made into QobjEvo.

argsdictionary, optional

Parameters to callback functions for time-dependent Hamiltonians and collapse operators.

optionsdict, optional

Options for the solver.

**kwargs

Extra parameters to use when creating the QobjEvo from a list format H. The most common are tlist and order for array-based time dependance. See QobjEvo for the full list.

Returns:
UQobj, list

Instance representing the propagator(s) \(U(t)\). Return a single Qobj when t is a number or a list when t is a list.

Notes

Unlike sesolve or mesolve, the output times in t are not used for array time dependent system. tlist must be passed as a keyword argument in those case. tlist and t can have different length and values.

propagator_steadystate(U: Qobj) Qobj[source]

Find the steady state for successive applications of the propagator \(U\).

Parameters:
UQobj

Operator representing the propagator.

Returns:
aQobj

Instance representing the steady-state density matrix.

class Propagator(
system: Qobj | QobjEvo | Solver,
*,
c_ops: QobjEvoLike | list[QobjEvoLike] = None,
args: dict[str, Any] = None,
options: dict[str, Any] = None,
memoize: int = 10,
tol: float = 1e-14,
)[source]

A generator of propagator for a system.

Usage:

U = Propagator(H, c_ops)

psi_t = U(t) @ psi_0

Save some previously computed propagator are stored to speed up subsequent computation. Changing args will erase these stored probagator.

Parameters:
systemQobj, QobjEvo, Solver

Possibly time-dependent system driving the evolution, either already packaged in a solver, such as SESolver or BRSolver, or the Liouvillian or Hamiltonian as a Qobj, QobjEvo. list of [Qobj, Coefficient] or callable that can be made into QobjEvo are also accepted.

Solvers that run non-deterministacilly, such as MCSolver, are not supported.

c_opslist, optional

List of Qobj or QobjEvo collapse operators.

argsdictionary, optional

Parameters to callback functions for time-dependent Hamiltonians and collapse operators.

optionsdict, optional

Options for the solver.

memoizeint, default: 10

Max number of propagator to save.

tolfloat, default: 1e-14

Absolute tolerance for the time. If a previous propagator was computed at a time within tolerance, that propagator will be returned.

Notes

The Propagator is not a QobjEvo so it cannot be used for operations with Qobj or QobjEvo. It can be made into a QobjEvo with

U = QobjEvo(Propagator(H))
__call__(t: float, t_start: float = 0, **args)[source]

Get the propagator from t_start to t.

Parameters:
tfloat

Time at which to compute the propagator.

t_start: float [0]
Time at which the propagator start such that:

psi[t] = U.prop(t, t_start) @ psi[t_start]

argsdict

Argument to pass to a time dependent Hamiltonian. Updating args take effect since t=0 and the new args will be used in future call.

inv(t: float, **args)[source]
Get the inverse of the propagator at t, such that

psi_0 = U.inv(t) @ psi_t

Parameters:
tfloat

Time at which to compute the propagator.

argsdict

Argument to pass to a time dependent Hamiltonian. Updating args take effect since t=0 and the new args will be used in future call.

Other dynamics functions

Correlation Functions

coherence_function_g1(
H,
state0,
taulist,
c_ops,
a_op,
solver='me',
args=None,
options=None,
)[source]

Calculate the normalized first-order quantum coherence function:

\[g^{(1)}(\tau) = \frac{\langle A^\dagger(\tau)A(0)\rangle} {\sqrt{\langle A^\dagger(\tau)A(\tau)\rangle \langle A^\dagger(0)A(0)\rangle}}\]

using the quantum regression theorem and the evolution solver indicated by the solver parameter.

Parameters:
HQobj, QobjEvo

System Hamiltonian, may be time-dependent for solver choice of me.

state0Qobj

Initial state density matrix \(\rho(t_0)\) or state vector \(\psi(t_0)\). If ‘state0’ is ‘None’, then the steady state will be used as the initial state. The ‘steady-state’ is only implemented if c_ops are provided and the Hamiltonian is constant.

taulistarray_like

List of times for \(\tau\). taulist must be positive and contain the element 0.

c_opslist of {Qobj, QobjEvo}

List of collapse operators

a_opQobj, QobjEvo

Operator A.

solverstr {‘me’, ‘es’}, default: ‘me’

Choice of solver, me for master-equation, and es for exponential series. es is equivalent to me with options={"method": "diag"}.

argsdict, optional

dictionary of parameters for time-dependent Hamiltonians

optionsdict, optional

Options for the solver.

Returns:
g1, G1tuple

The normalized and unnormalized second-order coherence function.

coherence_function_g2(
H,
state0,
taulist,
c_ops,
a_op,
solver='me',
args=None,
options=None,
)[source]

Calculate the normalized second-order quantum coherence function:

\[ g^{(2)}(\tau) = \frac{\langle A^\dagger(0)A^\dagger(\tau)A(\tau)A(0)\rangle} {\langle A^\dagger(\tau)A(\tau)\rangle \langle A^\dagger(0)A(0)\rangle}\]

using the quantum regression theorem and the evolution solver indicated by the solver parameter.

Parameters:
HQobj, QobjEvo

System Hamiltonian, may be time-dependent for solver choice of me.

state0Qobj

Initial state density matrix \(\rho(t_0)\) or state vector \(\psi(t_0)\). If ‘state0’ is ‘None’, then the steady state will be used as the initial state. The ‘steady-state’ is only implemented if c_ops are provided and the Hamiltonian is constant.

taulistarray_like

List of times for \(\tau\). taulist must be positive and contain the element 0.

c_opslist

List of collapse operators, may be time-dependent for solver choice of me.

a_opQobj

Operator A.

argsdict, optional

Dictionary of arguments to be passed to solver.

solverstr {‘me’, ‘es’}, default: ‘me’

Choice of solver, me for master-equation, and es for exponential series. es is equivalent to me with options={"method": "diag"}.

optionsdict, optional

Options for the solver.

Returns:
g2, G2tuple

The normalized and unnormalized second-order coherence function.

correlation_2op_1t(
H,
state0,
taulist,
c_ops,
a_op,
b_op,
solver='me',
reverse=False,
args=None,
options=None,
)[source]

Calculate the two-operator one-time correlation function: \(\left<A(\tau)B(0)\right>\) along one time axis using the quantum regression theorem and the evolution solver indicated by the solver parameter.

Parameters:
HQobj, QobjEvo

System Hamiltonian, may be time-dependent for solver choice of me.

state0Qobj

Initial state density matrix \(\rho(t_0)\) or state vector \(\psi(t_0)\). If ‘state0’ is ‘None’, then the steady state will be used as the initial state. The ‘steady-state’ is only implemented if c_ops are provided and the Hamiltonian is constant.

taulistarray_like

List of times for \(\tau\). taulist must be positive and contain the element 0.

c_opslist of {Qobj, QobjEvo}

List of collapse operators

a_opQobj, QobjEvo

Operator A.

b_opQobj, QobjEvo

Operator B.

reversebool, default: False

If True, calculate \(\left<A(t)B(t+\tau)\right>\) instead of \(\left<A(t+\tau)B(t)\right>\).

solverstr {‘me’, ‘es’}, default: ‘me’

Choice of solver, me for master-equation, and es for exponential series. es is equivalent to me with options={"method": "diag"}.

optionsdict, optional

Options for the solver.

Returns:
corr_vecndarray

An array of correlation values for the times specified by taulist.

See also

correlation_3op

Similar function supporting various solver types.

References

See, Gardiner, Quantum Noise, Section 5.2.

correlation_2op_2t(
H,
state0,
tlist,
taulist,
c_ops,
a_op,
b_op,
solver='me',
reverse=False,
args=None,
options=None,
)[source]

Calculate the two-operator two-time correlation function: \(\left<A(t+\tau)B(t)\right>\) along two time axes using the quantum regression theorem and the evolution solver indicated by the solver parameter.

Parameters:
HQobj, QobjEvo

System Hamiltonian, may be time-dependent for solver choice of me.

state0Qobj

Initial state density matrix \(\rho(t_0)\) or state vector \(\psi(t_0)\). If ‘state0’ is ‘None’, then the steady state will be used as the initial state. The ‘steady-state’ is only implemented if c_ops are provided and the Hamiltonian is constant.

tlistarray_like

List of times for \(t\). tlist must be positive and contain the element 0. When taking steady-steady correlations only one tlist value is necessary, i.e. when \(t \rightarrow \infty\). If tlist is None, tlist=[0] is assumed.

taulistarray_like

List of times for \(\tau\). taulist must be positive and contain the element 0.

c_opslist of {Qobj, QobjEvo}

List of collapse operators

a_opQobj, QobjEvo

Operator A.

b_opQobj, QobjEvo

Operator B.

reversebool, default: False

If True, calculate \(\left<A(t)B(t+\tau)\right>\) instead of \(\left<A(t+\tau)B(t)\right>\).

solverstr {‘me’, ‘es’}, default: ‘me’

Choice of solver, me for master-equation, and es for exponential series. es is equivalent to me with options={"method": "diag"}.

optionsdict, optional

Options for the solver.

Returns:
corr_matndarray

An 2-dimensional array (matrix) of correlation values for the times specified by tlist (first index) and taulist (second index).

See also

correlation_3op

Similar function supporting various solver types.

References

See, Gardiner, Quantum Noise, Section 5.2.

correlation_3op(solver, state0, tlist, taulist, A=None, B=None, C=None)[source]

Calculate the three-operator two-time correlation function:

\(\left<A(t)B(t+\tau)C(t)\right>\).

from a open system Solver.

Note: it is not possible to calculate a physically meaningful correlation where \(\tau<0\).

Parameters:
solverMESolver, BRSolver

Qutip solver for an open system.

state0Qobj

Initial state density matrix \(\rho(t_0)\) or state vector \(\psi(t_0)\).

tlistarray_like

List of times for \(t\). tlist must be positive and contain the element 0.

taulistarray_like

List of times for \(\tau\). taulist must be positive and contain the element 0.

A, B, CQobj, QobjEvo, optional, default=None

Operators A, B, C from the equation <A(t)B(t+\tau)C(t)> in the Schrodinger picture. They do not need to be all provided. For exemple, if A is not provided, <B(t+\tau)C(t)> is computed.

Returns:
corr_matarray

An 2-dimensional array (matrix) of correlation values for the times specified by tlist (first index) and taulist (second index). If tlist is None, then a 1-dimensional array of correlation values is returned instead.

correlation_3op_1t(
H,
state0,
taulist,
c_ops,
a_op,
b_op,
c_op,
solver='me',
args=None,
options=None,
)[source]

Calculate the three-operator two-time correlation function: \(\left<A(0)B(\tau)C(0)\right>\) along one time axis using the quantum regression theorem and the evolution solver indicated by the solver parameter.

Note: it is not possibly to calculate a physically meaningful correlation of this form where \(\tau<0\).

Parameters:
HQobj, QobjEvo

System Hamiltonian, may be time-dependent for solver choice of me.

state0Qobj

Initial state density matrix \(\rho(t_0)\) or state vector \(\psi(t_0)\). If ‘state0’ is ‘None’, then the steady state will be used as the initial state. The ‘steady-state’ is only implemented if c_ops are provided and the Hamiltonian is constant.

taulistarray_like

List of times for \(\tau\). taulist must be positive and contain the element 0.

c_opslist of {Qobj, QobjEvo}

List of collapse operators

a_opQobj, QobjEvo

Operator A.

b_opQobj, QobjEvo

Operator B.

c_opQobj, QobjEvo

Operator C.

solverstr {‘me’, ‘es’}, default: ‘me’

Choice of solver, me for master-equation, and es for exponential series. es is equivalent to me with options={"method": "diag"}.

optionsdict, optional

Options for the solver.

Returns:
corr_vecarray

An array of correlation values for the times specified by taulist.

See also

correlation_3op

Similar function supporting various solver types.

References

See, Gardiner, Quantum Noise, Section 5.2.

correlation_3op_2t(
H,
state0,
tlist,
taulist,
c_ops,
a_op,
b_op,
c_op,
solver='me',
args=None,
options=None,
)[source]

Calculate the three-operator two-time correlation function: \(\left<A(t)B(t+\tau)C(t)\right>\) along two time axes using the quantum regression theorem and the evolution solver indicated by the solver parameter.

Note: it is not possibly to calculate a physically meaningful correlation of this form where \(\tau<0\).

Parameters:
HQobj, QobjEvo

System Hamiltonian, may be time-dependent for solver choice of me.

state0Qobj

Initial state density matrix \(\rho(t_0)\) or state vector \(\psi(t_0)\). If ‘state0’ is ‘None’, then the steady state will be used as the initial state. The ‘steady-state’ is only implemented if c_ops are provided and the Hamiltonian is constant.

tlistarray_like

List of times for \(t\). tlist must be positive and contain the element 0. When taking steady-steady correlations only one tlist value is necessary, i.e. when \(t \rightarrow \infty\). If tlist is None, tlist=[0] is assumed.

taulistarray_like

List of times for \(\tau\). taulist must be positive and contain the element 0.

c_opslist of {Qobj, QobjEvo}

List of collapse operators

a_opQobj, QobjEvo

Operator A.

b_opQobj, QobjEvo

Operator B.

c_opQobj, QobjEvo

Operator C.

solverstr {‘me’, ‘es’}, default: ‘me’

Choice of solver, me for master-equation, and es for exponential series. es is equivalent to me with options={"method": "diag"}.

optionsdict, optional

Options for the solver. Only used with me solver.

Returns:
corr_matarray

An 2-dimensional array (matrix) of correlation values for the times specified by tlist (first index) and taulist (second index).

See also

correlation_3op

Similar function supporting various solver types.

References

See, Gardiner, Quantum Noise, Section 5.2.

spectrum(H, wlist, c_ops, a_op, b_op, solver='es')[source]

Calculate the spectrum of the correlation function \(\lim_{t \to \infty} \left<A(t+\tau)B(t)\right>\), i.e., the Fourier transform of the correlation function:

\[S(\omega) = \int_{-\infty}^{\infty} \lim_{t \to \infty} \left<A(t+\tau)B(t)\right> e^{-i\omega\tau} d\tau.\]

using the solver indicated by the solver parameter. Note: this spectrum is only defined for stationary statistics (uses steady state rho0)

Parameters:
Hqobj

system Hamiltonian.

wlistarray_like

List of frequencies for \(\omega\).

c_opslist

List of collapse operators.

a_opQobj

Operator A.

b_opQobj

Operator B.

solverstr, {‘es’, ‘pi’, ‘solve’}, default: ‘es’

Choice of solver, es for exponential series and pi for psuedo-inverse, solve for generic solver.

Returns:
spectrumarray

An array with spectrum \(S(\omega)\) for the frequencies specified in wlist.

spectrum_correlation_fft(tlist, y, inverse=False)[source]

Calculate the power spectrum corresponding to a two-time correlation function using FFT.

Parameters:
tlistarray_like

list/array of times \(t\) which the correlation function is given.

yarray_like

list/array of correlations corresponding to time delays \(t\).

inverse: bool, default: False

boolean parameter for using a positive exponent in the Fourier Transform instead. Default is False.

Returns:
w, Stuple

Returns an array of angular frequencies ‘w’ and the corresponding two-sided power spectrum ‘S(w)’.

Steady-state Solvers

pseudo_inverse(L, rhoss=None, w=None, method='splu', *, use_rcm=False, **kwargs)[source]

Compute the pseudo inverse for a Liouvillian superoperator, optionally given its steady state density matrix (which will be computed if not given).

Parameters:
LQobj

A Liouvillian superoperator for which to compute the pseudo inverse.

rhossQobj, optional

A steadystate density matrix as Qobj instance, for the Liouvillian superoperator L.

wdouble, optional

frequency at which to evaluate pseudo-inverse. Can be zero for dense systems and large sparse systems. Small sparse systems can fail for zero frequencies.

sparsebool, optional

Flag that indicate whether to use sparse or dense matrix methods when computing the pseudo inverse.

methodstr, optional

Method used to compte matrix inverse. Choice are ‘pinv’ to use scipy’s function of the same name, or a linear system solver. Default supported solver are:

  • “solve”, “lstsq” dense solver from numpy.linalg

  • “spsolve”, “gmres”, “lgmres”, “bicgstab”, “splu” sparse solver from scipy.sparse.linalg

  • “mkl_spsolve”, sparse solver by mkl.

Extension to qutip, such as qutip-tensorflow, can use come with their own solver. When L use these data backends, see the corresponding libraries linalg for available solver.

use_rcmbool, default: False

Use reverse Cuthill-Mckee reordering to minimize fill-in in the LU factorization of the Liouvillian.

kwargsdictionary

Additional keyword arguments for setting parameters for solver methods.

Returns:
RQobj

Returns a Qobj instance representing the pseudo inverse of L.

Notes

In general the inverse of a sparse matrix will be dense. If you are applying the inverse to a density matrix then it is better to cast the problem as an Ax=b type problem where the explicit calculation of the inverse is not required. See page 67 of “Electrons in nanostructures” C. Flindt, PhD Thesis available online: https://orbit.dtu.dk/en/publications/electrons-in-nanostructures-coherent-manipulation-and-counting-st

Note also that the definition of the pseudo-inverse herein is different from numpys pinv() alone, as it includes pre and post projection onto the subspace defined by the projector Q.

steadystate(A, c_ops=[], *, method='direct', solver=None, **kwargs)[source]

Calculates the steady state for quantum evolution subject to the supplied Hamiltonian or Liouvillian operator and (if given a Hamiltonian) a list of collapse operators.

If the user passes a Hamiltonian then it, along with the list of collapse operators, will be converted into a Liouvillian operator in Lindblad form.

Parameters:
AQobj

A Hamiltonian or Liouvillian operator.

c_op_listlist

A list of collapse operators.

methodstr, {“direct”, “eigen”, “svd”, “power”}, default: “direct”

The allowed methods are composed of 2 parts, the steadystate method:

  • “direct”: Solving L(rho_ss) = 0

  • “eigen” : Eigenvalue problem

  • “svd” : Singular value decomposition

  • “power” : Inverse-power method

  • “propagator” : Repeatedly applying the propagator

solverstr, optional

‘direct’ and ‘power’ methods only. Solver to use when solving the L(rho_ss) = 0 equation. Default supported solver are:

  • “solve”, “lstsq” dense solver from numpy.linalg

  • “spsolve”, “gmres”, “lgmres”, “bicgstab” sparse solver from scipy.sparse.linalg

  • “mkl_spsolve” sparse solver by mkl.

Extension to qutip, such as qutip-tensorflow, can use come with their own solver. When A and c_ops use these data backends, see the corresponding libraries linalg for available solver.

Extra options for these solver can be passed in **kw.

use_rcmbool, default: False

Use reverse Cuthill-Mckee reordering to minimize fill-in in the LU factorization of the Liouvillian. Used with ‘direct’ or ‘power’ method.

use_wbmbool, default: False

Use Weighted Bipartite Matching reordering to make the Liouvillian diagonally dominant. This is useful for iterative preconditioners only. Used with ‘direct’ or ‘power’ method.

weightfloat, optional

Sets the size of the elements used for adding the unity trace condition to the linear solvers. This is set to the average abs value of the Liouvillian elements if not specified by the user. Used with ‘direct’ method.

power_tolfloat, default: 1e-12

Tolerance for the solution when using the ‘power’ method.

power_maxiterint, default: 10

Maximum number of iteration to use when looking for a solution when using the ‘power’ method.

power_eps: double, default: 1e-15

Small weight used in the “power” method.

sparse: bool, default: True

Whether to use the sparse eigen solver with the “eigen” method (default sparse). With “direct” and “power” method, when the solver is not specified, it is used to set whether “solve” or “spsolve” is used as default solver.

rho: Qobj, default: None

Initial state for the “propagator” method.

propagator_T: float, default: 10

Initial time step for the propagator method. The time step is doubled each iteration.

propagator_tol: float, default: 1e-5

Tolerance for propagator method convergence. If the Hilbert distance between the states of a step is less than this tolerance, the state is considered to have converged to the steady state.

propagator_max_iter: int, default: 30

Maximum number of iterations until convergence. A RuntimeError is raised if the state did not converge.

**kwargs

Extra options to pass to the linear system solver. See the documentation of the used solver in numpy.linalg or scipy.sparse.linalg to see what extra arguments are supported.

Returns:
dmqobj

Steady state density matrix.

infodict, optional

Dictionary containing solver-specific information about the solution.

Notes

The SVD method works only for dense operators (i.e. small systems).

steadystate_floquet(
H_0,
c_ops,
Op_t,
w_d=1.0,
n_it=3,
sparse=False,
solver=None,
**kwargs,
)[source]
Calculates the effective steady state for a driven

system with a time-dependent cosinusoidal term:

\[\mathcal{\hat{H}}(t) = \hat{H}_0 + \mathcal{\hat{O}} \cos(\omega_d t)\]
Parameters:
H_0Qobj

A Hamiltonian or Liouvillian operator.

c_opslist

A list of collapse operators.

Op_tQobj

The the interaction operator which is multiplied by the cosine

w_dfloat, default: 1.0

The frequency of the drive

n_itint, default: 3

The number of iterations for the solver

sparsebool, default: False

Solve for the steady state using sparse algorithms.

solverstr, optional

Solver to use when solving the linear system. Default supported solver are:

  • “solve”, “lstsq” dense solver from numpy.linalg

  • “spsolve”, “gmres”, “lgmres”, “bicgstab” sparse solver from scipy.sparse.linalg

  • “mkl_spsolve” sparse solver by mkl.

Extensions to qutip, such as qutip-tensorflow, may provide their own solvers. When H_0 and c_ops use these data backends, see their documentation for the names and details of additional solvers they may provide.

**kwargs:

Extra options to pass to the linear system solver. See the documentation of the used solver in numpy.linalg or scipy.sparse.linalg to see what extra arguments are supported.

Returns:
dmqobj

Steady state density matrix.

Notes

See: Sze Meng Tan, https://painterlab.caltech.edu/wp-content/uploads/2019/06/qe_quantum_optics_toolbox.pdf, Section (16)

Scattering in Quantum Optical Systems

Photon scattering in quantum optical systems

This module includes a collection of functions for numerically computing photon scattering in driven arbitrary systems coupled to some configuration of output waveguides. The implementation of these functions closely follows the mathematical treatment given in K.A. Fischer, et. al., Scattering of Coherent Pulses from Quantum Optical Systems (2017, arXiv:1710.02875).

scattering_probability(
H,
psi0,
n_emissions,
c_ops,
tlist,
system_zero_state=None,
construct_effective_hamiltonian=True,
)[source]

Compute the integrated probability of scattering n photons in an arbitrary system. This function accepts a nonlinearly spaced array of times.

Parameters:
HQobj or list

System-waveguide(s) Hamiltonian or effective Hamiltonian in Qobj or list-callback format. If construct_effective_hamiltonian is not specified, an effective Hamiltonian is constructed from H and c_ops.

psi0Qobj

Initial state density matrix \(\rho(t_0)\) or state vector \(\psi(t_0)\).

n_emissionsint

Number of photons emitted by the system (into any combination of waveguides).

c_opslist

List of collapse operators for each waveguide; these are assumed to include spontaneous decay rates, e.g. \(\sigma = \sqrt \gamma \cdot a\).

tlistarray_like

List of times for \(\tau_i\). tlist should contain 0 and exceed the pulse duration / temporal region of interest; tlist need not be linearly spaced.

system_zero_stateQobj, optional

State representing zero excitations in the system. Defaults to basis(systemDims, 0).

construct_effective_hamiltonianbool, default: True

Whether an effective Hamiltonian should be constructed from H and c_ops: \(H_{eff} = H - \frac{i}{2} \sum_n \sigma_n^\dagger \sigma_n\) Default: True.

Returns:
scattering_probfloat

The probability of scattering n photons from the system over the time range specified.

temporal_basis_vector(waveguide_emission_indices, n_time_bins)[source]

Generate a temporal basis vector for emissions at specified time bins into specified waveguides.

Parameters:
waveguide_emission_indiceslist or tuple

List of indices where photon emission occurs for each waveguide, e.g. [[t1_wg1], [t1_wg2, t2_wg2], [], [t1_wg4, t2_wg4, t3_wg4]].

n_time_binsint

Number of time bins; the range over which each index can vary.

Returns:
temporal_basis_vectorQobj

A basis vector representing photon scattering at the specified indices. If there are W waveguides, T times, and N photon emissions, then the basis vector has dimensionality (W*T)^N.

temporal_scattered_state(
H,
psi0,
n_emissions,
c_ops,
tlist,
system_zero_state=None,
construct_effective_hamiltonian=True,
)[source]

Compute the scattered n-photon state projected onto the temporal basis.

Parameters:
HQobj or list

System-waveguide(s) Hamiltonian or effective Hamiltonian in Qobj or list-callback format. If construct_effective_hamiltonian is not specified, an effective Hamiltonian is constructed from H and c_ops.

psi0Qobj

Initial state density matrix \(\rho(t_0)\) or state vector \(\psi(t_0)\).

n_emissionsint

Number of photon emissions to calculate.

c_opslist

List of collapse operators for each waveguide; these are assumed to include spontaneous decay rates, e.g. \(\sigma = \sqrt \gamma \cdot a\)

tlistarray_like

List of times for \(\tau_i\). tlist should contain 0 and exceed the pulse duration / temporal region of interest.

system_zero_stateQobj, optional

State representing zero excitations in the system. Defaults to \(\psi(t_0)\)

construct_effective_hamiltonianbool, default: True

Whether an effective Hamiltonian should be constructed from H and c_ops: \(H_{eff} = H - \frac{i}{2} \sum_n \sigma_n^\dagger \sigma_n\) Default: True.

Returns:
phi_nQobj

The scattered bath state projected onto the temporal basis given by tlist. If there are W waveguides, T times, and N photon emissions, then the state is a tensor product state with dimensionality T^(W*N).