1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
|
{-# LANGUAGE QuasiQuotes #-}
{-# LANGUAGE TemplateHaskell #-}
{-# LANGUAGE MultiWayIf #-}
{-# LANGUAGE OverloadedStrings #-}
{-# LANGUAGE ScopedTypeVariables #-}
{-# LANGUAGE DeriveGeneric #-}
{-# LANGUAGE TypeOperators #-}
{-# LANGUAGE KindSignatures #-}
{-# LANGUAGE TypeSynonymInstances #-}
{-# LANGUAGE FlexibleInstances #-}
{-# LANGUAGE FlexibleContexts #-}
-----------------------------------------------------------------------------
-- |
-- Module : Numeric.Sundials.ARKode.ODE
-- Copyright : Dominic Steinitz 2018,
-- Novadiscovery 2018
-- License : BSD
-- Maintainer : Dominic Steinitz
-- Stability : provisional
--
-- Solution of ordinary differential equation (ODE) initial value problems.
-- See <https://computation.llnl.gov/projects/sundials/sundials-software> for more detail.
--
-- A simple example:
--
-- <<diagrams/brusselator.png#diagram=brusselator&height=400&width=500>>
--
-- @
-- import Numeric.Sundials.ARKode.ODE
-- import Numeric.LinearAlgebra
--
-- import Plots as P
-- import qualified Diagrams.Prelude as D
-- import Diagrams.Backend.Rasterific
--
-- brusselator :: Double -> [Double] -> [Double]
-- brusselator _t x = [ a - (w + 1) * u + v * u * u
-- , w * u - v * u * u
-- , (b - w) / eps - w * u
-- ]
-- where
-- a = 1.0
-- b = 3.5
-- eps = 5.0e-6
-- u = x !! 0
-- v = x !! 1
-- w = x !! 2
--
-- lSaxis :: [[Double]] -> P.Axis B D.V2 Double
-- lSaxis xs = P.r2Axis &~ do
-- let ts = xs!!0
-- us = xs!!1
-- vs = xs!!2
-- ws = xs!!3
-- P.linePlot' $ zip ts us
-- P.linePlot' $ zip ts vs
-- P.linePlot' $ zip ts ws
--
-- main = do
-- let res1 = odeSolve brusselator [1.2, 3.1, 3.0] (fromList [0.0, 0.1 .. 10.0])
-- renderRasterific "diagrams/brusselator.png"
-- (D.dims2D 500.0 500.0)
-- (renderAxis $ lSaxis $ [0.0, 0.1 .. 10.0]:(toLists $ tr res1))
-- @
--
-- With Sundials ARKode, it is possible to retrieve the Butcher tableau for the solver.
--
-- @
-- import Numeric.Sundials.ARKode.ODE
-- import Numeric.LinearAlgebra
--
-- import Data.List (intercalate)
--
-- import Text.PrettyPrint.HughesPJClass
--
--
-- butcherTableauTex :: ButcherTable -> String
-- butcherTableauTex (ButcherTable m c b b2) =
-- render $
-- vcat [ text ("\n\\begin{array}{c|" ++ (concat $ replicate n "c") ++ "}")
-- , us
-- , text "\\hline"
-- , text bs <+> text "\\\\"
-- , text b2s <+> text "\\\\"
-- , text "\\end{array}"
-- ]
-- where
-- n = rows m
-- rs = toLists m
-- ss = map (\r -> intercalate " & " $ map show r) rs
-- ts = zipWith (\i r -> show i ++ " & " ++ r) (toList c) ss
-- us = vcat $ map (\r -> text r <+> text "\\\\") ts
-- bs = " & " ++ (intercalate " & " $ map show $ toList b)
-- b2s = " & " ++ (intercalate " & " $ map show $ toList b2)
--
-- main :: IO ()
-- main = do
--
-- let res = butcherTable (SDIRK_2_1_2 undefined)
-- putStrLn $ show res
-- putStrLn $ butcherTableauTex res
--
-- let resA = butcherTable (KVAERNO_4_2_3 undefined)
-- putStrLn $ show resA
-- putStrLn $ butcherTableauTex resA
--
-- let resB = butcherTable (SDIRK_5_3_4 undefined)
-- putStrLn $ show resB
-- putStrLn $ butcherTableauTex resB
-- @
--
-- Using the code above from the examples gives
--
-- KVAERNO_4_2_3
--
-- \[
-- \begin{array}{c|cccc}
-- 0.0 & 0.0 & 0.0 & 0.0 & 0.0 \\
-- 0.871733043 & 0.4358665215 & 0.4358665215 & 0.0 & 0.0 \\
-- 1.0 & 0.490563388419108 & 7.3570090080892e-2 & 0.4358665215 & 0.0 \\
-- 1.0 & 0.308809969973036 & 1.490563388254106 & -1.235239879727145 & 0.4358665215 \\
-- \hline
-- & 0.308809969973036 & 1.490563388254106 & -1.235239879727145 & 0.4358665215 \\
-- & 0.490563388419108 & 7.3570090080892e-2 & 0.4358665215 & 0.0 \\
-- \end{array}
-- \]
--
-- SDIRK_2_1_2
--
-- \[
-- \begin{array}{c|cc}
-- 1.0 & 1.0 & 0.0 \\
-- 0.0 & -1.0 & 1.0 \\
-- \hline
-- & 0.5 & 0.5 \\
-- & 1.0 & 0.0 \\
-- \end{array}
-- \]
--
-- SDIRK_5_3_4
--
-- \[
-- \begin{array}{c|ccccc}
-- 0.25 & 0.25 & 0.0 & 0.0 & 0.0 & 0.0 \\
-- 0.75 & 0.5 & 0.25 & 0.0 & 0.0 & 0.0 \\
-- 0.55 & 0.34 & -4.0e-2 & 0.25 & 0.0 & 0.0 \\
-- 0.5 & 0.2727941176470588 & -5.036764705882353e-2 & 2.7573529411764705e-2 & 0.25 & 0.0 \\
-- 1.0 & 1.0416666666666667 & -1.0208333333333333 & 7.8125 & -7.083333333333333 & 0.25 \\
-- \hline
-- & 1.0416666666666667 & -1.0208333333333333 & 7.8125 & -7.083333333333333 & 0.25 \\
-- & 1.2291666666666667 & -0.17708333333333334 & 7.03125 & -7.083333333333333 & 0.0 \\
-- \end{array}
-- \]
-----------------------------------------------------------------------------
module Numeric.Sundials.ARKode.ODE ( odeSolve
, odeSolveV
, odeSolveVWith
, odeSolveVWith'
, ButcherTable(..)
, butcherTable
, ODEMethod(..)
, StepControl(..)
) where
import qualified Language.C.Inline as C
import qualified Language.C.Inline.Unsafe as CU
import Data.Monoid ((<>))
import Data.Maybe (isJust)
import Foreign.C.Types (CDouble, CInt, CLong)
import Foreign.Ptr (Ptr)
import Foreign.Storable (poke)
import qualified Data.Vector.Storable as V
import Data.Coerce (coerce)
import System.IO.Unsafe (unsafePerformIO)
import GHC.Generics (C1, Constructor, (:+:)(..), D1, Rep, Generic, M1(..),
from, conName)
import Numeric.LinearAlgebra.Devel (createVector)
import Numeric.LinearAlgebra.HMatrix (Vector, Matrix, toList, rows,
cols, toLists, size, reshape,
subVector, subMatrix, (><))
import Numeric.Sundials.ODEOpts (ODEOpts(..), Jacobian, SundialsDiagnostics(..))
import qualified Numeric.Sundials.Arkode as T
import Numeric.Sundials.Arkode (getDataFromContents, putDataInContents, arkSMax,
sDIRK_2_1_2,
bILLINGTON_3_3_2,
tRBDF2_3_3_2,
kVAERNO_4_2_3,
aRK324L2SA_DIRK_4_2_3,
cASH_5_2_4,
cASH_5_3_4,
sDIRK_5_3_4,
kVAERNO_5_3_4,
aRK436L2SA_DIRK_6_3_4,
kVAERNO_7_4_5,
aRK548L2SA_DIRK_8_4_5,
hEUN_EULER_2_1_2,
bOGACKI_SHAMPINE_4_2_3,
aRK324L2SA_ERK_4_2_3,
zONNEVELD_5_3_4,
aRK436L2SA_ERK_6_3_4,
sAYFY_ABURUB_6_3_4,
cASH_KARP_6_4_5,
fEHLBERG_6_4_5,
dORMAND_PRINCE_7_4_5,
aRK548L2SA_ERK_8_4_5,
vERNER_8_5_6,
fEHLBERG_13_7_8)
C.context (C.baseCtx <> C.vecCtx <> C.funCtx <> T.sunCtx)
C.include "<stdlib.h>"
C.include "<stdio.h>"
C.include "<math.h>"
C.include "<arkode/arkode.h>" -- prototypes for ARKODE fcts., consts.
C.include "<nvector/nvector_serial.h>" -- serial N_Vector types, fcts., macros
C.include "<sunmatrix/sunmatrix_dense.h>" -- access to dense SUNMatrix
C.include "<sunlinsol/sunlinsol_dense.h>" -- access to dense SUNLinearSolver
C.include "<arkode/arkode_direct.h>" -- access to ARKDls interface
C.include "<sundials/sundials_types.h>" -- definition of type realtype
C.include "<sundials/sundials_math.h>"
C.include "../../../helpers.h"
C.include "Numeric/Sundials/Arkode_hsc.h"
-- | Stepping functions
data ODEMethod = SDIRK_2_1_2 Jacobian
| SDIRK_2_1_2'
| BILLINGTON_3_3_2 Jacobian
| BILLINGTON_3_3_2'
| TRBDF2_3_3_2 Jacobian
| TRBDF2_3_3_2'
| KVAERNO_4_2_3 Jacobian
| KVAERNO_4_2_3'
| ARK324L2SA_DIRK_4_2_3 Jacobian
| ARK324L2SA_DIRK_4_2_3'
| CASH_5_2_4 Jacobian
| CASH_5_2_4'
| CASH_5_3_4 Jacobian
| CASH_5_3_4'
| SDIRK_5_3_4 Jacobian
| SDIRK_5_3_4'
| KVAERNO_5_3_4 Jacobian
| KVAERNO_5_3_4'
| ARK436L2SA_DIRK_6_3_4 Jacobian
| ARK436L2SA_DIRK_6_3_4'
| KVAERNO_7_4_5 Jacobian
| KVAERNO_7_4_5'
| ARK548L2SA_DIRK_8_4_5 Jacobian
| ARK548L2SA_DIRK_8_4_5'
| HEUN_EULER_2_1_2 Jacobian
| HEUN_EULER_2_1_2'
| BOGACKI_SHAMPINE_4_2_3 Jacobian
| BOGACKI_SHAMPINE_4_2_3'
| ARK324L2SA_ERK_4_2_3 Jacobian
| ARK324L2SA_ERK_4_2_3'
| ZONNEVELD_5_3_4 Jacobian
| ZONNEVELD_5_3_4'
| ARK436L2SA_ERK_6_3_4 Jacobian
| ARK436L2SA_ERK_6_3_4'
| SAYFY_ABURUB_6_3_4 Jacobian
| SAYFY_ABURUB_6_3_4'
| CASH_KARP_6_4_5 Jacobian
| CASH_KARP_6_4_5'
| FEHLBERG_6_4_5 Jacobian
| FEHLBERG_6_4_5'
| DORMAND_PRINCE_7_4_5 Jacobian
| DORMAND_PRINCE_7_4_5'
| ARK548L2SA_ERK_8_4_5 Jacobian
| ARK548L2SA_ERK_8_4_5'
| VERNER_8_5_6 Jacobian
| VERNER_8_5_6'
| FEHLBERG_13_7_8 Jacobian
| FEHLBERG_13_7_8'
deriving Generic
constrName :: (HasConstructor (Rep a), Generic a)=> a -> String
constrName = genericConstrName . from
class HasConstructor (f :: * -> *) where
genericConstrName :: f x -> String
instance HasConstructor f => HasConstructor (D1 c f) where
genericConstrName (M1 x) = genericConstrName x
instance (HasConstructor x, HasConstructor y) => HasConstructor (x :+: y) where
genericConstrName (L1 l) = genericConstrName l
genericConstrName (R1 r) = genericConstrName r
instance Constructor c => HasConstructor (C1 c f) where
genericConstrName x = conName x
instance Show ODEMethod where
show x = constrName x
-- FIXME: We can probably do better here with generics
getMethod :: ODEMethod -> Int
getMethod (SDIRK_2_1_2 _) = sDIRK_2_1_2
getMethod (SDIRK_2_1_2') = sDIRK_2_1_2
getMethod (BILLINGTON_3_3_2 _) = bILLINGTON_3_3_2
getMethod (BILLINGTON_3_3_2') = bILLINGTON_3_3_2
getMethod (TRBDF2_3_3_2 _) = tRBDF2_3_3_2
getMethod (TRBDF2_3_3_2') = tRBDF2_3_3_2
getMethod (KVAERNO_4_2_3 _) = kVAERNO_4_2_3
getMethod (KVAERNO_4_2_3') = kVAERNO_4_2_3
getMethod (ARK324L2SA_DIRK_4_2_3 _) = aRK324L2SA_DIRK_4_2_3
getMethod (ARK324L2SA_DIRK_4_2_3') = aRK324L2SA_DIRK_4_2_3
getMethod (CASH_5_2_4 _) = cASH_5_2_4
getMethod (CASH_5_2_4') = cASH_5_2_4
getMethod (CASH_5_3_4 _) = cASH_5_3_4
getMethod (CASH_5_3_4') = cASH_5_3_4
getMethod (SDIRK_5_3_4 _) = sDIRK_5_3_4
getMethod (SDIRK_5_3_4') = sDIRK_5_3_4
getMethod (KVAERNO_5_3_4 _) = kVAERNO_5_3_4
getMethod (KVAERNO_5_3_4') = kVAERNO_5_3_4
getMethod (ARK436L2SA_DIRK_6_3_4 _) = aRK436L2SA_DIRK_6_3_4
getMethod (ARK436L2SA_DIRK_6_3_4') = aRK436L2SA_DIRK_6_3_4
getMethod (KVAERNO_7_4_5 _) = kVAERNO_7_4_5
getMethod (KVAERNO_7_4_5') = kVAERNO_7_4_5
getMethod (ARK548L2SA_DIRK_8_4_5 _) = aRK548L2SA_DIRK_8_4_5
getMethod (ARK548L2SA_DIRK_8_4_5') = aRK548L2SA_DIRK_8_4_5
getMethod (HEUN_EULER_2_1_2 _) = hEUN_EULER_2_1_2
getMethod (HEUN_EULER_2_1_2') = hEUN_EULER_2_1_2
getMethod (BOGACKI_SHAMPINE_4_2_3 _) = bOGACKI_SHAMPINE_4_2_3
getMethod (BOGACKI_SHAMPINE_4_2_3') = bOGACKI_SHAMPINE_4_2_3
getMethod (ARK324L2SA_ERK_4_2_3 _) = aRK324L2SA_ERK_4_2_3
getMethod (ARK324L2SA_ERK_4_2_3') = aRK324L2SA_ERK_4_2_3
getMethod (ZONNEVELD_5_3_4 _) = zONNEVELD_5_3_4
getMethod (ZONNEVELD_5_3_4') = zONNEVELD_5_3_4
getMethod (ARK436L2SA_ERK_6_3_4 _) = aRK436L2SA_ERK_6_3_4
getMethod (ARK436L2SA_ERK_6_3_4') = aRK436L2SA_ERK_6_3_4
getMethod (SAYFY_ABURUB_6_3_4 _) = sAYFY_ABURUB_6_3_4
getMethod (SAYFY_ABURUB_6_3_4') = sAYFY_ABURUB_6_3_4
getMethod (CASH_KARP_6_4_5 _) = cASH_KARP_6_4_5
getMethod (CASH_KARP_6_4_5') = cASH_KARP_6_4_5
getMethod (FEHLBERG_6_4_5 _) = fEHLBERG_6_4_5
getMethod (FEHLBERG_6_4_5' ) = fEHLBERG_6_4_5
getMethod (DORMAND_PRINCE_7_4_5 _) = dORMAND_PRINCE_7_4_5
getMethod (DORMAND_PRINCE_7_4_5') = dORMAND_PRINCE_7_4_5
getMethod (ARK548L2SA_ERK_8_4_5 _) = aRK548L2SA_ERK_8_4_5
getMethod (ARK548L2SA_ERK_8_4_5') = aRK548L2SA_ERK_8_4_5
getMethod (VERNER_8_5_6 _) = vERNER_8_5_6
getMethod (VERNER_8_5_6') = vERNER_8_5_6
getMethod (FEHLBERG_13_7_8 _) = fEHLBERG_13_7_8
getMethod (FEHLBERG_13_7_8') = fEHLBERG_13_7_8
getJacobian :: ODEMethod -> Maybe Jacobian
getJacobian (SDIRK_2_1_2 j) = Just j
getJacobian (BILLINGTON_3_3_2 j) = Just j
getJacobian (TRBDF2_3_3_2 j) = Just j
getJacobian (KVAERNO_4_2_3 j) = Just j
getJacobian (ARK324L2SA_DIRK_4_2_3 j) = Just j
getJacobian (CASH_5_2_4 j) = Just j
getJacobian (CASH_5_3_4 j) = Just j
getJacobian (SDIRK_5_3_4 j) = Just j
getJacobian (KVAERNO_5_3_4 j) = Just j
getJacobian (ARK436L2SA_DIRK_6_3_4 j) = Just j
getJacobian (KVAERNO_7_4_5 j) = Just j
getJacobian (ARK548L2SA_DIRK_8_4_5 j) = Just j
getJacobian (HEUN_EULER_2_1_2 j) = Just j
getJacobian (BOGACKI_SHAMPINE_4_2_3 j) = Just j
getJacobian (ARK324L2SA_ERK_4_2_3 j) = Just j
getJacobian (ZONNEVELD_5_3_4 j) = Just j
getJacobian (ARK436L2SA_ERK_6_3_4 j) = Just j
getJacobian (SAYFY_ABURUB_6_3_4 j) = Just j
getJacobian (CASH_KARP_6_4_5 j) = Just j
getJacobian (FEHLBERG_6_4_5 j) = Just j
getJacobian (DORMAND_PRINCE_7_4_5 j) = Just j
getJacobian (ARK548L2SA_ERK_8_4_5 j) = Just j
getJacobian (VERNER_8_5_6 j) = Just j
getJacobian (FEHLBERG_13_7_8 j) = Just j
getJacobian _ = Nothing
-- | A version of 'odeSolveVWith' with reasonable default step control.
odeSolveV
:: ODEMethod
-> Maybe Double -- ^ initial step size - by default, ARKode
-- estimates the initial step size to be the
-- solution \(h\) of the equation
-- \(\|\frac{h^2\ddot{y}}{2}\| = 1\), where
-- \(\ddot{y}\) is an estimated value of the
-- second derivative of the solution at \(t_0\)
-> Double -- ^ absolute tolerance for the state vector
-> Double -- ^ relative tolerance for the state vector
-> (Double -> Vector Double -> Vector Double) -- ^ The RHS of the system \(\dot{y} = f(t,y)\)
-> Vector Double -- ^ initial conditions
-> Vector Double -- ^ desired solution times
-> Matrix Double -- ^ solution
odeSolveV meth hi epsAbs epsRel f y0 ts =
odeSolveVWith meth (X epsAbs epsRel) hi g y0 ts
where
g t x0 = coerce $ f t x0
-- | A version of 'odeSolveV' with reasonable default parameters and
-- system of equations defined using lists. FIXME: we should say
-- something about the fact we could use the Jacobian but don't for
-- compatibility with hmatrix-gsl.
odeSolve :: (Double -> [Double] -> [Double]) -- ^ The RHS of the system \(\dot{y} = f(t,y)\)
-> [Double] -- ^ initial conditions
-> Vector Double -- ^ desired solution times
-> Matrix Double -- ^ solution
odeSolve f y0 ts =
-- FIXME: These tolerances are different from the ones in GSL
odeSolveVWith SDIRK_5_3_4' (XX' 1.0e-6 1.0e-10 1 1) Nothing g (V.fromList y0) (V.fromList $ toList ts)
where
g t x0 = V.fromList $ f t (V.toList x0)
odeSolveVWith ::
ODEMethod
-> StepControl
-> Maybe Double -- ^ initial step size - by default, ARKode
-- estimates the initial step size to be the
-- solution \(h\) of the equation
-- \(\|\frac{h^2\ddot{y}}{2}\| = 1\), where
-- \(\ddot{y}\) is an estimated value of the second
-- derivative of the solution at \(t_0\)
-> (Double -> V.Vector Double -> V.Vector Double) -- ^ The RHS of the system \(\dot{y} = f(t,y)\)
-> V.Vector Double -- ^ Initial conditions
-> V.Vector Double -- ^ Desired solution times
-> Matrix Double -- ^ Error code or solution
odeSolveVWith method control initStepSize f y0 tt =
case odeSolveVWith' opts method control initStepSize f y0 tt of
Left (c, _v) -> error $ show c -- FIXME
Right (v, _d) -> v
where
opts = ODEOpts { maxNumSteps = 10000
, minStep = 1.0e-12
, relTol = error "relTol"
, absTols = error "absTol"
, initStep = error "initStep"
, maxFail = 10
}
odeSolveVWith' ::
ODEOpts
-> ODEMethod
-> StepControl
-> Maybe Double -- ^ initial step size - by default, ARKode
-- estimates the initial step size to be the
-- solution \(h\) of the equation
-- \(\|\frac{h^2\ddot{y}}{2}\| = 1\), where
-- \(\ddot{y}\) is an estimated value of the second
-- derivative of the solution at \(t_0\)
-> (Double -> V.Vector Double -> V.Vector Double) -- ^ The RHS of the system \(\dot{y} = f(t,y)\)
-> V.Vector Double -- ^ Initial conditions
-> V.Vector Double -- ^ Desired solution times
-> Either (Matrix Double, Int) (Matrix Double, SundialsDiagnostics) -- ^ Error code or solution
odeSolveVWith' opts method control initStepSize f y0 tt =
case solveOdeC (fromIntegral $ maxFail opts)
(fromIntegral $ maxNumSteps opts) (coerce $ minStep opts)
(fromIntegral $ getMethod method) (coerce initStepSize) jacH (scise control)
(coerce f) (coerce y0) (coerce tt) of
Left (v, c) -> Left (reshape l (coerce v), fromIntegral c)
Right (v, d) -> Right (reshape l (coerce v), d)
where
l = size y0
scise (X aTol rTol) = coerce (V.replicate l aTol, rTol)
scise (X' aTol rTol) = coerce (V.replicate l aTol, rTol)
scise (XX' aTol rTol yScale _yDotScale) = coerce (V.replicate l aTol, yScale * rTol)
-- FIXME; Should we check that the length of ss is correct?
scise (ScXX' aTol rTol yScale _yDotScale ss) = coerce (V.map (* aTol) ss, yScale * rTol)
jacH = fmap (\g t v -> matrixToSunMatrix $ g (coerce t) (coerce v)) $
getJacobian method
matrixToSunMatrix m = T.SunMatrix { T.rows = nr, T.cols = nc, T.vals = vs }
where
nr = fromIntegral $ rows m
nc = fromIntegral $ cols m
-- FIXME: efficiency
vs = V.fromList $ map coerce $ concat $ toLists m
solveOdeC ::
CInt ->
CLong ->
CDouble ->
CInt ->
Maybe CDouble ->
(Maybe (CDouble -> V.Vector CDouble -> T.SunMatrix)) ->
(V.Vector CDouble, CDouble) ->
(CDouble -> V.Vector CDouble -> V.Vector CDouble) -- ^ The RHS of the system \(\dot{y} = f(t,y)\)
-> V.Vector CDouble -- ^ Initial conditions
-> V.Vector CDouble -- ^ Desired solution times
-> Either (V.Vector CDouble, CInt) (V.Vector CDouble, SundialsDiagnostics) -- ^ Partial solution and error code or
-- solution and diagnostics
solveOdeC maxErrTestFails maxNumSteps_ minStep_ method initStepSize
jacH (aTols, rTol) fun f0 ts = unsafePerformIO $ do
let isInitStepSize :: CInt
isInitStepSize = fromIntegral $ fromEnum $ isJust initStepSize
ss :: CDouble
ss = case initStepSize of
-- It would be better to put an error message here but
-- inline-c seems to evaluate this even if it is never
-- used :(
Nothing -> 0.0
Just x -> x
let dim = V.length f0
nEq :: CLong
nEq = fromIntegral dim
nTs :: CInt
nTs = fromIntegral $ V.length ts
quasiMatrixRes <- createVector ((fromIntegral dim) * (fromIntegral nTs))
qMatMut <- V.thaw quasiMatrixRes
diagnostics :: V.Vector CLong <- createVector 10 -- FIXME
diagMut <- V.thaw diagnostics
-- We need the types that sundials expects. These are tied together
-- in 'CLangToHaskellTypes'. FIXME: The Haskell type is currently empty!
let funIO :: CDouble -> Ptr T.SunVector -> Ptr T.SunVector -> Ptr () -> IO CInt
funIO x y f _ptr = do
-- Convert the pointer we get from C (y) to a vector, and then
-- apply the user-supplied function.
fImm <- fun x <$> getDataFromContents dim y
-- Fill in the provided pointer with the resulting vector.
putDataInContents fImm dim f
-- FIXME: I don't understand what this comment means
-- Unsafe since the function will be called many times.
[CU.exp| int{ 0 } |]
let isJac :: CInt
isJac = fromIntegral $ fromEnum $ isJust jacH
jacIO :: CDouble -> Ptr T.SunVector -> Ptr T.SunVector -> Ptr T.SunMatrix ->
Ptr () -> Ptr T.SunVector -> Ptr T.SunVector -> Ptr T.SunVector ->
IO CInt
jacIO t y _fy jacS _ptr _tmp1 _tmp2 _tmp3 = do
case jacH of
Nothing -> error "Numeric.Sundials.ARKode.ODE: Jacobian not defined"
Just jacI -> do j <- jacI t <$> getDataFromContents dim y
poke jacS j
-- FIXME: I don't understand what this comment means
-- Unsafe since the function will be called many times.
[CU.exp| int{ 0 } |]
res <- [C.block| int {
/* general problem variables */
int flag; /* reusable error-checking flag */
int i, j; /* reusable loop indices */
N_Vector y = NULL; /* empty vector for storing solution */
N_Vector tv = NULL; /* empty vector for storing absolute tolerances */
SUNMatrix A = NULL; /* empty matrix for linear solver */
SUNLinearSolver LS = NULL; /* empty linear solver object */
void *arkode_mem = NULL; /* empty ARKode memory structure */
realtype t;
long int nst, nst_a, nfe, nfi, nsetups, nje, nfeLS, nni, ncfn, netf;
/* general problem parameters */
realtype T0 = RCONST(($vec-ptr:(double *ts))[0]); /* initial time */
sunindextype NEQ = $(sunindextype nEq); /* number of dependent vars. */
/* Initialize data structures */
y = N_VNew_Serial(NEQ); /* Create serial vector for solution */
if (check_flag((void *)y, "N_VNew_Serial", 0)) return 1;
/* Specify initial condition */
for (i = 0; i < NEQ; i++) {
NV_Ith_S(y,i) = ($vec-ptr:(double *f0))[i];
};
tv = N_VNew_Serial(NEQ); /* Create serial vector for absolute tolerances */
if (check_flag((void *)tv, "N_VNew_Serial", 0)) return 1;
/* Specify tolerances */
for (i = 0; i < NEQ; i++) {
NV_Ith_S(tv,i) = ($vec-ptr:(double *aTols))[i];
};
arkode_mem = ARKodeCreate(); /* Create the solver memory */
if (check_flag((void *)arkode_mem, "ARKodeCreate", 0)) return 1;
/* Call ARKodeInit to initialize the integrator memory and specify the */
/* right-hand side function in y'=f(t,y), the inital time T0, and */
/* the initial dependent variable vector y. Note: we treat the */
/* problem as fully implicit and set f_E to NULL and f_I to f. */
/* Here we use the C types defined in helpers.h which tie up with */
/* the Haskell types defined in CLangToHaskellTypes */
if ($(int method) < MIN_DIRK_NUM) {
flag = ARKodeInit(arkode_mem, $fun:(int (* funIO) (double t, SunVector y[], SunVector dydt[], void * params)), NULL, T0, y);
if (check_flag(&flag, "ARKodeInit", 1)) return 1;
} else {
flag = ARKodeInit(arkode_mem, NULL, $fun:(int (* funIO) (double t, SunVector y[], SunVector dydt[], void * params)), T0, y);
if (check_flag(&flag, "ARKodeInit", 1)) return 1;
}
flag = ARKodeSetMinStep(arkode_mem, $(double minStep_));
if (check_flag(&flag, "ARKodeSetMinStep", 1)) return 1;
flag = ARKodeSetMaxNumSteps(arkode_mem, $(long int maxNumSteps_));
if (check_flag(&flag, "ARKodeSetMaxNumSteps", 1)) return 1;
flag = ARKodeSetMaxErrTestFails(arkode_mem, $(int maxErrTestFails));
if (check_flag(&flag, "ARKodeSetMaxErrTestFails", 1)) return 1;
/* Set routines */
flag = ARKodeSVtolerances(arkode_mem, $(double rTol), tv);
if (check_flag(&flag, "ARKodeSVtolerances", 1)) return 1;
/* Initialize dense matrix data structure and solver */
A = SUNDenseMatrix(NEQ, NEQ);
if (check_flag((void *)A, "SUNDenseMatrix", 0)) return 1;
LS = SUNDenseLinearSolver(y, A);
if (check_flag((void *)LS, "SUNDenseLinearSolver", 0)) return 1;
/* Attach matrix and linear solver */
flag = ARKDlsSetLinearSolver(arkode_mem, LS, A);
if (check_flag(&flag, "ARKDlsSetLinearSolver", 1)) return 1;
/* Set the initial step size if there is one */
if ($(int isInitStepSize)) {
/* FIXME: We could check if the initial step size is 0 */
/* or even NaN and then throw an error */
flag = ARKodeSetInitStep(arkode_mem, $(double ss));
if (check_flag(&flag, "ARKodeSetInitStep", 1)) return 1;
}
/* Set the Jacobian if there is one */
if ($(int isJac)) {
flag = ARKDlsSetJacFn(arkode_mem, $fun:(int (* jacIO) (double t, SunVector y[], SunVector fy[], SunMatrix Jac[], void * params, SunVector tmp1[], SunVector tmp2[], SunVector tmp3[])));
if (check_flag(&flag, "ARKDlsSetJacFn", 1)) return 1;
}
/* Store initial conditions */
for (j = 0; j < NEQ; j++) {
($vec-ptr:(double *qMatMut))[0 * $(int nTs) + j] = NV_Ith_S(y,j);
}
/* Explicitly set the method */
if ($(int method) >= MIN_DIRK_NUM) {
flag = ARKodeSetIRKTableNum(arkode_mem, $(int method));
if (check_flag(&flag, "ARKodeSetIRKTableNum", 1)) return 1;
} else {
flag = ARKodeSetERKTableNum(arkode_mem, $(int method));
if (check_flag(&flag, "ARKodeSetERKTableNum", 1)) return 1;
}
/* Main time-stepping loop: calls ARKode to perform the integration */
/* Stops when the final time has been reached */
for (i = 1; i < $(int nTs); i++) {
flag = ARKode(arkode_mem, ($vec-ptr:(double *ts))[i], y, &t, ARK_NORMAL); /* call integrator */
if (check_flag(&flag, "ARKode solver failure, stopping integration", 1)) return 1;
/* Store the results for Haskell */
for (j = 0; j < NEQ; j++) {
($vec-ptr:(double *qMatMut))[i * NEQ + j] = NV_Ith_S(y,j);
}
}
/* Get some final statistics on how the solve progressed */
flag = ARKodeGetNumSteps(arkode_mem, &nst);
check_flag(&flag, "ARKodeGetNumSteps", 1);
($vec-ptr:(long int *diagMut))[0] = nst;
flag = ARKodeGetNumStepAttempts(arkode_mem, &nst_a);
check_flag(&flag, "ARKodeGetNumStepAttempts", 1);
($vec-ptr:(long int *diagMut))[1] = nst_a;
flag = ARKodeGetNumRhsEvals(arkode_mem, &nfe, &nfi);
check_flag(&flag, "ARKodeGetNumRhsEvals", 1);
($vec-ptr:(long int *diagMut))[2] = nfe;
($vec-ptr:(long int *diagMut))[3] = nfi;
flag = ARKodeGetNumLinSolvSetups(arkode_mem, &nsetups);
check_flag(&flag, "ARKodeGetNumLinSolvSetups", 1);
($vec-ptr:(long int *diagMut))[4] = nsetups;
flag = ARKodeGetNumErrTestFails(arkode_mem, &netf);
check_flag(&flag, "ARKodeGetNumErrTestFails", 1);
($vec-ptr:(long int *diagMut))[5] = netf;
flag = ARKodeGetNumNonlinSolvIters(arkode_mem, &nni);
check_flag(&flag, "ARKodeGetNumNonlinSolvIters", 1);
($vec-ptr:(long int *diagMut))[6] = nni;
flag = ARKodeGetNumNonlinSolvConvFails(arkode_mem, &ncfn);
check_flag(&flag, "ARKodeGetNumNonlinSolvConvFails", 1);
($vec-ptr:(long int *diagMut))[7] = ncfn;
flag = ARKDlsGetNumJacEvals(arkode_mem, &nje);
check_flag(&flag, "ARKDlsGetNumJacEvals", 1);
($vec-ptr:(long int *diagMut))[8] = ncfn;
flag = ARKDlsGetNumRhsEvals(arkode_mem, &nfeLS);
check_flag(&flag, "ARKDlsGetNumRhsEvals", 1);
($vec-ptr:(long int *diagMut))[9] = ncfn;
/* Clean up and return */
N_VDestroy(y); /* Free y vector */
N_VDestroy(tv); /* Free tv vector */
ARKodeFree(&arkode_mem); /* Free integrator memory */
SUNLinSolFree(LS); /* Free linear solver */
SUNMatDestroy(A); /* Free A matrix */
return flag;
} |]
preD <- V.freeze diagMut
let d = SundialsDiagnostics (fromIntegral $ preD V.!0)
(fromIntegral $ preD V.!1)
(fromIntegral $ preD V.!2)
(fromIntegral $ preD V.!3)
(fromIntegral $ preD V.!4)
(fromIntegral $ preD V.!5)
(fromIntegral $ preD V.!6)
(fromIntegral $ preD V.!7)
(fromIntegral $ preD V.!8)
(fromIntegral $ preD V.!9)
m <- V.freeze qMatMut
if res == 0
then do
return $ Right (m, d)
else do
return $ Left (m, res)
data ButcherTable = ButcherTable { am :: Matrix Double
, cv :: Vector Double
, bv :: Vector Double
, b2v :: Vector Double
}
deriving Show
data ButcherTable' a = ButcherTable' { am' :: V.Vector a
, cv' :: V.Vector a
, bv' :: V.Vector a
, b2v' :: V.Vector a
}
deriving Show
butcherTable :: ODEMethod -> ButcherTable
butcherTable method =
case getBT method of
Left c -> error $ show c -- FIXME
Right (ButcherTable' v w x y, sqp) ->
ButcherTable { am = subMatrix (0, 0) (s, s) $ (arkSMax >< arkSMax) (V.toList v)
, cv = subVector 0 s w
, bv = subVector 0 s x
, b2v = subVector 0 s y
}
where
s = fromIntegral $ sqp V.! 0
getBT :: ODEMethod -> Either Int (ButcherTable' Double, V.Vector Int)
getBT method = case getButcherTable method of
Left c ->
Left $ fromIntegral c
Right (ButcherTable' a b c d, sqp) ->
Right $ ( ButcherTable' (coerce a) (coerce b) (coerce c) (coerce d)
, V.map fromIntegral sqp )
getButcherTable :: ODEMethod
-> Either CInt (ButcherTable' CDouble, V.Vector CInt)
getButcherTable method = unsafePerformIO $ do
-- ARKode seems to want an ODE in order to set and then get the
-- Butcher tableau so here's one to keep it happy
let funI :: CDouble -> V.Vector CDouble -> V.Vector CDouble
funI _t ys = V.fromList [ ys V.! 0 ]
let funE :: CDouble -> V.Vector CDouble -> V.Vector CDouble
funE _t ys = V.fromList [ ys V.! 0 ]
f0 = V.fromList [ 1.0 ]
ts = V.fromList [ 0.0 ]
dim = V.length f0
nEq :: CLong
nEq = fromIntegral dim
mN :: CInt
mN = fromIntegral $ getMethod method
btSQP :: V.Vector CInt <- createVector 3
btSQPMut <- V.thaw btSQP
btAs :: V.Vector CDouble <- createVector (arkSMax * arkSMax)
btAsMut <- V.thaw btAs
btCs :: V.Vector CDouble <- createVector arkSMax
btBs :: V.Vector CDouble <- createVector arkSMax
btB2s :: V.Vector CDouble <- createVector arkSMax
btCsMut <- V.thaw btCs
btBsMut <- V.thaw btBs
btB2sMut <- V.thaw btB2s
let funIOI :: CDouble -> Ptr T.SunVector -> Ptr T.SunVector -> Ptr () -> IO CInt
funIOI x y f _ptr = do
fImm <- funI x <$> getDataFromContents dim y
putDataInContents fImm dim f
-- FIXME: I don't understand what this comment means
-- Unsafe since the function will be called many times.
[CU.exp| int{ 0 } |]
let funIOE :: CDouble -> Ptr T.SunVector -> Ptr T.SunVector -> Ptr () -> IO CInt
funIOE x y f _ptr = do
fImm <- funE x <$> getDataFromContents dim y
putDataInContents fImm dim f
-- FIXME: I don't understand what this comment means
-- Unsafe since the function will be called many times.
[CU.exp| int{ 0 } |]
res <- [C.block| int {
/* general problem variables */
int flag; /* reusable error-checking flag */
N_Vector y = NULL; /* empty vector for storing solution */
void *arkode_mem = NULL; /* empty ARKode memory structure */
int i, j; /* reusable loop indices */
/* general problem parameters */
realtype T0 = RCONST(($vec-ptr:(double *ts))[0]); /* initial time */
sunindextype NEQ = $(sunindextype nEq); /* number of dependent vars */
/* Initialize data structures */
y = N_VNew_Serial(NEQ); /* Create serial vector for solution */
if (check_flag((void *)y, "N_VNew_Serial", 0)) return 1;
/* Specify initial condition */
for (i = 0; i < NEQ; i++) {
NV_Ith_S(y,i) = ($vec-ptr:(double *f0))[i];
};
arkode_mem = ARKodeCreate(); /* Create the solver memory */
if (check_flag((void *)arkode_mem, "ARKodeCreate", 0)) return 1;
flag = ARKodeInit(arkode_mem, $fun:(int (* funIOE) (double t, SunVector y[], SunVector dydt[], void * params)), $fun:(int (* funIOI) (double t, SunVector y[], SunVector dydt[], void * params)), T0, y);
if (check_flag(&flag, "ARKodeInit", 1)) return 1;
if ($(int mN) >= MIN_DIRK_NUM) {
flag = ARKodeSetIRKTableNum(arkode_mem, $(int mN));
if (check_flag(&flag, "ARKodeSetIRKTableNum", 1)) return 1;
} else {
flag = ARKodeSetERKTableNum(arkode_mem, $(int mN));
if (check_flag(&flag, "ARKodeSetERKTableNum", 1)) return 1;
}
int s, q, p;
realtype *ai = (realtype *)malloc(ARK_S_MAX * ARK_S_MAX * sizeof(realtype));
realtype *ae = (realtype *)malloc(ARK_S_MAX * ARK_S_MAX * sizeof(realtype));
realtype *ci = (realtype *)malloc(ARK_S_MAX * sizeof(realtype));
realtype *ce = (realtype *)malloc(ARK_S_MAX * sizeof(realtype));
realtype *bi = (realtype *)malloc(ARK_S_MAX * sizeof(realtype));
realtype *be = (realtype *)malloc(ARK_S_MAX * sizeof(realtype));
realtype *b2i = (realtype *)malloc(ARK_S_MAX * sizeof(realtype));
realtype *b2e = (realtype *)malloc(ARK_S_MAX * sizeof(realtype));
flag = ARKodeGetCurrentButcherTables(arkode_mem, &s, &q, &p, ai, ae, ci, ce, bi, be, b2i, b2e);
if (check_flag(&flag, "ARKode", 1)) return 1;
$vec-ptr:(int *btSQPMut)[0] = s;
$vec-ptr:(int *btSQPMut)[1] = q;
$vec-ptr:(int *btSQPMut)[2] = p;
for (i = 0; i < s; i++) {
for (j = 0; j < s; j++) {
/* FIXME: double should be realtype */
($vec-ptr:(double *btAsMut))[i * ARK_S_MAX + j] = ai[i * ARK_S_MAX + j];
}
}
for (i = 0; i < s; i++) {
($vec-ptr:(double *btCsMut))[i] = ci[i];
($vec-ptr:(double *btBsMut))[i] = bi[i];
($vec-ptr:(double *btB2sMut))[i] = b2i[i];
}
/* Clean up and return */
N_VDestroy(y); /* Free y vector */
ARKodeFree(&arkode_mem); /* Free integrator memory */
return flag;
} |]
if res == 0
then do
x <- V.freeze btAsMut
y <- V.freeze btSQPMut
z <- V.freeze btCsMut
u <- V.freeze btBsMut
v <- V.freeze btB2sMut
return $ Right (ButcherTable' { am' = x, cv' = z, bv' = u, b2v' = v }, y)
else do
return $ Left res
-- | Adaptive step-size control
-- functions.
--
-- [GSL](https://www.gnu.org/software/gsl/doc/html/ode-initval.html#adaptive-step-size-control)
-- allows the user to control the step size adjustment using
-- \(D_i = \epsilon^{abs}s_i + \epsilon^{rel}(a_{y} |y_i| + a_{dy/dt} h |\dot{y}_i|)\) where
-- \(\epsilon^{abs}\) is the required absolute error, \(\epsilon^{rel}\)
-- is the required relative error, \(s_i\) is a vector of scaling
-- factors, \(a_{y}\) is a scaling factor for the solution \(y\) and
-- \(a_{dydt}\) is a scaling factor for the derivative of the solution \(dy/dt\).
--
-- [ARKode](https://computation.llnl.gov/projects/sundials/arkode)
-- allows the user to control the step size adjustment using
-- \(\eta^{rel}|y_i| + \eta^{abs}_i\). For compatibility with
-- [hmatrix-gsl](https://hackage.haskell.org/package/hmatrix-gsl),
-- tolerances for \(y\) and \(\dot{y}\) can be specified but the latter have no
-- effect.
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
| 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
| XX' Double Double Double Double -- ^ include both via relative tolerance
-- scaling factors \(a_y\), \(a_{{dy}/{dt}}\); in ARKode terms, the latter is ignored and \(\eta^{rel} = a_{y}\epsilon^{rel}\)
| 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}\)
|