xref: /petsc/src/dm/impls/plex/plexadapt.c (revision acfe6ead8bccb1a1822989ac988f66bbbf297487)
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
11f7bcbcbbSMatthew G. Knepley . vertexMetric - The metric to which the mesh is adapted, defined vertex-wise.
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 {
24f7bcbcbbSMatthew G. Knepley   const char        *bdName = "_boundary_";
25f7bcbcbbSMatthew G. Knepley   DM                 udm, cdm;
262ee120b0SMatthew G. Knepley   DMLabel            bdLabel = NULL, bdLabelFull;
27f7bcbcbbSMatthew G. Knepley   IS                 bdIS;
28f7bcbcbbSMatthew G. Knepley   PetscSection       coordSection;
29f7bcbcbbSMatthew G. Knepley   Vec                coordinates;
30f7bcbcbbSMatthew G. Knepley   const PetscScalar *coords, *met;
31f7bcbcbbSMatthew G. Knepley   const PetscInt    *bdFacesFull;
322ee120b0SMatthew G. Knepley   PetscInt          *bdFaces, *bdFaceIds;
33f7bcbcbbSMatthew G. Knepley   PetscReal         *x, *y, *z, *metric;
342ee120b0SMatthew G. Knepley   PetscInt          *cells;
352ee120b0SMatthew G. Knepley   PetscInt           dim, cStart, cEnd, numCells, c, coff, vStart, vEnd, numVertices, v;
362ee120b0SMatthew G. Knepley   PetscInt           off, maxConeSize, numBdFaces, f, bdSize;
37f7bcbcbbSMatthew G. Knepley   PetscBool          flg;
382ee120b0SMatthew G. Knepley #ifdef PETSC_HAVE_PRAGMATIC
392ee120b0SMatthew G. Knepley   DMLabel            bdLabelNew;
40*acfe6eadSToby Isaac   double            *coordsNew;
412ee120b0SMatthew G. Knepley   PetscInt          *bdTags;
422ee120b0SMatthew G. Knepley   PetscReal         *xNew[3] = {NULL, NULL, NULL};
432ee120b0SMatthew G. Knepley   PetscInt          *cellsNew;
442ee120b0SMatthew G. Knepley   PetscInt           d, numCellsNew, numVerticesNew;
452ee120b0SMatthew G. Knepley   PetscInt           numCornersNew, fStart, fEnd;
462ee120b0SMatthew G. Knepley #endif
47f7bcbcbbSMatthew G. Knepley   PetscErrorCode     ierr;
48f7bcbcbbSMatthew G. Knepley 
49f7bcbcbbSMatthew G. Knepley   PetscFunctionBegin;
50f7bcbcbbSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
51f7bcbcbbSMatthew G. Knepley   PetscValidHeaderSpecific(vertexMetric, VEC_CLASSID, 2);
52f7bcbcbbSMatthew G. Knepley   PetscValidCharPointer(bdLabelName, 3);
532e62ab5aSMatthew G. Knepley   PetscValidPointer(dmNew, 5);
54f7bcbcbbSMatthew G. Knepley   if (bdLabelName) {
55f7bcbcbbSMatthew G. Knepley     size_t len;
56f7bcbcbbSMatthew G. Knepley 
57f7bcbcbbSMatthew G. Knepley     ierr = PetscStrcmp(bdLabelName, bdName, &flg);CHKERRQ(ierr);
58f7bcbcbbSMatthew G. Knepley     if (flg) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName);
59f7bcbcbbSMatthew G. Knepley     ierr = PetscStrlen(bdLabelName, &len);CHKERRQ(ierr);
60f7bcbcbbSMatthew G. Knepley     if (len) {
61f7bcbcbbSMatthew G. Knepley       ierr = DMGetLabel(dm, bdLabelName, &bdLabel);CHKERRQ(ierr);
62f7bcbcbbSMatthew G. Knepley       if (!bdLabel) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Label \"%s\" does not exist in DM", bdLabelName);
63f7bcbcbbSMatthew G. Knepley     }
64f7bcbcbbSMatthew G. Knepley   }
65f7bcbcbbSMatthew G. Knepley   /* Get mesh information */
66f7bcbcbbSMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
67f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
68f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
69f7bcbcbbSMatthew G. Knepley   ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr);
70f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr);
71f7bcbcbbSMatthew G. Knepley   numCells    = cEnd - cStart;
72f7bcbcbbSMatthew G. Knepley   numVertices = vEnd - vStart;
73f7bcbcbbSMatthew G. Knepley   ierr = PetscCalloc5(numVertices, &x, numVertices, &y, numVertices, &z, numVertices*PetscSqr(dim), &metric, numCells*maxConeSize, &cells);CHKERRQ(ierr);
74f7bcbcbbSMatthew G. Knepley   for (c = 0, coff = 0; c < numCells; ++c) {
75f7bcbcbbSMatthew G. Knepley     const PetscInt *cone;
76f7bcbcbbSMatthew G. Knepley     PetscInt        coneSize, cl;
77f7bcbcbbSMatthew G. Knepley 
78f7bcbcbbSMatthew G. Knepley     ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr);
79f7bcbcbbSMatthew G. Knepley     ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr);
80f7bcbcbbSMatthew G. Knepley     for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart;
81f7bcbcbbSMatthew G. Knepley   }
82f7bcbcbbSMatthew G. Knepley   ierr = DMDestroy(&udm);CHKERRQ(ierr);
83f7bcbcbbSMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
84f7bcbcbbSMatthew G. Knepley   ierr = DMGetDefaultSection(cdm, &coordSection);CHKERRQ(ierr);
85f7bcbcbbSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
86f7bcbcbbSMatthew G. Knepley   ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
87f7bcbcbbSMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
88f7bcbcbbSMatthew G. Knepley     ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
895f2e139eSMatthew G. Knepley     x[v-vStart] = PetscRealPart(coords[off+0]);
905f2e139eSMatthew G. Knepley     if (dim > 1) y[v-vStart] = PetscRealPart(coords[off+1]);
915f2e139eSMatthew G. Knepley     if (dim > 2) z[v-vStart] = PetscRealPart(coords[off+2]);
92f7bcbcbbSMatthew G. Knepley   }
93f7bcbcbbSMatthew G. Knepley   ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
94f7bcbcbbSMatthew G. Knepley   /* Get boundary mesh */
95f7bcbcbbSMatthew G. Knepley   ierr = DMLabelCreate(bdName, &bdLabelFull);CHKERRQ(ierr);
96f7bcbcbbSMatthew G. Knepley   ierr = DMPlexMarkBoundaryFaces(dm, bdLabelFull);CHKERRQ(ierr);
97f7bcbcbbSMatthew G. Knepley   ierr = DMLabelGetStratumIS(bdLabelFull, 1, &bdIS);CHKERRQ(ierr);
98f7bcbcbbSMatthew G. Knepley   ierr = DMLabelGetStratumSize(bdLabelFull, 1, &numBdFaces);CHKERRQ(ierr);
99f7bcbcbbSMatthew G. Knepley   ierr = ISGetIndices(bdIS, &bdFacesFull);CHKERRQ(ierr);
100f7bcbcbbSMatthew G. Knepley   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
101f7bcbcbbSMatthew G. Knepley     PetscInt *closure = NULL;
102f7bcbcbbSMatthew G. Knepley     PetscInt  closureSize, cl;
103f7bcbcbbSMatthew G. Knepley 
104f7bcbcbbSMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
105f7bcbcbbSMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
106f7bcbcbbSMatthew G. Knepley       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
107f7bcbcbbSMatthew G. Knepley     }
108f7bcbcbbSMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
109f7bcbcbbSMatthew G. Knepley   }
110f7bcbcbbSMatthew G. Knepley   ierr = PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds);CHKERRQ(ierr);
111f7bcbcbbSMatthew G. Knepley   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
112f7bcbcbbSMatthew G. Knepley     PetscInt *closure = NULL;
113f7bcbcbbSMatthew G. Knepley     PetscInt  closureSize, cl;
114f7bcbcbbSMatthew G. Knepley 
115f7bcbcbbSMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
116f7bcbcbbSMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
117f7bcbcbbSMatthew G. Knepley       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart;
118f7bcbcbbSMatthew G. Knepley     }
119f7bcbcbbSMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
120f7bcbcbbSMatthew G. Knepley     if (bdLabel) {ierr = DMLabelGetValue(bdLabel, bdFacesFull[f], &bdFaceIds[f]);CHKERRQ(ierr);}
121f7bcbcbbSMatthew G. Knepley     else         {bdFaceIds[f] = 1;}
122f7bcbcbbSMatthew G. Knepley   }
123f7bcbcbbSMatthew G. Knepley   ierr = ISDestroy(&bdIS);CHKERRQ(ierr);
124f7bcbcbbSMatthew G. Knepley   ierr = DMLabelDestroy(&bdLabelFull);CHKERRQ(ierr);
125f7bcbcbbSMatthew G. Knepley   /* Get metric */
126f7bcbcbbSMatthew G. Knepley   ierr = VecGetArrayRead(vertexMetric, &met);CHKERRQ(ierr);
1275f2e139eSMatthew G. Knepley   for (v = 0; v < (vEnd-vStart)*PetscSqr(dim); ++v) metric[v] = PetscRealPart(met[v]);
128f7bcbcbbSMatthew G. Knepley   ierr = VecRestoreArrayRead(vertexMetric, &met);CHKERRQ(ierr);
129f7bcbcbbSMatthew G. Knepley   /* Create new mesh */
130f7bcbcbbSMatthew G. Knepley #ifdef PETSC_HAVE_PRAGMATIC
131f7bcbcbbSMatthew G. Knepley   switch (dim) {
132f7bcbcbbSMatthew G. Knepley   case 2:
133f7bcbcbbSMatthew G. Knepley     pragmatic_2d_init(&numVertices, &numCells, cells, x, y);break;
134f7bcbcbbSMatthew G. Knepley   case 3:
135f7bcbcbbSMatthew G. Knepley     pragmatic_3d_init(&numVertices, &numCells, cells, x, y, z);break;
136f7bcbcbbSMatthew G. Knepley   default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim);
137f7bcbcbbSMatthew G. Knepley   }
138f7bcbcbbSMatthew G. Knepley   pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds);
139f7bcbcbbSMatthew G. Knepley   pragmatic_set_metric(metric);
1402e62ab5aSMatthew G. Knepley   pragmatic_adapt(remeshBd ? 1 : 0);
141f7bcbcbbSMatthew G. Knepley   /* Read out mesh */
142f7bcbcbbSMatthew G. Knepley   pragmatic_get_info(&numVerticesNew, &numCellsNew);
143f7bcbcbbSMatthew G. Knepley   ierr = PetscMalloc1(numVerticesNew*dim, &coordsNew);CHKERRQ(ierr);
144f7bcbcbbSMatthew G. Knepley   switch (dim) {
145f7bcbcbbSMatthew G. Knepley   case 2:
146f7bcbcbbSMatthew G. Knepley     numCornersNew = 3;
147f7bcbcbbSMatthew G. Knepley     ierr = PetscMalloc2(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1]);CHKERRQ(ierr);
148f7bcbcbbSMatthew G. Knepley     pragmatic_get_coords_2d(xNew[0], xNew[1]);
149f7bcbcbbSMatthew G. Knepley     break;
150f7bcbcbbSMatthew G. Knepley   case 3:
151f7bcbcbbSMatthew G. Knepley     numCornersNew = 4;
152f7bcbcbbSMatthew G. Knepley     ierr = PetscMalloc3(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1], numVerticesNew, &xNew[2]);CHKERRQ(ierr);
153f7bcbcbbSMatthew G. Knepley     pragmatic_get_coords_3d(xNew[0], xNew[1], xNew[2]);
154f7bcbcbbSMatthew G. Knepley     break;
155f7bcbcbbSMatthew G. Knepley   default:
156f7bcbcbbSMatthew G. Knepley     SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim);
157f7bcbcbbSMatthew G. Knepley   }
158*acfe6eadSToby Isaac   for (v = 0; v < numVerticesNew; ++v) {for (d = 0; d < dim; ++d) coordsNew[v*dim+d] = (double) xNew[d][v];}
159f7bcbcbbSMatthew G. Knepley   ierr = PetscMalloc1(numCellsNew*(dim+1), &cellsNew);CHKERRQ(ierr);
160f7bcbcbbSMatthew G. Knepley   pragmatic_get_elements(cellsNew);
161f7bcbcbbSMatthew G. Knepley   ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject) dm), dim, numCellsNew, numVerticesNew, numCornersNew, PETSC_TRUE, cellsNew, dim, coordsNew, dmNew);CHKERRQ(ierr);
162f7bcbcbbSMatthew G. Knepley   /* Read out boundary label */
163f7bcbcbbSMatthew G. Knepley   pragmatic_get_boundaryTags(&bdTags);
164f7bcbcbbSMatthew G. Knepley   ierr = DMCreateLabel(*dmNew, bdLabel ? bdLabelName : bdName);CHKERRQ(ierr);
165f7bcbcbbSMatthew G. Knepley   ierr = DMGetLabel(*dmNew, bdLabel ? bdLabelName : bdName, &bdLabelNew);CHKERRQ(ierr);
166f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd);CHKERRQ(ierr);
167f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd);CHKERRQ(ierr);
168f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd);CHKERRQ(ierr);
169f7bcbcbbSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
170f7bcbcbbSMatthew G. Knepley     /* Only for simplicial meshes */
171f7bcbcbbSMatthew G. Knepley     coff = (c-cStart)*(dim+1);
172f7bcbcbbSMatthew G. Knepley     for (d = 0; d < dim+1; ++d) {
173f7bcbcbbSMatthew G. Knepley       if (bdTags[coff+d]) {
174f7bcbcbbSMatthew G. Knepley         const PetscInt *faces;
175f7bcbcbbSMatthew G. Knepley         PetscInt        numFaces, vertices[3], p;
176f7bcbcbbSMatthew G. Knepley 
177f7bcbcbbSMatthew G. Knepley         /* Mark face opposite to this vertex */
178f7bcbcbbSMatthew G. Knepley         for (p = 0, v = 0; p < dim+1; ++p) if (p != d) {vertices[v] = cellsNew[coff+p] + vStart; ++v;}
179f7bcbcbbSMatthew G. Knepley         ierr = DMPlexGetFullJoin(*dmNew, dim, vertices, &numFaces, &faces);CHKERRQ(ierr);
180f7bcbcbbSMatthew G. Knepley         if (numFaces != 1) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot find boundary face in new cell %D for vertex %D(%D) with tag %D", c, cellsNew[coff+d], d, bdTags[coff+d]);
181f7bcbcbbSMatthew G. Knepley         if (faces[0] < fStart || faces[0] >= fEnd) SETERRQ5(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Boundary point %D in new cell %D for vertex %D(%D) with tag %D is not a face", faces[0], c, cellsNew[coff+d], d, bdTags[coff+d]);
182f7bcbcbbSMatthew G. Knepley         ierr = DMLabelSetValue(bdLabelNew, faces[0], bdTags[coff+d]);CHKERRQ(ierr);
183f7bcbcbbSMatthew G. Knepley         ierr = DMPlexRestoreJoin(*dmNew, dim, vertices, &numFaces, &faces);CHKERRQ(ierr);
184f7bcbcbbSMatthew G. Knepley       }
185f7bcbcbbSMatthew G. Knepley     }
186f7bcbcbbSMatthew G. Knepley   }
187f7bcbcbbSMatthew G. Knepley   /* Cleanup */
188f7bcbcbbSMatthew G. Knepley   switch (dim) {
189f7bcbcbbSMatthew G. Knepley   case 2: ierr = PetscFree2(xNew[0], xNew[1]);CHKERRQ(ierr);break;
190f7bcbcbbSMatthew G. Knepley   case 3: ierr = PetscFree3(xNew[0], xNew[1], xNew[2]);CHKERRQ(ierr);break;
191f7bcbcbbSMatthew G. Knepley   }
192f7bcbcbbSMatthew G. Knepley   ierr = PetscFree(cellsNew);CHKERRQ(ierr);
193f7bcbcbbSMatthew G. Knepley   ierr = PetscFree5(x, y, z, metric, cells);CHKERRQ(ierr);
194f7bcbcbbSMatthew G. Knepley   ierr = PetscFree2(bdFaces, bdFaceIds);CHKERRQ(ierr);
195f7bcbcbbSMatthew G. Knepley   ierr = PetscFree(coordsNew);CHKERRQ(ierr);
196f7bcbcbbSMatthew G. Knepley   pragmatic_finalize();
197f7bcbcbbSMatthew G. Knepley #else
198f7bcbcbbSMatthew G. Knepley   SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Remeshing needs external package support.\nPlease reconfigure with --download-pragmatic.");
199f7bcbcbbSMatthew G. Knepley #endif
200f7bcbcbbSMatthew G. Knepley   PetscFunctionReturn(0);
201f7bcbcbbSMatthew G. Knepley }
2020bbe5d31Sbarral 
2030f7e110dSbarral /*@
2040f7e110dSbarral   DMPlexAdapt - Generates a mesh adapted to the specified metric field using the pragmatic library.
2050f7e110dSbarral 
2060f7e110dSbarral   Input Parameters:
2070f7e110dSbarral + dm - The DM object
2080f7e110dSbarral . metric - The metric to which the mesh is adapted, defined vertex-wise.
209379fadc5Sbarral - 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".
2100f7e110dSbarral 
2110f7e110dSbarral   Output Parameter:
212f7bcbcbbSMatthew G. Knepley . dmAdapt  - Pointer to the DM object containing the adapted mesh
2130f7e110dSbarral 
2140f7e110dSbarral   Level: advanced
2150f7e110dSbarral 
2160f7e110dSbarral .seealso: DMCoarsen(), DMRefine()
2170f7e110dSbarral @*/
218f7bcbcbbSMatthew G. Knepley PetscErrorCode DMPlexAdapt(DM dm, Vec metric, const char bdLabelName[], DM *dmAdapt)
2190bbe5d31Sbarral {
2202e62ab5aSMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex *) dm->data;
2210bbe5d31Sbarral   PetscErrorCode ierr;
2220bbe5d31Sbarral 
2230bbe5d31Sbarral   PetscFunctionBegin;
2242e62ab5aSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2252e62ab5aSMatthew G. Knepley   PetscValidHeaderSpecific(metric, VEC_CLASSID, 2);
2262e62ab5aSMatthew G. Knepley   PetscValidCharPointer(bdLabelName, 3);
2272e62ab5aSMatthew G. Knepley   PetscValidPointer(dmAdapt, 4);
2282e62ab5aSMatthew G. Knepley   ierr = DMPlexRemesh_Internal(dm, metric, bdLabelName, mesh->remeshBd, dmAdapt);CHKERRQ(ierr);
2290bbe5d31Sbarral   PetscFunctionReturn(0);
2300bbe5d31Sbarral }
231