Actual source code: ex11.c

petsc-3.9.4 2018-09-11
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> /* For SplitFaces() */

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

 46: static PetscFunctionList PhysicsList;

 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)(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];

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(PetscDS prob, Physics phys)
285: {
287:   const PetscInt inflowids[] = {100,200,300},outflowids[] = {101};

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

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

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

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

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

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

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

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

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

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

414: #if defined(PETSC_USE_COMPLEX)
415:   for (i = 0; i < 1+dim; i++) uLreal.vals[i] = PetscRealPart(xL[i]);
416:   for (i = 0; i < 1+dim; i++) uRreal.vals[i] = PetscRealPart(xR[i]);
417: #endif
418:   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"); */
419:   nn[0] = n[0];
420:   nn[1] = n[1];
421:   Normalize2Real(nn);
422:   SWFlux(phys,nn,uL,&(fL.swnode));
423:   SWFlux(phys,nn,uR,&(fR.swnode));
424:   cL    = PetscSqrtReal(sw->gravity*uL->h);
425:   cR    = PetscSqrtReal(sw->gravity*uR->h); /* gravity wave speed */
426:   speed = PetscMax(PetscAbsReal(Dot2Real(uL->uh,nn)/uL->h) + cL,PetscAbsReal(Dot2Real(uR->uh,nn)/uR->h) + cR);
427:   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);
428: }

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

435:   if (time != 0.0) SETERRQ1(mod->comm,PETSC_ERR_SUP,"No solution known for time %g",(double)time);
436:   dx[0] = x[0] - 1.5;
437:   dx[1] = x[1] - 1.0;
438:   r     = Norm2Real(dx);
439:   sigma = 0.5;
440:   u[0]  = 1 + 2*PetscExpReal(-PetscSqr(r)/(2*PetscSqr(sigma)));
441:   u[1]  = 0.0;
442:   u[2]  = 0.0;
443:   return(0);
444: }

446: static PetscErrorCode PhysicsFunctional_SW(Model mod,PetscReal time,const PetscReal *coord,const PetscScalar *xx,PetscReal *f,void *ctx)
447: {
448:   Physics      phys = (Physics)ctx;
449:   Physics_SW   *sw  = (Physics_SW*)phys->data;
450:   const SWNode *x   = (const SWNode*)xx;
451:   PetscReal  u[2];
452:   PetscReal    h;

455:   h = x->h;
456:   Scale2Real(1./x->h,x->uh,u);
457:   f[sw->functional.Height] = h;
458:   f[sw->functional.Speed]  = Norm2Real(u) + PetscSqrtReal(sw->gravity*h);
459:   f[sw->functional.Energy] = 0.5*(Dot2Real(x->uh,u) + sw->gravity*PetscSqr(h));
460:   return(0);
461: }

463: static PetscErrorCode SetUpBC_SW(PetscDS prob,Physics phys)
464: {
466:   const PetscInt wallids[] = {100,101,200,300};
468:   PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", "Face Sets", 0, 0, NULL, (void (*)(void)) PhysicsBoundary_SW_Wall, ALEN(wallids), wallids, phys);
469:   return(0);
470: }

472: static PetscErrorCode PhysicsCreate_SW(Model mod,Physics phys,PetscOptionItems *PetscOptionsObject)
473: {
474:   Physics_SW     *sw;

478:   phys->field_desc = PhysicsFields_SW;
479:   phys->riemann = (PetscRiemannFunc) PhysicsRiemann_SW;
480:   PetscNew(&sw);
481:   phys->data    = sw;
482:   mod->setupbc  = SetUpBC_SW;

484:   PetscOptionsHead(PetscOptionsObject,"SW options");
485:   {
486:     sw->gravity = 1.0;
487:     PetscOptionsReal("-sw_gravity","Gravitational constant","",sw->gravity,&sw->gravity,NULL);
488:   }
489:   PetscOptionsTail();
490:   phys->maxspeed = PetscSqrtReal(2.0*sw->gravity); /* Mach 1 for depth of 2 */

492:   ModelSolutionSetDefault(mod,PhysicsSolution_SW,phys);
493:   ModelFunctionalRegister(mod,"Height",&sw->functional.Height,PhysicsFunctional_SW,phys);
494:   ModelFunctionalRegister(mod,"Speed",&sw->functional.Speed,PhysicsFunctional_SW,phys);
495:   ModelFunctionalRegister(mod,"Energy",&sw->functional.Energy,PhysicsFunctional_SW,phys);

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

499:   return(0);
500: }

502: /******************* Euler Density Shock (EULER_IV_SHOCK,EULER_SS_SHOCK) ********************/
503: /* An initial-value and self-similar solutions of the compressible Euler equations */
504: /* Ravi Samtaney and D. I. Pullin */
505: /* Phys. Fluids 8, 2650 (1996); http://dx.doi.org/10.1063/1.869050 */
506: typedef enum {EULER_PAR_GAMMA,EULER_PAR_RHOR,EULER_PAR_AMACH,EULER_PAR_ITANA,EULER_PAR_SIZE} EulerParamIdx;
507: typedef enum {EULER_IV_SHOCK,EULER_SS_SHOCK,EULER_SHOCK_TUBE,EULER_LINEAR_WAVE} EulerType;
508: typedef struct {
509:   PetscReal r;
510:   PetscReal ru[DIM];
511:   PetscReal E;
512: } EulerNode;
513: typedef union {
514:   EulerNode eulernode;
515:   PetscReal vals[DIM+2];
516: } EulerNodeUnion;
517: typedef PetscErrorCode (*EquationOfState)(const PetscReal*, const EulerNode*, PetscReal*);
518: typedef struct {
519:   EulerType       type;
520:   PetscReal       pars[EULER_PAR_SIZE];
521:   EquationOfState sound;
522:   struct {
523:     PetscInt Density;
524:     PetscInt Momentum;
525:     PetscInt Energy;
526:     PetscInt Pressure;
527:     PetscInt Speed;
528:   } monitor;
529: } Physics_Euler;

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

533: /* initial condition */
534: int initLinearWave(EulerNode *ux, const PetscReal gamma, const PetscReal coord[], const PetscReal Lx);
535: static PetscErrorCode PhysicsSolution_Euler(Model mod, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx)
536: {
537:   PetscInt i;
538:   Physics         phys = (Physics)ctx;
539:   Physics_Euler   *eu  = (Physics_Euler*)phys->data;
540:   EulerNode       *uu  = (EulerNode*)u;
541:   PetscReal        p0,gamma,c;
543:   if (time != 0.0) SETERRQ1(mod->comm,PETSC_ERR_SUP,"No solution known for time %g",(double)time);

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

549:   if (eu->type==EULER_IV_SHOCK || eu->type==EULER_SS_SHOCK) {
550:     /******************* Euler Density Shock ********************/
551:     /* On initial-value and self-similar solutions of the compressible Euler equations */
552:     /* Ravi Samtaney and D. I. Pullin */
553:     /* Phys. Fluids 8, 2650 (1996); http://dx.doi.org/10.1063/1.869050 */
554:     /* initial conditions 1: left of shock, 0: left of discontinuity 2: right of discontinuity,  */
555:     p0 = 1.;
556:     if (x[0] < 0.0 + x[1]*eu->pars[EULER_PAR_ITANA]) {
557:       if (x[0] < mod->bounds[0]*0.5) { /* left of shock (1) */
558:         PetscReal amach,rho,press,gas1,p1;
559:         amach = eu->pars[EULER_PAR_AMACH];
560:         rho = 1.;
561:         press = p0;
562:         p1 = press*(1.0+2.0*gamma/(gamma+1.0)*(amach*amach-1.0));
563:         gas1 = (gamma-1.0)/(gamma+1.0);
564:         uu->r = rho*(p1/press+gas1)/(gas1*p1/press+1.0);
565:         uu->ru[0]   = ((uu->r - rho)*PetscSqrtReal(gamma*press/rho)*amach);
566:         uu->E = p1/(gamma-1.0) + .5/uu->r*uu->ru[0]*uu->ru[0];
567:       }
568:       else { /* left of discontinuity (0) */
569:         uu->r = 1.; /* rho = 1 */
570:         uu->E = p0/(gamma-1.0);
571:       }
572:     }
573:     else { /* right of discontinuity (2) */
574:       uu->r = eu->pars[EULER_PAR_RHOR];
575:       uu->E = p0/(gamma-1.0);
576:     }
577:   }
578:   else if (eu->type==EULER_SHOCK_TUBE) {
579:     /* 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. */
580:     if (x[0] < 0.0 ) {
581:       uu->r = 8.;
582:       uu->E = 10./(gamma-1.);
583:     }
584:     else {
585:       uu->r = 1.;
586:       uu->E = 1./(gamma-1.);
587:     }
588:   }
589:   else if (eu->type==EULER_LINEAR_WAVE) {
590:     initLinearWave( uu, gamma, x, mod->bounds[1] - mod->bounds[0]);
591:   }
592:   else SETERRQ1(mod->comm,PETSC_ERR_SUP,"Unknown type %d",eu->type);

594:   /* set phys->maxspeed: (mod->maxspeed = phys->maxspeed) in main; */
595:   eu->sound(&gamma,uu,&c);
596:   c = (uu->ru[0]/uu->r) + c;
597:   if (c > phys->maxspeed) phys->maxspeed = c;

599:   return(0);
600: }

602: static PetscErrorCode Pressure_PG(const PetscReal gamma,const EulerNode *x,PetscReal *p)
603: {
604:   PetscReal ru2;

607:   ru2  = DotDIMReal(x->ru,x->ru);
608:   (*p)=(x->E - 0.5*ru2/x->r)*(gamma - 1.0); /* (E - rho V^2/2)(gamma-1) = e rho (gamma-1) */
609:   return(0);
610: }

612: static PetscErrorCode SpeedOfSound_PG(const PetscReal *gamma, const EulerNode *x, PetscReal *c)
613: {
614:   PetscReal p;

617:   Pressure_PG(*gamma,x,&p);
618:   if (p<0.) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"negative pressure time %g -- NEED TO FIX!!!!!!",(double) p);
619:   /* pars[EULER_PAR_GAMMA] = heat capacity ratio */
620:   (*c)=PetscSqrtReal(*gamma * p / x->r);
621:   return(0);
622: }

624: /*
625:  * x = (rho,rho*(u_1),...,rho*e)^T
626:  * x_t+div(f_1(x))+...+div(f_DIM(x)) = 0
627:  *
628:  * f_i(x) = u_i*x+(0,0,...,p,...,p*u_i)^T
629:  *
630:  */
631: static PetscErrorCode EulerFlux(Physics phys,const PetscReal *n,const EulerNode *x,EulerNode *f)
632: {
633:   Physics_Euler *eu = (Physics_Euler*)phys->data;
634:   PetscReal     nu,p;
635:   PetscInt      i;

638:   Pressure_PG(eu->pars[EULER_PAR_GAMMA],x,&p);
639:   nu = DotDIMReal(x->ru,n);
640:   f->r = nu;   /* A rho u */
641:   nu /= x->r;  /* A u */
642:   for (i=0; i<DIM; i++) f->ru[i] = nu * x->ru[i] + n[i]*p;  /* r u^2 + p */
643:   f->E = nu * (x->E + p); /* u(e+p) */
644:   return(0);
645: }

647: /* PetscReal* => EulerNode* conversion */
648: static PetscErrorCode PhysicsBoundary_Euler_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *a_xI, PetscScalar *a_xG, void *ctx)
649: {
650:   PetscInt    i;
651:   const EulerNode *xI = (const EulerNode*)a_xI;
652:   EulerNode       *xG = (EulerNode*)a_xG;
653:   Physics         phys = (Physics)ctx;
654:   Physics_Euler   *eu  = (Physics_Euler*)phys->data;
656:   xG->r = xI->r;           /* ghost cell density - same */
657:   xG->E = xI->E;           /* ghost cell energy - same */
658:   if (n[1] != 0.) {        /* top and bottom */
659:     xG->ru[0] =  xI->ru[0]; /* copy tang to wall */
660:     xG->ru[1] = -xI->ru[1]; /* reflect perp to t/b wall */
661:   }
662:   else { /* sides */
663:     for (i=0; i<DIM; i++) xG->ru[i] = xI->ru[i]; /* copy */
664:   }
665:   if (eu->type == EULER_LINEAR_WAVE) { /* debug */
666: #if 0
667:     PetscPrintf(PETSC_COMM_WORLD,"%s coord=%g,%g\n",PETSC_FUNCTION_NAME,c[0],c[1]);
668: #endif
669:   }
670:   return(0);
671: }
672: int godunovflux( const PetscScalar *ul, const PetscScalar *ur, PetscScalar *flux, const PetscReal *nn, const int *ndim, const PetscReal *gamma);
673: /* PetscReal* => EulerNode* conversion */
674: static void PhysicsRiemann_Euler_Godunov( PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n,
675:                                           const PetscScalar *xL, const PetscScalar *xR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, Physics phys)
676: {
677:   Physics_Euler   *eu = (Physics_Euler*)phys->data;
678:   PetscReal       cL,cR,speed,velL,velR,nn[DIM],s2;
679:   PetscInt        i;
680:   PetscErrorCode  ierr;

683:   for (i=0,s2=0.; i<DIM; i++) {
684:     nn[i] = n[i];
685:     s2 += nn[i]*nn[i];
686:   }
687:   s2 = PetscSqrtReal(s2); /* |n|_2 = sum(n^2)^1/2 */
688:   for (i=0.; i<DIM; i++) nn[i] /= s2;
689:   if (0) { /* Rusanov */
690:     const EulerNode *uL = (const EulerNode*)xL,*uR = (const EulerNode*)xR;
691:     EulerNodeUnion  fL,fR;
692:     EulerFlux(phys,nn,uL,&(fL.eulernode));
693:     EulerFlux(phys,nn,uR,&(fR.eulernode));
694:     eu->sound(&eu->pars[EULER_PAR_GAMMA],uL,&cL);if (ierr) exit(13);
695:     eu->sound(&eu->pars[EULER_PAR_GAMMA],uR,&cR);if (ierr) exit(14);
696:     velL = DotDIMReal(uL->ru,nn)/uL->r;
697:     velR = DotDIMReal(uR->ru,nn)/uR->r;
698:     speed = PetscMax(velR + cR, velL + cL);
699:     for (i=0; i<2+dim; i++) flux[i] = 0.5*((fL.vals[i]+fR.vals[i]) + speed*(xL[i] - xR[i]))*s2;
700:   }
701:   else {
702:     int dim = DIM;
703:     /* int iwave =  */
704:     godunovflux(xL, xR, flux, nn, &dim, &eu->pars[EULER_PAR_GAMMA]);
705:     for (i=0; i<2+dim; i++) flux[i] *= s2;
706:   }
707:   PetscFunctionReturnVoid();
708: }

710: static PetscErrorCode PhysicsFunctional_Euler(Model mod,PetscReal time,const PetscReal *coord,const PetscScalar *xx,PetscReal *f,void *ctx)
711: {
712:   Physics         phys = (Physics)ctx;
713:   Physics_Euler   *eu  = (Physics_Euler*)phys->data;
714:   const EulerNode *x   = (const EulerNode*)xx;
715:   PetscReal       p;

718:   f[eu->monitor.Density]  = x->r;
719:   f[eu->monitor.Momentum] = NormDIM(x->ru);
720:   f[eu->monitor.Energy]   = x->E;
721:   f[eu->monitor.Speed]    = NormDIM(x->ru)/x->r;
722:   Pressure_PG(eu->pars[EULER_PAR_GAMMA], x, &p);
723:   f[eu->monitor.Pressure] = p;
724:   return(0);
725: }

727: static PetscErrorCode SetUpBC_Euler(PetscDS prob,Physics phys)
728: {
729:   PetscErrorCode  ierr;
730:   Physics_Euler   *eu = (Physics_Euler *) phys->data;
731:   if (eu->type == EULER_LINEAR_WAVE) {
732:     const PetscInt wallids[] = {100,101};
733:     PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", "Face Sets", 0, 0, NULL, (void (*)(void)) PhysicsBoundary_Euler_Wall, ALEN(wallids), wallids, phys);
734:   }
735:   else {
736:     const PetscInt wallids[] = {100,101,200,300};
737:     PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", "Face Sets", 0, 0, NULL, (void (*)(void)) PhysicsBoundary_Euler_Wall, ALEN(wallids), wallids, phys);
738:   }
739:   return(0);
740: }

742: static PetscErrorCode PhysicsCreate_Euler(Model mod,Physics phys,PetscOptionItems *PetscOptionsObject)
743: {
744:   Physics_Euler   *eu;
745:   PetscErrorCode  ierr;

748:   phys->field_desc = PhysicsFields_Euler;
749:   phys->riemann = (PetscRiemannFunc) PhysicsRiemann_Euler_Godunov;
750:   PetscNew(&eu);
751:   phys->data    = eu;
752:   mod->setupbc = SetUpBC_Euler;
753:   PetscOptionsHead(PetscOptionsObject,"Euler options");
754:   {
755:     PetscReal alpha;
756:     char type[64] = "linear_wave";
757:     PetscBool  is;
758:     mod->bcs[0] = mod->bcs[1] = mod->bcs[2] = DM_BOUNDARY_GHOSTED;
759:     eu->pars[EULER_PAR_GAMMA] = 1.4;
760:     eu->pars[EULER_PAR_AMACH] = 2.02;
761:     eu->pars[EULER_PAR_RHOR] = 3.0;
762:     eu->pars[EULER_PAR_ITANA] = 0.57735026918963; /* angle of Euler self similar (SS) shock */
763:     PetscOptionsReal("-eu_gamma","Heat capacity ratio","",eu->pars[EULER_PAR_GAMMA],&eu->pars[EULER_PAR_GAMMA],NULL);
764:     PetscOptionsReal("-eu_amach","Shock speed (Mach)","",eu->pars[EULER_PAR_AMACH],&eu->pars[EULER_PAR_AMACH],NULL);
765:     PetscOptionsReal("-eu_rho2","Density right of discontinuity","",eu->pars[EULER_PAR_RHOR],&eu->pars[EULER_PAR_RHOR],NULL);
766:     alpha = 60.;
767:     PetscOptionsReal("-eu_alpha","Angle of discontinuity","",alpha,&alpha,NULL);
768:     if (alpha<=0. || alpha>90.) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Alpha bust be > 0 and <= 90 (%g)",alpha);
769:     eu->pars[EULER_PAR_ITANA] = 1./PetscTanReal( alpha * PETSC_PI / 180.0 );
770:     PetscOptionsString("-eu_type","Type of Euler test","",type,type,sizeof(type),NULL);
771:     PetscStrcmp(type,"linear_wave", &is);
772:     if (is) {
773:       eu->type = EULER_LINEAR_WAVE;
774:       mod->bcs[0] = mod->bcs[1] = mod->bcs[2] = DM_BOUNDARY_PERIODIC;
775:       mod->bcs[1] = DM_BOUNDARY_GHOSTED; /* debug */
776:       PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"linear_wave");
777:     }
778:     else {
779:       if (DIM != 2) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DIM must be 2 unless linear wave test %s",type);
780:       PetscStrcmp(type,"iv_shock", &is);
781:       if (is) {
782:         eu->type = EULER_IV_SHOCK;
783:         PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"iv_shock");
784:       }
785:       else {
786:         PetscStrcmp(type,"ss_shock", &is);
787:         if (is) {
788:           eu->type = EULER_SS_SHOCK;
789:           PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"ss_shock");
790:         }
791:         else {
792:           PetscStrcmp(type,"shock_tube", &is);
793:           if (is) eu->type = EULER_SHOCK_TUBE;
794:           else SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Unknown Euler type %s",type);
795:           PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"shock_tube");
796:         }
797:       }
798:     }
799:   }
800:   PetscOptionsTail();
801:   eu->sound = SpeedOfSound_PG;
802:   phys->maxspeed = 0.; /* will get set in solution */
803:   ModelSolutionSetDefault(mod,PhysicsSolution_Euler,phys);
804:   ModelFunctionalRegister(mod,"Speed",&eu->monitor.Speed,PhysicsFunctional_Euler,phys);
805:   ModelFunctionalRegister(mod,"Energy",&eu->monitor.Energy,PhysicsFunctional_Euler,phys);
806:   ModelFunctionalRegister(mod,"Density",&eu->monitor.Density,PhysicsFunctional_Euler,phys);
807:   ModelFunctionalRegister(mod,"Momentum",&eu->monitor.Momentum,PhysicsFunctional_Euler,phys);
808:   ModelFunctionalRegister(mod,"Pressure",&eu->monitor.Pressure,PhysicsFunctional_Euler,phys);

810:   return(0);
811: }

813: static PetscErrorCode ErrorIndicator_Simple(PetscInt dim, PetscReal volume, PetscInt numComps, const PetscScalar u[], const PetscScalar grad[], PetscReal *error, void *ctx)
814: {
815:   PetscReal      err = 0.;
816:   PetscInt       i, j;

819:   for (i = 0; i < numComps; i++) {
820:     for (j = 0; j < dim; j++) {
821:       err += PetscSqr(PetscRealPart(grad[i * dim + j]));
822:     }
823:   }
824:   *error = volume * err;
825:   return(0);
826: }

828: PetscErrorCode ConstructCellBoundary(DM dm, User user)
829: {
830:   const char     *name   = "Cell Sets";
831:   const char     *bdname = "split faces";
832:   IS             regionIS, innerIS;
833:   const PetscInt *regions, *cells;
834:   PetscInt       numRegions, innerRegion, numCells, c;
835:   PetscInt       cStart, cEnd, cEndInterior, fStart, fEnd;
836:   PetscBool      hasLabel;

840:   DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
841:   DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);
842:   DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);

844:   DMHasLabel(dm, name, &hasLabel);
845:   if (!hasLabel) return(0);
846:   DMGetLabelSize(dm, name, &numRegions);
847:   if (numRegions != 2) return(0);
848:   /* Get the inner id */
849:   DMGetLabelIdIS(dm, name, &regionIS);
850:   ISGetIndices(regionIS, &regions);
851:   innerRegion = regions[0];
852:   ISRestoreIndices(regionIS, &regions);
853:   ISDestroy(&regionIS);
854:   /* Find the faces between cells in different regions, could call DMPlexCreateNeighborCSR() */
855:   DMGetStratumIS(dm, name, innerRegion, &innerIS);
856:   ISGetLocalSize(innerIS, &numCells);
857:   ISGetIndices(innerIS, &cells);
858:   DMCreateLabel(dm, bdname);
859:   for (c = 0; c < numCells; ++c) {
860:     const PetscInt cell = cells[c];
861:     const PetscInt *faces;
862:     PetscInt       numFaces, f;

864:     if ((cell < cStart) || (cell >= cEnd)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Got invalid point %d which is not a cell", cell);
865:     DMPlexGetConeSize(dm, cell, &numFaces);
866:     DMPlexGetCone(dm, cell, &faces);
867:     for (f = 0; f < numFaces; ++f) {
868:       const PetscInt face = faces[f];
869:       const PetscInt *neighbors;
870:       PetscInt       nC, regionA, regionB;

872:       if ((face < fStart) || (face >= fEnd)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Got invalid point %d which is not a face", face);
873:       DMPlexGetSupportSize(dm, face, &nC);
874:       if (nC != 2) continue;
875:       DMPlexGetSupport(dm, face, &neighbors);
876:       if ((neighbors[0] >= cEndInterior) || (neighbors[1] >= cEndInterior)) continue;
877:       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]);
878:       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]);
879:       DMGetLabelValue(dm, name, neighbors[0], &regionA);
880:       DMGetLabelValue(dm, name, neighbors[1], &regionB);
881:       if (regionA < 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid label %s: Cell %d has no value", name, neighbors[0]);
882:       if (regionB < 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid label %s: Cell %d has no value", name, neighbors[1]);
883:       if (regionA != regionB) {
884:         DMSetLabelValue(dm, bdname, faces[f], 1);
885:       }
886:     }
887:   }
888:   ISRestoreIndices(innerIS, &cells);
889:   ISDestroy(&innerIS);
890:   {
891:     DMLabel label;

893:     DMGetLabel(dm, bdname, &label);
894:     DMLabelView(label, PETSC_VIEWER_STDOUT_WORLD);
895:   }
896:   return(0);
897: }

899: /* Right now, I have just added duplicate faces, which see both cells. We can
900: - Add duplicate vertices and decouple the face cones
901: - Disconnect faces from cells across the rotation gap
902: */
903: PetscErrorCode SplitFaces(DM *dmSplit, const char labelName[], User user)
904: {
905:   DM             dm = *dmSplit, sdm;
906:   PetscSF        sfPoint, gsfPoint;
907:   PetscSection   coordSection, newCoordSection;
908:   Vec            coordinates;
909:   IS             idIS;
910:   const PetscInt *ids;
911:   PetscInt       *newpoints;
912:   PetscInt       dim, depth, maxConeSize, maxSupportSize, numLabels, numGhostCells;
913:   PetscInt       numFS, fs, pStart, pEnd, p, cEnd, cEndInterior, vStart, vEnd, v, fStart, fEnd, newf, d, l;
914:   PetscBool      hasLabel;

918:   DMHasLabel(dm, labelName, &hasLabel);
919:   if (!hasLabel) return(0);
920:   DMCreate(PetscObjectComm((PetscObject)dm), &sdm);
921:   DMSetType(sdm, DMPLEX);
922:   DMGetDimension(dm, &dim);
923:   DMSetDimension(sdm, dim);

925:   DMGetLabelIdIS(dm, labelName, &idIS);
926:   ISGetLocalSize(idIS, &numFS);
927:   ISGetIndices(idIS, &ids);

929:   user->numSplitFaces = 0;
930:   for (fs = 0; fs < numFS; ++fs) {
931:     PetscInt numBdFaces;

933:     DMGetStratumSize(dm, labelName, ids[fs], &numBdFaces);
934:     user->numSplitFaces += numBdFaces;
935:   }
936:   DMPlexGetChart(dm, &pStart, &pEnd);
937:   pEnd += user->numSplitFaces;
938:   DMPlexSetChart(sdm, pStart, pEnd);
939:   DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);
940:   DMPlexSetHybridBounds(sdm, cEndInterior, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE);
941:   DMPlexGetHeightStratum(dm, 0, NULL, &cEnd);
942:   numGhostCells = cEnd - cEndInterior;
943:   /* Set cone and support sizes */
944:   DMPlexGetDepth(dm, &depth);
945:   for (d = 0; d <= depth; ++d) {
946:     DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);
947:     for (p = pStart; p < pEnd; ++p) {
948:       PetscInt newp = p;
949:       PetscInt size;

951:       DMPlexGetConeSize(dm, p, &size);
952:       DMPlexSetConeSize(sdm, newp, size);
953:       DMPlexGetSupportSize(dm, p, &size);
954:       DMPlexSetSupportSize(sdm, newp, size);
955:     }
956:   }
957:   DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);
958:   for (fs = 0, newf = fEnd; fs < numFS; ++fs) {
959:     IS             faceIS;
960:     const PetscInt *faces;
961:     PetscInt       numFaces, f;

963:     DMGetStratumIS(dm, labelName, ids[fs], &faceIS);
964:     ISGetLocalSize(faceIS, &numFaces);
965:     ISGetIndices(faceIS, &faces);
966:     for (f = 0; f < numFaces; ++f, ++newf) {
967:       PetscInt size;

969:       /* Right now I think that both faces should see both cells */
970:       DMPlexGetConeSize(dm, faces[f], &size);
971:       DMPlexSetConeSize(sdm, newf, size);
972:       DMPlexGetSupportSize(dm, faces[f], &size);
973:       DMPlexSetSupportSize(sdm, newf, size);
974:     }
975:     ISRestoreIndices(faceIS, &faces);
976:     ISDestroy(&faceIS);
977:   }
978:   DMSetUp(sdm);
979:   /* Set cones and supports */
980:   DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);
981:   PetscMalloc1(PetscMax(maxConeSize, maxSupportSize), &newpoints);
982:   DMPlexGetChart(dm, &pStart, &pEnd);
983:   for (p = pStart; p < pEnd; ++p) {
984:     const PetscInt *points, *orientations;
985:     PetscInt       size, i, newp = p;

987:     DMPlexGetConeSize(dm, p, &size);
988:     DMPlexGetCone(dm, p, &points);
989:     DMPlexGetConeOrientation(dm, p, &orientations);
990:     for (i = 0; i < size; ++i) newpoints[i] = points[i];
991:     DMPlexSetCone(sdm, newp, newpoints);
992:     DMPlexSetConeOrientation(sdm, newp, orientations);
993:     DMPlexGetSupportSize(dm, p, &size);
994:     DMPlexGetSupport(dm, p, &points);
995:     for (i = 0; i < size; ++i) newpoints[i] = points[i];
996:     DMPlexSetSupport(sdm, newp, newpoints);
997:   }
998:   PetscFree(newpoints);
999:   for (fs = 0, newf = fEnd; fs < numFS; ++fs) {
1000:     IS             faceIS;
1001:     const PetscInt *faces;
1002:     PetscInt       numFaces, f;

1004:     DMGetStratumIS(dm, labelName, ids[fs], &faceIS);
1005:     ISGetLocalSize(faceIS, &numFaces);
1006:     ISGetIndices(faceIS, &faces);
1007:     for (f = 0; f < numFaces; ++f, ++newf) {
1008:       const PetscInt *points;

1010:       DMPlexGetCone(dm, faces[f], &points);
1011:       DMPlexSetCone(sdm, newf, points);
1012:       DMPlexGetSupport(dm, faces[f], &points);
1013:       DMPlexSetSupport(sdm, newf, points);
1014:     }
1015:     ISRestoreIndices(faceIS, &faces);
1016:     ISDestroy(&faceIS);
1017:   }
1018:   ISRestoreIndices(idIS, &ids);
1019:   ISDestroy(&idIS);
1020:   DMPlexStratify(sdm);
1021:   /* Convert coordinates */
1022:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
1023:   DMGetCoordinateSection(dm, &coordSection);
1024:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &newCoordSection);
1025:   PetscSectionSetNumFields(newCoordSection, 1);
1026:   PetscSectionSetFieldComponents(newCoordSection, 0, dim);
1027:   PetscSectionSetChart(newCoordSection, vStart, vEnd);
1028:   for (v = vStart; v < vEnd; ++v) {
1029:     PetscSectionSetDof(newCoordSection, v, dim);
1030:     PetscSectionSetFieldDof(newCoordSection, v, 0, dim);
1031:   }
1032:   PetscSectionSetUp(newCoordSection);
1033:   DMSetCoordinateSection(sdm, PETSC_DETERMINE, newCoordSection);
1034:   PetscSectionDestroy(&newCoordSection); /* relinquish our reference */
1035:   DMGetCoordinatesLocal(dm, &coordinates);
1036:   DMSetCoordinatesLocal(sdm, coordinates);
1037:   /* Convert labels */
1038:   DMGetNumLabels(dm, &numLabels);
1039:   for (l = 0; l < numLabels; ++l) {
1040:     const char *lname;
1041:     PetscBool  isDepth;

1043:     DMGetLabelName(dm, l, &lname);
1044:     PetscStrcmp(lname, "depth", &isDepth);
1045:     if (isDepth) continue;
1046:     DMCreateLabel(sdm, lname);
1047:     DMGetLabelIdIS(dm, lname, &idIS);
1048:     ISGetLocalSize(idIS, &numFS);
1049:     ISGetIndices(idIS, &ids);
1050:     for (fs = 0; fs < numFS; ++fs) {
1051:       IS             pointIS;
1052:       const PetscInt *points;
1053:       PetscInt       numPoints;

1055:       DMGetStratumIS(dm, lname, ids[fs], &pointIS);
1056:       ISGetLocalSize(pointIS, &numPoints);
1057:       ISGetIndices(pointIS, &points);
1058:       for (p = 0; p < numPoints; ++p) {
1059:         PetscInt newpoint = points[p];

1061:         DMSetLabelValue(sdm, lname, newpoint, ids[fs]);
1062:       }
1063:       ISRestoreIndices(pointIS, &points);
1064:       ISDestroy(&pointIS);
1065:     }
1066:     ISRestoreIndices(idIS, &ids);
1067:     ISDestroy(&idIS);
1068:   }
1069:   {
1070:     /* Convert pointSF */
1071:     const PetscSFNode *remotePoints;
1072:     PetscSFNode       *gremotePoints;
1073:     const PetscInt    *localPoints;
1074:     PetscInt          *glocalPoints,*newLocation,*newRemoteLocation;
1075:     PetscInt          numRoots, numLeaves;
1076:     PetscMPIInt       size;

1078:     MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size);
1079:     DMGetPointSF(dm, &sfPoint);
1080:     DMGetPointSF(sdm, &gsfPoint);
1081:     DMPlexGetChart(dm,&pStart,&pEnd);
1082:     PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);
1083:     if (numRoots >= 0) {
1084:       PetscMalloc2(numRoots,&newLocation,pEnd-pStart,&newRemoteLocation);
1085:       for (l=0; l<numRoots; l++) newLocation[l] = l; /* + (l >= cEnd ? numGhostCells : 0); */
1086:       PetscSFBcastBegin(sfPoint, MPIU_INT, newLocation, newRemoteLocation);
1087:       PetscSFBcastEnd(sfPoint, MPIU_INT, newLocation, newRemoteLocation);
1088:       PetscMalloc1(numLeaves,    &glocalPoints);
1089:       PetscMalloc1(numLeaves, &gremotePoints);
1090:       for (l = 0; l < numLeaves; ++l) {
1091:         glocalPoints[l]        = localPoints[l]; /* localPoints[l] >= cEnd ? localPoints[l] + numGhostCells : localPoints[l]; */
1092:         gremotePoints[l].rank  = remotePoints[l].rank;
1093:         gremotePoints[l].index = newRemoteLocation[localPoints[l]];
1094:       }
1095:       PetscFree2(newLocation,newRemoteLocation);
1096:       PetscSFSetGraph(gsfPoint, numRoots+numGhostCells, numLeaves, glocalPoints, PETSC_OWN_POINTER, gremotePoints, PETSC_OWN_POINTER);
1097:     }
1098:     DMDestroy(dmSplit);
1099:     *dmSplit = sdm;
1100:   }
1101:   return(0);
1102: }

1104: PetscErrorCode CreatePartitionVec(DM dm, DM *dmCell, Vec *partition)
1105: {
1106:   PetscSF        sfPoint;
1107:   PetscSection   coordSection;
1108:   Vec            coordinates;
1109:   PetscSection   sectionCell;
1110:   PetscScalar    *part;
1111:   PetscInt       cStart, cEnd, c;
1112:   PetscMPIInt    rank;

1116:   DMGetCoordinateSection(dm, &coordSection);
1117:   DMGetCoordinatesLocal(dm, &coordinates);
1118:   DMClone(dm, dmCell);
1119:   DMGetPointSF(dm, &sfPoint);
1120:   DMSetPointSF(*dmCell, sfPoint);
1121:   DMSetCoordinateSection(*dmCell, PETSC_DETERMINE, coordSection);
1122:   DMSetCoordinatesLocal(*dmCell, coordinates);
1123:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
1124:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &sectionCell);
1125:   DMPlexGetHeightStratum(*dmCell, 0, &cStart, &cEnd);
1126:   PetscSectionSetChart(sectionCell, cStart, cEnd);
1127:   for (c = cStart; c < cEnd; ++c) {
1128:     PetscSectionSetDof(sectionCell, c, 1);
1129:   }
1130:   PetscSectionSetUp(sectionCell);
1131:   DMSetDefaultSection(*dmCell, sectionCell);
1132:   PetscSectionDestroy(&sectionCell);
1133:   DMCreateLocalVector(*dmCell, partition);
1134:   PetscObjectSetName((PetscObject)*partition, "partition");
1135:   VecGetArray(*partition, &part);
1136:   for (c = cStart; c < cEnd; ++c) {
1137:     PetscScalar *p;

1139:     DMPlexPointLocalRef(*dmCell, c, part, &p);
1140:     p[0] = rank;
1141:   }
1142:   VecRestoreArray(*partition, &part);
1143:   return(0);
1144: }

1146: PetscErrorCode CreateMassMatrix(DM dm, Vec *massMatrix, User user)
1147: {
1148:   DM                dmMass, dmFace, dmCell, dmCoord;
1149:   PetscSection      coordSection;
1150:   Vec               coordinates, facegeom, cellgeom;
1151:   PetscSection      sectionMass;
1152:   PetscScalar       *m;
1153:   const PetscScalar *fgeom, *cgeom, *coords;
1154:   PetscInt          vStart, vEnd, v;
1155:   PetscErrorCode    ierr;

1158:   DMGetCoordinateSection(dm, &coordSection);
1159:   DMGetCoordinatesLocal(dm, &coordinates);
1160:   DMClone(dm, &dmMass);
1161:   DMSetCoordinateSection(dmMass, PETSC_DETERMINE, coordSection);
1162:   DMSetCoordinatesLocal(dmMass, coordinates);
1163:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &sectionMass);
1164:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
1165:   PetscSectionSetChart(sectionMass, vStart, vEnd);
1166:   for (v = vStart; v < vEnd; ++v) {
1167:     PetscInt numFaces;

1169:     DMPlexGetSupportSize(dmMass, v, &numFaces);
1170:     PetscSectionSetDof(sectionMass, v, numFaces*numFaces);
1171:   }
1172:   PetscSectionSetUp(sectionMass);
1173:   DMSetDefaultSection(dmMass, sectionMass);
1174:   PetscSectionDestroy(&sectionMass);
1175:   DMGetLocalVector(dmMass, massMatrix);
1176:   VecGetArray(*massMatrix, &m);
1177:   DMPlexTSGetGeometryFVM(dm, &facegeom, &cellgeom, NULL);
1178:   VecGetDM(facegeom, &dmFace);
1179:   VecGetArrayRead(facegeom, &fgeom);
1180:   VecGetDM(cellgeom, &dmCell);
1181:   VecGetArrayRead(cellgeom, &cgeom);
1182:   DMGetCoordinateDM(dm, &dmCoord);
1183:   VecGetArrayRead(coordinates, &coords);
1184:   for (v = vStart; v < vEnd; ++v) {
1185:     const PetscInt        *faces;
1186:     PetscFVFaceGeom       *fgA, *fgB, *cg;
1187:     PetscScalar           *vertex;
1188:     PetscInt               numFaces, sides[2], f, g;

1190:     DMPlexPointLocalRead(dmCoord, v, coords, &vertex);
1191:     DMPlexGetSupportSize(dmMass, v, &numFaces);
1192:     DMPlexGetSupport(dmMass, v, &faces);
1193:     for (f = 0; f < numFaces; ++f) {
1194:       sides[0] = faces[f];
1195:       DMPlexPointLocalRead(dmFace, faces[f], fgeom, &fgA);
1196:       for (g = 0; g < numFaces; ++g) {
1197:         const PetscInt *cells = NULL;
1198:         PetscReal      area   = 0.0;
1199:         PetscInt       numCells;

1201:         sides[1] = faces[g];
1202:         DMPlexPointLocalRead(dmFace, faces[g], fgeom, &fgB);
1203:         DMPlexGetJoin(dmMass, 2, sides, &numCells, &cells);
1204:         if (numCells != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Invalid join for faces");
1205:         DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cg);
1206:         area += PetscAbsScalar((vertex[0] - cg->centroid[0])*(fgA->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1])*(fgA->centroid[0] - cg->centroid[0]));
1207:         area += PetscAbsScalar((vertex[0] - cg->centroid[0])*(fgB->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1])*(fgB->centroid[0] - cg->centroid[0]));
1208:         m[f*numFaces+g] = Dot2Real(fgA->normal, fgB->normal)*area*0.5;
1209:         DMPlexRestoreJoin(dmMass, 2, sides, &numCells, &cells);
1210:       }
1211:     }
1212:   }
1213:   VecRestoreArrayRead(facegeom, &fgeom);
1214:   VecRestoreArrayRead(cellgeom, &cgeom);
1215:   VecRestoreArrayRead(coordinates, &coords);
1216:   VecRestoreArray(*massMatrix, &m);
1217:   DMDestroy(&dmMass);
1218:   return(0);
1219: }

1221: /* Behavior will be different for multi-physics or when using non-default boundary conditions */
1222: static PetscErrorCode ModelSolutionSetDefault(Model mod,SolutionFunction func,void *ctx)
1223: {
1225:   mod->solution    = func;
1226:   mod->solutionctx = ctx;
1227:   return(0);
1228: }

1230: static PetscErrorCode ModelFunctionalRegister(Model mod,const char *name,PetscInt *offset,FunctionalFunction func,void *ctx)
1231: {
1233:   FunctionalLink link,*ptr;
1234:   PetscInt       lastoffset = -1;

1237:   for (ptr=&mod->functionalRegistry; *ptr; ptr = &(*ptr)->next) lastoffset = (*ptr)->offset;
1238:   PetscNew(&link);
1239:   PetscStrallocpy(name,&link->name);
1240:   link->offset = lastoffset + 1;
1241:   link->func   = func;
1242:   link->ctx    = ctx;
1243:   link->next   = NULL;
1244:   *ptr         = link;
1245:   *offset      = link->offset;
1246:   return(0);
1247: }

1249: static PetscErrorCode ModelFunctionalSetFromOptions(Model mod,PetscOptionItems *PetscOptionsObject)
1250: {
1252:   PetscInt       i,j;
1253:   FunctionalLink link;
1254:   char           *names[256];

1257:   mod->numMonitored = ALEN(names);
1258:   PetscOptionsStringArray("-monitor","list of functionals to monitor","",names,&mod->numMonitored,NULL);
1259:   /* Create list of functionals that will be computed somehow */
1260:   PetscMalloc1(mod->numMonitored,&mod->functionalMonitored);
1261:   /* Create index of calls that we will have to make to compute these functionals (over-allocation in general). */
1262:   PetscMalloc1(mod->numMonitored,&mod->functionalCall);
1263:   mod->numCall = 0;
1264:   for (i=0; i<mod->numMonitored; i++) {
1265:     for (link=mod->functionalRegistry; link; link=link->next) {
1266:       PetscBool match;
1267:       PetscStrcasecmp(names[i],link->name,&match);
1268:       if (match) break;
1269:     }
1270:     if (!link) SETERRQ1(mod->comm,PETSC_ERR_USER,"No known functional '%s'",names[i]);
1271:     mod->functionalMonitored[i] = link;
1272:     for (j=0; j<i; j++) {
1273:       if (mod->functionalCall[j]->func == link->func && mod->functionalCall[j]->ctx == link->ctx) goto next_name;
1274:     }
1275:     mod->functionalCall[mod->numCall++] = link; /* Just points to the first link using the result. There may be more results. */
1276: next_name:
1277:     PetscFree(names[i]);
1278:   }

1280:   /* Find out the maximum index of any functional computed by a function we will be calling (even if we are not using it) */
1281:   mod->maxComputed = -1;
1282:   for (link=mod->functionalRegistry; link; link=link->next) {
1283:     for (i=0; i<mod->numCall; i++) {
1284:       FunctionalLink call = mod->functionalCall[i];
1285:       if (link->func == call->func && link->ctx == call->ctx) {
1286:         mod->maxComputed = PetscMax(mod->maxComputed,link->offset);
1287:       }
1288:     }
1289:   }
1290:   return(0);
1291: }

1293: static PetscErrorCode FunctionalLinkDestroy(FunctionalLink *link)
1294: {
1296:   FunctionalLink l,next;

1299:   if (!link) return(0);
1300:   l     = *link;
1301:   *link = NULL;
1302:   for (; l; l=next) {
1303:     next = l->next;
1304:     PetscFree(l->name);
1305:     PetscFree(l);
1306:   }
1307:   return(0);
1308: }

1310: /* put the solution callback into a functional callback */
1311: static PetscErrorCode SolutionFunctional(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *modctx)
1312: {
1313:   Model          mod;
1316:   mod  = (Model) modctx;
1317:   (*mod->solution)(mod, time, x, u, mod->solutionctx);
1318:   return(0);
1319: }

1321: PetscErrorCode SetInitialCondition(DM dm, Vec X, User user)
1322: {
1323:   PetscErrorCode     (*func[1]) (PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
1324:   void               *ctx[1];
1325:   Model              mod = user->model;
1326:   PetscErrorCode     ierr;

1329:   func[0] = SolutionFunctional;
1330:   ctx[0]  = (void *) mod;
1331:   DMProjectFunction(dm,0.0,func,ctx,INSERT_ALL_VALUES,X);
1332:   return(0);
1333: }

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

1340:   PetscViewerCreate(PetscObjectComm((PetscObject)dm), viewer);
1341:   PetscViewerSetType(*viewer, PETSCVIEWERVTK);
1342:   PetscViewerFileSetName(*viewer, filename);
1343:   return(0);
1344: }

1346: static PetscErrorCode MonitorVTK(TS ts,PetscInt stepnum,PetscReal time,Vec X,void *ctx)
1347: {
1348:   User           user = (User)ctx;
1349:   DM             dm;
1350:   Vec            cellgeom;
1351:   PetscViewer    viewer;
1352:   char           filename[PETSC_MAX_PATH_LEN],*ftable = NULL;
1353:   PetscReal      xnorm;
1354:   PetscInt       cEndInterior;

1358:   PetscObjectSetName((PetscObject) X, "u");
1359:   VecGetDM(X,&dm);
1360:   DMPlexTSGetGeometryFVM(dm, NULL, &cellgeom, NULL);
1361:   VecNorm(X,NORM_INFINITY,&xnorm);

1363:   if (stepnum >= 0) {
1364:     stepnum += user->monitorStepOffset;
1365:   }
1366:   if (stepnum >= 0) {           /* No summary for final time */
1367:     Model             mod = user->model;
1368:     PetscInt          c,cStart,cEnd,fcount,i;
1369:     size_t            ftableused,ftablealloc;
1370:     const PetscScalar *cgeom,*x;
1371:     DM                dmCell;
1372:     DMLabel           vtkLabel;
1373:     PetscReal         *fmin,*fmax,*fintegral,*ftmp;
1374:     fcount = mod->maxComputed+1;
1375:     PetscMalloc4(fcount,&fmin,fcount,&fmax,fcount,&fintegral,fcount,&ftmp);
1376:     for (i=0; i<fcount; i++) {
1377:       fmin[i]      = PETSC_MAX_REAL;
1378:       fmax[i]      = PETSC_MIN_REAL;
1379:       fintegral[i] = 0;
1380:     }
1381:     VecGetDM(cellgeom,&dmCell);
1382:     DMPlexGetHybridBounds(dmCell, &cEndInterior, NULL, NULL, NULL);
1383:     DMPlexGetHeightStratum(dmCell,0,&cStart,&cEnd);
1384:     VecGetArrayRead(cellgeom,&cgeom);
1385:     VecGetArrayRead(X,&x);
1386:     DMGetLabel(dm,"vtk",&vtkLabel);
1387:     for (c = cStart; c < cEndInterior; ++c) {
1388:       PetscFVCellGeom       *cg;
1389:       const PetscScalar     *cx    = NULL;
1390:       PetscInt              vtkVal = 0;

1392:       /* not that these two routines as currently implemented work for any dm with a
1393:        * defaultSection/defaultGlobalSection */
1394:       DMPlexPointLocalRead(dmCell,c,cgeom,&cg);
1395:       DMPlexPointGlobalRead(dm,c,x,&cx);
1396:       if (vtkLabel) {DMLabelGetValue(vtkLabel,c,&vtkVal);}
1397:       if (!vtkVal || !cx) continue;        /* ghost, or not a global cell */
1398:       for (i=0; i<mod->numCall; i++) {
1399:         FunctionalLink flink = mod->functionalCall[i];
1400:         (*flink->func)(mod,time,cg->centroid,cx,ftmp,flink->ctx);
1401:       }
1402:       for (i=0; i<fcount; i++) {
1403:         fmin[i]       = PetscMin(fmin[i],ftmp[i]);
1404:         fmax[i]       = PetscMax(fmax[i],ftmp[i]);
1405:         fintegral[i] += cg->volume * ftmp[i];
1406:       }
1407:     }
1408:     VecRestoreArrayRead(cellgeom,&cgeom);
1409:     VecRestoreArrayRead(X,&x);
1410:     MPI_Allreduce(MPI_IN_PLACE,fmin,fcount,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));
1411:     MPI_Allreduce(MPI_IN_PLACE,fmax,fcount,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));
1412:     MPI_Allreduce(MPI_IN_PLACE,fintegral,fcount,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)ts));

1414:     ftablealloc = fcount * 100;
1415:     ftableused  = 0;
1416:     PetscMalloc1(ftablealloc,&ftable);
1417:     for (i=0; i<mod->numMonitored; i++) {
1418:       size_t         countused;
1419:       char           buffer[256],*p;
1420:       FunctionalLink flink = mod->functionalMonitored[i];
1421:       PetscInt       id    = flink->offset;
1422:       if (i % 3) {
1423:         PetscMemcpy(buffer,"  ",2);
1424:         p    = buffer + 2;
1425:       } else if (i) {
1426:         char newline[] = "\n";
1427:         PetscMemcpy(buffer,newline,sizeof newline-1);
1428:         p    = buffer + sizeof newline - 1;
1429:       } else {
1430:         p = buffer;
1431:       }
1432:       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]);
1433:       countused += p - buffer;
1434:       if (countused > ftablealloc-ftableused-1) { /* reallocate */
1435:         char *ftablenew;
1436:         ftablealloc = 2*ftablealloc + countused;
1437:         PetscMalloc(ftablealloc,&ftablenew);
1438:         PetscMemcpy(ftablenew,ftable,ftableused);
1439:         PetscFree(ftable);
1440:         ftable = ftablenew;
1441:       }
1442:       PetscMemcpy(ftable+ftableused,buffer,countused);
1443:       ftableused += countused;
1444:       ftable[ftableused] = 0;
1445:     }
1446:     PetscFree4(fmin,fmax,fintegral,ftmp);

1448:     PetscPrintf(PetscObjectComm((PetscObject)ts),"% 3D  time %8.4g  |x| %8.4g  %s\n",stepnum,(double)time,(double)xnorm,ftable ? ftable : "");
1449:     PetscFree(ftable);
1450:   }
1451:   if (user->vtkInterval < 1) return(0);
1452:   if ((stepnum == -1) ^ (stepnum % user->vtkInterval == 0)) {
1453:     if (stepnum == -1) {        /* Final time is not multiple of normal time interval, write it anyway */
1454:       TSGetStepNumber(ts,&stepnum);
1455:     }
1456:     PetscSNPrintf(filename,sizeof filename,"%s-%03D.vtu",user->outputBasename,stepnum);
1457:     OutputVTK(dm,filename,&viewer);
1458:     VecView(X,viewer);
1459:     PetscViewerDestroy(&viewer);
1460:   }
1461:   return(0);
1462: }

1464: static PetscErrorCode initializeTS(DM dm, User user, TS *ts)
1465: {

1469:   TSCreate(PetscObjectComm((PetscObject)dm), ts);
1470:   TSSetType(*ts, TSSSP);
1471:   TSSetDM(*ts, dm);
1472:   if (user->vtkmon) {
1473:     TSMonitorSet(*ts,MonitorVTK,user,NULL);
1474:   }
1475:   DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, user);
1476:   DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, user);
1477:   TSSetMaxTime(*ts,2.0);
1478:   TSSetExactFinalTime(*ts,TS_EXACTFINALTIME_STEPOVER);
1479:   return(0);
1480: }

1482: static PetscErrorCode adaptToleranceFVM(PetscFV fvm, TS ts, Vec sol, VecTagger refineTag, VecTagger coarsenTag, User user, TS *tsNew, Vec *solNew)
1483: {
1484:   DM                dm, gradDM, plex, cellDM, adaptedDM = NULL;
1485:   Vec               cellGeom, faceGeom;
1486:   PetscBool         isForest, computeGradient;
1487:   Vec               grad, locGrad, locX, errVec;
1488:   PetscInt          cStart, cEnd, cEndInterior, c, dim, nRefine, nCoarsen;
1489:   PetscReal         minMaxInd[2] = {PETSC_MAX_REAL, PETSC_MIN_REAL}, minMaxIndGlobal[2], minInd, maxInd, time;
1490:   PetscScalar       *errArray;
1491:   const PetscScalar *pointVals;
1492:   const PetscScalar *pointGrads;
1493:   const PetscScalar *pointGeom;
1494:   DMLabel           adaptLabel = NULL;
1495:   IS                refineIS, coarsenIS;
1496:   PetscErrorCode    ierr;

1499:   TSGetTime(ts,&time);
1500:   VecGetDM(sol, &dm);
1501:   DMGetDimension(dm,&dim);
1502:   PetscFVGetComputeGradients(fvm,&computeGradient);
1503:   PetscFVSetComputeGradients(fvm,PETSC_TRUE);
1504:   DMIsForest(dm, &isForest);
1505:   DMConvert(dm, DMPLEX, &plex);
1506:   DMPlexGetDataFVM(plex, fvm, &cellGeom, &faceGeom, &gradDM);
1507:   DMCreateLocalVector(plex,&locX);
1508:   DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, 0.0, faceGeom, cellGeom, NULL);
1509:   DMGlobalToLocalBegin(plex, sol, INSERT_VALUES, locX);
1510:   DMGlobalToLocalEnd  (plex, sol, INSERT_VALUES, locX);
1511:   DMCreateGlobalVector(gradDM, &grad);
1512:   DMPlexReconstructGradientsFVM(plex, locX, grad);
1513:   DMCreateLocalVector(gradDM, &locGrad);
1514:   DMGlobalToLocalBegin(gradDM, grad, INSERT_VALUES, locGrad);
1515:   DMGlobalToLocalEnd(gradDM, grad, INSERT_VALUES, locGrad);
1516:   VecDestroy(&grad);
1517:   DMPlexGetHeightStratum(plex,0,&cStart,&cEnd);
1518:   DMPlexGetHybridBounds(plex,&cEndInterior,NULL,NULL,NULL);
1519:   cEnd = (cEndInterior < 0) ? cEnd : cEndInterior;
1520:   VecGetArrayRead(locGrad,&pointGrads);
1521:   VecGetArrayRead(cellGeom,&pointGeom);
1522:   VecGetArrayRead(locX,&pointVals);
1523:   VecGetDM(cellGeom,&cellDM);
1524:   DMLabelCreate("adapt",&adaptLabel);
1525:   VecCreateMPI(PetscObjectComm((PetscObject)plex),cEnd-cStart,PETSC_DETERMINE,&errVec);
1526:   VecSetUp(errVec);
1527:   VecGetArray(errVec,&errArray);
1528:   for (c = cStart; c < cEnd; c++) {
1529:     PetscReal             errInd = 0.;
1530:     PetscScalar           *pointGrad;
1531:     PetscScalar           *pointVal;
1532:     PetscFVCellGeom       *cg;

1534:     DMPlexPointLocalRead(gradDM,c,pointGrads,&pointGrad);
1535:     DMPlexPointLocalRead(cellDM,c,pointGeom,&cg);
1536:     DMPlexPointLocalRead(plex,c,pointVals,&pointVal);

1538:     (user->model->errorIndicator)(dim,cg->volume,user->model->physics->dof,pointVal,pointGrad,&errInd,user->model->errorCtx);
1539:     errArray[c-cStart] = errInd;
1540:     minMaxInd[0] = PetscMin(minMaxInd[0],errInd);
1541:     minMaxInd[1] = PetscMax(minMaxInd[1],errInd);
1542:   }
1543:   VecRestoreArray(errVec,&errArray);
1544:   VecRestoreArrayRead(locX,&pointVals);
1545:   VecRestoreArrayRead(cellGeom,&pointGeom);
1546:   VecRestoreArrayRead(locGrad,&pointGrads);
1547:   VecDestroy(&locGrad);
1548:   VecDestroy(&locX);
1549:   DMDestroy(&plex);

1551:   VecTaggerComputeIS(refineTag,errVec,&refineIS);
1552:   VecTaggerComputeIS(coarsenTag,errVec,&coarsenIS);
1553:   ISGetSize(refineIS,&nRefine);
1554:   ISGetSize(coarsenIS,&nCoarsen);
1555:   if (nRefine) {DMLabelSetStratumIS(adaptLabel,DM_ADAPT_REFINE,refineIS);}
1556:   if (nCoarsen) {DMLabelSetStratumIS(adaptLabel,DM_ADAPT_COARSEN,coarsenIS);}
1557:   ISDestroy(&coarsenIS);
1558:   ISDestroy(&refineIS);
1559:   VecDestroy(&errVec);

1561:   PetscFVSetComputeGradients(fvm,computeGradient);
1562:   minMaxInd[1] = -minMaxInd[1];
1563:   MPI_Allreduce(minMaxInd,minMaxIndGlobal,2,MPIU_REAL,MPI_MIN,PetscObjectComm((PetscObject)dm));
1564:   minInd = minMaxIndGlobal[0];
1565:   maxInd = -minMaxIndGlobal[1];
1566:   PetscInfo2(ts, "error indicator range (%E, %E)\n", minInd, maxInd);
1567:   if (nRefine || nCoarsen) { /* at least one cell is over the refinement threshold */
1568:     DMAdaptLabel(dm,adaptLabel,&adaptedDM);
1569:   }
1570:   DMLabelDestroy(&adaptLabel);
1571:   if (adaptedDM) {
1572:     PetscInfo2(ts, "Adapted mesh, marking %D cells for refinement, and %D cells for coarsening\n", nRefine, nCoarsen);
1573:     if (tsNew) {initializeTS(adaptedDM, user, tsNew);}
1574:     if (solNew) {
1575:       DMCreateGlobalVector(adaptedDM, solNew);
1576:       PetscObjectSetName((PetscObject) *solNew, "solution");
1577:       DMForestTransferVec(dm, sol, adaptedDM, *solNew, PETSC_TRUE, time);
1578:     }
1579:     if (isForest) {DMForestSetAdaptivityForest(adaptedDM,NULL);} /* clear internal references to the previous dm */
1580:     DMDestroy(&adaptedDM);
1581:   } else {
1582:     if (tsNew)  *tsNew  = NULL;
1583:     if (solNew) *solNew = NULL;
1584:   }
1585:   return(0);
1586: }

1588: int main(int argc, char **argv)
1589: {
1590:   MPI_Comm          comm;
1591:   PetscDS           prob;
1592:   PetscFV           fvm;
1593:   PetscLimiter      limiter = NULL, noneLimiter = NULL;
1594:   User              user;
1595:   Model             mod;
1596:   Physics           phys;
1597:   DM                dm;
1598:   PetscReal         ftime, cfl, dt, minRadius;
1599:   PetscInt          dim, nsteps;
1600:   TS                ts;
1601:   TSConvergedReason reason;
1602:   Vec               X;
1603:   PetscViewer       viewer;
1604:   PetscBool         simplex = PETSC_FALSE, vtkCellGeom, splitFaces, useAMR;
1605:   PetscInt          overlap, adaptInterval;
1606:   char              filename[PETSC_MAX_PATH_LEN] = "";
1607:   char              physname[256]  = "advect";
1608:   VecTagger         refineTag = NULL, coarsenTag = NULL;
1609:   PetscErrorCode    ierr;

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

1614:   PetscNew(&user);
1615:   PetscNew(&user->model);
1616:   PetscNew(&user->model->physics);
1617:   mod           = user->model;
1618:   phys          = mod->physics;
1619:   mod->comm     = comm;
1620:   useAMR        = PETSC_FALSE;
1621:   adaptInterval = 1;

1623:   /* Register physical models to be available on the command line */
1624:   PetscFunctionListAdd(&PhysicsList,"advect"          ,PhysicsCreate_Advect);
1625:   PetscFunctionListAdd(&PhysicsList,"sw"              ,PhysicsCreate_SW);
1626:   PetscFunctionListAdd(&PhysicsList,"euler"           ,PhysicsCreate_Euler);

1628:   PetscOptionsBegin(comm,NULL,"Unstructured Finite Volume Mesh Options","");
1629:   {
1630:     cfl  = 0.9 * 4; /* default SSPRKS2 with s=5 stages is stable for CFL number s-1 */
1631:     PetscOptionsReal("-ufv_cfl","CFL number per step","",cfl,&cfl,NULL);
1632:     PetscOptionsString("-f","Exodus.II filename to read","",filename,filename,sizeof(filename),NULL);
1633:     PetscOptionsBool("-simplex","Flag to use a simplex mesh","",simplex,&simplex,NULL);
1634:     splitFaces = PETSC_FALSE;
1635:     PetscOptionsBool("-ufv_split_faces","Split faces between cell sets","",splitFaces,&splitFaces,NULL);
1636:     overlap = 1;
1637:     PetscOptionsInt("-ufv_mesh_overlap","Number of cells to overlap partitions","",overlap,&overlap,NULL);
1638:     user->vtkInterval = 1;
1639:     PetscOptionsInt("-ufv_vtk_interval","VTK output interval (0 to disable)","",user->vtkInterval,&user->vtkInterval,NULL);
1640:     user->vtkmon = PETSC_TRUE;
1641:     PetscOptionsBool("-ufv_vtk_monitor","Use VTKMonitor routine","",user->vtkmon,&user->vtkmon,NULL);
1642:     vtkCellGeom = PETSC_FALSE;
1643:     PetscStrcpy(user->outputBasename, "ex11");
1644:     PetscOptionsString("-ufv_vtk_basename","VTK output basename","",user->outputBasename,user->outputBasename,PETSC_MAX_PATH_LEN,NULL);
1645:     PetscOptionsBool("-ufv_vtk_cellgeom","Write cell geometry (for debugging)","",vtkCellGeom,&vtkCellGeom,NULL);
1646:     PetscOptionsBool("-ufv_use_amr","use local adaptive mesh refinement","",useAMR,&useAMR,NULL);
1647:     PetscOptionsInt("-ufv_adapt_interval","time steps between AMR","",adaptInterval,&adaptInterval,NULL);
1648:   }
1649:   PetscOptionsEnd();

1651:   if (useAMR) {
1652:     VecTaggerBox refineBox, coarsenBox;

1654:     refineBox.min  = refineBox.max  = PETSC_MAX_REAL;
1655:     coarsenBox.min = coarsenBox.max = PETSC_MIN_REAL;

1657:     VecTaggerCreate(comm,&refineTag);
1658:     PetscObjectSetOptionsPrefix((PetscObject)refineTag,"refine_");
1659:     VecTaggerSetType(refineTag,VECTAGGERABSOLUTE);
1660:     VecTaggerAbsoluteSetBox(refineTag,&refineBox);
1661:     VecTaggerSetFromOptions(refineTag);
1662:     VecTaggerSetUp(refineTag);
1663:     PetscObjectViewFromOptions((PetscObject)refineTag,NULL,"-tag_view");

1665:     VecTaggerCreate(comm,&coarsenTag);
1666:     PetscObjectSetOptionsPrefix((PetscObject)coarsenTag,"coarsen_");
1667:     VecTaggerSetType(coarsenTag,VECTAGGERABSOLUTE);
1668:     VecTaggerAbsoluteSetBox(coarsenTag,&coarsenBox);
1669:     VecTaggerSetFromOptions(coarsenTag);
1670:     VecTaggerSetUp(coarsenTag);
1671:     PetscObjectViewFromOptions((PetscObject)coarsenTag,NULL,"-tag_view");
1672:   }

1674:   PetscOptionsBegin(comm,NULL,"Unstructured Finite Volume Physics Options","");
1675:   {
1676:     PetscErrorCode (*physcreate)(Model,Physics,PetscOptionItems*);
1677:     PetscOptionsFList("-physics","Physics module to solve","",PhysicsList,physname,physname,sizeof physname,NULL);
1678:     PetscFunctionListFind(PhysicsList,physname,&physcreate);
1679:     PetscMemzero(phys,sizeof(struct _n_Physics));
1680:     (*physcreate)(mod,phys,PetscOptionsObject);
1681:     /* Count number of fields and dofs */
1682:     for (phys->nfields=0,phys->dof=0; phys->field_desc[phys->nfields].name; phys->nfields++) phys->dof += phys->field_desc[phys->nfields].dof;
1683:     if (phys->dof <= 0) SETERRQ1(comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set dof",physname);
1684:     ModelFunctionalSetFromOptions(mod,PetscOptionsObject);
1685:   }
1686:   PetscOptionsEnd();

1688:   /* Create mesh */
1689:   {
1690:     size_t len,i;
1691:     for (i = 0; i < DIM; i++) { mod->bounds[2*i] = 0.; mod->bounds[2*i+1] = 1.;};
1692:     PetscStrlen(filename,&len);
1693:     dim = DIM;
1694:     if (!len) { /* a null name means just do a hex box */
1695:       PetscInt cells[3] = {1, 1, 1}; /* coarse mesh is one cell; refine from there */
1696:       PetscBool flg1, flg2, skew = PETSC_FALSE;
1697:       PetscInt nret1 = DIM;
1698:       PetscInt nret2 = 2*DIM;
1699:       PetscOptionsBegin(comm,NULL,"Rectangular mesh options","");
1700:       PetscOptionsIntArray("-grid_size","number of cells in each direction","",cells,&nret1,&flg1);
1701:       PetscOptionsRealArray("-grid_bounds","bounds of the mesh in each direction (i.e., x_min,x_max,y_min,y_max","",mod->bounds,&nret2,&flg2);
1702:       PetscOptionsBool("-grid_skew_60","Skew grid for 60 degree shock mesh","",skew,&skew,NULL);
1703:       PetscOptionsEnd();
1704:       if (flg1) {
1705:         dim = nret1;
1706:         if (dim != DIM) SETERRQ1(comm,PETSC_ERR_ARG_SIZ,"Dim wrong size %D in -grid_size",dim);
1707:       }
1708:       DMPlexCreateBoxMesh(comm, dim, simplex, cells, NULL, NULL, mod->bcs, PETSC_TRUE, &dm);
1709:       if (flg2) {
1710:         PetscInt dimEmbed, i;
1711:         PetscInt nCoords;
1712:         PetscScalar *coords;
1713:         Vec coordinates;

1715:         DMGetCoordinatesLocal(dm,&coordinates);
1716:         DMGetCoordinateDim(dm,&dimEmbed);
1717:         VecGetLocalSize(coordinates,&nCoords);
1718:         if (nCoords % dimEmbed) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Coordinate vector the wrong size");
1719:         VecGetArray(coordinates,&coords);
1720:         for (i = 0; i < nCoords; i += dimEmbed) {
1721:           PetscInt j;

1723:           PetscScalar *coord = &coords[i];
1724:           for (j = 0; j < dimEmbed; j++) {
1725:             coord[j] = mod->bounds[2 * j] + coord[j] * (mod->bounds[2 * j + 1] - mod->bounds[2 * j]);
1726:             if (dim==2 && cells[1]==1 && j==0 && skew) {
1727:               if (cells[0]==2 && i==8) {
1728:                 coord[j] = .57735026918963; /* hack to get 60 deg skewed mesh */
1729:               }
1730:               else if (cells[0]==3) {
1731:                 if(i==2 || i==10) coord[j] = mod->bounds[1]/4.;
1732:                 else if (i==4) coord[j] = mod->bounds[1]/2.;
1733:                 else if (i==12) coord[j] = 1.57735026918963*mod->bounds[1]/2.;
1734:               }
1735:             }
1736:           }
1737:         }
1738:         VecRestoreArray(coordinates,&coords);
1739:         DMSetCoordinatesLocal(dm,coordinates);
1740:       }
1741:     } else {
1742:       DMPlexCreateFromFile(comm, filename, PETSC_TRUE, &dm);
1743:     }
1744:   }
1745:   DMViewFromOptions(dm, NULL, "-orig_dm_view");
1746:   DMGetDimension(dm, &dim);

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

1751:   mod->errorIndicator = ErrorIndicator_Simple;

1753:   {
1754:     DM dmDist;

1756:     DMPlexSetAdjacencyUseCone(dm, PETSC_TRUE);
1757:     DMPlexSetAdjacencyUseClosure(dm, PETSC_FALSE);
1758:     DMPlexDistribute(dm, overlap, NULL, &dmDist);
1759:     if (dmDist) {
1760:       DMDestroy(&dm);
1761:       dm   = dmDist;
1762:     }
1763:   }

1765:   DMSetFromOptions(dm);

1767:   {
1768:     DM gdm;

1770:     DMPlexConstructGhostCells(dm, NULL, NULL, &gdm);
1771:     DMDestroy(&dm);
1772:     dm   = gdm;
1773:     DMViewFromOptions(dm, NULL, "-dm_view");
1774:   }
1775:   if (splitFaces) {ConstructCellBoundary(dm, user);}
1776:   SplitFaces(&dm, "split faces", user);

1778:   PetscDSCreate(PetscObjectComm((PetscObject)dm),&prob);
1779:   PetscFVCreate(comm, &fvm);
1780:   PetscFVSetFromOptions(fvm);
1781:   PetscFVSetNumComponents(fvm, phys->dof);
1782:   PetscFVSetSpatialDimension(fvm, dim);
1783:   PetscObjectSetName((PetscObject) fvm,"");
1784:   {
1785:     PetscInt f, dof;
1786:     for (f=0,dof=0; f < phys->nfields; f++) {
1787:       PetscInt newDof = phys->field_desc[f].dof;

1789:       if (newDof == 1) {
1790:         PetscFVSetComponentName(fvm,dof,phys->field_desc[f].name);
1791:       }
1792:       else {
1793:         PetscInt j;

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

1798:           PetscSNPrintf(compName,sizeof(compName),"%s_%d",phys->field_desc[f].name,j);
1799:           PetscFVSetComponentName(fvm,dof+j,compName);
1800:         }
1801:       }
1802:       dof += newDof;
1803:     }
1804:   }
1805:   /* FV is now structured with one field having all physics as components */
1806:   PetscDSAddDiscretization(prob, (PetscObject) fvm);
1807:   PetscDSSetRiemannSolver(prob, 0, user->model->physics->riemann);
1808:   PetscDSSetContext(prob, 0, user->model->physics);
1809:   (*mod->setupbc)(prob,phys);
1810:   PetscDSSetFromOptions(prob);
1811:   DMSetDS(dm,prob);
1812:   PetscDSDestroy(&prob);
1813:   {
1814:     char      convType[256];
1815:     PetscBool flg;

1817:     PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");
1818:     PetscOptionsFList("-dm_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg);
1819:     PetscOptionsEnd();
1820:     if (flg) {
1821:       DM dmConv;

1823:       DMConvert(dm,convType,&dmConv);
1824:       if (dmConv) {
1825:         DMViewFromOptions(dmConv, NULL, "-dm_conv_view");
1826:         DMDestroy(&dm);
1827:         dm   = dmConv;
1828:         DMSetFromOptions(dm);
1829:       }
1830:     }
1831:   }

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

1835:   DMCreateGlobalVector(dm, &X);
1836:   PetscObjectSetName((PetscObject) X, "solution");
1837:   SetInitialCondition(dm, X, user);
1838:   if (useAMR) {
1839:     PetscInt adaptIter;

1841:     /* use no limiting when reconstructing gradients for adaptivity */
1842:     PetscFVGetLimiter(fvm, &limiter);
1843:     PetscObjectReference((PetscObject) limiter);
1844:     PetscLimiterCreate(PetscObjectComm((PetscObject) fvm), &noneLimiter);
1845:     PetscLimiterSetType(noneLimiter, PETSCLIMITERNONE);

1847:     PetscFVSetLimiter(fvm, noneLimiter);
1848:     for (adaptIter = 0; ; ++adaptIter) {
1849:       PetscLogDouble bytes;
1850:       TS             tsNew = NULL;

1852:       PetscMemoryGetCurrentUsage(&bytes);
1853:       PetscInfo2(ts, "refinement loop %D: memory used %g\n", adaptIter, bytes);
1854:       DMViewFromOptions(dm, NULL, "-initial_dm_view");
1855:       VecViewFromOptions(X, NULL, "-initial_vec_view");
1856: #if 0
1857:       if (viewInitial) {
1858:         PetscViewer viewer;
1859:         char        buf[256];
1860:         PetscBool   isHDF5, isVTK;

1862:         PetscViewerCreate(comm,&viewer);
1863:         PetscViewerSetType(viewer,PETSCVIEWERVTK);
1864:         PetscViewerSetOptionsPrefix(viewer,"initial_");
1865:         PetscViewerSetFromOptions(viewer);
1866:         PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&isHDF5);
1867:         PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isVTK);
1868:         if (isHDF5) {
1869:           PetscSNPrintf(buf, 256, "ex11-initial-%d.h5", adaptIter);
1870:         } else if (isVTK) {
1871:           PetscSNPrintf(buf, 256, "ex11-initial-%d.vtu", adaptIter);
1872:           PetscViewerPushFormat(viewer,PETSC_VIEWER_VTK_VTU);
1873:         }
1874:         PetscViewerFileSetMode(viewer,FILE_MODE_WRITE);
1875:         PetscViewerFileSetName(viewer,buf);
1876:         if (isHDF5) {
1877:           DMView(dm,viewer);
1878:           PetscViewerFileSetMode(viewer,FILE_MODE_UPDATE);
1879:         }
1880:         VecView(X,viewer);
1881:         PetscViewerDestroy(&viewer);
1882:       }
1883: #endif

1885:       adaptToleranceFVM(fvm, ts, X, refineTag, coarsenTag, user, &tsNew, NULL);
1886:       if (!tsNew) {
1887:         break;
1888:       } else {
1889:         DMDestroy(&dm);
1890:         VecDestroy(&X);
1891:         TSDestroy(&ts);
1892:         ts   = tsNew;
1893:         TSGetDM(ts,&dm);
1894:         PetscObjectReference((PetscObject)dm);
1895:         DMCreateGlobalVector(dm,&X);
1896:         PetscObjectSetName((PetscObject) X, "solution");
1897:         SetInitialCondition(dm, X, user);
1898:       }
1899:     }
1900:     /* restore original limiter */
1901:     PetscFVSetLimiter(fvm, limiter);
1902:   }

1904:   if (vtkCellGeom) {
1905:     DM  dmCell;
1906:     Vec cellgeom, partition;

1908:     DMPlexTSGetGeometryFVM(dm, NULL, &cellgeom, NULL);
1909:     OutputVTK(dm, "ex11-cellgeom.vtk", &viewer);
1910:     VecView(cellgeom, viewer);
1911:     PetscViewerDestroy(&viewer);
1912:     CreatePartitionVec(dm, &dmCell, &partition);
1913:     OutputVTK(dmCell, "ex11-partition.vtk", &viewer);
1914:     VecView(partition, viewer);
1915:     PetscViewerDestroy(&viewer);
1916:     VecDestroy(&partition);
1917:     DMDestroy(&dmCell);
1918:   }

1920:   /* collect max maxspeed from all processes -- todo */
1921:   DMPlexTSGetGeometryFVM(dm, NULL, NULL, &minRadius);
1922:   MPI_Allreduce(&phys->maxspeed,&mod->maxspeed,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));
1923:   if (mod->maxspeed <= 0) SETERRQ1(comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set maxspeed",physname);
1924:   dt   = cfl * minRadius / mod->maxspeed;
1925:   TSSetTimeStep(ts,dt);
1926:   TSSetFromOptions(ts);
1927:   if (!useAMR) {
1928:     TSSolve(ts,X);
1929:     TSGetSolveTime(ts,&ftime);
1930:     TSGetStepNumber(ts,&nsteps);
1931:   } else {
1932:     PetscReal finalTime;
1933:     PetscInt  adaptIter;
1934:     TS        tsNew = NULL;
1935:     Vec       solNew = NULL;

1937:     TSGetMaxTime(ts,&finalTime);
1938:     TSSetMaxSteps(ts,adaptInterval);
1939:     TSSolve(ts,X);
1940:     TSGetSolveTime(ts,&ftime);
1941:     TSGetStepNumber(ts,&nsteps);
1942:     for (adaptIter = 0;ftime < finalTime;adaptIter++) {
1943:       PetscLogDouble bytes;

1945:       PetscMemoryGetCurrentUsage(&bytes);
1946:       PetscInfo2(ts, "AMR time step loop %D: memory used %g\n", adaptIter, bytes);
1947:       PetscFVSetLimiter(fvm,noneLimiter);
1948:       adaptToleranceFVM(fvm,ts,X,refineTag,coarsenTag,user,&tsNew,&solNew);
1949:       PetscFVSetLimiter(fvm,limiter);
1950:       if (tsNew) {
1951:         PetscInfo(ts, "AMR used\n");
1952:         DMDestroy(&dm);
1953:         VecDestroy(&X);
1954:         TSDestroy(&ts);
1955:         ts   = tsNew;
1956:         X    = solNew;
1957:         TSSetFromOptions(ts);
1958:         VecGetDM(X,&dm);
1959:         PetscObjectReference((PetscObject)dm);
1960:         DMPlexTSGetGeometryFVM(dm, NULL, NULL, &minRadius);
1961:         MPI_Allreduce(&phys->maxspeed,&mod->maxspeed,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));
1962:         if (mod->maxspeed <= 0) SETERRQ1(comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set maxspeed",physname);
1963:         dt   = cfl * minRadius / mod->maxspeed;
1964:         TSSetStepNumber(ts,nsteps);
1965:         TSSetTime(ts,ftime);
1966:         TSSetTimeStep(ts,dt);
1967:       } else {
1968:         PetscInfo(ts, "AMR not used\n");
1969:       }
1970:       user->monitorStepOffset = nsteps;
1971:       TSSetMaxSteps(ts,nsteps+adaptInterval);
1972:       TSSolve(ts,X);
1973:       TSGetSolveTime(ts,&ftime);
1974:       TSGetStepNumber(ts,&nsteps);
1975:     }
1976:   }
1977:   TSGetConvergedReason(ts,&reason);
1978:   PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,nsteps);
1979:   TSDestroy(&ts);

1981:   VecTaggerDestroy(&refineTag);
1982:   VecTaggerDestroy(&coarsenTag);
1983:   PetscFunctionListDestroy(&PhysicsList);
1984:   FunctionalLinkDestroy(&user->model->functionalRegistry);
1985:   PetscFree(user->model->functionalMonitored);
1986:   PetscFree(user->model->functionalCall);
1987:   PetscFree(user->model->physics->data);
1988:   PetscFree(user->model->physics);
1989:   PetscFree(user->model);
1990:   PetscFree(user);
1991:   VecDestroy(&X);
1992:   PetscLimiterDestroy(&limiter);
1993:   PetscLimiterDestroy(&noneLimiter);
1994:   PetscFVDestroy(&fvm);
1995:   DMDestroy(&dm);
1996:   PetscFinalize();
1997:   return ierr;
1998: }

2000: /* Godunov fluxs */
2001: PetscScalar cvmgp_(PetscScalar *a, PetscScalar *b, PetscScalar *test)
2002: {
2003:     /* System generated locals */
2004:     PetscScalar ret_val;

2006:     if (PetscRealPart(*test) > 0.) {
2007:         goto L10;
2008:     }
2009:     ret_val = *b;
2010:     return ret_val;
2011: L10:
2012:     ret_val = *a;
2013:     return ret_val;
2014: } /* cvmgp_ */

2016: PetscScalar cvmgm_(PetscScalar *a, PetscScalar *b, PetscScalar *test)
2017: {
2018:     /* System generated locals */
2019:     PetscScalar ret_val;

2021:     if (PetscRealPart(*test) < 0.) {
2022:         goto L10;
2023:     }
2024:     ret_val = *b;
2025:     return ret_val;
2026: L10:
2027:     ret_val = *a;
2028:     return ret_val;
2029: } /* cvmgm_ */

2031: int riem1mdt( PetscScalar *gaml, PetscScalar *gamr, PetscScalar *rl, PetscScalar *pl,
2032:               PetscScalar *uxl, PetscScalar *rr, PetscScalar *pr,
2033:               PetscScalar *uxr, PetscScalar *rstarl, PetscScalar *rstarr, PetscScalar *
2034:               pstar, PetscScalar *ustar)
2035: {
2036:     /* Initialized data */

2038:     static PetscScalar smallp = 1e-8;

2040:     /* System generated locals */
2041:     int i__1;
2042:     PetscScalar d__1, d__2;

2044:     /* Local variables */
2045:     static int i0;
2046:     static PetscScalar cl, cr, wl, zl, wr, zr, pst, durl, skpr1, skpr2;
2047:     static int iwave;
2048:     static PetscScalar gascl4, gascr4, cstarl, dpstar, cstarr;
2049:     /* static PetscScalar csqrl, csqrr, gascl1, gascl2, gascl3, gascr1, gascr2, gascr3; */
2050:     static int iterno;
2051:     static PetscScalar ustarl, ustarr, rarepr1, rarepr2;



2055:     /* gascl1 = *gaml - 1.; */
2056:     /* gascl2 = (*gaml + 1.) * .5; */
2057:     /* gascl3 = gascl2 / *gaml; */
2058:     gascl4 = 1. / (*gaml - 1.);

2060:     /* gascr1 = *gamr - 1.; */
2061:     /* gascr2 = (*gamr + 1.) * .5; */
2062:     /* gascr3 = gascr2 / *gamr; */
2063:     gascr4 = 1. / (*gamr - 1.);
2064:     iterno = 10;
2065: /*        find pstar: */
2066:     cl = PetscSqrtScalar(*gaml * *pl / *rl);
2067:     cr = PetscSqrtScalar(*gamr * *pr / *rr);
2068:     wl = *rl * cl;
2069:     wr = *rr * cr;
2070:     /* csqrl = wl * wl; */
2071:     /* csqrr = wr * wr; */
2072:     *pstar = (wl * *pr + wr * *pl) / (wl + wr);
2073:     *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp));
2074:     pst = *pl / *pr;
2075:     skpr1 = cr * (pst - 1.) * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
2076:     d__1 = (*gamr - 1.) / (*gamr * 2.);
2077:     rarepr2 = gascr4 * 2. * cr * (1. - PetscPowScalar(pst, d__1));
2078:     pst = *pr / *pl;
2079:     skpr2 = cl * (pst - 1.) * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
2080:     d__1 = (*gaml - 1.) / (*gaml * 2.);
2081:     rarepr1 = gascl4 * 2. * cl * (1. - PetscPowScalar(pst, d__1));
2082:     durl = *uxr - *uxl;
2083:     if (PetscRealPart(*pr) < PetscRealPart(*pl)) {
2084:         if (PetscRealPart(durl) >= PetscRealPart(rarepr1)) {
2085:             iwave = 100;
2086:         } else if (PetscRealPart(durl) <= PetscRealPart(-skpr1)) {
2087:             iwave = 300;
2088:         } else {
2089:             iwave = 400;
2090:         }
2091:     } else {
2092:         if (PetscRealPart(durl) >= PetscRealPart(rarepr2)) {
2093:             iwave = 100;
2094:         } else if (PetscRealPart(durl) <= PetscRealPart(-skpr2)) {
2095:             iwave = 300;
2096:         } else {
2097:             iwave = 200;
2098:         }
2099:     }
2100:     if (iwave == 100) {
2101: /*     1-wave: rarefaction wave, 3-wave: rarefaction wave */
2102: /*     case (100) */
2103:         i__1 = iterno;
2104:         for (i0 = 1; i0 <= i__1; ++i0) {
2105:             d__1 = *pstar / *pl;
2106:             d__2 = 1. / *gaml;
2107:             *rstarl = *rl * PetscPowScalar(d__1, d__2);
2108:             cstarl = PetscSqrtScalar(*gaml * *pstar / *rstarl);
2109:             ustarl = *uxl - gascl4 * 2. * (cstarl - cl);
2110:             zl = *rstarl * cstarl;
2111:             d__1 = *pstar / *pr;
2112:             d__2 = 1. / *gamr;
2113:             *rstarr = *rr * PetscPowScalar(d__1, d__2);
2114:             cstarr = PetscSqrtScalar(*gamr * *pstar / *rstarr);
2115:             ustarr = *uxr + gascr4 * 2. * (cstarr - cr);
2116:             zr = *rstarr * cstarr;
2117:             dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
2118:             *pstar -= dpstar;
2119:             *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp));
2120:             if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
2121: #if 0
2122:         break;
2123: #endif
2124:             }
2125:         }
2126: /*     1-wave: shock wave, 3-wave: rarefaction wave */
2127:     } else if (iwave == 200) {
2128: /*     case (200) */
2129:         i__1 = iterno;
2130:         for (i0 = 1; i0 <= i__1; ++i0) {
2131:             pst = *pstar / *pl;
2132:             ustarl = *uxl - (pst - 1.) * cl * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
2133:             zl = *pl / cl * PetscSqrtScalar(*gaml * 2. * (*gaml - 1. + (*gaml + 1.) * pst)) * (*gaml - 1. + (*gaml + 1.) * pst) / (*gaml * 3. - 1. + (*gaml + 1.) * pst);
2134:             d__1 = *pstar / *pr;
2135:             d__2 = 1. / *gamr;
2136:             *rstarr = *rr * PetscPowScalar(d__1, d__2);
2137:             cstarr = PetscSqrtScalar(*gamr * *pstar / *rstarr);
2138:             zr = *rstarr * cstarr;
2139:             ustarr = *uxr + gascr4 * 2. * (cstarr - cr);
2140:             dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
2141:             *pstar -= dpstar;
2142:             *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp));
2143:             if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
2144: #if 0
2145:         break;
2146: #endif
2147:             }
2148:         }
2149: /*     1-wave: shock wave, 3-wave: shock */
2150:     } else if (iwave == 300) {
2151: /*     case (300) */
2152:         i__1 = iterno;
2153:         for (i0 = 1; i0 <= i__1; ++i0) {
2154:             pst = *pstar / *pl;
2155:             ustarl = *uxl - (pst - 1.) * cl * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
2156:             zl = *pl / cl * PetscSqrtScalar(*gaml * 2. * (*gaml - 1. + (*gaml + 1.) * pst)) * (*gaml - 1. + (*gaml + 1.) * pst) / (*gaml * 3. - 1. + (*gaml + 1.) * pst);
2157:             pst = *pstar / *pr;
2158:             ustarr = *uxr + (pst - 1.) * cr * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
2159:             zr = *pr / cr * PetscSqrtScalar(*gamr * 2. * (*gamr - 1. + (*gamr + 1.) * pst)) * (*gamr - 1. + (*gamr + 1.) * pst) / (*gamr * 3. - 1. + (*gamr + 1.) * pst);
2160:             dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
2161:             *pstar -= dpstar;
2162:             *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp));
2163:             if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
2164: #if 0
2165:         break;
2166: #endif
2167:             }
2168:         }
2169: /*     1-wave: rarefaction wave, 3-wave: shock */
2170:     } else if (iwave == 400) {
2171: /*     case (400) */
2172:         i__1 = iterno;
2173:         for (i0 = 1; i0 <= i__1; ++i0) {
2174:             d__1 = *pstar / *pl;
2175:             d__2 = 1. / *gaml;
2176:             *rstarl = *rl * PetscPowScalar(d__1, d__2);
2177:             cstarl = PetscSqrtScalar(*gaml * *pstar / *rstarl);
2178:             ustarl = *uxl - gascl4 * 2. * (cstarl - cl);
2179:             zl = *rstarl * cstarl;
2180:             pst = *pstar / *pr;
2181:             ustarr = *uxr + (pst - 1.) * cr * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
2182:             zr = *pr / cr * PetscSqrtScalar(*gamr * 2. * (*gamr - 1. + (*gamr + 1.) * pst)) * (*gamr - 1. + (*gamr + 1.) * pst) / (*gamr * 3. - 1. + (*gamr + 1.) * pst);
2183:             dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
2184:             *pstar -= dpstar;
2185:             *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp));
2186:             if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
2187: #if 0
2188:               break;
2189: #endif
2190:             }
2191:         }
2192:     }

2194:     *ustar = (zl * ustarr + zr * ustarl) / (zl + zr);
2195:     if (PetscRealPart(*pstar) > PetscRealPart(*pl)) {
2196:         pst = *pstar / *pl;
2197:         *rstarl = ((*gaml + 1.) * pst + *gaml - 1.) / ((*gaml - 1.) * pst + *
2198:                 gaml + 1.) * *rl;
2199:     }
2200:     if (PetscRealPart(*pstar) > PetscRealPart(*pr)) {
2201:         pst = *pstar / *pr;
2202:         *rstarr = ((*gamr + 1.) * pst + *gamr - 1.) / ((*gamr - 1.) * pst + *
2203:                 gamr + 1.) * *rr;
2204:     }
2205:     return iwave;
2206: }

2208: PetscScalar sign(PetscScalar x)
2209: {
2210:     if(PetscRealPart(x) > 0) return 1.0;
2211:     if(PetscRealPart(x) < 0) return -1.0;
2212:     return 0.0;
2213: }
2214: /*        Riemann Solver */
2215: /* -------------------------------------------------------------------- */
2216: int riemannsolver(PetscScalar *xcen, PetscScalar *xp,
2217:                    PetscScalar *dtt, PetscScalar *rl, PetscScalar *uxl, PetscScalar *pl,
2218:                    PetscScalar *utl, PetscScalar *ubl, PetscScalar *gaml, PetscScalar *rho1l,
2219:                    PetscScalar *rr, PetscScalar *uxr, PetscScalar *pr, PetscScalar *utr,
2220:                    PetscScalar *ubr, PetscScalar *gamr, PetscScalar *rho1r, PetscScalar *rx,
2221:                    PetscScalar *uxm, PetscScalar *px, PetscScalar *utx, PetscScalar *ubx,
2222:                    PetscScalar *gam, PetscScalar *rho1)
2223: {
2224:     /* System generated locals */
2225:     PetscScalar d__1, d__2;

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

2232:     if (*rl == *rr && *pr == *pl && *uxl == *uxr && *gaml == *gamr) {
2233:         *rx = *rl;
2234:         *px = *pl;
2235:         *uxm = *uxl;
2236:         *gam = *gaml;
2237:         x2 = *xcen + *uxm * *dtt;

2239:         if (PetscRealPart(*xp) >= PetscRealPart(x2)) {
2240:             *utx = *utr;
2241:             *ubx = *ubr;
2242:             *rho1 = *rho1r;
2243:         } else {
2244:             *utx = *utl;
2245:             *ubx = *ubl;
2246:             *rho1 = *rho1l;
2247:         }
2248:         return 0;
2249:     }
2250:     iwave = riem1mdt(gaml, gamr, rl, pl, uxl, rr, pr, uxr, &rstarl, &rstarr, &pstar, &ustar);

2252:     x2 = *xcen + ustar * *dtt;
2253:     d__1 = *xp - x2;
2254:     sgn0 = sign(d__1);
2255: /*            x is in 3-wave if sgn0 = 1 */
2256: /*            x is in 1-wave if sgn0 = -1 */
2257:     r0 = cvmgm_(rl, rr, &sgn0);
2258:     p0 = cvmgm_(pl, pr, &sgn0);
2259:     u0 = cvmgm_(uxl, uxr, &sgn0);
2260:     *gam = cvmgm_(gaml, gamr, &sgn0);
2261:     gasc1 = *gam - 1.;
2262:     gasc2 = (*gam + 1.) * .5;
2263:     gasc3 = gasc2 / *gam;
2264:     gasc4 = 1. / (*gam - 1.);
2265:     c0 = PetscSqrtScalar(*gam * p0 / r0);
2266:     streng = pstar - p0;
2267:     w0 = *gam * r0 * p0 * (gasc3 * streng / p0 + 1.);
2268:     rstars = r0 / (1. - r0 * streng / w0);
2269:     d__1 = p0 / pstar;
2270:     d__2 = -1. / *gam;
2271:     rstarr = r0 * PetscPowScalar(d__1, d__2);
2272:     rstar = cvmgm_(&rstarr, &rstars, &streng);
2273:     w0 = PetscSqrtScalar(w0);
2274:     cstar = PetscSqrtScalar(*gam * pstar / rstar);
2275:     wsp0 = u0 + sgn0 * c0;
2276:     wspst = ustar + sgn0 * cstar;
2277:     ushock = ustar + sgn0 * w0 / rstar;
2278:     wspst = cvmgp_(&ushock, &wspst, &streng);
2279:     wsp0 = cvmgp_(&ushock, &wsp0, &streng);
2280:     x0 = *xcen + wsp0 * *dtt;
2281:     xstar = *xcen + wspst * *dtt;
2282: /*           using gas formula to evaluate rarefaction wave */
2283: /*            ri : reiman invariant */
2284:     ri = u0 - sgn0 * 2. * gasc4 * c0;
2285:     cx = sgn0 * .5 * gasc1 / gasc2 * ((*xp - *xcen) / *dtt - ri);
2286:     *uxm = ri + sgn0 * 2. * gasc4 * cx;
2287:     s = p0 / PetscPowScalar(r0, *gam);
2288:     d__1 = cx * cx / (*gam * s);
2289:     *rx = PetscPowScalar(d__1, gasc4);
2290:     *px = cx * cx * *rx / *gam;
2291:     d__1 = sgn0 * (x0 - *xp);
2292:     *rx = cvmgp_(rx, &r0, &d__1);
2293:     d__1 = sgn0 * (x0 - *xp);
2294:     *px = cvmgp_(px, &p0, &d__1);
2295:     d__1 = sgn0 * (x0 - *xp);
2296:     *uxm = cvmgp_(uxm, &u0, &d__1);
2297:     d__1 = sgn0 * (xstar - *xp);
2298:     *rx = cvmgm_(rx, &rstar, &d__1);
2299:     d__1 = sgn0 * (xstar - *xp);
2300:     *px = cvmgm_(px, &pstar, &d__1);
2301:     d__1 = sgn0 * (xstar - *xp);
2302:     *uxm = cvmgm_(uxm, &ustar, &d__1);
2303:     if (PetscRealPart(*xp) >= PetscRealPart(x2)) {
2304:         *utx = *utr;
2305:         *ubx = *ubr;
2306:         *rho1 = *rho1r;
2307:     } else {
2308:         *utx = *utl;
2309:         *ubx = *ubl;
2310:         *rho1 = *rho1l;
2311:     }
2312:     return iwave;
2313: }
2314: int godunovflux( const PetscScalar *ul, const PetscScalar *ur,
2315:                  PetscScalar *flux, const PetscReal *nn, const int *ndim,
2316:                  const PetscReal *gamma)
2317: {
2318:     /* System generated locals */
2319:   int i__1,iwave;
2320:     PetscScalar d__1, d__2, d__3;

2322:     /* Local variables */
2323:     static int k;
2324:     static PetscScalar bn[3], fn, ft, tg[3], pl, rl, pm, pr, rr, xp, ubl, ubm,
2325:             ubr, dtt, unm, tmp, utl, utm, uxl, utr, uxr, gaml, gamm, gamr,
2326:             xcen, rhom, rho1l, rho1m, rho1r;
2327:     /* Parameter adjustments */
2328:     --nn;
2329:     --flux;
2330:     --ur;
2331:     --ul;

2333:     /* Function Body */
2334:     xcen = 0.;
2335:     xp = 0.;
2336:     i__1 = *ndim;
2337:     for (k = 1; k <= i__1; ++k) {
2338:         tg[k - 1] = 0.;
2339:         bn[k - 1] = 0.;
2340:     }
2341:     dtt = 1.;
2342:     if (*ndim == 3) {
2343:         if (nn[1] == 0. && nn[2] == 0.) {
2344:             tg[0] = 1.;
2345:         } else {
2346:             tg[0] = -nn[2];
2347:             tg[1] = nn[1];
2348:         }
2349: /*           tmp=dsqrt(tg(1)**2+tg(2)**2) */
2350: /*           tg=tg/tmp */
2351:         bn[0] = -nn[3] * tg[1];
2352:         bn[1] = nn[3] * tg[0];
2353:         bn[2] = nn[1] * tg[1] - nn[2] * tg[0];
2354: /* Computing 2nd power */
2355:         d__1 = bn[0];
2356: /* Computing 2nd power */
2357:         d__2 = bn[1];
2358: /* Computing 2nd power */
2359:         d__3 = bn[2];
2360:         tmp = PetscSqrtScalar(d__1 * d__1 + d__2 * d__2 + d__3 * d__3);
2361:         i__1 = *ndim;
2362:         for (k = 1; k <= i__1; ++k) {
2363:             bn[k - 1] /= tmp;
2364:         }
2365:     } else if (*ndim == 2) {
2366:         tg[0] = -nn[2];
2367:         tg[1] = nn[1];
2368: /*           tmp=dsqrt(tg(1)**2+tg(2)**2) */
2369: /*           tg=tg/tmp */
2370:         bn[0] = 0.;
2371:         bn[1] = 0.;
2372:         bn[2] = 1.;
2373:     }
2374:     rl = ul[1];
2375:     rr = ur[1];
2376:     uxl = 0.;
2377:     uxr = 0.;
2378:     utl = 0.;
2379:     utr = 0.;
2380:     ubl = 0.;
2381:     ubr = 0.;
2382:     i__1 = *ndim;
2383:     for (k = 1; k <= i__1; ++k) {
2384:         uxl += ul[k + 1] * nn[k];
2385:         uxr += ur[k + 1] * nn[k];
2386:         utl += ul[k + 1] * tg[k - 1];
2387:         utr += ur[k + 1] * tg[k - 1];
2388:         ubl += ul[k + 1] * bn[k - 1];
2389:         ubr += ur[k + 1] * bn[k - 1];
2390:     }
2391:     uxl /= rl;
2392:     uxr /= rr;
2393:     utl /= rl;
2394:     utr /= rr;
2395:     ubl /= rl;
2396:     ubr /= rr;

2398:     gaml = *gamma;
2399:     gamr = *gamma;
2400: /* Computing 2nd power */
2401:     d__1 = uxl;
2402: /* Computing 2nd power */
2403:     d__2 = utl;
2404: /* Computing 2nd power */
2405:     d__3 = ubl;
2406:     pl = (*gamma - 1.) * (ul[*ndim + 2] - rl * .5 * (d__1 * d__1 + d__2 * d__2 + d__3 * d__3));
2407: /* Computing 2nd power */
2408:     d__1 = uxr;
2409: /* Computing 2nd power */
2410:     d__2 = utr;
2411: /* Computing 2nd power */
2412:     d__3 = ubr;
2413:     pr = (*gamma - 1.) * (ur[*ndim + 2] - rr * .5 * (d__1 * d__1 + d__2 * d__2 + d__3 * d__3));
2414:     rho1l = rl;
2415:     rho1r = rr;

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

2421:     flux[1] = rhom * unm;
2422:     fn = rhom * unm * unm + pm;
2423:     ft = rhom * unm * utm;
2424: /*           flux(2)=fn*nn(1)+ft*nn(2) */
2425: /*           flux(3)=fn*tg(1)+ft*tg(2) */
2426:     flux[2] = fn * nn[1] + ft * tg[0];
2427:     flux[3] = fn * nn[2] + ft * tg[1];
2428: /*           flux(2)=rhom*unm*(unm)+pm */
2429: /*           flux(3)=rhom*(unm)*utm */
2430:     if (*ndim == 3) {
2431:         flux[4] = rhom * unm * ubm;
2432:     }
2433:     flux[*ndim + 2] = (rhom * .5 * (unm * unm + utm * utm + ubm * ubm) + gamm / (gamm - 1.) * pm) * unm;
2434:     return iwave;
2435: } /* godunovflux_ */

2437: /* Subroutine to set up the initial conditions for the */
2438: /* Shock Interface interaction or linear wave (Ravi Samtaney,Mark Adams). */
2439: /* ----------------------------------------------------------------------- */
2440: int projecteqstate(PetscReal wc[], const PetscReal ueq[], PetscReal lv[][3])
2441: {
2442:   int j,k;
2443: /*      Wc=matmul(lv,Ueq) 3 vars */
2444:   for (k = 0; k < 3; ++k) {
2445:     wc[k] = 0.;
2446:     for (j = 0; j < 3; ++j) {
2447:       wc[k] += lv[k][j]*ueq[j];
2448:     }
2449:   }
2450:   return 0;
2451: }
2452: /* ----------------------------------------------------------------------- */
2453: int projecttoprim(PetscReal v[], const PetscReal wc[], PetscReal rv[][3])
2454: {
2455:   int k,j;
2456:   /*      V=matmul(rv,WC) 3 vars */
2457:   for (k = 0; k < 3; ++k) {
2458:     v[k] = 0.;
2459:     for (j = 0; j < 3; ++j) {
2460:       v[k] += rv[k][j]*wc[j];
2461:     }
2462:   }
2463:   return 0;
2464: }
2465: /* ---------------------------------------------------------------------- */
2466: int eigenvectors(PetscReal rv[][3], PetscReal lv[][3], const PetscReal ueq[], PetscReal gamma)
2467: {
2468:   int j,k;
2469:   PetscReal rho,csnd,p0;
2470:   /* PetscScalar u; */

2472:   for (k = 0; k < 3; ++k) for (j = 0; j < 3; ++j) { lv[k][j] = 0.; rv[k][j] = 0.; }
2473:   rho = ueq[0];
2474:   /* u = ueq[1]; */
2475:   p0 = ueq[2];
2476:   csnd = PetscSqrtReal(gamma * p0 / rho);
2477:   lv[0][1] = rho * .5;
2478:   lv[0][2] = -.5 / csnd;
2479:   lv[1][0] = csnd;
2480:   lv[1][2] = -1. / csnd;
2481:   lv[2][1] = rho * .5;
2482:   lv[2][2] = .5 / csnd;
2483:   rv[0][0] = -1. / csnd;
2484:   rv[1][0] = 1. / rho;
2485:   rv[2][0] = -csnd;
2486:   rv[0][1] = 1. / csnd;
2487:   rv[0][2] = 1. / csnd;
2488:   rv[1][2] = 1. / rho;
2489:   rv[2][2] = csnd;
2490:   return 0;
2491: }

2493: int initLinearWave(EulerNode *ux, const PetscReal gamma, const PetscReal coord[], const PetscReal Lx)
2494: {
2495:   PetscReal p0,u0,wcp[3],wc[3];
2496:   PetscReal lv[3][3];
2497:   PetscReal vp[3];
2498:   PetscReal rv[3][3];
2499:   PetscReal eps, ueq[3], rho0, twopi;

2501:   /* Function Body */
2502:   twopi = 2.*PETSC_PI;
2503:   eps = 1e-4; /* perturbation */
2504:   rho0 = 1e3;   /* density of water */
2505:   p0 = 101325.; /* init pressure of 1 atm (?) */
2506:   u0 = 0.;
2507:   ueq[0] = rho0;
2508:   ueq[1] = u0;
2509:   ueq[2] = p0;
2510:   /* Project initial state to characteristic variables */
2511:   eigenvectors(rv, lv, ueq, gamma);
2512:   projecteqstate(wc, ueq, lv);
2513:   wcp[0] = wc[0];
2514:   wcp[1] = wc[1];
2515:   wcp[2] = wc[2] + eps * PetscCosReal(coord[0] * 2. * twopi / Lx);
2516:   projecttoprim(vp, wcp, rv);
2517:   ux->r = vp[0]; /* density */
2518:   ux->ru[0] = vp[0] * vp[1]; /* x momentum */
2519:   ux->ru[1] = 0.;
2520: #if defined DIM > 2
2521:   if (dim>2) ux->ru[2] = 0.;
2522: #endif
2523:   /* E = rho * e + rho * v^2/2 = p/(gam-1) + rho*v^2/2 */
2524:   ux->E = vp[2]/(gamma - 1.) + 0.5*vp[0]*vp[1]*vp[1];
2525:   return 0;
2526: }

2528: /*TEST

2530:   # 2D Advection 0-10
2531:   test:
2532:     suffix: 0
2533:     requires: exodusii
2534:     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo

2536:   test:
2537:     suffix: 1
2538:     requires: exodusii
2539:     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo

2541:   test:
2542:     suffix: 2
2543:     requires: exodusii
2544:     nsize: 2
2545:     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo

2547:   test:
2548:     suffix: 3
2549:     requires: exodusii
2550:     nsize: 2
2551:     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo

2553:   test:
2554:     suffix: 4
2555:     requires: exodusii
2556:     nsize: 8
2557:     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo

2559:   test:
2560:     suffix: 5
2561:     requires: exodusii
2562:     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo -ts_type rosw -ts_adapt_reject_safety 1

2564:   test:
2565:     suffix: 6
2566:     requires: exodusii
2567:     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/squaremotor-30.exo -ufv_split_faces

2569:   test:
2570:     suffix: 7
2571:     requires: exodusii
2572:     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1

2574:   test:
2575:     suffix: 8
2576:     requires: exodusii
2577:     nsize: 2
2578:     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 2

2580:   test:
2581:     suffix: 9
2582:     requires: exodusii
2583:     nsize: 8
2584:     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 2

2586:   test:
2587:     suffix: 10
2588:     requires: exodusii
2589:     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo

2591:   # 2D Shallow water
2592:   test:
2593:     suffix: sw_0
2594:     requires: exodusii
2595:     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_final_time 1 -ts_ssp_type rks2 -ts_ssp_nstages 10 -monitor height,energy

2597:   # 2D Advection: p4est
2598:   test:
2599:     suffix: p4est_advec_2d
2600:     requires: p4est
2601:     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

2603:   # Advection in a box
2604:   test:
2605:     suffix: adv_2d_quad_0
2606:     args: -ufv_vtk_interval 0 -dm_refine 3 -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3

2608:   test:
2609:     suffix: adv_2d_quad_1
2610:     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
2611:     timeoutfactor: 3

2613:   test:
2614:     suffix: adv_2d_quad_p4est_0
2615:     requires: p4est
2616:     args: -ufv_vtk_interval 0 -dm_refine 5 -dm_type p4est -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3

2618:   test:
2619:     suffix: adv_2d_quad_p4est_1
2620:     requires: p4est
2621:     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
2622:     timeoutfactor: 3

2624:   test:
2625:     suffix: adv_2d_quad_p4est_adapt_0
2626:     requires: p4est
2627:     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_final_time 0.01
2628:     timeoutfactor: 3

2630:   test:
2631:     suffix: adv_2d_tri_0
2632:     requires: triangle
2633:     TODO: how did this ever get in master when there is no support for this
2634:     args: -ufv_vtk_interval 0 -simplex -dm_refine 3 -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3

2636:   test:
2637:     suffix: adv_2d_tri_1
2638:     requires: triangle
2639:     TODO: how did this ever get in master when there is no support for this
2640:     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
2641:   test:
2642:     suffix: adv_0
2643:     TODO: broken
2644:     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -bc_inflow 100,101,200 -bc_outflow 201

2646:   test:
2647:     suffix: shock_0
2648:     requires: p4est !single
2649:     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_final_time 0.5 -ts_ssp_type rks2 -ts_ssp_nstages 10 -ufv_vtk_basename ${wPETSC_DIR}/ex11
2650:     timeoutfactor: 3

2652:   # Test GLVis visualization of PetscFV fields
2653:   test:
2654:     suffix: glvis_adv_2d_tet
2655:     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

2657:   test:
2658:     suffix: glvis_adv_2d_quad
2659:     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

2661: TEST*/