From eab231907c3d11651e48e15510f4510fd6c77450 Mon Sep 17 00:00:00 2001 From: Alberto Ruiz Date: Sat, 25 Feb 2012 12:40:49 +0100 Subject: merge changelog --- lib/Numeric/GSL/ODE.hs | 47 ++++++++++++++++++++++++++++++++++++----------- 1 file changed, 36 insertions(+), 11 deletions(-) (limited to 'lib/Numeric/GSL/ODE.hs') diff --git a/lib/Numeric/GSL/ODE.hs b/lib/Numeric/GSL/ODE.hs index 2251acd..c243f4b 100644 --- a/lib/Numeric/GSL/ODE.hs +++ b/lib/Numeric/GSL/ODE.hs @@ -29,7 +29,7 @@ main = mplot (ts : toColumns sol)@ ----------------------------------------------------------------------------- module Numeric.GSL.ODE ( - odeSolve, odeSolveV, ODEMethod(..) + odeSolve, odeSolveV, ODEMethod(..), Jacobian ) where import Data.Packed.Internal @@ -41,18 +41,21 @@ import System.IO.Unsafe(unsafePerformIO) ------------------------------------------------------------------------- +type Jacobian = Double -> Vector Double -> Matrix Double + -- | Stepping functions data ODEMethod = 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 'RKf45'. + | 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 -- ^ Implicit 2nd order Runge-Kutta at Gaussian points. - | RK4imp -- ^ Implicit 4th order Runge-Kutta at Gaussian points. - | BSimp -- ^ Implicit Bulirsch-Stoer method of Bader and Deuflhard. This algorithm requires the Jacobian. - | Gear1 -- ^ M=1 implicit Gear method. - | Gear2 -- ^ M=2 implicit Gear method. - deriving (Enum,Eq,Show,Bounded) + | 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. + -- | A version of 'odeSolveV' with reasonable default parameters and system of equations defined using lists. odeSolve @@ -60,7 +63,7 @@ odeSolve -> [Double] -- ^ initial conditions -> Vector Double -- ^ desired solution times -> Matrix Double -- ^ solution -odeSolve xdot xi ts = odeSolveV RKf45 hi epsAbs epsRel (l2v xdot) Nothing (fromList xi) ts +odeSolve xdot xi ts = odeSolveV RKf45 hi epsAbs epsRel (l2v xdot) (fromList xi) ts where hi = (ts@>1 - ts@>0)/100 epsAbs = 1.49012e-08 epsRel = 1.49012e-08 @@ -73,11 +76,33 @@ odeSolveV -> 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 +odeSolveV RK2 = odeSolveV' 0 Nothing +odeSolveV RK4 = odeSolveV' 1 Nothing +odeSolveV RKf45 = odeSolveV' 2 Nothing +odeSolveV RKck = odeSolveV' 3 Nothing +odeSolveV RK8pd = odeSolveV' 4 Nothing +odeSolveV (RK2imp jac) = odeSolveV' 5 (Just jac) +odeSolveV (RK4imp jac) = odeSolveV' 6 (Just jac) +odeSolveV (BSimp jac) = odeSolveV' 7 (Just jac) +odeSolveV (RK1imp jac) = odeSolveV' 8 (Just jac) +odeSolveV MSAdams = odeSolveV' 9 Nothing +odeSolveV (MSBDF jac) = odeSolveV' 10 (Just jac) + + +odeSolveV' + :: CInt -> Maybe (Double -> Vector Double -> Matrix Double) -- ^ optional jacobian + -> 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 -odeSolveV method h epsAbs epsRel f mbjac xiv ts = unsafePerformIO $ do +odeSolveV' method mbjac h epsAbs epsRel f xiv ts = unsafePerformIO $ do let n = dim xiv fp <- mkDoubleVecVecfun (\t -> aux_vTov (checkdim1 n . f t)) jp <- case mbjac of @@ -86,7 +111,7 @@ odeSolveV method h epsAbs epsRel f mbjac xiv ts = unsafePerformIO $ do sol <- vec xiv $ \xiv' -> vec (checkTimes ts) $ \ts' -> createMIO (dim ts) n - (ode_c (fi (fromEnum method)) h epsAbs epsRel fp jp // xiv' // ts' ) + (ode_c (method) h epsAbs epsRel fp jp // xiv' // ts' ) "ode" freeHaskellFunPtr fp return sol -- cgit v1.2.3