-
Notifications
You must be signed in to change notification settings - Fork 132
/
arkode.c
3370 lines (2945 loc) · 109 KB
/
arkode.c
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
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
/*---------------------------------------------------------------
* Programmer(s): Daniel R. Reynolds @ SMU
*---------------------------------------------------------------
* SUNDIALS Copyright Start
* Copyright (c) 2002-2024, Lawrence Livermore National Security
* and Southern Methodist University.
* All rights reserved.
*
* See the top-level LICENSE and NOTICE files for details.
*
* SPDX-License-Identifier: BSD-3-Clause
* SUNDIALS Copyright End
*---------------------------------------------------------------
* This is the implementation file for the main ARKODE
* infrastructure. It is independent of the ARKODE time step
* module, nonlinear solver, linear solver and vector modules in
* use.
*--------------------------------------------------------------*/
/*===============================================================
Import Header Files
===============================================================*/
#include "arkode/arkode.h"
#include <stdarg.h>
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <sunadaptcontroller/sunadaptcontroller_soderlind.h>
#include <sundials/sundials_config.h>
#include <sundials/sundials_math.h>
#include <sundials/sundials_types.h>
#include "arkode_impl.h"
#include "arkode_interp_impl.h"
#include "sundials/priv/sundials_errors_impl.h"
#include "sundials/sundials_context.h"
#include "sundials/sundials_logger.h"
#include "sundials_utils.h"
/*===============================================================
EXPORTED FUNCTIONS
===============================================================*/
/*---------------------------------------------------------------
arkCreate:
arkCreate creates an internal memory block for a problem to
be solved by a time step module built on ARKODE. If successful,
arkCreate returns a pointer to the problem memory. If an
initialization error occurs, arkCreate prints an error message
to standard err and returns NULL.
---------------------------------------------------------------*/
ARKodeMem arkCreate(SUNContext sunctx)
{
int iret;
long int lenrw, leniw;
ARKodeMem ark_mem;
if (!sunctx)
{
arkProcessError(NULL, ARK_ILL_INPUT, __LINE__, __func__, __FILE__,
MSG_ARK_NULL_SUNCTX);
return (NULL);
}
ark_mem = NULL;
ark_mem = (ARKodeMem)malloc(sizeof(struct ARKodeMemRec));
if (ark_mem == NULL)
{
arkProcessError(NULL, ARK_MEM_FAIL, __LINE__, __func__, __FILE__,
MSG_ARK_ARKMEM_FAIL);
return (NULL);
}
/* Zero out ark_mem */
memset(ark_mem, 0, sizeof(struct ARKodeMemRec));
/* Set the context */
ark_mem->sunctx = sunctx;
/* Set uround */
ark_mem->uround = SUN_UNIT_ROUNDOFF;
/* Initialize time step module to NULL */
ark_mem->step_attachlinsol = NULL;
ark_mem->step_attachmasssol = NULL;
ark_mem->step_disablelsetup = NULL;
ark_mem->step_disablemsetup = NULL;
ark_mem->step_getlinmem = NULL;
ark_mem->step_getmassmem = NULL;
ark_mem->step_getimplicitrhs = NULL;
ark_mem->step_mmult = NULL;
ark_mem->step_getgammas = NULL;
ark_mem->step_init = NULL;
ark_mem->step_fullrhs = NULL;
ark_mem->step = NULL;
ark_mem->step_mem = NULL;
/* Initialize root finding variables */
ark_mem->root_mem = NULL;
/* Initialize inequality constraints variables */
ark_mem->constraintsSet = SUNFALSE;
ark_mem->constraints = NULL;
/* Initialize relaxation variables */
ark_mem->relax_enabled = SUNFALSE;
ark_mem->relax_mem = NULL;
/* Initialize lrw and liw */
ark_mem->lrw = 18;
ark_mem->liw = 41; /* fcn/data ptr, int, long int, sunindextype, sunbooleantype */
/* No mallocs have been done yet */
ark_mem->VabstolMallocDone = SUNFALSE;
ark_mem->VRabstolMallocDone = SUNFALSE;
ark_mem->MallocDone = SUNFALSE;
/* No user-supplied step postprocessing function yet */
ark_mem->ProcessStep = NULL;
ark_mem->ps_data = NULL;
/* No user-supplied stage postprocessing function yet */
ark_mem->ProcessStage = NULL;
/* No user_data pointer yet */
ark_mem->user_data = NULL;
/* Allocate step adaptivity structure and note storage */
ark_mem->hadapt_mem = arkAdaptInit();
if (ark_mem->hadapt_mem == NULL)
{
arkProcessError(NULL, ARK_MEM_FAIL, __LINE__, __func__, __FILE__,
"Allocation of step adaptivity structure failed");
return (NULL);
}
ark_mem->lrw += ARK_ADAPT_LRW;
ark_mem->liw += ARK_ADAPT_LIW;
/* Allocate default step controller (PID) and note storage */
ark_mem->hadapt_mem->hcontroller = SUNAdaptController_PID(sunctx);
if (ark_mem->hadapt_mem->hcontroller == NULL)
{
arkProcessError(NULL, ARK_MEM_FAIL, __LINE__, __func__, __FILE__,
"Allocation of step controller object failed");
return (NULL);
}
ark_mem->hadapt_mem->owncontroller = SUNTRUE;
(void)SUNAdaptController_Space(ark_mem->hadapt_mem->hcontroller, &lenrw,
&leniw);
ark_mem->lrw += lenrw;
ark_mem->liw += leniw;
/* Initialize the interpolation structure to NULL */
ark_mem->interp = NULL;
ark_mem->interp_type = -1;
/* Initially, rwt should point to ewt */
ark_mem->rwt_is_ewt = SUNTRUE;
/* Indicate that calling the full RHS function is not required, this flag is
updated to SUNTRUE by the interpolation module initialization function
and/or the stepper initialization function in arkInitialSetup */
ark_mem->call_fullrhs = SUNFALSE;
/* Indicate that the problem needs to be initialized */
ark_mem->initsetup = SUNTRUE;
ark_mem->init_type = FIRST_INIT;
ark_mem->firststage = SUNTRUE;
ark_mem->initialized = SUNFALSE;
/* Initial step size has not been determined yet */
ark_mem->h = ZERO;
ark_mem->h0u = ZERO;
/* Set default values for integrator optional inputs */
iret = arkSetDefaults(ark_mem);
if (iret != ARK_SUCCESS)
{
arkProcessError(NULL, 0, __LINE__, __func__, __FILE__,
"Error setting default solver options");
return (NULL);
}
/* Return pointer to ARKODE memory block */
return (ark_mem);
}
/*---------------------------------------------------------------
arkResize:
arkResize re-initializes ARKODE's memory for a problem with a
changing vector size. It is assumed that the problem dynamics
before and after the vector resize will be comparable, so that
all time-stepping heuristics prior to calling arkResize
remain valid after the call. If instead the dynamics should be
re-calibrated, the ARKODE memory structure should be deleted
with a call to *StepFree, and re-created with a call to
*StepCreate.
To aid in the vector-resize operation, the user can supply a
vector resize function, that will take as input an N_Vector with
the previous size, and return as output a corresponding vector
of the new size. If this function (of type ARKVecResizeFn) is
not supplied (i.e. is set to NULL), then all existing N_Vectors
will be destroyed and re-cloned from the input vector.
In the case that the dynamical time scale should be modified
slightly from the previous time scale, an input "hscale" is
allowed, that will re-scale the upcoming time step by the
specified factor. If a value <= 0 is specified, the default of
1.0 will be used.
Other arguments:
ark_mem Existing ARKODE memory data structure.
y0 The newly-sized solution vector, holding
the current dependent variable values.
t0 The current value of the independent
variable.
resize_data User-supplied data structure that will be
passed to the supplied resize function.
The return value is ARK_SUCCESS = 0 if no errors occurred, or
a negative value otherwise.
---------------------------------------------------------------*/
int arkResize(ARKodeMem ark_mem, N_Vector y0, sunrealtype hscale,
sunrealtype t0, ARKVecResizeFn resize, void* resize_data)
{
sunbooleantype resizeOK;
sunindextype lrw1, liw1, lrw_diff, liw_diff;
int retval;
/* Check ark_mem */
if (ark_mem == NULL)
{
arkProcessError(NULL, ARK_MEM_NULL, __LINE__, __func__, __FILE__,
MSG_ARK_NO_MEM);
return (ARK_MEM_NULL);
}
/* Check if ark_mem was allocated */
if (ark_mem->MallocDone == SUNFALSE)
{
arkProcessError(ark_mem, ARK_NO_MALLOC, __LINE__, __func__, __FILE__,
MSG_ARK_NO_MALLOC);
return (ARK_NO_MALLOC);
}
/* Check for legal input parameters */
if (y0 == NULL)
{
arkProcessError(ark_mem, ARK_ILL_INPUT, __LINE__, __func__, __FILE__,
MSG_ARK_NULL_Y0);
return (ARK_ILL_INPUT);
}
/* Copy the input parameters into ARKODE state */
ark_mem->tcur = t0;
ark_mem->tn = t0;
/* Update time-stepping parameters */
/* adjust upcoming step size depending on hscale */
if (hscale < ZERO) { hscale = ONE; }
if (hscale != ONE)
{
/* Encode hscale into ark_mem structure */
ark_mem->eta = hscale;
ark_mem->hprime *= hscale;
/* If next step would overtake tstop, adjust stepsize */
if (ark_mem->tstopset)
{
if ((ark_mem->tcur + ark_mem->hprime - ark_mem->tstop) * ark_mem->hprime >
ZERO)
{
ark_mem->hprime = (ark_mem->tstop - ark_mem->tcur) *
(ONE - FOUR * ark_mem->uround);
ark_mem->eta = ark_mem->hprime / ark_mem->h;
}
}
}
/* Determing change in vector sizes */
lrw1 = liw1 = 0;
if (y0->ops->nvspace != NULL) { N_VSpace(y0, &lrw1, &liw1); }
lrw_diff = lrw1 - ark_mem->lrw1;
liw_diff = liw1 - ark_mem->liw1;
ark_mem->lrw1 = lrw1;
ark_mem->liw1 = liw1;
/* Resize the solver vectors (using y0 as a template) */
resizeOK = arkResizeVectors(ark_mem, resize, resize_data, lrw_diff, liw_diff,
y0);
if (!resizeOK)
{
arkProcessError(ark_mem, ARK_MEM_FAIL, __LINE__, __func__, __FILE__,
"Unable to resize vector");
return (ARK_MEM_FAIL);
}
/* Resize the interpolation structure memory */
if (ark_mem->interp != NULL)
{
retval = arkInterpResize(ark_mem, ark_mem->interp, resize, resize_data,
lrw_diff, liw_diff, y0);
if (retval != ARK_SUCCESS)
{
arkProcessError(ark_mem, retval, __LINE__, __func__, __FILE__,
"Interpolation module resize failure");
return (retval);
}
}
/* Copy y0 into ark_yn to set the current solution */
N_VScale(ONE, y0, ark_mem->yn);
ark_mem->fn_is_current = SUNFALSE;
/* Disable constraints */
ark_mem->constraintsSet = SUNFALSE;
/* Indicate that problem needs to be initialized */
ark_mem->initsetup = SUNTRUE;
ark_mem->init_type = RESIZE_INIT;
ark_mem->firststage = SUNTRUE;
/* Problem has been successfully re-sized */
return (ARK_SUCCESS);
}
/*---------------------------------------------------------------
arkSStolerances, arkSVtolerances, arkWFtolerances:
These functions specify the integration tolerances. One of them
SHOULD be called before the first call to arkEvolve; otherwise
default values of reltol=1e-4 and abstol=1e-9 will be used,
which may be entirely incorrect for a specific problem.
arkSStolerances specifies scalar relative and absolute
tolerances.
arkSVtolerances specifies scalar relative tolerance and a
vector absolute tolerance (a potentially different absolute
tolerance for each vector component).
arkWFtolerances specifies a user-provides function (of type
ARKEwtFn) which will be called to set the error weight vector.
---------------------------------------------------------------*/
int arkSStolerances(ARKodeMem ark_mem, sunrealtype reltol, sunrealtype abstol)
{
/* Check inputs */
if (ark_mem == NULL)
{
arkProcessError(NULL, ARK_MEM_NULL, __LINE__, __func__, __FILE__,
MSG_ARK_NO_MEM);
return (ARK_MEM_NULL);
}
if (ark_mem->MallocDone == SUNFALSE)
{
arkProcessError(ark_mem, ARK_NO_MALLOC, __LINE__, __func__, __FILE__,
MSG_ARK_NO_MALLOC);
return (ARK_NO_MALLOC);
}
if (reltol < ZERO)
{
arkProcessError(ark_mem, ARK_ILL_INPUT, __LINE__, __func__, __FILE__,
MSG_ARK_BAD_RELTOL);
return (ARK_ILL_INPUT);
}
if (abstol < ZERO)
{
arkProcessError(ark_mem, ARK_ILL_INPUT, __LINE__, __func__, __FILE__,
MSG_ARK_BAD_ABSTOL);
return (ARK_ILL_INPUT);
}
/* Set flag indicating whether abstol == 0 */
ark_mem->atolmin0 = (abstol == ZERO);
/* Copy tolerances into memory */
ark_mem->reltol = reltol;
ark_mem->Sabstol = abstol;
ark_mem->itol = ARK_SS;
/* enforce use of arkEwtSetSS */
ark_mem->user_efun = SUNFALSE;
ark_mem->efun = arkEwtSetSS;
ark_mem->e_data = ark_mem;
return (ARK_SUCCESS);
}
int arkSVtolerances(ARKodeMem ark_mem, sunrealtype reltol, N_Vector abstol)
{
/* local variables */
sunrealtype abstolmin;
/* Check inputs */
if (ark_mem == NULL)
{
arkProcessError(NULL, ARK_MEM_NULL, __LINE__, __func__, __FILE__,
MSG_ARK_NO_MEM);
return (ARK_MEM_NULL);
}
if (ark_mem->MallocDone == SUNFALSE)
{
arkProcessError(ark_mem, ARK_NO_MALLOC, __LINE__, __func__, __FILE__,
MSG_ARK_NO_MALLOC);
return (ARK_NO_MALLOC);
}
if (reltol < ZERO)
{
arkProcessError(ark_mem, ARK_ILL_INPUT, __LINE__, __func__, __FILE__,
MSG_ARK_BAD_RELTOL);
return (ARK_ILL_INPUT);
}
if (abstol == NULL)
{
arkProcessError(ark_mem, ARK_ILL_INPUT, __LINE__, __func__, __FILE__,
MSG_ARK_NULL_ABSTOL);
return (ARK_ILL_INPUT);
}
if (abstol->ops->nvmin == NULL)
{
arkProcessError(ark_mem, ARK_ILL_INPUT, __LINE__, __func__, __FILE__,
"Missing N_VMin routine from N_Vector");
return (ARK_ILL_INPUT);
}
abstolmin = N_VMin(abstol);
if (abstolmin < ZERO)
{
arkProcessError(ark_mem, ARK_ILL_INPUT, __LINE__, __func__, __FILE__,
MSG_ARK_BAD_ABSTOL);
return (ARK_ILL_INPUT);
}
/* Set flag indicating whether min(abstol) == 0 */
ark_mem->atolmin0 = (abstolmin == ZERO);
/* Copy tolerances into memory */
if (!(ark_mem->VabstolMallocDone))
{
if (!arkAllocVec(ark_mem, ark_mem->ewt, &(ark_mem->Vabstol)))
{
arkProcessError(ark_mem, ARK_MEM_FAIL, __LINE__, __func__, __FILE__,
MSG_ARK_ARKMEM_FAIL);
return (ARK_ILL_INPUT);
}
ark_mem->VabstolMallocDone = SUNTRUE;
}
N_VScale(ONE, abstol, ark_mem->Vabstol);
ark_mem->reltol = reltol;
ark_mem->itol = ARK_SV;
/* enforce use of arkEwtSetSV */
ark_mem->user_efun = SUNFALSE;
ark_mem->efun = arkEwtSetSV;
ark_mem->e_data = ark_mem;
return (ARK_SUCCESS);
}
int arkWFtolerances(ARKodeMem ark_mem, ARKEwtFn efun)
{
if (ark_mem == NULL)
{
arkProcessError(NULL, ARK_MEM_NULL, __LINE__, __func__, __FILE__,
MSG_ARK_NO_MEM);
return (ARK_MEM_NULL);
}
if (ark_mem->MallocDone == SUNFALSE)
{
arkProcessError(ark_mem, ARK_NO_MALLOC, __LINE__, __func__, __FILE__,
MSG_ARK_NO_MALLOC);
return (ARK_NO_MALLOC);
}
/* Copy tolerance data into memory */
ark_mem->itol = ARK_WF;
ark_mem->user_efun = SUNTRUE;
ark_mem->efun = efun;
ark_mem->e_data = ark_mem->user_data;
return (ARK_SUCCESS);
}
/*---------------------------------------------------------------
arkResStolerance, arkResVtolerance, arkResFtolerance:
These functions specify the absolute residual tolerance.
Specification of the absolute residual tolerance is only
necessary for problems with non-identity mass matrices in which
the units of the solution vector y dramatically differ from the
units of the ODE right-hand side f(t,y). If this occurs, one
of these routines SHOULD be called before the first call to
ARKODE; otherwise the default value of rabstol=1e-9 will be
used, which may be entirely incorrect for a specific problem.
arkResStolerances specifies a scalar residual tolerance.
arkResVtolerances specifies a vector residual tolerance
(a potentially different absolute residual tolerance for
each vector component).
arkResFtolerances specifies a user-provides function (of
type ARKRwtFn) which will be called to set the residual
weight vector.
---------------------------------------------------------------*/
int arkResStolerance(ARKodeMem ark_mem, sunrealtype rabstol)
{
/* Check inputs */
if (ark_mem == NULL)
{
arkProcessError(NULL, ARK_MEM_NULL, __LINE__, __func__, __FILE__,
MSG_ARK_NO_MEM);
return (ARK_MEM_NULL);
}
if (ark_mem->MallocDone == SUNFALSE)
{
arkProcessError(ark_mem, ARK_NO_MALLOC, __LINE__, __func__, __FILE__,
MSG_ARK_NO_MALLOC);
return (ARK_NO_MALLOC);
}
if (rabstol < ZERO)
{
arkProcessError(ark_mem, ARK_ILL_INPUT, __LINE__, __func__, __FILE__,
MSG_ARK_BAD_RABSTOL);
return (ARK_ILL_INPUT);
}
/* Set flag indicating whether rabstol == 0 */
ark_mem->Ratolmin0 = (rabstol == ZERO);
/* Allocate space for rwt if necessary */
if (ark_mem->rwt_is_ewt)
{
ark_mem->rwt = NULL;
if (!arkAllocVec(ark_mem, ark_mem->ewt, &(ark_mem->rwt)))
{
arkProcessError(ark_mem, ARK_MEM_FAIL, __LINE__, __func__, __FILE__,
MSG_ARK_ARKMEM_FAIL);
return (ARK_ILL_INPUT);
}
ark_mem->rwt_is_ewt = SUNFALSE;
}
/* Copy tolerances into memory */
ark_mem->SRabstol = rabstol;
ark_mem->ritol = ARK_SS;
/* enforce use of arkRwtSet */
ark_mem->user_efun = SUNFALSE;
ark_mem->rfun = arkRwtSet;
ark_mem->r_data = ark_mem;
return (ARK_SUCCESS);
}
int arkResVtolerance(ARKodeMem ark_mem, N_Vector rabstol)
{
/* local variables */
sunrealtype rabstolmin;
/* Check inputs */
if (ark_mem == NULL)
{
arkProcessError(NULL, ARK_MEM_NULL, __LINE__, __func__, __FILE__,
MSG_ARK_NO_MEM);
return (ARK_MEM_NULL);
}
if (ark_mem->MallocDone == SUNFALSE)
{
arkProcessError(ark_mem, ARK_NO_MALLOC, __LINE__, __func__, __FILE__,
MSG_ARK_NO_MALLOC);
return (ARK_NO_MALLOC);
}
if (rabstol == NULL)
{
arkProcessError(ark_mem, ARK_NO_MALLOC, __LINE__, __func__, __FILE__,
MSG_ARK_NULL_RABSTOL);
return (ARK_NO_MALLOC);
}
if (rabstol->ops->nvmin == NULL)
{
arkProcessError(ark_mem, ARK_ILL_INPUT, __LINE__, __func__, __FILE__,
"Missing N_VMin routine from N_Vector");
return (ARK_ILL_INPUT);
}
rabstolmin = N_VMin(rabstol);
if (rabstolmin < ZERO)
{
arkProcessError(ark_mem, ARK_ILL_INPUT, __LINE__, __func__, __FILE__,
MSG_ARK_BAD_RABSTOL);
return (ARK_ILL_INPUT);
}
/* Set flag indicating whether min(abstol) == 0 */
ark_mem->Ratolmin0 = (rabstolmin == ZERO);
/* Allocate space for rwt if necessary */
if (ark_mem->rwt_is_ewt)
{
ark_mem->rwt = NULL;
if (!arkAllocVec(ark_mem, ark_mem->ewt, &(ark_mem->rwt)))
{
arkProcessError(ark_mem, ARK_MEM_FAIL, __LINE__, __func__, __FILE__,
MSG_ARK_ARKMEM_FAIL);
return (ARK_ILL_INPUT);
}
ark_mem->rwt_is_ewt = SUNFALSE;
}
/* Copy tolerances into memory */
if (!(ark_mem->VRabstolMallocDone))
{
if (!arkAllocVec(ark_mem, ark_mem->rwt, &(ark_mem->VRabstol)))
{
arkProcessError(ark_mem, ARK_MEM_FAIL, __LINE__, __func__, __FILE__,
MSG_ARK_ARKMEM_FAIL);
return (ARK_ILL_INPUT);
}
ark_mem->VRabstolMallocDone = SUNTRUE;
}
N_VScale(ONE, rabstol, ark_mem->VRabstol);
ark_mem->ritol = ARK_SV;
/* enforce use of arkRwtSet */
ark_mem->user_efun = SUNFALSE;
ark_mem->rfun = arkRwtSet;
ark_mem->r_data = ark_mem;
return (ARK_SUCCESS);
}
int arkResFtolerance(ARKodeMem ark_mem, ARKRwtFn rfun)
{
if (ark_mem == NULL)
{
arkProcessError(NULL, ARK_MEM_NULL, __LINE__, __func__, __FILE__,
MSG_ARK_NO_MEM);
return (ARK_MEM_NULL);
}
if (ark_mem->MallocDone == SUNFALSE)
{
arkProcessError(ark_mem, ARK_NO_MALLOC, __LINE__, __func__, __FILE__,
MSG_ARK_NO_MALLOC);
return (ARK_NO_MALLOC);
}
/* Allocate space for rwt if necessary */
if (ark_mem->rwt_is_ewt)
{
ark_mem->rwt = NULL;
if (!arkAllocVec(ark_mem, ark_mem->ewt, &(ark_mem->rwt)))
{
arkProcessError(ark_mem, ARK_MEM_FAIL, __LINE__, __func__, __FILE__,
MSG_ARK_ARKMEM_FAIL);
return (ARK_ILL_INPUT);
}
ark_mem->rwt_is_ewt = SUNFALSE;
}
/* Copy tolerance data into memory */
ark_mem->ritol = ARK_WF;
ark_mem->user_rfun = SUNTRUE;
ark_mem->rfun = rfun;
ark_mem->r_data = ark_mem->user_data;
return (ARK_SUCCESS);
}
/*---------------------------------------------------------------
arkEvolve:
This routine is the main driver of ARKODE-based integrators.
It integrates over a time interval defined by the user, by
calling the time step module to do internal time steps.
The first time that arkEvolve is called for a successfully
initialized problem, it computes a tentative initial step size.
arkEvolve supports two modes as specified by itask: ARK_NORMAL and
ARK_ONE_STEP. In the ARK_NORMAL mode, the solver steps until
it reaches or passes tout and then interpolates to obtain
y(tout). In the ARK_ONE_STEP mode, it takes one internal step
and returns. The behavior of both modes can be over-rided
through user-specification of ark_tstop (through the
*StepSetStopTime function), in which case if a solver step
would pass tstop, the step is shortened so that it stops at
exactly the specified stop time, and hence interpolation of
y(tout) is not required.
---------------------------------------------------------------*/
int arkEvolve(ARKodeMem ark_mem, sunrealtype tout, N_Vector yout,
sunrealtype* tret, int itask)
{
long int nstloc;
int retval, kflag, istate, ir;
int ewtsetOK;
sunrealtype troundoff, nrm;
sunbooleantype inactive_roots;
sunrealtype dsm;
int nflag, attempts, ncf, nef, constrfails;
int relax_fails;
/* Check and process inputs */
/* Check if ark_mem exists */
if (ark_mem == NULL)
{
arkProcessError(NULL, ARK_MEM_NULL, __LINE__, __func__, __FILE__,
MSG_ARK_NO_MEM);
return (ARK_MEM_NULL);
}
/* Check if ark_mem was allocated */
if (ark_mem->MallocDone == SUNFALSE)
{
arkProcessError(ark_mem, ARK_NO_MALLOC, __LINE__, __func__, __FILE__,
MSG_ARK_NO_MALLOC);
return (ARK_NO_MALLOC);
}
/* Check for yout != NULL */
if ((ark_mem->ycur = yout) == NULL)
{
arkProcessError(ark_mem, ARK_ILL_INPUT, __LINE__, __func__, __FILE__,
MSG_ARK_YOUT_NULL);
return (ARK_ILL_INPUT);
}
/* Check for tret != NULL */
if (tret == NULL)
{
arkProcessError(ark_mem, ARK_ILL_INPUT, __LINE__, __func__, __FILE__,
MSG_ARK_TRET_NULL);
return (ARK_ILL_INPUT);
}
/* Check for valid itask */
if ((itask != ARK_NORMAL) && (itask != ARK_ONE_STEP))
{
arkProcessError(ark_mem, ARK_ILL_INPUT, __LINE__, __func__, __FILE__,
MSG_ARK_BAD_ITASK);
return (ARK_ILL_INPUT);
}
/* store copy of itask if using root-finding */
if (ark_mem->root_mem != NULL)
{
if (itask == ARK_NORMAL) { ark_mem->root_mem->toutc = tout; }
ark_mem->root_mem->taskc = itask;
}
/* perform first-step-specific initializations:
- initialize tret values to initialization time
- perform initial integrator setup */
if (ark_mem->initsetup)
{
ark_mem->tretlast = *tret = ark_mem->tcur;
retval = arkInitialSetup(ark_mem, tout);
if (retval != ARK_SUCCESS) { return (retval); }
}
/* perform stopping tests */
if (!ark_mem->initsetup)
{
if (arkStopTests(ark_mem, tout, yout, tret, itask, &retval))
{
return (retval);
}
}
/*--------------------------------------------------
Looping point for successful internal steps
- update the ewt/rwt vectors for upcoming step
- check for errors (too many steps, too much
accuracy requested, step size too small)
- loop over attempts at a new step:
* try to take step (via time stepper module),
handle solver convergence or other failures
* perform constraint-handling (if selected)
* check temporal error
* if all of the above pass, complete step by
updating current time, solution, error &
stepsize history arrays.
- perform stop tests:
* check for root in last step taken
* check if tout was passed
* check if close to tstop
* check if in ONE_STEP mode (must return)
--------------------------------------------------*/
nstloc = 0;
for (;;)
{
ark_mem->next_h = ark_mem->h;
/* Reset and check ewt and rwt */
if (!ark_mem->initsetup)
{
ewtsetOK = ark_mem->efun(ark_mem->yn, ark_mem->ewt, ark_mem->e_data);
if (ewtsetOK != 0)
{
if (ark_mem->itol == ARK_WF)
{
arkProcessError(ark_mem, ARK_ILL_INPUT, __LINE__, __func__, __FILE__,
MSG_ARK_EWT_NOW_FAIL, ark_mem->tcur);
}
else
{
arkProcessError(ark_mem, ARK_ILL_INPUT, __LINE__, __func__, __FILE__,
MSG_ARK_EWT_NOW_BAD, ark_mem->tcur);
}
istate = ARK_ILL_INPUT;
ark_mem->tretlast = *tret = ark_mem->tcur;
N_VScale(ONE, ark_mem->yn, yout);
break;
}
if (!ark_mem->rwt_is_ewt)
{
ewtsetOK = ark_mem->rfun(ark_mem->yn, ark_mem->rwt, ark_mem->r_data);
if (ewtsetOK != 0)
{
if (ark_mem->itol == ARK_WF)
{
arkProcessError(ark_mem, ARK_ILL_INPUT, __LINE__, __func__,
__FILE__, MSG_ARK_RWT_NOW_FAIL, ark_mem->tcur);
}
else
{
arkProcessError(ark_mem, ARK_ILL_INPUT, __LINE__, __func__,
__FILE__, MSG_ARK_RWT_NOW_BAD, ark_mem->tcur);
}
istate = ARK_ILL_INPUT;
ark_mem->tretlast = *tret = ark_mem->tcur;
N_VScale(ONE, ark_mem->yn, yout);
break;
}
}
}
/* Check for too many steps */
if ((ark_mem->mxstep > 0) && (nstloc >= ark_mem->mxstep))
{
arkProcessError(ark_mem, ARK_TOO_MUCH_WORK, __LINE__, __func__, __FILE__,
MSG_ARK_MAX_STEPS, ark_mem->tcur);
istate = ARK_TOO_MUCH_WORK;
ark_mem->tretlast = *tret = ark_mem->tcur;
N_VScale(ONE, ark_mem->yn, yout);
break;
}
/* Check for too much accuracy requested */
nrm = N_VWrmsNorm(ark_mem->yn, ark_mem->ewt);
ark_mem->tolsf = ark_mem->uround * nrm;
if (ark_mem->tolsf > ONE)
{
arkProcessError(ark_mem, ARK_TOO_MUCH_ACC, __LINE__, __func__, __FILE__,
MSG_ARK_TOO_MUCH_ACC, ark_mem->tcur);
istate = ARK_TOO_MUCH_ACC;
ark_mem->tretlast = *tret = ark_mem->tcur;
N_VScale(ONE, ark_mem->yn, yout);
ark_mem->tolsf *= TWO;
break;
}
else { ark_mem->tolsf = ONE; }
/* Check for h below roundoff level in tn */
if (ark_mem->tcur + ark_mem->h == ark_mem->tcur)
{
ark_mem->nhnil++;
if (ark_mem->nhnil <= ark_mem->mxhnil)
{
arkProcessError(ark_mem, ARK_WARNING, __LINE__, __func__, __FILE__,
MSG_ARK_HNIL, ark_mem->tcur, ark_mem->h);
}
if (ark_mem->nhnil == ark_mem->mxhnil)
{
arkProcessError(ark_mem, ARK_WARNING, __LINE__, __func__, __FILE__,
MSG_ARK_HNIL_DONE);
}
}
/* Update parameter for upcoming step size */
if (ark_mem->hprime != ark_mem->h)
{
ark_mem->h = ark_mem->h * ark_mem->eta;
ark_mem->next_h = ark_mem->h;
}
if (ark_mem->fixedstep)
{
ark_mem->h = ark_mem->hin;
ark_mem->next_h = ark_mem->h;
/* patch for 'fixedstep' + 'tstop' use case:
limit fixed step size if step would overtake tstop */
if (ark_mem->tstopset)
{
if ((ark_mem->tcur + ark_mem->h - ark_mem->tstop) * ark_mem->h > ZERO)
{
ark_mem->h = (ark_mem->tstop - ark_mem->tcur) *
(ONE - FOUR * ark_mem->uround);
}
}
}
/* Looping point for step attempts */
dsm = ZERO;
attempts = ncf = nef = constrfails = ark_mem->last_kflag = 0;
relax_fails = 0;
nflag = FIRST_CALL;
for (;;)
{
/* increment attempt counters */
attempts++;
ark_mem->nst_attempts++;
#if SUNDIALS_LOGGING_LEVEL >= SUNDIALS_LOGGING_DEBUG
SUNLogger_QueueMsg(ARK_LOGGER, SUN_LOGLEVEL_DEBUG, "ARKODE::arkEvolve",
"start-step",
"step = %li, attempt = %i, h = %" RSYM
", tcur = %" RSYM,
ark_mem->nst, attempts, ark_mem->h, ark_mem->tcur);
#endif
/* Call time stepper module to attempt a step:
0 => step completed successfully
>0 => step encountered recoverable failure; reduce step if possible
<0 => step encountered unrecoverable failure */
kflag = ark_mem->step((void*)ark_mem, &dsm, &nflag);
if (kflag < 0) { break; }
/* handle solver convergence failures */
kflag = arkCheckConvergence(ark_mem, &nflag, &ncf);
if (kflag < 0) { break; }
/* Perform relaxation:
- computes relaxation parameter
- on success, updates ycur, h, and dsm
- on recoverable failure, updates eta and signals to retry step
- on fatal error, returns negative error flag */
if (ark_mem->relax_enabled && (kflag == ARK_SUCCESS))
{
kflag = arkRelax(ark_mem, &relax_fails, &dsm, &nflag);
if (kflag < 0) { break; }
}
/* perform constraint-handling (if selected, and if solver check passed) */
if (ark_mem->constraintsSet && (kflag == ARK_SUCCESS))
{
kflag = arkCheckConstraints(ark_mem, &constrfails, &nflag);
if (kflag < 0) { break; }
}
/* when fixed time-stepping is enabled, 'success' == successful stage solves
(checked in previous block), so just enforce no step size change */
if (ark_mem->fixedstep)
{
ark_mem->eta = ONE;
break;
}
/* check temporal error (if checks above passed) */
if (kflag == ARK_SUCCESS)
{
kflag = arkCheckTemporalError(ark_mem, &nflag, &nef, dsm);
if (kflag < 0) { break; }
}
/* if ignoring temporal error test result (XBraid) force step to pass */
if (ark_mem->force_pass)
{
ark_mem->last_kflag = kflag;
kflag = ARK_SUCCESS;
break;
}
/* break attempt loop on successful step */
if (kflag == ARK_SUCCESS) { break; }
/* unsuccessful step, if |h| = hmin, return ARK_ERR_FAILURE */
if (SUNRabs(ark_mem->h) <= ark_mem->hmin * ONEPSM)
{
return (ARK_ERR_FAILURE);
}
/* update h, hprime and next_h for next iteration */
ark_mem->h *= ark_mem->eta;
ark_mem->next_h = ark_mem->hprime = ark_mem->h;
} /* end looping for step attempts */
/* If step attempt loop succeeded, complete step (update current time, solution,
error stepsize history arrays; call user-supplied step postprocessing function)
(added stuff from arkStep_PrepareNextStep -- revisit) */
if (kflag == ARK_SUCCESS) { kflag = arkCompleteStep(ark_mem, dsm); }