xref: /petsc/src/ts/impls/rosw/rosw.c (revision 74778d6c3e36e5aa800405c1381bdd58dfd20219)
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 = TSROSWRA34PW2;
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   PetscInt  pinterp;            /* Interpolation order */
27   PetscReal *A;                 /* Propagation table, strictly lower triangular */
28   PetscReal *Gamma;             /* Stage table, lower triangular with nonzero diagonal */
29   PetscBool *GammaZeroDiag;     /* Diagonal entries that are zero in stage table Gamma, vector indicating explicit statages */
30   PetscReal *GammaExplicitCorr; /* Coefficients for correction terms needed for explicit stages in transformed variables*/
31   PetscReal *b;                 /* Step completion table */
32   PetscReal *bembed;            /* Step completion table for embedded method of order one less */
33   PetscReal *ASum;              /* Row sum of A */
34   PetscReal *GammaSum;          /* Row sum of Gamma, only needed for non-autonomous systems */
35   PetscReal *At;                /* Propagation table in transformed variables */
36   PetscReal *bt;                /* Step completion table in transformed variables */
37   PetscReal *bembedt;           /* Step completion table of order one less in transformed variables */
38   PetscReal *GammaInv;          /* Inverse of Gamma, used for transformed variables */
39   PetscReal ccfl;               /* Placeholder for CFL coefficient relative to forward Euler */
40   PetscReal *binterpt;          /* Dense output formula */
41 };
42 typedef struct _RosWTableauLink *RosWTableauLink;
43 struct _RosWTableauLink {
44   struct _RosWTableau tab;
45   RosWTableauLink next;
46 };
47 static RosWTableauLink RosWTableauList;
48 
49 typedef struct {
50   RosWTableau tableau;
51   Vec         *Y;               /* States computed during the step, used to complete the step */
52   Vec         Ydot;             /* Work vector holding Ydot during residual evaluation */
53   Vec         Ystage;           /* Work vector for the state value at each stage */
54   Vec         Zdot;             /* Ydot = Zdot + shift*Y */
55   Vec         Zstage;           /* Y = Zstage + Y */
56   Vec         VecSolPrev;       /* Work vector holding the solution from the previous step (used for interpolation)*/
57   PetscScalar *work;            /* Scalar work space of length number of stages, used to prepare VecMAXPY() */
58   PetscReal   shift;
59   PetscReal   stage_time;
60   PetscReal   stage_explicit;     /* Flag indicates that the current stage is explicit */
61   PetscBool   recompute_jacobian; /* Recompute the Jacobian at each stage, default is to freeze the Jacobian at the start of each step */
62   TSStepStatus status;
63 } TS_RosW;
64 
65 /*MC
66      TSROSWTHETA1 - One stage first order L-stable Rosenbrock-W scheme (aka theta method).
67 
68      Only an approximate Jacobian is needed.
69 
70      Level: intermediate
71 
72 .seealso: TSROSW
73 M*/
74 
75 /*MC
76      TSROSWTHETA2 - One stage second order A-stable Rosenbrock-W scheme (aka theta method).
77 
78      Only an approximate Jacobian is needed.
79 
80      Level: intermediate
81 
82 .seealso: TSROSW
83 M*/
84 
85 /*MC
86      TSROSW2M - Two stage second order L-stable Rosenbrock-W scheme.
87 
88      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2P.
89 
90      Level: intermediate
91 
92 .seealso: TSROSW
93 M*/
94 
95 /*MC
96      TSROSW2P - Two stage second order L-stable Rosenbrock-W scheme.
97 
98      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2M.
99 
100      Level: intermediate
101 
102 .seealso: TSROSW
103 M*/
104 
105 /*MC
106      TSROSWRA3PW - Three stage third order Rosenbrock-W scheme for PDAE of index 1.
107 
108      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.
109 
110      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.
111 
112      References:
113      Rang and Angermann, New Rosenbrock-W methods of order 3 for partial differential algebraic equations of index 1, 2005.
114 
115      Level: intermediate
116 
117 .seealso: TSROSW
118 M*/
119 
120 /*MC
121      TSROSWRA34PW2 - Four stage third order L-stable Rosenbrock-W scheme for PDAE of index 1.
122 
123      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.
124 
125      This is strongly A-stable with R(infty) = 0. The embedded method of order 2 is strongly A-stable with R(infty) = 0.48.
126 
127      References:
128      Rang and Angermann, New Rosenbrock-W methods of order 3 for partial differential algebraic equations of index 1, 2005.
129 
130      Level: intermediate
131 
132 .seealso: TSROSW
133 M*/
134 
135 /*MC
136      TSROSWRODAS3 - Four stage third order L-stable Rosenbrock scheme
137 
138      By default, the Jacobian is only recomputed once per step.
139 
140      Both the third order and embedded second order methods are stiffly accurate and L-stable.
141 
142      References:
143      Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997.
144 
145      Level: intermediate
146 
147 .seealso: TSROSW, TSROSWSANDU3
148 M*/
149 
150 /*MC
151      TSROSWSANDU3 - Three stage third order L-stable Rosenbrock scheme
152 
153      By default, the Jacobian is only recomputed once per step.
154 
155      The third order method is L-stable, but not stiffly accurate.
156      The second order embedded method is strongly A-stable with R(infty) = 0.5.
157      The internal stages are L-stable.
158      This method is called ROS3 in the paper.
159 
160      References:
161      Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997.
162 
163      Level: intermediate
164 
165 .seealso: TSROSW, TSROSWRODAS3
166 M*/
167 
168 /*MC
169      TSROSWASSP3P3S1C - A-stable Rosenbrock-W method with SSP explicit part, third order, three stages
170 
171      By default, the Jacobian is only recomputed once per step.
172 
173      A-stable SPP explicit order 3, 3 stages, CFL 1 (eff = 1/3)
174 
175      References:
176      Emil Constantinescu
177 
178      Level: intermediate
179 
180 .seealso: TSROSW, TSROSWLASSP3P4S2C, TSROSWLLSSP3P4S2C, SSP
181 M*/
182 
183 /*MC
184      TSROSWLASSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, three stages
185 
186      By default, the Jacobian is only recomputed once per step.
187 
188      L-stable (A-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2)
189 
190      References:
191      Emil Constantinescu
192 
193      Level: intermediate
194 
195 .seealso: TSROSW, TSROSWASSP3P3S1C, TSROSWLLSSP3P4S2C, TSSSP
196 M*/
197 
198 /*MC
199      TSROSWLLSSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, three stages
200 
201      By default, the Jacobian is only recomputed once per step.
202 
203      L-stable (L-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2)
204 
205      References:
206      Emil Constantinescu
207 
208      Level: intermediate
209 
210 .seealso: TSROSW, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSSSP
211 M*/
212 
213 #undef __FUNCT__
214 #define __FUNCT__ "TSRosWRegisterAll"
215 /*@C
216   TSRosWRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSRosW
217 
218   Not Collective, but should be called by all processes which will need the schemes to be registered
219 
220   Level: advanced
221 
222 .keywords: TS, TSRosW, register, all
223 
224 .seealso:  TSRosWRegisterDestroy()
225 @*/
226 PetscErrorCode TSRosWRegisterAll(void)
227 {
228   PetscErrorCode ierr;
229 
230   PetscFunctionBegin;
231   if (TSRosWRegisterAllCalled) PetscFunctionReturn(0);
232   TSRosWRegisterAllCalled = PETSC_TRUE;
233 
234   {
235     const PetscReal
236       A = 0,
237       Gamma = 1,
238       b = 1,
239       binterpt=1;
240 
241     ierr = TSRosWRegister(TSROSWTHETA1,1,1,&A,&Gamma,&b,PETSC_NULL,1,&binterpt);CHKERRQ(ierr);
242   }
243 
244   {
245     const PetscReal
246       A= 0,
247       Gamma = 0.5,
248       b = 1,
249       binterpt=1;
250     ierr = TSRosWRegister(TSROSWTHETA2,2,1,&A,&Gamma,&b,PETSC_NULL,1,&binterpt);CHKERRQ(ierr);
251   }
252 
253   {
254     const PetscReal g = 1. + 1./PetscSqrtReal(2.0);
255     const PetscReal
256       A[2][2] = {{0,0}, {1.,0}},
257       Gamma[2][2] = {{g,0}, {-2.*g,g}},
258       b[2] = {0.5,0.5},
259       b1[2] = {1.0,0.0};
260       PetscReal  binterpt[2][2];
261       binterpt[0][0]=g-1.0;
262       binterpt[1][0]=2.0-g;
263       binterpt[0][1]=g-1.5;
264       binterpt[1][1]=1.5-g;
265       ierr = TSRosWRegister(TSROSW2P,2,2,&A[0][0],&Gamma[0][0],b,b1,2,&binterpt[0][0]);CHKERRQ(ierr);
266   }
267   {
268     const PetscReal g = 1. - 1./PetscSqrtReal(2.0);
269     const PetscReal
270       A[2][2] = {{0,0}, {1.,0}},
271       Gamma[2][2] = {{g,0}, {-2.*g,g}},
272       b[2] = {0.5,0.5},
273       b1[2] = {1.0,0.0};
274       PetscReal  binterpt[2][2];
275       binterpt[0][0]=g-1.0;
276       binterpt[1][0]=2.0-g;
277       binterpt[0][1]=g-1.5;
278       binterpt[1][1]=1.5-g;
279     ierr = TSRosWRegister(TSROSW2M,2,2,&A[0][0],&Gamma[0][0],b,b1,2,&binterpt[0][0]);CHKERRQ(ierr);
280   }
281   {
282     const PetscReal g = 7.8867513459481287e-01;
283     PetscReal  binterpt[3][2];
284     const PetscReal
285       A[3][3] = {{0,0,0},
286                  {1.5773502691896257e+00,0,0},
287                  {0.5,0,0}},
288       Gamma[3][3] = {{g,0,0},
289                      {-1.5773502691896257e+00,g,0},
290                      {-6.7075317547305480e-01,-1.7075317547305482e-01,g}},
291       b[3] = {1.0566243270259355e-01,4.9038105676657971e-02,8.4529946162074843e-01},
292       b2[3] = {-1.7863279495408180e-01,1./3.,8.4529946162074843e-01};
293 
294       binterpt[0][0]=-0.8094010767585034;
295       binterpt[1][0]=-0.5;
296       binterpt[2][0]=2.3094010767585034;
297       binterpt[0][1]=0.9641016151377548;
298       binterpt[1][1]=0.5;
299       binterpt[2][1]=-1.4641016151377548;
300       ierr = TSRosWRegister(TSROSWRA3PW,3,3,&A[0][0],&Gamma[0][0],b,b2,2,&binterpt[0][0]);CHKERRQ(ierr);
301   }
302   {
303     PetscReal  binterpt[4][3];
304     const PetscReal g = 4.3586652150845900e-01;
305     const PetscReal
306       A[4][4] = {{0,0,0,0},
307                  {8.7173304301691801e-01,0,0,0},
308                  {8.4457060015369423e-01,-1.1299064236484185e-01,0,0},
309                  {0,0,1.,0}},
310       Gamma[4][4] = {{g,0,0,0},
311                      {-8.7173304301691801e-01,g,0,0},
312                      {-9.0338057013044082e-01,5.4180672388095326e-02,g,0},
313                      {2.4212380706095346e-01,-1.2232505839045147e+00,5.4526025533510214e-01,g}},
314         b[4] = {2.4212380706095346e-01,-1.2232505839045147e+00,1.5452602553351020e+00,4.3586652150845900e-01},
315           b2[4] = {3.7810903145819369e-01,-9.6042292212423178e-02,5.0000000000000000e-01,2.1793326075422950e-01};
316 
317           binterpt[0][0]=1.0564298455794094;
318           binterpt[1][0]=2.296429974281067;
319           binterpt[2][0]=-1.307599564525376;
320           binterpt[3][0]=-1.045260255335102;
321           binterpt[0][1]=-1.3864882699759573;
322           binterpt[1][1]=-8.262611700275677;
323           binterpt[2][1]=7.250979895056055;
324           binterpt[3][1]=2.398120075195581;
325           binterpt[0][2]=0.5721822314575016;
326           binterpt[1][2]=4.742931142090097;
327           binterpt[2][2]=-4.398120075195578;
328           binterpt[3][2]=-0.9169932983520199;
329 
330           ierr = TSRosWRegister(TSROSWRA34PW2,3,4,&A[0][0],&Gamma[0][0],b,b2,3,&binterpt[0][0]);CHKERRQ(ierr);
331   }
332   {
333     const PetscReal g = 0.5;
334     const PetscReal
335       A[4][4] = {{0,0,0,0},
336                  {0,0,0,0},
337                  {1.,0,0,0},
338                  {0.75,-0.25,0.5,0}},
339       Gamma[4][4] = {{g,0,0,0},
340                      {1.,g,0,0},
341                      {-0.25,-0.25,g,0},
342                      {1./12,1./12,-2./3,g}},
343       b[4] = {5./6,-1./6,-1./6,0.5},
344       b2[4] = {0.75,-0.25,0.5,0};
345     ierr = TSRosWRegister(TSROSWRODAS3,3,4,&A[0][0],&Gamma[0][0],b,b2,0,PETSC_NULL);CHKERRQ(ierr);
346   }
347   {
348     const PetscReal g = 0.43586652150845899941601945119356;
349     const PetscReal
350       A[3][3] = {{0,0,0},
351                  {g,0,0},
352                  {g,0,0}},
353       Gamma[3][3] = {{g,0,0},
354                      {-0.19294655696029095575009695436041,g,0},
355                      {0,1.74927148125794685173529749738960,g}},
356       b[3] = {-0.75457412385404315829818998646589,1.94100407061964420292840123379419,-0.18642994676560104463021124732829},
357       b2[3] = {-1.53358745784149585370766523913002,2.81745131148625772213931745457622,-0.28386385364476186843165221544619};
358 
359       PetscReal  binterpt[3][2];
360       binterpt[0][0]=3.793692883777660870425141387941;
361       binterpt[1][0]=-2.918692883777660870425141387941;
362       binterpt[2][0]=0.125;
363       binterpt[0][1]=-0.725741064379812106687651020584;
364       binterpt[1][1]=0.559074397713145440020984353917;
365       binterpt[2][1]=0.16666666666666666666666666666667;
366 
367       ierr = TSRosWRegister(TSROSWSANDU3,3,3,&A[0][0],&Gamma[0][0],b,b2,2,&binterpt[0][0]);CHKERRQ(ierr);
368   }
369   {
370     const PetscReal s3 = PetscSqrtReal(3.),g = (3.0+s3)/6.0;
371     const PetscReal
372       A[3][3] = {{0,0,0},
373                  {1,0,0},
374                  {0.25,0.25,0}},
375       Gamma[3][3] = {{0,0,0},
376                      {(-3.0-s3)/6.0,g,0},
377                      {(-3.0-s3)/24.0,(-3.0-s3)/8.0,g}},
378         b[3] = {1./6.,1./6.,2./3.},
379           b2[3] = {1./4.,1./4.,1./2.};
380     ierr = TSRosWRegister(TSROSWASSP3P3S1C,3,3,&A[0][0],&Gamma[0][0],b,b2,0,PETSC_NULL);CHKERRQ(ierr);
381   }
382 
383   {
384     const PetscReal
385       A[4][4] = {{0,0,0,0},
386                  {1./2.,0,0,0},
387                  {1./2.,1./2.,0,0},
388                  {1./6.,1./6.,1./6.,0}},
389       Gamma[4][4] = {{1./2.,0,0,0},
390                      {0.0,1./4.,0,0},
391                      {-2.,-2./3.,2./3.,0},
392                      {1./2.,5./36.,-2./9,0}},
393         b[4] = {1./6.,1./6.,1./6.,1./2.},
394         b2[4] = {1./8.,3./4.,1./8.,0};
395      ierr = TSRosWRegister(TSROSWLASSP3P4S2C,3,4,&A[0][0],&Gamma[0][0],b,b2,0,PETSC_NULL);CHKERRQ(ierr);
396   }
397 
398   {
399     const PetscReal
400       A[4][4] = {{0,0,0,0},
401                  {1./2.,0,0,0},
402                  {1./2.,1./2.,0,0},
403                  {1./6.,1./6.,1./6.,0}},
404       Gamma[4][4] = {{1./2.,0,0,0},
405                      {0.0,3./4.,0,0},
406                      {-2./3.,-23./9.,2./9.,0},
407                      {1./18.,65./108.,-2./27,0}},
408         b[4] = {1./6.,1./6.,1./6.,1./2.},
409         b2[4] = {3./16.,10./16.,3./16.,0};
410      ierr = TSRosWRegister(TSROSWLLSSP3P4S2C,3,4,&A[0][0],&Gamma[0][0],b,b2,0,PETSC_NULL);CHKERRQ(ierr);
411   }
412 
413  {
414    PetscReal A[4][4],Gamma[4][4],b[4],b2[4];
415    PetscReal  binterpt[4][3];
416 
417    Gamma[0][0]=0.4358665215084589994160194475295062513822671686978816;
418    Gamma[0][1]=0; Gamma[0][2]=0; Gamma[0][3]=0;
419    Gamma[1][0]=-1.997527830934941248426324674704153457289527280554476;
420    Gamma[1][1]=0.4358665215084589994160194475295062513822671686978816;
421    Gamma[1][2]=0; Gamma[1][3]=0;
422    Gamma[2][0]=-1.007948511795029620852002345345404191008352770119903;
423    Gamma[2][1]=-0.004648958462629345562774289390054679806993396798458131;
424    Gamma[2][2]=0.4358665215084589994160194475295062513822671686978816;
425    Gamma[2][3]=0;
426    Gamma[3][0]=-0.6685429734233467180451604600279552604364311322650783;
427    Gamma[3][1]=0.6056625986449338476089525334450053439525178740492984;
428    Gamma[3][2]=-0.9717899277217721234705114616271378792182450260943198;
429    Gamma[3][3]=0;
430 
431    A[0][0]=0; A[0][1]=0; A[0][2]=0; A[0][3]=0;
432    A[1][0]=0.8717330430169179988320388950590125027645343373957631;
433    A[1][1]=0; A[1][2]=0; A[1][3]=0;
434    A[2][0]=0.5275890119763004115618079766722914408876108660811028;
435    A[2][1]=0.07241098802369958843819203208518599088698057726988732;
436    A[2][2]=0; A[2][3]=0;
437    A[3][0]=0.3990960076760701320627260685975778145384666450351314;
438    A[3][1]=-0.4375576546135194437228463747348862825846903771419953;
439    A[3][2]=1.038461646937449311660120300601880176655352737312713;
440    A[3][3]=0;
441 
442    b[0]=0.1876410243467238251612921333138006734899663569186926;
443    b[1]=-0.5952974735769549480478230473706443582188442040780541;
444    b[2]=0.9717899277217721234705114616271378792182450260943198;
445    b[3]=0.4358665215084589994160194475295062513822671686978816;
446 
447    b2[0]=0.2147402862233891404862383521089097657790734483804460;
448    b2[1]=-0.4851622638849390928209050538171743017757490232519684;
449    b2[2]=0.8687250025203875511662123688667549217531982787600080;
450    b2[3]=0.4016969751411624011684543450940068201770721128357014;
451 
452    binterpt[0][0]=2.2565812720167954547104627844105;
453    binterpt[1][0]=1.349166413351089573796243820819;
454    binterpt[2][0]=-2.4695174540533503758652847586647;
455    binterpt[3][0]=-0.13623023131453465264142184656474;
456    binterpt[0][1]=-3.0826699111559187902922463354557;
457    binterpt[1][1]=-2.4689115685996042534544925650515;
458    binterpt[2][1]=5.7428279814696677152129332773553;
459    binterpt[3][1]=-0.19124650171414467146619437684812;
460    binterpt[0][2]=1.0137296634858471607430756831148;
461    binterpt[1][2]=0.52444768167155973161042570784064;
462    binterpt[2][2]=-2.3015205996945452158771370439586;
463    binterpt[3][2]=0.76334325453713832352363565300308;
464 
465    ierr = TSRosWRegister(TSROSWARK3,3,4,&A[0][0],&Gamma[0][0],b,b2,3,&binterpt[0][0]);CHKERRQ(ierr);
466   }
467 
468   PetscFunctionReturn(0);
469 }
470 
471 #undef __FUNCT__
472 #define __FUNCT__ "TSRosWRegisterDestroy"
473 /*@C
474    TSRosWRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister().
475 
476    Not Collective
477 
478    Level: advanced
479 
480 .keywords: TSRosW, register, destroy
481 .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegisterDynamic()
482 @*/
483 PetscErrorCode TSRosWRegisterDestroy(void)
484 {
485   PetscErrorCode ierr;
486   RosWTableauLink link;
487 
488   PetscFunctionBegin;
489   while ((link = RosWTableauList)) {
490     RosWTableau t = &link->tab;
491     RosWTableauList = link->next;
492     ierr = PetscFree5(t->A,t->Gamma,t->b,t->ASum,t->GammaSum);CHKERRQ(ierr);
493     ierr = PetscFree5(t->At,t->bt,t->GammaInv,t->GammaZeroDiag,t->GammaExplicitCorr);CHKERRQ(ierr);
494     ierr = PetscFree2(t->bembed,t->bembedt);CHKERRQ(ierr);
495     ierr = PetscFree(t->binterpt);CHKERRQ(ierr);
496     ierr = PetscFree(t->name);CHKERRQ(ierr);
497     ierr = PetscFree(link);CHKERRQ(ierr);
498   }
499   TSRosWRegisterAllCalled = PETSC_FALSE;
500   PetscFunctionReturn(0);
501 }
502 
503 #undef __FUNCT__
504 #define __FUNCT__ "TSRosWInitializePackage"
505 /*@C
506   TSRosWInitializePackage - This function initializes everything in the TSRosW package. It is called
507   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RosW()
508   when using static libraries.
509 
510   Input Parameter:
511   path - The dynamic library path, or PETSC_NULL
512 
513   Level: developer
514 
515 .keywords: TS, TSRosW, initialize, package
516 .seealso: PetscInitialize()
517 @*/
518 PetscErrorCode TSRosWInitializePackage(const char path[])
519 {
520   PetscErrorCode ierr;
521 
522   PetscFunctionBegin;
523   if (TSRosWPackageInitialized) PetscFunctionReturn(0);
524   TSRosWPackageInitialized = PETSC_TRUE;
525   ierr = TSRosWRegisterAll();CHKERRQ(ierr);
526   ierr = PetscRegisterFinalize(TSRosWFinalizePackage);CHKERRQ(ierr);
527   PetscFunctionReturn(0);
528 }
529 
530 #undef __FUNCT__
531 #define __FUNCT__ "TSRosWFinalizePackage"
532 /*@C
533   TSRosWFinalizePackage - This function destroys everything in the TSRosW package. It is
534   called from PetscFinalize().
535 
536   Level: developer
537 
538 .keywords: Petsc, destroy, package
539 .seealso: PetscFinalize()
540 @*/
541 PetscErrorCode TSRosWFinalizePackage(void)
542 {
543   PetscErrorCode ierr;
544 
545   PetscFunctionBegin;
546   TSRosWPackageInitialized = PETSC_FALSE;
547   ierr = TSRosWRegisterDestroy();CHKERRQ(ierr);
548   PetscFunctionReturn(0);
549 }
550 
551 #undef __FUNCT__
552 #define __FUNCT__ "TSRosWRegister"
553 /*@C
554    TSRosWRegister - register a Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
555 
556    Not Collective, but the same schemes should be registered on all processes on which they will be used
557 
558    Input Parameters:
559 +  name - identifier for method
560 .  order - approximation order of method
561 .  s - number of stages, this is the dimension of the matrices below
562 .  A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular
563 .  Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal
564 .  b - Step completion table (dimension s)
565 -  bembed - Step completion table for a scheme of order one less (dimension s, PETSC_NULL if no embedded scheme is available)
566 .  pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt
567 .  binterpt - Coefficients of the interpolation formula (dimension s*pinterp)
568 
569    Notes:
570    Several Rosenbrock W methods are provided, this function is only needed to create new methods.
571 
572    Level: advanced
573 
574 .keywords: TS, register
575 
576 .seealso: TSRosW
577 @*/
578 PetscErrorCode TSRosWRegister(const TSRosWType name,PetscInt order,PetscInt s,
579                               const PetscReal A[],const PetscReal Gamma[],const PetscReal b[],const PetscReal bembed[],
580                                  PetscInt pinterp,const PetscReal binterpt[])
581 {
582   PetscErrorCode ierr;
583   RosWTableauLink link;
584   RosWTableau t;
585   PetscInt i,j,k;
586   PetscScalar *GammaInv;
587 
588   PetscFunctionBegin;
589   PetscValidCharPointer(name,1);
590   PetscValidPointer(A,4);
591   PetscValidPointer(Gamma,5);
592   PetscValidPointer(b,6);
593   if (bembed) PetscValidPointer(bembed,7);
594 
595   ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
596   ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
597   t = &link->tab;
598   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
599   t->order = order;
600   t->s = s;
601   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);
602   ierr = PetscMalloc5(s*s,PetscReal,&t->At,s,PetscReal,&t->bt,s*s,PetscReal,&t->GammaInv,s,PetscBool,&t->GammaZeroDiag,s*s,PetscReal,&t->GammaExplicitCorr);CHKERRQ(ierr);
603   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
604   ierr = PetscMemcpy(t->Gamma,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr);
605   ierr = PetscMemcpy(t->GammaExplicitCorr,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr);
606   ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
607   if (bembed) {
608     ierr = PetscMalloc2(s,PetscReal,&t->bembed,s,PetscReal,&t->bembedt);CHKERRQ(ierr);
609     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
610   }
611   for (i=0; i<s; i++) {
612     t->ASum[i] = 0;
613     t->GammaSum[i] = 0;
614     for (j=0; j<s; j++) {
615       t->ASum[i] += A[i*s+j];
616       t->GammaSum[i] += Gamma[i*s+j];
617     }
618   }
619   ierr = PetscMalloc(s*s*sizeof(PetscScalar),&GammaInv);CHKERRQ(ierr); /* Need to use Scalar for inverse, then convert back to Real */
620   for (i=0; i<s*s; i++) GammaInv[i] = Gamma[i];
621   for (i=0; i<s; i++) {
622     if (Gamma[i*s+i] == 0.0) {
623       GammaInv[i*s+i] = 1.0;
624       t->GammaZeroDiag[i] = PETSC_TRUE;
625     } else {
626       t->GammaZeroDiag[i] = PETSC_FALSE;
627     }
628   }
629 
630   switch (s) {
631   case 1: GammaInv[0] = 1./GammaInv[0]; break;
632   case 2: ierr = Kernel_A_gets_inverse_A_2(GammaInv,0);CHKERRQ(ierr); break;
633   case 3: ierr = Kernel_A_gets_inverse_A_3(GammaInv,0);CHKERRQ(ierr); break;
634   case 4: ierr = Kernel_A_gets_inverse_A_4(GammaInv,0);CHKERRQ(ierr); break;
635   case 5: {
636     PetscInt ipvt5[5];
637     MatScalar work5[5*5];
638     ierr = Kernel_A_gets_inverse_A_5(GammaInv,ipvt5,work5,0);CHKERRQ(ierr); break;
639   }
640   case 6: ierr = Kernel_A_gets_inverse_A_6(GammaInv,0);CHKERRQ(ierr); break;
641   case 7: ierr = Kernel_A_gets_inverse_A_7(GammaInv,0);CHKERRQ(ierr); break;
642   default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for %D stages",s);
643   }
644   for (i=0; i<s*s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]);
645   ierr = PetscFree(GammaInv);CHKERRQ(ierr);
646 
647   for (i=0; i<s; i++) {
648     for (k=0; k<i+1; k++) {
649       t->GammaExplicitCorr[i*s+k]=(t->GammaExplicitCorr[i*s+k])*(t->GammaInv[k*s+k]);
650       for (j=k+1; j<i+1; j++) {
651         t->GammaExplicitCorr[i*s+k]+=(t->GammaExplicitCorr[i*s+j])*(t->GammaInv[j*s+k]);
652       }
653     }
654   }
655 
656   for (i=0; i<s; i++) {
657     for (j=0; j<s; j++) {
658       t->At[i*s+j] = 0;
659       for (k=0; k<s; k++) {
660         t->At[i*s+j] += t->A[i*s+k] * t->GammaInv[k*s+j];
661       }
662     }
663     t->bt[i] = 0;
664     for (j=0; j<s; j++) {
665       t->bt[i] += t->b[j] * t->GammaInv[j*s+i];
666     }
667     if (bembed) {
668       t->bembedt[i] = 0;
669       for (j=0; j<s; j++) {
670         t->bembedt[i] += t->bembed[j] * t->GammaInv[j*s+i];
671       }
672     }
673   }
674   t->ccfl = 1.0;                /* Fix this */
675 
676   t->pinterp = pinterp;
677   ierr = PetscMalloc(s*pinterp*sizeof(binterpt[0]),&t->binterpt);CHKERRQ(ierr);
678   ierr = PetscMemcpy(t->binterpt,binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr);
679   link->next = RosWTableauList;
680   RosWTableauList = link;
681   PetscFunctionReturn(0);
682 }
683 
684 #undef __FUNCT__
685 #define __FUNCT__ "TSEvaluateStep_RosW"
686 /*
687  The step completion formula is
688 
689  x1 = x0 + b^T Y
690 
691  where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been
692  updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write
693 
694  x1e = x0 + be^T Y
695      = x1 - b^T Y + be^T Y
696      = x1 + (be - b)^T Y
697 
698  so we can evaluate the method of different order even after the step has been optimistically completed.
699 */
700 static PetscErrorCode TSEvaluateStep_RosW(TS ts,PetscInt order,Vec X,PetscBool *done)
701 {
702   TS_RosW        *ros = (TS_RosW*)ts->data;
703   RosWTableau    tab  = ros->tableau;
704   PetscScalar    *w = ros->work;
705   PetscInt       i;
706   PetscErrorCode ierr;
707 
708   PetscFunctionBegin;
709   if (order == tab->order) {
710     if (ros->status == TS_STEP_INCOMPLETE) { /* Use standard completion formula */
711       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
712       for (i=0; i<tab->s; i++) w[i] = tab->bt[i];
713       ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr);
714     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
715     if (done) *done = PETSC_TRUE;
716     PetscFunctionReturn(0);
717   } else if (order == tab->order-1) {
718     if (!tab->bembedt) goto unavailable;
719     if (ros->status == TS_STEP_INCOMPLETE) { /* Use embedded completion formula */
720       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
721       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i];
722       ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr);
723     } else {                    /* Use rollback-and-recomplete formula (bembedt - bt) */
724       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i];
725       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
726       ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr);
727     }
728     if (done) *done = PETSC_TRUE;
729     PetscFunctionReturn(0);
730   }
731   unavailable:
732   if (done) *done = PETSC_FALSE;
733   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);
734   PetscFunctionReturn(0);
735 }
736 
737 #undef __FUNCT__
738 #define __FUNCT__ "TSStep_RosW"
739 static PetscErrorCode TSStep_RosW(TS ts)
740 {
741   TS_RosW         *ros = (TS_RosW*)ts->data;
742   RosWTableau     tab  = ros->tableau;
743   const PetscInt  s    = tab->s;
744   const PetscReal *At  = tab->At,*Gamma = tab->Gamma,*ASum = tab->ASum,*GammaInv = tab->GammaInv;
745   const PetscReal *GammaExplicitCorr = tab->GammaExplicitCorr;
746   const PetscBool *GammaZeroDiag = tab->GammaZeroDiag;
747   PetscScalar     *w   = ros->work;
748   Vec             *Y   = ros->Y,Ydot = ros->Ydot,Zdot = ros->Zdot,Zstage = ros->Zstage;
749   SNES            snes;
750   TSAdapt         adapt;
751   PetscInt        i,j,its,lits,reject,next_scheme;
752   PetscReal       next_time_step;
753   PetscBool       accept;
754   PetscErrorCode  ierr;
755   MatStructure    str;
756 
757   PetscFunctionBegin;
758   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
759   next_time_step = ts->time_step;
760   accept = PETSC_TRUE;
761   ros->status = TS_STEP_INCOMPLETE;
762 
763   for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) {
764     const PetscReal h = ts->time_step;
765     ierr = VecCopy(ts->vec_sol,ros->VecSolPrev);CHKERRQ(ierr);/*move this at the end*/
766     for (i=0; i<s; i++) {
767       ros->stage_time = ts->ptime + h*ASum[i];
768       if (GammaZeroDiag[i]) {
769         ros->stage_explicit = PETSC_TRUE;
770         ros->shift = 1./h;
771       } else {
772         ros->stage_explicit = PETSC_FALSE;
773         ros->shift = 1./(h*Gamma[i*s+i]);
774       }
775 
776       ierr = VecCopy(ts->vec_sol,Zstage);CHKERRQ(ierr);
777       for (j=0; j<i; j++) w[j] = At[i*s+j];
778       ierr = VecMAXPY(Zstage,i,w,Y);CHKERRQ(ierr);
779 
780       for (j=0; j<i; j++) w[j] = 1./h * GammaInv[i*s+j];
781       ierr = VecZeroEntries(Zdot);CHKERRQ(ierr);
782       ierr = VecMAXPY(Zdot,i,w,Y);CHKERRQ(ierr);
783 
784       /* Initial guess taken from last stage */
785       ierr = VecZeroEntries(Y[i]);CHKERRQ(ierr);
786 
787       if (!ros->stage_explicit) {
788         if (!ros->recompute_jacobian && !i) {
789           ierr = SNESSetLagJacobian(snes,-2);CHKERRQ(ierr); /* Recompute the Jacobian on this solve, but not again */
790         }
791         ierr = SNESSolve(snes,PETSC_NULL,Y[i]);CHKERRQ(ierr);
792         ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
793         ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
794         ts->nonlinear_its += its; ts->linear_its += lits;
795         ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
796         ierr = TSAdaptCheckStage(adapt,ts,&accept);CHKERRQ(ierr);
797         if (!accept) goto reject_step;
798       } else {
799         Mat J,Jp;
800         ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); /* Evaluate Y[i]=G(t,Ydot=0,Zstage) */
801         ierr = TSComputeIFunction(ts,ros->stage_time,Zstage,Ydot,Y[i],PETSC_FALSE);CHKERRQ(ierr);
802         ierr = VecScale(Y[i],-1.0);
803         ierr = VecAXPY(Y[i],-1.0,Zdot);CHKERRQ(ierr); /*Y[i]=F(Zstage)-Zdot[=GammaInv*Y]*/
804 
805         ierr = VecZeroEntries(Zstage);CHKERRQ(ierr); /* Zstage = GammaExplicitCorr[i,j] * Y[j] */
806         for (j=0; j<i; j++) w[j] = GammaExplicitCorr[i*s+j];
807         ierr = VecMAXPY(Zstage,i,w,Y);CHKERRQ(ierr);
808         /*Y[i] += Y[i] + Jac*Zstage[=Jac*GammaExplicitCorr[i,j] * Y[j]] */
809         str = SAME_NONZERO_PATTERN;
810         ierr = TSGetIJacobian(ts,&J,&Jp,PETSC_NULL,PETSC_NULL);
811         ierr = TSComputeIJacobian(ts,ros->stage_time,ts->vec_sol,Ydot,0,&J,&Jp,&str,PETSC_FALSE);CHKERRQ(ierr);
812         ierr = MatMult(J,Zstage,Zdot);
813 
814         ierr = VecAXPY(Y[i],-1.0,Zdot);CHKERRQ(ierr);
815         ierr = VecScale(Y[i],h);
816         ts->linear_its += 1;
817       }
818     }
819     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
820     ros->status = TS_STEP_PENDING;
821 
822     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
823     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
824     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
825     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
826     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
827     if (accept) {
828       /* ignore next_scheme for now */
829       ts->ptime += ts->time_step;
830       ts->time_step = next_time_step;
831       ts->steps++;
832       ros->status = TS_STEP_COMPLETE;
833       break;
834     } else {                    /* Roll back the current step */
835       for (i=0; i<s; i++) w[i] = -tab->bt[i];
836       ierr = VecMAXPY(ts->vec_sol,s,w,Y);CHKERRQ(ierr);
837       ts->time_step = next_time_step;
838       ros->status = TS_STEP_INCOMPLETE;
839     }
840     reject_step: continue;
841   }
842   if (ros->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED;
843   PetscFunctionReturn(0);
844 }
845 
846 #undef __FUNCT__
847 #define __FUNCT__ "TSInterpolate_RosW"
848 static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec X)
849 {
850   TS_RosW        *ros = (TS_RosW*)ts->data;
851   PetscInt s = ros->tableau->s,pinterp = ros->tableau->pinterp,i,j;
852   PetscReal h;
853   PetscReal tt,t;
854   PetscScalar *bt;
855   const PetscReal *Bt = ros->tableau->binterpt;
856   PetscErrorCode ierr;
857   const PetscReal *GammaInv = ros->tableau->GammaInv;
858   PetscScalar     *w   = ros->work;
859   Vec             *Y   = ros->Y;
860 
861   PetscFunctionBegin;
862   if (!Bt) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ros->tableau->name);
863 
864   switch (ros->status) {
865   case TS_STEP_INCOMPLETE:
866   case TS_STEP_PENDING:
867     h = ts->time_step;
868     t = (itime - ts->ptime)/h;
869     break;
870   case TS_STEP_COMPLETE:
871     h = ts->time_step_prev;
872     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
873     break;
874   default: SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Invalid TSStepStatus");
875   }
876   ierr = PetscMalloc(s*sizeof(bt[0]),&bt);CHKERRQ(ierr);
877   for (i=0; i<s; i++) bt[i] = 0;
878   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
879     for (i=0; i<s; i++) {
880       bt[i] += Bt[i*pinterp+j] * tt;
881     }
882   }
883 
884   /* y(t+tt*h) = y(t) + Sum bt(tt) * GammaInv * Ydot */
885   /*X<-0*/
886   ierr = VecZeroEntries(X);CHKERRQ(ierr);
887 
888   /*X<- Sum bt_i * GammaInv(i,1:i) * Y(1:i) */
889   for (j=0; j<s; j++)  w[j]=0;
890   for (j=0; j<s; j++) {
891     for (i=j; i<s; i++) {
892       w[j] +=  bt[i]*GammaInv[i*s+j];
893     }
894   }
895   ierr = VecMAXPY(X,i,w,Y);CHKERRQ(ierr);
896 
897   /*X<-y(t) + X*/
898   ierr = VecAXPY(X,1.0,ros->VecSolPrev);CHKERRQ(ierr);
899 
900   ierr = PetscFree(bt);CHKERRQ(ierr);
901 
902   PetscFunctionReturn(0);
903 }
904 
905 /*------------------------------------------------------------*/
906 #undef __FUNCT__
907 #define __FUNCT__ "TSReset_RosW"
908 static PetscErrorCode TSReset_RosW(TS ts)
909 {
910   TS_RosW        *ros = (TS_RosW*)ts->data;
911   PetscInt       s;
912   PetscErrorCode ierr;
913 
914   PetscFunctionBegin;
915   if (!ros->tableau) PetscFunctionReturn(0);
916    s = ros->tableau->s;
917   ierr = VecDestroyVecs(s,&ros->Y);CHKERRQ(ierr);
918   ierr = VecDestroy(&ros->Ydot);CHKERRQ(ierr);
919   ierr = VecDestroy(&ros->Ystage);CHKERRQ(ierr);
920   ierr = VecDestroy(&ros->Zdot);CHKERRQ(ierr);
921   ierr = VecDestroy(&ros->Zstage);CHKERRQ(ierr);
922   ierr = VecDestroy(&ros->VecSolPrev);CHKERRQ(ierr);
923   ierr = PetscFree(ros->work);CHKERRQ(ierr);
924   PetscFunctionReturn(0);
925 }
926 
927 #undef __FUNCT__
928 #define __FUNCT__ "TSDestroy_RosW"
929 static PetscErrorCode TSDestroy_RosW(TS ts)
930 {
931   PetscErrorCode  ierr;
932 
933   PetscFunctionBegin;
934   ierr = TSReset_RosW(ts);CHKERRQ(ierr);
935   ierr = PetscFree(ts->data);CHKERRQ(ierr);
936   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","",PETSC_NULL);CHKERRQ(ierr);
937   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","",PETSC_NULL);CHKERRQ(ierr);
938   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","",PETSC_NULL);CHKERRQ(ierr);
939   PetscFunctionReturn(0);
940 }
941 
942 /*
943   This defines the nonlinear equation that is to be solved with SNES
944   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
945 */
946 #undef __FUNCT__
947 #define __FUNCT__ "SNESTSFormFunction_RosW"
948 static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec X,Vec F,TS ts)
949 {
950   TS_RosW        *ros = (TS_RosW*)ts->data;
951   PetscErrorCode ierr;
952 
953   PetscFunctionBegin;
954   ierr = VecWAXPY(ros->Ydot,ros->shift,X,ros->Zdot);CHKERRQ(ierr); /* Ydot = shift*X + Zdot */
955   ierr = VecWAXPY(ros->Ystage,1.0,X,ros->Zstage);CHKERRQ(ierr);    /* Ystage = X + Zstage */
956   ierr = TSComputeIFunction(ts,ros->stage_time,ros->Ystage,ros->Ydot,F,PETSC_FALSE);CHKERRQ(ierr);
957   PetscFunctionReturn(0);
958 }
959 
960 #undef __FUNCT__
961 #define __FUNCT__ "SNESTSFormJacobian_RosW"
962 static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts)
963 {
964   TS_RosW        *ros = (TS_RosW*)ts->data;
965   PetscErrorCode ierr;
966 
967   PetscFunctionBegin;
968   /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */
969   ierr = TSComputeIJacobian(ts,ros->stage_time,ros->Ystage,ros->Ydot,ros->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr);
970   PetscFunctionReturn(0);
971 }
972 
973 #undef __FUNCT__
974 #define __FUNCT__ "TSSetUp_RosW"
975 static PetscErrorCode TSSetUp_RosW(TS ts)
976 {
977   TS_RosW        *ros = (TS_RosW*)ts->data;
978   RosWTableau    tab  = ros->tableau;
979   PetscInt       s    = tab->s;
980   PetscErrorCode ierr;
981 
982   PetscFunctionBegin;
983   if (!ros->tableau) {
984     ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);
985   }
986   ierr = VecDuplicateVecs(ts->vec_sol,s,&ros->Y);CHKERRQ(ierr);
987   ierr = VecDuplicate(ts->vec_sol,&ros->Ydot);CHKERRQ(ierr);
988   ierr = VecDuplicate(ts->vec_sol,&ros->Ystage);CHKERRQ(ierr);
989   ierr = VecDuplicate(ts->vec_sol,&ros->Zdot);CHKERRQ(ierr);
990   ierr = VecDuplicate(ts->vec_sol,&ros->Zstage);CHKERRQ(ierr);
991   ierr = VecDuplicate(ts->vec_sol,&ros->VecSolPrev);CHKERRQ(ierr);
992   ierr = PetscMalloc(s*sizeof(ros->work[0]),&ros->work);CHKERRQ(ierr);
993   PetscFunctionReturn(0);
994 }
995 /*------------------------------------------------------------*/
996 
997 #undef __FUNCT__
998 #define __FUNCT__ "TSSetFromOptions_RosW"
999 static PetscErrorCode TSSetFromOptions_RosW(TS ts)
1000 {
1001   TS_RosW        *ros = (TS_RosW*)ts->data;
1002   PetscErrorCode ierr;
1003   char           rostype[256];
1004 
1005   PetscFunctionBegin;
1006   ierr = PetscOptionsHead("RosW ODE solver options");CHKERRQ(ierr);
1007   {
1008     RosWTableauLink link;
1009     PetscInt count,choice;
1010     PetscBool flg;
1011     const char **namelist;
1012     SNES snes;
1013 
1014     ierr = PetscStrncpy(rostype,TSRosWDefault,sizeof rostype);CHKERRQ(ierr);
1015     for (link=RosWTableauList,count=0; link; link=link->next,count++) ;
1016     ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr);
1017     for (link=RosWTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
1018     ierr = PetscOptionsEList("-ts_rosw_type","Family of Rosenbrock-W method","TSRosWSetType",(const char*const*)namelist,count,rostype,&choice,&flg);CHKERRQ(ierr);
1019     ierr = TSRosWSetType(ts,flg ? namelist[choice] : rostype);CHKERRQ(ierr);
1020     ierr = PetscFree(namelist);CHKERRQ(ierr);
1021 
1022     ierr = PetscOptionsBool("-ts_rosw_recompute_jacobian","Recompute the Jacobian at each stage","TSRosWSetRecomputeJacobian",ros->recompute_jacobian,&ros->recompute_jacobian,PETSC_NULL);CHKERRQ(ierr);
1023 
1024     /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
1025     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1026     if (!((PetscObject)snes)->type_name) {
1027       ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);
1028     }
1029     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
1030   }
1031   ierr = PetscOptionsTail();CHKERRQ(ierr);
1032   PetscFunctionReturn(0);
1033 }
1034 
1035 #undef __FUNCT__
1036 #define __FUNCT__ "PetscFormatRealArray"
1037 static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
1038 {
1039   PetscErrorCode ierr;
1040   PetscInt i;
1041   size_t left,count;
1042   char *p;
1043 
1044   PetscFunctionBegin;
1045   for (i=0,p=buf,left=len; i<n; i++) {
1046     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
1047     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
1048     left -= count;
1049     p += count;
1050     *p++ = ' ';
1051   }
1052   p[i ? 0 : -1] = 0;
1053   PetscFunctionReturn(0);
1054 }
1055 
1056 #undef __FUNCT__
1057 #define __FUNCT__ "TSView_RosW"
1058 static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer)
1059 {
1060   TS_RosW        *ros = (TS_RosW*)ts->data;
1061   RosWTableau    tab  = ros->tableau;
1062   PetscBool      iascii;
1063   PetscErrorCode ierr;
1064 
1065   PetscFunctionBegin;
1066   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1067   if (iascii) {
1068     const TSRosWType rostype;
1069     PetscInt i;
1070     PetscReal abscissa[512];
1071     char buf[512];
1072     ierr = TSRosWGetType(ts,&rostype);CHKERRQ(ierr);
1073     ierr = PetscViewerASCIIPrintf(viewer,"  Rosenbrock-W %s\n",rostype);CHKERRQ(ierr);
1074     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ASum);CHKERRQ(ierr);
1075     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A       = %s\n",buf);CHKERRQ(ierr);
1076     for (i=0; i<tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i];
1077     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,abscissa);CHKERRQ(ierr);
1078     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A+Gamma = %s\n",buf);CHKERRQ(ierr);
1079   }
1080   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
1081   PetscFunctionReturn(0);
1082 }
1083 
1084 #undef __FUNCT__
1085 #define __FUNCT__ "TSRosWSetType"
1086 /*@C
1087   TSRosWSetType - Set the type of Rosenbrock-W scheme
1088 
1089   Logically collective
1090 
1091   Input Parameter:
1092 +  ts - timestepping context
1093 -  rostype - type of Rosenbrock-W scheme
1094 
1095   Level: beginner
1096 
1097 .seealso: TSRosWGetType(), TSROSW, TSROSW2M, TSROSW2P, TSROSWRA3PW, TSROSWRA34PW2, TSROSWRODAS3, TSROSWSANDU3, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSROSWLLSSP3P4S2C, TSROSWARK3
1098 @*/
1099 PetscErrorCode TSRosWSetType(TS ts,const TSRosWType rostype)
1100 {
1101   PetscErrorCode ierr;
1102 
1103   PetscFunctionBegin;
1104   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1105   ierr = PetscTryMethod(ts,"TSRosWSetType_C",(TS,const TSRosWType),(ts,rostype));CHKERRQ(ierr);
1106   PetscFunctionReturn(0);
1107 }
1108 
1109 #undef __FUNCT__
1110 #define __FUNCT__ "TSRosWGetType"
1111 /*@C
1112   TSRosWGetType - Get the type of Rosenbrock-W scheme
1113 
1114   Logically collective
1115 
1116   Input Parameter:
1117 .  ts - timestepping context
1118 
1119   Output Parameter:
1120 .  rostype - type of Rosenbrock-W scheme
1121 
1122   Level: intermediate
1123 
1124 .seealso: TSRosWGetType()
1125 @*/
1126 PetscErrorCode TSRosWGetType(TS ts,const TSRosWType *rostype)
1127 {
1128   PetscErrorCode ierr;
1129 
1130   PetscFunctionBegin;
1131   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1132   ierr = PetscUseMethod(ts,"TSRosWGetType_C",(TS,const TSRosWType*),(ts,rostype));CHKERRQ(ierr);
1133   PetscFunctionReturn(0);
1134 }
1135 
1136 #undef __FUNCT__
1137 #define __FUNCT__ "TSRosWSetRecomputeJacobian"
1138 /*@C
1139   TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step.
1140 
1141   Logically collective
1142 
1143   Input Parameter:
1144 +  ts - timestepping context
1145 -  flg - PETSC_TRUE to recompute the Jacobian at each stage
1146 
1147   Level: intermediate
1148 
1149 .seealso: TSRosWGetType()
1150 @*/
1151 PetscErrorCode TSRosWSetRecomputeJacobian(TS ts,PetscBool flg)
1152 {
1153   PetscErrorCode ierr;
1154 
1155   PetscFunctionBegin;
1156   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1157   ierr = PetscTryMethod(ts,"TSRosWSetRecomputeJacobian_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
1158   PetscFunctionReturn(0);
1159 }
1160 
1161 EXTERN_C_BEGIN
1162 #undef __FUNCT__
1163 #define __FUNCT__ "TSRosWGetType_RosW"
1164 PetscErrorCode  TSRosWGetType_RosW(TS ts,const TSRosWType *rostype)
1165 {
1166   TS_RosW        *ros = (TS_RosW*)ts->data;
1167   PetscErrorCode ierr;
1168 
1169   PetscFunctionBegin;
1170   if (!ros->tableau) {ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);}
1171   *rostype = ros->tableau->name;
1172   PetscFunctionReturn(0);
1173 }
1174 #undef __FUNCT__
1175 #define __FUNCT__ "TSRosWSetType_RosW"
1176 PetscErrorCode  TSRosWSetType_RosW(TS ts,const TSRosWType rostype)
1177 {
1178   TS_RosW         *ros = (TS_RosW*)ts->data;
1179   PetscErrorCode  ierr;
1180   PetscBool       match;
1181   RosWTableauLink link;
1182 
1183   PetscFunctionBegin;
1184   if (ros->tableau) {
1185     ierr = PetscStrcmp(ros->tableau->name,rostype,&match);CHKERRQ(ierr);
1186     if (match) PetscFunctionReturn(0);
1187   }
1188   for (link = RosWTableauList; link; link=link->next) {
1189     ierr = PetscStrcmp(link->tab.name,rostype,&match);CHKERRQ(ierr);
1190     if (match) {
1191       ierr = TSReset_RosW(ts);CHKERRQ(ierr);
1192       ros->tableau = &link->tab;
1193       PetscFunctionReturn(0);
1194     }
1195   }
1196   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rostype);
1197   PetscFunctionReturn(0);
1198 }
1199 
1200 #undef __FUNCT__
1201 #define __FUNCT__ "TSRosWSetRecomputeJacobian_RosW"
1202 PetscErrorCode  TSRosWSetRecomputeJacobian_RosW(TS ts,PetscBool flg)
1203 {
1204   TS_RosW *ros = (TS_RosW*)ts->data;
1205 
1206   PetscFunctionBegin;
1207   ros->recompute_jacobian = flg;
1208   PetscFunctionReturn(0);
1209 }
1210 EXTERN_C_END
1211 
1212 /* ------------------------------------------------------------ */
1213 /*MC
1214       TSROSW - ODE solver using Rosenbrock-W schemes
1215 
1216   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
1217   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
1218   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
1219 
1220   Notes:
1221   This method currently only works with autonomous ODE and DAE.
1222 
1223   Developer notes:
1224   Rosenbrock-W methods are typically specified for autonomous ODE
1225 
1226 $  xdot = f(x)
1227 
1228   by the stage equations
1229 
1230 $  k_i = h f(x_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j
1231 
1232   and step completion formula
1233 
1234 $  x_1 = x_0 + sum_j b_j k_j
1235 
1236   with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(x)
1237   and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner,
1238   we define new variables for the stage equations
1239 
1240 $  y_i = gamma_ij k_j
1241 
1242   The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define
1243 
1244 $  A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-i}
1245 
1246   to rewrite the method as
1247 
1248 $  [M/(h gamma_ii) - J] y_i = f(x_0 + sum_j a_ij y_j) + M sum_j (c_ij/h) y_j
1249 $  x_1 = x_0 + sum_j bt_j y_j
1250 
1251    where we have introduced the mass matrix M. Continue by defining
1252 
1253 $  ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j
1254 
1255    or, more compactly in tensor notation
1256 
1257 $  Ydot = 1/h (Gamma^{-1} \otimes I) Y .
1258 
1259    Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current
1260    stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the
1261    equation
1262 
1263 $  g(x_0 + sum_j a_ij y_j + y_i, ydot_i) = 0
1264 
1265    with initial guess y_i = 0.
1266 
1267   Level: beginner
1268 
1269 .seealso:  TSCreate(), TS, TSSetType(), TSRosWSetType(), TSRosWRegister()
1270 
1271 M*/
1272 EXTERN_C_BEGIN
1273 #undef __FUNCT__
1274 #define __FUNCT__ "TSCreate_RosW"
1275 PetscErrorCode  TSCreate_RosW(TS ts)
1276 {
1277   TS_RosW        *ros;
1278   PetscErrorCode ierr;
1279 
1280   PetscFunctionBegin;
1281 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
1282   ierr = TSRosWInitializePackage(PETSC_NULL);CHKERRQ(ierr);
1283 #endif
1284 
1285   ts->ops->reset          = TSReset_RosW;
1286   ts->ops->destroy        = TSDestroy_RosW;
1287   ts->ops->view           = TSView_RosW;
1288   ts->ops->setup          = TSSetUp_RosW;
1289   ts->ops->step           = TSStep_RosW;
1290   ts->ops->interpolate    = TSInterpolate_RosW;
1291   ts->ops->evaluatestep   = TSEvaluateStep_RosW;
1292   ts->ops->setfromoptions = TSSetFromOptions_RosW;
1293   ts->ops->snesfunction   = SNESTSFormFunction_RosW;
1294   ts->ops->snesjacobian   = SNESTSFormJacobian_RosW;
1295 
1296   ierr = PetscNewLog(ts,TS_RosW,&ros);CHKERRQ(ierr);
1297   ts->data = (void*)ros;
1298 
1299   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","TSRosWGetType_RosW",TSRosWGetType_RosW);CHKERRQ(ierr);
1300   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","TSRosWSetType_RosW",TSRosWSetType_RosW);CHKERRQ(ierr);
1301   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","TSRosWSetRecomputeJacobian_RosW",TSRosWSetRecomputeJacobian_RosW);CHKERRQ(ierr);
1302   PetscFunctionReturn(0);
1303 }
1304 EXTERN_C_END
1305