Actual source code: ex11.c

petsc-master 2016-09-28
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}};

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

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

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

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

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

203:       rad2 = 0.;
204:       for (i = 0; i < dim; i++) {
205:         comp2[i] = qp[i] * qp[i];
206:         rad2    += comp2[i];
207:       }

209:       wind[0] = -qp[1];
210:       wind[1] = qp[0];
211:       if (rad2 > 1.) {
212:         PetscInt  maxI = 0;
213:         PetscReal maxComp2 = comp2[0];

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

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

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

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

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

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

295:   /* Register "canned" boundary conditions and defaults for where to apply. */
296:   PetscDSAddBoundary(prob, PETSC_TRUE, "inflow",  "Face Sets", 0, 0, NULL, (void (*)()) PhysicsBoundary_Advect_Inflow,  ALEN(inflowids),   inflowids,  phys);
297:   PetscDSAddBoundary(prob, PETSC_FALSE, "outflow", "Face Sets", 0, 0, NULL, (void (*)()) PhysicsBoundary_Advect_Outflow, ALEN(outflowids), outflowids, phys);
298:   return(0);
299: }

303: static PetscErrorCode PhysicsCreate_Advect(Model mod,Physics phys,PetscOptionItems *PetscOptionsObject)
304: {
305:   Physics_Advect *advect;

309:   phys->field_desc = PhysicsFields_Advect;
310:   phys->riemann    = (PetscRiemannFunc)PhysicsRiemann_Advect;
311:   PetscNew(&advect);
312:   phys->data       = advect;
313:   mod->setupbc = SetUpBC_Advect;

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

356: /******************* Shallow Water ********************/
357: typedef struct {
358:   PetscReal gravity;
359:   PetscReal boundaryHeight;
360:   struct {
361:     PetscInt Height;
362:     PetscInt Speed;
363:     PetscInt Energy;
364:   } functional;
365: } Physics_SW;
366: typedef struct {
367:   PetscScalar vals[0];
368:   PetscScalar h;
369:   PetscScalar uh[DIM];
370: } SWNode;

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

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

388:   Scale2(1./x->h,x->uh,u);
389:   uhn  = Dot2(x->uh,n);
390:   f->h = uhn;
391:   for (i=0; i<DIM; i++) f->uh[i] = u[i] * uhn + sw->gravity * PetscSqr(x->h) * n[i];
392:   return(0);
393: }

397: static PetscErrorCode PhysicsBoundary_SW_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
398: {
400:   xG[0] = xI[0];
401:   xG[1] = -xI[1];
402:   xG[2] = -xI[2];
403:   return(0);
404: }

408: static void PhysicsRiemann_SW(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscScalar *flux, Physics phys)
409: {
410:   Physics_SW   *sw = (Physics_SW*)phys->data;
411:   PetscReal    cL,cR,speed,nn[DIM];
412:   const SWNode *uL = (const SWNode*)xL,*uR = (const SWNode*)xR;
413:   SWNode       fL,fR;
414:   PetscInt     i;

416:   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"); */
417:   nn[0] = n[0];
418:   nn[1] = n[1];
419:   Normalize2(nn);
420:   SWFlux(phys,nn,uL,&fL);
421:   SWFlux(phys,nn,uR,&fR);
422:   cL    = PetscSqrtReal(sw->gravity*PetscRealPart(uL->h));
423:   cR    = PetscSqrtReal(sw->gravity*PetscRealPart(uR->h)); /* gravity wave speed */
424:   speed = PetscMax(PetscAbsScalar(Dot2(uL->uh,nn)/uL->h) + cL,PetscAbsScalar(Dot2(uR->uh,nn)/uR->h) + cR);
425:   for (i=0; i<1+dim; i++) flux[i] = (0.5*(fL.vals[i] + fR.vals[i]) + 0.5*speed*(xL[i] - xR[i])) * Norm2(n);
426: }

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

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

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

457:   h = PetscRealPart(x->h);
458:   Scale2(1./x->h,x->uh,u);
459:   f[sw->functional.Height] = h;
460:   f[sw->functional.Speed]  = Norm2(u) + PetscSqrtReal(sw->gravity*h);
461:   f[sw->functional.Energy] = 0.5*(Dot2(x->uh,u) + sw->gravity*PetscSqr(h));
462:   return(0);
463: }

467: static PetscErrorCode SetUpBC_SW(PetscDS prob,Physics phys)
468: {
470:   const PetscInt wallids[] = {100,101,200,300};
472:   PetscDSAddBoundary(prob, PETSC_FALSE, "wall", "Face Sets", 0, 0, NULL, (void (*)()) PhysicsBoundary_SW_Wall, ALEN(wallids), wallids, phys);
473:   return(0);
474: }

478: static PetscErrorCode PhysicsCreate_SW(Model mod,Physics phys,PetscOptionItems *PetscOptionsObject)
479: {
480:   Physics_SW     *sw;

484:   phys->field_desc = PhysicsFields_SW;
485:   phys->riemann = (PetscRiemannFunc) PhysicsRiemann_SW;
486:   PetscNew(&sw);
487:   phys->data    = sw;
488:   mod->setupbc  = SetUpBC_SW;

490:   PetscOptionsHead(PetscOptionsObject,"SW options");
491:   {
492:     sw->gravity = 1.0;
493:     PetscOptionsReal("-sw_gravity","Gravitational constant","",sw->gravity,&sw->gravity,NULL);
494:   }
495:   PetscOptionsTail();
496:   phys->maxspeed = PetscSqrtReal(2.0*sw->gravity); /* Mach 1 for depth of 2 */

498:   ModelSolutionSetDefault(mod,PhysicsSolution_SW,phys);
499:   ModelFunctionalRegister(mod,"Height",&sw->functional.Height,PhysicsFunctional_SW,phys);
500:   ModelFunctionalRegister(mod,"Speed",&sw->functional.Speed,PhysicsFunctional_SW,phys);
501:   ModelFunctionalRegister(mod,"Energy",&sw->functional.Energy,PhysicsFunctional_SW,phys);

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

505:   return(0);
506: }

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

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

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

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

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

599:   // set phys->maxspeed: (mod->maxspeed = phys->maxspeed) in main;
600:   eu->sound(&gamma,uu,&c);
601:   c = PetscAbsScalar(uu->ru[0]/uu->r) + c;
602:   if (c > phys->maxspeed) phys->maxspeed = c;

604:   return(0);
605: }

609: static PetscErrorCode Pressure_PG(const PetscReal gamma,const EulerNode *x,PetscScalar *p)
610: {
611:   PetscScalar ru2;
613:   ru2  = DotDIM(x->ru,x->ru);
614:   (*p)=(x->E - 0.5*ru2/x->r)*(gamma - 1.0); /* (E - rho V^2/2)(gamma-1) = e rho (gamma-1) */
615:   return(0);
616: }

620: static PetscErrorCode SpeedOfSound_PG(const PetscReal *gamma, const EulerNode *x, PetscScalar *c)
621: {
622:   PetscScalar p;

625:   Pressure_PG(*gamma,x,&p);
626:   if (p<0.) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"negative pressure time %g -- NEED TO FIX!!!!!!",p);
627:   /* pars[EULER_PAR_GAMMA] = heat capacity ratio */
628:   (*c)=PetscSqrtScalar(*gamma * p / x->r);
629:   return(0);
630: }

634: /*
635:  * x = (rho,rho*(u_1),...,rho*e)^T
636:  * x_t+div(f_1(x))+...+div(f_DIM(x)) = 0
637:  *
638:  * f_i(x) = u_i*x+(0,0,...,p,...,p*u_i)^T
639:  *
640:  */
641: static PetscErrorCode EulerFlux(Physics phys,const PetscReal *n,const EulerNode *x,EulerNode *f)
642: {
643:   Physics_Euler *eu = (Physics_Euler*)phys->data;
644:   PetscScalar   nu,p;
645:   PetscInt      i;

648:   Pressure_PG(eu->pars[EULER_PAR_GAMMA],x,&p);
649:   nu = DotDIM(x->ru,n);
650:   f->r = nu;   /* A rho u */
651:   nu /= x->r;  /* A u */
652:   for (i=0; i<DIM; i++) f->ru[i] = nu * x->ru[i] + n[i]*p;  /* r u^2 + p */
653:   f->E = nu * (x->E + p); /* u(e+p) */
654:   return(0);
655: }

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

695:   for (i=0,s2=0.; i<DIM; i++) {
696:     nn[i] = n[i];
697:     s2 += nn[i]*nn[i];
698:   }
699:   s2 = PetscSqrtReal(s2); /* |n|_2 = sum(n^2)^1/2 */
700:   for (i=0.; i<DIM; i++) nn[i] /= s2;
701:   if (0) { /* Rusanov */
702:     const EulerNode *uL = (const EulerNode*)xL,*uR = (const EulerNode*)xR;
703:     EulerNode       fL,fR;
704:     EulerFlux(phys,nn,uL,&fL);
705:     EulerFlux(phys,nn,uR,&fR);
706:     eu->sound(&eu->pars[EULER_PAR_GAMMA],uL,&cL);if (ierr) exit(13);
707:     eu->sound(&eu->pars[EULER_PAR_GAMMA],uR,&cR);if (ierr) exit(14);
708:     velL = DotDIM(uL->ru,nn)/uL->r;
709:     velR = DotDIM(uR->ru,nn)/uR->r;
710:     speed = PetscMax(PetscAbsScalar(velR) + cR,PetscAbsScalar(velL) + cL);
711:     for (i=0; i<2+dim; i++) flux[i] = 0.5*((fL.vals[i]+fR.vals[i]) + speed*(xL[i] - xR[i]))*s2;
712:   }
713:   else {
714:     int dim = DIM;
715:     /* int iwave =  */
716:     godunovflux(xL, xR, flux, nn, &dim, &eu->pars[EULER_PAR_GAMMA]);
717:     for (i=0; i<2+dim; i++) flux[i] *= s2;
718:   }
719:   PetscFunctionReturnVoid();
720: }

724: static PetscErrorCode PhysicsFunctional_Euler(Model mod,PetscReal time,const PetscReal *coord,const PetscScalar *xx,PetscReal *f,void *ctx)
725: {
726:   Physics         phys = (Physics)ctx;
727:   Physics_Euler   *eu  = (Physics_Euler*)phys->data;
728:   const EulerNode *x   = (const EulerNode*)xx;
729:   PetscScalar     p;

732:   f[eu->monitor.Density]  = x->r;
733:   f[eu->monitor.Momentum] = NormDIM(x->ru);
734:   f[eu->monitor.Energy]   = x->E;
735:   f[eu->monitor.Speed]    = NormDIM(x->ru)/x->r;
736:   Pressure_PG(eu->pars[EULER_PAR_GAMMA], x, &p);
737:   f[eu->monitor.Pressure] = p;
738:   return(0);
739: }

743: static PetscErrorCode SetUpBC_Euler(PetscDS prob,Physics phys)
744: {
745:   PetscErrorCode  ierr;
746:   Physics_Euler   *eu = phys->data;
747:   if (eu->type == EULER_LINEAR_WAVE) {
748:     const PetscInt wallids[] = {100,101};
749:     PetscDSAddBoundary(prob, PETSC_FALSE, "wall", "Face Sets", 0, 0, NULL, (void (*)()) PhysicsBoundary_Euler_Wall, ALEN(wallids), wallids, phys);
750:   }
751:   else {
752:     const PetscInt wallids[] = {100,101,200,300};
753:     PetscDSAddBoundary(prob, PETSC_FALSE, "wall", "Face Sets", 0, 0, NULL, (void (*)()) PhysicsBoundary_Euler_Wall, ALEN(wallids), wallids, phys);
754:   }
755:   return(0);
756: }

760: static PetscErrorCode PhysicsCreate_Euler(Model mod,Physics phys,PetscOptionItems *PetscOptionsObject)
761: {
762:   Physics_Euler   *eu;
763:   PetscErrorCode  ierr;

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

828:   return(0);
829: }

833: static PetscErrorCode ErrorIndicator_Simple(PetscInt dim, PetscReal volume, PetscInt numComps, const PetscScalar u[], const PetscScalar grad[], PetscReal *error, void *ctx)
834: {
835:   PetscReal      err = 0.;
836:   PetscInt       i, j;

839:   for (i = 0; i < numComps; i++) {
840:     for (j = 0; j < dim; j++) {
841:       err += PetscSqr(PetscRealPart(grad[i * dim + j]));
842:     }
843:   }
844:   *error = volume * err;
845:   return(0);
846: }

850: PetscErrorCode ConstructCellBoundary(DM dm, User user)
851: {
852:   const char     *name   = "Cell Sets";
853:   const char     *bdname = "split faces";
854:   IS             regionIS, innerIS;
855:   const PetscInt *regions, *cells;
856:   PetscInt       numRegions, innerRegion, numCells, c;
857:   PetscInt       cStart, cEnd, cEndInterior, fStart, fEnd;
858:   PetscBool      hasLabel;

862:   DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
863:   DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);
864:   DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);

866:   DMHasLabel(dm, name, &hasLabel);
867:   if (!hasLabel) return(0);
868:   DMGetLabelSize(dm, name, &numRegions);
869:   if (numRegions != 2) return(0);
870:   /* Get the inner id */
871:   DMGetLabelIdIS(dm, name, &regionIS);
872:   ISGetIndices(regionIS, &regions);
873:   innerRegion = regions[0];
874:   ISRestoreIndices(regionIS, &regions);
875:   ISDestroy(&regionIS);
876:   /* Find the faces between cells in different regions, could call DMPlexCreateNeighborCSR() */
877:   DMGetStratumIS(dm, name, innerRegion, &innerIS);
878:   ISGetLocalSize(innerIS, &numCells);
879:   ISGetIndices(innerIS, &cells);
880:   DMCreateLabel(dm, bdname);
881:   for (c = 0; c < numCells; ++c) {
882:     const PetscInt cell = cells[c];
883:     const PetscInt *faces;
884:     PetscInt       numFaces, f;

886:     if ((cell < cStart) || (cell >= cEnd)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Got invalid point %d which is not a cell", cell);
887:     DMPlexGetConeSize(dm, cell, &numFaces);
888:     DMPlexGetCone(dm, cell, &faces);
889:     for (f = 0; f < numFaces; ++f) {
890:       const PetscInt face = faces[f];
891:       const PetscInt *neighbors;
892:       PetscInt       nC, regionA, regionB;

894:       if ((face < fStart) || (face >= fEnd)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Got invalid point %d which is not a face", face);
895:       DMPlexGetSupportSize(dm, face, &nC);
896:       if (nC != 2) continue;
897:       DMPlexGetSupport(dm, face, &neighbors);
898:       if ((neighbors[0] >= cEndInterior) || (neighbors[1] >= cEndInterior)) continue;
899:       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]);
900:       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]);
901:       DMGetLabelValue(dm, name, neighbors[0], &regionA);
902:       DMGetLabelValue(dm, name, neighbors[1], &regionB);
903:       if (regionA < 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid label %s: Cell %d has no value", name, neighbors[0]);
904:       if (regionB < 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid label %s: Cell %d has no value", name, neighbors[1]);
905:       if (regionA != regionB) {
906:         DMSetLabelValue(dm, bdname, faces[f], 1);
907:       }
908:     }
909:   }
910:   ISRestoreIndices(innerIS, &cells);
911:   ISDestroy(&innerIS);
912:   {
913:     DMLabel label;

915:     DMGetLabel(dm, bdname, &label);
916:     DMLabelView(label, PETSC_VIEWER_STDOUT_WORLD);
917:   }
918:   return(0);
919: }

923: /* Right now, I have just added duplicate faces, which see both cells. We can
924: - Add duplicate vertices and decouple the face cones
925: - Disconnect faces from cells across the rotation gap
926: */
927: PetscErrorCode SplitFaces(DM *dmSplit, const char labelName[], User user)
928: {
929:   DM             dm = *dmSplit, sdm;
930:   PetscSF        sfPoint, gsfPoint;
931:   PetscSection   coordSection, newCoordSection;
932:   Vec            coordinates;
933:   IS             idIS;
934:   const PetscInt *ids;
935:   PetscInt       *newpoints;
936:   PetscInt       dim, depth, maxConeSize, maxSupportSize, numLabels, numGhostCells;
937:   PetscInt       numFS, fs, pStart, pEnd, p, cEnd, cEndInterior, vStart, vEnd, v, fStart, fEnd, newf, d, l;
938:   PetscBool      hasLabel;

942:   DMHasLabel(dm, labelName, &hasLabel);
943:   if (!hasLabel) return(0);
944:   DMCreate(PetscObjectComm((PetscObject)dm), &sdm);
945:   DMSetType(sdm, DMPLEX);
946:   DMGetDimension(dm, &dim);
947:   DMSetDimension(sdm, dim);

949:   DMGetLabelIdIS(dm, labelName, &idIS);
950:   ISGetLocalSize(idIS, &numFS);
951:   ISGetIndices(idIS, &ids);

953:   user->numSplitFaces = 0;
954:   for (fs = 0; fs < numFS; ++fs) {
955:     PetscInt numBdFaces;

957:     DMGetStratumSize(dm, labelName, ids[fs], &numBdFaces);
958:     user->numSplitFaces += numBdFaces;
959:   }
960:   DMPlexGetChart(dm, &pStart, &pEnd);
961:   pEnd += user->numSplitFaces;
962:   DMPlexSetChart(sdm, pStart, pEnd);
963:   DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);
964:   DMPlexSetHybridBounds(sdm, cEndInterior, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE);
965:   DMPlexGetHeightStratum(dm, 0, NULL, &cEnd);
966:   numGhostCells = cEnd - cEndInterior;
967:   /* Set cone and support sizes */
968:   DMPlexGetDepth(dm, &depth);
969:   for (d = 0; d <= depth; ++d) {
970:     DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);
971:     for (p = pStart; p < pEnd; ++p) {
972:       PetscInt newp = p;
973:       PetscInt size;

975:       DMPlexGetConeSize(dm, p, &size);
976:       DMPlexSetConeSize(sdm, newp, size);
977:       DMPlexGetSupportSize(dm, p, &size);
978:       DMPlexSetSupportSize(sdm, newp, size);
979:     }
980:   }
981:   DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);
982:   for (fs = 0, newf = fEnd; fs < numFS; ++fs) {
983:     IS             faceIS;
984:     const PetscInt *faces;
985:     PetscInt       numFaces, f;

987:     DMGetStratumIS(dm, labelName, ids[fs], &faceIS);
988:     ISGetLocalSize(faceIS, &numFaces);
989:     ISGetIndices(faceIS, &faces);
990:     for (f = 0; f < numFaces; ++f, ++newf) {
991:       PetscInt size;

993:       /* Right now I think that both faces should see both cells */
994:       DMPlexGetConeSize(dm, faces[f], &size);
995:       DMPlexSetConeSize(sdm, newf, size);
996:       DMPlexGetSupportSize(dm, faces[f], &size);
997:       DMPlexSetSupportSize(sdm, newf, size);
998:     }
999:     ISRestoreIndices(faceIS, &faces);
1000:     ISDestroy(&faceIS);
1001:   }
1002:   DMSetUp(sdm);
1003:   /* Set cones and supports */
1004:   DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);
1005:   PetscMalloc1(PetscMax(maxConeSize, maxSupportSize), &newpoints);
1006:   DMPlexGetChart(dm, &pStart, &pEnd);
1007:   for (p = pStart; p < pEnd; ++p) {
1008:     const PetscInt *points, *orientations;
1009:     PetscInt       size, i, newp = p;

1011:     DMPlexGetConeSize(dm, p, &size);
1012:     DMPlexGetCone(dm, p, &points);
1013:     DMPlexGetConeOrientation(dm, p, &orientations);
1014:     for (i = 0; i < size; ++i) newpoints[i] = points[i];
1015:     DMPlexSetCone(sdm, newp, newpoints);
1016:     DMPlexSetConeOrientation(sdm, newp, orientations);
1017:     DMPlexGetSupportSize(dm, p, &size);
1018:     DMPlexGetSupport(dm, p, &points);
1019:     for (i = 0; i < size; ++i) newpoints[i] = points[i];
1020:     DMPlexSetSupport(sdm, newp, newpoints);
1021:   }
1022:   PetscFree(newpoints);
1023:   for (fs = 0, newf = fEnd; fs < numFS; ++fs) {
1024:     IS             faceIS;
1025:     const PetscInt *faces;
1026:     PetscInt       numFaces, f;

1028:     DMGetStratumIS(dm, labelName, ids[fs], &faceIS);
1029:     ISGetLocalSize(faceIS, &numFaces);
1030:     ISGetIndices(faceIS, &faces);
1031:     for (f = 0; f < numFaces; ++f, ++newf) {
1032:       const PetscInt *points;

1034:       DMPlexGetCone(dm, faces[f], &points);
1035:       DMPlexSetCone(sdm, newf, points);
1036:       DMPlexGetSupport(dm, faces[f], &points);
1037:       DMPlexSetSupport(sdm, newf, points);
1038:     }
1039:     ISRestoreIndices(faceIS, &faces);
1040:     ISDestroy(&faceIS);
1041:   }
1042:   ISRestoreIndices(idIS, &ids);
1043:   ISDestroy(&idIS);
1044:   DMPlexStratify(sdm);
1045:   /* Convert coordinates */
1046:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
1047:   DMGetCoordinateSection(dm, &coordSection);
1048:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &newCoordSection);
1049:   PetscSectionSetNumFields(newCoordSection, 1);
1050:   PetscSectionSetFieldComponents(newCoordSection, 0, dim);
1051:   PetscSectionSetChart(newCoordSection, vStart, vEnd);
1052:   for (v = vStart; v < vEnd; ++v) {
1053:     PetscSectionSetDof(newCoordSection, v, dim);
1054:     PetscSectionSetFieldDof(newCoordSection, v, 0, dim);
1055:   }
1056:   PetscSectionSetUp(newCoordSection);
1057:   DMSetCoordinateSection(sdm, PETSC_DETERMINE, newCoordSection);
1058:   PetscSectionDestroy(&newCoordSection); /* relinquish our reference */
1059:   DMGetCoordinatesLocal(dm, &coordinates);
1060:   DMSetCoordinatesLocal(sdm, coordinates);
1061:   /* Convert labels */
1062:   DMGetNumLabels(dm, &numLabels);
1063:   for (l = 0; l < numLabels; ++l) {
1064:     const char *lname;
1065:     PetscBool  isDepth;

1067:     DMGetLabelName(dm, l, &lname);
1068:     PetscStrcmp(lname, "depth", &isDepth);
1069:     if (isDepth) continue;
1070:     DMCreateLabel(sdm, lname);
1071:     DMGetLabelIdIS(dm, lname, &idIS);
1072:     ISGetLocalSize(idIS, &numFS);
1073:     ISGetIndices(idIS, &ids);
1074:     for (fs = 0; fs < numFS; ++fs) {
1075:       IS             pointIS;
1076:       const PetscInt *points;
1077:       PetscInt       numPoints;

1079:       DMGetStratumIS(dm, lname, ids[fs], &pointIS);
1080:       ISGetLocalSize(pointIS, &numPoints);
1081:       ISGetIndices(pointIS, &points);
1082:       for (p = 0; p < numPoints; ++p) {
1083:         PetscInt newpoint = points[p];

1085:         DMSetLabelValue(sdm, lname, newpoint, ids[fs]);
1086:       }
1087:       ISRestoreIndices(pointIS, &points);
1088:       ISDestroy(&pointIS);
1089:     }
1090:     ISRestoreIndices(idIS, &ids);
1091:     ISDestroy(&idIS);
1092:   }
1093:   /* Convert pointSF */
1094:   const PetscSFNode *remotePoints;
1095:   PetscSFNode       *gremotePoints;
1096:   const PetscInt    *localPoints;
1097:   PetscInt          *glocalPoints,*newLocation,*newRemoteLocation;
1098:   PetscInt          numRoots, numLeaves;
1099:   PetscMPIInt       numProcs;

1101:   MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);
1102:   DMGetPointSF(dm, &sfPoint);
1103:   DMGetPointSF(sdm, &gsfPoint);
1104:   DMPlexGetChart(dm,&pStart,&pEnd);
1105:   PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);
1106:   if (numRoots >= 0) {
1107:     PetscMalloc2(numRoots,&newLocation,pEnd-pStart,&newRemoteLocation);
1108:     for (l=0; l<numRoots; l++) newLocation[l] = l; /* + (l >= cEnd ? numGhostCells : 0); */
1109:     PetscSFBcastBegin(sfPoint, MPIU_INT, newLocation, newRemoteLocation);
1110:     PetscSFBcastEnd(sfPoint, MPIU_INT, newLocation, newRemoteLocation);
1111:     PetscMalloc1(numLeaves,    &glocalPoints);
1112:     PetscMalloc1(numLeaves, &gremotePoints);
1113:     for (l = 0; l < numLeaves; ++l) {
1114:       glocalPoints[l]        = localPoints[l]; /* localPoints[l] >= cEnd ? localPoints[l] + numGhostCells : localPoints[l]; */
1115:       gremotePoints[l].rank  = remotePoints[l].rank;
1116:       gremotePoints[l].index = newRemoteLocation[localPoints[l]];
1117:     }
1118:     PetscFree2(newLocation,newRemoteLocation);
1119:     PetscSFSetGraph(gsfPoint, numRoots+numGhostCells, numLeaves, glocalPoints, PETSC_OWN_POINTER, gremotePoints, PETSC_OWN_POINTER);
1120:   }
1121:   DMDestroy(dmSplit);
1122:   *dmSplit = sdm;
1123:   return(0);
1124: }

1128: PetscErrorCode CreatePartitionVec(DM dm, DM *dmCell, Vec *partition)
1129: {
1130:   PetscSF        sfPoint;
1131:   PetscSection   coordSection;
1132:   Vec            coordinates;
1133:   PetscSection   sectionCell;
1134:   PetscScalar    *part;
1135:   PetscInt       cStart, cEnd, c;
1136:   PetscMPIInt    rank;

1140:   DMGetCoordinateSection(dm, &coordSection);
1141:   DMGetCoordinatesLocal(dm, &coordinates);
1142:   DMClone(dm, dmCell);
1143:   DMGetPointSF(dm, &sfPoint);
1144:   DMSetPointSF(*dmCell, sfPoint);
1145:   DMSetCoordinateSection(*dmCell, PETSC_DETERMINE, coordSection);
1146:   DMSetCoordinatesLocal(*dmCell, coordinates);
1147:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
1148:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &sectionCell);
1149:   DMPlexGetHeightStratum(*dmCell, 0, &cStart, &cEnd);
1150:   PetscSectionSetChart(sectionCell, cStart, cEnd);
1151:   for (c = cStart; c < cEnd; ++c) {
1152:     PetscSectionSetDof(sectionCell, c, 1);
1153:   }
1154:   PetscSectionSetUp(sectionCell);
1155:   DMSetDefaultSection(*dmCell, sectionCell);
1156:   PetscSectionDestroy(&sectionCell);
1157:   DMCreateLocalVector(*dmCell, partition);
1158:   PetscObjectSetName((PetscObject)*partition, "partition");
1159:   VecGetArray(*partition, &part);
1160:   for (c = cStart; c < cEnd; ++c) {
1161:     PetscScalar *p;

1163:     DMPlexPointLocalRef(*dmCell, c, part, &p);
1164:     p[0] = rank;
1165:   }
1166:   VecRestoreArray(*partition, &part);
1167:   return(0);
1168: }

1172: PetscErrorCode CreateMassMatrix(DM dm, Vec *massMatrix, User user)
1173: {
1174:   DM                dmMass, dmFace, dmCell, dmCoord;
1175:   PetscSection      coordSection;
1176:   Vec               coordinates, facegeom, cellgeom;
1177:   PetscSection      sectionMass;
1178:   PetscScalar       *m;
1179:   const PetscScalar *fgeom, *cgeom, *coords;
1180:   PetscInt          vStart, vEnd, v;
1181:   PetscErrorCode    ierr;

1184:   DMGetCoordinateSection(dm, &coordSection);
1185:   DMGetCoordinatesLocal(dm, &coordinates);
1186:   DMClone(dm, &dmMass);
1187:   DMSetCoordinateSection(dmMass, PETSC_DETERMINE, coordSection);
1188:   DMSetCoordinatesLocal(dmMass, coordinates);
1189:   PetscSectionCreate(PetscObjectComm((PetscObject)dm), &sectionMass);
1190:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
1191:   PetscSectionSetChart(sectionMass, vStart, vEnd);
1192:   for (v = vStart; v < vEnd; ++v) {
1193:     PetscInt numFaces;

1195:     DMPlexGetSupportSize(dmMass, v, &numFaces);
1196:     PetscSectionSetDof(sectionMass, v, numFaces*numFaces);
1197:   }
1198:   PetscSectionSetUp(sectionMass);
1199:   DMSetDefaultSection(dmMass, sectionMass);
1200:   PetscSectionDestroy(&sectionMass);
1201:   DMGetLocalVector(dmMass, massMatrix);
1202:   VecGetArray(*massMatrix, &m);
1203:   DMPlexTSGetGeometryFVM(dm, &facegeom, &cellgeom, NULL);
1204:   VecGetDM(facegeom, &dmFace);
1205:   VecGetArrayRead(facegeom, &fgeom);
1206:   VecGetDM(cellgeom, &dmCell);
1207:   VecGetArrayRead(cellgeom, &cgeom);
1208:   DMGetCoordinateDM(dm, &dmCoord);
1209:   VecGetArrayRead(coordinates, &coords);
1210:   for (v = vStart; v < vEnd; ++v) {
1211:     const PetscInt        *faces;
1212:     const PetscFVFaceGeom *fgA, *fgB, *cg;
1213:     const PetscScalar     *vertex;
1214:     PetscInt               numFaces, sides[2], f, g;

1216:     DMPlexPointLocalRead(dmCoord, v, coords, &vertex);
1217:     DMPlexGetSupportSize(dmMass, v, &numFaces);
1218:     DMPlexGetSupport(dmMass, v, &faces);
1219:     for (f = 0; f < numFaces; ++f) {
1220:       sides[0] = faces[f];
1221:       DMPlexPointLocalRead(dmFace, faces[f], fgeom, &fgA);
1222:       for (g = 0; g < numFaces; ++g) {
1223:         const PetscInt *cells = NULL;;
1224:         PetscReal      area   = 0.0;
1225:         PetscInt       numCells;

1227:         sides[1] = faces[g];
1228:         DMPlexPointLocalRead(dmFace, faces[g], fgeom, &fgB);
1229:         DMPlexGetJoin(dmMass, 2, sides, &numCells, &cells);
1230:         if (numCells != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Invalid join for faces");
1231:         DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cg);
1232:         area += PetscAbsScalar((vertex[0] - cg->centroid[0])*(fgA->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1])*(fgA->centroid[0] - cg->centroid[0]));
1233:         area += PetscAbsScalar((vertex[0] - cg->centroid[0])*(fgB->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1])*(fgB->centroid[0] - cg->centroid[0]));
1234:         m[f*numFaces+g] = Dot2(fgA->normal, fgB->normal)*area*0.5;
1235:         DMPlexRestoreJoin(dmMass, 2, sides, &numCells, &cells);
1236:       }
1237:     }
1238:   }
1239:   VecRestoreArrayRead(facegeom, &fgeom);
1240:   VecRestoreArrayRead(cellgeom, &cgeom);
1241:   VecRestoreArrayRead(coordinates, &coords);
1242:   VecRestoreArray(*massMatrix, &m);
1243:   DMDestroy(&dmMass);
1244:   return(0);
1245: }

1249: /* Behavior will be different for multi-physics or when using non-default boundary conditions */
1250: static PetscErrorCode ModelSolutionSetDefault(Model mod,SolutionFunction func,void *ctx)
1251: {
1253:   mod->solution    = func;
1254:   mod->solutionctx = ctx;
1255:   return(0);
1256: }

1260: static PetscErrorCode ModelFunctionalRegister(Model mod,const char *name,PetscInt *offset,FunctionalFunction func,void *ctx)
1261: {
1263:   FunctionalLink link,*ptr;
1264:   PetscInt       lastoffset = -1;

1267:   for (ptr=&mod->functionalRegistry; *ptr; ptr = &(*ptr)->next) lastoffset = (*ptr)->offset;
1268:   PetscNew(&link);
1269:   PetscStrallocpy(name,&link->name);
1270:   link->offset = lastoffset + 1;
1271:   link->func   = func;
1272:   link->ctx    = ctx;
1273:   link->next   = NULL;
1274:   *ptr         = link;
1275:   *offset      = link->offset;
1276:   return(0);
1277: }

1281: static PetscErrorCode ModelFunctionalSetFromOptions(Model mod,PetscOptionItems *PetscOptionsObject)
1282: {
1284:   PetscInt       i,j;
1285:   FunctionalLink link;
1286:   char           *names[256];

1289:   mod->numMonitored = ALEN(names);
1290:   PetscOptionsStringArray("-monitor","list of functionals to monitor","",names,&mod->numMonitored,NULL);
1291:   /* Create list of functionals that will be computed somehow */
1292:   PetscMalloc1(mod->numMonitored,&mod->functionalMonitored);
1293:   /* Create index of calls that we will have to make to compute these functionals (over-allocation in general). */
1294:   PetscMalloc1(mod->numMonitored,&mod->functionalCall);
1295:   mod->numCall = 0;
1296:   for (i=0; i<mod->numMonitored; i++) {
1297:     for (link=mod->functionalRegistry; link; link=link->next) {
1298:       PetscBool match;
1299:       PetscStrcasecmp(names[i],link->name,&match);
1300:       if (match) break;
1301:     }
1302:     if (!link) SETERRQ1(mod->comm,PETSC_ERR_USER,"No known functional '%s'",names[i]);
1303:     mod->functionalMonitored[i] = link;
1304:     for (j=0; j<i; j++) {
1305:       if (mod->functionalCall[j]->func == link->func && mod->functionalCall[j]->ctx == link->ctx) goto next_name;
1306:     }
1307:     mod->functionalCall[mod->numCall++] = link; /* Just points to the first link using the result. There may be more results. */
1308: next_name:
1309:     PetscFree(names[i]);
1310:   }

1312:   /* Find out the maximum index of any functional computed by a function we will be calling (even if we are not using it) */
1313:   mod->maxComputed = -1;
1314:   for (link=mod->functionalRegistry; link; link=link->next) {
1315:     for (i=0; i<mod->numCall; i++) {
1316:       FunctionalLink call = mod->functionalCall[i];
1317:       if (link->func == call->func && link->ctx == call->ctx) {
1318:         mod->maxComputed = PetscMax(mod->maxComputed,link->offset);
1319:       }
1320:     }
1321:   }
1322:   return(0);
1323: }

1327: static PetscErrorCode FunctionalLinkDestroy(FunctionalLink *link)
1328: {
1330:   FunctionalLink l,next;

1333:   if (!link) return(0);
1334:   l     = *link;
1335:   *link = NULL;
1336:   for (; l; l=next) {
1337:     next = l->next;
1338:     PetscFree(l->name);
1339:     PetscFree(l);
1340:   }
1341:   return(0);
1342: }

1346: /* put the solution callback into a functional callback */
1347: static PetscErrorCode SolutionFunctional(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *modctx)
1348: {
1349:   Model          mod;
1352:   mod  = (Model) modctx;
1353:   (*mod->solution)(mod, time, x, u, mod->solutionctx);
1354:   return(0);
1355: }

1359: PetscErrorCode SetInitialCondition(DM dm, Vec X, User user)
1360: {
1361:   PetscErrorCode     (*func[1]) (PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
1362:   void               *ctx[1];
1363:   Model              mod = user->model;
1364:   PetscErrorCode     ierr;

1367:   func[0] = SolutionFunctional;
1368:   ctx[0]  = (void *) mod;
1369:   DMProjectFunction(dm,0.0,func,ctx,INSERT_ALL_VALUES,X);
1370:   return(0);
1371: }

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

1380:   PetscViewerCreate(PetscObjectComm((PetscObject)dm), viewer);
1381:   PetscViewerSetType(*viewer, PETSCVIEWERVTK);
1382:   PetscViewerFileSetName(*viewer, filename);
1383:   return(0);
1384: }

1388: static PetscErrorCode MonitorVTK(TS ts,PetscInt stepnum,PetscReal time,Vec X,void *ctx)
1389: {
1390:   User           user = (User)ctx;
1391:   DM             dm;
1392:   Vec            cellgeom;
1393:   PetscViewer    viewer;
1394:   char           filename[PETSC_MAX_PATH_LEN],*ftable = NULL;
1395:   PetscReal      xnorm;
1396:   PetscInt       cEndInterior;

1400:   PetscObjectSetName((PetscObject) X, "u");
1401:   VecGetDM(X,&dm);
1402:   DMPlexTSGetGeometryFVM(dm, NULL, &cellgeom, NULL);
1403:   VecNorm(X,NORM_INFINITY,&xnorm);

1405:   if (stepnum >= 0) {
1406:     stepnum += user->monitorStepOffset;
1407:   }
1408:   if (stepnum >= 0) {           /* No summary for final time */
1409:     Model             mod = user->model;
1410:     PetscInt          c,cStart,cEnd,fcount,i;
1411:     size_t            ftableused,ftablealloc;
1412:     const PetscScalar *cgeom,*x;
1413:     DM                dmCell;
1414:     DMLabel           vtkLabel;
1415:     PetscReal         *fmin,*fmax,*fintegral,*ftmp;
1416:     fcount = mod->maxComputed+1;
1417:     PetscMalloc4(fcount,&fmin,fcount,&fmax,fcount,&fintegral,fcount,&ftmp);
1418:     for (i=0; i<fcount; i++) {
1419:       fmin[i]      = PETSC_MAX_REAL;
1420:       fmax[i]      = PETSC_MIN_REAL;
1421:       fintegral[i] = 0;
1422:     }
1423:     VecGetDM(cellgeom,&dmCell);
1424:     DMPlexGetHybridBounds(dmCell, &cEndInterior, NULL, NULL, NULL);
1425:     DMPlexGetHeightStratum(dmCell,0,&cStart,&cEnd);
1426:     VecGetArrayRead(cellgeom,&cgeom);
1427:     VecGetArrayRead(X,&x);
1428:     DMGetLabel(dm,"vtk",&vtkLabel);
1429:     for (c = cStart; c < cEndInterior; ++c) {
1430:       const PetscFVCellGeom *cg;
1431:       const PetscScalar     *cx;
1432:       PetscInt              vtkVal = 0;

1434:       /* not that these two routines as currently implemented work for any dm with a
1435:        * defaultSection/defaultGlobalSection */
1436:       DMPlexPointLocalRead(dmCell,c,cgeom,&cg);
1437:       DMPlexPointGlobalRead(dm,c,x,&cx);
1438:       if (vtkLabel) {DMLabelGetValue(vtkLabel,c,&vtkVal);}
1439:       if (!vtkVal || !cx) continue;        /* ghost, or not a global cell */
1440:       for (i=0; i<mod->numCall; i++) {
1441:         FunctionalLink flink = mod->functionalCall[i];
1442:         (*flink->func)(mod,time,cg->centroid,cx,ftmp,flink->ctx);
1443:       }
1444:       for (i=0; i<fcount; i++) {
1445:         fmin[i]       = PetscMin(fmin[i],ftmp[i]);
1446:         fmax[i]       = PetscMax(fmax[i],ftmp[i]);
1447:         fintegral[i] += cg->volume * ftmp[i];
1448:       }
1449:     }
1450:     VecRestoreArrayRead(cellgeom,&cgeom);
1451:     VecRestoreArrayRead(X,&x);
1452:     MPI_Allreduce(MPI_IN_PLACE,fmin,fcount,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));
1453:     MPI_Allreduce(MPI_IN_PLACE,fmax,fcount,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));
1454:     MPI_Allreduce(MPI_IN_PLACE,fintegral,fcount,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)ts));

1456:     ftablealloc = fcount * 100;
1457:     ftableused  = 0;
1458:     PetscMalloc1(ftablealloc,&ftable);
1459:     for (i=0; i<mod->numMonitored; i++) {
1460:       size_t         countused;
1461:       char           buffer[256],*p;
1462:       FunctionalLink flink = mod->functionalMonitored[i];
1463:       PetscInt       id    = flink->offset;
1464:       if (i % 3) {
1465:         PetscMemcpy(buffer,"  ",2);
1466:         p    = buffer + 2;
1467:       } else if (i) {
1468:         char newline[] = "\n";
1469:         PetscMemcpy(buffer,newline,sizeof newline-1);
1470:         p    = buffer + sizeof newline - 1;
1471:       } else {
1472:         p = buffer;
1473:       }
1474:       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]);
1475:       countused += p - buffer;
1476:       if (countused > ftablealloc-ftableused-1) { /* reallocate */
1477:         char *ftablenew;
1478:         ftablealloc = 2*ftablealloc + countused;
1479:         PetscMalloc(ftablealloc,&ftablenew);
1480:         PetscMemcpy(ftablenew,ftable,ftableused);
1481:         PetscFree(ftable);
1482:         ftable = ftablenew;
1483:       }
1484:       PetscMemcpy(ftable+ftableused,buffer,countused);
1485:       ftableused += countused;
1486:       ftable[ftableused] = 0;
1487:     }
1488:     PetscFree4(fmin,fmax,fintegral,ftmp);

1490:     PetscPrintf(PetscObjectComm((PetscObject)ts),"% 3D  time %8.4g  |x| %8.4g  %s\n",stepnum,(double)time,(double)xnorm,ftable ? ftable : "");
1491:     PetscFree(ftable);
1492:   }
1493:   if (user->vtkInterval < 1) return(0);
1494:   if ((stepnum == -1) ^ (stepnum % user->vtkInterval == 0)) {
1495:     if (stepnum == -1) {        /* Final time is not multiple of normal time interval, write it anyway */
1496:       TSGetTimeStepNumber(ts,&stepnum);
1497:     }
1498:     PetscSNPrintf(filename,sizeof filename,"ex11-%03D.vtu",stepnum);
1499:     OutputVTK(dm,filename,&viewer);
1500:     VecView(X,viewer);
1501:     PetscViewerDestroy(&viewer);
1502:   }
1503:   return(0);
1504: }

1508: static PetscErrorCode initializeTS(DM dm, User user, TS *ts)
1509: {

1513:   TSCreate(PetscObjectComm((PetscObject)dm), ts);
1514:   TSSetType(*ts, TSSSP);
1515:   TSSetDM(*ts, dm);
1516:   TSMonitorSet(*ts,MonitorVTK,user,NULL);
1517:   DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, user);
1518:   DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, user);
1519:   TSSetDuration(*ts,1000,2.0);
1520:   TSSetExactFinalTime(*ts,TS_EXACTFINALTIME_STEPOVER);
1521:   return(0);
1522: }

1526: static PetscErrorCode adaptToleranceFVM(PetscFV fvm, TS ts, Vec sol, PetscReal refineTol, PetscReal coarsenTol, User user, TS *tsNew, Vec *solNew)
1527: {
1528:   DM                dm, gradDM, plex, cellDM, adaptedDM = NULL;
1529:   Vec               cellGeom, faceGeom;
1530:   PetscBool         isForest, computeGradient;
1531:   Vec               grad, locGrad, locX;
1532:   PetscInt          cStart, cEnd, cEndInterior, c, dim;
1533:   PetscReal         minMaxInd[2] = {PETSC_MAX_REAL, PETSC_MIN_REAL}, minMaxIndGlobal[2], minInd, maxInd, time;
1534:   const PetscScalar *pointVals;
1535:   const PetscScalar *pointGrads;
1536:   const PetscScalar *pointGeom;
1537:   DMLabel           adaptLabel = NULL;
1538:   PetscErrorCode    ierr;

1541:   TSGetTime(ts,&time);
1542:   VecGetDM(sol, &dm);
1543:   DMGetDimension(dm,&dim);
1544:   PetscFVGetComputeGradients(fvm,&computeGradient);
1545:   PetscFVSetComputeGradients(fvm,PETSC_TRUE);
1546:   DMIsForest(dm, &isForest);
1547:   DMConvert(dm, DMPLEX, &plex);
1548:   DMPlexGetDataFVM(plex, fvm, &cellGeom, &faceGeom, &gradDM);
1549:   DMCreateLocalVector(plex,&locX);
1550:   DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, 0.0, faceGeom, cellGeom, NULL);
1551:   DMGlobalToLocalBegin(plex, sol, INSERT_VALUES, locX);
1552:   DMGlobalToLocalEnd  (plex, sol, INSERT_VALUES, locX);
1553:   DMCreateGlobalVector(gradDM, &grad);
1554:   DMPlexReconstructGradientsFVM(plex, locX, grad);
1555:   DMCreateLocalVector(gradDM, &locGrad);
1556:   DMGlobalToLocalBegin(gradDM, grad, INSERT_VALUES, locGrad);
1557:   DMGlobalToLocalEnd(gradDM, grad, INSERT_VALUES, locGrad);
1558:   VecDestroy(&grad);
1559:   DMPlexGetHeightStratum(plex,0,&cStart,&cEnd);
1560:   DMPlexGetHybridBounds(plex,&cEndInterior,NULL,NULL,NULL);
1561:   cEnd = (cEndInterior < 0) ? cEnd : cEndInterior;
1562:   VecGetArrayRead(locGrad,&pointGrads);
1563:   VecGetArrayRead(cellGeom,&pointGeom);
1564:   VecGetArrayRead(locX,&pointVals);
1565:   VecGetDM(cellGeom,&cellDM);
1566:   DMLabelCreate("adapt",&adaptLabel);
1567:   for (c = cStart; c < cEnd; c++) {
1568:     PetscReal             errInd = 0.;
1569:     const PetscScalar     *pointGrad;
1570:     const PetscScalar     *pointVal;
1571:     const PetscFVCellGeom *cg;

1573:     DMPlexPointLocalRead(gradDM,c,pointGrads,&pointGrad);
1574:     DMPlexPointLocalRead(cellDM,c,pointGeom,&cg);
1575:     DMPlexPointLocalRead(plex,c,pointVals,&pointVal);

1577:     (user->model->errorIndicator)(dim,cg->volume,user->model->physics->dof,pointVal,pointGrad,&errInd,user->model->errorCtx);
1578:     minMaxInd[0] = PetscMin(minMaxInd[0],errInd);
1579:     minMaxInd[1] = PetscMax(minMaxInd[1],errInd);
1580:     if (errInd > refineTol)  {DMLabelSetValue(adaptLabel,c,DM_ADAPT_REFINE);}
1581:     if (errInd < coarsenTol) {DMLabelSetValue(adaptLabel,c,DM_ADAPT_COARSEN);}
1582:   }
1583:   VecRestoreArrayRead(locX,&pointVals);
1584:   VecRestoreArrayRead(cellGeom,&pointGeom);
1585:   VecRestoreArrayRead(locGrad,&pointGrads);
1586:   VecDestroy(&locGrad);
1587:   VecDestroy(&locX);
1588:   DMDestroy(&plex);
1589:   PetscFVSetComputeGradients(fvm,computeGradient);
1590:   minMaxInd[1] = -minMaxInd[1];
1591:   MPI_Allreduce(minMaxInd,minMaxIndGlobal,2,MPIU_REAL,MPI_MIN,PetscObjectComm((PetscObject)dm));
1592:   minInd = minMaxIndGlobal[0];
1593:   maxInd = -minMaxIndGlobal[1];
1594:   PetscInfo2(ts, "error indicator range (%E, %E)\n", minInd, maxInd);
1595:   if (maxInd > refineTol || minInd < coarsenTol) { /* at least one cell is over the refinement threshold */
1596:     DMAdaptLabel(dm,adaptLabel,&adaptedDM);
1597:   }
1598:   else if (maxInd < coarsenTol) { /* all cells are under the coarsening threshold */
1599:     DMCoarsen(dm,PetscObjectComm((PetscObject)dm),&adaptedDM);
1600:   }
1601:   else if (minInd > refineTol) { /* all cells are over the refinement threshold */
1602:     DMRefine(dm,PetscObjectComm((PetscObject)dm),&adaptedDM);
1603:   }
1604:   DMLabelDestroy(&adaptLabel);
1605:   if (adaptedDM) {
1606:     if (tsNew) {
1607:       initializeTS(adaptedDM,user,tsNew);
1608:     }
1609:     if (solNew) {
1610:       DMCreateGlobalVector(adaptedDM, solNew);
1611:       PetscObjectSetName((PetscObject) *solNew, "solution");
1612:       DMForestTransferVec(dm, sol, adaptedDM, *solNew, PETSC_TRUE, time);
1613:     }
1614:     if (isForest) {DMForestSetAdaptivityForest(adaptedDM,NULL);} /* clear internal references to the previous dm */
1615:     DMDestroy(&adaptedDM);
1616:   }
1617:   else {
1618:     if (tsNew) *tsNew = NULL;
1619:     if (solNew) *solNew = NULL;
1620:   }
1621:   return(0);
1622: }

1626: int main(int argc, char **argv)
1627: {
1628:   MPI_Comm          comm;
1629:   PetscDS           prob;
1630:   PetscFV           fvm;
1631:   PetscLimiter      limiter = NULL, noneLimiter = NULL;
1632:   User              user;
1633:   Model             mod;
1634:   Physics           phys;
1635:   DM                dm;
1636:   PetscReal         ftime, cfl, dt, minRadius;
1637:   PetscInt          dim, nsteps;
1638:   TS                ts;
1639:   TSConvergedReason reason;
1640:   Vec               X;
1641:   PetscViewer       viewer;
1642:   PetscBool         vtkCellGeom, splitFaces, useAMR, viewInitial;
1643:   PetscReal         refineTol, coarsenTol;
1644:   PetscInt          overlap, targetCells, adaptInterval;
1645:   char              filename[PETSC_MAX_PATH_LEN] = "sevenside.exo";
1646:   char              physname[256]  = "advect";
1647:   PetscErrorCode    ierr;

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

1652:   PetscNew(&user);
1653:   PetscNew(&user->model);
1654:   PetscNew(&user->model->physics);
1655:   mod           = user->model;
1656:   phys          = mod->physics;
1657:   mod->comm     = comm;
1658:   useAMR        = PETSC_FALSE;
1659:   viewInitial   = PETSC_FALSE;
1660:   refineTol     = PETSC_MAX_REAL;
1661:   coarsenTol    = 0.;
1662:   targetCells   = 0;
1663:   adaptInterval = 1;

1665:   /* Register physical models to be available on the command line */
1666:   PetscFunctionListAdd(&PhysicsList,"advect"          ,PhysicsCreate_Advect);
1667:   PetscFunctionListAdd(&PhysicsList,"sw"              ,PhysicsCreate_SW);
1668:   PetscFunctionListAdd(&PhysicsList,"euler"           ,PhysicsCreate_Euler);

1670:   PetscOptionsBegin(comm,NULL,"Unstructured Finite Volume Mesh Options","");
1671:   {
1672:     cfl  = 0.9 * 4; /* default SSPRKS2 with s=5 stages is stable for CFL number s-1 */
1673:     PetscOptionsReal("-ufv_cfl","CFL number per step","",cfl,&cfl,NULL);
1674:     PetscOptionsString("-f","Exodus.II filename to read","",filename,filename,sizeof(filename),NULL);
1675:     splitFaces = PETSC_FALSE;
1676:     PetscOptionsBool("-ufv_split_faces","Split faces between cell sets","",splitFaces,&splitFaces,NULL);
1677:     overlap = 1;
1678:     PetscOptionsInt("-ufv_mesh_overlap","Number of cells to overlap partitions","",overlap,&overlap,NULL);
1679:     user->vtkInterval = 1;
1680:     PetscOptionsInt("-ufv_vtk_interval","VTK output interval (0 to disable)","",user->vtkInterval,&user->vtkInterval,NULL);
1681:     vtkCellGeom = PETSC_FALSE;
1682:     PetscOptionsBool("-ufv_vtk_cellgeom","Write cell geometry (for debugging)","",vtkCellGeom,&vtkCellGeom,NULL);
1683:     PetscOptionsBool("-ufv_use_amr","use local adaptive mesh refinement","",useAMR,&useAMR,NULL);
1684:     PetscOptionsReal("-ufv_refine_tol","tolerance for refining cells in AMR","",refineTol,&refineTol,NULL);
1685:     PetscOptionsReal("-ufv_coarsen_tol","tolerance for coarsening cells in AMR","",coarsenTol,&coarsenTol,NULL);
1686:     PetscOptionsInt("-ufv_target_cells","target number of cells in AMR","",targetCells,&targetCells,NULL);
1687:     PetscOptionsBool("-ufv_view_initial_refinement","View initial conditions refinement history","",viewInitial,&viewInitial,NULL);
1688:     PetscOptionsInt("-ufv_adapt_interval","time steps between AMR","",adaptInterval,&adaptInterval,NULL);
1689:   }
1690:   PetscOptionsEnd();

1692:   PetscOptionsBegin(comm,NULL,"Unstructured Finite Volume Physics Options","");
1693:   {
1694:     PetscErrorCode (*physcreate)(Model,Physics,PetscOptionItems*);
1695:     PetscOptionsFList("-physics","Physics module to solve","",PhysicsList,physname,physname,sizeof physname,NULL);
1696:     PetscFunctionListFind(PhysicsList,physname,&physcreate);
1697:     PetscMemzero(phys,sizeof(struct _n_Physics));
1698:     (*physcreate)(mod,phys,PetscOptionsObject);
1699:     /* Count number of fields and dofs */
1700:     for (phys->nfields=0,phys->dof=0; phys->field_desc[phys->nfields].name; phys->nfields++) phys->dof += phys->field_desc[phys->nfields].dof;
1701:     if (phys->dof <= 0) SETERRQ1(comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set dof",physname);
1702:     ModelFunctionalSetFromOptions(mod,PetscOptionsObject);
1703:   }
1704:   PetscOptionsEnd();

1706:   {
1707:     size_t len,i;
1708:     for (i = 0; i < DIM; i++) { mod->bounds[2*i] = 0.; mod->bounds[2*i+1] = 1.;};
1709:     PetscStrlen(filename,&len);
1710:     dim = DIM;
1711:     if (!len) { /* a null name means just do a hex box */
1712:       PetscInt cells[3] = {1, 1, 1}; /* coarse mesh is one cell; refine from there */
1713:       PetscBool flg1, flg2, skew = PETSC_FALSE;
1714:       PetscInt nret1 = DIM;
1715:       PetscInt nret2 = 2*DIM;
1716:       PetscOptionsBegin(comm,NULL,"Rectangular mesh options","");
1717:       PetscOptionsIntArray("-grid_size","number of cells in each direction","",cells,&nret1,&flg1);
1718:       PetscOptionsRealArray("-grid_bounds","bounds of the mesh in each direction (i.e., x_min,x_max,y_min,y_max","",mod->bounds,&nret2,&flg2);
1719:       PetscOptionsBool("-grid_skew_60","Skew grid for 60 degree shock mesh","",skew,&skew,NULL);
1720:       PetscOptionsEnd();
1721:       if (flg1) {
1722:         dim = nret1;
1723:         if (dim != DIM) SETERRQ1(comm,PETSC_ERR_ARG_SIZ,"Dim wrong size %D in -grid_size",dim);
1724:       }
1725:       DMPlexCreateHexBoxMesh(comm, dim, cells, mod->bcs[0], mod->bcs[1], mod->bcs[2], &dm);
1726:       if (flg2) {
1727:         PetscInt dimEmbed, i;
1728:         PetscInt nCoords;
1729:         PetscScalar *coords;
1730:         Vec coordinates;

1732:         DMGetCoordinatesLocal(dm,&coordinates);
1733:         DMGetCoordinateDim(dm,&dimEmbed);
1734:         VecGetLocalSize(coordinates,&nCoords);
1735:         if (nCoords % dimEmbed) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Coordinate vector the wrong size");
1736:         VecGetArray(coordinates,&coords);
1737:         for (i = 0; i < nCoords; i += dimEmbed) {
1738:           PetscInt j;

1740:           PetscScalar *coord = &coords[i];
1741:           for (j = 0; j < dimEmbed; j++) {
1742:             coord[j] = mod->bounds[2 * j] + coord[j] * (mod->bounds[2 * j + 1] - mod->bounds[2 * j]);
1743:             if (dim==2 && cells[1]==1 && j==0 && skew) {
1744:               if (cells[0]==2 && i==8) {
1745:                 coord[j] = .57735026918963; /* hack to get 60 deg skewed mesh */
1746:               }
1747:               else if (cells[0]==3) {
1748:                 if(i==2 || i==10) coord[j] = mod->bounds[1]/4.;
1749:                 else if (i==4) coord[j] = mod->bounds[1]/2.;
1750:                 else if (i==12) coord[j] = 1.57735026918963*mod->bounds[1]/2.;
1751:               }
1752:             }
1753:           }
1754:         }
1755:         VecRestoreArray(coordinates,&coords);
1756:         DMSetCoordinatesLocal(dm,coordinates);
1757:       }
1758:     }
1759:     else {
1760:       DMPlexCreateFromFile(comm, filename, PETSC_TRUE, &dm);
1761:     }
1762:   }
1763:   DMViewFromOptions(dm, NULL, "-dm_view");
1764:   DMGetDimension(dm, &dim);

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

1769:   mod->errorIndicator = ErrorIndicator_Simple;

1771:   {
1772:     DM dmDist;

1774:     DMPlexSetAdjacencyUseCone(dm, PETSC_TRUE);
1775:     DMPlexSetAdjacencyUseClosure(dm, PETSC_FALSE);
1776:     DMPlexDistribute(dm, overlap, NULL, &dmDist);
1777:     if (dmDist) {
1778:       DMDestroy(&dm);
1779:       dm   = dmDist;
1780:     }
1781:   }

1783:   DMSetFromOptions(dm);

1785:   {
1786:     DM gdm;

1788:     DMPlexConstructGhostCells(dm, NULL, NULL, &gdm);
1789:     DMDestroy(&dm);
1790:     dm   = gdm;
1791:     DMViewFromOptions(dm, NULL, "-dm_view");
1792:   }
1793:   if (splitFaces) {ConstructCellBoundary(dm, user);}
1794:   SplitFaces(&dm, "split faces", user);

1796:   PetscDSCreate(PetscObjectComm((PetscObject)dm),&prob);
1797:   PetscFVCreate(comm, &fvm);
1798:   PetscFVSetFromOptions(fvm);
1799:   PetscFVSetNumComponents(fvm, phys->dof);
1800:   PetscFVSetSpatialDimension(fvm, dim);
1801:   PetscObjectSetName((PetscObject) fvm,"");
1802:   {
1803:     PetscInt f, dof;
1804:     for (f=0,dof=0; f < phys->nfields; f++) {
1805:       PetscInt newDof = phys->field_desc[f].dof;

1807:       if (newDof == 1) {
1808:         PetscFVSetComponentName(fvm,dof,phys->field_desc[f].name);
1809:       }
1810:       else {
1811:         PetscInt j;

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

1816:           PetscSNPrintf(compName,sizeof(compName),"%s_%d",phys->field_desc[f].name,j);
1817:           PetscFVSetComponentName(fvm,dof+j,compName);
1818:         }
1819:       }
1820:       dof += newDof;
1821:     }
1822:   }
1823:   /* FV is now structured with one field having all physics as components */
1824:   PetscDSAddDiscretization(prob, (PetscObject) fvm);
1825:   PetscDSSetRiemannSolver(prob, 0, user->model->physics->riemann);
1826:   PetscDSSetContext(prob, 0, user->model->physics);
1827:   (*mod->setupbc)(prob,phys);
1828:   PetscDSSetFromOptions(prob);
1829:   DMSetDS(dm,prob);
1830:   PetscDSDestroy(&prob);
1831:   {
1832:     char      convType[256];
1833:     PetscBool flg;

1835:     PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");
1836:     PetscOptionsFList("-dm_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg);
1837:     PetscOptionsEnd();
1838:     if (flg) {
1839:       DM dmConv;

1841:       DMConvert(dm,convType,&dmConv);
1842:       if (dmConv) {
1843:         DMViewFromOptions(dmConv, NULL, "-dm_conv_view");
1844:         DMDestroy(&dm);
1845:         dm   = dmConv;
1846:         DMSetFromOptions(dm);
1847:       }
1848:     }
1849:   }

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

1853:   DMCreateGlobalVector(dm, &X);
1854:   PetscObjectSetName((PetscObject) X, "solution");
1855:   SetInitialCondition(dm, X, user);
1856:   if (useAMR) {
1857:     /* use no limiting when reconstructing gradients for adaptivity */
1858:     PetscFVGetLimiter(fvm, &limiter);
1859:     PetscObjectReference((PetscObject)limiter);

1861:     PetscLimiterCreate(PetscObjectComm((PetscObject)fvm),&noneLimiter);
1862:     PetscLimiterSetType(noneLimiter,PETSCLIMITERNONE);

1864:     PetscFVSetLimiter(fvm,noneLimiter);
1865:     if (targetCells > 0) {
1866:       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"AMR with target cell count not implemented yet.");
1867:     } else {
1868:       PetscInt adaptIter;

1870:       for (adaptIter = 0;;adaptIter++) {
1871:         PetscLogDouble bytes;
1872:         TS             tsNew = NULL;

1874:         if (viewInitial) {
1875:           PetscViewer viewer;
1876:           char        buf[256];
1877:           PetscBool   isHDF5, isVTK;

1879:           PetscViewerCreate(comm,&viewer);
1880:           PetscViewerSetType(viewer,PETSCVIEWERVTK);
1881:           PetscViewerSetOptionsPrefix(viewer,"initial_");
1882:           PetscViewerSetFromOptions(viewer);
1883:           PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&isHDF5);
1884:           PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isVTK);
1885:           if (isHDF5) {
1886:             PetscSNPrintf(buf, 256, "ex11-initial-%d.h5", adaptIter);
1887:           } else if (isVTK) {
1888:             PetscSNPrintf(buf, 256, "ex11-initial-%d.vtu", adaptIter);
1889:             PetscViewerPushFormat(viewer,PETSC_VIEWER_VTK_VTU);
1890:           }
1891:           PetscViewerFileSetMode(viewer,FILE_MODE_WRITE);
1892:           PetscViewerFileSetName(viewer,buf);
1893:           if (isHDF5) {
1894:             DMView(dm,viewer);
1895:             PetscViewerFileSetMode(viewer,FILE_MODE_UPDATE);
1896:           }
1897:           VecView(X,viewer);
1898:           PetscViewerDestroy(&viewer);
1899:         }

1901:         PetscMemoryGetCurrentUsage(&bytes);
1902:         PetscInfo2(ts, "refinement loop %D: memory used %g\n", adaptIter, bytes);
1903:         adaptToleranceFVM(fvm, ts, X, refineTol, coarsenTol, user, &tsNew, NULL);
1904:         if (!tsNew) {
1905:           break;
1906:         } else {
1907:           DMDestroy(&dm);
1908:           VecDestroy(&X);
1909:           TSDestroy(&ts);
1910:           ts   = tsNew;
1911:           TSGetDM(ts,&dm);
1912:           PetscObjectReference((PetscObject)dm);
1913:           DMCreateGlobalVector(dm,&X);
1914:           PetscObjectSetName((PetscObject) X, "solution");
1915:           SetInitialCondition(dm, X, user);
1916:         }
1917:       }
1918:     }
1919:     /* restore original limiter */
1920:     PetscFVSetLimiter(fvm,limiter);
1921:   }

1923:   if (vtkCellGeom) {
1924:     DM  dmCell;
1925:     Vec cellgeom, partition;

1927:     DMPlexTSGetGeometryFVM(dm, NULL, &cellgeom, NULL);
1928:     OutputVTK(dm, "ex11-cellgeom.vtk", &viewer);
1929:     VecView(cellgeom, viewer);
1930:     PetscViewerDestroy(&viewer);
1931:     CreatePartitionVec(dm, &dmCell, &partition);
1932:     OutputVTK(dmCell, "ex11-partition.vtk", &viewer);
1933:     VecView(partition, viewer);
1934:     PetscViewerDestroy(&viewer);
1935:     VecDestroy(&partition);
1936:     DMDestroy(&dmCell);
1937:   }

1939:   /* collect max maxspeed from all processes -- todo */
1940:   DMPlexTSGetGeometryFVM(dm, NULL, NULL, &minRadius);
1941:   MPI_Allreduce(&phys->maxspeed,&mod->maxspeed,1,MPIU_REAL,MPI_MAX,PetscObjectComm((PetscObject)ts));
1942:   if (mod->maxspeed <= 0) SETERRQ1(comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set maxspeed",physname);
1943:   dt   = cfl * minRadius / mod->maxspeed;
1944:   TSSetInitialTimeStep(ts,0.0,dt);
1945:   TSSetFromOptions(ts);
1946:   if (!useAMR) {
1947:     TSSolve(ts,X);
1948:     TSGetSolveTime(ts,&ftime);
1949:     TSGetTimeStepNumber(ts,&nsteps);
1950:   } else {
1951:     PetscReal finalTime;
1952:     PetscInt  adaptIter;
1953:     TS        tsNew = NULL;
1954:     Vec       solNew = NULL;
1955:     PetscInt  incSteps;

1957:     TSGetDuration(ts,NULL,&finalTime);
1958:     TSSetDuration(ts,adaptInterval,finalTime);
1959:     TSSolve(ts,X);
1960:     TSGetSolveTime(ts,&ftime);
1961:     TSGetTimeStepNumber(ts,&nsteps);
1962:     for (adaptIter = 0;ftime < finalTime;adaptIter++) {
1963:       PetscLogDouble bytes;

1965:       PetscMemoryGetCurrentUsage(&bytes);
1966:       PetscInfo2(ts, "AMR time step loop %D: memory used %g\n", adaptIter, bytes);
1967:       PetscFVSetLimiter(fvm,noneLimiter);
1968:       if (targetCells > 0) {
1969:         SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"AMR with target cell count not implemented yet.");
1970:       }
1971:       else {
1972:         adaptToleranceFVM(fvm,ts,X,refineTol,coarsenTol,user,&tsNew,&solNew);
1973:       }
1974:       PetscFVSetLimiter(fvm,limiter);
1975:       if (tsNew) {
1976:         PetscInfo(ts, "AMR used\n");
1977:         DMDestroy(&dm);
1978:         VecDestroy(&X);
1979:         TSDestroy(&ts);
1980:         ts   = tsNew;
1981:         X    = solNew;
1982:         TSSetFromOptions(ts);
1983:         VecGetDM(X,&dm);
1984:         PetscObjectReference((PetscObject)dm);
1985:         DMPlexTSGetGeometryFVM(dm, NULL, NULL, &minRadius);
1986:         MPI_Allreduce(&phys->maxspeed,&mod->maxspeed,1,MPIU_REAL,MPI_MAX,PetscObjectComm((PetscObject)ts));
1987:         if (mod->maxspeed <= 0) SETERRQ1(comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set maxspeed",physname);
1988:         dt   = cfl * minRadius / mod->maxspeed;
1989:         TSSetInitialTimeStep(ts,ftime,dt);
1990:       }
1991:       else {
1992:         PetscInfo(ts, "AMR not used\n");
1993:       }
1994:       user->monitorStepOffset = nsteps;
1995:       TSSetDuration(ts,adaptInterval,finalTime);
1996:       TSSolve(ts,X);
1997:       TSGetSolveTime(ts,&ftime);
1998:       TSGetTimeStepNumber(ts,&incSteps);
1999:       nsteps += incSteps;
2000:     }
2001:   }
2002:   TSGetConvergedReason(ts,&reason);
2003:   PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,nsteps);
2004:   TSDestroy(&ts);

2006:   PetscFunctionListDestroy(&PhysicsList);
2007:   FunctionalLinkDestroy(&user->model->functionalRegistry);
2008:   PetscFree(user->model->functionalMonitored);
2009:   PetscFree(user->model->functionalCall);
2010:   PetscFree(user->model->physics->data);
2011:   PetscFree(user->model->physics);
2012:   PetscFree(user->model);
2013:   PetscFree(user);
2014:   VecDestroy(&X);
2015:   PetscLimiterDestroy(&limiter);
2016:   PetscLimiterDestroy(&noneLimiter);
2017:   PetscFVDestroy(&fvm);
2018:   DMDestroy(&dm);
2019:   PetscFinalize();
2020:   return ierr;
2021: }

2023: /* Godunov fluxs */
2024: PetscScalar cvmgp_(PetscScalar *a, PetscScalar *b, PetscScalar *test)
2025: {
2026:     /* System generated locals */
2027:     PetscScalar ret_val;

2029:     if (*test > 0.) {
2030:         goto L10;
2031:     }
2032:     ret_val = *b;
2033:     return ret_val;
2034: L10:
2035:     ret_val = *a;
2036:     return ret_val;
2037: } /* cvmgp_ */

2039: PetscScalar cvmgm_(PetscScalar *a, PetscScalar *b, PetscScalar *test)
2040: {
2041:     /* System generated locals */
2042:     PetscScalar ret_val;

2044:     if (*test < 0.) {
2045:         goto L10;
2046:     }
2047:     ret_val = *b;
2048:     return ret_val;
2049: L10:
2050:     ret_val = *a;
2051:     return ret_val;
2052: } /* cvmgm_ */

2056: int riem1mdt( PetscScalar *gaml, PetscScalar *gamr, PetscScalar *rl, PetscScalar *pl,
2057:               PetscScalar *uxl, PetscScalar *rr, PetscScalar *pr,
2058:               PetscScalar *uxr, PetscScalar *rstarl, PetscScalar *rstarr, PetscScalar *
2059:               pstar, PetscScalar *ustar)
2060: {
2061:     /* Initialized data */

2063:     static PetscScalar smallp = 1e-8;

2065:     /* System generated locals */
2066:     int i__1;
2067:     PetscScalar d__1, d__2;

2069:     /* Local variables */
2070:     static int i0;
2071:     static PetscScalar cl, cr, wl, zl, wr, zr, pst, durl, skpr1, skpr2;
2072:     static int iwave;
2073:     static PetscScalar csqrl, csqrr, gascl1, gascl2, gascl3, gascl4, gascr1, gascr2, gascr3, gascr4, cstarl, dpstar, cstarr;
2074:     static int iterno;
2075:     static PetscScalar ustarl, ustarr, rarepr1, rarepr2;



2079:     gascl1 = *gaml - 1.;
2080:     gascl2 = (*gaml + 1.) * .5;
2081:     gascl3 = gascl2 / *gaml;
2082:     gascl4 = 1.f / (*gaml - 1.);

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

2210:     *ustar = (zl * ustarr + zr * ustarl) / (zl + zr);
2211:     if (*pstar > *pl) {
2212:         pst = *pstar / *pl;
2213:         *rstarl = ((*gaml + 1.) * pst + *gaml - 1.) / ((*gaml - 1.) * pst + *
2214:                 gaml + 1.) * *rl;
2215:     }
2216:     if (*pstar > *pr) {
2217:         pst = *pstar / *pr;
2218:         *rstarr = ((*gamr + 1.) * pst + *gamr - 1.) / ((*gamr - 1.) * pst + *
2219:                 gamr + 1.) * *rr;
2220:     }
2221:     return iwave;
2222: }

2224: PetscScalar sign(PetscScalar x)
2225: {
2226:     if(x > 0) return 1.0;
2227:     if(x < 0) return -1.0;
2228:     return 0.0;
2229: }
2230: /*        Riemann Solver */
2233: /* -------------------------------------------------------------------- */
2234: int riemannsolver(PetscScalar *xcen, PetscScalar *xp,
2235:                    PetscScalar *dtt, PetscScalar *rl, PetscScalar *uxl, PetscScalar *pl,
2236:                    PetscScalar *utl, PetscScalar *ubl, PetscScalar *gaml, PetscScalar *rho1l,
2237:                    PetscScalar *rr, PetscScalar *uxr, PetscScalar *pr, PetscScalar *utr,
2238:                    PetscScalar *ubr, PetscScalar *gamr, PetscScalar *rho1r, PetscScalar *rx,
2239:                    PetscScalar *uxm, PetscScalar *px, PetscScalar *utx, PetscScalar *ubx,
2240:                    PetscScalar *gam, PetscScalar *rho1)
2241: {
2242:     /* System generated locals */
2243:     PetscScalar d__1, d__2;

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

2249:     if (*rl == *rr && *pr == *pl && *uxl == *uxr && *gaml == *gamr) {
2250:         *rx = *rl;
2251:         *px = *pl;
2252:         *uxm = *uxl;
2253:         *gam = *gaml;
2254:         x2 = *xcen + *uxm * *dtt;

2256:         if (*xp >= x2) {
2257:             *utx = *utr;
2258:             *ubx = *ubr;
2259:             *rho1 = *rho1r;
2260:         } else {
2261:             *utx = *utl;
2262:             *ubx = *ubl;
2263:             *rho1 = *rho1l;
2264:         }
2265:         return 0;
2266:     }
2267:     int iwave = riem1mdt(gaml, gamr, rl, pl, uxl, rr, pr, uxr, &rstarl, &rstarr, &pstar, &ustar);

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

2341:     /* Local variables */
2342:     static int k;
2343:     static PetscScalar bn[3], fn, ft, tg[3], pl, rl, pm, pr, rr, xp, ubl, ubm,
2344:             ubr, dtt, unm, tmp, utl, utm, uxl, utr, uxr, gaml, gamm, gamr,
2345:             xcen, rhom, rho1l, rho1m, rho1r;
2346:     /* Parameter adjustments */
2347:     --nn;
2348:     --flux;
2349:     --ur;
2350:     --ul;

2352:     /* Function Body */
2353:     xcen = 0.;
2354:     xp = 0.;
2355:     i__1 = *ndim;
2356:     for (k = 1; k <= i__1; ++k) {
2357:         tg[k - 1] = 0.;
2358:         bn[k - 1] = 0.;
2359:     }
2360:     dtt = 1.;
2361:     if (*ndim == 3) {
2362:         if (nn[1] == 0. && nn[2] == 0.) {
2363:             tg[0] = 1.;
2364:         } else {
2365:             tg[0] = -nn[2];
2366:             tg[1] = nn[1];
2367:         }
2368: /*           tmp=dsqrt(tg(1)**2+tg(2)**2) */
2369: /*           tg=tg/tmp */
2370:         bn[0] = -nn[3] * tg[1];
2371:         bn[1] = nn[3] * tg[0];
2372:         bn[2] = nn[1] * tg[1] - nn[2] * tg[0];
2373: /* Computing 2nd power */
2374:         d__1 = bn[0];
2375: /* Computing 2nd power */
2376:         d__2 = bn[1];
2377: /* Computing 2nd power */
2378:         d__3 = bn[2];
2379:         tmp = sqrt(d__1 * d__1 + d__2 * d__2 + d__3 * d__3);
2380:         i__1 = *ndim;
2381:         for (k = 1; k <= i__1; ++k) {
2382:             bn[k - 1] /= tmp;
2383:         }
2384:     } else if (*ndim == 2) {
2385:         tg[0] = -nn[2];
2386:         tg[1] = nn[1];
2387: /*           tmp=dsqrt(tg(1)**2+tg(2)**2) */
2388: /*           tg=tg/tmp */
2389:         bn[0] = 0.;
2390:         bn[1] = 0.;
2391:         bn[2] = 1.;
2392:     }
2393:     rl = ul[1];
2394:     rr = ur[1];
2395:     uxl = 0.;
2396:     uxr = 0.;
2397:     utl = 0.;
2398:     utr = 0.;
2399:     ubl = 0.;
2400:     ubr = 0.;
2401:     i__1 = *ndim;
2402:     for (k = 1; k <= i__1; ++k) {
2403:         uxl += ul[k + 1] * nn[k];
2404:         uxr += ur[k + 1] * nn[k];
2405:         utl += ul[k + 1] * tg[k - 1];
2406:         utr += ur[k + 1] * tg[k - 1];
2407:         ubl += ul[k + 1] * bn[k - 1];
2408:         ubr += ur[k + 1] * bn[k - 1];
2409:     }
2410:     uxl /= rl;
2411:     uxr /= rr;
2412:     utl /= rl;
2413:     utr /= rr;
2414:     ubl /= rl;
2415:     ubr /= rr;

2417:     gaml = *gamma;
2418:     gamr = *gamma;
2419: /* Computing 2nd power */
2420:     d__1 = uxl;
2421: /* Computing 2nd power */
2422:     d__2 = utl;
2423: /* Computing 2nd power */
2424:     d__3 = ubl;
2425:     pl = (*gamma - 1.) * (ul[*ndim + 2] - rl * .5 * (d__1 * d__1 + d__2 * d__2 + d__3 * d__3));
2426: /* Computing 2nd power */
2427:     d__1 = uxr;
2428: /* Computing 2nd power */
2429:     d__2 = utr;
2430: /* Computing 2nd power */
2431:     d__3 = ubr;
2432:     pr = (*gamma - 1.) * (ur[*ndim + 2] - rr * .5 * (d__1 * d__1 + d__2 * d__2 + d__3 * d__3));
2433:     rho1l = rl;
2434:     rho1r = rr;

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

2440:     flux[1] = rhom * unm;
2441:     fn = rhom * unm * unm + pm;
2442:     ft = rhom * unm * utm;
2443: /*           flux(2)=fn*nn(1)+ft*nn(2) */
2444: /*           flux(3)=fn*tg(1)+ft*tg(2) */
2445:     flux[2] = fn * nn[1] + ft * tg[0];
2446:     flux[3] = fn * nn[2] + ft * tg[1];
2447: /*           flux(2)=rhom*unm*(unm)+pm */
2448: /*           flux(3)=rhom*(unm)*utm */
2449:     if (*ndim == 3) {
2450:         flux[4] = rhom * unm * ubm;
2451:     }
2452:     flux[*ndim + 2] = (rhom * .5 * (unm * unm + utm * utm + ubm * ubm) + gamm / (gamm - 1.) * pm) * unm;
2453:     return iwave;
2454: } /* godunovflux_ */

2456: /* Subroutine to set up the initial conditions for the */
2457: /* Shock Interface interaction or linear wave (Ravi Samtaney,Mark Adams). */
2458: /* ----------------------------------------------------------------------- */
2459: int projecteqstate(PetscScalar wc[], const PetscScalar ueq[], const PetscScalar lv[][3])
2460: {
2461:   int j,k;
2462: /*      Wc=matmul(lv,Ueq) 3 vars */
2463:   for (k = 0; k < 3; ++k) {
2464:     wc[k] = 0.;
2465:     for (j = 0; j < 3; ++j) {
2466:       wc[k] += lv[k][j]*ueq[j];
2467:     }
2468:   }
2469:   return 0;
2470: }
2471: /* ----------------------------------------------------------------------- */
2472: int projecttoprim(PetscScalar v[], const PetscScalar wc[], PetscScalar rv[][3])
2473: {
2474:   int k,j;
2475:   /*      V=matmul(rv,WC) 3 vars */
2476:   for (k = 0; k < 3; ++k) {
2477:     v[k] = 0.;
2478:     for (j = 0; j < 3; ++j) {
2479:       v[k] += rv[k][j]*wc[j];
2480:     }
2481:   }
2482:   return 0;
2483: }
2484: /* ---------------------------------------------------------------------- */
2485: int eigenvectors(PetscScalar rv[][3], PetscScalar lv[][3], const PetscScalar ueq[], PetscScalar gamma)
2486: {
2487:   int j,k;
2488:   PetscScalar rho,csnd,p0,u;

2490:   for (k = 0; k < 3; ++k) for (j = 0; j < 3; ++j) { lv[k][j] = 0.; rv[k][j] = 0.; }
2491:   rho = ueq[0];
2492:   u = ueq[1];
2493:   p0 = ueq[2];
2494:   csnd = PetscSqrtScalar(gamma * p0 / rho);
2495:   lv[0][1] = rho * .5;
2496:   lv[0][2] = -.5 / csnd;
2497:   lv[1][0] = csnd;
2498:   lv[1][2] = -1. / csnd;
2499:   lv[2][1] = rho * .5;
2500:   lv[2][2] = .5 / csnd;
2501:   rv[0][0] = -1. / csnd;
2502:   rv[1][0] = 1. / rho;
2503:   rv[2][0] = -csnd;
2504:   rv[0][1] = 1. / csnd;
2505:   rv[0][2] = 1. / csnd;
2506:   rv[1][2] = 1. / rho;
2507:   rv[2][2] = csnd;
2508:   return 0;
2509: }
2512: int initLinearWave(EulerNode *ux, const PetscScalar gamma, const PetscReal coord[], const PetscReal Lx)
2513: {
2514:   PetscScalar p0,u0,wcp[3],wc[3];
2515:   PetscScalar lv[3][3];
2516:   PetscScalar vp[3];
2517:   PetscScalar rv[3][3];
2518:   PetscScalar eps, ueq[3], rho0, twopi;

2520:   /* Function Body */
2521:   twopi = 2.*M_PI;
2522:   eps = 1e-4; /* perturbation */
2523:   rho0 = 1e3;   /* density of water */
2524:   p0 = 101325.; /* init pressure of 1 atm (?) */
2525:   u0 = 0.;
2526:   ueq[0] = rho0;
2527:   ueq[1] = u0;
2528:   ueq[2] = p0;
2529:   /* Project initial state to characteristic variables */
2530:   eigenvectors(rv, lv, ueq, gamma);
2531:   projecteqstate(wc, ueq, lv);
2532:   wcp[0] = wc[0];
2533:   wcp[1] = wc[1];
2534:   wcp[2] = wc[2] + eps * cos(coord[0] * 2. * twopi / Lx);
2535:   projecttoprim(vp, wcp, rv);
2536:   ux->r = vp[0]; /* density */
2537:   ux->ru[0] = vp[0] * vp[1]; /* x momentum */
2538:   ux->ru[1] = 0.;
2539: #if defined DIM > 2
2540:   if (dim>2) ux->ru[2] = 0.;
2541: #endif
2542:   /* E = rho * e + rho * v^2/2 = p/(gam-1) + rho*v^2/2 */
2543:   ux->E = vp[2]/(gamma - 1.) + 0.5*vp[0]*vp[1]*vp[1];
2544:   return 0;
2545: }