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 6140e18c1SBarry Smith PetscFunctionList TSSSPList = 0; 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 18ef7bb5aaSJed Brown 1905d28066SBarry Smith static PetscErrorCode TSSSPGetWorkVectors(TS ts,PetscInt n,Vec **work) 20ef7bb5aaSJed Brown { 21ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 22ef7bb5aaSJed Brown PetscErrorCode ierr; 23ef7bb5aaSJed Brown 24ef7bb5aaSJed Brown PetscFunctionBegin; 25e32f2f54SBarry Smith if (ssp->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Work vectors already gotten"); 26ef7bb5aaSJed Brown if (ssp->nwork < n) { 27ef7bb5aaSJed Brown if (ssp->nwork > 0) { 28574deadeSBarry Smith ierr = VecDestroyVecs(ssp->nwork,&ssp->work);CHKERRQ(ierr); 29ef7bb5aaSJed Brown } 30ef7bb5aaSJed Brown ierr = VecDuplicateVecs(ts->vec_sol,n,&ssp->work);CHKERRQ(ierr); 31ef7bb5aaSJed Brown ssp->nwork = n; 32ef7bb5aaSJed Brown } 33ef7bb5aaSJed Brown *work = ssp->work; 34ef7bb5aaSJed Brown ssp->workout = PETSC_TRUE; 35ef7bb5aaSJed Brown PetscFunctionReturn(0); 36ef7bb5aaSJed Brown } 37ef7bb5aaSJed Brown 3805d28066SBarry Smith static PetscErrorCode TSSSPRestoreWorkVectors(TS ts,PetscInt n,Vec **work) 39ef7bb5aaSJed Brown { 40ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 41ef7bb5aaSJed Brown 42ef7bb5aaSJed Brown PetscFunctionBegin; 43e32f2f54SBarry Smith if (!ssp->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Work vectors have not been gotten"); 44e32f2f54SBarry Smith if (*work != ssp->work) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong work vectors checked out"); 45ef7bb5aaSJed Brown ssp->workout = PETSC_FALSE; 460298fd71SBarry Smith *work = NULL; 47ef7bb5aaSJed Brown PetscFunctionReturn(0); 48ef7bb5aaSJed Brown } 49ef7bb5aaSJed Brown 50815f1ad5SJed Brown /*MC 51815f1ad5SJed Brown TSSSPRKS2 - Optimal second order SSP Runge-Kutta method, low-storage, c_eff=(s-1)/s 52815f1ad5SJed Brown 53815f1ad5SJed Brown Pseudocode 2 of Ketcheson 2008 54815f1ad5SJed Brown 55b330ce4dSSatish Balay Level: beginner 56b330ce4dSSatish Balay 57815f1ad5SJed Brown .seealso: TSSSP, TSSSPSetType(), TSSSPSetNumStages() 58815f1ad5SJed Brown M*/ 5905d28066SBarry Smith static PetscErrorCode TSSSPStep_RK_2(TS ts,PetscReal t0,PetscReal dt,Vec sol) 60ef7bb5aaSJed Brown { 61ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 62ef7bb5aaSJed Brown Vec *work,F; 63ef7bb5aaSJed Brown PetscInt i,s; 64ef7bb5aaSJed Brown PetscErrorCode ierr; 65ef7bb5aaSJed Brown 66ef7bb5aaSJed Brown PetscFunctionBegin; 67ef7bb5aaSJed Brown s = ssp->nstages; 6805d28066SBarry Smith ierr = TSSSPGetWorkVectors(ts,2,&work);CHKERRQ(ierr); 69ef7bb5aaSJed Brown F = work[1]; 70ef7bb5aaSJed Brown ierr = VecCopy(sol,work[0]);CHKERRQ(ierr); 71ef7bb5aaSJed Brown for (i=0; i<s-1; i++) { 72b8123daeSJed Brown PetscReal stage_time = t0+dt*(i/(s-1.)); 73b8123daeSJed Brown ierr = TSPreStage(ts,stage_time);CHKERRQ(ierr); 74b8123daeSJed Brown ierr = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr); 75ef7bb5aaSJed Brown ierr = VecAXPY(work[0],dt/(s-1.),F);CHKERRQ(ierr); 76ef7bb5aaSJed Brown } 77ef7bb5aaSJed Brown ierr = TSComputeRHSFunction(ts,t0+dt,work[0],F);CHKERRQ(ierr); 78ef7bb5aaSJed Brown ierr = VecAXPBYPCZ(sol,(s-1.)/s,dt/s,1./s,work[0],F);CHKERRQ(ierr); 7905d28066SBarry Smith ierr = TSSSPRestoreWorkVectors(ts,2,&work);CHKERRQ(ierr); 80ef7bb5aaSJed Brown PetscFunctionReturn(0); 81ef7bb5aaSJed Brown } 82ef7bb5aaSJed Brown 83815f1ad5SJed Brown /*MC 84aaf9cf16SJed Brown TSSSPRKS3 - Optimal third order SSP Runge-Kutta, low-storage, c_eff=(PetscSqrtReal(s)-1)/PetscSqrtReal(s), where PetscSqrtReal(s) is an integer 85815f1ad5SJed Brown 86815f1ad5SJed Brown Pseudocode 2 of Ketcheson 2008 87815f1ad5SJed Brown 88b330ce4dSSatish Balay Level: beginner 89b330ce4dSSatish Balay 90815f1ad5SJed Brown .seealso: TSSSP, TSSSPSetType(), TSSSPSetNumStages() 91815f1ad5SJed Brown M*/ 9205d28066SBarry Smith static PetscErrorCode TSSSPStep_RK_3(TS ts,PetscReal t0,PetscReal dt,Vec sol) 93ef7bb5aaSJed Brown { 94ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 95ef7bb5aaSJed Brown Vec *work,F; 96ef7bb5aaSJed Brown PetscInt i,s,n,r; 97b8123daeSJed Brown PetscReal c,stage_time; 98ef7bb5aaSJed Brown PetscErrorCode ierr; 99ef7bb5aaSJed Brown 100ef7bb5aaSJed Brown PetscFunctionBegin; 101ef7bb5aaSJed Brown s = ssp->nstages; 102aaf9cf16SJed Brown n = (PetscInt)(PetscSqrtReal((PetscReal)s)+0.001); 103ef7bb5aaSJed Brown r = s-n; 104e32f2f54SBarry Smith if (n*n != s) SETERRQ1(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); 10505d28066SBarry Smith ierr = TSSSPGetWorkVectors(ts,3,&work);CHKERRQ(ierr); 106ef7bb5aaSJed Brown F = work[2]; 107ef7bb5aaSJed Brown ierr = VecCopy(sol,work[0]);CHKERRQ(ierr); 108ef7bb5aaSJed Brown for (i=0; i<(n-1)*(n-2)/2; 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; 111b8123daeSJed Brown ierr = TSPreStage(ts,stage_time);CHKERRQ(ierr); 112b8123daeSJed Brown ierr = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr); 113ef7bb5aaSJed Brown ierr = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr); 114ef7bb5aaSJed Brown } 115ef7bb5aaSJed Brown ierr = VecCopy(work[0],work[1]);CHKERRQ(ierr); 116ef7bb5aaSJed Brown for (; i<n*(n+1)/2-1; i++) { 117ef7bb5aaSJed Brown c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n); 118b8123daeSJed Brown stage_time = t0+c*dt; 119b8123daeSJed Brown ierr = TSPreStage(ts,stage_time);CHKERRQ(ierr); 120b8123daeSJed Brown ierr = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr); 121ef7bb5aaSJed Brown ierr = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr); 122ef7bb5aaSJed Brown } 123ef7bb5aaSJed Brown { 124ef7bb5aaSJed Brown c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n); 125b8123daeSJed Brown stage_time = t0+c*dt; 126b8123daeSJed Brown ierr = TSPreStage(ts,stage_time);CHKERRQ(ierr); 127b8123daeSJed Brown ierr = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr); 128ef7bb5aaSJed Brown ierr = VecAXPBYPCZ(work[0],1.*n/(2*n-1.),(n-1.)*dt/(r*(2*n-1)),(n-1.)/(2*n-1.),work[1],F);CHKERRQ(ierr); 129ef7bb5aaSJed Brown i++; 130ef7bb5aaSJed Brown } 131ef7bb5aaSJed Brown for (; i<s; i++) { 132ef7bb5aaSJed Brown c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n); 133b8123daeSJed Brown stage_time = t0+c*dt; 134b8123daeSJed Brown ierr = TSPreStage(ts,stage_time);CHKERRQ(ierr); 135b8123daeSJed Brown ierr = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr); 136ef7bb5aaSJed Brown ierr = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr); 137ef7bb5aaSJed Brown } 138ef7bb5aaSJed Brown ierr = VecCopy(work[0],sol);CHKERRQ(ierr); 13905d28066SBarry Smith ierr = TSSSPRestoreWorkVectors(ts,3,&work);CHKERRQ(ierr); 140ef7bb5aaSJed Brown PetscFunctionReturn(0); 141ef7bb5aaSJed Brown } 142ef7bb5aaSJed Brown 143815f1ad5SJed Brown /*MC 144b330ce4dSSatish Balay TSSSPRKS104 - Optimal fourth order SSP Runge-Kutta, low-storage (2N), c_eff=0.6 145815f1ad5SJed Brown 146815f1ad5SJed Brown SSPRK(10,4), Pseudocode 3 of Ketcheson 2008 147815f1ad5SJed Brown 148b330ce4dSSatish Balay Level: beginner 149b330ce4dSSatish Balay 150815f1ad5SJed Brown .seealso: TSSSP, TSSSPSetType() 151815f1ad5SJed Brown M*/ 15205d28066SBarry Smith static PetscErrorCode TSSSPStep_RK_10_4(TS ts,PetscReal t0,PetscReal dt,Vec sol) 153ef7bb5aaSJed Brown { 154ef7bb5aaSJed Brown const PetscReal c[10] = {0, 1./6, 2./6, 3./6, 4./6, 2./6, 3./6, 4./6, 5./6, 1}; 155ef7bb5aaSJed Brown Vec *work,F; 1565a586d82SBarry Smith PetscInt i; 157b8123daeSJed Brown PetscReal stage_time; 158ef7bb5aaSJed Brown PetscErrorCode ierr; 159ef7bb5aaSJed Brown 160ef7bb5aaSJed Brown PetscFunctionBegin; 16105d28066SBarry Smith ierr = TSSSPGetWorkVectors(ts,3,&work);CHKERRQ(ierr); 162ef7bb5aaSJed Brown F = work[2]; 163ef7bb5aaSJed Brown ierr = VecCopy(sol,work[0]);CHKERRQ(ierr); 164ef7bb5aaSJed Brown for (i=0; i<5; i++) { 165b8123daeSJed Brown stage_time = t0+c[i]*dt; 166b8123daeSJed Brown ierr = TSPreStage(ts,stage_time);CHKERRQ(ierr); 167b8123daeSJed Brown ierr = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr); 168ef7bb5aaSJed Brown ierr = VecAXPY(work[0],dt/6,F);CHKERRQ(ierr); 169ef7bb5aaSJed Brown } 170ef7bb5aaSJed Brown ierr = VecAXPBYPCZ(work[1],1./25,9./25,0,sol,work[0]);CHKERRQ(ierr); 171ef7bb5aaSJed Brown ierr = VecAXPBY(work[0],15,-5,work[1]);CHKERRQ(ierr); 172ef7bb5aaSJed Brown for (; i<9; i++) { 173b8123daeSJed Brown stage_time = t0+c[i]*dt; 174b8123daeSJed Brown ierr = TSPreStage(ts,stage_time);CHKERRQ(ierr); 175b8123daeSJed Brown ierr = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr); 176ef7bb5aaSJed Brown ierr = VecAXPY(work[0],dt/6,F);CHKERRQ(ierr); 177ef7bb5aaSJed Brown } 178b8123daeSJed Brown stage_time = t0+dt; 179b8123daeSJed Brown ierr = TSPreStage(ts,stage_time);CHKERRQ(ierr); 180b8123daeSJed Brown ierr = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr); 181ef7bb5aaSJed Brown ierr = VecAXPBYPCZ(work[1],3./5,dt/10,1,work[0],F);CHKERRQ(ierr); 182ef7bb5aaSJed Brown ierr = VecCopy(work[1],sol);CHKERRQ(ierr); 18305d28066SBarry Smith ierr = TSSSPRestoreWorkVectors(ts,3,&work);CHKERRQ(ierr); 184ef7bb5aaSJed Brown PetscFunctionReturn(0); 185ef7bb5aaSJed Brown } 186ef7bb5aaSJed Brown 187ef7bb5aaSJed Brown 188ef7bb5aaSJed Brown static PetscErrorCode TSSetUp_SSP(TS ts) 189ef7bb5aaSJed Brown { 1901566a47fSLisandro Dalcin PetscErrorCode ierr; 191ef7bb5aaSJed Brown 192ef7bb5aaSJed Brown PetscFunctionBegin; 1933dd83b38SBarry Smith ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr); 1941566a47fSLisandro Dalcin ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 1951566a47fSLisandro Dalcin ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 196ef7bb5aaSJed Brown PetscFunctionReturn(0); 197ef7bb5aaSJed Brown } 198ef7bb5aaSJed Brown 199193ac0bcSJed Brown static PetscErrorCode TSStep_SSP(TS ts) 200ef7bb5aaSJed Brown { 201ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 202ef7bb5aaSJed Brown Vec sol = ts->vec_sol; 2031566a47fSLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 2041566a47fSLisandro Dalcin PetscReal next_time_step = ts->time_step; 205ef7bb5aaSJed Brown PetscErrorCode ierr; 206ef7bb5aaSJed Brown 207ef7bb5aaSJed Brown PetscFunctionBegin; 208186e87acSLisandro Dalcin ierr = (*ssp->onestep)(ts,ts->ptime,ts->time_step,sol);CHKERRQ(ierr); 2091566a47fSLisandro Dalcin ierr = TSPostStage(ts,ts->ptime,0,&sol);CHKERRQ(ierr); 2101566a47fSLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,ts->ptime+ts->time_step,sol,&stageok);CHKERRQ(ierr); 2111566a47fSLisandro Dalcin if(!stageok) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);} 2121566a47fSLisandro Dalcin 2131566a47fSLisandro Dalcin ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 2141566a47fSLisandro Dalcin if (!accept) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);} 2151566a47fSLisandro Dalcin 216186e87acSLisandro Dalcin ts->ptime += ts->time_step; 2171566a47fSLisandro Dalcin ts->time_step = next_time_step; 218ef7bb5aaSJed Brown PetscFunctionReturn(0); 219ef7bb5aaSJed Brown } 220ef7bb5aaSJed Brown /*------------------------------------------------------------*/ 2211566a47fSLisandro Dalcin 222277b19d0SLisandro Dalcin static PetscErrorCode TSReset_SSP(TS ts) 223277b19d0SLisandro Dalcin { 224277b19d0SLisandro Dalcin TS_SSP *ssp = (TS_SSP*)ts->data; 225277b19d0SLisandro Dalcin PetscErrorCode ierr; 226277b19d0SLisandro Dalcin 227277b19d0SLisandro Dalcin PetscFunctionBegin; 228277b19d0SLisandro Dalcin if (ssp->work) {ierr = VecDestroyVecs(ssp->nwork,&ssp->work);CHKERRQ(ierr);} 229277b19d0SLisandro Dalcin ssp->nwork = 0; 230277b19d0SLisandro Dalcin ssp->workout = PETSC_FALSE; 231277b19d0SLisandro Dalcin PetscFunctionReturn(0); 232277b19d0SLisandro Dalcin } 233277b19d0SLisandro Dalcin 234ef7bb5aaSJed Brown static PetscErrorCode TSDestroy_SSP(TS ts) 235ef7bb5aaSJed Brown { 236815f1ad5SJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 237ef7bb5aaSJed Brown PetscErrorCode ierr; 238ef7bb5aaSJed Brown 239ef7bb5aaSJed Brown PetscFunctionBegin; 240277b19d0SLisandro Dalcin ierr = TSReset_SSP(ts);CHKERRQ(ierr); 2415164425aSJed Brown ierr = PetscFree(ssp->type_name);CHKERRQ(ierr); 242c31cb41cSBarry Smith ierr = PetscFree(ts->data);CHKERRQ(ierr); 243bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetType_C",NULL);CHKERRQ(ierr); 244bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetType_C",NULL);CHKERRQ(ierr); 245bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetNumStages_C",NULL);CHKERRQ(ierr); 246bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetNumStages_C",NULL);CHKERRQ(ierr); 247ef7bb5aaSJed Brown PetscFunctionReturn(0); 248ef7bb5aaSJed Brown } 249ef7bb5aaSJed Brown /*------------------------------------------------------------*/ 250ef7bb5aaSJed Brown 251815f1ad5SJed Brown /*@C 252815f1ad5SJed Brown TSSSPSetType - set the SSP time integration scheme to use 253815f1ad5SJed Brown 254815f1ad5SJed Brown Logically Collective 255815f1ad5SJed Brown 256815f1ad5SJed Brown Input Arguments: 257815f1ad5SJed Brown ts - time stepping object 258815f1ad5SJed Brown type - type of scheme to use 259815f1ad5SJed Brown 260815f1ad5SJed Brown Options Database Keys: 261815f1ad5SJed Brown -ts_ssp_type <rks2>: Type of SSP method (one of) rks2 rks3 rk104 262815f1ad5SJed Brown -ts_ssp_nstages <5>: Number of stages 263815f1ad5SJed Brown 264815f1ad5SJed Brown Level: beginner 265815f1ad5SJed Brown 266815f1ad5SJed Brown .seealso: TSSSP, TSSSPGetType(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104 267815f1ad5SJed Brown @*/ 26819fd82e9SBarry Smith PetscErrorCode TSSSPSetType(TS ts,TSSSPType type) 269815f1ad5SJed Brown { 270815f1ad5SJed Brown PetscErrorCode ierr; 271815f1ad5SJed Brown 272815f1ad5SJed Brown PetscFunctionBegin; 273815f1ad5SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 27419fd82e9SBarry Smith ierr = PetscTryMethod(ts,"TSSSPSetType_C",(TS,TSSSPType),(ts,type));CHKERRQ(ierr); 275815f1ad5SJed Brown PetscFunctionReturn(0); 276815f1ad5SJed Brown } 277815f1ad5SJed Brown 278815f1ad5SJed Brown /*@C 279815f1ad5SJed Brown TSSSPGetType - get the SSP time integration scheme 280815f1ad5SJed Brown 281815f1ad5SJed Brown Logically Collective 282815f1ad5SJed Brown 283815f1ad5SJed Brown Input Argument: 284815f1ad5SJed Brown ts - time stepping object 285815f1ad5SJed Brown 286815f1ad5SJed Brown Output Argument: 287815f1ad5SJed Brown type - type of scheme being used 288815f1ad5SJed Brown 289815f1ad5SJed Brown Level: beginner 290815f1ad5SJed Brown 291815f1ad5SJed Brown .seealso: TSSSP, TSSSPSettype(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104 292815f1ad5SJed Brown @*/ 29319fd82e9SBarry Smith PetscErrorCode TSSSPGetType(TS ts,TSSSPType *type) 294815f1ad5SJed Brown { 295815f1ad5SJed Brown PetscErrorCode ierr; 296815f1ad5SJed Brown 297815f1ad5SJed Brown PetscFunctionBegin; 298815f1ad5SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 299163d334eSBarry Smith ierr = PetscUseMethod(ts,"TSSSPGetType_C",(TS,TSSSPType*),(ts,type));CHKERRQ(ierr); 300815f1ad5SJed Brown PetscFunctionReturn(0); 301815f1ad5SJed Brown } 302815f1ad5SJed Brown 303815f1ad5SJed Brown /*@ 304815f1ad5SJed Brown TSSSPSetNumStages - set the number of stages to use with the SSP method 305815f1ad5SJed Brown 306815f1ad5SJed Brown Logically Collective 307815f1ad5SJed Brown 308815f1ad5SJed Brown Input Arguments: 309815f1ad5SJed Brown ts - time stepping object 310815f1ad5SJed Brown nstages - number of stages 311815f1ad5SJed Brown 312815f1ad5SJed Brown Options Database Keys: 313815f1ad5SJed Brown -ts_ssp_type <rks2>: NumStages of SSP method (one of) rks2 rks3 rk104 314815f1ad5SJed Brown -ts_ssp_nstages <5>: Number of stages 315815f1ad5SJed Brown 316815f1ad5SJed Brown Level: beginner 317815f1ad5SJed Brown 318815f1ad5SJed Brown .seealso: TSSSP, TSSSPGetNumStages(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104 319815f1ad5SJed Brown @*/ 320815f1ad5SJed Brown PetscErrorCode TSSSPSetNumStages(TS ts,PetscInt nstages) 321815f1ad5SJed Brown { 322815f1ad5SJed Brown PetscErrorCode ierr; 323815f1ad5SJed Brown 324815f1ad5SJed Brown PetscFunctionBegin; 325815f1ad5SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 326815f1ad5SJed Brown ierr = PetscTryMethod(ts,"TSSSPSetNumStages_C",(TS,PetscInt),(ts,nstages));CHKERRQ(ierr); 327815f1ad5SJed Brown PetscFunctionReturn(0); 328815f1ad5SJed Brown } 329815f1ad5SJed Brown 330815f1ad5SJed Brown /*@ 331815f1ad5SJed Brown TSSSPGetNumStages - get the number of stages in the SSP time integration scheme 332815f1ad5SJed Brown 333815f1ad5SJed Brown Logically Collective 334815f1ad5SJed Brown 335815f1ad5SJed Brown Input Argument: 336815f1ad5SJed Brown ts - time stepping object 337815f1ad5SJed Brown 338815f1ad5SJed Brown Output Argument: 339815f1ad5SJed Brown nstages - number of stages 340815f1ad5SJed Brown 341815f1ad5SJed Brown Level: beginner 342815f1ad5SJed Brown 343815f1ad5SJed Brown .seealso: TSSSP, TSSSPGetType(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104 344815f1ad5SJed Brown @*/ 345815f1ad5SJed Brown PetscErrorCode TSSSPGetNumStages(TS ts,PetscInt *nstages) 346815f1ad5SJed Brown { 347815f1ad5SJed Brown PetscErrorCode ierr; 348815f1ad5SJed Brown 349815f1ad5SJed Brown PetscFunctionBegin; 350815f1ad5SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 351163d334eSBarry Smith ierr = PetscUseMethod(ts,"TSSSPGetNumStages_C",(TS,PetscInt*),(ts,nstages));CHKERRQ(ierr); 352815f1ad5SJed Brown PetscFunctionReturn(0); 353815f1ad5SJed Brown } 354815f1ad5SJed Brown 355cc2e6a90SBarry Smith static PetscErrorCode TSSSPSetType_SSP(TS ts,TSSSPType type) 356ef7bb5aaSJed Brown { 357ef7bb5aaSJed Brown PetscErrorCode ierr,(*r)(TS,PetscReal,PetscReal,Vec); 358ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 359ef7bb5aaSJed Brown 360ef7bb5aaSJed Brown PetscFunctionBegin; 3611c9cd337SJed Brown ierr = PetscFunctionListFind(TSSSPList,type,&r);CHKERRQ(ierr); 362e32f2f54SBarry Smith if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TS_SSP type %s given",type); 363ef7bb5aaSJed Brown ssp->onestep = r; 3645164425aSJed Brown ierr = PetscFree(ssp->type_name);CHKERRQ(ierr); 3655164425aSJed Brown ierr = PetscStrallocpy(type,&ssp->type_name);CHKERRQ(ierr); 366*2ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 367ef7bb5aaSJed Brown PetscFunctionReturn(0); 368ef7bb5aaSJed Brown } 369cc2e6a90SBarry Smith static PetscErrorCode TSSSPGetType_SSP(TS ts,TSSSPType *type) 370815f1ad5SJed Brown { 371815f1ad5SJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 372815f1ad5SJed Brown 373815f1ad5SJed Brown PetscFunctionBegin; 3745164425aSJed Brown *type = ssp->type_name; 375815f1ad5SJed Brown PetscFunctionReturn(0); 376815f1ad5SJed Brown } 377cc2e6a90SBarry Smith static PetscErrorCode TSSSPSetNumStages_SSP(TS ts,PetscInt nstages) 378815f1ad5SJed Brown { 379815f1ad5SJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 380815f1ad5SJed Brown 381815f1ad5SJed Brown PetscFunctionBegin; 382815f1ad5SJed Brown ssp->nstages = nstages; 383815f1ad5SJed Brown PetscFunctionReturn(0); 384815f1ad5SJed Brown } 385cc2e6a90SBarry Smith static PetscErrorCode TSSSPGetNumStages_SSP(TS ts,PetscInt *nstages) 386815f1ad5SJed Brown { 387815f1ad5SJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 388815f1ad5SJed Brown 389815f1ad5SJed Brown PetscFunctionBegin; 390815f1ad5SJed Brown *nstages = ssp->nstages; 391815f1ad5SJed Brown PetscFunctionReturn(0); 392815f1ad5SJed Brown } 393ef7bb5aaSJed Brown 3944416b707SBarry Smith static PetscErrorCode TSSetFromOptions_SSP(PetscOptionItems *PetscOptionsObject,TS ts) 395ef7bb5aaSJed Brown { 396ef7bb5aaSJed Brown char tname[256] = TSSSPRKS2; 397ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 398ef7bb5aaSJed Brown PetscErrorCode ierr; 399ace3abfcSBarry Smith PetscBool flg; 400ef7bb5aaSJed Brown 401ef7bb5aaSJed Brown PetscFunctionBegin; 402e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SSP ODE solver options");CHKERRQ(ierr); 403ef7bb5aaSJed Brown { 404a264d7a6SBarry Smith ierr = PetscOptionsFList("-ts_ssp_type","Type of SSP method","TSSSPSetType",TSSSPList,tname,tname,sizeof(tname),&flg);CHKERRQ(ierr); 405ef7bb5aaSJed Brown if (flg) { 406ef7bb5aaSJed Brown ierr = TSSSPSetType(ts,tname);CHKERRQ(ierr); 407ef7bb5aaSJed Brown } 4080298fd71SBarry Smith ierr = PetscOptionsInt("-ts_ssp_nstages","Number of stages","TSSSPSetNumStages",ssp->nstages,&ssp->nstages,NULL);CHKERRQ(ierr); 409ef7bb5aaSJed Brown } 410ef7bb5aaSJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 411ef7bb5aaSJed Brown PetscFunctionReturn(0); 412ef7bb5aaSJed Brown } 413ef7bb5aaSJed Brown 414ef7bb5aaSJed Brown static PetscErrorCode TSView_SSP(TS ts,PetscViewer viewer) 415ef7bb5aaSJed Brown { 4161566a47fSLisandro Dalcin TS_SSP *ssp = (TS_SSP*)ts->data; 4171566a47fSLisandro Dalcin PetscBool ascii; 4181566a47fSLisandro Dalcin PetscErrorCode ierr; 4191566a47fSLisandro Dalcin 420ef7bb5aaSJed Brown PetscFunctionBegin; 4211566a47fSLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&ascii);CHKERRQ(ierr); 4221566a47fSLisandro Dalcin if (ascii) {ierr = PetscViewerASCIIPrintf(viewer," Scheme: %s\n",ssp->type_name);CHKERRQ(ierr);} 423ef7bb5aaSJed Brown PetscFunctionReturn(0); 424ef7bb5aaSJed Brown } 425ef7bb5aaSJed Brown 426ef7bb5aaSJed Brown /* ------------------------------------------------------------ */ 427ef7bb5aaSJed Brown 428ef7bb5aaSJed Brown /*MC 4298ab3e0fcSJed Brown TSSSP - Explicit strong stability preserving ODE solver 4308ab3e0fcSJed Brown 4318ab3e0fcSJed Brown Most hyperbolic conservation laws have exact solutions that are total variation diminishing (TVD) or total variation 4328ab3e0fcSJed Brown bounded (TVB) although these solutions often contain discontinuities. Spatial discretizations such as Godunov's 4338ab3e0fcSJed Brown scheme and high-resolution finite volume methods (TVD limiters, ENO/WENO) are designed to preserve these properties, 4348ab3e0fcSJed Brown but they are usually formulated using a forward Euler time discretization or by coupling the space and time 4358ab3e0fcSJed Brown discretization as in the classical Lax-Wendroff scheme. When the space and time discretization is coupled, it is very 4368ab3e0fcSJed Brown difficult to produce schemes with high temporal accuracy while preserving TVD properties. An alternative is the 4378ab3e0fcSJed Brown semidiscrete formulation where we choose a spatial discretization that is TVD with forward Euler and then choose a 4388ab3e0fcSJed Brown time discretization that preserves the TVD property. Such integrators are called strong stability preserving (SSP). 4398ab3e0fcSJed Brown 4408ab3e0fcSJed Brown Let c_eff be the minimum number of function evaluations required to step as far as one step of forward Euler while 4418ab3e0fcSJed Brown still being SSP. Some theoretical bounds 4428ab3e0fcSJed Brown 4438ab3e0fcSJed Brown 1. There are no explicit methods with c_eff > 1. 4440085e20eSJed Brown 4458ab3e0fcSJed Brown 2. There are no explicit methods beyond order 4 (for nonlinear problems) and c_eff > 0. 4460085e20eSJed Brown 4478ab3e0fcSJed Brown 3. There are no implicit methods with order greater than 1 and c_eff > 2. 4488ab3e0fcSJed Brown 4498ab3e0fcSJed Brown This integrator provides Runge-Kutta methods of order 2, 3, and 4 with maximal values of c_eff. More stages allows 4508ab3e0fcSJed Brown for larger values of c_eff which improves efficiency. These implementations are low-memory and only use 2 or 3 work 4518ab3e0fcSJed Brown vectors regardless of the total number of stages, so e.g. 25-stage 3rd order methods may be an excellent choice. 4528ab3e0fcSJed Brown 4538ab3e0fcSJed Brown Methods can be chosen with -ts_ssp_type {rks2,rks3,rk104} 4548ab3e0fcSJed Brown 4558ab3e0fcSJed Brown rks2: Second order methods with any number s>1 of stages. c_eff = (s-1)/s 4568ab3e0fcSJed Brown 4578ab3e0fcSJed Brown rks3: Third order methods with s=n^2 stages, n>1. c_eff = (s-n)/s 4588ab3e0fcSJed Brown 4598ab3e0fcSJed Brown rk104: A 10-stage fourth order method. c_eff = 0.6 460ef7bb5aaSJed Brown 461ef7bb5aaSJed Brown Level: beginner 462ef7bb5aaSJed Brown 4637b6bb2c6SJed Brown References: 46496a0c994SBarry Smith + 1. - Ketcheson, Highly efficient strong stability preserving Runge Kutta methods with low storage implementations, SISC, 2008. 46596a0c994SBarry Smith - 2. - Gottlieb, Ketcheson, and Shu, High order strong stability preserving time discretizations, J Scientific Computing, 2009. 4667b6bb2c6SJed Brown 467ef7bb5aaSJed Brown .seealso: TSCreate(), TS, TSSetType() 468ef7bb5aaSJed Brown 469ef7bb5aaSJed Brown M*/ 4708cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_SSP(TS ts) 471ef7bb5aaSJed Brown { 472ef7bb5aaSJed Brown TS_SSP *ssp; 473ef7bb5aaSJed Brown PetscErrorCode ierr; 474ef7bb5aaSJed Brown 475ef7bb5aaSJed Brown PetscFunctionBegin; 476787849ffSJed Brown ierr = TSSSPInitializePackage();CHKERRQ(ierr); 477ef7bb5aaSJed Brown 478ef7bb5aaSJed Brown ts->ops->setup = TSSetUp_SSP; 479ef7bb5aaSJed Brown ts->ops->step = TSStep_SSP; 480277b19d0SLisandro Dalcin ts->ops->reset = TSReset_SSP; 481ef7bb5aaSJed Brown ts->ops->destroy = TSDestroy_SSP; 482ef7bb5aaSJed Brown ts->ops->setfromoptions = TSSetFromOptions_SSP; 483ef7bb5aaSJed Brown ts->ops->view = TSView_SSP; 484ef7bb5aaSJed Brown 485b00a9115SJed Brown ierr = PetscNewLog(ts,&ssp);CHKERRQ(ierr); 486ef7bb5aaSJed Brown ts->data = (void*)ssp; 487ef7bb5aaSJed Brown 488bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetType_C",TSSSPGetType_SSP);CHKERRQ(ierr); 489bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetType_C",TSSSPSetType_SSP);CHKERRQ(ierr); 490bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetNumStages_C",TSSSPGetNumStages_SSP);CHKERRQ(ierr); 491bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetNumStages_C",TSSSPSetNumStages_SSP);CHKERRQ(ierr); 492815f1ad5SJed Brown 493ef7bb5aaSJed Brown ierr = TSSSPSetType(ts,TSSSPRKS2);CHKERRQ(ierr); 494ef7bb5aaSJed Brown ssp->nstages = 5; 495ef7bb5aaSJed Brown PetscFunctionReturn(0); 496ef7bb5aaSJed Brown } 497787849ffSJed Brown 498787849ffSJed Brown /*@C 499787849ffSJed Brown TSSSPInitializePackage - This function initializes everything in the TSSSP package. It is called 500787849ffSJed Brown from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_SSP() 501787849ffSJed Brown when using static libraries. 502787849ffSJed Brown 503787849ffSJed Brown Level: developer 504787849ffSJed Brown 505787849ffSJed Brown .keywords: TS, TSSSP, initialize, package 506787849ffSJed Brown .seealso: PetscInitialize() 507787849ffSJed Brown @*/ 508787849ffSJed Brown PetscErrorCode TSSSPInitializePackage(void) 509787849ffSJed Brown { 510787849ffSJed Brown PetscErrorCode ierr; 511787849ffSJed Brown 512787849ffSJed Brown PetscFunctionBegin; 513787849ffSJed Brown if (TSSSPPackageInitialized) PetscFunctionReturn(0); 514787849ffSJed Brown TSSSPPackageInitialized = PETSC_TRUE; 515787849ffSJed Brown ierr = PetscFunctionListAdd(&TSSSPList,TSSSPRKS2, TSSSPStep_RK_2);CHKERRQ(ierr); 516787849ffSJed Brown ierr = PetscFunctionListAdd(&TSSSPList,TSSSPRKS3, TSSSPStep_RK_3);CHKERRQ(ierr); 517787849ffSJed Brown ierr = PetscFunctionListAdd(&TSSSPList,TSSSPRK104,TSSSPStep_RK_10_4);CHKERRQ(ierr); 518787849ffSJed Brown ierr = PetscRegisterFinalize(TSSSPFinalizePackage);CHKERRQ(ierr); 519787849ffSJed Brown PetscFunctionReturn(0); 520787849ffSJed Brown } 521787849ffSJed Brown 522787849ffSJed Brown /*@C 523787849ffSJed Brown TSSSPFinalizePackage - This function destroys everything in the TSSSP package. It is 524787849ffSJed Brown called from PetscFinalize(). 525787849ffSJed Brown 526787849ffSJed Brown Level: developer 527787849ffSJed Brown 528787849ffSJed Brown .keywords: Petsc, destroy, package 529787849ffSJed Brown .seealso: PetscFinalize() 530787849ffSJed Brown @*/ 531787849ffSJed Brown PetscErrorCode TSSSPFinalizePackage(void) 532787849ffSJed Brown { 533787849ffSJed Brown PetscErrorCode ierr; 534787849ffSJed Brown 535787849ffSJed Brown PetscFunctionBegin; 536787849ffSJed Brown TSSSPPackageInitialized = PETSC_FALSE; 537787849ffSJed Brown ierr = PetscFunctionListDestroy(&TSSSPList);CHKERRQ(ierr); 538787849ffSJed Brown PetscFunctionReturn(0); 539787849ffSJed Brown } 540