xref: /petsc/src/dm/impls/plex/plexsfc.c (revision 4e2e9504eff61fd362b9eabd44dfe3883745e908)
13e72e933SJed Brown #include <petsc/private/dmpleximpl.h> /*I   "petscdmplex.h"   I*/
23e72e933SJed Brown #include <petscsf.h>
33e72e933SJed Brown #include <petsc/private/hashset.h>
43e72e933SJed Brown 
53e72e933SJed Brown typedef uint64_t ZCode;
63e72e933SJed Brown 
73e72e933SJed Brown PETSC_HASH_SET(ZSet, ZCode, PetscHash_UInt64, PetscHashEqual)
83e72e933SJed Brown 
93e72e933SJed Brown typedef struct {
103e72e933SJed Brown   PetscInt i, j, k;
113e72e933SJed Brown } Ijk;
123e72e933SJed Brown 
133e72e933SJed Brown typedef struct {
143e72e933SJed Brown   Ijk         eextent;
153e72e933SJed Brown   Ijk         vextent;
163e72e933SJed Brown   PetscMPIInt comm_size;
173e72e933SJed Brown   ZCode      *zstarts;
183e72e933SJed Brown } ZLayout;
193e72e933SJed Brown 
203e72e933SJed Brown static unsigned ZCodeSplit1(ZCode z)
213e72e933SJed Brown {
223e72e933SJed Brown   z = ((z & 01001001001001001) | ((z >> 2) & 02002002002002002) | ((z >> 4) & 04004004004004004));
233e72e933SJed Brown   z = (z | (z >> 6) | (z >> 12)) & 0000000777000000777;
243e72e933SJed Brown   z = (z | (z >> 18)) & 0777777;
253e72e933SJed Brown   return (unsigned)z;
263e72e933SJed Brown }
273e72e933SJed Brown 
283e72e933SJed Brown static ZCode ZEncode1(unsigned t)
293e72e933SJed Brown {
303e72e933SJed Brown   ZCode z = t;
313e72e933SJed Brown   z       = (z | (z << 18)) & 0777000000777;
323e72e933SJed Brown   z       = (z | (z << 6) | (z << 12)) & 07007007007007007;
333e72e933SJed Brown   z       = (z | (z << 2) | (z << 4)) & 0111111111111111111;
343e72e933SJed Brown   return z;
353e72e933SJed Brown }
363e72e933SJed Brown 
373e72e933SJed Brown static Ijk ZCodeSplit(ZCode z)
383e72e933SJed Brown {
393e72e933SJed Brown   Ijk c;
403e72e933SJed Brown   c.i = ZCodeSplit1(z >> 2);
413e72e933SJed Brown   c.j = ZCodeSplit1(z >> 1);
423e72e933SJed Brown   c.k = ZCodeSplit1(z >> 0);
433e72e933SJed Brown   return c;
443e72e933SJed Brown }
453e72e933SJed Brown 
463e72e933SJed Brown static ZCode ZEncode(Ijk c)
473e72e933SJed Brown {
483e72e933SJed Brown   ZCode z = (ZEncode1(c.i) << 2) | (ZEncode1(c.j) << 1) | ZEncode1(c.k);
493e72e933SJed Brown   return z;
503e72e933SJed Brown }
513e72e933SJed Brown 
523e72e933SJed Brown static PetscBool IjkActive(Ijk extent, Ijk l)
533e72e933SJed Brown {
543e72e933SJed Brown   if (l.i < extent.i && l.j < extent.j && l.k < extent.k) return PETSC_TRUE;
553e72e933SJed Brown   return PETSC_FALSE;
563e72e933SJed Brown }
573e72e933SJed Brown 
583e72e933SJed Brown // Since element/vertex box extents are typically not equal powers of 2, Z codes that lie within the domain are not contiguous.
59c77877e3SJed Brown static PetscErrorCode ZLayoutCreate(PetscMPIInt size, const PetscInt eextent[3], const PetscInt vextent[3], ZLayout *zlayout)
603e72e933SJed Brown {
613e72e933SJed Brown   ZLayout layout;
62c77877e3SJed Brown 
63c77877e3SJed Brown   PetscFunctionBegin;
643e72e933SJed Brown   layout.eextent.i = eextent[0];
653e72e933SJed Brown   layout.eextent.j = eextent[1];
663e72e933SJed Brown   layout.eextent.k = eextent[2];
673e72e933SJed Brown   layout.vextent.i = vextent[0];
683e72e933SJed Brown   layout.vextent.j = vextent[1];
693e72e933SJed Brown   layout.vextent.k = vextent[2];
703e72e933SJed Brown   layout.comm_size = size;
71c77877e3SJed Brown   PetscCall(PetscMalloc1(size + 1, &layout.zstarts));
723e72e933SJed Brown 
733e72e933SJed Brown   PetscInt total_elems = eextent[0] * eextent[1] * eextent[2];
743e72e933SJed Brown   ZCode    z           = 0;
753e72e933SJed Brown   layout.zstarts[0]    = 0;
763e72e933SJed Brown   for (PetscMPIInt r = 0; r < size; r++) {
773e72e933SJed Brown     PetscInt elems_needed = (total_elems / size) + (total_elems % size > r), count;
783e72e933SJed Brown     for (count = 0; count < elems_needed; z++) {
793e72e933SJed Brown       Ijk loc = ZCodeSplit(z);
803e72e933SJed Brown       if (IjkActive(layout.eextent, loc)) count++;
813e72e933SJed Brown     }
823e72e933SJed Brown     // Pick up any extra vertices in the Z ordering before the next rank's first owned element.
83c77877e3SJed Brown     //
84c77877e3SJed Brown     // TODO: This leads to poorly balanced vertices when eextent is a power of 2, since all the fringe vertices end up
85c77877e3SJed Brown     // on the last rank. A possible solution is to balance the Z-order vertices independently from the cells, which will
86c77877e3SJed Brown     // result in a lot of element closures being remote. We could finish marking boundary conditions, then do a round of
87c77877e3SJed Brown     // vertex ownership smoothing (which would reorder and redistribute vertices without touching element distribution).
88c77877e3SJed Brown     // Another would be to have an analytic ownership criteria for vertices in the fringe veextent - eextent. This would
89c77877e3SJed Brown     // complicate the job of identifying an owner and its offset.
903e72e933SJed Brown     for (; z <= ZEncode(layout.vextent); z++) {
913e72e933SJed Brown       Ijk loc = ZCodeSplit(z);
923e72e933SJed Brown       if (IjkActive(layout.eextent, loc)) break;
933e72e933SJed Brown     }
943e72e933SJed Brown     layout.zstarts[r + 1] = z;
953e72e933SJed Brown   }
96c77877e3SJed Brown   *zlayout = layout;
97c77877e3SJed Brown   PetscFunctionReturn(0);
983e72e933SJed Brown }
993e72e933SJed Brown 
100*4e2e9504SJed Brown static PetscInt ZLayoutElementsOnRank(const ZLayout *layout, PetscMPIInt rank)
101*4e2e9504SJed Brown {
102*4e2e9504SJed Brown   PetscInt remote_elem = 0;
103*4e2e9504SJed Brown   for (ZCode rz = layout->zstarts[rank]; rz < layout->zstarts[rank + 1]; rz++) {
104*4e2e9504SJed Brown     Ijk loc = ZCodeSplit(rz);
105*4e2e9504SJed Brown     if (IjkActive(layout->eextent, loc)) remote_elem++;
106*4e2e9504SJed Brown   }
107*4e2e9504SJed Brown   return remote_elem;
108*4e2e9504SJed Brown }
109*4e2e9504SJed Brown 
1103e72e933SJed Brown PetscInt ZCodeFind(ZCode key, PetscInt n, const ZCode X[])
1113e72e933SJed Brown {
1123e72e933SJed Brown   PetscInt lo = 0, hi = n;
1133e72e933SJed Brown 
1143e72e933SJed Brown   if (n == 0) return -1;
1153e72e933SJed Brown   while (hi - lo > 1) {
1163e72e933SJed Brown     PetscInt mid = lo + (hi - lo) / 2;
1173e72e933SJed Brown     if (key < X[mid]) hi = mid;
1183e72e933SJed Brown     else lo = mid;
1193e72e933SJed Brown   }
1203e72e933SJed Brown   return key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
1213e72e933SJed Brown }
1223e72e933SJed Brown 
123*4e2e9504SJed Brown static PetscErrorCode DMPlexCreateBoxMesh_Tensor_SFC_Periodicity_Private(DM dm, const ZLayout *layout, const ZCode *vert_z, PetscSegBuffer per_faces, PetscSegBuffer donor_face_closure, PetscSegBuffer my_donor_faces)
124*4e2e9504SJed Brown {
125*4e2e9504SJed Brown   MPI_Comm     comm;
126*4e2e9504SJed Brown   size_t       num_faces;
127*4e2e9504SJed Brown   PetscInt     dim, *faces, vStart, vEnd;
128*4e2e9504SJed Brown   PetscMPIInt  size;
129*4e2e9504SJed Brown   ZCode       *donor_verts, *donor_minz;
130*4e2e9504SJed Brown   PetscSFNode *leaf;
131*4e2e9504SJed Brown 
132*4e2e9504SJed Brown   PetscFunctionBegin;
133*4e2e9504SJed Brown   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
134*4e2e9504SJed Brown   PetscCallMPI(MPI_Comm_size(comm, &size));
135*4e2e9504SJed Brown   PetscCall(DMGetDimension(dm, &dim));
136*4e2e9504SJed Brown   const PetscInt csize = PetscPowInt(2, dim - 1);
137*4e2e9504SJed Brown   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
138*4e2e9504SJed Brown   PetscCall(PetscSegBufferGetSize(per_faces, &num_faces));
139*4e2e9504SJed Brown   PetscCall(PetscSegBufferExtractInPlace(per_faces, &faces));
140*4e2e9504SJed Brown   PetscCall(PetscSegBufferExtractInPlace(donor_face_closure, &donor_verts));
141*4e2e9504SJed Brown   PetscCall(PetscMalloc1(num_faces, &donor_minz));
142*4e2e9504SJed Brown   PetscCall(PetscMalloc1(num_faces, &leaf));
143*4e2e9504SJed Brown   for (PetscInt i = 0; i < (PetscInt)num_faces; i++) {
144*4e2e9504SJed Brown     ZCode minz = donor_verts[i * csize];
145*4e2e9504SJed Brown     for (PetscInt j = 1; j < csize; j++) minz = PetscMin(minz, donor_verts[i * csize + j]);
146*4e2e9504SJed Brown     donor_minz[i] = minz;
147*4e2e9504SJed Brown   }
148*4e2e9504SJed Brown   {
149*4e2e9504SJed Brown     PetscBool sorted;
150*4e2e9504SJed Brown     PetscCall(PetscSortedInt64(num_faces, (const PetscInt64 *)donor_minz, &sorted));
151*4e2e9504SJed Brown     PetscCheck(sorted, comm, PETSC_ERR_PLIB, "minz not sorted; periodicity in multiple dimensions not yet supported");
152*4e2e9504SJed Brown   }
153*4e2e9504SJed Brown   for (PetscInt i = 0; i < (PetscInt)num_faces;) {
154*4e2e9504SJed Brown     ZCode    z           = donor_minz[i];
155*4e2e9504SJed Brown     PetscInt remote_rank = ZCodeFind(z, size + 1, layout->zstarts), remote_count = 0;
156*4e2e9504SJed Brown     if (remote_rank < 0) remote_rank = -(remote_rank + 1) - 1;
157*4e2e9504SJed Brown     // Process all the vertices on this rank
158*4e2e9504SJed Brown     for (ZCode rz = layout->zstarts[remote_rank]; rz < layout->zstarts[remote_rank + 1]; rz++) {
159*4e2e9504SJed Brown       Ijk loc = ZCodeSplit(rz);
160*4e2e9504SJed Brown       if (rz == z) {
161*4e2e9504SJed Brown         leaf[i].rank  = remote_rank;
162*4e2e9504SJed Brown         leaf[i].index = remote_count;
163*4e2e9504SJed Brown         i++;
164*4e2e9504SJed Brown         if (i == (PetscInt)num_faces) break;
165*4e2e9504SJed Brown         z = donor_minz[i];
166*4e2e9504SJed Brown       }
167*4e2e9504SJed Brown       if (IjkActive(layout->vextent, loc)) remote_count++;
168*4e2e9504SJed Brown     }
169*4e2e9504SJed Brown   }
170*4e2e9504SJed Brown   PetscCall(PetscFree(donor_minz));
171*4e2e9504SJed Brown   PetscSF sfper;
172*4e2e9504SJed Brown   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), &sfper));
173*4e2e9504SJed Brown   PetscCall(PetscSFSetGraph(sfper, vEnd - vStart, num_faces, PETSC_NULL, PETSC_USE_POINTER, leaf, PETSC_USE_POINTER));
174*4e2e9504SJed Brown   const PetscInt *my_donor_degree;
175*4e2e9504SJed Brown   PetscCall(PetscSFComputeDegreeBegin(sfper, &my_donor_degree));
176*4e2e9504SJed Brown   PetscCall(PetscSFComputeDegreeEnd(sfper, &my_donor_degree));
177*4e2e9504SJed Brown   PetscInt num_multiroots = 0;
178*4e2e9504SJed Brown   for (PetscInt i = 0; i < vEnd - vStart; i++) {
179*4e2e9504SJed Brown     num_multiroots += my_donor_degree[i];
180*4e2e9504SJed Brown     if (my_donor_degree[i] == 0) continue;
181*4e2e9504SJed Brown     PetscAssert(my_donor_degree[i] == 1, comm, PETSC_ERR_SUP, "Local vertex has multiple faces");
182*4e2e9504SJed Brown   }
183*4e2e9504SJed Brown   PetscInt *my_donors, *donor_indices, *my_donor_indices;
184*4e2e9504SJed Brown   size_t    num_my_donors;
185*4e2e9504SJed Brown   PetscCall(PetscSegBufferGetSize(my_donor_faces, &num_my_donors));
186*4e2e9504SJed Brown   PetscCheck((PetscInt)num_my_donors == num_multiroots, PETSC_COMM_SELF, PETSC_ERR_SUP, "Donor request does not match expected donors");
187*4e2e9504SJed Brown   PetscCall(PetscSegBufferExtractInPlace(my_donor_faces, &my_donors));
188*4e2e9504SJed Brown   PetscCall(PetscMalloc1(vEnd - vStart, &my_donor_indices));
189*4e2e9504SJed Brown   for (PetscInt i = 0; i < (PetscInt)num_my_donors; i++) {
190*4e2e9504SJed Brown     PetscInt f = my_donors[i];
191*4e2e9504SJed Brown     PetscInt num_points, *points = NULL, minv = PETSC_MAX_INT;
192*4e2e9504SJed Brown     PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &num_points, &points));
193*4e2e9504SJed Brown     for (PetscInt j = 0; j < num_points; j++) {
194*4e2e9504SJed Brown       PetscInt p = points[2 * j];
195*4e2e9504SJed Brown       if (p < vStart || vEnd <= p) continue;
196*4e2e9504SJed Brown       minv = PetscMin(minv, p);
197*4e2e9504SJed Brown     }
198*4e2e9504SJed Brown     PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &num_points, &points));
199*4e2e9504SJed Brown     PetscAssert(my_donor_degree[minv - vStart] == 1, comm, PETSC_ERR_SUP, "Local vertex not requested");
200*4e2e9504SJed Brown     my_donor_indices[minv - vStart] = f;
201*4e2e9504SJed Brown   }
202*4e2e9504SJed Brown   PetscCall(PetscMalloc1(num_faces, &donor_indices));
203*4e2e9504SJed Brown   PetscCall(PetscSFBcastBegin(sfper, MPIU_INT, my_donor_indices, donor_indices, MPI_REPLACE));
204*4e2e9504SJed Brown   PetscCall(PetscSFBcastEnd(sfper, MPIU_INT, my_donor_indices, donor_indices, MPI_REPLACE));
205*4e2e9504SJed Brown   PetscCall(PetscFree(my_donor_indices));
206*4e2e9504SJed Brown   // Modify our leafs so they point to donor faces instead of donor minz. Additionally, give them indices as faces.
207*4e2e9504SJed Brown   for (PetscInt i = 0; i < (PetscInt)num_faces; i++) leaf[i].index = donor_indices[i];
208*4e2e9504SJed Brown   PetscCall(PetscFree(donor_indices));
209*4e2e9504SJed Brown   PetscInt pStart, pEnd;
210*4e2e9504SJed Brown   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
211*4e2e9504SJed Brown   PetscCall(PetscSFSetGraph(sfper, pEnd - pStart, num_faces, faces, PETSC_COPY_VALUES, leaf, PETSC_OWN_POINTER));
212*4e2e9504SJed Brown   PetscCall(PetscObjectSetName((PetscObject)sfper, "Periodic Faces"));
213*4e2e9504SJed Brown   PetscCall(PetscSFViewFromOptions(sfper, NULL, "-sfper_view"));
214*4e2e9504SJed Brown 
215*4e2e9504SJed Brown   PetscCall(PetscSegBufferDestroy(&per_faces));
216*4e2e9504SJed Brown   PetscCall(PetscSegBufferDestroy(&donor_face_closure));
217*4e2e9504SJed Brown   PetscCall(PetscSegBufferDestroy(&my_donor_faces));
218*4e2e9504SJed Brown   PetscCall(DMPlexSetPeriodicFaceSF(dm, sfper));
219*4e2e9504SJed Brown   PetscCall(PetscSFDestroy(&sfper));
220*4e2e9504SJed Brown   PetscFunctionReturn(0);
221*4e2e9504SJed Brown }
222*4e2e9504SJed Brown 
2233e72e933SJed Brown PetscErrorCode DMPlexCreateBoxMesh_Tensor_SFC_Internal(DM dm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate)
2243e72e933SJed Brown {
2253e72e933SJed Brown   PetscInt  eextent[3] = {1, 1, 1}, vextent[3] = {1, 1, 1};
2263e72e933SJed Brown   const Ijk closure_1[] = {
2273e72e933SJed Brown     {0, 0, 0},
2283e72e933SJed Brown     {1, 0, 0},
2293e72e933SJed Brown   };
2303e72e933SJed Brown   const Ijk closure_2[] = {
2313e72e933SJed Brown     {0, 0, 0},
2323e72e933SJed Brown     {1, 0, 0},
2333e72e933SJed Brown     {1, 1, 0},
2343e72e933SJed Brown     {0, 1, 0},
2353e72e933SJed Brown   };
2363e72e933SJed Brown   const Ijk closure_3[] = {
2373e72e933SJed Brown     {0, 0, 0},
2383e72e933SJed Brown     {0, 1, 0},
2393e72e933SJed Brown     {1, 1, 0},
2403e72e933SJed Brown     {1, 0, 0},
2413e72e933SJed Brown     {0, 0, 1},
2423e72e933SJed Brown     {1, 0, 1},
2433e72e933SJed Brown     {1, 1, 1},
2443e72e933SJed Brown     {0, 1, 1},
2453e72e933SJed Brown   };
2463431e603SJed Brown   const Ijk *const closure_dim[] = {NULL, closure_1, closure_2, closure_3};
2473431e603SJed Brown   // This must be kept consistent with DMPlexCreateCubeMesh_Internal
2483431e603SJed Brown   const PetscInt        face_marker_1[]   = {1, 2};
2493431e603SJed Brown   const PetscInt        face_marker_2[]   = {4, 2, 1, 3};
2503431e603SJed Brown   const PetscInt        face_marker_3[]   = {6, 5, 3, 4, 1, 2};
2513431e603SJed Brown   const PetscInt *const face_marker_dim[] = {NULL, face_marker_1, face_marker_2, face_marker_3};
2523e72e933SJed Brown 
2533e72e933SJed Brown   PetscFunctionBegin;
2543e72e933SJed Brown   PetscValidPointer(dm, 1);
2553e72e933SJed Brown   PetscValidLogicalCollectiveInt(dm, dim, 2);
2563e72e933SJed Brown   PetscCall(DMSetDimension(dm, dim));
2573e72e933SJed Brown   PetscMPIInt rank, size;
2583e72e933SJed Brown   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
2593e72e933SJed Brown   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
2603e72e933SJed Brown   for (PetscInt i = 0; i < dim; i++) {
2613e72e933SJed Brown     eextent[i] = faces[i];
2623e72e933SJed Brown     vextent[i] = faces[i] + 1;
2633e72e933SJed Brown   }
264c77877e3SJed Brown   ZLayout layout;
265c77877e3SJed Brown   PetscCall(ZLayoutCreate(size, eextent, vextent, &layout));
2663e72e933SJed Brown   PetscZSet vset; // set of all vertices in the closure of the owned elements
2673e72e933SJed Brown   PetscCall(PetscZSetCreate(&vset));
2683e72e933SJed Brown   PetscInt local_elems = 0;
2693e72e933SJed Brown   for (ZCode z = layout.zstarts[rank]; z < layout.zstarts[rank + 1]; z++) {
2703e72e933SJed Brown     Ijk loc = ZCodeSplit(z);
2713e72e933SJed Brown     if (IjkActive(layout.vextent, loc)) PetscZSetAdd(vset, z);
2723e72e933SJed Brown     if (IjkActive(layout.eextent, loc)) {
2733e72e933SJed Brown       local_elems++;
2743e72e933SJed Brown       // Add all neighboring vertices to set
2753e72e933SJed Brown       for (PetscInt n = 0; n < PetscPowInt(2, dim); n++) {
2763e72e933SJed Brown         Ijk   inc  = closure_dim[dim][n];
2773e72e933SJed Brown         Ijk   nloc = {loc.i + inc.i, loc.j + inc.j, loc.k + inc.k};
2783e72e933SJed Brown         ZCode v    = ZEncode(nloc);
2793e72e933SJed Brown         PetscZSetAdd(vset, v);
2803e72e933SJed Brown       }
2813e72e933SJed Brown     }
2823e72e933SJed Brown   }
2833e72e933SJed Brown   PetscInt local_verts, off = 0;
2843e72e933SJed Brown   ZCode   *vert_z;
2853e72e933SJed Brown   PetscCall(PetscZSetGetSize(vset, &local_verts));
2863e72e933SJed Brown   PetscCall(PetscMalloc1(local_verts, &vert_z));
2873e72e933SJed Brown   PetscCall(PetscZSetGetElems(vset, &off, vert_z));
2883e72e933SJed Brown   PetscCall(PetscZSetDestroy(&vset));
2893e72e933SJed Brown   // ZCode is unsigned for bitwise convenience, but highest bit should never be set, so can interpret as signed
2903e72e933SJed Brown   PetscCall(PetscSortInt64(local_verts, (PetscInt64 *)vert_z));
2913e72e933SJed Brown 
2923e72e933SJed Brown   PetscCall(DMPlexSetChart(dm, 0, local_elems + local_verts));
2933e72e933SJed Brown   for (PetscInt e = 0; e < local_elems; e++) PetscCall(DMPlexSetConeSize(dm, e, PetscPowInt(2, dim)));
2943e72e933SJed Brown   PetscCall(DMSetUp(dm));
2953e72e933SJed Brown   {
2963e72e933SJed Brown     PetscInt e = 0;
2973e72e933SJed Brown     for (ZCode z = layout.zstarts[rank]; z < layout.zstarts[rank + 1]; z++) {
2983e72e933SJed Brown       Ijk loc = ZCodeSplit(z);
2993e72e933SJed Brown       if (!IjkActive(layout.eextent, loc)) continue;
3003e72e933SJed Brown       PetscInt cone[8], orient[8] = {0};
3013e72e933SJed Brown       for (PetscInt n = 0; n < PetscPowInt(2, dim); n++) {
3023e72e933SJed Brown         Ijk      inc  = closure_dim[dim][n];
3033e72e933SJed Brown         Ijk      nloc = {loc.i + inc.i, loc.j + inc.j, loc.k + inc.k};
3043e72e933SJed Brown         ZCode    v    = ZEncode(nloc);
3053e72e933SJed Brown         PetscInt ci   = ZCodeFind(v, local_verts, vert_z);
3063e72e933SJed Brown         PetscAssert(ci >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Did not find neighbor vertex in set");
3073e72e933SJed Brown         cone[n] = local_elems + ci;
3083e72e933SJed Brown       }
3093e72e933SJed Brown       PetscCall(DMPlexSetCone(dm, e, cone));
3103e72e933SJed Brown       PetscCall(DMPlexSetConeOrientation(dm, e, orient));
3113e72e933SJed Brown       e++;
3123e72e933SJed Brown     }
3133e72e933SJed Brown   }
3143e72e933SJed Brown 
3153e72e933SJed Brown   if (0) {
3163e72e933SJed Brown     DMLabel depth;
3173e72e933SJed Brown     PetscCall(DMCreateLabel(dm, "depth"));
3183e72e933SJed Brown     PetscCall(DMPlexGetDepthLabel(dm, &depth));
3193e72e933SJed Brown     PetscCall(DMLabelSetStratumBounds(depth, 0, local_elems, local_elems + local_verts));
3203e72e933SJed Brown     PetscCall(DMLabelSetStratumBounds(depth, 1, 0, local_elems));
3213e72e933SJed Brown   } else {
3223e72e933SJed Brown     PetscCall(DMPlexSymmetrize(dm));
3233e72e933SJed Brown     PetscCall(DMPlexStratify(dm));
3243e72e933SJed Brown   }
3253e72e933SJed Brown   { // Create point SF
3263e72e933SJed Brown     PetscSF sf;
3273e72e933SJed Brown     PetscSFCreate(PetscObjectComm((PetscObject)dm), &sf);
3283e72e933SJed Brown     PetscInt owned_verts = ZCodeFind(layout.zstarts[rank + 1], local_verts, vert_z);
3293e72e933SJed Brown     if (owned_verts < 0) owned_verts = -(owned_verts + 1); // We don't care whether the key was found
3303e72e933SJed Brown     PetscInt     num_ghosts = local_verts - owned_verts;   // Due to sorting, owned vertices always come first
3313e72e933SJed Brown     PetscInt    *local_ghosts;
3323e72e933SJed Brown     PetscSFNode *ghosts;
3333e72e933SJed Brown     PetscCall(PetscMalloc1(num_ghosts, &local_ghosts));
3343e72e933SJed Brown     PetscCall(PetscMalloc1(num_ghosts, &ghosts));
3353e72e933SJed Brown     for (PetscInt i = 0; i < num_ghosts;) {
3363e72e933SJed Brown       ZCode    z           = vert_z[owned_verts + i];
3373e72e933SJed Brown       PetscInt remote_rank = ZCodeFind(z, size + 1, layout.zstarts), remote_count = 0;
3383e72e933SJed Brown       if (remote_rank < 0) remote_rank = -(remote_rank + 1) - 1;
3393e72e933SJed Brown       // We have a new remote rank; find all the ghost indices (which are contiguous in vert_z)
3403e72e933SJed Brown 
3413e72e933SJed Brown       // Count the elements on remote_rank
342*4e2e9504SJed Brown       PetscInt remote_elem = ZLayoutElementsOnRank(&layout, remote_rank);
3433e72e933SJed Brown 
3443e72e933SJed Brown       // Traverse vertices and make ghost links
3453e72e933SJed Brown       for (ZCode rz = layout.zstarts[remote_rank]; rz < layout.zstarts[remote_rank + 1]; rz++) {
3463e72e933SJed Brown         Ijk loc = ZCodeSplit(rz);
3473e72e933SJed Brown         if (rz == z) {
3483e72e933SJed Brown           local_ghosts[i] = local_elems + owned_verts + i;
3493e72e933SJed Brown           ghosts[i].rank  = remote_rank;
3503e72e933SJed Brown           ghosts[i].index = remote_elem + remote_count;
3513e72e933SJed Brown           i++;
3523e72e933SJed Brown           if (i == num_ghosts) break;
3533e72e933SJed Brown           z = vert_z[owned_verts + i];
3543e72e933SJed Brown         }
3553e72e933SJed Brown         if (IjkActive(layout.vextent, loc)) remote_count++;
3563e72e933SJed Brown       }
3573e72e933SJed Brown     }
3583e72e933SJed Brown     PetscCall(PetscSFSetGraph(sf, local_elems + local_verts, num_ghosts, local_ghosts, PETSC_OWN_POINTER, ghosts, PETSC_OWN_POINTER));
3593e72e933SJed Brown     PetscCall(PetscObjectSetName((PetscObject)sf, "SFC Point SF"));
3603e72e933SJed Brown     PetscCall(DMSetPointSF(dm, sf));
3613e72e933SJed Brown     PetscCall(PetscSFDestroy(&sf));
3623e72e933SJed Brown   }
3633e72e933SJed Brown   {
3643e72e933SJed Brown     Vec          coordinates;
3653e72e933SJed Brown     PetscScalar *coords;
3663e72e933SJed Brown     PetscSection coord_section;
3673e72e933SJed Brown     PetscInt     coord_size;
3683e72e933SJed Brown     PetscCall(DMGetCoordinateSection(dm, &coord_section));
3693e72e933SJed Brown     PetscCall(PetscSectionSetNumFields(coord_section, 1));
3703e72e933SJed Brown     PetscCall(PetscSectionSetFieldComponents(coord_section, 0, dim));
3713e72e933SJed Brown     PetscCall(PetscSectionSetChart(coord_section, local_elems, local_elems + local_verts));
3723e72e933SJed Brown     for (PetscInt v = 0; v < local_verts; v++) {
3733e72e933SJed Brown       PetscInt point = local_elems + v;
3743e72e933SJed Brown       PetscCall(PetscSectionSetDof(coord_section, point, dim));
3753e72e933SJed Brown       PetscCall(PetscSectionSetFieldDof(coord_section, point, 0, dim));
3763e72e933SJed Brown     }
3773e72e933SJed Brown     PetscCall(PetscSectionSetUp(coord_section));
3783e72e933SJed Brown     PetscCall(PetscSectionGetStorageSize(coord_section, &coord_size));
3793e72e933SJed Brown     PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates));
3803e72e933SJed Brown     PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
3813e72e933SJed Brown     PetscCall(VecSetSizes(coordinates, coord_size, PETSC_DETERMINE));
3823e72e933SJed Brown     PetscCall(VecSetBlockSize(coordinates, dim));
3833e72e933SJed Brown     PetscCall(VecSetType(coordinates, VECSTANDARD));
3843e72e933SJed Brown     PetscCall(VecGetArray(coordinates, &coords));
3853e72e933SJed Brown     for (PetscInt v = 0; v < local_verts; v++) {
3863e72e933SJed Brown       Ijk loc             = ZCodeSplit(vert_z[v]);
3873e72e933SJed Brown       coords[v * dim + 0] = lower[0] + loc.i * (upper[0] - lower[0]) / layout.eextent.i;
3883e72e933SJed Brown       if (dim > 1) coords[v * dim + 1] = lower[1] + loc.j * (upper[1] - lower[1]) / layout.eextent.j;
3893e72e933SJed Brown       if (dim > 2) coords[v * dim + 2] = lower[2] + loc.k * (upper[2] - lower[2]) / layout.eextent.k;
3903e72e933SJed Brown     }
3913e72e933SJed Brown     PetscCall(VecRestoreArray(coordinates, &coords));
3923e72e933SJed Brown     PetscCall(DMSetCoordinatesLocal(dm, coordinates));
3933e72e933SJed Brown     PetscCall(VecDestroy(&coordinates));
3943e72e933SJed Brown   }
3953e72e933SJed Brown   if (interpolate) {
3963431e603SJed Brown     PetscCall(DMPlexInterpolateInPlace_Internal(dm));
3973431e603SJed Brown 
3983431e603SJed Brown     DMLabel label;
3993431e603SJed Brown     PetscCall(DMCreateLabel(dm, "Face Sets"));
4003431e603SJed Brown     PetscCall(DMGetLabel(dm, "Face Sets", &label));
401*4e2e9504SJed Brown     PetscSegBuffer per_faces, donor_face_closure, my_donor_faces;
402*4e2e9504SJed Brown     PetscCall(PetscSegBufferCreate(sizeof(PetscInt), 64, &per_faces));
403*4e2e9504SJed Brown     PetscCall(PetscSegBufferCreate(sizeof(PetscInt), 64, &my_donor_faces));
404*4e2e9504SJed Brown     PetscCall(PetscSegBufferCreate(sizeof(ZCode), 64 * PetscPowInt(2, dim), &donor_face_closure));
4053431e603SJed Brown     PetscInt fStart, fEnd, vStart, vEnd;
4063431e603SJed Brown     PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
4073431e603SJed Brown     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
4083431e603SJed Brown     for (PetscInt f = fStart; f < fEnd; f++) {
409*4e2e9504SJed Brown       PetscInt npoints, *points = NULL, num_fverts = 0, fverts[8];
4103431e603SJed Brown       PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &npoints, &points));
4113431e603SJed Brown       PetscInt bc_count[6] = {0};
4123431e603SJed Brown       for (PetscInt i = 0; i < npoints; i++) {
4133431e603SJed Brown         PetscInt p = points[2 * i];
4143431e603SJed Brown         if (p < vStart || vEnd <= p) continue;
415*4e2e9504SJed Brown         fverts[num_fverts++] = p;
4163431e603SJed Brown         Ijk loc              = ZCodeSplit(vert_z[p - vStart]);
4173431e603SJed Brown         // Convention here matches DMPlexCreateCubeMesh_Internal
4183431e603SJed Brown         bc_count[0] += loc.i == 0;
4193431e603SJed Brown         bc_count[1] += loc.i == layout.vextent.i - 1;
4203431e603SJed Brown         bc_count[2] += loc.j == 0;
4213431e603SJed Brown         bc_count[3] += loc.j == layout.vextent.j - 1;
4223431e603SJed Brown         bc_count[4] += loc.k == 0;
4233431e603SJed Brown         bc_count[5] += loc.k == layout.vextent.k - 1;
4243e72e933SJed Brown       }
4253431e603SJed Brown       PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &npoints, &points));
4263431e603SJed Brown       for (PetscInt bc = 0, bc_match = 0; bc < 2 * dim; bc++) {
4273431e603SJed Brown         if (bc_count[bc] == PetscPowInt(2, dim - 1)) {
428*4e2e9504SJed Brown           if (periodicity[bc / 2] == DM_BOUNDARY_PERIODIC) {
429*4e2e9504SJed Brown             PetscInt *put;
430*4e2e9504SJed Brown             if (bc % 2 == 0) { // donor face; no label
431*4e2e9504SJed Brown               PetscCall(PetscSegBufferGet(my_donor_faces, 1, &put));
432*4e2e9504SJed Brown               *put = f;
433*4e2e9504SJed Brown             } else { // periodic face
434*4e2e9504SJed Brown               PetscCall(PetscSegBufferGet(per_faces, 1, &put));
435*4e2e9504SJed Brown               *put = f;
436*4e2e9504SJed Brown               ZCode *zput;
437*4e2e9504SJed Brown               PetscCall(PetscSegBufferGet(donor_face_closure, num_fverts, &zput));
438*4e2e9504SJed Brown               for (PetscInt i = 0; i < num_fverts; i++) {
439*4e2e9504SJed Brown                 Ijk loc = ZCodeSplit(vert_z[fverts[i] - vStart]);
440*4e2e9504SJed Brown                 switch (bc / 2) {
441*4e2e9504SJed Brown                 case 0:
442*4e2e9504SJed Brown                   loc.i = 0;
443*4e2e9504SJed Brown                   break;
444*4e2e9504SJed Brown                 case 1:
445*4e2e9504SJed Brown                   loc.j = 0;
446*4e2e9504SJed Brown                   break;
447*4e2e9504SJed Brown                 case 2:
448*4e2e9504SJed Brown                   loc.k = 0;
449*4e2e9504SJed Brown                   break;
450*4e2e9504SJed Brown                 }
451*4e2e9504SJed Brown                 *zput++ = ZEncode(loc);
452*4e2e9504SJed Brown               }
453*4e2e9504SJed Brown             }
454*4e2e9504SJed Brown             continue;
455*4e2e9504SJed Brown           }
4563431e603SJed Brown           PetscAssert(bc_match == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face matches multiple face sets");
4573431e603SJed Brown           PetscCall(DMLabelSetValue(label, f, face_marker_dim[dim][bc]));
4583431e603SJed Brown           bc_match++;
4593431e603SJed Brown         }
4603431e603SJed Brown       }
4613431e603SJed Brown     }
462*4e2e9504SJed Brown     PetscCall(DMPlexCreateBoxMesh_Tensor_SFC_Periodicity_Private(dm, &layout, vert_z, per_faces, donor_face_closure, my_donor_faces));
4633431e603SJed Brown     // Ensure that the Coordinate DM has our new boundary labels
4643431e603SJed Brown     DM cdm;
4653431e603SJed Brown     PetscCall(DMGetCoordinateDM(dm, &cdm));
4663431e603SJed Brown     PetscCall(DMCopyLabels(dm, cdm, PETSC_COPY_VALUES, PETSC_FALSE, DM_COPY_LABELS_FAIL));
467*4e2e9504SJed Brown     PetscSF sfper;
468*4e2e9504SJed Brown     PetscCall(DMPlexGetPeriodicFaceSF(dm, &sfper));
469*4e2e9504SJed Brown     PetscCall(DMPlexSetPeriodicFaceSF(cdm, sfper));
4703431e603SJed Brown   }
471*4e2e9504SJed Brown   PetscCall(PetscFree(layout.zstarts));
4723431e603SJed Brown   PetscCall(PetscFree(vert_z));
4733e72e933SJed Brown   PetscFunctionReturn(0);
4743e72e933SJed Brown }
475*4e2e9504SJed Brown 
476*4e2e9504SJed Brown /*@
477*4e2e9504SJed Brown   DMPlexSetPeriodicFaceSF - Express periodicity from an existing mesh
478*4e2e9504SJed Brown 
479*4e2e9504SJed Brown   Logically collective
480*4e2e9504SJed Brown 
481*4e2e9504SJed Brown   Input Parameters:
482*4e2e9504SJed Brown + dm - The `DMPLEX` on which to set periodicity
483*4e2e9504SJed Brown - face_sf - SF in which roots are (owned) donor faces and leaves are faces that must be matched to a (possibly remote) donor face.
484*4e2e9504SJed Brown 
485*4e2e9504SJed Brown   Level: advanced
486*4e2e9504SJed Brown 
487*4e2e9504SJed Brown   Notes:
488*4e2e9504SJed Brown 
489*4e2e9504SJed Brown   One can use `-dm_plex_box_sfc` to use this mode of periodicity, wherein the periodic points are distinct both globally
490*4e2e9504SJed Brown   and locally, but are paired when creating a global dof space.
491*4e2e9504SJed Brown 
492*4e2e9504SJed Brown .seealso: [](chapter_unstructured), `DMPLEX`, `DMGetGlobalSection()`, `DMPlexGetPeriodicFaceSF()`
493*4e2e9504SJed Brown @*/
494*4e2e9504SJed Brown PetscErrorCode DMPlexSetPeriodicFaceSF(DM dm, PetscSF face_sf)
495*4e2e9504SJed Brown {
496*4e2e9504SJed Brown   DM_Plex *plex = (DM_Plex *)dm->data;
497*4e2e9504SJed Brown   PetscFunctionBegin;
498*4e2e9504SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
499*4e2e9504SJed Brown   PetscCall(PetscObjectReference((PetscObject)face_sf));
500*4e2e9504SJed Brown   PetscCall(PetscSFDestroy(&plex->periodic.face_sf));
501*4e2e9504SJed Brown   plex->periodic.face_sf = face_sf;
502*4e2e9504SJed Brown   PetscFunctionReturn(0);
503*4e2e9504SJed Brown }
504*4e2e9504SJed Brown 
505*4e2e9504SJed Brown /*@
506*4e2e9504SJed Brown   DMPlexGetPeriodicFaceSF - Obtain periodicity for a mesh
507*4e2e9504SJed Brown 
508*4e2e9504SJed Brown   Logically collective
509*4e2e9504SJed Brown 
510*4e2e9504SJed Brown   Input Parameters:
511*4e2e9504SJed Brown . dm - The `DMPLEX` for which to obtain periodic relation
512*4e2e9504SJed Brown 
513*4e2e9504SJed Brown   Output Parameters:
514*4e2e9504SJed Brown . face_sf - SF in which roots are (owned) donor faces and leaves are faces that must be matched to a (possibly remote) donor face.
515*4e2e9504SJed Brown 
516*4e2e9504SJed Brown   Level: advanced
517*4e2e9504SJed Brown 
518*4e2e9504SJed Brown .seealso: [](chapter_unstructured), `DMPLEX`, `DMGetGlobalSection()`, `DMPlexSetPeriodicFaceSF()`
519*4e2e9504SJed Brown @*/
520*4e2e9504SJed Brown PetscErrorCode DMPlexGetPeriodicFaceSF(DM dm, PetscSF *face_sf)
521*4e2e9504SJed Brown {
522*4e2e9504SJed Brown   DM_Plex *plex = (DM_Plex *)dm->data;
523*4e2e9504SJed Brown   PetscFunctionBegin;
524*4e2e9504SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
525*4e2e9504SJed Brown   *face_sf = plex->periodic.face_sf;
526*4e2e9504SJed Brown   PetscFunctionReturn(0);
527*4e2e9504SJed Brown }
528