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