xref: /petsc/src/ts/impls/explicit/ssp/ssp.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
189371c9d4SSatish Balay static PetscErrorCode TSSSPGetWorkVectors(TS ts, PetscInt n, Vec **work) {
19ef7bb5aaSJed Brown   TS_SSP *ssp = (TS_SSP *)ts->data;
20ef7bb5aaSJed Brown 
21ef7bb5aaSJed Brown   PetscFunctionBegin;
225f80ce2aSJacob Faibussowitsch   PetscCheck(!ssp->workout, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Work vectors already gotten");
23ef7bb5aaSJed Brown   if (ssp->nwork < n) {
24*48a46eb9SPierre Jolivet     if (ssp->nwork > 0) PetscCall(VecDestroyVecs(ssp->nwork, &ssp->work));
259566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vec_sol, n, &ssp->work));
26ef7bb5aaSJed Brown     ssp->nwork = n;
27ef7bb5aaSJed Brown   }
28ef7bb5aaSJed Brown   *work        = ssp->work;
29ef7bb5aaSJed Brown   ssp->workout = PETSC_TRUE;
30ef7bb5aaSJed Brown   PetscFunctionReturn(0);
31ef7bb5aaSJed Brown }
32ef7bb5aaSJed Brown 
339371c9d4SSatish Balay static PetscErrorCode TSSSPRestoreWorkVectors(TS ts, PetscInt n, Vec **work) {
34ef7bb5aaSJed Brown   TS_SSP *ssp = (TS_SSP *)ts->data;
35ef7bb5aaSJed Brown 
36ef7bb5aaSJed Brown   PetscFunctionBegin;
373c633725SBarry Smith   PetscCheck(ssp->workout, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Work vectors have not been gotten");
3808401ef6SPierre Jolivet   PetscCheck(*work == ssp->work, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Wrong work vectors checked out");
39ef7bb5aaSJed Brown   ssp->workout = PETSC_FALSE;
400298fd71SBarry Smith   *work        = NULL;
41ef7bb5aaSJed Brown   PetscFunctionReturn(0);
42ef7bb5aaSJed Brown }
43ef7bb5aaSJed Brown 
44815f1ad5SJed Brown /*MC
45815f1ad5SJed Brown    TSSSPRKS2 - Optimal second order SSP Runge-Kutta method, low-storage, c_eff=(s-1)/s
46815f1ad5SJed Brown 
47815f1ad5SJed Brown    Pseudocode 2 of Ketcheson 2008
48815f1ad5SJed Brown 
49b330ce4dSSatish Balay    Level: beginner
50b330ce4dSSatish Balay 
51db781477SPatrick Sanan .seealso: `TSSSP`, `TSSSPSetType()`, `TSSSPSetNumStages()`
52815f1ad5SJed Brown M*/
539371c9d4SSatish Balay static PetscErrorCode TSSSPStep_RK_2(TS ts, PetscReal t0, PetscReal dt, Vec sol) {
54ef7bb5aaSJed Brown   TS_SSP  *ssp = (TS_SSP *)ts->data;
55ef7bb5aaSJed Brown   Vec     *work, F;
56ef7bb5aaSJed Brown   PetscInt i, s;
57ef7bb5aaSJed Brown 
58ef7bb5aaSJed Brown   PetscFunctionBegin;
59ef7bb5aaSJed Brown   s = ssp->nstages;
609566063dSJacob Faibussowitsch   PetscCall(TSSSPGetWorkVectors(ts, 2, &work));
61ef7bb5aaSJed Brown   F = work[1];
629566063dSJacob Faibussowitsch   PetscCall(VecCopy(sol, work[0]));
63ef7bb5aaSJed Brown   for (i = 0; i < s - 1; i++) {
64b8123daeSJed Brown     PetscReal stage_time = t0 + dt * (i / (s - 1.));
659566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, stage_time));
669566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F));
679566063dSJacob Faibussowitsch     PetscCall(VecAXPY(work[0], dt / (s - 1.), F));
68ef7bb5aaSJed Brown   }
699566063dSJacob Faibussowitsch   PetscCall(TSComputeRHSFunction(ts, t0 + dt, work[0], F));
709566063dSJacob Faibussowitsch   PetscCall(VecAXPBYPCZ(sol, (s - 1.) / s, dt / s, 1. / s, work[0], F));
719566063dSJacob Faibussowitsch   PetscCall(TSSSPRestoreWorkVectors(ts, 2, &work));
72ef7bb5aaSJed Brown   PetscFunctionReturn(0);
73ef7bb5aaSJed Brown }
74ef7bb5aaSJed Brown 
75815f1ad5SJed Brown /*MC
76aaf9cf16SJed Brown    TSSSPRKS3 - Optimal third order SSP Runge-Kutta, low-storage, c_eff=(PetscSqrtReal(s)-1)/PetscSqrtReal(s), where PetscSqrtReal(s) is an integer
77815f1ad5SJed Brown 
78815f1ad5SJed Brown    Pseudocode 2 of Ketcheson 2008
79815f1ad5SJed Brown 
80b330ce4dSSatish Balay    Level: beginner
81b330ce4dSSatish Balay 
82db781477SPatrick Sanan .seealso: `TSSSP`, `TSSSPSetType()`, `TSSSPSetNumStages()`
83815f1ad5SJed Brown M*/
849371c9d4SSatish Balay static PetscErrorCode TSSSPStep_RK_3(TS ts, PetscReal t0, PetscReal dt, Vec sol) {
85ef7bb5aaSJed Brown   TS_SSP   *ssp = (TS_SSP *)ts->data;
86ef7bb5aaSJed Brown   Vec      *work, F;
87ef7bb5aaSJed Brown   PetscInt  i, s, n, r;
88b8123daeSJed Brown   PetscReal c, stage_time;
89ef7bb5aaSJed Brown 
90ef7bb5aaSJed Brown   PetscFunctionBegin;
91ef7bb5aaSJed Brown   s = ssp->nstages;
92aaf9cf16SJed Brown   n = (PetscInt)(PetscSqrtReal((PetscReal)s) + 0.001);
93ef7bb5aaSJed Brown   r = s - n;
9463a3b9bcSJacob 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);
959566063dSJacob Faibussowitsch   PetscCall(TSSSPGetWorkVectors(ts, 3, &work));
96ef7bb5aaSJed Brown   F = work[2];
979566063dSJacob Faibussowitsch   PetscCall(VecCopy(sol, work[0]));
98ef7bb5aaSJed Brown   for (i = 0; i < (n - 1) * (n - 2) / 2; i++) {
99ef7bb5aaSJed Brown     c          = (i < n * (n + 1) / 2) ? 1. * i / (s - n) : (1. * i - n) / (s - n);
100b8123daeSJed Brown     stage_time = t0 + c * dt;
1019566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, stage_time));
1029566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F));
1039566063dSJacob Faibussowitsch     PetscCall(VecAXPY(work[0], dt / r, F));
104ef7bb5aaSJed Brown   }
1059566063dSJacob Faibussowitsch   PetscCall(VecCopy(work[0], work[1]));
106ef7bb5aaSJed Brown   for (; i < n * (n + 1) / 2 - 1; i++) {
107ef7bb5aaSJed Brown     c          = (i < n * (n + 1) / 2) ? 1. * i / (s - n) : (1. * i - n) / (s - n);
108b8123daeSJed Brown     stage_time = t0 + c * dt;
1099566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, stage_time));
1109566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F));
1119566063dSJacob Faibussowitsch     PetscCall(VecAXPY(work[0], dt / r, F));
112ef7bb5aaSJed Brown   }
113ef7bb5aaSJed Brown   {
114ef7bb5aaSJed Brown     c          = (i < n * (n + 1) / 2) ? 1. * i / (s - n) : (1. * i - n) / (s - n);
115b8123daeSJed Brown     stage_time = t0 + c * dt;
1169566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, stage_time));
1179566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F));
1189566063dSJacob 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));
119ef7bb5aaSJed Brown     i++;
120ef7bb5aaSJed Brown   }
121ef7bb5aaSJed Brown   for (; i < s; i++) {
122ef7bb5aaSJed Brown     c          = (i < n * (n + 1) / 2) ? 1. * i / (s - n) : (1. * i - n) / (s - n);
123b8123daeSJed Brown     stage_time = t0 + c * dt;
1249566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, stage_time));
1259566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F));
1269566063dSJacob Faibussowitsch     PetscCall(VecAXPY(work[0], dt / r, F));
127ef7bb5aaSJed Brown   }
1289566063dSJacob Faibussowitsch   PetscCall(VecCopy(work[0], sol));
1299566063dSJacob Faibussowitsch   PetscCall(TSSSPRestoreWorkVectors(ts, 3, &work));
130ef7bb5aaSJed Brown   PetscFunctionReturn(0);
131ef7bb5aaSJed Brown }
132ef7bb5aaSJed Brown 
133815f1ad5SJed Brown /*MC
134b330ce4dSSatish Balay    TSSSPRKS104 - Optimal fourth order SSP Runge-Kutta, low-storage (2N), c_eff=0.6
135815f1ad5SJed Brown 
136815f1ad5SJed Brown    SSPRK(10,4), Pseudocode 3 of Ketcheson 2008
137815f1ad5SJed Brown 
138b330ce4dSSatish Balay    Level: beginner
139b330ce4dSSatish Balay 
140db781477SPatrick Sanan .seealso: `TSSSP`, `TSSSPSetType()`
141815f1ad5SJed Brown M*/
1429371c9d4SSatish Balay static PetscErrorCode TSSSPStep_RK_10_4(TS ts, PetscReal t0, PetscReal dt, Vec sol) {
143ef7bb5aaSJed Brown   const PetscReal c[10] = {0, 1. / 6, 2. / 6, 3. / 6, 4. / 6, 2. / 6, 3. / 6, 4. / 6, 5. / 6, 1};
144ef7bb5aaSJed Brown   Vec            *work, F;
1455a586d82SBarry Smith   PetscInt        i;
146b8123daeSJed Brown   PetscReal       stage_time;
147ef7bb5aaSJed Brown 
148ef7bb5aaSJed Brown   PetscFunctionBegin;
1499566063dSJacob Faibussowitsch   PetscCall(TSSSPGetWorkVectors(ts, 3, &work));
150ef7bb5aaSJed Brown   F = work[2];
1519566063dSJacob Faibussowitsch   PetscCall(VecCopy(sol, work[0]));
152ef7bb5aaSJed Brown   for (i = 0; i < 5; i++) {
153b8123daeSJed Brown     stage_time = t0 + c[i] * dt;
1549566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, stage_time));
1559566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F));
1569566063dSJacob Faibussowitsch     PetscCall(VecAXPY(work[0], dt / 6, F));
157ef7bb5aaSJed Brown   }
1589566063dSJacob Faibussowitsch   PetscCall(VecAXPBYPCZ(work[1], 1. / 25, 9. / 25, 0, sol, work[0]));
1599566063dSJacob Faibussowitsch   PetscCall(VecAXPBY(work[0], 15, -5, work[1]));
160ef7bb5aaSJed Brown   for (; i < 9; i++) {
161b8123daeSJed Brown     stage_time = t0 + c[i] * dt;
1629566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, stage_time));
1639566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F));
1649566063dSJacob Faibussowitsch     PetscCall(VecAXPY(work[0], dt / 6, F));
165ef7bb5aaSJed Brown   }
166b8123daeSJed Brown   stage_time = t0 + dt;
1679566063dSJacob Faibussowitsch   PetscCall(TSPreStage(ts, stage_time));
1689566063dSJacob Faibussowitsch   PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F));
1699566063dSJacob Faibussowitsch   PetscCall(VecAXPBYPCZ(work[1], 3. / 5, dt / 10, 1, work[0], F));
1709566063dSJacob Faibussowitsch   PetscCall(VecCopy(work[1], sol));
1719566063dSJacob Faibussowitsch   PetscCall(TSSSPRestoreWorkVectors(ts, 3, &work));
172ef7bb5aaSJed Brown   PetscFunctionReturn(0);
173ef7bb5aaSJed Brown }
174ef7bb5aaSJed Brown 
1759371c9d4SSatish Balay static PetscErrorCode TSSetUp_SSP(TS ts) {
176ef7bb5aaSJed Brown   PetscFunctionBegin;
1779566063dSJacob Faibussowitsch   PetscCall(TSCheckImplicitTerm(ts));
1789566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &ts->adapt));
1799566063dSJacob Faibussowitsch   PetscCall(TSAdaptCandidatesClear(ts->adapt));
180ef7bb5aaSJed Brown   PetscFunctionReturn(0);
181ef7bb5aaSJed Brown }
182ef7bb5aaSJed Brown 
1839371c9d4SSatish Balay static PetscErrorCode TSStep_SSP(TS ts) {
184ef7bb5aaSJed Brown   TS_SSP   *ssp = (TS_SSP *)ts->data;
185ef7bb5aaSJed Brown   Vec       sol = ts->vec_sol;
1861566a47fSLisandro Dalcin   PetscBool stageok, accept = PETSC_TRUE;
1871566a47fSLisandro Dalcin   PetscReal next_time_step = ts->time_step;
188ef7bb5aaSJed Brown 
189ef7bb5aaSJed Brown   PetscFunctionBegin;
1909566063dSJacob Faibussowitsch   PetscCall((*ssp->onestep)(ts, ts->ptime, ts->time_step, sol));
1919566063dSJacob Faibussowitsch   PetscCall(TSPostStage(ts, ts->ptime, 0, &sol));
1929566063dSJacob Faibussowitsch   PetscCall(TSAdaptCheckStage(ts->adapt, ts, ts->ptime + ts->time_step, sol, &stageok));
1939371c9d4SSatish Balay   if (!stageok) {
1949371c9d4SSatish Balay     ts->reason = TS_DIVERGED_STEP_REJECTED;
1959371c9d4SSatish Balay     PetscFunctionReturn(0);
1969371c9d4SSatish Balay   }
1971566a47fSLisandro Dalcin 
1989566063dSJacob Faibussowitsch   PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
1999371c9d4SSatish Balay   if (!accept) {
2009371c9d4SSatish Balay     ts->reason = TS_DIVERGED_STEP_REJECTED;
2019371c9d4SSatish Balay     PetscFunctionReturn(0);
2029371c9d4SSatish Balay   }
2031566a47fSLisandro Dalcin 
204186e87acSLisandro Dalcin   ts->ptime += ts->time_step;
2051566a47fSLisandro Dalcin   ts->time_step = next_time_step;
206ef7bb5aaSJed Brown   PetscFunctionReturn(0);
207ef7bb5aaSJed Brown }
208ef7bb5aaSJed Brown /*------------------------------------------------------------*/
2091566a47fSLisandro Dalcin 
2109371c9d4SSatish Balay static PetscErrorCode TSReset_SSP(TS ts) {
211277b19d0SLisandro Dalcin   TS_SSP *ssp = (TS_SSP *)ts->data;
212277b19d0SLisandro Dalcin 
213277b19d0SLisandro Dalcin   PetscFunctionBegin;
2149566063dSJacob Faibussowitsch   if (ssp->work) PetscCall(VecDestroyVecs(ssp->nwork, &ssp->work));
215277b19d0SLisandro Dalcin   ssp->nwork   = 0;
216277b19d0SLisandro Dalcin   ssp->workout = PETSC_FALSE;
217277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
218277b19d0SLisandro Dalcin }
219277b19d0SLisandro Dalcin 
2209371c9d4SSatish Balay static PetscErrorCode TSDestroy_SSP(TS ts) {
221815f1ad5SJed Brown   TS_SSP *ssp = (TS_SSP *)ts->data;
222ef7bb5aaSJed Brown 
223ef7bb5aaSJed Brown   PetscFunctionBegin;
2249566063dSJacob Faibussowitsch   PetscCall(TSReset_SSP(ts));
2259566063dSJacob Faibussowitsch   PetscCall(PetscFree(ssp->type_name));
2269566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
2279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPGetType_C", NULL));
2289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPSetType_C", NULL));
2299566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPGetNumStages_C", NULL));
2309566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPSetNumStages_C", NULL));
231ef7bb5aaSJed Brown   PetscFunctionReturn(0);
232ef7bb5aaSJed Brown }
233ef7bb5aaSJed Brown /*------------------------------------------------------------*/
234ef7bb5aaSJed Brown 
235815f1ad5SJed Brown /*@C
236815f1ad5SJed Brown    TSSSPSetType - set the SSP time integration scheme to use
237815f1ad5SJed Brown 
238815f1ad5SJed Brown    Logically Collective
239815f1ad5SJed Brown 
2404165533cSJose E. Roman    Input Parameters:
2414165533cSJose E. Roman +  ts - time stepping object
2424165533cSJose E. Roman -  ssptype - type of scheme to use
243815f1ad5SJed Brown 
244815f1ad5SJed Brown    Options Database Keys:
245815f1ad5SJed Brown    -ts_ssp_type <rks2>               : Type of SSP method (one of) rks2 rks3 rk104
246a1f556f7SAidan Hamilton    -ts_ssp_nstages<rks2: 5, rks3: 9> : Number of stages
247815f1ad5SJed Brown 
248815f1ad5SJed Brown    Level: beginner
249815f1ad5SJed Brown 
250db781477SPatrick Sanan .seealso: `TSSSP`, `TSSSPGetType()`, `TSSSPSetNumStages()`, `TSSSPRKS2`, `TSSSPRKS3`, `TSSSPRK104`
251815f1ad5SJed Brown @*/
2529371c9d4SSatish Balay PetscErrorCode TSSSPSetType(TS ts, TSSSPType ssptype) {
253815f1ad5SJed Brown   PetscFunctionBegin;
254815f1ad5SJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
255b92453a8SLisandro Dalcin   PetscValidCharPointer(ssptype, 2);
256cac4c232SBarry Smith   PetscTryMethod(ts, "TSSSPSetType_C", (TS, TSSSPType), (ts, ssptype));
257815f1ad5SJed Brown   PetscFunctionReturn(0);
258815f1ad5SJed Brown }
259815f1ad5SJed Brown 
260815f1ad5SJed Brown /*@C
261815f1ad5SJed Brown    TSSSPGetType - get the SSP time integration scheme
262815f1ad5SJed Brown 
263815f1ad5SJed Brown    Logically Collective
264815f1ad5SJed Brown 
2654165533cSJose E. Roman    Input Parameter:
2664165533cSJose E. Roman .  ts - time stepping object
267815f1ad5SJed Brown 
2684165533cSJose E. Roman    Output Parameter:
2694165533cSJose E. Roman .  type - type of scheme being used
270815f1ad5SJed Brown 
271815f1ad5SJed Brown    Level: beginner
272815f1ad5SJed Brown 
273db781477SPatrick Sanan .seealso: `TSSSP`, `TSSSPSettype()`, `TSSSPSetNumStages()`, `TSSSPRKS2`, `TSSSPRKS3`, `TSSSPRK104`
274815f1ad5SJed Brown @*/
2759371c9d4SSatish Balay PetscErrorCode TSSSPGetType(TS ts, TSSSPType *type) {
276815f1ad5SJed Brown   PetscFunctionBegin;
277815f1ad5SJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
278cac4c232SBarry Smith   PetscUseMethod(ts, "TSSSPGetType_C", (TS, TSSSPType *), (ts, type));
279815f1ad5SJed Brown   PetscFunctionReturn(0);
280815f1ad5SJed Brown }
281815f1ad5SJed Brown 
282815f1ad5SJed Brown /*@
283a1f556f7SAidan Hamilton    TSSSPSetNumStages - set the number of stages to use with the SSP method. Must be called after
284a1f556f7SAidan Hamilton    TSSSPSetType().
285815f1ad5SJed Brown 
286815f1ad5SJed Brown    Logically Collective
287815f1ad5SJed Brown 
2884165533cSJose E. Roman    Input Parameters:
2894165533cSJose E. Roman +  ts - time stepping object
2904165533cSJose E. Roman -  nstages - number of stages
291815f1ad5SJed Brown 
292815f1ad5SJed Brown    Options Database Keys:
293a1f556f7SAidan Hamilton    -ts_ssp_type <rks2>               : Type of SSP method (one of) rks2 rks3 rk104
294a1f556f7SAidan Hamilton    -ts_ssp_nstages<rks2: 5, rks3: 9> : Number of stages
295815f1ad5SJed Brown 
296815f1ad5SJed Brown    Level: beginner
297815f1ad5SJed Brown 
298db781477SPatrick Sanan .seealso: `TSSSP`, `TSSSPGetNumStages()`, `TSSSPSetNumStages()`, `TSSSPRKS2`, `TSSSPRKS3`, `TSSSPRK104`
299815f1ad5SJed Brown @*/
3009371c9d4SSatish Balay PetscErrorCode TSSSPSetNumStages(TS ts, PetscInt nstages) {
301815f1ad5SJed Brown   PetscFunctionBegin;
302815f1ad5SJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
303cac4c232SBarry Smith   PetscTryMethod(ts, "TSSSPSetNumStages_C", (TS, PetscInt), (ts, nstages));
304815f1ad5SJed Brown   PetscFunctionReturn(0);
305815f1ad5SJed Brown }
306815f1ad5SJed Brown 
307815f1ad5SJed Brown /*@
308815f1ad5SJed Brown    TSSSPGetNumStages - get the number of stages in the SSP time integration scheme
309815f1ad5SJed Brown 
310815f1ad5SJed Brown    Logically Collective
311815f1ad5SJed Brown 
3124165533cSJose E. Roman    Input Parameter:
3134165533cSJose E. Roman .  ts - time stepping object
314815f1ad5SJed Brown 
3154165533cSJose E. Roman    Output Parameter:
3164165533cSJose E. Roman .  nstages - number of stages
317815f1ad5SJed Brown 
318815f1ad5SJed Brown    Level: beginner
319815f1ad5SJed Brown 
320db781477SPatrick Sanan .seealso: `TSSSP`, `TSSSPGetType()`, `TSSSPSetNumStages()`, `TSSSPRKS2`, `TSSSPRKS3`, `TSSSPRK104`
321815f1ad5SJed Brown @*/
3229371c9d4SSatish Balay PetscErrorCode TSSSPGetNumStages(TS ts, PetscInt *nstages) {
323815f1ad5SJed Brown   PetscFunctionBegin;
324815f1ad5SJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
325cac4c232SBarry Smith   PetscUseMethod(ts, "TSSSPGetNumStages_C", (TS, PetscInt *), (ts, nstages));
326815f1ad5SJed Brown   PetscFunctionReturn(0);
327815f1ad5SJed Brown }
328815f1ad5SJed Brown 
3299371c9d4SSatish Balay static PetscErrorCode TSSSPSetType_SSP(TS ts, TSSSPType type) {
330ef7bb5aaSJed Brown   TS_SSP *ssp = (TS_SSP *)ts->data;
3315f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(TS, PetscReal, PetscReal, Vec);
332a1f556f7SAidan Hamilton   PetscBool flag;
333ef7bb5aaSJed Brown 
334ef7bb5aaSJed Brown   PetscFunctionBegin;
3359566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(TSSSPList, type, &r));
3363c633725SBarry Smith   PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TS_SSP type %s given", type);
337ef7bb5aaSJed Brown   ssp->onestep = r;
3389566063dSJacob Faibussowitsch   PetscCall(PetscFree(ssp->type_name));
3399566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(type, &ssp->type_name));
3402ffb9264SLisandro Dalcin   ts->default_adapt_type = TSADAPTNONE;
341a1f556f7SAidan Hamilton   PetscCall(PetscStrcmp(type, TSSSPRKS2, &flag));
342a1f556f7SAidan Hamilton   if (flag) {
343a1f556f7SAidan Hamilton     ssp->nstages = 5;
344a1f556f7SAidan Hamilton   } else {
345a1f556f7SAidan Hamilton     PetscCall(PetscStrcmp(type, TSSSPRKS3, &flag));
346a1f556f7SAidan Hamilton     if (flag) {
347a1f556f7SAidan Hamilton       ssp->nstages = 9;
348a1f556f7SAidan Hamilton     } else {
349a1f556f7SAidan Hamilton       ssp->nstages = 5;
350a1f556f7SAidan Hamilton     }
351a1f556f7SAidan Hamilton   }
352ef7bb5aaSJed Brown   PetscFunctionReturn(0);
353ef7bb5aaSJed Brown }
3549371c9d4SSatish Balay static PetscErrorCode TSSSPGetType_SSP(TS ts, TSSSPType *type) {
355815f1ad5SJed Brown   TS_SSP *ssp = (TS_SSP *)ts->data;
356815f1ad5SJed Brown 
357815f1ad5SJed Brown   PetscFunctionBegin;
3585164425aSJed Brown   *type = ssp->type_name;
359815f1ad5SJed Brown   PetscFunctionReturn(0);
360815f1ad5SJed Brown }
3619371c9d4SSatish Balay static PetscErrorCode TSSSPSetNumStages_SSP(TS ts, PetscInt nstages) {
362815f1ad5SJed Brown   TS_SSP *ssp = (TS_SSP *)ts->data;
363815f1ad5SJed Brown 
364815f1ad5SJed Brown   PetscFunctionBegin;
365815f1ad5SJed Brown   ssp->nstages = nstages;
366815f1ad5SJed Brown   PetscFunctionReturn(0);
367815f1ad5SJed Brown }
3689371c9d4SSatish Balay static PetscErrorCode TSSSPGetNumStages_SSP(TS ts, PetscInt *nstages) {
369815f1ad5SJed Brown   TS_SSP *ssp = (TS_SSP *)ts->data;
370815f1ad5SJed Brown 
371815f1ad5SJed Brown   PetscFunctionBegin;
372815f1ad5SJed Brown   *nstages = ssp->nstages;
373815f1ad5SJed Brown   PetscFunctionReturn(0);
374815f1ad5SJed Brown }
375ef7bb5aaSJed Brown 
3769371c9d4SSatish Balay static PetscErrorCode TSSetFromOptions_SSP(TS ts, PetscOptionItems *PetscOptionsObject) {
377ef7bb5aaSJed Brown   char      tname[256] = TSSSPRKS2;
378ef7bb5aaSJed Brown   TS_SSP   *ssp        = (TS_SSP *)ts->data;
379ace3abfcSBarry Smith   PetscBool flg;
380ef7bb5aaSJed Brown 
381ef7bb5aaSJed Brown   PetscFunctionBegin;
382d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SSP ODE solver options");
383ef7bb5aaSJed Brown   {
3849566063dSJacob Faibussowitsch     PetscCall(PetscOptionsFList("-ts_ssp_type", "Type of SSP method", "TSSSPSetType", TSSSPList, tname, tname, sizeof(tname), &flg));
3851baa6e33SBarry Smith     if (flg) PetscCall(TSSSPSetType(ts, tname));
3869566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-ts_ssp_nstages", "Number of stages", "TSSSPSetNumStages", ssp->nstages, &ssp->nstages, NULL));
387ef7bb5aaSJed Brown   }
388d0609cedSBarry Smith   PetscOptionsHeadEnd();
389ef7bb5aaSJed Brown   PetscFunctionReturn(0);
390ef7bb5aaSJed Brown }
391ef7bb5aaSJed Brown 
3929371c9d4SSatish Balay static PetscErrorCode TSView_SSP(TS ts, PetscViewer viewer) {
3931566a47fSLisandro Dalcin   TS_SSP   *ssp = (TS_SSP *)ts->data;
3941566a47fSLisandro Dalcin   PetscBool ascii;
3951566a47fSLisandro Dalcin 
396ef7bb5aaSJed Brown   PetscFunctionBegin;
3979566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &ascii));
3989566063dSJacob Faibussowitsch   if (ascii) PetscCall(PetscViewerASCIIPrintf(viewer, "  Scheme: %s\n", ssp->type_name));
399ef7bb5aaSJed Brown   PetscFunctionReturn(0);
400ef7bb5aaSJed Brown }
401ef7bb5aaSJed Brown 
402ef7bb5aaSJed Brown /* ------------------------------------------------------------ */
403ef7bb5aaSJed Brown 
404ef7bb5aaSJed Brown /*MC
4058ab3e0fcSJed Brown       TSSSP - Explicit strong stability preserving ODE solver
4068ab3e0fcSJed Brown 
4078ab3e0fcSJed Brown   Most hyperbolic conservation laws have exact solutions that are total variation diminishing (TVD) or total variation
4088ab3e0fcSJed Brown   bounded (TVB) although these solutions often contain discontinuities.  Spatial discretizations such as Godunov's
4098ab3e0fcSJed Brown   scheme and high-resolution finite volume methods (TVD limiters, ENO/WENO) are designed to preserve these properties,
4108ab3e0fcSJed Brown   but they are usually formulated using a forward Euler time discretization or by coupling the space and time
4118ab3e0fcSJed Brown   discretization as in the classical Lax-Wendroff scheme.  When the space and time discretization is coupled, it is very
4128ab3e0fcSJed Brown   difficult to produce schemes with high temporal accuracy while preserving TVD properties.  An alternative is the
4138ab3e0fcSJed Brown   semidiscrete formulation where we choose a spatial discretization that is TVD with forward Euler and then choose a
4148ab3e0fcSJed Brown   time discretization that preserves the TVD property.  Such integrators are called strong stability preserving (SSP).
4158ab3e0fcSJed Brown 
4168ab3e0fcSJed Brown   Let c_eff be the minimum number of function evaluations required to step as far as one step of forward Euler while
4178ab3e0fcSJed Brown   still being SSP.  Some theoretical bounds
4188ab3e0fcSJed Brown 
4198ab3e0fcSJed Brown   1. There are no explicit methods with c_eff > 1.
4200085e20eSJed Brown 
4218ab3e0fcSJed Brown   2. There are no explicit methods beyond order 4 (for nonlinear problems) and c_eff > 0.
4220085e20eSJed Brown 
4238ab3e0fcSJed Brown   3. There are no implicit methods with order greater than 1 and c_eff > 2.
4248ab3e0fcSJed Brown 
4258ab3e0fcSJed Brown   This integrator provides Runge-Kutta methods of order 2, 3, and 4 with maximal values of c_eff.  More stages allows
4268ab3e0fcSJed Brown   for larger values of c_eff which improves efficiency.  These implementations are low-memory and only use 2 or 3 work
4278ab3e0fcSJed Brown   vectors regardless of the total number of stages, so e.g. 25-stage 3rd order methods may be an excellent choice.
4288ab3e0fcSJed Brown 
4298ab3e0fcSJed Brown   Methods can be chosen with -ts_ssp_type {rks2,rks3,rk104}
4308ab3e0fcSJed Brown 
4318ab3e0fcSJed Brown   rks2: Second order methods with any number s>1 of stages.  c_eff = (s-1)/s
4328ab3e0fcSJed Brown 
4338ab3e0fcSJed Brown   rks3: Third order methods with s=n^2 stages, n>1.  c_eff = (s-n)/s
4348ab3e0fcSJed Brown 
4358ab3e0fcSJed Brown   rk104: A 10-stage fourth order method.  c_eff = 0.6
436ef7bb5aaSJed Brown 
437ef7bb5aaSJed Brown   Level: beginner
438ef7bb5aaSJed Brown 
4397b6bb2c6SJed Brown   References:
440606c0280SSatish Balay +  * - Ketcheson, Highly efficient strong stability preserving Runge Kutta methods with low storage implementations, SISC, 2008.
441606c0280SSatish Balay -  * - Gottlieb, Ketcheson, and Shu, High order strong stability preserving time discretizations, J Scientific Computing, 2009.
4427b6bb2c6SJed Brown 
443db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()`
444ef7bb5aaSJed Brown 
445ef7bb5aaSJed Brown M*/
4469371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TSCreate_SSP(TS ts) {
447ef7bb5aaSJed Brown   TS_SSP *ssp;
448ef7bb5aaSJed Brown 
449ef7bb5aaSJed Brown   PetscFunctionBegin;
4509566063dSJacob Faibussowitsch   PetscCall(TSSSPInitializePackage());
451ef7bb5aaSJed Brown 
452ef7bb5aaSJed Brown   ts->ops->setup          = TSSetUp_SSP;
453ef7bb5aaSJed Brown   ts->ops->step           = TSStep_SSP;
454277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_SSP;
455ef7bb5aaSJed Brown   ts->ops->destroy        = TSDestroy_SSP;
456ef7bb5aaSJed Brown   ts->ops->setfromoptions = TSSetFromOptions_SSP;
457ef7bb5aaSJed Brown   ts->ops->view           = TSView_SSP;
458ef7bb5aaSJed Brown 
4599566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(ts, &ssp));
460ef7bb5aaSJed Brown   ts->data = (void *)ssp;
461ef7bb5aaSJed Brown 
4629566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPGetType_C", TSSSPGetType_SSP));
4639566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPSetType_C", TSSSPSetType_SSP));
4649566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPGetNumStages_C", TSSSPGetNumStages_SSP));
4659566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPSetNumStages_C", TSSSPSetNumStages_SSP));
466815f1ad5SJed Brown 
4679566063dSJacob Faibussowitsch   PetscCall(TSSSPSetType(ts, TSSSPRKS2));
468ef7bb5aaSJed Brown   PetscFunctionReturn(0);
469ef7bb5aaSJed Brown }
470787849ffSJed Brown 
471787849ffSJed Brown /*@C
472787849ffSJed Brown   TSSSPInitializePackage - This function initializes everything in the TSSSP package. It is called
4738a690491SBarry Smith   from TSInitializePackage().
474787849ffSJed Brown 
475787849ffSJed Brown   Level: developer
476787849ffSJed Brown 
477db781477SPatrick Sanan .seealso: `PetscInitialize()`
478787849ffSJed Brown @*/
4799371c9d4SSatish Balay PetscErrorCode TSSSPInitializePackage(void) {
480787849ffSJed Brown   PetscFunctionBegin;
481787849ffSJed Brown   if (TSSSPPackageInitialized) PetscFunctionReturn(0);
482787849ffSJed Brown   TSSSPPackageInitialized = PETSC_TRUE;
4839566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&TSSSPList, TSSSPRKS2, TSSSPStep_RK_2));
4849566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&TSSSPList, TSSSPRKS3, TSSSPStep_RK_3));
4859566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&TSSSPList, TSSSPRK104, TSSSPStep_RK_10_4));
4869566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(TSSSPFinalizePackage));
487787849ffSJed Brown   PetscFunctionReturn(0);
488787849ffSJed Brown }
489787849ffSJed Brown 
490787849ffSJed Brown /*@C
491787849ffSJed Brown   TSSSPFinalizePackage - This function destroys everything in the TSSSP package. It is
492787849ffSJed Brown   called from PetscFinalize().
493787849ffSJed Brown 
494787849ffSJed Brown   Level: developer
495787849ffSJed Brown 
496db781477SPatrick Sanan .seealso: `PetscFinalize()`
497787849ffSJed Brown @*/
4989371c9d4SSatish Balay PetscErrorCode TSSSPFinalizePackage(void) {
499787849ffSJed Brown   PetscFunctionBegin;
500787849ffSJed Brown   TSSSPPackageInitialized = PETSC_FALSE;
5019566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&TSSSPList));
502787849ffSJed Brown   PetscFunctionReturn(0);
503787849ffSJed Brown }
504