xref: /petsc/src/dm/impls/plex/plexsfc.c (revision 5f06a3dd92e34bdfe6112242961ce57bf1b830e2)
13e72e933SJed Brown #include <petsc/private/dmpleximpl.h> /*I   "petscdmplex.h"   I*/
23e72e933SJed Brown #include <petscsf.h>
33e72e933SJed Brown #include <petsc/private/hashset.h>
43e72e933SJed Brown 
53e72e933SJed Brown typedef uint64_t ZCode;
63e72e933SJed Brown 
73e72e933SJed Brown PETSC_HASH_SET(ZSet, ZCode, PetscHash_UInt64, PetscHashEqual)
83e72e933SJed Brown 
93e72e933SJed Brown typedef struct {
103e72e933SJed Brown   PetscInt i, j, k;
113e72e933SJed Brown } Ijk;
123e72e933SJed Brown 
133e72e933SJed Brown typedef struct {
143e72e933SJed Brown   Ijk         eextent;
153e72e933SJed Brown   Ijk         vextent;
163e72e933SJed Brown   PetscMPIInt comm_size;
173e72e933SJed Brown   ZCode      *zstarts;
183e72e933SJed Brown } ZLayout;
193e72e933SJed Brown 
203e72e933SJed Brown static unsigned ZCodeSplit1(ZCode z)
213e72e933SJed Brown {
223e72e933SJed Brown   z = ((z & 01001001001001001) | ((z >> 2) & 02002002002002002) | ((z >> 4) & 04004004004004004));
233e72e933SJed Brown   z = (z | (z >> 6) | (z >> 12)) & 0000000777000000777;
243e72e933SJed Brown   z = (z | (z >> 18)) & 0777777;
253e72e933SJed Brown   return (unsigned)z;
263e72e933SJed Brown }
273e72e933SJed Brown 
283e72e933SJed Brown static ZCode ZEncode1(unsigned t)
293e72e933SJed Brown {
303e72e933SJed Brown   ZCode z = t;
313e72e933SJed Brown   z       = (z | (z << 18)) & 0777000000777;
323e72e933SJed Brown   z       = (z | (z << 6) | (z << 12)) & 07007007007007007;
333e72e933SJed Brown   z       = (z | (z << 2) | (z << 4)) & 0111111111111111111;
343e72e933SJed Brown   return z;
353e72e933SJed Brown }
363e72e933SJed Brown 
373e72e933SJed Brown static Ijk ZCodeSplit(ZCode z)
383e72e933SJed Brown {
393e72e933SJed Brown   Ijk c;
403e72e933SJed Brown   c.i = ZCodeSplit1(z >> 2);
413e72e933SJed Brown   c.j = ZCodeSplit1(z >> 1);
423e72e933SJed Brown   c.k = ZCodeSplit1(z >> 0);
433e72e933SJed Brown   return c;
443e72e933SJed Brown }
453e72e933SJed Brown 
463e72e933SJed Brown static ZCode ZEncode(Ijk c)
473e72e933SJed Brown {
483e72e933SJed Brown   ZCode z = (ZEncode1(c.i) << 2) | (ZEncode1(c.j) << 1) | ZEncode1(c.k);
493e72e933SJed Brown   return z;
503e72e933SJed Brown }
513e72e933SJed Brown 
523e72e933SJed Brown static PetscBool IjkActive(Ijk extent, Ijk l)
533e72e933SJed Brown {
543e72e933SJed Brown   if (l.i < extent.i && l.j < extent.j && l.k < extent.k) return PETSC_TRUE;
553e72e933SJed Brown   return PETSC_FALSE;
563e72e933SJed Brown }
573e72e933SJed Brown 
583e72e933SJed Brown // Since element/vertex box extents are typically not equal powers of 2, Z codes that lie within the domain are not contiguous.
59*5f06a3ddSJed Brown static PetscErrorCode ZLayoutCreate(PetscMPIInt size, const PetscInt eextent[3], const PetscInt vextent[3], ZLayout *layout)
603e72e933SJed Brown {
61c77877e3SJed Brown   PetscFunctionBegin;
62*5f06a3ddSJed Brown   layout->eextent.i = eextent[0];
63*5f06a3ddSJed Brown   layout->eextent.j = eextent[1];
64*5f06a3ddSJed Brown   layout->eextent.k = eextent[2];
65*5f06a3ddSJed Brown   layout->vextent.i = vextent[0];
66*5f06a3ddSJed Brown   layout->vextent.j = vextent[1];
67*5f06a3ddSJed Brown   layout->vextent.k = vextent[2];
68*5f06a3ddSJed Brown   layout->comm_size = size;
69*5f06a3ddSJed Brown   layout->zstarts   = NULL;
70*5f06a3ddSJed Brown   PetscCall(PetscMalloc1(size + 1, &layout->zstarts));
713e72e933SJed Brown 
723e72e933SJed Brown   PetscInt total_elems = eextent[0] * eextent[1] * eextent[2];
733e72e933SJed Brown   ZCode    z           = 0;
74*5f06a3ddSJed Brown   layout->zstarts[0]   = 0;
753e72e933SJed Brown   for (PetscMPIInt r = 0; r < size; r++) {
763e72e933SJed Brown     PetscInt elems_needed = (total_elems / size) + (total_elems % size > r), count;
773e72e933SJed Brown     for (count = 0; count < elems_needed; z++) {
783e72e933SJed Brown       Ijk loc = ZCodeSplit(z);
79*5f06a3ddSJed Brown       if (IjkActive(layout->eextent, loc)) count++;
803e72e933SJed Brown     }
813e72e933SJed Brown     // Pick up any extra vertices in the Z ordering before the next rank's first owned element.
82c77877e3SJed Brown     //
83*5f06a3ddSJed Brown     // This leads to poorly balanced vertices when eextent is a power of 2, since all the fringe vertices end up
84c77877e3SJed Brown     // on the last rank. A possible solution is to balance the Z-order vertices independently from the cells, which will
85c77877e3SJed Brown     // result in a lot of element closures being remote. We could finish marking boundary conditions, then do a round of
86c77877e3SJed Brown     // vertex ownership smoothing (which would reorder and redistribute vertices without touching element distribution).
87c77877e3SJed Brown     // Another would be to have an analytic ownership criteria for vertices in the fringe veextent - eextent. This would
88c77877e3SJed Brown     // complicate the job of identifying an owner and its offset.
89*5f06a3ddSJed Brown     //
90*5f06a3ddSJed Brown     // The current recommended approach is to let `-dm_distribute 1` (default) resolve vertex ownership. This is
91*5f06a3ddSJed Brown     // *mandatory* with isoperiodicity (except in special cases) to remove standed vertices from local spaces. Here's
92*5f06a3ddSJed Brown     // the issue:
93*5f06a3ddSJed Brown     //
94*5f06a3ddSJed Brown     // Consider this partition on rank 0 (left) and rank 1.
95*5f06a3ddSJed Brown     //
96*5f06a3ddSJed Brown     //    4 --------  5 -- 14 --10 -- 21 --11
97*5f06a3ddSJed Brown     //                |          |          |
98*5f06a3ddSJed Brown     // 7 -- 16 --  8  |          |          |
99*5f06a3ddSJed Brown     // |           |  3 -------  7 -------  9
100*5f06a3ddSJed Brown     // |           |             |          |
101*5f06a3ddSJed Brown     // 4 --------  6 ------ 10   |          |
102*5f06a3ddSJed Brown     // |           |         |   6 -- 16 -- 8
103*5f06a3ddSJed Brown     // |           |         |
104*5f06a3ddSJed Brown     // 3 ---11---  5 --18--  9
105*5f06a3ddSJed Brown     //
106*5f06a3ddSJed Brown     // The periodic face SF looks like
107*5f06a3ddSJed Brown     // [0] Number of roots=21, leaves=1, remote ranks=1
108*5f06a3ddSJed Brown     // [0] 16 <- (0,11)
109*5f06a3ddSJed Brown     // [1] Number of roots=22, leaves=2, remote ranks=2
110*5f06a3ddSJed Brown     // [1] 14 <- (0,18)
111*5f06a3ddSJed Brown     // [1] 21 <- (1,16)
112*5f06a3ddSJed Brown     //
113*5f06a3ddSJed 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
114*5f06a3ddSJed Brown     // the point SF links to (1,4) and (1,5). Rank 1 learns about the periodic mapping of (1,5) while handling face
115*5f06a3ddSJed Brown     // (1,14), but never learns that vertex (1,4) has been mapped to (0,3) by face (0,16).
116*5f06a3ddSJed Brown     //
117*5f06a3ddSJed Brown     // We can relatively easily inform vertex (1,4) of this mapping, but it stays in rank 1's local space despite not
118*5f06a3ddSJed Brown     // being in the closure and thus not being contributed to. This would be mostly harmless except that some viewer
119*5f06a3ddSJed Brown     // routines expect all local points to be somehow significant. It is not easy to analytically remove the (1,4)
120*5f06a3ddSJed Brown     // vertex because the point SF and isoperiodic face SF would need to be updated to account for removal of the
121*5f06a3ddSJed Brown     // stranded vertices.
122*5f06a3ddSJed Brown     for (; z <= ZEncode(layout->vextent); z++) {
1233e72e933SJed Brown       Ijk loc = ZCodeSplit(z);
124*5f06a3ddSJed Brown       if (IjkActive(layout->eextent, loc)) break;
1253e72e933SJed Brown     }
126*5f06a3ddSJed Brown     layout->zstarts[r + 1] = z;
1273e72e933SJed Brown   }
128c77877e3SJed Brown   PetscFunctionReturn(0);
1293e72e933SJed Brown }
1303e72e933SJed Brown 
1314e2e9504SJed Brown static PetscInt ZLayoutElementsOnRank(const ZLayout *layout, PetscMPIInt rank)
1324e2e9504SJed Brown {
1334e2e9504SJed Brown   PetscInt remote_elem = 0;
1344e2e9504SJed Brown   for (ZCode rz = layout->zstarts[rank]; rz < layout->zstarts[rank + 1]; rz++) {
1354e2e9504SJed Brown     Ijk loc = ZCodeSplit(rz);
1364e2e9504SJed Brown     if (IjkActive(layout->eextent, loc)) remote_elem++;
1374e2e9504SJed Brown   }
1384e2e9504SJed Brown   return remote_elem;
1394e2e9504SJed Brown }
1404e2e9504SJed Brown 
1413e72e933SJed Brown PetscInt ZCodeFind(ZCode key, PetscInt n, const ZCode X[])
1423e72e933SJed Brown {
1433e72e933SJed Brown   PetscInt lo = 0, hi = n;
1443e72e933SJed Brown 
1453e72e933SJed Brown   if (n == 0) return -1;
1463e72e933SJed Brown   while (hi - lo > 1) {
1473e72e933SJed Brown     PetscInt mid = lo + (hi - lo) / 2;
1483e72e933SJed Brown     if (key < X[mid]) hi = mid;
1493e72e933SJed Brown     else lo = mid;
1503e72e933SJed Brown   }
1513e72e933SJed Brown   return key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
1523e72e933SJed Brown }
1533e72e933SJed Brown 
154*5f06a3ddSJed Brown static PetscErrorCode DMPlexCreateBoxMesh_Tensor_SFC_Periodicity_Private(DM dm, const ZLayout *layout, const ZCode *vert_z, PetscSegBuffer per_faces, const PetscReal *lower, const PetscReal *upper, const DMBoundaryType *periodicity, PetscSegBuffer donor_face_closure, PetscSegBuffer my_donor_faces)
1554e2e9504SJed Brown {
1564e2e9504SJed Brown   MPI_Comm     comm;
1574e2e9504SJed Brown   size_t       num_faces;
1584e2e9504SJed Brown   PetscInt     dim, *faces, vStart, vEnd;
1594e2e9504SJed Brown   PetscMPIInt  size;
1604e2e9504SJed Brown   ZCode       *donor_verts, *donor_minz;
1614e2e9504SJed Brown   PetscSFNode *leaf;
1624e2e9504SJed Brown 
1634e2e9504SJed Brown   PetscFunctionBegin;
1644e2e9504SJed Brown   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1654e2e9504SJed Brown   PetscCallMPI(MPI_Comm_size(comm, &size));
1664e2e9504SJed Brown   PetscCall(DMGetDimension(dm, &dim));
1674e2e9504SJed Brown   const PetscInt csize = PetscPowInt(2, dim - 1);
1684e2e9504SJed Brown   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1694e2e9504SJed Brown   PetscCall(PetscSegBufferGetSize(per_faces, &num_faces));
1704e2e9504SJed Brown   PetscCall(PetscSegBufferExtractInPlace(per_faces, &faces));
1714e2e9504SJed Brown   PetscCall(PetscSegBufferExtractInPlace(donor_face_closure, &donor_verts));
1724e2e9504SJed Brown   PetscCall(PetscMalloc1(num_faces, &donor_minz));
1734e2e9504SJed Brown   PetscCall(PetscMalloc1(num_faces, &leaf));
1744e2e9504SJed Brown   for (PetscInt i = 0; i < (PetscInt)num_faces; i++) {
1754e2e9504SJed Brown     ZCode minz = donor_verts[i * csize];
1764e2e9504SJed Brown     for (PetscInt j = 1; j < csize; j++) minz = PetscMin(minz, donor_verts[i * csize + j]);
1774e2e9504SJed Brown     donor_minz[i] = minz;
1784e2e9504SJed Brown   }
1794e2e9504SJed Brown   {
1804e2e9504SJed Brown     PetscBool sorted;
1814e2e9504SJed Brown     PetscCall(PetscSortedInt64(num_faces, (const PetscInt64 *)donor_minz, &sorted));
182*5f06a3ddSJed Brown     PetscCheck(sorted, PETSC_COMM_SELF, PETSC_ERR_PLIB, "minz not sorted; periodicity in multiple dimensions not yet supported");
1834e2e9504SJed Brown   }
1844e2e9504SJed Brown   for (PetscInt i = 0; i < (PetscInt)num_faces;) {
1854e2e9504SJed Brown     ZCode    z           = donor_minz[i];
1864e2e9504SJed Brown     PetscInt remote_rank = ZCodeFind(z, size + 1, layout->zstarts), remote_count = 0;
1874e2e9504SJed Brown     if (remote_rank < 0) remote_rank = -(remote_rank + 1) - 1;
1884e2e9504SJed Brown     // Process all the vertices on this rank
1894e2e9504SJed Brown     for (ZCode rz = layout->zstarts[remote_rank]; rz < layout->zstarts[remote_rank + 1]; rz++) {
1904e2e9504SJed Brown       Ijk loc = ZCodeSplit(rz);
1914e2e9504SJed Brown       if (rz == z) {
1924e2e9504SJed Brown         leaf[i].rank  = remote_rank;
1934e2e9504SJed Brown         leaf[i].index = remote_count;
1944e2e9504SJed Brown         i++;
1954e2e9504SJed Brown         if (i == (PetscInt)num_faces) break;
1964e2e9504SJed Brown         z = donor_minz[i];
1974e2e9504SJed Brown       }
1984e2e9504SJed Brown       if (IjkActive(layout->vextent, loc)) remote_count++;
1994e2e9504SJed Brown     }
2004e2e9504SJed Brown   }
2014e2e9504SJed Brown   PetscCall(PetscFree(donor_minz));
2024e2e9504SJed Brown   PetscSF sfper;
2034e2e9504SJed Brown   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), &sfper));
2044e2e9504SJed Brown   PetscCall(PetscSFSetGraph(sfper, vEnd - vStart, num_faces, PETSC_NULL, PETSC_USE_POINTER, leaf, PETSC_USE_POINTER));
2054e2e9504SJed Brown   const PetscInt *my_donor_degree;
2064e2e9504SJed Brown   PetscCall(PetscSFComputeDegreeBegin(sfper, &my_donor_degree));
2074e2e9504SJed Brown   PetscCall(PetscSFComputeDegreeEnd(sfper, &my_donor_degree));
2084e2e9504SJed Brown   PetscInt num_multiroots = 0;
2094e2e9504SJed Brown   for (PetscInt i = 0; i < vEnd - vStart; i++) {
2104e2e9504SJed Brown     num_multiroots += my_donor_degree[i];
2114e2e9504SJed Brown     if (my_donor_degree[i] == 0) continue;
212*5f06a3ddSJed Brown     PetscAssert(my_donor_degree[i] == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Local vertex has multiple faces");
2134e2e9504SJed Brown   }
2144e2e9504SJed Brown   PetscInt *my_donors, *donor_indices, *my_donor_indices;
2154e2e9504SJed Brown   size_t    num_my_donors;
2164e2e9504SJed Brown   PetscCall(PetscSegBufferGetSize(my_donor_faces, &num_my_donors));
2174e2e9504SJed Brown   PetscCheck((PetscInt)num_my_donors == num_multiroots, PETSC_COMM_SELF, PETSC_ERR_SUP, "Donor request does not match expected donors");
2184e2e9504SJed Brown   PetscCall(PetscSegBufferExtractInPlace(my_donor_faces, &my_donors));
2194e2e9504SJed Brown   PetscCall(PetscMalloc1(vEnd - vStart, &my_donor_indices));
2204e2e9504SJed Brown   for (PetscInt i = 0; i < (PetscInt)num_my_donors; i++) {
2214e2e9504SJed Brown     PetscInt f = my_donors[i];
2224e2e9504SJed Brown     PetscInt num_points, *points = NULL, minv = PETSC_MAX_INT;
2234e2e9504SJed Brown     PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &num_points, &points));
2244e2e9504SJed Brown     for (PetscInt j = 0; j < num_points; j++) {
2254e2e9504SJed Brown       PetscInt p = points[2 * j];
2264e2e9504SJed Brown       if (p < vStart || vEnd <= p) continue;
2274e2e9504SJed Brown       minv = PetscMin(minv, p);
2284e2e9504SJed Brown     }
2294e2e9504SJed Brown     PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &num_points, &points));
230*5f06a3ddSJed Brown     PetscAssert(my_donor_degree[minv - vStart] == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Local vertex not requested");
2314e2e9504SJed Brown     my_donor_indices[minv - vStart] = f;
2324e2e9504SJed Brown   }
2334e2e9504SJed Brown   PetscCall(PetscMalloc1(num_faces, &donor_indices));
2344e2e9504SJed Brown   PetscCall(PetscSFBcastBegin(sfper, MPIU_INT, my_donor_indices, donor_indices, MPI_REPLACE));
2354e2e9504SJed Brown   PetscCall(PetscSFBcastEnd(sfper, MPIU_INT, my_donor_indices, donor_indices, MPI_REPLACE));
2364e2e9504SJed Brown   PetscCall(PetscFree(my_donor_indices));
2374e2e9504SJed Brown   // Modify our leafs so they point to donor faces instead of donor minz. Additionally, give them indices as faces.
2384e2e9504SJed Brown   for (PetscInt i = 0; i < (PetscInt)num_faces; i++) leaf[i].index = donor_indices[i];
2394e2e9504SJed Brown   PetscCall(PetscFree(donor_indices));
2404e2e9504SJed Brown   PetscInt pStart, pEnd;
2414e2e9504SJed Brown   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
2424e2e9504SJed Brown   PetscCall(PetscSFSetGraph(sfper, pEnd - pStart, num_faces, faces, PETSC_COPY_VALUES, leaf, PETSC_OWN_POINTER));
243*5f06a3ddSJed Brown   PetscCall(PetscObjectSetName((PetscObject)sfper, "Z-order Isoperiodic Faces"));
2444e2e9504SJed Brown 
245*5f06a3ddSJed Brown   PetscCall(DMPlexSetIsoperiodicFaceSF(dm, sfper));
2466725e60dSJed Brown 
2476725e60dSJed Brown   PetscScalar t[4][4] = {{0}};
2486725e60dSJed Brown   t[0][0]             = 1;
2496725e60dSJed Brown   t[1][1]             = 1;
2506725e60dSJed Brown   t[2][2]             = 1;
2516725e60dSJed Brown   t[3][3]             = 1;
2526725e60dSJed Brown   for (PetscInt i = 0; i < dim; i++)
253*5f06a3ddSJed Brown     if (periodicity[i] == DM_BOUNDARY_PERIODIC) t[i][3] = upper[i] - lower[i];
254*5f06a3ddSJed Brown   PetscCall(DMPlexSetIsoperiodicFaceTransform(dm, &t[0][0]));
2554e2e9504SJed Brown   PetscCall(PetscSFDestroy(&sfper));
2564e2e9504SJed Brown   PetscFunctionReturn(0);
2574e2e9504SJed Brown }
2584e2e9504SJed Brown 
259*5f06a3ddSJed Brown // This is a DMGlobalToLocalHook that applies the affine offsets. When extended for rotated periodicity, it'll need to
260*5f06a3ddSJed Brown // apply a rotatonal transform and similar operations will be needed for fields (e.g., to rotate a velocity vector).
261*5f06a3ddSJed Brown // We use this crude approach here so we don't have to write new GPU kernels yet.
2626725e60dSJed Brown static PetscErrorCode DMCoordAddPeriodicOffsets_Private(DM dm, Vec g, InsertMode mode, Vec l, void *ctx)
2636725e60dSJed Brown {
2646725e60dSJed Brown   PetscFunctionBegin;
2656725e60dSJed Brown   PetscCall(VecScatterBegin(dm->periodic.affine_to_local, dm->periodic.affine, l, ADD_VALUES, SCATTER_FORWARD));
2666725e60dSJed Brown   PetscCall(VecScatterEnd(dm->periodic.affine_to_local, dm->periodic.affine, l, ADD_VALUES, SCATTER_FORWARD));
2676725e60dSJed Brown   PetscFunctionReturn(0);
2686725e60dSJed Brown }
2696725e60dSJed Brown 
270*5f06a3ddSJed Brown // Start with an SF for a positive depth (e.g., faces) and create a new SF for matched closure. The caller must ensure
271*5f06a3ddSJed Brown // that both the donor (root) face and the periodic (leaf) face have consistent orientation, meaning that their closures
272*5f06a3ddSJed Brown // are isomorphic. It may be useful/necessary for this restriction to be loosened.
2736725e60dSJed Brown //
274*5f06a3ddSJed Brown // Output Arguments:
275*5f06a3ddSJed Brown //
276*5f06a3ddSJed Brown // + closure_sf - augmented point SF (see `DMGetPointSF()`) that includes the faces and all points in its closure. This
277*5f06a3ddSJed Brown //   can be used to create a global section and section SF.
278*5f06a3ddSJed Brown // - is_points - index set for just the points in the closure of `face_sf`. These may be used to apply an affine
279*5f06a3ddSJed Brown //   transformation to periodic dofs; see DMPeriodicCoordinateSetUp_Internal().
280*5f06a3ddSJed Brown //
281*5f06a3ddSJed Brown static PetscErrorCode DMPlexCreateIsoperiodicPointSF_Private(DM dm, PetscSF face_sf, PetscSF *closure_sf, IS *is_points)
2826725e60dSJed Brown {
2836725e60dSJed Brown   MPI_Comm           comm;
2846725e60dSJed Brown   PetscInt           nroots, nleaves, npoints;
2856725e60dSJed Brown   const PetscInt    *filocal, *pilocal;
2866725e60dSJed Brown   const PetscSFNode *firemote, *piremote;
2876725e60dSJed Brown   PetscMPIInt        rank;
2886725e60dSJed Brown   PetscSF            point_sf;
2896725e60dSJed Brown 
2906725e60dSJed Brown   PetscFunctionBegin;
2916725e60dSJed Brown   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2926725e60dSJed Brown   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2936725e60dSJed Brown   PetscCall(PetscSFGetGraph(face_sf, &nroots, &nleaves, &filocal, &firemote));
2946725e60dSJed Brown   PetscCall(DMGetPointSF(dm, &point_sf)); // Point SF has remote points
2956725e60dSJed Brown   PetscCall(PetscSFGetGraph(point_sf, NULL, &npoints, &pilocal, &piremote));
2966725e60dSJed Brown   PetscInt *rootdata, *leafdata;
2976725e60dSJed Brown   PetscCall(PetscCalloc2(2 * nroots, &rootdata, 2 * nroots, &leafdata));
2986725e60dSJed Brown   for (PetscInt i = 0; i < nleaves; i++) {
2996725e60dSJed Brown     PetscInt point = filocal[i], cl_size, *closure = NULL;
3006725e60dSJed Brown     PetscCall(DMPlexGetTransitiveClosure(dm, point, PETSC_TRUE, &cl_size, &closure));
3016725e60dSJed Brown     leafdata[point] = cl_size - 1;
3026725e60dSJed Brown     PetscCall(DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, &cl_size, &closure));
3036725e60dSJed Brown   }
3046725e60dSJed Brown   PetscCall(PetscSFReduceBegin(face_sf, MPIU_INT, leafdata, rootdata + nroots, MPIU_SUM));
3056725e60dSJed Brown   PetscCall(PetscSFReduceEnd(face_sf, MPIU_INT, leafdata, rootdata + nroots, MPIU_SUM));
3066725e60dSJed Brown 
3076725e60dSJed Brown   PetscInt root_offset = 0;
3086725e60dSJed Brown   for (PetscInt p = 0; p < nroots; p++) {
3096725e60dSJed Brown     const PetscInt *donor_dof = rootdata + nroots;
3106725e60dSJed Brown     if (donor_dof[p] == 0) {
3116725e60dSJed Brown       rootdata[2 * p]     = -1;
3126725e60dSJed Brown       rootdata[2 * p + 1] = -1;
3136725e60dSJed Brown       continue;
3146725e60dSJed Brown     }
3156725e60dSJed Brown     PetscInt  cl_size;
3166725e60dSJed Brown     PetscInt *closure = NULL;
3176725e60dSJed Brown     PetscCall(DMPlexGetTransitiveClosure(dm, p, PETSC_TRUE, &cl_size, &closure));
3186725e60dSJed Brown     // cl_size - 1 = points not including self
319*5f06a3ddSJed Brown     PetscAssert(donor_dof[p] == cl_size - 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Reduced leaf cone sizes do not match root cone sizes");
3206725e60dSJed Brown     rootdata[2 * p]     = root_offset;
3216725e60dSJed Brown     rootdata[2 * p + 1] = cl_size - 1;
3226725e60dSJed Brown     root_offset += cl_size - 1;
3236725e60dSJed Brown     PetscCall(DMPlexRestoreTransitiveClosure(dm, p, PETSC_TRUE, &cl_size, &closure));
3246725e60dSJed Brown   }
3256725e60dSJed Brown   PetscCall(PetscSFBcastBegin(face_sf, MPIU_2INT, rootdata, leafdata, MPI_REPLACE));
3266725e60dSJed Brown   PetscCall(PetscSFBcastEnd(face_sf, MPIU_2INT, rootdata, leafdata, MPI_REPLACE));
3276725e60dSJed Brown   // Count how many leaves we need to communicate the closures
3286725e60dSJed Brown   PetscInt leaf_offset = 0;
3296725e60dSJed Brown   for (PetscInt i = 0; i < nleaves; i++) {
3306725e60dSJed Brown     PetscInt point = filocal[i];
3316725e60dSJed Brown     if (leafdata[2 * point + 1] < 0) continue;
3326725e60dSJed Brown     leaf_offset += leafdata[2 * point + 1];
3336725e60dSJed Brown   }
3346725e60dSJed Brown 
3356725e60dSJed Brown   PetscSFNode *closure_leaf;
3366725e60dSJed Brown   PetscCall(PetscMalloc1(leaf_offset, &closure_leaf));
3376725e60dSJed Brown   leaf_offset = 0;
3386725e60dSJed Brown   for (PetscInt i = 0; i < nleaves; i++) {
3396725e60dSJed Brown     PetscInt point   = filocal[i];
3406725e60dSJed Brown     PetscInt cl_size = leafdata[2 * point + 1];
3416725e60dSJed Brown     if (cl_size < 0) continue;
3426725e60dSJed Brown     for (PetscInt j = 0; j < cl_size; j++) {
3436725e60dSJed Brown       closure_leaf[leaf_offset].rank  = firemote[i].rank;
3446725e60dSJed Brown       closure_leaf[leaf_offset].index = leafdata[2 * point] + j;
3456725e60dSJed Brown       leaf_offset++;
3466725e60dSJed Brown     }
3476725e60dSJed Brown   }
3486725e60dSJed Brown 
3496725e60dSJed Brown   PetscSF sf_closure;
3506725e60dSJed Brown   PetscCall(PetscSFCreate(comm, &sf_closure));
3516725e60dSJed Brown   PetscCall(PetscSFSetGraph(sf_closure, root_offset, leaf_offset, NULL, PETSC_USE_POINTER, closure_leaf, PETSC_OWN_POINTER));
3526725e60dSJed Brown 
3536725e60dSJed Brown   // Pack root buffer with owner for every point in the root cones
3546725e60dSJed Brown   PetscSFNode *donor_closure;
3556725e60dSJed Brown   PetscCall(PetscCalloc1(root_offset, &donor_closure));
3566725e60dSJed Brown   root_offset = 0;
3576725e60dSJed Brown   for (PetscInt p = 0; p < nroots; p++) {
3586725e60dSJed Brown     if (rootdata[2 * p] < 0) continue;
3596725e60dSJed Brown     PetscInt  cl_size;
3606725e60dSJed Brown     PetscInt *closure = NULL;
3616725e60dSJed Brown     PetscCall(DMPlexGetTransitiveClosure(dm, p, PETSC_TRUE, &cl_size, &closure));
3626725e60dSJed Brown     for (PetscInt j = 1; j < cl_size; j++) {
3636725e60dSJed Brown       PetscInt c = closure[2 * j];
3646725e60dSJed Brown       if (pilocal) {
3656725e60dSJed Brown         PetscInt found = -1;
3666725e60dSJed Brown         if (npoints > 0) PetscCall(PetscFindInt(c, npoints, pilocal, &found));
3676725e60dSJed Brown         if (found >= 0) {
3686725e60dSJed Brown           donor_closure[root_offset++] = piremote[found];
3696725e60dSJed Brown           continue;
3706725e60dSJed Brown         }
3716725e60dSJed Brown       }
3726725e60dSJed Brown       // we own c
3736725e60dSJed Brown       donor_closure[root_offset].rank  = rank;
3746725e60dSJed Brown       donor_closure[root_offset].index = c;
3756725e60dSJed Brown       root_offset++;
3766725e60dSJed Brown     }
3776725e60dSJed Brown     PetscCall(DMPlexRestoreTransitiveClosure(dm, p, PETSC_TRUE, &cl_size, &closure));
3786725e60dSJed Brown   }
3796725e60dSJed Brown 
3806725e60dSJed Brown   PetscSFNode *leaf_donor_closure;
3816725e60dSJed Brown   PetscCall(PetscMalloc1(leaf_offset, &leaf_donor_closure));
3826725e60dSJed Brown   PetscCall(PetscSFBcastBegin(sf_closure, MPIU_2INT, donor_closure, leaf_donor_closure, MPI_REPLACE));
3836725e60dSJed Brown   PetscCall(PetscSFBcastEnd(sf_closure, MPIU_2INT, donor_closure, leaf_donor_closure, MPI_REPLACE));
3846725e60dSJed Brown   PetscCall(PetscSFDestroy(&sf_closure));
3856725e60dSJed Brown   PetscCall(PetscFree(donor_closure));
3866725e60dSJed Brown 
3876725e60dSJed Brown   PetscSFNode *new_iremote;
3886725e60dSJed Brown   PetscCall(PetscCalloc1(nroots, &new_iremote));
3896725e60dSJed Brown   for (PetscInt i = 0; i < nroots; i++) new_iremote[i].rank = -1;
3906725e60dSJed Brown   // Walk leaves and match vertices
3916725e60dSJed Brown   leaf_offset = 0;
3926725e60dSJed Brown   for (PetscInt i = 0; i < nleaves; i++) {
3936725e60dSJed Brown     PetscInt  point   = filocal[i], cl_size;
3946725e60dSJed Brown     PetscInt *closure = NULL;
3956725e60dSJed Brown     PetscCall(DMPlexGetTransitiveClosure(dm, point, PETSC_TRUE, &cl_size, &closure));
3966725e60dSJed Brown     for (PetscInt j = 1; j < cl_size; j++) { // TODO: should we send donor edge orientations so we can flip for consistency?
3976725e60dSJed Brown       PetscInt    c  = closure[2 * j];
3986725e60dSJed Brown       PetscSFNode lc = leaf_donor_closure[leaf_offset];
3996725e60dSJed Brown       // printf("[%d] face %d.%d: %d ?-- (%d,%d)\n", rank, point, j, c, lc.rank, lc.index);
4006725e60dSJed Brown       if (new_iremote[c].rank == -1) {
4016725e60dSJed Brown         new_iremote[c] = lc;
402*5f06a3ddSJed 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");
4036725e60dSJed Brown       leaf_offset++;
4046725e60dSJed Brown     }
4056725e60dSJed Brown     PetscCall(DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, &cl_size, &closure));
4066725e60dSJed Brown   }
4076725e60dSJed Brown   PetscCall(PetscFree(leaf_donor_closure));
4086725e60dSJed Brown 
4096725e60dSJed Brown   // Include face points in closure SF
4106725e60dSJed Brown   for (PetscInt i = 0; i < nleaves; i++) new_iremote[filocal[i]] = firemote[i];
4116725e60dSJed Brown   // consolidate leaves
4126725e60dSJed Brown   PetscInt num_new_leaves = 0;
4136725e60dSJed Brown   for (PetscInt i = 0; i < nroots; i++) {
4146725e60dSJed Brown     if (new_iremote[i].rank == -1) continue;
4156725e60dSJed Brown     new_iremote[num_new_leaves] = new_iremote[i];
4166725e60dSJed Brown     leafdata[num_new_leaves]    = i;
4176725e60dSJed Brown     num_new_leaves++;
4186725e60dSJed Brown   }
419*5f06a3ddSJed Brown   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, num_new_leaves, leafdata, PETSC_COPY_VALUES, is_points));
4206725e60dSJed Brown 
4216725e60dSJed Brown   PetscSF csf;
4226725e60dSJed Brown   PetscCall(PetscSFCreate(comm, &csf));
4236725e60dSJed Brown   PetscCall(PetscSFSetGraph(csf, nroots, num_new_leaves, leafdata, PETSC_COPY_VALUES, new_iremote, PETSC_COPY_VALUES));
4246725e60dSJed Brown   PetscCall(PetscFree(new_iremote)); // copy and delete because new_iremote is longer than it needs to be
4256725e60dSJed Brown   PetscCall(PetscFree2(rootdata, leafdata));
4266725e60dSJed Brown 
427*5f06a3ddSJed Brown   if (npoints < 0) { // empty point_sf
4286725e60dSJed Brown     *closure_sf = csf;
429*5f06a3ddSJed Brown   } else {
430*5f06a3ddSJed Brown     PetscCall(PetscSFMerge(point_sf, csf, closure_sf));
431*5f06a3ddSJed Brown     PetscCall(PetscSFDestroy(&csf));
432*5f06a3ddSJed Brown   }
433*5f06a3ddSJed Brown   PetscCall(PetscObjectSetName((PetscObject)*closure_sf, "Composed Periodic Points"));
4346725e60dSJed Brown   PetscFunctionReturn(0);
4356725e60dSJed Brown }
4366725e60dSJed Brown 
437*5f06a3ddSJed Brown static PetscErrorCode DMGetIsoperiodicPointSF_Plex(DM dm, PetscSF *sf)
4386725e60dSJed Brown {
4396725e60dSJed Brown   DM_Plex *plex = (DM_Plex *)dm->data;
4406725e60dSJed Brown 
4416725e60dSJed Brown   PetscFunctionBegin;
4426725e60dSJed Brown   if (!plex->periodic.composed_sf) {
4436725e60dSJed Brown     PetscSF face_sf = plex->periodic.face_sf;
4446725e60dSJed Brown 
445*5f06a3ddSJed Brown     PetscCall(DMPlexCreateIsoperiodicPointSF_Private(dm, face_sf, &plex->periodic.composed_sf, &plex->periodic.periodic_points));
4466725e60dSJed Brown   }
4476725e60dSJed Brown   if (sf) *sf = plex->periodic.composed_sf;
4486725e60dSJed Brown   PetscFunctionReturn(0);
4496725e60dSJed Brown }
4506725e60dSJed Brown 
451*5f06a3ddSJed Brown PetscErrorCode DMPlexMigrateIsoperiodicFaceSF_Internal(DM old_dm, DM dm, PetscSF sf_migration)
452*5f06a3ddSJed Brown {
453*5f06a3ddSJed Brown   DM_Plex    *plex = (DM_Plex *)old_dm->data;
454*5f06a3ddSJed Brown   PetscSF     sf_point;
455*5f06a3ddSJed Brown   PetscMPIInt rank;
456*5f06a3ddSJed Brown 
457*5f06a3ddSJed Brown   PetscFunctionBegin;
458*5f06a3ddSJed Brown   if (!plex->periodic.face_sf) PetscFunctionReturn(0);
459*5f06a3ddSJed Brown   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
460*5f06a3ddSJed Brown   PetscCall(DMGetPointSF(dm, &sf_point));
461*5f06a3ddSJed Brown   PetscInt           old_npoints, new_npoints, old_nleaf, new_nleaf, point_nleaf;
462*5f06a3ddSJed Brown   PetscSFNode       *new_leafdata, *rootdata, *leafdata;
463*5f06a3ddSJed Brown   const PetscInt    *old_local, *point_local;
464*5f06a3ddSJed Brown   const PetscSFNode *old_remote, *point_remote;
465*5f06a3ddSJed Brown   PetscCall(PetscSFGetGraph(plex->periodic.face_sf, &old_npoints, &old_nleaf, &old_local, &old_remote));
466*5f06a3ddSJed Brown   PetscCall(PetscSFGetGraph(sf_migration, NULL, &new_nleaf, NULL, NULL));
467*5f06a3ddSJed Brown   PetscCall(PetscSFGetGraph(sf_point, &new_npoints, &point_nleaf, &point_local, &point_remote));
468*5f06a3ddSJed Brown   PetscAssert(new_nleaf == new_npoints, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected migration leaf space to match new point root space");
469*5f06a3ddSJed Brown   PetscCall(PetscMalloc3(old_npoints, &rootdata, old_npoints, &leafdata, new_npoints, &new_leafdata));
470*5f06a3ddSJed Brown 
471*5f06a3ddSJed Brown   // Fill new_leafdata with new owners of all points
472*5f06a3ddSJed Brown   for (PetscInt i = 0; i < new_npoints; i++) {
473*5f06a3ddSJed Brown     new_leafdata[i].rank  = rank;
474*5f06a3ddSJed Brown     new_leafdata[i].index = i;
475*5f06a3ddSJed Brown   }
476*5f06a3ddSJed Brown   for (PetscInt i = 0; i < point_nleaf; i++) {
477*5f06a3ddSJed Brown     PetscInt j      = point_local[i];
478*5f06a3ddSJed Brown     new_leafdata[j] = point_remote[i];
479*5f06a3ddSJed Brown   }
480*5f06a3ddSJed Brown   // REPLACE is okay because every leaf agrees about the new owners
481*5f06a3ddSJed Brown   PetscCall(PetscSFReduceBegin(sf_migration, MPIU_2INT, new_leafdata, rootdata, MPI_REPLACE));
482*5f06a3ddSJed Brown   PetscCall(PetscSFReduceEnd(sf_migration, MPIU_2INT, new_leafdata, rootdata, MPI_REPLACE));
483*5f06a3ddSJed Brown   // rootdata now contains the new owners
484*5f06a3ddSJed Brown 
485*5f06a3ddSJed Brown   // Send to leaves of old space
486*5f06a3ddSJed Brown   for (PetscInt i = 0; i < old_npoints; i++) {
487*5f06a3ddSJed Brown     leafdata[i].rank  = -1;
488*5f06a3ddSJed Brown     leafdata[i].index = -1;
489*5f06a3ddSJed Brown   }
490*5f06a3ddSJed Brown   PetscCall(PetscSFBcastBegin(plex->periodic.face_sf, MPIU_2INT, rootdata, leafdata, MPI_REPLACE));
491*5f06a3ddSJed Brown   PetscCall(PetscSFBcastEnd(plex->periodic.face_sf, MPIU_2INT, rootdata, leafdata, MPI_REPLACE));
492*5f06a3ddSJed Brown 
493*5f06a3ddSJed Brown   // Send to new leaf space
494*5f06a3ddSJed Brown   PetscCall(PetscSFBcastBegin(sf_migration, MPIU_2INT, leafdata, new_leafdata, MPI_REPLACE));
495*5f06a3ddSJed Brown   PetscCall(PetscSFBcastEnd(sf_migration, MPIU_2INT, leafdata, new_leafdata, MPI_REPLACE));
496*5f06a3ddSJed Brown 
497*5f06a3ddSJed Brown   PetscInt     nface = 0, *new_local;
498*5f06a3ddSJed Brown   PetscSFNode *new_remote;
499*5f06a3ddSJed Brown   for (PetscInt i = 0; i < new_npoints; i++) nface += (new_leafdata[i].rank >= 0);
500*5f06a3ddSJed Brown   PetscCall(PetscMalloc1(nface, &new_local));
501*5f06a3ddSJed Brown   PetscCall(PetscMalloc1(nface, &new_remote));
502*5f06a3ddSJed Brown   nface = 0;
503*5f06a3ddSJed Brown   for (PetscInt i = 0; i < new_npoints; i++) {
504*5f06a3ddSJed Brown     if (new_leafdata[i].rank == -1) continue;
505*5f06a3ddSJed Brown     new_local[nface]  = i;
506*5f06a3ddSJed Brown     new_remote[nface] = new_leafdata[i];
507*5f06a3ddSJed Brown     nface++;
508*5f06a3ddSJed Brown   }
509*5f06a3ddSJed Brown   PetscCall(PetscFree3(rootdata, leafdata, new_leafdata));
510*5f06a3ddSJed Brown   PetscSF sf_face;
511*5f06a3ddSJed Brown   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), &sf_face));
512*5f06a3ddSJed Brown   PetscCall(PetscSFSetGraph(sf_face, new_npoints, nface, new_local, PETSC_OWN_POINTER, new_remote, PETSC_OWN_POINTER));
513*5f06a3ddSJed Brown   PetscCall(PetscObjectSetName((PetscObject)sf_face, "Migrated Isoperiodic Faces"));
514*5f06a3ddSJed Brown   PetscCall(DMPlexSetIsoperiodicFaceSF(dm, sf_face));
515*5f06a3ddSJed Brown   PetscCall(DMPlexSetIsoperiodicFaceTransform(dm, &plex->periodic.transform[0][0]));
516*5f06a3ddSJed Brown   PetscCall(PetscSFDestroy(&sf_face));
517*5f06a3ddSJed Brown   PetscFunctionReturn(0);
518*5f06a3ddSJed Brown }
519*5f06a3ddSJed Brown 
5206725e60dSJed Brown PetscErrorCode DMPeriodicCoordinateSetUp_Internal(DM dm)
5216725e60dSJed Brown {
5226725e60dSJed Brown   DM_Plex *plex = (DM_Plex *)dm->data;
5236725e60dSJed Brown   PetscFunctionBegin;
5246725e60dSJed Brown   if (!plex->periodic.face_sf) PetscFunctionReturn(0);
525*5f06a3ddSJed Brown   PetscCall(DMGetIsoperiodicPointSF_Plex(dm, NULL));
526*5f06a3ddSJed Brown   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMGetIsoperiodicPointSF_C", DMGetIsoperiodicPointSF_Plex));
5276725e60dSJed Brown 
5286725e60dSJed Brown   PetscInt dim;
5296725e60dSJed Brown   PetscCall(DMGetDimension(dm, &dim));
5306725e60dSJed Brown   size_t count;
5316725e60dSJed Brown   IS     isdof;
5326725e60dSJed Brown   {
5336725e60dSJed Brown     PetscInt        npoints;
5346725e60dSJed Brown     const PetscInt *points;
5356725e60dSJed Brown     IS              is = plex->periodic.periodic_points;
5366725e60dSJed Brown     PetscSegBuffer  seg;
5376725e60dSJed Brown     PetscSection    section;
5386725e60dSJed Brown     PetscCall(DMGetLocalSection(dm, &section));
5396725e60dSJed Brown     PetscCall(PetscSegBufferCreate(sizeof(PetscInt), 32, &seg));
5406725e60dSJed Brown     PetscCall(ISGetSize(is, &npoints));
5416725e60dSJed Brown     PetscCall(ISGetIndices(is, &points));
5426725e60dSJed Brown     for (PetscInt i = 0; i < npoints; i++) {
5436725e60dSJed Brown       PetscInt point = points[i], off, dof;
5446725e60dSJed Brown       PetscCall(PetscSectionGetOffset(section, point, &off));
5456725e60dSJed Brown       PetscCall(PetscSectionGetDof(section, point, &dof));
5466725e60dSJed Brown       PetscAssert(dof % dim == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unexpected dof %" PetscInt_FMT " not divisible by dimension %" PetscInt_FMT, dof, dim);
5476725e60dSJed Brown       for (PetscInt j = 0; j < dof / dim; j++) {
5486725e60dSJed Brown         PetscInt *slot;
5496725e60dSJed Brown         PetscCall(PetscSegBufferGetInts(seg, 1, &slot));
5506725e60dSJed Brown         *slot = off / dim;
5516725e60dSJed Brown       }
5526725e60dSJed Brown     }
5536725e60dSJed Brown     PetscInt *ind;
5546725e60dSJed Brown     PetscCall(PetscSegBufferGetSize(seg, &count));
5556725e60dSJed Brown     PetscCall(PetscSegBufferExtractAlloc(seg, &ind));
5566725e60dSJed Brown     PetscCall(PetscSegBufferDestroy(&seg));
5576725e60dSJed Brown     PetscCall(ISCreateBlock(PETSC_COMM_SELF, dim, count, ind, PETSC_OWN_POINTER, &isdof));
5586725e60dSJed Brown   }
5596725e60dSJed Brown   Vec        L, P;
5606725e60dSJed Brown   VecType    vec_type;
5616725e60dSJed Brown   VecScatter scatter;
5626725e60dSJed Brown   PetscCall(DMGetLocalVector(dm, &L));
5636725e60dSJed Brown   PetscCall(VecCreate(PETSC_COMM_SELF, &P));
5646725e60dSJed Brown   PetscCall(VecSetSizes(P, count * dim, count * dim));
5656725e60dSJed Brown   PetscCall(VecGetType(L, &vec_type));
5666725e60dSJed Brown   PetscCall(VecSetType(P, vec_type));
5676725e60dSJed Brown   PetscCall(VecScatterCreate(P, NULL, L, isdof, &scatter));
5686725e60dSJed Brown   PetscCall(DMRestoreLocalVector(dm, &L));
5696725e60dSJed Brown   PetscCall(ISDestroy(&isdof));
5706725e60dSJed Brown 
5716725e60dSJed Brown   {
5726725e60dSJed Brown     PetscScalar *x;
5736725e60dSJed Brown     PetscCall(VecGetArrayWrite(P, &x));
5746725e60dSJed Brown     for (PetscInt i = 0; i < (PetscInt)count; i++) {
5756725e60dSJed Brown       for (PetscInt j = 0; j < dim; j++) x[i * dim + j] = plex->periodic.transform[j][3];
5766725e60dSJed Brown     }
5776725e60dSJed Brown     PetscCall(VecRestoreArrayWrite(P, &x));
5786725e60dSJed Brown   }
5796725e60dSJed Brown 
5806725e60dSJed Brown   dm->periodic.affine_to_local = scatter;
5816725e60dSJed Brown   dm->periodic.affine          = P;
5826725e60dSJed Brown   PetscCall(DMGlobalToLocalHookAdd(dm, NULL, DMCoordAddPeriodicOffsets_Private, NULL));
5836725e60dSJed Brown   PetscFunctionReturn(0);
5846725e60dSJed Brown }
5856725e60dSJed Brown 
586*5f06a3ddSJed Brown // We'll just orient all the edges, though only periodic boundary edges need orientation
587*5f06a3ddSJed Brown static PetscErrorCode DMPlexOrientPositiveEdges_Private(DM dm)
588*5f06a3ddSJed Brown {
589*5f06a3ddSJed Brown   PetscInt dim, eStart, eEnd;
590*5f06a3ddSJed Brown   PetscFunctionBegin;
591*5f06a3ddSJed Brown   PetscCall(DMGetDimension(dm, &dim));
592*5f06a3ddSJed Brown   if (dim < 3) PetscFunctionReturn(0); // not necessary
593*5f06a3ddSJed Brown   PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd));
594*5f06a3ddSJed Brown   for (PetscInt e = eStart; e < eEnd; e++) {
595*5f06a3ddSJed Brown     const PetscInt *cone;
596*5f06a3ddSJed Brown     PetscCall(DMPlexGetCone(dm, e, &cone));
597*5f06a3ddSJed Brown     if (cone[0] > cone[1]) PetscCall(DMPlexOrientPoint(dm, e, -1));
598*5f06a3ddSJed Brown   }
599*5f06a3ddSJed Brown   PetscFunctionReturn(0);
600*5f06a3ddSJed Brown }
601*5f06a3ddSJed Brown 
6023e72e933SJed Brown PetscErrorCode DMPlexCreateBoxMesh_Tensor_SFC_Internal(DM dm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate)
6033e72e933SJed Brown {
6043e72e933SJed Brown   PetscInt  eextent[3] = {1, 1, 1}, vextent[3] = {1, 1, 1};
6053e72e933SJed Brown   const Ijk closure_1[] = {
6063e72e933SJed Brown     {0, 0, 0},
6073e72e933SJed Brown     {1, 0, 0},
6083e72e933SJed Brown   };
6093e72e933SJed Brown   const Ijk closure_2[] = {
6103e72e933SJed Brown     {0, 0, 0},
6113e72e933SJed Brown     {1, 0, 0},
6123e72e933SJed Brown     {1, 1, 0},
6133e72e933SJed Brown     {0, 1, 0},
6143e72e933SJed Brown   };
6153e72e933SJed Brown   const Ijk closure_3[] = {
6163e72e933SJed Brown     {0, 0, 0},
6173e72e933SJed Brown     {0, 1, 0},
6183e72e933SJed Brown     {1, 1, 0},
6193e72e933SJed Brown     {1, 0, 0},
6203e72e933SJed Brown     {0, 0, 1},
6213e72e933SJed Brown     {1, 0, 1},
6223e72e933SJed Brown     {1, 1, 1},
6233e72e933SJed Brown     {0, 1, 1},
6243e72e933SJed Brown   };
6253431e603SJed Brown   const Ijk *const closure_dim[] = {NULL, closure_1, closure_2, closure_3};
6263431e603SJed Brown   // This must be kept consistent with DMPlexCreateCubeMesh_Internal
6273431e603SJed Brown   const PetscInt        face_marker_1[]   = {1, 2};
6283431e603SJed Brown   const PetscInt        face_marker_2[]   = {4, 2, 1, 3};
6293431e603SJed Brown   const PetscInt        face_marker_3[]   = {6, 5, 3, 4, 1, 2};
6303431e603SJed Brown   const PetscInt *const face_marker_dim[] = {NULL, face_marker_1, face_marker_2, face_marker_3};
6316725e60dSJed Brown   // Orient faces so the normal is in the positive axis and the first vertex is the one closest to zero.
6326725e60dSJed Brown   // These orientations can be determined by examining cones of a reference quad and hex element.
6336725e60dSJed Brown   const PetscInt        face_orient_1[]   = {0, 0};
6346725e60dSJed Brown   const PetscInt        face_orient_2[]   = {-1, 0, 0, -1};
6356725e60dSJed Brown   const PetscInt        face_orient_3[]   = {-2, 0, -2, 1, -2, 0};
6366725e60dSJed Brown   const PetscInt *const face_orient_dim[] = {NULL, face_orient_1, face_orient_2, face_orient_3};
6373e72e933SJed Brown 
6383e72e933SJed Brown   PetscFunctionBegin;
6393e72e933SJed Brown   PetscValidPointer(dm, 1);
6403e72e933SJed Brown   PetscValidLogicalCollectiveInt(dm, dim, 2);
6413e72e933SJed Brown   PetscCall(DMSetDimension(dm, dim));
6423e72e933SJed Brown   PetscMPIInt rank, size;
6433e72e933SJed Brown   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
6443e72e933SJed Brown   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
6453e72e933SJed Brown   for (PetscInt i = 0; i < dim; i++) {
6463e72e933SJed Brown     eextent[i] = faces[i];
6473e72e933SJed Brown     vextent[i] = faces[i] + 1;
6483e72e933SJed Brown   }
649c77877e3SJed Brown   ZLayout layout;
650c77877e3SJed Brown   PetscCall(ZLayoutCreate(size, eextent, vextent, &layout));
6513e72e933SJed Brown   PetscZSet vset; // set of all vertices in the closure of the owned elements
6523e72e933SJed Brown   PetscCall(PetscZSetCreate(&vset));
6533e72e933SJed Brown   PetscInt local_elems = 0;
6543e72e933SJed Brown   for (ZCode z = layout.zstarts[rank]; z < layout.zstarts[rank + 1]; z++) {
6553e72e933SJed Brown     Ijk loc = ZCodeSplit(z);
6563e72e933SJed Brown     if (IjkActive(layout.vextent, loc)) PetscZSetAdd(vset, z);
6573e72e933SJed Brown     if (IjkActive(layout.eextent, loc)) {
6583e72e933SJed Brown       local_elems++;
6593e72e933SJed Brown       // Add all neighboring vertices to set
6603e72e933SJed Brown       for (PetscInt n = 0; n < PetscPowInt(2, dim); n++) {
6613e72e933SJed Brown         Ijk   inc  = closure_dim[dim][n];
6623e72e933SJed Brown         Ijk   nloc = {loc.i + inc.i, loc.j + inc.j, loc.k + inc.k};
6633e72e933SJed Brown         ZCode v    = ZEncode(nloc);
6643e72e933SJed Brown         PetscZSetAdd(vset, v);
6653e72e933SJed Brown       }
6663e72e933SJed Brown     }
6673e72e933SJed Brown   }
6683e72e933SJed Brown   PetscInt local_verts, off = 0;
6693e72e933SJed Brown   ZCode   *vert_z;
6703e72e933SJed Brown   PetscCall(PetscZSetGetSize(vset, &local_verts));
6713e72e933SJed Brown   PetscCall(PetscMalloc1(local_verts, &vert_z));
6723e72e933SJed Brown   PetscCall(PetscZSetGetElems(vset, &off, vert_z));
6733e72e933SJed Brown   PetscCall(PetscZSetDestroy(&vset));
6743e72e933SJed Brown   // ZCode is unsigned for bitwise convenience, but highest bit should never be set, so can interpret as signed
6753e72e933SJed Brown   PetscCall(PetscSortInt64(local_verts, (PetscInt64 *)vert_z));
6763e72e933SJed Brown 
6773e72e933SJed Brown   PetscCall(DMPlexSetChart(dm, 0, local_elems + local_verts));
6783e72e933SJed Brown   for (PetscInt e = 0; e < local_elems; e++) PetscCall(DMPlexSetConeSize(dm, e, PetscPowInt(2, dim)));
6793e72e933SJed Brown   PetscCall(DMSetUp(dm));
6803e72e933SJed Brown   {
6813e72e933SJed Brown     PetscInt e = 0;
6823e72e933SJed Brown     for (ZCode z = layout.zstarts[rank]; z < layout.zstarts[rank + 1]; z++) {
6833e72e933SJed Brown       Ijk loc = ZCodeSplit(z);
6843e72e933SJed Brown       if (!IjkActive(layout.eextent, loc)) continue;
6853e72e933SJed Brown       PetscInt cone[8], orient[8] = {0};
6863e72e933SJed Brown       for (PetscInt n = 0; n < PetscPowInt(2, dim); n++) {
6873e72e933SJed Brown         Ijk      inc  = closure_dim[dim][n];
6883e72e933SJed Brown         Ijk      nloc = {loc.i + inc.i, loc.j + inc.j, loc.k + inc.k};
6893e72e933SJed Brown         ZCode    v    = ZEncode(nloc);
6903e72e933SJed Brown         PetscInt ci   = ZCodeFind(v, local_verts, vert_z);
6913e72e933SJed Brown         PetscAssert(ci >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Did not find neighbor vertex in set");
6923e72e933SJed Brown         cone[n] = local_elems + ci;
6933e72e933SJed Brown       }
6943e72e933SJed Brown       PetscCall(DMPlexSetCone(dm, e, cone));
6953e72e933SJed Brown       PetscCall(DMPlexSetConeOrientation(dm, e, orient));
6963e72e933SJed Brown       e++;
6973e72e933SJed Brown     }
6983e72e933SJed Brown   }
6993e72e933SJed Brown 
7003e72e933SJed Brown   PetscCall(DMPlexSymmetrize(dm));
7013e72e933SJed Brown   PetscCall(DMPlexStratify(dm));
7026725e60dSJed Brown 
7033e72e933SJed Brown   { // Create point SF
7043e72e933SJed Brown     PetscSF sf;
7053e72e933SJed Brown     PetscSFCreate(PetscObjectComm((PetscObject)dm), &sf);
7063e72e933SJed Brown     PetscInt owned_verts = ZCodeFind(layout.zstarts[rank + 1], local_verts, vert_z);
7073e72e933SJed Brown     if (owned_verts < 0) owned_verts = -(owned_verts + 1); // We don't care whether the key was found
7083e72e933SJed Brown     PetscInt     num_ghosts = local_verts - owned_verts;   // Due to sorting, owned vertices always come first
7093e72e933SJed Brown     PetscInt    *local_ghosts;
7103e72e933SJed Brown     PetscSFNode *ghosts;
7113e72e933SJed Brown     PetscCall(PetscMalloc1(num_ghosts, &local_ghosts));
7123e72e933SJed Brown     PetscCall(PetscMalloc1(num_ghosts, &ghosts));
7133e72e933SJed Brown     for (PetscInt i = 0; i < num_ghosts;) {
7143e72e933SJed Brown       ZCode    z           = vert_z[owned_verts + i];
7153e72e933SJed Brown       PetscInt remote_rank = ZCodeFind(z, size + 1, layout.zstarts), remote_count = 0;
7163e72e933SJed Brown       if (remote_rank < 0) remote_rank = -(remote_rank + 1) - 1;
7173e72e933SJed Brown       // We have a new remote rank; find all the ghost indices (which are contiguous in vert_z)
7183e72e933SJed Brown 
7193e72e933SJed Brown       // Count the elements on remote_rank
7204e2e9504SJed Brown       PetscInt remote_elem = ZLayoutElementsOnRank(&layout, remote_rank);
7213e72e933SJed Brown 
7223e72e933SJed Brown       // Traverse vertices and make ghost links
7233e72e933SJed Brown       for (ZCode rz = layout.zstarts[remote_rank]; rz < layout.zstarts[remote_rank + 1]; rz++) {
7243e72e933SJed Brown         Ijk loc = ZCodeSplit(rz);
7253e72e933SJed Brown         if (rz == z) {
7263e72e933SJed Brown           local_ghosts[i] = local_elems + owned_verts + i;
7273e72e933SJed Brown           ghosts[i].rank  = remote_rank;
7283e72e933SJed Brown           ghosts[i].index = remote_elem + remote_count;
7293e72e933SJed Brown           i++;
7303e72e933SJed Brown           if (i == num_ghosts) break;
7313e72e933SJed Brown           z = vert_z[owned_verts + i];
7323e72e933SJed Brown         }
7333e72e933SJed Brown         if (IjkActive(layout.vextent, loc)) remote_count++;
7343e72e933SJed Brown       }
7353e72e933SJed Brown     }
7363e72e933SJed Brown     PetscCall(PetscSFSetGraph(sf, local_elems + local_verts, num_ghosts, local_ghosts, PETSC_OWN_POINTER, ghosts, PETSC_OWN_POINTER));
7373e72e933SJed Brown     PetscCall(PetscObjectSetName((PetscObject)sf, "SFC Point SF"));
7383e72e933SJed Brown     PetscCall(DMSetPointSF(dm, sf));
7393e72e933SJed Brown     PetscCall(PetscSFDestroy(&sf));
7403e72e933SJed Brown   }
7413e72e933SJed Brown   {
7423e72e933SJed Brown     Vec          coordinates;
7433e72e933SJed Brown     PetscScalar *coords;
7443e72e933SJed Brown     PetscSection coord_section;
7453e72e933SJed Brown     PetscInt     coord_size;
7463e72e933SJed Brown     PetscCall(DMGetCoordinateSection(dm, &coord_section));
7473e72e933SJed Brown     PetscCall(PetscSectionSetNumFields(coord_section, 1));
7483e72e933SJed Brown     PetscCall(PetscSectionSetFieldComponents(coord_section, 0, dim));
7493e72e933SJed Brown     PetscCall(PetscSectionSetChart(coord_section, local_elems, local_elems + local_verts));
7503e72e933SJed Brown     for (PetscInt v = 0; v < local_verts; v++) {
7513e72e933SJed Brown       PetscInt point = local_elems + v;
7523e72e933SJed Brown       PetscCall(PetscSectionSetDof(coord_section, point, dim));
7533e72e933SJed Brown       PetscCall(PetscSectionSetFieldDof(coord_section, point, 0, dim));
7543e72e933SJed Brown     }
7553e72e933SJed Brown     PetscCall(PetscSectionSetUp(coord_section));
7563e72e933SJed Brown     PetscCall(PetscSectionGetStorageSize(coord_section, &coord_size));
7573e72e933SJed Brown     PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates));
7583e72e933SJed Brown     PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
7593e72e933SJed Brown     PetscCall(VecSetSizes(coordinates, coord_size, PETSC_DETERMINE));
7603e72e933SJed Brown     PetscCall(VecSetBlockSize(coordinates, dim));
7613e72e933SJed Brown     PetscCall(VecSetType(coordinates, VECSTANDARD));
7623e72e933SJed Brown     PetscCall(VecGetArray(coordinates, &coords));
7633e72e933SJed Brown     for (PetscInt v = 0; v < local_verts; v++) {
7643e72e933SJed Brown       Ijk loc             = ZCodeSplit(vert_z[v]);
7653e72e933SJed Brown       coords[v * dim + 0] = lower[0] + loc.i * (upper[0] - lower[0]) / layout.eextent.i;
7663e72e933SJed Brown       if (dim > 1) coords[v * dim + 1] = lower[1] + loc.j * (upper[1] - lower[1]) / layout.eextent.j;
7673e72e933SJed Brown       if (dim > 2) coords[v * dim + 2] = lower[2] + loc.k * (upper[2] - lower[2]) / layout.eextent.k;
7683e72e933SJed Brown     }
7693e72e933SJed Brown     PetscCall(VecRestoreArray(coordinates, &coords));
7703e72e933SJed Brown     PetscCall(DMSetCoordinatesLocal(dm, coordinates));
7713e72e933SJed Brown     PetscCall(VecDestroy(&coordinates));
7723e72e933SJed Brown   }
7733e72e933SJed Brown   if (interpolate) {
7743431e603SJed Brown     PetscCall(DMPlexInterpolateInPlace_Internal(dm));
775*5f06a3ddSJed Brown     // It's currently necessary to orient the donor and periodic edges consistently. An easy way to ensure that is ot
776*5f06a3ddSJed Brown     // give all edges positive orientation. Since vertices are created in Z-order, all ranks will agree about the
777*5f06a3ddSJed Brown     // ordering cone[0] < cone[1]. This is overkill and it would be nice to remove this preparation and make
778*5f06a3ddSJed Brown     // DMPlexCreateIsoperiodicClosureSF_Private() more resilient, so it fixes any inconsistent orientations. That might
779*5f06a3ddSJed Brown     // be needed in a general CGNS reader, for example.
780*5f06a3ddSJed Brown     PetscCall(DMPlexOrientPositiveEdges_Private(dm));
7813431e603SJed Brown 
7823431e603SJed Brown     DMLabel label;
7833431e603SJed Brown     PetscCall(DMCreateLabel(dm, "Face Sets"));
7843431e603SJed Brown     PetscCall(DMGetLabel(dm, "Face Sets", &label));
7854e2e9504SJed Brown     PetscSegBuffer per_faces, donor_face_closure, my_donor_faces;
7864e2e9504SJed Brown     PetscCall(PetscSegBufferCreate(sizeof(PetscInt), 64, &per_faces));
7874e2e9504SJed Brown     PetscCall(PetscSegBufferCreate(sizeof(PetscInt), 64, &my_donor_faces));
7884e2e9504SJed Brown     PetscCall(PetscSegBufferCreate(sizeof(ZCode), 64 * PetscPowInt(2, dim), &donor_face_closure));
7893431e603SJed Brown     PetscInt fStart, fEnd, vStart, vEnd;
7903431e603SJed Brown     PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
7913431e603SJed Brown     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
7923431e603SJed Brown     for (PetscInt f = fStart; f < fEnd; f++) {
7934e2e9504SJed Brown       PetscInt npoints, *points = NULL, num_fverts = 0, fverts[8];
7943431e603SJed Brown       PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &npoints, &points));
7953431e603SJed Brown       PetscInt bc_count[6] = {0};
7963431e603SJed Brown       for (PetscInt i = 0; i < npoints; i++) {
7973431e603SJed Brown         PetscInt p = points[2 * i];
7983431e603SJed Brown         if (p < vStart || vEnd <= p) continue;
7994e2e9504SJed Brown         fverts[num_fverts++] = p;
8003431e603SJed Brown         Ijk loc              = ZCodeSplit(vert_z[p - vStart]);
8013431e603SJed Brown         // Convention here matches DMPlexCreateCubeMesh_Internal
8023431e603SJed Brown         bc_count[0] += loc.i == 0;
8033431e603SJed Brown         bc_count[1] += loc.i == layout.vextent.i - 1;
8043431e603SJed Brown         bc_count[2] += loc.j == 0;
8053431e603SJed Brown         bc_count[3] += loc.j == layout.vextent.j - 1;
8063431e603SJed Brown         bc_count[4] += loc.k == 0;
8073431e603SJed Brown         bc_count[5] += loc.k == layout.vextent.k - 1;
8083e72e933SJed Brown       }
8093431e603SJed Brown       PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &npoints, &points));
8103431e603SJed Brown       for (PetscInt bc = 0, bc_match = 0; bc < 2 * dim; bc++) {
8113431e603SJed Brown         if (bc_count[bc] == PetscPowInt(2, dim - 1)) {
8126725e60dSJed Brown           PetscCall(DMPlexOrientPoint(dm, f, face_orient_dim[dim][bc]));
8134e2e9504SJed Brown           if (periodicity[bc / 2] == DM_BOUNDARY_PERIODIC) {
8144e2e9504SJed Brown             PetscInt *put;
8154e2e9504SJed Brown             if (bc % 2 == 0) { // donor face; no label
8164e2e9504SJed Brown               PetscCall(PetscSegBufferGet(my_donor_faces, 1, &put));
8174e2e9504SJed Brown               *put = f;
8184e2e9504SJed Brown             } else { // periodic face
8194e2e9504SJed Brown               PetscCall(PetscSegBufferGet(per_faces, 1, &put));
8204e2e9504SJed Brown               *put = f;
8214e2e9504SJed Brown               ZCode *zput;
8224e2e9504SJed Brown               PetscCall(PetscSegBufferGet(donor_face_closure, num_fverts, &zput));
8234e2e9504SJed Brown               for (PetscInt i = 0; i < num_fverts; i++) {
8244e2e9504SJed Brown                 Ijk loc = ZCodeSplit(vert_z[fverts[i] - vStart]);
8254e2e9504SJed Brown                 switch (bc / 2) {
8264e2e9504SJed Brown                 case 0:
8274e2e9504SJed Brown                   loc.i = 0;
8284e2e9504SJed Brown                   break;
8294e2e9504SJed Brown                 case 1:
8304e2e9504SJed Brown                   loc.j = 0;
8314e2e9504SJed Brown                   break;
8324e2e9504SJed Brown                 case 2:
8334e2e9504SJed Brown                   loc.k = 0;
8344e2e9504SJed Brown                   break;
8354e2e9504SJed Brown                 }
8364e2e9504SJed Brown                 *zput++ = ZEncode(loc);
8374e2e9504SJed Brown               }
8384e2e9504SJed Brown             }
8394e2e9504SJed Brown             continue;
8404e2e9504SJed Brown           }
8413431e603SJed Brown           PetscAssert(bc_match == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face matches multiple face sets");
8423431e603SJed Brown           PetscCall(DMLabelSetValue(label, f, face_marker_dim[dim][bc]));
8433431e603SJed Brown           bc_match++;
8443431e603SJed Brown         }
8453431e603SJed Brown       }
8463431e603SJed Brown     }
8473431e603SJed Brown     // Ensure that the Coordinate DM has our new boundary labels
8483431e603SJed Brown     DM cdm;
8493431e603SJed Brown     PetscCall(DMGetCoordinateDM(dm, &cdm));
8503431e603SJed Brown     PetscCall(DMCopyLabels(dm, cdm, PETSC_COPY_VALUES, PETSC_FALSE, DM_COPY_LABELS_FAIL));
8516725e60dSJed Brown     if (periodicity[0] == DM_BOUNDARY_PERIODIC || (dim > 1 && periodicity[1] == DM_BOUNDARY_PERIODIC) || (dim > 2 && periodicity[2] == DM_BOUNDARY_PERIODIC)) {
852*5f06a3ddSJed Brown       PetscCall(DMPlexCreateBoxMesh_Tensor_SFC_Periodicity_Private(dm, &layout, vert_z, per_faces, lower, upper, periodicity, donor_face_closure, my_donor_faces));
8536725e60dSJed Brown     }
8546725e60dSJed Brown     PetscCall(PetscSegBufferDestroy(&per_faces));
8556725e60dSJed Brown     PetscCall(PetscSegBufferDestroy(&donor_face_closure));
8566725e60dSJed Brown     PetscCall(PetscSegBufferDestroy(&my_donor_faces));
8573431e603SJed Brown   }
8584e2e9504SJed Brown   PetscCall(PetscFree(layout.zstarts));
8593431e603SJed Brown   PetscCall(PetscFree(vert_z));
8603e72e933SJed Brown   PetscFunctionReturn(0);
8613e72e933SJed Brown }
8624e2e9504SJed Brown 
8634e2e9504SJed Brown /*@
864*5f06a3ddSJed Brown   DMPlexSetIsoperiodicFaceSF - Express periodicity from an existing mesh
8654e2e9504SJed Brown 
8664e2e9504SJed Brown   Logically collective
8674e2e9504SJed Brown 
8684e2e9504SJed Brown   Input Parameters:
8694e2e9504SJed Brown + dm - The `DMPLEX` on which to set periodicity
8704e2e9504SJed Brown - face_sf - SF in which roots are (owned) donor faces and leaves are faces that must be matched to a (possibly remote) donor face.
8714e2e9504SJed Brown 
8724e2e9504SJed Brown   Level: advanced
8734e2e9504SJed Brown 
8744e2e9504SJed Brown   Notes:
8754e2e9504SJed Brown 
8764e2e9504SJed Brown   One can use `-dm_plex_box_sfc` to use this mode of periodicity, wherein the periodic points are distinct both globally
8774e2e9504SJed Brown   and locally, but are paired when creating a global dof space.
8784e2e9504SJed Brown 
879*5f06a3ddSJed Brown .seealso: [](chapter_unstructured), `DMPLEX`, `DMGetGlobalSection()`, `DMPlexGetIsoperiodicFaceSF()`
8804e2e9504SJed Brown @*/
881*5f06a3ddSJed Brown PetscErrorCode DMPlexSetIsoperiodicFaceSF(DM dm, PetscSF face_sf)
8824e2e9504SJed Brown {
8834e2e9504SJed Brown   DM_Plex *plex = (DM_Plex *)dm->data;
8844e2e9504SJed Brown   PetscFunctionBegin;
8854e2e9504SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8864e2e9504SJed Brown   PetscCall(PetscObjectReference((PetscObject)face_sf));
8874e2e9504SJed Brown   PetscCall(PetscSFDestroy(&plex->periodic.face_sf));
8884e2e9504SJed Brown   plex->periodic.face_sf = face_sf;
889*5f06a3ddSJed Brown   if (face_sf) PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMGetIsoperiodicPointSF_C", DMGetIsoperiodicPointSF_Plex));
890*5f06a3ddSJed Brown 
891*5f06a3ddSJed Brown   DM cdm = dm->coordinates[0].dm; // Can't DMGetCoordinateDM because it automatically creates one
892*5f06a3ddSJed Brown   if (cdm) {
893*5f06a3ddSJed Brown     PetscCall(DMPlexSetIsoperiodicFaceSF(cdm, face_sf));
894*5f06a3ddSJed Brown     if (face_sf) cdm->periodic.setup = DMPeriodicCoordinateSetUp_Internal;
895*5f06a3ddSJed Brown   }
8964e2e9504SJed Brown   PetscFunctionReturn(0);
8974e2e9504SJed Brown }
8984e2e9504SJed Brown 
8994e2e9504SJed Brown /*@
900*5f06a3ddSJed Brown   DMPlexGetIsoperiodicFaceSF - Obtain periodicity for a mesh
9014e2e9504SJed Brown 
9024e2e9504SJed Brown   Logically collective
9034e2e9504SJed Brown 
9044e2e9504SJed Brown   Input Parameters:
9054e2e9504SJed Brown . dm - The `DMPLEX` for which to obtain periodic relation
9064e2e9504SJed Brown 
9074e2e9504SJed Brown   Output Parameters:
9084e2e9504SJed Brown . face_sf - SF in which roots are (owned) donor faces and leaves are faces that must be matched to a (possibly remote) donor face.
9094e2e9504SJed Brown 
9104e2e9504SJed Brown   Level: advanced
9114e2e9504SJed Brown 
912*5f06a3ddSJed Brown .seealso: [](chapter_unstructured), `DMPLEX`, `DMGetGlobalSection()`, `DMPlexSetIsoperiodicFaceSF()`
9134e2e9504SJed Brown @*/
914*5f06a3ddSJed Brown PetscErrorCode DMPlexGetIsoperiodicFaceSF(DM dm, PetscSF *face_sf)
9154e2e9504SJed Brown {
9164e2e9504SJed Brown   DM_Plex *plex = (DM_Plex *)dm->data;
9174e2e9504SJed Brown   PetscFunctionBegin;
9184e2e9504SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9194e2e9504SJed Brown   *face_sf = plex->periodic.face_sf;
9204e2e9504SJed Brown   PetscFunctionReturn(0);
9214e2e9504SJed Brown }
9226725e60dSJed Brown 
923*5f06a3ddSJed Brown /*@C
924*5f06a3ddSJed Brown   DMPlexSetIsoperiodicFaceTransform - set geometric transform from donor to periodic points
9256725e60dSJed Brown 
9266725e60dSJed Brown   Logically Collective
9276725e60dSJed Brown 
9286725e60dSJed Brown   Input Arguments:
929*5f06a3ddSJed Brown + dm - `DMPlex` that has been configured with `DMPlexSetIsoperiodicFaceSF()`
9306725e60dSJed Brown - t - 4x4 affine transformation basis.
9316725e60dSJed Brown 
932*5f06a3ddSJed Brown   Notes:
933*5f06a3ddSJed Brown   Affine transforms are 4x4 matrices in which the leading 3x3 block expresses a rotation (or identity for no rotation),
934*5f06a3ddSJed Brown   the last column contains a translation vector, and the bottom row is all zero except the last entry, which must always
935*5f06a3ddSJed Brown   be 1. This representation is common in geometric modeling and allows affine transformations to be composed using
936*5f06a3ddSJed Brown   simple matrix mulitplication.
937*5f06a3ddSJed Brown 
938*5f06a3ddSJed Brown   Although the interface accepts a general affine transform, only affine translation is supported at present.
939*5f06a3ddSJed Brown 
940*5f06a3ddSJed Brown   Developer Notes:
941*5f06a3ddSJed Brown   This interface should be replaced by making BasisTransform public, expanding it to support affine representations, and
942*5f06a3ddSJed Brown   adding GPU implementations to apply the G2L/L2G transforms.
9436725e60dSJed Brown @*/
944*5f06a3ddSJed Brown PetscErrorCode DMPlexSetIsoperiodicFaceTransform(DM dm, const PetscScalar t[])
9456725e60dSJed Brown {
9466725e60dSJed Brown   DM_Plex *plex = (DM_Plex *)dm->data;
9476725e60dSJed Brown   PetscFunctionBegin;
9486725e60dSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9496725e60dSJed Brown   for (PetscInt i = 0; i < 4; i++) {
9506725e60dSJed Brown     for (PetscInt j = 0; j < 4; j++) {
951*5f06a3ddSJed Brown       PetscCheck(i != j || t[i * 4 + j] == 1., PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Rotated transforms not supported");
952*5f06a3ddSJed Brown       plex->periodic.transform[i][j] = t[i * 4 + j];
9536725e60dSJed Brown     }
9546725e60dSJed Brown   }
9556725e60dSJed Brown   PetscFunctionReturn(0);
9566725e60dSJed Brown }
957