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