xref: /petsc/src/dm/impls/plex/plexadapt.c (revision 98921bda46e76d7aaed9e0138c5ff9d0ce93f355)
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   PetscErrorCode ierr;
70d1cd5e0SMatthew G. Knepley 
80d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
90d1cd5e0SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
100d1cd5e0SMatthew G. Knepley   refRatio = refRatio == PETSC_DEFAULT ? (PetscReal) ((PetscInt) 1 << dim) : refRatio;
110d1cd5e0SMatthew G. Knepley   for (c = cStart; c < cEnd; c++) {
120d1cd5e0SMatthew G. Knepley     PetscReal vol;
130d1cd5e0SMatthew G. Knepley     PetscInt  closureSize = 0, cl;
140d1cd5e0SMatthew G. Knepley     PetscInt *closure     = NULL;
150d1cd5e0SMatthew G. Knepley     PetscBool anyRefine   = PETSC_FALSE;
160d1cd5e0SMatthew G. Knepley     PetscBool anyCoarsen  = PETSC_FALSE;
170d1cd5e0SMatthew G. Knepley     PetscBool anyKeep     = PETSC_FALSE;
180d1cd5e0SMatthew G. Knepley 
190d1cd5e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, NULL, NULL);CHKERRQ(ierr);
200d1cd5e0SMatthew G. Knepley     maxVolumes[c - cStart] = vol;
210d1cd5e0SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
220d1cd5e0SMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
230d1cd5e0SMatthew G. Knepley       const PetscInt point = closure[cl];
240d1cd5e0SMatthew G. Knepley       PetscInt       refFlag;
250d1cd5e0SMatthew G. Knepley 
260d1cd5e0SMatthew G. Knepley       ierr = DMLabelGetValue(adaptLabel, point, &refFlag);CHKERRQ(ierr);
270d1cd5e0SMatthew G. Knepley       switch (refFlag) {
280d1cd5e0SMatthew G. Knepley       case DM_ADAPT_REFINE:
290d1cd5e0SMatthew G. Knepley         anyRefine  = PETSC_TRUE;break;
300d1cd5e0SMatthew G. Knepley       case DM_ADAPT_COARSEN:
310d1cd5e0SMatthew G. Knepley         anyCoarsen = PETSC_TRUE;break;
320d1cd5e0SMatthew G. Knepley       case DM_ADAPT_KEEP:
330d1cd5e0SMatthew G. Knepley         anyKeep    = PETSC_TRUE;break;
340d1cd5e0SMatthew G. Knepley       case DM_ADAPT_DETERMINE:
350d1cd5e0SMatthew G. Knepley         break;
360d1cd5e0SMatthew G. Knepley       default:
37*98921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "DMPlex does not support refinement flag %D", refFlag);
380d1cd5e0SMatthew G. Knepley       }
390d1cd5e0SMatthew G. Knepley       if (anyRefine) break;
400d1cd5e0SMatthew G. Knepley     }
410d1cd5e0SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
420d1cd5e0SMatthew G. Knepley     if (anyRefine) {
430d1cd5e0SMatthew G. Knepley       maxVolumes[c - cStart] = vol / refRatio;
440d1cd5e0SMatthew G. Knepley     } else if (anyKeep) {
450d1cd5e0SMatthew G. Knepley       maxVolumes[c - cStart] = vol;
460d1cd5e0SMatthew G. Knepley     } else if (anyCoarsen) {
470d1cd5e0SMatthew G. Knepley       maxVolumes[c - cStart] = vol * refRatio;
480d1cd5e0SMatthew G. Knepley     }
490d1cd5e0SMatthew G. Knepley   }
500d1cd5e0SMatthew G. Knepley   PetscFunctionReturn(0);
510d1cd5e0SMatthew G. Knepley }
520d1cd5e0SMatthew G. Knepley 
530d1cd5e0SMatthew G. Knepley static PetscErrorCode DMPlexLabelToMetricConstraint(DM dm, DMLabel adaptLabel, PetscInt cStart, PetscInt cEnd, PetscInt vStart, PetscInt vEnd, PetscReal refRatio, Vec *metricVec)
540d1cd5e0SMatthew G. Knepley {
550d1cd5e0SMatthew G. Knepley   DM              udm, coordDM;
560d1cd5e0SMatthew G. Knepley   PetscSection    coordSection;
570d1cd5e0SMatthew G. Knepley   Vec             coordinates, mb, mx;
580d1cd5e0SMatthew G. Knepley   Mat             A;
590d1cd5e0SMatthew G. Knepley   PetscScalar    *metric, *eqns;
600d1cd5e0SMatthew G. Knepley   const PetscReal coarseRatio = refRatio == PETSC_DEFAULT ? PetscSqr(0.5) : 1/refRatio;
610d1cd5e0SMatthew G. Knepley   PetscInt        dim, Nv, Neq, c, v;
620d1cd5e0SMatthew G. Knepley   PetscErrorCode  ierr;
630d1cd5e0SMatthew G. Knepley 
640d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
650d1cd5e0SMatthew G. Knepley   ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr);
660d1cd5e0SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
670d1cd5e0SMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &coordDM);CHKERRQ(ierr);
6892fd8e1eSJed Brown   ierr = DMGetLocalSection(coordDM, &coordSection);CHKERRQ(ierr);
690d1cd5e0SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
700d1cd5e0SMatthew G. Knepley   Nv   = vEnd - vStart;
710d1cd5e0SMatthew G. Knepley   ierr = VecCreateSeq(PETSC_COMM_SELF, Nv*PetscSqr(dim), metricVec);CHKERRQ(ierr);
720d1cd5e0SMatthew G. Knepley   ierr = VecGetArray(*metricVec, &metric);CHKERRQ(ierr);
730d1cd5e0SMatthew G. Knepley   Neq  = (dim*(dim+1))/2;
740d1cd5e0SMatthew G. Knepley   ierr = PetscMalloc1(PetscSqr(Neq), &eqns);CHKERRQ(ierr);
750d1cd5e0SMatthew G. Knepley   ierr = MatCreateSeqDense(PETSC_COMM_SELF, Neq, Neq, eqns, &A);CHKERRQ(ierr);
760d1cd5e0SMatthew G. Knepley   ierr = MatCreateVecs(A, &mx, &mb);CHKERRQ(ierr);
770d1cd5e0SMatthew G. Knepley   ierr = VecSet(mb, 1.0);CHKERRQ(ierr);
780d1cd5e0SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
790d1cd5e0SMatthew G. Knepley     const PetscScalar *sol;
800d1cd5e0SMatthew G. Knepley     PetscScalar       *cellCoords = NULL;
810d1cd5e0SMatthew G. Knepley     PetscReal          e[3], vol;
820d1cd5e0SMatthew G. Knepley     const PetscInt    *cone;
830d1cd5e0SMatthew G. Knepley     PetscInt           coneSize, cl, i, j, d, r;
840d1cd5e0SMatthew G. Knepley 
850d1cd5e0SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, c, NULL, &cellCoords);CHKERRQ(ierr);
860d1cd5e0SMatthew G. Knepley     /* Only works for simplices */
870d1cd5e0SMatthew G. Knepley     for (i = 0, r = 0; i < dim+1; ++i) {
880d1cd5e0SMatthew G. Knepley       for (j = 0; j < i; ++j, ++r) {
890d1cd5e0SMatthew G. Knepley         for (d = 0; d < dim; ++d) e[d] = PetscRealPart(cellCoords[i*dim+d] - cellCoords[j*dim+d]);
900d1cd5e0SMatthew G. Knepley         /* FORTRAN ORDERING */
910d1cd5e0SMatthew G. Knepley         switch (dim) {
920d1cd5e0SMatthew G. Knepley         case 2:
930d1cd5e0SMatthew G. Knepley           eqns[0*Neq+r] = PetscSqr(e[0]);
940d1cd5e0SMatthew G. Knepley           eqns[1*Neq+r] = 2.0*e[0]*e[1];
950d1cd5e0SMatthew G. Knepley           eqns[2*Neq+r] = PetscSqr(e[1]);
960d1cd5e0SMatthew G. Knepley           break;
970d1cd5e0SMatthew G. Knepley         case 3:
980d1cd5e0SMatthew G. Knepley           eqns[0*Neq+r] = PetscSqr(e[0]);
990d1cd5e0SMatthew G. Knepley           eqns[1*Neq+r] = 2.0*e[0]*e[1];
1000d1cd5e0SMatthew G. Knepley           eqns[2*Neq+r] = 2.0*e[0]*e[2];
1010d1cd5e0SMatthew G. Knepley           eqns[3*Neq+r] = PetscSqr(e[1]);
1020d1cd5e0SMatthew G. Knepley           eqns[4*Neq+r] = 2.0*e[1]*e[2];
1030d1cd5e0SMatthew G. Knepley           eqns[5*Neq+r] = PetscSqr(e[2]);
1040d1cd5e0SMatthew G. Knepley           break;
1050d1cd5e0SMatthew G. Knepley         }
1060d1cd5e0SMatthew G. Knepley       }
1070d1cd5e0SMatthew G. Knepley     }
1080d1cd5e0SMatthew G. Knepley     ierr = MatSetUnfactored(A);CHKERRQ(ierr);
1090d1cd5e0SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, c, NULL, &cellCoords);CHKERRQ(ierr);
1100d1cd5e0SMatthew G. Knepley     ierr = MatLUFactor(A, NULL, NULL, NULL);CHKERRQ(ierr);
1110d1cd5e0SMatthew G. Knepley     ierr = MatSolve(A, mb, mx);CHKERRQ(ierr);
1120d1cd5e0SMatthew G. Knepley     ierr = VecGetArrayRead(mx, &sol);CHKERRQ(ierr);
1130d1cd5e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, NULL, NULL);CHKERRQ(ierr);
1140d1cd5e0SMatthew G. Knepley     ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr);
1150d1cd5e0SMatthew G. Knepley     ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr);
1160d1cd5e0SMatthew G. Knepley     for (cl = 0; cl < coneSize; ++cl) {
1170d1cd5e0SMatthew G. Knepley       const PetscInt v = cone[cl] - vStart;
1180d1cd5e0SMatthew G. Knepley 
1190d1cd5e0SMatthew G. Knepley       if (dim == 2) {
1200d1cd5e0SMatthew G. Knepley         metric[v*4+0] += vol*coarseRatio*sol[0];
1210d1cd5e0SMatthew G. Knepley         metric[v*4+1] += vol*coarseRatio*sol[1];
1220d1cd5e0SMatthew G. Knepley         metric[v*4+2] += vol*coarseRatio*sol[1];
1230d1cd5e0SMatthew G. Knepley         metric[v*4+3] += vol*coarseRatio*sol[2];
1240d1cd5e0SMatthew G. Knepley       } else {
1250d1cd5e0SMatthew G. Knepley         metric[v*9+0] += vol*coarseRatio*sol[0];
1260d1cd5e0SMatthew G. Knepley         metric[v*9+1] += vol*coarseRatio*sol[1];
1270d1cd5e0SMatthew G. Knepley         metric[v*9+3] += vol*coarseRatio*sol[1];
1280d1cd5e0SMatthew G. Knepley         metric[v*9+2] += vol*coarseRatio*sol[2];
1290d1cd5e0SMatthew G. Knepley         metric[v*9+6] += vol*coarseRatio*sol[2];
1300d1cd5e0SMatthew G. Knepley         metric[v*9+4] += vol*coarseRatio*sol[3];
1310d1cd5e0SMatthew G. Knepley         metric[v*9+5] += vol*coarseRatio*sol[4];
1320d1cd5e0SMatthew G. Knepley         metric[v*9+7] += vol*coarseRatio*sol[4];
1330d1cd5e0SMatthew G. Knepley         metric[v*9+8] += vol*coarseRatio*sol[5];
1340d1cd5e0SMatthew G. Knepley       }
1350d1cd5e0SMatthew G. Knepley     }
1360d1cd5e0SMatthew G. Knepley     ierr = VecRestoreArrayRead(mx, &sol);CHKERRQ(ierr);
1370d1cd5e0SMatthew G. Knepley   }
1380d1cd5e0SMatthew G. Knepley   for (v = 0; v < Nv; ++v) {
1390d1cd5e0SMatthew G. Knepley     const PetscInt *support;
1400d1cd5e0SMatthew G. Knepley     PetscInt        supportSize, s;
1410d1cd5e0SMatthew G. Knepley     PetscReal       vol, totVol = 0.0;
1420d1cd5e0SMatthew G. Knepley 
1430d1cd5e0SMatthew G. Knepley     ierr = DMPlexGetSupport(udm, v+vStart, &support);CHKERRQ(ierr);
1440d1cd5e0SMatthew G. Knepley     ierr = DMPlexGetSupportSize(udm, v+vStart, &supportSize);CHKERRQ(ierr);
1450d1cd5e0SMatthew G. Knepley     for (s = 0; s < supportSize; ++s) {ierr = DMPlexComputeCellGeometryFVM(dm, support[s], &vol, NULL, NULL);CHKERRQ(ierr); totVol += vol;}
1460d1cd5e0SMatthew G. Knepley     for (s = 0; s < PetscSqr(dim); ++s) metric[v*PetscSqr(dim)+s] /= totVol;
1470d1cd5e0SMatthew G. Knepley   }
1480d1cd5e0SMatthew G. Knepley   ierr = PetscFree(eqns);CHKERRQ(ierr);
1490d1cd5e0SMatthew G. Knepley   ierr = VecRestoreArray(*metricVec, &metric);CHKERRQ(ierr);
1500d1cd5e0SMatthew G. Knepley   ierr = VecDestroy(&mx);CHKERRQ(ierr);
1510d1cd5e0SMatthew G. Knepley   ierr = VecDestroy(&mb);CHKERRQ(ierr);
1520d1cd5e0SMatthew G. Knepley   ierr = MatDestroy(&A);CHKERRQ(ierr);
1530d1cd5e0SMatthew G. Knepley   ierr = DMDestroy(&udm);CHKERRQ(ierr);
1540d1cd5e0SMatthew G. Knepley   PetscFunctionReturn(0);
1550d1cd5e0SMatthew G. Knepley }
1560d1cd5e0SMatthew G. Knepley 
1573a074057SBarry Smith /*
1583a074057SBarry Smith    Contains the list of registered DMPlexGenerators routines
1593a074057SBarry Smith */
1609fe9e680SJoe Wallwork PetscErrorCode DMPlexRefine_Internal(DM dm, PETSC_UNUSED Vec metric, DMLabel adaptLabel, PETSC_UNUSED DMLabel rgLabel, DM *dmRefined)
1610d1cd5e0SMatthew G. Knepley {
162c0517cd5SMatthew G. Knepley   DMGeneratorFunctionList fl;
163800b850cSBarry Smith   PetscErrorCode        (*refine)(DM,PetscReal*,DM*);
1649fe9e680SJoe Wallwork   PetscErrorCode        (*adapt)(DM,Vec,DMLabel,DMLabel,DM*);
1653f2a96e3SMatthew G. Knepley   PetscErrorCode        (*refinementFunc)(const PetscReal [], PetscReal *);
1663f2a96e3SMatthew G. Knepley   char                    genname[PETSC_MAX_PATH_LEN], *name = NULL;
1673f2a96e3SMatthew G. Knepley   PetscReal               refinementLimit;
1683a074057SBarry Smith   PetscReal              *maxVolumes;
1693f2a96e3SMatthew G. Knepley   PetscInt                dim, cStart, cEnd, c;
1703f2a96e3SMatthew G. Knepley   PetscBool               flg, flg2, localized;
1713f2a96e3SMatthew G. Knepley   PetscErrorCode          ierr;
1720d1cd5e0SMatthew G. Knepley 
1730d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
1740d1cd5e0SMatthew G. Knepley   ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr);
1750d1cd5e0SMatthew G. Knepley   ierr = DMPlexGetRefinementLimit(dm, &refinementLimit);CHKERRQ(ierr);
1760d1cd5e0SMatthew G. Knepley   ierr = DMPlexGetRefinementFunction(dm, &refinementFunc);CHKERRQ(ierr);
1770d1cd5e0SMatthew G. Knepley   if (refinementLimit == 0.0 && !refinementFunc && !adaptLabel) PetscFunctionReturn(0);
1780d1cd5e0SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1790d1cd5e0SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
180c0517cd5SMatthew G. Knepley   ierr = PetscOptionsGetString(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_adaptor", genname, sizeof(genname), &flg);CHKERRQ(ierr);
1810d1cd5e0SMatthew G. Knepley   if (flg) name = genname;
1823f2a96e3SMatthew G. Knepley   else {
183c0517cd5SMatthew G. Knepley     ierr = PetscOptionsGetString(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_generator", genname, sizeof(genname), &flg2);CHKERRQ(ierr);
1843f2a96e3SMatthew G. Knepley     if (flg2) name = genname;
1853f2a96e3SMatthew G. Knepley   }
1867d99616aSKarl Rupp 
187c0517cd5SMatthew G. Knepley   fl = DMGenerateList;
1883a074057SBarry Smith   if (name) {
1893a074057SBarry Smith     while (fl) {
1903a074057SBarry Smith       ierr = PetscStrcmp(fl->name,name,&flg);CHKERRQ(ierr);
1913a074057SBarry Smith       if (flg) {
1923a074057SBarry Smith         refine = fl->refine;
193c0517cd5SMatthew G. Knepley         adapt  = fl->adapt;
1943a074057SBarry Smith         goto gotit;
1953a074057SBarry Smith       }
1963a074057SBarry Smith       fl = fl->next;
1973a074057SBarry Smith     }
198*98921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Grid refiner %s not registered",name);
1993a074057SBarry Smith   } else {
2003a074057SBarry Smith     while (fl) {
201a12d352dSMatthew G. Knepley       if (fl->dim < 0 || dim-1 == fl->dim) {
2023a074057SBarry Smith         refine = fl->refine;
203c0517cd5SMatthew G. Knepley         adapt  = fl->adapt;
2043a074057SBarry Smith         goto gotit;
2053a074057SBarry Smith       }
2063a074057SBarry Smith       fl = fl->next;
2073a074057SBarry Smith     }
208*98921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No grid refiner of dimension %D registered",dim);
2093a074057SBarry Smith   }
2103a074057SBarry Smith 
2113f2a96e3SMatthew G. Knepley   gotit:
2123f2a96e3SMatthew G. Knepley   switch (dim) {
213a12d352dSMatthew G. Knepley     case 1:
2143a074057SBarry Smith     case 2:
2150d1cd5e0SMatthew G. Knepley     case 3:
2163f2a96e3SMatthew G. Knepley       if (adapt) {
2179fe9e680SJoe Wallwork         ierr = (*adapt)(dm, NULL, adaptLabel, NULL, dmRefined);CHKERRQ(ierr);
2183f2a96e3SMatthew G. Knepley       } else {
2190d1cd5e0SMatthew G. Knepley         ierr = PetscMalloc1(cEnd - cStart, &maxVolumes);CHKERRQ(ierr);
2200d1cd5e0SMatthew G. Knepley         if (adaptLabel) {
2210d1cd5e0SMatthew G. Knepley           ierr = DMPlexLabelToVolumeConstraint(dm, adaptLabel, cStart, cEnd, PETSC_DEFAULT, maxVolumes);CHKERRQ(ierr);
2220d1cd5e0SMatthew G. Knepley         } else if (refinementFunc) {
2230d1cd5e0SMatthew G. Knepley           for (c = cStart; c < cEnd; ++c) {
2240d1cd5e0SMatthew G. Knepley             PetscReal vol, centroid[3];
2250d1cd5e0SMatthew G. Knepley 
2260d1cd5e0SMatthew G. Knepley             ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL);CHKERRQ(ierr);
2270d1cd5e0SMatthew G. Knepley             ierr = (*refinementFunc)(centroid, &maxVolumes[c-cStart]);CHKERRQ(ierr);
2280d1cd5e0SMatthew G. Knepley           }
2290d1cd5e0SMatthew G. Knepley         } else {
2300d1cd5e0SMatthew G. Knepley           for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit;
2310d1cd5e0SMatthew G. Knepley         }
2323a074057SBarry Smith         ierr = (*refine)(dm, maxVolumes, dmRefined);CHKERRQ(ierr);
2330d1cd5e0SMatthew G. Knepley         ierr = PetscFree(maxVolumes);CHKERRQ(ierr);
2343f2a96e3SMatthew G. Knepley       }
2350d1cd5e0SMatthew G. Knepley       break;
236*98921bdaSJacob Faibussowitsch     default: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Mesh refinement in dimension %D is not supported.", dim);
2370d1cd5e0SMatthew G. Knepley   }
238a12d352dSMatthew G. Knepley   ((DM_Plex *) (*dmRefined)->data)->useHashLocation = ((DM_Plex *) dm->data)->useHashLocation;
23945480ffeSMatthew G. Knepley   ierr = DMCopyDisc(dm, *dmRefined);CHKERRQ(ierr);
2400d1cd5e0SMatthew G. Knepley   if (localized) {ierr = DMLocalizeCoordinates(*dmRefined);CHKERRQ(ierr);}
2410d1cd5e0SMatthew G. Knepley   PetscFunctionReturn(0);
2420d1cd5e0SMatthew G. Knepley }
2430d1cd5e0SMatthew G. Knepley 
2449fe9e680SJoe Wallwork PetscErrorCode DMPlexCoarsen_Internal(DM dm, PETSC_UNUSED Vec metric, DMLabel adaptLabel, PETSC_UNUSED DMLabel rgLabel, DM *dmCoarsened)
2450d1cd5e0SMatthew G. Knepley {
2460d1cd5e0SMatthew G. Knepley   Vec            metricVec;
2470d1cd5e0SMatthew G. Knepley   PetscInt       cStart, cEnd, vStart, vEnd;
2480d1cd5e0SMatthew G. Knepley   DMLabel        bdLabel = NULL;
2499fe9e680SJoe Wallwork   char           bdLabelName[PETSC_MAX_PATH_LEN], rgLabelName[PETSC_MAX_PATH_LEN];
2500d1cd5e0SMatthew G. Knepley   PetscBool      localized, flg;
2510d1cd5e0SMatthew G. Knepley   PetscErrorCode ierr;
2520d1cd5e0SMatthew G. Knepley 
2530d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
2540d1cd5e0SMatthew G. Knepley   ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr);
2550d1cd5e0SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2560d1cd5e0SMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
2570d1cd5e0SMatthew G. Knepley   ierr = DMPlexLabelToMetricConstraint(dm, adaptLabel, cStart, cEnd, vStart, vEnd, PETSC_DEFAULT, &metricVec);CHKERRQ(ierr);
258589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL, dm->hdr.prefix, "-dm_plex_coarsen_bd_label", bdLabelName, sizeof(bdLabelName), &flg);CHKERRQ(ierr);
2590d1cd5e0SMatthew G. Knepley   if (flg) {ierr = DMGetLabel(dm, bdLabelName, &bdLabel);CHKERRQ(ierr);}
2609fe9e680SJoe Wallwork   ierr = PetscOptionsGetString(NULL, dm->hdr.prefix, "-dm_plex_coarsen_rg_label", rgLabelName, sizeof(rgLabelName), &flg);CHKERRQ(ierr);
2619fe9e680SJoe Wallwork   if (flg) {ierr = DMGetLabel(dm, rgLabelName, &rgLabel);CHKERRQ(ierr);}
2629fe9e680SJoe Wallwork   ierr = DMAdaptMetric(dm, metricVec, bdLabel, rgLabel, dmCoarsened);CHKERRQ(ierr);
2630d1cd5e0SMatthew G. Knepley   ierr = VecDestroy(&metricVec);CHKERRQ(ierr);
264a12d352dSMatthew G. Knepley   ((DM_Plex *) (*dmCoarsened)->data)->useHashLocation = ((DM_Plex *) dm->data)->useHashLocation;
265a12d352dSMatthew G. Knepley   ierr = DMCopyDisc(dm, *dmCoarsened);CHKERRQ(ierr);
2660d1cd5e0SMatthew G. Knepley   if (localized) {ierr = DMLocalizeCoordinates(*dmCoarsened);CHKERRQ(ierr);}
2670d1cd5e0SMatthew G. Knepley   PetscFunctionReturn(0);
2680d1cd5e0SMatthew G. Knepley }
2690d1cd5e0SMatthew G. Knepley 
2709fe9e680SJoe Wallwork PetscErrorCode DMAdaptLabel_Plex(DM dm, PETSC_UNUSED Vec metric, DMLabel adaptLabel, PETSC_UNUSED DMLabel rgLabel, DM *dmAdapted)
2710d1cd5e0SMatthew G. Knepley {
2720d1cd5e0SMatthew G. Knepley   IS              flagIS;
2730d1cd5e0SMatthew G. Knepley   const PetscInt *flags;
2740d1cd5e0SMatthew G. Knepley   PetscInt        defFlag, minFlag, maxFlag, numFlags, f;
2750d1cd5e0SMatthew G. Knepley   PetscErrorCode  ierr;
2760d1cd5e0SMatthew G. Knepley 
2770d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
2780d1cd5e0SMatthew G. Knepley   ierr = DMLabelGetDefaultValue(adaptLabel, &defFlag);CHKERRQ(ierr);
2790d1cd5e0SMatthew G. Knepley   minFlag = defFlag;
2800d1cd5e0SMatthew G. Knepley   maxFlag = defFlag;
2810d1cd5e0SMatthew G. Knepley   ierr = DMLabelGetValueIS(adaptLabel, &flagIS);CHKERRQ(ierr);
2820d1cd5e0SMatthew G. Knepley   ierr = ISGetLocalSize(flagIS, &numFlags);CHKERRQ(ierr);
2830d1cd5e0SMatthew G. Knepley   ierr = ISGetIndices(flagIS, &flags);CHKERRQ(ierr);
2840d1cd5e0SMatthew G. Knepley   for (f = 0; f < numFlags; ++f) {
2850d1cd5e0SMatthew G. Knepley     const PetscInt flag = flags[f];
2860d1cd5e0SMatthew G. Knepley 
2870d1cd5e0SMatthew G. Knepley     minFlag = PetscMin(minFlag, flag);
2880d1cd5e0SMatthew G. Knepley     maxFlag = PetscMax(maxFlag, flag);
2890d1cd5e0SMatthew G. Knepley   }
2900d1cd5e0SMatthew G. Knepley   ierr = ISRestoreIndices(flagIS, &flags);CHKERRQ(ierr);
2910d1cd5e0SMatthew G. Knepley   ierr = ISDestroy(&flagIS);CHKERRQ(ierr);
2920d1cd5e0SMatthew G. Knepley   {
2930d1cd5e0SMatthew G. Knepley     PetscInt minMaxFlag[2], minMaxFlagGlobal[2];
2940d1cd5e0SMatthew G. Knepley 
2950d1cd5e0SMatthew G. Knepley     minMaxFlag[0] =  minFlag;
2960d1cd5e0SMatthew G. Knepley     minMaxFlag[1] = -maxFlag;
297ffc4695bSBarry Smith     ierr = MPI_Allreduce(minMaxFlag, minMaxFlagGlobal, 2, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm));CHKERRMPI(ierr);
2980d1cd5e0SMatthew G. Knepley     minFlag =  minMaxFlagGlobal[0];
2990d1cd5e0SMatthew G. Knepley     maxFlag = -minMaxFlagGlobal[1];
3000d1cd5e0SMatthew G. Knepley   }
3010d1cd5e0SMatthew G. Knepley   if (minFlag == maxFlag) {
3020d1cd5e0SMatthew G. Knepley     switch (minFlag) {
3030d1cd5e0SMatthew G. Knepley     case DM_ADAPT_DETERMINE:
3040d1cd5e0SMatthew G. Knepley       *dmAdapted = NULL;break;
3050d1cd5e0SMatthew G. Knepley     case DM_ADAPT_REFINE:
3060d1cd5e0SMatthew G. Knepley       ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE);CHKERRQ(ierr);
3070d1cd5e0SMatthew G. Knepley       ierr = DMRefine(dm, MPI_COMM_NULL, dmAdapted);CHKERRQ(ierr);break;
3080d1cd5e0SMatthew G. Knepley     case DM_ADAPT_COARSEN:
3090d1cd5e0SMatthew G. Knepley       ierr = DMCoarsen(dm, MPI_COMM_NULL, dmAdapted);CHKERRQ(ierr);break;
3100d1cd5e0SMatthew G. Knepley     default:
311*98921bdaSJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"DMPlex does not support refinement flag %D", minFlag);
3120d1cd5e0SMatthew G. Knepley     }
3130d1cd5e0SMatthew G. Knepley   } else {
3140d1cd5e0SMatthew G. Knepley     ierr = DMPlexSetRefinementUniform(dm, PETSC_FALSE);CHKERRQ(ierr);
3159fe9e680SJoe Wallwork     ierr = DMPlexRefine_Internal(dm, NULL, adaptLabel, NULL, dmAdapted);CHKERRQ(ierr);
3160d1cd5e0SMatthew G. Knepley   }
3170d1cd5e0SMatthew G. Knepley   PetscFunctionReturn(0);
3180d1cd5e0SMatthew G. Knepley }
319