xref: /petsc/src/dm/impls/plex/tests/ex60.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
100d473e3SJoe Wallwork static char help[] = "Test metric utils in the uniform, isotropic case.\n\n";
200d473e3SJoe Wallwork 
300d473e3SJoe Wallwork #include <petscdmplex.h>
400d473e3SJoe Wallwork 
500d473e3SJoe Wallwork static PetscErrorCode bowl(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
600d473e3SJoe Wallwork {
700d473e3SJoe Wallwork   PetscInt d;
800d473e3SJoe Wallwork 
900d473e3SJoe Wallwork   *u = 0.0;
1000d473e3SJoe Wallwork   for (d = 0; d < dim; d++) *u += 0.5*(x[d] - 0.5)*(x[d] - 0.5);
1100d473e3SJoe Wallwork 
1200d473e3SJoe Wallwork   return 0;
1300d473e3SJoe Wallwork }
1400d473e3SJoe Wallwork 
15d8b80fabSJoe Wallwork static PetscErrorCode CreateIndicator(DM dm, Vec *indicator, DM *dmIndi)
16d8b80fabSJoe Wallwork {
17d8b80fabSJoe Wallwork   MPI_Comm       comm;
18d8b80fabSJoe Wallwork   PetscErrorCode ierr;
19d8b80fabSJoe Wallwork   PetscFE        fe;
20d8b80fabSJoe Wallwork   PetscInt       dim;
21d8b80fabSJoe Wallwork 
22d8b80fabSJoe Wallwork   PetscFunctionBeginUser;
23d8b80fabSJoe Wallwork   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
24d8b80fabSJoe Wallwork   ierr = DMClone(dm, dmIndi);CHKERRQ(ierr);
25d8b80fabSJoe Wallwork   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
26d8b80fabSJoe Wallwork   ierr = PetscFECreateLagrange(comm, dim, 1, PETSC_TRUE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr);
27d8b80fabSJoe Wallwork   ierr = DMSetField(*dmIndi, 0, NULL, (PetscObject)fe);CHKERRQ(ierr);
28d8b80fabSJoe Wallwork   ierr = DMCreateDS(*dmIndi);CHKERRQ(ierr);
29d8b80fabSJoe Wallwork   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
30d8b80fabSJoe Wallwork   ierr = DMCreateLocalVector(*dmIndi, indicator);CHKERRQ(ierr);
31d8b80fabSJoe Wallwork   PetscFunctionReturn(0);
32d8b80fabSJoe Wallwork }
33d8b80fabSJoe Wallwork 
3400d473e3SJoe Wallwork int main(int argc, char **argv) {
3500d473e3SJoe Wallwork   DM              dm, dmDist, dmAdapt;
3621b3e102SJoe Wallwork   DMLabel         bdLabel = NULL, rgLabel = NULL;
3700d473e3SJoe Wallwork   MPI_Comm        comm;
3841473654SJoe Wallwork   PetscBool       uniform = PETSC_FALSE, isotropic = PETSC_FALSE, noTagging = PETSC_FALSE;
3900d473e3SJoe Wallwork   PetscErrorCode  ierr;
40f49e87b0SJoe Wallwork   PetscInt       *faces, dim = 3, d;
4100d473e3SJoe Wallwork   PetscReal       scaling = 1.0;
4200d473e3SJoe Wallwork   Vec             metric;
4300d473e3SJoe Wallwork 
4400d473e3SJoe Wallwork   /* Set up */
4500d473e3SJoe Wallwork   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
4600d473e3SJoe Wallwork   comm = PETSC_COMM_WORLD;
4700d473e3SJoe Wallwork   ierr = PetscOptionsBegin(comm, "", "Mesh adaptation options", "DMPLEX");CHKERRQ(ierr);
4800d473e3SJoe Wallwork   ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex60.c", dim, &dim, NULL, 2, 3);CHKERRQ(ierr);
4941473654SJoe Wallwork   ierr = PetscOptionsBool("-noTagging", "Should tag preservation testing be turned off?", "ex60.c", noTagging, &noTagging, NULL);CHKERRQ(ierr);
5000d473e3SJoe Wallwork   ierr = PetscOptionsEnd();
5100d473e3SJoe Wallwork 
5200d473e3SJoe Wallwork   /* Create box mesh */
5300d473e3SJoe Wallwork   ierr = PetscMalloc1(dim, &faces);CHKERRQ(ierr);
54f49e87b0SJoe Wallwork   for (d = 0; d < dim; ++d) faces[d] = 4;
5500d473e3SJoe Wallwork   ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_TRUE, faces, NULL, NULL, NULL, PETSC_TRUE, &dm);CHKERRQ(ierr);
5600d473e3SJoe Wallwork   ierr = PetscFree(faces);CHKERRQ(ierr);
57e5697243SJoe Wallwork   ierr = DMSetFromOptions(dm);CHKERRQ(ierr);
5800d473e3SJoe Wallwork 
5900d473e3SJoe Wallwork   /* Distribute mesh over processes */
6000d473e3SJoe Wallwork   ierr = DMPlexDistribute(dm, 0, NULL, &dmDist);CHKERRQ(ierr);
6100d473e3SJoe Wallwork   if (dmDist) {
6200d473e3SJoe Wallwork     ierr = DMDestroy(&dm);CHKERRQ(ierr);
6300d473e3SJoe Wallwork     dm = dmDist;
6400d473e3SJoe Wallwork   }
6500d473e3SJoe Wallwork   ierr = PetscObjectSetName((PetscObject) dm, "DM_init");CHKERRQ(ierr);
6600d473e3SJoe Wallwork   ierr = DMViewFromOptions(dm, NULL, "-initial_mesh_view");CHKERRQ(ierr);
6700d473e3SJoe Wallwork 
6841473654SJoe Wallwork   /* Set tags to be preserved */
6941473654SJoe Wallwork   if (!noTagging) {
7041473654SJoe Wallwork     DM                 cdm;
7141473654SJoe Wallwork     PetscInt           cStart, cEnd, c, fStart, fEnd, f, vStart, vEnd;
7241473654SJoe Wallwork     const PetscScalar *coords;
7341473654SJoe Wallwork     Vec                coordinates;
7441473654SJoe Wallwork 
7541473654SJoe Wallwork     /* Cell tags */
7641473654SJoe Wallwork     ierr = DMCreateLabel(dm, "Cell Sets");CHKERRQ(ierr);
7741473654SJoe Wallwork     ierr = DMGetLabel(dm, "Cell Sets", &rgLabel);CHKERRQ(ierr);
7841473654SJoe Wallwork     ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
7941473654SJoe Wallwork     for (c = cStart; c < cEnd; ++c) {
8041473654SJoe Wallwork       PetscReal centroid[3], volume, x;
8141473654SJoe Wallwork 
8241473654SJoe Wallwork       ierr = DMPlexComputeCellGeometryFVM(dm, c, &volume, centroid, NULL);CHKERRQ(ierr);
8341473654SJoe Wallwork       x = centroid[0];
8441473654SJoe Wallwork       if (x < 0.5) { ierr = DMLabelSetValue(rgLabel, c, 3);CHKERRQ(ierr); }
8541473654SJoe Wallwork       else         { ierr = DMLabelSetValue(rgLabel, c, 4);CHKERRQ(ierr); }
8641473654SJoe Wallwork     }
8741473654SJoe Wallwork 
8841473654SJoe Wallwork     /* Face tags */
8941473654SJoe Wallwork     ierr = DMCreateLabel(dm, "Face Sets");CHKERRQ(ierr);
9041473654SJoe Wallwork     ierr = DMGetLabel(dm, "Face Sets", &bdLabel);CHKERRQ(ierr);
9141473654SJoe Wallwork     ierr = DMPlexMarkBoundaryFaces(dm, 1, bdLabel);CHKERRQ(ierr);
9241473654SJoe Wallwork     ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
9341473654SJoe Wallwork     ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
9441473654SJoe Wallwork     ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
9541473654SJoe Wallwork     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
9641473654SJoe Wallwork     ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
9741473654SJoe Wallwork     for (f = fStart; f < fEnd; ++f) {
9841473654SJoe Wallwork       PetscBool flg = PETSC_TRUE;
9941473654SJoe Wallwork       PetscInt *closure = NULL, closureSize, cl;
10041473654SJoe Wallwork       PetscReal eps = 1.0e-08;
10141473654SJoe Wallwork 
10241473654SJoe Wallwork       ierr = DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
10341473654SJoe Wallwork       for (cl = 0; cl < closureSize*2; cl += 2) {
10441473654SJoe Wallwork         PetscInt   off = closure[cl];
10541473654SJoe Wallwork         PetscReal *x;
10641473654SJoe Wallwork 
10741473654SJoe Wallwork         if ((off < vStart) || (off >= vEnd)) continue;
10841473654SJoe Wallwork         ierr = DMPlexPointLocalRead(cdm, off, coords, &x);CHKERRQ(ierr);
10941473654SJoe Wallwork         if ((x[0] < 0.5 - eps) || (x[0] > 0.5 + eps)) flg = PETSC_FALSE;
11041473654SJoe Wallwork       }
11141473654SJoe Wallwork       if (flg) { ierr = DMLabelSetValue(bdLabel, f, 2);CHKERRQ(ierr); }
11241473654SJoe Wallwork       ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
11341473654SJoe Wallwork     }
11441473654SJoe Wallwork     ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
11541473654SJoe Wallwork   }
11641473654SJoe Wallwork 
11700d473e3SJoe Wallwork   /* Construct metric */
118d8b80fabSJoe Wallwork   ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr);
119d8b80fabSJoe Wallwork   ierr = DMPlexMetricIsUniform(dm, &uniform);CHKERRQ(ierr);
120d8b80fabSJoe Wallwork   ierr = DMPlexMetricIsIsotropic(dm, &isotropic);CHKERRQ(ierr);
12100d473e3SJoe Wallwork   if (uniform) {
122d8b80fabSJoe Wallwork     ierr = DMPlexMetricCreateUniform(dm, 0, scaling, &metric);CHKERRQ(ierr);
12300d473e3SJoe Wallwork   }
12400d473e3SJoe Wallwork   else {
12500d473e3SJoe Wallwork     DM  dmIndi;
12600d473e3SJoe Wallwork     Vec indicator;
12700d473e3SJoe Wallwork 
12800d473e3SJoe Wallwork     /* Construct "error indicator" */
129d8b80fabSJoe Wallwork     ierr = CreateIndicator(dm, &indicator, &dmIndi);CHKERRQ(ierr);
130e5697243SJoe Wallwork     if (isotropic) {
13100d473e3SJoe Wallwork 
13200d473e3SJoe Wallwork       /* Isotropic case: just specify unity */
13300d473e3SJoe Wallwork       ierr = VecSet(indicator, scaling);CHKERRQ(ierr);
13400d473e3SJoe Wallwork       ierr = DMPlexMetricCreateIsotropic(dm, 0, indicator, &metric);CHKERRQ(ierr);
13500d473e3SJoe Wallwork 
13600d473e3SJoe Wallwork     } else {
137d8b80fabSJoe Wallwork       PetscFE fe;
13800d473e3SJoe Wallwork 
13900d473e3SJoe Wallwork       /* 'Anisotropic' case: approximate the identity by recovering the Hessian of a parabola */
14000d473e3SJoe Wallwork       DM               dmGrad;
14100d473e3SJoe Wallwork       PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar*, void*) = {bowl};
14200d473e3SJoe Wallwork       Vec              gradient;
14300d473e3SJoe Wallwork 
14400d473e3SJoe Wallwork       /* Project the parabola into P1 space */
14500d473e3SJoe Wallwork       ierr = DMProjectFunctionLocal(dmIndi, 0.0, funcs, NULL, INSERT_ALL_VALUES, indicator);CHKERRQ(ierr);
14600d473e3SJoe Wallwork 
14700d473e3SJoe Wallwork       /* Approximate the gradient */
14800d473e3SJoe Wallwork       ierr = DMClone(dmIndi, &dmGrad);CHKERRQ(ierr);
14900d473e3SJoe Wallwork       ierr = PetscFECreateLagrange(comm, dim, dim, PETSC_TRUE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr);
1505767498bSJoe Wallwork       ierr = DMSetField(dmGrad, 0, NULL, (PetscObject)fe);CHKERRQ(ierr);
15100d473e3SJoe Wallwork       ierr = DMCreateDS(dmGrad);CHKERRQ(ierr);
15200d473e3SJoe Wallwork       ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
15300d473e3SJoe Wallwork       ierr = DMCreateLocalVector(dmGrad, &gradient);CHKERRQ(ierr);
15400d473e3SJoe Wallwork       ierr = DMPlexComputeGradientClementInterpolant(dmIndi, indicator, gradient);CHKERRQ(ierr);
15500d473e3SJoe Wallwork       ierr = VecViewFromOptions(gradient, NULL, "-adapt_gradient_view");CHKERRQ(ierr);
15600d473e3SJoe Wallwork 
15700d473e3SJoe Wallwork       /* Approximate the Hessian */
1585767498bSJoe Wallwork       ierr = DMPlexMetricCreate(dm, 0, &metric);CHKERRQ(ierr);
15900d473e3SJoe Wallwork       ierr = DMPlexComputeGradientClementInterpolant(dmGrad, gradient, metric);CHKERRQ(ierr);
16000d473e3SJoe Wallwork       ierr = VecViewFromOptions(metric, NULL, "-adapt_hessian_view");CHKERRQ(ierr);
16100d473e3SJoe Wallwork       ierr = VecDestroy(&gradient);CHKERRQ(ierr);
16200d473e3SJoe Wallwork       ierr = DMDestroy(&dmGrad);CHKERRQ(ierr);
16300d473e3SJoe Wallwork     }
16400d473e3SJoe Wallwork     ierr = VecDestroy(&indicator);CHKERRQ(ierr);
16500d473e3SJoe Wallwork     ierr = DMDestroy(&dmIndi);CHKERRQ(ierr);
16600d473e3SJoe Wallwork   }
16700d473e3SJoe Wallwork 
16800d473e3SJoe Wallwork   /* Test metric routines */
16900d473e3SJoe Wallwork   {
170f49e87b0SJoe Wallwork     DM        dmDet;
17100d473e3SJoe Wallwork     PetscReal errornorm, norm, tol = 1.0e-10, weights[2] = {0.8, 0.2};
172f49e87b0SJoe Wallwork     Vec       metric1, metric2, metricComb, determinant;
17300d473e3SJoe Wallwork     Vec       metrics[2];
17400d473e3SJoe Wallwork 
17500d473e3SJoe Wallwork     ierr = VecDuplicate(metric, &metric1);CHKERRQ(ierr);
17600d473e3SJoe Wallwork     ierr = VecSet(metric1, 0);CHKERRQ(ierr);
17700d473e3SJoe Wallwork     ierr = VecAXPY(metric1, 0.625, metric);CHKERRQ(ierr);
17800d473e3SJoe Wallwork     ierr = VecDuplicate(metric, &metric2);CHKERRQ(ierr);
17900d473e3SJoe Wallwork     ierr = VecSet(metric2, 0);CHKERRQ(ierr);
18000d473e3SJoe Wallwork     ierr = VecAXPY(metric2, 2.5, metric);CHKERRQ(ierr);
18100d473e3SJoe Wallwork     metrics[0] = metric1;
18200d473e3SJoe Wallwork     metrics[1] = metric2;
18300d473e3SJoe Wallwork 
18400d473e3SJoe Wallwork     /* Test metric average */
18500d473e3SJoe Wallwork     ierr = DMPlexMetricAverage(dm, 2, weights, metrics, &metricComb);CHKERRQ(ierr);
18600d473e3SJoe Wallwork     ierr = VecAXPY(metricComb, -1, metric);CHKERRQ(ierr);
18700d473e3SJoe Wallwork     ierr = VecNorm(metric, NORM_2, &norm);CHKERRQ(ierr);
18800d473e3SJoe Wallwork     ierr = VecNorm(metricComb, NORM_2, &errornorm);CHKERRQ(ierr);
18900d473e3SJoe Wallwork     errornorm /= norm;
1907cdd5544SJoe Wallwork     ierr = PetscPrintf(comm, "Metric average L2 error: %.4f%%\n", 100*errornorm);CHKERRQ(ierr);
191*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(errornorm > tol,comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric average test failed");
19200d473e3SJoe Wallwork     ierr = VecDestroy(&metricComb);CHKERRQ(ierr);
19300d473e3SJoe Wallwork 
19400d473e3SJoe Wallwork     /* Test metric intersection */
195e5697243SJoe Wallwork     if (isotropic) {
19600d473e3SJoe Wallwork       ierr = DMPlexMetricIntersection(dm, 2, metrics, &metricComb);CHKERRQ(ierr);
19700d473e3SJoe Wallwork       ierr = VecAXPY(metricComb, -1, metric1);CHKERRQ(ierr);
19800d473e3SJoe Wallwork       ierr = VecNorm(metricComb, NORM_2, &errornorm);CHKERRQ(ierr);
19900d473e3SJoe Wallwork       errornorm /= norm;
2007cdd5544SJoe Wallwork       ierr = PetscPrintf(comm, "Metric intersection L2 error: %.4f%%\n", 100*errornorm);CHKERRQ(ierr);
201*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(errornorm > tol,comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric intersection test failed");
20200d473e3SJoe Wallwork     }
203f49e87b0SJoe Wallwork     ierr = VecDestroy(&metric1);CHKERRQ(ierr);
20400d473e3SJoe Wallwork     ierr = VecDestroy(&metric2);CHKERRQ(ierr);
20500d473e3SJoe Wallwork     ierr = VecDestroy(&metricComb);CHKERRQ(ierr);
20600d473e3SJoe Wallwork 
20700d473e3SJoe Wallwork     /* Test metric SPD enforcement */
208f49e87b0SJoe Wallwork     ierr = DMPlexMetricEnforceSPD(dm, metric, PETSC_TRUE, PETSC_TRUE, &metric1, &determinant);CHKERRQ(ierr);
209e5697243SJoe Wallwork     if (isotropic) {
210f49e87b0SJoe Wallwork       Vec err;
211f49e87b0SJoe Wallwork 
212f49e87b0SJoe Wallwork       ierr = VecDuplicate(determinant, &err);CHKERRQ(ierr);
213f49e87b0SJoe Wallwork       ierr = VecSet(err, 1.0);CHKERRQ(ierr);
214f49e87b0SJoe Wallwork       ierr = VecNorm(err, NORM_2, &norm);CHKERRQ(ierr);
215f49e87b0SJoe Wallwork       ierr = VecAXPY(err, -1, determinant);CHKERRQ(ierr);
216f49e87b0SJoe Wallwork       ierr = VecNorm(err, NORM_2, &errornorm);CHKERRQ(ierr);
217f49e87b0SJoe Wallwork       ierr = VecDestroy(&err);CHKERRQ(ierr);
218f49e87b0SJoe Wallwork       errornorm /= norm;
2197cdd5544SJoe Wallwork       ierr = PetscPrintf(comm, "Metric determinant L2 error: %.4f%%\n", 100*errornorm);CHKERRQ(ierr);
220*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(errornorm > tol,comm, PETSC_ERR_ARG_OUTOFRANGE, "Determinant is not unit");
22100d473e3SJoe Wallwork       ierr = VecAXPY(metric1, -1, metric);CHKERRQ(ierr);
22200d473e3SJoe Wallwork       ierr = VecNorm(metric1, NORM_2, &errornorm);CHKERRQ(ierr);
22300d473e3SJoe Wallwork       errornorm /= norm;
2247cdd5544SJoe Wallwork       ierr = PetscPrintf(comm, "Metric SPD enforcement L2 error: %.4f%%\n", 100*errornorm);CHKERRQ(ierr);
225*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(errornorm > tol,comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric SPD enforcement test failed");
22600d473e3SJoe Wallwork     }
22700d473e3SJoe Wallwork     ierr = VecDestroy(&metric1);CHKERRQ(ierr);
228f49e87b0SJoe Wallwork     ierr = VecGetDM(determinant, &dmDet);CHKERRQ(ierr);
229f49e87b0SJoe Wallwork     ierr = VecDestroy(&determinant);CHKERRQ(ierr);
230f49e87b0SJoe Wallwork     ierr = DMDestroy(&dmDet);CHKERRQ(ierr);
23100d473e3SJoe Wallwork 
23200d473e3SJoe Wallwork     /* Test metric normalization */
233e5697243SJoe Wallwork     ierr = DMPlexMetricNormalize(dm, metric, PETSC_TRUE, PETSC_TRUE, &metric1);CHKERRQ(ierr);
234e5697243SJoe Wallwork     if (isotropic) {
235e5697243SJoe Wallwork       PetscReal target;
236e5697243SJoe Wallwork 
237e5697243SJoe Wallwork       ierr = DMPlexMetricGetTargetComplexity(dm, &target);CHKERRQ(ierr);
238e5697243SJoe Wallwork       scaling = PetscPowReal(target, 2.0/dim);
239d8b80fabSJoe Wallwork       if (uniform) {
24000d473e3SJoe Wallwork         ierr = DMPlexMetricCreateUniform(dm, 0, scaling, &metric2);CHKERRQ(ierr);
241d8b80fabSJoe Wallwork       } else {
242d8b80fabSJoe Wallwork         DM  dmIndi;
243d8b80fabSJoe Wallwork         Vec indicator;
244d8b80fabSJoe Wallwork 
245d8b80fabSJoe Wallwork         ierr = CreateIndicator(dm, &indicator, &dmIndi);CHKERRQ(ierr);
246d8b80fabSJoe Wallwork         ierr = VecSet(indicator, scaling);CHKERRQ(ierr);
247d8b80fabSJoe Wallwork         ierr = DMPlexMetricCreateIsotropic(dm, 0, indicator, &metric2);CHKERRQ(ierr);
248d8b80fabSJoe Wallwork         ierr = DMDestroy(&dmIndi);CHKERRQ(ierr);
249d8b80fabSJoe Wallwork         ierr = VecDestroy(&indicator);CHKERRQ(ierr);
250d8b80fabSJoe Wallwork       }
25100d473e3SJoe Wallwork       ierr = VecAXPY(metric2, -1, metric1);CHKERRQ(ierr);
25200d473e3SJoe Wallwork       ierr = VecNorm(metric2, NORM_2, &errornorm);CHKERRQ(ierr);
25300d473e3SJoe Wallwork       errornorm /= norm;
2547cdd5544SJoe Wallwork       ierr = PetscPrintf(comm, "Metric normalization L2 error: %.4f%%\n", 100*errornorm);CHKERRQ(ierr);
255*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(errornorm > tol,comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric normalization test failed");
25600d473e3SJoe Wallwork     }
25700d473e3SJoe Wallwork     ierr = VecCopy(metric1, metric);CHKERRQ(ierr);
25800d473e3SJoe Wallwork     ierr = VecDestroy(&metric2);CHKERRQ(ierr);
25900d473e3SJoe Wallwork     ierr = VecDestroy(&metric1);CHKERRQ(ierr);
26000d473e3SJoe Wallwork   }
26100d473e3SJoe Wallwork 
26200d473e3SJoe Wallwork   /* Adapt the mesh */
26321b3e102SJoe Wallwork   ierr = DMAdaptMetric(dm, metric, bdLabel, rgLabel, &dmAdapt);CHKERRQ(ierr);
26421b3e102SJoe Wallwork   ierr = DMDestroy(&dm);CHKERRQ(ierr);
26500d473e3SJoe Wallwork   ierr = PetscObjectSetName((PetscObject) dmAdapt, "DM_adapted");CHKERRQ(ierr);
26600d473e3SJoe Wallwork   ierr = VecDestroy(&metric);CHKERRQ(ierr);
26700d473e3SJoe Wallwork   ierr = DMViewFromOptions(dmAdapt, NULL, "-adapted_mesh_view");CHKERRQ(ierr);
26800d473e3SJoe Wallwork 
26941473654SJoe Wallwork   /* Test tag preservation */
27041473654SJoe Wallwork   if (!noTagging) {
27141473654SJoe Wallwork     PetscBool hasTag;
27241473654SJoe Wallwork     PetscInt  size;
27341473654SJoe Wallwork 
27441473654SJoe Wallwork     ierr = DMGetLabel(dmAdapt, "Face Sets", &bdLabel);CHKERRQ(ierr);
27541473654SJoe Wallwork     ierr = DMLabelHasStratum(bdLabel, 1, &hasTag);CHKERRQ(ierr);
276*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!hasTag,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have face tag 1");
27741473654SJoe Wallwork     ierr = DMLabelHasStratum(bdLabel, 2, &hasTag);CHKERRQ(ierr);
278*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!hasTag,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have face tag 2");
27941473654SJoe Wallwork     ierr = DMLabelGetNumValues(bdLabel, &size);CHKERRQ(ierr);
280*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(size != 2,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh has the wrong number of face tags (got %d, expected 2)", size);
28141473654SJoe Wallwork 
28241473654SJoe Wallwork     ierr = DMGetLabel(dmAdapt, "Cell Sets", &rgLabel);CHKERRQ(ierr);
28341473654SJoe Wallwork     ierr = DMLabelHasStratum(rgLabel, 3, &hasTag);CHKERRQ(ierr);
284*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!hasTag,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have cell tag 3");
28541473654SJoe Wallwork     ierr = DMLabelHasStratum(rgLabel, 4, &hasTag);CHKERRQ(ierr);
286*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!hasTag,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have cell tag 4");
28741473654SJoe Wallwork     ierr = DMLabelGetNumValues(rgLabel, &size);CHKERRQ(ierr);
288*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(size != 2,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh has the wrong number of cell tags (got %d, expected 2)", size);
28941473654SJoe Wallwork   }
29041473654SJoe Wallwork 
29100d473e3SJoe Wallwork   /* Clean up */
29200d473e3SJoe Wallwork   ierr = DMDestroy(&dmAdapt);CHKERRQ(ierr);
29300d473e3SJoe Wallwork   ierr = PetscFinalize();
29400d473e3SJoe Wallwork   return 0;
29500d473e3SJoe Wallwork }
29600d473e3SJoe Wallwork 
29700d473e3SJoe Wallwork /*TEST
29800d473e3SJoe Wallwork 
29980f7b1e7SJoe Wallwork   testset:
30080f7b1e7SJoe Wallwork     requires: pragmatic
30180f7b1e7SJoe Wallwork     args: -dm_plex_metric_target_complexity 100 -dm_adaptor pragmatic -noTagging -dim 2
30280f7b1e7SJoe Wallwork 
303dd0671eeSJoe Wallwork     test:
304dd0671eeSJoe Wallwork       suffix: uniform_2d_pragmatic
305d8b80fabSJoe Wallwork       args: -dm_plex_metric_uniform
30600d473e3SJoe Wallwork     test:
307dd0671eeSJoe Wallwork       suffix: iso_2d_pragmatic
308d8b80fabSJoe Wallwork       args: -dm_plex_metric_isotropic
30900d473e3SJoe Wallwork     test:
310dd0671eeSJoe Wallwork       suffix: hessian_2d_pragmatic
31180f7b1e7SJoe Wallwork 
31280f7b1e7SJoe Wallwork   testset:
31380f7b1e7SJoe Wallwork     requires: pragmatic tetgen
31480f7b1e7SJoe Wallwork     args: -dm_plex_metric_target_complexity 100 -dm_adaptor pragmatic -noTagging -dim 3
31580f7b1e7SJoe Wallwork 
31680f7b1e7SJoe Wallwork     test:
31780f7b1e7SJoe Wallwork       suffix: uniform_3d_pragmatic
318d8b80fabSJoe Wallwork       args: -dm_plex_metric_uniform -noTagging
31980f7b1e7SJoe Wallwork     test:
32080f7b1e7SJoe Wallwork       suffix: iso_3d_pragmatic
321d8b80fabSJoe Wallwork       args: -dm_plex_metric_isotropic -noTagging
32200d473e3SJoe Wallwork     test:
323dd0671eeSJoe Wallwork       suffix: hessian_3d_pragmatic
32480f7b1e7SJoe Wallwork 
32580f7b1e7SJoe Wallwork   testset:
32680f7b1e7SJoe Wallwork     requires: mmg
32780f7b1e7SJoe Wallwork     args: -dm_plex_metric_target_complexity 100 -dm_adaptor mmg -dim 2
32880f7b1e7SJoe Wallwork 
32900d473e3SJoe Wallwork     test:
330dd0671eeSJoe Wallwork       suffix: uniform_2d_mmg
331d8b80fabSJoe Wallwork       args: -dm_plex_metric_uniform
332dd0671eeSJoe Wallwork     test:
333dd0671eeSJoe Wallwork       suffix: iso_2d_mmg
334d8b80fabSJoe Wallwork       args: -dm_plex_metric_isotropic
335dd0671eeSJoe Wallwork     test:
336dd0671eeSJoe Wallwork       suffix: hessian_2d_mmg
33780f7b1e7SJoe Wallwork 
33880f7b1e7SJoe Wallwork   testset:
33980f7b1e7SJoe Wallwork     requires: mmg tetgen
34080f7b1e7SJoe Wallwork     args: -dm_plex_metric_target_complexity 100 -dm_adaptor mmg -dim 3
34180f7b1e7SJoe Wallwork 
34280f7b1e7SJoe Wallwork     test:
34380f7b1e7SJoe Wallwork       suffix: uniform_3d_mmg
344d8b80fabSJoe Wallwork       args: -dm_plex_metric_uniform
34580f7b1e7SJoe Wallwork     test:
34680f7b1e7SJoe Wallwork       suffix: iso_3d_mmg
347d8b80fabSJoe Wallwork       args: -dm_plex_metric_isotropic
348dd0671eeSJoe Wallwork     test:
349dd0671eeSJoe Wallwork       suffix: hessian_3d_mmg
35080f7b1e7SJoe Wallwork 
35180f7b1e7SJoe Wallwork   testset:
35280f7b1e7SJoe Wallwork     requires: parmmg tetgen
35380f7b1e7SJoe Wallwork     nsize: 2
35480f7b1e7SJoe Wallwork     args: -dm_plex_metric_target_complexity 100 -dm_adaptor parmmg -dim 3
35580f7b1e7SJoe Wallwork 
356dd0671eeSJoe Wallwork     test:
357dd0671eeSJoe Wallwork       suffix: uniform_3d_parmmg
358d8b80fabSJoe Wallwork       args: -dm_plex_metric_uniform
359dd0671eeSJoe Wallwork     test:
360dd0671eeSJoe Wallwork       suffix: iso_3d_parmmg
361d8b80fabSJoe Wallwork       args: -dm_plex_metric_isotropic
362dd0671eeSJoe Wallwork     test:
363dd0671eeSJoe Wallwork       suffix: hessian_3d_parmmg
36400d473e3SJoe Wallwork 
36500d473e3SJoe Wallwork TEST*/
366