Actual source code: plexadapt.c

petsc-3.8.3 2017-12-09
Report Typos and Errors
  1:  #include <petsc/private/dmpleximpl.h>
  2: #ifdef PETSC_HAVE_PRAGMATIC
  3: #include <pragmatic/cpragmatic.h>
  4: #endif

  6: static PetscErrorCode DMPlexLabelToVolumeConstraint(DM dm, DMLabel adaptLabel, PetscInt cStart, PetscInt cEnd, PetscReal refRatio, PetscReal maxVolumes[])
  7: {
  8:   PetscInt       dim, c;

 12:   DMGetDimension(dm, &dim);
 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;

 22:     DMPlexComputeCellGeometryFVM(dm, c, &vol, NULL, NULL);
 23:     maxVolumes[c - cStart] = vol;
 24:     DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
 25:     for (cl = 0; cl < closureSize*2; cl += 2) {
 26:       const PetscInt point = closure[cl];
 27:       PetscInt       refFlag;

 29:       DMLabelGetValue(adaptLabel, point, &refFlag);
 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:     DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
 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:   return(0);
 54: }

 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;

 68:   DMPlexUninterpolate(dm, &udm);
 69:   DMGetDimension(dm, &dim);
 70:   DMGetCoordinateDM(dm, &coordDM);
 71:   DMGetDefaultSection(coordDM, &coordSection);
 72:   DMGetCoordinatesLocal(dm, &coordinates);
 73:   Nv   = vEnd - vStart;
 74:   VecCreateSeq(PETSC_COMM_SELF, Nv*PetscSqr(dim), metricVec);
 75:   VecGetArray(*metricVec, &metric);
 76:   Neq  = (dim*(dim+1))/2;
 77:   PetscMalloc1(PetscSqr(Neq), &eqns);
 78:   MatCreateSeqDense(PETSC_COMM_SELF, Neq, Neq, eqns, &A);
 79:   MatCreateVecs(A, &mx, &mb);
 80:   VecSet(mb, 1.0);
 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;

 88:     DMPlexVecGetClosure(dm, coordSection, coordinates, c, NULL, &cellCoords);
 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:     MatSetUnfactored(A);
112:     DMPlexVecRestoreClosure(dm, coordSection, coordinates, c, NULL, &cellCoords);
113:     MatLUFactor(A, NULL, NULL, NULL);
114:     MatSolve(A, mb, mx);
115:     VecGetArrayRead(mx, &sol);
116:     DMPlexComputeCellGeometryFVM(dm, c, &vol, NULL, NULL);
117:     DMPlexGetCone(udm, c, &cone);
118:     DMPlexGetConeSize(udm, c, &coneSize);
119:     for (cl = 0; cl < coneSize; ++cl) {
120:       const PetscInt v = cone[cl] - vStart;

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:     VecRestoreArrayRead(mx, &sol);
140:   }
141:   for (v = 0; v < Nv; ++v) {
142:     const PetscInt *support;
143:     PetscInt        supportSize, s;
144:     PetscReal       vol, totVol = 0.0;

146:     DMPlexGetSupport(udm, v+vStart, &support);
147:     DMPlexGetSupportSize(udm, v+vStart, &supportSize);
148:     for (s = 0; s < supportSize; ++s) {DMPlexComputeCellGeometryFVM(dm, support[s], &vol, NULL, NULL); totVol += vol;}
149:     for (s = 0; s < PetscSqr(dim); ++s) metric[v*PetscSqr(dim)+s] /= totVol;
150:   }
151:   PetscFree(eqns);
152:   VecRestoreArray(*metricVec, &metric);
153:   VecDestroy(&mx);
154:   VecDestroy(&mb);
155:   MatDestroy(&A);
156:   DMDestroy(&udm);
157:   return(0);
158: }

160: PetscErrorCode DMPlexRefine_Internal(DM dm, DMLabel adaptLabel, DM *dmRefined)
161: {
162:   PetscErrorCode (*refinementFunc)(const PetscReal [], PetscReal *);
163:   PetscReal        refinementLimit;
164:   PetscInt         dim, cStart, cEnd;
165:   char             genname[1024], *name = NULL;
166:   PetscBool        isTriangle = PETSC_FALSE, isTetgen = PETSC_FALSE, isCTetgen = PETSC_FALSE, flg, localized;
167:   PetscErrorCode   ierr;

170:   DMGetCoordinatesLocalized(dm, &localized);
171:   DMPlexGetRefinementLimit(dm, &refinementLimit);
172:   DMPlexGetRefinementFunction(dm, &refinementFunc);
173:   if (refinementLimit == 0.0 && !refinementFunc && !adaptLabel) return(0);
174:   DMGetDimension(dm, &dim);
175:   DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
176:   PetscOptionsGetString(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_generator", genname, 1024, &flg);
177:   if (flg) name = genname;
178:   if (name) {
179:     PetscStrcmp(name, "triangle", &isTriangle);
180:     PetscStrcmp(name, "tetgen",   &isTetgen);
181:     PetscStrcmp(name, "ctetgen",  &isCTetgen);
182:   }
183:   switch (dim) {
184:   case 2:
185:     if (!name || isTriangle) {
186: #if defined(PETSC_HAVE_TRIANGLE)
187:       PetscReal *maxVolumes;
188:       PetscInt  c;

190:       PetscMalloc1(cEnd - cStart, &maxVolumes);
191:       if (adaptLabel) {
192:         DMPlexLabelToVolumeConstraint(dm, adaptLabel, cStart, cEnd, PETSC_DEFAULT, maxVolumes);
193:       } else if (refinementFunc) {
194:         for (c = cStart; c < cEnd; ++c) {
195:           PetscReal vol, centroid[3];
196:           PetscReal maxVol;

198:           DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL);
199:           (*refinementFunc)(centroid, &maxVol);
200:           maxVolumes[c - cStart] = (double) maxVol;
201:         }
202:       } else {
203:         for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit;
204:       }
205: #if !defined(PETSC_USE_REAL_DOUBLE)
206:       {
207:         double *mvols;
208:         PetscMalloc1(cEnd - cStart,&mvols);
209:         for (c = 0; c < cEnd-cStart; ++c) mvols[c] = (double)maxVolumes[c];
210:         DMPlexRefine_Triangle(dm, mvols, dmRefined);
211:         PetscFree(mvols);
212:       }
213: #else
214:       DMPlexRefine_Triangle(dm, maxVolumes, dmRefined);
215: #endif
216:       PetscFree(maxVolumes);
217: #else
218:       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Mesh refinement needs external package support.\nPlease reconfigure with --download-triangle.");
219: #endif
220:     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unknown 2D mesh generation package %s", name);
221:     break;
222:   case 3:
223:     if (!name || isCTetgen || isTetgen) {
224:       PetscReal *maxVolumes;
225:       PetscInt   c;

227:       PetscMalloc1(cEnd - cStart, &maxVolumes);
228:       if (adaptLabel) {
229:         DMPlexLabelToVolumeConstraint(dm, adaptLabel, cStart, cEnd, PETSC_DEFAULT, maxVolumes);
230:       } else if (refinementFunc) {
231:         for (c = cStart; c < cEnd; ++c) {
232:           PetscReal vol, centroid[3];

234:           DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL);
235:           (*refinementFunc)(centroid, &maxVolumes[c-cStart]);
236:         }
237:       } else {
238:         for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit;
239:       }
240:       if (!name) {
241: #if defined(PETSC_HAVE_CTETGEN)
242:         DMPlexRefine_CTetgen(dm, maxVolumes, dmRefined);
243: #elif defined(PETSC_HAVE_TETGEN)
244:         DMPlexRefine_Tetgen(dm, maxVolumes, dmRefined);
245: #else
246:         SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "External package CTetgen or Tetgen needed.\nPlease reconfigure with '--download-ctetgen' or '--with-clanguage=cxx --download-tetgen'.");
247: #endif
248:       } else if (isCTetgen) {
249: #if defined(PETSC_HAVE_CTETGEN)
250:         DMPlexRefine_CTetgen(dm, maxVolumes, dmRefined);
251: #else
252:         SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "CTetgen needs external package support.\nPlease reconfigure with --download-ctetgen.");
253: #endif
254:       } else {
255: #if defined(PETSC_HAVE_TETGEN)
256:         DMPlexRefine_Tetgen(dm, maxVolumes, dmRefined);
257: #else
258:         SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Tetgen needs external package support.\nPlease reconfigure with --download-tetgen.");
259: #endif
260:       }
261:       PetscFree(maxVolumes);
262:     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unknown 3D mesh generation package %s", name);
263:     break;
264:   default:
265:     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Mesh refinement in dimension %d is not supported.", dim);
266:   }
267:   DMCopyBoundary(dm, *dmRefined);
268:   if (localized) {DMLocalizeCoordinates(*dmRefined);}
269:   return(0);
270: }

272: PetscErrorCode DMPlexCoarsen_Internal(DM dm, DMLabel adaptLabel, DM *dmCoarsened)
273: {
274:   Vec            metricVec;
275:   PetscInt       cStart, cEnd, vStart, vEnd;
276:   DMLabel        bdLabel = NULL;
277:   char           bdLabelName[PETSC_MAX_PATH_LEN];
278:   PetscBool      localized, flg;

282:   DMGetCoordinatesLocalized(dm, &localized);
283:   DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
284:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
285:   DMPlexLabelToMetricConstraint(dm, adaptLabel, cStart, cEnd, vStart, vEnd, PETSC_DEFAULT, &metricVec);
286:   PetscOptionsGetString(NULL, dm->hdr.prefix, "-dm_plex_coarsen_bd_label", bdLabelName, PETSC_MAX_PATH_LEN-1, &flg);
287:   if (flg) {DMGetLabel(dm, bdLabelName, &bdLabel);}
288:   DMAdaptMetric_Plex(dm, metricVec, bdLabel, dmCoarsened);
289:   VecDestroy(&metricVec);
290:   if (localized) {DMLocalizeCoordinates(*dmCoarsened);}
291:   return(0);
292: }

294: PetscErrorCode DMAdaptLabel_Plex(DM dm, DMLabel adaptLabel, DM *dmAdapted)
295: {
296:   IS              flagIS;
297:   const PetscInt *flags;
298:   PetscInt        defFlag, minFlag, maxFlag, numFlags, f;
299:   PetscErrorCode  ierr;

302:   DMLabelGetDefaultValue(adaptLabel, &defFlag);
303:   minFlag = defFlag;
304:   maxFlag = defFlag;
305:   DMLabelGetValueIS(adaptLabel, &flagIS);
306:   ISGetLocalSize(flagIS, &numFlags);
307:   ISGetIndices(flagIS, &flags);
308:   for (f = 0; f < numFlags; ++f) {
309:     const PetscInt flag = flags[f];

311:     minFlag = PetscMin(minFlag, flag);
312:     maxFlag = PetscMax(maxFlag, flag);
313:   }
314:   ISRestoreIndices(flagIS, &flags);
315:   ISDestroy(&flagIS);
316:   {
317:     PetscInt minMaxFlag[2], minMaxFlagGlobal[2];

319:     minMaxFlag[0] =  minFlag;
320:     minMaxFlag[1] = -maxFlag;
321:     MPI_Allreduce(minMaxFlag, minMaxFlagGlobal, 2, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm));
322:     minFlag =  minMaxFlagGlobal[0];
323:     maxFlag = -minMaxFlagGlobal[1];
324:   }
325:   if (minFlag == maxFlag) {
326:     switch (minFlag) {
327:     case DM_ADAPT_DETERMINE:
328:       *dmAdapted = NULL;break;
329:     case DM_ADAPT_REFINE:
330:       DMPlexSetRefinementUniform(dm, PETSC_TRUE);
331:       DMRefine(dm, MPI_COMM_NULL, dmAdapted);break;
332:     case DM_ADAPT_COARSEN:
333:       DMCoarsen(dm, MPI_COMM_NULL, dmAdapted);break;
334:     default:
335:       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP,"DMPlex does not support refinement flag %D\n", minFlag);break;
336:     }
337:   } else {
338:     DMPlexSetRefinementUniform(dm, PETSC_FALSE);
339:     DMPlexRefine_Internal(dm, adaptLabel, dmAdapted);
340:   }
341:   return(0);
342: }

344: /*
345:   DMAdaptMetric_Plex - Generates a new mesh conforming to a metric field.

347:   Input Parameters:
348: + dm - The DM object
349: . vertexMetric - The metric to which the mesh is adapted, defined vertex-wise in a LOCAL vector
350: - bdLabel - Label for boundary tags which are preserved in dmNew, or NULL. Should not be named "_boundary_".

352:   Output Parameter:
353: . dmNew  - the new DM

355:   Level: advanced

357: .seealso: DMCoarsen(), DMRefine()
358: */
359: PetscErrorCode DMAdaptMetric_Plex(DM dm, Vec vertexMetric, DMLabel bdLabel, DM *dmNew)
360: {
361: #ifdef PETSC_HAVE_PRAGMATIC
362:   MPI_Comm           comm;
363:   const char        *bdName = "_boundary_";
364: #if 0
365:   DM                 odm = dm;
366: #endif
367:   DM                 udm, cdm;
368:   DMLabel            bdLabelFull;
369:   const char        *bdLabelName;
370:   IS                 bdIS, globalVertexNum;
371:   PetscSection       coordSection;
372:   Vec                coordinates;
373:   const PetscScalar *coords, *met;
374:   const PetscInt    *bdFacesFull, *gV;
375:   PetscInt          *bdFaces, *bdFaceIds, *l2gv;
376:   PetscReal         *x, *y, *z, *metric;
377:   PetscInt          *cells;
378:   PetscInt           dim, cStart, cEnd, numCells, c, coff, vStart, vEnd, numVertices, numLocVertices, v;
379:   PetscInt           off, maxConeSize, numBdFaces, f, bdSize;
380:   PetscBool          flg;
381:   DMLabel            bdLabelNew;
382:   double            *coordsNew;
383:   PetscInt          *bdTags;
384:   PetscReal         *xNew[3] = {NULL, NULL, NULL};
385:   PetscInt          *cellsNew;
386:   PetscInt           d, numCellsNew, numVerticesNew;
387:   PetscInt           numCornersNew, fStart, fEnd;
388:   PetscMPIInt        numProcs;
389:   PetscErrorCode     ierr;

392:   /* Check for FEM adjacency flags */
393:   PetscObjectGetComm((PetscObject) dm, &comm);
394:   MPI_Comm_size(comm, &numProcs);
395:   if (bdLabel) {
396:     DMLabelGetName(bdLabel, &bdLabelName);
397:     PetscStrcmp(bdLabelName, bdName, &flg);
398:     if (flg) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName);
399:   }
400:   /* Add overlap for Pragmatic */
401: #if 0
402:   /* Check for overlap by looking for cell in the SF */
403:   if (!overlapped) {
404:     DMPlexDistributeOverlap(odm, 1, NULL, &dm);
405:     if (!dm) {dm = odm; PetscObjectReference((PetscObject) dm);}
406:   }
407: #endif
408:   /* Get mesh information */
409:   DMGetDimension(dm, &dim);
410:   DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
411:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
412:   DMPlexUninterpolate(dm, &udm);
413:   DMPlexGetMaxSizes(udm, &maxConeSize, NULL);
414:   numCells    = cEnd - cStart;
415:   numVertices = vEnd - vStart;
416:   PetscCalloc5(numVertices, &x, numVertices, &y, numVertices, &z, numVertices*PetscSqr(dim), &metric, numCells*maxConeSize, &cells);
417:   for (c = 0, coff = 0; c < numCells; ++c) {
418:     const PetscInt *cone;
419:     PetscInt        coneSize, cl;

421:     DMPlexGetConeSize(udm, c, &coneSize);
422:     DMPlexGetCone(udm, c, &cone);
423:     for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart;
424:   }
425:   PetscCalloc1(numVertices, &l2gv);
426:   DMPlexGetVertexNumbering(udm, &globalVertexNum);
427:   ISGetIndices(globalVertexNum, &gV);
428:   for (v = 0, numLocVertices = 0; v < numVertices; ++v) {
429:     if (gV[v] >= 0) ++numLocVertices;
430:     l2gv[v] = gV[v] < 0 ? -(gV[v]+1) : gV[v];
431:   }
432:   ISRestoreIndices(globalVertexNum, &gV);
433:   DMDestroy(&udm);
434:   DMGetCoordinateDM(dm, &cdm);
435:   DMGetDefaultSection(cdm, &coordSection);
436:   DMGetCoordinatesLocal(dm, &coordinates);
437:   VecGetArrayRead(coordinates, &coords);
438:   for (v = vStart; v < vEnd; ++v) {
439:     PetscSectionGetOffset(coordSection, v, &off);
440:     x[v-vStart] = PetscRealPart(coords[off+0]);
441:     if (dim > 1) y[v-vStart] = PetscRealPart(coords[off+1]);
442:     if (dim > 2) z[v-vStart] = PetscRealPart(coords[off+2]);
443:   }
444:   VecRestoreArrayRead(coordinates, &coords);
445:   /* Get boundary mesh */
446:   DMLabelCreate(bdName, &bdLabelFull);
447:   DMPlexMarkBoundaryFaces(dm, bdLabelFull);
448:   DMLabelGetStratumIS(bdLabelFull, 1, &bdIS);
449:   DMLabelGetStratumSize(bdLabelFull, 1, &numBdFaces);
450:   ISGetIndices(bdIS, &bdFacesFull);
451:   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
452:     PetscInt *closure = NULL;
453:     PetscInt  closureSize, cl;

455:     DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);
456:     for (cl = 0; cl < closureSize*2; cl += 2) {
457:       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
458:     }
459:     DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);
460:   }
461:   PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds);
462:   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
463:     PetscInt *closure = NULL;
464:     PetscInt  closureSize, cl;

466:     DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);
467:     for (cl = 0; cl < closureSize*2; cl += 2) {
468:       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart;
469:     }
470:     DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);
471:     if (bdLabel) {DMLabelGetValue(bdLabel, bdFacesFull[f], &bdFaceIds[f]);}
472:     else         {bdFaceIds[f] = 1;}
473:   }
474:   ISDestroy(&bdIS);
475:   DMLabelDestroy(&bdLabelFull);
476:   /* Get metric */
477:   VecGetArrayRead(vertexMetric, &met);
478:   for (v = 0; v < (vEnd-vStart)*PetscSqr(dim); ++v) metric[v] = PetscRealPart(met[v]);
479:   VecRestoreArrayRead(vertexMetric, &met);
480: #if 0
481:   /* Destroy overlap mesh */
482:   DMDestroy(&dm);
483: #endif
484:   /* Create new mesh */
485:   switch (dim) {
486:   case 2:
487:     pragmatic_2d_mpi_init(&numVertices, &numCells, cells, x, y, l2gv, numLocVertices, comm);break;
488:   case 3:
489:     pragmatic_3d_mpi_init(&numVertices, &numCells, cells, x, y, z, l2gv, numLocVertices, comm);break;
490:   default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim);
491:   }
492:   pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds);
493:   pragmatic_set_metric(metric);
494:   pragmatic_adapt(((DM_Plex *) dm->data)->remeshBd ? 1 : 0);
495:   PetscFree(l2gv);
496:   /* Read out mesh */
497:   pragmatic_get_info_mpi(&numVerticesNew, &numCellsNew);
498:   PetscMalloc1(numVerticesNew*dim, &coordsNew);
499:   switch (dim) {
500:   case 2:
501:     numCornersNew = 3;
502:     PetscMalloc2(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1]);
503:     pragmatic_get_coords_2d_mpi(xNew[0], xNew[1]);
504:     break;
505:   case 3:
506:     numCornersNew = 4;
507:     PetscMalloc3(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1], numVerticesNew, &xNew[2]);
508:     pragmatic_get_coords_3d_mpi(xNew[0], xNew[1], xNew[2]);
509:     break;
510:   default:
511:     SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim);
512:   }
513:   for (v = 0; v < numVerticesNew; ++v) {for (d = 0; d < dim; ++d) coordsNew[v*dim+d] = (double) xNew[d][v];}
514:   PetscMalloc1(numCellsNew*(dim+1), &cellsNew);
515:   pragmatic_get_elements(cellsNew);
516:   DMPlexCreateFromCellListParallel(comm, dim, numCellsNew, numVerticesNew, numCornersNew, PETSC_TRUE, cellsNew, dim, coordsNew, NULL, dmNew);
517:   /* Read out boundary label */
518:   pragmatic_get_boundaryTags(&bdTags);
519:   DMCreateLabel(*dmNew, bdLabel ? bdLabelName : bdName);
520:   DMGetLabel(*dmNew, bdLabel ? bdLabelName : bdName, &bdLabelNew);
521:   DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd);
522:   DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd);
523:   DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd);
524:   for (c = cStart; c < cEnd; ++c) {
525:     /* Only for simplicial meshes */
526:     coff = (c-cStart)*(dim+1);
527:     /* d is the local cell number of the vertex opposite to the face we are marking */
528:     for (d = 0; d < dim+1; ++d) {
529:       if (bdTags[coff+d]) {
530:         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 */
531:         const PetscInt *cone;

533:         /* Mark face opposite to this vertex: This pattern is specified in DMPlexGetRawFaces_Internal() */
534:         DMPlexGetCone(*dmNew, c, &cone);
535:         DMLabelSetValue(bdLabelNew, cone[perm[dim][d]], bdTags[coff+d]);
536:       }
537:     }
538:   }
539:   /* Cleanup */
540:   switch (dim) {
541:   case 2: PetscFree2(xNew[0], xNew[1]);break;
542:   case 3: PetscFree3(xNew[0], xNew[1], xNew[2]);break;
543:   }
544:   PetscFree(cellsNew);
545:   PetscFree5(x, y, z, metric, cells);
546:   PetscFree2(bdFaces, bdFaceIds);
547:   PetscFree(coordsNew);
548:   pragmatic_finalize();
549:   return(0);
550: #else
552:   SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Remeshing needs external package support.\nPlease reconfigure with --download-pragmatic.");
553:   PetscFunctionReturn(PETSC_ERR_SUP);
554: #endif
555: }