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