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 47815f1ad5SJed Brown TSSSPRKS2 - Optimal second order SSP Runge-Kutta method, low-storage, c_eff=(s-1)/s 48815f1ad5SJed Brown 49815f1ad5SJed Brown Pseudocode 2 of Ketcheson 2008 50815f1ad5SJed Brown 51b330ce4dSSatish Balay Level: beginner 52b330ce4dSSatish Balay 53*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSSSP`, `TSSSPSetType()`, `TSSSPSetNumStages()` 54815f1ad5SJed Brown M*/ 55d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSSPStep_RK_2(TS ts, PetscReal t0, PetscReal dt, Vec sol) 56d71ae5a4SJacob Faibussowitsch { 57ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP *)ts->data; 58ef7bb5aaSJed Brown Vec *work, F; 59ef7bb5aaSJed Brown PetscInt i, s; 60ef7bb5aaSJed Brown 61ef7bb5aaSJed Brown PetscFunctionBegin; 62ef7bb5aaSJed Brown s = ssp->nstages; 639566063dSJacob Faibussowitsch PetscCall(TSSSPGetWorkVectors(ts, 2, &work)); 64ef7bb5aaSJed Brown F = work[1]; 659566063dSJacob Faibussowitsch PetscCall(VecCopy(sol, work[0])); 66ef7bb5aaSJed Brown for (i = 0; i < s - 1; i++) { 67b8123daeSJed Brown PetscReal stage_time = t0 + dt * (i / (s - 1.)); 689566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, stage_time)); 699566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F)); 709566063dSJacob Faibussowitsch PetscCall(VecAXPY(work[0], dt / (s - 1.), F)); 71ef7bb5aaSJed Brown } 729566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, t0 + dt, work[0], F)); 739566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(sol, (s - 1.) / s, dt / s, 1. / s, work[0], F)); 749566063dSJacob Faibussowitsch PetscCall(TSSSPRestoreWorkVectors(ts, 2, &work)); 753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 76ef7bb5aaSJed Brown } 77ef7bb5aaSJed Brown 78815f1ad5SJed Brown /*MC 79aaf9cf16SJed Brown TSSSPRKS3 - Optimal third order SSP Runge-Kutta, low-storage, c_eff=(PetscSqrtReal(s)-1)/PetscSqrtReal(s), where PetscSqrtReal(s) is an integer 80815f1ad5SJed Brown 81815f1ad5SJed Brown Pseudocode 2 of Ketcheson 2008 82815f1ad5SJed Brown 83b330ce4dSSatish Balay Level: beginner 84b330ce4dSSatish Balay 85*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSSSP`, `TSSSPSetType()`, `TSSSPSetNumStages()` 86815f1ad5SJed Brown M*/ 87d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSSPStep_RK_3(TS ts, PetscReal t0, PetscReal dt, Vec sol) 88d71ae5a4SJacob Faibussowitsch { 89ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP *)ts->data; 90ef7bb5aaSJed Brown Vec *work, F; 91ef7bb5aaSJed Brown PetscInt i, s, n, r; 92b8123daeSJed Brown PetscReal c, stage_time; 93ef7bb5aaSJed Brown 94ef7bb5aaSJed Brown PetscFunctionBegin; 95ef7bb5aaSJed Brown s = ssp->nstages; 96aaf9cf16SJed Brown n = (PetscInt)(PetscSqrtReal((PetscReal)s) + 0.001); 97ef7bb5aaSJed Brown r = s - n; 9863a3b9bcSJacob 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); 999566063dSJacob Faibussowitsch PetscCall(TSSSPGetWorkVectors(ts, 3, &work)); 100ef7bb5aaSJed Brown F = work[2]; 1019566063dSJacob Faibussowitsch PetscCall(VecCopy(sol, work[0])); 102ef7bb5aaSJed Brown for (i = 0; i < (n - 1) * (n - 2) / 2; i++) { 103ef7bb5aaSJed Brown c = (i < n * (n + 1) / 2) ? 1. * i / (s - n) : (1. * i - n) / (s - n); 104b8123daeSJed Brown stage_time = t0 + c * dt; 1059566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, stage_time)); 1069566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F)); 1079566063dSJacob Faibussowitsch PetscCall(VecAXPY(work[0], dt / r, F)); 108ef7bb5aaSJed Brown } 1099566063dSJacob Faibussowitsch PetscCall(VecCopy(work[0], work[1])); 110ef7bb5aaSJed Brown for (; i < n * (n + 1) / 2 - 1; i++) { 111ef7bb5aaSJed Brown c = (i < n * (n + 1) / 2) ? 1. * i / (s - n) : (1. * i - n) / (s - n); 112b8123daeSJed Brown stage_time = t0 + c * dt; 1139566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, stage_time)); 1149566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F)); 1159566063dSJacob Faibussowitsch PetscCall(VecAXPY(work[0], dt / r, F)); 116ef7bb5aaSJed Brown } 117ef7bb5aaSJed Brown { 118ef7bb5aaSJed Brown c = (i < n * (n + 1) / 2) ? 1. * i / (s - n) : (1. * i - n) / (s - n); 119b8123daeSJed Brown stage_time = t0 + c * dt; 1209566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, stage_time)); 1219566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F)); 1229566063dSJacob 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)); 123ef7bb5aaSJed Brown i++; 124ef7bb5aaSJed Brown } 125ef7bb5aaSJed Brown for (; i < s; i++) { 126ef7bb5aaSJed Brown c = (i < n * (n + 1) / 2) ? 1. * i / (s - n) : (1. * i - n) / (s - n); 127b8123daeSJed Brown stage_time = t0 + c * dt; 1289566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, stage_time)); 1299566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F)); 1309566063dSJacob Faibussowitsch PetscCall(VecAXPY(work[0], dt / r, F)); 131ef7bb5aaSJed Brown } 1329566063dSJacob Faibussowitsch PetscCall(VecCopy(work[0], sol)); 1339566063dSJacob Faibussowitsch PetscCall(TSSSPRestoreWorkVectors(ts, 3, &work)); 1343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 135ef7bb5aaSJed Brown } 136ef7bb5aaSJed Brown 137815f1ad5SJed Brown /*MC 138b330ce4dSSatish Balay TSSSPRKS104 - Optimal fourth order SSP Runge-Kutta, low-storage (2N), c_eff=0.6 139815f1ad5SJed Brown 140815f1ad5SJed Brown SSPRK(10,4), Pseudocode 3 of Ketcheson 2008 141815f1ad5SJed Brown 142b330ce4dSSatish Balay Level: beginner 143b330ce4dSSatish Balay 144*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSSSP`, `TSSSPSetType()` 145815f1ad5SJed Brown M*/ 146d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSSPStep_RK_10_4(TS ts, PetscReal t0, PetscReal dt, Vec sol) 147d71ae5a4SJacob Faibussowitsch { 148ef7bb5aaSJed Brown const PetscReal c[10] = {0, 1. / 6, 2. / 6, 3. / 6, 4. / 6, 2. / 6, 3. / 6, 4. / 6, 5. / 6, 1}; 149ef7bb5aaSJed Brown Vec *work, F; 1505a586d82SBarry Smith PetscInt i; 151b8123daeSJed Brown PetscReal stage_time; 152ef7bb5aaSJed Brown 153ef7bb5aaSJed Brown PetscFunctionBegin; 1549566063dSJacob Faibussowitsch PetscCall(TSSSPGetWorkVectors(ts, 3, &work)); 155ef7bb5aaSJed Brown F = work[2]; 1569566063dSJacob Faibussowitsch PetscCall(VecCopy(sol, work[0])); 157ef7bb5aaSJed Brown for (i = 0; i < 5; i++) { 158b8123daeSJed Brown stage_time = t0 + c[i] * dt; 1599566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, stage_time)); 1609566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F)); 1619566063dSJacob Faibussowitsch PetscCall(VecAXPY(work[0], dt / 6, F)); 162ef7bb5aaSJed Brown } 1639566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(work[1], 1. / 25, 9. / 25, 0, sol, work[0])); 1649566063dSJacob Faibussowitsch PetscCall(VecAXPBY(work[0], 15, -5, work[1])); 165ef7bb5aaSJed Brown for (; i < 9; i++) { 166b8123daeSJed Brown stage_time = t0 + c[i] * dt; 1679566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, stage_time)); 1689566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F)); 1699566063dSJacob Faibussowitsch PetscCall(VecAXPY(work[0], dt / 6, F)); 170ef7bb5aaSJed Brown } 171b8123daeSJed Brown stage_time = t0 + dt; 1729566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, stage_time)); 1739566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F)); 1749566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(work[1], 3. / 5, dt / 10, 1, work[0], F)); 1759566063dSJacob Faibussowitsch PetscCall(VecCopy(work[1], sol)); 1769566063dSJacob Faibussowitsch PetscCall(TSSSPRestoreWorkVectors(ts, 3, &work)); 1773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 178ef7bb5aaSJed Brown } 179ef7bb5aaSJed Brown 180d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_SSP(TS ts) 181d71ae5a4SJacob Faibussowitsch { 182ef7bb5aaSJed Brown PetscFunctionBegin; 1839566063dSJacob Faibussowitsch PetscCall(TSCheckImplicitTerm(ts)); 1849566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &ts->adapt)); 1859566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(ts->adapt)); 1863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 187ef7bb5aaSJed Brown } 188ef7bb5aaSJed Brown 189d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_SSP(TS ts) 190d71ae5a4SJacob Faibussowitsch { 191ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP *)ts->data; 192ef7bb5aaSJed Brown Vec sol = ts->vec_sol; 1931566a47fSLisandro Dalcin PetscBool stageok, accept = PETSC_TRUE; 1941566a47fSLisandro Dalcin PetscReal next_time_step = ts->time_step; 195ef7bb5aaSJed Brown 196ef7bb5aaSJed Brown PetscFunctionBegin; 1979566063dSJacob Faibussowitsch PetscCall((*ssp->onestep)(ts, ts->ptime, ts->time_step, sol)); 1989566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, ts->ptime, 0, &sol)); 1999566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt, ts, ts->ptime + ts->time_step, sol, &stageok)); 2009371c9d4SSatish Balay if (!stageok) { 2019371c9d4SSatish Balay ts->reason = TS_DIVERGED_STEP_REJECTED; 2023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2039371c9d4SSatish Balay } 2041566a47fSLisandro Dalcin 2059566063dSJacob Faibussowitsch PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept)); 2069371c9d4SSatish Balay if (!accept) { 2079371c9d4SSatish Balay ts->reason = TS_DIVERGED_STEP_REJECTED; 2083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2099371c9d4SSatish Balay } 2101566a47fSLisandro Dalcin 211186e87acSLisandro Dalcin ts->ptime += ts->time_step; 2121566a47fSLisandro Dalcin ts->time_step = next_time_step; 2133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 214ef7bb5aaSJed Brown } 215ef7bb5aaSJed Brown /*------------------------------------------------------------*/ 2161566a47fSLisandro Dalcin 217d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_SSP(TS ts) 218d71ae5a4SJacob Faibussowitsch { 219277b19d0SLisandro Dalcin TS_SSP *ssp = (TS_SSP *)ts->data; 220277b19d0SLisandro Dalcin 221277b19d0SLisandro Dalcin PetscFunctionBegin; 2229566063dSJacob Faibussowitsch if (ssp->work) PetscCall(VecDestroyVecs(ssp->nwork, &ssp->work)); 223277b19d0SLisandro Dalcin ssp->nwork = 0; 224277b19d0SLisandro Dalcin ssp->workout = PETSC_FALSE; 2253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 226277b19d0SLisandro Dalcin } 227277b19d0SLisandro Dalcin 228d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_SSP(TS ts) 229d71ae5a4SJacob Faibussowitsch { 230815f1ad5SJed Brown TS_SSP *ssp = (TS_SSP *)ts->data; 231ef7bb5aaSJed Brown 232ef7bb5aaSJed Brown PetscFunctionBegin; 2339566063dSJacob Faibussowitsch PetscCall(TSReset_SSP(ts)); 2349566063dSJacob Faibussowitsch PetscCall(PetscFree(ssp->type_name)); 2359566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 2369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPGetType_C", NULL)); 2379566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPSetType_C", NULL)); 2389566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPGetNumStages_C", NULL)); 2399566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPSetNumStages_C", NULL)); 2403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 241ef7bb5aaSJed Brown } 242ef7bb5aaSJed Brown /*------------------------------------------------------------*/ 243ef7bb5aaSJed Brown 244815f1ad5SJed Brown /*@C 245bcf0153eSBarry Smith TSSSPSetType - set the `TSSSP` time integration scheme to use 246815f1ad5SJed Brown 247815f1ad5SJed Brown Logically Collective 248815f1ad5SJed Brown 2494165533cSJose E. Roman Input Parameters: 2504165533cSJose E. Roman + ts - time stepping object 2514165533cSJose E. Roman - ssptype - type of scheme to use 252815f1ad5SJed Brown 253815f1ad5SJed Brown Options Database Keys: 254bcf0153eSBarry Smith + -ts_ssp_type <rks2> - Type of `TSSSP` method (one of) rks2 rks3 rk104 255bcf0153eSBarry Smith - -ts_ssp_nstages<rks2: 5, rks3: 9> - Number of stages 256815f1ad5SJed Brown 257815f1ad5SJed Brown Level: beginner 258815f1ad5SJed Brown 259*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSSSP`, `TSSSPGetType()`, `TSSSPSetNumStages()`, `TSSSPRKS2`, `TSSSPRKS3`, `TSSSPRK104` 260815f1ad5SJed Brown @*/ 261d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSSPSetType(TS ts, TSSSPType ssptype) 262d71ae5a4SJacob Faibussowitsch { 263815f1ad5SJed Brown PetscFunctionBegin; 264815f1ad5SJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 265b92453a8SLisandro Dalcin PetscValidCharPointer(ssptype, 2); 266cac4c232SBarry Smith PetscTryMethod(ts, "TSSSPSetType_C", (TS, TSSSPType), (ts, ssptype)); 2673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 268815f1ad5SJed Brown } 269815f1ad5SJed Brown 270815f1ad5SJed Brown /*@C 271bcf0153eSBarry Smith TSSSPGetType - get the `TSSSP` time integration scheme 272815f1ad5SJed Brown 273815f1ad5SJed Brown Logically Collective 274815f1ad5SJed Brown 2754165533cSJose E. Roman Input Parameter: 2764165533cSJose E. Roman . ts - time stepping object 277815f1ad5SJed Brown 2784165533cSJose E. Roman Output Parameter: 2794165533cSJose E. Roman . type - type of scheme being used 280815f1ad5SJed Brown 281815f1ad5SJed Brown Level: beginner 282815f1ad5SJed Brown 283*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSSSP`, `TSSSPSettype()`, `TSSSPSetNumStages()`, `TSSSPRKS2`, `TSSSPRKS3`, `TSSSPRK104` 284815f1ad5SJed Brown @*/ 285d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSSPGetType(TS ts, TSSSPType *type) 286d71ae5a4SJacob Faibussowitsch { 287815f1ad5SJed Brown PetscFunctionBegin; 288815f1ad5SJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 289cac4c232SBarry Smith PetscUseMethod(ts, "TSSSPGetType_C", (TS, TSSSPType *), (ts, type)); 2903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 291815f1ad5SJed Brown } 292815f1ad5SJed Brown 293815f1ad5SJed Brown /*@ 294bcf0153eSBarry Smith TSSSPSetNumStages - set the number of stages to use with the `TSSSP` method. Must be called after 295bcf0153eSBarry Smith `TSSSPSetType()`. 296815f1ad5SJed Brown 297815f1ad5SJed Brown Logically Collective 298815f1ad5SJed Brown 2994165533cSJose E. Roman Input Parameters: 3004165533cSJose E. Roman + ts - time stepping object 3014165533cSJose E. Roman - nstages - number of stages 302815f1ad5SJed Brown 303815f1ad5SJed Brown Options Database Keys: 304bcf0153eSBarry Smith + -ts_ssp_type <rks2> - Type of `TSSSP` method (one of) rks2 rks3 rk104 305bcf0153eSBarry Smith - -ts_ssp_nstages<rks2: 5, rks3: 9> - Number of stages 306815f1ad5SJed Brown 307815f1ad5SJed Brown Level: beginner 308815f1ad5SJed Brown 309*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSSSP`, `TSSSPGetNumStages()`, `TSSSPSetNumStages()`, `TSSSPRKS2`, `TSSSPRKS3`, `TSSSPRK104` 310815f1ad5SJed Brown @*/ 311d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSSPSetNumStages(TS ts, PetscInt nstages) 312d71ae5a4SJacob Faibussowitsch { 313815f1ad5SJed Brown PetscFunctionBegin; 314815f1ad5SJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 315cac4c232SBarry Smith PetscTryMethod(ts, "TSSSPSetNumStages_C", (TS, PetscInt), (ts, nstages)); 3163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 317815f1ad5SJed Brown } 318815f1ad5SJed Brown 319815f1ad5SJed Brown /*@ 320bcf0153eSBarry Smith TSSSPGetNumStages - get the number of stages in the `TSSSP` time integration scheme 321815f1ad5SJed Brown 322815f1ad5SJed Brown Logically Collective 323815f1ad5SJed Brown 3244165533cSJose E. Roman Input Parameter: 3254165533cSJose E. Roman . ts - time stepping object 326815f1ad5SJed Brown 3274165533cSJose E. Roman Output Parameter: 3284165533cSJose E. Roman . nstages - number of stages 329815f1ad5SJed Brown 330815f1ad5SJed Brown Level: beginner 331815f1ad5SJed Brown 332*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSSSP`, `TSSSPGetType()`, `TSSSPSetNumStages()`, `TSSSPRKS2`, `TSSSPRKS3`, `TSSSPRK104` 333815f1ad5SJed Brown @*/ 334d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSSPGetNumStages(TS ts, PetscInt *nstages) 335d71ae5a4SJacob Faibussowitsch { 336815f1ad5SJed Brown PetscFunctionBegin; 337815f1ad5SJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 338cac4c232SBarry Smith PetscUseMethod(ts, "TSSSPGetNumStages_C", (TS, PetscInt *), (ts, nstages)); 3393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 340815f1ad5SJed Brown } 341815f1ad5SJed Brown 342d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSSPSetType_SSP(TS ts, TSSSPType type) 343d71ae5a4SJacob Faibussowitsch { 344ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP *)ts->data; 3455f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(TS, PetscReal, PetscReal, Vec); 346a1f556f7SAidan Hamilton PetscBool flag; 347ef7bb5aaSJed Brown 348ef7bb5aaSJed Brown PetscFunctionBegin; 3499566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(TSSSPList, type, &r)); 3503c633725SBarry Smith PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TS_SSP type %s given", type); 351ef7bb5aaSJed Brown ssp->onestep = r; 3529566063dSJacob Faibussowitsch PetscCall(PetscFree(ssp->type_name)); 3539566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(type, &ssp->type_name)); 3542ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 355a1f556f7SAidan Hamilton PetscCall(PetscStrcmp(type, TSSSPRKS2, &flag)); 356a1f556f7SAidan Hamilton if (flag) { 357a1f556f7SAidan Hamilton ssp->nstages = 5; 358a1f556f7SAidan Hamilton } else { 359a1f556f7SAidan Hamilton PetscCall(PetscStrcmp(type, TSSSPRKS3, &flag)); 360a1f556f7SAidan Hamilton if (flag) { 361a1f556f7SAidan Hamilton ssp->nstages = 9; 362a1f556f7SAidan Hamilton } else { 363a1f556f7SAidan Hamilton ssp->nstages = 5; 364a1f556f7SAidan Hamilton } 365a1f556f7SAidan Hamilton } 3663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 367ef7bb5aaSJed Brown } 368d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSSPGetType_SSP(TS ts, TSSSPType *type) 369d71ae5a4SJacob Faibussowitsch { 370815f1ad5SJed Brown TS_SSP *ssp = (TS_SSP *)ts->data; 371815f1ad5SJed Brown 372815f1ad5SJed Brown PetscFunctionBegin; 3735164425aSJed Brown *type = ssp->type_name; 3743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 375815f1ad5SJed Brown } 376d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSSPSetNumStages_SSP(TS ts, PetscInt nstages) 377d71ae5a4SJacob Faibussowitsch { 378815f1ad5SJed Brown TS_SSP *ssp = (TS_SSP *)ts->data; 379815f1ad5SJed Brown 380815f1ad5SJed Brown PetscFunctionBegin; 381815f1ad5SJed Brown ssp->nstages = nstages; 3823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 383815f1ad5SJed Brown } 384d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSSPGetNumStages_SSP(TS ts, PetscInt *nstages) 385d71ae5a4SJacob Faibussowitsch { 386815f1ad5SJed Brown TS_SSP *ssp = (TS_SSP *)ts->data; 387815f1ad5SJed Brown 388815f1ad5SJed Brown PetscFunctionBegin; 389815f1ad5SJed Brown *nstages = ssp->nstages; 3903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 391815f1ad5SJed Brown } 392ef7bb5aaSJed Brown 393d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_SSP(TS ts, PetscOptionItems *PetscOptionsObject) 394d71ae5a4SJacob Faibussowitsch { 395ef7bb5aaSJed Brown char tname[256] = TSSSPRKS2; 396ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP *)ts->data; 397ace3abfcSBarry Smith PetscBool flg; 398ef7bb5aaSJed Brown 399ef7bb5aaSJed Brown PetscFunctionBegin; 400d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SSP ODE solver options"); 401ef7bb5aaSJed Brown { 4029566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-ts_ssp_type", "Type of SSP method", "TSSSPSetType", TSSSPList, tname, tname, sizeof(tname), &flg)); 4031baa6e33SBarry Smith if (flg) PetscCall(TSSSPSetType(ts, tname)); 4049566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_ssp_nstages", "Number of stages", "TSSSPSetNumStages", ssp->nstages, &ssp->nstages, NULL)); 405ef7bb5aaSJed Brown } 406d0609cedSBarry Smith PetscOptionsHeadEnd(); 4073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 408ef7bb5aaSJed Brown } 409ef7bb5aaSJed Brown 410d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_SSP(TS ts, PetscViewer viewer) 411d71ae5a4SJacob Faibussowitsch { 4121566a47fSLisandro Dalcin TS_SSP *ssp = (TS_SSP *)ts->data; 4131566a47fSLisandro Dalcin PetscBool ascii; 4141566a47fSLisandro Dalcin 415ef7bb5aaSJed Brown PetscFunctionBegin; 4169566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &ascii)); 4179566063dSJacob Faibussowitsch if (ascii) PetscCall(PetscViewerASCIIPrintf(viewer, " Scheme: %s\n", ssp->type_name)); 4183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 419ef7bb5aaSJed Brown } 420ef7bb5aaSJed Brown 421ef7bb5aaSJed Brown /* ------------------------------------------------------------ */ 422ef7bb5aaSJed Brown 423ef7bb5aaSJed Brown /*MC 4248ab3e0fcSJed Brown TSSSP - Explicit strong stability preserving ODE solver 4258ab3e0fcSJed Brown 4268ab3e0fcSJed Brown Most hyperbolic conservation laws have exact solutions that are total variation diminishing (TVD) or total variation 4278ab3e0fcSJed Brown bounded (TVB) although these solutions often contain discontinuities. Spatial discretizations such as Godunov's 4288ab3e0fcSJed Brown scheme and high-resolution finite volume methods (TVD limiters, ENO/WENO) are designed to preserve these properties, 4298ab3e0fcSJed Brown but they are usually formulated using a forward Euler time discretization or by coupling the space and time 4308ab3e0fcSJed Brown discretization as in the classical Lax-Wendroff scheme. When the space and time discretization is coupled, it is very 4318ab3e0fcSJed Brown difficult to produce schemes with high temporal accuracy while preserving TVD properties. An alternative is the 4328ab3e0fcSJed Brown semidiscrete formulation where we choose a spatial discretization that is TVD with forward Euler and then choose a 4338ab3e0fcSJed Brown time discretization that preserves the TVD property. Such integrators are called strong stability preserving (SSP). 4348ab3e0fcSJed Brown 4358ab3e0fcSJed Brown Let c_eff be the minimum number of function evaluations required to step as far as one step of forward Euler while 4368ab3e0fcSJed Brown still being SSP. Some theoretical bounds 4378ab3e0fcSJed Brown 4388ab3e0fcSJed Brown 1. There are no explicit methods with c_eff > 1. 4390085e20eSJed Brown 4408ab3e0fcSJed Brown 2. There are no explicit methods beyond order 4 (for nonlinear problems) and c_eff > 0. 4410085e20eSJed Brown 4428ab3e0fcSJed Brown 3. There are no implicit methods with order greater than 1 and c_eff > 2. 4438ab3e0fcSJed Brown 4448ab3e0fcSJed Brown This integrator provides Runge-Kutta methods of order 2, 3, and 4 with maximal values of c_eff. More stages allows 4458ab3e0fcSJed Brown for larger values of c_eff which improves efficiency. These implementations are low-memory and only use 2 or 3 work 4468ab3e0fcSJed Brown vectors regardless of the total number of stages, so e.g. 25-stage 3rd order methods may be an excellent choice. 4478ab3e0fcSJed Brown 4488ab3e0fcSJed Brown Methods can be chosen with -ts_ssp_type {rks2,rks3,rk104} 4498ab3e0fcSJed Brown 4508ab3e0fcSJed Brown rks2: Second order methods with any number s>1 of stages. c_eff = (s-1)/s 4518ab3e0fcSJed Brown 4528ab3e0fcSJed Brown rks3: Third order methods with s=n^2 stages, n>1. c_eff = (s-n)/s 4538ab3e0fcSJed Brown 4548ab3e0fcSJed Brown rk104: A 10-stage fourth order method. c_eff = 0.6 455ef7bb5aaSJed Brown 456ef7bb5aaSJed Brown Level: beginner 457ef7bb5aaSJed Brown 4587b6bb2c6SJed Brown References: 459606c0280SSatish Balay + * - Ketcheson, Highly efficient strong stability preserving Runge Kutta methods with low storage implementations, SISC, 2008. 460606c0280SSatish Balay - * - Gottlieb, Ketcheson, and Shu, High order strong stability preserving time discretizations, J Scientific Computing, 2009. 4617b6bb2c6SJed Brown 462*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()` 463ef7bb5aaSJed Brown 464ef7bb5aaSJed Brown M*/ 465d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_SSP(TS ts) 466d71ae5a4SJacob Faibussowitsch { 467ef7bb5aaSJed Brown TS_SSP *ssp; 468ef7bb5aaSJed Brown 469ef7bb5aaSJed Brown PetscFunctionBegin; 4709566063dSJacob Faibussowitsch PetscCall(TSSSPInitializePackage()); 471ef7bb5aaSJed Brown 472ef7bb5aaSJed Brown ts->ops->setup = TSSetUp_SSP; 473ef7bb5aaSJed Brown ts->ops->step = TSStep_SSP; 474277b19d0SLisandro Dalcin ts->ops->reset = TSReset_SSP; 475ef7bb5aaSJed Brown ts->ops->destroy = TSDestroy_SSP; 476ef7bb5aaSJed Brown ts->ops->setfromoptions = TSSetFromOptions_SSP; 477ef7bb5aaSJed Brown ts->ops->view = TSView_SSP; 478ef7bb5aaSJed Brown 4794dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&ssp)); 480ef7bb5aaSJed Brown ts->data = (void *)ssp; 481ef7bb5aaSJed Brown 4829566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPGetType_C", TSSSPGetType_SSP)); 4839566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPSetType_C", TSSSPSetType_SSP)); 4849566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPGetNumStages_C", TSSSPGetNumStages_SSP)); 4859566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPSetNumStages_C", TSSSPSetNumStages_SSP)); 486815f1ad5SJed Brown 4879566063dSJacob Faibussowitsch PetscCall(TSSSPSetType(ts, TSSSPRKS2)); 4883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 489ef7bb5aaSJed Brown } 490787849ffSJed Brown 491787849ffSJed Brown /*@C 492bcf0153eSBarry Smith TSSSPInitializePackage - This function initializes everything in the `TSSSP` package. It is called 493bcf0153eSBarry Smith from `TSInitializePackage()`. 494787849ffSJed Brown 495787849ffSJed Brown Level: developer 496787849ffSJed Brown 497*1cc06b55SBarry Smith .seealso: [](ch_ts), `PetscInitialize()`, `TSSSPFinalizePackage()`, `TSInitializePackage()` 498787849ffSJed Brown @*/ 499d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSSPInitializePackage(void) 500d71ae5a4SJacob Faibussowitsch { 501787849ffSJed Brown PetscFunctionBegin; 5023ba16761SJacob Faibussowitsch if (TSSSPPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS); 503787849ffSJed Brown TSSSPPackageInitialized = PETSC_TRUE; 5049566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&TSSSPList, TSSSPRKS2, TSSSPStep_RK_2)); 5059566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&TSSSPList, TSSSPRKS3, TSSSPStep_RK_3)); 5069566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&TSSSPList, TSSSPRK104, TSSSPStep_RK_10_4)); 5079566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSSSPFinalizePackage)); 5083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 509787849ffSJed Brown } 510787849ffSJed Brown 511787849ffSJed Brown /*@C 512bcf0153eSBarry Smith TSSSPFinalizePackage - This function destroys everything in the `TSSSP` package. It is 513bcf0153eSBarry Smith called from `PetscFinalize()`. 514787849ffSJed Brown 515787849ffSJed Brown Level: developer 516787849ffSJed Brown 517*1cc06b55SBarry Smith .seealso: [](ch_ts), `PetscFinalize()`, `TSSSPInitiallizePackage()`, `TSInitializePackage()` 518787849ffSJed Brown @*/ 519d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSSPFinalizePackage(void) 520d71ae5a4SJacob Faibussowitsch { 521787849ffSJed Brown PetscFunctionBegin; 522787849ffSJed Brown TSSSPPackageInitialized = PETSC_FALSE; 5239566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&TSSSPList)); 5243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 525787849ffSJed Brown } 526