xref: /petsc/src/ts/impls/rosw/rosw.c (revision daa037dfd3c3bec8dc8659548d2b20b07c1dc6de)
1 /*
2   Code for timestepping with Rosenbrock W methods
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   This method is designed to be linearly implicit on F and can use an approximate and lagged Jacobian.
11 
12 */
13 #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
14 #include <petscdm.h>
15 
16 #include <petsc/private/kernels/blockinvert.h>
17 
18 static TSRosWType TSRosWDefault = TSROSWRA34PW2;
19 static PetscBool  TSRosWRegisterAllCalled;
20 static PetscBool  TSRosWPackageInitialized;
21 
22 typedef struct _RosWTableau *RosWTableau;
23 struct _RosWTableau {
24   char      *name;
25   PetscInt  order;              /* Classical approximation order of the method */
26   PetscInt  s;                  /* Number of stages */
27   PetscInt  pinterp;            /* Interpolation order */
28   PetscReal *A;                 /* Propagation table, strictly lower triangular */
29   PetscReal *Gamma;             /* Stage table, lower triangular with nonzero diagonal */
30   PetscBool *GammaZeroDiag;     /* Diagonal entries that are zero in stage table Gamma, vector indicating explicit statages */
31   PetscReal *GammaExplicitCorr; /* Coefficients for correction terms needed for explicit stages in transformed variables*/
32   PetscReal *b;                 /* Step completion table */
33   PetscReal *bembed;            /* Step completion table for embedded method of order one less */
34   PetscReal *ASum;              /* Row sum of A */
35   PetscReal *GammaSum;          /* Row sum of Gamma, only needed for non-autonomous systems */
36   PetscReal *At;                /* Propagation table in transformed variables */
37   PetscReal *bt;                /* Step completion table in transformed variables */
38   PetscReal *bembedt;           /* Step completion table of order one less in transformed variables */
39   PetscReal *GammaInv;          /* Inverse of Gamma, used for transformed variables */
40   PetscReal ccfl;               /* Placeholder for CFL coefficient relative to forward Euler */
41   PetscReal *binterpt;          /* Dense output formula */
42 };
43 typedef struct _RosWTableauLink *RosWTableauLink;
44 struct _RosWTableauLink {
45   struct _RosWTableau tab;
46   RosWTableauLink     next;
47 };
48 static RosWTableauLink RosWTableauList;
49 
50 typedef struct {
51   RosWTableau  tableau;
52   Vec          *Y;               /* States computed during the step, used to complete the step */
53   Vec          Ydot;             /* Work vector holding Ydot during residual evaluation */
54   Vec          Ystage;           /* Work vector for the state value at each stage */
55   Vec          Zdot;             /* Ydot = Zdot + shift*Y */
56   Vec          Zstage;           /* Y = Zstage + Y */
57   Vec          vec_sol_prev;     /* Solution from the previous step (used for interpolation and rollback)*/
58   PetscScalar  *work;            /* Scalar work space of length number of stages, used to prepare VecMAXPY() */
59   PetscReal    scoeff;           /* shift = scoeff/dt */
60   PetscReal    stage_time;
61   PetscReal    stage_explicit;     /* Flag indicates that the current stage is explicit */
62   PetscBool    recompute_jacobian; /* Recompute the Jacobian at each stage, default is to freeze the Jacobian at the start of each step */
63   TSStepStatus status;
64 } TS_RosW;
65 
66 /*MC
67      TSROSWTHETA1 - One stage first order L-stable Rosenbrock-W scheme (aka theta method).
68 
69      Only an approximate Jacobian is needed.
70 
71      Level: intermediate
72 
73 .seealso: TSROSW
74 M*/
75 
76 /*MC
77      TSROSWTHETA2 - One stage second order A-stable Rosenbrock-W scheme (aka theta method).
78 
79      Only an approximate Jacobian is needed.
80 
81      Level: intermediate
82 
83 .seealso: TSROSW
84 M*/
85 
86 /*MC
87      TSROSW2M - Two stage second order L-stable Rosenbrock-W scheme.
88 
89      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2P.
90 
91      Level: intermediate
92 
93 .seealso: TSROSW
94 M*/
95 
96 /*MC
97      TSROSW2P - Two stage second order L-stable Rosenbrock-W scheme.
98 
99      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2M.
100 
101      Level: intermediate
102 
103 .seealso: TSROSW
104 M*/
105 
106 /*MC
107      TSROSWRA3PW - Three stage third order Rosenbrock-W scheme for PDAE of index 1.
108 
109      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.
110 
111      This is strongly A-stable with R(infty) = 0.73. The embedded method of order 2 is strongly A-stable with R(infty) = 0.73.
112 
113      References:
114 .  * - Rang and Angermann, New Rosenbrock W methods of order 3 for partial differential algebraic equations of index 1, 2005.
115 
116      Level: intermediate
117 
118 .seealso: TSROSW
119 M*/
120 
121 /*MC
122      TSROSWRA34PW2 - Four stage third order L-stable Rosenbrock-W scheme for PDAE of index 1.
123 
124      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.
125 
126      This is strongly A-stable with R(infty) = 0. The embedded method of order 2 is strongly A-stable with R(infty) = 0.48.
127 
128      References:
129 .  * - Rang and Angermann, New Rosenbrock W methods of order 3 for partial differential algebraic equations of index 1, 2005.
130 
131      Level: intermediate
132 
133 .seealso: TSROSW
134 M*/
135 
136 /*MC
137      TSROSWRODAS3 - Four stage third order L-stable Rosenbrock scheme
138 
139      By default, the Jacobian is only recomputed once per step.
140 
141      Both the third order and embedded second order methods are stiffly accurate and L-stable.
142 
143      References:
144 .  * - Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997.
145 
146      Level: intermediate
147 
148 .seealso: TSROSW, TSROSWSANDU3
149 M*/
150 
151 /*MC
152      TSROSWSANDU3 - Three stage third order L-stable Rosenbrock scheme
153 
154      By default, the Jacobian is only recomputed once per step.
155 
156      The third order method is L-stable, but not stiffly accurate.
157      The second order embedded method is strongly A-stable with R(infty) = 0.5.
158      The internal stages are L-stable.
159      This method is called ROS3 in the paper.
160 
161      References:
162 .  * - Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997.
163 
164      Level: intermediate
165 
166 .seealso: TSROSW, TSROSWRODAS3
167 M*/
168 
169 /*MC
170      TSROSWASSP3P3S1C - A-stable Rosenbrock-W method with SSP explicit part, third order, three stages
171 
172      By default, the Jacobian is only recomputed once per step.
173 
174      A-stable SPP explicit order 3, 3 stages, CFL 1 (eff = 1/3)
175 
176      References:
177 . * - Emil Constantinescu
178 
179      Level: intermediate
180 
181 .seealso: TSROSW, TSROSWLASSP3P4S2C, TSROSWLLSSP3P4S2C, SSP
182 M*/
183 
184 /*MC
185      TSROSWLASSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, four stages
186 
187      By default, the Jacobian is only recomputed once per step.
188 
189      L-stable (A-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2)
190 
191      References:
192 . * - Emil Constantinescu
193 
194      Level: intermediate
195 
196 .seealso: TSROSW, TSROSWASSP3P3S1C, TSROSWLLSSP3P4S2C, TSSSP
197 M*/
198 
199 /*MC
200      TSROSWLLSSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, four stages
201 
202      By default, the Jacobian is only recomputed once per step.
203 
204      L-stable (L-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2)
205 
206      References:
207 . * - Emil Constantinescu
208 
209      Level: intermediate
210 
211 .seealso: TSROSW, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSSSP
212 M*/
213 
214 /*MC
215      TSROSWGRK4T - four stage, fourth order Rosenbrock (not W) method from Kaps and Rentrop
216 
217      By default, the Jacobian is only recomputed once per step.
218 
219      A(89.3 degrees)-stable, |R(infty)| = 0.454.
220 
221      This method does not provide a dense output formula.
222 
223      References:
224 +   * -  Kaps and Rentrop, Generalized Runge Kutta methods of order four with stepsize control for stiff ordinary differential equations, 1979.
225 -   * -  Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2.
226 
227      Hairer's code ros4.f
228 
229      Level: intermediate
230 
231 .seealso: TSROSW, TSROSWSHAMP4, TSROSWVELDD4, TSROSW4L
232 M*/
233 
234 /*MC
235      TSROSWSHAMP4 - four stage, fourth order Rosenbrock (not W) method from Shampine
236 
237      By default, the Jacobian is only recomputed once per step.
238 
239      A-stable, |R(infty)| = 1/3.
240 
241      This method does not provide a dense output formula.
242 
243      References:
244 +   * -  Shampine, Implementation of Rosenbrock methods, 1982.
245 -   * -  Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2.
246 
247      Hairer's code ros4.f
248 
249      Level: intermediate
250 
251 .seealso: TSROSW, TSROSWGRK4T, TSROSWVELDD4, TSROSW4L
252 M*/
253 
254 /*MC
255      TSROSWVELDD4 - four stage, fourth order Rosenbrock (not W) method from van Veldhuizen
256 
257      By default, the Jacobian is only recomputed once per step.
258 
259      A(89.5 degrees)-stable, |R(infty)| = 0.24.
260 
261      This method does not provide a dense output formula.
262 
263      References:
264 +   * -  van Veldhuizen, D stability and Kaps Rentrop methods, 1984.
265 -   * -  Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2.
266 
267      Hairer's code ros4.f
268 
269      Level: intermediate
270 
271 .seealso: TSROSW, TSROSWGRK4T, TSROSWSHAMP4, TSROSW4L
272 M*/
273 
274 /*MC
275      TSROSW4L - four stage, fourth order Rosenbrock (not W) method
276 
277      By default, the Jacobian is only recomputed once per step.
278 
279      A-stable and L-stable
280 
281      This method does not provide a dense output formula.
282 
283      References:
284 .  * -   Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2.
285 
286      Hairer's code ros4.f
287 
288      Level: intermediate
289 
290 .seealso: TSROSW, TSROSWGRK4T, TSROSWSHAMP4, TSROSW4L
291 M*/
292 
293 /*@C
294   TSRosWRegisterAll - Registers all of the Rosenbrock-W methods in TSRosW
295 
296   Not Collective, but should be called by all processes which will need the schemes to be registered
297 
298   Level: advanced
299 
300 .seealso:  TSRosWRegisterDestroy()
301 @*/
302 PetscErrorCode TSRosWRegisterAll(void)
303 {
304   PetscFunctionBegin;
305   if (TSRosWRegisterAllCalled) PetscFunctionReturn(0);
306   TSRosWRegisterAllCalled = PETSC_TRUE;
307 
308   {
309     const PetscReal A = 0;
310     const PetscReal Gamma = 1;
311     const PetscReal b = 1;
312     const PetscReal binterpt=1;
313 
314     PetscCall(TSRosWRegister(TSROSWTHETA1,1,1,&A,&Gamma,&b,NULL,1,&binterpt));
315   }
316 
317   {
318     const PetscReal A = 0;
319     const PetscReal Gamma = 0.5;
320     const PetscReal b = 1;
321     const PetscReal binterpt=1;
322 
323     PetscCall(TSRosWRegister(TSROSWTHETA2,2,1,&A,&Gamma,&b,NULL,1,&binterpt));
324   }
325 
326   {
327     /*const PetscReal g = 1. + 1./PetscSqrtReal(2.0);   Direct evaluation: 1.707106781186547524401. Used for setting up arrays of values known at compile time below. */
328     const PetscReal
329       A[2][2]     = {{0,0}, {1.,0}},
330       Gamma[2][2] = {{1.707106781186547524401,0}, {-2.*1.707106781186547524401,1.707106781186547524401}},
331       b[2]        = {0.5,0.5},
332       b1[2]       = {1.0,0.0};
333     PetscReal binterpt[2][2];
334     binterpt[0][0] = 1.707106781186547524401 - 1.0;
335     binterpt[1][0] = 2.0 - 1.707106781186547524401;
336     binterpt[0][1] = 1.707106781186547524401 - 1.5;
337     binterpt[1][1] = 1.5 - 1.707106781186547524401;
338 
339     PetscCall(TSRosWRegister(TSROSW2P,2,2,&A[0][0],&Gamma[0][0],b,b1,2,&binterpt[0][0]));
340   }
341   {
342     /*const PetscReal g = 1. - 1./PetscSqrtReal(2.0);   Direct evaluation: 0.2928932188134524755992. Used for setting up arrays of values known at compile time below. */
343     const PetscReal
344       A[2][2]     = {{0,0}, {1.,0}},
345       Gamma[2][2] = {{0.2928932188134524755992,0}, {-2.*0.2928932188134524755992,0.2928932188134524755992}},
346       b[2]        = {0.5,0.5},
347       b1[2]       = {1.0,0.0};
348     PetscReal binterpt[2][2];
349     binterpt[0][0] = 0.2928932188134524755992 - 1.0;
350     binterpt[1][0] = 2.0 - 0.2928932188134524755992;
351     binterpt[0][1] = 0.2928932188134524755992 - 1.5;
352     binterpt[1][1] = 1.5 - 0.2928932188134524755992;
353 
354     PetscCall(TSRosWRegister(TSROSW2M,2,2,&A[0][0],&Gamma[0][0],b,b1,2,&binterpt[0][0]));
355   }
356   {
357     /*const PetscReal g = 7.8867513459481287e-01; Directly written in-place below */
358     PetscReal binterpt[3][2];
359     const PetscReal
360       A[3][3] = {{0,0,0},
361                  {1.5773502691896257e+00,0,0},
362                  {0.5,0,0}},
363       Gamma[3][3] = {{7.8867513459481287e-01,0,0},
364                      {-1.5773502691896257e+00,7.8867513459481287e-01,0},
365                      {-6.7075317547305480e-01,-1.7075317547305482e-01,7.8867513459481287e-01}},
366       b[3]  = {1.0566243270259355e-01,4.9038105676657971e-02,8.4529946162074843e-01},
367       b2[3] = {-1.7863279495408180e-01,1./3.,8.4529946162074843e-01};
368 
369       binterpt[0][0] = -0.8094010767585034;
370       binterpt[1][0] = -0.5;
371       binterpt[2][0] = 2.3094010767585034;
372       binterpt[0][1] = 0.9641016151377548;
373       binterpt[1][1] = 0.5;
374       binterpt[2][1] = -1.4641016151377548;
375 
376       PetscCall(TSRosWRegister(TSROSWRA3PW,3,3,&A[0][0],&Gamma[0][0],b,b2,2,&binterpt[0][0]));
377   }
378   {
379     PetscReal  binterpt[4][3];
380     /*const PetscReal g = 4.3586652150845900e-01; Directly written in-place below */
381     const PetscReal
382       A[4][4] = {{0,0,0,0},
383                  {8.7173304301691801e-01,0,0,0},
384                  {8.4457060015369423e-01,-1.1299064236484185e-01,0,0},
385                  {0,0,1.,0}},
386       Gamma[4][4] = {{4.3586652150845900e-01,0,0,0},
387                      {-8.7173304301691801e-01,4.3586652150845900e-01,0,0},
388                      {-9.0338057013044082e-01,5.4180672388095326e-02,4.3586652150845900e-01,0},
389                      {2.4212380706095346e-01,-1.2232505839045147e+00,5.4526025533510214e-01,4.3586652150845900e-01}},
390       b[4]  = {2.4212380706095346e-01,-1.2232505839045147e+00,1.5452602553351020e+00,4.3586652150845900e-01},
391       b2[4] = {3.7810903145819369e-01,-9.6042292212423178e-02,5.0000000000000000e-01,2.1793326075422950e-01};
392 
393     binterpt[0][0]=1.0564298455794094;
394     binterpt[1][0]=2.296429974281067;
395     binterpt[2][0]=-1.307599564525376;
396     binterpt[3][0]=-1.045260255335102;
397     binterpt[0][1]=-1.3864882699759573;
398     binterpt[1][1]=-8.262611700275677;
399     binterpt[2][1]=7.250979895056055;
400     binterpt[3][1]=2.398120075195581;
401     binterpt[0][2]=0.5721822314575016;
402     binterpt[1][2]=4.742931142090097;
403     binterpt[2][2]=-4.398120075195578;
404     binterpt[3][2]=-0.9169932983520199;
405 
406     PetscCall(TSRosWRegister(TSROSWRA34PW2,3,4,&A[0][0],&Gamma[0][0],b,b2,3,&binterpt[0][0]));
407   }
408   {
409     /* const PetscReal g = 0.5;       Directly written in-place below */
410     const PetscReal
411       A[4][4] = {{0,0,0,0},
412                  {0,0,0,0},
413                  {1.,0,0,0},
414                  {0.75,-0.25,0.5,0}},
415       Gamma[4][4] = {{0.5,0,0,0},
416                      {1.,0.5,0,0},
417                      {-0.25,-0.25,0.5,0},
418                      {1./12,1./12,-2./3,0.5}},
419       b[4]  = {5./6,-1./6,-1./6,0.5},
420       b2[4] = {0.75,-0.25,0.5,0};
421 
422     PetscCall(TSRosWRegister(TSROSWRODAS3,3,4,&A[0][0],&Gamma[0][0],b,b2,0,NULL));
423   }
424   {
425     /*const PetscReal g = 0.43586652150845899941601945119356;       Directly written in-place below */
426     const PetscReal
427       A[3][3] = {{0,0,0},
428                  {0.43586652150845899941601945119356,0,0},
429                  {0.43586652150845899941601945119356,0,0}},
430       Gamma[3][3] = {{0.43586652150845899941601945119356,0,0},
431                      {-0.19294655696029095575009695436041,0.43586652150845899941601945119356,0},
432                      {0,1.74927148125794685173529749738960,0.43586652150845899941601945119356}},
433       b[3]  = {-0.75457412385404315829818998646589,1.94100407061964420292840123379419,-0.18642994676560104463021124732829},
434       b2[3] = {-1.53358745784149585370766523913002,2.81745131148625772213931745457622,-0.28386385364476186843165221544619};
435 
436     PetscReal binterpt[3][2];
437     binterpt[0][0] = 3.793692883777660870425141387941;
438     binterpt[1][0] = -2.918692883777660870425141387941;
439     binterpt[2][0] = 0.125;
440     binterpt[0][1] = -0.725741064379812106687651020584;
441     binterpt[1][1] = 0.559074397713145440020984353917;
442     binterpt[2][1] = 0.16666666666666666666666666666667;
443 
444     PetscCall(TSRosWRegister(TSROSWSANDU3,3,3,&A[0][0],&Gamma[0][0],b,b2,2,&binterpt[0][0]));
445   }
446   {
447     /*const PetscReal s3 = PetscSqrtReal(3.),g = (3.0+s3)/6.0;
448      * Direct evaluation: s3 = 1.732050807568877293527;
449      *                     g = 0.7886751345948128822546;
450      * Values are directly inserted below to ensure availability at compile time (compiler warnings otherwise...) */
451     const PetscReal
452       A[3][3] = {{0,0,0},
453                  {1,0,0},
454                  {0.25,0.25,0}},
455       Gamma[3][3] = {{0,0,0},
456                      {(-3.0-1.732050807568877293527)/6.0,0.7886751345948128822546,0},
457                      {(-3.0-1.732050807568877293527)/24.0,(-3.0-1.732050807568877293527)/8.0,0.7886751345948128822546}},
458       b[3]  = {1./6.,1./6.,2./3.},
459       b2[3] = {1./4.,1./4.,1./2.};
460     PetscReal binterpt[3][2];
461 
462     binterpt[0][0]=0.089316397477040902157517886164709;
463     binterpt[1][0]=-0.91068360252295909784248211383529;
464     binterpt[2][0]=1.8213672050459181956849642276706;
465     binterpt[0][1]=0.077350269189625764509148780501957;
466     binterpt[1][1]=1.077350269189625764509148780502;
467     binterpt[2][1]=-1.1547005383792515290182975610039;
468 
469     PetscCall(TSRosWRegister(TSROSWASSP3P3S1C,3,3,&A[0][0],&Gamma[0][0],b,b2,2,&binterpt[0][0]));
470   }
471 
472   {
473     const PetscReal
474       A[4][4] = {{0,0,0,0},
475                  {1./2.,0,0,0},
476                  {1./2.,1./2.,0,0},
477                  {1./6.,1./6.,1./6.,0}},
478       Gamma[4][4] = {{1./2.,0,0,0},
479                      {0.0,1./4.,0,0},
480                      {-2.,-2./3.,2./3.,0},
481                      {1./2.,5./36.,-2./9,0}},
482       b[4]  = {1./6.,1./6.,1./6.,1./2.},
483       b2[4] = {1./8.,3./4.,1./8.,0};
484     PetscReal binterpt[4][3];
485 
486     binterpt[0][0]=6.25;
487     binterpt[1][0]=-30.25;
488     binterpt[2][0]=1.75;
489     binterpt[3][0]=23.25;
490     binterpt[0][1]=-9.75;
491     binterpt[1][1]=58.75;
492     binterpt[2][1]=-3.25;
493     binterpt[3][1]=-45.75;
494     binterpt[0][2]=3.6666666666666666666666666666667;
495     binterpt[1][2]=-28.333333333333333333333333333333;
496     binterpt[2][2]=1.6666666666666666666666666666667;
497     binterpt[3][2]=23.;
498 
499     PetscCall(TSRosWRegister(TSROSWLASSP3P4S2C,3,4,&A[0][0],&Gamma[0][0],b,b2,3,&binterpt[0][0]));
500   }
501 
502   {
503     const PetscReal
504       A[4][4] = {{0,0,0,0},
505                  {1./2.,0,0,0},
506                  {1./2.,1./2.,0,0},
507                  {1./6.,1./6.,1./6.,0}},
508       Gamma[4][4] = {{1./2.,0,0,0},
509                      {0.0,3./4.,0,0},
510                      {-2./3.,-23./9.,2./9.,0},
511                      {1./18.,65./108.,-2./27,0}},
512       b[4]  = {1./6.,1./6.,1./6.,1./2.},
513       b2[4] = {3./16.,10./16.,3./16.,0};
514     PetscReal binterpt[4][3];
515 
516     binterpt[0][0]=1.6911764705882352941176470588235;
517     binterpt[1][0]=3.6813725490196078431372549019608;
518     binterpt[2][0]=0.23039215686274509803921568627451;
519     binterpt[3][0]=-4.6029411764705882352941176470588;
520     binterpt[0][1]=-0.95588235294117647058823529411765;
521     binterpt[1][1]=-6.2401960784313725490196078431373;
522     binterpt[2][1]=-0.31862745098039215686274509803922;
523     binterpt[3][1]=7.5147058823529411764705882352941;
524     binterpt[0][2]=-0.56862745098039215686274509803922;
525     binterpt[1][2]=2.7254901960784313725490196078431;
526     binterpt[2][2]=0.25490196078431372549019607843137;
527     binterpt[3][2]=-2.4117647058823529411764705882353;
528 
529     PetscCall(TSRosWRegister(TSROSWLLSSP3P4S2C,3,4,&A[0][0],&Gamma[0][0],b,b2,3,&binterpt[0][0]));
530   }
531 
532   {
533     PetscReal A[4][4],Gamma[4][4],b[4],b2[4];
534     PetscReal binterpt[4][3];
535 
536     Gamma[0][0]=0.4358665215084589994160194475295062513822671686978816;
537     Gamma[0][1]=0; Gamma[0][2]=0; Gamma[0][3]=0;
538     Gamma[1][0]=-1.997527830934941248426324674704153457289527280554476;
539     Gamma[1][1]=0.4358665215084589994160194475295062513822671686978816;
540     Gamma[1][2]=0; Gamma[1][3]=0;
541     Gamma[2][0]=-1.007948511795029620852002345345404191008352770119903;
542     Gamma[2][1]=-0.004648958462629345562774289390054679806993396798458131;
543     Gamma[2][2]=0.4358665215084589994160194475295062513822671686978816;
544     Gamma[2][3]=0;
545     Gamma[3][0]=-0.6685429734233467180451604600279552604364311322650783;
546     Gamma[3][1]=0.6056625986449338476089525334450053439525178740492984;
547     Gamma[3][2]=-0.9717899277217721234705114616271378792182450260943198;
548     Gamma[3][3]=0;
549 
550     A[0][0]=0; A[0][1]=0; A[0][2]=0; A[0][3]=0;
551     A[1][0]=0.8717330430169179988320388950590125027645343373957631;
552     A[1][1]=0; A[1][2]=0; A[1][3]=0;
553     A[2][0]=0.5275890119763004115618079766722914408876108660811028;
554     A[2][1]=0.07241098802369958843819203208518599088698057726988732;
555     A[2][2]=0; A[2][3]=0;
556     A[3][0]=0.3990960076760701320627260685975778145384666450351314;
557     A[3][1]=-0.4375576546135194437228463747348862825846903771419953;
558     A[3][2]=1.038461646937449311660120300601880176655352737312713;
559     A[3][3]=0;
560 
561     b[0]=0.1876410243467238251612921333138006734899663569186926;
562     b[1]=-0.5952974735769549480478230473706443582188442040780541;
563     b[2]=0.9717899277217721234705114616271378792182450260943198;
564     b[3]=0.4358665215084589994160194475295062513822671686978816;
565 
566     b2[0]=0.2147402862233891404862383521089097657790734483804460;
567     b2[1]=-0.4851622638849390928209050538171743017757490232519684;
568     b2[2]=0.8687250025203875511662123688667549217531982787600080;
569     b2[3]=0.4016969751411624011684543450940068201770721128357014;
570 
571     binterpt[0][0]=2.2565812720167954547104627844105;
572     binterpt[1][0]=1.349166413351089573796243820819;
573     binterpt[2][0]=-2.4695174540533503758652847586647;
574     binterpt[3][0]=-0.13623023131453465264142184656474;
575     binterpt[0][1]=-3.0826699111559187902922463354557;
576     binterpt[1][1]=-2.4689115685996042534544925650515;
577     binterpt[2][1]=5.7428279814696677152129332773553;
578     binterpt[3][1]=-0.19124650171414467146619437684812;
579     binterpt[0][2]=1.0137296634858471607430756831148;
580     binterpt[1][2]=0.52444768167155973161042570784064;
581     binterpt[2][2]=-2.3015205996945452158771370439586;
582     binterpt[3][2]=0.76334325453713832352363565300308;
583 
584     PetscCall(TSRosWRegister(TSROSWARK3,3,4,&A[0][0],&Gamma[0][0],b,b2,3,&binterpt[0][0]));
585   }
586   PetscCall(TSRosWRegisterRos4(TSROSWGRK4T,0.231,PETSC_DEFAULT,PETSC_DEFAULT,0,-0.1282612945269037e+01));
587   PetscCall(TSRosWRegisterRos4(TSROSWSHAMP4,0.5,PETSC_DEFAULT,PETSC_DEFAULT,0,125./108.));
588   PetscCall(TSRosWRegisterRos4(TSROSWVELDD4,0.22570811482256823492,PETSC_DEFAULT,PETSC_DEFAULT,0,-1.355958941201148));
589   PetscCall(TSRosWRegisterRos4(TSROSW4L,0.57282,PETSC_DEFAULT,PETSC_DEFAULT,0,-1.093502252409163));
590   PetscFunctionReturn(0);
591 }
592 
593 /*@C
594    TSRosWRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister().
595 
596    Not Collective
597 
598    Level: advanced
599 
600 .seealso: TSRosWRegister(), TSRosWRegisterAll()
601 @*/
602 PetscErrorCode TSRosWRegisterDestroy(void)
603 {
604   RosWTableauLink link;
605 
606   PetscFunctionBegin;
607   while ((link = RosWTableauList)) {
608     RosWTableau t = &link->tab;
609     RosWTableauList = link->next;
610     PetscCall(PetscFree5(t->A,t->Gamma,t->b,t->ASum,t->GammaSum));
611     PetscCall(PetscFree5(t->At,t->bt,t->GammaInv,t->GammaZeroDiag,t->GammaExplicitCorr));
612     PetscCall(PetscFree2(t->bembed,t->bembedt));
613     PetscCall(PetscFree(t->binterpt));
614     PetscCall(PetscFree(t->name));
615     PetscCall(PetscFree(link));
616   }
617   TSRosWRegisterAllCalled = PETSC_FALSE;
618   PetscFunctionReturn(0);
619 }
620 
621 /*@C
622   TSRosWInitializePackage - This function initializes everything in the TSRosW package. It is called
623   from TSInitializePackage().
624 
625   Level: developer
626 
627 .seealso: PetscInitialize()
628 @*/
629 PetscErrorCode TSRosWInitializePackage(void)
630 {
631   PetscFunctionBegin;
632   if (TSRosWPackageInitialized) PetscFunctionReturn(0);
633   TSRosWPackageInitialized = PETSC_TRUE;
634   PetscCall(TSRosWRegisterAll());
635   PetscCall(PetscRegisterFinalize(TSRosWFinalizePackage));
636   PetscFunctionReturn(0);
637 }
638 
639 /*@C
640   TSRosWFinalizePackage - This function destroys everything in the TSRosW package. It is
641   called from PetscFinalize().
642 
643   Level: developer
644 
645 .seealso: PetscFinalize()
646 @*/
647 PetscErrorCode TSRosWFinalizePackage(void)
648 {
649   PetscFunctionBegin;
650   TSRosWPackageInitialized = PETSC_FALSE;
651   PetscCall(TSRosWRegisterDestroy());
652   PetscFunctionReturn(0);
653 }
654 
655 /*@C
656    TSRosWRegister - register a Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
657 
658    Not Collective, but the same schemes should be registered on all processes on which they will be used
659 
660    Input Parameters:
661 +  name - identifier for method
662 .  order - approximation order of method
663 .  s - number of stages, this is the dimension of the matrices below
664 .  A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular
665 .  Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal
666 .  b - Step completion table (dimension s)
667 .  bembed - Step completion table for a scheme of order one less (dimension s, NULL if no embedded scheme is available)
668 .  pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt
669 -  binterpt - Coefficients of the interpolation formula (dimension s*pinterp)
670 
671    Notes:
672    Several Rosenbrock W methods are provided, this function is only needed to create new methods.
673 
674    Level: advanced
675 
676 .seealso: TSRosW
677 @*/
678 PetscErrorCode TSRosWRegister(TSRosWType name,PetscInt order,PetscInt s,const PetscReal A[],const PetscReal Gamma[],const PetscReal b[],const PetscReal bembed[],
679                               PetscInt pinterp,const PetscReal binterpt[])
680 {
681   RosWTableauLink link;
682   RosWTableau     t;
683   PetscInt        i,j,k;
684   PetscScalar     *GammaInv;
685 
686   PetscFunctionBegin;
687   PetscValidCharPointer(name,1);
688   PetscValidRealPointer(A,4);
689   PetscValidRealPointer(Gamma,5);
690   PetscValidRealPointer(b,6);
691   if (bembed) PetscValidRealPointer(bembed,7);
692 
693   PetscCall(TSRosWInitializePackage());
694   PetscCall(PetscNew(&link));
695   t        = &link->tab;
696   PetscCall(PetscStrallocpy(name,&t->name));
697   t->order = order;
698   t->s     = s;
699   PetscCall(PetscMalloc5(s*s,&t->A,s*s,&t->Gamma,s,&t->b,s,&t->ASum,s,&t->GammaSum));
700   PetscCall(PetscMalloc5(s*s,&t->At,s,&t->bt,s*s,&t->GammaInv,s,&t->GammaZeroDiag,s*s,&t->GammaExplicitCorr));
701   PetscCall(PetscArraycpy(t->A,A,s*s));
702   PetscCall(PetscArraycpy(t->Gamma,Gamma,s*s));
703   PetscCall(PetscArraycpy(t->GammaExplicitCorr,Gamma,s*s));
704   PetscCall(PetscArraycpy(t->b,b,s));
705   if (bembed) {
706     PetscCall(PetscMalloc2(s,&t->bembed,s,&t->bembedt));
707     PetscCall(PetscArraycpy(t->bembed,bembed,s));
708   }
709   for (i=0; i<s; i++) {
710     t->ASum[i]     = 0;
711     t->GammaSum[i] = 0;
712     for (j=0; j<s; j++) {
713       t->ASum[i]     += A[i*s+j];
714       t->GammaSum[i] += Gamma[i*s+j];
715     }
716   }
717   PetscCall(PetscMalloc1(s*s,&GammaInv)); /* Need to use Scalar for inverse, then convert back to Real */
718   for (i=0; i<s*s; i++) GammaInv[i] = Gamma[i];
719   for (i=0; i<s; i++) {
720     if (Gamma[i*s+i] == 0.0) {
721       GammaInv[i*s+i] = 1.0;
722       t->GammaZeroDiag[i] = PETSC_TRUE;
723     } else {
724       t->GammaZeroDiag[i] = PETSC_FALSE;
725     }
726   }
727 
728   switch (s) {
729   case 1: GammaInv[0] = 1./GammaInv[0]; break;
730   case 2: PetscCall(PetscKernel_A_gets_inverse_A_2(GammaInv,0,PETSC_FALSE,NULL)); break;
731   case 3: PetscCall(PetscKernel_A_gets_inverse_A_3(GammaInv,0,PETSC_FALSE,NULL)); break;
732   case 4: PetscCall(PetscKernel_A_gets_inverse_A_4(GammaInv,0,PETSC_FALSE,NULL)); break;
733   case 5: {
734     PetscInt  ipvt5[5];
735     MatScalar work5[5*5];
736     PetscCall(PetscKernel_A_gets_inverse_A_5(GammaInv,ipvt5,work5,0,PETSC_FALSE,NULL)); break;
737   }
738   case 6: PetscCall(PetscKernel_A_gets_inverse_A_6(GammaInv,0,PETSC_FALSE,NULL)); break;
739   case 7: PetscCall(PetscKernel_A_gets_inverse_A_7(GammaInv,0,PETSC_FALSE,NULL)); break;
740   default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for %D stages",s);
741   }
742   for (i=0; i<s*s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]);
743   PetscCall(PetscFree(GammaInv));
744 
745   for (i=0; i<s; i++) {
746     for (k=0; k<i+1; k++) {
747       t->GammaExplicitCorr[i*s+k]=(t->GammaExplicitCorr[i*s+k])*(t->GammaInv[k*s+k]);
748       for (j=k+1; j<i+1; j++) {
749         t->GammaExplicitCorr[i*s+k]+=(t->GammaExplicitCorr[i*s+j])*(t->GammaInv[j*s+k]);
750       }
751     }
752   }
753 
754   for (i=0; i<s; i++) {
755     for (j=0; j<s; j++) {
756       t->At[i*s+j] = 0;
757       for (k=0; k<s; k++) {
758         t->At[i*s+j] += t->A[i*s+k] * t->GammaInv[k*s+j];
759       }
760     }
761     t->bt[i] = 0;
762     for (j=0; j<s; j++) {
763       t->bt[i] += t->b[j] * t->GammaInv[j*s+i];
764     }
765     if (bembed) {
766       t->bembedt[i] = 0;
767       for (j=0; j<s; j++) {
768         t->bembedt[i] += t->bembed[j] * t->GammaInv[j*s+i];
769       }
770     }
771   }
772   t->ccfl = 1.0;                /* Fix this */
773 
774   t->pinterp = pinterp;
775   PetscCall(PetscMalloc1(s*pinterp,&t->binterpt));
776   PetscCall(PetscArraycpy(t->binterpt,binterpt,s*pinterp));
777   link->next = RosWTableauList;
778   RosWTableauList = link;
779   PetscFunctionReturn(0);
780 }
781 
782 /*@C
783    TSRosWRegisterRos4 - register a fourth order Rosenbrock scheme by providing parameter choices
784 
785    Not Collective, but the same schemes should be registered on all processes on which they will be used
786 
787    Input Parameters:
788 +  name - identifier for method
789 .  gamma - leading coefficient (diagonal entry)
790 .  a2 - design parameter, see Table 7.2 of Hairer&Wanner
791 .  a3 - design parameter or PETSC_DEFAULT to satisfy one of the order five conditions (Eq 7.22)
792 .  b3 - design parameter, see Table 7.2 of Hairer&Wanner
793 .  beta43 - design parameter or PETSC_DEFAULT to use Equation 7.21 of Hairer&Wanner
794 -  e4 - design parameter for embedded method, see coefficient E4 in ros4.f code from Hairer
795 
796    Notes:
797    This routine encodes the design of fourth order Rosenbrock methods as described in Hairer and Wanner volume 2.
798    It is used here to implement several methods from the book and can be used to experiment with new methods.
799    It was written this way instead of by copying coefficients in order to provide better than double precision satisfaction of the order conditions.
800 
801    Level: developer
802 
803 .seealso: TSRosW, TSRosWRegister()
804 @*/
805 PetscErrorCode TSRosWRegisterRos4(TSRosWType name,PetscReal gamma,PetscReal a2,PetscReal a3,PetscReal b3,PetscReal e4)
806 {
807   /* Declare numeric constants so they can be quad precision without being truncated at double */
808   const PetscReal one = 1,two = 2,three = 3,four = 4,five = 5,six = 6,eight = 8,twelve = 12,twenty = 20,twentyfour = 24,
809     p32 = one/six - gamma + gamma*gamma,
810     p42 = one/eight - gamma/three,
811     p43 = one/twelve - gamma/three,
812     p44 = one/twentyfour - gamma/two + three/two*gamma*gamma - gamma*gamma*gamma,
813     p56 = one/twenty - gamma/four;
814   PetscReal   a4,a32,a42,a43,b1,b2,b4,beta2p,beta3p,beta4p,beta32,beta42,beta43,beta32beta2p,beta4jbetajp;
815   PetscReal   A[4][4],Gamma[4][4],b[4],bm[4];
816   PetscScalar M[3][3],rhs[3];
817 
818   PetscFunctionBegin;
819   /* Step 1: choose Gamma (input) */
820   /* Step 2: choose a2,a3,a4; b1,b2,b3,b4 to satisfy order conditions */
821   if (a3 == PETSC_DEFAULT) a3 = (one/five - a2/four)/(one/four - a2/three); /* Eq 7.22 */
822   a4 = a3;                                                  /* consequence of 7.20 */
823 
824   /* Solve order conditions 7.15a, 7.15c, 7.15e */
825   M[0][0] = one; M[0][1] = one;      M[0][2] = one;      /* 7.15a */
826   M[1][0] = 0.0; M[1][1] = a2*a2;    M[1][2] = a4*a4;    /* 7.15c */
827   M[2][0] = 0.0; M[2][1] = a2*a2*a2; M[2][2] = a4*a4*a4; /* 7.15e */
828   rhs[0]  = one - b3;
829   rhs[1]  = one/three - a3*a3*b3;
830   rhs[2]  = one/four - a3*a3*a3*b3;
831   PetscCall(PetscKernel_A_gets_inverse_A_3(&M[0][0],0,PETSC_FALSE,NULL));
832   b1      = PetscRealPart(M[0][0]*rhs[0] + M[0][1]*rhs[1] + M[0][2]*rhs[2]);
833   b2      = PetscRealPart(M[1][0]*rhs[0] + M[1][1]*rhs[1] + M[1][2]*rhs[2]);
834   b4      = PetscRealPart(M[2][0]*rhs[0] + M[2][1]*rhs[1] + M[2][2]*rhs[2]);
835 
836   /* Step 3 */
837   beta43       = (p56 - a2*p43) / (b4*a3*a3*(a3 - a2)); /* 7.21 */
838   beta32beta2p =  p44 / (b4*beta43);                    /* 7.15h */
839   beta4jbetajp = (p32 - b3*beta32beta2p) / b4;
840   M[0][0]      = b2;                                    M[0][1] = b3;                 M[0][2] = b4;
841   M[1][0]      = a4*a4*beta32beta2p-a3*a3*beta4jbetajp; M[1][1] = a2*a2*beta4jbetajp; M[1][2] = -a2*a2*beta32beta2p;
842   M[2][0]      = b4*beta43*a3*a3-p43;                   M[2][1] = -b4*beta43*a2*a2;   M[2][2] = 0;
843   rhs[0]       = one/two - gamma; rhs[1] = 0; rhs[2] = -a2*a2*p32;
844   PetscCall(PetscKernel_A_gets_inverse_A_3(&M[0][0],0,PETSC_FALSE,NULL));
845   beta2p       = PetscRealPart(M[0][0]*rhs[0] + M[0][1]*rhs[1] + M[0][2]*rhs[2]);
846   beta3p       = PetscRealPart(M[1][0]*rhs[0] + M[1][1]*rhs[1] + M[1][2]*rhs[2]);
847   beta4p       = PetscRealPart(M[2][0]*rhs[0] + M[2][1]*rhs[1] + M[2][2]*rhs[2]);
848 
849   /* Step 4: back-substitute */
850   beta32 = beta32beta2p / beta2p;
851   beta42 = (beta4jbetajp - beta43*beta3p) / beta2p;
852 
853   /* Step 5: 7.15f and 7.20, then 7.16 */
854   a43 = 0;
855   a32 = p42 / (b3*a3*beta2p + b4*a4*beta2p);
856   a42 = a32;
857 
858   A[0][0]     = 0;          A[0][1] = 0;   A[0][2] = 0;   A[0][3] = 0;
859   A[1][0]     = a2;         A[1][1] = 0;   A[1][2] = 0;   A[1][3] = 0;
860   A[2][0]     = a3-a32;     A[2][1] = a32; A[2][2] = 0;   A[2][3] = 0;
861   A[3][0]     = a4-a43-a42; A[3][1] = a42; A[3][2] = a43; A[3][3] = 0;
862   Gamma[0][0] = gamma;                        Gamma[0][1] = 0;              Gamma[0][2] = 0;              Gamma[0][3] = 0;
863   Gamma[1][0] = beta2p-A[1][0];               Gamma[1][1] = gamma;          Gamma[1][2] = 0;              Gamma[1][3] = 0;
864   Gamma[2][0] = beta3p-beta32-A[2][0];        Gamma[2][1] = beta32-A[2][1]; Gamma[2][2] = gamma;          Gamma[2][3] = 0;
865   Gamma[3][0] = beta4p-beta42-beta43-A[3][0]; Gamma[3][1] = beta42-A[3][1]; Gamma[3][2] = beta43-A[3][2]; Gamma[3][3] = gamma;
866   b[0] = b1; b[1] = b2; b[2] = b3; b[3] = b4;
867 
868   /* Construct embedded formula using given e4. We are solving Equation 7.18. */
869   bm[3] = b[3] - e4*gamma;                                          /* using definition of E4 */
870   bm[2] = (p32 - beta4jbetajp*bm[3]) / (beta32*beta2p);             /* fourth row of 7.18 */
871   bm[1] = (one/two - gamma - beta3p*bm[2] - beta4p*bm[3]) / beta2p; /* second row */
872   bm[0] = one - bm[1] - bm[2] - bm[3];                              /* first row */
873 
874   {
875     const PetscReal misfit = a2*a2*bm[1] + a3*a3*bm[2] + a4*a4*bm[3] - one/three;
876     PetscCheck(PetscAbs(misfit) <= PETSC_SMALL,PETSC_COMM_SELF,PETSC_ERR_SUP,"Assumptions violated, could not construct a third order embedded method");
877   }
878   PetscCall(TSRosWRegister(name,4,4,&A[0][0],&Gamma[0][0],b,bm,0,NULL));
879   PetscFunctionReturn(0);
880 }
881 
882 /*
883  The step completion formula is
884 
885  x1 = x0 + b^T Y
886 
887  where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been
888  updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write
889 
890  x1e = x0 + be^T Y
891      = x1 - b^T Y + be^T Y
892      = x1 + (be - b)^T Y
893 
894  so we can evaluate the method of different order even after the step has been optimistically completed.
895 */
896 static PetscErrorCode TSEvaluateStep_RosW(TS ts,PetscInt order,Vec U,PetscBool *done)
897 {
898   TS_RosW        *ros = (TS_RosW*)ts->data;
899   RosWTableau    tab  = ros->tableau;
900   PetscScalar    *w   = ros->work;
901   PetscInt       i;
902 
903   PetscFunctionBegin;
904   if (order == tab->order) {
905     if (ros->status == TS_STEP_INCOMPLETE) { /* Use standard completion formula */
906       PetscCall(VecCopy(ts->vec_sol,U));
907       for (i=0; i<tab->s; i++) w[i] = tab->bt[i];
908       PetscCall(VecMAXPY(U,tab->s,w,ros->Y));
909     } else PetscCall(VecCopy(ts->vec_sol,U));
910     if (done) *done = PETSC_TRUE;
911     PetscFunctionReturn(0);
912   } else if (order == tab->order-1) {
913     if (!tab->bembedt) goto unavailable;
914     if (ros->status == TS_STEP_INCOMPLETE) { /* Use embedded completion formula */
915       PetscCall(VecCopy(ts->vec_sol,U));
916       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i];
917       PetscCall(VecMAXPY(U,tab->s,w,ros->Y));
918     } else {                    /* Use rollback-and-recomplete formula (bembedt - bt) */
919       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i];
920       PetscCall(VecCopy(ts->vec_sol,U));
921       PetscCall(VecMAXPY(U,tab->s,w,ros->Y));
922     }
923     if (done) *done = PETSC_TRUE;
924     PetscFunctionReturn(0);
925   }
926   unavailable:
927   if (done) *done = PETSC_FALSE;
928   else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Rosenbrock-W '%s' of order %D cannot evaluate step at order %D. Consider using -ts_adapt_type none or a different method that has an embedded estimate.",tab->name,tab->order,order);
929   PetscFunctionReturn(0);
930 }
931 
932 static PetscErrorCode TSRollBack_RosW(TS ts)
933 {
934   TS_RosW        *ros = (TS_RosW*)ts->data;
935 
936   PetscFunctionBegin;
937   PetscCall(VecCopy(ros->vec_sol_prev,ts->vec_sol));
938   PetscFunctionReturn(0);
939 }
940 
941 static PetscErrorCode TSStep_RosW(TS ts)
942 {
943   TS_RosW         *ros = (TS_RosW*)ts->data;
944   RosWTableau     tab  = ros->tableau;
945   const PetscInt  s    = tab->s;
946   const PetscReal *At  = tab->At,*Gamma = tab->Gamma,*ASum = tab->ASum,*GammaInv = tab->GammaInv;
947   const PetscReal *GammaExplicitCorr = tab->GammaExplicitCorr;
948   const PetscBool *GammaZeroDiag = tab->GammaZeroDiag;
949   PetscScalar     *w   = ros->work;
950   Vec             *Y   = ros->Y,Ydot = ros->Ydot,Zdot = ros->Zdot,Zstage = ros->Zstage;
951   SNES            snes;
952   TSAdapt         adapt;
953   PetscInt        i,j,its,lits;
954   PetscInt        rejections = 0;
955   PetscBool       stageok,accept = PETSC_TRUE;
956   PetscReal       next_time_step = ts->time_step;
957   PetscInt        lag;
958 
959   PetscFunctionBegin;
960   if (!ts->steprollback) {
961     PetscCall(VecCopy(ts->vec_sol,ros->vec_sol_prev));
962   }
963 
964   ros->status = TS_STEP_INCOMPLETE;
965   while (!ts->reason && ros->status != TS_STEP_COMPLETE) {
966     const PetscReal h = ts->time_step;
967     for (i=0; i<s; i++) {
968       ros->stage_time = ts->ptime + h*ASum[i];
969       PetscCall(TSPreStage(ts,ros->stage_time));
970       if (GammaZeroDiag[i]) {
971         ros->stage_explicit = PETSC_TRUE;
972         ros->scoeff         = 1.;
973       } else {
974         ros->stage_explicit = PETSC_FALSE;
975         ros->scoeff         = 1./Gamma[i*s+i];
976       }
977 
978       PetscCall(VecCopy(ts->vec_sol,Zstage));
979       for (j=0; j<i; j++) w[j] = At[i*s+j];
980       PetscCall(VecMAXPY(Zstage,i,w,Y));
981 
982       for (j=0; j<i; j++) w[j] = 1./h * GammaInv[i*s+j];
983       PetscCall(VecZeroEntries(Zdot));
984       PetscCall(VecMAXPY(Zdot,i,w,Y));
985 
986       /* Initial guess taken from last stage */
987       PetscCall(VecZeroEntries(Y[i]));
988 
989       if (!ros->stage_explicit) {
990         PetscCall(TSGetSNES(ts,&snes));
991         if (!ros->recompute_jacobian && !i) {
992           PetscCall(SNESGetLagJacobian(snes,&lag));
993           if (lag == 1) {  /* use did not set a nontrival lag, so lag over all stages */
994             PetscCall(SNESSetLagJacobian(snes,-2)); /* Recompute the Jacobian on this solve, but not again for the rest of the stages */
995           }
996         }
997         PetscCall(SNESSolve(snes,NULL,Y[i]));
998         if (!ros->recompute_jacobian && i == s-1 && lag == 1) {
999           PetscCall(SNESSetLagJacobian(snes,lag)); /* Set lag back to 1 so we know user did not set it */
1000         }
1001         PetscCall(SNESGetIterationNumber(snes,&its));
1002         PetscCall(SNESGetLinearSolveIterations(snes,&lits));
1003         ts->snes_its += its; ts->ksp_its += lits;
1004       } else {
1005         Mat J,Jp;
1006         PetscCall(VecZeroEntries(Ydot)); /* Evaluate Y[i]=G(t,Ydot=0,Zstage) */
1007         PetscCall(TSComputeIFunction(ts,ros->stage_time,Zstage,Ydot,Y[i],PETSC_FALSE));
1008         PetscCall(VecScale(Y[i],-1.0));
1009         PetscCall(VecAXPY(Y[i],-1.0,Zdot)); /*Y[i] = F(Zstage)-Zdot[=GammaInv*Y]*/
1010 
1011         PetscCall(VecZeroEntries(Zstage)); /* Zstage = GammaExplicitCorr[i,j] * Y[j] */
1012         for (j=0; j<i; j++) w[j] = GammaExplicitCorr[i*s+j];
1013         PetscCall(VecMAXPY(Zstage,i,w,Y));
1014 
1015         /* Y[i] = Y[i] + Jac*Zstage[=Jac*GammaExplicitCorr[i,j] * Y[j]] */
1016         PetscCall(TSGetIJacobian(ts,&J,&Jp,NULL,NULL));
1017         PetscCall(TSComputeIJacobian(ts,ros->stage_time,ts->vec_sol,Ydot,0,J,Jp,PETSC_FALSE));
1018         PetscCall(MatMult(J,Zstage,Zdot));
1019         PetscCall(VecAXPY(Y[i],-1.0,Zdot));
1020         ts->ksp_its += 1;
1021 
1022         PetscCall(VecScale(Y[i],h));
1023       }
1024       PetscCall(TSPostStage(ts,ros->stage_time,i,Y));
1025       PetscCall(TSGetAdapt(ts,&adapt));
1026       PetscCall(TSAdaptCheckStage(adapt,ts,ros->stage_time,Y[i],&stageok));
1027       if (!stageok) goto reject_step;
1028     }
1029 
1030     ros->status = TS_STEP_INCOMPLETE;
1031     PetscCall(TSEvaluateStep_RosW(ts,tab->order,ts->vec_sol,NULL));
1032     ros->status = TS_STEP_PENDING;
1033     PetscCall(TSGetAdapt(ts,&adapt));
1034     PetscCall(TSAdaptCandidatesClear(adapt));
1035     PetscCall(TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE));
1036     PetscCall(TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept));
1037     ros->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
1038     if (!accept) { /* Roll back the current step */
1039       PetscCall(TSRollBack_RosW(ts));
1040       ts->time_step = next_time_step;
1041       goto reject_step;
1042     }
1043 
1044     ts->ptime += ts->time_step;
1045     ts->time_step = next_time_step;
1046     break;
1047 
1048   reject_step:
1049     ts->reject++; accept = PETSC_FALSE;
1050     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
1051       ts->reason = TS_DIVERGED_STEP_REJECTED;
1052       PetscCall(PetscInfo(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections));
1053     }
1054   }
1055   PetscFunctionReturn(0);
1056 }
1057 
1058 static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec U)
1059 {
1060   TS_RosW         *ros = (TS_RosW*)ts->data;
1061   PetscInt        s    = ros->tableau->s,pinterp = ros->tableau->pinterp,i,j;
1062   PetscReal       h;
1063   PetscReal       tt,t;
1064   PetscScalar     *bt;
1065   const PetscReal *Bt = ros->tableau->binterpt;
1066   const PetscReal *GammaInv = ros->tableau->GammaInv;
1067   PetscScalar     *w        = ros->work;
1068   Vec             *Y        = ros->Y;
1069 
1070   PetscFunctionBegin;
1071   PetscCheck(Bt,PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ros->tableau->name);
1072 
1073   switch (ros->status) {
1074   case TS_STEP_INCOMPLETE:
1075   case TS_STEP_PENDING:
1076     h = ts->time_step;
1077     t = (itime - ts->ptime)/h;
1078     break;
1079   case TS_STEP_COMPLETE:
1080     h = ts->ptime - ts->ptime_prev;
1081     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
1082     break;
1083   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
1084   }
1085   PetscCall(PetscMalloc1(s,&bt));
1086   for (i=0; i<s; i++) bt[i] = 0;
1087   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
1088     for (i=0; i<s; i++) {
1089       bt[i] += Bt[i*pinterp+j] * tt;
1090     }
1091   }
1092 
1093   /* y(t+tt*h) = y(t) + Sum bt(tt) * GammaInv * Ydot */
1094   /* U <- 0*/
1095   PetscCall(VecZeroEntries(U));
1096   /* U <- Sum bt_i * GammaInv(i,1:i) * Y(1:i) */
1097   for (j=0; j<s; j++) w[j] = 0;
1098   for (j=0; j<s; j++) {
1099     for (i=j; i<s; i++) {
1100       w[j] +=  bt[i]*GammaInv[i*s+j];
1101     }
1102   }
1103   PetscCall(VecMAXPY(U,i,w,Y));
1104   /* U <- y(t) + U */
1105   PetscCall(VecAXPY(U,1,ros->vec_sol_prev));
1106 
1107   PetscCall(PetscFree(bt));
1108   PetscFunctionReturn(0);
1109 }
1110 
1111 /*------------------------------------------------------------*/
1112 
1113 static PetscErrorCode TSRosWTableauReset(TS ts)
1114 {
1115   TS_RosW        *ros = (TS_RosW*)ts->data;
1116   RosWTableau    tab  = ros->tableau;
1117 
1118   PetscFunctionBegin;
1119   if (!tab) PetscFunctionReturn(0);
1120   PetscCall(VecDestroyVecs(tab->s,&ros->Y));
1121   PetscCall(PetscFree(ros->work));
1122   PetscFunctionReturn(0);
1123 }
1124 
1125 static PetscErrorCode TSReset_RosW(TS ts)
1126 {
1127   TS_RosW        *ros = (TS_RosW*)ts->data;
1128 
1129   PetscFunctionBegin;
1130   PetscCall(TSRosWTableauReset(ts));
1131   PetscCall(VecDestroy(&ros->Ydot));
1132   PetscCall(VecDestroy(&ros->Ystage));
1133   PetscCall(VecDestroy(&ros->Zdot));
1134   PetscCall(VecDestroy(&ros->Zstage));
1135   PetscCall(VecDestroy(&ros->vec_sol_prev));
1136   PetscFunctionReturn(0);
1137 }
1138 
1139 static PetscErrorCode TSRosWGetVecs(TS ts,DM dm,Vec *Ydot,Vec *Zdot,Vec *Ystage,Vec *Zstage)
1140 {
1141   TS_RosW        *rw = (TS_RosW*)ts->data;
1142 
1143   PetscFunctionBegin;
1144   if (Ydot) {
1145     if (dm && dm != ts->dm) {
1146       PetscCall(DMGetNamedGlobalVector(dm,"TSRosW_Ydot",Ydot));
1147     } else *Ydot = rw->Ydot;
1148   }
1149   if (Zdot) {
1150     if (dm && dm != ts->dm) {
1151       PetscCall(DMGetNamedGlobalVector(dm,"TSRosW_Zdot",Zdot));
1152     } else *Zdot = rw->Zdot;
1153   }
1154   if (Ystage) {
1155     if (dm && dm != ts->dm) {
1156       PetscCall(DMGetNamedGlobalVector(dm,"TSRosW_Ystage",Ystage));
1157     } else *Ystage = rw->Ystage;
1158   }
1159   if (Zstage) {
1160     if (dm && dm != ts->dm) {
1161       PetscCall(DMGetNamedGlobalVector(dm,"TSRosW_Zstage",Zstage));
1162     } else *Zstage = rw->Zstage;
1163   }
1164   PetscFunctionReturn(0);
1165 }
1166 
1167 static PetscErrorCode TSRosWRestoreVecs(TS ts,DM dm,Vec *Ydot,Vec *Zdot, Vec *Ystage, Vec *Zstage)
1168 {
1169   PetscFunctionBegin;
1170   if (Ydot) {
1171     if (dm && dm != ts->dm) {
1172       PetscCall(DMRestoreNamedGlobalVector(dm,"TSRosW_Ydot",Ydot));
1173     }
1174   }
1175   if (Zdot) {
1176     if (dm && dm != ts->dm) {
1177       PetscCall(DMRestoreNamedGlobalVector(dm,"TSRosW_Zdot",Zdot));
1178     }
1179   }
1180   if (Ystage) {
1181     if (dm && dm != ts->dm) {
1182       PetscCall(DMRestoreNamedGlobalVector(dm,"TSRosW_Ystage",Ystage));
1183     }
1184   }
1185   if (Zstage) {
1186     if (dm && dm != ts->dm) {
1187       PetscCall(DMRestoreNamedGlobalVector(dm,"TSRosW_Zstage",Zstage));
1188     }
1189   }
1190   PetscFunctionReturn(0);
1191 }
1192 
1193 static PetscErrorCode DMCoarsenHook_TSRosW(DM fine,DM coarse,void *ctx)
1194 {
1195   PetscFunctionBegin;
1196   PetscFunctionReturn(0);
1197 }
1198 
1199 static PetscErrorCode DMRestrictHook_TSRosW(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
1200 {
1201   TS             ts = (TS)ctx;
1202   Vec            Ydot,Zdot,Ystage,Zstage;
1203   Vec            Ydotc,Zdotc,Ystagec,Zstagec;
1204 
1205   PetscFunctionBegin;
1206   PetscCall(TSRosWGetVecs(ts,fine,&Ydot,&Ystage,&Zdot,&Zstage));
1207   PetscCall(TSRosWGetVecs(ts,coarse,&Ydotc,&Ystagec,&Zdotc,&Zstagec));
1208   PetscCall(MatRestrict(restrct,Ydot,Ydotc));
1209   PetscCall(VecPointwiseMult(Ydotc,rscale,Ydotc));
1210   PetscCall(MatRestrict(restrct,Ystage,Ystagec));
1211   PetscCall(VecPointwiseMult(Ystagec,rscale,Ystagec));
1212   PetscCall(MatRestrict(restrct,Zdot,Zdotc));
1213   PetscCall(VecPointwiseMult(Zdotc,rscale,Zdotc));
1214   PetscCall(MatRestrict(restrct,Zstage,Zstagec));
1215   PetscCall(VecPointwiseMult(Zstagec,rscale,Zstagec));
1216   PetscCall(TSRosWRestoreVecs(ts,fine,&Ydot,&Ystage,&Zdot,&Zstage));
1217   PetscCall(TSRosWRestoreVecs(ts,coarse,&Ydotc,&Ystagec,&Zdotc,&Zstagec));
1218   PetscFunctionReturn(0);
1219 }
1220 
1221 static PetscErrorCode DMSubDomainHook_TSRosW(DM fine,DM coarse,void *ctx)
1222 {
1223   PetscFunctionBegin;
1224   PetscFunctionReturn(0);
1225 }
1226 
1227 static PetscErrorCode DMSubDomainRestrictHook_TSRosW(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
1228 {
1229   TS             ts = (TS)ctx;
1230   Vec            Ydot,Zdot,Ystage,Zstage;
1231   Vec            Ydots,Zdots,Ystages,Zstages;
1232 
1233   PetscFunctionBegin;
1234   PetscCall(TSRosWGetVecs(ts,dm,&Ydot,&Ystage,&Zdot,&Zstage));
1235   PetscCall(TSRosWGetVecs(ts,subdm,&Ydots,&Ystages,&Zdots,&Zstages));
1236 
1237   PetscCall(VecScatterBegin(gscat,Ydot,Ydots,INSERT_VALUES,SCATTER_FORWARD));
1238   PetscCall(VecScatterEnd(gscat,Ydot,Ydots,INSERT_VALUES,SCATTER_FORWARD));
1239 
1240   PetscCall(VecScatterBegin(gscat,Ystage,Ystages,INSERT_VALUES,SCATTER_FORWARD));
1241   PetscCall(VecScatterEnd(gscat,Ystage,Ystages,INSERT_VALUES,SCATTER_FORWARD));
1242 
1243   PetscCall(VecScatterBegin(gscat,Zdot,Zdots,INSERT_VALUES,SCATTER_FORWARD));
1244   PetscCall(VecScatterEnd(gscat,Zdot,Zdots,INSERT_VALUES,SCATTER_FORWARD));
1245 
1246   PetscCall(VecScatterBegin(gscat,Zstage,Zstages,INSERT_VALUES,SCATTER_FORWARD));
1247   PetscCall(VecScatterEnd(gscat,Zstage,Zstages,INSERT_VALUES,SCATTER_FORWARD));
1248 
1249   PetscCall(TSRosWRestoreVecs(ts,dm,&Ydot,&Ystage,&Zdot,&Zstage));
1250   PetscCall(TSRosWRestoreVecs(ts,subdm,&Ydots,&Ystages,&Zdots,&Zstages));
1251   PetscFunctionReturn(0);
1252 }
1253 
1254 /*
1255   This defines the nonlinear equation that is to be solved with SNES
1256   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
1257 */
1258 static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec U,Vec F,TS ts)
1259 {
1260   TS_RosW        *ros = (TS_RosW*)ts->data;
1261   Vec            Ydot,Zdot,Ystage,Zstage;
1262   PetscReal      shift = ros->scoeff / ts->time_step;
1263   DM             dm,dmsave;
1264 
1265   PetscFunctionBegin;
1266   PetscCall(SNESGetDM(snes,&dm));
1267   PetscCall(TSRosWGetVecs(ts,dm,&Ydot,&Zdot,&Ystage,&Zstage));
1268   PetscCall(VecWAXPY(Ydot,shift,U,Zdot));    /* Ydot = shift*U + Zdot */
1269   PetscCall(VecWAXPY(Ystage,1.0,U,Zstage));  /* Ystage = U + Zstage */
1270   dmsave = ts->dm;
1271   ts->dm = dm;
1272   PetscCall(TSComputeIFunction(ts,ros->stage_time,Ystage,Ydot,F,PETSC_FALSE));
1273   ts->dm = dmsave;
1274   PetscCall(TSRosWRestoreVecs(ts,dm,&Ydot,&Zdot,&Ystage,&Zstage));
1275   PetscFunctionReturn(0);
1276 }
1277 
1278 static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec U,Mat A,Mat B,TS ts)
1279 {
1280   TS_RosW        *ros = (TS_RosW*)ts->data;
1281   Vec            Ydot,Zdot,Ystage,Zstage;
1282   PetscReal      shift = ros->scoeff / ts->time_step;
1283   DM             dm,dmsave;
1284 
1285   PetscFunctionBegin;
1286   /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */
1287   PetscCall(SNESGetDM(snes,&dm));
1288   PetscCall(TSRosWGetVecs(ts,dm,&Ydot,&Zdot,&Ystage,&Zstage));
1289   dmsave = ts->dm;
1290   ts->dm = dm;
1291   PetscCall(TSComputeIJacobian(ts,ros->stage_time,Ystage,Ydot,shift,A,B,PETSC_TRUE));
1292   ts->dm = dmsave;
1293   PetscCall(TSRosWRestoreVecs(ts,dm,&Ydot,&Zdot,&Ystage,&Zstage));
1294   PetscFunctionReturn(0);
1295 }
1296 
1297 static PetscErrorCode TSRosWTableauSetUp(TS ts)
1298 {
1299   TS_RosW        *ros = (TS_RosW*)ts->data;
1300   RosWTableau    tab  = ros->tableau;
1301 
1302   PetscFunctionBegin;
1303   PetscCall(VecDuplicateVecs(ts->vec_sol,tab->s,&ros->Y));
1304   PetscCall(PetscMalloc1(tab->s,&ros->work));
1305   PetscFunctionReturn(0);
1306 }
1307 
1308 static PetscErrorCode TSSetUp_RosW(TS ts)
1309 {
1310   TS_RosW        *ros = (TS_RosW*)ts->data;
1311   DM             dm;
1312   SNES           snes;
1313   TSRHSJacobian  rhsjacobian;
1314 
1315   PetscFunctionBegin;
1316   PetscCall(TSRosWTableauSetUp(ts));
1317   PetscCall(VecDuplicate(ts->vec_sol,&ros->Ydot));
1318   PetscCall(VecDuplicate(ts->vec_sol,&ros->Ystage));
1319   PetscCall(VecDuplicate(ts->vec_sol,&ros->Zdot));
1320   PetscCall(VecDuplicate(ts->vec_sol,&ros->Zstage));
1321   PetscCall(VecDuplicate(ts->vec_sol,&ros->vec_sol_prev));
1322   PetscCall(TSGetDM(ts,&dm));
1323   PetscCall(DMCoarsenHookAdd(dm,DMCoarsenHook_TSRosW,DMRestrictHook_TSRosW,ts));
1324   PetscCall(DMSubDomainHookAdd(dm,DMSubDomainHook_TSRosW,DMSubDomainRestrictHook_TSRosW,ts));
1325   /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
1326   PetscCall(TSGetSNES(ts,&snes));
1327   if (!((PetscObject)snes)->type_name) {
1328     PetscCall(SNESSetType(snes,SNESKSPONLY));
1329   }
1330   PetscCall(DMTSGetRHSJacobian(dm,&rhsjacobian,NULL));
1331   if (rhsjacobian == TSComputeRHSJacobianConstant) {
1332     Mat Amat,Pmat;
1333 
1334     /* Set the SNES matrix to be different from the RHS matrix because there is no way to reconstruct shift*M-J */
1335     PetscCall(SNESGetJacobian(snes,&Amat,&Pmat,NULL,NULL));
1336     if (Amat && Amat == ts->Arhs) {
1337       if (Amat == Pmat) {
1338         PetscCall(MatDuplicate(ts->Arhs,MAT_COPY_VALUES,&Amat));
1339         PetscCall(SNESSetJacobian(snes,Amat,Amat,NULL,NULL));
1340       } else {
1341         PetscCall(MatDuplicate(ts->Arhs,MAT_COPY_VALUES,&Amat));
1342         PetscCall(SNESSetJacobian(snes,Amat,NULL,NULL,NULL));
1343         if (Pmat && Pmat == ts->Brhs) {
1344           PetscCall(MatDuplicate(ts->Brhs,MAT_COPY_VALUES,&Pmat));
1345           PetscCall(SNESSetJacobian(snes,NULL,Pmat,NULL,NULL));
1346           PetscCall(MatDestroy(&Pmat));
1347         }
1348       }
1349       PetscCall(MatDestroy(&Amat));
1350     }
1351   }
1352   PetscFunctionReturn(0);
1353 }
1354 /*------------------------------------------------------------*/
1355 
1356 static PetscErrorCode TSSetFromOptions_RosW(PetscOptionItems *PetscOptionsObject,TS ts)
1357 {
1358   TS_RosW        *ros = (TS_RosW*)ts->data;
1359   SNES           snes;
1360 
1361   PetscFunctionBegin;
1362   PetscCall(PetscOptionsHead(PetscOptionsObject,"RosW ODE solver options"));
1363   {
1364     RosWTableauLink link;
1365     PetscInt        count,choice;
1366     PetscBool       flg;
1367     const char      **namelist;
1368 
1369     for (link=RosWTableauList,count=0; link; link=link->next,count++) ;
1370     PetscCall(PetscMalloc1(count,(char***)&namelist));
1371     for (link=RosWTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
1372     PetscCall(PetscOptionsEList("-ts_rosw_type","Family of Rosenbrock-W method","TSRosWSetType",(const char*const*)namelist,count,ros->tableau->name,&choice,&flg));
1373     if (flg) PetscCall(TSRosWSetType(ts,namelist[choice]));
1374     PetscCall(PetscFree(namelist));
1375 
1376     PetscCall(PetscOptionsBool("-ts_rosw_recompute_jacobian","Recompute the Jacobian at each stage","TSRosWSetRecomputeJacobian",ros->recompute_jacobian,&ros->recompute_jacobian,NULL));
1377   }
1378   PetscCall(PetscOptionsTail());
1379   /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
1380   PetscCall(TSGetSNES(ts,&snes));
1381   if (!((PetscObject)snes)->type_name) {
1382     PetscCall(SNESSetType(snes,SNESKSPONLY));
1383   }
1384   PetscFunctionReturn(0);
1385 }
1386 
1387 static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer)
1388 {
1389   TS_RosW        *ros = (TS_RosW*)ts->data;
1390   PetscBool      iascii;
1391 
1392   PetscFunctionBegin;
1393   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
1394   if (iascii) {
1395     RosWTableau tab  = ros->tableau;
1396     TSRosWType  rostype;
1397     char        buf[512];
1398     PetscInt    i;
1399     PetscReal   abscissa[512];
1400     PetscCall(TSRosWGetType(ts,&rostype));
1401     PetscCall(PetscViewerASCIIPrintf(viewer,"  Rosenbrock-W %s\n",rostype));
1402     PetscCall(PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->ASum));
1403     PetscCall(PetscViewerASCIIPrintf(viewer,"  Abscissa of A       = %s\n",buf));
1404     for (i=0; i<tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i];
1405     PetscCall(PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,abscissa));
1406     PetscCall(PetscViewerASCIIPrintf(viewer,"  Abscissa of A+Gamma = %s\n",buf));
1407   }
1408   PetscFunctionReturn(0);
1409 }
1410 
1411 static PetscErrorCode TSLoad_RosW(TS ts,PetscViewer viewer)
1412 {
1413   SNES           snes;
1414   TSAdapt        adapt;
1415 
1416   PetscFunctionBegin;
1417   PetscCall(TSGetAdapt(ts,&adapt));
1418   PetscCall(TSAdaptLoad(adapt,viewer));
1419   PetscCall(TSGetSNES(ts,&snes));
1420   PetscCall(SNESLoad(snes,viewer));
1421   /* function and Jacobian context for SNES when used with TS is always ts object */
1422   PetscCall(SNESSetFunction(snes,NULL,NULL,ts));
1423   PetscCall(SNESSetJacobian(snes,NULL,NULL,NULL,ts));
1424   PetscFunctionReturn(0);
1425 }
1426 
1427 /*@C
1428   TSRosWSetType - Set the type of Rosenbrock-W scheme
1429 
1430   Logically collective
1431 
1432   Input Parameters:
1433 +  ts - timestepping context
1434 -  roswtype - type of Rosenbrock-W scheme
1435 
1436   Level: beginner
1437 
1438 .seealso: TSRosWGetType(), TSROSW, TSROSW2M, TSROSW2P, TSROSWRA3PW, TSROSWRA34PW2, TSROSWRODAS3, TSROSWSANDU3, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSROSWLLSSP3P4S2C, TSROSWARK3
1439 @*/
1440 PetscErrorCode TSRosWSetType(TS ts,TSRosWType roswtype)
1441 {
1442   PetscFunctionBegin;
1443   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1444   PetscValidCharPointer(roswtype,2);
1445   PetscTryMethod(ts,"TSRosWSetType_C",(TS,TSRosWType),(ts,roswtype));
1446   PetscFunctionReturn(0);
1447 }
1448 
1449 /*@C
1450   TSRosWGetType - Get the type of Rosenbrock-W scheme
1451 
1452   Logically collective
1453 
1454   Input Parameter:
1455 .  ts - timestepping context
1456 
1457   Output Parameter:
1458 .  rostype - type of Rosenbrock-W scheme
1459 
1460   Level: intermediate
1461 
1462 .seealso: TSRosWGetType()
1463 @*/
1464 PetscErrorCode TSRosWGetType(TS ts,TSRosWType *rostype)
1465 {
1466   PetscFunctionBegin;
1467   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1468   PetscUseMethod(ts,"TSRosWGetType_C",(TS,TSRosWType*),(ts,rostype));
1469   PetscFunctionReturn(0);
1470 }
1471 
1472 /*@C
1473   TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step.
1474 
1475   Logically collective
1476 
1477   Input Parameters:
1478 +  ts - timestepping context
1479 -  flg - PETSC_TRUE to recompute the Jacobian at each stage
1480 
1481   Level: intermediate
1482 
1483 .seealso: TSRosWGetType()
1484 @*/
1485 PetscErrorCode TSRosWSetRecomputeJacobian(TS ts,PetscBool flg)
1486 {
1487   PetscFunctionBegin;
1488   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1489   PetscTryMethod(ts,"TSRosWSetRecomputeJacobian_C",(TS,PetscBool),(ts,flg));
1490   PetscFunctionReturn(0);
1491 }
1492 
1493 static PetscErrorCode  TSRosWGetType_RosW(TS ts,TSRosWType *rostype)
1494 {
1495   TS_RosW        *ros = (TS_RosW*)ts->data;
1496 
1497   PetscFunctionBegin;
1498   *rostype = ros->tableau->name;
1499   PetscFunctionReturn(0);
1500 }
1501 
1502 static PetscErrorCode  TSRosWSetType_RosW(TS ts,TSRosWType rostype)
1503 {
1504   TS_RosW         *ros = (TS_RosW*)ts->data;
1505   PetscBool       match;
1506   RosWTableauLink link;
1507 
1508   PetscFunctionBegin;
1509   if (ros->tableau) {
1510     PetscCall(PetscStrcmp(ros->tableau->name,rostype,&match));
1511     if (match) PetscFunctionReturn(0);
1512   }
1513   for (link = RosWTableauList; link; link=link->next) {
1514     PetscCall(PetscStrcmp(link->tab.name,rostype,&match));
1515     if (match) {
1516       if (ts->setupcalled) PetscCall(TSRosWTableauReset(ts));
1517       ros->tableau = &link->tab;
1518       if (ts->setupcalled) PetscCall(TSRosWTableauSetUp(ts));
1519       ts->default_adapt_type = ros->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
1520       PetscFunctionReturn(0);
1521     }
1522   }
1523   SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rostype);
1524 }
1525 
1526 static PetscErrorCode  TSRosWSetRecomputeJacobian_RosW(TS ts,PetscBool flg)
1527 {
1528   TS_RosW *ros = (TS_RosW*)ts->data;
1529 
1530   PetscFunctionBegin;
1531   ros->recompute_jacobian = flg;
1532   PetscFunctionReturn(0);
1533 }
1534 
1535 static PetscErrorCode TSDestroy_RosW(TS ts)
1536 {
1537   PetscFunctionBegin;
1538   PetscCall(TSReset_RosW(ts));
1539   if (ts->dm) {
1540     PetscCall(DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRosW,DMRestrictHook_TSRosW,ts));
1541     PetscCall(DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRosW,DMSubDomainRestrictHook_TSRosW,ts));
1542   }
1543   PetscCall(PetscFree(ts->data));
1544   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRosWGetType_C",NULL));
1545   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRosWSetType_C",NULL));
1546   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRosWSetRecomputeJacobian_C",NULL));
1547   PetscFunctionReturn(0);
1548 }
1549 
1550 /* ------------------------------------------------------------ */
1551 /*MC
1552       TSROSW - ODE solver using Rosenbrock-W schemes
1553 
1554   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
1555   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
1556   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
1557 
1558   Notes:
1559   This method currently only works with autonomous ODE and DAE.
1560 
1561   Consider trying TSARKIMEX if the stiff part is strongly nonlinear.
1562 
1563   Since this uses a single linear solve per time-step if you wish to lag the jacobian or preconditioner computation you must use also -snes_lag_jacobian_persists true or -snes_lag_jacobian_preconditioner true
1564 
1565   Developer Notes:
1566   Rosenbrock-W methods are typically specified for autonomous ODE
1567 
1568 $  udot = f(u)
1569 
1570   by the stage equations
1571 
1572 $  k_i = h f(u_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j
1573 
1574   and step completion formula
1575 
1576 $  u_1 = u_0 + sum_j b_j k_j
1577 
1578   with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(u)
1579   and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner,
1580   we define new variables for the stage equations
1581 
1582 $  y_i = gamma_ij k_j
1583 
1584   The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define
1585 
1586 $  A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-1}
1587 
1588   to rewrite the method as
1589 
1590 $  [M/(h gamma_ii) - J] y_i = f(u_0 + sum_j a_ij y_j) + M sum_j (c_ij/h) y_j
1591 $  u_1 = u_0 + sum_j bt_j y_j
1592 
1593    where we have introduced the mass matrix M. Continue by defining
1594 
1595 $  ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j
1596 
1597    or, more compactly in tensor notation
1598 
1599 $  Ydot = 1/h (Gamma^{-1} \otimes I) Y .
1600 
1601    Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current
1602    stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the
1603    equation
1604 
1605 $  g(u_0 + sum_j a_ij y_j + y_i, ydot_i) = 0
1606 
1607    with initial guess y_i = 0.
1608 
1609   Level: beginner
1610 
1611 .seealso:  TSCreate(), TS, TSSetType(), TSRosWSetType(), TSRosWRegister(), TSROSWTHETA1, TSROSWTHETA2, TSROSW2M, TSROSW2P, TSROSWRA3PW, TSROSWRA34PW2, TSROSWRODAS3,
1612            TSROSWSANDU3, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSROSWLLSSP3P4S2C, TSROSWGRK4T, TSROSWSHAMP4, TSROSWVELDD4, TSROSW4L
1613 M*/
1614 PETSC_EXTERN PetscErrorCode TSCreate_RosW(TS ts)
1615 {
1616   TS_RosW        *ros;
1617 
1618   PetscFunctionBegin;
1619   PetscCall(TSRosWInitializePackage());
1620 
1621   ts->ops->reset          = TSReset_RosW;
1622   ts->ops->destroy        = TSDestroy_RosW;
1623   ts->ops->view           = TSView_RosW;
1624   ts->ops->load           = TSLoad_RosW;
1625   ts->ops->setup          = TSSetUp_RosW;
1626   ts->ops->step           = TSStep_RosW;
1627   ts->ops->interpolate    = TSInterpolate_RosW;
1628   ts->ops->evaluatestep   = TSEvaluateStep_RosW;
1629   ts->ops->rollback       = TSRollBack_RosW;
1630   ts->ops->setfromoptions = TSSetFromOptions_RosW;
1631   ts->ops->snesfunction   = SNESTSFormFunction_RosW;
1632   ts->ops->snesjacobian   = SNESTSFormJacobian_RosW;
1633 
1634   ts->usessnes = PETSC_TRUE;
1635 
1636   PetscCall(PetscNewLog(ts,&ros));
1637   ts->data = (void*)ros;
1638 
1639   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRosWGetType_C",TSRosWGetType_RosW));
1640   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRosWSetType_C",TSRosWSetType_RosW));
1641   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRosWSetRecomputeJacobian_C",TSRosWSetRecomputeJacobian_RosW));
1642 
1643   PetscCall(TSRosWSetType(ts,TSRosWDefault));
1644   PetscFunctionReturn(0);
1645 }
1646