xref: /petsc/src/dm/impls/plex/tests/ex33.c (revision bef158480efac06de457f7a665168877ab3c2fd7)
1 static char help[] = "Tests for high order geometry\n\n";
2 
3 #include <petscdmplex.h>
4 #include <petscds.h>
5 
6 typedef enum {TRANSFORM_NONE, TRANSFORM_SHEAR, TRANSFORM_ANNULUS, TRANSFORM_SHELL} Transform;
7 const char * const TransformTypes[] = {"none", "shear", "annulus", "shell", "Mesh Transform", "TRANSFORM_", NULL};
8 
9 typedef struct {
10   DM          dm;
11   PetscInt    dim;              /* The topological mesh dimension */
12   PetscBool   simplex;          /* Use simplices or hexes */
13   char        filename[PETSC_MAX_PATH_LEN]; /* Import mesh from file */
14   Transform   meshTransform;    /* Transform for initial box mesh */
15   PetscReal   *transformDataReal; /* Parameters for mesh transform */
16   PetscScalar *transformData;   /* Parameters for mesh transform */
17   PetscReal   volume;           /* Analytical volume of the mesh */
18   PetscReal   tol;              /* Tolerance for volume check */
19 } AppCtx;
20 
21 PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
22 {
23   PetscInt       n = 0, i;
24   PetscErrorCode ierr;
25 
26   PetscFunctionBegin;
27   options->dim                = 2;
28   options->simplex            = PETSC_TRUE;
29   options->filename[0]        = '\0';
30   options->meshTransform      = TRANSFORM_NONE;
31   options->transformDataReal = NULL;
32   options->transformData      = NULL;
33   options->volume             = -1.0;
34   options->tol                = PETSC_SMALL;
35 
36   ierr = PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX");CHKERRQ(ierr);
37   ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex33.c", options->dim, &options->dim, NULL,1,3);CHKERRQ(ierr);
38   ierr = PetscOptionsBool("-simplex", "Use simplices if true, otherwise hexes", "ex33.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr);
39   ierr = PetscOptionsString("-filename", "The mesh file", "ex33.c", options->filename, options->filename, sizeof(options->filename), NULL);CHKERRQ(ierr);
40   ierr = PetscOptionsEnum("-mesh_transform", "Method to transform initial box mesh <none, shear, annulus, shell>", "ex33.c", TransformTypes, (PetscEnum) options->meshTransform, (PetscEnum *) &options->meshTransform, NULL);CHKERRQ(ierr);
41   switch (options->meshTransform) {
42     case TRANSFORM_NONE: break;
43     case TRANSFORM_SHEAR:
44       n = options->dim-1;
45       ierr = PetscMalloc1(n, &options->transformDataReal);CHKERRQ(ierr);
46       for (i = 0; i < n; ++i) options->transformDataReal[i] = 1.0;
47       ierr = PetscOptionsRealArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformDataReal, &n, NULL);CHKERRQ(ierr);
48       break;
49     case TRANSFORM_ANNULUS:
50       n = 2;
51       ierr = PetscMalloc1(n, &options->transformData);CHKERRQ(ierr);
52       options->transformData[0] = 1.0;
53       options->transformData[1] = 2.0;
54       ierr = PetscOptionsScalarArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformData, &n, NULL);CHKERRQ(ierr);
55       break;
56     case TRANSFORM_SHELL:
57       n = 2;
58       ierr = PetscMalloc1(n, &options->transformData);CHKERRQ(ierr);
59       options->transformData[0] = 1.0;
60       options->transformData[1] = 2.0;
61       ierr = PetscOptionsScalarArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformData, &n, NULL);CHKERRQ(ierr);
62       break;
63     default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Unknown mesh transform %D", options->meshTransform);
64   }
65   ierr = PetscOptionsReal("-volume", "The analytical volume of the mesh", "ex33.c", options->volume, &options->volume, NULL);CHKERRQ(ierr);
66   ierr = PetscOptionsReal("-tol", "The tolerance for the volume check", "ex33.c", options->tol, &options->tol, NULL);CHKERRQ(ierr);
67   ierr = PetscOptionsEnd();
68   PetscFunctionReturn(0);
69 }
70 
71 static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux,
72                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
73                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
74                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
75 {
76   const PetscInt Nc = uOff[1] - uOff[0];
77   PetscInt       c;
78 
79   for (c = 0; c < Nc; ++c) f0[c] = u[c];
80 }
81 
82 /*
83   We would like to map the unit square to a quarter of the annulus between circles of radius 1 and 2. We start by mapping the straight sections, which
84   will correspond to the top and bottom of our square. So
85 
86     (0,0)--(1,0)  ==>  (1,0)--(2,0)      Just a shift of (1,0)
87     (0,1)--(1,1)  ==>  (0,1)--(0,2)      Switch x and y
88 
89   So it looks like we want to map each layer in y to a ray, so x is the radius and y is the angle:
90 
91     (x, y)  ==>  (x+1, \pi/2 y)                           in (r', \theta') space
92             ==>  ((x+1) cos(\pi/2 y), (x+1) sin(\pi/2 y)) in (x', y') space
93 */
94 static void f0_annulus(PetscInt dim, PetscInt Nf, PetscInt NfAux,
95                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
96                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
97                        PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar xp[])
98 {
99   const PetscReal ri = PetscRealPart(constants[0]);
100   const PetscReal ro = PetscRealPart(constants[1]);
101 
102   xp[0] = (x[0] * (ro-ri) + ri) * PetscCosReal(0.5*PETSC_PI*x[1]);
103   xp[1] = (x[0] * (ro-ri) + ri) * PetscSinReal(0.5*PETSC_PI*x[1]);
104 }
105 
106 /*
107   We would like to map the unit cube to a hemisphere of the spherical shell between balls of radius 1 and 2. We want to map the bottom surface onto the
108   lower hemisphere and the upper surface onto the top, letting z be the radius.
109 
110     (x, y)  ==>  ((z+3)/2, \pi/2 (|x| or |y|), arctan y/x)                                                  in (r', \theta', \phi') space
111             ==>  ((z+3)/2 \cos(\theta') cos(\phi'), (z+3)/2 \cos(\theta') sin(\phi'), (z+3)/2 sin(\theta')) in (x', y', z') space
112 */
113 static void f0_shell(PetscInt dim, PetscInt Nf, PetscInt NfAux,
114                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
115                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
116                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar xp[])
117 {
118   const PetscReal pi4    = PETSC_PI/4.0;
119   const PetscReal ri     = PetscRealPart(constants[0]);
120   const PetscReal ro     = PetscRealPart(constants[1]);
121   const PetscReal rp     = (x[2]+1) * 0.5*(ro-ri) + ri;
122   const PetscReal phip   = PetscAtan2Real(x[1], x[0]);
123   const PetscReal thetap = 0.5*PETSC_PI * (1.0 - ((((phip <= pi4) && (phip >= -pi4)) || ((phip >= 3.0*pi4) || (phip <= -3.0*pi4))) ? PetscAbsReal(x[0]) : PetscAbsReal(x[1])));
124 
125   xp[0] = rp * PetscCosReal(thetap) * PetscCosReal(phip);
126   xp[1] = rp * PetscCosReal(thetap) * PetscSinReal(phip);
127   xp[2] = rp * PetscSinReal(thetap);
128 }
129 
130 static PetscErrorCode DMCreateCoordinateDisc(DM dm, PetscBool isSimplex)
131 {
132   DM              cdm;
133   PetscFE         fe;
134   PetscSpace      basis;
135   PetscDualSpace  dual;
136   PetscInt        dim, dE;
137   PetscErrorCode  ierr;
138 
139   PetscFunctionBegin;
140   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
141   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
142   ierr = DMGetCoordinateDim(dm, &dE);CHKERRQ(ierr);
143   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, dE, isSimplex, "geom_", -1, &fe);CHKERRQ(ierr);
144   ierr = DMSetField(cdm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr);
145   ierr = PetscFEGetBasisSpace(fe, &basis);CHKERRQ(ierr);
146   ierr = PetscFEGetDualSpace(fe, &dual);CHKERRQ(ierr);
147   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
148   ierr = DMCreateDS(cdm);CHKERRQ(ierr);
149 
150   {
151     DM      cdmLinear;
152     PetscFE feLinear;
153     Mat     In;
154     Vec     coordinates, coordinatesNew;
155 
156     /* Have to clear out previous coordinate structures */
157     ierr = DMGetCoordinates(dm, &coordinates);CHKERRQ(ierr);
158     ierr = DMSetCoordinateField(dm, NULL);CHKERRQ(ierr);
159     ierr = DMSetLocalSection(cdm, NULL);CHKERRQ(ierr);
160     ierr = DMSetGlobalSection(cdm, NULL);CHKERRQ(ierr);
161     ierr = DMSetSectionSF(cdm, NULL);CHKERRQ(ierr);
162     /* Project from vertex coordinates to chosen space */
163     ierr = DMClone(cdm, &cdmLinear);CHKERRQ(ierr);
164     ierr = PetscFECreateLagrange(PETSC_COMM_SELF, dim, dE, isSimplex, 1, -1, &feLinear);CHKERRQ(ierr);
165     ierr = DMSetField(cdmLinear, 0, NULL, (PetscObject) feLinear);CHKERRQ(ierr);
166     ierr = PetscFEDestroy(&feLinear);CHKERRQ(ierr);
167     ierr = DMCreateDS(cdmLinear);CHKERRQ(ierr);
168     ierr = DMCreateInterpolation(cdmLinear, cdm, &In, NULL);CHKERRQ(ierr);
169     ierr = DMCreateGlobalVector(cdm, &coordinatesNew);CHKERRQ(ierr);
170     ierr = MatInterpolate(In, coordinates, coordinatesNew);CHKERRQ(ierr);
171     ierr = MatDestroy(&In);CHKERRQ(ierr);
172     ierr = DMSetCoordinates(dm, coordinatesNew);CHKERRQ(ierr);
173     ierr = VecViewFromOptions(coordinatesNew, NULL, "-coords_view");CHKERRQ(ierr);
174     ierr = VecDestroy(&coordinatesNew);CHKERRQ(ierr);
175     ierr = DMDestroy(&cdmLinear);CHKERRQ(ierr);
176   }
177   PetscFunctionReturn(0);
178 }
179 
180 PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm)
181 {
182   PetscInt       dim      = ctx->dim;
183   PetscBool      simplex  = ctx->simplex;
184   const char    *filename = ctx->filename;
185   DM             cdm;
186   PetscDS        cds;
187   size_t         len;
188   PetscMPIInt    rank;
189   PetscErrorCode ierr;
190 
191   PetscFunctionBegin;
192   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
193   ierr = PetscStrlen(filename, &len);CHKERRQ(ierr);
194   if (len) {
195     ierr = DMPlexCreateFromFile(comm, filename, PETSC_FALSE, dm);CHKERRQ(ierr);
196     ierr = DMGetDimension(*dm, &dim);CHKERRQ(ierr);
197   } else {
198     ierr = DMPlexCreateBoxMesh(comm, dim, simplex, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
199   }
200   {
201     DM               pdm = NULL;
202     PetscPartitioner part;
203 
204     ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr);
205     ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
206     /* Distribute mesh over processes */
207     ierr = DMPlexDistribute(*dm, 0, NULL, &pdm);CHKERRQ(ierr);
208     if (pdm) {
209       ierr = DMDestroy(dm);CHKERRQ(ierr);
210       *dm  = pdm;
211     }
212   }
213   ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr);
214   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
215   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
216   ierr = DMCreateCoordinateDisc(*dm, simplex);CHKERRQ(ierr);
217   switch (ctx->meshTransform) {
218     case TRANSFORM_NONE:
219       ierr = DMPlexRemapGeometry(*dm, 0.0, identity);CHKERRQ(ierr);
220       break;
221     case TRANSFORM_SHEAR:
222       ierr = DMPlexShearGeometry(*dm, DM_X, ctx->transformDataReal);CHKERRQ(ierr);
223       break;
224     case TRANSFORM_ANNULUS:
225       ierr = DMGetCoordinateDM(*dm, &cdm);CHKERRQ(ierr);
226       ierr = DMGetDS(cdm, &cds);CHKERRQ(ierr);
227       ierr = PetscDSSetConstants(cds, 2, ctx->transformData);CHKERRQ(ierr);
228       ierr = DMPlexRemapGeometry(*dm, 0.0, f0_annulus);CHKERRQ(ierr);
229       break;
230     case TRANSFORM_SHELL:
231       ierr = DMGetCoordinateDM(*dm, &cdm);CHKERRQ(ierr);
232       ierr = DMGetDS(cdm, &cds);CHKERRQ(ierr);
233       ierr = PetscDSSetConstants(cds, 2, ctx->transformData);CHKERRQ(ierr);
234       ierr = DMPlexRemapGeometry(*dm, 0.0, f0_shell);CHKERRQ(ierr);
235       break;
236     default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Unknown mesh transform %D", ctx->meshTransform);
237   }
238   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
239   ctx->dm = *dm;
240   PetscFunctionReturn(0);
241 }
242 
243 static void volume(PetscInt dim, PetscInt Nf, PetscInt NfAux,
244                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
245                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
246                    PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar vol[])
247 {
248   vol[0] = 1.;
249 }
250 
251 static PetscErrorCode CreateDiscretization(DM dm, AppCtx *ctx)
252 {
253   PetscDS        ds;
254   PetscFE        fe;
255   PetscErrorCode ierr;
256 
257   PetscFunctionBeginUser;
258   ierr = PetscFECreateDefault(PETSC_COMM_SELF, ctx->dim, 1, ctx->simplex, NULL, PETSC_DETERMINE, &fe);CHKERRQ(ierr);
259   ierr = PetscFESetName(fe, "scalar");CHKERRQ(ierr);
260   ierr = DMAddField(dm, NULL, (PetscObject) fe);
261   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
262   ierr = DMCreateDS(dm);CHKERRQ(ierr);
263   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
264   ierr = PetscDSSetObjective(ds, 0, volume);CHKERRQ(ierr);
265   PetscFunctionReturn(0);
266 }
267 
268 static PetscErrorCode CheckVolume(DM dm, AppCtx *ctx)
269 {
270   Vec            u;
271   PetscScalar    result;
272   PetscReal      vol, tol = ctx->tol;
273   PetscErrorCode ierr;
274 
275   PetscFunctionBeginUser;
276   ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr);
277   ierr = DMPlexComputeIntegralFEM(dm, u, &result, ctx);CHKERRQ(ierr);
278   vol  = PetscRealPart(result);
279   ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr);
280   ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Volume: %g\n", (double) vol);CHKERRQ(ierr);
281   if (ctx->volume > 0.0 && PetscAbsReal(ctx->volume - vol) > tol) {
282     SETERRQ4(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "Calculated volume %g != %g actual volume (error %g > %g tol)", (double) vol, (double) ctx->volume, (double) PetscAbsReal(ctx->volume - vol), (double) tol);
283   }
284   PetscFunctionReturn(0);
285 }
286 
287 int main(int argc, char **argv)
288 {
289   AppCtx         user;
290   PetscErrorCode ierr;
291 
292   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
293   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
294   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &user.dm);CHKERRQ(ierr);
295   ierr = CreateDiscretization(user.dm, &user);CHKERRQ(ierr);
296   ierr = CheckVolume(user.dm, &user);CHKERRQ(ierr);
297   ierr = DMDestroy(&user.dm);CHKERRQ(ierr);
298   ierr = PetscFree(user.transformDataReal);CHKERRQ(ierr);
299   ierr = PetscFree(user.transformData);CHKERRQ(ierr);
300   ierr = PetscFinalize();
301   return ierr;
302 }
303 
304 /*TEST
305 
306   test:
307     suffix: square_0
308     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1.,-1. -dm_plex_box_upper 1.,1. -geom_petscspace_degree 1 -volume 4.
309 
310   test:
311     suffix: square_1
312     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1.,-1. -dm_plex_box_upper 1.,1. -geom_petscspace_degree 2 -volume 4.
313 
314   test:
315     suffix: square_2
316     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1.,-1. -dm_plex_box_upper 1.,1. -dm_refine 1 -geom_petscspace_degree 1 -volume 4.
317 
318   test:
319     suffix: square_3
320     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1.,-1. -dm_plex_box_upper 1.,1. -dm_refine 1 -geom_petscspace_degree 2 -volume 4.
321 
322   test:
323     suffix: cube_0
324     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -geom_petscspace_degree 1 -volume 8.
325 
326   test:
327     suffix: cube_1
328     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -geom_petscspace_degree 2 -volume 8.
329 
330   test:
331     suffix: cube_2
332     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -dm_refine 1 -geom_petscspace_degree 1 -volume 8.
333 
334   test:
335     suffix: cube_3
336     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -dm_refine 1 -geom_petscspace_degree 2 -volume 8.
337 
338   test:
339     suffix: shear_0
340     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1.,-1. -dm_plex_box_upper 1.,1. -geom_petscspace_degree 1 -mesh_transform shear -transform_data 3.0 -volume 4.
341 
342   test:
343     suffix: shear_1
344     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1.,-1. -dm_plex_box_upper 1.,1. -geom_petscspace_degree 2 -mesh_transform shear -transform_data 3.0 -volume 4.
345 
346   test:
347     suffix: shear_2
348     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1.,-1. -dm_plex_box_upper 1.,1. -dm_refine 1 -geom_petscspace_degree 1 -mesh_transform shear -transform_data 3.0 -volume 4.
349 
350   test:
351     suffix: shear_3
352     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1.,-1. -dm_plex_box_upper 1.,1. -dm_refine 1 -geom_petscspace_degree 2 -mesh_transform shear -transform_data 3.0 -volume 4.
353 
354   test:
355     suffix: shear_4
356     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -geom_petscspace_degree 1 -mesh_transform shear -transform_data 3.0 -volume 8.
357 
358   test:
359     suffix: shear_5
360     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -geom_petscspace_degree 2 -mesh_transform shear -transform_data 3.0 -volume 8.
361 
362   test:
363     suffix: shear_6
364     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1,-1. -dm_plex_box_upper 1.,1.,1. -dm_refine 1 -geom_petscspace_degree 1 -mesh_transform shear -transform_data 3.0,4.0 -volume 8.
365 
366   test:
367     suffix: shear_7
368     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -dm_refine 1 -geom_petscspace_degree 2 -mesh_transform shear -transform_data 3.0,4.0 -volume 8.
369 
370   test:
371     # Area: (a+b)/2 h = 3/\sqrt{2} (sqrt{2} - 1/\sqrt{2}) = 3/2
372     suffix: annulus_0
373     requires: double
374     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -geom_petscspace_degree 1 -mesh_transform annulus -volume 1.5
375 
376   test:
377     # Area: 3/4 \pi = 2.3562
378     suffix: annulus_1
379     requires: double
380     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_refine 3 -geom_petscspace_degree 1 -mesh_transform annulus -volume 2.35619449019235 -tol .016
381 
382   test:
383     # Area: 3/4 \pi = 2.3562
384     suffix: annulus_2
385     requires: double
386     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_refine 3 -geom_petscspace_degree 2 -mesh_transform annulus -volume 2.35619449019235 -tol .0038
387 
388   test:
389     # Area: 3/4 \pi = 2.3562
390     suffix: annulus_3
391     requires: double
392     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_refine 3 -geom_petscspace_degree 3 -mesh_transform annulus -volume 2.35619449019235 -tol 2.2e-6
393 
394   test:
395     # Area: 3/4 \pi = 2.3562
396     suffix: annulus_4
397     requires: double
398     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_refine 2 -geom_petscspace_degree 2 -petscfe_default_quadrature_order 2 -mesh_transform annulus -volume 2.35619449019235 -tol .00012
399 
400   test:
401     # Area: 3/4 \pi = 2.3562
402     suffix: annulus_5
403     requires: double
404     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_refine 2 -geom_petscspace_degree 3 -petscfe_default_quadrature_order 3 -mesh_transform annulus -volume 2.35619449019235 -tol 1.2e-7
405 
406   test:
407     suffix: shell_0
408     requires: double
409     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -dm_refine 1 -geom_petscspace_degree 1 -petscfe_default_quadrature_order 1 -mesh_transform shell -volume 5.633164922 -tol 1.0e-7
410 
411   test:
412     # Volume: 4/3 \pi (8 - 1)/2 = 14/3 \pi = 14.66076571675238
413     suffix: shell_1
414     requires: double
415     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -dm_refine 2 -geom_petscspace_degree 1 -petscfe_default_quadrature_order 1 -mesh_transform shell -volume 14.66076571675238 -tol 3.1
416 
417   test:
418     # Volume: 4/3 \pi (8 - 1)/2 = 14/3 \pi = 14.66076571675238
419     suffix: shell_2
420     requires: double
421     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -dm_refine 2 -geom_petscspace_degree 2 -petscfe_default_quadrature_order 2 -mesh_transform shell -volume 14.66076571675238 -tol .1
422 
423   test:
424     # Volume: 4/3 \pi (8 - 1)/2 = 14/3 \pi = 14.66076571675238
425     suffix: shell_3
426     requires: double
427     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -dm_refine 2 -geom_petscspace_degree 3 -petscfe_default_quadrature_order 3 -mesh_transform shell -volume 14.66076571675238 -tol .02
428 
429   test:
430     # Volume: 4/3 \pi (8 - 1)/2 = 14/3 \pi = 14.66076571675238
431     suffix: shell_4
432     requires: double
433     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -dm_refine 2 -geom_petscspace_degree 4 -petscfe_default_quadrature_order 4 -mesh_transform shell -volume 14.66076571675238 -tol .006
434 
435 TEST*/
436