xref: /petsc/src/ts/impls/arkimex/arkimex.c (revision 42ea6711cbf0098c0cd9833f8b285244145c85c5)
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. Spitheri, 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. Spitheri, 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,W = ark->Work,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 
710   PetscFunctionBegin;
711   if (ts->equation_type >= TS_EQ_IMPLICIT && tab->explicit_first_stage) {
712     PetscReal valid_time;
713     PetscBool isvalid;
714     ierr = PetscObjectComposedDataGetReal((PetscObject)ts->vec_sol,
715                                           explicit_stage_time_id,
716                                           valid_time,
717                                           isvalid);
718     CHKERRQ(ierr);
719     if (!isvalid || valid_time != ts->ptime) {
720       TS        ts_start;
721       SNES      snes_start;
722       DM        dm;
723       PetscReal atol;
724       Vec       vatol;
725       PetscReal rtol;
726       Vec       vrtol;
727 
728       ierr = TSCreate(PetscObjectComm((PetscObject)ts),&ts_start);CHKERRQ(ierr);
729       ierr = TSGetSNES(ts,&snes_start);CHKERRQ(ierr);
730       ierr = TSSetSNES(ts_start,snes_start);CHKERRQ(ierr);
731       ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
732       ierr = TSSetDM(ts_start,dm);CHKERRQ(ierr);
733 
734       ts_start->adapt=ts->adapt;
735       PetscObjectReference((PetscObject)ts_start->adapt);
736 
737       ierr = TSSetSolution(ts_start,ts->vec_sol);CHKERRQ(ierr);
738       ierr = TSSetTime(ts_start,ts->ptime);CHKERRQ(ierr);
739       ierr = TSSetDuration(ts_start,1,ts->ptime+ts->time_step);CHKERRQ(ierr);
740       ierr = TSSetTimeStep(ts_start,ts->time_step);CHKERRQ(ierr);
741       ierr = TSSetType(ts_start,TSARKIMEX);CHKERRQ(ierr);
742       ierr = TSARKIMEXSetFullyImplicit(ts_start,PETSC_TRUE);CHKERRQ(ierr);
743       ierr = TSARKIMEXSetType(ts_start,TSARKIMEX1BEE);CHKERRQ(ierr);
744       ierr = TSSetEquationType(ts_start,ts->equation_type);CHKERRQ(ierr);
745       ierr = TSGetTolerances(ts,&atol,&vatol,&rtol,&vrtol);CHKERRQ(ierr);
746       ierr = TSSetTolerances(ts_start,atol,vatol,rtol,vrtol);CHKERRQ(ierr);
747       ierr = TSSolve(ts_start,ts->vec_sol);CHKERRQ(ierr);
748       ierr = TSGetTime(ts_start,&ts->ptime);CHKERRQ(ierr);
749 
750       ts->time_step = ts_start->time_step;
751       ts->steps++;
752       ierr = VecCopy(((TS_ARKIMEX*)ts_start->data)->Ydot0,Ydot0);CHKERRQ(ierr);
753       ts_start->snes=NULL;
754       ierr = TSSetSNES(ts,snes_start);CHKERRQ(ierr);
755       ierr = SNESDestroy(&snes_start);CHKERRQ(ierr);
756       ierr = TSDestroy(&ts_start);CHKERRQ(ierr);
757     }
758   }
759 
760   ierr           = TSGetSNES(ts,&snes);CHKERRQ(ierr);
761   t              = ts->ptime;
762   next_time_step = ts->time_step;
763   accept         = PETSC_TRUE;
764   ark->status    = TS_STEP_INCOMPLETE;
765 
766 
767   for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) {
768     PetscReal h = ts->time_step;
769     ierr = TSPreStep(ts);CHKERRQ(ierr);
770     for (i=0; i<s; i++) {
771       ark->stage_time = t + h*ct[i];
772       if (At[i*s+i] == 0) {           /* This stage is explicit */
773         ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
774         for (j=0; j<i; j++) w[j] = h*At[i*s+j];
775         ierr = VecMAXPY(Y[i],i,w,YdotI);CHKERRQ(ierr);
776         for (j=0; j<i; j++) w[j] = h*A[i*s+j];
777         ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
778       } else {
779         ark->scoeff     = 1./At[i*s+i];
780         ierr            = TSPreStage(ts,ark->stage_time);CHKERRQ(ierr);
781         /* Affine part */
782         ierr = VecZeroEntries(W);CHKERRQ(ierr);
783         for (j=0; j<i; j++) w[j] = h*A[i*s+j];
784         ierr = VecMAXPY(W,i,w,YdotRHS);CHKERRQ(ierr);
785         ierr = VecScale(W, ark->scoeff/h);CHKERRQ(ierr);
786 
787         /* Ydot = shift*(Y-Z) */
788         ierr = VecCopy(ts->vec_sol,Z);CHKERRQ(ierr);
789         for (j=0; j<i; j++) w[j] = h*At[i*s+j];
790         ierr = VecMAXPY(Z,i,w,YdotI);CHKERRQ(ierr);
791 
792         if (init_guess_extrp && ark->prev_step_valid) {
793           /* Initial guess extrapolated from previous time step stage values */
794           ierr        = TSExtrapolate_ARKIMEX(ts,c[i],Y[i]);CHKERRQ(ierr);
795         } else {
796           /* Initial guess taken from last stage */
797           ierr        = VecCopy(i>0 ? Y[i-1] : ts->vec_sol,Y[i]);CHKERRQ(ierr);
798         }
799         ierr          = SNESSolve(snes,W,Y[i]);CHKERRQ(ierr);
800         ierr          = (ts->ops->snesfunction)(snes,Y[i],W,ts);CHKERRQ(ierr);
801         ierr          = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
802         ierr          = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
803         ts->snes_its += its; ts->ksp_its += lits;
804         ierr          = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
805         ierr          = TSAdaptCheckStage(adapt,ts,&accept);CHKERRQ(ierr);
806         if (!accept) {
807           /* We are likely rejecting the step because of solver or function domain problems so we should not attempt to
808            * use extrapolation to initialize the solves on the next attempt. */
809           ark->prev_step_valid = PETSC_FALSE;
810           goto reject_step;
811         }
812       }
813       ierr = TSPostStage(ts,ark->stage_time,i,Y); CHKERRQ(ierr);
814       if (ts->equation_type>=TS_EQ_IMPLICIT) {
815         if (i==0 && tab->explicit_first_stage) {
816           ierr = VecCopy(Ydot0,YdotI[0]);CHKERRQ(ierr);
817         } else {
818           ierr = VecAXPBYPCZ(YdotI[i],-ark->scoeff/h,ark->scoeff/h,0,Z,Y[i]);CHKERRQ(ierr); /* Ydot = shift*(X-Z) */
819         }
820       } else {
821         ierr = VecZeroEntries(Ydot);CHKERRQ(ierr);
822         ierr = TSComputeIFunction(ts,t+h*ct[i],Y[i],Ydot,YdotI[i],ark->imex);CHKERRQ(ierr);
823         ierr = VecScale(YdotI[i], -1.0);CHKERRQ(ierr);
824         if (ark->imex) {
825           ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
826         } else {
827           ierr = VecZeroEntries(YdotRHS[i]);CHKERRQ(ierr);
828         }
829       }
830     }
831     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
832     ark->status = TS_STEP_PENDING;
833 
834     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
835     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
836     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
837     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
838     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
839     if (accept) {
840       /* ignore next_scheme for now */
841       ts->ptime    += ts->time_step;
842       ts->time_step = next_time_step;
843       ts->steps++;
844       if (ts->equation_type>=TS_EQ_IMPLICIT) { /* save the initial slope for the next step*/
845         ierr = VecCopy(YdotI[s-1],Ydot0);CHKERRQ(ierr);
846       }
847       ark->status = TS_STEP_COMPLETE;
848       if (tab->explicit_first_stage) {
849         ierr = PetscObjectComposedDataSetReal((PetscObject)ts->vec_sol,explicit_stage_time_id,ts->ptime);CHKERRQ(ierr);
850       }
851       /* Save the Y, YdotI, YdotRHS for extrapolation initial guess */
852       if (ark->init_guess_extrp) {
853         for (i = 0; i<s; i++) {
854           ierr = VecCopy(Y[i],ark->Y_prev[i]);CHKERRQ(ierr);
855           ierr = VecCopy(YdotRHS[i],ark->YdotRHS_prev[i]);CHKERRQ(ierr);
856           ierr = VecCopy(YdotI[i],ark->YdotI_prev[i]);CHKERRQ(ierr);
857         }
858         ark->prev_step_valid = PETSC_TRUE;
859       }
860       break;
861     } else {                    /* Roll back the current step */
862       ts->ptime += next_time_step; /* This will be undone in rollback */
863       ark->status = TS_STEP_INCOMPLETE;
864       ierr = TSRollBack(ts);CHKERRQ(ierr);
865     }
866 reject_step: continue;
867   }
868   if (ark->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED;
869   PetscFunctionReturn(0);
870 }
871 
872 #undef __FUNCT__
873 #define __FUNCT__ "TSInterpolate_ARKIMEX"
874 static PetscErrorCode TSInterpolate_ARKIMEX(TS ts,PetscReal itime,Vec X)
875 {
876   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
877   PetscInt        s    = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j;
878   PetscReal       h;
879   PetscReal       tt,t;
880   PetscScalar     *bt,*b;
881   const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp;
882   PetscErrorCode  ierr;
883 
884   PetscFunctionBegin;
885   if (!Bt || !B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name);
886   switch (ark->status) {
887   case TS_STEP_INCOMPLETE:
888   case TS_STEP_PENDING:
889     h = ts->time_step;
890     t = (itime - ts->ptime)/h;
891     break;
892   case TS_STEP_COMPLETE:
893     h = ts->time_step_prev;
894     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
895     break;
896   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
897   }
898   ierr = PetscMalloc2(s,&bt,s,&b);CHKERRQ(ierr);
899   for (i=0; i<s; i++) bt[i] = b[i] = 0;
900   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
901     for (i=0; i<s; i++) {
902       bt[i] += h * Bt[i*pinterp+j] * tt;
903       b[i]  += h * B[i*pinterp+j] * tt;
904     }
905   }
906   ierr = VecCopy(ark->Y[0],X);CHKERRQ(ierr);
907   ierr = VecMAXPY(X,s,bt,ark->YdotI);CHKERRQ(ierr);
908   ierr = VecMAXPY(X,s,b,ark->YdotRHS);CHKERRQ(ierr);
909   ierr = PetscFree2(bt,b);CHKERRQ(ierr);
910   PetscFunctionReturn(0);
911 }
912 
913 #undef __FUNCT__
914 #define __FUNCT__ "TSExtrapolate_ARKIMEX"
915 static PetscErrorCode TSExtrapolate_ARKIMEX(TS ts,PetscReal c,Vec X)
916 {
917   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
918   PetscInt        s    = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j;
919   PetscReal       h;
920   PetscReal       tt,t;
921   PetscScalar     *bt,*b;
922   const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp;
923   PetscErrorCode  ierr;
924 
925   PetscFunctionBegin;
926   if (!Bt || !B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name);
927   t = 1.0 + (ts->time_step/ts->time_step_prev)*c;
928   h = ts->time_step;
929   ierr = PetscMalloc2(s,&bt,s,&b);CHKERRQ(ierr);
930   for (i=0; i<s; i++) bt[i] = b[i] = 0;
931   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
932     for (i=0; i<s; i++) {
933       bt[i] += h * Bt[i*pinterp+j] * tt;
934       b[i]  += h * B[i*pinterp+j] * tt;
935     }
936   }
937   if (!ark->prev_step_valid) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Stages from previous step have not been stored");
938   ierr = VecCopy(ark->Y_prev[0],X);CHKERRQ(ierr);
939   ierr = VecMAXPY(X,s,bt,ark->YdotI_prev);CHKERRQ(ierr);
940   ierr = VecMAXPY(X,s,b,ark->YdotRHS_prev);CHKERRQ(ierr);
941   ierr = PetscFree2(bt,b);CHKERRQ(ierr);
942   PetscFunctionReturn(0);
943 }
944 
945 /*------------------------------------------------------------*/
946 #undef __FUNCT__
947 #define __FUNCT__ "TSReset_ARKIMEX"
948 static PetscErrorCode TSReset_ARKIMEX(TS ts)
949 {
950   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
951   PetscInt       s;
952   PetscErrorCode ierr;
953 
954   PetscFunctionBegin;
955   if (!ark->tableau) PetscFunctionReturn(0);
956   s    = ark->tableau->s;
957   ierr = VecDestroyVecs(s,&ark->Y);CHKERRQ(ierr);
958   ierr = VecDestroyVecs(s,&ark->YdotI);CHKERRQ(ierr);
959   ierr = VecDestroyVecs(s,&ark->YdotRHS);CHKERRQ(ierr);
960   if (&ark->init_guess_extrp) {
961     ierr = VecDestroyVecs(s,&ark->Y_prev);CHKERRQ(ierr);
962     ierr = VecDestroyVecs(s,&ark->YdotI_prev);CHKERRQ(ierr);
963     ierr = VecDestroyVecs(s,&ark->YdotRHS_prev);CHKERRQ(ierr);
964   }
965   ierr = VecDestroy(&ark->Ydot);CHKERRQ(ierr);
966   ierr = VecDestroy(&ark->Work);CHKERRQ(ierr);
967   ierr = VecDestroy(&ark->Ydot0);CHKERRQ(ierr);
968   ierr = VecDestroy(&ark->Z);CHKERRQ(ierr);
969   ierr = PetscFree(ark->work);CHKERRQ(ierr);
970   PetscFunctionReturn(0);
971 }
972 
973 #undef __FUNCT__
974 #define __FUNCT__ "TSDestroy_ARKIMEX"
975 static PetscErrorCode TSDestroy_ARKIMEX(TS ts)
976 {
977   PetscErrorCode ierr;
978 
979   PetscFunctionBegin;
980   ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr);
981   ierr = PetscFree(ts->data);CHKERRQ(ierr);
982   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetType_C",NULL);CHKERRQ(ierr);
983   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetType_C",NULL);CHKERRQ(ierr);
984   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C",NULL);CHKERRQ(ierr);
985   PetscFunctionReturn(0);
986 }
987 
988 
989 #undef __FUNCT__
990 #define __FUNCT__ "TSARKIMEXGetVecs"
991 static PetscErrorCode TSARKIMEXGetVecs(TS ts,DM dm,Vec *Z,Vec *Ydot)
992 {
993   TS_ARKIMEX     *ax = (TS_ARKIMEX*)ts->data;
994   PetscErrorCode ierr;
995 
996   PetscFunctionBegin;
997   if (Z) {
998     if (dm && dm != ts->dm) {
999       ierr = DMGetNamedGlobalVector(dm,"TSARKIMEX_Z",Z);CHKERRQ(ierr);
1000     } else *Z = ax->Z;
1001   }
1002   if (Ydot) {
1003     if (dm && dm != ts->dm) {
1004       ierr = DMGetNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot);CHKERRQ(ierr);
1005     } else *Ydot = ax->Ydot;
1006   }
1007   PetscFunctionReturn(0);
1008 }
1009 
1010 
1011 #undef __FUNCT__
1012 #define __FUNCT__ "TSARKIMEXRestoreVecs"
1013 static PetscErrorCode TSARKIMEXRestoreVecs(TS ts,DM dm,Vec *Z,Vec *Ydot)
1014 {
1015   PetscErrorCode ierr;
1016 
1017   PetscFunctionBegin;
1018   if (Z) {
1019     if (dm && dm != ts->dm) {
1020       ierr = DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Z",Z);CHKERRQ(ierr);
1021     }
1022   }
1023   if (Ydot) {
1024     if (dm && dm != ts->dm) {
1025       ierr = DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot);CHKERRQ(ierr);
1026     }
1027   }
1028   PetscFunctionReturn(0);
1029 }
1030 
1031 /*
1032   This defines the nonlinear equation that is to be solved with SNES
1033   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
1034 */
1035 #undef __FUNCT__
1036 #define __FUNCT__ "SNESTSFormFunction_ARKIMEX"
1037 static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes,Vec X,Vec F,TS ts)
1038 {
1039   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1040   DM             dm,dmsave;
1041   Vec            Z,Ydot;
1042   PetscReal      shift = ark->scoeff / ts->time_step;
1043   PetscErrorCode ierr;
1044 
1045   PetscFunctionBegin;
1046   ierr   = SNESGetDM(snes,&dm);CHKERRQ(ierr);
1047   ierr   = TSARKIMEXGetVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr);
1048   ierr   = VecAXPBYPCZ(Ydot,-shift,shift,0,Z,X);CHKERRQ(ierr); /* Ydot = shift*(X-Z) */
1049   dmsave = ts->dm;
1050   ts->dm = dm;
1051 
1052   ierr = TSComputeIFunction(ts,ark->stage_time,X,Ydot,F,ark->imex);CHKERRQ(ierr);
1053 
1054   ts->dm = dmsave;
1055   ierr   = TSARKIMEXRestoreVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr);
1056   PetscFunctionReturn(0);
1057 }
1058 
1059 #undef __FUNCT__
1060 #define __FUNCT__ "SNESTSFormJacobian_ARKIMEX"
1061 static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes,Vec X,Mat A,Mat B,TS ts)
1062 {
1063   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1064   DM             dm,dmsave;
1065   Vec            Ydot;
1066   PetscReal      shift = ark->scoeff / ts->time_step;
1067   PetscErrorCode ierr;
1068 
1069   PetscFunctionBegin;
1070   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
1071   ierr = TSARKIMEXGetVecs(ts,dm,NULL,&Ydot);CHKERRQ(ierr);
1072   /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */
1073   dmsave = ts->dm;
1074   ts->dm = dm;
1075 
1076   ierr = TSComputeIJacobian(ts,ark->stage_time,X,Ydot,shift,A,B,ark->imex);CHKERRQ(ierr);
1077 
1078   ts->dm = dmsave;
1079   ierr   = TSARKIMEXRestoreVecs(ts,dm,NULL,&Ydot);CHKERRQ(ierr);
1080   PetscFunctionReturn(0);
1081 }
1082 
1083 #undef __FUNCT__
1084 #define __FUNCT__ "DMCoarsenHook_TSARKIMEX"
1085 static PetscErrorCode DMCoarsenHook_TSARKIMEX(DM fine,DM coarse,void *ctx)
1086 {
1087   PetscFunctionBegin;
1088   PetscFunctionReturn(0);
1089 }
1090 
1091 #undef __FUNCT__
1092 #define __FUNCT__ "DMRestrictHook_TSARKIMEX"
1093 static PetscErrorCode DMRestrictHook_TSARKIMEX(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
1094 {
1095   TS             ts = (TS)ctx;
1096   PetscErrorCode ierr;
1097   Vec            Z,Z_c;
1098 
1099   PetscFunctionBegin;
1100   ierr = TSARKIMEXGetVecs(ts,fine,&Z,NULL);CHKERRQ(ierr);
1101   ierr = TSARKIMEXGetVecs(ts,coarse,&Z_c,NULL);CHKERRQ(ierr);
1102   ierr = MatRestrict(restrct,Z,Z_c);CHKERRQ(ierr);
1103   ierr = VecPointwiseMult(Z_c,rscale,Z_c);CHKERRQ(ierr);
1104   ierr = TSARKIMEXRestoreVecs(ts,fine,&Z,NULL);CHKERRQ(ierr);
1105   ierr = TSARKIMEXRestoreVecs(ts,coarse,&Z_c,NULL);CHKERRQ(ierr);
1106   PetscFunctionReturn(0);
1107 }
1108 
1109 
1110 #undef __FUNCT__
1111 #define __FUNCT__ "DMSubDomainHook_TSARKIMEX"
1112 static PetscErrorCode DMSubDomainHook_TSARKIMEX(DM dm,DM subdm,void *ctx)
1113 {
1114   PetscFunctionBegin;
1115   PetscFunctionReturn(0);
1116 }
1117 
1118 #undef __FUNCT__
1119 #define __FUNCT__ "DMSubDomainRestrictHook_TSARKIMEX"
1120 static PetscErrorCode DMSubDomainRestrictHook_TSARKIMEX(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
1121 {
1122   TS             ts = (TS)ctx;
1123   PetscErrorCode ierr;
1124   Vec            Z,Z_c;
1125 
1126   PetscFunctionBegin;
1127   ierr = TSARKIMEXGetVecs(ts,dm,&Z,NULL);CHKERRQ(ierr);
1128   ierr = TSARKIMEXGetVecs(ts,subdm,&Z_c,NULL);CHKERRQ(ierr);
1129 
1130   ierr = VecScatterBegin(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1131   ierr = VecScatterEnd(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1132 
1133   ierr = TSARKIMEXRestoreVecs(ts,dm,&Z,NULL);CHKERRQ(ierr);
1134   ierr = TSARKIMEXRestoreVecs(ts,subdm,&Z_c,NULL);CHKERRQ(ierr);
1135   PetscFunctionReturn(0);
1136 }
1137 
1138 #undef __FUNCT__
1139 #define __FUNCT__ "TSSetUp_ARKIMEX"
1140 static PetscErrorCode TSSetUp_ARKIMEX(TS ts)
1141 {
1142   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1143   ARKTableau     tab;
1144   PetscInt       s;
1145   PetscErrorCode ierr;
1146   DM             dm;
1147 
1148   PetscFunctionBegin;
1149   if (!ark->tableau) {
1150     ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr);
1151   }
1152   tab  = ark->tableau;
1153   s    = tab->s;
1154   ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->Y);CHKERRQ(ierr);
1155   ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotI);CHKERRQ(ierr);
1156   ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotRHS);CHKERRQ(ierr);
1157   if (ark->init_guess_extrp) {
1158     ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->Y_prev);CHKERRQ(ierr);
1159     ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotI_prev);CHKERRQ(ierr);
1160     ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotRHS_prev);CHKERRQ(ierr);
1161   }
1162   ierr = VecDuplicate(ts->vec_sol,&ark->Ydot);CHKERRQ(ierr);
1163   ierr = VecDuplicate(ts->vec_sol,&ark->Work);CHKERRQ(ierr);
1164   ierr = VecDuplicate(ts->vec_sol,&ark->Ydot0);CHKERRQ(ierr);
1165   ierr = VecDuplicate(ts->vec_sol,&ark->Z);CHKERRQ(ierr);
1166   ierr = PetscMalloc1(s,&ark->work);CHKERRQ(ierr);
1167   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1168   if (dm) {
1169     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSARKIMEX,DMRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr);
1170     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSARKIMEX,DMSubDomainRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr);
1171   }
1172   PetscFunctionReturn(0);
1173 }
1174 /*------------------------------------------------------------*/
1175 
1176 #undef __FUNCT__
1177 #define __FUNCT__ "TSSetFromOptions_ARKIMEX"
1178 static PetscErrorCode TSSetFromOptions_ARKIMEX(TS ts)
1179 {
1180   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1181   PetscErrorCode ierr;
1182   char           arktype[256];
1183 
1184   PetscFunctionBegin;
1185   ierr = PetscOptionsHead("ARKIMEX ODE solver options");CHKERRQ(ierr);
1186   {
1187     ARKTableauLink link;
1188     PetscInt       count,choice;
1189     PetscBool      flg;
1190     const char     **namelist;
1191     ierr = PetscStrncpy(arktype,TSARKIMEXDefault,sizeof(arktype));CHKERRQ(ierr);
1192     for (link=ARKTableauList,count=0; link; link=link->next,count++) ;
1193     ierr = PetscMalloc1(count,&namelist);CHKERRQ(ierr);
1194     for (link=ARKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
1195     ierr      = PetscOptionsEList("-ts_arkimex_type","Family of ARK IMEX method","TSARKIMEXSetType",(const char*const*)namelist,count,arktype,&choice,&flg);CHKERRQ(ierr);
1196     ierr      = TSARKIMEXSetType(ts,flg ? namelist[choice] : arktype);CHKERRQ(ierr);
1197     ierr      = PetscFree(namelist);CHKERRQ(ierr);
1198     flg       = (PetscBool) !ark->imex;
1199     ierr      = PetscOptionsBool("-ts_arkimex_fully_implicit","Solve the problem fully implicitly","TSARKIMEXSetFullyImplicit",flg,&flg,NULL);CHKERRQ(ierr);
1200     ark->imex = (PetscBool) !flg;
1201     ark->init_guess_extrp = PETSC_FALSE;
1202     ierr      = PetscOptionsBool("-ts_arkimex_initial_guess_extrapolate","Extrapolate the initial guess for the stage solution from stage values of the previous time step","",ark->init_guess_extrp,&ark->init_guess_extrp,NULL);CHKERRQ(ierr);
1203     ierr      = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
1204   }
1205   ierr = PetscOptionsTail();CHKERRQ(ierr);
1206   PetscFunctionReturn(0);
1207 }
1208 
1209 #undef __FUNCT__
1210 #define __FUNCT__ "PetscFormatRealArray"
1211 static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
1212 {
1213   PetscErrorCode ierr;
1214   PetscInt       i;
1215   size_t         left,count;
1216   char           *p;
1217 
1218   PetscFunctionBegin;
1219   for (i=0,p=buf,left=len; i<n; i++) {
1220     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
1221     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
1222     left -= count;
1223     p    += count;
1224     *p++  = ' ';
1225   }
1226   p[i ? 0 : -1] = 0;
1227   PetscFunctionReturn(0);
1228 }
1229 
1230 #undef __FUNCT__
1231 #define __FUNCT__ "TSView_ARKIMEX"
1232 static PetscErrorCode TSView_ARKIMEX(TS ts,PetscViewer viewer)
1233 {
1234   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1235   ARKTableau     tab  = ark->tableau;
1236   PetscBool      iascii;
1237   PetscErrorCode ierr;
1238   TSAdapt        adapt;
1239 
1240   PetscFunctionBegin;
1241   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1242   if (iascii) {
1243     TSARKIMEXType arktype;
1244     char          buf[512];
1245     ierr = TSARKIMEXGetType(ts,&arktype);CHKERRQ(ierr);
1246     ierr = PetscViewerASCIIPrintf(viewer,"  ARK IMEX %s\n",arktype);CHKERRQ(ierr);
1247     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->ct);CHKERRQ(ierr);
1248     ierr = PetscViewerASCIIPrintf(viewer,"  Stiff abscissa       ct = %s\n",buf);CHKERRQ(ierr);
1249     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
1250     ierr = PetscViewerASCIIPrintf(viewer,"Stiffly accurate: %s\n",tab->stiffly_accurate ? "yes" : "no");CHKERRQ(ierr);
1251     ierr = PetscViewerASCIIPrintf(viewer,"Explicit first stage: %s\n",tab->explicit_first_stage ? "yes" : "no");CHKERRQ(ierr);
1252     ierr = PetscViewerASCIIPrintf(viewer,"FSAL property: %s\n",tab->FSAL_implicit ? "yes" : "no");CHKERRQ(ierr);
1253     ierr = PetscViewerASCIIPrintf(viewer,"  Nonstiff abscissa     c = %s\n",buf);CHKERRQ(ierr);
1254   }
1255   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
1256   ierr = TSAdaptView(adapt,viewer);CHKERRQ(ierr);
1257   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
1258   PetscFunctionReturn(0);
1259 }
1260 
1261 #undef __FUNCT__
1262 #define __FUNCT__ "TSLoad_ARKIMEX"
1263 static PetscErrorCode TSLoad_ARKIMEX(TS ts,PetscViewer viewer)
1264 {
1265   PetscErrorCode ierr;
1266   SNES           snes;
1267   TSAdapt        tsadapt;
1268 
1269   PetscFunctionBegin;
1270   ierr = TSGetAdapt(ts,&tsadapt);CHKERRQ(ierr);
1271   ierr = TSAdaptLoad(tsadapt,viewer);CHKERRQ(ierr);
1272   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1273   ierr = SNESLoad(snes,viewer);CHKERRQ(ierr);
1274   /* function and Jacobian context for SNES when used with TS is always ts object */
1275   ierr = SNESSetFunction(snes,NULL,NULL,ts);CHKERRQ(ierr);
1276   ierr = SNESSetJacobian(snes,NULL,NULL,NULL,ts);CHKERRQ(ierr);
1277   PetscFunctionReturn(0);
1278 }
1279 
1280 #undef __FUNCT__
1281 #define __FUNCT__ "TSARKIMEXSetType"
1282 /*@C
1283   TSARKIMEXSetType - Set the type of ARK IMEX scheme
1284 
1285   Logically collective
1286 
1287   Input Parameter:
1288 +  ts - timestepping context
1289 -  arktype - type of ARK-IMEX scheme
1290 
1291   Level: intermediate
1292 
1293 .seealso: TSARKIMEXGetType(), TSARKIMEX, TSARKIMEX2D, TSARKIMEX2E, TSARKIMEXPRSSP2, TSARKIMEX3, TSARKIMEXBPR3, TSARKIMEXARS443, TSARKIMEX4, TSARKIMEX5
1294 @*/
1295 PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType arktype)
1296 {
1297   PetscErrorCode ierr;
1298 
1299   PetscFunctionBegin;
1300   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1301   ierr = PetscTryMethod(ts,"TSARKIMEXSetType_C",(TS,TSARKIMEXType),(ts,arktype));CHKERRQ(ierr);
1302   PetscFunctionReturn(0);
1303 }
1304 
1305 #undef __FUNCT__
1306 #define __FUNCT__ "TSARKIMEXGetType"
1307 /*@C
1308   TSARKIMEXGetType - Get the type of ARK IMEX scheme
1309 
1310   Logically collective
1311 
1312   Input Parameter:
1313 .  ts - timestepping context
1314 
1315   Output Parameter:
1316 .  arktype - type of ARK-IMEX scheme
1317 
1318   Level: intermediate
1319 
1320 .seealso: TSARKIMEXGetType()
1321 @*/
1322 PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType *arktype)
1323 {
1324   PetscErrorCode ierr;
1325 
1326   PetscFunctionBegin;
1327   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1328   ierr = PetscUseMethod(ts,"TSARKIMEXGetType_C",(TS,TSARKIMEXType*),(ts,arktype));CHKERRQ(ierr);
1329   PetscFunctionReturn(0);
1330 }
1331 
1332 #undef __FUNCT__
1333 #define __FUNCT__ "TSARKIMEXSetFullyImplicit"
1334 /*@C
1335   TSARKIMEXSetFullyImplicit - Solve both parts of the equation implicitly
1336 
1337   Logically collective
1338 
1339   Input Parameter:
1340 +  ts - timestepping context
1341 -  flg - PETSC_TRUE for fully implicit
1342 
1343   Level: intermediate
1344 
1345 .seealso: TSARKIMEXGetType()
1346 @*/
1347 PetscErrorCode TSARKIMEXSetFullyImplicit(TS ts,PetscBool flg)
1348 {
1349   PetscErrorCode ierr;
1350 
1351   PetscFunctionBegin;
1352   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1353   ierr = PetscTryMethod(ts,"TSARKIMEXSetFullyImplicit_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
1354   PetscFunctionReturn(0);
1355 }
1356 
1357 #undef __FUNCT__
1358 #define __FUNCT__ "TSARKIMEXGetType_ARKIMEX"
1359 PetscErrorCode  TSARKIMEXGetType_ARKIMEX(TS ts,TSARKIMEXType *arktype)
1360 {
1361   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1362   PetscErrorCode ierr;
1363 
1364   PetscFunctionBegin;
1365   if (!ark->tableau) {
1366     ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr);
1367   }
1368   *arktype = ark->tableau->name;
1369   PetscFunctionReturn(0);
1370 }
1371 #undef __FUNCT__
1372 #define __FUNCT__ "TSARKIMEXSetType_ARKIMEX"
1373 PetscErrorCode  TSARKIMEXSetType_ARKIMEX(TS ts,TSARKIMEXType arktype)
1374 {
1375   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1376   PetscErrorCode ierr;
1377   PetscBool      match;
1378   ARKTableauLink link;
1379 
1380   PetscFunctionBegin;
1381   if (ark->tableau) {
1382     ierr = PetscStrcmp(ark->tableau->name,arktype,&match);CHKERRQ(ierr);
1383     if (match) PetscFunctionReturn(0);
1384   }
1385   for (link = ARKTableauList; link; link=link->next) {
1386     ierr = PetscStrcmp(link->tab.name,arktype,&match);CHKERRQ(ierr);
1387     if (match) {
1388       ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr);
1389       ark->tableau = &link->tab;
1390       PetscFunctionReturn(0);
1391     }
1392   }
1393   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",arktype);
1394   PetscFunctionReturn(0);
1395 }
1396 #undef __FUNCT__
1397 #define __FUNCT__ "TSARKIMEXSetFullyImplicit_ARKIMEX"
1398 PetscErrorCode  TSARKIMEXSetFullyImplicit_ARKIMEX(TS ts,PetscBool flg)
1399 {
1400   TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
1401 
1402   PetscFunctionBegin;
1403   ark->imex = (PetscBool)!flg;
1404   PetscFunctionReturn(0);
1405 }
1406 
1407 /* ------------------------------------------------------------ */
1408 /*MC
1409       TSARKIMEX - ODE and DAE solver using Additive Runge-Kutta IMEX schemes
1410 
1411   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
1412   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
1413   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
1414 
1415   Notes:
1416   The default is TSARKIMEX3, it can be changed with TSARKIMEXSetType() or -ts_arkimex_type
1417 
1418   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).
1419 
1420   Level: beginner
1421 
1422 .seealso:  TSCreate(), TS, TSSetType(), TSARKIMEXSetType(), TSARKIMEXGetType(), TSARKIMEXSetFullyImplicit(), TSARKIMEX2D, TTSARKIMEX2E, TSARKIMEX3,
1423            TSARKIMEX4, TSARKIMEX5, TSARKIMEXPRSSP2, TSARKIMEXBPR3, TSARKIMEXType, TSARKIMEXRegister()
1424 
1425 M*/
1426 #undef __FUNCT__
1427 #define __FUNCT__ "TSCreate_ARKIMEX"
1428 PETSC_EXTERN PetscErrorCode TSCreate_ARKIMEX(TS ts)
1429 {
1430   TS_ARKIMEX     *th;
1431   PetscErrorCode ierr;
1432 
1433   PetscFunctionBegin;
1434   ierr = TSARKIMEXInitializePackage();CHKERRQ(ierr);
1435 
1436   ts->ops->reset          = TSReset_ARKIMEX;
1437   ts->ops->destroy        = TSDestroy_ARKIMEX;
1438   ts->ops->view           = TSView_ARKIMEX;
1439   ts->ops->load           = TSLoad_ARKIMEX;
1440   ts->ops->setup          = TSSetUp_ARKIMEX;
1441   ts->ops->step           = TSStep_ARKIMEX;
1442   ts->ops->interpolate    = TSInterpolate_ARKIMEX;
1443   ts->ops->evaluatestep   = TSEvaluateStep_ARKIMEX;
1444   ts->ops->rollback       = TSRollBack_ARKIMEX;
1445   ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX;
1446   ts->ops->snesfunction   = SNESTSFormFunction_ARKIMEX;
1447   ts->ops->snesjacobian   = SNESTSFormJacobian_ARKIMEX;
1448 
1449   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
1450   ts->data = (void*)th;
1451   th->imex = PETSC_TRUE;
1452 
1453   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetType_C",TSARKIMEXGetType_ARKIMEX);CHKERRQ(ierr);
1454   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetType_C",TSARKIMEXSetType_ARKIMEX);CHKERRQ(ierr);
1455   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C",TSARKIMEXSetFullyImplicit_ARKIMEX);CHKERRQ(ierr);
1456   PetscFunctionReturn(0);
1457 }
1458