qom.solvers.differential module¶
Module to solve ordinary differential equations.
- class qom.solvers.differential.ODESolver(func, params: dict, cb_update=None)¶
Bases:
objectClass to solve ordinary differential equations using
scipy.integrate.Initializes
solver_methods,func,params,integratorandupdater.- Parameters:
func (callable) – Function returning the rate equations of the input variables, formatted as
func(t, v, c), wheretis the time at which the integration is performed,vis a list of variables andcis a list of constants. The output should match the dimension ofv.params (dict) – Parameters for the solver. Refer to Notes below for all available options.
cb_update (callable, optional) – Callback function to update status and progress, formatted as
cb_update(status, progress, reset), wherestatusis a string,progressis a float andresetis a boolean.
Notes
- The
paramsdictionary currently supports the following keys: key
value
‘show_progress’
(bool) option to display the progress of the integration. 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.ODESolver). 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.- Currently available Python-based methods are:
value
meaning
‘BDF’
backward-differentiation formulas.
‘DOP853’
explicit Runge-Kutta method of order 8(5, 3).
‘LSODA’
Adams/BDF method with automatic stiffness detection and switching.
‘Radau’
implicit Runge-Kutta of Radau IIA family of order 5.
‘RK23’
explicit Runge-Kutta method of order 3(2).
‘RK45’
explicit Runge-Kutta method of order 5(4) (fallback).
- Currently available FORTRAN-based are:
value
meaning
‘dop853’
explicit Runge-Kutta method of order 8(5, 3).
‘dopri5’
explicit Runge-Kutta method of order 5(4).
‘lsoda’
real-valued Adams/BDF method with automatic stiffness detection and switching.
‘vode’
real-valued implicit Adams/BDF methods.
‘zvode’
complex-valued implicit Adams/BDF methods.
- desc = 'Ordinary Differential Equations Solver'¶
Description of the solver.
- Type:
str
- name = 'ODESolver'¶
Name of the solver.
- Type:
str
- new_api_methods = ['BDF', 'DOP853', 'LSODA', 'Radau', 'RK23', 'RK45']¶
New Python-based methods availabile in
scipy.integrate.- Type:
list
- old_api_methods = ['dop853', 'dopri5', 'lsoda', 'vode', 'zvode']¶
Old FORTRAN-based methods availabile in
scipy.integrate.- Type:
list
- set_params(params)¶
Method to set the solver parameters.
- Parameters:
params (dict) – Parameters of the solver.
- solve(T, iv, c=None, func_c=None)¶
Method to obtain the solutions of the ODE at all times.
- Parameters:
T (numpy.ndarray) – Times at which the values are calculated.
iv (numpy.ndarray) – Initial values for the integration.
c (numpy.ndarray) – Constants of the integration.
func_c (callable, optional) – Function returning the time-dependent constants of the integration, formatted as
func_c(i), whereiis the i-th step of integration.
- Returns:
vs – Values of the variables at all times.
- Return type:
numpy.ndarray
- solve_new(T, iv, c)¶
Method to integrate with the new API methods.
- Parameters:
T (float) – Times at which the values are obtained.
iv (numpy.ndarray) – Initial values of the variables.
c (numpy.ndarray) – Constants of the integration.
- Returns:
vs – Values of the variables.
- Return type:
numpy.ndarray
- solver_defaults = {'ode_atol': 1e-12, 'ode_is_stiff': False, 'ode_method': 'RK45', 'ode_rtol': 1e-06, 'show_progress': False}¶
Default parameters of the solver.
- Type:
dict