odespy.odepack

class odespy.odepack.Lsoda(f, **kwargs)[source]

Bases: odespy.odepack.Odepack

A Python wrapper of the LSODA FORTRAN subroutine from ODEPACK. This subroutine automatically shifts between stiff (BDF) and nonstiff (Adams) methods, such that the user does not need to determine whether the problem is stiff or not (as is the case when using the LSODE subroutine and class Lsode). The integration always starts with the nonstiff method. Can generate the Jacobian, or apply a user-supplied Jacobian in dense or banded format.

Required input arguments:

Name Description
f Right-hand side f(u,t) defining the ODE.

Optional input arguments:

Name Description
f_args Extra positional arguments to f: f(u, t, *f_args, **f_kwargs). (default: ())
f_kwargs Extra keyword arguments to f: f(u, t, *f_args, **f_kwargs). (default: {})
complex_valued True if f is complex valued. (default: False)
disk_storage Indicates whether u is stored in memory or in file. If string, it is the filename; if False or “”, u is kept in memory; if True, a default filename tmp_odspy.dat is used. (default: False)
verbose Integer reflecting output of intermediate quantities. (default: 0)
u_exact Function of t returning exact solution. (default: None)
atol Absolute tolerance for solution. (default: 1e-08)
rtol Relative tolerance for solution. (default: 1e-06)
adams_or_bdf Method in vode or solvers in odepack: “adams” or “bdf”. (default: adams)
order Maximum order used by the integrator (<= 12 for “adams”, <= 5 for “bdf”). (default: 4)
nsteps Max no of internal solver steps per time step. (default: 1000)
first_step Suggested first time step size for an adaptive algorithm.
min_step Minimum step size for an adaptive algorithm.
max_step Maximum step size for an adaptive algorithm.
corrector_iter_method Corrector iteration method choice.
lrw Length of real work array.
liw Length of integer work array, similiar as <lrw>.
jac_banded Function for Jacobian in banded matrix form. Used in Lsode, Lsoda, Lsodar. jac_banded(u,t,ml,mu) returns df/du as an array of size neq times ml+mu+1.
jac_banded_f77

Fortran subroutine for jac_banded. This subroutine has the signature:

      subroutine jac_banded_f77
     1  (neq,t,u,ml, mu,pd,nrowpd)
Cf2py intent(hide) neq,ml,mu,nrowpd
Cf2py intent(out) pd
      integer neq,ml,mu,nrowpd
      double precision t,u,pd
      dimension u(neq),pd(nrowpd,neq)
      pd = ...
      return
      end
ml Number of lower non-zero diagonals in a banded Jacobian.
mu Number of upper non-zero diagonals in a banded Jacobian.
jac Jacobian of right-hand side function f (df/du). (default: None)
jac_f77

Fortran subroutine for jac. This subroutine has the signature:

      subroutine jac_f77
     1 (neq, t, u, ml, mu, pd, nrowpd)
Cf2py intent(hide) neq,ml,mu,nrowpd
Cf2py intent(out) pd
      integer neq,ml,mu,nrowpd
      double precision t,u,pd
      dimension u(neq),pd(neq,neq)
      pd = ...
      return
      end
max_ordn Maximum order in nonstiff methods.
max_ords Maximum order in stiff methods.
f_f77

Fortran subroutine for f. This subroutine has the signature:

      subroutine f_f77(neq,t,u,udot)
Cf2py intent(hide)   neq
Cf2py intent(out)    udot
      integer neq
      double precision t,u(neq),udot(neq)
      udot = ...
      return
      end

Methods

adjust_atol(u_current) Error tolerance tol(i) may become zero for some i
adjust_parameters() Properties for new parameters in this solver.
advance() This function intends to one step forward for linearly implicit solvers
check_conditional_parameters() This function is used to check whether conditional parameters are provided when specified condition fulfilled.
check_extra(**kwargs) A parameter may have a keyword extra_check for user-given functions that performs consistency checks on the parameter.
check_iaja() ia, ja, ic, jc are optional inputs to describe
check_input_range(**kwargs) Check whether all existing inputs are in right specified range.
check_input_types(**kwargs) Check whether all existing inputs are of right specified type.
check_liwlrw() If the lengths of work arrays are specified by users, check whether
check_pars() Some pairs of parameters should be input simutaneously.
check_tol() atol & rtol should be defined as scalars or vectors with length neq.
compile_string_functions(f, **kwargs) Compile functions which are supplied as Fortran strings.
constant_time_step() Check if self.t has a uniform partition.
expand_iwork(new_liw[, expand]) Extension module return an actually required length for integer work array when it is too short.
expand_rwork(new_lrw[, expand]) Length of real work array is smaller than actually required length.
func_wrapper() This function is defined to wrap user-defined functions with new forms of parameter-list, or wrap the returned values as numpy arrays.
get([parameter_name, print_info]) Return value of specified input parameters.
get_parameter_info([print_info]) Return a dictionary containing all properties of all legal parameters in current subclass (i.e., the parameters in self._parameters).
has_u_t_all() Return True if self.u_all and self.t_all, defined in
initialize() Import extension module _odesolver and check that it exists.
initialize_for_solve() In the long parameter-lists for solvers in ODEPACK, quite a few
new_stepnr() When Fortran code returns a status istate==-1, it indicates that
print_roots(jroot, t_current, u_current) Roots found at current T for some constraint functions.
set([strict]) Assign values to one or more parameters, specified as keyword arguments.
set_dummy_functions() Functions have to get dummy values before they are passed to extension
set_extra_args()
set_initial_condition(U0) Function set_initial_condition() is used to set initial value of
set_iopt()
set_iter_method()
set_iwork_rwork() Initialize arrays for optional inputs, and calculate the
set_jac()
set_liw_min()
set_lrw_min()
set_ydoti() ydoti is an array used in linearly implicit solvers.
solve(time_points[, terminate]) This function is involved for non-linearly implicit
switch_to(solver_target[, print_info]) Create a new solver instance which switch to another subclass with same values of common attributes.
tol_multiply(tolsf) This function is used to adjust tolerance parameters for Fortran part.
validate_data() Common validity check in Odepack.
__module__ = 'odespy.odepack'
adjust_parameters()[source]

Properties for new parameters in this solver.

quick_description = 'LSODA solver with stiff-nonstiff auto shift'
set_extra_args()[source]
set_iopt()[source]
set_iter_method()[source]
set_jac()[source]
set_liw_min()[source]
set_lrw_min()[source]
class odespy.odepack.Lsodar(f, **kwargs)[source]

Bases: odespy.odepack.Odepack

A Python wrapper of the LSODAR subroutine in ODEPACK. LSODAR is a variant of LSODE and differs from the latter in two ways.

Quote from the LSODAR source code documentation: “(a) It switches automatically between stiff and nonstiff methods. This means that the user does not have to determine whether the problem is stiff or not, and the solver will automatically choose the appropriate method. It always starts with the nonstiff method. (b) It finds the root of at least one of a set of constraint functions g(i) of the independent and dependent variables. It finds only those roots for which some g(i), as a function of t, changes sign in the interval of integration. It then returns the solution at the root, if that occurs sooner than the specified stop condition, and otherwise returns the solution according the specified stop condition.”

The mathematical problem reads

..:math::
u’ &= f(u, t),

g(u,t) &= 0.

Required input arguments:

Name Description
f Right-hand side f(u,t) defining the ODE.

Optional input arguments:

Name Description
f_args Extra positional arguments to f: f(u, t, *f_args, **f_kwargs). (default: ())
f_kwargs Extra keyword arguments to f: f(u, t, *f_args, **f_kwargs). (default: {})
complex_valued True if f is complex valued. (default: False)
disk_storage Indicates whether u is stored in memory or in file. If string, it is the filename; if False or “”, u is kept in memory; if True, a default filename tmp_odspy.dat is used. (default: False)
verbose Integer reflecting output of intermediate quantities. (default: 0)
u_exact Function of t returning exact solution. (default: None)
atol Absolute tolerance for solution. (default: 1e-08)
rtol Relative tolerance for solution. (default: 1e-06)
adams_or_bdf Method in vode or solvers in odepack: “adams” or “bdf”. (default: adams)
order Maximum order used by the integrator (<= 12 for “adams”, <= 5 for “bdf”). (default: 4)
nsteps Max no of internal solver steps per time step. (default: 1000)
first_step Suggested first time step size for an adaptive algorithm.
min_step Minimum step size for an adaptive algorithm.
max_step Maximum step size for an adaptive algorithm.
corrector_iter_method Corrector iteration method choice.
lrw Length of real work array.
liw Length of integer work array, similiar as <lrw>.
jac_banded Function for Jacobian in banded matrix form. Used in Lsode, Lsoda, Lsodar. jac_banded(u,t,ml,mu) returns df/du as an array of size neq times ml+mu+1.
jac_banded_f77

Fortran subroutine for jac_banded. This subroutine has the signature:

      subroutine jac_banded_f77
     1  (neq,t,u,ml, mu,pd,nrowpd)
Cf2py intent(hide) neq,ml,mu,nrowpd
Cf2py intent(out) pd
      integer neq,ml,mu,nrowpd
      double precision t,u,pd
      dimension u(neq),pd(nrowpd,neq)
      pd = ...
      return
      end
ml Number of lower non-zero diagonals in a banded Jacobian.
mu Number of upper non-zero diagonals in a banded Jacobian.
jac Jacobian of right-hand side function f (df/du). (default: None)
jac_f77

Fortran subroutine for jac. This subroutine has the signature:

      subroutine jac_f77
     1 (neq, t, u, ml, mu, pd, nrowpd)
Cf2py intent(hide) neq,ml,mu,nrowpd
Cf2py intent(out) pd
      integer neq,ml,mu,nrowpd
      double precision t,u,pd
      dimension u(neq),pd(neq,neq)
      pd = ...
      return
      end
max_ordn Maximum order in nonstiff methods.
max_ords Maximum order in stiff methods.
g Callable object to define constraint functions. g(u, t) returns a vector of the values of the constraints (left-hand sides in the constraint equations).
g_f77

Fortran subroutine for g. This subroutine has the signature:

      subroutine g_f77(neq, t, u, ng, groot)
Cf2py intent(hide) neq
Cf2py optional, intent(hide) ng
Cf2py intent(in) t, u
Cf2py intent(out) groot
      integer neq, ng
      double precision t, u, groot
      dimension u(neq), groot(ng)
      groot = ...
      return
      end
ng No of components in constraint function g.
f_f77

Fortran subroutine for f. This subroutine has the signature:

      subroutine f_f77(neq,t,u,udot)
Cf2py intent(hide)   neq
Cf2py intent(out)    udot
      integer neq
      double precision t,u(neq),udot(neq)
      udot = ...
      return
      end

Methods

adjust_atol(u_current) Error tolerance tol(i) may become zero for some i
adjust_parameters() Properties for new parameters in this solver.
advance() This function intends to one step forward for linearly implicit solvers
check_conditional_parameters() This function is used to check whether conditional parameters are provided when specified condition fulfilled.
check_extra(**kwargs) A parameter may have a keyword extra_check for user-given functions that performs consistency checks on the parameter.
check_iaja() ia, ja, ic, jc are optional inputs to describe
check_input_range(**kwargs) Check whether all existing inputs are in right specified range.
check_input_types(**kwargs) Check whether all existing inputs are of right specified type.
check_liwlrw() If the lengths of work arrays are specified by users, check whether
check_pars() Some pairs of parameters should be input simutaneously.
check_tol() atol & rtol should be defined as scalars or vectors with length neq.
compile_string_functions(f, **kwargs) Compile functions which are supplied as Fortran strings.
constant_time_step() Check if self.t has a uniform partition.
expand_iwork(new_liw[, expand]) Extension module return an actually required length for integer work array when it is too short.
expand_rwork(new_lrw[, expand]) Length of real work array is smaller than actually required length.
func_wrapper() This function is defined to wrap user-defined functions with new forms of parameter-list, or wrap the returned values as numpy arrays.
get([parameter_name, print_info]) Return value of specified input parameters.
get_parameter_info([print_info]) Return a dictionary containing all properties of all legal parameters in current subclass (i.e., the parameters in self._parameters).
has_u_t_all() Return True if self.u_all and self.t_all, defined in
initialize() Import extension module _odesolver and check that it exists.
initialize_for_solve() In the long parameter-lists for solvers in ODEPACK, quite a few
new_stepnr() When Fortran code returns a status istate==-1, it indicates that
print_roots(jroot, t_current, u_current) Roots found at current T for some constraint functions.
set([strict]) Assign values to one or more parameters, specified as keyword arguments.
set_dummy_functions() Functions have to get dummy values before they are passed to extension
set_extra_args()
set_initial_condition(U0) Function set_initial_condition() is used to set initial value of
set_iopt()
set_iter_method()
set_iwork_rwork() Initialize arrays for optional inputs, and calculate the
set_jac()
set_liw_min()
set_lrw_min()
set_ydoti() ydoti is an array used in linearly implicit solvers.
solve(time_points[, terminate])
switch_to(solver_target[, print_info]) Create a new solver instance which switch to another subclass with same values of common attributes.
tol_multiply(tolsf) This function is used to adjust tolerance parameters for Fortran part.
validate_data() Common validity check in Odepack.
__module__ = 'odespy.odepack'
adjust_parameters()[source]

Properties for new parameters in this solver.

quick_description = 'LSODAR method with stiff-nonstiff auto shift'
set_extra_args()[source]
set_iopt()[source]
set_iter_method()[source]
set_jac()[source]
set_liw_min()[source]
set_lrw_min()[source]
solve(time_points, terminate=None)[source]
class odespy.odepack.Lsode(f, **kwargs)[source]

Bases: odespy.odepack.Odepack

A Python wrapper of the LSODE (Livermore Solver for Ordinary Differential Equations) FORTRAN subroutine. Basic Solver in ODEPACK package from netib. Solves the initial-value problem for stiff or nonstiff systems of first-order ODE, \(u' = f(u,t)\) via Adams or BDF methods (for the nonstiff and stiff cases, respectively). Can generate the Jacobian, or apply a user-supplied Jacobian in dense or banded format.

Required input arguments:

Name Description
f Right-hand side f(u,t) defining the ODE.

Optional input arguments:

Name Description
f_args Extra positional arguments to f: f(u, t, *f_args, **f_kwargs). (default: ())
f_kwargs Extra keyword arguments to f: f(u, t, *f_args, **f_kwargs). (default: {})
complex_valued True if f is complex valued. (default: False)
disk_storage Indicates whether u is stored in memory or in file. If string, it is the filename; if False or “”, u is kept in memory; if True, a default filename tmp_odspy.dat is used. (default: False)
verbose Integer reflecting output of intermediate quantities. (default: 0)
u_exact Function of t returning exact solution. (default: None)
atol Absolute tolerance for solution. (default: 1e-08)
rtol Relative tolerance for solution. (default: 1e-06)
adams_or_bdf Method in vode or solvers in odepack: “adams” or “bdf”. (default: adams)
order Maximum order used by the integrator (<= 12 for “adams”, <= 5 for “bdf”). (default: 4)
nsteps Max no of internal solver steps per time step. (default: 1000)
first_step Suggested first time step size for an adaptive algorithm.
min_step Minimum step size for an adaptive algorithm.
max_step Maximum step size for an adaptive algorithm.
corrector_iter_method Corrector iteration method choice.
lrw Length of real work array.
liw Length of integer work array, similiar as <lrw>.
jac_banded Function for Jacobian in banded matrix form. Used in Lsode, Lsoda, Lsodar. jac_banded(u,t,ml,mu) returns df/du as an array of size neq times ml+mu+1.
jac_banded_f77

Fortran subroutine for jac_banded. This subroutine has the signature:

      subroutine jac_banded_f77
     1  (neq,t,u,ml, mu,pd,nrowpd)
Cf2py intent(hide) neq,ml,mu,nrowpd
Cf2py intent(out) pd
      integer neq,ml,mu,nrowpd
      double precision t,u,pd
      dimension u(neq),pd(nrowpd,neq)
      pd = ...
      return
      end
ml Number of lower non-zero diagonals in a banded Jacobian.
mu Number of upper non-zero diagonals in a banded Jacobian.
jac Jacobian of right-hand side function f (df/du). (default: None)
jac_f77

Fortran subroutine for jac. This subroutine has the signature:

      subroutine jac_f77
     1 (neq, t, u, ml, mu, pd, nrowpd)
Cf2py intent(hide) neq,ml,mu,nrowpd
Cf2py intent(out) pd
      integer neq,ml,mu,nrowpd
      double precision t,u,pd
      dimension u(neq),pd(neq,neq)
      pd = ...
      return
      end
order Maximum order used by the integrator (<= 12 for “adams”, <= 5 for “bdf”). (default: 4)
f_f77

Fortran subroutine for f. This subroutine has the signature:

      subroutine f_f77(neq,t,u,udot)
Cf2py intent(hide)   neq
Cf2py intent(out)    udot
      integer neq
      double precision t,u(neq),udot(neq)
      udot = ...
      return
      end

Methods

adjust_atol(u_current) Error tolerance tol(i) may become zero for some i
adjust_parameters() Properties for new parameters in this solver.
advance() This function intends to one step forward for linearly implicit solvers
check_conditional_parameters() This function is used to check whether conditional parameters are provided when specified condition fulfilled.
check_extra(**kwargs) A parameter may have a keyword extra_check for user-given functions that performs consistency checks on the parameter.
check_iaja() ia, ja, ic, jc are optional inputs to describe
check_input_range(**kwargs) Check whether all existing inputs are in right specified range.
check_input_types(**kwargs) Check whether all existing inputs are of right specified type.
check_liwlrw() If the lengths of work arrays are specified by users, check whether
check_pars() Some pairs of parameters should be input simutaneously.
check_tol() atol & rtol should be defined as scalars or vectors with length neq.
compile_string_functions(f, **kwargs) Compile functions which are supplied as Fortran strings.
constant_time_step() Check if self.t has a uniform partition.
expand_iwork(new_liw[, expand]) Extension module return an actually required length for integer work array when it is too short.
expand_rwork(new_lrw[, expand]) Length of real work array is smaller than actually required length.
func_wrapper() This function is defined to wrap user-defined functions with new forms of parameter-list, or wrap the returned values as numpy arrays.
get([parameter_name, print_info]) Return value of specified input parameters.
get_parameter_info([print_info]) Return a dictionary containing all properties of all legal parameters in current subclass (i.e., the parameters in self._parameters).
has_u_t_all() Return True if self.u_all and self.t_all, defined in
initialize() Import extension module _odesolver and check that it exists.
initialize_for_solve() In the long parameter-lists for solvers in ODEPACK, quite a few
new_stepnr() When Fortran code returns a status istate==-1, it indicates that
print_roots(jroot, t_current, u_current) Roots found at current T for some constraint functions.
set([strict]) Assign values to one or more parameters, specified as keyword arguments.
set_dummy_functions() Functions have to get dummy values before they are passed to extension
set_extra_args()
set_initial_condition(U0) Function set_initial_condition() is used to set initial value of
set_iopt()
set_iter_method()
set_iwork_rwork() Initialize arrays for optional inputs, and calculate the
set_jac()
set_liw_min()
set_lrw_min()
set_ydoti() ydoti is an array used in linearly implicit solvers.
solve(time_points[, terminate]) This function is involved for non-linearly implicit
switch_to(solver_target[, print_info]) Create a new solver instance which switch to another subclass with same values of common attributes.
tol_multiply(tolsf) This function is used to adjust tolerance parameters for Fortran part.
validate_data() Common validity check in Odepack.
__module__ = 'odespy.odepack'
adjust_parameters()[source]

Properties for new parameters in this solver.

quick_description = 'LSODE solver for a stiff or nonstiff system'
set_extra_args()[source]
set_iopt()[source]
set_iter_method()[source]
set_jac()[source]
set_liw_min()[source]
set_lrw_min()[source]
class odespy.odepack.Lsodes(f, **kwargs)[source]

Bases: odespy.odepack.Odepack

A Python wrapper of the LSODES subroutine from ODEPACK. LSODES is a variant of LSODE intended for problems where the Jacobian is provided via a sparse matrix data structure consisting of the parameters jac_column, ia, ja, and optionally jac_column_f77.

Required input arguments:

Name Description
f Right-hand side f(u,t) defining the ODE.

Optional input arguments:

Methods

adjust_atol(u_current) Error tolerance tol(i) may become zero for some i
adjust_parameters() Properties for new parameters in this solver.
advance() This function intends to one step forward for linearly implicit solvers
check_conditional_parameters() This function is used to check whether conditional parameters are provided when specified condition fulfilled.
check_extra(**kwargs) A parameter may have a keyword extra_check for user-given functions that performs consistency checks on the parameter.
check_iaja() ia, ja, ic, jc are optional inputs to describe
check_input_range(**kwargs) Check whether all existing inputs are in right specified range.
check_input_types(**kwargs) Check whether all existing inputs are of right specified type.
check_liwlrw() If the lengths of work arrays are specified by users, check whether
check_pars() Some pairs of parameters should be input simutaneously.
check_tol() atol & rtol should be defined as scalars or vectors with length neq.
compile_string_functions(f, **kwargs) Compile functions which are supplied as Fortran strings.
constant_time_step() Check if self.t has a uniform partition.
expand_iwork(new_liw[, expand]) Extension module return an actually required length for integer work array when it is too short.
expand_rwork(new_lrw[, expand]) Length of real work array is smaller than actually required length.
func_wrapper() This function is defined to wrap user-defined functions with new forms of parameter-list, or wrap the returned values as numpy arrays.
get([parameter_name, print_info]) Return value of specified input parameters.
get_parameter_info([print_info]) Return a dictionary containing all properties of all legal parameters in current subclass (i.e., the parameters in self._parameters).
has_u_t_all() Return True if self.u_all and self.t_all, defined in
initialize() Import extension module _odesolver and check that it exists.
initialize_for_solve() In the long parameter-lists for solvers in ODEPACK, quite a few
new_stepnr() When Fortran code returns a status istate==-1, it indicates that
print_roots(jroot, t_current, u_current) Roots found at current T for some constraint functions.
set([strict]) Assign values to one or more parameters, specified as keyword arguments.
set_dummy_functions() Functions have to get dummy values before they are passed to extension
set_extra_args() Setting for extra parameters of user-defined functions.
set_initial_condition(U0) Function set_initial_condition() is used to set initial value of
set_iopt()
set_iter_method()
set_iwork_rwork() Initialization of work arrays with caculated length and optional inputs.
set_jac() Set values for Jacobian matrix. In several solvers like Lsode,
set_liw_min() Calculate the necessary length of integer work arrays when it is not specified explicitly by users.
set_lrw_min() Calculate the necessary length of real work arrays for Fortran code.
set_ydoti() ydoti is an array used in linearly implicit solvers.
solve(time_points[, terminate]) This function is involved for non-linearly implicit
switch_to(solver_target[, print_info]) Create a new solver instance which switch to another subclass with same values of common attributes.
tol_multiply(tolsf) This function is used to adjust tolerance parameters for Fortran part.
validate_data() Common validity check in Odepack.
__module__ = 'odespy.odepack'
adjust_parameters()[source]

Properties for new parameters in this solver.

quick_description = 'LSODES solver for sparse Jacobians'
set_iopt()[source]
set_iter_method()[source]
set_iwork_rwork()[source]

Initialization of work arrays with caculated length and optional inputs. In ODEPACK, “iwork” & “rwork” should be initialized with the specific optional parameters in all the solvers. “liw” & “lrw” represented the length requirement of work arrays. Specially, in Dlsodes, ia & ja should be attached to iwork_in.

class odespy.odepack.Lsodi(**kwargs)[source]

Bases: odespy.odepack.Odepack

A Python wrapper of the LSODI subroutine from ODEPACK, targeting ODE problems of the form

\[\begin{split}A(u,t) u' & = g(u,t)\end{split}\]

where \(A(u,t)\) is a square matrix and \(g\) some function. If \(A\) is singular, this is a differential-algebraic system. The user input is in the form of a residual function \(r = g - As\) for some vector \(s\). The residual \(r\) is provided by the user function res (Python) or res_f77 (FORTRAN).

Required input arguments:

Name Description
   

Optional input arguments:

Name Description
f_args Extra positional arguments to f: f(u, t, *f_args, **f_kwargs). (default: ())
f_kwargs Extra keyword arguments to f: f(u, t, *f_args, **f_kwargs). (default: {})
complex_valued True if f is complex valued. (default: False)
disk_storage Indicates whether u is stored in memory or in file. If string, it is the filename; if False or “”, u is kept in memory; if True, a default filename tmp_odspy.dat is used. (default: False)
verbose Integer reflecting output of intermediate quantities. (default: 0)
u_exact Function of t returning exact solution. (default: None)
atol Absolute tolerance for solution. (default: 1e-08)
rtol Relative tolerance for solution. (default: 1e-06)
adams_or_bdf Method in vode or solvers in odepack: “adams” or “bdf”. (default: adams)
order Maximum order used by the integrator (<= 12 for “adams”, <= 5 for “bdf”). (default: 4)
nsteps Max no of internal solver steps per time step. (default: 1000)
first_step Suggested first time step size for an adaptive algorithm.
min_step Minimum step size for an adaptive algorithm.
max_step Maximum step size for an adaptive algorithm.
corrector_iter_method Corrector iteration method choice.
lrw Length of real work array.
liw Length of integer work array, similiar as <lrw>.
adda_lsodi Callable object to add the matrix A = A(u,t) to another matrix p stored in the same form as A. The signature is adda(u,t,p) and it returns a matrix p+A (square matrix as the Jacobian).
adda_banded_lsodi Callable object to add the banded matrix A = A(u,t) to another matrix stored P in the same form as A. For a banded matrix, A(i,j) is added to P(i-j+mu+1,j). The signature is adda(u,t,p,ml,mu) and it returns a banded matrix.
adda_lsodi_f77

Fortran subroutine for adda_lsodi. This subroutine has the signature:

      subroutine adda_lsodi_f77
     1  (neq, t, u, ml, mu, pd, nrowpd)
Cf2py intent(in, hide) neq, ml, mu
Cf2py intent(in, hide), depend(pd) nrowpd
Cf2py intent(in, out) pd
      integer neq, ml, mu, nrowpd
      double precision t, u, pd
      dimension u(neq), pd(nrowpd, neq)
      pd = ...
      return
      end
adda_banded_lsodi_f77

Fortran subroutine for adda_banded. This subroutine has the signature:

      subroutine adda_banded_lsodi_f77(neq, t,
     1                  u, ml, mu, pd, nrowpd)
Cf2py intent(in, hide) neq, ml, mu
Cf2py intent(in, hide), depend(pd) nrowpd
Cf2py intent(in, out) pd
      integer neq, ml, mu, nrowpd
      double precision t, u, pd
      dimension u(neq), pd(nrowpd, neq)
      pd = ...
      return
      end
jac_lsodi Callable object to define the full Jacobian matrix dr/du where r = g - A*s. The signature of this function is jac(u,t,s), returning a matrix like the other jac functions.
jac_banded_lsodi Callable object to define the banded Jacobian matrix dr/du where r = g - A*s. The signature is jac(u,t,s,ml,mu), where ml and mu are the lower and upper half-bandwidth of the banded matrix.
jac_lsodi_f77

Fortran subroutine for jac_lsodi. This subroutine has the signature:

      subroutine jac_lsodi_f77
     1  (neq, t, u, s, ml, mu, pd, nrowpd)
Cf2py intent(in, hide) neq, ml, mu, nrowpd
Cf2py intent(out) pd
      integer neq, ml, mu, nrowpd
      double precision t, u, pd, s
      dimension u(neq), s(neq), pd(nrowpd, neq)
      pd = ...
      return
      end
jac_banded_lsodi_f77

Fortran subroutine for jac_banded_lsodi. This subroutine has the signature:

      subroutine jac_banded_lsodi_f77
     1  (neq, t, u, s, ml, mu, pd, nrowpd)
Cf2py intent(in, hide) neq, ml, mu, nrowpd
Cf2py intent(out) pd
      integer neq, ml, mu, nrowpd
      double precision t, u, pd, s
      dimension u(neq), s(neq), pd(nrowpd, neq)
      pd = ...
      return
      end
res User-supplied function to calculate the residual vector, defined by r = g(u,t) - A(u,t) * s. Used in the linearly implicit solvers: Lsodi, Lsodis, Lsoibt. The res function has the signature res(u,t,s,ires) and returns the tuple (r,ires), where ires is an int. On input, ires indicates how ODEPACK would use the returned array “r”: ires=1 means the full residual exactly, ires=2 means that r is used only to compute the Jacobian dr/du by finite differences. res should set the flag ires if it encounters a halt condition or illegal input. Otherwise, it should not be reset. On output, the value 1 or -1 represents a normal return.
order Maximum order used by the integrator (<= 12 for “adams”, <= 5 for “bdf”). (default: 4)
ml Number of lower non-zero diagonals in a banded Jacobian.
mu Number of upper non-zero diagonals in a banded Jacobian.
ydoti Real array for the initial value of dy/dt. (default: [])
nsteps Max no of internal solver steps per time step. (default: 1000)
res_f77

Fortran subroutine for res. This subroutine has the signature:

     subroutine res_f77(neq, t, u, s, r, ires)
Cf2py intent(hide) neq
Cf2py intent(out) r
Cf2py intent(in,out) ires
     double precision t, u, s, r
     dimension u(neq, s(neq), r(neq)
     ...
     return
     end

Methods

adjust_atol(u_current) Error tolerance tol(i) may become zero for some i
adjust_parameters()
advance() This function intends to one step forward for linearly implicit solvers
check_conditional_parameters() This function is used to check whether conditional parameters are provided when specified condition fulfilled.
check_extra(**kwargs) A parameter may have a keyword extra_check for user-given functions that performs consistency checks on the parameter.
check_iaja() ia, ja, ic, jc are optional inputs to describe
check_input_range(**kwargs) Check whether all existing inputs are in right specified range.
check_input_types(**kwargs) Check whether all existing inputs are of right specified type.
check_liwlrw() If the lengths of work arrays are specified by users, check whether
check_pars() Some pairs of parameters should be input simutaneously.
check_tol() atol & rtol should be defined as scalars or vectors with length neq.
compile_string_functions(f, **kwargs) Compile functions which are supplied as Fortran strings.
constant_time_step() Check if self.t has a uniform partition.
expand_iwork(new_liw[, expand]) Extension module return an actually required length for integer work array when it is too short.
expand_rwork(new_lrw[, expand]) Length of real work array is smaller than actually required length.
func_wrapper() This function is defined to wrap user-defined functions with new forms of parameter-list, or wrap the returned values as numpy arrays.
get([parameter_name, print_info]) Return value of specified input parameters.
get_parameter_info([print_info]) Return a dictionary containing all properties of all legal parameters in current subclass (i.e., the parameters in self._parameters).
has_u_t_all() Return True if self.u_all and self.t_all, defined in
initialize() Import extension module _odesolver and check that it exists.
initialize_for_solve() In the long parameter-lists for solvers in ODEPACK, quite a few
new_stepnr() When Fortran code returns a status istate==-1, it indicates that
print_roots(jroot, t_current, u_current) Roots found at current T for some constraint functions.
set([strict]) Assign values to one or more parameters, specified as keyword arguments.
set_dummy_functions()
set_extra_args()
set_initial_condition(U0) Function set_initial_condition() is used to set initial value of
set_iopt()
set_iter_method()
set_iwork_rwork()
set_jac()
set_liw_min()
set_lrw_min()
set_ydoti() ydoti is an array used in linearly implicit solvers.
solve(time_points[, terminate])
switch_to(solver_target[, print_info]) Create a new solver instance which switch to another subclass with same values of common attributes.
tol_multiply(tolsf) This function is used to adjust tolerance parameters for Fortran part.
validate_data() Common validity check in Odepack.
__init__(**kwargs)[source]
__module__ = 'odespy.odepack'
adjust_parameters()[source]
quick_description = 'LSODI solver for linearly implicit systems'
set_dummy_functions()[source]
set_extra_args()[source]
set_iopt()[source]
set_iter_method()[source]
set_iwork_rwork()[source]
set_jac()[source]
set_liw_min()[source]
set_lrw_min()[source]
solve(time_points, terminate=None)[source]
class odespy.odepack.Lsodis(**kwargs)[source]

Bases: odespy.odepack.Odepack

A Python wrapper of the LSODIS subroutine from ODEPACK. This subroutine is a variant of LSODI, intended for stiff problems in which the matrix A and the Jacobian of the residual wrt the unknown functions have a sparse structure.

Required input arguments:

Name Description
res User-supplied function to calculate the residual vector, defined by r = g(u,t) - A(u,t) * s. Used in the linearly implicit solvers: Lsodi, Lsodis, Lsoibt. The res function has the signature res(u,t,s,ires) and returns the tuple (r,ires), where ires is an int. On input, ires indicates how ODEPACK would use the returned array “r”: ires=1 means the full residual exactly, ires=2 means that r is used only to compute the Jacobian dr/du by finite differences. res should set the flag ires if it encounters a halt condition or illegal input. Otherwise, it should not be reset. On output, the value 1 or -1 represents a normal return.
adda_lsodis Callable object to add j-th column of matrix A = A(u,t) to another matrix stored in sparse form. The signature is adda(u,t,j,ia,ja,p) and it returns a vector.

Optional input arguments:

Name Description
f_args Extra positional arguments to f: f(u, t, *f_args, **f_kwargs). (default: ())
f_kwargs Extra keyword arguments to f: f(u, t, *f_args, **f_kwargs). (default: {})
complex_valued True if f is complex valued. (default: False)
disk_storage Indicates whether u is stored in memory or in file. If string, it is the filename; if False or “”, u is kept in memory; if True, a default filename tmp_odspy.dat is used. (default: False)
verbose Integer reflecting output of intermediate quantities. (default: 0)
u_exact Function of t returning exact solution. (default: None)
atol Absolute tolerance for solution. (default: 1e-08)
rtol Relative tolerance for solution. (default: 1e-06)
adams_or_bdf Method in vode or solvers in odepack: “adams” or “bdf”. (default: adams)
order Maximum order used by the integrator (<= 12 for “adams”, <= 5 for “bdf”). (default: 4)
nsteps Max no of internal solver steps per time step. (default: 1000)
first_step Suggested first time step size for an adaptive algorithm.
min_step Minimum step size for an adaptive algorithm.
max_step Maximum step size for an adaptive algorithm.
corrector_iter_method Corrector iteration method choice.
lrw Length of real work array.
liw Length of integer work array, similiar as <lrw>.
jac_lsodis Callable object to supply the j-th column of the sparse Jacobian matrix dr/du where r = g - A*s. The signature is jac(u,t,s,j,ia,ja) and a vector is returned.
moss Method to obtain sparse structure of Jacobian.
ia Integer array with length neq+1 which contains starting locations in ja of the descriptions for columns 1...neq. ia(1) == 1. The last element ia[neq+1] should equal to the total number of nonzero locations assumed. For each column j = 1...neq, the values of the row index i in column j, where a nonzero element may occur, are given by i == ja(k) where ia(j) <=’ k < ia(j+1).
ja Integer array containing the row indices where the nonzero elements occur, in columnwise order. Describes the sparsity matrix structure together with ia. In Lsodes, ia and ja describe the structure of Jacobian matrix; while in Lsodis, ia and ja are used to describe the structure of matrix A.
ic Array which contains starting locations in jc.
jc Integer array which describes the sparsity Jacobian structure together with ic, like ia and ja. In Lsodis, ia and ja describe the sparse structure of matrix A, while ic and jc describe the sparse structure of Jacobian matrix.
ydoti Real array for the initial value of dy/dt. (default: [])
order Maximum order used by the integrator (<= 12 for “adams”, <= 5 for “bdf”). (default: 4)

Methods

adjust_atol(u_current) Error tolerance tol(i) may become zero for some i
adjust_parameters()
advance() This function intends to one step forward for linearly implicit solvers
check_conditional_parameters() This function is used to check whether conditional parameters are provided when specified condition fulfilled.
check_extra(**kwargs) A parameter may have a keyword extra_check for user-given functions that performs consistency checks on the parameter.
check_iaja() ia, ja, ic, jc are optional inputs to describe
check_input_range(**kwargs) Check whether all existing inputs are in right specified range.
check_input_types(**kwargs) Check whether all existing inputs are of right specified type.
check_liwlrw() If the lengths of work arrays are specified by users, check whether
check_pars() Some pairs of parameters should be input simutaneously.
check_tol() atol & rtol should be defined as scalars or vectors with length neq.
compile_string_functions(f, **kwargs) Compile functions which are supplied as Fortran strings.
constant_time_step() Check if self.t has a uniform partition.
expand_iwork(new_liw[, expand]) Extension module return an actually required length for integer work array when it is too short.
expand_rwork(new_lrw[, expand]) Length of real work array is smaller than actually required length.
func_wrapper() This function is defined to wrap user-defined functions with new forms of parameter-list, or wrap the returned values as numpy arrays.
get([parameter_name, print_info]) Return value of specified input parameters.
get_parameter_info([print_info]) Return a dictionary containing all properties of all legal parameters in current subclass (i.e., the parameters in self._parameters).
has_u_t_all() Return True if self.u_all and self.t_all, defined in
initialize() Import extension module _odesolver and check that it exists.
initialize_for_solve() In the long parameter-lists for solvers in ODEPACK, quite a few
new_stepnr() When Fortran code returns a status istate==-1, it indicates that
print_roots(jroot, t_current, u_current) Roots found at current T for some constraint functions.
set([strict]) Assign values to one or more parameters, specified as keyword arguments.
set_dummy_functions()
set_extra_args() Setting for extra parameters of user-defined functions.
set_initial_condition(U0) Function set_initial_condition() is used to set initial value of
set_iopt()
set_iter_method()
set_iwork_rwork() Initialization of work arrays with caculated length and optional inputs.
set_jac() Set values for Jacobian matrix. In several solvers like Lsode,
set_liw_min() Calculate the necessary length of integer work arrays when it is not specified explicitly by users.
set_lrw_min() Calculate the necessary length of real work arrays for Fortran code.
set_ydoti() ydoti is an array used in linearly implicit solvers.
solve(time_points[, terminate])
switch_to(solver_target[, print_info]) Create a new solver instance which switch to another subclass with same values of common attributes.
tol_multiply(tolsf) This function is used to adjust tolerance parameters for Fortran part.
validate_data() Common validity check in Odepack.
__init__(**kwargs)[source]
__module__ = 'odespy.odepack'
adjust_parameters()[source]
quick_description = 'LSODIS solver for linearly implicit sparse systems'
set_dummy_functions()[source]
set_iopt()[source]
set_iter_method()[source]
set_iwork_rwork()[source]

Initialization of work arrays with caculated length and optional inputs. In ODEPACK, “iwork” & “rwork” should be initialized with the specific optional parameters in all the solvers. “liw” & “lrw” represented the length requirement of work arrays. Specially, in Lsodis, (ia, ja, ic & jc) should be attached to iwork.

solve(time_points, terminate=None)[source]
class odespy.odepack.Lsoibt(**kwargs)[source]

Bases: odespy.odepack.Odepack

A Python wrapper of the LSOIBT subroutine from ODEPACK. This subroutine is a variant of LSODI for the case where the matrices \(A\), \(dg/du\), and \(d(As)/du\) are all block tridiagonal.

Required input arguments:

Name Description
res User-supplied function to calculate the residual vector, defined by r = g(u,t) - A(u,t) * s. Used in the linearly implicit solvers: Lsodi, Lsodis, Lsoibt. The res function has the signature res(u,t,s,ires) and returns the tuple (r,ires), where ires is an int. On input, ires indicates how ODEPACK would use the returned array “r”: ires=1 means the full residual exactly, ires=2 means that r is used only to compute the Jacobian dr/du by finite differences. res should set the flag ires if it encounters a halt condition or illegal input. Otherwise, it should not be reset. On output, the value 1 or -1 represents a normal return.
adda_lsoibt Callable object to add matrix A = A(u,t) to another matrix P, stored in block-tridiagonal form. The signature is adda(u,t,pa,pb,pc), which should return (pa,pb,pc) as described for the jac_lsoibt parameter.
mb Block size. Describe the block-tridiagonal form of matrix A together with nb.
nb Number of blocks in the main diagonal. In each of the nb block-rows of the matrix P (each consisting of mb consecutive rows), the nonzero elements are to lie in three consecutive mb by mb blocks. In block-rows 2 through nb-1, these are centered about the main diagonal. In block rows 1 and nb, they are the diagonal blocks and the two blocks adjacent to the diagonal block. (Thus block positions (1,3) and (nb,nb-2) can be nonzero.) Require: mb>=1, nb>=4, mb*nb==neq.

Optional input arguments:

Name Description
f_args Extra positional arguments to f: f(u, t, *f_args, **f_kwargs). (default: ())
f_kwargs Extra keyword arguments to f: f(u, t, *f_args, **f_kwargs). (default: {})
complex_valued True if f is complex valued. (default: False)
disk_storage Indicates whether u is stored in memory or in file. If string, it is the filename; if False or “”, u is kept in memory; if True, a default filename tmp_odspy.dat is used. (default: False)
verbose Integer reflecting output of intermediate quantities. (default: 0)
u_exact Function of t returning exact solution. (default: None)
atol Absolute tolerance for solution. (default: 1e-08)
rtol Relative tolerance for solution. (default: 1e-06)
adams_or_bdf Method in vode or solvers in odepack: “adams” or “bdf”. (default: adams)
order Maximum order used by the integrator (<= 12 for “adams”, <= 5 for “bdf”). (default: 4)
nsteps Max no of internal solver steps per time step. (default: 1000)
first_step Suggested first time step size for an adaptive algorithm.
min_step Minimum step size for an adaptive algorithm.
max_step Maximum step size for an adaptive algorithm.
corrector_iter_method Corrector iteration method choice.
lrw Length of real work array.
liw Length of integer work array, similiar as <lrw>.
jac_lsoibt

Callable object to supply the jth column of the Jacobian matrix dr/du where r = g - A*s, stored in block-tridiagonal form. The signature is jac(u,t,s), and the return value is a tuple (pa,pb,pc), where each of these arrays has size mb*mb*nb, mb and nb being two parameters that can be set (see the help line for the parameters in the doc string of Lsoibt). pa, pb, and pc are to be loaded with partial derivatives (elements of the Jacobian matrix) on output, in terms of the block-tridiagonal structure assumed. That is, load the diagonal blocks into pa, the superdiagonal blocks (and block (nb,nb-2)) into pb, and the subdiagonal blocks (and block (1,3)) into pc. The blocks in block-row k of dr/du are to be loaded into pa(,,k), pb(,,k), and pc(,,k). Thus the affect of this function should be:

pa(i,j,k) = ( (i.j) element of k-th diagonal
            block of dr/du)
pb(i,j,k) = ( (i,j) element of block (k,k+1)
            of dr/du, or block (nb,nb-2) if
            k == nb)                 pc(i,j,k) = ( (i,j) element of block (k,k-1)
            of dr/du, or block (1,3) if k==1).
ydoti Real array for the initial value of dy/dt. (default: [])
order Maximum order used by the integrator (<= 12 for “adams”, <= 5 for “bdf”). (default: 4)

Methods

adjust_atol(u_current) Error tolerance tol(i) may become zero for some i
adjust_parameters()
advance() This function intends to one step forward for linearly implicit solvers
check_conditional_parameters() This function is used to check whether conditional parameters are provided when specified condition fulfilled.
check_extra(**kwargs) A parameter may have a keyword extra_check for user-given functions that performs consistency checks on the parameter.
check_iaja() ia, ja, ic, jc are optional inputs to describe
check_input_range(**kwargs) Check whether all existing inputs are in right specified range.
check_input_types(**kwargs) Check whether all existing inputs are of right specified type.
check_liwlrw() If the lengths of work arrays are specified by users, check whether
check_pars() Some pairs of parameters should be input simutaneously.
check_tol() atol & rtol should be defined as scalars or vectors with length neq.
compile_string_functions(f, **kwargs) Compile functions which are supplied as Fortran strings.
constant_time_step() Check if self.t has a uniform partition.
expand_iwork(new_liw[, expand]) Extension module return an actually required length for integer work array when it is too short.
expand_rwork(new_lrw[, expand]) Length of real work array is smaller than actually required length.
func_wrapper() This function is defined to wrap user-defined functions with new forms of parameter-list, or wrap the returned values as numpy arrays.
get([parameter_name, print_info]) Return value of specified input parameters.
get_parameter_info([print_info]) Return a dictionary containing all properties of all legal parameters in current subclass (i.e., the parameters in self._parameters).
has_u_t_all() Return True if self.u_all and self.t_all, defined in
initialize() Import extension module _odesolver and check that it exists.
initialize_for_solve() In the long parameter-lists for solvers in ODEPACK, quite a few
new_stepnr() When Fortran code returns a status istate==-1, it indicates that
print_roots(jroot, t_current, u_current) Roots found at current T for some constraint functions.
set([strict]) Assign values to one or more parameters, specified as keyword arguments.
set_dummy_functions()
set_extra_args() Setting for extra parameters of user-defined functions.
set_initial_condition(U0) Function set_initial_condition() is used to set initial value of
set_iopt()
set_iter_method()
set_iwork_rwork()
set_jac() Set values for Jacobian matrix. In several solvers like Lsode,
set_liw_min()
set_lrw_min()
set_ydoti() ydoti is an array used in linearly implicit solvers.
solve(time_points[, terminate])
switch_to(solver_target[, print_info]) Create a new solver instance which switch to another subclass with same values of common attributes.
tol_multiply(tolsf) This function is used to adjust tolerance parameters for Fortran part.
validate_data()
__init__(**kwargs)[source]
__module__ = 'odespy.odepack'
adjust_parameters()[source]
quick_description = 'LSOIBIT solver for linearly implicit block tridiag systems'
set_dummy_functions()[source]
set_iopt()[source]
set_iter_method()[source]
set_iwork_rwork()[source]
set_liw_min()[source]
set_lrw_min()[source]
solve(time_points, terminate=None)[source]
validate_data()[source]
class odespy.odepack.Odepack(f, **kwargs)[source]

Bases: odespy.solvers.Solver

This is a superclass for wrapping for seven solvers in the Fortran package ODEPACK (available at the netlib repository: www.netlib.org/odepack).

Solvers for explicitly given systems. For each of the following solvers, it is assumed that the ODEs are given explicitly, so that the system can be written in the form du/dt = f(u,t), where u is a vector of dependent variables, and t is a scalar.

Name Description
Lsode

A wrapper of dlsode, the basic solver in ODEPACK for stiff and nonstiff systems of the form u’ = f.

In the stiff case, it treats the Jacobian matrix df/du as either a dense (full) or a banded matrix, and as either user-supplied or internally approximated by differences.

It uses Adams methods (predictor-corrector) in the nonstiff case, and Backward Differentiation Formula (BDF) methods (the Gear methods) in the stiff case. The linear systems that arise are solved by direct methods (LU factorization/backsolve).

Lsodes Solves systems u’ = f, and in the stiff case treats Jacobian matrix in general sparse form. It can determine the sparsity structure on its own, or optionally accepts this information from the user. It then uses parts of the Yale Sparse Matrix Package (YSMP) to solve the linear systems that arise, by a sparse (direct) LU factorization/backsolve method.
Lsoda Solves systems u’ = f, with a dense or banded Jacobian when the problem is stiff, but it automatically selects between nonstiff (Adams) and stiff (BDF) methods. It uses the nonstiff method initially, and dynamically monitors data in order to decide which method to use.
Lsodar A variant of Lsoda allowing for constraints. It solves u’ = with dense or banded Jacobian and automatic method selection, and at the same time, it solves g(u,t) = 0. This is often useful for finding stopping conditions, or for finding points at which a switch is to be made in the function f.

Solvers for linearly implicit systems. The following solvers treat systems in the linearly implicit form:

A(u,t) du/dt = g(u,t),

where A is a square matrix, i.e., with the derivative du/dt implicit, but linearly so. These solvers allow A to be singular, in which case the system is a differential-algebraic equation (DAE) system. In that case, the user must be very careful to supply a well-posed problem with consistent initial conditions.

Name Description
Lsodi Solves linearly implicit systems in which the matrices involved (A, dg/du, and d(A u’)/du) are all assumed to be either dense or banded.
Lsodibt Solves linearly implicit systems in which the matrices involved are all assumed to be block-tridiagonal. Linear systems are solved by the LU method.
Lsodis Solves linearly implicit systems in which the matrices involved are all assumed to be sparse. Either determines the sparsity structure or accepts it from the user, and uses parts of the Yale Sparse Matrix Package to solve the linear systems that arise, by a direct method.

Note: For large ODE systems the user is encouraged that users provide an f2py-compiled Fortran subroutine or a multi-line string Fortran code to define the ODE. This would help to improve efficiency.

Required input arguments:

Name Description
f Right-hand side f(u,t) defining the ODE.

Optional input arguments:

Name Description
f_args Extra positional arguments to f: f(u, t, *f_args, **f_kwargs). (default: ())
f_kwargs Extra keyword arguments to f: f(u, t, *f_args, **f_kwargs). (default: {})
complex_valued True if f is complex valued. (default: False)
disk_storage Indicates whether u is stored in memory or in file. If string, it is the filename; if False or “”, u is kept in memory; if True, a default filename tmp_odspy.dat is used. (default: False)
verbose Integer reflecting output of intermediate quantities. (default: 0)
u_exact Function of t returning exact solution. (default: None)
atol Absolute tolerance for solution. (default: 1e-08)
rtol Relative tolerance for solution. (default: 1e-06)
adams_or_bdf Method in vode or solvers in odepack: “adams” or “bdf”. (default: adams)
order Maximum order used by the integrator (<= 12 for “adams”, <= 5 for “bdf”). (default: 4)
nsteps Max no of internal solver steps per time step. (default: 1000)
first_step Suggested first time step size for an adaptive algorithm.
min_step Minimum step size for an adaptive algorithm.
max_step Maximum step size for an adaptive algorithm.
corrector_iter_method Corrector iteration method choice.
lrw Length of real work array.
liw Length of integer work array, similiar as <lrw>.

Methods

adjust_atol(u_current) Error tolerance tol(i) may become zero for some i
adjust_parameters() Special settings for properties of input parameters.
advance() This function intends to one step forward for linearly implicit solvers
check_conditional_parameters() This function is used to check whether conditional parameters are provided when specified condition fulfilled.
check_extra(**kwargs) A parameter may have a keyword extra_check for user-given functions that performs consistency checks on the parameter.
check_iaja() ia, ja, ic, jc are optional inputs to describe
check_input_range(**kwargs) Check whether all existing inputs are in right specified range.
check_input_types(**kwargs) Check whether all existing inputs are of right specified type.
check_liwlrw() If the lengths of work arrays are specified by users, check whether
check_pars() Some pairs of parameters should be input simutaneously.
check_tol() atol & rtol should be defined as scalars or vectors with length neq.
compile_string_functions(f, **kwargs) Compile functions which are supplied as Fortran strings.
constant_time_step() Check if self.t has a uniform partition.
expand_iwork(new_liw[, expand]) Extension module return an actually required length for integer work array when it is too short.
expand_rwork(new_lrw[, expand]) Length of real work array is smaller than actually required length.
func_wrapper() This function is defined to wrap user-defined functions with new forms of parameter-list, or wrap the returned values as numpy arrays.
get([parameter_name, print_info]) Return value of specified input parameters.
get_parameter_info([print_info]) Return a dictionary containing all properties of all legal parameters in current subclass (i.e., the parameters in self._parameters).
has_u_t_all() Return True if self.u_all and self.t_all, defined in
initialize() Import extension module _odesolver and check that it exists.
initialize_for_solve() In the long parameter-lists for solvers in ODEPACK, quite a few
new_stepnr() When Fortran code returns a status istate==-1, it indicates that
print_roots(jroot, t_current, u_current) Roots found at current T for some constraint functions.
set([strict]) Assign values to one or more parameters, specified as keyword arguments.
set_dummy_functions() Functions have to get dummy values before they are passed to extension
set_extra_args() Setting for extra parameters of user-defined functions.
set_initial_condition(U0) Function set_initial_condition() is used to set initial value of
set_iopt() Initialization for “ipot”, which is a flag to indicate whether optional
set_iter_method() Set proper values for method-choices when it is not specified
set_iwork_rwork() Initialize arrays for optional inputs, and calculate the
set_jac() Set values for Jacobian matrix. In several solvers like Lsode,
set_liw_min() Calculate the necessary length of integer work arrays when it is not specified explicitly by users.
set_lrw_min() Calculate the necessary length of real work arrays for Fortran code.
set_ydoti() ydoti is an array used in linearly implicit solvers.
solve(time_points[, terminate]) This function is involved for non-linearly implicit
switch_to(solver_target[, print_info]) Create a new solver instance which switch to another subclass with same values of common attributes.
tol_multiply(tolsf) This function is used to adjust tolerance parameters for Fortran part.
validate_data() Common validity check in Odepack.
__module__ = 'odespy.odepack'
adjust_atol(u_current)[source]

Error tolerance tol(i) may become zero for some i during integration, where:

tol = rtol(i) * u(i) + atol(i)

It indicates that pure absolute tolerance (atol(i)=0.0) was requested. In order to avoid possible divide-by-zero error, we find the indices of zero items and adjust atol to 1e-8.

adjust_parameters()[source]

Special settings for properties of input parameters.

advance()[source]

This function intends to one step forward for linearly implicit solvers (Lsodi, Lsodis, Lsoibt) in ODEPACK.

For these linearly implicit solvers, if extra wrappers are added in Fortran code, there are often memory errors. Besides, sometimes there are unavoidable errors caused by bugs in the Ubuntu/Linux libraries as libc. To make these solvers more reliable on all platforms, this function is used to call solvers in ODEPACK (dlsodi, dlsodis, dlsoibt) directly without any wrappers in Fortran. However, this would lead to efficiency lost with long work arrays as input parameters for Fortran code. In Lsodi, Lsodis and Lsoibt, Solver.solve() would be applied to get the desired solution, which will direct to this function to step forward.

check_iaja()[source]

ia, ja, ic, jc are optional inputs to describe arbitrary sparse structure of matrix.

ia and ja are used in dlsodes, dlsodis. ic and jc are used only in dlsodis.

There are special requirements for their values:

len(ia/ic) = neq + 1
(ia/ic)[0] = 1
(ia/ic)[-1] = 1 + len(ja/jc)
check_liwlrw()[source]

If the lengths of work arrays are specified by users, check whether they are greater than the required lengths of Fortran solvers.

check_pars()[source]

Some pairs of parameters should be input simutaneously.

check_tol()[source]

atol & rtol should be defined as scalars or vectors with length neq.

expand_iwork(new_liw, expand=False)[source]

Extension module return an actually required length for integer work array when it is too short. Then we could expand work array to required length to avoid this error.

expand_rwork(new_lrw, expand=False)[source]

Length of real work array is smaller than actually required length. Then we could expand work array to avoid this error.

func_wrapper()[source]

This function is defined to wrap user-defined functions with new forms of parameter-list, or wrap the returned values as numpy arrays.

Firstly, in odespy, all the user-supplied functions should have a parameter list starts with “u,t,...”. But in some special subclasses, (like solvers in ODEPACK), all the parameter lists of user-defined functions start with “t,u,...”. So we need this general function to wrap all these user-defined functions.

Secondly, in some user-defined functions, according to the different start indices in Fortran and Python, we need to make special wrapping for these uncompability. For an example, in user-defined function “jac_column”, column index is an internally valued parameter in Fortran code. In Python, it starts from 0 instead of 1 in Fortran. So we need to wrap the parameter list of user-defined “jac_column” from “u,t,j” to “t,u,j+1”. That is, define the Jacobian function as lambda t,u,j: jac_column(u,t,j-1).

Furthermore, the return value of user-defined functions need to be wrapped to Numpy arrays with great numerical features, e.g. vectorization and array slicing. In order to avoid unnecessary array copy by F2PY, it is always recommended to explicitly transform all Numpy arrays to Fortran order in Python code.

This functions is not intended for simple solvers. So it is not called automatically in current version. But for some complicated solvers as ones in ODEPACK, it is very useful and convenient.

Future developers can call this functions with appropriate locations and corresponding property-setting in adjust_parameters().

initialize()[source]

Import extension module _odesolver and check that it exists.

initialize_for_solve()[source]

In the long parameter-lists for solvers in ODEPACK, quite a few parameters can be set automatically with proper values depending on values of other parameters. In this function, all the parameters of this kind are initialized with proper values.

new_stepnr()[source]

When Fortran code returns a status istate==-1, it indicates that there are excessive amount of steps detected. Then we could try to increase nsteps to avoid this error.

print_roots(jroot, t_current, u_current)[source]

Roots found at current T for some constraint functions.

set_dummy_functions()[source]

Functions have to get dummy values before they are passed to extension module even if they are not involved in current solver.

set_extra_args()[source]

Setting for extra parameters of user-defined functions.

set_iopt()[source]

Initialization for “ipot”, which is a flag to indicate whether optional parameters are specified in work arrays.

set_iter_method()[source]

Set proper values for method-choices when it is not specified explicitly.

set_iwork_rwork()[source]

Initialize arrays for optional inputs, and calculate the required length of work arrays in Fortran code.

set_jac()[source]

Set values for Jacobian matrix. In several solvers like Lsode, Jacobian matrix could be supplied either in full form or in banded form.

This function intends to tell from which kind of Jacobian matrix is specified.

set_liw_min()[source]

Calculate the necessary length of integer work arrays when it is not specified explicitly by users. Different solvers have different formulas.

set_lrw_min()[source]

Calculate the necessary length of real work arrays for Fortran code. Different solvers have different formulas.

set_ydoti()[source]

ydoti is an array used in linearly implicit solvers. It has to be extended if its length is smaller than neq.

solve(time_points, terminate=None)[source]

This function is involved for non-linearly implicit solvers in ODEPACK, i.e., Lsode, Lsoda, Lsodar, and Lsodes.

tol_multiply(tolsf)[source]

This function is used to adjust tolerance parameters for Fortran part. When extension module returns a status “istate” as -2 or -3, it often indicates that there are excessive amount of steps detected. Then we could try to adjust tolerance settings with suggested factor to avoid this error.

validate_data()[source]

Common validity check in Odepack.

Previous topic

odespy.version