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