Actual source code: ex11.c

petsc-master 2017-06-24
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: };

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

283: static PetscErrorCode SetUpBC_Advect(PetscDS prob, Physics phys)
284: {
286:   const PetscInt inflowids[] = {100,200,300},outflowids[] = {101};

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

497:   return(0);
498: }

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

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

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

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

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

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

597:   return(0);
598: }

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

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

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

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

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

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

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

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

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

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

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

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

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

808:   return(0);
809: }

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

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

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

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

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

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

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

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

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

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

923:   DMGetLabelIdIS(dm, labelName, &idIS);
924:   ISGetLocalSize(idIS, &numFS);
925:   ISGetIndices(idIS, &ids);

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1291: static PetscErrorCode FunctionalLinkDestroy(FunctionalLink *link)
1292: {
1294:   FunctionalLink l,next;

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

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

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

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

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

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

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

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

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

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

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

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

1462: static PetscErrorCode initializeTS(DM dm, User user, TS *ts)
1463: {

1467:   TSCreate(PetscObjectComm((PetscObject)dm), ts);
1468:   TSSetType(*ts, TSSSP);
1469:   TSSetDM(*ts, dm);
1470:   TSMonitorSet(*ts,MonitorVTK,user,NULL);
1471:   DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, user);
1472:   DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, user);
1473:   TSSetDuration(*ts,1000,2.0);
1474:   TSSetExactFinalTime(*ts,TS_EXACTFINALTIME_STEPOVER);
1475:   return(0);
1476: }

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

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

1530:     DMPlexPointLocalRead(gradDM,c,pointGrads,&pointGrad);
1531:     DMPlexPointLocalRead(cellDM,c,pointGeom,&cg);
1532:     DMPlexPointLocalRead(plex,c,pointVals,&pointVal);

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

1547:   VecTaggerComputeIS(refineTag,errVec,&refineIS);
1548:   VecTaggerComputeIS(coarsenTag,errVec,&coarsenIS);
1549:   ISGetSize(refineIS,&nRefine);
1550:   ISGetSize(coarsenIS,&nCoarsen);
1551:   if (nRefine) {DMLabelSetStratumIS(adaptLabel,DM_ADAPT_REFINE,refineIS);}
1552:   if (nCoarsen) {DMLabelSetStratumIS(adaptLabel,DM_ADAPT_COARSEN,coarsenIS);}
1553:   ISDestroy(&coarsenIS);
1554:   ISDestroy(&refineIS);
1555:   VecDestroy(&errVec);

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

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

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

1610:   PetscNew(&user);
1611:   PetscNew(&user->model);
1612:   PetscNew(&user->model->physics);
1613:   mod           = user->model;
1614:   phys          = mod->physics;
1615:   mod->comm     = comm;
1616:   useAMR        = PETSC_FALSE;
1617:   adaptInterval = 1;

1619:   /* Register physical models to be available on the command line */
1620:   PetscFunctionListAdd(&PhysicsList,"advect"          ,PhysicsCreate_Advect);
1621:   PetscFunctionListAdd(&PhysicsList,"sw"              ,PhysicsCreate_SW);
1622:   PetscFunctionListAdd(&PhysicsList,"euler"           ,PhysicsCreate_Euler);

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

1645:   if (useAMR) {
1646:     VecTaggerBox refineBox, coarsenBox;

1648:     refineBox.min  = refineBox.max  = PETSC_MAX_REAL;
1649:     coarsenBox.min = coarsenBox.max = PETSC_MIN_REAL;

1651:     VecTaggerCreate(comm,&refineTag);
1652:     PetscObjectSetOptionsPrefix((PetscObject)refineTag,"refine_");
1653:     VecTaggerSetType(refineTag,VECTAGGERABSOLUTE);
1654:     VecTaggerAbsoluteSetBox(refineTag,&refineBox);
1655:     VecTaggerSetFromOptions(refineTag);
1656:     VecTaggerSetUp(refineTag);
1657:     PetscObjectViewFromOptions((PetscObject)refineTag,NULL,"-tag_view");

1659:     VecTaggerCreate(comm,&coarsenTag);
1660:     PetscObjectSetOptionsPrefix((PetscObject)coarsenTag,"coarsen_");
1661:     VecTaggerSetType(coarsenTag,VECTAGGERABSOLUTE);
1662:     VecTaggerAbsoluteSetBox(coarsenTag,&coarsenBox);
1663:     VecTaggerSetFromOptions(coarsenTag);
1664:     VecTaggerSetUp(coarsenTag);
1665:     PetscObjectViewFromOptions((PetscObject)coarsenTag,NULL,"-tag_view");
1666:   }

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

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

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

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

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

1749:   mod->errorIndicator = ErrorIndicator_Simple;

1751:   {
1752:     DM dmDist;

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

1763:   DMSetFromOptions(dm);

1765:   {
1766:     DM gdm;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1902:   if (vtkCellGeom) {
1903:     DM  dmCell;
1904:     Vec cellgeom, partition;

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

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

1936:     TSGetDuration(ts,NULL,&finalTime);
1937:     TSSetDuration(ts,adaptInterval,finalTime);
1938:     TSSolve(ts,X);
1939:     TSGetSolveTime(ts,&ftime);
1940:     TSGetTimeStepNumber(ts,&nsteps);
1941:     for (adaptIter = 0;ftime < finalTime;adaptIter++) {
1942:       PetscLogDouble bytes;

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

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

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

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

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

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

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

2036:     static PetscScalar smallp = 1e-8;

2038:     /* System generated locals */
2039:     int i__1;
2040:     PetscScalar d__1, d__2;

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



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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

2526: /*TEST

2528:   # 2D Advection 0-10
2529:   test:
2530:     suffix: 0
2531:     requires: exodusii
2532:     args: -ufv_vtk_interval 0 -f ${PETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo
2533:   test:
2534:     suffix: 1
2535:     requires: exodusii
2536:     args: -ufv_vtk_interval 0 -f ${PETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo
2537:   test:
2538:     suffix: 2
2539:     requires: exodusii
2540:     nsize: 2
2541:     args: -ufv_vtk_interval 0 -f ${PETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo
2542:   test:
2543:     suffix: 3
2544:     requires: exodusii
2545:     nsize: 2
2546:     args: -ufv_vtk_interval 0 -f ${PETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo
2547:   test:
2548:     suffix: 4
2549:     requires: exodusii
2550:     nsize: 8
2551:     args: -ufv_vtk_interval 0 -f ${PETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo
2552:   test:
2553:     suffix: 5
2554:     requires: exodusii
2555:     args: -ufv_vtk_interval 0 -f ${PETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo -ts_type rosw -ts_adapt_reject_safety 1
2556:   test:
2557:     suffix: 6
2558:     requires: exodusii
2559:     args: -ufv_vtk_interval 0 -f ${PETSC_DIR}/share/petsc/datafiles/meshes/squaremotor-30.exo -ufv_split_faces
2560:   test:
2561:     suffix: 7
2562:     requires: exodusii
2563:     args: -ufv_vtk_interval 0 -f ${PETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1
2564:   test:
2565:     suffix: 8
2566:     requires: exodusii
2567:     nsize: 2
2568:     args: -ufv_vtk_interval 0 -f ${PETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 2
2569:   test:
2570:     suffix: 9
2571:     requires: exodusii
2572:     nsize: 8
2573:     args: -ufv_vtk_interval 0 -f ${PETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 2
2574:   test:
2575:     suffix: 10
2576:     requires: exodusii
2577:     args: -ufv_vtk_interval 0 -f ${PETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo
2578:   # 2D Shallow water
2579:   test:
2580:     suffix: sw_0
2581:     requires: exodusii
2582:     args: -ufv_vtk_interval 0 -f ${PETSC_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
2583:   # 2D Advection: p4est
2584:   test:
2585:     suffix: p4est_advec_2d
2586:     requires: p4est
2587:     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
2588:   # Advection in a box
2589:   test:
2590:     suffix: adv_2d_quad_0
2591:     requires:  !mpiuni
2592:     args: -ufv_vtk_interval 0 -dm_refine 5 -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3
2593:   test:
2594:     suffix: adv_2d_quad_1
2595:     requires:  !mpiuni
2596:     args: -ufv_vtk_interval 0 -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
2597:   test:
2598:     suffix: adv_2d_quad_p4est_0
2599:     requires: p4est
2600:     args: -ufv_vtk_interval 0 -dm_refine 5 -dm_type p4est -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3
2601:   test:
2602:     suffix: adv_2d_quad_p4est_1
2603:     requires: p4est
2604:     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
2605:   test:
2606:     suffix: adv_2d_quad_p4est_adapt_0
2607:     requires: p4est
2608:     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
2609:   test:
2610:     suffix: adv_2d_tri_0
2611:     requires: triangle  !mpiuni
2612:     args: -ufv_vtk_interval 0 -simplex -dm_refine 3 -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3
2613:   test:
2614:     suffix: adv_2d_tri_1
2615:     requires: triangle  !mpiuni
2616:     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
2617:   test:
2618:     suffix: adv_0
2619:     TODO: broken
2620:     args: -ufv_vtk_interval 0 -f ${PETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -bc_inflow 100,101,200 -bc_outflow 201

2622: TEST*/