1 static char help[] = "Test metric utils in the uniform, isotropic case.\n\n"; 2 3 #include <petscdmplex.h> 4 5 static PetscErrorCode bowl(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 6 { 7 PetscInt d; 8 9 *u = 0.0; 10 for (d = 0; d < dim; d++) *u += 0.5*(x[d] - 0.5)*(x[d] - 0.5); 11 12 return 0; 13 } 14 15 static PetscErrorCode CreateIndicator(DM dm, Vec *indicator, DM *dmIndi) 16 { 17 MPI_Comm comm; 18 PetscFE fe; 19 PetscInt dim; 20 21 PetscFunctionBeginUser; 22 PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 23 PetscCall(DMClone(dm, dmIndi)); 24 PetscCall(DMGetDimension(dm, &dim)); 25 PetscCall(PetscFECreateLagrange(comm, dim, 1, PETSC_TRUE, 1, PETSC_DETERMINE, &fe)); 26 PetscCall(DMSetField(*dmIndi, 0, NULL, (PetscObject)fe)); 27 PetscCall(DMCreateDS(*dmIndi)); 28 PetscCall(PetscFEDestroy(&fe)); 29 PetscCall(DMCreateLocalVector(*dmIndi, indicator)); 30 PetscFunctionReturn(0); 31 } 32 33 int main(int argc, char **argv) { 34 DM dm, dmAdapt; 35 DMLabel bdLabel = NULL, rgLabel = NULL; 36 MPI_Comm comm; 37 PetscBool uniform = PETSC_FALSE, isotropic = PETSC_FALSE, noTagging = PETSC_FALSE; 38 PetscInt dim; 39 PetscReal scaling = 1.0; 40 Vec metric; 41 42 /* Set up */ 43 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 44 comm = PETSC_COMM_WORLD; 45 PetscOptionsBegin(comm, "", "Mesh adaptation options", "DMPLEX"); 46 PetscCall(PetscOptionsBool("-noTagging", "Should tag preservation testing be turned off?", "ex60.c", noTagging, &noTagging, NULL)); 47 PetscOptionsEnd(); 48 49 /* Create box mesh */ 50 PetscCall(DMCreate(comm, &dm)); 51 PetscCall(DMSetType(dm, DMPLEX)); 52 PetscCall(DMSetFromOptions(dm)); 53 PetscCall(PetscObjectSetName((PetscObject) dm, "DM_init")); 54 PetscCall(DMViewFromOptions(dm, NULL, "-initial_mesh_view")); 55 PetscCall(DMGetDimension(dm, &dim)); 56 57 /* Set tags to be preserved */ 58 if (!noTagging) { 59 DM cdm; 60 PetscInt cStart, cEnd, c, fStart, fEnd, f, vStart, vEnd; 61 const PetscScalar *coords; 62 Vec coordinates; 63 64 /* Cell tags */ 65 PetscCall(DMCreateLabel(dm, "Cell Sets")); 66 PetscCall(DMGetLabel(dm, "Cell Sets", &rgLabel)); 67 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 68 for (c = cStart; c < cEnd; ++c) { 69 PetscReal centroid[3], volume, x; 70 71 PetscCall(DMPlexComputeCellGeometryFVM(dm, c, &volume, centroid, NULL)); 72 x = centroid[0]; 73 if (x < 0.5) PetscCall(DMLabelSetValue(rgLabel, c, 3)); 74 else PetscCall(DMLabelSetValue(rgLabel, c, 4)); 75 } 76 77 /* Face tags */ 78 PetscCall(DMCreateLabel(dm, "Face Sets")); 79 PetscCall(DMGetLabel(dm, "Face Sets", &bdLabel)); 80 PetscCall(DMPlexMarkBoundaryFaces(dm, 1, bdLabel)); 81 PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd)); 82 PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 83 PetscCall(DMGetCoordinateDM(dm, &cdm)); 84 PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 85 PetscCall(VecGetArrayRead(coordinates, &coords)); 86 for (f = fStart; f < fEnd; ++f) { 87 PetscBool flg = PETSC_TRUE; 88 PetscInt *closure = NULL, closureSize, cl; 89 PetscReal eps = 1.0e-08; 90 91 PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure)); 92 for (cl = 0; cl < closureSize*2; cl += 2) { 93 PetscInt off = closure[cl]; 94 PetscReal *x; 95 96 if ((off < vStart) || (off >= vEnd)) continue; 97 PetscCall(DMPlexPointLocalRead(cdm, off, coords, &x)); 98 if ((x[0] < 0.5 - eps) || (x[0] > 0.5 + eps)) flg = PETSC_FALSE; 99 } 100 if (flg) PetscCall(DMLabelSetValue(bdLabel, f, 2)); 101 PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure)); 102 } 103 PetscCall(VecRestoreArrayRead(coordinates, &coords)); 104 } 105 106 /* Construct metric */ 107 PetscCall(DMPlexMetricSetFromOptions(dm)); 108 PetscCall(DMPlexMetricIsUniform(dm, &uniform)); 109 PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic)); 110 if (uniform) { 111 PetscCall(DMPlexMetricCreateUniform(dm, 0, scaling, &metric)); 112 } 113 else { 114 DM dmIndi; 115 Vec indicator; 116 117 /* Construct "error indicator" */ 118 PetscCall(CreateIndicator(dm, &indicator, &dmIndi)); 119 if (isotropic) { 120 121 /* Isotropic case: just specify unity */ 122 PetscCall(VecSet(indicator, scaling)); 123 PetscCall(DMPlexMetricCreateIsotropic(dm, 0, indicator, &metric)); 124 125 } else { 126 PetscFE fe; 127 128 /* 'Anisotropic' case: approximate the identity by recovering the Hessian of a parabola */ 129 DM dmGrad; 130 PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar*, void*) = {bowl}; 131 Vec gradient; 132 133 /* Project the parabola into P1 space */ 134 PetscCall(DMProjectFunctionLocal(dmIndi, 0.0, funcs, NULL, INSERT_ALL_VALUES, indicator)); 135 136 /* Approximate the gradient */ 137 PetscCall(DMClone(dmIndi, &dmGrad)); 138 PetscCall(PetscFECreateLagrange(comm, dim, dim, PETSC_TRUE, 1, PETSC_DETERMINE, &fe)); 139 PetscCall(DMSetField(dmGrad, 0, NULL, (PetscObject)fe)); 140 PetscCall(DMCreateDS(dmGrad)); 141 PetscCall(PetscFEDestroy(&fe)); 142 PetscCall(DMCreateLocalVector(dmGrad, &gradient)); 143 PetscCall(DMPlexComputeGradientClementInterpolant(dmIndi, indicator, gradient)); 144 PetscCall(VecViewFromOptions(gradient, NULL, "-adapt_gradient_view")); 145 146 /* Approximate the Hessian */ 147 PetscCall(DMPlexMetricCreate(dm, 0, &metric)); 148 PetscCall(DMPlexComputeGradientClementInterpolant(dmGrad, gradient, metric)); 149 PetscCall(VecViewFromOptions(metric, NULL, "-adapt_hessian_view")); 150 PetscCall(VecDestroy(&gradient)); 151 PetscCall(DMDestroy(&dmGrad)); 152 } 153 PetscCall(VecDestroy(&indicator)); 154 PetscCall(DMDestroy(&dmIndi)); 155 } 156 157 /* Test metric routines */ 158 { 159 DM dmDet; 160 PetscReal errornorm, norm, tol = 1.0e-10, weights[2] = {0.8, 0.2}; 161 Vec metric1, metric2, metricComb, determinant; 162 Vec metrics[2]; 163 164 PetscCall(VecDuplicate(metric, &metric1)); 165 PetscCall(VecSet(metric1, 0)); 166 PetscCall(VecAXPY(metric1, 0.625, metric)); 167 PetscCall(VecDuplicate(metric, &metric2)); 168 PetscCall(VecSet(metric2, 0)); 169 PetscCall(VecAXPY(metric2, 2.5, metric)); 170 metrics[0] = metric1; 171 metrics[1] = metric2; 172 173 /* Test metric average */ 174 PetscCall(DMPlexMetricAverage(dm, 2, weights, metrics, &metricComb)); 175 PetscCall(VecAXPY(metricComb, -1, metric)); 176 PetscCall(VecNorm(metric, NORM_2, &norm)); 177 PetscCall(VecNorm(metricComb, NORM_2, &errornorm)); 178 errornorm /= norm; 179 PetscCall(PetscPrintf(comm, "Metric average L2 error: %.4f%%\n", (double)(100*errornorm))); 180 PetscCheck(errornorm <= tol,comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric average test failed"); 181 PetscCall(VecDestroy(&metricComb)); 182 183 /* Test metric intersection */ 184 if (isotropic) { 185 PetscCall(DMPlexMetricIntersection(dm, 2, metrics, &metricComb)); 186 PetscCall(VecAXPY(metricComb, -1, metric1)); 187 PetscCall(VecNorm(metricComb, NORM_2, &errornorm)); 188 errornorm /= norm; 189 PetscCall(PetscPrintf(comm, "Metric intersection L2 error: %.4f%%\n", (double)(100*errornorm))); 190 PetscCheck(errornorm <= tol,comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric intersection test failed"); 191 } 192 PetscCall(VecDestroy(&metric1)); 193 PetscCall(VecDestroy(&metric2)); 194 PetscCall(VecDestroy(&metricComb)); 195 196 /* Test metric SPD enforcement */ 197 PetscCall(DMPlexMetricEnforceSPD(dm, metric, PETSC_TRUE, PETSC_TRUE, &metric1, &determinant)); 198 if (isotropic) { 199 Vec err; 200 201 PetscCall(VecDuplicate(determinant, &err)); 202 PetscCall(VecSet(err, 1.0)); 203 PetscCall(VecNorm(err, NORM_2, &norm)); 204 PetscCall(VecAXPY(err, -1, determinant)); 205 PetscCall(VecNorm(err, NORM_2, &errornorm)); 206 PetscCall(VecDestroy(&err)); 207 errornorm /= norm; 208 PetscCall(PetscPrintf(comm, "Metric determinant L2 error: %.4f%%\n", (double)(100*errornorm))); 209 PetscCheck(errornorm <= tol,comm, PETSC_ERR_ARG_OUTOFRANGE, "Determinant is not unit"); 210 PetscCall(VecAXPY(metric1, -1, metric)); 211 PetscCall(VecNorm(metric1, NORM_2, &errornorm)); 212 errornorm /= norm; 213 PetscCall(PetscPrintf(comm, "Metric SPD enforcement L2 error: %.4f%%\n", (double)(100*errornorm))); 214 PetscCheck(errornorm <= tol,comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric SPD enforcement test failed"); 215 } 216 PetscCall(VecDestroy(&metric1)); 217 PetscCall(VecGetDM(determinant, &dmDet)); 218 PetscCall(VecDestroy(&determinant)); 219 PetscCall(DMDestroy(&dmDet)); 220 221 /* Test metric normalization */ 222 PetscCall(DMPlexMetricNormalize(dm, metric, PETSC_TRUE, PETSC_TRUE, &metric1)); 223 if (isotropic) { 224 PetscReal target; 225 226 PetscCall(DMPlexMetricGetTargetComplexity(dm, &target)); 227 scaling = PetscPowReal(target, 2.0/dim); 228 if (uniform) { 229 PetscCall(DMPlexMetricCreateUniform(dm, 0, scaling, &metric2)); 230 } else { 231 DM dmIndi; 232 Vec indicator; 233 234 PetscCall(CreateIndicator(dm, &indicator, &dmIndi)); 235 PetscCall(VecSet(indicator, scaling)); 236 PetscCall(DMPlexMetricCreateIsotropic(dm, 0, indicator, &metric2)); 237 PetscCall(DMDestroy(&dmIndi)); 238 PetscCall(VecDestroy(&indicator)); 239 } 240 PetscCall(VecAXPY(metric2, -1, metric1)); 241 PetscCall(VecNorm(metric2, NORM_2, &errornorm)); 242 errornorm /= norm; 243 PetscCall(PetscPrintf(comm, "Metric normalization L2 error: %.4f%%\n", (double)(100*errornorm))); 244 PetscCheck(errornorm <= tol,comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric normalization test failed"); 245 } 246 PetscCall(VecCopy(metric1, metric)); 247 PetscCall(VecDestroy(&metric2)); 248 PetscCall(VecDestroy(&metric1)); 249 } 250 251 /* Adapt the mesh */ 252 PetscCall(DMAdaptMetric(dm, metric, bdLabel, rgLabel, &dmAdapt)); 253 PetscCall(DMDestroy(&dm)); 254 PetscCall(PetscObjectSetName((PetscObject) dmAdapt, "DM_adapted")); 255 PetscCall(VecDestroy(&metric)); 256 PetscCall(DMViewFromOptions(dmAdapt, NULL, "-adapted_mesh_view")); 257 258 /* Test tag preservation */ 259 if (!noTagging) { 260 PetscBool hasTag; 261 PetscInt size; 262 263 PetscCall(DMGetLabel(dmAdapt, "Face Sets", &bdLabel)); 264 PetscCall(DMLabelHasStratum(bdLabel, 1, &hasTag)); 265 PetscCheck(hasTag,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have face tag 1"); 266 PetscCall(DMLabelHasStratum(bdLabel, 2, &hasTag)); 267 PetscCheck(hasTag,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have face tag 2"); 268 PetscCall(DMLabelGetNumValues(bdLabel, &size)); 269 PetscCheck(size == 2,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh has the wrong number of face tags (got %" PetscInt_FMT ", expected 2)", size); 270 271 PetscCall(DMGetLabel(dmAdapt, "Cell Sets", &rgLabel)); 272 PetscCall(DMLabelHasStratum(rgLabel, 3, &hasTag)); 273 PetscCheck(hasTag,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have cell tag 3"); 274 PetscCall(DMLabelHasStratum(rgLabel, 4, &hasTag)); 275 PetscCheck(hasTag,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have cell tag 4"); 276 PetscCall(DMLabelGetNumValues(rgLabel, &size)); 277 PetscCheck(size == 2,comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh has the wrong number of cell tags (got %" PetscInt_FMT ", expected 2)", size); 278 } 279 280 /* Clean up */ 281 PetscCall(DMDestroy(&dmAdapt)); 282 PetscCall(PetscFinalize()); 283 return 0; 284 } 285 286 /*TEST 287 288 testset: 289 requires: pragmatic 290 args: -dm_plex_box_faces 4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor pragmatic -noTagging 291 292 test: 293 suffix: uniform_2d_pragmatic 294 args: -dm_plex_metric_uniform 295 test: 296 suffix: iso_2d_pragmatic 297 args: -dm_plex_metric_isotropic 298 test: 299 suffix: hessian_2d_pragmatic 300 301 testset: 302 requires: pragmatic tetgen 303 args: -dm_plex_dim 3 -dm_plex_box_faces 4,4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor pragmatic -noTagging 304 305 test: 306 suffix: uniform_3d_pragmatic 307 args: -dm_plex_metric_uniform -noTagging 308 test: 309 suffix: iso_3d_pragmatic 310 args: -dm_plex_metric_isotropic -noTagging 311 test: 312 suffix: hessian_3d_pragmatic 313 314 testset: 315 requires: mmg 316 args: -dm_plex_box_faces 4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor mmg 317 318 test: 319 suffix: uniform_2d_mmg 320 args: -dm_plex_metric_uniform 321 test: 322 suffix: iso_2d_mmg 323 args: -dm_plex_metric_isotropic 324 test: 325 suffix: hessian_2d_mmg 326 327 testset: 328 requires: mmg tetgen 329 args: -dm_plex_dim 3 -dm_plex_box_faces 4,4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor mmg 330 331 test: 332 suffix: uniform_3d_mmg 333 args: -dm_plex_metric_uniform 334 test: 335 suffix: iso_3d_mmg 336 args: -dm_plex_metric_isotropic 337 test: 338 suffix: hessian_3d_mmg 339 340 testset: 341 requires: parmmg tetgen 342 nsize: 2 343 args: -dm_plex_dim 3 -dm_plex_box_faces 4,4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor parmmg 344 345 test: 346 suffix: uniform_3d_parmmg 347 args: -dm_plex_metric_uniform 348 test: 349 suffix: iso_3d_parmmg 350 args: -dm_plex_metric_isotropic 351 test: 352 suffix: hessian_3d_parmmg 353 354 TEST*/ 355