xref: /petsc/src/ts/impls/explicit/ssp/ssp.c (revision 2f6eced2a19e978d64f27de66fbc6c26cc5c7934)
1 /*
2        Code for Timestepping with explicit SSP.
3 */
4 #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
5 
6 PetscFunctionList TSSSPList = 0;
7 static PetscBool TSSSPPackageInitialized;
8 
9 typedef struct {
10   PetscErrorCode (*onestep)(TS,PetscReal,PetscReal,Vec);
11   char           *type_name;
12   PetscInt       nstages;
13   Vec            *work;
14   PetscInt       nwork;
15   PetscBool      workout;
16 } TS_SSP;
17 
18 
19 static PetscErrorCode TSSSPGetWorkVectors(TS ts,PetscInt n,Vec **work)
20 {
21   TS_SSP         *ssp = (TS_SSP*)ts->data;
22   PetscErrorCode ierr;
23 
24   PetscFunctionBegin;
25   if (ssp->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Work vectors already gotten");
26   if (ssp->nwork < n) {
27     if (ssp->nwork > 0) {
28       ierr = VecDestroyVecs(ssp->nwork,&ssp->work);CHKERRQ(ierr);
29     }
30     ierr = VecDuplicateVecs(ts->vec_sol,n,&ssp->work);CHKERRQ(ierr);
31     ssp->nwork = n;
32   }
33   *work = ssp->work;
34   ssp->workout = PETSC_TRUE;
35   PetscFunctionReturn(0);
36 }
37 
38 static PetscErrorCode TSSSPRestoreWorkVectors(TS ts,PetscInt n,Vec **work)
39 {
40   TS_SSP *ssp = (TS_SSP*)ts->data;
41 
42   PetscFunctionBegin;
43   if (!ssp->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Work vectors have not been gotten");
44   if (*work != ssp->work) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong work vectors checked out");
45   ssp->workout = PETSC_FALSE;
46   *work = NULL;
47   PetscFunctionReturn(0);
48 }
49 
50 /*MC
51    TSSSPRKS2 - Optimal second order SSP Runge-Kutta method, low-storage, c_eff=(s-1)/s
52 
53    Pseudocode 2 of Ketcheson 2008
54 
55    Level: beginner
56 
57 .seealso: TSSSP, TSSSPSetType(), TSSSPSetNumStages()
58 M*/
59 static PetscErrorCode TSSSPStep_RK_2(TS ts,PetscReal t0,PetscReal dt,Vec sol)
60 {
61   TS_SSP         *ssp = (TS_SSP*)ts->data;
62   Vec            *work,F;
63   PetscInt       i,s;
64   PetscErrorCode ierr;
65 
66   PetscFunctionBegin;
67   s    = ssp->nstages;
68   ierr = TSSSPGetWorkVectors(ts,2,&work);CHKERRQ(ierr);
69   F    = work[1];
70   ierr = VecCopy(sol,work[0]);CHKERRQ(ierr);
71   for (i=0; i<s-1; i++) {
72     PetscReal stage_time = t0+dt*(i/(s-1.));
73     ierr = TSPreStage(ts,stage_time);CHKERRQ(ierr);
74     ierr = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr);
75     ierr = VecAXPY(work[0],dt/(s-1.),F);CHKERRQ(ierr);
76   }
77   ierr = TSComputeRHSFunction(ts,t0+dt,work[0],F);CHKERRQ(ierr);
78   ierr = VecAXPBYPCZ(sol,(s-1.)/s,dt/s,1./s,work[0],F);CHKERRQ(ierr);
79   ierr = TSSSPRestoreWorkVectors(ts,2,&work);CHKERRQ(ierr);
80   PetscFunctionReturn(0);
81 }
82 
83 /*MC
84    TSSSPRKS3 - Optimal third order SSP Runge-Kutta, low-storage, c_eff=(PetscSqrtReal(s)-1)/PetscSqrtReal(s), where PetscSqrtReal(s) is an integer
85 
86    Pseudocode 2 of Ketcheson 2008
87 
88    Level: beginner
89 
90 .seealso: TSSSP, TSSSPSetType(), TSSSPSetNumStages()
91 M*/
92 static PetscErrorCode TSSSPStep_RK_3(TS ts,PetscReal t0,PetscReal dt,Vec sol)
93 {
94   TS_SSP         *ssp = (TS_SSP*)ts->data;
95   Vec            *work,F;
96   PetscInt       i,s,n,r;
97   PetscReal      c,stage_time;
98   PetscErrorCode ierr;
99 
100   PetscFunctionBegin;
101   s = ssp->nstages;
102   n = (PetscInt)(PetscSqrtReal((PetscReal)s)+0.001);
103   r = s-n;
104   if (n*n != s) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for optimal third order schemes with %d stages, must be a square number at least 4",s);
105   ierr = TSSSPGetWorkVectors(ts,3,&work);CHKERRQ(ierr);
106   F    = work[2];
107   ierr = VecCopy(sol,work[0]);CHKERRQ(ierr);
108   for (i=0; i<(n-1)*(n-2)/2; i++) {
109     c          = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
110     stage_time = t0+c*dt;
111     ierr       = TSPreStage(ts,stage_time);CHKERRQ(ierr);
112     ierr       = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr);
113     ierr       = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr);
114   }
115   ierr = VecCopy(work[0],work[1]);CHKERRQ(ierr);
116   for (; i<n*(n+1)/2-1; i++) {
117     c          = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
118     stage_time = t0+c*dt;
119     ierr       = TSPreStage(ts,stage_time);CHKERRQ(ierr);
120     ierr       = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr);
121     ierr       = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr);
122   }
123   {
124     c          = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
125     stage_time = t0+c*dt;
126     ierr       = TSPreStage(ts,stage_time);CHKERRQ(ierr);
127     ierr       = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr);
128     ierr       = VecAXPBYPCZ(work[0],1.*n/(2*n-1.),(n-1.)*dt/(r*(2*n-1)),(n-1.)/(2*n-1.),work[1],F);CHKERRQ(ierr);
129     i++;
130   }
131   for (; i<s; i++) {
132     c          = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
133     stage_time = t0+c*dt;
134     ierr       = TSPreStage(ts,stage_time);CHKERRQ(ierr);
135     ierr       = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr);
136     ierr       = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr);
137   }
138   ierr = VecCopy(work[0],sol);CHKERRQ(ierr);
139   ierr = TSSSPRestoreWorkVectors(ts,3,&work);CHKERRQ(ierr);
140   PetscFunctionReturn(0);
141 }
142 
143 /*MC
144    TSSSPRKS104 - Optimal fourth order SSP Runge-Kutta, low-storage (2N), c_eff=0.6
145 
146    SSPRK(10,4), Pseudocode 3 of Ketcheson 2008
147 
148    Level: beginner
149 
150 .seealso: TSSSP, TSSSPSetType()
151 M*/
152 static PetscErrorCode TSSSPStep_RK_10_4(TS ts,PetscReal t0,PetscReal dt,Vec sol)
153 {
154   const PetscReal c[10] = {0, 1./6, 2./6, 3./6, 4./6, 2./6, 3./6, 4./6, 5./6, 1};
155   Vec             *work,F;
156   PetscInt        i;
157   PetscReal       stage_time;
158   PetscErrorCode  ierr;
159 
160   PetscFunctionBegin;
161   ierr = TSSSPGetWorkVectors(ts,3,&work);CHKERRQ(ierr);
162   F    = work[2];
163   ierr = VecCopy(sol,work[0]);CHKERRQ(ierr);
164   for (i=0; i<5; i++) {
165     stage_time = t0+c[i]*dt;
166     ierr       = TSPreStage(ts,stage_time);CHKERRQ(ierr);
167     ierr       = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr);
168     ierr       = VecAXPY(work[0],dt/6,F);CHKERRQ(ierr);
169   }
170   ierr = VecAXPBYPCZ(work[1],1./25,9./25,0,sol,work[0]);CHKERRQ(ierr);
171   ierr = VecAXPBY(work[0],15,-5,work[1]);CHKERRQ(ierr);
172   for (; i<9; i++) {
173     stage_time = t0+c[i]*dt;
174     ierr       = TSPreStage(ts,stage_time);CHKERRQ(ierr);
175     ierr       = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr);
176     ierr       = VecAXPY(work[0],dt/6,F);CHKERRQ(ierr);
177   }
178   stage_time = t0+dt;
179   ierr       = TSPreStage(ts,stage_time);CHKERRQ(ierr);
180   ierr       = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr);
181   ierr       = VecAXPBYPCZ(work[1],3./5,dt/10,1,work[0],F);CHKERRQ(ierr);
182   ierr       = VecCopy(work[1],sol);CHKERRQ(ierr);
183   ierr       = TSSSPRestoreWorkVectors(ts,3,&work);CHKERRQ(ierr);
184   PetscFunctionReturn(0);
185 }
186 
187 
188 static PetscErrorCode TSSetUp_SSP(TS ts)
189 {
190   PetscErrorCode ierr;
191 
192   PetscFunctionBegin;
193   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
194   ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr);
195   PetscFunctionReturn(0);
196 }
197 
198 static PetscErrorCode TSStep_SSP(TS ts)
199 {
200   TS_SSP         *ssp = (TS_SSP*)ts->data;
201   Vec            sol  = ts->vec_sol;
202   PetscBool      stageok,accept = PETSC_TRUE;
203   PetscReal      next_time_step = ts->time_step;
204   PetscErrorCode ierr;
205 
206   PetscFunctionBegin;
207   ierr = (*ssp->onestep)(ts,ts->ptime,ts->time_step,sol);CHKERRQ(ierr);
208   ierr = TSPostStage(ts,ts->ptime,0,&sol);CHKERRQ(ierr);
209   ierr = TSAdaptCheckStage(ts->adapt,ts,ts->ptime+ts->time_step,sol,&stageok);CHKERRQ(ierr);
210   if(!stageok) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);}
211 
212   ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
213   if (!accept) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);}
214 
215   ts->ptime += ts->time_step;
216   ts->time_step = next_time_step;
217   PetscFunctionReturn(0);
218 }
219 /*------------------------------------------------------------*/
220 
221 static PetscErrorCode TSReset_SSP(TS ts)
222 {
223   TS_SSP         *ssp = (TS_SSP*)ts->data;
224   PetscErrorCode ierr;
225 
226   PetscFunctionBegin;
227   if (ssp->work) {ierr = VecDestroyVecs(ssp->nwork,&ssp->work);CHKERRQ(ierr);}
228   ssp->nwork   = 0;
229   ssp->workout = PETSC_FALSE;
230   PetscFunctionReturn(0);
231 }
232 
233 static PetscErrorCode TSDestroy_SSP(TS ts)
234 {
235   TS_SSP         *ssp = (TS_SSP*)ts->data;
236   PetscErrorCode ierr;
237 
238   PetscFunctionBegin;
239   ierr = TSReset_SSP(ts);CHKERRQ(ierr);
240   ierr = PetscFree(ssp->type_name);CHKERRQ(ierr);
241   ierr = PetscFree(ts->data);CHKERRQ(ierr);
242   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetType_C",NULL);CHKERRQ(ierr);
243   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetType_C",NULL);CHKERRQ(ierr);
244   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetNumStages_C",NULL);CHKERRQ(ierr);
245   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetNumStages_C",NULL);CHKERRQ(ierr);
246   PetscFunctionReturn(0);
247 }
248 /*------------------------------------------------------------*/
249 
250 /*@C
251    TSSSPSetType - set the SSP time integration scheme to use
252 
253    Logically Collective
254 
255    Input Arguments:
256    ts - time stepping object
257    type - type of scheme to use
258 
259    Options Database Keys:
260    -ts_ssp_type <rks2>: Type of SSP method (one of) rks2 rks3 rk104
261    -ts_ssp_nstages <5>: Number of stages
262 
263    Level: beginner
264 
265 .seealso: TSSSP, TSSSPGetType(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104
266 @*/
267 PetscErrorCode TSSSPSetType(TS ts,TSSSPType type)
268 {
269   PetscErrorCode ierr;
270 
271   PetscFunctionBegin;
272   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
273   ierr = PetscTryMethod(ts,"TSSSPSetType_C",(TS,TSSSPType),(ts,type));CHKERRQ(ierr);
274   PetscFunctionReturn(0);
275 }
276 
277 /*@C
278    TSSSPGetType - get the SSP time integration scheme
279 
280    Logically Collective
281 
282    Input Argument:
283    ts - time stepping object
284 
285    Output Argument:
286    type - type of scheme being used
287 
288    Level: beginner
289 
290 .seealso: TSSSP, TSSSPSettype(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104
291 @*/
292 PetscErrorCode TSSSPGetType(TS ts,TSSSPType *type)
293 {
294   PetscErrorCode ierr;
295 
296   PetscFunctionBegin;
297   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
298   ierr = PetscUseMethod(ts,"TSSSPGetType_C",(TS,TSSSPType*),(ts,type));CHKERRQ(ierr);
299   PetscFunctionReturn(0);
300 }
301 
302 /*@
303    TSSSPSetNumStages - set the number of stages to use with the SSP method
304 
305    Logically Collective
306 
307    Input Arguments:
308    ts - time stepping object
309    nstages - number of stages
310 
311    Options Database Keys:
312    -ts_ssp_type <rks2>: NumStages of SSP method (one of) rks2 rks3 rk104
313    -ts_ssp_nstages <5>: Number of stages
314 
315    Level: beginner
316 
317 .seealso: TSSSP, TSSSPGetNumStages(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104
318 @*/
319 PetscErrorCode TSSSPSetNumStages(TS ts,PetscInt nstages)
320 {
321   PetscErrorCode ierr;
322 
323   PetscFunctionBegin;
324   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
325   ierr = PetscTryMethod(ts,"TSSSPSetNumStages_C",(TS,PetscInt),(ts,nstages));CHKERRQ(ierr);
326   PetscFunctionReturn(0);
327 }
328 
329 /*@
330    TSSSPGetNumStages - get the number of stages in the SSP time integration scheme
331 
332    Logically Collective
333 
334    Input Argument:
335    ts - time stepping object
336 
337    Output Argument:
338    nstages - number of stages
339 
340    Level: beginner
341 
342 .seealso: TSSSP, TSSSPGetType(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104
343 @*/
344 PetscErrorCode TSSSPGetNumStages(TS ts,PetscInt *nstages)
345 {
346   PetscErrorCode ierr;
347 
348   PetscFunctionBegin;
349   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
350   ierr = PetscUseMethod(ts,"TSSSPGetNumStages_C",(TS,PetscInt*),(ts,nstages));CHKERRQ(ierr);
351   PetscFunctionReturn(0);
352 }
353 
354 static PetscErrorCode TSSSPSetType_SSP(TS ts,TSSSPType type)
355 {
356   PetscErrorCode ierr,(*r)(TS,PetscReal,PetscReal,Vec);
357   TS_SSP         *ssp = (TS_SSP*)ts->data;
358 
359   PetscFunctionBegin;
360   ierr = PetscFunctionListFind(TSSSPList,type,&r);CHKERRQ(ierr);
361   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TS_SSP type %s given",type);
362   ssp->onestep = r;
363   ierr = PetscFree(ssp->type_name);CHKERRQ(ierr);
364   ierr = PetscStrallocpy(type,&ssp->type_name);CHKERRQ(ierr);
365   PetscFunctionReturn(0);
366 }
367 static PetscErrorCode TSSSPGetType_SSP(TS ts,TSSSPType *type)
368 {
369   TS_SSP *ssp = (TS_SSP*)ts->data;
370 
371   PetscFunctionBegin;
372   *type = ssp->type_name;
373   PetscFunctionReturn(0);
374 }
375 static PetscErrorCode TSSSPSetNumStages_SSP(TS ts,PetscInt nstages)
376 {
377   TS_SSP *ssp = (TS_SSP*)ts->data;
378 
379   PetscFunctionBegin;
380   ssp->nstages = nstages;
381   PetscFunctionReturn(0);
382 }
383 static PetscErrorCode TSSSPGetNumStages_SSP(TS ts,PetscInt *nstages)
384 {
385   TS_SSP *ssp = (TS_SSP*)ts->data;
386 
387   PetscFunctionBegin;
388   *nstages = ssp->nstages;
389   PetscFunctionReturn(0);
390 }
391 
392 static PetscErrorCode TSSetFromOptions_SSP(PetscOptionItems *PetscOptionsObject,TS ts)
393 {
394   char           tname[256] = TSSSPRKS2;
395   TS_SSP         *ssp       = (TS_SSP*)ts->data;
396   PetscErrorCode ierr;
397   PetscBool      flg;
398 
399   PetscFunctionBegin;
400   ierr = PetscOptionsHead(PetscOptionsObject,"SSP ODE solver options");CHKERRQ(ierr);
401   {
402     ierr = PetscOptionsFList("-ts_ssp_type","Type of SSP method","TSSSPSetType",TSSSPList,tname,tname,sizeof(tname),&flg);CHKERRQ(ierr);
403     if (flg) {
404       ierr = TSSSPSetType(ts,tname);CHKERRQ(ierr);
405     }
406     ierr = PetscOptionsInt("-ts_ssp_nstages","Number of stages","TSSSPSetNumStages",ssp->nstages,&ssp->nstages,NULL);CHKERRQ(ierr);
407   }
408   ierr = PetscOptionsTail();CHKERRQ(ierr);
409   PetscFunctionReturn(0);
410 }
411 
412 static PetscErrorCode TSView_SSP(TS ts,PetscViewer viewer)
413 {
414   TS_SSP         *ssp = (TS_SSP*)ts->data;
415   PetscBool      ascii;
416   PetscErrorCode ierr;
417 
418   PetscFunctionBegin;
419   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&ascii);CHKERRQ(ierr);
420   if (ascii) {ierr = PetscViewerASCIIPrintf(viewer,"  Scheme: %s\n",ssp->type_name);CHKERRQ(ierr);}
421   PetscFunctionReturn(0);
422 }
423 
424 /* ------------------------------------------------------------ */
425 
426 /*MC
427       TSSSP - Explicit strong stability preserving ODE solver
428 
429   Most hyperbolic conservation laws have exact solutions that are total variation diminishing (TVD) or total variation
430   bounded (TVB) although these solutions often contain discontinuities.  Spatial discretizations such as Godunov's
431   scheme and high-resolution finite volume methods (TVD limiters, ENO/WENO) are designed to preserve these properties,
432   but they are usually formulated using a forward Euler time discretization or by coupling the space and time
433   discretization as in the classical Lax-Wendroff scheme.  When the space and time discretization is coupled, it is very
434   difficult to produce schemes with high temporal accuracy while preserving TVD properties.  An alternative is the
435   semidiscrete formulation where we choose a spatial discretization that is TVD with forward Euler and then choose a
436   time discretization that preserves the TVD property.  Such integrators are called strong stability preserving (SSP).
437 
438   Let c_eff be the minimum number of function evaluations required to step as far as one step of forward Euler while
439   still being SSP.  Some theoretical bounds
440 
441   1. There are no explicit methods with c_eff > 1.
442 
443   2. There are no explicit methods beyond order 4 (for nonlinear problems) and c_eff > 0.
444 
445   3. There are no implicit methods with order greater than 1 and c_eff > 2.
446 
447   This integrator provides Runge-Kutta methods of order 2, 3, and 4 with maximal values of c_eff.  More stages allows
448   for larger values of c_eff which improves efficiency.  These implementations are low-memory and only use 2 or 3 work
449   vectors regardless of the total number of stages, so e.g. 25-stage 3rd order methods may be an excellent choice.
450 
451   Methods can be chosen with -ts_ssp_type {rks2,rks3,rk104}
452 
453   rks2: Second order methods with any number s>1 of stages.  c_eff = (s-1)/s
454 
455   rks3: Third order methods with s=n^2 stages, n>1.  c_eff = (s-n)/s
456 
457   rk104: A 10-stage fourth order method.  c_eff = 0.6
458 
459   Level: beginner
460 
461   References:
462 +  1. - Ketcheson, Highly efficient strong stability preserving Runge Kutta methods with low storage implementations, SISC, 2008.
463 -  2. - Gottlieb, Ketcheson, and Shu, High order strong stability preserving time discretizations, J Scientific Computing, 2009.
464 
465 .seealso:  TSCreate(), TS, TSSetType()
466 
467 M*/
468 PETSC_EXTERN PetscErrorCode TSCreate_SSP(TS ts)
469 {
470   TS_SSP         *ssp;
471   PetscErrorCode ierr;
472 
473   PetscFunctionBegin;
474   ierr = TSSSPInitializePackage();CHKERRQ(ierr);
475 
476   ts->ops->setup          = TSSetUp_SSP;
477   ts->ops->step           = TSStep_SSP;
478   ts->ops->reset          = TSReset_SSP;
479   ts->ops->destroy        = TSDestroy_SSP;
480   ts->ops->setfromoptions = TSSetFromOptions_SSP;
481   ts->ops->view           = TSView_SSP;
482 
483   ierr = PetscNewLog(ts,&ssp);CHKERRQ(ierr);
484   ts->data = (void*)ssp;
485 
486   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetType_C",TSSSPGetType_SSP);CHKERRQ(ierr);
487   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetType_C",TSSSPSetType_SSP);CHKERRQ(ierr);
488   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetNumStages_C",TSSSPGetNumStages_SSP);CHKERRQ(ierr);
489   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetNumStages_C",TSSSPSetNumStages_SSP);CHKERRQ(ierr);
490 
491   ierr = TSSSPSetType(ts,TSSSPRKS2);CHKERRQ(ierr);
492   ssp->nstages = 5;
493   PetscFunctionReturn(0);
494 }
495 
496 /*@C
497   TSSSPInitializePackage - This function initializes everything in the TSSSP package. It is called
498   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_SSP()
499   when using static libraries.
500 
501   Level: developer
502 
503 .keywords: TS, TSSSP, initialize, package
504 .seealso: PetscInitialize()
505 @*/
506 PetscErrorCode TSSSPInitializePackage(void)
507 {
508   PetscErrorCode ierr;
509 
510   PetscFunctionBegin;
511   if (TSSSPPackageInitialized) PetscFunctionReturn(0);
512   TSSSPPackageInitialized = PETSC_TRUE;
513   ierr = PetscFunctionListAdd(&TSSSPList,TSSSPRKS2, TSSSPStep_RK_2);CHKERRQ(ierr);
514   ierr = PetscFunctionListAdd(&TSSSPList,TSSSPRKS3, TSSSPStep_RK_3);CHKERRQ(ierr);
515   ierr = PetscFunctionListAdd(&TSSSPList,TSSSPRK104,TSSSPStep_RK_10_4);CHKERRQ(ierr);
516   ierr = PetscRegisterFinalize(TSSSPFinalizePackage);CHKERRQ(ierr);
517   PetscFunctionReturn(0);
518 }
519 
520 /*@C
521   TSSSPFinalizePackage - This function destroys everything in the TSSSP package. It is
522   called from PetscFinalize().
523 
524   Level: developer
525 
526 .keywords: Petsc, destroy, package
527 .seealso: PetscFinalize()
528 @*/
529 PetscErrorCode TSSSPFinalizePackage(void)
530 {
531   PetscErrorCode ierr;
532 
533   PetscFunctionBegin;
534   TSSSPPackageInitialized = PETSC_FALSE;
535   ierr = PetscFunctionListDestroy(&TSSSPList);CHKERRQ(ierr);
536   PetscFunctionReturn(0);
537 }
538