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