ExplicitEuler

This solver solves an explicit ordinary differential equation using the explicit Euler method.

We want to approximate the solution to the ordinary differential equation of the form,

\[\dot{y} = f(t,y), \quad y(t_0) = y_0 .\]

Using the explicit Euler method, the approximation is defined as follow,

\[y_{n+1} = y_n + hf(t_n,y_n)\]

with \(h\) being the step-size and \(y_n\) the previous solution to the equation.

Support

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

Usage

Import the solver together with the correct problem:

from assimulo.solvers import ExplicitEuler
from assimulo.problem import Explicit_Problem

Define the problem, such as:

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

y0 = [1.0]
t0 = 1.0

Create a problem instance:

mod = Explicit_Problem(rhs, y0, t0)

Note

For complex problems, it is recommended to check the available examples and the documentation in the problem class, Explicit_Problem. It is also recommended to define your problem as a subclass of Explicit_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 = ExplicitEuler(mod)

Modify (optionally) the solver parameters.

Parameters:

  • 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.
  • h Defines the step-size that is to be used by the solver.
  • 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.
  • 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.
  • verbosity This determines the level of the output.

Simulate the problem:

Information: