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