Actual source code: ex11.c

petsc-main 2021-04-20
Report Typos and Errors
  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>
 41: #include <petscsf.h>

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

 46: static PetscFunctionList PhysicsList, PhysicsRiemannList_SW;

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

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

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

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

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

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

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

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

116: PETSC_STATIC_INLINE PetscReal DotDIMReal(const PetscReal *x,const PetscReal *y)
117: {
118:   PetscInt  i;
119:   PetscReal prod=0.0;

121:   for (i=0; i<DIM; i++) prod += x[i]*y[i];
122:   return prod;
123: }
124: PETSC_STATIC_INLINE PetscReal NormDIM(const PetscReal *x) { return PetscSqrtReal(PetscAbsReal(DotDIMReal(x,x))); }

126: PETSC_STATIC_INLINE PetscReal Dot2Real(const PetscReal *x,const PetscReal *y) { return x[0]*y[0] + x[1]*y[1];}
127: PETSC_STATIC_INLINE PetscReal Norm2Real(const PetscReal *x) { return PetscSqrtReal(PetscAbsReal(Dot2Real(x,x)));}
128: PETSC_STATIC_INLINE void Normalize2Real(PetscReal *x) { PetscReal a = 1./Norm2Real(x); x[0] *= a; x[1] *= a; }
129: PETSC_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]; }
130: PETSC_STATIC_INLINE void Scale2Real(PetscReal a,const PetscReal *x,PetscReal *y) { y[0] = a*x[0]; y[1] = a*x[1]; }

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

138: typedef struct {
139:   PetscReal wind[DIM];
140: } Physics_Advect_Tilted;
141: typedef struct {
142:   PetscReal         center[DIM];
143:   PetscReal         radius;
144:   AdvectSolBumpType type;
145: } Physics_Advect_Bump;

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

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

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

168:   xG[0] = advect->inflowState;
169:   return(0);
170: }

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

179: 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)
180: {
181:   Physics_Advect *advect = (Physics_Advect*)phys->data;
182:   PetscReal      wind[DIM],wn;

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

199:       rad2 = 0.;
200:       for (i = 0; i < dim; i++) {
201:         comp2[i] = qp[i] * qp[i];
202:         rad2    += comp2[i];
203:       }

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

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

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

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

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

278:   PhysicsSolution_Advect(mod,time,x,yexact,phys);
279:   f[advect->functional.Solution] = PetscRealPart(y[0]);
280:   f[advect->functional.Error] = PetscAbsScalar(y[0]-yexact[0]);
281:   return(0);
282: }

284: static PetscErrorCode SetUpBC_Advect(DM dm, PetscDS prob, Physics phys)
285: {
287:   const PetscInt inflowids[] = {100,200,300},outflowids[] = {101};
288:   DMLabel        label;

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

298: static PetscErrorCode PhysicsCreate_Advect(Model mod,Physics phys,PetscOptionItems *PetscOptionsObject)
299: {
300:   Physics_Advect *advect;

304:   phys->field_desc = PhysicsFields_Advect;
305:   phys->riemann    = (PetscRiemannFunc)PhysicsRiemann_Advect;
306:   PetscNew(&advect);
307:   phys->data       = advect;
308:   mod->setupbc = SetUpBC_Advect;

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

351: /******************* Shallow Water ********************/
352: typedef struct {
353:   PetscReal gravity;
354:   PetscReal boundaryHeight;
355:   struct {
356:     PetscInt Height;
357:     PetscInt Speed;
358:     PetscInt Energy;
359:   } functional;
360: } Physics_SW;
361: typedef struct {
362:   PetscReal h;
363:   PetscReal uh[DIM];
364: } SWNode;
365: typedef union {
366:   SWNode    swnode;
367:   PetscReal vals[DIM+1];
368: } SWNodeUnion;

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

372: /*
373:  * h_t + div(uh) = 0
374:  * (uh)_t + div (u\otimes uh + g h^2 / 2 I) = 0
375:  *
376:  * */
377: static PetscErrorCode SWFlux(Physics phys,const PetscReal *n,const SWNode *x,SWNode *f)
378: {
379:   Physics_SW  *sw = (Physics_SW*)phys->data;
380:   PetscReal   uhn,u[DIM];
381:   PetscInt     i;

384:   Scale2Real(1./x->h,x->uh,u);
385:   uhn  = x->uh[0] * n[0] + x->uh[1] * n[1];
386:   f->h = uhn;
387:   for (i=0; i<DIM; i++) f->uh[i] = u[i] * uhn + sw->gravity * PetscSqr(x->h) * n[i];
388:   return(0);
389: }

391: static PetscErrorCode PhysicsBoundary_SW_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
392: {
394:   xG[0] = xI[0];
395:   xG[1] = -xI[1];
396:   xG[2] = -xI[2];
397:   return(0);
398: }

400: 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)
401: {
402:   Physics_SW *sw = (Physics_SW *) phys->data;
403:   PetscReal aL, aR;
404:   PetscReal nn[DIM];
405: #if !defined(PETSC_USE_COMPLEX)
406:   const SWNode *uL = (const SWNode *) xL, *uR = (const SWNode *) xR;
407: #else
408:   SWNodeUnion  uLreal, uRreal;
409:   const SWNode *uL = &uLreal.swnode;
410:   const SWNode *uR = &uRreal.swnode;
411: #endif
412:   SWNodeUnion fL, fR;
413:   PetscInt i;
414:   PetscReal zero = 0.;

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

455: 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)
456: {
457:   Physics_SW   *sw = (Physics_SW*)phys->data;
458:   PetscReal    cL,cR,speed;
459:   PetscReal    nn[DIM];
460: #if !defined(PETSC_USE_COMPLEX)
461:   const SWNode *uL = (const SWNode*)xL,*uR = (const SWNode*)xR;
462: #else
463:   SWNodeUnion  uLreal, uRreal;
464:   const SWNode *uL = &uLreal.swnode;
465:   const SWNode *uR = &uRreal.swnode;
466: #endif
467:   SWNodeUnion  fL,fR;
468:   PetscInt     i;
469:   PetscReal    zero=0.;

471: #if defined(PETSC_USE_COMPLEX)
472:   uLreal.swnode.h = 0; uRreal.swnode.h = 0;
473:   for (i = 0; i < 1+dim; i++) uLreal.vals[i] = PetscRealPart(xL[i]);
474:   for (i = 0; i < 1+dim; i++) uRreal.vals[i] = PetscRealPart(xR[i]);
475: #endif
476:   if (uL->h < 0 || uR->h < 0) {for (i=0; i<1+dim; i++) flux[i] = zero/zero; return;} /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reconstructed thickness is negative"); */
477:   nn[0] = n[0];
478:   nn[1] = n[1];
479:   Normalize2Real(nn);
480:   SWFlux(phys,nn,uL,&(fL.swnode));
481:   SWFlux(phys,nn,uR,&(fR.swnode));
482:   cL    = PetscSqrtReal(sw->gravity*uL->h);
483:   cR    = PetscSqrtReal(sw->gravity*uR->h); /* gravity wave speed */
484:   speed = PetscMax(PetscAbsReal(Dot2Real(uL->uh,nn)/uL->h) + cL,PetscAbsReal(Dot2Real(uR->uh,nn)/uR->h) + cR);
485:   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);
486: }

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

493:   if (time != 0.0) SETERRQ1(mod->comm,PETSC_ERR_SUP,"No solution known for time %g",(double)time);
494:   dx[0] = x[0] - 1.5;
495:   dx[1] = x[1] - 1.0;
496:   r     = Norm2Real(dx);
497:   sigma = 0.5;
498:   u[0]  = 1 + 2*PetscExpReal(-PetscSqr(r)/(2*PetscSqr(sigma)));
499:   u[1]  = 0.0;
500:   u[2]  = 0.0;
501:   return(0);
502: }

504: static PetscErrorCode PhysicsFunctional_SW(Model mod,PetscReal time,const PetscReal *coord,const PetscScalar *xx,PetscReal *f,void *ctx)
505: {
506:   Physics      phys = (Physics)ctx;
507:   Physics_SW   *sw  = (Physics_SW*)phys->data;
508:   const SWNode *x   = (const SWNode*)xx;
509:   PetscReal  u[2];
510:   PetscReal    h;

513:   h = x->h;
514:   Scale2Real(1./x->h,x->uh,u);
515:   f[sw->functional.Height] = h;
516:   f[sw->functional.Speed]  = Norm2Real(u) + PetscSqrtReal(sw->gravity*h);
517:   f[sw->functional.Energy] = 0.5*(Dot2Real(x->uh,u) + sw->gravity*PetscSqr(h));
518:   return(0);
519: }

521: static PetscErrorCode SetUpBC_SW(DM dm, PetscDS prob,Physics phys)
522: {
524:   const PetscInt wallids[] = {100,101,200,300};
525:   DMLabel        label;

528:   DMGetLabel(dm, "Face Sets", &label);
529:   PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", label, ALEN(wallids), wallids, 0, 0, NULL, (void (*)(void)) PhysicsBoundary_SW_Wall, NULL, phys, NULL);
530:   return(0);
531: }

533: static PetscErrorCode PhysicsCreate_SW(Model mod,Physics phys,PetscOptionItems *PetscOptionsObject)
534: {
535:   Physics_SW     *sw;
536:   char           sw_riemann[64] = "rusanov";

540:   phys->field_desc = PhysicsFields_SW;
541:   PetscNew(&sw);
542:   phys->data    = sw;
543:   mod->setupbc  = SetUpBC_SW;

545:   PetscFunctionListAdd(&PhysicsRiemannList_SW, "rusanov", PhysicsRiemann_SW_Rusanov);
546:   PetscFunctionListAdd(&PhysicsRiemannList_SW, "hll", PhysicsRiemann_SW_HLL);

548:   PetscOptionsHead(PetscOptionsObject,"SW options");
549:   {
550:     void (*PhysicsRiemann_SW)(PetscInt, PetscInt, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, PetscInt, const PetscScalar, PetscScalar *, Physics);
551:     sw->gravity = 1.0;
552:     PetscOptionsReal("-sw_gravity","Gravitational constant","",sw->gravity,&sw->gravity,NULL);
553:     PetscOptionsFList("-sw_riemann","Riemann solver","",PhysicsRiemannList_SW,sw_riemann,sw_riemann,sizeof sw_riemann,NULL);
554:     PetscFunctionListFind(PhysicsRiemannList_SW,sw_riemann,&PhysicsRiemann_SW);
555:     phys->riemann = (PetscRiemannFunc) PhysicsRiemann_SW;
556:   }
557:   PetscOptionsTail();
558:   phys->maxspeed = PetscSqrtReal(2.0*sw->gravity); /* Mach 1 for depth of 2 */

560:   ModelSolutionSetDefault(mod,PhysicsSolution_SW,phys);
561:   ModelFunctionalRegister(mod,"Height",&sw->functional.Height,PhysicsFunctional_SW,phys);
562:   ModelFunctionalRegister(mod,"Speed",&sw->functional.Speed,PhysicsFunctional_SW,phys);
563:   ModelFunctionalRegister(mod,"Energy",&sw->functional.Energy,PhysicsFunctional_SW,phys);

565:   mod->bcs[0] = mod->bcs[1] = mod->bcs[2] = DM_BOUNDARY_GHOSTED;

567:   return(0);
568: }

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

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

601: /* initial condition */
602: int initLinearWave(EulerNode *ux, const PetscReal gamma, const PetscReal coord[], const PetscReal Lx);
603: static PetscErrorCode PhysicsSolution_Euler(Model mod, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx)
604: {
605:   PetscInt i;
606:   Physics         phys = (Physics)ctx;
607:   Physics_Euler   *eu  = (Physics_Euler*)phys->data;
608:   EulerNode       *uu  = (EulerNode*)u;
609:   PetscReal        p0,gamma,c;
611:   if (time != 0.0) SETERRQ1(mod->comm,PETSC_ERR_SUP,"No solution known for time %g",(double)time);

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

617:   if (eu->type==EULER_IV_SHOCK || eu->type==EULER_SS_SHOCK) {
618:     /******************* Euler Density Shock ********************/
619:     /* On initial-value and self-similar solutions of the compressible Euler equations */
620:     /* Ravi Samtaney and D. I. Pullin */
621:     /* Phys. Fluids 8, 2650 (1996); http://dx.doi.org/10.1063/1.869050 */
622:     /* initial conditions 1: left of shock, 0: left of discontinuity 2: right of discontinuity,  */
623:     p0 = 1.;
624:     if (x[0] < 0.0 + x[1]*eu->pars[EULER_PAR_ITANA]) {
625:       if (x[0] < mod->bounds[0]*0.5) { /* left of shock (1) */
626:         PetscReal amach,rho,press,gas1,p1;
627:         amach = eu->pars[EULER_PAR_AMACH];
628:         rho = 1.;
629:         press = p0;
630:         p1 = press*(1.0+2.0*gamma/(gamma+1.0)*(amach*amach-1.0));
631:         gas1 = (gamma-1.0)/(gamma+1.0);
632:         uu->r = rho*(p1/press+gas1)/(gas1*p1/press+1.0);
633:         uu->ru[0]   = ((uu->r - rho)*PetscSqrtReal(gamma*press/rho)*amach);
634:         uu->E = p1/(gamma-1.0) + .5/uu->r*uu->ru[0]*uu->ru[0];
635:       }
636:       else { /* left of discontinuity (0) */
637:         uu->r = 1.; /* rho = 1 */
638:         uu->E = p0/(gamma-1.0);
639:       }
640:     }
641:     else { /* right of discontinuity (2) */
642:       uu->r = eu->pars[EULER_PAR_RHOR];
643:       uu->E = p0/(gamma-1.0);
644:     }
645:   }
646:   else if (eu->type==EULER_SHOCK_TUBE) {
647:     /* 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. */
648:     if (x[0] < 0.0) {
649:       uu->r = 8.;
650:       uu->E = 10./(gamma-1.);
651:     }
652:     else {
653:       uu->r = 1.;
654:       uu->E = 1./(gamma-1.);
655:     }
656:   }
657:   else if (eu->type==EULER_LINEAR_WAVE) {
658:     initLinearWave( uu, gamma, x, mod->bounds[1] - mod->bounds[0]);
659:   }
660:   else SETERRQ1(mod->comm,PETSC_ERR_SUP,"Unknown type %d",eu->type);

662:   /* set phys->maxspeed: (mod->maxspeed = phys->maxspeed) in main; */
663:   eu->sound(&gamma,uu,&c);
664:   c = (uu->ru[0]/uu->r) + c;
665:   if (c > phys->maxspeed) phys->maxspeed = c;

667:   return(0);
668: }

670: static PetscErrorCode Pressure_PG(const PetscReal gamma,const EulerNode *x,PetscReal *p)
671: {
672:   PetscReal ru2;

675:   ru2  = DotDIMReal(x->ru,x->ru);
676:   (*p)=(x->E - 0.5*ru2/x->r)*(gamma - 1.0); /* (E - rho V^2/2)(gamma-1) = e rho (gamma-1) */
677:   return(0);
678: }

680: static PetscErrorCode SpeedOfSound_PG(const PetscReal *gamma, const EulerNode *x, PetscReal *c)
681: {
682:   PetscReal p;

685:   Pressure_PG(*gamma,x,&p);
686:   if (p<0.) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"negative pressure time %g -- NEED TO FIX!!!!!!",(double) p);
687:   /* pars[EULER_PAR_GAMMA] = heat capacity ratio */
688:   (*c)=PetscSqrtReal(*gamma * p / x->r);
689:   return(0);
690: }

692: /*
693:  * x = (rho,rho*(u_1),...,rho*e)^T
694:  * x_t+div(f_1(x))+...+div(f_DIM(x)) = 0
695:  *
696:  * f_i(x) = u_i*x+(0,0,...,p,...,p*u_i)^T
697:  *
698:  */
699: static PetscErrorCode EulerFlux(Physics phys,const PetscReal *n,const EulerNode *x,EulerNode *f)
700: {
701:   Physics_Euler *eu = (Physics_Euler*)phys->data;
702:   PetscReal     nu,p;
703:   PetscInt      i;

706:   Pressure_PG(eu->pars[EULER_PAR_GAMMA],x,&p);
707:   nu = DotDIMReal(x->ru,n);
708:   f->r = nu;   /* A rho u */
709:   nu /= x->r;  /* A u */
710:   for (i=0; i<DIM; i++) f->ru[i] = nu * x->ru[i] + n[i]*p;  /* r u^2 + p */
711:   f->E = nu * (x->E + p); /* u(e+p) */
712:   return(0);
713: }

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

751:   for (i=0,s2=0.; i<DIM; i++) {
752:     nn[i] = n[i];
753:     s2 += nn[i]*nn[i];
754:   }
755:   s2 = PetscSqrtReal(s2); /* |n|_2 = sum(n^2)^1/2 */
756:   for (i=0.; i<DIM; i++) nn[i] /= s2;
757:   if (0) { /* Rusanov */
758:     const EulerNode *uL = (const EulerNode*)xL,*uR = (const EulerNode*)xR;
759:     EulerNodeUnion  fL,fR;
760:     EulerFlux(phys,nn,uL,&(fL.eulernode));
761:     EulerFlux(phys,nn,uR,&(fR.eulernode));
762:     eu->sound(&eu->pars[EULER_PAR_GAMMA],uL,&cL);if (ierr) exit(13);
763:     eu->sound(&eu->pars[EULER_PAR_GAMMA],uR,&cR);if (ierr) exit(14);
764:     velL = DotDIMReal(uL->ru,nn)/uL->r;
765:     velR = DotDIMReal(uR->ru,nn)/uR->r;
766:     speed = PetscMax(velR + cR, velL + cL);
767:     for (i=0; i<2+dim; i++) flux[i] = 0.5*((fL.vals[i]+fR.vals[i]) + speed*(xL[i] - xR[i]))*s2;
768:   }
769:   else {
770:     int dim = DIM;
771:     /* int iwave =  */
772:     godunovflux(xL, xR, flux, nn, &dim, &eu->pars[EULER_PAR_GAMMA]);
773:     for (i=0; i<2+dim; i++) flux[i] *= s2;
774:   }
775:   PetscFunctionReturnVoid();
776: }

778: static PetscErrorCode PhysicsFunctional_Euler(Model mod,PetscReal time,const PetscReal *coord,const PetscScalar *xx,PetscReal *f,void *ctx)
779: {
780:   Physics         phys = (Physics)ctx;
781:   Physics_Euler   *eu  = (Physics_Euler*)phys->data;
782:   const EulerNode *x   = (const EulerNode*)xx;
783:   PetscReal       p;

786:   f[eu->monitor.Density]  = x->r;
787:   f[eu->monitor.Momentum] = NormDIM(x->ru);
788:   f[eu->monitor.Energy]   = x->E;
789:   f[eu->monitor.Speed]    = NormDIM(x->ru)/x->r;
790:   Pressure_PG(eu->pars[EULER_PAR_GAMMA], x, &p);
791:   f[eu->monitor.Pressure] = p;
792:   return(0);
793: }

795: static PetscErrorCode SetUpBC_Euler(DM dm, PetscDS prob,Physics phys)
796: {
797:   PetscErrorCode  ierr;
798:   Physics_Euler   *eu = (Physics_Euler *) phys->data;
799:   DMLabel         label;

801:   DMGetLabel(dm, "Face Sets", &label);
802:   if (eu->type == EULER_LINEAR_WAVE) {
803:     const PetscInt wallids[] = {100,101};
804:     PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", label, ALEN(wallids), wallids, 0, 0, NULL, (void (*)(void)) PhysicsBoundary_Euler_Wall, NULL, phys, NULL);
805:   }
806:   else {
807:     const PetscInt wallids[] = {100,101,200,300};
808:     PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", label, ALEN(wallids), wallids, 0, 0, NULL, (void (*)(void)) PhysicsBoundary_Euler_Wall, NULL, phys, NULL);
809:   }
810:   return(0);
811: }

813: static PetscErrorCode PhysicsCreate_Euler(Model mod,Physics phys,PetscOptionItems *PetscOptionsObject)
814: {
815:   Physics_Euler   *eu;
816:   PetscErrorCode  ierr;

819:   phys->field_desc = PhysicsFields_Euler;
820:   phys->riemann = (PetscRiemannFunc) PhysicsRiemann_Euler_Godunov;
821:   PetscNew(&eu);
822:   phys->data    = eu;
823:   mod->setupbc = SetUpBC_Euler;
824:   PetscOptionsHead(PetscOptionsObject,"Euler options");
825:   {
826:     PetscReal alpha;
827:     char type[64] = "linear_wave";
828:     PetscBool  is;
829:     mod->bcs[0] = mod->bcs[1] = mod->bcs[2] = DM_BOUNDARY_GHOSTED;
830:     eu->pars[EULER_PAR_GAMMA] = 1.4;
831:     eu->pars[EULER_PAR_AMACH] = 2.02;
832:     eu->pars[EULER_PAR_RHOR] = 3.0;
833:     eu->pars[EULER_PAR_ITANA] = 0.57735026918963; /* angle of Euler self similar (SS) shock */
834:     PetscOptionsReal("-eu_gamma","Heat capacity ratio","",eu->pars[EULER_PAR_GAMMA],&eu->pars[EULER_PAR_GAMMA],NULL);
835:     PetscOptionsReal("-eu_amach","Shock speed (Mach)","",eu->pars[EULER_PAR_AMACH],&eu->pars[EULER_PAR_AMACH],NULL);
836:     PetscOptionsReal("-eu_rho2","Density right of discontinuity","",eu->pars[EULER_PAR_RHOR],&eu->pars[EULER_PAR_RHOR],NULL);
837:     alpha = 60.;
838:     PetscOptionsReal("-eu_alpha","Angle of discontinuity","",alpha,&alpha,NULL);
839:     if (alpha<=0. || alpha>90.) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Alpha bust be > 0 and <= 90 (%g)",alpha);
840:     eu->pars[EULER_PAR_ITANA] = 1./PetscTanReal( alpha * PETSC_PI / 180.0);
841:     PetscOptionsString("-eu_type","Type of Euler test","",type,type,sizeof(type),NULL);
842:     PetscStrcmp(type,"linear_wave", &is);
843:     if (is) {
844:       eu->type = EULER_LINEAR_WAVE;
845:       mod->bcs[0] = mod->bcs[1] = mod->bcs[2] = DM_BOUNDARY_PERIODIC;
846:       mod->bcs[1] = DM_BOUNDARY_GHOSTED; /* debug */
847:       PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"linear_wave");
848:     }
849:     else {
850:       if (DIM != 2) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DIM must be 2 unless linear wave test %s",type);
851:       PetscStrcmp(type,"iv_shock", &is);
852:       if (is) {
853:         eu->type = EULER_IV_SHOCK;
854:         PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"iv_shock");
855:       }
856:       else {
857:         PetscStrcmp(type,"ss_shock", &is);
858:         if (is) {
859:           eu->type = EULER_SS_SHOCK;
860:           PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"ss_shock");
861:         }
862:         else {
863:           PetscStrcmp(type,"shock_tube", &is);
864:           if (is) eu->type = EULER_SHOCK_TUBE;
865:           else SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Unknown Euler type %s",type);
866:           PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"shock_tube");
867:         }
868:       }
869:     }
870:   }
871:   PetscOptionsTail();
872:   eu->sound = SpeedOfSound_PG;
873:   phys->maxspeed = 0.; /* will get set in solution */
874:   ModelSolutionSetDefault(mod,PhysicsSolution_Euler,phys);
875:   ModelFunctionalRegister(mod,"Speed",&eu->monitor.Speed,PhysicsFunctional_Euler,phys);
876:   ModelFunctionalRegister(mod,"Energy",&eu->monitor.Energy,PhysicsFunctional_Euler,phys);
877:   ModelFunctionalRegister(mod,"Density",&eu->monitor.Density,PhysicsFunctional_Euler,phys);
878:   ModelFunctionalRegister(mod,"Momentum",&eu->monitor.Momentum,PhysicsFunctional_Euler,phys);
879:   ModelFunctionalRegister(mod,"Pressure",&eu->monitor.Pressure,PhysicsFunctional_Euler,phys);

881:   return(0);
882: }

884: static PetscErrorCode ErrorIndicator_Simple(PetscInt dim, PetscReal volume, PetscInt numComps, const PetscScalar u[], const PetscScalar grad[], PetscReal *error, void *ctx)
885: {
886:   PetscReal      err = 0.;
887:   PetscInt       i, j;

890:   for (i = 0; i < numComps; i++) {
891:     for (j = 0; j < dim; j++) {
892:       err += PetscSqr(PetscRealPart(grad[i * dim + j]));
893:     }
894:   }
895:   *error = volume * err;
896:   return(0);
897: }

899: PetscErrorCode ConstructCellBoundary(DM dm, User user)
900: {
901:   const char     *name   = "Cell Sets";
902:   const char     *bdname = "split faces";
903:   IS             regionIS, innerIS;
904:   const PetscInt *regions, *cells;
905:   PetscInt       numRegions, innerRegion, numCells, c;
906:   PetscInt       cStart, cEnd, cEndInterior, fStart, fEnd;
907:   PetscBool      hasLabel;

911:   DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
912:   DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);
913:   DMPlexGetGhostCellStratum(dm, &cEndInterior, NULL);

915:   DMHasLabel(dm, name, &hasLabel);
916:   if (!hasLabel) return(0);
917:   DMGetLabelSize(dm, name, &numRegions);
918:   if (numRegions != 2) return(0);
919:   /* Get the inner id */
920:   DMGetLabelIdIS(dm, name, &regionIS);
921:   ISGetIndices(regionIS, &regions);
922:   innerRegion = regions[0];
923:   ISRestoreIndices(regionIS, &regions);
924:   ISDestroy(&regionIS);
925:   /* Find the faces between cells in different regions, could call DMPlexCreateNeighborCSR() */
926:   DMGetStratumIS(dm, name, innerRegion, &innerIS);
927:   ISGetLocalSize(innerIS, &numCells);
928:   ISGetIndices(innerIS, &cells);
929:   DMCreateLabel(dm, bdname);
930:   for (c = 0; c < numCells; ++c) {
931:     const PetscInt cell = cells[c];
932:     const PetscInt *faces;
933:     PetscInt       numFaces, f;

935:     if ((cell < cStart) || (cell >= cEnd)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Got invalid point %d which is not a cell", cell);
936:     DMPlexGetConeSize(dm, cell, &numFaces);
937:     DMPlexGetCone(dm, cell, &faces);
938:     for (f = 0; f < numFaces; ++f) {
939:       const PetscInt face = faces[f];
940:       const PetscInt *neighbors;
941:       PetscInt       nC, regionA, regionB;

943:       if ((face < fStart) || (face >= fEnd)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Got invalid point %d which is not a face", face);
944:       DMPlexGetSupportSize(dm, face, &nC);
945:       if (nC != 2) continue;
946:       DMPlexGetSupport(dm, face, &neighbors);
947:       if ((neighbors[0] >= cEndInterior) || (neighbors[1] >= cEndInterior)) continue;
948:       if ((neighbors[0] < cStart) || (neighbors[0] >= cEnd)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Got invalid point %d which is not a cell", neighbors[0]);
949:       if ((neighbors[1] < cStart) || (neighbors[1] >= cEnd)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Got invalid point %d which is not a cell", neighbors[1]);
950:       DMGetLabelValue(dm, name, neighbors[0], &regionA);
951:       DMGetLabelValue(dm, name, neighbors[1], &regionB);
952:       if (regionA < 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid label %s: Cell %d has no value", name, neighbors[0]);
953:       if (regionB < 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid label %s: Cell %d has no value", name, neighbors[1]);
954:       if (regionA != regionB) {
955:         DMSetLabelValue(dm, bdname, faces[f], 1);
956:       }
957:     }
958:   }
959:   ISRestoreIndices(innerIS, &cells);
960:   ISDestroy(&innerIS);
961:   {
962:     DMLabel label;

964:     DMGetLabel(dm, bdname, &label);
965:     DMLabelView(label, PETSC_VIEWER_STDOUT_WORLD);
966:   }
967:   return(0);
968: }

970: /* Right now, I have just added duplicate faces, which see both cells. We can
971: - Add duplicate vertices and decouple the face cones
972: - Disconnect faces from cells across the rotation gap
973: */
974: PetscErrorCode SplitFaces(DM *dmSplit, const char labelName[], User user)
975: {
976:   DM             dm = *dmSplit, sdm;
977:   PetscSF        sfPoint, gsfPoint;
978:   PetscSection   coordSection, newCoordSection;
979:   Vec            coordinates;
980:   IS             idIS;
981:   const PetscInt *ids;
982:   PetscInt       *newpoints;
983:   PetscInt       dim, depth, maxConeSize, maxSupportSize, numLabels, numGhostCells;
984:   PetscInt       numFS, fs, pStart, pEnd, p, cEnd, cEndInterior, vStart, vEnd, v, fStart, fEnd, newf, d, l;
985:   PetscBool      hasLabel;

989:   DMHasLabel(dm, labelName, &hasLabel);
990:   if (!hasLabel) return(0);
991:   DMCreate(PetscObjectComm((PetscObject)dm), &sdm);
992:   DMSetType(sdm, DMPLEX);
993:   DMGetDimension(dm, &dim);
994:   DMSetDimension(sdm, dim);

996:   DMGetLabelIdIS(dm, labelName, &idIS);
997:   ISGetLocalSize(idIS, &numFS);
998:   ISGetIndices(idIS, &ids);

1000:   user->numSplitFaces = 0;
1001:   for (fs = 0; fs < numFS; ++fs) {
1002:     PetscInt numBdFaces;

1004:     DMGetStratumSize(dm, labelName, ids[fs], &numBdFaces);
1005:     user->numSplitFaces += numBdFaces;
1006:   }
1007:   DMPlexGetChart(dm, &pStart, &pEnd);
1008:   pEnd += user->numSplitFaces;
1009:   DMPlexSetChart(sdm, pStart, pEnd);
1010:   DMPlexGetGhostCellStratum(dm, &cEndInterior, NULL);
1011:   DMPlexGetHeightStratum(dm, 0, NULL, &cEnd);
1012:   numGhostCells = cEnd - cEndInterior;
1013:   /* Set cone and support sizes */
1014:   DMPlexGetDepth(dm, &depth);
1015:   for (d = 0; d <= depth; ++d) {
1016:     DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);
1017:     for (p = pStart; p < pEnd; ++p) {
1018:       PetscInt newp = p;
1019:       PetscInt size;

1021:       DMPlexGetConeSize(dm, p, &size);
1022:       DMPlexSetConeSize(sdm, newp, size);
1023:       DMPlexGetSupportSize(dm, p, &size);
1024:       DMPlexSetSupportSize(sdm, newp, size);
1025:     }
1026:   }
1027:   DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);
1028:   for (fs = 0, newf = fEnd; fs < numFS; ++fs) {
1029:     IS             faceIS;
1030:     const PetscInt *faces;
1031:     PetscInt       numFaces, f;

1033:     DMGetStratumIS(dm, labelName, ids[fs], &faceIS);
1034:     ISGetLocalSize(faceIS, &numFaces);
1035:     ISGetIndices(faceIS, &faces);
1036:     for (f = 0; f < numFaces; ++f, ++newf) {
1037:       PetscInt size;

1039:       /* Right now I think that both faces should see both cells */
1040:       DMPlexGetConeSize(dm, faces[f], &size);
1041:       DMPlexSetConeSize(sdm, newf, size);
1042:       DMPlexGetSupportSize(dm, faces[f], &size);
1043:       DMPlexSetSupportSize(sdm, newf, size);
1044:     }
1045:     ISRestoreIndices(faceIS, &faces);
1046:     ISDestroy(&faceIS);
1047:   }
1048:   DMSetUp(sdm);
1049:   /* Set cones and supports */
1050:   DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);
1051:   PetscMalloc1(PetscMax(maxConeSize, maxSupportSize), &newpoints);
1052:   DMPlexGetChart(dm, &pStart, &pEnd);
1053:   for (p = pStart; p < pEnd; ++p) {
1054:     const PetscInt *points, *orientations;
1055:     PetscInt       size, i, newp = p;

1057:     DMPlexGetConeSize(dm, p, &size);
1058:     DMPlexGetCone(dm, p, &points);
1059:     DMPlexGetConeOrientation(dm, p, &orientations);
1060:     for (i = 0; i < size; ++i) newpoints[i] = points[i];
1061:     DMPlexSetCone(sdm, newp, newpoints);
1062:     DMPlexSetConeOrientation(sdm, newp, orientations);
1063:     DMPlexGetSupportSize(dm, p, &size);
1064:     DMPlexGetSupport(dm, p, &points);
1065:     for (i = 0; i < size; ++i) newpoints[i] = points[i];
1066:     DMPlexSetSupport(sdm, newp, newpoints);
1067:   }
1068:   PetscFree(newpoints);
1069:   for (fs = 0, newf = fEnd; fs < numFS; ++fs) {
1070:     IS             faceIS;
1071:     const PetscInt *faces;
1072:     PetscInt       numFaces, f;

1074:     DMGetStratumIS(dm, labelName, ids[fs], &faceIS);
1075:     ISGetLocalSize(faceIS, &numFaces);
1076:     ISGetIndices(faceIS, &faces);
1077:     for (f = 0; f < numFaces; ++f, ++newf) {
1078:       const PetscInt *points;

1080:       DMPlexGetCone(dm, faces[f], &points);
1081:       DMPlexSetCone(sdm, newf, points);
1082:       DMPlexGetSupport(dm, faces[f], &points);
1083:       DMPlexSetSupport(sdm, newf, points);
1084:     }
1085:     ISRestoreIndices(faceIS, &faces);
1086:     ISDestroy(&faceIS);
1087:   }
1088:   ISRestoreIndices(idIS, &ids);
1089:   ISDestroy(&idIS);
1090:   DMPlexStratify(sdm);
1091:   /* Convert coordinates */
1092:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
1093:   DMGetCoordinateSection(dm, &coordSection);
1094:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &newCoordSection);
1095:   PetscSectionSetNumFields(newCoordSection, 1);
1096:   PetscSectionSetFieldComponents(newCoordSection, 0, dim);
1097:   PetscSectionSetChart(newCoordSection, vStart, vEnd);
1098:   for (v = vStart; v < vEnd; ++v) {
1099:     PetscSectionSetDof(newCoordSection, v, dim);
1100:     PetscSectionSetFieldDof(newCoordSection, v, 0, dim);
1101:   }
1102:   PetscSectionSetUp(newCoordSection);
1103:   DMSetCoordinateSection(sdm, PETSC_DETERMINE, newCoordSection);
1104:   PetscSectionDestroy(&newCoordSection); /* relinquish our reference */
1105:   DMGetCoordinatesLocal(dm, &coordinates);
1106:   DMSetCoordinatesLocal(sdm, coordinates);
1107:   /* Convert labels */
1108:   DMGetNumLabels(dm, &numLabels);
1109:   for (l = 0; l < numLabels; ++l) {
1110:     const char *lname;
1111:     PetscBool  isDepth, isDim;

1113:     DMGetLabelName(dm, l, &lname);
1114:     PetscStrcmp(lname, "depth", &isDepth);
1115:     if (isDepth) continue;
1116:     PetscStrcmp(lname, "dim", &isDim);
1117:     if (isDim) continue;
1118:     DMCreateLabel(sdm, lname);
1119:     DMGetLabelIdIS(dm, lname, &idIS);
1120:     ISGetLocalSize(idIS, &numFS);
1121:     ISGetIndices(idIS, &ids);
1122:     for (fs = 0; fs < numFS; ++fs) {
1123:       IS             pointIS;
1124:       const PetscInt *points;
1125:       PetscInt       numPoints;

1127:       DMGetStratumIS(dm, lname, ids[fs], &pointIS);
1128:       ISGetLocalSize(pointIS, &numPoints);
1129:       ISGetIndices(pointIS, &points);
1130:       for (p = 0; p < numPoints; ++p) {
1131:         PetscInt newpoint = points[p];

1133:         DMSetLabelValue(sdm, lname, newpoint, ids[fs]);
1134:       }
1135:       ISRestoreIndices(pointIS, &points);
1136:       ISDestroy(&pointIS);
1137:     }
1138:     ISRestoreIndices(idIS, &ids);
1139:     ISDestroy(&idIS);
1140:   }
1141:   {
1142:     /* Convert pointSF */
1143:     const PetscSFNode *remotePoints;
1144:     PetscSFNode       *gremotePoints;
1145:     const PetscInt    *localPoints;
1146:     PetscInt          *glocalPoints,*newLocation,*newRemoteLocation;
1147:     PetscInt          numRoots, numLeaves;
1148:     PetscMPIInt       size;

1150:     MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size);
1151:     DMGetPointSF(dm, &sfPoint);
1152:     DMGetPointSF(sdm, &gsfPoint);
1153:     DMPlexGetChart(dm,&pStart,&pEnd);
1154:     PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);
1155:     if (numRoots >= 0) {
1156:       PetscMalloc2(numRoots,&newLocation,pEnd-pStart,&newRemoteLocation);
1157:       for (l=0; l<numRoots; l++) newLocation[l] = l; /* + (l >= cEnd ? numGhostCells : 0); */
1158:       PetscSFBcastBegin(sfPoint, MPIU_INT, newLocation, newRemoteLocation,MPI_REPLACE);
1159:       PetscSFBcastEnd(sfPoint, MPIU_INT, newLocation, newRemoteLocation,MPI_REPLACE);
1160:       PetscMalloc1(numLeaves,    &glocalPoints);
1161:       PetscMalloc1(numLeaves, &gremotePoints);
1162:       for (l = 0; l < numLeaves; ++l) {
1163:         glocalPoints[l]        = localPoints[l]; /* localPoints[l] >= cEnd ? localPoints[l] + numGhostCells : localPoints[l]; */
1164:         gremotePoints[l].rank  = remotePoints[l].rank;
1165:         gremotePoints[l].index = newRemoteLocation[localPoints[l]];
1166:       }
1167:       PetscFree2(newLocation,newRemoteLocation);
1168:       PetscSFSetGraph(gsfPoint, numRoots+numGhostCells, numLeaves, glocalPoints, PETSC_OWN_POINTER, gremotePoints, PETSC_OWN_POINTER);
1169:     }
1170:     DMDestroy(dmSplit);
1171:     *dmSplit = sdm;
1172:   }
1173:   return(0);
1174: }

1176: PetscErrorCode CreatePartitionVec(DM dm, DM *dmCell, Vec *partition)
1177: {
1178:   PetscSF        sfPoint;
1179:   PetscSection   coordSection;
1180:   Vec            coordinates;
1181:   PetscSection   sectionCell;
1182:   PetscScalar    *part;
1183:   PetscInt       cStart, cEnd, c;
1184:   PetscMPIInt    rank;

1188:   DMGetCoordinateSection(dm, &coordSection);
1189:   DMGetCoordinatesLocal(dm, &coordinates);
1190:   DMClone(dm, dmCell);
1191:   DMGetPointSF(dm, &sfPoint);
1192:   DMSetPointSF(*dmCell, sfPoint);
1193:   DMSetCoordinateSection(*dmCell, PETSC_DETERMINE, coordSection);
1194:   DMSetCoordinatesLocal(*dmCell, coordinates);
1195:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
1196:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &sectionCell);
1197:   DMPlexGetHeightStratum(*dmCell, 0, &cStart, &cEnd);
1198:   PetscSectionSetChart(sectionCell, cStart, cEnd);
1199:   for (c = cStart; c < cEnd; ++c) {
1200:     PetscSectionSetDof(sectionCell, c, 1);
1201:   }
1202:   PetscSectionSetUp(sectionCell);
1203:   DMSetLocalSection(*dmCell, sectionCell);
1204:   PetscSectionDestroy(&sectionCell);
1205:   DMCreateLocalVector(*dmCell, partition);
1206:   PetscObjectSetName((PetscObject)*partition, "partition");
1207:   VecGetArray(*partition, &part);
1208:   for (c = cStart; c < cEnd; ++c) {
1209:     PetscScalar *p;

1211:     DMPlexPointLocalRef(*dmCell, c, part, &p);
1212:     p[0] = rank;
1213:   }
1214:   VecRestoreArray(*partition, &part);
1215:   return(0);
1216: }

1218: PetscErrorCode CreateMassMatrix(DM dm, Vec *massMatrix, User user)
1219: {
1220:   DM                plex, dmMass, dmFace, dmCell, dmCoord;
1221:   PetscSection      coordSection;
1222:   Vec               coordinates, facegeom, cellgeom;
1223:   PetscSection      sectionMass;
1224:   PetscScalar       *m;
1225:   const PetscScalar *fgeom, *cgeom, *coords;
1226:   PetscInt          vStart, vEnd, v;
1227:   PetscErrorCode    ierr;

1230:   DMConvert(dm, DMPLEX, &plex);
1231:   DMGetCoordinateSection(dm, &coordSection);
1232:   DMGetCoordinatesLocal(dm, &coordinates);
1233:   DMClone(dm, &dmMass);
1234:   DMSetCoordinateSection(dmMass, PETSC_DETERMINE, coordSection);
1235:   DMSetCoordinatesLocal(dmMass, coordinates);
1236:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &sectionMass);
1237:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
1238:   PetscSectionSetChart(sectionMass, vStart, vEnd);
1239:   for (v = vStart; v < vEnd; ++v) {
1240:     PetscInt numFaces;

1242:     DMPlexGetSupportSize(dmMass, v, &numFaces);
1243:     PetscSectionSetDof(sectionMass, v, numFaces*numFaces);
1244:   }
1245:   PetscSectionSetUp(sectionMass);
1246:   DMSetLocalSection(dmMass, sectionMass);
1247:   PetscSectionDestroy(&sectionMass);
1248:   DMGetLocalVector(dmMass, massMatrix);
1249:   VecGetArray(*massMatrix, &m);
1250:   DMPlexGetGeometryFVM(plex, &facegeom, &cellgeom, NULL);
1251:   VecGetDM(facegeom, &dmFace);
1252:   VecGetArrayRead(facegeom, &fgeom);
1253:   VecGetDM(cellgeom, &dmCell);
1254:   VecGetArrayRead(cellgeom, &cgeom);
1255:   DMGetCoordinateDM(dm, &dmCoord);
1256:   VecGetArrayRead(coordinates, &coords);
1257:   for (v = vStart; v < vEnd; ++v) {
1258:     const PetscInt        *faces;
1259:     PetscFVFaceGeom       *fgA, *fgB, *cg;
1260:     PetscScalar           *vertex;
1261:     PetscInt               numFaces, sides[2], f, g;

1263:     DMPlexPointLocalRead(dmCoord, v, coords, &vertex);
1264:     DMPlexGetSupportSize(dmMass, v, &numFaces);
1265:     DMPlexGetSupport(dmMass, v, &faces);
1266:     for (f = 0; f < numFaces; ++f) {
1267:       sides[0] = faces[f];
1268:       DMPlexPointLocalRead(dmFace, faces[f], fgeom, &fgA);
1269:       for (g = 0; g < numFaces; ++g) {
1270:         const PetscInt *cells = NULL;
1271:         PetscReal      area   = 0.0;
1272:         PetscInt       numCells;

1274:         sides[1] = faces[g];
1275:         DMPlexPointLocalRead(dmFace, faces[g], fgeom, &fgB);
1276:         DMPlexGetJoin(dmMass, 2, sides, &numCells, &cells);
1277:         if (numCells != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Invalid join for faces");
1278:         DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cg);
1279:         area += PetscAbsScalar((vertex[0] - cg->centroid[0])*(fgA->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1])*(fgA->centroid[0] - cg->centroid[0]));
1280:         area += PetscAbsScalar((vertex[0] - cg->centroid[0])*(fgB->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1])*(fgB->centroid[0] - cg->centroid[0]));
1281:         m[f*numFaces+g] = Dot2Real(fgA->normal, fgB->normal)*area*0.5;
1282:         DMPlexRestoreJoin(dmMass, 2, sides, &numCells, &cells);
1283:       }
1284:     }
1285:   }
1286:   VecRestoreArrayRead(facegeom, &fgeom);
1287:   VecRestoreArrayRead(cellgeom, &cgeom);
1288:   VecRestoreArrayRead(coordinates, &coords);
1289:   VecRestoreArray(*massMatrix, &m);
1290:   DMDestroy(&dmMass);
1291:   DMDestroy(&plex);
1292:   return(0);
1293: }

1295: /* Behavior will be different for multi-physics or when using non-default boundary conditions */
1296: static PetscErrorCode ModelSolutionSetDefault(Model mod,SolutionFunction func,void *ctx)
1297: {
1299:   mod->solution    = func;
1300:   mod->solutionctx = ctx;
1301:   return(0);
1302: }

1304: static PetscErrorCode ModelFunctionalRegister(Model mod,const char *name,PetscInt *offset,FunctionalFunction func,void *ctx)
1305: {
1307:   FunctionalLink link,*ptr;
1308:   PetscInt       lastoffset = -1;

1311:   for (ptr=&mod->functionalRegistry; *ptr; ptr = &(*ptr)->next) lastoffset = (*ptr)->offset;
1312:   PetscNew(&link);
1313:   PetscStrallocpy(name,&link->name);
1314:   link->offset = lastoffset + 1;
1315:   link->func   = func;
1316:   link->ctx    = ctx;
1317:   link->next   = NULL;
1318:   *ptr         = link;
1319:   *offset      = link->offset;
1320:   return(0);
1321: }

1323: static PetscErrorCode ModelFunctionalSetFromOptions(Model mod,PetscOptionItems *PetscOptionsObject)
1324: {
1326:   PetscInt       i,j;
1327:   FunctionalLink link;
1328:   char           *names[256];

1331:   mod->numMonitored = ALEN(names);
1332:   PetscOptionsStringArray("-monitor","list of functionals to monitor","",names,&mod->numMonitored,NULL);
1333:   /* Create list of functionals that will be computed somehow */
1334:   PetscMalloc1(mod->numMonitored,&mod->functionalMonitored);
1335:   /* Create index of calls that we will have to make to compute these functionals (over-allocation in general). */
1336:   PetscMalloc1(mod->numMonitored,&mod->functionalCall);
1337:   mod->numCall = 0;
1338:   for (i=0; i<mod->numMonitored; i++) {
1339:     for (link=mod->functionalRegistry; link; link=link->next) {
1340:       PetscBool match;
1341:       PetscStrcasecmp(names[i],link->name,&match);
1342:       if (match) break;
1343:     }
1344:     if (!link) SETERRQ1(mod->comm,PETSC_ERR_USER,"No known functional '%s'",names[i]);
1345:     mod->functionalMonitored[i] = link;
1346:     for (j=0; j<i; j++) {
1347:       if (mod->functionalCall[j]->func == link->func && mod->functionalCall[j]->ctx == link->ctx) goto next_name;
1348:     }
1349:     mod->functionalCall[mod->numCall++] = link; /* Just points to the first link using the result. There may be more results. */
1350: next_name:
1351:     PetscFree(names[i]);
1352:   }

1354:   /* Find out the maximum index of any functional computed by a function we will be calling (even if we are not using it) */
1355:   mod->maxComputed = -1;
1356:   for (link=mod->functionalRegistry; link; link=link->next) {
1357:     for (i=0; i<mod->numCall; i++) {
1358:       FunctionalLink call = mod->functionalCall[i];
1359:       if (link->func == call->func && link->ctx == call->ctx) {
1360:         mod->maxComputed = PetscMax(mod->maxComputed,link->offset);
1361:       }
1362:     }
1363:   }
1364:   return(0);
1365: }

1367: static PetscErrorCode FunctionalLinkDestroy(FunctionalLink *link)
1368: {
1370:   FunctionalLink l,next;

1373:   if (!link) return(0);
1374:   l     = *link;
1375:   *link = NULL;
1376:   for (; l; l=next) {
1377:     next = l->next;
1378:     PetscFree(l->name);
1379:     PetscFree(l);
1380:   }
1381:   return(0);
1382: }

1384: /* put the solution callback into a functional callback */
1385: static PetscErrorCode SolutionFunctional(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *modctx)
1386: {
1387:   Model          mod;
1390:   mod  = (Model) modctx;
1391:   (*mod->solution)(mod, time, x, u, mod->solutionctx);
1392:   return(0);
1393: }

1395: PetscErrorCode SetInitialCondition(DM dm, Vec X, User user)
1396: {
1397:   PetscErrorCode     (*func[1]) (PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
1398:   void               *ctx[1];
1399:   Model              mod = user->model;
1400:   PetscErrorCode     ierr;

1403:   func[0] = SolutionFunctional;
1404:   ctx[0]  = (void *) mod;
1405:   DMProjectFunction(dm,0.0,func,ctx,INSERT_ALL_VALUES,X);
1406:   return(0);
1407: }

1409: static PetscErrorCode OutputVTK(DM dm, const char *filename, PetscViewer *viewer)
1410: {

1414:   PetscViewerCreate(PetscObjectComm((PetscObject)dm), viewer);
1415:   PetscViewerSetType(*viewer, PETSCVIEWERVTK);
1416:   PetscViewerFileSetName(*viewer, filename);
1417:   return(0);
1418: }

1420: static PetscErrorCode MonitorVTK(TS ts,PetscInt stepnum,PetscReal time,Vec X,void *ctx)
1421: {
1422:   User           user = (User)ctx;
1423:   DM             dm, plex;
1424:   PetscViewer    viewer;
1425:   char           filename[PETSC_MAX_PATH_LEN],*ftable = NULL;
1426:   PetscReal      xnorm;

1430:   PetscObjectSetName((PetscObject) X, "u");
1431:   VecGetDM(X,&dm);
1432:   VecNorm(X,NORM_INFINITY,&xnorm);

1434:   if (stepnum >= 0) {
1435:     stepnum += user->monitorStepOffset;
1436:   }
1437:   if (stepnum >= 0) {           /* No summary for final time */
1438:     Model             mod = user->model;
1439:     Vec               cellgeom;
1440:     PetscInt          c,cStart,cEnd,fcount,i;
1441:     size_t            ftableused,ftablealloc;
1442:     const PetscScalar *cgeom,*x;
1443:     DM                dmCell;
1444:     DMLabel           vtkLabel;
1445:     PetscReal         *fmin,*fmax,*fintegral,*ftmp;

1447:     DMConvert(dm, DMPLEX, &plex);
1448:     DMPlexGetGeometryFVM(plex, NULL, &cellgeom, NULL);
1449:     fcount = mod->maxComputed+1;
1450:     PetscMalloc4(fcount,&fmin,fcount,&fmax,fcount,&fintegral,fcount,&ftmp);
1451:     for (i=0; i<fcount; i++) {
1452:       fmin[i]      = PETSC_MAX_REAL;
1453:       fmax[i]      = PETSC_MIN_REAL;
1454:       fintegral[i] = 0;
1455:     }
1456:     VecGetDM(cellgeom,&dmCell);
1457:     DMPlexGetSimplexOrBoxCells(dmCell,0,&cStart,&cEnd);
1458:     VecGetArrayRead(cellgeom,&cgeom);
1459:     VecGetArrayRead(X,&x);
1460:     DMGetLabel(dm,"vtk",&vtkLabel);
1461:     for (c = cStart; c < cEnd; ++c) {
1462:       PetscFVCellGeom       *cg;
1463:       const PetscScalar     *cx    = NULL;
1464:       PetscInt              vtkVal = 0;

1466:       /* not that these two routines as currently implemented work for any dm with a
1467:        * localSection/globalSection */
1468:       DMPlexPointLocalRead(dmCell,c,cgeom,&cg);
1469:       DMPlexPointGlobalRead(dm,c,x,&cx);
1470:       if (vtkLabel) {DMLabelGetValue(vtkLabel,c,&vtkVal);}
1471:       if (!vtkVal || !cx) continue;        /* ghost, or not a global cell */
1472:       for (i=0; i<mod->numCall; i++) {
1473:         FunctionalLink flink = mod->functionalCall[i];
1474:         (*flink->func)(mod,time,cg->centroid,cx,ftmp,flink->ctx);
1475:       }
1476:       for (i=0; i<fcount; i++) {
1477:         fmin[i]       = PetscMin(fmin[i],ftmp[i]);
1478:         fmax[i]       = PetscMax(fmax[i],ftmp[i]);
1479:         fintegral[i] += cg->volume * ftmp[i];
1480:       }
1481:     }
1482:     VecRestoreArrayRead(cellgeom,&cgeom);
1483:     VecRestoreArrayRead(X,&x);
1484:     DMDestroy(&plex);
1485:     MPI_Allreduce(MPI_IN_PLACE,fmin,fcount,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));
1486:     MPI_Allreduce(MPI_IN_PLACE,fmax,fcount,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));
1487:     MPI_Allreduce(MPI_IN_PLACE,fintegral,fcount,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)ts));

1489:     ftablealloc = fcount * 100;
1490:     ftableused  = 0;
1491:     PetscMalloc1(ftablealloc,&ftable);
1492:     for (i=0; i<mod->numMonitored; i++) {
1493:       size_t         countused;
1494:       char           buffer[256],*p;
1495:       FunctionalLink flink = mod->functionalMonitored[i];
1496:       PetscInt       id    = flink->offset;
1497:       if (i % 3) {
1498:         PetscArraycpy(buffer,"  ",2);
1499:         p    = buffer + 2;
1500:       } else if (i) {
1501:         char newline[] = "\n";
1502:         PetscMemcpy(buffer,newline,sizeof(newline)-1);
1503:         p    = buffer + sizeof(newline) - 1;
1504:       } else {
1505:         p = buffer;
1506:       }
1507:       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]);
1508:       countused--;
1509:       countused += p - buffer;
1510:       if (countused > ftablealloc-ftableused-1) { /* reallocate */
1511:         char *ftablenew;
1512:         ftablealloc = 2*ftablealloc + countused;
1513:         PetscMalloc(ftablealloc,&ftablenew);
1514:         PetscArraycpy(ftablenew,ftable,ftableused);
1515:         PetscFree(ftable);
1516:         ftable = ftablenew;
1517:       }
1518:       PetscArraycpy(ftable+ftableused,buffer,countused);
1519:       ftableused += countused;
1520:       ftable[ftableused] = 0;
1521:     }
1522:     PetscFree4(fmin,fmax,fintegral,ftmp);

1524:     PetscPrintf(PetscObjectComm((PetscObject)ts),"% 3D  time %8.4g  |x| %8.4g  %s\n",stepnum,(double)time,(double)xnorm,ftable ? ftable : "");
1525:     PetscFree(ftable);
1526:   }
1527:   if (user->vtkInterval < 1) return(0);
1528:   if ((stepnum == -1) ^ (stepnum % user->vtkInterval == 0)) {
1529:     if (stepnum == -1) {        /* Final time is not multiple of normal time interval, write it anyway */
1530:       TSGetStepNumber(ts,&stepnum);
1531:     }
1532:     PetscSNPrintf(filename,sizeof filename,"%s-%03D.vtu",user->outputBasename,stepnum);
1533:     OutputVTK(dm,filename,&viewer);
1534:     VecView(X,viewer);
1535:     PetscViewerDestroy(&viewer);
1536:   }
1537:   return(0);
1538: }

1540: static PetscErrorCode initializeTS(DM dm, User user, TS *ts)
1541: {

1545:   TSCreate(PetscObjectComm((PetscObject)dm), ts);
1546:   TSSetType(*ts, TSSSP);
1547:   TSSetDM(*ts, dm);
1548:   if (user->vtkmon) {
1549:     TSMonitorSet(*ts,MonitorVTK,user,NULL);
1550:   }
1551:   DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, user);
1552:   DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, user);
1553:   TSSetMaxTime(*ts,2.0);
1554:   TSSetExactFinalTime(*ts,TS_EXACTFINALTIME_STEPOVER);
1555:   return(0);
1556: }

1558: static PetscErrorCode adaptToleranceFVM(PetscFV fvm, TS ts, Vec sol, VecTagger refineTag, VecTagger coarsenTag, User user, TS *tsNew, Vec *solNew)
1559: {
1560:   DM                dm, gradDM, plex, cellDM, adaptedDM = NULL;
1561:   Vec               cellGeom, faceGeom;
1562:   PetscBool         isForest, computeGradient;
1563:   Vec               grad, locGrad, locX, errVec;
1564:   PetscInt          cStart, cEnd, c, dim, nRefine, nCoarsen;
1565:   PetscReal         minMaxInd[2] = {PETSC_MAX_REAL, PETSC_MIN_REAL}, minMaxIndGlobal[2], minInd, maxInd, time;
1566:   PetscScalar       *errArray;
1567:   const PetscScalar *pointVals;
1568:   const PetscScalar *pointGrads;
1569:   const PetscScalar *pointGeom;
1570:   DMLabel           adaptLabel = NULL;
1571:   IS                refineIS, coarsenIS;
1572:   PetscErrorCode    ierr;

1575:   TSGetTime(ts,&time);
1576:   VecGetDM(sol, &dm);
1577:   DMGetDimension(dm,&dim);
1578:   PetscFVGetComputeGradients(fvm,&computeGradient);
1579:   PetscFVSetComputeGradients(fvm,PETSC_TRUE);
1580:   DMIsForest(dm, &isForest);
1581:   DMConvert(dm, DMPLEX, &plex);
1582:   DMPlexGetDataFVM(plex, fvm, &cellGeom, &faceGeom, &gradDM);
1583:   DMCreateLocalVector(plex,&locX);
1584:   DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, 0.0, faceGeom, cellGeom, NULL);
1585:   DMGlobalToLocalBegin(plex, sol, INSERT_VALUES, locX);
1586:   DMGlobalToLocalEnd  (plex, sol, INSERT_VALUES, locX);
1587:   DMCreateGlobalVector(gradDM, &grad);
1588:   DMPlexReconstructGradientsFVM(plex, locX, grad);
1589:   DMCreateLocalVector(gradDM, &locGrad);
1590:   DMGlobalToLocalBegin(gradDM, grad, INSERT_VALUES, locGrad);
1591:   DMGlobalToLocalEnd(gradDM, grad, INSERT_VALUES, locGrad);
1592:   VecDestroy(&grad);
1593:   DMPlexGetSimplexOrBoxCells(plex,0,&cStart,&cEnd);
1594:   VecGetArrayRead(locGrad,&pointGrads);
1595:   VecGetArrayRead(cellGeom,&pointGeom);
1596:   VecGetArrayRead(locX,&pointVals);
1597:   VecGetDM(cellGeom,&cellDM);
1598:   DMLabelCreate(PETSC_COMM_SELF,"adapt",&adaptLabel);
1599:   VecCreateMPI(PetscObjectComm((PetscObject)plex),cEnd-cStart,PETSC_DETERMINE,&errVec);
1600:   VecSetUp(errVec);
1601:   VecGetArray(errVec,&errArray);
1602:   for (c = cStart; c < cEnd; c++) {
1603:     PetscReal             errInd = 0.;
1604:     PetscScalar           *pointGrad;
1605:     PetscScalar           *pointVal;
1606:     PetscFVCellGeom       *cg;

1608:     DMPlexPointLocalRead(gradDM,c,pointGrads,&pointGrad);
1609:     DMPlexPointLocalRead(cellDM,c,pointGeom,&cg);
1610:     DMPlexPointLocalRead(plex,c,pointVals,&pointVal);

1612:     (user->model->errorIndicator)(dim,cg->volume,user->model->physics->dof,pointVal,pointGrad,&errInd,user->model->errorCtx);
1613:     errArray[c-cStart] = errInd;
1614:     minMaxInd[0] = PetscMin(minMaxInd[0],errInd);
1615:     minMaxInd[1] = PetscMax(minMaxInd[1],errInd);
1616:   }
1617:   VecRestoreArray(errVec,&errArray);
1618:   VecRestoreArrayRead(locX,&pointVals);
1619:   VecRestoreArrayRead(cellGeom,&pointGeom);
1620:   VecRestoreArrayRead(locGrad,&pointGrads);
1621:   VecDestroy(&locGrad);
1622:   VecDestroy(&locX);
1623:   DMDestroy(&plex);

1625:   VecTaggerComputeIS(refineTag,errVec,&refineIS);
1626:   VecTaggerComputeIS(coarsenTag,errVec,&coarsenIS);
1627:   ISGetSize(refineIS,&nRefine);
1628:   ISGetSize(coarsenIS,&nCoarsen);
1629:   if (nRefine) {DMLabelSetStratumIS(adaptLabel,DM_ADAPT_REFINE,refineIS);}
1630:   if (nCoarsen) {DMLabelSetStratumIS(adaptLabel,DM_ADAPT_COARSEN,coarsenIS);}
1631:   ISDestroy(&coarsenIS);
1632:   ISDestroy(&refineIS);
1633:   VecDestroy(&errVec);

1635:   PetscFVSetComputeGradients(fvm,computeGradient);
1636:   minMaxInd[1] = -minMaxInd[1];
1637:   MPI_Allreduce(minMaxInd,minMaxIndGlobal,2,MPIU_REAL,MPI_MIN,PetscObjectComm((PetscObject)dm));
1638:   minInd = minMaxIndGlobal[0];
1639:   maxInd = -minMaxIndGlobal[1];
1640:   PetscInfo2(ts, "error indicator range (%E, %E)\n", minInd, maxInd);
1641:   if (nRefine || nCoarsen) { /* at least one cell is over the refinement threshold */
1642:     DMAdaptLabel(dm,adaptLabel,&adaptedDM);
1643:   }
1644:   DMLabelDestroy(&adaptLabel);
1645:   if (adaptedDM) {
1646:     PetscInfo2(ts, "Adapted mesh, marking %D cells for refinement, and %D cells for coarsening\n", nRefine, nCoarsen);
1647:     if (tsNew) {initializeTS(adaptedDM, user, tsNew);}
1648:     if (solNew) {
1649:       DMCreateGlobalVector(adaptedDM, solNew);
1650:       PetscObjectSetName((PetscObject) *solNew, "solution");
1651:       DMForestTransferVec(dm, sol, adaptedDM, *solNew, PETSC_TRUE, time);
1652:     }
1653:     if (isForest) {DMForestSetAdaptivityForest(adaptedDM,NULL);} /* clear internal references to the previous dm */
1654:     DMDestroy(&adaptedDM);
1655:   } else {
1656:     if (tsNew)  *tsNew  = NULL;
1657:     if (solNew) *solNew = NULL;
1658:   }
1659:   return(0);
1660: }

1662: int main(int argc, char **argv)
1663: {
1664:   MPI_Comm          comm;
1665:   PetscDS           prob;
1666:   PetscFV           fvm;
1667:   PetscLimiter      limiter = NULL, noneLimiter = NULL;
1668:   User              user;
1669:   Model             mod;
1670:   Physics           phys;
1671:   DM                dm, plex;
1672:   PetscReal         ftime, cfl, dt, minRadius;
1673:   PetscInt          dim, nsteps;
1674:   TS                ts;
1675:   TSConvergedReason reason;
1676:   Vec               X;
1677:   PetscViewer       viewer;
1678:   PetscBool         simplex = PETSC_FALSE, vtkCellGeom, splitFaces, useAMR;
1679:   PetscInt          overlap, adaptInterval;
1680:   char              filename[PETSC_MAX_PATH_LEN] = "";
1681:   char              physname[256]  = "advect";
1682:   VecTagger         refineTag = NULL, coarsenTag = NULL;
1683:   PetscErrorCode    ierr;

1685:   PetscInitialize(&argc, &argv, (char*) 0, help);if (ierr) return ierr;
1686:   comm = PETSC_COMM_WORLD;

1688:   PetscNew(&user);
1689:   PetscNew(&user->model);
1690:   PetscNew(&user->model->physics);
1691:   mod           = user->model;
1692:   phys          = mod->physics;
1693:   mod->comm     = comm;
1694:   useAMR        = PETSC_FALSE;
1695:   adaptInterval = 1;

1697:   /* Register physical models to be available on the command line */
1698:   PetscFunctionListAdd(&PhysicsList,"advect"          ,PhysicsCreate_Advect);
1699:   PetscFunctionListAdd(&PhysicsList,"sw"              ,PhysicsCreate_SW);
1700:   PetscFunctionListAdd(&PhysicsList,"euler"           ,PhysicsCreate_Euler);


1703:   PetscOptionsBegin(comm,NULL,"Unstructured Finite Volume Mesh Options","");
1704:   {
1705:     cfl  = 0.9 * 4; /* default SSPRKS2 with s=5 stages is stable for CFL number s-1 */
1706:     PetscOptionsReal("-ufv_cfl","CFL number per step","",cfl,&cfl,NULL);
1707:     PetscOptionsString("-f","Exodus.II filename to read","",filename,filename,sizeof(filename),NULL);
1708:     PetscOptionsBool("-simplex","Flag to use a simplex mesh","",simplex,&simplex,NULL);
1709:     splitFaces = PETSC_FALSE;
1710:     PetscOptionsBool("-ufv_split_faces","Split faces between cell sets","",splitFaces,&splitFaces,NULL);
1711:     overlap = 1;
1712:     PetscOptionsInt("-ufv_mesh_overlap","Number of cells to overlap partitions","",overlap,&overlap,NULL);
1713:     user->vtkInterval = 1;
1714:     PetscOptionsInt("-ufv_vtk_interval","VTK output interval (0 to disable)","",user->vtkInterval,&user->vtkInterval,NULL);
1715:     user->vtkmon = PETSC_TRUE;
1716:     PetscOptionsBool("-ufv_vtk_monitor","Use VTKMonitor routine","",user->vtkmon,&user->vtkmon,NULL);
1717:     vtkCellGeom = PETSC_FALSE;
1718:     PetscStrcpy(user->outputBasename, "ex11");
1719:     PetscOptionsString("-ufv_vtk_basename","VTK output basename","",user->outputBasename,user->outputBasename,sizeof(user->outputBasename),NULL);
1720:     PetscOptionsBool("-ufv_vtk_cellgeom","Write cell geometry (for debugging)","",vtkCellGeom,&vtkCellGeom,NULL);
1721:     PetscOptionsBool("-ufv_use_amr","use local adaptive mesh refinement","",useAMR,&useAMR,NULL);
1722:     PetscOptionsInt("-ufv_adapt_interval","time steps between AMR","",adaptInterval,&adaptInterval,NULL);
1723:   }
1724:   PetscOptionsEnd();

1726:   if (useAMR) {
1727:     VecTaggerBox refineBox, coarsenBox;

1729:     refineBox.min  = refineBox.max  = PETSC_MAX_REAL;
1730:     coarsenBox.min = coarsenBox.max = PETSC_MIN_REAL;

1732:     VecTaggerCreate(comm,&refineTag);
1733:     PetscObjectSetOptionsPrefix((PetscObject)refineTag,"refine_");
1734:     VecTaggerSetType(refineTag,VECTAGGERABSOLUTE);
1735:     VecTaggerAbsoluteSetBox(refineTag,&refineBox);
1736:     VecTaggerSetFromOptions(refineTag);
1737:     VecTaggerSetUp(refineTag);
1738:     PetscObjectViewFromOptions((PetscObject)refineTag,NULL,"-tag_view");

1740:     VecTaggerCreate(comm,&coarsenTag);
1741:     PetscObjectSetOptionsPrefix((PetscObject)coarsenTag,"coarsen_");
1742:     VecTaggerSetType(coarsenTag,VECTAGGERABSOLUTE);
1743:     VecTaggerAbsoluteSetBox(coarsenTag,&coarsenBox);
1744:     VecTaggerSetFromOptions(coarsenTag);
1745:     VecTaggerSetUp(coarsenTag);
1746:     PetscObjectViewFromOptions((PetscObject)coarsenTag,NULL,"-tag_view");
1747:   }

1749:   PetscOptionsBegin(comm,NULL,"Unstructured Finite Volume Physics Options","");
1750:   {
1751:     PetscErrorCode (*physcreate)(Model,Physics,PetscOptionItems*);
1752:     PetscOptionsFList("-physics","Physics module to solve","",PhysicsList,physname,physname,sizeof physname,NULL);
1753:     PetscFunctionListFind(PhysicsList,physname,&physcreate);
1754:     PetscMemzero(phys,sizeof(struct _n_Physics));
1755:     (*physcreate)(mod,phys,PetscOptionsObject);
1756:     /* Count number of fields and dofs */
1757:     for (phys->nfields=0,phys->dof=0; phys->field_desc[phys->nfields].name; phys->nfields++) phys->dof += phys->field_desc[phys->nfields].dof;
1758:     if (phys->dof <= 0) SETERRQ1(comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set dof",physname);
1759:     ModelFunctionalSetFromOptions(mod,PetscOptionsObject);
1760:   }
1761:   PetscOptionsEnd();

1763:   /* Create mesh */
1764:   {
1765:     size_t len,i;
1766:     for (i = 0; i < DIM; i++) { mod->bounds[2*i] = 0.; mod->bounds[2*i+1] = 1.;};
1767:     PetscStrlen(filename,&len);
1768:     dim = DIM;
1769:     if (!len) { /* a null name means just do a hex box */
1770:       PetscInt cells[3] = {1, 1, 1}; /* coarse mesh is one cell; refine from there */
1771:       PetscBool flg1, flg2, skew = PETSC_FALSE;
1772:       PetscInt nret1 = DIM;
1773:       PetscInt nret2 = 2*DIM;
1774:       PetscOptionsBegin(comm,NULL,"Rectangular mesh options","");
1775:       PetscOptionsIntArray("-grid_size","number of cells in each direction","",cells,&nret1,&flg1);
1776:       PetscOptionsRealArray("-grid_bounds","bounds of the mesh in each direction (i.e., x_min,x_max,y_min,y_max","",mod->bounds,&nret2,&flg2);
1777:       PetscOptionsBool("-grid_skew_60","Skew grid for 60 degree shock mesh","",skew,&skew,NULL);
1778:       PetscOptionsEnd();
1779:       if (flg1) {
1780:         dim = nret1;
1781:         if (dim != DIM) SETERRQ1(comm,PETSC_ERR_ARG_SIZ,"Dim wrong size %D in -grid_size",dim);
1782:       }
1783:       DMPlexCreateBoxMesh(comm, dim, simplex, cells, NULL, NULL, mod->bcs, PETSC_TRUE, &dm);
1784:       if (flg2) {
1785:         PetscInt dimEmbed, i;
1786:         PetscInt nCoords;
1787:         PetscScalar *coords;
1788:         Vec coordinates;

1790:         DMGetCoordinatesLocal(dm,&coordinates);
1791:         DMGetCoordinateDim(dm,&dimEmbed);
1792:         VecGetLocalSize(coordinates,&nCoords);
1793:         if (nCoords % dimEmbed) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Coordinate vector the wrong size");
1794:         VecGetArray(coordinates,&coords);
1795:         for (i = 0; i < nCoords; i += dimEmbed) {
1796:           PetscInt j;

1798:           PetscScalar *coord = &coords[i];
1799:           for (j = 0; j < dimEmbed; j++) {
1800:             coord[j] = mod->bounds[2 * j] + coord[j] * (mod->bounds[2 * j + 1] - mod->bounds[2 * j]);
1801:             if (dim==2 && cells[1]==1 && j==0 && skew) {
1802:               if (cells[0]==2 && i==8) {
1803:                 coord[j] = .57735026918963; /* hack to get 60 deg skewed mesh */
1804:               }
1805:               else if (cells[0]==3) {
1806:                 if (i==2 || i==10) coord[j] = mod->bounds[1]/4.;
1807:                 else if (i==4) coord[j] = mod->bounds[1]/2.;
1808:                 else if (i==12) coord[j] = 1.57735026918963*mod->bounds[1]/2.;
1809:               }
1810:             }
1811:           }
1812:         }
1813:         VecRestoreArray(coordinates,&coords);
1814:         DMSetCoordinatesLocal(dm,coordinates);
1815:       }
1816:     } else {
1817:       DMPlexCreateFromFile(comm, filename, PETSC_TRUE, &dm);
1818:     }
1819:   }
1820:   DMViewFromOptions(dm, NULL, "-orig_dm_view");
1821:   DMGetDimension(dm, &dim);

1823:   /* set up BCs, functions, tags */
1824:   DMCreateLabel(dm, "Face Sets");

1826:   mod->errorIndicator = ErrorIndicator_Simple;

1828:   {
1829:     DM dmDist;

1831:     DMSetBasicAdjacency(dm, PETSC_TRUE, PETSC_FALSE);
1832:     DMPlexDistribute(dm, overlap, NULL, &dmDist);
1833:     if (dmDist) {
1834:       DMDestroy(&dm);
1835:       dm   = dmDist;
1836:     }
1837:   }

1839:   DMSetFromOptions(dm);

1841:   {
1842:     DM gdm;

1844:     DMPlexConstructGhostCells(dm, NULL, NULL, &gdm);
1845:     DMDestroy(&dm);
1846:     dm   = gdm;
1847:     DMViewFromOptions(dm, NULL, "-dm_view");
1848:   }
1849:   if (splitFaces) {ConstructCellBoundary(dm, user);}
1850:   SplitFaces(&dm, "split faces", user);

1852:   PetscFVCreate(comm, &fvm);
1853:   PetscFVSetFromOptions(fvm);
1854:   PetscFVSetNumComponents(fvm, phys->dof);
1855:   PetscFVSetSpatialDimension(fvm, dim);
1856:   PetscObjectSetName((PetscObject) fvm,"");
1857:   {
1858:     PetscInt f, dof;
1859:     for (f=0,dof=0; f < phys->nfields; f++) {
1860:       PetscInt newDof = phys->field_desc[f].dof;

1862:       if (newDof == 1) {
1863:         PetscFVSetComponentName(fvm,dof,phys->field_desc[f].name);
1864:       }
1865:       else {
1866:         PetscInt j;

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

1871:           PetscSNPrintf(compName,sizeof(compName),"%s_%d",phys->field_desc[f].name,j);
1872:           PetscFVSetComponentName(fvm,dof+j,compName);
1873:         }
1874:       }
1875:       dof += newDof;
1876:     }
1877:   }
1878:   /* FV is now structured with one field having all physics as components */
1879:   DMAddField(dm, NULL, (PetscObject) fvm);
1880:   DMCreateDS(dm);
1881:   DMGetDS(dm, &prob);
1882:   PetscDSSetRiemannSolver(prob, 0, user->model->physics->riemann);
1883:   PetscDSSetContext(prob, 0, user->model->physics);
1884:   (*mod->setupbc)(dm, prob,phys);
1885:   PetscDSSetFromOptions(prob);
1886:   {
1887:     char      convType[256];
1888:     PetscBool flg;

1890:     PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");
1891:     PetscOptionsFList("-dm_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg);
1892:     PetscOptionsEnd();
1893:     if (flg) {
1894:       DM dmConv;

1896:       DMConvert(dm,convType,&dmConv);
1897:       if (dmConv) {
1898:         DMViewFromOptions(dmConv, NULL, "-dm_conv_view");
1899:         DMDestroy(&dm);
1900:         dm   = dmConv;
1901:         DMSetFromOptions(dm);
1902:       }
1903:     }
1904:   }

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

1908:   DMCreateGlobalVector(dm, &X);
1909:   PetscObjectSetName((PetscObject) X, "solution");
1910:   SetInitialCondition(dm, X, user);
1911:   if (useAMR) {
1912:     PetscInt adaptIter;

1914:     /* use no limiting when reconstructing gradients for adaptivity */
1915:     PetscFVGetLimiter(fvm, &limiter);
1916:     PetscObjectReference((PetscObject) limiter);
1917:     PetscLimiterCreate(PetscObjectComm((PetscObject) fvm), &noneLimiter);
1918:     PetscLimiterSetType(noneLimiter, PETSCLIMITERNONE);

1920:     PetscFVSetLimiter(fvm, noneLimiter);
1921:     for (adaptIter = 0; ; ++adaptIter) {
1922:       PetscLogDouble bytes;
1923:       TS             tsNew = NULL;

1925:       PetscMemoryGetCurrentUsage(&bytes);
1926:       PetscInfo2(ts, "refinement loop %D: memory used %g\n", adaptIter, bytes);
1927:       DMViewFromOptions(dm, NULL, "-initial_dm_view");
1928:       VecViewFromOptions(X, NULL, "-initial_vec_view");
1929: #if 0
1930:       if (viewInitial) {
1931:         PetscViewer viewer;
1932:         char        buf[256];
1933:         PetscBool   isHDF5, isVTK;

1935:         PetscViewerCreate(comm,&viewer);
1936:         PetscViewerSetType(viewer,PETSCVIEWERVTK);
1937:         PetscViewerSetOptionsPrefix(viewer,"initial_");
1938:         PetscViewerSetFromOptions(viewer);
1939:         PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&isHDF5);
1940:         PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isVTK);
1941:         if (isHDF5) {
1942:           PetscSNPrintf(buf, 256, "ex11-initial-%d.h5", adaptIter);
1943:         } else if (isVTK) {
1944:           PetscSNPrintf(buf, 256, "ex11-initial-%d.vtu", adaptIter);
1945:           PetscViewerPushFormat(viewer,PETSC_VIEWER_VTK_VTU);
1946:         }
1947:         PetscViewerFileSetMode(viewer,FILE_MODE_WRITE);
1948:         PetscViewerFileSetName(viewer,buf);
1949:         if (isHDF5) {
1950:           DMView(dm,viewer);
1951:           PetscViewerFileSetMode(viewer,FILE_MODE_UPDATE);
1952:         }
1953:         VecView(X,viewer);
1954:         PetscViewerDestroy(&viewer);
1955:       }
1956: #endif

1958:       adaptToleranceFVM(fvm, ts, X, refineTag, coarsenTag, user, &tsNew, NULL);
1959:       if (!tsNew) {
1960:         break;
1961:       } else {
1962:         DMDestroy(&dm);
1963:         VecDestroy(&X);
1964:         TSDestroy(&ts);
1965:         ts   = tsNew;
1966:         TSGetDM(ts,&dm);
1967:         PetscObjectReference((PetscObject)dm);
1968:         DMCreateGlobalVector(dm,&X);
1969:         PetscObjectSetName((PetscObject) X, "solution");
1970:         SetInitialCondition(dm, X, user);
1971:       }
1972:     }
1973:     /* restore original limiter */
1974:     PetscFVSetLimiter(fvm, limiter);
1975:   }

1977:   DMConvert(dm, DMPLEX, &plex);
1978:   if (vtkCellGeom) {
1979:     DM  dmCell;
1980:     Vec cellgeom, partition;

1982:     DMPlexGetGeometryFVM(plex, NULL, &cellgeom, NULL);
1983:     OutputVTK(dm, "ex11-cellgeom.vtk", &viewer);
1984:     VecView(cellgeom, viewer);
1985:     PetscViewerDestroy(&viewer);
1986:     CreatePartitionVec(dm, &dmCell, &partition);
1987:     OutputVTK(dmCell, "ex11-partition.vtk", &viewer);
1988:     VecView(partition, viewer);
1989:     PetscViewerDestroy(&viewer);
1990:     VecDestroy(&partition);
1991:     DMDestroy(&dmCell);
1992:   }
1993:   /* collect max maxspeed from all processes -- todo */
1994:   DMPlexGetGeometryFVM(plex, NULL, NULL, &minRadius);
1995:   DMDestroy(&plex);
1996:   MPI_Allreduce(&phys->maxspeed,&mod->maxspeed,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));
1997:   if (mod->maxspeed <= 0) SETERRQ1(comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set maxspeed",physname);
1998:   dt   = cfl * minRadius / mod->maxspeed;
1999:   TSSetTimeStep(ts,dt);
2000:   TSSetFromOptions(ts);
2001:   if (!useAMR) {
2002:     TSSolve(ts,X);
2003:     TSGetSolveTime(ts,&ftime);
2004:     TSGetStepNumber(ts,&nsteps);
2005:   } else {
2006:     PetscReal finalTime;
2007:     PetscInt  adaptIter;
2008:     TS        tsNew = NULL;
2009:     Vec       solNew = NULL;

2011:     TSGetMaxTime(ts,&finalTime);
2012:     TSSetMaxSteps(ts,adaptInterval);
2013:     TSSolve(ts,X);
2014:     TSGetSolveTime(ts,&ftime);
2015:     TSGetStepNumber(ts,&nsteps);
2016:     for (adaptIter = 0;ftime < finalTime;adaptIter++) {
2017:       PetscLogDouble bytes;

2019:       PetscMemoryGetCurrentUsage(&bytes);
2020:       PetscInfo2(ts, "AMR time step loop %D: memory used %g\n", adaptIter, bytes);
2021:       PetscFVSetLimiter(fvm,noneLimiter);
2022:       adaptToleranceFVM(fvm,ts,X,refineTag,coarsenTag,user,&tsNew,&solNew);
2023:       PetscFVSetLimiter(fvm,limiter);
2024:       if (tsNew) {
2025:         PetscInfo(ts, "AMR used\n");
2026:         DMDestroy(&dm);
2027:         VecDestroy(&X);
2028:         TSDestroy(&ts);
2029:         ts   = tsNew;
2030:         X    = solNew;
2031:         TSSetFromOptions(ts);
2032:         VecGetDM(X,&dm);
2033:         PetscObjectReference((PetscObject)dm);
2034:         DMConvert(dm, DMPLEX, &plex);
2035:         DMPlexGetGeometryFVM(dm, NULL, NULL, &minRadius);
2036:         DMDestroy(&plex);
2037:         MPI_Allreduce(&phys->maxspeed,&mod->maxspeed,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));
2038:         if (mod->maxspeed <= 0) SETERRQ1(comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set maxspeed",physname);
2039:         dt   = cfl * minRadius / mod->maxspeed;
2040:         TSSetStepNumber(ts,nsteps);
2041:         TSSetTime(ts,ftime);
2042:         TSSetTimeStep(ts,dt);
2043:       } else {
2044:         PetscInfo(ts, "AMR not used\n");
2045:       }
2046:       user->monitorStepOffset = nsteps;
2047:       TSSetMaxSteps(ts,nsteps+adaptInterval);
2048:       TSSolve(ts,X);
2049:       TSGetSolveTime(ts,&ftime);
2050:       TSGetStepNumber(ts,&nsteps);
2051:     }
2052:   }
2053:   TSGetConvergedReason(ts,&reason);
2054:   PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,nsteps);
2055:   TSDestroy(&ts);

2057:   VecTaggerDestroy(&refineTag);
2058:   VecTaggerDestroy(&coarsenTag);
2059:   PetscFunctionListDestroy(&PhysicsList);
2060:   PetscFunctionListDestroy(&PhysicsRiemannList_SW);
2061:   FunctionalLinkDestroy(&user->model->functionalRegistry);
2062:   PetscFree(user->model->functionalMonitored);
2063:   PetscFree(user->model->functionalCall);
2064:   PetscFree(user->model->physics->data);
2065:   PetscFree(user->model->physics);
2066:   PetscFree(user->model);
2067:   PetscFree(user);
2068:   VecDestroy(&X);
2069:   PetscLimiterDestroy(&limiter);
2070:   PetscLimiterDestroy(&noneLimiter);
2071:   PetscFVDestroy(&fvm);
2072:   DMDestroy(&dm);
2073:   PetscFinalize();
2074:   return ierr;
2075: }

2077: /* Godunov fluxs */
2078: PetscScalar cvmgp_(PetscScalar *a, PetscScalar *b, PetscScalar *test)
2079: {
2080:     /* System generated locals */
2081:     PetscScalar ret_val;

2083:     if (PetscRealPart(*test) > 0.) {
2084:         goto L10;
2085:     }
2086:     ret_val = *b;
2087:     return ret_val;
2088: L10:
2089:     ret_val = *a;
2090:     return ret_val;
2091: } /* cvmgp_ */

2093: PetscScalar cvmgm_(PetscScalar *a, PetscScalar *b, PetscScalar *test)
2094: {
2095:     /* System generated locals */
2096:     PetscScalar ret_val;

2098:     if (PetscRealPart(*test) < 0.) {
2099:         goto L10;
2100:     }
2101:     ret_val = *b;
2102:     return ret_val;
2103: L10:
2104:     ret_val = *a;
2105:     return ret_val;
2106: } /* cvmgm_ */

2108: int riem1mdt( PetscScalar *gaml, PetscScalar *gamr, PetscScalar *rl, PetscScalar *pl,
2109:               PetscScalar *uxl, PetscScalar *rr, PetscScalar *pr,
2110:               PetscScalar *uxr, PetscScalar *rstarl, PetscScalar *rstarr, PetscScalar *
2111:               pstar, PetscScalar *ustar)
2112: {
2113:     /* Initialized data */

2115:     static PetscScalar smallp = 1e-8;

2117:     /* System generated locals */
2118:     int i__1;
2119:     PetscScalar d__1, d__2;

2121:     /* Local variables */
2122:     static int i0;
2123:     static PetscScalar cl, cr, wl, zl, wr, zr, pst, durl, skpr1, skpr2;
2124:     static int iwave;
2125:     static PetscScalar gascl4, gascr4, cstarl, dpstar, cstarr;
2126:     /* static PetscScalar csqrl, csqrr, gascl1, gascl2, gascl3, gascr1, gascr2, gascr3; */
2127:     static int iterno;
2128:     static PetscScalar ustarl, ustarr, rarepr1, rarepr2;



2132:     /* gascl1 = *gaml - 1.; */
2133:     /* gascl2 = (*gaml + 1.) * .5; */
2134:     /* gascl3 = gascl2 / *gaml; */
2135:     gascl4 = 1. / (*gaml - 1.);

2137:     /* gascr1 = *gamr - 1.; */
2138:     /* gascr2 = (*gamr + 1.) * .5; */
2139:     /* gascr3 = gascr2 / *gamr; */
2140:     gascr4 = 1. / (*gamr - 1.);
2141:     iterno = 10;
2142: /*        find pstar: */
2143:     cl = PetscSqrtScalar(*gaml * *pl / *rl);
2144:     cr = PetscSqrtScalar(*gamr * *pr / *rr);
2145:     wl = *rl * cl;
2146:     wr = *rr * cr;
2147:     /* csqrl = wl * wl; */
2148:     /* csqrr = wr * wr; */
2149:     *pstar = (wl * *pr + wr * *pl) / (wl + wr);
2150:     *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp));
2151:     pst = *pl / *pr;
2152:     skpr1 = cr * (pst - 1.) * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
2153:     d__1 = (*gamr - 1.) / (*gamr * 2.);
2154:     rarepr2 = gascr4 * 2. * cr * (1. - PetscPowScalar(pst, d__1));
2155:     pst = *pr / *pl;
2156:     skpr2 = cl * (pst - 1.) * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
2157:     d__1 = (*gaml - 1.) / (*gaml * 2.);
2158:     rarepr1 = gascl4 * 2. * cl * (1. - PetscPowScalar(pst, d__1));
2159:     durl = *uxr - *uxl;
2160:     if (PetscRealPart(*pr) < PetscRealPart(*pl)) {
2161:         if (PetscRealPart(durl) >= PetscRealPart(rarepr1)) {
2162:             iwave = 100;
2163:         } else if (PetscRealPart(durl) <= PetscRealPart(-skpr1)) {
2164:             iwave = 300;
2165:         } else {
2166:             iwave = 400;
2167:         }
2168:     } else {
2169:         if (PetscRealPart(durl) >= PetscRealPart(rarepr2)) {
2170:             iwave = 100;
2171:         } else if (PetscRealPart(durl) <= PetscRealPart(-skpr2)) {
2172:             iwave = 300;
2173:         } else {
2174:             iwave = 200;
2175:         }
2176:     }
2177:     if (iwave == 100) {
2178: /*     1-wave: rarefaction wave, 3-wave: rarefaction wave */
2179: /*     case (100) */
2180:         i__1 = iterno;
2181:         for (i0 = 1; i0 <= i__1; ++i0) {
2182:             d__1 = *pstar / *pl;
2183:             d__2 = 1. / *gaml;
2184:             *rstarl = *rl * PetscPowScalar(d__1, d__2);
2185:             cstarl = PetscSqrtScalar(*gaml * *pstar / *rstarl);
2186:             ustarl = *uxl - gascl4 * 2. * (cstarl - cl);
2187:             zl = *rstarl * cstarl;
2188:             d__1 = *pstar / *pr;
2189:             d__2 = 1. / *gamr;
2190:             *rstarr = *rr * PetscPowScalar(d__1, d__2);
2191:             cstarr = PetscSqrtScalar(*gamr * *pstar / *rstarr);
2192:             ustarr = *uxr + gascr4 * 2. * (cstarr - cr);
2193:             zr = *rstarr * cstarr;
2194:             dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
2195:             *pstar -= dpstar;
2196:             *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp));
2197:             if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
2198: #if 0
2199:         break;
2200: #endif
2201:             }
2202:         }
2203: /*     1-wave: shock wave, 3-wave: rarefaction wave */
2204:     } else if (iwave == 200) {
2205: /*     case (200) */
2206:         i__1 = iterno;
2207:         for (i0 = 1; i0 <= i__1; ++i0) {
2208:             pst = *pstar / *pl;
2209:             ustarl = *uxl - (pst - 1.) * cl * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
2210:             zl = *pl / cl * PetscSqrtScalar(*gaml * 2. * (*gaml - 1. + (*gaml + 1.) * pst)) * (*gaml - 1. + (*gaml + 1.) * pst) / (*gaml * 3. - 1. + (*gaml + 1.) * pst);
2211:             d__1 = *pstar / *pr;
2212:             d__2 = 1. / *gamr;
2213:             *rstarr = *rr * PetscPowScalar(d__1, d__2);
2214:             cstarr = PetscSqrtScalar(*gamr * *pstar / *rstarr);
2215:             zr = *rstarr * cstarr;
2216:             ustarr = *uxr + gascr4 * 2. * (cstarr - cr);
2217:             dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
2218:             *pstar -= dpstar;
2219:             *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp));
2220:             if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
2221: #if 0
2222:         break;
2223: #endif
2224:             }
2225:         }
2226: /*     1-wave: shock wave, 3-wave: shock */
2227:     } else if (iwave == 300) {
2228: /*     case (300) */
2229:         i__1 = iterno;
2230:         for (i0 = 1; i0 <= i__1; ++i0) {
2231:             pst = *pstar / *pl;
2232:             ustarl = *uxl - (pst - 1.) * cl * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
2233:             zl = *pl / cl * PetscSqrtScalar(*gaml * 2. * (*gaml - 1. + (*gaml + 1.) * pst)) * (*gaml - 1. + (*gaml + 1.) * pst) / (*gaml * 3. - 1. + (*gaml + 1.) * pst);
2234:             pst = *pstar / *pr;
2235:             ustarr = *uxr + (pst - 1.) * cr * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
2236:             zr = *pr / cr * PetscSqrtScalar(*gamr * 2. * (*gamr - 1. + (*gamr + 1.) * pst)) * (*gamr - 1. + (*gamr + 1.) * pst) / (*gamr * 3. - 1. + (*gamr + 1.) * pst);
2237:             dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
2238:             *pstar -= dpstar;
2239:             *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp));
2240:             if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
2241: #if 0
2242:         break;
2243: #endif
2244:             }
2245:         }
2246: /*     1-wave: rarefaction wave, 3-wave: shock */
2247:     } else if (iwave == 400) {
2248: /*     case (400) */
2249:         i__1 = iterno;
2250:         for (i0 = 1; i0 <= i__1; ++i0) {
2251:             d__1 = *pstar / *pl;
2252:             d__2 = 1. / *gaml;
2253:             *rstarl = *rl * PetscPowScalar(d__1, d__2);
2254:             cstarl = PetscSqrtScalar(*gaml * *pstar / *rstarl);
2255:             ustarl = *uxl - gascl4 * 2. * (cstarl - cl);
2256:             zl = *rstarl * cstarl;
2257:             pst = *pstar / *pr;
2258:             ustarr = *uxr + (pst - 1.) * cr * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
2259:             zr = *pr / cr * PetscSqrtScalar(*gamr * 2. * (*gamr - 1. + (*gamr + 1.) * pst)) * (*gamr - 1. + (*gamr + 1.) * pst) / (*gamr * 3. - 1. + (*gamr + 1.) * pst);
2260:             dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
2261:             *pstar -= dpstar;
2262:             *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp));
2263:             if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
2264: #if 0
2265:               break;
2266: #endif
2267:             }
2268:         }
2269:     }

2271:     *ustar = (zl * ustarr + zr * ustarl) / (zl + zr);
2272:     if (PetscRealPart(*pstar) > PetscRealPart(*pl)) {
2273:         pst = *pstar / *pl;
2274:         *rstarl = ((*gaml + 1.) * pst + *gaml - 1.) / ((*gaml - 1.) * pst + *
2275:                 gaml + 1.) * *rl;
2276:     }
2277:     if (PetscRealPart(*pstar) > PetscRealPart(*pr)) {
2278:         pst = *pstar / *pr;
2279:         *rstarr = ((*gamr + 1.) * pst + *gamr - 1.) / ((*gamr - 1.) * pst + *
2280:                 gamr + 1.) * *rr;
2281:     }
2282:     return iwave;
2283: }

2285: PetscScalar sign(PetscScalar x)
2286: {
2287:     if (PetscRealPart(x) > 0) return 1.0;
2288:     if (PetscRealPart(x) < 0) return -1.0;
2289:     return 0.0;
2290: }
2291: /*        Riemann Solver */
2292: /* -------------------------------------------------------------------- */
2293: int riemannsolver(PetscScalar *xcen, PetscScalar *xp,
2294:                    PetscScalar *dtt, PetscScalar *rl, PetscScalar *uxl, PetscScalar *pl,
2295:                    PetscScalar *utl, PetscScalar *ubl, PetscScalar *gaml, PetscScalar *rho1l,
2296:                    PetscScalar *rr, PetscScalar *uxr, PetscScalar *pr, PetscScalar *utr,
2297:                    PetscScalar *ubr, PetscScalar *gamr, PetscScalar *rho1r, PetscScalar *rx,
2298:                    PetscScalar *uxm, PetscScalar *px, PetscScalar *utx, PetscScalar *ubx,
2299:                    PetscScalar *gam, PetscScalar *rho1)
2300: {
2301:     /* System generated locals */
2302:     PetscScalar d__1, d__2;

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

2309:     if (*rl == *rr && *pr == *pl && *uxl == *uxr && *gaml == *gamr) {
2310:         *rx = *rl;
2311:         *px = *pl;
2312:         *uxm = *uxl;
2313:         *gam = *gaml;
2314:         x2 = *xcen + *uxm * *dtt;

2316:         if (PetscRealPart(*xp) >= PetscRealPart(x2)) {
2317:             *utx = *utr;
2318:             *ubx = *ubr;
2319:             *rho1 = *rho1r;
2320:         } else {
2321:             *utx = *utl;
2322:             *ubx = *ubl;
2323:             *rho1 = *rho1l;
2324:         }
2325:         return 0;
2326:     }
2327:     iwave = riem1mdt(gaml, gamr, rl, pl, uxl, rr, pr, uxr, &rstarl, &rstarr, &pstar, &ustar);

2329:     x2 = *xcen + ustar * *dtt;
2330:     d__1 = *xp - x2;
2331:     sgn0 = sign(d__1);
2332: /*            x is in 3-wave if sgn0 = 1 */
2333: /*            x is in 1-wave if sgn0 = -1 */
2334:     r0 = cvmgm_(rl, rr, &sgn0);
2335:     p0 = cvmgm_(pl, pr, &sgn0);
2336:     u0 = cvmgm_(uxl, uxr, &sgn0);
2337:     *gam = cvmgm_(gaml, gamr, &sgn0);
2338:     gasc1 = *gam - 1.;
2339:     gasc2 = (*gam + 1.) * .5;
2340:     gasc3 = gasc2 / *gam;
2341:     gasc4 = 1. / (*gam - 1.);
2342:     c0 = PetscSqrtScalar(*gam * p0 / r0);
2343:     streng = pstar - p0;
2344:     w0 = *gam * r0 * p0 * (gasc3 * streng / p0 + 1.);
2345:     rstars = r0 / (1. - r0 * streng / w0);
2346:     d__1 = p0 / pstar;
2347:     d__2 = -1. / *gam;
2348:     rstarr = r0 * PetscPowScalar(d__1, d__2);
2349:     rstar = cvmgm_(&rstarr, &rstars, &streng);
2350:     w0 = PetscSqrtScalar(w0);
2351:     cstar = PetscSqrtScalar(*gam * pstar / rstar);
2352:     wsp0 = u0 + sgn0 * c0;
2353:     wspst = ustar + sgn0 * cstar;
2354:     ushock = ustar + sgn0 * w0 / rstar;
2355:     wspst = cvmgp_(&ushock, &wspst, &streng);
2356:     wsp0 = cvmgp_(&ushock, &wsp0, &streng);
2357:     x0 = *xcen + wsp0 * *dtt;
2358:     xstar = *xcen + wspst * *dtt;
2359: /*           using gas formula to evaluate rarefaction wave */
2360: /*            ri : reiman invariant */
2361:     ri = u0 - sgn0 * 2. * gasc4 * c0;
2362:     cx = sgn0 * .5 * gasc1 / gasc2 * ((*xp - *xcen) / *dtt - ri);
2363:     *uxm = ri + sgn0 * 2. * gasc4 * cx;
2364:     s = p0 / PetscPowScalar(r0, *gam);
2365:     d__1 = cx * cx / (*gam * s);
2366:     *rx = PetscPowScalar(d__1, gasc4);
2367:     *px = cx * cx * *rx / *gam;
2368:     d__1 = sgn0 * (x0 - *xp);
2369:     *rx = cvmgp_(rx, &r0, &d__1);
2370:     d__1 = sgn0 * (x0 - *xp);
2371:     *px = cvmgp_(px, &p0, &d__1);
2372:     d__1 = sgn0 * (x0 - *xp);
2373:     *uxm = cvmgp_(uxm, &u0, &d__1);
2374:     d__1 = sgn0 * (xstar - *xp);
2375:     *rx = cvmgm_(rx, &rstar, &d__1);
2376:     d__1 = sgn0 * (xstar - *xp);
2377:     *px = cvmgm_(px, &pstar, &d__1);
2378:     d__1 = sgn0 * (xstar - *xp);
2379:     *uxm = cvmgm_(uxm, &ustar, &d__1);
2380:     if (PetscRealPart(*xp) >= PetscRealPart(x2)) {
2381:         *utx = *utr;
2382:         *ubx = *ubr;
2383:         *rho1 = *rho1r;
2384:     } else {
2385:         *utx = *utl;
2386:         *ubx = *ubl;
2387:         *rho1 = *rho1l;
2388:     }
2389:     return iwave;
2390: }
2391: int godunovflux( const PetscScalar *ul, const PetscScalar *ur,
2392:                  PetscScalar *flux, const PetscReal *nn, const int *ndim,
2393:                  const PetscReal *gamma)
2394: {
2395:     /* System generated locals */
2396:   int i__1,iwave;
2397:     PetscScalar d__1, d__2, d__3;

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

2405:     /* Function Body */
2406:     xcen = 0.;
2407:     xp = 0.;
2408:     i__1 = *ndim;
2409:     for (k = 1; k <= i__1; ++k) {
2410:         tg[k - 1] = 0.;
2411:         bn[k - 1] = 0.;
2412:     }
2413:     dtt = 1.;
2414:     if (*ndim == 3) {
2415:         if (nn[0] == 0. && nn[1] == 0.) {
2416:             tg[0] = 1.;
2417:         } else {
2418:             tg[0] = -nn[1];
2419:             tg[1] = nn[0];
2420:         }
2421: /*           tmp=dsqrt(tg(1)**2+tg(2)**2) */
2422: /*           tg=tg/tmp */
2423:         bn[0] = -nn[2] * tg[1];
2424:         bn[1] = nn[2] * tg[0];
2425:         bn[2] = nn[0] * tg[1] - nn[1] * tg[0];
2426: /* Computing 2nd power */
2427:         d__1 = bn[0];
2428: /* Computing 2nd power */
2429:         d__2 = bn[1];
2430: /* Computing 2nd power */
2431:         d__3 = bn[2];
2432:         tmp = PetscSqrtScalar(d__1 * d__1 + d__2 * d__2 + d__3 * d__3);
2433:         i__1 = *ndim;
2434:         for (k = 1; k <= i__1; ++k) {
2435:             bn[k - 1] /= tmp;
2436:         }
2437:     } else if (*ndim == 2) {
2438:         tg[0] = -nn[1];
2439:         tg[1] = nn[0];
2440: /*           tmp=dsqrt(tg(1)**2+tg(2)**2) */
2441: /*           tg=tg/tmp */
2442:         bn[0] = 0.;
2443:         bn[1] = 0.;
2444:         bn[2] = 1.;
2445:     }
2446:     rl = ul[0];
2447:     rr = ur[0];
2448:     uxl = 0.;
2449:     uxr = 0.;
2450:     utl = 0.;
2451:     utr = 0.;
2452:     ubl = 0.;
2453:     ubr = 0.;
2454:     i__1 = *ndim;
2455:     for (k = 1; k <= i__1; ++k) {
2456:         uxl += ul[k] * nn[k-1];
2457:         uxr += ur[k] * nn[k-1];
2458:         utl += ul[k] * tg[k - 1];
2459:         utr += ur[k] * tg[k - 1];
2460:         ubl += ul[k] * bn[k - 1];
2461:         ubr += ur[k] * bn[k - 1];
2462:     }
2463:     uxl /= rl;
2464:     uxr /= rr;
2465:     utl /= rl;
2466:     utr /= rr;
2467:     ubl /= rl;
2468:     ubr /= rr;

2470:     gaml = *gamma;
2471:     gamr = *gamma;
2472: /* Computing 2nd power */
2473:     d__1 = uxl;
2474: /* Computing 2nd power */
2475:     d__2 = utl;
2476: /* Computing 2nd power */
2477:     d__3 = ubl;
2478:     pl = (*gamma - 1.) * (ul[*ndim + 1] - rl * .5 * (d__1 * d__1 + d__2 * d__2 + d__3 * d__3));
2479: /* Computing 2nd power */
2480:     d__1 = uxr;
2481: /* Computing 2nd power */
2482:     d__2 = utr;
2483: /* Computing 2nd power */
2484:     d__3 = ubr;
2485:     pr = (*gamma - 1.) * (ur[*ndim + 1] - rr * .5 * (d__1 * d__1 + d__2 * d__2 + d__3 * d__3));
2486:     rho1l = rl;
2487:     rho1r = rr;

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

2493:     flux[0] = rhom * unm;
2494:     fn = rhom * unm * unm + pm;
2495:     ft = rhom * unm * utm;
2496: /*           flux(2)=fn*nn(1)+ft*nn(2) */
2497: /*           flux(3)=fn*tg(1)+ft*tg(2) */
2498:     flux[1] = fn * nn[0] + ft * tg[0];
2499:     flux[2] = fn * nn[1] + ft * tg[1];
2500: /*           flux(2)=rhom*unm*(unm)+pm */
2501: /*           flux(3)=rhom*(unm)*utm */
2502:     if (*ndim == 3) {
2503:         flux[3] = rhom * unm * ubm;
2504:     }
2505:     flux[*ndim + 1] = (rhom * .5 * (unm * unm + utm * utm + ubm * ubm) + gamm / (gamm - 1.) * pm) * unm;
2506:     return iwave;
2507: } /* godunovflux_ */

2509: /* Subroutine to set up the initial conditions for the */
2510: /* Shock Interface interaction or linear wave (Ravi Samtaney,Mark Adams). */
2511: /* ----------------------------------------------------------------------- */
2512: int projecteqstate(PetscReal wc[], const PetscReal ueq[], PetscReal lv[][3])
2513: {
2514:   int j,k;
2515: /*      Wc=matmul(lv,Ueq) 3 vars */
2516:   for (k = 0; k < 3; ++k) {
2517:     wc[k] = 0.;
2518:     for (j = 0; j < 3; ++j) {
2519:       wc[k] += lv[k][j]*ueq[j];
2520:     }
2521:   }
2522:   return 0;
2523: }
2524: /* ----------------------------------------------------------------------- */
2525: int projecttoprim(PetscReal v[], const PetscReal wc[], PetscReal rv[][3])
2526: {
2527:   int k,j;
2528:   /*      V=matmul(rv,WC) 3 vars */
2529:   for (k = 0; k < 3; ++k) {
2530:     v[k] = 0.;
2531:     for (j = 0; j < 3; ++j) {
2532:       v[k] += rv[k][j]*wc[j];
2533:     }
2534:   }
2535:   return 0;
2536: }
2537: /* ---------------------------------------------------------------------- */
2538: int eigenvectors(PetscReal rv[][3], PetscReal lv[][3], const PetscReal ueq[], PetscReal gamma)
2539: {
2540:   int j,k;
2541:   PetscReal rho,csnd,p0;
2542:   /* PetscScalar u; */

2544:   for (k = 0; k < 3; ++k) for (j = 0; j < 3; ++j) { lv[k][j] = 0.; rv[k][j] = 0.; }
2545:   rho = ueq[0];
2546:   /* u = ueq[1]; */
2547:   p0 = ueq[2];
2548:   csnd = PetscSqrtReal(gamma * p0 / rho);
2549:   lv[0][1] = rho * .5;
2550:   lv[0][2] = -.5 / csnd;
2551:   lv[1][0] = csnd;
2552:   lv[1][2] = -1. / csnd;
2553:   lv[2][1] = rho * .5;
2554:   lv[2][2] = .5 / csnd;
2555:   rv[0][0] = -1. / csnd;
2556:   rv[1][0] = 1. / rho;
2557:   rv[2][0] = -csnd;
2558:   rv[0][1] = 1. / csnd;
2559:   rv[0][2] = 1. / csnd;
2560:   rv[1][2] = 1. / rho;
2561:   rv[2][2] = csnd;
2562:   return 0;
2563: }

2565: int initLinearWave(EulerNode *ux, const PetscReal gamma, const PetscReal coord[], const PetscReal Lx)
2566: {
2567:   PetscReal p0,u0,wcp[3],wc[3];
2568:   PetscReal lv[3][3];
2569:   PetscReal vp[3];
2570:   PetscReal rv[3][3];
2571:   PetscReal eps, ueq[3], rho0, twopi;

2573:   /* Function Body */
2574:   twopi = 2.*PETSC_PI;
2575:   eps = 1e-4; /* perturbation */
2576:   rho0 = 1e3;   /* density of water */
2577:   p0 = 101325.; /* init pressure of 1 atm (?) */
2578:   u0 = 0.;
2579:   ueq[0] = rho0;
2580:   ueq[1] = u0;
2581:   ueq[2] = p0;
2582:   /* Project initial state to characteristic variables */
2583:   eigenvectors(rv, lv, ueq, gamma);
2584:   projecteqstate(wc, ueq, lv);
2585:   wcp[0] = wc[0];
2586:   wcp[1] = wc[1];
2587:   wcp[2] = wc[2] + eps * PetscCosReal(coord[0] * 2. * twopi / Lx);
2588:   projecttoprim(vp, wcp, rv);
2589:   ux->r = vp[0]; /* density */
2590:   ux->ru[0] = vp[0] * vp[1]; /* x momentum */
2591:   ux->ru[1] = 0.;
2592: #if defined DIM > 2
2593:   if (dim>2) ux->ru[2] = 0.;
2594: #endif
2595:   /* E = rho * e + rho * v^2/2 = p/(gam-1) + rho*v^2/2 */
2596:   ux->E = vp[2]/(gamma - 1.) + 0.5*vp[0]*vp[1]*vp[1];
2597:   return 0;
2598: }

2600: /*TEST

2602:   # 2D Advection 0-10
2603:   test:
2604:     suffix: 0
2605:     requires: exodusii
2606:     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo

2608:   test:
2609:     suffix: 1
2610:     requires: exodusii
2611:     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo

2613:   test:
2614:     suffix: 2
2615:     requires: exodusii
2616:     nsize: 2
2617:     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo

2619:   test:
2620:     suffix: 3
2621:     requires: exodusii
2622:     nsize: 2
2623:     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo

2625:   test:
2626:     suffix: 4
2627:     requires: exodusii
2628:     nsize: 8
2629:     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo

2631:   test:
2632:     suffix: 5
2633:     requires: exodusii
2634:     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo -ts_type rosw -ts_adapt_reject_safety 1

2636:   test:
2637:     suffix: 6
2638:     requires: exodusii
2639:     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/squaremotor-30.exo -ufv_split_faces

2641:   test:
2642:     suffix: 7
2643:     requires: exodusii
2644:     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1

2646:   test:
2647:     suffix: 8
2648:     requires: exodusii
2649:     nsize: 2
2650:     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1

2652:   test:
2653:     suffix: 9
2654:     requires: exodusii
2655:     nsize: 8
2656:     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1

2658:   test:
2659:     suffix: 10
2660:     requires: exodusii
2661:     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo

2663:   # 2D Shallow water
2664:   test:
2665:     suffix: sw_0
2666:     requires: exodusii
2667:     args: -ufv_vtk_interval 0 -f ${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

2669:   test:
2670:     suffix: sw_hll
2671:     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 -grid_size 25,25 -sw_riemann hll

2673:   # 2D Advection: p4est
2674:   test:
2675:     suffix: p4est_advec_2d
2676:     requires: p4est
2677:     args: -ufv_vtk_interval 0 -f -dm_type p4est -dm_forest_minimum_refinement 1 -dm_forest_initial_refinement 2 -dm_p4est_refine_pattern hash -dm_forest_maximum_refinement 5

2679:   # Advection in a box
2680:   test:
2681:     suffix: adv_2d_quad_0
2682:     args: -ufv_vtk_interval 0 -dm_refine 3 -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3

2684:   test:
2685:     suffix: adv_2d_quad_1
2686:     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
2687:     timeoutfactor: 3

2689:   test:
2690:     suffix: adv_2d_quad_p4est_0
2691:     requires: p4est
2692:     args: -ufv_vtk_interval 0 -dm_refine 5 -dm_type p4est -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3

2694:   test:
2695:     suffix: adv_2d_quad_p4est_1
2696:     requires: p4est
2697:     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
2698:     timeoutfactor: 3

2700:   test:
2701:     suffix: adv_2d_quad_p4est_adapt_0
2702:     requires: p4est !__float128 #broken for quad precision
2703:     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
2704:     timeoutfactor: 3

2706:   test:
2707:     suffix: adv_2d_tri_0
2708:     requires: triangle
2709:     TODO: how did this ever get in main when there is no support for this
2710:     args: -ufv_vtk_interval 0 -simplex -dm_refine 3 -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3

2712:   test:
2713:     suffix: adv_2d_tri_1
2714:     requires: triangle
2715:     TODO: how did this ever get in main when there is no support for this
2716:     args: -ufv_vtk_interval 0 -simplex -dm_refine 5 -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

2718:   test:
2719:     suffix: adv_0
2720:     requires: exodusii
2721:     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -bc_inflow 100,101,200 -bc_outflow 201

2723:   test:
2724:     suffix: shock_0
2725:     requires: p4est !single !complex
2726:     args: -ufv_vtk_interval 0 -monitor density,energy -f -grid_size 2,1 -grid_bounds -1,1.,0.,1 -bc_wall 1,2,3,4 -dm_type p4est -dm_forest_partition_overlap 1 -dm_forest_maximum_refinement 6 -dm_forest_minimum_refinement 2 -dm_forest_initial_refinement 2 -ufv_use_amr -refine_vec_tagger_box 0.5,inf -coarsen_vec_tagger_box 0,1.e-2 -refine_tag_view -coarsen_tag_view -physics euler -eu_type iv_shock -ufv_cfl 10 -eu_alpha 60. -grid_skew_60 -eu_gamma 1.4 -eu_amach 2.02 -eu_rho2 3. -petscfv_type leastsquares -petsclimiter_type minmod -petscfv_compute_gradients 0 -ts_max_time 0.5 -ts_ssp_type rks2 -ts_ssp_nstages 10 -ufv_vtk_basename ${wPETSC_DIR}/ex11
2727:     timeoutfactor: 3

2729:   # Test GLVis visualization of PetscFV fields
2730:   test:
2731:     suffix: glvis_adv_2d_tet
2732:     args: -ufv_vtk_interval 0 -ts_monitor_solution glvis: -ts_max_steps 0 -ufv_vtk_monitor 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_gmsh_periodic 0

2734:   test:
2735:     suffix: glvis_adv_2d_quad
2736:     args: -ufv_vtk_interval 0 -ts_monitor_solution glvis: -ts_max_steps 0 -ufv_vtk_monitor 0 -dm_refine 5 -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3

2738:   test:
2739:     suffix: tut_1
2740:     requires: exodusii
2741:     nsize: 1
2742:     args: -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo

2744:   test:
2745:     suffix: tut_2
2746:     requires: exodusii
2747:     nsize: 1
2748:     args: -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo -ts_type rosw

2750:   test:
2751:     suffix: tut_3
2752:     requires: exodusii
2753:     nsize: 4
2754:     args: -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -monitor Error -advect_sol_type bump -petscfv_type leastsquares -petsclimiter_type sin

2756:   test:
2757:     suffix: tut_4
2758:     requires: exodusii
2759:     nsize: 4
2760:     args: -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -physics sw -monitor Height,Energy -petscfv_type leastsquares -petsclimiter_type minmod

2762: TEST*/