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