xref: /petsc/src/ts/impls/explicit/ssp/ssp.c (revision 4b91b6eae087708c6d6dc4e372cbc7b84f10f223)
1ef7bb5aaSJed Brown 
2ef7bb5aaSJed Brown /*
3ef7bb5aaSJed Brown        Code for Timestepping with explicit SSP.
4ef7bb5aaSJed Brown */
5c6db04a5SJed Brown #include <private/tsimpl.h>                /*I   "petscts.h"   I*/
6ef7bb5aaSJed Brown 
7ef7bb5aaSJed Brown PetscFList TSSSPList = 0;
8ef7bb5aaSJed Brown #define TSSSPType char*
9ef7bb5aaSJed Brown 
10ef7bb5aaSJed Brown #define TSSSPRKS2  "rks2"
11ef7bb5aaSJed Brown #define TSSSPRKS3  "rks3"
12ef7bb5aaSJed Brown #define TSSSPRK104 "rk104"
13ef7bb5aaSJed Brown 
14ef7bb5aaSJed Brown typedef struct {
15ef7bb5aaSJed Brown   PetscErrorCode (*onestep)(TS,PetscReal,PetscReal,Vec);
16ef7bb5aaSJed Brown   PetscInt nstages;
17ef7bb5aaSJed Brown   Vec xdot;
18ef7bb5aaSJed Brown   Vec *work;
19ef7bb5aaSJed Brown   PetscInt nwork;
20ace3abfcSBarry Smith   PetscBool  workout;
21ef7bb5aaSJed Brown } TS_SSP;
22ef7bb5aaSJed Brown 
23ef7bb5aaSJed Brown 
24ef7bb5aaSJed Brown #undef __FUNCT__
25ef7bb5aaSJed Brown #define __FUNCT__ "SSPGetWorkVectors"
26ef7bb5aaSJed Brown static PetscErrorCode SSPGetWorkVectors(TS ts,PetscInt n,Vec **work)
27ef7bb5aaSJed Brown {
28ef7bb5aaSJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
29ef7bb5aaSJed Brown   PetscErrorCode ierr;
30ef7bb5aaSJed Brown 
31ef7bb5aaSJed Brown   PetscFunctionBegin;
32e32f2f54SBarry Smith   if (ssp->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Work vectors already gotten");
33ef7bb5aaSJed Brown   if (ssp->nwork < n) {
34ef7bb5aaSJed Brown     if (ssp->nwork > 0) {
35574deadeSBarry Smith       ierr = VecDestroyVecs(ssp->nwork,&ssp->work);CHKERRQ(ierr);
36ef7bb5aaSJed Brown     }
37ef7bb5aaSJed Brown     ierr = VecDuplicateVecs(ts->vec_sol,n,&ssp->work);CHKERRQ(ierr);
38ef7bb5aaSJed Brown     ssp->nwork = n;
39ef7bb5aaSJed Brown   }
40ef7bb5aaSJed Brown   *work = ssp->work;
41ef7bb5aaSJed Brown   ssp->workout = PETSC_TRUE;
42ef7bb5aaSJed Brown   PetscFunctionReturn(0);
43ef7bb5aaSJed Brown }
44ef7bb5aaSJed Brown 
45ef7bb5aaSJed Brown #undef __FUNCT__
46ef7bb5aaSJed Brown #define __FUNCT__ "SSPRestoreWorkVectors"
47ef7bb5aaSJed Brown static PetscErrorCode SSPRestoreWorkVectors(TS ts,PetscInt n,Vec **work)
48ef7bb5aaSJed Brown {
49ef7bb5aaSJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
50ef7bb5aaSJed Brown 
51ef7bb5aaSJed Brown   PetscFunctionBegin;
52e32f2f54SBarry Smith   if (!ssp->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Work vectors have not been gotten");
53e32f2f54SBarry Smith   if (*work != ssp->work) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong work vectors checked out");
54ef7bb5aaSJed Brown   ssp->workout = PETSC_FALSE;
55ef7bb5aaSJed Brown   *work = PETSC_NULL;
56ef7bb5aaSJed Brown   PetscFunctionReturn(0);
57ef7bb5aaSJed Brown }
58ef7bb5aaSJed Brown 
59ef7bb5aaSJed Brown 
60ef7bb5aaSJed Brown #undef __FUNCT__
61ef7bb5aaSJed Brown #define __FUNCT__ "SSPStep_RK_2"
62ef7bb5aaSJed Brown /* Optimal second order SSP Runge-Kutta, low-storage, c_eff=(s-1)/s */
63ef7bb5aaSJed Brown /* Pseudocode 2 of Ketcheson 2008 */
64ef7bb5aaSJed Brown static PetscErrorCode SSPStep_RK_2(TS ts,PetscReal t0,PetscReal dt,Vec sol)
65ef7bb5aaSJed Brown {
66ef7bb5aaSJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
67ef7bb5aaSJed Brown   Vec *work,F;
68ef7bb5aaSJed Brown   PetscInt i,s;
69ef7bb5aaSJed Brown   PetscErrorCode ierr;
70ef7bb5aaSJed Brown 
71ef7bb5aaSJed Brown   PetscFunctionBegin;
72ef7bb5aaSJed Brown   s = ssp->nstages;
73ef7bb5aaSJed Brown   ierr = SSPGetWorkVectors(ts,2,&work);CHKERRQ(ierr);
74ef7bb5aaSJed Brown   F = work[1];
75ef7bb5aaSJed Brown   ierr = VecCopy(sol,work[0]);CHKERRQ(ierr);
76ef7bb5aaSJed Brown   for (i=0; i<s-1; i++) {
77ef7bb5aaSJed Brown     ierr = TSComputeRHSFunction(ts,t0+dt*(i/(s-1.)),work[0],F);CHKERRQ(ierr);
78ef7bb5aaSJed Brown     ierr = VecAXPY(work[0],dt/(s-1.),F);CHKERRQ(ierr);
79ef7bb5aaSJed Brown   }
80ef7bb5aaSJed Brown   ierr = TSComputeRHSFunction(ts,t0+dt,work[0],F);CHKERRQ(ierr);
81ef7bb5aaSJed Brown   ierr = VecAXPBYPCZ(sol,(s-1.)/s,dt/s,1./s,work[0],F);CHKERRQ(ierr);
82ef7bb5aaSJed Brown   ierr = SSPRestoreWorkVectors(ts,2,&work);CHKERRQ(ierr);
83ef7bb5aaSJed Brown   PetscFunctionReturn(0);
84ef7bb5aaSJed Brown }
85ef7bb5aaSJed Brown 
86ef7bb5aaSJed Brown #undef __FUNCT__
87ef7bb5aaSJed Brown #define __FUNCT__ "SSPStep_RK_3"
88ef7bb5aaSJed Brown /* Optimal third order SSP Runge-Kutta, low-storage, c_eff=(sqrt(s)-1)/sqrt(s), where sqrt(s) is an integer */
89ef7bb5aaSJed Brown /* Pseudocode 2 of Ketcheson 2008 */
90ef7bb5aaSJed Brown static PetscErrorCode SSPStep_RK_3(TS ts,PetscReal t0,PetscReal dt,Vec sol)
91ef7bb5aaSJed Brown {
92ef7bb5aaSJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
93ef7bb5aaSJed Brown   Vec *work,F;
94ef7bb5aaSJed Brown   PetscInt i,s,n,r;
95ef7bb5aaSJed Brown   PetscReal c;
96ef7bb5aaSJed Brown   PetscErrorCode ierr;
97ef7bb5aaSJed Brown 
98ef7bb5aaSJed Brown   PetscFunctionBegin;
99ef7bb5aaSJed Brown   s = ssp->nstages;
100fad8df86SSatish Balay   n = (PetscInt)(sqrt((PetscReal)s)+0.001);
101ef7bb5aaSJed Brown   r = s-n;
102e32f2f54SBarry Smith   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);
103ef7bb5aaSJed Brown   ierr = SSPGetWorkVectors(ts,3,&work);CHKERRQ(ierr);
104ef7bb5aaSJed Brown   F = work[2];
105ef7bb5aaSJed Brown   ierr = VecCopy(sol,work[0]);CHKERRQ(ierr);
106ef7bb5aaSJed Brown   for (i=0; i<(n-1)*(n-2)/2; i++) {
107ef7bb5aaSJed Brown     c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
108ef7bb5aaSJed Brown     ierr = TSComputeRHSFunction(ts,t0+c*dt,work[0],F);CHKERRQ(ierr);
109ef7bb5aaSJed Brown     ierr = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr);
110ef7bb5aaSJed Brown   }
111ef7bb5aaSJed Brown   ierr = VecCopy(work[0],work[1]);CHKERRQ(ierr);
112ef7bb5aaSJed Brown   for ( ; i<n*(n+1)/2-1; i++) {
113ef7bb5aaSJed Brown     c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
114ef7bb5aaSJed Brown     ierr = TSComputeRHSFunction(ts,t0+c*dt,work[0],F);CHKERRQ(ierr);
115ef7bb5aaSJed Brown     ierr = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr);
116ef7bb5aaSJed Brown   }
117ef7bb5aaSJed Brown   {
118ef7bb5aaSJed Brown     c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
119ef7bb5aaSJed Brown     ierr = TSComputeRHSFunction(ts,t0+c*dt,work[0],F);CHKERRQ(ierr);
120ef7bb5aaSJed 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);
121ef7bb5aaSJed Brown     i++;
122ef7bb5aaSJed Brown   }
123ef7bb5aaSJed Brown   for ( ; i<s; i++) {
124ef7bb5aaSJed Brown     c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
125ef7bb5aaSJed Brown     ierr = TSComputeRHSFunction(ts,t0+c*dt,work[0],F);CHKERRQ(ierr);
126ef7bb5aaSJed Brown     ierr = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr);
127ef7bb5aaSJed Brown   }
128ef7bb5aaSJed Brown   ierr = VecCopy(work[0],sol);CHKERRQ(ierr);
129ef7bb5aaSJed Brown   ierr = SSPRestoreWorkVectors(ts,3,&work);CHKERRQ(ierr);
130ef7bb5aaSJed Brown   PetscFunctionReturn(0);
131ef7bb5aaSJed Brown }
132ef7bb5aaSJed Brown 
133ef7bb5aaSJed Brown #undef __FUNCT__
134ef7bb5aaSJed Brown #define __FUNCT__ "SSPStep_RK_10_4"
135ef7bb5aaSJed Brown /* Optimal fourth order SSP Runge-Kutta, low-storage (2N), c_eff=0.6 */
136ef7bb5aaSJed Brown /* SSPRK(10,4), Pseudocode 3 of Ketcheson 2008 */
137ef7bb5aaSJed Brown static PetscErrorCode SSPStep_RK_10_4(TS ts,PetscReal t0,PetscReal dt,Vec sol)
138ef7bb5aaSJed Brown {
139ef7bb5aaSJed Brown   const PetscReal c[10] = {0, 1./6, 2./6, 3./6, 4./6, 2./6, 3./6, 4./6, 5./6, 1};
140ef7bb5aaSJed Brown   Vec *work,F;
1415a586d82SBarry Smith   PetscInt i;
142ef7bb5aaSJed Brown   PetscErrorCode ierr;
143ef7bb5aaSJed Brown 
144ef7bb5aaSJed Brown   PetscFunctionBegin;
145ef7bb5aaSJed Brown   ierr = SSPGetWorkVectors(ts,3,&work);CHKERRQ(ierr);
146ef7bb5aaSJed Brown   F = work[2];
147ef7bb5aaSJed Brown   ierr = VecCopy(sol,work[0]);CHKERRQ(ierr);
148ef7bb5aaSJed Brown   for (i=0; i<5; i++) {
149ef7bb5aaSJed Brown     ierr = TSComputeRHSFunction(ts,t0+c[i],work[0],F);CHKERRQ(ierr);
150ef7bb5aaSJed Brown     ierr = VecAXPY(work[0],dt/6,F);CHKERRQ(ierr);
151ef7bb5aaSJed Brown   }
152ef7bb5aaSJed Brown   ierr = VecAXPBYPCZ(work[1],1./25,9./25,0,sol,work[0]);CHKERRQ(ierr);
153ef7bb5aaSJed Brown   ierr = VecAXPBY(work[0],15,-5,work[1]);CHKERRQ(ierr);
154ef7bb5aaSJed Brown   for ( ; i<9; i++) {
155ef7bb5aaSJed Brown     ierr = TSComputeRHSFunction(ts,t0+c[i],work[0],F);CHKERRQ(ierr);
156ef7bb5aaSJed Brown     ierr = VecAXPY(work[0],dt/6,F);CHKERRQ(ierr);
157ef7bb5aaSJed Brown   }
158ef7bb5aaSJed Brown   ierr = TSComputeRHSFunction(ts,t0+dt,work[0],F);CHKERRQ(ierr);
159ef7bb5aaSJed Brown   ierr = VecAXPBYPCZ(work[1],3./5,dt/10,1,work[0],F);CHKERRQ(ierr);
160ef7bb5aaSJed Brown   ierr = VecCopy(work[1],sol);CHKERRQ(ierr);
161ef7bb5aaSJed Brown   ierr = SSPRestoreWorkVectors(ts,3,&work);CHKERRQ(ierr);
162ef7bb5aaSJed Brown   PetscFunctionReturn(0);
163ef7bb5aaSJed Brown }
164ef7bb5aaSJed Brown 
165ef7bb5aaSJed Brown 
166ef7bb5aaSJed Brown #undef __FUNCT__
167ef7bb5aaSJed Brown #define __FUNCT__ "TSSetUp_SSP"
168ef7bb5aaSJed Brown static PetscErrorCode TSSetUp_SSP(TS ts)
169ef7bb5aaSJed Brown {
170ef7bb5aaSJed Brown   /* TS_SSP       *ssp = (TS_SSP*)ts->data; */
171ef7bb5aaSJed Brown   /* PetscErrorCode ierr; */
172ef7bb5aaSJed Brown 
173ef7bb5aaSJed Brown   PetscFunctionBegin;
174ef7bb5aaSJed Brown   PetscFunctionReturn(0);
175ef7bb5aaSJed Brown }
176ef7bb5aaSJed Brown 
177ef7bb5aaSJed Brown #undef __FUNCT__
178ef7bb5aaSJed Brown #define __FUNCT__ "TSStep_SSP"
179ef7bb5aaSJed Brown static PetscErrorCode TSStep_SSP(TS ts,PetscInt *steps,PetscReal *ptime)
180ef7bb5aaSJed Brown {
181ef7bb5aaSJed Brown   TS_SSP        *ssp = (TS_SSP*)ts->data;
182ef7bb5aaSJed Brown   Vec            sol = ts->vec_sol;
183ef7bb5aaSJed Brown   PetscErrorCode ierr;
184ef7bb5aaSJed Brown   PetscInt       i,max_steps = ts->max_steps;
185ef7bb5aaSJed Brown 
186ef7bb5aaSJed Brown   PetscFunctionBegin;
187ef7bb5aaSJed Brown   *steps = -ts->steps;
188ef7bb5aaSJed Brown   ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
189ef7bb5aaSJed Brown 
190ef7bb5aaSJed Brown   for (i=0; i<max_steps; i++) {
191ef7bb5aaSJed Brown     PetscReal dt = ts->time_step;
192ef7bb5aaSJed Brown 
1933f2090d5SJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr);
194ef7bb5aaSJed Brown     ts->ptime += dt;
195ef7bb5aaSJed Brown     ierr = (*ssp->onestep)(ts,ts->ptime-dt,dt,sol);CHKERRQ(ierr);
196ef7bb5aaSJed Brown     ts->steps++;
1973f2090d5SJed Brown     ierr = TSPostStep(ts);CHKERRQ(ierr);
198ef7bb5aaSJed Brown     ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
199ef7bb5aaSJed Brown     if (ts->ptime > ts->max_time) break;
200ef7bb5aaSJed Brown   }
201ef7bb5aaSJed Brown 
202ef7bb5aaSJed Brown   *steps += ts->steps;
203ef7bb5aaSJed Brown   *ptime  = ts->ptime;
204ef7bb5aaSJed Brown   PetscFunctionReturn(0);
205ef7bb5aaSJed Brown }
206ef7bb5aaSJed Brown /*------------------------------------------------------------*/
207ef7bb5aaSJed Brown #undef __FUNCT__
208ef7bb5aaSJed Brown #define __FUNCT__ "TSDestroy_SSP"
209ef7bb5aaSJed Brown static PetscErrorCode TSDestroy_SSP(TS ts)
210ef7bb5aaSJed Brown {
211ef7bb5aaSJed Brown   TS_SSP         *ssp = (TS_SSP*)ts->data;
212ef7bb5aaSJed Brown   PetscErrorCode ierr;
213ef7bb5aaSJed Brown 
214ef7bb5aaSJed Brown   PetscFunctionBegin;
215574deadeSBarry Smith   if (ssp->work) {ierr = VecDestroyVecs(ssp->nwork,&ssp->work);CHKERRQ(ierr);}
216c31cb41cSBarry Smith   ierr = PetscFree(ts->data);CHKERRQ(ierr);
217ef7bb5aaSJed Brown   PetscFunctionReturn(0);
218ef7bb5aaSJed Brown }
219ef7bb5aaSJed Brown /*------------------------------------------------------------*/
220ef7bb5aaSJed Brown 
221ef7bb5aaSJed Brown #undef __FUNCT__
222ef7bb5aaSJed Brown #define __FUNCT__ "TSSSPSetType"
223ef7bb5aaSJed Brown static PetscErrorCode TSSSPSetType(TS ts,const TSSSPType type)
224ef7bb5aaSJed Brown {
225ef7bb5aaSJed Brown   PetscErrorCode ierr,(*r)(TS,PetscReal,PetscReal,Vec);
226ef7bb5aaSJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
227ef7bb5aaSJed Brown 
228ef7bb5aaSJed Brown   PetscFunctionBegin;
229*4b91b6eaSBarry Smith   ierr = PetscFListFind(TSSSPList,((PetscObject)ts)->comm,type,PETSC_TRUE,(void(**)(void))&r);CHKERRQ(ierr);
230e32f2f54SBarry Smith   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TS_SSP type %s given",type);
231ef7bb5aaSJed Brown   ssp->onestep = r;
232ef7bb5aaSJed Brown   PetscFunctionReturn(0);
233ef7bb5aaSJed Brown }
234ef7bb5aaSJed Brown 
235ef7bb5aaSJed Brown #undef __FUNCT__
236ef7bb5aaSJed Brown #define __FUNCT__ "TSSetFromOptions_SSP"
237ef7bb5aaSJed Brown static PetscErrorCode TSSetFromOptions_SSP(TS ts)
238ef7bb5aaSJed Brown {
239ef7bb5aaSJed Brown   char tname[256] = TSSSPRKS2;
240ef7bb5aaSJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
241ef7bb5aaSJed Brown   PetscErrorCode ierr;
242ace3abfcSBarry Smith   PetscBool  flg;
243ef7bb5aaSJed Brown 
244ef7bb5aaSJed Brown   PetscFunctionBegin;
245ef7bb5aaSJed Brown   ierr = PetscOptionsHead("SSP ODE solver options");CHKERRQ(ierr);
246ef7bb5aaSJed Brown   {
247ef7bb5aaSJed Brown     ierr = PetscOptionsList("-ts_ssp_type","Type of SSP method","TSSSPSetType",TSSSPList,tname,tname,sizeof(tname),&flg);CHKERRQ(ierr);
248ef7bb5aaSJed Brown     if (flg) {
249ef7bb5aaSJed Brown       ierr = TSSSPSetType(ts,tname);CHKERRQ(ierr);
250ef7bb5aaSJed Brown     }
251ef7bb5aaSJed Brown     ierr = PetscOptionsInt("-ts_ssp_nstages","Number of stages","TSSSPSetNumStages",ssp->nstages,&ssp->nstages,PETSC_NULL);CHKERRQ(ierr);
252ef7bb5aaSJed Brown   }
253ef7bb5aaSJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
254ef7bb5aaSJed Brown   PetscFunctionReturn(0);
255ef7bb5aaSJed Brown }
256ef7bb5aaSJed Brown 
257ef7bb5aaSJed Brown #undef __FUNCT__
258ef7bb5aaSJed Brown #define __FUNCT__ "TSView_SSP"
259ef7bb5aaSJed Brown static PetscErrorCode TSView_SSP(TS ts,PetscViewer viewer)
260ef7bb5aaSJed Brown {
261ef7bb5aaSJed Brown   PetscFunctionBegin;
262ef7bb5aaSJed Brown   PetscFunctionReturn(0);
263ef7bb5aaSJed Brown }
264ef7bb5aaSJed Brown 
265ef7bb5aaSJed Brown /* ------------------------------------------------------------ */
266ef7bb5aaSJed Brown 
267ef7bb5aaSJed Brown /*MC
2688ab3e0fcSJed Brown       TSSSP - Explicit strong stability preserving ODE solver
2698ab3e0fcSJed Brown 
2708ab3e0fcSJed Brown   Most hyperbolic conservation laws have exact solutions that are total variation diminishing (TVD) or total variation
2718ab3e0fcSJed Brown   bounded (TVB) although these solutions often contain discontinuities.  Spatial discretizations such as Godunov's
2728ab3e0fcSJed Brown   scheme and high-resolution finite volume methods (TVD limiters, ENO/WENO) are designed to preserve these properties,
2738ab3e0fcSJed Brown   but they are usually formulated using a forward Euler time discretization or by coupling the space and time
2748ab3e0fcSJed Brown   discretization as in the classical Lax-Wendroff scheme.  When the space and time discretization is coupled, it is very
2758ab3e0fcSJed Brown   difficult to produce schemes with high temporal accuracy while preserving TVD properties.  An alternative is the
2768ab3e0fcSJed Brown   semidiscrete formulation where we choose a spatial discretization that is TVD with forward Euler and then choose a
2778ab3e0fcSJed Brown   time discretization that preserves the TVD property.  Such integrators are called strong stability preserving (SSP).
2788ab3e0fcSJed Brown 
2798ab3e0fcSJed Brown   Let c_eff be the minimum number of function evaluations required to step as far as one step of forward Euler while
2808ab3e0fcSJed Brown   still being SSP.  Some theoretical bounds
2818ab3e0fcSJed Brown 
2828ab3e0fcSJed Brown   1. There are no explicit methods with c_eff > 1.
2830085e20eSJed Brown 
2848ab3e0fcSJed Brown   2. There are no explicit methods beyond order 4 (for nonlinear problems) and c_eff > 0.
2850085e20eSJed Brown 
2868ab3e0fcSJed Brown   3. There are no implicit methods with order greater than 1 and c_eff > 2.
2878ab3e0fcSJed Brown 
2888ab3e0fcSJed Brown   This integrator provides Runge-Kutta methods of order 2, 3, and 4 with maximal values of c_eff.  More stages allows
2898ab3e0fcSJed Brown   for larger values of c_eff which improves efficiency.  These implementations are low-memory and only use 2 or 3 work
2908ab3e0fcSJed Brown   vectors regardless of the total number of stages, so e.g. 25-stage 3rd order methods may be an excellent choice.
2918ab3e0fcSJed Brown 
2928ab3e0fcSJed Brown   Methods can be chosen with -ts_ssp_type {rks2,rks3,rk104}
2938ab3e0fcSJed Brown 
2948ab3e0fcSJed Brown   rks2: Second order methods with any number s>1 of stages.  c_eff = (s-1)/s
2958ab3e0fcSJed Brown 
2968ab3e0fcSJed Brown   rks3: Third order methods with s=n^2 stages, n>1.  c_eff = (s-n)/s
2978ab3e0fcSJed Brown 
2988ab3e0fcSJed Brown   rk104: A 10-stage fourth order method.  c_eff = 0.6
299ef7bb5aaSJed Brown 
300ef7bb5aaSJed Brown   Level: beginner
301ef7bb5aaSJed Brown 
302ef7bb5aaSJed Brown .seealso:  TSCreate(), TS, TSSetType()
303ef7bb5aaSJed Brown 
304ef7bb5aaSJed Brown M*/
305ef7bb5aaSJed Brown EXTERN_C_BEGIN
306ef7bb5aaSJed Brown #undef __FUNCT__
307ef7bb5aaSJed Brown #define __FUNCT__ "TSCreate_SSP"
3087087cfbeSBarry Smith PetscErrorCode  TSCreate_SSP(TS ts)
309ef7bb5aaSJed Brown {
310ef7bb5aaSJed Brown   TS_SSP       *ssp;
311ef7bb5aaSJed Brown   PetscErrorCode ierr;
312ef7bb5aaSJed Brown 
313ef7bb5aaSJed Brown   PetscFunctionBegin;
314ef7bb5aaSJed Brown   if (!TSSSPList) {
315ef7bb5aaSJed Brown     ierr = PetscFListAdd(&TSSSPList,TSSSPRKS2,  "SSPStep_RK_2",   (void(*)(void))SSPStep_RK_2);CHKERRQ(ierr);
316ef7bb5aaSJed Brown     ierr = PetscFListAdd(&TSSSPList,TSSSPRKS3,  "SSPStep_RK_3",   (void(*)(void))SSPStep_RK_3);CHKERRQ(ierr);
317ef7bb5aaSJed Brown     ierr = PetscFListAdd(&TSSSPList,TSSSPRK104, "SSPStep_RK_10_4",(void(*)(void))SSPStep_RK_10_4);CHKERRQ(ierr);
318ef7bb5aaSJed Brown   }
319ef7bb5aaSJed Brown 
320ef7bb5aaSJed Brown   ts->ops->setup           = TSSetUp_SSP;
321ef7bb5aaSJed Brown   ts->ops->step            = TSStep_SSP;
322ef7bb5aaSJed Brown   ts->ops->destroy         = TSDestroy_SSP;
323ef7bb5aaSJed Brown   ts->ops->setfromoptions  = TSSetFromOptions_SSP;
324ef7bb5aaSJed Brown   ts->ops->view            = TSView_SSP;
325ef7bb5aaSJed Brown 
326ef7bb5aaSJed Brown   ierr = PetscNewLog(ts,TS_SSP,&ssp);CHKERRQ(ierr);
327ef7bb5aaSJed Brown   ts->data = (void*)ssp;
328ef7bb5aaSJed Brown 
329ef7bb5aaSJed Brown   ierr = TSSSPSetType(ts,TSSSPRKS2);CHKERRQ(ierr);
330ef7bb5aaSJed Brown   ssp->nstages = 5;
331ef7bb5aaSJed Brown   PetscFunctionReturn(0);
332ef7bb5aaSJed Brown }
333ef7bb5aaSJed Brown EXTERN_C_END
334ef7bb5aaSJed Brown 
335ef7bb5aaSJed Brown 
336ef7bb5aaSJed Brown 
337ef7bb5aaSJed Brown 
338