Actual source code: ssp.c
1: #define PETSCTS_DLL
3: /*
4: Code for Timestepping with explicit SSP.
5: */
6: #include private/tsimpl.h
8: PetscFList TSSSPList = 0;
9: #define TSSSPType char*
11: #define TSSSPRKS2 "rks2"
12: #define TSSSPRKS3 "rks3"
13: #define TSSSPRK104 "rk104"
15: typedef struct {
16: PetscErrorCode (*onestep)(TS,PetscReal,PetscReal,Vec);
17: PetscInt nstages;
18: Vec xdot;
19: Vec *work;
20: PetscInt nwork;
21: PetscTruth workout;
22: } TS_SSP;
27: static PetscErrorCode SSPGetWorkVectors(TS ts,PetscInt n,Vec **work)
28: {
29: TS_SSP *ssp = (TS_SSP*)ts->data;
33: if (ssp->workout) SETERRQ(PETSC_ERR_PLIB,"Work vectors already gotten");
34: if (ssp->nwork < n) {
35: if (ssp->nwork > 0) {
36: VecDestroyVecs(ssp->work,ssp->nwork);
37: }
38: VecDuplicateVecs(ts->vec_sol,n,&ssp->work);
39: ssp->nwork = n;
40: }
41: *work = ssp->work;
42: ssp->workout = PETSC_TRUE;
43: return(0);
44: }
48: static PetscErrorCode SSPRestoreWorkVectors(TS ts,PetscInt n,Vec **work)
49: {
50: TS_SSP *ssp = (TS_SSP*)ts->data;
53: if (!ssp->workout) SETERRQ(PETSC_ERR_ORDER,"Work vectors have not been gotten");
54: if (*work != ssp->work) SETERRQ(PETSC_ERR_PLIB,"Wrong work vectors checked out");
55: ssp->workout = PETSC_FALSE;
56: *work = PETSC_NULL;
57: return(0);
58: }
63: /* Optimal second order SSP Runge-Kutta, low-storage, c_eff=(s-1)/s */
64: /* Pseudocode 2 of Ketcheson 2008 */
65: static PetscErrorCode SSPStep_RK_2(TS ts,PetscReal t0,PetscReal dt,Vec sol)
66: {
67: TS_SSP *ssp = (TS_SSP*)ts->data;
68: Vec *work,F;
69: PetscInt i,s;
73: s = ssp->nstages;
74: SSPGetWorkVectors(ts,2,&work);
75: F = work[1];
76: VecCopy(sol,work[0]);
77: for (i=0; i<s-1; i++) {
78: TSComputeRHSFunction(ts,t0+dt*(i/(s-1.)),work[0],F);
79: VecAXPY(work[0],dt/(s-1.),F);
80: }
81: TSComputeRHSFunction(ts,t0+dt,work[0],F);
82: VecAXPBYPCZ(sol,(s-1.)/s,dt/s,1./s,work[0],F);
83: SSPRestoreWorkVectors(ts,2,&work);
84: return(0);
85: }
89: /* Optimal third order SSP Runge-Kutta, low-storage, c_eff=(sqrt(s)-1)/sqrt(s), where sqrt(s) is an integer */
90: /* Pseudocode 2 of Ketcheson 2008 */
91: static PetscErrorCode SSPStep_RK_3(TS ts,PetscReal t0,PetscReal dt,Vec sol)
92: {
93: TS_SSP *ssp = (TS_SSP*)ts->data;
94: Vec *work,F;
95: PetscInt i,s,n,r;
96: PetscReal c;
100: s = ssp->nstages;
101: n = (PetscInt)(sqrt((PetscReal)s)+0.001);
102: r = s-n;
103: if (n*n != s) SETERRQ1(PETSC_ERR_SUP,"No support for optimal third order schemes with %d stages, must be a square number at least 4",s);
104: SSPGetWorkVectors(ts,3,&work);
105: F = work[2];
106: VecCopy(sol,work[0]);
107: for (i=0; i<(n-1)*(n-2)/2; i++) {
108: c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
109: TSComputeRHSFunction(ts,t0+c*dt,work[0],F);
110: VecAXPY(work[0],dt/r,F);
111: }
112: VecCopy(work[0],work[1]);
113: for ( ; i<n*(n+1)/2-1; i++) {
114: c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
115: TSComputeRHSFunction(ts,t0+c*dt,work[0],F);
116: VecAXPY(work[0],dt/r,F);
117: }
118: {
119: c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
120: TSComputeRHSFunction(ts,t0+c*dt,work[0],F);
121: VecAXPBYPCZ(work[0],1.*n/(2*n-1.),(n-1.)*dt/(r*(2*n-1)),(n-1.)/(2*n-1.),work[1],F);
122: i++;
123: }
124: for ( ; i<s; i++) {
125: c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
126: TSComputeRHSFunction(ts,t0+c*dt,work[0],F);
127: VecAXPY(work[0],dt/r,F);
128: }
129: VecCopy(work[0],sol);
130: SSPRestoreWorkVectors(ts,3,&work);
131: return(0);
132: }
136: /* Optimal fourth order SSP Runge-Kutta, low-storage (2N), c_eff=0.6 */
137: /* SSPRK(10,4), Pseudocode 3 of Ketcheson 2008 */
138: static PetscErrorCode SSPStep_RK_10_4(TS ts,PetscReal t0,PetscReal dt,Vec sol)
139: {
140: TS_SSP *ssp = (TS_SSP*)ts->data;
141: const PetscReal c[10] = {0, 1./6, 2./6, 3./6, 4./6, 2./6, 3./6, 4./6, 5./6, 1};
142: Vec *work,F;
143: PetscInt i,s;
147: s = ssp->nstages;
148: SSPGetWorkVectors(ts,3,&work);
149: F = work[2];
150: VecCopy(sol,work[0]);
151: for (i=0; i<5; i++) {
152: TSComputeRHSFunction(ts,t0+c[i],work[0],F);
153: VecAXPY(work[0],dt/6,F);
154: }
155: VecAXPBYPCZ(work[1],1./25,9./25,0,sol,work[0]);
156: VecAXPBY(work[0],15,-5,work[1]);
157: for ( ; i<9; i++) {
158: TSComputeRHSFunction(ts,t0+c[i],work[0],F);
159: VecAXPY(work[0],dt/6,F);
160: }
161: TSComputeRHSFunction(ts,t0+dt,work[0],F);
162: VecAXPBYPCZ(work[1],3./5,dt/10,1,work[0],F);
163: VecCopy(work[1],sol);
164: SSPRestoreWorkVectors(ts,3,&work);
165: return(0);
166: }
171: static PetscErrorCode TSSetUp_SSP(TS ts)
172: {
173: /* TS_SSP *ssp = (TS_SSP*)ts->data; */
177: return(0);
178: }
182: static PetscErrorCode TSStep_SSP(TS ts,PetscInt *steps,PetscReal *ptime)
183: {
184: TS_SSP *ssp = (TS_SSP*)ts->data;
185: Vec sol = ts->vec_sol;
187: PetscInt i,max_steps = ts->max_steps;
190: *steps = -ts->steps;
191: TSMonitor(ts,ts->steps,ts->ptime,sol);
193: for (i=0; i<max_steps; i++) {
194: PetscReal dt = ts->time_step;
196: TSPreStep(ts);
197: ts->ptime += dt;
198: (*ssp->onestep)(ts,ts->ptime-dt,dt,sol);
199: ts->steps++;
200: TSPostStep(ts);
201: TSMonitor(ts,ts->steps,ts->ptime,sol);
202: if (ts->ptime > ts->max_time) break;
203: }
205: *steps += ts->steps;
206: *ptime = ts->ptime;
207: return(0);
208: }
209: /*------------------------------------------------------------*/
212: static PetscErrorCode TSDestroy_SSP(TS ts)
213: {
214: TS_SSP *ssp = (TS_SSP*)ts->data;
218: if (ssp->work) {VecDestroyVecs(ssp->work,ssp->nwork);}
219: PetscFree(ssp);
220: return(0);
221: }
222: /*------------------------------------------------------------*/
226: static PetscErrorCode TSSSPSetType(TS ts,const TSSSPType type)
227: {
228: PetscErrorCode ierr,(*r)(TS,PetscReal,PetscReal,Vec);
229: TS_SSP *ssp = (TS_SSP*)ts->data;
232: PetscFListFind(TSSSPList,((PetscObject)ts)->comm,type,(void(**)(void))&r);
233: if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TS_SSP type %s given",type);
234: ssp->onestep = r;
235: return(0);
236: }
240: static PetscErrorCode TSSetFromOptions_SSP(TS ts)
241: {
242: char tname[256] = TSSSPRKS2;
243: TS_SSP *ssp = (TS_SSP*)ts->data;
245: PetscTruth flg;
248: PetscOptionsHead("SSP ODE solver options");
249: {
250: PetscOptionsList("-ts_ssp_type","Type of SSP method","TSSSPSetType",TSSSPList,tname,tname,sizeof(tname),&flg);
251: if (flg) {
252: TSSSPSetType(ts,tname);
253: }
254: PetscOptionsInt("-ts_ssp_nstages","Number of stages","TSSSPSetNumStages",ssp->nstages,&ssp->nstages,PETSC_NULL);
255: }
256: PetscOptionsTail();
257: return(0);
258: }
262: static PetscErrorCode TSView_SSP(TS ts,PetscViewer viewer)
263: {
265: return(0);
266: }
268: /* ------------------------------------------------------------ */
270: /*MC
271: TSSSP - Explicit strong stability preserving ODE solver
273: Most hyperbolic conservation laws have exact solutions that are total variation diminishing (TVD) or total variation
274: bounded (TVB) although these solutions often contain discontinuities. Spatial discretizations such as Godunov's
275: scheme and high-resolution finite volume methods (TVD limiters, ENO/WENO) are designed to preserve these properties,
276: but they are usually formulated using a forward Euler time discretization or by coupling the space and time
277: discretization as in the classical Lax-Wendroff scheme. When the space and time discretization is coupled, it is very
278: difficult to produce schemes with high temporal accuracy while preserving TVD properties. An alternative is the
279: semidiscrete formulation where we choose a spatial discretization that is TVD with forward Euler and then choose a
280: time discretization that preserves the TVD property. Such integrators are called strong stability preserving (SSP).
282: Let c_eff be the minimum number of function evaluations required to step as far as one step of forward Euler while
283: still being SSP. Some theoretical bounds
285: 1. There are no explicit methods with c_eff > 1.
287: 2. There are no explicit methods beyond order 4 (for nonlinear problems) and c_eff > 0.
289: 3. There are no implicit methods with order greater than 1 and c_eff > 2.
291: This integrator provides Runge-Kutta methods of order 2, 3, and 4 with maximal values of c_eff. More stages allows
292: for larger values of c_eff which improves efficiency. These implementations are low-memory and only use 2 or 3 work
293: vectors regardless of the total number of stages, so e.g. 25-stage 3rd order methods may be an excellent choice.
295: Methods can be chosen with -ts_ssp_type {rks2,rks3,rk104}
297: rks2: Second order methods with any number s>1 of stages. c_eff = (s-1)/s
299: rks3: Third order methods with s=n^2 stages, n>1. c_eff = (s-n)/s
301: rk104: A 10-stage fourth order method. c_eff = 0.6
303: Level: beginner
305: .seealso: TSCreate(), TS, TSSetType()
307: M*/
311: PetscErrorCode TSCreate_SSP(TS ts)
312: {
313: TS_SSP *ssp;
317: if (!TSSSPList) {
318: PetscFListAdd(&TSSSPList,TSSSPRKS2, "SSPStep_RK_2", (void(*)(void))SSPStep_RK_2);
319: PetscFListAdd(&TSSSPList,TSSSPRKS3, "SSPStep_RK_3", (void(*)(void))SSPStep_RK_3);
320: PetscFListAdd(&TSSSPList,TSSSPRK104, "SSPStep_RK_10_4",(void(*)(void))SSPStep_RK_10_4);
321: }
323: ts->ops->setup = TSSetUp_SSP;
324: ts->ops->step = TSStep_SSP;
325: ts->ops->destroy = TSDestroy_SSP;
326: ts->ops->setfromoptions = TSSetFromOptions_SSP;
327: ts->ops->view = TSView_SSP;
329: PetscNewLog(ts,TS_SSP,&ssp);
330: ts->data = (void*)ssp;
332: TSSSPSetType(ts,TSSSPRKS2);
333: ssp->nstages = 5;
334: return(0);
335: }