xref: /petsc/src/dm/impls/plex/plexreorder.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2af0996ceSBarry Smith #include <petsc/private/matorderimpl.h> /*I      "petscmat.h"      I*/
3695799ffSMatthew G. Knepley #include <petsc/private/dmlabelimpl.h>
48ed5f475SMatthew G. Knepley 
5d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCreateOrderingClosure_Static(DM dm, PetscInt numPoints, const PetscInt pperm[], PetscInt **clperm, PetscInt **invclperm)
6d71ae5a4SJacob Faibussowitsch {
78ed5f475SMatthew G. Knepley   PetscInt *perm, *iperm;
88ed5f475SMatthew G. Knepley   PetscInt  depth, d, pStart, pEnd, fStart, fMax, fEnd, p;
98ed5f475SMatthew G. Knepley 
108ed5f475SMatthew G. Knepley   PetscFunctionBegin;
119566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
129566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
139566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(pEnd - pStart, &perm));
149566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(pEnd - pStart, &iperm));
158ed5f475SMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) iperm[p] = -1;
168ed5f475SMatthew G. Knepley   for (d = depth; d > 0; --d) {
179566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, d, &pStart, &pEnd));
189566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, d - 1, &fStart, &fEnd));
198ed5f475SMatthew G. Knepley     fMax = fStart;
208ed5f475SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
218ed5f475SMatthew G. Knepley       const PetscInt *cone;
228ed5f475SMatthew G. Knepley       PetscInt        point, coneSize, c;
238ed5f475SMatthew G. Knepley 
248ed5f475SMatthew G. Knepley       if (d == depth) {
258ed5f475SMatthew G. Knepley         perm[p]         = pperm[p];
268ed5f475SMatthew G. Knepley         iperm[pperm[p]] = p;
278ed5f475SMatthew G. Knepley       }
288ed5f475SMatthew G. Knepley       point = perm[p];
299566063dSJacob Faibussowitsch       PetscCall(DMPlexGetConeSize(dm, point, &coneSize));
309566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(dm, point, &cone));
318ed5f475SMatthew G. Knepley       for (c = 0; c < coneSize; ++c) {
328ed5f475SMatthew G. Knepley         const PetscInt oldc = cone[c];
338ed5f475SMatthew G. Knepley         const PetscInt newc = iperm[oldc];
348ed5f475SMatthew G. Knepley 
358ed5f475SMatthew G. Knepley         if (newc < 0) {
368ed5f475SMatthew G. Knepley           perm[fMax]  = oldc;
378ed5f475SMatthew G. Knepley           iperm[oldc] = fMax++;
388ed5f475SMatthew G. Knepley         }
398ed5f475SMatthew G. Knepley       }
408ed5f475SMatthew G. Knepley     }
4163a3b9bcSJacob 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);
428ed5f475SMatthew G. Knepley   }
438ed5f475SMatthew G. Knepley   *clperm    = perm;
448ed5f475SMatthew G. Knepley   *invclperm = iperm;
45*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
468ed5f475SMatthew G. Knepley }
478ed5f475SMatthew G. Knepley 
488ed5f475SMatthew G. Knepley /*@
498ed5f475SMatthew G. Knepley   DMPlexGetOrdering - Calculate a reordering of the mesh
508ed5f475SMatthew G. Knepley 
51d083f849SBarry Smith   Collective on dm
528ed5f475SMatthew G. Knepley 
53d8d19677SJose E. Roman   Input Parameters:
548ed5f475SMatthew G. Knepley + dm - The DMPlex object
55c99f2573SMatthew G. Knepley . otype - type of reordering, one of the following:
568ed5f475SMatthew G. Knepley $     MATORDERINGNATURAL - Natural
578ed5f475SMatthew G. Knepley $     MATORDERINGND - Nested Dissection
588ed5f475SMatthew G. Knepley $     MATORDERING1WD - One-way Dissection
598ed5f475SMatthew G. Knepley $     MATORDERINGRCM - Reverse Cuthill-McKee
608ed5f475SMatthew G. Knepley $     MATORDERINGQMD - Quotient Minimum Degree
61c99f2573SMatthew G. Knepley - label - [Optional] Label used to segregate ordering into sets, or NULL
628ed5f475SMatthew G. Knepley 
638ed5f475SMatthew G. Knepley   Output Parameter:
64c99f2573SMatthew G. Knepley . perm - The point permutation as an IS, perm[old point number] = new point number
65c99f2573SMatthew G. Knepley 
66c99f2573SMatthew 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
67c99f2573SMatthew G. Knepley   has different types of cells, and then loop over each set of reordered cells for assembly.
688ed5f475SMatthew G. Knepley 
698ed5f475SMatthew G. Knepley   Level: intermediate
708ed5f475SMatthew G. Knepley 
71db781477SPatrick Sanan .seealso: `DMPlexPermute()`, `MatGetOrdering()`
728ed5f475SMatthew G. Knepley @*/
73d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetOrdering(DM dm, MatOrderingType otype, DMLabel label, IS *perm)
74d71ae5a4SJacob Faibussowitsch {
758ed5f475SMatthew G. Knepley   PetscInt  numCells = 0;
763c67d5adSSatish Balay   PetscInt *start = NULL, *adjacency = NULL, *cperm, *clperm = NULL, *invclperm = NULL, *mask, *xls, pStart, pEnd, c, i;
778ed5f475SMatthew G. Knepley 
788ed5f475SMatthew G. Knepley   PetscFunctionBegin;
798ed5f475SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
80064a246eSJacob Faibussowitsch   PetscValidPointer(perm, 4);
819566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateNeighborCSR(dm, 0, &numCells, &start, &adjacency));
829566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(numCells, &cperm, numCells, &mask, numCells * 2, &xls));
83d80c2872SMatthew G. Knepley   if (numCells) {
848ed5f475SMatthew G. Knepley     /* Shift for Fortran numbering */
858ed5f475SMatthew G. Knepley     for (i = 0; i < start[numCells]; ++i) ++adjacency[i];
868ed5f475SMatthew G. Knepley     for (i = 0; i <= numCells; ++i) ++start[i];
879566063dSJacob Faibussowitsch     PetscCall(SPARSEPACKgenrcm(&numCells, start, adjacency, cperm, mask, xls));
88d80c2872SMatthew G. Knepley   }
899566063dSJacob Faibussowitsch   PetscCall(PetscFree(start));
909566063dSJacob Faibussowitsch   PetscCall(PetscFree(adjacency));
918ed5f475SMatthew G. Knepley   /* Shift for Fortran numbering */
928ed5f475SMatthew G. Knepley   for (c = 0; c < numCells; ++c) --cperm[c];
93c99f2573SMatthew G. Knepley   /* Segregate */
94c99f2573SMatthew G. Knepley   if (label) {
95c99f2573SMatthew G. Knepley     IS              valueIS;
9650be1cf5SMatthew G. Knepley     const PetscInt *valuesTmp;
9750be1cf5SMatthew G. Knepley     PetscInt       *values;
98c99f2573SMatthew G. Knepley     PetscInt        numValues, numPoints = 0;
99c99f2573SMatthew G. Knepley     PetscInt       *sperm, *vsize, *voff, v;
100c99f2573SMatthew G. Knepley 
10150be1cf5SMatthew G. Knepley     // Can't directly sort the valueIS, since it is a view into the DMLabel
1029566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValueIS(label, &valueIS));
1039566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(valueIS, &numValues));
10450be1cf5SMatthew G. Knepley     PetscCall(ISGetIndices(valueIS, &valuesTmp));
10550be1cf5SMatthew G. Knepley     PetscCall(PetscCalloc4(numCells, &sperm, numValues, &values, numValues, &vsize, numValues + 1, &voff));
10650be1cf5SMatthew G. Knepley     PetscCall(PetscArraycpy(values, valuesTmp, numValues));
10750be1cf5SMatthew G. Knepley     PetscCall(PetscSortInt(numValues, values));
10850be1cf5SMatthew G. Knepley     PetscCall(ISRestoreIndices(valueIS, &valuesTmp));
10950be1cf5SMatthew G. Knepley     PetscCall(ISDestroy(&valueIS));
110c99f2573SMatthew G. Knepley     for (v = 0; v < numValues; ++v) {
1119566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumSize(label, values[v], &vsize[v]));
112c99f2573SMatthew G. Knepley       if (v < numValues - 1) voff[v + 2] += vsize[v] + voff[v + 1];
113c99f2573SMatthew G. Knepley       numPoints += vsize[v];
114c99f2573SMatthew G. Knepley     }
1156469bafaSMatthew G. Knepley     PetscCheck(numPoints == numCells, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label only covers %" PetscInt_FMT " cells != %" PetscInt_FMT " total cells", numPoints, numCells);
116c99f2573SMatthew G. Knepley     for (c = 0; c < numCells; ++c) {
117c99f2573SMatthew G. Knepley       const PetscInt oldc = cperm[c];
118c99f2573SMatthew G. Knepley       PetscInt       val, vloc;
119c99f2573SMatthew G. Knepley 
1209566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValue(label, oldc, &val));
12163a3b9bcSJacob Faibussowitsch       PetscCheck(val != -1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cell %" PetscInt_FMT " not present in label", oldc);
1229566063dSJacob Faibussowitsch       PetscCall(PetscFindInt(val, numValues, values, &vloc));
12363a3b9bcSJacob Faibussowitsch       PetscCheck(vloc >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Value %" PetscInt_FMT " not present label", val);
124c99f2573SMatthew G. Knepley       sperm[voff[vloc + 1]++] = oldc;
125c99f2573SMatthew G. Knepley     }
126ad540459SPierre Jolivet     for (v = 0; v < numValues; ++v) 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]);
1279566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(cperm, sperm, numCells));
12850be1cf5SMatthew G. Knepley     PetscCall(PetscFree4(sperm, values, vsize, voff));
129c99f2573SMatthew G. Knepley   }
1308ed5f475SMatthew G. Knepley   /* Construct closure */
1319566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateOrderingClosure_Static(dm, numCells, cperm, &clperm, &invclperm));
1329566063dSJacob Faibussowitsch   PetscCall(PetscFree3(cperm, mask, xls));
1339566063dSJacob Faibussowitsch   PetscCall(PetscFree(clperm));
1348ed5f475SMatthew G. Knepley   /* Invert permutation */
1359566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
1369566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), pEnd - pStart, invclperm, PETSC_OWN_POINTER, perm));
137*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1388ed5f475SMatthew G. Knepley }
1398ed5f475SMatthew G. Knepley 
1408ed5f475SMatthew G. Knepley /*@
1416913077dSMatthew G. Knepley   DMPlexGetOrdering1D - Reorder the vertices so that the mesh is in a line
1426913077dSMatthew G. Knepley 
1436913077dSMatthew G. Knepley   Collective on dm
1446913077dSMatthew G. Knepley 
1456913077dSMatthew G. Knepley   Input Parameter:
1466913077dSMatthew G. Knepley . dm - The DMPlex object
1476913077dSMatthew G. Knepley 
1486913077dSMatthew G. Knepley   Output Parameter:
1496913077dSMatthew G. Knepley . perm - The point permutation as an IS, perm[old point number] = new point number
1506913077dSMatthew G. Knepley 
1516913077dSMatthew G. Knepley   Level: intermediate
1526913077dSMatthew G. Knepley 
153db781477SPatrick Sanan .seealso: `DMPlexGetOrdering()`, `DMPlexPermute()`, `MatGetOrdering()`
1546913077dSMatthew G. Knepley @*/
155d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetOrdering1D(DM dm, IS *perm)
156d71ae5a4SJacob Faibussowitsch {
1576913077dSMatthew G. Knepley   PetscInt       *points;
1586913077dSMatthew G. Knepley   const PetscInt *support, *cone;
1596913077dSMatthew G. Knepley   PetscInt        dim, pStart, pEnd, cStart, cEnd, c, vStart, vEnd, v, suppSize, lastCell = 0;
1606913077dSMatthew G. Knepley 
1616913077dSMatthew G. Knepley   PetscFunctionBegin;
1629566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1636913077dSMatthew G. Knepley   PetscCheck(dim == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Input mesh must be one dimensional, not %" PetscInt_FMT, dim);
1649566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
1659566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1669566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1679566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(pEnd - pStart, &points));
1686913077dSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) points[c] = c;
1696913077dSMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) points[v] = v;
1706913077dSMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
1719566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupportSize(dm, v, &suppSize));
1729566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupport(dm, v, &support));
1739371c9d4SSatish Balay     if (suppSize == 1) {
1749371c9d4SSatish Balay       lastCell = support[0];
1759371c9d4SSatish Balay       break;
1769371c9d4SSatish Balay     }
1776913077dSMatthew G. Knepley   }
1786913077dSMatthew G. Knepley   if (v < vEnd) {
1796913077dSMatthew G. Knepley     PetscInt pos = cEnd;
1806913077dSMatthew G. Knepley 
1816913077dSMatthew G. Knepley     points[v] = pos++;
1826913077dSMatthew G. Knepley     while (lastCell >= cStart) {
1839566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(dm, lastCell, &cone));
1846913077dSMatthew G. Knepley       if (cone[0] == v) v = cone[1];
1856913077dSMatthew G. Knepley       else v = cone[0];
1869566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupport(dm, v, &support));
1879566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupportSize(dm, v, &suppSize));
1889371c9d4SSatish Balay       if (suppSize == 1) {
1899371c9d4SSatish Balay         lastCell = -1;
1909371c9d4SSatish Balay       } else {
1916913077dSMatthew G. Knepley         if (support[0] == lastCell) lastCell = support[1];
1926913077dSMatthew G. Knepley         else lastCell = support[0];
1936913077dSMatthew G. Knepley       }
1946913077dSMatthew G. Knepley       points[v] = pos++;
1956913077dSMatthew G. Knepley     }
1966913077dSMatthew G. Knepley     PetscCheck(pos == pEnd, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Last vertex was %" PetscInt_FMT ", not %" PetscInt_FMT, pos, pEnd);
1976913077dSMatthew G. Knepley   }
1989566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), pEnd - pStart, points, PETSC_OWN_POINTER, perm));
199*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2006913077dSMatthew G. Knepley }
2016913077dSMatthew G. Knepley 
202d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexRemapCoordinates_Private(IS perm, PetscSection cs, Vec coordinates, PetscSection *csNew, Vec *coordinatesNew)
203d71ae5a4SJacob Faibussowitsch {
2046858538eSMatthew G. Knepley   PetscScalar    *coords, *coordsNew;
2056858538eSMatthew G. Knepley   const PetscInt *pperm;
2066858538eSMatthew G. Knepley   PetscInt        pStart, pEnd, p;
2076858538eSMatthew G. Knepley   const char     *name;
2086858538eSMatthew G. Knepley 
2096858538eSMatthew G. Knepley   PetscFunctionBegin;
2106858538eSMatthew G. Knepley   PetscCall(PetscSectionPermute(cs, perm, csNew));
2116858538eSMatthew G. Knepley   PetscCall(VecDuplicate(coordinates, coordinatesNew));
2126858538eSMatthew G. Knepley   PetscCall(PetscObjectGetName((PetscObject)coordinates, &name));
2136858538eSMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)*coordinatesNew, name));
2146858538eSMatthew G. Knepley   PetscCall(VecGetArray(coordinates, &coords));
2156858538eSMatthew G. Knepley   PetscCall(VecGetArray(*coordinatesNew, &coordsNew));
2166858538eSMatthew G. Knepley   PetscCall(PetscSectionGetChart(*csNew, &pStart, &pEnd));
2176858538eSMatthew G. Knepley   PetscCall(ISGetIndices(perm, &pperm));
2186858538eSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
2196858538eSMatthew G. Knepley     PetscInt dof, off, offNew, d;
2206858538eSMatthew G. Knepley 
2216858538eSMatthew G. Knepley     PetscCall(PetscSectionGetDof(*csNew, p, &dof));
2226858538eSMatthew G. Knepley     PetscCall(PetscSectionGetOffset(cs, p, &off));
2236858538eSMatthew G. Knepley     PetscCall(PetscSectionGetOffset(*csNew, pperm[p], &offNew));
2246858538eSMatthew G. Knepley     for (d = 0; d < dof; ++d) coordsNew[offNew + d] = coords[off + d];
2256858538eSMatthew G. Knepley   }
2266858538eSMatthew G. Knepley   PetscCall(ISRestoreIndices(perm, &pperm));
2276858538eSMatthew G. Knepley   PetscCall(VecRestoreArray(coordinates, &coords));
2286858538eSMatthew G. Knepley   PetscCall(VecRestoreArray(*coordinatesNew, &coordsNew));
229*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2306858538eSMatthew G. Knepley }
2316858538eSMatthew G. Knepley 
2326913077dSMatthew G. Knepley /*@
2338ed5f475SMatthew G. Knepley   DMPlexPermute - Reorder the mesh according to the input permutation
2348ed5f475SMatthew G. Knepley 
235d083f849SBarry Smith   Collective on dm
2368ed5f475SMatthew G. Knepley 
237d8d19677SJose E. Roman   Input Parameters:
2388ed5f475SMatthew G. Knepley + dm - The DMPlex object
239c99f2573SMatthew G. Knepley - perm - The point permutation, perm[old point number] = new point number
2408ed5f475SMatthew G. Knepley 
2418ed5f475SMatthew G. Knepley   Output Parameter:
2428ed5f475SMatthew G. Knepley . pdm - The permuted DM
2438ed5f475SMatthew G. Knepley 
2448ed5f475SMatthew G. Knepley   Level: intermediate
2458ed5f475SMatthew G. Knepley 
246db781477SPatrick Sanan .seealso: `MatPermute()`
2478ed5f475SMatthew G. Knepley @*/
248d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexPermute(DM dm, IS perm, DM *pdm)
249d71ae5a4SJacob Faibussowitsch {
2508ed5f475SMatthew G. Knepley   DM_Plex    *plex = (DM_Plex *)dm->data, *plexNew;
2512f300c92SMatthew Knepley   PetscInt    dim, cdim;
2526b74b800SMatthew G. Knepley   const char *name;
2538ed5f475SMatthew G. Knepley 
2548ed5f475SMatthew G. Knepley   PetscFunctionBegin;
2558ed5f475SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2568ed5f475SMatthew G. Knepley   PetscValidHeaderSpecific(perm, IS_CLASSID, 2);
2578ed5f475SMatthew G. Knepley   PetscValidPointer(pdm, 3);
2589566063dSJacob Faibussowitsch   PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), pdm));
2599566063dSJacob Faibussowitsch   PetscCall(DMSetType(*pdm, DMPLEX));
2609566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)dm, &name));
2619566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*pdm, name));
2629566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
2639566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(*pdm, dim));
2649566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &cdim));
2659566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinateDim(*pdm, cdim));
2669566063dSJacob Faibussowitsch   PetscCall(DMCopyDisc(dm, *pdm));
2676b74b800SMatthew G. Knepley   if (dm->localSection) {
2686b74b800SMatthew G. Knepley     PetscSection section, sectionNew;
2696b74b800SMatthew G. Knepley 
2709566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dm, &section));
2719566063dSJacob Faibussowitsch     PetscCall(PetscSectionPermute(section, perm, &sectionNew));
2729566063dSJacob Faibussowitsch     PetscCall(DMSetLocalSection(*pdm, sectionNew));
2739566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&sectionNew));
274dd742cf1SMatthew G. Knepley   }
2758ed5f475SMatthew G. Knepley   plexNew = (DM_Plex *)(*pdm)->data;
2768ed5f475SMatthew G. Knepley   /* Ignore ltogmap, ltogmapb */
2771bb6d2a8SBarry Smith   /* Ignore sf, sectionSF */
2788ed5f475SMatthew G. Knepley   /* Ignore globalVertexNumbers, globalCellNumbers */
2798ed5f475SMatthew G. Knepley   /* Reorder labels */
2808ed5f475SMatthew G. Knepley   {
281a85475f2SMatthew G. Knepley     PetscInt numLabels, l;
282a85475f2SMatthew G. Knepley     DMLabel  label, labelNew;
2838ed5f475SMatthew G. Knepley 
2849566063dSJacob Faibussowitsch     PetscCall(DMGetNumLabels(dm, &numLabels));
2855d80c0bfSVaclav Hapla     for (l = 0; l < numLabels; ++l) {
2869566063dSJacob Faibussowitsch       PetscCall(DMGetLabelByNum(dm, l, &label));
2879566063dSJacob Faibussowitsch       PetscCall(DMLabelPermute(label, perm, &labelNew));
2889566063dSJacob Faibussowitsch       PetscCall(DMAddLabel(*pdm, labelNew));
2899566063dSJacob Faibussowitsch       PetscCall(DMLabelDestroy(&labelNew));
2908ed5f475SMatthew G. Knepley     }
2919566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(*pdm, "depth", &(*pdm)->depthLabel));
2929566063dSJacob Faibussowitsch     if (plex->subpointMap) PetscCall(DMLabelPermute(plex->subpointMap, perm, &plexNew->subpointMap));
2938ed5f475SMatthew G. Knepley   }
294695799ffSMatthew G. Knepley   if ((*pdm)->celltypeLabel) {
295695799ffSMatthew G. Knepley     DMLabel ctLabel;
296695799ffSMatthew G. Knepley 
297695799ffSMatthew G. Knepley     // Reset label for fast lookup
298695799ffSMatthew G. Knepley     PetscCall(DMPlexGetCellTypeLabel(*pdm, &ctLabel));
299695799ffSMatthew G. Knepley     PetscCall(DMLabelMakeAllInvalid_Internal(ctLabel));
300695799ffSMatthew G. Knepley   }
3018ed5f475SMatthew G. Knepley   /* Reorder topology */
3028ed5f475SMatthew G. Knepley   {
3038ed5f475SMatthew G. Knepley     const PetscInt *pperm;
3046302a7fbSVaclav Hapla     PetscInt        n, pStart, pEnd, p;
3058ed5f475SMatthew G. Knepley 
3069566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&plexNew->coneSection));
3079566063dSJacob Faibussowitsch     PetscCall(PetscSectionPermute(plex->coneSection, perm, &plexNew->coneSection));
3089566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetStorageSize(plexNew->coneSection, &n));
3099566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &plexNew->cones));
3109566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &plexNew->coneOrientations));
3119566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(perm, &pperm));
3129566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(plex->coneSection, &pStart, &pEnd));
3138ed5f475SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
3148ed5f475SMatthew G. Knepley       PetscInt dof, off, offNew, d;
3158ed5f475SMatthew G. Knepley 
3169566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(plexNew->coneSection, pperm[p], &dof));
3179566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(plex->coneSection, p, &off));
3189566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(plexNew->coneSection, pperm[p], &offNew));
3198ed5f475SMatthew G. Knepley       for (d = 0; d < dof; ++d) {
3208ed5f475SMatthew G. Knepley         plexNew->cones[offNew + d]            = pperm[plex->cones[off + d]];
3218ed5f475SMatthew G. Knepley         plexNew->coneOrientations[offNew + d] = plex->coneOrientations[off + d];
3228ed5f475SMatthew G. Knepley       }
3238ed5f475SMatthew G. Knepley     }
3249566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&plexNew->supportSection));
3259566063dSJacob Faibussowitsch     PetscCall(PetscSectionPermute(plex->supportSection, perm, &plexNew->supportSection));
3269566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetStorageSize(plexNew->supportSection, &n));
3279566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &plexNew->supports));
3289566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(plex->supportSection, &pStart, &pEnd));
3298ed5f475SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
3308ed5f475SMatthew G. Knepley       PetscInt dof, off, offNew, d;
3318ed5f475SMatthew G. Knepley 
3329566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(plexNew->supportSection, pperm[p], &dof));
3339566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(plex->supportSection, p, &off));
3349566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(plexNew->supportSection, pperm[p], &offNew));
335ad540459SPierre Jolivet       for (d = 0; d < dof; ++d) plexNew->supports[offNew + d] = pperm[plex->supports[off + d]];
3368ed5f475SMatthew G. Knepley     }
3379566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(perm, &pperm));
3388ed5f475SMatthew G. Knepley   }
339a6e0b375SMatthew G. Knepley   /* Remap coordinates */
340a6e0b375SMatthew G. Knepley   {
341a6e0b375SMatthew G. Knepley     DM           cdm, cdmNew;
3426858538eSMatthew G. Knepley     PetscSection cs, csNew;
343a6e0b375SMatthew G. Knepley     Vec          coordinates, coordinatesNew;
344a6e0b375SMatthew G. Knepley 
3456858538eSMatthew G. Knepley     PetscCall(DMGetCoordinateSection(dm, &cs));
3469566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
3476858538eSMatthew G. Knepley     PetscCall(DMPlexRemapCoordinates_Private(perm, cs, coordinates, &csNew, &coordinatesNew));
3486858538eSMatthew G. Knepley     PetscCall(DMSetCoordinateSection(*pdm, PETSC_DETERMINE, csNew));
3499566063dSJacob Faibussowitsch     PetscCall(DMSetCoordinatesLocal(*pdm, coordinatesNew));
3506858538eSMatthew G. Knepley     PetscCall(PetscSectionDestroy(&csNew));
3519566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&coordinatesNew));
3526858538eSMatthew G. Knepley 
3536858538eSMatthew G. Knepley     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
3546858538eSMatthew G. Knepley     if (cdm) {
3556858538eSMatthew G. Knepley       PetscCall(DMGetCoordinateDM(*pdm, &cdm));
3566858538eSMatthew G. Knepley       PetscCall(DMClone(cdm, &cdmNew));
3576858538eSMatthew G. Knepley       PetscCall(DMSetCellCoordinateDM(*pdm, cdmNew));
3586858538eSMatthew G. Knepley       PetscCall(DMDestroy(&cdmNew));
3596858538eSMatthew G. Knepley       PetscCall(DMGetCellCoordinateSection(dm, &cs));
3606858538eSMatthew G. Knepley       PetscCall(DMGetCellCoordinatesLocal(dm, &coordinates));
3616858538eSMatthew G. Knepley       PetscCall(DMPlexRemapCoordinates_Private(perm, cs, coordinates, &csNew, &coordinatesNew));
3626858538eSMatthew G. Knepley       PetscCall(DMSetCellCoordinateSection(*pdm, PETSC_DETERMINE, csNew));
3636858538eSMatthew G. Knepley       PetscCall(DMSetCellCoordinatesLocal(*pdm, coordinatesNew));
3646858538eSMatthew G. Knepley       PetscCall(PetscSectionDestroy(&csNew));
3656858538eSMatthew G. Knepley       PetscCall(VecDestroy(&coordinatesNew));
3666858538eSMatthew G. Knepley     }
367a6e0b375SMatthew G. Knepley   }
3685de52c6dSVaclav Hapla   PetscCall(DMPlexCopy_Internal(dm, PETSC_TRUE, PETSC_TRUE, *pdm));
369e66f2fa0SMatthew G. Knepley   (*pdm)->setupcalled = PETSC_TRUE;
370*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3718ed5f475SMatthew G. Knepley }
3726bc1bd01Sksagiyam 
373d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexReorderSetDefault_Plex(DM dm, DMPlexReorderDefaultFlag reorder)
374d71ae5a4SJacob Faibussowitsch {
3756bc1bd01Sksagiyam   DM_Plex *mesh = (DM_Plex *)dm->data;
3766bc1bd01Sksagiyam 
3776bc1bd01Sksagiyam   PetscFunctionBegin;
3786bc1bd01Sksagiyam   mesh->reorderDefault = reorder;
379*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3806bc1bd01Sksagiyam }
3816bc1bd01Sksagiyam 
3826bc1bd01Sksagiyam /*@
3836bc1bd01Sksagiyam   DMPlexReorderSetDefault - Set flag indicating whether the DM should be reordered by default
3846bc1bd01Sksagiyam 
3856bc1bd01Sksagiyam   Logically collective
3866bc1bd01Sksagiyam 
3876bc1bd01Sksagiyam   Input Parameters:
3886bc1bd01Sksagiyam + dm        - The DM
3896bc1bd01Sksagiyam - reorder   - Flag for reordering
3906bc1bd01Sksagiyam 
3916bc1bd01Sksagiyam   Level: intermediate
3926bc1bd01Sksagiyam 
3936bc1bd01Sksagiyam .seealso: `DMPlexReorderGetDefault()`
3946bc1bd01Sksagiyam @*/
395d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexReorderSetDefault(DM dm, DMPlexReorderDefaultFlag reorder)
396d71ae5a4SJacob Faibussowitsch {
3976bc1bd01Sksagiyam   PetscFunctionBegin;
3986bc1bd01Sksagiyam   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3996bc1bd01Sksagiyam   PetscTryMethod(dm, "DMPlexReorderSetDefault_C", (DM, DMPlexReorderDefaultFlag), (dm, reorder));
400*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4016bc1bd01Sksagiyam }
4026bc1bd01Sksagiyam 
403d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexReorderGetDefault_Plex(DM dm, DMPlexReorderDefaultFlag *reorder)
404d71ae5a4SJacob Faibussowitsch {
4056bc1bd01Sksagiyam   DM_Plex *mesh = (DM_Plex *)dm->data;
4066bc1bd01Sksagiyam 
4076bc1bd01Sksagiyam   PetscFunctionBegin;
4086bc1bd01Sksagiyam   *reorder = mesh->reorderDefault;
409*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4106bc1bd01Sksagiyam }
4116bc1bd01Sksagiyam 
4126bc1bd01Sksagiyam /*@
4136bc1bd01Sksagiyam   DMPlexReorderGetDefault - Get flag indicating whether the DM should be reordered by default
4146bc1bd01Sksagiyam 
4156bc1bd01Sksagiyam   Not collective
4166bc1bd01Sksagiyam 
4176bc1bd01Sksagiyam   Input Parameter:
4186bc1bd01Sksagiyam . dm      - The DM
4196bc1bd01Sksagiyam 
4206bc1bd01Sksagiyam   Output Parameter:
4216bc1bd01Sksagiyam . reorder - Flag for reordering
4226bc1bd01Sksagiyam 
4236bc1bd01Sksagiyam   Level: intermediate
4246bc1bd01Sksagiyam 
4256bc1bd01Sksagiyam .seealso: `DMPlexReorderSetDefault()`
4266bc1bd01Sksagiyam @*/
427d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexReorderGetDefault(DM dm, DMPlexReorderDefaultFlag *reorder)
428d71ae5a4SJacob Faibussowitsch {
4296bc1bd01Sksagiyam   PetscFunctionBegin;
4306bc1bd01Sksagiyam   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4316bc1bd01Sksagiyam   PetscValidPointer(reorder, 2);
4326bc1bd01Sksagiyam   PetscUseMethod(dm, "DMPlexReorderGetDefault_C", (DM, DMPlexReorderDefaultFlag *), (dm, reorder));
433*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4346bc1bd01Sksagiyam }
435