xref: /petsc/src/ts/tutorials/ex11.c (revision 4302210d98c281cf4e8bc0fee406cf222e373a4c)
1 static char help[] = "Second Order TVD Finite Volume Example.\n";
2 /*F
3 
4 We use a second order TVD finite volume method to evolve a system of PDEs. Our simple upwinded residual evaluation loops
5 over all mesh faces and uses a Riemann solver to produce the flux given the face geometry and cell values,
6 \begin{equation}
7   f_i = \mathrm{riemann}(\mathrm{phys}, p_\mathrm{centroid}, \hat n, x^L, x^R)
8 \end{equation}
9 and then update the cell values given the cell volume.
10 \begin{eqnarray}
11     f^L_i &-=& \frac{f_i}{vol^L} \\
12     f^R_i &+=& \frac{f_i}{vol^R}
13 \end{eqnarray}
14 
15 As an example, we can consider the shallow water wave equation,
16 \begin{eqnarray}
17      h_t + \nabla\cdot \left( uh                              \right) &=& 0 \\
18   (uh)_t + \nabla\cdot \left( u\otimes uh + \frac{g h^2}{2} I \right) &=& 0
19 \end{eqnarray}
20 where $h$ is wave height, $u$ is wave velocity, and $g$ is the acceleration due to gravity.
21 
22 A representative Riemann solver for the shallow water equations is given in the PhysicsRiemann_SW() function,
23 \begin{eqnarray}
24   f^{L,R}_h    &=& uh^{L,R} \cdot \hat n \\
25   f^{L,R}_{uh} &=& \frac{f^{L,R}_h}{h^{L,R}} uh^{L,R} + g (h^{L,R})^2 \hat n \\
26   c^{L,R}      &=& \sqrt{g h^{L,R}} \\
27   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) \\
28   f_i          &=& \frac{A_\mathrm{face}}{2} \left( f^L_i + f^R_i + s \left( x^L_i - x^R_i \right) \right)
29 \end{eqnarray}
30 where $c$ is the local gravity wave speed and $f_i$ is a Rusanov flux.
31 
32 The more sophisticated residual evaluation in RHSFunctionLocal_LS() uses a least-squares fit to a quadratic polynomial
33 over a neighborhood of the given element.
34 
35 The mesh is read in from an ExodusII file, usually generated by Cubit.
36 F*/
37 #include <petscdmplex.h>
38 #include <petscdmforest.h>
39 #include <petscds.h>
40 #include <petscts.h>
41 #include <petscsf.h> /* For SplitFaces() */
42 
43 #define DIM 2                   /* Geometric dimension */
44 #define ALEN(a) (sizeof(a)/sizeof((a)[0]))
45 
46 static PetscFunctionList PhysicsList, PhysicsRiemannList_SW;
47 
48 /* Represents continuum physical equations. */
49 typedef struct _n_Physics *Physics;
50 
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;
54 
55 /* 'User' implements a discretization of a continuous model. */
56 typedef struct _n_User *User;
57 typedef PetscErrorCode (*SolutionFunction)(Model,PetscReal,const PetscReal*,PetscScalar*,void*);
58 typedef PetscErrorCode (*SetUpBCFunction)(DM,PetscDS,Physics);
59 typedef PetscErrorCode (*FunctionalFunction)(Model,PetscReal,const PetscReal*,const PetscScalar*,PetscReal*,void*);
60 typedef PetscErrorCode (*SetupFields)(Physics,PetscSection);
61 static PetscErrorCode ModelSolutionSetDefault(Model,SolutionFunction,void*);
62 static PetscErrorCode ModelFunctionalRegister(Model,const char*,PetscInt*,FunctionalFunction,void*);
63 static PetscErrorCode OutputVTK(DM,const char*,PetscViewer*);
64 
65 struct FieldDescription {
66   const char *name;
67   PetscInt dof;
68 };
69 
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 };
78 
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 };
87 
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 };
106 
107 struct _n_User {
108   PetscInt numSplitFaces;
109   PetscInt vtkInterval;   /* For monitor */
110   char outputBasename[PETSC_MAX_PATH_LEN]; /* Basename for output files */
111   PetscInt monitorStepOffset;
112   Model    model;
113   PetscBool vtkmon;
114 };
115 
116 PETSC_STATIC_INLINE PetscReal DotDIMReal(const PetscReal *x,const PetscReal *y)
117 {
118   PetscInt  i;
119   PetscReal prod=0.0;
120 
121   for (i=0; i<DIM; i++) prod += x[i]*y[i];
122   return prod;
123 }
124 PETSC_STATIC_INLINE PetscReal NormDIM(const PetscReal *x) { return PetscSqrtReal(PetscAbsReal(DotDIMReal(x,x))); }
125 
126 PETSC_STATIC_INLINE PetscReal Dot2Real(const PetscReal *x,const PetscReal *y) { return x[0]*y[0] + x[1]*y[1];}
127 PETSC_STATIC_INLINE PetscReal Norm2Real(const PetscReal *x) { return PetscSqrtReal(PetscAbsReal(Dot2Real(x,x)));}
128 PETSC_STATIC_INLINE void Normalize2Real(PetscReal *x) { PetscReal a = 1./Norm2Real(x); x[0] *= a; x[1] *= a; }
129 PETSC_STATIC_INLINE void Waxpy2Real(PetscReal a,const PetscReal *x,const PetscReal *y,PetscReal *w) { w[0] = a*x[0] + y[0]; w[1] = a*x[1] + y[1]; }
130 PETSC_STATIC_INLINE void Scale2Real(PetscReal a,const PetscReal *x,PetscReal *y) { y[0] = a*x[0]; y[1] = a*x[1]; }
131 
132 /******************* Advect ********************/
133 typedef enum {ADVECT_SOL_TILTED,ADVECT_SOL_BUMP,ADVECT_SOL_BUMP_CAVITY} AdvectSolType;
134 static const char *const AdvectSolTypes[] = {"TILTED","BUMP","BUMP_CAVITY","AdvectSolType","ADVECT_SOL_",0};
135 typedef enum {ADVECT_SOL_BUMP_CONE,ADVECT_SOL_BUMP_COS} AdvectSolBumpType;
136 static const char *const AdvectSolBumpTypes[] = {"CONE","COS","AdvectSolBumpType","ADVECT_SOL_BUMP_",0};
137 
138 typedef struct {
139   PetscReal wind[DIM];
140 } Physics_Advect_Tilted;
141 typedef struct {
142   PetscReal         center[DIM];
143   PetscReal         radius;
144   AdvectSolBumpType type;
145 } Physics_Advect_Bump;
146 
147 typedef struct {
148   PetscReal     inflowState;
149   AdvectSolType soltype;
150   union {
151     Physics_Advect_Tilted tilted;
152     Physics_Advect_Bump   bump;
153   } sol;
154   struct {
155     PetscInt Solution;
156     PetscInt Error;
157   } functional;
158 } Physics_Advect;
159 
160 static const struct FieldDescription PhysicsFields_Advect[] = {{"U",1},{NULL,0}};
161 
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;
166 
167   PetscFunctionBeginUser;
168   xG[0] = advect->inflowState;
169   PetscFunctionReturn(0);
170 }
171 
172 static PetscErrorCode PhysicsBoundary_Advect_Outflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
173 {
174   PetscFunctionBeginUser;
175   xG[0] = xI[0];
176   PetscFunctionReturn(0);
177 }
178 
179 static void PhysicsRiemann_Advect(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, Physics phys)
180 {
181   Physics_Advect *advect = (Physics_Advect*)phys->data;
182   PetscReal      wind[DIM],wn;
183 
184   switch (advect->soltype) {
185   case ADVECT_SOL_TILTED: {
186     Physics_Advect_Tilted *tilted = &advect->sol.tilted;
187     wind[0] = tilted->wind[0];
188     wind[1] = tilted->wind[1];
189   } break;
190   case ADVECT_SOL_BUMP:
191     wind[0] = -qp[1];
192     wind[1] = qp[0];
193     break;
194   case ADVECT_SOL_BUMP_CAVITY:
195     {
196       PetscInt  i;
197       PetscReal comp2[3] = {0.,0.,0.}, rad2;
198 
199       rad2 = 0.;
200       for (i = 0; i < dim; i++) {
201         comp2[i] = qp[i] * qp[i];
202         rad2    += comp2[i];
203       }
204 
205       wind[0] = -qp[1];
206       wind[1] = qp[0];
207       if (rad2 > 1.) {
208         PetscInt  maxI = 0;
209         PetscReal maxComp2 = comp2[0];
210 
211         for (i = 1; i < dim; i++) {
212           if (comp2[i] > maxComp2) {
213             maxI     = i;
214             maxComp2 = comp2[i];
215           }
216         }
217         wind[maxI] = 0.;
218       }
219     }
220     break;
221   default:
222   {
223     PetscInt i;
224     for (i = 0; i < DIM; ++i) wind[i] = 0.0;
225   }
226   /* default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for solution type %s",AdvectSolBumpTypes[advect->soltype]); */
227   }
228   wn      = Dot2Real(wind, n);
229   flux[0] = (wn > 0 ? xL[0] : xR[0]) * wn;
230 }
231 
232 static PetscErrorCode PhysicsSolution_Advect(Model mod,PetscReal time,const PetscReal *x,PetscScalar *u,void *ctx)
233 {
234   Physics        phys    = (Physics)ctx;
235   Physics_Advect *advect = (Physics_Advect*)phys->data;
236 
237   PetscFunctionBeginUser;
238   switch (advect->soltype) {
239   case ADVECT_SOL_TILTED: {
240     PetscReal             x0[DIM];
241     Physics_Advect_Tilted *tilted = &advect->sol.tilted;
242     Waxpy2Real(-time,tilted->wind,x,x0);
243     if (x0[1] > 0) u[0] = 1.*x[0] + 3.*x[1];
244     else u[0] = advect->inflowState;
245   } break;
246   case ADVECT_SOL_BUMP_CAVITY:
247   case ADVECT_SOL_BUMP: {
248     Physics_Advect_Bump *bump = &advect->sol.bump;
249     PetscReal           x0[DIM],v[DIM],r,cost,sint;
250     cost  = PetscCosReal(time);
251     sint  = PetscSinReal(time);
252     x0[0] = cost*x[0] + sint*x[1];
253     x0[1] = -sint*x[0] + cost*x[1];
254     Waxpy2Real(-1,bump->center,x0,v);
255     r = Norm2Real(v);
256     switch (bump->type) {
257     case ADVECT_SOL_BUMP_CONE:
258       u[0] = PetscMax(1 - r/bump->radius,0);
259       break;
260     case ADVECT_SOL_BUMP_COS:
261       u[0] = 0.5 + 0.5*PetscCosReal(PetscMin(r/bump->radius,1)*PETSC_PI);
262       break;
263     }
264   } break;
265   default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown solution type");
266   }
267   PetscFunctionReturn(0);
268 }
269 
270 static PetscErrorCode PhysicsFunctional_Advect(Model mod,PetscReal time,const PetscReal *x,const PetscScalar *y,PetscReal *f,void *ctx)
271 {
272   Physics        phys    = (Physics)ctx;
273   Physics_Advect *advect = (Physics_Advect*)phys->data;
274   PetscScalar    yexact[1] = {0.0};
275   PetscErrorCode ierr;
276 
277   PetscFunctionBeginUser;
278   ierr = PhysicsSolution_Advect(mod,time,x,yexact,phys);CHKERRQ(ierr);
279   f[advect->functional.Solution] = PetscRealPart(y[0]);
280   f[advect->functional.Error] = PetscAbsScalar(y[0]-yexact[0]);
281   PetscFunctionReturn(0);
282 }
283 
284 static PetscErrorCode SetUpBC_Advect(DM dm, PetscDS prob, Physics phys)
285 {
286   PetscErrorCode ierr;
287   const PetscInt inflowids[] = {100,200,300},outflowids[] = {101};
288   DMLabel        label;
289 
290   PetscFunctionBeginUser;
291   /* Register "canned" boundary conditions and defaults for where to apply. */
292   ierr = DMGetLabel(dm, "Face Sets", &label);CHKERRQ(ierr);
293   ierr = PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "inflow",  label, ALEN(inflowids),  inflowids,  0, 0, NULL, (void (*)(void)) PhysicsBoundary_Advect_Inflow, NULL,  phys, NULL);CHKERRQ(ierr);
294   ierr = PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "outflow", label, ALEN(outflowids), outflowids, 0, 0, NULL, (void (*)(void)) PhysicsBoundary_Advect_Outflow, NULL, phys, NULL);CHKERRQ(ierr);
295   PetscFunctionReturn(0);
296 }
297 
298 static PetscErrorCode PhysicsCreate_Advect(Model mod,Physics phys,PetscOptionItems *PetscOptionsObject)
299 {
300   Physics_Advect *advect;
301   PetscErrorCode ierr;
302 
303   PetscFunctionBeginUser;
304   phys->field_desc = PhysicsFields_Advect;
305   phys->riemann    = (PetscRiemannFunc)PhysicsRiemann_Advect;
306   ierr = PetscNew(&advect);CHKERRQ(ierr);
307   phys->data       = advect;
308   mod->setupbc = SetUpBC_Advect;
309 
310   ierr = PetscOptionsHead(PetscOptionsObject,"Advect options");CHKERRQ(ierr);
311   {
312     PetscInt two = 2,dof = 1;
313     advect->soltype = ADVECT_SOL_TILTED;
314     ierr = PetscOptionsEnum("-advect_sol_type","solution type","",AdvectSolTypes,(PetscEnum)advect->soltype,(PetscEnum*)&advect->soltype,NULL);CHKERRQ(ierr);
315     switch (advect->soltype) {
316     case ADVECT_SOL_TILTED: {
317       Physics_Advect_Tilted *tilted = &advect->sol.tilted;
318       two = 2;
319       tilted->wind[0] = 0.0;
320       tilted->wind[1] = 1.0;
321       ierr = PetscOptionsRealArray("-advect_tilted_wind","background wind vx,vy","",tilted->wind,&two,NULL);CHKERRQ(ierr);
322       advect->inflowState = -2.0;
323       ierr = PetscOptionsRealArray("-advect_tilted_inflow","Inflow state","",&advect->inflowState,&dof,NULL);CHKERRQ(ierr);
324       phys->maxspeed = Norm2Real(tilted->wind);
325     } break;
326     case ADVECT_SOL_BUMP_CAVITY:
327     case ADVECT_SOL_BUMP: {
328       Physics_Advect_Bump *bump = &advect->sol.bump;
329       two = 2;
330       bump->center[0] = 2.;
331       bump->center[1] = 0.;
332       ierr = PetscOptionsRealArray("-advect_bump_center","location of center of bump x,y","",bump->center,&two,NULL);CHKERRQ(ierr);
333       bump->radius = 0.9;
334       ierr = PetscOptionsReal("-advect_bump_radius","radius of bump","",bump->radius,&bump->radius,NULL);CHKERRQ(ierr);
335       bump->type = ADVECT_SOL_BUMP_CONE;
336       ierr = PetscOptionsEnum("-advect_bump_type","type of bump","",AdvectSolBumpTypes,(PetscEnum)bump->type,(PetscEnum*)&bump->type,NULL);CHKERRQ(ierr);
337       phys->maxspeed = 3.;       /* radius of mesh, kludge */
338     } break;
339     }
340   }
341   ierr = PetscOptionsTail();CHKERRQ(ierr);
342   /* Initial/transient solution with default boundary conditions */
343   ierr = ModelSolutionSetDefault(mod,PhysicsSolution_Advect,phys);CHKERRQ(ierr);
344   /* Register "canned" functionals */
345   ierr = ModelFunctionalRegister(mod,"Solution",&advect->functional.Solution,PhysicsFunctional_Advect,phys);CHKERRQ(ierr);
346   ierr = ModelFunctionalRegister(mod,"Error",&advect->functional.Error,PhysicsFunctional_Advect,phys);CHKERRQ(ierr);
347   mod->bcs[0] = mod->bcs[1] = mod->bcs[2] = DM_BOUNDARY_GHOSTED;
348   PetscFunctionReturn(0);
349 }
350 
351 /******************* Shallow Water ********************/
352 typedef struct {
353   PetscReal gravity;
354   PetscReal boundaryHeight;
355   struct {
356     PetscInt Height;
357     PetscInt Speed;
358     PetscInt Energy;
359   } functional;
360 } Physics_SW;
361 typedef struct {
362   PetscReal h;
363   PetscReal uh[DIM];
364 } SWNode;
365 typedef union {
366   SWNode    swnode;
367   PetscReal vals[DIM+1];
368 } SWNodeUnion;
369 
370 static const struct FieldDescription PhysicsFields_SW[] = {{"Height",1},{"Momentum",DIM},{NULL,0}};
371 
372 /*
373  * h_t + div(uh) = 0
374  * (uh)_t + div (u\otimes uh + g h^2 / 2 I) = 0
375  *
376  * */
377 static PetscErrorCode SWFlux(Physics phys,const PetscReal *n,const SWNode *x,SWNode *f)
378 {
379   Physics_SW  *sw = (Physics_SW*)phys->data;
380   PetscReal   uhn,u[DIM];
381   PetscInt     i;
382 
383   PetscFunctionBeginUser;
384   Scale2Real(1./x->h,x->uh,u);
385   uhn  = x->uh[0] * n[0] + x->uh[1] * n[1];
386   f->h = uhn;
387   for (i=0; i<DIM; i++) f->uh[i] = u[i] * uhn + sw->gravity * PetscSqr(x->h) * n[i];
388   PetscFunctionReturn(0);
389 }
390 
391 static PetscErrorCode PhysicsBoundary_SW_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
392 {
393   PetscFunctionBeginUser;
394   xG[0] = xI[0];
395   xG[1] = -xI[1];
396   xG[2] = -xI[2];
397   PetscFunctionReturn(0);
398 }
399 
400 static void PhysicsRiemann_SW_HLL(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, Physics phys)
401 {
402   Physics_SW *sw = (Physics_SW *) phys->data;
403   PetscReal aL, aR;
404   PetscReal nn[DIM];
405 #if !defined(PETSC_USE_COMPLEX)
406   const SWNode *uL = (const SWNode *) xL, *uR = (const SWNode *) xR;
407 #else
408   SWNodeUnion  uLreal, uRreal;
409   const SWNode *uL = &uLreal.swnode;
410   const SWNode *uR = &uRreal.swnode;
411 #endif
412   SWNodeUnion fL, fR;
413   PetscInt i;
414   PetscReal zero = 0.;
415 
416 #if defined(PETSC_USE_COMPLEX)
417   uLreal.swnode.h = 0; uRreal.swnode.h = 0;
418   for (i = 0; i < 1+dim; i++) uLreal.vals[i] = PetscRealPart(xL[i]);
419   for (i = 0; i < 1+dim; i++) uRreal.vals[i] = PetscRealPart(xR[i]);
420 #endif
421   if (uL->h <= 0 || uR->h <= 0) {
422     for (i = 0; i < 1 + dim; i++) flux[i] = zero;
423     return;
424   } /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reconstructed thickness is negative"); */
425   nn[0] = n[0];
426   nn[1] = n[1];
427   Normalize2Real(nn);
428   SWFlux(phys, nn, uL, &(fL.swnode));
429   SWFlux(phys, nn, uR, &(fR.swnode));
430   /* gravity wave speed */
431   aL = PetscSqrtReal(sw->gravity * uL->h);
432   aR = PetscSqrtReal(sw->gravity * uR->h);
433   // Defining u_tilda and v_tilda as u and v
434   PetscReal u_L, u_R;
435   u_L = Dot2Real(uL->uh,nn)/uL->h;
436   u_R = Dot2Real(uR->uh,nn)/uR->h;
437   PetscReal sL, sR;
438   sL = PetscMin(u_L - aL, u_R - aR);
439   sR = PetscMax(u_L + aL, u_R + aR);
440   if (sL > zero) {
441     for (i = 0; i < dim + 1; i++) {
442       flux[i] = fL.vals[i] * Norm2Real(n);
443     }
444   } else if (sR < zero) {
445     for (i = 0; i < dim + 1; i++) {
446       flux[i] = fR.vals[i] * Norm2Real(n);
447     }
448   } else {
449     for (i = 0; i < dim + 1; i++) {
450       flux[i] = ((sR * fL.vals[i] - sL * fR.vals[i] + sR * sL * (xR[i] - xL[i])) / (sR - sL)) * Norm2Real(n);
451     }
452   }
453 }
454 
455 static void PhysicsRiemann_SW_Rusanov(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, Physics phys)
456 {
457   Physics_SW   *sw = (Physics_SW*)phys->data;
458   PetscReal    cL,cR,speed;
459   PetscReal    nn[DIM];
460 #if !defined(PETSC_USE_COMPLEX)
461   const SWNode *uL = (const SWNode*)xL,*uR = (const SWNode*)xR;
462 #else
463   SWNodeUnion  uLreal, uRreal;
464   const SWNode *uL = &uLreal.swnode;
465   const SWNode *uR = &uRreal.swnode;
466 #endif
467   SWNodeUnion  fL,fR;
468   PetscInt     i;
469   PetscReal    zero=0.;
470 
471 #if defined(PETSC_USE_COMPLEX)
472   uLreal.swnode.h = 0; uRreal.swnode.h = 0;
473   for (i = 0; i < 1+dim; i++) uLreal.vals[i] = PetscRealPart(xL[i]);
474   for (i = 0; i < 1+dim; i++) uRreal.vals[i] = PetscRealPart(xR[i]);
475 #endif
476   if (uL->h < 0 || uR->h < 0) {for (i=0; i<1+dim; i++) flux[i] = zero/zero; return;} /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reconstructed thickness is negative"); */
477   nn[0] = n[0];
478   nn[1] = n[1];
479   Normalize2Real(nn);
480   SWFlux(phys,nn,uL,&(fL.swnode));
481   SWFlux(phys,nn,uR,&(fR.swnode));
482   cL    = PetscSqrtReal(sw->gravity*uL->h);
483   cR    = PetscSqrtReal(sw->gravity*uR->h); /* gravity wave speed */
484   speed = PetscMax(PetscAbsReal(Dot2Real(uL->uh,nn)/uL->h) + cL,PetscAbsReal(Dot2Real(uR->uh,nn)/uR->h) + cR);
485   for (i=0; i<1+dim; i++) flux[i] = (0.5*(fL.vals[i] + fR.vals[i]) + 0.5*speed*(xL[i] - xR[i])) * Norm2Real(n);
486 }
487 
488 static PetscErrorCode PhysicsSolution_SW(Model mod,PetscReal time,const PetscReal *x,PetscScalar *u,void *ctx)
489 {
490   PetscReal dx[2],r,sigma;
491 
492   PetscFunctionBeginUser;
493   if (time != 0.0) SETERRQ1(mod->comm,PETSC_ERR_SUP,"No solution known for time %g",(double)time);
494   dx[0] = x[0] - 1.5;
495   dx[1] = x[1] - 1.0;
496   r     = Norm2Real(dx);
497   sigma = 0.5;
498   u[0]  = 1 + 2*PetscExpReal(-PetscSqr(r)/(2*PetscSqr(sigma)));
499   u[1]  = 0.0;
500   u[2]  = 0.0;
501   PetscFunctionReturn(0);
502 }
503 
504 static PetscErrorCode PhysicsFunctional_SW(Model mod,PetscReal time,const PetscReal *coord,const PetscScalar *xx,PetscReal *f,void *ctx)
505 {
506   Physics      phys = (Physics)ctx;
507   Physics_SW   *sw  = (Physics_SW*)phys->data;
508   const SWNode *x   = (const SWNode*)xx;
509   PetscReal  u[2];
510   PetscReal    h;
511 
512   PetscFunctionBeginUser;
513   h = x->h;
514   Scale2Real(1./x->h,x->uh,u);
515   f[sw->functional.Height] = h;
516   f[sw->functional.Speed]  = Norm2Real(u) + PetscSqrtReal(sw->gravity*h);
517   f[sw->functional.Energy] = 0.5*(Dot2Real(x->uh,u) + sw->gravity*PetscSqr(h));
518   PetscFunctionReturn(0);
519 }
520 
521 static PetscErrorCode SetUpBC_SW(DM dm, PetscDS prob,Physics phys)
522 {
523   PetscErrorCode ierr;
524   const PetscInt wallids[] = {100,101,200,300};
525   DMLabel        label;
526 
527   PetscFunctionBeginUser;
528   ierr = DMGetLabel(dm, "Face Sets", &label);CHKERRQ(ierr);
529   ierr = PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", label, ALEN(wallids), wallids, 0, 0, NULL, (void (*)(void)) PhysicsBoundary_SW_Wall, NULL, phys, NULL);CHKERRQ(ierr);
530   PetscFunctionReturn(0);
531 }
532 
533 static PetscErrorCode PhysicsCreate_SW(Model mod,Physics phys,PetscOptionItems *PetscOptionsObject)
534 {
535   Physics_SW     *sw;
536   char           sw_riemann[64] = "rusanov";
537   PetscErrorCode ierr;
538 
539   PetscFunctionBeginUser;
540   phys->field_desc = PhysicsFields_SW;
541   ierr          = PetscNew(&sw);CHKERRQ(ierr);
542   phys->data    = sw;
543   mod->setupbc  = SetUpBC_SW;
544 
545   PetscFunctionListAdd(&PhysicsRiemannList_SW, "rusanov", PhysicsRiemann_SW_Rusanov);
546   PetscFunctionListAdd(&PhysicsRiemannList_SW, "hll", PhysicsRiemann_SW_HLL);
547 
548   ierr          = PetscOptionsHead(PetscOptionsObject,"SW options");CHKERRQ(ierr);
549   {
550     void (*PhysicsRiemann_SW)(PetscInt, PetscInt, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, PetscInt, const PetscScalar, PetscScalar *, Physics);
551     sw->gravity = 1.0;
552     ierr = PetscOptionsReal("-sw_gravity","Gravitational constant","",sw->gravity,&sw->gravity,NULL);CHKERRQ(ierr);
553     ierr = PetscOptionsFList("-sw_riemann","Riemann solver","",PhysicsRiemannList_SW,sw_riemann,sw_riemann,sizeof sw_riemann,NULL);CHKERRQ(ierr);
554     ierr = PetscFunctionListFind(PhysicsRiemannList_SW,sw_riemann,&PhysicsRiemann_SW);CHKERRQ(ierr);
555     phys->riemann = (PetscRiemannFunc) PhysicsRiemann_SW;
556   }
557   ierr = PetscOptionsTail();CHKERRQ(ierr);
558   phys->maxspeed = PetscSqrtReal(2.0*sw->gravity); /* Mach 1 for depth of 2 */
559 
560   ierr = ModelSolutionSetDefault(mod,PhysicsSolution_SW,phys);CHKERRQ(ierr);
561   ierr = ModelFunctionalRegister(mod,"Height",&sw->functional.Height,PhysicsFunctional_SW,phys);CHKERRQ(ierr);
562   ierr = ModelFunctionalRegister(mod,"Speed",&sw->functional.Speed,PhysicsFunctional_SW,phys);CHKERRQ(ierr);
563   ierr = ModelFunctionalRegister(mod,"Energy",&sw->functional.Energy,PhysicsFunctional_SW,phys);CHKERRQ(ierr);
564 
565   mod->bcs[0] = mod->bcs[1] = mod->bcs[2] = DM_BOUNDARY_GHOSTED;
566 
567   PetscFunctionReturn(0);
568 }
569 
570 /******************* Euler Density Shock (EULER_IV_SHOCK,EULER_SS_SHOCK) ********************/
571 /* An initial-value and self-similar solutions of the compressible Euler equations */
572 /* Ravi Samtaney and D. I. Pullin */
573 /* Phys. Fluids 8, 2650 (1996); http://dx.doi.org/10.1063/1.869050 */
574 typedef enum {EULER_PAR_GAMMA,EULER_PAR_RHOR,EULER_PAR_AMACH,EULER_PAR_ITANA,EULER_PAR_SIZE} EulerParamIdx;
575 typedef enum {EULER_IV_SHOCK,EULER_SS_SHOCK,EULER_SHOCK_TUBE,EULER_LINEAR_WAVE} EulerType;
576 typedef struct {
577   PetscReal r;
578   PetscReal ru[DIM];
579   PetscReal E;
580 } EulerNode;
581 typedef union {
582   EulerNode eulernode;
583   PetscReal vals[DIM+2];
584 } EulerNodeUnion;
585 typedef PetscErrorCode (*EquationOfState)(const PetscReal*, const EulerNode*, PetscReal*);
586 typedef struct {
587   EulerType       type;
588   PetscReal       pars[EULER_PAR_SIZE];
589   EquationOfState sound;
590   struct {
591     PetscInt Density;
592     PetscInt Momentum;
593     PetscInt Energy;
594     PetscInt Pressure;
595     PetscInt Speed;
596   } monitor;
597 } Physics_Euler;
598 
599 static const struct FieldDescription PhysicsFields_Euler[] = {{"Density",1},{"Momentum",DIM},{"Energy",1},{NULL,0}};
600 
601 /* initial condition */
602 int initLinearWave(EulerNode *ux, const PetscReal gamma, const PetscReal coord[], const PetscReal Lx);
603 static PetscErrorCode PhysicsSolution_Euler(Model mod, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx)
604 {
605   PetscInt i;
606   Physics         phys = (Physics)ctx;
607   Physics_Euler   *eu  = (Physics_Euler*)phys->data;
608   EulerNode       *uu  = (EulerNode*)u;
609   PetscReal        p0,gamma,c;
610   PetscFunctionBeginUser;
611   if (time != 0.0) SETERRQ1(mod->comm,PETSC_ERR_SUP,"No solution known for time %g",(double)time);
612 
613   for (i=0; i<DIM; i++) uu->ru[i] = 0.0; /* zero out initial velocity */
614   /* set E and rho */
615   gamma = eu->pars[EULER_PAR_GAMMA];
616 
617   if (eu->type==EULER_IV_SHOCK || eu->type==EULER_SS_SHOCK) {
618     /******************* Euler Density Shock ********************/
619     /* On initial-value and self-similar solutions of the compressible Euler equations */
620     /* Ravi Samtaney and D. I. Pullin */
621     /* Phys. Fluids 8, 2650 (1996); http://dx.doi.org/10.1063/1.869050 */
622     /* initial conditions 1: left of shock, 0: left of discontinuity 2: right of discontinuity,  */
623     p0 = 1.;
624     if (x[0] < 0.0 + x[1]*eu->pars[EULER_PAR_ITANA]) {
625       if (x[0] < mod->bounds[0]*0.5) { /* left of shock (1) */
626         PetscReal amach,rho,press,gas1,p1;
627         amach = eu->pars[EULER_PAR_AMACH];
628         rho = 1.;
629         press = p0;
630         p1 = press*(1.0+2.0*gamma/(gamma+1.0)*(amach*amach-1.0));
631         gas1 = (gamma-1.0)/(gamma+1.0);
632         uu->r = rho*(p1/press+gas1)/(gas1*p1/press+1.0);
633         uu->ru[0]   = ((uu->r - rho)*PetscSqrtReal(gamma*press/rho)*amach);
634         uu->E = p1/(gamma-1.0) + .5/uu->r*uu->ru[0]*uu->ru[0];
635       }
636       else { /* left of discontinuity (0) */
637         uu->r = 1.; /* rho = 1 */
638         uu->E = p0/(gamma-1.0);
639       }
640     }
641     else { /* right of discontinuity (2) */
642       uu->r = eu->pars[EULER_PAR_RHOR];
643       uu->E = p0/(gamma-1.0);
644     }
645   }
646   else if (eu->type==EULER_SHOCK_TUBE) {
647     /* For (x<x0) set (rho,u,p)=(8,0,10) and for (x>x0) set (rho,u,p)=(1,0,1). Choose x0 to the midpoint of the domain in the x-direction. */
648     if (x[0] < 0.0) {
649       uu->r = 8.;
650       uu->E = 10./(gamma-1.);
651     }
652     else {
653       uu->r = 1.;
654       uu->E = 1./(gamma-1.);
655     }
656   }
657   else if (eu->type==EULER_LINEAR_WAVE) {
658     initLinearWave( uu, gamma, x, mod->bounds[1] - mod->bounds[0]);
659   }
660   else SETERRQ1(mod->comm,PETSC_ERR_SUP,"Unknown type %d",eu->type);
661 
662   /* set phys->maxspeed: (mod->maxspeed = phys->maxspeed) in main; */
663   eu->sound(&gamma,uu,&c);
664   c = (uu->ru[0]/uu->r) + c;
665   if (c > phys->maxspeed) phys->maxspeed = c;
666 
667   PetscFunctionReturn(0);
668 }
669 
670 static PetscErrorCode Pressure_PG(const PetscReal gamma,const EulerNode *x,PetscReal *p)
671 {
672   PetscReal ru2;
673 
674   PetscFunctionBeginUser;
675   ru2  = DotDIMReal(x->ru,x->ru);
676   (*p)=(x->E - 0.5*ru2/x->r)*(gamma - 1.0); /* (E - rho V^2/2)(gamma-1) = e rho (gamma-1) */
677   PetscFunctionReturn(0);
678 }
679 
680 static PetscErrorCode SpeedOfSound_PG(const PetscReal *gamma, const EulerNode *x, PetscReal *c)
681 {
682   PetscReal p;
683 
684   PetscFunctionBeginUser;
685   Pressure_PG(*gamma,x,&p);
686   if (p<0.) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"negative pressure time %g -- NEED TO FIX!!!!!!",(double) p);
687   /* pars[EULER_PAR_GAMMA] = heat capacity ratio */
688   (*c)=PetscSqrtReal(*gamma * p / x->r);
689   PetscFunctionReturn(0);
690 }
691 
692 /*
693  * x = (rho,rho*(u_1),...,rho*e)^T
694  * x_t+div(f_1(x))+...+div(f_DIM(x)) = 0
695  *
696  * f_i(x) = u_i*x+(0,0,...,p,...,p*u_i)^T
697  *
698  */
699 static PetscErrorCode EulerFlux(Physics phys,const PetscReal *n,const EulerNode *x,EulerNode *f)
700 {
701   Physics_Euler *eu = (Physics_Euler*)phys->data;
702   PetscReal     nu,p;
703   PetscInt      i;
704 
705   PetscFunctionBeginUser;
706   Pressure_PG(eu->pars[EULER_PAR_GAMMA],x,&p);
707   nu = DotDIMReal(x->ru,n);
708   f->r = nu;   /* A rho u */
709   nu /= x->r;  /* A u */
710   for (i=0; i<DIM; i++) f->ru[i] = nu * x->ru[i] + n[i]*p;  /* r u^2 + p */
711   f->E = nu * (x->E + p); /* u(e+p) */
712   PetscFunctionReturn(0);
713 }
714 
715 /* PetscReal* => EulerNode* conversion */
716 static PetscErrorCode PhysicsBoundary_Euler_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *a_xI, PetscScalar *a_xG, void *ctx)
717 {
718   PetscInt    i;
719   const EulerNode *xI = (const EulerNode*)a_xI;
720   EulerNode       *xG = (EulerNode*)a_xG;
721   Physics         phys = (Physics)ctx;
722   Physics_Euler   *eu  = (Physics_Euler*)phys->data;
723   PetscFunctionBeginUser;
724   xG->r = xI->r;           /* ghost cell density - same */
725   xG->E = xI->E;           /* ghost cell energy - same */
726   if (n[1] != 0.) {        /* top and bottom */
727     xG->ru[0] =  xI->ru[0]; /* copy tang to wall */
728     xG->ru[1] = -xI->ru[1]; /* reflect perp to t/b wall */
729   }
730   else { /* sides */
731     for (i=0; i<DIM; i++) xG->ru[i] = xI->ru[i]; /* copy */
732   }
733   if (eu->type == EULER_LINEAR_WAVE) { /* debug */
734 #if 0
735     PetscPrintf(PETSC_COMM_WORLD,"%s coord=%g,%g\n",PETSC_FUNCTION_NAME,c[0],c[1]);
736 #endif
737   }
738   PetscFunctionReturn(0);
739 }
740 int godunovflux( const PetscScalar *ul, const PetscScalar *ur, PetscScalar *flux, const PetscReal *nn, const int *ndim, const PetscReal *gamma);
741 /* PetscReal* => EulerNode* conversion */
742 static void PhysicsRiemann_Euler_Godunov( PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n,
743                                           const PetscScalar *xL, const PetscScalar *xR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, Physics phys)
744 {
745   Physics_Euler   *eu = (Physics_Euler*)phys->data;
746   PetscReal       cL,cR,speed,velL,velR,nn[DIM],s2;
747   PetscInt        i;
748   PetscErrorCode  ierr;
749   PetscFunctionBeginUser;
750 
751   for (i=0,s2=0.; i<DIM; i++) {
752     nn[i] = n[i];
753     s2 += nn[i]*nn[i];
754   }
755   s2 = PetscSqrtReal(s2); /* |n|_2 = sum(n^2)^1/2 */
756   for (i=0.; i<DIM; i++) nn[i] /= s2;
757   if (0) { /* Rusanov */
758     const EulerNode *uL = (const EulerNode*)xL,*uR = (const EulerNode*)xR;
759     EulerNodeUnion  fL,fR;
760     EulerFlux(phys,nn,uL,&(fL.eulernode));
761     EulerFlux(phys,nn,uR,&(fR.eulernode));
762     ierr = eu->sound(&eu->pars[EULER_PAR_GAMMA],uL,&cL);if (ierr) exit(13);
763     ierr = eu->sound(&eu->pars[EULER_PAR_GAMMA],uR,&cR);if (ierr) exit(14);
764     velL = DotDIMReal(uL->ru,nn)/uL->r;
765     velR = DotDIMReal(uR->ru,nn)/uR->r;
766     speed = PetscMax(velR + cR, velL + cL);
767     for (i=0; i<2+dim; i++) flux[i] = 0.5*((fL.vals[i]+fR.vals[i]) + speed*(xL[i] - xR[i]))*s2;
768   }
769   else {
770     int dim = DIM;
771     /* int iwave =  */
772     godunovflux(xL, xR, flux, nn, &dim, &eu->pars[EULER_PAR_GAMMA]);
773     for (i=0; i<2+dim; i++) flux[i] *= s2;
774   }
775   PetscFunctionReturnVoid();
776 }
777 
778 static PetscErrorCode PhysicsFunctional_Euler(Model mod,PetscReal time,const PetscReal *coord,const PetscScalar *xx,PetscReal *f,void *ctx)
779 {
780   Physics         phys = (Physics)ctx;
781   Physics_Euler   *eu  = (Physics_Euler*)phys->data;
782   const EulerNode *x   = (const EulerNode*)xx;
783   PetscReal       p;
784 
785   PetscFunctionBeginUser;
786   f[eu->monitor.Density]  = x->r;
787   f[eu->monitor.Momentum] = NormDIM(x->ru);
788   f[eu->monitor.Energy]   = x->E;
789   f[eu->monitor.Speed]    = NormDIM(x->ru)/x->r;
790   Pressure_PG(eu->pars[EULER_PAR_GAMMA], x, &p);
791   f[eu->monitor.Pressure] = p;
792   PetscFunctionReturn(0);
793 }
794 
795 static PetscErrorCode SetUpBC_Euler(DM dm, PetscDS prob,Physics phys)
796 {
797   PetscErrorCode  ierr;
798   Physics_Euler   *eu = (Physics_Euler *) phys->data;
799   DMLabel         label;
800 
801   ierr = DMGetLabel(dm, "Face Sets", &label);CHKERRQ(ierr);
802   if (eu->type == EULER_LINEAR_WAVE) {
803     const PetscInt wallids[] = {100,101};
804     ierr = PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", label, ALEN(wallids), wallids, 0, 0, NULL, (void (*)(void)) PhysicsBoundary_Euler_Wall, NULL, phys, NULL);CHKERRQ(ierr);
805   }
806   else {
807     const PetscInt wallids[] = {100,101,200,300};
808     ierr = PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", label, ALEN(wallids), wallids, 0, 0, NULL, (void (*)(void)) PhysicsBoundary_Euler_Wall, NULL, phys, NULL);CHKERRQ(ierr);
809   }
810   PetscFunctionReturn(0);
811 }
812 
813 static PetscErrorCode PhysicsCreate_Euler(Model mod,Physics phys,PetscOptionItems *PetscOptionsObject)
814 {
815   Physics_Euler   *eu;
816   PetscErrorCode  ierr;
817 
818   PetscFunctionBeginUser;
819   phys->field_desc = PhysicsFields_Euler;
820   phys->riemann = (PetscRiemannFunc) PhysicsRiemann_Euler_Godunov;
821   ierr = PetscNew(&eu);CHKERRQ(ierr);
822   phys->data    = eu;
823   mod->setupbc = SetUpBC_Euler;
824   ierr = PetscOptionsHead(PetscOptionsObject,"Euler options");CHKERRQ(ierr);
825   {
826     PetscReal alpha;
827     char type[64] = "linear_wave";
828     PetscBool  is;
829     mod->bcs[0] = mod->bcs[1] = mod->bcs[2] = DM_BOUNDARY_GHOSTED;
830     eu->pars[EULER_PAR_GAMMA] = 1.4;
831     eu->pars[EULER_PAR_AMACH] = 2.02;
832     eu->pars[EULER_PAR_RHOR] = 3.0;
833     eu->pars[EULER_PAR_ITANA] = 0.57735026918963; /* angle of Euler self similar (SS) shock */
834     ierr = PetscOptionsReal("-eu_gamma","Heat capacity ratio","",eu->pars[EULER_PAR_GAMMA],&eu->pars[EULER_PAR_GAMMA],NULL);CHKERRQ(ierr);
835     ierr = PetscOptionsReal("-eu_amach","Shock speed (Mach)","",eu->pars[EULER_PAR_AMACH],&eu->pars[EULER_PAR_AMACH],NULL);CHKERRQ(ierr);
836     ierr = PetscOptionsReal("-eu_rho2","Density right of discontinuity","",eu->pars[EULER_PAR_RHOR],&eu->pars[EULER_PAR_RHOR],NULL);CHKERRQ(ierr);
837     alpha = 60.;
838     ierr = PetscOptionsReal("-eu_alpha","Angle of discontinuity","",alpha,&alpha,NULL);CHKERRQ(ierr);
839     if (alpha<=0. || alpha>90.) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Alpha bust be > 0 and <= 90 (%g)",alpha);
840     eu->pars[EULER_PAR_ITANA] = 1./PetscTanReal( alpha * PETSC_PI / 180.0);
841     ierr = PetscOptionsString("-eu_type","Type of Euler test","",type,type,sizeof(type),NULL);CHKERRQ(ierr);
842     ierr = PetscStrcmp(type,"linear_wave", &is);CHKERRQ(ierr);
843     if (is) {
844       eu->type = EULER_LINEAR_WAVE;
845       mod->bcs[0] = mod->bcs[1] = mod->bcs[2] = DM_BOUNDARY_PERIODIC;
846       mod->bcs[1] = DM_BOUNDARY_GHOSTED; /* debug */
847       ierr = PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"linear_wave");CHKERRQ(ierr);
848     }
849     else {
850       if (DIM != 2) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DIM must be 2 unless linear wave test %s",type);
851       ierr = PetscStrcmp(type,"iv_shock", &is);CHKERRQ(ierr);
852       if (is) {
853         eu->type = EULER_IV_SHOCK;
854         ierr = PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"iv_shock");CHKERRQ(ierr);
855       }
856       else {
857         ierr = PetscStrcmp(type,"ss_shock", &is);CHKERRQ(ierr);
858         if (is) {
859           eu->type = EULER_SS_SHOCK;
860           ierr = PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"ss_shock");CHKERRQ(ierr);
861         }
862         else {
863           ierr = PetscStrcmp(type,"shock_tube", &is);CHKERRQ(ierr);
864           if (is) eu->type = EULER_SHOCK_TUBE;
865           else SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Unknown Euler type %s",type);
866           ierr = PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"shock_tube");CHKERRQ(ierr);
867         }
868       }
869     }
870   }
871   ierr = PetscOptionsTail();CHKERRQ(ierr);
872   eu->sound = SpeedOfSound_PG;
873   phys->maxspeed = 0.; /* will get set in solution */
874   ierr = ModelSolutionSetDefault(mod,PhysicsSolution_Euler,phys);CHKERRQ(ierr);
875   ierr = ModelFunctionalRegister(mod,"Speed",&eu->monitor.Speed,PhysicsFunctional_Euler,phys);CHKERRQ(ierr);
876   ierr = ModelFunctionalRegister(mod,"Energy",&eu->monitor.Energy,PhysicsFunctional_Euler,phys);CHKERRQ(ierr);
877   ierr = ModelFunctionalRegister(mod,"Density",&eu->monitor.Density,PhysicsFunctional_Euler,phys);CHKERRQ(ierr);
878   ierr = ModelFunctionalRegister(mod,"Momentum",&eu->monitor.Momentum,PhysicsFunctional_Euler,phys);CHKERRQ(ierr);
879   ierr = ModelFunctionalRegister(mod,"Pressure",&eu->monitor.Pressure,PhysicsFunctional_Euler,phys);CHKERRQ(ierr);
880 
881   PetscFunctionReturn(0);
882 }
883 
884 static PetscErrorCode ErrorIndicator_Simple(PetscInt dim, PetscReal volume, PetscInt numComps, const PetscScalar u[], const PetscScalar grad[], PetscReal *error, void *ctx)
885 {
886   PetscReal      err = 0.;
887   PetscInt       i, j;
888 
889   PetscFunctionBeginUser;
890   for (i = 0; i < numComps; i++) {
891     for (j = 0; j < dim; j++) {
892       err += PetscSqr(PetscRealPart(grad[i * dim + j]));
893     }
894   }
895   *error = volume * err;
896   PetscFunctionReturn(0);
897 }
898 
899 PetscErrorCode ConstructCellBoundary(DM dm, User user)
900 {
901   const char     *name   = "Cell Sets";
902   const char     *bdname = "split faces";
903   IS             regionIS, innerIS;
904   const PetscInt *regions, *cells;
905   PetscInt       numRegions, innerRegion, numCells, c;
906   PetscInt       cStart, cEnd, cEndInterior, fStart, fEnd;
907   PetscBool      hasLabel;
908   PetscErrorCode ierr;
909 
910   PetscFunctionBeginUser;
911   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
912   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
913   ierr = DMPlexGetGhostCellStratum(dm, &cEndInterior, NULL);CHKERRQ(ierr);
914 
915   ierr = DMHasLabel(dm, name, &hasLabel);CHKERRQ(ierr);
916   if (!hasLabel) PetscFunctionReturn(0);
917   ierr = DMGetLabelSize(dm, name, &numRegions);CHKERRQ(ierr);
918   if (numRegions != 2) PetscFunctionReturn(0);
919   /* Get the inner id */
920   ierr = DMGetLabelIdIS(dm, name, &regionIS);CHKERRQ(ierr);
921   ierr = ISGetIndices(regionIS, &regions);CHKERRQ(ierr);
922   innerRegion = regions[0];
923   ierr = ISRestoreIndices(regionIS, &regions);CHKERRQ(ierr);
924   ierr = ISDestroy(&regionIS);CHKERRQ(ierr);
925   /* Find the faces between cells in different regions, could call DMPlexCreateNeighborCSR() */
926   ierr = DMGetStratumIS(dm, name, innerRegion, &innerIS);CHKERRQ(ierr);
927   ierr = ISGetLocalSize(innerIS, &numCells);CHKERRQ(ierr);
928   ierr = ISGetIndices(innerIS, &cells);CHKERRQ(ierr);
929   ierr = DMCreateLabel(dm, bdname);CHKERRQ(ierr);
930   for (c = 0; c < numCells; ++c) {
931     const PetscInt cell = cells[c];
932     const PetscInt *faces;
933     PetscInt       numFaces, f;
934 
935     if ((cell < cStart) || (cell >= cEnd)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Got invalid point %d which is not a cell", cell);
936     ierr = DMPlexGetConeSize(dm, cell, &numFaces);CHKERRQ(ierr);
937     ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr);
938     for (f = 0; f < numFaces; ++f) {
939       const PetscInt face = faces[f];
940       const PetscInt *neighbors;
941       PetscInt       nC, regionA, regionB;
942 
943       if ((face < fStart) || (face >= fEnd)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Got invalid point %d which is not a face", face);
944       ierr = DMPlexGetSupportSize(dm, face, &nC);CHKERRQ(ierr);
945       if (nC != 2) continue;
946       ierr = DMPlexGetSupport(dm, face, &neighbors);CHKERRQ(ierr);
947       if ((neighbors[0] >= cEndInterior) || (neighbors[1] >= cEndInterior)) continue;
948       if ((neighbors[0] < cStart) || (neighbors[0] >= cEnd)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Got invalid point %d which is not a cell", neighbors[0]);
949       if ((neighbors[1] < cStart) || (neighbors[1] >= cEnd)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Got invalid point %d which is not a cell", neighbors[1]);
950       ierr = DMGetLabelValue(dm, name, neighbors[0], &regionA);CHKERRQ(ierr);
951       ierr = DMGetLabelValue(dm, name, neighbors[1], &regionB);CHKERRQ(ierr);
952       if (regionA < 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid label %s: Cell %d has no value", name, neighbors[0]);
953       if (regionB < 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid label %s: Cell %d has no value", name, neighbors[1]);
954       if (regionA != regionB) {
955         ierr = DMSetLabelValue(dm, bdname, faces[f], 1);CHKERRQ(ierr);
956       }
957     }
958   }
959   ierr = ISRestoreIndices(innerIS, &cells);CHKERRQ(ierr);
960   ierr = ISDestroy(&innerIS);CHKERRQ(ierr);
961   {
962     DMLabel label;
963 
964     ierr = DMGetLabel(dm, bdname, &label);CHKERRQ(ierr);
965     ierr = DMLabelView(label, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
966   }
967   PetscFunctionReturn(0);
968 }
969 
970 /* Right now, I have just added duplicate faces, which see both cells. We can
971 - Add duplicate vertices and decouple the face cones
972 - Disconnect faces from cells across the rotation gap
973 */
974 PetscErrorCode SplitFaces(DM *dmSplit, const char labelName[], User user)
975 {
976   DM             dm = *dmSplit, sdm;
977   PetscSF        sfPoint, gsfPoint;
978   PetscSection   coordSection, newCoordSection;
979   Vec            coordinates;
980   IS             idIS;
981   const PetscInt *ids;
982   PetscInt       *newpoints;
983   PetscInt       dim, depth, maxConeSize, maxSupportSize, numLabels, numGhostCells;
984   PetscInt       numFS, fs, pStart, pEnd, p, cEnd, cEndInterior, vStart, vEnd, v, fStart, fEnd, newf, d, l;
985   PetscBool      hasLabel;
986   PetscErrorCode ierr;
987 
988   PetscFunctionBeginUser;
989   ierr = DMHasLabel(dm, labelName, &hasLabel);CHKERRQ(ierr);
990   if (!hasLabel) PetscFunctionReturn(0);
991   ierr = DMCreate(PetscObjectComm((PetscObject)dm), &sdm);CHKERRQ(ierr);
992   ierr = DMSetType(sdm, DMPLEX);CHKERRQ(ierr);
993   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
994   ierr = DMSetDimension(sdm, dim);CHKERRQ(ierr);
995 
996   ierr = DMGetLabelIdIS(dm, labelName, &idIS);CHKERRQ(ierr);
997   ierr = ISGetLocalSize(idIS, &numFS);CHKERRQ(ierr);
998   ierr = ISGetIndices(idIS, &ids);CHKERRQ(ierr);
999 
1000   user->numSplitFaces = 0;
1001   for (fs = 0; fs < numFS; ++fs) {
1002     PetscInt numBdFaces;
1003 
1004     ierr = DMGetStratumSize(dm, labelName, ids[fs], &numBdFaces);CHKERRQ(ierr);
1005     user->numSplitFaces += numBdFaces;
1006   }
1007   ierr  = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
1008   pEnd += user->numSplitFaces;
1009   ierr  = DMPlexSetChart(sdm, pStart, pEnd);CHKERRQ(ierr);
1010   ierr  = DMPlexGetGhostCellStratum(dm, &cEndInterior, NULL);CHKERRQ(ierr);
1011   ierr  = DMPlexGetHeightStratum(dm, 0, NULL, &cEnd);CHKERRQ(ierr);
1012   numGhostCells = cEnd - cEndInterior;
1013   /* Set cone and support sizes */
1014   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
1015   for (d = 0; d <= depth; ++d) {
1016     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
1017     for (p = pStart; p < pEnd; ++p) {
1018       PetscInt newp = p;
1019       PetscInt size;
1020 
1021       ierr = DMPlexGetConeSize(dm, p, &size);CHKERRQ(ierr);
1022       ierr = DMPlexSetConeSize(sdm, newp, size);CHKERRQ(ierr);
1023       ierr = DMPlexGetSupportSize(dm, p, &size);CHKERRQ(ierr);
1024       ierr = DMPlexSetSupportSize(sdm, newp, size);CHKERRQ(ierr);
1025     }
1026   }
1027   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
1028   for (fs = 0, newf = fEnd; fs < numFS; ++fs) {
1029     IS             faceIS;
1030     const PetscInt *faces;
1031     PetscInt       numFaces, f;
1032 
1033     ierr = DMGetStratumIS(dm, labelName, ids[fs], &faceIS);CHKERRQ(ierr);
1034     ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr);
1035     ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr);
1036     for (f = 0; f < numFaces; ++f, ++newf) {
1037       PetscInt size;
1038 
1039       /* Right now I think that both faces should see both cells */
1040       ierr = DMPlexGetConeSize(dm, faces[f], &size);CHKERRQ(ierr);
1041       ierr = DMPlexSetConeSize(sdm, newf, size);CHKERRQ(ierr);
1042       ierr = DMPlexGetSupportSize(dm, faces[f], &size);CHKERRQ(ierr);
1043       ierr = DMPlexSetSupportSize(sdm, newf, size);CHKERRQ(ierr);
1044     }
1045     ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
1046     ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
1047   }
1048   ierr = DMSetUp(sdm);CHKERRQ(ierr);
1049   /* Set cones and supports */
1050   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr);
1051   ierr = PetscMalloc1(PetscMax(maxConeSize, maxSupportSize), &newpoints);CHKERRQ(ierr);
1052   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
1053   for (p = pStart; p < pEnd; ++p) {
1054     const PetscInt *points, *orientations;
1055     PetscInt       size, i, newp = p;
1056 
1057     ierr = DMPlexGetConeSize(dm, p, &size);CHKERRQ(ierr);
1058     ierr = DMPlexGetCone(dm, p, &points);CHKERRQ(ierr);
1059     ierr = DMPlexGetConeOrientation(dm, p, &orientations);CHKERRQ(ierr);
1060     for (i = 0; i < size; ++i) newpoints[i] = points[i];
1061     ierr = DMPlexSetCone(sdm, newp, newpoints);CHKERRQ(ierr);
1062     ierr = DMPlexSetConeOrientation(sdm, newp, orientations);CHKERRQ(ierr);
1063     ierr = DMPlexGetSupportSize(dm, p, &size);CHKERRQ(ierr);
1064     ierr = DMPlexGetSupport(dm, p, &points);CHKERRQ(ierr);
1065     for (i = 0; i < size; ++i) newpoints[i] = points[i];
1066     ierr = DMPlexSetSupport(sdm, newp, newpoints);CHKERRQ(ierr);
1067   }
1068   ierr = PetscFree(newpoints);CHKERRQ(ierr);
1069   for (fs = 0, newf = fEnd; fs < numFS; ++fs) {
1070     IS             faceIS;
1071     const PetscInt *faces;
1072     PetscInt       numFaces, f;
1073 
1074     ierr = DMGetStratumIS(dm, labelName, ids[fs], &faceIS);CHKERRQ(ierr);
1075     ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr);
1076     ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr);
1077     for (f = 0; f < numFaces; ++f, ++newf) {
1078       const PetscInt *points;
1079 
1080       ierr = DMPlexGetCone(dm, faces[f], &points);CHKERRQ(ierr);
1081       ierr = DMPlexSetCone(sdm, newf, points);CHKERRQ(ierr);
1082       ierr = DMPlexGetSupport(dm, faces[f], &points);CHKERRQ(ierr);
1083       ierr = DMPlexSetSupport(sdm, newf, points);CHKERRQ(ierr);
1084     }
1085     ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
1086     ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
1087   }
1088   ierr = ISRestoreIndices(idIS, &ids);CHKERRQ(ierr);
1089   ierr = ISDestroy(&idIS);CHKERRQ(ierr);
1090   ierr = DMPlexStratify(sdm);CHKERRQ(ierr);
1091   /* Convert coordinates */
1092   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
1093   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1094   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &newCoordSection);CHKERRQ(ierr);
1095   ierr = PetscSectionSetNumFields(newCoordSection, 1);CHKERRQ(ierr);
1096   ierr = PetscSectionSetFieldComponents(newCoordSection, 0, dim);CHKERRQ(ierr);
1097   ierr = PetscSectionSetChart(newCoordSection, vStart, vEnd);CHKERRQ(ierr);
1098   for (v = vStart; v < vEnd; ++v) {
1099     ierr = PetscSectionSetDof(newCoordSection, v, dim);CHKERRQ(ierr);
1100     ierr = PetscSectionSetFieldDof(newCoordSection, v, 0, dim);CHKERRQ(ierr);
1101   }
1102   ierr = PetscSectionSetUp(newCoordSection);CHKERRQ(ierr);
1103   ierr = DMSetCoordinateSection(sdm, PETSC_DETERMINE, newCoordSection);CHKERRQ(ierr);
1104   ierr = PetscSectionDestroy(&newCoordSection);CHKERRQ(ierr); /* relinquish our reference */
1105   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1106   ierr = DMSetCoordinatesLocal(sdm, coordinates);CHKERRQ(ierr);
1107   /* Convert labels */
1108   ierr = DMGetNumLabels(dm, &numLabels);CHKERRQ(ierr);
1109   for (l = 0; l < numLabels; ++l) {
1110     const char *lname;
1111     PetscBool  isDepth, isDim;
1112 
1113     ierr = DMGetLabelName(dm, l, &lname);CHKERRQ(ierr);
1114     ierr = PetscStrcmp(lname, "depth", &isDepth);CHKERRQ(ierr);
1115     if (isDepth) continue;
1116     ierr = PetscStrcmp(lname, "dim", &isDim);CHKERRQ(ierr);
1117     if (isDim) continue;
1118     ierr = DMCreateLabel(sdm, lname);CHKERRQ(ierr);
1119     ierr = DMGetLabelIdIS(dm, lname, &idIS);CHKERRQ(ierr);
1120     ierr = ISGetLocalSize(idIS, &numFS);CHKERRQ(ierr);
1121     ierr = ISGetIndices(idIS, &ids);CHKERRQ(ierr);
1122     for (fs = 0; fs < numFS; ++fs) {
1123       IS             pointIS;
1124       const PetscInt *points;
1125       PetscInt       numPoints;
1126 
1127       ierr = DMGetStratumIS(dm, lname, ids[fs], &pointIS);CHKERRQ(ierr);
1128       ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr);
1129       ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
1130       for (p = 0; p < numPoints; ++p) {
1131         PetscInt newpoint = points[p];
1132 
1133         ierr = DMSetLabelValue(sdm, lname, newpoint, ids[fs]);CHKERRQ(ierr);
1134       }
1135       ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
1136       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1137     }
1138     ierr = ISRestoreIndices(idIS, &ids);CHKERRQ(ierr);
1139     ierr = ISDestroy(&idIS);CHKERRQ(ierr);
1140   }
1141   {
1142     /* Convert pointSF */
1143     const PetscSFNode *remotePoints;
1144     PetscSFNode       *gremotePoints;
1145     const PetscInt    *localPoints;
1146     PetscInt          *glocalPoints,*newLocation,*newRemoteLocation;
1147     PetscInt          numRoots, numLeaves;
1148     PetscMPIInt       size;
1149 
1150     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size);CHKERRMPI(ierr);
1151     ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
1152     ierr = DMGetPointSF(sdm, &gsfPoint);CHKERRQ(ierr);
1153     ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
1154     ierr = PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr);
1155     if (numRoots >= 0) {
1156       ierr = PetscMalloc2(numRoots,&newLocation,pEnd-pStart,&newRemoteLocation);CHKERRQ(ierr);
1157       for (l=0; l<numRoots; l++) newLocation[l] = l; /* + (l >= cEnd ? numGhostCells : 0); */
1158       ierr = PetscSFBcastBegin(sfPoint, MPIU_INT, newLocation, newRemoteLocation,MPI_REPLACE);CHKERRQ(ierr);
1159       ierr = PetscSFBcastEnd(sfPoint, MPIU_INT, newLocation, newRemoteLocation,MPI_REPLACE);CHKERRQ(ierr);
1160       ierr = PetscMalloc1(numLeaves,    &glocalPoints);CHKERRQ(ierr);
1161       ierr = PetscMalloc1(numLeaves, &gremotePoints);CHKERRQ(ierr);
1162       for (l = 0; l < numLeaves; ++l) {
1163         glocalPoints[l]        = localPoints[l]; /* localPoints[l] >= cEnd ? localPoints[l] + numGhostCells : localPoints[l]; */
1164         gremotePoints[l].rank  = remotePoints[l].rank;
1165         gremotePoints[l].index = newRemoteLocation[localPoints[l]];
1166       }
1167       ierr = PetscFree2(newLocation,newRemoteLocation);CHKERRQ(ierr);
1168       ierr = PetscSFSetGraph(gsfPoint, numRoots+numGhostCells, numLeaves, glocalPoints, PETSC_OWN_POINTER, gremotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
1169     }
1170     ierr     = DMDestroy(dmSplit);CHKERRQ(ierr);
1171     *dmSplit = sdm;
1172   }
1173   PetscFunctionReturn(0);
1174 }
1175 
1176 PetscErrorCode CreatePartitionVec(DM dm, DM *dmCell, Vec *partition)
1177 {
1178   PetscSF        sfPoint;
1179   PetscSection   coordSection;
1180   Vec            coordinates;
1181   PetscSection   sectionCell;
1182   PetscScalar    *part;
1183   PetscInt       cStart, cEnd, c;
1184   PetscMPIInt    rank;
1185   PetscErrorCode ierr;
1186 
1187   PetscFunctionBeginUser;
1188   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1189   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1190   ierr = DMClone(dm, dmCell);CHKERRQ(ierr);
1191   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
1192   ierr = DMSetPointSF(*dmCell, sfPoint);CHKERRQ(ierr);
1193   ierr = DMSetCoordinateSection(*dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr);
1194   ierr = DMSetCoordinatesLocal(*dmCell, coordinates);CHKERRQ(ierr);
1195   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRMPI(ierr);
1196   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &sectionCell);CHKERRQ(ierr);
1197   ierr = DMPlexGetHeightStratum(*dmCell, 0, &cStart, &cEnd);CHKERRQ(ierr);
1198   ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr);
1199   for (c = cStart; c < cEnd; ++c) {
1200     ierr = PetscSectionSetDof(sectionCell, c, 1);CHKERRQ(ierr);
1201   }
1202   ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr);
1203   ierr = DMSetLocalSection(*dmCell, sectionCell);CHKERRQ(ierr);
1204   ierr = PetscSectionDestroy(&sectionCell);CHKERRQ(ierr);
1205   ierr = DMCreateLocalVector(*dmCell, partition);CHKERRQ(ierr);
1206   ierr = PetscObjectSetName((PetscObject)*partition, "partition");CHKERRQ(ierr);
1207   ierr = VecGetArray(*partition, &part);CHKERRQ(ierr);
1208   for (c = cStart; c < cEnd; ++c) {
1209     PetscScalar *p;
1210 
1211     ierr = DMPlexPointLocalRef(*dmCell, c, part, &p);CHKERRQ(ierr);
1212     p[0] = rank;
1213   }
1214   ierr = VecRestoreArray(*partition, &part);CHKERRQ(ierr);
1215   PetscFunctionReturn(0);
1216 }
1217 
1218 PetscErrorCode CreateMassMatrix(DM dm, Vec *massMatrix, User user)
1219 {
1220   DM                plex, dmMass, dmFace, dmCell, dmCoord;
1221   PetscSection      coordSection;
1222   Vec               coordinates, facegeom, cellgeom;
1223   PetscSection      sectionMass;
1224   PetscScalar       *m;
1225   const PetscScalar *fgeom, *cgeom, *coords;
1226   PetscInt          vStart, vEnd, v;
1227   PetscErrorCode    ierr;
1228 
1229   PetscFunctionBeginUser;
1230   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
1231   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1232   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1233   ierr = DMClone(dm, &dmMass);CHKERRQ(ierr);
1234   ierr = DMSetCoordinateSection(dmMass, PETSC_DETERMINE, coordSection);CHKERRQ(ierr);
1235   ierr = DMSetCoordinatesLocal(dmMass, coordinates);CHKERRQ(ierr);
1236   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &sectionMass);CHKERRQ(ierr);
1237   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
1238   ierr = PetscSectionSetChart(sectionMass, vStart, vEnd);CHKERRQ(ierr);
1239   for (v = vStart; v < vEnd; ++v) {
1240     PetscInt numFaces;
1241 
1242     ierr = DMPlexGetSupportSize(dmMass, v, &numFaces);CHKERRQ(ierr);
1243     ierr = PetscSectionSetDof(sectionMass, v, numFaces*numFaces);CHKERRQ(ierr);
1244   }
1245   ierr = PetscSectionSetUp(sectionMass);CHKERRQ(ierr);
1246   ierr = DMSetLocalSection(dmMass, sectionMass);CHKERRQ(ierr);
1247   ierr = PetscSectionDestroy(&sectionMass);CHKERRQ(ierr);
1248   ierr = DMGetLocalVector(dmMass, massMatrix);CHKERRQ(ierr);
1249   ierr = VecGetArray(*massMatrix, &m);CHKERRQ(ierr);
1250   ierr = DMPlexGetGeometryFVM(plex, &facegeom, &cellgeom, NULL);CHKERRQ(ierr);
1251   ierr = VecGetDM(facegeom, &dmFace);CHKERRQ(ierr);
1252   ierr = VecGetArrayRead(facegeom, &fgeom);CHKERRQ(ierr);
1253   ierr = VecGetDM(cellgeom, &dmCell);CHKERRQ(ierr);
1254   ierr = VecGetArrayRead(cellgeom, &cgeom);CHKERRQ(ierr);
1255   ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr);
1256   ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
1257   for (v = vStart; v < vEnd; ++v) {
1258     const PetscInt        *faces;
1259     PetscFVFaceGeom       *fgA, *fgB, *cg;
1260     PetscScalar           *vertex;
1261     PetscInt               numFaces, sides[2], f, g;
1262 
1263     ierr = DMPlexPointLocalRead(dmCoord, v, coords, &vertex);CHKERRQ(ierr);
1264     ierr = DMPlexGetSupportSize(dmMass, v, &numFaces);CHKERRQ(ierr);
1265     ierr = DMPlexGetSupport(dmMass, v, &faces);CHKERRQ(ierr);
1266     for (f = 0; f < numFaces; ++f) {
1267       sides[0] = faces[f];
1268       ierr = DMPlexPointLocalRead(dmFace, faces[f], fgeom, &fgA);CHKERRQ(ierr);
1269       for (g = 0; g < numFaces; ++g) {
1270         const PetscInt *cells = NULL;
1271         PetscReal      area   = 0.0;
1272         PetscInt       numCells;
1273 
1274         sides[1] = faces[g];
1275         ierr = DMPlexPointLocalRead(dmFace, faces[g], fgeom, &fgB);CHKERRQ(ierr);
1276         ierr = DMPlexGetJoin(dmMass, 2, sides, &numCells, &cells);CHKERRQ(ierr);
1277         if (numCells != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Invalid join for faces");
1278         ierr = DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cg);CHKERRQ(ierr);
1279         area += PetscAbsScalar((vertex[0] - cg->centroid[0])*(fgA->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1])*(fgA->centroid[0] - cg->centroid[0]));
1280         area += PetscAbsScalar((vertex[0] - cg->centroid[0])*(fgB->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1])*(fgB->centroid[0] - cg->centroid[0]));
1281         m[f*numFaces+g] = Dot2Real(fgA->normal, fgB->normal)*area*0.5;
1282         ierr = DMPlexRestoreJoin(dmMass, 2, sides, &numCells, &cells);CHKERRQ(ierr);
1283       }
1284     }
1285   }
1286   ierr = VecRestoreArrayRead(facegeom, &fgeom);CHKERRQ(ierr);
1287   ierr = VecRestoreArrayRead(cellgeom, &cgeom);CHKERRQ(ierr);
1288   ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
1289   ierr = VecRestoreArray(*massMatrix, &m);CHKERRQ(ierr);
1290   ierr = DMDestroy(&dmMass);CHKERRQ(ierr);
1291   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1292   PetscFunctionReturn(0);
1293 }
1294 
1295 /* Behavior will be different for multi-physics or when using non-default boundary conditions */
1296 static PetscErrorCode ModelSolutionSetDefault(Model mod,SolutionFunction func,void *ctx)
1297 {
1298   PetscFunctionBeginUser;
1299   mod->solution    = func;
1300   mod->solutionctx = ctx;
1301   PetscFunctionReturn(0);
1302 }
1303 
1304 static PetscErrorCode ModelFunctionalRegister(Model mod,const char *name,PetscInt *offset,FunctionalFunction func,void *ctx)
1305 {
1306   PetscErrorCode ierr;
1307   FunctionalLink link,*ptr;
1308   PetscInt       lastoffset = -1;
1309 
1310   PetscFunctionBeginUser;
1311   for (ptr=&mod->functionalRegistry; *ptr; ptr = &(*ptr)->next) lastoffset = (*ptr)->offset;
1312   ierr         = PetscNew(&link);CHKERRQ(ierr);
1313   ierr         = PetscStrallocpy(name,&link->name);CHKERRQ(ierr);
1314   link->offset = lastoffset + 1;
1315   link->func   = func;
1316   link->ctx    = ctx;
1317   link->next   = NULL;
1318   *ptr         = link;
1319   *offset      = link->offset;
1320   PetscFunctionReturn(0);
1321 }
1322 
1323 static PetscErrorCode ModelFunctionalSetFromOptions(Model mod,PetscOptionItems *PetscOptionsObject)
1324 {
1325   PetscErrorCode ierr;
1326   PetscInt       i,j;
1327   FunctionalLink link;
1328   char           *names[256];
1329 
1330   PetscFunctionBeginUser;
1331   mod->numMonitored = ALEN(names);
1332   ierr = PetscOptionsStringArray("-monitor","list of functionals to monitor","",names,&mod->numMonitored,NULL);CHKERRQ(ierr);
1333   /* Create list of functionals that will be computed somehow */
1334   ierr = PetscMalloc1(mod->numMonitored,&mod->functionalMonitored);CHKERRQ(ierr);
1335   /* Create index of calls that we will have to make to compute these functionals (over-allocation in general). */
1336   ierr = PetscMalloc1(mod->numMonitored,&mod->functionalCall);CHKERRQ(ierr);
1337   mod->numCall = 0;
1338   for (i=0; i<mod->numMonitored; i++) {
1339     for (link=mod->functionalRegistry; link; link=link->next) {
1340       PetscBool match;
1341       ierr = PetscStrcasecmp(names[i],link->name,&match);CHKERRQ(ierr);
1342       if (match) break;
1343     }
1344     if (!link) SETERRQ1(mod->comm,PETSC_ERR_USER,"No known functional '%s'",names[i]);
1345     mod->functionalMonitored[i] = link;
1346     for (j=0; j<i; j++) {
1347       if (mod->functionalCall[j]->func == link->func && mod->functionalCall[j]->ctx == link->ctx) goto next_name;
1348     }
1349     mod->functionalCall[mod->numCall++] = link; /* Just points to the first link using the result. There may be more results. */
1350 next_name:
1351     ierr = PetscFree(names[i]);CHKERRQ(ierr);
1352   }
1353 
1354   /* Find out the maximum index of any functional computed by a function we will be calling (even if we are not using it) */
1355   mod->maxComputed = -1;
1356   for (link=mod->functionalRegistry; link; link=link->next) {
1357     for (i=0; i<mod->numCall; i++) {
1358       FunctionalLink call = mod->functionalCall[i];
1359       if (link->func == call->func && link->ctx == call->ctx) {
1360         mod->maxComputed = PetscMax(mod->maxComputed,link->offset);
1361       }
1362     }
1363   }
1364   PetscFunctionReturn(0);
1365 }
1366 
1367 static PetscErrorCode FunctionalLinkDestroy(FunctionalLink *link)
1368 {
1369   PetscErrorCode ierr;
1370   FunctionalLink l,next;
1371 
1372   PetscFunctionBeginUser;
1373   if (!link) PetscFunctionReturn(0);
1374   l     = *link;
1375   *link = NULL;
1376   for (; l; l=next) {
1377     next = l->next;
1378     ierr = PetscFree(l->name);CHKERRQ(ierr);
1379     ierr = PetscFree(l);CHKERRQ(ierr);
1380   }
1381   PetscFunctionReturn(0);
1382 }
1383 
1384 /* put the solution callback into a functional callback */
1385 static PetscErrorCode SolutionFunctional(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *modctx)
1386 {
1387   Model          mod;
1388   PetscErrorCode ierr;
1389   PetscFunctionBegin;
1390   mod  = (Model) modctx;
1391   ierr = (*mod->solution)(mod, time, x, u, mod->solutionctx);CHKERRQ(ierr);
1392   PetscFunctionReturn(0);
1393 }
1394 
1395 PetscErrorCode SetInitialCondition(DM dm, Vec X, User user)
1396 {
1397   PetscErrorCode     (*func[1]) (PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
1398   void               *ctx[1];
1399   Model              mod = user->model;
1400   PetscErrorCode     ierr;
1401 
1402   PetscFunctionBeginUser;
1403   func[0] = SolutionFunctional;
1404   ctx[0]  = (void *) mod;
1405   ierr    = DMProjectFunction(dm,0.0,func,ctx,INSERT_ALL_VALUES,X);CHKERRQ(ierr);
1406   PetscFunctionReturn(0);
1407 }
1408 
1409 static PetscErrorCode OutputVTK(DM dm, const char *filename, PetscViewer *viewer)
1410 {
1411   PetscErrorCode ierr;
1412 
1413   PetscFunctionBeginUser;
1414   ierr = PetscViewerCreate(PetscObjectComm((PetscObject)dm), viewer);CHKERRQ(ierr);
1415   ierr = PetscViewerSetType(*viewer, PETSCVIEWERVTK);CHKERRQ(ierr);
1416   ierr = PetscViewerFileSetName(*viewer, filename);CHKERRQ(ierr);
1417   PetscFunctionReturn(0);
1418 }
1419 
1420 static PetscErrorCode MonitorVTK(TS ts,PetscInt stepnum,PetscReal time,Vec X,void *ctx)
1421 {
1422   User           user = (User)ctx;
1423   DM             dm, plex;
1424   PetscViewer    viewer;
1425   char           filename[PETSC_MAX_PATH_LEN],*ftable = NULL;
1426   PetscReal      xnorm;
1427   PetscErrorCode ierr;
1428 
1429   PetscFunctionBeginUser;
1430   ierr = PetscObjectSetName((PetscObject) X, "u");CHKERRQ(ierr);
1431   ierr = VecGetDM(X,&dm);CHKERRQ(ierr);
1432   ierr = VecNorm(X,NORM_INFINITY,&xnorm);CHKERRQ(ierr);
1433 
1434   if (stepnum >= 0) {
1435     stepnum += user->monitorStepOffset;
1436   }
1437   if (stepnum >= 0) {           /* No summary for final time */
1438     Model             mod = user->model;
1439     Vec               cellgeom;
1440     PetscInt          c,cStart,cEnd,fcount,i;
1441     size_t            ftableused,ftablealloc;
1442     const PetscScalar *cgeom,*x;
1443     DM                dmCell;
1444     DMLabel           vtkLabel;
1445     PetscReal         *fmin,*fmax,*fintegral,*ftmp;
1446 
1447     ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
1448     ierr = DMPlexGetGeometryFVM(plex, NULL, &cellgeom, NULL);CHKERRQ(ierr);
1449     fcount = mod->maxComputed+1;
1450     ierr   = PetscMalloc4(fcount,&fmin,fcount,&fmax,fcount,&fintegral,fcount,&ftmp);CHKERRQ(ierr);
1451     for (i=0; i<fcount; i++) {
1452       fmin[i]      = PETSC_MAX_REAL;
1453       fmax[i]      = PETSC_MIN_REAL;
1454       fintegral[i] = 0;
1455     }
1456     ierr = VecGetDM(cellgeom,&dmCell);CHKERRQ(ierr);
1457     ierr = DMPlexGetSimplexOrBoxCells(dmCell,0,&cStart,&cEnd);CHKERRQ(ierr);
1458     ierr = VecGetArrayRead(cellgeom,&cgeom);CHKERRQ(ierr);
1459     ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
1460     ierr = DMGetLabel(dm,"vtk",&vtkLabel);CHKERRQ(ierr);
1461     for (c = cStart; c < cEnd; ++c) {
1462       PetscFVCellGeom       *cg;
1463       const PetscScalar     *cx    = NULL;
1464       PetscInt              vtkVal = 0;
1465 
1466       /* not that these two routines as currently implemented work for any dm with a
1467        * localSection/globalSection */
1468       ierr = DMPlexPointLocalRead(dmCell,c,cgeom,&cg);CHKERRQ(ierr);
1469       ierr = DMPlexPointGlobalRead(dm,c,x,&cx);CHKERRQ(ierr);
1470       if (vtkLabel) {ierr = DMLabelGetValue(vtkLabel,c,&vtkVal);CHKERRQ(ierr);}
1471       if (!vtkVal || !cx) continue;        /* ghost, or not a global cell */
1472       for (i=0; i<mod->numCall; i++) {
1473         FunctionalLink flink = mod->functionalCall[i];
1474         ierr = (*flink->func)(mod,time,cg->centroid,cx,ftmp,flink->ctx);CHKERRQ(ierr);
1475       }
1476       for (i=0; i<fcount; i++) {
1477         fmin[i]       = PetscMin(fmin[i],ftmp[i]);
1478         fmax[i]       = PetscMax(fmax[i],ftmp[i]);
1479         fintegral[i] += cg->volume * ftmp[i];
1480       }
1481     }
1482     ierr = VecRestoreArrayRead(cellgeom,&cgeom);CHKERRQ(ierr);
1483     ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
1484     ierr = DMDestroy(&plex);CHKERRQ(ierr);
1485     ierr = MPI_Allreduce(MPI_IN_PLACE,fmin,fcount,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));CHKERRMPI(ierr);
1486     ierr = MPI_Allreduce(MPI_IN_PLACE,fmax,fcount,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));CHKERRMPI(ierr);
1487     ierr = MPI_Allreduce(MPI_IN_PLACE,fintegral,fcount,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)ts));CHKERRMPI(ierr);
1488 
1489     ftablealloc = fcount * 100;
1490     ftableused  = 0;
1491     ierr        = PetscMalloc1(ftablealloc,&ftable);CHKERRQ(ierr);
1492     for (i=0; i<mod->numMonitored; i++) {
1493       size_t         countused;
1494       char           buffer[256],*p;
1495       FunctionalLink flink = mod->functionalMonitored[i];
1496       PetscInt       id    = flink->offset;
1497       if (i % 3) {
1498         ierr = PetscArraycpy(buffer,"  ",2);CHKERRQ(ierr);
1499         p    = buffer + 2;
1500       } else if (i) {
1501         char newline[] = "\n";
1502         ierr = PetscMemcpy(buffer,newline,sizeof(newline)-1);CHKERRQ(ierr);
1503         p    = buffer + sizeof(newline) - 1;
1504       } else {
1505         p = buffer;
1506       }
1507       ierr = 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]);CHKERRQ(ierr);
1508       countused--;
1509       countused += p - buffer;
1510       if (countused > ftablealloc-ftableused-1) { /* reallocate */
1511         char *ftablenew;
1512         ftablealloc = 2*ftablealloc + countused;
1513         ierr = PetscMalloc(ftablealloc,&ftablenew);CHKERRQ(ierr);
1514         ierr = PetscArraycpy(ftablenew,ftable,ftableused);CHKERRQ(ierr);
1515         ierr = PetscFree(ftable);CHKERRQ(ierr);
1516         ftable = ftablenew;
1517       }
1518       ierr = PetscArraycpy(ftable+ftableused,buffer,countused);CHKERRQ(ierr);
1519       ftableused += countused;
1520       ftable[ftableused] = 0;
1521     }
1522     ierr = PetscFree4(fmin,fmax,fintegral,ftmp);CHKERRQ(ierr);
1523 
1524     ierr = PetscPrintf(PetscObjectComm((PetscObject)ts),"% 3D  time %8.4g  |x| %8.4g  %s\n",stepnum,(double)time,(double)xnorm,ftable ? ftable : "");CHKERRQ(ierr);
1525     ierr = PetscFree(ftable);CHKERRQ(ierr);
1526   }
1527   if (user->vtkInterval < 1) PetscFunctionReturn(0);
1528   if ((stepnum == -1) ^ (stepnum % user->vtkInterval == 0)) {
1529     if (stepnum == -1) {        /* Final time is not multiple of normal time interval, write it anyway */
1530       ierr = TSGetStepNumber(ts,&stepnum);CHKERRQ(ierr);
1531     }
1532     ierr = PetscSNPrintf(filename,sizeof filename,"%s-%03D.vtu",user->outputBasename,stepnum);CHKERRQ(ierr);
1533     ierr = OutputVTK(dm,filename,&viewer);CHKERRQ(ierr);
1534     ierr = VecView(X,viewer);CHKERRQ(ierr);
1535     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1536   }
1537   PetscFunctionReturn(0);
1538 }
1539 
1540 static PetscErrorCode initializeTS(DM dm, User user, TS *ts)
1541 {
1542   PetscErrorCode ierr;
1543 
1544   PetscFunctionBegin;
1545   ierr = TSCreate(PetscObjectComm((PetscObject)dm), ts);CHKERRQ(ierr);
1546   ierr = TSSetType(*ts, TSSSP);CHKERRQ(ierr);
1547   ierr = TSSetDM(*ts, dm);CHKERRQ(ierr);
1548   if (user->vtkmon) {
1549     ierr = TSMonitorSet(*ts,MonitorVTK,user,NULL);CHKERRQ(ierr);
1550   }
1551   ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, user);CHKERRQ(ierr);
1552   ierr = DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, user);CHKERRQ(ierr);
1553   ierr = TSSetMaxTime(*ts,2.0);CHKERRQ(ierr);
1554   ierr = TSSetExactFinalTime(*ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
1555   PetscFunctionReturn(0);
1556 }
1557 
1558 static PetscErrorCode adaptToleranceFVM(PetscFV fvm, TS ts, Vec sol, VecTagger refineTag, VecTagger coarsenTag, User user, TS *tsNew, Vec *solNew)
1559 {
1560   DM                dm, gradDM, plex, cellDM, adaptedDM = NULL;
1561   Vec               cellGeom, faceGeom;
1562   PetscBool         isForest, computeGradient;
1563   Vec               grad, locGrad, locX, errVec;
1564   PetscInt          cStart, cEnd, c, dim, nRefine, nCoarsen;
1565   PetscReal         minMaxInd[2] = {PETSC_MAX_REAL, PETSC_MIN_REAL}, minMaxIndGlobal[2], minInd, maxInd, time;
1566   PetscScalar       *errArray;
1567   const PetscScalar *pointVals;
1568   const PetscScalar *pointGrads;
1569   const PetscScalar *pointGeom;
1570   DMLabel           adaptLabel = NULL;
1571   IS                refineIS, coarsenIS;
1572   PetscErrorCode    ierr;
1573 
1574   PetscFunctionBegin;
1575   ierr = TSGetTime(ts,&time);CHKERRQ(ierr);
1576   ierr = VecGetDM(sol, &dm);CHKERRQ(ierr);
1577   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
1578   ierr = PetscFVGetComputeGradients(fvm,&computeGradient);CHKERRQ(ierr);
1579   ierr = PetscFVSetComputeGradients(fvm,PETSC_TRUE);CHKERRQ(ierr);
1580   ierr = DMIsForest(dm, &isForest);CHKERRQ(ierr);
1581   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
1582   ierr = DMPlexGetDataFVM(plex, fvm, &cellGeom, &faceGeom, &gradDM);CHKERRQ(ierr);
1583   ierr = DMCreateLocalVector(plex,&locX);CHKERRQ(ierr);
1584   ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, 0.0, faceGeom, cellGeom, NULL);CHKERRQ(ierr);
1585   ierr = DMGlobalToLocalBegin(plex, sol, INSERT_VALUES, locX);CHKERRQ(ierr);
1586   ierr = DMGlobalToLocalEnd  (plex, sol, INSERT_VALUES, locX);CHKERRQ(ierr);
1587   ierr = DMCreateGlobalVector(gradDM, &grad);CHKERRQ(ierr);
1588   ierr = DMPlexReconstructGradientsFVM(plex, locX, grad);CHKERRQ(ierr);
1589   ierr = DMCreateLocalVector(gradDM, &locGrad);CHKERRQ(ierr);
1590   ierr = DMGlobalToLocalBegin(gradDM, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr);
1591   ierr = DMGlobalToLocalEnd(gradDM, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr);
1592   ierr = VecDestroy(&grad);CHKERRQ(ierr);
1593   ierr = DMPlexGetSimplexOrBoxCells(plex,0,&cStart,&cEnd);CHKERRQ(ierr);
1594   ierr = VecGetArrayRead(locGrad,&pointGrads);CHKERRQ(ierr);
1595   ierr = VecGetArrayRead(cellGeom,&pointGeom);CHKERRQ(ierr);
1596   ierr = VecGetArrayRead(locX,&pointVals);CHKERRQ(ierr);
1597   ierr = VecGetDM(cellGeom,&cellDM);CHKERRQ(ierr);
1598   ierr = DMLabelCreate(PETSC_COMM_SELF,"adapt",&adaptLabel);CHKERRQ(ierr);
1599   ierr = VecCreateMPI(PetscObjectComm((PetscObject)plex),cEnd-cStart,PETSC_DETERMINE,&errVec);CHKERRQ(ierr);
1600   ierr = VecSetUp(errVec);CHKERRQ(ierr);
1601   ierr = VecGetArray(errVec,&errArray);CHKERRQ(ierr);
1602   for (c = cStart; c < cEnd; c++) {
1603     PetscReal             errInd = 0.;
1604     PetscScalar           *pointGrad;
1605     PetscScalar           *pointVal;
1606     PetscFVCellGeom       *cg;
1607 
1608     ierr = DMPlexPointLocalRead(gradDM,c,pointGrads,&pointGrad);CHKERRQ(ierr);
1609     ierr = DMPlexPointLocalRead(cellDM,c,pointGeom,&cg);CHKERRQ(ierr);
1610     ierr = DMPlexPointLocalRead(plex,c,pointVals,&pointVal);CHKERRQ(ierr);
1611 
1612     ierr = (user->model->errorIndicator)(dim,cg->volume,user->model->physics->dof,pointVal,pointGrad,&errInd,user->model->errorCtx);CHKERRQ(ierr);
1613     errArray[c-cStart] = errInd;
1614     minMaxInd[0] = PetscMin(minMaxInd[0],errInd);
1615     minMaxInd[1] = PetscMax(minMaxInd[1],errInd);
1616   }
1617   ierr = VecRestoreArray(errVec,&errArray);CHKERRQ(ierr);
1618   ierr = VecRestoreArrayRead(locX,&pointVals);CHKERRQ(ierr);
1619   ierr = VecRestoreArrayRead(cellGeom,&pointGeom);CHKERRQ(ierr);
1620   ierr = VecRestoreArrayRead(locGrad,&pointGrads);CHKERRQ(ierr);
1621   ierr = VecDestroy(&locGrad);CHKERRQ(ierr);
1622   ierr = VecDestroy(&locX);CHKERRQ(ierr);
1623   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1624 
1625   ierr = VecTaggerComputeIS(refineTag,errVec,&refineIS);CHKERRQ(ierr);
1626   ierr = VecTaggerComputeIS(coarsenTag,errVec,&coarsenIS);CHKERRQ(ierr);
1627   ierr = ISGetSize(refineIS,&nRefine);CHKERRQ(ierr);
1628   ierr = ISGetSize(coarsenIS,&nCoarsen);CHKERRQ(ierr);
1629   if (nRefine) {ierr = DMLabelSetStratumIS(adaptLabel,DM_ADAPT_REFINE,refineIS);CHKERRQ(ierr);}
1630   if (nCoarsen) {ierr = DMLabelSetStratumIS(adaptLabel,DM_ADAPT_COARSEN,coarsenIS);CHKERRQ(ierr);}
1631   ierr = ISDestroy(&coarsenIS);CHKERRQ(ierr);
1632   ierr = ISDestroy(&refineIS);CHKERRQ(ierr);
1633   ierr = VecDestroy(&errVec);CHKERRQ(ierr);
1634 
1635   ierr = PetscFVSetComputeGradients(fvm,computeGradient);CHKERRQ(ierr);
1636   minMaxInd[1] = -minMaxInd[1];
1637   ierr = MPI_Allreduce(minMaxInd,minMaxIndGlobal,2,MPIU_REAL,MPI_MIN,PetscObjectComm((PetscObject)dm));CHKERRMPI(ierr);
1638   minInd = minMaxIndGlobal[0];
1639   maxInd = -minMaxIndGlobal[1];
1640   ierr = PetscInfo2(ts, "error indicator range (%E, %E)\n", minInd, maxInd);CHKERRQ(ierr);
1641   if (nRefine || nCoarsen) { /* at least one cell is over the refinement threshold */
1642     ierr = DMAdaptLabel(dm,adaptLabel,&adaptedDM);CHKERRQ(ierr);
1643   }
1644   ierr = DMLabelDestroy(&adaptLabel);CHKERRQ(ierr);
1645   if (adaptedDM) {
1646     ierr = PetscInfo2(ts, "Adapted mesh, marking %D cells for refinement, and %D cells for coarsening\n", nRefine, nCoarsen);CHKERRQ(ierr);
1647     if (tsNew) {ierr = initializeTS(adaptedDM, user, tsNew);CHKERRQ(ierr);}
1648     if (solNew) {
1649       ierr = DMCreateGlobalVector(adaptedDM, solNew);CHKERRQ(ierr);
1650       ierr = PetscObjectSetName((PetscObject) *solNew, "solution");CHKERRQ(ierr);
1651       ierr = DMForestTransferVec(dm, sol, adaptedDM, *solNew, PETSC_TRUE, time);CHKERRQ(ierr);
1652     }
1653     if (isForest) {ierr = DMForestSetAdaptivityForest(adaptedDM,NULL);CHKERRQ(ierr);} /* clear internal references to the previous dm */
1654     ierr = DMDestroy(&adaptedDM);CHKERRQ(ierr);
1655   } else {
1656     if (tsNew)  *tsNew  = NULL;
1657     if (solNew) *solNew = NULL;
1658   }
1659   PetscFunctionReturn(0);
1660 }
1661 
1662 int main(int argc, char **argv)
1663 {
1664   MPI_Comm          comm;
1665   PetscDS           prob;
1666   PetscFV           fvm;
1667   PetscLimiter      limiter = NULL, noneLimiter = NULL;
1668   User              user;
1669   Model             mod;
1670   Physics           phys;
1671   DM                dm, plex;
1672   PetscReal         ftime, cfl, dt, minRadius;
1673   PetscInt          dim, nsteps;
1674   TS                ts;
1675   TSConvergedReason reason;
1676   Vec               X;
1677   PetscViewer       viewer;
1678   PetscBool         simplex = PETSC_FALSE, vtkCellGeom, splitFaces, useAMR;
1679   PetscInt          overlap, adaptInterval;
1680   char              filename[PETSC_MAX_PATH_LEN] = "";
1681   char              physname[256]  = "advect";
1682   VecTagger         refineTag = NULL, coarsenTag = NULL;
1683   PetscErrorCode    ierr;
1684 
1685   ierr = PetscInitialize(&argc, &argv, (char*) 0, help);if (ierr) return ierr;
1686   comm = PETSC_COMM_WORLD;
1687 
1688   ierr          = PetscNew(&user);CHKERRQ(ierr);
1689   ierr          = PetscNew(&user->model);CHKERRQ(ierr);
1690   ierr          = PetscNew(&user->model->physics);CHKERRQ(ierr);
1691   mod           = user->model;
1692   phys          = mod->physics;
1693   mod->comm     = comm;
1694   useAMR        = PETSC_FALSE;
1695   adaptInterval = 1;
1696 
1697   /* Register physical models to be available on the command line */
1698   ierr = PetscFunctionListAdd(&PhysicsList,"advect"          ,PhysicsCreate_Advect);CHKERRQ(ierr);
1699   ierr = PetscFunctionListAdd(&PhysicsList,"sw"              ,PhysicsCreate_SW);CHKERRQ(ierr);
1700   ierr = PetscFunctionListAdd(&PhysicsList,"euler"           ,PhysicsCreate_Euler);CHKERRQ(ierr);
1701 
1702 
1703   ierr = PetscOptionsBegin(comm,NULL,"Unstructured Finite Volume Mesh Options","");CHKERRQ(ierr);
1704   {
1705     cfl  = 0.9 * 4; /* default SSPRKS2 with s=5 stages is stable for CFL number s-1 */
1706     ierr = PetscOptionsReal("-ufv_cfl","CFL number per step","",cfl,&cfl,NULL);CHKERRQ(ierr);
1707     ierr = PetscOptionsString("-f","Exodus.II filename to read","",filename,filename,sizeof(filename),NULL);CHKERRQ(ierr);
1708     ierr = PetscOptionsBool("-simplex","Flag to use a simplex mesh","",simplex,&simplex,NULL);CHKERRQ(ierr);
1709     splitFaces = PETSC_FALSE;
1710     ierr = PetscOptionsBool("-ufv_split_faces","Split faces between cell sets","",splitFaces,&splitFaces,NULL);CHKERRQ(ierr);
1711     overlap = 1;
1712     ierr = PetscOptionsInt("-ufv_mesh_overlap","Number of cells to overlap partitions","",overlap,&overlap,NULL);CHKERRQ(ierr);
1713     user->vtkInterval = 1;
1714     ierr = PetscOptionsInt("-ufv_vtk_interval","VTK output interval (0 to disable)","",user->vtkInterval,&user->vtkInterval,NULL);CHKERRQ(ierr);
1715     user->vtkmon = PETSC_TRUE;
1716     ierr = PetscOptionsBool("-ufv_vtk_monitor","Use VTKMonitor routine","",user->vtkmon,&user->vtkmon,NULL);CHKERRQ(ierr);
1717     vtkCellGeom = PETSC_FALSE;
1718     ierr = PetscStrcpy(user->outputBasename, "ex11");CHKERRQ(ierr);
1719     ierr = PetscOptionsString("-ufv_vtk_basename","VTK output basename","",user->outputBasename,user->outputBasename,sizeof(user->outputBasename),NULL);CHKERRQ(ierr);
1720     ierr = PetscOptionsBool("-ufv_vtk_cellgeom","Write cell geometry (for debugging)","",vtkCellGeom,&vtkCellGeom,NULL);CHKERRQ(ierr);
1721     ierr = PetscOptionsBool("-ufv_use_amr","use local adaptive mesh refinement","",useAMR,&useAMR,NULL);CHKERRQ(ierr);
1722     ierr = PetscOptionsInt("-ufv_adapt_interval","time steps between AMR","",adaptInterval,&adaptInterval,NULL);CHKERRQ(ierr);
1723   }
1724   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1725 
1726   if (useAMR) {
1727     VecTaggerBox refineBox, coarsenBox;
1728 
1729     refineBox.min  = refineBox.max  = PETSC_MAX_REAL;
1730     coarsenBox.min = coarsenBox.max = PETSC_MIN_REAL;
1731 
1732     ierr = VecTaggerCreate(comm,&refineTag);CHKERRQ(ierr);
1733     ierr = PetscObjectSetOptionsPrefix((PetscObject)refineTag,"refine_");CHKERRQ(ierr);
1734     ierr = VecTaggerSetType(refineTag,VECTAGGERABSOLUTE);CHKERRQ(ierr);
1735     ierr = VecTaggerAbsoluteSetBox(refineTag,&refineBox);CHKERRQ(ierr);
1736     ierr = VecTaggerSetFromOptions(refineTag);CHKERRQ(ierr);
1737     ierr = VecTaggerSetUp(refineTag);CHKERRQ(ierr);
1738     ierr = PetscObjectViewFromOptions((PetscObject)refineTag,NULL,"-tag_view");CHKERRQ(ierr);
1739 
1740     ierr = VecTaggerCreate(comm,&coarsenTag);CHKERRQ(ierr);
1741     ierr = PetscObjectSetOptionsPrefix((PetscObject)coarsenTag,"coarsen_");CHKERRQ(ierr);
1742     ierr = VecTaggerSetType(coarsenTag,VECTAGGERABSOLUTE);CHKERRQ(ierr);
1743     ierr = VecTaggerAbsoluteSetBox(coarsenTag,&coarsenBox);CHKERRQ(ierr);
1744     ierr = VecTaggerSetFromOptions(coarsenTag);CHKERRQ(ierr);
1745     ierr = VecTaggerSetUp(coarsenTag);CHKERRQ(ierr);
1746     ierr = PetscObjectViewFromOptions((PetscObject)coarsenTag,NULL,"-tag_view");CHKERRQ(ierr);
1747   }
1748 
1749   ierr = PetscOptionsBegin(comm,NULL,"Unstructured Finite Volume Physics Options","");CHKERRQ(ierr);
1750   {
1751     PetscErrorCode (*physcreate)(Model,Physics,PetscOptionItems*);
1752     ierr = PetscOptionsFList("-physics","Physics module to solve","",PhysicsList,physname,physname,sizeof physname,NULL);CHKERRQ(ierr);
1753     ierr = PetscFunctionListFind(PhysicsList,physname,&physcreate);CHKERRQ(ierr);
1754     ierr = PetscMemzero(phys,sizeof(struct _n_Physics));CHKERRQ(ierr);
1755     ierr = (*physcreate)(mod,phys,PetscOptionsObject);CHKERRQ(ierr);
1756     /* Count number of fields and dofs */
1757     for (phys->nfields=0,phys->dof=0; phys->field_desc[phys->nfields].name; phys->nfields++) phys->dof += phys->field_desc[phys->nfields].dof;
1758     if (phys->dof <= 0) SETERRQ1(comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set dof",physname);
1759     ierr = ModelFunctionalSetFromOptions(mod,PetscOptionsObject);CHKERRQ(ierr);
1760   }
1761   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1762 
1763   /* Create mesh */
1764   {
1765     size_t len,i;
1766     for (i = 0; i < DIM; i++) { mod->bounds[2*i] = 0.; mod->bounds[2*i+1] = 1.;};
1767     ierr = PetscStrlen(filename,&len);CHKERRQ(ierr);
1768     dim = DIM;
1769     if (!len) { /* a null name means just do a hex box */
1770       PetscInt cells[3] = {1, 1, 1}; /* coarse mesh is one cell; refine from there */
1771       PetscBool flg1, flg2, skew = PETSC_FALSE;
1772       PetscInt nret1 = DIM;
1773       PetscInt nret2 = 2*DIM;
1774       ierr = PetscOptionsBegin(comm,NULL,"Rectangular mesh options","");CHKERRQ(ierr);
1775       ierr = PetscOptionsIntArray("-grid_size","number of cells in each direction","",cells,&nret1,&flg1);CHKERRQ(ierr);
1776       ierr = PetscOptionsRealArray("-grid_bounds","bounds of the mesh in each direction (i.e., x_min,x_max,y_min,y_max","",mod->bounds,&nret2,&flg2);CHKERRQ(ierr);
1777       ierr = PetscOptionsBool("-grid_skew_60","Skew grid for 60 degree shock mesh","",skew,&skew,NULL);CHKERRQ(ierr);
1778       ierr = PetscOptionsEnd();CHKERRQ(ierr);
1779       if (flg1) {
1780         dim = nret1;
1781         if (dim != DIM) SETERRQ1(comm,PETSC_ERR_ARG_SIZ,"Dim wrong size %D in -grid_size",dim);
1782       }
1783       ierr = DMPlexCreateBoxMesh(comm, dim, simplex, cells, NULL, NULL, mod->bcs, PETSC_TRUE, &dm);CHKERRQ(ierr);
1784       if (flg2) {
1785         PetscInt dimEmbed, i;
1786         PetscInt nCoords;
1787         PetscScalar *coords;
1788         Vec coordinates;
1789 
1790         ierr = DMGetCoordinatesLocal(dm,&coordinates);CHKERRQ(ierr);
1791         ierr = DMGetCoordinateDim(dm,&dimEmbed);CHKERRQ(ierr);
1792         ierr = VecGetLocalSize(coordinates,&nCoords);CHKERRQ(ierr);
1793         if (nCoords % dimEmbed) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Coordinate vector the wrong size");
1794         ierr = VecGetArray(coordinates,&coords);CHKERRQ(ierr);
1795         for (i = 0; i < nCoords; i += dimEmbed) {
1796           PetscInt j;
1797 
1798           PetscScalar *coord = &coords[i];
1799           for (j = 0; j < dimEmbed; j++) {
1800             coord[j] = mod->bounds[2 * j] + coord[j] * (mod->bounds[2 * j + 1] - mod->bounds[2 * j]);
1801             if (dim==2 && cells[1]==1 && j==0 && skew) {
1802               if (cells[0]==2 && i==8) {
1803                 coord[j] = .57735026918963; /* hack to get 60 deg skewed mesh */
1804               }
1805               else if (cells[0]==3) {
1806                 if (i==2 || i==10) coord[j] = mod->bounds[1]/4.;
1807                 else if (i==4) coord[j] = mod->bounds[1]/2.;
1808                 else if (i==12) coord[j] = 1.57735026918963*mod->bounds[1]/2.;
1809               }
1810             }
1811           }
1812         }
1813         ierr = VecRestoreArray(coordinates,&coords);CHKERRQ(ierr);
1814         ierr = DMSetCoordinatesLocal(dm,coordinates);CHKERRQ(ierr);
1815       }
1816     } else {
1817       ierr = DMPlexCreateFromFile(comm, filename, PETSC_TRUE, &dm);CHKERRQ(ierr);
1818     }
1819   }
1820   ierr = DMViewFromOptions(dm, NULL, "-orig_dm_view");CHKERRQ(ierr);
1821   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1822 
1823   /* set up BCs, functions, tags */
1824   ierr = DMCreateLabel(dm, "Face Sets");CHKERRQ(ierr);
1825 
1826   mod->errorIndicator = ErrorIndicator_Simple;
1827 
1828   {
1829     DM dmDist;
1830 
1831     ierr = DMSetBasicAdjacency(dm, PETSC_TRUE, PETSC_FALSE);CHKERRQ(ierr);
1832     ierr = DMPlexDistribute(dm, overlap, NULL, &dmDist);CHKERRQ(ierr);
1833     if (dmDist) {
1834       ierr = DMDestroy(&dm);CHKERRQ(ierr);
1835       dm   = dmDist;
1836     }
1837   }
1838 
1839   ierr = DMSetFromOptions(dm);CHKERRQ(ierr);
1840 
1841   {
1842     DM gdm;
1843 
1844     ierr = DMPlexConstructGhostCells(dm, NULL, NULL, &gdm);CHKERRQ(ierr);
1845     ierr = DMDestroy(&dm);CHKERRQ(ierr);
1846     dm   = gdm;
1847     ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
1848   }
1849   if (splitFaces) {ierr = ConstructCellBoundary(dm, user);CHKERRQ(ierr);}
1850   ierr = SplitFaces(&dm, "split faces", user);CHKERRQ(ierr);
1851 
1852   ierr = PetscFVCreate(comm, &fvm);CHKERRQ(ierr);
1853   ierr = PetscFVSetFromOptions(fvm);CHKERRQ(ierr);
1854   ierr = PetscFVSetNumComponents(fvm, phys->dof);CHKERRQ(ierr);
1855   ierr = PetscFVSetSpatialDimension(fvm, dim);CHKERRQ(ierr);
1856   ierr = PetscObjectSetName((PetscObject) fvm,"");CHKERRQ(ierr);
1857   {
1858     PetscInt f, dof;
1859     for (f=0,dof=0; f < phys->nfields; f++) {
1860       PetscInt newDof = phys->field_desc[f].dof;
1861 
1862       if (newDof == 1) {
1863         ierr = PetscFVSetComponentName(fvm,dof,phys->field_desc[f].name);CHKERRQ(ierr);
1864       }
1865       else {
1866         PetscInt j;
1867 
1868         for (j = 0; j < newDof; j++) {
1869           char     compName[256]  = "Unknown";
1870 
1871           ierr = PetscSNPrintf(compName,sizeof(compName),"%s_%d",phys->field_desc[f].name,j);CHKERRQ(ierr);
1872           ierr = PetscFVSetComponentName(fvm,dof+j,compName);CHKERRQ(ierr);
1873         }
1874       }
1875       dof += newDof;
1876     }
1877   }
1878   /* FV is now structured with one field having all physics as components */
1879   ierr = DMAddField(dm, NULL, (PetscObject) fvm);CHKERRQ(ierr);
1880   ierr = DMCreateDS(dm);CHKERRQ(ierr);
1881   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1882   ierr = PetscDSSetRiemannSolver(prob, 0, user->model->physics->riemann);CHKERRQ(ierr);
1883   ierr = PetscDSSetContext(prob, 0, user->model->physics);CHKERRQ(ierr);
1884   ierr = (*mod->setupbc)(dm, prob,phys);CHKERRQ(ierr);
1885   ierr = PetscDSSetFromOptions(prob);CHKERRQ(ierr);
1886   {
1887     char      convType[256];
1888     PetscBool flg;
1889 
1890     ierr = PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");CHKERRQ(ierr);
1891     ierr = PetscOptionsFList("-dm_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg);CHKERRQ(ierr);
1892     ierr = PetscOptionsEnd();CHKERRQ(ierr);
1893     if (flg) {
1894       DM dmConv;
1895 
1896       ierr = DMConvert(dm,convType,&dmConv);CHKERRQ(ierr);
1897       if (dmConv) {
1898         ierr = DMViewFromOptions(dmConv, NULL, "-dm_conv_view");CHKERRQ(ierr);
1899         ierr = DMDestroy(&dm);CHKERRQ(ierr);
1900         dm   = dmConv;
1901         ierr = DMSetFromOptions(dm);CHKERRQ(ierr);
1902       }
1903     }
1904   }
1905 
1906   ierr = initializeTS(dm, user, &ts);CHKERRQ(ierr);
1907 
1908   ierr = DMCreateGlobalVector(dm, &X);CHKERRQ(ierr);
1909   ierr = PetscObjectSetName((PetscObject) X, "solution");CHKERRQ(ierr);
1910   ierr = SetInitialCondition(dm, X, user);CHKERRQ(ierr);
1911   if (useAMR) {
1912     PetscInt adaptIter;
1913 
1914     /* use no limiting when reconstructing gradients for adaptivity */
1915     ierr = PetscFVGetLimiter(fvm, &limiter);CHKERRQ(ierr);
1916     ierr = PetscObjectReference((PetscObject) limiter);CHKERRQ(ierr);
1917     ierr = PetscLimiterCreate(PetscObjectComm((PetscObject) fvm), &noneLimiter);CHKERRQ(ierr);
1918     ierr = PetscLimiterSetType(noneLimiter, PETSCLIMITERNONE);CHKERRQ(ierr);
1919 
1920     ierr = PetscFVSetLimiter(fvm, noneLimiter);CHKERRQ(ierr);
1921     for (adaptIter = 0; ; ++adaptIter) {
1922       PetscLogDouble bytes;
1923       TS             tsNew = NULL;
1924 
1925       ierr = PetscMemoryGetCurrentUsage(&bytes);CHKERRQ(ierr);
1926       ierr = PetscInfo2(ts, "refinement loop %D: memory used %g\n", adaptIter, bytes);CHKERRQ(ierr);
1927       ierr = DMViewFromOptions(dm, NULL, "-initial_dm_view");CHKERRQ(ierr);
1928       ierr = VecViewFromOptions(X, NULL, "-initial_vec_view");CHKERRQ(ierr);
1929 #if 0
1930       if (viewInitial) {
1931         PetscViewer viewer;
1932         char        buf[256];
1933         PetscBool   isHDF5, isVTK;
1934 
1935         ierr = PetscViewerCreate(comm,&viewer);CHKERRQ(ierr);
1936         ierr = PetscViewerSetType(viewer,PETSCVIEWERVTK);CHKERRQ(ierr);
1937         ierr = PetscViewerSetOptionsPrefix(viewer,"initial_");CHKERRQ(ierr);
1938         ierr = PetscViewerSetFromOptions(viewer);CHKERRQ(ierr);
1939         ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&isHDF5);CHKERRQ(ierr);
1940         ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isVTK);CHKERRQ(ierr);
1941         if (isHDF5) {
1942           ierr = PetscSNPrintf(buf, 256, "ex11-initial-%d.h5", adaptIter);CHKERRQ(ierr);
1943         } else if (isVTK) {
1944           ierr = PetscSNPrintf(buf, 256, "ex11-initial-%d.vtu", adaptIter);CHKERRQ(ierr);
1945           ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_VTK_VTU);CHKERRQ(ierr);
1946         }
1947         ierr = PetscViewerFileSetMode(viewer,FILE_MODE_WRITE);CHKERRQ(ierr);
1948         ierr = PetscViewerFileSetName(viewer,buf);CHKERRQ(ierr);
1949         if (isHDF5) {
1950           ierr = DMView(dm,viewer);CHKERRQ(ierr);
1951           ierr = PetscViewerFileSetMode(viewer,FILE_MODE_UPDATE);CHKERRQ(ierr);
1952         }
1953         ierr = VecView(X,viewer);CHKERRQ(ierr);
1954         ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1955       }
1956 #endif
1957 
1958       ierr = adaptToleranceFVM(fvm, ts, X, refineTag, coarsenTag, user, &tsNew, NULL);CHKERRQ(ierr);
1959       if (!tsNew) {
1960         break;
1961       } else {
1962         ierr = DMDestroy(&dm);CHKERRQ(ierr);
1963         ierr = VecDestroy(&X);CHKERRQ(ierr);
1964         ierr = TSDestroy(&ts);CHKERRQ(ierr);
1965         ts   = tsNew;
1966         ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1967         ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
1968         ierr = DMCreateGlobalVector(dm,&X);CHKERRQ(ierr);
1969         ierr = PetscObjectSetName((PetscObject) X, "solution");CHKERRQ(ierr);
1970         ierr = SetInitialCondition(dm, X, user);CHKERRQ(ierr);
1971       }
1972     }
1973     /* restore original limiter */
1974     ierr = PetscFVSetLimiter(fvm, limiter);CHKERRQ(ierr);
1975   }
1976 
1977   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
1978   if (vtkCellGeom) {
1979     DM  dmCell;
1980     Vec cellgeom, partition;
1981 
1982     ierr = DMPlexGetGeometryFVM(plex, NULL, &cellgeom, NULL);CHKERRQ(ierr);
1983     ierr = OutputVTK(dm, "ex11-cellgeom.vtk", &viewer);CHKERRQ(ierr);
1984     ierr = VecView(cellgeom, viewer);CHKERRQ(ierr);
1985     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1986     ierr = CreatePartitionVec(dm, &dmCell, &partition);CHKERRQ(ierr);
1987     ierr = OutputVTK(dmCell, "ex11-partition.vtk", &viewer);CHKERRQ(ierr);
1988     ierr = VecView(partition, viewer);CHKERRQ(ierr);
1989     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1990     ierr = VecDestroy(&partition);CHKERRQ(ierr);
1991     ierr = DMDestroy(&dmCell);CHKERRQ(ierr);
1992   }
1993   /* collect max maxspeed from all processes -- todo */
1994   ierr = DMPlexGetGeometryFVM(plex, NULL, NULL, &minRadius);CHKERRQ(ierr);
1995   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1996   ierr = MPI_Allreduce(&phys->maxspeed,&mod->maxspeed,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));CHKERRMPI(ierr);
1997   if (mod->maxspeed <= 0) SETERRQ1(comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set maxspeed",physname);
1998   dt   = cfl * minRadius / mod->maxspeed;
1999   ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
2000   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
2001   if (!useAMR) {
2002     ierr = TSSolve(ts,X);CHKERRQ(ierr);
2003     ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
2004     ierr = TSGetStepNumber(ts,&nsteps);CHKERRQ(ierr);
2005   } else {
2006     PetscReal finalTime;
2007     PetscInt  adaptIter;
2008     TS        tsNew = NULL;
2009     Vec       solNew = NULL;
2010 
2011     ierr   = TSGetMaxTime(ts,&finalTime);CHKERRQ(ierr);
2012     ierr   = TSSetMaxSteps(ts,adaptInterval);CHKERRQ(ierr);
2013     ierr   = TSSolve(ts,X);CHKERRQ(ierr);
2014     ierr   = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
2015     ierr   = TSGetStepNumber(ts,&nsteps);CHKERRQ(ierr);
2016     for (adaptIter = 0;ftime < finalTime;adaptIter++) {
2017       PetscLogDouble bytes;
2018 
2019       ierr = PetscMemoryGetCurrentUsage(&bytes);CHKERRQ(ierr);
2020       ierr = PetscInfo2(ts, "AMR time step loop %D: memory used %g\n", adaptIter, bytes);CHKERRQ(ierr);
2021       ierr = PetscFVSetLimiter(fvm,noneLimiter);CHKERRQ(ierr);
2022       ierr = adaptToleranceFVM(fvm,ts,X,refineTag,coarsenTag,user,&tsNew,&solNew);CHKERRQ(ierr);
2023       ierr = PetscFVSetLimiter(fvm,limiter);CHKERRQ(ierr);
2024       if (tsNew) {
2025         ierr = PetscInfo(ts, "AMR used\n");CHKERRQ(ierr);
2026         ierr = DMDestroy(&dm);CHKERRQ(ierr);
2027         ierr = VecDestroy(&X);CHKERRQ(ierr);
2028         ierr = TSDestroy(&ts);CHKERRQ(ierr);
2029         ts   = tsNew;
2030         X    = solNew;
2031         ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
2032         ierr = VecGetDM(X,&dm);CHKERRQ(ierr);
2033         ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
2034         ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
2035         ierr = DMPlexGetGeometryFVM(dm, NULL, NULL, &minRadius);CHKERRQ(ierr);
2036         ierr = DMDestroy(&plex);CHKERRQ(ierr);
2037         ierr = MPI_Allreduce(&phys->maxspeed,&mod->maxspeed,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));CHKERRMPI(ierr);
2038         if (mod->maxspeed <= 0) SETERRQ1(comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set maxspeed",physname);
2039         dt   = cfl * minRadius / mod->maxspeed;
2040         ierr = TSSetStepNumber(ts,nsteps);CHKERRQ(ierr);
2041         ierr = TSSetTime(ts,ftime);CHKERRQ(ierr);
2042         ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
2043       } else {
2044         ierr = PetscInfo(ts, "AMR not used\n");CHKERRQ(ierr);
2045       }
2046       user->monitorStepOffset = nsteps;
2047       ierr = TSSetMaxSteps(ts,nsteps+adaptInterval);CHKERRQ(ierr);
2048       ierr = TSSolve(ts,X);CHKERRQ(ierr);
2049       ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
2050       ierr = TSGetStepNumber(ts,&nsteps);CHKERRQ(ierr);
2051     }
2052   }
2053   ierr = TSGetConvergedReason(ts,&reason);CHKERRQ(ierr);
2054   ierr = PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,nsteps);CHKERRQ(ierr);
2055   ierr = TSDestroy(&ts);CHKERRQ(ierr);
2056 
2057   ierr = VecTaggerDestroy(&refineTag);CHKERRQ(ierr);
2058   ierr = VecTaggerDestroy(&coarsenTag);CHKERRQ(ierr);
2059   ierr = PetscFunctionListDestroy(&PhysicsList);CHKERRQ(ierr);
2060   ierr = PetscFunctionListDestroy(&PhysicsRiemannList_SW);CHKERRQ(ierr);
2061   ierr = FunctionalLinkDestroy(&user->model->functionalRegistry);CHKERRQ(ierr);
2062   ierr = PetscFree(user->model->functionalMonitored);CHKERRQ(ierr);
2063   ierr = PetscFree(user->model->functionalCall);CHKERRQ(ierr);
2064   ierr = PetscFree(user->model->physics->data);CHKERRQ(ierr);
2065   ierr = PetscFree(user->model->physics);CHKERRQ(ierr);
2066   ierr = PetscFree(user->model);CHKERRQ(ierr);
2067   ierr = PetscFree(user);CHKERRQ(ierr);
2068   ierr = VecDestroy(&X);CHKERRQ(ierr);
2069   ierr = PetscLimiterDestroy(&limiter);CHKERRQ(ierr);
2070   ierr = PetscLimiterDestroy(&noneLimiter);CHKERRQ(ierr);
2071   ierr = PetscFVDestroy(&fvm);CHKERRQ(ierr);
2072   ierr = DMDestroy(&dm);CHKERRQ(ierr);
2073   ierr = PetscFinalize();
2074   return ierr;
2075 }
2076 
2077 /* Godunov fluxs */
2078 PetscScalar cvmgp_(PetscScalar *a, PetscScalar *b, PetscScalar *test)
2079 {
2080     /* System generated locals */
2081     PetscScalar ret_val;
2082 
2083     if (PetscRealPart(*test) > 0.) {
2084         goto L10;
2085     }
2086     ret_val = *b;
2087     return ret_val;
2088 L10:
2089     ret_val = *a;
2090     return ret_val;
2091 } /* cvmgp_ */
2092 
2093 PetscScalar cvmgm_(PetscScalar *a, PetscScalar *b, PetscScalar *test)
2094 {
2095     /* System generated locals */
2096     PetscScalar ret_val;
2097 
2098     if (PetscRealPart(*test) < 0.) {
2099         goto L10;
2100     }
2101     ret_val = *b;
2102     return ret_val;
2103 L10:
2104     ret_val = *a;
2105     return ret_val;
2106 } /* cvmgm_ */
2107 
2108 int riem1mdt( PetscScalar *gaml, PetscScalar *gamr, PetscScalar *rl, PetscScalar *pl,
2109               PetscScalar *uxl, PetscScalar *rr, PetscScalar *pr,
2110               PetscScalar *uxr, PetscScalar *rstarl, PetscScalar *rstarr, PetscScalar *
2111               pstar, PetscScalar *ustar)
2112 {
2113     /* Initialized data */
2114 
2115     static PetscScalar smallp = 1e-8;
2116 
2117     /* System generated locals */
2118     int i__1;
2119     PetscScalar d__1, d__2;
2120 
2121     /* Local variables */
2122     static int i0;
2123     static PetscScalar cl, cr, wl, zl, wr, zr, pst, durl, skpr1, skpr2;
2124     static int iwave;
2125     static PetscScalar gascl4, gascr4, cstarl, dpstar, cstarr;
2126     /* static PetscScalar csqrl, csqrr, gascl1, gascl2, gascl3, gascr1, gascr2, gascr3; */
2127     static int iterno;
2128     static PetscScalar ustarl, ustarr, rarepr1, rarepr2;
2129 
2130 
2131 
2132     /* gascl1 = *gaml - 1.; */
2133     /* gascl2 = (*gaml + 1.) * .5; */
2134     /* gascl3 = gascl2 / *gaml; */
2135     gascl4 = 1. / (*gaml - 1.);
2136 
2137     /* gascr1 = *gamr - 1.; */
2138     /* gascr2 = (*gamr + 1.) * .5; */
2139     /* gascr3 = gascr2 / *gamr; */
2140     gascr4 = 1. / (*gamr - 1.);
2141     iterno = 10;
2142 /*        find pstar: */
2143     cl = PetscSqrtScalar(*gaml * *pl / *rl);
2144     cr = PetscSqrtScalar(*gamr * *pr / *rr);
2145     wl = *rl * cl;
2146     wr = *rr * cr;
2147     /* csqrl = wl * wl; */
2148     /* csqrr = wr * wr; */
2149     *pstar = (wl * *pr + wr * *pl) / (wl + wr);
2150     *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp));
2151     pst = *pl / *pr;
2152     skpr1 = cr * (pst - 1.) * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
2153     d__1 = (*gamr - 1.) / (*gamr * 2.);
2154     rarepr2 = gascr4 * 2. * cr * (1. - PetscPowScalar(pst, d__1));
2155     pst = *pr / *pl;
2156     skpr2 = cl * (pst - 1.) * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
2157     d__1 = (*gaml - 1.) / (*gaml * 2.);
2158     rarepr1 = gascl4 * 2. * cl * (1. - PetscPowScalar(pst, d__1));
2159     durl = *uxr - *uxl;
2160     if (PetscRealPart(*pr) < PetscRealPart(*pl)) {
2161         if (PetscRealPart(durl) >= PetscRealPart(rarepr1)) {
2162             iwave = 100;
2163         } else if (PetscRealPart(durl) <= PetscRealPart(-skpr1)) {
2164             iwave = 300;
2165         } else {
2166             iwave = 400;
2167         }
2168     } else {
2169         if (PetscRealPart(durl) >= PetscRealPart(rarepr2)) {
2170             iwave = 100;
2171         } else if (PetscRealPart(durl) <= PetscRealPart(-skpr2)) {
2172             iwave = 300;
2173         } else {
2174             iwave = 200;
2175         }
2176     }
2177     if (iwave == 100) {
2178 /*     1-wave: rarefaction wave, 3-wave: rarefaction wave */
2179 /*     case (100) */
2180         i__1 = iterno;
2181         for (i0 = 1; i0 <= i__1; ++i0) {
2182             d__1 = *pstar / *pl;
2183             d__2 = 1. / *gaml;
2184             *rstarl = *rl * PetscPowScalar(d__1, d__2);
2185             cstarl = PetscSqrtScalar(*gaml * *pstar / *rstarl);
2186             ustarl = *uxl - gascl4 * 2. * (cstarl - cl);
2187             zl = *rstarl * cstarl;
2188             d__1 = *pstar / *pr;
2189             d__2 = 1. / *gamr;
2190             *rstarr = *rr * PetscPowScalar(d__1, d__2);
2191             cstarr = PetscSqrtScalar(*gamr * *pstar / *rstarr);
2192             ustarr = *uxr + gascr4 * 2. * (cstarr - cr);
2193             zr = *rstarr * cstarr;
2194             dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
2195             *pstar -= dpstar;
2196             *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp));
2197             if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
2198 #if 0
2199         break;
2200 #endif
2201             }
2202         }
2203 /*     1-wave: shock wave, 3-wave: rarefaction wave */
2204     } else if (iwave == 200) {
2205 /*     case (200) */
2206         i__1 = iterno;
2207         for (i0 = 1; i0 <= i__1; ++i0) {
2208             pst = *pstar / *pl;
2209             ustarl = *uxl - (pst - 1.) * cl * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
2210             zl = *pl / cl * PetscSqrtScalar(*gaml * 2. * (*gaml - 1. + (*gaml + 1.) * pst)) * (*gaml - 1. + (*gaml + 1.) * pst) / (*gaml * 3. - 1. + (*gaml + 1.) * pst);
2211             d__1 = *pstar / *pr;
2212             d__2 = 1. / *gamr;
2213             *rstarr = *rr * PetscPowScalar(d__1, d__2);
2214             cstarr = PetscSqrtScalar(*gamr * *pstar / *rstarr);
2215             zr = *rstarr * cstarr;
2216             ustarr = *uxr + gascr4 * 2. * (cstarr - cr);
2217             dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
2218             *pstar -= dpstar;
2219             *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp));
2220             if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
2221 #if 0
2222         break;
2223 #endif
2224             }
2225         }
2226 /*     1-wave: shock wave, 3-wave: shock */
2227     } else if (iwave == 300) {
2228 /*     case (300) */
2229         i__1 = iterno;
2230         for (i0 = 1; i0 <= i__1; ++i0) {
2231             pst = *pstar / *pl;
2232             ustarl = *uxl - (pst - 1.) * cl * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
2233             zl = *pl / cl * PetscSqrtScalar(*gaml * 2. * (*gaml - 1. + (*gaml + 1.) * pst)) * (*gaml - 1. + (*gaml + 1.) * pst) / (*gaml * 3. - 1. + (*gaml + 1.) * pst);
2234             pst = *pstar / *pr;
2235             ustarr = *uxr + (pst - 1.) * cr * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
2236             zr = *pr / cr * PetscSqrtScalar(*gamr * 2. * (*gamr - 1. + (*gamr + 1.) * pst)) * (*gamr - 1. + (*gamr + 1.) * pst) / (*gamr * 3. - 1. + (*gamr + 1.) * pst);
2237             dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
2238             *pstar -= dpstar;
2239             *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp));
2240             if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
2241 #if 0
2242         break;
2243 #endif
2244             }
2245         }
2246 /*     1-wave: rarefaction wave, 3-wave: shock */
2247     } else if (iwave == 400) {
2248 /*     case (400) */
2249         i__1 = iterno;
2250         for (i0 = 1; i0 <= i__1; ++i0) {
2251             d__1 = *pstar / *pl;
2252             d__2 = 1. / *gaml;
2253             *rstarl = *rl * PetscPowScalar(d__1, d__2);
2254             cstarl = PetscSqrtScalar(*gaml * *pstar / *rstarl);
2255             ustarl = *uxl - gascl4 * 2. * (cstarl - cl);
2256             zl = *rstarl * cstarl;
2257             pst = *pstar / *pr;
2258             ustarr = *uxr + (pst - 1.) * cr * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
2259             zr = *pr / cr * PetscSqrtScalar(*gamr * 2. * (*gamr - 1. + (*gamr + 1.) * pst)) * (*gamr - 1. + (*gamr + 1.) * pst) / (*gamr * 3. - 1. + (*gamr + 1.) * pst);
2260             dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
2261             *pstar -= dpstar;
2262             *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp));
2263             if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
2264 #if 0
2265               break;
2266 #endif
2267             }
2268         }
2269     }
2270 
2271     *ustar = (zl * ustarr + zr * ustarl) / (zl + zr);
2272     if (PetscRealPart(*pstar) > PetscRealPart(*pl)) {
2273         pst = *pstar / *pl;
2274         *rstarl = ((*gaml + 1.) * pst + *gaml - 1.) / ((*gaml - 1.) * pst + *
2275                 gaml + 1.) * *rl;
2276     }
2277     if (PetscRealPart(*pstar) > PetscRealPart(*pr)) {
2278         pst = *pstar / *pr;
2279         *rstarr = ((*gamr + 1.) * pst + *gamr - 1.) / ((*gamr - 1.) * pst + *
2280                 gamr + 1.) * *rr;
2281     }
2282     return iwave;
2283 }
2284 
2285 PetscScalar sign(PetscScalar x)
2286 {
2287     if (PetscRealPart(x) > 0) return 1.0;
2288     if (PetscRealPart(x) < 0) return -1.0;
2289     return 0.0;
2290 }
2291 /*        Riemann Solver */
2292 /* -------------------------------------------------------------------- */
2293 int riemannsolver(PetscScalar *xcen, PetscScalar *xp,
2294                    PetscScalar *dtt, PetscScalar *rl, PetscScalar *uxl, PetscScalar *pl,
2295                    PetscScalar *utl, PetscScalar *ubl, PetscScalar *gaml, PetscScalar *rho1l,
2296                    PetscScalar *rr, PetscScalar *uxr, PetscScalar *pr, PetscScalar *utr,
2297                    PetscScalar *ubr, PetscScalar *gamr, PetscScalar *rho1r, PetscScalar *rx,
2298                    PetscScalar *uxm, PetscScalar *px, PetscScalar *utx, PetscScalar *ubx,
2299                    PetscScalar *gam, PetscScalar *rho1)
2300 {
2301     /* System generated locals */
2302     PetscScalar d__1, d__2;
2303 
2304     /* Local variables */
2305     static PetscScalar s, c0, p0, r0, u0, w0, x0, x2, ri, cx, sgn0, wsp0, gasc1, gasc2, gasc3, gasc4;
2306     static PetscScalar cstar, pstar, rstar, ustar, xstar, wspst, ushock, streng, rstarl, rstarr, rstars;
2307     int iwave;
2308 
2309     if (*rl == *rr && *pr == *pl && *uxl == *uxr && *gaml == *gamr) {
2310         *rx = *rl;
2311         *px = *pl;
2312         *uxm = *uxl;
2313         *gam = *gaml;
2314         x2 = *xcen + *uxm * *dtt;
2315 
2316         if (PetscRealPart(*xp) >= PetscRealPart(x2)) {
2317             *utx = *utr;
2318             *ubx = *ubr;
2319             *rho1 = *rho1r;
2320         } else {
2321             *utx = *utl;
2322             *ubx = *ubl;
2323             *rho1 = *rho1l;
2324         }
2325         return 0;
2326     }
2327     iwave = riem1mdt(gaml, gamr, rl, pl, uxl, rr, pr, uxr, &rstarl, &rstarr, &pstar, &ustar);
2328 
2329     x2 = *xcen + ustar * *dtt;
2330     d__1 = *xp - x2;
2331     sgn0 = sign(d__1);
2332 /*            x is in 3-wave if sgn0 = 1 */
2333 /*            x is in 1-wave if sgn0 = -1 */
2334     r0 = cvmgm_(rl, rr, &sgn0);
2335     p0 = cvmgm_(pl, pr, &sgn0);
2336     u0 = cvmgm_(uxl, uxr, &sgn0);
2337     *gam = cvmgm_(gaml, gamr, &sgn0);
2338     gasc1 = *gam - 1.;
2339     gasc2 = (*gam + 1.) * .5;
2340     gasc3 = gasc2 / *gam;
2341     gasc4 = 1. / (*gam - 1.);
2342     c0 = PetscSqrtScalar(*gam * p0 / r0);
2343     streng = pstar - p0;
2344     w0 = *gam * r0 * p0 * (gasc3 * streng / p0 + 1.);
2345     rstars = r0 / (1. - r0 * streng / w0);
2346     d__1 = p0 / pstar;
2347     d__2 = -1. / *gam;
2348     rstarr = r0 * PetscPowScalar(d__1, d__2);
2349     rstar = cvmgm_(&rstarr, &rstars, &streng);
2350     w0 = PetscSqrtScalar(w0);
2351     cstar = PetscSqrtScalar(*gam * pstar / rstar);
2352     wsp0 = u0 + sgn0 * c0;
2353     wspst = ustar + sgn0 * cstar;
2354     ushock = ustar + sgn0 * w0 / rstar;
2355     wspst = cvmgp_(&ushock, &wspst, &streng);
2356     wsp0 = cvmgp_(&ushock, &wsp0, &streng);
2357     x0 = *xcen + wsp0 * *dtt;
2358     xstar = *xcen + wspst * *dtt;
2359 /*           using gas formula to evaluate rarefaction wave */
2360 /*            ri : reiman invariant */
2361     ri = u0 - sgn0 * 2. * gasc4 * c0;
2362     cx = sgn0 * .5 * gasc1 / gasc2 * ((*xp - *xcen) / *dtt - ri);
2363     *uxm = ri + sgn0 * 2. * gasc4 * cx;
2364     s = p0 / PetscPowScalar(r0, *gam);
2365     d__1 = cx * cx / (*gam * s);
2366     *rx = PetscPowScalar(d__1, gasc4);
2367     *px = cx * cx * *rx / *gam;
2368     d__1 = sgn0 * (x0 - *xp);
2369     *rx = cvmgp_(rx, &r0, &d__1);
2370     d__1 = sgn0 * (x0 - *xp);
2371     *px = cvmgp_(px, &p0, &d__1);
2372     d__1 = sgn0 * (x0 - *xp);
2373     *uxm = cvmgp_(uxm, &u0, &d__1);
2374     d__1 = sgn0 * (xstar - *xp);
2375     *rx = cvmgm_(rx, &rstar, &d__1);
2376     d__1 = sgn0 * (xstar - *xp);
2377     *px = cvmgm_(px, &pstar, &d__1);
2378     d__1 = sgn0 * (xstar - *xp);
2379     *uxm = cvmgm_(uxm, &ustar, &d__1);
2380     if (PetscRealPart(*xp) >= PetscRealPart(x2)) {
2381         *utx = *utr;
2382         *ubx = *ubr;
2383         *rho1 = *rho1r;
2384     } else {
2385         *utx = *utl;
2386         *ubx = *ubl;
2387         *rho1 = *rho1l;
2388     }
2389     return iwave;
2390 }
2391 int godunovflux( const PetscScalar *ul, const PetscScalar *ur,
2392                  PetscScalar *flux, const PetscReal *nn, const int *ndim,
2393                  const PetscReal *gamma)
2394 {
2395     /* System generated locals */
2396   int i__1,iwave;
2397     PetscScalar d__1, d__2, d__3;
2398 
2399     /* Local variables */
2400     static int k;
2401     static PetscScalar bn[3], fn, ft, tg[3], pl, rl, pm, pr, rr, xp, ubl, ubm,
2402             ubr, dtt, unm, tmp, utl, utm, uxl, utr, uxr, gaml, gamm, gamr,
2403             xcen, rhom, rho1l, rho1m, rho1r;
2404 
2405     /* Function Body */
2406     xcen = 0.;
2407     xp = 0.;
2408     i__1 = *ndim;
2409     for (k = 1; k <= i__1; ++k) {
2410         tg[k - 1] = 0.;
2411         bn[k - 1] = 0.;
2412     }
2413     dtt = 1.;
2414     if (*ndim == 3) {
2415         if (nn[0] == 0. && nn[1] == 0.) {
2416             tg[0] = 1.;
2417         } else {
2418             tg[0] = -nn[1];
2419             tg[1] = nn[0];
2420         }
2421 /*           tmp=dsqrt(tg(1)**2+tg(2)**2) */
2422 /*           tg=tg/tmp */
2423         bn[0] = -nn[2] * tg[1];
2424         bn[1] = nn[2] * tg[0];
2425         bn[2] = nn[0] * tg[1] - nn[1] * tg[0];
2426 /* Computing 2nd power */
2427         d__1 = bn[0];
2428 /* Computing 2nd power */
2429         d__2 = bn[1];
2430 /* Computing 2nd power */
2431         d__3 = bn[2];
2432         tmp = PetscSqrtScalar(d__1 * d__1 + d__2 * d__2 + d__3 * d__3);
2433         i__1 = *ndim;
2434         for (k = 1; k <= i__1; ++k) {
2435             bn[k - 1] /= tmp;
2436         }
2437     } else if (*ndim == 2) {
2438         tg[0] = -nn[1];
2439         tg[1] = nn[0];
2440 /*           tmp=dsqrt(tg(1)**2+tg(2)**2) */
2441 /*           tg=tg/tmp */
2442         bn[0] = 0.;
2443         bn[1] = 0.;
2444         bn[2] = 1.;
2445     }
2446     rl = ul[0];
2447     rr = ur[0];
2448     uxl = 0.;
2449     uxr = 0.;
2450     utl = 0.;
2451     utr = 0.;
2452     ubl = 0.;
2453     ubr = 0.;
2454     i__1 = *ndim;
2455     for (k = 1; k <= i__1; ++k) {
2456         uxl += ul[k] * nn[k-1];
2457         uxr += ur[k] * nn[k-1];
2458         utl += ul[k] * tg[k - 1];
2459         utr += ur[k] * tg[k - 1];
2460         ubl += ul[k] * bn[k - 1];
2461         ubr += ur[k] * bn[k - 1];
2462     }
2463     uxl /= rl;
2464     uxr /= rr;
2465     utl /= rl;
2466     utr /= rr;
2467     ubl /= rl;
2468     ubr /= rr;
2469 
2470     gaml = *gamma;
2471     gamr = *gamma;
2472 /* Computing 2nd power */
2473     d__1 = uxl;
2474 /* Computing 2nd power */
2475     d__2 = utl;
2476 /* Computing 2nd power */
2477     d__3 = ubl;
2478     pl = (*gamma - 1.) * (ul[*ndim + 1] - rl * .5 * (d__1 * d__1 + d__2 * d__2 + d__3 * d__3));
2479 /* Computing 2nd power */
2480     d__1 = uxr;
2481 /* Computing 2nd power */
2482     d__2 = utr;
2483 /* Computing 2nd power */
2484     d__3 = ubr;
2485     pr = (*gamma - 1.) * (ur[*ndim + 1] - rr * .5 * (d__1 * d__1 + d__2 * d__2 + d__3 * d__3));
2486     rho1l = rl;
2487     rho1r = rr;
2488 
2489     iwave = riemannsolver(&xcen, &xp, &dtt, &rl, &uxl, &pl, &utl, &ubl, &gaml, &
2490                           rho1l, &rr, &uxr, &pr, &utr, &ubr, &gamr, &rho1r, &rhom, &unm, &
2491                           pm, &utm, &ubm, &gamm, &rho1m);
2492 
2493     flux[0] = rhom * unm;
2494     fn = rhom * unm * unm + pm;
2495     ft = rhom * unm * utm;
2496 /*           flux(2)=fn*nn(1)+ft*nn(2) */
2497 /*           flux(3)=fn*tg(1)+ft*tg(2) */
2498     flux[1] = fn * nn[0] + ft * tg[0];
2499     flux[2] = fn * nn[1] + ft * tg[1];
2500 /*           flux(2)=rhom*unm*(unm)+pm */
2501 /*           flux(3)=rhom*(unm)*utm */
2502     if (*ndim == 3) {
2503         flux[3] = rhom * unm * ubm;
2504     }
2505     flux[*ndim + 1] = (rhom * .5 * (unm * unm + utm * utm + ubm * ubm) + gamm / (gamm - 1.) * pm) * unm;
2506     return iwave;
2507 } /* godunovflux_ */
2508 
2509 /* Subroutine to set up the initial conditions for the */
2510 /* Shock Interface interaction or linear wave (Ravi Samtaney,Mark Adams). */
2511 /* ----------------------------------------------------------------------- */
2512 int projecteqstate(PetscReal wc[], const PetscReal ueq[], PetscReal lv[][3])
2513 {
2514   int j,k;
2515 /*      Wc=matmul(lv,Ueq) 3 vars */
2516   for (k = 0; k < 3; ++k) {
2517     wc[k] = 0.;
2518     for (j = 0; j < 3; ++j) {
2519       wc[k] += lv[k][j]*ueq[j];
2520     }
2521   }
2522   return 0;
2523 }
2524 /* ----------------------------------------------------------------------- */
2525 int projecttoprim(PetscReal v[], const PetscReal wc[], PetscReal rv[][3])
2526 {
2527   int k,j;
2528   /*      V=matmul(rv,WC) 3 vars */
2529   for (k = 0; k < 3; ++k) {
2530     v[k] = 0.;
2531     for (j = 0; j < 3; ++j) {
2532       v[k] += rv[k][j]*wc[j];
2533     }
2534   }
2535   return 0;
2536 }
2537 /* ---------------------------------------------------------------------- */
2538 int eigenvectors(PetscReal rv[][3], PetscReal lv[][3], const PetscReal ueq[], PetscReal gamma)
2539 {
2540   int j,k;
2541   PetscReal rho,csnd,p0;
2542   /* PetscScalar u; */
2543 
2544   for (k = 0; k < 3; ++k) for (j = 0; j < 3; ++j) { lv[k][j] = 0.; rv[k][j] = 0.; }
2545   rho = ueq[0];
2546   /* u = ueq[1]; */
2547   p0 = ueq[2];
2548   csnd = PetscSqrtReal(gamma * p0 / rho);
2549   lv[0][1] = rho * .5;
2550   lv[0][2] = -.5 / csnd;
2551   lv[1][0] = csnd;
2552   lv[1][2] = -1. / csnd;
2553   lv[2][1] = rho * .5;
2554   lv[2][2] = .5 / csnd;
2555   rv[0][0] = -1. / csnd;
2556   rv[1][0] = 1. / rho;
2557   rv[2][0] = -csnd;
2558   rv[0][1] = 1. / csnd;
2559   rv[0][2] = 1. / csnd;
2560   rv[1][2] = 1. / rho;
2561   rv[2][2] = csnd;
2562   return 0;
2563 }
2564 
2565 int initLinearWave(EulerNode *ux, const PetscReal gamma, const PetscReal coord[], const PetscReal Lx)
2566 {
2567   PetscReal p0,u0,wcp[3],wc[3];
2568   PetscReal lv[3][3];
2569   PetscReal vp[3];
2570   PetscReal rv[3][3];
2571   PetscReal eps, ueq[3], rho0, twopi;
2572 
2573   /* Function Body */
2574   twopi = 2.*PETSC_PI;
2575   eps = 1e-4; /* perturbation */
2576   rho0 = 1e3;   /* density of water */
2577   p0 = 101325.; /* init pressure of 1 atm (?) */
2578   u0 = 0.;
2579   ueq[0] = rho0;
2580   ueq[1] = u0;
2581   ueq[2] = p0;
2582   /* Project initial state to characteristic variables */
2583   eigenvectors(rv, lv, ueq, gamma);
2584   projecteqstate(wc, ueq, lv);
2585   wcp[0] = wc[0];
2586   wcp[1] = wc[1];
2587   wcp[2] = wc[2] + eps * PetscCosReal(coord[0] * 2. * twopi / Lx);
2588   projecttoprim(vp, wcp, rv);
2589   ux->r = vp[0]; /* density */
2590   ux->ru[0] = vp[0] * vp[1]; /* x momentum */
2591   ux->ru[1] = 0.;
2592 #if defined DIM > 2
2593   if (dim>2) ux->ru[2] = 0.;
2594 #endif
2595   /* E = rho * e + rho * v^2/2 = p/(gam-1) + rho*v^2/2 */
2596   ux->E = vp[2]/(gamma - 1.) + 0.5*vp[0]*vp[1]*vp[1];
2597   return 0;
2598 }
2599 
2600 /*TEST
2601 
2602   # 2D Advection 0-10
2603   test:
2604     suffix: 0
2605     requires: exodusii
2606     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo
2607 
2608   test:
2609     suffix: 1
2610     requires: exodusii
2611     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo
2612 
2613   test:
2614     suffix: 2
2615     requires: exodusii
2616     nsize: 2
2617     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo
2618 
2619   test:
2620     suffix: 3
2621     requires: exodusii
2622     nsize: 2
2623     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo
2624 
2625   test:
2626     suffix: 4
2627     requires: exodusii
2628     nsize: 8
2629     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo
2630 
2631   test:
2632     suffix: 5
2633     requires: exodusii
2634     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo -ts_type rosw -ts_adapt_reject_safety 1
2635 
2636   test:
2637     suffix: 6
2638     requires: exodusii
2639     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/squaremotor-30.exo -ufv_split_faces
2640 
2641   test:
2642     suffix: 7
2643     requires: exodusii
2644     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1
2645 
2646   test:
2647     suffix: 8
2648     requires: exodusii
2649     nsize: 2
2650     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1
2651 
2652   test:
2653     suffix: 9
2654     requires: exodusii
2655     nsize: 8
2656     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1
2657 
2658   test:
2659     suffix: 10
2660     requires: exodusii
2661     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo
2662 
2663   # 2D Shallow water
2664   test:
2665     suffix: sw_0
2666     requires: exodusii
2667     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -bc_wall 100,101 -physics sw -ufv_cfl 5 -petscfv_type leastsquares -petsclimiter_type sin -ts_max_time 1 -ts_ssp_type rks2 -ts_ssp_nstages 10 -monitor height,energy
2668 
2669   test:
2670     suffix: sw_hll
2671     args: -ufv_vtk_interval 0 -bc_wall 1,2,3,4 -physics sw -ufv_cfl 3 -petscfv_type leastsquares -petsclimiter_type sin -ts_max_steps 5 -ts_ssp_type rks2 -ts_ssp_nstages 10 -monitor height,energy -grid_bounds 0,5,0,5 -grid_size 25,25 -sw_riemann hll
2672 
2673   # 2D Advection: p4est
2674   test:
2675     suffix: p4est_advec_2d
2676     requires: p4est
2677     args: -ufv_vtk_interval 0 -f -dm_type p4est -dm_forest_minimum_refinement 1 -dm_forest_initial_refinement 2 -dm_p4est_refine_pattern hash -dm_forest_maximum_refinement 5
2678 
2679   # Advection in a box
2680   test:
2681     suffix: adv_2d_quad_0
2682     args: -ufv_vtk_interval 0 -dm_refine 3 -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3
2683 
2684   test:
2685     suffix: adv_2d_quad_1
2686     args: -ufv_vtk_interval 0 -dm_refine 3 -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1
2687     timeoutfactor: 3
2688 
2689   test:
2690     suffix: adv_2d_quad_p4est_0
2691     requires: p4est
2692     args: -ufv_vtk_interval 0 -dm_refine 5 -dm_type p4est -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3
2693 
2694   test:
2695     suffix: adv_2d_quad_p4est_1
2696     requires: p4est
2697     args: -ufv_vtk_interval 0 -dm_refine 5 -dm_type p4est -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1
2698     timeoutfactor: 3
2699 
2700   test:
2701     suffix: adv_2d_quad_p4est_adapt_0
2702     requires: p4est !__float128 #broken for quad precision
2703     args: -ufv_vtk_interval 0 -dm_refine 3 -dm_type p4est -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1 -ufv_use_amr -refine_vec_tagger_box 0.005,inf -coarsen_vec_tagger_box 0,1.e-5 -petscfv_type leastsquares -ts_max_time 0.01
2704     timeoutfactor: 3
2705 
2706   test:
2707     suffix: adv_2d_tri_0
2708     requires: triangle
2709     TODO: how did this ever get in main when there is no support for this
2710     args: -ufv_vtk_interval 0 -simplex -dm_refine 3 -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3
2711 
2712   test:
2713     suffix: adv_2d_tri_1
2714     requires: triangle
2715     TODO: how did this ever get in main when there is no support for this
2716     args: -ufv_vtk_interval 0 -simplex -dm_refine 5 -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1
2717 
2718   test:
2719     suffix: adv_0
2720     requires: exodusii
2721     args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -bc_inflow 100,101,200 -bc_outflow 201
2722 
2723   test:
2724     suffix: shock_0
2725     requires: p4est !single !complex
2726     args: -ufv_vtk_interval 0 -monitor density,energy -f -grid_size 2,1 -grid_bounds -1,1.,0.,1 -bc_wall 1,2,3,4 -dm_type p4est -dm_forest_partition_overlap 1 -dm_forest_maximum_refinement 6 -dm_forest_minimum_refinement 2 -dm_forest_initial_refinement 2 -ufv_use_amr -refine_vec_tagger_box 0.5,inf -coarsen_vec_tagger_box 0,1.e-2 -refine_tag_view -coarsen_tag_view -physics euler -eu_type iv_shock -ufv_cfl 10 -eu_alpha 60. -grid_skew_60 -eu_gamma 1.4 -eu_amach 2.02 -eu_rho2 3. -petscfv_type leastsquares -petsclimiter_type minmod -petscfv_compute_gradients 0 -ts_max_time 0.5 -ts_ssp_type rks2 -ts_ssp_nstages 10 -ufv_vtk_basename ${wPETSC_DIR}/ex11
2727     timeoutfactor: 3
2728 
2729   # Test GLVis visualization of PetscFV fields
2730   test:
2731     suffix: glvis_adv_2d_tet
2732     args: -ufv_vtk_interval 0 -ts_monitor_solution glvis: -ts_max_steps 0 -ufv_vtk_monitor 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_gmsh_periodic 0
2733 
2734   test:
2735     suffix: glvis_adv_2d_quad
2736     args: -ufv_vtk_interval 0 -ts_monitor_solution glvis: -ts_max_steps 0 -ufv_vtk_monitor 0 -dm_refine 5 -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3
2737 
2738   test:
2739     suffix: tut_1
2740     requires: exodusii
2741     nsize: 1
2742     args: -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo
2743 
2744   test:
2745     suffix: tut_2
2746     requires: exodusii
2747     nsize: 1
2748     args: -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo -ts_type rosw
2749 
2750   test:
2751     suffix: tut_3
2752     requires: exodusii
2753     nsize: 4
2754     args: -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -monitor Error -advect_sol_type bump -petscfv_type leastsquares -petsclimiter_type sin
2755 
2756   test:
2757     suffix: tut_4
2758     requires: exodusii
2759     nsize: 4
2760     args: -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -physics sw -monitor Height,Energy -petscfv_type leastsquares -petsclimiter_type minmod
2761 
2762 TEST*/
2763