hmatrix-gsl-0.16.0.3: Numerical computation

Copyright(c) Alberto Ruiz 2010
LicenseGPL
MaintainerAlberto Ruiz
Stabilityprovisional
Safe HaskellNone
LanguageHaskell98

Numeric.GSL.ODE

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.ODE
import Numeric.LinearAlgebra
import Numeric.LinearAlgebra.Util(mplot)

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

ts = linspace 100 (0,20 :: Double)

sol = odeSolve xdot [10,0] ts

main = mplot (ts : toColumns sol)

Synopsis

Documentation

odeSolve Source

Arguments

:: (Double -> [Double] -> [Double])

xdot(t,x)

-> [Double]

initial conditions

-> Vector Double

desired solution times

-> Matrix Double

solution

A version of odeSolveV with reasonable default parameters and system of equations defined using lists.

odeSolveV Source

Arguments

:: ODEMethod 
-> Double

initial step size

-> Double

absolute tolerance for the state vector

-> Double

relative tolerance for the state vector

-> (Double -> Vector Double -> Vector Double)

xdot(t,x)

-> Vector Double

initial conditions

-> Vector Double

desired solution times

-> Matrix Double

solution

Evolution of the system with adaptive step-size control.

data ODEMethod Source

Stepping functions

Constructors

RK2

Embedded Runge-Kutta (2, 3) method.

RK4

4th order (classical) Runge-Kutta. The error estimate is obtained by halving the step-size. For more efficient estimate of the error, use the embedded methods.

RKf45

Embedded Runge-Kutta-Fehlberg (4, 5) method. This method is a good general-purpose integrator.

RKck

Embedded Runge-Kutta Cash-Karp (4, 5) method.

RK8pd

Embedded Runge-Kutta Prince-Dormand (8,9) method.

RK2imp Jacobian

Implicit 2nd order Runge-Kutta at Gaussian points.

RK4imp Jacobian

Implicit 4th order Runge-Kutta at Gaussian points.

BSimp Jacobian

Implicit Bulirsch-Stoer method of Bader and Deuflhard. The method is generally suitable for stiff problems.

RK1imp Jacobian

Implicit Gaussian first order Runge-Kutta. Also known as implicit Euler or backward Euler method. Error estimation is carried out by the step doubling method.

MSAdams

A variable-coefficient linear multistep Adams method in Nordsieck form. This stepper uses explicit Adams-Bashforth (predictor) and implicit Adams-Moulton (corrector) methods in P(EC)^m functional iteration mode. Method order varies dynamically between 1 and 12.

MSBDF Jacobian

A variable-coefficient linear multistep backward differentiation formula (BDF) method in Nordsieck form. This stepper uses the explicit BDF formula as predictor and implicit BDF formula as corrector. A modified Newton iteration method is used to solve the system of non-linear equations. Method order varies dynamically between 1 and 5. The method is generally suitable for stiff problems.