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