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; 95*50be1cf5SMatthew G. Knepley const PetscInt *valuesTmp; 96*50be1cf5SMatthew G. Knepley PetscInt *values; 97c99f2573SMatthew G. Knepley PetscInt numValues, numPoints = 0; 98c99f2573SMatthew G. Knepley PetscInt *sperm, *vsize, *voff, v; 99c99f2573SMatthew G. Knepley 100*50be1cf5SMatthew G. Knepley // Can't directly sort the valueIS, since it is a view into the DMLabel 1019566063dSJacob Faibussowitsch PetscCall(DMLabelGetValueIS(label, &valueIS)); 1029566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(valueIS, &numValues)); 103*50be1cf5SMatthew G. Knepley PetscCall(ISGetIndices(valueIS, &valuesTmp)); 104*50be1cf5SMatthew G. Knepley PetscCall(PetscCalloc4(numCells, &sperm, numValues, &values, numValues, &vsize, numValues+1,&voff)); 105*50be1cf5SMatthew G. Knepley PetscCall(PetscArraycpy(values, valuesTmp, numValues)); 106*50be1cf5SMatthew G. Knepley PetscCall(PetscSortInt(numValues, values)); 107*50be1cf5SMatthew G. Knepley PetscCall(ISRestoreIndices(valueIS, &valuesTmp)); 108*50be1cf5SMatthew G. Knepley PetscCall(ISDestroy(&valueIS)); 109c99f2573SMatthew G. Knepley for (v = 0; v < numValues; ++v) { 1109566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumSize(label, values[v], &vsize[v])); 111c99f2573SMatthew G. Knepley if (v < numValues-1) voff[v+2] += vsize[v] + voff[v+1]; 112c99f2573SMatthew G. Knepley numPoints += vsize[v]; 113c99f2573SMatthew G. Knepley } 1146469bafaSMatthew G. Knepley PetscCheck(numPoints == numCells, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label only covers %" PetscInt_FMT " cells != %" PetscInt_FMT " total cells", numPoints, numCells); 115c99f2573SMatthew G. Knepley for (c = 0; c < numCells; ++c) { 116c99f2573SMatthew G. Knepley const PetscInt oldc = cperm[c]; 117c99f2573SMatthew G. Knepley PetscInt val, vloc; 118c99f2573SMatthew G. Knepley 1199566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(label, oldc, &val)); 12063a3b9bcSJacob Faibussowitsch PetscCheck(val != -1,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cell %" PetscInt_FMT " not present in label", oldc); 1219566063dSJacob Faibussowitsch PetscCall(PetscFindInt(val, numValues, values, &vloc)); 12263a3b9bcSJacob Faibussowitsch PetscCheck(vloc >= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Value %" PetscInt_FMT " not present label", val); 123c99f2573SMatthew G. Knepley sperm[voff[vloc+1]++] = oldc; 124c99f2573SMatthew G. Knepley } 125c99f2573SMatthew G. Knepley for (v = 0; v < numValues; ++v) { 1261dca8a05SBarry 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]); 127c99f2573SMatthew G. Knepley } 1289566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(cperm, sperm, numCells)); 129*50be1cf5SMatthew G. Knepley PetscCall(PetscFree4(sperm, values, vsize, voff)); 130c99f2573SMatthew G. Knepley } 1318ed5f475SMatthew G. Knepley /* Construct closure */ 1329566063dSJacob Faibussowitsch PetscCall(DMPlexCreateOrderingClosure_Static(dm, numCells, cperm, &clperm, &invclperm)); 1339566063dSJacob Faibussowitsch PetscCall(PetscFree3(cperm,mask,xls)); 1349566063dSJacob Faibussowitsch PetscCall(PetscFree(clperm)); 1358ed5f475SMatthew G. Knepley /* Invert permutation */ 1369566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 1379566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject) dm), pEnd-pStart, invclperm, PETSC_OWN_POINTER, perm)); 1388ed5f475SMatthew G. Knepley PetscFunctionReturn(0); 1398ed5f475SMatthew G. Knepley } 1408ed5f475SMatthew G. Knepley 1418ed5f475SMatthew G. Knepley /*@ 1426913077dSMatthew G. Knepley DMPlexGetOrdering1D - Reorder the vertices so that the mesh is in a line 1436913077dSMatthew G. Knepley 1446913077dSMatthew G. Knepley Collective on dm 1456913077dSMatthew G. Knepley 1466913077dSMatthew G. Knepley Input Parameter: 1476913077dSMatthew G. Knepley . dm - The DMPlex object 1486913077dSMatthew G. Knepley 1496913077dSMatthew G. Knepley Output Parameter: 1506913077dSMatthew G. Knepley . perm - The point permutation as an IS, perm[old point number] = new point number 1516913077dSMatthew G. Knepley 1526913077dSMatthew G. Knepley Level: intermediate 1536913077dSMatthew G. Knepley 154db781477SPatrick Sanan .seealso: `DMPlexGetOrdering()`, `DMPlexPermute()`, `MatGetOrdering()` 1556913077dSMatthew G. Knepley @*/ 1566913077dSMatthew G. Knepley PetscErrorCode DMPlexGetOrdering1D(DM dm, IS *perm) 1576913077dSMatthew G. Knepley { 1586913077dSMatthew G. Knepley PetscInt *points; 1596913077dSMatthew G. Knepley const PetscInt *support, *cone; 1606913077dSMatthew G. Knepley PetscInt dim, pStart, pEnd, cStart, cEnd, c, vStart, vEnd, v, suppSize, lastCell = 0; 1616913077dSMatthew G. Knepley 1626913077dSMatthew G. Knepley PetscFunctionBegin; 1639566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1646913077dSMatthew G. Knepley PetscCheck(dim == 1, PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Input mesh must be one dimensional, not %" PetscInt_FMT, dim); 1659566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 1669566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 1679566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 1689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(pEnd-pStart, &points)); 1696913077dSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) points[c] = c; 1706913077dSMatthew G. Knepley for (v = vStart; v < vEnd; ++v) points[v] = v; 1716913077dSMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 1729566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, v, &suppSize)); 1739566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(dm, v, &support)); 1746913077dSMatthew G. Knepley if (suppSize == 1) {lastCell = support[0]; break;} 1756913077dSMatthew G. Knepley } 1766913077dSMatthew G. Knepley if (v < vEnd) { 1776913077dSMatthew G. Knepley PetscInt pos = cEnd; 1786913077dSMatthew G. Knepley 1796913077dSMatthew G. Knepley points[v] = pos++; 1806913077dSMatthew G. Knepley while (lastCell >= cStart) { 1819566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, lastCell, &cone)); 1826913077dSMatthew G. Knepley if (cone[0] == v) v = cone[1]; 1836913077dSMatthew G. Knepley else v = cone[0]; 1849566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(dm, v, &support)); 1859566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, v, &suppSize)); 1866913077dSMatthew G. Knepley if (suppSize == 1) {lastCell = -1;} 1876913077dSMatthew G. Knepley else { 1886913077dSMatthew G. Knepley if (support[0] == lastCell) lastCell = support[1]; 1896913077dSMatthew G. Knepley else lastCell = support[0]; 1906913077dSMatthew G. Knepley } 1916913077dSMatthew G. Knepley points[v] = pos++; 1926913077dSMatthew G. Knepley } 1936913077dSMatthew G. Knepley PetscCheck(pos == pEnd, PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Last vertex was %" PetscInt_FMT ", not %" PetscInt_FMT, pos, pEnd); 1946913077dSMatthew G. Knepley } 1959566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject) dm), pEnd-pStart, points, PETSC_OWN_POINTER, perm)); 1966913077dSMatthew G. Knepley PetscFunctionReturn(0); 1976913077dSMatthew G. Knepley } 1986913077dSMatthew G. Knepley 1996858538eSMatthew G. Knepley static PetscErrorCode DMPlexRemapCoordinates_Private(IS perm, PetscSection cs, Vec coordinates, PetscSection *csNew, Vec *coordinatesNew) 2006858538eSMatthew G. Knepley { 2016858538eSMatthew G. Knepley PetscScalar *coords, *coordsNew; 2026858538eSMatthew G. Knepley const PetscInt *pperm; 2036858538eSMatthew G. Knepley PetscInt pStart, pEnd, p; 2046858538eSMatthew G. Knepley const char *name; 2056858538eSMatthew G. Knepley 2066858538eSMatthew G. Knepley PetscFunctionBegin; 2076858538eSMatthew G. Knepley PetscCall(PetscSectionPermute(cs, perm, csNew)); 2086858538eSMatthew G. Knepley PetscCall(VecDuplicate(coordinates, coordinatesNew)); 2096858538eSMatthew G. Knepley PetscCall(PetscObjectGetName((PetscObject) coordinates, &name)); 2106858538eSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject) *coordinatesNew, name)); 2116858538eSMatthew G. Knepley PetscCall(VecGetArray(coordinates, &coords)); 2126858538eSMatthew G. Knepley PetscCall(VecGetArray(*coordinatesNew, &coordsNew)); 2136858538eSMatthew G. Knepley PetscCall(PetscSectionGetChart(*csNew, &pStart, &pEnd)); 2146858538eSMatthew G. Knepley PetscCall(ISGetIndices(perm, &pperm)); 2156858538eSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 2166858538eSMatthew G. Knepley PetscInt dof, off, offNew, d; 2176858538eSMatthew G. Knepley 2186858538eSMatthew G. Knepley PetscCall(PetscSectionGetDof(*csNew, p, &dof)); 2196858538eSMatthew G. Knepley PetscCall(PetscSectionGetOffset(cs, p, &off)); 2206858538eSMatthew G. Knepley PetscCall(PetscSectionGetOffset(*csNew, pperm[p], &offNew)); 2216858538eSMatthew G. Knepley for (d = 0; d < dof; ++d) coordsNew[offNew+d] = coords[off+d]; 2226858538eSMatthew G. Knepley } 2236858538eSMatthew G. Knepley PetscCall(ISRestoreIndices(perm, &pperm)); 2246858538eSMatthew G. Knepley PetscCall(VecRestoreArray(coordinates, &coords)); 2256858538eSMatthew G. Knepley PetscCall(VecRestoreArray(*coordinatesNew, &coordsNew)); 2266858538eSMatthew G. Knepley PetscFunctionReturn(0); 2276858538eSMatthew G. Knepley } 2286858538eSMatthew G. Knepley 2296913077dSMatthew G. Knepley /*@ 2308ed5f475SMatthew G. Knepley DMPlexPermute - Reorder the mesh according to the input permutation 2318ed5f475SMatthew G. Knepley 232d083f849SBarry Smith Collective on dm 2338ed5f475SMatthew G. Knepley 234d8d19677SJose E. Roman Input Parameters: 2358ed5f475SMatthew G. Knepley + dm - The DMPlex object 236c99f2573SMatthew G. Knepley - perm - The point permutation, perm[old point number] = new point number 2378ed5f475SMatthew G. Knepley 2388ed5f475SMatthew G. Knepley Output Parameter: 2398ed5f475SMatthew G. Knepley . pdm - The permuted DM 2408ed5f475SMatthew G. Knepley 2418ed5f475SMatthew G. Knepley Level: intermediate 2428ed5f475SMatthew G. Knepley 243db781477SPatrick Sanan .seealso: `MatPermute()` 2448ed5f475SMatthew G. Knepley @*/ 2458ed5f475SMatthew G. Knepley PetscErrorCode DMPlexPermute(DM dm, IS perm, DM *pdm) 2468ed5f475SMatthew G. Knepley { 2478ed5f475SMatthew G. Knepley DM_Plex *plex = (DM_Plex *) dm->data, *plexNew; 2482f300c92SMatthew Knepley PetscInt dim, cdim; 2496b74b800SMatthew G. Knepley const char *name; 2508ed5f475SMatthew G. Knepley 2518ed5f475SMatthew G. Knepley PetscFunctionBegin; 2528ed5f475SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2538ed5f475SMatthew G. Knepley PetscValidHeaderSpecific(perm, IS_CLASSID, 2); 2548ed5f475SMatthew G. Knepley PetscValidPointer(pdm, 3); 2559566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject) dm), pdm)); 2569566063dSJacob Faibussowitsch PetscCall(DMSetType(*pdm, DMPLEX)); 2579566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject) dm, &name)); 2589566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) *pdm, name)); 2599566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2609566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*pdm, dim)); 2619566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cdim)); 2629566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateDim(*pdm, cdim)); 2639566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, *pdm)); 2646b74b800SMatthew G. Knepley if (dm->localSection) { 2656b74b800SMatthew G. Knepley PetscSection section, sectionNew; 2666b74b800SMatthew G. Knepley 2679566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 2689566063dSJacob Faibussowitsch PetscCall(PetscSectionPermute(section, perm, §ionNew)); 2699566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(*pdm, sectionNew)); 2709566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(§ionNew)); 271dd742cf1SMatthew G. Knepley } 2728ed5f475SMatthew G. Knepley plexNew = (DM_Plex *) (*pdm)->data; 2738ed5f475SMatthew G. Knepley /* Ignore ltogmap, ltogmapb */ 2741bb6d2a8SBarry Smith /* Ignore sf, sectionSF */ 2758ed5f475SMatthew G. Knepley /* Ignore globalVertexNumbers, globalCellNumbers */ 2768ed5f475SMatthew G. Knepley /* Reorder labels */ 2778ed5f475SMatthew G. Knepley { 278a85475f2SMatthew G. Knepley PetscInt numLabels, l; 279a85475f2SMatthew G. Knepley DMLabel label, labelNew; 2808ed5f475SMatthew G. Knepley 2819566063dSJacob Faibussowitsch PetscCall(DMGetNumLabels(dm, &numLabels)); 2825d80c0bfSVaclav Hapla for (l = 0; l < numLabels; ++l) { 2839566063dSJacob Faibussowitsch PetscCall(DMGetLabelByNum(dm, l, &label)); 2849566063dSJacob Faibussowitsch PetscCall(DMLabelPermute(label, perm, &labelNew)); 2859566063dSJacob Faibussowitsch PetscCall(DMAddLabel(*pdm, labelNew)); 2869566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&labelNew)); 2878ed5f475SMatthew G. Knepley } 2889566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*pdm, "depth", &(*pdm)->depthLabel)); 2899566063dSJacob Faibussowitsch if (plex->subpointMap) PetscCall(DMLabelPermute(plex->subpointMap, perm, &plexNew->subpointMap)); 2908ed5f475SMatthew G. Knepley } 2918ed5f475SMatthew G. Knepley /* Reorder topology */ 2928ed5f475SMatthew G. Knepley { 2938ed5f475SMatthew G. Knepley const PetscInt *pperm; 2946302a7fbSVaclav Hapla PetscInt n, pStart, pEnd, p; 2958ed5f475SMatthew G. Knepley 2969566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&plexNew->coneSection)); 2979566063dSJacob Faibussowitsch PetscCall(PetscSectionPermute(plex->coneSection, perm, &plexNew->coneSection)); 2989566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(plexNew->coneSection, &n)); 2999566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &plexNew->cones)); 3009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &plexNew->coneOrientations)); 3019566063dSJacob Faibussowitsch PetscCall(ISGetIndices(perm, &pperm)); 3029566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(plex->coneSection, &pStart, &pEnd)); 3038ed5f475SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 3048ed5f475SMatthew G. Knepley PetscInt dof, off, offNew, d; 3058ed5f475SMatthew G. Knepley 3069566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(plexNew->coneSection, pperm[p], &dof)); 3079566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(plex->coneSection, p, &off)); 3089566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(plexNew->coneSection, pperm[p], &offNew)); 3098ed5f475SMatthew G. Knepley for (d = 0; d < dof; ++d) { 3108ed5f475SMatthew G. Knepley plexNew->cones[offNew+d] = pperm[plex->cones[off+d]]; 3118ed5f475SMatthew G. Knepley plexNew->coneOrientations[offNew+d] = plex->coneOrientations[off+d]; 3128ed5f475SMatthew G. Knepley } 3138ed5f475SMatthew G. Knepley } 3149566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&plexNew->supportSection)); 3159566063dSJacob Faibussowitsch PetscCall(PetscSectionPermute(plex->supportSection, perm, &plexNew->supportSection)); 3169566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(plexNew->supportSection, &n)); 3179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &plexNew->supports)); 3189566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(plex->supportSection, &pStart, &pEnd)); 3198ed5f475SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 3208ed5f475SMatthew G. Knepley PetscInt dof, off, offNew, d; 3218ed5f475SMatthew G. Knepley 3229566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(plexNew->supportSection, pperm[p], &dof)); 3239566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(plex->supportSection, p, &off)); 3249566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(plexNew->supportSection, pperm[p], &offNew)); 3258ed5f475SMatthew G. Knepley for (d = 0; d < dof; ++d) { 3268ed5f475SMatthew G. Knepley plexNew->supports[offNew+d] = pperm[plex->supports[off+d]]; 3278ed5f475SMatthew G. Knepley } 3288ed5f475SMatthew G. Knepley } 3299566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(perm, &pperm)); 3308ed5f475SMatthew G. Knepley } 331a6e0b375SMatthew G. Knepley /* Remap coordinates */ 332a6e0b375SMatthew G. Knepley { 333a6e0b375SMatthew G. Knepley DM cdm, cdmNew; 3346858538eSMatthew G. Knepley PetscSection cs, csNew; 335a6e0b375SMatthew G. Knepley Vec coordinates, coordinatesNew; 336a6e0b375SMatthew G. Knepley 3376858538eSMatthew G. Knepley PetscCall(DMGetCoordinateSection(dm, &cs)); 3389566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 3396858538eSMatthew G. Knepley PetscCall(DMPlexRemapCoordinates_Private(perm, cs, coordinates, &csNew, &coordinatesNew)); 3406858538eSMatthew G. Knepley PetscCall(DMSetCoordinateSection(*pdm, PETSC_DETERMINE, csNew)); 3419566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(*pdm, coordinatesNew)); 3426858538eSMatthew G. Knepley PetscCall(PetscSectionDestroy(&csNew)); 3439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coordinatesNew)); 3446858538eSMatthew G. Knepley 3456858538eSMatthew G. Knepley PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 3466858538eSMatthew G. Knepley if (cdm) { 3476858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDM(*pdm, &cdm)); 3486858538eSMatthew G. Knepley PetscCall(DMClone(cdm, &cdmNew)); 3496858538eSMatthew G. Knepley PetscCall(DMSetCellCoordinateDM(*pdm, cdmNew)); 3506858538eSMatthew G. Knepley PetscCall(DMDestroy(&cdmNew)); 3516858538eSMatthew G. Knepley PetscCall(DMGetCellCoordinateSection(dm, &cs)); 3526858538eSMatthew G. Knepley PetscCall(DMGetCellCoordinatesLocal(dm, &coordinates)); 3536858538eSMatthew G. Knepley PetscCall(DMPlexRemapCoordinates_Private(perm, cs, coordinates, &csNew, &coordinatesNew)); 3546858538eSMatthew G. Knepley PetscCall(DMSetCellCoordinateSection(*pdm, PETSC_DETERMINE, csNew)); 3556858538eSMatthew G. Knepley PetscCall(DMSetCellCoordinatesLocal(*pdm, coordinatesNew)); 3566858538eSMatthew G. Knepley PetscCall(PetscSectionDestroy(&csNew)); 3576858538eSMatthew G. Knepley PetscCall(VecDestroy(&coordinatesNew)); 3586858538eSMatthew G. Knepley } 359a6e0b375SMatthew G. Knepley } 3605de52c6dSVaclav Hapla PetscCall(DMPlexCopy_Internal(dm, PETSC_TRUE, PETSC_TRUE, *pdm)); 361e66f2fa0SMatthew G. Knepley (*pdm)->setupcalled = PETSC_TRUE; 3628ed5f475SMatthew G. Knepley PetscFunctionReturn(0); 3638ed5f475SMatthew G. Knepley } 3646bc1bd01Sksagiyam 3656bc1bd01Sksagiyam PetscErrorCode DMPlexReorderSetDefault_Plex(DM dm, DMPlexReorderDefaultFlag reorder) 3666bc1bd01Sksagiyam { 3676bc1bd01Sksagiyam DM_Plex *mesh = (DM_Plex *) dm->data; 3686bc1bd01Sksagiyam 3696bc1bd01Sksagiyam PetscFunctionBegin; 3706bc1bd01Sksagiyam mesh->reorderDefault = reorder; 3716bc1bd01Sksagiyam PetscFunctionReturn(0); 3726bc1bd01Sksagiyam } 3736bc1bd01Sksagiyam 3746bc1bd01Sksagiyam /*@ 3756bc1bd01Sksagiyam DMPlexReorderSetDefault - Set flag indicating whether the DM should be reordered by default 3766bc1bd01Sksagiyam 3776bc1bd01Sksagiyam Logically collective 3786bc1bd01Sksagiyam 3796bc1bd01Sksagiyam Input Parameters: 3806bc1bd01Sksagiyam + dm - The DM 3816bc1bd01Sksagiyam - reorder - Flag for reordering 3826bc1bd01Sksagiyam 3836bc1bd01Sksagiyam Level: intermediate 3846bc1bd01Sksagiyam 3856bc1bd01Sksagiyam .seealso: `DMPlexReorderGetDefault()` 3866bc1bd01Sksagiyam @*/ 3876bc1bd01Sksagiyam PetscErrorCode DMPlexReorderSetDefault(DM dm, DMPlexReorderDefaultFlag reorder) 3886bc1bd01Sksagiyam { 3896bc1bd01Sksagiyam PetscFunctionBegin; 3906bc1bd01Sksagiyam PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3916bc1bd01Sksagiyam PetscTryMethod(dm,"DMPlexReorderSetDefault_C",(DM,DMPlexReorderDefaultFlag),(dm,reorder)); 3926bc1bd01Sksagiyam PetscFunctionReturn(0); 3936bc1bd01Sksagiyam } 3946bc1bd01Sksagiyam 3956bc1bd01Sksagiyam PetscErrorCode DMPlexReorderGetDefault_Plex(DM dm, DMPlexReorderDefaultFlag *reorder) 3966bc1bd01Sksagiyam { 3976bc1bd01Sksagiyam DM_Plex *mesh = (DM_Plex *) dm->data; 3986bc1bd01Sksagiyam 3996bc1bd01Sksagiyam PetscFunctionBegin; 4006bc1bd01Sksagiyam *reorder = mesh->reorderDefault; 4016bc1bd01Sksagiyam PetscFunctionReturn(0); 4026bc1bd01Sksagiyam } 4036bc1bd01Sksagiyam 4046bc1bd01Sksagiyam /*@ 4056bc1bd01Sksagiyam DMPlexReorderGetDefault - Get flag indicating whether the DM should be reordered by default 4066bc1bd01Sksagiyam 4076bc1bd01Sksagiyam Not collective 4086bc1bd01Sksagiyam 4096bc1bd01Sksagiyam Input Parameter: 4106bc1bd01Sksagiyam . dm - The DM 4116bc1bd01Sksagiyam 4126bc1bd01Sksagiyam Output Parameter: 4136bc1bd01Sksagiyam . reorder - Flag for reordering 4146bc1bd01Sksagiyam 4156bc1bd01Sksagiyam Level: intermediate 4166bc1bd01Sksagiyam 4176bc1bd01Sksagiyam .seealso: `DMPlexReorderSetDefault()` 4186bc1bd01Sksagiyam @*/ 4196bc1bd01Sksagiyam PetscErrorCode DMPlexReorderGetDefault(DM dm, DMPlexReorderDefaultFlag *reorder) 4206bc1bd01Sksagiyam { 4216bc1bd01Sksagiyam PetscFunctionBegin; 4226bc1bd01Sksagiyam PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4236bc1bd01Sksagiyam PetscValidPointer(reorder, 2); 4246bc1bd01Sksagiyam PetscUseMethod(dm,"DMPlexReorderGetDefault_C",(DM,DMPlexReorderDefaultFlag*),(dm,reorder)); 4256bc1bd01Sksagiyam PetscFunctionReturn(0); 4266bc1bd01Sksagiyam } 427