xref: /petsc/src/ts/tutorials/ex18.c (revision f84028053d197aff3fdf060e86d88b3f1b0f0110)
1 static char help[] = "Hybrid Finite Element-Finite Volume Example.\n";
2 /*F
3   Here we are advecting a passive tracer in a harmonic velocity field, defined by
4 a forcing function $f$:
5 \begin{align}
6   -\Delta \mathbf{u} + f &= 0 \\
7   \frac{\partial\phi}{\partial t} + \nabla\cdot \phi \mathbf{u} &= 0
8 \end{align}
9 F*/
10 
11 #include <petscdmplex.h>
12 #include <petscds.h>
13 #include <petscts.h>
14 
15 #include <petsc/private/dmpleximpl.h> /* For DotD */
16 
17 typedef enum {VEL_ZERO, VEL_CONSTANT, VEL_HARMONIC, VEL_SHEAR} VelocityDistribution;
18 
19 typedef enum {ZERO, CONSTANT, GAUSSIAN, TILTED, DELTA} PorosityDistribution;
20 
21 static PetscErrorCode constant_u_2d(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
22 
23 /*
24   FunctionalFunc - Calculates the value of a functional of the solution at a point
25 
26   Input Parameters:
27 + dm   - The DM
28 . time - The TS time
29 . x    - The coordinates of the evaluation point
30 . u    - The field values at point x
31 - ctx  - A user context, or NULL
32 
33   Output Parameter:
34 . f    - The value of the functional at point x
35 
36 */
37 typedef PetscErrorCode (*FunctionalFunc)(DM, PetscReal, const PetscReal *, const PetscScalar *, PetscReal *, void *);
38 
39 typedef struct _n_Functional *Functional;
40 struct _n_Functional {
41   char          *name;
42   FunctionalFunc func;
43   void          *ctx;
44   PetscInt       offset;
45   Functional     next;
46 };
47 
48 typedef struct {
49   /* Problem definition */
50   PetscBool      useFV;             /* Use a finite volume scheme for advection */
51   PetscErrorCode (*initialGuess[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
52   VelocityDistribution velocityDist;
53   PorosityDistribution porosityDist;
54   PetscReal            inflowState;
55   PetscReal            source[3];
56   /* Monitoring */
57   PetscInt       numMonitorFuncs, maxMonitorFunc;
58   Functional    *monitorFuncs;
59   PetscInt       errorFunctional;
60   Functional     functionalRegistry;
61 } AppCtx;
62 
63 static  AppCtx *globalUser;
64 
65 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
66 {
67   const char    *velocityDist[4]  = {"zero", "constant", "harmonic", "shear"};
68   const char    *porosityDist[5]  = {"zero", "constant", "gaussian", "tilted", "delta"};
69   PetscInt       vd, pd, d;
70   PetscBool      flg;
71 
72   PetscFunctionBeginUser;
73   options->useFV        = PETSC_FALSE;
74   options->velocityDist = VEL_HARMONIC;
75   options->porosityDist = ZERO;
76   options->inflowState  = -2.0;
77   options->numMonitorFuncs = 0;
78   options->source[0]    = 0.5;
79   options->source[1]    = 0.5;
80   options->source[2]    = 0.5;
81 
82   PetscOptionsBegin(comm, "", "Magma Dynamics Options", "DMPLEX");
83   PetscCall(PetscOptionsBool("-use_fv", "Use the finite volume method for advection", "ex18.c", options->useFV, &options->useFV, NULL));
84   vd   = options->velocityDist;
85   PetscCall(PetscOptionsEList("-velocity_dist","Velocity distribution type","ex18.c",velocityDist,4,velocityDist[options->velocityDist],&vd,NULL));
86   options->velocityDist = (VelocityDistribution) vd;
87   pd   = options->porosityDist;
88   PetscCall(PetscOptionsEList("-porosity_dist","Initial porosity distribution type","ex18.c",porosityDist,5,porosityDist[options->porosityDist],&pd,NULL));
89   options->porosityDist = (PorosityDistribution) pd;
90   PetscCall(PetscOptionsReal("-inflow_state", "The inflow state", "ex18.c", options->inflowState, &options->inflowState, NULL));
91   d    = 2;
92   PetscCall(PetscOptionsRealArray("-source_loc", "The source location", "ex18.c", options->source, &d, &flg));
93   PetscCheck(!flg || d == 2,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Must give dim coordinates for the source location, not %" PetscInt_FMT, d);
94   PetscOptionsEnd();
95 
96   PetscFunctionReturn(0);
97 }
98 
99 static PetscErrorCode ProcessMonitorOptions(MPI_Comm comm, AppCtx *options)
100 {
101   Functional     func;
102   char          *names[256];
103   PetscInt       f;
104 
105   PetscFunctionBeginUser;
106   PetscOptionsBegin(comm, "", "Simulation Monitor Options", "DMPLEX");
107   options->numMonitorFuncs = PETSC_STATIC_ARRAY_LENGTH(names);
108   PetscCall(PetscOptionsStringArray("-monitor", "List of functionals to monitor", "", names, &options->numMonitorFuncs, NULL));
109   PetscCall(PetscMalloc1(options->numMonitorFuncs, &options->monitorFuncs));
110   for (f = 0; f < options->numMonitorFuncs; ++f) {
111     for (func = options->functionalRegistry; func; func = func->next) {
112       PetscBool match;
113 
114       PetscCall(PetscStrcasecmp(names[f], func->name, &match));
115       if (match) break;
116     }
117     PetscCheck(func,comm, PETSC_ERR_USER, "No known functional '%s'", names[f]);
118     options->monitorFuncs[f] = func;
119     /* Jed inserts a de-duplication of functionals here */
120     PetscCall(PetscFree(names[f]));
121   }
122   /* Find out the maximum index of any functional computed by a function we will be calling (even if we are not using it) */
123   options->maxMonitorFunc = -1;
124   for (func = options->functionalRegistry; func; func = func->next) {
125     for (f = 0; f < options->numMonitorFuncs; ++f) {
126       Functional call = options->monitorFuncs[f];
127 
128       if (func->func == call->func && func->ctx == call->ctx) options->maxMonitorFunc = PetscMax(options->maxMonitorFunc, func->offset);
129     }
130   }
131   PetscOptionsEnd();
132   PetscFunctionReturn(0);
133 }
134 
135 static PetscErrorCode FunctionalRegister(Functional *functionalRegistry, const char name[], PetscInt *offset, FunctionalFunc func, void *ctx)
136 {
137   Functional    *ptr, f;
138   PetscInt       lastoffset = -1;
139 
140   PetscFunctionBeginUser;
141   for (ptr = functionalRegistry; *ptr; ptr = &(*ptr)->next) lastoffset = (*ptr)->offset;
142   PetscCall(PetscNew(&f));
143   PetscCall(PetscStrallocpy(name, &f->name));
144   f->offset = lastoffset + 1;
145   f->func   = func;
146   f->ctx    = ctx;
147   f->next   = NULL;
148   *ptr      = f;
149   *offset   = f->offset;
150   PetscFunctionReturn(0);
151 }
152 
153 static PetscErrorCode FunctionalDestroy(Functional *link)
154 {
155   Functional     next, l;
156 
157   PetscFunctionBeginUser;
158   if (!link) PetscFunctionReturn(0);
159   l     = *link;
160   *link = NULL;
161   for (; l; l=next) {
162     next = l->next;
163     PetscCall(PetscFree(l->name));
164     PetscCall(PetscFree(l));
165   }
166   PetscFunctionReturn(0);
167 }
168 
169 static void f0_zero_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
170                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
171                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
172                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
173 {
174   PetscInt comp;
175   for (comp = 0; comp < dim; ++comp) f0[comp] = u[comp];
176 }
177 
178 static void f0_constant_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
179                           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
180                           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
181                           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
182 {
183   PetscScalar wind[3] = {0.0, 0.0, 0.0};
184   PetscInt    comp;
185 
186   constant_u_2d(dim, t, x, Nf, wind, NULL);
187   for (comp = 0; comp < dim && comp < 3; ++comp) f0[comp] = u[comp] - wind[comp];
188 }
189 
190 static void f1_constant_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
191                           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
192                           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
193                           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
194 {
195   PetscInt comp;
196   for (comp = 0; comp < dim*dim; ++comp) f1[comp] = 0.0;
197 }
198 
199 static void g0_constant_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
200                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
201                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
202                            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
203 {
204   PetscInt d;
205   for (d = 0; d < dim; ++d) g0[d*dim+d] = 1.0;
206 }
207 
208 static void g0_constant_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
209                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
210                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
211                            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
212 {
213   g0[0] = 1.0;
214 }
215 
216 static void f0_lap_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
217                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
218                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
219                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
220 {
221   PetscInt comp;
222   for (comp = 0; comp < dim; ++comp) f0[comp] = 4.0;
223 }
224 
225 static void f1_lap_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
226                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
227                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
228                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
229 {
230   PetscInt comp, d;
231   for (comp = 0; comp < dim; ++comp) {
232     for (d = 0; d < dim; ++d) {
233       f1[comp*dim+d] = u_x[comp*dim+d];
234     }
235   }
236 }
237 
238 static void f0_lap_periodic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
239                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
240                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
241                               PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
242 {
243   f0[0] = -PetscSinReal(2.0*PETSC_PI*x[0]);
244   f0[1] = 2.0*PETSC_PI*x[1]*PetscCosReal(2.0*PETSC_PI*x[0]);
245 }
246 
247 static void f0_lap_doubly_periodic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
248                                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
249                                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
250                                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
251 {
252   f0[0] = -2.0*PetscSinReal(2.0*PETSC_PI*x[0])*PetscCosReal(2.0*PETSC_PI*x[1]);
253   f0[1] =  2.0*PetscSinReal(2.0*PETSC_PI*x[1])*PetscCosReal(2.0*PETSC_PI*x[0]);
254 }
255 
256 void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
257            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
258            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
259            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
260 {
261   const PetscInt Ncomp = dim;
262   PetscInt       compI, d;
263 
264   for (compI = 0; compI < Ncomp; ++compI) {
265     for (d = 0; d < dim; ++d) {
266       g3[((compI*Ncomp+compI)*dim+d)*dim+d] = 1.0;
267     }
268   }
269 }
270 
271 /* \frac{\partial\phi}{\partial t} + \nabla\phi \cdot \mathbf{u} + \phi \nabla \cdot \mathbf{u} = 0 */
272 static void f0_advection(PetscInt dim, PetscInt Nf, PetscInt NfAux,
273                          const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
274                          const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
275                          PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
276 {
277   PetscInt d;
278   f0[0] = u_t[dim];
279   for (d = 0; d < dim; ++d) f0[0] += u[dim]*u_x[d*dim+d] + u_x[dim*dim+d]*u[d];
280 }
281 
282 static void f1_advection(PetscInt dim, PetscInt Nf, PetscInt NfAux,
283                          const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
284                          const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
285                          PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
286 {
287   PetscInt d;
288   for (d = 0; d < dim; ++d) f1[0] = 0.0;
289 }
290 
291 void g0_adv_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
292                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
293                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
294                PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
295 {
296   PetscInt d;
297   g0[0] = u_tShift;
298   for (d = 0; d < dim; ++d) g0[0] += u_x[d*dim+d];
299 }
300 
301 void g1_adv_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
302                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
303                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
304                PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
305 {
306   PetscInt d;
307   for (d = 0; d < dim; ++d) g1[d] = u[d];
308 }
309 
310 void g0_adv_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
311                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
312                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
313                PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
314 {
315   PetscInt d;
316   for (d = 0; d < dim; ++d) g0[0] += u_x[dim*dim+d];
317 }
318 
319 void g1_adv_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
320                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
321                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
322                PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
323 {
324   PetscInt d;
325   for (d = 0; d < dim; ++d) g1[d*dim+d] = u[dim];
326 }
327 
328 static void riemann_advection(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *uL, const PetscScalar *uR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, void *ctx)
329 {
330   PetscReal wind[3] = {0.0, 1.0, 0.0};
331   PetscReal wn = DMPlex_DotRealD_Internal(PetscMin(dim,3), wind, n);
332 
333   flux[0] = (wn > 0 ? uL[dim] : uR[dim]) * wn;
334 }
335 
336 static void riemann_coupled_advection(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *uL, const PetscScalar *uR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, void *ctx)
337 {
338   PetscReal wn = DMPlex_DotD_Internal(dim, uL, n);
339 
340 #if 1
341   flux[0] = (wn > 0 ? uL[dim] : uR[dim]) * wn;
342 #else
343   /* if (fabs(uL[0] - wind[0]) > 1.0e-7 || fabs(uL[1] - wind[1]) > 1.0e-7) PetscPrintf(PETSC_COMM_SELF, "wind (%g, %g) uL (%g, %g) uR (%g, %g)\n", wind[0], wind[1], uL[0], uL[1], uR[0], uR[1]); */
344   /* Smear it out */
345   flux[0] = 0.5*((uL[dim] + uR[dim]) + (uL[dim] - uR[dim])*tanh(1.0e5*wn)) * wn;
346 #endif
347 }
348 
349 static PetscErrorCode zero_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
350 {
351   u[0] = 0.0;
352   u[1] = 0.0;
353   return 0;
354 }
355 
356 static PetscErrorCode constant_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
357 {
358   u[0] = 0.0;
359   u[1] = 1.0;
360   return 0;
361 }
362 
363 /* Coordinates of the point which was at x at t = 0 */
364 static PetscErrorCode constant_x_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
365 {
366   const PetscReal t = *((PetscReal *) ctx);
367   u[0] = x[0];
368   u[1] = x[1] + t;
369 #if 0
370   PetscCall(DMLocalizeCoordinate(globalUser->dm, u, PETSC_FALSE, u));
371 #else
372   u[1] = u[1] - (int) PetscRealPart(u[1]);
373 #endif
374   return 0;
375 }
376 
377 /*
378   In 2D we use the exact solution:
379 
380     u   = x^2 + y^2
381     v   = 2 x^2 - 2xy
382     phi = h(x + y + (u + v) t)
383     f_x = f_y = 4
384 
385   so that
386 
387     -\Delta u + f = <-4, -4> + <4, 4> = 0
388     {\partial\phi}{\partial t} - \nabla\cdot \phi u = 0
389     h_t(x + y + (u + v) t) - u . grad phi - phi div u
390   = u h' + v h'              - u h_x - v h_y
391   = 0
392 
393 We will conserve phi since
394 
395     \nabla \cdot u = 2x - 2x = 0
396 
397 Also try h((x + ut)^2 + (y + vt)^2), so that
398 
399     h_t((x + ut)^2 + (y + vt)^2) - u . grad phi - phi div u
400   = 2 h' (u (x + ut) + v (y + vt)) - u h_x - v h_y
401   = 2 h' (u (x + ut) + v (y + vt)) - u h' 2 (x + u t) - v h' 2 (y + vt)
402   = 2 h' (u (x + ut) + v (y + vt)  - u (x + u t) - v (y + vt))
403   = 0
404 
405 */
406 static PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
407 {
408   u[0] = x[0]*x[0] + x[1]*x[1];
409   u[1] = 2.0*x[0]*x[0] - 2.0*x[0]*x[1];
410   return 0;
411 }
412 
413 /*
414   In 2D we use the exact, periodic solution:
415 
416     u   =  sin(2 pi x)/4 pi^2
417     v   = -y cos(2 pi x)/2 pi
418     phi = h(x + y + (u + v) t)
419     f_x = -sin(2 pi x)
420     f_y = 2 pi y cos(2 pi x)
421 
422   so that
423 
424     -\Delta u + f = <sin(2pi x),  -2pi y cos(2pi x)> + <-sin(2pi x), 2pi y cos(2pi x)> = 0
425 
426 We will conserve phi since
427 
428     \nabla \cdot u = cos(2pi x)/2pi - cos(2pi x)/2pi = 0
429 */
430 static PetscErrorCode periodic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
431 {
432   u[0] = PetscSinReal(2.0*PETSC_PI*x[0])/PetscSqr(2.0*PETSC_PI);
433   u[1] = -x[1]*PetscCosReal(2.0*PETSC_PI*x[0])/(2.0*PETSC_PI);
434   return 0;
435 }
436 
437 /*
438   In 2D we use the exact, doubly periodic solution:
439 
440     u   =  sin(2 pi x) cos(2 pi y)/4 pi^2
441     v   = -sin(2 pi y) cos(2 pi x)/4 pi^2
442     phi = h(x + y + (u + v) t)
443     f_x = -2sin(2 pi x) cos(2 pi y)
444     f_y =  2sin(2 pi y) cos(2 pi x)
445 
446   so that
447 
448     -\Delta u + f = <2 sin(2pi x) cos(2pi y),  -2 sin(2pi y) cos(2pi x)> + <-2 sin(2pi x) cos(2pi y), 2 sin(2pi y) cos(2pi x)> = 0
449 
450 We will conserve phi since
451 
452     \nabla \cdot u = cos(2pi x) cos(2pi y)/2pi - cos(2pi y) cos(2pi x)/2pi = 0
453 */
454 static PetscErrorCode doubly_periodic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
455 {
456   u[0] =  PetscSinReal(2.0*PETSC_PI*x[0])*PetscCosReal(2.0*PETSC_PI*x[1])/PetscSqr(2.0*PETSC_PI);
457   u[1] = -PetscSinReal(2.0*PETSC_PI*x[1])*PetscCosReal(2.0*PETSC_PI*x[0])/PetscSqr(2.0*PETSC_PI);
458   return 0;
459 }
460 
461 static PetscErrorCode shear_bc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
462 {
463   u[0] = x[1] - 0.5;
464   u[1] = 0.0;
465   return 0;
466 }
467 
468 static PetscErrorCode initialVelocity(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
469 {
470   PetscInt d;
471   for (d = 0; d < dim; ++d) u[d] = 0.0;
472   return 0;
473 }
474 
475 static PetscErrorCode zero_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
476 {
477   u[0] = 0.0;
478   return 0;
479 }
480 
481 static PetscErrorCode constant_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
482 {
483   u[0] = 1.0;
484   return 0;
485 }
486 
487 static PetscErrorCode delta_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
488 {
489   PetscReal   x0[2];
490   PetscScalar xn[2];
491 
492   x0[0] = globalUser->source[0];
493   x0[1] = globalUser->source[1];
494   constant_x_2d(dim, time, x0, Nf, xn, ctx);
495   {
496     const PetscReal xi  = x[0] - PetscRealPart(xn[0]);
497     const PetscReal eta = x[1] - PetscRealPart(xn[1]);
498     const PetscReal r2  = xi*xi + eta*eta;
499 
500     u[0] = r2 < 1.0e-7 ? 1.0 : 0.0;
501   }
502   return 0;
503 }
504 
505 /*
506   Gaussian blob, initially centered on (0.5, 0.5)
507 
508   xi = x(t) - x0, eta = y(t) - y0
509 
510 where x(t), y(t) are the integral curves of v(t),
511 
512   dx/dt . grad f = v . f
513 
514 Check: constant v(t) = {v0, w0}, x(t) = {x0 + v0 t, y0 + w0 t}
515 
516   v0 f_x + w0 f_y = v . f
517 */
518 static PetscErrorCode gaussian_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
519 {
520   const PetscReal x0[2] = {0.5, 0.5};
521   const PetscReal sigma = 1.0/6.0;
522   PetscScalar     xn[2];
523 
524   constant_x_2d(dim, time, x0, Nf, xn, ctx);
525   {
526     /* const PetscReal xi  = x[0] + (sin(2.0*PETSC_PI*x[0])/(4.0*PETSC_PI*PETSC_PI))*t - x0[0]; */
527     /* const PetscReal eta = x[1] + (-x[1]*cos(2.0*PETSC_PI*x[0])/(2.0*PETSC_PI))*t - x0[1]; */
528     const PetscReal xi  = x[0] - PetscRealPart(xn[0]);
529     const PetscReal eta = x[1] - PetscRealPart(xn[1]);
530     const PetscReal r2  = xi*xi + eta*eta;
531 
532     u[0] = PetscExpReal(-r2/(2.0*sigma*sigma))/(sigma*PetscSqrtReal(2.0*PETSC_PI));
533   }
534   return 0;
535 }
536 
537 static PetscErrorCode tilted_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
538 {
539   PetscReal       x0[3];
540   const PetscReal wind[3] = {0.0, 1.0, 0.0};
541   const PetscReal t       = *((PetscReal *) ctx);
542 
543   DMPlex_WaxpyD_Internal(2, -t, wind, x, x0);
544   if (x0[1] > 0) u[0] =  1.0*x[0] + 3.0*x[1];
545   else           u[0] = -2.0; /* Inflow state */
546   return 0;
547 }
548 
549 static PetscErrorCode tilted_phi_coupled_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
550 {
551   PetscReal       ur[3];
552   PetscReal       x0[3];
553   const PetscReal t = *((PetscReal *) ctx);
554 
555   ur[0] = PetscRealPart(u[0]); ur[1] = PetscRealPart(u[1]); ur[2] = PetscRealPart(u[2]);
556   DMPlex_WaxpyD_Internal(2, -t, ur, x, x0);
557   if (x0[1] > 0) u[0] =  1.0*x[0] + 3.0*x[1];
558   else           u[0] = -2.0; /* Inflow state */
559   return 0;
560 }
561 
562 static PetscErrorCode advect_inflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
563 {
564   AppCtx *user = (AppCtx *) ctx;
565 
566   PetscFunctionBeginUser;
567   xG[0] = user->inflowState;
568   PetscFunctionReturn(0);
569 }
570 
571 static PetscErrorCode advect_outflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
572 {
573   PetscFunctionBeginUser;
574   //xG[0] = xI[dim];
575   xG[0] = xI[2];
576   PetscFunctionReturn(0);
577 }
578 
579 static PetscErrorCode ExactSolution(DM dm, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx)
580 {
581   AppCtx        *user = (AppCtx *) ctx;
582   PetscInt       dim;
583 
584   PetscFunctionBeginUser;
585   PetscCall(DMGetDimension(dm, &dim));
586   switch (user->porosityDist) {
587   case TILTED:
588     if (user->velocityDist == VEL_ZERO) tilted_phi_2d(dim, time, x, 2, u, (void *) &time);
589     else                                tilted_phi_coupled_2d(dim, time, x, 2, u, (void *) &time);
590     break;
591   case GAUSSIAN:
592     gaussian_phi_2d(dim, time, x, 2, u, (void *) &time);
593     break;
594   case DELTA:
595     delta_phi_2d(dim, time, x, 2, u, (void *) &time);
596     break;
597   default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown solution type");
598   }
599   PetscFunctionReturn(0);
600 }
601 
602 static PetscErrorCode Functional_Error(DM dm, PetscReal time, const PetscReal *x, const PetscScalar *y, PetscReal *f, void *ctx)
603 {
604   AppCtx        *user = (AppCtx *) ctx;
605   PetscScalar    yexact[3]={0,0,0};
606 
607   PetscFunctionBeginUser;
608   PetscCall(ExactSolution(dm, time, x, yexact, ctx));
609   f[user->errorFunctional] = PetscAbsScalar(y[0] - yexact[0]);
610   PetscFunctionReturn(0);
611 }
612 
613 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
614 {
615   PetscFunctionBeginUser;
616   PetscCall(DMCreate(comm, dm));
617   PetscCall(DMSetType(*dm, DMPLEX));
618 #if 0
619   PetscBool       periodic = user->bd[0] == DM_BOUNDARY_PERIODIC || user->bd[0] == DM_BOUNDARY_TWIST || user->bd[1] == DM_BOUNDARY_PERIODIC || user->bd[1] == DM_BOUNDARY_TWIST ? PETSC_TRUE : PETSC_FALSE;
620   const PetscReal L[3]     = {1.0, 1.0, 1.0};
621   PetscReal       maxCell[3];
622   PetscInt        d;
623 
624   if (periodic) {for (d = 0; d < 3; ++d) maxCell[d] = 1.1*(L[d]/cells[d]); PetscCall(DMSetPeriodicity(*dm, PETSC_TRUE, maxCell, L, user->bd));}
625 #endif
626   PetscCall(DMSetFromOptions(*dm));
627   PetscCall(DMViewFromOptions(*dm, NULL, "-orig_dm_view"));
628   PetscFunctionReturn(0);
629 }
630 
631 static PetscErrorCode SetupBC(DM dm, AppCtx *user)
632 {
633   PetscErrorCode (*exactFuncs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
634   DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};
635   PetscDS        prob;
636   DMLabel        label;
637   PetscBool      check;
638   PetscInt       dim, n = 3;
639   const char    *prefix;
640 
641   PetscFunctionBeginUser;
642   PetscCall(PetscObjectGetOptionsPrefix((PetscObject) dm, &prefix));
643   PetscCall(PetscOptionsGetEnumArray(NULL, prefix, "-dm_plex_box_bd", DMBoundaryTypes, (PetscEnum *) bdt, &n, NULL));
644   PetscCall(DMGetDimension(dm, &dim));
645   /* Set initial guesses and exact solutions */
646   switch (dim) {
647     case 2:
648       user->initialGuess[0] = initialVelocity;
649       switch(user->porosityDist) {
650         case ZERO:     user->initialGuess[1] = zero_phi;break;
651         case CONSTANT: user->initialGuess[1] = constant_phi;break;
652         case GAUSSIAN: user->initialGuess[1] = gaussian_phi_2d;break;
653         case DELTA:    user->initialGuess[1] = delta_phi_2d;break;
654         case TILTED:
655         if (user->velocityDist == VEL_ZERO) user->initialGuess[1] = tilted_phi_2d;
656         else                                user->initialGuess[1] = tilted_phi_coupled_2d;
657         break;
658       }
659       break;
660     default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Dimension %" PetscInt_FMT " not supported", dim);
661   }
662   exactFuncs[0] = user->initialGuess[0];
663   exactFuncs[1] = user->initialGuess[1];
664   switch (dim) {
665     case 2:
666       switch (user->velocityDist) {
667         case VEL_ZERO:
668           exactFuncs[0] = zero_u_2d; break;
669         case VEL_CONSTANT:
670           exactFuncs[0] = constant_u_2d; break;
671         case VEL_HARMONIC:
672           switch (bdt[0]) {
673             case DM_BOUNDARY_PERIODIC:
674               switch (bdt[1]) {
675                 case DM_BOUNDARY_PERIODIC:
676                   exactFuncs[0] = doubly_periodic_u_2d; break;
677                 default:
678                   exactFuncs[0] = periodic_u_2d; break;
679               }
680               break;
681             default:
682               exactFuncs[0] = quadratic_u_2d; break;
683           }
684           break;
685         case VEL_SHEAR:
686           exactFuncs[0] = shear_bc; break;
687         default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %" PetscInt_FMT, dim);
688       }
689       break;
690     default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Dimension %" PetscInt_FMT " not supported", dim);
691   }
692   {
693     PetscBool isImplicit = PETSC_FALSE;
694 
695     PetscCall(PetscOptionsHasName(NULL,"", "-use_implicit", &isImplicit));
696     if (user->velocityDist == VEL_CONSTANT && !isImplicit) user->initialGuess[0] = exactFuncs[0];
697   }
698   PetscCall(PetscOptionsHasName(NULL,NULL, "-dmts_check", &check));
699   if (check) {
700     user->initialGuess[0] = exactFuncs[0];
701     user->initialGuess[1] = exactFuncs[1];
702   }
703   /* Set BC */
704   PetscCall(DMGetDS(dm, &prob));
705   PetscCall(DMGetLabel(dm, "marker", &label));
706   PetscCall(PetscDSSetExactSolution(prob, 0, exactFuncs[0], user));
707   PetscCall(PetscDSSetExactSolution(prob, 1, exactFuncs[1], user));
708   if (label) {
709     const PetscInt id = 1;
710 
711     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, user, NULL));
712   }
713   PetscCall(DMGetLabel(dm, "Face Sets", &label));
714   if (label && user->useFV) {
715     const PetscInt inflowids[] = {100,200,300}, outflowids[] = {101};
716 
717     PetscCall(DMAddBoundary(dm, DM_BC_NATURAL_RIEMANN, "inflow",  label,  PETSC_STATIC_ARRAY_LENGTH(inflowids),  inflowids, 1, 0, NULL, (void (*)(void)) advect_inflow, NULL, user, NULL));
718     PetscCall(DMAddBoundary(dm, DM_BC_NATURAL_RIEMANN, "outflow", label, PETSC_STATIC_ARRAY_LENGTH(outflowids), outflowids, 1, 0, NULL, (void (*)(void)) advect_outflow, NULL, user, NULL));
719   }
720   PetscFunctionReturn(0);
721 }
722 
723 static PetscErrorCode SetupProblem(DM dm, AppCtx *user)
724 {
725   DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};
726   PetscDS        prob;
727   PetscInt       n = 3;
728   const char    *prefix;
729 
730   PetscFunctionBeginUser;
731   PetscCall(PetscObjectGetOptionsPrefix((PetscObject) dm, &prefix));
732   PetscCall(PetscOptionsGetEnumArray(NULL, prefix, "-dm_plex_box_bd", DMBoundaryTypes, (PetscEnum *) bdt, &n, NULL));
733   PetscCall(DMGetDS(dm, &prob));
734   switch (user->velocityDist) {
735   case VEL_ZERO:
736     PetscCall(PetscDSSetResidual(prob, 0, f0_zero_u, f1_constant_u));
737     break;
738   case VEL_CONSTANT:
739     PetscCall(PetscDSSetResidual(prob, 0, f0_constant_u, f1_constant_u));
740     PetscCall(PetscDSSetJacobian(prob, 0, 0, g0_constant_uu, NULL, NULL, NULL));
741     PetscCall(PetscDSSetJacobian(prob, 1, 1, g0_constant_pp, NULL, NULL, NULL));
742     break;
743   case VEL_HARMONIC:
744     switch (bdt[0]) {
745     case DM_BOUNDARY_PERIODIC:
746       switch (bdt[1]) {
747       case DM_BOUNDARY_PERIODIC:
748         PetscCall(PetscDSSetResidual(prob, 0, f0_lap_doubly_periodic_u, f1_lap_u));
749         break;
750       default:
751         PetscCall(PetscDSSetResidual(prob, 0, f0_lap_periodic_u, f1_lap_u));
752         break;
753       }
754       break;
755     default:
756       PetscCall(PetscDSSetResidual(prob, 0, f0_lap_u, f1_lap_u));
757       break;
758     }
759     PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu));
760     break;
761   case VEL_SHEAR:
762     PetscCall(PetscDSSetResidual(prob, 0, f0_zero_u, f1_lap_u));
763     PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu));
764     break;
765   }
766   PetscCall(PetscDSSetResidual(prob, 1, f0_advection, f1_advection));
767   PetscCall(PetscDSSetJacobian(prob, 1, 1, g0_adv_pp, g1_adv_pp, NULL, NULL));
768   PetscCall(PetscDSSetJacobian(prob, 1, 0, g0_adv_pu, g1_adv_pu, NULL, NULL));
769   if (user->velocityDist == VEL_ZERO) PetscCall(PetscDSSetRiemannSolver(prob, 1, riemann_advection));
770   else                                PetscCall(PetscDSSetRiemannSolver(prob, 1, riemann_coupled_advection));
771 
772   PetscCall(FunctionalRegister(&user->functionalRegistry, "Error", &user->errorFunctional, Functional_Error, user));
773   PetscFunctionReturn(0);
774 }
775 
776 static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
777 {
778   DM              cdm = dm;
779   PetscQuadrature q;
780   PetscFE         fe[2];
781   PetscFV         fv;
782   MPI_Comm        comm;
783   PetscInt        dim;
784 
785   PetscFunctionBeginUser;
786   /* Create finite element */
787   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
788   PetscCall(DMGetDimension(dm, &dim));
789   PetscCall(PetscFECreateDefault(comm, dim, dim, PETSC_FALSE, "velocity_", PETSC_DEFAULT, &fe[0]));
790   PetscCall(PetscObjectSetName((PetscObject) fe[0], "velocity"));
791   PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_FALSE, "porosity_", PETSC_DEFAULT, &fe[1]));
792   PetscCall(PetscFECopyQuadrature(fe[0], fe[1]));
793   PetscCall(PetscObjectSetName((PetscObject) fe[1], "porosity"));
794 
795   PetscCall(PetscFVCreate(PetscObjectComm((PetscObject) dm), &fv));
796   PetscCall(PetscObjectSetName((PetscObject) fv, "porosity"));
797   PetscCall(PetscFVSetFromOptions(fv));
798   PetscCall(PetscFVSetNumComponents(fv, 1));
799   PetscCall(PetscFVSetSpatialDimension(fv, dim));
800   PetscCall(PetscFEGetQuadrature(fe[0], &q));
801   PetscCall(PetscFVSetQuadrature(fv, q));
802 
803   PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe[0]));
804   if (user->useFV) PetscCall(DMSetField(dm, 1, NULL, (PetscObject) fv));
805   else             PetscCall(DMSetField(dm, 1, NULL, (PetscObject) fe[1]));
806   PetscCall(DMCreateDS(dm));
807   PetscCall(SetupProblem(dm, user));
808 
809   /* Set discretization and boundary conditions for each mesh */
810   while (cdm) {
811     PetscCall(DMCopyDisc(dm, cdm));
812     PetscCall(DMGetCoarseDM(cdm, &cdm));
813     /* Coordinates were never localized for coarse meshes */
814     if (cdm) PetscCall(DMLocalizeCoordinates(cdm));
815   }
816   PetscCall(PetscFEDestroy(&fe[0]));
817   PetscCall(PetscFEDestroy(&fe[1]));
818   PetscCall(PetscFVDestroy(&fv));
819   PetscFunctionReturn(0);
820 }
821 
822 static PetscErrorCode CreateDM(MPI_Comm comm, AppCtx *user, DM *dm)
823 {
824   PetscFunctionBeginUser;
825   PetscCall(CreateMesh(comm, user, dm));
826   /* Handle refinement, etc. */
827   PetscCall(DMSetFromOptions(*dm));
828   /* Construct ghost cells */
829   if (user->useFV) {
830     DM gdm;
831 
832     PetscCall(DMPlexConstructGhostCells(*dm, NULL, NULL, &gdm));
833     PetscCall(DMDestroy(dm));
834     *dm  = gdm;
835   }
836   /* Localize coordinates */
837   PetscCall(DMLocalizeCoordinates(*dm));
838   PetscCall(PetscObjectSetName((PetscObject)(*dm),"Mesh"));
839   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
840   /* Setup problem */
841   PetscCall(SetupDiscretization(*dm, user));
842   /* Setup BC */
843   PetscCall(SetupBC(*dm, user));
844   PetscFunctionReturn(0);
845 }
846 
847 static PetscErrorCode SetInitialConditionFVM(DM dm, Vec X, PetscInt field, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *ctx)
848 {
849   PetscDS            prob;
850   DM                 dmCell;
851   Vec                cellgeom;
852   const PetscScalar *cgeom;
853   PetscScalar       *x;
854   PetscInt           dim, Nf, cStart, cEnd, c;
855 
856   PetscFunctionBeginUser;
857   PetscCall(DMGetDS(dm, &prob));
858   PetscCall(DMGetDimension(dm, &dim));
859   PetscCall(PetscDSGetNumFields(prob, &Nf));
860   PetscCall(DMPlexGetGeometryFVM(dm, NULL, &cellgeom, NULL));
861   PetscCall(VecGetDM(cellgeom, &dmCell));
862   PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd));
863   PetscCall(VecGetArrayRead(cellgeom, &cgeom));
864   PetscCall(VecGetArray(X, &x));
865   for (c = cStart; c < cEnd; ++c) {
866     PetscFVCellGeom       *cg;
867     PetscScalar           *xc;
868 
869     PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg));
870     PetscCall(DMPlexPointGlobalFieldRef(dm, c, field, x, &xc));
871     if (xc) PetscCall((*func)(dim, 0.0, cg->centroid, Nf, xc, ctx));
872   }
873   PetscCall(VecRestoreArrayRead(cellgeom, &cgeom));
874   PetscCall(VecRestoreArray(X, &x));
875   PetscFunctionReturn(0);
876 }
877 
878 static PetscErrorCode MonitorFunctionals(TS ts, PetscInt stepnum, PetscReal time, Vec X, void *ctx)
879 {
880   AppCtx            *user   = (AppCtx *) ctx;
881   char              *ftable = NULL;
882   DM                 dm;
883   PetscSection       s;
884   Vec                cellgeom;
885   const PetscScalar *x;
886   PetscScalar       *a;
887   PetscReal         *xnorms;
888   PetscInt           pStart, pEnd, p, Nf, f;
889 
890   PetscFunctionBeginUser;
891   PetscCall(VecViewFromOptions(X, (PetscObject) ts, "-view_solution"));
892   PetscCall(VecGetDM(X, &dm));
893   PetscCall(DMPlexGetGeometryFVM(dm, NULL, &cellgeom, NULL));
894   PetscCall(DMGetLocalSection(dm, &s));
895   PetscCall(PetscSectionGetNumFields(s, &Nf));
896   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
897   PetscCall(PetscCalloc1(Nf*2, &xnorms));
898   PetscCall(VecGetArrayRead(X, &x));
899   for (p = pStart; p < pEnd; ++p) {
900     for (f = 0; f < Nf; ++f) {
901       PetscInt dof, cdof, d;
902 
903       PetscCall(PetscSectionGetFieldDof(s, p, f, &dof));
904       PetscCall(PetscSectionGetFieldConstraintDof(s, p, f, &cdof));
905       PetscCall(DMPlexPointGlobalFieldRead(dm, p, f, x, &a));
906       /* TODO Use constrained indices here */
907       for (d = 0; d < dof-cdof; ++d) xnorms[f*2+0]  = PetscMax(xnorms[f*2+0], PetscAbsScalar(a[d]));
908       for (d = 0; d < dof-cdof; ++d) xnorms[f*2+1] += PetscAbsScalar(a[d]);
909     }
910   }
911   PetscCall(VecRestoreArrayRead(X, &x));
912   if (stepnum >= 0) { /* No summary for final time */
913     DM                 dmCell, *fdm;
914     Vec               *fv;
915     const PetscScalar *cgeom;
916     PetscScalar      **fx;
917     PetscReal         *fmin, *fmax, *fint, *ftmp, t;
918     PetscInt           cStart, cEnd, c, fcount, f, num;
919 
920     size_t             ftableused,ftablealloc;
921 
922     /* Functionals have indices after registering, this is an upper bound */
923     fcount = user->numMonitorFuncs;
924     PetscCall(PetscMalloc4(fcount,&fmin,fcount,&fmax,fcount,&fint,fcount,&ftmp));
925     PetscCall(PetscMalloc3(fcount,&fdm,fcount,&fv,fcount,&fx));
926     for (f = 0; f < fcount; ++f) {
927       PetscSection fs;
928       const char  *name = user->monitorFuncs[f]->name;
929 
930       fmin[f] = PETSC_MAX_REAL;
931       fmax[f] = PETSC_MIN_REAL;
932       fint[f] = 0;
933       /* Make monitor vecs */
934       PetscCall(DMClone(dm, &fdm[f]));
935       PetscCall(DMGetOutputSequenceNumber(dm, &num, &t));
936       PetscCall(DMSetOutputSequenceNumber(fdm[f], num, t));
937       PetscCall(PetscSectionClone(s, &fs));
938       PetscCall(PetscSectionSetFieldName(fs, 0, NULL));
939       PetscCall(PetscSectionSetFieldName(fs, 1, name));
940       PetscCall(DMSetLocalSection(fdm[f], fs));
941       PetscCall(PetscSectionDestroy(&fs));
942       PetscCall(DMGetGlobalVector(fdm[f], &fv[f]));
943       PetscCall(PetscObjectSetName((PetscObject) fv[f], name));
944       PetscCall(VecGetArray(fv[f], &fx[f]));
945     }
946     PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd));
947     PetscCall(VecGetDM(cellgeom, &dmCell));
948     PetscCall(VecGetArrayRead(cellgeom, &cgeom));
949     PetscCall(VecGetArrayRead(X, &x));
950     for (c = cStart; c < cEnd; ++c) {
951       PetscFVCellGeom *cg;
952       PetscScalar     *cx;
953 
954       PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg));
955       PetscCall(DMPlexPointGlobalFieldRead(dm, c, 1, x, &cx));
956       if (!cx) continue;        /* not a global cell */
957       for (f = 0;  f < user->numMonitorFuncs; ++f) {
958         Functional   func = user->monitorFuncs[f];
959         PetscScalar *fxc;
960 
961         PetscCall(DMPlexPointGlobalFieldRef(dm, c, 1, fx[f], &fxc));
962         /* I need to make it easier to get interpolated values here */
963         PetscCall((*func->func)(dm, time, cg->centroid, cx, ftmp, func->ctx));
964         fxc[0] = ftmp[user->monitorFuncs[f]->offset];
965       }
966       for (f = 0; f < fcount; ++f) {
967         fmin[f]  = PetscMin(fmin[f], ftmp[f]);
968         fmax[f]  = PetscMax(fmax[f], ftmp[f]);
969         fint[f] += cg->volume * ftmp[f];
970       }
971     }
972     PetscCall(VecRestoreArrayRead(cellgeom, &cgeom));
973     PetscCall(VecRestoreArrayRead(X, &x));
974     PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, fmin, fcount, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)ts)));
975     PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, fmax, fcount, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)ts)));
976     PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, fint, fcount, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)ts)));
977     /* Output functional data */
978     ftablealloc = fcount * 100;
979     ftableused  = 0;
980     PetscCall(PetscCalloc1(ftablealloc, &ftable));
981     for (f = 0; f < user->numMonitorFuncs; ++f) {
982       Functional func      = user->monitorFuncs[f];
983       PetscInt   id        = func->offset;
984       char       newline[] = "\n";
985       char       buffer[256], *p, *prefix;
986       size_t     countused, len;
987 
988       /* Create string with functional outputs */
989       if (f % 3) {
990         PetscCall(PetscArraycpy(buffer, "  ", 2));
991         p    = buffer + 2;
992       } else if (f) {
993         PetscCall(PetscArraycpy(buffer, newline, sizeof(newline)-1));
994         p    = buffer + sizeof(newline) - 1;
995       } else {
996         p = buffer;
997       }
998       PetscCall(PetscSNPrintfCount(p, sizeof buffer-(p-buffer), "%12s [%12.6g,%12.6g] int %12.6g", &countused, func->name, (double) fmin[id], (double) fmax[id], (double) fint[id]));
999       countused += p - buffer;
1000       /* reallocate */
1001       if (countused > ftablealloc-ftableused-1) {
1002         char *ftablenew;
1003 
1004         ftablealloc = 2*ftablealloc + countused;
1005         PetscCall(PetscMalloc1(ftablealloc, &ftablenew));
1006         PetscCall(PetscArraycpy(ftablenew, ftable, ftableused));
1007         PetscCall(PetscFree(ftable));
1008         ftable = ftablenew;
1009       }
1010       PetscCall(PetscArraycpy(ftable+ftableused, buffer, countused));
1011       ftableused += countused;
1012       ftable[ftableused] = 0;
1013       /* Output vecs */
1014       PetscCall(VecRestoreArray(fv[f], &fx[f]));
1015       PetscCall(PetscStrlen(func->name, &len));
1016       PetscCall(PetscMalloc1(len+2,&prefix));
1017       PetscCall(PetscStrcpy(prefix, func->name));
1018       PetscCall(PetscStrcat(prefix, "_"));
1019       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fv[f], prefix));
1020       PetscCall(VecViewFromOptions(fv[f], NULL, "-vec_view"));
1021       PetscCall(PetscFree(prefix));
1022       PetscCall(DMRestoreGlobalVector(fdm[f], &fv[f]));
1023       PetscCall(DMDestroy(&fdm[f]));
1024     }
1025     PetscCall(PetscFree4(fmin, fmax, fint, ftmp));
1026     PetscCall(PetscFree3(fdm, fv, fx));
1027     PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), "% 3" PetscInt_FMT "  time %8.4g  |x| (", stepnum, (double) time));
1028     for (f = 0; f < Nf; ++f) {
1029       if (f > 0) PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), ", "));
1030       PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), "%8.4g", (double) xnorms[f*2+0]));
1031     }
1032     PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), ") |x|_1 ("));
1033     for (f = 0; f < Nf; ++f) {
1034       if (f > 0) PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), ", "));
1035       PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), "%8.4g", (double) xnorms[f*2+1]));
1036     }
1037     PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), ")  %s\n", ftable ? ftable : ""));
1038     PetscCall(PetscFree(ftable));
1039   }
1040   PetscCall(PetscFree(xnorms));
1041   PetscFunctionReturn(0);
1042 }
1043 
1044 int main(int argc, char **argv)
1045 {
1046   MPI_Comm       comm;
1047   TS             ts;
1048   DM             dm;
1049   Vec            u;
1050   AppCtx         user;
1051   PetscReal      t0, t = 0.0;
1052   void          *ctxs[2];
1053 
1054   ctxs[0] = &t;
1055   ctxs[1] = &t;
1056   PetscCall(PetscInitialize(&argc, &argv, (char*) 0, help));
1057   comm = PETSC_COMM_WORLD;
1058   user.functionalRegistry = NULL;
1059   globalUser = &user;
1060   PetscCall(ProcessOptions(comm, &user));
1061   PetscCall(TSCreate(comm, &ts));
1062   PetscCall(TSSetType(ts, TSBEULER));
1063   PetscCall(CreateDM(comm, &user, &dm));
1064   PetscCall(TSSetDM(ts, dm));
1065   PetscCall(ProcessMonitorOptions(comm, &user));
1066 
1067   PetscCall(DMCreateGlobalVector(dm, &u));
1068   PetscCall(PetscObjectSetName((PetscObject) u, "solution"));
1069   if (user.useFV) {
1070     PetscBool isImplicit = PETSC_FALSE;
1071 
1072     PetscCall(PetscOptionsHasName(NULL,"", "-use_implicit", &isImplicit));
1073     if (isImplicit) {
1074       PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user));
1075       PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user));
1076     }
1077     PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user));
1078     PetscCall(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, &user));
1079   } else {
1080     PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user));
1081     PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user));
1082     PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user));
1083   }
1084   if (user.useFV) PetscCall(TSMonitorSet(ts, MonitorFunctionals, &user, NULL));
1085   PetscCall(TSSetMaxSteps(ts, 1));
1086   PetscCall(TSSetMaxTime(ts, 2.0));
1087   PetscCall(TSSetTimeStep(ts,0.01));
1088   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
1089   PetscCall(TSSetFromOptions(ts));
1090 
1091   PetscCall(DMProjectFunction(dm, 0.0, user.initialGuess, ctxs, INSERT_VALUES, u));
1092   if (user.useFV) PetscCall(SetInitialConditionFVM(dm, u, 1, user.initialGuess[1], ctxs[1]));
1093   PetscCall(VecViewFromOptions(u, NULL, "-init_vec_view"));
1094   PetscCall(TSGetTime(ts, &t));
1095   t0   = t;
1096   PetscCall(DMTSCheckFromOptions(ts, u));
1097   PetscCall(TSSolve(ts, u));
1098   PetscCall(TSGetTime(ts, &t));
1099   if (t > t0) PetscCall(DMTSCheckFromOptions(ts, u));
1100   PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
1101   {
1102     PetscReal ftime;
1103     PetscInt  nsteps;
1104     TSConvergedReason reason;
1105 
1106     PetscCall(TSGetSolveTime(ts, &ftime));
1107     PetscCall(TSGetStepNumber(ts, &nsteps));
1108     PetscCall(TSGetConvergedReason(ts, &reason));
1109     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %" PetscInt_FMT " steps\n", TSConvergedReasons[reason], (double) ftime, nsteps));
1110   }
1111 
1112   PetscCall(VecDestroy(&u));
1113   PetscCall(DMDestroy(&dm));
1114   PetscCall(TSDestroy(&ts));
1115   PetscCall(PetscFree(user.monitorFuncs));
1116   PetscCall(FunctionalDestroy(&user.functionalRegistry));
1117   PetscCall(PetscFinalize());
1118   return 0;
1119 }
1120 
1121 /*TEST
1122 
1123   testset:
1124     args: -dm_plex_simplex 0 -dm_plex_box_faces 3,3,3
1125 
1126     # 2D harmonic velocity, no porosity
1127     test:
1128       suffix: p1p1
1129       requires: !complex !single
1130       args: -velocity_petscspace_degree 1 -porosity_petscspace_degree 1 -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -pc_factor_shift_type nonzero -ts_monitor -snes_error_if_not_converged -ksp_error_if_not_converged -dmts_check
1131     test:
1132       suffix: p1p1_xper
1133       requires: !complex !single
1134       args: -dm_refine 1 -dm_plex_box_bd periodic,none -velocity_petscspace_degree 1 -porosity_petscspace_degree 1 -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -pc_type lu -pc_factor_shift_type nonzero -ksp_rtol 1.0e-8 -ts_monitor -snes_error_if_not_converged -ksp_error_if_not_converged -dmts_check
1135     test:
1136       suffix: p1p1_xper_ref
1137       requires: !complex !single
1138       args: -dm_refine 2 -dm_plex_box_bd periodic,none -velocity_petscspace_degree 1 -porosity_petscspace_degree 1 -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -pc_type lu -pc_factor_shift_type nonzero -ksp_rtol 1.0e-8 -ts_monitor -snes_error_if_not_converged -ksp_error_if_not_converged -dmts_check
1139     test:
1140       suffix: p1p1_xyper
1141       requires: !complex !single
1142       args: -dm_refine 1 -dm_plex_box_bd periodic,periodic -velocity_petscspace_degree 1 -porosity_petscspace_degree 1 -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -pc_type lu -pc_factor_shift_type nonzero -ksp_rtol 1.0e-8 -ts_monitor -snes_error_if_not_converged -ksp_error_if_not_converged -dmts_check
1143     test:
1144       suffix: p1p1_xyper_ref
1145       requires: !complex !single
1146       args: -dm_refine 2 -dm_plex_box_bd periodic,periodic -velocity_petscspace_degree 1 -porosity_petscspace_degree 1 -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -pc_type lu -pc_factor_shift_type nonzero -ksp_rtol 1.0e-8 -ts_monitor -snes_error_if_not_converged -ksp_error_if_not_converged -dmts_check
1147     test:
1148       suffix: p2p1
1149       requires: !complex !single
1150       args: -velocity_petscspace_degree 2 -porosity_petscspace_degree 1 -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -ts_monitor   -snes_error_if_not_converged -ksp_error_if_not_converged -dmts_check
1151     test:
1152       suffix: p2p1_xyper
1153       requires: !complex !single
1154       args: -dm_refine 1 -dm_plex_box_bd periodic,periodic -velocity_petscspace_degree 2 -porosity_petscspace_degree 1 -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -pc_type lu -pc_factor_shift_type nonzero -ksp_rtol 1.0e-8 -ts_monitor -snes_error_if_not_converged -ksp_error_if_not_converged -dmts_check
1155 
1156     test:
1157       suffix: adv_1
1158       requires: !complex !single
1159       args: -use_fv -velocity_dist zero -porosity_dist tilted -ts_type ssp -ts_max_time 2.0 -ts_max_steps 1000 -ts_dt 0.993392 -bc_inflow 1,2,4 -bc_outflow 3 -ts_view -dm_view
1160 
1161     test:
1162       suffix: adv_2
1163       requires: !complex
1164       TODO: broken memory corruption
1165       args: -use_fv -velocity_dist zero -porosity_dist tilted -ts_type beuler -ts_max_time 2.0 -ts_max_steps 1000 -ts_dt 0.993392 -bc_inflow 3,4 -bc_outflow 1,2 -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -ksp_max_it 100 -ts_view -dm_view -snes_converged_reason -ksp_converged_reason
1166 
1167     test:
1168       suffix: adv_3
1169       requires: !complex
1170       TODO: broken memory corruption
1171       args: -dm_plex_box_bd periodic,none -use_fv -velocity_dist zero -porosity_dist tilted -ts_type beuler -ts_max_time 2.0 -ts_max_steps 1000 -ts_dt 0.993392 -bc_inflow 3 -bc_outflow 1 -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -ksp_max_it 100 -ts_view -dm_view -snes_converged_reason
1172 
1173     test:
1174       suffix: adv_3_ex
1175       requires: !complex
1176       args: -dm_plex_box_bd periodic,none -use_fv -velocity_dist zero -porosity_dist tilted -ts_type ssp -ts_max_time 2.0 -ts_max_steps 1000 -ts_dt 0.1 -bc_inflow 3 -bc_outflow 1 -snes_fd_color -ksp_max_it 100 -ts_view -dm_view
1177 
1178     test:
1179       suffix: adv_4
1180       requires: !complex
1181       TODO: broken memory corruption
1182       args: -use_fv -velocity_dist zero -porosity_dist tilted -ts_type beuler -ts_max_time 2.0 -ts_max_steps 1000 -ts_dt 0.993392 -bc_inflow 3 -bc_outflow 1 -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -ksp_max_it 100 -ts_view -dm_view -snes_converged_reason
1183 
1184     # 2D Advection, box, delta
1185     test:
1186       suffix: adv_delta_yper_0
1187       requires: !complex
1188       TODO: broken
1189       args: -dm_plex_box_bd none,periodic -use_fv -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type euler -ts_max_time 5.0 -ts_max_steps 20 -ts_dt 0.333333 -bc_inflow 2 -bc_outflow 4 -ts_view -dm_view -monitor Error
1190 
1191     test:
1192       suffix: adv_delta_yper_1
1193       requires: !complex
1194       TODO: broken
1195       args: -dm_plex_box_bd none,periodic -use_fv -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type euler -ts_max_time 5.0 -ts_max_steps 40 -ts_dt 0.166666 -bc_inflow 2 -bc_outflow 4 -ts_view -dm_view -monitor Error -dm_refine 1 -source_loc 0.416666,0.416666
1196 
1197     test:
1198       suffix: adv_delta_yper_2
1199       requires: !complex
1200       TODO: broken
1201       args: -dm_plex_box_bd none,periodic -use_fv -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type euler -ts_max_time 5.0 -ts_max_steps 80 -ts_dt 0.083333 -bc_inflow 2 -bc_outflow 4 -ts_view -dm_view -monitor Error -dm_refine 2 -source_loc 0.458333,0.458333
1202 
1203     test:
1204       suffix: adv_delta_yper_fim_0
1205       requires: !complex
1206       TODO: broken
1207       args: -dm_plex_box_bd none,periodic -use_fv -use_implicit -velocity_petscspace_degree 0 -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type mimex -ts_max_time 5.0 -ts_max_steps 20 -ts_dt 0.333333 -bc_inflow 2 -bc_outflow 4 -ts_view -monitor Error -dm_view   -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -mat_coloring_greedy_symmetric 0 -snes_rtol 1.0e-7 -pc_type lu -snes_converged_reason
1208 
1209     test:
1210       suffix: adv_delta_yper_fim_1
1211       requires: !complex
1212       TODO: broken
1213       args: -dm_plex_box_bd none,periodic -use_fv -use_implicit -velocity_petscspace_degree 1 -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type mimex -ts_max_time 5.0 -ts_max_steps 20 -ts_dt 0.333333 -bc_inflow 2 -bc_outflow 4 -ts_view -monitor Error -dm_view   -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -mat_coloring_greedy_symmetric 0 -snes_rtol 1.0e-7 -pc_type lu -snes_converged_reason -snes_linesearch_type basic
1214 
1215     test:
1216       suffix: adv_delta_yper_fim_2
1217       requires: !complex
1218       TODO: broken
1219       args: -dm_plex_box_bd none,periodic -use_fv -use_implicit -velocity_petscspace_degree 2 -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type mimex -ts_max_time 5.0 -ts_max_steps 20 -ts_dt 0.333333 -bc_inflow 2 -bc_outflow 4 -ts_view -monitor Error -dm_view   -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -mat_coloring_greedy_symmetric 0 -snes_rtol 1.0e-7 -pc_type lu -snes_converged_reason -snes_linesearch_type basic
1220 
1221     test:
1222       suffix: adv_delta_yper_im_0
1223       requires: !complex
1224       TODO: broken
1225       args: -dm_plex_box_bd none,periodic -use_fv -use_implicit -velocity_petscspace_degree 0 -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type mimex -ts_mimex_version 0 -ts_max_time 5.0 -ts_max_steps 20 -ts_dt 0.333333 -bc_inflow 2 -bc_outflow 4 -ts_view -monitor Error -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -snes_rtol 1.0e-7 -pc_type lu -snes_converged_reason
1226 
1227     test:
1228       suffix: adv_delta_yper_im_1
1229       requires: !complex
1230       TODO: broken
1231       args: -dm_plex_box_bd none,periodic -use_fv -use_implicit -velocity_petscspace_degree 0 -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type mimex -ts_mimex_version 0 -ts_max_time 5.0 -ts_max_steps 40 -ts_dt 0.166666 -bc_inflow 2 -bc_outflow 4 -ts_view -monitor Error -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -snes_rtol 1.0e-7 -pc_type lu -snes_converged_reason -dm_refine 1 -source_loc 0.416666,0.416666
1232 
1233     test:
1234       suffix: adv_delta_yper_im_2
1235       requires: !complex
1236       TODO: broken
1237       args: -dm_plex_box_bd none,periodic -use_fv -use_implicit -velocity_petscspace_degree 0 -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type mimex -ts_mimex_version 0 -ts_max_time 5.0 -ts_max_steps 80 -ts_dt 0.083333 -bc_inflow 2 -bc_outflow 4 -ts_view -monitor Error -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -snes_rtol 1.0e-7 -pc_type lu -snes_converged_reason -dm_refine 2 -source_loc 0.458333,0.458333
1238 
1239     test:
1240       suffix: adv_delta_yper_im_3
1241       requires: !complex
1242       TODO: broken
1243       args: -dm_plex_box_bd none,periodic -use_fv -use_implicit -velocity_petscspace_degree 1 -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type mimex -ts_mimex_version 0 -ts_max_time 5.0 -ts_max_steps 20 -ts_dt 0.333333 -bc_inflow 2 -bc_outflow 4 -ts_view -monitor Error -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -snes_rtol 1.0e-7 -pc_type lu -snes_converged_reason
1244 
1245     #    I believe the nullspace is sin(pi y)
1246     test:
1247       suffix: adv_delta_yper_im_4
1248       requires: !complex
1249       TODO: broken
1250       args: -dm_plex_box_bd none,periodic -use_fv -use_implicit -velocity_petscspace_degree 1 -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type mimex -ts_mimex_version 0 -ts_max_time 5.0 -ts_max_steps 40 -ts_dt 0.166666 -bc_inflow 2 -bc_outflow 4 -ts_view -monitor Error -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -snes_rtol 1.0e-7 -pc_type lu -snes_converged_reason -dm_refine 1 -source_loc 0.416666,0.416666
1251 
1252     test:
1253       suffix: adv_delta_yper_im_5
1254       requires: !complex
1255       TODO: broken
1256       args: -dm_plex_box_bd none,periodic -use_fv -use_implicit -velocity_petscspace_degree 1 -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type mimex -ts_mimex_version 0 -ts_max_time 5.0 -ts_max_steps 80 -ts_dt 0.083333 -bc_inflow 2 -bc_outflow 4 -ts_view -monitor Error -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -snes_rtol 1.0e-7 -pc_type lu -snes_converged_reason -dm_refine 2 -source_loc 0.458333,0.458333
1257 
1258     test:
1259       suffix: adv_delta_yper_im_6
1260       requires: !complex
1261       TODO: broken
1262       args: -dm_plex_box_bd none,periodic -use_fv -use_implicit -velocity_petscspace_degree 2 -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type mimex -ts_max_time 5.0 -ts_max_steps 20 -ts_dt 0.333333 -bc_inflow 2 -bc_outflow 4 -ts_view -monitor Error -dm_view   -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -snes_rtol 1.0e-7 -pc_type svd -snes_converged_reason
1263     # 2D Advection, magma benchmark 1
1264 
1265     test:
1266       suffix: adv_delta_shear_im_0
1267       requires: !complex
1268       TODO: broken
1269       args: -dm_plex_box_bd periodic,none -dm_refine 2 -use_fv -use_implicit -velocity_petscspace_degree 1 -velocity_dist shear -porosity_dist   delta -inflow_state 0.0 -ts_type mimex -ts_max_time 5.0 -ts_max_steps 20 -ts_dt 0.333333 -bc_inflow 1,3 -ts_view -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -snes_rtol 1.0e-7 -pc_type lu -snes_converged_reason -source_loc 0.458333,0.708333
1270     # 2D Advection, box, gaussian
1271 
1272     test:
1273       suffix: adv_gauss
1274       requires: !complex
1275       TODO: broken
1276       args: -use_fv -velocity_dist constant -porosity_dist gaussian -inflow_state 0.0 -ts_type ssp -ts_max_time 2.0 -ts_max_steps 100 -ts_dt 0.01 -bc_inflow 1 -bc_outflow 3 -ts_view -dm_view
1277 
1278     test:
1279       suffix: adv_gauss_im
1280       requires: !complex
1281       TODO: broken
1282       args: -use_fv -use_implicit -velocity_dist constant -porosity_dist gaussian -inflow_state 0.0 -ts_type beuler -ts_max_time 2.0 -ts_max_steps   100 -ts_dt 0.01 -bc_inflow 1 -bc_outflow 3 -ts_view -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -snes_rtol 1.0e-7
1283 
1284     test:
1285       suffix: adv_gauss_im_1
1286       requires: !complex
1287       TODO: broken
1288       args: -use_fv -use_implicit -velocity_petscspace_degree 1 -velocity_dist constant -porosity_dist gaussian -inflow_state 0.0 -ts_type beuler -ts_max_time 2.0 -ts_max_steps 100 -ts_dt 0.01 -bc_inflow 1 -bc_outflow 3 -ts_view -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -snes_rtol 1.0e-7
1289 
1290     test:
1291       suffix: adv_gauss_im_2
1292       requires: !complex
1293       TODO: broken
1294       args: -use_fv -use_implicit -velocity_petscspace_degree 2 -velocity_dist constant -porosity_dist gaussian -inflow_state 0.0 -ts_type beuler -ts_max_time 2.0 -ts_max_steps 100 -ts_dt 0.01 -bc_inflow 1 -bc_outflow 3 -ts_view -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -snes_rtol 1.0e-7
1295 
1296     test:
1297       suffix: adv_gauss_corner
1298       requires: !complex
1299       TODO: broken
1300       args: -use_fv -velocity_dist constant -porosity_dist gaussian -inflow_state 0.0 -ts_type ssp -ts_max_time 2.0 -ts_max_steps 100 -ts_dt 0.01 -bc_inflow 1 -bc_outflow 2 -ts_view -dm_view
1301 
1302     # 2D Advection+Harmonic 12-
1303     test:
1304       suffix: adv_harm_0
1305       requires: !complex
1306       TODO: broken memory corruption
1307       args: -velocity_petscspace_degree 2 -use_fv -velocity_dist harmonic -porosity_dist gaussian -ts_type beuler -ts_max_time 2.0 -ts_max_steps   1000 -ts_dt 0.993392 -bc_inflow 1,2,4 -bc_outflow 3 -use_implicit -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -ksp_max_it   100 -ts_view -dm_view -snes_converged_reason -ksp_converged_reason -snes_monitor -dmts_check
1308 
1309   #   Must check that FV BCs propagate to coarse meshes
1310   #   Must check that FV BC ids propagate to coarse meshes
1311   #   Must check that FE+FV BCs work at the same time
1312   # 2D Advection, matching wind in ex11 8-11
1313   #   NOTE implicit solves are limited by accuracy of FD Jacobian
1314   test:
1315     suffix: adv_0
1316     requires: !complex !single exodusii
1317     args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo -use_fv -velocity_dist zero -porosity_dist tilted -ts_type ssp -ts_max_time 2.0 -ts_max_steps 1000 -ts_dt 0.993392 -ts_view -dm_view
1318 
1319   test:
1320     suffix: adv_0_im
1321     requires: !complex exodusii
1322     TODO: broken  memory corruption
1323     args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo -use_fv -use_implicit -velocity_dist zero -porosity_dist tilted -ts_type beuler -ts_max_time 2.0 -ts_max_steps 1000 -ts_dt 0.993392 -ts_view -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -pc_type lu
1324 
1325   test:
1326     suffix: adv_0_im_2
1327     requires: !complex exodusii
1328     TODO: broken
1329     args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo -use_fv -use_implicit -velocity_dist constant -porosity_dist tilted -ts_type beuler -ts_max_time 2.0 -ts_max_steps 1000 -ts_dt 0.993392 -ts_view -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -pc_type lu -snes_rtol 1.0e-7
1330 
1331   test:
1332     suffix: adv_0_im_3
1333     requires: !complex exodusii
1334     TODO: broken
1335     args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo -use_fv -use_implicit -velocity_petscspace_degree 1 -velocity_dist constant -porosity_dist tilted -ts_type beuler -ts_max_time 2.0 -ts_max_steps 1000 -ts_dt 0.993392 -ts_view -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -pc_type svd -snes_rtol 1.0e-7
1336 
1337   test:
1338     suffix: adv_0_im_4
1339     requires: !complex exodusii
1340     TODO: broken
1341     args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo -use_fv -use_implicit -velocity_petscspace_degree 2 -velocity_dist constant -porosity_dist tilted -ts_type beuler -ts_max_time 2.0 -ts_max_steps 1000 -ts_dt 0.993392 -ts_view -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -pc_type svd -snes_rtol 1.0e-7
1342   # 2D Advection, misc
1343 
1344 TEST*/
1345