pyodesys package

Straightforward numerical integration of ODE systems from Python.

The pyodesys package for enables intuitive solving of systems of first order differential equations numerically. The system may be symbolically defined.

pyodesys ties computer algebra systems (like SymPy and symengine), and numerical solvers (such as ODEPACK in SciPy, CVode in sundials, GSL or odeint) together.

Submodules

pyodesys.analysis module

pyodesys.convergence module

pyodesys.convergence.fit_factory(discard=1)[source]
pyodesys.convergence.integrate_tolerance_series(odesys, atols, rtols, x, y0, params=(), fit=<function <lambda> at 0x7f4c74e06dc0>, val=<function polyval at 0x7f4c72ee61f0>, \*\*kwargs)[source]
Parameters
odesysODESys
atolsarray_like

Positive, monotonically increasing 1D array.

rtolsarray_like

Positive, monotonically increasing 1D array.

xarray_like

Passed on to odesys.integrate for first set of tolerances. (subsequent calls will use xout from first integration).

y0array_like

Passed on to odesys.integrate.

paramsarray_like

Passed on to odesys.integrate.

fitcallable
valcallable
**kwargs:

Passed on to odesys.integrate.

Returns
result0Result
resultslist of Result instances
extradict

errest : 2D array of error estimates for result0.yout

pyodesys.core module

Core functionality for ODESys.

Note that it is possible to use new custom ODE integrators with pyodesys by providing a module with two functions named integrate_adaptive and integrate_predefined. See the pyodesys.integrators module for examples.

class pyodesys.core.ODESys(f, jac=None, dfdx=None, jtimes=None, first_step_cb=None, roots_cb=None, nroots=None, band=None, names=, param_names=, indep_name=None, description=None, dep_by_name=False, par_by_name=False, latex_names=, latex_param_names=, latex_indep_name=None, taken_names=None, pre_processors=None, post_processors=None, append_iv=False, autonomous_interface=None, to_arrays_callbacks=None, autonomous_exprs=None, _indep_autonomous_key=None, numpy=None, nnz=- 1, **kwargs)[source]

Bases: object

Object representing an ODE system.

ODESys provides unified interface to:

  • scipy.integarte.ode

  • pygslodeiv2

  • pyodeint

  • pycvodes

The numerical integration can be performed either in an adaptive() or predefined() mode. Where locations to report the solution is chosen by the stepper or the user respectively. For convenience in user code one may use integrate() which automatically chooses between the two based on the length of xout provided by the user.

Parameters
fcallback

first derivatives of dependent variables (y) with respect to dependent variable (x). Signature is any of:

  • rhs(x, y[:]) -> f[:]

  • rhs(x, y[:], p[:]) -> f[:]

  • rhs(x, y[:], p[:], backend=math) -> f[:].

jaccallback

Jacobian matrix (dfdy). Required for implicit methods. Signature should be either of:

  • jac(x, y[:]) -> J

  • jac(x, y[:], p[:]) -J.

If nnz < 0, J should be a 2D array-like object if nnz < 0 corresponding to a dense or banded jacobian (see also band). If nnz >= 0, J should be an instance of scipy.sparse.csc_matrix.

dfdxcallback

Signature dfdx(x, y[:], p[:]) -> out[:] (used by e.g. GSL)

jtimescallback

Jacobian-vector product (Jv). Signature is `jtimes(x, y[:], v[:]) -> Jv[:]` This is supported only by cvode.

first_step_cbcallback

Signature step1st(x, y[:], p[:]) -> dx0 (pass first_step==0 to use). This is available for cvode, odeint & gsl, but not for scipy.

roots_cbcallback

Signature roots_cb(x, y[:], p[:]=(), backend=math) -> discr[:].

nrootsint

Length of return vector from roots_cb.

bandtuple of 2 integers or None (default: None)

If jacobian is banded: number of sub- and super-diagonals

namesiterable of strings (default

Names of variables, used for referencing dependent variables by name and for labels in plots.

param_namesiterable of strings (default: None)

Names of the parameters, used for referencing parameters by name.

indep_namestr

Name of the independent variable

dep_by_namebool

When True integrate() expects a dictionary as input for y0.

par_by_namebool

When True integrate() expects a dictionary as input for params.

latex_namesiterable of strings (default

Names of variables in LaTeX format (e.g. for labels in plots).

latex_param_namesiterable of strings (default

Names of parameters in LaTeX format (e.g. for labels in plots).

latex_indep_namestr

LaTeX formatted name of independent variable.

taken_namesiterable of str

Names of dependent variables which are calculated in pre_processors

pre_processorsiterable of callables (optional)

signature: f(x1[:], y1[:], params1[:]) -> x2[:], y2[:], params2[:]. When modifying: insert at beginning.

post_processorsiterable of callables (optional)

signature: f(x2[:], y2[:, :], params2[:]) -> x1[:], y1[:, :], params1[:] When modifying: insert at end.

append_ivbool

See append_iv.

autonomous_interfacebool (optional)

If given, sets the autonomous_interface to indicate whether the system appears autonomous or not upon call to integrate().

autonomous_exprsbool

Describes whether the independent variable appears in the rhs expressions. If set to True the underlying solver is allowed to shift the independent variable during integration.

nnzint (default

Maximum number of non-zero entries in sparse jacobian. When non-negative, jacobian is assumed to be dense or banded.

Notes

Banded jacobians are supported by “scipy” and “cvode” integrators.

Examples

>>> odesys = ODESys(lambda x, y, p: p[0]*x + p[1]*y[0]*y[0])
>>> yout, info = odesys.predefined([1], [0, .2, .5], [2, 1])
>>> print(info['success'])
True
Attributes
f_cbcallback

For evaluating the vector of derivatives.

j_cbcallback or None

For evaluating the Jacobian matrix of f.

dfdx_cbcallback or None

For evaluating the second order derivatives.

jtimes_cbcallback or None

For evaluating the jacobian-vector product.

first_step_cbcallback or None

For calculating the first step based on x0, y0 & p.

roots_cbcallback
nrootsint
nnzint
namestuple of strings
param_namestuple of strings
descriptionstr
dep_by_namebool
par_by_namebool
latex_namestuple of str
latex_param_namestuple of str
pre_processorsiterable of callbacks
post_processorsiterable of callbacks
append_ivbool

If True params[:] passed to f_cb, jac_cb will contain initial values of y. Note that this happens after pre processors have been applied.

autonomous_interfacebool or None

Indicates whether the system appears autonomous upon call to integrate(). None indicates that it is unknown.

Methods

adaptive(self, y0, x0, xend[, params])

Integrate with integrator chosen output.

chained_parameter_variation(self, \*args, …)

See chained_parameter_variation().

integrate(self, x, y0[, params, atol, rtol])

Integrate the system of ordinary differential equations.

plot_phase_plane(self[, indices])

Plots a phase portrait from last integration.

plot_result(self, \*\*kwargs)

Plots the integrated dependent variables from last integration.

post_process(self, xout, yout, params)

Transforms internal values to output, used internally.

pre_process(self, xout, y0[, params])

Transforms input to internal values, used internally.

predefined(self, y0, xout[, params])

Integrate with user chosen output.

stiffness(self[, xyp, eigenvals_cb])

[DEPRECATED] Use Result.stiffness(), stiffness ration

to_arrays

adaptive(self, y0, x0, xend, params=(), \*\*kwargs)[source]

Integrate with integrator chosen output.

Parameters
integratorstr

See integrate().

y0array_like

See integrate().

x0float

Initial value of the independent variable.

xendfloat

Final value of the independent variable.

paramsarray_like

See integrate().

**kwargs :

See integrate().

Returns
Same as integrate()
chained_parameter_variation(self, \*args, \*\*kwargs)[source]

See chained_parameter_variation().

integrate(self, x, y0, params=(), atol=1e-08, rtol=1e-08, \*\*kwargs)[source]

Integrate the system of ordinary differential equations.

Solves the initial value problem (IVP).

Parameters
xarray_like or pair (start and final time) or float
if float:

make it a pair: (0, x)

if pair or length-2 array:

initial and final value of the independent variable

if array_like:

values of independent variable report at

y0array_like

Initial values at x[0] for the dependent variables.

paramsarray_like (default: tuple())

Value of parameters passed to user-supplied callbacks.

integratorstr or None
Name of integrator, one of:
  • ‘scipy’: _integrate_scipy()

  • ‘gsl’: _integrate_gsl()

  • ‘odeint’: _integrate_odeint()

  • ‘cvode’: _integrate_cvode()

See respective method for more information. If None: os.environ.get('PYODESYS_INTEGRATOR', 'scipy')

atolfloat

Absolute tolerance

rtolfloat

Relative tolerance

with_jacobianbool or None (default)

Whether to use the jacobian. When None the choice is done automatically (only used when required). This matters when jacobian is derived at runtime (high computational cost).

with_jtimesbool (default: False)

Whether to use the jacobian-vector product. This is only supported by cvode and only when linear_solver is one of: gmres’, ‘gmres_classic’, ‘bicgstab’, ‘tfqmr’. See the documentation for pycvodes for more information.

force_predefinedbool (default: False)

override behaviour of len(x) == 2 => adaptive()

**kwargs :

Additional keyword arguments for _integrate_$(integrator).

Returns
Length 3 tuple: (x, yout, info)

x : array of values of the independent variable yout : array of the dependent variable(s) for the different

values of x.

info : dict (‘nfev’ is guaranteed to be a key)

plot_phase_plane(self, indices=None, \*\*kwargs)[source]

Plots a phase portrait from last integration.

This method will be deprecated. Please use Result.plot_phase_plane(). See pyodesys.plotting.plot_phase_plane()

plot_result(self, \*\*kwargs)[source]

Plots the integrated dependent variables from last integration.

This method will be deprecated. Please use Result.plot(). See pyodesys.plotting.plot_result()

post_process(self, xout, yout, params)[source]

Transforms internal values to output, used internally.

pre_process(self, xout, y0, params=)[source]

Transforms input to internal values, used internally.

predefined(self, y0, xout, params=(), \*\*kwargs)[source]

Integrate with user chosen output.

Parameters
integratorstr

See integrate().

y0array_like

See integrate().

xoutarray_like
paramsarray_like

See integrate().

**kwargs:

See integrate()

Returns
Length 2 tuple(yout, info)

See integrate().

stiffness(self, xyp=None, eigenvals_cb=None)[source]

[DEPRECATED] Use Result.stiffness(), stiffness ration

Running stiffness ratio from last integration. Calculate sittness ratio, i.e. the ratio between the largest and smallest absolute eigenvalue of the jacobian matrix. The user may supply their own routine for calculating the eigenvalues, or they will be calculated from the SVD (singular value decomposition). Note that calculating the SVD for any but the smallest Jacobians may prove to be prohibitively expensive.

Parameters
xyplength 3 tuple (default: None)

internal_xout, internal_yout, internal_params, taken from last integration if not specified.

eigenvals_cbcallback (optional)

Signature (x, y, p) (internal variables), when not provided an internal routine will use self.j_cb and scipy.linalg.svd.

to_arrays(self, x, y, p, callbacks=None, reshape=True)[source]
class pyodesys.core.OdeSys(f, jac=None, dfdx=None, jtimes=None, first_step_cb=None, roots_cb=None, nroots=None, band=None, names=, param_names=, indep_name=None, description=None, dep_by_name=False, par_by_name=False, latex_names=, latex_param_names=, latex_indep_name=None, taken_names=None, pre_processors=None, post_processors=None, append_iv=False, autonomous_interface=None, to_arrays_callbacks=None, autonomous_exprs=None, _indep_autonomous_key=None, numpy=None, nnz=- 1, **kwargs)[source]

Bases: pyodesys.core.ODESys

DEPRECATED, use ODESys instead.

Methods

adaptive(self, y0, x0, xend[, params])

Integrate with integrator chosen output.

chained_parameter_variation(self, \*args, …)

See chained_parameter_variation().

integrate(self, x, y0[, params, atol, rtol])

Integrate the system of ordinary differential equations.

plot_phase_plane(self[, indices])

Plots a phase portrait from last integration.

plot_result(self, \*\*kwargs)

Plots the integrated dependent variables from last integration.

post_process(self, xout, yout, params)

Transforms internal values to output, used internally.

pre_process(self, xout, y0[, params])

Transforms input to internal values, used internally.

predefined(self, y0, xout[, params])

Integrate with user chosen output.

stiffness(self[, xyp, eigenvals_cb])

[DEPRECATED] Use Result.stiffness(), stiffness ration

to_arrays

exception pyodesys.core.RecoverableError[source]

Bases: Exception

pyodesys.core.chained_parameter_variation(subject, durations, y0, varied_params, default_params=None, integrate_kwargs=None, x0=None, npoints=1, numpy=None)[source]

Integrate an ODE-system for a serie of durations with some parameters changed in-between

Parameters
subjectfunction or ODESys instance

If a function: should have the signature of pyodesys.ODESys.integrate() (and resturn a pyodesys.results.Result object). If a ODESys instance: the integrate method will be used.

durationsiterable of floats

Spans of the independent variable.

y0dict or array_like
varied_paramsdict mapping parameter name (or index) to array_like

Each array_like need to be of same length as durations.

default_paramsdict or array_like

Default values for the parameters of the ODE system.

integrate_kwargsdict

Keyword arguments passed on to integrate.

x0float-like

First value of independent variable. default: 0.

npointsint

Number of points per sub-interval.

Examples

>>> odesys = ODESys(lambda t, y, p: [-p[0]*y[0]])
>>> int_kw = dict(integrator='cvode', method='adams', atol=1e-12, rtol=1e-12)
>>> kwargs = dict(default_params=[0], integrate_kwargs=int_kw)
>>> res = chained_parameter_variation(odesys, [2, 3], [42], {0: [.7, .1]}, **kwargs)
>>> mask1 = res.xout <= 2
>>> import numpy as np
>>> np.allclose(res.yout[mask1, 0], 42*np.exp(-.7*res.xout[mask1]))
True
>>> mask2 = 2 <= res.xout
>>> np.allclose(res.yout[mask2, 0], res.yout[mask2, 0][0]*np.exp(-.1*(res.xout[mask2] - res.xout[mask2][0])))
True
pyodesys.core.integrate_auto_switch(odes, kw, x, y0, params=(), \*\*kwargs)[source]

Auto-switching between formulations of ODE system.

In case one has a formulation of a system of ODEs which is preferential in the beginning of the integration, this function allows the user to run the integration with this system where it takes a user-specified maximum number of steps before switching to another formulation (unless final value of the independent variables has been reached). Number of systems used i returned as nsys in info dict.

Parameters
odesiterable of OdeSy instances
kwdict mapping kwarg to iterables of same legnth as odes
xarray_like
y0array_like
paramsarray_like
**kwargs:

See ODESys.integrate()

Notes

Plays particularly well with symbolic.TransformedSys.

pyodesys.core.integrate_chained(odes, kw, x, y0, params=(), \*\*kwargs)

Auto-switching between formulations of ODE system.

In case one has a formulation of a system of ODEs which is preferential in the beginning of the integration, this function allows the user to run the integration with this system where it takes a user-specified maximum number of steps before switching to another formulation (unless final value of the independent variables has been reached). Number of systems used i returned as nsys in info dict.

Parameters
odesiterable of OdeSy instances
kwdict mapping kwarg to iterables of same legnth as odes
xarray_like
y0array_like
paramsarray_like
**kwargs:

See ODESys.integrate()

Notes

Plays particularly well with symbolic.TransformedSys.

pyodesys.integrators module

This module is for demonstration purposes only and the integrators here are not meant for production use. Consider them provisional, i.e., API here may break without prior deprecation.

class pyodesys.integrators.EulerBackward_example_integrator[source]

Bases: object

Attributes
integrate_adaptive

Methods

integrate_predefined

integrate_adaptive = None
static integrate_predefined(rhs, jac, y0, xout, \*\*kwargs)[source]
with_jacobian = True
class pyodesys.integrators.EulerForward_example_integrator[source]

Bases: object

Attributes
integrate_adaptive

Methods

integrate_predefined

integrate_adaptive = None
static integrate_predefined(rhs, jac, y0, xout, \*\*kwargs)[source]
with_jacobian = False
class pyodesys.integrators.Midpoint_example_integrator[source]

Bases: object

Attributes
integrate_adaptive

Methods

integrate_predefined

integrate_adaptive = None
static integrate_predefined(rhs, jac, y0, xout, \*\*kwargs)[source]
with_jacobian = False
class pyodesys.integrators.RK4_example_integrator[source]

Bases: object

This is an example of how to implement a custom integrator. It uses fixed step size and is usually not useful for real problems.

Methods

integrate_adaptive

integrate_predefined

static integrate_adaptive(rhs, jac, y0, x0, xend, dx0, \*\*kwargs)[source]
static integrate_predefined(rhs, jac, y0, xout, \*\*kwargs)[source]
with_jacobian = False
class pyodesys.integrators.Trapezoidal_example_integrator[source]

Bases: object

Attributes
integrate_adaptive

Methods

integrate_predefined

integrate_adaptive = None
static integrate_predefined(rhs, jac, y0, xout, \*\*kwargs)[source]
with_jacobian = True

pyodesys.plotting module

pyodesys.plotting.info_vlines(ax, xout, info, vline_colors='maroon', 'purple', vline_keys='steps', 'rhs_xvals', 'jac_xvals', post_proc=None, alpha=None, fpes=None, every=None)[source]

Plot vertical lines in the background

Parameters
axaxes
xoutarray_like
infodict
vline_colorsiterable of str
vline_keysiterable of str

Choose from 'steps', 'rhs_xvals', 'jac_xvals', 'fe_underflow', 'fe_overflow', 'fe_invalid', 'fe_divbyzero'.

vline_post_proccallable
alphafloat
pyodesys.plotting.plot_phase_plane(x, y, indices=None, plot=None, names=None, ax=None, \*\*kwargs)[source]

Plot the phase portrait of two dependent variables

Parameters
x: array_like

Values of the independent variable.

y: array_like

Values of the dependent variables.

indices: pair of integers (default: None)

What dependent variable to plot for (None => (0, 1)).

plot: callable (default: None)

Uses matplotlib.pyplot.plot if None

names: iterable of strings

Labels for x and y axis.

**kwargs:

Keyword arguemtns passed to plot().

pyodesys.plotting.plot_result(x, y, indices=None, plot_kwargs_cb=None, ax=None, ls='-', '--', ':', '-.', c='tab:blue', 'tab:orange', 'tab:green', 'tab:red', 'tab:purple', 'tab:brown', 'tab:pink', 'tab:gray', 'tab:olive', 'tab:cyan', 'black', m='o', 'v', '8', 's', 'p', 'x', '+', 'd', 's', m_lim=- 1, lines=None, interpolate=None, interp_from_deriv=None, names=None, latex_names=None, xlabel=None, ylabel=None, xscale=None, yscale=None, legend=False, yerr=None, labels=None, tex_lbl_fmt='$%s$', fig_kw=None, xlim=None, ylim=None)[source]

Plot the depepndent variables vs. the independent variable

Parameters
xarray_like

Values of the independent variable.

yarray_like

Values of the independent variable. This must hold y.shape[0] == len(x), plot_results will draw y.shape[1] lines. If interpolate != None y is expected two be three dimensional, otherwise two dimensional.

indicesiterable of integers

What indices to plot (default: None => all).

plotcallback (default: None)

If None, use matplotlib.pyplot.plot.

plot_kwargs_cbcallback(int) -> dict

Keyword arguments for plot for each index (0:len(y)-1).

axAxes
lsiterable

Linestyles to cycle through (only used if plot and plot_kwargs_cb are both None).

citerable

Colors to cycle through (only used if plot and plot_kwargs_cb are both None).

miterable

Markers to cycle through (only used if plot and plot_kwargs_cb are both None and m_lim > 0).

m_limint (default: -1)

Upper limit (exclusive, number of points) for using markers instead of lines.

linesNone

default: draw between markers unless we are interpolating as well.

interpolatebool or int (default: None)

Density-multiplier for grid of independent variable when interpolating if True => 20. negative integer signifies log-spaced grid.

interp_from_derivcallback (default: None)

When None: scipy.interpolate.BPoly.from_derivatives

namesiterable of str
latex_namesiterable of str
labelsiterable of str

If None, use latex_names or names (in that order).

pyodesys.plotting.right_hand_ylabels(ax, labels)[source]

pyodesys.results module

class pyodesys.results.Result(xout, yout, params, info, odesys)[source]

Bases: object

Methods

at(self, x[, use_deriv, xdata, ydata, linear])

Returns interpolated result at a given time and an interpolation error-estimate

between(self, lower, upper[, xdata, ydata])

Get results inside span for independent variable

plot(self[, info_vlines_kw, between, deriv, …])

Plots the integrated dependent variables from last integration.

plot_phase_plane(self[, indices])

Plots a phase portrait from last integration.

stiffness(self[, xyp, eigenvals_cb])

Running stiffness ratio from last integration.

calc_invariant_violations

copy

extend_by_integration

named_dep

named_param

plot_invariant_violations

at(self, x, use_deriv=False, xdata=None, ydata=None, linear=False)[source]

Returns interpolated result at a given time and an interpolation error-estimate

By default interpolation is performed using cubic splines.

Parameters
xarray_like or float
use_derivbool

Calculate derivatives at spline knots for enhanced accuracy.

xdataarray
ydataarray
linearbool

Will use (cheaper) linear interpolation. Useful when x is an array. Error estimate will be None in this case.

Returns
interpolated_yarray
error_estim_yarray
between(self, lower, upper, xdata=None, ydata=None)[source]

Get results inside span for independent variable

calc_invariant_violations(self, xyp=None)[source]
copy(self)[source]
extend_by_integration(self, xend, params=None, odesys=None, autonomous=None, npoints=1, \*\*kwargs)[source]
named_dep(self, name)[source]
named_param(self, param_name)[source]
plot(self, info_vlines_kw=None, between=None, deriv=False, title_info=0, \*\*kwargs)[source]

Plots the integrated dependent variables from last integration.

Parameters
info_vlines_kwdict

Keyword arguments passed to plotting.info_vlines(), an empty dict will be used if True. Need to pass ax when given.

indicesiterable of int
betweenlength 2 tuple
derivbool

Plot derivatives (internal variables).

namesiterable of str
**kwargs:

See pyodesys.plotting.plot_result()

plot_invariant_violations(self, \*\*kwargs)[source]
plot_phase_plane(self, indices=None, \*\*kwargs)[source]

Plots a phase portrait from last integration.

Parameters
indicesiterable of int
namesiterable of str
**kwargs:

See pyodesys.plotting.plot_phase_plane()

stiffness(self, xyp=None, eigenvals_cb=None)[source]

Running stiffness ratio from last integration.

Calculate sittness ratio, i.e. the ratio between the largest and smallest absolute eigenvalue of the jacobian matrix. The user may supply their own routine for calculating the eigenvalues, or they will be calculated from the SVD (singular value decomposition). Note that calculating the SVD for any but the smallest Jacobians may prove to be prohibitively expensive.

Parameters
xyplength 3 tuple (default: None)

internal_xout, internal_yout, internal_params, taken from last integration if not specified.

eigenvals_cbcallback (optional)

Signature (x, y, p) (internal variables), when not provided an internal routine will use self.j_cb and scipy.linalg.svd.

pyodesys.symbolic module

This module contains a subclass of ODESys which allows the user to generate auxiliary expressions from a canonical set of symbolic expressions. Subclasses are also provided for dealing with variable transformations and partially solved systems.

class pyodesys.symbolic.PartiallySolvedSystem(original_system, analytic_factory, **kwargs)[source]

Bases: pyodesys.symbolic.SymbolicSys

Use analytic expressions for some dependent variables

Parameters
original_systemSymbolicSys
analytic_factorycallable

User provided callback for expressing analytic solutions to a set of dependent variables in original_system. The callback should have the signature: my_factory(x0, y0, p0, backend) -> dict, where the returned dictionary maps dependent variabels (from original_system.dep) to new expressions in remaining variables and initial conditions.

**kwargsdict

Keyword arguments passed onto SymbolicSys.

Examples

>>> odesys = SymbolicSys.from_callback(
...     lambda x, y, p: [
...         -p[0]*y[0],
...         p[0]*y[0] - p[1]*y[1]
...     ], 2, 2)
>>> dep0 = odesys.dep[0]
>>> partsys = PartiallySolvedSystem(odesys, lambda x0, y0, p0, be: {
...         dep0: y0[0]*be.exp(-p0[0]*(odesys.indep-x0))
...     })
>>> print(partsys.exprs)  
(_Dummy_29*p_0*exp(-p_0*(-_Dummy_28 + x)) - p_1*y_1,)
>>> y0, k = [3, 2], [3.5, 2.5]
>>> xout, yout, info = partsys.integrate([0, 1], y0, k, integrator='scipy')
>>> info['success'], yout.shape[1]
(True, 2)
Attributes
free_nameslist of str
analytic_exprslist of expressions
analytic_cbcallback
original_depdependent variable of original system

Methods

adaptive(self, y0, x0, xend[, params])

Integrate with integrator chosen output.

analytic_stiffness(self[, xyp])

Running stiffness ratio from last integration.

chained_parameter_variation(self, \*args, …)

See chained_parameter_variation().

from_callback(rhs[, ny, nparams, …])

Create an instance from a callback.

from_linear_invariants(ori_sys[, preferred])

Reformulates the ODE system in fewer variables.

from_other(ori, \*\*kwargs)

Creates a new instance with an existing one as a template.

from_other_new_params(ori, par_subs, new_pars)

Creates a new instance with an existing one as a template (with new parameters)

from_other_new_params_by_name(ori, par_subs)

Creates a new instance with an existing one as a template (with new parameters)

get_dfdx(self)

Calculates 2nd derivatives of self.exprs

get_dfdx_callback(self)

Generate a callback for evaluating derivative of self.exprs

get_f_ty_callback(self)

Generates a callback for evaluating self.exprs.

get_j_ty_callback(self)

Generates a callback for evaluating the jacobian.

get_jac(self)

Derives the jacobian from self.exprs and self.dep.

get_jtimes(self)

Derive the jacobian-vector product from self.exprs and self.dep

get_jtimes_callback(self)

Generate a callback fro evaluating the jacobian-vector product.

get_roots_callback(self)

Generate a callback for evaluating self.roots

integrate(self, \*args, \*\*kwargs)

Integrate the system of ordinary differential equations.

jacobian_singular(self)

Returns True if Jacobian is singular, else False.

plot_phase_plane(self[, indices])

Plots a phase portrait from last integration.

plot_result(self, \*\*kwargs)

Plots the integrated dependent variables from last integration.

post_process(self, xout, yout, params)

Transforms internal values to output, used internally.

pre_process(self, xout, y0[, params])

Transforms input to internal values, used internally.

predefined(self, y0, xout[, params])

Integrate with user chosen output.

stiffness(self[, xyp, eigenvals_cb])

[DEPRECATED] Use Result.stiffness(), stiffness ration

all_invariant_names

all_invariants

as_autonomous

get_first_step_callback

get_invariants_callback

to_arrays

classmethod from_linear_invariants(ori_sys, preferred=None, \*\*kwargs)[source]

Reformulates the ODE system in fewer variables.

Given linear invariant equations one can always reduce the number of dependent variables in the system by the rank of the matrix describing this linear system.

Parameters
ori_sysSymbolicSys instance
preferrediterable of preferred dependent variables

Due to numerical rounding it is preferable to choose the variables which are expected to be of the largest magnitude during integration.

**kwargs :

Keyword arguments passed on to constructor.

integrate(self, \*args, \*\*kwargs)[source]

Integrate the system of ordinary differential equations.

Solves the initial value problem (IVP).

Parameters
xarray_like or pair (start and final time) or float
if float:

make it a pair: (0, x)

if pair or length-2 array:

initial and final value of the independent variable

if array_like:

values of independent variable report at

y0array_like

Initial values at x[0] for the dependent variables.

paramsarray_like (default: tuple())

Value of parameters passed to user-supplied callbacks.

integratorstr or None
Name of integrator, one of:
  • ‘scipy’: _integrate_scipy()

  • ‘gsl’: _integrate_gsl()

  • ‘odeint’: _integrate_odeint()

  • ‘cvode’: _integrate_cvode()

See respective method for more information. If None: os.environ.get('PYODESYS_INTEGRATOR', 'scipy')

atolfloat

Absolute tolerance

rtolfloat

Relative tolerance

with_jacobianbool or None (default)

Whether to use the jacobian. When None the choice is done automatically (only used when required). This matters when jacobian is derived at runtime (high computational cost).

with_jtimesbool (default: False)

Whether to use the jacobian-vector product. This is only supported by cvode and only when linear_solver is one of: gmres’, ‘gmres_classic’, ‘bicgstab’, ‘tfqmr’. See the documentation for pycvodes for more information.

force_predefinedbool (default: False)

override behaviour of len(x) == 2 => adaptive()

**kwargs :

Additional keyword arguments for _integrate_$(integrator).

Returns
Length 3 tuple: (x, yout, info)

x : array of values of the independent variable yout : array of the dependent variable(s) for the different

values of x.

info : dict (‘nfev’ is guaranteed to be a key)

class pyodesys.symbolic.ScaledSys(dep_exprs, indep=None, dep_scaling=1, indep_scaling=1, params=, **kwargs)[source]

Bases: pyodesys.symbolic.TransformedSys

Transformed system where the variables have been scaled linearly.

Parameters
dep_exprsiterable of (symbol, expression)-pairs

see SymbolicSys

indepSymbol

see SymbolicSys

dep_scalingnumber (>0) or iterable of numbers

scaling of the dependent variables (default: 1)

indep_scalingnumber (>0)

scaling of the independent variable (default: 1)

params :

see SymbolicSys

**kwargs :

Keyword arguments passed onto TransformedSys.

Examples

>>> from sympy import Symbol
>>> x = Symbol('x')
>>> scaled1 = ScaledSys([(x, x*x)], dep_scaling=1000)
>>> scaled1.exprs == (x**2/1000,)
True
>>> scaled2 = ScaledSys([(x, x**3)], dep_scaling=1000)
>>> scaled2.exprs == (x**3/1000000,)
True
Attributes
linear_invariants
ny

Number of dependent variables in the system.

Methods

adaptive(self, y0, x0, xend[, params])

Integrate with integrator chosen output.

analytic_stiffness(self[, xyp])

Running stiffness ratio from last integration.

chained_parameter_variation(self, \*args, …)

See chained_parameter_variation().

from_callback(cb[, ny, nparams, …])

Create an instance from a callback.

from_other(ori, \*\*kwargs)

Creates a new instance with an existing one as a template.

from_other_new_params(ori, par_subs, new_pars)

Creates a new instance with an existing one as a template (with new parameters)

from_other_new_params_by_name(ori, par_subs)

Creates a new instance with an existing one as a template (with new parameters)

get_dfdx(self)

Calculates 2nd derivatives of self.exprs

get_dfdx_callback(self)

Generate a callback for evaluating derivative of self.exprs

get_f_ty_callback(self)

Generates a callback for evaluating self.exprs.

get_j_ty_callback(self)

Generates a callback for evaluating the jacobian.

get_jac(self)

Derives the jacobian from self.exprs and self.dep.

get_jtimes(self)

Derive the jacobian-vector product from self.exprs and self.dep

get_jtimes_callback(self)

Generate a callback fro evaluating the jacobian-vector product.

get_roots_callback(self)

Generate a callback for evaluating self.roots

integrate(self, x, y0[, params, atol, rtol])

Integrate the system of ordinary differential equations.

jacobian_singular(self)

Returns True if Jacobian is singular, else False.

plot_phase_plane(self[, indices])

Plots a phase portrait from last integration.

plot_result(self, \*\*kwargs)

Plots the integrated dependent variables from last integration.

post_process(self, xout, yout, params)

Transforms internal values to output, used internally.

pre_process(self, xout, y0[, params])

Transforms input to internal values, used internally.

predefined(self, y0, xout[, params])

Integrate with user chosen output.

stiffness(self[, xyp, eigenvals_cb])

[DEPRECATED] Use Result.stiffness(), stiffness ration

all_invariant_names

all_invariants

as_autonomous

get_first_step_callback

get_invariants_callback

to_arrays

classmethod from_callback(cb, ny=None, nparams=None, dep_scaling=1, indep_scaling=1, \*\*kwargs)[source]

Create an instance from a callback.

Analogous to SymbolicSys.from_callback().

Parameters
cbcallable

Signature rhs(x, y[:], p[:]) -> f[:]

nyint

length of y

nparamsint

length of p

dep_scalingnumber (>0) or iterable of numbers

scaling of the dependent variables (default: 1)

indep_scaling: number (>0)

scaling of the independent variable (default: 1)

**kwargs :

Keyword arguments passed onto ScaledSys.

Examples

>>> def f(x, y, p):
...     return [p[0]*y[0]**2]
>>> odesys = ScaledSys.from_callback(f, 1, 1, dep_scaling=10)
>>> odesys.exprs
(p_0*y_0**2/10,)
class pyodesys.symbolic.SymbolicSys(dep_exprs, indep=None, params=None, jac=True, dfdx=True, jtimes=False, first_step_expr=None, roots=None, backend=None, lower_bounds=None, upper_bounds=None, linear_invariants=None, nonlinear_invariants=None, linear_invariant_names=None, nonlinear_invariant_names=None, steady_state_root=False, init_indep=None, init_dep=None, sparse=False, **kwargs)[source]

Bases: pyodesys.core.ODESys

ODE System from symbolic expressions

Creates a ODESys instance from symbolic expressions. Jacobian and second derivatives are derived when needed.

Parameters
dep_exprsiterable of (symbol, expression)-pairs
indepSymbol

Independent variable (default: None => autonomous system).

paramsiterable of symbols

Problem parameters. If None: zero parameters assumed (violation of this will raise a ValueError), If True: params are deduced from (sorted) free_symbols.

jacImmutableMatrix or bool (default: True)
if True:

calculate jacobian from exprs

if False:

do not compute jacobian (use explicit steppers)

if instance of ImmutableMatrix:

user provided expressions for the jacobian

jtimesbool or iterable of (symbol, expression)-pairs
if True:

calculate jacobian-vector product from exprs

if False:

do not compute jacobian-vector product

if instance of iterable of (symbol, expression) pairs:

user-provided expressions for the jacobian-vector product, where the symbols are the elements of the vector being multiplied. Other than these symbols, the expressions should only involve the symbols in :attr`dep` and :attr`params`.

dfdxiterable of expressions

Derivatives of exprs with respect to :attr`indep`.

first_step_exprexpression

Closed form expression for calculating the first step. Be sure to pass first_step==0 to enable its use during integration. If not given, the solver default behavior will be invoked.

rootsiterable of expressions

Equations to look for root’s for during integration (currently available through cvode).

backendstr or sym.Backend

See documentation of sym.Backend.

lower_boundsarray_like

Convenience option setting magnitude constraint. (requires integrator with support for recoverable errors)

upper_boundsarray_like

Convenience option setting magnitude constraint. (requires integrator with support for recoverable errors)

linear_invariantsMatrix

Matrix specifing linear combinations of dependent variables that

nonlinear_invariantsiterable of expressions

Iterable collection of expressions of nonlinear invariants.

steady_state_rootbool or float

Generate an expressions for roots which is the sum the smaller of aboslute values of derivatives or relative derivatives subtracted by the value of steady_state_root (default 1e-10).

init_indepSymbol, True or None

When True construct using be.Symbol. See also init_indep.

init_deptuple of Symbols, True or None

When True construct using be.Symbol. See also init_dep.

sparsebool (default False)

When True, jacobian will be derived and stored (if jac = True) in compressed sparse column (CSC) format and the jacobian callback will return an instance of scipy.sparse.csc_matrix.

**kwargs:

See ODESys

Attributes
deptuple of symbols

Dependent variables.

exprstuple of expressions

Expressions for the derivatives of the dependent variables (dep) with respect to the independent variable (indep).

indepSymbol or None

Independent variable (None indicates autonomous system).

paramsiterable of symbols

Problem parameters.

first_step_exprexpression

Closed form expression for how to compute the first step.

rootsiterable of expressions or None

Roots to report for during integration.

backendmodule

Symbolic backend, e.g. ‘sympy’, ‘symengine’.

nyint

Number of dependent variables in the system.

init_indepSymbol, None

Symbol for initial value of independent variable (before pre processors).

init_deptuple of Symbols or None

Symbols for initial values of dependent variables (before pre processors).

Methods

adaptive(self, y0, x0, xend[, params])

Integrate with integrator chosen output.

analytic_stiffness(self[, xyp])

Running stiffness ratio from last integration.

chained_parameter_variation(self, \*args, …)

See chained_parameter_variation().

from_callback(rhs[, ny, nparams, …])

Create an instance from a callback.

from_other(ori, \*\*kwargs)

Creates a new instance with an existing one as a template.

from_other_new_params(ori, par_subs, new_pars)

Creates a new instance with an existing one as a template (with new parameters)

from_other_new_params_by_name(ori, par_subs)

Creates a new instance with an existing one as a template (with new parameters)

get_dfdx(self)

Calculates 2nd derivatives of self.exprs

get_dfdx_callback(self)

Generate a callback for evaluating derivative of self.exprs

get_f_ty_callback(self)

Generates a callback for evaluating self.exprs.

get_j_ty_callback(self)

Generates a callback for evaluating the jacobian.

get_jac(self)

Derives the jacobian from self.exprs and self.dep.

get_jtimes(self)

Derive the jacobian-vector product from self.exprs and self.dep

get_jtimes_callback(self)

Generate a callback fro evaluating the jacobian-vector product.

get_roots_callback(self)

Generate a callback for evaluating self.roots

integrate(self, x, y0[, params, atol, rtol])

Integrate the system of ordinary differential equations.

jacobian_singular(self)

Returns True if Jacobian is singular, else False.

plot_phase_plane(self[, indices])

Plots a phase portrait from last integration.

plot_result(self, \*\*kwargs)

Plots the integrated dependent variables from last integration.

post_process(self, xout, yout, params)

Transforms internal values to output, used internally.

pre_process(self, xout, y0[, params])

Transforms input to internal values, used internally.

predefined(self, y0, xout[, params])

Integrate with user chosen output.

stiffness(self[, xyp, eigenvals_cb])

[DEPRECATED] Use Result.stiffness(), stiffness ration

all_invariant_names

all_invariants

as_autonomous

get_first_step_callback

get_invariants_callback

to_arrays

all_invariant_names(self)[source]
all_invariants(self, linear_invariants=None, nonlinear_invariants=None, dep=None, backend=None)[source]
analytic_stiffness(self, xyp=None)[source]

Running stiffness ratio from last integration.

Calculate sittness ratio, i.e. the ratio between the largest and smallest absolute eigenvalue of the (analytic) jacobian matrix.

See ODESys.stiffness() for more info.

append_iv = True
as_autonomous(self, new_indep_name=None, new_latex_indep_name=None)[source]
classmethod from_callback(rhs, ny=None, nparams=None, first_step_factory=None, roots_cb=None, indep_name=None, \*\*kwargs)[source]

Create an instance from a callback.

Parameters
rhscallbable

Signature rhs(x, y[:], p[:], backend=math) -> f[:].

nyint

Length of y in rhs.

nparamsint

Length of p in rhs.

first_step_factorycallabble

Signature step1st(x, y[:], p[:]) -> dx0.

roots_cbcallable

Callback with signature roots(x, y[:], p[:], backend=math) -> r[:].

indep_namestr

Default ‘x’ if not already in names, otherwise indep0, or indep1, or …

dep_by_namebool

Make y passed to rhs a dict (keys from names) and convert its return value from dict to array.

par_by_namebool

Make p passed to rhs a dict (keys from param_names).

**kwargs :

Keyword arguments passed onto SymbolicSys.

Returns
An instance of SymbolicSys.

Examples

>>> def decay(x, y, p, backend=None):
...     rate = y['Po-210']*p[0]
...     return {'Po-210': -rate, 'Pb-206': rate}
...
>>> odesys = SymbolicSys.from_callback(decay, dep_by_name=True, names=('Po-210', 'Pb-206'), nparams=1)
>>> xout, yout, info = odesys.integrate([0, 138.4*24*3600], {'Po-210': 1.0, 'Pb-206': 0.0}, [5.798e-8])
>>> import numpy as np; np.allclose(yout[-1, :], [0.5, 0.5], rtol=1e-3, atol=1e-3)
True
classmethod from_other(ori, \*\*kwargs)[source]

Creates a new instance with an existing one as a template.

Parameters
oriSymbolicSys instance
**kwargs:

Keyword arguments used to create the new instance.

Returns
A new instance of the class.
classmethod from_other_new_params(ori, par_subs, new_pars, new_par_names=None, new_latex_par_names=None, \*\*kwargs)[source]

Creates a new instance with an existing one as a template (with new parameters)

Calls .from_other but first it replaces some parameters according to par_subs and (optionally) introduces new parameters given in new_pars.

Parameters
oriSymbolicSys instance
par_subsdict

Dictionary with substitutions (mapping symbols to new expressions) for parameters. Parameters appearing in this instance will be omitted in the new instance.

new_parsiterable (optional)

Iterable of symbols for new parameters.

new_par_namesiterable of str

Names of the new parameters given in new_pars.

new_latex_par_namesiterable of str

TeX formatted names of the new parameters given in new_pars.

**kwargs:

Keyword arguments passed to .from_other.

Returns
Intance of the class
extradict with keys:
  • recalc_params : f(t, y, p1) -> p0

classmethod from_other_new_params_by_name(ori, par_subs, new_par_names=(), \*\*kwargs)[source]

Creates a new instance with an existing one as a template (with new parameters)

Calls .from_other_new_params but first it creates the new instances from user provided callbacks generating the expressions the parameter substitutions.

Parameters
oriSymbolicSys instance
par_subsdict mapping str to f(t, y{}, p{}) -> expr

User provided callbacks for parameter names in ori.

new_par_namesiterable of str
**kwargs:

Keyword arguments passed to .from_other_new_params.

get_dfdx(self)[source]

Calculates 2nd derivatives of self.exprs

get_dfdx_callback(self)[source]

Generate a callback for evaluating derivative of self.exprs

get_f_ty_callback(self)[source]

Generates a callback for evaluating self.exprs.

get_first_step_callback(self)[source]
get_invariants_callback(self)[source]
get_j_ty_callback(self)[source]

Generates a callback for evaluating the jacobian.

get_jac(self)[source]

Derives the jacobian from self.exprs and self.dep.

get_jtimes(self)[source]

Derive the jacobian-vector product from self.exprs and self.dep

get_jtimes_callback(self)[source]

Generate a callback fro evaluating the jacobian-vector product.

get_roots_callback(self)[source]

Generate a callback for evaluating self.roots

jacobian_singular(self)[source]

Returns True if Jacobian is singular, else False.

property linear_invariants
property ny

Number of dependent variables in the system.

to_arrays(self, x, y, p, \*\*kwargs)[source]
class pyodesys.symbolic.TransformedSys(dep_exprs, indep=None, dep_transf=None, indep_transf=None, params=, exprs_process_cb=None, check_transforms=True, **kwargs)[source]

Bases: pyodesys.symbolic.SymbolicSys

SymbolicSys with variable transformations.

Parameters
dep_exprsiterable of pairs

see SymbolicSys

indepSymbol

see SymbolicSys

dep_transfiterable of (expression, expression) pairs

pairs of (forward, backward) transformations for the dependents variables

indep_transfpair of expressions

forward and backward transformation of the independent variable

params :

see SymbolicSys

exprs_process_cbcallbable

Post processing of the expressions (signature: f(exprs) -> exprs) for the derivatives of the dependent variables after transformation have been applied.

check_transformsbool

Passed as keyword argument check to util.transform_exprs_dep() and util.transform_exprs_indep().

**kwargs :

Keyword arguments passed onto SymbolicSys.

Attributes
linear_invariants
ny

Number of dependent variables in the system.

Methods

adaptive(self, y0, x0, xend[, params])

Integrate with integrator chosen output.

analytic_stiffness(self[, xyp])

Running stiffness ratio from last integration.

chained_parameter_variation(self, \*args, …)

See chained_parameter_variation().

from_callback(cb[, ny, nparams, …])

Create an instance from a callback.

from_other(ori, \*\*kwargs)

Creates a new instance with an existing one as a template.

from_other_new_params(ori, par_subs, new_pars)

Creates a new instance with an existing one as a template (with new parameters)

from_other_new_params_by_name(ori, par_subs)

Creates a new instance with an existing one as a template (with new parameters)

get_dfdx(self)

Calculates 2nd derivatives of self.exprs

get_dfdx_callback(self)

Generate a callback for evaluating derivative of self.exprs

get_f_ty_callback(self)

Generates a callback for evaluating self.exprs.

get_j_ty_callback(self)

Generates a callback for evaluating the jacobian.

get_jac(self)

Derives the jacobian from self.exprs and self.dep.

get_jtimes(self)

Derive the jacobian-vector product from self.exprs and self.dep

get_jtimes_callback(self)

Generate a callback fro evaluating the jacobian-vector product.

get_roots_callback(self)

Generate a callback for evaluating self.roots

integrate(self, x, y0[, params, atol, rtol])

Integrate the system of ordinary differential equations.

jacobian_singular(self)

Returns True if Jacobian is singular, else False.

plot_phase_plane(self[, indices])

Plots a phase portrait from last integration.

plot_result(self, \*\*kwargs)

Plots the integrated dependent variables from last integration.

post_process(self, xout, yout, params)

Transforms internal values to output, used internally.

pre_process(self, xout, y0[, params])

Transforms input to internal values, used internally.

predefined(self, y0, xout[, params])

Integrate with user chosen output.

stiffness(self[, xyp, eigenvals_cb])

[DEPRECATED] Use Result.stiffness(), stiffness ration

all_invariant_names

all_invariants

as_autonomous

get_first_step_callback

get_invariants_callback

to_arrays

classmethod from_callback(cb, ny=None, nparams=None, dep_transf_cbs=None, indep_transf_cbs=None, roots_cb=None, \*\*kwargs)[source]

Create an instance from a callback.

Analogous to SymbolicSys.from_callback().

Parameters
cbcallable

Signature rhs(x, y[:], p[:]) -> f[:]

nyint

length of y

nparamsint

length of p

dep_transf_cbsiterable of pairs callables

callables should have the signature f(yi) -> expression in yi

indep_transf_cbspair of callbacks

callables should have the signature f(x) -> expression in x

roots_cbcallable

Callback with signature roots(x, y[:], p[:], backend=math) -> r[:]. Callback should return untransformed roots.

**kwargs :

Keyword arguments passed onto TransformedSys.

pyodesys.symbolic.get_logexp(a=1, b=0, a2=None, b2=None, backend=None)[source]

Utility function for use with :func:symmetricsys.

Creates a pair of callbacks for logarithmic transformation (including scaling and shifting): u = ln(a*x + b).

Parameters
anumber

Scaling (forward).

bnumber

Shift (forward).

a2number

Scaling (backward).

b2number

Shift (backward).

Returns
Pair of callbacks.
pyodesys.symbolic.symmetricsys(dep_tr=None, indep_tr=None, SuperClass=<class 'pyodesys.symbolic.TransformedSys'>, \*\*kwargs)[source]

A factory function for creating symmetrically transformed systems.

Creates a new subclass which applies the same transformation for each dependent variable.

Parameters
dep_trpair of callables (default: None)

Forward and backward transformation callbacks to be applied to the dependent variables.

indep_trpair of callables (default: None)

Forward and backward transformation to be applied to the independent variable.

SuperClassclass
**kwargs :

Default keyword arguments for the TransformedSys subclass.

Returns
Subclass of SuperClass (by default TransformedSys).

Examples

>>> import sympy
>>> logexp = (sympy.log, sympy.exp)
>>> def psimp(exprs):
...     return [sympy.powsimp(expr.expand(), force=True) for expr in exprs]
...
>>> LogLogSys = symmetricsys(logexp, logexp, exprs_process_cb=psimp)
>>> mysys = LogLogSys.from_callback(lambda x, y, p: [-y[0], y[0] - y[1]], 2, 0)
>>> mysys.exprs
(-exp(x_0), -exp(x_0) + exp(x_0 + y_0 - y_1))

pyodesys.util module

class pyodesys.util.MissingImport(modname, exc)[source]

Bases: object

Methods

__call__(self, \*args, \*\*kwargs)

Call self as a function.

pyodesys.util.check_transforms(fw, bw, symbs)[source]

Verify validity of a pair of forward and backward transformations

Parameters
fw: expression

forward transformation

bw: expression

backward transformation

symbs: iterable of symbols

the variables that are transformed

pyodesys.util.import_(modname, \*args)[source]
pyodesys.util.merge_dicts(\*dicts)[source]

Merges dictionaries with incresing priority.

Parameters
*dicts: dictionaries
Returns
dict

Examples

>>> d1, d2 = {'a': 1, 'b': 2}, {'a': 2, 'c': 3}
>>> merge_dicts(d1, d2, {'a': 3}) == {'a': 3, 'b': 2, 'c': 3}
True
>>> d1 == {'a': 1, 'b': 2}
True
>>> from collections import defaultdict
>>> dd1 = defaultdict(lambda: 3, {'b': 4})
>>> dd2 = merge_dicts(dd1, {'c': 5}, {'c': 17})
>>> dd2['c'] - dd2['a'] - dd2['b'] == 10
True
pyodesys.util.pycvodes_double(cb)[source]
pyodesys.util.pycvodes_klu(cb)[source]
class pyodesys.util.requires(*reqs)[source]

Bases: object

Conditional skipping (on requirements) of tests in pytest

Examples

>>> @requires('numpy', 'scipy')
... def test_sqrt():
...     import numpy as np
...     assert np.sqrt(4) == 2
...     from scipy.special import zeta
...     assert zeta(2) < 2
...
>>> @requires('numpy>=1.9.0')
... def test_nanmedian():
...     import numpy as np
...     a = np.array([[10.0, 7, 4], [3, 2, 1]])
...     a[0, 1] = np.nan
...     assert np.nanmedian(a) == 3
...

Methods

__call__(self, cb)

Call self as a function.

eq(a, b, /)

Same as a == b.

ge(a, b, /)

Same as a >= b.

gt(a, b, /)

Same as a > b.

le(a, b, /)

Same as a <= b.

lt(a, b, /)

Same as a < b.

ne(a, b, /)

Same as a != b.

eq(a, b, /)

Same as a == b.

ge(a, b, /)

Same as a >= b.

gt(a, b, /)

Same as a > b.

le(a, b, /)

Same as a <= b.

lt(a, b, /)

Same as a < b.

ne(a, b, /)

Same as a != b.

pyodesys.util.stack_1d_on_left(x, y)[source]

Stack a 1D array on the left side of a 2D array

Parameters
x: 1D array
y: 2D array

Requirement: shape[0] == x.size

pyodesys.util.transform_exprs_dep(fw, bw, dep_exprs, check=True)[source]

Transform y[:] in dydx

Parameters
fw: expression

forward transformation

bw: expression

backward transformation

dep_exprs: iterable of (symbol, expression) pairs

pairs of (dependent variable, derivative expressions), i.e. (y, dydx) pairs

check: bool (default: True)

whether to verification of the analytic correctness should be performed

Returns
List of transformed expressions for dydx
pyodesys.util.transform_exprs_indep(fw, bw, dep_exprs, indep, check=True)[source]

Transform x in dydx

Parameters
fw: expression

forward transformation

bw: expression

backward transformation

dep_exprs: iterable of (symbol, expression) pairs

pairs of (dependent variable, derivative expressions)

check: bool (default: True)

whether to verification of the analytic correctness should be performed

Returns
List of transformed expressions for dydx