xref: /petsc/src/ts/impls/explicit/ssp/ssp.c (revision 9c334d8fdb557fc53fd345d68cbb3545b09ccab8)
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 #undef __FUNCT__
20 #define __FUNCT__ "TSSSPGetWorkVectors"
21 static PetscErrorCode TSSSPGetWorkVectors(TS ts,PetscInt n,Vec **work)
22 {
23   TS_SSP         *ssp = (TS_SSP*)ts->data;
24   PetscErrorCode ierr;
25 
26   PetscFunctionBegin;
27   if (ssp->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Work vectors already gotten");
28   if (ssp->nwork < n) {
29     if (ssp->nwork > 0) {
30       ierr = VecDestroyVecs(ssp->nwork,&ssp->work);CHKERRQ(ierr);
31     }
32     ierr = VecDuplicateVecs(ts->vec_sol,n,&ssp->work);CHKERRQ(ierr);
33     ssp->nwork = n;
34   }
35   *work = ssp->work;
36   ssp->workout = PETSC_TRUE;
37   PetscFunctionReturn(0);
38 }
39 
40 #undef __FUNCT__
41 #define __FUNCT__ "TSSSPRestoreWorkVectors"
42 static PetscErrorCode TSSSPRestoreWorkVectors(TS ts,PetscInt n,Vec **work)
43 {
44   TS_SSP *ssp = (TS_SSP*)ts->data;
45 
46   PetscFunctionBegin;
47   if (!ssp->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Work vectors have not been gotten");
48   if (*work != ssp->work) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong work vectors checked out");
49   ssp->workout = PETSC_FALSE;
50   *work = NULL;
51   PetscFunctionReturn(0);
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   PetscErrorCode ierr;
203 
204   PetscFunctionBegin;
205   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
206   ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr);
207   PetscFunctionReturn(0);
208 }
209 
210 #undef __FUNCT__
211 #define __FUNCT__ "TSStep_SSP"
212 static PetscErrorCode TSStep_SSP(TS ts)
213 {
214   TS_SSP         *ssp = (TS_SSP*)ts->data;
215   Vec            sol  = ts->vec_sol;
216   PetscBool      stageok,accept = PETSC_TRUE;
217   PetscReal      next_time_step = ts->time_step;
218   PetscErrorCode ierr;
219 
220   PetscFunctionBegin;
221   ierr = (*ssp->onestep)(ts,ts->ptime,ts->time_step,sol);CHKERRQ(ierr);
222   ierr = TSPostStage(ts,ts->ptime,0,&sol);CHKERRQ(ierr);
223   ierr = TSAdaptCheckStage(ts->adapt,ts,ts->ptime+ts->time_step,sol,&stageok);CHKERRQ(ierr);
224   if(!stageok) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);}
225 
226   ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
227   if (!accept) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);}
228 
229   ts->ptime += ts->time_step;
230   ts->time_step = next_time_step;
231   ts->steps++;
232   PetscFunctionReturn(0);
233 }
234 /*------------------------------------------------------------*/
235 
236 #undef __FUNCT__
237 #define __FUNCT__ "TSReset_SSP"
238 static PetscErrorCode TSReset_SSP(TS ts)
239 {
240   TS_SSP         *ssp = (TS_SSP*)ts->data;
241   PetscErrorCode ierr;
242 
243   PetscFunctionBegin;
244   if (ssp->work) {ierr = VecDestroyVecs(ssp->nwork,&ssp->work);CHKERRQ(ierr);}
245   ssp->nwork   = 0;
246   ssp->workout = PETSC_FALSE;
247   PetscFunctionReturn(0);
248 }
249 
250 #undef __FUNCT__
251 #define __FUNCT__ "TSDestroy_SSP"
252 static PetscErrorCode TSDestroy_SSP(TS ts)
253 {
254   TS_SSP         *ssp = (TS_SSP*)ts->data;
255   PetscErrorCode ierr;
256 
257   PetscFunctionBegin;
258   ierr = TSReset_SSP(ts);CHKERRQ(ierr);
259   ierr = PetscFree(ssp->type_name);CHKERRQ(ierr);
260   ierr = PetscFree(ts->data);CHKERRQ(ierr);
261   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetType_C",NULL);CHKERRQ(ierr);
262   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetType_C",NULL);CHKERRQ(ierr);
263   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetNumStages_C",NULL);CHKERRQ(ierr);
264   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetNumStages_C",NULL);CHKERRQ(ierr);
265   PetscFunctionReturn(0);
266 }
267 /*------------------------------------------------------------*/
268 
269 #undef __FUNCT__
270 #define __FUNCT__ "TSSSPSetType"
271 /*@C
272    TSSSPSetType - set the SSP time integration scheme to use
273 
274    Logically Collective
275 
276    Input Arguments:
277    ts - time stepping object
278    type - type of scheme to use
279 
280    Options Database Keys:
281    -ts_ssp_type <rks2>: Type of SSP method (one of) rks2 rks3 rk104
282    -ts_ssp_nstages <5>: Number of stages
283 
284    Level: beginner
285 
286 .seealso: TSSSP, TSSSPGetType(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104
287 @*/
288 PetscErrorCode TSSSPSetType(TS ts,TSSSPType type)
289 {
290   PetscErrorCode ierr;
291 
292   PetscFunctionBegin;
293   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
294   ierr = PetscTryMethod(ts,"TSSSPSetType_C",(TS,TSSSPType),(ts,type));CHKERRQ(ierr);
295   PetscFunctionReturn(0);
296 }
297 
298 #undef __FUNCT__
299 #define __FUNCT__ "TSSSPGetType"
300 /*@C
301    TSSSPGetType - get the SSP time integration scheme
302 
303    Logically Collective
304 
305    Input Argument:
306    ts - time stepping object
307 
308    Output Argument:
309    type - type of scheme being used
310 
311    Level: beginner
312 
313 .seealso: TSSSP, TSSSPSettype(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104
314 @*/
315 PetscErrorCode TSSSPGetType(TS ts,TSSSPType *type)
316 {
317   PetscErrorCode ierr;
318 
319   PetscFunctionBegin;
320   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
321   ierr = PetscUseMethod(ts,"TSSSPGetType_C",(TS,TSSSPType*),(ts,type));CHKERRQ(ierr);
322   PetscFunctionReturn(0);
323 }
324 
325 #undef __FUNCT__
326 #define __FUNCT__ "TSSSPSetNumStages"
327 /*@
328    TSSSPSetNumStages - set the number of stages to use with the SSP method
329 
330    Logically Collective
331 
332    Input Arguments:
333    ts - time stepping object
334    nstages - number of stages
335 
336    Options Database Keys:
337    -ts_ssp_type <rks2>: NumStages of SSP method (one of) rks2 rks3 rk104
338    -ts_ssp_nstages <5>: Number of stages
339 
340    Level: beginner
341 
342 .seealso: TSSSP, TSSSPGetNumStages(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104
343 @*/
344 PetscErrorCode TSSSPSetNumStages(TS ts,PetscInt nstages)
345 {
346   PetscErrorCode ierr;
347 
348   PetscFunctionBegin;
349   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
350   ierr = PetscTryMethod(ts,"TSSSPSetNumStages_C",(TS,PetscInt),(ts,nstages));CHKERRQ(ierr);
351   PetscFunctionReturn(0);
352 }
353 
354 #undef __FUNCT__
355 #define __FUNCT__ "TSSSPGetNumStages"
356 /*@
357    TSSSPGetNumStages - get the number of stages in the SSP time integration scheme
358 
359    Logically Collective
360 
361    Input Argument:
362    ts - time stepping object
363 
364    Output Argument:
365    nstages - number of stages
366 
367    Level: beginner
368 
369 .seealso: TSSSP, TSSSPGetType(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104
370 @*/
371 PetscErrorCode TSSSPGetNumStages(TS ts,PetscInt *nstages)
372 {
373   PetscErrorCode ierr;
374 
375   PetscFunctionBegin;
376   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
377   ierr = PetscUseMethod(ts,"TSSSPGetNumStages_C",(TS,PetscInt*),(ts,nstages));CHKERRQ(ierr);
378   PetscFunctionReturn(0);
379 }
380 
381 #undef __FUNCT__
382 #define __FUNCT__ "TSSSPSetType_SSP"
383 static PetscErrorCode TSSSPSetType_SSP(TS ts,TSSSPType type)
384 {
385   PetscErrorCode ierr,(*r)(TS,PetscReal,PetscReal,Vec);
386   TS_SSP         *ssp = (TS_SSP*)ts->data;
387 
388   PetscFunctionBegin;
389   ierr = PetscFunctionListFind(TSSSPList,type,&r);CHKERRQ(ierr);
390   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TS_SSP type %s given",type);
391   ssp->onestep = r;
392   ierr = PetscFree(ssp->type_name);CHKERRQ(ierr);
393   ierr = PetscStrallocpy(type,&ssp->type_name);CHKERRQ(ierr);
394   PetscFunctionReturn(0);
395 }
396 #undef __FUNCT__
397 #define __FUNCT__ "TSSSPGetType_SSP"
398 static PetscErrorCode TSSSPGetType_SSP(TS ts,TSSSPType *type)
399 {
400   TS_SSP *ssp = (TS_SSP*)ts->data;
401 
402   PetscFunctionBegin;
403   *type = ssp->type_name;
404   PetscFunctionReturn(0);
405 }
406 #undef __FUNCT__
407 #define __FUNCT__ "TSSSPSetNumStages_SSP"
408 static PetscErrorCode TSSSPSetNumStages_SSP(TS ts,PetscInt nstages)
409 {
410   TS_SSP *ssp = (TS_SSP*)ts->data;
411 
412   PetscFunctionBegin;
413   ssp->nstages = nstages;
414   PetscFunctionReturn(0);
415 }
416 #undef __FUNCT__
417 #define __FUNCT__ "TSSSPGetNumStages_SSP"
418 static PetscErrorCode TSSSPGetNumStages_SSP(TS ts,PetscInt *nstages)
419 {
420   TS_SSP *ssp = (TS_SSP*)ts->data;
421 
422   PetscFunctionBegin;
423   *nstages = ssp->nstages;
424   PetscFunctionReturn(0);
425 }
426 
427 #undef __FUNCT__
428 #define __FUNCT__ "TSSetFromOptions_SSP"
429 static PetscErrorCode TSSetFromOptions_SSP(PetscOptionItems *PetscOptionsObject,TS ts)
430 {
431   char           tname[256] = TSSSPRKS2;
432   TS_SSP         *ssp       = (TS_SSP*)ts->data;
433   PetscErrorCode ierr;
434   PetscBool      flg;
435 
436   PetscFunctionBegin;
437   ierr = PetscOptionsHead(PetscOptionsObject,"SSP ODE solver options");CHKERRQ(ierr);
438   {
439     ierr = PetscOptionsFList("-ts_ssp_type","Type of SSP method","TSSSPSetType",TSSSPList,tname,tname,sizeof(tname),&flg);CHKERRQ(ierr);
440     if (flg) {
441       ierr = TSSSPSetType(ts,tname);CHKERRQ(ierr);
442     }
443     ierr = PetscOptionsInt("-ts_ssp_nstages","Number of stages","TSSSPSetNumStages",ssp->nstages,&ssp->nstages,NULL);CHKERRQ(ierr);
444   }
445   ierr = PetscOptionsTail();CHKERRQ(ierr);
446   PetscFunctionReturn(0);
447 }
448 
449 #undef __FUNCT__
450 #define __FUNCT__ "TSView_SSP"
451 static PetscErrorCode TSView_SSP(TS ts,PetscViewer viewer)
452 {
453   TS_SSP         *ssp = (TS_SSP*)ts->data;
454   PetscBool      ascii;
455   PetscErrorCode ierr;
456 
457   PetscFunctionBegin;
458   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&ascii);CHKERRQ(ierr);
459   if (ascii) {ierr = PetscViewerASCIIPrintf(viewer,"  Scheme: %s\n",ssp->type_name);CHKERRQ(ierr);}
460   PetscFunctionReturn(0);
461 }
462 
463 /* ------------------------------------------------------------ */
464 
465 /*MC
466       TSSSP - Explicit strong stability preserving ODE solver
467 
468   Most hyperbolic conservation laws have exact solutions that are total variation diminishing (TVD) or total variation
469   bounded (TVB) although these solutions often contain discontinuities.  Spatial discretizations such as Godunov's
470   scheme and high-resolution finite volume methods (TVD limiters, ENO/WENO) are designed to preserve these properties,
471   but they are usually formulated using a forward Euler time discretization or by coupling the space and time
472   discretization as in the classical Lax-Wendroff scheme.  When the space and time discretization is coupled, it is very
473   difficult to produce schemes with high temporal accuracy while preserving TVD properties.  An alternative is the
474   semidiscrete formulation where we choose a spatial discretization that is TVD with forward Euler and then choose a
475   time discretization that preserves the TVD property.  Such integrators are called strong stability preserving (SSP).
476 
477   Let c_eff be the minimum number of function evaluations required to step as far as one step of forward Euler while
478   still being SSP.  Some theoretical bounds
479 
480   1. There are no explicit methods with c_eff > 1.
481 
482   2. There are no explicit methods beyond order 4 (for nonlinear problems) and c_eff > 0.
483 
484   3. There are no implicit methods with order greater than 1 and c_eff > 2.
485 
486   This integrator provides Runge-Kutta methods of order 2, 3, and 4 with maximal values of c_eff.  More stages allows
487   for larger values of c_eff which improves efficiency.  These implementations are low-memory and only use 2 or 3 work
488   vectors regardless of the total number of stages, so e.g. 25-stage 3rd order methods may be an excellent choice.
489 
490   Methods can be chosen with -ts_ssp_type {rks2,rks3,rk104}
491 
492   rks2: Second order methods with any number s>1 of stages.  c_eff = (s-1)/s
493 
494   rks3: Third order methods with s=n^2 stages, n>1.  c_eff = (s-n)/s
495 
496   rk104: A 10-stage fourth order method.  c_eff = 0.6
497 
498   Level: beginner
499 
500   References:
501 +  1. - Ketcheson, Highly efficient strong stability preserving Runge Kutta methods with low storage implementations, SISC, 2008.
502 -  2. - Gottlieb, Ketcheson, and Shu, High order strong stability preserving time discretizations, J Scientific Computing, 2009.
503 
504 .seealso:  TSCreate(), TS, TSSetType()
505 
506 M*/
507 #undef __FUNCT__
508 #define __FUNCT__ "TSCreate_SSP"
509 PETSC_EXTERN PetscErrorCode TSCreate_SSP(TS ts)
510 {
511   TS_SSP         *ssp;
512   PetscErrorCode ierr;
513 
514   PetscFunctionBegin;
515   ierr = TSSSPInitializePackage();CHKERRQ(ierr);
516 
517   ts->ops->setup          = TSSetUp_SSP;
518   ts->ops->step           = TSStep_SSP;
519   ts->ops->reset          = TSReset_SSP;
520   ts->ops->destroy        = TSDestroy_SSP;
521   ts->ops->setfromoptions = TSSetFromOptions_SSP;
522   ts->ops->view           = TSView_SSP;
523 
524   ierr = PetscNewLog(ts,&ssp);CHKERRQ(ierr);
525   ts->data = (void*)ssp;
526 
527   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetType_C",TSSSPGetType_SSP);CHKERRQ(ierr);
528   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetType_C",TSSSPSetType_SSP);CHKERRQ(ierr);
529   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetNumStages_C",TSSSPGetNumStages_SSP);CHKERRQ(ierr);
530   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetNumStages_C",TSSSPSetNumStages_SSP);CHKERRQ(ierr);
531 
532   ierr = TSSSPSetType(ts,TSSSPRKS2);CHKERRQ(ierr);
533   ssp->nstages = 5;
534   PetscFunctionReturn(0);
535 }
536 
537 #undef __FUNCT__
538 #define __FUNCT__ "TSSSPInitializePackage"
539 /*@C
540   TSSSPInitializePackage - This function initializes everything in the TSSSP package. It is called
541   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_SSP()
542   when using static libraries.
543 
544   Level: developer
545 
546 .keywords: TS, TSSSP, initialize, package
547 .seealso: PetscInitialize()
548 @*/
549 PetscErrorCode TSSSPInitializePackage(void)
550 {
551   PetscErrorCode ierr;
552 
553   PetscFunctionBegin;
554   if (TSSSPPackageInitialized) PetscFunctionReturn(0);
555   TSSSPPackageInitialized = PETSC_TRUE;
556   ierr = PetscFunctionListAdd(&TSSSPList,TSSSPRKS2, TSSSPStep_RK_2);CHKERRQ(ierr);
557   ierr = PetscFunctionListAdd(&TSSSPList,TSSSPRKS3, TSSSPStep_RK_3);CHKERRQ(ierr);
558   ierr = PetscFunctionListAdd(&TSSSPList,TSSSPRK104,TSSSPStep_RK_10_4);CHKERRQ(ierr);
559   ierr = PetscRegisterFinalize(TSSSPFinalizePackage);CHKERRQ(ierr);
560   PetscFunctionReturn(0);
561 }
562 
563 #undef __FUNCT__
564 #define __FUNCT__ "TSSSPFinalizePackage"
565 /*@C
566   TSSSPFinalizePackage - This function destroys everything in the TSSSP package. It is
567   called from PetscFinalize().
568 
569   Level: developer
570 
571 .keywords: Petsc, destroy, package
572 .seealso: PetscFinalize()
573 @*/
574 PetscErrorCode TSSSPFinalizePackage(void)
575 {
576   PetscErrorCode ierr;
577 
578   PetscFunctionBegin;
579   TSSSPPackageInitialized = PETSC_FALSE;
580   ierr = PetscFunctionListDestroy(&TSSSPList);CHKERRQ(ierr);
581   PetscFunctionReturn(0);
582 }
583