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