xref: /petsc/src/dm/impls/plex/tests/ex60.c (revision 00d473e3650901c295a590c8be836ff35aa7b0a9)
1*00d473e3SJoe Wallwork static char help[] = "Test metric utils in the uniform, isotropic case.\n\n";
2*00d473e3SJoe Wallwork 
3*00d473e3SJoe Wallwork #include <petscdmplex.h>
4*00d473e3SJoe Wallwork 
5*00d473e3SJoe Wallwork static PetscErrorCode bowl(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
6*00d473e3SJoe Wallwork {
7*00d473e3SJoe Wallwork   PetscInt d;
8*00d473e3SJoe Wallwork 
9*00d473e3SJoe Wallwork   *u = 0.0;
10*00d473e3SJoe Wallwork   for (d = 0; d < dim; d++) *u += 0.5*(x[d] - 0.5)*(x[d] - 0.5);
11*00d473e3SJoe Wallwork 
12*00d473e3SJoe Wallwork   return 0;
13*00d473e3SJoe Wallwork }
14*00d473e3SJoe Wallwork 
15*00d473e3SJoe Wallwork int main(int argc, char **argv) {
16*00d473e3SJoe Wallwork   DM              dm, dmDist, dmAdapt;
17*00d473e3SJoe Wallwork   DMLabel         bdLabel = NULL;
18*00d473e3SJoe Wallwork   DMPlexMetricCtx user;
19*00d473e3SJoe Wallwork   MPI_Comm        comm;
20*00d473e3SJoe Wallwork   PetscBool       uniform = PETSC_FALSE;
21*00d473e3SJoe Wallwork   PetscErrorCode  ierr;
22*00d473e3SJoe Wallwork   PetscInt       *faces, dim = 3, numEdges = 4, d;
23*00d473e3SJoe Wallwork   PetscReal       scaling = 1.0;
24*00d473e3SJoe Wallwork   Vec             metric;
25*00d473e3SJoe Wallwork 
26*00d473e3SJoe Wallwork   /* Default parameter values */
27*00d473e3SJoe Wallwork   user.p = 1.0;
28*00d473e3SJoe Wallwork   user.targetComplexity = 100.0;
29*00d473e3SJoe Wallwork   user.isotropic = PETSC_FALSE;
30*00d473e3SJoe Wallwork 
31*00d473e3SJoe Wallwork   /* Set up */
32*00d473e3SJoe Wallwork   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
33*00d473e3SJoe Wallwork   comm = PETSC_COMM_WORLD;
34*00d473e3SJoe Wallwork   ierr = PetscOptionsBegin(comm, "", "Mesh adaptation options", "DMPLEX");CHKERRQ(ierr);
35*00d473e3SJoe Wallwork   ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex60.c", dim, &dim, NULL, 2, 3);CHKERRQ(ierr);
36*00d473e3SJoe Wallwork   ierr = PetscOptionsInt("-num_edges", "Number of edges on each boundary of the initial mesh", "ex60.c", numEdges, &numEdges, NULL);CHKERRQ(ierr);
37*00d473e3SJoe Wallwork   ierr = PetscOptionsReal("-target", "Target metric complexity", "ex60.c", user.targetComplexity, &user.targetComplexity, NULL);CHKERRQ(ierr);
38*00d473e3SJoe Wallwork   ierr = PetscOptionsReal("-norm_order", "Metric normalization order", "ex60.c", user.p, &user.p, NULL);CHKERRQ(ierr);
39*00d473e3SJoe Wallwork   ierr = PetscOptionsBool("-uniform", "Should the metric be assumed uniform?", "ex60.c", uniform, &uniform, NULL);CHKERRQ(ierr);
40*00d473e3SJoe Wallwork   ierr = PetscOptionsBool("-isotropic", "Should the metric be assumed isotropic, or computed as a recovered Hessian?", "ex60.c", user.isotropic, &user.isotropic, NULL);CHKERRQ(ierr);
41*00d473e3SJoe Wallwork   ierr = PetscOptionsEnd();
42*00d473e3SJoe Wallwork 
43*00d473e3SJoe Wallwork   /* Create box mesh */
44*00d473e3SJoe Wallwork   ierr = PetscMalloc1(dim, &faces);CHKERRQ(ierr);
45*00d473e3SJoe Wallwork   for (d = 0; d < dim; ++d) faces[d] = numEdges;
46*00d473e3SJoe Wallwork   ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_TRUE, faces, NULL, NULL, NULL, PETSC_TRUE, &dm);CHKERRQ(ierr);
47*00d473e3SJoe Wallwork   ierr = PetscFree(faces);CHKERRQ(ierr);
48*00d473e3SJoe Wallwork 
49*00d473e3SJoe Wallwork   /* Distribute mesh over processes */
50*00d473e3SJoe Wallwork   ierr = DMPlexDistribute(dm, 0, NULL, &dmDist);CHKERRQ(ierr);
51*00d473e3SJoe Wallwork   if (dmDist) {
52*00d473e3SJoe Wallwork     ierr = DMDestroy(&dm);CHKERRQ(ierr);
53*00d473e3SJoe Wallwork     dm = dmDist;
54*00d473e3SJoe Wallwork   }
55*00d473e3SJoe Wallwork   ierr = DMSetApplicationContext(dm, (void**)&user);CHKERRQ(ierr);
56*00d473e3SJoe Wallwork   ierr = PetscObjectSetName((PetscObject) dm, "DM_init");CHKERRQ(ierr);
57*00d473e3SJoe Wallwork   ierr = DMViewFromOptions(dm, NULL, "-initial_mesh_view");CHKERRQ(ierr);
58*00d473e3SJoe Wallwork 
59*00d473e3SJoe Wallwork   /* Construct metric */
60*00d473e3SJoe Wallwork   if (uniform) {
61*00d473e3SJoe Wallwork     if (user.isotropic) { ierr = DMPlexMetricCreateUniform(dm, 0, scaling, &metric);CHKERRQ(ierr); }
62*00d473e3SJoe Wallwork     else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics not supported.");
63*00d473e3SJoe Wallwork   }
64*00d473e3SJoe Wallwork   else {
65*00d473e3SJoe Wallwork     DM      dmIndi;
66*00d473e3SJoe Wallwork     PetscFE fe;
67*00d473e3SJoe Wallwork     Vec     indicator;
68*00d473e3SJoe Wallwork 
69*00d473e3SJoe Wallwork     /* Construct "error indicator" */
70*00d473e3SJoe Wallwork     ierr = DMClone(dm, &dmIndi);CHKERRQ(ierr);
71*00d473e3SJoe Wallwork     ierr = PetscFECreateLagrange(comm, dim, 1, PETSC_TRUE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr);
72*00d473e3SJoe Wallwork     ierr = DMAddField(dmIndi, NULL, (PetscObject)fe);CHKERRQ(ierr);
73*00d473e3SJoe Wallwork     ierr = DMCreateDS(dmIndi);CHKERRQ(ierr);
74*00d473e3SJoe Wallwork     ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
75*00d473e3SJoe Wallwork     ierr = DMCreateLocalVector(dmIndi, &indicator);CHKERRQ(ierr);
76*00d473e3SJoe Wallwork     if (user.isotropic) {
77*00d473e3SJoe Wallwork 
78*00d473e3SJoe Wallwork       /* Isotropic case: just specify unity */
79*00d473e3SJoe Wallwork       ierr = VecSet(indicator, scaling);CHKERRQ(ierr);
80*00d473e3SJoe Wallwork       ierr = DMPlexMetricCreateIsotropic(dm, 0, indicator, &metric);CHKERRQ(ierr);
81*00d473e3SJoe Wallwork 
82*00d473e3SJoe Wallwork     } else {
83*00d473e3SJoe Wallwork 
84*00d473e3SJoe Wallwork       /* 'Anisotropic' case: approximate the identity by recovering the Hessian of a parabola */
85*00d473e3SJoe Wallwork       DM               dmGrad;
86*00d473e3SJoe Wallwork       PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar*, void*) = {bowl};
87*00d473e3SJoe Wallwork       PetscInt         pStart, pEnd, p;
88*00d473e3SJoe Wallwork       PetscSection     sec, secCoord;
89*00d473e3SJoe Wallwork       Vec              gradient;
90*00d473e3SJoe Wallwork 
91*00d473e3SJoe Wallwork       /* Project the parabola into P1 space */
92*00d473e3SJoe Wallwork       ierr = DMProjectFunctionLocal(dmIndi, 0.0, funcs, NULL, INSERT_ALL_VALUES, indicator);CHKERRQ(ierr);
93*00d473e3SJoe Wallwork 
94*00d473e3SJoe Wallwork       /* Approximate the gradient */
95*00d473e3SJoe Wallwork       ierr = DMClone(dmIndi, &dmGrad);CHKERRQ(ierr);
96*00d473e3SJoe Wallwork       ierr = DMGetCoordinateSection(dmGrad, &secCoord);CHKERRQ(ierr);
97*00d473e3SJoe Wallwork       ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dmGrad), &sec);CHKERRQ(ierr);
98*00d473e3SJoe Wallwork       ierr = PetscSectionSetNumFields(sec, 1);CHKERRQ(ierr);
99*00d473e3SJoe Wallwork       ierr = PetscSectionSetFieldComponents(sec, 0, dim);CHKERRQ(ierr);
100*00d473e3SJoe Wallwork       ierr = PetscSectionGetChart(secCoord, &pStart, &pEnd);CHKERRQ(ierr);
101*00d473e3SJoe Wallwork       ierr = PetscSectionSetChart(sec, pStart, pEnd);CHKERRQ(ierr);
102*00d473e3SJoe Wallwork       for (p = pStart; p < pEnd; ++p) {
103*00d473e3SJoe Wallwork         ierr = PetscSectionSetDof(sec, p, dim);CHKERRQ(ierr);
104*00d473e3SJoe Wallwork         ierr = PetscSectionSetFieldDof(sec, p, 0, dim);CHKERRQ(ierr);
105*00d473e3SJoe Wallwork       }
106*00d473e3SJoe Wallwork       ierr = PetscSectionSetUp(sec);CHKERRQ(ierr);
107*00d473e3SJoe Wallwork       ierr = DMSetLocalSection(dmGrad, sec);CHKERRQ(ierr);
108*00d473e3SJoe Wallwork       ierr = PetscSectionDestroy(&sec);CHKERRQ(ierr);
109*00d473e3SJoe Wallwork       ierr = PetscFECreateLagrange(comm, dim, dim, PETSC_TRUE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr);
110*00d473e3SJoe Wallwork       ierr = DMClearFields(dmGrad);CHKERRQ(ierr);
111*00d473e3SJoe Wallwork       ierr = DMAddField(dmGrad, NULL, (PetscObject)fe);CHKERRQ(ierr);
112*00d473e3SJoe Wallwork       ierr = DMCreateDS(dmGrad);CHKERRQ(ierr);
113*00d473e3SJoe Wallwork       ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
114*00d473e3SJoe Wallwork       ierr = DMCreateLocalVector(dmGrad, &gradient);CHKERRQ(ierr);
115*00d473e3SJoe Wallwork       ierr = DMPlexComputeGradientClementInterpolant(dmIndi, indicator, gradient);CHKERRQ(ierr);
116*00d473e3SJoe Wallwork       ierr = VecViewFromOptions(gradient, NULL, "-adapt_gradient_view");CHKERRQ(ierr);
117*00d473e3SJoe Wallwork 
118*00d473e3SJoe Wallwork       /* Approximate the Hessian */
119*00d473e3SJoe Wallwork       ierr = DMGetCoordinateSection(dm, &secCoord);CHKERRQ(ierr);
120*00d473e3SJoe Wallwork       ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sec);CHKERRQ(ierr);
121*00d473e3SJoe Wallwork       ierr = PetscSectionSetNumFields(sec, 1);CHKERRQ(ierr);
122*00d473e3SJoe Wallwork       ierr = PetscSectionSetFieldComponents(sec, 0, dim*dim);CHKERRQ(ierr);
123*00d473e3SJoe Wallwork       ierr = PetscSectionGetChart(secCoord, &pStart, &pEnd);CHKERRQ(ierr);
124*00d473e3SJoe Wallwork       ierr = PetscSectionSetChart(sec, pStart, pEnd);CHKERRQ(ierr);
125*00d473e3SJoe Wallwork       for (p = pStart; p < pEnd; ++p) {
126*00d473e3SJoe Wallwork         ierr = PetscSectionSetDof(sec, p, dim*dim);CHKERRQ(ierr);
127*00d473e3SJoe Wallwork         ierr = PetscSectionSetFieldDof(sec, p, 0, dim*dim);CHKERRQ(ierr);
128*00d473e3SJoe Wallwork       }
129*00d473e3SJoe Wallwork       ierr = PetscSectionSetUp(sec);CHKERRQ(ierr);
130*00d473e3SJoe Wallwork       ierr = DMSetLocalSection(dm, sec);CHKERRQ(ierr);
131*00d473e3SJoe Wallwork       ierr = PetscSectionDestroy(&sec);CHKERRQ(ierr);
132*00d473e3SJoe Wallwork       ierr = PetscFECreateLagrange(comm, dim, dim*dim, PETSC_TRUE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr);
133*00d473e3SJoe Wallwork       ierr = DMClearFields(dm);CHKERRQ(ierr);
134*00d473e3SJoe Wallwork       ierr = DMAddField(dm, NULL, (PetscObject)fe);CHKERRQ(ierr);
135*00d473e3SJoe Wallwork       ierr = DMCreateDS(dm);CHKERRQ(ierr);
136*00d473e3SJoe Wallwork       ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
137*00d473e3SJoe Wallwork       ierr = DMCreateLocalVector(dm, &metric);CHKERRQ(ierr);
138*00d473e3SJoe Wallwork       ierr = DMPlexComputeGradientClementInterpolant(dmGrad, gradient, metric);CHKERRQ(ierr);
139*00d473e3SJoe Wallwork       ierr = VecViewFromOptions(metric, NULL, "-adapt_hessian_view");CHKERRQ(ierr);
140*00d473e3SJoe Wallwork       ierr = VecDestroy(&gradient);CHKERRQ(ierr);
141*00d473e3SJoe Wallwork       ierr = DMDestroy(&dmGrad);CHKERRQ(ierr);
142*00d473e3SJoe Wallwork     }
143*00d473e3SJoe Wallwork     ierr = VecDestroy(&indicator);CHKERRQ(ierr);
144*00d473e3SJoe Wallwork     ierr = DMDestroy(&dmIndi);CHKERRQ(ierr);
145*00d473e3SJoe Wallwork   }
146*00d473e3SJoe Wallwork 
147*00d473e3SJoe Wallwork   /* Test metric routines */
148*00d473e3SJoe Wallwork   {
149*00d473e3SJoe Wallwork     PetscReal   errornorm, norm, tol = 1.0e-10, weights[2] = {0.8, 0.2};
150*00d473e3SJoe Wallwork     Vec         metric1, metric2, metricComb;
151*00d473e3SJoe Wallwork     Vec         metrics[2];
152*00d473e3SJoe Wallwork 
153*00d473e3SJoe Wallwork     ierr = VecDuplicate(metric, &metric1);CHKERRQ(ierr);
154*00d473e3SJoe Wallwork     ierr = VecSet(metric1, 0);CHKERRQ(ierr);
155*00d473e3SJoe Wallwork     ierr = VecAXPY(metric1, 0.625, metric);CHKERRQ(ierr);
156*00d473e3SJoe Wallwork     ierr = VecDuplicate(metric, &metric2);CHKERRQ(ierr);
157*00d473e3SJoe Wallwork     ierr = VecSet(metric2, 0);CHKERRQ(ierr);
158*00d473e3SJoe Wallwork     ierr = VecAXPY(metric2, 2.5, metric);CHKERRQ(ierr);
159*00d473e3SJoe Wallwork     metrics[0] = metric1;
160*00d473e3SJoe Wallwork     metrics[1] = metric2;
161*00d473e3SJoe Wallwork 
162*00d473e3SJoe Wallwork     /* Test metric average */
163*00d473e3SJoe Wallwork     ierr = DMPlexMetricAverage(dm, 2, weights, metrics, &metricComb);CHKERRQ(ierr);
164*00d473e3SJoe Wallwork     ierr = VecAXPY(metricComb, -1, metric);CHKERRQ(ierr);
165*00d473e3SJoe Wallwork     ierr = VecNorm(metric, NORM_2, &norm);CHKERRQ(ierr);
166*00d473e3SJoe Wallwork     ierr = VecNorm(metricComb, NORM_2, &errornorm);CHKERRQ(ierr);
167*00d473e3SJoe Wallwork     errornorm /= norm;
168*00d473e3SJoe Wallwork     if (errornorm > tol) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric average test failed (L2 error %f)", errornorm);
169*00d473e3SJoe Wallwork     ierr = VecDestroy(&metricComb);CHKERRQ(ierr);
170*00d473e3SJoe Wallwork 
171*00d473e3SJoe Wallwork     /* Test metric intersection */
172*00d473e3SJoe Wallwork     if (user.isotropic) {
173*00d473e3SJoe Wallwork       ierr = DMPlexMetricIntersection(dm, 2, metrics, &metricComb);CHKERRQ(ierr);
174*00d473e3SJoe Wallwork       ierr = VecAXPY(metricComb, -1, metric1);CHKERRQ(ierr);
175*00d473e3SJoe Wallwork       ierr = VecNorm(metricComb, NORM_2, &errornorm);CHKERRQ(ierr);
176*00d473e3SJoe Wallwork       errornorm /= norm;
177*00d473e3SJoe Wallwork       if (errornorm > tol) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric intersection test failed (L2 error %f)", errornorm);
178*00d473e3SJoe Wallwork     }
179*00d473e3SJoe Wallwork     ierr = VecDestroy(&metric2);CHKERRQ(ierr);
180*00d473e3SJoe Wallwork     ierr = VecDestroy(&metricComb);CHKERRQ(ierr);
181*00d473e3SJoe Wallwork     ierr = VecCopy(metric, metric1);CHKERRQ(ierr);
182*00d473e3SJoe Wallwork 
183*00d473e3SJoe Wallwork     /* Test metric SPD enforcement */
184*00d473e3SJoe Wallwork     ierr = DMPlexMetricEnforceSPD(dm, PETSC_TRUE, metric);CHKERRQ(ierr);
185*00d473e3SJoe Wallwork     if (user.isotropic) {
186*00d473e3SJoe Wallwork       ierr = VecAXPY(metric1, -1, metric);CHKERRQ(ierr);
187*00d473e3SJoe Wallwork       ierr = VecNorm(metric1, NORM_2, &errornorm);CHKERRQ(ierr);
188*00d473e3SJoe Wallwork       errornorm /= norm;
189*00d473e3SJoe Wallwork       if (errornorm > tol) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric SPD enforcement test failed (L2 error %f)", errornorm);
190*00d473e3SJoe Wallwork     }
191*00d473e3SJoe Wallwork     ierr = VecDestroy(&metric1);CHKERRQ(ierr);
192*00d473e3SJoe Wallwork 
193*00d473e3SJoe Wallwork     /* Test metric normalization */
194*00d473e3SJoe Wallwork     ierr = DMPlexMetricNormalize(dm, metric, PETSC_TRUE, &metric1);CHKERRQ(ierr);
195*00d473e3SJoe Wallwork     if (user.isotropic) {
196*00d473e3SJoe Wallwork       scaling = PetscPowReal(user.targetComplexity, 2.0/dim);
197*00d473e3SJoe Wallwork       ierr = DMPlexMetricCreateUniform(dm, 0, scaling, &metric2);CHKERRQ(ierr);
198*00d473e3SJoe Wallwork       ierr = VecAXPY(metric2, -1, metric1);CHKERRQ(ierr);
199*00d473e3SJoe Wallwork       ierr = VecNorm(metric2, NORM_2, &errornorm);CHKERRQ(ierr);
200*00d473e3SJoe Wallwork       errornorm /= norm;
201*00d473e3SJoe Wallwork       if (errornorm > tol) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric normalization test failed (L2 error %f)", errornorm);
202*00d473e3SJoe Wallwork     }
203*00d473e3SJoe Wallwork     ierr = VecCopy(metric1, metric);CHKERRQ(ierr);
204*00d473e3SJoe Wallwork     ierr = VecDestroy(&metric2);CHKERRQ(ierr);
205*00d473e3SJoe Wallwork     ierr = VecDestroy(&metric1);CHKERRQ(ierr);
206*00d473e3SJoe Wallwork   }
207*00d473e3SJoe Wallwork 
208*00d473e3SJoe Wallwork   /* Adapt the mesh */
209*00d473e3SJoe Wallwork   ierr = DMAdaptMetric(dm, metric, bdLabel, &dmAdapt);CHKERRQ(ierr);
210*00d473e3SJoe Wallwork   ierr = PetscObjectSetName((PetscObject) dmAdapt, "DM_adapted");CHKERRQ(ierr);
211*00d473e3SJoe Wallwork   ierr = VecDestroy(&metric);CHKERRQ(ierr);
212*00d473e3SJoe Wallwork   ierr = DMViewFromOptions(dmAdapt, NULL, "-adapted_mesh_view");CHKERRQ(ierr);
213*00d473e3SJoe Wallwork 
214*00d473e3SJoe Wallwork   /* Compare DMs */
215*00d473e3SJoe Wallwork   ierr = DMView(dm, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
216*00d473e3SJoe Wallwork   ierr = DMView(dmAdapt, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
217*00d473e3SJoe Wallwork 
218*00d473e3SJoe Wallwork   /* Clean up */
219*00d473e3SJoe Wallwork   ierr = DMDestroy(&dmAdapt);CHKERRQ(ierr);
220*00d473e3SJoe Wallwork   ierr = DMDestroy(&dm);CHKERRQ(ierr);
221*00d473e3SJoe Wallwork   ierr = PetscFinalize();
222*00d473e3SJoe Wallwork   return 0;
223*00d473e3SJoe Wallwork }
224*00d473e3SJoe Wallwork 
225*00d473e3SJoe Wallwork /*TEST
226*00d473e3SJoe Wallwork 
227*00d473e3SJoe Wallwork   build:
228*00d473e3SJoe Wallwork     requires: pragmatic
229*00d473e3SJoe Wallwork 
230*00d473e3SJoe Wallwork   test:
231*00d473e3SJoe Wallwork     suffix: uniform_2d
232*00d473e3SJoe Wallwork     args: -dim 2 -uniform -isotropic
233*00d473e3SJoe Wallwork   test:
234*00d473e3SJoe Wallwork     suffix: uniform_3d
235*00d473e3SJoe Wallwork     args: -dim 3 -uniform -isotropic
236*00d473e3SJoe Wallwork   test:
237*00d473e3SJoe Wallwork     suffix: iso_2d
238*00d473e3SJoe Wallwork     args: -dim 2 -isotropic
239*00d473e3SJoe Wallwork   test:
240*00d473e3SJoe Wallwork     suffix: iso_3d
241*00d473e3SJoe Wallwork     args: -dim 3 -isotropic
242*00d473e3SJoe Wallwork   test:
243*00d473e3SJoe Wallwork     suffix: hessian_2d
244*00d473e3SJoe Wallwork     args: -dim 2
245*00d473e3SJoe Wallwork   test:
246*00d473e3SJoe Wallwork     suffix: hessian_3d
247*00d473e3SJoe Wallwork     args: -dim 3
248*00d473e3SJoe Wallwork 
249*00d473e3SJoe Wallwork TEST*/
250