13a074057SBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 23a074057SBarry Smith 33316697fSBarry Smith #if !defined(ANSI_DECLARATORS) 43316697fSBarry Smith #define ANSI_DECLARATORS 53316697fSBarry Smith #endif 63a074057SBarry Smith #include <triangle.h> 73a074057SBarry Smith 8d71ae5a4SJacob Faibussowitsch static PetscErrorCode InitInput_Triangle(struct triangulateio *inputCtx) 9d71ae5a4SJacob Faibussowitsch { 103a074057SBarry Smith PetscFunctionBegin; 113a074057SBarry Smith inputCtx->numberofpoints = 0; 123a074057SBarry Smith inputCtx->numberofpointattributes = 0; 133a074057SBarry Smith inputCtx->pointlist = NULL; 143a074057SBarry Smith inputCtx->pointattributelist = NULL; 153a074057SBarry Smith inputCtx->pointmarkerlist = NULL; 163a074057SBarry Smith inputCtx->numberofsegments = 0; 173a074057SBarry Smith inputCtx->segmentlist = NULL; 183a074057SBarry Smith inputCtx->segmentmarkerlist = NULL; 193a074057SBarry Smith inputCtx->numberoftriangleattributes = 0; 203a074057SBarry Smith inputCtx->trianglelist = NULL; 213a074057SBarry Smith inputCtx->numberofholes = 0; 223a074057SBarry Smith inputCtx->holelist = NULL; 233a074057SBarry Smith inputCtx->numberofregions = 0; 243a074057SBarry Smith inputCtx->regionlist = NULL; 253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 263a074057SBarry Smith } 273a074057SBarry Smith 28d71ae5a4SJacob Faibussowitsch static PetscErrorCode InitOutput_Triangle(struct triangulateio *outputCtx) 29d71ae5a4SJacob Faibussowitsch { 303a074057SBarry Smith PetscFunctionBegin; 313a074057SBarry Smith outputCtx->numberofpoints = 0; 323a074057SBarry Smith outputCtx->pointlist = NULL; 333a074057SBarry Smith outputCtx->pointattributelist = NULL; 343a074057SBarry Smith outputCtx->pointmarkerlist = NULL; 353a074057SBarry Smith outputCtx->numberoftriangles = 0; 363a074057SBarry Smith outputCtx->trianglelist = NULL; 373a074057SBarry Smith outputCtx->triangleattributelist = NULL; 383a074057SBarry Smith outputCtx->neighborlist = NULL; 393a074057SBarry Smith outputCtx->segmentlist = NULL; 403a074057SBarry Smith outputCtx->segmentmarkerlist = NULL; 413a074057SBarry Smith outputCtx->numberofedges = 0; 423a074057SBarry Smith outputCtx->edgelist = NULL; 433a074057SBarry Smith outputCtx->edgemarkerlist = NULL; 443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 453a074057SBarry Smith } 463a074057SBarry Smith 47d71ae5a4SJacob Faibussowitsch static PetscErrorCode FiniOutput_Triangle(struct triangulateio *outputCtx) 48d71ae5a4SJacob Faibussowitsch { 493a074057SBarry Smith PetscFunctionBegin; 503a074057SBarry Smith free(outputCtx->pointlist); 513a074057SBarry Smith free(outputCtx->pointmarkerlist); 523a074057SBarry Smith free(outputCtx->segmentlist); 533a074057SBarry Smith free(outputCtx->segmentmarkerlist); 543a074057SBarry Smith free(outputCtx->edgelist); 553a074057SBarry Smith free(outputCtx->edgemarkerlist); 563a074057SBarry Smith free(outputCtx->trianglelist); 573a074057SBarry Smith free(outputCtx->neighborlist); 583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 593a074057SBarry Smith } 603a074057SBarry Smith 61d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMPlexGenerate_Triangle(DM boundary, PetscBool interpolate, DM *dm) 62d71ae5a4SJacob Faibussowitsch { 633a074057SBarry Smith MPI_Comm comm; 643a074057SBarry Smith DM_Plex *mesh = (DM_Plex *)boundary->data; 653a074057SBarry Smith PetscInt dim = 2; 663a074057SBarry Smith const PetscBool createConvexHull = PETSC_FALSE; 673a074057SBarry Smith const PetscBool constrained = PETSC_FALSE; 683a074057SBarry Smith const char *labelName = "marker"; 693a074057SBarry Smith const char *labelName2 = "Face Sets"; 703a074057SBarry Smith struct triangulateio in; 713a074057SBarry Smith struct triangulateio out; 723a074057SBarry Smith DMLabel label, label2; 733a074057SBarry Smith PetscInt vStart, vEnd, v, eStart, eEnd, e; 743a074057SBarry Smith PetscMPIInt rank; 753a074057SBarry Smith 763a074057SBarry Smith PetscFunctionBegin; 779566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)boundary, &comm)); 789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 799566063dSJacob Faibussowitsch PetscCall(InitInput_Triangle(&in)); 809566063dSJacob Faibussowitsch PetscCall(InitOutput_Triangle(&out)); 819566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd)); 829566063dSJacob Faibussowitsch PetscCall(DMGetLabel(boundary, labelName, &label)); 839566063dSJacob Faibussowitsch PetscCall(DMGetLabel(boundary, labelName2, &label2)); 843a074057SBarry Smith 853a074057SBarry Smith in.numberofpoints = vEnd - vStart; 863a074057SBarry Smith if (in.numberofpoints > 0) { 873a074057SBarry Smith PetscSection coordSection; 883a074057SBarry Smith Vec coordinates; 893a074057SBarry Smith PetscScalar *array; 903a074057SBarry Smith 919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(in.numberofpoints * dim, &in.pointlist)); 929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(in.numberofpoints, &in.pointmarkerlist)); 939566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(boundary, &coordinates)); 949566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(boundary, &coordSection)); 959566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &array)); 963a074057SBarry Smith for (v = vStart; v < vEnd; ++v) { 973a074057SBarry Smith const PetscInt idx = v - vStart; 98469e3fe5SMatthew G. Knepley PetscInt val, off, d; 993a074057SBarry Smith 1009566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordSection, v, &off)); 101ad540459SPierre Jolivet for (d = 0; d < dim; ++d) in.pointlist[idx * dim + d] = PetscRealPart(array[off + d]); 102469e3fe5SMatthew G. Knepley if (label) { 1039566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(label, v, &val)); 104469e3fe5SMatthew G. Knepley in.pointmarkerlist[idx] = val; 105469e3fe5SMatthew G. Knepley } 1063a074057SBarry Smith } 1079566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &array)); 1083a074057SBarry Smith } 1099566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(boundary, 0, &eStart, &eEnd)); 1103a074057SBarry Smith in.numberofsegments = eEnd - eStart; 1113a074057SBarry Smith if (in.numberofsegments > 0) { 1129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(in.numberofsegments * 2, &in.segmentlist)); 1139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(in.numberofsegments, &in.segmentmarkerlist)); 1143a074057SBarry Smith for (e = eStart; e < eEnd; ++e) { 1153a074057SBarry Smith const PetscInt idx = e - eStart; 1163a074057SBarry Smith const PetscInt *cone; 117469e3fe5SMatthew G. Knepley PetscInt val; 1183a074057SBarry Smith 1199566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(boundary, e, &cone)); 1203a074057SBarry Smith 1213a074057SBarry Smith in.segmentlist[idx * 2 + 0] = cone[0] - vStart; 1223a074057SBarry Smith in.segmentlist[idx * 2 + 1] = cone[1] - vStart; 1233a074057SBarry Smith 124469e3fe5SMatthew G. Knepley if (label) { 1259566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(label, e, &val)); 126469e3fe5SMatthew G. Knepley in.segmentmarkerlist[idx] = val; 127469e3fe5SMatthew G. Knepley } 1283a074057SBarry Smith } 1293a074057SBarry Smith } 1303a074057SBarry Smith #if 0 /* Do not currently support holes */ 1313a074057SBarry Smith PetscReal *holeCoords; 1323a074057SBarry Smith PetscInt h, d; 1333a074057SBarry Smith 1349566063dSJacob Faibussowitsch PetscCall(DMPlexGetHoles(boundary, &in.numberofholes, &holeCords)); 1353a074057SBarry Smith if (in.numberofholes > 0) { 1369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(in.numberofholes*dim, &in.holelist)); 1373a074057SBarry Smith for (h = 0; h < in.numberofholes; ++h) { 1383a074057SBarry Smith for (d = 0; d < dim; ++d) { 1393a074057SBarry Smith in.holelist[h*dim+d] = holeCoords[h*dim+d]; 1403a074057SBarry Smith } 1413a074057SBarry Smith } 1423a074057SBarry Smith } 1433a074057SBarry Smith #endif 144dd400576SPatrick Sanan if (rank == 0) { 1453a074057SBarry Smith char args[32]; 1463a074057SBarry Smith 1473a074057SBarry Smith /* Take away 'Q' for verbose output */ 148c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(args, "pqezQ", sizeof(args))); 149c6a7a370SJeremy L Thompson if (createConvexHull) PetscCall(PetscStrlcat(args, "c", sizeof(args))); 150c6a7a370SJeremy L Thompson if (constrained) PetscCall(PetscStrncpy(args, "zepDQ", sizeof(args))); 1519371c9d4SSatish Balay if (mesh->triangleOpts) { 1529371c9d4SSatish Balay triangulate(mesh->triangleOpts, &in, &out, NULL); 1539371c9d4SSatish Balay } else { 1549371c9d4SSatish Balay triangulate(args, &in, &out, NULL); 1559371c9d4SSatish Balay } 1563a074057SBarry Smith } 1579566063dSJacob Faibussowitsch PetscCall(PetscFree(in.pointlist)); 1589566063dSJacob Faibussowitsch PetscCall(PetscFree(in.pointmarkerlist)); 1599566063dSJacob Faibussowitsch PetscCall(PetscFree(in.segmentlist)); 1609566063dSJacob Faibussowitsch PetscCall(PetscFree(in.segmentmarkerlist)); 1619566063dSJacob Faibussowitsch PetscCall(PetscFree(in.holelist)); 1623a074057SBarry Smith 1633a074057SBarry Smith { 1643a074057SBarry Smith DMLabel glabel = NULL; 1653a074057SBarry Smith DMLabel glabel2 = NULL; 1663a074057SBarry Smith const PetscInt numCorners = 3; 1673a074057SBarry Smith const PetscInt numCells = out.numberoftriangles; 1683a074057SBarry Smith const PetscInt numVertices = out.numberofpoints; 169a4a685f2SJacob Faibussowitsch PetscInt *cells; 170a4a685f2SJacob Faibussowitsch PetscReal *meshCoords; 1713a074057SBarry Smith 172a4a685f2SJacob Faibussowitsch if (sizeof(PetscReal) == sizeof(out.pointlist[0])) { 173a4a685f2SJacob Faibussowitsch meshCoords = (PetscReal *)out.pointlist; 174a4a685f2SJacob Faibussowitsch } else { 175a4a685f2SJacob Faibussowitsch PetscInt i; 176a4a685f2SJacob Faibussowitsch 1779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim * numVertices, &meshCoords)); 178ad540459SPierre Jolivet for (i = 0; i < dim * numVertices; i++) meshCoords[i] = (PetscReal)out.pointlist[i]; 179a4a685f2SJacob Faibussowitsch } 180a4a685f2SJacob Faibussowitsch if (sizeof(PetscInt) == sizeof(out.trianglelist[0])) { 181a4a685f2SJacob Faibussowitsch cells = (PetscInt *)out.trianglelist; 182a4a685f2SJacob Faibussowitsch } else { 183a4a685f2SJacob Faibussowitsch PetscInt i; 184a4a685f2SJacob Faibussowitsch 1859566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numCells * numCorners, &cells)); 186ad540459SPierre Jolivet for (i = 0; i < numCells * numCorners; i++) cells[i] = (PetscInt)out.trianglelist[i]; 187a4a685f2SJacob Faibussowitsch } 1889566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm)); 18948a46eb9SPierre Jolivet if (sizeof(PetscReal) != sizeof(out.pointlist[0])) PetscCall(PetscFree(meshCoords)); 19048a46eb9SPierre Jolivet if (sizeof(PetscInt) != sizeof(out.trianglelist[0])) PetscCall(PetscFree(cells)); 191a4a685f2SJacob Faibussowitsch if (label) { 1929566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(*dm, labelName)); 1939566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dm, labelName, &glabel)); 194a4a685f2SJacob Faibussowitsch } 195a4a685f2SJacob Faibussowitsch if (label2) { 1969566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(*dm, labelName2)); 1979566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dm, labelName2, &glabel2)); 198a4a685f2SJacob Faibussowitsch } 1993a074057SBarry Smith /* Set labels */ 2003a074057SBarry Smith for (v = 0; v < numVertices; ++v) { 2013a074057SBarry Smith if (out.pointmarkerlist[v]) { 2029566063dSJacob Faibussowitsch if (glabel) PetscCall(DMLabelSetValue(glabel, v + numCells, out.pointmarkerlist[v])); 2033a074057SBarry Smith } 2043a074057SBarry Smith } 2053a074057SBarry Smith if (interpolate) { 2063a074057SBarry Smith for (e = 0; e < out.numberofedges; e++) { 2073a074057SBarry Smith if (out.edgemarkerlist[e]) { 2083a074057SBarry Smith const PetscInt vertices[2] = {out.edgelist[e * 2 + 0] + numCells, out.edgelist[e * 2 + 1] + numCells}; 2093a074057SBarry Smith const PetscInt *edges; 2103a074057SBarry Smith PetscInt numEdges; 2113a074057SBarry Smith 2129566063dSJacob Faibussowitsch PetscCall(DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges)); 21363a3b9bcSJacob Faibussowitsch PetscCheck(numEdges == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %" PetscInt_FMT, numEdges); 2149566063dSJacob Faibussowitsch if (glabel) PetscCall(DMLabelSetValue(glabel, edges[0], out.edgemarkerlist[e])); 2159566063dSJacob Faibussowitsch if (glabel2) PetscCall(DMLabelSetValue(glabel2, edges[0], out.edgemarkerlist[e])); 2169566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges)); 2173a074057SBarry Smith } 2183a074057SBarry Smith } 2193a074057SBarry Smith } 2209566063dSJacob Faibussowitsch PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_FALSE)); 2213a074057SBarry Smith } 2223a074057SBarry Smith #if 0 /* Do not currently support holes */ 2239566063dSJacob Faibussowitsch PetscCall(DMPlexCopyHoles(*dm, boundary)); 2243a074057SBarry Smith #endif 2259566063dSJacob Faibussowitsch PetscCall(FiniOutput_Triangle(&out)); 2263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2273a074057SBarry Smith } 2283a074057SBarry Smith 229d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMPlexRefine_Triangle(DM dm, PetscReal *inmaxVolumes, DM *dmRefined) 230d71ae5a4SJacob Faibussowitsch { 2313a074057SBarry Smith MPI_Comm comm; 2323a074057SBarry Smith PetscInt dim = 2; 2333a074057SBarry Smith const char *labelName = "marker"; 2343a074057SBarry Smith struct triangulateio in; 2353a074057SBarry Smith struct triangulateio out; 2363a074057SBarry Smith DMLabel label; 2378457c761SMatthew G. Knepley PetscInt vStart, vEnd, v, gcStart, cStart, cEnd, c, depth, depthGlobal; 2383a074057SBarry Smith PetscMPIInt rank; 239800b850cSBarry Smith double *maxVolumes; 2403a074057SBarry Smith 2413a074057SBarry Smith PetscFunctionBegin; 2429566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 2439566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2449566063dSJacob Faibussowitsch PetscCall(InitInput_Triangle(&in)); 2459566063dSJacob Faibussowitsch PetscCall(InitOutput_Triangle(&out)); 2469566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth)); 2471c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&depth, &depthGlobal, 1, MPIU_INT, MPI_MAX, comm)); 2489566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 2499566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, labelName, &label)); 2503a074057SBarry Smith 2513a074057SBarry Smith in.numberofpoints = vEnd - vStart; 2523a074057SBarry Smith if (in.numberofpoints > 0) { 2533a074057SBarry Smith PetscSection coordSection; 2543a074057SBarry Smith Vec coordinates; 2553a074057SBarry Smith PetscScalar *array; 2563a074057SBarry Smith 2579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(in.numberofpoints * dim, &in.pointlist)); 2589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(in.numberofpoints, &in.pointmarkerlist)); 2599566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 2609566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 2619566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &array)); 2623a074057SBarry Smith for (v = vStart; v < vEnd; ++v) { 2633a074057SBarry Smith const PetscInt idx = v - vStart; 264469e3fe5SMatthew G. Knepley PetscInt off, d, val; 2653a074057SBarry Smith 2669566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordSection, v, &off)); 267ad540459SPierre Jolivet for (d = 0; d < dim; ++d) in.pointlist[idx * dim + d] = PetscRealPart(array[off + d]); 268469e3fe5SMatthew G. Knepley if (label) { 2699566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(label, v, &val)); 270469e3fe5SMatthew G. Knepley in.pointmarkerlist[idx] = val; 271469e3fe5SMatthew G. Knepley } 2723a074057SBarry Smith } 2739566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &array)); 2743a074057SBarry Smith } 2759566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 276*2827ebadSStefano Zampini PetscCall(DMPlexGetCellTypeStratum(dm, DM_POLYTOPE_FV_GHOST, &gcStart, NULL)); 2778457c761SMatthew G. Knepley if (gcStart >= 0) cEnd = gcStart; 2783a074057SBarry Smith 2793a074057SBarry Smith in.numberofcorners = 3; 2803a074057SBarry Smith in.numberoftriangles = cEnd - cStart; 2813a074057SBarry Smith 282800b850cSBarry Smith #if !defined(PETSC_USE_REAL_DOUBLE) 2839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(cEnd - cStart, &maxVolumes)); 284800b850cSBarry Smith for (c = 0; c < cEnd - cStart; ++c) maxVolumes[c] = (double)inmaxVolumes[c]; 285800b850cSBarry Smith #else 286800b850cSBarry Smith maxVolumes = inmaxVolumes; 287800b850cSBarry Smith #endif 288800b850cSBarry Smith 2893a074057SBarry Smith in.trianglearealist = (double *)maxVolumes; 2903a074057SBarry Smith if (in.numberoftriangles > 0) { 2919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(in.numberoftriangles * in.numberofcorners, &in.trianglelist)); 2923a074057SBarry Smith for (c = cStart; c < cEnd; ++c) { 2933a074057SBarry Smith const PetscInt idx = c - cStart; 2943a074057SBarry Smith PetscInt *closure = NULL; 2953a074057SBarry Smith PetscInt closureSize; 2963a074057SBarry Smith 2979566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 29863a3b9bcSJacob Faibussowitsch PetscCheck(!(closureSize != 4) || !(closureSize != 7), comm, PETSC_ERR_ARG_WRONG, "Mesh has cell which is not a triangle, %" PetscInt_FMT " vertices in closure", closureSize); 299ad540459SPierre Jolivet for (v = 0; v < 3; ++v) in.trianglelist[idx * in.numberofcorners + v] = closure[(v + closureSize - 3) * 2] - vStart; 3009566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 3013a074057SBarry Smith } 3023a074057SBarry Smith } 3033a074057SBarry Smith /* TODO: Segment markers are missing on input */ 3043a074057SBarry Smith #if 0 /* Do not currently support holes */ 3053a074057SBarry Smith PetscReal *holeCoords; 3063a074057SBarry Smith PetscInt h, d; 3073a074057SBarry Smith 3089566063dSJacob Faibussowitsch PetscCall(DMPlexGetHoles(boundary, &in.numberofholes, &holeCords)); 3093a074057SBarry Smith if (in.numberofholes > 0) { 3109566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(in.numberofholes*dim, &in.holelist)); 3113a074057SBarry Smith for (h = 0; h < in.numberofholes; ++h) { 3123a074057SBarry Smith for (d = 0; d < dim; ++d) { 3133a074057SBarry Smith in.holelist[h*dim+d] = holeCoords[h*dim+d]; 3143a074057SBarry Smith } 3153a074057SBarry Smith } 3163a074057SBarry Smith } 3173a074057SBarry Smith #endif 318dd400576SPatrick Sanan if (rank == 0) { 3193a074057SBarry Smith char args[32]; 3203a074057SBarry Smith 3213a074057SBarry Smith /* Take away 'Q' for verbose output */ 322c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(args, "pqezQra", sizeof(args))); 3233a074057SBarry Smith triangulate(args, &in, &out, NULL); 3243a074057SBarry Smith } 3259566063dSJacob Faibussowitsch PetscCall(PetscFree(in.pointlist)); 3269566063dSJacob Faibussowitsch PetscCall(PetscFree(in.pointmarkerlist)); 3279566063dSJacob Faibussowitsch PetscCall(PetscFree(in.segmentlist)); 3289566063dSJacob Faibussowitsch PetscCall(PetscFree(in.segmentmarkerlist)); 3299566063dSJacob Faibussowitsch PetscCall(PetscFree(in.trianglelist)); 3303a074057SBarry Smith 3313a074057SBarry Smith { 3323a074057SBarry Smith DMLabel rlabel = NULL; 3333a074057SBarry Smith const PetscInt numCorners = 3; 3343a074057SBarry Smith const PetscInt numCells = out.numberoftriangles; 3353a074057SBarry Smith const PetscInt numVertices = out.numberofpoints; 336a4a685f2SJacob Faibussowitsch PetscInt *cells; 337a4a685f2SJacob Faibussowitsch PetscReal *meshCoords; 3383a074057SBarry Smith PetscBool interpolate = depthGlobal > 1 ? PETSC_TRUE : PETSC_FALSE; 3393a074057SBarry Smith 340a4a685f2SJacob Faibussowitsch if (sizeof(PetscReal) == sizeof(out.pointlist[0])) { 341a4a685f2SJacob Faibussowitsch meshCoords = (PetscReal *)out.pointlist; 342a4a685f2SJacob Faibussowitsch } else { 343a4a685f2SJacob Faibussowitsch PetscInt i; 344a4a685f2SJacob Faibussowitsch 3459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim * numVertices, &meshCoords)); 346ad540459SPierre Jolivet for (i = 0; i < dim * numVertices; i++) meshCoords[i] = (PetscReal)out.pointlist[i]; 347a4a685f2SJacob Faibussowitsch } 348a4a685f2SJacob Faibussowitsch if (sizeof(PetscInt) == sizeof(out.trianglelist[0])) { 349a4a685f2SJacob Faibussowitsch cells = (PetscInt *)out.trianglelist; 350a4a685f2SJacob Faibussowitsch } else { 351a4a685f2SJacob Faibussowitsch PetscInt i; 352a4a685f2SJacob Faibussowitsch 3539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numCells * numCorners, &cells)); 354ad540459SPierre Jolivet for (i = 0; i < numCells * numCorners; i++) cells[i] = (PetscInt)out.trianglelist[i]; 355a4a685f2SJacob Faibussowitsch } 356a4a685f2SJacob Faibussowitsch 3579566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined)); 358a4a685f2SJacob Faibussowitsch if (label) { 3599566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(*dmRefined, labelName)); 3609566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dmRefined, labelName, &rlabel)); 361a4a685f2SJacob Faibussowitsch } 36248a46eb9SPierre Jolivet if (sizeof(PetscReal) != sizeof(out.pointlist[0])) PetscCall(PetscFree(meshCoords)); 36348a46eb9SPierre Jolivet if (sizeof(PetscInt) != sizeof(out.trianglelist[0])) PetscCall(PetscFree(cells)); 3643a074057SBarry Smith /* Set labels */ 3653a074057SBarry Smith for (v = 0; v < numVertices; ++v) { 3663a074057SBarry Smith if (out.pointmarkerlist[v]) { 3679566063dSJacob Faibussowitsch if (rlabel) PetscCall(DMLabelSetValue(rlabel, v + numCells, out.pointmarkerlist[v])); 3683a074057SBarry Smith } 3693a074057SBarry Smith } 3703a074057SBarry Smith if (interpolate) { 3713a074057SBarry Smith PetscInt e; 3723a074057SBarry Smith 3733a074057SBarry Smith for (e = 0; e < out.numberofedges; e++) { 3743a074057SBarry Smith if (out.edgemarkerlist[e]) { 3753a074057SBarry Smith const PetscInt vertices[2] = {out.edgelist[e * 2 + 0] + numCells, out.edgelist[e * 2 + 1] + numCells}; 3763a074057SBarry Smith const PetscInt *edges; 3773a074057SBarry Smith PetscInt numEdges; 3783a074057SBarry Smith 3799566063dSJacob Faibussowitsch PetscCall(DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges)); 38063a3b9bcSJacob Faibussowitsch PetscCheck(numEdges == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %" PetscInt_FMT, numEdges); 3819566063dSJacob Faibussowitsch if (rlabel) PetscCall(DMLabelSetValue(rlabel, edges[0], out.edgemarkerlist[e])); 3829566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges)); 3833a074057SBarry Smith } 3843a074057SBarry Smith } 3853a074057SBarry Smith } 3869566063dSJacob Faibussowitsch PetscCall(DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE)); 3873a074057SBarry Smith } 3883a074057SBarry Smith #if 0 /* Do not currently support holes */ 3899566063dSJacob Faibussowitsch PetscCall(DMPlexCopyHoles(*dm, boundary)); 3903a074057SBarry Smith #endif 3919566063dSJacob Faibussowitsch PetscCall(FiniOutput_Triangle(&out)); 392800b850cSBarry Smith #if !defined(PETSC_USE_REAL_DOUBLE) 3939566063dSJacob Faibussowitsch PetscCall(PetscFree(maxVolumes)); 394800b850cSBarry Smith #endif 3953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3963a074057SBarry Smith } 397