summaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorDominic Steinitz <dominic@steinitz.org>2018-04-23 16:19:22 +0100
committerDominic Steinitz <dominic@steinitz.org>2018-04-23 16:19:22 +0100
commit9f571e009bc46c26334be8b6a635db1e1d5b0341 (patch)
tree38346a810d236bb2382bce55196c9f0f392f2e3e
parentdd96a98207dbafbf81b4a5f02613963cf5bd4b4c (diff)
Start of support for CVODE
-rw-r--r--packages/sundials/src/Numeric/Sundials/CVode/ODE.hs456
1 files changed, 456 insertions, 0 deletions
diff --git a/packages/sundials/src/Numeric/Sundials/CVode/ODE.hs b/packages/sundials/src/Numeric/Sundials/CVode/ODE.hs
new file mode 100644
index 0000000..f75d91f
--- /dev/null
+++ b/packages/sundials/src/Numeric/Sundials/CVode/ODE.hs
@@ -0,0 +1,456 @@
1{-# OPTIONS_GHC -Wall #-}
2
3{-# LANGUAGE QuasiQuotes #-}
4{-# LANGUAGE TemplateHaskell #-}
5{-# LANGUAGE MultiWayIf #-}
6{-# LANGUAGE OverloadedStrings #-}
7{-# LANGUAGE ScopedTypeVariables #-}
8
9-----------------------------------------------------------------------------
10-- |
11-- Module : Numeric.Sundials.CVode.ODE
12-- Copyright : Dominic Steinitz 2018,
13-- Novadiscovery 2018
14-- License : BSD
15-- Maintainer : Dominic Steinitz
16-- Stability : provisional
17--
18-- Solution of ordinary differential equation (ODE) initial value problems.
19--
20-- <https://computation.llnl.gov/projects/sundials/sundials-software>
21--
22-- A simple example:
23--
24-- <<diagrams/brusselator.png#diagram=brusselator&height=400&width=500>>
25--
26-- @
27-- import Numeric.Sundials.CVode.ODE
28-- import Numeric.LinearAlgebra
29--
30-- import Plots as P
31-- import qualified Diagrams.Prelude as D
32-- import Diagrams.Backend.Rasterific
33--
34-- brusselator :: Double -> [Double] -> [Double]
35-- brusselator _t x = [ a - (w + 1) * u + v * u * u
36-- , w * u - v * u * u
37-- , (b - w) / eps - w * u
38-- ]
39-- where
40-- a = 1.0
41-- b = 3.5
42-- eps = 5.0e-6
43-- u = x !! 0
44-- v = x !! 1
45-- w = x !! 2
46--
47-- lSaxis :: [[Double]] -> P.Axis B D.V2 Double
48-- lSaxis xs = P.r2Axis &~ do
49-- let ts = xs!!0
50-- us = xs!!1
51-- vs = xs!!2
52-- ws = xs!!3
53-- P.linePlot' $ zip ts us
54-- P.linePlot' $ zip ts vs
55-- P.linePlot' $ zip ts ws
56--
57-- main = do
58-- let res1 = odeSolve brusselator [1.2, 3.1, 3.0] (fromList [0.0, 0.1 .. 10.0])
59-- renderRasterific "diagrams/brusselator.png"
60-- (D.dims2D 500.0 500.0)
61-- (renderAxis $ lSaxis $ [0.0, 0.1 .. 10.0]:(toLists $ tr res1))
62-- @
63--
64-- KVAERNO_4_2_3
65--
66-- \[
67-- \begin{array}{c|cccc}
68-- 0.0 & 0.0 & 0.0 & 0.0 & 0.0 \\
69-- 0.871733043 & 0.4358665215 & 0.4358665215 & 0.0 & 0.0 \\
70-- 1.0 & 0.490563388419108 & 7.3570090080892e-2 & 0.4358665215 & 0.0 \\
71-- 1.0 & 0.308809969973036 & 1.490563388254106 & -1.235239879727145 & 0.4358665215 \\
72-- \hline
73-- & 0.308809969973036 & 1.490563388254106 & -1.235239879727145 & 0.4358665215 \\
74-- & 0.490563388419108 & 7.3570090080892e-2 & 0.4358665215 & 0.0 \\
75-- \end{array}
76-- \]
77--
78-- SDIRK_2_1_2
79--
80-- \[
81-- \begin{array}{c|cc}
82-- 1.0 & 1.0 & 0.0 \\
83-- 0.0 & -1.0 & 1.0 \\
84-- \hline
85-- & 0.5 & 0.5 \\
86-- & 1.0 & 0.0 \\
87-- \end{array}
88-- \]
89--
90-- SDIRK_5_3_4
91--
92-- \[
93-- \begin{array}{c|ccccc}
94-- 0.25 & 0.25 & 0.0 & 0.0 & 0.0 & 0.0 \\
95-- 0.75 & 0.5 & 0.25 & 0.0 & 0.0 & 0.0 \\
96-- 0.55 & 0.34 & -4.0e-2 & 0.25 & 0.0 & 0.0 \\
97-- 0.5 & 0.2727941176470588 & -5.036764705882353e-2 & 2.7573529411764705e-2 & 0.25 & 0.0 \\
98-- 1.0 & 1.0416666666666667 & -1.0208333333333333 & 7.8125 & -7.083333333333333 & 0.25 \\
99-- \hline
100-- & 1.0416666666666667 & -1.0208333333333333 & 7.8125 & -7.083333333333333 & 0.25 \\
101-- & 1.2291666666666667 & -0.17708333333333334 & 7.03125 & -7.083333333333333 & 0.0 \\
102-- \end{array}
103-- \]
104-----------------------------------------------------------------------------
105module Numeric.Sundials.CVode.ODE ( odeSolve
106 , odeSolveV
107 , odeSolveVWith
108 , odeSolveVWith'
109 , ODEMethod(..)
110 , StepControl(..)
111 , Jacobian
112 , SundialsDiagnostics(..)
113 ) where
114
115import qualified Language.C.Inline as C
116import qualified Language.C.Inline.Unsafe as CU
117
118import Data.Monoid ((<>))
119import Data.Maybe (isJust)
120
121import Foreign.C.Types
122import Foreign.Ptr (Ptr)
123import Foreign.ForeignPtr (newForeignPtr_)
124import Foreign.Storable (Storable)
125
126import qualified Data.Vector.Storable as V
127import qualified Data.Vector.Storable.Mutable as VM
128
129import Data.Coerce (coerce)
130import System.IO.Unsafe (unsafePerformIO)
131
132import Numeric.LinearAlgebra.Devel (createVector)
133
134import Numeric.LinearAlgebra.HMatrix (Vector, Matrix, toList, (><),
135 subMatrix, rows, cols, toLists,
136 size, subVector)
137
138import qualified Types as T
139import Arkode (cV_ADAMS, cV_BDF)
140import qualified Arkode as B
141
142
143C.context (C.baseCtx <> C.vecCtx <> C.funCtx <> T.sunCtx)
144
145C.include "<stdlib.h>"
146C.include "<stdio.h>"
147C.include "<math.h>"
148C.include "<cvode/cvode.h>" -- prototypes for CVODE fcts., consts.
149C.include "<nvector/nvector_serial.h>" -- serial N_Vector types, fcts., macros
150C.include "<sunmatrix/sunmatrix_dense.h>" -- access to dense SUNMatrix
151C.include "<sunlinsol/sunlinsol_dense.h>" -- access to dense SUNLinearSolver
152C.include "<cvode/cvode_direct.h>" -- access to CVDls interface
153C.include "<sundials/sundials_types.h>" -- definition of type realtype
154C.include "<sundials/sundials_math.h>"
155C.include "../../../helpers.h"
156C.include "Arkode_hsc.h"
157
158
159getDataFromContents :: Int -> Ptr T.SunVector -> IO (V.Vector CDouble)
160getDataFromContents len ptr = do
161 qtr <- B.getContentPtr ptr
162 rtr <- B.getData qtr
163 vectorFromC len rtr
164
165-- FIXME: Potentially an instance of Storable
166_getMatrixDataFromContents :: Ptr T.SunMatrix -> IO T.SunMatrix
167_getMatrixDataFromContents ptr = do
168 qtr <- B.getContentMatrixPtr ptr
169 rs <- B.getNRows qtr
170 cs <- B.getNCols qtr
171 rtr <- B.getMatrixData qtr
172 vs <- vectorFromC (fromIntegral $ rs * cs) rtr
173 return $ T.SunMatrix { T.rows = rs, T.cols = cs, T.vals = vs }
174
175putMatrixDataFromContents :: T.SunMatrix -> Ptr T.SunMatrix -> IO ()
176putMatrixDataFromContents mat ptr = do
177 let rs = T.rows mat
178 cs = T.cols mat
179 vs = T.vals mat
180 qtr <- B.getContentMatrixPtr ptr
181 B.putNRows rs qtr
182 B.putNCols cs qtr
183 rtr <- B.getMatrixData qtr
184 vectorToC vs (fromIntegral $ rs * cs) rtr
185-- FIXME: END
186
187putDataInContents :: Storable a => V.Vector a -> Int -> Ptr b -> IO ()
188putDataInContents vec len ptr = do
189 qtr <- B.getContentPtr ptr
190 rtr <- B.getData qtr
191 vectorToC vec len rtr
192
193-- Utils
194
195vectorFromC :: Storable a => Int -> Ptr a -> IO (V.Vector a)
196vectorFromC len ptr = do
197 ptr' <- newForeignPtr_ ptr
198 V.freeze $ VM.unsafeFromForeignPtr0 ptr' len
199
200vectorToC :: Storable a => V.Vector a -> Int -> Ptr a -> IO ()
201vectorToC vec len ptr = do
202 ptr' <- newForeignPtr_ ptr
203 V.copy (VM.unsafeFromForeignPtr0 ptr' len) vec
204
205data SundialsDiagnostics = SundialsDiagnostics {
206 aRKodeGetNumSteps :: Int
207 , aRKodeGetNumStepAttempts :: Int
208 , aRKodeGetNumRhsEvals_fe :: Int
209 , aRKodeGetNumRhsEvals_fi :: Int
210 , aRKodeGetNumLinSolvSetups :: Int
211 , aRKodeGetNumErrTestFails :: Int
212 , aRKodeGetNumNonlinSolvIters :: Int
213 , aRKodeGetNumNonlinSolvConvFails :: Int
214 , aRKDlsGetNumJacEvals :: Int
215 , aRKDlsGetNumRhsEvals :: Int
216 } deriving Show
217
218type Jacobian = Double -> Vector Double -> Matrix Double
219
220-- | Stepping functions
221data ODEMethod = ADAMS
222 | BDF
223
224getMethod :: ODEMethod -> Int
225getMethod (ADAMS) = cV_ADAMS
226getMethod (BDF) = cV_BDF
227
228getJacobian :: ODEMethod -> Maybe Jacobian
229getJacobian _ = Nothing
230
231-- | A version of 'odeSolveVWith' with reasonable default step control.
232odeSolveV
233 :: ODEMethod
234 -> Maybe Double -- ^ initial step size - by default, ARKode
235 -- estimates the initial step size to be the
236 -- solution \(h\) of the equation
237 -- \(\|\frac{h^2\ddot{y}}{2}\| = 1\), where
238 -- \(\ddot{y}\) is an estimated value of the
239 -- second derivative of the solution at \(t_0\)
240 -> Double -- ^ absolute tolerance for the state vector
241 -> Double -- ^ relative tolerance for the state vector
242 -> (Double -> Vector Double -> Vector Double) -- ^ The RHS of the system \(\dot{y} = f(t,y)\)
243 -> Vector Double -- ^ initial conditions
244 -> Vector Double -- ^ desired solution times
245 -> Matrix Double -- ^ solution
246odeSolveV meth hi epsAbs epsRel f y0 ts =
247 case odeSolveVWith meth (X epsAbs epsRel) hi g y0 ts of
248 Left c -> error $ show c -- FIXME
249 -- FIXME: Can we do better than using lists?
250 Right (v, d) -> (nR >< nC) (V.toList v)
251 where
252 us = toList ts
253 nR = length us
254 nC = size y0
255 g t x0 = coerce $ f t x0
256
257-- | A version of 'odeSolveV' with reasonable default parameters and
258-- system of equations defined using lists. FIXME: we should say
259-- something about the fact we could use the Jacobian but don't for
260-- compatibility with hmatrix-gsl.
261odeSolve :: (Double -> [Double] -> [Double]) -- ^ The RHS of the system \(\dot{y} = f(t,y)\)
262 -> [Double] -- ^ initial conditions
263 -> Vector Double -- ^ desired solution times
264 -> Matrix Double -- ^ solution
265odeSolve f y0 ts =
266 -- FIXME: These tolerances are different from the ones in GSL
267 case odeSolveVWith BDF (XX' 1.0e-6 1.0e-10 1 1) Nothing g (V.fromList y0) (V.fromList $ toList ts) of
268 Left c -> error $ show c -- FIXME
269 Right (v, d) -> (nR >< nC) (V.toList v)
270 where
271 us = toList ts
272 nR = length us
273 nC = length y0
274 g t x0 = V.fromList $ f t (V.toList x0)
275
276odeSolveVWith' ::
277 ODEMethod
278 -> StepControl
279 -> Maybe Double -- ^ initial step size - by default, ARKode
280 -- estimates the initial step size to be the
281 -- solution \(h\) of the equation
282 -- \(\|\frac{h^2\ddot{y}}{2}\| = 1\), where
283 -- \(\ddot{y}\) is an estimated value of the second
284 -- derivative of the solution at \(t_0\)
285 -> (Double -> V.Vector Double -> V.Vector Double) -- ^ The RHS of the system \(\dot{y} = f(t,y)\)
286 -> V.Vector Double -- ^ Initial conditions
287 -> V.Vector Double -- ^ Desired solution times
288 -> Matrix Double -- ^ Error code or solution
289odeSolveVWith' method control initStepSize f y0 tt =
290 case odeSolveVWith method control initStepSize f y0 tt of
291 Left c -> error $ show c -- FIXME
292 Right (v, _d) -> (nR >< nC) (V.toList v)
293 where
294 nR = V.length tt
295 nC = V.length y0
296
297odeSolveVWith ::
298 ODEMethod
299 -> StepControl
300 -> Maybe Double -- ^ initial step size - by default, ARKode
301 -- estimates the initial step size to be the
302 -- solution \(h\) of the equation
303 -- \(\|\frac{h^2\ddot{y}}{2}\| = 1\), where
304 -- \(\ddot{y}\) is an estimated value of the second
305 -- derivative of the solution at \(t_0\)
306 -> (Double -> V.Vector Double -> V.Vector Double) -- ^ The RHS of the system \(\dot{y} = f(t,y)\)
307 -> V.Vector Double -- ^ Initial conditions
308 -> V.Vector Double -- ^ Desired solution times
309 -> Either Int ((V.Vector Double), SundialsDiagnostics) -- ^ Error code or solution
310odeSolveVWith method control initStepSize f y0 tt =
311 case solveOdeC (fromIntegral $ getMethod method) (coerce initStepSize) jacH (scise control)
312 (coerce f) (coerce y0) (coerce tt) of
313 Left c -> Left $ fromIntegral c
314 Right (v, d) -> Right (coerce v, d)
315 where
316 l = size y0
317 scise (X absTol relTol) = coerce (V.replicate l absTol, relTol)
318 scise (X' absTol relTol) = coerce (V.replicate l absTol, relTol)
319 scise (XX' absTol relTol yScale _yDotScale) = coerce (V.replicate l absTol, yScale * relTol)
320 -- FIXME; Should we check that the length of ss is correct?
321 scise (ScXX' absTol relTol yScale _yDotScale ss) = coerce (V.map (* absTol) ss, yScale * relTol)
322 jacH = fmap (\g t v -> matrixToSunMatrix $ g (coerce t) (coerce v)) $
323 getJacobian method
324 matrixToSunMatrix m = T.SunMatrix { T.rows = nr, T.cols = nc, T.vals = vs }
325 where
326 nr = fromIntegral $ rows m
327 nc = fromIntegral $ cols m
328 -- FIXME: efficiency
329 vs = V.fromList $ map coerce $ concat $ toLists m
330
331solveOdeC ::
332 CInt ->
333 Maybe CDouble ->
334 (Maybe (CDouble -> V.Vector CDouble -> T.SunMatrix)) ->
335 (V.Vector CDouble, CDouble) ->
336 (CDouble -> V.Vector CDouble -> V.Vector CDouble) -- ^ The RHS of the system \(\dot{y} = f(t,y)\)
337 -> V.Vector CDouble -- ^ Initial conditions
338 -> V.Vector CDouble -- ^ Desired solution times
339 -> Either CInt ((V.Vector CDouble), SundialsDiagnostics) -- ^ Error code or solution
340solveOdeC method initStepSize jacH (absTols, relTol) fun f0 ts = unsafePerformIO $ do
341
342 let isInitStepSize :: CInt
343 isInitStepSize = fromIntegral $ fromEnum $ isJust initStepSize
344 ss :: CDouble
345 ss = case initStepSize of
346 -- It would be better to put an error message here but
347 -- inline-c seems to evaluate this even if it is never
348 -- used :(
349 Nothing -> 0.0
350 Just x -> x
351 let dim = V.length f0
352 nEq :: CLong
353 nEq = fromIntegral dim
354 nTs :: CInt
355 nTs = fromIntegral $ V.length ts
356 -- FIXME: fMut is not actually mutatated
357 fMut <- V.thaw f0
358 tMut <- V.thaw ts
359 -- FIXME: I believe this gets taken from the ghc heap and so should
360 -- be subject to garbage collection.
361 -- quasiMatrixRes <- createVector ((fromIntegral dim) * (fromIntegral nTs))
362 -- qMatMut <- V.thaw quasiMatrixRes
363 diagnostics :: V.Vector CLong <- createVector 10 -- FIXME
364 diagMut <- V.thaw diagnostics
365 -- We need the types that sundials expects. These are tied together
366 -- in 'Types'. FIXME: The Haskell type is currently empty!
367 let funIO :: CDouble -> Ptr T.SunVector -> Ptr T.SunVector -> Ptr () -> IO CInt
368 funIO x y f _ptr = do
369 -- Convert the pointer we get from C (y) to a vector, and then
370 -- apply the user-supplied function.
371 fImm <- fun x <$> getDataFromContents dim y
372 -- Fill in the provided pointer with the resulting vector.
373 putDataInContents fImm dim f
374 -- FIXME: I don't understand what this comment means
375 -- Unsafe since the function will be called many times.
376 [CU.exp| int{ 0 } |]
377 let isJac :: CInt
378 isJac = fromIntegral $ fromEnum $ isJust jacH
379 jacIO :: CDouble -> Ptr T.SunVector -> Ptr T.SunVector -> Ptr T.SunMatrix ->
380 Ptr () -> Ptr T.SunVector -> Ptr T.SunVector -> Ptr T.SunVector ->
381 IO CInt
382 jacIO t y _fy jacS _ptr _tmp1 _tmp2 _tmp3 = do
383 case jacH of
384 Nothing -> error "Numeric.Sundials.ARKode.ODE: Jacobian not defined"
385 Just jacI -> do j <- jacI t <$> getDataFromContents dim y
386 putMatrixDataFromContents j jacS
387 -- FIXME: I don't understand what this comment means
388 -- Unsafe since the function will be called many times.
389 [CU.exp| int{ 0 } |]
390
391 res <- [C.block| int {
392 /* general problem variables */
393
394 int flag; /* reusable error-checking flag */
395 int i, j; /* reusable loop indices */
396 N_Vector y = NULL; /* empty vector for storing solution */
397 void *cvode_mem = NULL; /* empty CVODE memory structure */
398
399 /* general problem parameters */
400
401 realtype T0 = RCONST(($vec-ptr:(double *ts))[0]); /* initial time */
402 sunindextype NEQ = $(sunindextype nEq); /* number of dependent vars. */
403
404 /* Initialize data structures */
405
406 y = N_VNew_Serial(NEQ); /* Create serial vector for solution */
407 if (check_flag((void *)y, "N_VNew_Serial", 0)) return 1;
408 /* Specify initial condition */
409 for (i = 0; i < NEQ; i++) {
410 NV_Ith_S(y,i) = ($vec-ptr:(double *f0))[i];
411 };
412
413 cvode_mem = CVodeCreate(CV_BDF, CV_NEWTON);
414 if (check_flag((void *)cvode_mem, "CVodeCreate", 0)) return(1);
415
416 /* Call CVodeInit to initialize the integrator memory and specify the
417 * user's right hand side function in y'=f(t,y), the inital time T0, and
418 * the initial dependent variable vector y. */
419 flag = CVodeInit(cvode_mem, $fun:(int (* funIO) (double t, SunVector y[], SunVector dydt[], void * params)), T0, y);
420 if (check_flag(&flag, "CVodeInit", 1)) return(1);
421
422 /* Clean up and return */
423
424 N_VDestroy(y); /* Free y vector */
425 CVodeFree(&cvode_mem); /* Free integrator memory */
426
427 return flag;
428 } |]
429 if res == 0
430 then do
431 return $ Left res
432 else do
433 return $ Left res
434
435-- | Adaptive step-size control
436-- functions.
437--
438-- [GSL](https://www.gnu.org/software/gsl/doc/html/ode-initval.html#adaptive-step-size-control)
439-- allows the user to control the step size adjustment using
440-- \(D_i = \epsilon^{abs}s_i + \epsilon^{rel}(a_{y} |y_i| + a_{dy/dt} h |\dot{y}_i|)\) where
441-- \(\epsilon^{abs}\) is the required absolute error, \(\epsilon^{rel}\)
442-- is the required relative error, \(s_i\) is a vector of scaling
443-- factors, \(a_{y}\) is a scaling factor for the solution \(y\) and
444-- \(a_{dydt}\) is a scaling factor for the derivative of the solution \(dy/dt\).
445--
446-- [ARKode](https://computation.llnl.gov/projects/sundials/arkode)
447-- allows the user to control the step size adjustment using
448-- \(\eta^{rel}|y_i| + \eta^{abs}_i\). For compatibility with
449-- [hmatrix-gsl](https://hackage.haskell.org/package/hmatrix-gsl),
450-- tolerances for \(y\) and \(\dot{y}\) can be specified but the latter have no
451-- effect.
452data StepControl = X Double Double -- ^ absolute and relative tolerance for \(y\); in GSL terms, \(a_{y} = 1\) and \(a_{dy/dt} = 0\); in ARKode terms, the \(\eta^{abs}_i\) are identical
453 | X' Double Double -- ^ absolute and relative tolerance for \(\dot{y}\); in GSL terms, \(a_{y} = 0\) and \(a_{dy/dt} = 1\); in ARKode terms, the latter is treated as the relative tolerance for \(y\) so this is the same as specifying 'X' which may be entirely incorrect for the given problem
454 | XX' Double Double Double Double -- ^ include both via relative tolerance
455 -- scaling factors \(a_y\), \(a_{{dy}/{dt}}\); in ARKode terms, the latter is ignored and \(\eta^{rel} = a_{y}\epsilon^{rel}\)
456 | ScXX' Double Double Double Double (Vector Double) -- ^ scale absolute tolerance of \(y_i\); in ARKode terms, \(a_{{dy}/{dt}}\) is ignored, \(\eta^{abs}_i = s_i \epsilon^{abs}\) and \(\eta^{rel} = a_{y}\epsilon^{rel}\)