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