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