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