xref: /petsc/src/ts/impls/explicit/ssp/ssp.c (revision 1d27aa22b2f6148b2c4e3f06a75e0638d6493e09)
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 
18d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSSPGetWorkVectors(TS ts, PetscInt n, Vec **work)
19d71ae5a4SJacob Faibussowitsch {
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) {
2548a46eb9SPierre Jolivet     if (ssp->nwork > 0) PetscCall(VecDestroyVecs(ssp->nwork, &ssp->work));
269566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vec_sol, n, &ssp->work));
27ef7bb5aaSJed Brown     ssp->nwork = n;
28ef7bb5aaSJed Brown   }
29ef7bb5aaSJed Brown   *work        = ssp->work;
30ef7bb5aaSJed Brown   ssp->workout = PETSC_TRUE;
313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
32ef7bb5aaSJed Brown }
33ef7bb5aaSJed Brown 
34d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSSPRestoreWorkVectors(TS ts, PetscInt n, Vec **work)
35d71ae5a4SJacob Faibussowitsch {
36ef7bb5aaSJed Brown   TS_SSP *ssp = (TS_SSP *)ts->data;
37ef7bb5aaSJed Brown 
38ef7bb5aaSJed Brown   PetscFunctionBegin;
393c633725SBarry Smith   PetscCheck(ssp->workout, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Work vectors have not been gotten");
4008401ef6SPierre Jolivet   PetscCheck(*work == ssp->work, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Wrong work vectors checked out");
41ef7bb5aaSJed Brown   ssp->workout = PETSC_FALSE;
420298fd71SBarry Smith   *work        = NULL;
433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
44ef7bb5aaSJed Brown }
45ef7bb5aaSJed Brown 
46815f1ad5SJed Brown /*MC
47*1d27aa22SBarry Smith    TSSSPRKS2 - Optimal second order SSP Runge-Kutta method, low-storage, c_eff=(s-1)/s. Pseudocode 2 of {cite}`ketcheson_2008`
48815f1ad5SJed Brown 
49b330ce4dSSatish Balay    Level: beginner
50b330ce4dSSatish Balay 
511cc06b55SBarry Smith .seealso: [](ch_ts), `TSSSP`, `TSSSPSetType()`, `TSSSPSetNumStages()`
52815f1ad5SJed Brown M*/
53d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSSPStep_RK_2(TS ts, PetscReal t0, PetscReal dt, Vec sol)
54d71ae5a4SJacob Faibussowitsch {
55ef7bb5aaSJed Brown   TS_SSP  *ssp = (TS_SSP *)ts->data;
56ef7bb5aaSJed Brown   Vec     *work, F;
57ef7bb5aaSJed Brown   PetscInt i, s;
58ef7bb5aaSJed Brown 
59ef7bb5aaSJed Brown   PetscFunctionBegin;
60ef7bb5aaSJed Brown   s = ssp->nstages;
619566063dSJacob Faibussowitsch   PetscCall(TSSSPGetWorkVectors(ts, 2, &work));
62ef7bb5aaSJed Brown   F = work[1];
639566063dSJacob Faibussowitsch   PetscCall(VecCopy(sol, work[0]));
64ef7bb5aaSJed Brown   for (i = 0; i < s - 1; i++) {
65b8123daeSJed Brown     PetscReal stage_time = t0 + dt * (i / (s - 1.));
669566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, stage_time));
679566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F));
689566063dSJacob Faibussowitsch     PetscCall(VecAXPY(work[0], dt / (s - 1.), F));
69ef7bb5aaSJed Brown   }
709566063dSJacob Faibussowitsch   PetscCall(TSComputeRHSFunction(ts, t0 + dt, work[0], F));
719566063dSJacob Faibussowitsch   PetscCall(VecAXPBYPCZ(sol, (s - 1.) / s, dt / s, 1. / s, work[0], F));
729566063dSJacob Faibussowitsch   PetscCall(TSSSPRestoreWorkVectors(ts, 2, &work));
733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
74ef7bb5aaSJed Brown }
75ef7bb5aaSJed Brown 
76815f1ad5SJed Brown /*MC
77*1d27aa22SBarry Smith    TSSSPRKS3 - Optimal third order SSP Runge-Kutta, low-storage, $c_eff=(PetscSqrtReal(s)-1)/PetscSqrtReal(s)$, where `PetscSqrtReal`(s) is an integer
78815f1ad5SJed Brown 
79*1d27aa22SBarry Smith    Pseudocode 2 of {cite}`ketcheson_2008`
80815f1ad5SJed Brown 
81b330ce4dSSatish Balay    Level: beginner
82b330ce4dSSatish Balay 
831cc06b55SBarry Smith .seealso: [](ch_ts), `TSSSP`, `TSSSPSetType()`, `TSSSPSetNumStages()`
84815f1ad5SJed Brown M*/
85d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSSPStep_RK_3(TS ts, PetscReal t0, PetscReal dt, Vec sol)
86d71ae5a4SJacob Faibussowitsch {
87ef7bb5aaSJed Brown   TS_SSP   *ssp = (TS_SSP *)ts->data;
88ef7bb5aaSJed Brown   Vec      *work, F;
89ef7bb5aaSJed Brown   PetscInt  i, s, n, r;
90b8123daeSJed Brown   PetscReal c, stage_time;
91ef7bb5aaSJed Brown 
92ef7bb5aaSJed Brown   PetscFunctionBegin;
93ef7bb5aaSJed Brown   s = ssp->nstages;
94aaf9cf16SJed Brown   n = (PetscInt)(PetscSqrtReal((PetscReal)s) + 0.001);
95ef7bb5aaSJed Brown   r = s - n;
9663a3b9bcSJacob 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);
979566063dSJacob Faibussowitsch   PetscCall(TSSSPGetWorkVectors(ts, 3, &work));
98ef7bb5aaSJed Brown   F = work[2];
999566063dSJacob Faibussowitsch   PetscCall(VecCopy(sol, work[0]));
100ef7bb5aaSJed Brown   for (i = 0; i < (n - 1) * (n - 2) / 2; i++) {
101ef7bb5aaSJed Brown     c          = (i < n * (n + 1) / 2) ? 1. * i / (s - n) : (1. * i - n) / (s - n);
102b8123daeSJed Brown     stage_time = t0 + c * dt;
1039566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, stage_time));
1049566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F));
1059566063dSJacob Faibussowitsch     PetscCall(VecAXPY(work[0], dt / r, F));
106ef7bb5aaSJed Brown   }
1079566063dSJacob Faibussowitsch   PetscCall(VecCopy(work[0], work[1]));
108ef7bb5aaSJed Brown   for (; i < n * (n + 1) / 2 - 1; i++) {
109ef7bb5aaSJed Brown     c          = (i < n * (n + 1) / 2) ? 1. * i / (s - n) : (1. * i - n) / (s - n);
110b8123daeSJed Brown     stage_time = t0 + c * dt;
1119566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, stage_time));
1129566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F));
1139566063dSJacob Faibussowitsch     PetscCall(VecAXPY(work[0], dt / r, F));
114ef7bb5aaSJed Brown   }
115ef7bb5aaSJed Brown   {
116ef7bb5aaSJed Brown     c          = (i < n * (n + 1) / 2) ? 1. * i / (s - n) : (1. * i - n) / (s - n);
117b8123daeSJed Brown     stage_time = t0 + c * dt;
1189566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, stage_time));
1199566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F));
1209566063dSJacob 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));
121ef7bb5aaSJed Brown     i++;
122ef7bb5aaSJed Brown   }
123ef7bb5aaSJed Brown   for (; i < s; i++) {
124ef7bb5aaSJed Brown     c          = (i < n * (n + 1) / 2) ? 1. * i / (s - n) : (1. * i - n) / (s - n);
125b8123daeSJed Brown     stage_time = t0 + c * dt;
1269566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, stage_time));
1279566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F));
1289566063dSJacob Faibussowitsch     PetscCall(VecAXPY(work[0], dt / r, F));
129ef7bb5aaSJed Brown   }
1309566063dSJacob Faibussowitsch   PetscCall(VecCopy(work[0], sol));
1319566063dSJacob Faibussowitsch   PetscCall(TSSSPRestoreWorkVectors(ts, 3, &work));
1323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
133ef7bb5aaSJed Brown }
134ef7bb5aaSJed Brown 
135815f1ad5SJed Brown /*MC
136b330ce4dSSatish Balay    TSSSPRKS104 - Optimal fourth order SSP Runge-Kutta, low-storage (2N), c_eff=0.6
137815f1ad5SJed Brown 
138*1d27aa22SBarry Smith    SSPRK(10,4), Pseudocode 3 of {cite}`ketcheson_2008`
139815f1ad5SJed Brown 
140b330ce4dSSatish Balay    Level: beginner
141b330ce4dSSatish Balay 
1421cc06b55SBarry Smith .seealso: [](ch_ts), `TSSSP`, `TSSSPSetType()`
143815f1ad5SJed Brown M*/
144d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSSPStep_RK_10_4(TS ts, PetscReal t0, PetscReal dt, Vec sol)
145d71ae5a4SJacob Faibussowitsch {
146ef7bb5aaSJed Brown   const PetscReal c[10] = {0, 1. / 6, 2. / 6, 3. / 6, 4. / 6, 2. / 6, 3. / 6, 4. / 6, 5. / 6, 1};
147ef7bb5aaSJed Brown   Vec            *work, F;
1485a586d82SBarry Smith   PetscInt        i;
149b8123daeSJed Brown   PetscReal       stage_time;
150ef7bb5aaSJed Brown 
151ef7bb5aaSJed Brown   PetscFunctionBegin;
1529566063dSJacob Faibussowitsch   PetscCall(TSSSPGetWorkVectors(ts, 3, &work));
153ef7bb5aaSJed Brown   F = work[2];
1549566063dSJacob Faibussowitsch   PetscCall(VecCopy(sol, work[0]));
155ef7bb5aaSJed Brown   for (i = 0; i < 5; i++) {
156b8123daeSJed Brown     stage_time = t0 + c[i] * dt;
1579566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, stage_time));
1589566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F));
1599566063dSJacob Faibussowitsch     PetscCall(VecAXPY(work[0], dt / 6, F));
160ef7bb5aaSJed Brown   }
1619566063dSJacob Faibussowitsch   PetscCall(VecAXPBYPCZ(work[1], 1. / 25, 9. / 25, 0, sol, work[0]));
1629566063dSJacob Faibussowitsch   PetscCall(VecAXPBY(work[0], 15, -5, work[1]));
163ef7bb5aaSJed Brown   for (; i < 9; i++) {
164b8123daeSJed Brown     stage_time = t0 + c[i] * dt;
1659566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, stage_time));
1669566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F));
1679566063dSJacob Faibussowitsch     PetscCall(VecAXPY(work[0], dt / 6, F));
168ef7bb5aaSJed Brown   }
169b8123daeSJed Brown   stage_time = t0 + dt;
1709566063dSJacob Faibussowitsch   PetscCall(TSPreStage(ts, stage_time));
1719566063dSJacob Faibussowitsch   PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F));
1729566063dSJacob Faibussowitsch   PetscCall(VecAXPBYPCZ(work[1], 3. / 5, dt / 10, 1, work[0], F));
1739566063dSJacob Faibussowitsch   PetscCall(VecCopy(work[1], sol));
1749566063dSJacob Faibussowitsch   PetscCall(TSSSPRestoreWorkVectors(ts, 3, &work));
1753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
176ef7bb5aaSJed Brown }
177ef7bb5aaSJed Brown 
178d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_SSP(TS ts)
179d71ae5a4SJacob Faibussowitsch {
180ef7bb5aaSJed Brown   PetscFunctionBegin;
1819566063dSJacob Faibussowitsch   PetscCall(TSCheckImplicitTerm(ts));
1829566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &ts->adapt));
1839566063dSJacob Faibussowitsch   PetscCall(TSAdaptCandidatesClear(ts->adapt));
1843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
185ef7bb5aaSJed Brown }
186ef7bb5aaSJed Brown 
187d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_SSP(TS ts)
188d71ae5a4SJacob Faibussowitsch {
189ef7bb5aaSJed Brown   TS_SSP   *ssp = (TS_SSP *)ts->data;
190ef7bb5aaSJed Brown   Vec       sol = ts->vec_sol;
1911566a47fSLisandro Dalcin   PetscBool stageok, accept = PETSC_TRUE;
1921566a47fSLisandro Dalcin   PetscReal next_time_step = ts->time_step;
193ef7bb5aaSJed Brown 
194ef7bb5aaSJed Brown   PetscFunctionBegin;
1959566063dSJacob Faibussowitsch   PetscCall((*ssp->onestep)(ts, ts->ptime, ts->time_step, sol));
1969566063dSJacob Faibussowitsch   PetscCall(TSPostStage(ts, ts->ptime, 0, &sol));
1979566063dSJacob Faibussowitsch   PetscCall(TSAdaptCheckStage(ts->adapt, ts, ts->ptime + ts->time_step, sol, &stageok));
1989371c9d4SSatish Balay   if (!stageok) {
1999371c9d4SSatish Balay     ts->reason = TS_DIVERGED_STEP_REJECTED;
2003ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2019371c9d4SSatish Balay   }
2021566a47fSLisandro Dalcin 
2039566063dSJacob Faibussowitsch   PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
2049371c9d4SSatish Balay   if (!accept) {
2059371c9d4SSatish Balay     ts->reason = TS_DIVERGED_STEP_REJECTED;
2063ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2079371c9d4SSatish Balay   }
2081566a47fSLisandro Dalcin 
209186e87acSLisandro Dalcin   ts->ptime += ts->time_step;
2101566a47fSLisandro Dalcin   ts->time_step = next_time_step;
2113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
212ef7bb5aaSJed Brown }
213ef7bb5aaSJed Brown /*------------------------------------------------------------*/
2141566a47fSLisandro Dalcin 
215d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_SSP(TS ts)
216d71ae5a4SJacob Faibussowitsch {
217277b19d0SLisandro Dalcin   TS_SSP *ssp = (TS_SSP *)ts->data;
218277b19d0SLisandro Dalcin 
219277b19d0SLisandro Dalcin   PetscFunctionBegin;
2209566063dSJacob Faibussowitsch   if (ssp->work) PetscCall(VecDestroyVecs(ssp->nwork, &ssp->work));
221277b19d0SLisandro Dalcin   ssp->nwork   = 0;
222277b19d0SLisandro Dalcin   ssp->workout = PETSC_FALSE;
2233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
224277b19d0SLisandro Dalcin }
225277b19d0SLisandro Dalcin 
226d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_SSP(TS ts)
227d71ae5a4SJacob Faibussowitsch {
228815f1ad5SJed Brown   TS_SSP *ssp = (TS_SSP *)ts->data;
229ef7bb5aaSJed Brown 
230ef7bb5aaSJed Brown   PetscFunctionBegin;
2319566063dSJacob Faibussowitsch   PetscCall(TSReset_SSP(ts));
2329566063dSJacob Faibussowitsch   PetscCall(PetscFree(ssp->type_name));
2339566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
2349566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPGetType_C", NULL));
2359566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPSetType_C", NULL));
2369566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPGetNumStages_C", NULL));
2379566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPSetNumStages_C", NULL));
2383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
239ef7bb5aaSJed Brown }
240ef7bb5aaSJed Brown /*------------------------------------------------------------*/
241ef7bb5aaSJed Brown 
242815f1ad5SJed Brown /*@C
243bcf0153eSBarry Smith   TSSSPSetType - set the `TSSSP` time integration scheme to use
244815f1ad5SJed Brown 
245815f1ad5SJed Brown   Logically Collective
246815f1ad5SJed Brown 
2474165533cSJose E. Roman   Input Parameters:
2484165533cSJose E. Roman + ts      - time stepping object
2494165533cSJose E. Roman - ssptype - type of scheme to use
250815f1ad5SJed Brown 
251815f1ad5SJed Brown   Options Database Keys:
252bcf0153eSBarry Smith + -ts_ssp_type <rks2>               - Type of `TSSSP` method (one of) rks2 rks3 rk104
253bcf0153eSBarry Smith - -ts_ssp_nstages<rks2: 5, rks3: 9> - Number of stages
254815f1ad5SJed Brown 
255815f1ad5SJed Brown   Level: beginner
256815f1ad5SJed Brown 
2571cc06b55SBarry Smith .seealso: [](ch_ts), `TSSSP`, `TSSSPGetType()`, `TSSSPSetNumStages()`, `TSSSPRKS2`, `TSSSPRKS3`, `TSSSPRK104`
258815f1ad5SJed Brown @*/
259d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSSPSetType(TS ts, TSSSPType ssptype)
260d71ae5a4SJacob Faibussowitsch {
261815f1ad5SJed Brown   PetscFunctionBegin;
262815f1ad5SJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2634f572ea9SToby Isaac   PetscAssertPointer(ssptype, 2);
264cac4c232SBarry Smith   PetscTryMethod(ts, "TSSSPSetType_C", (TS, TSSSPType), (ts, ssptype));
2653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
266815f1ad5SJed Brown }
267815f1ad5SJed Brown 
268815f1ad5SJed Brown /*@C
269bcf0153eSBarry Smith   TSSSPGetType - get the `TSSSP` time integration scheme
270815f1ad5SJed Brown 
271815f1ad5SJed Brown   Logically Collective
272815f1ad5SJed Brown 
2734165533cSJose E. Roman   Input Parameter:
2744165533cSJose E. Roman . ts - time stepping object
275815f1ad5SJed Brown 
2764165533cSJose E. Roman   Output Parameter:
2774165533cSJose E. Roman . type - type of scheme being used
278815f1ad5SJed Brown 
279815f1ad5SJed Brown   Level: beginner
280815f1ad5SJed Brown 
2811cc06b55SBarry Smith .seealso: [](ch_ts), `TSSSP`, `TSSSPSettype()`, `TSSSPSetNumStages()`, `TSSSPRKS2`, `TSSSPRKS3`, `TSSSPRK104`
282815f1ad5SJed Brown @*/
283d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSSPGetType(TS ts, TSSSPType *type)
284d71ae5a4SJacob Faibussowitsch {
285815f1ad5SJed Brown   PetscFunctionBegin;
286815f1ad5SJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
287cac4c232SBarry Smith   PetscUseMethod(ts, "TSSSPGetType_C", (TS, TSSSPType *), (ts, type));
2883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
289815f1ad5SJed Brown }
290815f1ad5SJed Brown 
291815f1ad5SJed Brown /*@
292bcf0153eSBarry Smith   TSSSPSetNumStages - set the number of stages to use with the `TSSSP` method. Must be called after
293bcf0153eSBarry Smith   `TSSSPSetType()`.
294815f1ad5SJed Brown 
295815f1ad5SJed Brown   Logically Collective
296815f1ad5SJed Brown 
2974165533cSJose E. Roman   Input Parameters:
2984165533cSJose E. Roman + ts      - time stepping object
2994165533cSJose E. Roman - nstages - number of stages
300815f1ad5SJed Brown 
301815f1ad5SJed Brown   Options Database Keys:
302bcf0153eSBarry Smith + -ts_ssp_type <rks2>               - Type of `TSSSP` method (one of) rks2 rks3 rk104
303bcf0153eSBarry Smith - -ts_ssp_nstages<rks2: 5, rks3: 9> - Number of stages
304815f1ad5SJed Brown 
305815f1ad5SJed Brown   Level: beginner
306815f1ad5SJed Brown 
30742747ad1SJacob Faibussowitsch .seealso: [](ch_ts), `TSSSP`, `TSSSPGetNumStages()`, `TSSSPRKS2`, `TSSSPRKS3`, `TSSSPRK104`
308815f1ad5SJed Brown @*/
309d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSSPSetNumStages(TS ts, PetscInt nstages)
310d71ae5a4SJacob Faibussowitsch {
311815f1ad5SJed Brown   PetscFunctionBegin;
312815f1ad5SJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
313cac4c232SBarry Smith   PetscTryMethod(ts, "TSSSPSetNumStages_C", (TS, PetscInt), (ts, nstages));
3143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
315815f1ad5SJed Brown }
316815f1ad5SJed Brown 
317815f1ad5SJed Brown /*@
318bcf0153eSBarry Smith   TSSSPGetNumStages - get the number of stages in the `TSSSP` time integration scheme
319815f1ad5SJed Brown 
320815f1ad5SJed Brown   Logically Collective
321815f1ad5SJed Brown 
3224165533cSJose E. Roman   Input Parameter:
3234165533cSJose E. Roman . ts - time stepping object
324815f1ad5SJed Brown 
3254165533cSJose E. Roman   Output Parameter:
3264165533cSJose E. Roman . nstages - number of stages
327815f1ad5SJed Brown 
328815f1ad5SJed Brown   Level: beginner
329815f1ad5SJed Brown 
3301cc06b55SBarry Smith .seealso: [](ch_ts), `TSSSP`, `TSSSPGetType()`, `TSSSPSetNumStages()`, `TSSSPRKS2`, `TSSSPRKS3`, `TSSSPRK104`
331815f1ad5SJed Brown @*/
332d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSSPGetNumStages(TS ts, PetscInt *nstages)
333d71ae5a4SJacob Faibussowitsch {
334815f1ad5SJed Brown   PetscFunctionBegin;
335815f1ad5SJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
336cac4c232SBarry Smith   PetscUseMethod(ts, "TSSSPGetNumStages_C", (TS, PetscInt *), (ts, nstages));
3373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
338815f1ad5SJed Brown }
339815f1ad5SJed Brown 
340d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSSPSetType_SSP(TS ts, TSSSPType type)
341d71ae5a4SJacob Faibussowitsch {
342ef7bb5aaSJed Brown   TS_SSP *ssp = (TS_SSP *)ts->data;
3435f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(TS, PetscReal, PetscReal, Vec);
344a1f556f7SAidan Hamilton   PetscBool flag;
345ef7bb5aaSJed Brown 
346ef7bb5aaSJed Brown   PetscFunctionBegin;
3479566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(TSSSPList, type, &r));
3486adde796SStefano Zampini   PetscCheck(r, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TS_SSP type %s given", type);
349ef7bb5aaSJed Brown   ssp->onestep = r;
3509566063dSJacob Faibussowitsch   PetscCall(PetscFree(ssp->type_name));
3519566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(type, &ssp->type_name));
3522ffb9264SLisandro Dalcin   ts->default_adapt_type = TSADAPTNONE;
353a1f556f7SAidan Hamilton   PetscCall(PetscStrcmp(type, TSSSPRKS2, &flag));
354a1f556f7SAidan Hamilton   if (flag) {
355a1f556f7SAidan Hamilton     ssp->nstages = 5;
356a1f556f7SAidan Hamilton   } else {
357a1f556f7SAidan Hamilton     PetscCall(PetscStrcmp(type, TSSSPRKS3, &flag));
358a1f556f7SAidan Hamilton     if (flag) {
359a1f556f7SAidan Hamilton       ssp->nstages = 9;
360a1f556f7SAidan Hamilton     } else {
361a1f556f7SAidan Hamilton       ssp->nstages = 5;
362a1f556f7SAidan Hamilton     }
363a1f556f7SAidan Hamilton   }
3643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
365ef7bb5aaSJed Brown }
366d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSSPGetType_SSP(TS ts, TSSSPType *type)
367d71ae5a4SJacob Faibussowitsch {
368815f1ad5SJed Brown   TS_SSP *ssp = (TS_SSP *)ts->data;
369815f1ad5SJed Brown 
370815f1ad5SJed Brown   PetscFunctionBegin;
3715164425aSJed Brown   *type = ssp->type_name;
3723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
373815f1ad5SJed Brown }
374d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSSPSetNumStages_SSP(TS ts, PetscInt nstages)
375d71ae5a4SJacob Faibussowitsch {
376815f1ad5SJed Brown   TS_SSP *ssp = (TS_SSP *)ts->data;
377815f1ad5SJed Brown 
378815f1ad5SJed Brown   PetscFunctionBegin;
379815f1ad5SJed Brown   ssp->nstages = nstages;
3803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
381815f1ad5SJed Brown }
382d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSSPGetNumStages_SSP(TS ts, PetscInt *nstages)
383d71ae5a4SJacob Faibussowitsch {
384815f1ad5SJed Brown   TS_SSP *ssp = (TS_SSP *)ts->data;
385815f1ad5SJed Brown 
386815f1ad5SJed Brown   PetscFunctionBegin;
387815f1ad5SJed Brown   *nstages = ssp->nstages;
3883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
389815f1ad5SJed Brown }
390ef7bb5aaSJed Brown 
391d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_SSP(TS ts, PetscOptionItems *PetscOptionsObject)
392d71ae5a4SJacob Faibussowitsch {
393ef7bb5aaSJed Brown   char      tname[256] = TSSSPRKS2;
394ef7bb5aaSJed Brown   TS_SSP   *ssp        = (TS_SSP *)ts->data;
395ace3abfcSBarry Smith   PetscBool flg;
396ef7bb5aaSJed Brown 
397ef7bb5aaSJed Brown   PetscFunctionBegin;
398d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SSP ODE solver options");
399ef7bb5aaSJed Brown   {
4009566063dSJacob Faibussowitsch     PetscCall(PetscOptionsFList("-ts_ssp_type", "Type of SSP method", "TSSSPSetType", TSSSPList, tname, tname, sizeof(tname), &flg));
4011baa6e33SBarry Smith     if (flg) PetscCall(TSSSPSetType(ts, tname));
4029566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-ts_ssp_nstages", "Number of stages", "TSSSPSetNumStages", ssp->nstages, &ssp->nstages, NULL));
403ef7bb5aaSJed Brown   }
404d0609cedSBarry Smith   PetscOptionsHeadEnd();
4053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
406ef7bb5aaSJed Brown }
407ef7bb5aaSJed Brown 
408d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_SSP(TS ts, PetscViewer viewer)
409d71ae5a4SJacob Faibussowitsch {
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));
4163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
417ef7bb5aaSJed Brown }
418ef7bb5aaSJed Brown 
419ef7bb5aaSJed Brown /* ------------------------------------------------------------ */
420ef7bb5aaSJed Brown 
421ef7bb5aaSJed Brown /*MC
422*1d27aa22SBarry Smith       TSSSP - Explicit strong stability preserving ODE solver {cite}`ketcheson_2008` {cite}`gottliebketchesonshu2009`
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 
4561cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`
457ef7bb5aaSJed Brown M*/
458d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_SSP(TS ts)
459d71ae5a4SJacob Faibussowitsch {
460ef7bb5aaSJed Brown   TS_SSP *ssp;
461ef7bb5aaSJed Brown 
462ef7bb5aaSJed Brown   PetscFunctionBegin;
4639566063dSJacob Faibussowitsch   PetscCall(TSSSPInitializePackage());
464ef7bb5aaSJed Brown 
465ef7bb5aaSJed Brown   ts->ops->setup          = TSSetUp_SSP;
466ef7bb5aaSJed Brown   ts->ops->step           = TSStep_SSP;
467277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_SSP;
468ef7bb5aaSJed Brown   ts->ops->destroy        = TSDestroy_SSP;
469ef7bb5aaSJed Brown   ts->ops->setfromoptions = TSSetFromOptions_SSP;
470ef7bb5aaSJed Brown   ts->ops->view           = TSView_SSP;
471ef7bb5aaSJed Brown 
4724dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&ssp));
473ef7bb5aaSJed Brown   ts->data = (void *)ssp;
474ef7bb5aaSJed Brown 
4759566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPGetType_C", TSSSPGetType_SSP));
4769566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPSetType_C", TSSSPSetType_SSP));
4779566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPGetNumStages_C", TSSSPGetNumStages_SSP));
4789566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPSetNumStages_C", TSSSPSetNumStages_SSP));
479815f1ad5SJed Brown 
4809566063dSJacob Faibussowitsch   PetscCall(TSSSPSetType(ts, TSSSPRKS2));
4813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
482ef7bb5aaSJed Brown }
483787849ffSJed Brown 
484787849ffSJed Brown /*@C
485bcf0153eSBarry Smith   TSSSPInitializePackage - This function initializes everything in the `TSSSP` package. It is called
486bcf0153eSBarry Smith   from `TSInitializePackage()`.
487787849ffSJed Brown 
488787849ffSJed Brown   Level: developer
489787849ffSJed Brown 
4901cc06b55SBarry Smith .seealso: [](ch_ts), `PetscInitialize()`, `TSSSPFinalizePackage()`, `TSInitializePackage()`
491787849ffSJed Brown @*/
492d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSSPInitializePackage(void)
493d71ae5a4SJacob Faibussowitsch {
494787849ffSJed Brown   PetscFunctionBegin;
4953ba16761SJacob Faibussowitsch   if (TSSSPPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
496787849ffSJed Brown   TSSSPPackageInitialized = PETSC_TRUE;
4979566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&TSSSPList, TSSSPRKS2, TSSSPStep_RK_2));
4989566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&TSSSPList, TSSSPRKS3, TSSSPStep_RK_3));
4999566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&TSSSPList, TSSSPRK104, TSSSPStep_RK_10_4));
5009566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(TSSSPFinalizePackage));
5013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
502787849ffSJed Brown }
503787849ffSJed Brown 
504787849ffSJed Brown /*@C
505bcf0153eSBarry Smith   TSSSPFinalizePackage - This function destroys everything in the `TSSSP` package. It is
506bcf0153eSBarry Smith   called from `PetscFinalize()`.
507787849ffSJed Brown 
508787849ffSJed Brown   Level: developer
509787849ffSJed Brown 
5101cc06b55SBarry Smith .seealso: [](ch_ts), `PetscFinalize()`, `TSSSPInitiallizePackage()`, `TSInitializePackage()`
511787849ffSJed Brown @*/
512d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSSPFinalizePackage(void)
513d71ae5a4SJacob Faibussowitsch {
514787849ffSJed Brown   PetscFunctionBegin;
515787849ffSJed Brown   TSSSPPackageInitialized = PETSC_FALSE;
5169566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&TSSSPList));
5173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
518787849ffSJed Brown }
519