xref: /petsc/src/dm/impls/plex/plexreorder.c (revision f1f2ae845fd5815e1ad61e968c2c51b97d3840d2)
1 #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2 #include <petsc/private/matorderimpl.h> /*I      "petscmat.h"      I*/
3 
4 static PetscErrorCode DMPlexCreateOrderingClosure_Static(DM dm, PetscInt numPoints, const PetscInt pperm[], PetscInt **clperm, PetscInt **invclperm)
5 {
6   PetscInt      *perm, *iperm;
7   PetscInt       depth, d, pStart, pEnd, fStart, fMax, fEnd, p;
8 
9   PetscFunctionBegin;
10   PetscCall(DMPlexGetDepth(dm, &depth));
11   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
12   PetscCall(PetscMalloc1(pEnd-pStart,&perm));
13   PetscCall(PetscMalloc1(pEnd-pStart,&iperm));
14   for (p = pStart; p < pEnd; ++p) iperm[p] = -1;
15   for (d = depth; d > 0; --d) {
16     PetscCall(DMPlexGetDepthStratum(dm, d,   &pStart, &pEnd));
17     PetscCall(DMPlexGetDepthStratum(dm, d-1, &fStart, &fEnd));
18     fMax = fStart;
19     for (p = pStart; p < pEnd; ++p) {
20       const PetscInt *cone;
21       PetscInt        point, coneSize, c;
22 
23       if (d == depth) {
24         perm[p]         = pperm[p];
25         iperm[pperm[p]] = p;
26       }
27       point = perm[p];
28       PetscCall(DMPlexGetConeSize(dm, point, &coneSize));
29       PetscCall(DMPlexGetCone(dm, point, &cone));
30       for (c = 0; c < coneSize; ++c) {
31         const PetscInt oldc = cone[c];
32         const PetscInt newc = iperm[oldc];
33 
34         if (newc < 0) {
35           perm[fMax]  = oldc;
36           iperm[oldc] = fMax++;
37         }
38       }
39     }
40     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);
41   }
42   *clperm    = perm;
43   *invclperm = iperm;
44   PetscFunctionReturn(0);
45 }
46 
47 /*@
48   DMPlexGetOrdering - Calculate a reordering of the mesh
49 
50   Collective on dm
51 
52   Input Parameters:
53 + dm - The DMPlex object
54 . otype - type of reordering, one of the following:
55 $     MATORDERINGNATURAL - Natural
56 $     MATORDERINGND - Nested Dissection
57 $     MATORDERING1WD - One-way Dissection
58 $     MATORDERINGRCM - Reverse Cuthill-McKee
59 $     MATORDERINGQMD - Quotient Minimum Degree
60 - label - [Optional] Label used to segregate ordering into sets, or NULL
61 
62   Output Parameter:
63 . perm - The point permutation as an IS, perm[old point number] = new point number
64 
65   Note: The label is used to group sets of points together by label value. This makes it easy to reorder a mesh which
66   has different types of cells, and then loop over each set of reordered cells for assembly.
67 
68   Level: intermediate
69 
70 .seealso: `DMPlexPermute()`, `MatGetOrdering()`
71 @*/
72 PetscErrorCode DMPlexGetOrdering(DM dm, MatOrderingType otype, DMLabel label, IS *perm)
73 {
74   PetscInt       numCells = 0;
75   PetscInt      *start = NULL, *adjacency = NULL, *cperm, *clperm = NULL, *invclperm = NULL, *mask, *xls, pStart, pEnd, c, i;
76 
77   PetscFunctionBegin;
78   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
79   PetscValidPointer(perm, 4);
80   PetscCall(DMPlexCreateNeighborCSR(dm, 0, &numCells, &start, &adjacency));
81   PetscCall(PetscMalloc3(numCells,&cperm,numCells,&mask,numCells*2,&xls));
82   if (numCells) {
83     /* Shift for Fortran numbering */
84     for (i = 0; i < start[numCells]; ++i) ++adjacency[i];
85     for (i = 0; i <= numCells; ++i)       ++start[i];
86     PetscCall(SPARSEPACKgenrcm(&numCells, start, adjacency, cperm, mask, xls));
87   }
88   PetscCall(PetscFree(start));
89   PetscCall(PetscFree(adjacency));
90   /* Shift for Fortran numbering */
91   for (c = 0; c < numCells; ++c) --cperm[c];
92   /* Segregate */
93   if (label) {
94     IS              valueIS;
95     const PetscInt *values;
96     PetscInt        numValues, numPoints = 0;
97     PetscInt       *sperm, *vsize, *voff, v;
98 
99     PetscCall(DMLabelGetValueIS(label, &valueIS));
100     PetscCall(ISSort(valueIS));
101     PetscCall(ISGetLocalSize(valueIS, &numValues));
102     PetscCall(ISGetIndices(valueIS, &values));
103     PetscCall(PetscCalloc3(numCells,&sperm,numValues,&vsize,numValues+1,&voff));
104     for (v = 0; v < numValues; ++v) {
105       PetscCall(DMLabelGetStratumSize(label, values[v], &vsize[v]));
106       if (v < numValues-1) voff[v+2] += vsize[v] + voff[v+1];
107       numPoints += vsize[v];
108     }
109     PetscCheck(numPoints == numCells, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label only covers %" PetscInt_FMT " cells != %" PetscInt_FMT " total cells", numPoints, numCells);
110     for (c = 0; c < numCells; ++c) {
111       const PetscInt oldc = cperm[c];
112       PetscInt       val, vloc;
113 
114       PetscCall(DMLabelGetValue(label, oldc, &val));
115       PetscCheck(val != -1,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cell %" PetscInt_FMT " not present in label", oldc);
116       PetscCall(PetscFindInt(val, numValues, values, &vloc));
117       PetscCheck(vloc >= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Value %" PetscInt_FMT " not present label", val);
118       sperm[voff[vloc+1]++] = oldc;
119     }
120     for (v = 0; v < numValues; ++v) {
121       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]);
122     }
123     PetscCall(ISRestoreIndices(valueIS, &values));
124     PetscCall(ISDestroy(&valueIS));
125     PetscCall(PetscArraycpy(cperm, sperm, numCells));
126     PetscCall(PetscFree3(sperm, vsize, voff));
127   }
128   /* Construct closure */
129   PetscCall(DMPlexCreateOrderingClosure_Static(dm, numCells, cperm, &clperm, &invclperm));
130   PetscCall(PetscFree3(cperm,mask,xls));
131   PetscCall(PetscFree(clperm));
132   /* Invert permutation */
133   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
134   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject) dm), pEnd-pStart, invclperm, PETSC_OWN_POINTER, perm));
135   PetscFunctionReturn(0);
136 }
137 
138 /*@
139   DMPlexGetOrdering1D - Reorder the vertices so that the mesh is in a line
140 
141   Collective on dm
142 
143   Input Parameter:
144 . dm - The DMPlex object
145 
146   Output Parameter:
147 . perm - The point permutation as an IS, perm[old point number] = new point number
148 
149   Level: intermediate
150 
151 .seealso: `DMPlexGetOrdering()`, `DMPlexPermute()`, `MatGetOrdering()`
152 @*/
153 PetscErrorCode DMPlexGetOrdering1D(DM dm, IS *perm)
154 {
155   PetscInt       *points;
156   const PetscInt *support, *cone;
157   PetscInt        dim, pStart, pEnd, cStart, cEnd, c, vStart, vEnd, v, suppSize, lastCell = 0;
158 
159   PetscFunctionBegin;
160   PetscCall(DMGetDimension(dm, &dim));
161   PetscCheck(dim == 1, PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Input mesh must be one dimensional, not %" PetscInt_FMT, dim);
162   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
163   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
164   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
165   PetscCall(PetscMalloc1(pEnd-pStart, &points));
166   for (c = cStart; c < cEnd; ++c) points[c] = c;
167   for (v = vStart; v < vEnd; ++v) points[v] = v;
168   for (v = vStart; v < vEnd; ++v) {
169     PetscCall(DMPlexGetSupportSize(dm, v, &suppSize));
170     PetscCall(DMPlexGetSupport(dm, v, &support));
171     if (suppSize == 1) {lastCell = support[0]; break;}
172   }
173   if (v < vEnd) {
174     PetscInt pos = cEnd;
175 
176     points[v] = pos++;
177     while (lastCell >= cStart) {
178       PetscCall(DMPlexGetCone(dm, lastCell, &cone));
179       if (cone[0] == v) v = cone[1];
180       else              v = cone[0];
181       PetscCall(DMPlexGetSupport(dm, v, &support));
182       PetscCall(DMPlexGetSupportSize(dm, v, &suppSize));
183       if (suppSize == 1) {lastCell = -1;}
184       else {
185         if (support[0] == lastCell) lastCell = support[1];
186         else                        lastCell = support[0];
187       }
188       points[v] = pos++;
189     }
190     PetscCheck(pos == pEnd, PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Last vertex was %" PetscInt_FMT ", not %" PetscInt_FMT, pos, pEnd);
191   }
192   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject) dm), pEnd-pStart, points, PETSC_OWN_POINTER, perm));
193   PetscFunctionReturn(0);
194 }
195 
196 static PetscErrorCode DMPlexRemapCoordinates_Private(IS perm, PetscSection cs, Vec coordinates, PetscSection *csNew, Vec *coordinatesNew)
197 {
198   PetscScalar    *coords, *coordsNew;
199   const PetscInt *pperm;
200   PetscInt        pStart, pEnd, p;
201   const char     *name;
202 
203   PetscFunctionBegin;
204   PetscCall(PetscSectionPermute(cs, perm, csNew));
205   PetscCall(VecDuplicate(coordinates, coordinatesNew));
206   PetscCall(PetscObjectGetName((PetscObject) coordinates, &name));
207   PetscCall(PetscObjectSetName((PetscObject) *coordinatesNew, name));
208   PetscCall(VecGetArray(coordinates, &coords));
209   PetscCall(VecGetArray(*coordinatesNew, &coordsNew));
210   PetscCall(PetscSectionGetChart(*csNew, &pStart, &pEnd));
211   PetscCall(ISGetIndices(perm, &pperm));
212   for (p = pStart; p < pEnd; ++p) {
213     PetscInt dof, off, offNew, d;
214 
215     PetscCall(PetscSectionGetDof(*csNew, p, &dof));
216     PetscCall(PetscSectionGetOffset(cs, p, &off));
217     PetscCall(PetscSectionGetOffset(*csNew, pperm[p], &offNew));
218     for (d = 0; d < dof; ++d) coordsNew[offNew+d] = coords[off+d];
219   }
220   PetscCall(ISRestoreIndices(perm, &pperm));
221   PetscCall(VecRestoreArray(coordinates, &coords));
222   PetscCall(VecRestoreArray(*coordinatesNew, &coordsNew));
223   PetscFunctionReturn(0);
224 }
225 
226 /*@
227   DMPlexPermute - Reorder the mesh according to the input permutation
228 
229   Collective on dm
230 
231   Input Parameters:
232 + dm - The DMPlex object
233 - perm - The point permutation, perm[old point number] = new point number
234 
235   Output Parameter:
236 . pdm - The permuted DM
237 
238   Level: intermediate
239 
240 .seealso: `MatPermute()`
241 @*/
242 PetscErrorCode DMPlexPermute(DM dm, IS perm, DM *pdm)
243 {
244   DM_Plex       *plex = (DM_Plex *) dm->data, *plexNew;
245   PetscInt       dim, cdim;
246   const char    *name;
247 
248   PetscFunctionBegin;
249   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
250   PetscValidHeaderSpecific(perm, IS_CLASSID, 2);
251   PetscValidPointer(pdm, 3);
252   PetscCall(DMCreate(PetscObjectComm((PetscObject) dm), pdm));
253   PetscCall(DMSetType(*pdm, DMPLEX));
254   PetscCall(PetscObjectGetName((PetscObject) dm, &name));
255   PetscCall(PetscObjectSetName((PetscObject) *pdm, name));
256   PetscCall(DMGetDimension(dm, &dim));
257   PetscCall(DMSetDimension(*pdm, dim));
258   PetscCall(DMGetCoordinateDim(dm, &cdim));
259   PetscCall(DMSetCoordinateDim(*pdm, cdim));
260   PetscCall(DMCopyDisc(dm, *pdm));
261   if (dm->localSection) {
262     PetscSection section, sectionNew;
263 
264     PetscCall(DMGetLocalSection(dm, &section));
265     PetscCall(PetscSectionPermute(section, perm, &sectionNew));
266     PetscCall(DMSetLocalSection(*pdm, sectionNew));
267     PetscCall(PetscSectionDestroy(&sectionNew));
268   }
269   plexNew = (DM_Plex *) (*pdm)->data;
270   /* Ignore ltogmap, ltogmapb */
271   /* Ignore sf, sectionSF */
272   /* Ignore globalVertexNumbers, globalCellNumbers */
273   /* Reorder labels */
274   {
275     PetscInt numLabels, l;
276     DMLabel  label, labelNew;
277 
278     PetscCall(DMGetNumLabels(dm, &numLabels));
279     for (l = 0; l < numLabels; ++l) {
280       PetscCall(DMGetLabelByNum(dm, l, &label));
281       PetscCall(DMLabelPermute(label, perm, &labelNew));
282       PetscCall(DMAddLabel(*pdm, labelNew));
283       PetscCall(DMLabelDestroy(&labelNew));
284     }
285     PetscCall(DMGetLabel(*pdm, "depth", &(*pdm)->depthLabel));
286     if (plex->subpointMap) PetscCall(DMLabelPermute(plex->subpointMap, perm, &plexNew->subpointMap));
287   }
288   /* Reorder topology */
289   {
290     const PetscInt *pperm;
291     PetscInt        n, pStart, pEnd, p;
292 
293     PetscCall(PetscSectionDestroy(&plexNew->coneSection));
294     PetscCall(PetscSectionPermute(plex->coneSection, perm, &plexNew->coneSection));
295     PetscCall(PetscSectionGetStorageSize(plexNew->coneSection, &n));
296     PetscCall(PetscMalloc1(n, &plexNew->cones));
297     PetscCall(PetscMalloc1(n, &plexNew->coneOrientations));
298     PetscCall(ISGetIndices(perm, &pperm));
299     PetscCall(PetscSectionGetChart(plex->coneSection, &pStart, &pEnd));
300     for (p = pStart; p < pEnd; ++p) {
301       PetscInt dof, off, offNew, d;
302 
303       PetscCall(PetscSectionGetDof(plexNew->coneSection, pperm[p], &dof));
304       PetscCall(PetscSectionGetOffset(plex->coneSection, p, &off));
305       PetscCall(PetscSectionGetOffset(plexNew->coneSection, pperm[p], &offNew));
306       for (d = 0; d < dof; ++d) {
307         plexNew->cones[offNew+d]            = pperm[plex->cones[off+d]];
308         plexNew->coneOrientations[offNew+d] = plex->coneOrientations[off+d];
309       }
310     }
311     PetscCall(PetscSectionDestroy(&plexNew->supportSection));
312     PetscCall(PetscSectionPermute(plex->supportSection, perm, &plexNew->supportSection));
313     PetscCall(PetscSectionGetStorageSize(plexNew->supportSection, &n));
314     PetscCall(PetscMalloc1(n, &plexNew->supports));
315     PetscCall(PetscSectionGetChart(plex->supportSection, &pStart, &pEnd));
316     for (p = pStart; p < pEnd; ++p) {
317       PetscInt dof, off, offNew, d;
318 
319       PetscCall(PetscSectionGetDof(plexNew->supportSection, pperm[p], &dof));
320       PetscCall(PetscSectionGetOffset(plex->supportSection, p, &off));
321       PetscCall(PetscSectionGetOffset(plexNew->supportSection, pperm[p], &offNew));
322       for (d = 0; d < dof; ++d) {
323         plexNew->supports[offNew+d] = pperm[plex->supports[off+d]];
324       }
325     }
326     PetscCall(ISRestoreIndices(perm, &pperm));
327   }
328   /* Remap coordinates */
329   {
330     DM           cdm, cdmNew;
331     PetscSection cs, csNew;
332     Vec          coordinates, coordinatesNew;
333 
334     PetscCall(DMGetCoordinateSection(dm, &cs));
335     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
336     PetscCall(DMPlexRemapCoordinates_Private(perm, cs, coordinates, &csNew, &coordinatesNew));
337     PetscCall(DMSetCoordinateSection(*pdm, PETSC_DETERMINE, csNew));
338     PetscCall(DMSetCoordinatesLocal(*pdm, coordinatesNew));
339     PetscCall(PetscSectionDestroy(&csNew));
340     PetscCall(VecDestroy(&coordinatesNew));
341 
342     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
343     if (cdm) {
344       PetscCall(DMGetCoordinateDM(*pdm, &cdm));
345       PetscCall(DMClone(cdm, &cdmNew));
346       PetscCall(DMSetCellCoordinateDM(*pdm, cdmNew));
347       PetscCall(DMDestroy(&cdmNew));
348       PetscCall(DMGetCellCoordinateSection(dm, &cs));
349       PetscCall(DMGetCellCoordinatesLocal(dm, &coordinates));
350       PetscCall(DMPlexRemapCoordinates_Private(perm, cs, coordinates, &csNew, &coordinatesNew));
351       PetscCall(DMSetCellCoordinateSection(*pdm, PETSC_DETERMINE, csNew));
352       PetscCall(DMSetCellCoordinatesLocal(*pdm, coordinatesNew));
353       PetscCall(PetscSectionDestroy(&csNew));
354       PetscCall(VecDestroy(&coordinatesNew));
355     }
356   }
357   PetscCall(DMPlexCopy_Internal(dm, PETSC_TRUE, PETSC_TRUE, *pdm));
358   (*pdm)->setupcalled = PETSC_TRUE;
359   PetscFunctionReturn(0);
360 }
361 
362 PetscErrorCode DMPlexReorderSetDefault_Plex(DM dm, DMPlexReorderDefaultFlag reorder)
363 {
364   DM_Plex *mesh = (DM_Plex *) dm->data;
365 
366   PetscFunctionBegin;
367   mesh->reorderDefault = reorder;
368   PetscFunctionReturn(0);
369 }
370 
371 /*@
372   DMPlexReorderSetDefault - Set flag indicating whether the DM should be reordered by default
373 
374   Logically collective
375 
376   Input Parameters:
377 + dm        - The DM
378 - reorder   - Flag for reordering
379 
380   Level: intermediate
381 
382 .seealso: `DMPlexReorderGetDefault()`
383 @*/
384 PetscErrorCode DMPlexReorderSetDefault(DM dm, DMPlexReorderDefaultFlag reorder)
385 {
386   PetscFunctionBegin;
387   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
388   PetscTryMethod(dm,"DMPlexReorderSetDefault_C",(DM,DMPlexReorderDefaultFlag),(dm,reorder));
389   PetscFunctionReturn(0);
390 }
391 
392 PetscErrorCode DMPlexReorderGetDefault_Plex(DM dm, DMPlexReorderDefaultFlag *reorder)
393 {
394   DM_Plex *mesh = (DM_Plex *) dm->data;
395 
396   PetscFunctionBegin;
397   *reorder = mesh->reorderDefault;
398   PetscFunctionReturn(0);
399 }
400 
401 /*@
402   DMPlexReorderGetDefault - Get flag indicating whether the DM should be reordered by default
403 
404   Not collective
405 
406   Input Parameter:
407 . dm      - The DM
408 
409   Output Parameter:
410 . reorder - Flag for reordering
411 
412   Level: intermediate
413 
414 .seealso: `DMPlexReorderSetDefault()`
415 @*/
416 PetscErrorCode DMPlexReorderGetDefault(DM dm, DMPlexReorderDefaultFlag *reorder)
417 {
418   PetscFunctionBegin;
419   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
420   PetscValidPointer(reorder, 2);
421   PetscUseMethod(dm,"DMPlexReorderGetDefault_C",(DM,DMPlexReorderDefaultFlag*),(dm,reorder));
422   PetscFunctionReturn(0);
423 }
424