Actual source code: ex11.c

petsc-master 2017-02-25
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:   PetscInt monitorStepOffset;
111:   Model    model;
112: };

114: PETSC_STATIC_INLINE PetscScalar DotDIM(const PetscScalar *x,const PetscScalar *y)
115: {
116:   PetscInt    i;
117:   PetscScalar prod=0.0;

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

124: PETSC_STATIC_INLINE PetscScalar Dot2(const PetscScalar *x,const PetscScalar *y) { return x[0]*y[0] + x[1]*y[1];}
125: PETSC_STATIC_INLINE PetscReal Norm2(const PetscScalar *x) { return PetscSqrtReal(PetscAbsScalar(Dot2(x,x)));}
126: PETSC_STATIC_INLINE void Normalize2(PetscScalar *x) { PetscReal a = 1./Norm2(x); x[0] *= a; x[1] *= a; }
127: PETSC_STATIC_INLINE void Waxpy2(PetscScalar a,const PetscScalar *x,const PetscScalar *y,PetscScalar *w) { w[0] = a*x[0] + y[0]; w[1] = a*x[1] + y[1]; }
128: PETSC_STATIC_INLINE void Scale2(PetscScalar a,const PetscScalar *x,PetscScalar *y) { y[0] = a*x[0]; y[1] = a*x[1]; }

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

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

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

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

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

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

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

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

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

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

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

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

225: static PetscErrorCode PhysicsSolution_Advect(Model mod,PetscReal time,const PetscReal *x,PetscScalar *u,void *ctx)
226: {
227:   Physics        phys    = (Physics)ctx;
228:   Physics_Advect *advect = (Physics_Advect*)phys->data;

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

263: static PetscErrorCode PhysicsFunctional_Advect(Model mod,PetscReal time,const PetscScalar *x,const PetscScalar *y,PetscReal *f,void *ctx)
264: {
265:   Physics        phys    = (Physics)ctx;
266:   Physics_Advect *advect = (Physics_Advect*)phys->data;
267:   PetscScalar    yexact[1];

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

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

283:   /* Register "canned" boundary conditions and defaults for where to apply. */
284:   PetscDSAddBoundary(prob, PETSC_TRUE, "inflow",  "Face Sets", 0, 0, NULL, (void (*)()) PhysicsBoundary_Advect_Inflow,  ALEN(inflowids),   inflowids,  phys);
285:   PetscDSAddBoundary(prob, PETSC_FALSE, "outflow", "Face Sets", 0, 0, NULL, (void (*)()) PhysicsBoundary_Advect_Outflow, ALEN(outflowids), outflowids, phys);
286:   return(0);
287: }

289: static PetscErrorCode PhysicsCreate_Advect(Model mod,Physics phys,PetscOptionItems *PetscOptionsObject)
290: {
291:   Physics_Advect *advect;

295:   phys->field_desc = PhysicsFields_Advect;
296:   phys->riemann    = (PetscRiemannFunc)PhysicsRiemann_Advect;
297:   PetscNew(&advect);
298:   phys->data       = advect;
299:   mod->setupbc = SetUpBC_Advect;

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

342: /******************* Shallow Water ********************/
343: typedef struct {
344:   PetscReal gravity;
345:   PetscReal boundaryHeight;
346:   struct {
347:     PetscInt Height;
348:     PetscInt Speed;
349:     PetscInt Energy;
350:   } functional;
351: } Physics_SW;
352: typedef struct {
353:   PetscScalar vals[0];
354:   PetscScalar h;
355:   PetscScalar uh[DIM];
356: } SWNode;

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

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

372:   Scale2(1./x->h,x->uh,u);
373:   uhn  = Dot2(x->uh,n);
374:   f->h = uhn;
375:   for (i=0; i<DIM; i++) f->uh[i] = u[i] * uhn + sw->gravity * PetscSqr(x->h) * n[i];
376:   return(0);
377: }

379: static PetscErrorCode PhysicsBoundary_SW_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
380: {
382:   xG[0] = xI[0];
383:   xG[1] = -xI[1];
384:   xG[2] = -xI[2];
385:   return(0);
386: }

388: static void PhysicsRiemann_SW(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscScalar *flux, Physics phys)
389: {
390:   Physics_SW   *sw = (Physics_SW*)phys->data;
391:   PetscReal    cL,cR,speed,nn[DIM];
392:   const SWNode *uL = (const SWNode*)xL,*uR = (const SWNode*)xR;
393:   SWNode       fL,fR;
394:   PetscInt     i;

396:   if (uL->h < 0 || uR->h < 0) {for (i=0; i<1+dim; i++) flux[i] = NAN; return;} /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reconstructed thickness is negative"); */
397:   nn[0] = n[0];
398:   nn[1] = n[1];
399:   Normalize2(nn);
400:   SWFlux(phys,nn,uL,&fL);
401:   SWFlux(phys,nn,uR,&fR);
402:   cL    = PetscSqrtReal(sw->gravity*PetscRealPart(uL->h));
403:   cR    = PetscSqrtReal(sw->gravity*PetscRealPart(uR->h)); /* gravity wave speed */
404:   speed = PetscMax(PetscAbsScalar(Dot2(uL->uh,nn)/uL->h) + cL,PetscAbsScalar(Dot2(uR->uh,nn)/uR->h) + cR);
405:   for (i=0; i<1+dim; i++) flux[i] = (0.5*(fL.vals[i] + fR.vals[i]) + 0.5*speed*(xL[i] - xR[i])) * Norm2(n);
406: }

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

413:   if (time != 0.0) SETERRQ1(mod->comm,PETSC_ERR_SUP,"No solution known for time %g",(double)time);
414:   dx[0] = x[0] - 1.5;
415:   dx[1] = x[1] - 1.0;
416:   r     = Norm2(dx);
417:   sigma = 0.5;
418:   u[0]  = 1 + 2*PetscExpScalar(-PetscSqr(r)/(2*PetscSqr(sigma)));
419:   u[1]  = 0.0;
420:   u[2]  = 0.0;
421:   return(0);
422: }

424: static PetscErrorCode PhysicsFunctional_SW(Model mod,PetscReal time,const PetscReal *coord,const PetscScalar *xx,PetscReal *f,void *ctx)
425: {
426:   Physics      phys = (Physics)ctx;
427:   Physics_SW   *sw  = (Physics_SW*)phys->data;
428:   const SWNode *x   = (const SWNode*)xx;
429:   PetscScalar  u[2];
430:   PetscReal    h;

433:   h = PetscRealPart(x->h);
434:   Scale2(1./x->h,x->uh,u);
435:   f[sw->functional.Height] = h;
436:   f[sw->functional.Speed]  = Norm2(u) + PetscSqrtReal(sw->gravity*h);
437:   f[sw->functional.Energy] = 0.5*(Dot2(x->uh,u) + sw->gravity*PetscSqr(h));
438:   return(0);
439: }

441: static PetscErrorCode SetUpBC_SW(PetscDS prob,Physics phys)
442: {
444:   const PetscInt wallids[] = {100,101,200,300};
446:   PetscDSAddBoundary(prob, PETSC_FALSE, "wall", "Face Sets", 0, 0, NULL, (void (*)()) PhysicsBoundary_SW_Wall, ALEN(wallids), wallids, phys);
447:   return(0);
448: }

450: static PetscErrorCode PhysicsCreate_SW(Model mod,Physics phys,PetscOptionItems *PetscOptionsObject)
451: {
452:   Physics_SW     *sw;

456:   phys->field_desc = PhysicsFields_SW;
457:   phys->riemann = (PetscRiemannFunc) PhysicsRiemann_SW;
458:   PetscNew(&sw);
459:   phys->data    = sw;
460:   mod->setupbc  = SetUpBC_SW;

462:   PetscOptionsHead(PetscOptionsObject,"SW options");
463:   {
464:     sw->gravity = 1.0;
465:     PetscOptionsReal("-sw_gravity","Gravitational constant","",sw->gravity,&sw->gravity,NULL);
466:   }
467:   PetscOptionsTail();
468:   phys->maxspeed = PetscSqrtReal(2.0*sw->gravity); /* Mach 1 for depth of 2 */

470:   ModelSolutionSetDefault(mod,PhysicsSolution_SW,phys);
471:   ModelFunctionalRegister(mod,"Height",&sw->functional.Height,PhysicsFunctional_SW,phys);
472:   ModelFunctionalRegister(mod,"Speed",&sw->functional.Speed,PhysicsFunctional_SW,phys);
473:   ModelFunctionalRegister(mod,"Energy",&sw->functional.Energy,PhysicsFunctional_SW,phys);

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

477:   return(0);
478: }

480: /******************* Euler Density Shock (EULER_IV_SHOCK,EULER_SS_SHOCK) ********************/
481: /* An initial-value and self-similar solutions of the compressible Euler equations */
482: /* Ravi Samtaney and D. I. Pullin */
483: /* Phys. Fluids 8, 2650 (1996); http://dx.doi.org/10.1063/1.869050 */
484: typedef enum {EULER_PAR_GAMMA,EULER_PAR_RHOR,EULER_PAR_AMACH,EULER_PAR_ITANA,EULER_PAR_SIZE} EulerParamIdx;
485: typedef enum {EULER_IV_SHOCK,EULER_SS_SHOCK,EULER_SHOCK_TUBE,EULER_LINEAR_WAVE} EulerType;
486: typedef struct {
487:   PetscScalar vals[0];
488:   PetscScalar r;
489:   PetscScalar ru[DIM];
490:   PetscScalar E;
491: } EulerNode;
492: typedef PetscErrorCode (*EquationOfState)(const PetscReal*, const EulerNode*, PetscScalar*);
493: typedef struct {
494:   EulerType       type;
495:   PetscReal       pars[EULER_PAR_SIZE];
496:   EquationOfState sound;
497:   struct {
498:     PetscInt Density;
499:     PetscInt Momentum;
500:     PetscInt Energy;
501:     PetscInt Pressure;
502:     PetscInt Speed;
503:   } monitor;
504: } Physics_Euler;

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

508: /* initial condition */
509: int initLinearWave(EulerNode *ux, const PetscScalar gamma, const PetscReal coord[], const PetscReal Lx);
510: static PetscErrorCode PhysicsSolution_Euler(Model mod, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx)
511: {
512:   PetscInt i;
513:   Physics         phys = (Physics)ctx;
514:   Physics_Euler   *eu  = (Physics_Euler*)phys->data;
515:   EulerNode       *uu  = (EulerNode*)u;
516:   PetscScalar p0,gamma,c;
518:   if (time != 0.0) SETERRQ1(mod->comm,PETSC_ERR_SUP,"No solution known for time %g",(double)time);

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

524:   if (eu->type==EULER_IV_SHOCK || eu->type==EULER_SS_SHOCK) {
525:     /******************* Euler Density Shock ********************/
526:     /* On initial-value and self-similar solutions of the compressible Euler equations */
527:     /* Ravi Samtaney and D. I. Pullin */
528:     /* Phys. Fluids 8, 2650 (1996); http://dx.doi.org/10.1063/1.869050 */
529:     /* initial conditions 1: left of shock, 0: left of discontinuity 2: right of discontinuity,  */
530:     p0 = 1.;
531:     if (x[0] < 0.0 + x[1]*eu->pars[EULER_PAR_ITANA]) {
532:       if (x[0] < mod->bounds[0]*0.5) { /* left of shock (1) */
533:         PetscScalar amach,rho,press,gas1,p1;
534:         amach = eu->pars[EULER_PAR_AMACH];
535:         rho = 1.;
536:         press = p0;
537:         p1 = press*(1.0+2.0*gamma/(gamma+1.0)*(amach*amach-1.0));
538:         gas1 = (gamma-1.0)/(gamma+1.0);
539:         uu->r = rho*(p1/press+gas1)/(gas1*p1/press+1.0);
540:         uu->ru[0]   = ((uu->r - rho)*sqrt(gamma*press/rho)*amach);
541:         uu->E = p1/(gamma-1.0) + .5/uu->r*uu->ru[0]*uu->ru[0];
542:       }
543:       else { /* left of discontinuity (0) */
544:         uu->r = 1.; /* rho = 1 */
545:         uu->E = p0/(gamma-1.0);
546:       }
547:     }
548:     else { /* right of discontinuity (2) */
549:       uu->r = eu->pars[EULER_PAR_RHOR];
550:       uu->E = p0/(gamma-1.0);
551:     }
552:   }
553:   else if (eu->type==EULER_SHOCK_TUBE) {
554:     /* 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. */
555:     if (x[0] < 0.0 ) {
556:       uu->r = 8.;
557:       uu->E = 10./(gamma-1.);
558:     }
559:     else {
560:       uu->r = 1.;
561:       uu->E = 1./(gamma-1.);
562:     }
563:   }
564:   else if (eu->type==EULER_LINEAR_WAVE) {
565:     initLinearWave( uu, gamma, x, mod->bounds[1] - mod->bounds[0]);
566:   }
567:   else SETERRQ1(mod->comm,PETSC_ERR_SUP,"Unknown type %d",eu->type);

569:   // set phys->maxspeed: (mod->maxspeed = phys->maxspeed) in main;
570:   eu->sound(&gamma,uu,&c);
571:   c = PetscAbsScalar(uu->ru[0]/uu->r) + c;
572:   if (c > phys->maxspeed) phys->maxspeed = c;

574:   return(0);
575: }

577: static PetscErrorCode Pressure_PG(const PetscReal gamma,const EulerNode *x,PetscScalar *p)
578: {
579:   PetscScalar ru2;
581:   ru2  = DotDIM(x->ru,x->ru);
582:   (*p)=(x->E - 0.5*ru2/x->r)*(gamma - 1.0); /* (E - rho V^2/2)(gamma-1) = e rho (gamma-1) */
583:   return(0);
584: }

586: static PetscErrorCode SpeedOfSound_PG(const PetscReal *gamma, const EulerNode *x, PetscScalar *c)
587: {
588:   PetscScalar p;

591:   Pressure_PG(*gamma,x,&p);
592:   if (p<0.) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"negative pressure time %g -- NEED TO FIX!!!!!!",p);
593:   /* pars[EULER_PAR_GAMMA] = heat capacity ratio */
594:   (*c)=PetscSqrtScalar(*gamma * p / x->r);
595:   return(0);
596: }

598: /*
599:  * x = (rho,rho*(u_1),...,rho*e)^T
600:  * x_t+div(f_1(x))+...+div(f_DIM(x)) = 0
601:  *
602:  * f_i(x) = u_i*x+(0,0,...,p,...,p*u_i)^T
603:  *
604:  */
605: static PetscErrorCode EulerFlux(Physics phys,const PetscReal *n,const EulerNode *x,EulerNode *f)
606: {
607:   Physics_Euler *eu = (Physics_Euler*)phys->data;
608:   PetscScalar   nu,p;
609:   PetscInt      i;

612:   Pressure_PG(eu->pars[EULER_PAR_GAMMA],x,&p);
613:   nu = DotDIM(x->ru,n);
614:   f->r = nu;   /* A rho u */
615:   nu /= x->r;  /* A u */
616:   for (i=0; i<DIM; i++) f->ru[i] = nu * x->ru[i] + n[i]*p;  /* r u^2 + p */
617:   f->E = nu * (x->E + p); /* u(e+p) */
618:   return(0);
619: }

621: /* PetscReal* => EulerNode* conversion */
622: static PetscErrorCode PhysicsBoundary_Euler_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *a_xI, PetscScalar *a_xG, void *ctx)
623: {
624:   PetscInt    i;
625:   const EulerNode *xI = (const EulerNode*)a_xI;
626:   EulerNode       *xG = (EulerNode*)a_xG;
627:   Physics         phys = (Physics)ctx;
628:   Physics_Euler   *eu  = (Physics_Euler*)phys->data;
630:   xG->r = xI->r;           /* ghost cell density - same */
631:   xG->E = xI->E;           /* ghost cell energy - same */
632:   if (n[1] != 0.) {        /* top and bottom */
633:     xG->ru[0] =  xI->ru[0]; /* copy tang to wall */
634:     xG->ru[1] = -xI->ru[1]; /* reflect perp to t/b wall */
635:   }
636:   else { /* sides */
637:     for (i=0; i<DIM; i++) xG->ru[i] = xI->ru[i]; /* copy */
638:   }
639:   if (eu->type == EULER_LINEAR_WAVE) { // debug
640:     // PetscPrintf(PETSC_COMM_WORLD,"%s coord=%g,%g\n",PETSC_FUNCTION_NAME,c[0],c[1]);
641:   }
642:   return(0);
643: }
644: int godunovflux( const PetscScalar *ul, const PetscScalar *ur, PetscScalar *flux, const PetscReal *nn, const int *ndim, const PetscReal *gamma);
645: /* PetscReal* => EulerNode* conversion */
646: static void PhysicsRiemann_Euler_Godunov( PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n,
647:                                           const PetscScalar *xL, const PetscScalar *xR, PetscScalar *flux, Physics phys)
648: {
649:   Physics_Euler   *eu = (Physics_Euler*)phys->data;
650:   PetscReal       cL,cR,speed,velL,velR,nn[DIM],s2;
651:   PetscInt        i;
652:   PetscErrorCode  ierr;

655:   for (i=0,s2=0.; i<DIM; i++) {
656:     nn[i] = n[i];
657:     s2 += nn[i]*nn[i];
658:   }
659:   s2 = PetscSqrtReal(s2); /* |n|_2 = sum(n^2)^1/2 */
660:   for (i=0.; i<DIM; i++) nn[i] /= s2;
661:   if (0) { /* Rusanov */
662:     const EulerNode *uL = (const EulerNode*)xL,*uR = (const EulerNode*)xR;
663:     EulerNode       fL,fR;
664:     EulerFlux(phys,nn,uL,&fL);
665:     EulerFlux(phys,nn,uR,&fR);
666:     eu->sound(&eu->pars[EULER_PAR_GAMMA],uL,&cL);if (ierr) exit(13);
667:     eu->sound(&eu->pars[EULER_PAR_GAMMA],uR,&cR);if (ierr) exit(14);
668:     velL = DotDIM(uL->ru,nn)/uL->r;
669:     velR = DotDIM(uR->ru,nn)/uR->r;
670:     speed = PetscMax(PetscAbsScalar(velR) + cR,PetscAbsScalar(velL) + cL);
671:     for (i=0; i<2+dim; i++) flux[i] = 0.5*((fL.vals[i]+fR.vals[i]) + speed*(xL[i] - xR[i]))*s2;
672:   }
673:   else {
674:     int dim = DIM;
675:     /* int iwave =  */
676:     godunovflux(xL, xR, flux, nn, &dim, &eu->pars[EULER_PAR_GAMMA]);
677:     for (i=0; i<2+dim; i++) flux[i] *= s2;
678:   }
679:   PetscFunctionReturnVoid();
680: }

682: static PetscErrorCode PhysicsFunctional_Euler(Model mod,PetscReal time,const PetscReal *coord,const PetscScalar *xx,PetscReal *f,void *ctx)
683: {
684:   Physics         phys = (Physics)ctx;
685:   Physics_Euler   *eu  = (Physics_Euler*)phys->data;
686:   const EulerNode *x   = (const EulerNode*)xx;
687:   PetscScalar     p;

690:   f[eu->monitor.Density]  = x->r;
691:   f[eu->monitor.Momentum] = NormDIM(x->ru);
692:   f[eu->monitor.Energy]   = x->E;
693:   f[eu->monitor.Speed]    = NormDIM(x->ru)/x->r;
694:   Pressure_PG(eu->pars[EULER_PAR_GAMMA], x, &p);
695:   f[eu->monitor.Pressure] = p;
696:   return(0);
697: }

699: static PetscErrorCode SetUpBC_Euler(PetscDS prob,Physics phys)
700: {
701:   PetscErrorCode  ierr;
702:   Physics_Euler   *eu = phys->data;
703:   if (eu->type == EULER_LINEAR_WAVE) {
704:     const PetscInt wallids[] = {100,101};
705:     PetscDSAddBoundary(prob, PETSC_FALSE, "wall", "Face Sets", 0, 0, NULL, (void (*)()) PhysicsBoundary_Euler_Wall, ALEN(wallids), wallids, phys);
706:   }
707:   else {
708:     const PetscInt wallids[] = {100,101,200,300};
709:     PetscDSAddBoundary(prob, PETSC_FALSE, "wall", "Face Sets", 0, 0, NULL, (void (*)()) PhysicsBoundary_Euler_Wall, ALEN(wallids), wallids, phys);
710:   }
711:   return(0);
712: }

714: static PetscErrorCode PhysicsCreate_Euler(Model mod,Physics phys,PetscOptionItems *PetscOptionsObject)
715: {
716:   Physics_Euler   *eu;
717:   PetscErrorCode  ierr;

720:   phys->field_desc = PhysicsFields_Euler;
721:   phys->riemann = (PetscRiemannFunc) PhysicsRiemann_Euler_Godunov;
722:   PetscNew(&eu);
723:   phys->data    = eu;
724:   mod->setupbc = SetUpBC_Euler;
725:   PetscOptionsHead(PetscOptionsObject,"Euler options");
726:   {
727:     PetscReal alpha;
728:     char type[64] = "linear_wave";
729:     PetscBool  is;
730:     mod->bcs[0] = mod->bcs[1] = mod->bcs[2] = DM_BOUNDARY_GHOSTED;
731:     eu->pars[EULER_PAR_GAMMA] = 1.4;
732:     eu->pars[EULER_PAR_AMACH] = 2.02;
733:     eu->pars[EULER_PAR_RHOR] = 3.0;
734:     eu->pars[EULER_PAR_ITANA] = 0.57735026918963; /* angle of Euler self similar (SS) shock */
735:     PetscOptionsReal("-eu_gamma","Heat capacity ratio","",eu->pars[EULER_PAR_GAMMA],&eu->pars[EULER_PAR_GAMMA],NULL);
736:     PetscOptionsReal("-eu_amach","Shock speed (Mach)","",eu->pars[EULER_PAR_AMACH],&eu->pars[EULER_PAR_AMACH],NULL);
737:     PetscOptionsReal("-eu_rho2","Density right of discontinuity","",eu->pars[EULER_PAR_RHOR],&eu->pars[EULER_PAR_RHOR],NULL);
738:     alpha = 60.;
739:     PetscOptionsReal("-eu_alpha","Angle of discontinuity","",alpha,&alpha,NULL);
740:     if (alpha<=0. || alpha>90.) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Alpha bust be > 0 and <= 90 (%g)",alpha);
741:     eu->pars[EULER_PAR_ITANA] = 1./tan ( alpha * M_PI / 180.0 );
742:     PetscOptionsString("-eu_type","Type of Euler test","",type,type,sizeof(type),NULL);
743:     PetscStrcmp(type,"linear_wave", &is);
744:     if (is) {
745:       eu->type = EULER_LINEAR_WAVE;
746:       mod->bcs[0] = mod->bcs[1] = mod->bcs[2] = DM_BOUNDARY_PERIODIC;
747:       mod->bcs[1] = DM_BOUNDARY_GHOSTED; /* debug */
748:       PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"linear_wave");
749:     }
750:     else {
751:       if (DIM != 2) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DIM must be 2 unless linear wave test %s",type);
752:       PetscStrcmp(type,"iv_shock", &is);
753:       if (is) {
754:         eu->type = EULER_IV_SHOCK;
755:         PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"iv_shock");
756:       }
757:       else {
758:         PetscStrcmp(type,"ss_shock", &is);
759:         if (is) {
760:           eu->type = EULER_SS_SHOCK;
761:           PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"ss_shock");
762:         }
763:         else {
764:           PetscStrcmp(type,"shock_tube", &is);
765:           if (is) eu->type = EULER_SHOCK_TUBE;
766:           else SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Unknown Euler type %s",type);
767:           PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"shock_tube");
768:         }
769:       }
770:     }
771:   }
772:   PetscOptionsTail();
773:   eu->sound = SpeedOfSound_PG;
774:   phys->maxspeed = 0.; /* will get set in solution */
775:   ModelSolutionSetDefault(mod,PhysicsSolution_Euler,phys);
776:   ModelFunctionalRegister(mod,"Speed",&eu->monitor.Speed,PhysicsFunctional_Euler,phys);
777:   ModelFunctionalRegister(mod,"Energy",&eu->monitor.Energy,PhysicsFunctional_Euler,phys);
778:   ModelFunctionalRegister(mod,"Density",&eu->monitor.Density,PhysicsFunctional_Euler,phys);
779:   ModelFunctionalRegister(mod,"Momentum",&eu->monitor.Momentum,PhysicsFunctional_Euler,phys);
780:   ModelFunctionalRegister(mod,"Pressure",&eu->monitor.Pressure,PhysicsFunctional_Euler,phys);

782:   return(0);
783: }

785: static PetscErrorCode ErrorIndicator_Simple(PetscInt dim, PetscReal volume, PetscInt numComps, const PetscScalar u[], const PetscScalar grad[], PetscReal *error, void *ctx)
786: {
787:   PetscReal      err = 0.;
788:   PetscInt       i, j;

791:   for (i = 0; i < numComps; i++) {
792:     for (j = 0; j < dim; j++) {
793:       err += PetscSqr(PetscRealPart(grad[i * dim + j]));
794:     }
795:   }
796:   *error = volume * err;
797:   return(0);
798: }

800: PetscErrorCode ConstructCellBoundary(DM dm, User user)
801: {
802:   const char     *name   = "Cell Sets";
803:   const char     *bdname = "split faces";
804:   IS             regionIS, innerIS;
805:   const PetscInt *regions, *cells;
806:   PetscInt       numRegions, innerRegion, numCells, c;
807:   PetscInt       cStart, cEnd, cEndInterior, fStart, fEnd;
808:   PetscBool      hasLabel;

812:   DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
813:   DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);
814:   DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);

816:   DMHasLabel(dm, name, &hasLabel);
817:   if (!hasLabel) return(0);
818:   DMGetLabelSize(dm, name, &numRegions);
819:   if (numRegions != 2) return(0);
820:   /* Get the inner id */
821:   DMGetLabelIdIS(dm, name, &regionIS);
822:   ISGetIndices(regionIS, &regions);
823:   innerRegion = regions[0];
824:   ISRestoreIndices(regionIS, &regions);
825:   ISDestroy(&regionIS);
826:   /* Find the faces between cells in different regions, could call DMPlexCreateNeighborCSR() */
827:   DMGetStratumIS(dm, name, innerRegion, &innerIS);
828:   ISGetLocalSize(innerIS, &numCells);
829:   ISGetIndices(innerIS, &cells);
830:   DMCreateLabel(dm, bdname);
831:   for (c = 0; c < numCells; ++c) {
832:     const PetscInt cell = cells[c];
833:     const PetscInt *faces;
834:     PetscInt       numFaces, f;

836:     if ((cell < cStart) || (cell >= cEnd)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Got invalid point %d which is not a cell", cell);
837:     DMPlexGetConeSize(dm, cell, &numFaces);
838:     DMPlexGetCone(dm, cell, &faces);
839:     for (f = 0; f < numFaces; ++f) {
840:       const PetscInt face = faces[f];
841:       const PetscInt *neighbors;
842:       PetscInt       nC, regionA, regionB;

844:       if ((face < fStart) || (face >= fEnd)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Got invalid point %d which is not a face", face);
845:       DMPlexGetSupportSize(dm, face, &nC);
846:       if (nC != 2) continue;
847:       DMPlexGetSupport(dm, face, &neighbors);
848:       if ((neighbors[0] >= cEndInterior) || (neighbors[1] >= cEndInterior)) continue;
849:       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]);
850:       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]);
851:       DMGetLabelValue(dm, name, neighbors[0], &regionA);
852:       DMGetLabelValue(dm, name, neighbors[1], &regionB);
853:       if (regionA < 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid label %s: Cell %d has no value", name, neighbors[0]);
854:       if (regionB < 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid label %s: Cell %d has no value", name, neighbors[1]);
855:       if (regionA != regionB) {
856:         DMSetLabelValue(dm, bdname, faces[f], 1);
857:       }
858:     }
859:   }
860:   ISRestoreIndices(innerIS, &cells);
861:   ISDestroy(&innerIS);
862:   {
863:     DMLabel label;

865:     DMGetLabel(dm, bdname, &label);
866:     DMLabelView(label, PETSC_VIEWER_STDOUT_WORLD);
867:   }
868:   return(0);
869: }

871: /* Right now, I have just added duplicate faces, which see both cells. We can
872: - Add duplicate vertices and decouple the face cones
873: - Disconnect faces from cells across the rotation gap
874: */
875: PetscErrorCode SplitFaces(DM *dmSplit, const char labelName[], User user)
876: {
877:   DM             dm = *dmSplit, sdm;
878:   PetscSF        sfPoint, gsfPoint;
879:   PetscSection   coordSection, newCoordSection;
880:   Vec            coordinates;
881:   IS             idIS;
882:   const PetscInt *ids;
883:   PetscInt       *newpoints;
884:   PetscInt       dim, depth, maxConeSize, maxSupportSize, numLabels, numGhostCells;
885:   PetscInt       numFS, fs, pStart, pEnd, p, cEnd, cEndInterior, vStart, vEnd, v, fStart, fEnd, newf, d, l;
886:   PetscBool      hasLabel;

890:   DMHasLabel(dm, labelName, &hasLabel);
891:   if (!hasLabel) return(0);
892:   DMCreate(PetscObjectComm((PetscObject)dm), &sdm);
893:   DMSetType(sdm, DMPLEX);
894:   DMGetDimension(dm, &dim);
895:   DMSetDimension(sdm, dim);

897:   DMGetLabelIdIS(dm, labelName, &idIS);
898:   ISGetLocalSize(idIS, &numFS);
899:   ISGetIndices(idIS, &ids);

901:   user->numSplitFaces = 0;
902:   for (fs = 0; fs < numFS; ++fs) {
903:     PetscInt numBdFaces;

905:     DMGetStratumSize(dm, labelName, ids[fs], &numBdFaces);
906:     user->numSplitFaces += numBdFaces;
907:   }
908:   DMPlexGetChart(dm, &pStart, &pEnd);
909:   pEnd += user->numSplitFaces;
910:   DMPlexSetChart(sdm, pStart, pEnd);
911:   DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);
912:   DMPlexSetHybridBounds(sdm, cEndInterior, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE);
913:   DMPlexGetHeightStratum(dm, 0, NULL, &cEnd);
914:   numGhostCells = cEnd - cEndInterior;
915:   /* Set cone and support sizes */
916:   DMPlexGetDepth(dm, &depth);
917:   for (d = 0; d <= depth; ++d) {
918:     DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);
919:     for (p = pStart; p < pEnd; ++p) {
920:       PetscInt newp = p;
921:       PetscInt size;

923:       DMPlexGetConeSize(dm, p, &size);
924:       DMPlexSetConeSize(sdm, newp, size);
925:       DMPlexGetSupportSize(dm, p, &size);
926:       DMPlexSetSupportSize(sdm, newp, size);
927:     }
928:   }
929:   DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);
930:   for (fs = 0, newf = fEnd; fs < numFS; ++fs) {
931:     IS             faceIS;
932:     const PetscInt *faces;
933:     PetscInt       numFaces, f;

935:     DMGetStratumIS(dm, labelName, ids[fs], &faceIS);
936:     ISGetLocalSize(faceIS, &numFaces);
937:     ISGetIndices(faceIS, &faces);
938:     for (f = 0; f < numFaces; ++f, ++newf) {
939:       PetscInt size;

941:       /* Right now I think that both faces should see both cells */
942:       DMPlexGetConeSize(dm, faces[f], &size);
943:       DMPlexSetConeSize(sdm, newf, size);
944:       DMPlexGetSupportSize(dm, faces[f], &size);
945:       DMPlexSetSupportSize(sdm, newf, size);
946:     }
947:     ISRestoreIndices(faceIS, &faces);
948:     ISDestroy(&faceIS);
949:   }
950:   DMSetUp(sdm);
951:   /* Set cones and supports */
952:   DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);
953:   PetscMalloc1(PetscMax(maxConeSize, maxSupportSize), &newpoints);
954:   DMPlexGetChart(dm, &pStart, &pEnd);
955:   for (p = pStart; p < pEnd; ++p) {
956:     const PetscInt *points, *orientations;
957:     PetscInt       size, i, newp = p;

959:     DMPlexGetConeSize(dm, p, &size);
960:     DMPlexGetCone(dm, p, &points);
961:     DMPlexGetConeOrientation(dm, p, &orientations);
962:     for (i = 0; i < size; ++i) newpoints[i] = points[i];
963:     DMPlexSetCone(sdm, newp, newpoints);
964:     DMPlexSetConeOrientation(sdm, newp, orientations);
965:     DMPlexGetSupportSize(dm, p, &size);
966:     DMPlexGetSupport(dm, p, &points);
967:     for (i = 0; i < size; ++i) newpoints[i] = points[i];
968:     DMPlexSetSupport(sdm, newp, newpoints);
969:   }
970:   PetscFree(newpoints);
971:   for (fs = 0, newf = fEnd; fs < numFS; ++fs) {
972:     IS             faceIS;
973:     const PetscInt *faces;
974:     PetscInt       numFaces, f;

976:     DMGetStratumIS(dm, labelName, ids[fs], &faceIS);
977:     ISGetLocalSize(faceIS, &numFaces);
978:     ISGetIndices(faceIS, &faces);
979:     for (f = 0; f < numFaces; ++f, ++newf) {
980:       const PetscInt *points;

982:       DMPlexGetCone(dm, faces[f], &points);
983:       DMPlexSetCone(sdm, newf, points);
984:       DMPlexGetSupport(dm, faces[f], &points);
985:       DMPlexSetSupport(sdm, newf, points);
986:     }
987:     ISRestoreIndices(faceIS, &faces);
988:     ISDestroy(&faceIS);
989:   }
990:   ISRestoreIndices(idIS, &ids);
991:   ISDestroy(&idIS);
992:   DMPlexStratify(sdm);
993:   /* Convert coordinates */
994:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
995:   DMGetCoordinateSection(dm, &coordSection);
996:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &newCoordSection);
997:   PetscSectionSetNumFields(newCoordSection, 1);
998:   PetscSectionSetFieldComponents(newCoordSection, 0, dim);
999:   PetscSectionSetChart(newCoordSection, vStart, vEnd);
1000:   for (v = vStart; v < vEnd; ++v) {
1001:     PetscSectionSetDof(newCoordSection, v, dim);
1002:     PetscSectionSetFieldDof(newCoordSection, v, 0, dim);
1003:   }
1004:   PetscSectionSetUp(newCoordSection);
1005:   DMSetCoordinateSection(sdm, PETSC_DETERMINE, newCoordSection);
1006:   PetscSectionDestroy(&newCoordSection); /* relinquish our reference */
1007:   DMGetCoordinatesLocal(dm, &coordinates);
1008:   DMSetCoordinatesLocal(sdm, coordinates);
1009:   /* Convert labels */
1010:   DMGetNumLabels(dm, &numLabels);
1011:   for (l = 0; l < numLabels; ++l) {
1012:     const char *lname;
1013:     PetscBool  isDepth;

1015:     DMGetLabelName(dm, l, &lname);
1016:     PetscStrcmp(lname, "depth", &isDepth);
1017:     if (isDepth) continue;
1018:     DMCreateLabel(sdm, lname);
1019:     DMGetLabelIdIS(dm, lname, &idIS);
1020:     ISGetLocalSize(idIS, &numFS);
1021:     ISGetIndices(idIS, &ids);
1022:     for (fs = 0; fs < numFS; ++fs) {
1023:       IS             pointIS;
1024:       const PetscInt *points;
1025:       PetscInt       numPoints;

1027:       DMGetStratumIS(dm, lname, ids[fs], &pointIS);
1028:       ISGetLocalSize(pointIS, &numPoints);
1029:       ISGetIndices(pointIS, &points);
1030:       for (p = 0; p < numPoints; ++p) {
1031:         PetscInt newpoint = points[p];

1033:         DMSetLabelValue(sdm, lname, newpoint, ids[fs]);
1034:       }
1035:       ISRestoreIndices(pointIS, &points);
1036:       ISDestroy(&pointIS);
1037:     }
1038:     ISRestoreIndices(idIS, &ids);
1039:     ISDestroy(&idIS);
1040:   }
1041:   /* Convert pointSF */
1042:   const PetscSFNode *remotePoints;
1043:   PetscSFNode       *gremotePoints;
1044:   const PetscInt    *localPoints;
1045:   PetscInt          *glocalPoints,*newLocation,*newRemoteLocation;
1046:   PetscInt          numRoots, numLeaves;
1047:   PetscMPIInt       size;

1049:   MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size);
1050:   DMGetPointSF(dm, &sfPoint);
1051:   DMGetPointSF(sdm, &gsfPoint);
1052:   DMPlexGetChart(dm,&pStart,&pEnd);
1053:   PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);
1054:   if (numRoots >= 0) {
1055:     PetscMalloc2(numRoots,&newLocation,pEnd-pStart,&newRemoteLocation);
1056:     for (l=0; l<numRoots; l++) newLocation[l] = l; /* + (l >= cEnd ? numGhostCells : 0); */
1057:     PetscSFBcastBegin(sfPoint, MPIU_INT, newLocation, newRemoteLocation);
1058:     PetscSFBcastEnd(sfPoint, MPIU_INT, newLocation, newRemoteLocation);
1059:     PetscMalloc1(numLeaves,    &glocalPoints);
1060:     PetscMalloc1(numLeaves, &gremotePoints);
1061:     for (l = 0; l < numLeaves; ++l) {
1062:       glocalPoints[l]        = localPoints[l]; /* localPoints[l] >= cEnd ? localPoints[l] + numGhostCells : localPoints[l]; */
1063:       gremotePoints[l].rank  = remotePoints[l].rank;
1064:       gremotePoints[l].index = newRemoteLocation[localPoints[l]];
1065:     }
1066:     PetscFree2(newLocation,newRemoteLocation);
1067:     PetscSFSetGraph(gsfPoint, numRoots+numGhostCells, numLeaves, glocalPoints, PETSC_OWN_POINTER, gremotePoints, PETSC_OWN_POINTER);
1068:   }
1069:   DMDestroy(dmSplit);
1070:   *dmSplit = sdm;
1071:   return(0);
1072: }

1074: PetscErrorCode CreatePartitionVec(DM dm, DM *dmCell, Vec *partition)
1075: {
1076:   PetscSF        sfPoint;
1077:   PetscSection   coordSection;
1078:   Vec            coordinates;
1079:   PetscSection   sectionCell;
1080:   PetscScalar    *part;
1081:   PetscInt       cStart, cEnd, c;
1082:   PetscMPIInt    rank;

1086:   DMGetCoordinateSection(dm, &coordSection);
1087:   DMGetCoordinatesLocal(dm, &coordinates);
1088:   DMClone(dm, dmCell);
1089:   DMGetPointSF(dm, &sfPoint);
1090:   DMSetPointSF(*dmCell, sfPoint);
1091:   DMSetCoordinateSection(*dmCell, PETSC_DETERMINE, coordSection);
1092:   DMSetCoordinatesLocal(*dmCell, coordinates);
1093:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
1094:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &sectionCell);
1095:   DMPlexGetHeightStratum(*dmCell, 0, &cStart, &cEnd);
1096:   PetscSectionSetChart(sectionCell, cStart, cEnd);
1097:   for (c = cStart; c < cEnd; ++c) {
1098:     PetscSectionSetDof(sectionCell, c, 1);
1099:   }
1100:   PetscSectionSetUp(sectionCell);
1101:   DMSetDefaultSection(*dmCell, sectionCell);
1102:   PetscSectionDestroy(&sectionCell);
1103:   DMCreateLocalVector(*dmCell, partition);
1104:   PetscObjectSetName((PetscObject)*partition, "partition");
1105:   VecGetArray(*partition, &part);
1106:   for (c = cStart; c < cEnd; ++c) {
1107:     PetscScalar *p;

1109:     DMPlexPointLocalRef(*dmCell, c, part, &p);
1110:     p[0] = rank;
1111:   }
1112:   VecRestoreArray(*partition, &part);
1113:   return(0);
1114: }

1116: PetscErrorCode CreateMassMatrix(DM dm, Vec *massMatrix, User user)
1117: {
1118:   DM                dmMass, dmFace, dmCell, dmCoord;
1119:   PetscSection      coordSection;
1120:   Vec               coordinates, facegeom, cellgeom;
1121:   PetscSection      sectionMass;
1122:   PetscScalar       *m;
1123:   const PetscScalar *fgeom, *cgeom, *coords;
1124:   PetscInt          vStart, vEnd, v;
1125:   PetscErrorCode    ierr;

1128:   DMGetCoordinateSection(dm, &coordSection);
1129:   DMGetCoordinatesLocal(dm, &coordinates);
1130:   DMClone(dm, &dmMass);
1131:   DMSetCoordinateSection(dmMass, PETSC_DETERMINE, coordSection);
1132:   DMSetCoordinatesLocal(dmMass, coordinates);
1133:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &sectionMass);
1134:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
1135:   PetscSectionSetChart(sectionMass, vStart, vEnd);
1136:   for (v = vStart; v < vEnd; ++v) {
1137:     PetscInt numFaces;

1139:     DMPlexGetSupportSize(dmMass, v, &numFaces);
1140:     PetscSectionSetDof(sectionMass, v, numFaces*numFaces);
1141:   }
1142:   PetscSectionSetUp(sectionMass);
1143:   DMSetDefaultSection(dmMass, sectionMass);
1144:   PetscSectionDestroy(&sectionMass);
1145:   DMGetLocalVector(dmMass, massMatrix);
1146:   VecGetArray(*massMatrix, &m);
1147:   DMPlexTSGetGeometryFVM(dm, &facegeom, &cellgeom, NULL);
1148:   VecGetDM(facegeom, &dmFace);
1149:   VecGetArrayRead(facegeom, &fgeom);
1150:   VecGetDM(cellgeom, &dmCell);
1151:   VecGetArrayRead(cellgeom, &cgeom);
1152:   DMGetCoordinateDM(dm, &dmCoord);
1153:   VecGetArrayRead(coordinates, &coords);
1154:   for (v = vStart; v < vEnd; ++v) {
1155:     const PetscInt        *faces;
1156:     const PetscFVFaceGeom *fgA, *fgB, *cg;
1157:     const PetscScalar     *vertex;
1158:     PetscInt               numFaces, sides[2], f, g;

1160:     DMPlexPointLocalRead(dmCoord, v, coords, &vertex);
1161:     DMPlexGetSupportSize(dmMass, v, &numFaces);
1162:     DMPlexGetSupport(dmMass, v, &faces);
1163:     for (f = 0; f < numFaces; ++f) {
1164:       sides[0] = faces[f];
1165:       DMPlexPointLocalRead(dmFace, faces[f], fgeom, &fgA);
1166:       for (g = 0; g < numFaces; ++g) {
1167:         const PetscInt *cells = NULL;;
1168:         PetscReal      area   = 0.0;
1169:         PetscInt       numCells;

1171:         sides[1] = faces[g];
1172:         DMPlexPointLocalRead(dmFace, faces[g], fgeom, &fgB);
1173:         DMPlexGetJoin(dmMass, 2, sides, &numCells, &cells);
1174:         if (numCells != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Invalid join for faces");
1175:         DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cg);
1176:         area += PetscAbsScalar((vertex[0] - cg->centroid[0])*(fgA->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1])*(fgA->centroid[0] - cg->centroid[0]));
1177:         area += PetscAbsScalar((vertex[0] - cg->centroid[0])*(fgB->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1])*(fgB->centroid[0] - cg->centroid[0]));
1178:         m[f*numFaces+g] = Dot2(fgA->normal, fgB->normal)*area*0.5;
1179:         DMPlexRestoreJoin(dmMass, 2, sides, &numCells, &cells);
1180:       }
1181:     }
1182:   }
1183:   VecRestoreArrayRead(facegeom, &fgeom);
1184:   VecRestoreArrayRead(cellgeom, &cgeom);
1185:   VecRestoreArrayRead(coordinates, &coords);
1186:   VecRestoreArray(*massMatrix, &m);
1187:   DMDestroy(&dmMass);
1188:   return(0);
1189: }

1191: /* Behavior will be different for multi-physics or when using non-default boundary conditions */
1192: static PetscErrorCode ModelSolutionSetDefault(Model mod,SolutionFunction func,void *ctx)
1193: {
1195:   mod->solution    = func;
1196:   mod->solutionctx = ctx;
1197:   return(0);
1198: }

1200: static PetscErrorCode ModelFunctionalRegister(Model mod,const char *name,PetscInt *offset,FunctionalFunction func,void *ctx)
1201: {
1203:   FunctionalLink link,*ptr;
1204:   PetscInt       lastoffset = -1;

1207:   for (ptr=&mod->functionalRegistry; *ptr; ptr = &(*ptr)->next) lastoffset = (*ptr)->offset;
1208:   PetscNew(&link);
1209:   PetscStrallocpy(name,&link->name);
1210:   link->offset = lastoffset + 1;
1211:   link->func   = func;
1212:   link->ctx    = ctx;
1213:   link->next   = NULL;
1214:   *ptr         = link;
1215:   *offset      = link->offset;
1216:   return(0);
1217: }

1219: static PetscErrorCode ModelFunctionalSetFromOptions(Model mod,PetscOptionItems *PetscOptionsObject)
1220: {
1222:   PetscInt       i,j;
1223:   FunctionalLink link;
1224:   char           *names[256];

1227:   mod->numMonitored = ALEN(names);
1228:   PetscOptionsStringArray("-monitor","list of functionals to monitor","",names,&mod->numMonitored,NULL);
1229:   /* Create list of functionals that will be computed somehow */
1230:   PetscMalloc1(mod->numMonitored,&mod->functionalMonitored);
1231:   /* Create index of calls that we will have to make to compute these functionals (over-allocation in general). */
1232:   PetscMalloc1(mod->numMonitored,&mod->functionalCall);
1233:   mod->numCall = 0;
1234:   for (i=0; i<mod->numMonitored; i++) {
1235:     for (link=mod->functionalRegistry; link; link=link->next) {
1236:       PetscBool match;
1237:       PetscStrcasecmp(names[i],link->name,&match);
1238:       if (match) break;
1239:     }
1240:     if (!link) SETERRQ1(mod->comm,PETSC_ERR_USER,"No known functional '%s'",names[i]);
1241:     mod->functionalMonitored[i] = link;
1242:     for (j=0; j<i; j++) {
1243:       if (mod->functionalCall[j]->func == link->func && mod->functionalCall[j]->ctx == link->ctx) goto next_name;
1244:     }
1245:     mod->functionalCall[mod->numCall++] = link; /* Just points to the first link using the result. There may be more results. */
1246: next_name:
1247:     PetscFree(names[i]);
1248:   }

1250:   /* Find out the maximum index of any functional computed by a function we will be calling (even if we are not using it) */
1251:   mod->maxComputed = -1;
1252:   for (link=mod->functionalRegistry; link; link=link->next) {
1253:     for (i=0; i<mod->numCall; i++) {
1254:       FunctionalLink call = mod->functionalCall[i];
1255:       if (link->func == call->func && link->ctx == call->ctx) {
1256:         mod->maxComputed = PetscMax(mod->maxComputed,link->offset);
1257:       }
1258:     }
1259:   }
1260:   return(0);
1261: }

1263: static PetscErrorCode FunctionalLinkDestroy(FunctionalLink *link)
1264: {
1266:   FunctionalLink l,next;

1269:   if (!link) return(0);
1270:   l     = *link;
1271:   *link = NULL;
1272:   for (; l; l=next) {
1273:     next = l->next;
1274:     PetscFree(l->name);
1275:     PetscFree(l);
1276:   }
1277:   return(0);
1278: }

1280: /* put the solution callback into a functional callback */
1281: static PetscErrorCode SolutionFunctional(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *modctx)
1282: {
1283:   Model          mod;
1286:   mod  = (Model) modctx;
1287:   (*mod->solution)(mod, time, x, u, mod->solutionctx);
1288:   return(0);
1289: }

1291: PetscErrorCode SetInitialCondition(DM dm, Vec X, User user)
1292: {
1293:   PetscErrorCode     (*func[1]) (PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
1294:   void               *ctx[1];
1295:   Model              mod = user->model;
1296:   PetscErrorCode     ierr;

1299:   func[0] = SolutionFunctional;
1300:   ctx[0]  = (void *) mod;
1301:   DMProjectFunction(dm,0.0,func,ctx,INSERT_ALL_VALUES,X);
1302:   return(0);
1303: }

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

1310:   PetscViewerCreate(PetscObjectComm((PetscObject)dm), viewer);
1311:   PetscViewerSetType(*viewer, PETSCVIEWERVTK);
1312:   PetscViewerFileSetName(*viewer, filename);
1313:   return(0);
1314: }

1316: static PetscErrorCode MonitorVTK(TS ts,PetscInt stepnum,PetscReal time,Vec X,void *ctx)
1317: {
1318:   User           user = (User)ctx;
1319:   DM             dm;
1320:   Vec            cellgeom;
1321:   PetscViewer    viewer;
1322:   char           filename[PETSC_MAX_PATH_LEN],*ftable = NULL;
1323:   PetscReal      xnorm;
1324:   PetscInt       cEndInterior;

1328:   PetscObjectSetName((PetscObject) X, "u");
1329:   VecGetDM(X,&dm);
1330:   DMPlexTSGetGeometryFVM(dm, NULL, &cellgeom, NULL);
1331:   VecNorm(X,NORM_INFINITY,&xnorm);

1333:   if (stepnum >= 0) {
1334:     stepnum += user->monitorStepOffset;
1335:   }
1336:   if (stepnum >= 0) {           /* No summary for final time */
1337:     Model             mod = user->model;
1338:     PetscInt          c,cStart,cEnd,fcount,i;
1339:     size_t            ftableused,ftablealloc;
1340:     const PetscScalar *cgeom,*x;
1341:     DM                dmCell;
1342:     DMLabel           vtkLabel;
1343:     PetscReal         *fmin,*fmax,*fintegral,*ftmp;
1344:     fcount = mod->maxComputed+1;
1345:     PetscMalloc4(fcount,&fmin,fcount,&fmax,fcount,&fintegral,fcount,&ftmp);
1346:     for (i=0; i<fcount; i++) {
1347:       fmin[i]      = PETSC_MAX_REAL;
1348:       fmax[i]      = PETSC_MIN_REAL;
1349:       fintegral[i] = 0;
1350:     }
1351:     VecGetDM(cellgeom,&dmCell);
1352:     DMPlexGetHybridBounds(dmCell, &cEndInterior, NULL, NULL, NULL);
1353:     DMPlexGetHeightStratum(dmCell,0,&cStart,&cEnd);
1354:     VecGetArrayRead(cellgeom,&cgeom);
1355:     VecGetArrayRead(X,&x);
1356:     DMGetLabel(dm,"vtk",&vtkLabel);
1357:     for (c = cStart; c < cEndInterior; ++c) {
1358:       const PetscFVCellGeom *cg;
1359:       const PetscScalar     *cx;
1360:       PetscInt              vtkVal = 0;

1362:       /* not that these two routines as currently implemented work for any dm with a
1363:        * defaultSection/defaultGlobalSection */
1364:       DMPlexPointLocalRead(dmCell,c,cgeom,&cg);
1365:       DMPlexPointGlobalRead(dm,c,x,&cx);
1366:       if (vtkLabel) {DMLabelGetValue(vtkLabel,c,&vtkVal);}
1367:       if (!vtkVal || !cx) continue;        /* ghost, or not a global cell */
1368:       for (i=0; i<mod->numCall; i++) {
1369:         FunctionalLink flink = mod->functionalCall[i];
1370:         (*flink->func)(mod,time,cg->centroid,cx,ftmp,flink->ctx);
1371:       }
1372:       for (i=0; i<fcount; i++) {
1373:         fmin[i]       = PetscMin(fmin[i],ftmp[i]);
1374:         fmax[i]       = PetscMax(fmax[i],ftmp[i]);
1375:         fintegral[i] += cg->volume * ftmp[i];
1376:       }
1377:     }
1378:     VecRestoreArrayRead(cellgeom,&cgeom);
1379:     VecRestoreArrayRead(X,&x);
1380:     MPI_Allreduce(MPI_IN_PLACE,fmin,fcount,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));
1381:     MPI_Allreduce(MPI_IN_PLACE,fmax,fcount,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));
1382:     MPI_Allreduce(MPI_IN_PLACE,fintegral,fcount,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)ts));

1384:     ftablealloc = fcount * 100;
1385:     ftableused  = 0;
1386:     PetscMalloc1(ftablealloc,&ftable);
1387:     for (i=0; i<mod->numMonitored; i++) {
1388:       size_t         countused;
1389:       char           buffer[256],*p;
1390:       FunctionalLink flink = mod->functionalMonitored[i];
1391:       PetscInt       id    = flink->offset;
1392:       if (i % 3) {
1393:         PetscMemcpy(buffer,"  ",2);
1394:         p    = buffer + 2;
1395:       } else if (i) {
1396:         char newline[] = "\n";
1397:         PetscMemcpy(buffer,newline,sizeof newline-1);
1398:         p    = buffer + sizeof newline - 1;
1399:       } else {
1400:         p = buffer;
1401:       }
1402:       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]);
1403:       countused += p - buffer;
1404:       if (countused > ftablealloc-ftableused-1) { /* reallocate */
1405:         char *ftablenew;
1406:         ftablealloc = 2*ftablealloc + countused;
1407:         PetscMalloc(ftablealloc,&ftablenew);
1408:         PetscMemcpy(ftablenew,ftable,ftableused);
1409:         PetscFree(ftable);
1410:         ftable = ftablenew;
1411:       }
1412:       PetscMemcpy(ftable+ftableused,buffer,countused);
1413:       ftableused += countused;
1414:       ftable[ftableused] = 0;
1415:     }
1416:     PetscFree4(fmin,fmax,fintegral,ftmp);

1418:     PetscPrintf(PetscObjectComm((PetscObject)ts),"% 3D  time %8.4g  |x| %8.4g  %s\n",stepnum,(double)time,(double)xnorm,ftable ? ftable : "");
1419:     PetscFree(ftable);
1420:   }
1421:   if (user->vtkInterval < 1) return(0);
1422:   if ((stepnum == -1) ^ (stepnum % user->vtkInterval == 0)) {
1423:     if (stepnum == -1) {        /* Final time is not multiple of normal time interval, write it anyway */
1424:       TSGetTimeStepNumber(ts,&stepnum);
1425:     }
1426:     PetscSNPrintf(filename,sizeof filename,"ex11-%03D.vtu",stepnum);
1427:     OutputVTK(dm,filename,&viewer);
1428:     VecView(X,viewer);
1429:     PetscViewerDestroy(&viewer);
1430:   }
1431:   return(0);
1432: }

1434: static PetscErrorCode initializeTS(DM dm, User user, TS *ts)
1435: {

1439:   TSCreate(PetscObjectComm((PetscObject)dm), ts);
1440:   TSSetType(*ts, TSSSP);
1441:   TSSetDM(*ts, dm);
1442:   TSMonitorSet(*ts,MonitorVTK,user,NULL);
1443:   DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, user);
1444:   DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, user);
1445:   TSSetDuration(*ts,1000,2.0);
1446:   TSSetExactFinalTime(*ts,TS_EXACTFINALTIME_STEPOVER);
1447:   return(0);
1448: }

1450: static PetscErrorCode adaptToleranceFVM(PetscFV fvm, TS ts, Vec sol, VecTagger refineTag, VecTagger coarsenTag, User user, TS *tsNew, Vec *solNew)
1451: {
1452:   DM                dm, gradDM, plex, cellDM, adaptedDM = NULL;
1453:   Vec               cellGeom, faceGeom;
1454:   PetscBool         isForest, computeGradient;
1455:   Vec               grad, locGrad, locX, errVec;
1456:   PetscInt          cStart, cEnd, cEndInterior, c, dim, nRefine, nCoarsen;
1457:   PetscReal         minMaxInd[2] = {PETSC_MAX_REAL, PETSC_MIN_REAL}, minMaxIndGlobal[2], minInd, maxInd, time;
1458:   PetscScalar       *errArray;
1459:   const PetscScalar *pointVals;
1460:   const PetscScalar *pointGrads;
1461:   const PetscScalar *pointGeom;
1462:   DMLabel           adaptLabel = NULL;
1463:   IS                refineIS, coarsenIS;
1464:   PetscErrorCode    ierr;

1467:   TSGetTime(ts,&time);
1468:   VecGetDM(sol, &dm);
1469:   DMGetDimension(dm,&dim);
1470:   PetscFVGetComputeGradients(fvm,&computeGradient);
1471:   PetscFVSetComputeGradients(fvm,PETSC_TRUE);
1472:   DMIsForest(dm, &isForest);
1473:   DMConvert(dm, DMPLEX, &plex);
1474:   DMPlexGetDataFVM(plex, fvm, &cellGeom, &faceGeom, &gradDM);
1475:   DMCreateLocalVector(plex,&locX);
1476:   DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, 0.0, faceGeom, cellGeom, NULL);
1477:   DMGlobalToLocalBegin(plex, sol, INSERT_VALUES, locX);
1478:   DMGlobalToLocalEnd  (plex, sol, INSERT_VALUES, locX);
1479:   DMCreateGlobalVector(gradDM, &grad);
1480:   DMPlexReconstructGradientsFVM(plex, locX, grad);
1481:   DMCreateLocalVector(gradDM, &locGrad);
1482:   DMGlobalToLocalBegin(gradDM, grad, INSERT_VALUES, locGrad);
1483:   DMGlobalToLocalEnd(gradDM, grad, INSERT_VALUES, locGrad);
1484:   VecDestroy(&grad);
1485:   DMPlexGetHeightStratum(plex,0,&cStart,&cEnd);
1486:   DMPlexGetHybridBounds(plex,&cEndInterior,NULL,NULL,NULL);
1487:   cEnd = (cEndInterior < 0) ? cEnd : cEndInterior;
1488:   VecGetArrayRead(locGrad,&pointGrads);
1489:   VecGetArrayRead(cellGeom,&pointGeom);
1490:   VecGetArrayRead(locX,&pointVals);
1491:   VecGetDM(cellGeom,&cellDM);
1492:   DMLabelCreate("adapt",&adaptLabel);
1493:   VecCreateMPI(PetscObjectComm((PetscObject)plex),cEnd-cStart,PETSC_DETERMINE,&errVec);
1494:   VecSetUp(errVec);
1495:   VecGetArray(errVec,&errArray);
1496:   for (c = cStart; c < cEnd; c++) {
1497:     PetscReal             errInd = 0.;
1498:     const PetscScalar     *pointGrad;
1499:     const PetscScalar     *pointVal;
1500:     const PetscFVCellGeom *cg;

1502:     DMPlexPointLocalRead(gradDM,c,pointGrads,&pointGrad);
1503:     DMPlexPointLocalRead(cellDM,c,pointGeom,&cg);
1504:     DMPlexPointLocalRead(plex,c,pointVals,&pointVal);

1506:     (user->model->errorIndicator)(dim,cg->volume,user->model->physics->dof,pointVal,pointGrad,&errInd,user->model->errorCtx);
1507:     errArray[c-cStart] = errInd;
1508:     minMaxInd[0] = PetscMin(minMaxInd[0],errInd);
1509:     minMaxInd[1] = PetscMax(minMaxInd[1],errInd);
1510:   }
1511:   VecRestoreArray(errVec,&errArray);
1512:   VecRestoreArrayRead(locX,&pointVals);
1513:   VecRestoreArrayRead(cellGeom,&pointGeom);
1514:   VecRestoreArrayRead(locGrad,&pointGrads);
1515:   VecDestroy(&locGrad);
1516:   VecDestroy(&locX);
1517:   DMDestroy(&plex);

1519:   VecTaggerComputeIS(refineTag,errVec,&refineIS);
1520:   VecTaggerComputeIS(coarsenTag,errVec,&coarsenIS);
1521:   ISGetSize(refineIS,&nRefine);
1522:   ISGetSize(coarsenIS,&nCoarsen);
1523:   if (nRefine) {DMLabelSetStratumIS(adaptLabel,DM_ADAPT_REFINE,refineIS);}
1524:   if (nCoarsen) {DMLabelSetStratumIS(adaptLabel,DM_ADAPT_COARSEN,coarsenIS);}
1525:   ISDestroy(&coarsenIS);
1526:   ISDestroy(&refineIS);
1527:   VecDestroy(&errVec);

1529:   PetscFVSetComputeGradients(fvm,computeGradient);
1530:   minMaxInd[1] = -minMaxInd[1];
1531:   MPI_Allreduce(minMaxInd,minMaxIndGlobal,2,MPIU_REAL,MPI_MIN,PetscObjectComm((PetscObject)dm));
1532:   minInd = minMaxIndGlobal[0];
1533:   maxInd = -minMaxIndGlobal[1];
1534:   PetscInfo2(ts, "error indicator range (%E, %E)\n", minInd, maxInd);
1535:   if (nRefine || nCoarsen) { /* at least one cell is over the refinement threshold */
1536:     DMAdaptLabel(dm,adaptLabel,&adaptedDM);
1537:   }
1538:   DMLabelDestroy(&adaptLabel);
1539:   if (adaptedDM) {
1540:     if (tsNew) {
1541:       initializeTS(adaptedDM,user,tsNew);
1542:     }
1543:     if (solNew) {
1544:       DMCreateGlobalVector(adaptedDM, solNew);
1545:       PetscObjectSetName((PetscObject) *solNew, "solution");
1546:       DMForestTransferVec(dm, sol, adaptedDM, *solNew, PETSC_TRUE, time);
1547:     }
1548:     if (isForest) {DMForestSetAdaptivityForest(adaptedDM,NULL);} /* clear internal references to the previous dm */
1549:     DMDestroy(&adaptedDM);
1550:   }
1551:   else {
1552:     if (tsNew) *tsNew = NULL;
1553:     if (solNew) *solNew = NULL;
1554:   }
1555:   return(0);
1556: }

1558: int main(int argc, char **argv)
1559: {
1560:   MPI_Comm          comm;
1561:   PetscDS           prob;
1562:   PetscFV           fvm;
1563:   PetscLimiter      limiter = NULL, noneLimiter = NULL;
1564:   User              user;
1565:   Model             mod;
1566:   Physics           phys;
1567:   DM                dm;
1568:   PetscReal         ftime, cfl, dt, minRadius;
1569:   PetscInt          dim, nsteps;
1570:   TS                ts;
1571:   TSConvergedReason reason;
1572:   Vec               X;
1573:   PetscViewer       viewer;
1574:   PetscBool         vtkCellGeom, splitFaces, useAMR, viewInitial;
1575:   PetscInt          overlap, adaptInterval;
1576:   char              filename[PETSC_MAX_PATH_LEN] = "sevenside.exo";
1577:   char              physname[256]  = "advect";
1578:   VecTagger         refineTag = NULL, coarsenTag = NULL;
1579:   PetscErrorCode    ierr;

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

1584:   PetscNew(&user);
1585:   PetscNew(&user->model);
1586:   PetscNew(&user->model->physics);
1587:   mod           = user->model;
1588:   phys          = mod->physics;
1589:   mod->comm     = comm;
1590:   useAMR        = PETSC_FALSE;
1591:   viewInitial   = PETSC_FALSE;
1592:   adaptInterval = 1;

1594:   /* Register physical models to be available on the command line */
1595:   PetscFunctionListAdd(&PhysicsList,"advect"          ,PhysicsCreate_Advect);
1596:   PetscFunctionListAdd(&PhysicsList,"sw"              ,PhysicsCreate_SW);
1597:   PetscFunctionListAdd(&PhysicsList,"euler"           ,PhysicsCreate_Euler);

1599:   PetscOptionsBegin(comm,NULL,"Unstructured Finite Volume Mesh Options","");
1600:   {
1601:     cfl  = 0.9 * 4; /* default SSPRKS2 with s=5 stages is stable for CFL number s-1 */
1602:     PetscOptionsReal("-ufv_cfl","CFL number per step","",cfl,&cfl,NULL);
1603:     PetscOptionsString("-f","Exodus.II filename to read","",filename,filename,sizeof(filename),NULL);
1604:     splitFaces = PETSC_FALSE;
1605:     PetscOptionsBool("-ufv_split_faces","Split faces between cell sets","",splitFaces,&splitFaces,NULL);
1606:     overlap = 1;
1607:     PetscOptionsInt("-ufv_mesh_overlap","Number of cells to overlap partitions","",overlap,&overlap,NULL);
1608:     user->vtkInterval = 1;
1609:     PetscOptionsInt("-ufv_vtk_interval","VTK output interval (0 to disable)","",user->vtkInterval,&user->vtkInterval,NULL);
1610:     vtkCellGeom = PETSC_FALSE;
1611:     PetscOptionsBool("-ufv_vtk_cellgeom","Write cell geometry (for debugging)","",vtkCellGeom,&vtkCellGeom,NULL);
1612:     PetscOptionsBool("-ufv_use_amr","use local adaptive mesh refinement","",useAMR,&useAMR,NULL);
1613:     PetscOptionsBool("-ufv_view_initial_refinement","View initial conditions refinement history","",viewInitial,&viewInitial,NULL);
1614:     PetscOptionsInt("-ufv_adapt_interval","time steps between AMR","",adaptInterval,&adaptInterval,NULL);
1615:   }
1616:   PetscOptionsEnd();

1618:   if (useAMR) {
1619:     VecTaggerBox refineBox, coarsenBox;

1621:     refineBox.min  = refineBox.max  = PETSC_MAX_REAL;
1622:     coarsenBox.min = coarsenBox.max = PETSC_MIN_REAL;

1624:     VecTaggerCreate(comm,&refineTag);
1625:     PetscObjectSetOptionsPrefix((PetscObject)refineTag,"refine_");
1626:     VecTaggerSetType(refineTag,VECTAGGERABSOLUTE);
1627:     VecTaggerAbsoluteSetBox(refineTag,&refineBox);
1628:     VecTaggerSetFromOptions(refineTag);
1629:     VecTaggerSetUp(refineTag);
1630:     PetscObjectViewFromOptions((PetscObject)refineTag,NULL,"-tag_view");

1632:     VecTaggerCreate(comm,&coarsenTag);
1633:     PetscObjectSetOptionsPrefix((PetscObject)coarsenTag,"coarsen_");
1634:     VecTaggerSetType(coarsenTag,VECTAGGERABSOLUTE);
1635:     VecTaggerAbsoluteSetBox(coarsenTag,&coarsenBox);
1636:     VecTaggerSetFromOptions(coarsenTag);
1637:     VecTaggerSetUp(coarsenTag);
1638:     PetscObjectViewFromOptions((PetscObject)coarsenTag,NULL,"-tag_view");
1639:   }

1641:   PetscOptionsBegin(comm,NULL,"Unstructured Finite Volume Physics Options","");
1642:   {
1643:     PetscErrorCode (*physcreate)(Model,Physics,PetscOptionItems*);
1644:     PetscOptionsFList("-physics","Physics module to solve","",PhysicsList,physname,physname,sizeof physname,NULL);
1645:     PetscFunctionListFind(PhysicsList,physname,&physcreate);
1646:     PetscMemzero(phys,sizeof(struct _n_Physics));
1647:     (*physcreate)(mod,phys,PetscOptionsObject);
1648:     /* Count number of fields and dofs */
1649:     for (phys->nfields=0,phys->dof=0; phys->field_desc[phys->nfields].name; phys->nfields++) phys->dof += phys->field_desc[phys->nfields].dof;
1650:     if (phys->dof <= 0) SETERRQ1(comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set dof",physname);
1651:     ModelFunctionalSetFromOptions(mod,PetscOptionsObject);
1652:   }
1653:   PetscOptionsEnd();

1655:   {
1656:     size_t len,i;
1657:     for (i = 0; i < DIM; i++) { mod->bounds[2*i] = 0.; mod->bounds[2*i+1] = 1.;};
1658:     PetscStrlen(filename,&len);
1659:     dim = DIM;
1660:     if (!len) { /* a null name means just do a hex box */
1661:       PetscInt cells[3] = {1, 1, 1}; /* coarse mesh is one cell; refine from there */
1662:       PetscBool flg1, flg2, skew = PETSC_FALSE;
1663:       PetscInt nret1 = DIM;
1664:       PetscInt nret2 = 2*DIM;
1665:       PetscOptionsBegin(comm,NULL,"Rectangular mesh options","");
1666:       PetscOptionsIntArray("-grid_size","number of cells in each direction","",cells,&nret1,&flg1);
1667:       PetscOptionsRealArray("-grid_bounds","bounds of the mesh in each direction (i.e., x_min,x_max,y_min,y_max","",mod->bounds,&nret2,&flg2);
1668:       PetscOptionsBool("-grid_skew_60","Skew grid for 60 degree shock mesh","",skew,&skew,NULL);
1669:       PetscOptionsEnd();
1670:       if (flg1) {
1671:         dim = nret1;
1672:         if (dim != DIM) SETERRQ1(comm,PETSC_ERR_ARG_SIZ,"Dim wrong size %D in -grid_size",dim);
1673:       }
1674:       DMPlexCreateHexBoxMesh(comm, dim, cells, mod->bcs[0], mod->bcs[1], mod->bcs[2], &dm);
1675:       if (flg2) {
1676:         PetscInt dimEmbed, i;
1677:         PetscInt nCoords;
1678:         PetscScalar *coords;
1679:         Vec coordinates;

1681:         DMGetCoordinatesLocal(dm,&coordinates);
1682:         DMGetCoordinateDim(dm,&dimEmbed);
1683:         VecGetLocalSize(coordinates,&nCoords);
1684:         if (nCoords % dimEmbed) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Coordinate vector the wrong size");
1685:         VecGetArray(coordinates,&coords);
1686:         for (i = 0; i < nCoords; i += dimEmbed) {
1687:           PetscInt j;

1689:           PetscScalar *coord = &coords[i];
1690:           for (j = 0; j < dimEmbed; j++) {
1691:             coord[j] = mod->bounds[2 * j] + coord[j] * (mod->bounds[2 * j + 1] - mod->bounds[2 * j]);
1692:             if (dim==2 && cells[1]==1 && j==0 && skew) {
1693:               if (cells[0]==2 && i==8) {
1694:                 coord[j] = .57735026918963; /* hack to get 60 deg skewed mesh */
1695:               }
1696:               else if (cells[0]==3) {
1697:                 if(i==2 || i==10) coord[j] = mod->bounds[1]/4.;
1698:                 else if (i==4) coord[j] = mod->bounds[1]/2.;
1699:                 else if (i==12) coord[j] = 1.57735026918963*mod->bounds[1]/2.;
1700:               }
1701:             }
1702:           }
1703:         }
1704:         VecRestoreArray(coordinates,&coords);
1705:         DMSetCoordinatesLocal(dm,coordinates);
1706:       }
1707:     }
1708:     else {
1709:       DMPlexCreateFromFile(comm, filename, PETSC_TRUE, &dm);
1710:     }
1711:   }
1712:   DMViewFromOptions(dm, NULL, "-dm_view");
1713:   DMGetDimension(dm, &dim);

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

1718:   mod->errorIndicator = ErrorIndicator_Simple;

1720:   {
1721:     DM dmDist;

1723:     DMPlexSetAdjacencyUseCone(dm, PETSC_TRUE);
1724:     DMPlexSetAdjacencyUseClosure(dm, PETSC_FALSE);
1725:     DMPlexDistribute(dm, overlap, NULL, &dmDist);
1726:     if (dmDist) {
1727:       DMDestroy(&dm);
1728:       dm   = dmDist;
1729:     }
1730:   }

1732:   DMSetFromOptions(dm);

1734:   {
1735:     DM gdm;

1737:     DMPlexConstructGhostCells(dm, NULL, NULL, &gdm);
1738:     DMDestroy(&dm);
1739:     dm   = gdm;
1740:     DMViewFromOptions(dm, NULL, "-dm_view");
1741:   }
1742:   if (splitFaces) {ConstructCellBoundary(dm, user);}
1743:   SplitFaces(&dm, "split faces", user);

1745:   PetscDSCreate(PetscObjectComm((PetscObject)dm),&prob);
1746:   PetscFVCreate(comm, &fvm);
1747:   PetscFVSetFromOptions(fvm);
1748:   PetscFVSetNumComponents(fvm, phys->dof);
1749:   PetscFVSetSpatialDimension(fvm, dim);
1750:   PetscObjectSetName((PetscObject) fvm,"");
1751:   {
1752:     PetscInt f, dof;
1753:     for (f=0,dof=0; f < phys->nfields; f++) {
1754:       PetscInt newDof = phys->field_desc[f].dof;

1756:       if (newDof == 1) {
1757:         PetscFVSetComponentName(fvm,dof,phys->field_desc[f].name);
1758:       }
1759:       else {
1760:         PetscInt j;

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

1765:           PetscSNPrintf(compName,sizeof(compName),"%s_%d",phys->field_desc[f].name,j);
1766:           PetscFVSetComponentName(fvm,dof+j,compName);
1767:         }
1768:       }
1769:       dof += newDof;
1770:     }
1771:   }
1772:   /* FV is now structured with one field having all physics as components */
1773:   PetscDSAddDiscretization(prob, (PetscObject) fvm);
1774:   PetscDSSetRiemannSolver(prob, 0, user->model->physics->riemann);
1775:   PetscDSSetContext(prob, 0, user->model->physics);
1776:   (*mod->setupbc)(prob,phys);
1777:   PetscDSSetFromOptions(prob);
1778:   DMSetDS(dm,prob);
1779:   PetscDSDestroy(&prob);
1780:   {
1781:     char      convType[256];
1782:     PetscBool flg;

1784:     PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");
1785:     PetscOptionsFList("-dm_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg);
1786:     PetscOptionsEnd();
1787:     if (flg) {
1788:       DM dmConv;

1790:       DMConvert(dm,convType,&dmConv);
1791:       if (dmConv) {
1792:         DMViewFromOptions(dmConv, NULL, "-dm_conv_view");
1793:         DMDestroy(&dm);
1794:         dm   = dmConv;
1795:         DMSetFromOptions(dm);
1796:       }
1797:     }
1798:   }

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

1802:   DMCreateGlobalVector(dm, &X);
1803:   PetscObjectSetName((PetscObject) X, "solution");
1804:   SetInitialCondition(dm, X, user);
1805:   if (useAMR) {

1807:     /* use no limiting when reconstructing gradients for adaptivity */
1808:     PetscFVGetLimiter(fvm, &limiter);
1809:     PetscObjectReference((PetscObject)limiter);

1811:     PetscLimiterCreate(PetscObjectComm((PetscObject)fvm),&noneLimiter);
1812:     PetscLimiterSetType(noneLimiter,PETSCLIMITERNONE);

1814:     PetscFVSetLimiter(fvm,noneLimiter);
1815:     {
1816:       PetscInt adaptIter;

1818:       for (adaptIter = 0;;adaptIter++) {
1819:         PetscLogDouble bytes;
1820:         TS             tsNew = NULL;

1822:         if (viewInitial) {
1823:           PetscViewer viewer;
1824:           char        buf[256];
1825:           PetscBool   isHDF5, isVTK;

1827:           PetscViewerCreate(comm,&viewer);
1828:           PetscViewerSetType(viewer,PETSCVIEWERVTK);
1829:           PetscViewerSetOptionsPrefix(viewer,"initial_");
1830:           PetscViewerSetFromOptions(viewer);
1831:           PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&isHDF5);
1832:           PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isVTK);
1833:           if (isHDF5) {
1834:             PetscSNPrintf(buf, 256, "ex11-initial-%d.h5", adaptIter);
1835:           } else if (isVTK) {
1836:             PetscSNPrintf(buf, 256, "ex11-initial-%d.vtu", adaptIter);
1837:             PetscViewerPushFormat(viewer,PETSC_VIEWER_VTK_VTU);
1838:           }
1839:           PetscViewerFileSetMode(viewer,FILE_MODE_WRITE);
1840:           PetscViewerFileSetName(viewer,buf);
1841:           if (isHDF5) {
1842:             DMView(dm,viewer);
1843:             PetscViewerFileSetMode(viewer,FILE_MODE_UPDATE);
1844:           }
1845:           VecView(X,viewer);
1846:           PetscViewerDestroy(&viewer);
1847:         }

1849:         PetscMemoryGetCurrentUsage(&bytes);
1850:         PetscInfo2(ts, "refinement loop %D: memory used %g\n", adaptIter, bytes);
1851:         adaptToleranceFVM(fvm, ts, X, refineTag, coarsenTag, user, &tsNew, NULL);
1852:         if (!tsNew) {
1853:           break;
1854:         } else {
1855:           DMDestroy(&dm);
1856:           VecDestroy(&X);
1857:           TSDestroy(&ts);
1858:           ts   = tsNew;
1859:           TSGetDM(ts,&dm);
1860:           PetscObjectReference((PetscObject)dm);
1861:           DMCreateGlobalVector(dm,&X);
1862:           PetscObjectSetName((PetscObject) X, "solution");
1863:           SetInitialCondition(dm, X, user);
1864:         }
1865:       }
1866:     }
1867:     /* restore original limiter */
1868:     PetscFVSetLimiter(fvm,limiter);
1869:   }

1871:   if (vtkCellGeom) {
1872:     DM  dmCell;
1873:     Vec cellgeom, partition;

1875:     DMPlexTSGetGeometryFVM(dm, NULL, &cellgeom, NULL);
1876:     OutputVTK(dm, "ex11-cellgeom.vtk", &viewer);
1877:     VecView(cellgeom, viewer);
1878:     PetscViewerDestroy(&viewer);
1879:     CreatePartitionVec(dm, &dmCell, &partition);
1880:     OutputVTK(dmCell, "ex11-partition.vtk", &viewer);
1881:     VecView(partition, viewer);
1882:     PetscViewerDestroy(&viewer);
1883:     VecDestroy(&partition);
1884:     DMDestroy(&dmCell);
1885:   }

1887:   /* collect max maxspeed from all processes -- todo */
1888:   DMPlexTSGetGeometryFVM(dm, NULL, NULL, &minRadius);
1889:   MPI_Allreduce(&phys->maxspeed,&mod->maxspeed,1,MPIU_REAL,MPI_MAX,PetscObjectComm((PetscObject)ts));
1890:   if (mod->maxspeed <= 0) SETERRQ1(comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set maxspeed",physname);
1891:   dt   = cfl * minRadius / mod->maxspeed;
1892:   TSSetInitialTimeStep(ts,0.0,dt);
1893:   TSSetFromOptions(ts);
1894:   if (!useAMR) {
1895:     TSSolve(ts,X);
1896:     TSGetSolveTime(ts,&ftime);
1897:     TSGetTimeStepNumber(ts,&nsteps);
1898:   } else {
1899:     PetscReal finalTime;
1900:     PetscInt  adaptIter;
1901:     TS        tsNew = NULL;
1902:     Vec       solNew = NULL;
1903:     PetscInt  incSteps;

1905:     TSGetDuration(ts,NULL,&finalTime);
1906:     TSSetDuration(ts,adaptInterval,finalTime);
1907:     TSSolve(ts,X);
1908:     TSGetSolveTime(ts,&ftime);
1909:     TSGetTimeStepNumber(ts,&nsteps);
1910:     for (adaptIter = 0;ftime < finalTime;adaptIter++) {
1911:       PetscLogDouble bytes;

1913:       PetscMemoryGetCurrentUsage(&bytes);
1914:       PetscInfo2(ts, "AMR time step loop %D: memory used %g\n", adaptIter, bytes);
1915:       PetscFVSetLimiter(fvm,noneLimiter);
1916:       adaptToleranceFVM(fvm,ts,X,refineTag,coarsenTag,user,&tsNew,&solNew);
1917:       PetscFVSetLimiter(fvm,limiter);
1918:       if (tsNew) {
1919:         PetscInfo(ts, "AMR used\n");
1920:         DMDestroy(&dm);
1921:         VecDestroy(&X);
1922:         TSDestroy(&ts);
1923:         ts   = tsNew;
1924:         X    = solNew;
1925:         TSSetFromOptions(ts);
1926:         VecGetDM(X,&dm);
1927:         PetscObjectReference((PetscObject)dm);
1928:         DMPlexTSGetGeometryFVM(dm, NULL, NULL, &minRadius);
1929:         MPI_Allreduce(&phys->maxspeed,&mod->maxspeed,1,MPIU_REAL,MPI_MAX,PetscObjectComm((PetscObject)ts));
1930:         if (mod->maxspeed <= 0) SETERRQ1(comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set maxspeed",physname);
1931:         dt   = cfl * minRadius / mod->maxspeed;
1932:         TSSetInitialTimeStep(ts,ftime,dt);
1933:       }
1934:       else {
1935:         PetscInfo(ts, "AMR not used\n");
1936:       }
1937:       user->monitorStepOffset = nsteps;
1938:       TSSetDuration(ts,adaptInterval,finalTime);
1939:       TSSolve(ts,X);
1940:       TSGetSolveTime(ts,&ftime);
1941:       TSGetTimeStepNumber(ts,&incSteps);
1942:       nsteps += incSteps;
1943:     }
1944:   }
1945:   TSGetConvergedReason(ts,&reason);
1946:   PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,nsteps);
1947:   TSDestroy(&ts);

1949:   VecTaggerDestroy(&refineTag);
1950:   VecTaggerDestroy(&coarsenTag);
1951:   PetscFunctionListDestroy(&PhysicsList);
1952:   FunctionalLinkDestroy(&user->model->functionalRegistry);
1953:   PetscFree(user->model->functionalMonitored);
1954:   PetscFree(user->model->functionalCall);
1955:   PetscFree(user->model->physics->data);
1956:   PetscFree(user->model->physics);
1957:   PetscFree(user->model);
1958:   PetscFree(user);
1959:   VecDestroy(&X);
1960:   PetscLimiterDestroy(&limiter);
1961:   PetscLimiterDestroy(&noneLimiter);
1962:   PetscFVDestroy(&fvm);
1963:   DMDestroy(&dm);
1964:   PetscFinalize();
1965:   return ierr;
1966: }

1968: /* Godunov fluxs */
1969: PetscScalar cvmgp_(PetscScalar *a, PetscScalar *b, PetscScalar *test)
1970: {
1971:     /* System generated locals */
1972:     PetscScalar ret_val;

1974:     if (*test > 0.) {
1975:         goto L10;
1976:     }
1977:     ret_val = *b;
1978:     return ret_val;
1979: L10:
1980:     ret_val = *a;
1981:     return ret_val;
1982: } /* cvmgp_ */

1984: PetscScalar cvmgm_(PetscScalar *a, PetscScalar *b, PetscScalar *test)
1985: {
1986:     /* System generated locals */
1987:     PetscScalar ret_val;

1989:     if (*test < 0.) {
1990:         goto L10;
1991:     }
1992:     ret_val = *b;
1993:     return ret_val;
1994: L10:
1995:     ret_val = *a;
1996:     return ret_val;
1997: } /* cvmgm_ */

1999: int riem1mdt( PetscScalar *gaml, PetscScalar *gamr, PetscScalar *rl, PetscScalar *pl,
2000:               PetscScalar *uxl, PetscScalar *rr, PetscScalar *pr,
2001:               PetscScalar *uxr, PetscScalar *rstarl, PetscScalar *rstarr, PetscScalar *
2002:               pstar, PetscScalar *ustar)
2003: {
2004:     /* Initialized data */

2006:     static PetscScalar smallp = 1e-8;

2008:     /* System generated locals */
2009:     int i__1;
2010:     PetscScalar d__1, d__2;

2012:     /* Local variables */
2013:     static int i0;
2014:     static PetscScalar cl, cr, wl, zl, wr, zr, pst, durl, skpr1, skpr2;
2015:     static int iwave;
2016:     static PetscScalar csqrl, csqrr, gascl1, gascl2, gascl3, gascl4, gascr1, gascr2, gascr3, gascr4, cstarl, dpstar, cstarr;
2017:     static int iterno;
2018:     static PetscScalar ustarl, ustarr, rarepr1, rarepr2;



2022:     gascl1 = *gaml - 1.;
2023:     gascl2 = (*gaml + 1.) * .5;
2024:     gascl3 = gascl2 / *gaml;
2025:     gascl4 = 1.f / (*gaml - 1.);

2027:     gascr1 = *gamr - 1.;
2028:     gascr2 = (*gamr + 1.) * .5;
2029:     gascr3 = gascr2 / *gamr;
2030:     gascr4 = 1. / (*gamr - 1.);
2031:     iterno = 10;
2032: /*        find pstar: */
2033:     cl = sqrt(*gaml * *pl / *rl);
2034:     cr = sqrt(*gamr * *pr / *rr);
2035:     wl = *rl * cl;
2036:     wr = *rr * cr;
2037:     csqrl = wl * wl;
2038:     csqrr = wr * wr;
2039:     *pstar = (wl * *pr + wr * *pl) / (wl + wr);
2040:     *pstar = PetscMax(*pstar,smallp);
2041:     pst = *pl / *pr;
2042:     skpr1 = cr * (pst - 1.) * sqrt(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
2043:     d__1 = (*gamr - 1.) / (*gamr * 2.);
2044:     rarepr2 = gascr4 * 2. * cr * (1. - PetscPowScalar(pst, d__1));
2045:     pst = *pr / *pl;
2046:     skpr2 = cl * (pst - 1.) * sqrt(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
2047:     d__1 = (*gaml - 1.) / (*gaml * 2.);
2048:     rarepr1 = gascl4 * 2. * cl * (1. - PetscPowScalar(pst, d__1));
2049:     durl = *uxr - *uxl;
2050:     if (*pr < *pl) {
2051:         if (durl >= rarepr1) {
2052:             iwave = 100;
2053:         } else if (durl <= -skpr1) {
2054:             iwave = 300;
2055:         } else {
2056:             iwave = 400;
2057:         }
2058:     } else {
2059:         if (durl >= rarepr2) {
2060:             iwave = 100;
2061:         } else if (durl <= -skpr2) {
2062:             iwave = 300;
2063:         } else {
2064:             iwave = 200;
2065:         }
2066:     }
2067:     if (iwave == 100) {
2068: /*     1-wave: rarefaction wave, 3-wave: rarefaction wave */
2069: /*     case (100) */
2070:         i__1 = iterno;
2071:         for (i0 = 1; i0 <= i__1; ++i0) {
2072:             d__1 = *pstar / *pl;
2073:             d__2 = 1. / *gaml;
2074:             *rstarl = *rl * PetscPowScalar(d__1, d__2);
2075:             cstarl = sqrt(*gaml * *pstar / *rstarl);
2076:             ustarl = *uxl - gascl4 * 2. * (cstarl - cl);
2077:             zl = *rstarl * cstarl;
2078:             d__1 = *pstar / *pr;
2079:             d__2 = 1. / *gamr;
2080:             *rstarr = *rr * PetscPowScalar(d__1, d__2);
2081:             cstarr = sqrt(*gamr * *pstar / *rstarr);
2082:             ustarr = *uxr + gascr4 * 2. * (cstarr - cr);
2083:             zr = *rstarr * cstarr;
2084:             dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
2085:             *pstar -= dpstar;
2086:             *pstar = PetscMax(*pstar,smallp);
2087:             if (PetscAbsScalar(dpstar) / *pstar <= 1e-8) {
2088:               //break;
2089:             }
2090:         }
2091: /*     1-wave: shock wave, 3-wave: rarefaction wave */
2092:     } else if (iwave == 200) {
2093: /*     case (200) */
2094:         i__1 = iterno;
2095:         for (i0 = 1; i0 <= i__1; ++i0) {
2096:             pst = *pstar / *pl;
2097:             ustarl = *uxl - (pst - 1.) * cl * sqrt(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
2098:             zl = *pl / cl * sqrt(*gaml * 2. * (*gaml - 1. + (*gaml + 1.) * pst)) * (*gaml - 1. + (*gaml + 1.) * pst) / (*gaml * 3. - 1. + (*gaml + 1.) * pst);
2099:             d__1 = *pstar / *pr;
2100:             d__2 = 1. / *gamr;
2101:             *rstarr = *rr * PetscPowScalar(d__1, d__2);
2102:             cstarr = sqrt(*gamr * *pstar / *rstarr);
2103:             zr = *rstarr * cstarr;
2104:             ustarr = *uxr + gascr4 * 2. * (cstarr - cr);
2105:             dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
2106:             *pstar -= dpstar;
2107:             *pstar = PetscMax(*pstar,smallp);
2108:             if (PetscAbsScalar(dpstar) / *pstar <= 1e-8) {
2109:               //break;
2110:             }
2111:         }
2112: /*     1-wave: shock wave, 3-wave: shock */
2113:     } else if (iwave == 300) {
2114: /*     case (300) */
2115:         i__1 = iterno;
2116:         for (i0 = 1; i0 <= i__1; ++i0) {
2117:             pst = *pstar / *pl;
2118:             ustarl = *uxl - (pst - 1.) * cl * sqrt(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
2119:             zl = *pl / cl * sqrt(*gaml * 2. * (*gaml - 1. + (*gaml + 1.) * pst)) * (*gaml - 1. + (*gaml + 1.) * pst) / (*gaml * 3. - 1. + (*gaml + 1.) * pst);
2120:             pst = *pstar / *pr;
2121:             ustarr = *uxr + (pst - 1.) * cr * sqrt(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
2122:             zr = *pr / cr * sqrt(*gamr * 2. * (*gamr - 1. + (*gamr + 1.) * pst)) * (*gamr - 1. + (*gamr + 1.) * pst) / (*gamr * 3. - 1. + (*gamr + 1.) * pst);
2123:             dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
2124:             *pstar -= dpstar;
2125:             *pstar = PetscMax(*pstar,smallp);
2126:             if (PetscAbsScalar(dpstar) / *pstar <= 1e-8) {
2127:               //break;
2128:             }
2129:         }
2130: /*     1-wave: rarefaction wave, 3-wave: shock */
2131:     } else if (iwave == 400) {
2132: /*     case (400) */
2133:         i__1 = iterno;
2134:         for (i0 = 1; i0 <= i__1; ++i0) {
2135:             d__1 = *pstar / *pl;
2136:             d__2 = 1. / *gaml;
2137:             *rstarl = *rl * PetscPowScalar(d__1, d__2);
2138:             cstarl = sqrt(*gaml * *pstar / *rstarl);
2139:             ustarl = *uxl - gascl4 * 2. * (cstarl - cl);
2140:             zl = *rstarl * cstarl;
2141:             pst = *pstar / *pr;
2142:             ustarr = *uxr + (pst - 1.) * cr * sqrt(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
2143:             zr = *pr / cr * sqrt(*gamr * 2. * (*gamr - 1. + (*gamr + 1.) * pst)) * (*gamr - 1. + (*gamr + 1.) * pst) / (*gamr * 3. - 1. + (*gamr + 1.) * pst);
2144:             dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
2145:             *pstar -= dpstar;
2146:             *pstar = PetscMax(*pstar,smallp);
2147:             if (PetscAbsScalar(dpstar) / *pstar <= 1e-8) {
2148:               //break;
2149:             }
2150:         }
2151:     }

2153:     *ustar = (zl * ustarr + zr * ustarl) / (zl + zr);
2154:     if (*pstar > *pl) {
2155:         pst = *pstar / *pl;
2156:         *rstarl = ((*gaml + 1.) * pst + *gaml - 1.) / ((*gaml - 1.) * pst + *
2157:                 gaml + 1.) * *rl;
2158:     }
2159:     if (*pstar > *pr) {
2160:         pst = *pstar / *pr;
2161:         *rstarr = ((*gamr + 1.) * pst + *gamr - 1.) / ((*gamr - 1.) * pst + *
2162:                 gamr + 1.) * *rr;
2163:     }
2164:     return iwave;
2165: }

2167: PetscScalar sign(PetscScalar x)
2168: {
2169:     if(x > 0) return 1.0;
2170:     if(x < 0) return -1.0;
2171:     return 0.0;
2172: }
2173: /*        Riemann Solver */
2174: /* -------------------------------------------------------------------- */
2175: int riemannsolver(PetscScalar *xcen, PetscScalar *xp,
2176:                    PetscScalar *dtt, PetscScalar *rl, PetscScalar *uxl, PetscScalar *pl,
2177:                    PetscScalar *utl, PetscScalar *ubl, PetscScalar *gaml, PetscScalar *rho1l,
2178:                    PetscScalar *rr, PetscScalar *uxr, PetscScalar *pr, PetscScalar *utr,
2179:                    PetscScalar *ubr, PetscScalar *gamr, PetscScalar *rho1r, PetscScalar *rx,
2180:                    PetscScalar *uxm, PetscScalar *px, PetscScalar *utx, PetscScalar *ubx,
2181:                    PetscScalar *gam, PetscScalar *rho1)
2182: {
2183:     /* System generated locals */
2184:     PetscScalar d__1, d__2;

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

2190:     if (*rl == *rr && *pr == *pl && *uxl == *uxr && *gaml == *gamr) {
2191:         *rx = *rl;
2192:         *px = *pl;
2193:         *uxm = *uxl;
2194:         *gam = *gaml;
2195:         x2 = *xcen + *uxm * *dtt;

2197:         if (*xp >= x2) {
2198:             *utx = *utr;
2199:             *ubx = *ubr;
2200:             *rho1 = *rho1r;
2201:         } else {
2202:             *utx = *utl;
2203:             *ubx = *ubl;
2204:             *rho1 = *rho1l;
2205:         }
2206:         return 0;
2207:     }
2208:     int iwave = riem1mdt(gaml, gamr, rl, pl, uxl, rr, pr, uxr, &rstarl, &rstarr, &pstar, &ustar);

2210:     x2 = *xcen + ustar * *dtt;
2211:     d__1 = *xp - x2;
2212:     sgn0 = sign(d__1);
2213: /*            x is in 3-wave if sgn0 = 1 */
2214: /*            x is in 1-wave if sgn0 = -1 */
2215:     r0 = cvmgm_(rl, rr, &sgn0);
2216:     p0 = cvmgm_(pl, pr, &sgn0);
2217:     u0 = cvmgm_(uxl, uxr, &sgn0);
2218:     *gam = cvmgm_(gaml, gamr, &sgn0);
2219:     gasc1 = *gam - 1.;
2220:     gasc2 = (*gam + 1.) * .5;
2221:     gasc3 = gasc2 / *gam;
2222:     gasc4 = 1. / (*gam - 1.);
2223:     c0 = sqrt(*gam * p0 / r0);
2224:     streng = pstar - p0;
2225:     w0 = *gam * r0 * p0 * (gasc3 * streng / p0 + 1.);
2226:     rstars = r0 / (1. - r0 * streng / w0);
2227:     d__1 = p0 / pstar;
2228:     d__2 = -1. / *gam;
2229:     rstarr = r0 * PetscPowScalar(d__1, d__2);
2230:     rstar = cvmgm_(&rstarr, &rstars, &streng);
2231:     w0 = sqrt(w0);
2232:     cstar = sqrt(*gam * pstar / rstar);
2233:     wsp0 = u0 + sgn0 * c0;
2234:     wspst = ustar + sgn0 * cstar;
2235:     ushock = ustar + sgn0 * w0 / rstar;
2236:     wspst = cvmgp_(&ushock, &wspst, &streng);
2237:     wsp0 = cvmgp_(&ushock, &wsp0, &streng);
2238:     x0 = *xcen + wsp0 * *dtt;
2239:     xstar = *xcen + wspst * *dtt;
2240: /*           using gas formula to evaluate rarefaction wave */
2241: /*            ri : reiman invariant */
2242:     ri = u0 - sgn0 * 2. * gasc4 * c0;
2243:     cx = sgn0 * .5 * gasc1 / gasc2 * ((*xp - *xcen) / *dtt - ri);
2244:     *uxm = ri + sgn0 * 2. * gasc4 * cx;
2245:     s = p0 / PetscPowScalar(r0, *gam);
2246:     d__1 = cx * cx / (*gam * s);
2247:     *rx = PetscPowScalar(d__1, gasc4);
2248:     *px = cx * cx * *rx / *gam;
2249:     d__1 = sgn0 * (x0 - *xp);
2250:     *rx = cvmgp_(rx, &r0, &d__1);
2251:     d__1 = sgn0 * (x0 - *xp);
2252:     *px = cvmgp_(px, &p0, &d__1);
2253:     d__1 = sgn0 * (x0 - *xp);
2254:     *uxm = cvmgp_(uxm, &u0, &d__1);
2255:     d__1 = sgn0 * (xstar - *xp);
2256:     *rx = cvmgm_(rx, &rstar, &d__1);
2257:     d__1 = sgn0 * (xstar - *xp);
2258:     *px = cvmgm_(px, &pstar, &d__1);
2259:     d__1 = sgn0 * (xstar - *xp);
2260:     *uxm = cvmgm_(uxm, &ustar, &d__1);
2261:     if (*xp >= x2) {
2262:         *utx = *utr;
2263:         *ubx = *ubr;
2264:         *rho1 = *rho1r;
2265:     } else {
2266:         *utx = *utl;
2267:         *ubx = *ubl;
2268:         *rho1 = *rho1l;
2269:     }
2270:     return iwave;
2271: }
2272: int godunovflux( const PetscScalar *ul, const PetscScalar *ur,
2273:                  PetscScalar *flux, const PetscScalar *nn, const int *ndim,
2274:                  const PetscScalar *gamma)
2275: {
2276:     /* System generated locals */
2277:   int i__1,iwave;
2278:     PetscScalar d__1, d__2, d__3;

2280:     /* Local variables */
2281:     static int k;
2282:     static PetscScalar bn[3], fn, ft, tg[3], pl, rl, pm, pr, rr, xp, ubl, ubm,
2283:             ubr, dtt, unm, tmp, utl, utm, uxl, utr, uxr, gaml, gamm, gamr,
2284:             xcen, rhom, rho1l, rho1m, rho1r;
2285:     /* Parameter adjustments */
2286:     --nn;
2287:     --flux;
2288:     --ur;
2289:     --ul;

2291:     /* Function Body */
2292:     xcen = 0.;
2293:     xp = 0.;
2294:     i__1 = *ndim;
2295:     for (k = 1; k <= i__1; ++k) {
2296:         tg[k - 1] = 0.;
2297:         bn[k - 1] = 0.;
2298:     }
2299:     dtt = 1.;
2300:     if (*ndim == 3) {
2301:         if (nn[1] == 0. && nn[2] == 0.) {
2302:             tg[0] = 1.;
2303:         } else {
2304:             tg[0] = -nn[2];
2305:             tg[1] = nn[1];
2306:         }
2307: /*           tmp=dsqrt(tg(1)**2+tg(2)**2) */
2308: /*           tg=tg/tmp */
2309:         bn[0] = -nn[3] * tg[1];
2310:         bn[1] = nn[3] * tg[0];
2311:         bn[2] = nn[1] * tg[1] - nn[2] * tg[0];
2312: /* Computing 2nd power */
2313:         d__1 = bn[0];
2314: /* Computing 2nd power */
2315:         d__2 = bn[1];
2316: /* Computing 2nd power */
2317:         d__3 = bn[2];
2318:         tmp = sqrt(d__1 * d__1 + d__2 * d__2 + d__3 * d__3);
2319:         i__1 = *ndim;
2320:         for (k = 1; k <= i__1; ++k) {
2321:             bn[k - 1] /= tmp;
2322:         }
2323:     } else if (*ndim == 2) {
2324:         tg[0] = -nn[2];
2325:         tg[1] = nn[1];
2326: /*           tmp=dsqrt(tg(1)**2+tg(2)**2) */
2327: /*           tg=tg/tmp */
2328:         bn[0] = 0.;
2329:         bn[1] = 0.;
2330:         bn[2] = 1.;
2331:     }
2332:     rl = ul[1];
2333:     rr = ur[1];
2334:     uxl = 0.;
2335:     uxr = 0.;
2336:     utl = 0.;
2337:     utr = 0.;
2338:     ubl = 0.;
2339:     ubr = 0.;
2340:     i__1 = *ndim;
2341:     for (k = 1; k <= i__1; ++k) {
2342:         uxl += ul[k + 1] * nn[k];
2343:         uxr += ur[k + 1] * nn[k];
2344:         utl += ul[k + 1] * tg[k - 1];
2345:         utr += ur[k + 1] * tg[k - 1];
2346:         ubl += ul[k + 1] * bn[k - 1];
2347:         ubr += ur[k + 1] * bn[k - 1];
2348:     }
2349:     uxl /= rl;
2350:     uxr /= rr;
2351:     utl /= rl;
2352:     utr /= rr;
2353:     ubl /= rl;
2354:     ubr /= rr;

2356:     gaml = *gamma;
2357:     gamr = *gamma;
2358: /* Computing 2nd power */
2359:     d__1 = uxl;
2360: /* Computing 2nd power */
2361:     d__2 = utl;
2362: /* Computing 2nd power */
2363:     d__3 = ubl;
2364:     pl = (*gamma - 1.) * (ul[*ndim + 2] - rl * .5 * (d__1 * d__1 + d__2 * d__2 + d__3 * d__3));
2365: /* Computing 2nd power */
2366:     d__1 = uxr;
2367: /* Computing 2nd power */
2368:     d__2 = utr;
2369: /* Computing 2nd power */
2370:     d__3 = ubr;
2371:     pr = (*gamma - 1.) * (ur[*ndim + 2] - rr * .5 * (d__1 * d__1 + d__2 * d__2 + d__3 * d__3));
2372:     rho1l = rl;
2373:     rho1r = rr;

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

2379:     flux[1] = rhom * unm;
2380:     fn = rhom * unm * unm + pm;
2381:     ft = rhom * unm * utm;
2382: /*           flux(2)=fn*nn(1)+ft*nn(2) */
2383: /*           flux(3)=fn*tg(1)+ft*tg(2) */
2384:     flux[2] = fn * nn[1] + ft * tg[0];
2385:     flux[3] = fn * nn[2] + ft * tg[1];
2386: /*           flux(2)=rhom*unm*(unm)+pm */
2387: /*           flux(3)=rhom*(unm)*utm */
2388:     if (*ndim == 3) {
2389:         flux[4] = rhom * unm * ubm;
2390:     }
2391:     flux[*ndim + 2] = (rhom * .5 * (unm * unm + utm * utm + ubm * ubm) + gamm / (gamm - 1.) * pm) * unm;
2392:     return iwave;
2393: } /* godunovflux_ */

2395: /* Subroutine to set up the initial conditions for the */
2396: /* Shock Interface interaction or linear wave (Ravi Samtaney,Mark Adams). */
2397: /* ----------------------------------------------------------------------- */
2398: int projecteqstate(PetscScalar wc[], const PetscScalar ueq[], const PetscScalar lv[][3])
2399: {
2400:   int j,k;
2401: /*      Wc=matmul(lv,Ueq) 3 vars */
2402:   for (k = 0; k < 3; ++k) {
2403:     wc[k] = 0.;
2404:     for (j = 0; j < 3; ++j) {
2405:       wc[k] += lv[k][j]*ueq[j];
2406:     }
2407:   }
2408:   return 0;
2409: }
2410: /* ----------------------------------------------------------------------- */
2411: int projecttoprim(PetscScalar v[], const PetscScalar wc[], PetscScalar rv[][3])
2412: {
2413:   int k,j;
2414:   /*      V=matmul(rv,WC) 3 vars */
2415:   for (k = 0; k < 3; ++k) {
2416:     v[k] = 0.;
2417:     for (j = 0; j < 3; ++j) {
2418:       v[k] += rv[k][j]*wc[j];
2419:     }
2420:   }
2421:   return 0;
2422: }
2423: /* ---------------------------------------------------------------------- */
2424: int eigenvectors(PetscScalar rv[][3], PetscScalar lv[][3], const PetscScalar ueq[], PetscScalar gamma)
2425: {
2426:   int j,k;
2427:   PetscScalar rho,csnd,p0,u;

2429:   for (k = 0; k < 3; ++k) for (j = 0; j < 3; ++j) { lv[k][j] = 0.; rv[k][j] = 0.; }
2430:   rho = ueq[0];
2431:   u = ueq[1];
2432:   p0 = ueq[2];
2433:   csnd = PetscSqrtScalar(gamma * p0 / rho);
2434:   lv[0][1] = rho * .5;
2435:   lv[0][2] = -.5 / csnd;
2436:   lv[1][0] = csnd;
2437:   lv[1][2] = -1. / csnd;
2438:   lv[2][1] = rho * .5;
2439:   lv[2][2] = .5 / csnd;
2440:   rv[0][0] = -1. / csnd;
2441:   rv[1][0] = 1. / rho;
2442:   rv[2][0] = -csnd;
2443:   rv[0][1] = 1. / csnd;
2444:   rv[0][2] = 1. / csnd;
2445:   rv[1][2] = 1. / rho;
2446:   rv[2][2] = csnd;
2447:   return 0;
2448: }
2449: int initLinearWave(EulerNode *ux, const PetscScalar gamma, const PetscReal coord[], const PetscReal Lx)
2450: {
2451:   PetscScalar p0,u0,wcp[3],wc[3];
2452:   PetscScalar lv[3][3];
2453:   PetscScalar vp[3];
2454:   PetscScalar rv[3][3];
2455:   PetscScalar eps, ueq[3], rho0, twopi;

2457:   /* Function Body */
2458:   twopi = 2.*M_PI;
2459:   eps = 1e-4; /* perturbation */
2460:   rho0 = 1e3;   /* density of water */
2461:   p0 = 101325.; /* init pressure of 1 atm (?) */
2462:   u0 = 0.;
2463:   ueq[0] = rho0;
2464:   ueq[1] = u0;
2465:   ueq[2] = p0;
2466:   /* Project initial state to characteristic variables */
2467:   eigenvectors(rv, lv, ueq, gamma);
2468:   projecteqstate(wc, ueq, lv);
2469:   wcp[0] = wc[0];
2470:   wcp[1] = wc[1];
2471:   wcp[2] = wc[2] + eps * cos(coord[0] * 2. * twopi / Lx);
2472:   projecttoprim(vp, wcp, rv);
2473:   ux->r = vp[0]; /* density */
2474:   ux->ru[0] = vp[0] * vp[1]; /* x momentum */
2475:   ux->ru[1] = 0.;
2476: #if defined DIM > 2
2477:   if (dim>2) ux->ru[2] = 0.;
2478: #endif
2479:   /* E = rho * e + rho * v^2/2 = p/(gam-1) + rho*v^2/2 */
2480:   ux->E = vp[2]/(gamma - 1.) + 0.5*vp[0]*vp[1]*vp[1];
2481:   return 0;
2482: }