xref: /petsc/src/dm/impls/plex/tests/ex3.c (revision 2827ebad7ef18d4ed3f6d7633d71fb685bb29d6f)
1c4762a1bSJed Brown static char help[] = "Check that a DM can accurately represent and interpolate functions of a given polynomial order\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdmplex.h>
4c4762a1bSJed Brown #include <petscdm.h>
5c4762a1bSJed Brown #include <petscdmda.h>
6c4762a1bSJed Brown #include <petscfe.h>
7c4762a1bSJed Brown #include <petscds.h>
8c4762a1bSJed Brown #include <petscksp.h>
9c4762a1bSJed Brown #include <petscsnes.h>
10c4762a1bSJed Brown 
11c4762a1bSJed Brown typedef struct {
12c4762a1bSJed Brown   /* Domain and mesh definition */
13c4762a1bSJed Brown   PetscBool useDA;           /* Flag DMDA tensor product mesh */
14c4762a1bSJed Brown   PetscBool shearCoords;     /* Flag for shear transform */
15c4762a1bSJed Brown   PetscBool nonaffineCoords; /* Flag for non-affine transform */
16c4762a1bSJed Brown   /* Element definition */
17c4762a1bSJed Brown   PetscInt qorder;        /* Order of the quadrature */
18c4762a1bSJed Brown   PetscInt numComponents; /* Number of field components */
19c4762a1bSJed Brown   PetscFE  fe;            /* The finite element */
20c4762a1bSJed Brown   /* Testing space */
21c4762a1bSJed Brown   PetscInt  porder;         /* Order of polynomials to test */
22c4762a1bSJed Brown   PetscBool convergence;    /* Test for order of convergence */
23c4762a1bSJed Brown   PetscBool convRefine;     /* Test for convergence using refinement, otherwise use coarsening */
24c4762a1bSJed Brown   PetscBool constraints;    /* Test local constraints */
25c4762a1bSJed Brown   PetscBool tree;           /* Test tree routines */
26c4762a1bSJed Brown   PetscBool testFEjacobian; /* Test finite element Jacobian assembly */
27c4762a1bSJed Brown   PetscBool testFVgrad;     /* Test finite difference gradient routine */
28c4762a1bSJed Brown   PetscBool testInjector;   /* Test finite element injection routines */
29c4762a1bSJed Brown   PetscInt  treeCell;       /* Cell to refine in tree test */
30c4762a1bSJed Brown   PetscReal constants[3];   /* Constant values for each dimension */
31c4762a1bSJed Brown } AppCtx;
32c4762a1bSJed Brown 
33c4762a1bSJed Brown /* u = 1 */
34d71ae5a4SJacob Faibussowitsch PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx)
35d71ae5a4SJacob Faibussowitsch {
36c4762a1bSJed Brown   AppCtx  *user = (AppCtx *)ctx;
37c4762a1bSJed Brown   PetscInt d;
3830602db0SMatthew G. Knepley   for (d = 0; d < dim; ++d) u[d] = user->constants[d];
393ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
40c4762a1bSJed Brown }
41d71ae5a4SJacob Faibussowitsch PetscErrorCode constantDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx)
42d71ae5a4SJacob Faibussowitsch {
43c4762a1bSJed Brown   PetscInt d;
4430602db0SMatthew G. Knepley   for (d = 0; d < dim; ++d) u[d] = 0.0;
453ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
46c4762a1bSJed Brown }
47c4762a1bSJed Brown 
48c4762a1bSJed Brown /* u = x */
49d71ae5a4SJacob Faibussowitsch PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx)
50d71ae5a4SJacob Faibussowitsch {
51c4762a1bSJed Brown   PetscInt d;
52c4762a1bSJed Brown   for (d = 0; d < dim; ++d) u[d] = coords[d];
533ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
54c4762a1bSJed Brown }
55d71ae5a4SJacob Faibussowitsch PetscErrorCode linearDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx)
56d71ae5a4SJacob Faibussowitsch {
57c4762a1bSJed Brown   PetscInt d, e;
58c4762a1bSJed Brown   for (d = 0; d < dim; ++d) {
59c4762a1bSJed Brown     u[d] = 0.0;
60c4762a1bSJed Brown     for (e = 0; e < dim; ++e) u[d] += (d == e ? 1.0 : 0.0) * n[e];
61c4762a1bSJed Brown   }
623ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
63c4762a1bSJed Brown }
64c4762a1bSJed Brown 
65c4762a1bSJed Brown /* u = x^2 or u = (x^2, xy) or u = (xy, yz, zx) */
66d71ae5a4SJacob Faibussowitsch PetscErrorCode quadratic(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx)
67d71ae5a4SJacob Faibussowitsch {
689371c9d4SSatish Balay   if (dim > 2) {
699371c9d4SSatish Balay     u[0] = coords[0] * coords[1];
709371c9d4SSatish Balay     u[1] = coords[1] * coords[2];
719371c9d4SSatish Balay     u[2] = coords[2] * coords[0];
729371c9d4SSatish Balay   } else if (dim > 1) {
739371c9d4SSatish Balay     u[0] = coords[0] * coords[0];
749371c9d4SSatish Balay     u[1] = coords[0] * coords[1];
759371c9d4SSatish Balay   } else if (dim > 0) {
769371c9d4SSatish Balay     u[0] = coords[0] * coords[0];
779371c9d4SSatish Balay   }
783ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
79c4762a1bSJed Brown }
80d71ae5a4SJacob Faibussowitsch PetscErrorCode quadraticDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx)
81d71ae5a4SJacob Faibussowitsch {
829371c9d4SSatish Balay   if (dim > 2) {
839371c9d4SSatish Balay     u[0] = coords[1] * n[0] + coords[0] * n[1];
849371c9d4SSatish Balay     u[1] = coords[2] * n[1] + coords[1] * n[2];
859371c9d4SSatish Balay     u[2] = coords[2] * n[0] + coords[0] * n[2];
869371c9d4SSatish Balay   } else if (dim > 1) {
879371c9d4SSatish Balay     u[0] = 2.0 * coords[0] * n[0];
889371c9d4SSatish Balay     u[1] = coords[1] * n[0] + coords[0] * n[1];
899371c9d4SSatish Balay   } else if (dim > 0) {
909371c9d4SSatish Balay     u[0] = 2.0 * coords[0] * n[0];
919371c9d4SSatish Balay   }
923ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
93c4762a1bSJed Brown }
94c4762a1bSJed Brown 
95c4762a1bSJed Brown /* u = x^3 or u = (x^3, x^2y) or u = (x^2y, y^2z, z^2x) */
96d71ae5a4SJacob Faibussowitsch PetscErrorCode cubic(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx)
97d71ae5a4SJacob Faibussowitsch {
989371c9d4SSatish Balay   if (dim > 2) {
999371c9d4SSatish Balay     u[0] = coords[0] * coords[0] * coords[1];
1009371c9d4SSatish Balay     u[1] = coords[1] * coords[1] * coords[2];
1019371c9d4SSatish Balay     u[2] = coords[2] * coords[2] * coords[0];
1029371c9d4SSatish Balay   } else if (dim > 1) {
1039371c9d4SSatish Balay     u[0] = coords[0] * coords[0] * coords[0];
1049371c9d4SSatish Balay     u[1] = coords[0] * coords[0] * coords[1];
1059371c9d4SSatish Balay   } else if (dim > 0) {
1069371c9d4SSatish Balay     u[0] = coords[0] * coords[0] * coords[0];
1079371c9d4SSatish Balay   }
1083ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
109c4762a1bSJed Brown }
110d71ae5a4SJacob Faibussowitsch PetscErrorCode cubicDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx)
111d71ae5a4SJacob Faibussowitsch {
1129371c9d4SSatish Balay   if (dim > 2) {
1139371c9d4SSatish Balay     u[0] = 2.0 * coords[0] * coords[1] * n[0] + coords[0] * coords[0] * n[1];
1149371c9d4SSatish Balay     u[1] = 2.0 * coords[1] * coords[2] * n[1] + coords[1] * coords[1] * n[2];
1159371c9d4SSatish Balay     u[2] = 2.0 * coords[2] * coords[0] * n[2] + coords[2] * coords[2] * n[0];
1169371c9d4SSatish Balay   } else if (dim > 1) {
1179371c9d4SSatish Balay     u[0] = 3.0 * coords[0] * coords[0] * n[0];
1189371c9d4SSatish Balay     u[1] = 2.0 * coords[0] * coords[1] * n[0] + coords[0] * coords[0] * n[1];
1199371c9d4SSatish Balay   } else if (dim > 0) {
1209371c9d4SSatish Balay     u[0] = 3.0 * coords[0] * coords[0] * n[0];
1219371c9d4SSatish Balay   }
1223ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
123c4762a1bSJed Brown }
124c4762a1bSJed Brown 
125c4762a1bSJed Brown /* u = tanh(x) */
126d71ae5a4SJacob Faibussowitsch PetscErrorCode trig(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx)
127d71ae5a4SJacob Faibussowitsch {
128c4762a1bSJed Brown   PetscInt d;
12930602db0SMatthew G. Knepley   for (d = 0; d < dim; ++d) u[d] = PetscTanhReal(coords[d] - 0.5);
1303ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
131c4762a1bSJed Brown }
132d71ae5a4SJacob Faibussowitsch PetscErrorCode trigDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx)
133d71ae5a4SJacob Faibussowitsch {
134c4762a1bSJed Brown   PetscInt d;
13530602db0SMatthew G. Knepley   for (d = 0; d < dim; ++d) u[d] = 1.0 / PetscSqr(PetscCoshReal(coords[d] - 0.5)) * n[d];
1363ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
137c4762a1bSJed Brown }
138c4762a1bSJed Brown 
139d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
140d71ae5a4SJacob Faibussowitsch {
141c4762a1bSJed Brown   PetscInt n = 3;
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   PetscFunctionBeginUser;
14430602db0SMatthew G. Knepley   options->useDA           = PETSC_FALSE;
145c4762a1bSJed Brown   options->shearCoords     = PETSC_FALSE;
146c4762a1bSJed Brown   options->nonaffineCoords = PETSC_FALSE;
147c4762a1bSJed Brown   options->qorder          = 0;
148c4762a1bSJed Brown   options->numComponents   = PETSC_DEFAULT;
149c4762a1bSJed Brown   options->porder          = 0;
150c4762a1bSJed Brown   options->convergence     = PETSC_FALSE;
151c4762a1bSJed Brown   options->convRefine      = PETSC_TRUE;
152c4762a1bSJed Brown   options->constraints     = PETSC_FALSE;
153c4762a1bSJed Brown   options->tree            = PETSC_FALSE;
154c4762a1bSJed Brown   options->treeCell        = 0;
155c4762a1bSJed Brown   options->testFEjacobian  = PETSC_FALSE;
156c4762a1bSJed Brown   options->testFVgrad      = PETSC_FALSE;
157c4762a1bSJed Brown   options->testInjector    = PETSC_FALSE;
158c4762a1bSJed Brown   options->constants[0]    = 1.0;
159c4762a1bSJed Brown   options->constants[1]    = 2.0;
160c4762a1bSJed Brown   options->constants[2]    = 3.0;
161c4762a1bSJed Brown 
162d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Projection Test Options", "DMPlex");
1639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-use_da", "Flag for DMDA mesh", "ex3.c", options->useDA, &options->useDA, NULL));
1649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-shear_coords", "Transform coordinates with a shear", "ex3.c", options->shearCoords, &options->shearCoords, NULL));
1659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-non_affine_coords", "Transform coordinates with a non-affine transform", "ex3.c", options->nonaffineCoords, &options->nonaffineCoords, NULL));
1669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-qorder", "The quadrature order", "ex3.c", options->qorder, &options->qorder, NULL, 0));
1679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-num_comp", "The number of field components", "ex3.c", options->numComponents, &options->numComponents, NULL, PETSC_DEFAULT));
1689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-porder", "The order of polynomials to test", "ex3.c", options->porder, &options->porder, NULL, 0));
1699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-convergence", "Check the convergence rate", "ex3.c", options->convergence, &options->convergence, NULL));
1709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-conv_refine", "Use refinement for the convergence rate", "ex3.c", options->convRefine, &options->convRefine, NULL));
1719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-constraints", "Test local constraints (serial only)", "ex3.c", options->constraints, &options->constraints, NULL));
1729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-tree", "Test tree routines", "ex3.c", options->tree, &options->tree, NULL));
1739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-tree_cell", "cell to refine in tree test", "ex3.c", options->treeCell, &options->treeCell, NULL, 0));
1749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_fe_jacobian", "Test finite element Jacobian assembly", "ex3.c", options->testFEjacobian, &options->testFEjacobian, NULL));
1759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_fv_grad", "Test finite volume gradient reconstruction", "ex3.c", options->testFVgrad, &options->testFVgrad, NULL));
1769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_injector", "Test finite element injection", "ex3.c", options->testInjector, &options->testInjector, NULL));
1779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRealArray("-constants", "Set the constant values", "ex3.c", options->constants, &n, NULL));
178d0609cedSBarry Smith   PetscOptionsEnd();
179c4762a1bSJed Brown 
1803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
181c4762a1bSJed Brown }
182c4762a1bSJed Brown 
183d71ae5a4SJacob Faibussowitsch static PetscErrorCode TransformCoordinates(DM dm, AppCtx *user)
184d71ae5a4SJacob Faibussowitsch {
185c4762a1bSJed Brown   PetscSection coordSection;
186c4762a1bSJed Brown   Vec          coordinates;
187c4762a1bSJed Brown   PetscScalar *coords;
188c4762a1bSJed Brown   PetscInt     vStart, vEnd, v;
189c4762a1bSJed Brown 
190c4762a1bSJed Brown   PetscFunctionBeginUser;
191c4762a1bSJed Brown   if (user->nonaffineCoords) {
192c4762a1bSJed Brown     /* x' = r^(1/p) (x/r), y' = r^(1/p) (y/r), z' = z */
1939566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateSection(dm, &coordSection));
1949566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
1959566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1969566063dSJacob Faibussowitsch     PetscCall(VecGetArray(coordinates, &coords));
197c4762a1bSJed Brown     for (v = vStart; v < vEnd; ++v) {
198c4762a1bSJed Brown       PetscInt  dof, off;
199c4762a1bSJed Brown       PetscReal p = 4.0, r;
200c4762a1bSJed Brown 
2019566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(coordSection, v, &dof));
2029566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
203c4762a1bSJed Brown       switch (dof) {
204c4762a1bSJed Brown       case 2:
205c4762a1bSJed Brown         r               = PetscSqr(PetscRealPart(coords[off + 0])) + PetscSqr(PetscRealPart(coords[off + 1]));
206c4762a1bSJed Brown         coords[off + 0] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p) / (2 * p)) * coords[off + 0];
207c4762a1bSJed Brown         coords[off + 1] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p) / (2 * p)) * coords[off + 1];
208c4762a1bSJed Brown         break;
209c4762a1bSJed Brown       case 3:
210c4762a1bSJed Brown         r               = PetscSqr(PetscRealPart(coords[off + 0])) + PetscSqr(PetscRealPart(coords[off + 1]));
211c4762a1bSJed Brown         coords[off + 0] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p) / (2 * p)) * coords[off + 0];
212c4762a1bSJed Brown         coords[off + 1] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p) / (2 * p)) * coords[off + 1];
213c4762a1bSJed Brown         coords[off + 2] = coords[off + 2];
214c4762a1bSJed Brown         break;
215c4762a1bSJed Brown       }
216c4762a1bSJed Brown     }
2179566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(coordinates, &coords));
218c4762a1bSJed Brown   }
219c4762a1bSJed Brown   if (user->shearCoords) {
220c4762a1bSJed Brown     /* x' = x + m y + m z, y' = y + m z,  z' = z */
2219566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateSection(dm, &coordSection));
2229566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
2239566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
2249566063dSJacob Faibussowitsch     PetscCall(VecGetArray(coordinates, &coords));
225c4762a1bSJed Brown     for (v = vStart; v < vEnd; ++v) {
226c4762a1bSJed Brown       PetscInt  dof, off;
227c4762a1bSJed Brown       PetscReal m = 1.0;
228c4762a1bSJed Brown 
2299566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(coordSection, v, &dof));
2309566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
231c4762a1bSJed Brown       switch (dof) {
232c4762a1bSJed Brown       case 2:
233c4762a1bSJed Brown         coords[off + 0] = coords[off + 0] + m * coords[off + 1];
234c4762a1bSJed Brown         coords[off + 1] = coords[off + 1];
235c4762a1bSJed Brown         break;
236c4762a1bSJed Brown       case 3:
237c4762a1bSJed Brown         coords[off + 0] = coords[off + 0] + m * coords[off + 1] + m * coords[off + 2];
238c4762a1bSJed Brown         coords[off + 1] = coords[off + 1] + m * coords[off + 2];
239c4762a1bSJed Brown         coords[off + 2] = coords[off + 2];
240c4762a1bSJed Brown         break;
241c4762a1bSJed Brown       }
242c4762a1bSJed Brown     }
2439566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(coordinates, &coords));
244c4762a1bSJed Brown   }
2453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
246c4762a1bSJed Brown }
247c4762a1bSJed Brown 
248d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
249d71ae5a4SJacob Faibussowitsch {
25030602db0SMatthew G. Knepley   PetscInt  dim = 2;
25130602db0SMatthew G. Knepley   PetscBool simplex;
252c4762a1bSJed Brown 
253c4762a1bSJed Brown   PetscFunctionBeginUser;
25430602db0SMatthew G. Knepley   if (user->useDA) {
25530602db0SMatthew G. Knepley     switch (dim) {
256c4762a1bSJed Brown     case 2:
2579566063dSJacob Faibussowitsch       PetscCall(DMDACreate2d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 2, 2, PETSC_DETERMINE, PETSC_DETERMINE, 1, 1, NULL, NULL, dm));
2589566063dSJacob Faibussowitsch       PetscCall(DMSetFromOptions(*dm));
2599566063dSJacob Faibussowitsch       PetscCall(DMSetUp(*dm));
2609566063dSJacob Faibussowitsch       PetscCall(DMDASetVertexCoordinates(*dm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
261c4762a1bSJed Brown       break;
262d71ae5a4SJacob Faibussowitsch     default:
263d71ae5a4SJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot create structured mesh of dimension %" PetscInt_FMT, dim);
264c4762a1bSJed Brown     }
2659566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)*dm, "Hexahedral Mesh"));
26630602db0SMatthew G. Knepley   } else {
2679566063dSJacob Faibussowitsch     PetscCall(DMCreate(comm, dm));
2689566063dSJacob Faibussowitsch     PetscCall(DMSetType(*dm, DMPLEX));
2699566063dSJacob Faibussowitsch     PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
2709566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(*dm));
271c4762a1bSJed Brown 
2729566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(*dm, &dim));
2739566063dSJacob Faibussowitsch     PetscCall(DMPlexIsSimplex(*dm, &simplex));
2749566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&simplex, 1, MPIU_BOOL, 0, comm));
275c4762a1bSJed Brown     if (user->tree) {
27630602db0SMatthew G. Knepley       DM refTree, ncdm = NULL;
277c4762a1bSJed Brown 
2789566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateDefaultReferenceTree(comm, dim, simplex, &refTree));
2799566063dSJacob Faibussowitsch       PetscCall(DMViewFromOptions(refTree, NULL, "-reftree_dm_view"));
2809566063dSJacob Faibussowitsch       PetscCall(DMPlexSetReferenceTree(*dm, refTree));
2819566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&refTree));
2829566063dSJacob Faibussowitsch       PetscCall(DMPlexTreeRefineCell(*dm, user->treeCell, &ncdm));
283c4762a1bSJed Brown       if (ncdm) {
2849566063dSJacob Faibussowitsch         PetscCall(DMDestroy(dm));
285c4762a1bSJed Brown         *dm = ncdm;
2869566063dSJacob Faibussowitsch         PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_FALSE));
287c4762a1bSJed Brown       }
2889566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, "tree_"));
2899566063dSJacob Faibussowitsch       PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
2909566063dSJacob Faibussowitsch       PetscCall(DMSetFromOptions(*dm));
2919566063dSJacob Faibussowitsch       PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
292c4762a1bSJed Brown     } else {
2939566063dSJacob Faibussowitsch       PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_TRUE));
294c4762a1bSJed Brown     }
2959566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, "dist_"));
2969566063dSJacob Faibussowitsch     PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
2979566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(*dm));
2989566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, NULL));
2999566063dSJacob Faibussowitsch     if (simplex) PetscCall(PetscObjectSetName((PetscObject)*dm, "Simplicial Mesh"));
3009566063dSJacob Faibussowitsch     else PetscCall(PetscObjectSetName((PetscObject)*dm, "Hexahedral Mesh"));
301c4762a1bSJed Brown   }
3029566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
3039566063dSJacob Faibussowitsch   PetscCall(TransformCoordinates(*dm, user));
3049566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
3053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
306c4762a1bSJed Brown }
307c4762a1bSJed Brown 
308d71ae5a4SJacob Faibussowitsch static void simple_mass(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
309d71ae5a4SJacob Faibussowitsch {
310c4762a1bSJed Brown   PetscInt d, e;
311ad540459SPierre Jolivet   for (d = 0, e = 0; d < dim; d++, e += dim + 1) g0[e] = 1.;
312c4762a1bSJed Brown }
313c4762a1bSJed Brown 
314c4762a1bSJed Brown /* < \nabla v, 1/2(\nabla u + {\nabla u}^T) > */
315d71ae5a4SJacob Faibussowitsch static void symmetric_gradient_inner_product(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar C[])
316d71ae5a4SJacob Faibussowitsch {
317c4762a1bSJed Brown   PetscInt compI, compJ, d, e;
318c4762a1bSJed Brown 
319c4762a1bSJed Brown   for (compI = 0; compI < dim; ++compI) {
320c4762a1bSJed Brown     for (compJ = 0; compJ < dim; ++compJ) {
321c4762a1bSJed Brown       for (d = 0; d < dim; ++d) {
322c4762a1bSJed Brown         for (e = 0; e < dim; e++) {
323c4762a1bSJed Brown           if (d == e && d == compI && d == compJ) {
324c4762a1bSJed Brown             C[((compI * dim + compJ) * dim + d) * dim + e] = 1.0;
325c4762a1bSJed Brown           } else if ((d == compJ && e == compI) || (d == e && compI == compJ)) {
326c4762a1bSJed Brown             C[((compI * dim + compJ) * dim + d) * dim + e] = 0.5;
327c4762a1bSJed Brown           } else {
328c4762a1bSJed Brown             C[((compI * dim + compJ) * dim + d) * dim + e] = 0.0;
329c4762a1bSJed Brown           }
330c4762a1bSJed Brown         }
331c4762a1bSJed Brown       }
332c4762a1bSJed Brown     }
333c4762a1bSJed Brown   }
334c4762a1bSJed Brown }
335c4762a1bSJed Brown 
336d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupSection(DM dm, AppCtx *user)
337d71ae5a4SJacob Faibussowitsch {
338c4762a1bSJed Brown   PetscFunctionBeginUser;
33930602db0SMatthew G. Knepley   if (user->constraints) {
340c4762a1bSJed Brown     /* test local constraints */
341c4762a1bSJed Brown     DM           coordDM;
342c4762a1bSJed Brown     PetscInt     fStart, fEnd, f, vStart, vEnd, v;
343c4762a1bSJed Brown     PetscInt     edgesx = 2, vertsx;
344c4762a1bSJed Brown     PetscInt     edgesy = 2, vertsy;
345c4762a1bSJed Brown     PetscMPIInt  size;
346c4762a1bSJed Brown     PetscInt     numConst;
347c4762a1bSJed Brown     PetscSection aSec;
348c4762a1bSJed Brown     PetscInt    *anchors;
349c4762a1bSJed Brown     PetscInt     offset;
350c4762a1bSJed Brown     IS           aIS;
351c4762a1bSJed Brown     MPI_Comm     comm = PetscObjectComm((PetscObject)dm);
352c4762a1bSJed Brown 
3539566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(comm, &size));
35408401ef6SPierre Jolivet     PetscCheck(size <= 1, comm, PETSC_ERR_SUP, "Local constraint test can only be performed in serial");
355c4762a1bSJed Brown 
356c4762a1bSJed Brown     /* we are going to test constraints by using them to enforce periodicity
357c4762a1bSJed Brown      * in one direction, and comparing to the existing method of enforcing
358c4762a1bSJed Brown      * periodicity */
359c4762a1bSJed Brown 
360c4762a1bSJed Brown     /* first create the coordinate section so that it does not clone the
361c4762a1bSJed Brown      * constraints */
3629566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dm, &coordDM));
363c4762a1bSJed Brown 
364c4762a1bSJed Brown     /* create the constrained-to-anchor section */
3659566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
3669566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 1, &fStart, &fEnd));
3679566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreate(PETSC_COMM_SELF, &aSec));
3689566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetChart(aSec, PetscMin(fStart, vStart), PetscMax(fEnd, vEnd)));
369c4762a1bSJed Brown 
370c4762a1bSJed Brown     /* define the constraints */
3719566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetInt(NULL, NULL, "-da_grid_x", &edgesx, NULL));
3729566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetInt(NULL, NULL, "-da_grid_y", &edgesy, NULL));
373c4762a1bSJed Brown     vertsx   = edgesx + 1;
374c4762a1bSJed Brown     vertsy   = edgesy + 1;
375c4762a1bSJed Brown     numConst = vertsy + edgesy;
3769566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numConst, &anchors));
377c4762a1bSJed Brown     offset = 0;
378c4762a1bSJed Brown     for (v = vStart + edgesx; v < vEnd; v += vertsx) {
3799566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetDof(aSec, v, 1));
380c4762a1bSJed Brown       anchors[offset++] = v - edgesx;
381c4762a1bSJed Brown     }
382c4762a1bSJed Brown     for (f = fStart + edgesx * vertsy + edgesx * edgesy; f < fEnd; f++) {
3839566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetDof(aSec, f, 1));
384c4762a1bSJed Brown       anchors[offset++] = f - edgesx * edgesy;
385c4762a1bSJed Brown     }
3869566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetUp(aSec));
3879566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numConst, anchors, PETSC_OWN_POINTER, &aIS));
388c4762a1bSJed Brown 
3899566063dSJacob Faibussowitsch     PetscCall(DMPlexSetAnchors(dm, aSec, aIS));
3909566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&aSec));
3919566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&aIS));
392c4762a1bSJed Brown   }
3939566063dSJacob Faibussowitsch   PetscCall(DMSetNumFields(dm, 1));
3949566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)user->fe));
3959566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
39630602db0SMatthew G. Knepley   if (user->constraints) {
397c4762a1bSJed Brown     /* test getting local constraint matrix that matches section */
398c4762a1bSJed Brown     PetscSection aSec;
399c4762a1bSJed Brown     IS           aIS;
400c4762a1bSJed Brown 
4019566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAnchors(dm, &aSec, &aIS));
402c4762a1bSJed Brown     if (aSec) {
403c4762a1bSJed Brown       PetscDS         ds;
404c4762a1bSJed Brown       PetscSection    cSec, section;
405c4762a1bSJed Brown       PetscInt        cStart, cEnd, c, numComp;
406c4762a1bSJed Brown       Mat             cMat, mass;
407c4762a1bSJed Brown       Vec             local;
408c4762a1bSJed Brown       const PetscInt *anchors;
409c4762a1bSJed Brown 
4109566063dSJacob Faibussowitsch       PetscCall(DMGetLocalSection(dm, &section));
411c4762a1bSJed Brown       /* this creates the matrix and preallocates the matrix structure: we
412c4762a1bSJed Brown        * just have to fill in the values */
4139566063dSJacob Faibussowitsch       PetscCall(DMGetDefaultConstraints(dm, &cSec, &cMat, NULL));
4149566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetChart(cSec, &cStart, &cEnd));
4159566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(aIS, &anchors));
4169566063dSJacob Faibussowitsch       PetscCall(PetscFEGetNumComponents(user->fe, &numComp));
417c4762a1bSJed Brown       for (c = cStart; c < cEnd; c++) {
418c4762a1bSJed Brown         PetscInt cDof;
419c4762a1bSJed Brown 
420c4762a1bSJed Brown         /* is this point constrained? (does it have an anchor?) */
4219566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(aSec, c, &cDof));
422c4762a1bSJed Brown         if (cDof) {
423c4762a1bSJed Brown           PetscInt cOff, a, aDof, aOff, j;
42463a3b9bcSJacob Faibussowitsch           PetscCheck(cDof == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Found %" PetscInt_FMT " anchor points: should be just one", cDof);
425c4762a1bSJed Brown 
426c4762a1bSJed Brown           /* find the anchor point */
4279566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetOffset(aSec, c, &cOff));
428c4762a1bSJed Brown           a = anchors[cOff];
429c4762a1bSJed Brown 
430c4762a1bSJed Brown           /* find the constrained dofs (row in constraint matrix) */
4319566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetDof(cSec, c, &cDof));
4329566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetOffset(cSec, c, &cOff));
433c4762a1bSJed Brown 
434c4762a1bSJed Brown           /* find the anchor dofs (column in constraint matrix) */
4359566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetDof(section, a, &aDof));
4369566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetOffset(section, a, &aOff));
437c4762a1bSJed Brown 
43863a3b9bcSJacob Faibussowitsch           PetscCheck(cDof == aDof, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point and anchor have different number of dofs: %" PetscInt_FMT ", %" PetscInt_FMT, cDof, aDof);
43963a3b9bcSJacob Faibussowitsch           PetscCheck(cDof % numComp == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point dofs not divisible by field components: %" PetscInt_FMT ", %" PetscInt_FMT, cDof, numComp);
440c4762a1bSJed Brown 
441c4762a1bSJed Brown           /* put in a simple equality constraint */
44248a46eb9SPierre Jolivet           for (j = 0; j < cDof; j++) PetscCall(MatSetValue(cMat, cOff + j, aOff + j, 1., INSERT_VALUES));
443c4762a1bSJed Brown         }
444c4762a1bSJed Brown       }
4459566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(cMat, MAT_FINAL_ASSEMBLY));
4469566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(cMat, MAT_FINAL_ASSEMBLY));
4479566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(aIS, &anchors));
448c4762a1bSJed Brown 
449c4762a1bSJed Brown       /* Now that we have constructed the constraint matrix, any FE matrix
450c4762a1bSJed Brown        * that we construct will apply the constraints during construction */
451c4762a1bSJed Brown 
4529566063dSJacob Faibussowitsch       PetscCall(DMCreateMatrix(dm, &mass));
453c4762a1bSJed Brown       /* get a dummy local variable to serve as the solution */
4549566063dSJacob Faibussowitsch       PetscCall(DMGetLocalVector(dm, &local));
4559566063dSJacob Faibussowitsch       PetscCall(DMGetDS(dm, &ds));
456c4762a1bSJed Brown       /* set the jacobian to be the mass matrix */
4579566063dSJacob Faibussowitsch       PetscCall(PetscDSSetJacobian(ds, 0, 0, simple_mass, NULL, NULL, NULL));
458c4762a1bSJed Brown       /* build the mass matrix */
4599566063dSJacob Faibussowitsch       PetscCall(DMPlexSNESComputeJacobianFEM(dm, local, mass, mass, NULL));
4609566063dSJacob Faibussowitsch       PetscCall(MatView(mass, PETSC_VIEWER_STDOUT_WORLD));
4619566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mass));
4629566063dSJacob Faibussowitsch       PetscCall(DMRestoreLocalVector(dm, &local));
463c4762a1bSJed Brown     }
464c4762a1bSJed Brown   }
4653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
466c4762a1bSJed Brown }
467c4762a1bSJed Brown 
468d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestFEJacobian(DM dm, AppCtx *user)
469d71ae5a4SJacob Faibussowitsch {
470c4762a1bSJed Brown   PetscFunctionBeginUser;
47130602db0SMatthew G. Knepley   if (!user->useDA) {
472c4762a1bSJed Brown     Vec          local;
473c4762a1bSJed Brown     const Vec   *vecs;
474c4762a1bSJed Brown     Mat          E;
475c4762a1bSJed Brown     MatNullSpace sp;
476c4762a1bSJed Brown     PetscBool    isNullSpace, hasConst;
47730602db0SMatthew G. Knepley     PetscInt     dim, n, i;
478c4762a1bSJed Brown     Vec          res = NULL, localX, localRes;
479c4762a1bSJed Brown     PetscDS      ds;
480c4762a1bSJed Brown 
4819566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
48263a3b9bcSJacob Faibussowitsch     PetscCheck(user->numComponents == dim, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "The number of components %" PetscInt_FMT " must be equal to the dimension %" PetscInt_FMT " for this test", user->numComponents, dim);
4839566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dm, &ds));
4849566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, symmetric_gradient_inner_product));
4859566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix(dm, &E));
4869566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &local));
4879566063dSJacob Faibussowitsch     PetscCall(DMPlexSNESComputeJacobianFEM(dm, local, E, E, NULL));
4889566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateRigidBody(dm, 0, &sp));
4899566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceGetVecs(sp, &hasConst, &n, &vecs));
4909566063dSJacob Faibussowitsch     if (n) PetscCall(VecDuplicate(vecs[0], &res));
4919566063dSJacob Faibussowitsch     PetscCall(DMCreateLocalVector(dm, &localX));
4929566063dSJacob Faibussowitsch     PetscCall(DMCreateLocalVector(dm, &localRes));
493c4762a1bSJed Brown     for (i = 0; i < n; i++) { /* also test via matrix-free Jacobian application */
494c4762a1bSJed Brown       PetscReal resNorm;
495c4762a1bSJed Brown 
4969566063dSJacob Faibussowitsch       PetscCall(VecSet(localRes, 0.));
4979566063dSJacob Faibussowitsch       PetscCall(VecSet(localX, 0.));
4989566063dSJacob Faibussowitsch       PetscCall(VecSet(local, 0.));
4999566063dSJacob Faibussowitsch       PetscCall(VecSet(res, 0.));
5009566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalBegin(dm, vecs[i], INSERT_VALUES, localX));
5019566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalEnd(dm, vecs[i], INSERT_VALUES, localX));
5029566063dSJacob Faibussowitsch       PetscCall(DMSNESComputeJacobianAction(dm, local, localX, localRes, NULL));
5039566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalBegin(dm, localRes, ADD_VALUES, res));
5049566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalEnd(dm, localRes, ADD_VALUES, res));
5059566063dSJacob Faibussowitsch       PetscCall(VecNorm(res, NORM_2, &resNorm));
50648a46eb9SPierre Jolivet       if (resNorm > PETSC_SMALL) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Symmetric gradient action null space vector %" PetscInt_FMT " residual: %E\n", i, (double)resNorm));
507c4762a1bSJed Brown     }
5089566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&localRes));
5099566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&localX));
5109566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&res));
5119566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceTest(sp, E, &isNullSpace));
512c4762a1bSJed Brown     if (isNullSpace) {
5139566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Symmetric gradient null space: PASS\n"));
514c4762a1bSJed Brown     } else {
5159566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Symmetric gradient null space: FAIL\n"));
516c4762a1bSJed Brown     }
5179566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceDestroy(&sp));
5189566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&E));
5199566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &local));
520c4762a1bSJed Brown   }
5213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
522c4762a1bSJed Brown }
523c4762a1bSJed Brown 
524d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestInjector(DM dm, AppCtx *user)
525d71ae5a4SJacob Faibussowitsch {
526c4762a1bSJed Brown   DM          refTree;
527c4762a1bSJed Brown   PetscMPIInt rank;
528c4762a1bSJed Brown 
529c4762a1bSJed Brown   PetscFunctionBegin;
5309566063dSJacob Faibussowitsch   PetscCall(DMPlexGetReferenceTree(dm, &refTree));
531c4762a1bSJed Brown   if (refTree) {
532c4762a1bSJed Brown     Mat inj;
533c4762a1bSJed Brown 
5349566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeInjectorReferenceTree(refTree, &inj));
5359566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)inj, "Reference Tree Injector"));
5369566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
53748a46eb9SPierre Jolivet     if (rank == 0) PetscCall(MatView(inj, PETSC_VIEWER_STDOUT_SELF));
5389566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&inj));
539c4762a1bSJed Brown   }
5403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
541c4762a1bSJed Brown }
542c4762a1bSJed Brown 
543d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestFVGrad(DM dm, AppCtx *user)
544d71ae5a4SJacob Faibussowitsch {
545c4762a1bSJed Brown   MPI_Comm           comm;
546c4762a1bSJed Brown   DM                 dmRedist, dmfv, dmgrad, dmCell, refTree;
547c4762a1bSJed Brown   PetscFV            fv;
54830602db0SMatthew G. Knepley   PetscInt           dim, nvecs, v, cStart, cEnd, cEndInterior;
549c4762a1bSJed Brown   PetscMPIInt        size;
550c4762a1bSJed Brown   Vec                cellgeom, grad, locGrad;
551c4762a1bSJed Brown   const PetscScalar *cgeom;
552c4762a1bSJed Brown   PetscReal          allVecMaxDiff = 0., fvTol = 100. * PETSC_MACHINE_EPSILON;
553c4762a1bSJed Brown 
554c4762a1bSJed Brown   PetscFunctionBeginUser;
555c4762a1bSJed Brown   comm = PetscObjectComm((PetscObject)dm);
5569566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
557c4762a1bSJed Brown   /* duplicate DM, give dup. a FV discretization */
5589566063dSJacob Faibussowitsch   PetscCall(DMSetBasicAdjacency(dm, PETSC_TRUE, PETSC_FALSE));
5599566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
560c4762a1bSJed Brown   dmRedist = NULL;
56148a46eb9SPierre Jolivet   if (size > 1) PetscCall(DMPlexDistributeOverlap(dm, 1, NULL, &dmRedist));
562c4762a1bSJed Brown   if (!dmRedist) {
563c4762a1bSJed Brown     dmRedist = dm;
5649566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)dmRedist));
565c4762a1bSJed Brown   }
5669566063dSJacob Faibussowitsch   PetscCall(PetscFVCreate(comm, &fv));
5679566063dSJacob Faibussowitsch   PetscCall(PetscFVSetType(fv, PETSCFVLEASTSQUARES));
5689566063dSJacob Faibussowitsch   PetscCall(PetscFVSetNumComponents(fv, user->numComponents));
5699566063dSJacob Faibussowitsch   PetscCall(PetscFVSetSpatialDimension(fv, dim));
5709566063dSJacob Faibussowitsch   PetscCall(PetscFVSetFromOptions(fv));
5719566063dSJacob Faibussowitsch   PetscCall(PetscFVSetUp(fv));
5729566063dSJacob Faibussowitsch   PetscCall(DMPlexConstructGhostCells(dmRedist, NULL, NULL, &dmfv));
5739566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmRedist));
5749566063dSJacob Faibussowitsch   PetscCall(DMSetNumFields(dmfv, 1));
5759566063dSJacob Faibussowitsch   PetscCall(DMSetField(dmfv, 0, NULL, (PetscObject)fv));
5769566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dmfv));
5779566063dSJacob Faibussowitsch   PetscCall(DMPlexGetReferenceTree(dm, &refTree));
5789566063dSJacob Faibussowitsch   if (refTree) PetscCall(DMCopyDisc(dmfv, refTree));
5799566063dSJacob Faibussowitsch   PetscCall(DMPlexGetGradientDM(dmfv, fv, &dmgrad));
5809566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dmfv, 0, &cStart, &cEnd));
58130602db0SMatthew G. Knepley   nvecs = dim * (dim + 1) / 2;
5829566063dSJacob Faibussowitsch   PetscCall(DMPlexGetGeometryFVM(dmfv, NULL, &cellgeom, NULL));
5839566063dSJacob Faibussowitsch   PetscCall(VecGetDM(cellgeom, &dmCell));
5849566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(cellgeom, &cgeom));
5859566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dmgrad, &grad));
5869566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dmgrad, &locGrad));
587*2827ebadSStefano Zampini   PetscCall(DMPlexGetCellTypeStratum(dmgrad, DM_POLYTOPE_FV_GHOST, &cEndInterior, NULL));
588c4762a1bSJed Brown   cEndInterior = (cEndInterior < 0) ? cEnd : cEndInterior;
589c4762a1bSJed Brown   for (v = 0; v < nvecs; v++) {
590c4762a1bSJed Brown     Vec                locX;
591c4762a1bSJed Brown     PetscInt           c;
592c4762a1bSJed Brown     PetscScalar        trueGrad[3][3] = {{0.}};
593c4762a1bSJed Brown     const PetscScalar *gradArray;
594c4762a1bSJed Brown     PetscReal          maxDiff, maxDiffGlob;
595c4762a1bSJed Brown 
5969566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dmfv, &locX));
597c4762a1bSJed Brown     /* get the local projection of the rigid body mode */
598c4762a1bSJed Brown     for (c = cStart; c < cEnd; c++) {
599c4762a1bSJed Brown       PetscFVCellGeom *cg;
600c4762a1bSJed Brown       PetscScalar      cx[3] = {0., 0., 0.};
601c4762a1bSJed Brown 
6029566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg));
60330602db0SMatthew G. Knepley       if (v < dim) {
604c4762a1bSJed Brown         cx[v] = 1.;
605c4762a1bSJed Brown       } else {
60630602db0SMatthew G. Knepley         PetscInt w = v - dim;
607c4762a1bSJed Brown 
60830602db0SMatthew G. Knepley         cx[(w + 1) % dim] = cg->centroid[(w + 2) % dim];
60930602db0SMatthew G. Knepley         cx[(w + 2) % dim] = -cg->centroid[(w + 1) % dim];
610c4762a1bSJed Brown       }
6119566063dSJacob Faibussowitsch       PetscCall(DMPlexVecSetClosure(dmfv, NULL, locX, c, cx, INSERT_ALL_VALUES));
612c4762a1bSJed Brown     }
613c4762a1bSJed Brown     /* TODO: this isn't in any header */
6149566063dSJacob Faibussowitsch     PetscCall(DMPlexReconstructGradientsFVM(dmfv, locX, grad));
6159566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(dmgrad, grad, INSERT_VALUES, locGrad));
6169566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(dmgrad, grad, INSERT_VALUES, locGrad));
6179566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(locGrad, &gradArray));
618c4762a1bSJed Brown     /* compare computed gradient to exact gradient */
61930602db0SMatthew G. Knepley     if (v >= dim) {
62030602db0SMatthew G. Knepley       PetscInt w = v - dim;
621c4762a1bSJed Brown 
62230602db0SMatthew G. Knepley       trueGrad[(w + 1) % dim][(w + 2) % dim] = 1.;
62330602db0SMatthew G. Knepley       trueGrad[(w + 2) % dim][(w + 1) % dim] = -1.;
624c4762a1bSJed Brown     }
625c4762a1bSJed Brown     maxDiff = 0.;
626c4762a1bSJed Brown     for (c = cStart; c < cEndInterior; c++) {
627c4762a1bSJed Brown       PetscScalar *compGrad;
628c4762a1bSJed Brown       PetscInt     i, j, k;
629c4762a1bSJed Brown       PetscReal    FrobDiff = 0.;
630c4762a1bSJed Brown 
6319566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRead(dmgrad, c, gradArray, &compGrad));
632c4762a1bSJed Brown 
63330602db0SMatthew G. Knepley       for (i = 0, k = 0; i < dim; i++) {
63430602db0SMatthew G. Knepley         for (j = 0; j < dim; j++, k++) {
635c4762a1bSJed Brown           PetscScalar diff = compGrad[k] - trueGrad[i][j];
636c4762a1bSJed Brown           FrobDiff += PetscRealPart(diff * PetscConj(diff));
637c4762a1bSJed Brown         }
638c4762a1bSJed Brown       }
639c4762a1bSJed Brown       FrobDiff = PetscSqrtReal(FrobDiff);
640c4762a1bSJed Brown       maxDiff  = PetscMax(maxDiff, FrobDiff);
641c4762a1bSJed Brown     }
642712fec58SPierre Jolivet     PetscCall(MPIU_Allreduce(&maxDiff, &maxDiffGlob, 1, MPIU_REAL, MPIU_MAX, comm));
643c4762a1bSJed Brown     allVecMaxDiff = PetscMax(allVecMaxDiff, maxDiffGlob);
6449566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(locGrad, &gradArray));
6459566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dmfv, &locX));
646c4762a1bSJed Brown   }
647c4762a1bSJed Brown   if (allVecMaxDiff < fvTol) {
6489566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Finite volume gradient reconstruction: PASS\n"));
649c4762a1bSJed Brown   } else {
65063a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Finite volume gradient reconstruction: FAIL at tolerance %g with max difference %g\n", (double)fvTol, (double)allVecMaxDiff));
651c4762a1bSJed Brown   }
6529566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dmgrad, &locGrad));
6539566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dmgrad, &grad));
6549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(cellgeom, &cgeom));
6559566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmfv));
6569566063dSJacob Faibussowitsch   PetscCall(PetscFVDestroy(&fv));
6573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
658c4762a1bSJed Brown }
659c4762a1bSJed Brown 
660d71ae5a4SJacob Faibussowitsch static PetscErrorCode ComputeError(DM dm, PetscErrorCode (**exactFuncs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), PetscErrorCode (**exactFuncDers)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), void **exactCtxs, PetscReal *error, PetscReal *errorDer, AppCtx *user)
661d71ae5a4SJacob Faibussowitsch {
662c4762a1bSJed Brown   Vec       u;
663c4762a1bSJed Brown   PetscReal n[3] = {1.0, 1.0, 1.0};
664c4762a1bSJed Brown 
665c4762a1bSJed Brown   PetscFunctionBeginUser;
6669566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &u));
667c4762a1bSJed Brown   /* Project function into FE function space */
6689566063dSJacob Faibussowitsch   PetscCall(DMProjectFunction(dm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, u));
6699566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-projection_view"));
670c4762a1bSJed Brown   /* Compare approximation to exact in L_2 */
6719566063dSJacob Faibussowitsch   PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, exactCtxs, u, error));
6729566063dSJacob Faibussowitsch   PetscCall(DMComputeL2GradientDiff(dm, 0.0, exactFuncDers, exactCtxs, u, n, errorDer));
6739566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &u));
6743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
675c4762a1bSJed Brown }
676c4762a1bSJed Brown 
677d71ae5a4SJacob Faibussowitsch static PetscErrorCode CheckFunctions(DM dm, PetscInt order, AppCtx *user)
678d71ae5a4SJacob Faibussowitsch {
679c4762a1bSJed Brown   PetscErrorCode (*exactFuncs[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
680c4762a1bSJed Brown   PetscErrorCode (*exactFuncDers[1])(PetscInt dim, PetscReal time, const PetscReal x[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx);
681c4762a1bSJed Brown   void     *exactCtxs[3];
682c4762a1bSJed Brown   MPI_Comm  comm;
683c4762a1bSJed Brown   PetscReal error, errorDer, tol = PETSC_SMALL;
684c4762a1bSJed Brown 
685c4762a1bSJed Brown   PetscFunctionBeginUser;
686c4762a1bSJed Brown   exactCtxs[0] = user;
687c4762a1bSJed Brown   exactCtxs[1] = user;
688c4762a1bSJed Brown   exactCtxs[2] = user;
6899566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
690c4762a1bSJed Brown   /* Setup functions to approximate */
691c4762a1bSJed Brown   switch (order) {
692c4762a1bSJed Brown   case 0:
693c4762a1bSJed Brown     exactFuncs[0]    = constant;
694c4762a1bSJed Brown     exactFuncDers[0] = constantDer;
695c4762a1bSJed Brown     break;
696c4762a1bSJed Brown   case 1:
697c4762a1bSJed Brown     exactFuncs[0]    = linear;
698c4762a1bSJed Brown     exactFuncDers[0] = linearDer;
699c4762a1bSJed Brown     break;
700c4762a1bSJed Brown   case 2:
701c4762a1bSJed Brown     exactFuncs[0]    = quadratic;
702c4762a1bSJed Brown     exactFuncDers[0] = quadraticDer;
703c4762a1bSJed Brown     break;
704c4762a1bSJed Brown   case 3:
705c4762a1bSJed Brown     exactFuncs[0]    = cubic;
706c4762a1bSJed Brown     exactFuncDers[0] = cubicDer;
707c4762a1bSJed Brown     break;
708d71ae5a4SJacob Faibussowitsch   default:
709d71ae5a4SJacob Faibussowitsch     SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Could not determine functions to test for order %" PetscInt_FMT, order);
710c4762a1bSJed Brown   }
7119566063dSJacob Faibussowitsch   PetscCall(ComputeError(dm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user));
712c4762a1bSJed Brown   /* Report result */
71363a3b9bcSJacob Faibussowitsch   if (error > tol) PetscCall(PetscPrintf(comm, "Function tests FAIL for order %" PetscInt_FMT " at tolerance %g error %g\n", order, (double)tol, (double)error));
71463a3b9bcSJacob Faibussowitsch   else PetscCall(PetscPrintf(comm, "Function tests pass for order %" PetscInt_FMT " at tolerance %g\n", order, (double)tol));
71563a3b9bcSJacob Faibussowitsch   if (errorDer > tol) PetscCall(PetscPrintf(comm, "Function tests FAIL for order %" PetscInt_FMT " derivatives at tolerance %g error %g\n", order, (double)tol, (double)errorDer));
71663a3b9bcSJacob Faibussowitsch   else PetscCall(PetscPrintf(comm, "Function tests pass for order %" PetscInt_FMT " derivatives at tolerance %g\n", order, (double)tol));
7173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
718c4762a1bSJed Brown }
719c4762a1bSJed Brown 
720d71ae5a4SJacob Faibussowitsch static PetscErrorCode CheckInterpolation(DM dm, PetscBool checkRestrict, PetscInt order, AppCtx *user)
721d71ae5a4SJacob Faibussowitsch {
722c4762a1bSJed Brown   PetscErrorCode (*exactFuncs[1])(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
723c4762a1bSJed Brown   PetscErrorCode (*exactFuncDers[1])(PetscInt, PetscReal, const PetscReal x[], const PetscReal n[], PetscInt, PetscScalar *u, void *ctx);
724c4762a1bSJed Brown   PetscReal n[3] = {1.0, 1.0, 1.0};
725c4762a1bSJed Brown   void     *exactCtxs[3];
726c4762a1bSJed Brown   DM        rdm, idm, fdm;
727c4762a1bSJed Brown   Mat       Interp;
728c4762a1bSJed Brown   Vec       iu, fu, scaling;
729c4762a1bSJed Brown   MPI_Comm  comm;
73030602db0SMatthew G. Knepley   PetscInt  dim;
731c4762a1bSJed Brown   PetscReal error, errorDer, tol = PETSC_SMALL;
732c4762a1bSJed Brown 
733c4762a1bSJed Brown   PetscFunctionBeginUser;
734c4762a1bSJed Brown   exactCtxs[0] = user;
735c4762a1bSJed Brown   exactCtxs[1] = user;
736c4762a1bSJed Brown   exactCtxs[2] = user;
7379566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
7389566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
7399566063dSJacob Faibussowitsch   PetscCall(DMRefine(dm, comm, &rdm));
7409566063dSJacob Faibussowitsch   PetscCall(DMSetCoarseDM(rdm, dm));
7419566063dSJacob Faibussowitsch   PetscCall(DMPlexSetRegularRefinement(rdm, user->convRefine));
74230602db0SMatthew G. Knepley   if (user->tree) {
743c4762a1bSJed Brown     DM refTree;
7449566063dSJacob Faibussowitsch     PetscCall(DMPlexGetReferenceTree(dm, &refTree));
7459566063dSJacob Faibussowitsch     PetscCall(DMPlexSetReferenceTree(rdm, refTree));
746c4762a1bSJed Brown   }
7479566063dSJacob Faibussowitsch   if (user->useDA) PetscCall(DMDASetVertexCoordinates(rdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
7489566063dSJacob Faibussowitsch   PetscCall(SetupSection(rdm, user));
749c4762a1bSJed Brown   /* Setup functions to approximate */
750c4762a1bSJed Brown   switch (order) {
751c4762a1bSJed Brown   case 0:
752c4762a1bSJed Brown     exactFuncs[0]    = constant;
753c4762a1bSJed Brown     exactFuncDers[0] = constantDer;
754c4762a1bSJed Brown     break;
755c4762a1bSJed Brown   case 1:
756c4762a1bSJed Brown     exactFuncs[0]    = linear;
757c4762a1bSJed Brown     exactFuncDers[0] = linearDer;
758c4762a1bSJed Brown     break;
759c4762a1bSJed Brown   case 2:
760c4762a1bSJed Brown     exactFuncs[0]    = quadratic;
761c4762a1bSJed Brown     exactFuncDers[0] = quadraticDer;
762c4762a1bSJed Brown     break;
763c4762a1bSJed Brown   case 3:
764c4762a1bSJed Brown     exactFuncs[0]    = cubic;
765c4762a1bSJed Brown     exactFuncDers[0] = cubicDer;
766c4762a1bSJed Brown     break;
767d71ae5a4SJacob Faibussowitsch   default:
768d71ae5a4SJacob Faibussowitsch     SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Could not determine functions to test for dimension %" PetscInt_FMT " order %" PetscInt_FMT, dim, order);
769c4762a1bSJed Brown   }
770c4762a1bSJed Brown   idm = checkRestrict ? rdm : dm;
771c4762a1bSJed Brown   fdm = checkRestrict ? dm : rdm;
7729566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(idm, &iu));
7739566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(fdm, &fu));
7749566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(dm, user));
7759566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(rdm, user));
7769566063dSJacob Faibussowitsch   PetscCall(DMCreateInterpolation(dm, rdm, &Interp, &scaling));
777c4762a1bSJed Brown   /* Project function into initial FE function space */
7789566063dSJacob Faibussowitsch   PetscCall(DMProjectFunction(idm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, iu));
779c4762a1bSJed Brown   /* Interpolate function into final FE function space */
7809371c9d4SSatish Balay   if (checkRestrict) {
7819371c9d4SSatish Balay     PetscCall(MatRestrict(Interp, iu, fu));
7829371c9d4SSatish Balay     PetscCall(VecPointwiseMult(fu, scaling, fu));
7839371c9d4SSatish Balay   } else PetscCall(MatInterpolate(Interp, iu, fu));
784c4762a1bSJed Brown   /* Compare approximation to exact in L_2 */
7859566063dSJacob Faibussowitsch   PetscCall(DMComputeL2Diff(fdm, 0.0, exactFuncs, exactCtxs, fu, &error));
7869566063dSJacob Faibussowitsch   PetscCall(DMComputeL2GradientDiff(fdm, 0.0, exactFuncDers, exactCtxs, fu, n, &errorDer));
787c4762a1bSJed Brown   /* Report result */
78863a3b9bcSJacob Faibussowitsch   if (error > tol) PetscCall(PetscPrintf(comm, "Interpolation tests FAIL for order %" PetscInt_FMT " at tolerance %g error %g\n", order, (double)tol, (double)error));
78963a3b9bcSJacob Faibussowitsch   else PetscCall(PetscPrintf(comm, "Interpolation tests pass for order %" PetscInt_FMT " at tolerance %g\n", order, (double)tol));
79063a3b9bcSJacob Faibussowitsch   if (errorDer > tol) PetscCall(PetscPrintf(comm, "Interpolation tests FAIL for order %" PetscInt_FMT " derivatives at tolerance %g error %g\n", order, (double)tol, (double)errorDer));
79163a3b9bcSJacob Faibussowitsch   else PetscCall(PetscPrintf(comm, "Interpolation tests pass for order %" PetscInt_FMT " derivatives at tolerance %g\n", order, (double)tol));
7929566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(idm, &iu));
7939566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(fdm, &fu));
7949566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Interp));
7959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&scaling));
7969566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&rdm));
7973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
798c4762a1bSJed Brown }
799c4762a1bSJed Brown 
800d71ae5a4SJacob Faibussowitsch static PetscErrorCode CheckConvergence(DM dm, PetscInt Nr, AppCtx *user)
801d71ae5a4SJacob Faibussowitsch {
802c4762a1bSJed Brown   DM odm = dm, rdm = NULL, cdm = NULL;
803c4762a1bSJed Brown   PetscErrorCode (*exactFuncs[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)                         = {trig};
804c4762a1bSJed Brown   PetscErrorCode (*exactFuncDers[1])(PetscInt dim, PetscReal time, const PetscReal x[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) = {trigDer};
805c4762a1bSJed Brown   void     *exactCtxs[3];
806c4762a1bSJed Brown   PetscInt  r, c, cStart, cEnd;
807c4762a1bSJed Brown   PetscReal errorOld, errorDerOld, error, errorDer, rel, len, lenOld;
808c4762a1bSJed Brown   double    p;
809c4762a1bSJed Brown 
810c4762a1bSJed Brown   PetscFunctionBeginUser;
8113ba16761SJacob Faibussowitsch   if (!user->convergence) PetscFunctionReturn(PETSC_SUCCESS);
812c4762a1bSJed Brown   exactCtxs[0] = user;
813c4762a1bSJed Brown   exactCtxs[1] = user;
814c4762a1bSJed Brown   exactCtxs[2] = user;
8159566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)odm));
816c4762a1bSJed Brown   if (!user->convRefine) {
817c4762a1bSJed Brown     for (r = 0; r < Nr; ++r) {
8189566063dSJacob Faibussowitsch       PetscCall(DMRefine(odm, PetscObjectComm((PetscObject)dm), &rdm));
8199566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&odm));
820c4762a1bSJed Brown       odm = rdm;
821c4762a1bSJed Brown     }
8229566063dSJacob Faibussowitsch     PetscCall(SetupSection(odm, user));
823c4762a1bSJed Brown   }
8249566063dSJacob Faibussowitsch   PetscCall(ComputeError(odm, exactFuncs, exactFuncDers, exactCtxs, &errorOld, &errorDerOld, user));
825c4762a1bSJed Brown   if (user->convRefine) {
826c4762a1bSJed Brown     for (r = 0; r < Nr; ++r) {
8279566063dSJacob Faibussowitsch       PetscCall(DMRefine(odm, PetscObjectComm((PetscObject)dm), &rdm));
8289566063dSJacob Faibussowitsch       if (user->useDA) PetscCall(DMDASetVertexCoordinates(rdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
8299566063dSJacob Faibussowitsch       PetscCall(SetupSection(rdm, user));
8309566063dSJacob Faibussowitsch       PetscCall(ComputeError(rdm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user));
831c4762a1bSJed Brown       p = PetscLog2Real(errorOld / error);
83263a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Function   convergence rate at refinement %" PetscInt_FMT ": %.2f\n", r, (double)p));
833c4762a1bSJed Brown       p = PetscLog2Real(errorDerOld / errorDer);
83463a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Derivative convergence rate at refinement %" PetscInt_FMT ": %.2f\n", r, (double)p));
8359566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&odm));
836c4762a1bSJed Brown       odm         = rdm;
837c4762a1bSJed Brown       errorOld    = error;
838c4762a1bSJed Brown       errorDerOld = errorDer;
839c4762a1bSJed Brown     }
840c4762a1bSJed Brown   } else {
8419566063dSJacob Faibussowitsch     /* PetscCall(ComputeLongestEdge(dm, &lenOld)); */
8429566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(odm, 0, &cStart, &cEnd));
843c4762a1bSJed Brown     lenOld = cEnd - cStart;
844c4762a1bSJed Brown     for (c = 0; c < Nr; ++c) {
8459566063dSJacob Faibussowitsch       PetscCall(DMCoarsen(odm, PetscObjectComm((PetscObject)dm), &cdm));
8469566063dSJacob Faibussowitsch       if (user->useDA) PetscCall(DMDASetVertexCoordinates(cdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
8479566063dSJacob Faibussowitsch       PetscCall(SetupSection(cdm, user));
8489566063dSJacob Faibussowitsch       PetscCall(ComputeError(cdm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user));
8499566063dSJacob Faibussowitsch       /* PetscCall(ComputeLongestEdge(cdm, &len)); */
8509566063dSJacob Faibussowitsch       PetscCall(DMPlexGetHeightStratum(cdm, 0, &cStart, &cEnd));
851c4762a1bSJed Brown       len = cEnd - cStart;
852c4762a1bSJed Brown       rel = error / errorOld;
853c4762a1bSJed Brown       p   = PetscLogReal(rel) / PetscLogReal(lenOld / len);
85463a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Function   convergence rate at coarsening %" PetscInt_FMT ": %.2f\n", c, (double)p));
855c4762a1bSJed Brown       rel = errorDer / errorDerOld;
856c4762a1bSJed Brown       p   = PetscLogReal(rel) / PetscLogReal(lenOld / len);
85763a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Derivative convergence rate at coarsening %" PetscInt_FMT ": %.2f\n", c, (double)p));
8589566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&odm));
859c4762a1bSJed Brown       odm         = cdm;
860c4762a1bSJed Brown       errorOld    = error;
861c4762a1bSJed Brown       errorDerOld = errorDer;
862c4762a1bSJed Brown       lenOld      = len;
863c4762a1bSJed Brown     }
864c4762a1bSJed Brown   }
8659566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&odm));
8663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
867c4762a1bSJed Brown }
868c4762a1bSJed Brown 
869d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
870d71ae5a4SJacob Faibussowitsch {
871c4762a1bSJed Brown   DM        dm;
872c4762a1bSJed Brown   AppCtx    user; /* user-defined work context */
87330602db0SMatthew G. Knepley   PetscInt  dim     = 2;
87430602db0SMatthew G. Knepley   PetscBool simplex = PETSC_FALSE;
875c4762a1bSJed Brown 
876327415f7SBarry Smith   PetscFunctionBeginUser;
8779566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
8789566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
8799566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
88030602db0SMatthew G. Knepley   if (!user.useDA) {
8819566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
8829566063dSJacob Faibussowitsch     PetscCall(DMPlexIsSimplex(dm, &simplex));
88330602db0SMatthew G. Knepley   }
8849566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetFromOptions(dm));
88530602db0SMatthew G. Knepley   user.numComponents = user.numComponents < 0 ? dim : user.numComponents;
8869566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_WORLD, dim, user.numComponents, simplex, NULL, user.qorder, &user.fe));
8879566063dSJacob Faibussowitsch   PetscCall(SetupSection(dm, &user));
8889566063dSJacob Faibussowitsch   if (user.testFEjacobian) PetscCall(TestFEJacobian(dm, &user));
8899566063dSJacob Faibussowitsch   if (user.testFVgrad) PetscCall(TestFVGrad(dm, &user));
8909566063dSJacob Faibussowitsch   if (user.testInjector) PetscCall(TestInjector(dm, &user));
8919566063dSJacob Faibussowitsch   PetscCall(CheckFunctions(dm, user.porder, &user));
892c4762a1bSJed Brown   {
893c4762a1bSJed Brown     PetscDualSpace dsp;
894c4762a1bSJed Brown     PetscInt       k;
895c4762a1bSJed Brown 
8969566063dSJacob Faibussowitsch     PetscCall(PetscFEGetDualSpace(user.fe, &dsp));
8979566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDeRahm(dsp, &k));
89830602db0SMatthew G. Knepley     if (dim == 2 && simplex == PETSC_TRUE && user.tree == PETSC_FALSE && k == 0) {
8999566063dSJacob Faibussowitsch       PetscCall(CheckInterpolation(dm, PETSC_FALSE, user.porder, &user));
9009566063dSJacob Faibussowitsch       PetscCall(CheckInterpolation(dm, PETSC_TRUE, user.porder, &user));
901c4762a1bSJed Brown     }
902c4762a1bSJed Brown   }
9039566063dSJacob Faibussowitsch   PetscCall(CheckConvergence(dm, 3, &user));
9049566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&user.fe));
9059566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
9069566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
907b122ec5aSJacob Faibussowitsch   return 0;
908c4762a1bSJed Brown }
909c4762a1bSJed Brown 
910c4762a1bSJed Brown /*TEST
911c4762a1bSJed Brown 
912c4762a1bSJed Brown   test:
913c4762a1bSJed Brown     suffix: 1
914c4762a1bSJed Brown     requires: triangle
915c4762a1bSJed Brown 
916c4762a1bSJed Brown   # 2D P_1 on a triangle
917c4762a1bSJed Brown   test:
918c4762a1bSJed Brown     suffix: p1_2d_0
919c4762a1bSJed Brown     requires: triangle
920c4762a1bSJed Brown     args: -petscspace_degree 1 -qorder 1 -convergence
921c4762a1bSJed Brown   test:
922c4762a1bSJed Brown     suffix: p1_2d_1
923c4762a1bSJed Brown     requires: triangle
924c4762a1bSJed Brown     args: -petscspace_degree 1 -qorder 1 -porder 1
925c4762a1bSJed Brown   test:
926c4762a1bSJed Brown     suffix: p1_2d_2
927c4762a1bSJed Brown     requires: triangle
928c4762a1bSJed Brown     args: -petscspace_degree 1 -qorder 1 -porder 2
929c4762a1bSJed Brown   test:
930c4762a1bSJed Brown     suffix: p1_2d_3
931e788613bSJoe Wallwork     requires: triangle mmg
932c4762a1bSJed Brown     args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0
933c4762a1bSJed Brown   test:
934c4762a1bSJed Brown     suffix: p1_2d_4
935e788613bSJoe Wallwork     requires: triangle mmg
936c4762a1bSJed Brown     args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0
937c4762a1bSJed Brown   test:
938c4762a1bSJed Brown     suffix: p1_2d_5
939e788613bSJoe Wallwork     requires: triangle mmg
940c4762a1bSJed Brown     args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0
941c4762a1bSJed Brown 
942c4762a1bSJed Brown   # 3D P_1 on a tetrahedron
943c4762a1bSJed Brown   test:
944c4762a1bSJed Brown     suffix: p1_3d_0
945c4762a1bSJed Brown     requires: ctetgen
94630602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -convergence
947c4762a1bSJed Brown   test:
948c4762a1bSJed Brown     suffix: p1_3d_1
949c4762a1bSJed Brown     requires: ctetgen
95030602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -porder 1
951c4762a1bSJed Brown   test:
952c4762a1bSJed Brown     suffix: p1_3d_2
953c4762a1bSJed Brown     requires: ctetgen
95430602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -porder 2
955c4762a1bSJed Brown   test:
956c4762a1bSJed Brown     suffix: p1_3d_3
957e788613bSJoe Wallwork     requires: ctetgen mmg
95830602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0
959c4762a1bSJed Brown   test:
960c4762a1bSJed Brown     suffix: p1_3d_4
961e788613bSJoe Wallwork     requires: ctetgen mmg
96230602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0
963c4762a1bSJed Brown   test:
964c4762a1bSJed Brown     suffix: p1_3d_5
965e788613bSJoe Wallwork     requires: ctetgen mmg
96630602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0
967c4762a1bSJed Brown 
968c4762a1bSJed Brown   # 2D P_2 on a triangle
969c4762a1bSJed Brown   test:
970c4762a1bSJed Brown     suffix: p2_2d_0
971c4762a1bSJed Brown     requires: triangle
972c4762a1bSJed Brown     args: -petscspace_degree 2 -qorder 2 -convergence
973c4762a1bSJed Brown   test:
974c4762a1bSJed Brown     suffix: p2_2d_1
975c4762a1bSJed Brown     requires: triangle
976c4762a1bSJed Brown     args: -petscspace_degree 2 -qorder 2 -porder 1
977c4762a1bSJed Brown   test:
978c4762a1bSJed Brown     suffix: p2_2d_2
979c4762a1bSJed Brown     requires: triangle
980c4762a1bSJed Brown     args: -petscspace_degree 2 -qorder 2 -porder 2
981c4762a1bSJed Brown   test:
982c4762a1bSJed Brown     suffix: p2_2d_3
983e788613bSJoe Wallwork     requires: triangle mmg
984c4762a1bSJed Brown     args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -convergence -conv_refine 0
985c4762a1bSJed Brown   test:
986c4762a1bSJed Brown     suffix: p2_2d_4
987e788613bSJoe Wallwork     requires: triangle mmg
988c4762a1bSJed Brown     args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 1 -conv_refine 0
989c4762a1bSJed Brown   test:
990c4762a1bSJed Brown     suffix: p2_2d_5
991e788613bSJoe Wallwork     requires: triangle mmg
992c4762a1bSJed Brown     args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 2 -conv_refine 0
993c4762a1bSJed Brown 
994c4762a1bSJed Brown   # 3D P_2 on a tetrahedron
995c4762a1bSJed Brown   test:
996c4762a1bSJed Brown     suffix: p2_3d_0
997c4762a1bSJed Brown     requires: ctetgen
99830602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -convergence
999c4762a1bSJed Brown   test:
1000c4762a1bSJed Brown     suffix: p2_3d_1
1001c4762a1bSJed Brown     requires: ctetgen
100230602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -porder 1
1003c4762a1bSJed Brown   test:
1004c4762a1bSJed Brown     suffix: p2_3d_2
1005c4762a1bSJed Brown     requires: ctetgen
100630602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -porder 2
1007c4762a1bSJed Brown   test:
1008c4762a1bSJed Brown     suffix: p2_3d_3
1009e788613bSJoe Wallwork     requires: ctetgen mmg
101030602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -convergence -conv_refine 0
1011c4762a1bSJed Brown   test:
1012c4762a1bSJed Brown     suffix: p2_3d_4
1013e788613bSJoe Wallwork     requires: ctetgen mmg
101430602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 1 -conv_refine 0
1015c4762a1bSJed Brown   test:
1016c4762a1bSJed Brown     suffix: p2_3d_5
1017e788613bSJoe Wallwork     requires: ctetgen mmg
101830602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 2 -conv_refine 0
1019c4762a1bSJed Brown 
1020c4762a1bSJed Brown   # 2D Q_1 on a quadrilaterial DA
1021c4762a1bSJed Brown   test:
1022c4762a1bSJed Brown     suffix: q1_2d_da_0
102399a192c5SJunchao Zhang     requires: broken
102430602db0SMatthew G. Knepley     args: -use_da 1 -petscspace_degree 1 -qorder 1 -convergence
1025c4762a1bSJed Brown   test:
1026c4762a1bSJed Brown     suffix: q1_2d_da_1
102799a192c5SJunchao Zhang     requires: broken
102830602db0SMatthew G. Knepley     args: -use_da 1 -petscspace_degree 1 -qorder 1 -porder 1
1029c4762a1bSJed Brown   test:
1030c4762a1bSJed Brown     suffix: q1_2d_da_2
103199a192c5SJunchao Zhang     requires: broken
103230602db0SMatthew G. Knepley     args: -use_da 1 -petscspace_degree 1 -qorder 1 -porder 2
1033c4762a1bSJed Brown 
1034c4762a1bSJed Brown   # 2D Q_1 on a quadrilaterial Plex
1035c4762a1bSJed Brown   test:
1036c4762a1bSJed Brown     suffix: q1_2d_plex_0
103730602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -convergence
1038c4762a1bSJed Brown   test:
1039c4762a1bSJed Brown     suffix: q1_2d_plex_1
104030602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 1
1041c4762a1bSJed Brown   test:
1042c4762a1bSJed Brown     suffix: q1_2d_plex_2
104330602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 2
1044c4762a1bSJed Brown   test:
1045c4762a1bSJed Brown     suffix: q1_2d_plex_3
104630602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 1 -shear_coords
1047c4762a1bSJed Brown   test:
1048c4762a1bSJed Brown     suffix: q1_2d_plex_4
104930602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 2 -shear_coords
1050c4762a1bSJed Brown   test:
1051c4762a1bSJed Brown     suffix: q1_2d_plex_5
105230602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 0 -non_affine_coords -convergence
1053c4762a1bSJed Brown   test:
1054c4762a1bSJed Brown     suffix: q1_2d_plex_6
105530602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 1 -non_affine_coords -convergence
1056c4762a1bSJed Brown   test:
1057c4762a1bSJed Brown     suffix: q1_2d_plex_7
105830602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 2 -non_affine_coords -convergence
1059c4762a1bSJed Brown 
1060c4762a1bSJed Brown   # 2D Q_2 on a quadrilaterial
1061c4762a1bSJed Brown   test:
1062c4762a1bSJed Brown     suffix: q2_2d_plex_0
106330602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -convergence
1064c4762a1bSJed Brown   test:
1065c4762a1bSJed Brown     suffix: q2_2d_plex_1
106630602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1
1067c4762a1bSJed Brown   test:
1068c4762a1bSJed Brown     suffix: q2_2d_plex_2
106930602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2
1070c4762a1bSJed Brown   test:
1071c4762a1bSJed Brown     suffix: q2_2d_plex_3
107230602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 -shear_coords
1073c4762a1bSJed Brown   test:
1074c4762a1bSJed Brown     suffix: q2_2d_plex_4
107530602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 -shear_coords
1076c4762a1bSJed Brown   test:
1077c4762a1bSJed Brown     suffix: q2_2d_plex_5
107830602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 0 -non_affine_coords -convergence
1079c4762a1bSJed Brown   test:
1080c4762a1bSJed Brown     suffix: q2_2d_plex_6
108130602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 1 -non_affine_coords -convergence
1082c4762a1bSJed Brown   test:
1083c4762a1bSJed Brown     suffix: q2_2d_plex_7
108430602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 2 -non_affine_coords -convergence
1085c4762a1bSJed Brown 
1086c4762a1bSJed Brown   # 2D P_3 on a triangle
1087c4762a1bSJed Brown   test:
1088c4762a1bSJed Brown     suffix: p3_2d_0
1089c4762a1bSJed Brown     requires: triangle !single
1090c4762a1bSJed Brown     args: -petscspace_degree 3 -qorder 3 -convergence
1091c4762a1bSJed Brown   test:
1092c4762a1bSJed Brown     suffix: p3_2d_1
1093c4762a1bSJed Brown     requires: triangle !single
1094c4762a1bSJed Brown     args: -petscspace_degree 3 -qorder 3 -porder 1
1095c4762a1bSJed Brown   test:
1096c4762a1bSJed Brown     suffix: p3_2d_2
1097c4762a1bSJed Brown     requires: triangle !single
1098c4762a1bSJed Brown     args: -petscspace_degree 3 -qorder 3 -porder 2
1099c4762a1bSJed Brown   test:
1100c4762a1bSJed Brown     suffix: p3_2d_3
1101c4762a1bSJed Brown     requires: triangle !single
1102c4762a1bSJed Brown     args: -petscspace_degree 3 -qorder 3 -porder 3
1103c4762a1bSJed Brown   test:
1104c4762a1bSJed Brown     suffix: p3_2d_4
1105e788613bSJoe Wallwork     requires: triangle mmg
1106c4762a1bSJed Brown     args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -convergence -conv_refine 0
1107c4762a1bSJed Brown   test:
1108c4762a1bSJed Brown     suffix: p3_2d_5
1109e788613bSJoe Wallwork     requires: triangle mmg
1110c4762a1bSJed Brown     args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -porder 1 -conv_refine 0
1111c4762a1bSJed Brown   test:
1112c4762a1bSJed Brown     suffix: p3_2d_6
1113e788613bSJoe Wallwork     requires: triangle mmg
1114c4762a1bSJed Brown     args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -porder 3 -conv_refine 0
1115c4762a1bSJed Brown 
1116c4762a1bSJed Brown   # 2D Q_3 on a quadrilaterial
1117c4762a1bSJed Brown   test:
1118c4762a1bSJed Brown     suffix: q3_2d_0
111999a192c5SJunchao Zhang     requires: !single
112030602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -convergence
1121c4762a1bSJed Brown   test:
1122c4762a1bSJed Brown     suffix: q3_2d_1
112399a192c5SJunchao Zhang     requires: !single
112430602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 1
1125c4762a1bSJed Brown   test:
1126c4762a1bSJed Brown     suffix: q3_2d_2
112799a192c5SJunchao Zhang     requires: !single
112830602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 2
1129c4762a1bSJed Brown   test:
1130c4762a1bSJed Brown     suffix: q3_2d_3
113199a192c5SJunchao Zhang     requires: !single
113230602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 3
1133c4762a1bSJed Brown 
1134c4762a1bSJed Brown   # 2D P_1disc on a triangle/quadrilateral
1135c4762a1bSJed Brown   test:
1136c4762a1bSJed Brown     suffix: p1d_2d_0
1137c4762a1bSJed Brown     requires: triangle
1138c4762a1bSJed Brown     args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -convergence
1139c4762a1bSJed Brown   test:
1140c4762a1bSJed Brown     suffix: p1d_2d_1
1141c4762a1bSJed Brown     requires: triangle
1142c4762a1bSJed Brown     args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 1
1143c4762a1bSJed Brown   test:
1144c4762a1bSJed Brown     suffix: p1d_2d_2
1145c4762a1bSJed Brown     requires: triangle
1146c4762a1bSJed Brown     args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 2
1147c4762a1bSJed Brown   test:
1148c4762a1bSJed Brown     suffix: p1d_2d_3
1149c4762a1bSJed Brown     requires: triangle
115030602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -convergence
1151c4762a1bSJed Brown     filter: sed  -e "s/convergence rate at refinement 0: 2/convergence rate at refinement 0: 1.9/g"
1152c4762a1bSJed Brown   test:
1153c4762a1bSJed Brown     suffix: p1d_2d_4
1154c4762a1bSJed Brown     requires: triangle
115530602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 1
1156c4762a1bSJed Brown   test:
1157c4762a1bSJed Brown     suffix: p1d_2d_5
1158c4762a1bSJed Brown     requires: triangle
115930602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 2
1160c4762a1bSJed Brown 
1161c4762a1bSJed Brown   # 2D BDM_1 on a triangle
1162c4762a1bSJed Brown   test:
1163c4762a1bSJed Brown     suffix: bdm1_2d_0
1164c4762a1bSJed Brown     requires: triangle
1165c4762a1bSJed Brown     args: -petscspace_degree 1 -petscdualspace_type bdm \
1166c4762a1bSJed Brown           -num_comp 2 -qorder 1 -convergence
1167c4762a1bSJed Brown   test:
1168c4762a1bSJed Brown     suffix: bdm1_2d_1
1169c4762a1bSJed Brown     requires: triangle
1170c4762a1bSJed Brown     args: -petscspace_degree 1 -petscdualspace_type bdm \
1171c4762a1bSJed Brown           -num_comp 2 -qorder 1 -porder 1
1172c4762a1bSJed Brown   test:
1173c4762a1bSJed Brown     suffix: bdm1_2d_2
1174c4762a1bSJed Brown     requires: triangle
1175c4762a1bSJed Brown     args: -petscspace_degree 1 -petscdualspace_type bdm \
1176c4762a1bSJed Brown           -num_comp 2 -qorder 1 -porder 2
1177c4762a1bSJed Brown 
1178c4762a1bSJed Brown   # 2D BDM_1 on a quadrilateral
1179c4762a1bSJed Brown   test:
1180c4762a1bSJed Brown     suffix: bdm1q_2d_0
1181c4762a1bSJed Brown     requires: triangle
1182c4762a1bSJed Brown     args: -petscspace_degree 1 -petscdualspace_type bdm \
11833f27d899SToby Isaac           -petscdualspace_lagrange_tensor 1 \
118430602db0SMatthew G. Knepley           -dm_plex_simplex 0 -num_comp 2 -qorder 1 -convergence
1185c4762a1bSJed Brown   test:
1186c4762a1bSJed Brown     suffix: bdm1q_2d_1
1187c4762a1bSJed Brown     requires: triangle
1188c4762a1bSJed Brown     args: -petscspace_degree 1 -petscdualspace_type bdm \
11893f27d899SToby Isaac           -petscdualspace_lagrange_tensor 1 \
119030602db0SMatthew G. Knepley           -dm_plex_simplex 0 -num_comp 2 -qorder 1 -porder 1
1191c4762a1bSJed Brown   test:
1192c4762a1bSJed Brown     suffix: bdm1q_2d_2
1193c4762a1bSJed Brown     requires: triangle
1194c4762a1bSJed Brown     args: -petscspace_degree 1 -petscdualspace_type bdm \
11953f27d899SToby Isaac           -petscdualspace_lagrange_tensor 1 \
119630602db0SMatthew G. Knepley           -dm_plex_simplex 0 -num_comp 2 -qorder 1 -porder 2
1197c4762a1bSJed Brown 
1198c4762a1bSJed Brown   # Test high order quadrature
1199c4762a1bSJed Brown   test:
1200c4762a1bSJed Brown     suffix: p1_quad_2
1201c4762a1bSJed Brown     requires: triangle
1202c4762a1bSJed Brown     args: -petscspace_degree 1 -qorder 2 -porder 1
1203c4762a1bSJed Brown   test:
1204c4762a1bSJed Brown     suffix: p1_quad_5
1205c4762a1bSJed Brown     requires: triangle
1206c4762a1bSJed Brown     args: -petscspace_degree 1 -qorder 5 -porder 1
1207c4762a1bSJed Brown   test:
1208c4762a1bSJed Brown     suffix: p2_quad_3
1209c4762a1bSJed Brown     requires: triangle
1210c4762a1bSJed Brown     args: -petscspace_degree 2 -qorder 3 -porder 2
1211c4762a1bSJed Brown   test:
1212c4762a1bSJed Brown     suffix: p2_quad_5
1213c4762a1bSJed Brown     requires: triangle
1214c4762a1bSJed Brown     args: -petscspace_degree 2 -qorder 5 -porder 2
1215c4762a1bSJed Brown   test:
1216c4762a1bSJed Brown     suffix: q1_quad_2
121730602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 2 -porder 1
1218c4762a1bSJed Brown   test:
1219c4762a1bSJed Brown     suffix: q1_quad_5
122030602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 5 -porder 1
1221c4762a1bSJed Brown   test:
1222c4762a1bSJed Brown     suffix: q2_quad_3
122330602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 3 -porder 1
1224c4762a1bSJed Brown   test:
1225c4762a1bSJed Brown     suffix: q2_quad_5
122630602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 5 -porder 1
1227c4762a1bSJed Brown 
1228c4762a1bSJed Brown   # Nonconforming tests
1229c4762a1bSJed Brown   test:
1230c4762a1bSJed Brown     suffix: constraints
123130602db0SMatthew G. Knepley     args: -dm_coord_space 0 -dm_plex_simplex 0 -petscspace_type tensor -petscspace_degree 1 -qorder 0 -constraints
1232c4762a1bSJed Brown   test:
1233c4762a1bSJed Brown     suffix: nonconforming_tensor_2
1234c4762a1bSJed Brown     nsize: 4
123530602db0SMatthew G. Knepley     args: -dist_dm_distribute -test_fe_jacobian -test_injector -petscpartitioner_type simple -tree -dm_plex_simplex 0 -dm_plex_max_projection_height 1 -petscspace_type tensor -petscspace_degree 2 -qorder 2 -dm_view ascii::ASCII_INFO_DETAIL
1236c4762a1bSJed Brown   test:
1237c4762a1bSJed Brown     suffix: nonconforming_tensor_3
1238c4762a1bSJed Brown     nsize: 4
123930602db0SMatthew G. Knepley     args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_plex_max_projection_height 2 -petscspace_type tensor -petscspace_degree 1 -qorder 1 -dm_view ascii::ASCII_INFO_DETAIL
1240c4762a1bSJed Brown   test:
1241c4762a1bSJed Brown     suffix: nonconforming_tensor_2_fv
1242c4762a1bSJed Brown     nsize: 4
124330602db0SMatthew G. Knepley     args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -dm_plex_simplex 0 -num_comp 2
1244c4762a1bSJed Brown   test:
1245c4762a1bSJed Brown     suffix: nonconforming_tensor_3_fv
1246c4762a1bSJed Brown     nsize: 4
124730602db0SMatthew G. Knepley     args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -num_comp 3
1248c4762a1bSJed Brown   test:
1249c4762a1bSJed Brown     suffix: nonconforming_tensor_2_hi
1250c4762a1bSJed Brown     requires: !single
1251c4762a1bSJed Brown     nsize: 4
125230602db0SMatthew G. Knepley     args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_simplex 0 -dm_plex_max_projection_height 1 -petscspace_type tensor -petscspace_degree 4 -qorder 4
1253c4762a1bSJed Brown   test:
1254c4762a1bSJed Brown     suffix: nonconforming_tensor_3_hi
1255c4762a1bSJed Brown     requires: !single skip
1256c4762a1bSJed Brown     nsize: 4
125730602db0SMatthew G. Knepley     args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_plex_max_projection_height 2 -petscspace_type tensor -petscspace_degree 4 -qorder 4
1258c4762a1bSJed Brown   test:
1259c4762a1bSJed Brown     suffix: nonconforming_simplex_2
1260c4762a1bSJed Brown     requires: triangle
1261c4762a1bSJed Brown     nsize: 4
126230602db0SMatthew G. Knepley     args: -dist_dm_distribute -test_fe_jacobian -test_injector -petscpartitioner_type simple -tree -dm_plex_max_projection_height 1 -petscspace_degree 2 -qorder 2 -dm_view ascii::ASCII_INFO_DETAIL
1263c4762a1bSJed Brown   test:
1264c4762a1bSJed Brown     suffix: nonconforming_simplex_2_hi
1265c4762a1bSJed Brown     requires: triangle !single
1266c4762a1bSJed Brown     nsize: 4
126730602db0SMatthew G. Knepley     args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_max_projection_height 1 -petscspace_degree 4 -qorder 4
1268c4762a1bSJed Brown   test:
1269c4762a1bSJed Brown     suffix: nonconforming_simplex_2_fv
1270c4762a1bSJed Brown     requires: triangle
1271c4762a1bSJed Brown     nsize: 4
127230602db0SMatthew G. Knepley     args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -num_comp 2
1273c4762a1bSJed Brown   test:
1274c4762a1bSJed Brown     suffix: nonconforming_simplex_3
1275c4762a1bSJed Brown     requires: ctetgen
1276c4762a1bSJed Brown     nsize: 4
127730602db0SMatthew G. Knepley     args: -dist_dm_distribute -test_fe_jacobian -test_injector -petscpartitioner_type simple -tree -dm_plex_dim 3 -dm_plex_max_projection_height 2 -petscspace_degree 2 -qorder 2 -dm_view ascii::ASCII_INFO_DETAIL
1278c4762a1bSJed Brown   test:
1279c4762a1bSJed Brown     suffix: nonconforming_simplex_3_hi
1280c4762a1bSJed Brown     requires: ctetgen skip
1281c4762a1bSJed Brown     nsize: 4
128230602db0SMatthew G. Knepley     args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_dim 3 -dm_plex_max_projection_height 2 -petscspace_degree 4 -qorder 4
1283c4762a1bSJed Brown   test:
1284c4762a1bSJed Brown     suffix: nonconforming_simplex_3_fv
1285c4762a1bSJed Brown     requires: ctetgen
1286c4762a1bSJed Brown     nsize: 4
128730602db0SMatthew G. Knepley     args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -dm_plex_dim 3 -num_comp 3
1288c4762a1bSJed Brown 
1289d21efd2eSMatthew G. Knepley   # 3D WXY on a triangular prism
1290d21efd2eSMatthew G. Knepley   test:
1291d21efd2eSMatthew G. Knepley     suffix: wxy_0
1292d21efd2eSMatthew G. Knepley     args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -qorder 2 -porder 0 \
1293e239af90SMatthew G. Knepley           -petscspace_type sum \
1294e239af90SMatthew G. Knepley           -petscspace_variables 3 \
1295e239af90SMatthew G. Knepley           -petscspace_components 3 \
1296e239af90SMatthew G. Knepley           -petscspace_sum_spaces 2 \
1297e239af90SMatthew G. Knepley           -petscspace_sum_concatenate false \
1298e239af90SMatthew G. Knepley           -sumcomp_0_petscspace_variables 3 \
1299e239af90SMatthew G. Knepley           -sumcomp_0_petscspace_components 3 \
1300e239af90SMatthew G. Knepley           -sumcomp_0_petscspace_degree 1 \
1301e239af90SMatthew G. Knepley           -sumcomp_1_petscspace_variables 3 \
1302e239af90SMatthew G. Knepley           -sumcomp_1_petscspace_components 3 \
1303e239af90SMatthew G. Knepley           -sumcomp_1_petscspace_type wxy \
1304e239af90SMatthew G. Knepley           -petscdualspace_refcell triangular_prism \
1305e239af90SMatthew G. Knepley           -petscdualspace_form_degree 0 \
1306e239af90SMatthew G. Knepley           -petscdualspace_order 1 \
1307e239af90SMatthew G. Knepley           -petscdualspace_components 3
1308d21efd2eSMatthew G. Knepley 
1309c4762a1bSJed Brown TEST*/
1310c4762a1bSJed Brown 
1311c4762a1bSJed Brown /*
1312c4762a1bSJed Brown    # 2D Q_2 on a quadrilaterial Plex
1313c4762a1bSJed Brown   test:
1314c4762a1bSJed Brown     suffix: q2_2d_plex_0
131530602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -convergence
1316c4762a1bSJed Brown   test:
1317c4762a1bSJed Brown     suffix: q2_2d_plex_1
131830602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1
1319c4762a1bSJed Brown   test:
1320c4762a1bSJed Brown     suffix: q2_2d_plex_2
132130602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2
1322c4762a1bSJed Brown   test:
1323c4762a1bSJed Brown     suffix: q2_2d_plex_3
132430602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 -shear_coords
1325c4762a1bSJed Brown   test:
1326c4762a1bSJed Brown     suffix: q2_2d_plex_4
132730602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 -shear_coords
1328c4762a1bSJed Brown   test:
1329c4762a1bSJed Brown     suffix: q2_2d_plex_5
133030602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 0 -non_affine_coords
1331c4762a1bSJed Brown   test:
1332c4762a1bSJed Brown     suffix: q2_2d_plex_6
133330602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 1 -non_affine_coords
1334c4762a1bSJed Brown   test:
1335c4762a1bSJed Brown     suffix: q2_2d_plex_7
133630602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 2 -non_affine_coords
1337c4762a1bSJed Brown 
1338c4762a1bSJed Brown   test:
1339c4762a1bSJed Brown     suffix: p1d_2d_6
1340e788613bSJoe Wallwork     requires: mmg
1341c4762a1bSJed Brown     args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0
1342c4762a1bSJed Brown   test:
1343c4762a1bSJed Brown     suffix: p1d_2d_7
1344e788613bSJoe Wallwork     requires: mmg
1345c4762a1bSJed Brown     args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0
1346c4762a1bSJed Brown   test:
1347c4762a1bSJed Brown     suffix: p1d_2d_8
1348e788613bSJoe Wallwork     requires: mmg
1349c4762a1bSJed Brown     args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0
1350c4762a1bSJed Brown */
1351