xref: /petsc/src/dm/impls/plex/plexadapt.c (revision 1e1ea65d8de51fde77ce8a787efbef25e407badc)
1 #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2 #if defined(PETSC_HAVE_PRAGMATIC)
3 #include <pragmatic/cpragmatic.h>
4 #endif
5 
6 static PetscErrorCode DMPlexLabelToVolumeConstraint(DM dm, DMLabel adaptLabel, PetscInt cStart, PetscInt cEnd, PetscReal refRatio, PetscReal maxVolumes[])
7 {
8   PetscInt       dim, c;
9   PetscErrorCode ierr;
10 
11   PetscFunctionBegin;
12   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
13   refRatio = refRatio == PETSC_DEFAULT ? (PetscReal) ((PetscInt) 1 << dim) : refRatio;
14   for (c = cStart; c < cEnd; c++) {
15     PetscReal vol;
16     PetscInt  closureSize = 0, cl;
17     PetscInt *closure     = NULL;
18     PetscBool anyRefine   = PETSC_FALSE;
19     PetscBool anyCoarsen  = PETSC_FALSE;
20     PetscBool anyKeep     = PETSC_FALSE;
21 
22     ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, NULL, NULL);CHKERRQ(ierr);
23     maxVolumes[c - cStart] = vol;
24     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
25     for (cl = 0; cl < closureSize*2; cl += 2) {
26       const PetscInt point = closure[cl];
27       PetscInt       refFlag;
28 
29       ierr = DMLabelGetValue(adaptLabel, point, &refFlag);CHKERRQ(ierr);
30       switch (refFlag) {
31       case DM_ADAPT_REFINE:
32         anyRefine  = PETSC_TRUE;break;
33       case DM_ADAPT_COARSEN:
34         anyCoarsen = PETSC_TRUE;break;
35       case DM_ADAPT_KEEP:
36         anyKeep    = PETSC_TRUE;break;
37       case DM_ADAPT_DETERMINE:
38         break;
39       default:
40         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP, "DMPlex does not support refinement flag %D\n", refFlag);
41       }
42       if (anyRefine) break;
43     }
44     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
45     if (anyRefine) {
46       maxVolumes[c - cStart] = vol / refRatio;
47     } else if (anyKeep) {
48       maxVolumes[c - cStart] = vol;
49     } else if (anyCoarsen) {
50       maxVolumes[c - cStart] = vol * refRatio;
51     }
52   }
53   PetscFunctionReturn(0);
54 }
55 
56 static PetscErrorCode DMPlexLabelToMetricConstraint(DM dm, DMLabel adaptLabel, PetscInt cStart, PetscInt cEnd, PetscInt vStart, PetscInt vEnd, PetscReal refRatio, Vec *metricVec)
57 {
58   DM              udm, coordDM;
59   PetscSection    coordSection;
60   Vec             coordinates, mb, mx;
61   Mat             A;
62   PetscScalar    *metric, *eqns;
63   const PetscReal coarseRatio = refRatio == PETSC_DEFAULT ? PetscSqr(0.5) : 1/refRatio;
64   PetscInt        dim, Nv, Neq, c, v;
65   PetscErrorCode  ierr;
66 
67   PetscFunctionBegin;
68   ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr);
69   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
70   ierr = DMGetCoordinateDM(dm, &coordDM);CHKERRQ(ierr);
71   ierr = DMGetLocalSection(coordDM, &coordSection);CHKERRQ(ierr);
72   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
73   Nv   = vEnd - vStart;
74   ierr = VecCreateSeq(PETSC_COMM_SELF, Nv*PetscSqr(dim), metricVec);CHKERRQ(ierr);
75   ierr = VecGetArray(*metricVec, &metric);CHKERRQ(ierr);
76   Neq  = (dim*(dim+1))/2;
77   ierr = PetscMalloc1(PetscSqr(Neq), &eqns);CHKERRQ(ierr);
78   ierr = MatCreateSeqDense(PETSC_COMM_SELF, Neq, Neq, eqns, &A);CHKERRQ(ierr);
79   ierr = MatCreateVecs(A, &mx, &mb);CHKERRQ(ierr);
80   ierr = VecSet(mb, 1.0);CHKERRQ(ierr);
81   for (c = cStart; c < cEnd; ++c) {
82     const PetscScalar *sol;
83     PetscScalar       *cellCoords = NULL;
84     PetscReal          e[3], vol;
85     const PetscInt    *cone;
86     PetscInt           coneSize, cl, i, j, d, r;
87 
88     ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, c, NULL, &cellCoords);CHKERRQ(ierr);
89     /* Only works for simplices */
90     for (i = 0, r = 0; i < dim+1; ++i) {
91       for (j = 0; j < i; ++j, ++r) {
92         for (d = 0; d < dim; ++d) e[d] = PetscRealPart(cellCoords[i*dim+d] - cellCoords[j*dim+d]);
93         /* FORTRAN ORDERING */
94         switch (dim) {
95         case 2:
96           eqns[0*Neq+r] = PetscSqr(e[0]);
97           eqns[1*Neq+r] = 2.0*e[0]*e[1];
98           eqns[2*Neq+r] = PetscSqr(e[1]);
99           break;
100         case 3:
101           eqns[0*Neq+r] = PetscSqr(e[0]);
102           eqns[1*Neq+r] = 2.0*e[0]*e[1];
103           eqns[2*Neq+r] = 2.0*e[0]*e[2];
104           eqns[3*Neq+r] = PetscSqr(e[1]);
105           eqns[4*Neq+r] = 2.0*e[1]*e[2];
106           eqns[5*Neq+r] = PetscSqr(e[2]);
107           break;
108         }
109       }
110     }
111     ierr = MatSetUnfactored(A);CHKERRQ(ierr);
112     ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, c, NULL, &cellCoords);CHKERRQ(ierr);
113     ierr = MatLUFactor(A, NULL, NULL, NULL);CHKERRQ(ierr);
114     ierr = MatSolve(A, mb, mx);CHKERRQ(ierr);
115     ierr = VecGetArrayRead(mx, &sol);CHKERRQ(ierr);
116     ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, NULL, NULL);CHKERRQ(ierr);
117     ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr);
118     ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr);
119     for (cl = 0; cl < coneSize; ++cl) {
120       const PetscInt v = cone[cl] - vStart;
121 
122       if (dim == 2) {
123         metric[v*4+0] += vol*coarseRatio*sol[0];
124         metric[v*4+1] += vol*coarseRatio*sol[1];
125         metric[v*4+2] += vol*coarseRatio*sol[1];
126         metric[v*4+3] += vol*coarseRatio*sol[2];
127       } else {
128         metric[v*9+0] += vol*coarseRatio*sol[0];
129         metric[v*9+1] += vol*coarseRatio*sol[1];
130         metric[v*9+3] += vol*coarseRatio*sol[1];
131         metric[v*9+2] += vol*coarseRatio*sol[2];
132         metric[v*9+6] += vol*coarseRatio*sol[2];
133         metric[v*9+4] += vol*coarseRatio*sol[3];
134         metric[v*9+5] += vol*coarseRatio*sol[4];
135         metric[v*9+7] += vol*coarseRatio*sol[4];
136         metric[v*9+8] += vol*coarseRatio*sol[5];
137       }
138     }
139     ierr = VecRestoreArrayRead(mx, &sol);CHKERRQ(ierr);
140   }
141   for (v = 0; v < Nv; ++v) {
142     const PetscInt *support;
143     PetscInt        supportSize, s;
144     PetscReal       vol, totVol = 0.0;
145 
146     ierr = DMPlexGetSupport(udm, v+vStart, &support);CHKERRQ(ierr);
147     ierr = DMPlexGetSupportSize(udm, v+vStart, &supportSize);CHKERRQ(ierr);
148     for (s = 0; s < supportSize; ++s) {ierr = DMPlexComputeCellGeometryFVM(dm, support[s], &vol, NULL, NULL);CHKERRQ(ierr); totVol += vol;}
149     for (s = 0; s < PetscSqr(dim); ++s) metric[v*PetscSqr(dim)+s] /= totVol;
150   }
151   ierr = PetscFree(eqns);CHKERRQ(ierr);
152   ierr = VecRestoreArray(*metricVec, &metric);CHKERRQ(ierr);
153   ierr = VecDestroy(&mx);CHKERRQ(ierr);
154   ierr = VecDestroy(&mb);CHKERRQ(ierr);
155   ierr = MatDestroy(&A);CHKERRQ(ierr);
156   ierr = DMDestroy(&udm);CHKERRQ(ierr);
157   PetscFunctionReturn(0);
158 }
159 
160 /*
161    Contains the list of registered DMPlexGenerators routines
162 */
163 extern PlexGeneratorFunctionList DMPlexGenerateList;
164 
165 PetscErrorCode DMPlexRefine_Internal(DM dm, DMLabel adaptLabel, DM *dmRefined)
166 {
167   PlexGeneratorFunctionList fl;
168   PetscErrorCode          (*refine)(DM,PetscReal*,DM*);
169   PetscErrorCode          (*adapt)(DM,DMLabel,DM*);
170   PetscErrorCode          (*refinementFunc)(const PetscReal [], PetscReal *);
171   char                      genname[PETSC_MAX_PATH_LEN], *name = NULL;
172   PetscReal                 refinementLimit;
173   PetscReal                *maxVolumes;
174   PetscInt                  dim, cStart, cEnd, c;
175   PetscBool                 flg, flg2, localized;
176   PetscErrorCode            ierr;
177 
178   PetscFunctionBegin;
179   ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr);
180   ierr = DMPlexGetRefinementLimit(dm, &refinementLimit);CHKERRQ(ierr);
181   ierr = DMPlexGetRefinementFunction(dm, &refinementFunc);CHKERRQ(ierr);
182   if (refinementLimit == 0.0 && !refinementFunc && !adaptLabel) PetscFunctionReturn(0);
183   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
184   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
185   ierr = PetscOptionsGetString(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_adaptor", genname, sizeof(genname), &flg);CHKERRQ(ierr);
186   if (flg) name = genname;
187   else {
188     ierr = PetscOptionsGetString(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_generator", genname, sizeof(genname), &flg2);CHKERRQ(ierr);
189     if (flg2) name = genname;
190   }
191 
192   fl = DMPlexGenerateList;
193   if (name) {
194     while (fl) {
195       ierr = PetscStrcmp(fl->name,name,&flg);CHKERRQ(ierr);
196       if (flg) {
197         refine = fl->refine;
198         adapt  = fl->adaptlabel;
199         goto gotit;
200       }
201       fl = fl->next;
202     }
203     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Grid refiner %s not registered",name);
204   } else {
205     while (fl) {
206       if (dim-1 == fl->dim) {
207         refine = fl->refine;
208         adapt  = fl->adaptlabel;
209         goto gotit;
210       }
211       fl = fl->next;
212     }
213     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No grid refiner of dimension %D registered",dim);
214   }
215 
216   gotit:
217   switch (dim) {
218     case 2:
219     case 3:
220       if (adapt) {
221         ierr = (*adapt)(dm, adaptLabel, dmRefined);CHKERRQ(ierr);
222       } else {
223         ierr = PetscMalloc1(cEnd - cStart, &maxVolumes);CHKERRQ(ierr);
224         if (adaptLabel) {
225           ierr = DMPlexLabelToVolumeConstraint(dm, adaptLabel, cStart, cEnd, PETSC_DEFAULT, maxVolumes);CHKERRQ(ierr);
226         } else if (refinementFunc) {
227           for (c = cStart; c < cEnd; ++c) {
228             PetscReal vol, centroid[3];
229 
230             ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL);CHKERRQ(ierr);
231             ierr = (*refinementFunc)(centroid, &maxVolumes[c-cStart]);CHKERRQ(ierr);
232           }
233         } else {
234           for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit;
235         }
236         ierr = (*refine)(dm, maxVolumes, dmRefined);CHKERRQ(ierr);
237         ierr = PetscFree(maxVolumes);CHKERRQ(ierr);
238       }
239       break;
240     default: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Mesh refinement in dimension %D is not supported.", dim);
241   }
242   ierr = DMCopyDisc(dm, *dmRefined);CHKERRQ(ierr);
243   if (localized) {ierr = DMLocalizeCoordinates(*dmRefined);CHKERRQ(ierr);}
244   PetscFunctionReturn(0);
245 }
246 
247 PetscErrorCode DMPlexCoarsen_Internal(DM dm, DMLabel adaptLabel, DM *dmCoarsened)
248 {
249   Vec            metricVec;
250   PetscInt       cStart, cEnd, vStart, vEnd;
251   DMLabel        bdLabel = NULL;
252   char           bdLabelName[PETSC_MAX_PATH_LEN];
253   PetscBool      localized, flg;
254   PetscErrorCode ierr;
255 
256   PetscFunctionBegin;
257   ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr);
258   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
259   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
260   ierr = DMPlexLabelToMetricConstraint(dm, adaptLabel, cStart, cEnd, vStart, vEnd, PETSC_DEFAULT, &metricVec);CHKERRQ(ierr);
261   ierr = PetscOptionsGetString(NULL, dm->hdr.prefix, "-dm_plex_coarsen_bd_label", bdLabelName, sizeof(bdLabelName), &flg);CHKERRQ(ierr);
262   if (flg) {ierr = DMGetLabel(dm, bdLabelName, &bdLabel);CHKERRQ(ierr);}
263   ierr = DMAdaptMetric_Plex(dm, metricVec, bdLabel, dmCoarsened);CHKERRQ(ierr);
264   ierr = VecDestroy(&metricVec);CHKERRQ(ierr);
265   if (localized) {ierr = DMLocalizeCoordinates(*dmCoarsened);CHKERRQ(ierr);}
266   PetscFunctionReturn(0);
267 }
268 
269 PetscErrorCode DMAdaptLabel_Plex(DM dm, DMLabel adaptLabel, DM *dmAdapted)
270 {
271   IS              flagIS;
272   const PetscInt *flags;
273   PetscInt        defFlag, minFlag, maxFlag, numFlags, f;
274   PetscErrorCode  ierr;
275 
276   PetscFunctionBegin;
277   ierr = DMLabelGetDefaultValue(adaptLabel, &defFlag);CHKERRQ(ierr);
278   minFlag = defFlag;
279   maxFlag = defFlag;
280   ierr = DMLabelGetValueIS(adaptLabel, &flagIS);CHKERRQ(ierr);
281   ierr = ISGetLocalSize(flagIS, &numFlags);CHKERRQ(ierr);
282   ierr = ISGetIndices(flagIS, &flags);CHKERRQ(ierr);
283   for (f = 0; f < numFlags; ++f) {
284     const PetscInt flag = flags[f];
285 
286     minFlag = PetscMin(minFlag, flag);
287     maxFlag = PetscMax(maxFlag, flag);
288   }
289   ierr = ISRestoreIndices(flagIS, &flags);CHKERRQ(ierr);
290   ierr = ISDestroy(&flagIS);CHKERRQ(ierr);
291   {
292     PetscInt minMaxFlag[2], minMaxFlagGlobal[2];
293 
294     minMaxFlag[0] =  minFlag;
295     minMaxFlag[1] = -maxFlag;
296     ierr = MPI_Allreduce(minMaxFlag, minMaxFlagGlobal, 2, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm));CHKERRMPI(ierr);
297     minFlag =  minMaxFlagGlobal[0];
298     maxFlag = -minMaxFlagGlobal[1];
299   }
300   if (minFlag == maxFlag) {
301     switch (minFlag) {
302     case DM_ADAPT_DETERMINE:
303       *dmAdapted = NULL;break;
304     case DM_ADAPT_REFINE:
305       ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE);CHKERRQ(ierr);
306       ierr = DMRefine(dm, MPI_COMM_NULL, dmAdapted);CHKERRQ(ierr);break;
307     case DM_ADAPT_COARSEN:
308       ierr = DMCoarsen(dm, MPI_COMM_NULL, dmAdapted);CHKERRQ(ierr);break;
309     default:
310       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP,"DMPlex does not support refinement flag %D\n", minFlag);
311     }
312   } else {
313     ierr = DMPlexSetRefinementUniform(dm, PETSC_FALSE);CHKERRQ(ierr);
314     ierr = DMPlexRefine_Internal(dm, adaptLabel, dmAdapted);CHKERRQ(ierr);
315   }
316   PetscFunctionReturn(0);
317 }
318 
319 /*
320   DMAdaptMetric_Plex - Generates a new mesh conforming to a metric field.
321 
322   Input Parameters:
323 + dm - The DM object
324 . vertexMetric - The metric to which the mesh is adapted, defined vertex-wise in a LOCAL vector
325 - bdLabel - Label for boundary tags which are preserved in dmNew, or NULL. Should not be named "_boundary_".
326 
327   Output Parameter:
328 . dmNew  - the new DM
329 
330   Level: advanced
331 
332 .seealso: DMCoarsen(), DMRefine()
333 */
334 PetscErrorCode DMAdaptMetric_Plex(DM dm, Vec vertexMetric, DMLabel bdLabel, DM *dmNew)
335 {
336 #if defined(PETSC_HAVE_PRAGMATIC)
337   MPI_Comm           comm;
338   const char        *bdName = "_boundary_";
339 #if 0
340   DM                 odm = dm;
341 #endif
342   DM                 udm, cdm;
343   DMLabel            bdLabelFull;
344   const char        *bdLabelName;
345   IS                 bdIS, globalVertexNum;
346   PetscSection       coordSection;
347   Vec                coordinates;
348   const PetscScalar *coords, *met;
349   const PetscInt    *bdFacesFull, *gV;
350   PetscInt          *bdFaces, *bdFaceIds, *l2gv;
351   PetscReal         *x, *y, *z, *metric;
352   PetscInt          *cells;
353   PetscInt           dim, cStart, cEnd, numCells, c, coff, vStart, vEnd, numVertices, numLocVertices, v;
354   PetscInt           off, maxConeSize, numBdFaces, f, bdSize;
355   PetscBool          flg;
356   DMLabel            bdLabelNew;
357   PetscReal         *coordsNew;
358   PetscInt          *bdTags;
359   PetscReal         *xNew[3] = {NULL, NULL, NULL};
360   PetscInt          *cellsNew;
361   PetscInt           d, numCellsNew, numVerticesNew;
362   PetscInt           numCornersNew, fStart, fEnd;
363   PetscMPIInt        numProcs;
364   PetscErrorCode     ierr;
365 
366   PetscFunctionBegin;
367   /* Check for FEM adjacency flags */
368   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
369   ierr = MPI_Comm_size(comm, &numProcs);CHKERRMPI(ierr);
370   if (bdLabel) {
371     ierr = PetscObjectGetName((PetscObject) bdLabel, &bdLabelName);CHKERRQ(ierr);
372     ierr = PetscStrcmp(bdLabelName, bdName, &flg);CHKERRQ(ierr);
373     if (flg) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName);
374   }
375   /* Add overlap for Pragmatic */
376 #if 0
377   /* Check for overlap by looking for cell in the SF */
378   if (!overlapped) {
379     ierr = DMPlexDistributeOverlap(odm, 1, NULL, &dm);CHKERRQ(ierr);
380     if (!dm) {dm = odm; ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);}
381   }
382 #endif
383   /* Get mesh information */
384   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
385   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
386   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
387   ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr);
388   ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr);
389   numCells    = cEnd - cStart;
390   numVertices = vEnd - vStart;
391   ierr = PetscCalloc5(numVertices, &x, numVertices, &y, numVertices, &z, numVertices*PetscSqr(dim), &metric, numCells*maxConeSize, &cells);CHKERRQ(ierr);
392   for (c = 0, coff = 0; c < numCells; ++c) {
393     const PetscInt *cone;
394     PetscInt        coneSize, cl;
395 
396     ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr);
397     ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr);
398     for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart;
399   }
400   ierr = PetscCalloc1(numVertices, &l2gv);CHKERRQ(ierr);
401   ierr = DMPlexGetVertexNumbering(udm, &globalVertexNum);CHKERRQ(ierr);
402   ierr = ISGetIndices(globalVertexNum, &gV);CHKERRQ(ierr);
403   for (v = 0, numLocVertices = 0; v < numVertices; ++v) {
404     if (gV[v] >= 0) ++numLocVertices;
405     l2gv[v] = gV[v] < 0 ? -(gV[v]+1) : gV[v];
406   }
407   ierr = ISRestoreIndices(globalVertexNum, &gV);CHKERRQ(ierr);
408   ierr = DMDestroy(&udm);CHKERRQ(ierr);
409   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
410   ierr = DMGetLocalSection(cdm, &coordSection);CHKERRQ(ierr);
411   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
412   ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
413   for (v = vStart; v < vEnd; ++v) {
414     ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
415     x[v-vStart] = PetscRealPart(coords[off+0]);
416     if (dim > 1) y[v-vStart] = PetscRealPart(coords[off+1]);
417     if (dim > 2) z[v-vStart] = PetscRealPart(coords[off+2]);
418   }
419   ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
420   /* Get boundary mesh */
421   ierr = DMLabelCreate(PETSC_COMM_SELF, bdName, &bdLabelFull);CHKERRQ(ierr);
422   ierr = DMPlexMarkBoundaryFaces(dm, 1, bdLabelFull);CHKERRQ(ierr);
423   ierr = DMLabelGetStratumIS(bdLabelFull, 1, &bdIS);CHKERRQ(ierr);
424   ierr = DMLabelGetStratumSize(bdLabelFull, 1, &numBdFaces);CHKERRQ(ierr);
425   ierr = ISGetIndices(bdIS, &bdFacesFull);CHKERRQ(ierr);
426   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
427     PetscInt *closure = NULL;
428     PetscInt  closureSize, cl;
429 
430     ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
431     for (cl = 0; cl < closureSize*2; cl += 2) {
432       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
433     }
434     ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
435   }
436   ierr = PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds);CHKERRQ(ierr);
437   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
438     PetscInt *closure = NULL;
439     PetscInt  closureSize, cl;
440 
441     ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
442     for (cl = 0; cl < closureSize*2; cl += 2) {
443       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart;
444     }
445     ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
446     if (bdLabel) {ierr = DMLabelGetValue(bdLabel, bdFacesFull[f], &bdFaceIds[f]);CHKERRQ(ierr);}
447     else         {bdFaceIds[f] = 1;}
448   }
449   ierr = ISDestroy(&bdIS);CHKERRQ(ierr);
450   ierr = DMLabelDestroy(&bdLabelFull);CHKERRQ(ierr);
451   /* Get metric */
452   ierr = VecViewFromOptions(vertexMetric, NULL, "-adapt_metric_view");CHKERRQ(ierr);
453   ierr = VecGetArrayRead(vertexMetric, &met);CHKERRQ(ierr);
454   for (v = 0; v < (vEnd-vStart)*PetscSqr(dim); ++v) metric[v] = PetscRealPart(met[v]);
455   ierr = VecRestoreArrayRead(vertexMetric, &met);CHKERRQ(ierr);
456 #if 0
457   /* Destroy overlap mesh */
458   ierr = DMDestroy(&dm);CHKERRQ(ierr);
459 #endif
460   /* Create new mesh */
461   switch (dim) {
462   case 2:
463     pragmatic_2d_mpi_init(&numVertices, &numCells, cells, x, y, l2gv, numLocVertices, comm);break;
464   case 3:
465     pragmatic_3d_mpi_init(&numVertices, &numCells, cells, x, y, z, l2gv, numLocVertices, comm);break;
466   default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %D", dim);
467   }
468   pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds);
469   pragmatic_set_metric(metric);
470   pragmatic_adapt(((DM_Plex *) dm->data)->remeshBd ? 1 : 0);
471   ierr = PetscFree(l2gv);CHKERRQ(ierr);
472   /* Read out mesh */
473   pragmatic_get_info_mpi(&numVerticesNew, &numCellsNew);
474   ierr = PetscMalloc1(numVerticesNew*dim, &coordsNew);CHKERRQ(ierr);
475   switch (dim) {
476   case 2:
477     numCornersNew = 3;
478     ierr = PetscMalloc2(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1]);CHKERRQ(ierr);
479     pragmatic_get_coords_2d_mpi(xNew[0], xNew[1]);
480     break;
481   case 3:
482     numCornersNew = 4;
483     ierr = PetscMalloc3(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1], numVerticesNew, &xNew[2]);CHKERRQ(ierr);
484     pragmatic_get_coords_3d_mpi(xNew[0], xNew[1], xNew[2]);
485     break;
486   default:
487     SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %D", dim);
488   }
489   for (v = 0; v < numVerticesNew; ++v) {for (d = 0; d < dim; ++d) coordsNew[v*dim+d] = xNew[d][v];}
490   ierr = PetscMalloc1(numCellsNew*(dim+1), &cellsNew);CHKERRQ(ierr);
491   pragmatic_get_elements(cellsNew);
492   ierr = DMPlexCreateFromCellListParallelPetsc(comm, dim, numCellsNew, numVerticesNew, PETSC_DECIDE, numCornersNew, PETSC_TRUE, cellsNew, dim, coordsNew, NULL, dmNew);CHKERRQ(ierr);
493   /* Read out boundary label */
494   pragmatic_get_boundaryTags(&bdTags);
495   ierr = DMCreateLabel(*dmNew, bdLabel ? bdLabelName : bdName);CHKERRQ(ierr);
496   ierr = DMGetLabel(*dmNew, bdLabel ? bdLabelName : bdName, &bdLabelNew);CHKERRQ(ierr);
497   ierr = DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd);CHKERRQ(ierr);
498   ierr = DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd);CHKERRQ(ierr);
499   ierr = DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd);CHKERRQ(ierr);
500   for (c = cStart; c < cEnd; ++c) {
501     /* Only for simplicial meshes */
502     coff = (c-cStart)*(dim+1);
503     /* d is the local cell number of the vertex opposite to the face we are marking */
504     for (d = 0; d < dim+1; ++d) {
505       if (bdTags[coff+d]) {
506         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 */
507         const PetscInt *cone;
508 
509         /* Mark face opposite to this vertex: This pattern is specified in DMPlexGetRawFaces_Internal() */
510         ierr = DMPlexGetCone(*dmNew, c, &cone);CHKERRQ(ierr);
511         ierr = DMLabelSetValue(bdLabelNew, cone[perm[dim][d]], bdTags[coff+d]);CHKERRQ(ierr);
512       }
513     }
514   }
515   /* Cleanup */
516   switch (dim) {
517   case 2: ierr = PetscFree2(xNew[0], xNew[1]);CHKERRQ(ierr);break;
518   case 3: ierr = PetscFree3(xNew[0], xNew[1], xNew[2]);CHKERRQ(ierr);break;
519   }
520   ierr = PetscFree(cellsNew);CHKERRQ(ierr);
521   ierr = PetscFree5(x, y, z, metric, cells);CHKERRQ(ierr);
522   ierr = PetscFree2(bdFaces, bdFaceIds);CHKERRQ(ierr);
523   ierr = PetscFree(coordsNew);CHKERRQ(ierr);
524   pragmatic_finalize();
525   PetscFunctionReturn(0);
526 #else
527   PetscFunctionBegin;
528   SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Remeshing needs external package support.\nPlease reconfigure with --download-pragmatic.");
529 #endif
530 }
531