xref: /petsc/src/dm/impls/plex/plexsfc.c (revision 53b735f80d8abba57081ee21f195942af5a4ca31)
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]
30bc3e2e66SJames Wright // Magic numbers taken from https://stackoverflow.com/a/18528775/7564988 (translated to octal)
313e72e933SJed Brown static unsigned ZCodeSplit1(ZCode z)
323e72e933SJed Brown {
33bc3e2e66SJames Wright   z &= 0111111111111111111111;
34bc3e2e66SJames Wright   z = (z | z >> 2) & 0103030303030303030303;
35bc3e2e66SJames Wright   z = (z | z >> 4) & 0100170017001700170017;
36bc3e2e66SJames Wright   z = (z | z >> 8) & 0370000037700000377;
37bc3e2e66SJames Wright   z = (z | z >> 16) & 0370000000000177777;
38bc3e2e66SJames 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;
47bc3e2e66SJames Wright   z &= 07777777;
48bc3e2e66SJames Wright   z = (z | z << 32) & 0370000000000177777;
49bc3e2e66SJames Wright   z = (z | z << 16) & 0370000037700000377;
50bc3e2e66SJames Wright   z = (z | z << 8) & 0100170017001700170017;
51bc3e2e66SJames Wright   z = (z | z << 4) & 0103030303030303030303;
52bc3e2e66SJames 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 
20585de0153SJames Wright static inline PetscBool IsPointInsideStratum(PetscInt point, PetscInt pStart, PetscInt pEnd)
20685de0153SJames Wright {
20785de0153SJames Wright   return (point >= pStart && point < pEnd) ? PETSC_TRUE : PETSC_FALSE;
20885de0153SJames Wright }
20985de0153SJames Wright 
2109ae3492bSJames 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])
2114e2e9504SJed Brown {
2124e2e9504SJed Brown   MPI_Comm    comm;
2139ae3492bSJames Wright   PetscInt    dim, vStart, vEnd;
2144e2e9504SJed Brown   PetscMPIInt size;
2159ae3492bSJames Wright   PetscSF     face_sfs[3];
2169ae3492bSJames Wright   PetscScalar transforms[3][4][4] = {{{0}}};
2174e2e9504SJed Brown 
2184e2e9504SJed Brown   PetscFunctionBegin;
2194e2e9504SJed Brown   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2204e2e9504SJed Brown   PetscCallMPI(MPI_Comm_size(comm, &size));
2214e2e9504SJed Brown   PetscCall(DMGetDimension(dm, &dim));
2224e2e9504SJed Brown   const PetscInt csize = PetscPowInt(2, dim - 1);
2234e2e9504SJed Brown   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
2249ae3492bSJames Wright 
2259ae3492bSJames Wright   PetscInt num_directions = 0;
2269ae3492bSJames Wright   for (PetscInt direction = 0; direction < dim; direction++) {
2276497c311SBarry Smith     PetscCount   num_faces;
2289ae3492bSJames Wright     PetscInt    *faces;
2299ae3492bSJames Wright     ZCode       *donor_verts, *donor_minz;
2309ae3492bSJames Wright     PetscSFNode *leaf;
2316497c311SBarry Smith     PetscCount   num_multiroots = 0;
2326497c311SBarry Smith     PetscInt     pStart, pEnd;
2336497c311SBarry Smith     PetscBool    sorted;
2346497c311SBarry Smith     PetscInt     inum_faces;
2359ae3492bSJames Wright 
2369ae3492bSJames Wright     if (periodicity[direction] != DM_BOUNDARY_PERIODIC) continue;
2379ae3492bSJames Wright     PetscCall(PetscSegBufferGetSize(per_faces[direction], &num_faces));
2389ae3492bSJames Wright     PetscCall(PetscSegBufferExtractInPlace(per_faces[direction], &faces));
2399ae3492bSJames Wright     PetscCall(PetscSegBufferExtractInPlace(donor_face_closure[direction], &donor_verts));
2404e2e9504SJed Brown     PetscCall(PetscMalloc1(num_faces, &donor_minz));
2414e2e9504SJed Brown     PetscCall(PetscMalloc1(num_faces, &leaf));
2426497c311SBarry Smith     for (PetscCount i = 0; i < num_faces; i++) {
2434e2e9504SJed Brown       ZCode minz = donor_verts[i * csize];
2446497c311SBarry Smith 
2454e2e9504SJed Brown       for (PetscInt j = 1; j < csize; j++) minz = PetscMin(minz, donor_verts[i * csize + j]);
2464e2e9504SJed Brown       donor_minz[i] = minz;
2474e2e9504SJed Brown     }
2486497c311SBarry Smith     PetscCall(PetscIntCast(num_faces, &inum_faces));
2496497c311SBarry Smith     PetscCall(PetscSortedInt64(inum_faces, (const PetscInt64 *)donor_minz, &sorted));
2509ae3492bSJames Wright     // If a donor vertex were chosen to broker multiple faces, we would have a logic error.
2519ae3492bSJames Wright     // Checking for sorting is a cheap check that there are no duplicates.
2529ae3492bSJames Wright     PetscCheck(sorted, PETSC_COMM_SELF, PETSC_ERR_PLIB, "minz not sorted; possible duplicates not checked");
2536497c311SBarry Smith     for (PetscCount i = 0; i < num_faces;) {
2544e2e9504SJed Brown       ZCode       z = donor_minz[i];
255835f2295SStefano Zampini       PetscMPIInt remote_rank, remote_count = 0;
2566497c311SBarry Smith 
257835f2295SStefano Zampini       PetscCall(PetscMPIIntCast(ZCodeFind(z, size + 1, layout->zstarts), &remote_rank));
2584e2e9504SJed Brown       if (remote_rank < 0) remote_rank = -(remote_rank + 1) - 1;
2594e2e9504SJed Brown       // Process all the vertices on this rank
2604e2e9504SJed Brown       for (ZCode rz = layout->zstarts[remote_rank]; rz < layout->zstarts[remote_rank + 1]; rz++) {
2614e2e9504SJed Brown         Ijk loc = ZCodeSplit(rz);
2626497c311SBarry Smith 
2634e2e9504SJed Brown         if (rz == z) {
2644e2e9504SJed Brown           leaf[i].rank  = remote_rank;
2654e2e9504SJed Brown           leaf[i].index = remote_count;
2664e2e9504SJed Brown           i++;
2676497c311SBarry Smith           if (i == num_faces) break;
2684e2e9504SJed Brown           z = donor_minz[i];
2694e2e9504SJed Brown         }
2704e2e9504SJed Brown         if (IjkActive(layout->vextent, loc)) remote_count++;
2714e2e9504SJed Brown       }
2724e2e9504SJed Brown     }
2734e2e9504SJed Brown     PetscCall(PetscFree(donor_minz));
2749ae3492bSJames Wright     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), &face_sfs[num_directions]));
2756497c311SBarry Smith     PetscCall(PetscSFSetGraph(face_sfs[num_directions], vEnd - vStart, inum_faces, NULL, PETSC_USE_POINTER, leaf, PETSC_USE_POINTER));
2764e2e9504SJed Brown     const PetscInt *my_donor_degree;
2779ae3492bSJames Wright     PetscCall(PetscSFComputeDegreeBegin(face_sfs[num_directions], &my_donor_degree));
2789ae3492bSJames Wright     PetscCall(PetscSFComputeDegreeEnd(face_sfs[num_directions], &my_donor_degree));
2796497c311SBarry Smith 
2804e2e9504SJed Brown     for (PetscInt i = 0; i < vEnd - vStart; i++) {
2814e2e9504SJed Brown       num_multiroots += my_donor_degree[i];
2824e2e9504SJed Brown       if (my_donor_degree[i] == 0) continue;
2835f06a3ddSJed Brown       PetscAssert(my_donor_degree[i] == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Local vertex has multiple faces");
2844e2e9504SJed Brown     }
2854e2e9504SJed Brown     PetscInt  *my_donors, *donor_indices, *my_donor_indices;
2866497c311SBarry Smith     PetscCount num_my_donors;
2876497c311SBarry Smith 
2889ae3492bSJames Wright     PetscCall(PetscSegBufferGetSize(my_donor_faces[direction], &num_my_donors));
289563af49cSJames 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);
2909ae3492bSJames Wright     PetscCall(PetscSegBufferExtractInPlace(my_donor_faces[direction], &my_donors));
2914e2e9504SJed Brown     PetscCall(PetscMalloc1(vEnd - vStart, &my_donor_indices));
2926497c311SBarry Smith     for (PetscCount i = 0; i < num_my_donors; i++) {
2934e2e9504SJed Brown       PetscInt f = my_donors[i];
2941690c2aeSBarry Smith       PetscInt num_points, *points = NULL, minv = PETSC_INT_MAX;
2956497c311SBarry Smith 
2964e2e9504SJed Brown       PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &num_points, &points));
2974e2e9504SJed Brown       for (PetscInt j = 0; j < num_points; j++) {
2984e2e9504SJed Brown         PetscInt p = points[2 * j];
29985de0153SJames Wright         if (!IsPointInsideStratum(p, vStart, vEnd)) continue;
3004e2e9504SJed Brown         minv = PetscMin(minv, p);
3014e2e9504SJed Brown       }
3024e2e9504SJed Brown       PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &num_points, &points));
3035f06a3ddSJed Brown       PetscAssert(my_donor_degree[minv - vStart] == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Local vertex not requested");
3044e2e9504SJed Brown       my_donor_indices[minv - vStart] = f;
3054e2e9504SJed Brown     }
3064e2e9504SJed Brown     PetscCall(PetscMalloc1(num_faces, &donor_indices));
3079ae3492bSJames Wright     PetscCall(PetscSFBcastBegin(face_sfs[num_directions], MPIU_INT, my_donor_indices, donor_indices, MPI_REPLACE));
3089ae3492bSJames Wright     PetscCall(PetscSFBcastEnd(face_sfs[num_directions], MPIU_INT, my_donor_indices, donor_indices, MPI_REPLACE));
3094e2e9504SJed Brown     PetscCall(PetscFree(my_donor_indices));
3104e2e9504SJed Brown     // Modify our leafs so they point to donor faces instead of donor minz. Additionally, give them indices as faces.
3116497c311SBarry Smith     for (PetscCount i = 0; i < num_faces; i++) leaf[i].index = donor_indices[i];
3124e2e9504SJed Brown     PetscCall(PetscFree(donor_indices));
3134e2e9504SJed Brown     PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
3146497c311SBarry Smith     PetscCall(PetscSFSetGraph(face_sfs[num_directions], pEnd - pStart, inum_faces, faces, PETSC_COPY_VALUES, leaf, PETSC_OWN_POINTER));
3159ae3492bSJames Wright     {
3169ae3492bSJames Wright       char face_sf_name[PETSC_MAX_PATH_LEN];
3179ae3492bSJames Wright       PetscCall(PetscSNPrintf(face_sf_name, sizeof face_sf_name, "Z-order Isoperiodic Faces #%" PetscInt_FMT, num_directions));
3189ae3492bSJames Wright       PetscCall(PetscObjectSetName((PetscObject)face_sfs[num_directions], face_sf_name));
3199ae3492bSJames Wright     }
3204e2e9504SJed Brown 
3219ae3492bSJames Wright     transforms[num_directions][0][0]         = 1;
3229ae3492bSJames Wright     transforms[num_directions][1][1]         = 1;
3239ae3492bSJames Wright     transforms[num_directions][2][2]         = 1;
3249ae3492bSJames Wright     transforms[num_directions][3][3]         = 1;
3259ae3492bSJames Wright     transforms[num_directions][direction][3] = upper[direction] - lower[direction];
3269ae3492bSJames Wright     num_directions++;
3279ae3492bSJames Wright   }
3286725e60dSJed Brown 
3291fca310dSJames Wright   PetscCall(DMPlexSetIsoperiodicFaceSF(dm, num_directions, face_sfs));
3301fca310dSJames Wright   PetscCall(DMPlexSetIsoperiodicFaceTransform(dm, num_directions, (PetscScalar *)transforms));
3319ae3492bSJames Wright 
3329ae3492bSJames Wright   for (PetscInt i = 0; i < num_directions; i++) PetscCall(PetscSFDestroy(&face_sfs[i]));
3333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3344e2e9504SJed Brown }
3354e2e9504SJed Brown 
3365f06a3ddSJed Brown // This is a DMGlobalToLocalHook that applies the affine offsets. When extended for rotated periodicity, it'll need to
3375f06a3ddSJed Brown // apply a rotatonal transform and similar operations will be needed for fields (e.g., to rotate a velocity vector).
3385f06a3ddSJed Brown // We use this crude approach here so we don't have to write new GPU kernels yet.
3396725e60dSJed Brown static PetscErrorCode DMCoordAddPeriodicOffsets_Private(DM dm, Vec g, InsertMode mode, Vec l, void *ctx)
3406725e60dSJed Brown {
3416725e60dSJed Brown   PetscFunctionBegin;
342ddedc8f6SJames Wright   // These `VecScatter`s should be merged to improve efficiency; the scatters cannot be overlapped.
343ddedc8f6SJames Wright   for (PetscInt i = 0; i < dm->periodic.num_affines; i++) {
344ddedc8f6SJames Wright     PetscCall(VecScatterBegin(dm->periodic.affine_to_local[i], dm->periodic.affine[i], l, ADD_VALUES, SCATTER_FORWARD));
345ddedc8f6SJames Wright     PetscCall(VecScatterEnd(dm->periodic.affine_to_local[i], dm->periodic.affine[i], l, ADD_VALUES, SCATTER_FORWARD));
346ddedc8f6SJames Wright   }
3473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3486725e60dSJed Brown }
3496725e60dSJed Brown 
350*53b735f8SJames Wright // Modify Vec based on the transformation of `point` for the given section and field
351*53b735f8SJames Wright static PetscErrorCode DMPlexOrientFieldPointVec(DM dm, PetscSection section, PetscInt field, Vec V, PetscInt point, PetscInt orientation)
352*53b735f8SJames Wright {
353*53b735f8SJames Wright   PetscScalar        *copy, *V_arr;
354*53b735f8SJames Wright   PetscInt            dof, off, point_ornt[2] = {point, orientation};
355*53b735f8SJames Wright   const PetscInt    **perms;
356*53b735f8SJames Wright   const PetscScalar **rots;
357*53b735f8SJames Wright 
358*53b735f8SJames Wright   PetscFunctionBeginUser;
359*53b735f8SJames Wright   PetscCall(PetscSectionGetDof(section, point, &dof));
360*53b735f8SJames Wright   PetscCall(PetscSectionGetOffset(section, point, &off));
361*53b735f8SJames Wright   PetscCall(VecGetArray(V, &V_arr));
362*53b735f8SJames Wright   PetscCall(DMGetWorkArray(dm, dof, MPIU_SCALAR, &copy));
363*53b735f8SJames Wright   PetscArraycpy(copy, &V_arr[off], dof);
364*53b735f8SJames Wright 
365*53b735f8SJames Wright   PetscCall(PetscSectionGetFieldPointSyms(section, field, 1, point_ornt, &perms, &rots));
366*53b735f8SJames Wright   for (PetscInt i = 0; i < dof; i++) {
367*53b735f8SJames Wright     if (perms[0]) V_arr[off + perms[0][i]] = copy[i];
368*53b735f8SJames Wright     if (rots[0]) V_arr[off + perms[0][i]] *= rots[0][i];
369*53b735f8SJames Wright   }
370*53b735f8SJames Wright 
371*53b735f8SJames Wright   PetscCall(PetscSectionRestoreFieldPointSyms(section, field, 1, point_ornt, &perms, &rots));
372*53b735f8SJames Wright   PetscCall(DMRestoreWorkArray(dm, dof, MPIU_SCALAR, &copy));
373*53b735f8SJames Wright   PetscCall(VecRestoreArray(V, &V_arr));
374*53b735f8SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
375*53b735f8SJames Wright }
376*53b735f8SJames Wright 
377*53b735f8SJames Wright // Reorient the point in the DMPlex while also applying necessary corrections to other structures (e.g. coordinates)
378*53b735f8SJames Wright static PetscErrorCode DMPlexOrientPointWithCorrections(DM dm, PetscInt point, PetscInt ornt)
379*53b735f8SJames Wright {
380*53b735f8SJames Wright   // TODO: Potential speed up if we early exit for ornt == 0 (i.e. if ornt is identity, we don't need to do anything)
381*53b735f8SJames Wright   PetscFunctionBeginUser;
382*53b735f8SJames Wright   PetscCall(DMPlexOrientPoint(dm, point, ornt));
383*53b735f8SJames Wright 
384*53b735f8SJames Wright   { // Correct coordinates based on new cone ordering
385*53b735f8SJames Wright     DM           cdm;
386*53b735f8SJames Wright     PetscSection csection;
387*53b735f8SJames Wright     Vec          coordinates;
388*53b735f8SJames Wright     PetscInt     pStart, pEnd;
389*53b735f8SJames Wright 
390*53b735f8SJames Wright     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
391*53b735f8SJames Wright     PetscCall(DMGetCoordinateDM(dm, &cdm));
392*53b735f8SJames Wright     PetscCall(DMGetLocalSection(cdm, &csection));
393*53b735f8SJames Wright     PetscCall(PetscSectionGetChart(csection, &pStart, &pEnd));
394*53b735f8SJames Wright     if (IsPointInsideStratum(point, pStart, pEnd)) PetscCall(DMPlexOrientFieldPointVec(cdm, csection, 0, coordinates, point, ornt));
395*53b735f8SJames Wright   }
396*53b735f8SJames Wright   // TODO: Correct sfNatural
397*53b735f8SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
398*53b735f8SJames Wright }
399*53b735f8SJames Wright 
40093defd67SJames Wright // Creates SF to communicate data from donor to periodic faces. The data can be different sizes per donor-periodic pair and is given in `point_sizes[]`
401*53b735f8SJames Wright static PetscErrorCode CreateDonorToPeriodicSF(DM dm, PetscSF face_sf, PetscInt pStart, PetscInt pEnd, const PetscInt point_sizes[], PetscInt *rootbuffersize, PetscInt *leafbuffersize, PetscBT *rootbt, PetscSF *sf_closure)
40293defd67SJames Wright {
40393defd67SJames Wright   MPI_Comm           comm;
40493defd67SJames Wright   PetscMPIInt        rank;
40593defd67SJames Wright   PetscInt           nroots, nleaves;
40693defd67SJames Wright   PetscInt          *rootdata, *leafdata;
40793defd67SJames Wright   const PetscInt    *filocal;
40893defd67SJames Wright   const PetscSFNode *firemote;
40993defd67SJames Wright 
41093defd67SJames Wright   PetscFunctionBeginUser;
41193defd67SJames Wright   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
41293defd67SJames Wright   PetscCallMPI(MPI_Comm_rank(comm, &rank));
41393defd67SJames Wright 
41493defd67SJames Wright   PetscCall(PetscSFGetGraph(face_sf, &nroots, &nleaves, &filocal, &firemote));
41593defd67SJames Wright   PetscCall(PetscCalloc2(2 * nroots, &rootdata, 2 * nroots, &leafdata));
41693defd67SJames Wright   for (PetscInt i = 0; i < nleaves; i++) {
41793defd67SJames Wright     PetscInt point = filocal[i];
41885de0153SJames Wright     PetscCheck(IsPointInsideStratum(point, pStart, pEnd), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %" PetscInt_FMT " in leaves exists outside of stratum [%" PetscInt_FMT ", %" PetscInt_FMT ")", point, pStart, pEnd);
41993defd67SJames Wright     leafdata[point] = point_sizes[point - pStart];
42093defd67SJames Wright   }
42193defd67SJames Wright   PetscCall(PetscSFReduceBegin(face_sf, MPIU_INT, leafdata, rootdata + nroots, MPIU_SUM));
42293defd67SJames Wright   PetscCall(PetscSFReduceEnd(face_sf, MPIU_INT, leafdata, rootdata + nroots, MPIU_SUM));
42393defd67SJames Wright 
42493defd67SJames Wright   PetscInt root_offset = 0;
42593defd67SJames Wright   PetscCall(PetscBTCreate(nroots, rootbt));
42693defd67SJames Wright   for (PetscInt p = 0; p < nroots; p++) {
42793defd67SJames Wright     const PetscInt *donor_dof = rootdata + nroots;
42893defd67SJames Wright     if (donor_dof[p] == 0) {
42993defd67SJames Wright       rootdata[2 * p]     = -1;
43093defd67SJames Wright       rootdata[2 * p + 1] = -1;
43193defd67SJames Wright       continue;
43293defd67SJames Wright     }
43393defd67SJames Wright     PetscCall(PetscBTSet(*rootbt, p));
43485de0153SJames Wright     PetscCheck(IsPointInsideStratum(p, pStart, pEnd), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %" PetscInt_FMT " in roots exists outside of stratum [%" PetscInt_FMT ", %" PetscInt_FMT ")", p, pStart, pEnd);
43593defd67SJames Wright     PetscInt p_size = point_sizes[p - pStart];
43693defd67SJames Wright     PetscCheck(donor_dof[p] == p_size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Reduced leaf data size (%" PetscInt_FMT ") does not match root data size (%" PetscInt_FMT ")", donor_dof[p], p_size);
43793defd67SJames Wright     rootdata[2 * p]     = root_offset;
43893defd67SJames Wright     rootdata[2 * p + 1] = p_size;
43993defd67SJames Wright     root_offset += p_size;
44093defd67SJames Wright   }
44193defd67SJames Wright   PetscCall(PetscSFBcastBegin(face_sf, MPIU_2INT, rootdata, leafdata, MPI_REPLACE));
44293defd67SJames Wright   PetscCall(PetscSFBcastEnd(face_sf, MPIU_2INT, rootdata, leafdata, MPI_REPLACE));
44393defd67SJames Wright   // Count how many leaves we need to communicate the closures
44493defd67SJames Wright   PetscInt leaf_offset = 0;
44593defd67SJames Wright   for (PetscInt i = 0; i < nleaves; i++) {
44693defd67SJames Wright     PetscInt point = filocal[i];
44793defd67SJames Wright     if (leafdata[2 * point + 1] < 0) continue;
44893defd67SJames Wright     leaf_offset += leafdata[2 * point + 1];
44993defd67SJames Wright   }
45093defd67SJames Wright 
45193defd67SJames Wright   PetscSFNode *closure_leaf;
45293defd67SJames Wright   PetscCall(PetscMalloc1(leaf_offset, &closure_leaf));
45393defd67SJames Wright   leaf_offset = 0;
45493defd67SJames Wright   for (PetscInt i = 0; i < nleaves; i++) {
45593defd67SJames Wright     PetscInt point   = filocal[i];
45693defd67SJames Wright     PetscInt cl_size = leafdata[2 * point + 1];
45793defd67SJames Wright     if (cl_size < 0) continue;
45893defd67SJames Wright     for (PetscInt j = 0; j < cl_size; j++) {
45993defd67SJames Wright       closure_leaf[leaf_offset].rank  = firemote[i].rank;
46093defd67SJames Wright       closure_leaf[leaf_offset].index = leafdata[2 * point] + j;
46193defd67SJames Wright       leaf_offset++;
46293defd67SJames Wright     }
46393defd67SJames Wright   }
46493defd67SJames Wright 
46593defd67SJames Wright   PetscCall(PetscSFCreate(comm, sf_closure));
46693defd67SJames Wright   PetscCall(PetscSFSetGraph(*sf_closure, root_offset, leaf_offset, NULL, PETSC_USE_POINTER, closure_leaf, PETSC_OWN_POINTER));
46793defd67SJames Wright   *rootbuffersize = root_offset;
46893defd67SJames Wright   *leafbuffersize = leaf_offset;
46993defd67SJames Wright   PetscCall(PetscFree2(rootdata, leafdata));
47093defd67SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
47193defd67SJames Wright }
47293defd67SJames Wright 
473*53b735f8SJames Wright // Determine if `key` is in `array`. `array` does NOT need to be sorted
474*53b735f8SJames Wright static inline PetscBool SearchIntArray(PetscInt key, PetscInt array_size, const PetscInt array[])
475*53b735f8SJames Wright {
476*53b735f8SJames Wright   for (PetscInt i = 0; i < array_size; i++)
477*53b735f8SJames Wright     if (array[i] == key) return PETSC_TRUE;
478*53b735f8SJames Wright   return PETSC_FALSE;
479*53b735f8SJames Wright }
480*53b735f8SJames Wright 
481*53b735f8SJames Wright // Translate a cone in periodic points to the cone in donor points based on the `periodic2donor` array
482*53b735f8SJames Wright static inline PetscErrorCode TranslateConeP2D(const PetscInt periodic_cone[], PetscInt cone_size, const PetscInt periodic2donor[], PetscInt p2d_count, PetscInt p2d_cone[])
483*53b735f8SJames Wright {
484*53b735f8SJames Wright   PetscFunctionBeginUser;
485*53b735f8SJames Wright   for (PetscInt p = 0; p < cone_size; p++) {
486*53b735f8SJames Wright     PetscInt p2d_index = -1;
487*53b735f8SJames Wright     for (PetscInt p2d = 0; p2d < p2d_count; p2d++) {
488*53b735f8SJames Wright       if (periodic2donor[p2d * 2] == periodic_cone[p]) p2d_index = p2d;
489*53b735f8SJames Wright     }
490*53b735f8SJames Wright     PetscCheck(p2d_index >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not find periodic point in periodic-to-donor array");
491*53b735f8SJames Wright     p2d_cone[p] = periodic2donor[2 * p2d_index + 1];
492*53b735f8SJames Wright   }
493*53b735f8SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
494*53b735f8SJames Wright }
495*53b735f8SJames Wright 
496*53b735f8SJames Wright // Corrects the cone order of periodic faces (and their transitive closure's cones) to match their donor face pair.
497*53b735f8SJames Wright //
498*53b735f8SJames Wright // This is done by:
499*53b735f8SJames Wright // 1. Communicating the donor's vertex coordinates and recursive cones (i.e. its own cone and those of it's constituent edges) to it's periodic pairs
500*53b735f8SJames Wright //    - The donor vertices have the isoperiodic transform applied to them such that they should match exactly
501*53b735f8SJames Wright // 2. Translating the periodic vertices into the donor vertices point IDs
502*53b735f8SJames Wright // 3. Translating the cone of each periodic point into the donor point IDs
503*53b735f8SJames Wright // 4. Comparing the periodic-to-donor cone to the donor cone for each point
504*53b735f8SJames Wright // 5. Apply the necessary transformation to the periodic cone to make it match the donor cone
505*53b735f8SJames Wright static PetscErrorCode DMPlexCorrectOrientationForIsoperiodic(DM dm)
506*53b735f8SJames Wright {
507*53b735f8SJames Wright   MPI_Comm        comm;
508*53b735f8SJames Wright   DM_Plex        *plex = (DM_Plex *)dm->data;
509*53b735f8SJames Wright   PetscInt        nroots, nleaves;
510*53b735f8SJames Wright   const PetscInt *filocal;
511*53b735f8SJames Wright   DM              cdm;
512*53b735f8SJames Wright   PetscSection    csection;
513*53b735f8SJames Wright   Vec             coordinates;
514*53b735f8SJames Wright   PetscInt        coords_field_id = 0;
515*53b735f8SJames Wright   PetscBool       debug_printing  = PETSC_FALSE;
516*53b735f8SJames Wright 
517*53b735f8SJames Wright   PetscFunctionBeginUser;
518*53b735f8SJames Wright   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
519*53b735f8SJames Wright   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
520*53b735f8SJames Wright   PetscCheck(coordinates, comm, PETSC_ERR_ARG_WRONGSTATE, "DM must have coordinates to setup isoperiodic");
521*53b735f8SJames Wright   PetscCall(DMGetCoordinateDM(dm, &cdm));
522*53b735f8SJames Wright   PetscCall(DMGetLocalSection(cdm, &csection));
523*53b735f8SJames Wright 
524*53b735f8SJames Wright   for (PetscInt f = 0; f < plex->periodic.num_face_sfs; f++) {
525*53b735f8SJames Wright     PetscSF face_sf                  = plex->periodic.face_sfs[f];
526*53b735f8SJames Wright     const PetscScalar(*transform)[4] = (const PetscScalar(*)[4])plex->periodic.transform[f];
527*53b735f8SJames Wright     PetscInt *face_vertices_size, *face_cones_size;
528*53b735f8SJames Wright     PetscInt  fStart, fEnd, vStart, vEnd, rootnumvert, leafnumvert, rootconesize, leafconesize, dim;
529*53b735f8SJames Wright     PetscSF   sf_vert_coords, sf_face_cones;
530*53b735f8SJames Wright     PetscBT   rootbt;
531*53b735f8SJames Wright 
532*53b735f8SJames Wright     PetscCall(DMGetCoordinateDim(dm, &dim));
533*53b735f8SJames Wright     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
534*53b735f8SJames Wright     PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
535*53b735f8SJames Wright     PetscCall(PetscCalloc2(fEnd - fStart, &face_vertices_size, fEnd - fStart, &face_cones_size));
536*53b735f8SJames Wright 
537*53b735f8SJames Wright     // Create SFs to communicate donor vertices and donor cones to periodic faces
538*53b735f8SJames Wright     for (PetscInt f = fStart, index = 0; f < fEnd; f++, index++) {
539*53b735f8SJames Wright       PetscInt cl_size, *closure = NULL, num_vertices = 0;
540*53b735f8SJames Wright       PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &cl_size, &closure));
541*53b735f8SJames Wright       for (PetscInt p = 0; p < cl_size; p++) {
542*53b735f8SJames Wright         PetscInt cl_point = closure[2 * p];
543*53b735f8SJames Wright         if (IsPointInsideStratum(cl_point, vStart, vEnd)) num_vertices++;
544*53b735f8SJames Wright         else {
545*53b735f8SJames Wright           PetscInt cone_size;
546*53b735f8SJames Wright           PetscCall(DMPlexGetConeSize(dm, cl_point, &cone_size));
547*53b735f8SJames Wright           face_cones_size[index] += cone_size + 2;
548*53b735f8SJames Wright         }
549*53b735f8SJames Wright       }
550*53b735f8SJames Wright       face_vertices_size[index] = num_vertices;
551*53b735f8SJames Wright       face_cones_size[index] += num_vertices;
552*53b735f8SJames Wright       PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &cl_size, &closure));
553*53b735f8SJames Wright     }
554*53b735f8SJames Wright     PetscCall(CreateDonorToPeriodicSF(dm, face_sf, fStart, fEnd, face_vertices_size, &rootnumvert, &leafnumvert, &rootbt, &sf_vert_coords));
555*53b735f8SJames Wright     PetscCall(PetscBTDestroy(&rootbt));
556*53b735f8SJames Wright     PetscCall(CreateDonorToPeriodicSF(dm, face_sf, fStart, fEnd, face_cones_size, &rootconesize, &leafconesize, &rootbt, &sf_face_cones));
557*53b735f8SJames Wright 
558*53b735f8SJames Wright     PetscCall(PetscSFGetGraph(face_sf, &nroots, &nleaves, &filocal, NULL));
559*53b735f8SJames Wright 
560*53b735f8SJames Wright     PetscReal *leaf_donor_coords;
561*53b735f8SJames Wright     PetscInt  *leaf_donor_cones;
562*53b735f8SJames Wright 
563*53b735f8SJames Wright     { // Communicate donor coords and cones to the periodic faces
564*53b735f8SJames Wright       PetscReal         *mydonor_vertices;
565*53b735f8SJames Wright       PetscInt          *mydonor_cones;
566*53b735f8SJames Wright       const PetscScalar *coords_arr;
567*53b735f8SJames Wright 
568*53b735f8SJames Wright       PetscCall(PetscCalloc2(rootnumvert * dim, &mydonor_vertices, rootconesize, &mydonor_cones));
569*53b735f8SJames Wright       PetscCall(VecGetArrayRead(coordinates, &coords_arr));
570*53b735f8SJames Wright       for (PetscInt donor_face = 0, donor_vert_offset = 0, donor_cone_offset = 0; donor_face < nroots; donor_face++) {
571*53b735f8SJames Wright         if (!PetscBTLookup(rootbt, donor_face)) continue;
572*53b735f8SJames Wright         PetscInt cl_size, *closure = NULL;
573*53b735f8SJames Wright 
574*53b735f8SJames Wright         PetscCall(DMPlexGetTransitiveClosure(dm, donor_face, PETSC_TRUE, &cl_size, &closure));
575*53b735f8SJames Wright         // Pack vertex coordinates
576*53b735f8SJames Wright         for (PetscInt p = 0; p < cl_size; p++) {
577*53b735f8SJames Wright           PetscInt cl_point = closure[2 * p], dof, offset;
578*53b735f8SJames Wright           if (!IsPointInsideStratum(cl_point, vStart, vEnd)) continue;
579*53b735f8SJames Wright           mydonor_cones[donor_cone_offset++] = cl_point;
580*53b735f8SJames Wright           PetscCall(PetscSectionGetFieldDof(csection, cl_point, coords_field_id, &dof));
581*53b735f8SJames Wright           PetscCall(PetscSectionGetFieldOffset(csection, cl_point, coords_field_id, &offset));
582*53b735f8SJames Wright           PetscAssert(dof == dim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %" PetscInt_FMT " has dof size %" PetscInt_FMT ", but should match dimension size %" PetscInt_FMT, cl_point, dof, dim);
583*53b735f8SJames Wright           // Apply isoperiodic transform to donor vertices such that corresponding periodic vertices should match exactly
584*53b735f8SJames Wright           for (PetscInt d = 0; d < dof; d++) mydonor_vertices[donor_vert_offset * dim + d] = PetscRealPart(coords_arr[offset + d]) + PetscRealPart(transform[d][3]);
585*53b735f8SJames Wright           donor_vert_offset++;
586*53b735f8SJames Wright         }
587*53b735f8SJames Wright         // Pack cones of face points (including face itself)
588*53b735f8SJames Wright         for (PetscInt p = 0; p < cl_size; p++) {
589*53b735f8SJames Wright           PetscInt        cl_point = closure[2 * p], cone_size, depth;
590*53b735f8SJames Wright           const PetscInt *cone;
591*53b735f8SJames Wright 
592*53b735f8SJames Wright           PetscCall(DMPlexGetConeSize(dm, cl_point, &cone_size));
593*53b735f8SJames Wright           PetscCall(DMPlexGetCone(dm, cl_point, &cone));
594*53b735f8SJames Wright           PetscCall(DMPlexGetPointDepth(dm, cl_point, &depth));
595*53b735f8SJames Wright           if (depth == 0) continue; // don't include vertex depth
596*53b735f8SJames Wright           mydonor_cones[donor_cone_offset++] = cone_size;
597*53b735f8SJames Wright           mydonor_cones[donor_cone_offset++] = cl_point;
598*53b735f8SJames Wright           PetscArraycpy(&mydonor_cones[donor_cone_offset], cone, cone_size);
599*53b735f8SJames Wright           donor_cone_offset += cone_size;
600*53b735f8SJames Wright         }
601*53b735f8SJames Wright         PetscCall(DMPlexRestoreTransitiveClosure(dm, donor_face, PETSC_TRUE, &cl_size, &closure));
602*53b735f8SJames Wright       }
603*53b735f8SJames Wright       PetscCall(VecRestoreArrayRead(coordinates, &coords_arr));
604*53b735f8SJames Wright       PetscCall(PetscBTDestroy(&rootbt));
605*53b735f8SJames Wright 
606*53b735f8SJames Wright       MPI_Datatype vertex_unit;
607*53b735f8SJames Wright       PetscCallMPI(MPI_Type_contiguous(dim, MPIU_REAL, &vertex_unit));
608*53b735f8SJames Wright       PetscCallMPI(MPI_Type_commit(&vertex_unit));
609*53b735f8SJames Wright       PetscCall(PetscMalloc2(leafnumvert * 3, &leaf_donor_coords, leafconesize, &leaf_donor_cones));
610*53b735f8SJames Wright       PetscCall(PetscSFBcastBegin(sf_vert_coords, vertex_unit, mydonor_vertices, leaf_donor_coords, MPI_REPLACE));
611*53b735f8SJames Wright       PetscCall(PetscSFBcastBegin(sf_face_cones, MPIU_INT, mydonor_cones, leaf_donor_cones, MPI_REPLACE));
612*53b735f8SJames Wright       PetscCall(PetscSFBcastEnd(sf_vert_coords, vertex_unit, mydonor_vertices, leaf_donor_coords, MPI_REPLACE));
613*53b735f8SJames Wright       PetscCall(PetscSFBcastEnd(sf_face_cones, MPIU_INT, mydonor_cones, leaf_donor_cones, MPI_REPLACE));
614*53b735f8SJames Wright       PetscCall(PetscSFDestroy(&sf_vert_coords));
615*53b735f8SJames Wright       PetscCall(PetscSFDestroy(&sf_face_cones));
616*53b735f8SJames Wright       PetscCallMPI(MPI_Type_free(&vertex_unit));
617*53b735f8SJames Wright       PetscCall(PetscFree2(mydonor_vertices, mydonor_cones));
618*53b735f8SJames Wright     }
619*53b735f8SJames Wright 
620*53b735f8SJames Wright     { // Determine periodic orientation w/rt donor vertices and reorient
621*53b735f8SJames Wright       PetscReal tol = PetscSqr(PETSC_MACHINE_EPSILON * 1e3);
622*53b735f8SJames Wright       PetscInt *periodic2donor, dm_depth, maxConeSize;
623*53b735f8SJames Wright       PetscInt  coords_offset = 0, cones_offset = 0;
624*53b735f8SJames Wright 
625*53b735f8SJames Wright       PetscCall(DMPlexGetDepth(dm, &dm_depth));
626*53b735f8SJames Wright       PetscCall(DMPlexGetMaxSizes(dm, &maxConeSize, NULL));
627*53b735f8SJames Wright       PetscCall(DMGetWorkArray(dm, 2 * PetscPowInt(maxConeSize, dm_depth - 1), MPIU_INT, &periodic2donor));
628*53b735f8SJames Wright 
629*53b735f8SJames Wright       // Translate the periodic face vertices into the donor vertices
630*53b735f8SJames Wright       // Translation stored in periodic2donor
631*53b735f8SJames Wright       for (PetscInt i = 0; i < nleaves; i++) {
632*53b735f8SJames Wright         PetscInt  periodic_face = filocal[i], cl_size, num_verts = face_vertices_size[periodic_face - fStart];
633*53b735f8SJames Wright         PetscInt  cones_size = face_cones_size[periodic_face - fStart], p2d_count = 0;
634*53b735f8SJames Wright         PetscInt *closure = NULL;
635*53b735f8SJames Wright 
636*53b735f8SJames Wright         PetscCall(DMPlexGetTransitiveClosure(dm, periodic_face, PETSC_TRUE, &cl_size, &closure));
637*53b735f8SJames Wright         for (PetscInt p = 0; p < cl_size; p++) {
638*53b735f8SJames Wright           PetscInt     cl_point = closure[2 * p], coords_size, donor_vertex = -1;
639*53b735f8SJames Wright           PetscScalar *coords = NULL;
640*53b735f8SJames Wright 
641*53b735f8SJames Wright           if (!IsPointInsideStratum(cl_point, vStart, vEnd)) continue;
642*53b735f8SJames Wright           PetscCall(DMPlexVecGetClosure(dm, csection, coordinates, cl_point, &coords_size, &coords));
643*53b735f8SJames Wright           PetscAssert(coords_size == dim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %" PetscInt_FMT " has dof size %" PetscInt_FMT ", but should match dimension size %" PetscInt_FMT, cl_point, coords_size, dim);
644*53b735f8SJames Wright 
645*53b735f8SJames Wright           for (PetscInt v = 0; v < num_verts; v++) {
646*53b735f8SJames Wright             PetscReal dist_sqr = 0;
647*53b735f8SJames Wright             for (PetscInt d = 0; d < coords_size; d++) dist_sqr += PetscSqr(PetscRealPart(coords[d]) - leaf_donor_coords[(v + coords_offset) * dim + d]);
648*53b735f8SJames Wright             if (dist_sqr < tol) {
649*53b735f8SJames Wright               donor_vertex = leaf_donor_cones[cones_offset + v];
650*53b735f8SJames Wright               break;
651*53b735f8SJames Wright             }
652*53b735f8SJames Wright           }
653*53b735f8SJames Wright           PetscCheck(donor_vertex >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Periodic face %" PetscInt_FMT " could not find matching donor vertex for vertex %" PetscInt_FMT, periodic_face, cl_point);
654*53b735f8SJames Wright           if (PetscDefined(USE_DEBUG)) {
655*53b735f8SJames Wright             for (PetscInt c = 0; c < p2d_count; c++) PetscCheck(periodic2donor[2 * c + 1] != donor_vertex, comm, PETSC_ERR_PLIB, "Found repeated cone_point in periodic_ordering");
656*53b735f8SJames Wright           }
657*53b735f8SJames Wright 
658*53b735f8SJames Wright           periodic2donor[2 * p2d_count + 0] = cl_point;
659*53b735f8SJames Wright           periodic2donor[2 * p2d_count + 1] = donor_vertex;
660*53b735f8SJames Wright           p2d_count++;
661*53b735f8SJames Wright           PetscCall(DMPlexVecRestoreClosure(dm, csection, coordinates, cl_point, &coords_size, &coords));
662*53b735f8SJames Wright         }
663*53b735f8SJames Wright         coords_offset += num_verts;
664*53b735f8SJames Wright         PetscCall(DMPlexRestoreTransitiveClosure(dm, periodic_face, PETSC_TRUE, &cl_size, &closure));
665*53b735f8SJames Wright 
666*53b735f8SJames Wright         { // Determine periodic orientation w/rt donor vertices and reorient
667*53b735f8SJames Wright           PetscInt      depth, *p2d_cone, face_is_array[1] = {periodic_face};
668*53b735f8SJames Wright           IS           *is_arrays, face_is;
669*53b735f8SJames Wright           PetscSection *section_arrays;
670*53b735f8SJames Wright           PetscInt     *donor_cone_array = &leaf_donor_cones[cones_offset + num_verts];
671*53b735f8SJames Wright 
672*53b735f8SJames Wright           PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 1, face_is_array, PETSC_USE_POINTER, &face_is));
673*53b735f8SJames Wright           PetscCall(DMPlexGetConeRecursive(dm, face_is, &depth, &is_arrays, &section_arrays));
674*53b735f8SJames Wright           PetscCall(ISDestroy(&face_is));
675*53b735f8SJames Wright           PetscCall(DMGetWorkArray(dm, maxConeSize, MPIU_INT, &p2d_cone));
676*53b735f8SJames Wright           for (PetscInt d = 0; d < depth - 1; d++) {
677*53b735f8SJames Wright             PetscInt        pStart, pEnd;
678*53b735f8SJames Wright             PetscSection    section = section_arrays[d];
679*53b735f8SJames Wright             const PetscInt *periodic_cone_arrays, *periodic_point_arrays;
680*53b735f8SJames Wright 
681*53b735f8SJames Wright             PetscCall(ISGetIndices(is_arrays[d], &periodic_cone_arrays));
682*53b735f8SJames Wright             PetscCall(ISGetIndices(is_arrays[d + 1], &periodic_point_arrays)); // Points at d+1 correspond to the cones at d
683*53b735f8SJames Wright             PetscCall(PetscSectionGetChart(section_arrays[d], &pStart, &pEnd));
684*53b735f8SJames Wright             for (PetscInt p = pStart; p < pEnd; p++) {
685*53b735f8SJames Wright               PetscInt periodic_cone_size, periodic_cone_offset, periodic_point = periodic_point_arrays[p];
686*53b735f8SJames Wright 
687*53b735f8SJames Wright               PetscCall(PetscSectionGetDof(section, p, &periodic_cone_size));
688*53b735f8SJames Wright               PetscCall(PetscSectionGetOffset(section, p, &periodic_cone_offset));
689*53b735f8SJames Wright               const PetscInt *periodic_cone = &periodic_cone_arrays[periodic_cone_offset];
690*53b735f8SJames Wright               PetscCall(TranslateConeP2D(periodic_cone, periodic_cone_size, periodic2donor, p2d_count, p2d_cone));
691*53b735f8SJames Wright 
692*53b735f8SJames Wright               // Find the donor cone that matches the periodic point's cone
693*53b735f8SJames Wright               PetscInt  donor_cone_offset = 0, donor_point = -1, *donor_cone = NULL;
694*53b735f8SJames Wright               PetscBool cone_matches = PETSC_FALSE;
695*53b735f8SJames Wright               while (donor_cone_offset < cones_size - num_verts) {
696*53b735f8SJames Wright                 PetscInt donor_cone_size = donor_cone_array[donor_cone_offset];
697*53b735f8SJames Wright                 donor_point              = donor_cone_array[donor_cone_offset + 1];
698*53b735f8SJames Wright                 donor_cone               = &donor_cone_array[donor_cone_offset + 2];
699*53b735f8SJames Wright 
700*53b735f8SJames Wright                 if (donor_cone_size != periodic_cone_size) goto next_cone;
701*53b735f8SJames Wright                 for (PetscInt c = 0; c < periodic_cone_size; c++) {
702*53b735f8SJames Wright                   cone_matches = SearchIntArray(donor_cone[c], periodic_cone_size, p2d_cone);
703*53b735f8SJames Wright                   if (!cone_matches) goto next_cone;
704*53b735f8SJames Wright                 }
705*53b735f8SJames Wright                 // Save the found donor cone's point to the translation array. These will be used for higher depth points.
706*53b735f8SJames Wright                 // i.e. we save the edge translations for when we look for face cones
707*53b735f8SJames Wright                 periodic2donor[2 * p2d_count + 0] = periodic_point;
708*53b735f8SJames Wright                 periodic2donor[2 * p2d_count + 1] = donor_point;
709*53b735f8SJames Wright                 p2d_count++;
710*53b735f8SJames Wright                 break;
711*53b735f8SJames Wright 
712*53b735f8SJames Wright               next_cone:
713*53b735f8SJames Wright                 donor_cone_offset += donor_cone_size + 2;
714*53b735f8SJames Wright               }
715*53b735f8SJames Wright               PetscCheck(donor_cone_offset < cones_size - num_verts, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not find donor cone equivalent to cone of periodic point %" PetscInt_FMT, periodic_point);
716*53b735f8SJames Wright 
717*53b735f8SJames Wright               { // Compare the donor cone with the translated periodic cone and reorient
718*53b735f8SJames Wright                 PetscInt       ornt;
719*53b735f8SJames Wright                 DMPolytopeType cell_type;
720*53b735f8SJames Wright                 PetscBool      found;
721*53b735f8SJames Wright                 PetscCall(DMPlexGetCellType(dm, periodic_point, &cell_type));
722*53b735f8SJames Wright                 PetscCall(DMPolytopeMatchOrientation(cell_type, donor_cone, p2d_cone, &ornt, &found));
723*53b735f8SJames Wright                 PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not find transformation between donor (%" PetscInt_FMT ") and periodic (%" PetscInt_FMT ") cone's", periodic_point, donor_point);
724*53b735f8SJames Wright                 if (debug_printing) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Reorienting point %" PetscInt_FMT " by %" PetscInt_FMT "\n", periodic_point, ornt));
725*53b735f8SJames Wright                 PetscCall(DMPlexOrientPointWithCorrections(dm, periodic_point, ornt));
726*53b735f8SJames Wright               }
727*53b735f8SJames Wright             }
728*53b735f8SJames Wright             PetscCall(ISRestoreIndices(is_arrays[d], &periodic_cone_arrays));
729*53b735f8SJames Wright             PetscCall(ISRestoreIndices(is_arrays[d + 1], &periodic_point_arrays));
730*53b735f8SJames Wright           }
731*53b735f8SJames Wright           PetscCall(DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &p2d_cone));
732*53b735f8SJames Wright           PetscCall(DMPlexRestoreConeRecursive(dm, face_is, &depth, &is_arrays, &section_arrays));
733*53b735f8SJames Wright         }
734*53b735f8SJames Wright 
735*53b735f8SJames Wright         PetscCall(DMPlexRestoreTransitiveClosure(dm, periodic_face, PETSC_TRUE, &cl_size, &closure));
736*53b735f8SJames Wright         cones_offset += cones_size;
737*53b735f8SJames Wright       }
738*53b735f8SJames Wright       PetscCall(DMRestoreWorkArray(dm, 2 * PetscPowInt(maxConeSize, dm_depth - 1), MPIU_INT, &periodic2donor));
739*53b735f8SJames Wright     }
740*53b735f8SJames Wright 
741*53b735f8SJames Wright     PetscCall(PetscFree2(leaf_donor_coords, leaf_donor_cones));
742*53b735f8SJames Wright     PetscCall(PetscFree2(face_vertices_size, face_cones_size));
743*53b735f8SJames Wright   }
744*53b735f8SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
745*53b735f8SJames Wright }
746*53b735f8SJames Wright 
747*53b735f8SJames Wright // Start with an SF for a positive depth (e.g., faces) and create a new SF for matched closure.
7486725e60dSJed Brown //
7495f06a3ddSJed Brown // Output Arguments:
7505f06a3ddSJed Brown //
7515f06a3ddSJed Brown // + closure_sf - augmented point SF (see `DMGetPointSF()`) that includes the faces and all points in its closure. This
7525f06a3ddSJed Brown //   can be used to create a global section and section SF.
753b83f62b0SJames 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
7545f06a3ddSJed Brown //   transformation to periodic dofs; see DMPeriodicCoordinateSetUp_Internal().
7555f06a3ddSJed Brown //
756b83f62b0SJames Wright static PetscErrorCode DMPlexCreateIsoperiodicPointSF_Private(DM dm, PetscInt num_face_sfs, PetscSF *face_sfs, PetscSF *closure_sf, IS **is_points)
7576725e60dSJed Brown {
7586725e60dSJed Brown   MPI_Comm           comm;
7596725e60dSJed Brown   PetscMPIInt        rank;
7606725e60dSJed Brown   PetscSF            point_sf;
761b83f62b0SJames Wright   PetscInt           nroots, nleaves;
762b83f62b0SJames Wright   const PetscInt    *filocal;
763b83f62b0SJames Wright   const PetscSFNode *firemote;
7646725e60dSJed Brown 
7656725e60dSJed Brown   PetscFunctionBegin;
7666725e60dSJed Brown   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
7676725e60dSJed Brown   PetscCallMPI(MPI_Comm_rank(comm, &rank));
7686725e60dSJed Brown   PetscCall(DMGetPointSF(dm, &point_sf)); // Point SF has remote points
769b83f62b0SJames Wright   PetscCall(PetscMalloc1(num_face_sfs, is_points));
770b83f62b0SJames Wright 
771*53b735f8SJames Wright   PetscCall(DMPlexCorrectOrientationForIsoperiodic(dm));
772*53b735f8SJames Wright 
773b83f62b0SJames Wright   for (PetscInt f = 0; f < num_face_sfs; f++) {
774b83f62b0SJames Wright     PetscSF   face_sf = face_sfs[f];
77593defd67SJames Wright     PetscInt *cl_sizes;
77693defd67SJames Wright     PetscInt  fStart, fEnd, rootbuffersize, leafbuffersize;
7776725e60dSJed Brown     PetscSF   sf_closure;
77893defd67SJames Wright     PetscBT   rootbt;
77993defd67SJames Wright 
78093defd67SJames Wright     PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
78193defd67SJames Wright     PetscCall(PetscMalloc1(fEnd - fStart, &cl_sizes));
78293defd67SJames Wright     for (PetscInt f = fStart, index = 0; f < fEnd; f++, index++) {
78393defd67SJames Wright       PetscInt cl_size, *closure = NULL;
78493defd67SJames Wright       PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &cl_size, &closure));
78593defd67SJames Wright       cl_sizes[index] = cl_size - 1;
78693defd67SJames Wright       PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &cl_size, &closure));
78793defd67SJames Wright     }
78893defd67SJames Wright 
78993defd67SJames Wright     PetscCall(CreateDonorToPeriodicSF(dm, face_sf, fStart, fEnd, cl_sizes, &rootbuffersize, &leafbuffersize, &rootbt, &sf_closure));
79093defd67SJames Wright     PetscCall(PetscFree(cl_sizes));
79193defd67SJames Wright     PetscCall(PetscSFGetGraph(face_sf, &nroots, &nleaves, &filocal, &firemote));
7926725e60dSJed Brown 
793b83f62b0SJames Wright     PetscSFNode *leaf_donor_closure;
794b83f62b0SJames Wright     { // Pack root buffer with owner for every point in the root cones
7956725e60dSJed Brown       PetscSFNode       *donor_closure;
796b83f62b0SJames Wright       const PetscInt    *pilocal;
797b83f62b0SJames Wright       const PetscSFNode *piremote;
798b83f62b0SJames Wright       PetscInt           npoints;
799b83f62b0SJames Wright 
800b83f62b0SJames Wright       PetscCall(PetscSFGetGraph(point_sf, NULL, &npoints, &pilocal, &piremote));
80193defd67SJames Wright       PetscCall(PetscCalloc1(rootbuffersize, &donor_closure));
80293defd67SJames Wright       for (PetscInt p = 0, root_offset = 0; p < nroots; p++) {
80393defd67SJames Wright         if (!PetscBTLookup(rootbt, p)) continue;
80493defd67SJames Wright         PetscInt cl_size, *closure = NULL;
8056725e60dSJed Brown         PetscCall(DMPlexGetTransitiveClosure(dm, p, PETSC_TRUE, &cl_size, &closure));
8066725e60dSJed Brown         for (PetscInt j = 1; j < cl_size; j++) {
8076725e60dSJed Brown           PetscInt c = closure[2 * j];
8086725e60dSJed Brown           if (pilocal) {
8096725e60dSJed Brown             PetscInt found = -1;
8106725e60dSJed Brown             if (npoints > 0) PetscCall(PetscFindInt(c, npoints, pilocal, &found));
8116725e60dSJed Brown             if (found >= 0) {
8126725e60dSJed Brown               donor_closure[root_offset++] = piremote[found];
8136725e60dSJed Brown               continue;
8146725e60dSJed Brown             }
8156725e60dSJed Brown           }
8166725e60dSJed Brown           // we own c
8176725e60dSJed Brown           donor_closure[root_offset].rank  = rank;
8186725e60dSJed Brown           donor_closure[root_offset].index = c;
8196725e60dSJed Brown           root_offset++;
8206725e60dSJed Brown         }
8216725e60dSJed Brown         PetscCall(DMPlexRestoreTransitiveClosure(dm, p, PETSC_TRUE, &cl_size, &closure));
8226725e60dSJed Brown       }
8236725e60dSJed Brown 
82493defd67SJames Wright       PetscCall(PetscMalloc1(leafbuffersize, &leaf_donor_closure));
8256497c311SBarry Smith       PetscCall(PetscSFBcastBegin(sf_closure, MPIU_SF_NODE, donor_closure, leaf_donor_closure, MPI_REPLACE));
8266497c311SBarry Smith       PetscCall(PetscSFBcastEnd(sf_closure, MPIU_SF_NODE, donor_closure, leaf_donor_closure, MPI_REPLACE));
8276725e60dSJed Brown       PetscCall(PetscSFDestroy(&sf_closure));
8286725e60dSJed Brown       PetscCall(PetscFree(donor_closure));
829b83f62b0SJames Wright     }
8306725e60dSJed Brown 
8316725e60dSJed Brown     PetscSFNode *new_iremote;
8326725e60dSJed Brown     PetscCall(PetscCalloc1(nroots, &new_iremote));
8336725e60dSJed Brown     for (PetscInt i = 0; i < nroots; i++) new_iremote[i].rank = -1;
8346725e60dSJed Brown     // Walk leaves and match vertices
83593defd67SJames Wright     for (PetscInt i = 0, leaf_offset = 0; i < nleaves; i++) {
8366725e60dSJed Brown       PetscInt  point   = filocal[i], cl_size;
8376725e60dSJed Brown       PetscInt *closure = NULL;
8386725e60dSJed Brown       PetscCall(DMPlexGetTransitiveClosure(dm, point, PETSC_TRUE, &cl_size, &closure));
839*53b735f8SJames Wright       for (PetscInt j = 1; j < cl_size; j++) {
8406725e60dSJed Brown         PetscInt    c  = closure[2 * j];
8416725e60dSJed Brown         PetscSFNode lc = leaf_donor_closure[leaf_offset];
8426725e60dSJed Brown         // printf("[%d] face %d.%d: %d ?-- (%d,%d)\n", rank, point, j, c, lc.rank, lc.index);
8436725e60dSJed Brown         if (new_iremote[c].rank == -1) {
8446725e60dSJed Brown           new_iremote[c] = lc;
8455f06a3ddSJed 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");
8466725e60dSJed Brown         leaf_offset++;
8476725e60dSJed Brown       }
8486725e60dSJed Brown       PetscCall(DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, &cl_size, &closure));
8496725e60dSJed Brown     }
8506725e60dSJed Brown     PetscCall(PetscFree(leaf_donor_closure));
8516725e60dSJed Brown 
8526725e60dSJed Brown     // Include face points in closure SF
8536725e60dSJed Brown     for (PetscInt i = 0; i < nleaves; i++) new_iremote[filocal[i]] = firemote[i];
8546725e60dSJed Brown     // consolidate leaves
85593defd67SJames Wright     PetscInt *leafdata;
85693defd67SJames Wright     PetscCall(PetscMalloc1(nroots, &leafdata));
8576725e60dSJed Brown     PetscInt num_new_leaves = 0;
8586725e60dSJed Brown     for (PetscInt i = 0; i < nroots; i++) {
8596725e60dSJed Brown       if (new_iremote[i].rank == -1) continue;
8606725e60dSJed Brown       new_iremote[num_new_leaves] = new_iremote[i];
8616725e60dSJed Brown       leafdata[num_new_leaves]    = i;
8626725e60dSJed Brown       num_new_leaves++;
8636725e60dSJed Brown     }
864b83f62b0SJames Wright     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, num_new_leaves, leafdata, PETSC_COPY_VALUES, &(*is_points)[f]));
8656725e60dSJed Brown 
8666725e60dSJed Brown     PetscSF csf;
8676725e60dSJed Brown     PetscCall(PetscSFCreate(comm, &csf));
8686725e60dSJed Brown     PetscCall(PetscSFSetGraph(csf, nroots, num_new_leaves, leafdata, PETSC_COPY_VALUES, new_iremote, PETSC_COPY_VALUES));
8696725e60dSJed Brown     PetscCall(PetscFree(new_iremote)); // copy and delete because new_iremote is longer than it needs to be
87093defd67SJames Wright     PetscCall(PetscFree(leafdata));
87193defd67SJames Wright     PetscCall(PetscBTDestroy(&rootbt));
8726725e60dSJed Brown 
873b83f62b0SJames Wright     PetscInt npoints;
874b83f62b0SJames Wright     PetscCall(PetscSFGetGraph(point_sf, NULL, &npoints, NULL, NULL));
8755f06a3ddSJed Brown     if (npoints < 0) { // empty point_sf
8766725e60dSJed Brown       *closure_sf = csf;
8775f06a3ddSJed Brown     } else {
8785f06a3ddSJed Brown       PetscCall(PetscSFMerge(point_sf, csf, closure_sf));
8795f06a3ddSJed Brown       PetscCall(PetscSFDestroy(&csf));
8805f06a3ddSJed Brown     }
881d8b4a066SPierre Jolivet     if (f > 0) PetscCall(PetscSFDestroy(&point_sf)); // Only destroy if point_sf is from previous calls to PetscSFMerge
882b83f62b0SJames Wright     point_sf = *closure_sf;                          // Use combined point + isoperiodic SF to define point ownership for further face_sf
883b83f62b0SJames Wright   }
8845f06a3ddSJed Brown   PetscCall(PetscObjectSetName((PetscObject)*closure_sf, "Composed Periodic Points"));
8853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8866725e60dSJed Brown }
8876725e60dSJed Brown 
8885f06a3ddSJed Brown static PetscErrorCode DMGetIsoperiodicPointSF_Plex(DM dm, PetscSF *sf)
8896725e60dSJed Brown {
8906725e60dSJed Brown   DM_Plex *plex = (DM_Plex *)dm->data;
8916725e60dSJed Brown 
8926725e60dSJed Brown   PetscFunctionBegin;
893b83f62b0SJames 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));
8946725e60dSJed Brown   if (sf) *sf = plex->periodic.composed_sf;
8953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8966725e60dSJed Brown }
8976725e60dSJed Brown 
8985f06a3ddSJed Brown PetscErrorCode DMPlexMigrateIsoperiodicFaceSF_Internal(DM old_dm, DM dm, PetscSF sf_migration)
8995f06a3ddSJed Brown {
9005f06a3ddSJed Brown   DM_Plex    *plex = (DM_Plex *)old_dm->data;
901a45b0079SJames Wright   PetscSF     sf_point, *new_face_sfs;
9025f06a3ddSJed Brown   PetscMPIInt rank;
9035f06a3ddSJed Brown 
9045f06a3ddSJed Brown   PetscFunctionBegin;
9051fca310dSJames Wright   if (!plex->periodic.face_sfs) PetscFunctionReturn(PETSC_SUCCESS);
9065f06a3ddSJed Brown   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
9075f06a3ddSJed Brown   PetscCall(DMGetPointSF(dm, &sf_point));
908a45b0079SJames Wright   PetscCall(PetscMalloc1(plex->periodic.num_face_sfs, &new_face_sfs));
909a45b0079SJames Wright 
910a45b0079SJames Wright   for (PetscInt f = 0; f < plex->periodic.num_face_sfs; f++) {
9115f06a3ddSJed Brown     PetscInt           old_npoints, new_npoints, old_nleaf, new_nleaf, point_nleaf;
9125f06a3ddSJed Brown     PetscSFNode       *new_leafdata, *rootdata, *leafdata;
9135f06a3ddSJed Brown     const PetscInt    *old_local, *point_local;
9145f06a3ddSJed Brown     const PetscSFNode *old_remote, *point_remote;
9156497c311SBarry Smith 
916a45b0079SJames Wright     PetscCall(PetscSFGetGraph(plex->periodic.face_sfs[f], &old_npoints, &old_nleaf, &old_local, &old_remote));
9175f06a3ddSJed Brown     PetscCall(PetscSFGetGraph(sf_migration, NULL, &new_nleaf, NULL, NULL));
9185f06a3ddSJed Brown     PetscCall(PetscSFGetGraph(sf_point, &new_npoints, &point_nleaf, &point_local, &point_remote));
9195f06a3ddSJed Brown     PetscAssert(new_nleaf == new_npoints, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected migration leaf space to match new point root space");
9205f06a3ddSJed Brown     PetscCall(PetscMalloc3(old_npoints, &rootdata, old_npoints, &leafdata, new_npoints, &new_leafdata));
9215f06a3ddSJed Brown 
9225f06a3ddSJed Brown     // Fill new_leafdata with new owners of all points
9235f06a3ddSJed Brown     for (PetscInt i = 0; i < new_npoints; i++) {
9245f06a3ddSJed Brown       new_leafdata[i].rank  = rank;
9255f06a3ddSJed Brown       new_leafdata[i].index = i;
9265f06a3ddSJed Brown     }
9275f06a3ddSJed Brown     for (PetscInt i = 0; i < point_nleaf; i++) {
9285f06a3ddSJed Brown       PetscInt j      = point_local[i];
9295f06a3ddSJed Brown       new_leafdata[j] = point_remote[i];
9305f06a3ddSJed Brown     }
9315f06a3ddSJed Brown     // REPLACE is okay because every leaf agrees about the new owners
9326497c311SBarry Smith     PetscCall(PetscSFReduceBegin(sf_migration, MPIU_SF_NODE, new_leafdata, rootdata, MPI_REPLACE));
9336497c311SBarry Smith     PetscCall(PetscSFReduceEnd(sf_migration, MPIU_SF_NODE, new_leafdata, rootdata, MPI_REPLACE));
9345f06a3ddSJed Brown     // rootdata now contains the new owners
9355f06a3ddSJed Brown 
9365f06a3ddSJed Brown     // Send to leaves of old space
9375f06a3ddSJed Brown     for (PetscInt i = 0; i < old_npoints; i++) {
9385f06a3ddSJed Brown       leafdata[i].rank  = -1;
9395f06a3ddSJed Brown       leafdata[i].index = -1;
9405f06a3ddSJed Brown     }
9416497c311SBarry Smith     PetscCall(PetscSFBcastBegin(plex->periodic.face_sfs[f], MPIU_SF_NODE, rootdata, leafdata, MPI_REPLACE));
9426497c311SBarry Smith     PetscCall(PetscSFBcastEnd(plex->periodic.face_sfs[f], MPIU_SF_NODE, rootdata, leafdata, MPI_REPLACE));
9435f06a3ddSJed Brown 
9445f06a3ddSJed Brown     // Send to new leaf space
9456497c311SBarry Smith     PetscCall(PetscSFBcastBegin(sf_migration, MPIU_SF_NODE, leafdata, new_leafdata, MPI_REPLACE));
9466497c311SBarry Smith     PetscCall(PetscSFBcastEnd(sf_migration, MPIU_SF_NODE, leafdata, new_leafdata, MPI_REPLACE));
9475f06a3ddSJed Brown 
9485f06a3ddSJed Brown     PetscInt     nface = 0, *new_local;
9495f06a3ddSJed Brown     PetscSFNode *new_remote;
9505f06a3ddSJed Brown     for (PetscInt i = 0; i < new_npoints; i++) nface += (new_leafdata[i].rank >= 0);
9515f06a3ddSJed Brown     PetscCall(PetscMalloc1(nface, &new_local));
9525f06a3ddSJed Brown     PetscCall(PetscMalloc1(nface, &new_remote));
9535f06a3ddSJed Brown     nface = 0;
9545f06a3ddSJed Brown     for (PetscInt i = 0; i < new_npoints; i++) {
9555f06a3ddSJed Brown       if (new_leafdata[i].rank == -1) continue;
9565f06a3ddSJed Brown       new_local[nface]  = i;
9575f06a3ddSJed Brown       new_remote[nface] = new_leafdata[i];
9585f06a3ddSJed Brown       nface++;
9595f06a3ddSJed Brown     }
9605f06a3ddSJed Brown     PetscCall(PetscFree3(rootdata, leafdata, new_leafdata));
961a45b0079SJames Wright     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), &new_face_sfs[f]));
962a45b0079SJames Wright     PetscCall(PetscSFSetGraph(new_face_sfs[f], new_npoints, nface, new_local, PETSC_OWN_POINTER, new_remote, PETSC_OWN_POINTER));
963a45b0079SJames Wright     {
964a45b0079SJames Wright       char new_face_sf_name[PETSC_MAX_PATH_LEN];
965a45b0079SJames Wright       PetscCall(PetscSNPrintf(new_face_sf_name, sizeof new_face_sf_name, "Migrated Isoperiodic Faces #%" PetscInt_FMT, f));
966a45b0079SJames Wright       PetscCall(PetscObjectSetName((PetscObject)new_face_sfs[f], new_face_sf_name));
967a45b0079SJames Wright     }
968a45b0079SJames Wright   }
969a45b0079SJames Wright 
970a45b0079SJames Wright   PetscCall(DMPlexSetIsoperiodicFaceSF(dm, plex->periodic.num_face_sfs, new_face_sfs));
971a45b0079SJames Wright   PetscCall(DMPlexSetIsoperiodicFaceTransform(dm, plex->periodic.num_face_sfs, (PetscScalar *)plex->periodic.transform));
972a45b0079SJames Wright   for (PetscInt f = 0; f < plex->periodic.num_face_sfs; f++) PetscCall(PetscSFDestroy(&new_face_sfs[f]));
973a45b0079SJames Wright   PetscCall(PetscFree(new_face_sfs));
9743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9755f06a3ddSJed Brown }
9765f06a3ddSJed Brown 
9776725e60dSJed Brown PetscErrorCode DMPeriodicCoordinateSetUp_Internal(DM dm)
9786725e60dSJed Brown {
9796725e60dSJed Brown   DM_Plex   *plex = (DM_Plex *)dm->data;
9806497c311SBarry Smith   PetscCount count;
981ddedc8f6SJames Wright   IS         isdof;
982ddedc8f6SJames Wright   PetscInt   dim;
9834d86920dSPierre Jolivet 
9846725e60dSJed Brown   PetscFunctionBegin;
9851fca310dSJames Wright   if (!plex->periodic.face_sfs) PetscFunctionReturn(PETSC_SUCCESS);
986*53b735f8SJames Wright   PetscCheck(plex->periodic.periodic_points, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Isoperiodic PointSF must be created before this function is called");
9876725e60dSJed Brown 
988*53b735f8SJames Wright   PetscCall(DMGetCoordinateDim(dm, &dim));
989ddedc8f6SJames Wright   dm->periodic.num_affines = plex->periodic.num_face_sfs;
990*53b735f8SJames Wright   PetscCall(PetscFree2(dm->periodic.affine_to_local, dm->periodic.affine));
991ddedc8f6SJames Wright   PetscCall(PetscMalloc2(dm->periodic.num_affines, &dm->periodic.affine_to_local, dm->periodic.num_affines, &dm->periodic.affine));
992ddedc8f6SJames Wright 
993ddedc8f6SJames Wright   for (PetscInt f = 0; f < plex->periodic.num_face_sfs; f++) {
9946497c311SBarry Smith     PetscInt        npoints, vsize, isize;
9956725e60dSJed Brown     const PetscInt *points;
996ddedc8f6SJames Wright     IS              is = plex->periodic.periodic_points[f];
9976725e60dSJed Brown     PetscSegBuffer  seg;
9986725e60dSJed Brown     PetscSection    section;
9996497c311SBarry Smith     PetscInt       *ind;
10006497c311SBarry Smith     Vec             L, P;
10016497c311SBarry Smith     VecType         vec_type;
10026497c311SBarry Smith     VecScatter      scatter;
10036497c311SBarry Smith     PetscScalar    *x;
10046497c311SBarry Smith 
10056725e60dSJed Brown     PetscCall(DMGetLocalSection(dm, &section));
10066725e60dSJed Brown     PetscCall(PetscSegBufferCreate(sizeof(PetscInt), 32, &seg));
10076725e60dSJed Brown     PetscCall(ISGetSize(is, &npoints));
10086725e60dSJed Brown     PetscCall(ISGetIndices(is, &points));
10096725e60dSJed Brown     for (PetscInt i = 0; i < npoints; i++) {
10106725e60dSJed Brown       PetscInt point = points[i], off, dof;
10116497c311SBarry Smith 
10126725e60dSJed Brown       PetscCall(PetscSectionGetOffset(section, point, &off));
10136725e60dSJed Brown       PetscCall(PetscSectionGetDof(section, point, &dof));
10146725e60dSJed Brown       PetscAssert(dof % dim == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unexpected dof %" PetscInt_FMT " not divisible by dimension %" PetscInt_FMT, dof, dim);
10156725e60dSJed Brown       for (PetscInt j = 0; j < dof / dim; j++) {
10166725e60dSJed Brown         PetscInt *slot;
10176497c311SBarry Smith 
10186725e60dSJed Brown         PetscCall(PetscSegBufferGetInts(seg, 1, &slot));
1019729bdd54SJed Brown         *slot = off / dim + j;
10206725e60dSJed Brown       }
10216725e60dSJed Brown     }
10226725e60dSJed Brown     PetscCall(PetscSegBufferGetSize(seg, &count));
10236725e60dSJed Brown     PetscCall(PetscSegBufferExtractAlloc(seg, &ind));
10246725e60dSJed Brown     PetscCall(PetscSegBufferDestroy(&seg));
10256497c311SBarry Smith     PetscCall(PetscIntCast(count, &isize));
10266497c311SBarry Smith     PetscCall(ISCreateBlock(PETSC_COMM_SELF, dim, isize, ind, PETSC_OWN_POINTER, &isdof));
10276497c311SBarry Smith 
10286497c311SBarry Smith     PetscCall(PetscIntCast(count * dim, &vsize));
10296725e60dSJed Brown     PetscCall(DMGetLocalVector(dm, &L));
10306725e60dSJed Brown     PetscCall(VecCreate(PETSC_COMM_SELF, &P));
10316497c311SBarry Smith     PetscCall(VecSetSizes(P, vsize, vsize));
10326725e60dSJed Brown     PetscCall(VecGetType(L, &vec_type));
10336725e60dSJed Brown     PetscCall(VecSetType(P, vec_type));
10346725e60dSJed Brown     PetscCall(VecScatterCreate(P, NULL, L, isdof, &scatter));
10356725e60dSJed Brown     PetscCall(DMRestoreLocalVector(dm, &L));
10366725e60dSJed Brown     PetscCall(ISDestroy(&isdof));
10376725e60dSJed Brown 
10386725e60dSJed Brown     PetscCall(VecGetArrayWrite(P, &x));
10396497c311SBarry Smith     for (PetscCount i = 0; i < count; i++) {
1040ddedc8f6SJames Wright       for (PetscInt j = 0; j < dim; j++) x[i * dim + j] = plex->periodic.transform[f][j][3];
10416725e60dSJed Brown     }
10426725e60dSJed Brown     PetscCall(VecRestoreArrayWrite(P, &x));
10436725e60dSJed Brown 
1044ddedc8f6SJames Wright     dm->periodic.affine_to_local[f] = scatter;
1045ddedc8f6SJames Wright     dm->periodic.affine[f]          = P;
1046ddedc8f6SJames Wright   }
10476725e60dSJed Brown   PetscCall(DMGlobalToLocalHookAdd(dm, NULL, DMCoordAddPeriodicOffsets_Private, NULL));
10483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10496725e60dSJed Brown }
10506725e60dSJed Brown 
10513e72e933SJed Brown PetscErrorCode DMPlexCreateBoxMesh_Tensor_SFC_Internal(DM dm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate)
10523e72e933SJed Brown {
10533e72e933SJed Brown   PetscInt  eextent[3] = {1, 1, 1}, vextent[3] = {1, 1, 1};
10543e72e933SJed Brown   const Ijk closure_1[] = {
10553e72e933SJed Brown     {0, 0, 0},
10563e72e933SJed Brown     {1, 0, 0},
10573e72e933SJed Brown   };
10583e72e933SJed Brown   const Ijk closure_2[] = {
10593e72e933SJed Brown     {0, 0, 0},
10603e72e933SJed Brown     {1, 0, 0},
10613e72e933SJed Brown     {1, 1, 0},
10623e72e933SJed Brown     {0, 1, 0},
10633e72e933SJed Brown   };
10643e72e933SJed Brown   const Ijk closure_3[] = {
10653e72e933SJed Brown     {0, 0, 0},
10663e72e933SJed Brown     {0, 1, 0},
10673e72e933SJed Brown     {1, 1, 0},
10683e72e933SJed Brown     {1, 0, 0},
10693e72e933SJed Brown     {0, 0, 1},
10703e72e933SJed Brown     {1, 0, 1},
10713e72e933SJed Brown     {1, 1, 1},
10723e72e933SJed Brown     {0, 1, 1},
10733e72e933SJed Brown   };
10743431e603SJed Brown   const Ijk *const closure_dim[] = {NULL, closure_1, closure_2, closure_3};
10753431e603SJed Brown   // This must be kept consistent with DMPlexCreateCubeMesh_Internal
10763431e603SJed Brown   const PetscInt        face_marker_1[]   = {1, 2};
10773431e603SJed Brown   const PetscInt        face_marker_2[]   = {4, 2, 1, 3};
10783431e603SJed Brown   const PetscInt        face_marker_3[]   = {6, 5, 3, 4, 1, 2};
10793431e603SJed Brown   const PetscInt *const face_marker_dim[] = {NULL, face_marker_1, face_marker_2, face_marker_3};
10806725e60dSJed Brown   // Orient faces so the normal is in the positive axis and the first vertex is the one closest to zero.
10816725e60dSJed Brown   // These orientations can be determined by examining cones of a reference quad and hex element.
10826725e60dSJed Brown   const PetscInt        face_orient_1[]   = {0, 0};
10836725e60dSJed Brown   const PetscInt        face_orient_2[]   = {-1, 0, 0, -1};
10846725e60dSJed Brown   const PetscInt        face_orient_3[]   = {-2, 0, -2, 1, -2, 0};
10856725e60dSJed Brown   const PetscInt *const face_orient_dim[] = {NULL, face_orient_1, face_orient_2, face_orient_3};
10863e72e933SJed Brown 
10873e72e933SJed Brown   PetscFunctionBegin;
1088ea7b3863SJames Wright   PetscCall(PetscLogEventBegin(DMPLEX_CreateBoxSFC, dm, 0, 0, 0));
10894f572ea9SToby Isaac   PetscAssertPointer(dm, 1);
10903e72e933SJed Brown   PetscValidLogicalCollectiveInt(dm, dim, 2);
10913e72e933SJed Brown   PetscCall(DMSetDimension(dm, dim));
10923e72e933SJed Brown   PetscMPIInt rank, size;
10933e72e933SJed Brown   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
10943e72e933SJed Brown   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
10953e72e933SJed Brown   for (PetscInt i = 0; i < dim; i++) {
10963e72e933SJed Brown     eextent[i] = faces[i];
10973e72e933SJed Brown     vextent[i] = faces[i] + 1;
10983e72e933SJed Brown   }
1099c77877e3SJed Brown   ZLayout layout;
1100c77877e3SJed Brown   PetscCall(ZLayoutCreate(size, eextent, vextent, &layout));
11013e72e933SJed Brown   PetscZSet vset; // set of all vertices in the closure of the owned elements
11023e72e933SJed Brown   PetscCall(PetscZSetCreate(&vset));
11033e72e933SJed Brown   PetscInt local_elems = 0;
11043e72e933SJed Brown   for (ZCode z = layout.zstarts[rank]; z < layout.zstarts[rank + 1]; z++) {
11053e72e933SJed Brown     Ijk loc = ZCodeSplit(z);
11063ba16761SJacob Faibussowitsch     if (IjkActive(layout.vextent, loc)) PetscCall(PetscZSetAdd(vset, z));
11076547ddbdSJed Brown     else {
11086547ddbdSJed Brown       z += ZStepOct(z);
11096547ddbdSJed Brown       continue;
11106547ddbdSJed Brown     }
11113e72e933SJed Brown     if (IjkActive(layout.eextent, loc)) {
11123e72e933SJed Brown       local_elems++;
11133e72e933SJed Brown       // Add all neighboring vertices to set
11143e72e933SJed Brown       for (PetscInt n = 0; n < PetscPowInt(2, dim); n++) {
11153e72e933SJed Brown         Ijk   inc  = closure_dim[dim][n];
11163e72e933SJed Brown         Ijk   nloc = {loc.i + inc.i, loc.j + inc.j, loc.k + inc.k};
11173e72e933SJed Brown         ZCode v    = ZEncode(nloc);
11183ba16761SJacob Faibussowitsch         PetscCall(PetscZSetAdd(vset, v));
11193e72e933SJed Brown       }
11203e72e933SJed Brown     }
11213e72e933SJed Brown   }
11223e72e933SJed Brown   PetscInt local_verts, off = 0;
11233e72e933SJed Brown   ZCode   *vert_z;
11243e72e933SJed Brown   PetscCall(PetscZSetGetSize(vset, &local_verts));
11253e72e933SJed Brown   PetscCall(PetscMalloc1(local_verts, &vert_z));
11263e72e933SJed Brown   PetscCall(PetscZSetGetElems(vset, &off, vert_z));
11273e72e933SJed Brown   PetscCall(PetscZSetDestroy(&vset));
11283e72e933SJed Brown   // ZCode is unsigned for bitwise convenience, but highest bit should never be set, so can interpret as signed
11293e72e933SJed Brown   PetscCall(PetscSortInt64(local_verts, (PetscInt64 *)vert_z));
11303e72e933SJed Brown 
11313e72e933SJed Brown   PetscCall(DMPlexSetChart(dm, 0, local_elems + local_verts));
11323e72e933SJed Brown   for (PetscInt e = 0; e < local_elems; e++) PetscCall(DMPlexSetConeSize(dm, e, PetscPowInt(2, dim)));
11333e72e933SJed Brown   PetscCall(DMSetUp(dm));
11343e72e933SJed Brown   {
11353e72e933SJed Brown     PetscInt e = 0;
11363e72e933SJed Brown     for (ZCode z = layout.zstarts[rank]; z < layout.zstarts[rank + 1]; z++) {
11373e72e933SJed Brown       Ijk loc = ZCodeSplit(z);
11386547ddbdSJed Brown       if (!IjkActive(layout.eextent, loc)) {
11396547ddbdSJed Brown         z += ZStepOct(z);
11406547ddbdSJed Brown         continue;
11416547ddbdSJed Brown       }
11423e72e933SJed Brown       PetscInt cone[8], orient[8] = {0};
11433e72e933SJed Brown       for (PetscInt n = 0; n < PetscPowInt(2, dim); n++) {
11443e72e933SJed Brown         Ijk      inc  = closure_dim[dim][n];
11453e72e933SJed Brown         Ijk      nloc = {loc.i + inc.i, loc.j + inc.j, loc.k + inc.k};
11463e72e933SJed Brown         ZCode    v    = ZEncode(nloc);
11473e72e933SJed Brown         PetscInt ci   = ZCodeFind(v, local_verts, vert_z);
11483e72e933SJed Brown         PetscAssert(ci >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Did not find neighbor vertex in set");
11493e72e933SJed Brown         cone[n] = local_elems + ci;
11503e72e933SJed Brown       }
11513e72e933SJed Brown       PetscCall(DMPlexSetCone(dm, e, cone));
11523e72e933SJed Brown       PetscCall(DMPlexSetConeOrientation(dm, e, orient));
11533e72e933SJed Brown       e++;
11543e72e933SJed Brown     }
11553e72e933SJed Brown   }
11563e72e933SJed Brown 
11573e72e933SJed Brown   PetscCall(DMPlexSymmetrize(dm));
11583e72e933SJed Brown   PetscCall(DMPlexStratify(dm));
11596725e60dSJed Brown 
11603e72e933SJed Brown   { // Create point SF
11613e72e933SJed Brown     PetscSF sf;
11623ba16761SJacob Faibussowitsch     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), &sf));
11633e72e933SJed Brown     PetscInt owned_verts = ZCodeFind(layout.zstarts[rank + 1], local_verts, vert_z);
11643e72e933SJed Brown     if (owned_verts < 0) owned_verts = -(owned_verts + 1); // We don't care whether the key was found
11653e72e933SJed Brown     PetscInt     num_ghosts = local_verts - owned_verts;   // Due to sorting, owned vertices always come first
11663e72e933SJed Brown     PetscInt    *local_ghosts;
11673e72e933SJed Brown     PetscSFNode *ghosts;
11683e72e933SJed Brown     PetscCall(PetscMalloc1(num_ghosts, &local_ghosts));
11693e72e933SJed Brown     PetscCall(PetscMalloc1(num_ghosts, &ghosts));
11703e72e933SJed Brown     for (PetscInt i = 0; i < num_ghosts;) {
11713e72e933SJed Brown       ZCode       z = vert_z[owned_verts + i];
1172835f2295SStefano Zampini       PetscMPIInt remote_rank, remote_count = 0;
1173835f2295SStefano Zampini 
1174835f2295SStefano Zampini       PetscCall(PetscMPIIntCast(ZCodeFind(z, size + 1, layout.zstarts), &remote_rank));
11753e72e933SJed Brown       if (remote_rank < 0) remote_rank = -(remote_rank + 1) - 1;
11763e72e933SJed Brown       // We have a new remote rank; find all the ghost indices (which are contiguous in vert_z)
11773e72e933SJed Brown 
11783e72e933SJed Brown       // Count the elements on remote_rank
11794e2e9504SJed Brown       PetscInt remote_elem = ZLayoutElementsOnRank(&layout, remote_rank);
11803e72e933SJed Brown 
11813e72e933SJed Brown       // Traverse vertices and make ghost links
11823e72e933SJed Brown       for (ZCode rz = layout.zstarts[remote_rank]; rz < layout.zstarts[remote_rank + 1]; rz++) {
11833e72e933SJed Brown         Ijk loc = ZCodeSplit(rz);
11843e72e933SJed Brown         if (rz == z) {
11853e72e933SJed Brown           local_ghosts[i] = local_elems + owned_verts + i;
11863e72e933SJed Brown           ghosts[i].rank  = remote_rank;
11873e72e933SJed Brown           ghosts[i].index = remote_elem + remote_count;
11883e72e933SJed Brown           i++;
11893e72e933SJed Brown           if (i == num_ghosts) break;
11903e72e933SJed Brown           z = vert_z[owned_verts + i];
11913e72e933SJed Brown         }
11923e72e933SJed Brown         if (IjkActive(layout.vextent, loc)) remote_count++;
11936547ddbdSJed Brown         else rz += ZStepOct(rz);
11943e72e933SJed Brown       }
11953e72e933SJed Brown     }
11963e72e933SJed Brown     PetscCall(PetscSFSetGraph(sf, local_elems + local_verts, num_ghosts, local_ghosts, PETSC_OWN_POINTER, ghosts, PETSC_OWN_POINTER));
11973e72e933SJed Brown     PetscCall(PetscObjectSetName((PetscObject)sf, "SFC Point SF"));
11983e72e933SJed Brown     PetscCall(DMSetPointSF(dm, sf));
11993e72e933SJed Brown     PetscCall(PetscSFDestroy(&sf));
12003e72e933SJed Brown   }
12013e72e933SJed Brown   {
12023e72e933SJed Brown     Vec          coordinates;
12033e72e933SJed Brown     PetscScalar *coords;
12043e72e933SJed Brown     PetscSection coord_section;
12053e72e933SJed Brown     PetscInt     coord_size;
12063e72e933SJed Brown     PetscCall(DMGetCoordinateSection(dm, &coord_section));
12073e72e933SJed Brown     PetscCall(PetscSectionSetNumFields(coord_section, 1));
12083e72e933SJed Brown     PetscCall(PetscSectionSetFieldComponents(coord_section, 0, dim));
12093e72e933SJed Brown     PetscCall(PetscSectionSetChart(coord_section, local_elems, local_elems + local_verts));
12103e72e933SJed Brown     for (PetscInt v = 0; v < local_verts; v++) {
12113e72e933SJed Brown       PetscInt point = local_elems + v;
12123e72e933SJed Brown       PetscCall(PetscSectionSetDof(coord_section, point, dim));
12133e72e933SJed Brown       PetscCall(PetscSectionSetFieldDof(coord_section, point, 0, dim));
12143e72e933SJed Brown     }
12153e72e933SJed Brown     PetscCall(PetscSectionSetUp(coord_section));
12163e72e933SJed Brown     PetscCall(PetscSectionGetStorageSize(coord_section, &coord_size));
12173e72e933SJed Brown     PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates));
12183e72e933SJed Brown     PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
12193e72e933SJed Brown     PetscCall(VecSetSizes(coordinates, coord_size, PETSC_DETERMINE));
12203e72e933SJed Brown     PetscCall(VecSetBlockSize(coordinates, dim));
12213e72e933SJed Brown     PetscCall(VecSetType(coordinates, VECSTANDARD));
12223e72e933SJed Brown     PetscCall(VecGetArray(coordinates, &coords));
12233e72e933SJed Brown     for (PetscInt v = 0; v < local_verts; v++) {
12243e72e933SJed Brown       Ijk loc             = ZCodeSplit(vert_z[v]);
12253e72e933SJed Brown       coords[v * dim + 0] = lower[0] + loc.i * (upper[0] - lower[0]) / layout.eextent.i;
12263e72e933SJed Brown       if (dim > 1) coords[v * dim + 1] = lower[1] + loc.j * (upper[1] - lower[1]) / layout.eextent.j;
12273e72e933SJed Brown       if (dim > 2) coords[v * dim + 2] = lower[2] + loc.k * (upper[2] - lower[2]) / layout.eextent.k;
12283e72e933SJed Brown     }
12293e72e933SJed Brown     PetscCall(VecRestoreArray(coordinates, &coords));
12303e72e933SJed Brown     PetscCall(DMSetCoordinatesLocal(dm, coordinates));
12313e72e933SJed Brown     PetscCall(VecDestroy(&coordinates));
12323e72e933SJed Brown   }
12333e72e933SJed Brown   if (interpolate) {
12343431e603SJed Brown     PetscCall(DMPlexInterpolateInPlace_Internal(dm));
12353431e603SJed Brown 
12363431e603SJed Brown     DMLabel label;
12373431e603SJed Brown     PetscCall(DMCreateLabel(dm, "Face Sets"));
12383431e603SJed Brown     PetscCall(DMGetLabel(dm, "Face Sets", &label));
12399ae3492bSJames Wright     PetscSegBuffer per_faces[3], donor_face_closure[3], my_donor_faces[3];
12409ae3492bSJames Wright     for (PetscInt i = 0; i < 3; i++) {
12419ae3492bSJames Wright       PetscCall(PetscSegBufferCreate(sizeof(PetscInt), 64, &per_faces[i]));
12429ae3492bSJames Wright       PetscCall(PetscSegBufferCreate(sizeof(PetscInt), 64, &my_donor_faces[i]));
12439ae3492bSJames Wright       PetscCall(PetscSegBufferCreate(sizeof(ZCode), 64 * PetscPowInt(2, dim), &donor_face_closure[i]));
12449ae3492bSJames Wright     }
12453431e603SJed Brown     PetscInt fStart, fEnd, vStart, vEnd;
12463431e603SJed Brown     PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
12473431e603SJed Brown     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
12483431e603SJed Brown     for (PetscInt f = fStart; f < fEnd; f++) {
12494e2e9504SJed Brown       PetscInt npoints, *points = NULL, num_fverts = 0, fverts[8];
12503431e603SJed Brown       PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &npoints, &points));
12513431e603SJed Brown       PetscInt bc_count[6] = {0};
12523431e603SJed Brown       for (PetscInt i = 0; i < npoints; i++) {
12533431e603SJed Brown         PetscInt p = points[2 * i];
125485de0153SJames Wright         if (!IsPointInsideStratum(p, vStart, vEnd)) continue;
12554e2e9504SJed Brown         fverts[num_fverts++] = p;
12563431e603SJed Brown         Ijk loc              = ZCodeSplit(vert_z[p - vStart]);
12573431e603SJed Brown         // Convention here matches DMPlexCreateCubeMesh_Internal
12583431e603SJed Brown         bc_count[0] += loc.i == 0;
12593431e603SJed Brown         bc_count[1] += loc.i == layout.vextent.i - 1;
12603431e603SJed Brown         bc_count[2] += loc.j == 0;
12613431e603SJed Brown         bc_count[3] += loc.j == layout.vextent.j - 1;
12623431e603SJed Brown         bc_count[4] += loc.k == 0;
12633431e603SJed Brown         bc_count[5] += loc.k == layout.vextent.k - 1;
12643e72e933SJed Brown       }
12653431e603SJed Brown       PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &npoints, &points));
12663431e603SJed Brown       for (PetscInt bc = 0, bc_match = 0; bc < 2 * dim; bc++) {
12673431e603SJed Brown         if (bc_count[bc] == PetscPowInt(2, dim - 1)) {
12686725e60dSJed Brown           PetscCall(DMPlexOrientPoint(dm, f, face_orient_dim[dim][bc]));
12694e2e9504SJed Brown           if (periodicity[bc / 2] == DM_BOUNDARY_PERIODIC) {
12704e2e9504SJed Brown             PetscInt *put;
12714e2e9504SJed Brown             if (bc % 2 == 0) { // donor face; no label
12729ae3492bSJames Wright               PetscCall(PetscSegBufferGet(my_donor_faces[bc / 2], 1, &put));
12734e2e9504SJed Brown               *put = f;
12744e2e9504SJed Brown             } else { // periodic face
12759ae3492bSJames Wright               PetscCall(PetscSegBufferGet(per_faces[bc / 2], 1, &put));
12764e2e9504SJed Brown               *put = f;
12774e2e9504SJed Brown               ZCode *zput;
12789ae3492bSJames Wright               PetscCall(PetscSegBufferGet(donor_face_closure[bc / 2], num_fverts, &zput));
12794e2e9504SJed Brown               for (PetscInt i = 0; i < num_fverts; i++) {
12804e2e9504SJed Brown                 Ijk loc = ZCodeSplit(vert_z[fverts[i] - vStart]);
12814e2e9504SJed Brown                 switch (bc / 2) {
12824e2e9504SJed Brown                 case 0:
12834e2e9504SJed Brown                   loc.i = 0;
12844e2e9504SJed Brown                   break;
12854e2e9504SJed Brown                 case 1:
12864e2e9504SJed Brown                   loc.j = 0;
12874e2e9504SJed Brown                   break;
12884e2e9504SJed Brown                 case 2:
12894e2e9504SJed Brown                   loc.k = 0;
12904e2e9504SJed Brown                   break;
12914e2e9504SJed Brown                 }
12924e2e9504SJed Brown                 *zput++ = ZEncode(loc);
12934e2e9504SJed Brown               }
12944e2e9504SJed Brown             }
12954e2e9504SJed Brown             continue;
12964e2e9504SJed Brown           }
12973431e603SJed Brown           PetscAssert(bc_match == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face matches multiple face sets");
12983431e603SJed Brown           PetscCall(DMLabelSetValue(label, f, face_marker_dim[dim][bc]));
12993431e603SJed Brown           bc_match++;
13003431e603SJed Brown         }
13013431e603SJed Brown       }
13023431e603SJed Brown     }
13033431e603SJed Brown     // Ensure that the Coordinate DM has our new boundary labels
13043431e603SJed Brown     DM cdm;
13053431e603SJed Brown     PetscCall(DMGetCoordinateDM(dm, &cdm));
13063431e603SJed Brown     PetscCall(DMCopyLabels(dm, cdm, PETSC_COPY_VALUES, PETSC_FALSE, DM_COPY_LABELS_FAIL));
13076725e60dSJed Brown     if (periodicity[0] == DM_BOUNDARY_PERIODIC || (dim > 1 && periodicity[1] == DM_BOUNDARY_PERIODIC) || (dim > 2 && periodicity[2] == DM_BOUNDARY_PERIODIC)) {
13085f06a3ddSJed Brown       PetscCall(DMPlexCreateBoxMesh_Tensor_SFC_Periodicity_Private(dm, &layout, vert_z, per_faces, lower, upper, periodicity, donor_face_closure, my_donor_faces));
13096725e60dSJed Brown     }
13109ae3492bSJames Wright     for (PetscInt i = 0; i < 3; i++) {
13119ae3492bSJames Wright       PetscCall(PetscSegBufferDestroy(&per_faces[i]));
13129ae3492bSJames Wright       PetscCall(PetscSegBufferDestroy(&donor_face_closure[i]));
13139ae3492bSJames Wright       PetscCall(PetscSegBufferDestroy(&my_donor_faces[i]));
13149ae3492bSJames Wright     }
13153431e603SJed Brown   }
13164e2e9504SJed Brown   PetscCall(PetscFree(layout.zstarts));
13173431e603SJed Brown   PetscCall(PetscFree(vert_z));
1318ea7b3863SJames Wright   PetscCall(PetscLogEventEnd(DMPLEX_CreateBoxSFC, dm, 0, 0, 0));
13193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13203e72e933SJed Brown }
13214e2e9504SJed Brown 
13224e2e9504SJed Brown /*@
13235f06a3ddSJed Brown   DMPlexSetIsoperiodicFaceSF - Express periodicity from an existing mesh
13244e2e9504SJed Brown 
132520f4b53cSBarry Smith   Logically Collective
13264e2e9504SJed Brown 
13274e2e9504SJed Brown   Input Parameters:
13284e2e9504SJed Brown + dm           - The `DMPLEX` on which to set periodicity
13291fca310dSJames Wright . num_face_sfs - Number of `PetscSF`s in `face_sfs`
13301fca310dSJames 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.
13314e2e9504SJed Brown 
13324e2e9504SJed Brown   Level: advanced
13334e2e9504SJed Brown 
133453c0d4aeSBarry Smith   Note:
13355dca41c3SJed Brown   One can use `-dm_plex_shape zbox` to use this mode of periodicity, wherein the periodic points are distinct both globally
13364e2e9504SJed Brown   and locally, but are paired when creating a global dof space.
13374e2e9504SJed Brown 
13381cc06b55SBarry Smith .seealso: [](ch_unstructured), `DMPLEX`, `DMGetGlobalSection()`, `DMPlexGetIsoperiodicFaceSF()`
13394e2e9504SJed Brown @*/
13401fca310dSJames Wright PetscErrorCode DMPlexSetIsoperiodicFaceSF(DM dm, PetscInt num_face_sfs, PetscSF *face_sfs)
13414e2e9504SJed Brown {
13424e2e9504SJed Brown   DM_Plex *plex = (DM_Plex *)dm->data;
13434d86920dSPierre Jolivet 
13444e2e9504SJed Brown   PetscFunctionBegin;
13454e2e9504SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1346d7d2d1d2SJames Wright   if (num_face_sfs) PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMGetIsoperiodicPointSF_C", DMGetIsoperiodicPointSF_Plex));
13472b4f33d9SJames Wright   else PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMGetIsoperiodicPointSF_C", NULL));
13482b4f33d9SJames Wright   if (num_face_sfs == plex->periodic.num_face_sfs && (num_face_sfs == 0 || face_sfs == plex->periodic.face_sfs)) PetscFunctionReturn(PETSC_SUCCESS);
13492b4f33d9SJames Wright   PetscCall(DMSetGlobalSection(dm, NULL));
13501fca310dSJames Wright 
13511fca310dSJames Wright   for (PetscInt i = 0; i < num_face_sfs; i++) PetscCall(PetscObjectReference((PetscObject)face_sfs[i]));
13521fca310dSJames Wright 
13531fca310dSJames Wright   if (plex->periodic.num_face_sfs > 0) {
13541fca310dSJames Wright     for (PetscInt i = 0; i < plex->periodic.num_face_sfs; i++) PetscCall(PetscSFDestroy(&plex->periodic.face_sfs[i]));
13551fca310dSJames Wright     PetscCall(PetscFree(plex->periodic.face_sfs));
13561fca310dSJames Wright   }
13571fca310dSJames Wright 
13581fca310dSJames Wright   plex->periodic.num_face_sfs = num_face_sfs;
13591fca310dSJames Wright   PetscCall(PetscCalloc1(num_face_sfs, &plex->periodic.face_sfs));
13601fca310dSJames Wright   for (PetscInt i = 0; i < num_face_sfs; i++) plex->periodic.face_sfs[i] = face_sfs[i];
13615f06a3ddSJed Brown 
13625f06a3ddSJed Brown   DM cdm = dm->coordinates[0].dm; // Can't DMGetCoordinateDM because it automatically creates one
13635f06a3ddSJed Brown   if (cdm) {
13641fca310dSJames Wright     PetscCall(DMPlexSetIsoperiodicFaceSF(cdm, num_face_sfs, face_sfs));
13651fca310dSJames Wright     if (face_sfs) cdm->periodic.setup = DMPeriodicCoordinateSetUp_Internal;
13665f06a3ddSJed Brown   }
13673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13684e2e9504SJed Brown }
13694e2e9504SJed Brown 
13701fca310dSJames Wright /*@C
13715f06a3ddSJed Brown   DMPlexGetIsoperiodicFaceSF - Obtain periodicity for a mesh
13724e2e9504SJed Brown 
137320f4b53cSBarry Smith   Logically Collective
13744e2e9504SJed Brown 
137520f4b53cSBarry Smith   Input Parameter:
13764e2e9504SJed Brown . dm - The `DMPLEX` for which to obtain periodic relation
13774e2e9504SJed Brown 
13789cde84edSJose E. Roman   Output Parameters:
13791fca310dSJames Wright + num_face_sfs - Number of `PetscSF`s in the array
13801fca310dSJames 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.
13814e2e9504SJed Brown 
13824e2e9504SJed Brown   Level: advanced
13834e2e9504SJed Brown 
13841cc06b55SBarry Smith .seealso: [](ch_unstructured), `DMPLEX`, `DMGetGlobalSection()`, `DMPlexSetIsoperiodicFaceSF()`
13854e2e9504SJed Brown @*/
13861fca310dSJames Wright PetscErrorCode DMPlexGetIsoperiodicFaceSF(DM dm, PetscInt *num_face_sfs, const PetscSF **face_sfs)
13874e2e9504SJed Brown {
13884e2e9504SJed Brown   DM_Plex *plex = (DM_Plex *)dm->data;
13894d86920dSPierre Jolivet 
13904e2e9504SJed Brown   PetscFunctionBegin;
13914e2e9504SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
13921fca310dSJames Wright   *face_sfs     = plex->periodic.face_sfs;
13931fca310dSJames Wright   *num_face_sfs = plex->periodic.num_face_sfs;
13943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13954e2e9504SJed Brown }
13966725e60dSJed Brown 
13975f06a3ddSJed Brown /*@C
13985f06a3ddSJed Brown   DMPlexSetIsoperiodicFaceTransform - set geometric transform from donor to periodic points
13996725e60dSJed Brown 
14006725e60dSJed Brown   Logically Collective
14016725e60dSJed Brown 
140260225df5SJacob Faibussowitsch   Input Parameters:
140353c0d4aeSBarry Smith + dm - `DMPLEX` that has been configured with `DMPlexSetIsoperiodicFaceSF()`
14041fca310dSJames Wright . n  - Number of transforms in array
14051fca310dSJames Wright - t  - Array of 4x4 affine transformation basis.
14066725e60dSJed Brown 
140753c0d4aeSBarry Smith   Level: advanced
140853c0d4aeSBarry Smith 
14095f06a3ddSJed Brown   Notes:
14105f06a3ddSJed Brown   Affine transforms are 4x4 matrices in which the leading 3x3 block expresses a rotation (or identity for no rotation),
14115f06a3ddSJed Brown   the last column contains a translation vector, and the bottom row is all zero except the last entry, which must always
14125f06a3ddSJed Brown   be 1. This representation is common in geometric modeling and allows affine transformations to be composed using
1413aaa8cc7dSPierre Jolivet   simple matrix multiplication.
14145f06a3ddSJed Brown 
14155f06a3ddSJed Brown   Although the interface accepts a general affine transform, only affine translation is supported at present.
14165f06a3ddSJed Brown 
141760225df5SJacob Faibussowitsch   Developer Notes:
14185f06a3ddSJed Brown   This interface should be replaced by making BasisTransform public, expanding it to support affine representations, and
14195f06a3ddSJed Brown   adding GPU implementations to apply the G2L/L2G transforms.
142053c0d4aeSBarry Smith 
14211cc06b55SBarry Smith .seealso: [](ch_unstructured), `DMPLEX`, `DMGetGlobalSection()`, `DMPlexSetIsoperiodicFaceSF()`
14226725e60dSJed Brown @*/
1423cc4c1da9SBarry Smith PetscErrorCode DMPlexSetIsoperiodicFaceTransform(DM dm, PetscInt n, const PetscScalar t[])
14246725e60dSJed Brown {
14256725e60dSJed Brown   DM_Plex *plex = (DM_Plex *)dm->data;
14264d86920dSPierre Jolivet 
14276725e60dSJed Brown   PetscFunctionBegin;
14286725e60dSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
14291fca310dSJames 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);
14301fca310dSJames Wright 
1431d7d2d1d2SJames Wright   PetscCall(PetscFree(plex->periodic.transform));
14321fca310dSJames Wright   PetscCall(PetscMalloc1(n, &plex->periodic.transform));
14331fca310dSJames Wright   for (PetscInt i = 0; i < n; i++) {
14346725e60dSJed Brown     for (PetscInt j = 0; j < 4; j++) {
14351fca310dSJames Wright       for (PetscInt k = 0; k < 4; k++) {
14361fca310dSJames Wright         PetscCheck(j != k || t[i * 16 + j * 4 + k] == 1., PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Rotated transforms not supported");
14371fca310dSJames Wright         plex->periodic.transform[i][j][k] = t[i * 16 + j * 4 + k];
14381fca310dSJames Wright       }
14396725e60dSJed Brown     }
14406725e60dSJed Brown   }
14413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14426725e60dSJed Brown }
1443