Package org.hipparchus.ode

This package provides classes to solve Ordinary Differential Equations problems.

This package solves Initial Value Problems of the form y'=f(t,y) with t0 and y(t0)=y0 known. The provided integrators compute an estimate of y(t) from t=t0 to t=t1. It is also possible to get thederivatives with respect to the initial state dy(t)/dy(t0) or the derivatives with respect to some ODE parameters dy(t)/dp.

All integrators provide dense output. This means that besides computing the state vector at discrete times, they also provide a cheap mean to get the state between the time steps. They do so through classes extending the ODEStateInterpolator abstract class, which are made available to the user at the end of each step.

All integrators handle multiple discrete events detection based on switching functions. This means that the integrator can be driven by user specified discrete events. The steps are shortened as needed to ensure the events occur at step boundaries (even if the integrator is a fixed-step integrator). When the events are triggered, integration can be stopped (this is called a G-stop facility), the state vector can be changed, or integration can simply go on. The latter case is useful to handle discontinuities in the differential equations gracefully and get accurate dense output even close to the discontinuity. See org.hipparchus.ode.events for more on how events are handled.

The user should describe his problem in his own classes (UserProblem in the diagram below) which should implement the OrdinaryDifferentialEquation interface. Then he should pass it to the integrator he prefers among all the classes that implement the ODEIntegrator interface.

The solution of the integration problem is provided by two means. The first one is aimed towards simple use: the state vector at the end of the integration process is copied in the y array of the ODEIntegrator.integrate method. The second one should be used when more in-depth information is needed throughout the integration process. The user can register an object implementing the ODEStepHandler interface or a StepNormalizer object wrapping a user-specified object implementing the ODEFixedStepHandler interface into the integrator before calling the ODEIntegrator.integrate method. The user object will be called appropriately during the integration process, allowing the user to process intermediate results. The default step handler does nothing.

DenseOutputModel is a special-purpose step handler that is able to store all steps and to provide transparent access to any intermediate result once the integration is over. An important feature of this class is that it implements the Serializable interface. This means that a complete continuous model of the integrated function throughout the integration range can be serialized and reused later (if stored into a persistent medium like a filesystem or a database) or elsewhere (if sent to another application). Only the result of the integration is stored, there is no reference to the integrated problem by itself.

Custom implementations can be developed for specific needs. As an example, if an application is to be completely driven by the integration process, then most of the application code will be run inside a step handler specific to this application.

Some integrators (the simple ones) use fixed steps that are set at creation time. The more efficient integrators use variable steps that are handled internally in order to control the integration error with respect to a specified accuracy (these integrators extend the AdaptiveStepsizeIntegrator abstract class). In this case, the step handler which is called after each successful step shows up the variable stepsize. The StepNormalizer class can be used to convert the variable stepsize into a fixed stepsize that can be handled by classes implementing the ODEFixedStepHandler interface. Adaptive stepsize integrators can automatically compute the initial stepsize by themselves, however the user can specify it if he prefers to retain full control over the integration or if the automatic guess is wrong.

Fixed Step Integrators
NameOrder
Euler1
Midpoint2
Classical Runge-Kutta4
Gill4
3/84
Luther6

Adaptive Stepsize Integrators
NameIntegration OrderError Estimation Order
Higham and Hall54
Dormand-Prince 5(4)54
Dormand-Prince 8(5,3)85 and 3
Gragg-Bulirsch-Stoervariable (up to 18 by default)variable
Adams-Bashforthvariablevariable
Adams-Moultonvariablevariable

In the table above, the Adams-Bashforth and Adams-Moulton integrators appear as variable-step ones. This is an extension to the classical algorithms using the Nordsieck vector representation.