xref: /petsc/src/ts/impls/rosw/rosw.c (revision 4489bed6aee10dbd846973bcb708f9fc25bfcf21) !
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   TSStepStatus status;
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[0][1]=0; Gamma[0][2]=0; Gamma[0][3]=0;
333    Gamma[1][0]=-1.997527830934941248426324674704153457289527280554476;
334    Gamma[1][1]=0.4358665215084589994160194475295062513822671686978816;
335    Gamma[1][2]=0; Gamma[1][3]=0;
336    Gamma[2][0]=-1.007948511795029620852002345345404191008352770119903;
337    Gamma[2][1]=-0.004648958462629345562774289390054679806993396798458131;
338    Gamma[2][2]=0.4358665215084589994160194475295062513822671686978816;
339    Gamma[2][3]=0;
340    Gamma[3][0]=-0.6685429734233467180451604600279552604364311322650783;
341    Gamma[3][1]=0.6056625986449338476089525334450053439525178740492984;
342    Gamma[3][2]=-0.9717899277217721234705114616271378792182450260943198;
343    Gamma[3][3]=0;
344 
345    A[0][0]=0; A[0][1]=0; A[0][2]=0; A[0][3]=0;
346    A[1][0]=0.8717330430169179988320388950590125027645343373957631;
347    A[1][1]=0; A[1][2]=0; A[1][3]=0;
348    A[2][0]=0.5275890119763004115618079766722914408876108660811028;
349    A[2][1]=0.07241098802369958843819203208518599088698057726988732;
350    A[2][2]=0; A[2][3]=0;
351    A[3][0]=0.3990960076760701320627260685975778145384666450351314;
352    A[3][1]=-0.4375576546135194437228463747348862825846903771419953;
353    A[3][2]=1.038461646937449311660120300601880176655352737312713;
354    A[3][3]=0;
355 
356    b[0]=0.1876410243467238251612921333138006734899663569186926;
357    b[1]=-0.5952974735769549480478230473706443582188442040780541;
358    b[2]=0.9717899277217721234705114616271378792182450260943198;
359    b[3]=0.4358665215084589994160194475295062513822671686978816;
360 
361    b2[0]=0.2147402862233891404862383521089097657790734483804460;
362    b2[1]=-0.4851622638849390928209050538171743017757490232519684;
363    b2[2]=0.8687250025203875511662123688667549217531982787600080;
364    b2[3]=0.4016969751411624011684543450940068201770721128357014;
365 
366    ierr = TSRosWRegister(TSROSWARK3,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
367   }
368 
369   PetscFunctionReturn(0);
370 }
371 
372 #undef __FUNCT__
373 #define __FUNCT__ "TSRosWRegisterDestroy"
374 /*@C
375    TSRosWRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister().
376 
377    Not Collective
378 
379    Level: advanced
380 
381 .keywords: TSRosW, register, destroy
382 .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegisterDynamic()
383 @*/
384 PetscErrorCode TSRosWRegisterDestroy(void)
385 {
386   PetscErrorCode ierr;
387   RosWTableauLink link;
388 
389   PetscFunctionBegin;
390   while ((link = RosWTableauList)) {
391     RosWTableau t = &link->tab;
392     RosWTableauList = link->next;
393     ierr = PetscFree5(t->A,t->Gamma,t->b,t->ASum,t->GammaSum);CHKERRQ(ierr);
394     ierr = PetscFree4(t->At,t->bt,t->GammaInv,t->GammaZeroDiag);CHKERRQ(ierr);
395     ierr = PetscFree2(t->bembed,t->bembedt);CHKERRQ(ierr);
396     ierr = PetscFree(t->name);CHKERRQ(ierr);
397     ierr = PetscFree(link);CHKERRQ(ierr);
398   }
399   TSRosWRegisterAllCalled = PETSC_FALSE;
400   PetscFunctionReturn(0);
401 }
402 
403 #undef __FUNCT__
404 #define __FUNCT__ "TSRosWInitializePackage"
405 /*@C
406   TSRosWInitializePackage - This function initializes everything in the TSRosW package. It is called
407   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RosW()
408   when using static libraries.
409 
410   Input Parameter:
411   path - The dynamic library path, or PETSC_NULL
412 
413   Level: developer
414 
415 .keywords: TS, TSRosW, initialize, package
416 .seealso: PetscInitialize()
417 @*/
418 PetscErrorCode TSRosWInitializePackage(const char path[])
419 {
420   PetscErrorCode ierr;
421 
422   PetscFunctionBegin;
423   if (TSRosWPackageInitialized) PetscFunctionReturn(0);
424   TSRosWPackageInitialized = PETSC_TRUE;
425   ierr = TSRosWRegisterAll();CHKERRQ(ierr);
426   ierr = PetscRegisterFinalize(TSRosWFinalizePackage);CHKERRQ(ierr);
427   PetscFunctionReturn(0);
428 }
429 
430 #undef __FUNCT__
431 #define __FUNCT__ "TSRosWFinalizePackage"
432 /*@C
433   TSRosWFinalizePackage - This function destroys everything in the TSRosW package. It is
434   called from PetscFinalize().
435 
436   Level: developer
437 
438 .keywords: Petsc, destroy, package
439 .seealso: PetscFinalize()
440 @*/
441 PetscErrorCode TSRosWFinalizePackage(void)
442 {
443   PetscErrorCode ierr;
444 
445   PetscFunctionBegin;
446   TSRosWPackageInitialized = PETSC_FALSE;
447   ierr = TSRosWRegisterDestroy();CHKERRQ(ierr);
448   PetscFunctionReturn(0);
449 }
450 
451 #undef __FUNCT__
452 #define __FUNCT__ "TSRosWRegister"
453 /*@C
454    TSRosWRegister - register a Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
455 
456    Not Collective, but the same schemes should be registered on all processes on which they will be used
457 
458    Input Parameters:
459 +  name - identifier for method
460 .  order - approximation order of method
461 .  s - number of stages, this is the dimension of the matrices below
462 .  A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular
463 .  Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal
464 .  b - Step completion table (dimension s)
465 -  bembed - Step completion table for a scheme of order one less (dimension s, PETSC_NULL if no embedded scheme is available)
466 
467    Notes:
468    Several Rosenbrock W methods are provided, this function is only needed to create new methods.
469 
470    Level: advanced
471 
472 .keywords: TS, register
473 
474 .seealso: TSRosW
475 @*/
476 PetscErrorCode TSRosWRegister(const TSRosWType name,PetscInt order,PetscInt s,
477                               const PetscReal A[],const PetscReal Gamma[],const PetscReal b[],const PetscReal bembed[])
478 {
479   PetscErrorCode ierr;
480   RosWTableauLink link;
481   RosWTableau t;
482   PetscInt i,j,k;
483   PetscScalar *GammaInv;
484 
485   PetscFunctionBegin;
486   PetscValidCharPointer(name,1);
487   PetscValidPointer(A,4);
488   PetscValidPointer(Gamma,5);
489   PetscValidPointer(b,6);
490   if (bembed) PetscValidPointer(bembed,7);
491 
492   ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
493   ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
494   t = &link->tab;
495   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
496   t->order = order;
497   t->s = s;
498   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);
499   ierr = PetscMalloc4(s*s,PetscReal,&t->At,s,PetscReal,&t->bt,s*s,PetscReal,&t->GammaInv,s,PetscBool,&t->GammaZeroDiag);CHKERRQ(ierr);
500   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
501   ierr = PetscMemcpy(t->Gamma,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr);
502   ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
503   if (bembed) {
504     ierr = PetscMalloc2(s,PetscReal,&t->bembed,s,PetscReal,&t->bembedt);CHKERRQ(ierr);
505     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
506   }
507   for (i=0; i<s; i++) {
508     t->ASum[i] = 0;
509     t->GammaSum[i] = 0;
510     for (j=0; j<s; j++) {
511       t->ASum[i] += A[i*s+j];
512       t->GammaSum[i] += Gamma[i*s+j];
513     }
514   }
515   ierr = PetscMalloc(s*s*sizeof(PetscScalar),&GammaInv);CHKERRQ(ierr); /* Need to use Scalar for inverse, then convert back to Real */
516   for (i=0; i<s*s; i++) GammaInv[i] = Gamma[i];
517   for (i=0; i<s; i++) {
518     if (Gamma[i*s+i] == 0.0) {
519       GammaInv[i*s+i] = 1.0;
520       t->GammaZeroDiag[i] = PETSC_TRUE;
521     } else {
522       t->GammaZeroDiag[i] = PETSC_FALSE;
523     }
524   }
525 
526   switch (s) {
527   case 1: GammaInv[0] = 1./GammaInv[0]; break;
528   case 2: ierr = Kernel_A_gets_inverse_A_2(GammaInv,0);CHKERRQ(ierr); break;
529   case 3: ierr = Kernel_A_gets_inverse_A_3(GammaInv,0);CHKERRQ(ierr); break;
530   case 4: ierr = Kernel_A_gets_inverse_A_4(GammaInv,0);CHKERRQ(ierr); break;
531   case 5: {
532     PetscInt ipvt5[5];
533     MatScalar work5[5*5];
534     ierr = Kernel_A_gets_inverse_A_5(GammaInv,ipvt5,work5,0);CHKERRQ(ierr); break;
535   }
536   case 6: ierr = Kernel_A_gets_inverse_A_6(GammaInv,0);CHKERRQ(ierr); break;
537   case 7: ierr = Kernel_A_gets_inverse_A_7(GammaInv,0);CHKERRQ(ierr); break;
538   default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for %D stages",s);
539   }
540   for (i=0; i<s*s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]);
541   ierr = PetscFree(GammaInv);CHKERRQ(ierr);
542   for (i=0; i<s; i++) {
543     for (j=0; j<s; j++) {
544       t->At[i*s+j] = 0;
545       for (k=0; k<s; k++) {
546         t->At[i*s+j] += t->A[i*s+k] * t->GammaInv[k*s+j];
547       }
548     }
549     t->bt[i] = 0;
550     for (j=0; j<s; j++) {
551       t->bt[i] += t->b[j] * t->GammaInv[j*s+i];
552     }
553     if (bembed) {
554       t->bembedt[i] = 0;
555       for (j=0; j<s; j++) {
556         t->bembedt[i] += t->bembed[j] * t->GammaInv[j*s+i];
557       }
558     }
559   }
560   t->ccfl = 1.0;                /* Fix this */
561 
562   link->next = RosWTableauList;
563   RosWTableauList = link;
564   PetscFunctionReturn(0);
565 }
566 
567 #undef __FUNCT__
568 #define __FUNCT__ "TSEvaluateStep_RosW"
569 /*
570  The step completion formula is
571 
572  x1 = x0 + b^T Y
573 
574  where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been
575  updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write
576 
577  x1e = x0 + be^T Y
578      = x1 - b^T Y + be^T Y
579      = x1 + (be - b)^T Y
580 
581  so we can evaluate the method of different order even after the step has been optimistically completed.
582 */
583 static PetscErrorCode TSEvaluateStep_RosW(TS ts,PetscInt order,Vec X,PetscBool *done)
584 {
585   TS_RosW        *ros = (TS_RosW*)ts->data;
586   RosWTableau    tab  = ros->tableau;
587   PetscScalar    *w = ros->work;
588   PetscInt       i;
589   PetscErrorCode ierr;
590 
591   PetscFunctionBegin;
592   if (order == tab->order) {
593     if (ros->status == TS_STEP_INCOMPLETE) { /* Use standard completion formula */
594       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
595       for (i=0; i<tab->s; i++) w[i] = tab->bt[i];
596       ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr);
597     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
598     if (done) *done = PETSC_TRUE;
599     PetscFunctionReturn(0);
600   } else if (order == tab->order-1) {
601     if (!tab->bembedt) goto unavailable;
602     if (ros->status == TS_STEP_INCOMPLETE) { /* Use embedded completion formula */
603       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
604       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i];
605       ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr);
606     } else {                    /* Use rollback-and-recomplete formula (bembedt - bt) */
607       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i];
608       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
609       ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr);
610     }
611     if (done) *done = PETSC_TRUE;
612     PetscFunctionReturn(0);
613   }
614   unavailable:
615   if (done) *done = PETSC_FALSE;
616   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);
617   PetscFunctionReturn(0);
618 }
619 
620 #undef __FUNCT__
621 #define __FUNCT__ "TSStep_RosW"
622 static PetscErrorCode TSStep_RosW(TS ts)
623 {
624   TS_RosW         *ros = (TS_RosW*)ts->data;
625   RosWTableau     tab  = ros->tableau;
626   const PetscInt  s    = tab->s;
627   const PetscReal *At  = tab->At,*Gamma = tab->Gamma,*ASum = tab->ASum,*GammaInv = tab->GammaInv;
628   const PetscBool *GammaZeroDiag = tab->GammaZeroDiag;
629   PetscScalar     *w   = ros->work;
630   Vec             *Y   = ros->Y,Ydot = ros->Ydot,Zdot = ros->Zdot,Zstage = ros->Zstage;
631   SNES            snes;
632   TSAdapt         adapt;
633   PetscInt        i,j,its,lits,reject,next_scheme;
634   PetscReal       next_time_step;
635   PetscBool       accept;
636   PetscErrorCode  ierr;
637 
638   PetscFunctionBegin;
639   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
640   next_time_step = ts->time_step;
641   accept = PETSC_TRUE;
642   ros->status = TS_STEP_INCOMPLETE;
643 
644   for (reject=0; reject<ts->max_reject; reject++,ts->reject++) {
645     const PetscReal h = ts->time_step;
646     for (i=0; i<s; i++) {
647       ros->stage_time = ts->ptime + h*ASum[i];
648       if (GammaZeroDiag[i]) {
649         ros->stage_explicit = PETSC_TRUE;
650         ros->shift = 1./h;
651       } else {
652         ros->stage_explicit = PETSC_FALSE;
653         ros->shift = 1./(h*Gamma[i*s+i]);
654       }
655 
656       ierr = VecCopy(ts->vec_sol,Zstage);CHKERRQ(ierr);
657       for (j=0; j<i; j++) w[j] = At[i*s+j];
658       ierr = VecMAXPY(Zstage,i,w,Y);CHKERRQ(ierr);
659 
660       for (j=0; j<i; j++) w[j] = 1./h * GammaInv[i*s+j];
661       ierr = VecZeroEntries(Zdot);CHKERRQ(ierr);
662       ierr = VecMAXPY(Zdot,i,w,Y);CHKERRQ(ierr);
663 
664       /* Initial guess taken from last stage */
665       ierr = VecZeroEntries(Y[i]);CHKERRQ(ierr);
666 
667       if (!ros->stage_explicit) {
668         if (!ros->recompute_jacobian && !i) {
669           ierr = SNESSetLagJacobian(snes,-2);CHKERRQ(ierr); /* Recompute the Jacobian on this solve, but not again */
670         }
671         ierr = SNESSolve(snes,PETSC_NULL,Y[i]);CHKERRQ(ierr);
672         ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
673         ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
674         ts->nonlinear_its += its; ts->linear_its += lits;
675       } else {
676         ierr = VecWAXPY(Ydot,1,ts->vec_sol,Zdot);CHKERRQ(ierr); /* Ydot = x0 + Zdot */
677         ierr = TSComputeIFunction(ts,ros->stage_time,ros->Ystage,Ydot,Zdot,PETSC_FALSE);CHKERRQ(ierr);
678         ierr = VecWAXPY(ros->Ystage,1.0,Zdot,ros->Zstage);CHKERRQ(ierr);    /* Ystage = F + Zstage */
679         ts->linear_its += 1;
680       }
681     }
682     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
683     ros->status = TS_STEP_PENDING;
684 
685     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
686     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
687     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
688     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
689     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
690     if (accept) {
691       /* ignore next_scheme for now */
692       ts->ptime += ts->time_step;
693       ts->time_step = next_time_step;
694       ts->steps++;
695       ros->status = TS_STEP_COMPLETE;
696       break;
697     } else {                    /* Roll back the current step */
698       for (i=0; i<s; i++) w[i] = -tab->bt[i];
699       ierr = VecMAXPY(ts->vec_sol,s,w,Y);CHKERRQ(ierr);
700       ts->time_step = next_time_step;
701       ros->status = TS_STEP_INCOMPLETE;
702     }
703   }
704   PetscFunctionReturn(0);
705 }
706 
707 #undef __FUNCT__
708 #define __FUNCT__ "TSInterpolate_RosW"
709 static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec X)
710 {
711   TS_RosW        *ros = (TS_RosW*)ts->data;
712 
713   PetscFunctionBegin;
714   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ros->tableau->name);
715   PetscFunctionReturn(0);
716 }
717 
718 /*------------------------------------------------------------*/
719 #undef __FUNCT__
720 #define __FUNCT__ "TSReset_RosW"
721 static PetscErrorCode TSReset_RosW(TS ts)
722 {
723   TS_RosW        *ros = (TS_RosW*)ts->data;
724   PetscInt       s;
725   PetscErrorCode ierr;
726 
727   PetscFunctionBegin;
728   if (!ros->tableau) PetscFunctionReturn(0);
729    s = ros->tableau->s;
730   ierr = VecDestroyVecs(s,&ros->Y);CHKERRQ(ierr);
731   ierr = VecDestroy(&ros->Ydot);CHKERRQ(ierr);
732   ierr = VecDestroy(&ros->Ystage);CHKERRQ(ierr);
733   ierr = VecDestroy(&ros->Zdot);CHKERRQ(ierr);
734   ierr = VecDestroy(&ros->Zstage);CHKERRQ(ierr);
735   ierr = PetscFree(ros->work);CHKERRQ(ierr);
736   PetscFunctionReturn(0);
737 }
738 
739 #undef __FUNCT__
740 #define __FUNCT__ "TSDestroy_RosW"
741 static PetscErrorCode TSDestroy_RosW(TS ts)
742 {
743   PetscErrorCode  ierr;
744 
745   PetscFunctionBegin;
746   ierr = TSReset_RosW(ts);CHKERRQ(ierr);
747   ierr = PetscFree(ts->data);CHKERRQ(ierr);
748   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","",PETSC_NULL);CHKERRQ(ierr);
749   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","",PETSC_NULL);CHKERRQ(ierr);
750   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","",PETSC_NULL);CHKERRQ(ierr);
751   PetscFunctionReturn(0);
752 }
753 
754 /*
755   This defines the nonlinear equation that is to be solved with SNES
756   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
757 */
758 #undef __FUNCT__
759 #define __FUNCT__ "SNESTSFormFunction_RosW"
760 static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec X,Vec F,TS ts)
761 {
762   TS_RosW        *ros = (TS_RosW*)ts->data;
763   PetscErrorCode ierr;
764 
765   PetscFunctionBegin;
766   ierr = VecWAXPY(ros->Ydot,ros->shift,X,ros->Zdot);CHKERRQ(ierr); /* Ydot = shift*X + Zdot */
767   ierr = VecWAXPY(ros->Ystage,1.0,X,ros->Zstage);CHKERRQ(ierr);    /* Ystage = X + Zstage */
768   ierr = TSComputeIFunction(ts,ros->stage_time,ros->Ystage,ros->Ydot,F,PETSC_FALSE);CHKERRQ(ierr);
769   PetscFunctionReturn(0);
770 }
771 
772 #undef __FUNCT__
773 #define __FUNCT__ "SNESTSFormJacobian_RosW"
774 static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts)
775 {
776   TS_RosW        *ros = (TS_RosW*)ts->data;
777   PetscErrorCode ierr;
778 
779   PetscFunctionBegin;
780   /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */
781   ierr = TSComputeIJacobian(ts,ros->stage_time,ros->Ystage,ros->Ydot,ros->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr);
782   PetscFunctionReturn(0);
783 }
784 
785 #undef __FUNCT__
786 #define __FUNCT__ "TSSetUp_RosW"
787 static PetscErrorCode TSSetUp_RosW(TS ts)
788 {
789   TS_RosW        *ros = (TS_RosW*)ts->data;
790   RosWTableau    tab  = ros->tableau;
791   PetscInt       s    = tab->s;
792   PetscErrorCode ierr;
793 
794   PetscFunctionBegin;
795   if (!ros->tableau) {
796     ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);
797   }
798   ierr = VecDuplicateVecs(ts->vec_sol,s,&ros->Y);CHKERRQ(ierr);
799   ierr = VecDuplicate(ts->vec_sol,&ros->Ydot);CHKERRQ(ierr);
800   ierr = VecDuplicate(ts->vec_sol,&ros->Ystage);CHKERRQ(ierr);
801   ierr = VecDuplicate(ts->vec_sol,&ros->Zdot);CHKERRQ(ierr);
802   ierr = VecDuplicate(ts->vec_sol,&ros->Zstage);CHKERRQ(ierr);
803   ierr = PetscMalloc(s*sizeof(ros->work[0]),&ros->work);CHKERRQ(ierr);
804   PetscFunctionReturn(0);
805 }
806 /*------------------------------------------------------------*/
807 
808 #undef __FUNCT__
809 #define __FUNCT__ "TSSetFromOptions_RosW"
810 static PetscErrorCode TSSetFromOptions_RosW(TS ts)
811 {
812   TS_RosW        *ros = (TS_RosW*)ts->data;
813   PetscErrorCode ierr;
814   char           rostype[256];
815 
816   PetscFunctionBegin;
817   ierr = PetscOptionsHead("RosW ODE solver options");CHKERRQ(ierr);
818   {
819     RosWTableauLink link;
820     PetscInt count,choice;
821     PetscBool flg;
822     const char **namelist;
823     SNES snes;
824 
825     ierr = PetscStrncpy(rostype,TSRosWDefault,sizeof rostype);CHKERRQ(ierr);
826     for (link=RosWTableauList,count=0; link; link=link->next,count++) ;
827     ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr);
828     for (link=RosWTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
829     ierr = PetscOptionsEList("-ts_rosw_type","Family of Rosenbrock-W method","TSRosWSetType",(const char*const*)namelist,count,rostype,&choice,&flg);CHKERRQ(ierr);
830     ierr = TSRosWSetType(ts,flg ? namelist[choice] : rostype);CHKERRQ(ierr);
831     ierr = PetscFree(namelist);CHKERRQ(ierr);
832 
833     ierr = PetscOptionsBool("-ts_rosw_recompute_jacobian","Recompute the Jacobian at each stage","TSRosWSetRecomputeJacobian",ros->recompute_jacobian,&ros->recompute_jacobian,PETSC_NULL);CHKERRQ(ierr);
834 
835     /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
836     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
837     if (!((PetscObject)snes)->type_name) {
838       ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);
839     }
840     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
841   }
842   ierr = PetscOptionsTail();CHKERRQ(ierr);
843   PetscFunctionReturn(0);
844 }
845 
846 #undef __FUNCT__
847 #define __FUNCT__ "PetscFormatRealArray"
848 static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
849 {
850   PetscErrorCode ierr;
851   PetscInt i;
852   size_t left,count;
853   char *p;
854 
855   PetscFunctionBegin;
856   for (i=0,p=buf,left=len; i<n; i++) {
857     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
858     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
859     left -= count;
860     p += count;
861     *p++ = ' ';
862   }
863   p[i ? 0 : -1] = 0;
864   PetscFunctionReturn(0);
865 }
866 
867 #undef __FUNCT__
868 #define __FUNCT__ "TSView_RosW"
869 static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer)
870 {
871   TS_RosW        *ros = (TS_RosW*)ts->data;
872   RosWTableau    tab  = ros->tableau;
873   PetscBool      iascii;
874   PetscErrorCode ierr;
875 
876   PetscFunctionBegin;
877   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
878   if (iascii) {
879     const TSRosWType rostype;
880     PetscInt i;
881     PetscReal abscissa[512];
882     char buf[512];
883     ierr = TSRosWGetType(ts,&rostype);CHKERRQ(ierr);
884     ierr = PetscViewerASCIIPrintf(viewer,"  Rosenbrock-W %s\n",rostype);CHKERRQ(ierr);
885     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ASum);CHKERRQ(ierr);
886     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A       = %s\n",buf);CHKERRQ(ierr);
887     for (i=0; i<tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i];
888     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,abscissa);CHKERRQ(ierr);
889     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A+Gamma = %s\n",buf);CHKERRQ(ierr);
890   }
891   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
892   PetscFunctionReturn(0);
893 }
894 
895 #undef __FUNCT__
896 #define __FUNCT__ "TSRosWSetType"
897 /*@C
898   TSRosWSetType - Set the type of Rosenbrock-W scheme
899 
900   Logically collective
901 
902   Input Parameter:
903 +  ts - timestepping context
904 -  rostype - type of Rosenbrock-W scheme
905 
906   Level: intermediate
907 
908 .seealso: TSRosWGetType()
909 @*/
910 PetscErrorCode TSRosWSetType(TS ts,const TSRosWType rostype)
911 {
912   PetscErrorCode ierr;
913 
914   PetscFunctionBegin;
915   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
916   ierr = PetscTryMethod(ts,"TSRosWSetType_C",(TS,const TSRosWType),(ts,rostype));CHKERRQ(ierr);
917   PetscFunctionReturn(0);
918 }
919 
920 #undef __FUNCT__
921 #define __FUNCT__ "TSRosWGetType"
922 /*@C
923   TSRosWGetType - Get the type of Rosenbrock-W scheme
924 
925   Logically collective
926 
927   Input Parameter:
928 .  ts - timestepping context
929 
930   Output Parameter:
931 .  rostype - type of Rosenbrock-W scheme
932 
933   Level: intermediate
934 
935 .seealso: TSRosWGetType()
936 @*/
937 PetscErrorCode TSRosWGetType(TS ts,const TSRosWType *rostype)
938 {
939   PetscErrorCode ierr;
940 
941   PetscFunctionBegin;
942   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
943   ierr = PetscUseMethod(ts,"TSRosWGetType_C",(TS,const TSRosWType*),(ts,rostype));CHKERRQ(ierr);
944   PetscFunctionReturn(0);
945 }
946 
947 #undef __FUNCT__
948 #define __FUNCT__ "TSRosWSetRecomputeJacobian"
949 /*@C
950   TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step.
951 
952   Logically collective
953 
954   Input Parameter:
955 +  ts - timestepping context
956 -  flg - PETSC_TRUE to recompute the Jacobian at each stage
957 
958   Level: intermediate
959 
960 .seealso: TSRosWGetType()
961 @*/
962 PetscErrorCode TSRosWSetRecomputeJacobian(TS ts,PetscBool flg)
963 {
964   PetscErrorCode ierr;
965 
966   PetscFunctionBegin;
967   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
968   ierr = PetscTryMethod(ts,"TSRosWSetRecomputeJacobian_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
969   PetscFunctionReturn(0);
970 }
971 
972 EXTERN_C_BEGIN
973 #undef __FUNCT__
974 #define __FUNCT__ "TSRosWGetType_RosW"
975 PetscErrorCode  TSRosWGetType_RosW(TS ts,const TSRosWType *rostype)
976 {
977   TS_RosW        *ros = (TS_RosW*)ts->data;
978   PetscErrorCode ierr;
979 
980   PetscFunctionBegin;
981   if (!ros->tableau) {ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);}
982   *rostype = ros->tableau->name;
983   PetscFunctionReturn(0);
984 }
985 #undef __FUNCT__
986 #define __FUNCT__ "TSRosWSetType_RosW"
987 PetscErrorCode  TSRosWSetType_RosW(TS ts,const TSRosWType rostype)
988 {
989   TS_RosW         *ros = (TS_RosW*)ts->data;
990   PetscErrorCode  ierr;
991   PetscBool       match;
992   RosWTableauLink link;
993 
994   PetscFunctionBegin;
995   if (ros->tableau) {
996     ierr = PetscStrcmp(ros->tableau->name,rostype,&match);CHKERRQ(ierr);
997     if (match) PetscFunctionReturn(0);
998   }
999   for (link = RosWTableauList; link; link=link->next) {
1000     ierr = PetscStrcmp(link->tab.name,rostype,&match);CHKERRQ(ierr);
1001     if (match) {
1002       ierr = TSReset_RosW(ts);CHKERRQ(ierr);
1003       ros->tableau = &link->tab;
1004       PetscFunctionReturn(0);
1005     }
1006   }
1007   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rostype);
1008   PetscFunctionReturn(0);
1009 }
1010 
1011 #undef __FUNCT__
1012 #define __FUNCT__ "TSRosWSetRecomputeJacobian_RosW"
1013 PetscErrorCode  TSRosWSetRecomputeJacobian_RosW(TS ts,PetscBool flg)
1014 {
1015   TS_RosW *ros = (TS_RosW*)ts->data;
1016 
1017   PetscFunctionBegin;
1018   ros->recompute_jacobian = flg;
1019   PetscFunctionReturn(0);
1020 }
1021 EXTERN_C_END
1022 
1023 /* ------------------------------------------------------------ */
1024 /*MC
1025       TSRosW - ODE solver using Rosenbrock-W schemes
1026 
1027   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
1028   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
1029   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
1030 
1031   Notes:
1032   This method currently only works with autonomous ODE and DAE.
1033 
1034   Developer notes:
1035   Rosenbrock-W methods are typically specified for autonomous ODE
1036 
1037 $  xdot = f(x)
1038 
1039   by the stage equations
1040 
1041 $  k_i = h f(x_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j
1042 
1043   and step completion formula
1044 
1045 $  x_1 = x_0 + sum_j b_j k_j
1046 
1047   with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(x)
1048   and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner,
1049   we define new variables for the stage equations
1050 
1051 $  y_i = gamma_ij k_j
1052 
1053   The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define
1054 
1055 $  A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-i}
1056 
1057   to rewrite the method as
1058 
1059 $  [M/(h gamma_ii) - J] y_i = f(x_0 + sum_j a_ij y_j) + M sum_j (c_ij/h) y_j
1060 $  x_1 = x_0 + sum_j bt_j y_j
1061 
1062    where we have introduced the mass matrix M. Continue by defining
1063 
1064 $  ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j
1065 
1066    or, more compactly in tensor notation
1067 
1068 $  Ydot = 1/h (Gamma^{-1} \otimes I) Y .
1069 
1070    Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current
1071    stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the
1072    equation
1073 
1074 $  g(x_0 + sum_j a_ij y_j + y_i, ydot_i) = 0
1075 
1076    with initial guess y_i = 0.
1077 
1078   Level: beginner
1079 
1080 .seealso:  TSCreate(), TS, TSSetType(), TSRosWRegister()
1081 
1082 M*/
1083 EXTERN_C_BEGIN
1084 #undef __FUNCT__
1085 #define __FUNCT__ "TSCreate_RosW"
1086 PetscErrorCode  TSCreate_RosW(TS ts)
1087 {
1088   TS_RosW        *ros;
1089   PetscErrorCode ierr;
1090 
1091   PetscFunctionBegin;
1092 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
1093   ierr = TSRosWInitializePackage(PETSC_NULL);CHKERRQ(ierr);
1094 #endif
1095 
1096   ts->ops->reset          = TSReset_RosW;
1097   ts->ops->destroy        = TSDestroy_RosW;
1098   ts->ops->view           = TSView_RosW;
1099   ts->ops->setup          = TSSetUp_RosW;
1100   ts->ops->step           = TSStep_RosW;
1101   ts->ops->interpolate    = TSInterpolate_RosW;
1102   ts->ops->evaluatestep   = TSEvaluateStep_RosW;
1103   ts->ops->setfromoptions = TSSetFromOptions_RosW;
1104   ts->ops->snesfunction   = SNESTSFormFunction_RosW;
1105   ts->ops->snesjacobian   = SNESTSFormJacobian_RosW;
1106 
1107   ierr = PetscNewLog(ts,TS_RosW,&ros);CHKERRQ(ierr);
1108   ts->data = (void*)ros;
1109 
1110   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","TSRosWGetType_RosW",TSRosWGetType_RosW);CHKERRQ(ierr);
1111   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","TSRosWSetType_RosW",TSRosWSetType_RosW);CHKERRQ(ierr);
1112   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","TSRosWSetRecomputeJacobian_RosW",TSRosWSetRecomputeJacobian_RosW);CHKERRQ(ierr);
1113   PetscFunctionReturn(0);
1114 }
1115 EXTERN_C_END
1116