Lower Level API

scikits.odes

scikits.odes is a scikit offering a cython wrapper around some extra ode/dae solvers, so they can mature outside of scipy.

It offers wrappers around the following solvers from SUNDIALS
  • CVODE

  • IDA

It additionally offers wrappers around
  • ddaspk <http://www.netlib.org/ode/ddaspk.f> (included)

  • lsodi <http://www.netlib.org/ode/lsodi.f> (included)

scikits.odes.ode

First-order ODE solver

User-friendly interface to various numerical integrators for solving a system of first order ODEs with prescribed initial conditions:

\[ \begin{align}\begin{aligned}\frac{dy(t)}{dt} = f(t, y(t)), \quad y(t_0)=y_0\\y(t_0)[i] = y_0[i], i = 0, ..., \mathrm{len}(y_0) - 1\end{aligned}\end{align} \]

\(f(t,y)\) is the right hand side function and returns a vector of size \(\mathrm{len}(y_0)\).

class scikits.odes.ode.OdeBase(Rfn, **options)[source]

The interface which ODE solvers must implement.

Parameters
  • Rfn (function) – A function which computes the required derivatives. The signature should be func(t, y, y_dot, *args, **kwargs). Note that *args and **kwargs handling are solver dependent.

  • options (mapping) – Additional options for initialization, solver dependent

init_step(t0, y0)[source]

Initializes the solver and allocates memory.

Parameters
  • t0 (number) – initial time

  • y0 (array) – initial condition for y (can be list or numpy array)

Returns

  • old_api is False (namedtuple) – namedtuple with the following attributes

    Field

    Meaning

    flag

    An integer flag (StatusEnum)

    values

    Named tuple with entries t and y

    errors

    Named tuple with entries t and y

    roots

    Named tuple with entries t and y

    tstop

    Named tuple with entries t and y

    message

    String with message in case of an error

  • old_api is True (tuple) – tuple with the following elements in order

    Field

    Meaning

    flag

    boolean status of the computation (successful or error occured)

    t_out

    inititial time

set_options(**options)[source]

Set specific options for the solver.

Calling set_options a second time, normally resets the solver.

solve(tspan, y0)[source]

Runs the solver.

Parameters
  • tspan (array (or similar)) – a list of times at which the computed value will be returned. Must contain the start time as first entry.

  • y0 (array (or similar)) – a list of initial values

Returns

  • old_api is False (namedtuple) – namedtuple with the following attributes

    Field

    Meaning

    flag

    An integer flag (StatusEnum)

    values

    Named tuple with entries t and y

    errors

    Named tuple with entries t and y

    roots

    Named tuple with entries t and y

    tstop

    Named tuple with entries t and y

    message

    String with message in case of an error

  • old_api is True (tuple) – tuple with the following elements in order

    Field

    Meaning

    flag

    indicating return status of the solver

    t

    numpy array of times at which the computations were successful

    y

    numpy array of values corresponding to times t (values of y[i, :] ~ t[i])

    t_err

    float or None - if recoverable error occured (for example reached maximum number of allowed iterations), this is the time at which it happened

    y_err

    numpy array of values corresponding to time t_err

step(t, y_retn=None)[source]

Method for calling successive next step of the ODE solver to allow more precise control over the solver. The ‘init_step’ method has to be called before the ‘step’ method.

A step is done towards time t, and output at t returned. This time can be higher or lower than the previous time. If option ‘one_step_compute’==True, and the solver supports it, only one internal solver step is done in the direction of t starting at the current step.

If old_api=True, the old behavior is used: if t>0.0 then integration is performed until this time and results at this time are returned in y_retn if t<0.0 only one internal step is perfomed towards time abs(t) and results after this one time step are returned

Parameters
  • t (number) –

  • y_retn (numpy vector (ndim = 1)) – in which the computed value will be stored (needs to be preallocated). If None y_retn is not used.

Returns

  • old_api is False (namedtuple) – namedtuple with the following attributes

    Field

    Meaning

    flag

    An integer flag (StatusEnum)

    values

    Named tuple with entries t and y

    errors

    Named tuple with entries t and y

    roots

    Named tuple with entries t and y

    tstop

    Named tuple with entries t and y

    message

    String with message in case of an error

  • old_api is True (tuple) – tuple with the following elements in order

    Field

    Meaning

    flag

    status of the computation (successful or error occured)

    t_out

    time, where the solver stopped (when no error occured, t_out == t)

scikits.odes.dae

First-order DAE solver

User-friendly interface to various numerical integrators for solving an algebraic system of first order ODEs with prescribed initial conditions:

\[ \begin{align}\begin{aligned}A \frac{dy(t)}{dt} = f(t,y(t)),\\y(t=0)[i] = y0[i],\\\frac{d y(t=0)}{dt}[i] = yprime0[i],\end{aligned}\end{align} \]

where \(i = 0, ..., len(y0) - 1\); \(A\) is a (possibly singular) matrix of size \(i × i\); and \(f(t,y)\) is a vector of size \(i\) or more generally, equations of the form

\[G(t,y,y') = 0\]
class scikits.odes.dae.DaeBase(Rfn, **options)[source]

The interface which DAE solvers must implement.

Parameters
  • Rfn – residual function

  • options (mapping) – Additional options for initialization, solver dependent

init_step(t0, y0, yp0, y_ic0_retn=None, yp_ic0_retn=None)[source]

Initializes the solver and allocates memory.

Parameters
  • t0 (number) – initial time

  • y0 (list/array) – initial condition for y

  • yp0 (list/array) – initial condition for yp

  • y_ic0 (numpy array) – (optional) returns the calculated consistent initial condition for y

  • yp_ic0 (numpy array) – (optional) returns the calculated consistent initial condition for y derivated.

Returns

  • old_api is False (namedtuple) – namedtuple with the following attributes

    Field

    Meaning

    flag

    An integer flag (StatusEnumXXX)

    values

    Named tuple with entries t and y and ydot. y will correspond to y_retn value and ydot to yp_retn!

    errors

    Named tuple with entries t and y and ydot

    roots

    Named tuple with entries t and y and ydot

    tstop

    Named tuple with entries t and y and ydot

    message

    String with message in case of an error

  • old_api is True (tuple) – tuple with the following elements in order

    Field

    Meaning

    flag

    status of the computation (successful or error occured)

    t_out

    time, where the solver stopped (when no error occured, t_out == t)

set_options(**options)[source]

Set specific options for the solver.

Calling set_options a second time, normally resets the solver.

solve(tspan, y0, yp0)[source]

Runs the solver.

Parameters
  • tspan (list/array) – A list of times at which the computed value will be returned. Must contain the start time as first entry.

  • y0 (list/array) – List of initial values

  • yp0 (list/array) – List of initial values of derivatives

Returns

  • old_api is False (namedtuple) – namedtuple with the following attributes

    Field

    Meaning

    flag

    An integer flag (StatusEnumXXX)

    values

    Named tuple with entries array t and array y and array ydot. y will correspond to y_retn value and ydot to yp_retn!

    errors

    Named tuple with entries t and y and ydot of error

    roots

    Named tuple with entries array t and array y and array ydot

    tstop

    Named tuple with entries array t and array y and array ydot

    message

    String with message in case of an error

  • old_api is True (tuple) – tuple with the following elements in order

    Field

    Meaning

    flag

    indicating return status of the solver

    t

    numpy array of times at which the computations were successful

    y

    numpy array of values corresponding to times t (values of y[i, :] ~ t[i])

    yp

    numpy array of derivatives corresponding to times t (values of yp[i, :] ~ t[i])

    t_err

    float or None - if recoverable error occured (for example reached maximum number of allowed iterations), this is the time at which it happened

    y_err

    numpy array of values corresponding to time t_err

    yp_err

    numpy array of derivatives corresponding to time t_err

step(t, y_retn=None, yp_retn=None)[source]

Method for calling successive next step of the IDA solver to allow more precise control over the IDA solver. The ‘init_step’ method has to be called before the ‘step’ method.

A step is done towards time t, and output at t returned. This time can be higher or lower than the previous time. If option ‘one_step_compute’==True, and the solver supports it, only one internal solver step is done in the direction of t starting at the current step.

If old_api=True, the old behavior is used: if t>0.0 then integration is performed until this time and results at this time are returned in y_retn; else if if t<0.0 only one internal step is perfomed towards time abs(t) and results after this one time step are returned.

Parameters
  • t (number) –

  • y_retn (numpy array (ndim = 1) or None.) – (Needs to be preallocated) If not None, will be filled with y at time t. If None y_retn is not used.

  • yp_retn (numpy array (ndim = 1) or None.) – (Needs to be preallocated) If not None, will be filled with derivatives of y at time t. If None yp_retn is not used.

Returns

  • old_api is False (namedtuple) – namedtuple with the following attributes

    Field

    Meaning

    flag

    An integer flag (StatusEnumXXX)

    values

    Named tuple with entries t and y and ydot. y will correspond to y_retn value and ydot to yp_retn!

    errors

    Named tuple with entries t and y and ydot

    roots

    Named tuple with entries t and y and ydot

    tstop

    Named tuple with entries t and y and ydot

    message

    String with message in case of an error

  • old_api is True (tuple) – tuple with the following elements in order

    Field

    Meaning

    flag

    status of the computation (successful or error occured)

    t_out

    time, where the solver stopped (when no error occured, t_out == t)

scikits.odes.dopri5

Making scipy ode solvers available via the ode API

dopri5

This is an explicit runge-kutta method of order (4)5 due to Dormand & Prince (with stepsize control and dense output). The API of this solver is as the other scikit.odes ODE solvers

Authors:

  1. Hairer and G. Wanner Universite de Geneve, Dept. de Mathematiques CH-1211 Geneve 24, Switzerland e-mail: ernst.hairer@math.unige.ch, gerhard.wanner@math.unige.ch

This code is described in [HNW93].

This integrator accepts the following options:

atol : float or sequence absolute tolerance for solution, default 1e-12 rtol : float or sequence relative tolerance for solution, default 1e-6 nsteps : int Maximum number of (internally defined) steps allowed during one call to the solver. Default=500 first_step : float max_step : float safety : float Safety factor on new step selection (default 0.9) ifactor : float dfactor : float Maximum factor to increase/decrease step size by in one step beta : float Beta parameter for stabilised step size control. verbosity : int Switch for printing messages (< 0 for no messages).

References [HNW93] (1, 2) E. Hairer, S.P. Norsett and G. Wanner, Solving Ordinary Differential Equations i. Nonstiff Problems. 2nd edition. Springer Series in Computational Mathematics, Springer-Verlag (1993)

dop853

This is an explicit runge-kutta method of order 8(5,3) due to Dormand & Prince (with stepsize control and dense output). Options and references the same as “dopri5”.

exception scikits.odes.dopri5.DOPSolveException(soln)[source]

Base class for exceptions raised by DOP.validate_flags.

exception scikits.odes.dopri5.DOPSolveFailed(soln)[source]

DOP.solve failed to reach endpoint

class scikits.odes.dopri5.SolverReturn(flag, values, errors, roots, tstop, message)
property errors

Alias for field number 2

property flag

Alias for field number 0

property message

Alias for field number 5

property roots

Alias for field number 3

property tstop

Alias for field number 4

property values

Alias for field number 1

class scikits.odes.dopri5.SolverVariables(t, y)
property t

Alias for field number 0

property y

Alias for field number 1

class scikits.odes.dopri5.StatusEnumDOP[source]

An enumeration.

class scikits.odes.dopri5.dop853(Rfn, **options)[source]
class scikits.odes.dopri5.dopri5(Rfn, **options)[source]
init_step(t0, y0)[source]

Initializes the solver and allocates memory.

Parameters
  • - initial time (t0) –

  • - initial condition for y (can be list or numpy array) (y0) –

Returns

  • if old_api – not supported

  • if old_api False

    A named tuple, with entries:

    flag = An integer flag (StatusEnumDop) values = Named tuple with entries t and y and ydot. y will

    correspond to y_retn value and ydot to yp_retn!

    errors = Named tuple with entries t_err and y_err roots = Named tuple with entries t_roots and y_roots tstop = Named tuple with entries t_stop and y_tstop message= String with message in case of an error

set_options(**options)[source]

Set specific options for the solver.

Calling set_options a second time, normally resets the solver.

solve(tspan, y0)[source]

Runs the solver.

Parameters
  • - an list/array of times at which the computed value will be (tspan) – returned. Must contain the start time as first entry..

  • - list/numpy array of initial values (y0) –

Returns

  • if old_api – Not supported

  • if old_api False

    A named tuple, with entries:

    flag = An integer flag values = Named tuple with entries t and y errors = Named tuple with entries t and y roots = Named tuple with entries t and y tstop = Named tuple with entries t and y message= String with message in case of an error

step(t, y_retn=None)[source]

Method for calling successive next step of the ODE solver to allow more precise control over the solver. The ‘init_step’ method has to be called before the ‘step’ method.

Parameters
  • - A step is done towards time t, and output at t returned. (t) –

    This time can be higher or lower than the previous time. If option ‘one_step_compute’==True, and the solver supports it, only one internal solver step is done in the direction of t starting at the current step.

    If old_api=True, the old behavior is used:
    if t>0.0 then integration is performed until this time

    and results at this time are returned in y_retn

    if t<0.0 only one internal step is perfomed towards time abs(t)

    and results after this one time step are returned

  • - numpy vector (ndim = 1) in which the computed (y_retn) – value will be stored (needs to be preallocated). If None y_retn is not used.

Returns

  • if old_api – not supported

  • if old_api False

    A named tuple, with entries:

    flag = An integer flag (StatusEnumDOP) values = Named tuple with entries t and y. y will

    correspond to y_retn value

    errors = Named tuple with entries t_err and y_err roots = Named tuple with entries t_roots and y_roots tstop = Named tuple with entries t_stop and y_tstop message= String with message in case of an error

validate_flags(soln)[source]

Validates the flag returned by dopri.solve.

Validation happens using the following scheme:
  • failures (flag < 0) raise DOPSolveFailed or a subclass of it;

  • otherwise, return an instance of SolverReturn.

scikits.odes.ddaspkint

ddaspk

Solver developed 1989 to 1996, with some corrections from 2000 - Fortran

This code solves a system of differential/algebraic equations of the form G(t,y,y’) = 0 , using a combination of Backward Differentiation Formula (BDF) methods and a choice of two linear system solution methods: direct (dense or band) or Krylov (iterative). Krylov is not supported from within scikits.odes. In order to support it, a new interface should be created ddaspk_krylov, with a different signature, reflecting the changes needed.

Source: http://www.netlib.org/ode/ddaspk.f

On construction the function calculating the residual (res) must be given and optionally also the function calculating the jacobian (jac). res has the signature: res(x, y, yprime, return_residual) with

x : independent variable, eg the time, float y : array of n unknowns in x yprime : dy/dx array of n unknowns in x, dimension = dim(y) return_residual: array that must be updated with the value of the residuals,

so G(t,y,y’). The dimension is equal to dim(y)

return value: Not needed. However, for use with other solvers, consider

returning 0 on success.

jac has the signature jac(x, y, yprime, cj, return_jac) with

x : independent variable, eg the time, float y : array of n unknowns in x yprime : dy/dx array of n unknowns in x, dimension = dim(y) cj : internal variable of ddaspk algorithm you can use, don’t change it!

cj can be ignored, or used to rescale constraint equations in the system

return_jactwo dimensional array of the Jacobian, as per the Jacobian
definition of the ddaspk solver. This means it should be or a

full nxn shaped array in general (n=dim(y)), or a banded shaped array as per the definition of lband/uband. Jac is optional and should be set with the jacfn option keyword. Note that jac is defined as

dres(i)/dy(j) + cj*dres(i)/dyprime(j)

return value: Not needed. However, for use with other solvers, consider

returning 0 on success.

This integrator accepts the following parameters in the initializer or set_options method of the dae class:

  • rfnresidual function, see above for signature. This option need not be

    set, as rfn will be set during initialization. If the residual function of initialization needs to be reset, this option can be used

  • jacfn : jacobian function, see above for signature. Default is None.

  • atol : float or sequence of length i absolute tolerance for solution

  • rtol : float or sequence of length i relative tolerance for solution

  • lband : None or int

  • uband : None or int Jacobian band width, jac[i,j] != 0 for i-lband <= j <= i+uband. Setting these requires your jac routine to return the jacobian in packed format, jac_packed[i-j+lband, j] = jac[i,j].

  • tstop : None or float. If given, tstop is a critical time point beyond which no integration occurs

  • max_steps : int Maximum number of (internally defined) steps allowed during one call to the solver.

  • first_step : float Set initial stepsize. DAE solver can suffer on the first step, set this to circumvent this.

  • max_step_size : float Limits for the step sizes used by the integrator. This overrides the internal value

  • order : int Maximum order used by the integrator, >=1, <= 5 for BDF. 5 is the default

  • enforce_nonnegativity: bool Enforce the nonnegativity of Y during integration Note: best is to run code first with no condition

  • compute_initcond: None or ‘yp0’ or ‘y0’ DDASPK may be able to compute the initial conditions if you do not know them precisely. If y0, then y0 will be calculated If yp0, then the differential variables will be used to solve for the

    algebraic variables and the derivative of the differential variables. Which are the algebraic variables must be indicated with algebraic_var method

  • exclude_algvar_from_error: bool To determine solution, do not take the algebraic variables error control into account. Default=False

  • constraint_init: bool Enforce the constraint checks of Y during initial condition computation Note: try first with no constraints

  • constraint_type: if constraint_init, give an integer array with for every unknown the specific condition to check:

    1: y0[i] >= 0 2: y0[i] > 0

    -1: y0[i] <= 0 -2: y0[i] < 0

    0: y0[i] not constrained

    Alternatively, pass only one integer that then applies to all unknowns

  • algebraic_vars_idx: an array or None (= default)
    Description:

    If given problem is of type DAE, some items of the residual vector returned by the ‘resfn’ have to be treated as algebraic variables. These are denoted by the position (index) in the residual vector. All these indexes have to be specified in the ‘algebraic_vars_idx’ array.

class scikits.odes.ddaspkint.ddaspk(resfn, **options)[source]
init_step(t0, y0, yp0, y_ic0_retn=None, yp_ic0_retn=None)[source]

See dae.DaeBase

set_options(**options)[source]

See dae.DaeBase

solve(tspan, y0, yp0, hook_fn=None)[source]

See dae.DaeBase

step(t, y_retn, yp_retn=None)[source]

See dae.DaeBase

scikits.odes.lsodiint

lsodi

Solver developed during the 1980s, this is the version from 1987 - Fortran

Integrator for linearly implicit systems of first-order odes. lsodi provides the same methods as vode (adams and bdf). This integrator accepts the following parameters in initializer or set_options method of the dae class:

  • rfnresidual function, see below for signature. This option need not be

    set, as rfn will be set during initialization. If the residual function of initialization needs to be reset, this option can be used

  • jacfn : jacobian function, see below for signature. Default is None.

  • adda_func = function, required, see below

  • atol=float|seq

  • rtol=float

  • lband=None|int

  • rband=None|int

  • method=’adams’|’bdf’

  • with_jacobian=0|1

  • max_steps = int

  • max_step_size = float

  • (first|min|max)_step = float

  • tstop = None|float

  • order = int # <=12 for adams, <=5 for bdf

  • compute_initcond = None|’yp0’

Details:

Implicit integration requires two functions do be defined to integrate

d y(t)[i]

a(t,y) * ——— = g(t,y)[i]

dt

where a(t,y) is a square matrix.

  • rfn computes the residual which the user provides. rfn has the signature: rfn(x, y, yprime, return_residual)

with

x : independent variable, eg the time, float y : array of n unknowns in x yprime : dy/dx array of n unknowns in x, dimension = dim(y) return_residual: array that must be updated with the value of the residuals,

so G(t,y,y’). The dimension is equal to dim(y)

return value: Not needed. However, for use with other solvers, consider

returning 0 on success.

It calculates something like
def res(t, y, yprime, r):

r(t,y,s)=g(t,y)-a(t,y)*yprime

  • adda_func must modify the provided matrix p and is a required option. If banded storage is used,

ml and mu provide the upper and lower diagonals, see the lsodi.f source for full documentation. adda_func is passed to the dae via set_options with the adda_func keyword argument. Schematically, it is

def adda(t, y, ml, mu, p, nrowp):

p += a(t,y) return p

Note: if your residual is a yprime - g, then adda must substract a, see def of rfn !

  • jacfn is not a supported option for lsodi.

  • compute_initcond: None or ‘yp0’ LSODI may be able to compute the initial conditions if you do not know them precisely and the problem is not algebraic at t=0. If yp0, then the differential variables (y of the ode system at time 0)

    will be used to solve for the derivatives of the differential variables, so yp0 will be calculated.

Source: http://www.netlib.org/ode/lsodi.f

class scikits.odes.lsodiint.lsodi(resfn, **options)[source]
init_step(t0, y0, yp0, y_ic0_retn=None, yp_ic0_retn=None)[source]

See dae.DaeBase

set_options(**options)[source]

See dae.DaeBase

solve(tspan, y0, yp0, hook_fn=None)[source]

See dae.DaeBase

step(t, y_retn, yp_retn=None)[source]

See dae.DaeBase

scikits.odes.sundials

exception scikits.odes.sundials.CVODESolveException(soln)[source]

Base class for exceptions raised by CVODE.validate_flags.

exception scikits.odes.sundials.CVODESolveFailed(soln)[source]

CVODE.solve failed to reach endpoint

exception scikits.odes.sundials.CVODESolveFoundRoot(soln)[source]

CVODE.solve found a root

exception scikits.odes.sundials.CVODESolveReachedTSTOP(soln)[source]

CVODE.solve reached the endpoint specified by tstop.

exception scikits.odes.sundials.IDASolveException(soln)[source]

Base class for exceptions raised by IDA.validate_flags.

exception scikits.odes.sundials.IDASolveFailed(soln)[source]

IDA.solve failed to reach endpoint

exception scikits.odes.sundials.IDASolveFoundRoot(soln)[source]

IDA.solve found a root

exception scikits.odes.sundials.IDASolveReachedTSTOP(soln)[source]

IDA.solve reached the endpoint specified by tstop.

scikits.odes.sundials.cvode

class scikits.odes.sundials.cvode.CV_ContinuationFunction

Simple wrapper for functions called when ROOT or TSTOP are returned.

evaluate(self, DTYPE_t t, ndarray y, CVODE solver) → int
class scikits.odes.sundials.cvode.CV_JacRhsFunction

Prototype for jacobian function.

Note that evaluate must return a integer, 0 for success, positive for recoverable failure, negative for unrecoverable failure (as per CVODE documentation).

evaluate(self, DTYPE_t t, ndarray y, ndarray fy, ndarray J) → int
Returns the Jacobi matrix of the right hand side function, as

d(rhs)/d y

(for dense the full matrix, for band only bands). Result has to be stored in the variable J, which is preallocated to the corresponding size.

This is a generic class, you should subclass is for the problem specific purposes.

class scikits.odes.sundials.cvode.CV_JacTimesSetupFunction

Prototype for jacobian times setup function.

Note that evaluate must return a integer, 0 for success, non-zero for error (as per CVODE documentation), with >0 a recoverable error (step is retried).

evaluate(self, DTYPE_t t, ndarray y, ndarray fy, userdata=None) → int

This function calculates the product of the Jacobian with a given vector v. Use the userdata object to expose Jacobian related data to the solve function.

This is a generic class, you should subclass it for the problem specific purposes.

class scikits.odes.sundials.cvode.CV_JacTimesVecFunction

Prototype for jacobian times vector function.

Note that evaluate must return a integer, 0 for success, non-zero for error (as per CVODE documentation).

evaluate(self, ndarray v, ndarray Jv, DTYPE_t t, ndarray y, userdata=None) → int

This function calculates the product of the Jacobian with a given vector v. Use the userdata object to expose Jacobian related data to the solve function.

This is a generic class, you should subclass it for the problem specific purposes.

class scikits.odes.sundials.cvode.CV_PrecSetupFunction

Prototype for preconditioning setup function.

Note that evaluate must return a integer, 0 for success, positive for recoverable failure, negative for unrecoverable failure (as per CVODE documentation).

evaluate(self, DTYPE_t t, ndarray y, bool jok, jcurPtr, DTYPE_t gamma, userdata=None) → int

This function preprocesses and/or evaluates Jacobian-related data needed by the preconditioner. Use the userdata object to expose the preprocessed data to the solve function.

This is a generic class, you should subclass it for the problem specific purposes.

class scikits.odes.sundials.cvode.CV_PrecSolveFunction

Prototype for precondititioning solution function.

Note that evaluate must return a integer, 0 for success, positive for recoverable failure, negative for unrecoverable failure (as per CVODE documentation).

evaluate(self, DTYPE_t t, ndarray y, ndarray r, ndarray z, DTYPE_t gamma, DTYPE_t delta, int lr, userdata=None) → int

This function solves the preconditioned system P*z = r, where P may be either a left or right preconditioner matrix. Here P should approximate (at least crudely) the Newton matrix M = I − gamma*J, where J is the Jacobian of the system. If preconditioning is done on both sides, the product of the two preconditioner matrices should approximate M.

This is a generic class, you should subclass it for the problem specific purposes.

class scikits.odes.sundials.cvode.CV_RhsFunction

Prototype for rhs function.

Note that evaluate must return a integer, 0 for success, positive for recoverable failure, negative for unrecoverable failure (as per CVODE documentation).

evaluate(self, DTYPE_t t, ndarray y, ndarray ydot, userdata=None) → int
class scikits.odes.sundials.cvode.CV_RootFunction

Prototype for root function.

Note that evaluate must return a integer, 0 for success, non-zero if error (as per CVODE documentation).

evaluate(self, DTYPE_t t, ndarray y, ndarray g, userdata=None) → int
class scikits.odes.sundials.cvode.CV_WrapJacRhsFunction
evaluate(self, DTYPE_t t, ndarray y, ndarray fy, ndarray J) → int

Returns the Jacobi matrix (for dense the full matrix, for band only bands. Result has to be stored in the variable J, which is preallocated to the corresponding size.

set_jacfn(self, jacfn)

Set some jacobian equations as a JacRhsFunction executable class.

class scikits.odes.sundials.cvode.CV_WrapJacTimesSetupFunction
evaluate(self, DTYPE_t t, ndarray y, ndarray fy, userdata=None) → int
set_jac_times_setupfn(self, jac_times_setupfn)

Set some CV_JacTimesSetupFn executable class.

class scikits.odes.sundials.cvode.CV_WrapJacTimesVecFunction
evaluate(self, ndarray v, ndarray Jv, DTYPE_t t, ndarray y, userdata=None) → int
set_jac_times_vecfn(self, jac_times_vecfn)

Set some CV_JacTimesVecFn executable class.

class scikits.odes.sundials.cvode.CV_WrapPrecSetupFunction
evaluate(self, DTYPE_t t, ndarray y, bool jok, jcurPtr, DTYPE_t gamma, userdata=None) → int
set_prec_setupfn(self, prec_setupfn)

set a precondititioning setup method as a CV_PrecSetupFunction executable class

class scikits.odes.sundials.cvode.CV_WrapPrecSolveFunction
evaluate(self, DTYPE_t t, ndarray y, ndarray r, ndarray z, DTYPE_t gamma, DTYPE_t delta, int lr, userdata=None) → int
set_prec_solvefn(self, prec_solvefn)

set a precondititioning solve method as a CV_PrecSolveFunction executable class

class scikits.odes.sundials.cvode.CV_WrapRhsFunction
evaluate(self, DTYPE_t t, ndarray y, ndarray ydot, userdata=None) → int
set_rhsfn(self, rhsfn)

set some rhs equations as a RhsFunction executable class

with_userdata

‘int’

Type

with_userdata

class scikits.odes.sundials.cvode.CV_WrapRootFunction
evaluate(self, DTYPE_t t, ndarray y, ndarray g, userdata=None) → int
set_rootfn(self, rootfn)

set root-ing condition(equations) as a RootFunction executable class

class scikits.odes.sundials.cvode.StatusEnum

An enumeration.

scikits.odes.sundials.cvode.no_continue_fn(t, y, solver)

scikits.odes.sundials.ida

class scikits.odes.sundials.ida.IDA_ContinuationFunction

Simple wrapper for functions called when ROOT or TSTOP are returned.

evaluate(self, DTYPE_t t, ndarray y, ndarray yp, IDA solver) → int
class scikits.odes.sundials.ida.IDA_JacRhsFunction

Prototype for jacobian function.

Note that evaluate must return a integer, 0 for success, positive for recoverable failure, negative for unrecoverable failure (as per IDA documentation).

evaluate(self, DTYPE_t t, ndarray y, ndarray ydot, ndarray residual, DTYPE_t cj, ndarray J) → int
Returns the Jacobi matrix of the residual function, as

d(res)/d y + cj d(res)/d ydot

(for dense the full matrix, for band only bands). Result has to be stored in the variable J, which is preallocated to the corresponding size.

This is a generic class, you should subclass is for the problem specific purposes.”

class scikits.odes.sundials.ida.IDA_JacTimesSetupFunction

Prototype for jacobian times setup function.

Note that evaluate must return a integer, 0 for success, non-zero for error (as per CVODE documentation), with >0 a recoverable error (step is retried).

evaluate(self, DTYPE_t tt, ndarray yy, ndarray yp, ndarray rr, DTYPE_t cj, userdata=None) → int

This function calculates the product of the Jacobian with a given vector v. Use the userdata object to expose Jacobian related data to the solve function.

This is a generic class, you should subclass it for the problem specific purposes.

class scikits.odes.sundials.ida.IDA_JacTimesVecFunction

Prototype for jacobian times vector function.

Note that evaluate must return a integer, 0 for success, non-zero for error (as per IDA documentation).

evaluate(self, DTYPE_t t, ndarray yy, ndarray yp, ndarray rr, ndarray v, ndarray Jv, DTYPE_t cj, userdata=None) → int

This function calculates the product of the Jacobian with a given vector v. Use the userdata object to expose Jacobian related data to the solve function.

This is a generic class, you should subclass it for the problem specific purposes.

class scikits.odes.sundials.ida.IDA_PrecSetupFunction

Prototype for preconditioning setup function.

Note that evaluate must return a integer, 0 for success, positive for recoverable failure, negative for unrecoverable failure (as per CVODE documentation).

evaluate(self, DTYPE_t t, ndarray y, ndarray yp, ndarray rr, DTYPE_t cj, userdata=None) → int

This function preprocesses and/or evaluates Jacobian-related data needed by the preconditioner. Use the userdata object to expose the preprocessed data to the solve function.

This is a generic class, you should subclass it for the problem specific purposes.

class scikits.odes.sundials.ida.IDA_PrecSolveFunction

Prototype for precondititioning solution function.

Note that evaluate must return a integer, 0 for success, positive for recoverable failure, negative for unrecoverable failure (as per CVODE documentation).

evaluate(self, DTYPE_t t, ndarray y, ndarray yp, ndarray r, ndarray rvec, ndarray z, DTYPE_t cj, DTYPE_t delta, userdata=None) → int

This function solves the preconditioned system P*z = r, where P may be either a left or right preconditioner matrix. Here P should approximate (at least crudely) the Newton matrix M = I − gamma*J, where J is the Jacobian of the system. If preconditioning is done on both sides, the product of the two preconditioner matrices should approximate M.

This is a generic class, you should subclass it for the problem specific purposes.

class scikits.odes.sundials.ida.IDA_RhsFunction

Prototype for rhs function.

Note that evaluate must return a integer, 0 for success, positive for recoverable failure, negative for unrecoverable failure (as per IDA documentation).

evaluate(self, DTYPE_t t, ndarray y, ndarray ydot, ndarray result, userdata=None) → int
class scikits.odes.sundials.ida.IDA_RootFunction

Prototype for root function.

Note that evaluate must return a integer, 0 for success, non-zero for error (as per IDA documentation).

evaluate(self, DTYPE_t t, ndarray y, ndarray ydot, ndarray g, userdata=None) → int
class scikits.odes.sundials.ida.IDA_WrapJacRhsFunction
evaluate(self, DTYPE_t t, ndarray y, ndarray ydot, ndarray residual, DTYPE_t cj, ndarray J) → int

Returns the Jacobi matrix (for dense the full matrix, for band only bands. Result has to be stored in the variable J, which is preallocated to the corresponding size.

set_jacfn(self, jacfn)

Set some jacobian equations as a JacResFunction executable class.

class scikits.odes.sundials.ida.IDA_WrapJacTimesSetupFunction
evaluate(self, DTYPE_t tt, ndarray yy, ndarray yp, ndarray rr, DTYPE_t cj, userdata=None) → int
set_jac_times_setupfn(self, jac_times_setupfn)

Set some IDA_JacTimesSetupFn executable class.

class scikits.odes.sundials.ida.IDA_WrapJacTimesVecFunction
evaluate(self, DTYPE_t t, ndarray yy, ndarray yp, ndarray rr, ndarray v, ndarray Jv, DTYPE_t cj, userdata=None) → int
set_jac_times_vecfn(self, jac_times_vecfn)

Set some IDA_JacTimesVecFn executable class.

class scikits.odes.sundials.ida.IDA_WrapPrecSetupFunction
evaluate(self, DTYPE_t t, ndarray y, ndarray yp, ndarray rr, DTYPE_t cj, userdata=None) → int
set_prec_setupfn(self, prec_setupfn)

set a precondititioning setup method as a IDA_PrecSetupFunction executable class

class scikits.odes.sundials.ida.IDA_WrapPrecSolveFunction
evaluate(self, DTYPE_t t, ndarray y, ndarray yp, ndarray r, ndarray rvec, ndarray z, DTYPE_t cj, DTYPE_t delta, userdata=None) → int
set_prec_solvefn(self, prec_solvefn)

set a precondititioning solve method as a IDA_PrecSolveFunction executable class

class scikits.odes.sundials.ida.IDA_WrapRhsFunction
evaluate(self, DTYPE_t t, ndarray y, ndarray ydot, ndarray result, userdata=None) → int
set_resfn(self, resfn)

set some residual equations as a ResFunction executable class

class scikits.odes.sundials.ida.IDA_WrapRootFunction
evaluate(self, DTYPE_t t, ndarray y, ndarray ydot, ndarray g, userdata=None) → int
set_rootfn(self, rootfn)

set root-ing condition(equations) as a RootFunction executable class

class scikits.odes.sundials.ida.StatusEnumIDA

An enumeration.

scikits.odes.sundials.ida.no_continue_fn(t, y, yp, solver)