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