10bbe5d31Sbarral #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 20bbe5d31Sbarral 3d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexLabelToVolumeConstraint(DM dm, DMLabel adaptLabel, PetscInt cStart, PetscInt cEnd, PetscReal refRatio, PetscReal maxVolumes[]) 4d71ae5a4SJacob Faibussowitsch { 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) { 27d71ae5a4SJacob Faibussowitsch case DM_ADAPT_REFINE: 28d71ae5a4SJacob Faibussowitsch anyRefine = PETSC_TRUE; 29d71ae5a4SJacob Faibussowitsch break; 30d71ae5a4SJacob Faibussowitsch case DM_ADAPT_COARSEN: 31d71ae5a4SJacob Faibussowitsch anyCoarsen = PETSC_TRUE; 32d71ae5a4SJacob Faibussowitsch break; 33d71ae5a4SJacob Faibussowitsch case DM_ADAPT_KEEP: 34d71ae5a4SJacob Faibussowitsch anyKeep = PETSC_TRUE; 35d71ae5a4SJacob Faibussowitsch break; 36d71ae5a4SJacob Faibussowitsch case DM_ADAPT_DETERMINE: 37d71ae5a4SJacob Faibussowitsch break; 38d71ae5a4SJacob Faibussowitsch default: 39d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "DMPlex does not support refinement flag %" PetscInt_FMT, refFlag); 400d1cd5e0SMatthew G. Knepley } 410d1cd5e0SMatthew G. Knepley if (anyRefine) break; 420d1cd5e0SMatthew G. Knepley } 439566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 440d1cd5e0SMatthew G. Knepley if (anyRefine) { 450d1cd5e0SMatthew G. Knepley maxVolumes[c - cStart] = vol / refRatio; 460d1cd5e0SMatthew G. Knepley } else if (anyKeep) { 470d1cd5e0SMatthew G. Knepley maxVolumes[c - cStart] = vol; 480d1cd5e0SMatthew G. Knepley } else if (anyCoarsen) { 490d1cd5e0SMatthew G. Knepley maxVolumes[c - cStart] = vol * refRatio; 500d1cd5e0SMatthew G. Knepley } 510d1cd5e0SMatthew G. Knepley } 52*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 530d1cd5e0SMatthew G. Knepley } 540d1cd5e0SMatthew G. Knepley 55d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexLabelToMetricConstraint(DM dm, DMLabel adaptLabel, PetscInt cStart, PetscInt cEnd, PetscInt vStart, PetscInt vEnd, PetscReal refRatio, Vec *metricVec) 56d71ae5a4SJacob Faibussowitsch { 570d1cd5e0SMatthew G. Knepley DM udm, coordDM; 580d1cd5e0SMatthew G. Knepley PetscSection coordSection; 590d1cd5e0SMatthew G. Knepley Vec coordinates, mb, mx; 600d1cd5e0SMatthew G. Knepley Mat A; 610d1cd5e0SMatthew G. Knepley PetscScalar *metric, *eqns; 620d1cd5e0SMatthew G. Knepley const PetscReal coarseRatio = refRatio == PETSC_DEFAULT ? PetscSqr(0.5) : 1 / refRatio; 630d1cd5e0SMatthew G. Knepley PetscInt dim, Nv, Neq, c, v; 640d1cd5e0SMatthew G. Knepley 650d1cd5e0SMatthew G. Knepley PetscFunctionBegin; 669566063dSJacob Faibussowitsch PetscCall(DMPlexUninterpolate(dm, &udm)); 679566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 689566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &coordDM)); 699566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(coordDM, &coordSection)); 709566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 710d1cd5e0SMatthew G. Knepley Nv = vEnd - vStart; 729566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, Nv * PetscSqr(dim), metricVec)); 739566063dSJacob Faibussowitsch PetscCall(VecGetArray(*metricVec, &metric)); 740d1cd5e0SMatthew G. Knepley Neq = (dim * (dim + 1)) / 2; 759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(PetscSqr(Neq), &eqns)); 769566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, Neq, Neq, eqns, &A)); 779566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &mx, &mb)); 789566063dSJacob Faibussowitsch PetscCall(VecSet(mb, 1.0)); 790d1cd5e0SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 800d1cd5e0SMatthew G. Knepley const PetscScalar *sol; 810d1cd5e0SMatthew G. Knepley PetscScalar *cellCoords = NULL; 820d1cd5e0SMatthew G. Knepley PetscReal e[3], vol; 830d1cd5e0SMatthew G. Knepley const PetscInt *cone; 840d1cd5e0SMatthew G. Knepley PetscInt coneSize, cl, i, j, d, r; 850d1cd5e0SMatthew G. Knepley 869566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, coordSection, coordinates, c, NULL, &cellCoords)); 870d1cd5e0SMatthew G. Knepley /* Only works for simplices */ 880d1cd5e0SMatthew G. Knepley for (i = 0, r = 0; i < dim + 1; ++i) { 890d1cd5e0SMatthew G. Knepley for (j = 0; j < i; ++j, ++r) { 900d1cd5e0SMatthew G. Knepley for (d = 0; d < dim; ++d) e[d] = PetscRealPart(cellCoords[i * dim + d] - cellCoords[j * dim + d]); 910d1cd5e0SMatthew G. Knepley /* FORTRAN ORDERING */ 920d1cd5e0SMatthew G. Knepley switch (dim) { 930d1cd5e0SMatthew G. Knepley case 2: 940d1cd5e0SMatthew G. Knepley eqns[0 * Neq + r] = PetscSqr(e[0]); 950d1cd5e0SMatthew G. Knepley eqns[1 * Neq + r] = 2.0 * e[0] * e[1]; 960d1cd5e0SMatthew G. Knepley eqns[2 * Neq + r] = PetscSqr(e[1]); 970d1cd5e0SMatthew G. Knepley break; 980d1cd5e0SMatthew G. Knepley case 3: 990d1cd5e0SMatthew G. Knepley eqns[0 * Neq + r] = PetscSqr(e[0]); 1000d1cd5e0SMatthew G. Knepley eqns[1 * Neq + r] = 2.0 * e[0] * e[1]; 1010d1cd5e0SMatthew G. Knepley eqns[2 * Neq + r] = 2.0 * e[0] * e[2]; 1020d1cd5e0SMatthew G. Knepley eqns[3 * Neq + r] = PetscSqr(e[1]); 1030d1cd5e0SMatthew G. Knepley eqns[4 * Neq + r] = 2.0 * e[1] * e[2]; 1040d1cd5e0SMatthew G. Knepley eqns[5 * Neq + r] = PetscSqr(e[2]); 1050d1cd5e0SMatthew G. Knepley break; 1060d1cd5e0SMatthew G. Knepley } 1070d1cd5e0SMatthew G. Knepley } 1080d1cd5e0SMatthew G. Knepley } 1099566063dSJacob Faibussowitsch PetscCall(MatSetUnfactored(A)); 1109566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, coordSection, coordinates, c, NULL, &cellCoords)); 1119566063dSJacob Faibussowitsch PetscCall(MatLUFactor(A, NULL, NULL, NULL)); 1129566063dSJacob Faibussowitsch PetscCall(MatSolve(A, mb, mx)); 1139566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(mx, &sol)); 1149566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dm, c, &vol, NULL, NULL)); 1159566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(udm, c, &cone)); 1169566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(udm, c, &coneSize)); 1170d1cd5e0SMatthew G. Knepley for (cl = 0; cl < coneSize; ++cl) { 1180d1cd5e0SMatthew G. Knepley const PetscInt v = cone[cl] - vStart; 1190d1cd5e0SMatthew G. Knepley 1200d1cd5e0SMatthew G. Knepley if (dim == 2) { 1210d1cd5e0SMatthew G. Knepley metric[v * 4 + 0] += vol * coarseRatio * sol[0]; 1220d1cd5e0SMatthew G. Knepley metric[v * 4 + 1] += vol * coarseRatio * sol[1]; 1230d1cd5e0SMatthew G. Knepley metric[v * 4 + 2] += vol * coarseRatio * sol[1]; 1240d1cd5e0SMatthew G. Knepley metric[v * 4 + 3] += vol * coarseRatio * sol[2]; 1250d1cd5e0SMatthew G. Knepley } else { 1260d1cd5e0SMatthew G. Knepley metric[v * 9 + 0] += vol * coarseRatio * sol[0]; 1270d1cd5e0SMatthew G. Knepley metric[v * 9 + 1] += vol * coarseRatio * sol[1]; 1280d1cd5e0SMatthew G. Knepley metric[v * 9 + 3] += vol * coarseRatio * sol[1]; 1290d1cd5e0SMatthew G. Knepley metric[v * 9 + 2] += vol * coarseRatio * sol[2]; 1300d1cd5e0SMatthew G. Knepley metric[v * 9 + 6] += vol * coarseRatio * sol[2]; 1310d1cd5e0SMatthew G. Knepley metric[v * 9 + 4] += vol * coarseRatio * sol[3]; 1320d1cd5e0SMatthew G. Knepley metric[v * 9 + 5] += vol * coarseRatio * sol[4]; 1330d1cd5e0SMatthew G. Knepley metric[v * 9 + 7] += vol * coarseRatio * sol[4]; 1340d1cd5e0SMatthew G. Knepley metric[v * 9 + 8] += vol * coarseRatio * sol[5]; 1350d1cd5e0SMatthew G. Knepley } 1360d1cd5e0SMatthew G. Knepley } 1379566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(mx, &sol)); 1380d1cd5e0SMatthew G. Knepley } 1390d1cd5e0SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 1400d1cd5e0SMatthew G. Knepley const PetscInt *support; 1410d1cd5e0SMatthew G. Knepley PetscInt supportSize, s; 1420d1cd5e0SMatthew G. Knepley PetscReal vol, totVol = 0.0; 1430d1cd5e0SMatthew G. Knepley 1449566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(udm, v + vStart, &support)); 1459566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(udm, v + vStart, &supportSize)); 1469371c9d4SSatish Balay for (s = 0; s < supportSize; ++s) { 1479371c9d4SSatish Balay PetscCall(DMPlexComputeCellGeometryFVM(dm, support[s], &vol, NULL, NULL)); 1489371c9d4SSatish Balay totVol += vol; 1499371c9d4SSatish Balay } 1500d1cd5e0SMatthew G. Knepley for (s = 0; s < PetscSqr(dim); ++s) metric[v * PetscSqr(dim) + s] /= totVol; 1510d1cd5e0SMatthew G. Knepley } 1529566063dSJacob Faibussowitsch PetscCall(PetscFree(eqns)); 1539566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(*metricVec, &metric)); 1549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mx)); 1559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mb)); 1569566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1579566063dSJacob Faibussowitsch PetscCall(DMDestroy(&udm)); 158*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1590d1cd5e0SMatthew G. Knepley } 1600d1cd5e0SMatthew G. Knepley 1613a074057SBarry Smith /* 1623a074057SBarry Smith Contains the list of registered DMPlexGenerators routines 1633a074057SBarry Smith */ 164d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexRefine_Internal(DM dm, PETSC_UNUSED Vec metric, DMLabel adaptLabel, PETSC_UNUSED DMLabel rgLabel, DM *dmRefined) 165d71ae5a4SJacob Faibussowitsch { 166c0517cd5SMatthew G. Knepley DMGeneratorFunctionList fl; 167800b850cSBarry Smith PetscErrorCode (*refine)(DM, PetscReal *, DM *); 1689fe9e680SJoe Wallwork PetscErrorCode (*adapt)(DM, Vec, DMLabel, DMLabel, DM *); 1693f2a96e3SMatthew G. Knepley PetscErrorCode (*refinementFunc)(const PetscReal[], PetscReal *); 1703f2a96e3SMatthew G. Knepley char genname[PETSC_MAX_PATH_LEN], *name = NULL; 1713f2a96e3SMatthew G. Knepley PetscReal refinementLimit; 1723a074057SBarry Smith PetscReal *maxVolumes; 1733f2a96e3SMatthew G. Knepley PetscInt dim, cStart, cEnd, c; 1743f2a96e3SMatthew G. Knepley PetscBool flg, flg2, localized; 1750d1cd5e0SMatthew G. Knepley 1760d1cd5e0SMatthew G. Knepley PetscFunctionBegin; 1779566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocalized(dm, &localized)); 1789566063dSJacob Faibussowitsch PetscCall(DMPlexGetRefinementLimit(dm, &refinementLimit)); 1799566063dSJacob Faibussowitsch PetscCall(DMPlexGetRefinementFunction(dm, &refinementFunc)); 180*3ba16761SJacob Faibussowitsch if (refinementLimit == 0.0 && !refinementFunc && !adaptLabel) PetscFunctionReturn(PETSC_SUCCESS); 1819566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1829566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 1839566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_adaptor", genname, sizeof(genname), &flg)); 1840d1cd5e0SMatthew G. Knepley if (flg) name = genname; 1853f2a96e3SMatthew G. Knepley else { 1869566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_generator", genname, sizeof(genname), &flg2)); 1873f2a96e3SMatthew G. Knepley if (flg2) name = genname; 1883f2a96e3SMatthew G. Knepley } 1897d99616aSKarl Rupp 190c0517cd5SMatthew G. Knepley fl = DMGenerateList; 1913a074057SBarry Smith if (name) { 1923a074057SBarry Smith while (fl) { 1939566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(fl->name, name, &flg)); 1943a074057SBarry Smith if (flg) { 1953a074057SBarry Smith refine = fl->refine; 196c0517cd5SMatthew G. Knepley adapt = fl->adapt; 1973a074057SBarry Smith goto gotit; 1983a074057SBarry Smith } 1993a074057SBarry Smith fl = fl->next; 2003a074057SBarry Smith } 20198921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Grid refiner %s not registered", name); 2023a074057SBarry Smith } else { 2033a074057SBarry Smith while (fl) { 204a12d352dSMatthew G. Knepley if (fl->dim < 0 || dim - 1 == fl->dim) { 2053a074057SBarry Smith refine = fl->refine; 206c0517cd5SMatthew G. Knepley adapt = fl->adapt; 2073a074057SBarry Smith goto gotit; 2083a074057SBarry Smith } 2093a074057SBarry Smith fl = fl->next; 2103a074057SBarry Smith } 21163a3b9bcSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No grid refiner of dimension %" PetscInt_FMT " registered", dim); 2123a074057SBarry Smith } 2133a074057SBarry Smith 2143f2a96e3SMatthew G. Knepley gotit: 2153f2a96e3SMatthew G. Knepley switch (dim) { 216a12d352dSMatthew G. Knepley case 1: 2173a074057SBarry Smith case 2: 2180d1cd5e0SMatthew G. Knepley case 3: 2193f2a96e3SMatthew G. Knepley if (adapt) { 2209566063dSJacob Faibussowitsch PetscCall((*adapt)(dm, NULL, adaptLabel, NULL, dmRefined)); 2213f2a96e3SMatthew G. Knepley } else { 2229566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(cEnd - cStart, &maxVolumes)); 2230d1cd5e0SMatthew G. Knepley if (adaptLabel) { 2249566063dSJacob Faibussowitsch PetscCall(DMPlexLabelToVolumeConstraint(dm, adaptLabel, cStart, cEnd, PETSC_DEFAULT, maxVolumes)); 2250d1cd5e0SMatthew G. Knepley } else if (refinementFunc) { 2260d1cd5e0SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 2270d1cd5e0SMatthew G. Knepley PetscReal vol, centroid[3]; 2280d1cd5e0SMatthew G. Knepley 2299566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL)); 2309566063dSJacob Faibussowitsch PetscCall((*refinementFunc)(centroid, &maxVolumes[c - cStart])); 2310d1cd5e0SMatthew G. Knepley } 2320d1cd5e0SMatthew G. Knepley } else { 2330d1cd5e0SMatthew G. Knepley for (c = 0; c < cEnd - cStart; ++c) maxVolumes[c] = refinementLimit; 2340d1cd5e0SMatthew G. Knepley } 2359566063dSJacob Faibussowitsch PetscCall((*refine)(dm, maxVolumes, dmRefined)); 2369566063dSJacob Faibussowitsch PetscCall(PetscFree(maxVolumes)); 2373f2a96e3SMatthew G. Knepley } 2380d1cd5e0SMatthew G. Knepley break; 239d71ae5a4SJacob Faibussowitsch default: 240d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Mesh refinement in dimension %" PetscInt_FMT " is not supported.", dim); 2410d1cd5e0SMatthew G. Knepley } 2429566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, *dmRefined)); 2435de52c6dSVaclav Hapla PetscCall(DMPlexCopy_Internal(dm, PETSC_TRUE, PETSC_TRUE, *dmRefined)); 2449566063dSJacob Faibussowitsch if (localized) PetscCall(DMLocalizeCoordinates(*dmRefined)); 245*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2460d1cd5e0SMatthew G. Knepley } 2470d1cd5e0SMatthew G. Knepley 248d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCoarsen_Internal(DM dm, PETSC_UNUSED Vec metric, DMLabel adaptLabel, PETSC_UNUSED DMLabel rgLabel, DM *dmCoarsened) 249d71ae5a4SJacob Faibussowitsch { 2500d1cd5e0SMatthew G. Knepley Vec metricVec; 2510d1cd5e0SMatthew G. Knepley PetscInt cStart, cEnd, vStart, vEnd; 2520d1cd5e0SMatthew G. Knepley DMLabel bdLabel = NULL; 2539fe9e680SJoe Wallwork char bdLabelName[PETSC_MAX_PATH_LEN], rgLabelName[PETSC_MAX_PATH_LEN]; 2540d1cd5e0SMatthew G. Knepley PetscBool localized, flg; 2550d1cd5e0SMatthew G. Knepley 2560d1cd5e0SMatthew G. Knepley PetscFunctionBegin; 2579566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocalized(dm, &localized)); 2589566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 2599566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 2609566063dSJacob Faibussowitsch PetscCall(DMPlexLabelToMetricConstraint(dm, adaptLabel, cStart, cEnd, vStart, vEnd, PETSC_DEFAULT, &metricVec)); 2619566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, dm->hdr.prefix, "-dm_plex_coarsen_bd_label", bdLabelName, sizeof(bdLabelName), &flg)); 2629566063dSJacob Faibussowitsch if (flg) PetscCall(DMGetLabel(dm, bdLabelName, &bdLabel)); 2639566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, dm->hdr.prefix, "-dm_plex_coarsen_rg_label", rgLabelName, sizeof(rgLabelName), &flg)); 2649566063dSJacob Faibussowitsch if (flg) PetscCall(DMGetLabel(dm, rgLabelName, &rgLabel)); 2659566063dSJacob Faibussowitsch PetscCall(DMAdaptMetric(dm, metricVec, bdLabel, rgLabel, dmCoarsened)); 2669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&metricVec)); 2679566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, *dmCoarsened)); 2685de52c6dSVaclav Hapla PetscCall(DMPlexCopy_Internal(dm, PETSC_TRUE, PETSC_TRUE, *dmCoarsened)); 2699566063dSJacob Faibussowitsch if (localized) PetscCall(DMLocalizeCoordinates(*dmCoarsened)); 270*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2710d1cd5e0SMatthew G. Knepley } 2720d1cd5e0SMatthew G. Knepley 273d71ae5a4SJacob Faibussowitsch PetscErrorCode DMAdaptLabel_Plex(DM dm, PETSC_UNUSED Vec metric, DMLabel adaptLabel, PETSC_UNUSED DMLabel rgLabel, DM *dmAdapted) 274d71ae5a4SJacob Faibussowitsch { 2750d1cd5e0SMatthew G. Knepley IS flagIS; 2760d1cd5e0SMatthew G. Knepley const PetscInt *flags; 2770d1cd5e0SMatthew G. Knepley PetscInt defFlag, minFlag, maxFlag, numFlags, f; 2780d1cd5e0SMatthew G. Knepley 2790d1cd5e0SMatthew G. Knepley PetscFunctionBegin; 2809566063dSJacob Faibussowitsch PetscCall(DMLabelGetDefaultValue(adaptLabel, &defFlag)); 2810d1cd5e0SMatthew G. Knepley minFlag = defFlag; 2820d1cd5e0SMatthew G. Knepley maxFlag = defFlag; 2839566063dSJacob Faibussowitsch PetscCall(DMLabelGetValueIS(adaptLabel, &flagIS)); 2849566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(flagIS, &numFlags)); 2859566063dSJacob Faibussowitsch PetscCall(ISGetIndices(flagIS, &flags)); 2860d1cd5e0SMatthew G. Knepley for (f = 0; f < numFlags; ++f) { 2870d1cd5e0SMatthew G. Knepley const PetscInt flag = flags[f]; 2880d1cd5e0SMatthew G. Knepley 2890d1cd5e0SMatthew G. Knepley minFlag = PetscMin(minFlag, flag); 2900d1cd5e0SMatthew G. Knepley maxFlag = PetscMax(maxFlag, flag); 2910d1cd5e0SMatthew G. Knepley } 2929566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(flagIS, &flags)); 2939566063dSJacob Faibussowitsch PetscCall(ISDestroy(&flagIS)); 2940d1cd5e0SMatthew G. Knepley { 2950d1cd5e0SMatthew G. Knepley PetscInt minMaxFlag[2], minMaxFlagGlobal[2]; 2960d1cd5e0SMatthew G. Knepley 2970d1cd5e0SMatthew G. Knepley minMaxFlag[0] = minFlag; 2980d1cd5e0SMatthew G. Knepley minMaxFlag[1] = -maxFlag; 2999566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(minMaxFlag, minMaxFlagGlobal, 2, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm))); 3000d1cd5e0SMatthew G. Knepley minFlag = minMaxFlagGlobal[0]; 3010d1cd5e0SMatthew G. Knepley maxFlag = -minMaxFlagGlobal[1]; 3020d1cd5e0SMatthew G. Knepley } 3030d1cd5e0SMatthew G. Knepley if (minFlag == maxFlag) { 3040d1cd5e0SMatthew G. Knepley switch (minFlag) { 305d71ae5a4SJacob Faibussowitsch case DM_ADAPT_DETERMINE: 306d71ae5a4SJacob Faibussowitsch *dmAdapted = NULL; 307d71ae5a4SJacob Faibussowitsch break; 3080d1cd5e0SMatthew G. Knepley case DM_ADAPT_REFINE: 3099566063dSJacob Faibussowitsch PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE)); 3109371c9d4SSatish Balay PetscCall(DMRefine(dm, MPI_COMM_NULL, dmAdapted)); 3119371c9d4SSatish Balay break; 312d71ae5a4SJacob Faibussowitsch case DM_ADAPT_COARSEN: 313d71ae5a4SJacob Faibussowitsch PetscCall(DMCoarsen(dm, MPI_COMM_NULL, dmAdapted)); 314d71ae5a4SJacob Faibussowitsch break; 315d71ae5a4SJacob Faibussowitsch default: 316d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "DMPlex does not support refinement flag %" PetscInt_FMT, minFlag); 3170d1cd5e0SMatthew G. Knepley } 3180d1cd5e0SMatthew G. Knepley } else { 3199566063dSJacob Faibussowitsch PetscCall(DMPlexSetRefinementUniform(dm, PETSC_FALSE)); 3209566063dSJacob Faibussowitsch PetscCall(DMPlexRefine_Internal(dm, NULL, adaptLabel, NULL, dmAdapted)); 3210d1cd5e0SMatthew G. Knepley } 322*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3230d1cd5e0SMatthew G. Knepley } 324