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"); 422c71b3e2SJacob Faibussowitsch PetscCheckFalse(*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 55815f1ad5SJed Brown .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 87815f1ad5SJed Brown .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; 1005f80ce2aSJacob Faibussowitsch PetscCheck(n*n == s,PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for optimal third order schemes with %d stages, must be a square number at least 4",s); 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 146815f1ad5SJed Brown .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 251815f1ad5SJed Brown -ts_ssp_nstages <5>: Number of stages 252815f1ad5SJed Brown 253815f1ad5SJed Brown Level: beginner 254815f1ad5SJed Brown 255815f1ad5SJed Brown .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); 262*cac4c232SBarry 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 279815f1ad5SJed Brown .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); 285*cac4c232SBarry Smith PetscUseMethod(ts,"TSSSPGetType_C",(TS,TSSSPType*),(ts,type)); 286815f1ad5SJed Brown PetscFunctionReturn(0); 287815f1ad5SJed Brown } 288815f1ad5SJed Brown 289815f1ad5SJed Brown /*@ 290815f1ad5SJed Brown TSSSPSetNumStages - set the number of stages to use with the SSP method 291815f1ad5SJed Brown 292815f1ad5SJed Brown Logically Collective 293815f1ad5SJed Brown 2944165533cSJose E. Roman Input Parameters: 2954165533cSJose E. Roman + ts - time stepping object 2964165533cSJose E. Roman - nstages - number of stages 297815f1ad5SJed Brown 298815f1ad5SJed Brown Options Database Keys: 299815f1ad5SJed Brown -ts_ssp_type <rks2>: NumStages of SSP method (one of) rks2 rks3 rk104 300815f1ad5SJed Brown -ts_ssp_nstages <5>: Number of stages 301815f1ad5SJed Brown 302815f1ad5SJed Brown Level: beginner 303815f1ad5SJed Brown 304815f1ad5SJed Brown .seealso: TSSSP, TSSSPGetNumStages(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104 305815f1ad5SJed Brown @*/ 306815f1ad5SJed Brown PetscErrorCode TSSSPSetNumStages(TS ts,PetscInt nstages) 307815f1ad5SJed Brown { 308815f1ad5SJed Brown PetscFunctionBegin; 309815f1ad5SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 310*cac4c232SBarry Smith PetscTryMethod(ts,"TSSSPSetNumStages_C",(TS,PetscInt),(ts,nstages)); 311815f1ad5SJed Brown PetscFunctionReturn(0); 312815f1ad5SJed Brown } 313815f1ad5SJed Brown 314815f1ad5SJed Brown /*@ 315815f1ad5SJed Brown TSSSPGetNumStages - get the number of stages in the SSP time integration scheme 316815f1ad5SJed Brown 317815f1ad5SJed Brown Logically Collective 318815f1ad5SJed Brown 3194165533cSJose E. Roman Input Parameter: 3204165533cSJose E. Roman . ts - time stepping object 321815f1ad5SJed Brown 3224165533cSJose E. Roman Output Parameter: 3234165533cSJose E. Roman . nstages - number of stages 324815f1ad5SJed Brown 325815f1ad5SJed Brown Level: beginner 326815f1ad5SJed Brown 327815f1ad5SJed Brown .seealso: TSSSP, TSSSPGetType(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104 328815f1ad5SJed Brown @*/ 329815f1ad5SJed Brown PetscErrorCode TSSSPGetNumStages(TS ts,PetscInt *nstages) 330815f1ad5SJed Brown { 331815f1ad5SJed Brown PetscFunctionBegin; 332815f1ad5SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 333*cac4c232SBarry Smith PetscUseMethod(ts,"TSSSPGetNumStages_C",(TS,PetscInt*),(ts,nstages)); 334815f1ad5SJed Brown PetscFunctionReturn(0); 335815f1ad5SJed Brown } 336815f1ad5SJed Brown 337cc2e6a90SBarry Smith static PetscErrorCode TSSSPSetType_SSP(TS ts,TSSSPType type) 338ef7bb5aaSJed Brown { 339ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 3405f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(TS,PetscReal,PetscReal,Vec); 341ef7bb5aaSJed Brown 342ef7bb5aaSJed Brown PetscFunctionBegin; 3439566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(TSSSPList,type,&r)); 3443c633725SBarry Smith PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TS_SSP type %s given",type); 345ef7bb5aaSJed Brown ssp->onestep = r; 3469566063dSJacob Faibussowitsch PetscCall(PetscFree(ssp->type_name)); 3479566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(type,&ssp->type_name)); 3482ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 349ef7bb5aaSJed Brown PetscFunctionReturn(0); 350ef7bb5aaSJed Brown } 351cc2e6a90SBarry Smith static PetscErrorCode TSSSPGetType_SSP(TS ts,TSSSPType *type) 352815f1ad5SJed Brown { 353815f1ad5SJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 354815f1ad5SJed Brown 355815f1ad5SJed Brown PetscFunctionBegin; 3565164425aSJed Brown *type = ssp->type_name; 357815f1ad5SJed Brown PetscFunctionReturn(0); 358815f1ad5SJed Brown } 359cc2e6a90SBarry Smith static PetscErrorCode TSSSPSetNumStages_SSP(TS ts,PetscInt nstages) 360815f1ad5SJed Brown { 361815f1ad5SJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 362815f1ad5SJed Brown 363815f1ad5SJed Brown PetscFunctionBegin; 364815f1ad5SJed Brown ssp->nstages = nstages; 365815f1ad5SJed Brown PetscFunctionReturn(0); 366815f1ad5SJed Brown } 367cc2e6a90SBarry Smith static PetscErrorCode TSSSPGetNumStages_SSP(TS ts,PetscInt *nstages) 368815f1ad5SJed Brown { 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 3764416b707SBarry Smith static PetscErrorCode TSSetFromOptions_SSP(PetscOptionItems *PetscOptionsObject,TS ts) 377ef7bb5aaSJed Brown { 378ef7bb5aaSJed Brown char tname[256] = TSSSPRKS2; 379ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 380ace3abfcSBarry Smith PetscBool flg; 381ef7bb5aaSJed Brown 382ef7bb5aaSJed Brown PetscFunctionBegin; 3839566063dSJacob Faibussowitsch PetscCall(PetscOptionsHead(PetscOptionsObject,"SSP ODE solver options")); 384ef7bb5aaSJed Brown { 3859566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-ts_ssp_type","Type of SSP method","TSSSPSetType",TSSSPList,tname,tname,sizeof(tname),&flg)); 386ef7bb5aaSJed Brown if (flg) { 3879566063dSJacob Faibussowitsch PetscCall(TSSSPSetType(ts,tname)); 388ef7bb5aaSJed Brown } 3899566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_ssp_nstages","Number of stages","TSSSPSetNumStages",ssp->nstages,&ssp->nstages,NULL)); 390ef7bb5aaSJed Brown } 3919566063dSJacob Faibussowitsch PetscCall(PetscOptionsTail()); 392ef7bb5aaSJed Brown PetscFunctionReturn(0); 393ef7bb5aaSJed Brown } 394ef7bb5aaSJed Brown 395ef7bb5aaSJed Brown static PetscErrorCode TSView_SSP(TS ts,PetscViewer viewer) 396ef7bb5aaSJed Brown { 3971566a47fSLisandro Dalcin TS_SSP *ssp = (TS_SSP*)ts->data; 3981566a47fSLisandro Dalcin PetscBool ascii; 3991566a47fSLisandro Dalcin 400ef7bb5aaSJed Brown PetscFunctionBegin; 4019566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&ascii)); 4029566063dSJacob Faibussowitsch if (ascii) PetscCall(PetscViewerASCIIPrintf(viewer," Scheme: %s\n",ssp->type_name)); 403ef7bb5aaSJed Brown PetscFunctionReturn(0); 404ef7bb5aaSJed Brown } 405ef7bb5aaSJed Brown 406ef7bb5aaSJed Brown /* ------------------------------------------------------------ */ 407ef7bb5aaSJed Brown 408ef7bb5aaSJed Brown /*MC 4098ab3e0fcSJed Brown TSSSP - Explicit strong stability preserving ODE solver 4108ab3e0fcSJed Brown 4118ab3e0fcSJed Brown Most hyperbolic conservation laws have exact solutions that are total variation diminishing (TVD) or total variation 4128ab3e0fcSJed Brown bounded (TVB) although these solutions often contain discontinuities. Spatial discretizations such as Godunov's 4138ab3e0fcSJed Brown scheme and high-resolution finite volume methods (TVD limiters, ENO/WENO) are designed to preserve these properties, 4148ab3e0fcSJed Brown but they are usually formulated using a forward Euler time discretization or by coupling the space and time 4158ab3e0fcSJed Brown discretization as in the classical Lax-Wendroff scheme. When the space and time discretization is coupled, it is very 4168ab3e0fcSJed Brown difficult to produce schemes with high temporal accuracy while preserving TVD properties. An alternative is the 4178ab3e0fcSJed Brown semidiscrete formulation where we choose a spatial discretization that is TVD with forward Euler and then choose a 4188ab3e0fcSJed Brown time discretization that preserves the TVD property. Such integrators are called strong stability preserving (SSP). 4198ab3e0fcSJed Brown 4208ab3e0fcSJed Brown Let c_eff be the minimum number of function evaluations required to step as far as one step of forward Euler while 4218ab3e0fcSJed Brown still being SSP. Some theoretical bounds 4228ab3e0fcSJed Brown 4238ab3e0fcSJed Brown 1. There are no explicit methods with c_eff > 1. 4240085e20eSJed Brown 4258ab3e0fcSJed Brown 2. There are no explicit methods beyond order 4 (for nonlinear problems) and c_eff > 0. 4260085e20eSJed Brown 4278ab3e0fcSJed Brown 3. There are no implicit methods with order greater than 1 and c_eff > 2. 4288ab3e0fcSJed Brown 4298ab3e0fcSJed Brown This integrator provides Runge-Kutta methods of order 2, 3, and 4 with maximal values of c_eff. More stages allows 4308ab3e0fcSJed Brown for larger values of c_eff which improves efficiency. These implementations are low-memory and only use 2 or 3 work 4318ab3e0fcSJed Brown vectors regardless of the total number of stages, so e.g. 25-stage 3rd order methods may be an excellent choice. 4328ab3e0fcSJed Brown 4338ab3e0fcSJed Brown Methods can be chosen with -ts_ssp_type {rks2,rks3,rk104} 4348ab3e0fcSJed Brown 4358ab3e0fcSJed Brown rks2: Second order methods with any number s>1 of stages. c_eff = (s-1)/s 4368ab3e0fcSJed Brown 4378ab3e0fcSJed Brown rks3: Third order methods with s=n^2 stages, n>1. c_eff = (s-n)/s 4388ab3e0fcSJed Brown 4398ab3e0fcSJed Brown rk104: A 10-stage fourth order method. c_eff = 0.6 440ef7bb5aaSJed Brown 441ef7bb5aaSJed Brown Level: beginner 442ef7bb5aaSJed Brown 4437b6bb2c6SJed Brown References: 444606c0280SSatish Balay + * - Ketcheson, Highly efficient strong stability preserving Runge Kutta methods with low storage implementations, SISC, 2008. 445606c0280SSatish Balay - * - Gottlieb, Ketcheson, and Shu, High order strong stability preserving time discretizations, J Scientific Computing, 2009. 4467b6bb2c6SJed Brown 447ef7bb5aaSJed Brown .seealso: TSCreate(), TS, TSSetType() 448ef7bb5aaSJed Brown 449ef7bb5aaSJed Brown M*/ 4508cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_SSP(TS ts) 451ef7bb5aaSJed Brown { 452ef7bb5aaSJed Brown TS_SSP *ssp; 453ef7bb5aaSJed Brown 454ef7bb5aaSJed Brown PetscFunctionBegin; 4559566063dSJacob Faibussowitsch PetscCall(TSSSPInitializePackage()); 456ef7bb5aaSJed Brown 457ef7bb5aaSJed Brown ts->ops->setup = TSSetUp_SSP; 458ef7bb5aaSJed Brown ts->ops->step = TSStep_SSP; 459277b19d0SLisandro Dalcin ts->ops->reset = TSReset_SSP; 460ef7bb5aaSJed Brown ts->ops->destroy = TSDestroy_SSP; 461ef7bb5aaSJed Brown ts->ops->setfromoptions = TSSetFromOptions_SSP; 462ef7bb5aaSJed Brown ts->ops->view = TSView_SSP; 463ef7bb5aaSJed Brown 4649566063dSJacob Faibussowitsch PetscCall(PetscNewLog(ts,&ssp)); 465ef7bb5aaSJed Brown ts->data = (void*)ssp; 466ef7bb5aaSJed Brown 4679566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetType_C",TSSSPGetType_SSP)); 4689566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetType_C",TSSSPSetType_SSP)); 4699566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetNumStages_C",TSSSPGetNumStages_SSP)); 4709566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetNumStages_C",TSSSPSetNumStages_SSP)); 471815f1ad5SJed Brown 4729566063dSJacob Faibussowitsch PetscCall(TSSSPSetType(ts,TSSSPRKS2)); 473ef7bb5aaSJed Brown ssp->nstages = 5; 474ef7bb5aaSJed Brown PetscFunctionReturn(0); 475ef7bb5aaSJed Brown } 476787849ffSJed Brown 477787849ffSJed Brown /*@C 478787849ffSJed Brown TSSSPInitializePackage - This function initializes everything in the TSSSP package. It is called 4798a690491SBarry Smith from TSInitializePackage(). 480787849ffSJed Brown 481787849ffSJed Brown Level: developer 482787849ffSJed Brown 483787849ffSJed Brown .seealso: PetscInitialize() 484787849ffSJed Brown @*/ 485787849ffSJed Brown PetscErrorCode TSSSPInitializePackage(void) 486787849ffSJed Brown { 487787849ffSJed Brown PetscFunctionBegin; 488787849ffSJed Brown if (TSSSPPackageInitialized) PetscFunctionReturn(0); 489787849ffSJed Brown TSSSPPackageInitialized = PETSC_TRUE; 4909566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&TSSSPList,TSSSPRKS2, TSSSPStep_RK_2)); 4919566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&TSSSPList,TSSSPRKS3, TSSSPStep_RK_3)); 4929566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&TSSSPList,TSSSPRK104,TSSSPStep_RK_10_4)); 4939566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSSSPFinalizePackage)); 494787849ffSJed Brown PetscFunctionReturn(0); 495787849ffSJed Brown } 496787849ffSJed Brown 497787849ffSJed Brown /*@C 498787849ffSJed Brown TSSSPFinalizePackage - This function destroys everything in the TSSSP package. It is 499787849ffSJed Brown called from PetscFinalize(). 500787849ffSJed Brown 501787849ffSJed Brown Level: developer 502787849ffSJed Brown 503787849ffSJed Brown .seealso: PetscFinalize() 504787849ffSJed Brown @*/ 505787849ffSJed Brown PetscErrorCode TSSSPFinalizePackage(void) 506787849ffSJed Brown { 507787849ffSJed Brown PetscFunctionBegin; 508787849ffSJed Brown TSSSPPackageInitialized = PETSC_FALSE; 5099566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&TSSSPList)); 510787849ffSJed Brown PetscFunctionReturn(0); 511787849ffSJed Brown } 512