xref: /petsc/src/dm/impls/plex/plexadapt.c (revision cc29c4979a13630694839f32bcf893e20c172eb2)
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_";
26a2e1bd45SMatthew G. Knepley   DM                 odm = dm, udm, cdm;
272ee120b0SMatthew G. Knepley   DMLabel            bdLabel = NULL, bdLabelFull;
281838b4a6SMatthew G. Knepley   IS                 bdIS, globalVertexNum;
29f7bcbcbbSMatthew G. Knepley   PetscSection       coordSection;
30f7bcbcbbSMatthew G. Knepley   Vec                coordinates;
31f7bcbcbbSMatthew G. Knepley   const PetscScalar *coords, *met;
321838b4a6SMatthew G. Knepley   const PetscInt    *bdFacesFull, *gV;
331838b4a6SMatthew G. Knepley   PetscInt          *bdFaces, *bdFaceIds, *l2gv;
34f7bcbcbbSMatthew G. Knepley   PetscReal         *x, *y, *z, *metric;
352ee120b0SMatthew G. Knepley   PetscInt          *cells;
361838b4a6SMatthew G. Knepley   PetscInt           dim, cStart, cEnd, numCells, c, coff, vStart, vEnd, numVertices, numLocVertices, v;
372ee120b0SMatthew G. Knepley   PetscInt           off, maxConeSize, numBdFaces, f, bdSize;
38f7bcbcbbSMatthew G. Knepley   PetscBool          flg;
392ee120b0SMatthew G. Knepley #ifdef PETSC_HAVE_PRAGMATIC
402ee120b0SMatthew G. Knepley   DMLabel            bdLabelNew;
412ee120b0SMatthew G. Knepley   PetscScalar       *coordsNew;
422ee120b0SMatthew G. Knepley   PetscInt          *bdTags;
432ee120b0SMatthew G. Knepley   PetscReal         *xNew[3] = {NULL, NULL, NULL};
442ee120b0SMatthew G. Knepley   PetscInt          *cellsNew;
452ee120b0SMatthew G. Knepley   PetscInt           d, numCellsNew, numVerticesNew;
462ee120b0SMatthew G. Knepley   PetscInt           numCornersNew, fStart, fEnd;
472ee120b0SMatthew G. Knepley #endif
48bc1d7e78SMatthew G. Knepley   PetscMPIInt        numProcs;
49f7bcbcbbSMatthew G. Knepley   PetscErrorCode     ierr;
50f7bcbcbbSMatthew G. Knepley 
51f7bcbcbbSMatthew G. Knepley   PetscFunctionBegin;
52f7bcbcbbSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
53f7bcbcbbSMatthew G. Knepley   PetscValidHeaderSpecific(vertexMetric, VEC_CLASSID, 2);
54f7bcbcbbSMatthew G. Knepley   PetscValidCharPointer(bdLabelName, 3);
552e62ab5aSMatthew G. Knepley   PetscValidPointer(dmNew, 5);
56bc1d7e78SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
57bc1d7e78SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
58f7bcbcbbSMatthew G. Knepley   if (bdLabelName) {
59f7bcbcbbSMatthew G. Knepley     size_t len;
60f7bcbcbbSMatthew G. Knepley 
61f7bcbcbbSMatthew G. Knepley     ierr = PetscStrcmp(bdLabelName, bdName, &flg);CHKERRQ(ierr);
62bc1d7e78SMatthew G. Knepley     if (flg) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName);
63f7bcbcbbSMatthew G. Knepley     ierr = PetscStrlen(bdLabelName, &len);CHKERRQ(ierr);
64f7bcbcbbSMatthew G. Knepley     if (len) {
65f7bcbcbbSMatthew G. Knepley       ierr = DMGetLabel(dm, bdLabelName, &bdLabel);CHKERRQ(ierr);
66bc1d7e78SMatthew G. Knepley       if (!bdLabel) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Label \"%s\" does not exist in DM", bdLabelName);
67f7bcbcbbSMatthew G. Knepley     }
68f7bcbcbbSMatthew G. Knepley   }
69a2e1bd45SMatthew G. Knepley   /* Add overlap for Pragmatic */
70a2e1bd45SMatthew G. Knepley   ierr = DMPlexDistributeOverlap(odm, 1, NULL, &dm);CHKERRQ(ierr);
71a2e1bd45SMatthew G. Knepley   if (!dm) {dm = odm; ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);}
72f7bcbcbbSMatthew G. Knepley   /* Get mesh information */
73f7bcbcbbSMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
74f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
75f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
76f7bcbcbbSMatthew G. Knepley   ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr);
77f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr);
78f7bcbcbbSMatthew G. Knepley   numCells    = cEnd - cStart;
79f7bcbcbbSMatthew G. Knepley   numVertices = vEnd - vStart;
80f7bcbcbbSMatthew G. Knepley   ierr = PetscCalloc5(numVertices, &x, numVertices, &y, numVertices, &z, numVertices*PetscSqr(dim), &metric, numCells*maxConeSize, &cells);CHKERRQ(ierr);
81f7bcbcbbSMatthew G. Knepley   for (c = 0, coff = 0; c < numCells; ++c) {
82f7bcbcbbSMatthew G. Knepley     const PetscInt *cone;
83f7bcbcbbSMatthew G. Knepley     PetscInt        coneSize, cl;
84f7bcbcbbSMatthew G. Knepley 
85f7bcbcbbSMatthew G. Knepley     ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr);
86f7bcbcbbSMatthew G. Knepley     ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr);
87f7bcbcbbSMatthew G. Knepley     for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart;
88f7bcbcbbSMatthew G. Knepley   }
891838b4a6SMatthew G. Knepley   ierr = PetscCalloc1(numVertices, &l2gv);CHKERRQ(ierr);
901838b4a6SMatthew G. Knepley   ierr = DMPlexGetVertexNumbering(udm, &globalVertexNum);CHKERRQ(ierr);
911838b4a6SMatthew G. Knepley   ierr = ISGetIndices(globalVertexNum, &gV);CHKERRQ(ierr);
921838b4a6SMatthew G. Knepley   for (v = 0, numLocVertices = 0; v < numVertices; ++v) {
931838b4a6SMatthew G. Knepley     if (gV[v] >= 0) ++numLocVertices;
941838b4a6SMatthew G. Knepley     l2gv[v] = gV[v] < 0 ? -(gV[v]+1) : gV[v];
951838b4a6SMatthew G. Knepley   }
961838b4a6SMatthew G. Knepley   ierr = ISRestoreIndices(globalVertexNum, &gV);CHKERRQ(ierr);
97f7bcbcbbSMatthew G. Knepley   ierr = DMDestroy(&udm);CHKERRQ(ierr);
98f7bcbcbbSMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
99f7bcbcbbSMatthew G. Knepley   ierr = DMGetDefaultSection(cdm, &coordSection);CHKERRQ(ierr);
100f7bcbcbbSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
101f7bcbcbbSMatthew G. Knepley   ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
102f7bcbcbbSMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
103f7bcbcbbSMatthew G. Knepley     ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
1045f2e139eSMatthew G. Knepley     x[v-vStart] = PetscRealPart(coords[off+0]);
1055f2e139eSMatthew G. Knepley     if (dim > 1) y[v-vStart] = PetscRealPart(coords[off+1]);
1065f2e139eSMatthew G. Knepley     if (dim > 2) z[v-vStart] = PetscRealPart(coords[off+2]);
107f7bcbcbbSMatthew G. Knepley   }
108f7bcbcbbSMatthew G. Knepley   ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
109f7bcbcbbSMatthew G. Knepley   /* Get boundary mesh */
110f7bcbcbbSMatthew G. Knepley   ierr = DMLabelCreate(bdName, &bdLabelFull);CHKERRQ(ierr);
111f7bcbcbbSMatthew G. Knepley   ierr = DMPlexMarkBoundaryFaces(dm, bdLabelFull);CHKERRQ(ierr);
112f7bcbcbbSMatthew G. Knepley   ierr = DMLabelGetStratumIS(bdLabelFull, 1, &bdIS);CHKERRQ(ierr);
113f7bcbcbbSMatthew G. Knepley   ierr = DMLabelGetStratumSize(bdLabelFull, 1, &numBdFaces);CHKERRQ(ierr);
114f7bcbcbbSMatthew G. Knepley   ierr = ISGetIndices(bdIS, &bdFacesFull);CHKERRQ(ierr);
115f7bcbcbbSMatthew G. Knepley   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
116f7bcbcbbSMatthew G. Knepley     PetscInt *closure = NULL;
117f7bcbcbbSMatthew G. Knepley     PetscInt  closureSize, cl;
118f7bcbcbbSMatthew G. Knepley 
119f7bcbcbbSMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
120f7bcbcbbSMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
121f7bcbcbbSMatthew G. Knepley       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
122f7bcbcbbSMatthew G. Knepley     }
123f7bcbcbbSMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
124f7bcbcbbSMatthew G. Knepley   }
125f7bcbcbbSMatthew G. Knepley   ierr = PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds);CHKERRQ(ierr);
126f7bcbcbbSMatthew G. Knepley   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
127f7bcbcbbSMatthew G. Knepley     PetscInt *closure = NULL;
128f7bcbcbbSMatthew G. Knepley     PetscInt  closureSize, cl;
129f7bcbcbbSMatthew G. Knepley 
130f7bcbcbbSMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
131f7bcbcbbSMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
132f7bcbcbbSMatthew G. Knepley       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart;
133f7bcbcbbSMatthew G. Knepley     }
134f7bcbcbbSMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
135f7bcbcbbSMatthew G. Knepley     if (bdLabel) {ierr = DMLabelGetValue(bdLabel, bdFacesFull[f], &bdFaceIds[f]);CHKERRQ(ierr);}
136f7bcbcbbSMatthew G. Knepley     else         {bdFaceIds[f] = 1;}
137f7bcbcbbSMatthew G. Knepley   }
138f7bcbcbbSMatthew G. Knepley   ierr = ISDestroy(&bdIS);CHKERRQ(ierr);
139f7bcbcbbSMatthew G. Knepley   ierr = DMLabelDestroy(&bdLabelFull);CHKERRQ(ierr);
140f7bcbcbbSMatthew G. Knepley   /* Get metric */
141f7bcbcbbSMatthew G. Knepley   ierr = VecGetArrayRead(vertexMetric, &met);CHKERRQ(ierr);
1425f2e139eSMatthew G. Knepley   for (v = 0; v < (vEnd-vStart)*PetscSqr(dim); ++v) metric[v] = PetscRealPart(met[v]);
143f7bcbcbbSMatthew G. Knepley   ierr = VecRestoreArrayRead(vertexMetric, &met);CHKERRQ(ierr);
144a2e1bd45SMatthew G. Knepley   /* Destroy overlap mesh */
145a2e1bd45SMatthew G. Knepley   ierr = DMDestroy(&dm);CHKERRQ(ierr);
146f7bcbcbbSMatthew G. Knepley   /* Create new mesh */
147f7bcbcbbSMatthew G. Knepley #ifdef PETSC_HAVE_PRAGMATIC
148f7bcbcbbSMatthew G. Knepley   switch (dim) {
149f7bcbcbbSMatthew G. Knepley   case 2:
150*cc29c497SMatthew G. Knepley     pragmatic_2d_mpi_init(&numVertices, &numCells, cells, x, y, l2gv, numLocVertices, comm);break;
151f7bcbcbbSMatthew G. Knepley   case 3:
152*cc29c497SMatthew G. Knepley     pragmatic_3d_mpi_init(&numVertices, &numCells, cells, x, y, z, l2gv, numLocVertices, comm);break;
1531838b4a6SMatthew G. Knepley   default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim);
154f7bcbcbbSMatthew G. Knepley   }
155f7bcbcbbSMatthew G. Knepley   pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds);
156f7bcbcbbSMatthew G. Knepley   pragmatic_set_metric(metric);
1572e62ab5aSMatthew G. Knepley   pragmatic_adapt(remeshBd ? 1 : 0);
1581838b4a6SMatthew G. Knepley   ierr = PetscFree(l2gv);CHKERRQ(ierr);
159f7bcbcbbSMatthew G. Knepley   /* Read out mesh */
160f7bcbcbbSMatthew G. Knepley   pragmatic_get_info(&numVerticesNew, &numCellsNew);
161f7bcbcbbSMatthew G. Knepley   ierr = PetscMalloc1(numVerticesNew*dim, &coordsNew);CHKERRQ(ierr);
162f7bcbcbbSMatthew G. Knepley   switch (dim) {
163f7bcbcbbSMatthew G. Knepley   case 2:
164f7bcbcbbSMatthew G. Knepley     numCornersNew = 3;
165f7bcbcbbSMatthew G. Knepley     ierr = PetscMalloc2(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1]);CHKERRQ(ierr);
166f7bcbcbbSMatthew G. Knepley     pragmatic_get_coords_2d(xNew[0], xNew[1]);
167f7bcbcbbSMatthew G. Knepley     break;
168f7bcbcbbSMatthew G. Knepley   case 3:
169f7bcbcbbSMatthew G. Knepley     numCornersNew = 4;
170f7bcbcbbSMatthew G. Knepley     ierr = PetscMalloc3(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1], numVerticesNew, &xNew[2]);CHKERRQ(ierr);
171f7bcbcbbSMatthew G. Knepley     pragmatic_get_coords_3d(xNew[0], xNew[1], xNew[2]);
172f7bcbcbbSMatthew G. Knepley     break;
173f7bcbcbbSMatthew G. Knepley   default:
174bc1d7e78SMatthew G. Knepley     SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim);
175f7bcbcbbSMatthew G. Knepley   }
176f7bcbcbbSMatthew G. Knepley   for (v = 0; v < numVerticesNew; ++v) {for (d = 0; d < dim; ++d) coordsNew[v*dim+d] = xNew[d][v];}
177f7bcbcbbSMatthew G. Knepley   ierr = PetscMalloc1(numCellsNew*(dim+1), &cellsNew);CHKERRQ(ierr);
178f7bcbcbbSMatthew G. Knepley   pragmatic_get_elements(cellsNew);
17951ff0a1cSMatthew G. Knepley   ierr = DMPlexCreateFromCellListParallel(comm, dim, numCellsNew, numVerticesNew, numCornersNew, PETSC_TRUE, cellsNew, dim, coordsNew, NULL, dmNew);CHKERRQ(ierr);
180f7bcbcbbSMatthew G. Knepley   /* Read out boundary label */
181f7bcbcbbSMatthew G. Knepley   pragmatic_get_boundaryTags(&bdTags);
182f7bcbcbbSMatthew G. Knepley   ierr = DMCreateLabel(*dmNew, bdLabel ? bdLabelName : bdName);CHKERRQ(ierr);
183f7bcbcbbSMatthew G. Knepley   ierr = DMGetLabel(*dmNew, bdLabel ? bdLabelName : bdName, &bdLabelNew);CHKERRQ(ierr);
184f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd);CHKERRQ(ierr);
185f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd);CHKERRQ(ierr);
186f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd);CHKERRQ(ierr);
187f7bcbcbbSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
188f7bcbcbbSMatthew G. Knepley     /* Only for simplicial meshes */
189f7bcbcbbSMatthew G. Knepley     coff = (c-cStart)*(dim+1);
190f7bcbcbbSMatthew G. Knepley     for (d = 0; d < dim+1; ++d) {
191f7bcbcbbSMatthew G. Knepley       if (bdTags[coff+d]) {
192f7bcbcbbSMatthew G. Knepley         const PetscInt *faces;
193f7bcbcbbSMatthew G. Knepley         PetscInt        numFaces, vertices[3], p;
194f7bcbcbbSMatthew G. Knepley 
195f7bcbcbbSMatthew G. Knepley         /* Mark face opposite to this vertex */
196f7bcbcbbSMatthew G. Knepley         for (p = 0, v = 0; p < dim+1; ++p) if (p != d) {vertices[v] = cellsNew[coff+p] + vStart; ++v;}
197f7bcbcbbSMatthew G. Knepley         ierr = DMPlexGetFullJoin(*dmNew, dim, vertices, &numFaces, &faces);CHKERRQ(ierr);
198f7bcbcbbSMatthew 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]);
199f7bcbcbbSMatthew 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]);
200f7bcbcbbSMatthew G. Knepley         ierr = DMLabelSetValue(bdLabelNew, faces[0], bdTags[coff+d]);CHKERRQ(ierr);
201f7bcbcbbSMatthew G. Knepley         ierr = DMPlexRestoreJoin(*dmNew, dim, vertices, &numFaces, &faces);CHKERRQ(ierr);
202f7bcbcbbSMatthew G. Knepley       }
203f7bcbcbbSMatthew G. Knepley     }
204f7bcbcbbSMatthew G. Knepley   }
205f7bcbcbbSMatthew G. Knepley   /* Cleanup */
206f7bcbcbbSMatthew G. Knepley   switch (dim) {
207f7bcbcbbSMatthew G. Knepley   case 2: ierr = PetscFree2(xNew[0], xNew[1]);CHKERRQ(ierr);break;
208f7bcbcbbSMatthew G. Knepley   case 3: ierr = PetscFree3(xNew[0], xNew[1], xNew[2]);CHKERRQ(ierr);break;
209f7bcbcbbSMatthew G. Knepley   }
210f7bcbcbbSMatthew G. Knepley   ierr = PetscFree(cellsNew);CHKERRQ(ierr);
211f7bcbcbbSMatthew G. Knepley   ierr = PetscFree5(x, y, z, metric, cells);CHKERRQ(ierr);
212f7bcbcbbSMatthew G. Knepley   ierr = PetscFree2(bdFaces, bdFaceIds);CHKERRQ(ierr);
213f7bcbcbbSMatthew G. Knepley   ierr = PetscFree(coordsNew);CHKERRQ(ierr);
214f7bcbcbbSMatthew G. Knepley   pragmatic_finalize();
215f7bcbcbbSMatthew G. Knepley #else
216bc1d7e78SMatthew G. Knepley   SETERRQ(comm, PETSC_ERR_SUP, "Remeshing needs external package support.\nPlease reconfigure with --download-pragmatic.");
217f7bcbcbbSMatthew G. Knepley #endif
218f7bcbcbbSMatthew G. Knepley   PetscFunctionReturn(0);
219f7bcbcbbSMatthew G. Knepley }
2200bbe5d31Sbarral 
2210f7e110dSbarral /*@
2220f7e110dSbarral   DMPlexAdapt - Generates a mesh adapted to the specified metric field using the pragmatic library.
2230f7e110dSbarral 
2240f7e110dSbarral   Input Parameters:
2250f7e110dSbarral + dm - The DM object
2260f7e110dSbarral . metric - The metric to which the mesh is adapted, defined vertex-wise.
227379fadc5Sbarral - 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".
2280f7e110dSbarral 
2290f7e110dSbarral   Output Parameter:
230f7bcbcbbSMatthew G. Knepley . dmAdapt  - Pointer to the DM object containing the adapted mesh
2310f7e110dSbarral 
2320f7e110dSbarral   Level: advanced
2330f7e110dSbarral 
2340f7e110dSbarral .seealso: DMCoarsen(), DMRefine()
2350f7e110dSbarral @*/
236f7bcbcbbSMatthew G. Knepley PetscErrorCode DMPlexAdapt(DM dm, Vec metric, const char bdLabelName[], DM *dmAdapt)
2370bbe5d31Sbarral {
2382e62ab5aSMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex *) dm->data;
2390bbe5d31Sbarral   PetscErrorCode ierr;
2400bbe5d31Sbarral 
2410bbe5d31Sbarral   PetscFunctionBegin;
2422e62ab5aSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2432e62ab5aSMatthew G. Knepley   PetscValidHeaderSpecific(metric, VEC_CLASSID, 2);
2442e62ab5aSMatthew G. Knepley   PetscValidCharPointer(bdLabelName, 3);
2452e62ab5aSMatthew G. Knepley   PetscValidPointer(dmAdapt, 4);
2462e62ab5aSMatthew G. Knepley   ierr = DMPlexRemesh_Internal(dm, metric, bdLabelName, mesh->remeshBd, dmAdapt);CHKERRQ(ierr);
2470bbe5d31Sbarral   PetscFunctionReturn(0);
2480bbe5d31Sbarral }
249