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