xref: /petsc/src/ts/impls/explicit/ssp/ssp.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
1ef7bb5aaSJed Brown /*
2ef7bb5aaSJed Brown        Code for Timestepping with explicit SSP.
3ef7bb5aaSJed Brown */
4af0996ceSBarry Smith #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
5ef7bb5aaSJed Brown 
6c793f718SLisandro Dalcin PetscFunctionList TSSSPList = NULL;
7787849ffSJed Brown static PetscBool TSSSPPackageInitialized;
8ef7bb5aaSJed Brown 
9ef7bb5aaSJed Brown typedef struct {
10ef7bb5aaSJed Brown   PetscErrorCode (*onestep)(TS,PetscReal,PetscReal,Vec);
115164425aSJed Brown   char           *type_name;
12ef7bb5aaSJed Brown   PetscInt       nstages;
13ef7bb5aaSJed Brown   Vec            *work;
14ef7bb5aaSJed Brown   PetscInt       nwork;
15ace3abfcSBarry Smith   PetscBool      workout;
16ef7bb5aaSJed Brown } TS_SSP;
17ef7bb5aaSJed Brown 
1805d28066SBarry Smith static PetscErrorCode TSSSPGetWorkVectors(TS ts,PetscInt n,Vec **work)
19ef7bb5aaSJed Brown {
20ef7bb5aaSJed Brown   TS_SSP         *ssp = (TS_SSP*)ts->data;
21ef7bb5aaSJed Brown   PetscErrorCode ierr;
22ef7bb5aaSJed Brown 
23ef7bb5aaSJed Brown   PetscFunctionBegin;
24*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(ssp->workout,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Work vectors already gotten");
25ef7bb5aaSJed Brown   if (ssp->nwork < n) {
26ef7bb5aaSJed Brown     if (ssp->nwork > 0) {
27574deadeSBarry Smith       ierr = VecDestroyVecs(ssp->nwork,&ssp->work);CHKERRQ(ierr);
28ef7bb5aaSJed Brown     }
29ef7bb5aaSJed Brown     ierr = VecDuplicateVecs(ts->vec_sol,n,&ssp->work);CHKERRQ(ierr);
30ef7bb5aaSJed Brown     ssp->nwork = n;
31ef7bb5aaSJed Brown   }
32ef7bb5aaSJed Brown   *work = ssp->work;
33ef7bb5aaSJed Brown   ssp->workout = PETSC_TRUE;
34ef7bb5aaSJed Brown   PetscFunctionReturn(0);
35ef7bb5aaSJed Brown }
36ef7bb5aaSJed Brown 
3705d28066SBarry Smith static PetscErrorCode TSSSPRestoreWorkVectors(TS ts,PetscInt n,Vec **work)
38ef7bb5aaSJed Brown {
39ef7bb5aaSJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
40ef7bb5aaSJed Brown 
41ef7bb5aaSJed Brown   PetscFunctionBegin;
42*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!ssp->workout,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Work vectors have not been gotten");
43*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(*work != ssp->work,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong work vectors checked out");
44ef7bb5aaSJed Brown   ssp->workout = PETSC_FALSE;
450298fd71SBarry Smith   *work = NULL;
46ef7bb5aaSJed Brown   PetscFunctionReturn(0);
47ef7bb5aaSJed Brown }
48ef7bb5aaSJed Brown 
49815f1ad5SJed Brown /*MC
50815f1ad5SJed Brown    TSSSPRKS2 - Optimal second order SSP Runge-Kutta method, low-storage, c_eff=(s-1)/s
51815f1ad5SJed Brown 
52815f1ad5SJed Brown    Pseudocode 2 of Ketcheson 2008
53815f1ad5SJed Brown 
54b330ce4dSSatish Balay    Level: beginner
55b330ce4dSSatish Balay 
56815f1ad5SJed Brown .seealso: TSSSP, TSSSPSetType(), TSSSPSetNumStages()
57815f1ad5SJed Brown M*/
5805d28066SBarry Smith static PetscErrorCode TSSSPStep_RK_2(TS ts,PetscReal t0,PetscReal dt,Vec sol)
59ef7bb5aaSJed Brown {
60ef7bb5aaSJed Brown   TS_SSP         *ssp = (TS_SSP*)ts->data;
61ef7bb5aaSJed Brown   Vec            *work,F;
62ef7bb5aaSJed Brown   PetscInt       i,s;
63ef7bb5aaSJed Brown   PetscErrorCode ierr;
64ef7bb5aaSJed Brown 
65ef7bb5aaSJed Brown   PetscFunctionBegin;
66ef7bb5aaSJed Brown   s    = ssp->nstages;
6705d28066SBarry Smith   ierr = TSSSPGetWorkVectors(ts,2,&work);CHKERRQ(ierr);
68ef7bb5aaSJed Brown   F    = work[1];
69ef7bb5aaSJed Brown   ierr = VecCopy(sol,work[0]);CHKERRQ(ierr);
70ef7bb5aaSJed Brown   for (i=0; i<s-1; i++) {
71b8123daeSJed Brown     PetscReal stage_time = t0+dt*(i/(s-1.));
72b8123daeSJed Brown     ierr = TSPreStage(ts,stage_time);CHKERRQ(ierr);
73b8123daeSJed Brown     ierr = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr);
74ef7bb5aaSJed Brown     ierr = VecAXPY(work[0],dt/(s-1.),F);CHKERRQ(ierr);
75ef7bb5aaSJed Brown   }
76ef7bb5aaSJed Brown   ierr = TSComputeRHSFunction(ts,t0+dt,work[0],F);CHKERRQ(ierr);
77ef7bb5aaSJed Brown   ierr = VecAXPBYPCZ(sol,(s-1.)/s,dt/s,1./s,work[0],F);CHKERRQ(ierr);
7805d28066SBarry Smith   ierr = TSSSPRestoreWorkVectors(ts,2,&work);CHKERRQ(ierr);
79ef7bb5aaSJed Brown   PetscFunctionReturn(0);
80ef7bb5aaSJed Brown }
81ef7bb5aaSJed Brown 
82815f1ad5SJed Brown /*MC
83aaf9cf16SJed Brown    TSSSPRKS3 - Optimal third order SSP Runge-Kutta, low-storage, c_eff=(PetscSqrtReal(s)-1)/PetscSqrtReal(s), where PetscSqrtReal(s) is an integer
84815f1ad5SJed Brown 
85815f1ad5SJed Brown    Pseudocode 2 of Ketcheson 2008
86815f1ad5SJed Brown 
87b330ce4dSSatish Balay    Level: beginner
88b330ce4dSSatish Balay 
89815f1ad5SJed Brown .seealso: TSSSP, TSSSPSetType(), TSSSPSetNumStages()
90815f1ad5SJed Brown M*/
9105d28066SBarry Smith static PetscErrorCode TSSSPStep_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;
96b8123daeSJed Brown   PetscReal      c,stage_time;
97ef7bb5aaSJed Brown   PetscErrorCode ierr;
98ef7bb5aaSJed Brown 
99ef7bb5aaSJed Brown   PetscFunctionBegin;
100ef7bb5aaSJed Brown   s = ssp->nstages;
101aaf9cf16SJed Brown   n = (PetscInt)(PetscSqrtReal((PetscReal)s)+0.001);
102ef7bb5aaSJed Brown   r = s-n;
103*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(n*n != s,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);
10405d28066SBarry Smith   ierr = TSSSPGetWorkVectors(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);
109b8123daeSJed Brown     stage_time = t0+c*dt;
110b8123daeSJed Brown     ierr       = TSPreStage(ts,stage_time);CHKERRQ(ierr);
111b8123daeSJed Brown     ierr       = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr);
112ef7bb5aaSJed Brown     ierr       = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr);
113ef7bb5aaSJed Brown   }
114ef7bb5aaSJed Brown   ierr = VecCopy(work[0],work[1]);CHKERRQ(ierr);
115ef7bb5aaSJed Brown   for (; i<n*(n+1)/2-1; i++) {
116ef7bb5aaSJed Brown     c          = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
117b8123daeSJed Brown     stage_time = t0+c*dt;
118b8123daeSJed Brown     ierr       = TSPreStage(ts,stage_time);CHKERRQ(ierr);
119b8123daeSJed Brown     ierr       = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr);
120ef7bb5aaSJed Brown     ierr       = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr);
121ef7bb5aaSJed Brown   }
122ef7bb5aaSJed Brown   {
123ef7bb5aaSJed Brown     c          = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
124b8123daeSJed Brown     stage_time = t0+c*dt;
125b8123daeSJed Brown     ierr       = TSPreStage(ts,stage_time);CHKERRQ(ierr);
126b8123daeSJed Brown     ierr       = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr);
127ef7bb5aaSJed 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);
128ef7bb5aaSJed Brown     i++;
129ef7bb5aaSJed Brown   }
130ef7bb5aaSJed Brown   for (; i<s; i++) {
131ef7bb5aaSJed Brown     c          = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
132b8123daeSJed Brown     stage_time = t0+c*dt;
133b8123daeSJed Brown     ierr       = TSPreStage(ts,stage_time);CHKERRQ(ierr);
134b8123daeSJed Brown     ierr       = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr);
135ef7bb5aaSJed Brown     ierr       = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr);
136ef7bb5aaSJed Brown   }
137ef7bb5aaSJed Brown   ierr = VecCopy(work[0],sol);CHKERRQ(ierr);
13805d28066SBarry Smith   ierr = TSSSPRestoreWorkVectors(ts,3,&work);CHKERRQ(ierr);
139ef7bb5aaSJed Brown   PetscFunctionReturn(0);
140ef7bb5aaSJed Brown }
141ef7bb5aaSJed Brown 
142815f1ad5SJed Brown /*MC
143b330ce4dSSatish Balay    TSSSPRKS104 - Optimal fourth order SSP Runge-Kutta, low-storage (2N), c_eff=0.6
144815f1ad5SJed Brown 
145815f1ad5SJed Brown    SSPRK(10,4), Pseudocode 3 of Ketcheson 2008
146815f1ad5SJed Brown 
147b330ce4dSSatish Balay    Level: beginner
148b330ce4dSSatish Balay 
149815f1ad5SJed Brown .seealso: TSSSP, TSSSPSetType()
150815f1ad5SJed Brown M*/
15105d28066SBarry Smith static PetscErrorCode TSSSPStep_RK_10_4(TS ts,PetscReal t0,PetscReal dt,Vec sol)
152ef7bb5aaSJed Brown {
153ef7bb5aaSJed Brown   const PetscReal c[10] = {0, 1./6, 2./6, 3./6, 4./6, 2./6, 3./6, 4./6, 5./6, 1};
154ef7bb5aaSJed Brown   Vec             *work,F;
1555a586d82SBarry Smith   PetscInt        i;
156b8123daeSJed Brown   PetscReal       stage_time;
157ef7bb5aaSJed Brown   PetscErrorCode  ierr;
158ef7bb5aaSJed Brown 
159ef7bb5aaSJed Brown   PetscFunctionBegin;
16005d28066SBarry Smith   ierr = TSSSPGetWorkVectors(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++) {
164b8123daeSJed Brown     stage_time = t0+c[i]*dt;
165b8123daeSJed Brown     ierr       = TSPreStage(ts,stage_time);CHKERRQ(ierr);
166b8123daeSJed Brown     ierr       = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr);
167ef7bb5aaSJed Brown     ierr       = VecAXPY(work[0],dt/6,F);CHKERRQ(ierr);
168ef7bb5aaSJed Brown   }
169ef7bb5aaSJed Brown   ierr = VecAXPBYPCZ(work[1],1./25,9./25,0,sol,work[0]);CHKERRQ(ierr);
170ef7bb5aaSJed Brown   ierr = VecAXPBY(work[0],15,-5,work[1]);CHKERRQ(ierr);
171ef7bb5aaSJed Brown   for (; i<9; i++) {
172b8123daeSJed Brown     stage_time = t0+c[i]*dt;
173b8123daeSJed Brown     ierr       = TSPreStage(ts,stage_time);CHKERRQ(ierr);
174b8123daeSJed Brown     ierr       = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr);
175ef7bb5aaSJed Brown     ierr       = VecAXPY(work[0],dt/6,F);CHKERRQ(ierr);
176ef7bb5aaSJed Brown   }
177b8123daeSJed Brown   stage_time = t0+dt;
178b8123daeSJed Brown   ierr       = TSPreStage(ts,stage_time);CHKERRQ(ierr);
179b8123daeSJed Brown   ierr       = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr);
180ef7bb5aaSJed Brown   ierr       = VecAXPBYPCZ(work[1],3./5,dt/10,1,work[0],F);CHKERRQ(ierr);
181ef7bb5aaSJed Brown   ierr       = VecCopy(work[1],sol);CHKERRQ(ierr);
18205d28066SBarry Smith   ierr       = TSSSPRestoreWorkVectors(ts,3,&work);CHKERRQ(ierr);
183ef7bb5aaSJed Brown   PetscFunctionReturn(0);
184ef7bb5aaSJed Brown }
185ef7bb5aaSJed Brown 
186ef7bb5aaSJed Brown static PetscErrorCode TSSetUp_SSP(TS ts)
187ef7bb5aaSJed Brown {
1881566a47fSLisandro Dalcin   PetscErrorCode ierr;
189ef7bb5aaSJed Brown 
190ef7bb5aaSJed Brown   PetscFunctionBegin;
1913dd83b38SBarry Smith   ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr);
1921566a47fSLisandro Dalcin   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
1931566a47fSLisandro Dalcin   ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr);
194ef7bb5aaSJed Brown   PetscFunctionReturn(0);
195ef7bb5aaSJed Brown }
196ef7bb5aaSJed Brown 
197193ac0bcSJed Brown static PetscErrorCode TSStep_SSP(TS ts)
198ef7bb5aaSJed Brown {
199ef7bb5aaSJed Brown   TS_SSP         *ssp = (TS_SSP*)ts->data;
200ef7bb5aaSJed Brown   Vec            sol  = ts->vec_sol;
2011566a47fSLisandro Dalcin   PetscBool      stageok,accept = PETSC_TRUE;
2021566a47fSLisandro Dalcin   PetscReal      next_time_step = ts->time_step;
203ef7bb5aaSJed Brown   PetscErrorCode ierr;
204ef7bb5aaSJed Brown 
205ef7bb5aaSJed Brown   PetscFunctionBegin;
206186e87acSLisandro Dalcin   ierr = (*ssp->onestep)(ts,ts->ptime,ts->time_step,sol);CHKERRQ(ierr);
2071566a47fSLisandro Dalcin   ierr = TSPostStage(ts,ts->ptime,0,&sol);CHKERRQ(ierr);
2081566a47fSLisandro Dalcin   ierr = TSAdaptCheckStage(ts->adapt,ts,ts->ptime+ts->time_step,sol,&stageok);CHKERRQ(ierr);
2091566a47fSLisandro Dalcin   if (!stageok) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);}
2101566a47fSLisandro Dalcin 
2111566a47fSLisandro Dalcin   ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
2121566a47fSLisandro Dalcin   if (!accept) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);}
2131566a47fSLisandro Dalcin 
214186e87acSLisandro Dalcin   ts->ptime += ts->time_step;
2151566a47fSLisandro Dalcin   ts->time_step = next_time_step;
216ef7bb5aaSJed Brown   PetscFunctionReturn(0);
217ef7bb5aaSJed Brown }
218ef7bb5aaSJed Brown /*------------------------------------------------------------*/
2191566a47fSLisandro Dalcin 
220277b19d0SLisandro Dalcin static PetscErrorCode TSReset_SSP(TS ts)
221277b19d0SLisandro Dalcin {
222277b19d0SLisandro Dalcin   TS_SSP         *ssp = (TS_SSP*)ts->data;
223277b19d0SLisandro Dalcin   PetscErrorCode ierr;
224277b19d0SLisandro Dalcin 
225277b19d0SLisandro Dalcin   PetscFunctionBegin;
226277b19d0SLisandro Dalcin   if (ssp->work) {ierr = VecDestroyVecs(ssp->nwork,&ssp->work);CHKERRQ(ierr);}
227277b19d0SLisandro Dalcin   ssp->nwork   = 0;
228277b19d0SLisandro Dalcin   ssp->workout = PETSC_FALSE;
229277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
230277b19d0SLisandro Dalcin }
231277b19d0SLisandro Dalcin 
232ef7bb5aaSJed Brown static PetscErrorCode TSDestroy_SSP(TS ts)
233ef7bb5aaSJed Brown {
234815f1ad5SJed Brown   TS_SSP         *ssp = (TS_SSP*)ts->data;
235ef7bb5aaSJed Brown   PetscErrorCode ierr;
236ef7bb5aaSJed Brown 
237ef7bb5aaSJed Brown   PetscFunctionBegin;
238277b19d0SLisandro Dalcin   ierr = TSReset_SSP(ts);CHKERRQ(ierr);
2395164425aSJed Brown   ierr = PetscFree(ssp->type_name);CHKERRQ(ierr);
240c31cb41cSBarry Smith   ierr = PetscFree(ts->data);CHKERRQ(ierr);
241bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetType_C",NULL);CHKERRQ(ierr);
242bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetType_C",NULL);CHKERRQ(ierr);
243bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetNumStages_C",NULL);CHKERRQ(ierr);
244bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetNumStages_C",NULL);CHKERRQ(ierr);
245ef7bb5aaSJed Brown   PetscFunctionReturn(0);
246ef7bb5aaSJed Brown }
247ef7bb5aaSJed Brown /*------------------------------------------------------------*/
248ef7bb5aaSJed Brown 
249815f1ad5SJed Brown /*@C
250815f1ad5SJed Brown    TSSSPSetType - set the SSP time integration scheme to use
251815f1ad5SJed Brown 
252815f1ad5SJed Brown    Logically Collective
253815f1ad5SJed Brown 
2544165533cSJose E. Roman    Input Parameters:
2554165533cSJose E. Roman +  ts - time stepping object
2564165533cSJose E. Roman -  ssptype - type of scheme to use
257815f1ad5SJed Brown 
258815f1ad5SJed Brown    Options Database Keys:
259815f1ad5SJed Brown    -ts_ssp_type <rks2>: Type of SSP method (one of) rks2 rks3 rk104
260815f1ad5SJed Brown    -ts_ssp_nstages <5>: Number of stages
261815f1ad5SJed Brown 
262815f1ad5SJed Brown    Level: beginner
263815f1ad5SJed Brown 
264815f1ad5SJed Brown .seealso: TSSSP, TSSSPGetType(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104
265815f1ad5SJed Brown @*/
266b92453a8SLisandro Dalcin PetscErrorCode TSSSPSetType(TS ts,TSSSPType ssptype)
267815f1ad5SJed Brown {
268815f1ad5SJed Brown   PetscErrorCode ierr;
269815f1ad5SJed Brown 
270815f1ad5SJed Brown   PetscFunctionBegin;
271815f1ad5SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
272b92453a8SLisandro Dalcin   PetscValidCharPointer(ssptype,2);
273b92453a8SLisandro Dalcin   ierr = PetscTryMethod(ts,"TSSSPSetType_C",(TS,TSSSPType),(ts,ssptype));CHKERRQ(ierr);
274815f1ad5SJed Brown   PetscFunctionReturn(0);
275815f1ad5SJed Brown }
276815f1ad5SJed Brown 
277815f1ad5SJed Brown /*@C
278815f1ad5SJed Brown    TSSSPGetType - get the SSP time integration scheme
279815f1ad5SJed Brown 
280815f1ad5SJed Brown    Logically Collective
281815f1ad5SJed Brown 
2824165533cSJose E. Roman    Input Parameter:
2834165533cSJose E. Roman .  ts - time stepping object
284815f1ad5SJed Brown 
2854165533cSJose E. Roman    Output Parameter:
2864165533cSJose E. Roman .  type - type of scheme being used
287815f1ad5SJed Brown 
288815f1ad5SJed Brown    Level: beginner
289815f1ad5SJed Brown 
290815f1ad5SJed Brown .seealso: TSSSP, TSSSPSettype(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104
291815f1ad5SJed Brown @*/
29219fd82e9SBarry Smith PetscErrorCode TSSSPGetType(TS ts,TSSSPType *type)
293815f1ad5SJed Brown {
294815f1ad5SJed Brown   PetscErrorCode ierr;
295815f1ad5SJed Brown 
296815f1ad5SJed Brown   PetscFunctionBegin;
297815f1ad5SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
298163d334eSBarry Smith   ierr = PetscUseMethod(ts,"TSSSPGetType_C",(TS,TSSSPType*),(ts,type));CHKERRQ(ierr);
299815f1ad5SJed Brown   PetscFunctionReturn(0);
300815f1ad5SJed Brown }
301815f1ad5SJed Brown 
302815f1ad5SJed Brown /*@
303815f1ad5SJed Brown    TSSSPSetNumStages - set the number of stages to use with the SSP method
304815f1ad5SJed Brown 
305815f1ad5SJed Brown    Logically Collective
306815f1ad5SJed Brown 
3074165533cSJose E. Roman    Input Parameters:
3084165533cSJose E. Roman +  ts - time stepping object
3094165533cSJose E. Roman -  nstages - number of stages
310815f1ad5SJed Brown 
311815f1ad5SJed Brown    Options Database Keys:
312815f1ad5SJed Brown    -ts_ssp_type <rks2>: NumStages of SSP method (one of) rks2 rks3 rk104
313815f1ad5SJed Brown    -ts_ssp_nstages <5>: Number of stages
314815f1ad5SJed Brown 
315815f1ad5SJed Brown    Level: beginner
316815f1ad5SJed Brown 
317815f1ad5SJed Brown .seealso: TSSSP, TSSSPGetNumStages(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104
318815f1ad5SJed Brown @*/
319815f1ad5SJed Brown PetscErrorCode TSSSPSetNumStages(TS ts,PetscInt nstages)
320815f1ad5SJed Brown {
321815f1ad5SJed Brown   PetscErrorCode ierr;
322815f1ad5SJed Brown 
323815f1ad5SJed Brown   PetscFunctionBegin;
324815f1ad5SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
325815f1ad5SJed Brown   ierr = PetscTryMethod(ts,"TSSSPSetNumStages_C",(TS,PetscInt),(ts,nstages));CHKERRQ(ierr);
326815f1ad5SJed Brown   PetscFunctionReturn(0);
327815f1ad5SJed Brown }
328815f1ad5SJed Brown 
329815f1ad5SJed Brown /*@
330815f1ad5SJed Brown    TSSSPGetNumStages - get the number of stages in the SSP time integration scheme
331815f1ad5SJed Brown 
332815f1ad5SJed Brown    Logically Collective
333815f1ad5SJed Brown 
3344165533cSJose E. Roman    Input Parameter:
3354165533cSJose E. Roman .  ts - time stepping object
336815f1ad5SJed Brown 
3374165533cSJose E. Roman    Output Parameter:
3384165533cSJose E. Roman .  nstages - number of stages
339815f1ad5SJed Brown 
340815f1ad5SJed Brown    Level: beginner
341815f1ad5SJed Brown 
342815f1ad5SJed Brown .seealso: TSSSP, TSSSPGetType(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104
343815f1ad5SJed Brown @*/
344815f1ad5SJed Brown PetscErrorCode TSSSPGetNumStages(TS ts,PetscInt *nstages)
345815f1ad5SJed Brown {
346815f1ad5SJed Brown   PetscErrorCode ierr;
347815f1ad5SJed Brown 
348815f1ad5SJed Brown   PetscFunctionBegin;
349815f1ad5SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
350163d334eSBarry Smith   ierr = PetscUseMethod(ts,"TSSSPGetNumStages_C",(TS,PetscInt*),(ts,nstages));CHKERRQ(ierr);
351815f1ad5SJed Brown   PetscFunctionReturn(0);
352815f1ad5SJed Brown }
353815f1ad5SJed Brown 
354cc2e6a90SBarry Smith static PetscErrorCode TSSSPSetType_SSP(TS ts,TSSSPType type)
355ef7bb5aaSJed Brown {
356ef7bb5aaSJed Brown   PetscErrorCode ierr,(*r)(TS,PetscReal,PetscReal,Vec);
357ef7bb5aaSJed Brown   TS_SSP         *ssp = (TS_SSP*)ts->data;
358ef7bb5aaSJed Brown 
359ef7bb5aaSJed Brown   PetscFunctionBegin;
3601c9cd337SJed Brown   ierr = PetscFunctionListFind(TSSSPList,type,&r);CHKERRQ(ierr);
361*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TS_SSP type %s given",type);
362ef7bb5aaSJed Brown   ssp->onestep = r;
3635164425aSJed Brown   ierr = PetscFree(ssp->type_name);CHKERRQ(ierr);
3645164425aSJed Brown   ierr = PetscStrallocpy(type,&ssp->type_name);CHKERRQ(ierr);
3652ffb9264SLisandro Dalcin   ts->default_adapt_type = TSADAPTNONE;
366ef7bb5aaSJed Brown   PetscFunctionReturn(0);
367ef7bb5aaSJed Brown }
368cc2e6a90SBarry Smith static PetscErrorCode TSSSPGetType_SSP(TS ts,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 }
376cc2e6a90SBarry Smith static PetscErrorCode TSSSPSetNumStages_SSP(TS ts,PetscInt nstages)
377815f1ad5SJed Brown {
378815f1ad5SJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
379815f1ad5SJed Brown 
380815f1ad5SJed Brown   PetscFunctionBegin;
381815f1ad5SJed Brown   ssp->nstages = nstages;
382815f1ad5SJed Brown   PetscFunctionReturn(0);
383815f1ad5SJed Brown }
384cc2e6a90SBarry Smith static PetscErrorCode TSSSPGetNumStages_SSP(TS ts,PetscInt *nstages)
385815f1ad5SJed Brown {
386815f1ad5SJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
387815f1ad5SJed Brown 
388815f1ad5SJed Brown   PetscFunctionBegin;
389815f1ad5SJed Brown   *nstages = ssp->nstages;
390815f1ad5SJed Brown   PetscFunctionReturn(0);
391815f1ad5SJed Brown }
392ef7bb5aaSJed Brown 
3934416b707SBarry Smith static PetscErrorCode TSSetFromOptions_SSP(PetscOptionItems *PetscOptionsObject,TS ts)
394ef7bb5aaSJed Brown {
395ef7bb5aaSJed Brown   char           tname[256] = TSSSPRKS2;
396ef7bb5aaSJed Brown   TS_SSP         *ssp       = (TS_SSP*)ts->data;
397ef7bb5aaSJed Brown   PetscErrorCode ierr;
398ace3abfcSBarry Smith   PetscBool      flg;
399ef7bb5aaSJed Brown 
400ef7bb5aaSJed Brown   PetscFunctionBegin;
401e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SSP ODE solver options");CHKERRQ(ierr);
402ef7bb5aaSJed Brown   {
403a264d7a6SBarry Smith     ierr = PetscOptionsFList("-ts_ssp_type","Type of SSP method","TSSSPSetType",TSSSPList,tname,tname,sizeof(tname),&flg);CHKERRQ(ierr);
404ef7bb5aaSJed Brown     if (flg) {
405ef7bb5aaSJed Brown       ierr = TSSSPSetType(ts,tname);CHKERRQ(ierr);
406ef7bb5aaSJed Brown     }
4070298fd71SBarry Smith     ierr = PetscOptionsInt("-ts_ssp_nstages","Number of stages","TSSSPSetNumStages",ssp->nstages,&ssp->nstages,NULL);CHKERRQ(ierr);
408ef7bb5aaSJed Brown   }
409ef7bb5aaSJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
410ef7bb5aaSJed Brown   PetscFunctionReturn(0);
411ef7bb5aaSJed Brown }
412ef7bb5aaSJed Brown 
413ef7bb5aaSJed Brown static PetscErrorCode TSView_SSP(TS ts,PetscViewer viewer)
414ef7bb5aaSJed Brown {
4151566a47fSLisandro Dalcin   TS_SSP         *ssp = (TS_SSP*)ts->data;
4161566a47fSLisandro Dalcin   PetscBool      ascii;
4171566a47fSLisandro Dalcin   PetscErrorCode ierr;
4181566a47fSLisandro Dalcin 
419ef7bb5aaSJed Brown   PetscFunctionBegin;
4201566a47fSLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&ascii);CHKERRQ(ierr);
4211566a47fSLisandro Dalcin   if (ascii) {ierr = PetscViewerASCIIPrintf(viewer,"  Scheme: %s\n",ssp->type_name);CHKERRQ(ierr);}
422ef7bb5aaSJed Brown   PetscFunctionReturn(0);
423ef7bb5aaSJed Brown }
424ef7bb5aaSJed Brown 
425ef7bb5aaSJed Brown /* ------------------------------------------------------------ */
426ef7bb5aaSJed Brown 
427ef7bb5aaSJed Brown /*MC
4288ab3e0fcSJed Brown       TSSSP - Explicit strong stability preserving ODE solver
4298ab3e0fcSJed Brown 
4308ab3e0fcSJed Brown   Most hyperbolic conservation laws have exact solutions that are total variation diminishing (TVD) or total variation
4318ab3e0fcSJed Brown   bounded (TVB) although these solutions often contain discontinuities.  Spatial discretizations such as Godunov's
4328ab3e0fcSJed Brown   scheme and high-resolution finite volume methods (TVD limiters, ENO/WENO) are designed to preserve these properties,
4338ab3e0fcSJed Brown   but they are usually formulated using a forward Euler time discretization or by coupling the space and time
4348ab3e0fcSJed Brown   discretization as in the classical Lax-Wendroff scheme.  When the space and time discretization is coupled, it is very
4358ab3e0fcSJed Brown   difficult to produce schemes with high temporal accuracy while preserving TVD properties.  An alternative is the
4368ab3e0fcSJed Brown   semidiscrete formulation where we choose a spatial discretization that is TVD with forward Euler and then choose a
4378ab3e0fcSJed Brown   time discretization that preserves the TVD property.  Such integrators are called strong stability preserving (SSP).
4388ab3e0fcSJed Brown 
4398ab3e0fcSJed Brown   Let c_eff be the minimum number of function evaluations required to step as far as one step of forward Euler while
4408ab3e0fcSJed Brown   still being SSP.  Some theoretical bounds
4418ab3e0fcSJed Brown 
4428ab3e0fcSJed Brown   1. There are no explicit methods with c_eff > 1.
4430085e20eSJed Brown 
4448ab3e0fcSJed Brown   2. There are no explicit methods beyond order 4 (for nonlinear problems) and c_eff > 0.
4450085e20eSJed Brown 
4468ab3e0fcSJed Brown   3. There are no implicit methods with order greater than 1 and c_eff > 2.
4478ab3e0fcSJed Brown 
4488ab3e0fcSJed Brown   This integrator provides Runge-Kutta methods of order 2, 3, and 4 with maximal values of c_eff.  More stages allows
4498ab3e0fcSJed Brown   for larger values of c_eff which improves efficiency.  These implementations are low-memory and only use 2 or 3 work
4508ab3e0fcSJed Brown   vectors regardless of the total number of stages, so e.g. 25-stage 3rd order methods may be an excellent choice.
4518ab3e0fcSJed Brown 
4528ab3e0fcSJed Brown   Methods can be chosen with -ts_ssp_type {rks2,rks3,rk104}
4538ab3e0fcSJed Brown 
4548ab3e0fcSJed Brown   rks2: Second order methods with any number s>1 of stages.  c_eff = (s-1)/s
4558ab3e0fcSJed Brown 
4568ab3e0fcSJed Brown   rks3: Third order methods with s=n^2 stages, n>1.  c_eff = (s-n)/s
4578ab3e0fcSJed Brown 
4588ab3e0fcSJed Brown   rk104: A 10-stage fourth order method.  c_eff = 0.6
459ef7bb5aaSJed Brown 
460ef7bb5aaSJed Brown   Level: beginner
461ef7bb5aaSJed Brown 
4627b6bb2c6SJed Brown   References:
46396a0c994SBarry Smith +  1. - Ketcheson, Highly efficient strong stability preserving Runge Kutta methods with low storage implementations, SISC, 2008.
46496a0c994SBarry Smith -  2. - Gottlieb, Ketcheson, and Shu, High order strong stability preserving time discretizations, J Scientific Computing, 2009.
4657b6bb2c6SJed Brown 
466ef7bb5aaSJed Brown .seealso:  TSCreate(), TS, TSSetType()
467ef7bb5aaSJed Brown 
468ef7bb5aaSJed Brown M*/
4698cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_SSP(TS ts)
470ef7bb5aaSJed Brown {
471ef7bb5aaSJed Brown   TS_SSP         *ssp;
472ef7bb5aaSJed Brown   PetscErrorCode ierr;
473ef7bb5aaSJed Brown 
474ef7bb5aaSJed Brown   PetscFunctionBegin;
475787849ffSJed Brown   ierr = TSSSPInitializePackage();CHKERRQ(ierr);
476ef7bb5aaSJed Brown 
477ef7bb5aaSJed Brown   ts->ops->setup          = TSSetUp_SSP;
478ef7bb5aaSJed Brown   ts->ops->step           = TSStep_SSP;
479277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_SSP;
480ef7bb5aaSJed Brown   ts->ops->destroy        = TSDestroy_SSP;
481ef7bb5aaSJed Brown   ts->ops->setfromoptions = TSSetFromOptions_SSP;
482ef7bb5aaSJed Brown   ts->ops->view           = TSView_SSP;
483ef7bb5aaSJed Brown 
484b00a9115SJed Brown   ierr = PetscNewLog(ts,&ssp);CHKERRQ(ierr);
485ef7bb5aaSJed Brown   ts->data = (void*)ssp;
486ef7bb5aaSJed Brown 
487bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetType_C",TSSSPGetType_SSP);CHKERRQ(ierr);
488bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetType_C",TSSSPSetType_SSP);CHKERRQ(ierr);
489bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetNumStages_C",TSSSPGetNumStages_SSP);CHKERRQ(ierr);
490bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetNumStages_C",TSSSPSetNumStages_SSP);CHKERRQ(ierr);
491815f1ad5SJed Brown 
492ef7bb5aaSJed Brown   ierr = TSSSPSetType(ts,TSSSPRKS2);CHKERRQ(ierr);
493ef7bb5aaSJed Brown   ssp->nstages = 5;
494ef7bb5aaSJed Brown   PetscFunctionReturn(0);
495ef7bb5aaSJed Brown }
496787849ffSJed Brown 
497787849ffSJed Brown /*@C
498787849ffSJed Brown   TSSSPInitializePackage - This function initializes everything in the TSSSP package. It is called
4998a690491SBarry Smith   from TSInitializePackage().
500787849ffSJed Brown 
501787849ffSJed Brown   Level: developer
502787849ffSJed Brown 
503787849ffSJed Brown .seealso: PetscInitialize()
504787849ffSJed Brown @*/
505787849ffSJed Brown PetscErrorCode TSSSPInitializePackage(void)
506787849ffSJed Brown {
507787849ffSJed Brown   PetscErrorCode ierr;
508787849ffSJed Brown 
509787849ffSJed Brown   PetscFunctionBegin;
510787849ffSJed Brown   if (TSSSPPackageInitialized) PetscFunctionReturn(0);
511787849ffSJed Brown   TSSSPPackageInitialized = PETSC_TRUE;
512787849ffSJed Brown   ierr = PetscFunctionListAdd(&TSSSPList,TSSSPRKS2, TSSSPStep_RK_2);CHKERRQ(ierr);
513787849ffSJed Brown   ierr = PetscFunctionListAdd(&TSSSPList,TSSSPRKS3, TSSSPStep_RK_3);CHKERRQ(ierr);
514787849ffSJed Brown   ierr = PetscFunctionListAdd(&TSSSPList,TSSSPRK104,TSSSPStep_RK_10_4);CHKERRQ(ierr);
515787849ffSJed Brown   ierr = PetscRegisterFinalize(TSSSPFinalizePackage);CHKERRQ(ierr);
516787849ffSJed Brown   PetscFunctionReturn(0);
517787849ffSJed Brown }
518787849ffSJed Brown 
519787849ffSJed Brown /*@C
520787849ffSJed Brown   TSSSPFinalizePackage - This function destroys everything in the TSSSP package. It is
521787849ffSJed Brown   called from PetscFinalize().
522787849ffSJed Brown 
523787849ffSJed Brown   Level: developer
524787849ffSJed Brown 
525787849ffSJed Brown .seealso: PetscFinalize()
526787849ffSJed Brown @*/
527787849ffSJed Brown PetscErrorCode TSSSPFinalizePackage(void)
528787849ffSJed Brown {
529787849ffSJed Brown   PetscErrorCode ierr;
530787849ffSJed Brown 
531787849ffSJed Brown   PetscFunctionBegin;
532787849ffSJed Brown   TSSSPPackageInitialized = PETSC_FALSE;
533787849ffSJed Brown   ierr = PetscFunctionListDestroy(&TSSSPList);CHKERRQ(ierr);
534787849ffSJed Brown   PetscFunctionReturn(0);
535787849ffSJed Brown }
536