qom.solvers.stochastic module¶
Module to solve stochastic equations of motion.
- class qom.solvers.stochastic.MCQTSolver(system, params: dict, num_trajs: int, cb_update=None, parallel: bool = False, p_index: int = 0, p_start: float = None)¶
Bases:
object
Class to solve for the expectation values of operators using the Monte-Carlo quantum trajectories method.
Initializes
system
,num_trajs
,T
andupdater
.- Parameters:
system (
qom.systems.*
) – Instance of the system. Requires predefined system methodsget_ops_collapse
,get_ops_expect
,get_H_0
andget_ivc
. For time-dependent Hamiltonians, the methodget_H_t
should also be defined. Refer to Notes below for their implementations.params (dict) –
- Parameters for the solver. Refer to Notes below for all available options. Required options are:
key
value
’t_min’
(float) minimum time at which integration starts.
’t_max’
(float) maximum time at which integration stops.
’t_dim’
(int) number of values from
't_max'
to't_min'
, both inclusive.
num_trajs (int) – Number of trajectories.
cb_update (callable, optional) – Callback function to update status and progress, formatted as
cb_update(status, progress, reset)
, wherestatus
is a string,progress
is a float andreset
is a boolean.parallel (bool, default=False) – Option to format outputs when running in parallel.
p_index (int, default=0) – Index of the process.
p_start (float, optional) – Time at which the process was started. If not provided, the value is initialized to current time.
Notes
- The following system methods follow a strict formatting:
method
returns
get_H_0
the time-independent part of the Hamiltonian, formatted as
get_H_0(c)
, wherec
is an array of the derived constants and controls of the system.get_H_t
the time-dependent part of the Hamiltonian, formatted as
get_H_t(c, t)
, wherec
is an array of the derived constants and controls of the system andt
is the time at which the Hamiltonian is obtained.get_ops_collapse
the collapse operators of the system. It follows the same formatting as
get_H_0
.get_ops_expect
the operators of the system for which the expectation values are to be obtained. It follows the same formatting as
get_H_0
.get_ivc
the initial state vector and the derived constants and controls, formatted as
get_ivc()
.- The
params
dictionary currently supports the following keys: key
value
‘show_progress’
(bool) option to display the progress of the solver. Default is
False
.‘ode_method’
(str) method used to solve the ODEs. Available options are
'BDF'
,'DOP853'
,'LSODA'
,'Radau'
,'RK23'
,'RK45'
(fallback),'dop853'
,'dopri5'
,'lsoda'
,'vode'
and'zvode'
. Refer toqom.solvers.differential.ODESolver
for details of each method. Default is'RK45'
.‘ode_is_stiff’
(bool) option to select whether the integration is a stiff problem or a non-stiff one. Default is
False
.‘ode_atol’
(float) absolute tolerance of the integrator. Default is
1e-12
.‘ode_rtol’
(float) relative tolerance of the integrator. Default is
1e-6
.‘t_min’
(float) minimum time at which integration starts. Default is
0.0
.‘t_max’
(float) maximum time at which integration stops. Default is
1000.0
.‘t_dim’
(int) number of values from
't_max'
to't_min'
, both inclusive. Default is10001
.
- desc = 'Monte-Carlo Quantum Trajectories Solver'¶
Description of the solver.
- Type:
str
- name = 'MCQTSolver'¶
Name of the solver.
- Type:
str
- set_params(params)¶
Method to set the solver parameters.
- Parameters:
params (dict) – Parameters of the solver.
- solve()¶
Method to solve for the quantum trajectories.
Sets the results with keys
'times'
,'trajs'
,'expects'
and'runtimes'
for the time steps, individual trajectories of expectation values, pooled expectation values and execution times of each trajectory respectively.
- solver_defaults = {'ode_atol': 1e-08, 'ode_is_stiff': False, 'ode_method': 'zvode', 'ode_rtol': 1e-06, 'show_progress': False, 't_dim': 10001, 't_max': 1000.0, 't_min': 0.0}¶
Default parameters of the solver.
- Type:
dict