xref: /petsc/src/dm/impls/plex/tests/ex3.c (revision 4446c3cd466bac3b2df10d1a17c9bf818279221a)
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>
10*4446c3cdSksagiyam #include <petscsf.h>
11c4762a1bSJed Brown 
12c4762a1bSJed Brown typedef struct {
13c4762a1bSJed Brown   /* Domain and mesh definition */
14c4762a1bSJed Brown   PetscBool useDA;           /* Flag DMDA tensor product mesh */
15c4762a1bSJed Brown   PetscBool shearCoords;     /* Flag for shear transform */
16c4762a1bSJed Brown   PetscBool nonaffineCoords; /* Flag for non-affine transform */
17c4762a1bSJed Brown   /* Element definition */
18c4762a1bSJed Brown   PetscInt qorder;        /* Order of the quadrature */
19c4762a1bSJed Brown   PetscInt numComponents; /* Number of field components */
20c4762a1bSJed Brown   PetscFE  fe;            /* The finite element */
21c4762a1bSJed Brown   /* Testing space */
22c4762a1bSJed Brown   PetscInt  porder;         /* Order of polynomials to test */
23c4762a1bSJed Brown   PetscBool convergence;    /* Test for order of convergence */
24c4762a1bSJed Brown   PetscBool convRefine;     /* Test for convergence using refinement, otherwise use coarsening */
25c4762a1bSJed Brown   PetscBool constraints;    /* Test local constraints */
26c4762a1bSJed Brown   PetscBool tree;           /* Test tree routines */
27c4762a1bSJed Brown   PetscBool testFEjacobian; /* Test finite element Jacobian assembly */
28c4762a1bSJed Brown   PetscBool testFVgrad;     /* Test finite difference gradient routine */
29c4762a1bSJed Brown   PetscBool testInjector;   /* Test finite element injection routines */
30c4762a1bSJed Brown   PetscInt  treeCell;       /* Cell to refine in tree test */
31c4762a1bSJed Brown   PetscReal constants[3];   /* Constant values for each dimension */
32c4762a1bSJed Brown } AppCtx;
33c4762a1bSJed Brown 
34c4762a1bSJed Brown /* u = 1 */
35d71ae5a4SJacob Faibussowitsch PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx)
36d71ae5a4SJacob Faibussowitsch {
37c4762a1bSJed Brown   AppCtx  *user = (AppCtx *)ctx;
38c4762a1bSJed Brown   PetscInt d;
3930602db0SMatthew G. Knepley   for (d = 0; d < dim; ++d) u[d] = user->constants[d];
403ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
41c4762a1bSJed Brown }
42d71ae5a4SJacob Faibussowitsch PetscErrorCode constantDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx)
43d71ae5a4SJacob Faibussowitsch {
44c4762a1bSJed Brown   PetscInt d;
4530602db0SMatthew G. Knepley   for (d = 0; d < dim; ++d) u[d] = 0.0;
463ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
47c4762a1bSJed Brown }
48c4762a1bSJed Brown 
49c4762a1bSJed Brown /* u = x */
50d71ae5a4SJacob Faibussowitsch PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx)
51d71ae5a4SJacob Faibussowitsch {
52c4762a1bSJed Brown   PetscInt d;
53c4762a1bSJed Brown   for (d = 0; d < dim; ++d) u[d] = coords[d];
543ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
55c4762a1bSJed Brown }
56d71ae5a4SJacob Faibussowitsch PetscErrorCode linearDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx)
57d71ae5a4SJacob Faibussowitsch {
58c4762a1bSJed Brown   PetscInt d, e;
59c4762a1bSJed Brown   for (d = 0; d < dim; ++d) {
60c4762a1bSJed Brown     u[d] = 0.0;
61c4762a1bSJed Brown     for (e = 0; e < dim; ++e) u[d] += (d == e ? 1.0 : 0.0) * n[e];
62c4762a1bSJed Brown   }
633ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
64c4762a1bSJed Brown }
65c4762a1bSJed Brown 
66c4762a1bSJed Brown /* u = x^2 or u = (x^2, xy) or u = (xy, yz, zx) */
67d71ae5a4SJacob Faibussowitsch PetscErrorCode quadratic(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx)
68d71ae5a4SJacob Faibussowitsch {
699371c9d4SSatish Balay   if (dim > 2) {
709371c9d4SSatish Balay     u[0] = coords[0] * coords[1];
719371c9d4SSatish Balay     u[1] = coords[1] * coords[2];
729371c9d4SSatish Balay     u[2] = coords[2] * coords[0];
739371c9d4SSatish Balay   } else if (dim > 1) {
749371c9d4SSatish Balay     u[0] = coords[0] * coords[0];
759371c9d4SSatish Balay     u[1] = coords[0] * coords[1];
769371c9d4SSatish Balay   } else if (dim > 0) {
779371c9d4SSatish Balay     u[0] = coords[0] * coords[0];
789371c9d4SSatish Balay   }
793ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
80c4762a1bSJed Brown }
81d71ae5a4SJacob Faibussowitsch PetscErrorCode quadraticDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx)
82d71ae5a4SJacob Faibussowitsch {
839371c9d4SSatish Balay   if (dim > 2) {
849371c9d4SSatish Balay     u[0] = coords[1] * n[0] + coords[0] * n[1];
859371c9d4SSatish Balay     u[1] = coords[2] * n[1] + coords[1] * n[2];
869371c9d4SSatish Balay     u[2] = coords[2] * n[0] + coords[0] * n[2];
879371c9d4SSatish Balay   } else if (dim > 1) {
889371c9d4SSatish Balay     u[0] = 2.0 * coords[0] * n[0];
899371c9d4SSatish Balay     u[1] = coords[1] * n[0] + coords[0] * n[1];
909371c9d4SSatish Balay   } else if (dim > 0) {
919371c9d4SSatish Balay     u[0] = 2.0 * coords[0] * n[0];
929371c9d4SSatish Balay   }
933ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
94c4762a1bSJed Brown }
95c4762a1bSJed Brown 
96c4762a1bSJed Brown /* u = x^3 or u = (x^3, x^2y) or u = (x^2y, y^2z, z^2x) */
97d71ae5a4SJacob Faibussowitsch PetscErrorCode cubic(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx)
98d71ae5a4SJacob Faibussowitsch {
999371c9d4SSatish Balay   if (dim > 2) {
1009371c9d4SSatish Balay     u[0] = coords[0] * coords[0] * coords[1];
1019371c9d4SSatish Balay     u[1] = coords[1] * coords[1] * coords[2];
1029371c9d4SSatish Balay     u[2] = coords[2] * coords[2] * coords[0];
1039371c9d4SSatish Balay   } else if (dim > 1) {
1049371c9d4SSatish Balay     u[0] = coords[0] * coords[0] * coords[0];
1059371c9d4SSatish Balay     u[1] = coords[0] * coords[0] * coords[1];
1069371c9d4SSatish Balay   } else if (dim > 0) {
1079371c9d4SSatish Balay     u[0] = coords[0] * coords[0] * coords[0];
1089371c9d4SSatish Balay   }
1093ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
110c4762a1bSJed Brown }
111d71ae5a4SJacob Faibussowitsch PetscErrorCode cubicDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx)
112d71ae5a4SJacob Faibussowitsch {
1139371c9d4SSatish Balay   if (dim > 2) {
1149371c9d4SSatish Balay     u[0] = 2.0 * coords[0] * coords[1] * n[0] + coords[0] * coords[0] * n[1];
1159371c9d4SSatish Balay     u[1] = 2.0 * coords[1] * coords[2] * n[1] + coords[1] * coords[1] * n[2];
1169371c9d4SSatish Balay     u[2] = 2.0 * coords[2] * coords[0] * n[2] + coords[2] * coords[2] * n[0];
1179371c9d4SSatish Balay   } else if (dim > 1) {
1189371c9d4SSatish Balay     u[0] = 3.0 * coords[0] * coords[0] * n[0];
1199371c9d4SSatish Balay     u[1] = 2.0 * coords[0] * coords[1] * n[0] + coords[0] * coords[0] * n[1];
1209371c9d4SSatish Balay   } else if (dim > 0) {
1219371c9d4SSatish Balay     u[0] = 3.0 * coords[0] * coords[0] * n[0];
1229371c9d4SSatish Balay   }
1233ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
124c4762a1bSJed Brown }
125c4762a1bSJed Brown 
126c4762a1bSJed Brown /* u = tanh(x) */
127d71ae5a4SJacob Faibussowitsch PetscErrorCode trig(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx)
128d71ae5a4SJacob Faibussowitsch {
129c4762a1bSJed Brown   PetscInt d;
13030602db0SMatthew G. Knepley   for (d = 0; d < dim; ++d) u[d] = PetscTanhReal(coords[d] - 0.5);
1313ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
132c4762a1bSJed Brown }
133d71ae5a4SJacob Faibussowitsch PetscErrorCode trigDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx)
134d71ae5a4SJacob Faibussowitsch {
135c4762a1bSJed Brown   PetscInt d;
13630602db0SMatthew G. Knepley   for (d = 0; d < dim; ++d) u[d] = 1.0 / PetscSqr(PetscCoshReal(coords[d] - 0.5)) * n[d];
1373ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
138c4762a1bSJed Brown }
139c4762a1bSJed Brown 
140d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
141d71ae5a4SJacob Faibussowitsch {
142c4762a1bSJed Brown   PetscInt n = 3;
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   PetscFunctionBeginUser;
14530602db0SMatthew G. Knepley   options->useDA           = PETSC_FALSE;
146c4762a1bSJed Brown   options->shearCoords     = PETSC_FALSE;
147c4762a1bSJed Brown   options->nonaffineCoords = PETSC_FALSE;
148c4762a1bSJed Brown   options->qorder          = 0;
149c4762a1bSJed Brown   options->numComponents   = PETSC_DEFAULT;
150c4762a1bSJed Brown   options->porder          = 0;
151c4762a1bSJed Brown   options->convergence     = PETSC_FALSE;
152c4762a1bSJed Brown   options->convRefine      = PETSC_TRUE;
153c4762a1bSJed Brown   options->constraints     = PETSC_FALSE;
154c4762a1bSJed Brown   options->tree            = PETSC_FALSE;
155c4762a1bSJed Brown   options->treeCell        = 0;
156c4762a1bSJed Brown   options->testFEjacobian  = PETSC_FALSE;
157c4762a1bSJed Brown   options->testFVgrad      = PETSC_FALSE;
158c4762a1bSJed Brown   options->testInjector    = PETSC_FALSE;
159c4762a1bSJed Brown   options->constants[0]    = 1.0;
160c4762a1bSJed Brown   options->constants[1]    = 2.0;
161c4762a1bSJed Brown   options->constants[2]    = 3.0;
162c4762a1bSJed Brown 
163d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Projection Test Options", "DMPlex");
1649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-use_da", "Flag for DMDA mesh", "ex3.c", options->useDA, &options->useDA, NULL));
1659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-shear_coords", "Transform coordinates with a shear", "ex3.c", options->shearCoords, &options->shearCoords, NULL));
1669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-non_affine_coords", "Transform coordinates with a non-affine transform", "ex3.c", options->nonaffineCoords, &options->nonaffineCoords, NULL));
1679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-qorder", "The quadrature order", "ex3.c", options->qorder, &options->qorder, NULL, 0));
1689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-num_comp", "The number of field components", "ex3.c", options->numComponents, &options->numComponents, NULL, PETSC_DEFAULT));
1699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-porder", "The order of polynomials to test", "ex3.c", options->porder, &options->porder, NULL, 0));
1709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-convergence", "Check the convergence rate", "ex3.c", options->convergence, &options->convergence, NULL));
1719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-conv_refine", "Use refinement for the convergence rate", "ex3.c", options->convRefine, &options->convRefine, NULL));
1729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-constraints", "Test local constraints (serial only)", "ex3.c", options->constraints, &options->constraints, NULL));
1739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-tree", "Test tree routines", "ex3.c", options->tree, &options->tree, NULL));
1749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-tree_cell", "cell to refine in tree test", "ex3.c", options->treeCell, &options->treeCell, NULL, 0));
1759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_fe_jacobian", "Test finite element Jacobian assembly", "ex3.c", options->testFEjacobian, &options->testFEjacobian, NULL));
1769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_fv_grad", "Test finite volume gradient reconstruction", "ex3.c", options->testFVgrad, &options->testFVgrad, NULL));
1779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_injector", "Test finite element injection", "ex3.c", options->testInjector, &options->testInjector, NULL));
1789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRealArray("-constants", "Set the constant values", "ex3.c", options->constants, &n, NULL));
179d0609cedSBarry Smith   PetscOptionsEnd();
180c4762a1bSJed Brown 
1813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
182c4762a1bSJed Brown }
183c4762a1bSJed Brown 
184d71ae5a4SJacob Faibussowitsch static PetscErrorCode TransformCoordinates(DM dm, AppCtx *user)
185d71ae5a4SJacob Faibussowitsch {
186c4762a1bSJed Brown   PetscSection coordSection;
187c4762a1bSJed Brown   Vec          coordinates;
188c4762a1bSJed Brown   PetscScalar *coords;
189c4762a1bSJed Brown   PetscInt     vStart, vEnd, v;
190c4762a1bSJed Brown 
191c4762a1bSJed Brown   PetscFunctionBeginUser;
192c4762a1bSJed Brown   if (user->nonaffineCoords) {
193c4762a1bSJed Brown     /* x' = r^(1/p) (x/r), y' = r^(1/p) (y/r), z' = z */
1949566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateSection(dm, &coordSection));
1959566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
1969566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1979566063dSJacob Faibussowitsch     PetscCall(VecGetArray(coordinates, &coords));
198c4762a1bSJed Brown     for (v = vStart; v < vEnd; ++v) {
199c4762a1bSJed Brown       PetscInt  dof, off;
200c4762a1bSJed Brown       PetscReal p = 4.0, r;
201c4762a1bSJed Brown 
2029566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(coordSection, v, &dof));
2039566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
204c4762a1bSJed Brown       switch (dof) {
205c4762a1bSJed Brown       case 2:
206c4762a1bSJed Brown         r               = PetscSqr(PetscRealPart(coords[off + 0])) + PetscSqr(PetscRealPart(coords[off + 1]));
207c4762a1bSJed Brown         coords[off + 0] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p) / (2 * p)) * coords[off + 0];
208c4762a1bSJed Brown         coords[off + 1] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p) / (2 * p)) * coords[off + 1];
209c4762a1bSJed Brown         break;
210c4762a1bSJed Brown       case 3:
211c4762a1bSJed Brown         r               = PetscSqr(PetscRealPart(coords[off + 0])) + PetscSqr(PetscRealPart(coords[off + 1]));
212c4762a1bSJed Brown         coords[off + 0] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p) / (2 * p)) * coords[off + 0];
213c4762a1bSJed Brown         coords[off + 1] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p) / (2 * p)) * coords[off + 1];
214c4762a1bSJed Brown         coords[off + 2] = coords[off + 2];
215c4762a1bSJed Brown         break;
216c4762a1bSJed Brown       }
217c4762a1bSJed Brown     }
2189566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(coordinates, &coords));
219c4762a1bSJed Brown   }
220c4762a1bSJed Brown   if (user->shearCoords) {
221c4762a1bSJed Brown     /* x' = x + m y + m z, y' = y + m z,  z' = z */
2229566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateSection(dm, &coordSection));
2239566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
2249566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
2259566063dSJacob Faibussowitsch     PetscCall(VecGetArray(coordinates, &coords));
226c4762a1bSJed Brown     for (v = vStart; v < vEnd; ++v) {
227c4762a1bSJed Brown       PetscInt  dof, off;
228c4762a1bSJed Brown       PetscReal m = 1.0;
229c4762a1bSJed Brown 
2309566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(coordSection, v, &dof));
2319566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
232c4762a1bSJed Brown       switch (dof) {
233c4762a1bSJed Brown       case 2:
234c4762a1bSJed Brown         coords[off + 0] = coords[off + 0] + m * coords[off + 1];
235c4762a1bSJed Brown         coords[off + 1] = coords[off + 1];
236c4762a1bSJed Brown         break;
237c4762a1bSJed Brown       case 3:
238c4762a1bSJed Brown         coords[off + 0] = coords[off + 0] + m * coords[off + 1] + m * coords[off + 2];
239c4762a1bSJed Brown         coords[off + 1] = coords[off + 1] + m * coords[off + 2];
240c4762a1bSJed Brown         coords[off + 2] = coords[off + 2];
241c4762a1bSJed Brown         break;
242c4762a1bSJed Brown       }
243c4762a1bSJed Brown     }
2449566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(coordinates, &coords));
245c4762a1bSJed Brown   }
2463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
247c4762a1bSJed Brown }
248c4762a1bSJed Brown 
249d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
250d71ae5a4SJacob Faibussowitsch {
25130602db0SMatthew G. Knepley   PetscInt  dim = 2;
25230602db0SMatthew G. Knepley   PetscBool simplex;
253c4762a1bSJed Brown 
254c4762a1bSJed Brown   PetscFunctionBeginUser;
25530602db0SMatthew G. Knepley   if (user->useDA) {
25630602db0SMatthew G. Knepley     switch (dim) {
257c4762a1bSJed Brown     case 2:
2589566063dSJacob Faibussowitsch       PetscCall(DMDACreate2d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 2, 2, PETSC_DETERMINE, PETSC_DETERMINE, 1, 1, NULL, NULL, dm));
2599566063dSJacob Faibussowitsch       PetscCall(DMSetFromOptions(*dm));
2609566063dSJacob Faibussowitsch       PetscCall(DMSetUp(*dm));
2619566063dSJacob Faibussowitsch       PetscCall(DMDASetVertexCoordinates(*dm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
262c4762a1bSJed Brown       break;
263d71ae5a4SJacob Faibussowitsch     default:
264d71ae5a4SJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot create structured mesh of dimension %" PetscInt_FMT, dim);
265c4762a1bSJed Brown     }
2669566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)*dm, "Hexahedral Mesh"));
26730602db0SMatthew G. Knepley   } else {
2689566063dSJacob Faibussowitsch     PetscCall(DMCreate(comm, dm));
2699566063dSJacob Faibussowitsch     PetscCall(DMSetType(*dm, DMPLEX));
2709566063dSJacob Faibussowitsch     PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
2719566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(*dm));
272c4762a1bSJed Brown 
2739566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(*dm, &dim));
2749566063dSJacob Faibussowitsch     PetscCall(DMPlexIsSimplex(*dm, &simplex));
2759566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&simplex, 1, MPIU_BOOL, 0, comm));
276c4762a1bSJed Brown     if (user->tree) {
27730602db0SMatthew G. Knepley       DM refTree, ncdm = NULL;
278c4762a1bSJed Brown 
2799566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateDefaultReferenceTree(comm, dim, simplex, &refTree));
2809566063dSJacob Faibussowitsch       PetscCall(DMViewFromOptions(refTree, NULL, "-reftree_dm_view"));
2819566063dSJacob Faibussowitsch       PetscCall(DMPlexSetReferenceTree(*dm, refTree));
2829566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&refTree));
2839566063dSJacob Faibussowitsch       PetscCall(DMPlexTreeRefineCell(*dm, user->treeCell, &ncdm));
284c4762a1bSJed Brown       if (ncdm) {
2859566063dSJacob Faibussowitsch         PetscCall(DMDestroy(dm));
286c4762a1bSJed Brown         *dm = ncdm;
2879566063dSJacob Faibussowitsch         PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_FALSE));
288c4762a1bSJed Brown       }
2899566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, "tree_"));
2909566063dSJacob Faibussowitsch       PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
2919566063dSJacob Faibussowitsch       PetscCall(DMSetFromOptions(*dm));
2929566063dSJacob Faibussowitsch       PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
293c4762a1bSJed Brown     } else {
2949566063dSJacob Faibussowitsch       PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_TRUE));
295c4762a1bSJed Brown     }
2969566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, "dist_"));
2979566063dSJacob Faibussowitsch     PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
2989566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(*dm));
2999566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, NULL));
3009566063dSJacob Faibussowitsch     if (simplex) PetscCall(PetscObjectSetName((PetscObject)*dm, "Simplicial Mesh"));
3019566063dSJacob Faibussowitsch     else PetscCall(PetscObjectSetName((PetscObject)*dm, "Hexahedral Mesh"));
302c4762a1bSJed Brown   }
3039566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
3049566063dSJacob Faibussowitsch   PetscCall(TransformCoordinates(*dm, user));
3059566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
3063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
307c4762a1bSJed Brown }
308c4762a1bSJed Brown 
309d71ae5a4SJacob 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[])
310d71ae5a4SJacob Faibussowitsch {
311c4762a1bSJed Brown   PetscInt d, e;
312ad540459SPierre Jolivet   for (d = 0, e = 0; d < dim; d++, e += dim + 1) g0[e] = 1.;
313c4762a1bSJed Brown }
314c4762a1bSJed Brown 
315c4762a1bSJed Brown /* < \nabla v, 1/2(\nabla u + {\nabla u}^T) > */
316d71ae5a4SJacob 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[])
317d71ae5a4SJacob Faibussowitsch {
318c4762a1bSJed Brown   PetscInt compI, compJ, d, e;
319c4762a1bSJed Brown 
320c4762a1bSJed Brown   for (compI = 0; compI < dim; ++compI) {
321c4762a1bSJed Brown     for (compJ = 0; compJ < dim; ++compJ) {
322c4762a1bSJed Brown       for (d = 0; d < dim; ++d) {
323c4762a1bSJed Brown         for (e = 0; e < dim; e++) {
324c4762a1bSJed Brown           if (d == e && d == compI && d == compJ) {
325c4762a1bSJed Brown             C[((compI * dim + compJ) * dim + d) * dim + e] = 1.0;
326c4762a1bSJed Brown           } else if ((d == compJ && e == compI) || (d == e && compI == compJ)) {
327c4762a1bSJed Brown             C[((compI * dim + compJ) * dim + d) * dim + e] = 0.5;
328c4762a1bSJed Brown           } else {
329c4762a1bSJed Brown             C[((compI * dim + compJ) * dim + d) * dim + e] = 0.0;
330c4762a1bSJed Brown           }
331c4762a1bSJed Brown         }
332c4762a1bSJed Brown       }
333c4762a1bSJed Brown     }
334c4762a1bSJed Brown   }
335c4762a1bSJed Brown }
336c4762a1bSJed Brown 
337d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupSection(DM dm, AppCtx *user)
338d71ae5a4SJacob Faibussowitsch {
339c4762a1bSJed Brown   PetscFunctionBeginUser;
34030602db0SMatthew G. Knepley   if (user->constraints) {
341c4762a1bSJed Brown     /* test local constraints */
342c4762a1bSJed Brown     DM           coordDM;
343c4762a1bSJed Brown     PetscInt     fStart, fEnd, f, vStart, vEnd, v;
344c4762a1bSJed Brown     PetscInt     edgesx = 2, vertsx;
345c4762a1bSJed Brown     PetscInt     edgesy = 2, vertsy;
346c4762a1bSJed Brown     PetscMPIInt  size;
347c4762a1bSJed Brown     PetscInt     numConst;
348c4762a1bSJed Brown     PetscSection aSec;
349c4762a1bSJed Brown     PetscInt    *anchors;
350c4762a1bSJed Brown     PetscInt     offset;
351c4762a1bSJed Brown     IS           aIS;
352c4762a1bSJed Brown     MPI_Comm     comm = PetscObjectComm((PetscObject)dm);
353c4762a1bSJed Brown 
3549566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(comm, &size));
35508401ef6SPierre Jolivet     PetscCheck(size <= 1, comm, PETSC_ERR_SUP, "Local constraint test can only be performed in serial");
356c4762a1bSJed Brown 
357c4762a1bSJed Brown     /* we are going to test constraints by using them to enforce periodicity
358c4762a1bSJed Brown      * in one direction, and comparing to the existing method of enforcing
359c4762a1bSJed Brown      * periodicity */
360c4762a1bSJed Brown 
361c4762a1bSJed Brown     /* first create the coordinate section so that it does not clone the
362c4762a1bSJed Brown      * constraints */
3639566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dm, &coordDM));
364c4762a1bSJed Brown 
365c4762a1bSJed Brown     /* create the constrained-to-anchor section */
3669566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
3679566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 1, &fStart, &fEnd));
3689566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreate(PETSC_COMM_SELF, &aSec));
3699566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetChart(aSec, PetscMin(fStart, vStart), PetscMax(fEnd, vEnd)));
370c4762a1bSJed Brown 
371c4762a1bSJed Brown     /* define the constraints */
3729566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetInt(NULL, NULL, "-da_grid_x", &edgesx, NULL));
3739566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetInt(NULL, NULL, "-da_grid_y", &edgesy, NULL));
374c4762a1bSJed Brown     vertsx   = edgesx + 1;
375c4762a1bSJed Brown     vertsy   = edgesy + 1;
376c4762a1bSJed Brown     numConst = vertsy + edgesy;
3779566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numConst, &anchors));
378c4762a1bSJed Brown     offset = 0;
379c4762a1bSJed Brown     for (v = vStart + edgesx; v < vEnd; v += vertsx) {
3809566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetDof(aSec, v, 1));
381c4762a1bSJed Brown       anchors[offset++] = v - edgesx;
382c4762a1bSJed Brown     }
383c4762a1bSJed Brown     for (f = fStart + edgesx * vertsy + edgesx * edgesy; f < fEnd; f++) {
3849566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetDof(aSec, f, 1));
385c4762a1bSJed Brown       anchors[offset++] = f - edgesx * edgesy;
386c4762a1bSJed Brown     }
3879566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetUp(aSec));
3889566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numConst, anchors, PETSC_OWN_POINTER, &aIS));
389c4762a1bSJed Brown 
3909566063dSJacob Faibussowitsch     PetscCall(DMPlexSetAnchors(dm, aSec, aIS));
3919566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&aSec));
3929566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&aIS));
393c4762a1bSJed Brown   }
3949566063dSJacob Faibussowitsch   PetscCall(DMSetNumFields(dm, 1));
3959566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)user->fe));
3969566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
39730602db0SMatthew G. Knepley   if (user->constraints) {
398c4762a1bSJed Brown     /* test getting local constraint matrix that matches section */
399c4762a1bSJed Brown     PetscSection aSec;
400c4762a1bSJed Brown     IS           aIS;
401c4762a1bSJed Brown 
4029566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAnchors(dm, &aSec, &aIS));
403c4762a1bSJed Brown     if (aSec) {
404c4762a1bSJed Brown       PetscDS         ds;
405c4762a1bSJed Brown       PetscSection    cSec, section;
406c4762a1bSJed Brown       PetscInt        cStart, cEnd, c, numComp;
407c4762a1bSJed Brown       Mat             cMat, mass;
408c4762a1bSJed Brown       Vec             local;
409c4762a1bSJed Brown       const PetscInt *anchors;
410c4762a1bSJed Brown 
4119566063dSJacob Faibussowitsch       PetscCall(DMGetLocalSection(dm, &section));
412c4762a1bSJed Brown       /* this creates the matrix and preallocates the matrix structure: we
413c4762a1bSJed Brown        * just have to fill in the values */
4149566063dSJacob Faibussowitsch       PetscCall(DMGetDefaultConstraints(dm, &cSec, &cMat, NULL));
4159566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetChart(cSec, &cStart, &cEnd));
4169566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(aIS, &anchors));
4179566063dSJacob Faibussowitsch       PetscCall(PetscFEGetNumComponents(user->fe, &numComp));
418c4762a1bSJed Brown       for (c = cStart; c < cEnd; c++) {
419c4762a1bSJed Brown         PetscInt cDof;
420c4762a1bSJed Brown 
421c4762a1bSJed Brown         /* is this point constrained? (does it have an anchor?) */
4229566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(aSec, c, &cDof));
423c4762a1bSJed Brown         if (cDof) {
424c4762a1bSJed Brown           PetscInt cOff, a, aDof, aOff, j;
42563a3b9bcSJacob Faibussowitsch           PetscCheck(cDof == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Found %" PetscInt_FMT " anchor points: should be just one", cDof);
426c4762a1bSJed Brown 
427c4762a1bSJed Brown           /* find the anchor point */
4289566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetOffset(aSec, c, &cOff));
429c4762a1bSJed Brown           a = anchors[cOff];
430c4762a1bSJed Brown 
431c4762a1bSJed Brown           /* find the constrained dofs (row in constraint matrix) */
4329566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetDof(cSec, c, &cDof));
4339566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetOffset(cSec, c, &cOff));
434c4762a1bSJed Brown 
435c4762a1bSJed Brown           /* find the anchor dofs (column in constraint matrix) */
4369566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetDof(section, a, &aDof));
4379566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetOffset(section, a, &aOff));
438c4762a1bSJed Brown 
43963a3b9bcSJacob Faibussowitsch           PetscCheck(cDof == aDof, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point and anchor have different number of dofs: %" PetscInt_FMT ", %" PetscInt_FMT, cDof, aDof);
44063a3b9bcSJacob Faibussowitsch           PetscCheck(cDof % numComp == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point dofs not divisible by field components: %" PetscInt_FMT ", %" PetscInt_FMT, cDof, numComp);
441c4762a1bSJed Brown 
442c4762a1bSJed Brown           /* put in a simple equality constraint */
44348a46eb9SPierre Jolivet           for (j = 0; j < cDof; j++) PetscCall(MatSetValue(cMat, cOff + j, aOff + j, 1., INSERT_VALUES));
444c4762a1bSJed Brown         }
445c4762a1bSJed Brown       }
4469566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(cMat, MAT_FINAL_ASSEMBLY));
4479566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(cMat, MAT_FINAL_ASSEMBLY));
4489566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(aIS, &anchors));
449c4762a1bSJed Brown 
450c4762a1bSJed Brown       /* Now that we have constructed the constraint matrix, any FE matrix
451c4762a1bSJed Brown        * that we construct will apply the constraints during construction */
452c4762a1bSJed Brown 
4539566063dSJacob Faibussowitsch       PetscCall(DMCreateMatrix(dm, &mass));
454c4762a1bSJed Brown       /* get a dummy local variable to serve as the solution */
4559566063dSJacob Faibussowitsch       PetscCall(DMGetLocalVector(dm, &local));
4569566063dSJacob Faibussowitsch       PetscCall(DMGetDS(dm, &ds));
457c4762a1bSJed Brown       /* set the jacobian to be the mass matrix */
4589566063dSJacob Faibussowitsch       PetscCall(PetscDSSetJacobian(ds, 0, 0, simple_mass, NULL, NULL, NULL));
459c4762a1bSJed Brown       /* build the mass matrix */
4609566063dSJacob Faibussowitsch       PetscCall(DMPlexSNESComputeJacobianFEM(dm, local, mass, mass, NULL));
4619566063dSJacob Faibussowitsch       PetscCall(MatView(mass, PETSC_VIEWER_STDOUT_WORLD));
4629566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mass));
4639566063dSJacob Faibussowitsch       PetscCall(DMRestoreLocalVector(dm, &local));
464c4762a1bSJed Brown     }
465c4762a1bSJed Brown   }
4663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
467c4762a1bSJed Brown }
468c4762a1bSJed Brown 
469d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestFEJacobian(DM dm, AppCtx *user)
470d71ae5a4SJacob Faibussowitsch {
471c4762a1bSJed Brown   PetscFunctionBeginUser;
47230602db0SMatthew G. Knepley   if (!user->useDA) {
473c4762a1bSJed Brown     Vec          local;
474c4762a1bSJed Brown     const Vec   *vecs;
475c4762a1bSJed Brown     Mat          E;
476c4762a1bSJed Brown     MatNullSpace sp;
477c4762a1bSJed Brown     PetscBool    isNullSpace, hasConst;
47830602db0SMatthew G. Knepley     PetscInt     dim, n, i;
479c4762a1bSJed Brown     Vec          res = NULL, localX, localRes;
480c4762a1bSJed Brown     PetscDS      ds;
481c4762a1bSJed Brown 
4829566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
48363a3b9bcSJacob 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);
4849566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dm, &ds));
4859566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, symmetric_gradient_inner_product));
4869566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix(dm, &E));
4879566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &local));
4889566063dSJacob Faibussowitsch     PetscCall(DMPlexSNESComputeJacobianFEM(dm, local, E, E, NULL));
4899566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateRigidBody(dm, 0, &sp));
4909566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceGetVecs(sp, &hasConst, &n, &vecs));
4919566063dSJacob Faibussowitsch     if (n) PetscCall(VecDuplicate(vecs[0], &res));
4929566063dSJacob Faibussowitsch     PetscCall(DMCreateLocalVector(dm, &localX));
4939566063dSJacob Faibussowitsch     PetscCall(DMCreateLocalVector(dm, &localRes));
494c4762a1bSJed Brown     for (i = 0; i < n; i++) { /* also test via matrix-free Jacobian application */
495c4762a1bSJed Brown       PetscReal resNorm;
496c4762a1bSJed Brown 
4979566063dSJacob Faibussowitsch       PetscCall(VecSet(localRes, 0.));
4989566063dSJacob Faibussowitsch       PetscCall(VecSet(localX, 0.));
4999566063dSJacob Faibussowitsch       PetscCall(VecSet(local, 0.));
5009566063dSJacob Faibussowitsch       PetscCall(VecSet(res, 0.));
5019566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalBegin(dm, vecs[i], INSERT_VALUES, localX));
5029566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalEnd(dm, vecs[i], INSERT_VALUES, localX));
5039566063dSJacob Faibussowitsch       PetscCall(DMSNESComputeJacobianAction(dm, local, localX, localRes, NULL));
5049566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalBegin(dm, localRes, ADD_VALUES, res));
5059566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalEnd(dm, localRes, ADD_VALUES, res));
5069566063dSJacob Faibussowitsch       PetscCall(VecNorm(res, NORM_2, &resNorm));
50748a46eb9SPierre Jolivet       if (resNorm > PETSC_SMALL) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Symmetric gradient action null space vector %" PetscInt_FMT " residual: %E\n", i, (double)resNorm));
508c4762a1bSJed Brown     }
5099566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&localRes));
5109566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&localX));
5119566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&res));
5129566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceTest(sp, E, &isNullSpace));
513c4762a1bSJed Brown     if (isNullSpace) {
5149566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Symmetric gradient null space: PASS\n"));
515c4762a1bSJed Brown     } else {
5169566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Symmetric gradient null space: FAIL\n"));
517c4762a1bSJed Brown     }
5189566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceDestroy(&sp));
5199566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&E));
5209566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &local));
521c4762a1bSJed Brown   }
5223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
523c4762a1bSJed Brown }
524c4762a1bSJed Brown 
525d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestInjector(DM dm, AppCtx *user)
526d71ae5a4SJacob Faibussowitsch {
527c4762a1bSJed Brown   DM          refTree;
528c4762a1bSJed Brown   PetscMPIInt rank;
529c4762a1bSJed Brown 
530c4762a1bSJed Brown   PetscFunctionBegin;
5319566063dSJacob Faibussowitsch   PetscCall(DMPlexGetReferenceTree(dm, &refTree));
532c4762a1bSJed Brown   if (refTree) {
533c4762a1bSJed Brown     Mat inj;
534c4762a1bSJed Brown 
5359566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeInjectorReferenceTree(refTree, &inj));
5369566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)inj, "Reference Tree Injector"));
5379566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
53848a46eb9SPierre Jolivet     if (rank == 0) PetscCall(MatView(inj, PETSC_VIEWER_STDOUT_SELF));
5399566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&inj));
540c4762a1bSJed Brown   }
5413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
542c4762a1bSJed Brown }
543c4762a1bSJed Brown 
544d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestFVGrad(DM dm, AppCtx *user)
545d71ae5a4SJacob Faibussowitsch {
546c4762a1bSJed Brown   MPI_Comm           comm;
547c4762a1bSJed Brown   DM                 dmRedist, dmfv, dmgrad, dmCell, refTree;
548c4762a1bSJed Brown   PetscFV            fv;
54930602db0SMatthew G. Knepley   PetscInt           dim, nvecs, v, cStart, cEnd, cEndInterior;
550c4762a1bSJed Brown   PetscMPIInt        size;
551c4762a1bSJed Brown   Vec                cellgeom, grad, locGrad;
552c4762a1bSJed Brown   const PetscScalar *cgeom;
553c4762a1bSJed Brown   PetscReal          allVecMaxDiff = 0., fvTol = 100. * PETSC_MACHINE_EPSILON;
554c4762a1bSJed Brown 
555c4762a1bSJed Brown   PetscFunctionBeginUser;
556c4762a1bSJed Brown   comm = PetscObjectComm((PetscObject)dm);
5579566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
558c4762a1bSJed Brown   /* duplicate DM, give dup. a FV discretization */
5599566063dSJacob Faibussowitsch   PetscCall(DMSetBasicAdjacency(dm, PETSC_TRUE, PETSC_FALSE));
5609566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
561c4762a1bSJed Brown   dmRedist = NULL;
56248a46eb9SPierre Jolivet   if (size > 1) PetscCall(DMPlexDistributeOverlap(dm, 1, NULL, &dmRedist));
563c4762a1bSJed Brown   if (!dmRedist) {
564c4762a1bSJed Brown     dmRedist = dm;
5659566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)dmRedist));
566c4762a1bSJed Brown   }
5679566063dSJacob Faibussowitsch   PetscCall(PetscFVCreate(comm, &fv));
5689566063dSJacob Faibussowitsch   PetscCall(PetscFVSetType(fv, PETSCFVLEASTSQUARES));
5699566063dSJacob Faibussowitsch   PetscCall(PetscFVSetNumComponents(fv, user->numComponents));
5709566063dSJacob Faibussowitsch   PetscCall(PetscFVSetSpatialDimension(fv, dim));
5719566063dSJacob Faibussowitsch   PetscCall(PetscFVSetFromOptions(fv));
5729566063dSJacob Faibussowitsch   PetscCall(PetscFVSetUp(fv));
573*4446c3cdSksagiyam   {
574*4446c3cdSksagiyam     PetscSF pointSF;
575*4446c3cdSksagiyam     DMLabel label;
576*4446c3cdSksagiyam 
577*4446c3cdSksagiyam     PetscCall(DMCreateLabel(dmRedist, "Face Sets"));
578*4446c3cdSksagiyam     PetscCall(DMGetLabel(dmRedist, "Face Sets", &label));
579*4446c3cdSksagiyam     PetscCall(DMGetPointSF(dmRedist, &pointSF));
580*4446c3cdSksagiyam     PetscCall(PetscObjectReference((PetscObject)pointSF));
581*4446c3cdSksagiyam     PetscCall(DMSetPointSF(dmRedist, NULL));
582*4446c3cdSksagiyam     PetscCall(DMPlexMarkBoundaryFaces(dmRedist, 1, label));
583*4446c3cdSksagiyam     PetscCall(DMSetPointSF(dmRedist, pointSF));
584*4446c3cdSksagiyam     PetscCall(PetscSFDestroy(&pointSF));
585*4446c3cdSksagiyam   }
5869566063dSJacob Faibussowitsch   PetscCall(DMPlexConstructGhostCells(dmRedist, NULL, NULL, &dmfv));
5879566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmRedist));
5889566063dSJacob Faibussowitsch   PetscCall(DMSetNumFields(dmfv, 1));
5899566063dSJacob Faibussowitsch   PetscCall(DMSetField(dmfv, 0, NULL, (PetscObject)fv));
5909566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dmfv));
5919566063dSJacob Faibussowitsch   PetscCall(DMPlexGetReferenceTree(dm, &refTree));
5929566063dSJacob Faibussowitsch   if (refTree) PetscCall(DMCopyDisc(dmfv, refTree));
5939566063dSJacob Faibussowitsch   PetscCall(DMPlexGetGradientDM(dmfv, fv, &dmgrad));
5949566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dmfv, 0, &cStart, &cEnd));
59530602db0SMatthew G. Knepley   nvecs = dim * (dim + 1) / 2;
5969566063dSJacob Faibussowitsch   PetscCall(DMPlexGetGeometryFVM(dmfv, NULL, &cellgeom, NULL));
5979566063dSJacob Faibussowitsch   PetscCall(VecGetDM(cellgeom, &dmCell));
5989566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(cellgeom, &cgeom));
5999566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dmgrad, &grad));
6009566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dmgrad, &locGrad));
6012827ebadSStefano Zampini   PetscCall(DMPlexGetCellTypeStratum(dmgrad, DM_POLYTOPE_FV_GHOST, &cEndInterior, NULL));
602c4762a1bSJed Brown   cEndInterior = (cEndInterior < 0) ? cEnd : cEndInterior;
603c4762a1bSJed Brown   for (v = 0; v < nvecs; v++) {
604c4762a1bSJed Brown     Vec                locX;
605c4762a1bSJed Brown     PetscInt           c;
606c4762a1bSJed Brown     PetscScalar        trueGrad[3][3] = {{0.}};
607c4762a1bSJed Brown     const PetscScalar *gradArray;
608c4762a1bSJed Brown     PetscReal          maxDiff, maxDiffGlob;
609c4762a1bSJed Brown 
6109566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dmfv, &locX));
611c4762a1bSJed Brown     /* get the local projection of the rigid body mode */
612c4762a1bSJed Brown     for (c = cStart; c < cEnd; c++) {
613c4762a1bSJed Brown       PetscFVCellGeom *cg;
614c4762a1bSJed Brown       PetscScalar      cx[3] = {0., 0., 0.};
615c4762a1bSJed Brown 
6169566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg));
61730602db0SMatthew G. Knepley       if (v < dim) {
618c4762a1bSJed Brown         cx[v] = 1.;
619c4762a1bSJed Brown       } else {
62030602db0SMatthew G. Knepley         PetscInt w = v - dim;
621c4762a1bSJed Brown 
62230602db0SMatthew G. Knepley         cx[(w + 1) % dim] = cg->centroid[(w + 2) % dim];
62330602db0SMatthew G. Knepley         cx[(w + 2) % dim] = -cg->centroid[(w + 1) % dim];
624c4762a1bSJed Brown       }
6259566063dSJacob Faibussowitsch       PetscCall(DMPlexVecSetClosure(dmfv, NULL, locX, c, cx, INSERT_ALL_VALUES));
626c4762a1bSJed Brown     }
627c4762a1bSJed Brown     /* TODO: this isn't in any header */
6289566063dSJacob Faibussowitsch     PetscCall(DMPlexReconstructGradientsFVM(dmfv, locX, grad));
6299566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(dmgrad, grad, INSERT_VALUES, locGrad));
6309566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(dmgrad, grad, INSERT_VALUES, locGrad));
6319566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(locGrad, &gradArray));
632c4762a1bSJed Brown     /* compare computed gradient to exact gradient */
63330602db0SMatthew G. Knepley     if (v >= dim) {
63430602db0SMatthew G. Knepley       PetscInt w = v - dim;
635c4762a1bSJed Brown 
63630602db0SMatthew G. Knepley       trueGrad[(w + 1) % dim][(w + 2) % dim] = 1.;
63730602db0SMatthew G. Knepley       trueGrad[(w + 2) % dim][(w + 1) % dim] = -1.;
638c4762a1bSJed Brown     }
639c4762a1bSJed Brown     maxDiff = 0.;
640c4762a1bSJed Brown     for (c = cStart; c < cEndInterior; c++) {
641c4762a1bSJed Brown       PetscScalar *compGrad;
642c4762a1bSJed Brown       PetscInt     i, j, k;
643c4762a1bSJed Brown       PetscReal    FrobDiff = 0.;
644c4762a1bSJed Brown 
6459566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRead(dmgrad, c, gradArray, &compGrad));
646c4762a1bSJed Brown 
64730602db0SMatthew G. Knepley       for (i = 0, k = 0; i < dim; i++) {
64830602db0SMatthew G. Knepley         for (j = 0; j < dim; j++, k++) {
649c4762a1bSJed Brown           PetscScalar diff = compGrad[k] - trueGrad[i][j];
650c4762a1bSJed Brown           FrobDiff += PetscRealPart(diff * PetscConj(diff));
651c4762a1bSJed Brown         }
652c4762a1bSJed Brown       }
653c4762a1bSJed Brown       FrobDiff = PetscSqrtReal(FrobDiff);
654c4762a1bSJed Brown       maxDiff  = PetscMax(maxDiff, FrobDiff);
655c4762a1bSJed Brown     }
656712fec58SPierre Jolivet     PetscCall(MPIU_Allreduce(&maxDiff, &maxDiffGlob, 1, MPIU_REAL, MPIU_MAX, comm));
657c4762a1bSJed Brown     allVecMaxDiff = PetscMax(allVecMaxDiff, maxDiffGlob);
6589566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(locGrad, &gradArray));
6599566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dmfv, &locX));
660c4762a1bSJed Brown   }
661c4762a1bSJed Brown   if (allVecMaxDiff < fvTol) {
6629566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Finite volume gradient reconstruction: PASS\n"));
663c4762a1bSJed Brown   } else {
66463a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Finite volume gradient reconstruction: FAIL at tolerance %g with max difference %g\n", (double)fvTol, (double)allVecMaxDiff));
665c4762a1bSJed Brown   }
6669566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dmgrad, &locGrad));
6679566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dmgrad, &grad));
6689566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(cellgeom, &cgeom));
6699566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmfv));
6709566063dSJacob Faibussowitsch   PetscCall(PetscFVDestroy(&fv));
6713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
672c4762a1bSJed Brown }
673c4762a1bSJed Brown 
674d71ae5a4SJacob 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)
675d71ae5a4SJacob Faibussowitsch {
676c4762a1bSJed Brown   Vec       u;
677c4762a1bSJed Brown   PetscReal n[3] = {1.0, 1.0, 1.0};
678c4762a1bSJed Brown 
679c4762a1bSJed Brown   PetscFunctionBeginUser;
6809566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &u));
681c4762a1bSJed Brown   /* Project function into FE function space */
6829566063dSJacob Faibussowitsch   PetscCall(DMProjectFunction(dm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, u));
6839566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-projection_view"));
684c4762a1bSJed Brown   /* Compare approximation to exact in L_2 */
6859566063dSJacob Faibussowitsch   PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, exactCtxs, u, error));
6869566063dSJacob Faibussowitsch   PetscCall(DMComputeL2GradientDiff(dm, 0.0, exactFuncDers, exactCtxs, u, n, errorDer));
6879566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &u));
6883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
689c4762a1bSJed Brown }
690c4762a1bSJed Brown 
691d71ae5a4SJacob Faibussowitsch static PetscErrorCode CheckFunctions(DM dm, PetscInt order, AppCtx *user)
692d71ae5a4SJacob Faibussowitsch {
693c4762a1bSJed Brown   PetscErrorCode (*exactFuncs[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
694c4762a1bSJed Brown   PetscErrorCode (*exactFuncDers[1])(PetscInt dim, PetscReal time, const PetscReal x[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx);
695c4762a1bSJed Brown   void     *exactCtxs[3];
696c4762a1bSJed Brown   MPI_Comm  comm;
697c4762a1bSJed Brown   PetscReal error, errorDer, tol = PETSC_SMALL;
698c4762a1bSJed Brown 
699c4762a1bSJed Brown   PetscFunctionBeginUser;
700c4762a1bSJed Brown   exactCtxs[0] = user;
701c4762a1bSJed Brown   exactCtxs[1] = user;
702c4762a1bSJed Brown   exactCtxs[2] = user;
7039566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
704c4762a1bSJed Brown   /* Setup functions to approximate */
705c4762a1bSJed Brown   switch (order) {
706c4762a1bSJed Brown   case 0:
707c4762a1bSJed Brown     exactFuncs[0]    = constant;
708c4762a1bSJed Brown     exactFuncDers[0] = constantDer;
709c4762a1bSJed Brown     break;
710c4762a1bSJed Brown   case 1:
711c4762a1bSJed Brown     exactFuncs[0]    = linear;
712c4762a1bSJed Brown     exactFuncDers[0] = linearDer;
713c4762a1bSJed Brown     break;
714c4762a1bSJed Brown   case 2:
715c4762a1bSJed Brown     exactFuncs[0]    = quadratic;
716c4762a1bSJed Brown     exactFuncDers[0] = quadraticDer;
717c4762a1bSJed Brown     break;
718c4762a1bSJed Brown   case 3:
719c4762a1bSJed Brown     exactFuncs[0]    = cubic;
720c4762a1bSJed Brown     exactFuncDers[0] = cubicDer;
721c4762a1bSJed Brown     break;
722d71ae5a4SJacob Faibussowitsch   default:
723d71ae5a4SJacob Faibussowitsch     SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Could not determine functions to test for order %" PetscInt_FMT, order);
724c4762a1bSJed Brown   }
7259566063dSJacob Faibussowitsch   PetscCall(ComputeError(dm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user));
726c4762a1bSJed Brown   /* Report result */
72763a3b9bcSJacob 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));
72863a3b9bcSJacob Faibussowitsch   else PetscCall(PetscPrintf(comm, "Function tests pass for order %" PetscInt_FMT " at tolerance %g\n", order, (double)tol));
72963a3b9bcSJacob 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));
73063a3b9bcSJacob Faibussowitsch   else PetscCall(PetscPrintf(comm, "Function tests pass for order %" PetscInt_FMT " derivatives at tolerance %g\n", order, (double)tol));
7313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
732c4762a1bSJed Brown }
733c4762a1bSJed Brown 
734d71ae5a4SJacob Faibussowitsch static PetscErrorCode CheckInterpolation(DM dm, PetscBool checkRestrict, PetscInt order, AppCtx *user)
735d71ae5a4SJacob Faibussowitsch {
736c4762a1bSJed Brown   PetscErrorCode (*exactFuncs[1])(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
737c4762a1bSJed Brown   PetscErrorCode (*exactFuncDers[1])(PetscInt, PetscReal, const PetscReal x[], const PetscReal n[], PetscInt, PetscScalar *u, void *ctx);
738c4762a1bSJed Brown   PetscReal n[3] = {1.0, 1.0, 1.0};
739c4762a1bSJed Brown   void     *exactCtxs[3];
740c4762a1bSJed Brown   DM        rdm, idm, fdm;
741c4762a1bSJed Brown   Mat       Interp;
742c4762a1bSJed Brown   Vec       iu, fu, scaling;
743c4762a1bSJed Brown   MPI_Comm  comm;
74430602db0SMatthew G. Knepley   PetscInt  dim;
745c4762a1bSJed Brown   PetscReal error, errorDer, tol = PETSC_SMALL;
746c4762a1bSJed Brown 
747c4762a1bSJed Brown   PetscFunctionBeginUser;
748c4762a1bSJed Brown   exactCtxs[0] = user;
749c4762a1bSJed Brown   exactCtxs[1] = user;
750c4762a1bSJed Brown   exactCtxs[2] = user;
7519566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
7529566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
7539566063dSJacob Faibussowitsch   PetscCall(DMRefine(dm, comm, &rdm));
7549566063dSJacob Faibussowitsch   PetscCall(DMSetCoarseDM(rdm, dm));
7559566063dSJacob Faibussowitsch   PetscCall(DMPlexSetRegularRefinement(rdm, user->convRefine));
75630602db0SMatthew G. Knepley   if (user->tree) {
757c4762a1bSJed Brown     DM refTree;
7589566063dSJacob Faibussowitsch     PetscCall(DMPlexGetReferenceTree(dm, &refTree));
7599566063dSJacob Faibussowitsch     PetscCall(DMPlexSetReferenceTree(rdm, refTree));
760c4762a1bSJed Brown   }
7619566063dSJacob Faibussowitsch   if (user->useDA) PetscCall(DMDASetVertexCoordinates(rdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
7629566063dSJacob Faibussowitsch   PetscCall(SetupSection(rdm, user));
763c4762a1bSJed Brown   /* Setup functions to approximate */
764c4762a1bSJed Brown   switch (order) {
765c4762a1bSJed Brown   case 0:
766c4762a1bSJed Brown     exactFuncs[0]    = constant;
767c4762a1bSJed Brown     exactFuncDers[0] = constantDer;
768c4762a1bSJed Brown     break;
769c4762a1bSJed Brown   case 1:
770c4762a1bSJed Brown     exactFuncs[0]    = linear;
771c4762a1bSJed Brown     exactFuncDers[0] = linearDer;
772c4762a1bSJed Brown     break;
773c4762a1bSJed Brown   case 2:
774c4762a1bSJed Brown     exactFuncs[0]    = quadratic;
775c4762a1bSJed Brown     exactFuncDers[0] = quadraticDer;
776c4762a1bSJed Brown     break;
777c4762a1bSJed Brown   case 3:
778c4762a1bSJed Brown     exactFuncs[0]    = cubic;
779c4762a1bSJed Brown     exactFuncDers[0] = cubicDer;
780c4762a1bSJed Brown     break;
781d71ae5a4SJacob Faibussowitsch   default:
782d71ae5a4SJacob Faibussowitsch     SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Could not determine functions to test for dimension %" PetscInt_FMT " order %" PetscInt_FMT, dim, order);
783c4762a1bSJed Brown   }
784c4762a1bSJed Brown   idm = checkRestrict ? rdm : dm;
785c4762a1bSJed Brown   fdm = checkRestrict ? dm : rdm;
7869566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(idm, &iu));
7879566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(fdm, &fu));
7889566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(dm, user));
7899566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(rdm, user));
7909566063dSJacob Faibussowitsch   PetscCall(DMCreateInterpolation(dm, rdm, &Interp, &scaling));
791c4762a1bSJed Brown   /* Project function into initial FE function space */
7929566063dSJacob Faibussowitsch   PetscCall(DMProjectFunction(idm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, iu));
793c4762a1bSJed Brown   /* Interpolate function into final FE function space */
7949371c9d4SSatish Balay   if (checkRestrict) {
7959371c9d4SSatish Balay     PetscCall(MatRestrict(Interp, iu, fu));
7969371c9d4SSatish Balay     PetscCall(VecPointwiseMult(fu, scaling, fu));
7979371c9d4SSatish Balay   } else PetscCall(MatInterpolate(Interp, iu, fu));
798c4762a1bSJed Brown   /* Compare approximation to exact in L_2 */
7999566063dSJacob Faibussowitsch   PetscCall(DMComputeL2Diff(fdm, 0.0, exactFuncs, exactCtxs, fu, &error));
8009566063dSJacob Faibussowitsch   PetscCall(DMComputeL2GradientDiff(fdm, 0.0, exactFuncDers, exactCtxs, fu, n, &errorDer));
801c4762a1bSJed Brown   /* Report result */
80263a3b9bcSJacob 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));
80363a3b9bcSJacob Faibussowitsch   else PetscCall(PetscPrintf(comm, "Interpolation tests pass for order %" PetscInt_FMT " at tolerance %g\n", order, (double)tol));
80463a3b9bcSJacob 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));
80563a3b9bcSJacob Faibussowitsch   else PetscCall(PetscPrintf(comm, "Interpolation tests pass for order %" PetscInt_FMT " derivatives at tolerance %g\n", order, (double)tol));
8069566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(idm, &iu));
8079566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(fdm, &fu));
8089566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Interp));
8099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&scaling));
8109566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&rdm));
8113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
812c4762a1bSJed Brown }
813c4762a1bSJed Brown 
814d71ae5a4SJacob Faibussowitsch static PetscErrorCode CheckConvergence(DM dm, PetscInt Nr, AppCtx *user)
815d71ae5a4SJacob Faibussowitsch {
816c4762a1bSJed Brown   DM odm = dm, rdm = NULL, cdm = NULL;
817c4762a1bSJed Brown   PetscErrorCode (*exactFuncs[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)                         = {trig};
818c4762a1bSJed Brown   PetscErrorCode (*exactFuncDers[1])(PetscInt dim, PetscReal time, const PetscReal x[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) = {trigDer};
819c4762a1bSJed Brown   void     *exactCtxs[3];
820c4762a1bSJed Brown   PetscInt  r, c, cStart, cEnd;
821c4762a1bSJed Brown   PetscReal errorOld, errorDerOld, error, errorDer, rel, len, lenOld;
822c4762a1bSJed Brown   double    p;
823c4762a1bSJed Brown 
824c4762a1bSJed Brown   PetscFunctionBeginUser;
8253ba16761SJacob Faibussowitsch   if (!user->convergence) PetscFunctionReturn(PETSC_SUCCESS);
826c4762a1bSJed Brown   exactCtxs[0] = user;
827c4762a1bSJed Brown   exactCtxs[1] = user;
828c4762a1bSJed Brown   exactCtxs[2] = user;
8299566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)odm));
830c4762a1bSJed Brown   if (!user->convRefine) {
831c4762a1bSJed Brown     for (r = 0; r < Nr; ++r) {
8329566063dSJacob Faibussowitsch       PetscCall(DMRefine(odm, PetscObjectComm((PetscObject)dm), &rdm));
8339566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&odm));
834c4762a1bSJed Brown       odm = rdm;
835c4762a1bSJed Brown     }
8369566063dSJacob Faibussowitsch     PetscCall(SetupSection(odm, user));
837c4762a1bSJed Brown   }
8389566063dSJacob Faibussowitsch   PetscCall(ComputeError(odm, exactFuncs, exactFuncDers, exactCtxs, &errorOld, &errorDerOld, user));
839c4762a1bSJed Brown   if (user->convRefine) {
840c4762a1bSJed Brown     for (r = 0; r < Nr; ++r) {
8419566063dSJacob Faibussowitsch       PetscCall(DMRefine(odm, PetscObjectComm((PetscObject)dm), &rdm));
8429566063dSJacob Faibussowitsch       if (user->useDA) PetscCall(DMDASetVertexCoordinates(rdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
8439566063dSJacob Faibussowitsch       PetscCall(SetupSection(rdm, user));
8449566063dSJacob Faibussowitsch       PetscCall(ComputeError(rdm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user));
845c4762a1bSJed Brown       p = PetscLog2Real(errorOld / error);
84663a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Function   convergence rate at refinement %" PetscInt_FMT ": %.2f\n", r, (double)p));
847c4762a1bSJed Brown       p = PetscLog2Real(errorDerOld / errorDer);
84863a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Derivative convergence rate at refinement %" PetscInt_FMT ": %.2f\n", r, (double)p));
8499566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&odm));
850c4762a1bSJed Brown       odm         = rdm;
851c4762a1bSJed Brown       errorOld    = error;
852c4762a1bSJed Brown       errorDerOld = errorDer;
853c4762a1bSJed Brown     }
854c4762a1bSJed Brown   } else {
8559566063dSJacob Faibussowitsch     /* PetscCall(ComputeLongestEdge(dm, &lenOld)); */
8569566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(odm, 0, &cStart, &cEnd));
857c4762a1bSJed Brown     lenOld = cEnd - cStart;
858c4762a1bSJed Brown     for (c = 0; c < Nr; ++c) {
8599566063dSJacob Faibussowitsch       PetscCall(DMCoarsen(odm, PetscObjectComm((PetscObject)dm), &cdm));
8609566063dSJacob Faibussowitsch       if (user->useDA) PetscCall(DMDASetVertexCoordinates(cdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
8619566063dSJacob Faibussowitsch       PetscCall(SetupSection(cdm, user));
8629566063dSJacob Faibussowitsch       PetscCall(ComputeError(cdm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user));
8639566063dSJacob Faibussowitsch       /* PetscCall(ComputeLongestEdge(cdm, &len)); */
8649566063dSJacob Faibussowitsch       PetscCall(DMPlexGetHeightStratum(cdm, 0, &cStart, &cEnd));
865c4762a1bSJed Brown       len = cEnd - cStart;
866c4762a1bSJed Brown       rel = error / errorOld;
867c4762a1bSJed Brown       p   = PetscLogReal(rel) / PetscLogReal(lenOld / len);
86863a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Function   convergence rate at coarsening %" PetscInt_FMT ": %.2f\n", c, (double)p));
869c4762a1bSJed Brown       rel = errorDer / errorDerOld;
870c4762a1bSJed Brown       p   = PetscLogReal(rel) / PetscLogReal(lenOld / len);
87163a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Derivative convergence rate at coarsening %" PetscInt_FMT ": %.2f\n", c, (double)p));
8729566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&odm));
873c4762a1bSJed Brown       odm         = cdm;
874c4762a1bSJed Brown       errorOld    = error;
875c4762a1bSJed Brown       errorDerOld = errorDer;
876c4762a1bSJed Brown       lenOld      = len;
877c4762a1bSJed Brown     }
878c4762a1bSJed Brown   }
8799566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&odm));
8803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
881c4762a1bSJed Brown }
882c4762a1bSJed Brown 
883d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
884d71ae5a4SJacob Faibussowitsch {
885c4762a1bSJed Brown   DM        dm;
886c4762a1bSJed Brown   AppCtx    user; /* user-defined work context */
88730602db0SMatthew G. Knepley   PetscInt  dim     = 2;
88830602db0SMatthew G. Knepley   PetscBool simplex = PETSC_FALSE;
889c4762a1bSJed Brown 
890327415f7SBarry Smith   PetscFunctionBeginUser;
8919566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
8929566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
8939566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
89430602db0SMatthew G. Knepley   if (!user.useDA) {
8959566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
8969566063dSJacob Faibussowitsch     PetscCall(DMPlexIsSimplex(dm, &simplex));
89730602db0SMatthew G. Knepley   }
8989566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetFromOptions(dm));
89930602db0SMatthew G. Knepley   user.numComponents = user.numComponents < 0 ? dim : user.numComponents;
9009566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_WORLD, dim, user.numComponents, simplex, NULL, user.qorder, &user.fe));
9019566063dSJacob Faibussowitsch   PetscCall(SetupSection(dm, &user));
9029566063dSJacob Faibussowitsch   if (user.testFEjacobian) PetscCall(TestFEJacobian(dm, &user));
9039566063dSJacob Faibussowitsch   if (user.testFVgrad) PetscCall(TestFVGrad(dm, &user));
9049566063dSJacob Faibussowitsch   if (user.testInjector) PetscCall(TestInjector(dm, &user));
9059566063dSJacob Faibussowitsch   PetscCall(CheckFunctions(dm, user.porder, &user));
906c4762a1bSJed Brown   {
907c4762a1bSJed Brown     PetscDualSpace dsp;
908c4762a1bSJed Brown     PetscInt       k;
909c4762a1bSJed Brown 
9109566063dSJacob Faibussowitsch     PetscCall(PetscFEGetDualSpace(user.fe, &dsp));
9119566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDeRahm(dsp, &k));
91230602db0SMatthew G. Knepley     if (dim == 2 && simplex == PETSC_TRUE && user.tree == PETSC_FALSE && k == 0) {
9139566063dSJacob Faibussowitsch       PetscCall(CheckInterpolation(dm, PETSC_FALSE, user.porder, &user));
9149566063dSJacob Faibussowitsch       PetscCall(CheckInterpolation(dm, PETSC_TRUE, user.porder, &user));
915c4762a1bSJed Brown     }
916c4762a1bSJed Brown   }
9179566063dSJacob Faibussowitsch   PetscCall(CheckConvergence(dm, 3, &user));
9189566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&user.fe));
9199566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
9209566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
921b122ec5aSJacob Faibussowitsch   return 0;
922c4762a1bSJed Brown }
923c4762a1bSJed Brown 
924c4762a1bSJed Brown /*TEST
925c4762a1bSJed Brown 
926c4762a1bSJed Brown   test:
927c4762a1bSJed Brown     suffix: 1
928c4762a1bSJed Brown     requires: triangle
929c4762a1bSJed Brown 
930c4762a1bSJed Brown   # 2D P_1 on a triangle
931c4762a1bSJed Brown   test:
932c4762a1bSJed Brown     suffix: p1_2d_0
933c4762a1bSJed Brown     requires: triangle
934c4762a1bSJed Brown     args: -petscspace_degree 1 -qorder 1 -convergence
935c4762a1bSJed Brown   test:
936c4762a1bSJed Brown     suffix: p1_2d_1
937c4762a1bSJed Brown     requires: triangle
938c4762a1bSJed Brown     args: -petscspace_degree 1 -qorder 1 -porder 1
939c4762a1bSJed Brown   test:
940c4762a1bSJed Brown     suffix: p1_2d_2
941c4762a1bSJed Brown     requires: triangle
942c4762a1bSJed Brown     args: -petscspace_degree 1 -qorder 1 -porder 2
943c4762a1bSJed Brown   test:
944c4762a1bSJed Brown     suffix: p1_2d_3
945e788613bSJoe Wallwork     requires: triangle mmg
946c4762a1bSJed Brown     args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0
947c4762a1bSJed Brown   test:
948c4762a1bSJed Brown     suffix: p1_2d_4
949e788613bSJoe Wallwork     requires: triangle mmg
950c4762a1bSJed Brown     args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0
951c4762a1bSJed Brown   test:
952c4762a1bSJed Brown     suffix: p1_2d_5
953e788613bSJoe Wallwork     requires: triangle mmg
954c4762a1bSJed Brown     args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0
955c4762a1bSJed Brown 
956c4762a1bSJed Brown   # 3D P_1 on a tetrahedron
957c4762a1bSJed Brown   test:
958c4762a1bSJed Brown     suffix: p1_3d_0
959c4762a1bSJed Brown     requires: ctetgen
96030602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -convergence
961c4762a1bSJed Brown   test:
962c4762a1bSJed Brown     suffix: p1_3d_1
963c4762a1bSJed Brown     requires: ctetgen
96430602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -porder 1
965c4762a1bSJed Brown   test:
966c4762a1bSJed Brown     suffix: p1_3d_2
967c4762a1bSJed Brown     requires: ctetgen
96830602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -porder 2
969c4762a1bSJed Brown   test:
970c4762a1bSJed Brown     suffix: p1_3d_3
971e788613bSJoe Wallwork     requires: ctetgen mmg
97230602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0
973c4762a1bSJed Brown   test:
974c4762a1bSJed Brown     suffix: p1_3d_4
975e788613bSJoe Wallwork     requires: ctetgen mmg
97630602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0
977c4762a1bSJed Brown   test:
978c4762a1bSJed Brown     suffix: p1_3d_5
979e788613bSJoe Wallwork     requires: ctetgen mmg
98030602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0
981c4762a1bSJed Brown 
982c4762a1bSJed Brown   # 2D P_2 on a triangle
983c4762a1bSJed Brown   test:
984c4762a1bSJed Brown     suffix: p2_2d_0
985c4762a1bSJed Brown     requires: triangle
986c4762a1bSJed Brown     args: -petscspace_degree 2 -qorder 2 -convergence
987c4762a1bSJed Brown   test:
988c4762a1bSJed Brown     suffix: p2_2d_1
989c4762a1bSJed Brown     requires: triangle
990c4762a1bSJed Brown     args: -petscspace_degree 2 -qorder 2 -porder 1
991c4762a1bSJed Brown   test:
992c4762a1bSJed Brown     suffix: p2_2d_2
993c4762a1bSJed Brown     requires: triangle
994c4762a1bSJed Brown     args: -petscspace_degree 2 -qorder 2 -porder 2
995c4762a1bSJed Brown   test:
996c4762a1bSJed Brown     suffix: p2_2d_3
997e788613bSJoe Wallwork     requires: triangle mmg
998c4762a1bSJed Brown     args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -convergence -conv_refine 0
999c4762a1bSJed Brown   test:
1000c4762a1bSJed Brown     suffix: p2_2d_4
1001e788613bSJoe Wallwork     requires: triangle mmg
1002c4762a1bSJed Brown     args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 1 -conv_refine 0
1003c4762a1bSJed Brown   test:
1004c4762a1bSJed Brown     suffix: p2_2d_5
1005e788613bSJoe Wallwork     requires: triangle mmg
1006c4762a1bSJed Brown     args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 2 -conv_refine 0
1007c4762a1bSJed Brown 
1008c4762a1bSJed Brown   # 3D P_2 on a tetrahedron
1009c4762a1bSJed Brown   test:
1010c4762a1bSJed Brown     suffix: p2_3d_0
1011c4762a1bSJed Brown     requires: ctetgen
101230602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -convergence
1013c4762a1bSJed Brown   test:
1014c4762a1bSJed Brown     suffix: p2_3d_1
1015c4762a1bSJed Brown     requires: ctetgen
101630602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -porder 1
1017c4762a1bSJed Brown   test:
1018c4762a1bSJed Brown     suffix: p2_3d_2
1019c4762a1bSJed Brown     requires: ctetgen
102030602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -porder 2
1021c4762a1bSJed Brown   test:
1022c4762a1bSJed Brown     suffix: p2_3d_3
1023e788613bSJoe Wallwork     requires: ctetgen mmg
102430602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -convergence -conv_refine 0
1025c4762a1bSJed Brown   test:
1026c4762a1bSJed Brown     suffix: p2_3d_4
1027e788613bSJoe Wallwork     requires: ctetgen mmg
102830602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 1 -conv_refine 0
1029c4762a1bSJed Brown   test:
1030c4762a1bSJed Brown     suffix: p2_3d_5
1031e788613bSJoe Wallwork     requires: ctetgen mmg
103230602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 2 -conv_refine 0
1033c4762a1bSJed Brown 
1034c4762a1bSJed Brown   # 2D Q_1 on a quadrilaterial DA
1035c4762a1bSJed Brown   test:
1036c4762a1bSJed Brown     suffix: q1_2d_da_0
103799a192c5SJunchao Zhang     requires: broken
103830602db0SMatthew G. Knepley     args: -use_da 1 -petscspace_degree 1 -qorder 1 -convergence
1039c4762a1bSJed Brown   test:
1040c4762a1bSJed Brown     suffix: q1_2d_da_1
104199a192c5SJunchao Zhang     requires: broken
104230602db0SMatthew G. Knepley     args: -use_da 1 -petscspace_degree 1 -qorder 1 -porder 1
1043c4762a1bSJed Brown   test:
1044c4762a1bSJed Brown     suffix: q1_2d_da_2
104599a192c5SJunchao Zhang     requires: broken
104630602db0SMatthew G. Knepley     args: -use_da 1 -petscspace_degree 1 -qorder 1 -porder 2
1047c4762a1bSJed Brown 
1048c4762a1bSJed Brown   # 2D Q_1 on a quadrilaterial Plex
1049c4762a1bSJed Brown   test:
1050c4762a1bSJed Brown     suffix: q1_2d_plex_0
105130602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -convergence
1052c4762a1bSJed Brown   test:
1053c4762a1bSJed Brown     suffix: q1_2d_plex_1
105430602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 1
1055c4762a1bSJed Brown   test:
1056c4762a1bSJed Brown     suffix: q1_2d_plex_2
105730602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 2
1058c4762a1bSJed Brown   test:
1059c4762a1bSJed Brown     suffix: q1_2d_plex_3
106030602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 1 -shear_coords
1061c4762a1bSJed Brown   test:
1062c4762a1bSJed Brown     suffix: q1_2d_plex_4
106330602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 2 -shear_coords
1064c4762a1bSJed Brown   test:
1065c4762a1bSJed Brown     suffix: q1_2d_plex_5
106630602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 0 -non_affine_coords -convergence
1067c4762a1bSJed Brown   test:
1068c4762a1bSJed Brown     suffix: q1_2d_plex_6
106930602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 1 -non_affine_coords -convergence
1070c4762a1bSJed Brown   test:
1071c4762a1bSJed Brown     suffix: q1_2d_plex_7
107230602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 2 -non_affine_coords -convergence
1073c4762a1bSJed Brown 
1074c4762a1bSJed Brown   # 2D Q_2 on a quadrilaterial
1075c4762a1bSJed Brown   test:
1076c4762a1bSJed Brown     suffix: q2_2d_plex_0
107730602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -convergence
1078c4762a1bSJed Brown   test:
1079c4762a1bSJed Brown     suffix: q2_2d_plex_1
108030602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1
1081c4762a1bSJed Brown   test:
1082c4762a1bSJed Brown     suffix: q2_2d_plex_2
108330602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2
1084c4762a1bSJed Brown   test:
1085c4762a1bSJed Brown     suffix: q2_2d_plex_3
108630602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 -shear_coords
1087c4762a1bSJed Brown   test:
1088c4762a1bSJed Brown     suffix: q2_2d_plex_4
108930602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 -shear_coords
1090c4762a1bSJed Brown   test:
1091c4762a1bSJed Brown     suffix: q2_2d_plex_5
109230602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 0 -non_affine_coords -convergence
1093c4762a1bSJed Brown   test:
1094c4762a1bSJed Brown     suffix: q2_2d_plex_6
109530602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 1 -non_affine_coords -convergence
1096c4762a1bSJed Brown   test:
1097c4762a1bSJed Brown     suffix: q2_2d_plex_7
109830602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 2 -non_affine_coords -convergence
1099c4762a1bSJed Brown 
1100c4762a1bSJed Brown   # 2D P_3 on a triangle
1101c4762a1bSJed Brown   test:
1102c4762a1bSJed Brown     suffix: p3_2d_0
1103c4762a1bSJed Brown     requires: triangle !single
1104c4762a1bSJed Brown     args: -petscspace_degree 3 -qorder 3 -convergence
1105c4762a1bSJed Brown   test:
1106c4762a1bSJed Brown     suffix: p3_2d_1
1107c4762a1bSJed Brown     requires: triangle !single
1108c4762a1bSJed Brown     args: -petscspace_degree 3 -qorder 3 -porder 1
1109c4762a1bSJed Brown   test:
1110c4762a1bSJed Brown     suffix: p3_2d_2
1111c4762a1bSJed Brown     requires: triangle !single
1112c4762a1bSJed Brown     args: -petscspace_degree 3 -qorder 3 -porder 2
1113c4762a1bSJed Brown   test:
1114c4762a1bSJed Brown     suffix: p3_2d_3
1115c4762a1bSJed Brown     requires: triangle !single
1116c4762a1bSJed Brown     args: -petscspace_degree 3 -qorder 3 -porder 3
1117c4762a1bSJed Brown   test:
1118c4762a1bSJed Brown     suffix: p3_2d_4
1119e788613bSJoe Wallwork     requires: triangle mmg
1120c4762a1bSJed Brown     args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -convergence -conv_refine 0
1121c4762a1bSJed Brown   test:
1122c4762a1bSJed Brown     suffix: p3_2d_5
1123e788613bSJoe Wallwork     requires: triangle mmg
1124c4762a1bSJed Brown     args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -porder 1 -conv_refine 0
1125c4762a1bSJed Brown   test:
1126c4762a1bSJed Brown     suffix: p3_2d_6
1127e788613bSJoe Wallwork     requires: triangle mmg
1128c4762a1bSJed Brown     args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -porder 3 -conv_refine 0
1129c4762a1bSJed Brown 
1130c4762a1bSJed Brown   # 2D Q_3 on a quadrilaterial
1131c4762a1bSJed Brown   test:
1132c4762a1bSJed Brown     suffix: q3_2d_0
113399a192c5SJunchao Zhang     requires: !single
113430602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -convergence
1135c4762a1bSJed Brown   test:
1136c4762a1bSJed Brown     suffix: q3_2d_1
113799a192c5SJunchao Zhang     requires: !single
113830602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 1
1139c4762a1bSJed Brown   test:
1140c4762a1bSJed Brown     suffix: q3_2d_2
114199a192c5SJunchao Zhang     requires: !single
114230602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 2
1143c4762a1bSJed Brown   test:
1144c4762a1bSJed Brown     suffix: q3_2d_3
114599a192c5SJunchao Zhang     requires: !single
114630602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 3
1147c4762a1bSJed Brown 
1148c4762a1bSJed Brown   # 2D P_1disc on a triangle/quadrilateral
1149c4762a1bSJed Brown   test:
1150c4762a1bSJed Brown     suffix: p1d_2d_0
1151c4762a1bSJed Brown     requires: triangle
1152c4762a1bSJed Brown     args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -convergence
1153c4762a1bSJed Brown   test:
1154c4762a1bSJed Brown     suffix: p1d_2d_1
1155c4762a1bSJed Brown     requires: triangle
1156c4762a1bSJed Brown     args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 1
1157c4762a1bSJed Brown   test:
1158c4762a1bSJed Brown     suffix: p1d_2d_2
1159c4762a1bSJed Brown     requires: triangle
1160c4762a1bSJed Brown     args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 2
1161c4762a1bSJed Brown   test:
1162c4762a1bSJed Brown     suffix: p1d_2d_3
1163c4762a1bSJed Brown     requires: triangle
116430602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -convergence
1165c4762a1bSJed Brown     filter: sed  -e "s/convergence rate at refinement 0: 2/convergence rate at refinement 0: 1.9/g"
1166c4762a1bSJed Brown   test:
1167c4762a1bSJed Brown     suffix: p1d_2d_4
1168c4762a1bSJed Brown     requires: triangle
116930602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 1
1170c4762a1bSJed Brown   test:
1171c4762a1bSJed Brown     suffix: p1d_2d_5
1172c4762a1bSJed Brown     requires: triangle
117330602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 2
1174c4762a1bSJed Brown 
1175c4762a1bSJed Brown   # 2D BDM_1 on a triangle
1176c4762a1bSJed Brown   test:
1177c4762a1bSJed Brown     suffix: bdm1_2d_0
1178c4762a1bSJed Brown     requires: triangle
1179c4762a1bSJed Brown     args: -petscspace_degree 1 -petscdualspace_type bdm \
1180c4762a1bSJed Brown           -num_comp 2 -qorder 1 -convergence
1181c4762a1bSJed Brown   test:
1182c4762a1bSJed Brown     suffix: bdm1_2d_1
1183c4762a1bSJed Brown     requires: triangle
1184c4762a1bSJed Brown     args: -petscspace_degree 1 -petscdualspace_type bdm \
1185c4762a1bSJed Brown           -num_comp 2 -qorder 1 -porder 1
1186c4762a1bSJed Brown   test:
1187c4762a1bSJed Brown     suffix: bdm1_2d_2
1188c4762a1bSJed Brown     requires: triangle
1189c4762a1bSJed Brown     args: -petscspace_degree 1 -petscdualspace_type bdm \
1190c4762a1bSJed Brown           -num_comp 2 -qorder 1 -porder 2
1191c4762a1bSJed Brown 
1192c4762a1bSJed Brown   # 2D BDM_1 on a quadrilateral
1193c4762a1bSJed Brown   test:
1194c4762a1bSJed Brown     suffix: bdm1q_2d_0
1195c4762a1bSJed Brown     requires: triangle
1196c4762a1bSJed Brown     args: -petscspace_degree 1 -petscdualspace_type bdm \
11973f27d899SToby Isaac           -petscdualspace_lagrange_tensor 1 \
119830602db0SMatthew G. Knepley           -dm_plex_simplex 0 -num_comp 2 -qorder 1 -convergence
1199c4762a1bSJed Brown   test:
1200c4762a1bSJed Brown     suffix: bdm1q_2d_1
1201c4762a1bSJed Brown     requires: triangle
1202c4762a1bSJed Brown     args: -petscspace_degree 1 -petscdualspace_type bdm \
12033f27d899SToby Isaac           -petscdualspace_lagrange_tensor 1 \
120430602db0SMatthew G. Knepley           -dm_plex_simplex 0 -num_comp 2 -qorder 1 -porder 1
1205c4762a1bSJed Brown   test:
1206c4762a1bSJed Brown     suffix: bdm1q_2d_2
1207c4762a1bSJed Brown     requires: triangle
1208c4762a1bSJed Brown     args: -petscspace_degree 1 -petscdualspace_type bdm \
12093f27d899SToby Isaac           -petscdualspace_lagrange_tensor 1 \
121030602db0SMatthew G. Knepley           -dm_plex_simplex 0 -num_comp 2 -qorder 1 -porder 2
1211c4762a1bSJed Brown 
1212c4762a1bSJed Brown   # Test high order quadrature
1213c4762a1bSJed Brown   test:
1214c4762a1bSJed Brown     suffix: p1_quad_2
1215c4762a1bSJed Brown     requires: triangle
1216c4762a1bSJed Brown     args: -petscspace_degree 1 -qorder 2 -porder 1
1217c4762a1bSJed Brown   test:
1218c4762a1bSJed Brown     suffix: p1_quad_5
1219c4762a1bSJed Brown     requires: triangle
1220c4762a1bSJed Brown     args: -petscspace_degree 1 -qorder 5 -porder 1
1221c4762a1bSJed Brown   test:
1222c4762a1bSJed Brown     suffix: p2_quad_3
1223c4762a1bSJed Brown     requires: triangle
1224c4762a1bSJed Brown     args: -petscspace_degree 2 -qorder 3 -porder 2
1225c4762a1bSJed Brown   test:
1226c4762a1bSJed Brown     suffix: p2_quad_5
1227c4762a1bSJed Brown     requires: triangle
1228c4762a1bSJed Brown     args: -petscspace_degree 2 -qorder 5 -porder 2
1229c4762a1bSJed Brown   test:
1230c4762a1bSJed Brown     suffix: q1_quad_2
123130602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 2 -porder 1
1232c4762a1bSJed Brown   test:
1233c4762a1bSJed Brown     suffix: q1_quad_5
123430602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 5 -porder 1
1235c4762a1bSJed Brown   test:
1236c4762a1bSJed Brown     suffix: q2_quad_3
123730602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 3 -porder 1
1238c4762a1bSJed Brown   test:
1239c4762a1bSJed Brown     suffix: q2_quad_5
124030602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 5 -porder 1
1241c4762a1bSJed Brown 
1242c4762a1bSJed Brown   # Nonconforming tests
1243c4762a1bSJed Brown   test:
1244c4762a1bSJed Brown     suffix: constraints
124530602db0SMatthew G. Knepley     args: -dm_coord_space 0 -dm_plex_simplex 0 -petscspace_type tensor -petscspace_degree 1 -qorder 0 -constraints
1246c4762a1bSJed Brown   test:
1247c4762a1bSJed Brown     suffix: nonconforming_tensor_2
1248c4762a1bSJed Brown     nsize: 4
124930602db0SMatthew 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
1250c4762a1bSJed Brown   test:
1251c4762a1bSJed Brown     suffix: nonconforming_tensor_3
1252c4762a1bSJed Brown     nsize: 4
125330602db0SMatthew 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
1254c4762a1bSJed Brown   test:
1255c4762a1bSJed Brown     suffix: nonconforming_tensor_2_fv
1256c4762a1bSJed Brown     nsize: 4
125730602db0SMatthew 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
1258c4762a1bSJed Brown   test:
1259c4762a1bSJed Brown     suffix: nonconforming_tensor_3_fv
1260c4762a1bSJed Brown     nsize: 4
126130602db0SMatthew 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
1262c4762a1bSJed Brown   test:
1263c4762a1bSJed Brown     suffix: nonconforming_tensor_2_hi
1264c4762a1bSJed Brown     requires: !single
1265c4762a1bSJed Brown     nsize: 4
126630602db0SMatthew 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
1267c4762a1bSJed Brown   test:
1268c4762a1bSJed Brown     suffix: nonconforming_tensor_3_hi
1269c4762a1bSJed Brown     requires: !single skip
1270c4762a1bSJed Brown     nsize: 4
127130602db0SMatthew 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
1272c4762a1bSJed Brown   test:
1273c4762a1bSJed Brown     suffix: nonconforming_simplex_2
1274c4762a1bSJed Brown     requires: triangle
1275c4762a1bSJed Brown     nsize: 4
127630602db0SMatthew 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
1277c4762a1bSJed Brown   test:
1278c4762a1bSJed Brown     suffix: nonconforming_simplex_2_hi
1279c4762a1bSJed Brown     requires: triangle !single
1280c4762a1bSJed Brown     nsize: 4
128130602db0SMatthew G. Knepley     args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_max_projection_height 1 -petscspace_degree 4 -qorder 4
1282c4762a1bSJed Brown   test:
1283c4762a1bSJed Brown     suffix: nonconforming_simplex_2_fv
1284c4762a1bSJed Brown     requires: triangle
1285c4762a1bSJed Brown     nsize: 4
128630602db0SMatthew G. Knepley     args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -num_comp 2
1287c4762a1bSJed Brown   test:
1288c4762a1bSJed Brown     suffix: nonconforming_simplex_3
1289c4762a1bSJed Brown     requires: ctetgen
1290c4762a1bSJed Brown     nsize: 4
129130602db0SMatthew 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
1292c4762a1bSJed Brown   test:
1293c4762a1bSJed Brown     suffix: nonconforming_simplex_3_hi
1294c4762a1bSJed Brown     requires: ctetgen skip
1295c4762a1bSJed Brown     nsize: 4
129630602db0SMatthew 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
1297c4762a1bSJed Brown   test:
1298c4762a1bSJed Brown     suffix: nonconforming_simplex_3_fv
1299c4762a1bSJed Brown     requires: ctetgen
1300c4762a1bSJed Brown     nsize: 4
130130602db0SMatthew 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
1302c4762a1bSJed Brown 
1303d21efd2eSMatthew G. Knepley   # 3D WXY on a triangular prism
1304d21efd2eSMatthew G. Knepley   test:
1305d21efd2eSMatthew G. Knepley     suffix: wxy_0
1306d21efd2eSMatthew G. Knepley     args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -qorder 2 -porder 0 \
1307e239af90SMatthew G. Knepley           -petscspace_type sum \
1308e239af90SMatthew G. Knepley           -petscspace_variables 3 \
1309e239af90SMatthew G. Knepley           -petscspace_components 3 \
1310e239af90SMatthew G. Knepley           -petscspace_sum_spaces 2 \
1311e239af90SMatthew G. Knepley           -petscspace_sum_concatenate false \
1312e239af90SMatthew G. Knepley           -sumcomp_0_petscspace_variables 3 \
1313e239af90SMatthew G. Knepley           -sumcomp_0_petscspace_components 3 \
1314e239af90SMatthew G. Knepley           -sumcomp_0_petscspace_degree 1 \
1315e239af90SMatthew G. Knepley           -sumcomp_1_petscspace_variables 3 \
1316e239af90SMatthew G. Knepley           -sumcomp_1_petscspace_components 3 \
1317e239af90SMatthew G. Knepley           -sumcomp_1_petscspace_type wxy \
1318e239af90SMatthew G. Knepley           -petscdualspace_refcell triangular_prism \
1319e239af90SMatthew G. Knepley           -petscdualspace_form_degree 0 \
1320e239af90SMatthew G. Knepley           -petscdualspace_order 1 \
1321e239af90SMatthew G. Knepley           -petscdualspace_components 3
1322d21efd2eSMatthew G. Knepley 
1323c4762a1bSJed Brown TEST*/
1324c4762a1bSJed Brown 
1325c4762a1bSJed Brown /*
1326c4762a1bSJed Brown    # 2D Q_2 on a quadrilaterial Plex
1327c4762a1bSJed Brown   test:
1328c4762a1bSJed Brown     suffix: q2_2d_plex_0
132930602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -convergence
1330c4762a1bSJed Brown   test:
1331c4762a1bSJed Brown     suffix: q2_2d_plex_1
133230602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1
1333c4762a1bSJed Brown   test:
1334c4762a1bSJed Brown     suffix: q2_2d_plex_2
133530602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2
1336c4762a1bSJed Brown   test:
1337c4762a1bSJed Brown     suffix: q2_2d_plex_3
133830602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 -shear_coords
1339c4762a1bSJed Brown   test:
1340c4762a1bSJed Brown     suffix: q2_2d_plex_4
134130602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 -shear_coords
1342c4762a1bSJed Brown   test:
1343c4762a1bSJed Brown     suffix: q2_2d_plex_5
134430602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 0 -non_affine_coords
1345c4762a1bSJed Brown   test:
1346c4762a1bSJed Brown     suffix: q2_2d_plex_6
134730602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 1 -non_affine_coords
1348c4762a1bSJed Brown   test:
1349c4762a1bSJed Brown     suffix: q2_2d_plex_7
135030602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 2 -non_affine_coords
1351c4762a1bSJed Brown 
1352c4762a1bSJed Brown   test:
1353c4762a1bSJed Brown     suffix: p1d_2d_6
1354e788613bSJoe Wallwork     requires: mmg
1355c4762a1bSJed Brown     args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0
1356c4762a1bSJed Brown   test:
1357c4762a1bSJed Brown     suffix: p1d_2d_7
1358e788613bSJoe Wallwork     requires: mmg
1359c4762a1bSJed Brown     args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0
1360c4762a1bSJed Brown   test:
1361c4762a1bSJed Brown     suffix: p1d_2d_8
1362e788613bSJoe Wallwork     requires: mmg
1363c4762a1bSJed Brown     args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0
1364c4762a1bSJed Brown */
1365