10bbe5d31Sbarral #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 20bbe5d31Sbarral #ifdef PETSC_HAVE_PRAGMATIC 30bbe5d31Sbarral #include <pragmatic/cpragmatic.h> 40bbe5d31Sbarral #endif 50bbe5d31Sbarral 6f7bcbcbbSMatthew G. Knepley /* 7f7bcbcbbSMatthew G. Knepley DMPlexRemesh_Internal - Generates a new mesh conforming to a metric field. 80bbe5d31Sbarral 9f7bcbcbbSMatthew G. Knepley Input Parameters: 10f7bcbcbbSMatthew G. Knepley + dm - The DM object 1151ff0a1cSMatthew G. Knepley . vertexMetric - The metric to which the mesh is adapted, defined vertex-wise in a LOCAL vector 122e62ab5aSMatthew G. Knepley . remeshBd - Flag to allow boundary changes 13f7bcbcbbSMatthew G. Knepley - bdLabelName - Label name for boundary tags which are preserved in dmNew, or NULL. Should not be "_boundary_". 14f7bcbcbbSMatthew G. Knepley 15f7bcbcbbSMatthew G. Knepley Output Parameter: 16f7bcbcbbSMatthew G. Knepley . dmNew - the new DM 17f7bcbcbbSMatthew G. Knepley 18f7bcbcbbSMatthew G. Knepley Level: advanced 19f7bcbcbbSMatthew G. Knepley 20f7bcbcbbSMatthew G. Knepley .seealso: DMCoarsen(), DMRefine() 21f7bcbcbbSMatthew G. Knepley */ 222e62ab5aSMatthew G. Knepley PetscErrorCode DMPlexRemesh_Internal(DM dm, Vec vertexMetric, const char bdLabelName[], PetscBool remeshBd, DM *dmNew) 23f7bcbcbbSMatthew G. Knepley { 24bc1d7e78SMatthew G. Knepley MPI_Comm comm; 25f7bcbcbbSMatthew G. Knepley const char *bdName = "_boundary_"; 26*2b5f8d01SMatthew G. Knepley #if 0 27*2b5f8d01SMatthew G. Knepley DM odm = dm; 28*2b5f8d01SMatthew G. Knepley #endif 29*2b5f8d01SMatthew G. Knepley DM udm, cdm; 302ee120b0SMatthew G. Knepley DMLabel bdLabel = NULL, bdLabelFull; 311838b4a6SMatthew G. Knepley IS bdIS, globalVertexNum; 32f7bcbcbbSMatthew G. Knepley PetscSection coordSection; 33f7bcbcbbSMatthew G. Knepley Vec coordinates; 34f7bcbcbbSMatthew G. Knepley const PetscScalar *coords, *met; 351838b4a6SMatthew G. Knepley const PetscInt *bdFacesFull, *gV; 361838b4a6SMatthew G. Knepley PetscInt *bdFaces, *bdFaceIds, *l2gv; 37f7bcbcbbSMatthew G. Knepley PetscReal *x, *y, *z, *metric; 382ee120b0SMatthew G. Knepley PetscInt *cells; 391838b4a6SMatthew G. Knepley PetscInt dim, cStart, cEnd, numCells, c, coff, vStart, vEnd, numVertices, numLocVertices, v; 402ee120b0SMatthew G. Knepley PetscInt off, maxConeSize, numBdFaces, f, bdSize; 41f7bcbcbbSMatthew G. Knepley PetscBool flg; 422ee120b0SMatthew G. Knepley #ifdef PETSC_HAVE_PRAGMATIC 432ee120b0SMatthew G. Knepley DMLabel bdLabelNew; 442ee120b0SMatthew G. Knepley PetscScalar *coordsNew; 452ee120b0SMatthew G. Knepley PetscInt *bdTags; 462ee120b0SMatthew G. Knepley PetscReal *xNew[3] = {NULL, NULL, NULL}; 472ee120b0SMatthew G. Knepley PetscInt *cellsNew; 482ee120b0SMatthew G. Knepley PetscInt d, numCellsNew, numVerticesNew; 492ee120b0SMatthew G. Knepley PetscInt numCornersNew, fStart, fEnd; 502ee120b0SMatthew G. Knepley #endif 51bc1d7e78SMatthew G. Knepley PetscMPIInt numProcs; 52f7bcbcbbSMatthew G. Knepley PetscErrorCode ierr; 53f7bcbcbbSMatthew G. Knepley 54f7bcbcbbSMatthew G. Knepley PetscFunctionBegin; 55bbc23918SMatthew G. Knepley /* Check for FEM adjacency flags */ 56f7bcbcbbSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 57f7bcbcbbSMatthew G. Knepley PetscValidHeaderSpecific(vertexMetric, VEC_CLASSID, 2); 58f7bcbcbbSMatthew G. Knepley PetscValidCharPointer(bdLabelName, 3); 592e62ab5aSMatthew G. Knepley PetscValidPointer(dmNew, 5); 60bc1d7e78SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 61bc1d7e78SMatthew G. Knepley ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 62f7bcbcbbSMatthew G. Knepley if (bdLabelName) { 63f7bcbcbbSMatthew G. Knepley size_t len; 64f7bcbcbbSMatthew G. Knepley 65f7bcbcbbSMatthew G. Knepley ierr = PetscStrcmp(bdLabelName, bdName, &flg);CHKERRQ(ierr); 66bc1d7e78SMatthew G. Knepley if (flg) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName); 67f7bcbcbbSMatthew G. Knepley ierr = PetscStrlen(bdLabelName, &len);CHKERRQ(ierr); 68f7bcbcbbSMatthew G. Knepley if (len) { 69f7bcbcbbSMatthew G. Knepley ierr = DMGetLabel(dm, bdLabelName, &bdLabel);CHKERRQ(ierr); 70bc1d7e78SMatthew G. Knepley if (!bdLabel) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Label \"%s\" does not exist in DM", bdLabelName); 71f7bcbcbbSMatthew G. Knepley } 72f7bcbcbbSMatthew G. Knepley } 73a2e1bd45SMatthew G. Knepley /* Add overlap for Pragmatic */ 74bbc23918SMatthew G. Knepley #if 0 75bbc23918SMatthew G. Knepley /* Check for overlap by looking for cell in the SF */ 76bbc23918SMatthew G. Knepley if (!overlapped) { 77a2e1bd45SMatthew G. Knepley ierr = DMPlexDistributeOverlap(odm, 1, NULL, &dm);CHKERRQ(ierr); 78a2e1bd45SMatthew G. Knepley if (!dm) {dm = odm; ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);} 79bbc23918SMatthew G. Knepley } 80bbc23918SMatthew G. Knepley #endif 81f7bcbcbbSMatthew G. Knepley /* Get mesh information */ 82f7bcbcbbSMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 83f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 84f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 85f7bcbcbbSMatthew G. Knepley ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr); 86f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr); 87f7bcbcbbSMatthew G. Knepley numCells = cEnd - cStart; 88f7bcbcbbSMatthew G. Knepley numVertices = vEnd - vStart; 89f7bcbcbbSMatthew G. Knepley ierr = PetscCalloc5(numVertices, &x, numVertices, &y, numVertices, &z, numVertices*PetscSqr(dim), &metric, numCells*maxConeSize, &cells);CHKERRQ(ierr); 90f7bcbcbbSMatthew G. Knepley for (c = 0, coff = 0; c < numCells; ++c) { 91f7bcbcbbSMatthew G. Knepley const PetscInt *cone; 92f7bcbcbbSMatthew G. Knepley PetscInt coneSize, cl; 93f7bcbcbbSMatthew G. Knepley 94f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr); 95f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr); 96f7bcbcbbSMatthew G. Knepley for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart; 97f7bcbcbbSMatthew G. Knepley } 981838b4a6SMatthew G. Knepley ierr = PetscCalloc1(numVertices, &l2gv);CHKERRQ(ierr); 991838b4a6SMatthew G. Knepley ierr = DMPlexGetVertexNumbering(udm, &globalVertexNum);CHKERRQ(ierr); 1001838b4a6SMatthew G. Knepley ierr = ISGetIndices(globalVertexNum, &gV);CHKERRQ(ierr); 1011838b4a6SMatthew G. Knepley for (v = 0, numLocVertices = 0; v < numVertices; ++v) { 1021838b4a6SMatthew G. Knepley if (gV[v] >= 0) ++numLocVertices; 1031838b4a6SMatthew G. Knepley l2gv[v] = gV[v] < 0 ? -(gV[v]+1) : gV[v]; 1041838b4a6SMatthew G. Knepley } 1051838b4a6SMatthew G. Knepley ierr = ISRestoreIndices(globalVertexNum, &gV);CHKERRQ(ierr); 106f7bcbcbbSMatthew G. Knepley ierr = DMDestroy(&udm);CHKERRQ(ierr); 107f7bcbcbbSMatthew G. Knepley ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 108f7bcbcbbSMatthew G. Knepley ierr = DMGetDefaultSection(cdm, &coordSection);CHKERRQ(ierr); 109f7bcbcbbSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 110f7bcbcbbSMatthew G. Knepley ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 111f7bcbcbbSMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 112f7bcbcbbSMatthew G. Knepley ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 1135f2e139eSMatthew G. Knepley x[v-vStart] = PetscRealPart(coords[off+0]); 1145f2e139eSMatthew G. Knepley if (dim > 1) y[v-vStart] = PetscRealPart(coords[off+1]); 1155f2e139eSMatthew G. Knepley if (dim > 2) z[v-vStart] = PetscRealPart(coords[off+2]); 116f7bcbcbbSMatthew G. Knepley } 117f7bcbcbbSMatthew G. Knepley ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 118f7bcbcbbSMatthew G. Knepley /* Get boundary mesh */ 119f7bcbcbbSMatthew G. Knepley ierr = DMLabelCreate(bdName, &bdLabelFull);CHKERRQ(ierr); 120f7bcbcbbSMatthew G. Knepley ierr = DMPlexMarkBoundaryFaces(dm, bdLabelFull);CHKERRQ(ierr); 121f7bcbcbbSMatthew G. Knepley ierr = DMLabelGetStratumIS(bdLabelFull, 1, &bdIS);CHKERRQ(ierr); 122f7bcbcbbSMatthew G. Knepley ierr = DMLabelGetStratumSize(bdLabelFull, 1, &numBdFaces);CHKERRQ(ierr); 123f7bcbcbbSMatthew G. Knepley ierr = ISGetIndices(bdIS, &bdFacesFull);CHKERRQ(ierr); 124f7bcbcbbSMatthew G. Knepley for (f = 0, bdSize = 0; f < numBdFaces; ++f) { 125f7bcbcbbSMatthew G. Knepley PetscInt *closure = NULL; 126f7bcbcbbSMatthew G. Knepley PetscInt closureSize, cl; 127f7bcbcbbSMatthew G. Knepley 128f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 129f7bcbcbbSMatthew G. Knepley for (cl = 0; cl < closureSize*2; cl += 2) { 130f7bcbcbbSMatthew G. Knepley if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize; 131f7bcbcbbSMatthew G. Knepley } 132f7bcbcbbSMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 133f7bcbcbbSMatthew G. Knepley } 134f7bcbcbbSMatthew G. Knepley ierr = PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds);CHKERRQ(ierr); 135f7bcbcbbSMatthew G. Knepley for (f = 0, bdSize = 0; f < numBdFaces; ++f) { 136f7bcbcbbSMatthew G. Knepley PetscInt *closure = NULL; 137f7bcbcbbSMatthew G. Knepley PetscInt closureSize, cl; 138f7bcbcbbSMatthew G. Knepley 139f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 140f7bcbcbbSMatthew G. Knepley for (cl = 0; cl < closureSize*2; cl += 2) { 141f7bcbcbbSMatthew G. Knepley if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart; 142f7bcbcbbSMatthew G. Knepley } 143f7bcbcbbSMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 144f7bcbcbbSMatthew G. Knepley if (bdLabel) {ierr = DMLabelGetValue(bdLabel, bdFacesFull[f], &bdFaceIds[f]);CHKERRQ(ierr);} 145f7bcbcbbSMatthew G. Knepley else {bdFaceIds[f] = 1;} 146f7bcbcbbSMatthew G. Knepley } 147f7bcbcbbSMatthew G. Knepley ierr = ISDestroy(&bdIS);CHKERRQ(ierr); 148f7bcbcbbSMatthew G. Knepley ierr = DMLabelDestroy(&bdLabelFull);CHKERRQ(ierr); 149f7bcbcbbSMatthew G. Knepley /* Get metric */ 150f7bcbcbbSMatthew G. Knepley ierr = VecGetArrayRead(vertexMetric, &met);CHKERRQ(ierr); 1515f2e139eSMatthew G. Knepley for (v = 0; v < (vEnd-vStart)*PetscSqr(dim); ++v) metric[v] = PetscRealPart(met[v]); 152f7bcbcbbSMatthew G. Knepley ierr = VecRestoreArrayRead(vertexMetric, &met);CHKERRQ(ierr); 153bbc23918SMatthew G. Knepley #if 0 154a2e1bd45SMatthew G. Knepley /* Destroy overlap mesh */ 155a2e1bd45SMatthew G. Knepley ierr = DMDestroy(&dm);CHKERRQ(ierr); 156bbc23918SMatthew G. Knepley #endif 157f7bcbcbbSMatthew G. Knepley /* Create new mesh */ 158f7bcbcbbSMatthew G. Knepley #ifdef PETSC_HAVE_PRAGMATIC 159f7bcbcbbSMatthew G. Knepley switch (dim) { 160f7bcbcbbSMatthew G. Knepley case 2: 161cc29c497SMatthew G. Knepley pragmatic_2d_mpi_init(&numVertices, &numCells, cells, x, y, l2gv, numLocVertices, comm);break; 162f7bcbcbbSMatthew G. Knepley case 3: 163cc29c497SMatthew G. Knepley pragmatic_3d_mpi_init(&numVertices, &numCells, cells, x, y, z, l2gv, numLocVertices, comm);break; 1641838b4a6SMatthew G. Knepley default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim); 165f7bcbcbbSMatthew G. Knepley } 166f7bcbcbbSMatthew G. Knepley pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds); 167f7bcbcbbSMatthew G. Knepley pragmatic_set_metric(metric); 1682e62ab5aSMatthew G. Knepley pragmatic_adapt(remeshBd ? 1 : 0); 1691838b4a6SMatthew G. Knepley ierr = PetscFree(l2gv);CHKERRQ(ierr); 170f7bcbcbbSMatthew G. Knepley /* Read out mesh */ 1718713909fSMatthew G. Knepley pragmatic_get_info_mpi(&numVerticesNew, &numCellsNew); 172f7bcbcbbSMatthew G. Knepley ierr = PetscMalloc1(numVerticesNew*dim, &coordsNew);CHKERRQ(ierr); 173f7bcbcbbSMatthew G. Knepley switch (dim) { 174f7bcbcbbSMatthew G. Knepley case 2: 175f7bcbcbbSMatthew G. Knepley numCornersNew = 3; 176f7bcbcbbSMatthew G. Knepley ierr = PetscMalloc2(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1]);CHKERRQ(ierr); 1778713909fSMatthew G. Knepley pragmatic_get_coords_2d_mpi(xNew[0], xNew[1]); 178f7bcbcbbSMatthew G. Knepley break; 179f7bcbcbbSMatthew G. Knepley case 3: 180f7bcbcbbSMatthew G. Knepley numCornersNew = 4; 181f7bcbcbbSMatthew G. Knepley ierr = PetscMalloc3(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1], numVerticesNew, &xNew[2]);CHKERRQ(ierr); 1828713909fSMatthew G. Knepley pragmatic_get_coords_3d_mpi(xNew[0], xNew[1], xNew[2]); 183f7bcbcbbSMatthew G. Knepley break; 184f7bcbcbbSMatthew G. Knepley default: 185bc1d7e78SMatthew G. Knepley SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim); 186f7bcbcbbSMatthew G. Knepley } 187f7bcbcbbSMatthew G. Knepley for (v = 0; v < numVerticesNew; ++v) {for (d = 0; d < dim; ++d) coordsNew[v*dim+d] = xNew[d][v];} 188f7bcbcbbSMatthew G. Knepley ierr = PetscMalloc1(numCellsNew*(dim+1), &cellsNew);CHKERRQ(ierr); 189f7bcbcbbSMatthew G. Knepley pragmatic_get_elements(cellsNew); 19051ff0a1cSMatthew G. Knepley ierr = DMPlexCreateFromCellListParallel(comm, dim, numCellsNew, numVerticesNew, numCornersNew, PETSC_TRUE, cellsNew, dim, coordsNew, NULL, dmNew);CHKERRQ(ierr); 191f7bcbcbbSMatthew G. Knepley /* Read out boundary label */ 192f7bcbcbbSMatthew G. Knepley pragmatic_get_boundaryTags(&bdTags); 193f7bcbcbbSMatthew G. Knepley ierr = DMCreateLabel(*dmNew, bdLabel ? bdLabelName : bdName);CHKERRQ(ierr); 194f7bcbcbbSMatthew G. Knepley ierr = DMGetLabel(*dmNew, bdLabel ? bdLabelName : bdName, &bdLabelNew);CHKERRQ(ierr); 195f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd);CHKERRQ(ierr); 196f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd);CHKERRQ(ierr); 197f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd);CHKERRQ(ierr); 198f7bcbcbbSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 199f7bcbcbbSMatthew G. Knepley /* Only for simplicial meshes */ 200f7bcbcbbSMatthew G. Knepley coff = (c-cStart)*(dim+1); 2018713909fSMatthew G. Knepley /* d is the local cell number of the vertex opposite to the face we are marking */ 202f7bcbcbbSMatthew G. Knepley for (d = 0; d < dim+1; ++d) { 203f7bcbcbbSMatthew G. Knepley if (bdTags[coff+d]) { 2048713909fSMatthew G. Knepley const PetscInt perm[4][4] = {{-1, -1, -1, -1}, {-1, -1, -1, -1}, {1, 2, 0, -1}, {3, 2, 1, 0}}; /* perm[d] = face opposite */ 2058713909fSMatthew G. Knepley const PetscInt *cone; 206f7bcbcbbSMatthew G. Knepley 2078713909fSMatthew G. Knepley /* Mark face opposite to this vertex: This pattern is specified in DMPlexGetRawFaces_Internal() */ 2088713909fSMatthew G. Knepley ierr = DMPlexGetCone(*dmNew, c, &cone);CHKERRQ(ierr); 2098713909fSMatthew G. Knepley ierr = DMLabelSetValue(bdLabelNew, cone[perm[dim][d]], bdTags[coff+d]);CHKERRQ(ierr); 210f7bcbcbbSMatthew G. Knepley } 211f7bcbcbbSMatthew G. Knepley } 212f7bcbcbbSMatthew G. Knepley } 213f7bcbcbbSMatthew G. Knepley /* Cleanup */ 214f7bcbcbbSMatthew G. Knepley switch (dim) { 215f7bcbcbbSMatthew G. Knepley case 2: ierr = PetscFree2(xNew[0], xNew[1]);CHKERRQ(ierr);break; 216f7bcbcbbSMatthew G. Knepley case 3: ierr = PetscFree3(xNew[0], xNew[1], xNew[2]);CHKERRQ(ierr);break; 217f7bcbcbbSMatthew G. Knepley } 218f7bcbcbbSMatthew G. Knepley ierr = PetscFree(cellsNew);CHKERRQ(ierr); 219f7bcbcbbSMatthew G. Knepley ierr = PetscFree5(x, y, z, metric, cells);CHKERRQ(ierr); 220f7bcbcbbSMatthew G. Knepley ierr = PetscFree2(bdFaces, bdFaceIds);CHKERRQ(ierr); 221f7bcbcbbSMatthew G. Knepley ierr = PetscFree(coordsNew);CHKERRQ(ierr); 222f7bcbcbbSMatthew G. Knepley pragmatic_finalize(); 223f7bcbcbbSMatthew G. Knepley #else 224bc1d7e78SMatthew G. Knepley SETERRQ(comm, PETSC_ERR_SUP, "Remeshing needs external package support.\nPlease reconfigure with --download-pragmatic."); 225f7bcbcbbSMatthew G. Knepley #endif 226f7bcbcbbSMatthew G. Knepley PetscFunctionReturn(0); 227f7bcbcbbSMatthew G. Knepley } 2280bbe5d31Sbarral 2290f7e110dSbarral /*@ 2300f7e110dSbarral DMPlexAdapt - Generates a mesh adapted to the specified metric field using the pragmatic library. 2310f7e110dSbarral 2320f7e110dSbarral Input Parameters: 2330f7e110dSbarral + dm - The DM object 2340f7e110dSbarral . metric - The metric to which the mesh is adapted, defined vertex-wise. 235379fadc5Sbarral - bdyLabelName - Label name for boundary tags. These will be preserved in the output mesh. bdyLabelName should be "" (empty string) if there is no such label, and should be different from "boundary". 2360f7e110dSbarral 2370f7e110dSbarral Output Parameter: 238f7bcbcbbSMatthew G. Knepley . dmAdapt - Pointer to the DM object containing the adapted mesh 2390f7e110dSbarral 2400f7e110dSbarral Level: advanced 2410f7e110dSbarral 2420f7e110dSbarral .seealso: DMCoarsen(), DMRefine() 2430f7e110dSbarral @*/ 244f7bcbcbbSMatthew G. Knepley PetscErrorCode DMPlexAdapt(DM dm, Vec metric, const char bdLabelName[], DM *dmAdapt) 2450bbe5d31Sbarral { 2462e62ab5aSMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 2470bbe5d31Sbarral PetscErrorCode ierr; 2480bbe5d31Sbarral 2490bbe5d31Sbarral PetscFunctionBegin; 2502e62ab5aSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2512e62ab5aSMatthew G. Knepley PetscValidHeaderSpecific(metric, VEC_CLASSID, 2); 2522e62ab5aSMatthew G. Knepley PetscValidCharPointer(bdLabelName, 3); 2532e62ab5aSMatthew G. Knepley PetscValidPointer(dmAdapt, 4); 2542e62ab5aSMatthew G. Knepley ierr = DMPlexRemesh_Internal(dm, metric, bdLabelName, mesh->remeshBd, dmAdapt);CHKERRQ(ierr); 2550bbe5d31Sbarral PetscFunctionReturn(0); 2560bbe5d31Sbarral } 257