xref: /petsc/src/dm/impls/plex/adaptors/parmmg/parmmgadapt.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
1 #include "../mmgcommon.h" /*I      "petscdmplex.h"   I*/
2 #include <parmmg/libparmmg.h>
3 
4 PETSC_EXTERN PetscErrorCode DMAdaptMetric_ParMmg_Plex(DM dm, Vec vertexMetric, DMLabel bdLabel, DMLabel rgLabel, DM *dmNew) {
5   MPI_Comm           comm;
6   const char        *bdName = "_boundary_";
7   const char        *rgName = "_regions_";
8   DM                 udm, cdm;
9   DMLabel            bdLabelNew, rgLabelNew;
10   const char        *bdLabelName, *rgLabelName;
11   IS                 globalVertexNum;
12   PetscSection       coordSection;
13   Vec                coordinates;
14   PetscSF            sf;
15   const PetscScalar *coords, *met;
16   PetscReal         *vertices, *metric, *verticesNew, *verticesNewLoc, gradationFactor, hausdorffNumber;
17   PetscInt          *cells, *cellsNew, *cellTags, *cellTagsNew, *verTags, *verTagsNew;
18   PetscInt          *bdFaces, *faceTags, *facesNew, *faceTagsNew;
19   PetscInt          *corners, *requiredCells, *requiredVer, *ridges, *requiredFaces;
20   PetscInt           cStart, cEnd, c, numCells, fStart, fEnd, f, numFaceTags, vStart, vEnd, v, numVertices;
21   PetscInt           dim, off, coff, maxConeSize, bdSize, i, j, k, Neq, verbosity, numIter;
22   PetscInt          *numVerInterfaces, *ngbRanks, *verNgbRank, *interfaces_lv, *interfaces_gv, *intOffset;
23   PetscInt           niranks, nrranks, numNgbRanks, numVerNgbRanksTotal, count, sliceSize, p, r, n, lv, gv;
24   PetscInt          *gv_new, *owners, *verticesNewSorted, pStart, pEnd;
25   PetscInt           numCellsNew, numVerticesNew, numCornersNew, numFacesNew, numVerticesNewLoc;
26   const PetscInt    *gV, *ioffset, *irootloc, *roffset, *rmine, *rremote;
27   PetscBool          flg = PETSC_FALSE, noInsert, noSwap, noMove, noSurf, isotropic, uniform;
28   const PetscMPIInt *iranks, *rranks;
29   PetscMPIInt        numProcs, rank;
30   PMMG_pParMesh      parmesh = NULL;
31 
32   PetscFunctionBegin;
33   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
34   PetscCallMPI(MPI_Comm_size(comm, &numProcs));
35   PetscCallMPI(MPI_Comm_rank(comm, &rank));
36   if (bdLabel) {
37     PetscCall(PetscObjectGetName((PetscObject)bdLabel, &bdLabelName));
38     PetscCall(PetscStrcmp(bdLabelName, bdName, &flg));
39     PetscCheck(!flg, comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName);
40   }
41   if (rgLabel) {
42     PetscCall(PetscObjectGetName((PetscObject)rgLabel, &rgLabelName));
43     PetscCall(PetscStrcmp(rgLabelName, rgName, &flg));
44     PetscCheck(!flg, comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for element tags", rgLabelName);
45   }
46 
47   /* Get mesh information */
48   PetscCall(DMGetDimension(dm, &dim));
49   PetscCheck(dim == 3, comm, PETSC_ERR_ARG_OUTOFRANGE, "ParMmg only works in 3D.");
50   Neq = (dim * (dim + 1)) / 2;
51   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
52   PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
53   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
54   PetscCall(DMPlexUninterpolate(dm, &udm));
55   PetscCall(DMPlexGetMaxSizes(udm, &maxConeSize, NULL));
56   numCells    = cEnd - cStart;
57   numVertices = vEnd - vStart;
58 
59   /* Get cell offsets */
60   PetscCall(PetscMalloc1(numCells * maxConeSize, &cells));
61   for (c = 0, coff = 0; c < numCells; ++c) {
62     const PetscInt *cone;
63     PetscInt        coneSize, cl;
64 
65     PetscCall(DMPlexGetConeSize(udm, c, &coneSize));
66     PetscCall(DMPlexGetCone(udm, c, &cone));
67     for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart + 1;
68   }
69 
70   /* Get vertex coordinate array */
71   PetscCall(DMGetCoordinateDM(dm, &cdm));
72   PetscCall(DMGetLocalSection(cdm, &coordSection));
73   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
74   PetscCall(VecGetArrayRead(coordinates, &coords));
75   PetscCall(PetscMalloc2(numVertices * Neq, &metric, dim * numVertices, &vertices));
76   for (v = 0; v < vEnd - vStart; ++v) {
77     PetscCall(PetscSectionGetOffset(coordSection, v + vStart, &off));
78     for (i = 0; i < dim; ++i) vertices[dim * v + i] = PetscRealPart(coords[off + i]);
79   }
80   PetscCall(VecRestoreArrayRead(coordinates, &coords));
81 
82   /* Get face tags */
83   if (!bdLabel) {
84     flg = PETSC_TRUE;
85     PetscCall(DMLabelCreate(PETSC_COMM_SELF, bdName, &bdLabel));
86     PetscCall(DMPlexMarkBoundaryFaces(dm, 1, bdLabel));
87   }
88   PetscCall(DMLabelGetBounds(bdLabel, &pStart, &pEnd));
89   for (f = pStart, bdSize = 0, numFaceTags = 0; f < pEnd; ++f) {
90     PetscBool hasPoint;
91     PetscInt *closure = NULL, closureSize, cl;
92 
93     PetscCall(DMLabelHasPoint(bdLabel, f, &hasPoint));
94     if ((!hasPoint) || (f < fStart) || (f >= fEnd)) continue;
95     numFaceTags++;
96 
97     PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
98     for (cl = 0; cl < closureSize * 2; cl += 2) {
99       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
100     }
101     PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
102   }
103   PetscCall(PetscMalloc2(bdSize, &bdFaces, numFaceTags, &faceTags));
104   for (f = pStart, bdSize = 0, numFaceTags = 0; f < pEnd; ++f) {
105     PetscBool hasPoint;
106     PetscInt *closure = NULL, closureSize, cl;
107 
108     PetscCall(DMLabelHasPoint(bdLabel, f, &hasPoint));
109     if ((!hasPoint) || (f < fStart) || (f >= fEnd)) continue;
110 
111     PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
112     for (cl = 0; cl < closureSize * 2; cl += 2) {
113       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart + 1;
114     }
115     PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
116     PetscCall(DMLabelGetValue(bdLabel, f, &faceTags[numFaceTags++]));
117   }
118 
119   /* Get cell tags */
120   PetscCall(PetscCalloc2(numVertices, &verTags, numCells, &cellTags));
121   if (rgLabel) {
122     for (c = cStart; c < cEnd; ++c) PetscCall(DMLabelGetValue(rgLabel, c, &cellTags[c]));
123   }
124 
125   /* Get metric */
126   PetscCall(VecViewFromOptions(vertexMetric, NULL, "-adapt_metric_view"));
127   PetscCall(VecGetArrayRead(vertexMetric, &met));
128   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
129   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
130   for (v = 0; v < (vEnd - vStart); ++v) {
131     for (i = 0, k = 0; i < dim; ++i) {
132       for (j = i; j < dim; ++j, ++k) {
133         if (isotropic) {
134           if (i == j) {
135             if (uniform) metric[Neq * v + k] = PetscRealPart(met[0]);
136             else metric[Neq * v + k] = PetscRealPart(met[v]);
137           } else metric[Neq * v + k] = 0.0;
138         } else metric[Neq * v + k] = PetscRealPart(met[dim * dim * v + dim * i + j]);
139       }
140     }
141   }
142   PetscCall(VecRestoreArrayRead(vertexMetric, &met));
143 
144   /* Build ParMMG communicators: the list of vertices between two partitions  */
145   niranks = nrranks = 0;
146   numNgbRanks       = 0;
147   if (numProcs > 1) {
148     PetscCall(DMGetPointSF(dm, &sf));
149     PetscCall(PetscSFSetUp(sf));
150     PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc));
151     PetscCall(PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote));
152     PetscCall(PetscCalloc1(numProcs, &numVerInterfaces));
153 
154     /* Count number of roots associated with each leaf */
155     for (r = 0; r < niranks; ++r) {
156       for (i = ioffset[r], count = 0; i < ioffset[r + 1]; ++i) {
157         if (irootloc[i] >= vStart && irootloc[i] < vEnd) count++;
158       }
159       numVerInterfaces[iranks[r]] += count;
160     }
161 
162     /* Count number of leaves associated with each root */
163     for (r = 0; r < nrranks; ++r) {
164       for (i = roffset[r], count = 0; i < roffset[r + 1]; ++i) {
165         if (rmine[i] >= vStart && rmine[i] < vEnd) count++;
166       }
167       numVerInterfaces[rranks[r]] += count;
168     }
169 
170     /* Count global number of ranks */
171     for (p = 0; p < numProcs; ++p) {
172       if (numVerInterfaces[p]) numNgbRanks++;
173     }
174 
175     /* Provide numbers of vertex interfaces */
176     PetscCall(PetscMalloc2(numNgbRanks, &ngbRanks, numNgbRanks, &verNgbRank));
177     for (p = 0, n = 0; p < numProcs; ++p) {
178       if (numVerInterfaces[p]) {
179         ngbRanks[n]   = p;
180         verNgbRank[n] = numVerInterfaces[p];
181         n++;
182       }
183     }
184     numVerNgbRanksTotal = 0;
185     for (p = 0; p < numNgbRanks; ++p) numVerNgbRanksTotal += verNgbRank[p];
186 
187     /* For each neighbor, fill in interface arrays */
188     PetscCall(PetscMalloc3(numVerNgbRanksTotal, &interfaces_lv, numVerNgbRanksTotal, &interfaces_gv, numNgbRanks + 1, &intOffset));
189     intOffset[0] = 0;
190     for (p = 0, r = 0, i = 0; p < numNgbRanks; ++p) {
191       intOffset[p + 1] = intOffset[p];
192 
193       /* Leaf case */
194       if (iranks && iranks[i] == ngbRanks[p]) {
195         /* Add the right slice of irootloc at the right place */
196         sliceSize = ioffset[i + 1] - ioffset[i];
197         for (j = 0, count = 0; j < sliceSize; ++j) {
198           PetscCheck(ioffset[i] + j < ioffset[niranks], comm, PETSC_ERR_ARG_OUTOFRANGE, "Leaf index %" PetscInt_FMT " out of range (expected < %" PetscInt_FMT ")", ioffset[i] + j, ioffset[niranks]);
199           v = irootloc[ioffset[i] + j];
200           if (v >= vStart && v < vEnd) {
201             PetscCheck(intOffset[p + 1] + count < numVerNgbRanksTotal, comm, PETSC_ERR_ARG_OUTOFRANGE, "Leaf interface index %" PetscInt_FMT " out of range (expected < %" PetscInt_FMT ")", intOffset[p + 1] + count, numVerNgbRanksTotal);
202             interfaces_lv[intOffset[p + 1] + count] = v - vStart;
203             count++;
204           }
205         }
206         intOffset[p + 1] += count;
207         i++;
208       }
209 
210       /* Root case */
211       if (rranks && rranks[r] == ngbRanks[p]) {
212         /* Add the right slice of rmine at the right place */
213         sliceSize = roffset[r + 1] - roffset[r];
214         for (j = 0, count = 0; j < sliceSize; ++j) {
215           PetscCheck(roffset[r] + j < roffset[nrranks], comm, PETSC_ERR_ARG_OUTOFRANGE, "Root index %" PetscInt_FMT " out of range (expected < %" PetscInt_FMT ")", roffset[r] + j, roffset[nrranks]);
216           v = rmine[roffset[r] + j];
217           if (v >= vStart && v < vEnd) {
218             PetscCheck(intOffset[p + 1] + count < numVerNgbRanksTotal, comm, PETSC_ERR_ARG_OUTOFRANGE, "Root interface index %" PetscInt_FMT " out of range (expected < %" PetscInt_FMT ")", intOffset[p + 1] + count, numVerNgbRanksTotal);
219             interfaces_lv[intOffset[p + 1] + count] = v - vStart;
220             count++;
221           }
222         }
223         intOffset[p + 1] += count;
224         r++;
225       }
226 
227       /* Check validity of offsets */
228       PetscCheck(intOffset[p + 1] == intOffset[p] + verNgbRank[p], comm, PETSC_ERR_ARG_OUTOFRANGE, "Missing offsets (expected %" PetscInt_FMT ", got %" PetscInt_FMT ")", intOffset[p] + verNgbRank[p], intOffset[p + 1]);
229     }
230     PetscCall(DMPlexGetVertexNumbering(udm, &globalVertexNum));
231     PetscCall(ISGetIndices(globalVertexNum, &gV));
232     for (i = 0; i < numVerNgbRanksTotal; ++i) {
233       v                = gV[interfaces_lv[i]];
234       interfaces_gv[i] = v < 0 ? -v - 1 : v;
235       interfaces_lv[i] += 1;
236       interfaces_gv[i] += 1;
237     }
238     PetscCall(ISRestoreIndices(globalVertexNum, &gV));
239     PetscCall(PetscFree(numVerInterfaces));
240   }
241   PetscCall(DMDestroy(&udm));
242 
243   /* Send the data to ParMmg and remesh */
244   PetscCall(DMPlexMetricNoInsertion(dm, &noInsert));
245   PetscCall(DMPlexMetricNoSwapping(dm, &noSwap));
246   PetscCall(DMPlexMetricNoMovement(dm, &noMove));
247   PetscCall(DMPlexMetricNoSurf(dm, &noSurf));
248   PetscCall(DMPlexMetricGetVerbosity(dm, &verbosity));
249   PetscCall(DMPlexMetricGetNumIterations(dm, &numIter));
250   PetscCall(DMPlexMetricGetGradationFactor(dm, &gradationFactor));
251   PetscCall(DMPlexMetricGetHausdorffNumber(dm, &hausdorffNumber));
252   PetscCallMMG_NONSTANDARD(PMMG_Init_parMesh(PMMG_ARG_start, PMMG_ARG_ppParMesh, &parmesh, PMMG_ARG_pMesh, PMMG_ARG_pMet, PMMG_ARG_dim, 3, PMMG_ARG_MPIComm, comm, PMMG_ARG_end));
253   PetscCallMMG_NONSTANDARD(PMMG_Set_meshSize(parmesh, numVertices, numCells, 0, numFaceTags, 0, 0));
254   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_APImode, PMMG_APIDISTRIB_nodes));
255   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_noinsert, noInsert));
256   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_noswap, noSwap));
257   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_nomove, noMove));
258   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_nosurf, noSurf));
259   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_verbose, verbosity));
260   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_globalNum, 1));
261   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_niter, numIter));
262   PetscCallMMG_NONSTANDARD(PMMG_Set_dparameter(parmesh, PMMG_DPARAM_hgrad, gradationFactor));
263   PetscCallMMG_NONSTANDARD(PMMG_Set_dparameter(parmesh, PMMG_DPARAM_hausd, hausdorffNumber));
264   PetscCallMMG_NONSTANDARD(PMMG_Set_vertices(parmesh, vertices, verTags));
265   PetscCallMMG_NONSTANDARD(PMMG_Set_tetrahedra(parmesh, cells, cellTags));
266   PetscCallMMG_NONSTANDARD(PMMG_Set_triangles(parmesh, bdFaces, faceTags));
267   PetscCallMMG_NONSTANDARD(PMMG_Set_metSize(parmesh, MMG5_Vertex, numVertices, MMG5_Tensor));
268   PetscCallMMG_NONSTANDARD(PMMG_Set_tensorMets(parmesh, metric));
269   PetscCallMMG_NONSTANDARD(PMMG_Set_numberOfNodeCommunicators(parmesh, numNgbRanks));
270   for (c = 0; c < numNgbRanks; ++c) {
271     PetscCallMMG_NONSTANDARD(PMMG_Set_ithNodeCommunicatorSize(parmesh, c, ngbRanks[c], intOffset[c + 1] - intOffset[c]));
272     PetscCallMMG_NONSTANDARD(PMMG_Set_ithNodeCommunicator_nodes(parmesh, c, &interfaces_lv[intOffset[c]], &interfaces_gv[intOffset[c]], 1));
273   }
274   PetscCallMMG(PMMG_parmmglib_distributed(parmesh));
275   PetscCall(PetscFree(cells));
276   PetscCall(PetscFree2(metric, vertices));
277   PetscCall(PetscFree2(bdFaces, faceTags));
278   PetscCall(PetscFree2(verTags, cellTags));
279   if (numProcs > 1) {
280     PetscCall(PetscFree2(ngbRanks, verNgbRank));
281     PetscCall(PetscFree3(interfaces_lv, interfaces_gv, intOffset));
282   }
283 
284   /* Retrieve mesh from Mmg */
285   numCornersNew = 4;
286   PetscCallMMG_NONSTANDARD(PMMG_Get_meshSize(parmesh, &numVerticesNew, &numCellsNew, 0, &numFacesNew, 0, 0));
287   PetscCall(PetscMalloc4(dim * numVerticesNew, &verticesNew, numVerticesNew, &verTagsNew, numVerticesNew, &corners, numVerticesNew, &requiredVer));
288   PetscCall(PetscMalloc3((dim + 1) * numCellsNew, &cellsNew, numCellsNew, &cellTagsNew, numCellsNew, &requiredCells));
289   PetscCall(PetscMalloc4(dim * numFacesNew, &facesNew, numFacesNew, &faceTagsNew, numFacesNew, &ridges, numFacesNew, &requiredFaces));
290   PetscCallMMG_NONSTANDARD(PMMG_Get_vertices(parmesh, verticesNew, verTagsNew, corners, requiredVer));
291   PetscCallMMG_NONSTANDARD(PMMG_Get_tetrahedra(parmesh, cellsNew, cellTagsNew, requiredCells));
292   PetscCallMMG_NONSTANDARD(PMMG_Get_triangles(parmesh, facesNew, faceTagsNew, requiredFaces));
293   PetscCall(PetscMalloc2(numVerticesNew, &owners, numVerticesNew, &gv_new));
294   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_globalNum, 1));
295   PetscCallMMG_NONSTANDARD(PMMG_Get_verticesGloNum(parmesh, gv_new, owners));
296   for (i = 0; i < dim * numFacesNew; ++i) facesNew[i] -= 1;
297   for (i = 0; i < (dim + 1) * numCellsNew; ++i) cellsNew[i] = gv_new[cellsNew[i] - 1] - 1;
298   for (i = 0, numVerticesNewLoc = 0; i < numVerticesNew; ++i) {
299     if (owners[i] == rank) numVerticesNewLoc++;
300   }
301   PetscCall(PetscMalloc2(numVerticesNewLoc * dim, &verticesNewLoc, numVerticesNew, &verticesNewSorted));
302   for (i = 0, c = 0; i < numVerticesNew; i++) {
303     if (owners[i] == rank) {
304       for (j = 0; j < dim; ++j) verticesNewLoc[dim * c + j] = verticesNew[dim * i + j];
305       c++;
306     }
307   }
308 
309   /* Reorder for consistency with DMPlex */
310   for (i = 0; i < numCellsNew; ++i) PetscCall(DMPlexInvertCell(DM_POLYTOPE_TETRAHEDRON, &cellsNew[4 * i]));
311 
312   /* Create new plex */
313   PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, dim, numCellsNew, numVerticesNewLoc, PETSC_DECIDE, numCornersNew, PETSC_TRUE, cellsNew, dim, verticesNewLoc, NULL, &verticesNewSorted, dmNew));
314   PetscCallMMG_NONSTANDARD(PMMG_Free_all(PMMG_ARG_start, PMMG_ARG_ppParMesh, &parmesh, PMMG_ARG_end));
315   PetscCall(PetscFree4(verticesNew, verTagsNew, corners, requiredVer));
316 
317   /* Get adapted mesh information */
318   PetscCall(DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd));
319   PetscCall(DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd));
320   PetscCall(DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd));
321 
322   /* Rebuild boundary label */
323   PetscCall(DMCreateLabel(*dmNew, flg ? bdName : bdLabelName));
324   PetscCall(DMGetLabel(*dmNew, flg ? bdName : bdLabelName, &bdLabelNew));
325   for (i = 0; i < numFacesNew; i++) {
326     PetscBool       hasTag = PETSC_FALSE;
327     PetscInt        numCoveredPoints, numFaces = 0, facePoints[3];
328     const PetscInt *coveredPoints = NULL;
329 
330     for (j = 0; j < dim; ++j) {
331       lv = facesNew[i * dim + j];
332       gv = gv_new[lv] - 1;
333       PetscCall(PetscFindInt(gv, numVerticesNew, verticesNewSorted, &lv));
334       facePoints[j] = lv + vStart;
335     }
336     PetscCall(DMPlexGetFullJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints));
337     for (j = 0; j < numCoveredPoints; ++j) {
338       if (coveredPoints[j] >= fStart && coveredPoints[j] < fEnd) {
339         numFaces++;
340         f = j;
341       }
342     }
343     PetscCheck(numFaces == 1, comm, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt_FMT " vertices cannot define more than 1 facet (%" PetscInt_FMT ")", dim, numFaces);
344     PetscCall(DMLabelHasStratum(bdLabel, faceTagsNew[i], &hasTag));
345     if (hasTag) PetscCall(DMLabelSetValue(bdLabelNew, coveredPoints[f], faceTagsNew[i]));
346     PetscCall(DMPlexRestoreJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints));
347   }
348   PetscCall(PetscFree4(facesNew, faceTagsNew, ridges, requiredFaces));
349   PetscCall(PetscFree2(owners, gv_new));
350   PetscCall(PetscFree2(verticesNewLoc, verticesNewSorted));
351   if (flg) PetscCall(DMLabelDestroy(&bdLabel));
352 
353   /* Rebuild cell labels */
354   PetscCall(DMCreateLabel(*dmNew, rgLabel ? rgLabelName : rgName));
355   PetscCall(DMGetLabel(*dmNew, rgLabel ? rgLabelName : rgName, &rgLabelNew));
356   for (c = cStart; c < cEnd; ++c) PetscCall(DMLabelSetValue(rgLabelNew, c, cellTagsNew[c - cStart]));
357   PetscCall(PetscFree3(cellsNew, cellTagsNew, requiredCells));
358 
359   PetscFunctionReturn(0);
360 }
361