xref: /petsc/src/dm/impls/plex/plexadapt.c (revision 5b6bfdb9644f185dbf5e5a09b808ec241507e1e7)
1 #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2 #ifdef 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);break;
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 = DMGetDefaultSection(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 PetscFunctionList DMPlexGenerateList;
164 
165 struct _n_PetscFunctionList {
166   PetscErrorCode    (*generate)(DM, PetscBool, DM*);
167   PetscErrorCode    (*refine)(DM,PetscReal*, DM*);
168   char              *name;               /* string to identify routine */
169   PetscInt          dim;
170   PetscFunctionList next;                /* next pointer */
171 };
172 
173 PetscErrorCode DMPlexRefine_Internal(DM dm, DMLabel adaptLabel, DM *dmRefined)
174 {
175   PetscErrorCode    (*refinementFunc)(const PetscReal [], PetscReal *);
176   PetscReal         refinementLimit;
177   PetscInt          dim, cStart, cEnd;
178   char              genname[1024], *name = NULL;
179   PetscBool         flg, localized;
180   PetscErrorCode    ierr;
181   PetscErrorCode    (*refine)(DM,PetscReal*,DM*);
182   PetscFunctionList fl;
183   PetscReal         *maxVolumes;
184   PetscInt          c;
185 
186   PetscFunctionBegin;
187   ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr);
188   ierr = DMPlexGetRefinementLimit(dm, &refinementLimit);CHKERRQ(ierr);
189   ierr = DMPlexGetRefinementFunction(dm, &refinementFunc);CHKERRQ(ierr);
190   if (refinementLimit == 0.0 && !refinementFunc && !adaptLabel) PetscFunctionReturn(0);
191   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
192   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
193   ierr = PetscOptionsGetString(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_generator", genname, 1024, &flg);CHKERRQ(ierr);
194   if (flg) name = genname;
195 
196   fl = DMPlexGenerateList;
197   if (name) {
198     while (fl) {
199       ierr = PetscStrcmp(fl->name,name,&flg);CHKERRQ(ierr);
200       if (flg) {
201         refine = fl->refine;
202         goto gotit;
203       }
204       fl = fl->next;
205     }
206     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Grid refiner %g not registered",name);CHKERRQ(ierr);
207   } else {
208     while (fl) {
209       if (dim-1 == fl->dim) {
210         refine = fl->refine;
211         goto gotit;
212       }
213       fl = fl->next;
214     }
215     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No grid refiner of dimension %D registered",dim);CHKERRQ(ierr);
216   }
217 
218   gotit: switch (dim) {
219   case 2:
220       ierr = PetscMalloc1(cEnd - cStart, &maxVolumes);CHKERRQ(ierr);
221       if (adaptLabel) {
222         ierr = DMPlexLabelToVolumeConstraint(dm, adaptLabel, cStart, cEnd, PETSC_DEFAULT, maxVolumes);CHKERRQ(ierr);
223       } else if (refinementFunc) {
224         for (c = cStart; c < cEnd; ++c) {
225           PetscReal vol, centroid[3];
226           PetscReal maxVol;
227 
228           ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL);CHKERRQ(ierr);
229           ierr = (*refinementFunc)(centroid, &maxVol);CHKERRQ(ierr);
230           maxVolumes[c - cStart] = (double) maxVol;
231         }
232       } else {
233         for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit;
234       }
235       ierr = (*refine)(dm, maxVolumes, dmRefined);CHKERRQ(ierr);
236       ierr = PetscFree(maxVolumes);CHKERRQ(ierr);
237     break;
238   case 3:
239       ierr = PetscMalloc1(cEnd - cStart, &maxVolumes);CHKERRQ(ierr);
240       if (adaptLabel) {
241         ierr = DMPlexLabelToVolumeConstraint(dm, adaptLabel, cStart, cEnd, PETSC_DEFAULT, maxVolumes);CHKERRQ(ierr);
242       } else if (refinementFunc) {
243         for (c = cStart; c < cEnd; ++c) {
244           PetscReal vol, centroid[3];
245 
246           ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL);CHKERRQ(ierr);
247           ierr = (*refinementFunc)(centroid, &maxVolumes[c-cStart]);CHKERRQ(ierr);
248         }
249       } else {
250         for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit;
251       }
252       ierr = (*refine)(dm, maxVolumes, dmRefined);CHKERRQ(ierr);
253       ierr = PetscFree(maxVolumes);CHKERRQ(ierr);
254     break;
255   default:
256     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Mesh refinement in dimension %d is not supported.", dim);
257   }
258   ierr = DMCopyBoundary(dm, *dmRefined);CHKERRQ(ierr);
259   if (localized) {ierr = DMLocalizeCoordinates(*dmRefined);CHKERRQ(ierr);}
260   PetscFunctionReturn(0);
261 }
262 
263 PetscErrorCode DMPlexCoarsen_Internal(DM dm, DMLabel adaptLabel, DM *dmCoarsened)
264 {
265   Vec            metricVec;
266   PetscInt       cStart, cEnd, vStart, vEnd;
267   DMLabel        bdLabel = NULL;
268   char           bdLabelName[PETSC_MAX_PATH_LEN];
269   PetscBool      localized, flg;
270   PetscErrorCode ierr;
271 
272   PetscFunctionBegin;
273   ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr);
274   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
275   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
276   ierr = DMPlexLabelToMetricConstraint(dm, adaptLabel, cStart, cEnd, vStart, vEnd, PETSC_DEFAULT, &metricVec);CHKERRQ(ierr);
277   ierr = PetscOptionsGetString(NULL, dm->hdr.prefix, "-dm_plex_coarsen_bd_label", bdLabelName, PETSC_MAX_PATH_LEN-1, &flg);CHKERRQ(ierr);
278   if (flg) {ierr = DMGetLabel(dm, bdLabelName, &bdLabel);CHKERRQ(ierr);}
279   ierr = DMAdaptMetric_Plex(dm, metricVec, bdLabel, dmCoarsened);CHKERRQ(ierr);
280   ierr = VecDestroy(&metricVec);CHKERRQ(ierr);
281   if (localized) {ierr = DMLocalizeCoordinates(*dmCoarsened);CHKERRQ(ierr);}
282   PetscFunctionReturn(0);
283 }
284 
285 PetscErrorCode DMAdaptLabel_Plex(DM dm, DMLabel adaptLabel, DM *dmAdapted)
286 {
287   IS              flagIS;
288   const PetscInt *flags;
289   PetscInt        defFlag, minFlag, maxFlag, numFlags, f;
290   PetscErrorCode  ierr;
291 
292   PetscFunctionBegin;
293   ierr = DMLabelGetDefaultValue(adaptLabel, &defFlag);CHKERRQ(ierr);
294   minFlag = defFlag;
295   maxFlag = defFlag;
296   ierr = DMLabelGetValueIS(adaptLabel, &flagIS);CHKERRQ(ierr);
297   ierr = ISGetLocalSize(flagIS, &numFlags);CHKERRQ(ierr);
298   ierr = ISGetIndices(flagIS, &flags);CHKERRQ(ierr);
299   for (f = 0; f < numFlags; ++f) {
300     const PetscInt flag = flags[f];
301 
302     minFlag = PetscMin(minFlag, flag);
303     maxFlag = PetscMax(maxFlag, flag);
304   }
305   ierr = ISRestoreIndices(flagIS, &flags);CHKERRQ(ierr);
306   ierr = ISDestroy(&flagIS);CHKERRQ(ierr);
307   {
308     PetscInt minMaxFlag[2], minMaxFlagGlobal[2];
309 
310     minMaxFlag[0] =  minFlag;
311     minMaxFlag[1] = -maxFlag;
312     ierr = MPI_Allreduce(minMaxFlag, minMaxFlagGlobal, 2, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
313     minFlag =  minMaxFlagGlobal[0];
314     maxFlag = -minMaxFlagGlobal[1];
315   }
316   if (minFlag == maxFlag) {
317     switch (minFlag) {
318     case DM_ADAPT_DETERMINE:
319       *dmAdapted = NULL;break;
320     case DM_ADAPT_REFINE:
321       ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE);CHKERRQ(ierr);
322       ierr = DMRefine(dm, MPI_COMM_NULL, dmAdapted);CHKERRQ(ierr);break;
323     case DM_ADAPT_COARSEN:
324       ierr = DMCoarsen(dm, MPI_COMM_NULL, dmAdapted);CHKERRQ(ierr);break;
325     default:
326       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP,"DMPlex does not support refinement flag %D\n", minFlag);break;
327     }
328   } else {
329     ierr = DMPlexSetRefinementUniform(dm, PETSC_FALSE);CHKERRQ(ierr);
330     ierr = DMPlexRefine_Internal(dm, adaptLabel, dmAdapted);CHKERRQ(ierr);
331   }
332   PetscFunctionReturn(0);
333 }
334 
335 /*
336   DMAdaptMetric_Plex - Generates a new mesh conforming to a metric field.
337 
338   Input Parameters:
339 + dm - The DM object
340 . vertexMetric - The metric to which the mesh is adapted, defined vertex-wise in a LOCAL vector
341 - bdLabel - Label for boundary tags which are preserved in dmNew, or NULL. Should not be named "_boundary_".
342 
343   Output Parameter:
344 . dmNew  - the new DM
345 
346   Level: advanced
347 
348 .seealso: DMCoarsen(), DMRefine()
349 */
350 PetscErrorCode DMAdaptMetric_Plex(DM dm, Vec vertexMetric, DMLabel bdLabel, DM *dmNew)
351 {
352 #ifdef PETSC_HAVE_PRAGMATIC
353   MPI_Comm           comm;
354   const char        *bdName = "_boundary_";
355 #if 0
356   DM                 odm = dm;
357 #endif
358   DM                 udm, cdm;
359   DMLabel            bdLabelFull;
360   const char        *bdLabelName;
361   IS                 bdIS, globalVertexNum;
362   PetscSection       coordSection;
363   Vec                coordinates;
364   const PetscScalar *coords, *met;
365   const PetscInt    *bdFacesFull, *gV;
366   PetscInt          *bdFaces, *bdFaceIds, *l2gv;
367   PetscReal         *x, *y, *z, *metric;
368   PetscInt          *cells;
369   PetscInt           dim, cStart, cEnd, numCells, c, coff, vStart, vEnd, numVertices, numLocVertices, v;
370   PetscInt           off, maxConeSize, numBdFaces, f, bdSize;
371   PetscBool          flg;
372   DMLabel            bdLabelNew;
373   double            *coordsNew;
374   PetscInt          *bdTags;
375   PetscReal         *xNew[3] = {NULL, NULL, NULL};
376   PetscInt          *cellsNew;
377   PetscInt           d, numCellsNew, numVerticesNew;
378   PetscInt           numCornersNew, fStart, fEnd;
379   PetscMPIInt        numProcs;
380   PetscErrorCode     ierr;
381 
382   PetscFunctionBegin;
383   /* Check for FEM adjacency flags */
384   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
385   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
386   if (bdLabel) {
387     ierr = DMLabelGetName(bdLabel, &bdLabelName);CHKERRQ(ierr);
388     ierr = PetscStrcmp(bdLabelName, bdName, &flg);CHKERRQ(ierr);
389     if (flg) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName);
390   }
391   /* Add overlap for Pragmatic */
392 #if 0
393   /* Check for overlap by looking for cell in the SF */
394   if (!overlapped) {
395     ierr = DMPlexDistributeOverlap(odm, 1, NULL, &dm);CHKERRQ(ierr);
396     if (!dm) {dm = odm; ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);}
397   }
398 #endif
399   /* Get mesh information */
400   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
401   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
402   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
403   ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr);
404   ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr);
405   numCells    = cEnd - cStart;
406   numVertices = vEnd - vStart;
407   ierr = PetscCalloc5(numVertices, &x, numVertices, &y, numVertices, &z, numVertices*PetscSqr(dim), &metric, numCells*maxConeSize, &cells);CHKERRQ(ierr);
408   for (c = 0, coff = 0; c < numCells; ++c) {
409     const PetscInt *cone;
410     PetscInt        coneSize, cl;
411 
412     ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr);
413     ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr);
414     for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart;
415   }
416   ierr = PetscCalloc1(numVertices, &l2gv);CHKERRQ(ierr);
417   ierr = DMPlexGetVertexNumbering(udm, &globalVertexNum);CHKERRQ(ierr);
418   ierr = ISGetIndices(globalVertexNum, &gV);CHKERRQ(ierr);
419   for (v = 0, numLocVertices = 0; v < numVertices; ++v) {
420     if (gV[v] >= 0) ++numLocVertices;
421     l2gv[v] = gV[v] < 0 ? -(gV[v]+1) : gV[v];
422   }
423   ierr = ISRestoreIndices(globalVertexNum, &gV);CHKERRQ(ierr);
424   ierr = DMDestroy(&udm);CHKERRQ(ierr);
425   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
426   ierr = DMGetDefaultSection(cdm, &coordSection);CHKERRQ(ierr);
427   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
428   ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
429   for (v = vStart; v < vEnd; ++v) {
430     ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
431     x[v-vStart] = PetscRealPart(coords[off+0]);
432     if (dim > 1) y[v-vStart] = PetscRealPart(coords[off+1]);
433     if (dim > 2) z[v-vStart] = PetscRealPart(coords[off+2]);
434   }
435   ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
436   /* Get boundary mesh */
437   ierr = DMLabelCreate(bdName, &bdLabelFull);CHKERRQ(ierr);
438   ierr = DMPlexMarkBoundaryFaces(dm, 1, bdLabelFull);CHKERRQ(ierr);
439   ierr = DMLabelGetStratumIS(bdLabelFull, 1, &bdIS);CHKERRQ(ierr);
440   ierr = DMLabelGetStratumSize(bdLabelFull, 1, &numBdFaces);CHKERRQ(ierr);
441   ierr = ISGetIndices(bdIS, &bdFacesFull);CHKERRQ(ierr);
442   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
443     PetscInt *closure = NULL;
444     PetscInt  closureSize, cl;
445 
446     ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
447     for (cl = 0; cl < closureSize*2; cl += 2) {
448       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
449     }
450     ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
451   }
452   ierr = PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds);CHKERRQ(ierr);
453   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
454     PetscInt *closure = NULL;
455     PetscInt  closureSize, cl;
456 
457     ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
458     for (cl = 0; cl < closureSize*2; cl += 2) {
459       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart;
460     }
461     ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
462     if (bdLabel) {ierr = DMLabelGetValue(bdLabel, bdFacesFull[f], &bdFaceIds[f]);CHKERRQ(ierr);}
463     else         {bdFaceIds[f] = 1;}
464   }
465   ierr = ISDestroy(&bdIS);CHKERRQ(ierr);
466   ierr = DMLabelDestroy(&bdLabelFull);CHKERRQ(ierr);
467   /* Get metric */
468   ierr = VecViewFromOptions(vertexMetric, NULL, "-adapt_metric_view");CHKERRQ(ierr);
469   ierr = VecGetArrayRead(vertexMetric, &met);CHKERRQ(ierr);
470   for (v = 0; v < (vEnd-vStart)*PetscSqr(dim); ++v) metric[v] = PetscRealPart(met[v]);
471   ierr = VecRestoreArrayRead(vertexMetric, &met);CHKERRQ(ierr);
472 #if 0
473   /* Destroy overlap mesh */
474   ierr = DMDestroy(&dm);CHKERRQ(ierr);
475 #endif
476   /* Create new mesh */
477   switch (dim) {
478   case 2:
479     pragmatic_2d_mpi_init(&numVertices, &numCells, cells, x, y, l2gv, numLocVertices, comm);break;
480   case 3:
481     pragmatic_3d_mpi_init(&numVertices, &numCells, cells, x, y, z, l2gv, numLocVertices, comm);break;
482   default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim);
483   }
484   pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds);
485   pragmatic_set_metric(metric);
486   pragmatic_adapt(((DM_Plex *) dm->data)->remeshBd ? 1 : 0);
487   ierr = PetscFree(l2gv);CHKERRQ(ierr);
488   /* Read out mesh */
489   pragmatic_get_info_mpi(&numVerticesNew, &numCellsNew);
490   ierr = PetscMalloc1(numVerticesNew*dim, &coordsNew);CHKERRQ(ierr);
491   switch (dim) {
492   case 2:
493     numCornersNew = 3;
494     ierr = PetscMalloc2(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1]);CHKERRQ(ierr);
495     pragmatic_get_coords_2d_mpi(xNew[0], xNew[1]);
496     break;
497   case 3:
498     numCornersNew = 4;
499     ierr = PetscMalloc3(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1], numVerticesNew, &xNew[2]);CHKERRQ(ierr);
500     pragmatic_get_coords_3d_mpi(xNew[0], xNew[1], xNew[2]);
501     break;
502   default:
503     SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim);
504   }
505   for (v = 0; v < numVerticesNew; ++v) {for (d = 0; d < dim; ++d) coordsNew[v*dim+d] = (double) xNew[d][v];}
506   ierr = PetscMalloc1(numCellsNew*(dim+1), &cellsNew);CHKERRQ(ierr);
507   pragmatic_get_elements(cellsNew);
508   ierr = DMPlexCreateFromCellListParallel(comm, dim, numCellsNew, numVerticesNew, numCornersNew, PETSC_TRUE, cellsNew, dim, coordsNew, NULL, dmNew);CHKERRQ(ierr);
509   /* Read out boundary label */
510   pragmatic_get_boundaryTags(&bdTags);
511   ierr = DMCreateLabel(*dmNew, bdLabel ? bdLabelName : bdName);CHKERRQ(ierr);
512   ierr = DMGetLabel(*dmNew, bdLabel ? bdLabelName : bdName, &bdLabelNew);CHKERRQ(ierr);
513   ierr = DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd);CHKERRQ(ierr);
514   ierr = DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd);CHKERRQ(ierr);
515   ierr = DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd);CHKERRQ(ierr);
516   for (c = cStart; c < cEnd; ++c) {
517     /* Only for simplicial meshes */
518     coff = (c-cStart)*(dim+1);
519     /* d is the local cell number of the vertex opposite to the face we are marking */
520     for (d = 0; d < dim+1; ++d) {
521       if (bdTags[coff+d]) {
522         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 */
523         const PetscInt *cone;
524 
525         /* Mark face opposite to this vertex: This pattern is specified in DMPlexGetRawFaces_Internal() */
526         ierr = DMPlexGetCone(*dmNew, c, &cone);CHKERRQ(ierr);
527         ierr = DMLabelSetValue(bdLabelNew, cone[perm[dim][d]], bdTags[coff+d]);CHKERRQ(ierr);
528       }
529     }
530   }
531   /* Cleanup */
532   switch (dim) {
533   case 2: ierr = PetscFree2(xNew[0], xNew[1]);CHKERRQ(ierr);break;
534   case 3: ierr = PetscFree3(xNew[0], xNew[1], xNew[2]);CHKERRQ(ierr);break;
535   }
536   ierr = PetscFree(cellsNew);CHKERRQ(ierr);
537   ierr = PetscFree5(x, y, z, metric, cells);CHKERRQ(ierr);
538   ierr = PetscFree2(bdFaces, bdFaceIds);CHKERRQ(ierr);
539   ierr = PetscFree(coordsNew);CHKERRQ(ierr);
540   pragmatic_finalize();
541   PetscFunctionReturn(0);
542 #else
543   PetscFunctionBegin;
544   SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Remeshing needs external package support.\nPlease reconfigure with --download-pragmatic.");
545   PetscFunctionReturn(0);
546 #endif
547 }
548