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