xref: /petsc/src/ts/impls/explicit/ssp/ssp.c (revision b330ce4d4dd0241bb385d3087e2abb900b4ce54a)
1ef7bb5aaSJed Brown /*
2ef7bb5aaSJed Brown        Code for Timestepping with explicit SSP.
3ef7bb5aaSJed Brown */
4c6db04a5SJed Brown #include <private/tsimpl.h>                /*I   "petscts.h"   I*/
5ef7bb5aaSJed Brown 
6ef7bb5aaSJed Brown PetscFList TSSSPList = 0;
7ef7bb5aaSJed Brown 
8ef7bb5aaSJed Brown typedef struct {
9ef7bb5aaSJed Brown   PetscErrorCode (*onestep)(TS,PetscReal,PetscReal,Vec);
105164425aSJed Brown   char *type_name;
11ef7bb5aaSJed Brown   PetscInt nstages;
12ef7bb5aaSJed Brown   Vec *work;
13ef7bb5aaSJed Brown   PetscInt nwork;
14ace3abfcSBarry Smith   PetscBool  workout;
15ef7bb5aaSJed Brown } TS_SSP;
16ef7bb5aaSJed Brown 
17ef7bb5aaSJed Brown 
18ef7bb5aaSJed Brown #undef __FUNCT__
19ef7bb5aaSJed Brown #define __FUNCT__ "SSPGetWorkVectors"
20ef7bb5aaSJed Brown static PetscErrorCode SSPGetWorkVectors(TS ts,PetscInt n,Vec **work)
21ef7bb5aaSJed Brown {
22ef7bb5aaSJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
23ef7bb5aaSJed Brown   PetscErrorCode ierr;
24ef7bb5aaSJed Brown 
25ef7bb5aaSJed Brown   PetscFunctionBegin;
26e32f2f54SBarry Smith   if (ssp->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Work vectors already gotten");
27ef7bb5aaSJed Brown   if (ssp->nwork < n) {
28ef7bb5aaSJed Brown     if (ssp->nwork > 0) {
29574deadeSBarry Smith       ierr = VecDestroyVecs(ssp->nwork,&ssp->work);CHKERRQ(ierr);
30ef7bb5aaSJed Brown     }
31ef7bb5aaSJed Brown     ierr = VecDuplicateVecs(ts->vec_sol,n,&ssp->work);CHKERRQ(ierr);
32ef7bb5aaSJed Brown     ssp->nwork = n;
33ef7bb5aaSJed Brown   }
34ef7bb5aaSJed Brown   *work = ssp->work;
35ef7bb5aaSJed Brown   ssp->workout = PETSC_TRUE;
36ef7bb5aaSJed Brown   PetscFunctionReturn(0);
37ef7bb5aaSJed Brown }
38ef7bb5aaSJed Brown 
39ef7bb5aaSJed Brown #undef __FUNCT__
40ef7bb5aaSJed Brown #define __FUNCT__ "SSPRestoreWorkVectors"
41ef7bb5aaSJed Brown static PetscErrorCode SSPRestoreWorkVectors(TS ts,PetscInt n,Vec **work)
42ef7bb5aaSJed Brown {
43ef7bb5aaSJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
44ef7bb5aaSJed Brown 
45ef7bb5aaSJed Brown   PetscFunctionBegin;
46e32f2f54SBarry Smith   if (!ssp->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Work vectors have not been gotten");
47e32f2f54SBarry Smith   if (*work != ssp->work) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong work vectors checked out");
48ef7bb5aaSJed Brown   ssp->workout = PETSC_FALSE;
49ef7bb5aaSJed Brown   *work = PETSC_NULL;
50ef7bb5aaSJed Brown   PetscFunctionReturn(0);
51ef7bb5aaSJed Brown }
52ef7bb5aaSJed Brown 
53ef7bb5aaSJed Brown 
54ef7bb5aaSJed Brown #undef __FUNCT__
55ef7bb5aaSJed Brown #define __FUNCT__ "SSPStep_RK_2"
56815f1ad5SJed Brown /*MC
57815f1ad5SJed Brown    TSSSPRKS2 - Optimal second order SSP Runge-Kutta method, low-storage, c_eff=(s-1)/s
58815f1ad5SJed Brown 
59815f1ad5SJed Brown    Pseudocode 2 of Ketcheson 2008
60815f1ad5SJed Brown 
61*b330ce4dSSatish Balay    Level: beginner
62*b330ce4dSSatish Balay 
63815f1ad5SJed Brown .seealso: TSSSP, TSSSPSetType(), TSSSPSetNumStages()
64815f1ad5SJed Brown M*/
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"
89815f1ad5SJed Brown /*MC
90815f1ad5SJed Brown    TSSSPRKS3 - Optimal third order SSP Runge-Kutta, low-storage, c_eff=(sqrt(s)-1)/sqrt(s), where sqrt(s) is an integer
91815f1ad5SJed Brown 
92815f1ad5SJed Brown    Pseudocode 2 of Ketcheson 2008
93815f1ad5SJed Brown 
94*b330ce4dSSatish Balay    Level: beginner
95*b330ce4dSSatish Balay 
96815f1ad5SJed Brown .seealso: TSSSP, TSSSPSetType(), TSSSPSetNumStages()
97815f1ad5SJed Brown M*/
98ef7bb5aaSJed Brown static PetscErrorCode SSPStep_RK_3(TS ts,PetscReal t0,PetscReal dt,Vec sol)
99ef7bb5aaSJed Brown {
100ef7bb5aaSJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
101ef7bb5aaSJed Brown   Vec *work,F;
102ef7bb5aaSJed Brown   PetscInt i,s,n,r;
103ef7bb5aaSJed Brown   PetscReal c;
104ef7bb5aaSJed Brown   PetscErrorCode ierr;
105ef7bb5aaSJed Brown 
106ef7bb5aaSJed Brown   PetscFunctionBegin;
107ef7bb5aaSJed Brown   s = ssp->nstages;
108fad8df86SSatish Balay   n = (PetscInt)(sqrt((PetscReal)s)+0.001);
109ef7bb5aaSJed Brown   r = s-n;
110e32f2f54SBarry 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);
111ef7bb5aaSJed Brown   ierr = SSPGetWorkVectors(ts,3,&work);CHKERRQ(ierr);
112ef7bb5aaSJed Brown   F = work[2];
113ef7bb5aaSJed Brown   ierr = VecCopy(sol,work[0]);CHKERRQ(ierr);
114ef7bb5aaSJed Brown   for (i=0; i<(n-1)*(n-2)/2; i++) {
115ef7bb5aaSJed Brown     c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
116ef7bb5aaSJed Brown     ierr = TSComputeRHSFunction(ts,t0+c*dt,work[0],F);CHKERRQ(ierr);
117ef7bb5aaSJed Brown     ierr = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr);
118ef7bb5aaSJed Brown   }
119ef7bb5aaSJed Brown   ierr = VecCopy(work[0],work[1]);CHKERRQ(ierr);
120ef7bb5aaSJed Brown   for ( ; i<n*(n+1)/2-1; i++) {
121ef7bb5aaSJed Brown     c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
122ef7bb5aaSJed Brown     ierr = TSComputeRHSFunction(ts,t0+c*dt,work[0],F);CHKERRQ(ierr);
123ef7bb5aaSJed Brown     ierr = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr);
124ef7bb5aaSJed Brown   }
125ef7bb5aaSJed Brown   {
126ef7bb5aaSJed Brown     c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
127ef7bb5aaSJed Brown     ierr = TSComputeRHSFunction(ts,t0+c*dt,work[0],F);CHKERRQ(ierr);
128ef7bb5aaSJed 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);
129ef7bb5aaSJed Brown     i++;
130ef7bb5aaSJed Brown   }
131ef7bb5aaSJed Brown   for ( ; i<s; i++) {
132ef7bb5aaSJed Brown     c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
133ef7bb5aaSJed Brown     ierr = TSComputeRHSFunction(ts,t0+c*dt,work[0],F);CHKERRQ(ierr);
134ef7bb5aaSJed Brown     ierr = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr);
135ef7bb5aaSJed Brown   }
136ef7bb5aaSJed Brown   ierr = VecCopy(work[0],sol);CHKERRQ(ierr);
137ef7bb5aaSJed Brown   ierr = SSPRestoreWorkVectors(ts,3,&work);CHKERRQ(ierr);
138ef7bb5aaSJed Brown   PetscFunctionReturn(0);
139ef7bb5aaSJed Brown }
140ef7bb5aaSJed Brown 
141ef7bb5aaSJed Brown #undef __FUNCT__
142ef7bb5aaSJed Brown #define __FUNCT__ "SSPStep_RK_10_4"
143815f1ad5SJed Brown /*MC
144*b330ce4dSSatish Balay    TSSSPRKS104 - Optimal fourth order SSP Runge-Kutta, low-storage (2N), c_eff=0.6
145815f1ad5SJed Brown 
146815f1ad5SJed Brown    SSPRK(10,4), Pseudocode 3 of Ketcheson 2008
147815f1ad5SJed Brown 
148*b330ce4dSSatish Balay    Level: beginner
149*b330ce4dSSatish Balay 
150815f1ad5SJed Brown .seealso: TSSSP, TSSSPSetType()
151815f1ad5SJed Brown M*/
152ef7bb5aaSJed Brown static PetscErrorCode SSPStep_RK_10_4(TS ts,PetscReal t0,PetscReal dt,Vec sol)
153ef7bb5aaSJed Brown {
154ef7bb5aaSJed Brown   const PetscReal c[10] = {0, 1./6, 2./6, 3./6, 4./6, 2./6, 3./6, 4./6, 5./6, 1};
155ef7bb5aaSJed Brown   Vec *work,F;
1565a586d82SBarry Smith   PetscInt i;
157ef7bb5aaSJed Brown   PetscErrorCode ierr;
158ef7bb5aaSJed Brown 
159ef7bb5aaSJed Brown   PetscFunctionBegin;
160ef7bb5aaSJed Brown   ierr = SSPGetWorkVectors(ts,3,&work);CHKERRQ(ierr);
161ef7bb5aaSJed Brown   F = work[2];
162ef7bb5aaSJed Brown   ierr = VecCopy(sol,work[0]);CHKERRQ(ierr);
163ef7bb5aaSJed Brown   for (i=0; i<5; i++) {
164ef7bb5aaSJed Brown     ierr = TSComputeRHSFunction(ts,t0+c[i],work[0],F);CHKERRQ(ierr);
165ef7bb5aaSJed Brown     ierr = VecAXPY(work[0],dt/6,F);CHKERRQ(ierr);
166ef7bb5aaSJed Brown   }
167ef7bb5aaSJed Brown   ierr = VecAXPBYPCZ(work[1],1./25,9./25,0,sol,work[0]);CHKERRQ(ierr);
168ef7bb5aaSJed Brown   ierr = VecAXPBY(work[0],15,-5,work[1]);CHKERRQ(ierr);
169ef7bb5aaSJed Brown   for ( ; i<9; i++) {
170ef7bb5aaSJed Brown     ierr = TSComputeRHSFunction(ts,t0+c[i],work[0],F);CHKERRQ(ierr);
171ef7bb5aaSJed Brown     ierr = VecAXPY(work[0],dt/6,F);CHKERRQ(ierr);
172ef7bb5aaSJed Brown   }
173ef7bb5aaSJed Brown   ierr = TSComputeRHSFunction(ts,t0+dt,work[0],F);CHKERRQ(ierr);
174ef7bb5aaSJed Brown   ierr = VecAXPBYPCZ(work[1],3./5,dt/10,1,work[0],F);CHKERRQ(ierr);
175ef7bb5aaSJed Brown   ierr = VecCopy(work[1],sol);CHKERRQ(ierr);
176ef7bb5aaSJed Brown   ierr = SSPRestoreWorkVectors(ts,3,&work);CHKERRQ(ierr);
177ef7bb5aaSJed Brown   PetscFunctionReturn(0);
178ef7bb5aaSJed Brown }
179ef7bb5aaSJed Brown 
180ef7bb5aaSJed Brown 
181ef7bb5aaSJed Brown #undef __FUNCT__
182ef7bb5aaSJed Brown #define __FUNCT__ "TSSetUp_SSP"
183ef7bb5aaSJed Brown static PetscErrorCode TSSetUp_SSP(TS ts)
184ef7bb5aaSJed Brown {
185ef7bb5aaSJed Brown 
186ef7bb5aaSJed Brown   PetscFunctionBegin;
187ef7bb5aaSJed Brown   PetscFunctionReturn(0);
188ef7bb5aaSJed Brown }
189ef7bb5aaSJed Brown 
190ef7bb5aaSJed Brown #undef __FUNCT__
191ef7bb5aaSJed Brown #define __FUNCT__ "TSStep_SSP"
192193ac0bcSJed Brown static PetscErrorCode TSStep_SSP(TS ts)
193ef7bb5aaSJed Brown {
194ef7bb5aaSJed Brown   TS_SSP        *ssp = (TS_SSP*)ts->data;
195ef7bb5aaSJed Brown   Vec            sol = ts->vec_sol;
196ef7bb5aaSJed Brown   PetscErrorCode ierr;
197ef7bb5aaSJed Brown 
198ef7bb5aaSJed Brown   PetscFunctionBegin;
199186e87acSLisandro Dalcin   ierr = (*ssp->onestep)(ts,ts->ptime,ts->time_step,sol);CHKERRQ(ierr);
200186e87acSLisandro Dalcin   ts->ptime += ts->time_step;
201ef7bb5aaSJed Brown   ts->steps++;
202ef7bb5aaSJed Brown   PetscFunctionReturn(0);
203ef7bb5aaSJed Brown }
204ef7bb5aaSJed Brown /*------------------------------------------------------------*/
205ef7bb5aaSJed Brown #undef __FUNCT__
206277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_SSP"
207277b19d0SLisandro Dalcin static PetscErrorCode TSReset_SSP(TS ts)
208277b19d0SLisandro Dalcin {
209277b19d0SLisandro Dalcin   TS_SSP         *ssp = (TS_SSP*)ts->data;
210277b19d0SLisandro Dalcin   PetscErrorCode ierr;
211277b19d0SLisandro Dalcin 
212277b19d0SLisandro Dalcin   PetscFunctionBegin;
213277b19d0SLisandro Dalcin   if (ssp->work) {ierr = VecDestroyVecs(ssp->nwork,&ssp->work);CHKERRQ(ierr);}
214277b19d0SLisandro Dalcin   ssp->nwork = 0;
215277b19d0SLisandro Dalcin   ssp->workout = PETSC_FALSE;
216277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
217277b19d0SLisandro Dalcin }
218277b19d0SLisandro Dalcin 
219277b19d0SLisandro Dalcin #undef __FUNCT__
220ef7bb5aaSJed Brown #define __FUNCT__ "TSDestroy_SSP"
221ef7bb5aaSJed Brown static PetscErrorCode TSDestroy_SSP(TS ts)
222ef7bb5aaSJed Brown {
223815f1ad5SJed Brown   TS_SSP         *ssp = (TS_SSP*)ts->data;
224ef7bb5aaSJed Brown   PetscErrorCode ierr;
225ef7bb5aaSJed Brown 
226ef7bb5aaSJed Brown   PetscFunctionBegin;
227277b19d0SLisandro Dalcin   ierr = TSReset_SSP(ts);CHKERRQ(ierr);
2285164425aSJed Brown   ierr = PetscFree(ssp->type_name);CHKERRQ(ierr);
229c31cb41cSBarry Smith   ierr = PetscFree(ts->data);CHKERRQ(ierr);
230815f1ad5SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSSPGetType_C","",PETSC_NULL);CHKERRQ(ierr);
231815f1ad5SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSSPSetType_C","",PETSC_NULL);CHKERRQ(ierr);
232815f1ad5SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSSPGetNumStages_C","",PETSC_NULL);CHKERRQ(ierr);
233815f1ad5SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSSPSetNumStages_C","",PETSC_NULL);CHKERRQ(ierr);
234ef7bb5aaSJed Brown   PetscFunctionReturn(0);
235ef7bb5aaSJed Brown }
236ef7bb5aaSJed Brown /*------------------------------------------------------------*/
237ef7bb5aaSJed Brown 
238ef7bb5aaSJed Brown #undef __FUNCT__
239ef7bb5aaSJed Brown #define __FUNCT__ "TSSSPSetType"
240815f1ad5SJed Brown /*@C
241815f1ad5SJed Brown    TSSSPSetType - set the SSP time integration scheme to use
242815f1ad5SJed Brown 
243815f1ad5SJed Brown    Logically Collective
244815f1ad5SJed Brown 
245815f1ad5SJed Brown    Input Arguments:
246815f1ad5SJed Brown    ts - time stepping object
247815f1ad5SJed Brown    type - type of scheme to use
248815f1ad5SJed Brown 
249815f1ad5SJed Brown    Options Database Keys:
250815f1ad5SJed Brown    -ts_ssp_type <rks2>: Type of SSP method (one of) rks2 rks3 rk104
251815f1ad5SJed Brown    -ts_ssp_nstages <5>: Number of stages
252815f1ad5SJed Brown 
253815f1ad5SJed Brown    Level: beginner
254815f1ad5SJed Brown 
255815f1ad5SJed Brown .seealso: TSSSP, TSSSPGetType(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104
256815f1ad5SJed Brown @*/
257815f1ad5SJed Brown PetscErrorCode TSSSPSetType(TS ts,const TSSSPType type)
258815f1ad5SJed Brown {
259815f1ad5SJed Brown   PetscErrorCode ierr;
260815f1ad5SJed Brown 
261815f1ad5SJed Brown   PetscFunctionBegin;
262815f1ad5SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
263815f1ad5SJed Brown   ierr = PetscTryMethod(ts,"TSSSPSetType_C",(TS,const TSSSPType),(ts,type));CHKERRQ(ierr);
264815f1ad5SJed Brown   PetscFunctionReturn(0);
265815f1ad5SJed Brown }
266815f1ad5SJed Brown 
267815f1ad5SJed Brown #undef __FUNCT__
268815f1ad5SJed Brown #define __FUNCT__ "TSSSPGetType"
269815f1ad5SJed Brown /*@C
270815f1ad5SJed Brown    TSSSPGetType - get the SSP time integration scheme
271815f1ad5SJed Brown 
272815f1ad5SJed Brown    Logically Collective
273815f1ad5SJed Brown 
274815f1ad5SJed Brown    Input Argument:
275815f1ad5SJed Brown    ts - time stepping object
276815f1ad5SJed Brown 
277815f1ad5SJed Brown    Output Argument:
278815f1ad5SJed Brown    type - type of scheme being used
279815f1ad5SJed Brown 
280815f1ad5SJed Brown    Level: beginner
281815f1ad5SJed Brown 
282815f1ad5SJed Brown .seealso: TSSSP, TSSSPSettype(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104
283815f1ad5SJed Brown @*/
284815f1ad5SJed Brown PetscErrorCode TSSSPGetType(TS ts,const TSSSPType *type)
285815f1ad5SJed Brown {
286815f1ad5SJed Brown   PetscErrorCode ierr;
287815f1ad5SJed Brown 
288815f1ad5SJed Brown   PetscFunctionBegin;
289815f1ad5SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
290815f1ad5SJed Brown   ierr = PetscTryMethod(ts,"TSSSPGetType_C",(TS,const TSSSPType*),(ts,type));CHKERRQ(ierr);
291815f1ad5SJed Brown   PetscFunctionReturn(0);
292815f1ad5SJed Brown }
293815f1ad5SJed Brown 
294815f1ad5SJed Brown #undef __FUNCT__
295815f1ad5SJed Brown #define __FUNCT__ "TSSSPSetNumStages"
296815f1ad5SJed Brown /*@
297815f1ad5SJed Brown    TSSSPSetNumStages - set the number of stages to use with the SSP method
298815f1ad5SJed Brown 
299815f1ad5SJed Brown    Logically Collective
300815f1ad5SJed Brown 
301815f1ad5SJed Brown    Input Arguments:
302815f1ad5SJed Brown    ts - time stepping object
303815f1ad5SJed Brown    nstages - number of stages
304815f1ad5SJed Brown 
305815f1ad5SJed Brown    Options Database Keys:
306815f1ad5SJed Brown    -ts_ssp_type <rks2>: NumStages of SSP method (one of) rks2 rks3 rk104
307815f1ad5SJed Brown    -ts_ssp_nstages <5>: Number of stages
308815f1ad5SJed Brown 
309815f1ad5SJed Brown    Level: beginner
310815f1ad5SJed Brown 
311815f1ad5SJed Brown .seealso: TSSSP, TSSSPGetNumStages(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104
312815f1ad5SJed Brown @*/
313815f1ad5SJed Brown PetscErrorCode TSSSPSetNumStages(TS ts,PetscInt nstages)
314815f1ad5SJed Brown {
315815f1ad5SJed Brown   PetscErrorCode ierr;
316815f1ad5SJed Brown 
317815f1ad5SJed Brown   PetscFunctionBegin;
318815f1ad5SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
319815f1ad5SJed Brown   ierr = PetscTryMethod(ts,"TSSSPSetNumStages_C",(TS,PetscInt),(ts,nstages));CHKERRQ(ierr);
320815f1ad5SJed Brown   PetscFunctionReturn(0);
321815f1ad5SJed Brown }
322815f1ad5SJed Brown 
323815f1ad5SJed Brown #undef __FUNCT__
324815f1ad5SJed Brown #define __FUNCT__ "TSSSPGetNumStages"
325815f1ad5SJed Brown /*@
326815f1ad5SJed Brown    TSSSPGetNumStages - get the number of stages in the SSP time integration scheme
327815f1ad5SJed Brown 
328815f1ad5SJed Brown    Logically Collective
329815f1ad5SJed Brown 
330815f1ad5SJed Brown    Input Argument:
331815f1ad5SJed Brown    ts - time stepping object
332815f1ad5SJed Brown 
333815f1ad5SJed Brown    Output Argument:
334815f1ad5SJed Brown    nstages - number of stages
335815f1ad5SJed Brown 
336815f1ad5SJed Brown    Level: beginner
337815f1ad5SJed Brown 
338815f1ad5SJed Brown .seealso: TSSSP, TSSSPGetType(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104
339815f1ad5SJed Brown @*/
340815f1ad5SJed Brown PetscErrorCode TSSSPGetNumStages(TS ts,PetscInt *nstages)
341815f1ad5SJed Brown {
342815f1ad5SJed Brown   PetscErrorCode ierr;
343815f1ad5SJed Brown 
344815f1ad5SJed Brown   PetscFunctionBegin;
345815f1ad5SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
346815f1ad5SJed Brown   ierr = PetscTryMethod(ts,"TSSSPGetNumStages_C",(TS,PetscInt*),(ts,nstages));CHKERRQ(ierr);
347815f1ad5SJed Brown   PetscFunctionReturn(0);
348815f1ad5SJed Brown }
349815f1ad5SJed Brown 
350815f1ad5SJed Brown EXTERN_C_BEGIN
351815f1ad5SJed Brown #undef __FUNCT__
352815f1ad5SJed Brown #define __FUNCT__ "TSSSPSetType_SSP"
353815f1ad5SJed Brown PetscErrorCode TSSSPSetType_SSP(TS ts,const TSSSPType type)
354ef7bb5aaSJed Brown {
355ef7bb5aaSJed Brown   PetscErrorCode ierr,(*r)(TS,PetscReal,PetscReal,Vec);
356ef7bb5aaSJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
357ef7bb5aaSJed Brown 
358ef7bb5aaSJed Brown   PetscFunctionBegin;
3594b91b6eaSBarry Smith   ierr = PetscFListFind(TSSSPList,((PetscObject)ts)->comm,type,PETSC_TRUE,(void(**)(void))&r);CHKERRQ(ierr);
360e32f2f54SBarry Smith   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TS_SSP type %s given",type);
361ef7bb5aaSJed Brown   ssp->onestep = r;
3625164425aSJed Brown   ierr = PetscFree(ssp->type_name);CHKERRQ(ierr);
3635164425aSJed Brown   ierr = PetscStrallocpy(type,&ssp->type_name);CHKERRQ(ierr);
364ef7bb5aaSJed Brown   PetscFunctionReturn(0);
365ef7bb5aaSJed Brown }
366815f1ad5SJed Brown #undef __FUNCT__
367815f1ad5SJed Brown #define __FUNCT__ "TSSSPGetType_SSP"
368815f1ad5SJed Brown PetscErrorCode TSSSPGetType_SSP(TS ts,const TSSSPType *type)
369815f1ad5SJed Brown {
370815f1ad5SJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
371815f1ad5SJed Brown 
372815f1ad5SJed Brown   PetscFunctionBegin;
3735164425aSJed Brown   *type = ssp->type_name;
374815f1ad5SJed Brown   PetscFunctionReturn(0);
375815f1ad5SJed Brown }
376815f1ad5SJed Brown #undef __FUNCT__
377815f1ad5SJed Brown #define __FUNCT__ "TSSSPSetNumStages_SSP"
378815f1ad5SJed Brown PetscErrorCode TSSSPSetNumStages_SSP(TS ts,PetscInt nstages)
379815f1ad5SJed Brown {
380815f1ad5SJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
381815f1ad5SJed Brown 
382815f1ad5SJed Brown   PetscFunctionBegin;
383815f1ad5SJed Brown   ssp->nstages = nstages;
384815f1ad5SJed Brown   PetscFunctionReturn(0);
385815f1ad5SJed Brown }
386815f1ad5SJed Brown #undef __FUNCT__
387815f1ad5SJed Brown #define __FUNCT__ "TSSSPGetNumStages_SSP"
388815f1ad5SJed Brown PetscErrorCode TSSSPGetNumStages_SSP(TS ts,PetscInt *nstages)
389815f1ad5SJed Brown {
390815f1ad5SJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
391815f1ad5SJed Brown 
392815f1ad5SJed Brown   PetscFunctionBegin;
393815f1ad5SJed Brown   *nstages = ssp->nstages;
394815f1ad5SJed Brown   PetscFunctionReturn(0);
395815f1ad5SJed Brown }
3965164425aSJed Brown EXTERN_C_END
397ef7bb5aaSJed Brown 
398ef7bb5aaSJed Brown #undef __FUNCT__
399ef7bb5aaSJed Brown #define __FUNCT__ "TSSetFromOptions_SSP"
400ef7bb5aaSJed Brown static PetscErrorCode TSSetFromOptions_SSP(TS ts)
401ef7bb5aaSJed Brown {
402ef7bb5aaSJed Brown   char tname[256] = TSSSPRKS2;
403ef7bb5aaSJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
404ef7bb5aaSJed Brown   PetscErrorCode ierr;
405ace3abfcSBarry Smith   PetscBool  flg;
406ef7bb5aaSJed Brown 
407ef7bb5aaSJed Brown   PetscFunctionBegin;
408ef7bb5aaSJed Brown   ierr = PetscOptionsHead("SSP ODE solver options");CHKERRQ(ierr);
409ef7bb5aaSJed Brown   {
410ef7bb5aaSJed Brown     ierr = PetscOptionsList("-ts_ssp_type","Type of SSP method","TSSSPSetType",TSSSPList,tname,tname,sizeof(tname),&flg);CHKERRQ(ierr);
411ef7bb5aaSJed Brown     if (flg) {
412ef7bb5aaSJed Brown       ierr = TSSSPSetType(ts,tname);CHKERRQ(ierr);
413ef7bb5aaSJed Brown     }
414ef7bb5aaSJed Brown     ierr = PetscOptionsInt("-ts_ssp_nstages","Number of stages","TSSSPSetNumStages",ssp->nstages,&ssp->nstages,PETSC_NULL);CHKERRQ(ierr);
415ef7bb5aaSJed Brown   }
416ef7bb5aaSJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
417ef7bb5aaSJed Brown   PetscFunctionReturn(0);
418ef7bb5aaSJed Brown }
419ef7bb5aaSJed Brown 
420ef7bb5aaSJed Brown #undef __FUNCT__
421ef7bb5aaSJed Brown #define __FUNCT__ "TSView_SSP"
422ef7bb5aaSJed Brown static PetscErrorCode TSView_SSP(TS ts,PetscViewer viewer)
423ef7bb5aaSJed Brown {
424ef7bb5aaSJed Brown   PetscFunctionBegin;
425ef7bb5aaSJed Brown   PetscFunctionReturn(0);
426ef7bb5aaSJed Brown }
427ef7bb5aaSJed Brown 
428ef7bb5aaSJed Brown /* ------------------------------------------------------------ */
429ef7bb5aaSJed Brown 
430ef7bb5aaSJed Brown /*MC
4318ab3e0fcSJed Brown       TSSSP - Explicit strong stability preserving ODE solver
4328ab3e0fcSJed Brown 
4338ab3e0fcSJed Brown   Most hyperbolic conservation laws have exact solutions that are total variation diminishing (TVD) or total variation
4348ab3e0fcSJed Brown   bounded (TVB) although these solutions often contain discontinuities.  Spatial discretizations such as Godunov's
4358ab3e0fcSJed Brown   scheme and high-resolution finite volume methods (TVD limiters, ENO/WENO) are designed to preserve these properties,
4368ab3e0fcSJed Brown   but they are usually formulated using a forward Euler time discretization or by coupling the space and time
4378ab3e0fcSJed Brown   discretization as in the classical Lax-Wendroff scheme.  When the space and time discretization is coupled, it is very
4388ab3e0fcSJed Brown   difficult to produce schemes with high temporal accuracy while preserving TVD properties.  An alternative is the
4398ab3e0fcSJed Brown   semidiscrete formulation where we choose a spatial discretization that is TVD with forward Euler and then choose a
4408ab3e0fcSJed Brown   time discretization that preserves the TVD property.  Such integrators are called strong stability preserving (SSP).
4418ab3e0fcSJed Brown 
4428ab3e0fcSJed Brown   Let c_eff be the minimum number of function evaluations required to step as far as one step of forward Euler while
4438ab3e0fcSJed Brown   still being SSP.  Some theoretical bounds
4448ab3e0fcSJed Brown 
4458ab3e0fcSJed Brown   1. There are no explicit methods with c_eff > 1.
4460085e20eSJed Brown 
4478ab3e0fcSJed Brown   2. There are no explicit methods beyond order 4 (for nonlinear problems) and c_eff > 0.
4480085e20eSJed Brown 
4498ab3e0fcSJed Brown   3. There are no implicit methods with order greater than 1 and c_eff > 2.
4508ab3e0fcSJed Brown 
4518ab3e0fcSJed Brown   This integrator provides Runge-Kutta methods of order 2, 3, and 4 with maximal values of c_eff.  More stages allows
4528ab3e0fcSJed Brown   for larger values of c_eff which improves efficiency.  These implementations are low-memory and only use 2 or 3 work
4538ab3e0fcSJed Brown   vectors regardless of the total number of stages, so e.g. 25-stage 3rd order methods may be an excellent choice.
4548ab3e0fcSJed Brown 
4558ab3e0fcSJed Brown   Methods can be chosen with -ts_ssp_type {rks2,rks3,rk104}
4568ab3e0fcSJed Brown 
4578ab3e0fcSJed Brown   rks2: Second order methods with any number s>1 of stages.  c_eff = (s-1)/s
4588ab3e0fcSJed Brown 
4598ab3e0fcSJed Brown   rks3: Third order methods with s=n^2 stages, n>1.  c_eff = (s-n)/s
4608ab3e0fcSJed Brown 
4618ab3e0fcSJed Brown   rk104: A 10-stage fourth order method.  c_eff = 0.6
462ef7bb5aaSJed Brown 
463ef7bb5aaSJed Brown   Level: beginner
464ef7bb5aaSJed Brown 
4657b6bb2c6SJed Brown   References:
4667b6bb2c6SJed Brown   Ketcheson, Highly efficient strong stability preserving Runge-Kutta methods with low-storage implementations, SISC, 2008.
4677b6bb2c6SJed Brown 
4687b6bb2c6SJed Brown   Gottlieb, Ketcheson, and Shu, High order strong stability preserving time discretizations, J Scientific Computing, 2009.
4697b6bb2c6SJed Brown 
470ef7bb5aaSJed Brown .seealso:  TSCreate(), TS, TSSetType()
471ef7bb5aaSJed Brown 
472ef7bb5aaSJed Brown M*/
473ef7bb5aaSJed Brown EXTERN_C_BEGIN
474ef7bb5aaSJed Brown #undef __FUNCT__
475ef7bb5aaSJed Brown #define __FUNCT__ "TSCreate_SSP"
4767087cfbeSBarry Smith PetscErrorCode  TSCreate_SSP(TS ts)
477ef7bb5aaSJed Brown {
478ef7bb5aaSJed Brown   TS_SSP       *ssp;
479ef7bb5aaSJed Brown   PetscErrorCode ierr;
480ef7bb5aaSJed Brown 
481ef7bb5aaSJed Brown   PetscFunctionBegin;
482ef7bb5aaSJed Brown   if (!TSSSPList) {
483ef7bb5aaSJed Brown     ierr = PetscFListAdd(&TSSSPList,TSSSPRKS2,  "SSPStep_RK_2",   (void(*)(void))SSPStep_RK_2);CHKERRQ(ierr);
484ef7bb5aaSJed Brown     ierr = PetscFListAdd(&TSSSPList,TSSSPRKS3,  "SSPStep_RK_3",   (void(*)(void))SSPStep_RK_3);CHKERRQ(ierr);
485ef7bb5aaSJed Brown     ierr = PetscFListAdd(&TSSSPList,TSSSPRK104, "SSPStep_RK_10_4",(void(*)(void))SSPStep_RK_10_4);CHKERRQ(ierr);
486ef7bb5aaSJed Brown   }
487ef7bb5aaSJed Brown 
488ef7bb5aaSJed Brown   ts->ops->setup           = TSSetUp_SSP;
489ef7bb5aaSJed Brown   ts->ops->step            = TSStep_SSP;
490277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_SSP;
491ef7bb5aaSJed Brown   ts->ops->destroy         = TSDestroy_SSP;
492ef7bb5aaSJed Brown   ts->ops->setfromoptions  = TSSetFromOptions_SSP;
493ef7bb5aaSJed Brown   ts->ops->view            = TSView_SSP;
494ef7bb5aaSJed Brown 
495ef7bb5aaSJed Brown   ierr = PetscNewLog(ts,TS_SSP,&ssp);CHKERRQ(ierr);
496ef7bb5aaSJed Brown   ts->data = (void*)ssp;
497ef7bb5aaSJed Brown 
498815f1ad5SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSSPGetType_C","TSSSPGetType_SSP",TSSSPGetType_SSP);CHKERRQ(ierr);
499815f1ad5SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSSPSetType_C","TSSSPSetType_SSP",TSSSPSetType_SSP);CHKERRQ(ierr);
500815f1ad5SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSSPGetNumStages_C","TSSSPGetNumStages_SSP",TSSSPGetNumStages_SSP);CHKERRQ(ierr);
501815f1ad5SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSSPSetNumStages_C","TSSSPSetNumStages_SSP",TSSSPSetNumStages_SSP);CHKERRQ(ierr);
502815f1ad5SJed Brown 
503ef7bb5aaSJed Brown   ierr = TSSSPSetType(ts,TSSSPRKS2);CHKERRQ(ierr);
504ef7bb5aaSJed Brown   ssp->nstages = 5;
505ef7bb5aaSJed Brown   PetscFunctionReturn(0);
506ef7bb5aaSJed Brown }
507ef7bb5aaSJed Brown EXTERN_C_END
508