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