1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2af0996ceSBarry Smith #include <petsc/private/matorderimpl.h> /*I "petscmat.h" I*/ 38ed5f475SMatthew G. Knepley 4a1fd77bcSMatthew G. Knepley static PetscErrorCode DMPlexCreateOrderingClosure_Static(DM dm, PetscInt numPoints, const PetscInt pperm[], PetscInt **clperm, PetscInt **invclperm) 58ed5f475SMatthew G. Knepley { 68ed5f475SMatthew G. Knepley PetscInt *perm, *iperm; 78ed5f475SMatthew G. Knepley PetscInt depth, d, pStart, pEnd, fStart, fMax, fEnd, p; 88ed5f475SMatthew G. Knepley 98ed5f475SMatthew G. Knepley PetscFunctionBegin; 109566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth)); 119566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(pEnd-pStart,&perm)); 139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(pEnd-pStart,&iperm)); 148ed5f475SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) iperm[p] = -1; 158ed5f475SMatthew G. Knepley for (d = depth; d > 0; --d) { 169566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, d, &pStart, &pEnd)); 179566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, d-1, &fStart, &fEnd)); 188ed5f475SMatthew G. Knepley fMax = fStart; 198ed5f475SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 208ed5f475SMatthew G. Knepley const PetscInt *cone; 218ed5f475SMatthew G. Knepley PetscInt point, coneSize, c; 228ed5f475SMatthew G. Knepley 238ed5f475SMatthew G. Knepley if (d == depth) { 248ed5f475SMatthew G. Knepley perm[p] = pperm[p]; 258ed5f475SMatthew G. Knepley iperm[pperm[p]] = p; 268ed5f475SMatthew G. Knepley } 278ed5f475SMatthew G. Knepley point = perm[p]; 289566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, point, &coneSize)); 299566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, point, &cone)); 308ed5f475SMatthew G. Knepley for (c = 0; c < coneSize; ++c) { 318ed5f475SMatthew G. Knepley const PetscInt oldc = cone[c]; 328ed5f475SMatthew G. Knepley const PetscInt newc = iperm[oldc]; 338ed5f475SMatthew G. Knepley 348ed5f475SMatthew G. Knepley if (newc < 0) { 358ed5f475SMatthew G. Knepley perm[fMax] = oldc; 368ed5f475SMatthew G. Knepley iperm[oldc] = fMax++; 378ed5f475SMatthew G. Knepley } 388ed5f475SMatthew G. Knepley } 398ed5f475SMatthew G. Knepley } 4063a3b9bcSJacob Faibussowitsch PetscCheck(fMax == fEnd,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of depth %" PetscInt_FMT " faces %" PetscInt_FMT " does not match permuted number %" PetscInt_FMT, d, fEnd-fStart, fMax-fStart); 418ed5f475SMatthew G. Knepley } 428ed5f475SMatthew G. Knepley *clperm = perm; 438ed5f475SMatthew G. Knepley *invclperm = iperm; 448ed5f475SMatthew G. Knepley PetscFunctionReturn(0); 458ed5f475SMatthew G. Knepley } 468ed5f475SMatthew G. Knepley 478ed5f475SMatthew G. Knepley /*@ 488ed5f475SMatthew G. Knepley DMPlexGetOrdering - Calculate a reordering of the mesh 498ed5f475SMatthew G. Knepley 50d083f849SBarry Smith Collective on dm 518ed5f475SMatthew G. Knepley 52d8d19677SJose E. Roman Input Parameters: 538ed5f475SMatthew G. Knepley + dm - The DMPlex object 54c99f2573SMatthew G. Knepley . otype - type of reordering, one of the following: 558ed5f475SMatthew G. Knepley $ MATORDERINGNATURAL - Natural 568ed5f475SMatthew G. Knepley $ MATORDERINGND - Nested Dissection 578ed5f475SMatthew G. Knepley $ MATORDERING1WD - One-way Dissection 588ed5f475SMatthew G. Knepley $ MATORDERINGRCM - Reverse Cuthill-McKee 598ed5f475SMatthew G. Knepley $ MATORDERINGQMD - Quotient Minimum Degree 60c99f2573SMatthew G. Knepley - label - [Optional] Label used to segregate ordering into sets, or NULL 618ed5f475SMatthew G. Knepley 628ed5f475SMatthew G. Knepley Output Parameter: 63c99f2573SMatthew G. Knepley . perm - The point permutation as an IS, perm[old point number] = new point number 64c99f2573SMatthew G. Knepley 65c99f2573SMatthew G. Knepley Note: The label is used to group sets of points together by label value. This makes it easy to reorder a mesh which 66c99f2573SMatthew G. Knepley has different types of cells, and then loop over each set of reordered cells for assembly. 678ed5f475SMatthew G. Knepley 688ed5f475SMatthew G. Knepley Level: intermediate 698ed5f475SMatthew G. Knepley 70db781477SPatrick Sanan .seealso: `DMPlexPermute()`, `MatGetOrdering()` 718ed5f475SMatthew G. Knepley @*/ 72c99f2573SMatthew G. Knepley PetscErrorCode DMPlexGetOrdering(DM dm, MatOrderingType otype, DMLabel label, IS *perm) 738ed5f475SMatthew G. Knepley { 748ed5f475SMatthew G. Knepley PetscInt numCells = 0; 753c67d5adSSatish Balay PetscInt *start = NULL, *adjacency = NULL, *cperm, *clperm = NULL, *invclperm = NULL, *mask, *xls, pStart, pEnd, c, i; 768ed5f475SMatthew G. Knepley 778ed5f475SMatthew G. Knepley PetscFunctionBegin; 788ed5f475SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 79064a246eSJacob Faibussowitsch PetscValidPointer(perm, 4); 809566063dSJacob Faibussowitsch PetscCall(DMPlexCreateNeighborCSR(dm, 0, &numCells, &start, &adjacency)); 819566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(numCells,&cperm,numCells,&mask,numCells*2,&xls)); 82d80c2872SMatthew G. Knepley if (numCells) { 838ed5f475SMatthew G. Knepley /* Shift for Fortran numbering */ 848ed5f475SMatthew G. Knepley for (i = 0; i < start[numCells]; ++i) ++adjacency[i]; 858ed5f475SMatthew G. Knepley for (i = 0; i <= numCells; ++i) ++start[i]; 869566063dSJacob Faibussowitsch PetscCall(SPARSEPACKgenrcm(&numCells, start, adjacency, cperm, mask, xls)); 87d80c2872SMatthew G. Knepley } 889566063dSJacob Faibussowitsch PetscCall(PetscFree(start)); 899566063dSJacob Faibussowitsch PetscCall(PetscFree(adjacency)); 908ed5f475SMatthew G. Knepley /* Shift for Fortran numbering */ 918ed5f475SMatthew G. Knepley for (c = 0; c < numCells; ++c) --cperm[c]; 92c99f2573SMatthew G. Knepley /* Segregate */ 93c99f2573SMatthew G. Knepley if (label) { 94c99f2573SMatthew G. Knepley IS valueIS; 95c99f2573SMatthew G. Knepley const PetscInt *values; 96c99f2573SMatthew G. Knepley PetscInt numValues, numPoints = 0; 97c99f2573SMatthew G. Knepley PetscInt *sperm, *vsize, *voff, v; 98c99f2573SMatthew G. Knepley 999566063dSJacob Faibussowitsch PetscCall(DMLabelGetValueIS(label, &valueIS)); 1009566063dSJacob Faibussowitsch PetscCall(ISSort(valueIS)); 1019566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(valueIS, &numValues)); 1029566063dSJacob Faibussowitsch PetscCall(ISGetIndices(valueIS, &values)); 1039566063dSJacob Faibussowitsch PetscCall(PetscCalloc3(numCells,&sperm,numValues,&vsize,numValues+1,&voff)); 104c99f2573SMatthew G. Knepley for (v = 0; v < numValues; ++v) { 1059566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumSize(label, values[v], &vsize[v])); 106c99f2573SMatthew G. Knepley if (v < numValues-1) voff[v+2] += vsize[v] + voff[v+1]; 107c99f2573SMatthew G. Knepley numPoints += vsize[v]; 108c99f2573SMatthew G. Knepley } 10963a3b9bcSJacob Faibussowitsch PetscCheck(numPoints == numCells,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label only covers %" PetscInt_FMT " cells < %" PetscInt_FMT " total", numPoints, numCells); 110c99f2573SMatthew G. Knepley for (c = 0; c < numCells; ++c) { 111c99f2573SMatthew G. Knepley const PetscInt oldc = cperm[c]; 112c99f2573SMatthew G. Knepley PetscInt val, vloc; 113c99f2573SMatthew G. Knepley 1149566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(label, oldc, &val)); 11563a3b9bcSJacob Faibussowitsch PetscCheck(val != -1,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cell %" PetscInt_FMT " not present in label", oldc); 1169566063dSJacob Faibussowitsch PetscCall(PetscFindInt(val, numValues, values, &vloc)); 11763a3b9bcSJacob Faibussowitsch PetscCheck(vloc >= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Value %" PetscInt_FMT " not present label", val); 118c99f2573SMatthew G. Knepley sperm[voff[vloc+1]++] = oldc; 119c99f2573SMatthew G. Knepley } 120c99f2573SMatthew G. Knepley for (v = 0; v < numValues; ++v) { 1211dca8a05SBarry Smith PetscCheck(voff[v+1] - voff[v] == vsize[v],PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of %" PetscInt_FMT " values found is %" PetscInt_FMT " != %" PetscInt_FMT, values[v], voff[v+1] - voff[v], vsize[v]); 122c99f2573SMatthew G. Knepley } 1239566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(valueIS, &values)); 1249566063dSJacob Faibussowitsch PetscCall(ISDestroy(&valueIS)); 1259566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(cperm, sperm, numCells)); 1269566063dSJacob Faibussowitsch PetscCall(PetscFree3(sperm, vsize, voff)); 127c99f2573SMatthew G. Knepley } 1288ed5f475SMatthew G. Knepley /* Construct closure */ 1299566063dSJacob Faibussowitsch PetscCall(DMPlexCreateOrderingClosure_Static(dm, numCells, cperm, &clperm, &invclperm)); 1309566063dSJacob Faibussowitsch PetscCall(PetscFree3(cperm,mask,xls)); 1319566063dSJacob Faibussowitsch PetscCall(PetscFree(clperm)); 1328ed5f475SMatthew G. Knepley /* Invert permutation */ 1339566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 1349566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject) dm), pEnd-pStart, invclperm, PETSC_OWN_POINTER, perm)); 1358ed5f475SMatthew G. Knepley PetscFunctionReturn(0); 1368ed5f475SMatthew G. Knepley } 1378ed5f475SMatthew G. Knepley 1388ed5f475SMatthew G. Knepley /*@ 1396913077dSMatthew G. Knepley DMPlexGetOrdering1D - Reorder the vertices so that the mesh is in a line 1406913077dSMatthew G. Knepley 1416913077dSMatthew G. Knepley Collective on dm 1426913077dSMatthew G. Knepley 1436913077dSMatthew G. Knepley Input Parameter: 1446913077dSMatthew G. Knepley . dm - The DMPlex object 1456913077dSMatthew G. Knepley 1466913077dSMatthew G. Knepley Output Parameter: 1476913077dSMatthew G. Knepley . perm - The point permutation as an IS, perm[old point number] = new point number 1486913077dSMatthew G. Knepley 1496913077dSMatthew G. Knepley Level: intermediate 1506913077dSMatthew G. Knepley 151db781477SPatrick Sanan .seealso: `DMPlexGetOrdering()`, `DMPlexPermute()`, `MatGetOrdering()` 1526913077dSMatthew G. Knepley @*/ 1536913077dSMatthew G. Knepley PetscErrorCode DMPlexGetOrdering1D(DM dm, IS *perm) 1546913077dSMatthew G. Knepley { 1556913077dSMatthew G. Knepley PetscInt *points; 1566913077dSMatthew G. Knepley const PetscInt *support, *cone; 1576913077dSMatthew G. Knepley PetscInt dim, pStart, pEnd, cStart, cEnd, c, vStart, vEnd, v, suppSize, lastCell = 0; 1586913077dSMatthew G. Knepley 1596913077dSMatthew G. Knepley PetscFunctionBegin; 1609566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1616913077dSMatthew G. Knepley PetscCheck(dim == 1, PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Input mesh must be one dimensional, not %" PetscInt_FMT, dim); 1629566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 1639566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 1649566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 1659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(pEnd-pStart, &points)); 1666913077dSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) points[c] = c; 1676913077dSMatthew G. Knepley for (v = vStart; v < vEnd; ++v) points[v] = v; 1686913077dSMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 1699566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, v, &suppSize)); 1709566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(dm, v, &support)); 1716913077dSMatthew G. Knepley if (suppSize == 1) {lastCell = support[0]; break;} 1726913077dSMatthew G. Knepley } 1736913077dSMatthew G. Knepley if (v < vEnd) { 1746913077dSMatthew G. Knepley PetscInt pos = cEnd; 1756913077dSMatthew G. Knepley 1766913077dSMatthew G. Knepley points[v] = pos++; 1776913077dSMatthew G. Knepley while (lastCell >= cStart) { 1789566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, lastCell, &cone)); 1796913077dSMatthew G. Knepley if (cone[0] == v) v = cone[1]; 1806913077dSMatthew G. Knepley else v = cone[0]; 1819566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(dm, v, &support)); 1829566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, v, &suppSize)); 1836913077dSMatthew G. Knepley if (suppSize == 1) {lastCell = -1;} 1846913077dSMatthew G. Knepley else { 1856913077dSMatthew G. Knepley if (support[0] == lastCell) lastCell = support[1]; 1866913077dSMatthew G. Knepley else lastCell = support[0]; 1876913077dSMatthew G. Knepley } 1886913077dSMatthew G. Knepley points[v] = pos++; 1896913077dSMatthew G. Knepley } 1906913077dSMatthew G. Knepley PetscCheck(pos == pEnd, PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Last vertex was %" PetscInt_FMT ", not %" PetscInt_FMT, pos, pEnd); 1916913077dSMatthew G. Knepley } 1929566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject) dm), pEnd-pStart, points, PETSC_OWN_POINTER, perm)); 1936913077dSMatthew G. Knepley PetscFunctionReturn(0); 1946913077dSMatthew G. Knepley } 1956913077dSMatthew G. Knepley 1966913077dSMatthew G. Knepley /*@ 1978ed5f475SMatthew G. Knepley DMPlexPermute - Reorder the mesh according to the input permutation 1988ed5f475SMatthew G. Knepley 199d083f849SBarry Smith Collective on dm 2008ed5f475SMatthew G. Knepley 201d8d19677SJose E. Roman Input Parameters: 2028ed5f475SMatthew G. Knepley + dm - The DMPlex object 203c99f2573SMatthew G. Knepley - perm - The point permutation, perm[old point number] = new point number 2048ed5f475SMatthew G. Knepley 2058ed5f475SMatthew G. Knepley Output Parameter: 2068ed5f475SMatthew G. Knepley . pdm - The permuted DM 2078ed5f475SMatthew G. Knepley 2088ed5f475SMatthew G. Knepley Level: intermediate 2098ed5f475SMatthew G. Knepley 210db781477SPatrick Sanan .seealso: `MatPermute()` 2118ed5f475SMatthew G. Knepley @*/ 2128ed5f475SMatthew G. Knepley PetscErrorCode DMPlexPermute(DM dm, IS perm, DM *pdm) 2138ed5f475SMatthew G. Knepley { 2148ed5f475SMatthew G. Knepley DM_Plex *plex = (DM_Plex *) dm->data, *plexNew; 2152f300c92SMatthew Knepley PetscInt dim, cdim; 2166b74b800SMatthew G. Knepley const char *name; 2178ed5f475SMatthew G. Knepley 2188ed5f475SMatthew G. Knepley PetscFunctionBegin; 2198ed5f475SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2208ed5f475SMatthew G. Knepley PetscValidHeaderSpecific(perm, IS_CLASSID, 2); 2218ed5f475SMatthew G. Knepley PetscValidPointer(pdm, 3); 2229566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject) dm), pdm)); 2239566063dSJacob Faibussowitsch PetscCall(DMSetType(*pdm, DMPLEX)); 2249566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject) dm, &name)); 2259566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) *pdm, name)); 2269566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2279566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*pdm, dim)); 2289566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cdim)); 2299566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateDim(*pdm, cdim)); 2309566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, *pdm)); 2316b74b800SMatthew G. Knepley if (dm->localSection) { 2326b74b800SMatthew G. Knepley PetscSection section, sectionNew; 2336b74b800SMatthew G. Knepley 2349566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 2359566063dSJacob Faibussowitsch PetscCall(PetscSectionPermute(section, perm, §ionNew)); 2369566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(*pdm, sectionNew)); 2379566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(§ionNew)); 238dd742cf1SMatthew G. Knepley } 2398ed5f475SMatthew G. Knepley plexNew = (DM_Plex *) (*pdm)->data; 2408ed5f475SMatthew G. Knepley /* Ignore ltogmap, ltogmapb */ 2411bb6d2a8SBarry Smith /* Ignore sf, sectionSF */ 2428ed5f475SMatthew G. Knepley /* Ignore globalVertexNumbers, globalCellNumbers */ 2438ed5f475SMatthew G. Knepley /* Reorder labels */ 2448ed5f475SMatthew G. Knepley { 245a85475f2SMatthew G. Knepley PetscInt numLabels, l; 246a85475f2SMatthew G. Knepley DMLabel label, labelNew; 2478ed5f475SMatthew G. Knepley 2489566063dSJacob Faibussowitsch PetscCall(DMGetNumLabels(dm, &numLabels)); 2495d80c0bfSVaclav Hapla for (l = 0; l < numLabels; ++l) { 2509566063dSJacob Faibussowitsch PetscCall(DMGetLabelByNum(dm, l, &label)); 2519566063dSJacob Faibussowitsch PetscCall(DMLabelPermute(label, perm, &labelNew)); 2529566063dSJacob Faibussowitsch PetscCall(DMAddLabel(*pdm, labelNew)); 2539566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&labelNew)); 2548ed5f475SMatthew G. Knepley } 2559566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*pdm, "depth", &(*pdm)->depthLabel)); 2569566063dSJacob Faibussowitsch if (plex->subpointMap) PetscCall(DMLabelPermute(plex->subpointMap, perm, &plexNew->subpointMap)); 2578ed5f475SMatthew G. Knepley } 2588ed5f475SMatthew G. Knepley /* Reorder topology */ 2598ed5f475SMatthew G. Knepley { 2608ed5f475SMatthew G. Knepley const PetscInt *pperm; 2616302a7fbSVaclav Hapla PetscInt n, pStart, pEnd, p; 2628ed5f475SMatthew G. Knepley 2639566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&plexNew->coneSection)); 2649566063dSJacob Faibussowitsch PetscCall(PetscSectionPermute(plex->coneSection, perm, &plexNew->coneSection)); 2659566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(plexNew->coneSection, &n)); 2669566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &plexNew->cones)); 2679566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &plexNew->coneOrientations)); 2689566063dSJacob Faibussowitsch PetscCall(ISGetIndices(perm, &pperm)); 2699566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(plex->coneSection, &pStart, &pEnd)); 2708ed5f475SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 2718ed5f475SMatthew G. Knepley PetscInt dof, off, offNew, d; 2728ed5f475SMatthew G. Knepley 2739566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(plexNew->coneSection, pperm[p], &dof)); 2749566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(plex->coneSection, p, &off)); 2759566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(plexNew->coneSection, pperm[p], &offNew)); 2768ed5f475SMatthew G. Knepley for (d = 0; d < dof; ++d) { 2778ed5f475SMatthew G. Knepley plexNew->cones[offNew+d] = pperm[plex->cones[off+d]]; 2788ed5f475SMatthew G. Knepley plexNew->coneOrientations[offNew+d] = plex->coneOrientations[off+d]; 2798ed5f475SMatthew G. Knepley } 2808ed5f475SMatthew G. Knepley } 2819566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&plexNew->supportSection)); 2829566063dSJacob Faibussowitsch PetscCall(PetscSectionPermute(plex->supportSection, perm, &plexNew->supportSection)); 2839566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(plexNew->supportSection, &n)); 2849566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &plexNew->supports)); 2859566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(plex->supportSection, &pStart, &pEnd)); 2868ed5f475SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 2878ed5f475SMatthew G. Knepley PetscInt dof, off, offNew, d; 2888ed5f475SMatthew G. Knepley 2899566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(plexNew->supportSection, pperm[p], &dof)); 2909566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(plex->supportSection, p, &off)); 2919566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(plexNew->supportSection, pperm[p], &offNew)); 2928ed5f475SMatthew G. Knepley for (d = 0; d < dof; ++d) { 2938ed5f475SMatthew G. Knepley plexNew->supports[offNew+d] = pperm[plex->supports[off+d]]; 2948ed5f475SMatthew G. Knepley } 2958ed5f475SMatthew G. Knepley } 2969566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(perm, &pperm)); 2978ed5f475SMatthew G. Knepley } 298a6e0b375SMatthew G. Knepley /* Remap coordinates */ 299a6e0b375SMatthew G. Knepley { 300a6e0b375SMatthew G. Knepley DM cdm, cdmNew; 301a6e0b375SMatthew G. Knepley PetscSection csection, csectionNew; 302a6e0b375SMatthew G. Knepley Vec coordinates, coordinatesNew; 303a6e0b375SMatthew G. Knepley PetscScalar *coords, *coordsNew; 304a6e0b375SMatthew G. Knepley const PetscInt *pperm; 305a6e0b375SMatthew G. Knepley PetscInt pStart, pEnd, p; 306a6e0b375SMatthew G. Knepley const char *name; 307a6e0b375SMatthew G. Knepley 3089566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 3099566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(cdm, &csection)); 3109566063dSJacob Faibussowitsch PetscCall(PetscSectionPermute(csection, perm, &csectionNew)); 3119566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 3129566063dSJacob Faibussowitsch PetscCall(VecDuplicate(coordinates, &coordinatesNew)); 3139566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)coordinates,&name)); 3149566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)coordinatesNew,name)); 3159566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords)); 3169566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinatesNew, &coordsNew)); 3179566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(csectionNew, &pStart, &pEnd)); 3189566063dSJacob Faibussowitsch PetscCall(ISGetIndices(perm, &pperm)); 319a6e0b375SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 320a6e0b375SMatthew G. Knepley PetscInt dof, off, offNew, d; 321a6e0b375SMatthew G. Knepley 3229566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(csectionNew, p, &dof)); 3239566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(csection, p, &off)); 3249566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(csectionNew, pperm[p], &offNew)); 325a6e0b375SMatthew G. Knepley for (d = 0; d < dof; ++d) coordsNew[offNew+d] = coords[off+d]; 326a6e0b375SMatthew G. Knepley } 3279566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(perm, &pperm)); 3289566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 3299566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinatesNew, &coordsNew)); 3309566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(*pdm, &cdmNew)); 3319566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(cdmNew, csectionNew)); 3329566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(*pdm, coordinatesNew)); 3339566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&csectionNew)); 3349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coordinatesNew)); 335a6e0b375SMatthew G. Knepley } 336*5de52c6dSVaclav Hapla PetscCall(DMPlexCopy_Internal(dm, PETSC_TRUE, PETSC_TRUE, *pdm)); 337e66f2fa0SMatthew G. Knepley (*pdm)->setupcalled = PETSC_TRUE; 3388ed5f475SMatthew G. Knepley PetscFunctionReturn(0); 3398ed5f475SMatthew G. Knepley } 340