xref: /petsc/src/dm/impls/plex/plexadapt.c (revision 5de52c6d2b8d6de382140bd9fae367e1f02f2f13)
10bbe5d31Sbarral #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
20bbe5d31Sbarral 
30d1cd5e0SMatthew G. Knepley static PetscErrorCode DMPlexLabelToVolumeConstraint(DM dm, DMLabel adaptLabel, PetscInt cStart, PetscInt cEnd, PetscReal refRatio, PetscReal maxVolumes[])
40d1cd5e0SMatthew G. Knepley {
50d1cd5e0SMatthew G. Knepley   PetscInt       dim, c;
60d1cd5e0SMatthew G. Knepley 
70d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
89566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
90d1cd5e0SMatthew G. Knepley   refRatio = refRatio == PETSC_DEFAULT ? (PetscReal) ((PetscInt) 1 << dim) : refRatio;
100d1cd5e0SMatthew G. Knepley   for (c = cStart; c < cEnd; c++) {
110d1cd5e0SMatthew G. Knepley     PetscReal vol;
120d1cd5e0SMatthew G. Knepley     PetscInt  closureSize = 0, cl;
130d1cd5e0SMatthew G. Knepley     PetscInt *closure     = NULL;
140d1cd5e0SMatthew G. Knepley     PetscBool anyRefine   = PETSC_FALSE;
150d1cd5e0SMatthew G. Knepley     PetscBool anyCoarsen  = PETSC_FALSE;
160d1cd5e0SMatthew G. Knepley     PetscBool anyKeep     = PETSC_FALSE;
170d1cd5e0SMatthew G. Knepley 
189566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeCellGeometryFVM(dm, c, &vol, NULL, NULL));
190d1cd5e0SMatthew G. Knepley     maxVolumes[c - cStart] = vol;
209566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
210d1cd5e0SMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
220d1cd5e0SMatthew G. Knepley       const PetscInt point = closure[cl];
230d1cd5e0SMatthew G. Knepley       PetscInt       refFlag;
240d1cd5e0SMatthew G. Knepley 
259566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValue(adaptLabel, point, &refFlag));
260d1cd5e0SMatthew G. Knepley       switch (refFlag) {
270d1cd5e0SMatthew G. Knepley       case DM_ADAPT_REFINE:
280d1cd5e0SMatthew G. Knepley         anyRefine  = PETSC_TRUE;break;
290d1cd5e0SMatthew G. Knepley       case DM_ADAPT_COARSEN:
300d1cd5e0SMatthew G. Knepley         anyCoarsen = PETSC_TRUE;break;
310d1cd5e0SMatthew G. Knepley       case DM_ADAPT_KEEP:
320d1cd5e0SMatthew G. Knepley         anyKeep    = PETSC_TRUE;break;
330d1cd5e0SMatthew G. Knepley       case DM_ADAPT_DETERMINE:
340d1cd5e0SMatthew G. Knepley         break;
350d1cd5e0SMatthew G. Knepley       default:
3663a3b9bcSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "DMPlex does not support refinement flag %" PetscInt_FMT, refFlag);
370d1cd5e0SMatthew G. Knepley       }
380d1cd5e0SMatthew G. Knepley       if (anyRefine) break;
390d1cd5e0SMatthew G. Knepley     }
409566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
410d1cd5e0SMatthew G. Knepley     if (anyRefine) {
420d1cd5e0SMatthew G. Knepley       maxVolumes[c - cStart] = vol / refRatio;
430d1cd5e0SMatthew G. Knepley     } else if (anyKeep) {
440d1cd5e0SMatthew G. Knepley       maxVolumes[c - cStart] = vol;
450d1cd5e0SMatthew G. Knepley     } else if (anyCoarsen) {
460d1cd5e0SMatthew G. Knepley       maxVolumes[c - cStart] = vol * refRatio;
470d1cd5e0SMatthew G. Knepley     }
480d1cd5e0SMatthew G. Knepley   }
490d1cd5e0SMatthew G. Knepley   PetscFunctionReturn(0);
500d1cd5e0SMatthew G. Knepley }
510d1cd5e0SMatthew G. Knepley 
520d1cd5e0SMatthew G. Knepley static PetscErrorCode DMPlexLabelToMetricConstraint(DM dm, DMLabel adaptLabel, PetscInt cStart, PetscInt cEnd, PetscInt vStart, PetscInt vEnd, PetscReal refRatio, Vec *metricVec)
530d1cd5e0SMatthew G. Knepley {
540d1cd5e0SMatthew G. Knepley   DM              udm, coordDM;
550d1cd5e0SMatthew G. Knepley   PetscSection    coordSection;
560d1cd5e0SMatthew G. Knepley   Vec             coordinates, mb, mx;
570d1cd5e0SMatthew G. Knepley   Mat             A;
580d1cd5e0SMatthew G. Knepley   PetscScalar    *metric, *eqns;
590d1cd5e0SMatthew G. Knepley   const PetscReal coarseRatio = refRatio == PETSC_DEFAULT ? PetscSqr(0.5) : 1/refRatio;
600d1cd5e0SMatthew G. Knepley   PetscInt        dim, Nv, Neq, c, v;
610d1cd5e0SMatthew G. Knepley 
620d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
639566063dSJacob Faibussowitsch   PetscCall(DMPlexUninterpolate(dm, &udm));
649566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
659566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &coordDM));
669566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(coordDM, &coordSection));
679566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
680d1cd5e0SMatthew G. Knepley   Nv   = vEnd - vStart;
699566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, Nv*PetscSqr(dim), metricVec));
709566063dSJacob Faibussowitsch   PetscCall(VecGetArray(*metricVec, &metric));
710d1cd5e0SMatthew G. Knepley   Neq  = (dim*(dim+1))/2;
729566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(PetscSqr(Neq), &eqns));
739566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, Neq, Neq, eqns, &A));
749566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &mx, &mb));
759566063dSJacob Faibussowitsch   PetscCall(VecSet(mb, 1.0));
760d1cd5e0SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
770d1cd5e0SMatthew G. Knepley     const PetscScalar *sol;
780d1cd5e0SMatthew G. Knepley     PetscScalar       *cellCoords = NULL;
790d1cd5e0SMatthew G. Knepley     PetscReal          e[3], vol;
800d1cd5e0SMatthew G. Knepley     const PetscInt    *cone;
810d1cd5e0SMatthew G. Knepley     PetscInt           coneSize, cl, i, j, d, r;
820d1cd5e0SMatthew G. Knepley 
839566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dm, coordSection, coordinates, c, NULL, &cellCoords));
840d1cd5e0SMatthew G. Knepley     /* Only works for simplices */
850d1cd5e0SMatthew G. Knepley     for (i = 0, r = 0; i < dim+1; ++i) {
860d1cd5e0SMatthew G. Knepley       for (j = 0; j < i; ++j, ++r) {
870d1cd5e0SMatthew G. Knepley         for (d = 0; d < dim; ++d) e[d] = PetscRealPart(cellCoords[i*dim+d] - cellCoords[j*dim+d]);
880d1cd5e0SMatthew G. Knepley         /* FORTRAN ORDERING */
890d1cd5e0SMatthew G. Knepley         switch (dim) {
900d1cd5e0SMatthew G. Knepley         case 2:
910d1cd5e0SMatthew G. Knepley           eqns[0*Neq+r] = PetscSqr(e[0]);
920d1cd5e0SMatthew G. Knepley           eqns[1*Neq+r] = 2.0*e[0]*e[1];
930d1cd5e0SMatthew G. Knepley           eqns[2*Neq+r] = PetscSqr(e[1]);
940d1cd5e0SMatthew G. Knepley           break;
950d1cd5e0SMatthew G. Knepley         case 3:
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] = 2.0*e[0]*e[2];
990d1cd5e0SMatthew G. Knepley           eqns[3*Neq+r] = PetscSqr(e[1]);
1000d1cd5e0SMatthew G. Knepley           eqns[4*Neq+r] = 2.0*e[1]*e[2];
1010d1cd5e0SMatthew G. Knepley           eqns[5*Neq+r] = PetscSqr(e[2]);
1020d1cd5e0SMatthew G. Knepley           break;
1030d1cd5e0SMatthew G. Knepley         }
1040d1cd5e0SMatthew G. Knepley       }
1050d1cd5e0SMatthew G. Knepley     }
1069566063dSJacob Faibussowitsch     PetscCall(MatSetUnfactored(A));
1079566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dm, coordSection, coordinates, c, NULL, &cellCoords));
1089566063dSJacob Faibussowitsch     PetscCall(MatLUFactor(A, NULL, NULL, NULL));
1099566063dSJacob Faibussowitsch     PetscCall(MatSolve(A, mb, mx));
1109566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(mx, &sol));
1119566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeCellGeometryFVM(dm, c, &vol, NULL, NULL));
1129566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCone(udm, c, &cone));
1139566063dSJacob Faibussowitsch     PetscCall(DMPlexGetConeSize(udm, c, &coneSize));
1140d1cd5e0SMatthew G. Knepley     for (cl = 0; cl < coneSize; ++cl) {
1150d1cd5e0SMatthew G. Knepley       const PetscInt v = cone[cl] - vStart;
1160d1cd5e0SMatthew G. Knepley 
1170d1cd5e0SMatthew G. Knepley       if (dim == 2) {
1180d1cd5e0SMatthew G. Knepley         metric[v*4+0] += vol*coarseRatio*sol[0];
1190d1cd5e0SMatthew G. Knepley         metric[v*4+1] += vol*coarseRatio*sol[1];
1200d1cd5e0SMatthew G. Knepley         metric[v*4+2] += vol*coarseRatio*sol[1];
1210d1cd5e0SMatthew G. Knepley         metric[v*4+3] += vol*coarseRatio*sol[2];
1220d1cd5e0SMatthew G. Knepley       } else {
1230d1cd5e0SMatthew G. Knepley         metric[v*9+0] += vol*coarseRatio*sol[0];
1240d1cd5e0SMatthew G. Knepley         metric[v*9+1] += vol*coarseRatio*sol[1];
1250d1cd5e0SMatthew G. Knepley         metric[v*9+3] += vol*coarseRatio*sol[1];
1260d1cd5e0SMatthew G. Knepley         metric[v*9+2] += vol*coarseRatio*sol[2];
1270d1cd5e0SMatthew G. Knepley         metric[v*9+6] += vol*coarseRatio*sol[2];
1280d1cd5e0SMatthew G. Knepley         metric[v*9+4] += vol*coarseRatio*sol[3];
1290d1cd5e0SMatthew G. Knepley         metric[v*9+5] += vol*coarseRatio*sol[4];
1300d1cd5e0SMatthew G. Knepley         metric[v*9+7] += vol*coarseRatio*sol[4];
1310d1cd5e0SMatthew G. Knepley         metric[v*9+8] += vol*coarseRatio*sol[5];
1320d1cd5e0SMatthew G. Knepley       }
1330d1cd5e0SMatthew G. Knepley     }
1349566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(mx, &sol));
1350d1cd5e0SMatthew G. Knepley   }
1360d1cd5e0SMatthew G. Knepley   for (v = 0; v < Nv; ++v) {
1370d1cd5e0SMatthew G. Knepley     const PetscInt *support;
1380d1cd5e0SMatthew G. Knepley     PetscInt        supportSize, s;
1390d1cd5e0SMatthew G. Knepley     PetscReal       vol, totVol = 0.0;
1400d1cd5e0SMatthew G. Knepley 
1419566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupport(udm, v+vStart, &support));
1429566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupportSize(udm, v+vStart, &supportSize));
1439566063dSJacob Faibussowitsch     for (s = 0; s < supportSize; ++s) {PetscCall(DMPlexComputeCellGeometryFVM(dm, support[s], &vol, NULL, NULL)); totVol += vol;}
1440d1cd5e0SMatthew G. Knepley     for (s = 0; s < PetscSqr(dim); ++s) metric[v*PetscSqr(dim)+s] /= totVol;
1450d1cd5e0SMatthew G. Knepley   }
1469566063dSJacob Faibussowitsch   PetscCall(PetscFree(eqns));
1479566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(*metricVec, &metric));
1489566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mx));
1499566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mb));
1509566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1519566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&udm));
1520d1cd5e0SMatthew G. Knepley   PetscFunctionReturn(0);
1530d1cd5e0SMatthew G. Knepley }
1540d1cd5e0SMatthew G. Knepley 
1553a074057SBarry Smith /*
1563a074057SBarry Smith    Contains the list of registered DMPlexGenerators routines
1573a074057SBarry Smith */
1589fe9e680SJoe Wallwork PetscErrorCode DMPlexRefine_Internal(DM dm, PETSC_UNUSED Vec metric, DMLabel adaptLabel, PETSC_UNUSED DMLabel rgLabel, DM *dmRefined)
1590d1cd5e0SMatthew G. Knepley {
160c0517cd5SMatthew G. Knepley   DMGeneratorFunctionList fl;
161800b850cSBarry Smith   PetscErrorCode        (*refine)(DM,PetscReal*,DM*);
1629fe9e680SJoe Wallwork   PetscErrorCode        (*adapt)(DM,Vec,DMLabel,DMLabel,DM*);
1633f2a96e3SMatthew G. Knepley   PetscErrorCode        (*refinementFunc)(const PetscReal [], PetscReal *);
1643f2a96e3SMatthew G. Knepley   char                    genname[PETSC_MAX_PATH_LEN], *name = NULL;
1653f2a96e3SMatthew G. Knepley   PetscReal               refinementLimit;
1663a074057SBarry Smith   PetscReal              *maxVolumes;
1673f2a96e3SMatthew G. Knepley   PetscInt                dim, cStart, cEnd, c;
1683f2a96e3SMatthew G. Knepley   PetscBool               flg, flg2, localized;
1690d1cd5e0SMatthew G. Knepley 
1700d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
1719566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocalized(dm, &localized));
1729566063dSJacob Faibussowitsch   PetscCall(DMPlexGetRefinementLimit(dm, &refinementLimit));
1739566063dSJacob Faibussowitsch   PetscCall(DMPlexGetRefinementFunction(dm, &refinementFunc));
1740d1cd5e0SMatthew G. Knepley   if (refinementLimit == 0.0 && !refinementFunc && !adaptLabel) PetscFunctionReturn(0);
1759566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1769566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_adaptor", genname, sizeof(genname), &flg));
1780d1cd5e0SMatthew G. Knepley   if (flg) name = genname;
1793f2a96e3SMatthew G. Knepley   else {
1809566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_generator", genname, sizeof(genname), &flg2));
1813f2a96e3SMatthew G. Knepley     if (flg2) name = genname;
1823f2a96e3SMatthew G. Knepley   }
1837d99616aSKarl Rupp 
184c0517cd5SMatthew G. Knepley   fl = DMGenerateList;
1853a074057SBarry Smith   if (name) {
1863a074057SBarry Smith     while (fl) {
1879566063dSJacob Faibussowitsch       PetscCall(PetscStrcmp(fl->name,name,&flg));
1883a074057SBarry Smith       if (flg) {
1893a074057SBarry Smith         refine = fl->refine;
190c0517cd5SMatthew G. Knepley         adapt  = fl->adapt;
1913a074057SBarry Smith         goto gotit;
1923a074057SBarry Smith       }
1933a074057SBarry Smith       fl = fl->next;
1943a074057SBarry Smith     }
19598921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Grid refiner %s not registered",name);
1963a074057SBarry Smith   } else {
1973a074057SBarry Smith     while (fl) {
198a12d352dSMatthew G. Knepley       if (fl->dim < 0 || dim-1 == fl->dim) {
1993a074057SBarry Smith         refine = fl->refine;
200c0517cd5SMatthew G. Knepley         adapt  = fl->adapt;
2013a074057SBarry Smith         goto gotit;
2023a074057SBarry Smith       }
2033a074057SBarry Smith       fl = fl->next;
2043a074057SBarry Smith     }
20563a3b9bcSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No grid refiner of dimension %" PetscInt_FMT " registered",dim);
2063a074057SBarry Smith   }
2073a074057SBarry Smith 
2083f2a96e3SMatthew G. Knepley   gotit:
2093f2a96e3SMatthew G. Knepley   switch (dim) {
210a12d352dSMatthew G. Knepley     case 1:
2113a074057SBarry Smith     case 2:
2120d1cd5e0SMatthew G. Knepley     case 3:
2133f2a96e3SMatthew G. Knepley       if (adapt) {
2149566063dSJacob Faibussowitsch         PetscCall((*adapt)(dm, NULL, adaptLabel, NULL, dmRefined));
2153f2a96e3SMatthew G. Knepley       } else {
2169566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(cEnd - cStart, &maxVolumes));
2170d1cd5e0SMatthew G. Knepley         if (adaptLabel) {
2189566063dSJacob Faibussowitsch           PetscCall(DMPlexLabelToVolumeConstraint(dm, adaptLabel, cStart, cEnd, PETSC_DEFAULT, maxVolumes));
2190d1cd5e0SMatthew G. Knepley         } else if (refinementFunc) {
2200d1cd5e0SMatthew G. Knepley           for (c = cStart; c < cEnd; ++c) {
2210d1cd5e0SMatthew G. Knepley             PetscReal vol, centroid[3];
2220d1cd5e0SMatthew G. Knepley 
2239566063dSJacob Faibussowitsch             PetscCall(DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL));
2249566063dSJacob Faibussowitsch             PetscCall((*refinementFunc)(centroid, &maxVolumes[c-cStart]));
2250d1cd5e0SMatthew G. Knepley           }
2260d1cd5e0SMatthew G. Knepley         } else {
2270d1cd5e0SMatthew G. Knepley           for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit;
2280d1cd5e0SMatthew G. Knepley         }
2299566063dSJacob Faibussowitsch         PetscCall((*refine)(dm, maxVolumes, dmRefined));
2309566063dSJacob Faibussowitsch         PetscCall(PetscFree(maxVolumes));
2313f2a96e3SMatthew G. Knepley       }
2320d1cd5e0SMatthew G. Knepley       break;
23363a3b9bcSJacob Faibussowitsch     default: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Mesh refinement in dimension %" PetscInt_FMT " is not supported.", dim);
2340d1cd5e0SMatthew G. Knepley   }
2359566063dSJacob Faibussowitsch   PetscCall(DMCopyDisc(dm, *dmRefined));
236*5de52c6dSVaclav Hapla   PetscCall(DMPlexCopy_Internal(dm, PETSC_TRUE, PETSC_TRUE, *dmRefined));
2379566063dSJacob Faibussowitsch   if (localized) PetscCall(DMLocalizeCoordinates(*dmRefined));
2380d1cd5e0SMatthew G. Knepley   PetscFunctionReturn(0);
2390d1cd5e0SMatthew G. Knepley }
2400d1cd5e0SMatthew G. Knepley 
2419fe9e680SJoe Wallwork PetscErrorCode DMPlexCoarsen_Internal(DM dm, PETSC_UNUSED Vec metric, DMLabel adaptLabel, PETSC_UNUSED DMLabel rgLabel, DM *dmCoarsened)
2420d1cd5e0SMatthew G. Knepley {
2430d1cd5e0SMatthew G. Knepley   Vec            metricVec;
2440d1cd5e0SMatthew G. Knepley   PetscInt       cStart, cEnd, vStart, vEnd;
2450d1cd5e0SMatthew G. Knepley   DMLabel        bdLabel = NULL;
2469fe9e680SJoe Wallwork   char           bdLabelName[PETSC_MAX_PATH_LEN], rgLabelName[PETSC_MAX_PATH_LEN];
2470d1cd5e0SMatthew G. Knepley   PetscBool      localized, flg;
2480d1cd5e0SMatthew G. Knepley 
2490d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
2509566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocalized(dm, &localized));
2519566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
2529566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
2539566063dSJacob Faibussowitsch   PetscCall(DMPlexLabelToMetricConstraint(dm, adaptLabel, cStart, cEnd, vStart, vEnd, PETSC_DEFAULT, &metricVec));
2549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, dm->hdr.prefix, "-dm_plex_coarsen_bd_label", bdLabelName, sizeof(bdLabelName), &flg));
2559566063dSJacob Faibussowitsch   if (flg) PetscCall(DMGetLabel(dm, bdLabelName, &bdLabel));
2569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, dm->hdr.prefix, "-dm_plex_coarsen_rg_label", rgLabelName, sizeof(rgLabelName), &flg));
2579566063dSJacob Faibussowitsch   if (flg) PetscCall(DMGetLabel(dm, rgLabelName, &rgLabel));
2589566063dSJacob Faibussowitsch   PetscCall(DMAdaptMetric(dm, metricVec, bdLabel, rgLabel, dmCoarsened));
2599566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&metricVec));
2609566063dSJacob Faibussowitsch   PetscCall(DMCopyDisc(dm, *dmCoarsened));
261*5de52c6dSVaclav Hapla   PetscCall(DMPlexCopy_Internal(dm, PETSC_TRUE, PETSC_TRUE, *dmCoarsened));
2629566063dSJacob Faibussowitsch   if (localized) PetscCall(DMLocalizeCoordinates(*dmCoarsened));
2630d1cd5e0SMatthew G. Knepley   PetscFunctionReturn(0);
2640d1cd5e0SMatthew G. Knepley }
2650d1cd5e0SMatthew G. Knepley 
2669fe9e680SJoe Wallwork PetscErrorCode DMAdaptLabel_Plex(DM dm, PETSC_UNUSED Vec metric, DMLabel adaptLabel, PETSC_UNUSED DMLabel rgLabel, DM *dmAdapted)
2670d1cd5e0SMatthew G. Knepley {
2680d1cd5e0SMatthew G. Knepley   IS              flagIS;
2690d1cd5e0SMatthew G. Knepley   const PetscInt *flags;
2700d1cd5e0SMatthew G. Knepley   PetscInt        defFlag, minFlag, maxFlag, numFlags, f;
2710d1cd5e0SMatthew G. Knepley 
2720d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
2739566063dSJacob Faibussowitsch   PetscCall(DMLabelGetDefaultValue(adaptLabel, &defFlag));
2740d1cd5e0SMatthew G. Knepley   minFlag = defFlag;
2750d1cd5e0SMatthew G. Knepley   maxFlag = defFlag;
2769566063dSJacob Faibussowitsch   PetscCall(DMLabelGetValueIS(adaptLabel, &flagIS));
2779566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(flagIS, &numFlags));
2789566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(flagIS, &flags));
2790d1cd5e0SMatthew G. Knepley   for (f = 0; f < numFlags; ++f) {
2800d1cd5e0SMatthew G. Knepley     const PetscInt flag = flags[f];
2810d1cd5e0SMatthew G. Knepley 
2820d1cd5e0SMatthew G. Knepley     minFlag = PetscMin(minFlag, flag);
2830d1cd5e0SMatthew G. Knepley     maxFlag = PetscMax(maxFlag, flag);
2840d1cd5e0SMatthew G. Knepley   }
2859566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(flagIS, &flags));
2869566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&flagIS));
2870d1cd5e0SMatthew G. Knepley   {
2880d1cd5e0SMatthew G. Knepley     PetscInt minMaxFlag[2], minMaxFlagGlobal[2];
2890d1cd5e0SMatthew G. Knepley 
2900d1cd5e0SMatthew G. Knepley     minMaxFlag[0] =  minFlag;
2910d1cd5e0SMatthew G. Knepley     minMaxFlag[1] = -maxFlag;
2929566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allreduce(minMaxFlag, minMaxFlagGlobal, 2, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm)));
2930d1cd5e0SMatthew G. Knepley     minFlag =  minMaxFlagGlobal[0];
2940d1cd5e0SMatthew G. Knepley     maxFlag = -minMaxFlagGlobal[1];
2950d1cd5e0SMatthew G. Knepley   }
2960d1cd5e0SMatthew G. Knepley   if (minFlag == maxFlag) {
2970d1cd5e0SMatthew G. Knepley     switch (minFlag) {
2980d1cd5e0SMatthew G. Knepley     case DM_ADAPT_DETERMINE:
2990d1cd5e0SMatthew G. Knepley       *dmAdapted = NULL;break;
3000d1cd5e0SMatthew G. Knepley     case DM_ADAPT_REFINE:
3019566063dSJacob Faibussowitsch       PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE));
3029566063dSJacob Faibussowitsch       PetscCall(DMRefine(dm, MPI_COMM_NULL, dmAdapted));break;
3030d1cd5e0SMatthew G. Knepley     case DM_ADAPT_COARSEN:
3049566063dSJacob Faibussowitsch       PetscCall(DMCoarsen(dm, MPI_COMM_NULL, dmAdapted));break;
3050d1cd5e0SMatthew G. Knepley     default:
30663a3b9bcSJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"DMPlex does not support refinement flag %" PetscInt_FMT, minFlag);
3070d1cd5e0SMatthew G. Knepley     }
3080d1cd5e0SMatthew G. Knepley   } else {
3099566063dSJacob Faibussowitsch     PetscCall(DMPlexSetRefinementUniform(dm, PETSC_FALSE));
3109566063dSJacob Faibussowitsch     PetscCall(DMPlexRefine_Internal(dm, NULL, adaptLabel, NULL, dmAdapted));
3110d1cd5e0SMatthew G. Knepley   }
3120d1cd5e0SMatthew G. Knepley   PetscFunctionReturn(0);
3130d1cd5e0SMatthew G. Knepley }
314