xref: /petsc/src/dm/impls/plex/tests/ex60.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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;
229566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
239566063dSJacob Faibussowitsch   PetscCall(DMClone(dm, dmIndi));
249566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
259566063dSJacob Faibussowitsch   PetscCall(PetscFECreateLagrange(comm, dim, 1, PETSC_TRUE, 1, PETSC_DETERMINE, &fe));
269566063dSJacob Faibussowitsch   PetscCall(DMSetField(*dmIndi, 0, NULL, (PetscObject)fe));
279566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(*dmIndi));
289566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
299566063dSJacob Faibussowitsch   PetscCall(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;
38e600fa54SMatthew G. Knepley   PetscInt        dim;
3900d473e3SJoe Wallwork   PetscReal       scaling = 1.0;
4000d473e3SJoe Wallwork   Vec             metric;
4100d473e3SJoe Wallwork 
4200d473e3SJoe Wallwork   /* Set up */
43*327415f7SBarry Smith   PetscFunctionBeginUser;
449566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
4500d473e3SJoe Wallwork   comm = PETSC_COMM_WORLD;
46d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Mesh adaptation options", "DMPLEX");
479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-noTagging", "Should tag preservation testing be turned off?", "ex60.c", noTagging, &noTagging, NULL));
48d0609cedSBarry Smith   PetscOptionsEnd();
4900d473e3SJoe Wallwork 
5000d473e3SJoe Wallwork   /* Create box mesh */
519566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, &dm));
529566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
539566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
549566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) dm, "DM_init"));
559566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-initial_mesh_view"));
569566063dSJacob Faibussowitsch   PetscCall(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 */
668fb5bd83SMatthew G. Knepley     PetscCall(DMGetCoordinatesLocalSetUp(dm));
679566063dSJacob Faibussowitsch     PetscCall(DMCreateLabel(dm, "Cell Sets"));
689566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(dm, "Cell Sets", &rgLabel));
699566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
7041473654SJoe Wallwork     for (c = cStart; c < cEnd; ++c) {
7141473654SJoe Wallwork       PetscReal centroid[3], volume, x;
7241473654SJoe Wallwork 
739566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFVM(dm, c, &volume, centroid, NULL));
7441473654SJoe Wallwork       x = centroid[0];
759566063dSJacob Faibussowitsch       if (x < 0.5) PetscCall(DMLabelSetValue(rgLabel, c, 3));
769566063dSJacob Faibussowitsch       else         PetscCall(DMLabelSetValue(rgLabel, c, 4));
7741473654SJoe Wallwork     }
7841473654SJoe Wallwork 
7941473654SJoe Wallwork     /* Face tags */
809566063dSJacob Faibussowitsch     PetscCall(DMCreateLabel(dm, "Face Sets"));
819566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(dm, "Face Sets", &bdLabel));
829566063dSJacob Faibussowitsch     PetscCall(DMPlexMarkBoundaryFaces(dm, 1, bdLabel));
839566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
849566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
859566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dm, &cdm));
869566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
879566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(coordinates, &coords));
8841473654SJoe Wallwork     for (f = fStart; f < fEnd; ++f) {
8941473654SJoe Wallwork       PetscBool flg = PETSC_TRUE;
9041473654SJoe Wallwork       PetscInt *closure = NULL, closureSize, cl;
9141473654SJoe Wallwork       PetscReal eps = 1.0e-08;
9241473654SJoe Wallwork 
939566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
9441473654SJoe Wallwork       for (cl = 0; cl < closureSize*2; cl += 2) {
9541473654SJoe Wallwork         PetscInt   off = closure[cl];
9641473654SJoe Wallwork         PetscReal *x;
9741473654SJoe Wallwork 
9841473654SJoe Wallwork         if ((off < vStart) || (off >= vEnd)) continue;
999566063dSJacob Faibussowitsch         PetscCall(DMPlexPointLocalRead(cdm, off, coords, &x));
10041473654SJoe Wallwork         if ((x[0] < 0.5 - eps) || (x[0] > 0.5 + eps)) flg = PETSC_FALSE;
10141473654SJoe Wallwork       }
1029566063dSJacob Faibussowitsch       if (flg) PetscCall(DMLabelSetValue(bdLabel, f, 2));
1039566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
10441473654SJoe Wallwork     }
1059566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(coordinates, &coords));
10641473654SJoe Wallwork   }
10741473654SJoe Wallwork 
10800d473e3SJoe Wallwork   /* Construct metric */
1099566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetFromOptions(dm));
1109566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1119566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
11200d473e3SJoe Wallwork   if (uniform) {
1139566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricCreateUniform(dm, 0, scaling, &metric));
11400d473e3SJoe Wallwork   }
11500d473e3SJoe Wallwork   else {
11600d473e3SJoe Wallwork     DM  dmIndi;
11700d473e3SJoe Wallwork     Vec indicator;
11800d473e3SJoe Wallwork 
11900d473e3SJoe Wallwork     /* Construct "error indicator" */
1209566063dSJacob Faibussowitsch     PetscCall(CreateIndicator(dm, &indicator, &dmIndi));
121e5697243SJoe Wallwork     if (isotropic) {
12200d473e3SJoe Wallwork 
12300d473e3SJoe Wallwork       /* Isotropic case: just specify unity */
1249566063dSJacob Faibussowitsch       PetscCall(VecSet(indicator, scaling));
1259566063dSJacob Faibussowitsch       PetscCall(DMPlexMetricCreateIsotropic(dm, 0, indicator, &metric));
12600d473e3SJoe Wallwork 
12700d473e3SJoe Wallwork     } else {
128d8b80fabSJoe Wallwork       PetscFE fe;
12900d473e3SJoe Wallwork 
13000d473e3SJoe Wallwork       /* 'Anisotropic' case: approximate the identity by recovering the Hessian of a parabola */
13100d473e3SJoe Wallwork       DM               dmGrad;
13200d473e3SJoe Wallwork       PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar*, void*) = {bowl};
13300d473e3SJoe Wallwork       Vec              gradient;
13400d473e3SJoe Wallwork 
13500d473e3SJoe Wallwork       /* Project the parabola into P1 space */
1369566063dSJacob Faibussowitsch       PetscCall(DMProjectFunctionLocal(dmIndi, 0.0, funcs, NULL, INSERT_ALL_VALUES, indicator));
13700d473e3SJoe Wallwork 
13800d473e3SJoe Wallwork       /* Approximate the gradient */
1399566063dSJacob Faibussowitsch       PetscCall(DMClone(dmIndi, &dmGrad));
1409566063dSJacob Faibussowitsch       PetscCall(PetscFECreateLagrange(comm, dim, dim, PETSC_TRUE, 1, PETSC_DETERMINE, &fe));
1419566063dSJacob Faibussowitsch       PetscCall(DMSetField(dmGrad, 0, NULL, (PetscObject)fe));
1429566063dSJacob Faibussowitsch       PetscCall(DMCreateDS(dmGrad));
1439566063dSJacob Faibussowitsch       PetscCall(PetscFEDestroy(&fe));
1449566063dSJacob Faibussowitsch       PetscCall(DMCreateLocalVector(dmGrad, &gradient));
1459566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeGradientClementInterpolant(dmIndi, indicator, gradient));
1469566063dSJacob Faibussowitsch       PetscCall(VecViewFromOptions(gradient, NULL, "-adapt_gradient_view"));
14700d473e3SJoe Wallwork 
14800d473e3SJoe Wallwork       /* Approximate the Hessian */
1499566063dSJacob Faibussowitsch       PetscCall(DMPlexMetricCreate(dm, 0, &metric));
1509566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeGradientClementInterpolant(dmGrad, gradient, metric));
1519566063dSJacob Faibussowitsch       PetscCall(VecViewFromOptions(metric, NULL, "-adapt_hessian_view"));
1529566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&gradient));
1539566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&dmGrad));
15400d473e3SJoe Wallwork     }
1559566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&indicator));
1569566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dmIndi));
15700d473e3SJoe Wallwork   }
15800d473e3SJoe Wallwork 
15900d473e3SJoe Wallwork   /* Test metric routines */
16000d473e3SJoe Wallwork   {
161f49e87b0SJoe Wallwork     DM        dmDet;
16200d473e3SJoe Wallwork     PetscReal errornorm, norm, tol = 1.0e-10, weights[2] = {0.8, 0.2};
163f49e87b0SJoe Wallwork     Vec       metric1, metric2, metricComb, determinant;
16400d473e3SJoe Wallwork     Vec       metrics[2];
16500d473e3SJoe Wallwork 
1669566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(metric, &metric1));
1679566063dSJacob Faibussowitsch     PetscCall(VecSet(metric1, 0));
1689566063dSJacob Faibussowitsch     PetscCall(VecAXPY(metric1, 0.625, metric));
1699566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(metric, &metric2));
1709566063dSJacob Faibussowitsch     PetscCall(VecSet(metric2, 0));
1719566063dSJacob Faibussowitsch     PetscCall(VecAXPY(metric2, 2.5, metric));
17200d473e3SJoe Wallwork     metrics[0] = metric1;
17300d473e3SJoe Wallwork     metrics[1] = metric2;
17400d473e3SJoe Wallwork 
17500d473e3SJoe Wallwork     /* Test metric average */
17640a2a1afSJoe Wallwork     PetscCall(DMPlexMetricCreate(dm, 0, &metricComb));
17740a2a1afSJoe Wallwork     PetscCall(DMPlexMetricAverage(dm, 2, weights, metrics, metricComb));
1789566063dSJacob Faibussowitsch     PetscCall(VecAXPY(metricComb, -1, metric));
1799566063dSJacob Faibussowitsch     PetscCall(VecNorm(metric, NORM_2, &norm));
1809566063dSJacob Faibussowitsch     PetscCall(VecNorm(metricComb, NORM_2, &errornorm));
18100d473e3SJoe Wallwork     errornorm /= norm;
18263a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Metric average L2 error: %.4f%%\n", (double)(100*errornorm)));
18340a2a1afSJoe Wallwork     PetscCheck(errornorm < tol, comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric average test failed");
18400d473e3SJoe Wallwork 
18500d473e3SJoe Wallwork     /* Test metric intersection */
18661a4b293SJoe Wallwork     PetscCall(DMPlexMetricDeterminantCreate(dm, 0, &determinant, &dmDet));
18761a4b293SJoe Wallwork     if (!isotropic) {
18861a4b293SJoe Wallwork       PetscCall(DMPlexMetricEnforceSPD(dm, metrics[0], PETSC_FALSE, PETSC_FALSE, metricComb, determinant));
18961a4b293SJoe Wallwork       PetscCall(VecCopy(metricComb, metrics[0]));
19061a4b293SJoe Wallwork       PetscCall(DMPlexMetricEnforceSPD(dm, metrics[1], PETSC_FALSE, PETSC_FALSE, metricComb, determinant));
19161a4b293SJoe Wallwork       PetscCall(VecCopy(metricComb, metrics[1]));
19261a4b293SJoe Wallwork     }
1936f71502aSJoe Wallwork     PetscCall(DMPlexMetricIntersection(dm, 2, metrics, metricComb));
1946f71502aSJoe Wallwork     PetscCall(VecAXPY(metricComb, -1, metric2));
1959566063dSJacob Faibussowitsch     PetscCall(VecNorm(metricComb, NORM_2, &errornorm));
19600d473e3SJoe Wallwork     errornorm /= norm;
19763a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Metric intersection L2 error: %.4f%%\n", (double)(100*errornorm)));
19840a2a1afSJoe Wallwork     PetscCheck(errornorm < tol, comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric intersection test failed");
1999566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&metric2));
2009566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&metricComb));
20100d473e3SJoe Wallwork 
20200d473e3SJoe Wallwork     /* Test metric SPD enforcement */
20340a2a1afSJoe Wallwork     PetscCall(DMPlexMetricEnforceSPD(dm, metric, PETSC_TRUE, PETSC_TRUE, metric1, determinant));
204e5697243SJoe Wallwork     if (isotropic) {
205f49e87b0SJoe Wallwork       Vec err;
206f49e87b0SJoe Wallwork 
2079566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(determinant, &err));
2089566063dSJacob Faibussowitsch       PetscCall(VecSet(err, 1.0));
2099566063dSJacob Faibussowitsch       PetscCall(VecNorm(err, NORM_2, &norm));
2109566063dSJacob Faibussowitsch       PetscCall(VecAXPY(err, -1, determinant));
2119566063dSJacob Faibussowitsch       PetscCall(VecNorm(err, NORM_2, &errornorm));
2129566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&err));
213f49e87b0SJoe Wallwork       errornorm /= norm;
21463a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Metric determinant L2 error: %.4f%%\n", (double)(100*errornorm)));
21540a2a1afSJoe Wallwork       PetscCheck(errornorm < tol, comm, PETSC_ERR_ARG_OUTOFRANGE, "Determinant is not unit");
2169566063dSJacob Faibussowitsch       PetscCall(VecAXPY(metric1, -1, metric));
2179566063dSJacob Faibussowitsch       PetscCall(VecNorm(metric1, NORM_2, &errornorm));
21800d473e3SJoe Wallwork       errornorm /= norm;
21963a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Metric SPD enforcement L2 error: %.4f%%\n", (double)(100*errornorm)));
22040a2a1afSJoe Wallwork       PetscCheck(errornorm < tol, comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric SPD enforcement test failed");
22100d473e3SJoe Wallwork     }
22200d473e3SJoe Wallwork 
22300d473e3SJoe Wallwork     /* Test metric normalization */
22440a2a1afSJoe Wallwork     PetscCall(DMPlexMetricNormalize(dm, metric, PETSC_TRUE, PETSC_TRUE, metric1, determinant));
225e5697243SJoe Wallwork     if (isotropic) {
226e5697243SJoe Wallwork       PetscReal target;
227e5697243SJoe Wallwork 
2289566063dSJacob Faibussowitsch       PetscCall(DMPlexMetricGetTargetComplexity(dm, &target));
229e5697243SJoe Wallwork       scaling = PetscPowReal(target, 2.0/dim);
230d8b80fabSJoe Wallwork       if (uniform) {
2319566063dSJacob Faibussowitsch         PetscCall(DMPlexMetricCreateUniform(dm, 0, scaling, &metric2));
232d8b80fabSJoe Wallwork       } else {
233d8b80fabSJoe Wallwork         DM  dmIndi;
234d8b80fabSJoe Wallwork         Vec indicator;
235d8b80fabSJoe Wallwork 
2369566063dSJacob Faibussowitsch         PetscCall(CreateIndicator(dm, &indicator, &dmIndi));
2379566063dSJacob Faibussowitsch         PetscCall(VecSet(indicator, scaling));
2389566063dSJacob Faibussowitsch         PetscCall(DMPlexMetricCreateIsotropic(dm, 0, indicator, &metric2));
2399566063dSJacob Faibussowitsch         PetscCall(DMDestroy(&dmIndi));
2409566063dSJacob Faibussowitsch         PetscCall(VecDestroy(&indicator));
241d8b80fabSJoe Wallwork       }
2429566063dSJacob Faibussowitsch       PetscCall(VecAXPY(metric2, -1, metric1));
2439566063dSJacob Faibussowitsch       PetscCall(VecNorm(metric2, NORM_2, &errornorm));
24400d473e3SJoe Wallwork       errornorm /= norm;
24563a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Metric normalization L2 error: %.4f%%\n", (double)(100*errornorm)));
24640a2a1afSJoe Wallwork       PetscCheck(errornorm < tol, comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric normalization test failed");
24700d473e3SJoe Wallwork     }
24840a2a1afSJoe Wallwork     PetscCall(VecDestroy(&determinant));
24940a2a1afSJoe Wallwork     PetscCall(DMDestroy(&dmDet));
2509566063dSJacob Faibussowitsch     PetscCall(VecCopy(metric1, metric));
2519566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&metric2));
2529566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&metric1));
25300d473e3SJoe Wallwork   }
25400d473e3SJoe Wallwork 
25500d473e3SJoe Wallwork   /* Adapt the mesh */
2569566063dSJacob Faibussowitsch   PetscCall(DMAdaptMetric(dm, metric, bdLabel, rgLabel, &dmAdapt));
2579566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
2589566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) dmAdapt, "DM_adapted"));
2599566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&metric));
2609566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dmAdapt, NULL, "-adapted_mesh_view"));
26100d473e3SJoe Wallwork 
26241473654SJoe Wallwork   /* Test tag preservation */
26341473654SJoe Wallwork   if (!noTagging) {
26441473654SJoe Wallwork     PetscBool hasTag;
26541473654SJoe Wallwork     PetscInt  size;
26641473654SJoe Wallwork 
2679566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(dmAdapt, "Face Sets", &bdLabel));
2689566063dSJacob Faibussowitsch     PetscCall(DMLabelHasStratum(bdLabel, 1, &hasTag));
26928b400f6SJacob Faibussowitsch     PetscCheck(hasTag, comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have face tag 1");
2709566063dSJacob Faibussowitsch     PetscCall(DMLabelHasStratum(bdLabel, 2, &hasTag));
27128b400f6SJacob Faibussowitsch     PetscCheck(hasTag, comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have face tag 2");
2729566063dSJacob Faibussowitsch     PetscCall(DMLabelGetNumValues(bdLabel, &size));
27363a3b9bcSJacob Faibussowitsch     PetscCheck(size == 2, comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh has the wrong number of face tags (got %" PetscInt_FMT ", expected 2)", size);
27441473654SJoe Wallwork 
2759566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(dmAdapt, "Cell Sets", &rgLabel));
2769566063dSJacob Faibussowitsch     PetscCall(DMLabelHasStratum(rgLabel, 3, &hasTag));
27728b400f6SJacob Faibussowitsch     PetscCheck(hasTag, comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have cell tag 3");
2789566063dSJacob Faibussowitsch     PetscCall(DMLabelHasStratum(rgLabel, 4, &hasTag));
27928b400f6SJacob Faibussowitsch     PetscCheck(hasTag, comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have cell tag 4");
2809566063dSJacob Faibussowitsch     PetscCall(DMLabelGetNumValues(rgLabel, &size));
28163a3b9bcSJacob Faibussowitsch     PetscCheck(size == 2, comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh has the wrong number of cell tags (got %" PetscInt_FMT ", expected 2)", size);
28241473654SJoe Wallwork   }
28341473654SJoe Wallwork 
28400d473e3SJoe Wallwork   /* Clean up */
2859566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmAdapt));
2869566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
28700d473e3SJoe Wallwork   return 0;
28800d473e3SJoe Wallwork }
28900d473e3SJoe Wallwork 
29000d473e3SJoe Wallwork /*TEST
29100d473e3SJoe Wallwork 
29280f7b1e7SJoe Wallwork   testset:
29380f7b1e7SJoe Wallwork     requires: pragmatic
294e600fa54SMatthew G. Knepley     args: -dm_plex_box_faces 4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor pragmatic -noTagging
29580f7b1e7SJoe Wallwork 
296dd0671eeSJoe Wallwork     test:
297dd0671eeSJoe Wallwork       suffix: uniform_2d_pragmatic
298d8b80fabSJoe Wallwork       args: -dm_plex_metric_uniform
29900d473e3SJoe Wallwork     test:
300dd0671eeSJoe Wallwork       suffix: iso_2d_pragmatic
301d8b80fabSJoe Wallwork       args: -dm_plex_metric_isotropic
30200d473e3SJoe Wallwork     test:
303dd0671eeSJoe Wallwork       suffix: hessian_2d_pragmatic
30480f7b1e7SJoe Wallwork 
30580f7b1e7SJoe Wallwork   testset:
30680f7b1e7SJoe Wallwork     requires: pragmatic tetgen
307e600fa54SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_box_faces 4,4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor pragmatic -noTagging
30880f7b1e7SJoe Wallwork 
30980f7b1e7SJoe Wallwork     test:
31080f7b1e7SJoe Wallwork       suffix: uniform_3d_pragmatic
311d8b80fabSJoe Wallwork       args: -dm_plex_metric_uniform -noTagging
31280f7b1e7SJoe Wallwork     test:
31380f7b1e7SJoe Wallwork       suffix: iso_3d_pragmatic
314d8b80fabSJoe Wallwork       args: -dm_plex_metric_isotropic -noTagging
31500d473e3SJoe Wallwork     test:
316dd0671eeSJoe Wallwork       suffix: hessian_3d_pragmatic
31780f7b1e7SJoe Wallwork 
31880f7b1e7SJoe Wallwork   testset:
31980f7b1e7SJoe Wallwork     requires: mmg
320e600fa54SMatthew G. Knepley     args: -dm_plex_box_faces 4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor mmg
32180f7b1e7SJoe Wallwork 
32200d473e3SJoe Wallwork     test:
323dd0671eeSJoe Wallwork       suffix: uniform_2d_mmg
324d8b80fabSJoe Wallwork       args: -dm_plex_metric_uniform
325dd0671eeSJoe Wallwork     test:
326dd0671eeSJoe Wallwork       suffix: iso_2d_mmg
327d8b80fabSJoe Wallwork       args: -dm_plex_metric_isotropic
328dd0671eeSJoe Wallwork     test:
329dd0671eeSJoe Wallwork       suffix: hessian_2d_mmg
33080f7b1e7SJoe Wallwork 
33180f7b1e7SJoe Wallwork   testset:
33280f7b1e7SJoe Wallwork     requires: mmg tetgen
333e600fa54SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_box_faces 4,4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor mmg
33480f7b1e7SJoe Wallwork 
33580f7b1e7SJoe Wallwork     test:
33680f7b1e7SJoe Wallwork       suffix: uniform_3d_mmg
337d8b80fabSJoe Wallwork       args: -dm_plex_metric_uniform
33880f7b1e7SJoe Wallwork     test:
33980f7b1e7SJoe Wallwork       suffix: iso_3d_mmg
340d8b80fabSJoe Wallwork       args: -dm_plex_metric_isotropic
341dd0671eeSJoe Wallwork     test:
342dd0671eeSJoe Wallwork       suffix: hessian_3d_mmg
34380f7b1e7SJoe Wallwork 
34480f7b1e7SJoe Wallwork   testset:
34580f7b1e7SJoe Wallwork     requires: parmmg tetgen
34680f7b1e7SJoe Wallwork     nsize: 2
347e600fa54SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_box_faces 4,4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor parmmg
34880f7b1e7SJoe Wallwork 
349dd0671eeSJoe Wallwork     test:
350dd0671eeSJoe Wallwork       suffix: uniform_3d_parmmg
351d8b80fabSJoe Wallwork       args: -dm_plex_metric_uniform
352dd0671eeSJoe Wallwork     test:
353dd0671eeSJoe Wallwork       suffix: iso_3d_parmmg
354d8b80fabSJoe Wallwork       args: -dm_plex_metric_isotropic
355dd0671eeSJoe Wallwork     test:
356dd0671eeSJoe Wallwork       suffix: hessian_3d_parmmg
35700d473e3SJoe Wallwork 
35800d473e3SJoe Wallwork TEST*/
359