xref: /petsc/src/ts/impls/explicit/ssp/ssp.c (revision ef7bb5aa64da98453fa5fb4571451d73f643854f)
1*ef7bb5aaSJed Brown #define PETSCTS_DLL
2*ef7bb5aaSJed Brown 
3*ef7bb5aaSJed Brown /*
4*ef7bb5aaSJed Brown        Code for Timestepping with explicit SSP.
5*ef7bb5aaSJed Brown */
6*ef7bb5aaSJed Brown #include "private/tsimpl.h"                /*I   "petscts.h"   I*/
7*ef7bb5aaSJed Brown 
8*ef7bb5aaSJed Brown PetscFList TSSSPList = 0;
9*ef7bb5aaSJed Brown #define TSSSPType char*
10*ef7bb5aaSJed Brown 
11*ef7bb5aaSJed Brown #define TSSSPRKS2  "rks2"
12*ef7bb5aaSJed Brown #define TSSSPRKS3  "rks3"
13*ef7bb5aaSJed Brown #define TSSSPRK104 "rk104"
14*ef7bb5aaSJed Brown 
15*ef7bb5aaSJed Brown typedef struct {
16*ef7bb5aaSJed Brown   PetscErrorCode (*onestep)(TS,PetscReal,PetscReal,Vec);
17*ef7bb5aaSJed Brown   PetscInt nstages;
18*ef7bb5aaSJed Brown   Vec xdot;
19*ef7bb5aaSJed Brown   Vec *work;
20*ef7bb5aaSJed Brown   PetscInt nwork;
21*ef7bb5aaSJed Brown   PetscTruth workout;
22*ef7bb5aaSJed Brown } TS_SSP;
23*ef7bb5aaSJed Brown 
24*ef7bb5aaSJed Brown 
25*ef7bb5aaSJed Brown #undef __FUNCT__
26*ef7bb5aaSJed Brown #define __FUNCT__ "SSPGetWorkVectors"
27*ef7bb5aaSJed Brown static PetscErrorCode SSPGetWorkVectors(TS ts,PetscInt n,Vec **work)
28*ef7bb5aaSJed Brown {
29*ef7bb5aaSJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
30*ef7bb5aaSJed Brown   PetscErrorCode ierr;
31*ef7bb5aaSJed Brown 
32*ef7bb5aaSJed Brown   PetscFunctionBegin;
33*ef7bb5aaSJed Brown   if (ssp->workout) SETERRQ(PETSC_ERR_PLIB,"Work vectors already gotten");
34*ef7bb5aaSJed Brown   if (ssp->nwork < n) {
35*ef7bb5aaSJed Brown     if (ssp->nwork > 0) {
36*ef7bb5aaSJed Brown       ierr = VecDestroyVecs(ssp->work,ssp->nwork);CHKERRQ(ierr);
37*ef7bb5aaSJed Brown     }
38*ef7bb5aaSJed Brown     ierr = VecDuplicateVecs(ts->vec_sol,n,&ssp->work);CHKERRQ(ierr);
39*ef7bb5aaSJed Brown     ssp->nwork = n;
40*ef7bb5aaSJed Brown   }
41*ef7bb5aaSJed Brown   *work = ssp->work;
42*ef7bb5aaSJed Brown   ssp->workout = PETSC_TRUE;
43*ef7bb5aaSJed Brown   PetscFunctionReturn(0);
44*ef7bb5aaSJed Brown }
45*ef7bb5aaSJed Brown 
46*ef7bb5aaSJed Brown #undef __FUNCT__
47*ef7bb5aaSJed Brown #define __FUNCT__ "SSPRestoreWorkVectors"
48*ef7bb5aaSJed Brown static PetscErrorCode SSPRestoreWorkVectors(TS ts,PetscInt n,Vec **work)
49*ef7bb5aaSJed Brown {
50*ef7bb5aaSJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
51*ef7bb5aaSJed Brown 
52*ef7bb5aaSJed Brown   PetscFunctionBegin;
53*ef7bb5aaSJed Brown   if (!ssp->workout) SETERRQ(PETSC_ERR_ORDER,"Work vectors have not been gotten");
54*ef7bb5aaSJed Brown   if (*work != ssp->work) SETERRQ(PETSC_ERR_PLIB,"Wrong work vectors checked out");
55*ef7bb5aaSJed Brown   ssp->workout = PETSC_FALSE;
56*ef7bb5aaSJed Brown   *work = PETSC_NULL;
57*ef7bb5aaSJed Brown   PetscFunctionReturn(0);
58*ef7bb5aaSJed Brown }
59*ef7bb5aaSJed Brown 
60*ef7bb5aaSJed Brown 
61*ef7bb5aaSJed Brown #undef __FUNCT__
62*ef7bb5aaSJed Brown #define __FUNCT__ "SSPStep_RK_2"
63*ef7bb5aaSJed Brown /* Optimal second order SSP Runge-Kutta, low-storage, c_eff=(s-1)/s */
64*ef7bb5aaSJed Brown /* Pseudocode 2 of Ketcheson 2008 */
65*ef7bb5aaSJed Brown static PetscErrorCode SSPStep_RK_2(TS ts,PetscReal t0,PetscReal dt,Vec sol)
66*ef7bb5aaSJed Brown {
67*ef7bb5aaSJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
68*ef7bb5aaSJed Brown   Vec *work,F;
69*ef7bb5aaSJed Brown   PetscInt i,s;
70*ef7bb5aaSJed Brown   PetscErrorCode ierr;
71*ef7bb5aaSJed Brown 
72*ef7bb5aaSJed Brown   PetscFunctionBegin;
73*ef7bb5aaSJed Brown   s = ssp->nstages;
74*ef7bb5aaSJed Brown   ierr = SSPGetWorkVectors(ts,2,&work);CHKERRQ(ierr);
75*ef7bb5aaSJed Brown   F = work[1];
76*ef7bb5aaSJed Brown   ierr = VecCopy(sol,work[0]);CHKERRQ(ierr);
77*ef7bb5aaSJed Brown   for (i=0; i<s-1; i++) {
78*ef7bb5aaSJed Brown     ierr = TSComputeRHSFunction(ts,t0+dt*(i/(s-1.)),work[0],F);CHKERRQ(ierr);
79*ef7bb5aaSJed Brown     ierr = VecAXPY(work[0],dt/(s-1.),F);CHKERRQ(ierr);
80*ef7bb5aaSJed Brown   }
81*ef7bb5aaSJed Brown   ierr = TSComputeRHSFunction(ts,t0+dt,work[0],F);CHKERRQ(ierr);
82*ef7bb5aaSJed Brown   ierr = VecAXPBYPCZ(sol,(s-1.)/s,dt/s,1./s,work[0],F);CHKERRQ(ierr);
83*ef7bb5aaSJed Brown   ierr = SSPRestoreWorkVectors(ts,2,&work);CHKERRQ(ierr);
84*ef7bb5aaSJed Brown   PetscFunctionReturn(0);
85*ef7bb5aaSJed Brown }
86*ef7bb5aaSJed Brown 
87*ef7bb5aaSJed Brown #undef __FUNCT__
88*ef7bb5aaSJed Brown #define __FUNCT__ "SSPStep_RK_3"
89*ef7bb5aaSJed Brown /* Optimal third order SSP Runge-Kutta, low-storage, c_eff=(sqrt(s)-1)/sqrt(s), where sqrt(s) is an integer */
90*ef7bb5aaSJed Brown /* Pseudocode 2 of Ketcheson 2008 */
91*ef7bb5aaSJed Brown static PetscErrorCode SSPStep_RK_3(TS ts,PetscReal t0,PetscReal dt,Vec sol)
92*ef7bb5aaSJed Brown {
93*ef7bb5aaSJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
94*ef7bb5aaSJed Brown   Vec *work,F;
95*ef7bb5aaSJed Brown   PetscInt i,s,n,r;
96*ef7bb5aaSJed Brown   PetscReal c;
97*ef7bb5aaSJed Brown   PetscErrorCode ierr;
98*ef7bb5aaSJed Brown 
99*ef7bb5aaSJed Brown   PetscFunctionBegin;
100*ef7bb5aaSJed Brown   s = ssp->nstages;
101*ef7bb5aaSJed Brown   n = floor(sqrt(s));
102*ef7bb5aaSJed Brown   r = s-n;
103*ef7bb5aaSJed Brown   if (n*n != s) SETERRQ1(PETSC_ERR_SUP,"No support for optimal third order schemes with %d stages, must be a square number at least 4",s);
104*ef7bb5aaSJed Brown   ierr = SSPGetWorkVectors(ts,3,&work);CHKERRQ(ierr);
105*ef7bb5aaSJed Brown   F = work[2];
106*ef7bb5aaSJed Brown   ierr = VecCopy(sol,work[0]);CHKERRQ(ierr);
107*ef7bb5aaSJed Brown   for (i=0; i<(n-1)*(n-2)/2; i++) {
108*ef7bb5aaSJed Brown     c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
109*ef7bb5aaSJed Brown     ierr = TSComputeRHSFunction(ts,t0+c*dt,work[0],F);CHKERRQ(ierr);
110*ef7bb5aaSJed Brown     ierr = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr);
111*ef7bb5aaSJed Brown   }
112*ef7bb5aaSJed Brown   ierr = VecCopy(work[0],work[1]);CHKERRQ(ierr);
113*ef7bb5aaSJed Brown   for ( ; i<n*(n+1)/2-1; i++) {
114*ef7bb5aaSJed Brown     c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
115*ef7bb5aaSJed Brown     ierr = TSComputeRHSFunction(ts,t0+c*dt,work[0],F);CHKERRQ(ierr);
116*ef7bb5aaSJed Brown     ierr = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr);
117*ef7bb5aaSJed Brown   }
118*ef7bb5aaSJed Brown   {
119*ef7bb5aaSJed Brown     c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
120*ef7bb5aaSJed Brown     ierr = TSComputeRHSFunction(ts,t0+c*dt,work[0],F);CHKERRQ(ierr);
121*ef7bb5aaSJed Brown     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);
122*ef7bb5aaSJed Brown     i++;
123*ef7bb5aaSJed Brown   }
124*ef7bb5aaSJed Brown   for ( ; i<s; i++) {
125*ef7bb5aaSJed Brown     c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
126*ef7bb5aaSJed Brown     ierr = TSComputeRHSFunction(ts,t0+c*dt,work[0],F);CHKERRQ(ierr);
127*ef7bb5aaSJed Brown     ierr = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr);
128*ef7bb5aaSJed Brown   }
129*ef7bb5aaSJed Brown   ierr = VecCopy(work[0],sol);CHKERRQ(ierr);
130*ef7bb5aaSJed Brown   ierr = SSPRestoreWorkVectors(ts,3,&work);CHKERRQ(ierr);
131*ef7bb5aaSJed Brown   PetscFunctionReturn(0);
132*ef7bb5aaSJed Brown }
133*ef7bb5aaSJed Brown 
134*ef7bb5aaSJed Brown #undef __FUNCT__
135*ef7bb5aaSJed Brown #define __FUNCT__ "SSPStep_RK_10_4"
136*ef7bb5aaSJed Brown /* Optimal fourth order SSP Runge-Kutta, low-storage (2N), c_eff=0.6 */
137*ef7bb5aaSJed Brown /* SSPRK(10,4), Pseudocode 3 of Ketcheson 2008 */
138*ef7bb5aaSJed Brown static PetscErrorCode SSPStep_RK_10_4(TS ts,PetscReal t0,PetscReal dt,Vec sol)
139*ef7bb5aaSJed Brown {
140*ef7bb5aaSJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
141*ef7bb5aaSJed Brown   const PetscReal c[10] = {0, 1./6, 2./6, 3./6, 4./6, 2./6, 3./6, 4./6, 5./6, 1};
142*ef7bb5aaSJed Brown   Vec *work,F;
143*ef7bb5aaSJed Brown   PetscInt i,s;
144*ef7bb5aaSJed Brown   PetscErrorCode ierr;
145*ef7bb5aaSJed Brown 
146*ef7bb5aaSJed Brown   PetscFunctionBegin;
147*ef7bb5aaSJed Brown   s = ssp->nstages;
148*ef7bb5aaSJed Brown   ierr = SSPGetWorkVectors(ts,3,&work);CHKERRQ(ierr);
149*ef7bb5aaSJed Brown   F = work[2];
150*ef7bb5aaSJed Brown   ierr = VecCopy(sol,work[0]);CHKERRQ(ierr);
151*ef7bb5aaSJed Brown   for (i=0; i<5; i++) {
152*ef7bb5aaSJed Brown     ierr = TSComputeRHSFunction(ts,t0+c[i],work[0],F);CHKERRQ(ierr);
153*ef7bb5aaSJed Brown     ierr = VecAXPY(work[0],dt/6,F);CHKERRQ(ierr);
154*ef7bb5aaSJed Brown   }
155*ef7bb5aaSJed Brown   ierr = VecAXPBYPCZ(work[1],1./25,9./25,0,sol,work[0]);CHKERRQ(ierr);
156*ef7bb5aaSJed Brown   ierr = VecAXPBY(work[0],15,-5,work[1]);CHKERRQ(ierr);
157*ef7bb5aaSJed Brown   for ( ; i<9; i++) {
158*ef7bb5aaSJed Brown     ierr = TSComputeRHSFunction(ts,t0+c[i],work[0],F);CHKERRQ(ierr);
159*ef7bb5aaSJed Brown     ierr = VecAXPY(work[0],dt/6,F);CHKERRQ(ierr);
160*ef7bb5aaSJed Brown   }
161*ef7bb5aaSJed Brown   ierr = TSComputeRHSFunction(ts,t0+dt,work[0],F);CHKERRQ(ierr);
162*ef7bb5aaSJed Brown   ierr = VecAXPBYPCZ(work[1],3./5,dt/10,1,work[0],F);CHKERRQ(ierr);
163*ef7bb5aaSJed Brown   ierr = VecCopy(work[1],sol);CHKERRQ(ierr);
164*ef7bb5aaSJed Brown   ierr = SSPRestoreWorkVectors(ts,3,&work);CHKERRQ(ierr);
165*ef7bb5aaSJed Brown   PetscFunctionReturn(0);
166*ef7bb5aaSJed Brown }
167*ef7bb5aaSJed Brown 
168*ef7bb5aaSJed Brown 
169*ef7bb5aaSJed Brown #undef __FUNCT__
170*ef7bb5aaSJed Brown #define __FUNCT__ "TSSetUp_SSP"
171*ef7bb5aaSJed Brown static PetscErrorCode TSSetUp_SSP(TS ts)
172*ef7bb5aaSJed Brown {
173*ef7bb5aaSJed Brown   /* TS_SSP       *ssp = (TS_SSP*)ts->data; */
174*ef7bb5aaSJed Brown   /* PetscErrorCode ierr; */
175*ef7bb5aaSJed Brown 
176*ef7bb5aaSJed Brown   PetscFunctionBegin;
177*ef7bb5aaSJed Brown   PetscFunctionReturn(0);
178*ef7bb5aaSJed Brown }
179*ef7bb5aaSJed Brown 
180*ef7bb5aaSJed Brown #undef __FUNCT__
181*ef7bb5aaSJed Brown #define __FUNCT__ "TSStep_SSP"
182*ef7bb5aaSJed Brown static PetscErrorCode TSStep_SSP(TS ts,PetscInt *steps,PetscReal *ptime)
183*ef7bb5aaSJed Brown {
184*ef7bb5aaSJed Brown   TS_SSP        *ssp = (TS_SSP*)ts->data;
185*ef7bb5aaSJed Brown   Vec            sol = ts->vec_sol;
186*ef7bb5aaSJed Brown   PetscErrorCode ierr;
187*ef7bb5aaSJed Brown   PetscInt       i,max_steps = ts->max_steps;
188*ef7bb5aaSJed Brown 
189*ef7bb5aaSJed Brown   PetscFunctionBegin;
190*ef7bb5aaSJed Brown   *steps = -ts->steps;
191*ef7bb5aaSJed Brown   ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
192*ef7bb5aaSJed Brown 
193*ef7bb5aaSJed Brown   for (i=0; i<max_steps; i++) {
194*ef7bb5aaSJed Brown     PetscReal dt = ts->time_step;
195*ef7bb5aaSJed Brown 
196*ef7bb5aaSJed Brown     ts->ptime += dt;
197*ef7bb5aaSJed Brown     ierr = (*ssp->onestep)(ts,ts->ptime-dt,dt,sol);CHKERRQ(ierr);
198*ef7bb5aaSJed Brown     ts->steps++;
199*ef7bb5aaSJed Brown     ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
200*ef7bb5aaSJed Brown     if (ts->ptime > ts->max_time) break;
201*ef7bb5aaSJed Brown   }
202*ef7bb5aaSJed Brown 
203*ef7bb5aaSJed Brown   *steps += ts->steps;
204*ef7bb5aaSJed Brown   *ptime  = ts->ptime;
205*ef7bb5aaSJed Brown   PetscFunctionReturn(0);
206*ef7bb5aaSJed Brown }
207*ef7bb5aaSJed Brown /*------------------------------------------------------------*/
208*ef7bb5aaSJed Brown #undef __FUNCT__
209*ef7bb5aaSJed Brown #define __FUNCT__ "TSDestroy_SSP"
210*ef7bb5aaSJed Brown static PetscErrorCode TSDestroy_SSP(TS ts)
211*ef7bb5aaSJed Brown {
212*ef7bb5aaSJed Brown   TS_SSP       *ssp = (TS_SSP*)ts->data;
213*ef7bb5aaSJed Brown   PetscErrorCode ierr;
214*ef7bb5aaSJed Brown 
215*ef7bb5aaSJed Brown   PetscFunctionBegin;
216*ef7bb5aaSJed Brown   if (ssp->work) {ierr = VecDestroyVecs(ssp->work,ssp->nwork);CHKERRQ(ierr);}
217*ef7bb5aaSJed Brown   ierr = PetscFree(ssp);CHKERRQ(ierr);
218*ef7bb5aaSJed Brown   PetscFunctionReturn(0);
219*ef7bb5aaSJed Brown }
220*ef7bb5aaSJed Brown /*------------------------------------------------------------*/
221*ef7bb5aaSJed Brown 
222*ef7bb5aaSJed Brown #undef __FUNCT__
223*ef7bb5aaSJed Brown #define __FUNCT__ "TSSSPSetType"
224*ef7bb5aaSJed Brown static PetscErrorCode TSSSPSetType(TS ts,const TSSSPType type)
225*ef7bb5aaSJed Brown {
226*ef7bb5aaSJed Brown   PetscErrorCode ierr,(*r)(TS,PetscReal,PetscReal,Vec);
227*ef7bb5aaSJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
228*ef7bb5aaSJed Brown 
229*ef7bb5aaSJed Brown   PetscFunctionBegin;
230*ef7bb5aaSJed Brown   ierr = PetscFListFind(TSSSPList,((PetscObject)ts)->comm,type,(void(**)(void))&r);CHKERRQ(ierr);
231*ef7bb5aaSJed Brown   if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TS_SSP type %s given",type);
232*ef7bb5aaSJed Brown   ssp->onestep = r;
233*ef7bb5aaSJed Brown   PetscFunctionReturn(0);
234*ef7bb5aaSJed Brown }
235*ef7bb5aaSJed Brown 
236*ef7bb5aaSJed Brown #undef __FUNCT__
237*ef7bb5aaSJed Brown #define __FUNCT__ "TSSetFromOptions_SSP"
238*ef7bb5aaSJed Brown static PetscErrorCode TSSetFromOptions_SSP(TS ts)
239*ef7bb5aaSJed Brown {
240*ef7bb5aaSJed Brown   char tname[256] = TSSSPRKS2;
241*ef7bb5aaSJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
242*ef7bb5aaSJed Brown   PetscErrorCode ierr;
243*ef7bb5aaSJed Brown   PetscTruth flg;
244*ef7bb5aaSJed Brown 
245*ef7bb5aaSJed Brown   PetscFunctionBegin;
246*ef7bb5aaSJed Brown   ierr = PetscOptionsHead("SSP ODE solver options");CHKERRQ(ierr);
247*ef7bb5aaSJed Brown   {
248*ef7bb5aaSJed Brown     ierr = PetscOptionsList("-ts_ssp_type","Type of SSP method","TSSSPSetType",TSSSPList,tname,tname,sizeof(tname),&flg);CHKERRQ(ierr);
249*ef7bb5aaSJed Brown     if (flg) {
250*ef7bb5aaSJed Brown       ierr = TSSSPSetType(ts,tname);CHKERRQ(ierr);
251*ef7bb5aaSJed Brown     }
252*ef7bb5aaSJed Brown     ierr = PetscOptionsInt("-ts_ssp_nstages","Number of stages","TSSSPSetNumStages",ssp->nstages,&ssp->nstages,PETSC_NULL);CHKERRQ(ierr);
253*ef7bb5aaSJed Brown   }
254*ef7bb5aaSJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
255*ef7bb5aaSJed Brown   PetscFunctionReturn(0);
256*ef7bb5aaSJed Brown }
257*ef7bb5aaSJed Brown 
258*ef7bb5aaSJed Brown #undef __FUNCT__
259*ef7bb5aaSJed Brown #define __FUNCT__ "TSView_SSP"
260*ef7bb5aaSJed Brown static PetscErrorCode TSView_SSP(TS ts,PetscViewer viewer)
261*ef7bb5aaSJed Brown {
262*ef7bb5aaSJed Brown   PetscFunctionBegin;
263*ef7bb5aaSJed Brown   PetscFunctionReturn(0);
264*ef7bb5aaSJed Brown }
265*ef7bb5aaSJed Brown 
266*ef7bb5aaSJed Brown /* ------------------------------------------------------------ */
267*ef7bb5aaSJed Brown 
268*ef7bb5aaSJed Brown /*MC
269*ef7bb5aaSJed Brown       TSSSP - ODE solver using the explicit SSP method
270*ef7bb5aaSJed Brown 
271*ef7bb5aaSJed Brown   Level: beginner
272*ef7bb5aaSJed Brown 
273*ef7bb5aaSJed Brown .seealso:  TSCreate(), TS, TSSetType()
274*ef7bb5aaSJed Brown 
275*ef7bb5aaSJed Brown M*/
276*ef7bb5aaSJed Brown EXTERN_C_BEGIN
277*ef7bb5aaSJed Brown #undef __FUNCT__
278*ef7bb5aaSJed Brown #define __FUNCT__ "TSCreate_SSP"
279*ef7bb5aaSJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSCreate_SSP(TS ts)
280*ef7bb5aaSJed Brown {
281*ef7bb5aaSJed Brown   TS_SSP       *ssp;
282*ef7bb5aaSJed Brown   PetscErrorCode ierr;
283*ef7bb5aaSJed Brown 
284*ef7bb5aaSJed Brown   PetscFunctionBegin;
285*ef7bb5aaSJed Brown   if (!TSSSPList) {
286*ef7bb5aaSJed Brown     ierr = PetscFListAdd(&TSSSPList,TSSSPRKS2,  "SSPStep_RK_2",   (void(*)(void))SSPStep_RK_2);CHKERRQ(ierr);
287*ef7bb5aaSJed Brown     ierr = PetscFListAdd(&TSSSPList,TSSSPRKS3,  "SSPStep_RK_3",   (void(*)(void))SSPStep_RK_3);CHKERRQ(ierr);
288*ef7bb5aaSJed Brown     ierr = PetscFListAdd(&TSSSPList,TSSSPRK104, "SSPStep_RK_10_4",(void(*)(void))SSPStep_RK_10_4);CHKERRQ(ierr);
289*ef7bb5aaSJed Brown   }
290*ef7bb5aaSJed Brown 
291*ef7bb5aaSJed Brown   ts->ops->setup           = TSSetUp_SSP;
292*ef7bb5aaSJed Brown   ts->ops->step            = TSStep_SSP;
293*ef7bb5aaSJed Brown   ts->ops->destroy         = TSDestroy_SSP;
294*ef7bb5aaSJed Brown   ts->ops->setfromoptions  = TSSetFromOptions_SSP;
295*ef7bb5aaSJed Brown   ts->ops->view            = TSView_SSP;
296*ef7bb5aaSJed Brown 
297*ef7bb5aaSJed Brown   ierr = PetscNewLog(ts,TS_SSP,&ssp);CHKERRQ(ierr);
298*ef7bb5aaSJed Brown   ts->data = (void*)ssp;
299*ef7bb5aaSJed Brown 
300*ef7bb5aaSJed Brown   ierr = TSSSPSetType(ts,TSSSPRKS2);CHKERRQ(ierr);
301*ef7bb5aaSJed Brown   ssp->nstages = 5;
302*ef7bb5aaSJed Brown   PetscFunctionReturn(0);
303*ef7bb5aaSJed Brown }
304*ef7bb5aaSJed Brown EXTERN_C_END
305*ef7bb5aaSJed Brown 
306*ef7bb5aaSJed Brown 
307*ef7bb5aaSJed Brown 
308*ef7bb5aaSJed Brown 
309