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