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