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