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 19ef7bb5aaSJed Brown #undef __FUNCT__ 2005d28066SBarry Smith #define __FUNCT__ "TSSSPGetWorkVectors" 2105d28066SBarry Smith static PetscErrorCode TSSSPGetWorkVectors(TS ts,PetscInt n,Vec **work) 22ef7bb5aaSJed Brown { 23ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 24ef7bb5aaSJed Brown PetscErrorCode ierr; 25ef7bb5aaSJed Brown 26ef7bb5aaSJed Brown PetscFunctionBegin; 27e32f2f54SBarry Smith if (ssp->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Work vectors already gotten"); 28ef7bb5aaSJed Brown if (ssp->nwork < n) { 29ef7bb5aaSJed Brown if (ssp->nwork > 0) { 30574deadeSBarry Smith ierr = VecDestroyVecs(ssp->nwork,&ssp->work);CHKERRQ(ierr); 31ef7bb5aaSJed Brown } 32ef7bb5aaSJed Brown ierr = VecDuplicateVecs(ts->vec_sol,n,&ssp->work);CHKERRQ(ierr); 33ef7bb5aaSJed Brown ssp->nwork = n; 34ef7bb5aaSJed Brown } 35ef7bb5aaSJed Brown *work = ssp->work; 36ef7bb5aaSJed Brown ssp->workout = PETSC_TRUE; 37ef7bb5aaSJed Brown PetscFunctionReturn(0); 38ef7bb5aaSJed Brown } 39ef7bb5aaSJed Brown 40ef7bb5aaSJed Brown #undef __FUNCT__ 4105d28066SBarry Smith #define __FUNCT__ "TSSSPRestoreWorkVectors" 4205d28066SBarry Smith static PetscErrorCode TSSSPRestoreWorkVectors(TS ts,PetscInt n,Vec **work) 43ef7bb5aaSJed Brown { 44ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 45ef7bb5aaSJed Brown 46ef7bb5aaSJed Brown PetscFunctionBegin; 47e32f2f54SBarry Smith if (!ssp->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Work vectors have not been gotten"); 48e32f2f54SBarry Smith if (*work != ssp->work) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong work vectors checked out"); 49ef7bb5aaSJed Brown ssp->workout = PETSC_FALSE; 500298fd71SBarry Smith *work = NULL; 51ef7bb5aaSJed Brown PetscFunctionReturn(0); 52ef7bb5aaSJed Brown } 53ef7bb5aaSJed Brown 54ef7bb5aaSJed Brown #undef __FUNCT__ 5505d28066SBarry Smith #define __FUNCT__ "TSSSPStep_RK_2" 56815f1ad5SJed Brown /*MC 57815f1ad5SJed Brown TSSSPRKS2 - Optimal second order SSP Runge-Kutta method, low-storage, c_eff=(s-1)/s 58815f1ad5SJed Brown 59815f1ad5SJed Brown Pseudocode 2 of Ketcheson 2008 60815f1ad5SJed Brown 61b330ce4dSSatish Balay Level: beginner 62b330ce4dSSatish Balay 63815f1ad5SJed Brown .seealso: TSSSP, TSSSPSetType(), TSSSPSetNumStages() 64815f1ad5SJed Brown M*/ 6505d28066SBarry Smith static PetscErrorCode TSSSPStep_RK_2(TS ts,PetscReal t0,PetscReal dt,Vec sol) 66ef7bb5aaSJed Brown { 67ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 68ef7bb5aaSJed Brown Vec *work,F; 69ef7bb5aaSJed Brown PetscInt i,s; 70ef7bb5aaSJed Brown PetscErrorCode ierr; 71ef7bb5aaSJed Brown 72ef7bb5aaSJed Brown PetscFunctionBegin; 73ef7bb5aaSJed Brown s = ssp->nstages; 7405d28066SBarry Smith ierr = TSSSPGetWorkVectors(ts,2,&work);CHKERRQ(ierr); 75ef7bb5aaSJed Brown F = work[1]; 76ef7bb5aaSJed Brown ierr = VecCopy(sol,work[0]);CHKERRQ(ierr); 77ef7bb5aaSJed Brown for (i=0; i<s-1; i++) { 78b8123daeSJed Brown PetscReal stage_time = t0+dt*(i/(s-1.)); 79b8123daeSJed Brown ierr = TSPreStage(ts,stage_time);CHKERRQ(ierr); 80b8123daeSJed Brown ierr = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr); 81ef7bb5aaSJed Brown ierr = VecAXPY(work[0],dt/(s-1.),F);CHKERRQ(ierr); 82ef7bb5aaSJed Brown } 83ef7bb5aaSJed Brown ierr = TSComputeRHSFunction(ts,t0+dt,work[0],F);CHKERRQ(ierr); 84ef7bb5aaSJed Brown ierr = VecAXPBYPCZ(sol,(s-1.)/s,dt/s,1./s,work[0],F);CHKERRQ(ierr); 8505d28066SBarry Smith ierr = TSSSPRestoreWorkVectors(ts,2,&work);CHKERRQ(ierr); 86ef7bb5aaSJed Brown PetscFunctionReturn(0); 87ef7bb5aaSJed Brown } 88ef7bb5aaSJed Brown 89ef7bb5aaSJed Brown #undef __FUNCT__ 9005d28066SBarry Smith #define __FUNCT__ "TSSSPStep_RK_3" 91815f1ad5SJed Brown /*MC 92aaf9cf16SJed Brown TSSSPRKS3 - Optimal third order SSP Runge-Kutta, low-storage, c_eff=(PetscSqrtReal(s)-1)/PetscSqrtReal(s), where PetscSqrtReal(s) is an integer 93815f1ad5SJed Brown 94815f1ad5SJed Brown Pseudocode 2 of Ketcheson 2008 95815f1ad5SJed Brown 96b330ce4dSSatish Balay Level: beginner 97b330ce4dSSatish Balay 98815f1ad5SJed Brown .seealso: TSSSP, TSSSPSetType(), TSSSPSetNumStages() 99815f1ad5SJed Brown M*/ 10005d28066SBarry Smith static PetscErrorCode TSSSPStep_RK_3(TS ts,PetscReal t0,PetscReal dt,Vec sol) 101ef7bb5aaSJed Brown { 102ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 103ef7bb5aaSJed Brown Vec *work,F; 104ef7bb5aaSJed Brown PetscInt i,s,n,r; 105b8123daeSJed Brown PetscReal c,stage_time; 106ef7bb5aaSJed Brown PetscErrorCode ierr; 107ef7bb5aaSJed Brown 108ef7bb5aaSJed Brown PetscFunctionBegin; 109ef7bb5aaSJed Brown s = ssp->nstages; 110aaf9cf16SJed Brown n = (PetscInt)(PetscSqrtReal((PetscReal)s)+0.001); 111ef7bb5aaSJed Brown r = s-n; 112e32f2f54SBarry 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); 11305d28066SBarry Smith ierr = TSSSPGetWorkVectors(ts,3,&work);CHKERRQ(ierr); 114ef7bb5aaSJed Brown F = work[2]; 115ef7bb5aaSJed Brown ierr = VecCopy(sol,work[0]);CHKERRQ(ierr); 116ef7bb5aaSJed Brown for (i=0; i<(n-1)*(n-2)/2; 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 ierr = VecCopy(work[0],work[1]);CHKERRQ(ierr); 124ef7bb5aaSJed Brown for (; i<n*(n+1)/2-1; i++) { 125ef7bb5aaSJed Brown c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n); 126b8123daeSJed Brown stage_time = t0+c*dt; 127b8123daeSJed Brown ierr = TSPreStage(ts,stage_time);CHKERRQ(ierr); 128b8123daeSJed Brown ierr = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr); 129ef7bb5aaSJed Brown ierr = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr); 130ef7bb5aaSJed Brown } 131ef7bb5aaSJed Brown { 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 = 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); 137ef7bb5aaSJed Brown i++; 138ef7bb5aaSJed Brown } 139ef7bb5aaSJed Brown for (; i<s; i++) { 140ef7bb5aaSJed Brown c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n); 141b8123daeSJed Brown stage_time = t0+c*dt; 142b8123daeSJed Brown ierr = TSPreStage(ts,stage_time);CHKERRQ(ierr); 143b8123daeSJed Brown ierr = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr); 144ef7bb5aaSJed Brown ierr = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr); 145ef7bb5aaSJed Brown } 146ef7bb5aaSJed Brown ierr = VecCopy(work[0],sol);CHKERRQ(ierr); 14705d28066SBarry Smith ierr = TSSSPRestoreWorkVectors(ts,3,&work);CHKERRQ(ierr); 148ef7bb5aaSJed Brown PetscFunctionReturn(0); 149ef7bb5aaSJed Brown } 150ef7bb5aaSJed Brown 151ef7bb5aaSJed Brown #undef __FUNCT__ 15205d28066SBarry Smith #define __FUNCT__ "TSSSPStep_RK_10_4" 153815f1ad5SJed Brown /*MC 154b330ce4dSSatish Balay TSSSPRKS104 - Optimal fourth order SSP Runge-Kutta, low-storage (2N), c_eff=0.6 155815f1ad5SJed Brown 156815f1ad5SJed Brown SSPRK(10,4), Pseudocode 3 of Ketcheson 2008 157815f1ad5SJed Brown 158b330ce4dSSatish Balay Level: beginner 159b330ce4dSSatish Balay 160815f1ad5SJed Brown .seealso: TSSSP, TSSSPSetType() 161815f1ad5SJed Brown M*/ 16205d28066SBarry Smith static PetscErrorCode TSSSPStep_RK_10_4(TS ts,PetscReal t0,PetscReal dt,Vec sol) 163ef7bb5aaSJed Brown { 164ef7bb5aaSJed Brown const PetscReal c[10] = {0, 1./6, 2./6, 3./6, 4./6, 2./6, 3./6, 4./6, 5./6, 1}; 165ef7bb5aaSJed Brown Vec *work,F; 1665a586d82SBarry Smith PetscInt i; 167b8123daeSJed Brown PetscReal stage_time; 168ef7bb5aaSJed Brown PetscErrorCode ierr; 169ef7bb5aaSJed Brown 170ef7bb5aaSJed Brown PetscFunctionBegin; 17105d28066SBarry Smith ierr = TSSSPGetWorkVectors(ts,3,&work);CHKERRQ(ierr); 172ef7bb5aaSJed Brown F = work[2]; 173ef7bb5aaSJed Brown ierr = VecCopy(sol,work[0]);CHKERRQ(ierr); 174ef7bb5aaSJed Brown for (i=0; i<5; i++) { 175b8123daeSJed Brown stage_time = t0+c[i]*dt; 176b8123daeSJed Brown ierr = TSPreStage(ts,stage_time);CHKERRQ(ierr); 177b8123daeSJed Brown ierr = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr); 178ef7bb5aaSJed Brown ierr = VecAXPY(work[0],dt/6,F);CHKERRQ(ierr); 179ef7bb5aaSJed Brown } 180ef7bb5aaSJed Brown ierr = VecAXPBYPCZ(work[1],1./25,9./25,0,sol,work[0]);CHKERRQ(ierr); 181ef7bb5aaSJed Brown ierr = VecAXPBY(work[0],15,-5,work[1]);CHKERRQ(ierr); 182ef7bb5aaSJed Brown for (; i<9; i++) { 183b8123daeSJed Brown stage_time = t0+c[i]*dt; 184b8123daeSJed Brown ierr = TSPreStage(ts,stage_time);CHKERRQ(ierr); 185b8123daeSJed Brown ierr = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr); 186ef7bb5aaSJed Brown ierr = VecAXPY(work[0],dt/6,F);CHKERRQ(ierr); 187ef7bb5aaSJed Brown } 188b8123daeSJed Brown stage_time = t0+dt; 189b8123daeSJed Brown ierr = TSPreStage(ts,stage_time);CHKERRQ(ierr); 190b8123daeSJed Brown ierr = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr); 191ef7bb5aaSJed Brown ierr = VecAXPBYPCZ(work[1],3./5,dt/10,1,work[0],F);CHKERRQ(ierr); 192ef7bb5aaSJed Brown ierr = VecCopy(work[1],sol);CHKERRQ(ierr); 19305d28066SBarry Smith ierr = TSSSPRestoreWorkVectors(ts,3,&work);CHKERRQ(ierr); 194ef7bb5aaSJed Brown PetscFunctionReturn(0); 195ef7bb5aaSJed Brown } 196ef7bb5aaSJed Brown 197ef7bb5aaSJed Brown 198ef7bb5aaSJed Brown #undef __FUNCT__ 199ef7bb5aaSJed Brown #define __FUNCT__ "TSSetUp_SSP" 200ef7bb5aaSJed Brown static PetscErrorCode TSSetUp_SSP(TS ts) 201ef7bb5aaSJed Brown { 2021566a47fSLisandro Dalcin PetscErrorCode ierr; 203ef7bb5aaSJed Brown 204ef7bb5aaSJed Brown PetscFunctionBegin; 205*3dd83b38SBarry Smith ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr); 2061566a47fSLisandro Dalcin ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 2071566a47fSLisandro Dalcin ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 208ef7bb5aaSJed Brown PetscFunctionReturn(0); 209ef7bb5aaSJed Brown } 210ef7bb5aaSJed Brown 211ef7bb5aaSJed Brown #undef __FUNCT__ 212ef7bb5aaSJed Brown #define __FUNCT__ "TSStep_SSP" 213193ac0bcSJed Brown static PetscErrorCode TSStep_SSP(TS ts) 214ef7bb5aaSJed Brown { 215ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 216ef7bb5aaSJed Brown Vec sol = ts->vec_sol; 2171566a47fSLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 2181566a47fSLisandro Dalcin PetscReal next_time_step = ts->time_step; 219ef7bb5aaSJed Brown PetscErrorCode ierr; 220ef7bb5aaSJed Brown 221ef7bb5aaSJed Brown PetscFunctionBegin; 222186e87acSLisandro Dalcin ierr = (*ssp->onestep)(ts,ts->ptime,ts->time_step,sol);CHKERRQ(ierr); 2231566a47fSLisandro Dalcin ierr = TSPostStage(ts,ts->ptime,0,&sol);CHKERRQ(ierr); 2241566a47fSLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,ts->ptime+ts->time_step,sol,&stageok);CHKERRQ(ierr); 2251566a47fSLisandro Dalcin if(!stageok) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);} 2261566a47fSLisandro Dalcin 2271566a47fSLisandro Dalcin ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 2281566a47fSLisandro Dalcin if (!accept) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);} 2291566a47fSLisandro Dalcin 230186e87acSLisandro Dalcin ts->ptime += ts->time_step; 2311566a47fSLisandro Dalcin ts->time_step = next_time_step; 232ef7bb5aaSJed Brown PetscFunctionReturn(0); 233ef7bb5aaSJed Brown } 234ef7bb5aaSJed Brown /*------------------------------------------------------------*/ 2351566a47fSLisandro Dalcin 236ef7bb5aaSJed Brown #undef __FUNCT__ 237277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_SSP" 238277b19d0SLisandro Dalcin static PetscErrorCode TSReset_SSP(TS ts) 239277b19d0SLisandro Dalcin { 240277b19d0SLisandro Dalcin TS_SSP *ssp = (TS_SSP*)ts->data; 241277b19d0SLisandro Dalcin PetscErrorCode ierr; 242277b19d0SLisandro Dalcin 243277b19d0SLisandro Dalcin PetscFunctionBegin; 244277b19d0SLisandro Dalcin if (ssp->work) {ierr = VecDestroyVecs(ssp->nwork,&ssp->work);CHKERRQ(ierr);} 245277b19d0SLisandro Dalcin ssp->nwork = 0; 246277b19d0SLisandro Dalcin ssp->workout = PETSC_FALSE; 247277b19d0SLisandro Dalcin PetscFunctionReturn(0); 248277b19d0SLisandro Dalcin } 249277b19d0SLisandro Dalcin 250277b19d0SLisandro Dalcin #undef __FUNCT__ 251ef7bb5aaSJed Brown #define __FUNCT__ "TSDestroy_SSP" 252ef7bb5aaSJed Brown static PetscErrorCode TSDestroy_SSP(TS ts) 253ef7bb5aaSJed Brown { 254815f1ad5SJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 255ef7bb5aaSJed Brown PetscErrorCode ierr; 256ef7bb5aaSJed Brown 257ef7bb5aaSJed Brown PetscFunctionBegin; 258277b19d0SLisandro Dalcin ierr = TSReset_SSP(ts);CHKERRQ(ierr); 2595164425aSJed Brown ierr = PetscFree(ssp->type_name);CHKERRQ(ierr); 260c31cb41cSBarry Smith ierr = PetscFree(ts->data);CHKERRQ(ierr); 261bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetType_C",NULL);CHKERRQ(ierr); 262bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetType_C",NULL);CHKERRQ(ierr); 263bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetNumStages_C",NULL);CHKERRQ(ierr); 264bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetNumStages_C",NULL);CHKERRQ(ierr); 265ef7bb5aaSJed Brown PetscFunctionReturn(0); 266ef7bb5aaSJed Brown } 267ef7bb5aaSJed Brown /*------------------------------------------------------------*/ 268ef7bb5aaSJed Brown 269ef7bb5aaSJed Brown #undef __FUNCT__ 270ef7bb5aaSJed Brown #define __FUNCT__ "TSSSPSetType" 271815f1ad5SJed Brown /*@C 272815f1ad5SJed Brown TSSSPSetType - set the SSP time integration scheme to use 273815f1ad5SJed Brown 274815f1ad5SJed Brown Logically Collective 275815f1ad5SJed Brown 276815f1ad5SJed Brown Input Arguments: 277815f1ad5SJed Brown ts - time stepping object 278815f1ad5SJed Brown type - type of scheme to use 279815f1ad5SJed Brown 280815f1ad5SJed Brown Options Database Keys: 281815f1ad5SJed Brown -ts_ssp_type <rks2>: Type of SSP method (one of) rks2 rks3 rk104 282815f1ad5SJed Brown -ts_ssp_nstages <5>: Number of stages 283815f1ad5SJed Brown 284815f1ad5SJed Brown Level: beginner 285815f1ad5SJed Brown 286815f1ad5SJed Brown .seealso: TSSSP, TSSSPGetType(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104 287815f1ad5SJed Brown @*/ 28819fd82e9SBarry Smith PetscErrorCode TSSSPSetType(TS ts,TSSSPType type) 289815f1ad5SJed Brown { 290815f1ad5SJed Brown PetscErrorCode ierr; 291815f1ad5SJed Brown 292815f1ad5SJed Brown PetscFunctionBegin; 293815f1ad5SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 29419fd82e9SBarry Smith ierr = PetscTryMethod(ts,"TSSSPSetType_C",(TS,TSSSPType),(ts,type));CHKERRQ(ierr); 295815f1ad5SJed Brown PetscFunctionReturn(0); 296815f1ad5SJed Brown } 297815f1ad5SJed Brown 298815f1ad5SJed Brown #undef __FUNCT__ 299815f1ad5SJed Brown #define __FUNCT__ "TSSSPGetType" 300815f1ad5SJed Brown /*@C 301815f1ad5SJed Brown TSSSPGetType - get the SSP time integration scheme 302815f1ad5SJed Brown 303815f1ad5SJed Brown Logically Collective 304815f1ad5SJed Brown 305815f1ad5SJed Brown Input Argument: 306815f1ad5SJed Brown ts - time stepping object 307815f1ad5SJed Brown 308815f1ad5SJed Brown Output Argument: 309815f1ad5SJed Brown type - type of scheme being used 310815f1ad5SJed Brown 311815f1ad5SJed Brown Level: beginner 312815f1ad5SJed Brown 313815f1ad5SJed Brown .seealso: TSSSP, TSSSPSettype(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104 314815f1ad5SJed Brown @*/ 31519fd82e9SBarry Smith PetscErrorCode TSSSPGetType(TS ts,TSSSPType *type) 316815f1ad5SJed Brown { 317815f1ad5SJed Brown PetscErrorCode ierr; 318815f1ad5SJed Brown 319815f1ad5SJed Brown PetscFunctionBegin; 320815f1ad5SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 321163d334eSBarry Smith ierr = PetscUseMethod(ts,"TSSSPGetType_C",(TS,TSSSPType*),(ts,type));CHKERRQ(ierr); 322815f1ad5SJed Brown PetscFunctionReturn(0); 323815f1ad5SJed Brown } 324815f1ad5SJed Brown 325815f1ad5SJed Brown #undef __FUNCT__ 326815f1ad5SJed Brown #define __FUNCT__ "TSSSPSetNumStages" 327815f1ad5SJed Brown /*@ 328815f1ad5SJed Brown TSSSPSetNumStages - set the number of stages to use with the SSP method 329815f1ad5SJed Brown 330815f1ad5SJed Brown Logically Collective 331815f1ad5SJed Brown 332815f1ad5SJed Brown Input Arguments: 333815f1ad5SJed Brown ts - time stepping object 334815f1ad5SJed Brown nstages - number of stages 335815f1ad5SJed Brown 336815f1ad5SJed Brown Options Database Keys: 337815f1ad5SJed Brown -ts_ssp_type <rks2>: NumStages of SSP method (one of) rks2 rks3 rk104 338815f1ad5SJed Brown -ts_ssp_nstages <5>: Number of stages 339815f1ad5SJed Brown 340815f1ad5SJed Brown Level: beginner 341815f1ad5SJed Brown 342815f1ad5SJed Brown .seealso: TSSSP, TSSSPGetNumStages(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104 343815f1ad5SJed Brown @*/ 344815f1ad5SJed Brown PetscErrorCode TSSSPSetNumStages(TS ts,PetscInt nstages) 345815f1ad5SJed Brown { 346815f1ad5SJed Brown PetscErrorCode ierr; 347815f1ad5SJed Brown 348815f1ad5SJed Brown PetscFunctionBegin; 349815f1ad5SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 350815f1ad5SJed Brown ierr = PetscTryMethod(ts,"TSSSPSetNumStages_C",(TS,PetscInt),(ts,nstages));CHKERRQ(ierr); 351815f1ad5SJed Brown PetscFunctionReturn(0); 352815f1ad5SJed Brown } 353815f1ad5SJed Brown 354815f1ad5SJed Brown #undef __FUNCT__ 355815f1ad5SJed Brown #define __FUNCT__ "TSSSPGetNumStages" 356815f1ad5SJed Brown /*@ 357815f1ad5SJed Brown TSSSPGetNumStages - get the number of stages in the SSP time integration scheme 358815f1ad5SJed Brown 359815f1ad5SJed Brown Logically Collective 360815f1ad5SJed Brown 361815f1ad5SJed Brown Input Argument: 362815f1ad5SJed Brown ts - time stepping object 363815f1ad5SJed Brown 364815f1ad5SJed Brown Output Argument: 365815f1ad5SJed Brown nstages - number of stages 366815f1ad5SJed Brown 367815f1ad5SJed Brown Level: beginner 368815f1ad5SJed Brown 369815f1ad5SJed Brown .seealso: TSSSP, TSSSPGetType(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104 370815f1ad5SJed Brown @*/ 371815f1ad5SJed Brown PetscErrorCode TSSSPGetNumStages(TS ts,PetscInt *nstages) 372815f1ad5SJed Brown { 373815f1ad5SJed Brown PetscErrorCode ierr; 374815f1ad5SJed Brown 375815f1ad5SJed Brown PetscFunctionBegin; 376815f1ad5SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 377163d334eSBarry Smith ierr = PetscUseMethod(ts,"TSSSPGetNumStages_C",(TS,PetscInt*),(ts,nstages));CHKERRQ(ierr); 378815f1ad5SJed Brown PetscFunctionReturn(0); 379815f1ad5SJed Brown } 380815f1ad5SJed Brown 381815f1ad5SJed Brown #undef __FUNCT__ 382815f1ad5SJed Brown #define __FUNCT__ "TSSSPSetType_SSP" 383cc2e6a90SBarry Smith static PetscErrorCode TSSSPSetType_SSP(TS ts,TSSSPType type) 384ef7bb5aaSJed Brown { 385ef7bb5aaSJed Brown PetscErrorCode ierr,(*r)(TS,PetscReal,PetscReal,Vec); 386ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 387ef7bb5aaSJed Brown 388ef7bb5aaSJed Brown PetscFunctionBegin; 3891c9cd337SJed Brown ierr = PetscFunctionListFind(TSSSPList,type,&r);CHKERRQ(ierr); 390e32f2f54SBarry Smith if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TS_SSP type %s given",type); 391ef7bb5aaSJed Brown ssp->onestep = r; 3925164425aSJed Brown ierr = PetscFree(ssp->type_name);CHKERRQ(ierr); 3935164425aSJed Brown ierr = PetscStrallocpy(type,&ssp->type_name);CHKERRQ(ierr); 394ef7bb5aaSJed Brown PetscFunctionReturn(0); 395ef7bb5aaSJed Brown } 396815f1ad5SJed Brown #undef __FUNCT__ 397815f1ad5SJed Brown #define __FUNCT__ "TSSSPGetType_SSP" 398cc2e6a90SBarry Smith static PetscErrorCode TSSSPGetType_SSP(TS ts,TSSSPType *type) 399815f1ad5SJed Brown { 400815f1ad5SJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 401815f1ad5SJed Brown 402815f1ad5SJed Brown PetscFunctionBegin; 4035164425aSJed Brown *type = ssp->type_name; 404815f1ad5SJed Brown PetscFunctionReturn(0); 405815f1ad5SJed Brown } 406815f1ad5SJed Brown #undef __FUNCT__ 407815f1ad5SJed Brown #define __FUNCT__ "TSSSPSetNumStages_SSP" 408cc2e6a90SBarry Smith static PetscErrorCode TSSSPSetNumStages_SSP(TS ts,PetscInt nstages) 409815f1ad5SJed Brown { 410815f1ad5SJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 411815f1ad5SJed Brown 412815f1ad5SJed Brown PetscFunctionBegin; 413815f1ad5SJed Brown ssp->nstages = nstages; 414815f1ad5SJed Brown PetscFunctionReturn(0); 415815f1ad5SJed Brown } 416815f1ad5SJed Brown #undef __FUNCT__ 417815f1ad5SJed Brown #define __FUNCT__ "TSSSPGetNumStages_SSP" 418cc2e6a90SBarry Smith static PetscErrorCode TSSSPGetNumStages_SSP(TS ts,PetscInt *nstages) 419815f1ad5SJed Brown { 420815f1ad5SJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 421815f1ad5SJed Brown 422815f1ad5SJed Brown PetscFunctionBegin; 423815f1ad5SJed Brown *nstages = ssp->nstages; 424815f1ad5SJed Brown PetscFunctionReturn(0); 425815f1ad5SJed Brown } 426ef7bb5aaSJed Brown 427ef7bb5aaSJed Brown #undef __FUNCT__ 428ef7bb5aaSJed Brown #define __FUNCT__ "TSSetFromOptions_SSP" 4294416b707SBarry Smith static PetscErrorCode TSSetFromOptions_SSP(PetscOptionItems *PetscOptionsObject,TS ts) 430ef7bb5aaSJed Brown { 431ef7bb5aaSJed Brown char tname[256] = TSSSPRKS2; 432ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 433ef7bb5aaSJed Brown PetscErrorCode ierr; 434ace3abfcSBarry Smith PetscBool flg; 435ef7bb5aaSJed Brown 436ef7bb5aaSJed Brown PetscFunctionBegin; 437e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SSP ODE solver options");CHKERRQ(ierr); 438ef7bb5aaSJed Brown { 439a264d7a6SBarry Smith ierr = PetscOptionsFList("-ts_ssp_type","Type of SSP method","TSSSPSetType",TSSSPList,tname,tname,sizeof(tname),&flg);CHKERRQ(ierr); 440ef7bb5aaSJed Brown if (flg) { 441ef7bb5aaSJed Brown ierr = TSSSPSetType(ts,tname);CHKERRQ(ierr); 442ef7bb5aaSJed Brown } 4430298fd71SBarry Smith ierr = PetscOptionsInt("-ts_ssp_nstages","Number of stages","TSSSPSetNumStages",ssp->nstages,&ssp->nstages,NULL);CHKERRQ(ierr); 444ef7bb5aaSJed Brown } 445ef7bb5aaSJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 446ef7bb5aaSJed Brown PetscFunctionReturn(0); 447ef7bb5aaSJed Brown } 448ef7bb5aaSJed Brown 449ef7bb5aaSJed Brown #undef __FUNCT__ 450ef7bb5aaSJed Brown #define __FUNCT__ "TSView_SSP" 451ef7bb5aaSJed Brown static PetscErrorCode TSView_SSP(TS ts,PetscViewer viewer) 452ef7bb5aaSJed Brown { 4531566a47fSLisandro Dalcin TS_SSP *ssp = (TS_SSP*)ts->data; 4541566a47fSLisandro Dalcin PetscBool ascii; 4551566a47fSLisandro Dalcin PetscErrorCode ierr; 4561566a47fSLisandro Dalcin 457ef7bb5aaSJed Brown PetscFunctionBegin; 4581566a47fSLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&ascii);CHKERRQ(ierr); 4591566a47fSLisandro Dalcin if (ascii) {ierr = PetscViewerASCIIPrintf(viewer," Scheme: %s\n",ssp->type_name);CHKERRQ(ierr);} 460ef7bb5aaSJed Brown PetscFunctionReturn(0); 461ef7bb5aaSJed Brown } 462ef7bb5aaSJed Brown 463ef7bb5aaSJed Brown /* ------------------------------------------------------------ */ 464ef7bb5aaSJed Brown 465ef7bb5aaSJed Brown /*MC 4668ab3e0fcSJed Brown TSSSP - Explicit strong stability preserving ODE solver 4678ab3e0fcSJed Brown 4688ab3e0fcSJed Brown Most hyperbolic conservation laws have exact solutions that are total variation diminishing (TVD) or total variation 4698ab3e0fcSJed Brown bounded (TVB) although these solutions often contain discontinuities. Spatial discretizations such as Godunov's 4708ab3e0fcSJed Brown scheme and high-resolution finite volume methods (TVD limiters, ENO/WENO) are designed to preserve these properties, 4718ab3e0fcSJed Brown but they are usually formulated using a forward Euler time discretization or by coupling the space and time 4728ab3e0fcSJed Brown discretization as in the classical Lax-Wendroff scheme. When the space and time discretization is coupled, it is very 4738ab3e0fcSJed Brown difficult to produce schemes with high temporal accuracy while preserving TVD properties. An alternative is the 4748ab3e0fcSJed Brown semidiscrete formulation where we choose a spatial discretization that is TVD with forward Euler and then choose a 4758ab3e0fcSJed Brown time discretization that preserves the TVD property. Such integrators are called strong stability preserving (SSP). 4768ab3e0fcSJed Brown 4778ab3e0fcSJed Brown Let c_eff be the minimum number of function evaluations required to step as far as one step of forward Euler while 4788ab3e0fcSJed Brown still being SSP. Some theoretical bounds 4798ab3e0fcSJed Brown 4808ab3e0fcSJed Brown 1. There are no explicit methods with c_eff > 1. 4810085e20eSJed Brown 4828ab3e0fcSJed Brown 2. There are no explicit methods beyond order 4 (for nonlinear problems) and c_eff > 0. 4830085e20eSJed Brown 4848ab3e0fcSJed Brown 3. There are no implicit methods with order greater than 1 and c_eff > 2. 4858ab3e0fcSJed Brown 4868ab3e0fcSJed Brown This integrator provides Runge-Kutta methods of order 2, 3, and 4 with maximal values of c_eff. More stages allows 4878ab3e0fcSJed Brown for larger values of c_eff which improves efficiency. These implementations are low-memory and only use 2 or 3 work 4888ab3e0fcSJed Brown vectors regardless of the total number of stages, so e.g. 25-stage 3rd order methods may be an excellent choice. 4898ab3e0fcSJed Brown 4908ab3e0fcSJed Brown Methods can be chosen with -ts_ssp_type {rks2,rks3,rk104} 4918ab3e0fcSJed Brown 4928ab3e0fcSJed Brown rks2: Second order methods with any number s>1 of stages. c_eff = (s-1)/s 4938ab3e0fcSJed Brown 4948ab3e0fcSJed Brown rks3: Third order methods with s=n^2 stages, n>1. c_eff = (s-n)/s 4958ab3e0fcSJed Brown 4968ab3e0fcSJed Brown rk104: A 10-stage fourth order method. c_eff = 0.6 497ef7bb5aaSJed Brown 498ef7bb5aaSJed Brown Level: beginner 499ef7bb5aaSJed Brown 5007b6bb2c6SJed Brown References: 50196a0c994SBarry Smith + 1. - Ketcheson, Highly efficient strong stability preserving Runge Kutta methods with low storage implementations, SISC, 2008. 50296a0c994SBarry Smith - 2. - Gottlieb, Ketcheson, and Shu, High order strong stability preserving time discretizations, J Scientific Computing, 2009. 5037b6bb2c6SJed Brown 504ef7bb5aaSJed Brown .seealso: TSCreate(), TS, TSSetType() 505ef7bb5aaSJed Brown 506ef7bb5aaSJed Brown M*/ 507ef7bb5aaSJed Brown #undef __FUNCT__ 508ef7bb5aaSJed Brown #define __FUNCT__ "TSCreate_SSP" 5098cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_SSP(TS ts) 510ef7bb5aaSJed Brown { 511ef7bb5aaSJed Brown TS_SSP *ssp; 512ef7bb5aaSJed Brown PetscErrorCode ierr; 513ef7bb5aaSJed Brown 514ef7bb5aaSJed Brown PetscFunctionBegin; 515787849ffSJed Brown ierr = TSSSPInitializePackage();CHKERRQ(ierr); 516ef7bb5aaSJed Brown 517ef7bb5aaSJed Brown ts->ops->setup = TSSetUp_SSP; 518ef7bb5aaSJed Brown ts->ops->step = TSStep_SSP; 519277b19d0SLisandro Dalcin ts->ops->reset = TSReset_SSP; 520ef7bb5aaSJed Brown ts->ops->destroy = TSDestroy_SSP; 521ef7bb5aaSJed Brown ts->ops->setfromoptions = TSSetFromOptions_SSP; 522ef7bb5aaSJed Brown ts->ops->view = TSView_SSP; 523ef7bb5aaSJed Brown 524b00a9115SJed Brown ierr = PetscNewLog(ts,&ssp);CHKERRQ(ierr); 525ef7bb5aaSJed Brown ts->data = (void*)ssp; 526ef7bb5aaSJed Brown 527bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetType_C",TSSSPGetType_SSP);CHKERRQ(ierr); 528bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetType_C",TSSSPSetType_SSP);CHKERRQ(ierr); 529bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetNumStages_C",TSSSPGetNumStages_SSP);CHKERRQ(ierr); 530bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetNumStages_C",TSSSPSetNumStages_SSP);CHKERRQ(ierr); 531815f1ad5SJed Brown 532ef7bb5aaSJed Brown ierr = TSSSPSetType(ts,TSSSPRKS2);CHKERRQ(ierr); 533ef7bb5aaSJed Brown ssp->nstages = 5; 534ef7bb5aaSJed Brown PetscFunctionReturn(0); 535ef7bb5aaSJed Brown } 536787849ffSJed Brown 537787849ffSJed Brown #undef __FUNCT__ 538787849ffSJed Brown #define __FUNCT__ "TSSSPInitializePackage" 539787849ffSJed Brown /*@C 540787849ffSJed Brown TSSSPInitializePackage - This function initializes everything in the TSSSP package. It is called 541787849ffSJed Brown from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_SSP() 542787849ffSJed Brown when using static libraries. 543787849ffSJed Brown 544787849ffSJed Brown Level: developer 545787849ffSJed Brown 546787849ffSJed Brown .keywords: TS, TSSSP, initialize, package 547787849ffSJed Brown .seealso: PetscInitialize() 548787849ffSJed Brown @*/ 549787849ffSJed Brown PetscErrorCode TSSSPInitializePackage(void) 550787849ffSJed Brown { 551787849ffSJed Brown PetscErrorCode ierr; 552787849ffSJed Brown 553787849ffSJed Brown PetscFunctionBegin; 554787849ffSJed Brown if (TSSSPPackageInitialized) PetscFunctionReturn(0); 555787849ffSJed Brown TSSSPPackageInitialized = PETSC_TRUE; 556787849ffSJed Brown ierr = PetscFunctionListAdd(&TSSSPList,TSSSPRKS2, TSSSPStep_RK_2);CHKERRQ(ierr); 557787849ffSJed Brown ierr = PetscFunctionListAdd(&TSSSPList,TSSSPRKS3, TSSSPStep_RK_3);CHKERRQ(ierr); 558787849ffSJed Brown ierr = PetscFunctionListAdd(&TSSSPList,TSSSPRK104,TSSSPStep_RK_10_4);CHKERRQ(ierr); 559787849ffSJed Brown ierr = PetscRegisterFinalize(TSSSPFinalizePackage);CHKERRQ(ierr); 560787849ffSJed Brown PetscFunctionReturn(0); 561787849ffSJed Brown } 562787849ffSJed Brown 563787849ffSJed Brown #undef __FUNCT__ 564787849ffSJed Brown #define __FUNCT__ "TSSSPFinalizePackage" 565787849ffSJed Brown /*@C 566787849ffSJed Brown TSSSPFinalizePackage - This function destroys everything in the TSSSP package. It is 567787849ffSJed Brown called from PetscFinalize(). 568787849ffSJed Brown 569787849ffSJed Brown Level: developer 570787849ffSJed Brown 571787849ffSJed Brown .keywords: Petsc, destroy, package 572787849ffSJed Brown .seealso: PetscFinalize() 573787849ffSJed Brown @*/ 574787849ffSJed Brown PetscErrorCode TSSSPFinalizePackage(void) 575787849ffSJed Brown { 576787849ffSJed Brown PetscErrorCode ierr; 577787849ffSJed Brown 578787849ffSJed Brown PetscFunctionBegin; 579787849ffSJed Brown TSSSPPackageInitialized = PETSC_FALSE; 580787849ffSJed Brown ierr = PetscFunctionListDestroy(&TSSSPList);CHKERRQ(ierr); 581787849ffSJed Brown PetscFunctionReturn(0); 582787849ffSJed Brown } 583