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