xref: /petsc/src/ts/impls/rosw/rosw.c (revision 7d4bf2def842f67c84f7ff7d6f813bbb229a818f)
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[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->step_taken) {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
594     else {
595       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
596       for (i=0; i<tab->s; i++) w[i] = tab->bt[i];
597       ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr);
598     }
599     if (done) *done = PETSC_TRUE;
600     PetscFunctionReturn(0);
601   } else if (order == tab->order-1) {
602     if (!tab->bembedt) goto unavailable;
603     if (ros->step_taken) {
604       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i];
605       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
606       ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr);
607     } else {
608       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
609       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i];
610       ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr);
611     }
612     if (done) *done = PETSC_TRUE;
613     PetscFunctionReturn(0);
614   }
615   unavailable:
616   if (done) *done = PETSC_FALSE;
617   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);
618   PetscFunctionReturn(0);
619 }
620 
621 #undef __FUNCT__
622 #define __FUNCT__ "TSStep_RosW"
623 static PetscErrorCode TSStep_RosW(TS ts)
624 {
625   TS_RosW         *ros = (TS_RosW*)ts->data;
626   RosWTableau     tab  = ros->tableau;
627   const PetscInt  s    = tab->s;
628   const PetscReal *At  = tab->At,*Gamma = tab->Gamma,*ASum = tab->ASum,*GammaInv = tab->GammaInv;
629   const PetscBool *GammaZeroDiag = tab->GammaZeroDiag;
630   PetscScalar     *w   = ros->work;
631   Vec             *Y   = ros->Y,Ydot = ros->Ydot,Zdot = ros->Zdot,Zstage = ros->Zstage;
632   SNES            snes;
633   TSAdapt         adapt;
634   PetscInt        i,j,its,lits,reject,next_scheme;
635   PetscReal       next_time_step;
636   PetscBool       accept;
637   PetscErrorCode  ierr;
638 
639   PetscFunctionBegin;
640   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
641   next_time_step = ts->time_step;
642   accept = PETSC_TRUE;
643   ros->step_taken = PETSC_FALSE;
644 
645   for (reject=0; reject<ts->max_reject; reject++,ts->reject++) {
646     const PetscReal h = ts->time_step;
647     for (i=0; i<s; i++) {
648       ros->stage_time = ts->ptime + h*ASum[i];
649       if (GammaZeroDiag[i]) {
650         ros->stage_explicit = PETSC_TRUE;
651         ros->shift = 1./h;
652       } else {
653         ros->stage_explicit = PETSC_FALSE;
654         ros->shift = 1./(h*Gamma[i*s+i]);
655       }
656 
657       ierr = VecCopy(ts->vec_sol,Zstage);CHKERRQ(ierr);
658       for (j=0; j<i; j++) w[j] = At[i*s+j];
659       ierr = VecMAXPY(Zstage,i,w,Y);CHKERRQ(ierr);
660 
661       for (j=0; j<i; j++) w[j] = 1./h * GammaInv[i*s+j];
662       ierr = VecZeroEntries(Zdot);CHKERRQ(ierr);
663       ierr = VecMAXPY(Zdot,i,w,Y);CHKERRQ(ierr);
664 
665       /* Initial guess taken from last stage */
666       ierr = VecZeroEntries(Y[i]);CHKERRQ(ierr);
667 
668       if (!ros->stage_explicit) {
669         if (!ros->recompute_jacobian && !i) {
670           ierr = SNESSetLagJacobian(snes,-2);CHKERRQ(ierr); /* Recompute the Jacobian on this solve, but not again */
671         }
672         ierr = SNESSolve(snes,PETSC_NULL,Y[i]);CHKERRQ(ierr);
673         ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
674         ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
675         ts->nonlinear_its += its; ts->linear_its += lits;
676       } else {
677         ierr = VecWAXPY(Ydot,1,ts->vec_sol,Zdot);CHKERRQ(ierr); /* Ydot = x0 + Zdot */
678         ierr = TSComputeIFunction(ts,ros->stage_time,ros->Ystage,Ydot,Zdot,PETSC_FALSE);CHKERRQ(ierr);
679         ierr = VecWAXPY(ros->Ystage,1.0,Zdot,ros->Zstage);CHKERRQ(ierr);    /* Ystage = F + Zstage */
680         ts->linear_its += 1;
681       }
682     }
683     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
684     ros->step_taken = PETSC_TRUE;
685 
686     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
687     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
688     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
689     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
690     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
691     if (accept) {
692       /* ignore next_scheme for now */
693       ts->ptime += ts->time_step;
694       ts->time_step = next_time_step;
695       ts->steps++;
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->step_taken = PETSC_FALSE;
702     }
703   }
704 
705   PetscFunctionReturn(0);
706 }
707 
708 #undef __FUNCT__
709 #define __FUNCT__ "TSInterpolate_RosW"
710 static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec X)
711 {
712   TS_RosW        *ros = (TS_RosW*)ts->data;
713 
714   PetscFunctionBegin;
715   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ros->tableau->name);
716   PetscFunctionReturn(0);
717 }
718 
719 /*------------------------------------------------------------*/
720 #undef __FUNCT__
721 #define __FUNCT__ "TSReset_RosW"
722 static PetscErrorCode TSReset_RosW(TS ts)
723 {
724   TS_RosW        *ros = (TS_RosW*)ts->data;
725   PetscInt       s;
726   PetscErrorCode ierr;
727 
728   PetscFunctionBegin;
729   if (!ros->tableau) PetscFunctionReturn(0);
730    s = ros->tableau->s;
731   ierr = VecDestroyVecs(s,&ros->Y);CHKERRQ(ierr);
732   ierr = VecDestroy(&ros->Ydot);CHKERRQ(ierr);
733   ierr = VecDestroy(&ros->Ystage);CHKERRQ(ierr);
734   ierr = VecDestroy(&ros->Zdot);CHKERRQ(ierr);
735   ierr = VecDestroy(&ros->Zstage);CHKERRQ(ierr);
736   ierr = PetscFree(ros->work);CHKERRQ(ierr);
737   PetscFunctionReturn(0);
738 }
739 
740 #undef __FUNCT__
741 #define __FUNCT__ "TSDestroy_RosW"
742 static PetscErrorCode TSDestroy_RosW(TS ts)
743 {
744   PetscErrorCode  ierr;
745 
746   PetscFunctionBegin;
747   ierr = TSReset_RosW(ts);CHKERRQ(ierr);
748   ierr = PetscFree(ts->data);CHKERRQ(ierr);
749   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","",PETSC_NULL);CHKERRQ(ierr);
750   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","",PETSC_NULL);CHKERRQ(ierr);
751   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","",PETSC_NULL);CHKERRQ(ierr);
752   PetscFunctionReturn(0);
753 }
754 
755 /*
756   This defines the nonlinear equation that is to be solved with SNES
757   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
758 */
759 #undef __FUNCT__
760 #define __FUNCT__ "SNESTSFormFunction_RosW"
761 static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec X,Vec F,TS ts)
762 {
763   TS_RosW        *ros = (TS_RosW*)ts->data;
764   PetscErrorCode ierr;
765 
766   PetscFunctionBegin;
767   ierr = VecWAXPY(ros->Ydot,ros->shift,X,ros->Zdot);CHKERRQ(ierr); /* Ydot = shift*X + Zdot */
768   ierr = VecWAXPY(ros->Ystage,1.0,X,ros->Zstage);CHKERRQ(ierr);    /* Ystage = X + Zstage */
769   ierr = TSComputeIFunction(ts,ros->stage_time,ros->Ystage,ros->Ydot,F,PETSC_FALSE);CHKERRQ(ierr);
770   PetscFunctionReturn(0);
771 }
772 
773 #undef __FUNCT__
774 #define __FUNCT__ "SNESTSFormJacobian_RosW"
775 static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts)
776 {
777   TS_RosW        *ros = (TS_RosW*)ts->data;
778   PetscErrorCode ierr;
779 
780   PetscFunctionBegin;
781   /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */
782   ierr = TSComputeIJacobian(ts,ros->stage_time,ros->Ystage,ros->Ydot,ros->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr);
783   PetscFunctionReturn(0);
784 }
785 
786 #undef __FUNCT__
787 #define __FUNCT__ "TSSetUp_RosW"
788 static PetscErrorCode TSSetUp_RosW(TS ts)
789 {
790   TS_RosW        *ros = (TS_RosW*)ts->data;
791   RosWTableau    tab  = ros->tableau;
792   PetscInt       s    = tab->s;
793   PetscErrorCode ierr;
794 
795   PetscFunctionBegin;
796   if (!ros->tableau) {
797     ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);
798   }
799   ierr = VecDuplicateVecs(ts->vec_sol,s,&ros->Y);CHKERRQ(ierr);
800   ierr = VecDuplicate(ts->vec_sol,&ros->Ydot);CHKERRQ(ierr);
801   ierr = VecDuplicate(ts->vec_sol,&ros->Ystage);CHKERRQ(ierr);
802   ierr = VecDuplicate(ts->vec_sol,&ros->Zdot);CHKERRQ(ierr);
803   ierr = VecDuplicate(ts->vec_sol,&ros->Zstage);CHKERRQ(ierr);
804   ierr = PetscMalloc(s*sizeof(ros->work[0]),&ros->work);CHKERRQ(ierr);
805   PetscFunctionReturn(0);
806 }
807 /*------------------------------------------------------------*/
808 
809 #undef __FUNCT__
810 #define __FUNCT__ "TSSetFromOptions_RosW"
811 static PetscErrorCode TSSetFromOptions_RosW(TS ts)
812 {
813   TS_RosW        *ros = (TS_RosW*)ts->data;
814   PetscErrorCode ierr;
815   char           rostype[256];
816 
817   PetscFunctionBegin;
818   ierr = PetscOptionsHead("RosW ODE solver options");CHKERRQ(ierr);
819   {
820     RosWTableauLink link;
821     PetscInt count,choice;
822     PetscBool flg;
823     const char **namelist;
824     SNES snes;
825 
826     ierr = PetscStrncpy(rostype,TSRosWDefault,sizeof rostype);CHKERRQ(ierr);
827     for (link=RosWTableauList,count=0; link; link=link->next,count++) ;
828     ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr);
829     for (link=RosWTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
830     ierr = PetscOptionsEList("-ts_rosw_type","Family of Rosenbrock-W method","TSRosWSetType",(const char*const*)namelist,count,rostype,&choice,&flg);CHKERRQ(ierr);
831     ierr = TSRosWSetType(ts,flg ? namelist[choice] : rostype);CHKERRQ(ierr);
832     ierr = PetscFree(namelist);CHKERRQ(ierr);
833 
834     ierr = PetscOptionsBool("-ts_rosw_recompute_jacobian","Recompute the Jacobian at each stage","TSRosWSetRecomputeJacobian",ros->recompute_jacobian,&ros->recompute_jacobian,PETSC_NULL);CHKERRQ(ierr);
835 
836     /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
837     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
838     if (!((PetscObject)snes)->type_name) {
839       ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);
840     }
841     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
842   }
843   ierr = PetscOptionsTail();CHKERRQ(ierr);
844   PetscFunctionReturn(0);
845 }
846 
847 #undef __FUNCT__
848 #define __FUNCT__ "PetscFormatRealArray"
849 static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
850 {
851   PetscErrorCode ierr;
852   PetscInt i;
853   size_t left,count;
854   char *p;
855 
856   PetscFunctionBegin;
857   for (i=0,p=buf,left=len; i<n; i++) {
858     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
859     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
860     left -= count;
861     p += count;
862     *p++ = ' ';
863   }
864   p[i ? 0 : -1] = 0;
865   PetscFunctionReturn(0);
866 }
867 
868 #undef __FUNCT__
869 #define __FUNCT__ "TSView_RosW"
870 static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer)
871 {
872   TS_RosW        *ros = (TS_RosW*)ts->data;
873   RosWTableau    tab  = ros->tableau;
874   PetscBool      iascii;
875   PetscErrorCode ierr;
876 
877   PetscFunctionBegin;
878   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
879   if (iascii) {
880     const TSRosWType rostype;
881     PetscInt i;
882     PetscReal abscissa[512];
883     char buf[512];
884     ierr = TSRosWGetType(ts,&rostype);CHKERRQ(ierr);
885     ierr = PetscViewerASCIIPrintf(viewer,"  Rosenbrock-W %s\n",rostype);CHKERRQ(ierr);
886     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ASum);CHKERRQ(ierr);
887     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A       = %s\n",buf);CHKERRQ(ierr);
888     for (i=0; i<tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i];
889     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,abscissa);CHKERRQ(ierr);
890     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A+Gamma = %s\n",buf);CHKERRQ(ierr);
891   }
892   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
893   PetscFunctionReturn(0);
894 }
895 
896 #undef __FUNCT__
897 #define __FUNCT__ "TSRosWSetType"
898 /*@C
899   TSRosWSetType - Set the type of Rosenbrock-W scheme
900 
901   Logically collective
902 
903   Input Parameter:
904 +  ts - timestepping context
905 -  rostype - type of Rosenbrock-W scheme
906 
907   Level: intermediate
908 
909 .seealso: TSRosWGetType()
910 @*/
911 PetscErrorCode TSRosWSetType(TS ts,const TSRosWType rostype)
912 {
913   PetscErrorCode ierr;
914 
915   PetscFunctionBegin;
916   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
917   ierr = PetscTryMethod(ts,"TSRosWSetType_C",(TS,const TSRosWType),(ts,rostype));CHKERRQ(ierr);
918   PetscFunctionReturn(0);
919 }
920 
921 #undef __FUNCT__
922 #define __FUNCT__ "TSRosWGetType"
923 /*@C
924   TSRosWGetType - Get the type of Rosenbrock-W scheme
925 
926   Logically collective
927 
928   Input Parameter:
929 .  ts - timestepping context
930 
931   Output Parameter:
932 .  rostype - type of Rosenbrock-W scheme
933 
934   Level: intermediate
935 
936 .seealso: TSRosWGetType()
937 @*/
938 PetscErrorCode TSRosWGetType(TS ts,const TSRosWType *rostype)
939 {
940   PetscErrorCode ierr;
941 
942   PetscFunctionBegin;
943   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
944   ierr = PetscUseMethod(ts,"TSRosWGetType_C",(TS,const TSRosWType*),(ts,rostype));CHKERRQ(ierr);
945   PetscFunctionReturn(0);
946 }
947 
948 #undef __FUNCT__
949 #define __FUNCT__ "TSRosWSetRecomputeJacobian"
950 /*@C
951   TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step.
952 
953   Logically collective
954 
955   Input Parameter:
956 +  ts - timestepping context
957 -  flg - PETSC_TRUE to recompute the Jacobian at each stage
958 
959   Level: intermediate
960 
961 .seealso: TSRosWGetType()
962 @*/
963 PetscErrorCode TSRosWSetRecomputeJacobian(TS ts,PetscBool flg)
964 {
965   PetscErrorCode ierr;
966 
967   PetscFunctionBegin;
968   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
969   ierr = PetscTryMethod(ts,"TSRosWSetRecomputeJacobian_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
970   PetscFunctionReturn(0);
971 }
972 
973 EXTERN_C_BEGIN
974 #undef __FUNCT__
975 #define __FUNCT__ "TSRosWGetType_RosW"
976 PetscErrorCode  TSRosWGetType_RosW(TS ts,const TSRosWType *rostype)
977 {
978   TS_RosW        *ros = (TS_RosW*)ts->data;
979   PetscErrorCode ierr;
980 
981   PetscFunctionBegin;
982   if (!ros->tableau) {ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);}
983   *rostype = ros->tableau->name;
984   PetscFunctionReturn(0);
985 }
986 #undef __FUNCT__
987 #define __FUNCT__ "TSRosWSetType_RosW"
988 PetscErrorCode  TSRosWSetType_RosW(TS ts,const TSRosWType rostype)
989 {
990   TS_RosW         *ros = (TS_RosW*)ts->data;
991   PetscErrorCode  ierr;
992   PetscBool       match;
993   RosWTableauLink link;
994 
995   PetscFunctionBegin;
996   if (ros->tableau) {
997     ierr = PetscStrcmp(ros->tableau->name,rostype,&match);CHKERRQ(ierr);
998     if (match) PetscFunctionReturn(0);
999   }
1000   for (link = RosWTableauList; link; link=link->next) {
1001     ierr = PetscStrcmp(link->tab.name,rostype,&match);CHKERRQ(ierr);
1002     if (match) {
1003       ierr = TSReset_RosW(ts);CHKERRQ(ierr);
1004       ros->tableau = &link->tab;
1005       PetscFunctionReturn(0);
1006     }
1007   }
1008   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rostype);
1009   PetscFunctionReturn(0);
1010 }
1011 
1012 #undef __FUNCT__
1013 #define __FUNCT__ "TSRosWSetRecomputeJacobian_RosW"
1014 PetscErrorCode  TSRosWSetRecomputeJacobian_RosW(TS ts,PetscBool flg)
1015 {
1016   TS_RosW *ros = (TS_RosW*)ts->data;
1017 
1018   PetscFunctionBegin;
1019   ros->recompute_jacobian = flg;
1020   PetscFunctionReturn(0);
1021 }
1022 EXTERN_C_END
1023 
1024 /* ------------------------------------------------------------ */
1025 /*MC
1026       TSRosW - ODE solver using Rosenbrock-W schemes
1027 
1028   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
1029   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
1030   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
1031 
1032   Notes:
1033   This method currently only works with autonomous ODE and DAE.
1034 
1035   Developer notes:
1036   Rosenbrock-W methods are typically specified for autonomous ODE
1037 
1038 $  xdot = f(x)
1039 
1040   by the stage equations
1041 
1042 $  k_i = h f(x_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j
1043 
1044   and step completion formula
1045 
1046 $  x_1 = x_0 + sum_j b_j k_j
1047 
1048   with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(x)
1049   and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner,
1050   we define new variables for the stage equations
1051 
1052 $  y_i = gamma_ij k_j
1053 
1054   The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define
1055 
1056 $  A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-i}
1057 
1058   to rewrite the method as
1059 
1060 $  [M/(h gamma_ii) - J] y_i = f(x_0 + sum_j a_ij y_j) + M sum_j (c_ij/h) y_j
1061 $  x_1 = x_0 + sum_j bt_j y_j
1062 
1063    where we have introduced the mass matrix M. Continue by defining
1064 
1065 $  ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j
1066 
1067    or, more compactly in tensor notation
1068 
1069 $  Ydot = 1/h (Gamma^{-1} \otimes I) Y .
1070 
1071    Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current
1072    stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the
1073    equation
1074 
1075 $  g(x_0 + sum_j a_ij y_j + y_i, ydot_i) = 0
1076 
1077    with initial guess y_i = 0.
1078 
1079   Level: beginner
1080 
1081 .seealso:  TSCreate(), TS, TSSetType(), TSRosWRegister()
1082 
1083 M*/
1084 EXTERN_C_BEGIN
1085 #undef __FUNCT__
1086 #define __FUNCT__ "TSCreate_RosW"
1087 PetscErrorCode  TSCreate_RosW(TS ts)
1088 {
1089   TS_RosW        *ros;
1090   PetscErrorCode ierr;
1091 
1092   PetscFunctionBegin;
1093 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
1094   ierr = TSRosWInitializePackage(PETSC_NULL);CHKERRQ(ierr);
1095 #endif
1096 
1097   ts->ops->reset          = TSReset_RosW;
1098   ts->ops->destroy        = TSDestroy_RosW;
1099   ts->ops->view           = TSView_RosW;
1100   ts->ops->setup          = TSSetUp_RosW;
1101   ts->ops->step           = TSStep_RosW;
1102   ts->ops->interpolate    = TSInterpolate_RosW;
1103   ts->ops->evaluatestep   = TSEvaluateStep_RosW;
1104   ts->ops->setfromoptions = TSSetFromOptions_RosW;
1105   ts->ops->snesfunction   = SNESTSFormFunction_RosW;
1106   ts->ops->snesjacobian   = SNESTSFormJacobian_RosW;
1107 
1108   ierr = PetscNewLog(ts,TS_RosW,&ros);CHKERRQ(ierr);
1109   ts->data = (void*)ros;
1110 
1111   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","TSRosWGetType_RosW",TSRosWGetType_RosW);CHKERRQ(ierr);
1112   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","TSRosWSetType_RosW",TSRosWSetType_RosW);CHKERRQ(ierr);
1113   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","TSRosWSetRecomputeJacobian_RosW",TSRosWSetRecomputeJacobian_RosW);CHKERRQ(ierr);
1114   PetscFunctionReturn(0);
1115 }
1116 EXTERN_C_END
1117