1ef7bb5aaSJed Brown /* 2ef7bb5aaSJed Brown Code for Timestepping with explicit SSP. 3ef7bb5aaSJed Brown */ 4c6db04a5SJed Brown #include <private/tsimpl.h> /*I "petscts.h" I*/ 5ef7bb5aaSJed Brown 6ef7bb5aaSJed Brown PetscFList TSSSPList = 0; 7ef7bb5aaSJed Brown #define TSSSPType char* 8ef7bb5aaSJed Brown 9ef7bb5aaSJed Brown #define TSSSPRKS2 "rks2" 10ef7bb5aaSJed Brown #define TSSSPRKS3 "rks3" 11ef7bb5aaSJed Brown #define TSSSPRK104 "rk104" 12ef7bb5aaSJed Brown 13ef7bb5aaSJed Brown typedef struct { 14ef7bb5aaSJed Brown PetscErrorCode (*onestep)(TS,PetscReal,PetscReal,Vec); 15ef7bb5aaSJed Brown PetscInt nstages; 16ef7bb5aaSJed Brown Vec *work; 17ef7bb5aaSJed Brown PetscInt nwork; 18ace3abfcSBarry Smith PetscBool workout; 19ef7bb5aaSJed Brown } TS_SSP; 20ef7bb5aaSJed Brown 21ef7bb5aaSJed Brown 22ef7bb5aaSJed Brown #undef __FUNCT__ 23ef7bb5aaSJed Brown #define __FUNCT__ "SSPGetWorkVectors" 24ef7bb5aaSJed Brown static PetscErrorCode SSPGetWorkVectors(TS ts,PetscInt n,Vec **work) 25ef7bb5aaSJed Brown { 26ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 27ef7bb5aaSJed Brown PetscErrorCode ierr; 28ef7bb5aaSJed Brown 29ef7bb5aaSJed Brown PetscFunctionBegin; 30e32f2f54SBarry Smith if (ssp->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Work vectors already gotten"); 31ef7bb5aaSJed Brown if (ssp->nwork < n) { 32ef7bb5aaSJed Brown if (ssp->nwork > 0) { 33574deadeSBarry Smith ierr = VecDestroyVecs(ssp->nwork,&ssp->work);CHKERRQ(ierr); 34ef7bb5aaSJed Brown } 35ef7bb5aaSJed Brown ierr = VecDuplicateVecs(ts->vec_sol,n,&ssp->work);CHKERRQ(ierr); 36ef7bb5aaSJed Brown ssp->nwork = n; 37ef7bb5aaSJed Brown } 38ef7bb5aaSJed Brown *work = ssp->work; 39ef7bb5aaSJed Brown ssp->workout = PETSC_TRUE; 40ef7bb5aaSJed Brown PetscFunctionReturn(0); 41ef7bb5aaSJed Brown } 42ef7bb5aaSJed Brown 43ef7bb5aaSJed Brown #undef __FUNCT__ 44ef7bb5aaSJed Brown #define __FUNCT__ "SSPRestoreWorkVectors" 45ef7bb5aaSJed Brown static PetscErrorCode SSPRestoreWorkVectors(TS ts,PetscInt n,Vec **work) 46ef7bb5aaSJed Brown { 47ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 48ef7bb5aaSJed Brown 49ef7bb5aaSJed Brown PetscFunctionBegin; 50e32f2f54SBarry Smith if (!ssp->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Work vectors have not been gotten"); 51e32f2f54SBarry Smith if (*work != ssp->work) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong work vectors checked out"); 52ef7bb5aaSJed Brown ssp->workout = PETSC_FALSE; 53ef7bb5aaSJed Brown *work = PETSC_NULL; 54ef7bb5aaSJed Brown PetscFunctionReturn(0); 55ef7bb5aaSJed Brown } 56ef7bb5aaSJed Brown 57ef7bb5aaSJed Brown 58ef7bb5aaSJed Brown #undef __FUNCT__ 59ef7bb5aaSJed Brown #define __FUNCT__ "SSPStep_RK_2" 60ef7bb5aaSJed Brown /* Optimal second order SSP Runge-Kutta, low-storage, c_eff=(s-1)/s */ 61ef7bb5aaSJed Brown /* Pseudocode 2 of Ketcheson 2008 */ 62ef7bb5aaSJed Brown static PetscErrorCode SSPStep_RK_2(TS ts,PetscReal t0,PetscReal dt,Vec sol) 63ef7bb5aaSJed Brown { 64ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 65ef7bb5aaSJed Brown Vec *work,F; 66ef7bb5aaSJed Brown PetscInt i,s; 67ef7bb5aaSJed Brown PetscErrorCode ierr; 68ef7bb5aaSJed Brown 69ef7bb5aaSJed Brown PetscFunctionBegin; 70ef7bb5aaSJed Brown s = ssp->nstages; 71ef7bb5aaSJed Brown ierr = SSPGetWorkVectors(ts,2,&work);CHKERRQ(ierr); 72ef7bb5aaSJed Brown F = work[1]; 73ef7bb5aaSJed Brown ierr = VecCopy(sol,work[0]);CHKERRQ(ierr); 74ef7bb5aaSJed Brown for (i=0; i<s-1; i++) { 75ef7bb5aaSJed Brown ierr = TSComputeRHSFunction(ts,t0+dt*(i/(s-1.)),work[0],F);CHKERRQ(ierr); 76ef7bb5aaSJed Brown ierr = VecAXPY(work[0],dt/(s-1.),F);CHKERRQ(ierr); 77ef7bb5aaSJed Brown } 78ef7bb5aaSJed Brown ierr = TSComputeRHSFunction(ts,t0+dt,work[0],F);CHKERRQ(ierr); 79ef7bb5aaSJed Brown ierr = VecAXPBYPCZ(sol,(s-1.)/s,dt/s,1./s,work[0],F);CHKERRQ(ierr); 80ef7bb5aaSJed Brown ierr = SSPRestoreWorkVectors(ts,2,&work);CHKERRQ(ierr); 81ef7bb5aaSJed Brown PetscFunctionReturn(0); 82ef7bb5aaSJed Brown } 83ef7bb5aaSJed Brown 84ef7bb5aaSJed Brown #undef __FUNCT__ 85ef7bb5aaSJed Brown #define __FUNCT__ "SSPStep_RK_3" 86ef7bb5aaSJed Brown /* Optimal third order SSP Runge-Kutta, low-storage, c_eff=(sqrt(s)-1)/sqrt(s), where sqrt(s) is an integer */ 87ef7bb5aaSJed Brown /* Pseudocode 2 of Ketcheson 2008 */ 88ef7bb5aaSJed Brown static PetscErrorCode SSPStep_RK_3(TS ts,PetscReal t0,PetscReal dt,Vec sol) 89ef7bb5aaSJed Brown { 90ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 91ef7bb5aaSJed Brown Vec *work,F; 92ef7bb5aaSJed Brown PetscInt i,s,n,r; 93ef7bb5aaSJed Brown PetscReal c; 94ef7bb5aaSJed Brown PetscErrorCode ierr; 95ef7bb5aaSJed Brown 96ef7bb5aaSJed Brown PetscFunctionBegin; 97ef7bb5aaSJed Brown s = ssp->nstages; 98fad8df86SSatish Balay n = (PetscInt)(sqrt((PetscReal)s)+0.001); 99ef7bb5aaSJed Brown r = s-n; 100e32f2f54SBarry 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); 101ef7bb5aaSJed Brown ierr = SSPGetWorkVectors(ts,3,&work);CHKERRQ(ierr); 102ef7bb5aaSJed Brown F = work[2]; 103ef7bb5aaSJed Brown ierr = VecCopy(sol,work[0]);CHKERRQ(ierr); 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); 106ef7bb5aaSJed Brown ierr = TSComputeRHSFunction(ts,t0+c*dt,work[0],F);CHKERRQ(ierr); 107ef7bb5aaSJed Brown ierr = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr); 108ef7bb5aaSJed Brown } 109ef7bb5aaSJed Brown ierr = VecCopy(work[0],work[1]);CHKERRQ(ierr); 110ef7bb5aaSJed Brown for ( ; i<n*(n+1)/2-1; i++) { 111ef7bb5aaSJed Brown c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n); 112ef7bb5aaSJed Brown ierr = TSComputeRHSFunction(ts,t0+c*dt,work[0],F);CHKERRQ(ierr); 113ef7bb5aaSJed Brown ierr = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr); 114ef7bb5aaSJed Brown } 115ef7bb5aaSJed Brown { 116ef7bb5aaSJed Brown c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n); 117ef7bb5aaSJed Brown ierr = TSComputeRHSFunction(ts,t0+c*dt,work[0],F);CHKERRQ(ierr); 118ef7bb5aaSJed 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); 119ef7bb5aaSJed Brown i++; 120ef7bb5aaSJed Brown } 121ef7bb5aaSJed Brown for ( ; i<s; i++) { 122ef7bb5aaSJed Brown c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n); 123ef7bb5aaSJed Brown ierr = TSComputeRHSFunction(ts,t0+c*dt,work[0],F);CHKERRQ(ierr); 124ef7bb5aaSJed Brown ierr = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr); 125ef7bb5aaSJed Brown } 126ef7bb5aaSJed Brown ierr = VecCopy(work[0],sol);CHKERRQ(ierr); 127ef7bb5aaSJed Brown ierr = SSPRestoreWorkVectors(ts,3,&work);CHKERRQ(ierr); 128ef7bb5aaSJed Brown PetscFunctionReturn(0); 129ef7bb5aaSJed Brown } 130ef7bb5aaSJed Brown 131ef7bb5aaSJed Brown #undef __FUNCT__ 132ef7bb5aaSJed Brown #define __FUNCT__ "SSPStep_RK_10_4" 133ef7bb5aaSJed Brown /* Optimal fourth order SSP Runge-Kutta, low-storage (2N), c_eff=0.6 */ 134ef7bb5aaSJed Brown /* SSPRK(10,4), Pseudocode 3 of Ketcheson 2008 */ 135ef7bb5aaSJed Brown static PetscErrorCode SSPStep_RK_10_4(TS ts,PetscReal t0,PetscReal dt,Vec sol) 136ef7bb5aaSJed Brown { 137ef7bb5aaSJed Brown const PetscReal c[10] = {0, 1./6, 2./6, 3./6, 4./6, 2./6, 3./6, 4./6, 5./6, 1}; 138ef7bb5aaSJed Brown Vec *work,F; 1395a586d82SBarry Smith PetscInt i; 140ef7bb5aaSJed Brown PetscErrorCode ierr; 141ef7bb5aaSJed Brown 142ef7bb5aaSJed Brown PetscFunctionBegin; 143ef7bb5aaSJed Brown ierr = SSPGetWorkVectors(ts,3,&work);CHKERRQ(ierr); 144ef7bb5aaSJed Brown F = work[2]; 145ef7bb5aaSJed Brown ierr = VecCopy(sol,work[0]);CHKERRQ(ierr); 146ef7bb5aaSJed Brown for (i=0; i<5; i++) { 147ef7bb5aaSJed Brown ierr = TSComputeRHSFunction(ts,t0+c[i],work[0],F);CHKERRQ(ierr); 148ef7bb5aaSJed Brown ierr = VecAXPY(work[0],dt/6,F);CHKERRQ(ierr); 149ef7bb5aaSJed Brown } 150ef7bb5aaSJed Brown ierr = VecAXPBYPCZ(work[1],1./25,9./25,0,sol,work[0]);CHKERRQ(ierr); 151ef7bb5aaSJed Brown ierr = VecAXPBY(work[0],15,-5,work[1]);CHKERRQ(ierr); 152ef7bb5aaSJed Brown for ( ; i<9; i++) { 153ef7bb5aaSJed Brown ierr = TSComputeRHSFunction(ts,t0+c[i],work[0],F);CHKERRQ(ierr); 154ef7bb5aaSJed Brown ierr = VecAXPY(work[0],dt/6,F);CHKERRQ(ierr); 155ef7bb5aaSJed Brown } 156ef7bb5aaSJed Brown ierr = TSComputeRHSFunction(ts,t0+dt,work[0],F);CHKERRQ(ierr); 157ef7bb5aaSJed Brown ierr = VecAXPBYPCZ(work[1],3./5,dt/10,1,work[0],F);CHKERRQ(ierr); 158ef7bb5aaSJed Brown ierr = VecCopy(work[1],sol);CHKERRQ(ierr); 159ef7bb5aaSJed Brown ierr = SSPRestoreWorkVectors(ts,3,&work);CHKERRQ(ierr); 160ef7bb5aaSJed Brown PetscFunctionReturn(0); 161ef7bb5aaSJed Brown } 162ef7bb5aaSJed Brown 163ef7bb5aaSJed Brown 164ef7bb5aaSJed Brown #undef __FUNCT__ 165ef7bb5aaSJed Brown #define __FUNCT__ "TSSetUp_SSP" 166ef7bb5aaSJed Brown static PetscErrorCode TSSetUp_SSP(TS ts) 167ef7bb5aaSJed Brown { 168ef7bb5aaSJed Brown /* TS_SSP *ssp = (TS_SSP*)ts->data; */ 169ef7bb5aaSJed Brown /* PetscErrorCode ierr; */ 170ef7bb5aaSJed Brown 171ef7bb5aaSJed Brown PetscFunctionBegin; 172ef7bb5aaSJed Brown PetscFunctionReturn(0); 173ef7bb5aaSJed Brown } 174ef7bb5aaSJed Brown 175ef7bb5aaSJed Brown #undef __FUNCT__ 176ef7bb5aaSJed Brown #define __FUNCT__ "TSStep_SSP" 177*193ac0bcSJed Brown static PetscErrorCode TSStep_SSP(TS ts) 178ef7bb5aaSJed Brown { 179ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 180ef7bb5aaSJed Brown Vec sol = ts->vec_sol; 181ef7bb5aaSJed Brown PetscErrorCode ierr; 182ef7bb5aaSJed Brown 183ef7bb5aaSJed Brown PetscFunctionBegin; 184186e87acSLisandro Dalcin ierr = (*ssp->onestep)(ts,ts->ptime,ts->time_step,sol);CHKERRQ(ierr); 185186e87acSLisandro Dalcin ts->ptime += ts->time_step; 186ef7bb5aaSJed Brown ts->steps++; 187ef7bb5aaSJed Brown PetscFunctionReturn(0); 188ef7bb5aaSJed Brown } 189ef7bb5aaSJed Brown /*------------------------------------------------------------*/ 190ef7bb5aaSJed Brown #undef __FUNCT__ 191277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_SSP" 192277b19d0SLisandro Dalcin static PetscErrorCode TSReset_SSP(TS ts) 193277b19d0SLisandro Dalcin { 194277b19d0SLisandro Dalcin TS_SSP *ssp = (TS_SSP*)ts->data; 195277b19d0SLisandro Dalcin PetscErrorCode ierr; 196277b19d0SLisandro Dalcin 197277b19d0SLisandro Dalcin PetscFunctionBegin; 198277b19d0SLisandro Dalcin if (ssp->work) {ierr = VecDestroyVecs(ssp->nwork,&ssp->work);CHKERRQ(ierr);} 199277b19d0SLisandro Dalcin ssp->nwork = 0; 200277b19d0SLisandro Dalcin ssp->workout = PETSC_FALSE; 201277b19d0SLisandro Dalcin PetscFunctionReturn(0); 202277b19d0SLisandro Dalcin } 203277b19d0SLisandro Dalcin 204277b19d0SLisandro Dalcin #undef __FUNCT__ 205ef7bb5aaSJed Brown #define __FUNCT__ "TSDestroy_SSP" 206ef7bb5aaSJed Brown static PetscErrorCode TSDestroy_SSP(TS ts) 207ef7bb5aaSJed Brown { 208ef7bb5aaSJed Brown PetscErrorCode ierr; 209ef7bb5aaSJed Brown 210ef7bb5aaSJed Brown PetscFunctionBegin; 211277b19d0SLisandro Dalcin ierr = TSReset_SSP(ts);CHKERRQ(ierr); 212c31cb41cSBarry Smith ierr = PetscFree(ts->data);CHKERRQ(ierr); 213ef7bb5aaSJed Brown PetscFunctionReturn(0); 214ef7bb5aaSJed Brown } 215ef7bb5aaSJed Brown /*------------------------------------------------------------*/ 216ef7bb5aaSJed Brown 217ef7bb5aaSJed Brown #undef __FUNCT__ 218ef7bb5aaSJed Brown #define __FUNCT__ "TSSSPSetType" 219ef7bb5aaSJed Brown static PetscErrorCode TSSSPSetType(TS ts,const TSSSPType type) 220ef7bb5aaSJed Brown { 221ef7bb5aaSJed Brown PetscErrorCode ierr,(*r)(TS,PetscReal,PetscReal,Vec); 222ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 223ef7bb5aaSJed Brown 224ef7bb5aaSJed Brown PetscFunctionBegin; 2254b91b6eaSBarry Smith ierr = PetscFListFind(TSSSPList,((PetscObject)ts)->comm,type,PETSC_TRUE,(void(**)(void))&r);CHKERRQ(ierr); 226e32f2f54SBarry Smith if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TS_SSP type %s given",type); 227ef7bb5aaSJed Brown ssp->onestep = r; 228ef7bb5aaSJed Brown PetscFunctionReturn(0); 229ef7bb5aaSJed Brown } 230ef7bb5aaSJed Brown 231ef7bb5aaSJed Brown #undef __FUNCT__ 232ef7bb5aaSJed Brown #define __FUNCT__ "TSSetFromOptions_SSP" 233ef7bb5aaSJed Brown static PetscErrorCode TSSetFromOptions_SSP(TS ts) 234ef7bb5aaSJed Brown { 235ef7bb5aaSJed Brown char tname[256] = TSSSPRKS2; 236ef7bb5aaSJed Brown TS_SSP *ssp = (TS_SSP*)ts->data; 237ef7bb5aaSJed Brown PetscErrorCode ierr; 238ace3abfcSBarry Smith PetscBool flg; 239ef7bb5aaSJed Brown 240ef7bb5aaSJed Brown PetscFunctionBegin; 241ef7bb5aaSJed Brown ierr = PetscOptionsHead("SSP ODE solver options");CHKERRQ(ierr); 242ef7bb5aaSJed Brown { 243ef7bb5aaSJed Brown ierr = PetscOptionsList("-ts_ssp_type","Type of SSP method","TSSSPSetType",TSSSPList,tname,tname,sizeof(tname),&flg);CHKERRQ(ierr); 244ef7bb5aaSJed Brown if (flg) { 245ef7bb5aaSJed Brown ierr = TSSSPSetType(ts,tname);CHKERRQ(ierr); 246ef7bb5aaSJed Brown } 247ef7bb5aaSJed Brown ierr = PetscOptionsInt("-ts_ssp_nstages","Number of stages","TSSSPSetNumStages",ssp->nstages,&ssp->nstages,PETSC_NULL);CHKERRQ(ierr); 248ef7bb5aaSJed Brown } 249ef7bb5aaSJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 250ef7bb5aaSJed Brown PetscFunctionReturn(0); 251ef7bb5aaSJed Brown } 252ef7bb5aaSJed Brown 253ef7bb5aaSJed Brown #undef __FUNCT__ 254ef7bb5aaSJed Brown #define __FUNCT__ "TSView_SSP" 255ef7bb5aaSJed Brown static PetscErrorCode TSView_SSP(TS ts,PetscViewer viewer) 256ef7bb5aaSJed Brown { 257ef7bb5aaSJed Brown PetscFunctionBegin; 258ef7bb5aaSJed Brown PetscFunctionReturn(0); 259ef7bb5aaSJed Brown } 260ef7bb5aaSJed Brown 261ef7bb5aaSJed Brown /* ------------------------------------------------------------ */ 262ef7bb5aaSJed Brown 263ef7bb5aaSJed Brown /*MC 2648ab3e0fcSJed Brown TSSSP - Explicit strong stability preserving ODE solver 2658ab3e0fcSJed Brown 2668ab3e0fcSJed Brown Most hyperbolic conservation laws have exact solutions that are total variation diminishing (TVD) or total variation 2678ab3e0fcSJed Brown bounded (TVB) although these solutions often contain discontinuities. Spatial discretizations such as Godunov's 2688ab3e0fcSJed Brown scheme and high-resolution finite volume methods (TVD limiters, ENO/WENO) are designed to preserve these properties, 2698ab3e0fcSJed Brown but they are usually formulated using a forward Euler time discretization or by coupling the space and time 2708ab3e0fcSJed Brown discretization as in the classical Lax-Wendroff scheme. When the space and time discretization is coupled, it is very 2718ab3e0fcSJed Brown difficult to produce schemes with high temporal accuracy while preserving TVD properties. An alternative is the 2728ab3e0fcSJed Brown semidiscrete formulation where we choose a spatial discretization that is TVD with forward Euler and then choose a 2738ab3e0fcSJed Brown time discretization that preserves the TVD property. Such integrators are called strong stability preserving (SSP). 2748ab3e0fcSJed Brown 2758ab3e0fcSJed Brown Let c_eff be the minimum number of function evaluations required to step as far as one step of forward Euler while 2768ab3e0fcSJed Brown still being SSP. Some theoretical bounds 2778ab3e0fcSJed Brown 2788ab3e0fcSJed Brown 1. There are no explicit methods with c_eff > 1. 2790085e20eSJed Brown 2808ab3e0fcSJed Brown 2. There are no explicit methods beyond order 4 (for nonlinear problems) and c_eff > 0. 2810085e20eSJed Brown 2828ab3e0fcSJed Brown 3. There are no implicit methods with order greater than 1 and c_eff > 2. 2838ab3e0fcSJed Brown 2848ab3e0fcSJed Brown This integrator provides Runge-Kutta methods of order 2, 3, and 4 with maximal values of c_eff. More stages allows 2858ab3e0fcSJed Brown for larger values of c_eff which improves efficiency. These implementations are low-memory and only use 2 or 3 work 2868ab3e0fcSJed Brown vectors regardless of the total number of stages, so e.g. 25-stage 3rd order methods may be an excellent choice. 2878ab3e0fcSJed Brown 2888ab3e0fcSJed Brown Methods can be chosen with -ts_ssp_type {rks2,rks3,rk104} 2898ab3e0fcSJed Brown 2908ab3e0fcSJed Brown rks2: Second order methods with any number s>1 of stages. c_eff = (s-1)/s 2918ab3e0fcSJed Brown 2928ab3e0fcSJed Brown rks3: Third order methods with s=n^2 stages, n>1. c_eff = (s-n)/s 2938ab3e0fcSJed Brown 2948ab3e0fcSJed Brown rk104: A 10-stage fourth order method. c_eff = 0.6 295ef7bb5aaSJed Brown 296ef7bb5aaSJed Brown Level: beginner 297ef7bb5aaSJed Brown 2987b6bb2c6SJed Brown References: 2997b6bb2c6SJed Brown Ketcheson, Highly efficient strong stability preserving Runge-Kutta methods with low-storage implementations, SISC, 2008. 3007b6bb2c6SJed Brown 3017b6bb2c6SJed Brown Gottlieb, Ketcheson, and Shu, High order strong stability preserving time discretizations, J Scientific Computing, 2009. 3027b6bb2c6SJed Brown 303ef7bb5aaSJed Brown .seealso: TSCreate(), TS, TSSetType() 304ef7bb5aaSJed Brown 305ef7bb5aaSJed Brown M*/ 306ef7bb5aaSJed Brown EXTERN_C_BEGIN 307ef7bb5aaSJed Brown #undef __FUNCT__ 308ef7bb5aaSJed Brown #define __FUNCT__ "TSCreate_SSP" 3097087cfbeSBarry Smith PetscErrorCode TSCreate_SSP(TS ts) 310ef7bb5aaSJed Brown { 311ef7bb5aaSJed Brown TS_SSP *ssp; 312ef7bb5aaSJed Brown PetscErrorCode ierr; 313ef7bb5aaSJed Brown 314ef7bb5aaSJed Brown PetscFunctionBegin; 315ef7bb5aaSJed Brown if (!TSSSPList) { 316ef7bb5aaSJed Brown ierr = PetscFListAdd(&TSSSPList,TSSSPRKS2, "SSPStep_RK_2", (void(*)(void))SSPStep_RK_2);CHKERRQ(ierr); 317ef7bb5aaSJed Brown ierr = PetscFListAdd(&TSSSPList,TSSSPRKS3, "SSPStep_RK_3", (void(*)(void))SSPStep_RK_3);CHKERRQ(ierr); 318ef7bb5aaSJed Brown ierr = PetscFListAdd(&TSSSPList,TSSSPRK104, "SSPStep_RK_10_4",(void(*)(void))SSPStep_RK_10_4);CHKERRQ(ierr); 319ef7bb5aaSJed Brown } 320ef7bb5aaSJed Brown 321ef7bb5aaSJed Brown ts->ops->setup = TSSetUp_SSP; 322ef7bb5aaSJed Brown ts->ops->step = TSStep_SSP; 323277b19d0SLisandro Dalcin ts->ops->reset = TSReset_SSP; 324ef7bb5aaSJed Brown ts->ops->destroy = TSDestroy_SSP; 325ef7bb5aaSJed Brown ts->ops->setfromoptions = TSSetFromOptions_SSP; 326ef7bb5aaSJed Brown ts->ops->view = TSView_SSP; 327ef7bb5aaSJed Brown 328ef7bb5aaSJed Brown ierr = PetscNewLog(ts,TS_SSP,&ssp);CHKERRQ(ierr); 329ef7bb5aaSJed Brown ts->data = (void*)ssp; 330ef7bb5aaSJed Brown 331ef7bb5aaSJed Brown ierr = TSSSPSetType(ts,TSSSPRKS2);CHKERRQ(ierr); 332ef7bb5aaSJed Brown ssp->nstages = 5; 333ef7bb5aaSJed Brown PetscFunctionReturn(0); 334ef7bb5aaSJed Brown } 335ef7bb5aaSJed Brown EXTERN_C_END 336