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