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 PetscErrorCode ierr; 19d8b80fabSJoe Wallwork PetscFE fe; 20d8b80fabSJoe Wallwork PetscInt dim; 21d8b80fabSJoe Wallwork 22d8b80fabSJoe Wallwork PetscFunctionBeginUser; 23d8b80fabSJoe Wallwork ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 24d8b80fabSJoe Wallwork ierr = DMClone(dm, dmIndi);CHKERRQ(ierr); 25d8b80fabSJoe Wallwork ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 26d8b80fabSJoe Wallwork ierr = PetscFECreateLagrange(comm, dim, 1, PETSC_TRUE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr); 27d8b80fabSJoe Wallwork ierr = DMSetField(*dmIndi, 0, NULL, (PetscObject)fe);CHKERRQ(ierr); 28d8b80fabSJoe Wallwork ierr = DMCreateDS(*dmIndi);CHKERRQ(ierr); 29d8b80fabSJoe Wallwork ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 30d8b80fabSJoe Wallwork ierr = DMCreateLocalVector(*dmIndi, indicator);CHKERRQ(ierr); 31d8b80fabSJoe Wallwork PetscFunctionReturn(0); 32d8b80fabSJoe Wallwork } 33d8b80fabSJoe Wallwork 3400d473e3SJoe Wallwork int main(int argc, char **argv) { 35*e600fa54SMatthew G. Knepley DM dm, dmAdapt; 3621b3e102SJoe Wallwork DMLabel bdLabel = NULL, rgLabel = NULL; 3700d473e3SJoe Wallwork MPI_Comm comm; 3841473654SJoe Wallwork PetscBool uniform = PETSC_FALSE, isotropic = PETSC_FALSE, noTagging = PETSC_FALSE; 3900d473e3SJoe Wallwork PetscErrorCode ierr; 40*e600fa54SMatthew G. Knepley PetscInt dim; 4100d473e3SJoe Wallwork PetscReal scaling = 1.0; 4200d473e3SJoe Wallwork Vec metric; 4300d473e3SJoe Wallwork 4400d473e3SJoe Wallwork /* Set up */ 4500d473e3SJoe Wallwork ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 4600d473e3SJoe Wallwork comm = PETSC_COMM_WORLD; 4700d473e3SJoe Wallwork ierr = PetscOptionsBegin(comm, "", "Mesh adaptation options", "DMPLEX");CHKERRQ(ierr); 4841473654SJoe Wallwork ierr = PetscOptionsBool("-noTagging", "Should tag preservation testing be turned off?", "ex60.c", noTagging, &noTagging, NULL);CHKERRQ(ierr); 4900d473e3SJoe Wallwork ierr = PetscOptionsEnd(); 5000d473e3SJoe Wallwork 5100d473e3SJoe Wallwork /* Create box mesh */ 52*e600fa54SMatthew G. Knepley ierr = DMCreate(comm, &dm);CHKERRQ(ierr); 53*e600fa54SMatthew G. Knepley ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr); 54e5697243SJoe Wallwork ierr = DMSetFromOptions(dm);CHKERRQ(ierr); 5500d473e3SJoe Wallwork ierr = PetscObjectSetName((PetscObject) dm, "DM_init");CHKERRQ(ierr); 5600d473e3SJoe Wallwork ierr = DMViewFromOptions(dm, NULL, "-initial_mesh_view");CHKERRQ(ierr); 57*e600fa54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 5800d473e3SJoe Wallwork 5941473654SJoe Wallwork /* Set tags to be preserved */ 6041473654SJoe Wallwork if (!noTagging) { 6141473654SJoe Wallwork DM cdm; 6241473654SJoe Wallwork PetscInt cStart, cEnd, c, fStart, fEnd, f, vStart, vEnd; 6341473654SJoe Wallwork const PetscScalar *coords; 6441473654SJoe Wallwork Vec coordinates; 6541473654SJoe Wallwork 6641473654SJoe Wallwork /* Cell tags */ 6741473654SJoe Wallwork ierr = DMCreateLabel(dm, "Cell Sets");CHKERRQ(ierr); 6841473654SJoe Wallwork ierr = DMGetLabel(dm, "Cell Sets", &rgLabel);CHKERRQ(ierr); 6941473654SJoe Wallwork ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 7041473654SJoe Wallwork for (c = cStart; c < cEnd; ++c) { 7141473654SJoe Wallwork PetscReal centroid[3], volume, x; 7241473654SJoe Wallwork 7341473654SJoe Wallwork ierr = DMPlexComputeCellGeometryFVM(dm, c, &volume, centroid, NULL);CHKERRQ(ierr); 7441473654SJoe Wallwork x = centroid[0]; 7541473654SJoe Wallwork if (x < 0.5) { ierr = DMLabelSetValue(rgLabel, c, 3);CHKERRQ(ierr); } 7641473654SJoe Wallwork else { ierr = DMLabelSetValue(rgLabel, c, 4);CHKERRQ(ierr); } 7741473654SJoe Wallwork } 7841473654SJoe Wallwork 7941473654SJoe Wallwork /* Face tags */ 8041473654SJoe Wallwork ierr = DMCreateLabel(dm, "Face Sets");CHKERRQ(ierr); 8141473654SJoe Wallwork ierr = DMGetLabel(dm, "Face Sets", &bdLabel);CHKERRQ(ierr); 8241473654SJoe Wallwork ierr = DMPlexMarkBoundaryFaces(dm, 1, bdLabel);CHKERRQ(ierr); 8341473654SJoe Wallwork ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 8441473654SJoe Wallwork ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 8541473654SJoe Wallwork ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 8641473654SJoe Wallwork ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 8741473654SJoe Wallwork ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 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 9341473654SJoe Wallwork ierr = DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 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; 9941473654SJoe Wallwork ierr = DMPlexPointLocalRead(cdm, off, coords, &x);CHKERRQ(ierr); 10041473654SJoe Wallwork if ((x[0] < 0.5 - eps) || (x[0] > 0.5 + eps)) flg = PETSC_FALSE; 10141473654SJoe Wallwork } 10241473654SJoe Wallwork if (flg) { ierr = DMLabelSetValue(bdLabel, f, 2);CHKERRQ(ierr); } 10341473654SJoe Wallwork ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 10441473654SJoe Wallwork } 10541473654SJoe Wallwork ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 10641473654SJoe Wallwork } 10741473654SJoe Wallwork 10800d473e3SJoe Wallwork /* Construct metric */ 109d8b80fabSJoe Wallwork ierr = DMPlexMetricSetFromOptions(dm);CHKERRQ(ierr); 110d8b80fabSJoe Wallwork ierr = DMPlexMetricIsUniform(dm, &uniform);CHKERRQ(ierr); 111d8b80fabSJoe Wallwork ierr = DMPlexMetricIsIsotropic(dm, &isotropic);CHKERRQ(ierr); 11200d473e3SJoe Wallwork if (uniform) { 113d8b80fabSJoe Wallwork ierr = DMPlexMetricCreateUniform(dm, 0, scaling, &metric);CHKERRQ(ierr); 11400d473e3SJoe Wallwork } 11500d473e3SJoe Wallwork else { 11600d473e3SJoe Wallwork DM dmIndi; 11700d473e3SJoe Wallwork Vec indicator; 11800d473e3SJoe Wallwork 11900d473e3SJoe Wallwork /* Construct "error indicator" */ 120d8b80fabSJoe Wallwork ierr = CreateIndicator(dm, &indicator, &dmIndi);CHKERRQ(ierr); 121e5697243SJoe Wallwork if (isotropic) { 12200d473e3SJoe Wallwork 12300d473e3SJoe Wallwork /* Isotropic case: just specify unity */ 12400d473e3SJoe Wallwork ierr = VecSet(indicator, scaling);CHKERRQ(ierr); 12500d473e3SJoe Wallwork ierr = DMPlexMetricCreateIsotropic(dm, 0, indicator, &metric);CHKERRQ(ierr); 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 */ 13600d473e3SJoe Wallwork ierr = DMProjectFunctionLocal(dmIndi, 0.0, funcs, NULL, INSERT_ALL_VALUES, indicator);CHKERRQ(ierr); 13700d473e3SJoe Wallwork 13800d473e3SJoe Wallwork /* Approximate the gradient */ 13900d473e3SJoe Wallwork ierr = DMClone(dmIndi, &dmGrad);CHKERRQ(ierr); 14000d473e3SJoe Wallwork ierr = PetscFECreateLagrange(comm, dim, dim, PETSC_TRUE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr); 1415767498bSJoe Wallwork ierr = DMSetField(dmGrad, 0, NULL, (PetscObject)fe);CHKERRQ(ierr); 14200d473e3SJoe Wallwork ierr = DMCreateDS(dmGrad);CHKERRQ(ierr); 14300d473e3SJoe Wallwork ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 14400d473e3SJoe Wallwork ierr = DMCreateLocalVector(dmGrad, &gradient);CHKERRQ(ierr); 14500d473e3SJoe Wallwork ierr = DMPlexComputeGradientClementInterpolant(dmIndi, indicator, gradient);CHKERRQ(ierr); 14600d473e3SJoe Wallwork ierr = VecViewFromOptions(gradient, NULL, "-adapt_gradient_view");CHKERRQ(ierr); 14700d473e3SJoe Wallwork 14800d473e3SJoe Wallwork /* Approximate the Hessian */ 1495767498bSJoe Wallwork ierr = DMPlexMetricCreate(dm, 0, &metric);CHKERRQ(ierr); 15000d473e3SJoe Wallwork ierr = DMPlexComputeGradientClementInterpolant(dmGrad, gradient, metric);CHKERRQ(ierr); 15100d473e3SJoe Wallwork ierr = VecViewFromOptions(metric, NULL, "-adapt_hessian_view");CHKERRQ(ierr); 15200d473e3SJoe Wallwork ierr = VecDestroy(&gradient);CHKERRQ(ierr); 15300d473e3SJoe Wallwork ierr = DMDestroy(&dmGrad);CHKERRQ(ierr); 15400d473e3SJoe Wallwork } 15500d473e3SJoe Wallwork ierr = VecDestroy(&indicator);CHKERRQ(ierr); 15600d473e3SJoe Wallwork ierr = DMDestroy(&dmIndi);CHKERRQ(ierr); 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 16600d473e3SJoe Wallwork ierr = VecDuplicate(metric, &metric1);CHKERRQ(ierr); 16700d473e3SJoe Wallwork ierr = VecSet(metric1, 0);CHKERRQ(ierr); 16800d473e3SJoe Wallwork ierr = VecAXPY(metric1, 0.625, metric);CHKERRQ(ierr); 16900d473e3SJoe Wallwork ierr = VecDuplicate(metric, &metric2);CHKERRQ(ierr); 17000d473e3SJoe Wallwork ierr = VecSet(metric2, 0);CHKERRQ(ierr); 17100d473e3SJoe Wallwork ierr = VecAXPY(metric2, 2.5, metric);CHKERRQ(ierr); 17200d473e3SJoe Wallwork metrics[0] = metric1; 17300d473e3SJoe Wallwork metrics[1] = metric2; 17400d473e3SJoe Wallwork 17500d473e3SJoe Wallwork /* Test metric average */ 17600d473e3SJoe Wallwork ierr = DMPlexMetricAverage(dm, 2, weights, metrics, &metricComb);CHKERRQ(ierr); 17700d473e3SJoe Wallwork ierr = VecAXPY(metricComb, -1, metric);CHKERRQ(ierr); 17800d473e3SJoe Wallwork ierr = VecNorm(metric, NORM_2, &norm);CHKERRQ(ierr); 17900d473e3SJoe Wallwork ierr = VecNorm(metricComb, NORM_2, &errornorm);CHKERRQ(ierr); 18000d473e3SJoe Wallwork errornorm /= norm; 1817cdd5544SJoe Wallwork ierr = PetscPrintf(comm, "Metric average L2 error: %.4f%%\n", 100*errornorm);CHKERRQ(ierr); 1822c71b3e2SJacob Faibussowitsch PetscCheckFalse(errornorm > tol,comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric average test failed"); 18300d473e3SJoe Wallwork ierr = VecDestroy(&metricComb);CHKERRQ(ierr); 18400d473e3SJoe Wallwork 18500d473e3SJoe Wallwork /* Test metric intersection */ 186e5697243SJoe Wallwork if (isotropic) { 18700d473e3SJoe Wallwork ierr = DMPlexMetricIntersection(dm, 2, metrics, &metricComb);CHKERRQ(ierr); 18800d473e3SJoe Wallwork ierr = VecAXPY(metricComb, -1, metric1);CHKERRQ(ierr); 18900d473e3SJoe Wallwork ierr = VecNorm(metricComb, NORM_2, &errornorm);CHKERRQ(ierr); 19000d473e3SJoe Wallwork errornorm /= norm; 1917cdd5544SJoe Wallwork ierr = PetscPrintf(comm, "Metric intersection L2 error: %.4f%%\n", 100*errornorm);CHKERRQ(ierr); 1922c71b3e2SJacob Faibussowitsch PetscCheckFalse(errornorm > tol,comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric intersection test failed"); 19300d473e3SJoe Wallwork } 194f49e87b0SJoe Wallwork ierr = VecDestroy(&metric1);CHKERRQ(ierr); 19500d473e3SJoe Wallwork ierr = VecDestroy(&metric2);CHKERRQ(ierr); 19600d473e3SJoe Wallwork ierr = VecDestroy(&metricComb);CHKERRQ(ierr); 19700d473e3SJoe Wallwork 19800d473e3SJoe Wallwork /* Test metric SPD enforcement */ 199f49e87b0SJoe Wallwork ierr = DMPlexMetricEnforceSPD(dm, metric, PETSC_TRUE, PETSC_TRUE, &metric1, &determinant);CHKERRQ(ierr); 200e5697243SJoe Wallwork if (isotropic) { 201f49e87b0SJoe Wallwork Vec err; 202f49e87b0SJoe Wallwork 203f49e87b0SJoe Wallwork ierr = VecDuplicate(determinant, &err);CHKERRQ(ierr); 204f49e87b0SJoe Wallwork ierr = VecSet(err, 1.0);CHKERRQ(ierr); 205f49e87b0SJoe Wallwork ierr = VecNorm(err, NORM_2, &norm);CHKERRQ(ierr); 206f49e87b0SJoe Wallwork ierr = VecAXPY(err, -1, determinant);CHKERRQ(ierr); 207f49e87b0SJoe Wallwork ierr = VecNorm(err, NORM_2, &errornorm);CHKERRQ(ierr); 208f49e87b0SJoe Wallwork ierr = VecDestroy(&err);CHKERRQ(ierr); 209f49e87b0SJoe Wallwork errornorm /= norm; 2107cdd5544SJoe Wallwork ierr = PetscPrintf(comm, "Metric determinant L2 error: %.4f%%\n", 100*errornorm);CHKERRQ(ierr); 2112c71b3e2SJacob Faibussowitsch PetscCheckFalse(errornorm > tol,comm, PETSC_ERR_ARG_OUTOFRANGE, "Determinant is not unit"); 21200d473e3SJoe Wallwork ierr = VecAXPY(metric1, -1, metric);CHKERRQ(ierr); 21300d473e3SJoe Wallwork ierr = VecNorm(metric1, NORM_2, &errornorm);CHKERRQ(ierr); 21400d473e3SJoe Wallwork errornorm /= norm; 2157cdd5544SJoe Wallwork ierr = PetscPrintf(comm, "Metric SPD enforcement L2 error: %.4f%%\n", 100*errornorm);CHKERRQ(ierr); 2162c71b3e2SJacob Faibussowitsch PetscCheckFalse(errornorm > tol,comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric SPD enforcement test failed"); 21700d473e3SJoe Wallwork } 21800d473e3SJoe Wallwork ierr = VecDestroy(&metric1);CHKERRQ(ierr); 219f49e87b0SJoe Wallwork ierr = VecGetDM(determinant, &dmDet);CHKERRQ(ierr); 220f49e87b0SJoe Wallwork ierr = VecDestroy(&determinant);CHKERRQ(ierr); 221f49e87b0SJoe Wallwork ierr = DMDestroy(&dmDet);CHKERRQ(ierr); 22200d473e3SJoe Wallwork 22300d473e3SJoe Wallwork /* Test metric normalization */ 224e5697243SJoe Wallwork ierr = DMPlexMetricNormalize(dm, metric, PETSC_TRUE, PETSC_TRUE, &metric1);CHKERRQ(ierr); 225e5697243SJoe Wallwork if (isotropic) { 226e5697243SJoe Wallwork PetscReal target; 227e5697243SJoe Wallwork 228e5697243SJoe Wallwork ierr = DMPlexMetricGetTargetComplexity(dm, &target);CHKERRQ(ierr); 229e5697243SJoe Wallwork scaling = PetscPowReal(target, 2.0/dim); 230d8b80fabSJoe Wallwork if (uniform) { 23100d473e3SJoe Wallwork ierr = DMPlexMetricCreateUniform(dm, 0, scaling, &metric2);CHKERRQ(ierr); 232d8b80fabSJoe Wallwork } else { 233d8b80fabSJoe Wallwork DM dmIndi; 234d8b80fabSJoe Wallwork Vec indicator; 235d8b80fabSJoe Wallwork 236d8b80fabSJoe Wallwork ierr = CreateIndicator(dm, &indicator, &dmIndi);CHKERRQ(ierr); 237d8b80fabSJoe Wallwork ierr = VecSet(indicator, scaling);CHKERRQ(ierr); 238d8b80fabSJoe Wallwork ierr = DMPlexMetricCreateIsotropic(dm, 0, indicator, &metric2);CHKERRQ(ierr); 239d8b80fabSJoe Wallwork ierr = DMDestroy(&dmIndi);CHKERRQ(ierr); 240d8b80fabSJoe Wallwork ierr = VecDestroy(&indicator);CHKERRQ(ierr); 241d8b80fabSJoe Wallwork } 24200d473e3SJoe Wallwork ierr = VecAXPY(metric2, -1, metric1);CHKERRQ(ierr); 24300d473e3SJoe Wallwork ierr = VecNorm(metric2, NORM_2, &errornorm);CHKERRQ(ierr); 24400d473e3SJoe Wallwork errornorm /= norm; 2457cdd5544SJoe Wallwork ierr = PetscPrintf(comm, "Metric normalization L2 error: %.4f%%\n", 100*errornorm);CHKERRQ(ierr); 2462c71b3e2SJacob Faibussowitsch PetscCheckFalse(errornorm > tol,comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric normalization test failed"); 24700d473e3SJoe Wallwork } 24800d473e3SJoe Wallwork ierr = VecCopy(metric1, metric);CHKERRQ(ierr); 24900d473e3SJoe Wallwork ierr = VecDestroy(&metric2);CHKERRQ(ierr); 25000d473e3SJoe Wallwork ierr = VecDestroy(&metric1);CHKERRQ(ierr); 25100d473e3SJoe Wallwork } 25200d473e3SJoe Wallwork 25300d473e3SJoe Wallwork /* Adapt the mesh */ 25421b3e102SJoe Wallwork ierr = DMAdaptMetric(dm, metric, bdLabel, rgLabel, &dmAdapt);CHKERRQ(ierr); 25521b3e102SJoe Wallwork ierr = DMDestroy(&dm);CHKERRQ(ierr); 25600d473e3SJoe Wallwork ierr = PetscObjectSetName((PetscObject) dmAdapt, "DM_adapted");CHKERRQ(ierr); 25700d473e3SJoe Wallwork ierr = VecDestroy(&metric);CHKERRQ(ierr); 25800d473e3SJoe Wallwork ierr = DMViewFromOptions(dmAdapt, NULL, "-adapted_mesh_view");CHKERRQ(ierr); 25900d473e3SJoe Wallwork 26041473654SJoe Wallwork /* Test tag preservation */ 26141473654SJoe Wallwork if (!noTagging) { 26241473654SJoe Wallwork PetscBool hasTag; 26341473654SJoe Wallwork PetscInt size; 26441473654SJoe Wallwork 26541473654SJoe Wallwork ierr = DMGetLabel(dmAdapt, "Face Sets", &bdLabel);CHKERRQ(ierr); 26641473654SJoe Wallwork ierr = DMLabelHasStratum(bdLabel, 1, &hasTag);CHKERRQ(ierr); 2672c71b3e2SJacob Faibussowitsch PetscCheckFalse(!hasTag,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have face tag 1"); 26841473654SJoe Wallwork ierr = DMLabelHasStratum(bdLabel, 2, &hasTag);CHKERRQ(ierr); 2692c71b3e2SJacob Faibussowitsch PetscCheckFalse(!hasTag,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have face tag 2"); 27041473654SJoe Wallwork ierr = DMLabelGetNumValues(bdLabel, &size);CHKERRQ(ierr); 2712c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 2,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh has the wrong number of face tags (got %d, expected 2)", size); 27241473654SJoe Wallwork 27341473654SJoe Wallwork ierr = DMGetLabel(dmAdapt, "Cell Sets", &rgLabel);CHKERRQ(ierr); 27441473654SJoe Wallwork ierr = DMLabelHasStratum(rgLabel, 3, &hasTag);CHKERRQ(ierr); 2752c71b3e2SJacob Faibussowitsch PetscCheckFalse(!hasTag,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have cell tag 3"); 27641473654SJoe Wallwork ierr = DMLabelHasStratum(rgLabel, 4, &hasTag);CHKERRQ(ierr); 2772c71b3e2SJacob Faibussowitsch PetscCheckFalse(!hasTag,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have cell tag 4"); 27841473654SJoe Wallwork ierr = DMLabelGetNumValues(rgLabel, &size);CHKERRQ(ierr); 2792c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 2,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh has the wrong number of cell tags (got %d, expected 2)", size); 28041473654SJoe Wallwork } 28141473654SJoe Wallwork 28200d473e3SJoe Wallwork /* Clean up */ 28300d473e3SJoe Wallwork ierr = DMDestroy(&dmAdapt);CHKERRQ(ierr); 28400d473e3SJoe Wallwork ierr = PetscFinalize(); 28500d473e3SJoe Wallwork return 0; 28600d473e3SJoe Wallwork } 28700d473e3SJoe Wallwork 28800d473e3SJoe Wallwork /*TEST 28900d473e3SJoe Wallwork 29080f7b1e7SJoe Wallwork testset: 29180f7b1e7SJoe Wallwork requires: pragmatic 292*e600fa54SMatthew G. Knepley args: -dm_plex_box_faces 4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor pragmatic -noTagging 29380f7b1e7SJoe Wallwork 294dd0671eeSJoe Wallwork test: 295dd0671eeSJoe Wallwork suffix: uniform_2d_pragmatic 296d8b80fabSJoe Wallwork args: -dm_plex_metric_uniform 29700d473e3SJoe Wallwork test: 298dd0671eeSJoe Wallwork suffix: iso_2d_pragmatic 299d8b80fabSJoe Wallwork args: -dm_plex_metric_isotropic 30000d473e3SJoe Wallwork test: 301dd0671eeSJoe Wallwork suffix: hessian_2d_pragmatic 30280f7b1e7SJoe Wallwork 30380f7b1e7SJoe Wallwork testset: 30480f7b1e7SJoe Wallwork requires: pragmatic tetgen 305*e600fa54SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 4,4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor pragmatic -noTagging 30680f7b1e7SJoe Wallwork 30780f7b1e7SJoe Wallwork test: 30880f7b1e7SJoe Wallwork suffix: uniform_3d_pragmatic 309d8b80fabSJoe Wallwork args: -dm_plex_metric_uniform -noTagging 31080f7b1e7SJoe Wallwork test: 31180f7b1e7SJoe Wallwork suffix: iso_3d_pragmatic 312d8b80fabSJoe Wallwork args: -dm_plex_metric_isotropic -noTagging 31300d473e3SJoe Wallwork test: 314dd0671eeSJoe Wallwork suffix: hessian_3d_pragmatic 31580f7b1e7SJoe Wallwork 31680f7b1e7SJoe Wallwork testset: 31780f7b1e7SJoe Wallwork requires: mmg 318*e600fa54SMatthew G. Knepley args: -dm_plex_box_faces 4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor mmg 31980f7b1e7SJoe Wallwork 32000d473e3SJoe Wallwork test: 321dd0671eeSJoe Wallwork suffix: uniform_2d_mmg 322d8b80fabSJoe Wallwork args: -dm_plex_metric_uniform 323dd0671eeSJoe Wallwork test: 324dd0671eeSJoe Wallwork suffix: iso_2d_mmg 325d8b80fabSJoe Wallwork args: -dm_plex_metric_isotropic 326dd0671eeSJoe Wallwork test: 327dd0671eeSJoe Wallwork suffix: hessian_2d_mmg 32880f7b1e7SJoe Wallwork 32980f7b1e7SJoe Wallwork testset: 33080f7b1e7SJoe Wallwork requires: mmg tetgen 331*e600fa54SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 4,4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor mmg 33280f7b1e7SJoe Wallwork 33380f7b1e7SJoe Wallwork test: 33480f7b1e7SJoe Wallwork suffix: uniform_3d_mmg 335d8b80fabSJoe Wallwork args: -dm_plex_metric_uniform 33680f7b1e7SJoe Wallwork test: 33780f7b1e7SJoe Wallwork suffix: iso_3d_mmg 338d8b80fabSJoe Wallwork args: -dm_plex_metric_isotropic 339dd0671eeSJoe Wallwork test: 340dd0671eeSJoe Wallwork suffix: hessian_3d_mmg 34180f7b1e7SJoe Wallwork 34280f7b1e7SJoe Wallwork testset: 34380f7b1e7SJoe Wallwork requires: parmmg tetgen 34480f7b1e7SJoe Wallwork nsize: 2 345*e600fa54SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 4,4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor parmmg 34680f7b1e7SJoe Wallwork 347dd0671eeSJoe Wallwork test: 348dd0671eeSJoe Wallwork suffix: uniform_3d_parmmg 349d8b80fabSJoe Wallwork args: -dm_plex_metric_uniform 350dd0671eeSJoe Wallwork test: 351dd0671eeSJoe Wallwork suffix: iso_3d_parmmg 352d8b80fabSJoe Wallwork args: -dm_plex_metric_isotropic 353dd0671eeSJoe Wallwork test: 354dd0671eeSJoe Wallwork suffix: hessian_3d_parmmg 35500d473e3SJoe Wallwork 35600d473e3SJoe Wallwork TEST*/ 357