100d473e3SJoe Wallwork static char help[] = "Test metric utils in the uniform, isotropic case.\n\n"; 200d473e3SJoe Wallwork 300d473e3SJoe Wallwork #include <petscdmplex.h> 400d473e3SJoe Wallwork 5*9371c9d4SSatish Balay static PetscErrorCode bowl(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) { 600d473e3SJoe Wallwork PetscInt d; 700d473e3SJoe Wallwork 800d473e3SJoe Wallwork *u = 0.0; 900d473e3SJoe Wallwork for (d = 0; d < dim; d++) *u += 0.5 * (x[d] - 0.5) * (x[d] - 0.5); 1000d473e3SJoe Wallwork 1100d473e3SJoe Wallwork return 0; 1200d473e3SJoe Wallwork } 1300d473e3SJoe Wallwork 14*9371c9d4SSatish Balay static PetscErrorCode CreateIndicator(DM dm, Vec *indicator, DM *dmIndi) { 15d8b80fabSJoe Wallwork MPI_Comm comm; 16d8b80fabSJoe Wallwork PetscFE fe; 17d8b80fabSJoe Wallwork PetscInt dim; 18d8b80fabSJoe Wallwork 19d8b80fabSJoe Wallwork PetscFunctionBeginUser; 209566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 219566063dSJacob Faibussowitsch PetscCall(DMClone(dm, dmIndi)); 229566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 239566063dSJacob Faibussowitsch PetscCall(PetscFECreateLagrange(comm, dim, 1, PETSC_TRUE, 1, PETSC_DETERMINE, &fe)); 249566063dSJacob Faibussowitsch PetscCall(DMSetField(*dmIndi, 0, NULL, (PetscObject)fe)); 259566063dSJacob Faibussowitsch PetscCall(DMCreateDS(*dmIndi)); 269566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 279566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(*dmIndi, indicator)); 28d8b80fabSJoe Wallwork PetscFunctionReturn(0); 29d8b80fabSJoe Wallwork } 30d8b80fabSJoe Wallwork 3100d473e3SJoe Wallwork int main(int argc, char **argv) { 32e600fa54SMatthew G. Knepley DM dm, dmAdapt; 3321b3e102SJoe Wallwork DMLabel bdLabel = NULL, rgLabel = NULL; 3400d473e3SJoe Wallwork MPI_Comm comm; 3541473654SJoe Wallwork PetscBool uniform = PETSC_FALSE, isotropic = PETSC_FALSE, noTagging = PETSC_FALSE; 36e600fa54SMatthew G. Knepley PetscInt dim; 3700d473e3SJoe Wallwork PetscReal scaling = 1.0; 3800d473e3SJoe Wallwork Vec metric; 3900d473e3SJoe Wallwork 4000d473e3SJoe Wallwork /* Set up */ 41327415f7SBarry Smith PetscFunctionBeginUser; 429566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 4300d473e3SJoe Wallwork comm = PETSC_COMM_WORLD; 44d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Mesh adaptation options", "DMPLEX"); 459566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-noTagging", "Should tag preservation testing be turned off?", "ex60.c", noTagging, &noTagging, NULL)); 46d0609cedSBarry Smith PetscOptionsEnd(); 4700d473e3SJoe Wallwork 4800d473e3SJoe Wallwork /* Create box mesh */ 499566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, &dm)); 509566063dSJacob Faibussowitsch PetscCall(DMSetType(dm, DMPLEX)); 519566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 529566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)dm, "DM_init")); 539566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-initial_mesh_view")); 549566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 5500d473e3SJoe Wallwork 5641473654SJoe Wallwork /* Set tags to be preserved */ 5741473654SJoe Wallwork if (!noTagging) { 5841473654SJoe Wallwork DM cdm; 5941473654SJoe Wallwork PetscInt cStart, cEnd, c, fStart, fEnd, f, vStart, vEnd; 6041473654SJoe Wallwork const PetscScalar *coords; 6141473654SJoe Wallwork Vec coordinates; 6241473654SJoe Wallwork 6341473654SJoe Wallwork /* Cell tags */ 648fb5bd83SMatthew G. Knepley PetscCall(DMGetCoordinatesLocalSetUp(dm)); 659566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "Cell Sets")); 669566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Cell Sets", &rgLabel)); 679566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 6841473654SJoe Wallwork for (c = cStart; c < cEnd; ++c) { 6941473654SJoe Wallwork PetscReal centroid[3], volume, x; 7041473654SJoe Wallwork 719566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dm, c, &volume, centroid, NULL)); 7241473654SJoe Wallwork x = centroid[0]; 739566063dSJacob Faibussowitsch if (x < 0.5) PetscCall(DMLabelSetValue(rgLabel, c, 3)); 749566063dSJacob Faibussowitsch else PetscCall(DMLabelSetValue(rgLabel, c, 4)); 7541473654SJoe Wallwork } 7641473654SJoe Wallwork 7741473654SJoe Wallwork /* Face tags */ 789566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "Face Sets")); 799566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Face Sets", &bdLabel)); 809566063dSJacob Faibussowitsch PetscCall(DMPlexMarkBoundaryFaces(dm, 1, bdLabel)); 819566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd)); 829566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 839566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 849566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 859566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coordinates, &coords)); 8641473654SJoe Wallwork for (f = fStart; f < fEnd; ++f) { 8741473654SJoe Wallwork PetscBool flg = PETSC_TRUE; 8841473654SJoe Wallwork PetscInt *closure = NULL, closureSize, cl; 8941473654SJoe Wallwork PetscReal eps = 1.0e-08; 9041473654SJoe Wallwork 919566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure)); 9241473654SJoe Wallwork for (cl = 0; cl < closureSize * 2; cl += 2) { 9341473654SJoe Wallwork PetscInt off = closure[cl]; 9441473654SJoe Wallwork PetscReal *x; 9541473654SJoe Wallwork 9641473654SJoe Wallwork if ((off < vStart) || (off >= vEnd)) continue; 979566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRead(cdm, off, coords, &x)); 9841473654SJoe Wallwork if ((x[0] < 0.5 - eps) || (x[0] > 0.5 + eps)) flg = PETSC_FALSE; 9941473654SJoe Wallwork } 1009566063dSJacob Faibussowitsch if (flg) PetscCall(DMLabelSetValue(bdLabel, f, 2)); 1019566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure)); 10241473654SJoe Wallwork } 1039566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coordinates, &coords)); 10441473654SJoe Wallwork } 10541473654SJoe Wallwork 10600d473e3SJoe Wallwork /* Construct metric */ 1079566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 1089566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIsUniform(dm, &uniform)); 1099566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic)); 11000d473e3SJoe Wallwork if (uniform) { 1119566063dSJacob Faibussowitsch PetscCall(DMPlexMetricCreateUniform(dm, 0, scaling, &metric)); 112*9371c9d4SSatish Balay } else { 11300d473e3SJoe Wallwork DM dmIndi; 11400d473e3SJoe Wallwork Vec indicator; 11500d473e3SJoe Wallwork 11600d473e3SJoe Wallwork /* Construct "error indicator" */ 1179566063dSJacob Faibussowitsch PetscCall(CreateIndicator(dm, &indicator, &dmIndi)); 118e5697243SJoe Wallwork if (isotropic) { 11900d473e3SJoe Wallwork /* Isotropic case: just specify unity */ 1209566063dSJacob Faibussowitsch PetscCall(VecSet(indicator, scaling)); 1219566063dSJacob Faibussowitsch PetscCall(DMPlexMetricCreateIsotropic(dm, 0, indicator, &metric)); 12200d473e3SJoe Wallwork 12300d473e3SJoe Wallwork } else { 124d8b80fabSJoe Wallwork PetscFE fe; 12500d473e3SJoe Wallwork 12600d473e3SJoe Wallwork /* 'Anisotropic' case: approximate the identity by recovering the Hessian of a parabola */ 12700d473e3SJoe Wallwork DM dmGrad; 12800d473e3SJoe Wallwork PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {bowl}; 12900d473e3SJoe Wallwork Vec gradient; 13000d473e3SJoe Wallwork 13100d473e3SJoe Wallwork /* Project the parabola into P1 space */ 1329566063dSJacob Faibussowitsch PetscCall(DMProjectFunctionLocal(dmIndi, 0.0, funcs, NULL, INSERT_ALL_VALUES, indicator)); 13300d473e3SJoe Wallwork 13400d473e3SJoe Wallwork /* Approximate the gradient */ 1359566063dSJacob Faibussowitsch PetscCall(DMClone(dmIndi, &dmGrad)); 1369566063dSJacob Faibussowitsch PetscCall(PetscFECreateLagrange(comm, dim, dim, PETSC_TRUE, 1, PETSC_DETERMINE, &fe)); 1379566063dSJacob Faibussowitsch PetscCall(DMSetField(dmGrad, 0, NULL, (PetscObject)fe)); 1389566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dmGrad)); 1399566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 1409566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dmGrad, &gradient)); 1419566063dSJacob Faibussowitsch PetscCall(DMPlexComputeGradientClementInterpolant(dmIndi, indicator, gradient)); 1429566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(gradient, NULL, "-adapt_gradient_view")); 14300d473e3SJoe Wallwork 14400d473e3SJoe Wallwork /* Approximate the Hessian */ 1459566063dSJacob Faibussowitsch PetscCall(DMPlexMetricCreate(dm, 0, &metric)); 1469566063dSJacob Faibussowitsch PetscCall(DMPlexComputeGradientClementInterpolant(dmGrad, gradient, metric)); 1479566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(metric, NULL, "-adapt_hessian_view")); 1489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&gradient)); 1499566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmGrad)); 15000d473e3SJoe Wallwork } 1519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&indicator)); 1529566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmIndi)); 15300d473e3SJoe Wallwork } 15400d473e3SJoe Wallwork 15500d473e3SJoe Wallwork /* Test metric routines */ 15600d473e3SJoe Wallwork { 157f49e87b0SJoe Wallwork DM dmDet; 15800d473e3SJoe Wallwork PetscReal errornorm, norm, tol = 1.0e-10, weights[2] = {0.8, 0.2}; 159f49e87b0SJoe Wallwork Vec metric1, metric2, metricComb, determinant; 16000d473e3SJoe Wallwork Vec metrics[2]; 16100d473e3SJoe Wallwork 1629566063dSJacob Faibussowitsch PetscCall(VecDuplicate(metric, &metric1)); 1639566063dSJacob Faibussowitsch PetscCall(VecSet(metric1, 0)); 1649566063dSJacob Faibussowitsch PetscCall(VecAXPY(metric1, 0.625, metric)); 1659566063dSJacob Faibussowitsch PetscCall(VecDuplicate(metric, &metric2)); 1669566063dSJacob Faibussowitsch PetscCall(VecSet(metric2, 0)); 1679566063dSJacob Faibussowitsch PetscCall(VecAXPY(metric2, 2.5, metric)); 16800d473e3SJoe Wallwork metrics[0] = metric1; 16900d473e3SJoe Wallwork metrics[1] = metric2; 17000d473e3SJoe Wallwork 17100d473e3SJoe Wallwork /* Test metric average */ 17240a2a1afSJoe Wallwork PetscCall(DMPlexMetricCreate(dm, 0, &metricComb)); 17340a2a1afSJoe Wallwork PetscCall(DMPlexMetricAverage(dm, 2, weights, metrics, metricComb)); 1749566063dSJacob Faibussowitsch PetscCall(VecAXPY(metricComb, -1, metric)); 1759566063dSJacob Faibussowitsch PetscCall(VecNorm(metric, NORM_2, &norm)); 1769566063dSJacob Faibussowitsch PetscCall(VecNorm(metricComb, NORM_2, &errornorm)); 17700d473e3SJoe Wallwork errornorm /= norm; 17863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Metric average L2 error: %.4f%%\n", (double)(100 * errornorm))); 17940a2a1afSJoe Wallwork PetscCheck(errornorm < tol, comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric average test failed"); 18000d473e3SJoe Wallwork 18100d473e3SJoe Wallwork /* Test metric intersection */ 18261a4b293SJoe Wallwork PetscCall(DMPlexMetricDeterminantCreate(dm, 0, &determinant, &dmDet)); 18361a4b293SJoe Wallwork if (!isotropic) { 18461a4b293SJoe Wallwork PetscCall(DMPlexMetricEnforceSPD(dm, metrics[0], PETSC_FALSE, PETSC_FALSE, metricComb, determinant)); 18561a4b293SJoe Wallwork PetscCall(VecCopy(metricComb, metrics[0])); 18661a4b293SJoe Wallwork PetscCall(DMPlexMetricEnforceSPD(dm, metrics[1], PETSC_FALSE, PETSC_FALSE, metricComb, determinant)); 18761a4b293SJoe Wallwork PetscCall(VecCopy(metricComb, metrics[1])); 18861a4b293SJoe Wallwork } 1896f71502aSJoe Wallwork PetscCall(DMPlexMetricIntersection(dm, 2, metrics, metricComb)); 1906f71502aSJoe Wallwork PetscCall(VecAXPY(metricComb, -1, metric2)); 1919566063dSJacob Faibussowitsch PetscCall(VecNorm(metricComb, NORM_2, &errornorm)); 19200d473e3SJoe Wallwork errornorm /= norm; 19363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Metric intersection L2 error: %.4f%%\n", (double)(100 * errornorm))); 19440a2a1afSJoe Wallwork PetscCheck(errornorm < tol, comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric intersection test failed"); 1959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&metric2)); 1969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&metricComb)); 19700d473e3SJoe Wallwork 19800d473e3SJoe Wallwork /* Test metric SPD enforcement */ 19940a2a1afSJoe Wallwork PetscCall(DMPlexMetricEnforceSPD(dm, metric, PETSC_TRUE, PETSC_TRUE, metric1, determinant)); 200e5697243SJoe Wallwork if (isotropic) { 201f49e87b0SJoe Wallwork Vec err; 202f49e87b0SJoe Wallwork 2039566063dSJacob Faibussowitsch PetscCall(VecDuplicate(determinant, &err)); 2049566063dSJacob Faibussowitsch PetscCall(VecSet(err, 1.0)); 2059566063dSJacob Faibussowitsch PetscCall(VecNorm(err, NORM_2, &norm)); 2069566063dSJacob Faibussowitsch PetscCall(VecAXPY(err, -1, determinant)); 2079566063dSJacob Faibussowitsch PetscCall(VecNorm(err, NORM_2, &errornorm)); 2089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&err)); 209f49e87b0SJoe Wallwork errornorm /= norm; 21063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Metric determinant L2 error: %.4f%%\n", (double)(100 * errornorm))); 21140a2a1afSJoe Wallwork PetscCheck(errornorm < tol, comm, PETSC_ERR_ARG_OUTOFRANGE, "Determinant is not unit"); 2129566063dSJacob Faibussowitsch PetscCall(VecAXPY(metric1, -1, metric)); 2139566063dSJacob Faibussowitsch PetscCall(VecNorm(metric1, NORM_2, &errornorm)); 21400d473e3SJoe Wallwork errornorm /= norm; 21563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Metric SPD enforcement L2 error: %.4f%%\n", (double)(100 * errornorm))); 21640a2a1afSJoe Wallwork PetscCheck(errornorm < tol, comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric SPD enforcement test failed"); 21700d473e3SJoe Wallwork } 21800d473e3SJoe Wallwork 21900d473e3SJoe Wallwork /* Test metric normalization */ 22040a2a1afSJoe Wallwork PetscCall(DMPlexMetricNormalize(dm, metric, PETSC_TRUE, PETSC_TRUE, metric1, determinant)); 221e5697243SJoe Wallwork if (isotropic) { 222e5697243SJoe Wallwork PetscReal target; 223e5697243SJoe Wallwork 2249566063dSJacob Faibussowitsch PetscCall(DMPlexMetricGetTargetComplexity(dm, &target)); 225e5697243SJoe Wallwork scaling = PetscPowReal(target, 2.0 / dim); 226d8b80fabSJoe Wallwork if (uniform) { 2279566063dSJacob Faibussowitsch PetscCall(DMPlexMetricCreateUniform(dm, 0, scaling, &metric2)); 228d8b80fabSJoe Wallwork } else { 229d8b80fabSJoe Wallwork DM dmIndi; 230d8b80fabSJoe Wallwork Vec indicator; 231d8b80fabSJoe Wallwork 2329566063dSJacob Faibussowitsch PetscCall(CreateIndicator(dm, &indicator, &dmIndi)); 2339566063dSJacob Faibussowitsch PetscCall(VecSet(indicator, scaling)); 2349566063dSJacob Faibussowitsch PetscCall(DMPlexMetricCreateIsotropic(dm, 0, indicator, &metric2)); 2359566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmIndi)); 2369566063dSJacob Faibussowitsch PetscCall(VecDestroy(&indicator)); 237d8b80fabSJoe Wallwork } 2389566063dSJacob Faibussowitsch PetscCall(VecAXPY(metric2, -1, metric1)); 2399566063dSJacob Faibussowitsch PetscCall(VecNorm(metric2, NORM_2, &errornorm)); 24000d473e3SJoe Wallwork errornorm /= norm; 24163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Metric normalization L2 error: %.4f%%\n", (double)(100 * errornorm))); 24240a2a1afSJoe Wallwork PetscCheck(errornorm < tol, comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric normalization test failed"); 24300d473e3SJoe Wallwork } 24440a2a1afSJoe Wallwork PetscCall(VecDestroy(&determinant)); 24540a2a1afSJoe Wallwork PetscCall(DMDestroy(&dmDet)); 2469566063dSJacob Faibussowitsch PetscCall(VecCopy(metric1, metric)); 2479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&metric2)); 2489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&metric1)); 24900d473e3SJoe Wallwork } 25000d473e3SJoe Wallwork 25100d473e3SJoe Wallwork /* Adapt the mesh */ 2529566063dSJacob Faibussowitsch PetscCall(DMAdaptMetric(dm, metric, bdLabel, rgLabel, &dmAdapt)); 2539566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 2549566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)dmAdapt, "DM_adapted")); 2559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&metric)); 2569566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dmAdapt, NULL, "-adapted_mesh_view")); 25700d473e3SJoe Wallwork 25841473654SJoe Wallwork /* Test tag preservation */ 25941473654SJoe Wallwork if (!noTagging) { 26041473654SJoe Wallwork PetscBool hasTag; 26141473654SJoe Wallwork PetscInt size; 26241473654SJoe Wallwork 2639566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dmAdapt, "Face Sets", &bdLabel)); 2649566063dSJacob Faibussowitsch PetscCall(DMLabelHasStratum(bdLabel, 1, &hasTag)); 26528b400f6SJacob Faibussowitsch PetscCheck(hasTag, comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have face tag 1"); 2669566063dSJacob Faibussowitsch PetscCall(DMLabelHasStratum(bdLabel, 2, &hasTag)); 26728b400f6SJacob Faibussowitsch PetscCheck(hasTag, comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have face tag 2"); 2689566063dSJacob Faibussowitsch PetscCall(DMLabelGetNumValues(bdLabel, &size)); 26963a3b9bcSJacob Faibussowitsch PetscCheck(size == 2, comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh has the wrong number of face tags (got %" PetscInt_FMT ", expected 2)", size); 27041473654SJoe Wallwork 2719566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dmAdapt, "Cell Sets", &rgLabel)); 2729566063dSJacob Faibussowitsch PetscCall(DMLabelHasStratum(rgLabel, 3, &hasTag)); 27328b400f6SJacob Faibussowitsch PetscCheck(hasTag, comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have cell tag 3"); 2749566063dSJacob Faibussowitsch PetscCall(DMLabelHasStratum(rgLabel, 4, &hasTag)); 27528b400f6SJacob Faibussowitsch PetscCheck(hasTag, comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have cell tag 4"); 2769566063dSJacob Faibussowitsch PetscCall(DMLabelGetNumValues(rgLabel, &size)); 27763a3b9bcSJacob Faibussowitsch PetscCheck(size == 2, comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh has the wrong number of cell tags (got %" PetscInt_FMT ", expected 2)", size); 27841473654SJoe Wallwork } 27941473654SJoe Wallwork 28000d473e3SJoe Wallwork /* Clean up */ 2819566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmAdapt)); 2829566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 28300d473e3SJoe Wallwork return 0; 28400d473e3SJoe Wallwork } 28500d473e3SJoe Wallwork 28600d473e3SJoe Wallwork /*TEST 28700d473e3SJoe Wallwork 28880f7b1e7SJoe Wallwork testset: 28980f7b1e7SJoe Wallwork requires: pragmatic 290e600fa54SMatthew G. Knepley args: -dm_plex_box_faces 4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor pragmatic -noTagging 29180f7b1e7SJoe Wallwork 292dd0671eeSJoe Wallwork test: 293dd0671eeSJoe Wallwork suffix: uniform_2d_pragmatic 294d8b80fabSJoe Wallwork args: -dm_plex_metric_uniform 29500d473e3SJoe Wallwork test: 296dd0671eeSJoe Wallwork suffix: iso_2d_pragmatic 297d8b80fabSJoe Wallwork args: -dm_plex_metric_isotropic 29800d473e3SJoe Wallwork test: 299dd0671eeSJoe Wallwork suffix: hessian_2d_pragmatic 30080f7b1e7SJoe Wallwork 30180f7b1e7SJoe Wallwork testset: 30280f7b1e7SJoe Wallwork requires: pragmatic tetgen 303e600fa54SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 4,4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor pragmatic -noTagging 30480f7b1e7SJoe Wallwork 30580f7b1e7SJoe Wallwork test: 30680f7b1e7SJoe Wallwork suffix: uniform_3d_pragmatic 307d8b80fabSJoe Wallwork args: -dm_plex_metric_uniform -noTagging 30880f7b1e7SJoe Wallwork test: 30980f7b1e7SJoe Wallwork suffix: iso_3d_pragmatic 310d8b80fabSJoe Wallwork args: -dm_plex_metric_isotropic -noTagging 31100d473e3SJoe Wallwork test: 312dd0671eeSJoe Wallwork suffix: hessian_3d_pragmatic 31380f7b1e7SJoe Wallwork 31480f7b1e7SJoe Wallwork testset: 31580f7b1e7SJoe Wallwork requires: mmg 316e600fa54SMatthew G. Knepley args: -dm_plex_box_faces 4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor mmg 31780f7b1e7SJoe Wallwork 31800d473e3SJoe Wallwork test: 319dd0671eeSJoe Wallwork suffix: uniform_2d_mmg 320d8b80fabSJoe Wallwork args: -dm_plex_metric_uniform 321dd0671eeSJoe Wallwork test: 322dd0671eeSJoe Wallwork suffix: iso_2d_mmg 323d8b80fabSJoe Wallwork args: -dm_plex_metric_isotropic 324dd0671eeSJoe Wallwork test: 325dd0671eeSJoe Wallwork suffix: hessian_2d_mmg 32680f7b1e7SJoe Wallwork 32780f7b1e7SJoe Wallwork testset: 32880f7b1e7SJoe Wallwork requires: mmg tetgen 329e600fa54SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 4,4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor mmg 33080f7b1e7SJoe Wallwork 33180f7b1e7SJoe Wallwork test: 33280f7b1e7SJoe Wallwork suffix: uniform_3d_mmg 333d8b80fabSJoe Wallwork args: -dm_plex_metric_uniform 33480f7b1e7SJoe Wallwork test: 33580f7b1e7SJoe Wallwork suffix: iso_3d_mmg 336d8b80fabSJoe Wallwork args: -dm_plex_metric_isotropic 337dd0671eeSJoe Wallwork test: 338dd0671eeSJoe Wallwork suffix: hessian_3d_mmg 33980f7b1e7SJoe Wallwork 34080f7b1e7SJoe Wallwork testset: 34180f7b1e7SJoe Wallwork requires: parmmg tetgen 34280f7b1e7SJoe Wallwork nsize: 2 343e600fa54SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 4,4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor parmmg 34480f7b1e7SJoe Wallwork 345dd0671eeSJoe Wallwork test: 346dd0671eeSJoe Wallwork suffix: uniform_3d_parmmg 347d8b80fabSJoe Wallwork args: -dm_plex_metric_uniform 348dd0671eeSJoe Wallwork test: 349dd0671eeSJoe Wallwork suffix: iso_3d_parmmg 350d8b80fabSJoe Wallwork args: -dm_plex_metric_isotropic 351dd0671eeSJoe Wallwork test: 352dd0671eeSJoe Wallwork suffix: hessian_3d_parmmg 35300d473e3SJoe Wallwork 35400d473e3SJoe Wallwork TEST*/ 355