10bbe5d31Sbarral #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 20bbe5d31Sbarral 3*9371c9d4SSatish Balay static PetscErrorCode DMPlexLabelToVolumeConstraint(DM dm, DMLabel adaptLabel, PetscInt cStart, PetscInt cEnd, PetscReal refRatio, PetscReal maxVolumes[]) { 40d1cd5e0SMatthew G. Knepley PetscInt dim, c; 50d1cd5e0SMatthew G. Knepley 60d1cd5e0SMatthew G. Knepley PetscFunctionBegin; 79566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 80d1cd5e0SMatthew G. Knepley refRatio = refRatio == PETSC_DEFAULT ? (PetscReal)((PetscInt)1 << dim) : refRatio; 90d1cd5e0SMatthew G. Knepley for (c = cStart; c < cEnd; c++) { 100d1cd5e0SMatthew G. Knepley PetscReal vol; 110d1cd5e0SMatthew G. Knepley PetscInt closureSize = 0, cl; 120d1cd5e0SMatthew G. Knepley PetscInt *closure = NULL; 130d1cd5e0SMatthew G. Knepley PetscBool anyRefine = PETSC_FALSE; 140d1cd5e0SMatthew G. Knepley PetscBool anyCoarsen = PETSC_FALSE; 150d1cd5e0SMatthew G. Knepley PetscBool anyKeep = PETSC_FALSE; 160d1cd5e0SMatthew G. Knepley 179566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dm, c, &vol, NULL, NULL)); 180d1cd5e0SMatthew G. Knepley maxVolumes[c - cStart] = vol; 199566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 200d1cd5e0SMatthew G. Knepley for (cl = 0; cl < closureSize * 2; cl += 2) { 210d1cd5e0SMatthew G. Knepley const PetscInt point = closure[cl]; 220d1cd5e0SMatthew G. Knepley PetscInt refFlag; 230d1cd5e0SMatthew G. Knepley 249566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(adaptLabel, point, &refFlag)); 250d1cd5e0SMatthew G. Knepley switch (refFlag) { 26*9371c9d4SSatish Balay case DM_ADAPT_REFINE: anyRefine = PETSC_TRUE; break; 27*9371c9d4SSatish Balay case DM_ADAPT_COARSEN: anyCoarsen = PETSC_TRUE; break; 28*9371c9d4SSatish Balay case DM_ADAPT_KEEP: anyKeep = PETSC_TRUE; break; 29*9371c9d4SSatish Balay case DM_ADAPT_DETERMINE: break; 30*9371c9d4SSatish Balay default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "DMPlex does not support refinement flag %" PetscInt_FMT, refFlag); 310d1cd5e0SMatthew G. Knepley } 320d1cd5e0SMatthew G. Knepley if (anyRefine) break; 330d1cd5e0SMatthew G. Knepley } 349566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 350d1cd5e0SMatthew G. Knepley if (anyRefine) { 360d1cd5e0SMatthew G. Knepley maxVolumes[c - cStart] = vol / refRatio; 370d1cd5e0SMatthew G. Knepley } else if (anyKeep) { 380d1cd5e0SMatthew G. Knepley maxVolumes[c - cStart] = vol; 390d1cd5e0SMatthew G. Knepley } else if (anyCoarsen) { 400d1cd5e0SMatthew G. Knepley maxVolumes[c - cStart] = vol * refRatio; 410d1cd5e0SMatthew G. Knepley } 420d1cd5e0SMatthew G. Knepley } 430d1cd5e0SMatthew G. Knepley PetscFunctionReturn(0); 440d1cd5e0SMatthew G. Knepley } 450d1cd5e0SMatthew G. Knepley 46*9371c9d4SSatish Balay static PetscErrorCode DMPlexLabelToMetricConstraint(DM dm, DMLabel adaptLabel, PetscInt cStart, PetscInt cEnd, PetscInt vStart, PetscInt vEnd, PetscReal refRatio, Vec *metricVec) { 470d1cd5e0SMatthew G. Knepley DM udm, coordDM; 480d1cd5e0SMatthew G. Knepley PetscSection coordSection; 490d1cd5e0SMatthew G. Knepley Vec coordinates, mb, mx; 500d1cd5e0SMatthew G. Knepley Mat A; 510d1cd5e0SMatthew G. Knepley PetscScalar *metric, *eqns; 520d1cd5e0SMatthew G. Knepley const PetscReal coarseRatio = refRatio == PETSC_DEFAULT ? PetscSqr(0.5) : 1 / refRatio; 530d1cd5e0SMatthew G. Knepley PetscInt dim, Nv, Neq, c, v; 540d1cd5e0SMatthew G. Knepley 550d1cd5e0SMatthew G. Knepley PetscFunctionBegin; 569566063dSJacob Faibussowitsch PetscCall(DMPlexUninterpolate(dm, &udm)); 579566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 589566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &coordDM)); 599566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(coordDM, &coordSection)); 609566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 610d1cd5e0SMatthew G. Knepley Nv = vEnd - vStart; 629566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, Nv * PetscSqr(dim), metricVec)); 639566063dSJacob Faibussowitsch PetscCall(VecGetArray(*metricVec, &metric)); 640d1cd5e0SMatthew G. Knepley Neq = (dim * (dim + 1)) / 2; 659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(PetscSqr(Neq), &eqns)); 669566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, Neq, Neq, eqns, &A)); 679566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &mx, &mb)); 689566063dSJacob Faibussowitsch PetscCall(VecSet(mb, 1.0)); 690d1cd5e0SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 700d1cd5e0SMatthew G. Knepley const PetscScalar *sol; 710d1cd5e0SMatthew G. Knepley PetscScalar *cellCoords = NULL; 720d1cd5e0SMatthew G. Knepley PetscReal e[3], vol; 730d1cd5e0SMatthew G. Knepley const PetscInt *cone; 740d1cd5e0SMatthew G. Knepley PetscInt coneSize, cl, i, j, d, r; 750d1cd5e0SMatthew G. Knepley 769566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, coordSection, coordinates, c, NULL, &cellCoords)); 770d1cd5e0SMatthew G. Knepley /* Only works for simplices */ 780d1cd5e0SMatthew G. Knepley for (i = 0, r = 0; i < dim + 1; ++i) { 790d1cd5e0SMatthew G. Knepley for (j = 0; j < i; ++j, ++r) { 800d1cd5e0SMatthew G. Knepley for (d = 0; d < dim; ++d) e[d] = PetscRealPart(cellCoords[i * dim + d] - cellCoords[j * dim + d]); 810d1cd5e0SMatthew G. Knepley /* FORTRAN ORDERING */ 820d1cd5e0SMatthew G. Knepley switch (dim) { 830d1cd5e0SMatthew G. Knepley case 2: 840d1cd5e0SMatthew G. Knepley eqns[0 * Neq + r] = PetscSqr(e[0]); 850d1cd5e0SMatthew G. Knepley eqns[1 * Neq + r] = 2.0 * e[0] * e[1]; 860d1cd5e0SMatthew G. Knepley eqns[2 * Neq + r] = PetscSqr(e[1]); 870d1cd5e0SMatthew G. Knepley break; 880d1cd5e0SMatthew G. Knepley case 3: 890d1cd5e0SMatthew G. Knepley eqns[0 * Neq + r] = PetscSqr(e[0]); 900d1cd5e0SMatthew G. Knepley eqns[1 * Neq + r] = 2.0 * e[0] * e[1]; 910d1cd5e0SMatthew G. Knepley eqns[2 * Neq + r] = 2.0 * e[0] * e[2]; 920d1cd5e0SMatthew G. Knepley eqns[3 * Neq + r] = PetscSqr(e[1]); 930d1cd5e0SMatthew G. Knepley eqns[4 * Neq + r] = 2.0 * e[1] * e[2]; 940d1cd5e0SMatthew G. Knepley eqns[5 * Neq + r] = PetscSqr(e[2]); 950d1cd5e0SMatthew G. Knepley break; 960d1cd5e0SMatthew G. Knepley } 970d1cd5e0SMatthew G. Knepley } 980d1cd5e0SMatthew G. Knepley } 999566063dSJacob Faibussowitsch PetscCall(MatSetUnfactored(A)); 1009566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, coordSection, coordinates, c, NULL, &cellCoords)); 1019566063dSJacob Faibussowitsch PetscCall(MatLUFactor(A, NULL, NULL, NULL)); 1029566063dSJacob Faibussowitsch PetscCall(MatSolve(A, mb, mx)); 1039566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(mx, &sol)); 1049566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dm, c, &vol, NULL, NULL)); 1059566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(udm, c, &cone)); 1069566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(udm, c, &coneSize)); 1070d1cd5e0SMatthew G. Knepley for (cl = 0; cl < coneSize; ++cl) { 1080d1cd5e0SMatthew G. Knepley const PetscInt v = cone[cl] - vStart; 1090d1cd5e0SMatthew G. Knepley 1100d1cd5e0SMatthew G. Knepley if (dim == 2) { 1110d1cd5e0SMatthew G. Knepley metric[v * 4 + 0] += vol * coarseRatio * sol[0]; 1120d1cd5e0SMatthew G. Knepley metric[v * 4 + 1] += vol * coarseRatio * sol[1]; 1130d1cd5e0SMatthew G. Knepley metric[v * 4 + 2] += vol * coarseRatio * sol[1]; 1140d1cd5e0SMatthew G. Knepley metric[v * 4 + 3] += vol * coarseRatio * sol[2]; 1150d1cd5e0SMatthew G. Knepley } else { 1160d1cd5e0SMatthew G. Knepley metric[v * 9 + 0] += vol * coarseRatio * sol[0]; 1170d1cd5e0SMatthew G. Knepley metric[v * 9 + 1] += vol * coarseRatio * sol[1]; 1180d1cd5e0SMatthew G. Knepley metric[v * 9 + 3] += vol * coarseRatio * sol[1]; 1190d1cd5e0SMatthew G. Knepley metric[v * 9 + 2] += vol * coarseRatio * sol[2]; 1200d1cd5e0SMatthew G. Knepley metric[v * 9 + 6] += vol * coarseRatio * sol[2]; 1210d1cd5e0SMatthew G. Knepley metric[v * 9 + 4] += vol * coarseRatio * sol[3]; 1220d1cd5e0SMatthew G. Knepley metric[v * 9 + 5] += vol * coarseRatio * sol[4]; 1230d1cd5e0SMatthew G. Knepley metric[v * 9 + 7] += vol * coarseRatio * sol[4]; 1240d1cd5e0SMatthew G. Knepley metric[v * 9 + 8] += vol * coarseRatio * sol[5]; 1250d1cd5e0SMatthew G. Knepley } 1260d1cd5e0SMatthew G. Knepley } 1279566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(mx, &sol)); 1280d1cd5e0SMatthew G. Knepley } 1290d1cd5e0SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 1300d1cd5e0SMatthew G. Knepley const PetscInt *support; 1310d1cd5e0SMatthew G. Knepley PetscInt supportSize, s; 1320d1cd5e0SMatthew G. Knepley PetscReal vol, totVol = 0.0; 1330d1cd5e0SMatthew G. Knepley 1349566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(udm, v + vStart, &support)); 1359566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(udm, v + vStart, &supportSize)); 136*9371c9d4SSatish Balay for (s = 0; s < supportSize; ++s) { 137*9371c9d4SSatish Balay PetscCall(DMPlexComputeCellGeometryFVM(dm, support[s], &vol, NULL, NULL)); 138*9371c9d4SSatish Balay totVol += vol; 139*9371c9d4SSatish Balay } 1400d1cd5e0SMatthew G. Knepley for (s = 0; s < PetscSqr(dim); ++s) metric[v * PetscSqr(dim) + s] /= totVol; 1410d1cd5e0SMatthew G. Knepley } 1429566063dSJacob Faibussowitsch PetscCall(PetscFree(eqns)); 1439566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(*metricVec, &metric)); 1449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mx)); 1459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mb)); 1469566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1479566063dSJacob Faibussowitsch PetscCall(DMDestroy(&udm)); 1480d1cd5e0SMatthew G. Knepley PetscFunctionReturn(0); 1490d1cd5e0SMatthew G. Knepley } 1500d1cd5e0SMatthew G. Knepley 1513a074057SBarry Smith /* 1523a074057SBarry Smith Contains the list of registered DMPlexGenerators routines 1533a074057SBarry Smith */ 154*9371c9d4SSatish Balay PetscErrorCode DMPlexRefine_Internal(DM dm, PETSC_UNUSED Vec metric, DMLabel adaptLabel, PETSC_UNUSED DMLabel rgLabel, DM *dmRefined) { 155c0517cd5SMatthew G. Knepley DMGeneratorFunctionList fl; 156800b850cSBarry Smith PetscErrorCode (*refine)(DM, PetscReal *, DM *); 1579fe9e680SJoe Wallwork PetscErrorCode (*adapt)(DM, Vec, DMLabel, DMLabel, DM *); 1583f2a96e3SMatthew G. Knepley PetscErrorCode (*refinementFunc)(const PetscReal[], PetscReal *); 1593f2a96e3SMatthew G. Knepley char genname[PETSC_MAX_PATH_LEN], *name = NULL; 1603f2a96e3SMatthew G. Knepley PetscReal refinementLimit; 1613a074057SBarry Smith PetscReal *maxVolumes; 1623f2a96e3SMatthew G. Knepley PetscInt dim, cStart, cEnd, c; 1633f2a96e3SMatthew G. Knepley PetscBool flg, flg2, localized; 1640d1cd5e0SMatthew G. Knepley 1650d1cd5e0SMatthew G. Knepley PetscFunctionBegin; 1669566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocalized(dm, &localized)); 1679566063dSJacob Faibussowitsch PetscCall(DMPlexGetRefinementLimit(dm, &refinementLimit)); 1689566063dSJacob Faibussowitsch PetscCall(DMPlexGetRefinementFunction(dm, &refinementFunc)); 1690d1cd5e0SMatthew G. Knepley if (refinementLimit == 0.0 && !refinementFunc && !adaptLabel) PetscFunctionReturn(0); 1709566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1719566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 1729566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_adaptor", genname, sizeof(genname), &flg)); 1730d1cd5e0SMatthew G. Knepley if (flg) name = genname; 1743f2a96e3SMatthew G. Knepley else { 1759566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_generator", genname, sizeof(genname), &flg2)); 1763f2a96e3SMatthew G. Knepley if (flg2) name = genname; 1773f2a96e3SMatthew G. Knepley } 1787d99616aSKarl Rupp 179c0517cd5SMatthew G. Knepley fl = DMGenerateList; 1803a074057SBarry Smith if (name) { 1813a074057SBarry Smith while (fl) { 1829566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(fl->name, name, &flg)); 1833a074057SBarry Smith if (flg) { 1843a074057SBarry Smith refine = fl->refine; 185c0517cd5SMatthew G. Knepley adapt = fl->adapt; 1863a074057SBarry Smith goto gotit; 1873a074057SBarry Smith } 1883a074057SBarry Smith fl = fl->next; 1893a074057SBarry Smith } 19098921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Grid refiner %s not registered", name); 1913a074057SBarry Smith } else { 1923a074057SBarry Smith while (fl) { 193a12d352dSMatthew G. Knepley if (fl->dim < 0 || dim - 1 == fl->dim) { 1943a074057SBarry Smith refine = fl->refine; 195c0517cd5SMatthew G. Knepley adapt = fl->adapt; 1963a074057SBarry Smith goto gotit; 1973a074057SBarry Smith } 1983a074057SBarry Smith fl = fl->next; 1993a074057SBarry Smith } 20063a3b9bcSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No grid refiner of dimension %" PetscInt_FMT " registered", dim); 2013a074057SBarry Smith } 2023a074057SBarry Smith 2033f2a96e3SMatthew G. Knepley gotit: 2043f2a96e3SMatthew G. Knepley switch (dim) { 205a12d352dSMatthew G. Knepley case 1: 2063a074057SBarry Smith case 2: 2070d1cd5e0SMatthew G. Knepley case 3: 2083f2a96e3SMatthew G. Knepley if (adapt) { 2099566063dSJacob Faibussowitsch PetscCall((*adapt)(dm, NULL, adaptLabel, NULL, dmRefined)); 2103f2a96e3SMatthew G. Knepley } else { 2119566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(cEnd - cStart, &maxVolumes)); 2120d1cd5e0SMatthew G. Knepley if (adaptLabel) { 2139566063dSJacob Faibussowitsch PetscCall(DMPlexLabelToVolumeConstraint(dm, adaptLabel, cStart, cEnd, PETSC_DEFAULT, maxVolumes)); 2140d1cd5e0SMatthew G. Knepley } else if (refinementFunc) { 2150d1cd5e0SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 2160d1cd5e0SMatthew G. Knepley PetscReal vol, centroid[3]; 2170d1cd5e0SMatthew G. Knepley 2189566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL)); 2199566063dSJacob Faibussowitsch PetscCall((*refinementFunc)(centroid, &maxVolumes[c - cStart])); 2200d1cd5e0SMatthew G. Knepley } 2210d1cd5e0SMatthew G. Knepley } else { 2220d1cd5e0SMatthew G. Knepley for (c = 0; c < cEnd - cStart; ++c) maxVolumes[c] = refinementLimit; 2230d1cd5e0SMatthew G. Knepley } 2249566063dSJacob Faibussowitsch PetscCall((*refine)(dm, maxVolumes, dmRefined)); 2259566063dSJacob Faibussowitsch PetscCall(PetscFree(maxVolumes)); 2263f2a96e3SMatthew G. Knepley } 2270d1cd5e0SMatthew G. Knepley break; 22863a3b9bcSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Mesh refinement in dimension %" PetscInt_FMT " is not supported.", dim); 2290d1cd5e0SMatthew G. Knepley } 2309566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, *dmRefined)); 2315de52c6dSVaclav Hapla PetscCall(DMPlexCopy_Internal(dm, PETSC_TRUE, PETSC_TRUE, *dmRefined)); 2329566063dSJacob Faibussowitsch if (localized) PetscCall(DMLocalizeCoordinates(*dmRefined)); 2330d1cd5e0SMatthew G. Knepley PetscFunctionReturn(0); 2340d1cd5e0SMatthew G. Knepley } 2350d1cd5e0SMatthew G. Knepley 236*9371c9d4SSatish Balay PetscErrorCode DMPlexCoarsen_Internal(DM dm, PETSC_UNUSED Vec metric, DMLabel adaptLabel, PETSC_UNUSED DMLabel rgLabel, DM *dmCoarsened) { 2370d1cd5e0SMatthew G. Knepley Vec metricVec; 2380d1cd5e0SMatthew G. Knepley PetscInt cStart, cEnd, vStart, vEnd; 2390d1cd5e0SMatthew G. Knepley DMLabel bdLabel = NULL; 2409fe9e680SJoe Wallwork char bdLabelName[PETSC_MAX_PATH_LEN], rgLabelName[PETSC_MAX_PATH_LEN]; 2410d1cd5e0SMatthew G. Knepley PetscBool localized, flg; 2420d1cd5e0SMatthew G. Knepley 2430d1cd5e0SMatthew G. Knepley PetscFunctionBegin; 2449566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocalized(dm, &localized)); 2459566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 2469566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 2479566063dSJacob Faibussowitsch PetscCall(DMPlexLabelToMetricConstraint(dm, adaptLabel, cStart, cEnd, vStart, vEnd, PETSC_DEFAULT, &metricVec)); 2489566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, dm->hdr.prefix, "-dm_plex_coarsen_bd_label", bdLabelName, sizeof(bdLabelName), &flg)); 2499566063dSJacob Faibussowitsch if (flg) PetscCall(DMGetLabel(dm, bdLabelName, &bdLabel)); 2509566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, dm->hdr.prefix, "-dm_plex_coarsen_rg_label", rgLabelName, sizeof(rgLabelName), &flg)); 2519566063dSJacob Faibussowitsch if (flg) PetscCall(DMGetLabel(dm, rgLabelName, &rgLabel)); 2529566063dSJacob Faibussowitsch PetscCall(DMAdaptMetric(dm, metricVec, bdLabel, rgLabel, dmCoarsened)); 2539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&metricVec)); 2549566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, *dmCoarsened)); 2555de52c6dSVaclav Hapla PetscCall(DMPlexCopy_Internal(dm, PETSC_TRUE, PETSC_TRUE, *dmCoarsened)); 2569566063dSJacob Faibussowitsch if (localized) PetscCall(DMLocalizeCoordinates(*dmCoarsened)); 2570d1cd5e0SMatthew G. Knepley PetscFunctionReturn(0); 2580d1cd5e0SMatthew G. Knepley } 2590d1cd5e0SMatthew G. Knepley 260*9371c9d4SSatish Balay PetscErrorCode DMAdaptLabel_Plex(DM dm, PETSC_UNUSED Vec metric, DMLabel adaptLabel, PETSC_UNUSED DMLabel rgLabel, DM *dmAdapted) { 2610d1cd5e0SMatthew G. Knepley IS flagIS; 2620d1cd5e0SMatthew G. Knepley const PetscInt *flags; 2630d1cd5e0SMatthew G. Knepley PetscInt defFlag, minFlag, maxFlag, numFlags, f; 2640d1cd5e0SMatthew G. Knepley 2650d1cd5e0SMatthew G. Knepley PetscFunctionBegin; 2669566063dSJacob Faibussowitsch PetscCall(DMLabelGetDefaultValue(adaptLabel, &defFlag)); 2670d1cd5e0SMatthew G. Knepley minFlag = defFlag; 2680d1cd5e0SMatthew G. Knepley maxFlag = defFlag; 2699566063dSJacob Faibussowitsch PetscCall(DMLabelGetValueIS(adaptLabel, &flagIS)); 2709566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(flagIS, &numFlags)); 2719566063dSJacob Faibussowitsch PetscCall(ISGetIndices(flagIS, &flags)); 2720d1cd5e0SMatthew G. Knepley for (f = 0; f < numFlags; ++f) { 2730d1cd5e0SMatthew G. Knepley const PetscInt flag = flags[f]; 2740d1cd5e0SMatthew G. Knepley 2750d1cd5e0SMatthew G. Knepley minFlag = PetscMin(minFlag, flag); 2760d1cd5e0SMatthew G. Knepley maxFlag = PetscMax(maxFlag, flag); 2770d1cd5e0SMatthew G. Knepley } 2789566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(flagIS, &flags)); 2799566063dSJacob Faibussowitsch PetscCall(ISDestroy(&flagIS)); 2800d1cd5e0SMatthew G. Knepley { 2810d1cd5e0SMatthew G. Knepley PetscInt minMaxFlag[2], minMaxFlagGlobal[2]; 2820d1cd5e0SMatthew G. Knepley 2830d1cd5e0SMatthew G. Knepley minMaxFlag[0] = minFlag; 2840d1cd5e0SMatthew G. Knepley minMaxFlag[1] = -maxFlag; 2859566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(minMaxFlag, minMaxFlagGlobal, 2, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm))); 2860d1cd5e0SMatthew G. Knepley minFlag = minMaxFlagGlobal[0]; 2870d1cd5e0SMatthew G. Knepley maxFlag = -minMaxFlagGlobal[1]; 2880d1cd5e0SMatthew G. Knepley } 2890d1cd5e0SMatthew G. Knepley if (minFlag == maxFlag) { 2900d1cd5e0SMatthew G. Knepley switch (minFlag) { 291*9371c9d4SSatish Balay case DM_ADAPT_DETERMINE: *dmAdapted = NULL; break; 2920d1cd5e0SMatthew G. Knepley case DM_ADAPT_REFINE: 2939566063dSJacob Faibussowitsch PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE)); 294*9371c9d4SSatish Balay PetscCall(DMRefine(dm, MPI_COMM_NULL, dmAdapted)); 295*9371c9d4SSatish Balay break; 296*9371c9d4SSatish Balay case DM_ADAPT_COARSEN: PetscCall(DMCoarsen(dm, MPI_COMM_NULL, dmAdapted)); break; 297*9371c9d4SSatish Balay default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "DMPlex does not support refinement flag %" PetscInt_FMT, minFlag); 2980d1cd5e0SMatthew G. Knepley } 2990d1cd5e0SMatthew G. Knepley } else { 3009566063dSJacob Faibussowitsch PetscCall(DMPlexSetRefinementUniform(dm, PETSC_FALSE)); 3019566063dSJacob Faibussowitsch PetscCall(DMPlexRefine_Internal(dm, NULL, adaptLabel, NULL, dmAdapted)); 3020d1cd5e0SMatthew G. Knepley } 3030d1cd5e0SMatthew G. Knepley PetscFunctionReturn(0); 3040d1cd5e0SMatthew G. Knepley } 305