10bbe5d31Sbarral #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2497880caSRichard Tran Mills #if defined(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: 40b458e8f1SJose E. Roman SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP, "DMPlex does not support refinement flag %D\n", refFlag); 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); 7192fd8e1eSJed Brown ierr = DMGetLocalSection(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 1603a074057SBarry Smith /* 1613a074057SBarry Smith Contains the list of registered DMPlexGenerators routines 1623a074057SBarry Smith */ 163*3f2a96e3SMatthew G. Knepley extern PlexGeneratorFunctionList DMPlexGenerateList; 1643a074057SBarry Smith 1650d1cd5e0SMatthew G. Knepley PetscErrorCode DMPlexRefine_Internal(DM dm, DMLabel adaptLabel, DM *dmRefined) 1660d1cd5e0SMatthew G. Knepley { 167*3f2a96e3SMatthew G. Knepley PlexGeneratorFunctionList fl; 168800b850cSBarry Smith PetscErrorCode (*refine)(DM,PetscReal*,DM*); 169*3f2a96e3SMatthew G. Knepley PetscErrorCode (*adapt)(DM,DMLabel,DM*); 170*3f2a96e3SMatthew G. Knepley PetscErrorCode (*refinementFunc)(const PetscReal [], PetscReal *); 171*3f2a96e3SMatthew G. Knepley char genname[PETSC_MAX_PATH_LEN], *name = NULL; 172*3f2a96e3SMatthew G. Knepley PetscReal refinementLimit; 1733a074057SBarry Smith PetscReal *maxVolumes; 174*3f2a96e3SMatthew G. Knepley PetscInt dim, cStart, cEnd, c; 175*3f2a96e3SMatthew G. Knepley PetscBool flg, flg2, localized; 176*3f2a96e3SMatthew G. Knepley PetscErrorCode ierr; 1770d1cd5e0SMatthew G. Knepley 1780d1cd5e0SMatthew G. Knepley PetscFunctionBegin; 1790d1cd5e0SMatthew G. Knepley ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr); 1800d1cd5e0SMatthew G. Knepley ierr = DMPlexGetRefinementLimit(dm, &refinementLimit);CHKERRQ(ierr); 1810d1cd5e0SMatthew G. Knepley ierr = DMPlexGetRefinementFunction(dm, &refinementFunc);CHKERRQ(ierr); 1820d1cd5e0SMatthew G. Knepley if (refinementLimit == 0.0 && !refinementFunc && !adaptLabel) PetscFunctionReturn(0); 1830d1cd5e0SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1840d1cd5e0SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 185*3f2a96e3SMatthew G. Knepley ierr = PetscOptionsGetString(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_adaptor", genname, sizeof(genname), &flg);CHKERRQ(ierr); 1860d1cd5e0SMatthew G. Knepley if (flg) name = genname; 187*3f2a96e3SMatthew G. Knepley else { 188*3f2a96e3SMatthew G. Knepley ierr = PetscOptionsGetString(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_generator", genname, sizeof(genname), &flg2);CHKERRQ(ierr); 189*3f2a96e3SMatthew G. Knepley if (flg2) name = genname; 190*3f2a96e3SMatthew G. Knepley } 1917d99616aSKarl Rupp 1923a074057SBarry Smith fl = DMPlexGenerateList; 1933a074057SBarry Smith if (name) { 1943a074057SBarry Smith while (fl) { 1953a074057SBarry Smith ierr = PetscStrcmp(fl->name,name,&flg);CHKERRQ(ierr); 1963a074057SBarry Smith if (flg) { 1973a074057SBarry Smith refine = fl->refine; 198*3f2a96e3SMatthew G. Knepley adapt = fl->adaptlabel; 1993a074057SBarry Smith goto gotit; 2003a074057SBarry Smith } 2013a074057SBarry Smith fl = fl->next; 2023a074057SBarry Smith } 20389c010cfSBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Grid refiner %s not registered",name); 2043a074057SBarry Smith } else { 2053a074057SBarry Smith while (fl) { 2063a074057SBarry Smith if (dim-1 == fl->dim) { 2073a074057SBarry Smith refine = fl->refine; 208*3f2a96e3SMatthew G. Knepley adapt = fl->adaptlabel; 2093a074057SBarry Smith goto gotit; 2103a074057SBarry Smith } 2113a074057SBarry Smith fl = fl->next; 2123a074057SBarry Smith } 213367003a6SStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No grid refiner of dimension %D registered",dim); 2143a074057SBarry Smith } 2153a074057SBarry Smith 216*3f2a96e3SMatthew G. Knepley gotit: 217*3f2a96e3SMatthew G. Knepley switch (dim) { 2183a074057SBarry Smith case 2: 2190d1cd5e0SMatthew G. Knepley case 3: 220*3f2a96e3SMatthew G. Knepley if (adapt) { 221*3f2a96e3SMatthew G. Knepley ierr = (*adapt)(dm, adaptLabel, dmRefined);CHKERRQ(ierr); 222*3f2a96e3SMatthew G. Knepley } else { 2230d1cd5e0SMatthew G. Knepley ierr = PetscMalloc1(cEnd - cStart, &maxVolumes);CHKERRQ(ierr); 2240d1cd5e0SMatthew G. Knepley if (adaptLabel) { 2250d1cd5e0SMatthew G. Knepley ierr = DMPlexLabelToVolumeConstraint(dm, adaptLabel, cStart, cEnd, PETSC_DEFAULT, maxVolumes);CHKERRQ(ierr); 2260d1cd5e0SMatthew G. Knepley } else if (refinementFunc) { 2270d1cd5e0SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 2280d1cd5e0SMatthew G. Knepley PetscReal vol, centroid[3]; 2290d1cd5e0SMatthew G. Knepley 2300d1cd5e0SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL);CHKERRQ(ierr); 2310d1cd5e0SMatthew G. Knepley ierr = (*refinementFunc)(centroid, &maxVolumes[c-cStart]);CHKERRQ(ierr); 2320d1cd5e0SMatthew G. Knepley } 2330d1cd5e0SMatthew G. Knepley } else { 2340d1cd5e0SMatthew G. Knepley for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit; 2350d1cd5e0SMatthew G. Knepley } 2363a074057SBarry Smith ierr = (*refine)(dm, maxVolumes, dmRefined);CHKERRQ(ierr); 2370d1cd5e0SMatthew G. Knepley ierr = PetscFree(maxVolumes);CHKERRQ(ierr); 238*3f2a96e3SMatthew G. Knepley } 2390d1cd5e0SMatthew G. Knepley break; 240*3f2a96e3SMatthew G. Knepley default: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Mesh refinement in dimension %D is not supported.", dim); 2410d1cd5e0SMatthew G. Knepley } 2420d1cd5e0SMatthew G. Knepley ierr = DMCopyBoundary(dm, *dmRefined);CHKERRQ(ierr); 2430d1cd5e0SMatthew G. Knepley if (localized) {ierr = DMLocalizeCoordinates(*dmRefined);CHKERRQ(ierr);} 2440d1cd5e0SMatthew G. Knepley PetscFunctionReturn(0); 2450d1cd5e0SMatthew G. Knepley } 2460d1cd5e0SMatthew G. Knepley 2470d1cd5e0SMatthew G. Knepley PetscErrorCode DMPlexCoarsen_Internal(DM dm, DMLabel adaptLabel, DM *dmCoarsened) 2480d1cd5e0SMatthew G. Knepley { 2490d1cd5e0SMatthew G. Knepley Vec metricVec; 2500d1cd5e0SMatthew G. Knepley PetscInt cStart, cEnd, vStart, vEnd; 2510d1cd5e0SMatthew G. Knepley DMLabel bdLabel = NULL; 2520d1cd5e0SMatthew G. Knepley char bdLabelName[PETSC_MAX_PATH_LEN]; 2530d1cd5e0SMatthew G. Knepley PetscBool localized, flg; 2540d1cd5e0SMatthew G. Knepley PetscErrorCode ierr; 2550d1cd5e0SMatthew G. Knepley 2560d1cd5e0SMatthew G. Knepley PetscFunctionBegin; 2570d1cd5e0SMatthew G. Knepley ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr); 2580d1cd5e0SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 2590d1cd5e0SMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 2600d1cd5e0SMatthew G. Knepley ierr = DMPlexLabelToMetricConstraint(dm, adaptLabel, cStart, cEnd, vStart, vEnd, PETSC_DEFAULT, &metricVec);CHKERRQ(ierr); 261589a23caSBarry Smith ierr = PetscOptionsGetString(NULL, dm->hdr.prefix, "-dm_plex_coarsen_bd_label", bdLabelName, sizeof(bdLabelName), &flg);CHKERRQ(ierr); 2620d1cd5e0SMatthew G. Knepley if (flg) {ierr = DMGetLabel(dm, bdLabelName, &bdLabel);CHKERRQ(ierr);} 2630d1cd5e0SMatthew G. Knepley ierr = DMAdaptMetric_Plex(dm, metricVec, bdLabel, dmCoarsened);CHKERRQ(ierr); 2640d1cd5e0SMatthew G. Knepley ierr = VecDestroy(&metricVec);CHKERRQ(ierr); 2650d1cd5e0SMatthew G. Knepley if (localized) {ierr = DMLocalizeCoordinates(*dmCoarsened);CHKERRQ(ierr);} 2660d1cd5e0SMatthew G. Knepley PetscFunctionReturn(0); 2670d1cd5e0SMatthew G. Knepley } 2680d1cd5e0SMatthew G. Knepley 2690d1cd5e0SMatthew G. Knepley PetscErrorCode DMAdaptLabel_Plex(DM dm, DMLabel adaptLabel, DM *dmAdapted) 2700d1cd5e0SMatthew G. Knepley { 2710d1cd5e0SMatthew G. Knepley IS flagIS; 2720d1cd5e0SMatthew G. Knepley const PetscInt *flags; 2730d1cd5e0SMatthew G. Knepley PetscInt defFlag, minFlag, maxFlag, numFlags, f; 2740d1cd5e0SMatthew G. Knepley PetscErrorCode ierr; 2750d1cd5e0SMatthew G. Knepley 2760d1cd5e0SMatthew G. Knepley PetscFunctionBegin; 2770d1cd5e0SMatthew G. Knepley ierr = DMLabelGetDefaultValue(adaptLabel, &defFlag);CHKERRQ(ierr); 2780d1cd5e0SMatthew G. Knepley minFlag = defFlag; 2790d1cd5e0SMatthew G. Knepley maxFlag = defFlag; 2800d1cd5e0SMatthew G. Knepley ierr = DMLabelGetValueIS(adaptLabel, &flagIS);CHKERRQ(ierr); 2810d1cd5e0SMatthew G. Knepley ierr = ISGetLocalSize(flagIS, &numFlags);CHKERRQ(ierr); 2820d1cd5e0SMatthew G. Knepley ierr = ISGetIndices(flagIS, &flags);CHKERRQ(ierr); 2830d1cd5e0SMatthew G. Knepley for (f = 0; f < numFlags; ++f) { 2840d1cd5e0SMatthew G. Knepley const PetscInt flag = flags[f]; 2850d1cd5e0SMatthew G. Knepley 2860d1cd5e0SMatthew G. Knepley minFlag = PetscMin(minFlag, flag); 2870d1cd5e0SMatthew G. Knepley maxFlag = PetscMax(maxFlag, flag); 2880d1cd5e0SMatthew G. Knepley } 2890d1cd5e0SMatthew G. Knepley ierr = ISRestoreIndices(flagIS, &flags);CHKERRQ(ierr); 2900d1cd5e0SMatthew G. Knepley ierr = ISDestroy(&flagIS);CHKERRQ(ierr); 2910d1cd5e0SMatthew G. Knepley { 2920d1cd5e0SMatthew G. Knepley PetscInt minMaxFlag[2], minMaxFlagGlobal[2]; 2930d1cd5e0SMatthew G. Knepley 2940d1cd5e0SMatthew G. Knepley minMaxFlag[0] = minFlag; 2950d1cd5e0SMatthew G. Knepley minMaxFlag[1] = -maxFlag; 2960d1cd5e0SMatthew G. Knepley ierr = MPI_Allreduce(minMaxFlag, minMaxFlagGlobal, 2, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 2970d1cd5e0SMatthew G. Knepley minFlag = minMaxFlagGlobal[0]; 2980d1cd5e0SMatthew G. Knepley maxFlag = -minMaxFlagGlobal[1]; 2990d1cd5e0SMatthew G. Knepley } 3000d1cd5e0SMatthew G. Knepley if (minFlag == maxFlag) { 3010d1cd5e0SMatthew G. Knepley switch (minFlag) { 3020d1cd5e0SMatthew G. Knepley case DM_ADAPT_DETERMINE: 3030d1cd5e0SMatthew G. Knepley *dmAdapted = NULL;break; 3040d1cd5e0SMatthew G. Knepley case DM_ADAPT_REFINE: 3050d1cd5e0SMatthew G. Knepley ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE);CHKERRQ(ierr); 3060d1cd5e0SMatthew G. Knepley ierr = DMRefine(dm, MPI_COMM_NULL, dmAdapted);CHKERRQ(ierr);break; 3070d1cd5e0SMatthew G. Knepley case DM_ADAPT_COARSEN: 3080d1cd5e0SMatthew G. Knepley ierr = DMCoarsen(dm, MPI_COMM_NULL, dmAdapted);CHKERRQ(ierr);break; 3090d1cd5e0SMatthew G. Knepley default: 310b458e8f1SJose E. Roman SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP,"DMPlex does not support refinement flag %D\n", minFlag); 3110d1cd5e0SMatthew G. Knepley } 3120d1cd5e0SMatthew G. Knepley } else { 3130d1cd5e0SMatthew G. Knepley ierr = DMPlexSetRefinementUniform(dm, PETSC_FALSE);CHKERRQ(ierr); 3140d1cd5e0SMatthew G. Knepley ierr = DMPlexRefine_Internal(dm, adaptLabel, dmAdapted);CHKERRQ(ierr); 3150d1cd5e0SMatthew G. Knepley } 3160d1cd5e0SMatthew G. Knepley PetscFunctionReturn(0); 3170d1cd5e0SMatthew G. Knepley } 3180d1cd5e0SMatthew G. Knepley 319f7bcbcbbSMatthew G. Knepley /* 3200d1cd5e0SMatthew G. Knepley DMAdaptMetric_Plex - Generates a new mesh conforming to a metric field. 3210bbe5d31Sbarral 322f7bcbcbbSMatthew G. Knepley Input Parameters: 323f7bcbcbbSMatthew G. Knepley + dm - The DM object 32451ff0a1cSMatthew G. Knepley . vertexMetric - The metric to which the mesh is adapted, defined vertex-wise in a LOCAL vector 3250d1cd5e0SMatthew G. Knepley - bdLabel - Label for boundary tags which are preserved in dmNew, or NULL. Should not be named "_boundary_". 326f7bcbcbbSMatthew G. Knepley 327f7bcbcbbSMatthew G. Knepley Output Parameter: 328f7bcbcbbSMatthew G. Knepley . dmNew - the new DM 329f7bcbcbbSMatthew G. Knepley 330f7bcbcbbSMatthew G. Knepley Level: advanced 331f7bcbcbbSMatthew G. Knepley 332f7bcbcbbSMatthew G. Knepley .seealso: DMCoarsen(), DMRefine() 333f7bcbcbbSMatthew G. Knepley */ 3340d1cd5e0SMatthew G. Knepley PetscErrorCode DMAdaptMetric_Plex(DM dm, Vec vertexMetric, DMLabel bdLabel, DM *dmNew) 335f7bcbcbbSMatthew G. Knepley { 336497880caSRichard Tran Mills #if defined(PETSC_HAVE_PRAGMATIC) 337bc1d7e78SMatthew G. Knepley MPI_Comm comm; 338f7bcbcbbSMatthew G. Knepley const char *bdName = "_boundary_"; 3392b5f8d01SMatthew G. Knepley #if 0 3402b5f8d01SMatthew G. Knepley DM odm = dm; 3412b5f8d01SMatthew G. Knepley #endif 342f7bcbcbbSMatthew G. Knepley DM udm, cdm; 3430d1cd5e0SMatthew G. Knepley DMLabel bdLabelFull; 3440d1cd5e0SMatthew G. Knepley const char *bdLabelName; 3451838b4a6SMatthew G. Knepley IS bdIS, globalVertexNum; 346f7bcbcbbSMatthew G. Knepley PetscSection coordSection; 347f7bcbcbbSMatthew G. Knepley Vec coordinates; 348f7bcbcbbSMatthew G. Knepley const PetscScalar *coords, *met; 3491838b4a6SMatthew G. Knepley const PetscInt *bdFacesFull, *gV; 3501838b4a6SMatthew G. Knepley PetscInt *bdFaces, *bdFaceIds, *l2gv; 351f7bcbcbbSMatthew G. Knepley PetscReal *x, *y, *z, *metric; 3522ee120b0SMatthew G. Knepley PetscInt *cells; 3531838b4a6SMatthew G. Knepley PetscInt dim, cStart, cEnd, numCells, c, coff, vStart, vEnd, numVertices, numLocVertices, v; 3542ee120b0SMatthew G. Knepley PetscInt off, maxConeSize, numBdFaces, f, bdSize; 355f7bcbcbbSMatthew G. Knepley PetscBool flg; 3562ee120b0SMatthew G. Knepley DMLabel bdLabelNew; 357a4a685f2SJacob Faibussowitsch PetscReal *coordsNew; 3582ee120b0SMatthew G. Knepley PetscInt *bdTags; 3592ee120b0SMatthew G. Knepley PetscReal *xNew[3] = {NULL, NULL, NULL}; 3602ee120b0SMatthew G. Knepley PetscInt *cellsNew; 3612ee120b0SMatthew G. Knepley PetscInt d, numCellsNew, numVerticesNew; 3622ee120b0SMatthew G. Knepley PetscInt numCornersNew, fStart, fEnd; 363bc1d7e78SMatthew G. Knepley PetscMPIInt numProcs; 364f7bcbcbbSMatthew G. Knepley PetscErrorCode ierr; 365f7bcbcbbSMatthew G. Knepley 366f7bcbcbbSMatthew G. Knepley PetscFunctionBegin; 367bbc23918SMatthew G. Knepley /* Check for FEM adjacency flags */ 368bc1d7e78SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 369bc1d7e78SMatthew G. Knepley ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 3700d1cd5e0SMatthew G. Knepley if (bdLabel) { 371d67d17b1SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) bdLabel, &bdLabelName);CHKERRQ(ierr); 372f7bcbcbbSMatthew G. Knepley ierr = PetscStrcmp(bdLabelName, bdName, &flg);CHKERRQ(ierr); 373bc1d7e78SMatthew G. Knepley if (flg) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName); 374f7bcbcbbSMatthew G. Knepley } 375a2e1bd45SMatthew G. Knepley /* Add overlap for Pragmatic */ 376bbc23918SMatthew G. Knepley #if 0 377bbc23918SMatthew G. Knepley /* Check for overlap by looking for cell in the SF */ 378bbc23918SMatthew G. Knepley if (!overlapped) { 379a2e1bd45SMatthew G. Knepley ierr = DMPlexDistributeOverlap(odm, 1, NULL, &dm);CHKERRQ(ierr); 380a2e1bd45SMatthew G. Knepley if (!dm) {dm = odm; ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);} 381bbc23918SMatthew G. Knepley } 382bbc23918SMatthew G. Knepley #endif 383f7bcbcbbSMatthew G. Knepley /* Get mesh information */ 384f7bcbcbbSMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 385f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 386f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 387f7bcbcbbSMatthew G. Knepley ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr); 388f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr); 389f7bcbcbbSMatthew G. Knepley numCells = cEnd - cStart; 390f7bcbcbbSMatthew G. Knepley numVertices = vEnd - vStart; 391f7bcbcbbSMatthew G. Knepley ierr = PetscCalloc5(numVertices, &x, numVertices, &y, numVertices, &z, numVertices*PetscSqr(dim), &metric, numCells*maxConeSize, &cells);CHKERRQ(ierr); 392f7bcbcbbSMatthew G. Knepley for (c = 0, coff = 0; c < numCells; ++c) { 393f7bcbcbbSMatthew G. Knepley const PetscInt *cone; 394f7bcbcbbSMatthew G. Knepley PetscInt coneSize, cl; 395f7bcbcbbSMatthew G. Knepley 396f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr); 397f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr); 398f7bcbcbbSMatthew G. Knepley for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart; 399f7bcbcbbSMatthew G. Knepley } 4001838b4a6SMatthew G. Knepley ierr = PetscCalloc1(numVertices, &l2gv);CHKERRQ(ierr); 4011838b4a6SMatthew G. Knepley ierr = DMPlexGetVertexNumbering(udm, &globalVertexNum);CHKERRQ(ierr); 4021838b4a6SMatthew G. Knepley ierr = ISGetIndices(globalVertexNum, &gV);CHKERRQ(ierr); 4031838b4a6SMatthew G. Knepley for (v = 0, numLocVertices = 0; v < numVertices; ++v) { 4041838b4a6SMatthew G. Knepley if (gV[v] >= 0) ++numLocVertices; 4051838b4a6SMatthew G. Knepley l2gv[v] = gV[v] < 0 ? -(gV[v]+1) : gV[v]; 4061838b4a6SMatthew G. Knepley } 4071838b4a6SMatthew G. Knepley ierr = ISRestoreIndices(globalVertexNum, &gV);CHKERRQ(ierr); 408f7bcbcbbSMatthew G. Knepley ierr = DMDestroy(&udm);CHKERRQ(ierr); 409f7bcbcbbSMatthew G. Knepley ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 41092fd8e1eSJed Brown ierr = DMGetLocalSection(cdm, &coordSection);CHKERRQ(ierr); 411f7bcbcbbSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 412f7bcbcbbSMatthew G. Knepley ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 413f7bcbcbbSMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 414f7bcbcbbSMatthew G. Knepley ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 4155f2e139eSMatthew G. Knepley x[v-vStart] = PetscRealPart(coords[off+0]); 4165f2e139eSMatthew G. Knepley if (dim > 1) y[v-vStart] = PetscRealPart(coords[off+1]); 4175f2e139eSMatthew G. Knepley if (dim > 2) z[v-vStart] = PetscRealPart(coords[off+2]); 418f7bcbcbbSMatthew G. Knepley } 419f7bcbcbbSMatthew G. Knepley ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 420f7bcbcbbSMatthew G. Knepley /* Get boundary mesh */ 421d67d17b1SMatthew G. Knepley ierr = DMLabelCreate(PETSC_COMM_SELF, bdName, &bdLabelFull);CHKERRQ(ierr); 422e752be1aSMatthew G. Knepley ierr = DMPlexMarkBoundaryFaces(dm, 1, bdLabelFull);CHKERRQ(ierr); 423f7bcbcbbSMatthew G. Knepley ierr = DMLabelGetStratumIS(bdLabelFull, 1, &bdIS);CHKERRQ(ierr); 424f7bcbcbbSMatthew G. Knepley ierr = DMLabelGetStratumSize(bdLabelFull, 1, &numBdFaces);CHKERRQ(ierr); 425f7bcbcbbSMatthew G. Knepley ierr = ISGetIndices(bdIS, &bdFacesFull);CHKERRQ(ierr); 426f7bcbcbbSMatthew G. Knepley for (f = 0, bdSize = 0; f < numBdFaces; ++f) { 427f7bcbcbbSMatthew G. Knepley PetscInt *closure = NULL; 428f7bcbcbbSMatthew G. Knepley PetscInt closureSize, cl; 429f7bcbcbbSMatthew G. Knepley 430f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 431f7bcbcbbSMatthew G. Knepley for (cl = 0; cl < closureSize*2; cl += 2) { 432f7bcbcbbSMatthew G. Knepley if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize; 433f7bcbcbbSMatthew G. Knepley } 434f7bcbcbbSMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 435f7bcbcbbSMatthew G. Knepley } 436f7bcbcbbSMatthew G. Knepley ierr = PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds);CHKERRQ(ierr); 437f7bcbcbbSMatthew G. Knepley for (f = 0, bdSize = 0; f < numBdFaces; ++f) { 438f7bcbcbbSMatthew G. Knepley PetscInt *closure = NULL; 439f7bcbcbbSMatthew G. Knepley PetscInt closureSize, cl; 440f7bcbcbbSMatthew G. Knepley 441f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 442f7bcbcbbSMatthew G. Knepley for (cl = 0; cl < closureSize*2; cl += 2) { 443f7bcbcbbSMatthew G. Knepley if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart; 444f7bcbcbbSMatthew G. Knepley } 445f7bcbcbbSMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 446f7bcbcbbSMatthew G. Knepley if (bdLabel) {ierr = DMLabelGetValue(bdLabel, bdFacesFull[f], &bdFaceIds[f]);CHKERRQ(ierr);} 447f7bcbcbbSMatthew G. Knepley else {bdFaceIds[f] = 1;} 448f7bcbcbbSMatthew G. Knepley } 449f7bcbcbbSMatthew G. Knepley ierr = ISDestroy(&bdIS);CHKERRQ(ierr); 450f7bcbcbbSMatthew G. Knepley ierr = DMLabelDestroy(&bdLabelFull);CHKERRQ(ierr); 451f7bcbcbbSMatthew G. Knepley /* Get metric */ 45204996bafSMatthew G. Knepley ierr = VecViewFromOptions(vertexMetric, NULL, "-adapt_metric_view");CHKERRQ(ierr); 453f7bcbcbbSMatthew G. Knepley ierr = VecGetArrayRead(vertexMetric, &met);CHKERRQ(ierr); 4545f2e139eSMatthew G. Knepley for (v = 0; v < (vEnd-vStart)*PetscSqr(dim); ++v) metric[v] = PetscRealPart(met[v]); 455f7bcbcbbSMatthew G. Knepley ierr = VecRestoreArrayRead(vertexMetric, &met);CHKERRQ(ierr); 456bbc23918SMatthew G. Knepley #if 0 457a2e1bd45SMatthew G. Knepley /* Destroy overlap mesh */ 458a2e1bd45SMatthew G. Knepley ierr = DMDestroy(&dm);CHKERRQ(ierr); 459bbc23918SMatthew G. Knepley #endif 460f7bcbcbbSMatthew G. Knepley /* Create new mesh */ 461f7bcbcbbSMatthew G. Knepley switch (dim) { 462f7bcbcbbSMatthew G. Knepley case 2: 463cc29c497SMatthew G. Knepley pragmatic_2d_mpi_init(&numVertices, &numCells, cells, x, y, l2gv, numLocVertices, comm);break; 464f7bcbcbbSMatthew G. Knepley case 3: 465cc29c497SMatthew G. Knepley pragmatic_3d_mpi_init(&numVertices, &numCells, cells, x, y, z, l2gv, numLocVertices, comm);break; 466a4a685f2SJacob Faibussowitsch default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %D", dim); 467f7bcbcbbSMatthew G. Knepley } 468f7bcbcbbSMatthew G. Knepley pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds); 469f7bcbcbbSMatthew G. Knepley pragmatic_set_metric(metric); 4700d1cd5e0SMatthew G. Knepley pragmatic_adapt(((DM_Plex *) dm->data)->remeshBd ? 1 : 0); 4711838b4a6SMatthew G. Knepley ierr = PetscFree(l2gv);CHKERRQ(ierr); 472f7bcbcbbSMatthew G. Knepley /* Read out mesh */ 4738713909fSMatthew G. Knepley pragmatic_get_info_mpi(&numVerticesNew, &numCellsNew); 474f7bcbcbbSMatthew G. Knepley ierr = PetscMalloc1(numVerticesNew*dim, &coordsNew);CHKERRQ(ierr); 475f7bcbcbbSMatthew G. Knepley switch (dim) { 476f7bcbcbbSMatthew G. Knepley case 2: 477f7bcbcbbSMatthew G. Knepley numCornersNew = 3; 478f7bcbcbbSMatthew G. Knepley ierr = PetscMalloc2(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1]);CHKERRQ(ierr); 4798713909fSMatthew G. Knepley pragmatic_get_coords_2d_mpi(xNew[0], xNew[1]); 480f7bcbcbbSMatthew G. Knepley break; 481f7bcbcbbSMatthew G. Knepley case 3: 482f7bcbcbbSMatthew G. Knepley numCornersNew = 4; 483f7bcbcbbSMatthew G. Knepley ierr = PetscMalloc3(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1], numVerticesNew, &xNew[2]);CHKERRQ(ierr); 4848713909fSMatthew G. Knepley pragmatic_get_coords_3d_mpi(xNew[0], xNew[1], xNew[2]); 485f7bcbcbbSMatthew G. Knepley break; 486f7bcbcbbSMatthew G. Knepley default: 487a4a685f2SJacob Faibussowitsch SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %D", dim); 488f7bcbcbbSMatthew G. Knepley } 489a4a685f2SJacob Faibussowitsch for (v = 0; v < numVerticesNew; ++v) {for (d = 0; d < dim; ++d) coordsNew[v*dim+d] = xNew[d][v];} 490f7bcbcbbSMatthew G. Knepley ierr = PetscMalloc1(numCellsNew*(dim+1), &cellsNew);CHKERRQ(ierr); 491f7bcbcbbSMatthew G. Knepley pragmatic_get_elements(cellsNew); 49225b6865aSVaclav Hapla ierr = DMPlexCreateFromCellListParallelPetsc(comm, dim, numCellsNew, numVerticesNew, PETSC_DECIDE, numCornersNew, PETSC_TRUE, cellsNew, dim, coordsNew, NULL, dmNew);CHKERRQ(ierr); 493f7bcbcbbSMatthew G. Knepley /* Read out boundary label */ 494f7bcbcbbSMatthew G. Knepley pragmatic_get_boundaryTags(&bdTags); 495f7bcbcbbSMatthew G. Knepley ierr = DMCreateLabel(*dmNew, bdLabel ? bdLabelName : bdName);CHKERRQ(ierr); 496f7bcbcbbSMatthew G. Knepley ierr = DMGetLabel(*dmNew, bdLabel ? bdLabelName : bdName, &bdLabelNew);CHKERRQ(ierr); 497f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd);CHKERRQ(ierr); 498f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd);CHKERRQ(ierr); 499f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd);CHKERRQ(ierr); 500f7bcbcbbSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 501f7bcbcbbSMatthew G. Knepley /* Only for simplicial meshes */ 502f7bcbcbbSMatthew G. Knepley coff = (c-cStart)*(dim+1); 5038713909fSMatthew G. Knepley /* d is the local cell number of the vertex opposite to the face we are marking */ 504f7bcbcbbSMatthew G. Knepley for (d = 0; d < dim+1; ++d) { 505f7bcbcbbSMatthew G. Knepley if (bdTags[coff+d]) { 5068713909fSMatthew 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 */ 5078713909fSMatthew G. Knepley const PetscInt *cone; 508f7bcbcbbSMatthew G. Knepley 5098713909fSMatthew G. Knepley /* Mark face opposite to this vertex: This pattern is specified in DMPlexGetRawFaces_Internal() */ 5108713909fSMatthew G. Knepley ierr = DMPlexGetCone(*dmNew, c, &cone);CHKERRQ(ierr); 5118713909fSMatthew G. Knepley ierr = DMLabelSetValue(bdLabelNew, cone[perm[dim][d]], bdTags[coff+d]);CHKERRQ(ierr); 512f7bcbcbbSMatthew G. Knepley } 513f7bcbcbbSMatthew G. Knepley } 514f7bcbcbbSMatthew G. Knepley } 515f7bcbcbbSMatthew G. Knepley /* Cleanup */ 516f7bcbcbbSMatthew G. Knepley switch (dim) { 517f7bcbcbbSMatthew G. Knepley case 2: ierr = PetscFree2(xNew[0], xNew[1]);CHKERRQ(ierr);break; 518f7bcbcbbSMatthew G. Knepley case 3: ierr = PetscFree3(xNew[0], xNew[1], xNew[2]);CHKERRQ(ierr);break; 519f7bcbcbbSMatthew G. Knepley } 520f7bcbcbbSMatthew G. Knepley ierr = PetscFree(cellsNew);CHKERRQ(ierr); 521f7bcbcbbSMatthew G. Knepley ierr = PetscFree5(x, y, z, metric, cells);CHKERRQ(ierr); 522f7bcbcbbSMatthew G. Knepley ierr = PetscFree2(bdFaces, bdFaceIds);CHKERRQ(ierr); 523f7bcbcbbSMatthew G. Knepley ierr = PetscFree(coordsNew);CHKERRQ(ierr); 524f7bcbcbbSMatthew G. Knepley pragmatic_finalize(); 525f7bcbcbbSMatthew G. Knepley PetscFunctionReturn(0); 5266f25b0d8SLisandro Dalcin #else 5276f25b0d8SLisandro Dalcin PetscFunctionBegin; 5286f25b0d8SLisandro Dalcin SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Remeshing needs external package support.\nPlease reconfigure with --download-pragmatic."); 5296f25b0d8SLisandro Dalcin #endif 530f7bcbcbbSMatthew G. Knepley } 531