xref: /petsc/src/dm/impls/plex/tests/ex3.c (revision a69119a591a03a9d906b29c0a4e9802e4d7c9795)
1 static char help[] = "Check that a DM can accurately represent and interpolate functions of a given polynomial order\n\n";
2 
3 #include <petscdmplex.h>
4 #include <petscdm.h>
5 #include <petscdmda.h>
6 #include <petscfe.h>
7 #include <petscds.h>
8 #include <petscksp.h>
9 #include <petscsnes.h>
10 
11 typedef struct {
12   /* Domain and mesh definition */
13   PetscBool useDA;           /* Flag DMDA tensor product mesh */
14   PetscBool shearCoords;     /* Flag for shear transform */
15   PetscBool nonaffineCoords; /* Flag for non-affine transform */
16   /* Element definition */
17   PetscInt  qorder;        /* Order of the quadrature */
18   PetscInt  numComponents; /* Number of field components */
19   PetscFE   fe;            /* The finite element */
20   /* Testing space */
21   PetscInt  porder;         /* Order of polynomials to test */
22   PetscBool convergence;    /* Test for order of convergence */
23   PetscBool convRefine;     /* Test for convergence using refinement, otherwise use coarsening */
24   PetscBool constraints;    /* Test local constraints */
25   PetscBool tree;           /* Test tree routines */
26   PetscBool testFEjacobian; /* Test finite element Jacobian assembly */
27   PetscBool testFVgrad;     /* Test finite difference gradient routine */
28   PetscBool testInjector;   /* Test finite element injection routines */
29   PetscInt  treeCell;       /* Cell to refine in tree test */
30   PetscReal constants[3];   /* Constant values for each dimension */
31 } AppCtx;
32 
33 /* u = 1 */
34 PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx) {
35   AppCtx  *user = (AppCtx *)ctx;
36   PetscInt d;
37   for (d = 0; d < dim; ++d) u[d] = user->constants[d];
38   return 0;
39 }
40 PetscErrorCode constantDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) {
41   PetscInt d;
42   for (d = 0; d < dim; ++d) u[d] = 0.0;
43   return 0;
44 }
45 
46 /* u = x */
47 PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx) {
48   PetscInt d;
49   for (d = 0; d < dim; ++d) u[d] = coords[d];
50   return 0;
51 }
52 PetscErrorCode linearDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) {
53   PetscInt d, e;
54   for (d = 0; d < dim; ++d) {
55     u[d] = 0.0;
56     for (e = 0; e < dim; ++e) u[d] += (d == e ? 1.0 : 0.0) * n[e];
57   }
58   return 0;
59 }
60 
61 /* u = x^2 or u = (x^2, xy) or u = (xy, yz, zx) */
62 PetscErrorCode quadratic(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx) {
63   if (dim > 2) {
64     u[0] = coords[0] * coords[1];
65     u[1] = coords[1] * coords[2];
66     u[2] = coords[2] * coords[0];
67   } else if (dim > 1) {
68     u[0] = coords[0] * coords[0];
69     u[1] = coords[0] * coords[1];
70   } else if (dim > 0) {
71     u[0] = coords[0] * coords[0];
72   }
73   return 0;
74 }
75 PetscErrorCode quadraticDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) {
76   if (dim > 2) {
77     u[0] = coords[1] * n[0] + coords[0] * n[1];
78     u[1] = coords[2] * n[1] + coords[1] * n[2];
79     u[2] = coords[2] * n[0] + coords[0] * n[2];
80   } else if (dim > 1) {
81     u[0] = 2.0 * coords[0] * n[0];
82     u[1] = coords[1] * n[0] + coords[0] * n[1];
83   } else if (dim > 0) {
84     u[0] = 2.0 * coords[0] * n[0];
85   }
86   return 0;
87 }
88 
89 /* u = x^3 or u = (x^3, x^2y) or u = (x^2y, y^2z, z^2x) */
90 PetscErrorCode cubic(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx) {
91   if (dim > 2) {
92     u[0] = coords[0] * coords[0] * coords[1];
93     u[1] = coords[1] * coords[1] * coords[2];
94     u[2] = coords[2] * coords[2] * coords[0];
95   } else if (dim > 1) {
96     u[0] = coords[0] * coords[0] * coords[0];
97     u[1] = coords[0] * coords[0] * coords[1];
98   } else if (dim > 0) {
99     u[0] = coords[0] * coords[0] * coords[0];
100   }
101   return 0;
102 }
103 PetscErrorCode cubicDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) {
104   if (dim > 2) {
105     u[0] = 2.0 * coords[0] * coords[1] * n[0] + coords[0] * coords[0] * n[1];
106     u[1] = 2.0 * coords[1] * coords[2] * n[1] + coords[1] * coords[1] * n[2];
107     u[2] = 2.0 * coords[2] * coords[0] * n[2] + coords[2] * coords[2] * n[0];
108   } else if (dim > 1) {
109     u[0] = 3.0 * coords[0] * coords[0] * n[0];
110     u[1] = 2.0 * coords[0] * coords[1] * n[0] + coords[0] * coords[0] * n[1];
111   } else if (dim > 0) {
112     u[0] = 3.0 * coords[0] * coords[0] * n[0];
113   }
114   return 0;
115 }
116 
117 /* u = tanh(x) */
118 PetscErrorCode trig(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx) {
119   PetscInt d;
120   for (d = 0; d < dim; ++d) u[d] = PetscTanhReal(coords[d] - 0.5);
121   return 0;
122 }
123 PetscErrorCode trigDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) {
124   PetscInt d;
125   for (d = 0; d < dim; ++d) u[d] = 1.0 / PetscSqr(PetscCoshReal(coords[d] - 0.5)) * n[d];
126   return 0;
127 }
128 
129 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
130   PetscInt n = 3;
131 
132   PetscFunctionBeginUser;
133   options->useDA           = PETSC_FALSE;
134   options->shearCoords     = PETSC_FALSE;
135   options->nonaffineCoords = PETSC_FALSE;
136   options->qorder          = 0;
137   options->numComponents   = PETSC_DEFAULT;
138   options->porder          = 0;
139   options->convergence     = PETSC_FALSE;
140   options->convRefine      = PETSC_TRUE;
141   options->constraints     = PETSC_FALSE;
142   options->tree            = PETSC_FALSE;
143   options->treeCell        = 0;
144   options->testFEjacobian  = PETSC_FALSE;
145   options->testFVgrad      = PETSC_FALSE;
146   options->testInjector    = PETSC_FALSE;
147   options->constants[0]    = 1.0;
148   options->constants[1]    = 2.0;
149   options->constants[2]    = 3.0;
150 
151   PetscOptionsBegin(comm, "", "Projection Test Options", "DMPlex");
152   PetscCall(PetscOptionsBool("-use_da", "Flag for DMDA mesh", "ex3.c", options->useDA, &options->useDA, NULL));
153   PetscCall(PetscOptionsBool("-shear_coords", "Transform coordinates with a shear", "ex3.c", options->shearCoords, &options->shearCoords, NULL));
154   PetscCall(PetscOptionsBool("-non_affine_coords", "Transform coordinates with a non-affine transform", "ex3.c", options->nonaffineCoords, &options->nonaffineCoords, NULL));
155   PetscCall(PetscOptionsBoundedInt("-qorder", "The quadrature order", "ex3.c", options->qorder, &options->qorder, NULL, 0));
156   PetscCall(PetscOptionsBoundedInt("-num_comp", "The number of field components", "ex3.c", options->numComponents, &options->numComponents, NULL, PETSC_DEFAULT));
157   PetscCall(PetscOptionsBoundedInt("-porder", "The order of polynomials to test", "ex3.c", options->porder, &options->porder, NULL, 0));
158   PetscCall(PetscOptionsBool("-convergence", "Check the convergence rate", "ex3.c", options->convergence, &options->convergence, NULL));
159   PetscCall(PetscOptionsBool("-conv_refine", "Use refinement for the convergence rate", "ex3.c", options->convRefine, &options->convRefine, NULL));
160   PetscCall(PetscOptionsBool("-constraints", "Test local constraints (serial only)", "ex3.c", options->constraints, &options->constraints, NULL));
161   PetscCall(PetscOptionsBool("-tree", "Test tree routines", "ex3.c", options->tree, &options->tree, NULL));
162   PetscCall(PetscOptionsBoundedInt("-tree_cell", "cell to refine in tree test", "ex3.c", options->treeCell, &options->treeCell, NULL, 0));
163   PetscCall(PetscOptionsBool("-test_fe_jacobian", "Test finite element Jacobian assembly", "ex3.c", options->testFEjacobian, &options->testFEjacobian, NULL));
164   PetscCall(PetscOptionsBool("-test_fv_grad", "Test finite volume gradient reconstruction", "ex3.c", options->testFVgrad, &options->testFVgrad, NULL));
165   PetscCall(PetscOptionsBool("-test_injector", "Test finite element injection", "ex3.c", options->testInjector, &options->testInjector, NULL));
166   PetscCall(PetscOptionsRealArray("-constants", "Set the constant values", "ex3.c", options->constants, &n, NULL));
167   PetscOptionsEnd();
168 
169   PetscFunctionReturn(0);
170 }
171 
172 static PetscErrorCode TransformCoordinates(DM dm, AppCtx *user) {
173   PetscSection coordSection;
174   Vec          coordinates;
175   PetscScalar *coords;
176   PetscInt     vStart, vEnd, v;
177 
178   PetscFunctionBeginUser;
179   if (user->nonaffineCoords) {
180     /* x' = r^(1/p) (x/r), y' = r^(1/p) (y/r), z' = z */
181     PetscCall(DMGetCoordinateSection(dm, &coordSection));
182     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
183     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
184     PetscCall(VecGetArray(coordinates, &coords));
185     for (v = vStart; v < vEnd; ++v) {
186       PetscInt  dof, off;
187       PetscReal p = 4.0, r;
188 
189       PetscCall(PetscSectionGetDof(coordSection, v, &dof));
190       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
191       switch (dof) {
192       case 2:
193         r               = PetscSqr(PetscRealPart(coords[off + 0])) + PetscSqr(PetscRealPart(coords[off + 1]));
194         coords[off + 0] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p) / (2 * p)) * coords[off + 0];
195         coords[off + 1] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p) / (2 * p)) * coords[off + 1];
196         break;
197       case 3:
198         r               = PetscSqr(PetscRealPart(coords[off + 0])) + PetscSqr(PetscRealPart(coords[off + 1]));
199         coords[off + 0] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p) / (2 * p)) * coords[off + 0];
200         coords[off + 1] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p) / (2 * p)) * coords[off + 1];
201         coords[off + 2] = coords[off + 2];
202         break;
203       }
204     }
205     PetscCall(VecRestoreArray(coordinates, &coords));
206   }
207   if (user->shearCoords) {
208     /* x' = x + m y + m z, y' = y + m z,  z' = z */
209     PetscCall(DMGetCoordinateSection(dm, &coordSection));
210     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
211     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
212     PetscCall(VecGetArray(coordinates, &coords));
213     for (v = vStart; v < vEnd; ++v) {
214       PetscInt  dof, off;
215       PetscReal m = 1.0;
216 
217       PetscCall(PetscSectionGetDof(coordSection, v, &dof));
218       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
219       switch (dof) {
220       case 2:
221         coords[off + 0] = coords[off + 0] + m * coords[off + 1];
222         coords[off + 1] = coords[off + 1];
223         break;
224       case 3:
225         coords[off + 0] = coords[off + 0] + m * coords[off + 1] + m * coords[off + 2];
226         coords[off + 1] = coords[off + 1] + m * coords[off + 2];
227         coords[off + 2] = coords[off + 2];
228         break;
229       }
230     }
231     PetscCall(VecRestoreArray(coordinates, &coords));
232   }
233   PetscFunctionReturn(0);
234 }
235 
236 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) {
237   PetscInt  dim = 2;
238   PetscBool simplex;
239 
240   PetscFunctionBeginUser;
241   if (user->useDA) {
242     switch (dim) {
243     case 2:
244       PetscCall(DMDACreate2d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 2, 2, PETSC_DETERMINE, PETSC_DETERMINE, 1, 1, NULL, NULL, dm));
245       PetscCall(DMSetFromOptions(*dm));
246       PetscCall(DMSetUp(*dm));
247       PetscCall(DMDASetVertexCoordinates(*dm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
248       break;
249     default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot create structured mesh of dimension %" PetscInt_FMT, dim);
250     }
251     PetscCall(PetscObjectSetName((PetscObject)*dm, "Hexahedral Mesh"));
252   } else {
253     PetscCall(DMCreate(comm, dm));
254     PetscCall(DMSetType(*dm, DMPLEX));
255     PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
256     PetscCall(DMSetFromOptions(*dm));
257 
258     PetscCall(DMGetDimension(*dm, &dim));
259     PetscCall(DMPlexIsSimplex(*dm, &simplex));
260     PetscCallMPI(MPI_Bcast(&simplex, 1, MPIU_BOOL, 0, comm));
261     if (user->tree) {
262       DM refTree, ncdm = NULL;
263 
264       PetscCall(DMPlexCreateDefaultReferenceTree(comm, dim, simplex, &refTree));
265       PetscCall(DMViewFromOptions(refTree, NULL, "-reftree_dm_view"));
266       PetscCall(DMPlexSetReferenceTree(*dm, refTree));
267       PetscCall(DMDestroy(&refTree));
268       PetscCall(DMPlexTreeRefineCell(*dm, user->treeCell, &ncdm));
269       if (ncdm) {
270         PetscCall(DMDestroy(dm));
271         *dm = ncdm;
272         PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_FALSE));
273       }
274       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, "tree_"));
275       PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
276       PetscCall(DMSetFromOptions(*dm));
277       PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
278     } else {
279       PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_TRUE));
280     }
281     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, "dist_"));
282     PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
283     PetscCall(DMSetFromOptions(*dm));
284     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, NULL));
285     if (simplex) PetscCall(PetscObjectSetName((PetscObject)*dm, "Simplicial Mesh"));
286     else PetscCall(PetscObjectSetName((PetscObject)*dm, "Hexahedral Mesh"));
287   }
288   PetscCall(DMSetFromOptions(*dm));
289   PetscCall(TransformCoordinates(*dm, user));
290   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
291   PetscFunctionReturn(0);
292 }
293 
294 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[]) {
295   PetscInt d, e;
296   for (d = 0, e = 0; d < dim; d++, e += dim + 1) g0[e] = 1.;
297 }
298 
299 /* < \nabla v, 1/2(\nabla u + {\nabla u}^T) > */
300 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[]) {
301   PetscInt compI, compJ, d, e;
302 
303   for (compI = 0; compI < dim; ++compI) {
304     for (compJ = 0; compJ < dim; ++compJ) {
305       for (d = 0; d < dim; ++d) {
306         for (e = 0; e < dim; e++) {
307           if (d == e && d == compI && d == compJ) {
308             C[((compI * dim + compJ) * dim + d) * dim + e] = 1.0;
309           } else if ((d == compJ && e == compI) || (d == e && compI == compJ)) {
310             C[((compI * dim + compJ) * dim + d) * dim + e] = 0.5;
311           } else {
312             C[((compI * dim + compJ) * dim + d) * dim + e] = 0.0;
313           }
314         }
315       }
316     }
317   }
318 }
319 
320 static PetscErrorCode SetupSection(DM dm, AppCtx *user) {
321   PetscFunctionBeginUser;
322   if (user->constraints) {
323     /* test local constraints */
324     DM           coordDM;
325     PetscInt     fStart, fEnd, f, vStart, vEnd, v;
326     PetscInt     edgesx = 2, vertsx;
327     PetscInt     edgesy = 2, vertsy;
328     PetscMPIInt  size;
329     PetscInt     numConst;
330     PetscSection aSec;
331     PetscInt    *anchors;
332     PetscInt     offset;
333     IS           aIS;
334     MPI_Comm     comm = PetscObjectComm((PetscObject)dm);
335 
336     PetscCallMPI(MPI_Comm_size(comm, &size));
337     PetscCheck(size <= 1, comm, PETSC_ERR_SUP, "Local constraint test can only be performed in serial");
338 
339     /* we are going to test constraints by using them to enforce periodicity
340      * in one direction, and comparing to the existing method of enforcing
341      * periodicity */
342 
343     /* first create the coordinate section so that it does not clone the
344      * constraints */
345     PetscCall(DMGetCoordinateDM(dm, &coordDM));
346 
347     /* create the constrained-to-anchor section */
348     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
349     PetscCall(DMPlexGetDepthStratum(dm, 1, &fStart, &fEnd));
350     PetscCall(PetscSectionCreate(PETSC_COMM_SELF, &aSec));
351     PetscCall(PetscSectionSetChart(aSec, PetscMin(fStart, vStart), PetscMax(fEnd, vEnd)));
352 
353     /* define the constraints */
354     PetscCall(PetscOptionsGetInt(NULL, NULL, "-da_grid_x", &edgesx, NULL));
355     PetscCall(PetscOptionsGetInt(NULL, NULL, "-da_grid_y", &edgesy, NULL));
356     vertsx   = edgesx + 1;
357     vertsy   = edgesy + 1;
358     numConst = vertsy + edgesy;
359     PetscCall(PetscMalloc1(numConst, &anchors));
360     offset = 0;
361     for (v = vStart + edgesx; v < vEnd; v += vertsx) {
362       PetscCall(PetscSectionSetDof(aSec, v, 1));
363       anchors[offset++] = v - edgesx;
364     }
365     for (f = fStart + edgesx * vertsy + edgesx * edgesy; f < fEnd; f++) {
366       PetscCall(PetscSectionSetDof(aSec, f, 1));
367       anchors[offset++] = f - edgesx * edgesy;
368     }
369     PetscCall(PetscSectionSetUp(aSec));
370     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numConst, anchors, PETSC_OWN_POINTER, &aIS));
371 
372     PetscCall(DMPlexSetAnchors(dm, aSec, aIS));
373     PetscCall(PetscSectionDestroy(&aSec));
374     PetscCall(ISDestroy(&aIS));
375   }
376   PetscCall(DMSetNumFields(dm, 1));
377   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)user->fe));
378   PetscCall(DMCreateDS(dm));
379   if (user->constraints) {
380     /* test getting local constraint matrix that matches section */
381     PetscSection aSec;
382     IS           aIS;
383 
384     PetscCall(DMPlexGetAnchors(dm, &aSec, &aIS));
385     if (aSec) {
386       PetscDS         ds;
387       PetscSection    cSec, section;
388       PetscInt        cStart, cEnd, c, numComp;
389       Mat             cMat, mass;
390       Vec             local;
391       const PetscInt *anchors;
392 
393       PetscCall(DMGetLocalSection(dm, &section));
394       /* this creates the matrix and preallocates the matrix structure: we
395        * just have to fill in the values */
396       PetscCall(DMGetDefaultConstraints(dm, &cSec, &cMat, NULL));
397       PetscCall(PetscSectionGetChart(cSec, &cStart, &cEnd));
398       PetscCall(ISGetIndices(aIS, &anchors));
399       PetscCall(PetscFEGetNumComponents(user->fe, &numComp));
400       for (c = cStart; c < cEnd; c++) {
401         PetscInt cDof;
402 
403         /* is this point constrained? (does it have an anchor?) */
404         PetscCall(PetscSectionGetDof(aSec, c, &cDof));
405         if (cDof) {
406           PetscInt cOff, a, aDof, aOff, j;
407           PetscCheck(cDof == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Found %" PetscInt_FMT " anchor points: should be just one", cDof);
408 
409           /* find the anchor point */
410           PetscCall(PetscSectionGetOffset(aSec, c, &cOff));
411           a = anchors[cOff];
412 
413           /* find the constrained dofs (row in constraint matrix) */
414           PetscCall(PetscSectionGetDof(cSec, c, &cDof));
415           PetscCall(PetscSectionGetOffset(cSec, c, &cOff));
416 
417           /* find the anchor dofs (column in constraint matrix) */
418           PetscCall(PetscSectionGetDof(section, a, &aDof));
419           PetscCall(PetscSectionGetOffset(section, a, &aOff));
420 
421           PetscCheck(cDof == aDof, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point and anchor have different number of dofs: %" PetscInt_FMT ", %" PetscInt_FMT, cDof, aDof);
422           PetscCheck(cDof % numComp == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point dofs not divisible by field components: %" PetscInt_FMT ", %" PetscInt_FMT, cDof, numComp);
423 
424           /* put in a simple equality constraint */
425           for (j = 0; j < cDof; j++) PetscCall(MatSetValue(cMat, cOff + j, aOff + j, 1., INSERT_VALUES));
426         }
427       }
428       PetscCall(MatAssemblyBegin(cMat, MAT_FINAL_ASSEMBLY));
429       PetscCall(MatAssemblyEnd(cMat, MAT_FINAL_ASSEMBLY));
430       PetscCall(ISRestoreIndices(aIS, &anchors));
431 
432       /* Now that we have constructed the constraint matrix, any FE matrix
433        * that we construct will apply the constraints during construction */
434 
435       PetscCall(DMCreateMatrix(dm, &mass));
436       /* get a dummy local variable to serve as the solution */
437       PetscCall(DMGetLocalVector(dm, &local));
438       PetscCall(DMGetDS(dm, &ds));
439       /* set the jacobian to be the mass matrix */
440       PetscCall(PetscDSSetJacobian(ds, 0, 0, simple_mass, NULL, NULL, NULL));
441       /* build the mass matrix */
442       PetscCall(DMPlexSNESComputeJacobianFEM(dm, local, mass, mass, NULL));
443       PetscCall(MatView(mass, PETSC_VIEWER_STDOUT_WORLD));
444       PetscCall(MatDestroy(&mass));
445       PetscCall(DMRestoreLocalVector(dm, &local));
446     }
447   }
448   PetscFunctionReturn(0);
449 }
450 
451 static PetscErrorCode TestFEJacobian(DM dm, AppCtx *user) {
452   PetscFunctionBeginUser;
453   if (!user->useDA) {
454     Vec          local;
455     const Vec   *vecs;
456     Mat          E;
457     MatNullSpace sp;
458     PetscBool    isNullSpace, hasConst;
459     PetscInt     dim, n, i;
460     Vec          res = NULL, localX, localRes;
461     PetscDS      ds;
462 
463     PetscCall(DMGetDimension(dm, &dim));
464     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);
465     PetscCall(DMGetDS(dm, &ds));
466     PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, symmetric_gradient_inner_product));
467     PetscCall(DMCreateMatrix(dm, &E));
468     PetscCall(DMGetLocalVector(dm, &local));
469     PetscCall(DMPlexSNESComputeJacobianFEM(dm, local, E, E, NULL));
470     PetscCall(DMPlexCreateRigidBody(dm, 0, &sp));
471     PetscCall(MatNullSpaceGetVecs(sp, &hasConst, &n, &vecs));
472     if (n) PetscCall(VecDuplicate(vecs[0], &res));
473     PetscCall(DMCreateLocalVector(dm, &localX));
474     PetscCall(DMCreateLocalVector(dm, &localRes));
475     for (i = 0; i < n; i++) { /* also test via matrix-free Jacobian application */
476       PetscReal resNorm;
477 
478       PetscCall(VecSet(localRes, 0.));
479       PetscCall(VecSet(localX, 0.));
480       PetscCall(VecSet(local, 0.));
481       PetscCall(VecSet(res, 0.));
482       PetscCall(DMGlobalToLocalBegin(dm, vecs[i], INSERT_VALUES, localX));
483       PetscCall(DMGlobalToLocalEnd(dm, vecs[i], INSERT_VALUES, localX));
484       PetscCall(DMSNESComputeJacobianAction(dm, local, localX, localRes, NULL));
485       PetscCall(DMLocalToGlobalBegin(dm, localRes, ADD_VALUES, res));
486       PetscCall(DMLocalToGlobalEnd(dm, localRes, ADD_VALUES, res));
487       PetscCall(VecNorm(res, NORM_2, &resNorm));
488       if (resNorm > PETSC_SMALL) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Symmetric gradient action null space vector %" PetscInt_FMT " residual: %E\n", i, (double)resNorm));
489     }
490     PetscCall(VecDestroy(&localRes));
491     PetscCall(VecDestroy(&localX));
492     PetscCall(VecDestroy(&res));
493     PetscCall(MatNullSpaceTest(sp, E, &isNullSpace));
494     if (isNullSpace) {
495       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Symmetric gradient null space: PASS\n"));
496     } else {
497       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Symmetric gradient null space: FAIL\n"));
498     }
499     PetscCall(MatNullSpaceDestroy(&sp));
500     PetscCall(MatDestroy(&E));
501     PetscCall(DMRestoreLocalVector(dm, &local));
502   }
503   PetscFunctionReturn(0);
504 }
505 
506 static PetscErrorCode TestInjector(DM dm, AppCtx *user) {
507   DM          refTree;
508   PetscMPIInt rank;
509 
510   PetscFunctionBegin;
511   PetscCall(DMPlexGetReferenceTree(dm, &refTree));
512   if (refTree) {
513     Mat inj;
514 
515     PetscCall(DMPlexComputeInjectorReferenceTree(refTree, &inj));
516     PetscCall(PetscObjectSetName((PetscObject)inj, "Reference Tree Injector"));
517     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
518     if (rank == 0) PetscCall(MatView(inj, PETSC_VIEWER_STDOUT_SELF));
519     PetscCall(MatDestroy(&inj));
520   }
521   PetscFunctionReturn(0);
522 }
523 
524 static PetscErrorCode TestFVGrad(DM dm, AppCtx *user) {
525   MPI_Comm           comm;
526   DM                 dmRedist, dmfv, dmgrad, dmCell, refTree;
527   PetscFV            fv;
528   PetscInt           dim, nvecs, v, cStart, cEnd, cEndInterior;
529   PetscMPIInt        size;
530   Vec                cellgeom, grad, locGrad;
531   const PetscScalar *cgeom;
532   PetscReal          allVecMaxDiff = 0., fvTol = 100. * PETSC_MACHINE_EPSILON;
533 
534   PetscFunctionBeginUser;
535   comm = PetscObjectComm((PetscObject)dm);
536   PetscCall(DMGetDimension(dm, &dim));
537   /* duplicate DM, give dup. a FV discretization */
538   PetscCall(DMSetBasicAdjacency(dm, PETSC_TRUE, PETSC_FALSE));
539   PetscCallMPI(MPI_Comm_size(comm, &size));
540   dmRedist = NULL;
541   if (size > 1) PetscCall(DMPlexDistributeOverlap(dm, 1, NULL, &dmRedist));
542   if (!dmRedist) {
543     dmRedist = dm;
544     PetscCall(PetscObjectReference((PetscObject)dmRedist));
545   }
546   PetscCall(PetscFVCreate(comm, &fv));
547   PetscCall(PetscFVSetType(fv, PETSCFVLEASTSQUARES));
548   PetscCall(PetscFVSetNumComponents(fv, user->numComponents));
549   PetscCall(PetscFVSetSpatialDimension(fv, dim));
550   PetscCall(PetscFVSetFromOptions(fv));
551   PetscCall(PetscFVSetUp(fv));
552   PetscCall(DMPlexConstructGhostCells(dmRedist, NULL, NULL, &dmfv));
553   PetscCall(DMDestroy(&dmRedist));
554   PetscCall(DMSetNumFields(dmfv, 1));
555   PetscCall(DMSetField(dmfv, 0, NULL, (PetscObject)fv));
556   PetscCall(DMCreateDS(dmfv));
557   PetscCall(DMPlexGetReferenceTree(dm, &refTree));
558   if (refTree) PetscCall(DMCopyDisc(dmfv, refTree));
559   PetscCall(DMPlexGetGradientDM(dmfv, fv, &dmgrad));
560   PetscCall(DMPlexGetHeightStratum(dmfv, 0, &cStart, &cEnd));
561   nvecs = dim * (dim + 1) / 2;
562   PetscCall(DMPlexGetGeometryFVM(dmfv, NULL, &cellgeom, NULL));
563   PetscCall(VecGetDM(cellgeom, &dmCell));
564   PetscCall(VecGetArrayRead(cellgeom, &cgeom));
565   PetscCall(DMGetGlobalVector(dmgrad, &grad));
566   PetscCall(DMGetLocalVector(dmgrad, &locGrad));
567   PetscCall(DMPlexGetGhostCellStratum(dmgrad, &cEndInterior, NULL));
568   cEndInterior = (cEndInterior < 0) ? cEnd : cEndInterior;
569   for (v = 0; v < nvecs; v++) {
570     Vec                locX;
571     PetscInt           c;
572     PetscScalar        trueGrad[3][3] = {{0.}};
573     const PetscScalar *gradArray;
574     PetscReal          maxDiff, maxDiffGlob;
575 
576     PetscCall(DMGetLocalVector(dmfv, &locX));
577     /* get the local projection of the rigid body mode */
578     for (c = cStart; c < cEnd; c++) {
579       PetscFVCellGeom *cg;
580       PetscScalar      cx[3] = {0., 0., 0.};
581 
582       PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg));
583       if (v < dim) {
584         cx[v] = 1.;
585       } else {
586         PetscInt w = v - dim;
587 
588         cx[(w + 1) % dim] = cg->centroid[(w + 2) % dim];
589         cx[(w + 2) % dim] = -cg->centroid[(w + 1) % dim];
590       }
591       PetscCall(DMPlexVecSetClosure(dmfv, NULL, locX, c, cx, INSERT_ALL_VALUES));
592     }
593     /* TODO: this isn't in any header */
594     PetscCall(DMPlexReconstructGradientsFVM(dmfv, locX, grad));
595     PetscCall(DMGlobalToLocalBegin(dmgrad, grad, INSERT_VALUES, locGrad));
596     PetscCall(DMGlobalToLocalEnd(dmgrad, grad, INSERT_VALUES, locGrad));
597     PetscCall(VecGetArrayRead(locGrad, &gradArray));
598     /* compare computed gradient to exact gradient */
599     if (v >= dim) {
600       PetscInt w = v - dim;
601 
602       trueGrad[(w + 1) % dim][(w + 2) % dim] = 1.;
603       trueGrad[(w + 2) % dim][(w + 1) % dim] = -1.;
604     }
605     maxDiff = 0.;
606     for (c = cStart; c < cEndInterior; c++) {
607       PetscScalar *compGrad;
608       PetscInt     i, j, k;
609       PetscReal    FrobDiff = 0.;
610 
611       PetscCall(DMPlexPointLocalRead(dmgrad, c, gradArray, &compGrad));
612 
613       for (i = 0, k = 0; i < dim; i++) {
614         for (j = 0; j < dim; j++, k++) {
615           PetscScalar diff = compGrad[k] - trueGrad[i][j];
616           FrobDiff += PetscRealPart(diff * PetscConj(diff));
617         }
618       }
619       FrobDiff = PetscSqrtReal(FrobDiff);
620       maxDiff  = PetscMax(maxDiff, FrobDiff);
621     }
622     PetscCallMPI(MPI_Allreduce(&maxDiff, &maxDiffGlob, 1, MPIU_REAL, MPIU_MAX, comm));
623     allVecMaxDiff = PetscMax(allVecMaxDiff, maxDiffGlob);
624     PetscCall(VecRestoreArrayRead(locGrad, &gradArray));
625     PetscCall(DMRestoreLocalVector(dmfv, &locX));
626   }
627   if (allVecMaxDiff < fvTol) {
628     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Finite volume gradient reconstruction: PASS\n"));
629   } else {
630     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Finite volume gradient reconstruction: FAIL at tolerance %g with max difference %g\n", (double)fvTol, (double)allVecMaxDiff));
631   }
632   PetscCall(DMRestoreLocalVector(dmgrad, &locGrad));
633   PetscCall(DMRestoreGlobalVector(dmgrad, &grad));
634   PetscCall(VecRestoreArrayRead(cellgeom, &cgeom));
635   PetscCall(DMDestroy(&dmfv));
636   PetscCall(PetscFVDestroy(&fv));
637   PetscFunctionReturn(0);
638 }
639 
640 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) {
641   Vec       u;
642   PetscReal n[3] = {1.0, 1.0, 1.0};
643 
644   PetscFunctionBeginUser;
645   PetscCall(DMGetGlobalVector(dm, &u));
646   /* Project function into FE function space */
647   PetscCall(DMProjectFunction(dm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, u));
648   PetscCall(VecViewFromOptions(u, NULL, "-projection_view"));
649   /* Compare approximation to exact in L_2 */
650   PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, exactCtxs, u, error));
651   PetscCall(DMComputeL2GradientDiff(dm, 0.0, exactFuncDers, exactCtxs, u, n, errorDer));
652   PetscCall(DMRestoreGlobalVector(dm, &u));
653   PetscFunctionReturn(0);
654 }
655 
656 static PetscErrorCode CheckFunctions(DM dm, PetscInt order, AppCtx *user) {
657   PetscErrorCode (*exactFuncs[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
658   PetscErrorCode (*exactFuncDers[1])(PetscInt dim, PetscReal time, const PetscReal x[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx);
659   void     *exactCtxs[3];
660   MPI_Comm  comm;
661   PetscReal error, errorDer, tol = PETSC_SMALL;
662 
663   PetscFunctionBeginUser;
664   exactCtxs[0] = user;
665   exactCtxs[1] = user;
666   exactCtxs[2] = user;
667   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
668   /* Setup functions to approximate */
669   switch (order) {
670   case 0:
671     exactFuncs[0]    = constant;
672     exactFuncDers[0] = constantDer;
673     break;
674   case 1:
675     exactFuncs[0]    = linear;
676     exactFuncDers[0] = linearDer;
677     break;
678   case 2:
679     exactFuncs[0]    = quadratic;
680     exactFuncDers[0] = quadraticDer;
681     break;
682   case 3:
683     exactFuncs[0]    = cubic;
684     exactFuncDers[0] = cubicDer;
685     break;
686   default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Could not determine functions to test for order %" PetscInt_FMT, order);
687   }
688   PetscCall(ComputeError(dm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user));
689   /* Report result */
690   if (error > tol) PetscCall(PetscPrintf(comm, "Function tests FAIL for order %" PetscInt_FMT " at tolerance %g error %g\n", order, (double)tol, (double)error));
691   else PetscCall(PetscPrintf(comm, "Function tests pass for order %" PetscInt_FMT " at tolerance %g\n", order, (double)tol));
692   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));
693   else PetscCall(PetscPrintf(comm, "Function tests pass for order %" PetscInt_FMT " derivatives at tolerance %g\n", order, (double)tol));
694   PetscFunctionReturn(0);
695 }
696 
697 static PetscErrorCode CheckInterpolation(DM dm, PetscBool checkRestrict, PetscInt order, AppCtx *user) {
698   PetscErrorCode (*exactFuncs[1])(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
699   PetscErrorCode (*exactFuncDers[1])(PetscInt, PetscReal, const PetscReal x[], const PetscReal n[], PetscInt, PetscScalar *u, void *ctx);
700   PetscReal n[3] = {1.0, 1.0, 1.0};
701   void     *exactCtxs[3];
702   DM        rdm, idm, fdm;
703   Mat       Interp;
704   Vec       iu, fu, scaling;
705   MPI_Comm  comm;
706   PetscInt  dim;
707   PetscReal error, errorDer, tol = PETSC_SMALL;
708 
709   PetscFunctionBeginUser;
710   exactCtxs[0] = user;
711   exactCtxs[1] = user;
712   exactCtxs[2] = user;
713   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
714   PetscCall(DMGetDimension(dm, &dim));
715   PetscCall(DMRefine(dm, comm, &rdm));
716   PetscCall(DMSetCoarseDM(rdm, dm));
717   PetscCall(DMPlexSetRegularRefinement(rdm, user->convRefine));
718   if (user->tree) {
719     DM refTree;
720     PetscCall(DMPlexGetReferenceTree(dm, &refTree));
721     PetscCall(DMPlexSetReferenceTree(rdm, refTree));
722   }
723   if (user->useDA) PetscCall(DMDASetVertexCoordinates(rdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
724   PetscCall(SetupSection(rdm, user));
725   /* Setup functions to approximate */
726   switch (order) {
727   case 0:
728     exactFuncs[0]    = constant;
729     exactFuncDers[0] = constantDer;
730     break;
731   case 1:
732     exactFuncs[0]    = linear;
733     exactFuncDers[0] = linearDer;
734     break;
735   case 2:
736     exactFuncs[0]    = quadratic;
737     exactFuncDers[0] = quadraticDer;
738     break;
739   case 3:
740     exactFuncs[0]    = cubic;
741     exactFuncDers[0] = cubicDer;
742     break;
743   default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Could not determine functions to test for dimension %" PetscInt_FMT " order %" PetscInt_FMT, dim, order);
744   }
745   idm = checkRestrict ? rdm : dm;
746   fdm = checkRestrict ? dm : rdm;
747   PetscCall(DMGetGlobalVector(idm, &iu));
748   PetscCall(DMGetGlobalVector(fdm, &fu));
749   PetscCall(DMSetApplicationContext(dm, user));
750   PetscCall(DMSetApplicationContext(rdm, user));
751   PetscCall(DMCreateInterpolation(dm, rdm, &Interp, &scaling));
752   /* Project function into initial FE function space */
753   PetscCall(DMProjectFunction(idm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, iu));
754   /* Interpolate function into final FE function space */
755   if (checkRestrict) {
756     PetscCall(MatRestrict(Interp, iu, fu));
757     PetscCall(VecPointwiseMult(fu, scaling, fu));
758   } else PetscCall(MatInterpolate(Interp, iu, fu));
759   /* Compare approximation to exact in L_2 */
760   PetscCall(DMComputeL2Diff(fdm, 0.0, exactFuncs, exactCtxs, fu, &error));
761   PetscCall(DMComputeL2GradientDiff(fdm, 0.0, exactFuncDers, exactCtxs, fu, n, &errorDer));
762   /* Report result */
763   if (error > tol) PetscCall(PetscPrintf(comm, "Interpolation tests FAIL for order %" PetscInt_FMT " at tolerance %g error %g\n", order, (double)tol, (double)error));
764   else PetscCall(PetscPrintf(comm, "Interpolation tests pass for order %" PetscInt_FMT " at tolerance %g\n", order, (double)tol));
765   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));
766   else PetscCall(PetscPrintf(comm, "Interpolation tests pass for order %" PetscInt_FMT " derivatives at tolerance %g\n", order, (double)tol));
767   PetscCall(DMRestoreGlobalVector(idm, &iu));
768   PetscCall(DMRestoreGlobalVector(fdm, &fu));
769   PetscCall(MatDestroy(&Interp));
770   PetscCall(VecDestroy(&scaling));
771   PetscCall(DMDestroy(&rdm));
772   PetscFunctionReturn(0);
773 }
774 
775 static PetscErrorCode CheckConvergence(DM dm, PetscInt Nr, AppCtx *user) {
776   DM odm = dm, rdm = NULL, cdm = NULL;
777   PetscErrorCode (*exactFuncs[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)                         = {trig};
778   PetscErrorCode (*exactFuncDers[1])(PetscInt dim, PetscReal time, const PetscReal x[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) = {trigDer};
779   void     *exactCtxs[3];
780   PetscInt  r, c, cStart, cEnd;
781   PetscReal errorOld, errorDerOld, error, errorDer, rel, len, lenOld;
782   double    p;
783 
784   PetscFunctionBeginUser;
785   if (!user->convergence) PetscFunctionReturn(0);
786   exactCtxs[0] = user;
787   exactCtxs[1] = user;
788   exactCtxs[2] = user;
789   PetscCall(PetscObjectReference((PetscObject)odm));
790   if (!user->convRefine) {
791     for (r = 0; r < Nr; ++r) {
792       PetscCall(DMRefine(odm, PetscObjectComm((PetscObject)dm), &rdm));
793       PetscCall(DMDestroy(&odm));
794       odm = rdm;
795     }
796     PetscCall(SetupSection(odm, user));
797   }
798   PetscCall(ComputeError(odm, exactFuncs, exactFuncDers, exactCtxs, &errorOld, &errorDerOld, user));
799   if (user->convRefine) {
800     for (r = 0; r < Nr; ++r) {
801       PetscCall(DMRefine(odm, PetscObjectComm((PetscObject)dm), &rdm));
802       if (user->useDA) PetscCall(DMDASetVertexCoordinates(rdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
803       PetscCall(SetupSection(rdm, user));
804       PetscCall(ComputeError(rdm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user));
805       p = PetscLog2Real(errorOld / error);
806       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Function   convergence rate at refinement %" PetscInt_FMT ": %.2f\n", r, (double)p));
807       p = PetscLog2Real(errorDerOld / errorDer);
808       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Derivative convergence rate at refinement %" PetscInt_FMT ": %.2f\n", r, (double)p));
809       PetscCall(DMDestroy(&odm));
810       odm         = rdm;
811       errorOld    = error;
812       errorDerOld = errorDer;
813     }
814   } else {
815     /* PetscCall(ComputeLongestEdge(dm, &lenOld)); */
816     PetscCall(DMPlexGetHeightStratum(odm, 0, &cStart, &cEnd));
817     lenOld = cEnd - cStart;
818     for (c = 0; c < Nr; ++c) {
819       PetscCall(DMCoarsen(odm, PetscObjectComm((PetscObject)dm), &cdm));
820       if (user->useDA) PetscCall(DMDASetVertexCoordinates(cdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
821       PetscCall(SetupSection(cdm, user));
822       PetscCall(ComputeError(cdm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user));
823       /* PetscCall(ComputeLongestEdge(cdm, &len)); */
824       PetscCall(DMPlexGetHeightStratum(cdm, 0, &cStart, &cEnd));
825       len = cEnd - cStart;
826       rel = error / errorOld;
827       p   = PetscLogReal(rel) / PetscLogReal(lenOld / len);
828       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Function   convergence rate at coarsening %" PetscInt_FMT ": %.2f\n", c, (double)p));
829       rel = errorDer / errorDerOld;
830       p   = PetscLogReal(rel) / PetscLogReal(lenOld / len);
831       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Derivative convergence rate at coarsening %" PetscInt_FMT ": %.2f\n", c, (double)p));
832       PetscCall(DMDestroy(&odm));
833       odm         = cdm;
834       errorOld    = error;
835       errorDerOld = errorDer;
836       lenOld      = len;
837     }
838   }
839   PetscCall(DMDestroy(&odm));
840   PetscFunctionReturn(0);
841 }
842 
843 int main(int argc, char **argv) {
844   DM        dm;
845   AppCtx    user; /* user-defined work context */
846   PetscInt  dim     = 2;
847   PetscBool simplex = PETSC_FALSE;
848 
849   PetscFunctionBeginUser;
850   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
851   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
852   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
853   if (!user.useDA) {
854     PetscCall(DMGetDimension(dm, &dim));
855     PetscCall(DMPlexIsSimplex(dm, &simplex));
856   }
857   PetscCall(DMPlexMetricSetFromOptions(dm));
858   user.numComponents = user.numComponents < 0 ? dim : user.numComponents;
859   PetscCall(PetscFECreateDefault(PETSC_COMM_WORLD, dim, user.numComponents, simplex, NULL, user.qorder, &user.fe));
860   PetscCall(SetupSection(dm, &user));
861   if (user.testFEjacobian) PetscCall(TestFEJacobian(dm, &user));
862   if (user.testFVgrad) PetscCall(TestFVGrad(dm, &user));
863   if (user.testInjector) PetscCall(TestInjector(dm, &user));
864   PetscCall(CheckFunctions(dm, user.porder, &user));
865   {
866     PetscDualSpace dsp;
867     PetscInt       k;
868 
869     PetscCall(PetscFEGetDualSpace(user.fe, &dsp));
870     PetscCall(PetscDualSpaceGetDeRahm(dsp, &k));
871     if (dim == 2 && simplex == PETSC_TRUE && user.tree == PETSC_FALSE && k == 0) {
872       PetscCall(CheckInterpolation(dm, PETSC_FALSE, user.porder, &user));
873       PetscCall(CheckInterpolation(dm, PETSC_TRUE, user.porder, &user));
874     }
875   }
876   PetscCall(CheckConvergence(dm, 3, &user));
877   PetscCall(PetscFEDestroy(&user.fe));
878   PetscCall(DMDestroy(&dm));
879   PetscCall(PetscFinalize());
880   return 0;
881 }
882 
883 /*TEST
884 
885   test:
886     suffix: 1
887     requires: triangle
888 
889   # 2D P_1 on a triangle
890   test:
891     suffix: p1_2d_0
892     requires: triangle
893     args: -petscspace_degree 1 -qorder 1 -convergence
894   test:
895     suffix: p1_2d_1
896     requires: triangle
897     args: -petscspace_degree 1 -qorder 1 -porder 1
898   test:
899     suffix: p1_2d_2
900     requires: triangle
901     args: -petscspace_degree 1 -qorder 1 -porder 2
902   test:
903     suffix: p1_2d_3
904     requires: triangle mmg
905     args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0
906   test:
907     suffix: p1_2d_4
908     requires: triangle mmg
909     args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0
910   test:
911     suffix: p1_2d_5
912     requires: triangle mmg
913     args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0
914 
915   # 3D P_1 on a tetrahedron
916   test:
917     suffix: p1_3d_0
918     requires: ctetgen
919     args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -convergence
920   test:
921     suffix: p1_3d_1
922     requires: ctetgen
923     args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -porder 1
924   test:
925     suffix: p1_3d_2
926     requires: ctetgen
927     args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -porder 2
928   test:
929     suffix: p1_3d_3
930     requires: ctetgen mmg
931     args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0
932   test:
933     suffix: p1_3d_4
934     requires: ctetgen mmg
935     args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0
936   test:
937     suffix: p1_3d_5
938     requires: ctetgen mmg
939     args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0
940 
941   # 2D P_2 on a triangle
942   test:
943     suffix: p2_2d_0
944     requires: triangle
945     args: -petscspace_degree 2 -qorder 2 -convergence
946   test:
947     suffix: p2_2d_1
948     requires: triangle
949     args: -petscspace_degree 2 -qorder 2 -porder 1
950   test:
951     suffix: p2_2d_2
952     requires: triangle
953     args: -petscspace_degree 2 -qorder 2 -porder 2
954   test:
955     suffix: p2_2d_3
956     requires: triangle mmg
957     args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -convergence -conv_refine 0
958   test:
959     suffix: p2_2d_4
960     requires: triangle mmg
961     args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 1 -conv_refine 0
962   test:
963     suffix: p2_2d_5
964     requires: triangle mmg
965     args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 2 -conv_refine 0
966 
967   # 3D P_2 on a tetrahedron
968   test:
969     suffix: p2_3d_0
970     requires: ctetgen
971     args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -convergence
972   test:
973     suffix: p2_3d_1
974     requires: ctetgen
975     args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -porder 1
976   test:
977     suffix: p2_3d_2
978     requires: ctetgen
979     args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -porder 2
980   test:
981     suffix: p2_3d_3
982     requires: ctetgen mmg
983     args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -convergence -conv_refine 0
984   test:
985     suffix: p2_3d_4
986     requires: ctetgen mmg
987     args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 1 -conv_refine 0
988   test:
989     suffix: p2_3d_5
990     requires: ctetgen mmg
991     args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 2 -conv_refine 0
992 
993   # 2D Q_1 on a quadrilaterial DA
994   test:
995     suffix: q1_2d_da_0
996     requires: broken
997     args: -use_da 1 -petscspace_degree 1 -qorder 1 -convergence
998   test:
999     suffix: q1_2d_da_1
1000     requires: broken
1001     args: -use_da 1 -petscspace_degree 1 -qorder 1 -porder 1
1002   test:
1003     suffix: q1_2d_da_2
1004     requires: broken
1005     args: -use_da 1 -petscspace_degree 1 -qorder 1 -porder 2
1006 
1007   # 2D Q_1 on a quadrilaterial Plex
1008   test:
1009     suffix: q1_2d_plex_0
1010     args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -convergence
1011   test:
1012     suffix: q1_2d_plex_1
1013     args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 1
1014   test:
1015     suffix: q1_2d_plex_2
1016     args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 2
1017   test:
1018     suffix: q1_2d_plex_3
1019     args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 1 -shear_coords
1020   test:
1021     suffix: q1_2d_plex_4
1022     args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 2 -shear_coords
1023   test:
1024     suffix: q1_2d_plex_5
1025     args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 0 -non_affine_coords -convergence
1026   test:
1027     suffix: q1_2d_plex_6
1028     args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 1 -non_affine_coords -convergence
1029   test:
1030     suffix: q1_2d_plex_7
1031     args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 2 -non_affine_coords -convergence
1032 
1033   # 2D Q_2 on a quadrilaterial
1034   test:
1035     suffix: q2_2d_plex_0
1036     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -convergence
1037   test:
1038     suffix: q2_2d_plex_1
1039     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1
1040   test:
1041     suffix: q2_2d_plex_2
1042     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2
1043   test:
1044     suffix: q2_2d_plex_3
1045     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 -shear_coords
1046   test:
1047     suffix: q2_2d_plex_4
1048     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 -shear_coords
1049   test:
1050     suffix: q2_2d_plex_5
1051     args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 0 -non_affine_coords -convergence
1052   test:
1053     suffix: q2_2d_plex_6
1054     args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 1 -non_affine_coords -convergence
1055   test:
1056     suffix: q2_2d_plex_7
1057     args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 2 -non_affine_coords -convergence
1058 
1059   # 2D P_3 on a triangle
1060   test:
1061     suffix: p3_2d_0
1062     requires: triangle !single
1063     args: -petscspace_degree 3 -qorder 3 -convergence
1064   test:
1065     suffix: p3_2d_1
1066     requires: triangle !single
1067     args: -petscspace_degree 3 -qorder 3 -porder 1
1068   test:
1069     suffix: p3_2d_2
1070     requires: triangle !single
1071     args: -petscspace_degree 3 -qorder 3 -porder 2
1072   test:
1073     suffix: p3_2d_3
1074     requires: triangle !single
1075     args: -petscspace_degree 3 -qorder 3 -porder 3
1076   test:
1077     suffix: p3_2d_4
1078     requires: triangle mmg
1079     args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -convergence -conv_refine 0
1080   test:
1081     suffix: p3_2d_5
1082     requires: triangle mmg
1083     args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -porder 1 -conv_refine 0
1084   test:
1085     suffix: p3_2d_6
1086     requires: triangle mmg
1087     args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -porder 3 -conv_refine 0
1088 
1089   # 2D Q_3 on a quadrilaterial
1090   test:
1091     suffix: q3_2d_0
1092     requires: !single
1093     args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -convergence
1094   test:
1095     suffix: q3_2d_1
1096     requires: !single
1097     args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 1
1098   test:
1099     suffix: q3_2d_2
1100     requires: !single
1101     args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 2
1102   test:
1103     suffix: q3_2d_3
1104     requires: !single
1105     args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 3
1106 
1107   # 2D P_1disc on a triangle/quadrilateral
1108   test:
1109     suffix: p1d_2d_0
1110     requires: triangle
1111     args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -convergence
1112   test:
1113     suffix: p1d_2d_1
1114     requires: triangle
1115     args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 1
1116   test:
1117     suffix: p1d_2d_2
1118     requires: triangle
1119     args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 2
1120   test:
1121     suffix: p1d_2d_3
1122     requires: triangle
1123     args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -convergence
1124     filter: sed  -e "s/convergence rate at refinement 0: 2/convergence rate at refinement 0: 1.9/g"
1125   test:
1126     suffix: p1d_2d_4
1127     requires: triangle
1128     args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 1
1129   test:
1130     suffix: p1d_2d_5
1131     requires: triangle
1132     args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 2
1133 
1134   # 2D BDM_1 on a triangle
1135   test:
1136     suffix: bdm1_2d_0
1137     requires: triangle
1138     args: -petscspace_degree 1 -petscdualspace_type bdm \
1139           -num_comp 2 -qorder 1 -convergence
1140   test:
1141     suffix: bdm1_2d_1
1142     requires: triangle
1143     args: -petscspace_degree 1 -petscdualspace_type bdm \
1144           -num_comp 2 -qorder 1 -porder 1
1145   test:
1146     suffix: bdm1_2d_2
1147     requires: triangle
1148     args: -petscspace_degree 1 -petscdualspace_type bdm \
1149           -num_comp 2 -qorder 1 -porder 2
1150 
1151   # 2D BDM_1 on a quadrilateral
1152   test:
1153     suffix: bdm1q_2d_0
1154     requires: triangle
1155     args: -petscspace_degree 1 -petscdualspace_type bdm \
1156           -petscdualspace_lagrange_tensor 1 \
1157           -dm_plex_simplex 0 -num_comp 2 -qorder 1 -convergence
1158   test:
1159     suffix: bdm1q_2d_1
1160     requires: triangle
1161     args: -petscspace_degree 1 -petscdualspace_type bdm \
1162           -petscdualspace_lagrange_tensor 1 \
1163           -dm_plex_simplex 0 -num_comp 2 -qorder 1 -porder 1
1164   test:
1165     suffix: bdm1q_2d_2
1166     requires: triangle
1167     args: -petscspace_degree 1 -petscdualspace_type bdm \
1168           -petscdualspace_lagrange_tensor 1 \
1169           -dm_plex_simplex 0 -num_comp 2 -qorder 1 -porder 2
1170 
1171   # Test high order quadrature
1172   test:
1173     suffix: p1_quad_2
1174     requires: triangle
1175     args: -petscspace_degree 1 -qorder 2 -porder 1
1176   test:
1177     suffix: p1_quad_5
1178     requires: triangle
1179     args: -petscspace_degree 1 -qorder 5 -porder 1
1180   test:
1181     suffix: p2_quad_3
1182     requires: triangle
1183     args: -petscspace_degree 2 -qorder 3 -porder 2
1184   test:
1185     suffix: p2_quad_5
1186     requires: triangle
1187     args: -petscspace_degree 2 -qorder 5 -porder 2
1188   test:
1189     suffix: q1_quad_2
1190     args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 2 -porder 1
1191   test:
1192     suffix: q1_quad_5
1193     args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 5 -porder 1
1194   test:
1195     suffix: q2_quad_3
1196     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 3 -porder 1
1197   test:
1198     suffix: q2_quad_5
1199     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 5 -porder 1
1200 
1201   # Nonconforming tests
1202   test:
1203     suffix: constraints
1204     args: -dm_coord_space 0 -dm_plex_simplex 0 -petscspace_type tensor -petscspace_degree 1 -qorder 0 -constraints
1205   test:
1206     suffix: nonconforming_tensor_2
1207     nsize: 4
1208     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
1209   test:
1210     suffix: nonconforming_tensor_3
1211     nsize: 4
1212     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
1213   test:
1214     suffix: nonconforming_tensor_2_fv
1215     nsize: 4
1216     args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -dm_plex_simplex 0 -num_comp 2
1217   test:
1218     suffix: nonconforming_tensor_3_fv
1219     nsize: 4
1220     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
1221   test:
1222     suffix: nonconforming_tensor_2_hi
1223     requires: !single
1224     nsize: 4
1225     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
1226   test:
1227     suffix: nonconforming_tensor_3_hi
1228     requires: !single skip
1229     nsize: 4
1230     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
1231   test:
1232     suffix: nonconforming_simplex_2
1233     requires: triangle
1234     nsize: 4
1235     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
1236   test:
1237     suffix: nonconforming_simplex_2_hi
1238     requires: triangle !single
1239     nsize: 4
1240     args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_max_projection_height 1 -petscspace_degree 4 -qorder 4
1241   test:
1242     suffix: nonconforming_simplex_2_fv
1243     requires: triangle
1244     nsize: 4
1245     args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -num_comp 2
1246   test:
1247     suffix: nonconforming_simplex_3
1248     requires: ctetgen
1249     nsize: 4
1250     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
1251   test:
1252     suffix: nonconforming_simplex_3_hi
1253     requires: ctetgen skip
1254     nsize: 4
1255     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
1256   test:
1257     suffix: nonconforming_simplex_3_fv
1258     requires: ctetgen
1259     nsize: 4
1260     args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -dm_plex_dim 3 -num_comp 3
1261 
1262   # 3D WXY on a triangular prism
1263   test:
1264     suffix: wxy_0
1265     args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -qorder 2 -porder 0 \
1266           -petscspace_type sum \
1267           -petscspace_variables 3 \
1268           -petscspace_components 3 \
1269           -petscspace_sum_spaces 2 \
1270           -petscspace_sum_concatenate false \
1271           -sumcomp_0_petscspace_variables 3 \
1272           -sumcomp_0_petscspace_components 3 \
1273           -sumcomp_0_petscspace_degree 1 \
1274           -sumcomp_1_petscspace_variables 3 \
1275           -sumcomp_1_petscspace_components 3 \
1276           -sumcomp_1_petscspace_type wxy \
1277           -petscdualspace_refcell triangular_prism \
1278           -petscdualspace_form_degree 0 \
1279           -petscdualspace_order 1 \
1280           -petscdualspace_components 3
1281 
1282 TEST*/
1283 
1284 /*
1285    # 2D Q_2 on a quadrilaterial Plex
1286   test:
1287     suffix: q2_2d_plex_0
1288     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -convergence
1289   test:
1290     suffix: q2_2d_plex_1
1291     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1
1292   test:
1293     suffix: q2_2d_plex_2
1294     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2
1295   test:
1296     suffix: q2_2d_plex_3
1297     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 -shear_coords
1298   test:
1299     suffix: q2_2d_plex_4
1300     args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 -shear_coords
1301   test:
1302     suffix: q2_2d_plex_5
1303     args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 0 -non_affine_coords
1304   test:
1305     suffix: q2_2d_plex_6
1306     args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 1 -non_affine_coords
1307   test:
1308     suffix: q2_2d_plex_7
1309     args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 2 -non_affine_coords
1310 
1311   test:
1312     suffix: p1d_2d_6
1313     requires: mmg
1314     args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0
1315   test:
1316     suffix: p1d_2d_7
1317     requires: mmg
1318     args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0
1319   test:
1320     suffix: p1d_2d_8
1321     requires: mmg
1322     args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0
1323 */
1324