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