Actual source code: ex11.c

  1: static char help[] = "Second Order TVD Finite Volume Example.\n";

We use a second order TVD finite volume method to evolve a system of PDEs. Our simple upwinded residual evaluation loops
over all mesh faces and uses a Riemann solver to produce the flux given the face geometry and cell values,
\begin{equation}
f_i = \mathrm{riemann}(\mathrm{phys}, p_\mathrm{centroid}, \hat n, x^L, x^R)
\end{equation}
and then update the cell values given the cell volume.
\begin{eqnarray}
f^L_i &-=& \frac{f_i}{vol^L} \\
f^R_i &+=& \frac{f_i}{vol^R}
\end{eqnarray}

As an example, we can consider the shallow water wave equation,
\begin{eqnarray}
h_t + \nabla\cdot \left( uh \right) &=& 0 \\
(uh)_t + \nabla\cdot \left( u\otimes uh + \frac{g h^2}{2} I \right) &=& 0
\end{eqnarray}
where $h$ is wave height, $u$ is wave velocity, and $g$ is the acceleration due to gravity.

A representative Riemann solver for the shallow water equations is given in the PhysicsRiemann_SW() function,
\begin{eqnarray}
f^{L,R}_h &=& uh^{L,R} \cdot \hat n \\
f^{L,R}_{uh} &=& \frac{f^{L,R}_h}{h^{L,R}} uh^{L,R} + g (h^{L,R})^2 \hat n \\
c^{L,R} &=& \sqrt{g h^{L,R}} \\
s &=& \max\left( \left|\frac{uh^L \cdot \hat n}{h^L}\right| + c^L, \left|\frac{uh^R \cdot \hat n}{h^R}\right| + c^R \right) \\
f_i &=& \frac{A_\mathrm{face}}{2} \left( f^L_i + f^R_i + s \left( x^L_i - x^R_i \right) \right)
\end{eqnarray}
where $c$ is the local gravity wave speed and $f_i$ is a Rusanov flux.

The more sophisticated residual evaluation in RHSFunctionLocal_LS() uses a least-squares fit to a quadratic polynomial
over a neighborhood of the given element.

The mesh is read in from an ExodusII file, usually generated by Cubit.
 37: #include <petscdmplex.h>
 38: #include <petscdmforest.h>
 39: #include <petscds.h>
 40: #include <petscts.h>

 42: #define DIM 2                   /* Geometric dimension */
 43: #define ALEN(a) (sizeof(a)/sizeof((a)[0]))

 45: static PetscFunctionList PhysicsList, PhysicsRiemannList_SW;

 47: /* Represents continuum physical equations. */
 48: typedef struct _n_Physics *Physics;

 50: /* Physical model includes boundary conditions, initial conditions, and functionals of interest. It is
 51:  * discretization-independent, but its members depend on the scenario being solved. */
 52: typedef struct _n_Model *Model;

 54: /* 'User' implements a discretization of a continuous model. */
 55: typedef struct _n_User *User;
 56: typedef PetscErrorCode (*SolutionFunction)(Model,PetscReal,const PetscReal*,PetscScalar*,void*);
 57: typedef PetscErrorCode (*SetUpBCFunction)(DM,PetscDS,Physics);
 58: typedef PetscErrorCode (*FunctionalFunction)(Model,PetscReal,const PetscReal*,const PetscScalar*,PetscReal*,void*);
 59: typedef PetscErrorCode (*SetupFields)(Physics,PetscSection);
 60: static PetscErrorCode ModelSolutionSetDefault(Model,SolutionFunction,void*);
 61: static PetscErrorCode ModelFunctionalRegister(Model,const char*,PetscInt*,FunctionalFunction,void*);
 62: static PetscErrorCode OutputVTK(DM,const char*,PetscViewer*);

 64: struct FieldDescription {
 65:   const char *name;
 66:   PetscInt dof;
 67: };

 69: typedef struct _n_FunctionalLink *FunctionalLink;
 70: struct _n_FunctionalLink {
 71:   char               *name;
 72:   FunctionalFunction func;
 73:   void               *ctx;
 74:   PetscInt           offset;
 75:   FunctionalLink     next;
 76: };

 78: struct _n_Physics {
 79:   PetscRiemannFunc riemann;
 80:   PetscInt         dof;          /* number of degrees of freedom per cell */
 81:   PetscReal        maxspeed;     /* kludge to pick initial time step, need to add monitoring and step control */
 82:   void             *data;
 83:   PetscInt         nfields;
 84:   const struct FieldDescription *field_desc;
 85: };

 87: struct _n_Model {
 88:   MPI_Comm         comm;        /* Does not do collective communicaton, but some error conditions can be collective */
 89:   Physics          physics;
 90:   FunctionalLink   functionalRegistry;
 91:   PetscInt         maxComputed;
 92:   PetscInt         numMonitored;
 93:   FunctionalLink   *functionalMonitored;
 94:   PetscInt         numCall;
 95:   FunctionalLink   *functionalCall;
 96:   SolutionFunction solution;
 97:   SetUpBCFunction  setupbc;
 98:   void             *solutionctx;
 99:   PetscReal        maxspeed;    /* estimate of global maximum speed (for CFL calculation) */
100:   PetscReal        bounds[2*DIM];
101:   PetscErrorCode   (*errorIndicator)(PetscInt, PetscReal, PetscInt, const PetscScalar[], const PetscScalar[], PetscReal *, void *);
102:   void             *errorCtx;
103: };

105: struct _n_User {
106:   PetscInt vtkInterval;   /* For monitor */
107:   char outputBasename[PETSC_MAX_PATH_LEN]; /* Basename for output files */
108:   PetscInt monitorStepOffset;
109:   Model    model;
110:   PetscBool vtkmon;
111: };

113: static inline PetscReal DotDIMReal(const PetscReal *x,const PetscReal *y)
114: {
115:   PetscInt  i;
116:   PetscReal prod=0.0;

118:   for (i=0; i<DIM; i++) prod += x[i]*y[i];
119:   return prod;
120: }
121: static inline PetscReal NormDIM(const PetscReal *x) { return PetscSqrtReal(PetscAbsReal(DotDIMReal(x,x))); }

123: static inline PetscReal Dot2Real(const PetscReal *x,const PetscReal *y) { return x[0]*y[0] + x[1]*y[1];}
124: static inline PetscReal Norm2Real(const PetscReal *x) { return PetscSqrtReal(PetscAbsReal(Dot2Real(x,x)));}
125: static inline void Normalize2Real(PetscReal *x) { PetscReal a = 1./Norm2Real(x); x[0] *= a; x[1] *= a; }
126: static inline void Waxpy2Real(PetscReal a,const PetscReal *x,const PetscReal *y,PetscReal *w) { w[0] = a*x[0] + y[0]; w[1] = a*x[1] + y[1]; }
127: static inline void Scale2Real(PetscReal a,const PetscReal *x,PetscReal *y) { y[0] = a*x[0]; y[1] = a*x[1]; }

129: /******************* Advect ********************/
130: typedef enum {ADVECT_SOL_TILTED,ADVECT_SOL_BUMP,ADVECT_SOL_BUMP_CAVITY} AdvectSolType;
131: static const char *const AdvectSolTypes[] = {"TILTED","BUMP","BUMP_CAVITY","AdvectSolType","ADVECT_SOL_",0};
132: typedef enum {ADVECT_SOL_BUMP_CONE,ADVECT_SOL_BUMP_COS} AdvectSolBumpType;
133: static const char *const AdvectSolBumpTypes[] = {"CONE","COS","AdvectSolBumpType","ADVECT_SOL_BUMP_",0};

135: typedef struct {
136:   PetscReal wind[DIM];
137: } Physics_Advect_Tilted;
138: typedef struct {
139:   PetscReal         center[DIM];
140:   PetscReal         radius;
141:   AdvectSolBumpType type;
142: } Physics_Advect_Bump;

144: typedef struct {
145:   PetscReal     inflowState;
146:   AdvectSolType soltype;
147:   union {
148:     Physics_Advect_Tilted tilted;
149:     Physics_Advect_Bump   bump;
150:   } sol;
151:   struct {
152:     PetscInt Solution;
153:     PetscInt Error;
154:   } functional;
155: } Physics_Advect;

157: static const struct FieldDescription PhysicsFields_Advect[] = {{"U",1},{NULL,0}};

159: static PetscErrorCode PhysicsBoundary_Advect_Inflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
160: {
161:   Physics        phys    = (Physics)ctx;
162:   Physics_Advect *advect = (Physics_Advect*)phys->data;

165:   xG[0] = advect->inflowState;
166:   return 0;
167: }

169: static PetscErrorCode PhysicsBoundary_Advect_Outflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
170: {
172:   xG[0] = xI[0];
173:   return 0;
174: }

176: static void PhysicsRiemann_Advect(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, Physics phys)
177: {
178:   Physics_Advect *advect = (Physics_Advect*)phys->data;
179:   PetscReal      wind[DIM],wn;

181:   switch (advect->soltype) {
182:   case ADVECT_SOL_TILTED: {
183:     Physics_Advect_Tilted *tilted = &advect->sol.tilted;
184:     wind[0] = tilted->wind[0];
185:     wind[1] = tilted->wind[1];
186:   } break;
187:   case ADVECT_SOL_BUMP:
188:     wind[0] = -qp[1];
189:     wind[1] = qp[0];
190:     break;
191:   case ADVECT_SOL_BUMP_CAVITY:
192:     {
193:       PetscInt  i;
194:       PetscReal comp2[3] = {0.,0.,0.}, rad2;

196:       rad2 = 0.;
197:       for (i = 0; i < dim; i++) {
198:         comp2[i] = qp[i] * qp[i];
199:         rad2    += comp2[i];
200:       }

202:       wind[0] = -qp[1];
203:       wind[1] = qp[0];
204:       if (rad2 > 1.) {
205:         PetscInt  maxI = 0;
206:         PetscReal maxComp2 = comp2[0];

208:         for (i = 1; i < dim; i++) {
209:           if (comp2[i] > maxComp2) {
210:             maxI     = i;
211:             maxComp2 = comp2[i];
212:           }
213:         }
214:         wind[maxI] = 0.;
215:       }
216:     }
217:     break;
218:   default:
219:   {
220:     PetscInt i;
221:     for (i = 0; i < DIM; ++i) wind[i] = 0.0;
222:   }
223:   /* default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for solution type %s",AdvectSolBumpTypes[advect->soltype]); */
224:   }
225:   wn      = Dot2Real(wind, n);
226:   flux[0] = (wn > 0 ? xL[0] : xR[0]) * wn;
227: }

229: static PetscErrorCode PhysicsSolution_Advect(Model mod,PetscReal time,const PetscReal *x,PetscScalar *u,void *ctx)
230: {
231:   Physics        phys    = (Physics)ctx;
232:   Physics_Advect *advect = (Physics_Advect*)phys->data;

235:   switch (advect->soltype) {
236:   case ADVECT_SOL_TILTED: {
237:     PetscReal             x0[DIM];
238:     Physics_Advect_Tilted *tilted = &advect->sol.tilted;
239:     Waxpy2Real(-time,tilted->wind,x,x0);
240:     if (x0[1] > 0) u[0] = 1.*x[0] + 3.*x[1];
241:     else u[0] = advect->inflowState;
242:   } break;
243:   case ADVECT_SOL_BUMP_CAVITY:
244:   case ADVECT_SOL_BUMP: {
245:     Physics_Advect_Bump *bump = &advect->sol.bump;
246:     PetscReal           x0[DIM],v[DIM],r,cost,sint;
247:     cost  = PetscCosReal(time);
248:     sint  = PetscSinReal(time);
249:     x0[0] = cost*x[0] + sint*x[1];
250:     x0[1] = -sint*x[0] + cost*x[1];
251:     Waxpy2Real(-1,bump->center,x0,v);
252:     r = Norm2Real(v);
253:     switch (bump->type) {
254:     case ADVECT_SOL_BUMP_CONE:
255:       u[0] = PetscMax(1 - r/bump->radius,0);
256:       break;
257:     case ADVECT_SOL_BUMP_COS:
258:       u[0] = 0.5 + 0.5*PetscCosReal(PetscMin(r/bump->radius,1)*PETSC_PI);
259:       break;
260:     }
261:   } break;
262:   default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown solution type");
263:   }
264:   return 0;
265: }

267: static PetscErrorCode PhysicsFunctional_Advect(Model mod,PetscReal time,const PetscReal *x,const PetscScalar *y,PetscReal *f,void *ctx)
268: {
269:   Physics        phys    = (Physics)ctx;
270:   Physics_Advect *advect = (Physics_Advect*)phys->data;
271:   PetscScalar    yexact[1] = {0.0};

274:   PhysicsSolution_Advect(mod,time,x,yexact,phys);
275:   f[advect->functional.Solution] = PetscRealPart(y[0]);
276:   f[advect->functional.Error] = PetscAbsScalar(y[0]-yexact[0]);
277:   return 0;
278: }

280: static PetscErrorCode SetUpBC_Advect(DM dm, PetscDS prob, Physics phys)
281: {
282:   const PetscInt inflowids[] = {100,200,300},outflowids[] = {101};
283:   DMLabel        label;

286:   /* Register "canned" boundary conditions and defaults for where to apply. */
287:   DMGetLabel(dm, "Face Sets", &label);
288:   PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "inflow",  label, ALEN(inflowids),  inflowids,  0, 0, NULL, (void (*)(void)) PhysicsBoundary_Advect_Inflow, NULL,  phys, NULL);
289:   PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "outflow", label, ALEN(outflowids), outflowids, 0, 0, NULL, (void (*)(void)) PhysicsBoundary_Advect_Outflow, NULL, phys, NULL);
290:   return 0;
291: }

293: static PetscErrorCode PhysicsCreate_Advect(Model mod,Physics phys,PetscOptionItems *PetscOptionsObject)
294: {
295:   Physics_Advect *advect;

298:   phys->field_desc = PhysicsFields_Advect;
299:   phys->riemann    = (PetscRiemannFunc)PhysicsRiemann_Advect;
300:   PetscNew(&advect);
301:   phys->data       = advect;
302:   mod->setupbc = SetUpBC_Advect;

304:   PetscOptionsHead(PetscOptionsObject,"Advect options");
305:   {
306:     PetscInt two = 2,dof = 1;
307:     advect->soltype = ADVECT_SOL_TILTED;
308:     PetscOptionsEnum("-advect_sol_type","solution type","",AdvectSolTypes,(PetscEnum)advect->soltype,(PetscEnum*)&advect->soltype,NULL);
309:     switch (advect->soltype) {
310:     case ADVECT_SOL_TILTED: {
311:       Physics_Advect_Tilted *tilted = &advect->sol.tilted;
312:       two = 2;
313:       tilted->wind[0] = 0.0;
314:       tilted->wind[1] = 1.0;
315:       PetscOptionsRealArray("-advect_tilted_wind","background wind vx,vy","",tilted->wind,&two,NULL);
316:       advect->inflowState = -2.0;
317:       PetscOptionsRealArray("-advect_tilted_inflow","Inflow state","",&advect->inflowState,&dof,NULL);
318:       phys->maxspeed = Norm2Real(tilted->wind);
319:     } break;
320:     case ADVECT_SOL_BUMP_CAVITY:
321:     case ADVECT_SOL_BUMP: {
322:       Physics_Advect_Bump *bump = &advect->sol.bump;
323:       two = 2;
324:       bump->center[0] = 2.;
325:       bump->center[1] = 0.;
326:       PetscOptionsRealArray("-advect_bump_center","location of center of bump x,y","",bump->center,&two,NULL);
327:       bump->radius = 0.9;
328:       PetscOptionsReal("-advect_bump_radius","radius of bump","",bump->radius,&bump->radius,NULL);
329:       bump->type = ADVECT_SOL_BUMP_CONE;
330:       PetscOptionsEnum("-advect_bump_type","type of bump","",AdvectSolBumpTypes,(PetscEnum)bump->type,(PetscEnum*)&bump->type,NULL);
331:       phys->maxspeed = 3.;       /* radius of mesh, kludge */
332:     } break;
333:     }
334:   }
335:   PetscOptionsTail();
336:   /* Initial/transient solution with default boundary conditions */
337:   ModelSolutionSetDefault(mod,PhysicsSolution_Advect,phys);
338:   /* Register "canned" functionals */
339:   ModelFunctionalRegister(mod,"Solution",&advect->functional.Solution,PhysicsFunctional_Advect,phys);
340:   ModelFunctionalRegister(mod,"Error",&advect->functional.Error,PhysicsFunctional_Advect,phys);
341:   return 0;
342: }

344: /******************* Shallow Water ********************/
345: typedef struct {
346:   PetscReal gravity;
347:   PetscReal boundaryHeight;
348:   struct {
349:     PetscInt Height;
350:     PetscInt Speed;
351:     PetscInt Energy;
352:   } functional;
353: } Physics_SW;
354: typedef struct {
355:   PetscReal h;
356:   PetscReal uh[DIM];
357: } SWNode;
358: typedef union {
359:   SWNode    swnode;
360:   PetscReal vals[DIM+1];
361: } SWNodeUnion;

363: static const struct FieldDescription PhysicsFields_SW[] = {{"Height",1},{"Momentum",DIM},{NULL,0}};

365: /*
366:  * h_t + div(uh) = 0
367:  * (uh)_t + div (u\otimes uh + g h^2 / 2 I) = 0
368:  *
369:  * */
370: static PetscErrorCode SWFlux(Physics phys,const PetscReal *n,const SWNode *x,SWNode *f)
371: {
372:   Physics_SW  *sw = (Physics_SW*)phys->data;
373:   PetscReal   uhn,u[DIM];
374:   PetscInt     i;

377:   Scale2Real(1./x->h,x->uh,u);
378:   uhn  = x->uh[0] * n[0] + x->uh[1] * n[1];
379:   f->h = uhn;
380:   for (i=0; i<DIM; i++) f->uh[i] = u[i] * uhn + sw->gravity * PetscSqr(x->h) * n[i];
381:   return 0;
382: }

384: static PetscErrorCode PhysicsBoundary_SW_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
385: {
387:   xG[0] = xI[0];
388:   xG[1] = -xI[1];
389:   xG[2] = -xI[2];
390:   return 0;
391: }

393: static void PhysicsRiemann_SW_HLL(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, Physics phys)
394: {
395:   Physics_SW *sw = (Physics_SW *) phys->data;
396:   PetscReal aL, aR;
397:   PetscReal nn[DIM];
398: #if !defined(PETSC_USE_COMPLEX)
399:   const SWNode *uL = (const SWNode *) xL, *uR = (const SWNode *) xR;
400: #else
401:   SWNodeUnion  uLreal, uRreal;
402:   const SWNode *uL = &uLreal.swnode;
403:   const SWNode *uR = &uRreal.swnode;
404: #endif
405:   SWNodeUnion fL, fR;
406:   PetscInt i;
407:   PetscReal zero = 0.;

409: #if defined(PETSC_USE_COMPLEX)
410:   uLreal.swnode.h = 0; uRreal.swnode.h = 0;
411:   for (i = 0; i < 1+dim; i++) uLreal.vals[i] = PetscRealPart(xL[i]);
412:   for (i = 0; i < 1+dim; i++) uRreal.vals[i] = PetscRealPart(xR[i]);
413: #endif
414:   if (uL->h <= 0 || uR->h <= 0) {
415:     for (i = 0; i < 1 + dim; i++) flux[i] = zero;
416:     return;
417:   } /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reconstructed thickness is negative"); */
418:   nn[0] = n[0];
419:   nn[1] = n[1];
420:   Normalize2Real(nn);
421:   SWFlux(phys, nn, uL, &(fL.swnode));
422:   SWFlux(phys, nn, uR, &(fR.swnode));
423:   /* gravity wave speed */
424:   aL = PetscSqrtReal(sw->gravity * uL->h);
425:   aR = PetscSqrtReal(sw->gravity * uR->h);
426:   // Defining u_tilda and v_tilda as u and v
427:   PetscReal u_L, u_R;
428:   u_L = Dot2Real(uL->uh,nn)/uL->h;
429:   u_R = Dot2Real(uR->uh,nn)/uR->h;
430:   PetscReal sL, sR;
431:   sL = PetscMin(u_L - aL, u_R - aR);
432:   sR = PetscMax(u_L + aL, u_R + aR);
433:   if (sL > zero) {
434:     for (i = 0; i < dim + 1; i++) {
435:       flux[i] = fL.vals[i] * Norm2Real(n);
436:     }
437:   } else if (sR < zero) {
438:     for (i = 0; i < dim + 1; i++) {
439:       flux[i] = fR.vals[i] * Norm2Real(n);
440:     }
441:   } else {
442:     for (i = 0; i < dim + 1; i++) {
443:       flux[i] = ((sR * fL.vals[i] - sL * fR.vals[i] + sR * sL * (xR[i] - xL[i])) / (sR - sL)) * Norm2Real(n);
444:     }
445:   }
446: }

448: static void PhysicsRiemann_SW_Rusanov(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, Physics phys)
449: {
450:   Physics_SW   *sw = (Physics_SW*)phys->data;
451:   PetscReal    cL,cR,speed;
452:   PetscReal    nn[DIM];
453: #if !defined(PETSC_USE_COMPLEX)
454:   const SWNode *uL = (const SWNode*)xL,*uR = (const SWNode*)xR;
455: #else
456:   SWNodeUnion  uLreal, uRreal;
457:   const SWNode *uL = &uLreal.swnode;
458:   const SWNode *uR = &uRreal.swnode;
459: #endif
460:   SWNodeUnion  fL,fR;
461:   PetscInt     i;
462:   PetscReal    zero=0.;

464: #if defined(PETSC_USE_COMPLEX)
465:   uLreal.swnode.h = 0; uRreal.swnode.h = 0;
466:   for (i = 0; i < 1+dim; i++) uLreal.vals[i] = PetscRealPart(xL[i]);
467:   for (i = 0; i < 1+dim; i++) uRreal.vals[i] = PetscRealPart(xR[i]);
468: #endif
469:   if (uL->h < 0 || uR->h < 0) {for (i=0; i<1+dim; i++) flux[i] = zero/zero; return;} /* reconstructed thickness is negative */
470:   nn[0] = n[0];
471:   nn[1] = n[1];
472:   Normalize2Real(nn);
473:   SWFlux(phys,nn,uL,&(fL.swnode));
474:   SWFlux(phys,nn,uR,&(fR.swnode));
475:   cL    = PetscSqrtReal(sw->gravity*uL->h);
476:   cR    = PetscSqrtReal(sw->gravity*uR->h); /* gravity wave speed */
477:   speed = PetscMax(PetscAbsReal(Dot2Real(uL->uh,nn)/uL->h) + cL,PetscAbsReal(Dot2Real(uR->uh,nn)/uR->h) + cR);
478:   for (i=0; i<1+dim; i++) flux[i] = (0.5*(fL.vals[i] + fR.vals[i]) + 0.5*speed*(xL[i] - xR[i])) * Norm2Real(n);
479: }

481: static PetscErrorCode PhysicsSolution_SW(Model mod,PetscReal time,const PetscReal *x,PetscScalar *u,void *ctx)
482: {
483:   PetscReal dx[2],r,sigma;

487:   dx[0] = x[0] - 1.5;
488:   dx[1] = x[1] - 1.0;
489:   r     = Norm2Real(dx);
490:   sigma = 0.5;
491:   u[0]  = 1 + 2*PetscExpReal(-PetscSqr(r)/(2*PetscSqr(sigma)));
492:   u[1]  = 0.0;
493:   u[2]  = 0.0;
494:   return 0;
495: }

497: static PetscErrorCode PhysicsFunctional_SW(Model mod,PetscReal time,const PetscReal *coord,const PetscScalar *xx,PetscReal *f,void *ctx)
498: {
499:   Physics      phys = (Physics)ctx;
500:   Physics_SW   *sw  = (Physics_SW*)phys->data;
501:   const SWNode *x   = (const SWNode*)xx;
502:   PetscReal  u[2];
503:   PetscReal    h;

506:   h = x->h;
507:   Scale2Real(1./x->h,x->uh,u);
508:   f[sw->functional.Height] = h;
509:   f[sw->functional.Speed]  = Norm2Real(u) + PetscSqrtReal(sw->gravity*h);
510:   f[sw->functional.Energy] = 0.5*(Dot2Real(x->uh,u) + sw->gravity*PetscSqr(h));
511:   return 0;
512: }

514: static PetscErrorCode SetUpBC_SW(DM dm, PetscDS prob,Physics phys)
515: {
516:   const PetscInt wallids[] = {100,101,200,300};
517:   DMLabel        label;

520:   DMGetLabel(dm, "Face Sets", &label);
521:   PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", label, ALEN(wallids), wallids, 0, 0, NULL, (void (*)(void)) PhysicsBoundary_SW_Wall, NULL, phys, NULL);
522:   return 0;
523: }

525: static PetscErrorCode PhysicsCreate_SW(Model mod,Physics phys,PetscOptionItems *PetscOptionsObject)
526: {
527:   Physics_SW     *sw;
528:   char           sw_riemann[64] = "rusanov";

531:   phys->field_desc = PhysicsFields_SW;
532:   PetscNew(&sw);
533:   phys->data    = sw;
534:   mod->setupbc  = SetUpBC_SW;

536:   PetscFunctionListAdd(&PhysicsRiemannList_SW, "rusanov", PhysicsRiemann_SW_Rusanov);
537:   PetscFunctionListAdd(&PhysicsRiemannList_SW, "hll", PhysicsRiemann_SW_HLL);

539:   PetscOptionsHead(PetscOptionsObject,"SW options");
540:   {
541:     void (*PhysicsRiemann_SW)(PetscInt, PetscInt, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, PetscInt, const PetscScalar, PetscScalar *, Physics);
542:     sw->gravity = 1.0;
543:     PetscOptionsReal("-sw_gravity","Gravitational constant","",sw->gravity,&sw->gravity,NULL);
544:     PetscOptionsFList("-sw_riemann","Riemann solver","",PhysicsRiemannList_SW,sw_riemann,sw_riemann,sizeof sw_riemann,NULL);
545:     PetscFunctionListFind(PhysicsRiemannList_SW,sw_riemann,&PhysicsRiemann_SW);
546:     phys->riemann = (PetscRiemannFunc) PhysicsRiemann_SW;
547:   }
548:   PetscOptionsTail();
549:   phys->maxspeed = PetscSqrtReal(2.0*sw->gravity); /* Mach 1 for depth of 2 */

551:   ModelSolutionSetDefault(mod,PhysicsSolution_SW,phys);
552:   ModelFunctionalRegister(mod,"Height",&sw->functional.Height,PhysicsFunctional_SW,phys);
553:   ModelFunctionalRegister(mod,"Speed",&sw->functional.Speed,PhysicsFunctional_SW,phys);
554:   ModelFunctionalRegister(mod,"Energy",&sw->functional.Energy,PhysicsFunctional_SW,phys);

556:   return 0;
557: }

559: /******************* Euler Density Shock (EULER_IV_SHOCK,EULER_SS_SHOCK) ********************/
560: /* An initial-value and self-similar solutions of the compressible Euler equations */
561: /* Ravi Samtaney and D. I. Pullin */
562: /* Phys. Fluids 8, 2650 (1996); http://dx.doi.org/10.1063/1.869050 */
563: typedef enum {EULER_PAR_GAMMA,EULER_PAR_RHOR,EULER_PAR_AMACH,EULER_PAR_ITANA,EULER_PAR_SIZE} EulerParamIdx;
564: typedef enum {EULER_IV_SHOCK,EULER_SS_SHOCK,EULER_SHOCK_TUBE,EULER_LINEAR_WAVE} EulerType;
565: typedef struct {
566:   PetscReal r;
567:   PetscReal ru[DIM];
568:   PetscReal E;
569: } EulerNode;
570: typedef union {
571:   EulerNode eulernode;
572:   PetscReal vals[DIM+2];
573: } EulerNodeUnion;
574: typedef PetscErrorCode (*EquationOfState)(const PetscReal*, const EulerNode*, PetscReal*);
575: typedef struct {
576:   EulerType       type;
577:   PetscReal       pars[EULER_PAR_SIZE];
578:   EquationOfState sound;
579:   struct {
580:     PetscInt Density;
581:     PetscInt Momentum;
582:     PetscInt Energy;
583:     PetscInt Pressure;
584:     PetscInt Speed;
585:   } monitor;
586: } Physics_Euler;

588: static const struct FieldDescription PhysicsFields_Euler[] = {{"Density",1},{"Momentum",DIM},{"Energy",1},{NULL,0}};

590: /* initial condition */
591: int initLinearWave(EulerNode *ux, const PetscReal gamma, const PetscReal coord[], const PetscReal Lx);
592: static PetscErrorCode PhysicsSolution_Euler(Model mod, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx)
593: {
594:   PetscInt i;
595:   Physics         phys = (Physics)ctx;
596:   Physics_Euler   *eu  = (Physics_Euler*)phys->data;
597:   EulerNode       *uu  = (EulerNode*)u;
598:   PetscReal        p0,gamma,c;

602:   for (i=0; i<DIM; i++) uu->ru[i] = 0.0; /* zero out initial velocity */
603:   /* set E and rho */
604:   gamma = eu->pars[EULER_PAR_GAMMA];

606:   if (eu->type==EULER_IV_SHOCK || eu->type==EULER_SS_SHOCK) {
607:     /******************* Euler Density Shock ********************/
608:     /* On initial-value and self-similar solutions of the compressible Euler equations */
609:     /* Ravi Samtaney and D. I. Pullin */
610:     /* Phys. Fluids 8, 2650 (1996); http://dx.doi.org/10.1063/1.869050 */
611:     /* initial conditions 1: left of shock, 0: left of discontinuity 2: right of discontinuity,  */
612:     p0 = 1.;
613:     if (x[0] < 0.0 + x[1]*eu->pars[EULER_PAR_ITANA]) {
614:       if (x[0] < mod->bounds[0]*0.5) { /* left of shock (1) */
615:         PetscReal amach,rho,press,gas1,p1;
616:         amach = eu->pars[EULER_PAR_AMACH];
617:         rho = 1.;
618:         press = p0;
619:         p1 = press*(1.0+2.0*gamma/(gamma+1.0)*(amach*amach-1.0));
620:         gas1 = (gamma-1.0)/(gamma+1.0);
621:         uu->r = rho*(p1/press+gas1)/(gas1*p1/press+1.0);
622:         uu->ru[0]   = ((uu->r - rho)*PetscSqrtReal(gamma*press/rho)*amach);
623:         uu->E = p1/(gamma-1.0) + .5/uu->r*uu->ru[0]*uu->ru[0];
624:       }
625:       else { /* left of discontinuity (0) */
626:         uu->r = 1.; /* rho = 1 */
627:         uu->E = p0/(gamma-1.0);
628:       }
629:     }
630:     else { /* right of discontinuity (2) */
631:       uu->r = eu->pars[EULER_PAR_RHOR];
632:       uu->E = p0/(gamma-1.0);
633:     }
634:   }
635:   else if (eu->type==EULER_SHOCK_TUBE) {
636:     /* For (x<x0) set (rho,u,p)=(8,0,10) and for (x>x0) set (rho,u,p)=(1,0,1). Choose x0 to the midpoint of the domain in the x-direction. */
637:     if (x[0] < 0.0) {
638:       uu->r = 8.;
639:       uu->E = 10./(gamma-1.);
640:     }
641:     else {
642:       uu->r = 1.;
643:       uu->E = 1./(gamma-1.);
644:     }
645:   }
646:   else if (eu->type==EULER_LINEAR_WAVE) {
647:     initLinearWave( uu, gamma, x, mod->bounds[1] - mod->bounds[0]);
648:   }
649:   else SETERRQ(mod->comm,PETSC_ERR_SUP,"Unknown type %d",eu->type);

651:   /* set phys->maxspeed: (mod->maxspeed = phys->maxspeed) in main; */
652:   eu->sound(&gamma,uu,&c);
653:   c = (uu->ru[0]/uu->r) + c;
654:   if (c > phys->maxspeed) phys->maxspeed = c;

656:   return 0;
657: }

659: static PetscErrorCode Pressure_PG(const PetscReal gamma,const EulerNode *x,PetscReal *p)
660: {
661:   PetscReal ru2;

664:   ru2  = DotDIMReal(x->ru,x->ru);
665:   (*p)=(x->E - 0.5*ru2/x->r)*(gamma - 1.0); /* (E - rho V^2/2)(gamma-1) = e rho (gamma-1) */
666:   return 0;
667: }

669: static PetscErrorCode SpeedOfSound_PG(const PetscReal *gamma, const EulerNode *x, PetscReal *c)
670: {
671:   PetscReal p;

674:   Pressure_PG(*gamma,x,&p);
676:   /* pars[EULER_PAR_GAMMA] = heat capacity ratio */
677:   (*c)=PetscSqrtReal(*gamma * p / x->r);
678:   return 0;
679: }

681: /*
682:  * x = (rho,rho*(u_1),...,rho*e)^T
683:  * x_t+div(f_1(x))+...+div(f_DIM(x)) = 0
684:  *
685:  * f_i(x) = u_i*x+(0,0,...,p,...,p*u_i)^T
686:  *
687:  */
688: static PetscErrorCode EulerFlux(Physics phys,const PetscReal *n,const EulerNode *x,EulerNode *f)
689: {
690:   Physics_Euler *eu = (Physics_Euler*)phys->data;
691:   PetscReal     nu,p;
692:   PetscInt      i;

695:   Pressure_PG(eu->pars[EULER_PAR_GAMMA],x,&p);
696:   nu = DotDIMReal(x->ru,n);
697:   f->r = nu;   /* A rho u */
698:   nu /= x->r;  /* A u */
699:   for (i=0; i<DIM; i++) f->ru[i] = nu * x->ru[i] + n[i]*p;  /* r u^2 + p */
700:   f->E = nu * (x->E + p); /* u(e+p) */
701:   return 0;
702: }

704: /* PetscReal* => EulerNode* conversion */
705: static PetscErrorCode PhysicsBoundary_Euler_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *a_xI, PetscScalar *a_xG, void *ctx)
706: {
707:   PetscInt    i;
708:   const EulerNode *xI = (const EulerNode*)a_xI;
709:   EulerNode       *xG = (EulerNode*)a_xG;
710:   Physics         phys = (Physics)ctx;
711:   Physics_Euler   *eu  = (Physics_Euler*)phys->data;
713:   xG->r = xI->r;           /* ghost cell density - same */
714:   xG->E = xI->E;           /* ghost cell energy - same */
715:   if (n[1] != 0.) {        /* top and bottom */
716:     xG->ru[0] =  xI->ru[0]; /* copy tang to wall */
717:     xG->ru[1] = -xI->ru[1]; /* reflect perp to t/b wall */
718:   }
719:   else { /* sides */
720:     for (i=0; i<DIM; i++) xG->ru[i] = xI->ru[i]; /* copy */
721:   }
722:   if (eu->type == EULER_LINEAR_WAVE) { /* debug */
723: #if 0
724:     PetscPrintf(PETSC_COMM_WORLD,"%s coord=%g,%g\n",PETSC_FUNCTION_NAME,c[0],c[1]);
725: #endif
726:   }
727:   return 0;
728: }
729: int godunovflux( const PetscScalar *ul, const PetscScalar *ur, PetscScalar *flux, const PetscReal *nn, const int *ndim, const PetscReal *gamma);
730: /* PetscReal* => EulerNode* conversion */
731: static void PhysicsRiemann_Euler_Godunov( PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n,
732:                                           const PetscScalar *xL, const PetscScalar *xR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, Physics phys)
733: {
734:   Physics_Euler   *eu = (Physics_Euler*)phys->data;
735:   PetscReal       cL,cR,speed,velL,velR,nn[DIM],s2;
736:   PetscInt        i;
737:   PetscErrorCode  ierr;

740:   for (i=0,s2=0.; i<DIM; i++) {
741:     nn[i] = n[i];
742:     s2 += nn[i]*nn[i];
743:   }
744:   s2 = PetscSqrtReal(s2); /* |n|_2 = sum(n^2)^1/2 */
745:   for (i=0.; i<DIM; i++) nn[i] /= s2;
746:   if (0) { /* Rusanov */
747:     const EulerNode *uL = (const EulerNode*)xL,*uR = (const EulerNode*)xR;
748:     EulerNodeUnion  fL,fR;
749:     EulerFlux(phys,nn,uL,&(fL.eulernode));
750:     EulerFlux(phys,nn,uR,&(fR.eulernode));
751:     eu->sound(&eu->pars[EULER_PAR_GAMMA],uL,&cL);if (ierr) exit(13);
752:     eu->sound(&eu->pars[EULER_PAR_GAMMA],uR,&cR);if (ierr) exit(14);
753:     velL = DotDIMReal(uL->ru,nn)/uL->r;
754:     velR = DotDIMReal(uR->ru,nn)/uR->r;
755:     speed = PetscMax(velR + cR, velL + cL);
756:     for (i=0; i<2+dim; i++) flux[i] = 0.5*((fL.vals[i]+fR.vals[i]) + speed*(xL[i] - xR[i]))*s2;
757:   }
758:   else {
759:     int dim = DIM;
760:     /* int iwave =  */
761:     godunovflux(xL, xR, flux, nn, &dim, &eu->pars[EULER_PAR_GAMMA]);
762:     for (i=0; i<2+dim; i++) flux[i] *= s2;
763:   }
764:   return;
765: }

767: static PetscErrorCode PhysicsFunctional_Euler(Model mod,PetscReal time,const PetscReal *coord,const PetscScalar *xx,PetscReal *f,void *ctx)
768: {
769:   Physics         phys = (Physics)ctx;
770:   Physics_Euler   *eu  = (Physics_Euler*)phys->data;
771:   const EulerNode *x   = (const EulerNode*)xx;
772:   PetscReal       p;

775:   f[eu->monitor.Density]  = x->r;
776:   f[eu->monitor.Momentum] = NormDIM(x->ru);
777:   f[eu->monitor.Energy]   = x->E;
778:   f[eu->monitor.Speed]    = NormDIM(x->ru)/x->r;
779:   Pressure_PG(eu->pars[EULER_PAR_GAMMA], x, &p);
780:   f[eu->monitor.Pressure] = p;
781:   return 0;
782: }

784: static PetscErrorCode SetUpBC_Euler(DM dm, PetscDS prob,Physics phys)
785: {
786:   Physics_Euler   *eu = (Physics_Euler *) phys->data;
787:   DMLabel         label;

790:   DMGetLabel(dm, "Face Sets", &label);
791:   if (eu->type == EULER_LINEAR_WAVE) {
792:     const PetscInt wallids[] = {100,101};
793:     PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", label, ALEN(wallids), wallids, 0, 0, NULL, (void (*)(void)) PhysicsBoundary_Euler_Wall, NULL, phys, NULL);
794:   }
795:   else {
796:     const PetscInt wallids[] = {100,101,200,300};
797:     PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", label, ALEN(wallids), wallids, 0, 0, NULL, (void (*)(void)) PhysicsBoundary_Euler_Wall, NULL, phys, NULL);
798:   }
799:   return 0;
800: }

802: static PetscErrorCode PhysicsCreate_Euler(Model mod,Physics phys,PetscOptionItems *PetscOptionsObject)
803: {
804:   Physics_Euler   *eu;

807:   phys->field_desc = PhysicsFields_Euler;
808:   phys->riemann = (PetscRiemannFunc) PhysicsRiemann_Euler_Godunov;
809:   PetscNew(&eu);
810:   phys->data    = eu;
811:   mod->setupbc = SetUpBC_Euler;
812:   PetscOptionsHead(PetscOptionsObject,"Euler options");
813:   {
814:     PetscReal alpha;
815:     char type[64] = "linear_wave";
816:     PetscBool  is;
817:     eu->pars[EULER_PAR_GAMMA] = 1.4;
818:     eu->pars[EULER_PAR_AMACH] = 2.02;
819:     eu->pars[EULER_PAR_RHOR] = 3.0;
820:     eu->pars[EULER_PAR_ITANA] = 0.57735026918963; /* angle of Euler self similar (SS) shock */
821:     PetscOptionsReal("-eu_gamma","Heat capacity ratio","",eu->pars[EULER_PAR_GAMMA],&eu->pars[EULER_PAR_GAMMA],NULL);
822:     PetscOptionsReal("-eu_amach","Shock speed (Mach)","",eu->pars[EULER_PAR_AMACH],&eu->pars[EULER_PAR_AMACH],NULL);
823:     PetscOptionsReal("-eu_rho2","Density right of discontinuity","",eu->pars[EULER_PAR_RHOR],&eu->pars[EULER_PAR_RHOR],NULL);
824:     alpha = 60.;
825:     PetscOptionsReal("-eu_alpha","Angle of discontinuity","",alpha,&alpha,NULL);
827:     eu->pars[EULER_PAR_ITANA] = 1./PetscTanReal( alpha * PETSC_PI / 180.0);
828:     PetscOptionsString("-eu_type","Type of Euler test","",type,type,sizeof(type),NULL);
829:     PetscStrcmp(type,"linear_wave", &is);
830:     if (is) {
831:       /* Remember this should be periodic */
832:       eu->type = EULER_LINEAR_WAVE;
833:       PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"linear_wave");
834:     }
835:     else {
837:       PetscStrcmp(type,"iv_shock", &is);
838:       if (is) {
839:         eu->type = EULER_IV_SHOCK;
840:         PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"iv_shock");
841:       }
842:       else {
843:         PetscStrcmp(type,"ss_shock", &is);
844:         if (is) {
845:           eu->type = EULER_SS_SHOCK;
846:           PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"ss_shock");
847:         }
848:         else {
849:           PetscStrcmp(type,"shock_tube", &is);
850:           if (is) eu->type = EULER_SHOCK_TUBE;
851:           else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Unknown Euler type %s",type);
852:           PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"shock_tube");
853:         }
854:       }
855:     }
856:   }
857:   PetscOptionsTail();
858:   eu->sound = SpeedOfSound_PG;
859:   phys->maxspeed = 0.; /* will get set in solution */
860:   ModelSolutionSetDefault(mod,PhysicsSolution_Euler,phys);
861:   ModelFunctionalRegister(mod,"Speed",&eu->monitor.Speed,PhysicsFunctional_Euler,phys);
862:   ModelFunctionalRegister(mod,"Energy",&eu->monitor.Energy,PhysicsFunctional_Euler,phys);
863:   ModelFunctionalRegister(mod,"Density",&eu->monitor.Density,PhysicsFunctional_Euler,phys);
864:   ModelFunctionalRegister(mod,"Momentum",&eu->monitor.Momentum,PhysicsFunctional_Euler,phys);
865:   ModelFunctionalRegister(mod,"Pressure",&eu->monitor.Pressure,PhysicsFunctional_Euler,phys);

867:   return 0;
868: }

870: static PetscErrorCode ErrorIndicator_Simple(PetscInt dim, PetscReal volume, PetscInt numComps, const PetscScalar u[], const PetscScalar grad[], PetscReal *error, void *ctx)
871: {
872:   PetscReal      err = 0.;
873:   PetscInt       i, j;

876:   for (i = 0; i < numComps; i++) {
877:     for (j = 0; j < dim; j++) {
878:       err += PetscSqr(PetscRealPart(grad[i * dim + j]));
879:     }
880:   }
881:   *error = volume * err;
882:   return 0;
883: }

885: PetscErrorCode CreatePartitionVec(DM dm, DM *dmCell, Vec *partition)
886: {
887:   PetscSF        sfPoint;
888:   PetscSection   coordSection;
889:   Vec            coordinates;
890:   PetscSection   sectionCell;
891:   PetscScalar    *part;
892:   PetscInt       cStart, cEnd, c;
893:   PetscMPIInt    rank;

896:   DMGetCoordinateSection(dm, &coordSection);
897:   DMGetCoordinatesLocal(dm, &coordinates);
898:   DMClone(dm, dmCell);
899:   DMGetPointSF(dm, &sfPoint);
900:   DMSetPointSF(*dmCell, sfPoint);
901:   DMSetCoordinateSection(*dmCell, PETSC_DETERMINE, coordSection);
902:   DMSetCoordinatesLocal(*dmCell, coordinates);
903:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
904:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &sectionCell);
905:   DMPlexGetHeightStratum(*dmCell, 0, &cStart, &cEnd);
906:   PetscSectionSetChart(sectionCell, cStart, cEnd);
907:   for (c = cStart; c < cEnd; ++c) {
908:     PetscSectionSetDof(sectionCell, c, 1);
909:   }
910:   PetscSectionSetUp(sectionCell);
911:   DMSetLocalSection(*dmCell, sectionCell);
912:   PetscSectionDestroy(&sectionCell);
913:   DMCreateLocalVector(*dmCell, partition);
914:   PetscObjectSetName((PetscObject)*partition, "partition");
915:   VecGetArray(*partition, &part);
916:   for (c = cStart; c < cEnd; ++c) {
917:     PetscScalar *p;

919:     DMPlexPointLocalRef(*dmCell, c, part, &p);
920:     p[0] = rank;
921:   }
922:   VecRestoreArray(*partition, &part);
923:   return 0;
924: }

926: PetscErrorCode CreateMassMatrix(DM dm, Vec *massMatrix, User user)
927: {
928:   DM                plex, dmMass, dmFace, dmCell, dmCoord;
929:   PetscSection      coordSection;
930:   Vec               coordinates, facegeom, cellgeom;
931:   PetscSection      sectionMass;
932:   PetscScalar       *m;
933:   const PetscScalar *fgeom, *cgeom, *coords;
934:   PetscInt          vStart, vEnd, v;

937:   DMConvert(dm, DMPLEX, &plex);
938:   DMGetCoordinateSection(dm, &coordSection);
939:   DMGetCoordinatesLocal(dm, &coordinates);
940:   DMClone(dm, &dmMass);
941:   DMSetCoordinateSection(dmMass, PETSC_DETERMINE, coordSection);
942:   DMSetCoordinatesLocal(dmMass, coordinates);
943:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &sectionMass);
944:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
945:   PetscSectionSetChart(sectionMass, vStart, vEnd);
946:   for (v = vStart; v < vEnd; ++v) {
947:     PetscInt numFaces;

949:     DMPlexGetSupportSize(dmMass, v, &numFaces);
950:     PetscSectionSetDof(sectionMass, v, numFaces*numFaces);
951:   }
952:   PetscSectionSetUp(sectionMass);
953:   DMSetLocalSection(dmMass, sectionMass);
954:   PetscSectionDestroy(&sectionMass);
955:   DMGetLocalVector(dmMass, massMatrix);
956:   VecGetArray(*massMatrix, &m);
957:   DMPlexGetGeometryFVM(plex, &facegeom, &cellgeom, NULL);
958:   VecGetDM(facegeom, &dmFace);
959:   VecGetArrayRead(facegeom, &fgeom);
960:   VecGetDM(cellgeom, &dmCell);
961:   VecGetArrayRead(cellgeom, &cgeom);
962:   DMGetCoordinateDM(dm, &dmCoord);
963:   VecGetArrayRead(coordinates, &coords);
964:   for (v = vStart; v < vEnd; ++v) {
965:     const PetscInt        *faces;
966:     PetscFVFaceGeom       *fgA, *fgB, *cg;
967:     PetscScalar           *vertex;
968:     PetscInt               numFaces, sides[2], f, g;

970:     DMPlexPointLocalRead(dmCoord, v, coords, &vertex);
971:     DMPlexGetSupportSize(dmMass, v, &numFaces);
972:     DMPlexGetSupport(dmMass, v, &faces);
973:     for (f = 0; f < numFaces; ++f) {
974:       sides[0] = faces[f];
975:       DMPlexPointLocalRead(dmFace, faces[f], fgeom, &fgA);
976:       for (g = 0; g < numFaces; ++g) {
977:         const PetscInt *cells = NULL;
978:         PetscReal      area   = 0.0;
979:         PetscInt       numCells;

981:         sides[1] = faces[g];
982:         DMPlexPointLocalRead(dmFace, faces[g], fgeom, &fgB);
983:         DMPlexGetJoin(dmMass, 2, sides, &numCells, &cells);
985:         DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cg);
986:         area += PetscAbsScalar((vertex[0] - cg->centroid[0])*(fgA->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1])*(fgA->centroid[0] - cg->centroid[0]));
987:         area += PetscAbsScalar((vertex[0] - cg->centroid[0])*(fgB->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1])*(fgB->centroid[0] - cg->centroid[0]));
988:         m[f*numFaces+g] = Dot2Real(fgA->normal, fgB->normal)*area*0.5;
989:         DMPlexRestoreJoin(dmMass, 2, sides, &numCells, &cells);
990:       }
991:     }
992:   }
993:   VecRestoreArrayRead(facegeom, &fgeom);
994:   VecRestoreArrayRead(cellgeom, &cgeom);
995:   VecRestoreArrayRead(coordinates, &coords);
996:   VecRestoreArray(*massMatrix, &m);
997:   DMDestroy(&dmMass);
998:   DMDestroy(&plex);
999:   return 0;
1000: }

1002: /* Behavior will be different for multi-physics or when using non-default boundary conditions */
1003: static PetscErrorCode ModelSolutionSetDefault(Model mod,SolutionFunction func,void *ctx)
1004: {
1006:   mod->solution    = func;
1007:   mod->solutionctx = ctx;
1008:   return 0;
1009: }

1011: static PetscErrorCode ModelFunctionalRegister(Model mod,const char *name,PetscInt *offset,FunctionalFunction func,void *ctx)
1012: {
1013:   FunctionalLink link,*ptr;
1014:   PetscInt       lastoffset = -1;

1017:   for (ptr=&mod->functionalRegistry; *ptr; ptr = &(*ptr)->next) lastoffset = (*ptr)->offset;
1018:   PetscNew(&link);
1019:   PetscStrallocpy(name,&link->name);
1020:   link->offset = lastoffset + 1;
1021:   link->func   = func;
1022:   link->ctx    = ctx;
1023:   link->next   = NULL;
1024:   *ptr         = link;
1025:   *offset      = link->offset;
1026:   return 0;
1027: }

1029: static PetscErrorCode ModelFunctionalSetFromOptions(Model mod,PetscOptionItems *PetscOptionsObject)
1030: {
1031:   PetscInt       i,j;
1032:   FunctionalLink link;
1033:   char           *names[256];

1036:   mod->numMonitored = ALEN(names);
1037:   PetscOptionsStringArray("-monitor","list of functionals to monitor","",names,&mod->numMonitored,NULL);
1038:   /* Create list of functionals that will be computed somehow */
1039:   PetscMalloc1(mod->numMonitored,&mod->functionalMonitored);
1040:   /* Create index of calls that we will have to make to compute these functionals (over-allocation in general). */
1041:   PetscMalloc1(mod->numMonitored,&mod->functionalCall);
1042:   mod->numCall = 0;
1043:   for (i=0; i<mod->numMonitored; i++) {
1044:     for (link=mod->functionalRegistry; link; link=link->next) {
1045:       PetscBool match;
1046:       PetscStrcasecmp(names[i],link->name,&match);
1047:       if (match) break;
1048:     }
1050:     mod->functionalMonitored[i] = link;
1051:     for (j=0; j<i; j++) {
1052:       if (mod->functionalCall[j]->func == link->func && mod->functionalCall[j]->ctx == link->ctx) goto next_name;
1053:     }
1054:     mod->functionalCall[mod->numCall++] = link; /* Just points to the first link using the result. There may be more results. */
1055: next_name:
1056:     PetscFree(names[i]);
1057:   }

1059:   /* Find out the maximum index of any functional computed by a function we will be calling (even if we are not using it) */
1060:   mod->maxComputed = -1;
1061:   for (link=mod->functionalRegistry; link; link=link->next) {
1062:     for (i=0; i<mod->numCall; i++) {
1063:       FunctionalLink call = mod->functionalCall[i];
1064:       if (link->func == call->func && link->ctx == call->ctx) {
1065:         mod->maxComputed = PetscMax(mod->maxComputed,link->offset);
1066:       }
1067:     }
1068:   }
1069:   return 0;
1070: }

1072: static PetscErrorCode FunctionalLinkDestroy(FunctionalLink *link)
1073: {
1074:   FunctionalLink l,next;

1077:   if (!link) return 0;
1078:   l     = *link;
1079:   *link = NULL;
1080:   for (; l; l=next) {
1081:     next = l->next;
1082:     PetscFree(l->name);
1083:     PetscFree(l);
1084:   }
1085:   return 0;
1086: }

1088: /* put the solution callback into a functional callback */
1089: static PetscErrorCode SolutionFunctional(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *modctx)
1090: {
1091:   Model          mod;
1092:   mod  = (Model) modctx;
1093:   (*mod->solution)(mod, time, x, u, mod->solutionctx);
1094:   return 0;
1095: }

1097: PetscErrorCode SetInitialCondition(DM dm, Vec X, User user)
1098: {
1099:   PetscErrorCode     (*func[1]) (PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
1100:   void               *ctx[1];
1101:   Model              mod = user->model;

1104:   func[0] = SolutionFunctional;
1105:   ctx[0]  = (void *) mod;
1106:   DMProjectFunction(dm,0.0,func,ctx,INSERT_ALL_VALUES,X);
1107:   return 0;
1108: }

1110: static PetscErrorCode OutputVTK(DM dm, const char *filename, PetscViewer *viewer)
1111: {
1113:   PetscViewerCreate(PetscObjectComm((PetscObject)dm), viewer);
1114:   PetscViewerSetType(*viewer, PETSCVIEWERVTK);
1115:   PetscViewerFileSetName(*viewer, filename);
1116:   return 0;
1117: }

1119: static PetscErrorCode MonitorVTK(TS ts,PetscInt stepnum,PetscReal time,Vec X,void *ctx)
1120: {
1121:   User           user = (User)ctx;
1122:   DM             dm, plex;
1123:   PetscViewer    viewer;
1124:   char           filename[PETSC_MAX_PATH_LEN],*ftable = NULL;
1125:   PetscReal      xnorm;

1128:   PetscObjectSetName((PetscObject) X, "u");
1129:   VecGetDM(X,&dm);
1130:   VecNorm(X,NORM_INFINITY,&xnorm);

1132:   if (stepnum >= 0) {
1133:     stepnum += user->monitorStepOffset;
1134:   }
1135:   if (stepnum >= 0) {           /* No summary for final time */
1136:     Model             mod = user->model;
1137:     Vec               cellgeom;
1138:     PetscInt          c,cStart,cEnd,fcount,i;
1139:     size_t            ftableused,ftablealloc;
1140:     const PetscScalar *cgeom,*x;
1141:     DM                dmCell;
1142:     DMLabel           vtkLabel;
1143:     PetscReal         *fmin,*fmax,*fintegral,*ftmp;

1145:     DMConvert(dm, DMPLEX, &plex);
1146:     DMPlexGetGeometryFVM(plex, NULL, &cellgeom, NULL);
1147:     fcount = mod->maxComputed+1;
1148:     PetscMalloc4(fcount,&fmin,fcount,&fmax,fcount,&fintegral,fcount,&ftmp);
1149:     for (i=0; i<fcount; i++) {
1150:       fmin[i]      = PETSC_MAX_REAL;
1151:       fmax[i]      = PETSC_MIN_REAL;
1152:       fintegral[i] = 0;
1153:     }
1154:     VecGetDM(cellgeom,&dmCell);
1155:     DMPlexGetSimplexOrBoxCells(dmCell,0,&cStart,&cEnd);
1156:     VecGetArrayRead(cellgeom,&cgeom);
1157:     VecGetArrayRead(X,&x);
1158:     DMGetLabel(dm,"vtk",&vtkLabel);
1159:     for (c = cStart; c < cEnd; ++c) {
1160:       PetscFVCellGeom       *cg;
1161:       const PetscScalar     *cx    = NULL;
1162:       PetscInt              vtkVal = 0;

1164:       /* not that these two routines as currently implemented work for any dm with a
1165:        * localSection/globalSection */
1166:       DMPlexPointLocalRead(dmCell,c,cgeom,&cg);
1167:       DMPlexPointGlobalRead(dm,c,x,&cx);
1168:       if (vtkLabel) DMLabelGetValue(vtkLabel,c,&vtkVal);
1169:       if (!vtkVal || !cx) continue;        /* ghost, or not a global cell */
1170:       for (i=0; i<mod->numCall; i++) {
1171:         FunctionalLink flink = mod->functionalCall[i];
1172:         (*flink->func)(mod,time,cg->centroid,cx,ftmp,flink->ctx);
1173:       }
1174:       for (i=0; i<fcount; i++) {
1175:         fmin[i]       = PetscMin(fmin[i],ftmp[i]);
1176:         fmax[i]       = PetscMax(fmax[i],ftmp[i]);
1177:         fintegral[i] += cg->volume * ftmp[i];
1178:       }
1179:     }
1180:     VecRestoreArrayRead(cellgeom,&cgeom);
1181:     VecRestoreArrayRead(X,&x);
1182:     DMDestroy(&plex);
1183:     MPI_Allreduce(MPI_IN_PLACE,fmin,fcount,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));
1184:     MPI_Allreduce(MPI_IN_PLACE,fmax,fcount,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));
1185:     MPI_Allreduce(MPI_IN_PLACE,fintegral,fcount,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)ts));

1187:     ftablealloc = fcount * 100;
1188:     ftableused  = 0;
1189:     PetscMalloc1(ftablealloc,&ftable);
1190:     for (i=0; i<mod->numMonitored; i++) {
1191:       size_t         countused;
1192:       char           buffer[256],*p;
1193:       FunctionalLink flink = mod->functionalMonitored[i];
1194:       PetscInt       id    = flink->offset;
1195:       if (i % 3) {
1196:         PetscArraycpy(buffer,"  ",2);
1197:         p    = buffer + 2;
1198:       } else if (i) {
1199:         char newline[] = "\n";
1200:         PetscMemcpy(buffer,newline,sizeof(newline)-1);
1201:         p    = buffer + sizeof(newline) - 1;
1202:       } else {
1203:         p = buffer;
1204:       }
1205:       PetscSNPrintfCount(p,sizeof buffer-(p-buffer),"%12s [%10.7g,%10.7g] int %10.7g",&countused,flink->name,(double)fmin[id],(double)fmax[id],(double)fintegral[id]);
1206:       countused--;
1207:       countused += p - buffer;
1208:       if (countused > ftablealloc-ftableused-1) { /* reallocate */
1209:         char *ftablenew;
1210:         ftablealloc = 2*ftablealloc + countused;
1211:         PetscMalloc(ftablealloc,&ftablenew);
1212:         PetscArraycpy(ftablenew,ftable,ftableused);
1213:         PetscFree(ftable);
1214:         ftable = ftablenew;
1215:       }
1216:       PetscArraycpy(ftable+ftableused,buffer,countused);
1217:       ftableused += countused;
1218:       ftable[ftableused] = 0;
1219:     }
1220:     PetscFree4(fmin,fmax,fintegral,ftmp);

1222:     PetscPrintf(PetscObjectComm((PetscObject)ts),"% 3D  time %8.4g  |x| %8.4g  %s\n",stepnum,(double)time,(double)xnorm,ftable ? ftable : "");
1223:     PetscFree(ftable);
1224:   }
1225:   if (user->vtkInterval < 1) return 0;
1226:   if ((stepnum == -1) ^ (stepnum % user->vtkInterval == 0)) {
1227:     if (stepnum == -1) {        /* Final time is not multiple of normal time interval, write it anyway */
1228:       TSGetStepNumber(ts,&stepnum);
1229:     }
1230:     PetscSNPrintf(filename,sizeof filename,"%s-%03D.vtu",user->outputBasename,stepnum);
1231:     OutputVTK(dm,filename,&viewer);
1232:     VecView(X,viewer);
1233:     PetscViewerDestroy(&viewer);
1234:   }
1235:   return 0;
1236: }

1238: static PetscErrorCode initializeTS(DM dm, User user, TS *ts)
1239: {
1240:   TSCreate(PetscObjectComm((PetscObject)dm), ts);
1241:   TSSetType(*ts, TSSSP);
1242:   TSSetDM(*ts, dm);
1243:   if (user->vtkmon) {
1244:     TSMonitorSet(*ts,MonitorVTK,user,NULL);
1245:   }
1246:   DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, user);
1247:   DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, user);
1248:   TSSetMaxTime(*ts,2.0);
1249:   TSSetExactFinalTime(*ts,TS_EXACTFINALTIME_STEPOVER);
1250:   return 0;
1251: }

1253: static PetscErrorCode adaptToleranceFVM(PetscFV fvm, TS ts, Vec sol, VecTagger refineTag, VecTagger coarsenTag, User user, TS *tsNew, Vec *solNew)
1254: {
1255:   DM                dm, gradDM, plex, cellDM, adaptedDM = NULL;
1256:   Vec               cellGeom, faceGeom;
1257:   PetscBool         isForest, computeGradient;
1258:   Vec               grad, locGrad, locX, errVec;
1259:   PetscInt          cStart, cEnd, c, dim, nRefine, nCoarsen;
1260:   PetscReal         minMaxInd[2] = {PETSC_MAX_REAL, PETSC_MIN_REAL}, minMaxIndGlobal[2], minInd, maxInd, time;
1261:   PetscScalar       *errArray;
1262:   const PetscScalar *pointVals;
1263:   const PetscScalar *pointGrads;
1264:   const PetscScalar *pointGeom;
1265:   DMLabel           adaptLabel = NULL;
1266:   IS                refineIS, coarsenIS;

1268:   TSGetTime(ts,&time);
1269:   VecGetDM(sol, &dm);
1270:   DMGetDimension(dm,&dim);
1271:   PetscFVGetComputeGradients(fvm,&computeGradient);
1272:   PetscFVSetComputeGradients(fvm,PETSC_TRUE);
1273:   DMIsForest(dm, &isForest);
1274:   DMConvert(dm, DMPLEX, &plex);
1275:   DMPlexGetDataFVM(plex, fvm, &cellGeom, &faceGeom, &gradDM);
1276:   DMCreateLocalVector(plex,&locX);
1277:   DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, 0.0, faceGeom, cellGeom, NULL);
1278:   DMGlobalToLocalBegin(plex, sol, INSERT_VALUES, locX);
1279:   DMGlobalToLocalEnd  (plex, sol, INSERT_VALUES, locX);
1280:   DMCreateGlobalVector(gradDM, &grad);
1281:   DMPlexReconstructGradientsFVM(plex, locX, grad);
1282:   DMCreateLocalVector(gradDM, &locGrad);
1283:   DMGlobalToLocalBegin(gradDM, grad, INSERT_VALUES, locGrad);
1284:   DMGlobalToLocalEnd(gradDM, grad, INSERT_VALUES, locGrad);
1285:   VecDestroy(&grad);
1286:   DMPlexGetSimplexOrBoxCells(plex,0,&cStart,&cEnd);
1287:   VecGetArrayRead(locGrad,&pointGrads);
1288:   VecGetArrayRead(cellGeom,&pointGeom);
1289:   VecGetArrayRead(locX,&pointVals);
1290:   VecGetDM(cellGeom,&cellDM);
1291:   DMLabelCreate(PETSC_COMM_SELF,"adapt",&adaptLabel);
1292:   VecCreateMPI(PetscObjectComm((PetscObject)plex),cEnd-cStart,PETSC_DETERMINE,&errVec);
1293:   VecSetUp(errVec);
1294:   VecGetArray(errVec,&errArray);
1295:   for (c = cStart; c < cEnd; c++) {
1296:     PetscReal             errInd = 0.;
1297:     PetscScalar           *pointGrad;
1298:     PetscScalar           *pointVal;
1299:     PetscFVCellGeom       *cg;

1301:     DMPlexPointLocalRead(gradDM,c,pointGrads,&pointGrad);
1302:     DMPlexPointLocalRead(cellDM,c,pointGeom,&cg);
1303:     DMPlexPointLocalRead(plex,c,pointVals,&pointVal);

1305:     (user->model->errorIndicator)(dim,cg->volume,user->model->physics->dof,pointVal,pointGrad,&errInd,user->model->errorCtx);
1306:     errArray[c-cStart] = errInd;
1307:     minMaxInd[0] = PetscMin(minMaxInd[0],errInd);
1308:     minMaxInd[1] = PetscMax(minMaxInd[1],errInd);
1309:   }
1310:   VecRestoreArray(errVec,&errArray);
1311:   VecRestoreArrayRead(locX,&pointVals);
1312:   VecRestoreArrayRead(cellGeom,&pointGeom);
1313:   VecRestoreArrayRead(locGrad,&pointGrads);
1314:   VecDestroy(&locGrad);
1315:   VecDestroy(&locX);
1316:   DMDestroy(&plex);

1318:   VecTaggerComputeIS(refineTag,errVec,&refineIS,NULL);
1319:   VecTaggerComputeIS(coarsenTag,errVec,&coarsenIS,NULL);
1320:   ISGetSize(refineIS,&nRefine);
1321:   ISGetSize(coarsenIS,&nCoarsen);
1322:   if (nRefine) DMLabelSetStratumIS(adaptLabel,DM_ADAPT_REFINE,refineIS);
1323:   if (nCoarsen) DMLabelSetStratumIS(adaptLabel,DM_ADAPT_COARSEN,coarsenIS);
1324:   ISDestroy(&coarsenIS);
1325:   ISDestroy(&refineIS);
1326:   VecDestroy(&errVec);

1328:   PetscFVSetComputeGradients(fvm,computeGradient);
1329:   minMaxInd[1] = -minMaxInd[1];
1330:   MPI_Allreduce(minMaxInd,minMaxIndGlobal,2,MPIU_REAL,MPI_MIN,PetscObjectComm((PetscObject)dm));
1331:   minInd = minMaxIndGlobal[0];
1332:   maxInd = -minMaxIndGlobal[1];
1333:   PetscInfo(ts, "error indicator range (%E, %E)\n", minInd, maxInd);
1334:   if (nRefine || nCoarsen) { /* at least one cell is over the refinement threshold */
1335:     DMAdaptLabel(dm,adaptLabel,&adaptedDM);
1336:   }
1337:   DMLabelDestroy(&adaptLabel);
1338:   if (adaptedDM) {
1339:     PetscInfo(ts, "Adapted mesh, marking %D cells for refinement, and %D cells for coarsening\n", nRefine, nCoarsen);
1340:     if (tsNew) initializeTS(adaptedDM, user, tsNew);
1341:     if (solNew) {
1342:       DMCreateGlobalVector(adaptedDM, solNew);
1343:       PetscObjectSetName((PetscObject) *solNew, "solution");
1344:       DMForestTransferVec(dm, sol, adaptedDM, *solNew, PETSC_TRUE, time);
1345:     }
1346:     if (isForest) DMForestSetAdaptivityForest(adaptedDM,NULL); /* clear internal references to the previous dm */
1347:     DMDestroy(&adaptedDM);
1348:   } else {
1349:     if (tsNew)  *tsNew  = NULL;
1350:     if (solNew) *solNew = NULL;
1351:   }
1352:   return 0;
1353: }

1355: int main(int argc, char **argv)
1356: {
1357:   MPI_Comm          comm;
1358:   PetscDS           prob;
1359:   PetscFV           fvm;
1360:   PetscLimiter      limiter = NULL, noneLimiter = NULL;
1361:   User              user;
1362:   Model             mod;
1363:   Physics           phys;
1364:   DM                dm, plex;
1365:   PetscReal         ftime, cfl, dt, minRadius;
1366:   PetscInt          dim, nsteps;
1367:   TS                ts;
1368:   TSConvergedReason reason;
1369:   Vec               X;
1370:   PetscViewer       viewer;
1371:   PetscBool         vtkCellGeom, useAMR;
1372:   PetscInt          adaptInterval;
1373:   char              physname[256]  = "advect";
1374:   VecTagger         refineTag = NULL, coarsenTag = NULL;
1375:   PetscErrorCode    ierr;

1377:   PetscInitialize(&argc, &argv, (char*) 0, help);
1378:   comm = PETSC_COMM_WORLD;

1380:   PetscNew(&user);
1381:   PetscNew(&user->model);
1382:   PetscNew(&user->model->physics);
1383:   mod           = user->model;
1384:   phys          = mod->physics;
1385:   mod->comm     = comm;
1386:   useAMR        = PETSC_FALSE;
1387:   adaptInterval = 1;

1389:   /* Register physical models to be available on the command line */
1390:   PetscFunctionListAdd(&PhysicsList,"advect"          ,PhysicsCreate_Advect);
1391:   PetscFunctionListAdd(&PhysicsList,"sw"              ,PhysicsCreate_SW);
1392:   PetscFunctionListAdd(&PhysicsList,"euler"           ,PhysicsCreate_Euler);

1394:   PetscOptionsBegin(comm,NULL,"Unstructured Finite Volume Mesh Options","");
1395:   {
1396:     cfl  = 0.9 * 4; /* default SSPRKS2 with s=5 stages is stable for CFL number s-1 */
1397:     PetscOptionsReal("-ufv_cfl","CFL number per step","",cfl,&cfl,NULL);
1398:     user->vtkInterval = 1;
1399:     PetscOptionsInt("-ufv_vtk_interval","VTK output interval (0 to disable)","",user->vtkInterval,&user->vtkInterval,NULL);
1400:     user->vtkmon = PETSC_TRUE;
1401:     PetscOptionsBool("-ufv_vtk_monitor","Use VTKMonitor routine","",user->vtkmon,&user->vtkmon,NULL);
1402:     vtkCellGeom = PETSC_FALSE;
1403:     PetscStrcpy(user->outputBasename, "ex11");
1404:     PetscOptionsString("-ufv_vtk_basename","VTK output basename","",user->outputBasename,user->outputBasename,sizeof(user->outputBasename),NULL);
1405:     PetscOptionsBool("-ufv_vtk_cellgeom","Write cell geometry (for debugging)","",vtkCellGeom,&vtkCellGeom,NULL);
1406:     PetscOptionsBool("-ufv_use_amr","use local adaptive mesh refinement","",useAMR,&useAMR,NULL);
1407:     PetscOptionsInt("-ufv_adapt_interval","time steps between AMR","",adaptInterval,&adaptInterval,NULL);
1408:   }
1409:   PetscOptionsEnd();

1411:   if (useAMR) {
1412:     VecTaggerBox refineBox, coarsenBox;

1414:     refineBox.min  = refineBox.max  = PETSC_MAX_REAL;
1415:     coarsenBox.min = coarsenBox.max = PETSC_MIN_REAL;

1417:     VecTaggerCreate(comm,&refineTag);
1418:     PetscObjectSetOptionsPrefix((PetscObject)refineTag,"refine_");
1419:     VecTaggerSetType(refineTag,VECTAGGERABSOLUTE);
1420:     VecTaggerAbsoluteSetBox(refineTag,&refineBox);
1421:     VecTaggerSetFromOptions(refineTag);
1422:     VecTaggerSetUp(refineTag);
1423:     PetscObjectViewFromOptions((PetscObject)refineTag,NULL,"-tag_view");

1425:     VecTaggerCreate(comm,&coarsenTag);
1426:     PetscObjectSetOptionsPrefix((PetscObject)coarsenTag,"coarsen_");
1427:     VecTaggerSetType(coarsenTag,VECTAGGERABSOLUTE);
1428:     VecTaggerAbsoluteSetBox(coarsenTag,&coarsenBox);
1429:     VecTaggerSetFromOptions(coarsenTag);
1430:     VecTaggerSetUp(coarsenTag);
1431:     PetscObjectViewFromOptions((PetscObject)coarsenTag,NULL,"-tag_view");
1432:   }

1434:   PetscOptionsBegin(comm,NULL,"Unstructured Finite Volume Physics Options","");
1435:   {
1436:     PetscErrorCode (*physcreate)(Model,Physics,PetscOptionItems*);
1437:     PetscOptionsFList("-physics","Physics module to solve","",PhysicsList,physname,physname,sizeof physname,NULL);
1438:     PetscFunctionListFind(PhysicsList,physname,&physcreate);
1439:     PetscMemzero(phys,sizeof(struct _n_Physics));
1440:     (*physcreate)(mod,phys,PetscOptionsObject);
1441:     /* Count number of fields and dofs */
1442:     for (phys->nfields=0,phys->dof=0; phys->field_desc[phys->nfields].name; phys->nfields++) phys->dof += phys->field_desc[phys->nfields].dof;
1444:     ModelFunctionalSetFromOptions(mod,PetscOptionsObject);
1445:   }
1446:   PetscOptionsEnd();

1448:   /* Create mesh */
1449:   {
1450:     PetscInt i;

1452:     DMCreate(comm, &dm);
1453:     DMSetType(dm, DMPLEX);
1454:     DMSetFromOptions(dm);
1455:     for (i = 0; i < DIM; i++) { mod->bounds[2*i] = 0.; mod->bounds[2*i+1] = 1.;};
1456:     dim = DIM;
1457:     { /* a null name means just do a hex box */
1458:       PetscInt  cells[3] = {1, 1, 1}, n = 3;
1459:       PetscBool flg2, skew = PETSC_FALSE;
1460:       PetscInt nret2 = 2*DIM;
1461:       PetscOptionsBegin(comm,NULL,"Rectangular mesh options","");
1462:       PetscOptionsRealArray("-grid_bounds","bounds of the mesh in each direction (i.e., x_min,x_max,y_min,y_max","",mod->bounds,&nret2,&flg2);
1463:       PetscOptionsBool("-grid_skew_60","Skew grid for 60 degree shock mesh","",skew,&skew,NULL);
1464:       PetscOptionsIntArray("-dm_plex_box_faces", "Number of faces along each dimension", "", cells, &n, NULL);
1465:       PetscOptionsEnd();
1466:       /* TODO Rewrite this with Mark, and remove grid_bounds at that time */
1467:       if (flg2) {
1468:         PetscInt dimEmbed, i;
1469:         PetscInt nCoords;
1470:         PetscScalar *coords;
1471:         Vec coordinates;

1473:         DMGetCoordinatesLocal(dm,&coordinates);
1474:         DMGetCoordinateDim(dm,&dimEmbed);
1475:         VecGetLocalSize(coordinates,&nCoords);
1477:         VecGetArray(coordinates,&coords);
1478:         for (i = 0; i < nCoords; i += dimEmbed) {
1479:           PetscInt j;

1481:           PetscScalar *coord = &coords[i];
1482:           for (j = 0; j < dimEmbed; j++) {
1483:             coord[j] = mod->bounds[2 * j] + coord[j] * (mod->bounds[2 * j + 1] - mod->bounds[2 * j]);
1484:             if (dim==2 && cells[1]==1 && j==0 && skew) {
1485:               if (cells[0] == 2 && i == 8) {
1486:                 coord[j] = .57735026918963; /* hack to get 60 deg skewed mesh */
1487:               } else if (cells[0] == 3) {
1488:                 if (i==2 || i==10) coord[j] = mod->bounds[1]/4.;
1489:                 else if (i==4) coord[j] = mod->bounds[1]/2.;
1490:                 else if (i==12) coord[j] = 1.57735026918963*mod->bounds[1]/2.;
1491:               }
1492:             }
1493:           }
1494:         }
1495:         VecRestoreArray(coordinates,&coords);
1496:         DMSetCoordinatesLocal(dm,coordinates);
1497:       }
1498:     }
1499:   }
1500:   DMViewFromOptions(dm, NULL, "-orig_dm_view");
1501:   DMGetDimension(dm, &dim);

1503:   /* set up BCs, functions, tags */
1504:   DMCreateLabel(dm, "Face Sets");
1505:   mod->errorIndicator = ErrorIndicator_Simple;

1507:   {
1508:     DM gdm;

1510:     DMPlexConstructGhostCells(dm, NULL, NULL, &gdm);
1511:     DMDestroy(&dm);
1512:     dm   = gdm;
1513:     DMViewFromOptions(dm, NULL, "-dm_view");
1514:   }

1516:   PetscFVCreate(comm, &fvm);
1517:   PetscFVSetFromOptions(fvm);
1518:   PetscFVSetNumComponents(fvm, phys->dof);
1519:   PetscFVSetSpatialDimension(fvm, dim);
1520:   PetscObjectSetName((PetscObject) fvm,"");
1521:   {
1522:     PetscInt f, dof;
1523:     for (f=0,dof=0; f < phys->nfields; f++) {
1524:       PetscInt newDof = phys->field_desc[f].dof;

1526:       if (newDof == 1) {
1527:         PetscFVSetComponentName(fvm,dof,phys->field_desc[f].name);
1528:       }
1529:       else {
1530:         PetscInt j;

1532:         for (j = 0; j < newDof; j++) {
1533:           char     compName[256]  = "Unknown";

1535:           PetscSNPrintf(compName,sizeof(compName),"%s_%d",phys->field_desc[f].name,j);
1536:           PetscFVSetComponentName(fvm,dof+j,compName);
1537:         }
1538:       }
1539:       dof += newDof;
1540:     }
1541:   }
1542:   /* FV is now structured with one field having all physics as components */
1543:   DMAddField(dm, NULL, (PetscObject) fvm);
1544:   DMCreateDS(dm);
1545:   DMGetDS(dm, &prob);
1546:   PetscDSSetRiemannSolver(prob, 0, user->model->physics->riemann);
1547:   PetscDSSetContext(prob, 0, user->model->physics);
1548:   (*mod->setupbc)(dm, prob,phys);
1549:   PetscDSSetFromOptions(prob);
1550:   {
1551:     char      convType[256];
1552:     PetscBool flg;

1554:     PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");
1555:     PetscOptionsFList("-dm_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg);
1556:     PetscOptionsEnd();
1557:     if (flg) {
1558:       DM dmConv;

1560:       DMConvert(dm,convType,&dmConv);
1561:       if (dmConv) {
1562:         DMViewFromOptions(dmConv, NULL, "-dm_conv_view");
1563:         DMDestroy(&dm);
1564:         dm   = dmConv;
1565:         DMSetFromOptions(dm);
1566:       }
1567:     }
1568:   }

1570:   initializeTS(dm, user, &ts);

1572:   DMCreateGlobalVector(dm, &X);
1573:   PetscObjectSetName((PetscObject) X, "solution");
1574:   SetInitialCondition(dm, X, user);
1575:   if (useAMR) {
1576:     PetscInt adaptIter;

1578:     /* use no limiting when reconstructing gradients for adaptivity */
1579:     PetscFVGetLimiter(fvm, &limiter);
1580:     PetscObjectReference((PetscObject) limiter);
1581:     PetscLimiterCreate(PetscObjectComm((PetscObject) fvm), &noneLimiter);
1582:     PetscLimiterSetType(noneLimiter, PETSCLIMITERNONE);

1584:     PetscFVSetLimiter(fvm, noneLimiter);
1585:     for (adaptIter = 0; ; ++adaptIter) {
1586:       PetscLogDouble bytes;
1587:       TS             tsNew = NULL;

1589:       PetscMemoryGetCurrentUsage(&bytes);
1590:       PetscInfo(ts, "refinement loop %D: memory used %g\n", adaptIter, bytes);
1591:       DMViewFromOptions(dm, NULL, "-initial_dm_view");
1592:       VecViewFromOptions(X, NULL, "-initial_vec_view");
1593: #if 0
1594:       if (viewInitial) {
1595:         PetscViewer viewer;
1596:         char        buf[256];
1597:         PetscBool   isHDF5, isVTK;

1599:         PetscViewerCreate(comm,&viewer);
1600:         PetscViewerSetType(viewer,PETSCVIEWERVTK);
1601:         PetscViewerSetOptionsPrefix(viewer,"initial_");
1602:         PetscViewerSetFromOptions(viewer);
1603:         PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&isHDF5);
1604:         PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isVTK);
1605:         if (isHDF5) {
1606:           PetscSNPrintf(buf, 256, "ex11-initial-%d.h5", adaptIter);
1607:         } else if (isVTK) {
1608:           PetscSNPrintf(buf, 256, "ex11-initial-%d.vtu", adaptIter);
1609:           PetscViewerPushFormat(viewer,PETSC_VIEWER_VTK_VTU);
1610:         }
1611:         PetscViewerFileSetMode(viewer,FILE_MODE_WRITE);
1612:         PetscViewerFileSetName(viewer,buf);
1613:         if (isHDF5) {
1614:           DMView(dm,viewer);
1615:           PetscViewerFileSetMode(viewer,FILE_MODE_UPDATE);
1616:         }
1617:         VecView(X,viewer);
1618:         PetscViewerDestroy(&viewer);
1619:       }
1620: #endif

1622:       adaptToleranceFVM(fvm, ts, X, refineTag, coarsenTag, user, &tsNew, NULL);
1623:       if (!tsNew) {
1624:         break;
1625:       } else {
1626:         DMDestroy(&dm);
1627:         VecDestroy(&X);
1628:         TSDestroy(&ts);
1629:         ts   = tsNew;
1630:         TSGetDM(ts,&dm);
1631:         PetscObjectReference((PetscObject)dm);
1632:         DMCreateGlobalVector(dm,&X);
1633:         PetscObjectSetName((PetscObject) X, "solution");
1634:         SetInitialCondition(dm, X, user);
1635:       }
1636:     }
1637:     /* restore original limiter */
1638:     PetscFVSetLimiter(fvm, limiter);
1639:   }

1641:   DMConvert(dm, DMPLEX, &plex);
1642:   if (vtkCellGeom) {
1643:     DM  dmCell;
1644:     Vec cellgeom, partition;

1646:     DMPlexGetGeometryFVM(plex, NULL, &cellgeom, NULL);
1647:     OutputVTK(dm, "ex11-cellgeom.vtk", &viewer);
1648:     VecView(cellgeom, viewer);
1649:     PetscViewerDestroy(&viewer);
1650:     CreatePartitionVec(dm, &dmCell, &partition);
1651:     OutputVTK(dmCell, "ex11-partition.vtk", &viewer);
1652:     VecView(partition, viewer);
1653:     PetscViewerDestroy(&viewer);
1654:     VecDestroy(&partition);
1655:     DMDestroy(&dmCell);
1656:   }
1657:   /* collect max maxspeed from all processes -- todo */
1658:   DMPlexGetGeometryFVM(plex, NULL, NULL, &minRadius);
1659:   DMDestroy(&plex);
1660:   MPI_Allreduce(&phys->maxspeed,&mod->maxspeed,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));
1662:   dt   = cfl * minRadius / mod->maxspeed;
1663:   TSSetTimeStep(ts,dt);
1664:   TSSetFromOptions(ts);
1665:   if (!useAMR) {
1666:     TSSolve(ts,X);
1667:     TSGetSolveTime(ts,&ftime);
1668:     TSGetStepNumber(ts,&nsteps);
1669:   } else {
1670:     PetscReal finalTime;
1671:     PetscInt  adaptIter;
1672:     TS        tsNew = NULL;
1673:     Vec       solNew = NULL;

1675:     TSGetMaxTime(ts,&finalTime);
1676:     TSSetMaxSteps(ts,adaptInterval);
1677:     TSSolve(ts,X);
1678:     TSGetSolveTime(ts,&ftime);
1679:     TSGetStepNumber(ts,&nsteps);
1680:     for (adaptIter = 0;ftime < finalTime;adaptIter++) {
1681:       PetscLogDouble bytes;

1683:       PetscMemoryGetCurrentUsage(&bytes);
1684:       PetscInfo(ts, "AMR time step loop %D: memory used %g\n", adaptIter, bytes);
1685:       PetscFVSetLimiter(fvm,noneLimiter);
1686:       adaptToleranceFVM(fvm,ts,X,refineTag,coarsenTag,user,&tsNew,&solNew);
1687:       PetscFVSetLimiter(fvm,limiter);
1688:       if (tsNew) {
1689:         PetscInfo(ts, "AMR used\n");
1690:         DMDestroy(&dm);
1691:         VecDestroy(&X);
1692:         TSDestroy(&ts);
1693:         ts   = tsNew;
1694:         X    = solNew;
1695:         TSSetFromOptions(ts);
1696:         VecGetDM(X,&dm);
1697:         PetscObjectReference((PetscObject)dm);
1698:         DMConvert(dm, DMPLEX, &plex);
1699:         DMPlexGetGeometryFVM(dm, NULL, NULL, &minRadius);
1700:         DMDestroy(&plex);
1701:         MPI_Allreduce(&phys->maxspeed,&mod->maxspeed,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));
1703:         dt   = cfl * minRadius / mod->maxspeed;
1704:         TSSetStepNumber(ts,nsteps);
1705:         TSSetTime(ts,ftime);
1706:         TSSetTimeStep(ts,dt);
1707:       } else {
1708:         PetscInfo(ts, "AMR not used\n");
1709:       }
1710:       user->monitorStepOffset = nsteps;
1711:       TSSetMaxSteps(ts,nsteps+adaptInterval);
1712:       TSSolve(ts,X);
1713:       TSGetSolveTime(ts,&ftime);
1714:       TSGetStepNumber(ts,&nsteps);
1715:     }
1716:   }
1717:   TSGetConvergedReason(ts,&reason);
1718:   PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,nsteps);
1719:   TSDestroy(&ts);

1721:   VecTaggerDestroy(&refineTag);
1722:   VecTaggerDestroy(&coarsenTag);
1723:   PetscFunctionListDestroy(&PhysicsList);
1724:   PetscFunctionListDestroy(&PhysicsRiemannList_SW);
1725:   FunctionalLinkDestroy(&user->model->functionalRegistry);
1726:   PetscFree(user->model->functionalMonitored);
1727:   PetscFree(user->model->functionalCall);
1728:   PetscFree(user->model->physics->data);
1729:   PetscFree(user->model->physics);
1730:   PetscFree(user->model);
1731:   PetscFree(user);
1732:   VecDestroy(&X);
1733:   PetscLimiterDestroy(&limiter);
1734:   PetscLimiterDestroy(&noneLimiter);
1735:   PetscFVDestroy(&fvm);
1736:   DMDestroy(&dm);
1737:   PetscFinalize();
1738:   return 0;
1739: }

1741: /* Godunov fluxs */
1742: PetscScalar cvmgp_(PetscScalar *a, PetscScalar *b, PetscScalar *test)
1743: {
1744:     /* System generated locals */
1745:     PetscScalar ret_val;

1747:     if (PetscRealPart(*test) > 0.) {
1748:         goto L10;
1749:     }
1750:     ret_val = *b;
1751:     return ret_val;
1752: L10:
1753:     ret_val = *a;
1754:     return ret_val;
1755: } /* cvmgp_ */

1757: PetscScalar cvmgm_(PetscScalar *a, PetscScalar *b, PetscScalar *test)
1758: {
1759:     /* System generated locals */
1760:     PetscScalar ret_val;

1762:     if (PetscRealPart(*test) < 0.) {
1763:         goto L10;
1764:     }
1765:     ret_val = *b;
1766:     return ret_val;
1767: L10:
1768:     ret_val = *a;
1769:     return ret_val;
1770: } /* cvmgm_ */

1772: int riem1mdt( PetscScalar *gaml, PetscScalar *gamr, PetscScalar *rl, PetscScalar *pl,
1773:               PetscScalar *uxl, PetscScalar *rr, PetscScalar *pr,
1774:               PetscScalar *uxr, PetscScalar *rstarl, PetscScalar *rstarr, PetscScalar *
1775:               pstar, PetscScalar *ustar)
1776: {
1777:     /* Initialized data */

1779:     static PetscScalar smallp = 1e-8;

1781:     /* System generated locals */
1782:     int i__1;
1783:     PetscScalar d__1, d__2;

1785:     /* Local variables */
1786:     static int i0;
1787:     static PetscScalar cl, cr, wl, zl, wr, zr, pst, durl, skpr1, skpr2;
1788:     static int iwave;
1789:     static PetscScalar gascl4, gascr4, cstarl, dpstar, cstarr;
1790:     /* static PetscScalar csqrl, csqrr, gascl1, gascl2, gascl3, gascr1, gascr2, gascr3; */
1791:     static int iterno;
1792:     static PetscScalar ustarl, ustarr, rarepr1, rarepr2;

1794:     /* gascl1 = *gaml - 1.; */
1795:     /* gascl2 = (*gaml + 1.) * .5; */
1796:     /* gascl3 = gascl2 / *gaml; */
1797:     gascl4 = 1. / (*gaml - 1.);

1799:     /* gascr1 = *gamr - 1.; */
1800:     /* gascr2 = (*gamr + 1.) * .5; */
1801:     /* gascr3 = gascr2 / *gamr; */
1802:     gascr4 = 1. / (*gamr - 1.);
1803:     iterno = 10;
1804: /*        find pstar: */
1805:     cl = PetscSqrtScalar(*gaml * *pl / *rl);
1806:     cr = PetscSqrtScalar(*gamr * *pr / *rr);
1807:     wl = *rl * cl;
1808:     wr = *rr * cr;
1809:     /* csqrl = wl * wl; */
1810:     /* csqrr = wr * wr; */
1811:     *pstar = (wl * *pr + wr * *pl) / (wl + wr);
1812:     *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp));
1813:     pst = *pl / *pr;
1814:     skpr1 = cr * (pst - 1.) * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
1815:     d__1 = (*gamr - 1.) / (*gamr * 2.);
1816:     rarepr2 = gascr4 * 2. * cr * (1. - PetscPowScalar(pst, d__1));
1817:     pst = *pr / *pl;
1818:     skpr2 = cl * (pst - 1.) * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
1819:     d__1 = (*gaml - 1.) / (*gaml * 2.);
1820:     rarepr1 = gascl4 * 2. * cl * (1. - PetscPowScalar(pst, d__1));
1821:     durl = *uxr - *uxl;
1822:     if (PetscRealPart(*pr) < PetscRealPart(*pl)) {
1823:         if (PetscRealPart(durl) >= PetscRealPart(rarepr1)) {
1824:             iwave = 100;
1825:         } else if (PetscRealPart(durl) <= PetscRealPart(-skpr1)) {
1826:             iwave = 300;
1827:         } else {
1828:             iwave = 400;
1829:         }
1830:     } else {
1831:         if (PetscRealPart(durl) >= PetscRealPart(rarepr2)) {
1832:             iwave = 100;
1833:         } else if (PetscRealPart(durl) <= PetscRealPart(-skpr2)) {
1834:             iwave = 300;
1835:         } else {
1836:             iwave = 200;
1837:         }
1838:     }
1839:     if (iwave == 100) {
1840: /*     1-wave: rarefaction wave, 3-wave: rarefaction wave */
1841: /*     case (100) */
1842:         i__1 = iterno;
1843:         for (i0 = 1; i0 <= i__1; ++i0) {
1844:             d__1 = *pstar / *pl;
1845:             d__2 = 1. / *gaml;
1846:             *rstarl = *rl * PetscPowScalar(d__1, d__2);
1847:             cstarl = PetscSqrtScalar(*gaml * *pstar / *rstarl);
1848:             ustarl = *uxl - gascl4 * 2. * (cstarl - cl);
1849:             zl = *rstarl * cstarl;
1850:             d__1 = *pstar / *pr;
1851:             d__2 = 1. / *gamr;
1852:             *rstarr = *rr * PetscPowScalar(d__1, d__2);
1853:             cstarr = PetscSqrtScalar(*gamr * *pstar / *rstarr);
1854:             ustarr = *uxr + gascr4 * 2. * (cstarr - cr);
1855:             zr = *rstarr * cstarr;
1856:             dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
1857:             *pstar -= dpstar;
1858:             *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp));
1859:             if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
1860: #if 0
1861:         break;
1862: #endif
1863:             }
1864:         }
1865: /*     1-wave: shock wave, 3-wave: rarefaction wave */
1866:     } else if (iwave == 200) {
1867: /*     case (200) */
1868:         i__1 = iterno;
1869:         for (i0 = 1; i0 <= i__1; ++i0) {
1870:             pst = *pstar / *pl;
1871:             ustarl = *uxl - (pst - 1.) * cl * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
1872:             zl = *pl / cl * PetscSqrtScalar(*gaml * 2. * (*gaml - 1. + (*gaml + 1.) * pst)) * (*gaml - 1. + (*gaml + 1.) * pst) / (*gaml * 3. - 1. + (*gaml + 1.) * pst);
1873:             d__1 = *pstar / *pr;
1874:             d__2 = 1. / *gamr;
1875:             *rstarr = *rr * PetscPowScalar(d__1, d__2);
1876:             cstarr = PetscSqrtScalar(*gamr * *pstar / *rstarr);
1877:             zr = *rstarr * cstarr;
1878:             ustarr = *uxr + gascr4 * 2. * (cstarr - cr);
1879:             dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
1880:             *pstar -= dpstar;
1881:             *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp));
1882:             if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
1883: #if 0
1884:         break;
1885: #endif
1886:             }
1887:         }
1888: /*     1-wave: shock wave, 3-wave: shock */
1889:     } else if (iwave == 300) {
1890: /*     case (300) */
1891:         i__1 = iterno;
1892:         for (i0 = 1; i0 <= i__1; ++i0) {
1893:             pst = *pstar / *pl;
1894:             ustarl = *uxl - (pst - 1.) * cl * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
1895:             zl = *pl / cl * PetscSqrtScalar(*gaml * 2. * (*gaml - 1. + (*gaml + 1.) * pst)) * (*gaml - 1. + (*gaml + 1.) * pst) / (*gaml * 3. - 1. + (*gaml + 1.) * pst);
1896:             pst = *pstar / *pr;
1897:             ustarr = *uxr + (pst - 1.) * cr * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
1898:             zr = *pr / cr * PetscSqrtScalar(*gamr * 2. * (*gamr - 1. + (*gamr + 1.) * pst)) * (*gamr - 1. + (*gamr + 1.) * pst) / (*gamr * 3. - 1. + (*gamr + 1.) * pst);
1899:             dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
1900:             *pstar -= dpstar;
1901:             *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp));
1902:             if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
1903: #if 0
1904:         break;
1905: #endif
1906:             }
1907:         }
1908: /*     1-wave: rarefaction wave, 3-wave: shock */
1909:     } else if (iwave == 400) {
1910: /*     case (400) */
1911:         i__1 = iterno;
1912:         for (i0 = 1; i0 <= i__1; ++i0) {
1913:             d__1 = *pstar / *pl;
1914:             d__2 = 1. / *gaml;
1915:             *rstarl = *rl * PetscPowScalar(d__1, d__2);
1916:             cstarl = PetscSqrtScalar(*gaml * *pstar / *rstarl);
1917:             ustarl = *uxl - gascl4 * 2. * (cstarl - cl);
1918:             zl = *rstarl * cstarl;
1919:             pst = *pstar / *pr;
1920:             ustarr = *uxr + (pst - 1.) * cr * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
1921:             zr = *pr / cr * PetscSqrtScalar(*gamr * 2. * (*gamr - 1. + (*gamr + 1.) * pst)) * (*gamr - 1. + (*gamr + 1.) * pst) / (*gamr * 3. - 1. + (*gamr + 1.) * pst);
1922:             dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
1923:             *pstar -= dpstar;
1924:             *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp));
1925:             if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
1926: #if 0
1927:               break;
1928: #endif
1929:             }
1930:         }
1931:     }

1933:     *ustar = (zl * ustarr + zr * ustarl) / (zl + zr);
1934:     if (PetscRealPart(*pstar) > PetscRealPart(*pl)) {
1935:         pst = *pstar / *pl;
1936:         *rstarl = ((*gaml + 1.) * pst + *gaml - 1.) / ((*gaml - 1.) * pst + *
1937:                 gaml + 1.) * *rl;
1938:     }
1939:     if (PetscRealPart(*pstar) > PetscRealPart(*pr)) {
1940:         pst = *pstar / *pr;
1941:         *rstarr = ((*gamr + 1.) * pst + *gamr - 1.) / ((*gamr - 1.) * pst + *
1942:                 gamr + 1.) * *rr;
1943:     }
1944:     return iwave;
1945: }

1947: PetscScalar sign(PetscScalar x)
1948: {
1949:     if (PetscRealPart(x) > 0) return 1.0;
1950:     if (PetscRealPart(x) < 0) return -1.0;
1951:     return 0.0;
1952: }
1953: /*        Riemann Solver */
1954: /* -------------------------------------------------------------------- */
1955: int riemannsolver(PetscScalar *xcen, PetscScalar *xp,
1956:                    PetscScalar *dtt, PetscScalar *rl, PetscScalar *uxl, PetscScalar *pl,
1957:                    PetscScalar *utl, PetscScalar *ubl, PetscScalar *gaml, PetscScalar *rho1l,
1958:                    PetscScalar *rr, PetscScalar *uxr, PetscScalar *pr, PetscScalar *utr,
1959:                    PetscScalar *ubr, PetscScalar *gamr, PetscScalar *rho1r, PetscScalar *rx,
1960:                    PetscScalar *uxm, PetscScalar *px, PetscScalar *utx, PetscScalar *ubx,
1961:                    PetscScalar *gam, PetscScalar *rho1)
1962: {
1963:     /* System generated locals */
1964:     PetscScalar d__1, d__2;

1966:     /* Local variables */
1967:     static PetscScalar s, c0, p0, r0, u0, w0, x0, x2, ri, cx, sgn0, wsp0, gasc1, gasc2, gasc3, gasc4;
1968:     static PetscScalar cstar, pstar, rstar, ustar, xstar, wspst, ushock, streng, rstarl, rstarr, rstars;
1969:     int iwave;

1971:     if (*rl == *rr && *pr == *pl && *uxl == *uxr && *gaml == *gamr) {
1972:         *rx = *rl;
1973:         *px = *pl;
1974:         *uxm = *uxl;
1975:         *gam = *gaml;
1976:         x2 = *xcen + *uxm * *dtt;

1978:         if (PetscRealPart(*xp) >= PetscRealPart(x2)) {
1979:             *utx = *utr;
1980:             *ubx = *ubr;
1981:             *rho1 = *rho1r;
1982:         } else {
1983:             *utx = *utl;
1984:             *ubx = *ubl;
1985:             *rho1 = *rho1l;
1986:         }
1987:         return 0;
1988:     }
1989:     iwave = riem1mdt(gaml, gamr, rl, pl, uxl, rr, pr, uxr, &rstarl, &rstarr, &pstar, &ustar);

1991:     x2 = *xcen + ustar * *dtt;
1992:     d__1 = *xp - x2;
1993:     sgn0 = sign(d__1);
1994: /*            x is in 3-wave if sgn0 = 1 */
1995: /*            x is in 1-wave if sgn0 = -1 */
1996:     r0 = cvmgm_(rl, rr, &sgn0);
1997:     p0 = cvmgm_(pl, pr, &sgn0);
1998:     u0 = cvmgm_(uxl, uxr, &sgn0);
1999:     *gam = cvmgm_(gaml, gamr, &sgn0);
2000:     gasc1 = *gam - 1.;
2001:     gasc2 = (*gam + 1.) * .5;
2002:     gasc3 = gasc2 / *gam;
2003:     gasc4 = 1. / (*gam - 1.);
2004:     c0 = PetscSqrtScalar(*gam * p0 / r0);
2005:     streng = pstar - p0;
2006:     w0 = *gam * r0 * p0 * (gasc3 * streng / p0 + 1.);
2007:     rstars = r0 / (1. - r0 * streng / w0);
2008:     d__1 = p0 / pstar;
2009:     d__2 = -1. / *gam;
2010:     rstarr = r0 * PetscPowScalar(d__1, d__2);
2011:     rstar = cvmgm_(&rstarr, &rstars, &streng);
2012:     w0 = PetscSqrtScalar(w0);
2013:     cstar = PetscSqrtScalar(*gam * pstar / rstar);
2014:     wsp0 = u0 + sgn0 * c0;
2015:     wspst = ustar + sgn0 * cstar;
2016:     ushock = ustar + sgn0 * w0 / rstar;
2017:     wspst = cvmgp_(&ushock, &wspst, &streng);
2018:     wsp0 = cvmgp_(&ushock, &wsp0, &streng);
2019:     x0 = *xcen + wsp0 * *dtt;
2020:     xstar = *xcen + wspst * *dtt;
2021: /*           using gas formula to evaluate rarefaction wave */
2022: /*            ri : reiman invariant */
2023:     ri = u0 - sgn0 * 2. * gasc4 * c0;
2024:     cx = sgn0 * .5 * gasc1 / gasc2 * ((*xp - *xcen) / *dtt - ri);
2025:     *uxm = ri + sgn0 * 2. * gasc4 * cx;
2026:     s = p0 / PetscPowScalar(r0, *gam);
2027:     d__1 = cx * cx / (*gam * s);
2028:     *rx = PetscPowScalar(d__1, gasc4);
2029:     *px = cx * cx * *rx / *gam;
2030:     d__1 = sgn0 * (x0 - *xp);
2031:     *rx = cvmgp_(rx, &r0, &d__1);
2032:     d__1 = sgn0 * (x0 - *xp);
2033:     *px = cvmgp_(px, &p0, &d__1);
2034:     d__1 = sgn0 * (x0 - *xp);
2035:     *uxm = cvmgp_(uxm, &u0, &d__1);
2036:     d__1 = sgn0 * (xstar - *xp);
2037:     *rx = cvmgm_(rx, &rstar, &d__1);
2038:     d__1 = sgn0 * (xstar - *xp);
2039:     *px = cvmgm_(px, &pstar, &d__1);
2040:     d__1 = sgn0 * (xstar - *xp);
2041:     *uxm = cvmgm_(uxm, &ustar, &d__1);
2042:     if (PetscRealPart(*xp) >= PetscRealPart(x2)) {
2043:         *utx = *utr;
2044:         *ubx = *ubr;
2045:         *rho1 = *rho1r;
2046:     } else {
2047:         *utx = *utl;
2048:         *ubx = *ubl;
2049:         *rho1 = *rho1l;
2050:     }
2051:     return iwave;
2052: }
2053: int godunovflux( const PetscScalar *ul, const PetscScalar *ur,
2054:                  PetscScalar *flux, const PetscReal *nn, const int *ndim,
2055:                  const PetscReal *gamma)
2056: {
2057:     /* System generated locals */
2058:   int i__1,iwave;
2059:     PetscScalar d__1, d__2, d__3;

2061:     /* Local variables */
2062:     static int k;
2063:     static PetscScalar bn[3], fn, ft, tg[3], pl, rl, pm, pr, rr, xp, ubl, ubm,
2064:             ubr, dtt, unm, tmp, utl, utm, uxl, utr, uxr, gaml, gamm, gamr,
2065:             xcen, rhom, rho1l, rho1m, rho1r;

2067:     /* Function Body */
2068:     xcen = 0.;
2069:     xp = 0.;
2070:     i__1 = *ndim;
2071:     for (k = 1; k <= i__1; ++k) {
2072:         tg[k - 1] = 0.;
2073:         bn[k - 1] = 0.;
2074:     }
2075:     dtt = 1.;
2076:     if (*ndim == 3) {
2077:         if (nn[0] == 0. && nn[1] == 0.) {
2078:             tg[0] = 1.;
2079:         } else {
2080:             tg[0] = -nn[1];
2081:             tg[1] = nn[0];
2082:         }
2083: /*           tmp=dsqrt(tg(1)**2+tg(2)**2) */
2084: /*           tg=tg/tmp */
2085:         bn[0] = -nn[2] * tg[1];
2086:         bn[1] = nn[2] * tg[0];
2087:         bn[2] = nn[0] * tg[1] - nn[1] * tg[0];
2088: /* Computing 2nd power */
2089:         d__1 = bn[0];
2090: /* Computing 2nd power */
2091:         d__2 = bn[1];
2092: /* Computing 2nd power */
2093:         d__3 = bn[2];
2094:         tmp = PetscSqrtScalar(d__1 * d__1 + d__2 * d__2 + d__3 * d__3);
2095:         i__1 = *ndim;
2096:         for (k = 1; k <= i__1; ++k) {
2097:             bn[k - 1] /= tmp;
2098:         }
2099:     } else if (*ndim == 2) {
2100:         tg[0] = -nn[1];
2101:         tg[1] = nn[0];
2102: /*           tmp=dsqrt(tg(1)**2+tg(2)**2) */
2103: /*           tg=tg/tmp */
2104:         bn[0] = 0.;
2105:         bn[1] = 0.;
2106:         bn[2] = 1.;
2107:     }
2108:     rl = ul[0];
2109:     rr = ur[0];
2110:     uxl = 0.;
2111:     uxr = 0.;
2112:     utl = 0.;
2113:     utr = 0.;
2114:     ubl = 0.;
2115:     ubr = 0.;
2116:     i__1 = *ndim;
2117:     for (k = 1; k <= i__1; ++k) {
2118:         uxl += ul[k] * nn[k-1];
2119:         uxr += ur[k] * nn[k-1];
2120:         utl += ul[k] * tg[k - 1];
2121:         utr += ur[k] * tg[k - 1];
2122:         ubl += ul[k] * bn[k - 1];
2123:         ubr += ur[k] * bn[k - 1];
2124:     }
2125:     uxl /= rl;
2126:     uxr /= rr;
2127:     utl /= rl;
2128:     utr /= rr;
2129:     ubl /= rl;
2130:     ubr /= rr;

2132:     gaml = *gamma;
2133:     gamr = *gamma;
2134: /* Computing 2nd power */
2135:     d__1 = uxl;
2136: /* Computing 2nd power */
2137:     d__2 = utl;
2138: /* Computing 2nd power */
2139:     d__3 = ubl;
2140:     pl = (*gamma - 1.) * (ul[*ndim + 1] - rl * .5 * (d__1 * d__1 + d__2 * d__2 + d__3 * d__3));
2141: /* Computing 2nd power */
2142:     d__1 = uxr;
2143: /* Computing 2nd power */
2144:     d__2 = utr;
2145: /* Computing 2nd power */
2146:     d__3 = ubr;
2147:     pr = (*gamma - 1.) * (ur[*ndim + 1] - rr * .5 * (d__1 * d__1 + d__2 * d__2 + d__3 * d__3));
2148:     rho1l = rl;
2149:     rho1r = rr;

2151:     iwave = riemannsolver(&xcen, &xp, &dtt, &rl, &uxl, &pl, &utl, &ubl, &gaml, &
2152:                           rho1l, &rr, &uxr, &pr, &utr, &ubr, &gamr, &rho1r, &rhom, &unm, &
2153:                           pm, &utm, &ubm, &gamm, &rho1m);

2155:     flux[0] = rhom * unm;
2156:     fn = rhom * unm * unm + pm;
2157:     ft = rhom * unm * utm;
2158: /*           flux(2)=fn*nn(1)+ft*nn(2) */
2159: /*           flux(3)=fn*tg(1)+ft*tg(2) */
2160:     flux[1] = fn * nn[0] + ft * tg[0];
2161:     flux[2] = fn * nn[1] + ft * tg[1];
2162: /*           flux(2)=rhom*unm*(unm)+pm */
2163: /*           flux(3)=rhom*(unm)*utm */
2164:     if (*ndim == 3) {
2165:         flux[3] = rhom * unm * ubm;
2166:     }
2167:     flux[*ndim + 1] = (rhom * .5 * (unm * unm + utm * utm + ubm * ubm) + gamm / (gamm - 1.) * pm) * unm;
2168:     return iwave;
2169: } /* godunovflux_ */

2171: /* Subroutine to set up the initial conditions for the */
2172: /* Shock Interface interaction or linear wave (Ravi Samtaney,Mark Adams). */
2173: /* ----------------------------------------------------------------------- */
2174: int projecteqstate(PetscReal wc[], const PetscReal ueq[], PetscReal lv[][3])
2175: {
2176:   int j,k;
2177: /*      Wc=matmul(lv,Ueq) 3 vars */
2178:   for (k = 0; k < 3; ++k) {
2179:     wc[k] = 0.;
2180:     for (j = 0; j < 3; ++j) {
2181:       wc[k] += lv[k][j]*ueq[j];
2182:     }
2183:   }
2184:   return 0;
2185: }
2186: /* ----------------------------------------------------------------------- */
2187: int projecttoprim(PetscReal v[], const PetscReal wc[], PetscReal rv[][3])
2188: {
2189:   int k,j;
2190:   /*      V=matmul(rv,WC) 3 vars */
2191:   for (k = 0; k < 3; ++k) {
2192:     v[k] = 0.;
2193:     for (j = 0; j < 3; ++j) {
2194:       v[k] += rv[k][j]*wc[j];
2195:     }
2196:   }
2197:   return 0;
2198: }
2199: /* ---------------------------------------------------------------------- */
2200: int eigenvectors(PetscReal rv[][3], PetscReal lv[][3], const PetscReal ueq[], PetscReal gamma)
2201: {
2202:   int j,k;
2203:   PetscReal rho,csnd,p0;
2204:   /* PetscScalar u; */

2206:   for (k = 0; k < 3; ++k) for (j = 0; j < 3; ++j) { lv[k][j] = 0.; rv[k][j] = 0.; }
2207:   rho = ueq[0];
2208:   /* u = ueq[1]; */
2209:   p0 = ueq[2];
2210:   csnd = PetscSqrtReal(gamma * p0 / rho);
2211:   lv[0][1] = rho * .5;
2212:   lv[0][2] = -.5 / csnd;
2213:   lv[1][0] = csnd;
2214:   lv[1][2] = -1. / csnd;
2215:   lv[2][1] = rho * .5;
2216:   lv[2][2] = .5 / csnd;
2217:   rv[0][0] = -1. / csnd;
2218:   rv[1][0] = 1. / rho;
2219:   rv[2][0] = -csnd;
2220:   rv[0][1] = 1. / csnd;
2221:   rv[0][2] = 1. / csnd;
2222:   rv[1][2] = 1. / rho;
2223:   rv[2][2] = csnd;
2224:   return 0;
2225: }

2227: int initLinearWave(EulerNode *ux, const PetscReal gamma, const PetscReal coord[], const PetscReal Lx)
2228: {
2229:   PetscReal p0,u0,wcp[3],wc[3];
2230:   PetscReal lv[3][3];
2231:   PetscReal vp[3];
2232:   PetscReal rv[3][3];
2233:   PetscReal eps, ueq[3], rho0, twopi;

2235:   /* Function Body */
2236:   twopi = 2.*PETSC_PI;
2237:   eps = 1e-4; /* perturbation */
2238:   rho0 = 1e3;   /* density of water */
2239:   p0 = 101325.; /* init pressure of 1 atm (?) */
2240:   u0 = 0.;
2241:   ueq[0] = rho0;
2242:   ueq[1] = u0;
2243:   ueq[2] = p0;
2244:   /* Project initial state to characteristic variables */
2245:   eigenvectors(rv, lv, ueq, gamma);
2246:   projecteqstate(wc, ueq, lv);
2247:   wcp[0] = wc[0];
2248:   wcp[1] = wc[1];
2249:   wcp[2] = wc[2] + eps * PetscCosReal(coord[0] * 2. * twopi / Lx);
2250:   projecttoprim(vp, wcp, rv);
2251:   ux->r = vp[0]; /* density */
2252:   ux->ru[0] = vp[0] * vp[1]; /* x momentum */
2253:   ux->ru[1] = 0.;
2254: #if defined DIM > 2
2255:   if (dim>2) ux->ru[2] = 0.;
2256: #endif
2257:   /* E = rho * e + rho * v^2/2 = p/(gam-1) + rho*v^2/2 */
2258:   ux->E = vp[2]/(gamma - 1.) + 0.5*vp[0]*vp[1]*vp[1];
2259:   return 0;
2260: }

2262: /*TEST

2264:   testset:
2265:     args: -dm_plex_adj_cone -dm_plex_adj_closure 0

2267:     test:
2268:       suffix: adv_2d_tri_0
2269:       requires: triangle
2270:       TODO: how did this ever get in main when there is no support for this
2271:       args: -ufv_vtk_interval 0 -simplex -dm_refine 3 -dm_plex_faces 1,1 -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3

2273:     test:
2274:       suffix: adv_2d_tri_1
2275:       requires: triangle
2276:       TODO: how did this ever get in main when there is no support for this
2277:       args: -ufv_vtk_interval 0 -simplex -dm_refine 5 -dm_plex_faces 1,1 -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1

2279:     test:
2280:       suffix: tut_1
2281:       requires: exodusii
2282:       nsize: 1
2283:       args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo

2285:     test:
2286:       suffix: tut_2
2287:       requires: exodusii
2288:       nsize: 1
2289:       args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo -ts_type rosw

2291:     test:
2292:       suffix: tut_3
2293:       requires: exodusii
2294:       nsize: 4
2295:       args: -dm_distribute_overlap 1 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -monitor Error -advect_sol_type bump -petscfv_type leastsquares -petsclimiter_type sin

2297:     test:
2298:       suffix: tut_4
2299:       requires: exodusii
2300:       nsize: 4
2301:       args: -dm_distribute_overlap 1 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -physics sw -monitor Height,Energy -petscfv_type leastsquares -petsclimiter_type minmod

2303:   testset:
2304:     args: -dm_plex_adj_cone -dm_plex_adj_closure 0 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1

2306:     # 2D Advection 0-10
2307:     test:
2308:       suffix: 0
2309:       requires: exodusii
2310:       args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo

2312:     test:
2313:       suffix: 1
2314:       requires: exodusii
2315:       args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo

2317:     test:
2318:       suffix: 2
2319:       requires: exodusii
2320:       nsize: 2
2321:       args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo

2323:     test:
2324:       suffix: 3
2325:       requires: exodusii
2326:       nsize: 2
2327:       args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo

2329:     test:
2330:       suffix: 4
2331:       requires: exodusii
2332:       nsize: 8
2333:       args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo

2335:     test:
2336:       suffix: 5
2337:       requires: exodusii
2338:       args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo -ts_type rosw -ts_adapt_reject_safety 1

2340:     test:
2341:       suffix: 7
2342:       requires: exodusii
2343:       args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1

2345:     test:
2346:       suffix: 8
2347:       requires: exodusii
2348:       nsize: 2
2349:       args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1

2351:     test:
2352:       suffix: 9
2353:       requires: exodusii
2354:       nsize: 8
2355:       args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1

2357:     test:
2358:       suffix: 10
2359:       requires: exodusii
2360:       args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo

2362:     # 2D Shallow water
2363:     test:
2364:       suffix: sw_0
2365:       requires: exodusii
2366:       args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -bc_wall 100,101 -physics sw -ufv_cfl 5 -petscfv_type   leastsquares -petsclimiter_type sin -ts_max_time 1 -ts_ssp_type rks2 -ts_ssp_nstages 10 -monitor height,energy

2368:     test:
2369:       suffix: sw_hll
2370:       args: -ufv_vtk_interval 0 -bc_wall 1,2,3,4 -physics sw -ufv_cfl 3 -petscfv_type leastsquares -petsclimiter_type sin -ts_max_steps 5 -ts_ssp_type rks2 -ts_ssp_nstages 10 -monitor height,energy -grid_bounds 0,5,0,5 -dm_plex_box_faces 25,25 -sw_riemann hll

2372:     # 2D Advection: p4est
2373:     test:
2374:       suffix: p4est_advec_2d
2375:       requires: p4est
2376:       args: -ufv_vtk_interval 0 -dm_type p4est -dm_forest_minimum_refinement 1 -dm_forest_initial_refinement 2 -dm_p4est_refine_pattern hash   -dm_forest_maximum_refinement 5

2378:     # Advection in a box
2379:     test:
2380:       suffix: adv_2d_quad_0
2381:       args: -ufv_vtk_interval 0 -dm_refine 3 -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3

2383:     test:
2384:       suffix: adv_2d_quad_1
2385:       args: -ufv_vtk_interval 0 -dm_refine 3 -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1
2386:       timeoutfactor: 3

2388:     test:
2389:       suffix: adv_2d_quad_p4est_0
2390:       requires: p4est
2391:       args: -ufv_vtk_interval 0 -dm_refine 5 -dm_type p4est -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3

2393:     test:
2394:       suffix: adv_2d_quad_p4est_1
2395:       requires: p4est
2396:       args: -ufv_vtk_interval 0 -dm_refine 5 -dm_type p4est -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow   3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1
2397:       timeoutfactor: 3

2399:     test:
2400:       suffix: adv_2d_quad_p4est_adapt_0
2401:       requires: p4est !__float128 #broken for quad precision
2402:       args: -ufv_vtk_interval 0 -dm_refine 3 -dm_type p4est -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow   3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1 -ufv_use_amr -refine_vec_tagger_box 0.005,inf -coarsen_vec_tagger_box   0,1.e-5 -petscfv_type leastsquares -ts_max_time 0.01
2403:       timeoutfactor: 3

2405:     test:
2406:       suffix: adv_0
2407:       requires: exodusii
2408:       args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -bc_inflow 100,101,200 -bc_outflow 201

2410:     test:
2411:       suffix: shock_0
2412:       requires: p4est !single !complex
2413:       args: -dm_plex_box_faces 2,1 -grid_bounds -1,1.,0.,1 -grid_skew_60 \
2414:       -dm_type p4est -dm_forest_partition_overlap 1 -dm_forest_maximum_refinement 6 -dm_forest_minimum_refinement 2 -dm_forest_initial_refinement 2 \
2415:       -ufv_use_amr -refine_vec_tagger_box 0.5,inf -coarsen_vec_tagger_box 0,1.e-2 -refine_tag_view -coarsen_tag_view \
2416:       -bc_wall 1,2,3,4 -physics euler -eu_type iv_shock -ufv_cfl 10 -eu_alpha 60. -eu_gamma 1.4 -eu_amach 2.02 -eu_rho2 3. \
2417:       -petscfv_type leastsquares -petsclimiter_type minmod -petscfv_compute_gradients 0 \
2418:       -ts_max_time 0.5 -ts_ssp_type rks2 -ts_ssp_nstages 10 \
2419:       -ufv_vtk_basename ${wPETSC_DIR}/ex11 -ufv_vtk_interval 0 -monitor density,energy
2420:       timeoutfactor: 3

2422:     # Test GLVis visualization of PetscFV fields
2423:     test:
2424:       suffix: glvis_adv_2d_tet
2425:       args: -ufv_vtk_interval 0 -ufv_vtk_monitor 0 \
2426:             -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_gmsh_periodic 0 \
2427:             -ts_monitor_solution glvis: -ts_max_steps 0

2429:     test:
2430:       suffix: glvis_adv_2d_quad
2431:       args: -ufv_vtk_interval 0 -ufv_vtk_monitor 0 -bc_inflow 1,2,4 -bc_outflow 3 \
2432:             -dm_refine 5 -dm_plex_separate_marker \
2433:             -ts_monitor_solution glvis: -ts_max_steps 0

2435: TEST*/