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