xref: /petsc/src/dm/impls/plex/plexreorder.c (revision 1ebf93c6b7d760d592de6ebe6cdc0debaa3caf75)
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 = PetscMalloc1(pEnd-pStart,&perm);CHKERRQ(ierr);
16   ierr = PetscMalloc1(pEnd-pStart,&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 - label - [Optional] Label used to segregate ordering into sets, or NULL
66 
67 
68   Output Parameter:
69 . perm - The point permutation as an IS, perm[old point number] = new point number
70 
71   Note: The label is used to group sets of points together by label value. This makes it easy to reorder a mesh which
72   has different types of cells, and then loop over each set of reordered cells for assembly.
73 
74   Level: intermediate
75 
76 .keywords: mesh
77 .seealso: MatGetOrdering()
78 @*/
79 PetscErrorCode DMPlexGetOrdering(DM dm, MatOrderingType otype, DMLabel label, IS *perm)
80 {
81   PetscInt       numCells = 0;
82   PetscInt      *start = NULL, *adjacency = NULL, *cperm, *clperm, *invclperm, *mask, *xls, pStart, pEnd, c, i;
83   PetscErrorCode ierr;
84 
85   PetscFunctionBegin;
86   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
87   PetscValidPointer(perm, 3);
88   ierr = DMPlexCreateNeighborCSR(dm, 0, &numCells, &start, &adjacency);CHKERRQ(ierr);
89   ierr = PetscMalloc3(numCells,&cperm,numCells,&mask,numCells*2,&xls);CHKERRQ(ierr);
90   if (numCells) {
91     /* Shift for Fortran numbering */
92     for (i = 0; i < start[numCells]; ++i) ++adjacency[i];
93     for (i = 0; i <= numCells; ++i)       ++start[i];
94     ierr = SPARSEPACKgenrcm(&numCells, start, adjacency, cperm, mask, xls);CHKERRQ(ierr);
95   }
96   ierr = PetscFree(start);CHKERRQ(ierr);
97   ierr = PetscFree(adjacency);CHKERRQ(ierr);
98   /* Shift for Fortran numbering */
99   for (c = 0; c < numCells; ++c) --cperm[c];
100   /* Segregate */
101   if (label) {
102     IS              valueIS;
103     const PetscInt *values;
104     PetscInt        numValues, numPoints = 0;
105     PetscInt       *sperm, *vsize, *voff, v;
106 
107     ierr = DMLabelGetValueIS(label, &valueIS);CHKERRQ(ierr);
108     ierr = ISGetLocalSize(valueIS, &numValues);CHKERRQ(ierr);
109     ierr = ISGetIndices(valueIS, &values);CHKERRQ(ierr);
110     ierr = PetscCalloc3(numCells,&sperm,numValues,&vsize,numValues+1,&voff);CHKERRQ(ierr);
111     for (v = 0; v < numValues; ++v) {
112       ierr = DMLabelGetStratumSize(label, values[v], &vsize[v]);CHKERRQ(ierr);
113       if (v < numValues-1) voff[v+2] += vsize[v] + voff[v+1];
114       numPoints += vsize[v];
115     }
116     if (numPoints != numCells) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label only covers %D cells < %D total", numPoints, numCells);
117     for (c = 0; c < numCells; ++c) {
118       const PetscInt oldc = cperm[c];
119       PetscInt       val, vloc;
120 
121       ierr = DMLabelGetValue(label, oldc, &val);CHKERRQ(ierr);
122       if (val == -1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cell %D not present in label", oldc);
123       ierr = PetscFindInt(val, numValues, values, &vloc);CHKERRQ(ierr);
124       if (vloc < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Value %D not present label", val);
125       sperm[voff[vloc+1]++] = oldc;
126     }
127     for (v = 0; v < numValues; ++v) {
128       if (voff[v+1] - voff[v] != vsize[v]) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of %D values found is %D != %D", values[v], voff[v+1] - voff[v], vsize[v]);
129     }
130     ierr = ISRestoreIndices(valueIS, &values);CHKERRQ(ierr);
131     ierr = ISDestroy(&valueIS);CHKERRQ(ierr);
132     ierr = PetscMemcpy(cperm, sperm, numCells * sizeof(PetscInt));CHKERRQ(ierr);
133     ierr = PetscFree3(sperm, vsize, voff);CHKERRQ(ierr);
134   }
135   /* Construct closure */
136   ierr = DMPlexCreateOrderingClosure_Static(dm, numCells, cperm, &clperm, &invclperm);CHKERRQ(ierr);
137   ierr = PetscFree3(cperm,mask,xls);CHKERRQ(ierr);
138   ierr = PetscFree(clperm);CHKERRQ(ierr);
139   /* Invert permutation */
140   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
141   ierr = ISCreateGeneral(PetscObjectComm((PetscObject) dm), pEnd-pStart, invclperm, PETSC_OWN_POINTER, perm);CHKERRQ(ierr);
142   PetscFunctionReturn(0);
143 }
144 
145 #undef __FUNCT__
146 #define __FUNCT__ "DMPlexPermute"
147 /*@
148   DMPlexPermute - Reorder the mesh according to the input permutation
149 
150   Collective on DM
151 
152   Input Parameter:
153 + dm - The DMPlex object
154 - perm - The point permutation, perm[old point number] = new point number
155 
156   Output Parameter:
157 . pdm - The permuted DM
158 
159   Level: intermediate
160 
161 .keywords: mesh
162 .seealso: MatPermute()
163 @*/
164 PetscErrorCode DMPlexPermute(DM dm, IS perm, DM *pdm)
165 {
166   DM_Plex       *plex = (DM_Plex *) dm->data, *plexNew;
167   PetscSection   section, sectionNew;
168   PetscInt       dim;
169   PetscErrorCode ierr;
170 
171   PetscFunctionBegin;
172   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
173   PetscValidHeaderSpecific(perm, IS_CLASSID, 2);
174   PetscValidPointer(pdm, 3);
175   ierr = DMCreate(PetscObjectComm((PetscObject) dm), pdm);CHKERRQ(ierr);
176   ierr = DMSetType(*pdm, DMPLEX);CHKERRQ(ierr);
177   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
178   ierr = DMSetDimension(*pdm, dim);CHKERRQ(ierr);
179   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
180   if (section) {
181     ierr = PetscSectionPermute(section, perm, &sectionNew);CHKERRQ(ierr);
182     ierr = DMSetDefaultSection(*pdm, sectionNew);CHKERRQ(ierr);
183     ierr = PetscSectionDestroy(&sectionNew);CHKERRQ(ierr);
184   }
185   plexNew = (DM_Plex *) (*pdm)->data;
186   /* Ignore ltogmap, ltogmapb */
187   /* Ignore sf, defaultSF */
188   /* Ignore globalVertexNumbers, globalCellNumbers */
189   /* Remap coordinates */
190   {
191     DM              cdm, cdmNew;
192     PetscSection    csection, csectionNew;
193     Vec             coordinates, coordinatesNew;
194     PetscScalar    *coords, *coordsNew;
195     const PetscInt *pperm;
196     PetscInt        pStart, pEnd, p;
197     const char     *name;
198 
199     ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
200     ierr = DMGetDefaultSection(cdm, &csection);CHKERRQ(ierr);
201     ierr = PetscSectionPermute(csection, perm, &csectionNew);CHKERRQ(ierr);
202     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
203     ierr = VecDuplicate(coordinates, &coordinatesNew);CHKERRQ(ierr);
204     ierr = PetscObjectGetName((PetscObject)coordinates,&name);CHKERRQ(ierr);
205     ierr = PetscObjectSetName((PetscObject)coordinatesNew,name);CHKERRQ(ierr);
206     ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
207     ierr = VecGetArray(coordinatesNew, &coordsNew);CHKERRQ(ierr);
208     ierr = PetscSectionGetChart(csectionNew, &pStart, &pEnd);CHKERRQ(ierr);
209     ierr = ISGetIndices(perm, &pperm);CHKERRQ(ierr);
210     for (p = pStart; p < pEnd; ++p) {
211       PetscInt dof, off, offNew, d;
212 
213       ierr = PetscSectionGetDof(csectionNew, p, &dof);CHKERRQ(ierr);
214       ierr = PetscSectionGetOffset(csection, p, &off);CHKERRQ(ierr);
215       ierr = PetscSectionGetOffset(csectionNew, pperm[p], &offNew);CHKERRQ(ierr);
216       for (d = 0; d < dof; ++d) coordsNew[offNew+d] = coords[off+d];
217     }
218     ierr = ISRestoreIndices(perm, &pperm);CHKERRQ(ierr);
219     ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
220     ierr = VecRestoreArray(coordinatesNew, &coordsNew);CHKERRQ(ierr);
221     ierr = DMGetCoordinateDM(*pdm, &cdmNew);CHKERRQ(ierr);
222     ierr = DMSetDefaultSection(cdmNew, csectionNew);CHKERRQ(ierr);
223     ierr = DMSetCoordinatesLocal(*pdm, coordinatesNew);CHKERRQ(ierr);
224     ierr = PetscSectionDestroy(&csectionNew);CHKERRQ(ierr);
225     ierr = VecDestroy(&coordinatesNew);CHKERRQ(ierr);
226   }
227   /* Reorder labels */
228   {
229     PetscInt numLabels, l;
230     DMLabel  label, labelNew;
231 
232     ierr = DMGetNumLabels(dm, &numLabels);CHKERRQ(ierr);
233     for (l = numLabels-1; l >= 0; --l) {
234       ierr = DMGetLabelByNum(dm, l, &label);CHKERRQ(ierr);
235       ierr = DMLabelPermute(label, perm, &labelNew);CHKERRQ(ierr);
236       ierr = DMAddLabel(*pdm, labelNew);CHKERRQ(ierr);
237     }
238     if (plex->subpointMap) {ierr = DMLabelPermute(plex->subpointMap, perm, &plexNew->subpointMap);CHKERRQ(ierr);}
239   }
240   /* Reorder topology */
241   {
242     const PetscInt *pperm;
243     PetscInt        maxConeSize, maxSupportSize, n, pStart, pEnd, p;
244 
245     ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr);
246     plexNew->maxConeSize    = maxConeSize;
247     plexNew->maxSupportSize = maxSupportSize;
248     ierr = PetscSectionDestroy(&plexNew->coneSection);CHKERRQ(ierr);
249     ierr = PetscSectionPermute(plex->coneSection, perm, &plexNew->coneSection);CHKERRQ(ierr);
250     ierr = PetscSectionGetStorageSize(plexNew->coneSection, &n);CHKERRQ(ierr);
251     ierr = PetscMalloc1(n, &plexNew->cones);CHKERRQ(ierr);
252     ierr = PetscMalloc1(n, &plexNew->coneOrientations);CHKERRQ(ierr);
253     ierr = ISGetIndices(perm, &pperm);CHKERRQ(ierr);
254     ierr = PetscSectionGetChart(plex->coneSection, &pStart, &pEnd);CHKERRQ(ierr);
255     for (p = pStart; p < pEnd; ++p) {
256       PetscInt dof, off, offNew, d;
257 
258       ierr = PetscSectionGetDof(plexNew->coneSection, pperm[p], &dof);CHKERRQ(ierr);
259       ierr = PetscSectionGetOffset(plex->coneSection, p, &off);CHKERRQ(ierr);
260       ierr = PetscSectionGetOffset(plexNew->coneSection, pperm[p], &offNew);CHKERRQ(ierr);
261       for (d = 0; d < dof; ++d) {
262         plexNew->cones[offNew+d]            = pperm[plex->cones[off+d]];
263         plexNew->coneOrientations[offNew+d] = plex->coneOrientations[off+d];
264       }
265     }
266     ierr = PetscSectionDestroy(&plexNew->supportSection);CHKERRQ(ierr);
267     ierr = PetscSectionPermute(plex->supportSection, perm, &plexNew->supportSection);CHKERRQ(ierr);
268     ierr = PetscSectionGetStorageSize(plexNew->supportSection, &n);CHKERRQ(ierr);
269     ierr = PetscMalloc1(n, &plexNew->supports);CHKERRQ(ierr);
270     ierr = PetscSectionGetChart(plex->supportSection, &pStart, &pEnd);CHKERRQ(ierr);
271     for (p = pStart; p < pEnd; ++p) {
272       PetscInt dof, off, offNew, d;
273 
274       ierr = PetscSectionGetDof(plexNew->supportSection, pperm[p], &dof);CHKERRQ(ierr);
275       ierr = PetscSectionGetOffset(plex->supportSection, p, &off);CHKERRQ(ierr);
276       ierr = PetscSectionGetOffset(plexNew->supportSection, pperm[p], &offNew);CHKERRQ(ierr);
277       for (d = 0; d < dof; ++d) {
278         plexNew->supports[offNew+d] = pperm[plex->supports[off+d]];
279       }
280     }
281     ierr = ISRestoreIndices(perm, &pperm);CHKERRQ(ierr);
282   }
283   PetscFunctionReturn(0);
284 }
285