xref: /petsc/src/dm/impls/plex/plexsfc.c (revision bc3e2e66c7b74b02889f07219b9c373dff6e4370)
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 
20376e073fSJames Wright // ***** Overview of ZCode *******
21376e073fSJames Wright // The SFC uses integer indexing for each dimension and encodes them into a single integer by interleaving the bits of each index.
22376e073fSJames Wright // This is known as Morton encoding, and is refered to as ZCode in this code.
23376e073fSJames Wright // So for index i having bits [i2,i1,i0], and similar for indexes j and k, the ZCode (Morton number) would be:
24376e073fSJames Wright //    [k2,j2,i2,k1,j1,i1,k0,j0,i0]
25376e073fSJames Wright // This encoding allows for easier traversal of the SFC structure (see https://en.wikipedia.org/wiki/Z-order_curve and `ZStepOct()`).
26376e073fSJames Wright // `ZEncode()` is used to go from indices to ZCode, while `ZCodeSplit()` goes from ZCode back to indices.
27376e073fSJames Wright 
28376e073fSJames Wright // Decodes the leading interleaved index from a ZCode
29376e073fSJames Wright // e.g. [k2,j2,i2,k1,j1,i1,k0,j0,i0] -> [i2,i1,i0]
30*bc3e2e66SJames Wright // Magic numbers taken from https://stackoverflow.com/a/18528775/7564988 (translated to octal)
313e72e933SJed Brown static unsigned ZCodeSplit1(ZCode z)
323e72e933SJed Brown {
33*bc3e2e66SJames Wright   z &= 0111111111111111111111;
34*bc3e2e66SJames Wright   z = (z | z >> 2) & 0103030303030303030303;
35*bc3e2e66SJames Wright   z = (z | z >> 4) & 0100170017001700170017;
36*bc3e2e66SJames Wright   z = (z | z >> 8) & 0370000037700000377;
37*bc3e2e66SJames Wright   z = (z | z >> 16) & 0370000000000177777;
38*bc3e2e66SJames Wright   z = (z | z >> 32) & 07777777;
393e72e933SJed Brown   return (unsigned)z;
403e72e933SJed Brown }
413e72e933SJed Brown 
42376e073fSJames Wright // Encodes the leading interleaved index from a ZCode
43376e073fSJames Wright // e.g. [i2,i1,i0] -> [0,0,i2,0,0,i1,0,0,i0]
443e72e933SJed Brown static ZCode ZEncode1(unsigned t)
453e72e933SJed Brown {
463e72e933SJed Brown   ZCode z = t;
47*bc3e2e66SJames Wright   z &= 07777777;
48*bc3e2e66SJames Wright   z = (z | z << 32) & 0370000000000177777;
49*bc3e2e66SJames Wright   z = (z | z << 16) & 0370000037700000377;
50*bc3e2e66SJames Wright   z = (z | z << 8) & 0100170017001700170017;
51*bc3e2e66SJames Wright   z = (z | z << 4) & 0103030303030303030303;
52*bc3e2e66SJames Wright   z = (z | z << 2) & 0111111111111111111111;
533e72e933SJed Brown   return z;
543e72e933SJed Brown }
553e72e933SJed Brown 
56376e073fSJames Wright // Decodes i j k indices from a ZCode.
57376e073fSJames Wright // Uses `ZCodeSplit1()` by shifting ZCode so that the leading index is the desired one to decode
583e72e933SJed Brown static Ijk ZCodeSplit(ZCode z)
593e72e933SJed Brown {
603e72e933SJed Brown   Ijk c;
613e72e933SJed Brown   c.i = ZCodeSplit1(z >> 2);
623e72e933SJed Brown   c.j = ZCodeSplit1(z >> 1);
633e72e933SJed Brown   c.k = ZCodeSplit1(z >> 0);
643e72e933SJed Brown   return c;
653e72e933SJed Brown }
663e72e933SJed Brown 
67376e073fSJames Wright // Encodes i j k indices to a ZCode.
68376e073fSJames Wright // Uses `ZCodeEncode1()` by shifting resulting ZCode to the appropriate bit placement
693e72e933SJed Brown static ZCode ZEncode(Ijk c)
703e72e933SJed Brown {
716497c311SBarry Smith   ZCode z = (ZEncode1((unsigned int)c.i) << 2) | (ZEncode1((unsigned int)c.j) << 1) | ZEncode1((unsigned int)c.k);
723e72e933SJed Brown   return z;
733e72e933SJed Brown }
743e72e933SJed Brown 
753e72e933SJed Brown static PetscBool IjkActive(Ijk extent, Ijk l)
763e72e933SJed Brown {
773e72e933SJed Brown   if (l.i < extent.i && l.j < extent.j && l.k < extent.k) return PETSC_TRUE;
783e72e933SJed Brown   return PETSC_FALSE;
793e72e933SJed Brown }
803e72e933SJed Brown 
816547ddbdSJed Brown // If z is not the base of an octet (last 3 bits 0), return 0.
826547ddbdSJed Brown //
836547ddbdSJed Brown // If z is the base of an octet, we recursively grow to the biggest structured octet. This is typically useful when a z
846547ddbdSJed Brown // is outside the domain and we wish to skip a (possibly recursively large) octet to find our next interesting point.
856547ddbdSJed Brown static ZCode ZStepOct(ZCode z)
866547ddbdSJed Brown {
876547ddbdSJed Brown   if (PetscUnlikely(z == 0)) return 0; // Infinite loop below if z == 0
886547ddbdSJed Brown   ZCode step = 07;
896547ddbdSJed Brown   for (; (z & step) == 0; step = (step << 3) | 07) { }
906547ddbdSJed Brown   return step >> 3;
916547ddbdSJed Brown }
926547ddbdSJed Brown 
933e72e933SJed Brown // Since element/vertex box extents are typically not equal powers of 2, Z codes that lie within the domain are not contiguous.
945f06a3ddSJed Brown static PetscErrorCode ZLayoutCreate(PetscMPIInt size, const PetscInt eextent[3], const PetscInt vextent[3], ZLayout *layout)
953e72e933SJed Brown {
96c77877e3SJed Brown   PetscFunctionBegin;
975f06a3ddSJed Brown   layout->eextent.i = eextent[0];
985f06a3ddSJed Brown   layout->eextent.j = eextent[1];
995f06a3ddSJed Brown   layout->eextent.k = eextent[2];
1005f06a3ddSJed Brown   layout->vextent.i = vextent[0];
1015f06a3ddSJed Brown   layout->vextent.j = vextent[1];
1025f06a3ddSJed Brown   layout->vextent.k = vextent[2];
1035f06a3ddSJed Brown   layout->comm_size = size;
1045f06a3ddSJed Brown   layout->zstarts   = NULL;
1055f06a3ddSJed Brown   PetscCall(PetscMalloc1(size + 1, &layout->zstarts));
1063e72e933SJed Brown 
1073e72e933SJed Brown   PetscInt total_elems = eextent[0] * eextent[1] * eextent[2];
1083e72e933SJed Brown   ZCode    z           = 0;
1095f06a3ddSJed Brown   layout->zstarts[0]   = 0;
1106547ddbdSJed Brown   // This loop traverses all vertices in the global domain, so is worth making fast. We use ZStepBound
1113e72e933SJed Brown   for (PetscMPIInt r = 0; r < size; r++) {
1123e72e933SJed Brown     PetscInt elems_needed = (total_elems / size) + (total_elems % size > r), count;
1133e72e933SJed Brown     for (count = 0; count < elems_needed; z++) {
1146547ddbdSJed Brown       ZCode skip = ZStepOct(z); // optimistically attempt a longer step
1156547ddbdSJed Brown       for (ZCode s = skip;; s >>= 3) {
1166547ddbdSJed Brown         Ijk trial = ZCodeSplit(z + s);
1176547ddbdSJed Brown         if (IjkActive(layout->eextent, trial)) {
1186547ddbdSJed Brown           while (count + s + 1 > (ZCode)elems_needed) s >>= 3; // Shrink the octet
1196547ddbdSJed Brown           count += s + 1;
1206547ddbdSJed Brown           z += s;
1216547ddbdSJed Brown           break;
1226547ddbdSJed Brown         }
1236547ddbdSJed Brown         if (s == 0) { // the whole skip octet is inactive
1246547ddbdSJed Brown           z += skip;
1256547ddbdSJed Brown           break;
1266547ddbdSJed Brown         }
1276547ddbdSJed Brown       }
1283e72e933SJed Brown     }
1293e72e933SJed Brown     // Pick up any extra vertices in the Z ordering before the next rank's first owned element.
130c77877e3SJed Brown     //
1315f06a3ddSJed Brown     // This leads to poorly balanced vertices when eextent is a power of 2, since all the fringe vertices end up
132c77877e3SJed Brown     // on the last rank. A possible solution is to balance the Z-order vertices independently from the cells, which will
133c77877e3SJed Brown     // result in a lot of element closures being remote. We could finish marking boundary conditions, then do a round of
134c77877e3SJed Brown     // vertex ownership smoothing (which would reorder and redistribute vertices without touching element distribution).
135c77877e3SJed Brown     // Another would be to have an analytic ownership criteria for vertices in the fringe veextent - eextent. This would
136c77877e3SJed Brown     // complicate the job of identifying an owner and its offset.
1375f06a3ddSJed Brown     //
1385f06a3ddSJed Brown     // The current recommended approach is to let `-dm_distribute 1` (default) resolve vertex ownership. This is
1395f06a3ddSJed Brown     // *mandatory* with isoperiodicity (except in special cases) to remove standed vertices from local spaces. Here's
1405f06a3ddSJed Brown     // the issue:
1415f06a3ddSJed Brown     //
1425f06a3ddSJed Brown     // Consider this partition on rank 0 (left) and rank 1.
1435f06a3ddSJed Brown     //
1445f06a3ddSJed Brown     //    4 --------  5 -- 14 --10 -- 21 --11
1455f06a3ddSJed Brown     //                |          |          |
1465f06a3ddSJed Brown     // 7 -- 16 --  8  |          |          |
1475f06a3ddSJed Brown     // |           |  3 -------  7 -------  9
1485f06a3ddSJed Brown     // |           |             |          |
1495f06a3ddSJed Brown     // 4 --------  6 ------ 10   |          |
1505f06a3ddSJed Brown     // |           |         |   6 -- 16 -- 8
1515f06a3ddSJed Brown     // |           |         |
1525f06a3ddSJed Brown     // 3 ---11---  5 --18--  9
1535f06a3ddSJed Brown     //
1545f06a3ddSJed Brown     // The periodic face SF looks like
1555f06a3ddSJed Brown     // [0] Number of roots=21, leaves=1, remote ranks=1
1565f06a3ddSJed Brown     // [0] 16 <- (0,11)
1575f06a3ddSJed Brown     // [1] Number of roots=22, leaves=2, remote ranks=2
1585f06a3ddSJed Brown     // [1] 14 <- (0,18)
1595f06a3ddSJed Brown     // [1] 21 <- (1,16)
1605f06a3ddSJed Brown     //
1615f06a3ddSJed 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
1625f06a3ddSJed Brown     // the point SF links to (1,4) and (1,5). Rank 1 learns about the periodic mapping of (1,5) while handling face
1635f06a3ddSJed Brown     // (1,14), but never learns that vertex (1,4) has been mapped to (0,3) by face (0,16).
1645f06a3ddSJed Brown     //
1655f06a3ddSJed Brown     // We can relatively easily inform vertex (1,4) of this mapping, but it stays in rank 1's local space despite not
1665f06a3ddSJed Brown     // being in the closure and thus not being contributed to. This would be mostly harmless except that some viewer
1675f06a3ddSJed Brown     // routines expect all local points to be somehow significant. It is not easy to analytically remove the (1,4)
1685f06a3ddSJed Brown     // vertex because the point SF and isoperiodic face SF would need to be updated to account for removal of the
1695f06a3ddSJed Brown     // stranded vertices.
1705f06a3ddSJed Brown     for (; z <= ZEncode(layout->vextent); z++) {
1713e72e933SJed Brown       Ijk loc = ZCodeSplit(z);
1725f06a3ddSJed Brown       if (IjkActive(layout->eextent, loc)) break;
1736547ddbdSJed Brown       z += ZStepOct(z);
1743e72e933SJed Brown     }
1755f06a3ddSJed Brown     layout->zstarts[r + 1] = z;
1763e72e933SJed Brown   }
1776547ddbdSJed Brown   layout->zstarts[size] = ZEncode(layout->vextent);
1783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1793e72e933SJed Brown }
1803e72e933SJed Brown 
1814e2e9504SJed Brown static PetscInt ZLayoutElementsOnRank(const ZLayout *layout, PetscMPIInt rank)
1824e2e9504SJed Brown {
1834e2e9504SJed Brown   PetscInt remote_elem = 0;
1844e2e9504SJed Brown   for (ZCode rz = layout->zstarts[rank]; rz < layout->zstarts[rank + 1]; rz++) {
1854e2e9504SJed Brown     Ijk loc = ZCodeSplit(rz);
1864e2e9504SJed Brown     if (IjkActive(layout->eextent, loc)) remote_elem++;
1876547ddbdSJed Brown     else rz += ZStepOct(rz);
1884e2e9504SJed Brown   }
1894e2e9504SJed Brown   return remote_elem;
1904e2e9504SJed Brown }
1914e2e9504SJed Brown 
19266976f2fSJacob Faibussowitsch static PetscInt ZCodeFind(ZCode key, PetscInt n, const ZCode X[])
1933e72e933SJed Brown {
1943e72e933SJed Brown   PetscInt lo = 0, hi = n;
1953e72e933SJed Brown 
1963e72e933SJed Brown   if (n == 0) return -1;
1973e72e933SJed Brown   while (hi - lo > 1) {
1983e72e933SJed Brown     PetscInt mid = lo + (hi - lo) / 2;
1993e72e933SJed Brown     if (key < X[mid]) hi = mid;
2003e72e933SJed Brown     else lo = mid;
2013e72e933SJed Brown   }
2023e72e933SJed Brown   return key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
2033e72e933SJed Brown }
2043e72e933SJed Brown 
2059ae3492bSJames 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])
2064e2e9504SJed Brown {
2074e2e9504SJed Brown   MPI_Comm    comm;
2089ae3492bSJames Wright   PetscInt    dim, vStart, vEnd;
2094e2e9504SJed Brown   PetscMPIInt size;
2109ae3492bSJames Wright   PetscSF     face_sfs[3];
2119ae3492bSJames Wright   PetscScalar transforms[3][4][4] = {{{0}}};
2124e2e9504SJed Brown 
2134e2e9504SJed Brown   PetscFunctionBegin;
2144e2e9504SJed Brown   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2154e2e9504SJed Brown   PetscCallMPI(MPI_Comm_size(comm, &size));
2164e2e9504SJed Brown   PetscCall(DMGetDimension(dm, &dim));
2174e2e9504SJed Brown   const PetscInt csize = PetscPowInt(2, dim - 1);
2184e2e9504SJed Brown   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
2199ae3492bSJames Wright 
2209ae3492bSJames Wright   PetscInt num_directions = 0;
2219ae3492bSJames Wright   for (PetscInt direction = 0; direction < dim; direction++) {
2226497c311SBarry Smith     PetscCount   num_faces;
2239ae3492bSJames Wright     PetscInt    *faces;
2249ae3492bSJames Wright     ZCode       *donor_verts, *donor_minz;
2259ae3492bSJames Wright     PetscSFNode *leaf;
2266497c311SBarry Smith     PetscCount   num_multiroots = 0;
2276497c311SBarry Smith     PetscInt     pStart, pEnd;
2286497c311SBarry Smith     PetscBool    sorted;
2296497c311SBarry Smith     PetscInt     inum_faces;
2309ae3492bSJames Wright 
2319ae3492bSJames Wright     if (periodicity[direction] != DM_BOUNDARY_PERIODIC) continue;
2329ae3492bSJames Wright     PetscCall(PetscSegBufferGetSize(per_faces[direction], &num_faces));
2339ae3492bSJames Wright     PetscCall(PetscSegBufferExtractInPlace(per_faces[direction], &faces));
2349ae3492bSJames Wright     PetscCall(PetscSegBufferExtractInPlace(donor_face_closure[direction], &donor_verts));
2354e2e9504SJed Brown     PetscCall(PetscMalloc1(num_faces, &donor_minz));
2364e2e9504SJed Brown     PetscCall(PetscMalloc1(num_faces, &leaf));
2376497c311SBarry Smith     for (PetscCount i = 0; i < num_faces; i++) {
2384e2e9504SJed Brown       ZCode minz = donor_verts[i * csize];
2396497c311SBarry Smith 
2404e2e9504SJed Brown       for (PetscInt j = 1; j < csize; j++) minz = PetscMin(minz, donor_verts[i * csize + j]);
2414e2e9504SJed Brown       donor_minz[i] = minz;
2424e2e9504SJed Brown     }
2436497c311SBarry Smith     PetscCall(PetscIntCast(num_faces, &inum_faces));
2446497c311SBarry Smith     PetscCall(PetscSortedInt64(inum_faces, (const PetscInt64 *)donor_minz, &sorted));
2459ae3492bSJames Wright     // If a donor vertex were chosen to broker multiple faces, we would have a logic error.
2469ae3492bSJames Wright     // Checking for sorting is a cheap check that there are no duplicates.
2479ae3492bSJames Wright     PetscCheck(sorted, PETSC_COMM_SELF, PETSC_ERR_PLIB, "minz not sorted; possible duplicates not checked");
2486497c311SBarry Smith     for (PetscCount i = 0; i < num_faces;) {
2494e2e9504SJed Brown       ZCode       z           = donor_minz[i];
2506497c311SBarry Smith       PetscMPIInt remote_rank = (PetscMPIInt)ZCodeFind(z, size + 1, layout->zstarts), remote_count = 0;
2516497c311SBarry Smith 
2524e2e9504SJed Brown       if (remote_rank < 0) remote_rank = -(remote_rank + 1) - 1;
2534e2e9504SJed Brown       // Process all the vertices on this rank
2544e2e9504SJed Brown       for (ZCode rz = layout->zstarts[remote_rank]; rz < layout->zstarts[remote_rank + 1]; rz++) {
2554e2e9504SJed Brown         Ijk loc = ZCodeSplit(rz);
2566497c311SBarry Smith 
2574e2e9504SJed Brown         if (rz == z) {
2584e2e9504SJed Brown           leaf[i].rank  = remote_rank;
2594e2e9504SJed Brown           leaf[i].index = remote_count;
2604e2e9504SJed Brown           i++;
2616497c311SBarry Smith           if (i == num_faces) break;
2624e2e9504SJed Brown           z = donor_minz[i];
2634e2e9504SJed Brown         }
2644e2e9504SJed Brown         if (IjkActive(layout->vextent, loc)) remote_count++;
2654e2e9504SJed Brown       }
2664e2e9504SJed Brown     }
2674e2e9504SJed Brown     PetscCall(PetscFree(donor_minz));
2689ae3492bSJames Wright     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), &face_sfs[num_directions]));
2696497c311SBarry Smith     PetscCall(PetscSFSetGraph(face_sfs[num_directions], vEnd - vStart, inum_faces, NULL, PETSC_USE_POINTER, leaf, PETSC_USE_POINTER));
2704e2e9504SJed Brown     const PetscInt *my_donor_degree;
2719ae3492bSJames Wright     PetscCall(PetscSFComputeDegreeBegin(face_sfs[num_directions], &my_donor_degree));
2729ae3492bSJames Wright     PetscCall(PetscSFComputeDegreeEnd(face_sfs[num_directions], &my_donor_degree));
2736497c311SBarry Smith 
2744e2e9504SJed Brown     for (PetscInt i = 0; i < vEnd - vStart; i++) {
2754e2e9504SJed Brown       num_multiroots += my_donor_degree[i];
2764e2e9504SJed Brown       if (my_donor_degree[i] == 0) continue;
2775f06a3ddSJed Brown       PetscAssert(my_donor_degree[i] == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Local vertex has multiple faces");
2784e2e9504SJed Brown     }
2794e2e9504SJed Brown     PetscInt  *my_donors, *donor_indices, *my_donor_indices;
2806497c311SBarry Smith     PetscCount num_my_donors;
2816497c311SBarry Smith 
2829ae3492bSJames Wright     PetscCall(PetscSegBufferGetSize(my_donor_faces[direction], &num_my_donors));
283563af49cSJames 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);
2849ae3492bSJames Wright     PetscCall(PetscSegBufferExtractInPlace(my_donor_faces[direction], &my_donors));
2854e2e9504SJed Brown     PetscCall(PetscMalloc1(vEnd - vStart, &my_donor_indices));
2866497c311SBarry Smith     for (PetscCount i = 0; i < num_my_donors; i++) {
2874e2e9504SJed Brown       PetscInt f = my_donors[i];
2881690c2aeSBarry Smith       PetscInt num_points, *points = NULL, minv = PETSC_INT_MAX;
2896497c311SBarry Smith 
2904e2e9504SJed Brown       PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &num_points, &points));
2914e2e9504SJed Brown       for (PetscInt j = 0; j < num_points; j++) {
2924e2e9504SJed Brown         PetscInt p = points[2 * j];
2934e2e9504SJed Brown         if (p < vStart || vEnd <= p) continue;
2944e2e9504SJed Brown         minv = PetscMin(minv, p);
2954e2e9504SJed Brown       }
2964e2e9504SJed Brown       PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &num_points, &points));
2975f06a3ddSJed Brown       PetscAssert(my_donor_degree[minv - vStart] == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Local vertex not requested");
2984e2e9504SJed Brown       my_donor_indices[minv - vStart] = f;
2994e2e9504SJed Brown     }
3004e2e9504SJed Brown     PetscCall(PetscMalloc1(num_faces, &donor_indices));
3019ae3492bSJames Wright     PetscCall(PetscSFBcastBegin(face_sfs[num_directions], MPIU_INT, my_donor_indices, donor_indices, MPI_REPLACE));
3029ae3492bSJames Wright     PetscCall(PetscSFBcastEnd(face_sfs[num_directions], MPIU_INT, my_donor_indices, donor_indices, MPI_REPLACE));
3034e2e9504SJed Brown     PetscCall(PetscFree(my_donor_indices));
3044e2e9504SJed Brown     // Modify our leafs so they point to donor faces instead of donor minz. Additionally, give them indices as faces.
3056497c311SBarry Smith     for (PetscCount i = 0; i < num_faces; i++) leaf[i].index = donor_indices[i];
3064e2e9504SJed Brown     PetscCall(PetscFree(donor_indices));
3074e2e9504SJed Brown     PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
3086497c311SBarry Smith     PetscCall(PetscSFSetGraph(face_sfs[num_directions], pEnd - pStart, inum_faces, faces, PETSC_COPY_VALUES, leaf, PETSC_OWN_POINTER));
3099ae3492bSJames Wright     {
3109ae3492bSJames Wright       char face_sf_name[PETSC_MAX_PATH_LEN];
3119ae3492bSJames Wright       PetscCall(PetscSNPrintf(face_sf_name, sizeof face_sf_name, "Z-order Isoperiodic Faces #%" PetscInt_FMT, num_directions));
3129ae3492bSJames Wright       PetscCall(PetscObjectSetName((PetscObject)face_sfs[num_directions], face_sf_name));
3139ae3492bSJames Wright     }
3144e2e9504SJed Brown 
3159ae3492bSJames Wright     transforms[num_directions][0][0]         = 1;
3169ae3492bSJames Wright     transforms[num_directions][1][1]         = 1;
3179ae3492bSJames Wright     transforms[num_directions][2][2]         = 1;
3189ae3492bSJames Wright     transforms[num_directions][3][3]         = 1;
3199ae3492bSJames Wright     transforms[num_directions][direction][3] = upper[direction] - lower[direction];
3209ae3492bSJames Wright     num_directions++;
3219ae3492bSJames Wright   }
3226725e60dSJed Brown 
3231fca310dSJames Wright   PetscCall(DMPlexSetIsoperiodicFaceSF(dm, num_directions, face_sfs));
3241fca310dSJames Wright   PetscCall(DMPlexSetIsoperiodicFaceTransform(dm, num_directions, (PetscScalar *)transforms));
3259ae3492bSJames Wright 
3269ae3492bSJames Wright   for (PetscInt i = 0; i < num_directions; i++) PetscCall(PetscSFDestroy(&face_sfs[i]));
3273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3284e2e9504SJed Brown }
3294e2e9504SJed Brown 
3305f06a3ddSJed Brown // This is a DMGlobalToLocalHook that applies the affine offsets. When extended for rotated periodicity, it'll need to
3315f06a3ddSJed Brown // apply a rotatonal transform and similar operations will be needed for fields (e.g., to rotate a velocity vector).
3325f06a3ddSJed Brown // We use this crude approach here so we don't have to write new GPU kernels yet.
3336725e60dSJed Brown static PetscErrorCode DMCoordAddPeriodicOffsets_Private(DM dm, Vec g, InsertMode mode, Vec l, void *ctx)
3346725e60dSJed Brown {
3356725e60dSJed Brown   PetscFunctionBegin;
336ddedc8f6SJames Wright   // These `VecScatter`s should be merged to improve efficiency; the scatters cannot be overlapped.
337ddedc8f6SJames Wright   for (PetscInt i = 0; i < dm->periodic.num_affines; i++) {
338ddedc8f6SJames Wright     PetscCall(VecScatterBegin(dm->periodic.affine_to_local[i], dm->periodic.affine[i], l, ADD_VALUES, SCATTER_FORWARD));
339ddedc8f6SJames Wright     PetscCall(VecScatterEnd(dm->periodic.affine_to_local[i], dm->periodic.affine[i], l, ADD_VALUES, SCATTER_FORWARD));
340ddedc8f6SJames Wright   }
3413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3426725e60dSJed Brown }
3436725e60dSJed Brown 
3445f06a3ddSJed Brown // Start with an SF for a positive depth (e.g., faces) and create a new SF for matched closure. The caller must ensure
3455f06a3ddSJed Brown // that both the donor (root) face and the periodic (leaf) face have consistent orientation, meaning that their closures
3465f06a3ddSJed Brown // are isomorphic. It may be useful/necessary for this restriction to be loosened.
3476725e60dSJed Brown //
3485f06a3ddSJed Brown // Output Arguments:
3495f06a3ddSJed Brown //
3505f06a3ddSJed Brown // + closure_sf - augmented point SF (see `DMGetPointSF()`) that includes the faces and all points in its closure. This
3515f06a3ddSJed Brown //   can be used to create a global section and section SF.
352b83f62b0SJames 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
3535f06a3ddSJed Brown //   transformation to periodic dofs; see DMPeriodicCoordinateSetUp_Internal().
3545f06a3ddSJed Brown //
355b83f62b0SJames Wright static PetscErrorCode DMPlexCreateIsoperiodicPointSF_Private(DM dm, PetscInt num_face_sfs, PetscSF *face_sfs, PetscSF *closure_sf, IS **is_points)
3566725e60dSJed Brown {
3576725e60dSJed Brown   MPI_Comm           comm;
3586725e60dSJed Brown   PetscMPIInt        rank;
3596725e60dSJed Brown   PetscSF            point_sf;
360b83f62b0SJames Wright   PetscInt           nroots, nleaves;
361b83f62b0SJames Wright   const PetscInt    *filocal;
362b83f62b0SJames Wright   const PetscSFNode *firemote;
3636725e60dSJed Brown 
3646725e60dSJed Brown   PetscFunctionBegin;
3656725e60dSJed Brown   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
3666725e60dSJed Brown   PetscCallMPI(MPI_Comm_rank(comm, &rank));
3676725e60dSJed Brown   PetscCall(DMGetPointSF(dm, &point_sf)); // Point SF has remote points
368b83f62b0SJames Wright   PetscCall(PetscMalloc1(num_face_sfs, is_points));
369b83f62b0SJames Wright 
370b83f62b0SJames Wright   for (PetscInt f = 0; f < num_face_sfs; f++) {
371b83f62b0SJames Wright     PetscSF   face_sf = face_sfs[f];
3726725e60dSJed Brown     PetscInt *rootdata, *leafdata;
373b83f62b0SJames Wright 
374b83f62b0SJames Wright     PetscCall(PetscSFGetGraph(face_sf, &nroots, &nleaves, &filocal, &firemote));
3756725e60dSJed Brown     PetscCall(PetscCalloc2(2 * nroots, &rootdata, 2 * nroots, &leafdata));
3766725e60dSJed Brown     for (PetscInt i = 0; i < nleaves; i++) {
3776725e60dSJed Brown       PetscInt point = filocal[i], cl_size, *closure = NULL;
3786725e60dSJed Brown       PetscCall(DMPlexGetTransitiveClosure(dm, point, PETSC_TRUE, &cl_size, &closure));
3796725e60dSJed Brown       leafdata[point] = cl_size - 1;
3806725e60dSJed Brown       PetscCall(DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, &cl_size, &closure));
3816725e60dSJed Brown     }
3826725e60dSJed Brown     PetscCall(PetscSFReduceBegin(face_sf, MPIU_INT, leafdata, rootdata + nroots, MPIU_SUM));
3836725e60dSJed Brown     PetscCall(PetscSFReduceEnd(face_sf, MPIU_INT, leafdata, rootdata + nroots, MPIU_SUM));
3846725e60dSJed Brown 
3856725e60dSJed Brown     PetscInt root_offset = 0;
3866725e60dSJed Brown     for (PetscInt p = 0; p < nroots; p++) {
3876725e60dSJed Brown       const PetscInt *donor_dof = rootdata + nroots;
3886725e60dSJed Brown       if (donor_dof[p] == 0) {
3896725e60dSJed Brown         rootdata[2 * p]     = -1;
3906725e60dSJed Brown         rootdata[2 * p + 1] = -1;
3916725e60dSJed Brown         continue;
3926725e60dSJed Brown       }
3936725e60dSJed Brown       PetscInt  cl_size;
3946725e60dSJed Brown       PetscInt *closure = NULL;
3956725e60dSJed Brown       PetscCall(DMPlexGetTransitiveClosure(dm, p, PETSC_TRUE, &cl_size, &closure));
3966725e60dSJed Brown       // cl_size - 1 = points not including self
3975f06a3ddSJed Brown       PetscAssert(donor_dof[p] == cl_size - 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Reduced leaf cone sizes do not match root cone sizes");
3986725e60dSJed Brown       rootdata[2 * p]     = root_offset;
3996725e60dSJed Brown       rootdata[2 * p + 1] = cl_size - 1;
4006725e60dSJed Brown       root_offset += cl_size - 1;
4016725e60dSJed Brown       PetscCall(DMPlexRestoreTransitiveClosure(dm, p, PETSC_TRUE, &cl_size, &closure));
4026725e60dSJed Brown     }
4036725e60dSJed Brown     PetscCall(PetscSFBcastBegin(face_sf, MPIU_2INT, rootdata, leafdata, MPI_REPLACE));
4046725e60dSJed Brown     PetscCall(PetscSFBcastEnd(face_sf, MPIU_2INT, rootdata, leafdata, MPI_REPLACE));
4056725e60dSJed Brown     // Count how many leaves we need to communicate the closures
4066725e60dSJed Brown     PetscInt leaf_offset = 0;
4076725e60dSJed Brown     for (PetscInt i = 0; i < nleaves; i++) {
4086725e60dSJed Brown       PetscInt point = filocal[i];
4096725e60dSJed Brown       if (leafdata[2 * point + 1] < 0) continue;
4106725e60dSJed Brown       leaf_offset += leafdata[2 * point + 1];
4116725e60dSJed Brown     }
4126725e60dSJed Brown 
4136725e60dSJed Brown     PetscSFNode *closure_leaf;
4146725e60dSJed Brown     PetscCall(PetscMalloc1(leaf_offset, &closure_leaf));
4156725e60dSJed Brown     leaf_offset = 0;
4166725e60dSJed Brown     for (PetscInt i = 0; i < nleaves; i++) {
4176725e60dSJed Brown       PetscInt point   = filocal[i];
4186725e60dSJed Brown       PetscInt cl_size = leafdata[2 * point + 1];
4196725e60dSJed Brown       if (cl_size < 0) continue;
4206725e60dSJed Brown       for (PetscInt j = 0; j < cl_size; j++) {
4216725e60dSJed Brown         closure_leaf[leaf_offset].rank  = firemote[i].rank;
4226725e60dSJed Brown         closure_leaf[leaf_offset].index = leafdata[2 * point] + j;
4236725e60dSJed Brown         leaf_offset++;
4246725e60dSJed Brown       }
4256725e60dSJed Brown     }
4266725e60dSJed Brown 
4276725e60dSJed Brown     PetscSF sf_closure;
4286725e60dSJed Brown     PetscCall(PetscSFCreate(comm, &sf_closure));
4296725e60dSJed Brown     PetscCall(PetscSFSetGraph(sf_closure, root_offset, leaf_offset, NULL, PETSC_USE_POINTER, closure_leaf, PETSC_OWN_POINTER));
4306725e60dSJed Brown 
431b83f62b0SJames Wright     PetscSFNode *leaf_donor_closure;
432b83f62b0SJames Wright     { // Pack root buffer with owner for every point in the root cones
4336725e60dSJed Brown       PetscSFNode       *donor_closure;
434b83f62b0SJames Wright       const PetscInt    *pilocal;
435b83f62b0SJames Wright       const PetscSFNode *piremote;
436b83f62b0SJames Wright       PetscInt           npoints;
437b83f62b0SJames Wright 
438b83f62b0SJames Wright       PetscCall(PetscSFGetGraph(point_sf, NULL, &npoints, &pilocal, &piremote));
4396725e60dSJed Brown       PetscCall(PetscCalloc1(root_offset, &donor_closure));
4406725e60dSJed Brown       root_offset = 0;
4416725e60dSJed Brown       for (PetscInt p = 0; p < nroots; p++) {
4426725e60dSJed Brown         if (rootdata[2 * p] < 0) continue;
4436725e60dSJed Brown         PetscInt  cl_size;
4446725e60dSJed Brown         PetscInt *closure = NULL;
4456725e60dSJed Brown         PetscCall(DMPlexGetTransitiveClosure(dm, p, PETSC_TRUE, &cl_size, &closure));
4466725e60dSJed Brown         for (PetscInt j = 1; j < cl_size; j++) {
4476725e60dSJed Brown           PetscInt c = closure[2 * j];
4486725e60dSJed Brown           if (pilocal) {
4496725e60dSJed Brown             PetscInt found = -1;
4506725e60dSJed Brown             if (npoints > 0) PetscCall(PetscFindInt(c, npoints, pilocal, &found));
4516725e60dSJed Brown             if (found >= 0) {
4526725e60dSJed Brown               donor_closure[root_offset++] = piremote[found];
4536725e60dSJed Brown               continue;
4546725e60dSJed Brown             }
4556725e60dSJed Brown           }
4566725e60dSJed Brown           // we own c
4576725e60dSJed Brown           donor_closure[root_offset].rank  = rank;
4586725e60dSJed Brown           donor_closure[root_offset].index = c;
4596725e60dSJed Brown           root_offset++;
4606725e60dSJed Brown         }
4616725e60dSJed Brown         PetscCall(DMPlexRestoreTransitiveClosure(dm, p, PETSC_TRUE, &cl_size, &closure));
4626725e60dSJed Brown       }
4636725e60dSJed Brown 
4646725e60dSJed Brown       PetscCall(PetscMalloc1(leaf_offset, &leaf_donor_closure));
4656497c311SBarry Smith       PetscCall(PetscSFBcastBegin(sf_closure, MPIU_SF_NODE, donor_closure, leaf_donor_closure, MPI_REPLACE));
4666497c311SBarry Smith       PetscCall(PetscSFBcastEnd(sf_closure, MPIU_SF_NODE, donor_closure, leaf_donor_closure, MPI_REPLACE));
4676725e60dSJed Brown       PetscCall(PetscSFDestroy(&sf_closure));
4686725e60dSJed Brown       PetscCall(PetscFree(donor_closure));
469b83f62b0SJames Wright     }
4706725e60dSJed Brown 
4716725e60dSJed Brown     PetscSFNode *new_iremote;
4726725e60dSJed Brown     PetscCall(PetscCalloc1(nroots, &new_iremote));
4736725e60dSJed Brown     for (PetscInt i = 0; i < nroots; i++) new_iremote[i].rank = -1;
4746725e60dSJed Brown     // Walk leaves and match vertices
4756725e60dSJed Brown     leaf_offset = 0;
4766725e60dSJed Brown     for (PetscInt i = 0; i < nleaves; i++) {
4776725e60dSJed Brown       PetscInt  point   = filocal[i], cl_size;
4786725e60dSJed Brown       PetscInt *closure = NULL;
4796725e60dSJed Brown       PetscCall(DMPlexGetTransitiveClosure(dm, point, PETSC_TRUE, &cl_size, &closure));
4806725e60dSJed Brown       for (PetscInt j = 1; j < cl_size; j++) { // TODO: should we send donor edge orientations so we can flip for consistency?
4816725e60dSJed Brown         PetscInt    c  = closure[2 * j];
4826725e60dSJed Brown         PetscSFNode lc = leaf_donor_closure[leaf_offset];
4836725e60dSJed Brown         // printf("[%d] face %d.%d: %d ?-- (%d,%d)\n", rank, point, j, c, lc.rank, lc.index);
4846725e60dSJed Brown         if (new_iremote[c].rank == -1) {
4856725e60dSJed Brown           new_iremote[c] = lc;
4865f06a3ddSJed 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");
4876725e60dSJed Brown         leaf_offset++;
4886725e60dSJed Brown       }
4896725e60dSJed Brown       PetscCall(DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, &cl_size, &closure));
4906725e60dSJed Brown     }
4916725e60dSJed Brown     PetscCall(PetscFree(leaf_donor_closure));
4926725e60dSJed Brown 
4936725e60dSJed Brown     // Include face points in closure SF
4946725e60dSJed Brown     for (PetscInt i = 0; i < nleaves; i++) new_iremote[filocal[i]] = firemote[i];
4956725e60dSJed Brown     // consolidate leaves
4966725e60dSJed Brown     PetscInt num_new_leaves = 0;
4976725e60dSJed Brown     for (PetscInt i = 0; i < nroots; i++) {
4986725e60dSJed Brown       if (new_iremote[i].rank == -1) continue;
4996725e60dSJed Brown       new_iremote[num_new_leaves] = new_iremote[i];
5006725e60dSJed Brown       leafdata[num_new_leaves]    = i;
5016725e60dSJed Brown       num_new_leaves++;
5026725e60dSJed Brown     }
503b83f62b0SJames Wright     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, num_new_leaves, leafdata, PETSC_COPY_VALUES, &(*is_points)[f]));
5046725e60dSJed Brown 
5056725e60dSJed Brown     PetscSF csf;
5066725e60dSJed Brown     PetscCall(PetscSFCreate(comm, &csf));
5076725e60dSJed Brown     PetscCall(PetscSFSetGraph(csf, nroots, num_new_leaves, leafdata, PETSC_COPY_VALUES, new_iremote, PETSC_COPY_VALUES));
5086725e60dSJed Brown     PetscCall(PetscFree(new_iremote)); // copy and delete because new_iremote is longer than it needs to be
5096725e60dSJed Brown     PetscCall(PetscFree2(rootdata, leafdata));
5106725e60dSJed Brown 
511b83f62b0SJames Wright     PetscInt npoints;
512b83f62b0SJames Wright     PetscCall(PetscSFGetGraph(point_sf, NULL, &npoints, NULL, NULL));
5135f06a3ddSJed Brown     if (npoints < 0) { // empty point_sf
5146725e60dSJed Brown       *closure_sf = csf;
5155f06a3ddSJed Brown     } else {
5165f06a3ddSJed Brown       PetscCall(PetscSFMerge(point_sf, csf, closure_sf));
5175f06a3ddSJed Brown       PetscCall(PetscSFDestroy(&csf));
5185f06a3ddSJed Brown     }
519d8b4a066SPierre Jolivet     if (f > 0) PetscCall(PetscSFDestroy(&point_sf)); // Only destroy if point_sf is from previous calls to PetscSFMerge
520b83f62b0SJames Wright     point_sf = *closure_sf;                          // Use combined point + isoperiodic SF to define point ownership for further face_sf
521b83f62b0SJames Wright   }
5225f06a3ddSJed Brown   PetscCall(PetscObjectSetName((PetscObject)*closure_sf, "Composed Periodic Points"));
5233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5246725e60dSJed Brown }
5256725e60dSJed Brown 
5265f06a3ddSJed Brown static PetscErrorCode DMGetIsoperiodicPointSF_Plex(DM dm, PetscSF *sf)
5276725e60dSJed Brown {
5286725e60dSJed Brown   DM_Plex *plex = (DM_Plex *)dm->data;
5296725e60dSJed Brown 
5306725e60dSJed Brown   PetscFunctionBegin;
531b83f62b0SJames 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));
5326725e60dSJed Brown   if (sf) *sf = plex->periodic.composed_sf;
5333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5346725e60dSJed Brown }
5356725e60dSJed Brown 
5365f06a3ddSJed Brown PetscErrorCode DMPlexMigrateIsoperiodicFaceSF_Internal(DM old_dm, DM dm, PetscSF sf_migration)
5375f06a3ddSJed Brown {
5385f06a3ddSJed Brown   DM_Plex    *plex = (DM_Plex *)old_dm->data;
539a45b0079SJames Wright   PetscSF     sf_point, *new_face_sfs;
5405f06a3ddSJed Brown   PetscMPIInt rank;
5415f06a3ddSJed Brown 
5425f06a3ddSJed Brown   PetscFunctionBegin;
5431fca310dSJames Wright   if (!plex->periodic.face_sfs) PetscFunctionReturn(PETSC_SUCCESS);
5445f06a3ddSJed Brown   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
5455f06a3ddSJed Brown   PetscCall(DMGetPointSF(dm, &sf_point));
546a45b0079SJames Wright   PetscCall(PetscMalloc1(plex->periodic.num_face_sfs, &new_face_sfs));
547a45b0079SJames Wright 
548a45b0079SJames Wright   for (PetscInt f = 0; f < plex->periodic.num_face_sfs; f++) {
5495f06a3ddSJed Brown     PetscInt           old_npoints, new_npoints, old_nleaf, new_nleaf, point_nleaf;
5505f06a3ddSJed Brown     PetscSFNode       *new_leafdata, *rootdata, *leafdata;
5515f06a3ddSJed Brown     const PetscInt    *old_local, *point_local;
5525f06a3ddSJed Brown     const PetscSFNode *old_remote, *point_remote;
5536497c311SBarry Smith 
554a45b0079SJames Wright     PetscCall(PetscSFGetGraph(plex->periodic.face_sfs[f], &old_npoints, &old_nleaf, &old_local, &old_remote));
5555f06a3ddSJed Brown     PetscCall(PetscSFGetGraph(sf_migration, NULL, &new_nleaf, NULL, NULL));
5565f06a3ddSJed Brown     PetscCall(PetscSFGetGraph(sf_point, &new_npoints, &point_nleaf, &point_local, &point_remote));
5575f06a3ddSJed Brown     PetscAssert(new_nleaf == new_npoints, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected migration leaf space to match new point root space");
5585f06a3ddSJed Brown     PetscCall(PetscMalloc3(old_npoints, &rootdata, old_npoints, &leafdata, new_npoints, &new_leafdata));
5595f06a3ddSJed Brown 
5605f06a3ddSJed Brown     // Fill new_leafdata with new owners of all points
5615f06a3ddSJed Brown     for (PetscInt i = 0; i < new_npoints; i++) {
5625f06a3ddSJed Brown       new_leafdata[i].rank  = rank;
5635f06a3ddSJed Brown       new_leafdata[i].index = i;
5645f06a3ddSJed Brown     }
5655f06a3ddSJed Brown     for (PetscInt i = 0; i < point_nleaf; i++) {
5665f06a3ddSJed Brown       PetscInt j      = point_local[i];
5675f06a3ddSJed Brown       new_leafdata[j] = point_remote[i];
5685f06a3ddSJed Brown     }
5695f06a3ddSJed Brown     // REPLACE is okay because every leaf agrees about the new owners
5706497c311SBarry Smith     PetscCall(PetscSFReduceBegin(sf_migration, MPIU_SF_NODE, new_leafdata, rootdata, MPI_REPLACE));
5716497c311SBarry Smith     PetscCall(PetscSFReduceEnd(sf_migration, MPIU_SF_NODE, new_leafdata, rootdata, MPI_REPLACE));
5725f06a3ddSJed Brown     // rootdata now contains the new owners
5735f06a3ddSJed Brown 
5745f06a3ddSJed Brown     // Send to leaves of old space
5755f06a3ddSJed Brown     for (PetscInt i = 0; i < old_npoints; i++) {
5765f06a3ddSJed Brown       leafdata[i].rank  = -1;
5775f06a3ddSJed Brown       leafdata[i].index = -1;
5785f06a3ddSJed Brown     }
5796497c311SBarry Smith     PetscCall(PetscSFBcastBegin(plex->periodic.face_sfs[f], MPIU_SF_NODE, rootdata, leafdata, MPI_REPLACE));
5806497c311SBarry Smith     PetscCall(PetscSFBcastEnd(plex->periodic.face_sfs[f], MPIU_SF_NODE, rootdata, leafdata, MPI_REPLACE));
5815f06a3ddSJed Brown 
5825f06a3ddSJed Brown     // Send to new leaf space
5836497c311SBarry Smith     PetscCall(PetscSFBcastBegin(sf_migration, MPIU_SF_NODE, leafdata, new_leafdata, MPI_REPLACE));
5846497c311SBarry Smith     PetscCall(PetscSFBcastEnd(sf_migration, MPIU_SF_NODE, leafdata, new_leafdata, MPI_REPLACE));
5855f06a3ddSJed Brown 
5865f06a3ddSJed Brown     PetscInt     nface = 0, *new_local;
5875f06a3ddSJed Brown     PetscSFNode *new_remote;
5885f06a3ddSJed Brown     for (PetscInt i = 0; i < new_npoints; i++) nface += (new_leafdata[i].rank >= 0);
5895f06a3ddSJed Brown     PetscCall(PetscMalloc1(nface, &new_local));
5905f06a3ddSJed Brown     PetscCall(PetscMalloc1(nface, &new_remote));
5915f06a3ddSJed Brown     nface = 0;
5925f06a3ddSJed Brown     for (PetscInt i = 0; i < new_npoints; i++) {
5935f06a3ddSJed Brown       if (new_leafdata[i].rank == -1) continue;
5945f06a3ddSJed Brown       new_local[nface]  = i;
5955f06a3ddSJed Brown       new_remote[nface] = new_leafdata[i];
5965f06a3ddSJed Brown       nface++;
5975f06a3ddSJed Brown     }
5985f06a3ddSJed Brown     PetscCall(PetscFree3(rootdata, leafdata, new_leafdata));
599a45b0079SJames Wright     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), &new_face_sfs[f]));
600a45b0079SJames Wright     PetscCall(PetscSFSetGraph(new_face_sfs[f], new_npoints, nface, new_local, PETSC_OWN_POINTER, new_remote, PETSC_OWN_POINTER));
601a45b0079SJames Wright     {
602a45b0079SJames Wright       char new_face_sf_name[PETSC_MAX_PATH_LEN];
603a45b0079SJames Wright       PetscCall(PetscSNPrintf(new_face_sf_name, sizeof new_face_sf_name, "Migrated Isoperiodic Faces #%" PetscInt_FMT, f));
604a45b0079SJames Wright       PetscCall(PetscObjectSetName((PetscObject)new_face_sfs[f], new_face_sf_name));
605a45b0079SJames Wright     }
606a45b0079SJames Wright   }
607a45b0079SJames Wright 
608a45b0079SJames Wright   PetscCall(DMPlexSetIsoperiodicFaceSF(dm, plex->periodic.num_face_sfs, new_face_sfs));
609a45b0079SJames Wright   PetscCall(DMPlexSetIsoperiodicFaceTransform(dm, plex->periodic.num_face_sfs, (PetscScalar *)plex->periodic.transform));
610a45b0079SJames Wright   for (PetscInt f = 0; f < plex->periodic.num_face_sfs; f++) PetscCall(PetscSFDestroy(&new_face_sfs[f]));
611a45b0079SJames Wright   PetscCall(PetscFree(new_face_sfs));
6123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6135f06a3ddSJed Brown }
6145f06a3ddSJed Brown 
6156725e60dSJed Brown PetscErrorCode DMPeriodicCoordinateSetUp_Internal(DM dm)
6166725e60dSJed Brown {
6176725e60dSJed Brown   DM_Plex   *plex = (DM_Plex *)dm->data;
6186497c311SBarry Smith   PetscCount count;
619ddedc8f6SJames Wright   IS         isdof;
620ddedc8f6SJames Wright   PetscInt   dim;
6214d86920dSPierre Jolivet 
6226725e60dSJed Brown   PetscFunctionBegin;
6231fca310dSJames Wright   if (!plex->periodic.face_sfs) PetscFunctionReturn(PETSC_SUCCESS);
6245f06a3ddSJed Brown   PetscCall(DMGetIsoperiodicPointSF_Plex(dm, NULL));
6255f06a3ddSJed Brown   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMGetIsoperiodicPointSF_C", DMGetIsoperiodicPointSF_Plex));
6266725e60dSJed Brown 
6276725e60dSJed Brown   PetscCall(DMGetDimension(dm, &dim));
628ddedc8f6SJames Wright   dm->periodic.num_affines = plex->periodic.num_face_sfs;
629ddedc8f6SJames Wright   PetscCall(PetscMalloc2(dm->periodic.num_affines, &dm->periodic.affine_to_local, dm->periodic.num_affines, &dm->periodic.affine));
630ddedc8f6SJames Wright 
631ddedc8f6SJames Wright   for (PetscInt f = 0; f < plex->periodic.num_face_sfs; f++) {
6326497c311SBarry Smith     PetscInt        npoints, vsize, isize;
6336725e60dSJed Brown     const PetscInt *points;
634ddedc8f6SJames Wright     IS              is = plex->periodic.periodic_points[f];
6356725e60dSJed Brown     PetscSegBuffer  seg;
6366725e60dSJed Brown     PetscSection    section;
6376497c311SBarry Smith     PetscInt       *ind;
6386497c311SBarry Smith     Vec             L, P;
6396497c311SBarry Smith     VecType         vec_type;
6406497c311SBarry Smith     VecScatter      scatter;
6416497c311SBarry Smith     PetscScalar    *x;
6426497c311SBarry Smith 
6436725e60dSJed Brown     PetscCall(DMGetLocalSection(dm, &section));
6446725e60dSJed Brown     PetscCall(PetscSegBufferCreate(sizeof(PetscInt), 32, &seg));
6456725e60dSJed Brown     PetscCall(ISGetSize(is, &npoints));
6466725e60dSJed Brown     PetscCall(ISGetIndices(is, &points));
6476725e60dSJed Brown     for (PetscInt i = 0; i < npoints; i++) {
6486725e60dSJed Brown       PetscInt point = points[i], off, dof;
6496497c311SBarry Smith 
6506725e60dSJed Brown       PetscCall(PetscSectionGetOffset(section, point, &off));
6516725e60dSJed Brown       PetscCall(PetscSectionGetDof(section, point, &dof));
6526725e60dSJed Brown       PetscAssert(dof % dim == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unexpected dof %" PetscInt_FMT " not divisible by dimension %" PetscInt_FMT, dof, dim);
6536725e60dSJed Brown       for (PetscInt j = 0; j < dof / dim; j++) {
6546725e60dSJed Brown         PetscInt *slot;
6556497c311SBarry Smith 
6566725e60dSJed Brown         PetscCall(PetscSegBufferGetInts(seg, 1, &slot));
657729bdd54SJed Brown         *slot = off / dim + j;
6586725e60dSJed Brown       }
6596725e60dSJed Brown     }
6606725e60dSJed Brown     PetscCall(PetscSegBufferGetSize(seg, &count));
6616725e60dSJed Brown     PetscCall(PetscSegBufferExtractAlloc(seg, &ind));
6626725e60dSJed Brown     PetscCall(PetscSegBufferDestroy(&seg));
6636497c311SBarry Smith     PetscCall(PetscIntCast(count, &isize));
6646497c311SBarry Smith     PetscCall(ISCreateBlock(PETSC_COMM_SELF, dim, isize, ind, PETSC_OWN_POINTER, &isdof));
6656497c311SBarry Smith 
6666497c311SBarry Smith     PetscCall(PetscIntCast(count * dim, &vsize));
6676725e60dSJed Brown     PetscCall(DMGetLocalVector(dm, &L));
6686725e60dSJed Brown     PetscCall(VecCreate(PETSC_COMM_SELF, &P));
6696497c311SBarry Smith     PetscCall(VecSetSizes(P, vsize, vsize));
6706725e60dSJed Brown     PetscCall(VecGetType(L, &vec_type));
6716725e60dSJed Brown     PetscCall(VecSetType(P, vec_type));
6726725e60dSJed Brown     PetscCall(VecScatterCreate(P, NULL, L, isdof, &scatter));
6736725e60dSJed Brown     PetscCall(DMRestoreLocalVector(dm, &L));
6746725e60dSJed Brown     PetscCall(ISDestroy(&isdof));
6756725e60dSJed Brown 
6766725e60dSJed Brown     PetscCall(VecGetArrayWrite(P, &x));
6776497c311SBarry Smith     for (PetscCount i = 0; i < count; i++) {
678ddedc8f6SJames Wright       for (PetscInt j = 0; j < dim; j++) x[i * dim + j] = plex->periodic.transform[f][j][3];
6796725e60dSJed Brown     }
6806725e60dSJed Brown     PetscCall(VecRestoreArrayWrite(P, &x));
6816725e60dSJed Brown 
682ddedc8f6SJames Wright     dm->periodic.affine_to_local[f] = scatter;
683ddedc8f6SJames Wright     dm->periodic.affine[f]          = P;
684ddedc8f6SJames Wright   }
6856725e60dSJed Brown   PetscCall(DMGlobalToLocalHookAdd(dm, NULL, DMCoordAddPeriodicOffsets_Private, NULL));
6863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6876725e60dSJed Brown }
6886725e60dSJed Brown 
6895f06a3ddSJed Brown // We'll just orient all the edges, though only periodic boundary edges need orientation
6905f06a3ddSJed Brown static PetscErrorCode DMPlexOrientPositiveEdges_Private(DM dm)
6915f06a3ddSJed Brown {
6925f06a3ddSJed Brown   PetscInt dim, eStart, eEnd;
6934d86920dSPierre Jolivet 
6945f06a3ddSJed Brown   PetscFunctionBegin;
6955f06a3ddSJed Brown   PetscCall(DMGetDimension(dm, &dim));
6963ba16761SJacob Faibussowitsch   if (dim < 3) PetscFunctionReturn(PETSC_SUCCESS); // not necessary
6975f06a3ddSJed Brown   PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd));
6985f06a3ddSJed Brown   for (PetscInt e = eStart; e < eEnd; e++) {
6995f06a3ddSJed Brown     const PetscInt *cone;
7005f06a3ddSJed Brown     PetscCall(DMPlexGetCone(dm, e, &cone));
7015f06a3ddSJed Brown     if (cone[0] > cone[1]) PetscCall(DMPlexOrientPoint(dm, e, -1));
7025f06a3ddSJed Brown   }
7033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7045f06a3ddSJed Brown }
7055f06a3ddSJed Brown 
7063e72e933SJed Brown PetscErrorCode DMPlexCreateBoxMesh_Tensor_SFC_Internal(DM dm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate)
7073e72e933SJed Brown {
7083e72e933SJed Brown   PetscInt  eextent[3] = {1, 1, 1}, vextent[3] = {1, 1, 1};
7093e72e933SJed Brown   const Ijk closure_1[] = {
7103e72e933SJed Brown     {0, 0, 0},
7113e72e933SJed Brown     {1, 0, 0},
7123e72e933SJed Brown   };
7133e72e933SJed Brown   const Ijk closure_2[] = {
7143e72e933SJed Brown     {0, 0, 0},
7153e72e933SJed Brown     {1, 0, 0},
7163e72e933SJed Brown     {1, 1, 0},
7173e72e933SJed Brown     {0, 1, 0},
7183e72e933SJed Brown   };
7193e72e933SJed Brown   const Ijk closure_3[] = {
7203e72e933SJed Brown     {0, 0, 0},
7213e72e933SJed Brown     {0, 1, 0},
7223e72e933SJed Brown     {1, 1, 0},
7233e72e933SJed Brown     {1, 0, 0},
7243e72e933SJed Brown     {0, 0, 1},
7253e72e933SJed Brown     {1, 0, 1},
7263e72e933SJed Brown     {1, 1, 1},
7273e72e933SJed Brown     {0, 1, 1},
7283e72e933SJed Brown   };
7293431e603SJed Brown   const Ijk *const closure_dim[] = {NULL, closure_1, closure_2, closure_3};
7303431e603SJed Brown   // This must be kept consistent with DMPlexCreateCubeMesh_Internal
7313431e603SJed Brown   const PetscInt        face_marker_1[]   = {1, 2};
7323431e603SJed Brown   const PetscInt        face_marker_2[]   = {4, 2, 1, 3};
7333431e603SJed Brown   const PetscInt        face_marker_3[]   = {6, 5, 3, 4, 1, 2};
7343431e603SJed Brown   const PetscInt *const face_marker_dim[] = {NULL, face_marker_1, face_marker_2, face_marker_3};
7356725e60dSJed Brown   // Orient faces so the normal is in the positive axis and the first vertex is the one closest to zero.
7366725e60dSJed Brown   // These orientations can be determined by examining cones of a reference quad and hex element.
7376725e60dSJed Brown   const PetscInt        face_orient_1[]   = {0, 0};
7386725e60dSJed Brown   const PetscInt        face_orient_2[]   = {-1, 0, 0, -1};
7396725e60dSJed Brown   const PetscInt        face_orient_3[]   = {-2, 0, -2, 1, -2, 0};
7406725e60dSJed Brown   const PetscInt *const face_orient_dim[] = {NULL, face_orient_1, face_orient_2, face_orient_3};
7413e72e933SJed Brown 
7423e72e933SJed Brown   PetscFunctionBegin;
743ea7b3863SJames Wright   PetscCall(PetscLogEventBegin(DMPLEX_CreateBoxSFC, dm, 0, 0, 0));
7444f572ea9SToby Isaac   PetscAssertPointer(dm, 1);
7453e72e933SJed Brown   PetscValidLogicalCollectiveInt(dm, dim, 2);
7463e72e933SJed Brown   PetscCall(DMSetDimension(dm, dim));
7473e72e933SJed Brown   PetscMPIInt rank, size;
7483e72e933SJed Brown   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
7493e72e933SJed Brown   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
7503e72e933SJed Brown   for (PetscInt i = 0; i < dim; i++) {
7513e72e933SJed Brown     eextent[i] = faces[i];
7523e72e933SJed Brown     vextent[i] = faces[i] + 1;
7533e72e933SJed Brown   }
754c77877e3SJed Brown   ZLayout layout;
755c77877e3SJed Brown   PetscCall(ZLayoutCreate(size, eextent, vextent, &layout));
7563e72e933SJed Brown   PetscZSet vset; // set of all vertices in the closure of the owned elements
7573e72e933SJed Brown   PetscCall(PetscZSetCreate(&vset));
7583e72e933SJed Brown   PetscInt local_elems = 0;
7593e72e933SJed Brown   for (ZCode z = layout.zstarts[rank]; z < layout.zstarts[rank + 1]; z++) {
7603e72e933SJed Brown     Ijk loc = ZCodeSplit(z);
7613ba16761SJacob Faibussowitsch     if (IjkActive(layout.vextent, loc)) PetscCall(PetscZSetAdd(vset, z));
7626547ddbdSJed Brown     else {
7636547ddbdSJed Brown       z += ZStepOct(z);
7646547ddbdSJed Brown       continue;
7656547ddbdSJed Brown     }
7663e72e933SJed Brown     if (IjkActive(layout.eextent, loc)) {
7673e72e933SJed Brown       local_elems++;
7683e72e933SJed Brown       // Add all neighboring vertices to set
7693e72e933SJed Brown       for (PetscInt n = 0; n < PetscPowInt(2, dim); n++) {
7703e72e933SJed Brown         Ijk   inc  = closure_dim[dim][n];
7713e72e933SJed Brown         Ijk   nloc = {loc.i + inc.i, loc.j + inc.j, loc.k + inc.k};
7723e72e933SJed Brown         ZCode v    = ZEncode(nloc);
7733ba16761SJacob Faibussowitsch         PetscCall(PetscZSetAdd(vset, v));
7743e72e933SJed Brown       }
7753e72e933SJed Brown     }
7763e72e933SJed Brown   }
7773e72e933SJed Brown   PetscInt local_verts, off = 0;
7783e72e933SJed Brown   ZCode   *vert_z;
7793e72e933SJed Brown   PetscCall(PetscZSetGetSize(vset, &local_verts));
7803e72e933SJed Brown   PetscCall(PetscMalloc1(local_verts, &vert_z));
7813e72e933SJed Brown   PetscCall(PetscZSetGetElems(vset, &off, vert_z));
7823e72e933SJed Brown   PetscCall(PetscZSetDestroy(&vset));
7833e72e933SJed Brown   // ZCode is unsigned for bitwise convenience, but highest bit should never be set, so can interpret as signed
7843e72e933SJed Brown   PetscCall(PetscSortInt64(local_verts, (PetscInt64 *)vert_z));
7853e72e933SJed Brown 
7863e72e933SJed Brown   PetscCall(DMPlexSetChart(dm, 0, local_elems + local_verts));
7873e72e933SJed Brown   for (PetscInt e = 0; e < local_elems; e++) PetscCall(DMPlexSetConeSize(dm, e, PetscPowInt(2, dim)));
7883e72e933SJed Brown   PetscCall(DMSetUp(dm));
7893e72e933SJed Brown   {
7903e72e933SJed Brown     PetscInt e = 0;
7913e72e933SJed Brown     for (ZCode z = layout.zstarts[rank]; z < layout.zstarts[rank + 1]; z++) {
7923e72e933SJed Brown       Ijk loc = ZCodeSplit(z);
7936547ddbdSJed Brown       if (!IjkActive(layout.eextent, loc)) {
7946547ddbdSJed Brown         z += ZStepOct(z);
7956547ddbdSJed Brown         continue;
7966547ddbdSJed Brown       }
7973e72e933SJed Brown       PetscInt cone[8], orient[8] = {0};
7983e72e933SJed Brown       for (PetscInt n = 0; n < PetscPowInt(2, dim); n++) {
7993e72e933SJed Brown         Ijk      inc  = closure_dim[dim][n];
8003e72e933SJed Brown         Ijk      nloc = {loc.i + inc.i, loc.j + inc.j, loc.k + inc.k};
8013e72e933SJed Brown         ZCode    v    = ZEncode(nloc);
8023e72e933SJed Brown         PetscInt ci   = ZCodeFind(v, local_verts, vert_z);
8033e72e933SJed Brown         PetscAssert(ci >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Did not find neighbor vertex in set");
8043e72e933SJed Brown         cone[n] = local_elems + ci;
8053e72e933SJed Brown       }
8063e72e933SJed Brown       PetscCall(DMPlexSetCone(dm, e, cone));
8073e72e933SJed Brown       PetscCall(DMPlexSetConeOrientation(dm, e, orient));
8083e72e933SJed Brown       e++;
8093e72e933SJed Brown     }
8103e72e933SJed Brown   }
8113e72e933SJed Brown 
8123e72e933SJed Brown   PetscCall(DMPlexSymmetrize(dm));
8133e72e933SJed Brown   PetscCall(DMPlexStratify(dm));
8146725e60dSJed Brown 
8153e72e933SJed Brown   { // Create point SF
8163e72e933SJed Brown     PetscSF sf;
8173ba16761SJacob Faibussowitsch     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), &sf));
8183e72e933SJed Brown     PetscInt owned_verts = ZCodeFind(layout.zstarts[rank + 1], local_verts, vert_z);
8193e72e933SJed Brown     if (owned_verts < 0) owned_verts = -(owned_verts + 1); // We don't care whether the key was found
8203e72e933SJed Brown     PetscInt     num_ghosts = local_verts - owned_verts;   // Due to sorting, owned vertices always come first
8213e72e933SJed Brown     PetscInt    *local_ghosts;
8223e72e933SJed Brown     PetscSFNode *ghosts;
8233e72e933SJed Brown     PetscCall(PetscMalloc1(num_ghosts, &local_ghosts));
8243e72e933SJed Brown     PetscCall(PetscMalloc1(num_ghosts, &ghosts));
8253e72e933SJed Brown     for (PetscInt i = 0; i < num_ghosts;) {
8263e72e933SJed Brown       ZCode       z           = vert_z[owned_verts + i];
8276497c311SBarry Smith       PetscMPIInt remote_rank = (PetscMPIInt)ZCodeFind(z, size + 1, layout.zstarts), remote_count = 0;
8283e72e933SJed Brown       if (remote_rank < 0) remote_rank = -(remote_rank + 1) - 1;
8293e72e933SJed Brown       // We have a new remote rank; find all the ghost indices (which are contiguous in vert_z)
8303e72e933SJed Brown 
8313e72e933SJed Brown       // Count the elements on remote_rank
8324e2e9504SJed Brown       PetscInt remote_elem = ZLayoutElementsOnRank(&layout, remote_rank);
8333e72e933SJed Brown 
8343e72e933SJed Brown       // Traverse vertices and make ghost links
8353e72e933SJed Brown       for (ZCode rz = layout.zstarts[remote_rank]; rz < layout.zstarts[remote_rank + 1]; rz++) {
8363e72e933SJed Brown         Ijk loc = ZCodeSplit(rz);
8373e72e933SJed Brown         if (rz == z) {
8383e72e933SJed Brown           local_ghosts[i] = local_elems + owned_verts + i;
8393e72e933SJed Brown           ghosts[i].rank  = remote_rank;
8403e72e933SJed Brown           ghosts[i].index = remote_elem + remote_count;
8413e72e933SJed Brown           i++;
8423e72e933SJed Brown           if (i == num_ghosts) break;
8433e72e933SJed Brown           z = vert_z[owned_verts + i];
8443e72e933SJed Brown         }
8453e72e933SJed Brown         if (IjkActive(layout.vextent, loc)) remote_count++;
8466547ddbdSJed Brown         else rz += ZStepOct(rz);
8473e72e933SJed Brown       }
8483e72e933SJed Brown     }
8493e72e933SJed Brown     PetscCall(PetscSFSetGraph(sf, local_elems + local_verts, num_ghosts, local_ghosts, PETSC_OWN_POINTER, ghosts, PETSC_OWN_POINTER));
8503e72e933SJed Brown     PetscCall(PetscObjectSetName((PetscObject)sf, "SFC Point SF"));
8513e72e933SJed Brown     PetscCall(DMSetPointSF(dm, sf));
8523e72e933SJed Brown     PetscCall(PetscSFDestroy(&sf));
8533e72e933SJed Brown   }
8543e72e933SJed Brown   {
8553e72e933SJed Brown     Vec          coordinates;
8563e72e933SJed Brown     PetscScalar *coords;
8573e72e933SJed Brown     PetscSection coord_section;
8583e72e933SJed Brown     PetscInt     coord_size;
8593e72e933SJed Brown     PetscCall(DMGetCoordinateSection(dm, &coord_section));
8603e72e933SJed Brown     PetscCall(PetscSectionSetNumFields(coord_section, 1));
8613e72e933SJed Brown     PetscCall(PetscSectionSetFieldComponents(coord_section, 0, dim));
8623e72e933SJed Brown     PetscCall(PetscSectionSetChart(coord_section, local_elems, local_elems + local_verts));
8633e72e933SJed Brown     for (PetscInt v = 0; v < local_verts; v++) {
8643e72e933SJed Brown       PetscInt point = local_elems + v;
8653e72e933SJed Brown       PetscCall(PetscSectionSetDof(coord_section, point, dim));
8663e72e933SJed Brown       PetscCall(PetscSectionSetFieldDof(coord_section, point, 0, dim));
8673e72e933SJed Brown     }
8683e72e933SJed Brown     PetscCall(PetscSectionSetUp(coord_section));
8693e72e933SJed Brown     PetscCall(PetscSectionGetStorageSize(coord_section, &coord_size));
8703e72e933SJed Brown     PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates));
8713e72e933SJed Brown     PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
8723e72e933SJed Brown     PetscCall(VecSetSizes(coordinates, coord_size, PETSC_DETERMINE));
8733e72e933SJed Brown     PetscCall(VecSetBlockSize(coordinates, dim));
8743e72e933SJed Brown     PetscCall(VecSetType(coordinates, VECSTANDARD));
8753e72e933SJed Brown     PetscCall(VecGetArray(coordinates, &coords));
8763e72e933SJed Brown     for (PetscInt v = 0; v < local_verts; v++) {
8773e72e933SJed Brown       Ijk loc             = ZCodeSplit(vert_z[v]);
8783e72e933SJed Brown       coords[v * dim + 0] = lower[0] + loc.i * (upper[0] - lower[0]) / layout.eextent.i;
8793e72e933SJed Brown       if (dim > 1) coords[v * dim + 1] = lower[1] + loc.j * (upper[1] - lower[1]) / layout.eextent.j;
8803e72e933SJed Brown       if (dim > 2) coords[v * dim + 2] = lower[2] + loc.k * (upper[2] - lower[2]) / layout.eextent.k;
8813e72e933SJed Brown     }
8823e72e933SJed Brown     PetscCall(VecRestoreArray(coordinates, &coords));
8833e72e933SJed Brown     PetscCall(DMSetCoordinatesLocal(dm, coordinates));
8843e72e933SJed Brown     PetscCall(VecDestroy(&coordinates));
8853e72e933SJed Brown   }
8863e72e933SJed Brown   if (interpolate) {
8873431e603SJed Brown     PetscCall(DMPlexInterpolateInPlace_Internal(dm));
8885f06a3ddSJed Brown     // It's currently necessary to orient the donor and periodic edges consistently. An easy way to ensure that is ot
8895f06a3ddSJed Brown     // give all edges positive orientation. Since vertices are created in Z-order, all ranks will agree about the
8905f06a3ddSJed Brown     // ordering cone[0] < cone[1]. This is overkill and it would be nice to remove this preparation and make
8915f06a3ddSJed Brown     // DMPlexCreateIsoperiodicClosureSF_Private() more resilient, so it fixes any inconsistent orientations. That might
8925f06a3ddSJed Brown     // be needed in a general CGNS reader, for example.
8935f06a3ddSJed Brown     PetscCall(DMPlexOrientPositiveEdges_Private(dm));
8943431e603SJed Brown 
8953431e603SJed Brown     DMLabel label;
8963431e603SJed Brown     PetscCall(DMCreateLabel(dm, "Face Sets"));
8973431e603SJed Brown     PetscCall(DMGetLabel(dm, "Face Sets", &label));
8989ae3492bSJames Wright     PetscSegBuffer per_faces[3], donor_face_closure[3], my_donor_faces[3];
8999ae3492bSJames Wright     for (PetscInt i = 0; i < 3; i++) {
9009ae3492bSJames Wright       PetscCall(PetscSegBufferCreate(sizeof(PetscInt), 64, &per_faces[i]));
9019ae3492bSJames Wright       PetscCall(PetscSegBufferCreate(sizeof(PetscInt), 64, &my_donor_faces[i]));
9029ae3492bSJames Wright       PetscCall(PetscSegBufferCreate(sizeof(ZCode), 64 * PetscPowInt(2, dim), &donor_face_closure[i]));
9039ae3492bSJames Wright     }
9043431e603SJed Brown     PetscInt fStart, fEnd, vStart, vEnd;
9053431e603SJed Brown     PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
9063431e603SJed Brown     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
9073431e603SJed Brown     for (PetscInt f = fStart; f < fEnd; f++) {
9084e2e9504SJed Brown       PetscInt npoints, *points = NULL, num_fverts = 0, fverts[8];
9093431e603SJed Brown       PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &npoints, &points));
9103431e603SJed Brown       PetscInt bc_count[6] = {0};
9113431e603SJed Brown       for (PetscInt i = 0; i < npoints; i++) {
9123431e603SJed Brown         PetscInt p = points[2 * i];
9133431e603SJed Brown         if (p < vStart || vEnd <= p) continue;
9144e2e9504SJed Brown         fverts[num_fverts++] = p;
9153431e603SJed Brown         Ijk loc              = ZCodeSplit(vert_z[p - vStart]);
9163431e603SJed Brown         // Convention here matches DMPlexCreateCubeMesh_Internal
9173431e603SJed Brown         bc_count[0] += loc.i == 0;
9183431e603SJed Brown         bc_count[1] += loc.i == layout.vextent.i - 1;
9193431e603SJed Brown         bc_count[2] += loc.j == 0;
9203431e603SJed Brown         bc_count[3] += loc.j == layout.vextent.j - 1;
9213431e603SJed Brown         bc_count[4] += loc.k == 0;
9223431e603SJed Brown         bc_count[5] += loc.k == layout.vextent.k - 1;
9233e72e933SJed Brown       }
9243431e603SJed Brown       PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &npoints, &points));
9253431e603SJed Brown       for (PetscInt bc = 0, bc_match = 0; bc < 2 * dim; bc++) {
9263431e603SJed Brown         if (bc_count[bc] == PetscPowInt(2, dim - 1)) {
9276725e60dSJed Brown           PetscCall(DMPlexOrientPoint(dm, f, face_orient_dim[dim][bc]));
9284e2e9504SJed Brown           if (periodicity[bc / 2] == DM_BOUNDARY_PERIODIC) {
9294e2e9504SJed Brown             PetscInt *put;
9304e2e9504SJed Brown             if (bc % 2 == 0) { // donor face; no label
9319ae3492bSJames Wright               PetscCall(PetscSegBufferGet(my_donor_faces[bc / 2], 1, &put));
9324e2e9504SJed Brown               *put = f;
9334e2e9504SJed Brown             } else { // periodic face
9349ae3492bSJames Wright               PetscCall(PetscSegBufferGet(per_faces[bc / 2], 1, &put));
9354e2e9504SJed Brown               *put = f;
9364e2e9504SJed Brown               ZCode *zput;
9379ae3492bSJames Wright               PetscCall(PetscSegBufferGet(donor_face_closure[bc / 2], num_fverts, &zput));
9384e2e9504SJed Brown               for (PetscInt i = 0; i < num_fverts; i++) {
9394e2e9504SJed Brown                 Ijk loc = ZCodeSplit(vert_z[fverts[i] - vStart]);
9404e2e9504SJed Brown                 switch (bc / 2) {
9414e2e9504SJed Brown                 case 0:
9424e2e9504SJed Brown                   loc.i = 0;
9434e2e9504SJed Brown                   break;
9444e2e9504SJed Brown                 case 1:
9454e2e9504SJed Brown                   loc.j = 0;
9464e2e9504SJed Brown                   break;
9474e2e9504SJed Brown                 case 2:
9484e2e9504SJed Brown                   loc.k = 0;
9494e2e9504SJed Brown                   break;
9504e2e9504SJed Brown                 }
9514e2e9504SJed Brown                 *zput++ = ZEncode(loc);
9524e2e9504SJed Brown               }
9534e2e9504SJed Brown             }
9544e2e9504SJed Brown             continue;
9554e2e9504SJed Brown           }
9563431e603SJed Brown           PetscAssert(bc_match == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face matches multiple face sets");
9573431e603SJed Brown           PetscCall(DMLabelSetValue(label, f, face_marker_dim[dim][bc]));
9583431e603SJed Brown           bc_match++;
9593431e603SJed Brown         }
9603431e603SJed Brown       }
9613431e603SJed Brown     }
9623431e603SJed Brown     // Ensure that the Coordinate DM has our new boundary labels
9633431e603SJed Brown     DM cdm;
9643431e603SJed Brown     PetscCall(DMGetCoordinateDM(dm, &cdm));
9653431e603SJed Brown     PetscCall(DMCopyLabels(dm, cdm, PETSC_COPY_VALUES, PETSC_FALSE, DM_COPY_LABELS_FAIL));
9666725e60dSJed Brown     if (periodicity[0] == DM_BOUNDARY_PERIODIC || (dim > 1 && periodicity[1] == DM_BOUNDARY_PERIODIC) || (dim > 2 && periodicity[2] == DM_BOUNDARY_PERIODIC)) {
9675f06a3ddSJed Brown       PetscCall(DMPlexCreateBoxMesh_Tensor_SFC_Periodicity_Private(dm, &layout, vert_z, per_faces, lower, upper, periodicity, donor_face_closure, my_donor_faces));
9686725e60dSJed Brown     }
9699ae3492bSJames Wright     for (PetscInt i = 0; i < 3; i++) {
9709ae3492bSJames Wright       PetscCall(PetscSegBufferDestroy(&per_faces[i]));
9719ae3492bSJames Wright       PetscCall(PetscSegBufferDestroy(&donor_face_closure[i]));
9729ae3492bSJames Wright       PetscCall(PetscSegBufferDestroy(&my_donor_faces[i]));
9739ae3492bSJames Wright     }
9743431e603SJed Brown   }
9754e2e9504SJed Brown   PetscCall(PetscFree(layout.zstarts));
9763431e603SJed Brown   PetscCall(PetscFree(vert_z));
977ea7b3863SJames Wright   PetscCall(PetscLogEventEnd(DMPLEX_CreateBoxSFC, dm, 0, 0, 0));
9783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9793e72e933SJed Brown }
9804e2e9504SJed Brown 
9814e2e9504SJed Brown /*@
9825f06a3ddSJed Brown   DMPlexSetIsoperiodicFaceSF - Express periodicity from an existing mesh
9834e2e9504SJed Brown 
98420f4b53cSBarry Smith   Logically Collective
9854e2e9504SJed Brown 
9864e2e9504SJed Brown   Input Parameters:
9874e2e9504SJed Brown + dm           - The `DMPLEX` on which to set periodicity
9881fca310dSJames Wright . num_face_sfs - Number of `PetscSF`s in `face_sfs`
9891fca310dSJames 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.
9904e2e9504SJed Brown 
9914e2e9504SJed Brown   Level: advanced
9924e2e9504SJed Brown 
99353c0d4aeSBarry Smith   Note:
9945dca41c3SJed Brown   One can use `-dm_plex_shape zbox` to use this mode of periodicity, wherein the periodic points are distinct both globally
9954e2e9504SJed Brown   and locally, but are paired when creating a global dof space.
9964e2e9504SJed Brown 
9971cc06b55SBarry Smith .seealso: [](ch_unstructured), `DMPLEX`, `DMGetGlobalSection()`, `DMPlexGetIsoperiodicFaceSF()`
9984e2e9504SJed Brown @*/
9991fca310dSJames Wright PetscErrorCode DMPlexSetIsoperiodicFaceSF(DM dm, PetscInt num_face_sfs, PetscSF *face_sfs)
10004e2e9504SJed Brown {
10014e2e9504SJed Brown   DM_Plex *plex = (DM_Plex *)dm->data;
10024d86920dSPierre Jolivet 
10034e2e9504SJed Brown   PetscFunctionBegin;
10044e2e9504SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10051fca310dSJames Wright   if (face_sfs) PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMGetIsoperiodicPointSF_C", DMGetIsoperiodicPointSF_Plex));
10061fca310dSJames Wright   if (face_sfs == plex->periodic.face_sfs && num_face_sfs == plex->periodic.num_face_sfs) PetscFunctionReturn(PETSC_SUCCESS);
10071fca310dSJames Wright 
10081fca310dSJames Wright   for (PetscInt i = 0; i < num_face_sfs; i++) PetscCall(PetscObjectReference((PetscObject)face_sfs[i]));
10091fca310dSJames Wright 
10101fca310dSJames Wright   if (plex->periodic.num_face_sfs > 0) {
10111fca310dSJames Wright     for (PetscInt i = 0; i < plex->periodic.num_face_sfs; i++) PetscCall(PetscSFDestroy(&plex->periodic.face_sfs[i]));
10121fca310dSJames Wright     PetscCall(PetscFree(plex->periodic.face_sfs));
10131fca310dSJames Wright   }
10141fca310dSJames Wright 
10151fca310dSJames Wright   plex->periodic.num_face_sfs = num_face_sfs;
10161fca310dSJames Wright   PetscCall(PetscCalloc1(num_face_sfs, &plex->periodic.face_sfs));
10171fca310dSJames Wright   for (PetscInt i = 0; i < num_face_sfs; i++) plex->periodic.face_sfs[i] = face_sfs[i];
10185f06a3ddSJed Brown 
10195f06a3ddSJed Brown   DM cdm = dm->coordinates[0].dm; // Can't DMGetCoordinateDM because it automatically creates one
10205f06a3ddSJed Brown   if (cdm) {
10211fca310dSJames Wright     PetscCall(DMPlexSetIsoperiodicFaceSF(cdm, num_face_sfs, face_sfs));
10221fca310dSJames Wright     if (face_sfs) cdm->periodic.setup = DMPeriodicCoordinateSetUp_Internal;
10235f06a3ddSJed Brown   }
10243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10254e2e9504SJed Brown }
10264e2e9504SJed Brown 
10271fca310dSJames Wright /*@C
10285f06a3ddSJed Brown   DMPlexGetIsoperiodicFaceSF - Obtain periodicity for a mesh
10294e2e9504SJed Brown 
103020f4b53cSBarry Smith   Logically Collective
10314e2e9504SJed Brown 
103220f4b53cSBarry Smith   Input Parameter:
10334e2e9504SJed Brown . dm - The `DMPLEX` for which to obtain periodic relation
10344e2e9504SJed Brown 
10359cde84edSJose E. Roman   Output Parameters:
10361fca310dSJames Wright + num_face_sfs - Number of `PetscSF`s in the array
10371fca310dSJames 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.
10384e2e9504SJed Brown 
10394e2e9504SJed Brown   Level: advanced
10404e2e9504SJed Brown 
10411cc06b55SBarry Smith .seealso: [](ch_unstructured), `DMPLEX`, `DMGetGlobalSection()`, `DMPlexSetIsoperiodicFaceSF()`
10424e2e9504SJed Brown @*/
10431fca310dSJames Wright PetscErrorCode DMPlexGetIsoperiodicFaceSF(DM dm, PetscInt *num_face_sfs, const PetscSF **face_sfs)
10444e2e9504SJed Brown {
10454e2e9504SJed Brown   DM_Plex *plex = (DM_Plex *)dm->data;
10464d86920dSPierre Jolivet 
10474e2e9504SJed Brown   PetscFunctionBegin;
10484e2e9504SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10491fca310dSJames Wright   *face_sfs     = plex->periodic.face_sfs;
10501fca310dSJames Wright   *num_face_sfs = plex->periodic.num_face_sfs;
10513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10524e2e9504SJed Brown }
10536725e60dSJed Brown 
10545f06a3ddSJed Brown /*@C
10555f06a3ddSJed Brown   DMPlexSetIsoperiodicFaceTransform - set geometric transform from donor to periodic points
10566725e60dSJed Brown 
10576725e60dSJed Brown   Logically Collective
10586725e60dSJed Brown 
105960225df5SJacob Faibussowitsch   Input Parameters:
106053c0d4aeSBarry Smith + dm - `DMPLEX` that has been configured with `DMPlexSetIsoperiodicFaceSF()`
10611fca310dSJames Wright . n  - Number of transforms in array
10621fca310dSJames Wright - t  - Array of 4x4 affine transformation basis.
10636725e60dSJed Brown 
106453c0d4aeSBarry Smith   Level: advanced
106553c0d4aeSBarry Smith 
10665f06a3ddSJed Brown   Notes:
10675f06a3ddSJed Brown   Affine transforms are 4x4 matrices in which the leading 3x3 block expresses a rotation (or identity for no rotation),
10685f06a3ddSJed Brown   the last column contains a translation vector, and the bottom row is all zero except the last entry, which must always
10695f06a3ddSJed Brown   be 1. This representation is common in geometric modeling and allows affine transformations to be composed using
1070aaa8cc7dSPierre Jolivet   simple matrix multiplication.
10715f06a3ddSJed Brown 
10725f06a3ddSJed Brown   Although the interface accepts a general affine transform, only affine translation is supported at present.
10735f06a3ddSJed Brown 
107460225df5SJacob Faibussowitsch   Developer Notes:
10755f06a3ddSJed Brown   This interface should be replaced by making BasisTransform public, expanding it to support affine representations, and
10765f06a3ddSJed Brown   adding GPU implementations to apply the G2L/L2G transforms.
107753c0d4aeSBarry Smith 
10781cc06b55SBarry Smith .seealso: [](ch_unstructured), `DMPLEX`, `DMGetGlobalSection()`, `DMPlexSetIsoperiodicFaceSF()`
10796725e60dSJed Brown @*/
1080cc4c1da9SBarry Smith PetscErrorCode DMPlexSetIsoperiodicFaceTransform(DM dm, PetscInt n, const PetscScalar t[])
10816725e60dSJed Brown {
10826725e60dSJed Brown   DM_Plex *plex = (DM_Plex *)dm->data;
10834d86920dSPierre Jolivet 
10846725e60dSJed Brown   PetscFunctionBegin;
10856725e60dSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10861fca310dSJames 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);
10871fca310dSJames Wright 
10881fca310dSJames Wright   PetscCall(PetscMalloc1(n, &plex->periodic.transform));
10891fca310dSJames Wright   for (PetscInt i = 0; i < n; i++) {
10906725e60dSJed Brown     for (PetscInt j = 0; j < 4; j++) {
10911fca310dSJames Wright       for (PetscInt k = 0; k < 4; k++) {
10921fca310dSJames Wright         PetscCheck(j != k || t[i * 16 + j * 4 + k] == 1., PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Rotated transforms not supported");
10931fca310dSJames Wright         plex->periodic.transform[i][j][k] = t[i * 16 + j * 4 + k];
10941fca310dSJames Wright       }
10956725e60dSJed Brown     }
10966725e60dSJed Brown   }
10973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10986725e60dSJed Brown }
1099