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