Actual source code: ex9.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  1: static const char help[] = "1D periodic Finite Volume solver in slope-limiter form with semidiscrete time stepping.\n"
  2:   "Solves scalar and vector problems, choose the physical model with -physics\n"
  3:   "  advection   - Constant coefficient scalar advection\n"
  4:   "                u_t       + (a*u)_x               = 0\n"
  5:   "  burgers     - Burgers equation\n"
  6:   "                u_t       + (u^2/2)_x             = 0\n"
  7:   "  traffic     - Traffic equation\n"
  8:   "                u_t       + (u*(1-u))_x           = 0\n"
  9:   "  acoustics   - Acoustic wave propagation\n"
 10:   "                u_t       + (c*z*v)_x             = 0\n"
 11:   "                v_t       + (c/z*u)_x             = 0\n"
 12:   "  isogas      - Isothermal gas dynamics\n"
 13:   "                rho_t     + (rho*u)_x             = 0\n"
 14:   "                (rho*u)_t + (rho*u^2 + c^2*rho)_x = 0\n"
 15:   "  shallow     - Shallow water equations\n"
 16:   "                h_t       + (h*u)_x               = 0\n"
 17:   "                (h*u)_t   + (h*u^2 + g*h^2/2)_x   = 0\n"
 18:   "Some of these physical models have multiple Riemann solvers, select these with -physics_xxx_riemann\n"
 19:   "  exact       - Exact Riemann solver which usually needs to perform a Newton iteration to connect\n"
 20:   "                the states across shocks and rarefactions\n"
 21:   "  roe         - Linearized scheme, usually with an entropy fix inside sonic rarefactions\n"
 22:   "The systems provide a choice of reconstructions with -physics_xxx_reconstruct\n"
 23:   "  characteristic - Limit the characteristic variables, this is usually preferred (default)\n"
 24:   "  conservative   - Limit the conservative variables directly, can cause undesired interaction of waves\n\n"
 25:   "A variety of limiters for high-resolution TVD limiters are available with -limit\n"
 26:   "  upwind,minmod,superbee,mc,vanleer,vanalbada,koren,cada-torillhon (last two are nominally third order)\n"
 27:   "  and non-TVD schemes lax-wendroff,beam-warming,fromm\n\n"
 28:   "To preserve the TVD property, one should time step with a strong stability preserving method.\n"
 29:   "The optimal high order explicit Runge-Kutta methods in TSSSP are recommended for non-stiff problems.\n\n"
 30:   "Several initial conditions can be chosen with -initial N\n\n"
 31:   "The problem size should be set with -da_grid_x M\n\n";

 33: #include <petscts.h>
 34: #include <petscdm.h>
 35: #include <petscdmda.h>
 36: #include <petscdraw.h>

 38: #include <petsc/private/kernels/blockinvert.h> /* For the Kernel_*_gets_* stuff for BAIJ */

 40: PETSC_STATIC_INLINE PetscReal Sgn(PetscReal a) { return (a<0) ? -1 : 1; }
 41: PETSC_STATIC_INLINE PetscReal Abs(PetscReal a) { return (a<0) ? 0 : a; }
 42: PETSC_STATIC_INLINE PetscReal Sqr(PetscReal a) { return a*a; }
 43: PETSC_STATIC_INLINE PetscReal MaxAbs(PetscReal a,PetscReal b) { return (PetscAbs(a) > PetscAbs(b)) ? a : b; }
 44: PETSC_UNUSED PETSC_STATIC_INLINE PetscReal MinAbs(PetscReal a,PetscReal b) { return (PetscAbs(a) < PetscAbs(b)) ? a : b; }
 45: PETSC_STATIC_INLINE PetscReal MinMod2(PetscReal a,PetscReal b) { return (a*b<0) ? 0 : Sgn(a)*PetscMin(PetscAbs(a),PetscAbs(b)); }
 46: PETSC_STATIC_INLINE PetscReal MaxMod2(PetscReal a,PetscReal b) { return (a*b<0) ? 0 : Sgn(a)*PetscMax(PetscAbs(a),PetscAbs(b)); }
 47: PETSC_STATIC_INLINE PetscReal MinMod3(PetscReal a,PetscReal b,PetscReal c) {return (a*b<0 || a*c<0) ? 0 : Sgn(a)*PetscMin(PetscAbs(a),PetscMin(PetscAbs(b),PetscAbs(c))); }

 49: PETSC_STATIC_INLINE PetscReal RangeMod(PetscReal a,PetscReal xmin,PetscReal xmax) { PetscReal range = xmax-xmin; return xmin +PetscFmodReal(range+PetscFmodReal(a,range),range); }


 52: /* ----------------------- Lots of limiters, these could go in a separate library ------------------------- */
 53: typedef struct _LimitInfo {
 54:   PetscReal hx;
 55:   PetscInt  m;
 56: } *LimitInfo;
 57: static void Limit_Upwind(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
 58: {
 59:   PetscInt i;
 60:   for (i=0; i<info->m; i++) lmt[i] = 0;
 61: }
 62: static void Limit_LaxWendroff(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
 63: {
 64:   PetscInt i;
 65:   for (i=0; i<info->m; i++) lmt[i] = jR[i];
 66: }
 67: static void Limit_BeamWarming(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
 68: {
 69:   PetscInt i;
 70:   for (i=0; i<info->m; i++) lmt[i] = jL[i];
 71: }
 72: static void Limit_Fromm(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
 73: {
 74:   PetscInt i;
 75:   for (i=0; i<info->m; i++) lmt[i] = 0.5*(jL[i] + jR[i]);
 76: }
 77: static void Limit_Minmod(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
 78: {
 79:   PetscInt i;
 80:   for (i=0; i<info->m; i++) lmt[i] = MinMod2(jL[i],jR[i]);
 81: }
 82: static void Limit_Superbee(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
 83: {
 84:   PetscInt i;
 85:   for (i=0; i<info->m; i++) lmt[i] = MaxMod2(MinMod2(jL[i],2*jR[i]),MinMod2(2*jL[i],jR[i]));
 86: }
 87: static void Limit_MC(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
 88: {
 89:   PetscInt i;
 90:   for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i],0.5*(jL[i]+jR[i]),2*jR[i]);
 91: }
 92: static void Limit_VanLeer(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
 93: { /* phi = (t + abs(t)) / (1 + abs(t)) */
 94:   PetscInt i;
 95:   for (i=0; i<info->m; i++) lmt[i] = (jL[i]*Abs(jR[i]) + Abs(jL[i])*jR[i]) / (Abs(jL[i]) + Abs(jR[i]) + 1e-15);
 96: }
 97: static void Limit_VanAlbada(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) /* differentiable */
 98: { /* phi = (t + t^2) / (1 + t^2) */
 99:   PetscInt i;
100:   for (i=0; i<info->m; i++) lmt[i] = (jL[i]*Sqr(jR[i]) + Sqr(jL[i])*jR[i]) / (Sqr(jL[i]) + Sqr(jR[i]) + 1e-15);
101: }
102: static void Limit_VanAlbadaTVD(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
103: { /* phi = (t + t^2) / (1 + t^2) */
104:   PetscInt i;
105:   for (i=0; i<info->m; i++) lmt[i] = (jL[i]*jR[i]<0) ? 0 : (jL[i]*Sqr(jR[i]) + Sqr(jL[i])*jR[i]) / (Sqr(jL[i]) + Sqr(jR[i]) + 1e-15);
106: }
107: static void Limit_Koren(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) /* differentiable */
108: { /* phi = (t + 2*t^2) / (2 - t + 2*t^2) */
109:   PetscInt i;
110:   for (i=0; i<info->m; i++) lmt[i] = ((jL[i]*Sqr(jR[i]) + 2*Sqr(jL[i])*jR[i])/(2*Sqr(jL[i]) - jL[i]*jR[i] + 2*Sqr(jR[i]) + 1e-15));
111: }
112: static void Limit_KorenSym(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) /* differentiable */
113: { /* Symmetric version of above */
114:   PetscInt i;
115:   for (i=0; i<info->m; i++) lmt[i] = (1.5*(jL[i]*Sqr(jR[i]) + Sqr(jL[i])*jR[i])/(2*Sqr(jL[i]) - jL[i]*jR[i] + 2*Sqr(jR[i]) + 1e-15));
116: }
117: static void Limit_Koren3(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
118: { /* Eq 11 of Cada-Torrilhon 2009 */
119:   PetscInt i;
120:   for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i],(jL[i]+2*jR[i])/3,2*jR[i]);
121: }
122: static PetscReal CadaTorrilhonPhiHatR_Eq13(PetscReal L,PetscReal R)
123: {
124:   return PetscMax(0,PetscMin((L+2*R)/3,PetscMax(-0.5*L,PetscMin(2*L,PetscMin((L+2*R)/3,1.6*R)))));
125: }
126: static void Limit_CadaTorrilhon2(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
127: { /* Cada-Torrilhon 2009, Eq 13 */
128:   PetscInt i;
129:   for (i=0; i<info->m; i++) lmt[i] = CadaTorrilhonPhiHatR_Eq13(jL[i],jR[i]);
130: }
131: static void Limit_CadaTorrilhon3R(PetscReal r,LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
132: { /* Cada-Torrilhon 2009, Eq 22 */
133:   /* They recommend 0.001 < r < 1, but larger values are more accurate in smooth regions */
134:   const PetscReal eps = 1e-7,hx = info->hx;
135:   PetscInt        i;
136:   for (i=0; i<info->m; i++) {
137:     const PetscReal eta = (Sqr(jL[i]) + Sqr(jR[i])) / Sqr(r*hx);
138:     lmt[i] = ((eta < 1-eps) ? (jL[i] + 2*jR[i]) / 3 : ((eta > 1+eps) ? CadaTorrilhonPhiHatR_Eq13(jL[i],jR[i]) : 0.5*((1-(eta-1)/eps)*(jL[i]+2*jR[i])/3 + (1+(eta+1)/eps)*CadaTorrilhonPhiHatR_Eq13(jL[i],jR[i]))));
139:   }
140: }
141: static void Limit_CadaTorrilhon3R0p1(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
142: {
143:   Limit_CadaTorrilhon3R(0.1,info,jL,jR,lmt);
144: }
145: static void Limit_CadaTorrilhon3R1(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
146: {
147:   Limit_CadaTorrilhon3R(1,info,jL,jR,lmt);
148: }
149: static void Limit_CadaTorrilhon3R10(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
150: {
151:   Limit_CadaTorrilhon3R(10,info,jL,jR,lmt);
152: }
153: static void Limit_CadaTorrilhon3R100(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
154: {
155:   Limit_CadaTorrilhon3R(100,info,jL,jR,lmt);
156: }


159: /* --------------------------------- Finite Volume data structures ----------------------------------- */

161: typedef enum {FVBC_PERIODIC, FVBC_OUTFLOW} FVBCType;
162: static const char *FVBCTypes[] = {"PERIODIC","OUTFLOW","FVBCType","FVBC_",0};
163: typedef PetscErrorCode (*RiemannFunction)(void*,PetscInt,const PetscScalar*,const PetscScalar*,PetscScalar*,PetscReal*);
164: typedef PetscErrorCode (*ReconstructFunction)(void*,PetscInt,const PetscScalar*,PetscScalar*,PetscScalar*,PetscReal*);

166: typedef struct {
167:   PetscErrorCode      (*sample)(void*,PetscInt,FVBCType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal*);
168:   RiemannFunction     riemann;
169:   ReconstructFunction characteristic;
170:   PetscErrorCode      (*destroy)(void*);
171:   void                *user;
172:   PetscInt            dof;
173:   char                *fieldname[16];
174: } PhysicsCtx;

176: typedef struct {
177:   void        (*limit)(LimitInfo,const PetscScalar*,const PetscScalar*,PetscScalar*);
178:   PhysicsCtx  physics;
179:   MPI_Comm    comm;
180:   char        prefix[256];

182:   /* Local work arrays */
183:   PetscScalar *R,*Rinv;         /* Characteristic basis, and it's inverse.  COLUMN-MAJOR */
184:   PetscScalar *cjmpLR;          /* Jumps at left and right edge of cell, in characteristic basis, len=2*dof */
185:   PetscScalar *cslope;          /* Limited slope, written in characteristic basis */
186:   PetscScalar *uLR;             /* Solution at left and right of interface, conservative variables, len=2*dof */
187:   PetscScalar *flux;            /* Flux across interface */
188:   PetscReal   *speeds;          /* Speeds of each wave */

190:   PetscReal   cfl_idt;            /* Max allowable value of 1/Delta t */
191:   PetscReal   cfl;
192:   PetscReal   xmin,xmax;
193:   PetscInt    initial;
194:   PetscBool   exact;
195:   FVBCType    bctype;
196: } FVCtx;

200: PetscErrorCode RiemannListAdd(PetscFunctionList *flist,const char *name,RiemannFunction rsolve)
201: {

205:   PetscFunctionListAdd(flist,name,rsolve);
206:   return(0);
207: }

211: PetscErrorCode RiemannListFind(PetscFunctionList flist,const char *name,RiemannFunction *rsolve)
212: {

216:   PetscFunctionListFind(flist,name,rsolve);
217:   if (!*rsolve) SETERRQ1(PETSC_COMM_SELF,1,"Riemann solver \"%s\" could not be found",name);
218:   return(0);
219: }

223: PetscErrorCode ReconstructListAdd(PetscFunctionList *flist,const char *name,ReconstructFunction r)
224: {

228:   PetscFunctionListAdd(flist,name,r);
229:   return(0);
230: }

234: PetscErrorCode ReconstructListFind(PetscFunctionList flist,const char *name,ReconstructFunction *r)
235: {

239:   PetscFunctionListFind(flist,name,r);
240:   if (!*r) SETERRQ1(PETSC_COMM_SELF,1,"Reconstruction \"%s\" could not be found",name);
241:   return(0);
242: }

244: /* --------------------------------- Physics ----------------------------------- */
245: /**
246: * Each physical model consists of Riemann solver and a function to determine the basis to use for reconstruction.  These
247: * are set with the PhysicsCreate_XXX function which allocates private storage and sets these methods as well as the
248: * number of fields and their names, and a function to deallocate private storage.
249: **/

251: /* First a few functions useful to several different physics */
254: static PetscErrorCode PhysicsCharacteristic_Conservative(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds)
255: {
256:   PetscInt i,j;

259:   for (i=0; i<m; i++) {
260:     for (j=0; j<m; j++) Xi[i*m+j] = X[i*m+j] = (PetscScalar)(i==j);
261:     speeds[i] = PETSC_MAX_REAL; /* Indicates invalid */
262:   }
263:   return(0);
264: }

268: static PetscErrorCode PhysicsDestroy_SimpleFree(void *vctx)
269: {

273:   PetscFree(vctx);
274:   return(0);
275: }



279: /* --------------------------------- Advection ----------------------------------- */

281: typedef struct {
282:   PetscReal a;                  /* advective velocity */
283: } AdvectCtx;

287: static PetscErrorCode PhysicsRiemann_Advect(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
288: {
289:   AdvectCtx *ctx = (AdvectCtx*)vctx;
290:   PetscReal speed;

293:   speed     = ctx->a;
294:   flux[0]   = PetscMax(0,speed)*uL[0] + PetscMin(0,speed)*uR[0];
295:   *maxspeed = speed;
296:   return(0);
297: }

301: static PetscErrorCode PhysicsCharacteristic_Advect(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds)
302: {
303:   AdvectCtx *ctx = (AdvectCtx*)vctx;

306:   X[0]      = 1.;
307:   Xi[0]     = 1.;
308:   speeds[0] = ctx->a;
309:   return(0);
310: }

314: static PetscErrorCode PhysicsSample_Advect(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u)
315: {
316:   AdvectCtx *ctx = (AdvectCtx*)vctx;
317:   PetscReal a    = ctx->a,x0;

320:   switch (bctype) {
321:     case FVBC_OUTFLOW:   x0 = x-a*t; break;
322:     case FVBC_PERIODIC: x0 = RangeMod(x-a*t,xmin,xmax); break;
323:     default: SETERRQ(PETSC_COMM_SELF,1,"unknown BCType");
324:   }
325:   switch (initial) {
326:     case 0: u[0] = (x0 < 0) ? 1 : -1; break;
327:     case 1: u[0] = (x0 < 0) ? -1 : 1; break;
328:     case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break;
329:     case 3: u[0] = PetscSinReal(2*PETSC_PI*x0); break;
330:     case 4: u[0] = PetscAbs(x0); break;
331:     case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2*PETSC_PI*x0)); break;
332:     case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2-x0 : 0)); break;
333:     default: SETERRQ(PETSC_COMM_SELF,1,"unknown initial condition");
334:   }
335:   return(0);
336: }

340: static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx)
341: {
343:   AdvectCtx      *user;

346:   PetscNew(&user);
347:   ctx->physics.sample         = PhysicsSample_Advect;
348:   ctx->physics.riemann        = PhysicsRiemann_Advect;
349:   ctx->physics.characteristic = PhysicsCharacteristic_Advect;
350:   ctx->physics.destroy        = PhysicsDestroy_SimpleFree;
351:   ctx->physics.user           = user;
352:   ctx->physics.dof            = 1;
353:   PetscStrallocpy("u",&ctx->physics.fieldname[0]);
354:   user->a = 1;
355:   PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for advection","");
356:   {
357:     PetscOptionsReal("-physics_advect_a","Speed","",user->a,&user->a,NULL);
358:   }
359:   PetscOptionsEnd();
360:   return(0);
361: }

363: /* --------------------------------- Burgers ----------------------------------- */

365: typedef struct {
366:   PetscReal lxf_speed;
367: } BurgersCtx;

371: static PetscErrorCode PhysicsSample_Burgers(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u)
372: {
374:   if (bctype == FVBC_PERIODIC && t > 0) SETERRQ(PETSC_COMM_SELF,1,"Exact solution not implemented for periodic");
375:   switch (initial) {
376:     case 0: u[0] = (x < 0) ? 1 : -1; break;
377:     case 1:
378:       if       (x < -t) u[0] = -1;
379:       else if  (x < t)  u[0] = x/t;
380:       else              u[0] = 1;
381:       break;
382:     case 2:
383:       if      (x < 0)       u[0] = 0;
384:       else if (x <= t)      u[0] = x/t;
385:       else if (x < 1+0.5*t) u[0] = 1;
386:       else                  u[0] = 0;
387:       break;
388:     case 3:
389:       if       (x < 0.2*t) u[0] = 0.2;
390:       else if  (x < t) u[0] = x/t;
391:       else             u[0] = 1;
392:       break;
393:     case 4:
394:       if (t > 0) SETERRQ(PETSC_COMM_SELF,1,"Only initial condition available");
395:       u[0] = 0.7 + 0.3*PetscSinReal(2*PETSC_PI*((x-xmin)/(xmax-xmin)));
396:       break;
397:     case 5:                     /* Pure shock solution */
398:       if (x < 0.5*t) u[0] = 1;
399:       else u[0] = 0;
400:       break;
401:     default: SETERRQ(PETSC_COMM_SELF,1,"unknown initial condition");
402:   }
403:   return(0);
404: }

408: static PetscErrorCode PhysicsRiemann_Burgers_Exact(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
409: {
411:   if (uL[0] < uR[0]) {          /* rarefaction */
412:     flux[0] = (uL[0]*uR[0] < 0)
413:       ? 0                       /* sonic rarefaction */
414:       : 0.5*PetscMin(PetscSqr(uL[0]),PetscSqr(uR[0]));
415:   } else {                      /* shock */
416:     flux[0] = 0.5*PetscMax(PetscSqr(uL[0]),PetscSqr(uR[0]));
417:   }
418:   *maxspeed = (PetscAbs(uL[0]) > PetscAbs(uR[0])) ? uL[0] : uR[0];
419:   return(0);
420: }

424: static PetscErrorCode PhysicsRiemann_Burgers_Roe(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
425: {
426:   PetscReal speed;

429:   speed   = 0.5*(uL[0] + uR[0]);
430:   flux[0] = 0.25*(PetscSqr(uL[0]) + PetscSqr(uR[0])) - 0.5*PetscAbs(speed)*(uR[0]-uL[0]);
431:   if (uL[0] <= 0 && 0 <= uR[0]) flux[0] = 0; /* Entropy fix for sonic rarefaction */
432:   *maxspeed = speed;
433:   return(0);
434: }

438: static PetscErrorCode PhysicsRiemann_Burgers_LxF(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
439: {
440:   PetscReal   c;
441:   PetscScalar fL,fR;

444:   c         = ((BurgersCtx*)vctx)->lxf_speed;
445:   fL        = 0.5*PetscSqr(uL[0]);
446:   fR        = 0.5*PetscSqr(uR[0]);
447:   flux[0]   = 0.5*(fL + fR) - 0.5*c*(uR[0] - uL[0]);
448:   *maxspeed = c;
449:   return(0);
450: }

454: static PetscErrorCode PhysicsRiemann_Burgers_Rusanov(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
455: {
456:   PetscReal   c;
457:   PetscScalar fL,fR;

460:   c         = PetscMax(PetscAbs(uL[0]),PetscAbs(uR[0]));
461:   fL        = 0.5*PetscSqr(uL[0]);
462:   fR        = 0.5*PetscSqr(uR[0]);
463:   flux[0]   = 0.5*(fL + fR) - 0.5*c*(uR[0] - uL[0]);
464:   *maxspeed = c;
465:   return(0);
466: }

470: static PetscErrorCode PhysicsCreate_Burgers(FVCtx *ctx)
471: {
472:   BurgersCtx        *user;
473:   PetscErrorCode    ierr;
474:   RiemannFunction   r;
475:   PetscFunctionList rlist      = 0;
476:   char              rname[256] = "exact";

479:   PetscNew(&user);

481:   ctx->physics.sample         = PhysicsSample_Burgers;
482:   ctx->physics.characteristic = PhysicsCharacteristic_Conservative;
483:   ctx->physics.destroy        = PhysicsDestroy_SimpleFree;
484:   ctx->physics.user           = user;
485:   ctx->physics.dof            = 1;

487:   PetscStrallocpy("u",&ctx->physics.fieldname[0]);
488:   RiemannListAdd(&rlist,"exact",  PhysicsRiemann_Burgers_Exact);
489:   RiemannListAdd(&rlist,"roe",    PhysicsRiemann_Burgers_Roe);
490:   RiemannListAdd(&rlist,"lxf",    PhysicsRiemann_Burgers_LxF);
491:   RiemannListAdd(&rlist,"rusanov",PhysicsRiemann_Burgers_Rusanov);
492:   PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for advection","");
493:   {
494:     PetscOptionsFList("-physics_burgers_riemann","Riemann solver","",rlist,rname,rname,sizeof(rname),NULL);
495:   }
496:   PetscOptionsEnd();
497:   RiemannListFind(rlist,rname,&r);
498:   PetscFunctionListDestroy(&rlist);
499:   ctx->physics.riemann = r;

501:   /* *
502:   * Hack to deal with LxF in semi-discrete form
503:   * max speed is 1 for the basic initial conditions (where |u| <= 1)
504:   * */
505:   if (r == PhysicsRiemann_Burgers_LxF) user->lxf_speed = 1;
506:   return(0);
507: }

509: /* --------------------------------- Traffic ----------------------------------- */

511: typedef struct {
512:   PetscReal lxf_speed;
513:   PetscReal a;
514: } TrafficCtx;

516: PETSC_STATIC_INLINE PetscScalar TrafficFlux(PetscScalar a,PetscScalar u) { return a*u*(1-u); }

520: static PetscErrorCode PhysicsSample_Traffic(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u)
521: {
522:   PetscReal a = ((TrafficCtx*)vctx)->a;

525:   if (bctype == FVBC_PERIODIC && t > 0) SETERRQ(PETSC_COMM_SELF,1,"Exact solution not implemented for periodic");
526:   switch (initial) {
527:     case 0:
528:       u[0] = (-a*t < x) ? 2 : 0; break;
529:     case 1:
530:       if      (x < PetscMin(2*a*t,0.5+a*t)) u[0] = -1;
531:       else if (x < 1)                       u[0] = 0;
532:       else                                  u[0] = 1;
533:       break;
534:     case 2:
535:       if (t > 0) SETERRQ(PETSC_COMM_SELF,1,"Only initial condition available");
536:       u[0] = 0.7 + 0.3*PetscSinReal(2*PETSC_PI*((x-xmin)/(xmax-xmin)));
537:       break;
538:     default: SETERRQ(PETSC_COMM_SELF,1,"unknown initial condition");
539:   }
540:   return(0);
541: }

545: static PetscErrorCode PhysicsRiemann_Traffic_Exact(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
546: {
547:   PetscReal a = ((TrafficCtx*)vctx)->a;

550:   if (uL[0] < uR[0]) {
551:     flux[0] = PetscMin(TrafficFlux(a,uL[0]),TrafficFlux(a,uR[0]));
552:   } else {
553:     flux[0] = (uR[0] < 0.5 && 0.5 < uL[0]) ? TrafficFlux(a,0.5) : PetscMax(TrafficFlux(a,uL[0]),TrafficFlux(a,uR[0]));
554:   }
555:   *maxspeed = a*MaxAbs(1-2*uL[0],1-2*uR[0]);
556:   return(0);
557: }

561: static PetscErrorCode PhysicsRiemann_Traffic_Roe(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
562: {
563:   PetscReal a = ((TrafficCtx*)vctx)->a;
564:   PetscReal speed;

567:   speed = a*(1 - (uL[0] + uR[0]));
568:   flux[0] = 0.5*(TrafficFlux(a,uL[0]) + TrafficFlux(a,uR[0])) - 0.5*PetscAbs(speed)*(uR[0]-uL[0]);
569:   *maxspeed = speed;
570:   return(0);
571: }

575: static PetscErrorCode PhysicsRiemann_Traffic_LxF(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
576: {
577:   TrafficCtx *phys = (TrafficCtx*)vctx;
578:   PetscReal  a     = phys->a;
579:   PetscReal  speed;

582:   speed     = a*(1 - (uL[0] + uR[0]));
583:   flux[0]   = 0.5*(TrafficFlux(a,uL[0]) + TrafficFlux(a,uR[0])) - 0.5*phys->lxf_speed*(uR[0]-uL[0]);
584:   *maxspeed = speed;
585:   return(0);
586: }

590: static PetscErrorCode PhysicsRiemann_Traffic_Rusanov(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
591: {
592:   PetscReal a = ((TrafficCtx*)vctx)->a;
593:   PetscReal speed;

596:   speed     = a*PetscMax(PetscAbs(1-2*uL[0]),PetscAbs(1-2*uR[0]));
597:   flux[0]   = 0.5*(TrafficFlux(a,uL[0]) + TrafficFlux(a,uR[0])) - 0.5*speed*(uR[0]-uL[0]);
598:   *maxspeed = speed;
599:   return(0);
600: }

604: static PetscErrorCode PhysicsCreate_Traffic(FVCtx *ctx)
605: {
606:   PetscErrorCode    ierr;
607:   TrafficCtx        *user;
608:   RiemannFunction   r;
609:   PetscFunctionList rlist      = 0;
610:   char              rname[256] = "exact";

613:   PetscNew(&user);
614:   ctx->physics.sample         = PhysicsSample_Traffic;
615:   ctx->physics.characteristic = PhysicsCharacteristic_Conservative;
616:   ctx->physics.destroy        = PhysicsDestroy_SimpleFree;
617:   ctx->physics.user           = user;
618:   ctx->physics.dof            = 1;

620:   PetscStrallocpy("density",&ctx->physics.fieldname[0]);
621:   user->a = 0.5;
622:   RiemannListAdd(&rlist,"exact",  PhysicsRiemann_Traffic_Exact);
623:   RiemannListAdd(&rlist,"roe",    PhysicsRiemann_Traffic_Roe);
624:   RiemannListAdd(&rlist,"lxf",    PhysicsRiemann_Traffic_LxF);
625:   RiemannListAdd(&rlist,"rusanov",PhysicsRiemann_Traffic_Rusanov);
626:   PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for Traffic","");
627:     PetscOptionsReal("-physics_traffic_a","Flux = a*u*(1-u)","",user->a,&user->a,NULL);
628:     PetscOptionsFList("-physics_traffic_riemann","Riemann solver","",rlist,rname,rname,sizeof(rname),NULL);
629:   PetscOptionsEnd();

631:   RiemannListFind(rlist,rname,&r);
632:   PetscFunctionListDestroy(&rlist);

634:   ctx->physics.riemann = r;

636:   /* *
637:   * Hack to deal with LxF in semi-discrete form
638:   * max speed is 3*a for the basic initial conditions (-1 <= u <= 2)
639:   * */
640:   if (r == PhysicsRiemann_Traffic_LxF) user->lxf_speed = 3*user->a;
641:   return(0);
642: }

644: /* --------------------------------- Linear Acoustics ----------------------------------- */

646: /* Flux: u_t + (A u)_x
647:  * z = sqrt(rho*bulk), c = sqrt(rho/bulk)
648:  * Spectral decomposition: A = R * D * Rinv
649:  * [    cz] = [-z   z] [-c    ] [-1/2z  1/2]
650:  * [c/z   ] = [ 1   1] [     c] [ 1/2z  1/2]
651:  *
652:  * We decompose this into the left-traveling waves Al = R * D^- Rinv
653:  * and the right-traveling waves Ar = R * D^+ * Rinv
654:  * Multiplying out these expressions produces the following two matrices
655:  */

657: typedef struct {
658:   PetscReal c;                  /* speed of sound: c = sqrt(bulk/rho) */
659:   PetscReal z;                  /* impedence: z = sqrt(rho*bulk) */
660: } AcousticsCtx;

662: PETSC_UNUSED PETSC_STATIC_INLINE void AcousticsFlux(AcousticsCtx *ctx,const PetscScalar *u,PetscScalar *f)
663: {
664:   f[0] = ctx->c*ctx->z*u[1];
665:   f[1] = ctx->c/ctx->z*u[0];
666: }

670: static PetscErrorCode PhysicsCharacteristic_Acoustics(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds)
671: {
672:   AcousticsCtx *phys = (AcousticsCtx*)vctx;
673:   PetscReal    z     = phys->z,c = phys->c;

676:   X[0*2+0]  = -z;
677:   X[0*2+1]  = z;
678:   X[1*2+0]  = 1;
679:   X[1*2+1]  = 1;
680:   Xi[0*2+0] = -1./(2*z);
681:   Xi[0*2+1] = 1./2;
682:   Xi[1*2+0] = 1./(2*z);
683:   Xi[1*2+1] = 1./2;
684:   speeds[0] = -c;
685:   speeds[1] = c;
686:   return(0);
687: }

691: static PetscErrorCode PhysicsSample_Acoustics_Initial(AcousticsCtx *phys,PetscInt initial,PetscReal xmin,PetscReal xmax,PetscReal x,PetscReal *u)
692: {
694:   switch (initial) {
695:   case 0:
696:     u[0] = (PetscAbs((x - xmin)/(xmax - xmin) - 0.2) < 0.1) ? 1 : 0.5;
697:     u[1] = (PetscAbs((x - xmin)/(xmax - xmin) - 0.7) < 0.1) ? 1 : -0.5;
698:     break;
699:   case 1:
700:     u[0] = PetscCosReal(3 * 2*PETSC_PI*x/(xmax-xmin));
701:     u[1] = PetscExpReal(-PetscSqr(x - (xmax + xmin)/2) / (2*PetscSqr(0.2*(xmax - xmin)))) - 0.5;
702:     break;
703:   default: SETERRQ(PETSC_COMM_SELF,1,"unknown initial condition");
704:   }
705:   return(0);
706: }

710: static PetscErrorCode PhysicsSample_Acoustics(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u)
711: {
712:   AcousticsCtx   *phys = (AcousticsCtx*)vctx;
713:   PetscReal      c     = phys->c;
714:   PetscReal      x0a,x0b,u0a[2],u0b[2],tmp[2];
715:   PetscReal      X[2][2],Xi[2][2],dummy[2];

719:   switch (bctype) {
720:   case FVBC_OUTFLOW:
721:     x0a = x+c*t;
722:     x0b = x-c*t;
723:     break;
724:   case FVBC_PERIODIC:
725:     x0a = RangeMod(x+c*t,xmin,xmax);
726:     x0b = RangeMod(x-c*t,xmin,xmax);
727:     break;
728:   default: SETERRQ(PETSC_COMM_SELF,1,"unknown BCType");
729:   }
730:   PhysicsSample_Acoustics_Initial(phys,initial,xmin,xmax,x0a,u0a);
731:   PhysicsSample_Acoustics_Initial(phys,initial,xmin,xmax,x0b,u0b);
732:   PhysicsCharacteristic_Acoustics(vctx,2,u,&X[0][0],&Xi[0][0],dummy);
733:   tmp[0] = Xi[0][0]*u0a[0] + Xi[0][1]*u0a[1];
734:   tmp[1] = Xi[1][0]*u0b[0] + Xi[1][1]*u0b[1];
735:   u[0]   = X[0][0]*tmp[0] + X[0][1]*tmp[1];
736:   u[1]   = X[1][0]*tmp[0] + X[1][1]*tmp[1];
737:   return(0);
738: }

742: static PetscErrorCode PhysicsRiemann_Acoustics_Exact(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
743: {
744:   AcousticsCtx *phys = (AcousticsCtx*)vctx;
745:   PetscReal    c     = phys->c,z = phys->z;
746:   PetscReal
747:     Al[2][2] = {{-c/2     , c*z/2  },
748:                 {c/(2*z)  , -c/2   }}, /* Left traveling waves */
749:     Ar[2][2] = {{c/2      , c*z/2  },
750:                 {c/(2*z)  , c/2    }}; /* Right traveling waves */

753:   flux[0]   = Al[0][0]*uR[0] + Al[0][1]*uR[1] + Ar[0][0]*uL[0] + Ar[0][1]*uL[1];
754:   flux[1]   = Al[1][0]*uR[0] + Al[1][1]*uR[1] + Ar[1][0]*uL[0] + Ar[1][1]*uL[1];
755:   *maxspeed = c;
756:   return(0);
757: }

761: static PetscErrorCode PhysicsCreate_Acoustics(FVCtx *ctx)
762: {
763:   PetscErrorCode    ierr;
764:   AcousticsCtx      *user;
765:   PetscFunctionList rlist      = 0,rclist = 0;
766:   char              rname[256] = "exact",rcname[256] = "characteristic";

769:   PetscNew(&user);
770:   ctx->physics.sample         = PhysicsSample_Acoustics;
771:   ctx->physics.destroy        = PhysicsDestroy_SimpleFree;
772:   ctx->physics.user           = user;
773:   ctx->physics.dof            = 2;

775:   PetscStrallocpy("u",&ctx->physics.fieldname[0]);
776:   PetscStrallocpy("v",&ctx->physics.fieldname[1]);

778:   user->c = 1;
779:   user->z = 1;

781:   RiemannListAdd(&rlist,"exact",  PhysicsRiemann_Acoustics_Exact);
782:   ReconstructListAdd(&rclist,"characteristic",PhysicsCharacteristic_Acoustics);
783:   ReconstructListAdd(&rclist,"conservative",PhysicsCharacteristic_Conservative);
784:   PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for linear Acoustics","");
785:   {
786:     PetscOptionsReal("-physics_acoustics_c","c = sqrt(bulk/rho)","",user->c,&user->c,NULL);
787:     PetscOptionsReal("-physics_acoustics_z","z = sqrt(bulk*rho)","",user->z,&user->z,NULL);
788:     PetscOptionsFList("-physics_acoustics_riemann","Riemann solver","",rlist,rname,rname,sizeof(rname),NULL);
789:     PetscOptionsFList("-physics_acoustics_reconstruct","Reconstruction","",rclist,rcname,rcname,sizeof(rcname),NULL);
790:   }
791:   PetscOptionsEnd();
792:   RiemannListFind(rlist,rname,&ctx->physics.riemann);
793:   ReconstructListFind(rclist,rcname,&ctx->physics.characteristic);
794:   PetscFunctionListDestroy(&rlist);
795:   PetscFunctionListDestroy(&rclist);
796:   return(0);
797: }

799: /* --------------------------------- Isothermal Gas Dynamics ----------------------------------- */

801: typedef struct {
802:   PetscReal acoustic_speed;
803: } IsoGasCtx;

805: PETSC_STATIC_INLINE void IsoGasFlux(PetscReal c,const PetscScalar *u,PetscScalar *f)
806: {
807:   f[0] = u[1];
808:   f[1] = PetscSqr(u[1])/u[0] + c*c*u[0];
809: }

813: static PetscErrorCode PhysicsSample_IsoGas(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u)
814: {
816:   if (t > 0) SETERRQ(PETSC_COMM_SELF,1,"Exact solutions not implemented for t > 0");
817:   switch (initial) {
818:     case 0:
819:       u[0] = (x < 0) ? 1 : 0.5;
820:       u[1] = (x < 0) ? 1 : 0.7;
821:       break;
822:     case 1:
823:       u[0] = 1+0.5*PetscSinReal(2*PETSC_PI*x);
824:       u[1] = 1*u[0];
825:       break;
826:     default: SETERRQ(PETSC_COMM_SELF,1,"unknown initial condition");
827:   }
828:   return(0);
829: }

833: static PetscErrorCode PhysicsRiemann_IsoGas_Roe(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
834: {
835:   IsoGasCtx   *phys = (IsoGasCtx*)vctx;
836:   PetscReal   c     = phys->acoustic_speed;
837:   PetscScalar ubar,du[2],a[2],fL[2],fR[2],lam[2],ustar[2],R[2][2];
838:   PetscInt    i;

841:   ubar = (uL[1]/PetscSqrtScalar(uL[0]) + uR[1]/PetscSqrtScalar(uR[0])) / (PetscSqrtScalar(uL[0]) + PetscSqrtScalar(uR[0]));
842:   /* write fluxuations in characteristic basis */
843:   du[0] = uR[0] - uL[0];
844:   du[1] = uR[1] - uL[1];
845:   a[0]  = (1/(2*c)) * ((ubar + c)*du[0] - du[1]);
846:   a[1]  = (1/(2*c)) * ((-ubar + c)*du[0] + du[1]);
847:   /* wave speeds */
848:   lam[0] = ubar - c;
849:   lam[1] = ubar + c;
850:   /* Right eigenvectors */
851:   R[0][0] = 1; R[0][1] = ubar-c;
852:   R[1][0] = 1; R[1][1] = ubar+c;
853:   /* Compute state in star region (between the 1-wave and 2-wave) */
854:   for (i=0; i<2; i++) ustar[i] = uL[i] + a[0]*R[0][i];
855:   if (uL[1]/uL[0] < c && c < ustar[1]/ustar[0]) { /* 1-wave is sonic rarefaction */
856:     PetscScalar ufan[2];
857:     ufan[0] = uL[0]*PetscExpScalar(uL[1]/(uL[0]*c) - 1);
858:     ufan[1] = c*ufan[0];
859:     IsoGasFlux(c,ufan,flux);
860:   } else if (ustar[1]/ustar[0] < -c && -c < uR[1]/uR[0]) { /* 2-wave is sonic rarefaction */
861:     PetscScalar ufan[2];
862:     ufan[0] = uR[0]*PetscExpScalar(-uR[1]/(uR[0]*c) - 1);
863:     ufan[1] = -c*ufan[0];
864:     IsoGasFlux(c,ufan,flux);
865:   } else {                      /* Centered form */
866:     IsoGasFlux(c,uL,fL);
867:     IsoGasFlux(c,uR,fR);
868:     for (i=0; i<2; i++) {
869:       PetscScalar absdu = PetscAbsScalar(lam[0])*a[0]*R[0][i] + PetscAbsScalar(lam[1])*a[1]*R[1][i];
870:       flux[i] = 0.5*(fL[i]+fR[i]) - 0.5*absdu;
871:     }
872:   }
873:   *maxspeed = MaxAbs(lam[0],lam[1]);
874:   return(0);
875: }

879: static PetscErrorCode PhysicsRiemann_IsoGas_Exact(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
880: {
881:   IsoGasCtx                   *phys = (IsoGasCtx*)vctx;
882:   PetscReal                   c     = phys->acoustic_speed;
883:   PetscScalar                 ustar[2];
884:   struct {PetscScalar rho,u;} L = {uL[0],uL[1]/uL[0]},R = {uR[0],uR[1]/uR[0]},star;
885:   PetscInt                    i;
886:   PetscErrorCode              ierr;

889:   if (!(L.rho > 0 && R.rho > 0)) SETERRQ(PETSC_COMM_SELF,1,"Reconstructed density is negative");
890:   {
891:     /* Solve for star state */
892:     PetscScalar res,tmp,rho = 0.5*(L.rho + R.rho); /* initial guess */
893:     for (i=0; i<20; i++) {
894:       PetscScalar fr,fl,dfr,dfl;
895:       fl = (L.rho < rho)
896:         ? (rho-L.rho)/PetscSqrtScalar(L.rho*rho)       /* shock */
897:         : PetscLogScalar(rho) - PetscLogScalar(L.rho); /* rarefaction */
898:       fr = (R.rho < rho)
899:         ? (rho-R.rho)/PetscSqrtScalar(R.rho*rho)       /* shock */
900:         : PetscLogScalar(rho) - PetscLogScalar(R.rho); /* rarefaction */
901:       res = R.u-L.u + c*(fr+fl);
902:       PetscIsInfOrNanScalar(res);
903:       if (PetscAbsScalar(res) < 1e-10) {
904:         star.rho = rho;
905:         star.u   = L.u - c*fl;
906:         goto converged;
907:       }
908:       dfl = (L.rho < rho) ? 1/PetscSqrtScalar(L.rho*rho)*(1 - 0.5*(rho-L.rho)/rho) : 1/rho;
909:       dfr = (R.rho < rho) ? 1/PetscSqrtScalar(R.rho*rho)*(1 - 0.5*(rho-R.rho)/rho) : 1/rho;
910:       tmp = rho - res/(c*(dfr+dfl));
911:       if (tmp <= 0) rho /= 2;   /* Guard against Newton shooting off to a negative density */
912:       else rho = tmp;
913:       if (!((rho > 0) && PetscIsNormalScalar(rho))) SETERRQ1(PETSC_COMM_SELF,1,"non-normal iterate rho=%g",(double)PetscRealPart(rho));
914:     }
915:     SETERRQ1(PETSC_COMM_SELF,1,"Newton iteration for star.rho diverged after %D iterations",i);
916:   }
917: converged:
918:   if (L.u-c < 0 && 0 < star.u-c) { /* 1-wave is sonic rarefaction */
919:     PetscScalar ufan[2];
920:     ufan[0] = L.rho*PetscExpScalar(L.u/c - 1);
921:     ufan[1] = c*ufan[0];
922:     IsoGasFlux(c,ufan,flux);
923:   } else if (star.u+c < 0 && 0 < R.u+c) { /* 2-wave is sonic rarefaction */
924:     PetscScalar ufan[2];
925:     ufan[0] = R.rho*PetscExpScalar(-R.u/c - 1);
926:     ufan[1] = -c*ufan[0];
927:     IsoGasFlux(c,ufan,flux);
928:   } else if ((L.rho >= star.rho && L.u-c >= 0) || (L.rho < star.rho && (star.rho*star.u-L.rho*L.u)/(star.rho-L.rho) > 0)) {
929:     /* 1-wave is supersonic rarefaction, or supersonic shock */
930:     IsoGasFlux(c,uL,flux);
931:   } else if ((star.rho <= R.rho && R.u+c <= 0) || (star.rho > R.rho && (R.rho*R.u-star.rho*star.u)/(R.rho-star.rho) < 0)) {
932:     /* 2-wave is supersonic rarefaction or supersonic shock */
933:     IsoGasFlux(c,uR,flux);
934:   } else {
935:     ustar[0] = star.rho;
936:     ustar[1] = star.rho*star.u;
937:     IsoGasFlux(c,ustar,flux);
938:   }
939:   *maxspeed = MaxAbs(MaxAbs(star.u-c,star.u+c),MaxAbs(L.u-c,R.u+c));
940:   return(0);
941: }

945: static PetscErrorCode PhysicsRiemann_IsoGas_Rusanov(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
946: {
947:   IsoGasCtx                   *phys = (IsoGasCtx*)vctx;
948:   PetscScalar                 c = phys->acoustic_speed,fL[2],fR[2],s;
949:   struct {PetscScalar rho,u;} L = {uL[0],uL[1]/uL[0]},R = {uR[0],uR[1]/uR[0]};

952:   if (!(L.rho > 0 && R.rho > 0)) SETERRQ(PETSC_COMM_SELF,1,"Reconstructed density is negative");
953:   IsoGasFlux(c,uL,fL);
954:   IsoGasFlux(c,uR,fR);
955:   s         = PetscMax(PetscAbs(L.u),PetscAbs(R.u))+c;
956:   flux[0]   = 0.5*(fL[0] + fR[0]) + 0.5*s*(uL[0] - uR[0]);
957:   flux[1]   = 0.5*(fL[1] + fR[1]) + 0.5*s*(uL[1] - uR[1]);
958:   *maxspeed = s;
959:   return(0);
960: }

964: static PetscErrorCode PhysicsCharacteristic_IsoGas(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds)
965: {
966:   IsoGasCtx      *phys = (IsoGasCtx*)vctx;
967:   PetscReal      c     = phys->acoustic_speed;

971:   speeds[0] = u[1]/u[0] - c;
972:   speeds[1] = u[1]/u[0] + c;
973:   X[0*2+0]  = 1;
974:   X[0*2+1]  = speeds[0];
975:   X[1*2+0]  = 1;
976:   X[1*2+1]  = speeds[1];
977:   PetscMemcpy(Xi,X,4*sizeof(X[0]));
978:   PetscKernel_A_gets_inverse_A_2(Xi,0,PETSC_FALSE,NULL);
979:   return(0);
980: }

984: static PetscErrorCode PhysicsCreate_IsoGas(FVCtx *ctx)
985: {
986:   PetscErrorCode    ierr;
987:   IsoGasCtx         *user;
988:   PetscFunctionList rlist = 0,rclist = 0;
989:   char              rname[256] = "exact",rcname[256] = "characteristic";

992:   PetscNew(&user);
993:   ctx->physics.sample         = PhysicsSample_IsoGas;
994:   ctx->physics.destroy        = PhysicsDestroy_SimpleFree;
995:   ctx->physics.user           = user;
996:   ctx->physics.dof            = 2;

998:   PetscStrallocpy("density",&ctx->physics.fieldname[0]);
999:   PetscStrallocpy("momentum",&ctx->physics.fieldname[1]);

1001:   user->acoustic_speed = 1;

1003:   RiemannListAdd(&rlist,"exact",  PhysicsRiemann_IsoGas_Exact);
1004:   RiemannListAdd(&rlist,"roe",    PhysicsRiemann_IsoGas_Roe);
1005:   RiemannListAdd(&rlist,"rusanov",PhysicsRiemann_IsoGas_Rusanov);
1006:   ReconstructListAdd(&rclist,"characteristic",PhysicsCharacteristic_IsoGas);
1007:   ReconstructListAdd(&rclist,"conservative",PhysicsCharacteristic_Conservative);
1008:   PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for IsoGas","");
1009:     PetscOptionsReal("-physics_isogas_acoustic_speed","Acoustic speed","",user->acoustic_speed,&user->acoustic_speed,NULL);
1010:     PetscOptionsFList("-physics_isogas_riemann","Riemann solver","",rlist,rname,rname,sizeof(rname),NULL);
1011:     PetscOptionsFList("-physics_isogas_reconstruct","Reconstruction","",rclist,rcname,rcname,sizeof(rcname),NULL);
1012:   PetscOptionsEnd();
1013:   RiemannListFind(rlist,rname,&ctx->physics.riemann);
1014:   ReconstructListFind(rclist,rcname,&ctx->physics.characteristic);
1015:   PetscFunctionListDestroy(&rlist);
1016:   PetscFunctionListDestroy(&rclist);
1017:   return(0);
1018: }

1020: /* --------------------------------- Shallow Water ----------------------------------- */
1021: typedef struct {
1022:   PetscReal gravity;
1023: } ShallowCtx;

1025: PETSC_STATIC_INLINE void ShallowFlux(ShallowCtx *phys,const PetscScalar *u,PetscScalar *f)
1026: {
1027:   f[0] = u[1];
1028:   f[1] = PetscSqr(u[1])/u[0] + 0.5*phys->gravity*PetscSqr(u[0]);
1029: }

1033: static PetscErrorCode PhysicsRiemann_Shallow_Exact(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
1034: {
1035:   ShallowCtx                *phys = (ShallowCtx*)vctx;
1036:   PetscScalar               g    = phys->gravity,ustar[2],cL,cR,c,cstar;
1037:   struct {PetscScalar h,u;} L = {uL[0],uL[1]/uL[0]},R = {uR[0],uR[1]/uR[0]},star;
1038:   PetscInt                  i;
1039:   PetscErrorCode            ierr;

1042:   if (!(L.h > 0 && R.h > 0)) SETERRQ(PETSC_COMM_SELF,1,"Reconstructed thickness is negative");
1043:   cL = PetscSqrtScalar(g*L.h);
1044:   cR = PetscSqrtScalar(g*R.h);
1045:   c  = PetscMax(cL,cR);
1046:   {
1047:     /* Solve for star state */
1048:     const PetscInt maxits = 50;
1049:     PetscScalar tmp,res,res0=0,h0,h = 0.5*(L.h + R.h); /* initial guess */
1050:     h0 = h;
1051:     for (i=0; i<maxits; i++) {
1052:       PetscScalar fr,fl,dfr,dfl;
1053:       fl = (L.h < h)
1054:         ? PetscSqrtScalar(0.5*g*(h*h - L.h*L.h)*(1/L.h - 1/h)) /* shock */
1055:         : 2*PetscSqrtScalar(g*h) - 2*PetscSqrtScalar(g*L.h);   /* rarefaction */
1056:       fr = (R.h < h)
1057:         ? PetscSqrtScalar(0.5*g*(h*h - R.h*R.h)*(1/R.h - 1/h)) /* shock */
1058:         : 2*PetscSqrtScalar(g*h) - 2*PetscSqrtScalar(g*R.h);   /* rarefaction */
1059:       res = R.u - L.u + fr + fl;
1060:       PetscIsInfOrNanScalar(res);
1061:       if (PetscAbsScalar(res) < 1e-8 || (i > 0 && PetscAbsScalar(h-h0) < 1e-8)) {
1062:         star.h = h;
1063:         star.u = L.u - fl;
1064:         goto converged;
1065:       } else if (i > 0 && PetscAbsScalar(res) >= PetscAbsScalar(res0)) {        /* Line search */
1066:         h = 0.8*h0 + 0.2*h;
1067:         continue;
1068:       }
1069:       /* Accept the last step and take another */
1070:       res0 = res;
1071:       h0 = h;
1072:       dfl = (L.h < h) ? 0.5/fl*0.5*g*(-L.h*L.h/(h*h) - 1 + 2*h/L.h) : PetscSqrtScalar(g/h);
1073:       dfr = (R.h < h) ? 0.5/fr*0.5*g*(-R.h*R.h/(h*h) - 1 + 2*h/R.h) : PetscSqrtScalar(g/h);
1074:       tmp = h - res/(dfr+dfl);
1075:       if (tmp <= 0) h /= 2;   /* Guard against Newton shooting off to a negative thickness */
1076:       else h = tmp;
1077:       if (!((h > 0) && PetscIsNormalScalar(h))) SETERRQ1(PETSC_COMM_SELF,1,"non-normal iterate h=%g",(double)h);
1078:     }
1079:     SETERRQ1(PETSC_COMM_SELF,1,"Newton iteration for star.h diverged after %D iterations",i);
1080:   }
1081: converged:
1082:   cstar = PetscSqrtScalar(g*star.h);
1083:   if (L.u-cL < 0 && 0 < star.u-cstar) { /* 1-wave is sonic rarefaction */
1084:     PetscScalar ufan[2];
1085:     ufan[0] = 1/g*PetscSqr(L.u/3 + 2./3*cL);
1086:     ufan[1] = PetscSqrtScalar(g*ufan[0])*ufan[0];
1087:     ShallowFlux(phys,ufan,flux);
1088:   } else if (star.u+cstar < 0 && 0 < R.u+cR) { /* 2-wave is sonic rarefaction */
1089:     PetscScalar ufan[2];
1090:     ufan[0] = 1/g*PetscSqr(R.u/3 - 2./3*cR);
1091:     ufan[1] = -PetscSqrtScalar(g*ufan[0])*ufan[0];
1092:     ShallowFlux(phys,ufan,flux);
1093:   } else if ((L.h >= star.h && L.u-c >= 0) || (L.h<star.h && (star.h*star.u-L.h*L.u)/(star.h-L.h) > 0)) {
1094:     /* 1-wave is right-travelling shock (supersonic) */
1095:     ShallowFlux(phys,uL,flux);
1096:   } else if ((star.h <= R.h && R.u+c <= 0) || (star.h>R.h && (R.h*R.u-star.h*star.h)/(R.h-star.h) < 0)) {
1097:     /* 2-wave is left-travelling shock (supersonic) */
1098:     ShallowFlux(phys,uR,flux);
1099:   } else {
1100:     ustar[0] = star.h;
1101:     ustar[1] = star.h*star.u;
1102:     ShallowFlux(phys,ustar,flux);
1103:   }
1104:   *maxspeed = MaxAbs(MaxAbs(star.u-cstar,star.u+cstar),MaxAbs(L.u-cL,R.u+cR));
1105:   return(0);
1106: }

1110: static PetscErrorCode PhysicsRiemann_Shallow_Rusanov(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
1111: {
1112:   ShallowCtx                *phys = (ShallowCtx*)vctx;
1113:   PetscScalar               g = phys->gravity,fL[2],fR[2],s;
1114:   struct {PetscScalar h,u;} L = {uL[0],uL[1]/uL[0]},R = {uR[0],uR[1]/uR[0]};

1117:   if (!(L.h > 0 && R.h > 0)) SETERRQ(PETSC_COMM_SELF,1,"Reconstructed thickness is negative");
1118:   ShallowFlux(phys,uL,fL);
1119:   ShallowFlux(phys,uR,fR);
1120:   s         = PetscMax(PetscAbs(L.u)+PetscSqrtScalar(g*L.h),PetscAbs(R.u)+PetscSqrtScalar(g*R.h));
1121:   flux[0]   = 0.5*(fL[0] + fR[0]) + 0.5*s*(uL[0] - uR[0]);
1122:   flux[1]   = 0.5*(fL[1] + fR[1]) + 0.5*s*(uL[1] - uR[1]);
1123:   *maxspeed = s;
1124:   return(0);
1125: }

1129: static PetscErrorCode PhysicsCharacteristic_Shallow(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds)
1130: {
1131:   ShallowCtx     *phys = (ShallowCtx*)vctx;
1132:   PetscReal      c;

1136:   c         = PetscSqrtScalar(u[0]*phys->gravity);
1137:   speeds[0] = u[1]/u[0] - c;
1138:   speeds[1] = u[1]/u[0] + c;
1139:   X[0*2+0]  = 1;
1140:   X[0*2+1]  = speeds[0];
1141:   X[1*2+0]  = 1;
1142:   X[1*2+1]  = speeds[1];
1143:   PetscMemcpy(Xi,X,4*sizeof(X[0]));
1144:   PetscKernel_A_gets_inverse_A_2(Xi,0,PETSC_FALSE,NULL);
1145:   return(0);
1146: }

1150: static PetscErrorCode PhysicsCreate_Shallow(FVCtx *ctx)
1151: {
1152:   PetscErrorCode    ierr;
1153:   ShallowCtx        *user;
1154:   PetscFunctionList rlist = 0,rclist = 0;
1155:   char              rname[256] = "exact",rcname[256] = "characteristic";

1158:   PetscNew(&user);
1159:   /* Shallow water and Isothermal Gas dynamics are similar so we reuse initial conditions for now */
1160:   ctx->physics.sample         = PhysicsSample_IsoGas;
1161:   ctx->physics.destroy        = PhysicsDestroy_SimpleFree;
1162:   ctx->physics.user           = user;
1163:   ctx->physics.dof            = 2;

1165:   PetscStrallocpy("density",&ctx->physics.fieldname[0]);
1166:   PetscStrallocpy("momentum",&ctx->physics.fieldname[1]);

1168:   user->gravity = 1;

1170:   RiemannListAdd(&rlist,"exact",  PhysicsRiemann_Shallow_Exact);
1171:   RiemannListAdd(&rlist,"rusanov",PhysicsRiemann_Shallow_Rusanov);
1172:   ReconstructListAdd(&rclist,"characteristic",PhysicsCharacteristic_Shallow);
1173:   ReconstructListAdd(&rclist,"conservative",PhysicsCharacteristic_Conservative);
1174:   PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for Shallow","");
1175:     PetscOptionsReal("-physics_shallow_gravity","Gravity","",user->gravity,&user->gravity,NULL);
1176:     PetscOptionsFList("-physics_shallow_riemann","Riemann solver","",rlist,rname,rname,sizeof(rname),NULL);
1177:     PetscOptionsFList("-physics_shallow_reconstruct","Reconstruction","",rclist,rcname,rcname,sizeof(rcname),NULL);
1178:   PetscOptionsEnd();
1179:   RiemannListFind(rlist,rname,&ctx->physics.riemann);
1180:   ReconstructListFind(rclist,rcname,&ctx->physics.characteristic);
1181:   PetscFunctionListDestroy(&rlist);
1182:   PetscFunctionListDestroy(&rclist);
1183:   return(0);
1184: }

1186: /* --------------------------------- Finite Volume Solver ----------------------------------- */

1190: static PetscErrorCode FVRHSFunction(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
1191: {
1192:   FVCtx             *ctx = (FVCtx*)vctx;
1193:   PetscErrorCode    ierr;
1194:   PetscInt          i,j,k,Mx,dof,xs,xm;
1195:   PetscReal         hx,cfl_idt = 0;
1196:   PetscScalar       *x,*f,*slope;
1197:   Vec               Xloc;
1198:   DM                da;

1201:   TSGetDM(ts,&da);
1202:   DMGetLocalVector(da,&Xloc);
1203:   DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);
1204:   hx   = (ctx->xmax - ctx->xmin)/Mx;
1205:   DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);
1206:   DMGlobalToLocalEnd  (da,X,INSERT_VALUES,Xloc);

1208:   VecZeroEntries(F);

1210:   DMDAVecGetArray(da,Xloc,&x);
1211:   DMDAVecGetArray(da,F,&f);
1212:   DMDAGetArray(da,PETSC_TRUE,&slope);

1214:   DMDAGetCorners(da,&xs,0,0,&xm,0,0);

1216:   if (ctx->bctype == FVBC_OUTFLOW) {
1217:     for (i=xs-2; i<0; i++) {
1218:       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
1219:     }
1220:     for (i=Mx; i<xs+xm+2; i++) {
1221:       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
1222:     }
1223:   }
1224:   for (i=xs-1; i<xs+xm+1; i++) {
1225:     struct _LimitInfo info;
1226:     PetscScalar       *cjmpL,*cjmpR;
1227:     /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
1228:     (*ctx->physics.characteristic)(ctx->physics.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);
1229:     /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
1230:     PetscMemzero(ctx->cjmpLR,2*dof*sizeof(ctx->cjmpLR[0]));
1231:     cjmpL = &ctx->cjmpLR[0];
1232:     cjmpR = &ctx->cjmpLR[dof];
1233:     for (j=0; j<dof; j++) {
1234:       PetscScalar jmpL,jmpR;
1235:       jmpL = x[(i+0)*dof+j] - x[(i-1)*dof+j];
1236:       jmpR = x[(i+1)*dof+j] - x[(i+0)*dof+j];
1237:       for (k=0; k<dof; k++) {
1238:         cjmpL[k] += ctx->Rinv[k+j*dof] * jmpL;
1239:         cjmpR[k] += ctx->Rinv[k+j*dof] * jmpR;
1240:       }
1241:     }
1242:     /* Apply limiter to the left and right characteristic jumps */
1243:     info.m  = dof;
1244:     info.hx = hx;
1245:     (*ctx->limit)(&info,cjmpL,cjmpR,ctx->cslope);
1246:     for (j=0; j<dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */
1247:     for (j=0; j<dof; j++) {
1248:       PetscScalar tmp = 0;
1249:       for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof] * ctx->cslope[k];
1250:       slope[i*dof+j] = tmp;
1251:     }
1252:   }

1254:   for (i=xs; i<xs+xm+1; i++) {
1255:     PetscReal   maxspeed;
1256:     PetscScalar *uL,*uR;
1257:     uL = &ctx->uLR[0];
1258:     uR = &ctx->uLR[dof];
1259:     for (j=0; j<dof; j++) {
1260:       uL[j] = x[(i-1)*dof+j] + slope[(i-1)*dof+j]*hx/2;
1261:       uR[j] = x[(i-0)*dof+j] - slope[(i-0)*dof+j]*hx/2;
1262:     }
1263:     (*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed);
1264:     cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */

1266:     if (i > xs) {
1267:       for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hx;
1268:     }
1269:     if (i < xs+xm) {
1270:       for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hx;
1271:     }
1272:   }

1274:   DMDAVecRestoreArray(da,Xloc,&x);
1275:   DMDAVecRestoreArray(da,F,&f);
1276:   DMDARestoreArray(da,PETSC_TRUE,&slope);
1277:   DMRestoreLocalVector(da,&Xloc);

1279:   MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)da));
1280:   if (0) {
1281:     /* We need to a way to inform the TS of a CFL constraint, this is a debugging fragment */
1282:     PetscReal dt,tnow;
1283:     TSGetTimeStep(ts,&dt);
1284:     TSGetTime(ts,&tnow);
1285:     if (dt > 0.5/ctx->cfl_idt) {
1286:       if (1) {
1287:         PetscPrintf(ctx->comm,"Stability constraint exceeded at t=%g, dt %g > %g\n",(double)tnow,(double)dt,(double)(0.5/ctx->cfl_idt));
1288:       } else SETERRQ2(PETSC_COMM_SELF,1,"Stability constraint exceeded, %g > %g",(double)dt,(double)(ctx->cfl/ctx->cfl_idt));
1289:     }
1290:   }
1291:   return(0);
1292: }

1296: static PetscErrorCode SmallMatMultADB(PetscScalar *C,PetscInt bs,const PetscScalar *A,const PetscReal *D,const PetscScalar *B)
1297: {
1298:   PetscInt i,j,k;

1301:   for (i=0; i<bs; i++) {
1302:     for (j=0; j<bs; j++) {
1303:       PetscScalar tmp = 0;
1304:       for (k=0; k<bs; k++) tmp += A[i*bs+k] * D[k] * B[k*bs+j];
1305:       C[i*bs+j] = tmp;
1306:     }
1307:   }
1308:   return(0);
1309: }


1314: static PetscErrorCode FVIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat A,Mat B,void *vctx)
1315: {
1316:   FVCtx             *ctx = (FVCtx*)vctx;
1317:   PetscErrorCode    ierr;
1318:   PetscInt          i,j,dof = ctx->physics.dof;
1319:   PetscScalar       *J;
1320:   const PetscScalar *x;
1321:   PetscReal         hx;
1322:   DM                da;
1323:   DMDALocalInfo     dainfo;

1326:   TSGetDM(ts,&da);
1327:   DMDAVecGetArrayRead(da,X,(void*)&x);
1328:   DMDAGetLocalInfo(da,&dainfo);
1329:   hx   = (ctx->xmax - ctx->xmin)/dainfo.mx;
1330:   PetscMalloc1(dof*dof,&J);
1331:   for (i=dainfo.xs; i<dainfo.xs+dainfo.xm; i++) {
1332:     (*ctx->physics.characteristic)(ctx->physics.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);
1333:     for (j=0; j<dof; j++) ctx->speeds[j] = PetscAbs(ctx->speeds[j]);
1334:     SmallMatMultADB(J,dof,ctx->R,ctx->speeds,ctx->Rinv);
1335:     for (j=0; j<dof*dof; j++) J[j] = J[j]/hx + shift*(j/dof == j%dof);
1336:     MatSetValuesBlocked(B,1,&i,1,&i,J,INSERT_VALUES);
1337:   }
1338:   PetscFree(J);
1339:   DMDAVecRestoreArrayRead(da,X,(void*)&x);

1341:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
1342:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
1343:   if (A != B) {
1344:     MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
1345:     MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
1346:   }
1347:   return(0);
1348: }

1352: static PetscErrorCode FVSample(FVCtx *ctx,DM da,PetscReal time,Vec U)
1353: {
1355:   PetscScalar    *u,*uj;
1356:   PetscInt       i,j,k,dof,xs,xm,Mx;

1359:   if (!ctx->physics.sample) SETERRQ(PETSC_COMM_SELF,1,"Physics has not provided a sampling function");
1360:   DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);
1361:   DMDAGetCorners(da,&xs,0,0,&xm,0,0);
1362:   DMDAVecGetArray(da,U,&u);
1363:   PetscMalloc1(dof,&uj);
1364:   for (i=xs; i<xs+xm; i++) {
1365:     const PetscReal h = (ctx->xmax-ctx->xmin)/Mx,xi = ctx->xmin+h/2+i*h;
1366:     const PetscInt  N = 200;
1367:     /* Integrate over cell i using trapezoid rule with N points. */
1368:     for (k=0; k<dof; k++) u[i*dof+k] = 0;
1369:     for (j=0; j<N+1; j++) {
1370:       PetscScalar xj = xi+h*(j-N/2)/(PetscReal)N;
1371:       (*ctx->physics.sample)(ctx->physics.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);
1372:       for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N;
1373:     }
1374:   }
1375:   DMDAVecRestoreArray(da,U,&u);
1376:   PetscFree(uj);
1377:   return(0);
1378: }

1382: static PetscErrorCode SolutionStatsView(DM da,Vec X,PetscViewer viewer)
1383: {
1384:   PetscErrorCode    ierr;
1385:   PetscReal         xmin,xmax;
1386:   PetscScalar       sum,tvsum,tvgsum;
1387:   const PetscScalar *x;
1388:   PetscInt          imin,imax,Mx,i,j,xs,xm,dof;
1389:   Vec               Xloc;
1390:   PetscBool         iascii;

1393:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1394:   if (iascii) {
1395:     /* PETSc lacks a function to compute total variation norm (difficult in multiple dimensions), we do it here */
1396:     DMGetLocalVector(da,&Xloc);
1397:     DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);
1398:     DMGlobalToLocalEnd  (da,X,INSERT_VALUES,Xloc);
1399:     DMDAVecGetArrayRead(da,Xloc,(void*)&x);
1400:     DMDAGetCorners(da,&xs,0,0,&xm,0,0);
1401:     DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);
1402:     tvsum = 0;
1403:     for (i=xs; i<xs+xm; i++) {
1404:       for (j=0; j<dof; j++) tvsum += PetscAbsScalar(x[i*dof+j] - x[(i-1)*dof+j]);
1405:     }
1406:     MPI_Allreduce(&tvsum,&tvgsum,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)da));
1407:     DMDAVecRestoreArrayRead(da,Xloc,(void*)&x);
1408:     DMRestoreLocalVector(da,&Xloc);

1410:     VecMin(X,&imin,&xmin);
1411:     VecMax(X,&imax,&xmax);
1412:     VecSum(X,&sum);
1413:     PetscViewerASCIIPrintf(viewer,"Solution range [%8.5f,%8.5f] with extrema at %D and %D, mean %8.5f, ||x||_TV %8.5f\n",(double)xmin,(double)xmax,imin,imax,(double)(sum/Mx),(double)(tvgsum/Mx));
1414:   } else SETERRQ(PETSC_COMM_SELF,1,"Viewer type not supported");
1415:   return(0);
1416: }

1420: static PetscErrorCode SolutionErrorNorms(FVCtx *ctx,DM da,PetscReal t,Vec X,PetscReal *nrm1,PetscReal *nrmsup)
1421: {
1423:   Vec            Y;
1424:   PetscInt       Mx;

1427:   VecGetSize(X,&Mx);
1428:   VecDuplicate(X,&Y);
1429:   FVSample(ctx,da,t,Y);
1430:   VecAYPX(Y,-1,X);
1431:   VecNorm(Y,NORM_1,nrm1);
1432:   VecNorm(Y,NORM_INFINITY,nrmsup);
1433:   *nrm1 /= Mx;
1434:   VecDestroy(&Y);
1435:   return(0);
1436: }

1440: int main(int argc,char *argv[])
1441: {
1442:   char              lname[256] = "mc",physname[256] = "advect",final_fname[256] = "solution.m";
1443:   PetscFunctionList limiters   = 0,physics = 0;
1444:   MPI_Comm          comm;
1445:   TS                ts;
1446:   DM                da;
1447:   Vec               X,X0,R;
1448:   Mat               B;
1449:   FVCtx             ctx;
1450:   PetscInt          i,dof,xs,xm,Mx,draw = 0;
1451:   PetscBool         view_final = PETSC_FALSE;
1452:   PetscReal         ptime;
1453:   PetscErrorCode    ierr;

1455:   PetscInitialize(&argc,&argv,0,help);
1456:   comm = PETSC_COMM_WORLD;
1457:   PetscMemzero(&ctx,sizeof(ctx));

1459:   /* Register limiters to be available on the command line */
1460:   PetscFunctionListAdd(&limiters,"upwind"              ,Limit_Upwind);
1461:   PetscFunctionListAdd(&limiters,"lax-wendroff"        ,Limit_LaxWendroff);
1462:   PetscFunctionListAdd(&limiters,"beam-warming"        ,Limit_BeamWarming);
1463:   PetscFunctionListAdd(&limiters,"fromm"               ,Limit_Fromm);
1464:   PetscFunctionListAdd(&limiters,"minmod"              ,Limit_Minmod);
1465:   PetscFunctionListAdd(&limiters,"superbee"            ,Limit_Superbee);
1466:   PetscFunctionListAdd(&limiters,"mc"                  ,Limit_MC);
1467:   PetscFunctionListAdd(&limiters,"vanleer"             ,Limit_VanLeer);
1468:   PetscFunctionListAdd(&limiters,"vanalbada"           ,Limit_VanAlbada);
1469:   PetscFunctionListAdd(&limiters,"vanalbadatvd"        ,Limit_VanAlbadaTVD);
1470:   PetscFunctionListAdd(&limiters,"koren"               ,Limit_Koren);
1471:   PetscFunctionListAdd(&limiters,"korensym"            ,Limit_KorenSym);
1472:   PetscFunctionListAdd(&limiters,"koren3"              ,Limit_Koren3);
1473:   PetscFunctionListAdd(&limiters,"cada-torrilhon2"     ,Limit_CadaTorrilhon2);
1474:   PetscFunctionListAdd(&limiters,"cada-torrilhon3-r0p1",Limit_CadaTorrilhon3R0p1);
1475:   PetscFunctionListAdd(&limiters,"cada-torrilhon3-r1"  ,Limit_CadaTorrilhon3R1);
1476:   PetscFunctionListAdd(&limiters,"cada-torrilhon3-r10" ,Limit_CadaTorrilhon3R10);
1477:   PetscFunctionListAdd(&limiters,"cada-torrilhon3-r100",Limit_CadaTorrilhon3R100);

1479:   /* Register physical models to be available on the command line */
1480:   PetscFunctionListAdd(&physics,"advect"          ,PhysicsCreate_Advect);
1481:   PetscFunctionListAdd(&physics,"burgers"         ,PhysicsCreate_Burgers);
1482:   PetscFunctionListAdd(&physics,"traffic"         ,PhysicsCreate_Traffic);
1483:   PetscFunctionListAdd(&physics,"acoustics"       ,PhysicsCreate_Acoustics);
1484:   PetscFunctionListAdd(&physics,"isogas"          ,PhysicsCreate_IsoGas);
1485:   PetscFunctionListAdd(&physics,"shallow"         ,PhysicsCreate_Shallow);

1487:   ctx.comm = comm;
1488:   ctx.cfl  = 0.9; ctx.bctype = FVBC_PERIODIC;
1489:   ctx.xmin = -1; ctx.xmax = 1;
1490:   PetscOptionsBegin(comm,NULL,"Finite Volume solver options","");
1491:     PetscOptionsReal("-xmin","X min","",ctx.xmin,&ctx.xmin,NULL);
1492:     PetscOptionsReal("-xmax","X max","",ctx.xmax,&ctx.xmax,NULL);
1493:     PetscOptionsFList("-limit","Name of flux limiter to use","",limiters,lname,lname,sizeof(lname),NULL);
1494:     PetscOptionsFList("-physics","Name of physics (Riemann solver and characteristics) to use","",physics,physname,physname,sizeof(physname),NULL);
1495:     PetscOptionsInt("-draw","Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)","",draw,&draw,NULL);
1496:     PetscOptionsString("-view_final","Write final solution in ASCII MATLAB format to given file name","",final_fname,final_fname,sizeof(final_fname),&view_final);
1497:     PetscOptionsInt("-initial","Initial condition (depends on the physics)","",ctx.initial,&ctx.initial,NULL);
1498:     PetscOptionsBool("-exact","Compare errors with exact solution","",ctx.exact,&ctx.exact,NULL);
1499:     PetscOptionsReal("-cfl","CFL number to time step at","",ctx.cfl,&ctx.cfl,NULL);
1500:     PetscOptionsEnum("-bc_type","Boundary condition","",FVBCTypes,(PetscEnum)ctx.bctype,(PetscEnum*)&ctx.bctype,NULL);
1501:   PetscOptionsEnd();

1503:   /* Choose the limiter from the list of registered limiters */
1504:   PetscFunctionListFind(limiters,lname,&ctx.limit);
1505:   if (!ctx.limit) SETERRQ1(PETSC_COMM_SELF,1,"Limiter '%s' not found",lname);

1507:   /* Choose the physics from the list of registered models */
1508:   {
1509:     PetscErrorCode (*r)(FVCtx*);
1510:     PetscFunctionListFind(physics,physname,&r);
1511:     if (!r) SETERRQ1(PETSC_COMM_SELF,1,"Physics '%s' not found",physname);
1512:     /* Create the physics, will set the number of fields and their names */
1513:     (*r)(&ctx);
1514:   }

1516:   /* Create a DMDA to manage the parallel grid */
1517:   DMDACreate1d(comm,DM_BOUNDARY_PERIODIC,-50,ctx.physics.dof,2,NULL,&da);
1518:   /* Inform the DMDA of the field names provided by the physics. */
1519:   /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */
1520:   for (i=0; i<ctx.physics.dof; i++) {
1521:     DMDASetFieldName(da,i,ctx.physics.fieldname[i]);
1522:   }
1523:   DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);
1524:   DMDAGetCorners(da,&xs,0,0,&xm,0,0);

1526:   /* Set coordinates of cell centers */
1527:   DMDASetUniformCoordinates(da,ctx.xmin+0.5*(ctx.xmax-ctx.xmin)/Mx,ctx.xmax+0.5*(ctx.xmax-ctx.xmin)/Mx,0,0,0,0);

1529:   /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */
1530:   PetscMalloc4(dof*dof,&ctx.R,dof*dof,&ctx.Rinv,2*dof,&ctx.cjmpLR,1*dof,&ctx.cslope);
1531:   PetscMalloc3(2*dof,&ctx.uLR,dof,&ctx.flux,dof,&ctx.speeds);

1533:   /* Create a vector to store the solution and to save the initial state */
1534:   DMCreateGlobalVector(da,&X);
1535:   VecDuplicate(X,&X0);
1536:   VecDuplicate(X,&R);

1538:   DMCreateMatrix(da,&B);

1540:   /* Create a time-stepping object */
1541:   TSCreate(comm,&ts);
1542:   TSSetDM(ts,da);
1543:   TSSetRHSFunction(ts,R,FVRHSFunction,&ctx);
1544:   TSSetIJacobian(ts,B,B,FVIJacobian,&ctx);
1545:   TSSetType(ts,TSSSP);
1546:   TSSetDuration(ts,1000,10);
1547:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
1548: 
1549:   /* Compute initial conditions and starting time step */
1550:   FVSample(&ctx,da,0,X0);
1551:   FVRHSFunction(ts,0,X0,X,(void*)&ctx); /* Initial function evaluation, only used to determine max speed */
1552:   VecCopy(X0,X);                        /* The function value was not used so we set X=X0 again */
1553:   TSSetInitialTimeStep(ts,0,ctx.cfl/ctx.cfl_idt);
1554:   TSSetFromOptions(ts); /* Take runtime options */
1555:   SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD);
1556:   {
1557:     PetscReal nrm1,nrmsup;
1558:     PetscInt  steps;

1560:     TSSolve(ts,X);
1561:     TSGetSolveTime(ts,&ptime);
1562:     TSGetTimeStepNumber(ts,&steps);

1564:     PetscPrintf(comm,"Final time %8.5f, steps %D\n",(double)ptime,steps);
1565:     if (ctx.exact) {
1566:       SolutionErrorNorms(&ctx,da,ptime,X,&nrm1,&nrmsup);
1567:       PetscPrintf(comm,"Error ||x-x_e||_1 %8.4e  ||x-x_e||_sup %8.4e\n",(double)nrm1,(double)nrmsup);
1568:     }
1569:   }

1571:   SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD);
1572:   if (draw & 0x1) {VecView(X0,PETSC_VIEWER_DRAW_WORLD);}
1573:   if (draw & 0x2) {VecView(X,PETSC_VIEWER_DRAW_WORLD);}
1574:   if (draw & 0x4) {
1575:     Vec Y;
1576:     VecDuplicate(X,&Y);
1577:     FVSample(&ctx,da,ptime,Y);
1578:     VecAYPX(Y,-1,X);
1579:     VecView(Y,PETSC_VIEWER_DRAW_WORLD);
1580:     VecDestroy(&Y);
1581:   }

1583:   if (view_final) {
1584:     PetscViewer viewer;
1585:     PetscViewerASCIIOpen(PETSC_COMM_WORLD,final_fname,&viewer);
1586:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
1587:     VecView(X,viewer);
1588:     PetscViewerPopFormat(viewer);
1589:     PetscViewerDestroy(&viewer);
1590:   }

1592:   /* Clean up */
1593:   (*ctx.physics.destroy)(ctx.physics.user);
1594:   for (i=0; i<ctx.physics.dof; i++) {PetscFree(ctx.physics.fieldname[i]);}
1595:   PetscFree4(ctx.R,ctx.Rinv,ctx.cjmpLR,ctx.cslope);
1596:   PetscFree3(ctx.uLR,ctx.flux,ctx.speeds);
1597:   VecDestroy(&X);
1598:   VecDestroy(&X0);
1599:   VecDestroy(&R);
1600:   MatDestroy(&B);
1601:   DMDestroy(&da);
1602:   TSDestroy(&ts);
1603:   PetscFunctionListDestroy(&limiters);
1604:   PetscFunctionListDestroy(&physics);
1605:   PetscFinalize();
1606:   return 0;
1607: }