xref: /petsc/src/ts/impls/rosw/rosw.c (revision f4aed992dd926c4e715d02cabab6ea04d956d070) !
1 /*
2   Code for timestepping with Rosenbrock W methods
3 
4   Notes:
5   The general system is written as
6 
7   G(t,X,Xdot) = F(t,X)
8 
9   where G represents the stiff part of the physics and F represents the non-stiff part.
10   This method is designed to be linearly implicit on G and can use an approximate and lagged Jacobian.
11 
12 */
13 #include <private/tsimpl.h>                /*I   "petscts.h"   I*/
14 
15 #include <../src/mat/blockinvert.h>
16 
17 static const 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,*binterp; /* 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   PetscScalar *work;            /* Scalar work space of length number of stages, used to prepare VecMAXPY() */
57   PetscReal   shift;
58   PetscReal   stage_time;
59   PetscReal   stage_explicit;     /* Flag indicates that the current stage is explicit */
60   PetscBool   recompute_jacobian; /* Recompute the Jacobian at each stage, default is to freeze the Jacobian at the start of each step */
61   TSStepStatus status;
62 } TS_RosW;
63 
64 /*MC
65      TSROSWTHETA1 - One stage first order L-stable Rosenbrock-W scheme (aka theta method).
66 
67      Only an approximate Jacobian is needed.
68 
69      Level: intermediate
70 
71 .seealso: TSROSW
72 M*/
73 
74 /*MC
75      TSROSWTHETA2 - One stage second order A-stable Rosenbrock-W scheme (aka theta method).
76 
77      Only an approximate Jacobian is needed.
78 
79      Level: intermediate
80 
81 .seealso: TSROSW
82 M*/
83 
84 /*MC
85      TSROSW2M - Two stage second order L-stable Rosenbrock-W scheme.
86 
87      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2P.
88 
89      Level: intermediate
90 
91 .seealso: TSROSW
92 M*/
93 
94 /*MC
95      TSROSW2P - Two stage second order L-stable Rosenbrock-W scheme.
96 
97      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2M.
98 
99      Level: intermediate
100 
101 .seealso: TSROSW
102 M*/
103 
104 /*MC
105      TSROSWRA3PW - Three stage third order Rosenbrock-W scheme for PDAE of index 1.
106 
107      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.
108 
109      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.
110 
111      References:
112      Rang and Angermann, New Rosenbrock-W methods of order 3 for partial differential algebraic equations of index 1, 2005.
113 
114      Level: intermediate
115 
116 .seealso: TSROSW
117 M*/
118 
119 /*MC
120      TSROSWRA34PW2 - Four stage third order L-stable Rosenbrock-W scheme for PDAE of index 1.
121 
122      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.
123 
124      This is strongly A-stable with R(infty) = 0. The embedded method of order 2 is strongly A-stable with R(infty) = 0.48.
125 
126      References:
127      Rang and Angermann, New Rosenbrock-W methods of order 3 for partial differential algebraic equations of index 1, 2005.
128 
129      Level: intermediate
130 
131 .seealso: TSROSW
132 M*/
133 
134 /*MC
135      TSROSWRODAS3 - Four stage third order L-stable Rosenbrock scheme
136 
137      By default, the Jacobian is only recomputed once per step.
138 
139      Both the third order and embedded second order methods are stiffly accurate and L-stable.
140 
141      References:
142      Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997.
143 
144      Level: intermediate
145 
146 .seealso: TSROSW, TSROSWSANDU3
147 M*/
148 
149 /*MC
150      TSROSWSANDU3 - Three stage third order L-stable Rosenbrock scheme
151 
152      By default, the Jacobian is only recomputed once per step.
153 
154      The third order method is L-stable, but not stiffly accurate.
155      The second order embedded method is strongly A-stable with R(infty) = 0.5.
156      The internal stages are L-stable.
157      This method is called ROS3 in the paper.
158 
159      References:
160      Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997.
161 
162      Level: intermediate
163 
164 .seealso: TSROSW, TSROSWRODAS3
165 M*/
166 
167 /*MC
168      TSROSWASSP3P3S1C - A-stable Rosenbrock-W method with SSP explicit part, third order, three stages
169 
170      By default, the Jacobian is only recomputed once per step.
171 
172      A-stable SPP explicit order 3, 3 stages, CFL 1 (eff = 1/3)
173 
174      References:
175      Emil Constantinescu
176 
177      Level: intermediate
178 
179 .seealso: TSROSW, TSROSWLASSP3P4S2C, TSROSWLLSSP3P4S2C, SSP
180 M*/
181 
182 /*MC
183      TSROSWLASSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, three stages
184 
185      By default, the Jacobian is only recomputed once per step.
186 
187      L-stable (A-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2)
188 
189      References:
190      Emil Constantinescu
191 
192      Level: intermediate
193 
194 .seealso: TSROSW, TSROSWASSP3P3S1C, TSROSWLLSSP3P4S2C, TSSSP
195 M*/
196 
197 /*MC
198      TSROSWLLSSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, three stages
199 
200      By default, the Jacobian is only recomputed once per step.
201 
202      L-stable (L-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2)
203 
204      References:
205      Emil Constantinescu
206 
207      Level: intermediate
208 
209 .seealso: TSROSW, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSSSP
210 M*/
211 
212 #undef __FUNCT__
213 #define __FUNCT__ "TSRosWRegisterAll"
214 /*@C
215   TSRosWRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSRosW
216 
217   Not Collective, but should be called by all processes which will need the schemes to be registered
218 
219   Level: advanced
220 
221 .keywords: TS, TSRosW, register, all
222 
223 .seealso:  TSRosWRegisterDestroy()
224 @*/
225 PetscErrorCode TSRosWRegisterAll(void)
226 {
227   PetscErrorCode ierr;
228 
229   PetscFunctionBegin;
230   if (TSRosWRegisterAllCalled) PetscFunctionReturn(0);
231   TSRosWRegisterAllCalled = PETSC_TRUE;
232 
233   {
234     const PetscReal
235       A = 0,
236       Gamma = 1,
237       b = 1;
238     ierr = TSRosWRegister(TSROSWTHETA1,1,1,&A,&Gamma,&b,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
239   }
240 
241   {
242     const PetscReal
243       A= 0,
244       Gamma = 0.5,
245       b = 1;
246     ierr = TSRosWRegister(TSROSWTHETA2,2,1,&A,&Gamma,&b,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
247   }
248 
249   {
250     const PetscReal g = 1. + 1./PetscSqrtReal(2.0);
251     const PetscReal
252       A[2][2] = {{0,0}, {1.,0}},
253       Gamma[2][2] = {{g,0}, {-2.*g,g}},
254       b[2] = {0.5,0.5},
255       b1[2] = {1.0,0.0};
256     ierr = TSRosWRegister(TSROSW2P,2,2,&A[0][0],&Gamma[0][0],b,b1,0,PETSC_NULL);CHKERRQ(ierr);
257   }
258   {
259     const PetscReal g = 1. - 1./PetscSqrtReal(2.0);
260     const PetscReal
261       A[2][2] = {{0,0}, {1.,0}},
262       Gamma[2][2] = {{g,0}, {-2.*g,g}},
263       b[2] = {0.5,0.5},
264       b1[2] = {1.0,0.0};
265     ierr = TSRosWRegister(TSROSW2M,2,2,&A[0][0],&Gamma[0][0],b,b1,0,PETSC_NULL);CHKERRQ(ierr);
266   }
267   {
268     const PetscReal g = 7.8867513459481287e-01;
269     const PetscReal
270       A[3][3] = {{0,0,0},
271                  {1.5773502691896257e+00,0,0},
272                  {0.5,0,0}},
273       Gamma[3][3] = {{g,0,0},
274                      {-1.5773502691896257e+00,g,0},
275                      {-6.7075317547305480e-01,-1.7075317547305482e-01,g}},
276       b[3] = {1.0566243270259355e-01,4.9038105676657971e-02,8.4529946162074843e-01},
277       b2[3] = {-1.7863279495408180e-01,1./3.,8.4529946162074843e-01};
278     ierr = TSRosWRegister(TSROSWRA3PW,3,3,&A[0][0],&Gamma[0][0],b,b2,0,PETSC_NULL);CHKERRQ(ierr);
279   }
280   {
281     const PetscReal g = 4.3586652150845900e-01;
282     const PetscReal
283       A[4][4] = {{0,0,0,0},
284                  {8.7173304301691801e-01,0,0,0},
285                  {8.4457060015369423e-01,-1.1299064236484185e-01,0,0},
286                  {0,0,1.,0}},
287       Gamma[4][4] = {{g,0,0,0},
288                      {-8.7173304301691801e-01,g,0,0},
289                      {-9.0338057013044082e-01,5.4180672388095326e-02,g,0},
290                      {2.4212380706095346e-01,-1.2232505839045147e+00,5.4526025533510214e-01,g}},
291         b[4] = {2.4212380706095346e-01,-1.2232505839045147e+00,1.5452602553351020e+00,4.3586652150845900e-01},
292           b2[4] = {3.7810903145819369e-01,-9.6042292212423178e-02,5.0000000000000000e-01,2.1793326075422950e-01},
293             binterpt[4][3]={{1.0564298455794094,-1.3864882699759573,0.5721822314575016},{2.296429974281067,-8.262611700275677,4.742931142090097},{-1.307599564525376,7.250979895056055,-4.398120075195578},{-1.045260255335102,2.398120075195581,-0.9169932983520199}};
294       ierr = TSRosWRegister(TSROSWRA34PW2,3,4,&A[0][0],&Gamma[0][0],b,b2,3,binterpt[0]);CHKERRQ(ierr);
295   }
296   {
297     const PetscReal g = 0.5;
298     const PetscReal
299       A[4][4] = {{0,0,0,0},
300                  {0,0,0,0},
301                  {1.,0,0,0},
302                  {0.75,-0.25,0.5,0}},
303       Gamma[4][4] = {{g,0,0,0},
304                      {1.,g,0,0},
305                      {-0.25,-0.25,g,0},
306                      {1./12,1./12,-2./3,g}},
307       b[4] = {5./6,-1./6,-1./6,0.5},
308       b2[4] = {0.75,-0.25,0.5,0};
309     ierr = TSRosWRegister(TSROSWRODAS3,3,4,&A[0][0],&Gamma[0][0],b,b2,0,PETSC_NULL);CHKERRQ(ierr);
310   }
311   {
312     const PetscReal g = 0.43586652150845899941601945119356;
313     const PetscReal
314       A[3][3] = {{0,0,0},
315                  {g,0,0},
316                  {g,0,0}},
317       Gamma[3][3] = {{g,0,0},
318                      {-0.19294655696029095575009695436041,g,0},
319                      {0,1.74927148125794685173529749738960,g}},
320       b[3] = {-0.75457412385404315829818998646589,1.94100407061964420292840123379419,-0.18642994676560104463021124732829},
321       b2[3] = {-1.53358745784149585370766523913002,2.81745131148625772213931745457622,-0.28386385364476186843165221544619};
322     ierr = TSRosWRegister(TSROSWSANDU3,3,3,&A[0][0],&Gamma[0][0],b,b2,0,PETSC_NULL);CHKERRQ(ierr);
323   }
324   {
325     const PetscReal s3 = PetscSqrtReal(3.),g = (3.0+s3)/6.0;
326     const PetscReal
327       A[3][3] = {{0,0,0},
328                  {1,0,0},
329                  {0.25,0.25,0}},
330       Gamma[3][3] = {{0,0,0},
331                      {(-3.0-s3)/6.0,g,0},
332                      {(-3.0-s3)/24.0,(-3.0-s3)/8.0,g}},
333         b[3] = {1./6.,1./6.,2./3.},
334           b2[3] = {1./4.,1./4.,1./2.};
335     ierr = TSRosWRegister(TSROSWASSP3P3S1C,3,3,&A[0][0],&Gamma[0][0],b,b2,0,PETSC_NULL);CHKERRQ(ierr);
336   }
337 
338   {
339     const PetscReal
340       A[4][4] = {{0,0,0,0},
341                  {1./2.,0,0,0},
342                  {1./2.,1./2.,0,0},
343                  {1./6.,1./6.,1./6.,0}},
344       Gamma[4][4] = {{1./2.,0,0,0},
345                      {0.0,1./4.,0,0},
346                      {-2.,-2./3.,2./3.,0},
347                      {1./2.,5./36.,-2./9,0}},
348         b[4] = {1./6.,1./6.,1./6.,1./2.},
349         b2[4] = {1./8.,3./4.,1./8.,0};
350      ierr = TSRosWRegister(TSROSWLASSP3P4S2C,3,4,&A[0][0],&Gamma[0][0],b,b2,0,PETSC_NULL);CHKERRQ(ierr);
351   }
352 
353   {
354     const PetscReal
355       A[4][4] = {{0,0,0,0},
356                  {1./2.,0,0,0},
357                  {1./2.,1./2.,0,0},
358                  {1./6.,1./6.,1./6.,0}},
359       Gamma[4][4] = {{1./2.,0,0,0},
360                      {0.0,3./4.,0,0},
361                      {-2./3.,-23./9.,2./9.,0},
362                      {1./18.,65./108.,-2./27,0}},
363         b[4] = {1./6.,1./6.,1./6.,1./2.},
364         b2[4] = {3./16.,10./16.,3./16.,0};
365      ierr = TSRosWRegister(TSROSWLLSSP3P4S2C,3,4,&A[0][0],&Gamma[0][0],b,b2,0,PETSC_NULL);CHKERRQ(ierr);
366   }
367 
368  {
369    PetscReal A[4][4],Gamma[4][4],b[4],b2[4];
370 
371    Gamma[0][0]=0.4358665215084589994160194475295062513822671686978816;
372    Gamma[0][1]=0; Gamma[0][2]=0; Gamma[0][3]=0;
373    Gamma[1][0]=-1.997527830934941248426324674704153457289527280554476;
374    Gamma[1][1]=0.4358665215084589994160194475295062513822671686978816;
375    Gamma[1][2]=0; Gamma[1][3]=0;
376    Gamma[2][0]=-1.007948511795029620852002345345404191008352770119903;
377    Gamma[2][1]=-0.004648958462629345562774289390054679806993396798458131;
378    Gamma[2][2]=0.4358665215084589994160194475295062513822671686978816;
379    Gamma[2][3]=0;
380    Gamma[3][0]=-0.6685429734233467180451604600279552604364311322650783;
381    Gamma[3][1]=0.6056625986449338476089525334450053439525178740492984;
382    Gamma[3][2]=-0.9717899277217721234705114616271378792182450260943198;
383    Gamma[3][3]=0;
384 
385    A[0][0]=0; A[0][1]=0; A[0][2]=0; A[0][3]=0;
386    A[1][0]=0.8717330430169179988320388950590125027645343373957631;
387    A[1][1]=0; A[1][2]=0; A[1][3]=0;
388    A[2][0]=0.5275890119763004115618079766722914408876108660811028;
389    A[2][1]=0.07241098802369958843819203208518599088698057726988732;
390    A[2][2]=0; A[2][3]=0;
391    A[3][0]=0.3990960076760701320627260685975778145384666450351314;
392    A[3][1]=-0.4375576546135194437228463747348862825846903771419953;
393    A[3][2]=1.038461646937449311660120300601880176655352737312713;
394    A[3][3]=0;
395 
396    b[0]=0.1876410243467238251612921333138006734899663569186926;
397    b[1]=-0.5952974735769549480478230473706443582188442040780541;
398    b[2]=0.9717899277217721234705114616271378792182450260943198;
399    b[3]=0.4358665215084589994160194475295062513822671686978816;
400 
401    b2[0]=0.2147402862233891404862383521089097657790734483804460;
402    b2[1]=-0.4851622638849390928209050538171743017757490232519684;
403    b2[2]=0.8687250025203875511662123688667549217531982787600080;
404    b2[3]=0.4016969751411624011684543450940068201770721128357014;
405 
406    PetscReal binterpt[4][3]={{2.2565812720167954547104627844105,-3.0826699111559187902922463354557,1.0137296634858471607430756831148},{1.349166413351089573796243820819,-2.4689115685996042534544925650515,0.52444768167155973161042570784064},{-2.4695174540533503758652847586647,5.7428279814696677152129332773553,-2.3015205996945452158771370439586},{-0.13623023131453465264142184656474,-0.19124650171414467146619437684812,0.76334325453713832352363565300308}};
407 
408    ierr = TSRosWRegister(TSROSWARK3,3,4,&A[0][0],&Gamma[0][0],b,b2,3,binterpt);CHKERRQ(ierr);
409   }
410 
411   PetscFunctionReturn(0);
412 }
413 
414 #undef __FUNCT__
415 #define __FUNCT__ "TSRosWRegisterDestroy"
416 /*@C
417    TSRosWRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister().
418 
419    Not Collective
420 
421    Level: advanced
422 
423 .keywords: TSRosW, register, destroy
424 .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegisterDynamic()
425 @*/
426 PetscErrorCode TSRosWRegisterDestroy(void)
427 {
428   PetscErrorCode ierr;
429   RosWTableauLink link;
430 
431   PetscFunctionBegin;
432   while ((link = RosWTableauList)) {
433     RosWTableau t = &link->tab;
434     RosWTableauList = link->next;
435     ierr = PetscFree5(t->A,t->Gamma,t->b,t->ASum,t->GammaSum);CHKERRQ(ierr);
436     ierr = PetscFree5(t->At,t->bt,t->GammaInv,t->GammaZeroDiag,t->GammaExplicitCorr);CHKERRQ(ierr);
437     ierr = PetscFree2(t->bembed,t->bembedt);CHKERRQ(ierr);
438     ierr = PetscFree(t->binterpt);CHKERRQ(ierr);
439     ierr = PetscFree(t->name);CHKERRQ(ierr);
440     ierr = PetscFree(link);CHKERRQ(ierr);
441   }
442   TSRosWRegisterAllCalled = PETSC_FALSE;
443   PetscFunctionReturn(0);
444 }
445 
446 #undef __FUNCT__
447 #define __FUNCT__ "TSRosWInitializePackage"
448 /*@C
449   TSRosWInitializePackage - This function initializes everything in the TSRosW package. It is called
450   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RosW()
451   when using static libraries.
452 
453   Input Parameter:
454   path - The dynamic library path, or PETSC_NULL
455 
456   Level: developer
457 
458 .keywords: TS, TSRosW, initialize, package
459 .seealso: PetscInitialize()
460 @*/
461 PetscErrorCode TSRosWInitializePackage(const char path[])
462 {
463   PetscErrorCode ierr;
464 
465   PetscFunctionBegin;
466   if (TSRosWPackageInitialized) PetscFunctionReturn(0);
467   TSRosWPackageInitialized = PETSC_TRUE;
468   ierr = TSRosWRegisterAll();CHKERRQ(ierr);
469   ierr = PetscRegisterFinalize(TSRosWFinalizePackage);CHKERRQ(ierr);
470   PetscFunctionReturn(0);
471 }
472 
473 #undef __FUNCT__
474 #define __FUNCT__ "TSRosWFinalizePackage"
475 /*@C
476   TSRosWFinalizePackage - This function destroys everything in the TSRosW package. It is
477   called from PetscFinalize().
478 
479   Level: developer
480 
481 .keywords: Petsc, destroy, package
482 .seealso: PetscFinalize()
483 @*/
484 PetscErrorCode TSRosWFinalizePackage(void)
485 {
486   PetscErrorCode ierr;
487 
488   PetscFunctionBegin;
489   TSRosWPackageInitialized = PETSC_FALSE;
490   ierr = TSRosWRegisterDestroy();CHKERRQ(ierr);
491   PetscFunctionReturn(0);
492 }
493 
494 #undef __FUNCT__
495 #define __FUNCT__ "TSRosWRegister"
496 /*@C
497    TSRosWRegister - register a Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
498 
499    Not Collective, but the same schemes should be registered on all processes on which they will be used
500 
501    Input Parameters:
502 +  name - identifier for method
503 .  order - approximation order of method
504 .  s - number of stages, this is the dimension of the matrices below
505 .  A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular
506 .  Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal
507 .  b - Step completion table (dimension s)
508 -  bembed - Step completion table for a scheme of order one less (dimension s, PETSC_NULL if no embedded scheme is available)
509 .  pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt
510 .  binterpt - Coefficients of the interpolation formula (dimension s*pinterp)
511 
512    Notes:
513    Several Rosenbrock W methods are provided, this function is only needed to create new methods.
514 
515    Level: advanced
516 
517 .keywords: TS, register
518 
519 .seealso: TSRosW
520 @*/
521 PetscErrorCode TSRosWRegister(const TSRosWType name,PetscInt order,PetscInt s,
522                               const PetscReal A[],const PetscReal Gamma[],const PetscReal b[],const PetscReal bembed[],
523                                  PetscInt pinterp,const PetscReal binterpt[])
524 {
525   PetscErrorCode ierr;
526   RosWTableauLink link;
527   RosWTableau t;
528   PetscInt i,j,k;
529   PetscScalar *GammaInv;
530 
531   PetscFunctionBegin;
532   PetscValidCharPointer(name,1);
533   PetscValidPointer(A,4);
534   PetscValidPointer(Gamma,5);
535   PetscValidPointer(b,6);
536   if (bembed) PetscValidPointer(bembed,7);
537 
538   ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
539   ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
540   t = &link->tab;
541   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
542   t->order = order;
543   t->s = s;
544   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);
545   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);
546   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
547   ierr = PetscMemcpy(t->Gamma,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr);
548   ierr = PetscMemcpy(t->GammaExplicitCorr,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr);
549   ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
550   if (bembed) {
551     ierr = PetscMalloc2(s,PetscReal,&t->bembed,s,PetscReal,&t->bembedt);CHKERRQ(ierr);
552     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
553   }
554   for (i=0; i<s; i++) {
555     t->ASum[i] = 0;
556     t->GammaSum[i] = 0;
557     for (j=0; j<s; j++) {
558       t->ASum[i] += A[i*s+j];
559       t->GammaSum[i] += Gamma[i*s+j];
560     }
561   }
562   ierr = PetscMalloc(s*s*sizeof(PetscScalar),&GammaInv);CHKERRQ(ierr); /* Need to use Scalar for inverse, then convert back to Real */
563   for (i=0; i<s*s; i++) GammaInv[i] = Gamma[i];
564   for (i=0; i<s; i++) {
565     if (Gamma[i*s+i] == 0.0) {
566       GammaInv[i*s+i] = 1.0;
567       t->GammaZeroDiag[i] = PETSC_TRUE;
568     } else {
569       t->GammaZeroDiag[i] = PETSC_FALSE;
570     }
571   }
572 
573   switch (s) {
574   case 1: GammaInv[0] = 1./GammaInv[0]; break;
575   case 2: ierr = Kernel_A_gets_inverse_A_2(GammaInv,0);CHKERRQ(ierr); break;
576   case 3: ierr = Kernel_A_gets_inverse_A_3(GammaInv,0);CHKERRQ(ierr); break;
577   case 4: ierr = Kernel_A_gets_inverse_A_4(GammaInv,0);CHKERRQ(ierr); break;
578   case 5: {
579     PetscInt ipvt5[5];
580     MatScalar work5[5*5];
581     ierr = Kernel_A_gets_inverse_A_5(GammaInv,ipvt5,work5,0);CHKERRQ(ierr); break;
582   }
583   case 6: ierr = Kernel_A_gets_inverse_A_6(GammaInv,0);CHKERRQ(ierr); break;
584   case 7: ierr = Kernel_A_gets_inverse_A_7(GammaInv,0);CHKERRQ(ierr); break;
585   default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for %D stages",s);
586   }
587   for (i=0; i<s*s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]);
588   ierr = PetscFree(GammaInv);CHKERRQ(ierr);
589 
590   for (i=0; i<s; i++) {
591     for (k=0; k<i+1; k++) {
592       t->GammaExplicitCorr[i*s+k]=(t->GammaExplicitCorr[i*s+k])*(t->GammaInv[k*s+k]);
593       for (j=k+1; j<i+1; j++) {
594         t->GammaExplicitCorr[i*s+k]+=(t->GammaExplicitCorr[i*s+j])*(t->GammaInv[j*s+k]);
595       }
596     }
597   }
598 
599   for (i=0; i<s; i++) {
600     for (j=0; j<s; j++) {
601       t->At[i*s+j] = 0;
602       for (k=0; k<s; k++) {
603         t->At[i*s+j] += t->A[i*s+k] * t->GammaInv[k*s+j];
604       }
605     }
606     t->bt[i] = 0;
607     for (j=0; j<s; j++) {
608       t->bt[i] += t->b[j] * t->GammaInv[j*s+i];
609     }
610     if (bembed) {
611       t->bembedt[i] = 0;
612       for (j=0; j<s; j++) {
613         t->bembedt[i] += t->bembed[j] * t->GammaInv[j*s+i];
614       }
615     }
616   }
617   t->ccfl = 1.0;                /* Fix this */
618 
619   t->pinterp = pinterp;
620   ierr = PetscMalloc(s*pinterp,&t->binterpt);CHKERRQ(ierr);
621 
622   link->next = RosWTableauList;
623   RosWTableauList = link;
624   PetscFunctionReturn(0);
625 }
626 
627 #undef __FUNCT__
628 #define __FUNCT__ "TSEvaluateStep_RosW"
629 /*
630  The step completion formula is
631 
632  x1 = x0 + b^T Y
633 
634  where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been
635  updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write
636 
637  x1e = x0 + be^T Y
638      = x1 - b^T Y + be^T Y
639      = x1 + (be - b)^T Y
640 
641  so we can evaluate the method of different order even after the step has been optimistically completed.
642 */
643 static PetscErrorCode TSEvaluateStep_RosW(TS ts,PetscInt order,Vec X,PetscBool *done)
644 {
645   TS_RosW        *ros = (TS_RosW*)ts->data;
646   RosWTableau    tab  = ros->tableau;
647   PetscScalar    *w = ros->work;
648   PetscInt       i;
649   PetscErrorCode ierr;
650 
651   PetscFunctionBegin;
652   if (order == tab->order) {
653     if (ros->status == TS_STEP_INCOMPLETE) { /* Use standard completion formula */
654       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
655       for (i=0; i<tab->s; i++) w[i] = tab->bt[i];
656       ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr);
657     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
658     if (done) *done = PETSC_TRUE;
659     PetscFunctionReturn(0);
660   } else if (order == tab->order-1) {
661     if (!tab->bembedt) goto unavailable;
662     if (ros->status == TS_STEP_INCOMPLETE) { /* Use embedded completion formula */
663       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
664       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i];
665       ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr);
666     } else {                    /* Use rollback-and-recomplete formula (bembedt - bt) */
667       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i];
668       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
669       ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr);
670     }
671     if (done) *done = PETSC_TRUE;
672     PetscFunctionReturn(0);
673   }
674   unavailable:
675   if (done) *done = PETSC_FALSE;
676   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);
677   PetscFunctionReturn(0);
678 }
679 
680 #undef __FUNCT__
681 #define __FUNCT__ "TSStep_RosW"
682 static PetscErrorCode TSStep_RosW(TS ts)
683 {
684   TS_RosW         *ros = (TS_RosW*)ts->data;
685   RosWTableau     tab  = ros->tableau;
686   const PetscInt  s    = tab->s;
687   const PetscReal *At  = tab->At,*Gamma = tab->Gamma,*ASum = tab->ASum,*GammaInv = tab->GammaInv;
688   const PetscReal *GammaExplicitCorr = tab->GammaExplicitCorr;
689   const PetscBool *GammaZeroDiag = tab->GammaZeroDiag;
690   PetscScalar     *w   = ros->work;
691   Vec             *Y   = ros->Y,Ydot = ros->Ydot,Zdot = ros->Zdot,Zstage = ros->Zstage;
692   SNES            snes;
693   TSAdapt         adapt;
694   PetscInt        i,j,its,lits,reject,next_scheme;
695   PetscReal       next_time_step;
696   PetscBool       accept;
697   PetscErrorCode  ierr;
698   MatStructure    str;
699 
700   PetscFunctionBegin;
701   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
702   next_time_step = ts->time_step;
703   accept = PETSC_TRUE;
704   ros->status = TS_STEP_INCOMPLETE;
705 
706   for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) {
707     const PetscReal h = ts->time_step;
708     for (i=0; i<s; i++) {
709       ros->stage_time = ts->ptime + h*ASum[i];
710       if (GammaZeroDiag[i]) {
711         ros->stage_explicit = PETSC_TRUE;
712         ros->shift = 1./h;
713       } else {
714         ros->stage_explicit = PETSC_FALSE;
715         ros->shift = 1./(h*Gamma[i*s+i]);
716       }
717 
718       ierr = VecCopy(ts->vec_sol,Zstage);CHKERRQ(ierr);
719       for (j=0; j<i; j++) w[j] = At[i*s+j];
720       ierr = VecMAXPY(Zstage,i,w,Y);CHKERRQ(ierr);
721 
722       for (j=0; j<i; j++) w[j] = 1./h * GammaInv[i*s+j];
723       ierr = VecZeroEntries(Zdot);CHKERRQ(ierr);
724       ierr = VecMAXPY(Zdot,i,w,Y);CHKERRQ(ierr);
725 
726       /* Initial guess taken from last stage */
727       ierr = VecZeroEntries(Y[i]);CHKERRQ(ierr);
728 
729       if (!ros->stage_explicit) {
730         if (!ros->recompute_jacobian && !i) {
731           ierr = SNESSetLagJacobian(snes,-2);CHKERRQ(ierr); /* Recompute the Jacobian on this solve, but not again */
732         }
733         ierr = SNESSolve(snes,PETSC_NULL,Y[i]);CHKERRQ(ierr);
734         ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
735         ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
736         ts->nonlinear_its += its; ts->linear_its += lits;
737         ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
738         ierr = TSAdaptCheckStage(adapt,ts,&accept);CHKERRQ(ierr);
739         if (!accept) goto reject_step;
740       } else {
741         Mat J,Jp;
742         ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); /* Evaluate Y[i]=G(t,Ydot=0,Zstage) */
743         ierr = TSComputeIFunction(ts,ros->stage_time,Zstage,Ydot,Y[i],PETSC_FALSE);CHKERRQ(ierr);
744         ierr = VecScale(Y[i],-1.0);
745         ierr = VecAXPY(Y[i],-1.0,Zdot);CHKERRQ(ierr); /*Y[i]=F(Zstage)-Zdot[=GammaInv*Y]*/
746 
747         ierr = VecZeroEntries(Zstage);CHKERRQ(ierr); /* Zstage = GammaExplicitCorr[i,j] * Y[j] */
748         for (j=0; j<i; j++) w[j] = GammaExplicitCorr[i*s+j];
749         ierr = VecMAXPY(Zstage,i,w,Y);CHKERRQ(ierr);
750         /*Y[i] += Y[i] + Jac*Zstage[=Jac*GammaExplicitCorr[i,j] * Y[j]] */
751         str = SAME_NONZERO_PATTERN;
752         ierr = TSGetIJacobian(ts,&J,&Jp,PETSC_NULL,PETSC_NULL);
753         ierr = TSComputeIJacobian(ts,ros->stage_time,ts->vec_sol,Ydot,0,&J,&Jp,&str,PETSC_FALSE);CHKERRQ(ierr);
754         ierr = MatMult(J,Zstage,Zdot);
755 
756         ierr = VecAXPY(Y[i],-1.0,Zdot);CHKERRQ(ierr);
757         ierr = VecScale(Y[i],h);
758         ts->linear_its += 1;
759       }
760     }
761     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
762     ros->status = TS_STEP_PENDING;
763 
764     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
765     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
766     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
767     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
768     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
769     if (accept) {
770       /* ignore next_scheme for now */
771       ts->ptime += ts->time_step;
772       ts->time_step = next_time_step;
773       ts->steps++;
774       ros->status = TS_STEP_COMPLETE;
775       break;
776     } else {                    /* Roll back the current step */
777       for (i=0; i<s; i++) w[i] = -tab->bt[i];
778       ierr = VecMAXPY(ts->vec_sol,s,w,Y);CHKERRQ(ierr);
779       ts->time_step = next_time_step;
780       ros->status = TS_STEP_INCOMPLETE;
781     }
782     reject_step: continue;
783   }
784   if (ros->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED;
785   PetscFunctionReturn(0);
786 }
787 
788 #undef __FUNCT__
789 #define __FUNCT__ "TSInterpolate_RosW"
790 static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec X)
791 {
792   TS_RosW        *ros = (TS_RosW*)ts->data;
793   PetscInt s = ros->tableau->s,pinterp = ros->tableau->pinterp,i,j;
794   PetscReal h;
795   PetscReal tt,t;
796   PetscScalar *bt;
797   const PetscReal *Bt = ros->tableau->binterpt;
798   PetscErrorCode ierr;
799   const PetscReal *GammaInv = ros->tableau->GammaInv;
800   PetscScalar     *w   = ros->work;
801   Vec             *Y   = ros->Y;
802 
803   PetscFunctionBegin;
804   //  SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ros->tableau->name);
805   if (!Bt) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ros->tableau->name);
806 
807   switch (ros->status) {
808   case TS_STEP_INCOMPLETE:
809   case TS_STEP_PENDING:
810     h = ts->time_step;
811     t = (itime - ts->ptime)/h;
812     break;
813   case TS_STEP_COMPLETE:
814     h = ts->time_step_prev;
815     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
816     break;
817   default: SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Invalid TSStepStatus");
818   }
819   ierr = PetscMalloc(s,&bt);CHKERRQ(ierr);
820   for (i=0; i<s; i++) bt[i] = 0;
821   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
822     for (i=0; i<s; i++) {
823       bt[i] += h * Bt[i*pinterp+j] * tt;
824     }
825   }
826 
827   /* y(t+tt*h) = y(t) + Sum bt(tt) * GammaInv * Ydot */
828   /*X<-0*/
829   ierr = VecZeroEntries(X);CHKERRQ(ierr);
830 
831   /*X<- Sum bt_i * GammaInv(i,1:i) * Y(1:i) */
832   for (i=0; i<s; i++) {
833     for (j=0; j<=i; j++) w[j] =  bt[i]*GammaInv[i*s+j];
834     ierr = VecMAXPY(X,i,w,Y);CHKERRQ(ierr);
835   }
836 
837   /*X<-y(t) + X*/
838   ierr = VecAXPY(X,1.0,ts->vec_sol);CHKERRQ(ierr);
839 
840   ierr = PetscFree(bt);CHKERRQ(ierr);
841 
842   PetscFunctionReturn(0);
843 }
844 
845 /*------------------------------------------------------------*/
846 #undef __FUNCT__
847 #define __FUNCT__ "TSReset_RosW"
848 static PetscErrorCode TSReset_RosW(TS ts)
849 {
850   TS_RosW        *ros = (TS_RosW*)ts->data;
851   PetscInt       s;
852   PetscErrorCode ierr;
853 
854   PetscFunctionBegin;
855   if (!ros->tableau) PetscFunctionReturn(0);
856    s = ros->tableau->s;
857   ierr = VecDestroyVecs(s,&ros->Y);CHKERRQ(ierr);
858   ierr = VecDestroy(&ros->Ydot);CHKERRQ(ierr);
859   ierr = VecDestroy(&ros->Ystage);CHKERRQ(ierr);
860   ierr = VecDestroy(&ros->Zdot);CHKERRQ(ierr);
861   ierr = VecDestroy(&ros->Zstage);CHKERRQ(ierr);
862   ierr = PetscFree(ros->work);CHKERRQ(ierr);
863   PetscFunctionReturn(0);
864 }
865 
866 #undef __FUNCT__
867 #define __FUNCT__ "TSDestroy_RosW"
868 static PetscErrorCode TSDestroy_RosW(TS ts)
869 {
870   PetscErrorCode  ierr;
871 
872   PetscFunctionBegin;
873   ierr = TSReset_RosW(ts);CHKERRQ(ierr);
874   ierr = PetscFree(ts->data);CHKERRQ(ierr);
875   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","",PETSC_NULL);CHKERRQ(ierr);
876   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","",PETSC_NULL);CHKERRQ(ierr);
877   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","",PETSC_NULL);CHKERRQ(ierr);
878   PetscFunctionReturn(0);
879 }
880 
881 /*
882   This defines the nonlinear equation that is to be solved with SNES
883   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
884 */
885 #undef __FUNCT__
886 #define __FUNCT__ "SNESTSFormFunction_RosW"
887 static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec X,Vec F,TS ts)
888 {
889   TS_RosW        *ros = (TS_RosW*)ts->data;
890   PetscErrorCode ierr;
891 
892   PetscFunctionBegin;
893   ierr = VecWAXPY(ros->Ydot,ros->shift,X,ros->Zdot);CHKERRQ(ierr); /* Ydot = shift*X + Zdot */
894   ierr = VecWAXPY(ros->Ystage,1.0,X,ros->Zstage);CHKERRQ(ierr);    /* Ystage = X + Zstage */
895   ierr = TSComputeIFunction(ts,ros->stage_time,ros->Ystage,ros->Ydot,F,PETSC_FALSE);CHKERRQ(ierr);
896   PetscFunctionReturn(0);
897 }
898 
899 #undef __FUNCT__
900 #define __FUNCT__ "SNESTSFormJacobian_RosW"
901 static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts)
902 {
903   TS_RosW        *ros = (TS_RosW*)ts->data;
904   PetscErrorCode ierr;
905 
906   PetscFunctionBegin;
907   /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */
908   ierr = TSComputeIJacobian(ts,ros->stage_time,ros->Ystage,ros->Ydot,ros->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr);
909   PetscFunctionReturn(0);
910 }
911 
912 #undef __FUNCT__
913 #define __FUNCT__ "TSSetUp_RosW"
914 static PetscErrorCode TSSetUp_RosW(TS ts)
915 {
916   TS_RosW        *ros = (TS_RosW*)ts->data;
917   RosWTableau    tab  = ros->tableau;
918   PetscInt       s    = tab->s;
919   PetscErrorCode ierr;
920 
921   PetscFunctionBegin;
922   if (!ros->tableau) {
923     ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);
924   }
925   ierr = VecDuplicateVecs(ts->vec_sol,s,&ros->Y);CHKERRQ(ierr);
926   ierr = VecDuplicate(ts->vec_sol,&ros->Ydot);CHKERRQ(ierr);
927   ierr = VecDuplicate(ts->vec_sol,&ros->Ystage);CHKERRQ(ierr);
928   ierr = VecDuplicate(ts->vec_sol,&ros->Zdot);CHKERRQ(ierr);
929   ierr = VecDuplicate(ts->vec_sol,&ros->Zstage);CHKERRQ(ierr);
930   ierr = PetscMalloc(s*sizeof(ros->work[0]),&ros->work);CHKERRQ(ierr);
931   PetscFunctionReturn(0);
932 }
933 /*------------------------------------------------------------*/
934 
935 #undef __FUNCT__
936 #define __FUNCT__ "TSSetFromOptions_RosW"
937 static PetscErrorCode TSSetFromOptions_RosW(TS ts)
938 {
939   TS_RosW        *ros = (TS_RosW*)ts->data;
940   PetscErrorCode ierr;
941   char           rostype[256];
942 
943   PetscFunctionBegin;
944   ierr = PetscOptionsHead("RosW ODE solver options");CHKERRQ(ierr);
945   {
946     RosWTableauLink link;
947     PetscInt count,choice;
948     PetscBool flg;
949     const char **namelist;
950     SNES snes;
951 
952     ierr = PetscStrncpy(rostype,TSRosWDefault,sizeof rostype);CHKERRQ(ierr);
953     for (link=RosWTableauList,count=0; link; link=link->next,count++) ;
954     ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr);
955     for (link=RosWTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
956     ierr = PetscOptionsEList("-ts_rosw_type","Family of Rosenbrock-W method","TSRosWSetType",(const char*const*)namelist,count,rostype,&choice,&flg);CHKERRQ(ierr);
957     ierr = TSRosWSetType(ts,flg ? namelist[choice] : rostype);CHKERRQ(ierr);
958     ierr = PetscFree(namelist);CHKERRQ(ierr);
959 
960     ierr = PetscOptionsBool("-ts_rosw_recompute_jacobian","Recompute the Jacobian at each stage","TSRosWSetRecomputeJacobian",ros->recompute_jacobian,&ros->recompute_jacobian,PETSC_NULL);CHKERRQ(ierr);
961 
962     /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
963     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
964     if (!((PetscObject)snes)->type_name) {
965       ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);
966     }
967     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
968   }
969   ierr = PetscOptionsTail();CHKERRQ(ierr);
970   PetscFunctionReturn(0);
971 }
972 
973 #undef __FUNCT__
974 #define __FUNCT__ "PetscFormatRealArray"
975 static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
976 {
977   PetscErrorCode ierr;
978   PetscInt i;
979   size_t left,count;
980   char *p;
981 
982   PetscFunctionBegin;
983   for (i=0,p=buf,left=len; i<n; i++) {
984     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
985     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
986     left -= count;
987     p += count;
988     *p++ = ' ';
989   }
990   p[i ? 0 : -1] = 0;
991   PetscFunctionReturn(0);
992 }
993 
994 #undef __FUNCT__
995 #define __FUNCT__ "TSView_RosW"
996 static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer)
997 {
998   TS_RosW        *ros = (TS_RosW*)ts->data;
999   RosWTableau    tab  = ros->tableau;
1000   PetscBool      iascii;
1001   PetscErrorCode ierr;
1002 
1003   PetscFunctionBegin;
1004   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1005   if (iascii) {
1006     const TSRosWType rostype;
1007     PetscInt i;
1008     PetscReal abscissa[512];
1009     char buf[512];
1010     ierr = TSRosWGetType(ts,&rostype);CHKERRQ(ierr);
1011     ierr = PetscViewerASCIIPrintf(viewer,"  Rosenbrock-W %s\n",rostype);CHKERRQ(ierr);
1012     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ASum);CHKERRQ(ierr);
1013     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A       = %s\n",buf);CHKERRQ(ierr);
1014     for (i=0; i<tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i];
1015     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,abscissa);CHKERRQ(ierr);
1016     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A+Gamma = %s\n",buf);CHKERRQ(ierr);
1017   }
1018   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
1019   PetscFunctionReturn(0);
1020 }
1021 
1022 #undef __FUNCT__
1023 #define __FUNCT__ "TSRosWSetType"
1024 /*@C
1025   TSRosWSetType - Set the type of Rosenbrock-W scheme
1026 
1027   Logically collective
1028 
1029   Input Parameter:
1030 +  ts - timestepping context
1031 -  rostype - type of Rosenbrock-W scheme
1032 
1033   Level: beginner
1034 
1035 .seealso: TSRosWGetType(), TSROSW, TSROSW2M, TSROSW2P, TSROSWRA3PW, TSROSWRA34PW2, TSROSWRODAS3, TSROSWSANDU3, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSROSWLLSSP3P4S2C, TSROSWARK3
1036 @*/
1037 PetscErrorCode TSRosWSetType(TS ts,const TSRosWType rostype)
1038 {
1039   PetscErrorCode ierr;
1040 
1041   PetscFunctionBegin;
1042   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1043   ierr = PetscTryMethod(ts,"TSRosWSetType_C",(TS,const TSRosWType),(ts,rostype));CHKERRQ(ierr);
1044   PetscFunctionReturn(0);
1045 }
1046 
1047 #undef __FUNCT__
1048 #define __FUNCT__ "TSRosWGetType"
1049 /*@C
1050   TSRosWGetType - Get the type of Rosenbrock-W scheme
1051 
1052   Logically collective
1053 
1054   Input Parameter:
1055 .  ts - timestepping context
1056 
1057   Output Parameter:
1058 .  rostype - type of Rosenbrock-W scheme
1059 
1060   Level: intermediate
1061 
1062 .seealso: TSRosWGetType()
1063 @*/
1064 PetscErrorCode TSRosWGetType(TS ts,const TSRosWType *rostype)
1065 {
1066   PetscErrorCode ierr;
1067 
1068   PetscFunctionBegin;
1069   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1070   ierr = PetscUseMethod(ts,"TSRosWGetType_C",(TS,const TSRosWType*),(ts,rostype));CHKERRQ(ierr);
1071   PetscFunctionReturn(0);
1072 }
1073 
1074 #undef __FUNCT__
1075 #define __FUNCT__ "TSRosWSetRecomputeJacobian"
1076 /*@C
1077   TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step.
1078 
1079   Logically collective
1080 
1081   Input Parameter:
1082 +  ts - timestepping context
1083 -  flg - PETSC_TRUE to recompute the Jacobian at each stage
1084 
1085   Level: intermediate
1086 
1087 .seealso: TSRosWGetType()
1088 @*/
1089 PetscErrorCode TSRosWSetRecomputeJacobian(TS ts,PetscBool flg)
1090 {
1091   PetscErrorCode ierr;
1092 
1093   PetscFunctionBegin;
1094   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1095   ierr = PetscTryMethod(ts,"TSRosWSetRecomputeJacobian_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
1096   PetscFunctionReturn(0);
1097 }
1098 
1099 EXTERN_C_BEGIN
1100 #undef __FUNCT__
1101 #define __FUNCT__ "TSRosWGetType_RosW"
1102 PetscErrorCode  TSRosWGetType_RosW(TS ts,const TSRosWType *rostype)
1103 {
1104   TS_RosW        *ros = (TS_RosW*)ts->data;
1105   PetscErrorCode ierr;
1106 
1107   PetscFunctionBegin;
1108   if (!ros->tableau) {ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);}
1109   *rostype = ros->tableau->name;
1110   PetscFunctionReturn(0);
1111 }
1112 #undef __FUNCT__
1113 #define __FUNCT__ "TSRosWSetType_RosW"
1114 PetscErrorCode  TSRosWSetType_RosW(TS ts,const TSRosWType rostype)
1115 {
1116   TS_RosW         *ros = (TS_RosW*)ts->data;
1117   PetscErrorCode  ierr;
1118   PetscBool       match;
1119   RosWTableauLink link;
1120 
1121   PetscFunctionBegin;
1122   if (ros->tableau) {
1123     ierr = PetscStrcmp(ros->tableau->name,rostype,&match);CHKERRQ(ierr);
1124     if (match) PetscFunctionReturn(0);
1125   }
1126   for (link = RosWTableauList; link; link=link->next) {
1127     ierr = PetscStrcmp(link->tab.name,rostype,&match);CHKERRQ(ierr);
1128     if (match) {
1129       ierr = TSReset_RosW(ts);CHKERRQ(ierr);
1130       ros->tableau = &link->tab;
1131       PetscFunctionReturn(0);
1132     }
1133   }
1134   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rostype);
1135   PetscFunctionReturn(0);
1136 }
1137 
1138 #undef __FUNCT__
1139 #define __FUNCT__ "TSRosWSetRecomputeJacobian_RosW"
1140 PetscErrorCode  TSRosWSetRecomputeJacobian_RosW(TS ts,PetscBool flg)
1141 {
1142   TS_RosW *ros = (TS_RosW*)ts->data;
1143 
1144   PetscFunctionBegin;
1145   ros->recompute_jacobian = flg;
1146   PetscFunctionReturn(0);
1147 }
1148 EXTERN_C_END
1149 
1150 /* ------------------------------------------------------------ */
1151 /*MC
1152       TSROSW - ODE solver using Rosenbrock-W schemes
1153 
1154   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
1155   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
1156   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
1157 
1158   Notes:
1159   This method currently only works with autonomous ODE and DAE.
1160 
1161   Developer notes:
1162   Rosenbrock-W methods are typically specified for autonomous ODE
1163 
1164 $  xdot = f(x)
1165 
1166   by the stage equations
1167 
1168 $  k_i = h f(x_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j
1169 
1170   and step completion formula
1171 
1172 $  x_1 = x_0 + sum_j b_j k_j
1173 
1174   with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(x)
1175   and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner,
1176   we define new variables for the stage equations
1177 
1178 $  y_i = gamma_ij k_j
1179 
1180   The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define
1181 
1182 $  A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-i}
1183 
1184   to rewrite the method as
1185 
1186 $  [M/(h gamma_ii) - J] y_i = f(x_0 + sum_j a_ij y_j) + M sum_j (c_ij/h) y_j
1187 $  x_1 = x_0 + sum_j bt_j y_j
1188 
1189    where we have introduced the mass matrix M. Continue by defining
1190 
1191 $  ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j
1192 
1193    or, more compactly in tensor notation
1194 
1195 $  Ydot = 1/h (Gamma^{-1} \otimes I) Y .
1196 
1197    Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current
1198    stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the
1199    equation
1200 
1201 $  g(x_0 + sum_j a_ij y_j + y_i, ydot_i) = 0
1202 
1203    with initial guess y_i = 0.
1204 
1205   Level: beginner
1206 
1207 .seealso:  TSCreate(), TS, TSSetType(), TSRosWSetType(), TSRosWRegister()
1208 
1209 M*/
1210 EXTERN_C_BEGIN
1211 #undef __FUNCT__
1212 #define __FUNCT__ "TSCreate_RosW"
1213 PetscErrorCode  TSCreate_RosW(TS ts)
1214 {
1215   TS_RosW        *ros;
1216   PetscErrorCode ierr;
1217 
1218   PetscFunctionBegin;
1219 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
1220   ierr = TSRosWInitializePackage(PETSC_NULL);CHKERRQ(ierr);
1221 #endif
1222 
1223   ts->ops->reset          = TSReset_RosW;
1224   ts->ops->destroy        = TSDestroy_RosW;
1225   ts->ops->view           = TSView_RosW;
1226   ts->ops->setup          = TSSetUp_RosW;
1227   ts->ops->step           = TSStep_RosW;
1228   ts->ops->interpolate    = TSInterpolate_RosW;
1229   ts->ops->evaluatestep   = TSEvaluateStep_RosW;
1230   ts->ops->setfromoptions = TSSetFromOptions_RosW;
1231   ts->ops->snesfunction   = SNESTSFormFunction_RosW;
1232   ts->ops->snesjacobian   = SNESTSFormJacobian_RosW;
1233 
1234   ierr = PetscNewLog(ts,TS_RosW,&ros);CHKERRQ(ierr);
1235   ts->data = (void*)ros;
1236 
1237   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","TSRosWGetType_RosW",TSRosWGetType_RosW);CHKERRQ(ierr);
1238   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","TSRosWSetType_RosW",TSRosWSetType_RosW);CHKERRQ(ierr);
1239   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","TSRosWSetRecomputeJacobian_RosW",TSRosWSetRecomputeJacobian_RosW);CHKERRQ(ierr);
1240   PetscFunctionReturn(0);
1241 }
1242 EXTERN_C_END
1243