diff options
Diffstat (limited to 'packages/sundials/src/Numeric/Sundials/CVode')
-rw-r--r-- | packages/sundials/src/Numeric/Sundials/CVode/ODE.hs | 476 |
1 files changed, 476 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..a6f185e --- /dev/null +++ b/packages/sundials/src/Numeric/Sundials/CVode/ODE.hs | |||
@@ -0,0 +1,476 @@ | |||
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 | ----------------------------------------------------------------------------- | ||
65 | module Numeric.Sundials.CVode.ODE ( odeSolve | ||
66 | , odeSolveV | ||
67 | , odeSolveVWith | ||
68 | , odeSolveVWith' | ||
69 | , ODEMethod(..) | ||
70 | , StepControl(..) | ||
71 | ) where | ||
72 | |||
73 | import qualified Language.C.Inline as C | ||
74 | import qualified Language.C.Inline.Unsafe as CU | ||
75 | |||
76 | import Data.Monoid ((<>)) | ||
77 | import Data.Maybe (isJust) | ||
78 | |||
79 | import Foreign.C.Types (CDouble, CInt, CLong) | ||
80 | import Foreign.Ptr (Ptr) | ||
81 | import Foreign.Storable (poke) | ||
82 | |||
83 | import qualified Data.Vector.Storable as V | ||
84 | |||
85 | import Data.Coerce (coerce) | ||
86 | import System.IO.Unsafe (unsafePerformIO) | ||
87 | |||
88 | import Numeric.LinearAlgebra.Devel (createVector) | ||
89 | |||
90 | import Numeric.LinearAlgebra.HMatrix (Vector, Matrix, toList, rows, | ||
91 | cols, toLists, size, reshape) | ||
92 | |||
93 | import Numeric.Sundials.Arkode (cV_ADAMS, cV_BDF, | ||
94 | getDataFromContents, putDataInContents) | ||
95 | import qualified Numeric.Sundials.Arkode as T | ||
96 | import Numeric.Sundials.ODEOpts (ODEOpts(..), Jacobian, SundialsDiagnostics(..)) | ||
97 | |||
98 | |||
99 | C.context (C.baseCtx <> C.vecCtx <> C.funCtx <> T.sunCtx) | ||
100 | |||
101 | C.include "<stdlib.h>" | ||
102 | C.include "<stdio.h>" | ||
103 | C.include "<math.h>" | ||
104 | C.include "<cvode/cvode.h>" -- prototypes for CVODE fcts., consts. | ||
105 | C.include "<nvector/nvector_serial.h>" -- serial N_Vector types, fcts., macros | ||
106 | C.include "<sunmatrix/sunmatrix_dense.h>" -- access to dense SUNMatrix | ||
107 | C.include "<sunlinsol/sunlinsol_dense.h>" -- access to dense SUNLinearSolver | ||
108 | C.include "<cvode/cvode_direct.h>" -- access to CVDls interface | ||
109 | C.include "<sundials/sundials_types.h>" -- definition of type realtype | ||
110 | C.include "<sundials/sundials_math.h>" | ||
111 | C.include "../../../helpers.h" | ||
112 | C.include "Numeric/Sundials/Arkode_hsc.h" | ||
113 | |||
114 | |||
115 | -- | Stepping functions | ||
116 | data ODEMethod = ADAMS | ||
117 | | BDF | ||
118 | |||
119 | getMethod :: ODEMethod -> Int | ||
120 | getMethod (ADAMS) = cV_ADAMS | ||
121 | getMethod (BDF) = cV_BDF | ||
122 | |||
123 | getJacobian :: ODEMethod -> Maybe Jacobian | ||
124 | getJacobian _ = Nothing | ||
125 | |||
126 | -- | A version of 'odeSolveVWith' with reasonable default step control. | ||
127 | odeSolveV | ||
128 | :: ODEMethod | ||
129 | -> Maybe Double -- ^ initial step size - by default, CVode | ||
130 | -- estimates the initial step size to be the | ||
131 | -- solution \(h\) of the equation | ||
132 | -- \(\|\frac{h^2\ddot{y}}{2}\| = 1\), where | ||
133 | -- \(\ddot{y}\) is an estimated value of the | ||
134 | -- second derivative of the solution at \(t_0\) | ||
135 | -> Double -- ^ absolute tolerance for the state vector | ||
136 | -> Double -- ^ relative tolerance for the state vector | ||
137 | -> (Double -> Vector Double -> Vector Double) -- ^ The RHS of the system \(\dot{y} = f(t,y)\) | ||
138 | -> Vector Double -- ^ initial conditions | ||
139 | -> Vector Double -- ^ desired solution times | ||
140 | -> Matrix Double -- ^ solution | ||
141 | odeSolveV meth hi epsAbs epsRel f y0 ts = | ||
142 | odeSolveVWith meth (X epsAbs epsRel) hi g y0 ts | ||
143 | where | ||
144 | g t x0 = coerce $ f t x0 | ||
145 | |||
146 | -- | A version of 'odeSolveV' with reasonable default parameters and | ||
147 | -- system of equations defined using lists. FIXME: we should say | ||
148 | -- something about the fact we could use the Jacobian but don't for | ||
149 | -- compatibility with hmatrix-gsl. | ||
150 | odeSolve :: (Double -> [Double] -> [Double]) -- ^ The RHS of the system \(\dot{y} = f(t,y)\) | ||
151 | -> [Double] -- ^ initial conditions | ||
152 | -> Vector Double -- ^ desired solution times | ||
153 | -> Matrix Double -- ^ solution | ||
154 | odeSolve f y0 ts = | ||
155 | -- FIXME: These tolerances are different from the ones in GSL | ||
156 | odeSolveVWith BDF (XX' 1.0e-6 1.0e-10 1 1) Nothing g (V.fromList y0) (V.fromList $ toList ts) | ||
157 | where | ||
158 | g t x0 = V.fromList $ f t (V.toList x0) | ||
159 | |||
160 | odeSolveVWith :: | ||
161 | ODEMethod | ||
162 | -> StepControl | ||
163 | -> Maybe Double -- ^ initial step size - by default, CVode | ||
164 | -- estimates the initial step size to be the | ||
165 | -- solution \(h\) of the equation | ||
166 | -- \(\|\frac{h^2\ddot{y}}{2}\| = 1\), where | ||
167 | -- \(\ddot{y}\) is an estimated value of the second | ||
168 | -- derivative of the solution at \(t_0\) | ||
169 | -> (Double -> V.Vector Double -> V.Vector Double) -- ^ The RHS of the system \(\dot{y} = f(t,y)\) | ||
170 | -> V.Vector Double -- ^ Initial conditions | ||
171 | -> V.Vector Double -- ^ Desired solution times | ||
172 | -> Matrix Double -- ^ Error code or solution | ||
173 | odeSolveVWith method control initStepSize f y0 tt = | ||
174 | case odeSolveVWith' opts method control initStepSize f y0 tt of | ||
175 | Left c -> error $ show c -- FIXME | ||
176 | Right (v, _d) -> v | ||
177 | where | ||
178 | opts = ODEOpts { maxNumSteps = 10000 | ||
179 | , minStep = 1.0e-12 | ||
180 | , relTol = error "relTol" | ||
181 | , absTols = error "absTol" | ||
182 | , initStep = error "initStep" | ||
183 | , maxFail = 10 | ||
184 | } | ||
185 | |||
186 | odeSolveVWith' :: | ||
187 | ODEOpts | ||
188 | -> ODEMethod | ||
189 | -> StepControl | ||
190 | -> Maybe Double -- ^ initial step size - by default, CVode | ||
191 | -- estimates the initial step size to be the | ||
192 | -- solution \(h\) of the equation | ||
193 | -- \(\|\frac{h^2\ddot{y}}{2}\| = 1\), where | ||
194 | -- \(\ddot{y}\) is an estimated value of the second | ||
195 | -- derivative of the solution at \(t_0\) | ||
196 | -> (Double -> V.Vector Double -> V.Vector Double) -- ^ The RHS of the system \(\dot{y} = f(t,y)\) | ||
197 | -> V.Vector Double -- ^ Initial conditions | ||
198 | -> V.Vector Double -- ^ Desired solution times | ||
199 | -> Either Int (Matrix Double, SundialsDiagnostics) -- ^ Error code or solution | ||
200 | odeSolveVWith' opts method control initStepSize f y0 tt = | ||
201 | case solveOdeC (fromIntegral $ maxFail opts) | ||
202 | (fromIntegral $ maxNumSteps opts) (coerce $ minStep opts) | ||
203 | (fromIntegral $ getMethod method) (coerce initStepSize) jacH (scise control) | ||
204 | (coerce f) (coerce y0) (coerce tt) of | ||
205 | Left c -> Left $ fromIntegral c | ||
206 | Right (v, d) -> Right (reshape l (coerce v), d) | ||
207 | where | ||
208 | l = size y0 | ||
209 | scise (X aTol rTol) = coerce (V.replicate l aTol, rTol) | ||
210 | scise (X' aTol rTol) = coerce (V.replicate l aTol, rTol) | ||
211 | scise (XX' aTol rTol yScale _yDotScale) = coerce (V.replicate l aTol, yScale * rTol) | ||
212 | -- FIXME; Should we check that the length of ss is correct? | ||
213 | scise (ScXX' aTol rTol yScale _yDotScale ss) = coerce (V.map (* aTol) ss, yScale * rTol) | ||
214 | jacH = fmap (\g t v -> matrixToSunMatrix $ g (coerce t) (coerce v)) $ | ||
215 | getJacobian method | ||
216 | matrixToSunMatrix m = T.SunMatrix { T.rows = nr, T.cols = nc, T.vals = vs } | ||
217 | where | ||
218 | nr = fromIntegral $ rows m | ||
219 | nc = fromIntegral $ cols m | ||
220 | -- FIXME: efficiency | ||
221 | vs = V.fromList $ map coerce $ concat $ toLists m | ||
222 | |||
223 | solveOdeC :: | ||
224 | CInt -> | ||
225 | CLong -> | ||
226 | CDouble -> | ||
227 | CInt -> | ||
228 | Maybe CDouble -> | ||
229 | (Maybe (CDouble -> V.Vector CDouble -> T.SunMatrix)) -> | ||
230 | (V.Vector CDouble, CDouble) -> | ||
231 | (CDouble -> V.Vector CDouble -> V.Vector CDouble) -- ^ The RHS of the system \(\dot{y} = f(t,y)\) | ||
232 | -> V.Vector CDouble -- ^ Initial conditions | ||
233 | -> V.Vector CDouble -- ^ Desired solution times | ||
234 | -> Either CInt ((V.Vector CDouble), SundialsDiagnostics) -- ^ Error code or solution | ||
235 | solveOdeC maxErrTestFails maxNumSteps_ minStep_ method initStepSize | ||
236 | jacH (aTols, rTol) fun f0 ts = | ||
237 | unsafePerformIO $ do | ||
238 | |||
239 | let isInitStepSize :: CInt | ||
240 | isInitStepSize = fromIntegral $ fromEnum $ isJust initStepSize | ||
241 | ss :: CDouble | ||
242 | ss = case initStepSize of | ||
243 | -- It would be better to put an error message here but | ||
244 | -- inline-c seems to evaluate this even if it is never | ||
245 | -- used :( | ||
246 | Nothing -> 0.0 | ||
247 | Just x -> x | ||
248 | |||
249 | let dim = V.length f0 | ||
250 | nEq :: CLong | ||
251 | nEq = fromIntegral dim | ||
252 | nTs :: CInt | ||
253 | nTs = fromIntegral $ V.length ts | ||
254 | quasiMatrixRes <- createVector ((fromIntegral dim) * (fromIntegral nTs)) | ||
255 | qMatMut <- V.thaw quasiMatrixRes | ||
256 | diagnostics :: V.Vector CLong <- createVector 10 -- FIXME | ||
257 | diagMut <- V.thaw diagnostics | ||
258 | -- We need the types that sundials expects. These are tied together | ||
259 | -- in 'CLangToHaskellTypes'. FIXME: The Haskell type is currently empty! | ||
260 | let funIO :: CDouble -> Ptr T.SunVector -> Ptr T.SunVector -> Ptr () -> IO CInt | ||
261 | funIO x y f _ptr = do | ||
262 | -- Convert the pointer we get from C (y) to a vector, and then | ||
263 | -- apply the user-supplied function. | ||
264 | fImm <- fun x <$> getDataFromContents dim y | ||
265 | -- Fill in the provided pointer with the resulting vector. | ||
266 | putDataInContents fImm dim f | ||
267 | -- FIXME: I don't understand what this comment means | ||
268 | -- Unsafe since the function will be called many times. | ||
269 | [CU.exp| int{ 0 } |] | ||
270 | let isJac :: CInt | ||
271 | isJac = fromIntegral $ fromEnum $ isJust jacH | ||
272 | jacIO :: CDouble -> Ptr T.SunVector -> Ptr T.SunVector -> Ptr T.SunMatrix -> | ||
273 | Ptr () -> Ptr T.SunVector -> Ptr T.SunVector -> Ptr T.SunVector -> | ||
274 | IO CInt | ||
275 | jacIO t y _fy jacS _ptr _tmp1 _tmp2 _tmp3 = do | ||
276 | case jacH of | ||
277 | Nothing -> error "Numeric.Sundials.CVode.ODE: Jacobian not defined" | ||
278 | Just jacI -> do j <- jacI t <$> getDataFromContents dim y | ||
279 | poke jacS j | ||
280 | -- FIXME: I don't understand what this comment means | ||
281 | -- Unsafe since the function will be called many times. | ||
282 | [CU.exp| int{ 0 } |] | ||
283 | |||
284 | res <- [C.block| int { | ||
285 | /* general problem variables */ | ||
286 | |||
287 | int flag; /* reusable error-checking flag */ | ||
288 | int i, j; /* reusable loop indices */ | ||
289 | N_Vector y = NULL; /* empty vector for storing solution */ | ||
290 | N_Vector tv = NULL; /* empty vector for storing absolute tolerances */ | ||
291 | |||
292 | SUNMatrix A = NULL; /* empty matrix for linear solver */ | ||
293 | SUNLinearSolver LS = NULL; /* empty linear solver object */ | ||
294 | void *cvode_mem = NULL; /* empty CVODE memory structure */ | ||
295 | realtype t; | ||
296 | long int nst, nfe, nsetups, nje, nfeLS, nni, ncfn, netf, nge; | ||
297 | |||
298 | /* general problem parameters */ | ||
299 | |||
300 | realtype T0 = RCONST(($vec-ptr:(double *ts))[0]); /* initial time */ | ||
301 | sunindextype NEQ = $(sunindextype nEq); /* number of dependent vars. */ | ||
302 | |||
303 | /* Initialize data structures */ | ||
304 | |||
305 | y = N_VNew_Serial(NEQ); /* Create serial vector for solution */ | ||
306 | if (check_flag((void *)y, "N_VNew_Serial", 0)) return 1; | ||
307 | /* Specify initial condition */ | ||
308 | for (i = 0; i < NEQ; i++) { | ||
309 | NV_Ith_S(y,i) = ($vec-ptr:(double *f0))[i]; | ||
310 | }; | ||
311 | |||
312 | cvode_mem = CVodeCreate($(int method), CV_NEWTON); | ||
313 | if (check_flag((void *)cvode_mem, "CVodeCreate", 0)) return(1); | ||
314 | |||
315 | /* Call CVodeInit to initialize the integrator memory and specify the | ||
316 | * user's right hand side function in y'=f(t,y), the inital time T0, and | ||
317 | * the initial dependent variable vector y. */ | ||
318 | flag = CVodeInit(cvode_mem, $fun:(int (* funIO) (double t, SunVector y[], SunVector dydt[], void * params)), T0, y); | ||
319 | if (check_flag(&flag, "CVodeInit", 1)) return(1); | ||
320 | |||
321 | tv = N_VNew_Serial(NEQ); /* Create serial vector for absolute tolerances */ | ||
322 | if (check_flag((void *)tv, "N_VNew_Serial", 0)) return 1; | ||
323 | /* Specify tolerances */ | ||
324 | for (i = 0; i < NEQ; i++) { | ||
325 | NV_Ith_S(tv,i) = ($vec-ptr:(double *aTols))[i]; | ||
326 | }; | ||
327 | |||
328 | flag = CVodeSetMinStep(cvode_mem, $(double minStep_)); | ||
329 | if (check_flag(&flag, "CVodeSetMinStep", 1)) return 1; | ||
330 | flag = CVodeSetMaxNumSteps(cvode_mem, $(long int maxNumSteps_)); | ||
331 | if (check_flag(&flag, "CVodeSetMaxNumSteps", 1)) return 1; | ||
332 | flag = CVodeSetMaxErrTestFails(cvode_mem, $(int maxErrTestFails)); | ||
333 | if (check_flag(&flag, "CVodeSetMaxErrTestFails", 1)) return 1; | ||
334 | |||
335 | /* Call CVodeSVtolerances to specify the scalar relative tolerance | ||
336 | * and vector absolute tolerances */ | ||
337 | flag = CVodeSVtolerances(cvode_mem, $(double rTol), tv); | ||
338 | if (check_flag(&flag, "CVodeSVtolerances", 1)) return(1); | ||
339 | |||
340 | /* Initialize dense matrix data structure and solver */ | ||
341 | A = SUNDenseMatrix(NEQ, NEQ); | ||
342 | if (check_flag((void *)A, "SUNDenseMatrix", 0)) return 1; | ||
343 | LS = SUNDenseLinearSolver(y, A); | ||
344 | if (check_flag((void *)LS, "SUNDenseLinearSolver", 0)) return 1; | ||
345 | |||
346 | /* Attach matrix and linear solver */ | ||
347 | flag = CVDlsSetLinearSolver(cvode_mem, LS, A); | ||
348 | if (check_flag(&flag, "CVDlsSetLinearSolver", 1)) return 1; | ||
349 | |||
350 | /* Set the initial step size if there is one */ | ||
351 | if ($(int isInitStepSize)) { | ||
352 | /* FIXME: We could check if the initial step size is 0 */ | ||
353 | /* or even NaN and then throw an error */ | ||
354 | flag = CVodeSetInitStep(cvode_mem, $(double ss)); | ||
355 | if (check_flag(&flag, "CVodeSetInitStep", 1)) return 1; | ||
356 | } | ||
357 | |||
358 | /* Set the Jacobian if there is one */ | ||
359 | if ($(int isJac)) { | ||
360 | flag = CVDlsSetJacFn(cvode_mem, $fun:(int (* jacIO) (double t, SunVector y[], SunVector fy[], SunMatrix Jac[], void * params, SunVector tmp1[], SunVector tmp2[], SunVector tmp3[]))); | ||
361 | if (check_flag(&flag, "CVDlsSetJacFn", 1)) return 1; | ||
362 | } | ||
363 | |||
364 | /* Store initial conditions */ | ||
365 | for (j = 0; j < NEQ; j++) { | ||
366 | ($vec-ptr:(double *qMatMut))[0 * $(int nTs) + j] = NV_Ith_S(y,j); | ||
367 | } | ||
368 | |||
369 | /* Main time-stepping loop: calls CVode to perform the integration */ | ||
370 | /* Stops when the final time has been reached */ | ||
371 | for (i = 1; i < $(int nTs); i++) { | ||
372 | |||
373 | flag = CVode(cvode_mem, ($vec-ptr:(double *ts))[i], y, &t, CV_NORMAL); /* call integrator */ | ||
374 | if (check_flag(&flag, "CVode", 1)) break; | ||
375 | |||
376 | /* Store the results for Haskell */ | ||
377 | for (j = 0; j < NEQ; j++) { | ||
378 | ($vec-ptr:(double *qMatMut))[i * NEQ + j] = NV_Ith_S(y,j); | ||
379 | } | ||
380 | |||
381 | /* unsuccessful solve: break */ | ||
382 | if (flag < 0) { | ||
383 | fprintf(stderr,"Solver failure, stopping integration\n"); | ||
384 | break; | ||
385 | } | ||
386 | } | ||
387 | |||
388 | /* Get some final statistics on how the solve progressed */ | ||
389 | |||
390 | flag = CVodeGetNumSteps(cvode_mem, &nst); | ||
391 | check_flag(&flag, "CVodeGetNumSteps", 1); | ||
392 | ($vec-ptr:(long int *diagMut))[0] = nst; | ||
393 | |||
394 | /* FIXME */ | ||
395 | ($vec-ptr:(long int *diagMut))[1] = 0; | ||
396 | |||
397 | flag = CVodeGetNumRhsEvals(cvode_mem, &nfe); | ||
398 | check_flag(&flag, "CVodeGetNumRhsEvals", 1); | ||
399 | ($vec-ptr:(long int *diagMut))[2] = nfe; | ||
400 | /* FIXME */ | ||
401 | ($vec-ptr:(long int *diagMut))[3] = 0; | ||
402 | |||
403 | flag = CVodeGetNumLinSolvSetups(cvode_mem, &nsetups); | ||
404 | check_flag(&flag, "CVodeGetNumLinSolvSetups", 1); | ||
405 | ($vec-ptr:(long int *diagMut))[4] = nsetups; | ||
406 | |||
407 | flag = CVodeGetNumErrTestFails(cvode_mem, &netf); | ||
408 | check_flag(&flag, "CVodeGetNumErrTestFails", 1); | ||
409 | ($vec-ptr:(long int *diagMut))[5] = netf; | ||
410 | |||
411 | flag = CVodeGetNumNonlinSolvIters(cvode_mem, &nni); | ||
412 | check_flag(&flag, "CVodeGetNumNonlinSolvIters", 1); | ||
413 | ($vec-ptr:(long int *diagMut))[6] = nni; | ||
414 | |||
415 | flag = CVodeGetNumNonlinSolvConvFails(cvode_mem, &ncfn); | ||
416 | check_flag(&flag, "CVodeGetNumNonlinSolvConvFails", 1); | ||
417 | ($vec-ptr:(long int *diagMut))[7] = ncfn; | ||
418 | |||
419 | flag = CVDlsGetNumJacEvals(cvode_mem, &nje); | ||
420 | check_flag(&flag, "CVDlsGetNumJacEvals", 1); | ||
421 | ($vec-ptr:(long int *diagMut))[8] = ncfn; | ||
422 | |||
423 | flag = CVDlsGetNumRhsEvals(cvode_mem, &nfeLS); | ||
424 | check_flag(&flag, "CVDlsGetNumRhsEvals", 1); | ||
425 | ($vec-ptr:(long int *diagMut))[9] = ncfn; | ||
426 | |||
427 | /* Clean up and return */ | ||
428 | |||
429 | N_VDestroy(y); /* Free y vector */ | ||
430 | N_VDestroy(tv); /* Free tv vector */ | ||
431 | CVodeFree(&cvode_mem); /* Free integrator memory */ | ||
432 | SUNLinSolFree(LS); /* Free linear solver */ | ||
433 | SUNMatDestroy(A); /* Free A matrix */ | ||
434 | |||
435 | return flag; | ||
436 | } |] | ||
437 | if res == 0 | ||
438 | then do | ||
439 | preD <- V.freeze diagMut | ||
440 | let d = SundialsDiagnostics (fromIntegral $ preD V.!0) | ||
441 | (fromIntegral $ preD V.!1) | ||
442 | (fromIntegral $ preD V.!2) | ||
443 | (fromIntegral $ preD V.!3) | ||
444 | (fromIntegral $ preD V.!4) | ||
445 | (fromIntegral $ preD V.!5) | ||
446 | (fromIntegral $ preD V.!6) | ||
447 | (fromIntegral $ preD V.!7) | ||
448 | (fromIntegral $ preD V.!8) | ||
449 | (fromIntegral $ preD V.!9) | ||
450 | m <- V.freeze qMatMut | ||
451 | return $ Right (m, d) | ||
452 | else do | ||
453 | return $ Left res | ||
454 | |||
455 | -- | Adaptive step-size control | ||
456 | -- functions. | ||
457 | -- | ||
458 | -- [GSL](https://www.gnu.org/software/gsl/doc/html/ode-initval.html#adaptive-step-size-control) | ||
459 | -- allows the user to control the step size adjustment using | ||
460 | -- \(D_i = \epsilon^{abs}s_i + \epsilon^{rel}(a_{y} |y_i| + a_{dy/dt} h |\dot{y}_i|)\) where | ||
461 | -- \(\epsilon^{abs}\) is the required absolute error, \(\epsilon^{rel}\) | ||
462 | -- is the required relative error, \(s_i\) is a vector of scaling | ||
463 | -- factors, \(a_{y}\) is a scaling factor for the solution \(y\) and | ||
464 | -- \(a_{dydt}\) is a scaling factor for the derivative of the solution \(dy/dt\). | ||
465 | -- | ||
466 | -- [ARKode](https://computation.llnl.gov/projects/sundials/arkode) | ||
467 | -- allows the user to control the step size adjustment using | ||
468 | -- \(\eta^{rel}|y_i| + \eta^{abs}_i\). For compatibility with | ||
469 | -- [hmatrix-gsl](https://hackage.haskell.org/package/hmatrix-gsl), | ||
470 | -- tolerances for \(y\) and \(\dot{y}\) can be specified but the latter have no | ||
471 | -- effect. | ||
472 | data 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 | ||
473 | | 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 | ||
474 | | XX' Double Double Double Double -- ^ include both via relative tolerance | ||
475 | -- scaling factors \(a_y\), \(a_{{dy}/{dt}}\); in ARKode terms, the latter is ignored and \(\eta^{rel} = a_{y}\epsilon^{rel}\) | ||
476 | | 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}\) | ||