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