Actual source code: ssp.c

  1: /*
  2:        Code for Timestepping with explicit SSP.
  3: */
  4: #include <private/tsimpl.h>                /*I   "petscts.h"   I*/

  6: PetscFList TSSSPList = 0;

  8: typedef struct {
  9:   PetscErrorCode (*onestep)(TS,PetscReal,PetscReal,Vec);
 10:   char *type_name;
 11:   PetscInt nstages;
 12:   Vec *work;
 13:   PetscInt nwork;
 14:   PetscBool  workout;
 15: } TS_SSP;


 20: static PetscErrorCode SSPGetWorkVectors(TS ts,PetscInt n,Vec **work)
 21: {
 22:   TS_SSP *ssp = (TS_SSP*)ts->data;

 26:   if (ssp->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Work vectors already gotten");
 27:   if (ssp->nwork < n) {
 28:     if (ssp->nwork > 0) {
 29:       VecDestroyVecs(ssp->nwork,&ssp->work);
 30:     }
 31:     VecDuplicateVecs(ts->vec_sol,n,&ssp->work);
 32:     ssp->nwork = n;
 33:   }
 34:   *work = ssp->work;
 35:   ssp->workout = PETSC_TRUE;
 36:   return(0);
 37: }

 41: static PetscErrorCode SSPRestoreWorkVectors(TS ts,PetscInt n,Vec **work)
 42: {
 43:   TS_SSP *ssp = (TS_SSP*)ts->data;

 46:   if (!ssp->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Work vectors have not been gotten");
 47:   if (*work != ssp->work) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong work vectors checked out");
 48:   ssp->workout = PETSC_FALSE;
 49:   *work = PETSC_NULL;
 50:   return(0);
 51: }


 56: /*MC
 57:    TSSSPRKS2 - Optimal second order SSP Runge-Kutta method, low-storage, c_eff=(s-1)/s

 59:    Pseudocode 2 of Ketcheson 2008

 61:    Level: beginner

 63: .seealso: TSSSP, TSSSPSetType(), TSSSPSetNumStages()
 64: M*/
 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: /*MC
 90:    TSSSPRKS3 - Optimal third order SSP Runge-Kutta, low-storage, c_eff=(sqrt(s)-1)/sqrt(s), where sqrt(s) is an integer

 92:    Pseudocode 2 of Ketcheson 2008

 94:    Level: beginner

 96: .seealso: TSSSP, TSSSPSetType(), TSSSPSetNumStages()
 97: M*/
 98: static PetscErrorCode SSPStep_RK_3(TS ts,PetscReal t0,PetscReal dt,Vec sol)
 99: {
100:   TS_SSP *ssp = (TS_SSP*)ts->data;
101:   Vec *work,F;
102:   PetscInt i,s,n,r;
103:   PetscReal c;

107:   s = ssp->nstages;
108:   n = (PetscInt)(sqrt((PetscReal)s)+0.001);
109:   r = s-n;
110:   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);
111:   SSPGetWorkVectors(ts,3,&work);
112:   F = work[2];
113:   VecCopy(sol,work[0]);
114:   for (i=0; i<(n-1)*(n-2)/2; i++) {
115:     c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
116:     TSComputeRHSFunction(ts,t0+c*dt,work[0],F);
117:     VecAXPY(work[0],dt/r,F);
118:   }
119:   VecCopy(work[0],work[1]);
120:   for ( ; i<n*(n+1)/2-1; i++) {
121:     c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
122:     TSComputeRHSFunction(ts,t0+c*dt,work[0],F);
123:     VecAXPY(work[0],dt/r,F);
124:   }
125:   {
126:     c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
127:     TSComputeRHSFunction(ts,t0+c*dt,work[0],F);
128:     VecAXPBYPCZ(work[0],1.*n/(2*n-1.),(n-1.)*dt/(r*(2*n-1)),(n-1.)/(2*n-1.),work[1],F);
129:     i++;
130:   }
131:   for ( ; i<s; i++) {
132:     c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n);
133:     TSComputeRHSFunction(ts,t0+c*dt,work[0],F);
134:     VecAXPY(work[0],dt/r,F);
135:   }
136:   VecCopy(work[0],sol);
137:   SSPRestoreWorkVectors(ts,3,&work);
138:   return(0);
139: }

143: /*MC
144:    TSSSPRKS104 - Optimal fourth order SSP Runge-Kutta, low-storage (2N), c_eff=0.6

146:    SSPRK(10,4), Pseudocode 3 of Ketcheson 2008

148:    Level: beginner

150: .seealso: TSSSP, TSSSPSetType()
151: M*/
152: static PetscErrorCode SSPStep_RK_10_4(TS ts,PetscReal t0,PetscReal dt,Vec sol)
153: {
154:   const PetscReal c[10] = {0, 1./6, 2./6, 3./6, 4./6, 2./6, 3./6, 4./6, 5./6, 1};
155:   Vec *work,F;
156:   PetscInt i;

160:   SSPGetWorkVectors(ts,3,&work);
161:   F = work[2];
162:   VecCopy(sol,work[0]);
163:   for (i=0; i<5; i++) {
164:     TSComputeRHSFunction(ts,t0+c[i],work[0],F);
165:     VecAXPY(work[0],dt/6,F);
166:   }
167:   VecAXPBYPCZ(work[1],1./25,9./25,0,sol,work[0]);
168:   VecAXPBY(work[0],15,-5,work[1]);
169:   for ( ; i<9; i++) {
170:     TSComputeRHSFunction(ts,t0+c[i],work[0],F);
171:     VecAXPY(work[0],dt/6,F);
172:   }
173:   TSComputeRHSFunction(ts,t0+dt,work[0],F);
174:   VecAXPBYPCZ(work[1],3./5,dt/10,1,work[0],F);
175:   VecCopy(work[1],sol);
176:   SSPRestoreWorkVectors(ts,3,&work);
177:   return(0);
178: }


183: static PetscErrorCode TSSetUp_SSP(TS ts)
184: {

187:   return(0);
188: }

192: static PetscErrorCode TSStep_SSP(TS ts)
193: {
194:   TS_SSP        *ssp = (TS_SSP*)ts->data;
195:   Vec            sol = ts->vec_sol;

199:   (*ssp->onestep)(ts,ts->ptime,ts->time_step,sol);
200:   ts->ptime += ts->time_step;
201:   ts->steps++;
202:   return(0);
203: }
204: /*------------------------------------------------------------*/
207: static PetscErrorCode TSReset_SSP(TS ts)
208: {
209:   TS_SSP         *ssp = (TS_SSP*)ts->data;

213:   if (ssp->work) {VecDestroyVecs(ssp->nwork,&ssp->work);}
214:   ssp->nwork = 0;
215:   ssp->workout = PETSC_FALSE;
216:   return(0);
217: }

221: static PetscErrorCode TSDestroy_SSP(TS ts)
222: {
223:   TS_SSP         *ssp = (TS_SSP*)ts->data;

227:   TSReset_SSP(ts);
228:   PetscFree(ssp->type_name);
229:   PetscFree(ts->data);
230:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSSPGetType_C","",PETSC_NULL);
231:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSSPSetType_C","",PETSC_NULL);
232:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSSPGetNumStages_C","",PETSC_NULL);
233:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSSPSetNumStages_C","",PETSC_NULL);
234:   return(0);
235: }
236: /*------------------------------------------------------------*/

240: /*@C
241:    TSSSPSetType - set the SSP time integration scheme to use

243:    Logically Collective

245:    Input Arguments:
246:    ts - time stepping object
247:    type - type of scheme to use

249:    Options Database Keys:
250:    -ts_ssp_type <rks2>: Type of SSP method (one of) rks2 rks3 rk104
251:    -ts_ssp_nstages <5>: Number of stages

253:    Level: beginner

255: .seealso: TSSSP, TSSSPGetType(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104
256: @*/
257: PetscErrorCode TSSSPSetType(TS ts,const TSSSPType type)
258: {

263:   PetscTryMethod(ts,"TSSSPSetType_C",(TS,const TSSSPType),(ts,type));
264:   return(0);
265: }

269: /*@C
270:    TSSSPGetType - get the SSP time integration scheme

272:    Logically Collective

274:    Input Argument:
275:    ts - time stepping object

277:    Output Argument:
278:    type - type of scheme being used

280:    Level: beginner

282: .seealso: TSSSP, TSSSPSettype(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104
283: @*/
284: PetscErrorCode TSSSPGetType(TS ts,const TSSSPType *type)
285: {

290:   PetscTryMethod(ts,"TSSSPGetType_C",(TS,const TSSSPType*),(ts,type));
291:   return(0);
292: }

296: /*@
297:    TSSSPSetNumStages - set the number of stages to use with the SSP method

299:    Logically Collective

301:    Input Arguments:
302:    ts - time stepping object
303:    nstages - number of stages

305:    Options Database Keys:
306:    -ts_ssp_type <rks2>: NumStages of SSP method (one of) rks2 rks3 rk104
307:    -ts_ssp_nstages <5>: Number of stages

309:    Level: beginner

311: .seealso: TSSSP, TSSSPGetNumStages(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104
312: @*/
313: PetscErrorCode TSSSPSetNumStages(TS ts,PetscInt nstages)
314: {

319:   PetscTryMethod(ts,"TSSSPSetNumStages_C",(TS,PetscInt),(ts,nstages));
320:   return(0);
321: }

325: /*@
326:    TSSSPGetNumStages - get the number of stages in the SSP time integration scheme

328:    Logically Collective

330:    Input Argument:
331:    ts - time stepping object

333:    Output Argument:
334:    nstages - number of stages

336:    Level: beginner

338: .seealso: TSSSP, TSSSPGetType(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104
339: @*/
340: PetscErrorCode TSSSPGetNumStages(TS ts,PetscInt *nstages)
341: {

346:   PetscTryMethod(ts,"TSSSPGetNumStages_C",(TS,PetscInt*),(ts,nstages));
347:   return(0);
348: }

353: PetscErrorCode TSSSPSetType_SSP(TS ts,const TSSSPType type)
354: {
355:   PetscErrorCode ierr,(*r)(TS,PetscReal,PetscReal,Vec);
356:   TS_SSP *ssp = (TS_SSP*)ts->data;

359:   PetscFListFind(TSSSPList,((PetscObject)ts)->comm,type,PETSC_TRUE,(PetscVoidStarFunction)&r);
360:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TS_SSP type %s given",type);
361:   ssp->onestep = r;
362:   PetscFree(ssp->type_name);
363:   PetscStrallocpy(type,&ssp->type_name);
364:   return(0);
365: }
368: PetscErrorCode TSSSPGetType_SSP(TS ts,const TSSSPType *type)
369: {
370:   TS_SSP *ssp = (TS_SSP*)ts->data;

373:   *type = ssp->type_name;
374:   return(0);
375: }
378: PetscErrorCode TSSSPSetNumStages_SSP(TS ts,PetscInt nstages)
379: {
380:   TS_SSP *ssp = (TS_SSP*)ts->data;

383:   ssp->nstages = nstages;
384:   return(0);
385: }
388: PetscErrorCode TSSSPGetNumStages_SSP(TS ts,PetscInt *nstages)
389: {
390:   TS_SSP *ssp = (TS_SSP*)ts->data;

393:   *nstages = ssp->nstages;
394:   return(0);
395: }

400: static PetscErrorCode TSSetFromOptions_SSP(TS ts)
401: {
402:   char tname[256] = TSSSPRKS2;
403:   TS_SSP *ssp = (TS_SSP*)ts->data;
405:   PetscBool  flg;

408:   PetscOptionsHead("SSP ODE solver options");
409:   {
410:     PetscOptionsList("-ts_ssp_type","Type of SSP method","TSSSPSetType",TSSSPList,tname,tname,sizeof(tname),&flg);
411:     if (flg) {
412:       TSSSPSetType(ts,tname);
413:     }
414:     PetscOptionsInt("-ts_ssp_nstages","Number of stages","TSSSPSetNumStages",ssp->nstages,&ssp->nstages,PETSC_NULL);
415:   }
416:   PetscOptionsTail();
417:   return(0);
418: }

422: static PetscErrorCode TSView_SSP(TS ts,PetscViewer viewer)
423: {
425:   return(0);
426: }

428: /* ------------------------------------------------------------ */

430: /*MC
431:       TSSSP - Explicit strong stability preserving ODE solver

433:   Most hyperbolic conservation laws have exact solutions that are total variation diminishing (TVD) or total variation
434:   bounded (TVB) although these solutions often contain discontinuities.  Spatial discretizations such as Godunov's
435:   scheme and high-resolution finite volume methods (TVD limiters, ENO/WENO) are designed to preserve these properties,
436:   but they are usually formulated using a forward Euler time discretization or by coupling the space and time
437:   discretization as in the classical Lax-Wendroff scheme.  When the space and time discretization is coupled, it is very
438:   difficult to produce schemes with high temporal accuracy while preserving TVD properties.  An alternative is the
439:   semidiscrete formulation where we choose a spatial discretization that is TVD with forward Euler and then choose a
440:   time discretization that preserves the TVD property.  Such integrators are called strong stability preserving (SSP).

442:   Let c_eff be the minimum number of function evaluations required to step as far as one step of forward Euler while
443:   still being SSP.  Some theoretical bounds

445:   1. There are no explicit methods with c_eff > 1.

447:   2. There are no explicit methods beyond order 4 (for nonlinear problems) and c_eff > 0.

449:   3. There are no implicit methods with order greater than 1 and c_eff > 2.

451:   This integrator provides Runge-Kutta methods of order 2, 3, and 4 with maximal values of c_eff.  More stages allows
452:   for larger values of c_eff which improves efficiency.  These implementations are low-memory and only use 2 or 3 work
453:   vectors regardless of the total number of stages, so e.g. 25-stage 3rd order methods may be an excellent choice.

455:   Methods can be chosen with -ts_ssp_type {rks2,rks3,rk104}

457:   rks2: Second order methods with any number s>1 of stages.  c_eff = (s-1)/s

459:   rks3: Third order methods with s=n^2 stages, n>1.  c_eff = (s-n)/s

461:   rk104: A 10-stage fourth order method.  c_eff = 0.6

463:   Level: beginner

465:   References:
466:   Ketcheson, Highly efficient strong stability preserving Runge-Kutta methods with low-storage implementations, SISC, 2008.

468:   Gottlieb, Ketcheson, and Shu, High order strong stability preserving time discretizations, J Scientific Computing, 2009.

470: .seealso:  TSCreate(), TS, TSSetType()

472: M*/
476: PetscErrorCode  TSCreate_SSP(TS ts)
477: {
478:   TS_SSP       *ssp;

482:   if (!TSSSPList) {
483:     PetscFListAdd(&TSSSPList,TSSSPRKS2,  "SSPStep_RK_2",   (void(*)(void))SSPStep_RK_2);
484:     PetscFListAdd(&TSSSPList,TSSSPRKS3,  "SSPStep_RK_3",   (void(*)(void))SSPStep_RK_3);
485:     PetscFListAdd(&TSSSPList,TSSSPRK104, "SSPStep_RK_10_4",(void(*)(void))SSPStep_RK_10_4);
486:   }

488:   ts->ops->setup           = TSSetUp_SSP;
489:   ts->ops->step            = TSStep_SSP;
490:   ts->ops->reset           = TSReset_SSP;
491:   ts->ops->destroy         = TSDestroy_SSP;
492:   ts->ops->setfromoptions  = TSSetFromOptions_SSP;
493:   ts->ops->view            = TSView_SSP;

495:   PetscNewLog(ts,TS_SSP,&ssp);
496:   ts->data = (void*)ssp;

498:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSSPGetType_C","TSSSPGetType_SSP",TSSSPGetType_SSP);
499:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSSPSetType_C","TSSSPSetType_SSP",TSSSPSetType_SSP);
500:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSSPGetNumStages_C","TSSSPGetNumStages_SSP",TSSSPGetNumStages_SSP);
501:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSSPSetNumStages_C","TSSSPSetNumStages_SSP",TSSSPSetNumStages_SSP);

503:   TSSSPSetType(ts,TSSSPRKS2);
504:   ssp->nstages = 5;
505:   return(0);
506: }