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