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