xref: /petsc/src/ts/impls/explicit/ssp/ssp.c (revision a1f556f7e8f03e34584e3a8c859fc94ee1fe1aad)
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 
22ef7bb5aaSJed Brown   PetscFunctionBegin;
235f80ce2aSJacob Faibussowitsch   PetscCheck(!ssp->workout,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Work vectors already gotten");
24ef7bb5aaSJed Brown   if (ssp->nwork < n) {
25ef7bb5aaSJed Brown     if (ssp->nwork > 0) {
269566063dSJacob Faibussowitsch       PetscCall(VecDestroyVecs(ssp->nwork,&ssp->work));
27ef7bb5aaSJed Brown     }
289566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vec_sol,n,&ssp->work));
29ef7bb5aaSJed Brown     ssp->nwork = n;
30ef7bb5aaSJed Brown   }
31ef7bb5aaSJed Brown   *work = ssp->work;
32ef7bb5aaSJed Brown   ssp->workout = PETSC_TRUE;
33ef7bb5aaSJed Brown   PetscFunctionReturn(0);
34ef7bb5aaSJed Brown }
35ef7bb5aaSJed Brown 
3605d28066SBarry Smith static PetscErrorCode TSSSPRestoreWorkVectors(TS ts,PetscInt n,Vec **work)
37ef7bb5aaSJed Brown {
38ef7bb5aaSJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
39ef7bb5aaSJed Brown 
40ef7bb5aaSJed Brown   PetscFunctionBegin;
413c633725SBarry Smith   PetscCheck(ssp->workout,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Work vectors have not been gotten");
4208401ef6SPierre Jolivet   PetscCheck(*work == ssp->work,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong work vectors checked out");
43ef7bb5aaSJed Brown   ssp->workout = PETSC_FALSE;
440298fd71SBarry Smith   *work = NULL;
45ef7bb5aaSJed Brown   PetscFunctionReturn(0);
46ef7bb5aaSJed Brown }
47ef7bb5aaSJed Brown 
48815f1ad5SJed Brown /*MC
49815f1ad5SJed Brown    TSSSPRKS2 - Optimal second order SSP Runge-Kutta method, low-storage, c_eff=(s-1)/s
50815f1ad5SJed Brown 
51815f1ad5SJed Brown    Pseudocode 2 of Ketcheson 2008
52815f1ad5SJed Brown 
53b330ce4dSSatish Balay    Level: beginner
54b330ce4dSSatish Balay 
55db781477SPatrick Sanan .seealso: `TSSSP`, `TSSSPSetType()`, `TSSSPSetNumStages()`
56815f1ad5SJed Brown M*/
5705d28066SBarry Smith static PetscErrorCode TSSSPStep_RK_2(TS ts,PetscReal t0,PetscReal dt,Vec sol)
58ef7bb5aaSJed Brown {
59ef7bb5aaSJed Brown   TS_SSP         *ssp = (TS_SSP*)ts->data;
60ef7bb5aaSJed Brown   Vec            *work,F;
61ef7bb5aaSJed Brown   PetscInt       i,s;
62ef7bb5aaSJed Brown 
63ef7bb5aaSJed Brown   PetscFunctionBegin;
64ef7bb5aaSJed Brown   s    = ssp->nstages;
659566063dSJacob Faibussowitsch   PetscCall(TSSSPGetWorkVectors(ts,2,&work));
66ef7bb5aaSJed Brown   F    = work[1];
679566063dSJacob Faibussowitsch   PetscCall(VecCopy(sol,work[0]));
68ef7bb5aaSJed Brown   for (i=0; i<s-1; i++) {
69b8123daeSJed Brown     PetscReal stage_time = t0+dt*(i/(s-1.));
709566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts,stage_time));
719566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(ts,stage_time,work[0],F));
729566063dSJacob Faibussowitsch     PetscCall(VecAXPY(work[0],dt/(s-1.),F));
73ef7bb5aaSJed Brown   }
749566063dSJacob Faibussowitsch   PetscCall(TSComputeRHSFunction(ts,t0+dt,work[0],F));
759566063dSJacob Faibussowitsch   PetscCall(VecAXPBYPCZ(sol,(s-1.)/s,dt/s,1./s,work[0],F));
769566063dSJacob Faibussowitsch   PetscCall(TSSSPRestoreWorkVectors(ts,2,&work));
77ef7bb5aaSJed Brown   PetscFunctionReturn(0);
78ef7bb5aaSJed Brown }
79ef7bb5aaSJed Brown 
80815f1ad5SJed Brown /*MC
81aaf9cf16SJed Brown    TSSSPRKS3 - Optimal third order SSP Runge-Kutta, low-storage, c_eff=(PetscSqrtReal(s)-1)/PetscSqrtReal(s), where PetscSqrtReal(s) is an integer
82815f1ad5SJed Brown 
83815f1ad5SJed Brown    Pseudocode 2 of Ketcheson 2008
84815f1ad5SJed Brown 
85b330ce4dSSatish Balay    Level: beginner
86b330ce4dSSatish Balay 
87db781477SPatrick Sanan .seealso: `TSSSP`, `TSSSPSetType()`, `TSSSPSetNumStages()`
88815f1ad5SJed Brown M*/
8905d28066SBarry Smith static PetscErrorCode TSSSPStep_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;
94b8123daeSJed Brown   PetscReal      c,stage_time;
95ef7bb5aaSJed Brown 
96ef7bb5aaSJed Brown   PetscFunctionBegin;
97ef7bb5aaSJed Brown   s = ssp->nstages;
98aaf9cf16SJed Brown   n = (PetscInt)(PetscSqrtReal((PetscReal)s)+0.001);
99ef7bb5aaSJed Brown   r = s-n;
10063a3b9bcSJacob Faibussowitsch   PetscCheck(n*n == s,PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for optimal third order schemes with %" PetscInt_FMT " stages, must be a square number at least 4",s);
1019566063dSJacob Faibussowitsch   PetscCall(TSSSPGetWorkVectors(ts,3,&work));
102ef7bb5aaSJed Brown   F    = work[2];
1039566063dSJacob Faibussowitsch   PetscCall(VecCopy(sol,work[0]));
104ef7bb5aaSJed Brown   for (i=0; i<(n-1)*(n-2)/2; i++) {
105ef7bb5aaSJed Brown     c          = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
106b8123daeSJed Brown     stage_time = t0+c*dt;
1079566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts,stage_time));
1089566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(ts,stage_time,work[0],F));
1099566063dSJacob Faibussowitsch     PetscCall(VecAXPY(work[0],dt/r,F));
110ef7bb5aaSJed Brown   }
1119566063dSJacob Faibussowitsch   PetscCall(VecCopy(work[0],work[1]));
112ef7bb5aaSJed Brown   for (; i<n*(n+1)/2-1; i++) {
113ef7bb5aaSJed Brown     c          = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
114b8123daeSJed Brown     stage_time = t0+c*dt;
1159566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts,stage_time));
1169566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(ts,stage_time,work[0],F));
1179566063dSJacob Faibussowitsch     PetscCall(VecAXPY(work[0],dt/r,F));
118ef7bb5aaSJed Brown   }
119ef7bb5aaSJed Brown   {
120ef7bb5aaSJed Brown     c          = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
121b8123daeSJed Brown     stage_time = t0+c*dt;
1229566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts,stage_time));
1239566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(ts,stage_time,work[0],F));
1249566063dSJacob Faibussowitsch     PetscCall(VecAXPBYPCZ(work[0],1.*n/(2*n-1.),(n-1.)*dt/(r*(2*n-1)),(n-1.)/(2*n-1.),work[1],F));
125ef7bb5aaSJed Brown     i++;
126ef7bb5aaSJed Brown   }
127ef7bb5aaSJed Brown   for (; i<s; i++) {
128ef7bb5aaSJed Brown     c          = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
129b8123daeSJed Brown     stage_time = t0+c*dt;
1309566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts,stage_time));
1319566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(ts,stage_time,work[0],F));
1329566063dSJacob Faibussowitsch     PetscCall(VecAXPY(work[0],dt/r,F));
133ef7bb5aaSJed Brown   }
1349566063dSJacob Faibussowitsch   PetscCall(VecCopy(work[0],sol));
1359566063dSJacob Faibussowitsch   PetscCall(TSSSPRestoreWorkVectors(ts,3,&work));
136ef7bb5aaSJed Brown   PetscFunctionReturn(0);
137ef7bb5aaSJed Brown }
138ef7bb5aaSJed Brown 
139815f1ad5SJed Brown /*MC
140b330ce4dSSatish Balay    TSSSPRKS104 - Optimal fourth order SSP Runge-Kutta, low-storage (2N), c_eff=0.6
141815f1ad5SJed Brown 
142815f1ad5SJed Brown    SSPRK(10,4), Pseudocode 3 of Ketcheson 2008
143815f1ad5SJed Brown 
144b330ce4dSSatish Balay    Level: beginner
145b330ce4dSSatish Balay 
146db781477SPatrick Sanan .seealso: `TSSSP`, `TSSSPSetType()`
147815f1ad5SJed Brown M*/
14805d28066SBarry Smith static PetscErrorCode TSSSPStep_RK_10_4(TS ts,PetscReal t0,PetscReal dt,Vec sol)
149ef7bb5aaSJed Brown {
150ef7bb5aaSJed Brown   const PetscReal c[10] = {0, 1./6, 2./6, 3./6, 4./6, 2./6, 3./6, 4./6, 5./6, 1};
151ef7bb5aaSJed Brown   Vec             *work,F;
1525a586d82SBarry Smith   PetscInt        i;
153b8123daeSJed Brown   PetscReal       stage_time;
154ef7bb5aaSJed Brown 
155ef7bb5aaSJed Brown   PetscFunctionBegin;
1569566063dSJacob Faibussowitsch   PetscCall(TSSSPGetWorkVectors(ts,3,&work));
157ef7bb5aaSJed Brown   F    = work[2];
1589566063dSJacob Faibussowitsch   PetscCall(VecCopy(sol,work[0]));
159ef7bb5aaSJed Brown   for (i=0; i<5; i++) {
160b8123daeSJed Brown     stage_time = t0+c[i]*dt;
1619566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts,stage_time));
1629566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(ts,stage_time,work[0],F));
1639566063dSJacob Faibussowitsch     PetscCall(VecAXPY(work[0],dt/6,F));
164ef7bb5aaSJed Brown   }
1659566063dSJacob Faibussowitsch   PetscCall(VecAXPBYPCZ(work[1],1./25,9./25,0,sol,work[0]));
1669566063dSJacob Faibussowitsch   PetscCall(VecAXPBY(work[0],15,-5,work[1]));
167ef7bb5aaSJed Brown   for (; i<9; i++) {
168b8123daeSJed Brown     stage_time = t0+c[i]*dt;
1699566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts,stage_time));
1709566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(ts,stage_time,work[0],F));
1719566063dSJacob Faibussowitsch     PetscCall(VecAXPY(work[0],dt/6,F));
172ef7bb5aaSJed Brown   }
173b8123daeSJed Brown   stage_time = t0+dt;
1749566063dSJacob Faibussowitsch   PetscCall(TSPreStage(ts,stage_time));
1759566063dSJacob Faibussowitsch   PetscCall(TSComputeRHSFunction(ts,stage_time,work[0],F));
1769566063dSJacob Faibussowitsch   PetscCall(VecAXPBYPCZ(work[1],3./5,dt/10,1,work[0],F));
1779566063dSJacob Faibussowitsch   PetscCall(VecCopy(work[1],sol));
1789566063dSJacob Faibussowitsch   PetscCall(TSSSPRestoreWorkVectors(ts,3,&work));
179ef7bb5aaSJed Brown   PetscFunctionReturn(0);
180ef7bb5aaSJed Brown }
181ef7bb5aaSJed Brown 
182ef7bb5aaSJed Brown static PetscErrorCode TSSetUp_SSP(TS ts)
183ef7bb5aaSJed Brown {
184ef7bb5aaSJed Brown   PetscFunctionBegin;
1859566063dSJacob Faibussowitsch   PetscCall(TSCheckImplicitTerm(ts));
1869566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts,&ts->adapt));
1879566063dSJacob Faibussowitsch   PetscCall(TSAdaptCandidatesClear(ts->adapt));
188ef7bb5aaSJed Brown   PetscFunctionReturn(0);
189ef7bb5aaSJed Brown }
190ef7bb5aaSJed Brown 
191193ac0bcSJed Brown static PetscErrorCode TSStep_SSP(TS ts)
192ef7bb5aaSJed Brown {
193ef7bb5aaSJed Brown   TS_SSP         *ssp = (TS_SSP*)ts->data;
194ef7bb5aaSJed Brown   Vec            sol  = ts->vec_sol;
1951566a47fSLisandro Dalcin   PetscBool      stageok,accept = PETSC_TRUE;
1961566a47fSLisandro Dalcin   PetscReal      next_time_step = ts->time_step;
197ef7bb5aaSJed Brown 
198ef7bb5aaSJed Brown   PetscFunctionBegin;
1999566063dSJacob Faibussowitsch   PetscCall((*ssp->onestep)(ts,ts->ptime,ts->time_step,sol));
2009566063dSJacob Faibussowitsch   PetscCall(TSPostStage(ts,ts->ptime,0,&sol));
2019566063dSJacob Faibussowitsch   PetscCall(TSAdaptCheckStage(ts->adapt,ts,ts->ptime+ts->time_step,sol,&stageok));
2021566a47fSLisandro Dalcin   if (!stageok) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);}
2031566a47fSLisandro Dalcin 
2049566063dSJacob Faibussowitsch   PetscCall(TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept));
2051566a47fSLisandro Dalcin   if (!accept) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);}
2061566a47fSLisandro Dalcin 
207186e87acSLisandro Dalcin   ts->ptime += ts->time_step;
2081566a47fSLisandro Dalcin   ts->time_step = next_time_step;
209ef7bb5aaSJed Brown   PetscFunctionReturn(0);
210ef7bb5aaSJed Brown }
211ef7bb5aaSJed Brown /*------------------------------------------------------------*/
2121566a47fSLisandro Dalcin 
213277b19d0SLisandro Dalcin static PetscErrorCode TSReset_SSP(TS ts)
214277b19d0SLisandro Dalcin {
215277b19d0SLisandro Dalcin   TS_SSP         *ssp = (TS_SSP*)ts->data;
216277b19d0SLisandro Dalcin 
217277b19d0SLisandro Dalcin   PetscFunctionBegin;
2189566063dSJacob Faibussowitsch   if (ssp->work) PetscCall(VecDestroyVecs(ssp->nwork,&ssp->work));
219277b19d0SLisandro Dalcin   ssp->nwork   = 0;
220277b19d0SLisandro Dalcin   ssp->workout = PETSC_FALSE;
221277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
222277b19d0SLisandro Dalcin }
223277b19d0SLisandro Dalcin 
224ef7bb5aaSJed Brown static PetscErrorCode TSDestroy_SSP(TS ts)
225ef7bb5aaSJed Brown {
226815f1ad5SJed Brown   TS_SSP         *ssp = (TS_SSP*)ts->data;
227ef7bb5aaSJed Brown 
228ef7bb5aaSJed Brown   PetscFunctionBegin;
2299566063dSJacob Faibussowitsch   PetscCall(TSReset_SSP(ts));
2309566063dSJacob Faibussowitsch   PetscCall(PetscFree(ssp->type_name));
2319566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
2329566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetType_C",NULL));
2339566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetType_C",NULL));
2349566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetNumStages_C",NULL));
2359566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetNumStages_C",NULL));
236ef7bb5aaSJed Brown   PetscFunctionReturn(0);
237ef7bb5aaSJed Brown }
238ef7bb5aaSJed Brown /*------------------------------------------------------------*/
239ef7bb5aaSJed Brown 
240815f1ad5SJed Brown /*@C
241815f1ad5SJed Brown    TSSSPSetType - set the SSP time integration scheme to use
242815f1ad5SJed Brown 
243815f1ad5SJed Brown    Logically Collective
244815f1ad5SJed Brown 
2454165533cSJose E. Roman    Input Parameters:
2464165533cSJose E. Roman +  ts - time stepping object
2474165533cSJose E. Roman -  ssptype - 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
251*a1f556f7SAidan Hamilton    -ts_ssp_nstages<rks2: 5, rks3: 9> : Number of stages
252815f1ad5SJed Brown 
253815f1ad5SJed Brown    Level: beginner
254815f1ad5SJed Brown 
255db781477SPatrick Sanan .seealso: `TSSSP`, `TSSSPGetType()`, `TSSSPSetNumStages()`, `TSSSPRKS2`, `TSSSPRKS3`, `TSSSPRK104`
256815f1ad5SJed Brown @*/
257b92453a8SLisandro Dalcin PetscErrorCode TSSSPSetType(TS ts,TSSSPType ssptype)
258815f1ad5SJed Brown {
259815f1ad5SJed Brown   PetscFunctionBegin;
260815f1ad5SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
261b92453a8SLisandro Dalcin   PetscValidCharPointer(ssptype,2);
262cac4c232SBarry Smith   PetscTryMethod(ts,"TSSSPSetType_C",(TS,TSSSPType),(ts,ssptype));
263815f1ad5SJed Brown   PetscFunctionReturn(0);
264815f1ad5SJed Brown }
265815f1ad5SJed Brown 
266815f1ad5SJed Brown /*@C
267815f1ad5SJed Brown    TSSSPGetType - get the SSP time integration scheme
268815f1ad5SJed Brown 
269815f1ad5SJed Brown    Logically Collective
270815f1ad5SJed Brown 
2714165533cSJose E. Roman    Input Parameter:
2724165533cSJose E. Roman .  ts - time stepping object
273815f1ad5SJed Brown 
2744165533cSJose E. Roman    Output Parameter:
2754165533cSJose E. Roman .  type - type of scheme being used
276815f1ad5SJed Brown 
277815f1ad5SJed Brown    Level: beginner
278815f1ad5SJed Brown 
279db781477SPatrick Sanan .seealso: `TSSSP`, `TSSSPSettype()`, `TSSSPSetNumStages()`, `TSSSPRKS2`, `TSSSPRKS3`, `TSSSPRK104`
280815f1ad5SJed Brown @*/
28119fd82e9SBarry Smith PetscErrorCode TSSSPGetType(TS ts,TSSSPType *type)
282815f1ad5SJed Brown {
283815f1ad5SJed Brown   PetscFunctionBegin;
284815f1ad5SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
285cac4c232SBarry Smith   PetscUseMethod(ts,"TSSSPGetType_C",(TS,TSSSPType*),(ts,type));
286815f1ad5SJed Brown   PetscFunctionReturn(0);
287815f1ad5SJed Brown }
288815f1ad5SJed Brown 
289815f1ad5SJed Brown /*@
290*a1f556f7SAidan Hamilton    TSSSPSetNumStages - set the number of stages to use with the SSP method. Must be called after
291*a1f556f7SAidan Hamilton    TSSSPSetType().
292815f1ad5SJed Brown 
293815f1ad5SJed Brown    Logically Collective
294815f1ad5SJed Brown 
2954165533cSJose E. Roman    Input Parameters:
2964165533cSJose E. Roman +  ts - time stepping object
2974165533cSJose E. Roman -  nstages - number of stages
298815f1ad5SJed Brown 
299815f1ad5SJed Brown    Options Database Keys:
300*a1f556f7SAidan Hamilton    -ts_ssp_type <rks2>               : Type of SSP method (one of) rks2 rks3 rk104
301*a1f556f7SAidan Hamilton    -ts_ssp_nstages<rks2: 5, rks3: 9> : Number of stages
302815f1ad5SJed Brown 
303815f1ad5SJed Brown    Level: beginner
304815f1ad5SJed Brown 
305db781477SPatrick Sanan .seealso: `TSSSP`, `TSSSPGetNumStages()`, `TSSSPSetNumStages()`, `TSSSPRKS2`, `TSSSPRKS3`, `TSSSPRK104`
306815f1ad5SJed Brown @*/
307815f1ad5SJed Brown PetscErrorCode TSSSPSetNumStages(TS ts,PetscInt nstages)
308815f1ad5SJed Brown {
309815f1ad5SJed Brown   PetscFunctionBegin;
310815f1ad5SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
311cac4c232SBarry Smith   PetscTryMethod(ts,"TSSSPSetNumStages_C",(TS,PetscInt),(ts,nstages));
312815f1ad5SJed Brown   PetscFunctionReturn(0);
313815f1ad5SJed Brown }
314815f1ad5SJed Brown 
315815f1ad5SJed Brown /*@
316815f1ad5SJed Brown    TSSSPGetNumStages - get the number of stages in the SSP time integration scheme
317815f1ad5SJed Brown 
318815f1ad5SJed Brown    Logically Collective
319815f1ad5SJed Brown 
3204165533cSJose E. Roman    Input Parameter:
3214165533cSJose E. Roman .  ts - time stepping object
322815f1ad5SJed Brown 
3234165533cSJose E. Roman    Output Parameter:
3244165533cSJose E. Roman .  nstages - number of stages
325815f1ad5SJed Brown 
326815f1ad5SJed Brown    Level: beginner
327815f1ad5SJed Brown 
328db781477SPatrick Sanan .seealso: `TSSSP`, `TSSSPGetType()`, `TSSSPSetNumStages()`, `TSSSPRKS2`, `TSSSPRKS3`, `TSSSPRK104`
329815f1ad5SJed Brown @*/
330815f1ad5SJed Brown PetscErrorCode TSSSPGetNumStages(TS ts,PetscInt *nstages)
331815f1ad5SJed Brown {
332815f1ad5SJed Brown   PetscFunctionBegin;
333815f1ad5SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
334cac4c232SBarry Smith   PetscUseMethod(ts,"TSSSPGetNumStages_C",(TS,PetscInt*),(ts,nstages));
335815f1ad5SJed Brown   PetscFunctionReturn(0);
336815f1ad5SJed Brown }
337815f1ad5SJed Brown 
338cc2e6a90SBarry Smith static PetscErrorCode TSSSPSetType_SSP(TS ts,TSSSPType type)
339ef7bb5aaSJed Brown {
340ef7bb5aaSJed Brown   TS_SSP         *ssp = (TS_SSP*)ts->data;
3415f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(TS,PetscReal,PetscReal,Vec);
342*a1f556f7SAidan Hamilton   PetscBool      flag;
343ef7bb5aaSJed Brown 
344ef7bb5aaSJed Brown   PetscFunctionBegin;
3459566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(TSSSPList,type,&r));
3463c633725SBarry Smith   PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TS_SSP type %s given",type);
347ef7bb5aaSJed Brown   ssp->onestep = r;
3489566063dSJacob Faibussowitsch   PetscCall(PetscFree(ssp->type_name));
3499566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(type,&ssp->type_name));
3502ffb9264SLisandro Dalcin   ts->default_adapt_type = TSADAPTNONE;
351*a1f556f7SAidan Hamilton   PetscCall(PetscStrcmp(type,TSSSPRKS2,&flag));
352*a1f556f7SAidan Hamilton   if (flag) {
353*a1f556f7SAidan Hamilton     ssp->nstages = 5;
354*a1f556f7SAidan Hamilton   } else {
355*a1f556f7SAidan Hamilton     PetscCall(PetscStrcmp(type,TSSSPRKS3,&flag));
356*a1f556f7SAidan Hamilton     if (flag) {
357*a1f556f7SAidan Hamilton       ssp->nstages = 9;
358*a1f556f7SAidan Hamilton     } else {
359*a1f556f7SAidan Hamilton       ssp->nstages = 5;
360*a1f556f7SAidan Hamilton     }
361*a1f556f7SAidan Hamilton   }
362ef7bb5aaSJed Brown   PetscFunctionReturn(0);
363ef7bb5aaSJed Brown }
364cc2e6a90SBarry Smith static PetscErrorCode TSSSPGetType_SSP(TS ts,TSSSPType *type)
365815f1ad5SJed Brown {
366815f1ad5SJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
367815f1ad5SJed Brown 
368815f1ad5SJed Brown   PetscFunctionBegin;
3695164425aSJed Brown   *type = ssp->type_name;
370815f1ad5SJed Brown   PetscFunctionReturn(0);
371815f1ad5SJed Brown }
372cc2e6a90SBarry Smith static PetscErrorCode TSSSPSetNumStages_SSP(TS ts,PetscInt nstages)
373815f1ad5SJed Brown {
374815f1ad5SJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
375815f1ad5SJed Brown 
376815f1ad5SJed Brown   PetscFunctionBegin;
377815f1ad5SJed Brown   ssp->nstages = nstages;
378815f1ad5SJed Brown   PetscFunctionReturn(0);
379815f1ad5SJed Brown }
380cc2e6a90SBarry Smith static PetscErrorCode TSSSPGetNumStages_SSP(TS ts,PetscInt *nstages)
381815f1ad5SJed Brown {
382815f1ad5SJed Brown   TS_SSP *ssp = (TS_SSP*)ts->data;
383815f1ad5SJed Brown 
384815f1ad5SJed Brown   PetscFunctionBegin;
385815f1ad5SJed Brown   *nstages = ssp->nstages;
386815f1ad5SJed Brown   PetscFunctionReturn(0);
387815f1ad5SJed Brown }
388ef7bb5aaSJed Brown 
3894416b707SBarry Smith static PetscErrorCode TSSetFromOptions_SSP(PetscOptionItems *PetscOptionsObject,TS ts)
390ef7bb5aaSJed Brown {
391ef7bb5aaSJed Brown   char           tname[256] = TSSSPRKS2;
392ef7bb5aaSJed Brown   TS_SSP         *ssp       = (TS_SSP*)ts->data;
393ace3abfcSBarry Smith   PetscBool      flg;
394ef7bb5aaSJed Brown 
395ef7bb5aaSJed Brown   PetscFunctionBegin;
396d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"SSP ODE solver options");
397ef7bb5aaSJed Brown   {
3989566063dSJacob Faibussowitsch     PetscCall(PetscOptionsFList("-ts_ssp_type","Type of SSP method","TSSSPSetType",TSSSPList,tname,tname,sizeof(tname),&flg));
399ef7bb5aaSJed Brown     if (flg) {
4009566063dSJacob Faibussowitsch       PetscCall(TSSSPSetType(ts,tname));
401ef7bb5aaSJed Brown     }
4029566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-ts_ssp_nstages","Number of stages","TSSSPSetNumStages",ssp->nstages,&ssp->nstages,NULL));
403ef7bb5aaSJed Brown   }
404d0609cedSBarry Smith   PetscOptionsHeadEnd();
405ef7bb5aaSJed Brown   PetscFunctionReturn(0);
406ef7bb5aaSJed Brown }
407ef7bb5aaSJed Brown 
408ef7bb5aaSJed Brown static PetscErrorCode TSView_SSP(TS ts,PetscViewer viewer)
409ef7bb5aaSJed Brown {
4101566a47fSLisandro Dalcin   TS_SSP         *ssp = (TS_SSP*)ts->data;
4111566a47fSLisandro Dalcin   PetscBool      ascii;
4121566a47fSLisandro Dalcin 
413ef7bb5aaSJed Brown   PetscFunctionBegin;
4149566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&ascii));
4159566063dSJacob Faibussowitsch   if (ascii) PetscCall(PetscViewerASCIIPrintf(viewer,"  Scheme: %s\n",ssp->type_name));
416ef7bb5aaSJed Brown   PetscFunctionReturn(0);
417ef7bb5aaSJed Brown }
418ef7bb5aaSJed Brown 
419ef7bb5aaSJed Brown /* ------------------------------------------------------------ */
420ef7bb5aaSJed Brown 
421ef7bb5aaSJed Brown /*MC
4228ab3e0fcSJed Brown       TSSSP - Explicit strong stability preserving ODE solver
4238ab3e0fcSJed Brown 
4248ab3e0fcSJed Brown   Most hyperbolic conservation laws have exact solutions that are total variation diminishing (TVD) or total variation
4258ab3e0fcSJed Brown   bounded (TVB) although these solutions often contain discontinuities.  Spatial discretizations such as Godunov's
4268ab3e0fcSJed Brown   scheme and high-resolution finite volume methods (TVD limiters, ENO/WENO) are designed to preserve these properties,
4278ab3e0fcSJed Brown   but they are usually formulated using a forward Euler time discretization or by coupling the space and time
4288ab3e0fcSJed Brown   discretization as in the classical Lax-Wendroff scheme.  When the space and time discretization is coupled, it is very
4298ab3e0fcSJed Brown   difficult to produce schemes with high temporal accuracy while preserving TVD properties.  An alternative is the
4308ab3e0fcSJed Brown   semidiscrete formulation where we choose a spatial discretization that is TVD with forward Euler and then choose a
4318ab3e0fcSJed Brown   time discretization that preserves the TVD property.  Such integrators are called strong stability preserving (SSP).
4328ab3e0fcSJed Brown 
4338ab3e0fcSJed Brown   Let c_eff be the minimum number of function evaluations required to step as far as one step of forward Euler while
4348ab3e0fcSJed Brown   still being SSP.  Some theoretical bounds
4358ab3e0fcSJed Brown 
4368ab3e0fcSJed Brown   1. There are no explicit methods with c_eff > 1.
4370085e20eSJed Brown 
4388ab3e0fcSJed Brown   2. There are no explicit methods beyond order 4 (for nonlinear problems) and c_eff > 0.
4390085e20eSJed Brown 
4408ab3e0fcSJed Brown   3. There are no implicit methods with order greater than 1 and c_eff > 2.
4418ab3e0fcSJed Brown 
4428ab3e0fcSJed Brown   This integrator provides Runge-Kutta methods of order 2, 3, and 4 with maximal values of c_eff.  More stages allows
4438ab3e0fcSJed Brown   for larger values of c_eff which improves efficiency.  These implementations are low-memory and only use 2 or 3 work
4448ab3e0fcSJed Brown   vectors regardless of the total number of stages, so e.g. 25-stage 3rd order methods may be an excellent choice.
4458ab3e0fcSJed Brown 
4468ab3e0fcSJed Brown   Methods can be chosen with -ts_ssp_type {rks2,rks3,rk104}
4478ab3e0fcSJed Brown 
4488ab3e0fcSJed Brown   rks2: Second order methods with any number s>1 of stages.  c_eff = (s-1)/s
4498ab3e0fcSJed Brown 
4508ab3e0fcSJed Brown   rks3: Third order methods with s=n^2 stages, n>1.  c_eff = (s-n)/s
4518ab3e0fcSJed Brown 
4528ab3e0fcSJed Brown   rk104: A 10-stage fourth order method.  c_eff = 0.6
453ef7bb5aaSJed Brown 
454ef7bb5aaSJed Brown   Level: beginner
455ef7bb5aaSJed Brown 
4567b6bb2c6SJed Brown   References:
457606c0280SSatish Balay +  * - Ketcheson, Highly efficient strong stability preserving Runge Kutta methods with low storage implementations, SISC, 2008.
458606c0280SSatish Balay -  * - Gottlieb, Ketcheson, and Shu, High order strong stability preserving time discretizations, J Scientific Computing, 2009.
4597b6bb2c6SJed Brown 
460db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()`
461ef7bb5aaSJed Brown 
462ef7bb5aaSJed Brown M*/
4638cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_SSP(TS ts)
464ef7bb5aaSJed Brown {
465ef7bb5aaSJed Brown   TS_SSP         *ssp;
466ef7bb5aaSJed Brown 
467ef7bb5aaSJed Brown   PetscFunctionBegin;
4689566063dSJacob Faibussowitsch   PetscCall(TSSSPInitializePackage());
469ef7bb5aaSJed Brown 
470ef7bb5aaSJed Brown   ts->ops->setup          = TSSetUp_SSP;
471ef7bb5aaSJed Brown   ts->ops->step           = TSStep_SSP;
472277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_SSP;
473ef7bb5aaSJed Brown   ts->ops->destroy        = TSDestroy_SSP;
474ef7bb5aaSJed Brown   ts->ops->setfromoptions = TSSetFromOptions_SSP;
475ef7bb5aaSJed Brown   ts->ops->view           = TSView_SSP;
476ef7bb5aaSJed Brown 
4779566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(ts,&ssp));
478ef7bb5aaSJed Brown   ts->data = (void*)ssp;
479ef7bb5aaSJed Brown 
4809566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetType_C",TSSSPGetType_SSP));
4819566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetType_C",TSSSPSetType_SSP));
4829566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetNumStages_C",TSSSPGetNumStages_SSP));
4839566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetNumStages_C",TSSSPSetNumStages_SSP));
484815f1ad5SJed Brown 
4859566063dSJacob Faibussowitsch   PetscCall(TSSSPSetType(ts,TSSSPRKS2));
486ef7bb5aaSJed Brown   PetscFunctionReturn(0);
487ef7bb5aaSJed Brown }
488787849ffSJed Brown 
489787849ffSJed Brown /*@C
490787849ffSJed Brown   TSSSPInitializePackage - This function initializes everything in the TSSSP package. It is called
4918a690491SBarry Smith   from TSInitializePackage().
492787849ffSJed Brown 
493787849ffSJed Brown   Level: developer
494787849ffSJed Brown 
495db781477SPatrick Sanan .seealso: `PetscInitialize()`
496787849ffSJed Brown @*/
497787849ffSJed Brown PetscErrorCode TSSSPInitializePackage(void)
498787849ffSJed Brown {
499787849ffSJed Brown   PetscFunctionBegin;
500787849ffSJed Brown   if (TSSSPPackageInitialized) PetscFunctionReturn(0);
501787849ffSJed Brown   TSSSPPackageInitialized = PETSC_TRUE;
5029566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&TSSSPList,TSSSPRKS2, TSSSPStep_RK_2));
5039566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&TSSSPList,TSSSPRKS3, TSSSPStep_RK_3));
5049566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&TSSSPList,TSSSPRK104,TSSSPStep_RK_10_4));
5059566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(TSSSPFinalizePackage));
506787849ffSJed Brown   PetscFunctionReturn(0);
507787849ffSJed Brown }
508787849ffSJed Brown 
509787849ffSJed Brown /*@C
510787849ffSJed Brown   TSSSPFinalizePackage - This function destroys everything in the TSSSP package. It is
511787849ffSJed Brown   called from PetscFinalize().
512787849ffSJed Brown 
513787849ffSJed Brown   Level: developer
514787849ffSJed Brown 
515db781477SPatrick Sanan .seealso: `PetscFinalize()`
516787849ffSJed Brown @*/
517787849ffSJed Brown PetscErrorCode TSSSPFinalizePackage(void)
518787849ffSJed Brown {
519787849ffSJed Brown   PetscFunctionBegin;
520787849ffSJed Brown   TSSSPPackageInitialized = PETSC_FALSE;
5219566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&TSSSPList));
522787849ffSJed Brown   PetscFunctionReturn(0);
523787849ffSJed Brown }
524