hmatrix-0.9.3.0: Linear algebra and numerical computationSource codeContentsIndex
Numeric.GSL.ODE
Portabilityuses ffi
Stabilityprovisional
MaintainerAlberto Ruiz (aruiz at um dot es)
Description

Solution of ordinary differential equation (ODE) initial value problems.

http://www.gnu.org/software/gsl/manual/html_node/Ordinary-Differential-Equations.html

A simple example:

import Numeric.GSL
import Numeric.LinearAlgebra
import Graphics.Plot

xdot t [x,v] = [v, -0.95*x - 0.1*v]

ts = linspace 100 (0,20)

sol = odeSolve xdot [10,0] ts

main = mplot (ts : toColumns sol)
Synopsis
odeSolve :: (Double -> [Double] -> [Double]) -> [Double] -> Vector Double -> Matrix Double
odeSolveV :: ODEMethod -> Double -> Double -> Double -> (Double -> Vector Double -> Vector Double) -> Maybe (Double -> Vector Double -> Matrix Double) -> Vector Double -> Vector Double -> Matrix Double
data ODEMethod
= RK2
| RK4
| RKf45
| RKck
| RK8pd
| RK2imp
| RK4imp
| BSimp
| Gear1
| Gear2
Documentation
odeSolveSource
:: Double -> [Double] -> [Double]xdot(t,x)
-> [Double]initial conditions
-> Vector Doubledesired solution times
-> Matrix Doublesolution
A version of odeSolveV with reasonable default parameters and system of equations defined using lists.
odeSolveVSource
:: ODEMethod
-> Doubleinitial step size
-> Doubleabsolute tolerance for the state vector
-> Doublerelative tolerance for the state vector
-> Double -> Vector Double -> Vector Doublexdot(t,x)
-> Maybe (Double -> Vector Double -> Matrix Double)optional jacobian
-> Vector Doubleinitial conditions
-> Vector Doubledesired solution times
-> Matrix Doublesolution
Evolution of the system with adaptive step-size control.
data ODEMethod Source
Stepping functions
Constructors
RK2Embedded Runge-Kutta (2, 3) method.
RK44th order (classical) Runge-Kutta. The error estimate is obtained by halving the step-size. For more efficient estimate of the error, use RKf45.
RKf45Embedded Runge-Kutta-Fehlberg (4, 5) method. This method is a good general-purpose integrator.
RKckEmbedded Runge-Kutta Cash-Karp (4, 5) method.
RK8pdEmbedded Runge-Kutta Prince-Dormand (8,9) method.
RK2impImplicit 2nd order Runge-Kutta at Gaussian points.
RK4impImplicit 4th order Runge-Kutta at Gaussian points.
BSimpImplicit Bulirsch-Stoer method of Bader and Deuflhard. This algorithm requires the Jacobian.
Gear1M=1 implicit Gear method.
Gear2M=2 implicit Gear method.
show/hide Instances
Produced by Haddock version 2.6.1