xref: /petsc/src/dm/impls/plex/tests/ex33.c (revision 4302210d98c281cf4e8bc0fee406cf222e373a4c)
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   char        filename[PETSC_MAX_PATH_LEN]; /* Import mesh from file */
11   Transform   meshTransform;      /* Transform for initial box mesh */
12   PetscReal   *transformDataReal; /* Parameters for mesh transform */
13   PetscScalar *transformData;     /* Parameters for mesh transform */
14   PetscReal   volume;             /* Analytical volume of the mesh */
15   PetscReal   tol;                /* Tolerance for volume check */
16 } AppCtx;
17 
18 PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
19 {
20   PetscInt       n = 0, i;
21   PetscErrorCode ierr;
22 
23   PetscFunctionBegin;
24   options->filename[0]       = '\0';
25   options->meshTransform     = TRANSFORM_NONE;
26   options->transformDataReal = NULL;
27   options->transformData     = NULL;
28   options->volume            = -1.0;
29   options->tol               = PETSC_SMALL;
30 
31   ierr = PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX");CHKERRQ(ierr);
32   ierr = PetscOptionsString("-filename", "The mesh file", "ex33.c", options->filename, options->filename, sizeof(options->filename), NULL);CHKERRQ(ierr);
33   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);
34   switch (options->meshTransform) {
35     case TRANSFORM_NONE: break;
36     case TRANSFORM_SHEAR:
37       n = 2;
38       ierr = PetscMalloc1(n, &options->transformDataReal);CHKERRQ(ierr);
39       for (i = 0; i < n; ++i) options->transformDataReal[i] = 1.0;
40       ierr = PetscOptionsRealArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformDataReal, &n, NULL);CHKERRQ(ierr);
41       break;
42     case TRANSFORM_ANNULUS:
43       n = 2;
44       ierr = PetscMalloc1(n, &options->transformData);CHKERRQ(ierr);
45       options->transformData[0] = 1.0;
46       options->transformData[1] = 2.0;
47       ierr = PetscOptionsScalarArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformData, &n, NULL);CHKERRQ(ierr);
48       break;
49     case TRANSFORM_SHELL:
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     default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Unknown mesh transform %D", options->meshTransform);
57   }
58   ierr = PetscOptionsReal("-volume", "The analytical volume of the mesh", "ex33.c", options->volume, &options->volume, NULL);CHKERRQ(ierr);
59   ierr = PetscOptionsReal("-tol", "The tolerance for the volume check", "ex33.c", options->tol, &options->tol, NULL);CHKERRQ(ierr);
60   ierr = PetscOptionsEnd();
61   PetscFunctionReturn(0);
62 }
63 
64 static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux,
65                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
66                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
67                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
68 {
69   const PetscInt Nc = uOff[1] - uOff[0];
70   PetscInt       c;
71 
72   for (c = 0; c < Nc; ++c) f0[c] = u[c];
73 }
74 
75 /*
76   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
77   will correspond to the top and bottom of our square. So
78 
79     (0,0)--(1,0)  ==>  (1,0)--(2,0)      Just a shift of (1,0)
80     (0,1)--(1,1)  ==>  (0,1)--(0,2)      Switch x and y
81 
82   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:
83 
84     (x, y)  ==>  (x+1, \pi/2 y)                           in (r', \theta') space
85             ==>  ((x+1) cos(\pi/2 y), (x+1) sin(\pi/2 y)) in (x', y') space
86 */
87 static void f0_annulus(PetscInt dim, PetscInt Nf, PetscInt NfAux,
88                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
89                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
90                        PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar xp[])
91 {
92   const PetscReal ri = PetscRealPart(constants[0]);
93   const PetscReal ro = PetscRealPart(constants[1]);
94 
95   xp[0] = (x[0] * (ro-ri) + ri) * PetscCosReal(0.5*PETSC_PI*x[1]);
96   xp[1] = (x[0] * (ro-ri) + ri) * PetscSinReal(0.5*PETSC_PI*x[1]);
97 }
98 
99 /*
100   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
101   lower hemisphere and the upper surface onto the top, letting z be the radius.
102 
103     (x, y)  ==>  ((z+3)/2, \pi/2 (|x| or |y|), arctan y/x)                                                  in (r', \theta', \phi') space
104             ==>  ((z+3)/2 \cos(\theta') cos(\phi'), (z+3)/2 \cos(\theta') sin(\phi'), (z+3)/2 sin(\theta')) in (x', y', z') space
105 */
106 static void f0_shell(PetscInt dim, PetscInt Nf, PetscInt NfAux,
107                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
108                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
109                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar xp[])
110 {
111   const PetscReal pi4    = PETSC_PI/4.0;
112   const PetscReal ri     = PetscRealPart(constants[0]);
113   const PetscReal ro     = PetscRealPart(constants[1]);
114   const PetscReal rp     = (x[2]+1) * 0.5*(ro-ri) + ri;
115   const PetscReal phip   = PetscAtan2Real(x[1], x[0]);
116   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])));
117 
118   xp[0] = rp * PetscCosReal(thetap) * PetscCosReal(phip);
119   xp[1] = rp * PetscCosReal(thetap) * PetscSinReal(phip);
120   xp[2] = rp * PetscSinReal(thetap);
121 }
122 
123 static PetscErrorCode DMCreateCoordinateDisc(DM dm)
124 {
125   DM             cdm;
126   PetscFE        fe;
127   DMPolytopeType ct;
128   PetscInt       dim, dE, cStart;
129   PetscBool      simplex;
130   PetscErrorCode ierr;
131 
132   PetscFunctionBegin;
133   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
134   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
135   ierr = DMGetCoordinateDim(dm, &dE);CHKERRQ(ierr);
136   ierr = DMPlexGetHeightStratum(cdm, 0, &cStart, NULL);CHKERRQ(ierr);
137   ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr);
138   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
139   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, dE, simplex, "geom_", -1, &fe);CHKERRQ(ierr);
140   ierr = DMProjectCoordinates(dm, fe);CHKERRQ(ierr);
141   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
142   PetscFunctionReturn(0);
143 }
144 
145 PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm)
146 {
147   const char    *filename = ctx->filename;
148   DM             cdm;
149   PetscDS        cds;
150   size_t         len;
151   PetscMPIInt    rank;
152   PetscErrorCode ierr;
153 
154   PetscFunctionBegin;
155   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
156   ierr = PetscStrlen(filename, &len);CHKERRQ(ierr);
157   if (len) {ierr = DMPlexCreateFromFile(comm, filename, PETSC_TRUE, dm);CHKERRQ(ierr);}
158   else     {ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_FALSE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);}
159   ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr);
160   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
161   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
162   if (!len) {ierr = DMCreateCoordinateDisc(*dm);CHKERRQ(ierr);}
163   switch (ctx->meshTransform) {
164     case TRANSFORM_NONE:
165       ierr = DMPlexRemapGeometry(*dm, 0.0, identity);CHKERRQ(ierr);
166       break;
167     case TRANSFORM_SHEAR:
168       ierr = DMPlexShearGeometry(*dm, DM_X, ctx->transformDataReal);CHKERRQ(ierr);
169       break;
170     case TRANSFORM_ANNULUS:
171       ierr = DMGetCoordinateDM(*dm, &cdm);CHKERRQ(ierr);
172       ierr = DMGetDS(cdm, &cds);CHKERRQ(ierr);
173       ierr = PetscDSSetConstants(cds, 2, ctx->transformData);CHKERRQ(ierr);
174       ierr = DMPlexRemapGeometry(*dm, 0.0, f0_annulus);CHKERRQ(ierr);
175       break;
176     case TRANSFORM_SHELL:
177       ierr = DMGetCoordinateDM(*dm, &cdm);CHKERRQ(ierr);
178       ierr = DMGetDS(cdm, &cds);CHKERRQ(ierr);
179       ierr = PetscDSSetConstants(cds, 2, ctx->transformData);CHKERRQ(ierr);
180       ierr = DMPlexRemapGeometry(*dm, 0.0, f0_shell);CHKERRQ(ierr);
181       break;
182     default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Unknown mesh transform %D", ctx->meshTransform);
183   }
184   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
185   PetscFunctionReturn(0);
186 }
187 
188 static void volume(PetscInt dim, PetscInt Nf, PetscInt NfAux,
189                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
190                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
191                    PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar vol[])
192 {
193   vol[0] = 1.;
194 }
195 
196 static PetscErrorCode CreateDiscretization(DM dm, AppCtx *ctx)
197 {
198   PetscDS        ds;
199   PetscFE        fe;
200   DMPolytopeType ct;
201   PetscInt       dim, cStart;
202   PetscBool      simplex;
203   PetscErrorCode ierr;
204 
205   PetscFunctionBeginUser;
206   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
207   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr);
208   ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr);
209   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
210   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, PETSC_DETERMINE, &fe);CHKERRQ(ierr);
211   ierr = PetscFESetName(fe, "scalar");CHKERRQ(ierr);
212   ierr = DMAddField(dm, NULL, (PetscObject) fe);
213   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
214   ierr = DMCreateDS(dm);CHKERRQ(ierr);
215   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
216   ierr = PetscDSSetObjective(ds, 0, volume);CHKERRQ(ierr);
217   PetscFunctionReturn(0);
218 }
219 
220 static PetscErrorCode CheckVolume(DM dm, AppCtx *ctx)
221 {
222   Vec            u;
223   PetscScalar    result;
224   PetscReal      vol, tol = ctx->tol;
225   PetscErrorCode ierr;
226 
227   PetscFunctionBeginUser;
228   ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr);
229   ierr = DMPlexComputeIntegralFEM(dm, u, &result, ctx);CHKERRQ(ierr);
230   vol  = PetscRealPart(result);
231   ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr);
232   ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Volume: %g\n", (double) vol);CHKERRQ(ierr);
233   if (ctx->volume > 0.0 && PetscAbsReal(ctx->volume - vol) > tol) {
234     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);
235   }
236   PetscFunctionReturn(0);
237 }
238 
239 int main(int argc, char **argv)
240 {
241   DM             dm;
242   AppCtx         user;
243   PetscErrorCode ierr;
244 
245   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
246   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
247   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
248   ierr = CreateDiscretization(dm, &user);CHKERRQ(ierr);
249   ierr = CheckVolume(dm, &user);CHKERRQ(ierr);
250   ierr = DMDestroy(&dm);CHKERRQ(ierr);
251   ierr = PetscFree(user.transformDataReal);CHKERRQ(ierr);
252   ierr = PetscFree(user.transformData);CHKERRQ(ierr);
253   ierr = PetscFinalize();
254   return ierr;
255 }
256 
257 /*TEST
258 
259   test:
260     suffix: square_0
261     args: -dm_plex_box_faces 1,1 -dm_plex_box_lower -1.,-1. -dm_plex_box_upper 1.,1. -geom_petscspace_degree 1 -volume 4.
262 
263   test:
264     suffix: square_1
265     args: -dm_plex_box_faces 1,1 -dm_plex_box_lower -1.,-1. -dm_plex_box_upper 1.,1. -geom_petscspace_degree 2 -volume 4.
266 
267   test:
268     suffix: square_2
269     args: -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.
270 
271   test:
272     suffix: square_3
273     args: -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.
274 
275   test:
276     suffix: cube_0
277     args: -dm_plex_box_dim 3 -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.
278 
279   test:
280     suffix: cube_1
281     args: -dm_plex_box_dim 3 -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.
282 
283   test:
284     suffix: cube_2
285     args: -dm_plex_box_dim 3 -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.
286 
287   test:
288     suffix: cube_3
289     args: -dm_plex_box_dim 3 -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.
290 
291   test:
292     suffix: shear_0
293     args: -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.
294 
295   test:
296     suffix: shear_1
297     args: -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.
298 
299   test:
300     suffix: shear_2
301     args: -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.
302 
303   test:
304     suffix: shear_3
305     args: -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.
306 
307   test:
308     suffix: shear_4
309     args: -dm_plex_box_dim 3 -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.
310 
311   test:
312     suffix: shear_5
313     args: -dm_plex_box_dim 3 -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.
314 
315   test:
316     suffix: shear_6
317     args: -dm_plex_box_dim 3 -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.
318 
319   test:
320     suffix: shear_7
321     args: -dm_plex_box_dim 3 -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.
322 
323   test:
324     # Area: (a+b)/2 h = 3/\sqrt{2} (sqrt{2} - 1/\sqrt{2}) = 3/2
325     suffix: annulus_0
326     requires: double
327     args: -dm_plex_box_faces 1,1 -geom_petscspace_degree 1 -mesh_transform annulus -volume 1.5
328 
329   test:
330     # Area: 3/4 \pi = 2.3562
331     suffix: annulus_1
332     requires: double
333     args: -dm_plex_box_faces 1,1 -dm_refine 3 -geom_petscspace_degree 1 -mesh_transform annulus -volume 2.35619449019235 -tol .016
334 
335   test:
336     # Area: 3/4 \pi = 2.3562
337     suffix: annulus_2
338     requires: double
339     args: -dm_plex_box_faces 1,1 -dm_refine 3 -geom_petscspace_degree 2 -mesh_transform annulus -volume 2.35619449019235 -tol .0038
340 
341   test:
342     # Area: 3/4 \pi = 2.3562
343     suffix: annulus_3
344     requires: double
345     args: -dm_plex_box_faces 1,1 -dm_refine 3 -geom_petscspace_degree 3 -mesh_transform annulus -volume 2.35619449019235 -tol 2.2e-6
346 
347   test:
348     # Area: 3/4 \pi = 2.3562
349     suffix: annulus_4
350     requires: double
351     args: -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
352 
353   test:
354     # Area: 3/4 \pi = 2.3562
355     suffix: annulus_5
356     requires: double
357     args: -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
358 
359   test:
360     suffix: shell_0
361     requires: double
362     args: -dm_plex_box_dim 3 -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
363 
364   test:
365     # Volume: 4/3 \pi (8 - 1)/2 = 14/3 \pi = 14.66076571675238
366     suffix: shell_1
367     requires: double
368     args: -dm_plex_box_dim 3 -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
369 
370   test:
371     # Volume: 4/3 \pi (8 - 1)/2 = 14/3 \pi = 14.66076571675238
372     suffix: shell_2
373     requires: double
374     args: -dm_plex_box_dim 3 -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
375 
376   test:
377     # Volume: 4/3 \pi (8 - 1)/2 = 14/3 \pi = 14.66076571675238
378     suffix: shell_3
379     requires: double
380     args: -dm_plex_box_dim 3 -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
381 
382   test:
383     # Volume: 4/3 \pi (8 - 1)/2 = 14/3 \pi = 14.66076571675238
384     suffix: shell_4
385     requires: double
386     args: -dm_plex_box_dim 3 -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
387 
388   test:
389     # Volume: 1.0
390     suffix: gmsh_q2
391     requires: double
392     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/quads-q2.msh -dm_plex_gmsh_project -volume 1.0 -tol 1e-6
393 
394   test:
395     # Volume: 1.0
396     suffix: gmsh_q3
397     requires: double
398     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/quads-q3.msh -dm_plex_gmsh_project -volume 1.0 -tol 1e-6
399 
400 TEST*/
401