ODASSL

Modified version of DASSL for solving overdetermined systems of (singularily) implicit ODEs. The main difference to DASSL is in the corrector iteration part.

ODASSL ad-ons : FUEHRER, CLAUS
DEUTSCHE FORSCHUNGSANSTALT
FUER LUFT- UND RAUMFAHRT (DLR)
INST. DYNAMIC DER FLUGSYSTEME
D-8031 WESSLING  (F.R.G)

Based on DASSL version dated to 900103 by:

DASSL-Author:  PETZOLD, LINDA
APPLIED MATHEMATICS DIVISION 8331
SANDIA NATIONAL LABORATORIES
LIVERMORE, CA.    94550

Support

  • State events (root funtions) : False
  • Step events (completed step) : True
  • Time events : True

Usage

Import the solver together with the correct problem:

from assimulo.solvers import ODASSL
from assimulo.problem import Overdetermined_Problem

Define the problem, such as:

def res(t, y, yd): #Note that y and yd are 1-D numpy arrays.
    res = [yd[0]-1.0, y[0]-1.0]
    return N.array([res]) #Note that the return must be numpy array, NOT a scalar.

y0  = [1.0]
yd0 = [1.0]
t0  = 1.0

Create a problem instance:

mod = Overdetermined_Problem(res, y0, yd0, t0)

Note

For complex problems, it is recommended to check the available examples and the documentation in the problem class, Overdetermined_Problem. It is also recommended to define your problem as a subclass of Overdetermined_Problem.

Warning

When subclassing from a problem class, the function for calculating the right-hand-side (for ODEs) must be named rhs and in the case with a residual function (for DAEs) it must be named res.

Create a solver instance:

sim = ODASSL(mod)

Modify (optionally) the solver parameters.

Parameters:

  • atol Defines the absolute tolerance(s) that is to be used by the solver.
  • backward Specifies if the simulation is done in reverse time.
  • clock_step Specifies if the elapsed time of an integrator step should be timed or not.
  • display_progress This option actives output during the integration in terms of that the current integration is periodically printed to the stdout.
  • inith This determines the initial step-size to be used in the integration.
  • maxh Defines the maximal step-size that is to be used by the solver.
  • maxord Defines the maximal order that is to be used by the solver.
  • maxsteps .
  • num_threads This options specifies the number of threads to be used for those solvers that supports it.
  • report_continuously This options specifies if the solver should report the solution continuously after steps.
  • rtol Defines the relative tolerance that is to be used by the solver.
  • safe .
  • store_event_points This options specifies if the solver should save additional points at the events, \(t_e^-, t_e^+\).
  • time_limit This option can be used to limit the time of an integration.
  • usejac This sets the option to use the user defined Jacobian.
  • verbosity This determines the level of the output.

Simulate the problem:

Information: