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