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