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