xref: /petsc/src/dm/impls/plex/plexsfc.c (revision ea7b3863513b7120da5bbc46fbc9df9d66f981c3) !
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 
20563af49cSJames Wright // Magic numbers taken from https://stackoverflow.com/a/18528775/7564988
213e72e933SJed Brown static unsigned ZCodeSplit1(ZCode z)
223e72e933SJed Brown {
23563af49cSJames Wright   z &= 0x1249249249249249;
24563af49cSJames Wright   z = (z | z >> 2) & 0x10c30c30c30c30c3;
25563af49cSJames Wright   z = (z | z >> 4) & 0x100f00f00f00f00f;
26563af49cSJames Wright   z = (z | z >> 8) & 0x1f0000ff0000ff;
27563af49cSJames Wright   z = (z | z >> 16) & 0x1f00000000ffff;
28563af49cSJames Wright   z = (z | z >> 32) & 0x1fffff;
293e72e933SJed Brown   return (unsigned)z;
303e72e933SJed Brown }
313e72e933SJed Brown 
323e72e933SJed Brown static ZCode ZEncode1(unsigned t)
333e72e933SJed Brown {
343e72e933SJed Brown   ZCode z = t;
35563af49cSJames Wright   z &= 0x1fffff;
36563af49cSJames Wright   z = (z | z << 32) & 0x1f00000000ffff;
37563af49cSJames Wright   z = (z | z << 16) & 0x1f0000ff0000ff;
38563af49cSJames Wright   z = (z | z << 8) & 0x100f00f00f00f00f;
39563af49cSJames Wright   z = (z | z << 4) & 0x10c30c30c30c30c3;
40563af49cSJames Wright   z = (z | z << 2) & 0x1249249249249249;
413e72e933SJed Brown   return z;
423e72e933SJed Brown }
433e72e933SJed Brown 
443e72e933SJed Brown static Ijk ZCodeSplit(ZCode z)
453e72e933SJed Brown {
463e72e933SJed Brown   Ijk c;
473e72e933SJed Brown   c.i = ZCodeSplit1(z >> 2);
483e72e933SJed Brown   c.j = ZCodeSplit1(z >> 1);
493e72e933SJed Brown   c.k = ZCodeSplit1(z >> 0);
503e72e933SJed Brown   return c;
513e72e933SJed Brown }
523e72e933SJed Brown 
533e72e933SJed Brown static ZCode ZEncode(Ijk c)
543e72e933SJed Brown {
556497c311SBarry Smith   ZCode z = (ZEncode1((unsigned int)c.i) << 2) | (ZEncode1((unsigned int)c.j) << 1) | ZEncode1((unsigned int)c.k);
563e72e933SJed Brown   return z;
573e72e933SJed Brown }
583e72e933SJed Brown 
593e72e933SJed Brown static PetscBool IjkActive(Ijk extent, Ijk l)
603e72e933SJed Brown {
613e72e933SJed Brown   if (l.i < extent.i && l.j < extent.j && l.k < extent.k) return PETSC_TRUE;
623e72e933SJed Brown   return PETSC_FALSE;
633e72e933SJed Brown }
643e72e933SJed Brown 
656547ddbdSJed Brown // If z is not the base of an octet (last 3 bits 0), return 0.
666547ddbdSJed Brown //
676547ddbdSJed Brown // If z is the base of an octet, we recursively grow to the biggest structured octet. This is typically useful when a z
686547ddbdSJed Brown // is outside the domain and we wish to skip a (possibly recursively large) octet to find our next interesting point.
696547ddbdSJed Brown static ZCode ZStepOct(ZCode z)
706547ddbdSJed Brown {
716547ddbdSJed Brown   if (PetscUnlikely(z == 0)) return 0; // Infinite loop below if z == 0
726547ddbdSJed Brown   ZCode step = 07;
736547ddbdSJed Brown   for (; (z & step) == 0; step = (step << 3) | 07) { }
746547ddbdSJed Brown   return step >> 3;
756547ddbdSJed Brown }
766547ddbdSJed Brown 
773e72e933SJed Brown // Since element/vertex box extents are typically not equal powers of 2, Z codes that lie within the domain are not contiguous.
785f06a3ddSJed Brown static PetscErrorCode ZLayoutCreate(PetscMPIInt size, const PetscInt eextent[3], const PetscInt vextent[3], ZLayout *layout)
793e72e933SJed Brown {
80c77877e3SJed Brown   PetscFunctionBegin;
815f06a3ddSJed Brown   layout->eextent.i = eextent[0];
825f06a3ddSJed Brown   layout->eextent.j = eextent[1];
835f06a3ddSJed Brown   layout->eextent.k = eextent[2];
845f06a3ddSJed Brown   layout->vextent.i = vextent[0];
855f06a3ddSJed Brown   layout->vextent.j = vextent[1];
865f06a3ddSJed Brown   layout->vextent.k = vextent[2];
875f06a3ddSJed Brown   layout->comm_size = size;
885f06a3ddSJed Brown   layout->zstarts   = NULL;
895f06a3ddSJed Brown   PetscCall(PetscMalloc1(size + 1, &layout->zstarts));
903e72e933SJed Brown 
913e72e933SJed Brown   PetscInt total_elems = eextent[0] * eextent[1] * eextent[2];
923e72e933SJed Brown   ZCode    z           = 0;
935f06a3ddSJed Brown   layout->zstarts[0]   = 0;
946547ddbdSJed Brown   // This loop traverses all vertices in the global domain, so is worth making fast. We use ZStepBound
953e72e933SJed Brown   for (PetscMPIInt r = 0; r < size; r++) {
963e72e933SJed Brown     PetscInt elems_needed = (total_elems / size) + (total_elems % size > r), count;
973e72e933SJed Brown     for (count = 0; count < elems_needed; z++) {
986547ddbdSJed Brown       ZCode skip = ZStepOct(z); // optimistically attempt a longer step
996547ddbdSJed Brown       for (ZCode s = skip;; s >>= 3) {
1006547ddbdSJed Brown         Ijk trial = ZCodeSplit(z + s);
1016547ddbdSJed Brown         if (IjkActive(layout->eextent, trial)) {
1026547ddbdSJed Brown           while (count + s + 1 > (ZCode)elems_needed) s >>= 3; // Shrink the octet
1036547ddbdSJed Brown           count += s + 1;
1046547ddbdSJed Brown           z += s;
1056547ddbdSJed Brown           break;
1066547ddbdSJed Brown         }
1076547ddbdSJed Brown         if (s == 0) { // the whole skip octet is inactive
1086547ddbdSJed Brown           z += skip;
1096547ddbdSJed Brown           break;
1106547ddbdSJed Brown         }
1116547ddbdSJed Brown       }
1123e72e933SJed Brown     }
1133e72e933SJed Brown     // Pick up any extra vertices in the Z ordering before the next rank's first owned element.
114c77877e3SJed Brown     //
1155f06a3ddSJed Brown     // This leads to poorly balanced vertices when eextent is a power of 2, since all the fringe vertices end up
116c77877e3SJed Brown     // on the last rank. A possible solution is to balance the Z-order vertices independently from the cells, which will
117c77877e3SJed Brown     // result in a lot of element closures being remote. We could finish marking boundary conditions, then do a round of
118c77877e3SJed Brown     // vertex ownership smoothing (which would reorder and redistribute vertices without touching element distribution).
119c77877e3SJed Brown     // Another would be to have an analytic ownership criteria for vertices in the fringe veextent - eextent. This would
120c77877e3SJed Brown     // complicate the job of identifying an owner and its offset.
1215f06a3ddSJed Brown     //
1225f06a3ddSJed Brown     // The current recommended approach is to let `-dm_distribute 1` (default) resolve vertex ownership. This is
1235f06a3ddSJed Brown     // *mandatory* with isoperiodicity (except in special cases) to remove standed vertices from local spaces. Here's
1245f06a3ddSJed Brown     // the issue:
1255f06a3ddSJed Brown     //
1265f06a3ddSJed Brown     // Consider this partition on rank 0 (left) and rank 1.
1275f06a3ddSJed Brown     //
1285f06a3ddSJed Brown     //    4 --------  5 -- 14 --10 -- 21 --11
1295f06a3ddSJed Brown     //                |          |          |
1305f06a3ddSJed Brown     // 7 -- 16 --  8  |          |          |
1315f06a3ddSJed Brown     // |           |  3 -------  7 -------  9
1325f06a3ddSJed Brown     // |           |             |          |
1335f06a3ddSJed Brown     // 4 --------  6 ------ 10   |          |
1345f06a3ddSJed Brown     // |           |         |   6 -- 16 -- 8
1355f06a3ddSJed Brown     // |           |         |
1365f06a3ddSJed Brown     // 3 ---11---  5 --18--  9
1375f06a3ddSJed Brown     //
1385f06a3ddSJed Brown     // The periodic face SF looks like
1395f06a3ddSJed Brown     // [0] Number of roots=21, leaves=1, remote ranks=1
1405f06a3ddSJed Brown     // [0] 16 <- (0,11)
1415f06a3ddSJed Brown     // [1] Number of roots=22, leaves=2, remote ranks=2
1425f06a3ddSJed Brown     // [1] 14 <- (0,18)
1435f06a3ddSJed Brown     // [1] 21 <- (1,16)
1445f06a3ddSJed Brown     //
1455f06a3ddSJed Brown     // In handling face (0,16), rank 0 learns that (0,7) and (0,8) map to (0,3) and (0,5) respectively, thus we won't use
1465f06a3ddSJed Brown     // the point SF links to (1,4) and (1,5). Rank 1 learns about the periodic mapping of (1,5) while handling face
1475f06a3ddSJed Brown     // (1,14), but never learns that vertex (1,4) has been mapped to (0,3) by face (0,16).
1485f06a3ddSJed Brown     //
1495f06a3ddSJed Brown     // We can relatively easily inform vertex (1,4) of this mapping, but it stays in rank 1's local space despite not
1505f06a3ddSJed Brown     // being in the closure and thus not being contributed to. This would be mostly harmless except that some viewer
1515f06a3ddSJed Brown     // routines expect all local points to be somehow significant. It is not easy to analytically remove the (1,4)
1525f06a3ddSJed Brown     // vertex because the point SF and isoperiodic face SF would need to be updated to account for removal of the
1535f06a3ddSJed Brown     // stranded vertices.
1545f06a3ddSJed Brown     for (; z <= ZEncode(layout->vextent); z++) {
1553e72e933SJed Brown       Ijk loc = ZCodeSplit(z);
1565f06a3ddSJed Brown       if (IjkActive(layout->eextent, loc)) break;
1576547ddbdSJed Brown       z += ZStepOct(z);
1583e72e933SJed Brown     }
1595f06a3ddSJed Brown     layout->zstarts[r + 1] = z;
1603e72e933SJed Brown   }
1616547ddbdSJed Brown   layout->zstarts[size] = ZEncode(layout->vextent);
1623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1633e72e933SJed Brown }
1643e72e933SJed Brown 
1654e2e9504SJed Brown static PetscInt ZLayoutElementsOnRank(const ZLayout *layout, PetscMPIInt rank)
1664e2e9504SJed Brown {
1674e2e9504SJed Brown   PetscInt remote_elem = 0;
1684e2e9504SJed Brown   for (ZCode rz = layout->zstarts[rank]; rz < layout->zstarts[rank + 1]; rz++) {
1694e2e9504SJed Brown     Ijk loc = ZCodeSplit(rz);
1704e2e9504SJed Brown     if (IjkActive(layout->eextent, loc)) remote_elem++;
1716547ddbdSJed Brown     else rz += ZStepOct(rz);
1724e2e9504SJed Brown   }
1734e2e9504SJed Brown   return remote_elem;
1744e2e9504SJed Brown }
1754e2e9504SJed Brown 
17666976f2fSJacob Faibussowitsch static PetscInt ZCodeFind(ZCode key, PetscInt n, const ZCode X[])
1773e72e933SJed Brown {
1783e72e933SJed Brown   PetscInt lo = 0, hi = n;
1793e72e933SJed Brown 
1803e72e933SJed Brown   if (n == 0) return -1;
1813e72e933SJed Brown   while (hi - lo > 1) {
1823e72e933SJed Brown     PetscInt mid = lo + (hi - lo) / 2;
1833e72e933SJed Brown     if (key < X[mid]) hi = mid;
1843e72e933SJed Brown     else lo = mid;
1853e72e933SJed Brown   }
1863e72e933SJed Brown   return key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
1873e72e933SJed Brown }
1883e72e933SJed Brown 
1899ae3492bSJames Wright static PetscErrorCode DMPlexCreateBoxMesh_Tensor_SFC_Periodicity_Private(DM dm, const ZLayout *layout, const ZCode *vert_z, PetscSegBuffer per_faces[3], const PetscReal *lower, const PetscReal *upper, const DMBoundaryType *periodicity, PetscSegBuffer donor_face_closure[3], PetscSegBuffer my_donor_faces[3])
1904e2e9504SJed Brown {
1914e2e9504SJed Brown   MPI_Comm    comm;
1929ae3492bSJames Wright   PetscInt    dim, vStart, vEnd;
1934e2e9504SJed Brown   PetscMPIInt size;
1949ae3492bSJames Wright   PetscSF     face_sfs[3];
1959ae3492bSJames Wright   PetscScalar transforms[3][4][4] = {{{0}}};
1964e2e9504SJed Brown 
1974e2e9504SJed Brown   PetscFunctionBegin;
1984e2e9504SJed Brown   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1994e2e9504SJed Brown   PetscCallMPI(MPI_Comm_size(comm, &size));
2004e2e9504SJed Brown   PetscCall(DMGetDimension(dm, &dim));
2014e2e9504SJed Brown   const PetscInt csize = PetscPowInt(2, dim - 1);
2024e2e9504SJed Brown   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
2039ae3492bSJames Wright 
2049ae3492bSJames Wright   PetscInt num_directions = 0;
2059ae3492bSJames Wright   for (PetscInt direction = 0; direction < dim; direction++) {
2066497c311SBarry Smith     PetscCount   num_faces;
2079ae3492bSJames Wright     PetscInt    *faces;
2089ae3492bSJames Wright     ZCode       *donor_verts, *donor_minz;
2099ae3492bSJames Wright     PetscSFNode *leaf;
2106497c311SBarry Smith     PetscCount   num_multiroots = 0;
2116497c311SBarry Smith     PetscInt     pStart, pEnd;
2126497c311SBarry Smith     PetscBool    sorted;
2136497c311SBarry Smith     PetscInt     inum_faces;
2149ae3492bSJames Wright 
2159ae3492bSJames Wright     if (periodicity[direction] != DM_BOUNDARY_PERIODIC) continue;
2169ae3492bSJames Wright     PetscCall(PetscSegBufferGetSize(per_faces[direction], &num_faces));
2179ae3492bSJames Wright     PetscCall(PetscSegBufferExtractInPlace(per_faces[direction], &faces));
2189ae3492bSJames Wright     PetscCall(PetscSegBufferExtractInPlace(donor_face_closure[direction], &donor_verts));
2194e2e9504SJed Brown     PetscCall(PetscMalloc1(num_faces, &donor_minz));
2204e2e9504SJed Brown     PetscCall(PetscMalloc1(num_faces, &leaf));
2216497c311SBarry Smith     for (PetscCount i = 0; i < num_faces; i++) {
2224e2e9504SJed Brown       ZCode minz = donor_verts[i * csize];
2236497c311SBarry Smith 
2244e2e9504SJed Brown       for (PetscInt j = 1; j < csize; j++) minz = PetscMin(minz, donor_verts[i * csize + j]);
2254e2e9504SJed Brown       donor_minz[i] = minz;
2264e2e9504SJed Brown     }
2276497c311SBarry Smith     PetscCall(PetscIntCast(num_faces, &inum_faces));
2286497c311SBarry Smith     PetscCall(PetscSortedInt64(inum_faces, (const PetscInt64 *)donor_minz, &sorted));
2299ae3492bSJames Wright     // If a donor vertex were chosen to broker multiple faces, we would have a logic error.
2309ae3492bSJames Wright     // Checking for sorting is a cheap check that there are no duplicates.
2319ae3492bSJames Wright     PetscCheck(sorted, PETSC_COMM_SELF, PETSC_ERR_PLIB, "minz not sorted; possible duplicates not checked");
2326497c311SBarry Smith     for (PetscCount i = 0; i < num_faces;) {
2334e2e9504SJed Brown       ZCode       z           = donor_minz[i];
2346497c311SBarry Smith       PetscMPIInt remote_rank = (PetscMPIInt)ZCodeFind(z, size + 1, layout->zstarts), remote_count = 0;
2356497c311SBarry Smith 
2364e2e9504SJed Brown       if (remote_rank < 0) remote_rank = -(remote_rank + 1) - 1;
2374e2e9504SJed Brown       // Process all the vertices on this rank
2384e2e9504SJed Brown       for (ZCode rz = layout->zstarts[remote_rank]; rz < layout->zstarts[remote_rank + 1]; rz++) {
2394e2e9504SJed Brown         Ijk loc = ZCodeSplit(rz);
2406497c311SBarry Smith 
2414e2e9504SJed Brown         if (rz == z) {
2424e2e9504SJed Brown           leaf[i].rank  = remote_rank;
2434e2e9504SJed Brown           leaf[i].index = remote_count;
2444e2e9504SJed Brown           i++;
2456497c311SBarry Smith           if (i == num_faces) break;
2464e2e9504SJed Brown           z = donor_minz[i];
2474e2e9504SJed Brown         }
2484e2e9504SJed Brown         if (IjkActive(layout->vextent, loc)) remote_count++;
2494e2e9504SJed Brown       }
2504e2e9504SJed Brown     }
2514e2e9504SJed Brown     PetscCall(PetscFree(donor_minz));
2529ae3492bSJames Wright     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), &face_sfs[num_directions]));
2536497c311SBarry Smith     PetscCall(PetscSFSetGraph(face_sfs[num_directions], vEnd - vStart, inum_faces, NULL, PETSC_USE_POINTER, leaf, PETSC_USE_POINTER));
2544e2e9504SJed Brown     const PetscInt *my_donor_degree;
2559ae3492bSJames Wright     PetscCall(PetscSFComputeDegreeBegin(face_sfs[num_directions], &my_donor_degree));
2569ae3492bSJames Wright     PetscCall(PetscSFComputeDegreeEnd(face_sfs[num_directions], &my_donor_degree));
2576497c311SBarry Smith 
2584e2e9504SJed Brown     for (PetscInt i = 0; i < vEnd - vStart; i++) {
2594e2e9504SJed Brown       num_multiroots += my_donor_degree[i];
2604e2e9504SJed Brown       if (my_donor_degree[i] == 0) continue;
2615f06a3ddSJed Brown       PetscAssert(my_donor_degree[i] == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Local vertex has multiple faces");
2624e2e9504SJed Brown     }
2634e2e9504SJed Brown     PetscInt  *my_donors, *donor_indices, *my_donor_indices;
2646497c311SBarry Smith     PetscCount num_my_donors;
2656497c311SBarry Smith 
2669ae3492bSJames Wright     PetscCall(PetscSegBufferGetSize(my_donor_faces[direction], &num_my_donors));
267563af49cSJames Wright     PetscCheck(num_my_donors == num_multiroots, PETSC_COMM_SELF, PETSC_ERR_SUP, "Donor request (%" PetscCount_FMT ") does not match expected donors (%" PetscCount_FMT ")", num_my_donors, num_multiroots);
2689ae3492bSJames Wright     PetscCall(PetscSegBufferExtractInPlace(my_donor_faces[direction], &my_donors));
2694e2e9504SJed Brown     PetscCall(PetscMalloc1(vEnd - vStart, &my_donor_indices));
2706497c311SBarry Smith     for (PetscCount i = 0; i < num_my_donors; i++) {
2714e2e9504SJed Brown       PetscInt f = my_donors[i];
2721690c2aeSBarry Smith       PetscInt num_points, *points = NULL, minv = PETSC_INT_MAX;
2736497c311SBarry Smith 
2744e2e9504SJed Brown       PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &num_points, &points));
2754e2e9504SJed Brown       for (PetscInt j = 0; j < num_points; j++) {
2764e2e9504SJed Brown         PetscInt p = points[2 * j];
2774e2e9504SJed Brown         if (p < vStart || vEnd <= p) continue;
2784e2e9504SJed Brown         minv = PetscMin(minv, p);
2794e2e9504SJed Brown       }
2804e2e9504SJed Brown       PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &num_points, &points));
2815f06a3ddSJed Brown       PetscAssert(my_donor_degree[minv - vStart] == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Local vertex not requested");
2824e2e9504SJed Brown       my_donor_indices[minv - vStart] = f;
2834e2e9504SJed Brown     }
2844e2e9504SJed Brown     PetscCall(PetscMalloc1(num_faces, &donor_indices));
2859ae3492bSJames Wright     PetscCall(PetscSFBcastBegin(face_sfs[num_directions], MPIU_INT, my_donor_indices, donor_indices, MPI_REPLACE));
2869ae3492bSJames Wright     PetscCall(PetscSFBcastEnd(face_sfs[num_directions], MPIU_INT, my_donor_indices, donor_indices, MPI_REPLACE));
2874e2e9504SJed Brown     PetscCall(PetscFree(my_donor_indices));
2884e2e9504SJed Brown     // Modify our leafs so they point to donor faces instead of donor minz. Additionally, give them indices as faces.
2896497c311SBarry Smith     for (PetscCount i = 0; i < num_faces; i++) leaf[i].index = donor_indices[i];
2904e2e9504SJed Brown     PetscCall(PetscFree(donor_indices));
2914e2e9504SJed Brown     PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
2926497c311SBarry Smith     PetscCall(PetscSFSetGraph(face_sfs[num_directions], pEnd - pStart, inum_faces, faces, PETSC_COPY_VALUES, leaf, PETSC_OWN_POINTER));
2939ae3492bSJames Wright     {
2949ae3492bSJames Wright       char face_sf_name[PETSC_MAX_PATH_LEN];
2959ae3492bSJames Wright       PetscCall(PetscSNPrintf(face_sf_name, sizeof face_sf_name, "Z-order Isoperiodic Faces #%" PetscInt_FMT, num_directions));
2969ae3492bSJames Wright       PetscCall(PetscObjectSetName((PetscObject)face_sfs[num_directions], face_sf_name));
2979ae3492bSJames Wright     }
2984e2e9504SJed Brown 
2999ae3492bSJames Wright     transforms[num_directions][0][0]         = 1;
3009ae3492bSJames Wright     transforms[num_directions][1][1]         = 1;
3019ae3492bSJames Wright     transforms[num_directions][2][2]         = 1;
3029ae3492bSJames Wright     transforms[num_directions][3][3]         = 1;
3039ae3492bSJames Wright     transforms[num_directions][direction][3] = upper[direction] - lower[direction];
3049ae3492bSJames Wright     num_directions++;
3059ae3492bSJames Wright   }
3066725e60dSJed Brown 
3071fca310dSJames Wright   PetscCall(DMPlexSetIsoperiodicFaceSF(dm, num_directions, face_sfs));
3081fca310dSJames Wright   PetscCall(DMPlexSetIsoperiodicFaceTransform(dm, num_directions, (PetscScalar *)transforms));
3099ae3492bSJames Wright 
3109ae3492bSJames Wright   for (PetscInt i = 0; i < num_directions; i++) PetscCall(PetscSFDestroy(&face_sfs[i]));
3113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3124e2e9504SJed Brown }
3134e2e9504SJed Brown 
3145f06a3ddSJed Brown // This is a DMGlobalToLocalHook that applies the affine offsets. When extended for rotated periodicity, it'll need to
3155f06a3ddSJed Brown // apply a rotatonal transform and similar operations will be needed for fields (e.g., to rotate a velocity vector).
3165f06a3ddSJed Brown // We use this crude approach here so we don't have to write new GPU kernels yet.
3176725e60dSJed Brown static PetscErrorCode DMCoordAddPeriodicOffsets_Private(DM dm, Vec g, InsertMode mode, Vec l, void *ctx)
3186725e60dSJed Brown {
3196725e60dSJed Brown   PetscFunctionBegin;
320ddedc8f6SJames Wright   // These `VecScatter`s should be merged to improve efficiency; the scatters cannot be overlapped.
321ddedc8f6SJames Wright   for (PetscInt i = 0; i < dm->periodic.num_affines; i++) {
322ddedc8f6SJames Wright     PetscCall(VecScatterBegin(dm->periodic.affine_to_local[i], dm->periodic.affine[i], l, ADD_VALUES, SCATTER_FORWARD));
323ddedc8f6SJames Wright     PetscCall(VecScatterEnd(dm->periodic.affine_to_local[i], dm->periodic.affine[i], l, ADD_VALUES, SCATTER_FORWARD));
324ddedc8f6SJames Wright   }
3253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3266725e60dSJed Brown }
3276725e60dSJed Brown 
3285f06a3ddSJed Brown // Start with an SF for a positive depth (e.g., faces) and create a new SF for matched closure. The caller must ensure
3295f06a3ddSJed Brown // that both the donor (root) face and the periodic (leaf) face have consistent orientation, meaning that their closures
3305f06a3ddSJed Brown // are isomorphic. It may be useful/necessary for this restriction to be loosened.
3316725e60dSJed Brown //
3325f06a3ddSJed Brown // Output Arguments:
3335f06a3ddSJed Brown //
3345f06a3ddSJed Brown // + closure_sf - augmented point SF (see `DMGetPointSF()`) that includes the faces and all points in its closure. This
3355f06a3ddSJed Brown //   can be used to create a global section and section SF.
336b83f62b0SJames Wright // - is_points - array of index sets for just the points in the closure of `face_sf`. These may be used to apply an affine
3375f06a3ddSJed Brown //   transformation to periodic dofs; see DMPeriodicCoordinateSetUp_Internal().
3385f06a3ddSJed Brown //
339b83f62b0SJames Wright static PetscErrorCode DMPlexCreateIsoperiodicPointSF_Private(DM dm, PetscInt num_face_sfs, PetscSF *face_sfs, PetscSF *closure_sf, IS **is_points)
3406725e60dSJed Brown {
3416725e60dSJed Brown   MPI_Comm           comm;
3426725e60dSJed Brown   PetscMPIInt        rank;
3436725e60dSJed Brown   PetscSF            point_sf;
344b83f62b0SJames Wright   PetscInt           nroots, nleaves;
345b83f62b0SJames Wright   const PetscInt    *filocal;
346b83f62b0SJames Wright   const PetscSFNode *firemote;
3476725e60dSJed Brown 
3486725e60dSJed Brown   PetscFunctionBegin;
3496725e60dSJed Brown   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
3506725e60dSJed Brown   PetscCallMPI(MPI_Comm_rank(comm, &rank));
3516725e60dSJed Brown   PetscCall(DMGetPointSF(dm, &point_sf)); // Point SF has remote points
352b83f62b0SJames Wright   PetscCall(PetscMalloc1(num_face_sfs, is_points));
353b83f62b0SJames Wright 
354b83f62b0SJames Wright   for (PetscInt f = 0; f < num_face_sfs; f++) {
355b83f62b0SJames Wright     PetscSF   face_sf = face_sfs[f];
3566725e60dSJed Brown     PetscInt *rootdata, *leafdata;
357b83f62b0SJames Wright 
358b83f62b0SJames Wright     PetscCall(PetscSFGetGraph(face_sf, &nroots, &nleaves, &filocal, &firemote));
3596725e60dSJed Brown     PetscCall(PetscCalloc2(2 * nroots, &rootdata, 2 * nroots, &leafdata));
3606725e60dSJed Brown     for (PetscInt i = 0; i < nleaves; i++) {
3616725e60dSJed Brown       PetscInt point = filocal[i], cl_size, *closure = NULL;
3626725e60dSJed Brown       PetscCall(DMPlexGetTransitiveClosure(dm, point, PETSC_TRUE, &cl_size, &closure));
3636725e60dSJed Brown       leafdata[point] = cl_size - 1;
3646725e60dSJed Brown       PetscCall(DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, &cl_size, &closure));
3656725e60dSJed Brown     }
3666725e60dSJed Brown     PetscCall(PetscSFReduceBegin(face_sf, MPIU_INT, leafdata, rootdata + nroots, MPIU_SUM));
3676725e60dSJed Brown     PetscCall(PetscSFReduceEnd(face_sf, MPIU_INT, leafdata, rootdata + nroots, MPIU_SUM));
3686725e60dSJed Brown 
3696725e60dSJed Brown     PetscInt root_offset = 0;
3706725e60dSJed Brown     for (PetscInt p = 0; p < nroots; p++) {
3716725e60dSJed Brown       const PetscInt *donor_dof = rootdata + nroots;
3726725e60dSJed Brown       if (donor_dof[p] == 0) {
3736725e60dSJed Brown         rootdata[2 * p]     = -1;
3746725e60dSJed Brown         rootdata[2 * p + 1] = -1;
3756725e60dSJed Brown         continue;
3766725e60dSJed Brown       }
3776725e60dSJed Brown       PetscInt  cl_size;
3786725e60dSJed Brown       PetscInt *closure = NULL;
3796725e60dSJed Brown       PetscCall(DMPlexGetTransitiveClosure(dm, p, PETSC_TRUE, &cl_size, &closure));
3806725e60dSJed Brown       // cl_size - 1 = points not including self
3815f06a3ddSJed Brown       PetscAssert(donor_dof[p] == cl_size - 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Reduced leaf cone sizes do not match root cone sizes");
3826725e60dSJed Brown       rootdata[2 * p]     = root_offset;
3836725e60dSJed Brown       rootdata[2 * p + 1] = cl_size - 1;
3846725e60dSJed Brown       root_offset += cl_size - 1;
3856725e60dSJed Brown       PetscCall(DMPlexRestoreTransitiveClosure(dm, p, PETSC_TRUE, &cl_size, &closure));
3866725e60dSJed Brown     }
3876725e60dSJed Brown     PetscCall(PetscSFBcastBegin(face_sf, MPIU_2INT, rootdata, leafdata, MPI_REPLACE));
3886725e60dSJed Brown     PetscCall(PetscSFBcastEnd(face_sf, MPIU_2INT, rootdata, leafdata, MPI_REPLACE));
3896725e60dSJed Brown     // Count how many leaves we need to communicate the closures
3906725e60dSJed Brown     PetscInt leaf_offset = 0;
3916725e60dSJed Brown     for (PetscInt i = 0; i < nleaves; i++) {
3926725e60dSJed Brown       PetscInt point = filocal[i];
3936725e60dSJed Brown       if (leafdata[2 * point + 1] < 0) continue;
3946725e60dSJed Brown       leaf_offset += leafdata[2 * point + 1];
3956725e60dSJed Brown     }
3966725e60dSJed Brown 
3976725e60dSJed Brown     PetscSFNode *closure_leaf;
3986725e60dSJed Brown     PetscCall(PetscMalloc1(leaf_offset, &closure_leaf));
3996725e60dSJed Brown     leaf_offset = 0;
4006725e60dSJed Brown     for (PetscInt i = 0; i < nleaves; i++) {
4016725e60dSJed Brown       PetscInt point   = filocal[i];
4026725e60dSJed Brown       PetscInt cl_size = leafdata[2 * point + 1];
4036725e60dSJed Brown       if (cl_size < 0) continue;
4046725e60dSJed Brown       for (PetscInt j = 0; j < cl_size; j++) {
4056725e60dSJed Brown         closure_leaf[leaf_offset].rank  = firemote[i].rank;
4066725e60dSJed Brown         closure_leaf[leaf_offset].index = leafdata[2 * point] + j;
4076725e60dSJed Brown         leaf_offset++;
4086725e60dSJed Brown       }
4096725e60dSJed Brown     }
4106725e60dSJed Brown 
4116725e60dSJed Brown     PetscSF sf_closure;
4126725e60dSJed Brown     PetscCall(PetscSFCreate(comm, &sf_closure));
4136725e60dSJed Brown     PetscCall(PetscSFSetGraph(sf_closure, root_offset, leaf_offset, NULL, PETSC_USE_POINTER, closure_leaf, PETSC_OWN_POINTER));
4146725e60dSJed Brown 
415b83f62b0SJames Wright     PetscSFNode *leaf_donor_closure;
416b83f62b0SJames Wright     { // Pack root buffer with owner for every point in the root cones
4176725e60dSJed Brown       PetscSFNode       *donor_closure;
418b83f62b0SJames Wright       const PetscInt    *pilocal;
419b83f62b0SJames Wright       const PetscSFNode *piremote;
420b83f62b0SJames Wright       PetscInt           npoints;
421b83f62b0SJames Wright 
422b83f62b0SJames Wright       PetscCall(PetscSFGetGraph(point_sf, NULL, &npoints, &pilocal, &piremote));
4236725e60dSJed Brown       PetscCall(PetscCalloc1(root_offset, &donor_closure));
4246725e60dSJed Brown       root_offset = 0;
4256725e60dSJed Brown       for (PetscInt p = 0; p < nroots; p++) {
4266725e60dSJed Brown         if (rootdata[2 * p] < 0) continue;
4276725e60dSJed Brown         PetscInt  cl_size;
4286725e60dSJed Brown         PetscInt *closure = NULL;
4296725e60dSJed Brown         PetscCall(DMPlexGetTransitiveClosure(dm, p, PETSC_TRUE, &cl_size, &closure));
4306725e60dSJed Brown         for (PetscInt j = 1; j < cl_size; j++) {
4316725e60dSJed Brown           PetscInt c = closure[2 * j];
4326725e60dSJed Brown           if (pilocal) {
4336725e60dSJed Brown             PetscInt found = -1;
4346725e60dSJed Brown             if (npoints > 0) PetscCall(PetscFindInt(c, npoints, pilocal, &found));
4356725e60dSJed Brown             if (found >= 0) {
4366725e60dSJed Brown               donor_closure[root_offset++] = piremote[found];
4376725e60dSJed Brown               continue;
4386725e60dSJed Brown             }
4396725e60dSJed Brown           }
4406725e60dSJed Brown           // we own c
4416725e60dSJed Brown           donor_closure[root_offset].rank  = rank;
4426725e60dSJed Brown           donor_closure[root_offset].index = c;
4436725e60dSJed Brown           root_offset++;
4446725e60dSJed Brown         }
4456725e60dSJed Brown         PetscCall(DMPlexRestoreTransitiveClosure(dm, p, PETSC_TRUE, &cl_size, &closure));
4466725e60dSJed Brown       }
4476725e60dSJed Brown 
4486725e60dSJed Brown       PetscCall(PetscMalloc1(leaf_offset, &leaf_donor_closure));
4496497c311SBarry Smith       PetscCall(PetscSFBcastBegin(sf_closure, MPIU_SF_NODE, donor_closure, leaf_donor_closure, MPI_REPLACE));
4506497c311SBarry Smith       PetscCall(PetscSFBcastEnd(sf_closure, MPIU_SF_NODE, donor_closure, leaf_donor_closure, MPI_REPLACE));
4516725e60dSJed Brown       PetscCall(PetscSFDestroy(&sf_closure));
4526725e60dSJed Brown       PetscCall(PetscFree(donor_closure));
453b83f62b0SJames Wright     }
4546725e60dSJed Brown 
4556725e60dSJed Brown     PetscSFNode *new_iremote;
4566725e60dSJed Brown     PetscCall(PetscCalloc1(nroots, &new_iremote));
4576725e60dSJed Brown     for (PetscInt i = 0; i < nroots; i++) new_iremote[i].rank = -1;
4586725e60dSJed Brown     // Walk leaves and match vertices
4596725e60dSJed Brown     leaf_offset = 0;
4606725e60dSJed Brown     for (PetscInt i = 0; i < nleaves; i++) {
4616725e60dSJed Brown       PetscInt  point   = filocal[i], cl_size;
4626725e60dSJed Brown       PetscInt *closure = NULL;
4636725e60dSJed Brown       PetscCall(DMPlexGetTransitiveClosure(dm, point, PETSC_TRUE, &cl_size, &closure));
4646725e60dSJed Brown       for (PetscInt j = 1; j < cl_size; j++) { // TODO: should we send donor edge orientations so we can flip for consistency?
4656725e60dSJed Brown         PetscInt    c  = closure[2 * j];
4666725e60dSJed Brown         PetscSFNode lc = leaf_donor_closure[leaf_offset];
4676725e60dSJed Brown         // printf("[%d] face %d.%d: %d ?-- (%d,%d)\n", rank, point, j, c, lc.rank, lc.index);
4686725e60dSJed Brown         if (new_iremote[c].rank == -1) {
4696725e60dSJed Brown           new_iremote[c] = lc;
4705f06a3ddSJed Brown         } else PetscCheck(new_iremote[c].rank == lc.rank && new_iremote[c].index == lc.index, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Mismatched cone ordering between faces");
4716725e60dSJed Brown         leaf_offset++;
4726725e60dSJed Brown       }
4736725e60dSJed Brown       PetscCall(DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, &cl_size, &closure));
4746725e60dSJed Brown     }
4756725e60dSJed Brown     PetscCall(PetscFree(leaf_donor_closure));
4766725e60dSJed Brown 
4776725e60dSJed Brown     // Include face points in closure SF
4786725e60dSJed Brown     for (PetscInt i = 0; i < nleaves; i++) new_iremote[filocal[i]] = firemote[i];
4796725e60dSJed Brown     // consolidate leaves
4806725e60dSJed Brown     PetscInt num_new_leaves = 0;
4816725e60dSJed Brown     for (PetscInt i = 0; i < nroots; i++) {
4826725e60dSJed Brown       if (new_iremote[i].rank == -1) continue;
4836725e60dSJed Brown       new_iremote[num_new_leaves] = new_iremote[i];
4846725e60dSJed Brown       leafdata[num_new_leaves]    = i;
4856725e60dSJed Brown       num_new_leaves++;
4866725e60dSJed Brown     }
487b83f62b0SJames Wright     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, num_new_leaves, leafdata, PETSC_COPY_VALUES, &(*is_points)[f]));
4886725e60dSJed Brown 
4896725e60dSJed Brown     PetscSF csf;
4906725e60dSJed Brown     PetscCall(PetscSFCreate(comm, &csf));
4916725e60dSJed Brown     PetscCall(PetscSFSetGraph(csf, nroots, num_new_leaves, leafdata, PETSC_COPY_VALUES, new_iremote, PETSC_COPY_VALUES));
4926725e60dSJed Brown     PetscCall(PetscFree(new_iremote)); // copy and delete because new_iremote is longer than it needs to be
4936725e60dSJed Brown     PetscCall(PetscFree2(rootdata, leafdata));
4946725e60dSJed Brown 
495b83f62b0SJames Wright     PetscInt npoints;
496b83f62b0SJames Wright     PetscCall(PetscSFGetGraph(point_sf, NULL, &npoints, NULL, NULL));
4975f06a3ddSJed Brown     if (npoints < 0) { // empty point_sf
4986725e60dSJed Brown       *closure_sf = csf;
4995f06a3ddSJed Brown     } else {
5005f06a3ddSJed Brown       PetscCall(PetscSFMerge(point_sf, csf, closure_sf));
5015f06a3ddSJed Brown       PetscCall(PetscSFDestroy(&csf));
5025f06a3ddSJed Brown     }
503d8b4a066SPierre Jolivet     if (f > 0) PetscCall(PetscSFDestroy(&point_sf)); // Only destroy if point_sf is from previous calls to PetscSFMerge
504b83f62b0SJames Wright     point_sf = *closure_sf;                          // Use combined point + isoperiodic SF to define point ownership for further face_sf
505b83f62b0SJames Wright   }
5065f06a3ddSJed Brown   PetscCall(PetscObjectSetName((PetscObject)*closure_sf, "Composed Periodic Points"));
5073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5086725e60dSJed Brown }
5096725e60dSJed Brown 
5105f06a3ddSJed Brown static PetscErrorCode DMGetIsoperiodicPointSF_Plex(DM dm, PetscSF *sf)
5116725e60dSJed Brown {
5126725e60dSJed Brown   DM_Plex *plex = (DM_Plex *)dm->data;
5136725e60dSJed Brown 
5146725e60dSJed Brown   PetscFunctionBegin;
515b83f62b0SJames Wright   if (!plex->periodic.composed_sf) PetscCall(DMPlexCreateIsoperiodicPointSF_Private(dm, plex->periodic.num_face_sfs, plex->periodic.face_sfs, &plex->periodic.composed_sf, &plex->periodic.periodic_points));
5166725e60dSJed Brown   if (sf) *sf = plex->periodic.composed_sf;
5173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5186725e60dSJed Brown }
5196725e60dSJed Brown 
5205f06a3ddSJed Brown PetscErrorCode DMPlexMigrateIsoperiodicFaceSF_Internal(DM old_dm, DM dm, PetscSF sf_migration)
5215f06a3ddSJed Brown {
5225f06a3ddSJed Brown   DM_Plex    *plex = (DM_Plex *)old_dm->data;
523a45b0079SJames Wright   PetscSF     sf_point, *new_face_sfs;
5245f06a3ddSJed Brown   PetscMPIInt rank;
5255f06a3ddSJed Brown 
5265f06a3ddSJed Brown   PetscFunctionBegin;
5271fca310dSJames Wright   if (!plex->periodic.face_sfs) PetscFunctionReturn(PETSC_SUCCESS);
5285f06a3ddSJed Brown   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
5295f06a3ddSJed Brown   PetscCall(DMGetPointSF(dm, &sf_point));
530a45b0079SJames Wright   PetscCall(PetscMalloc1(plex->periodic.num_face_sfs, &new_face_sfs));
531a45b0079SJames Wright 
532a45b0079SJames Wright   for (PetscInt f = 0; f < plex->periodic.num_face_sfs; f++) {
5335f06a3ddSJed Brown     PetscInt           old_npoints, new_npoints, old_nleaf, new_nleaf, point_nleaf;
5345f06a3ddSJed Brown     PetscSFNode       *new_leafdata, *rootdata, *leafdata;
5355f06a3ddSJed Brown     const PetscInt    *old_local, *point_local;
5365f06a3ddSJed Brown     const PetscSFNode *old_remote, *point_remote;
5376497c311SBarry Smith 
538a45b0079SJames Wright     PetscCall(PetscSFGetGraph(plex->periodic.face_sfs[f], &old_npoints, &old_nleaf, &old_local, &old_remote));
5395f06a3ddSJed Brown     PetscCall(PetscSFGetGraph(sf_migration, NULL, &new_nleaf, NULL, NULL));
5405f06a3ddSJed Brown     PetscCall(PetscSFGetGraph(sf_point, &new_npoints, &point_nleaf, &point_local, &point_remote));
5415f06a3ddSJed Brown     PetscAssert(new_nleaf == new_npoints, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected migration leaf space to match new point root space");
5425f06a3ddSJed Brown     PetscCall(PetscMalloc3(old_npoints, &rootdata, old_npoints, &leafdata, new_npoints, &new_leafdata));
5435f06a3ddSJed Brown 
5445f06a3ddSJed Brown     // Fill new_leafdata with new owners of all points
5455f06a3ddSJed Brown     for (PetscInt i = 0; i < new_npoints; i++) {
5465f06a3ddSJed Brown       new_leafdata[i].rank  = rank;
5475f06a3ddSJed Brown       new_leafdata[i].index = i;
5485f06a3ddSJed Brown     }
5495f06a3ddSJed Brown     for (PetscInt i = 0; i < point_nleaf; i++) {
5505f06a3ddSJed Brown       PetscInt j      = point_local[i];
5515f06a3ddSJed Brown       new_leafdata[j] = point_remote[i];
5525f06a3ddSJed Brown     }
5535f06a3ddSJed Brown     // REPLACE is okay because every leaf agrees about the new owners
5546497c311SBarry Smith     PetscCall(PetscSFReduceBegin(sf_migration, MPIU_SF_NODE, new_leafdata, rootdata, MPI_REPLACE));
5556497c311SBarry Smith     PetscCall(PetscSFReduceEnd(sf_migration, MPIU_SF_NODE, new_leafdata, rootdata, MPI_REPLACE));
5565f06a3ddSJed Brown     // rootdata now contains the new owners
5575f06a3ddSJed Brown 
5585f06a3ddSJed Brown     // Send to leaves of old space
5595f06a3ddSJed Brown     for (PetscInt i = 0; i < old_npoints; i++) {
5605f06a3ddSJed Brown       leafdata[i].rank  = -1;
5615f06a3ddSJed Brown       leafdata[i].index = -1;
5625f06a3ddSJed Brown     }
5636497c311SBarry Smith     PetscCall(PetscSFBcastBegin(plex->periodic.face_sfs[f], MPIU_SF_NODE, rootdata, leafdata, MPI_REPLACE));
5646497c311SBarry Smith     PetscCall(PetscSFBcastEnd(plex->periodic.face_sfs[f], MPIU_SF_NODE, rootdata, leafdata, MPI_REPLACE));
5655f06a3ddSJed Brown 
5665f06a3ddSJed Brown     // Send to new leaf space
5676497c311SBarry Smith     PetscCall(PetscSFBcastBegin(sf_migration, MPIU_SF_NODE, leafdata, new_leafdata, MPI_REPLACE));
5686497c311SBarry Smith     PetscCall(PetscSFBcastEnd(sf_migration, MPIU_SF_NODE, leafdata, new_leafdata, MPI_REPLACE));
5695f06a3ddSJed Brown 
5705f06a3ddSJed Brown     PetscInt     nface = 0, *new_local;
5715f06a3ddSJed Brown     PetscSFNode *new_remote;
5725f06a3ddSJed Brown     for (PetscInt i = 0; i < new_npoints; i++) nface += (new_leafdata[i].rank >= 0);
5735f06a3ddSJed Brown     PetscCall(PetscMalloc1(nface, &new_local));
5745f06a3ddSJed Brown     PetscCall(PetscMalloc1(nface, &new_remote));
5755f06a3ddSJed Brown     nface = 0;
5765f06a3ddSJed Brown     for (PetscInt i = 0; i < new_npoints; i++) {
5775f06a3ddSJed Brown       if (new_leafdata[i].rank == -1) continue;
5785f06a3ddSJed Brown       new_local[nface]  = i;
5795f06a3ddSJed Brown       new_remote[nface] = new_leafdata[i];
5805f06a3ddSJed Brown       nface++;
5815f06a3ddSJed Brown     }
5825f06a3ddSJed Brown     PetscCall(PetscFree3(rootdata, leafdata, new_leafdata));
583a45b0079SJames Wright     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), &new_face_sfs[f]));
584a45b0079SJames Wright     PetscCall(PetscSFSetGraph(new_face_sfs[f], new_npoints, nface, new_local, PETSC_OWN_POINTER, new_remote, PETSC_OWN_POINTER));
585a45b0079SJames Wright     {
586a45b0079SJames Wright       char new_face_sf_name[PETSC_MAX_PATH_LEN];
587a45b0079SJames Wright       PetscCall(PetscSNPrintf(new_face_sf_name, sizeof new_face_sf_name, "Migrated Isoperiodic Faces #%" PetscInt_FMT, f));
588a45b0079SJames Wright       PetscCall(PetscObjectSetName((PetscObject)new_face_sfs[f], new_face_sf_name));
589a45b0079SJames Wright     }
590a45b0079SJames Wright   }
591a45b0079SJames Wright 
592a45b0079SJames Wright   PetscCall(DMPlexSetIsoperiodicFaceSF(dm, plex->periodic.num_face_sfs, new_face_sfs));
593a45b0079SJames Wright   PetscCall(DMPlexSetIsoperiodicFaceTransform(dm, plex->periodic.num_face_sfs, (PetscScalar *)plex->periodic.transform));
594a45b0079SJames Wright   for (PetscInt f = 0; f < plex->periodic.num_face_sfs; f++) PetscCall(PetscSFDestroy(&new_face_sfs[f]));
595a45b0079SJames Wright   PetscCall(PetscFree(new_face_sfs));
5963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5975f06a3ddSJed Brown }
5985f06a3ddSJed Brown 
5996725e60dSJed Brown PetscErrorCode DMPeriodicCoordinateSetUp_Internal(DM dm)
6006725e60dSJed Brown {
6016725e60dSJed Brown   DM_Plex   *plex = (DM_Plex *)dm->data;
6026497c311SBarry Smith   PetscCount count;
603ddedc8f6SJames Wright   IS         isdof;
604ddedc8f6SJames Wright   PetscInt   dim;
6054d86920dSPierre Jolivet 
6066725e60dSJed Brown   PetscFunctionBegin;
6071fca310dSJames Wright   if (!plex->periodic.face_sfs) PetscFunctionReturn(PETSC_SUCCESS);
6085f06a3ddSJed Brown   PetscCall(DMGetIsoperiodicPointSF_Plex(dm, NULL));
6095f06a3ddSJed Brown   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMGetIsoperiodicPointSF_C", DMGetIsoperiodicPointSF_Plex));
6106725e60dSJed Brown 
6116725e60dSJed Brown   PetscCall(DMGetDimension(dm, &dim));
612ddedc8f6SJames Wright   dm->periodic.num_affines = plex->periodic.num_face_sfs;
613ddedc8f6SJames Wright   PetscCall(PetscMalloc2(dm->periodic.num_affines, &dm->periodic.affine_to_local, dm->periodic.num_affines, &dm->periodic.affine));
614ddedc8f6SJames Wright 
615ddedc8f6SJames Wright   for (PetscInt f = 0; f < plex->periodic.num_face_sfs; f++) {
6166497c311SBarry Smith     PetscInt        npoints, vsize, isize;
6176725e60dSJed Brown     const PetscInt *points;
618ddedc8f6SJames Wright     IS              is = plex->periodic.periodic_points[f];
6196725e60dSJed Brown     PetscSegBuffer  seg;
6206725e60dSJed Brown     PetscSection    section;
6216497c311SBarry Smith     PetscInt       *ind;
6226497c311SBarry Smith     Vec             L, P;
6236497c311SBarry Smith     VecType         vec_type;
6246497c311SBarry Smith     VecScatter      scatter;
6256497c311SBarry Smith     PetscScalar    *x;
6266497c311SBarry Smith 
6276725e60dSJed Brown     PetscCall(DMGetLocalSection(dm, &section));
6286725e60dSJed Brown     PetscCall(PetscSegBufferCreate(sizeof(PetscInt), 32, &seg));
6296725e60dSJed Brown     PetscCall(ISGetSize(is, &npoints));
6306725e60dSJed Brown     PetscCall(ISGetIndices(is, &points));
6316725e60dSJed Brown     for (PetscInt i = 0; i < npoints; i++) {
6326725e60dSJed Brown       PetscInt point = points[i], off, dof;
6336497c311SBarry Smith 
6346725e60dSJed Brown       PetscCall(PetscSectionGetOffset(section, point, &off));
6356725e60dSJed Brown       PetscCall(PetscSectionGetDof(section, point, &dof));
6366725e60dSJed Brown       PetscAssert(dof % dim == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unexpected dof %" PetscInt_FMT " not divisible by dimension %" PetscInt_FMT, dof, dim);
6376725e60dSJed Brown       for (PetscInt j = 0; j < dof / dim; j++) {
6386725e60dSJed Brown         PetscInt *slot;
6396497c311SBarry Smith 
6406725e60dSJed Brown         PetscCall(PetscSegBufferGetInts(seg, 1, &slot));
641729bdd54SJed Brown         *slot = off / dim + j;
6426725e60dSJed Brown       }
6436725e60dSJed Brown     }
6446725e60dSJed Brown     PetscCall(PetscSegBufferGetSize(seg, &count));
6456725e60dSJed Brown     PetscCall(PetscSegBufferExtractAlloc(seg, &ind));
6466725e60dSJed Brown     PetscCall(PetscSegBufferDestroy(&seg));
6476497c311SBarry Smith     PetscCall(PetscIntCast(count, &isize));
6486497c311SBarry Smith     PetscCall(ISCreateBlock(PETSC_COMM_SELF, dim, isize, ind, PETSC_OWN_POINTER, &isdof));
6496497c311SBarry Smith 
6506497c311SBarry Smith     PetscCall(PetscIntCast(count * dim, &vsize));
6516725e60dSJed Brown     PetscCall(DMGetLocalVector(dm, &L));
6526725e60dSJed Brown     PetscCall(VecCreate(PETSC_COMM_SELF, &P));
6536497c311SBarry Smith     PetscCall(VecSetSizes(P, vsize, vsize));
6546725e60dSJed Brown     PetscCall(VecGetType(L, &vec_type));
6556725e60dSJed Brown     PetscCall(VecSetType(P, vec_type));
6566725e60dSJed Brown     PetscCall(VecScatterCreate(P, NULL, L, isdof, &scatter));
6576725e60dSJed Brown     PetscCall(DMRestoreLocalVector(dm, &L));
6586725e60dSJed Brown     PetscCall(ISDestroy(&isdof));
6596725e60dSJed Brown 
6606725e60dSJed Brown     PetscCall(VecGetArrayWrite(P, &x));
6616497c311SBarry Smith     for (PetscCount i = 0; i < count; i++) {
662ddedc8f6SJames Wright       for (PetscInt j = 0; j < dim; j++) x[i * dim + j] = plex->periodic.transform[f][j][3];
6636725e60dSJed Brown     }
6646725e60dSJed Brown     PetscCall(VecRestoreArrayWrite(P, &x));
6656725e60dSJed Brown 
666ddedc8f6SJames Wright     dm->periodic.affine_to_local[f] = scatter;
667ddedc8f6SJames Wright     dm->periodic.affine[f]          = P;
668ddedc8f6SJames Wright   }
6696725e60dSJed Brown   PetscCall(DMGlobalToLocalHookAdd(dm, NULL, DMCoordAddPeriodicOffsets_Private, NULL));
6703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6716725e60dSJed Brown }
6726725e60dSJed Brown 
6735f06a3ddSJed Brown // We'll just orient all the edges, though only periodic boundary edges need orientation
6745f06a3ddSJed Brown static PetscErrorCode DMPlexOrientPositiveEdges_Private(DM dm)
6755f06a3ddSJed Brown {
6765f06a3ddSJed Brown   PetscInt dim, eStart, eEnd;
6774d86920dSPierre Jolivet 
6785f06a3ddSJed Brown   PetscFunctionBegin;
6795f06a3ddSJed Brown   PetscCall(DMGetDimension(dm, &dim));
6803ba16761SJacob Faibussowitsch   if (dim < 3) PetscFunctionReturn(PETSC_SUCCESS); // not necessary
6815f06a3ddSJed Brown   PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd));
6825f06a3ddSJed Brown   for (PetscInt e = eStart; e < eEnd; e++) {
6835f06a3ddSJed Brown     const PetscInt *cone;
6845f06a3ddSJed Brown     PetscCall(DMPlexGetCone(dm, e, &cone));
6855f06a3ddSJed Brown     if (cone[0] > cone[1]) PetscCall(DMPlexOrientPoint(dm, e, -1));
6865f06a3ddSJed Brown   }
6873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6885f06a3ddSJed Brown }
6895f06a3ddSJed Brown 
6903e72e933SJed Brown PetscErrorCode DMPlexCreateBoxMesh_Tensor_SFC_Internal(DM dm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate)
6913e72e933SJed Brown {
6923e72e933SJed Brown   PetscInt  eextent[3] = {1, 1, 1}, vextent[3] = {1, 1, 1};
6933e72e933SJed Brown   const Ijk closure_1[] = {
6943e72e933SJed Brown     {0, 0, 0},
6953e72e933SJed Brown     {1, 0, 0},
6963e72e933SJed Brown   };
6973e72e933SJed Brown   const Ijk closure_2[] = {
6983e72e933SJed Brown     {0, 0, 0},
6993e72e933SJed Brown     {1, 0, 0},
7003e72e933SJed Brown     {1, 1, 0},
7013e72e933SJed Brown     {0, 1, 0},
7023e72e933SJed Brown   };
7033e72e933SJed Brown   const Ijk closure_3[] = {
7043e72e933SJed Brown     {0, 0, 0},
7053e72e933SJed Brown     {0, 1, 0},
7063e72e933SJed Brown     {1, 1, 0},
7073e72e933SJed Brown     {1, 0, 0},
7083e72e933SJed Brown     {0, 0, 1},
7093e72e933SJed Brown     {1, 0, 1},
7103e72e933SJed Brown     {1, 1, 1},
7113e72e933SJed Brown     {0, 1, 1},
7123e72e933SJed Brown   };
7133431e603SJed Brown   const Ijk *const closure_dim[] = {NULL, closure_1, closure_2, closure_3};
7143431e603SJed Brown   // This must be kept consistent with DMPlexCreateCubeMesh_Internal
7153431e603SJed Brown   const PetscInt        face_marker_1[]   = {1, 2};
7163431e603SJed Brown   const PetscInt        face_marker_2[]   = {4, 2, 1, 3};
7173431e603SJed Brown   const PetscInt        face_marker_3[]   = {6, 5, 3, 4, 1, 2};
7183431e603SJed Brown   const PetscInt *const face_marker_dim[] = {NULL, face_marker_1, face_marker_2, face_marker_3};
7196725e60dSJed Brown   // Orient faces so the normal is in the positive axis and the first vertex is the one closest to zero.
7206725e60dSJed Brown   // These orientations can be determined by examining cones of a reference quad and hex element.
7216725e60dSJed Brown   const PetscInt        face_orient_1[]   = {0, 0};
7226725e60dSJed Brown   const PetscInt        face_orient_2[]   = {-1, 0, 0, -1};
7236725e60dSJed Brown   const PetscInt        face_orient_3[]   = {-2, 0, -2, 1, -2, 0};
7246725e60dSJed Brown   const PetscInt *const face_orient_dim[] = {NULL, face_orient_1, face_orient_2, face_orient_3};
7253e72e933SJed Brown 
7263e72e933SJed Brown   PetscFunctionBegin;
727*ea7b3863SJames Wright   PetscCall(PetscLogEventBegin(DMPLEX_CreateBoxSFC, dm, 0, 0, 0));
7284f572ea9SToby Isaac   PetscAssertPointer(dm, 1);
7293e72e933SJed Brown   PetscValidLogicalCollectiveInt(dm, dim, 2);
7303e72e933SJed Brown   PetscCall(DMSetDimension(dm, dim));
7313e72e933SJed Brown   PetscMPIInt rank, size;
7323e72e933SJed Brown   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
7333e72e933SJed Brown   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
7343e72e933SJed Brown   for (PetscInt i = 0; i < dim; i++) {
7353e72e933SJed Brown     eextent[i] = faces[i];
7363e72e933SJed Brown     vextent[i] = faces[i] + 1;
7373e72e933SJed Brown   }
738c77877e3SJed Brown   ZLayout layout;
739c77877e3SJed Brown   PetscCall(ZLayoutCreate(size, eextent, vextent, &layout));
7403e72e933SJed Brown   PetscZSet vset; // set of all vertices in the closure of the owned elements
7413e72e933SJed Brown   PetscCall(PetscZSetCreate(&vset));
7423e72e933SJed Brown   PetscInt local_elems = 0;
7433e72e933SJed Brown   for (ZCode z = layout.zstarts[rank]; z < layout.zstarts[rank + 1]; z++) {
7443e72e933SJed Brown     Ijk loc = ZCodeSplit(z);
7453ba16761SJacob Faibussowitsch     if (IjkActive(layout.vextent, loc)) PetscCall(PetscZSetAdd(vset, z));
7466547ddbdSJed Brown     else {
7476547ddbdSJed Brown       z += ZStepOct(z);
7486547ddbdSJed Brown       continue;
7496547ddbdSJed Brown     }
7503e72e933SJed Brown     if (IjkActive(layout.eextent, loc)) {
7513e72e933SJed Brown       local_elems++;
7523e72e933SJed Brown       // Add all neighboring vertices to set
7533e72e933SJed Brown       for (PetscInt n = 0; n < PetscPowInt(2, dim); n++) {
7543e72e933SJed Brown         Ijk   inc  = closure_dim[dim][n];
7553e72e933SJed Brown         Ijk   nloc = {loc.i + inc.i, loc.j + inc.j, loc.k + inc.k};
7563e72e933SJed Brown         ZCode v    = ZEncode(nloc);
7573ba16761SJacob Faibussowitsch         PetscCall(PetscZSetAdd(vset, v));
7583e72e933SJed Brown       }
7593e72e933SJed Brown     }
7603e72e933SJed Brown   }
7613e72e933SJed Brown   PetscInt local_verts, off = 0;
7623e72e933SJed Brown   ZCode   *vert_z;
7633e72e933SJed Brown   PetscCall(PetscZSetGetSize(vset, &local_verts));
7643e72e933SJed Brown   PetscCall(PetscMalloc1(local_verts, &vert_z));
7653e72e933SJed Brown   PetscCall(PetscZSetGetElems(vset, &off, vert_z));
7663e72e933SJed Brown   PetscCall(PetscZSetDestroy(&vset));
7673e72e933SJed Brown   // ZCode is unsigned for bitwise convenience, but highest bit should never be set, so can interpret as signed
7683e72e933SJed Brown   PetscCall(PetscSortInt64(local_verts, (PetscInt64 *)vert_z));
7693e72e933SJed Brown 
7703e72e933SJed Brown   PetscCall(DMPlexSetChart(dm, 0, local_elems + local_verts));
7713e72e933SJed Brown   for (PetscInt e = 0; e < local_elems; e++) PetscCall(DMPlexSetConeSize(dm, e, PetscPowInt(2, dim)));
7723e72e933SJed Brown   PetscCall(DMSetUp(dm));
7733e72e933SJed Brown   {
7743e72e933SJed Brown     PetscInt e = 0;
7753e72e933SJed Brown     for (ZCode z = layout.zstarts[rank]; z < layout.zstarts[rank + 1]; z++) {
7763e72e933SJed Brown       Ijk loc = ZCodeSplit(z);
7776547ddbdSJed Brown       if (!IjkActive(layout.eextent, loc)) {
7786547ddbdSJed Brown         z += ZStepOct(z);
7796547ddbdSJed Brown         continue;
7806547ddbdSJed Brown       }
7813e72e933SJed Brown       PetscInt cone[8], orient[8] = {0};
7823e72e933SJed Brown       for (PetscInt n = 0; n < PetscPowInt(2, dim); n++) {
7833e72e933SJed Brown         Ijk      inc  = closure_dim[dim][n];
7843e72e933SJed Brown         Ijk      nloc = {loc.i + inc.i, loc.j + inc.j, loc.k + inc.k};
7853e72e933SJed Brown         ZCode    v    = ZEncode(nloc);
7863e72e933SJed Brown         PetscInt ci   = ZCodeFind(v, local_verts, vert_z);
7873e72e933SJed Brown         PetscAssert(ci >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Did not find neighbor vertex in set");
7883e72e933SJed Brown         cone[n] = local_elems + ci;
7893e72e933SJed Brown       }
7903e72e933SJed Brown       PetscCall(DMPlexSetCone(dm, e, cone));
7913e72e933SJed Brown       PetscCall(DMPlexSetConeOrientation(dm, e, orient));
7923e72e933SJed Brown       e++;
7933e72e933SJed Brown     }
7943e72e933SJed Brown   }
7953e72e933SJed Brown 
7963e72e933SJed Brown   PetscCall(DMPlexSymmetrize(dm));
7973e72e933SJed Brown   PetscCall(DMPlexStratify(dm));
7986725e60dSJed Brown 
7993e72e933SJed Brown   { // Create point SF
8003e72e933SJed Brown     PetscSF sf;
8013ba16761SJacob Faibussowitsch     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), &sf));
8023e72e933SJed Brown     PetscInt owned_verts = ZCodeFind(layout.zstarts[rank + 1], local_verts, vert_z);
8033e72e933SJed Brown     if (owned_verts < 0) owned_verts = -(owned_verts + 1); // We don't care whether the key was found
8043e72e933SJed Brown     PetscInt     num_ghosts = local_verts - owned_verts;   // Due to sorting, owned vertices always come first
8053e72e933SJed Brown     PetscInt    *local_ghosts;
8063e72e933SJed Brown     PetscSFNode *ghosts;
8073e72e933SJed Brown     PetscCall(PetscMalloc1(num_ghosts, &local_ghosts));
8083e72e933SJed Brown     PetscCall(PetscMalloc1(num_ghosts, &ghosts));
8093e72e933SJed Brown     for (PetscInt i = 0; i < num_ghosts;) {
8103e72e933SJed Brown       ZCode       z           = vert_z[owned_verts + i];
8116497c311SBarry Smith       PetscMPIInt remote_rank = (PetscMPIInt)ZCodeFind(z, size + 1, layout.zstarts), remote_count = 0;
8123e72e933SJed Brown       if (remote_rank < 0) remote_rank = -(remote_rank + 1) - 1;
8133e72e933SJed Brown       // We have a new remote rank; find all the ghost indices (which are contiguous in vert_z)
8143e72e933SJed Brown 
8153e72e933SJed Brown       // Count the elements on remote_rank
8164e2e9504SJed Brown       PetscInt remote_elem = ZLayoutElementsOnRank(&layout, remote_rank);
8173e72e933SJed Brown 
8183e72e933SJed Brown       // Traverse vertices and make ghost links
8193e72e933SJed Brown       for (ZCode rz = layout.zstarts[remote_rank]; rz < layout.zstarts[remote_rank + 1]; rz++) {
8203e72e933SJed Brown         Ijk loc = ZCodeSplit(rz);
8213e72e933SJed Brown         if (rz == z) {
8223e72e933SJed Brown           local_ghosts[i] = local_elems + owned_verts + i;
8233e72e933SJed Brown           ghosts[i].rank  = remote_rank;
8243e72e933SJed Brown           ghosts[i].index = remote_elem + remote_count;
8253e72e933SJed Brown           i++;
8263e72e933SJed Brown           if (i == num_ghosts) break;
8273e72e933SJed Brown           z = vert_z[owned_verts + i];
8283e72e933SJed Brown         }
8293e72e933SJed Brown         if (IjkActive(layout.vextent, loc)) remote_count++;
8306547ddbdSJed Brown         else rz += ZStepOct(rz);
8313e72e933SJed Brown       }
8323e72e933SJed Brown     }
8333e72e933SJed Brown     PetscCall(PetscSFSetGraph(sf, local_elems + local_verts, num_ghosts, local_ghosts, PETSC_OWN_POINTER, ghosts, PETSC_OWN_POINTER));
8343e72e933SJed Brown     PetscCall(PetscObjectSetName((PetscObject)sf, "SFC Point SF"));
8353e72e933SJed Brown     PetscCall(DMSetPointSF(dm, sf));
8363e72e933SJed Brown     PetscCall(PetscSFDestroy(&sf));
8373e72e933SJed Brown   }
8383e72e933SJed Brown   {
8393e72e933SJed Brown     Vec          coordinates;
8403e72e933SJed Brown     PetscScalar *coords;
8413e72e933SJed Brown     PetscSection coord_section;
8423e72e933SJed Brown     PetscInt     coord_size;
8433e72e933SJed Brown     PetscCall(DMGetCoordinateSection(dm, &coord_section));
8443e72e933SJed Brown     PetscCall(PetscSectionSetNumFields(coord_section, 1));
8453e72e933SJed Brown     PetscCall(PetscSectionSetFieldComponents(coord_section, 0, dim));
8463e72e933SJed Brown     PetscCall(PetscSectionSetChart(coord_section, local_elems, local_elems + local_verts));
8473e72e933SJed Brown     for (PetscInt v = 0; v < local_verts; v++) {
8483e72e933SJed Brown       PetscInt point = local_elems + v;
8493e72e933SJed Brown       PetscCall(PetscSectionSetDof(coord_section, point, dim));
8503e72e933SJed Brown       PetscCall(PetscSectionSetFieldDof(coord_section, point, 0, dim));
8513e72e933SJed Brown     }
8523e72e933SJed Brown     PetscCall(PetscSectionSetUp(coord_section));
8533e72e933SJed Brown     PetscCall(PetscSectionGetStorageSize(coord_section, &coord_size));
8543e72e933SJed Brown     PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates));
8553e72e933SJed Brown     PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
8563e72e933SJed Brown     PetscCall(VecSetSizes(coordinates, coord_size, PETSC_DETERMINE));
8573e72e933SJed Brown     PetscCall(VecSetBlockSize(coordinates, dim));
8583e72e933SJed Brown     PetscCall(VecSetType(coordinates, VECSTANDARD));
8593e72e933SJed Brown     PetscCall(VecGetArray(coordinates, &coords));
8603e72e933SJed Brown     for (PetscInt v = 0; v < local_verts; v++) {
8613e72e933SJed Brown       Ijk loc             = ZCodeSplit(vert_z[v]);
8623e72e933SJed Brown       coords[v * dim + 0] = lower[0] + loc.i * (upper[0] - lower[0]) / layout.eextent.i;
8633e72e933SJed Brown       if (dim > 1) coords[v * dim + 1] = lower[1] + loc.j * (upper[1] - lower[1]) / layout.eextent.j;
8643e72e933SJed Brown       if (dim > 2) coords[v * dim + 2] = lower[2] + loc.k * (upper[2] - lower[2]) / layout.eextent.k;
8653e72e933SJed Brown     }
8663e72e933SJed Brown     PetscCall(VecRestoreArray(coordinates, &coords));
8673e72e933SJed Brown     PetscCall(DMSetCoordinatesLocal(dm, coordinates));
8683e72e933SJed Brown     PetscCall(VecDestroy(&coordinates));
8693e72e933SJed Brown   }
8703e72e933SJed Brown   if (interpolate) {
8713431e603SJed Brown     PetscCall(DMPlexInterpolateInPlace_Internal(dm));
8725f06a3ddSJed Brown     // It's currently necessary to orient the donor and periodic edges consistently. An easy way to ensure that is ot
8735f06a3ddSJed Brown     // give all edges positive orientation. Since vertices are created in Z-order, all ranks will agree about the
8745f06a3ddSJed Brown     // ordering cone[0] < cone[1]. This is overkill and it would be nice to remove this preparation and make
8755f06a3ddSJed Brown     // DMPlexCreateIsoperiodicClosureSF_Private() more resilient, so it fixes any inconsistent orientations. That might
8765f06a3ddSJed Brown     // be needed in a general CGNS reader, for example.
8775f06a3ddSJed Brown     PetscCall(DMPlexOrientPositiveEdges_Private(dm));
8783431e603SJed Brown 
8793431e603SJed Brown     DMLabel label;
8803431e603SJed Brown     PetscCall(DMCreateLabel(dm, "Face Sets"));
8813431e603SJed Brown     PetscCall(DMGetLabel(dm, "Face Sets", &label));
8829ae3492bSJames Wright     PetscSegBuffer per_faces[3], donor_face_closure[3], my_donor_faces[3];
8839ae3492bSJames Wright     for (PetscInt i = 0; i < 3; i++) {
8849ae3492bSJames Wright       PetscCall(PetscSegBufferCreate(sizeof(PetscInt), 64, &per_faces[i]));
8859ae3492bSJames Wright       PetscCall(PetscSegBufferCreate(sizeof(PetscInt), 64, &my_donor_faces[i]));
8869ae3492bSJames Wright       PetscCall(PetscSegBufferCreate(sizeof(ZCode), 64 * PetscPowInt(2, dim), &donor_face_closure[i]));
8879ae3492bSJames Wright     }
8883431e603SJed Brown     PetscInt fStart, fEnd, vStart, vEnd;
8893431e603SJed Brown     PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
8903431e603SJed Brown     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
8913431e603SJed Brown     for (PetscInt f = fStart; f < fEnd; f++) {
8924e2e9504SJed Brown       PetscInt npoints, *points = NULL, num_fverts = 0, fverts[8];
8933431e603SJed Brown       PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &npoints, &points));
8943431e603SJed Brown       PetscInt bc_count[6] = {0};
8953431e603SJed Brown       for (PetscInt i = 0; i < npoints; i++) {
8963431e603SJed Brown         PetscInt p = points[2 * i];
8973431e603SJed Brown         if (p < vStart || vEnd <= p) continue;
8984e2e9504SJed Brown         fverts[num_fverts++] = p;
8993431e603SJed Brown         Ijk loc              = ZCodeSplit(vert_z[p - vStart]);
9003431e603SJed Brown         // Convention here matches DMPlexCreateCubeMesh_Internal
9013431e603SJed Brown         bc_count[0] += loc.i == 0;
9023431e603SJed Brown         bc_count[1] += loc.i == layout.vextent.i - 1;
9033431e603SJed Brown         bc_count[2] += loc.j == 0;
9043431e603SJed Brown         bc_count[3] += loc.j == layout.vextent.j - 1;
9053431e603SJed Brown         bc_count[4] += loc.k == 0;
9063431e603SJed Brown         bc_count[5] += loc.k == layout.vextent.k - 1;
9073e72e933SJed Brown       }
9083431e603SJed Brown       PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &npoints, &points));
9093431e603SJed Brown       for (PetscInt bc = 0, bc_match = 0; bc < 2 * dim; bc++) {
9103431e603SJed Brown         if (bc_count[bc] == PetscPowInt(2, dim - 1)) {
9116725e60dSJed Brown           PetscCall(DMPlexOrientPoint(dm, f, face_orient_dim[dim][bc]));
9124e2e9504SJed Brown           if (periodicity[bc / 2] == DM_BOUNDARY_PERIODIC) {
9134e2e9504SJed Brown             PetscInt *put;
9144e2e9504SJed Brown             if (bc % 2 == 0) { // donor face; no label
9159ae3492bSJames Wright               PetscCall(PetscSegBufferGet(my_donor_faces[bc / 2], 1, &put));
9164e2e9504SJed Brown               *put = f;
9174e2e9504SJed Brown             } else { // periodic face
9189ae3492bSJames Wright               PetscCall(PetscSegBufferGet(per_faces[bc / 2], 1, &put));
9194e2e9504SJed Brown               *put = f;
9204e2e9504SJed Brown               ZCode *zput;
9219ae3492bSJames Wright               PetscCall(PetscSegBufferGet(donor_face_closure[bc / 2], num_fverts, &zput));
9224e2e9504SJed Brown               for (PetscInt i = 0; i < num_fverts; i++) {
9234e2e9504SJed Brown                 Ijk loc = ZCodeSplit(vert_z[fverts[i] - vStart]);
9244e2e9504SJed Brown                 switch (bc / 2) {
9254e2e9504SJed Brown                 case 0:
9264e2e9504SJed Brown                   loc.i = 0;
9274e2e9504SJed Brown                   break;
9284e2e9504SJed Brown                 case 1:
9294e2e9504SJed Brown                   loc.j = 0;
9304e2e9504SJed Brown                   break;
9314e2e9504SJed Brown                 case 2:
9324e2e9504SJed Brown                   loc.k = 0;
9334e2e9504SJed Brown                   break;
9344e2e9504SJed Brown                 }
9354e2e9504SJed Brown                 *zput++ = ZEncode(loc);
9364e2e9504SJed Brown               }
9374e2e9504SJed Brown             }
9384e2e9504SJed Brown             continue;
9394e2e9504SJed Brown           }
9403431e603SJed Brown           PetscAssert(bc_match == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face matches multiple face sets");
9413431e603SJed Brown           PetscCall(DMLabelSetValue(label, f, face_marker_dim[dim][bc]));
9423431e603SJed Brown           bc_match++;
9433431e603SJed Brown         }
9443431e603SJed Brown       }
9453431e603SJed Brown     }
9463431e603SJed Brown     // Ensure that the Coordinate DM has our new boundary labels
9473431e603SJed Brown     DM cdm;
9483431e603SJed Brown     PetscCall(DMGetCoordinateDM(dm, &cdm));
9493431e603SJed Brown     PetscCall(DMCopyLabels(dm, cdm, PETSC_COPY_VALUES, PETSC_FALSE, DM_COPY_LABELS_FAIL));
9506725e60dSJed Brown     if (periodicity[0] == DM_BOUNDARY_PERIODIC || (dim > 1 && periodicity[1] == DM_BOUNDARY_PERIODIC) || (dim > 2 && periodicity[2] == DM_BOUNDARY_PERIODIC)) {
9515f06a3ddSJed Brown       PetscCall(DMPlexCreateBoxMesh_Tensor_SFC_Periodicity_Private(dm, &layout, vert_z, per_faces, lower, upper, periodicity, donor_face_closure, my_donor_faces));
9526725e60dSJed Brown     }
9539ae3492bSJames Wright     for (PetscInt i = 0; i < 3; i++) {
9549ae3492bSJames Wright       PetscCall(PetscSegBufferDestroy(&per_faces[i]));
9559ae3492bSJames Wright       PetscCall(PetscSegBufferDestroy(&donor_face_closure[i]));
9569ae3492bSJames Wright       PetscCall(PetscSegBufferDestroy(&my_donor_faces[i]));
9579ae3492bSJames Wright     }
9583431e603SJed Brown   }
9594e2e9504SJed Brown   PetscCall(PetscFree(layout.zstarts));
9603431e603SJed Brown   PetscCall(PetscFree(vert_z));
961*ea7b3863SJames Wright   PetscCall(PetscLogEventEnd(DMPLEX_CreateBoxSFC, dm, 0, 0, 0));
9623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9633e72e933SJed Brown }
9644e2e9504SJed Brown 
9654e2e9504SJed Brown /*@
9665f06a3ddSJed Brown   DMPlexSetIsoperiodicFaceSF - Express periodicity from an existing mesh
9674e2e9504SJed Brown 
96820f4b53cSBarry Smith   Logically Collective
9694e2e9504SJed Brown 
9704e2e9504SJed Brown   Input Parameters:
9714e2e9504SJed Brown + dm           - The `DMPLEX` on which to set periodicity
9721fca310dSJames Wright . num_face_sfs - Number of `PetscSF`s in `face_sfs`
9731fca310dSJames Wright - face_sfs     - Array of `PetscSF` in which roots are (owned) donor faces and leaves are faces that must be matched to a (possibly remote) donor face.
9744e2e9504SJed Brown 
9754e2e9504SJed Brown   Level: advanced
9764e2e9504SJed Brown 
97753c0d4aeSBarry Smith   Note:
9785dca41c3SJed Brown   One can use `-dm_plex_shape zbox` to use this mode of periodicity, wherein the periodic points are distinct both globally
9794e2e9504SJed Brown   and locally, but are paired when creating a global dof space.
9804e2e9504SJed Brown 
9811cc06b55SBarry Smith .seealso: [](ch_unstructured), `DMPLEX`, `DMGetGlobalSection()`, `DMPlexGetIsoperiodicFaceSF()`
9824e2e9504SJed Brown @*/
9831fca310dSJames Wright PetscErrorCode DMPlexSetIsoperiodicFaceSF(DM dm, PetscInt num_face_sfs, PetscSF *face_sfs)
9844e2e9504SJed Brown {
9854e2e9504SJed Brown   DM_Plex *plex = (DM_Plex *)dm->data;
9864d86920dSPierre Jolivet 
9874e2e9504SJed Brown   PetscFunctionBegin;
9884e2e9504SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9891fca310dSJames Wright   if (face_sfs) PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMGetIsoperiodicPointSF_C", DMGetIsoperiodicPointSF_Plex));
9901fca310dSJames Wright   if (face_sfs == plex->periodic.face_sfs && num_face_sfs == plex->periodic.num_face_sfs) PetscFunctionReturn(PETSC_SUCCESS);
9911fca310dSJames Wright 
9921fca310dSJames Wright   for (PetscInt i = 0; i < num_face_sfs; i++) PetscCall(PetscObjectReference((PetscObject)face_sfs[i]));
9931fca310dSJames Wright 
9941fca310dSJames Wright   if (plex->periodic.num_face_sfs > 0) {
9951fca310dSJames Wright     for (PetscInt i = 0; i < plex->periodic.num_face_sfs; i++) PetscCall(PetscSFDestroy(&plex->periodic.face_sfs[i]));
9961fca310dSJames Wright     PetscCall(PetscFree(plex->periodic.face_sfs));
9971fca310dSJames Wright   }
9981fca310dSJames Wright 
9991fca310dSJames Wright   plex->periodic.num_face_sfs = num_face_sfs;
10001fca310dSJames Wright   PetscCall(PetscCalloc1(num_face_sfs, &plex->periodic.face_sfs));
10011fca310dSJames Wright   for (PetscInt i = 0; i < num_face_sfs; i++) plex->periodic.face_sfs[i] = face_sfs[i];
10025f06a3ddSJed Brown 
10035f06a3ddSJed Brown   DM cdm = dm->coordinates[0].dm; // Can't DMGetCoordinateDM because it automatically creates one
10045f06a3ddSJed Brown   if (cdm) {
10051fca310dSJames Wright     PetscCall(DMPlexSetIsoperiodicFaceSF(cdm, num_face_sfs, face_sfs));
10061fca310dSJames Wright     if (face_sfs) cdm->periodic.setup = DMPeriodicCoordinateSetUp_Internal;
10075f06a3ddSJed Brown   }
10083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10094e2e9504SJed Brown }
10104e2e9504SJed Brown 
10111fca310dSJames Wright /*@C
10125f06a3ddSJed Brown   DMPlexGetIsoperiodicFaceSF - Obtain periodicity for a mesh
10134e2e9504SJed Brown 
101420f4b53cSBarry Smith   Logically Collective
10154e2e9504SJed Brown 
101620f4b53cSBarry Smith   Input Parameter:
10174e2e9504SJed Brown . dm - The `DMPLEX` for which to obtain periodic relation
10184e2e9504SJed Brown 
10199cde84edSJose E. Roman   Output Parameters:
10201fca310dSJames Wright + num_face_sfs - Number of `PetscSF`s in the array
10211fca310dSJames Wright - face_sfs     - Array of `PetscSF` in which roots are (owned) donor faces and leaves are faces that must be matched to a (possibly remote) donor face.
10224e2e9504SJed Brown 
10234e2e9504SJed Brown   Level: advanced
10244e2e9504SJed Brown 
10251cc06b55SBarry Smith .seealso: [](ch_unstructured), `DMPLEX`, `DMGetGlobalSection()`, `DMPlexSetIsoperiodicFaceSF()`
10264e2e9504SJed Brown @*/
10271fca310dSJames Wright PetscErrorCode DMPlexGetIsoperiodicFaceSF(DM dm, PetscInt *num_face_sfs, const PetscSF **face_sfs)
10284e2e9504SJed Brown {
10294e2e9504SJed Brown   DM_Plex *plex = (DM_Plex *)dm->data;
10304d86920dSPierre Jolivet 
10314e2e9504SJed Brown   PetscFunctionBegin;
10324e2e9504SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10331fca310dSJames Wright   *face_sfs     = plex->periodic.face_sfs;
10341fca310dSJames Wright   *num_face_sfs = plex->periodic.num_face_sfs;
10353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10364e2e9504SJed Brown }
10376725e60dSJed Brown 
10385f06a3ddSJed Brown /*@C
10395f06a3ddSJed Brown   DMPlexSetIsoperiodicFaceTransform - set geometric transform from donor to periodic points
10406725e60dSJed Brown 
10416725e60dSJed Brown   Logically Collective
10426725e60dSJed Brown 
104360225df5SJacob Faibussowitsch   Input Parameters:
104453c0d4aeSBarry Smith + dm - `DMPLEX` that has been configured with `DMPlexSetIsoperiodicFaceSF()`
10451fca310dSJames Wright . n  - Number of transforms in array
10461fca310dSJames Wright - t  - Array of 4x4 affine transformation basis.
10476725e60dSJed Brown 
104853c0d4aeSBarry Smith   Level: advanced
104953c0d4aeSBarry Smith 
10505f06a3ddSJed Brown   Notes:
10515f06a3ddSJed Brown   Affine transforms are 4x4 matrices in which the leading 3x3 block expresses a rotation (or identity for no rotation),
10525f06a3ddSJed Brown   the last column contains a translation vector, and the bottom row is all zero except the last entry, which must always
10535f06a3ddSJed Brown   be 1. This representation is common in geometric modeling and allows affine transformations to be composed using
1054aaa8cc7dSPierre Jolivet   simple matrix multiplication.
10555f06a3ddSJed Brown 
10565f06a3ddSJed Brown   Although the interface accepts a general affine transform, only affine translation is supported at present.
10575f06a3ddSJed Brown 
105860225df5SJacob Faibussowitsch   Developer Notes:
10595f06a3ddSJed Brown   This interface should be replaced by making BasisTransform public, expanding it to support affine representations, and
10605f06a3ddSJed Brown   adding GPU implementations to apply the G2L/L2G transforms.
106153c0d4aeSBarry Smith 
10621cc06b55SBarry Smith .seealso: [](ch_unstructured), `DMPLEX`, `DMGetGlobalSection()`, `DMPlexSetIsoperiodicFaceSF()`
10636725e60dSJed Brown @*/
1064cc4c1da9SBarry Smith PetscErrorCode DMPlexSetIsoperiodicFaceTransform(DM dm, PetscInt n, const PetscScalar t[])
10656725e60dSJed Brown {
10666725e60dSJed Brown   DM_Plex *plex = (DM_Plex *)dm->data;
10674d86920dSPierre Jolivet 
10686725e60dSJed Brown   PetscFunctionBegin;
10696725e60dSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10701fca310dSJames Wright   PetscCheck(n == plex->periodic.num_face_sfs, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Number of transforms (%" PetscInt_FMT ") must equal number of isoperiodc face SFs (%" PetscInt_FMT ")", n, plex->periodic.num_face_sfs);
10711fca310dSJames Wright 
10721fca310dSJames Wright   PetscCall(PetscMalloc1(n, &plex->periodic.transform));
10731fca310dSJames Wright   for (PetscInt i = 0; i < n; i++) {
10746725e60dSJed Brown     for (PetscInt j = 0; j < 4; j++) {
10751fca310dSJames Wright       for (PetscInt k = 0; k < 4; k++) {
10761fca310dSJames Wright         PetscCheck(j != k || t[i * 16 + j * 4 + k] == 1., PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Rotated transforms not supported");
10771fca310dSJames Wright         plex->periodic.transform[i][j][k] = t[i * 16 + j * 4 + k];
10781fca310dSJames Wright       }
10796725e60dSJed Brown     }
10806725e60dSJed Brown   }
10813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10826725e60dSJed Brown }
1083