1ef7bb5aaSJed Brown /* 2ef7bb5aaSJed Brown Code for Timestepping with explicit SSP. 3ef7bb5aaSJed Brown */ 4b45d2f2cSJed Brown #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 5ef7bb5aaSJed Brown 6140e18c1SBarry Smith PetscFunctionList TSSSPList = 0; 7ef7bb5aaSJed Brown 8ef7bb5aaSJed Brown typedef struct { 9ef7bb5aaSJed Brown PetscErrorCode (*onestep)(TS,PetscReal,PetscReal,Vec); 105164425aSJed Brown char *type_name; 11ef7bb5aaSJed Brown PetscInt nstages; 12ef7bb5aaSJed Brown Vec *work; 13ef7bb5aaSJed Brown PetscInt nwork; 14ace3abfcSBarry Smith PetscBool workout; 15ef7bb5aaSJed Brown } TS_SSP; 16ef7bb5aaSJed Brown 17ef7bb5aaSJed Brown 18ef7bb5aaSJed Brown #undef __FUNCT__ 1905d28066SBarry Smith #define __FUNCT__ "TSSSPGetWorkVectors" 2005d28066SBarry Smith static PetscErrorCode TSSSPGetWorkVectors(TS ts,PetscInt n,Vec **work) 21ef7bb5aaSJed Brown { 22ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 23ef7bb5aaSJed Brown PetscErrorCode ierr; 24ef7bb5aaSJed Brown 25ef7bb5aaSJed Brown PetscFunctionBegin; 26e32f2f54SBarry Smith if (ssp->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Work vectors already gotten"); 27ef7bb5aaSJed Brown if (ssp->nwork < n) { 28ef7bb5aaSJed Brown if (ssp->nwork > 0) { 29574deadeSBarry Smith ierr = VecDestroyVecs(ssp->nwork,&ssp->work);CHKERRQ(ierr); 30ef7bb5aaSJed Brown } 31ef7bb5aaSJed Brown ierr = VecDuplicateVecs(ts->vec_sol,n,&ssp->work);CHKERRQ(ierr); 32ef7bb5aaSJed Brown ssp->nwork = n; 33ef7bb5aaSJed Brown } 34ef7bb5aaSJed Brown *work = ssp->work; 35ef7bb5aaSJed Brown ssp->workout = PETSC_TRUE; 36ef7bb5aaSJed Brown PetscFunctionReturn(0); 37ef7bb5aaSJed Brown } 38ef7bb5aaSJed Brown 39ef7bb5aaSJed Brown #undef __FUNCT__ 4005d28066SBarry Smith #define __FUNCT__ "TSSSPRestoreWorkVectors" 4105d28066SBarry Smith static PetscErrorCode TSSSPRestoreWorkVectors(TS ts,PetscInt n,Vec **work) 42ef7bb5aaSJed Brown { 43ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 44ef7bb5aaSJed Brown 45ef7bb5aaSJed Brown PetscFunctionBegin; 46e32f2f54SBarry Smith if (!ssp->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Work vectors have not been gotten"); 47e32f2f54SBarry Smith if (*work != ssp->work) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong work vectors checked out"); 48ef7bb5aaSJed Brown ssp->workout = PETSC_FALSE; 490298fd71SBarry Smith *work = NULL; 50ef7bb5aaSJed Brown PetscFunctionReturn(0); 51ef7bb5aaSJed Brown } 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 { 202ef7bb5aaSJed Brown 203ef7bb5aaSJed Brown PetscFunctionBegin; 204ef7bb5aaSJed Brown PetscFunctionReturn(0); 205ef7bb5aaSJed Brown } 206ef7bb5aaSJed Brown 207ef7bb5aaSJed Brown #undef __FUNCT__ 208ef7bb5aaSJed Brown #define __FUNCT__ "TSStep_SSP" 209193ac0bcSJed Brown static PetscErrorCode TSStep_SSP(TS ts) 210ef7bb5aaSJed Brown { 211ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 212ef7bb5aaSJed Brown Vec sol = ts->vec_sol; 213ef7bb5aaSJed Brown PetscErrorCode ierr; 214ef7bb5aaSJed Brown 215ef7bb5aaSJed Brown PetscFunctionBegin; 216b8123daeSJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); 217186e87acSLisandro Dalcin ierr = (*ssp->onestep)(ts,ts->ptime,ts->time_step,sol);CHKERRQ(ierr); 218186e87acSLisandro Dalcin ts->ptime += ts->time_step; 219ef7bb5aaSJed Brown ts->steps++; 220ef7bb5aaSJed Brown PetscFunctionReturn(0); 221ef7bb5aaSJed Brown } 222ef7bb5aaSJed Brown /*------------------------------------------------------------*/ 223ef7bb5aaSJed Brown #undef __FUNCT__ 224277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_SSP" 225277b19d0SLisandro Dalcin static PetscErrorCode TSReset_SSP(TS ts) 226277b19d0SLisandro Dalcin { 227277b19d0SLisandro Dalcin TS_SSP *ssp = (TS_SSP*)ts->data; 228277b19d0SLisandro Dalcin PetscErrorCode ierr; 229277b19d0SLisandro Dalcin 230277b19d0SLisandro Dalcin PetscFunctionBegin; 231277b19d0SLisandro Dalcin if (ssp->work) {ierr = VecDestroyVecs(ssp->nwork,&ssp->work);CHKERRQ(ierr);} 232277b19d0SLisandro Dalcin ssp->nwork = 0; 233277b19d0SLisandro Dalcin ssp->workout = PETSC_FALSE; 234277b19d0SLisandro Dalcin PetscFunctionReturn(0); 235277b19d0SLisandro Dalcin } 236277b19d0SLisandro Dalcin 237277b19d0SLisandro Dalcin #undef __FUNCT__ 238ef7bb5aaSJed Brown #define __FUNCT__ "TSDestroy_SSP" 239ef7bb5aaSJed Brown static PetscErrorCode TSDestroy_SSP(TS ts) 240ef7bb5aaSJed Brown { 241815f1ad5SJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 242ef7bb5aaSJed Brown PetscErrorCode ierr; 243ef7bb5aaSJed Brown 244ef7bb5aaSJed Brown PetscFunctionBegin; 245277b19d0SLisandro Dalcin ierr = TSReset_SSP(ts);CHKERRQ(ierr); 2465164425aSJed Brown ierr = PetscFree(ssp->type_name);CHKERRQ(ierr); 247c31cb41cSBarry Smith ierr = PetscFree(ts->data);CHKERRQ(ierr); 248bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetType_C",NULL);CHKERRQ(ierr); 249bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetType_C",NULL);CHKERRQ(ierr); 250bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetNumStages_C",NULL);CHKERRQ(ierr); 251bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetNumStages_C",NULL);CHKERRQ(ierr); 252ef7bb5aaSJed Brown PetscFunctionReturn(0); 253ef7bb5aaSJed Brown } 254ef7bb5aaSJed Brown /*------------------------------------------------------------*/ 255ef7bb5aaSJed Brown 256ef7bb5aaSJed Brown #undef __FUNCT__ 257ef7bb5aaSJed Brown #define __FUNCT__ "TSSSPSetType" 258815f1ad5SJed Brown /*@C 259815f1ad5SJed Brown TSSSPSetType - set the SSP time integration scheme to use 260815f1ad5SJed Brown 261815f1ad5SJed Brown Logically Collective 262815f1ad5SJed Brown 263815f1ad5SJed Brown Input Arguments: 264815f1ad5SJed Brown ts - time stepping object 265815f1ad5SJed Brown type - type of scheme to use 266815f1ad5SJed Brown 267815f1ad5SJed Brown Options Database Keys: 268815f1ad5SJed Brown -ts_ssp_type <rks2>: Type of SSP method (one of) rks2 rks3 rk104 269815f1ad5SJed Brown -ts_ssp_nstages <5>: Number of stages 270815f1ad5SJed Brown 271815f1ad5SJed Brown Level: beginner 272815f1ad5SJed Brown 273815f1ad5SJed Brown .seealso: TSSSP, TSSSPGetType(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104 274815f1ad5SJed Brown @*/ 27519fd82e9SBarry Smith PetscErrorCode TSSSPSetType(TS ts,TSSSPType type) 276815f1ad5SJed Brown { 277815f1ad5SJed Brown PetscErrorCode ierr; 278815f1ad5SJed Brown 279815f1ad5SJed Brown PetscFunctionBegin; 280815f1ad5SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 28119fd82e9SBarry Smith ierr = PetscTryMethod(ts,"TSSSPSetType_C",(TS,TSSSPType),(ts,type));CHKERRQ(ierr); 282815f1ad5SJed Brown PetscFunctionReturn(0); 283815f1ad5SJed Brown } 284815f1ad5SJed Brown 285815f1ad5SJed Brown #undef __FUNCT__ 286815f1ad5SJed Brown #define __FUNCT__ "TSSSPGetType" 287815f1ad5SJed Brown /*@C 288815f1ad5SJed Brown TSSSPGetType - get the SSP time integration scheme 289815f1ad5SJed Brown 290815f1ad5SJed Brown Logically Collective 291815f1ad5SJed Brown 292815f1ad5SJed Brown Input Argument: 293815f1ad5SJed Brown ts - time stepping object 294815f1ad5SJed Brown 295815f1ad5SJed Brown Output Argument: 296815f1ad5SJed Brown type - type of scheme being used 297815f1ad5SJed Brown 298815f1ad5SJed Brown Level: beginner 299815f1ad5SJed Brown 300815f1ad5SJed Brown .seealso: TSSSP, TSSSPSettype(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104 301815f1ad5SJed Brown @*/ 30219fd82e9SBarry Smith PetscErrorCode TSSSPGetType(TS ts,TSSSPType *type) 303815f1ad5SJed Brown { 304815f1ad5SJed Brown PetscErrorCode ierr; 305815f1ad5SJed Brown 306815f1ad5SJed Brown PetscFunctionBegin; 307815f1ad5SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 30819fd82e9SBarry Smith ierr = PetscTryMethod(ts,"TSSSPGetType_C",(TS,TSSSPType*),(ts,type));CHKERRQ(ierr); 309815f1ad5SJed Brown PetscFunctionReturn(0); 310815f1ad5SJed Brown } 311815f1ad5SJed Brown 312815f1ad5SJed Brown #undef __FUNCT__ 313815f1ad5SJed Brown #define __FUNCT__ "TSSSPSetNumStages" 314815f1ad5SJed Brown /*@ 315815f1ad5SJed Brown TSSSPSetNumStages - set the number of stages to use with the SSP method 316815f1ad5SJed Brown 317815f1ad5SJed Brown Logically Collective 318815f1ad5SJed Brown 319815f1ad5SJed Brown Input Arguments: 320815f1ad5SJed Brown ts - time stepping object 321815f1ad5SJed Brown nstages - number of stages 322815f1ad5SJed Brown 323815f1ad5SJed Brown Options Database Keys: 324815f1ad5SJed Brown -ts_ssp_type <rks2>: NumStages of SSP method (one of) rks2 rks3 rk104 325815f1ad5SJed Brown -ts_ssp_nstages <5>: Number of stages 326815f1ad5SJed Brown 327815f1ad5SJed Brown Level: beginner 328815f1ad5SJed Brown 329815f1ad5SJed Brown .seealso: TSSSP, TSSSPGetNumStages(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104 330815f1ad5SJed Brown @*/ 331815f1ad5SJed Brown PetscErrorCode TSSSPSetNumStages(TS ts,PetscInt nstages) 332815f1ad5SJed Brown { 333815f1ad5SJed Brown PetscErrorCode ierr; 334815f1ad5SJed Brown 335815f1ad5SJed Brown PetscFunctionBegin; 336815f1ad5SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 337815f1ad5SJed Brown ierr = PetscTryMethod(ts,"TSSSPSetNumStages_C",(TS,PetscInt),(ts,nstages));CHKERRQ(ierr); 338815f1ad5SJed Brown PetscFunctionReturn(0); 339815f1ad5SJed Brown } 340815f1ad5SJed Brown 341815f1ad5SJed Brown #undef __FUNCT__ 342815f1ad5SJed Brown #define __FUNCT__ "TSSSPGetNumStages" 343815f1ad5SJed Brown /*@ 344815f1ad5SJed Brown TSSSPGetNumStages - get the number of stages in the SSP time integration scheme 345815f1ad5SJed Brown 346815f1ad5SJed Brown Logically Collective 347815f1ad5SJed Brown 348815f1ad5SJed Brown Input Argument: 349815f1ad5SJed Brown ts - time stepping object 350815f1ad5SJed Brown 351815f1ad5SJed Brown Output Argument: 352815f1ad5SJed Brown nstages - number of stages 353815f1ad5SJed Brown 354815f1ad5SJed Brown Level: beginner 355815f1ad5SJed Brown 356815f1ad5SJed Brown .seealso: TSSSP, TSSSPGetType(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104 357815f1ad5SJed Brown @*/ 358815f1ad5SJed Brown PetscErrorCode TSSSPGetNumStages(TS ts,PetscInt *nstages) 359815f1ad5SJed Brown { 360815f1ad5SJed Brown PetscErrorCode ierr; 361815f1ad5SJed Brown 362815f1ad5SJed Brown PetscFunctionBegin; 363815f1ad5SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 364815f1ad5SJed Brown ierr = PetscTryMethod(ts,"TSSSPGetNumStages_C",(TS,PetscInt*),(ts,nstages));CHKERRQ(ierr); 365815f1ad5SJed Brown PetscFunctionReturn(0); 366815f1ad5SJed Brown } 367815f1ad5SJed Brown 368815f1ad5SJed Brown #undef __FUNCT__ 369815f1ad5SJed Brown #define __FUNCT__ "TSSSPSetType_SSP" 3708cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSSSPSetType_SSP(TS ts,TSSSPType type) 371ef7bb5aaSJed Brown { 372ef7bb5aaSJed Brown PetscErrorCode ierr,(*r)(TS,PetscReal,PetscReal,Vec); 373ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 374ef7bb5aaSJed Brown 375ef7bb5aaSJed Brown PetscFunctionBegin; 376*1c9cd337SJed Brown ierr = PetscFunctionListFind(TSSSPList,type,&r);CHKERRQ(ierr); 377e32f2f54SBarry Smith if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TS_SSP type %s given",type); 378ef7bb5aaSJed Brown ssp->onestep = r; 3795164425aSJed Brown ierr = PetscFree(ssp->type_name);CHKERRQ(ierr); 3805164425aSJed Brown ierr = PetscStrallocpy(type,&ssp->type_name);CHKERRQ(ierr); 381ef7bb5aaSJed Brown PetscFunctionReturn(0); 382ef7bb5aaSJed Brown } 383815f1ad5SJed Brown #undef __FUNCT__ 384815f1ad5SJed Brown #define __FUNCT__ "TSSSPGetType_SSP" 38519fd82e9SBarry Smith PetscErrorCode TSSSPGetType_SSP(TS ts,TSSSPType *type) 386815f1ad5SJed Brown { 387815f1ad5SJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 388815f1ad5SJed Brown 389815f1ad5SJed Brown PetscFunctionBegin; 3905164425aSJed Brown *type = ssp->type_name; 391815f1ad5SJed Brown PetscFunctionReturn(0); 392815f1ad5SJed Brown } 393815f1ad5SJed Brown #undef __FUNCT__ 394815f1ad5SJed Brown #define __FUNCT__ "TSSSPSetNumStages_SSP" 395815f1ad5SJed Brown PetscErrorCode TSSSPSetNumStages_SSP(TS ts,PetscInt nstages) 396815f1ad5SJed Brown { 397815f1ad5SJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 398815f1ad5SJed Brown 399815f1ad5SJed Brown PetscFunctionBegin; 400815f1ad5SJed Brown ssp->nstages = nstages; 401815f1ad5SJed Brown PetscFunctionReturn(0); 402815f1ad5SJed Brown } 403815f1ad5SJed Brown #undef __FUNCT__ 404815f1ad5SJed Brown #define __FUNCT__ "TSSSPGetNumStages_SSP" 405815f1ad5SJed Brown PetscErrorCode TSSSPGetNumStages_SSP(TS ts,PetscInt *nstages) 406815f1ad5SJed Brown { 407815f1ad5SJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 408815f1ad5SJed Brown 409815f1ad5SJed Brown PetscFunctionBegin; 410815f1ad5SJed Brown *nstages = ssp->nstages; 411815f1ad5SJed Brown PetscFunctionReturn(0); 412815f1ad5SJed Brown } 413ef7bb5aaSJed Brown 414ef7bb5aaSJed Brown #undef __FUNCT__ 415ef7bb5aaSJed Brown #define __FUNCT__ "TSSetFromOptions_SSP" 416ef7bb5aaSJed Brown static PetscErrorCode TSSetFromOptions_SSP(TS ts) 417ef7bb5aaSJed Brown { 418ef7bb5aaSJed Brown char tname[256] = TSSSPRKS2; 419ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 420ef7bb5aaSJed Brown PetscErrorCode ierr; 421ace3abfcSBarry Smith PetscBool flg; 422ef7bb5aaSJed Brown 423ef7bb5aaSJed Brown PetscFunctionBegin; 424ef7bb5aaSJed Brown ierr = PetscOptionsHead("SSP ODE solver options");CHKERRQ(ierr); 425ef7bb5aaSJed Brown { 426ef7bb5aaSJed Brown ierr = PetscOptionsList("-ts_ssp_type","Type of SSP method","TSSSPSetType",TSSSPList,tname,tname,sizeof(tname),&flg);CHKERRQ(ierr); 427ef7bb5aaSJed Brown if (flg) { 428ef7bb5aaSJed Brown ierr = TSSSPSetType(ts,tname);CHKERRQ(ierr); 429ef7bb5aaSJed Brown } 4300298fd71SBarry Smith ierr = PetscOptionsInt("-ts_ssp_nstages","Number of stages","TSSSPSetNumStages",ssp->nstages,&ssp->nstages,NULL);CHKERRQ(ierr); 431ef7bb5aaSJed Brown } 432ef7bb5aaSJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 433ef7bb5aaSJed Brown PetscFunctionReturn(0); 434ef7bb5aaSJed Brown } 435ef7bb5aaSJed Brown 436ef7bb5aaSJed Brown #undef __FUNCT__ 437ef7bb5aaSJed Brown #define __FUNCT__ "TSView_SSP" 438ef7bb5aaSJed Brown static PetscErrorCode TSView_SSP(TS ts,PetscViewer viewer) 439ef7bb5aaSJed Brown { 440ef7bb5aaSJed Brown PetscFunctionBegin; 441ef7bb5aaSJed Brown PetscFunctionReturn(0); 442ef7bb5aaSJed Brown } 443ef7bb5aaSJed Brown 444ef7bb5aaSJed Brown /* ------------------------------------------------------------ */ 445ef7bb5aaSJed Brown 446ef7bb5aaSJed Brown /*MC 4478ab3e0fcSJed Brown TSSSP - Explicit strong stability preserving ODE solver 4488ab3e0fcSJed Brown 4498ab3e0fcSJed Brown Most hyperbolic conservation laws have exact solutions that are total variation diminishing (TVD) or total variation 4508ab3e0fcSJed Brown bounded (TVB) although these solutions often contain discontinuities. Spatial discretizations such as Godunov's 4518ab3e0fcSJed Brown scheme and high-resolution finite volume methods (TVD limiters, ENO/WENO) are designed to preserve these properties, 4528ab3e0fcSJed Brown but they are usually formulated using a forward Euler time discretization or by coupling the space and time 4538ab3e0fcSJed Brown discretization as in the classical Lax-Wendroff scheme. When the space and time discretization is coupled, it is very 4548ab3e0fcSJed Brown difficult to produce schemes with high temporal accuracy while preserving TVD properties. An alternative is the 4558ab3e0fcSJed Brown semidiscrete formulation where we choose a spatial discretization that is TVD with forward Euler and then choose a 4568ab3e0fcSJed Brown time discretization that preserves the TVD property. Such integrators are called strong stability preserving (SSP). 4578ab3e0fcSJed Brown 4588ab3e0fcSJed Brown Let c_eff be the minimum number of function evaluations required to step as far as one step of forward Euler while 4598ab3e0fcSJed Brown still being SSP. Some theoretical bounds 4608ab3e0fcSJed Brown 4618ab3e0fcSJed Brown 1. There are no explicit methods with c_eff > 1. 4620085e20eSJed Brown 4638ab3e0fcSJed Brown 2. There are no explicit methods beyond order 4 (for nonlinear problems) and c_eff > 0. 4640085e20eSJed Brown 4658ab3e0fcSJed Brown 3. There are no implicit methods with order greater than 1 and c_eff > 2. 4668ab3e0fcSJed Brown 4678ab3e0fcSJed Brown This integrator provides Runge-Kutta methods of order 2, 3, and 4 with maximal values of c_eff. More stages allows 4688ab3e0fcSJed Brown for larger values of c_eff which improves efficiency. These implementations are low-memory and only use 2 or 3 work 4698ab3e0fcSJed Brown vectors regardless of the total number of stages, so e.g. 25-stage 3rd order methods may be an excellent choice. 4708ab3e0fcSJed Brown 4718ab3e0fcSJed Brown Methods can be chosen with -ts_ssp_type {rks2,rks3,rk104} 4728ab3e0fcSJed Brown 4738ab3e0fcSJed Brown rks2: Second order methods with any number s>1 of stages. c_eff = (s-1)/s 4748ab3e0fcSJed Brown 4758ab3e0fcSJed Brown rks3: Third order methods with s=n^2 stages, n>1. c_eff = (s-n)/s 4768ab3e0fcSJed Brown 4778ab3e0fcSJed Brown rk104: A 10-stage fourth order method. c_eff = 0.6 478ef7bb5aaSJed Brown 479ef7bb5aaSJed Brown Level: beginner 480ef7bb5aaSJed Brown 4817b6bb2c6SJed Brown References: 4827b6bb2c6SJed Brown Ketcheson, Highly efficient strong stability preserving Runge-Kutta methods with low-storage implementations, SISC, 2008. 4837b6bb2c6SJed Brown 4847b6bb2c6SJed Brown Gottlieb, Ketcheson, and Shu, High order strong stability preserving time discretizations, J Scientific Computing, 2009. 4857b6bb2c6SJed Brown 486ef7bb5aaSJed Brown .seealso: TSCreate(), TS, TSSetType() 487ef7bb5aaSJed Brown 488ef7bb5aaSJed Brown M*/ 489ef7bb5aaSJed Brown #undef __FUNCT__ 490ef7bb5aaSJed Brown #define __FUNCT__ "TSCreate_SSP" 4918cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_SSP(TS ts) 492ef7bb5aaSJed Brown { 493ef7bb5aaSJed Brown TS_SSP *ssp; 494ef7bb5aaSJed Brown PetscErrorCode ierr; 495ef7bb5aaSJed Brown 496ef7bb5aaSJed Brown PetscFunctionBegin; 497ef7bb5aaSJed Brown if (!TSSSPList) { 498bdf89e91SBarry Smith ierr = PetscFunctionListAdd(&TSSSPList,TSSSPRKS2, (void(*)(void))TSSSPStep_RK_2);CHKERRQ(ierr); 499bdf89e91SBarry Smith ierr = PetscFunctionListAdd(&TSSSPList,TSSSPRKS3, (void(*)(void))TSSSPStep_RK_3);CHKERRQ(ierr); 500bdf89e91SBarry Smith ierr = PetscFunctionListAdd(&TSSSPList,TSSSPRK104, (void(*)(void))TSSSPStep_RK_10_4);CHKERRQ(ierr); 501ef7bb5aaSJed Brown } 502ef7bb5aaSJed Brown 503ef7bb5aaSJed Brown ts->ops->setup = TSSetUp_SSP; 504ef7bb5aaSJed Brown ts->ops->step = TSStep_SSP; 505277b19d0SLisandro Dalcin ts->ops->reset = TSReset_SSP; 506ef7bb5aaSJed Brown ts->ops->destroy = TSDestroy_SSP; 507ef7bb5aaSJed Brown ts->ops->setfromoptions = TSSetFromOptions_SSP; 508ef7bb5aaSJed Brown ts->ops->view = TSView_SSP; 509ef7bb5aaSJed Brown 510ef7bb5aaSJed Brown ierr = PetscNewLog(ts,TS_SSP,&ssp);CHKERRQ(ierr); 511ef7bb5aaSJed Brown ts->data = (void*)ssp; 512ef7bb5aaSJed Brown 513bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetType_C",TSSSPGetType_SSP);CHKERRQ(ierr); 514bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetType_C",TSSSPSetType_SSP);CHKERRQ(ierr); 515bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetNumStages_C",TSSSPGetNumStages_SSP);CHKERRQ(ierr); 516bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetNumStages_C",TSSSPSetNumStages_SSP);CHKERRQ(ierr); 517815f1ad5SJed Brown 518ef7bb5aaSJed Brown ierr = TSSSPSetType(ts,TSSSPRKS2);CHKERRQ(ierr); 519ef7bb5aaSJed Brown ssp->nstages = 5; 520ef7bb5aaSJed Brown PetscFunctionReturn(0); 521ef7bb5aaSJed Brown } 522