xref: /petsc/src/dm/impls/plex/plexreorder.c (revision 5520554388890bd89a1c1cf7870aedf4e71d512f)
1 #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2 #include <petsc/private/matorderimpl.h> /*I      "petscmat.h"      I*/
3 #include <petsc/private/dmlabelimpl.h>
4 
5 static PetscErrorCode DMPlexCreateOrderingClosure_Static(DM dm, PetscInt numPoints, const PetscInt pperm[], PetscInt **clperm, PetscInt **invclperm) {
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   PetscInt  numCells = 0;
74   PetscInt *start = NULL, *adjacency = NULL, *cperm, *clperm = NULL, *invclperm = NULL, *mask, *xls, pStart, pEnd, c, i;
75 
76   PetscFunctionBegin;
77   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
78   PetscValidPointer(perm, 4);
79   PetscCall(DMPlexCreateNeighborCSR(dm, 0, &numCells, &start, &adjacency));
80   PetscCall(PetscMalloc3(numCells, &cperm, numCells, &mask, numCells * 2, &xls));
81   if (numCells) {
82     /* Shift for Fortran numbering */
83     for (i = 0; i < start[numCells]; ++i) ++adjacency[i];
84     for (i = 0; i <= numCells; ++i) ++start[i];
85     PetscCall(SPARSEPACKgenrcm(&numCells, start, adjacency, cperm, mask, xls));
86   }
87   PetscCall(PetscFree(start));
88   PetscCall(PetscFree(adjacency));
89   /* Shift for Fortran numbering */
90   for (c = 0; c < numCells; ++c) --cperm[c];
91   /* Segregate */
92   if (label) {
93     IS              valueIS;
94     const PetscInt *valuesTmp;
95     PetscInt       *values;
96     PetscInt        numValues, numPoints = 0;
97     PetscInt       *sperm, *vsize, *voff, v;
98 
99     // Can't directly sort the valueIS, since it is a view into the DMLabel
100     PetscCall(DMLabelGetValueIS(label, &valueIS));
101     PetscCall(ISGetLocalSize(valueIS, &numValues));
102     PetscCall(ISGetIndices(valueIS, &valuesTmp));
103     PetscCall(PetscCalloc4(numCells, &sperm, numValues, &values, numValues, &vsize, numValues + 1, &voff));
104     PetscCall(PetscArraycpy(values, valuesTmp, numValues));
105     PetscCall(PetscSortInt(numValues, values));
106     PetscCall(ISRestoreIndices(valueIS, &valuesTmp));
107     PetscCall(ISDestroy(&valueIS));
108     for (v = 0; v < numValues; ++v) {
109       PetscCall(DMLabelGetStratumSize(label, values[v], &vsize[v]));
110       if (v < numValues - 1) voff[v + 2] += vsize[v] + voff[v + 1];
111       numPoints += vsize[v];
112     }
113     PetscCheck(numPoints == numCells, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label only covers %" PetscInt_FMT " cells != %" PetscInt_FMT " total cells", numPoints, numCells);
114     for (c = 0; c < numCells; ++c) {
115       const PetscInt oldc = cperm[c];
116       PetscInt       val, vloc;
117 
118       PetscCall(DMLabelGetValue(label, oldc, &val));
119       PetscCheck(val != -1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cell %" PetscInt_FMT " not present in label", oldc);
120       PetscCall(PetscFindInt(val, numValues, values, &vloc));
121       PetscCheck(vloc >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Value %" PetscInt_FMT " not present label", val);
122       sperm[voff[vloc + 1]++] = oldc;
123     }
124     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]);
125     PetscCall(PetscArraycpy(cperm, sperm, numCells));
126     PetscCall(PetscFree4(sperm, values, 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   PetscInt       *points;
155   const PetscInt *support, *cone;
156   PetscInt        dim, pStart, pEnd, cStart, cEnd, c, vStart, vEnd, v, suppSize, lastCell = 0;
157 
158   PetscFunctionBegin;
159   PetscCall(DMGetDimension(dm, &dim));
160   PetscCheck(dim == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Input mesh must be one dimensional, not %" PetscInt_FMT, dim);
161   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
162   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
163   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
164   PetscCall(PetscMalloc1(pEnd - pStart, &points));
165   for (c = cStart; c < cEnd; ++c) points[c] = c;
166   for (v = vStart; v < vEnd; ++v) points[v] = v;
167   for (v = vStart; v < vEnd; ++v) {
168     PetscCall(DMPlexGetSupportSize(dm, v, &suppSize));
169     PetscCall(DMPlexGetSupport(dm, v, &support));
170     if (suppSize == 1) {
171       lastCell = support[0];
172       break;
173     }
174   }
175   if (v < vEnd) {
176     PetscInt pos = cEnd;
177 
178     points[v] = pos++;
179     while (lastCell >= cStart) {
180       PetscCall(DMPlexGetCone(dm, lastCell, &cone));
181       if (cone[0] == v) v = cone[1];
182       else v = cone[0];
183       PetscCall(DMPlexGetSupport(dm, v, &support));
184       PetscCall(DMPlexGetSupportSize(dm, v, &suppSize));
185       if (suppSize == 1) {
186         lastCell = -1;
187       } else {
188         if (support[0] == lastCell) lastCell = support[1];
189         else lastCell = support[0];
190       }
191       points[v] = pos++;
192     }
193     PetscCheck(pos == pEnd, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Last vertex was %" PetscInt_FMT ", not %" PetscInt_FMT, pos, pEnd);
194   }
195   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), pEnd - pStart, points, PETSC_OWN_POINTER, perm));
196   PetscFunctionReturn(0);
197 }
198 
199 static PetscErrorCode DMPlexRemapCoordinates_Private(IS perm, PetscSection cs, Vec coordinates, PetscSection *csNew, Vec *coordinatesNew) {
200   PetscScalar    *coords, *coordsNew;
201   const PetscInt *pperm;
202   PetscInt        pStart, pEnd, p;
203   const char     *name;
204 
205   PetscFunctionBegin;
206   PetscCall(PetscSectionPermute(cs, perm, csNew));
207   PetscCall(VecDuplicate(coordinates, coordinatesNew));
208   PetscCall(PetscObjectGetName((PetscObject)coordinates, &name));
209   PetscCall(PetscObjectSetName((PetscObject)*coordinatesNew, name));
210   PetscCall(VecGetArray(coordinates, &coords));
211   PetscCall(VecGetArray(*coordinatesNew, &coordsNew));
212   PetscCall(PetscSectionGetChart(*csNew, &pStart, &pEnd));
213   PetscCall(ISGetIndices(perm, &pperm));
214   for (p = pStart; p < pEnd; ++p) {
215     PetscInt dof, off, offNew, d;
216 
217     PetscCall(PetscSectionGetDof(*csNew, p, &dof));
218     PetscCall(PetscSectionGetOffset(cs, p, &off));
219     PetscCall(PetscSectionGetOffset(*csNew, pperm[p], &offNew));
220     for (d = 0; d < dof; ++d) coordsNew[offNew + d] = coords[off + d];
221   }
222   PetscCall(ISRestoreIndices(perm, &pperm));
223   PetscCall(VecRestoreArray(coordinates, &coords));
224   PetscCall(VecRestoreArray(*coordinatesNew, &coordsNew));
225   PetscFunctionReturn(0);
226 }
227 
228 /*@
229   DMPlexPermute - Reorder the mesh according to the input permutation
230 
231   Collective on dm
232 
233   Input Parameters:
234 + dm - The DMPlex object
235 - perm - The point permutation, perm[old point number] = new point number
236 
237   Output Parameter:
238 . pdm - The permuted DM
239 
240   Level: intermediate
241 
242 .seealso: `MatPermute()`
243 @*/
244 PetscErrorCode DMPlexPermute(DM dm, IS perm, DM *pdm) {
245   DM_Plex    *plex = (DM_Plex *)dm->data, *plexNew;
246   PetscInt    dim, cdim;
247   const char *name;
248 
249   PetscFunctionBegin;
250   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
251   PetscValidHeaderSpecific(perm, IS_CLASSID, 2);
252   PetscValidPointer(pdm, 3);
253   PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), pdm));
254   PetscCall(DMSetType(*pdm, DMPLEX));
255   PetscCall(PetscObjectGetName((PetscObject)dm, &name));
256   PetscCall(PetscObjectSetName((PetscObject)*pdm, name));
257   PetscCall(DMGetDimension(dm, &dim));
258   PetscCall(DMSetDimension(*pdm, dim));
259   PetscCall(DMGetCoordinateDim(dm, &cdim));
260   PetscCall(DMSetCoordinateDim(*pdm, cdim));
261   PetscCall(DMCopyDisc(dm, *pdm));
262   if (dm->localSection) {
263     PetscSection section, sectionNew;
264 
265     PetscCall(DMGetLocalSection(dm, &section));
266     PetscCall(PetscSectionPermute(section, perm, &sectionNew));
267     PetscCall(DMSetLocalSection(*pdm, sectionNew));
268     PetscCall(PetscSectionDestroy(&sectionNew));
269   }
270   plexNew = (DM_Plex *)(*pdm)->data;
271   /* Ignore ltogmap, ltogmapb */
272   /* Ignore sf, sectionSF */
273   /* Ignore globalVertexNumbers, globalCellNumbers */
274   /* Reorder labels */
275   {
276     PetscInt numLabels, l;
277     DMLabel  label, labelNew;
278 
279     PetscCall(DMGetNumLabels(dm, &numLabels));
280     for (l = 0; l < numLabels; ++l) {
281       PetscCall(DMGetLabelByNum(dm, l, &label));
282       PetscCall(DMLabelPermute(label, perm, &labelNew));
283       PetscCall(DMAddLabel(*pdm, labelNew));
284       PetscCall(DMLabelDestroy(&labelNew));
285     }
286     PetscCall(DMGetLabel(*pdm, "depth", &(*pdm)->depthLabel));
287     if (plex->subpointMap) PetscCall(DMLabelPermute(plex->subpointMap, perm, &plexNew->subpointMap));
288   }
289   if ((*pdm)->celltypeLabel) {
290     DMLabel ctLabel;
291 
292     // Reset label for fast lookup
293     PetscCall(DMPlexGetCellTypeLabel(*pdm, &ctLabel));
294     PetscCall(DMLabelMakeAllInvalid_Internal(ctLabel));
295   }
296   /* Reorder topology */
297   {
298     const PetscInt *pperm;
299     PetscInt        n, pStart, pEnd, p;
300 
301     PetscCall(PetscSectionDestroy(&plexNew->coneSection));
302     PetscCall(PetscSectionPermute(plex->coneSection, perm, &plexNew->coneSection));
303     PetscCall(PetscSectionGetStorageSize(plexNew->coneSection, &n));
304     PetscCall(PetscMalloc1(n, &plexNew->cones));
305     PetscCall(PetscMalloc1(n, &plexNew->coneOrientations));
306     PetscCall(ISGetIndices(perm, &pperm));
307     PetscCall(PetscSectionGetChart(plex->coneSection, &pStart, &pEnd));
308     for (p = pStart; p < pEnd; ++p) {
309       PetscInt dof, off, offNew, d;
310 
311       PetscCall(PetscSectionGetDof(plexNew->coneSection, pperm[p], &dof));
312       PetscCall(PetscSectionGetOffset(plex->coneSection, p, &off));
313       PetscCall(PetscSectionGetOffset(plexNew->coneSection, pperm[p], &offNew));
314       for (d = 0; d < dof; ++d) {
315         plexNew->cones[offNew + d]            = pperm[plex->cones[off + d]];
316         plexNew->coneOrientations[offNew + d] = plex->coneOrientations[off + d];
317       }
318     }
319     PetscCall(PetscSectionDestroy(&plexNew->supportSection));
320     PetscCall(PetscSectionPermute(plex->supportSection, perm, &plexNew->supportSection));
321     PetscCall(PetscSectionGetStorageSize(plexNew->supportSection, &n));
322     PetscCall(PetscMalloc1(n, &plexNew->supports));
323     PetscCall(PetscSectionGetChart(plex->supportSection, &pStart, &pEnd));
324     for (p = pStart; p < pEnd; ++p) {
325       PetscInt dof, off, offNew, d;
326 
327       PetscCall(PetscSectionGetDof(plexNew->supportSection, pperm[p], &dof));
328       PetscCall(PetscSectionGetOffset(plex->supportSection, p, &off));
329       PetscCall(PetscSectionGetOffset(plexNew->supportSection, pperm[p], &offNew));
330       for (d = 0; d < dof; ++d) plexNew->supports[offNew + d] = pperm[plex->supports[off + d]];
331     }
332     PetscCall(ISRestoreIndices(perm, &pperm));
333   }
334   /* Remap coordinates */
335   {
336     DM           cdm, cdmNew;
337     PetscSection cs, csNew;
338     Vec          coordinates, coordinatesNew;
339 
340     PetscCall(DMGetCoordinateSection(dm, &cs));
341     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
342     PetscCall(DMPlexRemapCoordinates_Private(perm, cs, coordinates, &csNew, &coordinatesNew));
343     PetscCall(DMSetCoordinateSection(*pdm, PETSC_DETERMINE, csNew));
344     PetscCall(DMSetCoordinatesLocal(*pdm, coordinatesNew));
345     PetscCall(PetscSectionDestroy(&csNew));
346     PetscCall(VecDestroy(&coordinatesNew));
347 
348     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
349     if (cdm) {
350       PetscCall(DMGetCoordinateDM(*pdm, &cdm));
351       PetscCall(DMClone(cdm, &cdmNew));
352       PetscCall(DMSetCellCoordinateDM(*pdm, cdmNew));
353       PetscCall(DMDestroy(&cdmNew));
354       PetscCall(DMGetCellCoordinateSection(dm, &cs));
355       PetscCall(DMGetCellCoordinatesLocal(dm, &coordinates));
356       PetscCall(DMPlexRemapCoordinates_Private(perm, cs, coordinates, &csNew, &coordinatesNew));
357       PetscCall(DMSetCellCoordinateSection(*pdm, PETSC_DETERMINE, csNew));
358       PetscCall(DMSetCellCoordinatesLocal(*pdm, coordinatesNew));
359       PetscCall(PetscSectionDestroy(&csNew));
360       PetscCall(VecDestroy(&coordinatesNew));
361     }
362   }
363   PetscCall(DMPlexCopy_Internal(dm, PETSC_TRUE, PETSC_TRUE, *pdm));
364   (*pdm)->setupcalled = PETSC_TRUE;
365   PetscFunctionReturn(0);
366 }
367 
368 PetscErrorCode DMPlexReorderSetDefault_Plex(DM dm, DMPlexReorderDefaultFlag reorder) {
369   DM_Plex *mesh = (DM_Plex *)dm->data;
370 
371   PetscFunctionBegin;
372   mesh->reorderDefault = reorder;
373   PetscFunctionReturn(0);
374 }
375 
376 /*@
377   DMPlexReorderSetDefault - Set flag indicating whether the DM should be reordered by default
378 
379   Logically collective
380 
381   Input Parameters:
382 + dm        - The DM
383 - reorder   - Flag for reordering
384 
385   Level: intermediate
386 
387 .seealso: `DMPlexReorderGetDefault()`
388 @*/
389 PetscErrorCode DMPlexReorderSetDefault(DM dm, DMPlexReorderDefaultFlag reorder) {
390   PetscFunctionBegin;
391   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
392   PetscTryMethod(dm, "DMPlexReorderSetDefault_C", (DM, DMPlexReorderDefaultFlag), (dm, reorder));
393   PetscFunctionReturn(0);
394 }
395 
396 PetscErrorCode DMPlexReorderGetDefault_Plex(DM dm, DMPlexReorderDefaultFlag *reorder) {
397   DM_Plex *mesh = (DM_Plex *)dm->data;
398 
399   PetscFunctionBegin;
400   *reorder = mesh->reorderDefault;
401   PetscFunctionReturn(0);
402 }
403 
404 /*@
405   DMPlexReorderGetDefault - Get flag indicating whether the DM should be reordered by default
406 
407   Not collective
408 
409   Input Parameter:
410 . dm      - The DM
411 
412   Output Parameter:
413 . reorder - Flag for reordering
414 
415   Level: intermediate
416 
417 .seealso: `DMPlexReorderSetDefault()`
418 @*/
419 PetscErrorCode DMPlexReorderGetDefault(DM dm, DMPlexReorderDefaultFlag *reorder) {
420   PetscFunctionBegin;
421   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
422   PetscValidPointer(reorder, 2);
423   PetscUseMethod(dm, "DMPlexReorderGetDefault_C", (DM, DMPlexReorderDefaultFlag *), (dm, reorder));
424   PetscFunctionReturn(0);
425 }
426