xref: /petsc/src/dm/impls/plex/plexadapt.c (revision 3cbbc5ff037ae03fd9010ee4a8e4d7c21ef7cb71)
10bbe5d31Sbarral #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
20bbe5d31Sbarral #ifdef PETSC_HAVE_PRAGMATIC
30bbe5d31Sbarral #include <pragmatic/cpragmatic.h>
40bbe5d31Sbarral #endif
50bbe5d31Sbarral 
60d1cd5e0SMatthew G. Knepley static PetscErrorCode DMPlexLabelToVolumeConstraint(DM dm, DMLabel adaptLabel, PetscInt cStart, PetscInt cEnd, PetscReal refRatio, PetscReal maxVolumes[])
70d1cd5e0SMatthew G. Knepley {
80d1cd5e0SMatthew G. Knepley   PetscInt       dim, c;
90d1cd5e0SMatthew G. Knepley   PetscErrorCode ierr;
100d1cd5e0SMatthew G. Knepley 
110d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
120d1cd5e0SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
130d1cd5e0SMatthew G. Knepley   refRatio = refRatio == PETSC_DEFAULT ? (PetscReal) ((PetscInt) 1 << dim) : refRatio;
140d1cd5e0SMatthew G. Knepley   for (c = cStart; c < cEnd; c++) {
150d1cd5e0SMatthew G. Knepley     PetscReal vol;
160d1cd5e0SMatthew G. Knepley     PetscInt  closureSize = 0, cl;
170d1cd5e0SMatthew G. Knepley     PetscInt *closure     = NULL;
180d1cd5e0SMatthew G. Knepley     PetscBool anyRefine   = PETSC_FALSE;
190d1cd5e0SMatthew G. Knepley     PetscBool anyCoarsen  = PETSC_FALSE;
200d1cd5e0SMatthew G. Knepley     PetscBool anyKeep     = PETSC_FALSE;
210d1cd5e0SMatthew G. Knepley 
220d1cd5e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, NULL, NULL);CHKERRQ(ierr);
230d1cd5e0SMatthew G. Knepley     maxVolumes[c - cStart] = vol;
240d1cd5e0SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
250d1cd5e0SMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
260d1cd5e0SMatthew G. Knepley       const PetscInt point = closure[cl];
270d1cd5e0SMatthew G. Knepley       PetscInt       refFlag;
280d1cd5e0SMatthew G. Knepley 
290d1cd5e0SMatthew G. Knepley       ierr = DMLabelGetValue(adaptLabel, point, &refFlag);CHKERRQ(ierr);
300d1cd5e0SMatthew G. Knepley       switch (refFlag) {
310d1cd5e0SMatthew G. Knepley       case DM_ADAPT_REFINE:
320d1cd5e0SMatthew G. Knepley         anyRefine  = PETSC_TRUE;break;
330d1cd5e0SMatthew G. Knepley       case DM_ADAPT_COARSEN:
340d1cd5e0SMatthew G. Knepley         anyCoarsen = PETSC_TRUE;break;
350d1cd5e0SMatthew G. Knepley       case DM_ADAPT_KEEP:
360d1cd5e0SMatthew G. Knepley         anyKeep    = PETSC_TRUE;break;
370d1cd5e0SMatthew G. Knepley       case DM_ADAPT_DETERMINE:
380d1cd5e0SMatthew G. Knepley         break;
390d1cd5e0SMatthew G. Knepley       default:
400d1cd5e0SMatthew G. Knepley         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP, "DMPlex does not support refinement flag %D\n", refFlag);break;
410d1cd5e0SMatthew G. Knepley       }
420d1cd5e0SMatthew G. Knepley       if (anyRefine) break;
430d1cd5e0SMatthew G. Knepley     }
440d1cd5e0SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
450d1cd5e0SMatthew G. Knepley     if (anyRefine) {
460d1cd5e0SMatthew G. Knepley       maxVolumes[c - cStart] = vol / refRatio;
470d1cd5e0SMatthew G. Knepley     } else if (anyKeep) {
480d1cd5e0SMatthew G. Knepley       maxVolumes[c - cStart] = vol;
490d1cd5e0SMatthew G. Knepley     } else if (anyCoarsen) {
500d1cd5e0SMatthew G. Knepley       maxVolumes[c - cStart] = vol * refRatio;
510d1cd5e0SMatthew G. Knepley     }
520d1cd5e0SMatthew G. Knepley   }
530d1cd5e0SMatthew G. Knepley   PetscFunctionReturn(0);
540d1cd5e0SMatthew G. Knepley }
550d1cd5e0SMatthew G. Knepley 
560d1cd5e0SMatthew G. Knepley static PetscErrorCode DMPlexLabelToMetricConstraint(DM dm, DMLabel adaptLabel, PetscInt cStart, PetscInt cEnd, PetscInt vStart, PetscInt vEnd, PetscReal refRatio, Vec *metricVec)
570d1cd5e0SMatthew G. Knepley {
580d1cd5e0SMatthew G. Knepley   DM              udm, coordDM;
590d1cd5e0SMatthew G. Knepley   PetscSection    coordSection;
600d1cd5e0SMatthew G. Knepley   Vec             coordinates, mb, mx;
610d1cd5e0SMatthew G. Knepley   Mat             A;
620d1cd5e0SMatthew G. Knepley   PetscScalar    *metric, *eqns;
630d1cd5e0SMatthew G. Knepley   const PetscReal coarseRatio = refRatio == PETSC_DEFAULT ? PetscSqr(0.5) : 1/refRatio;
640d1cd5e0SMatthew G. Knepley   PetscInt        dim, Nv, Neq, c, v;
650d1cd5e0SMatthew G. Knepley   PetscErrorCode  ierr;
660d1cd5e0SMatthew G. Knepley 
670d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
680d1cd5e0SMatthew G. Knepley   ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr);
690d1cd5e0SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
700d1cd5e0SMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &coordDM);CHKERRQ(ierr);
710d1cd5e0SMatthew G. Knepley   ierr = DMGetDefaultSection(coordDM, &coordSection);CHKERRQ(ierr);
720d1cd5e0SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
730d1cd5e0SMatthew G. Knepley   Nv   = vEnd - vStart;
740d1cd5e0SMatthew G. Knepley   ierr = VecCreateSeq(PETSC_COMM_SELF, Nv*PetscSqr(dim), metricVec);CHKERRQ(ierr);
750d1cd5e0SMatthew G. Knepley   ierr = VecGetArray(*metricVec, &metric);CHKERRQ(ierr);
760d1cd5e0SMatthew G. Knepley   Neq  = (dim*(dim+1))/2;
770d1cd5e0SMatthew G. Knepley   ierr = PetscMalloc1(PetscSqr(Neq), &eqns);CHKERRQ(ierr);
780d1cd5e0SMatthew G. Knepley   ierr = MatCreateSeqDense(PETSC_COMM_SELF, Neq, Neq, eqns, &A);CHKERRQ(ierr);
790d1cd5e0SMatthew G. Knepley   ierr = MatCreateVecs(A, &mx, &mb);CHKERRQ(ierr);
800d1cd5e0SMatthew G. Knepley   ierr = VecSet(mb, 1.0);CHKERRQ(ierr);
810d1cd5e0SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
820d1cd5e0SMatthew G. Knepley     const PetscScalar *sol;
830d1cd5e0SMatthew G. Knepley     PetscScalar       *cellCoords = NULL;
840d1cd5e0SMatthew G. Knepley     PetscReal          e[3], vol;
850d1cd5e0SMatthew G. Knepley     const PetscInt    *cone;
860d1cd5e0SMatthew G. Knepley     PetscInt           coneSize, cl, i, j, d, r;
870d1cd5e0SMatthew G. Knepley 
880d1cd5e0SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, c, NULL, &cellCoords);CHKERRQ(ierr);
890d1cd5e0SMatthew G. Knepley     /* Only works for simplices */
900d1cd5e0SMatthew G. Knepley     for (i = 0, r = 0; i < dim+1; ++i) {
910d1cd5e0SMatthew G. Knepley       for (j = 0; j < i; ++j, ++r) {
920d1cd5e0SMatthew G. Knepley         for (d = 0; d < dim; ++d) e[d] = PetscRealPart(cellCoords[i*dim+d] - cellCoords[j*dim+d]);
930d1cd5e0SMatthew G. Knepley         /* FORTRAN ORDERING */
940d1cd5e0SMatthew G. Knepley         switch (dim) {
950d1cd5e0SMatthew G. Knepley         case 2:
960d1cd5e0SMatthew G. Knepley           eqns[0*Neq+r] = PetscSqr(e[0]);
970d1cd5e0SMatthew G. Knepley           eqns[1*Neq+r] = 2.0*e[0]*e[1];
980d1cd5e0SMatthew G. Knepley           eqns[2*Neq+r] = PetscSqr(e[1]);
990d1cd5e0SMatthew G. Knepley           break;
1000d1cd5e0SMatthew G. Knepley         case 3:
1010d1cd5e0SMatthew G. Knepley           eqns[0*Neq+r] = PetscSqr(e[0]);
1020d1cd5e0SMatthew G. Knepley           eqns[1*Neq+r] = 2.0*e[0]*e[1];
1030d1cd5e0SMatthew G. Knepley           eqns[2*Neq+r] = 2.0*e[0]*e[2];
1040d1cd5e0SMatthew G. Knepley           eqns[3*Neq+r] = PetscSqr(e[1]);
1050d1cd5e0SMatthew G. Knepley           eqns[4*Neq+r] = 2.0*e[1]*e[2];
1060d1cd5e0SMatthew G. Knepley           eqns[5*Neq+r] = PetscSqr(e[2]);
1070d1cd5e0SMatthew G. Knepley           break;
1080d1cd5e0SMatthew G. Knepley         }
1090d1cd5e0SMatthew G. Knepley       }
1100d1cd5e0SMatthew G. Knepley     }
1110d1cd5e0SMatthew G. Knepley     ierr = MatSetUnfactored(A);CHKERRQ(ierr);
1120d1cd5e0SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, c, NULL, &cellCoords);CHKERRQ(ierr);
1130d1cd5e0SMatthew G. Knepley     ierr = MatLUFactor(A, NULL, NULL, NULL);CHKERRQ(ierr);
1140d1cd5e0SMatthew G. Knepley     ierr = MatSolve(A, mb, mx);CHKERRQ(ierr);
1150d1cd5e0SMatthew G. Knepley     ierr = VecGetArrayRead(mx, &sol);CHKERRQ(ierr);
1160d1cd5e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, NULL, NULL);CHKERRQ(ierr);
1170d1cd5e0SMatthew G. Knepley     ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr);
1180d1cd5e0SMatthew G. Knepley     ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr);
1190d1cd5e0SMatthew G. Knepley     for (cl = 0; cl < coneSize; ++cl) {
1200d1cd5e0SMatthew G. Knepley       const PetscInt v = cone[cl] - vStart;
1210d1cd5e0SMatthew G. Knepley 
1220d1cd5e0SMatthew G. Knepley       if (dim == 2) {
1230d1cd5e0SMatthew G. Knepley         metric[v*4+0] += vol*coarseRatio*sol[0];
1240d1cd5e0SMatthew G. Knepley         metric[v*4+1] += vol*coarseRatio*sol[1];
1250d1cd5e0SMatthew G. Knepley         metric[v*4+2] += vol*coarseRatio*sol[1];
1260d1cd5e0SMatthew G. Knepley         metric[v*4+3] += vol*coarseRatio*sol[2];
1270d1cd5e0SMatthew G. Knepley       } else {
1280d1cd5e0SMatthew G. Knepley         metric[v*9+0] += vol*coarseRatio*sol[0];
1290d1cd5e0SMatthew G. Knepley         metric[v*9+1] += vol*coarseRatio*sol[1];
1300d1cd5e0SMatthew G. Knepley         metric[v*9+3] += vol*coarseRatio*sol[1];
1310d1cd5e0SMatthew G. Knepley         metric[v*9+2] += vol*coarseRatio*sol[2];
1320d1cd5e0SMatthew G. Knepley         metric[v*9+6] += vol*coarseRatio*sol[2];
1330d1cd5e0SMatthew G. Knepley         metric[v*9+4] += vol*coarseRatio*sol[3];
1340d1cd5e0SMatthew G. Knepley         metric[v*9+5] += vol*coarseRatio*sol[4];
1350d1cd5e0SMatthew G. Knepley         metric[v*9+7] += vol*coarseRatio*sol[4];
1360d1cd5e0SMatthew G. Knepley         metric[v*9+8] += vol*coarseRatio*sol[5];
1370d1cd5e0SMatthew G. Knepley       }
1380d1cd5e0SMatthew G. Knepley     }
1390d1cd5e0SMatthew G. Knepley     ierr = VecRestoreArrayRead(mx, &sol);CHKERRQ(ierr);
1400d1cd5e0SMatthew G. Knepley   }
1410d1cd5e0SMatthew G. Knepley   for (v = 0; v < Nv; ++v) {
1420d1cd5e0SMatthew G. Knepley     const PetscInt *support;
1430d1cd5e0SMatthew G. Knepley     PetscInt        supportSize, s;
1440d1cd5e0SMatthew G. Knepley     PetscReal       vol, totVol = 0.0;
1450d1cd5e0SMatthew G. Knepley 
1460d1cd5e0SMatthew G. Knepley     ierr = DMPlexGetSupport(udm, v+vStart, &support);CHKERRQ(ierr);
1470d1cd5e0SMatthew G. Knepley     ierr = DMPlexGetSupportSize(udm, v+vStart, &supportSize);CHKERRQ(ierr);
1480d1cd5e0SMatthew G. Knepley     for (s = 0; s < supportSize; ++s) {ierr = DMPlexComputeCellGeometryFVM(dm, support[s], &vol, NULL, NULL);CHKERRQ(ierr); totVol += vol;}
1490d1cd5e0SMatthew G. Knepley     for (s = 0; s < PetscSqr(dim); ++s) metric[v*PetscSqr(dim)+s] /= totVol;
1500d1cd5e0SMatthew G. Knepley   }
1510d1cd5e0SMatthew G. Knepley   ierr = PetscFree(eqns);CHKERRQ(ierr);
1520d1cd5e0SMatthew G. Knepley   ierr = VecRestoreArray(*metricVec, &metric);CHKERRQ(ierr);
1530d1cd5e0SMatthew G. Knepley   ierr = VecDestroy(&mx);CHKERRQ(ierr);
1540d1cd5e0SMatthew G. Knepley   ierr = VecDestroy(&mb);CHKERRQ(ierr);
1550d1cd5e0SMatthew G. Knepley   ierr = MatDestroy(&A);CHKERRQ(ierr);
1560d1cd5e0SMatthew G. Knepley   ierr = DMDestroy(&udm);CHKERRQ(ierr);
1570d1cd5e0SMatthew G. Knepley   PetscFunctionReturn(0);
1580d1cd5e0SMatthew G. Knepley }
1590d1cd5e0SMatthew G. Knepley 
1600d1cd5e0SMatthew G. Knepley PetscErrorCode DMPlexRefine_Internal(DM dm, DMLabel adaptLabel, DM *dmRefined)
1610d1cd5e0SMatthew G. Knepley {
1620d1cd5e0SMatthew G. Knepley   PetscErrorCode (*refinementFunc)(const PetscReal [], PetscReal *);
1630d1cd5e0SMatthew G. Knepley   PetscReal        refinementLimit;
1640d1cd5e0SMatthew G. Knepley   PetscInt         dim, cStart, cEnd;
1650d1cd5e0SMatthew G. Knepley   char             genname[1024], *name = NULL;
1660d1cd5e0SMatthew G. Knepley   PetscBool        isTriangle = PETSC_FALSE, isTetgen = PETSC_FALSE, isCTetgen = PETSC_FALSE, flg, localized;
1670d1cd5e0SMatthew G. Knepley   PetscErrorCode   ierr;
1680d1cd5e0SMatthew G. Knepley 
1690d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
1700d1cd5e0SMatthew G. Knepley   ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr);
1710d1cd5e0SMatthew G. Knepley   ierr = DMPlexGetRefinementLimit(dm, &refinementLimit);CHKERRQ(ierr);
1720d1cd5e0SMatthew G. Knepley   ierr = DMPlexGetRefinementFunction(dm, &refinementFunc);CHKERRQ(ierr);
1730d1cd5e0SMatthew G. Knepley   if (refinementLimit == 0.0 && !refinementFunc && !adaptLabel) PetscFunctionReturn(0);
1740d1cd5e0SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1750d1cd5e0SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1760d1cd5e0SMatthew G. Knepley   ierr = PetscOptionsGetString(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_generator", genname, 1024, &flg);CHKERRQ(ierr);
1770d1cd5e0SMatthew G. Knepley   if (flg) name = genname;
1780d1cd5e0SMatthew G. Knepley   if (name) {
1790d1cd5e0SMatthew G. Knepley     ierr = PetscStrcmp(name, "triangle", &isTriangle);CHKERRQ(ierr);
1800d1cd5e0SMatthew G. Knepley     ierr = PetscStrcmp(name, "tetgen",   &isTetgen);CHKERRQ(ierr);
1810d1cd5e0SMatthew G. Knepley     ierr = PetscStrcmp(name, "ctetgen",  &isCTetgen);CHKERRQ(ierr);
1820d1cd5e0SMatthew G. Knepley   }
1830d1cd5e0SMatthew G. Knepley   switch (dim) {
1840d1cd5e0SMatthew G. Knepley   case 2:
1850d1cd5e0SMatthew G. Knepley     if (!name || isTriangle) {
1860d1cd5e0SMatthew G. Knepley #if defined(PETSC_HAVE_TRIANGLE)
1870d1cd5e0SMatthew G. Knepley       PetscReal *maxVolumes;
1880d1cd5e0SMatthew G. Knepley       PetscInt  c;
1890d1cd5e0SMatthew G. Knepley 
1900d1cd5e0SMatthew G. Knepley       ierr = PetscMalloc1(cEnd - cStart, &maxVolumes);CHKERRQ(ierr);
1910d1cd5e0SMatthew G. Knepley       if (adaptLabel) {
1920d1cd5e0SMatthew G. Knepley         ierr = DMPlexLabelToVolumeConstraint(dm, adaptLabel, cStart, cEnd, PETSC_DEFAULT, maxVolumes);CHKERRQ(ierr);
1930d1cd5e0SMatthew G. Knepley       } else if (refinementFunc) {
1940d1cd5e0SMatthew G. Knepley         for (c = cStart; c < cEnd; ++c) {
1950d1cd5e0SMatthew G. Knepley           PetscReal vol, centroid[3];
1960d1cd5e0SMatthew G. Knepley           PetscReal maxVol;
1970d1cd5e0SMatthew G. Knepley 
1980d1cd5e0SMatthew G. Knepley           ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL);CHKERRQ(ierr);
1990d1cd5e0SMatthew G. Knepley           ierr = (*refinementFunc)(centroid, &maxVol);CHKERRQ(ierr);
2000d1cd5e0SMatthew G. Knepley           maxVolumes[c - cStart] = (double) maxVol;
2010d1cd5e0SMatthew G. Knepley         }
2020d1cd5e0SMatthew G. Knepley       } else {
2030d1cd5e0SMatthew G. Knepley         for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit;
2040d1cd5e0SMatthew G. Knepley       }
2050d1cd5e0SMatthew G. Knepley #if !defined(PETSC_USE_REAL_DOUBLE)
2060d1cd5e0SMatthew G. Knepley       {
2070d1cd5e0SMatthew G. Knepley         double *mvols;
2080d1cd5e0SMatthew G. Knepley         ierr = PetscMalloc1(cEnd - cStart,&mvols);CHKERRQ(ierr);
2090d1cd5e0SMatthew G. Knepley         for (c = 0; c < cEnd-cStart; ++c) mvols[c] = (double)maxVolumes[c];
2100d1cd5e0SMatthew G. Knepley         ierr = DMPlexRefine_Triangle(dm, mvols, dmRefined);CHKERRQ(ierr);
2110d1cd5e0SMatthew G. Knepley         ierr = PetscFree(mvols);CHKERRQ(ierr);
2120d1cd5e0SMatthew G. Knepley       }
2130d1cd5e0SMatthew G. Knepley #else
2140d1cd5e0SMatthew G. Knepley       ierr = DMPlexRefine_Triangle(dm, maxVolumes, dmRefined);CHKERRQ(ierr);
2150d1cd5e0SMatthew G. Knepley #endif
2160d1cd5e0SMatthew G. Knepley       ierr = PetscFree(maxVolumes);CHKERRQ(ierr);
2170d1cd5e0SMatthew G. Knepley #else
2180d1cd5e0SMatthew G. Knepley       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Mesh refinement needs external package support.\nPlease reconfigure with --download-triangle.");
2190d1cd5e0SMatthew G. Knepley #endif
2200d1cd5e0SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unknown 2D mesh generation package %s", name);
2210d1cd5e0SMatthew G. Knepley     break;
2220d1cd5e0SMatthew G. Knepley   case 3:
2230d1cd5e0SMatthew G. Knepley     if (!name || isCTetgen || isTetgen) {
2240d1cd5e0SMatthew G. Knepley       PetscReal *maxVolumes;
2250d1cd5e0SMatthew G. Knepley       PetscInt   c;
2260d1cd5e0SMatthew G. Knepley 
2270d1cd5e0SMatthew G. Knepley       ierr = PetscMalloc1(cEnd - cStart, &maxVolumes);CHKERRQ(ierr);
2280d1cd5e0SMatthew G. Knepley       if (adaptLabel) {
2290d1cd5e0SMatthew G. Knepley         ierr = DMPlexLabelToVolumeConstraint(dm, adaptLabel, cStart, cEnd, PETSC_DEFAULT, maxVolumes);CHKERRQ(ierr);
2300d1cd5e0SMatthew G. Knepley       } else if (refinementFunc) {
2310d1cd5e0SMatthew G. Knepley         for (c = cStart; c < cEnd; ++c) {
2320d1cd5e0SMatthew G. Knepley           PetscReal vol, centroid[3];
2330d1cd5e0SMatthew G. Knepley 
2340d1cd5e0SMatthew G. Knepley           ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL);CHKERRQ(ierr);
2350d1cd5e0SMatthew G. Knepley           ierr = (*refinementFunc)(centroid, &maxVolumes[c-cStart]);CHKERRQ(ierr);
2360d1cd5e0SMatthew G. Knepley         }
2370d1cd5e0SMatthew G. Knepley       } else {
2380d1cd5e0SMatthew G. Knepley         for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit;
2390d1cd5e0SMatthew G. Knepley       }
2400d1cd5e0SMatthew G. Knepley       if (!name) {
2410d1cd5e0SMatthew G. Knepley #if defined(PETSC_HAVE_CTETGEN)
2420d1cd5e0SMatthew G. Knepley         ierr = DMPlexRefine_CTetgen(dm, maxVolumes, dmRefined);CHKERRQ(ierr);
2430d1cd5e0SMatthew G. Knepley #elif defined(PETSC_HAVE_TETGEN)
2440d1cd5e0SMatthew G. Knepley         ierr = DMPlexRefine_Tetgen(dm, maxVolumes, dmRefined);CHKERRQ(ierr);
2450d1cd5e0SMatthew G. Knepley #else
2460d1cd5e0SMatthew G. Knepley         SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "External package CTetgen or Tetgen needed.\nPlease reconfigure with '--download-ctetgen' or '--with-clanguage=cxx --download-tetgen'.");
2470d1cd5e0SMatthew G. Knepley #endif
2480d1cd5e0SMatthew G. Knepley       } else if (isCTetgen) {
2490d1cd5e0SMatthew G. Knepley #if defined(PETSC_HAVE_CTETGEN)
2500d1cd5e0SMatthew G. Knepley         ierr = DMPlexRefine_CTetgen(dm, maxVolumes, dmRefined);CHKERRQ(ierr);
2510d1cd5e0SMatthew G. Knepley #else
2520d1cd5e0SMatthew G. Knepley         SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "CTetgen needs external package support.\nPlease reconfigure with --download-ctetgen.");
2530d1cd5e0SMatthew G. Knepley #endif
2540d1cd5e0SMatthew G. Knepley       } else {
2550d1cd5e0SMatthew G. Knepley #if defined(PETSC_HAVE_TETGEN)
2560d1cd5e0SMatthew G. Knepley         ierr = DMPlexRefine_Tetgen(dm, maxVolumes, dmRefined);CHKERRQ(ierr);
2570d1cd5e0SMatthew G. Knepley #else
2580d1cd5e0SMatthew G. Knepley         SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Tetgen needs external package support.\nPlease reconfigure with --download-tetgen.");
2590d1cd5e0SMatthew G. Knepley #endif
2600d1cd5e0SMatthew G. Knepley       }
2610d1cd5e0SMatthew G. Knepley       ierr = PetscFree(maxVolumes);CHKERRQ(ierr);
2620d1cd5e0SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unknown 3D mesh generation package %s", name);
2630d1cd5e0SMatthew G. Knepley     break;
2640d1cd5e0SMatthew G. Knepley   default:
2650d1cd5e0SMatthew G. Knepley     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Mesh refinement in dimension %d is not supported.", dim);
2660d1cd5e0SMatthew G. Knepley   }
2670d1cd5e0SMatthew G. Knepley   ierr = DMCopyBoundary(dm, *dmRefined);CHKERRQ(ierr);
2680d1cd5e0SMatthew G. Knepley   if (localized) {ierr = DMLocalizeCoordinates(*dmRefined);CHKERRQ(ierr);}
2690d1cd5e0SMatthew G. Knepley   PetscFunctionReturn(0);
2700d1cd5e0SMatthew G. Knepley }
2710d1cd5e0SMatthew G. Knepley 
2720d1cd5e0SMatthew G. Knepley PetscErrorCode DMPlexCoarsen_Internal(DM dm, DMLabel adaptLabel, DM *dmCoarsened)
2730d1cd5e0SMatthew G. Knepley {
2740d1cd5e0SMatthew G. Knepley   Vec            metricVec;
2750d1cd5e0SMatthew G. Knepley   PetscInt       cStart, cEnd, vStart, vEnd;
2760d1cd5e0SMatthew G. Knepley   DMLabel        bdLabel = NULL;
2770d1cd5e0SMatthew G. Knepley   char           bdLabelName[PETSC_MAX_PATH_LEN];
2780d1cd5e0SMatthew G. Knepley   PetscBool      localized, flg;
2790d1cd5e0SMatthew G. Knepley   PetscErrorCode ierr;
2800d1cd5e0SMatthew G. Knepley 
2810d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
2820d1cd5e0SMatthew G. Knepley   ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr);
2830d1cd5e0SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2840d1cd5e0SMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
2850d1cd5e0SMatthew G. Knepley   ierr = DMPlexLabelToMetricConstraint(dm, adaptLabel, cStart, cEnd, vStart, vEnd, PETSC_DEFAULT, &metricVec);CHKERRQ(ierr);
2860d1cd5e0SMatthew G. Knepley   ierr = PetscOptionsGetString(NULL, dm->hdr.prefix, "-dm_plex_coarsen_bd_label", bdLabelName, PETSC_MAX_PATH_LEN-1, &flg);CHKERRQ(ierr);
2870d1cd5e0SMatthew G. Knepley   if (flg) {ierr = DMGetLabel(dm, bdLabelName, &bdLabel);CHKERRQ(ierr);}
2880d1cd5e0SMatthew G. Knepley   ierr = DMAdaptMetric_Plex(dm, metricVec, bdLabel, dmCoarsened);CHKERRQ(ierr);
2890d1cd5e0SMatthew G. Knepley   ierr = VecDestroy(&metricVec);CHKERRQ(ierr);
2900d1cd5e0SMatthew G. Knepley   if (localized) {ierr = DMLocalizeCoordinates(*dmCoarsened);CHKERRQ(ierr);}
2910d1cd5e0SMatthew G. Knepley   PetscFunctionReturn(0);
2920d1cd5e0SMatthew G. Knepley }
2930d1cd5e0SMatthew G. Knepley 
2940d1cd5e0SMatthew G. Knepley PetscErrorCode DMAdaptLabel_Plex(DM dm, DMLabel adaptLabel, DM *dmAdapted)
2950d1cd5e0SMatthew G. Knepley {
2960d1cd5e0SMatthew G. Knepley   IS              flagIS;
2970d1cd5e0SMatthew G. Knepley   const PetscInt *flags;
2980d1cd5e0SMatthew G. Knepley   PetscInt        defFlag, minFlag, maxFlag, numFlags, f;
2990d1cd5e0SMatthew G. Knepley   PetscErrorCode  ierr;
3000d1cd5e0SMatthew G. Knepley 
3010d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
3020d1cd5e0SMatthew G. Knepley   ierr = DMLabelGetDefaultValue(adaptLabel, &defFlag);CHKERRQ(ierr);
3030d1cd5e0SMatthew G. Knepley   minFlag = defFlag;
3040d1cd5e0SMatthew G. Knepley   maxFlag = defFlag;
3050d1cd5e0SMatthew G. Knepley   ierr = DMLabelGetValueIS(adaptLabel, &flagIS);CHKERRQ(ierr);
3060d1cd5e0SMatthew G. Knepley   ierr = ISGetLocalSize(flagIS, &numFlags);CHKERRQ(ierr);
3070d1cd5e0SMatthew G. Knepley   ierr = ISGetIndices(flagIS, &flags);CHKERRQ(ierr);
3080d1cd5e0SMatthew G. Knepley   for (f = 0; f < numFlags; ++f) {
3090d1cd5e0SMatthew G. Knepley     const PetscInt flag = flags[f];
3100d1cd5e0SMatthew G. Knepley 
3110d1cd5e0SMatthew G. Knepley     minFlag = PetscMin(minFlag, flag);
3120d1cd5e0SMatthew G. Knepley     maxFlag = PetscMax(maxFlag, flag);
3130d1cd5e0SMatthew G. Knepley   }
3140d1cd5e0SMatthew G. Knepley   ierr = ISRestoreIndices(flagIS, &flags);CHKERRQ(ierr);
3150d1cd5e0SMatthew G. Knepley   ierr = ISDestroy(&flagIS);CHKERRQ(ierr);
3160d1cd5e0SMatthew G. Knepley   {
3170d1cd5e0SMatthew G. Knepley     PetscInt minMaxFlag[2], minMaxFlagGlobal[2];
3180d1cd5e0SMatthew G. Knepley 
3190d1cd5e0SMatthew G. Knepley     minMaxFlag[0] =  minFlag;
3200d1cd5e0SMatthew G. Knepley     minMaxFlag[1] = -maxFlag;
3210d1cd5e0SMatthew G. Knepley     ierr = MPI_Allreduce(minMaxFlag, minMaxFlagGlobal, 2, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
3220d1cd5e0SMatthew G. Knepley     minFlag =  minMaxFlagGlobal[0];
3230d1cd5e0SMatthew G. Knepley     maxFlag = -minMaxFlagGlobal[1];
3240d1cd5e0SMatthew G. Knepley   }
3250d1cd5e0SMatthew G. Knepley   if (minFlag == maxFlag) {
3260d1cd5e0SMatthew G. Knepley     switch (minFlag) {
3270d1cd5e0SMatthew G. Knepley     case DM_ADAPT_DETERMINE:
3280d1cd5e0SMatthew G. Knepley       *dmAdapted = NULL;break;
3290d1cd5e0SMatthew G. Knepley     case DM_ADAPT_REFINE:
3300d1cd5e0SMatthew G. Knepley       ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE);CHKERRQ(ierr);
3310d1cd5e0SMatthew G. Knepley       ierr = DMRefine(dm, MPI_COMM_NULL, dmAdapted);CHKERRQ(ierr);break;
3320d1cd5e0SMatthew G. Knepley     case DM_ADAPT_COARSEN:
3330d1cd5e0SMatthew G. Knepley       ierr = DMCoarsen(dm, MPI_COMM_NULL, dmAdapted);CHKERRQ(ierr);break;
3340d1cd5e0SMatthew G. Knepley     default:
3350d1cd5e0SMatthew G. Knepley       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP,"DMPlex does not support refinement flag %D\n", minFlag);break;
3360d1cd5e0SMatthew G. Knepley     }
3370d1cd5e0SMatthew G. Knepley   } else {
3380d1cd5e0SMatthew G. Knepley     ierr = DMPlexSetRefinementUniform(dm, PETSC_FALSE);CHKERRQ(ierr);
3390d1cd5e0SMatthew G. Knepley     ierr = DMPlexRefine_Internal(dm, adaptLabel, dmAdapted);CHKERRQ(ierr);
3400d1cd5e0SMatthew G. Knepley   }
3410d1cd5e0SMatthew G. Knepley   PetscFunctionReturn(0);
3420d1cd5e0SMatthew G. Knepley }
3430d1cd5e0SMatthew G. Knepley 
344f7bcbcbbSMatthew G. Knepley /*
3450d1cd5e0SMatthew G. Knepley   DMAdaptMetric_Plex - Generates a new mesh conforming to a metric field.
3460bbe5d31Sbarral 
347f7bcbcbbSMatthew G. Knepley   Input Parameters:
348f7bcbcbbSMatthew G. Knepley + dm - The DM object
34951ff0a1cSMatthew G. Knepley . vertexMetric - The metric to which the mesh is adapted, defined vertex-wise in a LOCAL vector
3500d1cd5e0SMatthew G. Knepley - bdLabel - Label for boundary tags which are preserved in dmNew, or NULL. Should not be named "_boundary_".
351f7bcbcbbSMatthew G. Knepley 
352f7bcbcbbSMatthew G. Knepley   Output Parameter:
353f7bcbcbbSMatthew G. Knepley . dmNew  - the new DM
354f7bcbcbbSMatthew G. Knepley 
355f7bcbcbbSMatthew G. Knepley   Level: advanced
356f7bcbcbbSMatthew G. Knepley 
357f7bcbcbbSMatthew G. Knepley .seealso: DMCoarsen(), DMRefine()
358f7bcbcbbSMatthew G. Knepley */
3590d1cd5e0SMatthew G. Knepley PetscErrorCode DMAdaptMetric_Plex(DM dm, Vec vertexMetric, DMLabel bdLabel, DM *dmNew)
360f7bcbcbbSMatthew G. Knepley {
3616f25b0d8SLisandro Dalcin #ifdef PETSC_HAVE_PRAGMATIC
362bc1d7e78SMatthew G. Knepley   MPI_Comm           comm;
363f7bcbcbbSMatthew G. Knepley   const char        *bdName = "_boundary_";
3642b5f8d01SMatthew G. Knepley #if 0
3652b5f8d01SMatthew G. Knepley   DM                 odm = dm;
3662b5f8d01SMatthew G. Knepley #endif
367f7bcbcbbSMatthew G. Knepley   DM                 udm, cdm;
3680d1cd5e0SMatthew G. Knepley   DMLabel            bdLabelFull;
3690d1cd5e0SMatthew G. Knepley   const char        *bdLabelName;
3701838b4a6SMatthew G. Knepley   IS                 bdIS, globalVertexNum;
371f7bcbcbbSMatthew G. Knepley   PetscSection       coordSection;
372f7bcbcbbSMatthew G. Knepley   Vec                coordinates;
373f7bcbcbbSMatthew G. Knepley   const PetscScalar *coords, *met;
3741838b4a6SMatthew G. Knepley   const PetscInt    *bdFacesFull, *gV;
3751838b4a6SMatthew G. Knepley   PetscInt          *bdFaces, *bdFaceIds, *l2gv;
376f7bcbcbbSMatthew G. Knepley   PetscReal         *x, *y, *z, *metric;
3772ee120b0SMatthew G. Knepley   PetscInt          *cells;
3781838b4a6SMatthew G. Knepley   PetscInt           dim, cStart, cEnd, numCells, c, coff, vStart, vEnd, numVertices, numLocVertices, v;
3792ee120b0SMatthew G. Knepley   PetscInt           off, maxConeSize, numBdFaces, f, bdSize;
380f7bcbcbbSMatthew G. Knepley   PetscBool          flg;
3812ee120b0SMatthew G. Knepley   DMLabel            bdLabelNew;
382acfe6eadSToby Isaac   double            *coordsNew;
3832ee120b0SMatthew G. Knepley   PetscInt          *bdTags;
3842ee120b0SMatthew G. Knepley   PetscReal         *xNew[3] = {NULL, NULL, NULL};
3852ee120b0SMatthew G. Knepley   PetscInt          *cellsNew;
3862ee120b0SMatthew G. Knepley   PetscInt           d, numCellsNew, numVerticesNew;
3872ee120b0SMatthew G. Knepley   PetscInt           numCornersNew, fStart, fEnd;
388bc1d7e78SMatthew G. Knepley   PetscMPIInt        numProcs;
389f7bcbcbbSMatthew G. Knepley   PetscErrorCode     ierr;
390f7bcbcbbSMatthew G. Knepley 
391f7bcbcbbSMatthew G. Knepley   PetscFunctionBegin;
392bbc23918SMatthew G. Knepley   /* Check for FEM adjacency flags */
393bc1d7e78SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
394bc1d7e78SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
3950d1cd5e0SMatthew G. Knepley   if (bdLabel) {
3960d1cd5e0SMatthew G. Knepley     ierr = DMLabelGetName(bdLabel, &bdLabelName);CHKERRQ(ierr);
397f7bcbcbbSMatthew G. Knepley     ierr = PetscStrcmp(bdLabelName, bdName, &flg);CHKERRQ(ierr);
398bc1d7e78SMatthew G. Knepley     if (flg) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName);
399f7bcbcbbSMatthew G. Knepley   }
400a2e1bd45SMatthew G. Knepley   /* Add overlap for Pragmatic */
401bbc23918SMatthew G. Knepley #if 0
402bbc23918SMatthew G. Knepley   /* Check for overlap by looking for cell in the SF */
403bbc23918SMatthew G. Knepley   if (!overlapped) {
404a2e1bd45SMatthew G. Knepley     ierr = DMPlexDistributeOverlap(odm, 1, NULL, &dm);CHKERRQ(ierr);
405a2e1bd45SMatthew G. Knepley     if (!dm) {dm = odm; ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);}
406bbc23918SMatthew G. Knepley   }
407bbc23918SMatthew G. Knepley #endif
408f7bcbcbbSMatthew G. Knepley   /* Get mesh information */
409f7bcbcbbSMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
410f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
411f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
412f7bcbcbbSMatthew G. Knepley   ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr);
413f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr);
414f7bcbcbbSMatthew G. Knepley   numCells    = cEnd - cStart;
415f7bcbcbbSMatthew G. Knepley   numVertices = vEnd - vStart;
416f7bcbcbbSMatthew G. Knepley   ierr = PetscCalloc5(numVertices, &x, numVertices, &y, numVertices, &z, numVertices*PetscSqr(dim), &metric, numCells*maxConeSize, &cells);CHKERRQ(ierr);
417f7bcbcbbSMatthew G. Knepley   for (c = 0, coff = 0; c < numCells; ++c) {
418f7bcbcbbSMatthew G. Knepley     const PetscInt *cone;
419f7bcbcbbSMatthew G. Knepley     PetscInt        coneSize, cl;
420f7bcbcbbSMatthew G. Knepley 
421f7bcbcbbSMatthew G. Knepley     ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr);
422f7bcbcbbSMatthew G. Knepley     ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr);
423f7bcbcbbSMatthew G. Knepley     for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart;
424f7bcbcbbSMatthew G. Knepley   }
4251838b4a6SMatthew G. Knepley   ierr = PetscCalloc1(numVertices, &l2gv);CHKERRQ(ierr);
4261838b4a6SMatthew G. Knepley   ierr = DMPlexGetVertexNumbering(udm, &globalVertexNum);CHKERRQ(ierr);
4271838b4a6SMatthew G. Knepley   ierr = ISGetIndices(globalVertexNum, &gV);CHKERRQ(ierr);
4281838b4a6SMatthew G. Knepley   for (v = 0, numLocVertices = 0; v < numVertices; ++v) {
4291838b4a6SMatthew G. Knepley     if (gV[v] >= 0) ++numLocVertices;
4301838b4a6SMatthew G. Knepley     l2gv[v] = gV[v] < 0 ? -(gV[v]+1) : gV[v];
4311838b4a6SMatthew G. Knepley   }
4321838b4a6SMatthew G. Knepley   ierr = ISRestoreIndices(globalVertexNum, &gV);CHKERRQ(ierr);
433f7bcbcbbSMatthew G. Knepley   ierr = DMDestroy(&udm);CHKERRQ(ierr);
434f7bcbcbbSMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
435f7bcbcbbSMatthew G. Knepley   ierr = DMGetDefaultSection(cdm, &coordSection);CHKERRQ(ierr);
436f7bcbcbbSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
437f7bcbcbbSMatthew G. Knepley   ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
438f7bcbcbbSMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
439f7bcbcbbSMatthew G. Knepley     ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
4405f2e139eSMatthew G. Knepley     x[v-vStart] = PetscRealPart(coords[off+0]);
4415f2e139eSMatthew G. Knepley     if (dim > 1) y[v-vStart] = PetscRealPart(coords[off+1]);
4425f2e139eSMatthew G. Knepley     if (dim > 2) z[v-vStart] = PetscRealPart(coords[off+2]);
443f7bcbcbbSMatthew G. Knepley   }
444f7bcbcbbSMatthew G. Knepley   ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
445f7bcbcbbSMatthew G. Knepley   /* Get boundary mesh */
446f7bcbcbbSMatthew G. Knepley   ierr = DMLabelCreate(bdName, &bdLabelFull);CHKERRQ(ierr);
447f7bcbcbbSMatthew G. Knepley   ierr = DMPlexMarkBoundaryFaces(dm, bdLabelFull);CHKERRQ(ierr);
448f7bcbcbbSMatthew G. Knepley   ierr = DMLabelGetStratumIS(bdLabelFull, 1, &bdIS);CHKERRQ(ierr);
449f7bcbcbbSMatthew G. Knepley   ierr = DMLabelGetStratumSize(bdLabelFull, 1, &numBdFaces);CHKERRQ(ierr);
450f7bcbcbbSMatthew G. Knepley   ierr = ISGetIndices(bdIS, &bdFacesFull);CHKERRQ(ierr);
451f7bcbcbbSMatthew G. Knepley   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
452f7bcbcbbSMatthew G. Knepley     PetscInt *closure = NULL;
453f7bcbcbbSMatthew G. Knepley     PetscInt  closureSize, cl;
454f7bcbcbbSMatthew G. Knepley 
455f7bcbcbbSMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
456f7bcbcbbSMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
457f7bcbcbbSMatthew G. Knepley       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
458f7bcbcbbSMatthew G. Knepley     }
459f7bcbcbbSMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
460f7bcbcbbSMatthew G. Knepley   }
461f7bcbcbbSMatthew G. Knepley   ierr = PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds);CHKERRQ(ierr);
462f7bcbcbbSMatthew G. Knepley   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
463f7bcbcbbSMatthew G. Knepley     PetscInt *closure = NULL;
464f7bcbcbbSMatthew G. Knepley     PetscInt  closureSize, cl;
465f7bcbcbbSMatthew G. Knepley 
466f7bcbcbbSMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
467f7bcbcbbSMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
468f7bcbcbbSMatthew G. Knepley       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart;
469f7bcbcbbSMatthew G. Knepley     }
470f7bcbcbbSMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
471f7bcbcbbSMatthew G. Knepley     if (bdLabel) {ierr = DMLabelGetValue(bdLabel, bdFacesFull[f], &bdFaceIds[f]);CHKERRQ(ierr);}
472f7bcbcbbSMatthew G. Knepley     else         {bdFaceIds[f] = 1;}
473f7bcbcbbSMatthew G. Knepley   }
474f7bcbcbbSMatthew G. Knepley   ierr = ISDestroy(&bdIS);CHKERRQ(ierr);
475f7bcbcbbSMatthew G. Knepley   ierr = DMLabelDestroy(&bdLabelFull);CHKERRQ(ierr);
476f7bcbcbbSMatthew G. Knepley   /* Get metric */
47704996bafSMatthew G. Knepley   ierr = VecViewFromOptions(vertexMetric, NULL, "-adapt_metric_view");CHKERRQ(ierr);
478f7bcbcbbSMatthew G. Knepley   ierr = VecGetArrayRead(vertexMetric, &met);CHKERRQ(ierr);
4795f2e139eSMatthew G. Knepley   for (v = 0; v < (vEnd-vStart)*PetscSqr(dim); ++v) metric[v] = PetscRealPart(met[v]);
480f7bcbcbbSMatthew G. Knepley   ierr = VecRestoreArrayRead(vertexMetric, &met);CHKERRQ(ierr);
481bbc23918SMatthew G. Knepley #if 0
482a2e1bd45SMatthew G. Knepley   /* Destroy overlap mesh */
483a2e1bd45SMatthew G. Knepley   ierr = DMDestroy(&dm);CHKERRQ(ierr);
484bbc23918SMatthew G. Knepley #endif
485f7bcbcbbSMatthew G. Knepley   /* Create new mesh */
486f7bcbcbbSMatthew G. Knepley   switch (dim) {
487f7bcbcbbSMatthew G. Knepley   case 2:
488cc29c497SMatthew G. Knepley     pragmatic_2d_mpi_init(&numVertices, &numCells, cells, x, y, l2gv, numLocVertices, comm);break;
489f7bcbcbbSMatthew G. Knepley   case 3:
490cc29c497SMatthew G. Knepley     pragmatic_3d_mpi_init(&numVertices, &numCells, cells, x, y, z, l2gv, numLocVertices, comm);break;
4911838b4a6SMatthew G. Knepley   default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim);
492f7bcbcbbSMatthew G. Knepley   }
493f7bcbcbbSMatthew G. Knepley   pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds);
494f7bcbcbbSMatthew G. Knepley   pragmatic_set_metric(metric);
4950d1cd5e0SMatthew G. Knepley   pragmatic_adapt(((DM_Plex *) dm->data)->remeshBd ? 1 : 0);
4961838b4a6SMatthew G. Knepley   ierr = PetscFree(l2gv);CHKERRQ(ierr);
497f7bcbcbbSMatthew G. Knepley   /* Read out mesh */
4988713909fSMatthew G. Knepley   pragmatic_get_info_mpi(&numVerticesNew, &numCellsNew);
499f7bcbcbbSMatthew G. Knepley   ierr = PetscMalloc1(numVerticesNew*dim, &coordsNew);CHKERRQ(ierr);
500f7bcbcbbSMatthew G. Knepley   switch (dim) {
501f7bcbcbbSMatthew G. Knepley   case 2:
502f7bcbcbbSMatthew G. Knepley     numCornersNew = 3;
503f7bcbcbbSMatthew G. Knepley     ierr = PetscMalloc2(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1]);CHKERRQ(ierr);
5048713909fSMatthew G. Knepley     pragmatic_get_coords_2d_mpi(xNew[0], xNew[1]);
505f7bcbcbbSMatthew G. Knepley     break;
506f7bcbcbbSMatthew G. Knepley   case 3:
507f7bcbcbbSMatthew G. Knepley     numCornersNew = 4;
508f7bcbcbbSMatthew G. Knepley     ierr = PetscMalloc3(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1], numVerticesNew, &xNew[2]);CHKERRQ(ierr);
5098713909fSMatthew G. Knepley     pragmatic_get_coords_3d_mpi(xNew[0], xNew[1], xNew[2]);
510f7bcbcbbSMatthew G. Knepley     break;
511f7bcbcbbSMatthew G. Knepley   default:
512bc1d7e78SMatthew G. Knepley     SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim);
513f7bcbcbbSMatthew G. Knepley   }
514acfe6eadSToby Isaac   for (v = 0; v < numVerticesNew; ++v) {for (d = 0; d < dim; ++d) coordsNew[v*dim+d] = (double) xNew[d][v];}
515f7bcbcbbSMatthew G. Knepley   ierr = PetscMalloc1(numCellsNew*(dim+1), &cellsNew);CHKERRQ(ierr);
516f7bcbcbbSMatthew G. Knepley   pragmatic_get_elements(cellsNew);
51751ff0a1cSMatthew G. Knepley   ierr = DMPlexCreateFromCellListParallel(comm, dim, numCellsNew, numVerticesNew, numCornersNew, PETSC_TRUE, cellsNew, dim, coordsNew, NULL, dmNew);CHKERRQ(ierr);
518f7bcbcbbSMatthew G. Knepley   /* Read out boundary label */
519f7bcbcbbSMatthew G. Knepley   pragmatic_get_boundaryTags(&bdTags);
520f7bcbcbbSMatthew G. Knepley   ierr = DMCreateLabel(*dmNew, bdLabel ? bdLabelName : bdName);CHKERRQ(ierr);
521f7bcbcbbSMatthew G. Knepley   ierr = DMGetLabel(*dmNew, bdLabel ? bdLabelName : bdName, &bdLabelNew);CHKERRQ(ierr);
522f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd);CHKERRQ(ierr);
523f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd);CHKERRQ(ierr);
524f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd);CHKERRQ(ierr);
525f7bcbcbbSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
526f7bcbcbbSMatthew G. Knepley     /* Only for simplicial meshes */
527f7bcbcbbSMatthew G. Knepley     coff = (c-cStart)*(dim+1);
5288713909fSMatthew G. Knepley     /* d is the local cell number of the vertex opposite to the face we are marking */
529f7bcbcbbSMatthew G. Knepley     for (d = 0; d < dim+1; ++d) {
530f7bcbcbbSMatthew G. Knepley       if (bdTags[coff+d]) {
5318713909fSMatthew G. Knepley         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 */
5328713909fSMatthew G. Knepley         const PetscInt *cone;
533f7bcbcbbSMatthew G. Knepley 
5348713909fSMatthew G. Knepley         /* Mark face opposite to this vertex: This pattern is specified in DMPlexGetRawFaces_Internal() */
5358713909fSMatthew G. Knepley         ierr = DMPlexGetCone(*dmNew, c, &cone);CHKERRQ(ierr);
5368713909fSMatthew G. Knepley         ierr = DMLabelSetValue(bdLabelNew, cone[perm[dim][d]], bdTags[coff+d]);CHKERRQ(ierr);
537f7bcbcbbSMatthew G. Knepley       }
538f7bcbcbbSMatthew G. Knepley     }
539f7bcbcbbSMatthew G. Knepley   }
540f7bcbcbbSMatthew G. Knepley   /* Cleanup */
541f7bcbcbbSMatthew G. Knepley   switch (dim) {
542f7bcbcbbSMatthew G. Knepley   case 2: ierr = PetscFree2(xNew[0], xNew[1]);CHKERRQ(ierr);break;
543f7bcbcbbSMatthew G. Knepley   case 3: ierr = PetscFree3(xNew[0], xNew[1], xNew[2]);CHKERRQ(ierr);break;
544f7bcbcbbSMatthew G. Knepley   }
545f7bcbcbbSMatthew G. Knepley   ierr = PetscFree(cellsNew);CHKERRQ(ierr);
546f7bcbcbbSMatthew G. Knepley   ierr = PetscFree5(x, y, z, metric, cells);CHKERRQ(ierr);
547f7bcbcbbSMatthew G. Knepley   ierr = PetscFree2(bdFaces, bdFaceIds);CHKERRQ(ierr);
548f7bcbcbbSMatthew G. Knepley   ierr = PetscFree(coordsNew);CHKERRQ(ierr);
549f7bcbcbbSMatthew G. Knepley   pragmatic_finalize();
550f7bcbcbbSMatthew G. Knepley   PetscFunctionReturn(0);
5516f25b0d8SLisandro Dalcin #else
5526f25b0d8SLisandro Dalcin   PetscFunctionBegin;
5536f25b0d8SLisandro Dalcin   SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Remeshing needs external package support.\nPlease reconfigure with --download-pragmatic.");
554*3cbbc5ffSBarry Smith   PetscFunctionReturn(0);
5556f25b0d8SLisandro Dalcin #endif
556f7bcbcbbSMatthew G. Knepley }
557