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 */ 439566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 4400d473e3SJoe Wallwork comm = PETSC_COMM_WORLD; 45d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Mesh adaptation options", "DMPLEX"); 469566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-noTagging", "Should tag preservation testing be turned off?", "ex60.c", noTagging, &noTagging, NULL)); 47d0609cedSBarry Smith PetscOptionsEnd(); 4800d473e3SJoe Wallwork 4900d473e3SJoe Wallwork /* Create box mesh */ 509566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, &dm)); 519566063dSJacob Faibussowitsch PetscCall(DMSetType(dm, DMPLEX)); 529566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 539566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) dm, "DM_init")); 549566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-initial_mesh_view")); 559566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 5600d473e3SJoe Wallwork 5741473654SJoe Wallwork /* Set tags to be preserved */ 5841473654SJoe Wallwork if (!noTagging) { 5941473654SJoe Wallwork DM cdm; 6041473654SJoe Wallwork PetscInt cStart, cEnd, c, fStart, fEnd, f, vStart, vEnd; 6141473654SJoe Wallwork const PetscScalar *coords; 6241473654SJoe Wallwork Vec coordinates; 6341473654SJoe Wallwork 6441473654SJoe Wallwork /* Cell tags */ 65*8fb5bd83SMatthew G. Knepley PetscCall(DMGetCoordinatesLocalSetUp(dm)); 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 */ 17540a2a1afSJoe Wallwork PetscCall(DMPlexMetricCreate(dm, 0, &metricComb)); 17640a2a1afSJoe Wallwork PetscCall(DMPlexMetricAverage(dm, 2, weights, metrics, metricComb)); 1779566063dSJacob Faibussowitsch PetscCall(VecAXPY(metricComb, -1, metric)); 1789566063dSJacob Faibussowitsch PetscCall(VecNorm(metric, NORM_2, &norm)); 1799566063dSJacob Faibussowitsch PetscCall(VecNorm(metricComb, NORM_2, &errornorm)); 18000d473e3SJoe Wallwork errornorm /= norm; 18163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Metric average L2 error: %.4f%%\n", (double)(100*errornorm))); 18240a2a1afSJoe Wallwork PetscCheck(errornorm < tol, comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric average test failed"); 18300d473e3SJoe Wallwork 18400d473e3SJoe Wallwork /* Test metric intersection */ 185e5697243SJoe Wallwork if (isotropic) { 1866f71502aSJoe Wallwork PetscCall(DMPlexMetricIntersection(dm, 2, metrics, metricComb)); 1876f71502aSJoe Wallwork PetscCall(VecAXPY(metricComb, -1, metric2)); 1889566063dSJacob Faibussowitsch PetscCall(VecNorm(metricComb, NORM_2, &errornorm)); 18900d473e3SJoe Wallwork errornorm /= norm; 19063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Metric intersection L2 error: %.4f%%\n", (double)(100*errornorm))); 19140a2a1afSJoe Wallwork PetscCheck(errornorm < tol, comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric intersection test failed"); 19200d473e3SJoe Wallwork } 1939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&metric2)); 1949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&metricComb)); 19500d473e3SJoe Wallwork 19600d473e3SJoe Wallwork /* Test metric SPD enforcement */ 19740a2a1afSJoe Wallwork PetscCall(DMPlexMetricDeterminantCreate(dm, 0, &determinant, &dmDet)); 19840a2a1afSJoe Wallwork 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; 20963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Metric determinant L2 error: %.4f%%\n", (double)(100*errornorm))); 21040a2a1afSJoe Wallwork 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; 21463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Metric SPD enforcement L2 error: %.4f%%\n", (double)(100*errornorm))); 21540a2a1afSJoe Wallwork PetscCheck(errornorm < tol, comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric SPD enforcement test failed"); 21600d473e3SJoe Wallwork } 21700d473e3SJoe Wallwork 21800d473e3SJoe Wallwork /* Test metric normalization */ 21940a2a1afSJoe Wallwork PetscCall(DMPlexMetricNormalize(dm, metric, PETSC_TRUE, PETSC_TRUE, metric1, determinant)); 220e5697243SJoe Wallwork if (isotropic) { 221e5697243SJoe Wallwork PetscReal target; 222e5697243SJoe Wallwork 2239566063dSJacob Faibussowitsch PetscCall(DMPlexMetricGetTargetComplexity(dm, &target)); 224e5697243SJoe Wallwork scaling = PetscPowReal(target, 2.0/dim); 225d8b80fabSJoe Wallwork if (uniform) { 2269566063dSJacob Faibussowitsch PetscCall(DMPlexMetricCreateUniform(dm, 0, scaling, &metric2)); 227d8b80fabSJoe Wallwork } else { 228d8b80fabSJoe Wallwork DM dmIndi; 229d8b80fabSJoe Wallwork Vec indicator; 230d8b80fabSJoe Wallwork 2319566063dSJacob Faibussowitsch PetscCall(CreateIndicator(dm, &indicator, &dmIndi)); 2329566063dSJacob Faibussowitsch PetscCall(VecSet(indicator, scaling)); 2339566063dSJacob Faibussowitsch PetscCall(DMPlexMetricCreateIsotropic(dm, 0, indicator, &metric2)); 2349566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmIndi)); 2359566063dSJacob Faibussowitsch PetscCall(VecDestroy(&indicator)); 236d8b80fabSJoe Wallwork } 2379566063dSJacob Faibussowitsch PetscCall(VecAXPY(metric2, -1, metric1)); 2389566063dSJacob Faibussowitsch PetscCall(VecNorm(metric2, NORM_2, &errornorm)); 23900d473e3SJoe Wallwork errornorm /= norm; 24063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Metric normalization L2 error: %.4f%%\n", (double)(100*errornorm))); 24140a2a1afSJoe Wallwork PetscCheck(errornorm < tol, comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric normalization test failed"); 24200d473e3SJoe Wallwork } 24340a2a1afSJoe Wallwork PetscCall(VecDestroy(&determinant)); 24440a2a1afSJoe Wallwork PetscCall(DMDestroy(&dmDet)); 2459566063dSJacob Faibussowitsch PetscCall(VecCopy(metric1, metric)); 2469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&metric2)); 2479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&metric1)); 24800d473e3SJoe Wallwork } 24900d473e3SJoe Wallwork 25000d473e3SJoe Wallwork /* Adapt the mesh */ 2519566063dSJacob Faibussowitsch PetscCall(DMAdaptMetric(dm, metric, bdLabel, rgLabel, &dmAdapt)); 2529566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 2539566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) dmAdapt, "DM_adapted")); 2549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&metric)); 2559566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dmAdapt, NULL, "-adapted_mesh_view")); 25600d473e3SJoe Wallwork 25741473654SJoe Wallwork /* Test tag preservation */ 25841473654SJoe Wallwork if (!noTagging) { 25941473654SJoe Wallwork PetscBool hasTag; 26041473654SJoe Wallwork PetscInt size; 26141473654SJoe Wallwork 2629566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dmAdapt, "Face Sets", &bdLabel)); 2639566063dSJacob Faibussowitsch PetscCall(DMLabelHasStratum(bdLabel, 1, &hasTag)); 26428b400f6SJacob Faibussowitsch PetscCheck(hasTag, comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have face tag 1"); 2659566063dSJacob Faibussowitsch PetscCall(DMLabelHasStratum(bdLabel, 2, &hasTag)); 26628b400f6SJacob Faibussowitsch PetscCheck(hasTag, comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have face tag 2"); 2679566063dSJacob Faibussowitsch PetscCall(DMLabelGetNumValues(bdLabel, &size)); 26863a3b9bcSJacob Faibussowitsch PetscCheck(size == 2, comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh has the wrong number of face tags (got %" PetscInt_FMT ", expected 2)", size); 26941473654SJoe Wallwork 2709566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dmAdapt, "Cell Sets", &rgLabel)); 2719566063dSJacob Faibussowitsch PetscCall(DMLabelHasStratum(rgLabel, 3, &hasTag)); 27228b400f6SJacob Faibussowitsch PetscCheck(hasTag, comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have cell tag 3"); 2739566063dSJacob Faibussowitsch PetscCall(DMLabelHasStratum(rgLabel, 4, &hasTag)); 27428b400f6SJacob Faibussowitsch PetscCheck(hasTag, comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have cell tag 4"); 2759566063dSJacob Faibussowitsch PetscCall(DMLabelGetNumValues(rgLabel, &size)); 27663a3b9bcSJacob Faibussowitsch PetscCheck(size == 2, comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh has the wrong number of cell tags (got %" PetscInt_FMT ", expected 2)", size); 27741473654SJoe Wallwork } 27841473654SJoe Wallwork 27900d473e3SJoe Wallwork /* Clean up */ 2809566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmAdapt)); 2819566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 28200d473e3SJoe Wallwork return 0; 28300d473e3SJoe Wallwork } 28400d473e3SJoe Wallwork 28500d473e3SJoe Wallwork /*TEST 28600d473e3SJoe Wallwork 28780f7b1e7SJoe Wallwork testset: 28880f7b1e7SJoe Wallwork requires: pragmatic 289e600fa54SMatthew G. Knepley args: -dm_plex_box_faces 4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor pragmatic -noTagging 29080f7b1e7SJoe Wallwork 291dd0671eeSJoe Wallwork test: 292dd0671eeSJoe Wallwork suffix: uniform_2d_pragmatic 293d8b80fabSJoe Wallwork args: -dm_plex_metric_uniform 29400d473e3SJoe Wallwork test: 295dd0671eeSJoe Wallwork suffix: iso_2d_pragmatic 296d8b80fabSJoe Wallwork args: -dm_plex_metric_isotropic 29700d473e3SJoe Wallwork test: 298dd0671eeSJoe Wallwork suffix: hessian_2d_pragmatic 29980f7b1e7SJoe Wallwork 30080f7b1e7SJoe Wallwork testset: 30180f7b1e7SJoe Wallwork requires: pragmatic tetgen 302e600fa54SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 4,4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor pragmatic -noTagging 30380f7b1e7SJoe Wallwork 30480f7b1e7SJoe Wallwork test: 30580f7b1e7SJoe Wallwork suffix: uniform_3d_pragmatic 306d8b80fabSJoe Wallwork args: -dm_plex_metric_uniform -noTagging 30780f7b1e7SJoe Wallwork test: 30880f7b1e7SJoe Wallwork suffix: iso_3d_pragmatic 309d8b80fabSJoe Wallwork args: -dm_plex_metric_isotropic -noTagging 31000d473e3SJoe Wallwork test: 311dd0671eeSJoe Wallwork suffix: hessian_3d_pragmatic 31280f7b1e7SJoe Wallwork 31380f7b1e7SJoe Wallwork testset: 31480f7b1e7SJoe Wallwork requires: mmg 315e600fa54SMatthew G. Knepley args: -dm_plex_box_faces 4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor mmg 31680f7b1e7SJoe Wallwork 31700d473e3SJoe Wallwork test: 318dd0671eeSJoe Wallwork suffix: uniform_2d_mmg 319d8b80fabSJoe Wallwork args: -dm_plex_metric_uniform 320dd0671eeSJoe Wallwork test: 321dd0671eeSJoe Wallwork suffix: iso_2d_mmg 322d8b80fabSJoe Wallwork args: -dm_plex_metric_isotropic 323dd0671eeSJoe Wallwork test: 324dd0671eeSJoe Wallwork suffix: hessian_2d_mmg 32580f7b1e7SJoe Wallwork 32680f7b1e7SJoe Wallwork testset: 32780f7b1e7SJoe Wallwork requires: mmg tetgen 328e600fa54SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 4,4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor mmg 32980f7b1e7SJoe Wallwork 33080f7b1e7SJoe Wallwork test: 33180f7b1e7SJoe Wallwork suffix: uniform_3d_mmg 332d8b80fabSJoe Wallwork args: -dm_plex_metric_uniform 33380f7b1e7SJoe Wallwork test: 33480f7b1e7SJoe Wallwork suffix: iso_3d_mmg 335d8b80fabSJoe Wallwork args: -dm_plex_metric_isotropic 336dd0671eeSJoe Wallwork test: 337dd0671eeSJoe Wallwork suffix: hessian_3d_mmg 33880f7b1e7SJoe Wallwork 33980f7b1e7SJoe Wallwork testset: 34080f7b1e7SJoe Wallwork requires: parmmg tetgen 34180f7b1e7SJoe Wallwork nsize: 2 342e600fa54SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 4,4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor parmmg 34380f7b1e7SJoe Wallwork 344dd0671eeSJoe Wallwork test: 345dd0671eeSJoe Wallwork suffix: uniform_3d_parmmg 346d8b80fabSJoe Wallwork args: -dm_plex_metric_uniform 347dd0671eeSJoe Wallwork test: 348dd0671eeSJoe Wallwork suffix: iso_3d_parmmg 349d8b80fabSJoe Wallwork args: -dm_plex_metric_isotropic 350dd0671eeSJoe Wallwork test: 351dd0671eeSJoe Wallwork suffix: hessian_3d_parmmg 35200d473e3SJoe Wallwork 35300d473e3SJoe Wallwork TEST*/ 354