1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 3 PetscErrorCode DMPlexInvertCell_Internal(PetscInt dim, PetscInt numCorners, PetscInt cone[]) 4 { 5 int tmpc; 6 7 PetscFunctionBegin; 8 if (dim != 3) PetscFunctionReturn(0); 9 switch (numCorners) { 10 case 4: 11 tmpc = cone[0]; 12 cone[0] = cone[1]; 13 cone[1] = tmpc; 14 break; 15 case 8: 16 tmpc = cone[1]; 17 cone[1] = cone[3]; 18 cone[3] = tmpc; 19 break; 20 default: break; 21 } 22 PetscFunctionReturn(0); 23 } 24 25 /*@C 26 DMPlexInvertCell - This flips tetrahedron and hexahedron orientation since Plex stores them internally with outward normals. Other cells are left untouched. 27 28 Input Parameters: 29 + numCorners - The number of vertices in a cell 30 - cone - The incoming cone 31 32 Output Parameter: 33 . cone - The inverted cone (in-place) 34 35 Level: developer 36 37 .seealso: DMPlexGenerate() 38 @*/ 39 PetscErrorCode DMPlexInvertCell(PetscInt dim, PetscInt numCorners, int cone[]) 40 { 41 int tmpc; 42 43 PetscFunctionBegin; 44 if (dim != 3) PetscFunctionReturn(0); 45 switch (numCorners) { 46 case 4: 47 tmpc = cone[0]; 48 cone[0] = cone[1]; 49 cone[1] = tmpc; 50 break; 51 case 8: 52 tmpc = cone[1]; 53 cone[1] = cone[3]; 54 cone[3] = tmpc; 55 break; 56 default: break; 57 } 58 PetscFunctionReturn(0); 59 } 60 61 /* This is to fix the tetrahedron orientation from TetGen */ 62 PETSC_UNUSED static PetscErrorCode DMPlexInvertCells_Internal(PetscInt dim, PetscInt numCells, PetscInt numCorners, int cells[]) 63 { 64 PetscInt bound = numCells*numCorners, coff; 65 PetscErrorCode ierr; 66 67 PetscFunctionBegin; 68 for (coff = 0; coff < bound; coff += numCorners) { 69 ierr = DMPlexInvertCell(dim, numCorners, &cells[coff]);CHKERRQ(ierr); 70 } 71 PetscFunctionReturn(0); 72 } 73 74 /*@ 75 DMPlexTriangleSetOptions - Set the options used for the Triangle mesh generator 76 77 Not Collective 78 79 Inputs Parameters: 80 + dm - The DMPlex object 81 - opts - The command line options 82 83 Level: developer 84 85 .keywords: mesh, points 86 .seealso: DMPlexTetgenSetOptions(), DMPlexGenerate() 87 @*/ 88 PetscErrorCode DMPlexTriangleSetOptions(DM dm, const char *opts) 89 { 90 DM_Plex *mesh = (DM_Plex*) dm->data; 91 PetscErrorCode ierr; 92 93 PetscFunctionBegin; 94 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 95 PetscValidPointer(opts, 2); 96 ierr = PetscFree(mesh->triangleOpts);CHKERRQ(ierr); 97 ierr = PetscStrallocpy(opts, &mesh->triangleOpts);CHKERRQ(ierr); 98 PetscFunctionReturn(0); 99 } 100 101 /*@ 102 DMPlexTetgenSetOptions - Set the options used for the Tetgen mesh generator 103 104 Not Collective 105 106 Inputs Parameters: 107 + dm - The DMPlex object 108 - opts - The command line options 109 110 Level: developer 111 112 .keywords: mesh, points 113 .seealso: DMPlexTriangleSetOptions(), DMPlexGenerate() 114 @*/ 115 PetscErrorCode DMPlexTetgenSetOptions(DM dm, const char *opts) 116 { 117 DM_Plex *mesh = (DM_Plex*) dm->data; 118 PetscErrorCode ierr; 119 120 PetscFunctionBegin; 121 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 122 PetscValidPointer(opts, 2); 123 ierr = PetscFree(mesh->tetgenOpts);CHKERRQ(ierr); 124 ierr = PetscStrallocpy(opts, &mesh->tetgenOpts);CHKERRQ(ierr); 125 PetscFunctionReturn(0); 126 } 127 128 #if defined(PETSC_HAVE_TRIANGLE) 129 #include <triangle.h> 130 131 static PetscErrorCode InitInput_Triangle(struct triangulateio *inputCtx) 132 { 133 PetscFunctionBegin; 134 inputCtx->numberofpoints = 0; 135 inputCtx->numberofpointattributes = 0; 136 inputCtx->pointlist = NULL; 137 inputCtx->pointattributelist = NULL; 138 inputCtx->pointmarkerlist = NULL; 139 inputCtx->numberofsegments = 0; 140 inputCtx->segmentlist = NULL; 141 inputCtx->segmentmarkerlist = NULL; 142 inputCtx->numberoftriangleattributes = 0; 143 inputCtx->trianglelist = NULL; 144 inputCtx->numberofholes = 0; 145 inputCtx->holelist = NULL; 146 inputCtx->numberofregions = 0; 147 inputCtx->regionlist = NULL; 148 PetscFunctionReturn(0); 149 } 150 151 static PetscErrorCode InitOutput_Triangle(struct triangulateio *outputCtx) 152 { 153 PetscFunctionBegin; 154 outputCtx->numberofpoints = 0; 155 outputCtx->pointlist = NULL; 156 outputCtx->pointattributelist = NULL; 157 outputCtx->pointmarkerlist = NULL; 158 outputCtx->numberoftriangles = 0; 159 outputCtx->trianglelist = NULL; 160 outputCtx->triangleattributelist = NULL; 161 outputCtx->neighborlist = NULL; 162 outputCtx->segmentlist = NULL; 163 outputCtx->segmentmarkerlist = NULL; 164 outputCtx->numberofedges = 0; 165 outputCtx->edgelist = NULL; 166 outputCtx->edgemarkerlist = NULL; 167 PetscFunctionReturn(0); 168 } 169 170 static PetscErrorCode FiniOutput_Triangle(struct triangulateio *outputCtx) 171 { 172 PetscFunctionBegin; 173 free(outputCtx->pointlist); 174 free(outputCtx->pointmarkerlist); 175 free(outputCtx->segmentlist); 176 free(outputCtx->segmentmarkerlist); 177 free(outputCtx->edgelist); 178 free(outputCtx->edgemarkerlist); 179 free(outputCtx->trianglelist); 180 free(outputCtx->neighborlist); 181 PetscFunctionReturn(0); 182 } 183 184 static PetscErrorCode DMPlexGenerate_Triangle(DM boundary, PetscBool interpolate, DM *dm) 185 { 186 MPI_Comm comm; 187 DM_Plex *mesh = (DM_Plex *) boundary->data; 188 PetscInt dim = 2; 189 const PetscBool createConvexHull = PETSC_FALSE; 190 const PetscBool constrained = PETSC_FALSE; 191 const char *labelName = "marker"; 192 struct triangulateio in; 193 struct triangulateio out; 194 DMLabel label; 195 PetscInt vStart, vEnd, v, eStart, eEnd, e; 196 PetscMPIInt rank; 197 PetscErrorCode ierr; 198 199 PetscFunctionBegin; 200 ierr = PetscObjectGetComm((PetscObject)boundary,&comm);CHKERRQ(ierr); 201 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 202 ierr = InitInput_Triangle(&in);CHKERRQ(ierr); 203 ierr = InitOutput_Triangle(&out);CHKERRQ(ierr); 204 ierr = DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd);CHKERRQ(ierr); 205 ierr = DMGetLabel(boundary, labelName, &label);CHKERRQ(ierr); 206 207 in.numberofpoints = vEnd - vStart; 208 if (in.numberofpoints > 0) { 209 PetscSection coordSection; 210 Vec coordinates; 211 PetscScalar *array; 212 213 ierr = PetscMalloc1(in.numberofpoints*dim, &in.pointlist);CHKERRQ(ierr); 214 ierr = PetscMalloc1(in.numberofpoints, &in.pointmarkerlist);CHKERRQ(ierr); 215 ierr = DMGetCoordinatesLocal(boundary, &coordinates);CHKERRQ(ierr); 216 ierr = DMGetCoordinateSection(boundary, &coordSection);CHKERRQ(ierr); 217 ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr); 218 for (v = vStart; v < vEnd; ++v) { 219 const PetscInt idx = v - vStart; 220 PetscInt off, d; 221 222 ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 223 for (d = 0; d < dim; ++d) { 224 in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]); 225 } 226 if (label) {ierr = DMLabelGetValue(label, v, &in.pointmarkerlist[idx]);CHKERRQ(ierr);} 227 } 228 ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr); 229 } 230 ierr = DMPlexGetHeightStratum(boundary, 0, &eStart, &eEnd);CHKERRQ(ierr); 231 in.numberofsegments = eEnd - eStart; 232 if (in.numberofsegments > 0) { 233 ierr = PetscMalloc1(in.numberofsegments*2, &in.segmentlist);CHKERRQ(ierr); 234 ierr = PetscMalloc1(in.numberofsegments, &in.segmentmarkerlist);CHKERRQ(ierr); 235 for (e = eStart; e < eEnd; ++e) { 236 const PetscInt idx = e - eStart; 237 const PetscInt *cone; 238 239 ierr = DMPlexGetCone(boundary, e, &cone);CHKERRQ(ierr); 240 241 in.segmentlist[idx*2+0] = cone[0] - vStart; 242 in.segmentlist[idx*2+1] = cone[1] - vStart; 243 244 if (label) {ierr = DMLabelGetValue(label, e, &in.segmentmarkerlist[idx]);CHKERRQ(ierr);} 245 } 246 } 247 #if 0 /* Do not currently support holes */ 248 PetscReal *holeCoords; 249 PetscInt h, d; 250 251 ierr = DMPlexGetHoles(boundary, &in.numberofholes, &holeCords);CHKERRQ(ierr); 252 if (in.numberofholes > 0) { 253 ierr = PetscMalloc1(in.numberofholes*dim, &in.holelist);CHKERRQ(ierr); 254 for (h = 0; h < in.numberofholes; ++h) { 255 for (d = 0; d < dim; ++d) { 256 in.holelist[h*dim+d] = holeCoords[h*dim+d]; 257 } 258 } 259 } 260 #endif 261 if (!rank) { 262 char args[32]; 263 264 /* Take away 'Q' for verbose output */ 265 ierr = PetscStrcpy(args, "pqezQ");CHKERRQ(ierr); 266 if (createConvexHull) {ierr = PetscStrcat(args, "c");CHKERRQ(ierr);} 267 if (constrained) {ierr = PetscStrcpy(args, "zepDQ");CHKERRQ(ierr);} 268 if (mesh->triangleOpts) {triangulate(mesh->triangleOpts, &in, &out, NULL);} 269 else {triangulate(args, &in, &out, NULL);} 270 } 271 ierr = PetscFree(in.pointlist);CHKERRQ(ierr); 272 ierr = PetscFree(in.pointmarkerlist);CHKERRQ(ierr); 273 ierr = PetscFree(in.segmentlist);CHKERRQ(ierr); 274 ierr = PetscFree(in.segmentmarkerlist);CHKERRQ(ierr); 275 ierr = PetscFree(in.holelist);CHKERRQ(ierr); 276 277 { 278 DMLabel glabel = NULL; 279 const PetscInt numCorners = 3; 280 const PetscInt numCells = out.numberoftriangles; 281 const PetscInt numVertices = out.numberofpoints; 282 const int *cells = out.trianglelist; 283 const double *meshCoords = out.pointlist; 284 285 ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm);CHKERRQ(ierr); 286 if (label) {ierr = DMCreateLabel(*dm, labelName); ierr = DMGetLabel(*dm, labelName, &glabel);} 287 /* Set labels */ 288 for (v = 0; v < numVertices; ++v) { 289 if (out.pointmarkerlist[v]) { 290 if (glabel) {ierr = DMLabelSetValue(glabel, v+numCells, out.pointmarkerlist[v]);CHKERRQ(ierr);} 291 } 292 } 293 if (interpolate) { 294 for (e = 0; e < out.numberofedges; e++) { 295 if (out.edgemarkerlist[e]) { 296 const PetscInt vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells}; 297 const PetscInt *edges; 298 PetscInt numEdges; 299 300 ierr = DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 301 if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges); 302 if (glabel) {ierr = DMLabelSetValue(glabel, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);} 303 ierr = DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 304 } 305 } 306 } 307 ierr = DMPlexSetRefinementUniform(*dm, PETSC_FALSE);CHKERRQ(ierr); 308 } 309 #if 0 /* Do not currently support holes */ 310 ierr = DMPlexCopyHoles(*dm, boundary);CHKERRQ(ierr); 311 #endif 312 ierr = FiniOutput_Triangle(&out);CHKERRQ(ierr); 313 PetscFunctionReturn(0); 314 } 315 316 static PetscErrorCode DMPlexRefine_Triangle(DM dm, double *maxVolumes, DM *dmRefined) 317 { 318 MPI_Comm comm; 319 PetscInt dim = 2; 320 const char *labelName = "marker"; 321 struct triangulateio in; 322 struct triangulateio out; 323 DMLabel label; 324 PetscInt vStart, vEnd, v, cStart, cEnd, c, depth, depthGlobal; 325 PetscMPIInt rank; 326 PetscErrorCode ierr; 327 328 PetscFunctionBegin; 329 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 330 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 331 ierr = InitInput_Triangle(&in);CHKERRQ(ierr); 332 ierr = InitOutput_Triangle(&out);CHKERRQ(ierr); 333 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 334 ierr = MPIU_Allreduce(&depth, &depthGlobal, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr); 335 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 336 ierr = DMGetLabel(dm, labelName, &label);CHKERRQ(ierr); 337 338 in.numberofpoints = vEnd - vStart; 339 if (in.numberofpoints > 0) { 340 PetscSection coordSection; 341 Vec coordinates; 342 PetscScalar *array; 343 344 ierr = PetscMalloc1(in.numberofpoints*dim, &in.pointlist);CHKERRQ(ierr); 345 ierr = PetscMalloc1(in.numberofpoints, &in.pointmarkerlist);CHKERRQ(ierr); 346 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 347 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 348 ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr); 349 for (v = vStart; v < vEnd; ++v) { 350 const PetscInt idx = v - vStart; 351 PetscInt off, d; 352 353 ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 354 for (d = 0; d < dim; ++d) { 355 in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]); 356 } 357 if (label) {ierr = DMLabelGetValue(label, v, &in.pointmarkerlist[idx]);CHKERRQ(ierr);} 358 } 359 ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr); 360 } 361 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 362 363 in.numberofcorners = 3; 364 in.numberoftriangles = cEnd - cStart; 365 366 in.trianglearealist = (double*) maxVolumes; 367 if (in.numberoftriangles > 0) { 368 ierr = PetscMalloc1(in.numberoftriangles*in.numberofcorners, &in.trianglelist);CHKERRQ(ierr); 369 for (c = cStart; c < cEnd; ++c) { 370 const PetscInt idx = c - cStart; 371 PetscInt *closure = NULL; 372 PetscInt closureSize; 373 374 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 375 if ((closureSize != 4) && (closureSize != 7)) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Mesh has cell which is not a triangle, %D vertices in closure", closureSize); 376 for (v = 0; v < 3; ++v) { 377 in.trianglelist[idx*in.numberofcorners + v] = closure[(v+closureSize-3)*2] - vStart; 378 } 379 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 380 } 381 } 382 /* TODO: Segment markers are missing on input */ 383 #if 0 /* Do not currently support holes */ 384 PetscReal *holeCoords; 385 PetscInt h, d; 386 387 ierr = DMPlexGetHoles(boundary, &in.numberofholes, &holeCords);CHKERRQ(ierr); 388 if (in.numberofholes > 0) { 389 ierr = PetscMalloc1(in.numberofholes*dim, &in.holelist);CHKERRQ(ierr); 390 for (h = 0; h < in.numberofholes; ++h) { 391 for (d = 0; d < dim; ++d) { 392 in.holelist[h*dim+d] = holeCoords[h*dim+d]; 393 } 394 } 395 } 396 #endif 397 if (!rank) { 398 char args[32]; 399 400 /* Take away 'Q' for verbose output */ 401 ierr = PetscStrcpy(args, "pqezQra");CHKERRQ(ierr); 402 triangulate(args, &in, &out, NULL); 403 } 404 ierr = PetscFree(in.pointlist);CHKERRQ(ierr); 405 ierr = PetscFree(in.pointmarkerlist);CHKERRQ(ierr); 406 ierr = PetscFree(in.segmentlist);CHKERRQ(ierr); 407 ierr = PetscFree(in.segmentmarkerlist);CHKERRQ(ierr); 408 ierr = PetscFree(in.trianglelist);CHKERRQ(ierr); 409 410 { 411 DMLabel rlabel = NULL; 412 const PetscInt numCorners = 3; 413 const PetscInt numCells = out.numberoftriangles; 414 const PetscInt numVertices = out.numberofpoints; 415 const int *cells = out.trianglelist; 416 const double *meshCoords = out.pointlist; 417 PetscBool interpolate = depthGlobal > 1 ? PETSC_TRUE : PETSC_FALSE; 418 419 ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined);CHKERRQ(ierr); 420 if (label) {ierr = DMCreateLabel(*dmRefined, labelName); ierr = DMGetLabel(*dmRefined, labelName, &rlabel);} 421 /* Set labels */ 422 for (v = 0; v < numVertices; ++v) { 423 if (out.pointmarkerlist[v]) { 424 if (rlabel) {ierr = DMLabelSetValue(rlabel, v+numCells, out.pointmarkerlist[v]);CHKERRQ(ierr);} 425 } 426 } 427 if (interpolate) { 428 PetscInt e; 429 430 for (e = 0; e < out.numberofedges; e++) { 431 if (out.edgemarkerlist[e]) { 432 const PetscInt vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells}; 433 const PetscInt *edges; 434 PetscInt numEdges; 435 436 ierr = DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 437 if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges); 438 if (rlabel) {ierr = DMLabelSetValue(rlabel, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);} 439 ierr = DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 440 } 441 } 442 } 443 ierr = DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE);CHKERRQ(ierr); 444 } 445 #if 0 /* Do not currently support holes */ 446 ierr = DMPlexCopyHoles(*dm, boundary);CHKERRQ(ierr); 447 #endif 448 ierr = FiniOutput_Triangle(&out);CHKERRQ(ierr); 449 PetscFunctionReturn(0); 450 } 451 #endif /* PETSC_HAVE_TRIANGLE */ 452 453 #if defined(PETSC_HAVE_TETGEN) 454 #include <tetgen.h> 455 static PetscErrorCode DMPlexGenerate_Tetgen(DM boundary, PetscBool interpolate, DM *dm) 456 { 457 MPI_Comm comm; 458 DM_Plex *mesh = (DM_Plex *) boundary->data; 459 const PetscInt dim = 3; 460 const char *labelName = "marker"; 461 ::tetgenio in; 462 ::tetgenio out; 463 DMLabel label; 464 PetscInt vStart, vEnd, v, fStart, fEnd, f; 465 PetscMPIInt rank; 466 PetscErrorCode ierr; 467 468 PetscFunctionBegin; 469 ierr = PetscObjectGetComm((PetscObject)boundary,&comm);CHKERRQ(ierr); 470 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 471 ierr = DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd);CHKERRQ(ierr); 472 ierr = DMGetLabel(boundary, labelName, &label);CHKERRQ(ierr); 473 474 in.numberofpoints = vEnd - vStart; 475 if (in.numberofpoints > 0) { 476 PetscSection coordSection; 477 Vec coordinates; 478 PetscScalar *array; 479 480 in.pointlist = new double[in.numberofpoints*dim]; 481 in.pointmarkerlist = new int[in.numberofpoints]; 482 483 ierr = DMGetCoordinatesLocal(boundary, &coordinates);CHKERRQ(ierr); 484 ierr = DMGetCoordinateSection(boundary, &coordSection);CHKERRQ(ierr); 485 ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr); 486 for (v = vStart; v < vEnd; ++v) { 487 const PetscInt idx = v - vStart; 488 PetscInt off, d; 489 490 ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 491 for (d = 0; d < dim; ++d) in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]); 492 if (label) {ierr = DMLabelGetValue(label, v, &in.pointmarkerlist[idx]);CHKERRQ(ierr);} 493 } 494 ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr); 495 } 496 ierr = DMPlexGetHeightStratum(boundary, 0, &fStart, &fEnd);CHKERRQ(ierr); 497 498 in.numberoffacets = fEnd - fStart; 499 if (in.numberoffacets > 0) { 500 in.facetlist = new tetgenio::facet[in.numberoffacets]; 501 in.facetmarkerlist = new int[in.numberoffacets]; 502 for (f = fStart; f < fEnd; ++f) { 503 const PetscInt idx = f - fStart; 504 PetscInt *points = NULL, numPoints, p, numVertices = 0, v; 505 506 in.facetlist[idx].numberofpolygons = 1; 507 in.facetlist[idx].polygonlist = new tetgenio::polygon[in.facetlist[idx].numberofpolygons]; 508 in.facetlist[idx].numberofholes = 0; 509 in.facetlist[idx].holelist = NULL; 510 511 ierr = DMPlexGetTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr); 512 for (p = 0; p < numPoints*2; p += 2) { 513 const PetscInt point = points[p]; 514 if ((point >= vStart) && (point < vEnd)) points[numVertices++] = point; 515 } 516 517 tetgenio::polygon *poly = in.facetlist[idx].polygonlist; 518 poly->numberofvertices = numVertices; 519 poly->vertexlist = new int[poly->numberofvertices]; 520 for (v = 0; v < numVertices; ++v) { 521 const PetscInt vIdx = points[v] - vStart; 522 poly->vertexlist[v] = vIdx; 523 } 524 if (label) {ierr = DMLabelGetValue(label, f, &in.facetmarkerlist[idx]);CHKERRQ(ierr);} 525 ierr = DMPlexRestoreTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr); 526 } 527 } 528 if (!rank) { 529 char args[32]; 530 531 /* Take away 'Q' for verbose output */ 532 ierr = PetscStrcpy(args, "pqezQ");CHKERRQ(ierr); 533 if (mesh->tetgenOpts) {::tetrahedralize(mesh->tetgenOpts, &in, &out);} 534 else {::tetrahedralize(args, &in, &out);} 535 } 536 { 537 DMLabel glabel = NULL; 538 const PetscInt numCorners = 4; 539 const PetscInt numCells = out.numberoftetrahedra; 540 const PetscInt numVertices = out.numberofpoints; 541 const double *meshCoords = out.pointlist; 542 int *cells = out.tetrahedronlist; 543 544 ierr = DMPlexInvertCells_Internal(dim, numCells, numCorners, cells);CHKERRQ(ierr); 545 ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm);CHKERRQ(ierr); 546 if (label) {ierr = DMCreateLabel(*dm, labelName); ierr = DMGetLabel(*dm, labelName, &glabel);} 547 /* Set labels */ 548 for (v = 0; v < numVertices; ++v) { 549 if (out.pointmarkerlist[v]) { 550 if (glabel) {ierr = DMLabelSetValue(glabel, v+numCells, out.pointmarkerlist[v]);CHKERRQ(ierr);} 551 } 552 } 553 if (interpolate) { 554 #if 0 555 PetscInt e; 556 557 /* This check is never actually executed for ctetgen (which never returns edgemarkers) and seems to be broken for 558 * tetgen */ 559 for (e = 0; e < out.numberofedges; e++) { 560 if (out.edgemarkerlist[e]) { 561 const PetscInt vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells}; 562 const PetscInt *edges; 563 PetscInt numEdges; 564 565 ierr = DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 566 if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges); 567 if (glabel) {ierr = DMLabelSetValue(glabel, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);} 568 ierr = DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 569 } 570 } 571 #endif 572 for (f = 0; f < out.numberoftrifaces; f++) { 573 if (out.trifacemarkerlist[f]) { 574 const PetscInt vertices[3] = {out.trifacelist[f*3+0]+numCells, out.trifacelist[f*3+1]+numCells, out.trifacelist[f*3+2]+numCells}; 575 const PetscInt *faces; 576 PetscInt numFaces; 577 578 ierr = DMPlexGetFullJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr); 579 if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces); 580 if (glabel) {ierr = DMLabelSetValue(glabel, faces[0], out.trifacemarkerlist[f]);CHKERRQ(ierr);} 581 ierr = DMPlexRestoreJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr); 582 } 583 } 584 } 585 ierr = DMPlexSetRefinementUniform(*dm, PETSC_FALSE);CHKERRQ(ierr); 586 } 587 PetscFunctionReturn(0); 588 } 589 590 static PetscErrorCode DMPlexRefine_Tetgen(DM dm, double *maxVolumes, DM *dmRefined) 591 { 592 MPI_Comm comm; 593 const PetscInt dim = 3; 594 const char *labelName = "marker"; 595 ::tetgenio in; 596 ::tetgenio out; 597 DMLabel label; 598 PetscInt vStart, vEnd, v, cStart, cEnd, c, depth, depthGlobal; 599 PetscMPIInt rank; 600 PetscErrorCode ierr; 601 602 PetscFunctionBegin; 603 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 604 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 605 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 606 ierr = MPIU_Allreduce(&depth, &depthGlobal, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr); 607 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 608 ierr = DMGetLabel(dm, labelName, &label);CHKERRQ(ierr); 609 610 in.numberofpoints = vEnd - vStart; 611 if (in.numberofpoints > 0) { 612 PetscSection coordSection; 613 Vec coordinates; 614 PetscScalar *array; 615 616 in.pointlist = new double[in.numberofpoints*dim]; 617 in.pointmarkerlist = new int[in.numberofpoints]; 618 619 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 620 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 621 ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr); 622 for (v = vStart; v < vEnd; ++v) { 623 const PetscInt idx = v - vStart; 624 PetscInt off, d; 625 626 ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 627 for (d = 0; d < dim; ++d) in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]); 628 if (label) {ierr = DMLabelGetValue(label, v, &in.pointmarkerlist[idx]);CHKERRQ(ierr);} 629 } 630 ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr); 631 } 632 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 633 634 in.numberofcorners = 4; 635 in.numberoftetrahedra = cEnd - cStart; 636 in.tetrahedronvolumelist = (double*) maxVolumes; 637 if (in.numberoftetrahedra > 0) { 638 in.tetrahedronlist = new int[in.numberoftetrahedra*in.numberofcorners]; 639 for (c = cStart; c < cEnd; ++c) { 640 const PetscInt idx = c - cStart; 641 PetscInt *closure = NULL; 642 PetscInt closureSize; 643 644 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 645 if ((closureSize != 5) && (closureSize != 15)) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Mesh has cell which is not a tetrahedron, %D vertices in closure", closureSize); 646 for (v = 0; v < 4; ++v) { 647 in.tetrahedronlist[idx*in.numberofcorners + v] = closure[(v+closureSize-4)*2] - vStart; 648 } 649 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 650 } 651 } 652 /* TODO: Put in boundary faces with markers */ 653 if (!rank) { 654 char args[32]; 655 656 /* Take away 'Q' for verbose output */ 657 /*ierr = PetscStrcpy(args, "qezQra");CHKERRQ(ierr); */ 658 ierr = PetscStrcpy(args, "qezraVVVV");CHKERRQ(ierr); 659 ::tetrahedralize(args, &in, &out); 660 } 661 in.tetrahedronvolumelist = NULL; 662 663 { 664 DMLabel rlabel = NULL; 665 const PetscInt numCorners = 4; 666 const PetscInt numCells = out.numberoftetrahedra; 667 const PetscInt numVertices = out.numberofpoints; 668 const double *meshCoords = out.pointlist; 669 int *cells = out.tetrahedronlist; 670 671 PetscBool interpolate = depthGlobal > 1 ? PETSC_TRUE : PETSC_FALSE; 672 673 ierr = DMPlexInvertCells_Internal(dim, numCells, numCorners, cells);CHKERRQ(ierr); 674 ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined);CHKERRQ(ierr); 675 if (label) {ierr = DMCreateLabel(*dmRefined, labelName); ierr = DMGetLabel(*dmRefined, labelName, &rlabel);} 676 /* Set labels */ 677 for (v = 0; v < numVertices; ++v) { 678 if (out.pointmarkerlist[v]) { 679 if (rlabel) {ierr = DMLabelSetValue(rlabel, v+numCells, out.pointmarkerlist[v]);CHKERRQ(ierr);} 680 } 681 } 682 if (interpolate) { 683 PetscInt e, f; 684 685 for (e = 0; e < out.numberofedges; e++) { 686 if (out.edgemarkerlist[e]) { 687 const PetscInt vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells}; 688 const PetscInt *edges; 689 PetscInt numEdges; 690 691 ierr = DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 692 if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges); 693 if (rlabel) {ierr = DMLabelSetValue(rlabel, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);} 694 ierr = DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 695 } 696 } 697 for (f = 0; f < out.numberoftrifaces; f++) { 698 if (out.trifacemarkerlist[f]) { 699 const PetscInt vertices[3] = {out.trifacelist[f*3+0]+numCells, out.trifacelist[f*3+1]+numCells, out.trifacelist[f*3+2]+numCells}; 700 const PetscInt *faces; 701 PetscInt numFaces; 702 703 ierr = DMPlexGetJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr); 704 if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces); 705 if (rlabel) {ierr = DMLabelSetValue(rlabel, faces[0], out.trifacemarkerlist[f]);CHKERRQ(ierr);} 706 ierr = DMPlexRestoreJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr); 707 } 708 } 709 } 710 ierr = DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE);CHKERRQ(ierr); 711 } 712 PetscFunctionReturn(0); 713 } 714 #endif /* PETSC_HAVE_TETGEN */ 715 716 #if defined(PETSC_HAVE_CTETGEN) 717 #include <ctetgen.h> 718 719 static PetscErrorCode DMPlexGenerate_CTetgen(DM boundary, PetscBool interpolate, DM *dm) 720 { 721 MPI_Comm comm; 722 const PetscInt dim = 3; 723 const char *labelName = "marker"; 724 PLC *in, *out; 725 DMLabel label; 726 PetscInt verbose = 0, vStart, vEnd, v, fStart, fEnd, f; 727 PetscMPIInt rank; 728 PetscErrorCode ierr; 729 730 PetscFunctionBegin; 731 ierr = PetscObjectGetComm((PetscObject)boundary,&comm);CHKERRQ(ierr); 732 ierr = PetscOptionsGetInt(NULL,((PetscObject) boundary)->prefix, "-ctetgen_verbose", &verbose, NULL);CHKERRQ(ierr); 733 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 734 ierr = DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd);CHKERRQ(ierr); 735 ierr = DMGetLabel(boundary, labelName, &label);CHKERRQ(ierr); 736 ierr = PLCCreate(&in);CHKERRQ(ierr); 737 ierr = PLCCreate(&out);CHKERRQ(ierr); 738 739 in->numberofpoints = vEnd - vStart; 740 if (in->numberofpoints > 0) { 741 PetscSection coordSection; 742 Vec coordinates; 743 PetscScalar *array; 744 745 ierr = PetscMalloc1(in->numberofpoints*dim, &in->pointlist);CHKERRQ(ierr); 746 ierr = PetscMalloc1(in->numberofpoints, &in->pointmarkerlist);CHKERRQ(ierr); 747 ierr = DMGetCoordinatesLocal(boundary, &coordinates);CHKERRQ(ierr); 748 ierr = DMGetCoordinateSection(boundary, &coordSection);CHKERRQ(ierr); 749 ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr); 750 for (v = vStart; v < vEnd; ++v) { 751 const PetscInt idx = v - vStart; 752 PetscInt off, d, m = -1; 753 754 ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 755 for (d = 0; d < dim; ++d) { 756 in->pointlist[idx*dim + d] = PetscRealPart(array[off+d]); 757 } 758 if (label) {ierr = DMLabelGetValue(label, v, &m);CHKERRQ(ierr);} 759 760 in->pointmarkerlist[idx] = (int) m; 761 } 762 ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr); 763 } 764 ierr = DMPlexGetHeightStratum(boundary, 0, &fStart, &fEnd);CHKERRQ(ierr); 765 766 in->numberoffacets = fEnd - fStart; 767 if (in->numberoffacets > 0) { 768 ierr = PetscMalloc1(in->numberoffacets, &in->facetlist);CHKERRQ(ierr); 769 ierr = PetscMalloc1(in->numberoffacets, &in->facetmarkerlist);CHKERRQ(ierr); 770 for (f = fStart; f < fEnd; ++f) { 771 const PetscInt idx = f - fStart; 772 PetscInt *points = NULL, numPoints, p, numVertices = 0, v, m = -1; 773 polygon *poly; 774 775 in->facetlist[idx].numberofpolygons = 1; 776 777 ierr = PetscMalloc1(in->facetlist[idx].numberofpolygons, &in->facetlist[idx].polygonlist);CHKERRQ(ierr); 778 779 in->facetlist[idx].numberofholes = 0; 780 in->facetlist[idx].holelist = NULL; 781 782 ierr = DMPlexGetTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr); 783 for (p = 0; p < numPoints*2; p += 2) { 784 const PetscInt point = points[p]; 785 if ((point >= vStart) && (point < vEnd)) points[numVertices++] = point; 786 } 787 788 poly = in->facetlist[idx].polygonlist; 789 poly->numberofvertices = numVertices; 790 ierr = PetscMalloc1(poly->numberofvertices, &poly->vertexlist);CHKERRQ(ierr); 791 for (v = 0; v < numVertices; ++v) { 792 const PetscInt vIdx = points[v] - vStart; 793 poly->vertexlist[v] = vIdx; 794 } 795 if (label) {ierr = DMLabelGetValue(label, f, &m);CHKERRQ(ierr);} 796 in->facetmarkerlist[idx] = (int) m; 797 ierr = DMPlexRestoreTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr); 798 } 799 } 800 if (!rank) { 801 TetGenOpts t; 802 803 ierr = TetGenOptsInitialize(&t);CHKERRQ(ierr); 804 t.in = boundary; /* Should go away */ 805 t.plc = 1; 806 t.quality = 1; 807 t.edgesout = 1; 808 t.zeroindex = 1; 809 t.quiet = 1; 810 t.verbose = verbose; 811 ierr = TetGenCheckOpts(&t);CHKERRQ(ierr); 812 ierr = TetGenTetrahedralize(&t, in, out);CHKERRQ(ierr); 813 } 814 { 815 DMLabel glabel = NULL; 816 const PetscInt numCorners = 4; 817 const PetscInt numCells = out->numberoftetrahedra; 818 const PetscInt numVertices = out->numberofpoints; 819 double *meshCoords; 820 int *cells = out->tetrahedronlist; 821 822 if (sizeof (PetscReal) == sizeof (double)) { 823 meshCoords = (double *) out->pointlist; 824 } 825 else { 826 PetscInt i; 827 828 ierr = PetscMalloc1(3 * numVertices,&meshCoords);CHKERRQ(ierr); 829 for (i = 0; i < 3 * numVertices; i++) { 830 meshCoords[i] = (PetscReal) out->pointlist[i]; 831 } 832 } 833 834 ierr = DMPlexInvertCells_Internal(dim, numCells, numCorners, cells);CHKERRQ(ierr); 835 ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm);CHKERRQ(ierr); 836 if (sizeof (PetscReal) != sizeof (double)) { 837 ierr = PetscFree(meshCoords);CHKERRQ(ierr); 838 } 839 if (label) {ierr = DMCreateLabel(*dm, labelName); ierr = DMGetLabel(*dm, labelName, &glabel);} 840 /* Set labels */ 841 for (v = 0; v < numVertices; ++v) { 842 if (out->pointmarkerlist[v]) { 843 if (glabel) {ierr = DMLabelSetValue(glabel, v+numCells, out->pointmarkerlist[v]);CHKERRQ(ierr);} 844 } 845 } 846 if (interpolate) { 847 PetscInt e; 848 849 for (e = 0; e < out->numberofedges; e++) { 850 if (out->edgemarkerlist[e]) { 851 const PetscInt vertices[2] = {out->edgelist[e*2+0]+numCells, out->edgelist[e*2+1]+numCells}; 852 const PetscInt *edges; 853 PetscInt numEdges; 854 855 ierr = DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 856 if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges); 857 if (glabel) {ierr = DMLabelSetValue(glabel, edges[0], out->edgemarkerlist[e]);CHKERRQ(ierr);} 858 ierr = DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 859 } 860 } 861 for (f = 0; f < out->numberoftrifaces; f++) { 862 if (out->trifacemarkerlist[f]) { 863 const PetscInt vertices[3] = {out->trifacelist[f*3+0]+numCells, out->trifacelist[f*3+1]+numCells, out->trifacelist[f*3+2]+numCells}; 864 const PetscInt *faces; 865 PetscInt numFaces; 866 867 ierr = DMPlexGetFullJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr); 868 if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces); 869 if (glabel) {ierr = DMLabelSetValue(glabel, faces[0], out->trifacemarkerlist[f]);CHKERRQ(ierr);} 870 ierr = DMPlexRestoreJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr); 871 } 872 } 873 } 874 ierr = DMPlexSetRefinementUniform(*dm, PETSC_FALSE);CHKERRQ(ierr); 875 } 876 877 ierr = PLCDestroy(&in);CHKERRQ(ierr); 878 ierr = PLCDestroy(&out);CHKERRQ(ierr); 879 PetscFunctionReturn(0); 880 } 881 882 static PetscErrorCode DMPlexRefine_CTetgen(DM dm, PetscReal *maxVolumes, DM *dmRefined) 883 { 884 MPI_Comm comm; 885 const PetscInt dim = 3; 886 const char *labelName = "marker"; 887 PLC *in, *out; 888 DMLabel label; 889 PetscInt verbose = 0, vStart, vEnd, v, cStart, cEnd, c, depth, depthGlobal; 890 PetscMPIInt rank; 891 PetscErrorCode ierr; 892 893 PetscFunctionBegin; 894 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 895 ierr = PetscOptionsGetInt(NULL,((PetscObject) dm)->prefix, "-ctetgen_verbose", &verbose, NULL);CHKERRQ(ierr); 896 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 897 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 898 ierr = MPIU_Allreduce(&depth, &depthGlobal, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr); 899 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 900 ierr = DMGetLabel(dm, labelName, &label);CHKERRQ(ierr); 901 ierr = PLCCreate(&in);CHKERRQ(ierr); 902 ierr = PLCCreate(&out);CHKERRQ(ierr); 903 904 in->numberofpoints = vEnd - vStart; 905 if (in->numberofpoints > 0) { 906 PetscSection coordSection; 907 Vec coordinates; 908 PetscScalar *array; 909 910 ierr = PetscMalloc1(in->numberofpoints*dim, &in->pointlist);CHKERRQ(ierr); 911 ierr = PetscMalloc1(in->numberofpoints, &in->pointmarkerlist);CHKERRQ(ierr); 912 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 913 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 914 ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr); 915 for (v = vStart; v < vEnd; ++v) { 916 const PetscInt idx = v - vStart; 917 PetscInt off, d, m = -1; 918 919 ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 920 for (d = 0; d < dim; ++d) { 921 in->pointlist[idx*dim + d] = PetscRealPart(array[off+d]); 922 } 923 if (label) {ierr = DMLabelGetValue(label, v, &m);CHKERRQ(ierr);} 924 925 in->pointmarkerlist[idx] = (int) m; 926 } 927 ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr); 928 } 929 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 930 931 in->numberofcorners = 4; 932 in->numberoftetrahedra = cEnd - cStart; 933 in->tetrahedronvolumelist = maxVolumes; 934 if (in->numberoftetrahedra > 0) { 935 ierr = PetscMalloc1(in->numberoftetrahedra*in->numberofcorners, &in->tetrahedronlist);CHKERRQ(ierr); 936 for (c = cStart; c < cEnd; ++c) { 937 const PetscInt idx = c - cStart; 938 PetscInt *closure = NULL; 939 PetscInt closureSize; 940 941 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 942 if ((closureSize != 5) && (closureSize != 15)) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Mesh has cell which is not a tetrahedron, %D vertices in closure", closureSize); 943 for (v = 0; v < 4; ++v) { 944 in->tetrahedronlist[idx*in->numberofcorners + v] = closure[(v+closureSize-4)*2] - vStart; 945 } 946 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 947 } 948 } 949 if (!rank) { 950 TetGenOpts t; 951 952 ierr = TetGenOptsInitialize(&t);CHKERRQ(ierr); 953 954 t.in = dm; /* Should go away */ 955 t.refine = 1; 956 t.varvolume = 1; 957 t.quality = 1; 958 t.edgesout = 1; 959 t.zeroindex = 1; 960 t.quiet = 1; 961 t.verbose = verbose; /* Change this */ 962 963 ierr = TetGenCheckOpts(&t);CHKERRQ(ierr); 964 ierr = TetGenTetrahedralize(&t, in, out);CHKERRQ(ierr); 965 } 966 { 967 DMLabel rlabel = NULL; 968 const PetscInt numCorners = 4; 969 const PetscInt numCells = out->numberoftetrahedra; 970 const PetscInt numVertices = out->numberofpoints; 971 double *meshCoords; 972 int *cells = out->tetrahedronlist; 973 PetscBool interpolate = depthGlobal > 1 ? PETSC_TRUE : PETSC_FALSE; 974 975 if (sizeof (PetscReal) == sizeof (double)) { 976 meshCoords = (double *) out->pointlist; 977 } 978 else { 979 PetscInt i; 980 981 ierr = PetscMalloc1(3 * numVertices,&meshCoords);CHKERRQ(ierr); 982 for (i = 0; i < 3 * numVertices; i++) { 983 meshCoords[i] = (PetscReal) out->pointlist[i]; 984 } 985 } 986 987 ierr = DMPlexInvertCells_Internal(dim, numCells, numCorners, cells);CHKERRQ(ierr); 988 ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined);CHKERRQ(ierr); 989 if (sizeof (PetscReal) != sizeof (double)) { 990 ierr = PetscFree(meshCoords);CHKERRQ(ierr); 991 } 992 if (label) {ierr = DMCreateLabel(*dmRefined, labelName); ierr = DMGetLabel(*dmRefined, labelName, &rlabel);} 993 /* Set labels */ 994 for (v = 0; v < numVertices; ++v) { 995 if (out->pointmarkerlist[v]) { 996 if (rlabel) {ierr = DMLabelSetValue(rlabel, v+numCells, out->pointmarkerlist[v]);CHKERRQ(ierr);} 997 } 998 } 999 if (interpolate) { 1000 PetscInt e, f; 1001 1002 for (e = 0; e < out->numberofedges; e++) { 1003 if (out->edgemarkerlist[e]) { 1004 const PetscInt vertices[2] = {out->edgelist[e*2+0]+numCells, out->edgelist[e*2+1]+numCells}; 1005 const PetscInt *edges; 1006 PetscInt numEdges; 1007 1008 ierr = DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 1009 if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges); 1010 if (rlabel) {ierr = DMLabelSetValue(rlabel, edges[0], out->edgemarkerlist[e]);CHKERRQ(ierr);} 1011 ierr = DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 1012 } 1013 } 1014 for (f = 0; f < out->numberoftrifaces; f++) { 1015 if (out->trifacemarkerlist[f]) { 1016 const PetscInt vertices[3] = {out->trifacelist[f*3+0]+numCells, out->trifacelist[f*3+1]+numCells, out->trifacelist[f*3+2]+numCells}; 1017 const PetscInt *faces; 1018 PetscInt numFaces; 1019 1020 ierr = DMPlexGetFullJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr); 1021 if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces); 1022 if (rlabel) {ierr = DMLabelSetValue(rlabel, faces[0], out->trifacemarkerlist[f]);CHKERRQ(ierr);} 1023 ierr = DMPlexRestoreJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr); 1024 } 1025 } 1026 } 1027 ierr = DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE);CHKERRQ(ierr); 1028 } 1029 ierr = PLCDestroy(&in);CHKERRQ(ierr); 1030 ierr = PLCDestroy(&out);CHKERRQ(ierr); 1031 PetscFunctionReturn(0); 1032 } 1033 #endif /* PETSC_HAVE_CTETGEN */ 1034 1035 /*@C 1036 DMPlexGenerate - Generates a mesh. 1037 1038 Not Collective 1039 1040 Input Parameters: 1041 + boundary - The DMPlex boundary object 1042 . name - The mesh generation package name 1043 - interpolate - Flag to create intermediate mesh elements 1044 1045 Output Parameter: 1046 . mesh - The DMPlex object 1047 1048 Level: intermediate 1049 1050 .keywords: mesh, elements 1051 .seealso: DMPlexCreate(), DMRefine() 1052 @*/ 1053 PetscErrorCode DMPlexGenerate(DM boundary, const char name[], PetscBool interpolate, DM *mesh) 1054 { 1055 PetscInt dim; 1056 char genname[1024]; 1057 PetscBool isTriangle = PETSC_FALSE, isTetgen = PETSC_FALSE, isCTetgen = PETSC_FALSE, flg; 1058 PetscErrorCode ierr; 1059 1060 PetscFunctionBegin; 1061 PetscValidHeaderSpecific(boundary, DM_CLASSID, 1); 1062 PetscValidLogicalCollectiveBool(boundary, interpolate, 2); 1063 ierr = DMGetDimension(boundary, &dim);CHKERRQ(ierr); 1064 ierr = PetscOptionsGetString(((PetscObject) boundary)->options,((PetscObject) boundary)->prefix, "-dm_plex_generator", genname, 1024, &flg);CHKERRQ(ierr); 1065 if (flg) name = genname; 1066 if (name) { 1067 ierr = PetscStrcmp(name, "triangle", &isTriangle);CHKERRQ(ierr); 1068 ierr = PetscStrcmp(name, "tetgen", &isTetgen);CHKERRQ(ierr); 1069 ierr = PetscStrcmp(name, "ctetgen", &isCTetgen);CHKERRQ(ierr); 1070 } 1071 switch (dim) { 1072 case 1: 1073 if (!name || isTriangle) { 1074 #if defined(PETSC_HAVE_TRIANGLE) 1075 ierr = DMPlexGenerate_Triangle(boundary, interpolate, mesh);CHKERRQ(ierr); 1076 #else 1077 SETERRQ(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "Mesh generation needs external package support.\nPlease reconfigure with --download-triangle."); 1078 #endif 1079 } else SETERRQ1(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "Unknown 2D mesh generation package %s", name); 1080 break; 1081 case 2: 1082 if (!name) { 1083 #if defined(PETSC_HAVE_CTETGEN) 1084 ierr = DMPlexGenerate_CTetgen(boundary, interpolate, mesh);CHKERRQ(ierr); 1085 #elif defined(PETSC_HAVE_TETGEN) 1086 ierr = DMPlexGenerate_Tetgen(boundary, interpolate, mesh);CHKERRQ(ierr); 1087 #else 1088 SETERRQ(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "External package CTetgen or Tetgen needed.\nPlease reconfigure with --download-ctetgen or --download-tetgen."); 1089 #endif 1090 } else if (isCTetgen) { 1091 #if defined(PETSC_HAVE_CTETGEN) 1092 ierr = DMPlexGenerate_CTetgen(boundary, interpolate, mesh);CHKERRQ(ierr); 1093 #else 1094 SETERRQ(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "CTetgen needs external package support.\nPlease reconfigure with --download-ctetgen."); 1095 #endif 1096 } else if (isTetgen) { 1097 #if defined(PETSC_HAVE_TETGEN) 1098 ierr = DMPlexGenerate_Tetgen(boundary, interpolate, mesh);CHKERRQ(ierr); 1099 #else 1100 SETERRQ(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "Tetgen needs external package support.\nPlease reconfigure with --download-tetgen."); 1101 #endif 1102 } else SETERRQ1(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "Unknown 3D mesh generation package %s", name); 1103 break; 1104 default: 1105 SETERRQ1(PetscObjectComm((PetscObject)boundary), PETSC_ERR_SUP, "Mesh generation for a dimension %d boundary is not supported.", dim); 1106 } 1107 PetscFunctionReturn(0); 1108 } 1109 1110 #if defined(PETSC_HAVE_TRIANGLE) || defined(PETSC_HAVE_CTETGEN) || defined(PETSC_HAVE_TETGEN) 1111 static PetscErrorCode DMRefine_Plex_Label(DM dm, DMLabel adaptLabel, PetscInt cStart, PetscInt cEnd, PetscReal maxVolumes[]) 1112 { 1113 PetscInt dim, c; 1114 PetscReal refRatio; 1115 PetscErrorCode ierr; 1116 1117 PetscFunctionBegin; 1118 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 1119 refRatio = (PetscReal) ((PetscInt) 1 << dim); 1120 for (c = cStart; c < cEnd; c++) { 1121 PetscReal vol; 1122 PetscInt i, closureSize = 0; 1123 PetscInt *closure = NULL; 1124 PetscBool anyRefine = PETSC_FALSE; 1125 PetscBool anyCoarsen = PETSC_FALSE; 1126 PetscBool anyKeep = PETSC_FALSE; 1127 1128 ierr = DMPlexComputeCellGeometryFVM(dm,c,&vol,NULL,NULL);CHKERRQ(ierr); 1129 maxVolumes[c - cStart] = vol; 1130 ierr = DMPlexGetTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 1131 for (i = 0; i < closureSize; i++) { 1132 PetscInt point = closure[2 * i], refFlag; 1133 1134 ierr = DMLabelGetValue(adaptLabel,point,&refFlag);CHKERRQ(ierr); 1135 switch (refFlag) { 1136 case DM_ADAPT_REFINE: 1137 anyRefine = PETSC_TRUE; 1138 break; 1139 case DM_ADAPT_COARSEN: 1140 anyCoarsen = PETSC_TRUE; 1141 break; 1142 case DM_ADAPT_KEEP: 1143 anyKeep = PETSC_TRUE; 1144 break; 1145 case DM_ADAPT_DETERMINE: 1146 break; 1147 default: 1148 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"DMPlex does not support refinement flag %D\n",refFlag); 1149 break; 1150 } 1151 if (anyRefine) break; 1152 } 1153 ierr = DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 1154 if (anyRefine) { 1155 maxVolumes[c - cStart] = vol / refRatio; 1156 } else if (anyKeep) { 1157 maxVolumes[c - cStart] = vol; 1158 } else if (anyCoarsen) { 1159 maxVolumes[c - cStart] = vol * refRatio; 1160 } 1161 } 1162 PetscFunctionReturn(0); 1163 } 1164 #endif 1165 1166 static PetscErrorCode DMRefine_Plex_Internal(DM dm, MPI_Comm comm, DMLabel adaptLabel, DM *dmRefined) 1167 { 1168 PetscErrorCode (*refinementFunc)(const PetscReal [], PetscReal *); 1169 PetscReal refinementLimit; 1170 PetscInt dim, cStart, cEnd; 1171 char genname[1024], *name = NULL; 1172 PetscBool isUniform, isTriangle = PETSC_FALSE, isTetgen = PETSC_FALSE, isCTetgen = PETSC_FALSE, flg, localized; 1173 PetscErrorCode ierr; 1174 1175 PetscFunctionBegin; 1176 ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr); 1177 ierr = DMPlexGetRefinementUniform(dm, &isUniform);CHKERRQ(ierr); 1178 if (isUniform) { 1179 CellRefiner cellRefiner; 1180 1181 ierr = DMPlexGetCellRefiner_Internal(dm, &cellRefiner);CHKERRQ(ierr); 1182 ierr = DMPlexRefineUniform_Internal(dm, cellRefiner, dmRefined);CHKERRQ(ierr); 1183 ierr = DMCopyBoundary(dm, *dmRefined);CHKERRQ(ierr); 1184 if (localized) { 1185 ierr = DMLocalizeCoordinates(*dmRefined);CHKERRQ(ierr); 1186 } 1187 PetscFunctionReturn(0); 1188 } 1189 ierr = DMPlexGetRefinementLimit(dm, &refinementLimit);CHKERRQ(ierr); 1190 ierr = DMPlexGetRefinementFunction(dm, &refinementFunc);CHKERRQ(ierr); 1191 if (refinementLimit == 0.0 && !refinementFunc) PetscFunctionReturn(0); 1192 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1193 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1194 ierr = PetscOptionsGetString(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_generator", genname, 1024, &flg);CHKERRQ(ierr); 1195 if (flg) name = genname; 1196 if (name) { 1197 ierr = PetscStrcmp(name, "triangle", &isTriangle);CHKERRQ(ierr); 1198 ierr = PetscStrcmp(name, "tetgen", &isTetgen);CHKERRQ(ierr); 1199 ierr = PetscStrcmp(name, "ctetgen", &isCTetgen);CHKERRQ(ierr); 1200 } 1201 switch (dim) { 1202 case 2: 1203 if (!name || isTriangle) { 1204 #if defined(PETSC_HAVE_TRIANGLE) 1205 PetscReal *maxVolumes; 1206 PetscInt c; 1207 1208 ierr = PetscMalloc1(cEnd - cStart, &maxVolumes);CHKERRQ(ierr); 1209 if (adaptLabel) { 1210 ierr = DMRefine_Plex_Label(dm,adaptLabel,cStart,cEnd,maxVolumes);CHKERRQ(ierr); 1211 } else if (refinementFunc) { 1212 for (c = cStart; c < cEnd; ++c) { 1213 PetscReal vol, centroid[3]; 1214 PetscReal maxVol; 1215 1216 ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL);CHKERRQ(ierr); 1217 ierr = (*refinementFunc)(centroid, &maxVol);CHKERRQ(ierr); 1218 maxVolumes[c - cStart] = (double) maxVol; 1219 } 1220 } else { 1221 for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit; 1222 } 1223 #if !defined(PETSC_USE_REAL_DOUBLE) 1224 { 1225 double *mvols; 1226 ierr = PetscMalloc1(cEnd - cStart,&mvols);CHKERRQ(ierr); 1227 for (c = 0; c < cEnd-cStart; ++c) mvols[c] = (double)maxVolumes[c]; 1228 ierr = DMPlexRefine_Triangle(dm, mvols, dmRefined);CHKERRQ(ierr); 1229 ierr = PetscFree(mvols);CHKERRQ(ierr); 1230 } 1231 #else 1232 ierr = DMPlexRefine_Triangle(dm, maxVolumes, dmRefined);CHKERRQ(ierr); 1233 #endif 1234 ierr = PetscFree(maxVolumes);CHKERRQ(ierr); 1235 #else 1236 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Mesh refinement needs external package support.\nPlease reconfigure with --download-triangle."); 1237 #endif 1238 } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unknown 2D mesh generation package %s", name); 1239 break; 1240 case 3: 1241 if (!name || isCTetgen) { 1242 #if defined(PETSC_HAVE_CTETGEN) 1243 PetscReal *maxVolumes; 1244 PetscInt c; 1245 1246 ierr = PetscMalloc1(cEnd - cStart, &maxVolumes);CHKERRQ(ierr); 1247 if (adaptLabel) { 1248 ierr = DMRefine_Plex_Label(dm,adaptLabel,cStart,cEnd,maxVolumes);CHKERRQ(ierr); 1249 } else if (refinementFunc) { 1250 for (c = cStart; c < cEnd; ++c) { 1251 PetscReal vol, centroid[3]; 1252 1253 ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL);CHKERRQ(ierr); 1254 ierr = (*refinementFunc)(centroid, &maxVolumes[c-cStart]);CHKERRQ(ierr); 1255 } 1256 } else { 1257 for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit; 1258 } 1259 ierr = DMPlexRefine_CTetgen(dm, maxVolumes, dmRefined);CHKERRQ(ierr); 1260 #else 1261 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "CTetgen needs external package support.\nPlease reconfigure with --download-ctetgen."); 1262 #endif 1263 } else if (isTetgen) { 1264 #if defined(PETSC_HAVE_TETGEN) 1265 double *maxVolumes; 1266 PetscInt c; 1267 1268 ierr = PetscMalloc1(cEnd - cStart, &maxVolumes);CHKERRQ(ierr); 1269 if (adaptLabel) { 1270 ierr = DMRefine_Plex_Label(dm,adaptLabel,cStart,cEnd,maxVolumes);CHKERRQ(ierr); 1271 } else if (refinementFunc) { 1272 for (c = cStart; c < cEnd; ++c) { 1273 PetscReal vol, centroid[3]; 1274 1275 ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL);CHKERRQ(ierr); 1276 ierr = (*refinementFunc)(centroid, &maxVolumes[c-cStart]);CHKERRQ(ierr); 1277 } 1278 } else { 1279 for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit; 1280 } 1281 ierr = DMPlexRefine_Tetgen(dm, maxVolumes, dmRefined);CHKERRQ(ierr); 1282 #else 1283 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Tetgen needs external package support.\nPlease reconfigure with --with-c-language=cxx --download-tetgen."); 1284 #endif 1285 } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unknown 3D mesh generation package %s", name); 1286 break; 1287 default: 1288 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Mesh refinement in dimension %d is not supported.", dim); 1289 } 1290 ierr = DMCopyBoundary(dm, *dmRefined);CHKERRQ(ierr); 1291 if (localized) { 1292 ierr = DMLocalizeCoordinates(*dmRefined);CHKERRQ(ierr); 1293 } 1294 PetscFunctionReturn(0); 1295 } 1296 1297 PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *dmRefined) 1298 { 1299 PetscErrorCode ierr; 1300 1301 PetscFunctionBegin; 1302 ierr = DMRefine_Plex_Internal(dm,comm,NULL,dmRefined);CHKERRQ(ierr); 1303 PetscFunctionReturn(0); 1304 } 1305 1306 PetscErrorCode DMAdaptLabel_Plex(DM dm, DMLabel adaptLabel, DM *dmAdapted) 1307 { 1308 PetscInt defFlag, minFlag, maxFlag, numFlags, i; 1309 const PetscInt *flags; 1310 IS flagIS; 1311 PetscErrorCode ierr; 1312 1313 PetscFunctionBegin; 1314 ierr = DMLabelGetDefaultValue(adaptLabel,&defFlag);CHKERRQ(ierr); 1315 minFlag = defFlag; 1316 maxFlag = defFlag; 1317 ierr = DMLabelGetValueIS(adaptLabel,&flagIS);CHKERRQ(ierr); 1318 ierr = ISGetLocalSize(flagIS,&numFlags);CHKERRQ(ierr); 1319 ierr = ISGetIndices(flagIS,&flags);CHKERRQ(ierr); 1320 for (i = 0; i < numFlags; i++) { 1321 PetscInt flag = flags[i]; 1322 1323 minFlag = PetscMin(minFlag,flag); 1324 maxFlag = PetscMax(maxFlag,flag); 1325 } 1326 ierr = ISRestoreIndices(flagIS,&flags);CHKERRQ(ierr); 1327 ierr = ISDestroy(&flagIS);CHKERRQ(ierr); 1328 { 1329 PetscInt minMaxFlag[2], minMaxFlagGlobal[2]; 1330 1331 minMaxFlag[0] = minFlag; 1332 minMaxFlag[1] = -maxFlag; 1333 ierr = MPI_Allreduce(minMaxFlag,minMaxFlagGlobal,2,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 1334 minFlag = minMaxFlagGlobal[0]; 1335 maxFlag = -minMaxFlagGlobal[1]; 1336 } 1337 if (minFlag == maxFlag) { 1338 switch (minFlag) { 1339 case DM_ADAPT_DETERMINE: 1340 *dmAdapted = NULL; 1341 break; 1342 case DM_ADAPT_REFINE: 1343 ierr = DMPlexSetRefinementUniform(dm,PETSC_TRUE);CHKERRQ(ierr); 1344 ierr = DMRefine(dm,MPI_COMM_NULL,dmAdapted);CHKERRQ(ierr); 1345 break; 1346 case DM_ADAPT_COARSEN: 1347 ierr = DMCoarsen(dm,MPI_COMM_NULL,dmAdapted);CHKERRQ(ierr); 1348 break; 1349 default: 1350 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"DMPlex does not support refinement flag %D\n",minFlag); 1351 break; 1352 } 1353 } else { 1354 ierr = DMPlexSetRefinementUniform(dm,PETSC_FALSE);CHKERRQ(ierr); 1355 ierr = DMRefine_Plex_Internal(dm,MPI_COMM_NULL,adaptLabel,dmAdapted);CHKERRQ(ierr); 1356 } 1357 PetscFunctionReturn(0); 1358 } 1359 1360 PetscErrorCode DMRefineHierarchy_Plex(DM dm, PetscInt nlevels, DM dmRefined[]) 1361 { 1362 DM cdm = dm; 1363 PetscInt r; 1364 PetscBool isUniform, localized; 1365 PetscErrorCode ierr; 1366 1367 PetscFunctionBegin; 1368 ierr = DMPlexGetRefinementUniform(dm, &isUniform);CHKERRQ(ierr); 1369 ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr); 1370 if (!isUniform) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Non-uniform refinement is incompatible with the hierarchy"); 1371 for (r = 0; r < nlevels; ++r) { 1372 CellRefiner cellRefiner; 1373 1374 ierr = DMPlexGetCellRefiner_Internal(cdm, &cellRefiner);CHKERRQ(ierr); 1375 ierr = DMPlexRefineUniform_Internal(cdm, cellRefiner, &dmRefined[r]);CHKERRQ(ierr); 1376 ierr = DMCopyBoundary(cdm, dmRefined[r]);CHKERRQ(ierr); 1377 if (localized) { 1378 ierr = DMLocalizeCoordinates(dmRefined[r]);CHKERRQ(ierr); 1379 } 1380 ierr = DMSetCoarseDM(dmRefined[r], cdm);CHKERRQ(ierr); 1381 ierr = DMPlexSetRegularRefinement(dmRefined[r], PETSC_TRUE);CHKERRQ(ierr); 1382 cdm = dmRefined[r]; 1383 } 1384 PetscFunctionReturn(0); 1385 } 1386