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