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