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