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