xref: /petsc/src/dm/impls/plex/plexreorder.c (revision dd5a1fb8a1b1dc1df5cd02569a30d8a384bbbc95)
18ed5f475SMatthew G. Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
28ed5f475SMatthew G. Knepley #include <petsc-private/matorderimpl.h> /*I      "petscmat.h"      I*/
38ed5f475SMatthew G. Knepley 
48ed5f475SMatthew G. Knepley #undef __FUNCT__
58ed5f475SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateOrderingClosure_Static"
68ed5f475SMatthew G. Knepley PetscErrorCode DMPlexCreateOrderingClosure_Static(DM dm, PetscInt numPoints, const PetscInt pperm[], PetscInt **clperm, PetscInt **invclperm)
78ed5f475SMatthew G. Knepley {
88ed5f475SMatthew G. Knepley   PetscInt      *perm, *iperm;
98ed5f475SMatthew G. Knepley   PetscInt       depth, d, pStart, pEnd, fStart, fMax, fEnd, p;
108ed5f475SMatthew G. Knepley   PetscErrorCode ierr;
118ed5f475SMatthew G. Knepley 
128ed5f475SMatthew G. Knepley   PetscFunctionBegin;
138ed5f475SMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
148ed5f475SMatthew G. Knepley   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
158ed5f475SMatthew G. Knepley   ierr = PetscMalloc2(pEnd-pStart,PetscInt,&perm,pEnd-pStart,PetscInt,&iperm);CHKERRQ(ierr);
168ed5f475SMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) iperm[p] = -1;
178ed5f475SMatthew G. Knepley   for (d = depth; d > 0; --d) {
188ed5f475SMatthew G. Knepley     ierr = DMPlexGetDepthStratum(dm, d,   &pStart, &pEnd);CHKERRQ(ierr);
198ed5f475SMatthew G. Knepley     ierr = DMPlexGetDepthStratum(dm, d-1, &fStart, &fEnd);CHKERRQ(ierr);
208ed5f475SMatthew G. Knepley     fMax = fStart;
218ed5f475SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
228ed5f475SMatthew G. Knepley       const PetscInt *cone;
238ed5f475SMatthew G. Knepley       PetscInt        point, coneSize, c;
248ed5f475SMatthew G. Knepley 
258ed5f475SMatthew G. Knepley       if (d == depth) {
268ed5f475SMatthew G. Knepley         perm[p]         = pperm[p];
278ed5f475SMatthew G. Knepley         iperm[pperm[p]] = p;
288ed5f475SMatthew G. Knepley       }
298ed5f475SMatthew G. Knepley       point = perm[p];
308ed5f475SMatthew G. Knepley       ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr);
318ed5f475SMatthew G. Knepley       ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr);
328ed5f475SMatthew G. Knepley       for (c = 0; c < coneSize; ++c) {
338ed5f475SMatthew G. Knepley         const PetscInt oldc = cone[c];
348ed5f475SMatthew G. Knepley         const PetscInt newc = iperm[oldc];
358ed5f475SMatthew G. Knepley 
368ed5f475SMatthew G. Knepley         if (newc < 0) {
378ed5f475SMatthew G. Knepley           perm[fMax]  = oldc;
388ed5f475SMatthew G. Knepley           iperm[oldc] = fMax++;
398ed5f475SMatthew G. Knepley         }
408ed5f475SMatthew G. Knepley       }
418ed5f475SMatthew G. Knepley     }
428ed5f475SMatthew G. Knepley     if (fMax != fEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of depth %d faces %d does not match permuted nubmer %d", d, fEnd-fStart, fMax-fStart);
438ed5f475SMatthew G. Knepley   }
448ed5f475SMatthew G. Knepley   *clperm    = perm;
458ed5f475SMatthew G. Knepley   *invclperm = iperm;
468ed5f475SMatthew G. Knepley   PetscFunctionReturn(0);
478ed5f475SMatthew G. Knepley }
488ed5f475SMatthew G. Knepley 
498ed5f475SMatthew G. Knepley #undef __FUNCT__
508ed5f475SMatthew G. Knepley #define __FUNCT__ "DMPlexGetOrdering"
518ed5f475SMatthew G. Knepley /*@
528ed5f475SMatthew G. Knepley   DMPlexGetOrdering - Calculate a reordering of the mesh
538ed5f475SMatthew G. Knepley 
548ed5f475SMatthew G. Knepley   Collective on DM
558ed5f475SMatthew G. Knepley 
568ed5f475SMatthew G. Knepley   Input Parameter:
578ed5f475SMatthew G. Knepley + dm - The DMPlex object
588ed5f475SMatthew G. Knepley - otype - type of reordering, one of the following:
598ed5f475SMatthew G. Knepley $     MATORDERINGNATURAL - Natural
608ed5f475SMatthew G. Knepley $     MATORDERINGND - Nested Dissection
618ed5f475SMatthew G. Knepley $     MATORDERING1WD - One-way Dissection
628ed5f475SMatthew G. Knepley $     MATORDERINGRCM - Reverse Cuthill-McKee
638ed5f475SMatthew G. Knepley $     MATORDERINGQMD - Quotient Minimum Degree
648ed5f475SMatthew G. Knepley 
658ed5f475SMatthew G. Knepley 
668ed5f475SMatthew G. Knepley   Output Parameter:
678ed5f475SMatthew G. Knepley . perm - The point permutation as an IS
688ed5f475SMatthew G. Knepley 
698ed5f475SMatthew G. Knepley   Level: intermediate
708ed5f475SMatthew G. Knepley 
718ed5f475SMatthew G. Knepley .keywords: mesh
728ed5f475SMatthew G. Knepley .seealso: MatGetOrdering()
738ed5f475SMatthew G. Knepley @*/
748ed5f475SMatthew G. Knepley PetscErrorCode DMPlexGetOrdering(DM dm, MatOrderingType otype, IS *perm)
758ed5f475SMatthew G. Knepley {
768ed5f475SMatthew G. Knepley   PetscInt       numCells = 0;
778ed5f475SMatthew G. Knepley   PetscInt      *start = NULL, *adjacency = NULL, *cperm, *clperm, *invclperm, *mask, *xls, pStart, pEnd, c, i;
788ed5f475SMatthew G. Knepley   PetscErrorCode ierr;
798ed5f475SMatthew G. Knepley 
808ed5f475SMatthew G. Knepley   PetscFunctionBegin;
818ed5f475SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
828ed5f475SMatthew G. Knepley   PetscValidPointer(perm, 3);
838ed5f475SMatthew G. Knepley   ierr = DMPlexCreateNeighborCSR(dm, 0, &numCells, &start, &adjacency);CHKERRQ(ierr);
848ed5f475SMatthew G. Knepley   ierr = PetscMalloc3(numCells,PetscInt,&cperm,numCells,PetscInt,&mask,numCells*2,PetscInt,&xls);CHKERRQ(ierr);
858ed5f475SMatthew G. Knepley   /* Shift for Fortran numbering */
868ed5f475SMatthew G. Knepley   for (i = 0; i < start[numCells]; ++i) ++adjacency[i];
878ed5f475SMatthew G. Knepley   for (i = 0; i <= numCells; ++i)       ++start[i];
888ed5f475SMatthew G. Knepley   ierr = SPARSEPACKgenrcm(&numCells, start, adjacency, cperm, mask, xls);CHKERRQ(ierr);
898ed5f475SMatthew G. Knepley   ierr = PetscFree(start);CHKERRQ(ierr);
908ed5f475SMatthew G. Knepley   ierr = PetscFree(adjacency);CHKERRQ(ierr);
918ed5f475SMatthew G. Knepley   /* Shift for Fortran numbering */
928ed5f475SMatthew G. Knepley   for (c = 0; c < numCells; ++c) --cperm[c];
938ed5f475SMatthew G. Knepley   /* Construct closure */
948ed5f475SMatthew G. Knepley   ierr = DMPlexCreateOrderingClosure_Static(dm, numCells, cperm, &clperm, &invclperm);
958ed5f475SMatthew G. Knepley   ierr = PetscFree3(cperm,mask,xls);CHKERRQ(ierr);
968ed5f475SMatthew G. Knepley   ierr = PetscFree(clperm);CHKERRQ(ierr);
978ed5f475SMatthew G. Knepley   /* Invert permutation */
988ed5f475SMatthew G. Knepley   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
998ed5f475SMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject) dm), pEnd-pStart, invclperm, PETSC_OWN_POINTER, perm);CHKERRQ(ierr);
1008ed5f475SMatthew G. Knepley   PetscFunctionReturn(0);
1018ed5f475SMatthew G. Knepley }
1028ed5f475SMatthew G. Knepley 
1038ed5f475SMatthew G. Knepley #undef __FUNCT__
1048ed5f475SMatthew G. Knepley #define __FUNCT__ "DMPlexPermute"
1058ed5f475SMatthew G. Knepley /*@
1068ed5f475SMatthew G. Knepley   DMPlexPermute - Reorder the mesh according to the input permutation
1078ed5f475SMatthew G. Knepley 
1088ed5f475SMatthew G. Knepley   Collective on DM
1098ed5f475SMatthew G. Knepley 
1108ed5f475SMatthew G. Knepley   Input Parameter:
1118ed5f475SMatthew G. Knepley + dm - The DMPlex object
1128ed5f475SMatthew G. Knepley - perm - The point permutation
1138ed5f475SMatthew G. Knepley 
1148ed5f475SMatthew G. Knepley   Output Parameter:
1158ed5f475SMatthew G. Knepley . pdm - The permuted DM
1168ed5f475SMatthew G. Knepley 
1178ed5f475SMatthew G. Knepley   Level: intermediate
1188ed5f475SMatthew G. Knepley 
1198ed5f475SMatthew G. Knepley .keywords: mesh
1208ed5f475SMatthew G. Knepley .seealso: MatPermute()
1218ed5f475SMatthew G. Knepley @*/
1228ed5f475SMatthew G. Knepley PetscErrorCode DMPlexPermute(DM dm, IS perm, DM *pdm)
1238ed5f475SMatthew G. Knepley {
1248ed5f475SMatthew G. Knepley   DM_Plex       *plex = (DM_Plex *) dm->data, *plexNew;
1258ed5f475SMatthew G. Knepley   PetscSection   section, sectionNew;
1268ed5f475SMatthew G. Knepley   PetscInt       dim;
1278ed5f475SMatthew G. Knepley   PetscErrorCode ierr;
1288ed5f475SMatthew G. Knepley 
1298ed5f475SMatthew G. Knepley   PetscFunctionBegin;
1308ed5f475SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1318ed5f475SMatthew G. Knepley   PetscValidHeaderSpecific(perm, IS_CLASSID, 2);
1328ed5f475SMatthew G. Knepley   PetscValidPointer(pdm, 3);
1338ed5f475SMatthew G. Knepley   ierr = DMCreate(PetscObjectComm((PetscObject) dm), pdm);CHKERRQ(ierr);
1348ed5f475SMatthew G. Knepley   ierr = DMSetType(*pdm, DMPLEX);CHKERRQ(ierr);
1358ed5f475SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1368ed5f475SMatthew G. Knepley   ierr = DMPlexSetDimension(*pdm, dim);CHKERRQ(ierr);
1378ed5f475SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
138dd742cf1SMatthew G. Knepley   if (section) {
1398ed5f475SMatthew G. Knepley     ierr = PetscSectionPermute(section, perm, &sectionNew);CHKERRQ(ierr);
1408ed5f475SMatthew G. Knepley     ierr = DMSetDefaultSection(*pdm, sectionNew);CHKERRQ(ierr);
1418ed5f475SMatthew G. Knepley     ierr = PetscSectionDestroy(&sectionNew);CHKERRQ(ierr);
142dd742cf1SMatthew G. Knepley   }
1438ed5f475SMatthew G. Knepley   plexNew = (DM_Plex *) (*pdm)->data;
1448ed5f475SMatthew G. Knepley   /* Ignore ltogmap, ltogmapb */
1458ed5f475SMatthew G. Knepley   /* Ignore sf, defaultSF */
1468ed5f475SMatthew G. Knepley   /* Ignore globalVertexNumbers, globalCellNumbers */
1478ed5f475SMatthew G. Knepley   /* Remap coordinates */
1488ed5f475SMatthew G. Knepley   {
1498ed5f475SMatthew G. Knepley     DM              cdm, cdmNew;
1508ed5f475SMatthew G. Knepley     PetscSection    csection, csectionNew;
1518ed5f475SMatthew G. Knepley     Vec             coordinates, coordinatesNew;
1528ed5f475SMatthew G. Knepley     PetscScalar    *coords, *coordsNew;
153*dd5a1fb8SMatthew G. Knepley     const PetscInt *pperm;
1548ed5f475SMatthew G. Knepley     PetscInt        pStart, pEnd, p;
1558ed5f475SMatthew G. Knepley 
1568ed5f475SMatthew G. Knepley     ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
1578ed5f475SMatthew G. Knepley     ierr = DMGetDefaultSection(cdm, &csection);CHKERRQ(ierr);
1588ed5f475SMatthew G. Knepley     ierr = PetscSectionPermute(csection, perm, &csectionNew);CHKERRQ(ierr);
1598ed5f475SMatthew G. Knepley     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1608ed5f475SMatthew G. Knepley     ierr = VecDuplicate(coordinates, &coordinatesNew);CHKERRQ(ierr);
1618ed5f475SMatthew G. Knepley     ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
1628ed5f475SMatthew G. Knepley     ierr = VecGetArray(coordinatesNew, &coordsNew);CHKERRQ(ierr);
1638ed5f475SMatthew G. Knepley     ierr = PetscSectionGetChart(csectionNew, &pStart, &pEnd);CHKERRQ(ierr);
164*dd5a1fb8SMatthew G. Knepley     ierr = ISGetIndices(perm, &pperm);CHKERRQ(ierr);
1658ed5f475SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
1668ed5f475SMatthew G. Knepley       PetscInt dof, off, offNew, d;
1678ed5f475SMatthew G. Knepley 
1688ed5f475SMatthew G. Knepley       ierr = PetscSectionGetDof(csectionNew, p, &dof);CHKERRQ(ierr);
1698ed5f475SMatthew G. Knepley       ierr = PetscSectionGetOffset(csection, p, &off);CHKERRQ(ierr);
170*dd5a1fb8SMatthew G. Knepley       ierr = PetscSectionGetOffset(csectionNew, pperm[p], &offNew);CHKERRQ(ierr);
1718ed5f475SMatthew G. Knepley       for (d = 0; d < dof; ++d) coordsNew[offNew+d] = coords[off+d];
1728ed5f475SMatthew G. Knepley     }
173*dd5a1fb8SMatthew G. Knepley     ierr = ISRestoreIndices(perm, &pperm);CHKERRQ(ierr);
1748ed5f475SMatthew G. Knepley     ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
1758ed5f475SMatthew G. Knepley     ierr = VecRestoreArray(coordinatesNew, &coordsNew);CHKERRQ(ierr);
1768ed5f475SMatthew G. Knepley     ierr = DMGetCoordinateDM(*pdm, &cdmNew);CHKERRQ(ierr);
1778ed5f475SMatthew G. Knepley     ierr = DMSetDefaultSection(cdmNew, csectionNew);CHKERRQ(ierr);
1788ed5f475SMatthew G. Knepley     ierr = DMSetCoordinatesLocal(*pdm, coordinatesNew);CHKERRQ(ierr);
1798ed5f475SMatthew G. Knepley     ierr = PetscSectionDestroy(&csectionNew);CHKERRQ(ierr);
1808ed5f475SMatthew G. Knepley     ierr = VecDestroy(&coordinatesNew);CHKERRQ(ierr);
1818ed5f475SMatthew G. Knepley   }
1828ed5f475SMatthew G. Knepley   /* Reorder labels */
1838ed5f475SMatthew G. Knepley   {
184211f4af5SMatthew G. Knepley     DMLabel label = plex->labels, labelNew = NULL;
1858ed5f475SMatthew G. Knepley 
1868ed5f475SMatthew G. Knepley     while (label) {
1878ed5f475SMatthew G. Knepley       if (!plexNew->labels) {
1888ed5f475SMatthew G. Knepley         ierr = DMLabelPermute(label, perm, &plexNew->labels);CHKERRQ(ierr);
1898ed5f475SMatthew G. Knepley         labelNew = plexNew->labels;
1908ed5f475SMatthew G. Knepley       } else {
1918ed5f475SMatthew G. Knepley         ierr = DMLabelPermute(label, perm, &labelNew->next);CHKERRQ(ierr);
19283808cfdSMatthew G. Knepley         labelNew = labelNew->next;
1938ed5f475SMatthew G. Knepley       }
1948ed5f475SMatthew G. Knepley       label = label->next;
1958ed5f475SMatthew G. Knepley     }
1968ed5f475SMatthew G. Knepley     if (plex->subpointMap) {ierr = DMLabelPermute(plex->subpointMap, perm, &plexNew->subpointMap);CHKERRQ(ierr);}
1978ed5f475SMatthew G. Knepley   }
1988ed5f475SMatthew G. Knepley   /* Reorder topology */
1998ed5f475SMatthew G. Knepley   {
2008ed5f475SMatthew G. Knepley     const PetscInt *pperm;
2018ed5f475SMatthew G. Knepley     PetscInt        maxConeSize, maxSupportSize, n, pStart, pEnd, p;
2028ed5f475SMatthew G. Knepley 
2038ed5f475SMatthew G. Knepley     ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr);
2048ed5f475SMatthew G. Knepley     plexNew->maxConeSize    = maxConeSize;
2058ed5f475SMatthew G. Knepley     plexNew->maxSupportSize = maxSupportSize;
2068ed5f475SMatthew G. Knepley     ierr = PetscSectionDestroy(&plexNew->coneSection);CHKERRQ(ierr);
2078ed5f475SMatthew G. Knepley     ierr = PetscSectionPermute(plex->coneSection, perm, &plexNew->coneSection);CHKERRQ(ierr);
2088ed5f475SMatthew G. Knepley     ierr = PetscSectionGetStorageSize(plexNew->coneSection, &n);CHKERRQ(ierr);
2098ed5f475SMatthew G. Knepley     ierr = PetscMalloc(n * sizeof(PetscInt), &plexNew->cones);CHKERRQ(ierr);
2108ed5f475SMatthew G. Knepley     ierr = PetscMalloc(n * sizeof(PetscInt), &plexNew->coneOrientations);CHKERRQ(ierr);
2118ed5f475SMatthew G. Knepley     ierr = ISGetIndices(perm, &pperm);CHKERRQ(ierr);
2128ed5f475SMatthew G. Knepley     ierr = PetscSectionGetChart(plex->coneSection, &pStart, &pEnd);CHKERRQ(ierr);
2138ed5f475SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
2148ed5f475SMatthew G. Knepley       PetscInt dof, off, offNew, d;
2158ed5f475SMatthew G. Knepley 
2168ed5f475SMatthew G. Knepley       ierr = PetscSectionGetDof(plexNew->coneSection, pperm[p], &dof);CHKERRQ(ierr);
2178ed5f475SMatthew G. Knepley       ierr = PetscSectionGetOffset(plex->coneSection, p, &off);CHKERRQ(ierr);
2188ed5f475SMatthew G. Knepley       ierr = PetscSectionGetOffset(plexNew->coneSection, pperm[p], &offNew);CHKERRQ(ierr);
2198ed5f475SMatthew G. Knepley       for (d = 0; d < dof; ++d) {
2208ed5f475SMatthew G. Knepley         plexNew->cones[offNew+d]            = pperm[plex->cones[off+d]];
2218ed5f475SMatthew G. Knepley         plexNew->coneOrientations[offNew+d] = plex->coneOrientations[off+d];
2228ed5f475SMatthew G. Knepley       }
2238ed5f475SMatthew G. Knepley     }
2248ed5f475SMatthew G. Knepley     ierr = PetscSectionDestroy(&plexNew->supportSection);CHKERRQ(ierr);
2258ed5f475SMatthew G. Knepley     ierr = PetscSectionPermute(plex->supportSection, perm, &plexNew->supportSection);CHKERRQ(ierr);
2268ed5f475SMatthew G. Knepley     ierr = PetscSectionGetStorageSize(plexNew->supportSection, &n);CHKERRQ(ierr);
2278ed5f475SMatthew G. Knepley     ierr = PetscMalloc(n * sizeof(PetscInt), &plexNew->supports);CHKERRQ(ierr);
2288ed5f475SMatthew G. Knepley     ierr = PetscSectionGetChart(plex->supportSection, &pStart, &pEnd);CHKERRQ(ierr);
2298ed5f475SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
2308ed5f475SMatthew G. Knepley       PetscInt dof, off, offNew, d;
2318ed5f475SMatthew G. Knepley 
2328ed5f475SMatthew G. Knepley       ierr = PetscSectionGetDof(plexNew->supportSection, pperm[p], &dof);CHKERRQ(ierr);
2338ed5f475SMatthew G. Knepley       ierr = PetscSectionGetOffset(plex->supportSection, p, &off);CHKERRQ(ierr);
2348ed5f475SMatthew G. Knepley       ierr = PetscSectionGetOffset(plexNew->supportSection, pperm[p], &offNew);CHKERRQ(ierr);
2358ed5f475SMatthew G. Knepley       for (d = 0; d < dof; ++d) {
2368ed5f475SMatthew G. Knepley         plexNew->supports[offNew+d] = pperm[plex->supports[off+d]];
2378ed5f475SMatthew G. Knepley       }
2388ed5f475SMatthew G. Knepley     }
2398ed5f475SMatthew G. Knepley     ierr = ISRestoreIndices(perm, &pperm);CHKERRQ(ierr);
2408ed5f475SMatthew G. Knepley   }
2418ed5f475SMatthew G. Knepley   PetscFunctionReturn(0);
2428ed5f475SMatthew G. Knepley }
243