xref: /petsc/src/dm/impls/plex/tests/ex33.c (revision 4e278199b78715991f5c71ebbd945c1489263e6c)
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_FLARE, TRANSFORM_ANNULUS, TRANSFORM_SHELL} Transform;
7 const char * const TransformTypes[] = {"none", "shear", "flare", "annulus", "shell", "Mesh Transform", "TRANSFORM_", NULL};
8 
9 typedef struct {
10   PetscBool   coordSpace;         /* Flag to create coordinate space */
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->coordSpace        = PETSC_TRUE;
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 = PetscOptionsBool("-coord_space", "Flag to create a coordinate space", "ex33.c", options->coordSpace, &options->coordSpace, 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_FLARE:
43       n = 4;
44       ierr = PetscMalloc1(n, &options->transformData);CHKERRQ(ierr);
45       for (i = 0; i < n; ++i) options->transformData[i] = 1.0;
46       options->transformData[0] = (PetscScalar) 0;
47       ierr = PetscOptionsScalarArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformData, &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 /* Flare applies the transformation, assuming we fix x_f,
83 
84    x_i = x_i * alpha_i x_f
85 */
86 static void f0_flare(PetscInt dim, PetscInt Nf, PetscInt NfAux,
87                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
88                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
89                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar coords[])
90 {
91   const PetscInt Nc = uOff[1]-uOff[0];
92   const PetscInt cf = (PetscInt) PetscRealPart(constants[0]);
93   PetscInt       c;
94 
95   for (c = 0; c < Nc; ++c) {
96     coords[c] = u[c] * (c == cf ? 1.0 : constants[c+1]*u[cf]);
97   }
98 }
99 
100 /*
101   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
102   will correspond to the top and bottom of our square. So
103 
104     (0,0)--(1,0)  ==>  (1,0)--(2,0)      Just a shift of (1,0)
105     (0,1)--(1,1)  ==>  (0,1)--(0,2)      Switch x and y
106 
107   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:
108 
109     (x, y)  ==>  (x+1, \pi/2 y)                           in (r', \theta') space
110             ==>  ((x+1) cos(\pi/2 y), (x+1) sin(\pi/2 y)) in (x', y') space
111 */
112 static void f0_annulus(PetscInt dim, PetscInt Nf, PetscInt NfAux,
113                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
114                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
115                        PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar xp[])
116 {
117   const PetscReal ri = PetscRealPart(constants[0]);
118   const PetscReal ro = PetscRealPart(constants[1]);
119 
120   xp[0] = (x[0] * (ro-ri) + ri) * PetscCosReal(0.5*PETSC_PI*x[1]);
121   xp[1] = (x[0] * (ro-ri) + ri) * PetscSinReal(0.5*PETSC_PI*x[1]);
122 }
123 
124 /*
125   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
126   lower hemisphere and the upper surface onto the top, letting z be the radius.
127 
128     (x, y)  ==>  ((z+3)/2, \pi/2 (|x| or |y|), arctan y/x)                                                  in (r', \theta', \phi') space
129             ==>  ((z+3)/2 \cos(\theta') cos(\phi'), (z+3)/2 \cos(\theta') sin(\phi'), (z+3)/2 sin(\theta')) in (x', y', z') space
130 */
131 static void f0_shell(PetscInt dim, PetscInt Nf, PetscInt NfAux,
132                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
133                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
134                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar xp[])
135 {
136   const PetscReal pi4    = PETSC_PI/4.0;
137   const PetscReal ri     = PetscRealPart(constants[0]);
138   const PetscReal ro     = PetscRealPart(constants[1]);
139   const PetscReal rp     = (x[2]+1) * 0.5*(ro-ri) + ri;
140   const PetscReal phip   = PetscAtan2Real(x[1], x[0]);
141   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])));
142 
143   xp[0] = rp * PetscCosReal(thetap) * PetscCosReal(phip);
144   xp[1] = rp * PetscCosReal(thetap) * PetscSinReal(phip);
145   xp[2] = rp * PetscSinReal(thetap);
146 }
147 
148 static PetscErrorCode DMCreateCoordinateDisc(DM dm)
149 {
150   DM             cdm;
151   PetscFE        fe;
152   DMPolytopeType ct;
153   PetscInt       dim, dE, cStart;
154   PetscBool      simplex;
155   PetscErrorCode ierr;
156 
157   PetscFunctionBegin;
158   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
159   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
160   ierr = DMGetCoordinateDim(dm, &dE);CHKERRQ(ierr);
161   ierr = DMPlexGetHeightStratum(cdm, 0, &cStart, NULL);CHKERRQ(ierr);
162   ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr);
163   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
164   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, dE, simplex, "dm_coord_", -1, &fe);CHKERRQ(ierr);
165   ierr = DMProjectCoordinates(dm, fe);CHKERRQ(ierr);
166   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
167   PetscFunctionReturn(0);
168 }
169 
170 PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm)
171 {
172   DM             cdm;
173   PetscDS        cds;
174   PetscErrorCode ierr;
175 
176   PetscFunctionBegin;
177   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
178   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
179   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
180 
181   if (ctx->coordSpace) {ierr = DMCreateCoordinateDisc(*dm);CHKERRQ(ierr);}
182   switch (ctx->meshTransform) {
183     case TRANSFORM_NONE:
184       ierr = DMPlexRemapGeometry(*dm, 0.0, identity);CHKERRQ(ierr);
185       break;
186     case TRANSFORM_SHEAR:
187       ierr = DMPlexShearGeometry(*dm, DM_X, ctx->transformDataReal);CHKERRQ(ierr);
188       break;
189     case TRANSFORM_FLARE:
190       ierr = DMGetCoordinateDM(*dm, &cdm);CHKERRQ(ierr);
191       ierr = DMGetDS(cdm, &cds);CHKERRQ(ierr);
192       ierr = PetscDSSetConstants(cds, 4, ctx->transformData);CHKERRQ(ierr);
193       ierr = DMPlexRemapGeometry(*dm, 0.0, f0_flare);CHKERRQ(ierr);
194       break;
195     case TRANSFORM_ANNULUS:
196       ierr = DMGetCoordinateDM(*dm, &cdm);CHKERRQ(ierr);
197       ierr = DMGetDS(cdm, &cds);CHKERRQ(ierr);
198       ierr = PetscDSSetConstants(cds, 2, ctx->transformData);CHKERRQ(ierr);
199       ierr = DMPlexRemapGeometry(*dm, 0.0, f0_annulus);CHKERRQ(ierr);
200       break;
201     case TRANSFORM_SHELL:
202       ierr = DMGetCoordinateDM(*dm, &cdm);CHKERRQ(ierr);
203       ierr = DMGetDS(cdm, &cds);CHKERRQ(ierr);
204       ierr = PetscDSSetConstants(cds, 2, ctx->transformData);CHKERRQ(ierr);
205       ierr = DMPlexRemapGeometry(*dm, 0.0, f0_shell);CHKERRQ(ierr);
206       break;
207     default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Unknown mesh transform %D", ctx->meshTransform);
208   }
209   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
210   PetscFunctionReturn(0);
211 }
212 
213 static void volume(PetscInt dim, PetscInt Nf, PetscInt NfAux,
214                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
215                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
216                    PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar vol[])
217 {
218   vol[0] = 1.;
219 }
220 
221 static PetscErrorCode CreateDiscretization(DM dm, AppCtx *ctx)
222 {
223   PetscDS        ds;
224   PetscFE        fe;
225   DMPolytopeType ct;
226   PetscInt       dim, cStart;
227   PetscBool      simplex;
228   PetscErrorCode ierr;
229 
230   PetscFunctionBeginUser;
231   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
232   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr);
233   ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr);
234   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
235   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, PETSC_DETERMINE, &fe);CHKERRQ(ierr);
236   ierr = PetscFESetName(fe, "scalar");CHKERRQ(ierr);
237   ierr = DMAddField(dm, NULL, (PetscObject) fe);
238   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
239   ierr = DMCreateDS(dm);CHKERRQ(ierr);
240   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
241   ierr = PetscDSSetObjective(ds, 0, volume);CHKERRQ(ierr);
242   PetscFunctionReturn(0);
243 }
244 
245 static PetscErrorCode CheckVolume(DM dm, AppCtx *ctx)
246 {
247   Vec            u;
248   PetscScalar    result;
249   PetscReal      vol, tol = ctx->tol;
250   PetscErrorCode ierr;
251 
252   PetscFunctionBeginUser;
253   ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr);
254   ierr = DMPlexComputeIntegralFEM(dm, u, &result, ctx);CHKERRQ(ierr);
255   vol  = PetscRealPart(result);
256   ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr);
257   ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Volume: %g\n", (double) vol);CHKERRQ(ierr);
258   if (ctx->volume > 0.0 && PetscAbsReal(ctx->volume - vol) > tol) {
259     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);
260   }
261   PetscFunctionReturn(0);
262 }
263 
264 int main(int argc, char **argv)
265 {
266   DM             dm;
267   AppCtx         user;
268   PetscErrorCode ierr;
269 
270   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
271   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
272   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
273   ierr = CreateDiscretization(dm, &user);CHKERRQ(ierr);
274   ierr = CheckVolume(dm, &user);CHKERRQ(ierr);
275   ierr = DMDestroy(&dm);CHKERRQ(ierr);
276   ierr = PetscFree(user.transformDataReal);CHKERRQ(ierr);
277   ierr = PetscFree(user.transformData);CHKERRQ(ierr);
278   ierr = PetscFinalize();
279   return ierr;
280 }
281 
282 /*TEST
283 
284   testset:
285     args: -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1.,-1. -dm_plex_box_upper 1.,1. -volume 4. -dm_coord_space 0
286 
287     test:
288       suffix: square_0
289       args: -dm_coord_petscspace_degree 1
290 
291     test:
292       suffix: square_1
293       args: -dm_coord_petscspace_degree 2
294 
295     test:
296       suffix: square_2
297       args: -dm_refine 1 -dm_coord_petscspace_degree 1
298 
299     test:
300       suffix: square_3
301       args: -dm_refine 1 -dm_coord_petscspace_degree 2
302 
303   testset:
304     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -volume 8. -dm_coord_space 0
305 
306     test:
307       suffix: cube_0
308       args: -dm_coord_petscspace_degree 1
309 
310     test:
311       suffix: cube_1
312       args: -dm_coord_petscspace_degree 2
313 
314     test:
315       suffix: cube_2
316       args: -dm_refine 1 -dm_coord_petscspace_degree 1
317 
318     test:
319       suffix: cube_3
320       args: -dm_refine 1 -dm_coord_petscspace_degree 2
321 
322   testset:
323     args: -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1.,-1. -dm_plex_box_upper 1.,1. -volume 4. -dm_coord_space 0
324 
325     test:
326       suffix: shear_0
327       args: -dm_coord_petscspace_degree 1 -mesh_transform shear -transform_data 3.0
328 
329     test:
330       suffix: shear_1
331       args: -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0
332 
333     test:
334       suffix: shear_2
335       args: -dm_refine 1 -dm_coord_petscspace_degree 1 -mesh_transform shear -transform_data 3.0
336 
337     test:
338       suffix: shear_3
339       args: -dm_refine 1 -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0
340 
341   testset:
342     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -volume 8. -dm_coord_space 0
343 
344     test:
345       suffix: shear_4
346       args: -dm_coord_petscspace_degree 1 -mesh_transform shear -transform_data 3.0
347 
348     test:
349       suffix: shear_5
350       args: -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0
351 
352     test:
353       suffix: shear_6
354       args: -dm_refine 1 -dm_coord_petscspace_degree 1 -mesh_transform shear -transform_data 3.0,4.0
355 
356     test:
357       suffix: shear_7
358       args: -dm_refine 1 -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0,4.0
359 
360   testset:
361     args: -dm_coord_space 0 -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower 1.,-1. -dm_plex_box_upper 3.,1. \
362           -dm_coord_petscspace_degree 1 -mesh_transform flare -volume 8.
363 
364     test:
365       suffix: flare_0
366       args:
367 
368     test:
369       suffix: flare_1
370       args: -dm_refine 2
371 
372   testset:
373     # Area: 3/4 \pi = 2.3562
374     args: -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -mesh_transform annulus -volume 2.35619449019235 -dm_coord_space 0
375 
376     test:
377       # Area: (a+b)/2 h = 3/\sqrt{2} (sqrt{2} - 1/\sqrt{2}) = 3/2
378       suffix: annulus_0
379       requires: double
380       args: -dm_coord_petscspace_degree 1 -volume 1.5
381 
382     test:
383       suffix: annulus_1
384       requires: double
385       args: -dm_refine 3 -dm_coord_petscspace_degree 1 -tol .016
386 
387     test:
388       suffix: annulus_2
389       requires: double
390       args: -dm_refine 3 -dm_coord_petscspace_degree 2 -tol .0038
391 
392     test:
393       suffix: annulus_3
394       requires: double
395       args: -dm_refine 3 -dm_coord_petscspace_degree 3 -tol 2.2e-6
396 
397     test:
398       suffix: annulus_4
399       requires: double
400       args: -dm_refine 2 -dm_coord_petscspace_degree 2 -petscfe_default_quadrature_order 2 -tol .00012
401 
402     test:
403       suffix: annulus_5
404       requires: double
405       args: -dm_refine 2 -dm_coord_petscspace_degree 3 -petscfe_default_quadrature_order 3 -tol 1.2e-7
406 
407   testset:
408     # Volume: 4/3 \pi (8 - 1)/2 = 14/3 \pi = 14.66076571675238
409     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -mesh_transform shell -volume 14.66076571675238 -dm_coord_space 0
410 
411     test:
412       suffix: shell_0
413       requires: double
414       args: -dm_refine 1 -dm_coord_petscspace_degree 1 -petscfe_default_quadrature_order 1 -volume 5.633164922 -tol 1.0e-7
415 
416     test:
417       suffix: shell_1
418       requires: double
419       args: -dm_refine 2 -dm_coord_petscspace_degree 1 -petscfe_default_quadrature_order 1 -tol 3.1
420 
421     test:
422       suffix: shell_2
423       requires: double
424       args: -dm_refine 2 -dm_coord_petscspace_degree 2 -petscfe_default_quadrature_order 2 -tol .1
425 
426     test:
427       suffix: shell_3
428       requires: double
429       args: -dm_refine 2 -dm_coord_petscspace_degree 3 -petscfe_default_quadrature_order 3 -tol .02
430 
431     test:
432       suffix: shell_4
433       requires: double
434       args: -dm_refine 2 -dm_coord_petscspace_degree 4 -petscfe_default_quadrature_order 4 -tol .006
435 
436   test:
437     # Volume: 1.0
438     suffix: gmsh_q2
439     requires: double
440     args: -coord_space 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/quads-q2.msh -dm_plex_gmsh_project -volume 1.0 -tol 1e-6
441 
442   test:
443     # Volume: 1.0
444     suffix: gmsh_q3
445     requires: double
446     args: -coord_space 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/quads-q3.msh -dm_plex_gmsh_project -volume 1.0 -tol 1e-6
447 
448 TEST*/
449