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