xref: /petsc/src/dm/impls/plex/tests/ex60.c (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
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;
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 */
449566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
4500d473e3SJoe Wallwork   comm = PETSC_COMM_WORLD;
469566063dSJacob Faibussowitsch   ierr = PetscOptionsBegin(comm, "", "Mesh adaptation options", "DMPLEX");PetscCall(ierr);
479566063dSJacob Faibussowitsch   PetscCall(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 */
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 */
669566063dSJacob Faibussowitsch     PetscCall(DMCreateLabel(dm, "Cell Sets"));
679566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(dm, "Cell Sets", &rgLabel));
689566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
6941473654SJoe Wallwork     for (c = cStart; c < cEnd; ++c) {
7041473654SJoe Wallwork       PetscReal centroid[3], volume, x;
7141473654SJoe Wallwork 
729566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFVM(dm, c, &volume, centroid, NULL));
7341473654SJoe Wallwork       x = centroid[0];
749566063dSJacob Faibussowitsch       if (x < 0.5) PetscCall(DMLabelSetValue(rgLabel, c, 3));
759566063dSJacob Faibussowitsch       else         PetscCall(DMLabelSetValue(rgLabel, c, 4));
7641473654SJoe Wallwork     }
7741473654SJoe Wallwork 
7841473654SJoe Wallwork     /* Face tags */
799566063dSJacob Faibussowitsch     PetscCall(DMCreateLabel(dm, "Face Sets"));
809566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(dm, "Face Sets", &bdLabel));
819566063dSJacob Faibussowitsch     PetscCall(DMPlexMarkBoundaryFaces(dm, 1, bdLabel));
829566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
839566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
849566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dm, &cdm));
859566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
869566063dSJacob Faibussowitsch     PetscCall(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 
929566063dSJacob Faibussowitsch       PetscCall(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;
989566063dSJacob Faibussowitsch         PetscCall(DMPlexPointLocalRead(cdm, off, coords, &x));
9941473654SJoe Wallwork         if ((x[0] < 0.5 - eps) || (x[0] > 0.5 + eps)) flg = PETSC_FALSE;
10041473654SJoe Wallwork       }
1019566063dSJacob Faibussowitsch       if (flg) PetscCall(DMLabelSetValue(bdLabel, f, 2));
1029566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
10341473654SJoe Wallwork     }
1049566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(coordinates, &coords));
10541473654SJoe Wallwork   }
10641473654SJoe Wallwork 
10700d473e3SJoe Wallwork   /* Construct metric */
1089566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricSetFromOptions(dm));
1099566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1109566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
11100d473e3SJoe Wallwork   if (uniform) {
1129566063dSJacob Faibussowitsch     PetscCall(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" */
1199566063dSJacob Faibussowitsch     PetscCall(CreateIndicator(dm, &indicator, &dmIndi));
120e5697243SJoe Wallwork     if (isotropic) {
12100d473e3SJoe Wallwork 
12200d473e3SJoe Wallwork       /* Isotropic case: just specify unity */
1239566063dSJacob Faibussowitsch       PetscCall(VecSet(indicator, scaling));
1249566063dSJacob Faibussowitsch       PetscCall(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 */
1359566063dSJacob Faibussowitsch       PetscCall(DMProjectFunctionLocal(dmIndi, 0.0, funcs, NULL, INSERT_ALL_VALUES, indicator));
13600d473e3SJoe Wallwork 
13700d473e3SJoe Wallwork       /* Approximate the gradient */
1389566063dSJacob Faibussowitsch       PetscCall(DMClone(dmIndi, &dmGrad));
1399566063dSJacob Faibussowitsch       PetscCall(PetscFECreateLagrange(comm, dim, dim, PETSC_TRUE, 1, PETSC_DETERMINE, &fe));
1409566063dSJacob Faibussowitsch       PetscCall(DMSetField(dmGrad, 0, NULL, (PetscObject)fe));
1419566063dSJacob Faibussowitsch       PetscCall(DMCreateDS(dmGrad));
1429566063dSJacob Faibussowitsch       PetscCall(PetscFEDestroy(&fe));
1439566063dSJacob Faibussowitsch       PetscCall(DMCreateLocalVector(dmGrad, &gradient));
1449566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeGradientClementInterpolant(dmIndi, indicator, gradient));
1459566063dSJacob Faibussowitsch       PetscCall(VecViewFromOptions(gradient, NULL, "-adapt_gradient_view"));
14600d473e3SJoe Wallwork 
14700d473e3SJoe Wallwork       /* Approximate the Hessian */
1489566063dSJacob Faibussowitsch       PetscCall(DMPlexMetricCreate(dm, 0, &metric));
1499566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeGradientClementInterpolant(dmGrad, gradient, metric));
1509566063dSJacob Faibussowitsch       PetscCall(VecViewFromOptions(metric, NULL, "-adapt_hessian_view"));
1519566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&gradient));
1529566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&dmGrad));
15300d473e3SJoe Wallwork     }
1549566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&indicator));
1559566063dSJacob Faibussowitsch     PetscCall(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 
1659566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(metric, &metric1));
1669566063dSJacob Faibussowitsch     PetscCall(VecSet(metric1, 0));
1679566063dSJacob Faibussowitsch     PetscCall(VecAXPY(metric1, 0.625, metric));
1689566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(metric, &metric2));
1699566063dSJacob Faibussowitsch     PetscCall(VecSet(metric2, 0));
1709566063dSJacob Faibussowitsch     PetscCall(VecAXPY(metric2, 2.5, metric));
17100d473e3SJoe Wallwork     metrics[0] = metric1;
17200d473e3SJoe Wallwork     metrics[1] = metric2;
17300d473e3SJoe Wallwork 
17400d473e3SJoe Wallwork     /* Test metric average */
1759566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricAverage(dm, 2, weights, metrics, &metricComb));
1769566063dSJacob Faibussowitsch     PetscCall(VecAXPY(metricComb, -1, metric));
1779566063dSJacob Faibussowitsch     PetscCall(VecNorm(metric, NORM_2, &norm));
1789566063dSJacob Faibussowitsch     PetscCall(VecNorm(metricComb, NORM_2, &errornorm));
17900d473e3SJoe Wallwork     errornorm /= norm;
1809566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Metric average L2 error: %.4f%%\n", 100*errornorm));
181*08401ef6SPierre Jolivet     PetscCheck(errornorm <= tol,comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric average test failed");
1829566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&metricComb));
18300d473e3SJoe Wallwork 
18400d473e3SJoe Wallwork     /* Test metric intersection */
185e5697243SJoe Wallwork     if (isotropic) {
1869566063dSJacob Faibussowitsch       PetscCall(DMPlexMetricIntersection(dm, 2, metrics, &metricComb));
1879566063dSJacob Faibussowitsch       PetscCall(VecAXPY(metricComb, -1, metric1));
1889566063dSJacob Faibussowitsch       PetscCall(VecNorm(metricComb, NORM_2, &errornorm));
18900d473e3SJoe Wallwork       errornorm /= norm;
1909566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Metric intersection L2 error: %.4f%%\n", 100*errornorm));
191*08401ef6SPierre Jolivet       PetscCheck(errornorm <= tol,comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric intersection test failed");
19200d473e3SJoe Wallwork     }
1939566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&metric1));
1949566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&metric2));
1959566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&metricComb));
19600d473e3SJoe Wallwork 
19700d473e3SJoe Wallwork     /* Test metric SPD enforcement */
1989566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricEnforceSPD(dm, metric, PETSC_TRUE, PETSC_TRUE, &metric1, &determinant));
199e5697243SJoe Wallwork     if (isotropic) {
200f49e87b0SJoe Wallwork       Vec err;
201f49e87b0SJoe Wallwork 
2029566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(determinant, &err));
2039566063dSJacob Faibussowitsch       PetscCall(VecSet(err, 1.0));
2049566063dSJacob Faibussowitsch       PetscCall(VecNorm(err, NORM_2, &norm));
2059566063dSJacob Faibussowitsch       PetscCall(VecAXPY(err, -1, determinant));
2069566063dSJacob Faibussowitsch       PetscCall(VecNorm(err, NORM_2, &errornorm));
2079566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&err));
208f49e87b0SJoe Wallwork       errornorm /= norm;
2099566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Metric determinant L2 error: %.4f%%\n", 100*errornorm));
210*08401ef6SPierre Jolivet       PetscCheck(errornorm <= tol,comm, PETSC_ERR_ARG_OUTOFRANGE, "Determinant is not unit");
2119566063dSJacob Faibussowitsch       PetscCall(VecAXPY(metric1, -1, metric));
2129566063dSJacob Faibussowitsch       PetscCall(VecNorm(metric1, NORM_2, &errornorm));
21300d473e3SJoe Wallwork       errornorm /= norm;
2149566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Metric SPD enforcement L2 error: %.4f%%\n", 100*errornorm));
215*08401ef6SPierre Jolivet       PetscCheck(errornorm <= tol,comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric SPD enforcement test failed");
21600d473e3SJoe Wallwork     }
2179566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&metric1));
2189566063dSJacob Faibussowitsch     PetscCall(VecGetDM(determinant, &dmDet));
2199566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&determinant));
2209566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dmDet));
22100d473e3SJoe Wallwork 
22200d473e3SJoe Wallwork     /* Test metric normalization */
2239566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricNormalize(dm, metric, PETSC_TRUE, PETSC_TRUE, &metric1));
224e5697243SJoe Wallwork     if (isotropic) {
225e5697243SJoe Wallwork       PetscReal target;
226e5697243SJoe Wallwork 
2279566063dSJacob Faibussowitsch       PetscCall(DMPlexMetricGetTargetComplexity(dm, &target));
228e5697243SJoe Wallwork       scaling = PetscPowReal(target, 2.0/dim);
229d8b80fabSJoe Wallwork       if (uniform) {
2309566063dSJacob Faibussowitsch         PetscCall(DMPlexMetricCreateUniform(dm, 0, scaling, &metric2));
231d8b80fabSJoe Wallwork       } else {
232d8b80fabSJoe Wallwork         DM  dmIndi;
233d8b80fabSJoe Wallwork         Vec indicator;
234d8b80fabSJoe Wallwork 
2359566063dSJacob Faibussowitsch         PetscCall(CreateIndicator(dm, &indicator, &dmIndi));
2369566063dSJacob Faibussowitsch         PetscCall(VecSet(indicator, scaling));
2379566063dSJacob Faibussowitsch         PetscCall(DMPlexMetricCreateIsotropic(dm, 0, indicator, &metric2));
2389566063dSJacob Faibussowitsch         PetscCall(DMDestroy(&dmIndi));
2399566063dSJacob Faibussowitsch         PetscCall(VecDestroy(&indicator));
240d8b80fabSJoe Wallwork       }
2419566063dSJacob Faibussowitsch       PetscCall(VecAXPY(metric2, -1, metric1));
2429566063dSJacob Faibussowitsch       PetscCall(VecNorm(metric2, NORM_2, &errornorm));
24300d473e3SJoe Wallwork       errornorm /= norm;
2449566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Metric normalization L2 error: %.4f%%\n", 100*errornorm));
245*08401ef6SPierre Jolivet       PetscCheck(errornorm <= tol,comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric normalization test failed");
24600d473e3SJoe Wallwork     }
2479566063dSJacob Faibussowitsch     PetscCall(VecCopy(metric1, metric));
2489566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&metric2));
2499566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&metric1));
25000d473e3SJoe Wallwork   }
25100d473e3SJoe Wallwork 
25200d473e3SJoe Wallwork   /* Adapt the mesh */
2539566063dSJacob Faibussowitsch   PetscCall(DMAdaptMetric(dm, metric, bdLabel, rgLabel, &dmAdapt));
2549566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
2559566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) dmAdapt, "DM_adapted"));
2569566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&metric));
2579566063dSJacob Faibussowitsch   PetscCall(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 
2649566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(dmAdapt, "Face Sets", &bdLabel));
2659566063dSJacob Faibussowitsch     PetscCall(DMLabelHasStratum(bdLabel, 1, &hasTag));
26628b400f6SJacob Faibussowitsch     PetscCheck(hasTag,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have face tag 1");
2679566063dSJacob Faibussowitsch     PetscCall(DMLabelHasStratum(bdLabel, 2, &hasTag));
26828b400f6SJacob Faibussowitsch     PetscCheck(hasTag,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have face tag 2");
2699566063dSJacob Faibussowitsch     PetscCall(DMLabelGetNumValues(bdLabel, &size));
270*08401ef6SPierre Jolivet     PetscCheck(size == 2,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh has the wrong number of face tags (got %d, expected 2)", size);
27141473654SJoe Wallwork 
2729566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(dmAdapt, "Cell Sets", &rgLabel));
2739566063dSJacob Faibussowitsch     PetscCall(DMLabelHasStratum(rgLabel, 3, &hasTag));
27428b400f6SJacob Faibussowitsch     PetscCheck(hasTag,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have cell tag 3");
2759566063dSJacob Faibussowitsch     PetscCall(DMLabelHasStratum(rgLabel, 4, &hasTag));
27628b400f6SJacob Faibussowitsch     PetscCheck(hasTag,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have cell tag 4");
2779566063dSJacob Faibussowitsch     PetscCall(DMLabelGetNumValues(rgLabel, &size));
278*08401ef6SPierre Jolivet     PetscCheck(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 */
2829566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmAdapt));
2839566063dSJacob Faibussowitsch   PetscCall(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