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