xref: /petsc/src/dm/impls/plex/plexreorder.c (revision 563ffbab739d96c559dd28f4ad9be3f5bd74c24d)
1 #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2 #include <petsc-private/matorderimpl.h> /*I      "petscmat.h"      I*/
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "DMPlexCreateOrderingClosure_Static"
6 PetscErrorCode DMPlexCreateOrderingClosure_Static(DM dm, PetscInt numPoints, const PetscInt pperm[], PetscInt **clperm, PetscInt **invclperm)
7 {
8   PetscInt      *perm, *iperm;
9   PetscInt       depth, d, pStart, pEnd, fStart, fMax, fEnd, p;
10   PetscErrorCode ierr;
11 
12   PetscFunctionBegin;
13   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
14   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
15   ierr = PetscMalloc((pEnd-pStart) * sizeof(PetscInt), &perm);CHKERRQ(ierr);
16   ierr = PetscMalloc((pEnd-pStart) * sizeof(PetscInt) ,&iperm);CHKERRQ(ierr);
17   for (p = pStart; p < pEnd; ++p) iperm[p] = -1;
18   for (d = depth; d > 0; --d) {
19     ierr = DMPlexGetDepthStratum(dm, d,   &pStart, &pEnd);CHKERRQ(ierr);
20     ierr = DMPlexGetDepthStratum(dm, d-1, &fStart, &fEnd);CHKERRQ(ierr);
21     fMax = fStart;
22     for (p = pStart; p < pEnd; ++p) {
23       const PetscInt *cone;
24       PetscInt        point, coneSize, c;
25 
26       if (d == depth) {
27         perm[p]         = pperm[p];
28         iperm[pperm[p]] = p;
29       }
30       point = perm[p];
31       ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr);
32       ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr);
33       for (c = 0; c < coneSize; ++c) {
34         const PetscInt oldc = cone[c];
35         const PetscInt newc = iperm[oldc];
36 
37         if (newc < 0) {
38           perm[fMax]  = oldc;
39           iperm[oldc] = fMax++;
40         }
41       }
42     }
43     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);
44   }
45   *clperm    = perm;
46   *invclperm = iperm;
47   PetscFunctionReturn(0);
48 }
49 
50 #undef __FUNCT__
51 #define __FUNCT__ "DMPlexGetOrdering"
52 /*@
53   DMPlexGetOrdering - Calculate a reordering of the mesh
54 
55   Collective on DM
56 
57   Input Parameter:
58 + dm - The DMPlex object
59 - otype - type of reordering, one of the following:
60 $     MATORDERINGNATURAL - Natural
61 $     MATORDERINGND - Nested Dissection
62 $     MATORDERING1WD - One-way Dissection
63 $     MATORDERINGRCM - Reverse Cuthill-McKee
64 $     MATORDERINGQMD - Quotient Minimum Degree
65 
66 
67   Output Parameter:
68 . perm - The point permutation as an IS
69 
70   Level: intermediate
71 
72 .keywords: mesh
73 .seealso: MatGetOrdering()
74 @*/
75 PetscErrorCode DMPlexGetOrdering(DM dm, MatOrderingType otype, IS *perm)
76 {
77   PetscInt       numCells = 0;
78   PetscInt      *start = NULL, *adjacency = NULL, *cperm, *clperm, *invclperm, *mask, *xls, pStart, pEnd, c, i;
79   PetscErrorCode ierr;
80 
81   PetscFunctionBegin;
82   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
83   PetscValidPointer(perm, 3);
84   ierr = DMPlexCreateNeighborCSR(dm, 0, &numCells, &start, &adjacency);CHKERRQ(ierr);
85   ierr = PetscMalloc3(numCells,PetscInt,&cperm,numCells,PetscInt,&mask,numCells*2,PetscInt,&xls);CHKERRQ(ierr);
86   /* Shift for Fortran numbering */
87   for (i = 0; i < start[numCells]; ++i) ++adjacency[i];
88   for (i = 0; i <= numCells; ++i)       ++start[i];
89   ierr = SPARSEPACKgenrcm(&numCells, start, adjacency, cperm, mask, xls);CHKERRQ(ierr);
90   ierr = PetscFree(start);CHKERRQ(ierr);
91   ierr = PetscFree(adjacency);CHKERRQ(ierr);
92   /* Shift for Fortran numbering */
93   for (c = 0; c < numCells; ++c) --cperm[c];
94   /* Construct closure */
95   ierr = DMPlexCreateOrderingClosure_Static(dm, numCells, cperm, &clperm, &invclperm);
96   ierr = PetscFree3(cperm,mask,xls);CHKERRQ(ierr);
97   ierr = PetscFree(clperm);CHKERRQ(ierr);
98   /* Invert permutation */
99   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
100   ierr = ISCreateGeneral(PetscObjectComm((PetscObject) dm), pEnd-pStart, invclperm, PETSC_OWN_POINTER, perm);CHKERRQ(ierr);
101   PetscFunctionReturn(0);
102 }
103 
104 #undef __FUNCT__
105 #define __FUNCT__ "DMPlexPermute"
106 /*@
107   DMPlexPermute - Reorder the mesh according to the input permutation
108 
109   Collective on DM
110 
111   Input Parameter:
112 + dm - The DMPlex object
113 - perm - The point permutation
114 
115   Output Parameter:
116 . pdm - The permuted DM
117 
118   Level: intermediate
119 
120 .keywords: mesh
121 .seealso: MatPermute()
122 @*/
123 PetscErrorCode DMPlexPermute(DM dm, IS perm, DM *pdm)
124 {
125   DM_Plex       *plex = (DM_Plex *) dm->data, *plexNew;
126   PetscSection   section, sectionNew;
127   PetscInt       dim;
128   PetscErrorCode ierr;
129 
130   PetscFunctionBegin;
131   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
132   PetscValidHeaderSpecific(perm, IS_CLASSID, 2);
133   PetscValidPointer(pdm, 3);
134   ierr = DMCreate(PetscObjectComm((PetscObject) dm), pdm);CHKERRQ(ierr);
135   ierr = DMSetType(*pdm, DMPLEX);CHKERRQ(ierr);
136   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
137   ierr = DMPlexSetDimension(*pdm, dim);CHKERRQ(ierr);
138   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
139   if (section) {
140     ierr = PetscSectionPermute(section, perm, &sectionNew);CHKERRQ(ierr);
141     ierr = DMSetDefaultSection(*pdm, sectionNew);CHKERRQ(ierr);
142     ierr = PetscSectionDestroy(&sectionNew);CHKERRQ(ierr);
143   }
144   plexNew = (DM_Plex *) (*pdm)->data;
145   /* Ignore ltogmap, ltogmapb */
146   /* Ignore sf, defaultSF */
147   /* Ignore globalVertexNumbers, globalCellNumbers */
148   /* Remap coordinates */
149   {
150     DM           cdm, cdmNew;
151     PetscSection csection, csectionNew;
152     Vec          coordinates, coordinatesNew;
153     PetscScalar *coords, *coordsNew;
154     PetscInt     pStart, pEnd, p;
155 
156     ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
157     ierr = DMGetDefaultSection(cdm, &csection);CHKERRQ(ierr);
158     ierr = PetscSectionPermute(csection, perm, &csectionNew);CHKERRQ(ierr);
159     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
160     ierr = VecDuplicate(coordinates, &coordinatesNew);CHKERRQ(ierr);
161     ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
162     ierr = VecGetArray(coordinatesNew, &coordsNew);CHKERRQ(ierr);
163     ierr = PetscSectionGetChart(csectionNew, &pStart, &pEnd);CHKERRQ(ierr);
164     for (p = pStart; p < pEnd; ++p) {
165       PetscInt dof, off, offNew, d;
166 
167       ierr = PetscSectionGetDof(csectionNew, p, &dof);CHKERRQ(ierr);
168       ierr = PetscSectionGetOffset(csection, p, &off);CHKERRQ(ierr);
169       ierr = PetscSectionGetOffset(csectionNew, p, &offNew);CHKERRQ(ierr);
170       for (d = 0; d < dof; ++d) coordsNew[offNew+d] = coords[off+d];
171     }
172     ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
173     ierr = VecRestoreArray(coordinatesNew, &coordsNew);CHKERRQ(ierr);
174     ierr = DMGetCoordinateDM(*pdm, &cdmNew);CHKERRQ(ierr);
175     ierr = DMSetDefaultSection(cdmNew, csectionNew);CHKERRQ(ierr);
176     ierr = DMSetCoordinatesLocal(*pdm, coordinatesNew);CHKERRQ(ierr);
177     ierr = PetscSectionDestroy(&csectionNew);CHKERRQ(ierr);
178     ierr = VecDestroy(&coordinatesNew);CHKERRQ(ierr);
179   }
180   /* Reorder labels */
181   {
182     DMLabel label = plex->labels, labelNew = NULL;
183 
184     while (label) {
185       if (!plexNew->labels) {
186         ierr = DMLabelPermute(label, perm, &plexNew->labels);CHKERRQ(ierr);
187         labelNew = plexNew->labels;
188       } else {
189         ierr = DMLabelPermute(label, perm, &labelNew->next);CHKERRQ(ierr);
190         labelNew = labelNew->next;
191       }
192       label = label->next;
193     }
194     if (plex->subpointMap) {ierr = DMLabelPermute(plex->subpointMap, perm, &plexNew->subpointMap);CHKERRQ(ierr);}
195   }
196   /* Reorder topology */
197   {
198     const PetscInt *pperm;
199     PetscInt        maxConeSize, maxSupportSize, n, pStart, pEnd, p;
200 
201     ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr);
202     plexNew->maxConeSize    = maxConeSize;
203     plexNew->maxSupportSize = maxSupportSize;
204     ierr = PetscSectionDestroy(&plexNew->coneSection);CHKERRQ(ierr);
205     ierr = PetscSectionPermute(plex->coneSection, perm, &plexNew->coneSection);CHKERRQ(ierr);
206     ierr = PetscSectionGetStorageSize(plexNew->coneSection, &n);CHKERRQ(ierr);
207     ierr = PetscMalloc(n * sizeof(PetscInt), &plexNew->cones);CHKERRQ(ierr);
208     ierr = PetscMalloc(n * sizeof(PetscInt), &plexNew->coneOrientations);CHKERRQ(ierr);
209     ierr = ISGetIndices(perm, &pperm);CHKERRQ(ierr);
210     ierr = PetscSectionGetChart(plex->coneSection, &pStart, &pEnd);CHKERRQ(ierr);
211     for (p = pStart; p < pEnd; ++p) {
212       PetscInt dof, off, offNew, d;
213 
214       ierr = PetscSectionGetDof(plexNew->coneSection, pperm[p], &dof);CHKERRQ(ierr);
215       ierr = PetscSectionGetOffset(plex->coneSection, p, &off);CHKERRQ(ierr);
216       ierr = PetscSectionGetOffset(plexNew->coneSection, pperm[p], &offNew);CHKERRQ(ierr);
217       for (d = 0; d < dof; ++d) {
218         plexNew->cones[offNew+d]            = pperm[plex->cones[off+d]];
219         plexNew->coneOrientations[offNew+d] = plex->coneOrientations[off+d];
220       }
221     }
222     ierr = PetscSectionDestroy(&plexNew->supportSection);CHKERRQ(ierr);
223     ierr = PetscSectionPermute(plex->supportSection, perm, &plexNew->supportSection);CHKERRQ(ierr);
224     ierr = PetscSectionGetStorageSize(plexNew->supportSection, &n);CHKERRQ(ierr);
225     ierr = PetscMalloc(n * sizeof(PetscInt), &plexNew->supports);CHKERRQ(ierr);
226     ierr = PetscSectionGetChart(plex->supportSection, &pStart, &pEnd);CHKERRQ(ierr);
227     for (p = pStart; p < pEnd; ++p) {
228       PetscInt dof, off, offNew, d;
229 
230       ierr = PetscSectionGetDof(plexNew->supportSection, pperm[p], &dof);CHKERRQ(ierr);
231       ierr = PetscSectionGetOffset(plex->supportSection, p, &off);CHKERRQ(ierr);
232       ierr = PetscSectionGetOffset(plexNew->supportSection, pperm[p], &offNew);CHKERRQ(ierr);
233       for (d = 0; d < dof; ++d) {
234         plexNew->supports[offNew+d] = pperm[plex->supports[off+d]];
235       }
236     }
237     ierr = ISRestoreIndices(perm, &pperm);CHKERRQ(ierr);
238   }
239   PetscFunctionReturn(0);
240 }
241