1*3a074057SBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2*3a074057SBarry Smith 3*3a074057SBarry Smith #include <tetgen.h> 4*3a074057SBarry Smith 5*3a074057SBarry Smith /* This is to fix the tetrahedron orientation from TetGen */ 6*3a074057SBarry Smith static PetscErrorCode DMPlexInvertCells_Internal(PetscInt dim, PetscInt numCells, PetscInt numCorners, int cells[]) 7*3a074057SBarry Smith { 8*3a074057SBarry Smith PetscInt bound = numCells*numCorners, coff; 9*3a074057SBarry Smith PetscErrorCode ierr; 10*3a074057SBarry Smith 11*3a074057SBarry Smith PetscFunctionBegin; 12*3a074057SBarry Smith for (coff = 0; coff < bound; coff += numCorners) { 13*3a074057SBarry Smith ierr = DMPlexInvertCell(dim, numCorners, &cells[coff]);CHKERRQ(ierr); 14*3a074057SBarry Smith } 15*3a074057SBarry Smith PetscFunctionReturn(0); 16*3a074057SBarry Smith } 17*3a074057SBarry Smith 18*3a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexGenerate_Tetgen(DM boundary, PetscBool interpolate, DM *dm) 19*3a074057SBarry Smith { 20*3a074057SBarry Smith MPI_Comm comm; 21*3a074057SBarry Smith DM_Plex *mesh = (DM_Plex *) boundary->data; 22*3a074057SBarry Smith const PetscInt dim = 3; 23*3a074057SBarry Smith const char *labelName = "marker"; 24*3a074057SBarry Smith ::tetgenio in; 25*3a074057SBarry Smith ::tetgenio out; 26*3a074057SBarry Smith DMLabel label; 27*3a074057SBarry Smith PetscInt vStart, vEnd, v, fStart, fEnd, f; 28*3a074057SBarry Smith PetscMPIInt rank; 29*3a074057SBarry Smith PetscErrorCode ierr; 30*3a074057SBarry Smith 31*3a074057SBarry Smith PetscFunctionBegin; 32*3a074057SBarry Smith ierr = PetscObjectGetComm((PetscObject)boundary,&comm);CHKERRQ(ierr); 33*3a074057SBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 34*3a074057SBarry Smith ierr = DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd);CHKERRQ(ierr); 35*3a074057SBarry Smith ierr = DMGetLabel(boundary, labelName, &label);CHKERRQ(ierr); 36*3a074057SBarry Smith 37*3a074057SBarry Smith in.numberofpoints = vEnd - vStart; 38*3a074057SBarry Smith if (in.numberofpoints > 0) { 39*3a074057SBarry Smith PetscSection coordSection; 40*3a074057SBarry Smith Vec coordinates; 41*3a074057SBarry Smith PetscScalar *array; 42*3a074057SBarry Smith 43*3a074057SBarry Smith in.pointlist = new double[in.numberofpoints*dim]; 44*3a074057SBarry Smith in.pointmarkerlist = new int[in.numberofpoints]; 45*3a074057SBarry Smith 46*3a074057SBarry Smith ierr = DMGetCoordinatesLocal(boundary, &coordinates);CHKERRQ(ierr); 47*3a074057SBarry Smith ierr = DMGetCoordinateSection(boundary, &coordSection);CHKERRQ(ierr); 48*3a074057SBarry Smith ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr); 49*3a074057SBarry Smith for (v = vStart; v < vEnd; ++v) { 50*3a074057SBarry Smith const PetscInt idx = v - vStart; 51*3a074057SBarry Smith PetscInt off, d; 52*3a074057SBarry Smith 53*3a074057SBarry Smith ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 54*3a074057SBarry Smith for (d = 0; d < dim; ++d) in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]); 55*3a074057SBarry Smith if (label) { 56*3a074057SBarry Smith PetscInt val; 57*3a074057SBarry Smith 58*3a074057SBarry Smith ierr = DMLabelGetValue(label, v, &val);CHKERRQ(ierr); 59*3a074057SBarry Smith in.pointmarkerlist[idx] = (int) val; 60*3a074057SBarry Smith } 61*3a074057SBarry Smith } 62*3a074057SBarry Smith ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr); 63*3a074057SBarry Smith } 64*3a074057SBarry Smith ierr = DMPlexGetHeightStratum(boundary, 0, &fStart, &fEnd);CHKERRQ(ierr); 65*3a074057SBarry Smith 66*3a074057SBarry Smith in.numberoffacets = fEnd - fStart; 67*3a074057SBarry Smith if (in.numberoffacets > 0) { 68*3a074057SBarry Smith in.facetlist = new tetgenio::facet[in.numberoffacets]; 69*3a074057SBarry Smith in.facetmarkerlist = new int[in.numberoffacets]; 70*3a074057SBarry Smith for (f = fStart; f < fEnd; ++f) { 71*3a074057SBarry Smith const PetscInt idx = f - fStart; 72*3a074057SBarry Smith PetscInt *points = NULL, numPoints, p, numVertices = 0, v; 73*3a074057SBarry Smith 74*3a074057SBarry Smith in.facetlist[idx].numberofpolygons = 1; 75*3a074057SBarry Smith in.facetlist[idx].polygonlist = new tetgenio::polygon[in.facetlist[idx].numberofpolygons]; 76*3a074057SBarry Smith in.facetlist[idx].numberofholes = 0; 77*3a074057SBarry Smith in.facetlist[idx].holelist = NULL; 78*3a074057SBarry Smith 79*3a074057SBarry Smith ierr = DMPlexGetTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr); 80*3a074057SBarry Smith for (p = 0; p < numPoints*2; p += 2) { 81*3a074057SBarry Smith const PetscInt point = points[p]; 82*3a074057SBarry Smith if ((point >= vStart) && (point < vEnd)) points[numVertices++] = point; 83*3a074057SBarry Smith } 84*3a074057SBarry Smith 85*3a074057SBarry Smith tetgenio::polygon *poly = in.facetlist[idx].polygonlist; 86*3a074057SBarry Smith poly->numberofvertices = numVertices; 87*3a074057SBarry Smith poly->vertexlist = new int[poly->numberofvertices]; 88*3a074057SBarry Smith for (v = 0; v < numVertices; ++v) { 89*3a074057SBarry Smith const PetscInt vIdx = points[v] - vStart; 90*3a074057SBarry Smith poly->vertexlist[v] = vIdx; 91*3a074057SBarry Smith } 92*3a074057SBarry Smith if (label) { 93*3a074057SBarry Smith PetscInt val; 94*3a074057SBarry Smith 95*3a074057SBarry Smith ierr = DMLabelGetValue(label, f, &val);CHKERRQ(ierr); 96*3a074057SBarry Smith in.facetmarkerlist[idx] = (int) val; 97*3a074057SBarry Smith } 98*3a074057SBarry Smith ierr = DMPlexRestoreTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr); 99*3a074057SBarry Smith } 100*3a074057SBarry Smith } 101*3a074057SBarry Smith if (!rank) { 102*3a074057SBarry Smith char args[32]; 103*3a074057SBarry Smith 104*3a074057SBarry Smith /* Take away 'Q' for verbose output */ 105*3a074057SBarry Smith ierr = PetscStrcpy(args, "pqezQ");CHKERRQ(ierr); 106*3a074057SBarry Smith if (mesh->tetgenOpts) {::tetrahedralize(mesh->tetgenOpts, &in, &out);} 107*3a074057SBarry Smith else {::tetrahedralize(args, &in, &out);} 108*3a074057SBarry Smith } 109*3a074057SBarry Smith { 110*3a074057SBarry Smith DMLabel glabel = NULL; 111*3a074057SBarry Smith const PetscInt numCorners = 4; 112*3a074057SBarry Smith const PetscInt numCells = out.numberoftetrahedra; 113*3a074057SBarry Smith const PetscInt numVertices = out.numberofpoints; 114*3a074057SBarry Smith const double *meshCoords = out.pointlist; 115*3a074057SBarry Smith int *cells = out.tetrahedronlist; 116*3a074057SBarry Smith 117*3a074057SBarry Smith ierr = DMPlexInvertCells_Internal(dim, numCells, numCorners, cells);CHKERRQ(ierr); 118*3a074057SBarry Smith ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm);CHKERRQ(ierr); 119*3a074057SBarry Smith if (label) {ierr = DMCreateLabel(*dm, labelName); ierr = DMGetLabel(*dm, labelName, &glabel);} 120*3a074057SBarry Smith /* Set labels */ 121*3a074057SBarry Smith for (v = 0; v < numVertices; ++v) { 122*3a074057SBarry Smith if (out.pointmarkerlist[v]) { 123*3a074057SBarry Smith if (glabel) {ierr = DMLabelSetValue(glabel, v+numCells, out.pointmarkerlist[v]);CHKERRQ(ierr);} 124*3a074057SBarry Smith } 125*3a074057SBarry Smith } 126*3a074057SBarry Smith if (interpolate) { 127*3a074057SBarry Smith #if 0 128*3a074057SBarry Smith PetscInt e; 129*3a074057SBarry Smith 130*3a074057SBarry Smith /* This check is never actually executed for ctetgen (which never returns edgemarkers) and seems to be broken for 131*3a074057SBarry Smith * tetgen */ 132*3a074057SBarry Smith for (e = 0; e < out.numberofedges; e++) { 133*3a074057SBarry Smith if (out.edgemarkerlist[e]) { 134*3a074057SBarry Smith const PetscInt vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells}; 135*3a074057SBarry Smith const PetscInt *edges; 136*3a074057SBarry Smith PetscInt numEdges; 137*3a074057SBarry Smith 138*3a074057SBarry Smith ierr = DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 139*3a074057SBarry Smith if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges); 140*3a074057SBarry Smith if (glabel) {ierr = DMLabelSetValue(glabel, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);} 141*3a074057SBarry Smith ierr = DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 142*3a074057SBarry Smith } 143*3a074057SBarry Smith } 144*3a074057SBarry Smith #endif 145*3a074057SBarry Smith for (f = 0; f < out.numberoftrifaces; f++) { 146*3a074057SBarry Smith if (out.trifacemarkerlist[f]) { 147*3a074057SBarry Smith const PetscInt vertices[3] = {out.trifacelist[f*3+0]+numCells, out.trifacelist[f*3+1]+numCells, out.trifacelist[f*3+2]+numCells}; 148*3a074057SBarry Smith const PetscInt *faces; 149*3a074057SBarry Smith PetscInt numFaces; 150*3a074057SBarry Smith 151*3a074057SBarry Smith ierr = DMPlexGetFullJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr); 152*3a074057SBarry Smith if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces); 153*3a074057SBarry Smith if (glabel) {ierr = DMLabelSetValue(glabel, faces[0], out.trifacemarkerlist[f]);CHKERRQ(ierr);} 154*3a074057SBarry Smith ierr = DMPlexRestoreJoin(*dm, 3, vertices, &numFaces, &faces);CHKERRQ(ierr); 155*3a074057SBarry Smith } 156*3a074057SBarry Smith } 157*3a074057SBarry Smith } 158*3a074057SBarry Smith ierr = DMPlexSetRefinementUniform(*dm, PETSC_FALSE);CHKERRQ(ierr); 159*3a074057SBarry Smith } 160*3a074057SBarry Smith PetscFunctionReturn(0); 161*3a074057SBarry Smith } 162*3a074057SBarry Smith 163*3a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexRefine_Tetgen(DM dm, double *maxVolumes, DM *dmRefined) 164*3a074057SBarry Smith { 165*3a074057SBarry Smith MPI_Comm comm; 166*3a074057SBarry Smith const PetscInt dim = 3; 167*3a074057SBarry Smith const char *labelName = "marker"; 168*3a074057SBarry Smith ::tetgenio in; 169*3a074057SBarry Smith ::tetgenio out; 170*3a074057SBarry Smith DMLabel label; 171*3a074057SBarry Smith PetscInt vStart, vEnd, v, cStart, cEnd, c, depth, depthGlobal; 172*3a074057SBarry Smith PetscMPIInt rank; 173*3a074057SBarry Smith PetscErrorCode ierr; 174*3a074057SBarry Smith 175*3a074057SBarry Smith PetscFunctionBegin; 176*3a074057SBarry Smith ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 177*3a074057SBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 178*3a074057SBarry Smith ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 179*3a074057SBarry Smith ierr = MPIU_Allreduce(&depth, &depthGlobal, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr); 180*3a074057SBarry Smith ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 181*3a074057SBarry Smith ierr = DMGetLabel(dm, labelName, &label);CHKERRQ(ierr); 182*3a074057SBarry Smith 183*3a074057SBarry Smith in.numberofpoints = vEnd - vStart; 184*3a074057SBarry Smith if (in.numberofpoints > 0) { 185*3a074057SBarry Smith PetscSection coordSection; 186*3a074057SBarry Smith Vec coordinates; 187*3a074057SBarry Smith PetscScalar *array; 188*3a074057SBarry Smith 189*3a074057SBarry Smith in.pointlist = new double[in.numberofpoints*dim]; 190*3a074057SBarry Smith in.pointmarkerlist = new int[in.numberofpoints]; 191*3a074057SBarry Smith 192*3a074057SBarry Smith ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 193*3a074057SBarry Smith ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 194*3a074057SBarry Smith ierr = VecGetArray(coordinates, &array);CHKERRQ(ierr); 195*3a074057SBarry Smith for (v = vStart; v < vEnd; ++v) { 196*3a074057SBarry Smith const PetscInt idx = v - vStart; 197*3a074057SBarry Smith PetscInt off, d; 198*3a074057SBarry Smith 199*3a074057SBarry Smith ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 200*3a074057SBarry Smith for (d = 0; d < dim; ++d) in.pointlist[idx*dim + d] = PetscRealPart(array[off+d]); 201*3a074057SBarry Smith if (label) { 202*3a074057SBarry Smith PetscInt val; 203*3a074057SBarry Smith 204*3a074057SBarry Smith ierr = DMLabelGetValue(label, v, &val);CHKERRQ(ierr); 205*3a074057SBarry Smith in.pointmarkerlist[idx] = (int) val; 206*3a074057SBarry Smith } 207*3a074057SBarry Smith } 208*3a074057SBarry Smith ierr = VecRestoreArray(coordinates, &array);CHKERRQ(ierr); 209*3a074057SBarry Smith } 210*3a074057SBarry Smith ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 211*3a074057SBarry Smith 212*3a074057SBarry Smith in.numberofcorners = 4; 213*3a074057SBarry Smith in.numberoftetrahedra = cEnd - cStart; 214*3a074057SBarry Smith in.tetrahedronvolumelist = (double*) maxVolumes; 215*3a074057SBarry Smith if (in.numberoftetrahedra > 0) { 216*3a074057SBarry Smith in.tetrahedronlist = new int[in.numberoftetrahedra*in.numberofcorners]; 217*3a074057SBarry Smith for (c = cStart; c < cEnd; ++c) { 218*3a074057SBarry Smith const PetscInt idx = c - cStart; 219*3a074057SBarry Smith PetscInt *closure = NULL; 220*3a074057SBarry Smith PetscInt closureSize; 221*3a074057SBarry Smith 222*3a074057SBarry Smith ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 223*3a074057SBarry Smith if ((closureSize != 5) && (closureSize != 15)) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Mesh has cell which is not a tetrahedron, %D vertices in closure", closureSize); 224*3a074057SBarry Smith for (v = 0; v < 4; ++v) { 225*3a074057SBarry Smith in.tetrahedronlist[idx*in.numberofcorners + v] = closure[(v+closureSize-4)*2] - vStart; 226*3a074057SBarry Smith } 227*3a074057SBarry Smith ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 228*3a074057SBarry Smith } 229*3a074057SBarry Smith } 230*3a074057SBarry Smith /* TODO: Put in boundary faces with markers */ 231*3a074057SBarry Smith if (!rank) { 232*3a074057SBarry Smith char args[32]; 233*3a074057SBarry Smith 234*3a074057SBarry Smith #if 1 235*3a074057SBarry Smith /* Take away 'Q' for verbose output */ 236*3a074057SBarry Smith ierr = PetscStrcpy(args, "qezQra");CHKERRQ(ierr); 237*3a074057SBarry Smith #else 238*3a074057SBarry Smith ierr = PetscStrcpy(args, "qezraVVVV");CHKERRQ(ierr); 239*3a074057SBarry Smith #endif 240*3a074057SBarry Smith ::tetrahedralize(args, &in, &out); 241*3a074057SBarry Smith } 242*3a074057SBarry Smith in.tetrahedronvolumelist = NULL; 243*3a074057SBarry Smith 244*3a074057SBarry Smith { 245*3a074057SBarry Smith DMLabel rlabel = NULL; 246*3a074057SBarry Smith const PetscInt numCorners = 4; 247*3a074057SBarry Smith const PetscInt numCells = out.numberoftetrahedra; 248*3a074057SBarry Smith const PetscInt numVertices = out.numberofpoints; 249*3a074057SBarry Smith const double *meshCoords = out.pointlist; 250*3a074057SBarry Smith int *cells = out.tetrahedronlist; 251*3a074057SBarry Smith 252*3a074057SBarry Smith PetscBool interpolate = depthGlobal > 1 ? PETSC_TRUE : PETSC_FALSE; 253*3a074057SBarry Smith 254*3a074057SBarry Smith ierr = DMPlexInvertCells_Internal(dim, numCells, numCorners, cells);CHKERRQ(ierr); 255*3a074057SBarry Smith ierr = DMPlexCreateFromCellList(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined);CHKERRQ(ierr); 256*3a074057SBarry Smith if (label) {ierr = DMCreateLabel(*dmRefined, labelName); ierr = DMGetLabel(*dmRefined, labelName, &rlabel);} 257*3a074057SBarry Smith /* Set labels */ 258*3a074057SBarry Smith for (v = 0; v < numVertices; ++v) { 259*3a074057SBarry Smith if (out.pointmarkerlist[v]) { 260*3a074057SBarry Smith if (rlabel) {ierr = DMLabelSetValue(rlabel, v+numCells, out.pointmarkerlist[v]);CHKERRQ(ierr);} 261*3a074057SBarry Smith } 262*3a074057SBarry Smith } 263*3a074057SBarry Smith if (interpolate) { 264*3a074057SBarry Smith PetscInt f; 265*3a074057SBarry Smith #if 0 266*3a074057SBarry Smith PetscInt e; 267*3a074057SBarry Smith 268*3a074057SBarry Smith for (e = 0; e < out.numberofedges; e++) { 269*3a074057SBarry Smith if (out.edgemarkerlist[e]) { 270*3a074057SBarry Smith const PetscInt vertices[2] = {out.edgelist[e*2+0]+numCells, out.edgelist[e*2+1]+numCells}; 271*3a074057SBarry Smith const PetscInt *edges; 272*3a074057SBarry Smith PetscInt numEdges; 273*3a074057SBarry Smith 274*3a074057SBarry Smith ierr = DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 275*3a074057SBarry Smith if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges); 276*3a074057SBarry Smith if (rlabel) {ierr = DMLabelSetValue(rlabel, edges[0], out.edgemarkerlist[e]);CHKERRQ(ierr);} 277*3a074057SBarry Smith ierr = DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges);CHKERRQ(ierr); 278*3a074057SBarry Smith } 279*3a074057SBarry Smith } 280*3a074057SBarry Smith #endif 281*3a074057SBarry Smith for (f = 0; f < out.numberoftrifaces; f++) { 282*3a074057SBarry Smith if (out.trifacemarkerlist[f]) { 283*3a074057SBarry Smith const PetscInt vertices[3] = {out.trifacelist[f*3+0]+numCells, out.trifacelist[f*3+1]+numCells, out.trifacelist[f*3+2]+numCells}; 284*3a074057SBarry Smith const PetscInt *faces; 285*3a074057SBarry Smith PetscInt numFaces; 286*3a074057SBarry Smith 287*3a074057SBarry Smith ierr = DMPlexGetFullJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr); 288*3a074057SBarry Smith if (numFaces != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces); 289*3a074057SBarry Smith if (rlabel) {ierr = DMLabelSetValue(rlabel, faces[0], out.trifacemarkerlist[f]);CHKERRQ(ierr);} 290*3a074057SBarry Smith ierr = DMPlexRestoreJoin(*dmRefined, 3, vertices, &numFaces, &faces);CHKERRQ(ierr); 291*3a074057SBarry Smith } 292*3a074057SBarry Smith } 293*3a074057SBarry Smith } 294*3a074057SBarry Smith ierr = DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE);CHKERRQ(ierr); 295*3a074057SBarry Smith } 296*3a074057SBarry Smith PetscFunctionReturn(0); 297*3a074057SBarry Smith } 298