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