xref: /petsc/src/ts/impls/arkimex/arkimex.c (revision c3b5f7ba6bc5ce25a01a67bb37ba5d34b02bbbd7)
1 /*
2   Code for timestepping with additive Runge-Kutta IMEX method
3 
4   Notes:
5   The general system is written as
6 
7   F(t,U,Udot) = G(t,U)
8 
9   where F represents the stiff part of the physics and G represents the non-stiff part.
10 
11 */
12 #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
13 #include <petscdm.h>
14 
15 static TSARKIMEXType  TSARKIMEXDefault = TSARKIMEX3;
16 static PetscBool      TSARKIMEXRegisterAllCalled;
17 static PetscBool      TSARKIMEXPackageInitialized;
18 static PetscErrorCode TSExtrapolate_ARKIMEX(TS,PetscReal,Vec);
19 
20 typedef struct _ARKTableau *ARKTableau;
21 struct _ARKTableau {
22   char      *name;
23   PetscInt  order;                /* Classical approximation order of the method */
24   PetscInt  s;                    /* Number of stages */
25   PetscBool stiffly_accurate;     /* The implicit part is stiffly accurate*/
26   PetscBool FSAL_implicit;        /* The implicit part is FSAL*/
27   PetscBool explicit_first_stage; /* The implicit part has an explicit first stage*/
28   PetscInt  pinterp;              /* Interpolation order */
29   PetscReal *At,*bt,*ct;          /* Stiff tableau */
30   PetscReal *A,*b,*c;             /* Non-stiff tableau */
31   PetscReal *bembedt,*bembed;     /* Embedded formula of order one less (order-1) */
32   PetscReal *binterpt,*binterp;   /* Dense output formula */
33   PetscReal ccfl;                 /* Placeholder for CFL coefficient relative to forward Euler */
34 };
35 typedef struct _ARKTableauLink *ARKTableauLink;
36 struct _ARKTableauLink {
37   struct _ARKTableau tab;
38   ARKTableauLink     next;
39 };
40 static ARKTableauLink ARKTableauList;
41 
42 typedef struct {
43   ARKTableau   tableau;
44   Vec          *Y;               /* States computed during the step */
45   Vec          *YdotI;           /* Time derivatives for the stiff part */
46   Vec          *YdotRHS;         /* Function evaluations for the non-stiff part */
47   Vec          *Y_prev;          /* States computed during the previous time step */
48   Vec          *YdotI_prev;      /* Time derivatives for the stiff part for the previous time step*/
49   Vec          *YdotRHS_prev;    /* Function evaluations for the non-stiff part for the previous time step*/
50   Vec          Ydot0;            /* Holds the slope from the previous step in FSAL case */
51   Vec          Ydot;             /* Work vector holding Ydot during residual evaluation */
52   Vec          Z;                /* Ydot = shift(Y-Z) */
53   PetscScalar  *work;            /* Scalar work */
54   PetscReal    scoeff;           /* shift = scoeff/dt */
55   PetscReal    stage_time;
56   PetscBool    imex;
57   PetscBool    extrapolate;      /* Extrapolate initial guess from previous time-step stage values */
58   TSStepStatus status;
59 } TS_ARKIMEX;
60 /*MC
61      TSARKIMEXARS122 - Second order ARK IMEX scheme.
62 
63      This method has one explicit stage and one implicit stage.
64 
65      Options Database:
66 .      -ts_arkimex_type ars122 - set arkimex type to ars122
67 
68      References:
69 .    * -  U. Ascher, S. Ruuth, R. J. Spiteri, Implicit explicit Runge Kutta methods for time dependent Partial Differential Equations. Appl. Numer. Math. 25, (1997).
70 
71      Level: advanced
72 
73 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
74 M*/
75 /*MC
76      TSARKIMEXA2 - Second order ARK IMEX scheme with A-stable implicit part.
77 
78      This method has an explicit stage and one implicit stage, and has an A-stable implicit scheme. This method was provided by Emil Constantinescu.
79 
80      Options Database:
81 .      -ts_arkimex_type a2 - set arkimex type to a2
82 
83      Level: advanced
84 
85 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
86 M*/
87 /*MC
88      TSARKIMEXL2 - Second order ARK IMEX scheme with L-stable implicit part.
89 
90      This method has two implicit stages, and L-stable implicit scheme.
91 
92      Options Database:
93 .      -ts_arkimex_type l2 - set arkimex type to l2
94 
95     References:
96 .   * -  L. Pareschi, G. Russo, Implicit Explicit Runge Kutta schemes and applications to hyperbolic systems with relaxations. Journal of Scientific Computing Volume: 25, Issue: 1, October, 2005.
97 
98      Level: advanced
99 
100 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
101 M*/
102 /*MC
103      TSARKIMEX1BEE - First order backward Euler represented as an ARK IMEX scheme with extrapolation as error estimator. This is a 3-stage method.
104 
105      This method is aimed at starting the integration of implicit DAEs when explicit first-stage ARK methods are used.
106 
107      Options Database:
108 .      -ts_arkimex_type 1bee - set arkimex type to 1bee
109 
110      Level: advanced
111 
112 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
113 M*/
114 /*MC
115      TSARKIMEX2C - Second order ARK IMEX scheme with L-stable implicit part.
116 
117      This method has one explicit stage and two implicit stages. The implicit part is the same as in TSARKIMEX2D and TSARKIMEX2E, but the explicit part has a larger stability region on the negative real axis. This method was provided by Emil Constantinescu.
118 
119      Options Database:
120 .      -ts_arkimex_type 2c - set arkimex type to 2c
121 
122      Level: advanced
123 
124 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
125 M*/
126 /*MC
127      TSARKIMEX2D - Second order ARK IMEX scheme with L-stable implicit part.
128 
129      This method has one explicit stage and two implicit stages. The stability function is independent of the explicit part in the infinity limit of the implict component. This method was provided by Emil Constantinescu.
130 
131      Options Database:
132 .      -ts_arkimex_type 2d - set arkimex type to 2d
133 
134      Level: advanced
135 
136 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
137 M*/
138 /*MC
139      TSARKIMEX2E - Second order ARK IMEX scheme with L-stable implicit part.
140 
141      This method has one explicit stage and two implicit stages. It is is an optimal method developed by Emil Constantinescu.
142 
143      Options Database:
144 .      -ts_arkimex_type 2e - set arkimex type to 2e
145 
146     Level: advanced
147 
148 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
149 M*/
150 /*MC
151      TSARKIMEXPRSSP2 - Second order SSP ARK IMEX scheme.
152 
153      This method has three implicit stages.
154 
155      References:
156 .    * -  L. Pareschi, G. Russo, Implicit Explicit Runge Kutta schemes and applications to hyperbolic systems with relaxations. Journal of Scientific Computing Volume: 25, Issue: 1, October, 2005.
157 
158      This method is referred to as SSP2-(3,3,2) in https://arxiv.org/abs/1110.4375
159 
160      Options Database:
161 .      -ts_arkimex_type prssp2 - set arkimex type to prssp2
162 
163      Level: advanced
164 
165 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
166 M*/
167 /*MC
168      TSARKIMEX3 - Third order ARK IMEX scheme with L-stable implicit part.
169 
170      This method has one explicit stage and three implicit stages.
171 
172      Options Database:
173 .      -ts_arkimex_type 3 - set arkimex type to 3
174 
175      References:
176 .    * -  Kennedy and Carpenter 2003.
177 
178      Level: advanced
179 
180 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
181 M*/
182 /*MC
183      TSARKIMEXARS443 - Third order ARK IMEX scheme.
184 
185      This method has one explicit stage and four implicit stages.
186 
187      Options Database:
188 .      -ts_arkimex_type ars443 - set arkimex type to ars443
189 
190      References:
191 +    * -  U. Ascher, S. Ruuth, R. J. Spiteri, Implicit explicit Runge Kutta methods for time dependent Partial Differential Equations. Appl. Numer. Math. 25, (1997).
192 -    * -  This method is referred to as ARS(4,4,3) in https://arxiv.org/abs/1110.4375
193 
194      Level: advanced
195 
196 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
197 M*/
198 /*MC
199      TSARKIMEXBPR3 - Third order ARK IMEX scheme.
200 
201      This method has one explicit stage and four implicit stages.
202 
203      Options Database:
204 .      -ts_arkimex_type bpr3 - set arkimex type to bpr3
205 
206      References:
207 .    * - This method is referred to as ARK3 in https://arxiv.org/abs/1110.4375
208 
209      Level: advanced
210 
211 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
212 M*/
213 /*MC
214      TSARKIMEX4 - Fourth order ARK IMEX scheme with L-stable implicit part.
215 
216      This method has one explicit stage and four implicit stages.
217 
218      Options Database:
219 .      -ts_arkimex_type 4 - set arkimex type to4
220 
221      References:
222 .    * - Kennedy and Carpenter 2003.
223 
224      Level: advanced
225 
226 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
227 M*/
228 /*MC
229      TSARKIMEX5 - Fifth order ARK IMEX scheme with L-stable implicit part.
230 
231      This method has one explicit stage and five implicit stages.
232 
233      Options Database:
234 .      -ts_arkimex_type 5 - set arkimex type to 5
235 
236      References:
237 .    * - Kennedy and Carpenter 2003.
238 
239      Level: advanced
240 
241 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
242 M*/
243 
244 /*@C
245   TSARKIMEXRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSARKIMEX
246 
247   Not Collective, but should be called by all processes which will need the schemes to be registered
248 
249   Level: advanced
250 
251 .seealso:  TSARKIMEXRegisterDestroy()
252 @*/
253 PetscErrorCode TSARKIMEXRegisterAll(void)
254 {
255   PetscFunctionBegin;
256   if (TSARKIMEXRegisterAllCalled) PetscFunctionReturn(0);
257   TSARKIMEXRegisterAllCalled = PETSC_TRUE;
258 
259   {
260     const PetscReal
261       A[3][3] = {{0.0,0.0,0.0},
262                  {0.0,0.0,0.0},
263                  {0.0,0.5,0.0}},
264       At[3][3] = {{1.0,0.0,0.0},
265                   {0.0,0.5,0.0},
266                   {0.0,0.5,0.5}},
267       b[3]       = {0.0,0.5,0.5},
268       bembedt[3] = {1.0,0.0,0.0};
269     PetscCall(TSARKIMEXRegister(TSARKIMEX1BEE,2,3,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL));
270   }
271   {
272     const PetscReal
273       A[2][2] = {{0.0,0.0},
274                  {0.5,0.0}},
275       At[2][2] = {{0.0,0.0},
276                   {0.0,0.5}},
277       b[2]       = {0.0,1.0},
278       bembedt[2] = {0.5,0.5};
279     /* binterpt[2][2] = {{1.0,-1.0},{0.0,1.0}};  second order dense output has poor stability properties and hence it is not currently in use*/
280     PetscCall(TSARKIMEXRegister(TSARKIMEXARS122,2,2,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL));
281   }
282   {
283     const PetscReal
284       A[2][2] = {{0.0,0.0},
285                  {1.0,0.0}},
286       At[2][2] = {{0.0,0.0},
287                   {0.5,0.5}},
288       b[2]       = {0.5,0.5},
289       bembedt[2] = {0.0,1.0};
290     /* binterpt[2][2] = {{1.0,-0.5},{0.0,0.5}}  second order dense output has poor stability properties and hence it is not currently in use*/
291     PetscCall(TSARKIMEXRegister(TSARKIMEXA2,2,2,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL));
292   }
293   {
294     /* const PetscReal us2 = 1.0-1.0/PetscSqrtReal((PetscReal)2.0);    Direct evaluation: 0.2928932188134524755992. Used below to ensure all values are available at compile time   */
295     const PetscReal
296       A[2][2] = {{0.0,0.0},
297                  {1.0,0.0}},
298       At[2][2] = {{0.2928932188134524755992,0.0},
299                   {1.0-2.0*0.2928932188134524755992,0.2928932188134524755992}},
300       b[2]       = {0.5,0.5},
301       bembedt[2] = {0.0,1.0},
302       binterpt[2][2] = {{  (0.2928932188134524755992-1.0)/(2.0*0.2928932188134524755992-1.0),-1/(2.0*(1.0-2.0*0.2928932188134524755992))},
303                         {1-(0.2928932188134524755992-1.0)/(2.0*0.2928932188134524755992-1.0),-1/(2.0*(1.0-2.0*0.2928932188134524755992))}},
304       binterp[2][2] = {{1.0,-0.5},{0.0,0.5}};
305     PetscCall(TSARKIMEXRegister(TSARKIMEXL2,2,2,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,2,binterpt[0],binterp[0]));
306   }
307   {
308     /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0),  Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time   */
309     const PetscReal
310       A[3][3] = {{0,0,0},
311                  {2-1.414213562373095048802,0,0},
312                  {0.5,0.5,0}},
313       At[3][3] = {{0,0,0},
314                   {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0},
315                   {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}},
316       bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)},
317       binterpt[3][2] = {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
318                         {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
319                         {1.0-1.414213562373095048802,1.0/1.414213562373095048802}};
320     PetscCall(TSARKIMEXRegister(TSARKIMEX2C,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL));
321   }
322   {
323     /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0),  Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time   */
324     const PetscReal
325       A[3][3] = {{0,0,0},
326                  {2-1.414213562373095048802,0,0},
327                  {0.75,0.25,0}},
328       At[3][3] = {{0,0,0},
329                   {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0},
330                   {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}},
331       bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)},
332       binterpt[3][2] =  {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
333                          {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
334                          {1.0-1.414213562373095048802,1.0/1.414213562373095048802}};
335     PetscCall(TSARKIMEXRegister(TSARKIMEX2D,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL));
336   }
337   {                             /* Optimal for linear implicit part */
338     /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0),  Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time   */
339     const PetscReal
340       A[3][3] = {{0,0,0},
341                  {2-1.414213562373095048802,0,0},
342                  {(3-2*1.414213562373095048802)/6,(3+2*1.414213562373095048802)/6,0}},
343       At[3][3] = {{0,0,0},
344                   {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0},
345                   {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}},
346       bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)},
347       binterpt[3][2] =  {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
348                          {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
349                          {1.0-1.414213562373095048802,1.0/1.414213562373095048802}};
350     PetscCall(TSARKIMEXRegister(TSARKIMEX2E,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL));
351   }
352   {                             /* Optimal for linear implicit part */
353     const PetscReal
354       A[3][3] = {{0,0,0},
355                  {0.5,0,0},
356                  {0.5,0.5,0}},
357       At[3][3] = {{0.25,0,0},
358                   {0,0.25,0},
359                   {1./3,1./3,1./3}};
360     PetscCall(TSARKIMEXRegister(TSARKIMEXPRSSP2,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,NULL,NULL,0,NULL,NULL));
361   }
362   {
363     const PetscReal
364       A[4][4] = {{0,0,0,0},
365                  {1767732205903./2027836641118.,0,0,0},
366                  {5535828885825./10492691773637.,788022342437./10882634858940.,0,0},
367                  {6485989280629./16251701735622.,-4246266847089./9704473918619.,10755448449292./10357097424841.,0}},
368       At[4][4] = {{0,0,0,0},
369                   {1767732205903./4055673282236.,1767732205903./4055673282236.,0,0},
370                   {2746238789719./10658868560708.,-640167445237./6845629431997.,1767732205903./4055673282236.,0},
371                   {1471266399579./7840856788654.,-4482444167858./7529755066697.,11266239266428./11593286722821.,1767732205903./4055673282236.}},
372       bembedt[4]     = {2756255671327./12835298489170.,-10771552573575./22201958757719.,9247589265047./10645013368117.,2193209047091./5459859503100.},
373       binterpt[4][2] = {{4655552711362./22874653954995., -215264564351./13552729205753.},
374                         {-18682724506714./9892148508045.,17870216137069./13817060693119.},
375                         {34259539580243./13192909600954.,-28141676662227./17317692491321.},
376                         {584795268549./6622622206610.,   2508943948391./7218656332882.}};
377     PetscCall(TSARKIMEXRegister(TSARKIMEX3,3,4,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL));
378   }
379   {
380     const PetscReal
381       A[5][5] = {{0,0,0,0,0},
382                  {1./2,0,0,0,0},
383                  {11./18,1./18,0,0,0},
384                  {5./6,-5./6,.5,0,0},
385                  {1./4,7./4,3./4,-7./4,0}},
386       At[5][5] = {{0,0,0,0,0},
387                   {0,1./2,0,0,0},
388                   {0,1./6,1./2,0,0},
389                   {0,-1./2,1./2,1./2,0},
390                   {0,3./2,-3./2,1./2,1./2}},
391     *bembedt = NULL;
392     PetscCall(TSARKIMEXRegister(TSARKIMEXARS443,3,5,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,0,NULL,NULL));
393   }
394   {
395     const PetscReal
396       A[5][5] = {{0,0,0,0,0},
397                  {1,0,0,0,0},
398                  {4./9,2./9,0,0,0},
399                  {1./4,0,3./4,0,0},
400                  {1./4,0,3./5,0,0}},
401       At[5][5] = {{0,0,0,0,0},
402                   {.5,.5,0,0,0},
403                   {5./18,-1./9,.5,0,0},
404                   {.5,0,0,.5,0},
405                   {.25,0,.75,-.5,.5}},
406     *bembedt = NULL;
407     PetscCall(TSARKIMEXRegister(TSARKIMEXBPR3,3,5,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,0,NULL,NULL));
408   }
409   {
410     const PetscReal
411       A[6][6] = {{0,0,0,0,0,0},
412                  {1./2,0,0,0,0,0},
413                  {13861./62500.,6889./62500.,0,0,0,0},
414                  {-116923316275./2393684061468.,-2731218467317./15368042101831.,9408046702089./11113171139209.,0,0,0},
415                  {-451086348788./2902428689909.,-2682348792572./7519795681897.,12662868775082./11960479115383.,3355817975965./11060851509271.,0,0},
416                  {647845179188./3216320057751.,73281519250./8382639484533.,552539513391./3454668386233.,3354512671639./8306763924573.,4040./17871.,0}},
417       At[6][6] = {{0,0,0,0,0,0},
418                   {1./4,1./4,0,0,0,0},
419                   {8611./62500.,-1743./31250.,1./4,0,0,0},
420                   {5012029./34652500.,-654441./2922500.,174375./388108.,1./4,0,0},
421                   {15267082809./155376265600.,-71443401./120774400.,730878875./902184768.,2285395./8070912.,1./4,0},
422                   {82889./524892.,0,15625./83664.,69875./102672.,-2260./8211,1./4}},
423       bembedt[6]     = {4586570599./29645900160.,0,178811875./945068544.,814220225./1159782912.,-3700637./11593932.,61727./225920.},
424       binterpt[6][3] = {{6943876665148./7220017795957.,-54480133./30881146.,6818779379841./7100303317025.},
425                         {0,0,0},
426                         {7640104374378./9702883013639.,-11436875./14766696.,2173542590792./12501825683035.},
427                         {-20649996744609./7521556579894.,174696575./18121608.,-31592104683404./5083833661969.},
428                         {8854892464581./2390941311638.,-12120380./966161.,61146701046299./7138195549469.},
429                         {-11397109935349./6675773540249.,3843./706.,-17219254887155./4939391667607.}};
430     PetscCall(TSARKIMEXRegister(TSARKIMEX4,4,6,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,3,binterpt[0],NULL));
431   }
432   {
433     const PetscReal
434       A[8][8] = {{0,0,0,0,0,0,0,0},
435                  {41./100,0,0,0,0,0,0,0},
436                  {367902744464./2072280473677.,677623207551./8224143866563.,0,0,0,0,0,0},
437                  {1268023523408./10340822734521.,0,1029933939417./13636558850479.,0,0,0,0,0},
438                  {14463281900351./6315353703477.,0,66114435211212./5879490589093.,-54053170152839./4284798021562.,0,0,0,0},
439                  {14090043504691./34967701212078.,0,15191511035443./11219624916014.,-18461159152457./12425892160975.,-281667163811./9011619295870.,0,0,0},
440                  {19230459214898./13134317526959.,0,21275331358303./2942455364971.,-38145345988419./4862620318723.,-1./8,-1./8,0,0},
441                  {-19977161125411./11928030595625.,0,-40795976796054./6384907823539.,177454434618887./12078138498510.,782672205425./8267701900261.,-69563011059811./9646580694205.,7356628210526./4942186776405.,0}},
442       At[8][8] = {{0,0,0,0,0,0,0,0},
443                   {41./200.,41./200.,0,0,0,0,0,0},
444                   {41./400.,-567603406766./11931857230679.,41./200.,0,0,0,0,0},
445                   {683785636431./9252920307686.,0,-110385047103./1367015193373.,41./200.,0,0,0,0},
446                   {3016520224154./10081342136671.,0,30586259806659./12414158314087.,-22760509404356./11113319521817.,41./200.,0,0,0},
447                   {218866479029./1489978393911.,0,638256894668./5436446318841.,-1179710474555./5321154724896.,-60928119172./8023461067671.,41./200.,0,0},
448                   {1020004230633./5715676835656.,0,25762820946817./25263940353407.,-2161375909145./9755907335909.,-211217309593./5846859502534.,-4269925059573./7827059040749.,41./200,0},
449                   {-872700587467./9133579230613.,0,0,22348218063261./9555858737531.,-1143369518992./8141816002931.,-39379526789629./19018526304540.,32727382324388./42900044865799.,41./200.}},
450       bembedt[8]     = {-975461918565./9796059967033.,0,0,78070527104295./32432590147079.,-548382580838./3424219808633.,-33438840321285./15594753105479.,3629800801594./4656183773603.,4035322873751./18575991585200.},
451       binterpt[8][3] = {{-17674230611817./10670229744614.,  43486358583215./12773830924787., -9257016797708./5021505065439.},
452                         {0,  0, 0                            },
453                         {0,  0, 0                            },
454                         {65168852399939./7868540260826.,  -91478233927265./11067650958493., 26096422576131./11239449250142.},
455                         {15494834004392./5936557850923.,  -79368583304911./10890268929626., 92396832856987./20362823103730.},
456                         {-99329723586156./26959484932159.,  -12239297817655./9152339842473., 30029262896817./10175596800299.},
457                         {-19024464361622./5461577185407.,  115839755401235./10719374521269., -26136350496073./3983972220547.},
458                         {-6511271360970./6095937251113.,  5843115559534./2180450260947., -5289405421727./3760307252460. }};
459     PetscCall(TSARKIMEXRegister(TSARKIMEX5,5,8,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,3,binterpt[0],NULL));
460   }
461   PetscFunctionReturn(0);
462 }
463 
464 /*@C
465    TSARKIMEXRegisterDestroy - Frees the list of schemes that were registered by TSARKIMEXRegister().
466 
467    Not Collective
468 
469    Level: advanced
470 
471 .seealso: TSARKIMEXRegister(), TSARKIMEXRegisterAll()
472 @*/
473 PetscErrorCode TSARKIMEXRegisterDestroy(void)
474 {
475   ARKTableauLink link;
476 
477   PetscFunctionBegin;
478   while ((link = ARKTableauList)) {
479     ARKTableau t = &link->tab;
480     ARKTableauList = link->next;
481     PetscCall(PetscFree6(t->At,t->bt,t->ct,t->A,t->b,t->c));
482     PetscCall(PetscFree2(t->bembedt,t->bembed));
483     PetscCall(PetscFree2(t->binterpt,t->binterp));
484     PetscCall(PetscFree(t->name));
485     PetscCall(PetscFree(link));
486   }
487   TSARKIMEXRegisterAllCalled = PETSC_FALSE;
488   PetscFunctionReturn(0);
489 }
490 
491 /*@C
492   TSARKIMEXInitializePackage - This function initializes everything in the TSARKIMEX package. It is called
493   from TSInitializePackage().
494 
495   Level: developer
496 
497 .seealso: PetscInitialize()
498 @*/
499 PetscErrorCode TSARKIMEXInitializePackage(void)
500 {
501   PetscFunctionBegin;
502   if (TSARKIMEXPackageInitialized) PetscFunctionReturn(0);
503   TSARKIMEXPackageInitialized = PETSC_TRUE;
504   PetscCall(TSARKIMEXRegisterAll());
505   PetscCall(PetscRegisterFinalize(TSARKIMEXFinalizePackage));
506   PetscFunctionReturn(0);
507 }
508 
509 /*@C
510   TSARKIMEXFinalizePackage - This function destroys everything in the TSARKIMEX package. It is
511   called from PetscFinalize().
512 
513   Level: developer
514 
515 .seealso: PetscFinalize()
516 @*/
517 PetscErrorCode TSARKIMEXFinalizePackage(void)
518 {
519   PetscFunctionBegin;
520   TSARKIMEXPackageInitialized = PETSC_FALSE;
521   PetscCall(TSARKIMEXRegisterDestroy());
522   PetscFunctionReturn(0);
523 }
524 
525 /*@C
526    TSARKIMEXRegister - register an ARK IMEX scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
527 
528    Not Collective, but the same schemes should be registered on all processes on which they will be used
529 
530    Input Parameters:
531 +  name - identifier for method
532 .  order - approximation order of method
533 .  s - number of stages, this is the dimension of the matrices below
534 .  At - Butcher table of stage coefficients for stiff part (dimension s*s, row-major)
535 .  bt - Butcher table for completing the stiff part of the step (dimension s; NULL to use the last row of At)
536 .  ct - Abscissa of each stiff stage (dimension s, NULL to use row sums of At)
537 .  A - Non-stiff stage coefficients (dimension s*s, row-major)
538 .  b - Non-stiff step completion table (dimension s; NULL to use last row of At)
539 .  c - Non-stiff abscissa (dimension s; NULL to use row sums of A)
540 .  bembedt - Stiff part of completion table for embedded method (dimension s; NULL if not available)
541 .  bembed - Non-stiff part of completion table for embedded method (dimension s; NULL to use bembedt if provided)
542 .  pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt and binterp
543 .  binterpt - Coefficients of the interpolation formula for the stiff part (dimension s*pinterp)
544 -  binterp - Coefficients of the interpolation formula for the non-stiff part (dimension s*pinterp; NULL to reuse binterpt)
545 
546    Notes:
547    Several ARK IMEX methods are provided, this function is only needed to create new methods.
548 
549    Level: advanced
550 
551 .seealso: TSARKIMEX
552 @*/
553 PetscErrorCode TSARKIMEXRegister(TSARKIMEXType name,PetscInt order,PetscInt s,
554                                  const PetscReal At[],const PetscReal bt[],const PetscReal ct[],
555                                  const PetscReal A[],const PetscReal b[],const PetscReal c[],
556                                  const PetscReal bembedt[],const PetscReal bembed[],
557                                  PetscInt pinterp,const PetscReal binterpt[],const PetscReal binterp[])
558 {
559   ARKTableauLink link;
560   ARKTableau     t;
561   PetscInt       i,j;
562 
563   PetscFunctionBegin;
564   PetscCall(TSARKIMEXInitializePackage());
565   PetscCall(PetscNew(&link));
566   t        = &link->tab;
567   PetscCall(PetscStrallocpy(name,&t->name));
568   t->order = order;
569   t->s     = s;
570   PetscCall(PetscMalloc6(s*s,&t->At,s,&t->bt,s,&t->ct,s*s,&t->A,s,&t->b,s,&t->c));
571   PetscCall(PetscArraycpy(t->At,At,s*s));
572   PetscCall(PetscArraycpy(t->A,A,s*s));
573   if (bt) PetscCall(PetscArraycpy(t->bt,bt,s));
574   else for (i=0; i<s; i++) t->bt[i] = At[(s-1)*s+i];
575   if (b)  PetscCall(PetscArraycpy(t->b,b,s));
576   else for (i=0; i<s; i++) t->b[i] = t->bt[i];
577   if (ct) PetscCall(PetscArraycpy(t->ct,ct,s));
578   else for (i=0; i<s; i++) for (j=0,t->ct[i]=0; j<s; j++) t->ct[i] += At[i*s+j];
579   if (c)  PetscCall(PetscArraycpy(t->c,c,s));
580   else for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j];
581   t->stiffly_accurate = PETSC_TRUE;
582   for (i=0; i<s; i++) if (t->At[(s-1)*s+i] != t->bt[i]) t->stiffly_accurate = PETSC_FALSE;
583   t->explicit_first_stage = PETSC_TRUE;
584   for (i=0; i<s; i++) if (t->At[i] != 0.0) t->explicit_first_stage = PETSC_FALSE;
585   /*def of FSAL can be made more precise*/
586   t->FSAL_implicit = (PetscBool)(t->explicit_first_stage && t->stiffly_accurate);
587   if (bembedt) {
588     PetscCall(PetscMalloc2(s,&t->bembedt,s,&t->bembed));
589     PetscCall(PetscArraycpy(t->bembedt,bembedt,s));
590     PetscCall(PetscArraycpy(t->bembed,bembed ? bembed : bembedt,s));
591   }
592 
593   t->pinterp     = pinterp;
594   PetscCall(PetscMalloc2(s*pinterp,&t->binterpt,s*pinterp,&t->binterp));
595   PetscCall(PetscArraycpy(t->binterpt,binterpt,s*pinterp));
596   PetscCall(PetscArraycpy(t->binterp,binterp ? binterp : binterpt,s*pinterp));
597   link->next     = ARKTableauList;
598   ARKTableauList = link;
599   PetscFunctionReturn(0);
600 }
601 
602 /*
603  The step completion formula is
604 
605  x1 = x0 - h bt^T YdotI + h b^T YdotRHS
606 
607  This function can be called before or after ts->vec_sol has been updated.
608  Suppose we have a completion formula (bt,b) and an embedded formula (bet,be) of different order.
609  We can write
610 
611  x1e = x0 - h bet^T YdotI + h be^T YdotRHS
612      = x1 + h bt^T YdotI - h b^T YdotRHS - h bet^T YdotI + h be^T YdotRHS
613      = x1 - h (bet - bt)^T YdotI + h (be - b)^T YdotRHS
614 
615  so we can evaluate the method with different order even after the step has been optimistically completed.
616 */
617 static PetscErrorCode TSEvaluateStep_ARKIMEX(TS ts,PetscInt order,Vec X,PetscBool *done)
618 {
619   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
620   ARKTableau     tab  = ark->tableau;
621   PetscScalar    *w   = ark->work;
622   PetscReal      h;
623   PetscInt       s = tab->s,j;
624 
625   PetscFunctionBegin;
626   switch (ark->status) {
627   case TS_STEP_INCOMPLETE:
628   case TS_STEP_PENDING:
629     h = ts->time_step; break;
630   case TS_STEP_COMPLETE:
631     h = ts->ptime - ts->ptime_prev; break;
632   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
633   }
634   if (order == tab->order) {
635     if (ark->status == TS_STEP_INCOMPLETE) {
636       if (!ark->imex && tab->stiffly_accurate) { /* Only the stiffly accurate implicit formula is used */
637         PetscCall(VecCopy(ark->Y[s-1],X));
638       } else { /* Use the standard completion formula (bt,b) */
639         PetscCall(VecCopy(ts->vec_sol,X));
640         for (j=0; j<s; j++) w[j] = h*tab->bt[j];
641         PetscCall(VecMAXPY(X,s,w,ark->YdotI));
642         if (ark->imex) { /* Method is IMEX, complete the explicit formula */
643           for (j=0; j<s; j++) w[j] = h*tab->b[j];
644           PetscCall(VecMAXPY(X,s,w,ark->YdotRHS));
645         }
646       }
647     } else PetscCall(VecCopy(ts->vec_sol,X));
648     if (done) *done = PETSC_TRUE;
649     PetscFunctionReturn(0);
650   } else if (order == tab->order-1) {
651     if (!tab->bembedt) goto unavailable;
652     if (ark->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (bet,be) */
653       PetscCall(VecCopy(ts->vec_sol,X));
654       for (j=0; j<s; j++) w[j] = h*tab->bembedt[j];
655       PetscCall(VecMAXPY(X,s,w,ark->YdotI));
656       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
657       PetscCall(VecMAXPY(X,s,w,ark->YdotRHS));
658     } else { /* Rollback and re-complete using (bet-be,be-b) */
659       PetscCall(VecCopy(ts->vec_sol,X));
660       for (j=0; j<s; j++) w[j] = h*(tab->bembedt[j] - tab->bt[j]);
661       PetscCall(VecMAXPY(X,tab->s,w,ark->YdotI));
662       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
663       PetscCall(VecMAXPY(X,s,w,ark->YdotRHS));
664     }
665     if (done) *done = PETSC_TRUE;
666     PetscFunctionReturn(0);
667   }
668 unavailable:
669   if (done) *done = PETSC_FALSE;
670   else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"ARKIMEX '%s' of order %D cannot evaluate step at order %D. Consider using -ts_adapt_type none or a different method that has an embedded estimate.",tab->name,tab->order,order);
671   PetscFunctionReturn(0);
672 }
673 
674 static PetscErrorCode TSARKIMEXTestMassIdentity(TS ts,PetscBool *id)
675 {
676   Vec            Udot,Y1,Y2;
677   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
678   PetscReal      norm;
679 
680   PetscFunctionBegin;
681   PetscCall(VecDuplicate(ts->vec_sol,&Udot));
682   PetscCall(VecDuplicate(ts->vec_sol,&Y1));
683   PetscCall(VecDuplicate(ts->vec_sol,&Y2));
684   PetscCall(TSComputeIFunction(ts,ts->ptime,ts->vec_sol,Udot,Y1,ark->imex));
685   PetscCall(VecSetRandom(Udot,NULL));
686   PetscCall(TSComputeIFunction(ts,ts->ptime,ts->vec_sol,Udot,Y2,ark->imex));
687   PetscCall(VecAXPY(Y2,-1.0,Y1));
688   PetscCall(VecAXPY(Y2,-1.0,Udot));
689   PetscCall(VecNorm(Y2,NORM_2,&norm));
690   if (norm < 100.0*PETSC_MACHINE_EPSILON) {
691     *id = PETSC_TRUE;
692   } else {
693     *id = PETSC_FALSE;
694     PetscCall(PetscInfo((PetscObject)ts,"IFunction(Udot = random) - IFunction(Udot = 0) is not near Udot, %g, suspect mass matrix implied in IFunction() is not the identity as required\n",(double)norm));
695   }
696   PetscCall(VecDestroy(&Udot));
697   PetscCall(VecDestroy(&Y1));
698   PetscCall(VecDestroy(&Y2));
699   PetscFunctionReturn(0);
700 }
701 
702 static PetscErrorCode TSRollBack_ARKIMEX(TS ts)
703 {
704   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
705   ARKTableau      tab  = ark->tableau;
706   const PetscInt  s    = tab->s;
707   const PetscReal *bt  = tab->bt,*b = tab->b;
708   PetscScalar     *w   = ark->work;
709   Vec             *YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS;
710   PetscInt        j;
711   PetscReal       h;
712 
713   PetscFunctionBegin;
714   switch (ark->status) {
715   case TS_STEP_INCOMPLETE:
716   case TS_STEP_PENDING:
717     h = ts->time_step; break;
718   case TS_STEP_COMPLETE:
719     h = ts->ptime - ts->ptime_prev; break;
720   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
721   }
722   for (j=0; j<s; j++) w[j] = -h*bt[j];
723   PetscCall(VecMAXPY(ts->vec_sol,s,w,YdotI));
724   for (j=0; j<s; j++) w[j] = -h*b[j];
725   PetscCall(VecMAXPY(ts->vec_sol,s,w,YdotRHS));
726   PetscFunctionReturn(0);
727 }
728 
729 static PetscErrorCode TSStep_ARKIMEX(TS ts)
730 {
731   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
732   ARKTableau      tab  = ark->tableau;
733   const PetscInt  s    = tab->s;
734   const PetscReal *At  = tab->At,*A = tab->A,*ct = tab->ct,*c = tab->c;
735   PetscScalar     *w   = ark->work;
736   Vec             *Y   = ark->Y,*YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS,Ydot = ark->Ydot,Ydot0 = ark->Ydot0,Z = ark->Z;
737   PetscBool       extrapolate = ark->extrapolate;
738   TSAdapt         adapt;
739   SNES            snes;
740   PetscInt        i,j,its,lits;
741   PetscInt        rejections = 0;
742   PetscBool       stageok,accept = PETSC_TRUE;
743   PetscReal       next_time_step = ts->time_step;
744 
745   PetscFunctionBegin;
746   if (ark->extrapolate && !ark->Y_prev) {
747     PetscCall(VecDuplicateVecs(ts->vec_sol,tab->s,&ark->Y_prev));
748     PetscCall(VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotI_prev));
749     PetscCall(VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotRHS_prev));
750   }
751 
752   if (!ts->steprollback) {
753     if (ts->equation_type >= TS_EQ_IMPLICIT) { /* Save the initial slope for the next step */
754       PetscCall(VecCopy(YdotI[s-1],Ydot0));
755     }
756     if (ark->extrapolate && !ts->steprestart) { /* Save the Y, YdotI, YdotRHS for extrapolation initial guess */
757       for (i = 0; i<s; i++) {
758         PetscCall(VecCopy(Y[i],ark->Y_prev[i]));
759         PetscCall(VecCopy(YdotRHS[i],ark->YdotRHS_prev[i]));
760         PetscCall(VecCopy(YdotI[i],ark->YdotI_prev[i]));
761       }
762     }
763   }
764 
765   if (ts->equation_type >= TS_EQ_IMPLICIT && tab->explicit_first_stage && ts->steprestart) {
766     TS ts_start;
767     if (PetscDefined(USE_DEBUG)) {
768       PetscBool id = PETSC_FALSE;
769       PetscCall(TSARKIMEXTestMassIdentity(ts,&id));
770       PetscCheck(id,PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_INCOMP,"This scheme requires an identity mass matrix, however the TSIFunction you provide does not utilize an identity mass matrix");
771     }
772     PetscCall(TSClone(ts,&ts_start));
773     PetscCall(TSSetSolution(ts_start,ts->vec_sol));
774     PetscCall(TSSetTime(ts_start,ts->ptime));
775     PetscCall(TSSetMaxSteps(ts_start,ts->steps+1));
776     PetscCall(TSSetMaxTime(ts_start,ts->ptime+ts->time_step));
777     PetscCall(TSSetExactFinalTime(ts_start,TS_EXACTFINALTIME_STEPOVER));
778     PetscCall(TSSetTimeStep(ts_start,ts->time_step));
779     PetscCall(TSSetType(ts_start,TSARKIMEX));
780     PetscCall(TSARKIMEXSetFullyImplicit(ts_start,PETSC_TRUE));
781     PetscCall(TSARKIMEXSetType(ts_start,TSARKIMEX1BEE));
782 
783     PetscCall(TSRestartStep(ts_start));
784     PetscCall(TSSolve(ts_start,ts->vec_sol));
785     PetscCall(TSGetTime(ts_start,&ts->ptime));
786     PetscCall(TSGetTimeStep(ts_start,&ts->time_step));
787 
788     { /* Save the initial slope for the next step */
789       TS_ARKIMEX *ark_start = (TS_ARKIMEX*)ts_start->data;
790       PetscCall(VecCopy(ark_start->YdotI[ark_start->tableau->s-1],Ydot0));
791     }
792     ts->steps++;
793 
794     /* Set the correct TS in SNES */
795     /* We'll try to bypass this by changing the method on the fly */
796     {
797       PetscCall(TSGetSNES(ts,&snes));
798       PetscCall(TSSetSNES(ts,snes));
799     }
800     PetscCall(TSDestroy(&ts_start));
801   }
802 
803   ark->status = TS_STEP_INCOMPLETE;
804   while (!ts->reason && ark->status != TS_STEP_COMPLETE) {
805     PetscReal t = ts->ptime;
806     PetscReal h = ts->time_step;
807     for (i=0; i<s; i++) {
808       ark->stage_time = t + h*ct[i];
809       PetscCall(TSPreStage(ts,ark->stage_time));
810       if (At[i*s+i] == 0) { /* This stage is explicit */
811         PetscCheck(i == 0 || ts->equation_type < TS_EQ_IMPLICIT,PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Explicit stages other than the first one are not supported for implicit problems");
812         PetscCall(VecCopy(ts->vec_sol,Y[i]));
813         for (j=0; j<i; j++) w[j] = h*At[i*s+j];
814         PetscCall(VecMAXPY(Y[i],i,w,YdotI));
815         for (j=0; j<i; j++) w[j] = h*A[i*s+j];
816         PetscCall(VecMAXPY(Y[i],i,w,YdotRHS));
817       } else {
818         ark->scoeff = 1./At[i*s+i];
819         /* Ydot = shift*(Y-Z) */
820         PetscCall(VecCopy(ts->vec_sol,Z));
821         for (j=0; j<i; j++) w[j] = h*At[i*s+j];
822         PetscCall(VecMAXPY(Z,i,w,YdotI));
823         for (j=0; j<i; j++) w[j] = h*A[i*s+j];
824         PetscCall(VecMAXPY(Z,i,w,YdotRHS));
825         if (extrapolate && !ts->steprestart) {
826           /* Initial guess extrapolated from previous time step stage values */
827           PetscCall(TSExtrapolate_ARKIMEX(ts,c[i],Y[i]));
828         } else {
829           /* Initial guess taken from last stage */
830           PetscCall(VecCopy(i>0 ? Y[i-1] : ts->vec_sol,Y[i]));
831         }
832         PetscCall(TSGetSNES(ts,&snes));
833         PetscCall(SNESSolve(snes,NULL,Y[i]));
834         PetscCall(SNESGetIterationNumber(snes,&its));
835         PetscCall(SNESGetLinearSolveIterations(snes,&lits));
836         ts->snes_its += its; ts->ksp_its += lits;
837         PetscCall(TSGetAdapt(ts,&adapt));
838         PetscCall(TSAdaptCheckStage(adapt,ts,ark->stage_time,Y[i],&stageok));
839         if (!stageok) {
840           /* We are likely rejecting the step because of solver or function domain problems so we should not attempt to
841            * use extrapolation to initialize the solves on the next attempt. */
842           extrapolate = PETSC_FALSE;
843           goto reject_step;
844         }
845       }
846       if (ts->equation_type >= TS_EQ_IMPLICIT) {
847         if (i==0 && tab->explicit_first_stage) {
848           PetscCheck(tab->stiffly_accurate,PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s is not stiffly accurate and therefore explicit-first stage methods cannot be used if the equation is implicit because the slope cannot be evaluated",ark->tableau->name);
849           PetscCall(VecCopy(Ydot0,YdotI[0]));                                      /* YdotI = YdotI(tn-1) */
850         } else {
851           PetscCall(VecAXPBYPCZ(YdotI[i],-ark->scoeff/h,ark->scoeff/h,0,Z,Y[i]));  /* YdotI = shift*(X-Z) */
852         }
853       } else {
854         if (i==0 && tab->explicit_first_stage) {
855           PetscCall(VecZeroEntries(Ydot));
856           PetscCall(TSComputeIFunction(ts,t+h*ct[i],Y[i],Ydot,YdotI[i],ark->imex));/* YdotI = -G(t,Y,0)   */
857           PetscCall(VecScale(YdotI[i],-1.0));
858         } else {
859           PetscCall(VecAXPBYPCZ(YdotI[i],-ark->scoeff/h,ark->scoeff/h,0,Z,Y[i]));  /* YdotI = shift*(X-Z) */
860         }
861         if (ark->imex) {
862           PetscCall(TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]));
863         } else {
864           PetscCall(VecZeroEntries(YdotRHS[i]));
865         }
866       }
867       PetscCall(TSPostStage(ts,ark->stage_time,i,Y));
868     }
869 
870     ark->status = TS_STEP_INCOMPLETE;
871     PetscCall(TSEvaluateStep_ARKIMEX(ts,tab->order,ts->vec_sol,NULL));
872     ark->status = TS_STEP_PENDING;
873     PetscCall(TSGetAdapt(ts,&adapt));
874     PetscCall(TSAdaptCandidatesClear(adapt));
875     PetscCall(TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE));
876     PetscCall(TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept));
877     ark->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
878     if (!accept) { /* Roll back the current step */
879       PetscCall(TSRollBack_ARKIMEX(ts));
880       ts->time_step = next_time_step;
881       goto reject_step;
882     }
883 
884     ts->ptime += ts->time_step;
885     ts->time_step = next_time_step;
886     break;
887 
888   reject_step:
889     ts->reject++; accept = PETSC_FALSE;
890     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
891       ts->reason = TS_DIVERGED_STEP_REJECTED;
892       PetscCall(PetscInfo(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections));
893     }
894   }
895   PetscFunctionReturn(0);
896 }
897 
898 static PetscErrorCode TSInterpolate_ARKIMEX(TS ts,PetscReal itime,Vec X)
899 {
900   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
901   PetscInt        s    = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j;
902   PetscReal       h;
903   PetscReal       tt,t;
904   PetscScalar     *bt,*b;
905   const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp;
906 
907   PetscFunctionBegin;
908   PetscCheck(Bt && B,PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name);
909   switch (ark->status) {
910   case TS_STEP_INCOMPLETE:
911   case TS_STEP_PENDING:
912     h = ts->time_step;
913     t = (itime - ts->ptime)/h;
914     break;
915   case TS_STEP_COMPLETE:
916     h = ts->ptime - ts->ptime_prev;
917     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
918     break;
919   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
920   }
921   PetscCall(PetscMalloc2(s,&bt,s,&b));
922   for (i=0; i<s; i++) bt[i] = b[i] = 0;
923   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
924     for (i=0; i<s; i++) {
925       bt[i] += h * Bt[i*pinterp+j] * tt;
926       b[i]  += h * B[i*pinterp+j] * tt;
927     }
928   }
929   PetscCall(VecCopy(ark->Y[0],X));
930   PetscCall(VecMAXPY(X,s,bt,ark->YdotI));
931   PetscCall(VecMAXPY(X,s,b,ark->YdotRHS));
932   PetscCall(PetscFree2(bt,b));
933   PetscFunctionReturn(0);
934 }
935 
936 static PetscErrorCode TSExtrapolate_ARKIMEX(TS ts,PetscReal c,Vec X)
937 {
938   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
939   PetscInt        s = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j;
940   PetscReal       h,h_prev,t,tt;
941   PetscScalar     *bt,*b;
942   const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp;
943 
944   PetscFunctionBegin;
945   PetscCheck(Bt && B,PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name);
946   PetscCall(PetscCalloc2(s,&bt,s,&b));
947   h = ts->time_step;
948   h_prev = ts->ptime - ts->ptime_prev;
949   t = 1 + h/h_prev*c;
950   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
951     for (i=0; i<s; i++) {
952       bt[i] += h * Bt[i*pinterp+j] * tt;
953       b[i]  += h * B[i*pinterp+j] * tt;
954     }
955   }
956   PetscCheck(ark->Y_prev,PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Stages from previous step have not been stored");
957   PetscCall(VecCopy(ark->Y_prev[0],X));
958   PetscCall(VecMAXPY(X,s,bt,ark->YdotI_prev));
959   PetscCall(VecMAXPY(X,s,b,ark->YdotRHS_prev));
960   PetscCall(PetscFree2(bt,b));
961   PetscFunctionReturn(0);
962 }
963 
964 /*------------------------------------------------------------*/
965 
966 static PetscErrorCode TSARKIMEXTableauReset(TS ts)
967 {
968   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
969   ARKTableau     tab  = ark->tableau;
970 
971   PetscFunctionBegin;
972   if (!tab) PetscFunctionReturn(0);
973   PetscCall(PetscFree(ark->work));
974   PetscCall(VecDestroyVecs(tab->s,&ark->Y));
975   PetscCall(VecDestroyVecs(tab->s,&ark->YdotI));
976   PetscCall(VecDestroyVecs(tab->s,&ark->YdotRHS));
977   PetscCall(VecDestroyVecs(tab->s,&ark->Y_prev));
978   PetscCall(VecDestroyVecs(tab->s,&ark->YdotI_prev));
979   PetscCall(VecDestroyVecs(tab->s,&ark->YdotRHS_prev));
980   PetscFunctionReturn(0);
981 }
982 
983 static PetscErrorCode TSReset_ARKIMEX(TS ts)
984 {
985   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
986 
987   PetscFunctionBegin;
988   PetscCall(TSARKIMEXTableauReset(ts));
989   PetscCall(VecDestroy(&ark->Ydot));
990   PetscCall(VecDestroy(&ark->Ydot0));
991   PetscCall(VecDestroy(&ark->Z));
992   PetscFunctionReturn(0);
993 }
994 
995 static PetscErrorCode TSARKIMEXGetVecs(TS ts,DM dm,Vec *Z,Vec *Ydot)
996 {
997   TS_ARKIMEX     *ax = (TS_ARKIMEX*)ts->data;
998 
999   PetscFunctionBegin;
1000   if (Z) {
1001     if (dm && dm != ts->dm) {
1002       PetscCall(DMGetNamedGlobalVector(dm,"TSARKIMEX_Z",Z));
1003     } else *Z = ax->Z;
1004   }
1005   if (Ydot) {
1006     if (dm && dm != ts->dm) {
1007       PetscCall(DMGetNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot));
1008     } else *Ydot = ax->Ydot;
1009   }
1010   PetscFunctionReturn(0);
1011 }
1012 
1013 static PetscErrorCode TSARKIMEXRestoreVecs(TS ts,DM dm,Vec *Z,Vec *Ydot)
1014 {
1015   PetscFunctionBegin;
1016   if (Z) {
1017     if (dm && dm != ts->dm) {
1018       PetscCall(DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Z",Z));
1019     }
1020   }
1021   if (Ydot) {
1022     if (dm && dm != ts->dm) {
1023       PetscCall(DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot));
1024     }
1025   }
1026   PetscFunctionReturn(0);
1027 }
1028 
1029 /*
1030   This defines the nonlinear equation that is to be solved with SNES
1031   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
1032 */
1033 static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes,Vec X,Vec F,TS ts)
1034 {
1035   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1036   DM             dm,dmsave;
1037   Vec            Z,Ydot;
1038   PetscReal      shift = ark->scoeff / ts->time_step;
1039 
1040   PetscFunctionBegin;
1041   PetscCall(SNESGetDM(snes,&dm));
1042   PetscCall(TSARKIMEXGetVecs(ts,dm,&Z,&Ydot));
1043   PetscCall(VecAXPBYPCZ(Ydot,-shift,shift,0,Z,X)); /* Ydot = shift*(X-Z) */
1044   dmsave = ts->dm;
1045   ts->dm = dm;
1046 
1047   PetscCall(TSComputeIFunction(ts,ark->stage_time,X,Ydot,F,ark->imex));
1048 
1049   ts->dm = dmsave;
1050   PetscCall(TSARKIMEXRestoreVecs(ts,dm,&Z,&Ydot));
1051   PetscFunctionReturn(0);
1052 }
1053 
1054 static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes,Vec X,Mat A,Mat B,TS ts)
1055 {
1056   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1057   DM             dm,dmsave;
1058   Vec            Ydot;
1059   PetscReal      shift = ark->scoeff / ts->time_step;
1060 
1061   PetscFunctionBegin;
1062   PetscCall(SNESGetDM(snes,&dm));
1063   PetscCall(TSARKIMEXGetVecs(ts,dm,NULL,&Ydot));
1064   /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */
1065   dmsave = ts->dm;
1066   ts->dm = dm;
1067 
1068   PetscCall(TSComputeIJacobian(ts,ark->stage_time,X,Ydot,shift,A,B,ark->imex));
1069 
1070   ts->dm = dmsave;
1071   PetscCall(TSARKIMEXRestoreVecs(ts,dm,NULL,&Ydot));
1072   PetscFunctionReturn(0);
1073 }
1074 
1075 static PetscErrorCode DMCoarsenHook_TSARKIMEX(DM fine,DM coarse,void *ctx)
1076 {
1077   PetscFunctionBegin;
1078   PetscFunctionReturn(0);
1079 }
1080 
1081 static PetscErrorCode DMRestrictHook_TSARKIMEX(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
1082 {
1083   TS             ts = (TS)ctx;
1084   Vec            Z,Z_c;
1085 
1086   PetscFunctionBegin;
1087   PetscCall(TSARKIMEXGetVecs(ts,fine,&Z,NULL));
1088   PetscCall(TSARKIMEXGetVecs(ts,coarse,&Z_c,NULL));
1089   PetscCall(MatRestrict(restrct,Z,Z_c));
1090   PetscCall(VecPointwiseMult(Z_c,rscale,Z_c));
1091   PetscCall(TSARKIMEXRestoreVecs(ts,fine,&Z,NULL));
1092   PetscCall(TSARKIMEXRestoreVecs(ts,coarse,&Z_c,NULL));
1093   PetscFunctionReturn(0);
1094 }
1095 
1096 static PetscErrorCode DMSubDomainHook_TSARKIMEX(DM dm,DM subdm,void *ctx)
1097 {
1098   PetscFunctionBegin;
1099   PetscFunctionReturn(0);
1100 }
1101 
1102 static PetscErrorCode DMSubDomainRestrictHook_TSARKIMEX(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
1103 {
1104   TS             ts = (TS)ctx;
1105   Vec            Z,Z_c;
1106 
1107   PetscFunctionBegin;
1108   PetscCall(TSARKIMEXGetVecs(ts,dm,&Z,NULL));
1109   PetscCall(TSARKIMEXGetVecs(ts,subdm,&Z_c,NULL));
1110 
1111   PetscCall(VecScatterBegin(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD));
1112   PetscCall(VecScatterEnd(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD));
1113 
1114   PetscCall(TSARKIMEXRestoreVecs(ts,dm,&Z,NULL));
1115   PetscCall(TSARKIMEXRestoreVecs(ts,subdm,&Z_c,NULL));
1116   PetscFunctionReturn(0);
1117 }
1118 
1119 static PetscErrorCode TSARKIMEXTableauSetUp(TS ts)
1120 {
1121   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1122   ARKTableau     tab  = ark->tableau;
1123 
1124   PetscFunctionBegin;
1125   PetscCall(PetscMalloc1(tab->s,&ark->work));
1126   PetscCall(VecDuplicateVecs(ts->vec_sol,tab->s,&ark->Y));
1127   PetscCall(VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotI));
1128   PetscCall(VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotRHS));
1129   if (ark->extrapolate) {
1130     PetscCall(VecDuplicateVecs(ts->vec_sol,tab->s,&ark->Y_prev));
1131     PetscCall(VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotI_prev));
1132     PetscCall(VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotRHS_prev));
1133   }
1134   PetscFunctionReturn(0);
1135 }
1136 
1137 static PetscErrorCode TSSetUp_ARKIMEX(TS ts)
1138 {
1139   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1140   DM             dm;
1141   SNES           snes;
1142 
1143   PetscFunctionBegin;
1144   PetscCall(TSARKIMEXTableauSetUp(ts));
1145   PetscCall(VecDuplicate(ts->vec_sol,&ark->Ydot));
1146   PetscCall(VecDuplicate(ts->vec_sol,&ark->Ydot0));
1147   PetscCall(VecDuplicate(ts->vec_sol,&ark->Z));
1148   PetscCall(TSGetDM(ts,&dm));
1149   PetscCall(DMCoarsenHookAdd(dm,DMCoarsenHook_TSARKIMEX,DMRestrictHook_TSARKIMEX,ts));
1150   PetscCall(DMSubDomainHookAdd(dm,DMSubDomainHook_TSARKIMEX,DMSubDomainRestrictHook_TSARKIMEX,ts));
1151   PetscCall(TSGetSNES(ts,&snes));
1152   PetscFunctionReturn(0);
1153 }
1154 /*------------------------------------------------------------*/
1155 
1156 static PetscErrorCode TSSetFromOptions_ARKIMEX(PetscOptionItems *PetscOptionsObject,TS ts)
1157 {
1158   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1159 
1160   PetscFunctionBegin;
1161   PetscCall(PetscOptionsHead(PetscOptionsObject,"ARKIMEX ODE solver options"));
1162   {
1163     ARKTableauLink link;
1164     PetscInt       count,choice;
1165     PetscBool      flg;
1166     const char     **namelist;
1167     for (link=ARKTableauList,count=0; link; link=link->next,count++) ;
1168     PetscCall(PetscMalloc1(count,(char***)&namelist));
1169     for (link=ARKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
1170     PetscCall(PetscOptionsEList("-ts_arkimex_type","Family of ARK IMEX method","TSARKIMEXSetType",(const char*const*)namelist,count,ark->tableau->name,&choice,&flg));
1171     if (flg) PetscCall(TSARKIMEXSetType(ts,namelist[choice]));
1172     PetscCall(PetscFree(namelist));
1173 
1174     flg  = (PetscBool) !ark->imex;
1175     PetscCall(PetscOptionsBool("-ts_arkimex_fully_implicit","Solve the problem fully implicitly","TSARKIMEXSetFullyImplicit",flg,&flg,NULL));
1176     ark->imex = (PetscBool) !flg;
1177     PetscCall(PetscOptionsBool("-ts_arkimex_initial_guess_extrapolate","Extrapolate the initial guess for the stage solution from stage values of the previous time step","",ark->extrapolate,&ark->extrapolate,NULL));
1178   }
1179   PetscCall(PetscOptionsTail());
1180   PetscFunctionReturn(0);
1181 }
1182 
1183 static PetscErrorCode TSView_ARKIMEX(TS ts,PetscViewer viewer)
1184 {
1185   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1186   PetscBool      iascii;
1187 
1188   PetscFunctionBegin;
1189   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
1190   if (iascii) {
1191     ARKTableau    tab = ark->tableau;
1192     TSARKIMEXType arktype;
1193     char          buf[512];
1194     PetscBool     flg;
1195 
1196     PetscCall(TSARKIMEXGetType(ts,&arktype));
1197     PetscCall(TSARKIMEXGetFullyImplicit(ts,&flg));
1198     PetscCall(PetscViewerASCIIPrintf(viewer,"  ARK IMEX %s\n",arktype));
1199     PetscCall(PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->ct));
1200     PetscCall(PetscViewerASCIIPrintf(viewer,"  Stiff abscissa       ct = %s\n",buf));
1201     PetscCall(PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c));
1202     PetscCall(PetscViewerASCIIPrintf(viewer,"Fully implicit: %s\n",flg ? "yes" : "no"));
1203     PetscCall(PetscViewerASCIIPrintf(viewer,"Stiffly accurate: %s\n",tab->stiffly_accurate ? "yes" : "no"));
1204     PetscCall(PetscViewerASCIIPrintf(viewer,"Explicit first stage: %s\n",tab->explicit_first_stage ? "yes" : "no"));
1205     PetscCall(PetscViewerASCIIPrintf(viewer,"FSAL property: %s\n",tab->FSAL_implicit ? "yes" : "no"));
1206     PetscCall(PetscViewerASCIIPrintf(viewer,"  Nonstiff abscissa     c = %s\n",buf));
1207   }
1208   PetscFunctionReturn(0);
1209 }
1210 
1211 static PetscErrorCode TSLoad_ARKIMEX(TS ts,PetscViewer viewer)
1212 {
1213   SNES           snes;
1214   TSAdapt        adapt;
1215 
1216   PetscFunctionBegin;
1217   PetscCall(TSGetAdapt(ts,&adapt));
1218   PetscCall(TSAdaptLoad(adapt,viewer));
1219   PetscCall(TSGetSNES(ts,&snes));
1220   PetscCall(SNESLoad(snes,viewer));
1221   /* function and Jacobian context for SNES when used with TS is always ts object */
1222   PetscCall(SNESSetFunction(snes,NULL,NULL,ts));
1223   PetscCall(SNESSetJacobian(snes,NULL,NULL,NULL,ts));
1224   PetscFunctionReturn(0);
1225 }
1226 
1227 /*@C
1228   TSARKIMEXSetType - Set the type of ARK IMEX scheme
1229 
1230   Logically collective
1231 
1232   Input Parameters:
1233 +  ts - timestepping context
1234 -  arktype - type of ARK-IMEX scheme
1235 
1236   Options Database:
1237 .  -ts_arkimex_type <1bee,a2,l2,ars122,2c,2d,2e,prssp2,3,bpr3,ars443,4,5> - set ARK IMEX scheme type
1238 
1239   Level: intermediate
1240 
1241 .seealso: TSARKIMEXGetType(), TSARKIMEX, TSARKIMEXType, TSARKIMEX1BEE, TSARKIMEXA2, TSARKIMEXL2, TSARKIMEXARS122, TSARKIMEX2C, TSARKIMEX2D, TSARKIMEX2E, TSARKIMEXPRSSP2,
1242           TSARKIMEX3, TSARKIMEXBPR3, TSARKIMEXARS443, TSARKIMEX4, TSARKIMEX5
1243 @*/
1244 PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType arktype)
1245 {
1246   PetscFunctionBegin;
1247   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1248   PetscValidCharPointer(arktype,2);
1249   PetscTryMethod(ts,"TSARKIMEXSetType_C",(TS,TSARKIMEXType),(ts,arktype));
1250   PetscFunctionReturn(0);
1251 }
1252 
1253 /*@C
1254   TSARKIMEXGetType - Get the type of ARK IMEX scheme
1255 
1256   Logically collective
1257 
1258   Input Parameter:
1259 .  ts - timestepping context
1260 
1261   Output Parameter:
1262 .  arktype - type of ARK-IMEX scheme
1263 
1264   Level: intermediate
1265 
1266 .seealso: TSARKIMEXGetType()
1267 @*/
1268 PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType *arktype)
1269 {
1270   PetscFunctionBegin;
1271   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1272   PetscUseMethod(ts,"TSARKIMEXGetType_C",(TS,TSARKIMEXType*),(ts,arktype));
1273   PetscFunctionReturn(0);
1274 }
1275 
1276 /*@
1277   TSARKIMEXSetFullyImplicit - Solve both parts of the equation implicitly
1278 
1279   Logically collective
1280 
1281   Input Parameters:
1282 +  ts - timestepping context
1283 -  flg - PETSC_TRUE for fully implicit
1284 
1285   Level: intermediate
1286 
1287 .seealso: TSARKIMEXGetType(), TSARKIMEXGetFullyImplicit()
1288 @*/
1289 PetscErrorCode TSARKIMEXSetFullyImplicit(TS ts,PetscBool flg)
1290 {
1291   PetscFunctionBegin;
1292   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1293   PetscValidLogicalCollectiveBool(ts,flg,2);
1294   PetscTryMethod(ts,"TSARKIMEXSetFullyImplicit_C",(TS,PetscBool),(ts,flg));
1295   PetscFunctionReturn(0);
1296 }
1297 
1298 /*@
1299   TSARKIMEXGetFullyImplicit - Inquires if both parts of the equation are solved implicitly
1300 
1301   Logically collective
1302 
1303   Input Parameter:
1304 .  ts - timestepping context
1305 
1306   Output Parameter:
1307 .  flg - PETSC_TRUE for fully implicit
1308 
1309   Level: intermediate
1310 
1311 .seealso: TSARKIMEXGetType(), TSARKIMEXSetFullyImplicit()
1312 @*/
1313 PetscErrorCode TSARKIMEXGetFullyImplicit(TS ts,PetscBool *flg)
1314 {
1315   PetscFunctionBegin;
1316   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1317   PetscValidBoolPointer(flg,2);
1318   PetscUseMethod(ts,"TSARKIMEXGetFullyImplicit_C",(TS,PetscBool*),(ts,flg));
1319   PetscFunctionReturn(0);
1320 }
1321 
1322 static PetscErrorCode  TSARKIMEXGetType_ARKIMEX(TS ts,TSARKIMEXType *arktype)
1323 {
1324   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1325 
1326   PetscFunctionBegin;
1327   *arktype = ark->tableau->name;
1328   PetscFunctionReturn(0);
1329 }
1330 static PetscErrorCode  TSARKIMEXSetType_ARKIMEX(TS ts,TSARKIMEXType arktype)
1331 {
1332   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1333   PetscBool      match;
1334   ARKTableauLink link;
1335 
1336   PetscFunctionBegin;
1337   if (ark->tableau) {
1338     PetscCall(PetscStrcmp(ark->tableau->name,arktype,&match));
1339     if (match) PetscFunctionReturn(0);
1340   }
1341   for (link = ARKTableauList; link; link=link->next) {
1342     PetscCall(PetscStrcmp(link->tab.name,arktype,&match));
1343     if (match) {
1344       if (ts->setupcalled) PetscCall(TSARKIMEXTableauReset(ts));
1345       ark->tableau = &link->tab;
1346       if (ts->setupcalled) PetscCall(TSARKIMEXTableauSetUp(ts));
1347       ts->default_adapt_type = ark->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
1348       PetscFunctionReturn(0);
1349     }
1350   }
1351   SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",arktype);
1352 }
1353 
1354 static PetscErrorCode  TSARKIMEXSetFullyImplicit_ARKIMEX(TS ts,PetscBool flg)
1355 {
1356   TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
1357 
1358   PetscFunctionBegin;
1359   ark->imex = (PetscBool)!flg;
1360   PetscFunctionReturn(0);
1361 }
1362 
1363 static PetscErrorCode  TSARKIMEXGetFullyImplicit_ARKIMEX(TS ts,PetscBool *flg)
1364 {
1365   TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
1366 
1367   PetscFunctionBegin;
1368   *flg = (PetscBool)!ark->imex;
1369   PetscFunctionReturn(0);
1370 }
1371 
1372 static PetscErrorCode TSDestroy_ARKIMEX(TS ts)
1373 {
1374   PetscFunctionBegin;
1375   PetscCall(TSReset_ARKIMEX(ts));
1376   if (ts->dm) {
1377     PetscCall(DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSARKIMEX,DMRestrictHook_TSARKIMEX,ts));
1378     PetscCall(DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSARKIMEX,DMSubDomainRestrictHook_TSARKIMEX,ts));
1379   }
1380   PetscCall(PetscFree(ts->data));
1381   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetType_C",NULL));
1382   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetType_C",NULL));
1383   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C",NULL));
1384   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C",NULL));
1385   PetscFunctionReturn(0);
1386 }
1387 
1388 /* ------------------------------------------------------------ */
1389 /*MC
1390       TSARKIMEX - ODE and DAE solver using additive Runge-Kutta IMEX schemes
1391 
1392   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
1393   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
1394   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
1395 
1396   Notes:
1397   The default is TSARKIMEX3, it can be changed with TSARKIMEXSetType() or -ts_arkimex_type
1398 
1399   If the equation is implicit or a DAE, then TSSetEquationType() needs to be set accordingly. Refer to the manual for further information.
1400 
1401   Methods with an explicit stage can only be used with ODE in which the stiff part G(t,X,Xdot) has the form Xdot + Ghat(t,X).
1402 
1403   Consider trying TSROSW if the stiff part is linear or weakly nonlinear.
1404 
1405   Level: beginner
1406 
1407 .seealso:  TSCreate(), TS, TSSetType(), TSARKIMEXSetType(), TSARKIMEXGetType(), TSARKIMEXSetFullyImplicit(), TSARKIMEXGetFullyImplicit(),
1408            TSARKIMEX1BEE, TSARKIMEX2C, TSARKIMEX2D, TSARKIMEX2E, TSARKIMEX3, TSARKIMEXL2, TSARKIMEXA2, TSARKIMEXARS122,
1409            TSARKIMEX4, TSARKIMEX5, TSARKIMEXPRSSP2, TSARKIMEXARS443, TSARKIMEXBPR3, TSARKIMEXType, TSARKIMEXRegister()
1410 
1411 M*/
1412 PETSC_EXTERN PetscErrorCode TSCreate_ARKIMEX(TS ts)
1413 {
1414   TS_ARKIMEX     *th;
1415 
1416   PetscFunctionBegin;
1417   PetscCall(TSARKIMEXInitializePackage());
1418 
1419   ts->ops->reset          = TSReset_ARKIMEX;
1420   ts->ops->destroy        = TSDestroy_ARKIMEX;
1421   ts->ops->view           = TSView_ARKIMEX;
1422   ts->ops->load           = TSLoad_ARKIMEX;
1423   ts->ops->setup          = TSSetUp_ARKIMEX;
1424   ts->ops->step           = TSStep_ARKIMEX;
1425   ts->ops->interpolate    = TSInterpolate_ARKIMEX;
1426   ts->ops->evaluatestep   = TSEvaluateStep_ARKIMEX;
1427   ts->ops->rollback       = TSRollBack_ARKIMEX;
1428   ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX;
1429   ts->ops->snesfunction   = SNESTSFormFunction_ARKIMEX;
1430   ts->ops->snesjacobian   = SNESTSFormJacobian_ARKIMEX;
1431 
1432   ts->usessnes = PETSC_TRUE;
1433 
1434   PetscCall(PetscNewLog(ts,&th));
1435   ts->data = (void*)th;
1436   th->imex = PETSC_TRUE;
1437 
1438   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetType_C",TSARKIMEXGetType_ARKIMEX));
1439   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetType_C",TSARKIMEXSetType_ARKIMEX));
1440   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C",TSARKIMEXSetFullyImplicit_ARKIMEX));
1441   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetFullyImplicit_C",TSARKIMEXGetFullyImplicit_ARKIMEX));
1442 
1443   PetscCall(TSARKIMEXSetType(ts,TSARKIMEXDefault));
1444   PetscFunctionReturn(0);
1445 }
1446