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