xref: /petsc/src/dm/impls/plex/tests/ex60.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
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   PetscFE        fe;
19d8b80fabSJoe Wallwork   PetscInt       dim;
20d8b80fabSJoe Wallwork 
21d8b80fabSJoe Wallwork   PetscFunctionBeginUser;
225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject)dm, &comm));
235f80ce2aSJacob Faibussowitsch   CHKERRQ(DMClone(dm, dmIndi));
245f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFECreateLagrange(comm, dim, 1, PETSC_TRUE, 1, PETSC_DETERMINE, &fe));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetField(*dmIndi, 0, NULL, (PetscObject)fe));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateDS(*dmIndi));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEDestroy(&fe));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateLocalVector(*dmIndi, indicator));
30d8b80fabSJoe Wallwork   PetscFunctionReturn(0);
31d8b80fabSJoe Wallwork }
32d8b80fabSJoe Wallwork 
3300d473e3SJoe Wallwork int main(int argc, char **argv) {
34e600fa54SMatthew G. Knepley   DM              dm, dmAdapt;
3521b3e102SJoe Wallwork   DMLabel         bdLabel = NULL, rgLabel = NULL;
3600d473e3SJoe Wallwork   MPI_Comm        comm;
3741473654SJoe Wallwork   PetscBool       uniform = PETSC_FALSE, isotropic = PETSC_FALSE, noTagging = PETSC_FALSE;
3800d473e3SJoe Wallwork   PetscErrorCode  ierr;
39e600fa54SMatthew G. Knepley   PetscInt        dim;
4000d473e3SJoe Wallwork   PetscReal       scaling = 1.0;
4100d473e3SJoe Wallwork   Vec             metric;
4200d473e3SJoe Wallwork 
4300d473e3SJoe Wallwork   /* Set up */
4400d473e3SJoe Wallwork   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
4500d473e3SJoe Wallwork   comm = PETSC_COMM_WORLD;
4600d473e3SJoe Wallwork   ierr = PetscOptionsBegin(comm, "", "Mesh adaptation options", "DMPLEX");CHKERRQ(ierr);
475f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-noTagging", "Should tag preservation testing be turned off?", "ex60.c", noTagging, &noTagging, NULL));
4800d473e3SJoe Wallwork   ierr = PetscOptionsEnd();
4900d473e3SJoe Wallwork 
5000d473e3SJoe Wallwork   /* Create box mesh */
515f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(comm, &dm));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(dm, DMPLEX));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(dm));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) dm, "DM_init"));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(dm, NULL, "-initial_mesh_view"));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
5700d473e3SJoe Wallwork 
5841473654SJoe Wallwork   /* Set tags to be preserved */
5941473654SJoe Wallwork   if (!noTagging) {
6041473654SJoe Wallwork     DM                 cdm;
6141473654SJoe Wallwork     PetscInt           cStart, cEnd, c, fStart, fEnd, f, vStart, vEnd;
6241473654SJoe Wallwork     const PetscScalar *coords;
6341473654SJoe Wallwork     Vec                coordinates;
6441473654SJoe Wallwork 
6541473654SJoe Wallwork     /* Cell tags */
665f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCreateLabel(dm, "Cell Sets"));
675f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetLabel(dm, "Cell Sets", &rgLabel));
685f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
6941473654SJoe Wallwork     for (c = cStart; c < cEnd; ++c) {
7041473654SJoe Wallwork       PetscReal centroid[3], volume, x;
7141473654SJoe Wallwork 
725f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexComputeCellGeometryFVM(dm, c, &volume, centroid, NULL));
7341473654SJoe Wallwork       x = centroid[0];
745f80ce2aSJacob Faibussowitsch       if (x < 0.5) CHKERRQ(DMLabelSetValue(rgLabel, c, 3));
755f80ce2aSJacob Faibussowitsch       else         CHKERRQ(DMLabelSetValue(rgLabel, c, 4));
7641473654SJoe Wallwork     }
7741473654SJoe Wallwork 
7841473654SJoe Wallwork     /* Face tags */
795f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCreateLabel(dm, "Face Sets"));
805f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetLabel(dm, "Face Sets", &bdLabel));
815f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexMarkBoundaryFaces(dm, 1, bdLabel));
825f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
835f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
845f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetCoordinateDM(dm, &cdm));
855f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetCoordinatesLocal(dm, &coordinates));
865f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArrayRead(coordinates, &coords));
8741473654SJoe Wallwork     for (f = fStart; f < fEnd; ++f) {
8841473654SJoe Wallwork       PetscBool flg = PETSC_TRUE;
8941473654SJoe Wallwork       PetscInt *closure = NULL, closureSize, cl;
9041473654SJoe Wallwork       PetscReal eps = 1.0e-08;
9141473654SJoe Wallwork 
925f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
9341473654SJoe Wallwork       for (cl = 0; cl < closureSize*2; cl += 2) {
9441473654SJoe Wallwork         PetscInt   off = closure[cl];
9541473654SJoe Wallwork         PetscReal *x;
9641473654SJoe Wallwork 
9741473654SJoe Wallwork         if ((off < vStart) || (off >= vEnd)) continue;
985f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexPointLocalRead(cdm, off, coords, &x));
9941473654SJoe Wallwork         if ((x[0] < 0.5 - eps) || (x[0] > 0.5 + eps)) flg = PETSC_FALSE;
10041473654SJoe Wallwork       }
1015f80ce2aSJacob Faibussowitsch       if (flg) CHKERRQ(DMLabelSetValue(bdLabel, f, 2));
1025f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
10341473654SJoe Wallwork     }
1045f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArrayRead(coordinates, &coords));
10541473654SJoe Wallwork   }
10641473654SJoe Wallwork 
10700d473e3SJoe Wallwork   /* Construct metric */
1085f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexMetricSetFromOptions(dm));
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexMetricIsUniform(dm, &uniform));
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexMetricIsIsotropic(dm, &isotropic));
11100d473e3SJoe Wallwork   if (uniform) {
1125f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexMetricCreateUniform(dm, 0, scaling, &metric));
11300d473e3SJoe Wallwork   }
11400d473e3SJoe Wallwork   else {
11500d473e3SJoe Wallwork     DM  dmIndi;
11600d473e3SJoe Wallwork     Vec indicator;
11700d473e3SJoe Wallwork 
11800d473e3SJoe Wallwork     /* Construct "error indicator" */
1195f80ce2aSJacob Faibussowitsch     CHKERRQ(CreateIndicator(dm, &indicator, &dmIndi));
120e5697243SJoe Wallwork     if (isotropic) {
12100d473e3SJoe Wallwork 
12200d473e3SJoe Wallwork       /* Isotropic case: just specify unity */
1235f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSet(indicator, scaling));
1245f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexMetricCreateIsotropic(dm, 0, indicator, &metric));
12500d473e3SJoe Wallwork 
12600d473e3SJoe Wallwork     } else {
127d8b80fabSJoe Wallwork       PetscFE fe;
12800d473e3SJoe Wallwork 
12900d473e3SJoe Wallwork       /* 'Anisotropic' case: approximate the identity by recovering the Hessian of a parabola */
13000d473e3SJoe Wallwork       DM               dmGrad;
13100d473e3SJoe Wallwork       PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar*, void*) = {bowl};
13200d473e3SJoe Wallwork       Vec              gradient;
13300d473e3SJoe Wallwork 
13400d473e3SJoe Wallwork       /* Project the parabola into P1 space */
1355f80ce2aSJacob Faibussowitsch       CHKERRQ(DMProjectFunctionLocal(dmIndi, 0.0, funcs, NULL, INSERT_ALL_VALUES, indicator));
13600d473e3SJoe Wallwork 
13700d473e3SJoe Wallwork       /* Approximate the gradient */
1385f80ce2aSJacob Faibussowitsch       CHKERRQ(DMClone(dmIndi, &dmGrad));
1395f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFECreateLagrange(comm, dim, dim, PETSC_TRUE, 1, PETSC_DETERMINE, &fe));
1405f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSetField(dmGrad, 0, NULL, (PetscObject)fe));
1415f80ce2aSJacob Faibussowitsch       CHKERRQ(DMCreateDS(dmGrad));
1425f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFEDestroy(&fe));
1435f80ce2aSJacob Faibussowitsch       CHKERRQ(DMCreateLocalVector(dmGrad, &gradient));
1445f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexComputeGradientClementInterpolant(dmIndi, indicator, gradient));
1455f80ce2aSJacob Faibussowitsch       CHKERRQ(VecViewFromOptions(gradient, NULL, "-adapt_gradient_view"));
14600d473e3SJoe Wallwork 
14700d473e3SJoe Wallwork       /* Approximate the Hessian */
1485f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexMetricCreate(dm, 0, &metric));
1495f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexComputeGradientClementInterpolant(dmGrad, gradient, metric));
1505f80ce2aSJacob Faibussowitsch       CHKERRQ(VecViewFromOptions(metric, NULL, "-adapt_hessian_view"));
1515f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&gradient));
1525f80ce2aSJacob Faibussowitsch       CHKERRQ(DMDestroy(&dmGrad));
15300d473e3SJoe Wallwork     }
1545f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&indicator));
1555f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDestroy(&dmIndi));
15600d473e3SJoe Wallwork   }
15700d473e3SJoe Wallwork 
15800d473e3SJoe Wallwork   /* Test metric routines */
15900d473e3SJoe Wallwork   {
160f49e87b0SJoe Wallwork     DM        dmDet;
16100d473e3SJoe Wallwork     PetscReal errornorm, norm, tol = 1.0e-10, weights[2] = {0.8, 0.2};
162f49e87b0SJoe Wallwork     Vec       metric1, metric2, metricComb, determinant;
16300d473e3SJoe Wallwork     Vec       metrics[2];
16400d473e3SJoe Wallwork 
1655f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(metric, &metric1));
1665f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(metric1, 0));
1675f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(metric1, 0.625, metric));
1685f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(metric, &metric2));
1695f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(metric2, 0));
1705f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(metric2, 2.5, metric));
17100d473e3SJoe Wallwork     metrics[0] = metric1;
17200d473e3SJoe Wallwork     metrics[1] = metric2;
17300d473e3SJoe Wallwork 
17400d473e3SJoe Wallwork     /* Test metric average */
1755f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexMetricAverage(dm, 2, weights, metrics, &metricComb));
1765f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(metricComb, -1, metric));
1775f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(metric, NORM_2, &norm));
1785f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(metricComb, NORM_2, &errornorm));
17900d473e3SJoe Wallwork     errornorm /= norm;
1805f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm, "Metric average L2 error: %.4f%%\n", 100*errornorm));
1812c71b3e2SJacob Faibussowitsch     PetscCheckFalse(errornorm > tol,comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric average test failed");
1825f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&metricComb));
18300d473e3SJoe Wallwork 
18400d473e3SJoe Wallwork     /* Test metric intersection */
185e5697243SJoe Wallwork     if (isotropic) {
1865f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexMetricIntersection(dm, 2, metrics, &metricComb));
1875f80ce2aSJacob Faibussowitsch       CHKERRQ(VecAXPY(metricComb, -1, metric1));
1885f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(metricComb, NORM_2, &errornorm));
18900d473e3SJoe Wallwork       errornorm /= norm;
1905f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(comm, "Metric intersection L2 error: %.4f%%\n", 100*errornorm));
1912c71b3e2SJacob Faibussowitsch       PetscCheckFalse(errornorm > tol,comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric intersection test failed");
19200d473e3SJoe Wallwork     }
1935f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&metric1));
1945f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&metric2));
1955f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&metricComb));
19600d473e3SJoe Wallwork 
19700d473e3SJoe Wallwork     /* Test metric SPD enforcement */
1985f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexMetricEnforceSPD(dm, metric, PETSC_TRUE, PETSC_TRUE, &metric1, &determinant));
199e5697243SJoe Wallwork     if (isotropic) {
200f49e87b0SJoe Wallwork       Vec err;
201f49e87b0SJoe Wallwork 
2025f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDuplicate(determinant, &err));
2035f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSet(err, 1.0));
2045f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(err, NORM_2, &norm));
2055f80ce2aSJacob Faibussowitsch       CHKERRQ(VecAXPY(err, -1, determinant));
2065f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(err, NORM_2, &errornorm));
2075f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&err));
208f49e87b0SJoe Wallwork       errornorm /= norm;
2095f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(comm, "Metric determinant L2 error: %.4f%%\n", 100*errornorm));
2102c71b3e2SJacob Faibussowitsch       PetscCheckFalse(errornorm > tol,comm, PETSC_ERR_ARG_OUTOFRANGE, "Determinant is not unit");
2115f80ce2aSJacob Faibussowitsch       CHKERRQ(VecAXPY(metric1, -1, metric));
2125f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(metric1, NORM_2, &errornorm));
21300d473e3SJoe Wallwork       errornorm /= norm;
2145f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(comm, "Metric SPD enforcement L2 error: %.4f%%\n", 100*errornorm));
2152c71b3e2SJacob Faibussowitsch       PetscCheckFalse(errornorm > tol,comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric SPD enforcement test failed");
21600d473e3SJoe Wallwork     }
2175f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&metric1));
2185f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetDM(determinant, &dmDet));
2195f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&determinant));
2205f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDestroy(&dmDet));
22100d473e3SJoe Wallwork 
22200d473e3SJoe Wallwork     /* Test metric normalization */
2235f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexMetricNormalize(dm, metric, PETSC_TRUE, PETSC_TRUE, &metric1));
224e5697243SJoe Wallwork     if (isotropic) {
225e5697243SJoe Wallwork       PetscReal target;
226e5697243SJoe Wallwork 
2275f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexMetricGetTargetComplexity(dm, &target));
228e5697243SJoe Wallwork       scaling = PetscPowReal(target, 2.0/dim);
229d8b80fabSJoe Wallwork       if (uniform) {
2305f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexMetricCreateUniform(dm, 0, scaling, &metric2));
231d8b80fabSJoe Wallwork       } else {
232d8b80fabSJoe Wallwork         DM  dmIndi;
233d8b80fabSJoe Wallwork         Vec indicator;
234d8b80fabSJoe Wallwork 
2355f80ce2aSJacob Faibussowitsch         CHKERRQ(CreateIndicator(dm, &indicator, &dmIndi));
2365f80ce2aSJacob Faibussowitsch         CHKERRQ(VecSet(indicator, scaling));
2375f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexMetricCreateIsotropic(dm, 0, indicator, &metric2));
2385f80ce2aSJacob Faibussowitsch         CHKERRQ(DMDestroy(&dmIndi));
2395f80ce2aSJacob Faibussowitsch         CHKERRQ(VecDestroy(&indicator));
240d8b80fabSJoe Wallwork       }
2415f80ce2aSJacob Faibussowitsch       CHKERRQ(VecAXPY(metric2, -1, metric1));
2425f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(metric2, NORM_2, &errornorm));
24300d473e3SJoe Wallwork       errornorm /= norm;
2445f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(comm, "Metric normalization L2 error: %.4f%%\n", 100*errornorm));
2452c71b3e2SJacob Faibussowitsch       PetscCheckFalse(errornorm > tol,comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric normalization test failed");
24600d473e3SJoe Wallwork     }
2475f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(metric1, metric));
2485f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&metric2));
2495f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&metric1));
25000d473e3SJoe Wallwork   }
25100d473e3SJoe Wallwork 
25200d473e3SJoe Wallwork   /* Adapt the mesh */
2535f80ce2aSJacob Faibussowitsch   CHKERRQ(DMAdaptMetric(dm, metric, bdLabel, rgLabel, &dmAdapt));
2545f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
2555f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) dmAdapt, "DM_adapted"));
2565f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&metric));
2575f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(dmAdapt, NULL, "-adapted_mesh_view"));
25800d473e3SJoe Wallwork 
25941473654SJoe Wallwork   /* Test tag preservation */
26041473654SJoe Wallwork   if (!noTagging) {
26141473654SJoe Wallwork     PetscBool hasTag;
26241473654SJoe Wallwork     PetscInt  size;
26341473654SJoe Wallwork 
2645f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetLabel(dmAdapt, "Face Sets", &bdLabel));
2655f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLabelHasStratum(bdLabel, 1, &hasTag));
266*28b400f6SJacob Faibussowitsch     PetscCheck(hasTag,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have face tag 1");
2675f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLabelHasStratum(bdLabel, 2, &hasTag));
268*28b400f6SJacob Faibussowitsch     PetscCheck(hasTag,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have face tag 2");
2695f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLabelGetNumValues(bdLabel, &size));
2702c71b3e2SJacob Faibussowitsch     PetscCheckFalse(size != 2,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh has the wrong number of face tags (got %d, expected 2)", size);
27141473654SJoe Wallwork 
2725f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetLabel(dmAdapt, "Cell Sets", &rgLabel));
2735f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLabelHasStratum(rgLabel, 3, &hasTag));
274*28b400f6SJacob Faibussowitsch     PetscCheck(hasTag,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have cell tag 3");
2755f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLabelHasStratum(rgLabel, 4, &hasTag));
276*28b400f6SJacob Faibussowitsch     PetscCheck(hasTag,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have cell tag 4");
2775f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLabelGetNumValues(rgLabel, &size));
2782c71b3e2SJacob Faibussowitsch     PetscCheckFalse(size != 2,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh has the wrong number of cell tags (got %d, expected 2)", size);
27941473654SJoe Wallwork   }
28041473654SJoe Wallwork 
28100d473e3SJoe Wallwork   /* Clean up */
2825f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dmAdapt));
28300d473e3SJoe Wallwork   ierr = PetscFinalize();
28400d473e3SJoe Wallwork   return 0;
28500d473e3SJoe Wallwork }
28600d473e3SJoe Wallwork 
28700d473e3SJoe Wallwork /*TEST
28800d473e3SJoe Wallwork 
28980f7b1e7SJoe Wallwork   testset:
29080f7b1e7SJoe Wallwork     requires: pragmatic
291e600fa54SMatthew G. Knepley     args: -dm_plex_box_faces 4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor pragmatic -noTagging
29280f7b1e7SJoe Wallwork 
293dd0671eeSJoe Wallwork     test:
294dd0671eeSJoe Wallwork       suffix: uniform_2d_pragmatic
295d8b80fabSJoe Wallwork       args: -dm_plex_metric_uniform
29600d473e3SJoe Wallwork     test:
297dd0671eeSJoe Wallwork       suffix: iso_2d_pragmatic
298d8b80fabSJoe Wallwork       args: -dm_plex_metric_isotropic
29900d473e3SJoe Wallwork     test:
300dd0671eeSJoe Wallwork       suffix: hessian_2d_pragmatic
30180f7b1e7SJoe Wallwork 
30280f7b1e7SJoe Wallwork   testset:
30380f7b1e7SJoe Wallwork     requires: pragmatic tetgen
304e600fa54SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_box_faces 4,4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor pragmatic -noTagging
30580f7b1e7SJoe Wallwork 
30680f7b1e7SJoe Wallwork     test:
30780f7b1e7SJoe Wallwork       suffix: uniform_3d_pragmatic
308d8b80fabSJoe Wallwork       args: -dm_plex_metric_uniform -noTagging
30980f7b1e7SJoe Wallwork     test:
31080f7b1e7SJoe Wallwork       suffix: iso_3d_pragmatic
311d8b80fabSJoe Wallwork       args: -dm_plex_metric_isotropic -noTagging
31200d473e3SJoe Wallwork     test:
313dd0671eeSJoe Wallwork       suffix: hessian_3d_pragmatic
31480f7b1e7SJoe Wallwork 
31580f7b1e7SJoe Wallwork   testset:
31680f7b1e7SJoe Wallwork     requires: mmg
317e600fa54SMatthew G. Knepley     args: -dm_plex_box_faces 4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor mmg
31880f7b1e7SJoe Wallwork 
31900d473e3SJoe Wallwork     test:
320dd0671eeSJoe Wallwork       suffix: uniform_2d_mmg
321d8b80fabSJoe Wallwork       args: -dm_plex_metric_uniform
322dd0671eeSJoe Wallwork     test:
323dd0671eeSJoe Wallwork       suffix: iso_2d_mmg
324d8b80fabSJoe Wallwork       args: -dm_plex_metric_isotropic
325dd0671eeSJoe Wallwork     test:
326dd0671eeSJoe Wallwork       suffix: hessian_2d_mmg
32780f7b1e7SJoe Wallwork 
32880f7b1e7SJoe Wallwork   testset:
32980f7b1e7SJoe Wallwork     requires: mmg tetgen
330e600fa54SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_box_faces 4,4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor mmg
33180f7b1e7SJoe Wallwork 
33280f7b1e7SJoe Wallwork     test:
33380f7b1e7SJoe Wallwork       suffix: uniform_3d_mmg
334d8b80fabSJoe Wallwork       args: -dm_plex_metric_uniform
33580f7b1e7SJoe Wallwork     test:
33680f7b1e7SJoe Wallwork       suffix: iso_3d_mmg
337d8b80fabSJoe Wallwork       args: -dm_plex_metric_isotropic
338dd0671eeSJoe Wallwork     test:
339dd0671eeSJoe Wallwork       suffix: hessian_3d_mmg
34080f7b1e7SJoe Wallwork 
34180f7b1e7SJoe Wallwork   testset:
34280f7b1e7SJoe Wallwork     requires: parmmg tetgen
34380f7b1e7SJoe Wallwork     nsize: 2
344e600fa54SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_box_faces 4,4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor parmmg
34580f7b1e7SJoe Wallwork 
346dd0671eeSJoe Wallwork     test:
347dd0671eeSJoe Wallwork       suffix: uniform_3d_parmmg
348d8b80fabSJoe Wallwork       args: -dm_plex_metric_uniform
349dd0671eeSJoe Wallwork     test:
350dd0671eeSJoe Wallwork       suffix: iso_3d_parmmg
351d8b80fabSJoe Wallwork       args: -dm_plex_metric_isotropic
352dd0671eeSJoe Wallwork     test:
353dd0671eeSJoe Wallwork       suffix: hessian_3d_parmmg
35400d473e3SJoe Wallwork 
35500d473e3SJoe Wallwork TEST*/
356