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