xref: /petsc/src/dm/impls/plex/tests/ex73.c (revision a2424aed8403ea1ba9ad665c8e808ee45e27c88b)
1*a2424aedSMatthew G. Knepley static char help[] = "Tests for Gauss' Law\n\n";
2*a2424aedSMatthew G. Knepley 
3*a2424aedSMatthew G. Knepley /* We want to check the weak version of Gauss' Law, namely that
4*a2424aedSMatthew G. Knepley 
5*a2424aedSMatthew G. Knepley   \int_\Omega v div q - \int_\Gamma v (q \cdot n) = 0
6*a2424aedSMatthew G. Knepley 
7*a2424aedSMatthew G. Knepley */
8*a2424aedSMatthew G. Knepley 
9*a2424aedSMatthew G. Knepley #include <petscdmplex.h>
10*a2424aedSMatthew G. Knepley #include <petscsnes.h>
11*a2424aedSMatthew G. Knepley #include <petscds.h>
12*a2424aedSMatthew G. Knepley 
13*a2424aedSMatthew G. Knepley typedef struct {
14*a2424aedSMatthew G. Knepley   PetscInt  degree;  // The degree of the discretization
15*a2424aedSMatthew G. Knepley   PetscBool divFree; // True if the solution is divergence-free
16*a2424aedSMatthew G. Knepley } AppCtx;
17*a2424aedSMatthew G. Knepley 
18*a2424aedSMatthew G. Knepley static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
19*a2424aedSMatthew G. Knepley {
20*a2424aedSMatthew G. Knepley   u[0] = 0.0;
21*a2424aedSMatthew G. Knepley   return PETSC_SUCCESS;
22*a2424aedSMatthew G. Knepley }
23*a2424aedSMatthew G. Knepley 
24*a2424aedSMatthew G. Knepley // div <x^d y^{d-1}, -x^{d-1} y^d> = 0
25*a2424aedSMatthew G. Knepley static void solenoidal_2d(PetscInt n, const PetscReal x[], PetscScalar u[])
26*a2424aedSMatthew G. Knepley {
27*a2424aedSMatthew G. Knepley   u[0] = PetscPowRealInt(x[0], n) * PetscPowRealInt(x[1], n - 1);
28*a2424aedSMatthew G. Knepley   u[1] = -PetscPowRealInt(x[0], n - 1) * PetscPowRealInt(x[1], n);
29*a2424aedSMatthew G. Knepley }
30*a2424aedSMatthew G. Knepley // div <x^d y^{d-1} z^{d-1}, -2 x^{d-1} y^d z^{d-1}, x^{d-1} y^{d-1} z^d> = 0
31*a2424aedSMatthew G. Knepley static void solenoidal_3d(PetscInt n, const PetscReal x[], PetscScalar u[])
32*a2424aedSMatthew G. Knepley {
33*a2424aedSMatthew G. Knepley   u[0] = PetscPowRealInt(x[0], n) * PetscPowRealInt(x[1], n - 1) * PetscPowRealInt(x[2], n - 1);
34*a2424aedSMatthew G. Knepley   u[1] = -2. * PetscPowRealInt(x[0], n - 1) * PetscPowRealInt(x[1], n) * PetscPowRealInt(x[2], n - 1);
35*a2424aedSMatthew G. Knepley   u[2] = PetscPowRealInt(x[0], n - 1) * PetscPowRealInt(x[1], n - 1) * PetscPowRealInt(x[2], n);
36*a2424aedSMatthew G. Knepley }
37*a2424aedSMatthew G. Knepley 
38*a2424aedSMatthew G. Knepley static PetscErrorCode solenoidal_totaldeg_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
39*a2424aedSMatthew G. Knepley {
40*a2424aedSMatthew G. Knepley   const PetscInt deg = *(PetscInt *)ctx;
41*a2424aedSMatthew G. Knepley   const PetscInt n   = deg / 2 + deg % 2;
42*a2424aedSMatthew G. Knepley 
43*a2424aedSMatthew G. Knepley   solenoidal_2d(n, x, u);
44*a2424aedSMatthew G. Knepley   return PETSC_SUCCESS;
45*a2424aedSMatthew G. Knepley }
46*a2424aedSMatthew G. Knepley 
47*a2424aedSMatthew G. Knepley static PetscErrorCode solenoidal_totaldeg_3d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
48*a2424aedSMatthew G. Knepley {
49*a2424aedSMatthew G. Knepley   const PetscInt deg = *(PetscInt *)ctx;
50*a2424aedSMatthew G. Knepley   const PetscInt n   = deg / 3 + (deg % 3 ? 1 : 0);
51*a2424aedSMatthew G. Knepley 
52*a2424aedSMatthew G. Knepley   solenoidal_3d(n, x, u);
53*a2424aedSMatthew G. Knepley   return PETSC_SUCCESS;
54*a2424aedSMatthew G. Knepley }
55*a2424aedSMatthew G. Knepley 
56*a2424aedSMatthew G. Knepley // This is in P_n^{-}
57*a2424aedSMatthew G. Knepley static PetscErrorCode source_totaldeg(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
58*a2424aedSMatthew G. Knepley {
59*a2424aedSMatthew G. Knepley   const PetscInt n = *(PetscInt *)ctx;
60*a2424aedSMatthew G. Knepley 
61*a2424aedSMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) u[d] = PetscPowRealInt(x[d], n + 1);
62*a2424aedSMatthew G. Knepley   return PETSC_SUCCESS;
63*a2424aedSMatthew G. Knepley }
64*a2424aedSMatthew G. Knepley 
65*a2424aedSMatthew G. Knepley static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
66*a2424aedSMatthew G. Knepley {
67*a2424aedSMatthew G. Knepley   const PetscInt deg = (PetscInt)PetscRealPart(constants[0]);
68*a2424aedSMatthew G. Knepley   PetscScalar    p[3];
69*a2424aedSMatthew G. Knepley 
70*a2424aedSMatthew G. Knepley   if (dim == 2) PetscCallVoid(solenoidal_totaldeg_2d(dim, t, x, uOff[1] - uOff[0], p, (void *)&deg));
71*a2424aedSMatthew G. Knepley   else PetscCallVoid(solenoidal_totaldeg_3d(dim, t, x, uOff[1] - uOff[0], p, (void *)&deg));
72*a2424aedSMatthew G. Knepley   for (PetscInt c = 0; c < dim; ++c) f0[c] = -u[c] + p[c];
73*a2424aedSMatthew G. Knepley }
74*a2424aedSMatthew G. Knepley 
75*a2424aedSMatthew G. Knepley static void zero_bd(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
76*a2424aedSMatthew G. Knepley {
77*a2424aedSMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) f0[0] = 0.;
78*a2424aedSMatthew G. Knepley }
79*a2424aedSMatthew G. Knepley 
80*a2424aedSMatthew G. Knepley static void flux(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
81*a2424aedSMatthew G. Knepley {
82*a2424aedSMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) f0[0] -= u[d] * n[d];
83*a2424aedSMatthew G. Knepley }
84*a2424aedSMatthew G. Knepley 
85*a2424aedSMatthew G. Knepley static void divergence(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
86*a2424aedSMatthew G. Knepley {
87*a2424aedSMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) f0[0] += u_x[d * dim + d];
88*a2424aedSMatthew G. Knepley }
89*a2424aedSMatthew G. Knepley 
90*a2424aedSMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm)
91*a2424aedSMatthew G. Knepley {
92*a2424aedSMatthew G. Knepley   PetscFunctionBegin;
93*a2424aedSMatthew G. Knepley   PetscCall(DMCreate(comm, dm));
94*a2424aedSMatthew G. Knepley   PetscCall(DMSetType(*dm, DMPLEX));
95*a2424aedSMatthew G. Knepley   PetscCall(DMSetFromOptions(*dm));
96*a2424aedSMatthew G. Knepley   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
97*a2424aedSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
98*a2424aedSMatthew G. Knepley }
99*a2424aedSMatthew G. Knepley 
100*a2424aedSMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
101*a2424aedSMatthew G. Knepley {
102*a2424aedSMatthew G. Knepley   options->degree = -1;
103*a2424aedSMatthew G. Knepley 
104*a2424aedSMatthew G. Knepley   PetscFunctionBeginUser;
105*a2424aedSMatthew G. Knepley   PetscOptionsBegin(comm, "", "Gauss' Law Test Options", "DMPLEX");
106*a2424aedSMatthew G. Knepley   PetscOptionsEnd();
107*a2424aedSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
108*a2424aedSMatthew G. Knepley }
109*a2424aedSMatthew G. Knepley 
110*a2424aedSMatthew G. Knepley static PetscErrorCode CreateDiscretization(DM dm, AppCtx *user)
111*a2424aedSMatthew G. Knepley {
112*a2424aedSMatthew G. Knepley   PetscFE         feq, fep;
113*a2424aedSMatthew G. Knepley   PetscSpace      sp;
114*a2424aedSMatthew G. Knepley   PetscQuadrature quad, fquad;
115*a2424aedSMatthew G. Knepley   PetscDS         ds;
116*a2424aedSMatthew G. Knepley   DMLabel         label;
117*a2424aedSMatthew G. Knepley   DMPolytopeType  ct;
118*a2424aedSMatthew G. Knepley   PetscInt        dim, cStart, minDeg, maxDeg;
119*a2424aedSMatthew G. Knepley   PetscBool       isTrimmed, isSum;
120*a2424aedSMatthew G. Knepley 
121*a2424aedSMatthew G. Knepley   PetscFunctionBeginUser;
122*a2424aedSMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
123*a2424aedSMatthew G. Knepley   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
124*a2424aedSMatthew G. Knepley   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
125*a2424aedSMatthew G. Knepley   PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, dim, ct, "field_", -1, &feq));
126*a2424aedSMatthew G. Knepley   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)feq));
127*a2424aedSMatthew G. Knepley   PetscCall(PetscFEGetQuadrature(feq, &quad));
128*a2424aedSMatthew G. Knepley   PetscCall(PetscFEGetFaceQuadrature(feq, &fquad));
129*a2424aedSMatthew G. Knepley   PetscCall(PetscFEGetBasisSpace(feq, &sp));
130*a2424aedSMatthew G. Knepley   PetscCall(PetscSpaceGetDegree(sp, &minDeg, &maxDeg));
131*a2424aedSMatthew G. Knepley   PetscCall(PetscObjectTypeCompare((PetscObject)sp, PETSCSPACEPTRIMMED, &isTrimmed));
132*a2424aedSMatthew G. Knepley   PetscCall(PetscObjectTypeCompare((PetscObject)sp, PETSCSPACESUM, &isSum));
133*a2424aedSMatthew G. Knepley   if (isSum) {
134*a2424aedSMatthew G. Knepley     PetscSpace subsp, xsp, ysp;
135*a2424aedSMatthew G. Knepley     PetscInt   xdeg, ydeg;
136*a2424aedSMatthew G. Knepley     PetscBool  isTensor;
137*a2424aedSMatthew G. Knepley 
138*a2424aedSMatthew G. Knepley     PetscCall(PetscSpaceSumGetSubspace(sp, 0, &subsp));
139*a2424aedSMatthew G. Knepley     PetscCall(PetscObjectTypeCompare((PetscObject)subsp, PETSCSPACETENSOR, &isTensor));
140*a2424aedSMatthew G. Knepley     if (isTensor) {
141*a2424aedSMatthew G. Knepley       PetscCall(PetscSpaceTensorGetSubspace(subsp, 0, &xsp));
142*a2424aedSMatthew G. Knepley       PetscCall(PetscSpaceTensorGetSubspace(subsp, 1, &ysp));
143*a2424aedSMatthew G. Knepley       PetscCall(PetscSpaceGetDegree(xsp, &xdeg, NULL));
144*a2424aedSMatthew G. Knepley       PetscCall(PetscSpaceGetDegree(ysp, &ydeg, NULL));
145*a2424aedSMatthew G. Knepley       isTrimmed = xdeg != ydeg ? PETSC_TRUE : PETSC_FALSE;
146*a2424aedSMatthew G. Knepley     }
147*a2424aedSMatthew G. Knepley   }
148*a2424aedSMatthew G. Knepley   user->degree = minDeg;
149*a2424aedSMatthew G. Knepley   if (isTrimmed) user->divFree = PETSC_FALSE;
150*a2424aedSMatthew G. Knepley   else user->divFree = PETSC_TRUE;
151*a2424aedSMatthew G. Knepley   PetscCheck(!user->divFree || user->degree, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Degree 0 solution not available");
152*a2424aedSMatthew G. Knepley   PetscCall(PetscFEDestroy(&feq));
153*a2424aedSMatthew G. Knepley   PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, 1, ct, "pot_", -1, &fep));
154*a2424aedSMatthew G. Knepley   PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fep));
155*a2424aedSMatthew G. Knepley   PetscCall(PetscFESetQuadrature(fep, quad));
156*a2424aedSMatthew G. Knepley   PetscCall(PetscFESetFaceQuadrature(fep, fquad));
157*a2424aedSMatthew G. Knepley   PetscCall(PetscFEDestroy(&fep));
158*a2424aedSMatthew G. Knepley   PetscCall(DMCreateDS(dm));
159*a2424aedSMatthew G. Knepley 
160*a2424aedSMatthew G. Knepley   PetscCall(DMGetDS(dm, &ds));
161*a2424aedSMatthew G. Knepley   PetscCall(PetscDSSetResidual(ds, 0, identity, NULL));
162*a2424aedSMatthew G. Knepley   PetscCall(PetscDSSetResidual(ds, 1, divergence, NULL));
163*a2424aedSMatthew G. Knepley   if (user->divFree) {
164*a2424aedSMatthew G. Knepley     if (dim == 2) PetscCall(PetscDSSetExactSolution(ds, 0, solenoidal_totaldeg_2d, &user->degree));
165*a2424aedSMatthew G. Knepley     else PetscCall(PetscDSSetExactSolution(ds, 0, solenoidal_totaldeg_3d, &user->degree));
166*a2424aedSMatthew G. Knepley   } else {
167*a2424aedSMatthew G. Knepley     PetscCall(PetscDSSetExactSolution(ds, 0, source_totaldeg, &user->degree));
168*a2424aedSMatthew G. Knepley   }
169*a2424aedSMatthew G. Knepley   PetscCall(PetscDSSetExactSolution(ds, 1, zero, &user->degree));
170*a2424aedSMatthew G. Knepley   PetscCall(DMGetLabel(dm, "marker", &label));
171*a2424aedSMatthew G. Knepley 
172*a2424aedSMatthew G. Knepley   // TODO Can we also test the boundary residual integration?
173*a2424aedSMatthew G. Knepley   //PetscWeakForm wf;
174*a2424aedSMatthew G. Knepley   //PetscInt      bd, id = 1;
175*a2424aedSMatthew G. Knepley   //PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "boundary", label, 1, &id, 1, 0, NULL, NULL, NULL, user, &bd));
176*a2424aedSMatthew G. Knepley   //PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
177*a2424aedSMatthew G. Knepley   //PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 1, 0, 0, flux, 0, NULL));
178*a2424aedSMatthew G. Knepley 
179*a2424aedSMatthew G. Knepley   {
180*a2424aedSMatthew G. Knepley     PetscScalar constants[1];
181*a2424aedSMatthew G. Knepley 
182*a2424aedSMatthew G. Knepley     constants[0] = user->degree;
183*a2424aedSMatthew G. Knepley     PetscCall(PetscDSSetConstants(ds, 1, constants));
184*a2424aedSMatthew G. Knepley   }
185*a2424aedSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
186*a2424aedSMatthew G. Knepley }
187*a2424aedSMatthew G. Knepley 
188*a2424aedSMatthew G. Knepley int main(int argc, char **argv)
189*a2424aedSMatthew G. Knepley {
190*a2424aedSMatthew G. Knepley   DM          dm;
191*a2424aedSMatthew G. Knepley   SNES        snes;
192*a2424aedSMatthew G. Knepley   Vec         u;
193*a2424aedSMatthew G. Knepley   PetscReal   error[2], residual;
194*a2424aedSMatthew G. Knepley   PetscScalar source[2], outflow[2];
195*a2424aedSMatthew G. Knepley   AppCtx      user;
196*a2424aedSMatthew G. Knepley 
197*a2424aedSMatthew G. Knepley   PetscFunctionBeginUser;
198*a2424aedSMatthew G. Knepley   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
199*a2424aedSMatthew G. Knepley   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
200*a2424aedSMatthew G. Knepley   PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm));
201*a2424aedSMatthew G. Knepley   PetscCall(CreateDiscretization(dm, &user));
202*a2424aedSMatthew G. Knepley   PetscCall(DMGetGlobalVector(dm, &u));
203*a2424aedSMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)u, "solution"));
204*a2424aedSMatthew G. Knepley   PetscCall(DMComputeExactSolution(dm, 0., u, NULL));
205*a2424aedSMatthew G. Knepley 
206*a2424aedSMatthew G. Knepley   PetscCall(SNESCreate(PetscObjectComm((PetscObject)dm), &snes));
207*a2424aedSMatthew G. Knepley   PetscCall(SNESSetDM(snes, dm));
208*a2424aedSMatthew G. Knepley   PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, &user));
209*a2424aedSMatthew G. Knepley   PetscCall(SNESSetFromOptions(snes));
210*a2424aedSMatthew G. Knepley   PetscCall(DMSNESCheckDiscretization(snes, dm, 0., u, PETSC_DETERMINE, error));
211*a2424aedSMatthew G. Knepley   PetscCheck(PetscAbsReal(error[0]) < PETSC_SMALL, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Exact solution does not fit into FEM space: %g should be zero", (double)error[0]);
212*a2424aedSMatthew G. Knepley   if (user.divFree) {
213*a2424aedSMatthew G. Knepley     PetscCall(DMSNESCheckResidual(snes, dm, u, PETSC_DETERMINE, &residual));
214*a2424aedSMatthew G. Knepley     PetscCheck(PetscAbsReal(residual) < PETSC_SMALL, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Exact solution is not divergence-free: %g should be zero", (double)residual);
215*a2424aedSMatthew G. Knepley   } else {
216*a2424aedSMatthew G. Knepley     PetscDS ds;
217*a2424aedSMatthew G. Knepley 
218*a2424aedSMatthew G. Knepley     PetscCall(DMGetDS(dm, &ds));
219*a2424aedSMatthew G. Knepley     PetscCall(PetscDSSetObjective(ds, 1, divergence));
220*a2424aedSMatthew G. Knepley     PetscCall(DMPlexComputeIntegralFEM(dm, u, source, &user));
221*a2424aedSMatthew G. Knepley   }
222*a2424aedSMatthew G. Knepley   PetscCall(SNESDestroy(&snes));
223*a2424aedSMatthew G. Knepley 
224*a2424aedSMatthew G. Knepley   PetscBdPointFunc funcs[] = {zero_bd, flux};
225*a2424aedSMatthew G. Knepley   DMLabel          label;
226*a2424aedSMatthew G. Knepley   PetscInt         id = 1;
227*a2424aedSMatthew G. Knepley 
228*a2424aedSMatthew G. Knepley   PetscCall(DMGetLabel(dm, "marker", &label));
229*a2424aedSMatthew G. Knepley   PetscCall(DMPlexComputeBdIntegral(dm, u, label, 1, &id, funcs, outflow, &user));
230*a2424aedSMatthew G. Knepley   if (user.divFree) PetscCheck(PetscAbsScalar(outflow[1]) < PETSC_SMALL, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Outflow %g should be zero for a divergence-free field", (double)PetscRealPart(outflow[1]));
231*a2424aedSMatthew G. Knepley   else PetscCheck(PetscAbsScalar(source[1] + outflow[1]) < PETSC_SMALL, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Outflow %g should oppose source %g", (double)PetscRealPart(outflow[1]), (double)PetscRealPart(source[1]));
232*a2424aedSMatthew G. Knepley 
233*a2424aedSMatthew G. Knepley   PetscCall(DMRestoreGlobalVector(dm, &u));
234*a2424aedSMatthew G. Knepley   PetscCall(DMDestroy(&dm));
235*a2424aedSMatthew G. Knepley   PetscCall(PetscFinalize());
236*a2424aedSMatthew G. Knepley   return 0;
237*a2424aedSMatthew G. Knepley }
238*a2424aedSMatthew G. Knepley 
239*a2424aedSMatthew G. Knepley /*TEST
240*a2424aedSMatthew G. Knepley 
241*a2424aedSMatthew G. Knepley   testset:
242*a2424aedSMatthew G. Knepley     suffix: p
243*a2424aedSMatthew G. Knepley     requires: triangle ctetgen
244*a2424aedSMatthew G. Knepley     args: -dm_plex_dim {{2 3}} -dm_plex_box_faces 2,2,2
245*a2424aedSMatthew G. Knepley 
246*a2424aedSMatthew G. Knepley     test:
247*a2424aedSMatthew G. Knepley       suffix: 1
248*a2424aedSMatthew G. Knepley       args: -field_petscspace_degree 1 -pot_petscspace_degree 1
249*a2424aedSMatthew G. Knepley     test:
250*a2424aedSMatthew G. Knepley       suffix: 2
251*a2424aedSMatthew G. Knepley       args: -field_petscspace_degree 2 -pot_petscspace_degree 2
252*a2424aedSMatthew G. Knepley     test:
253*a2424aedSMatthew G. Knepley       suffix: 3
254*a2424aedSMatthew G. Knepley       args: -field_petscspace_degree 3 -pot_petscspace_degree 3
255*a2424aedSMatthew G. Knepley     test:
256*a2424aedSMatthew G. Knepley       suffix: 4
257*a2424aedSMatthew G. Knepley       args: -field_petscspace_degree 4 -pot_petscspace_degree 4
258*a2424aedSMatthew G. Knepley 
259*a2424aedSMatthew G. Knepley   testset:
260*a2424aedSMatthew G. Knepley     suffix: q
261*a2424aedSMatthew G. Knepley     args: -dm_plex_dim {{2 3}} -dm_plex_simplex 0 -dm_plex_box_faces 2,2
262*a2424aedSMatthew G. Knepley 
263*a2424aedSMatthew G. Knepley     test:
264*a2424aedSMatthew G. Knepley       suffix: 1
265*a2424aedSMatthew G. Knepley       args: -field_petscspace_degree 1 -pot_petscspace_degree 1
266*a2424aedSMatthew G. Knepley     test:
267*a2424aedSMatthew G. Knepley       suffix: 2
268*a2424aedSMatthew G. Knepley       args: -field_petscspace_degree 2 -pot_petscspace_degree 2
269*a2424aedSMatthew G. Knepley     test:
270*a2424aedSMatthew G. Knepley       suffix: 3
271*a2424aedSMatthew G. Knepley       args: -field_petscspace_degree 3 -pot_petscspace_degree 3
272*a2424aedSMatthew G. Knepley     test:
273*a2424aedSMatthew G. Knepley       suffix: 4
274*a2424aedSMatthew G. Knepley       args: -field_petscspace_degree 4 -pot_petscspace_degree 4
275*a2424aedSMatthew G. Knepley 
276*a2424aedSMatthew G. Knepley   testset:
277*a2424aedSMatthew G. Knepley     suffix: bdm
278*a2424aedSMatthew G. Knepley     requires: triangle ctetgen
279*a2424aedSMatthew G. Knepley     args: -dm_plex_dim 2 -dm_plex_box_faces 2,2
280*a2424aedSMatthew G. Knepley 
281*a2424aedSMatthew G. Knepley     test:
282*a2424aedSMatthew G. Knepley       suffix: 1
283*a2424aedSMatthew G. Knepley       args: -pot_petscspace_degree 0 -pot_petscdualspace_lagrange_continuity 0 \
284*a2424aedSMatthew G. Knepley             -field_petscspace_degree 1 -field_petscdualspace_type bdm \
285*a2424aedSMatthew G. Knepley             -field_petscfe_default_quadrature_order 2
286*a2424aedSMatthew G. Knepley 
287*a2424aedSMatthew G. Knepley   testset:
288*a2424aedSMatthew G. Knepley     suffix: rt
289*a2424aedSMatthew G. Knepley     requires: triangle ctetgen
290*a2424aedSMatthew G. Knepley     args: -dm_plex_dim 2 -dm_plex_box_faces 2,2
291*a2424aedSMatthew G. Knepley 
292*a2424aedSMatthew G. Knepley     test:
293*a2424aedSMatthew G. Knepley       suffix: 1
294*a2424aedSMatthew G. Knepley       args: -pot_petscspace_degree 0 -pot_petscdualspace_lagrange_continuity 0 \
295*a2424aedSMatthew G. Knepley             -field_petscspace_type ptrimmed \
296*a2424aedSMatthew G. Knepley             -field_petscspace_components 2 \
297*a2424aedSMatthew G. Knepley             -field_petscspace_ptrimmed_form_degree -1 \
298*a2424aedSMatthew G. Knepley             -field_petscdualspace_order 1 \
299*a2424aedSMatthew G. Knepley             -field_petscdualspace_form_degree -1 \
300*a2424aedSMatthew G. Knepley             -field_petscdualspace_lagrange_trimmed true \
301*a2424aedSMatthew G. Knepley             -field_petscfe_default_quadrature_order 2
302*a2424aedSMatthew G. Knepley 
303*a2424aedSMatthew G. Knepley   testset:
304*a2424aedSMatthew G. Knepley     suffix: rtq
305*a2424aedSMatthew G. Knepley     requires: triangle ctetgen
306*a2424aedSMatthew G. Knepley     args: -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 2,2
307*a2424aedSMatthew G. Knepley 
308*a2424aedSMatthew G. Knepley     test:
309*a2424aedSMatthew G. Knepley       suffix: 1
310*a2424aedSMatthew G. Knepley       args: -pot_petscspace_degree 0 -pot_petscdualspace_lagrange_continuity 0 \
311*a2424aedSMatthew G. Knepley             -field_petscspace_degree 1 \
312*a2424aedSMatthew G. Knepley             -field_petscspace_type sum \
313*a2424aedSMatthew G. Knepley             -field_petscspace_variables 2 \
314*a2424aedSMatthew G. Knepley             -field_petscspace_components 2 \
315*a2424aedSMatthew G. Knepley             -field_petscspace_sum_spaces 2 \
316*a2424aedSMatthew G. Knepley             -field_petscspace_sum_concatenate true \
317*a2424aedSMatthew G. Knepley             -field_sumcomp_0_petscspace_variables 2 \
318*a2424aedSMatthew G. Knepley             -field_sumcomp_0_petscspace_type tensor \
319*a2424aedSMatthew G. Knepley             -field_sumcomp_0_petscspace_tensor_spaces 2 \
320*a2424aedSMatthew G. Knepley             -field_sumcomp_0_petscspace_tensor_uniform false \
321*a2424aedSMatthew G. Knepley             -field_sumcomp_0_tensorcomp_0_petscspace_degree 1 \
322*a2424aedSMatthew G. Knepley             -field_sumcomp_0_tensorcomp_1_petscspace_degree 0 \
323*a2424aedSMatthew G. Knepley             -field_sumcomp_1_petscspace_variables 2 \
324*a2424aedSMatthew G. Knepley             -field_sumcomp_1_petscspace_type tensor \
325*a2424aedSMatthew G. Knepley             -field_sumcomp_1_petscspace_tensor_spaces 2 \
326*a2424aedSMatthew G. Knepley             -field_sumcomp_1_petscspace_tensor_uniform false \
327*a2424aedSMatthew G. Knepley             -field_sumcomp_1_tensorcomp_0_petscspace_degree 0 \
328*a2424aedSMatthew G. Knepley             -field_sumcomp_1_tensorcomp_1_petscspace_degree 1 \
329*a2424aedSMatthew G. Knepley             -field_petscdualspace_order 1 \
330*a2424aedSMatthew G. Knepley             -field_petscdualspace_form_degree -1 \
331*a2424aedSMatthew G. Knepley             -field_petscdualspace_lagrange_trimmed true \
332*a2424aedSMatthew G. Knepley             -field_petscfe_default_quadrature_order 2
333*a2424aedSMatthew G. Knepley 
334*a2424aedSMatthew G. Knepley TEST*/
335