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*92a26154SJames Wright // Modify index array based on the transformation of `point` for the given section and field 351*92a26154SJames Wright // Used for correcting the sfNatural based on point reorientation 352*92a26154SJames Wright static PetscErrorCode DMPlexOrientFieldPointIndex(DM dm, PetscSection section, PetscInt field, PetscInt array_size, PetscInt array[], PetscInt point, PetscInt orientation) 353*92a26154SJames Wright { 354*92a26154SJames Wright PetscInt *copy; 355*92a26154SJames Wright PetscInt dof, off, point_ornt[2] = {point, orientation}; 356*92a26154SJames Wright const PetscInt **perms; 357*92a26154SJames Wright 358*92a26154SJames Wright PetscFunctionBeginUser; 359*92a26154SJames Wright PetscCall(PetscSectionGetDof(section, point, &dof)); 360*92a26154SJames Wright PetscCall(PetscSectionGetOffset(section, point, &off)); 361*92a26154SJames Wright PetscCheck(off + dof <= array_size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Section indices exceed index array bounds"); 362*92a26154SJames Wright PetscCall(DMGetWorkArray(dm, dof, MPIU_INT, ©)); 363*92a26154SJames Wright PetscArraycpy(copy, &array[off], dof); 364*92a26154SJames Wright 365*92a26154SJames Wright PetscCall(PetscSectionGetFieldPointSyms(section, field, 1, point_ornt, &perms, NULL)); 366*92a26154SJames Wright for (PetscInt i = 0; i < dof; i++) { 367*92a26154SJames Wright if (perms[0]) array[off + perms[0][i]] = copy[i]; 368*92a26154SJames Wright } 369*92a26154SJames Wright 370*92a26154SJames Wright PetscCall(PetscSectionRestoreFieldPointSyms(section, field, 1, point_ornt, &perms, NULL)); 371*92a26154SJames Wright PetscCall(DMRestoreWorkArray(dm, dof, MPIU_INT, ©)); 372*92a26154SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 373*92a26154SJames Wright } 374*92a26154SJames Wright 37553b735f8SJames Wright // Modify Vec based on the transformation of `point` for the given section and field 37653b735f8SJames Wright static PetscErrorCode DMPlexOrientFieldPointVec(DM dm, PetscSection section, PetscInt field, Vec V, PetscInt point, PetscInt orientation) 37753b735f8SJames Wright { 37853b735f8SJames Wright PetscScalar *copy, *V_arr; 37953b735f8SJames Wright PetscInt dof, off, point_ornt[2] = {point, orientation}; 38053b735f8SJames Wright const PetscInt **perms; 38153b735f8SJames Wright const PetscScalar **rots; 38253b735f8SJames Wright 38353b735f8SJames Wright PetscFunctionBeginUser; 38453b735f8SJames Wright PetscCall(PetscSectionGetDof(section, point, &dof)); 38553b735f8SJames Wright PetscCall(PetscSectionGetOffset(section, point, &off)); 38653b735f8SJames Wright PetscCall(VecGetArray(V, &V_arr)); 38753b735f8SJames Wright PetscCall(DMGetWorkArray(dm, dof, MPIU_SCALAR, ©)); 38853b735f8SJames Wright PetscArraycpy(copy, &V_arr[off], dof); 38953b735f8SJames Wright 39053b735f8SJames Wright PetscCall(PetscSectionGetFieldPointSyms(section, field, 1, point_ornt, &perms, &rots)); 39153b735f8SJames Wright for (PetscInt i = 0; i < dof; i++) { 39253b735f8SJames Wright if (perms[0]) V_arr[off + perms[0][i]] = copy[i]; 39353b735f8SJames Wright if (rots[0]) V_arr[off + perms[0][i]] *= rots[0][i]; 39453b735f8SJames Wright } 39553b735f8SJames Wright 39653b735f8SJames Wright PetscCall(PetscSectionRestoreFieldPointSyms(section, field, 1, point_ornt, &perms, &rots)); 39753b735f8SJames Wright PetscCall(DMRestoreWorkArray(dm, dof, MPIU_SCALAR, ©)); 39853b735f8SJames Wright PetscCall(VecRestoreArray(V, &V_arr)); 39953b735f8SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 40053b735f8SJames Wright } 40153b735f8SJames Wright 40253b735f8SJames Wright // Reorient the point in the DMPlex while also applying necessary corrections to other structures (e.g. coordinates) 403*92a26154SJames Wright static PetscErrorCode DMPlexOrientPointWithCorrections(DM dm, PetscInt point, PetscInt ornt, PetscSection perm_section, PetscInt perm_length, PetscInt perm[]) 40453b735f8SJames Wright { 40553b735f8SJames 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) 40653b735f8SJames Wright PetscFunctionBeginUser; 40753b735f8SJames Wright PetscCall(DMPlexOrientPoint(dm, point, ornt)); 40853b735f8SJames Wright 40953b735f8SJames Wright { // Correct coordinates based on new cone ordering 41053b735f8SJames Wright DM cdm; 41153b735f8SJames Wright PetscSection csection; 41253b735f8SJames Wright Vec coordinates; 41353b735f8SJames Wright PetscInt pStart, pEnd; 41453b735f8SJames Wright 41553b735f8SJames Wright PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 41653b735f8SJames Wright PetscCall(DMGetCoordinateDM(dm, &cdm)); 41753b735f8SJames Wright PetscCall(DMGetLocalSection(cdm, &csection)); 41853b735f8SJames Wright PetscCall(PetscSectionGetChart(csection, &pStart, &pEnd)); 41953b735f8SJames Wright if (IsPointInsideStratum(point, pStart, pEnd)) PetscCall(DMPlexOrientFieldPointVec(cdm, csection, 0, coordinates, point, ornt)); 42053b735f8SJames Wright } 421*92a26154SJames Wright 422*92a26154SJames Wright if (perm_section) { 423*92a26154SJames Wright PetscInt num_fields; 424*92a26154SJames Wright PetscCall(PetscSectionGetNumFields(perm_section, &num_fields)); 425*92a26154SJames Wright for (PetscInt f = 0; f < num_fields; f++) PetscCall(DMPlexOrientFieldPointIndex(dm, perm_section, f, perm_length, perm, point, ornt)); 426*92a26154SJames Wright } 42753b735f8SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 42853b735f8SJames Wright } 42953b735f8SJames Wright 43093defd67SJames 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[]` 43153b735f8SJames 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) 43293defd67SJames Wright { 43393defd67SJames Wright MPI_Comm comm; 43493defd67SJames Wright PetscMPIInt rank; 43593defd67SJames Wright PetscInt nroots, nleaves; 43693defd67SJames Wright PetscInt *rootdata, *leafdata; 43793defd67SJames Wright const PetscInt *filocal; 43893defd67SJames Wright const PetscSFNode *firemote; 43993defd67SJames Wright 44093defd67SJames Wright PetscFunctionBeginUser; 44193defd67SJames Wright PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 44293defd67SJames Wright PetscCallMPI(MPI_Comm_rank(comm, &rank)); 44393defd67SJames Wright 44493defd67SJames Wright PetscCall(PetscSFGetGraph(face_sf, &nroots, &nleaves, &filocal, &firemote)); 44593defd67SJames Wright PetscCall(PetscCalloc2(2 * nroots, &rootdata, 2 * nroots, &leafdata)); 44693defd67SJames Wright for (PetscInt i = 0; i < nleaves; i++) { 44793defd67SJames Wright PetscInt point = filocal[i]; 44885de0153SJames 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); 44993defd67SJames Wright leafdata[point] = point_sizes[point - pStart]; 45093defd67SJames Wright } 45193defd67SJames Wright PetscCall(PetscSFReduceBegin(face_sf, MPIU_INT, leafdata, rootdata + nroots, MPIU_SUM)); 45293defd67SJames Wright PetscCall(PetscSFReduceEnd(face_sf, MPIU_INT, leafdata, rootdata + nroots, MPIU_SUM)); 45393defd67SJames Wright 45493defd67SJames Wright PetscInt root_offset = 0; 45593defd67SJames Wright PetscCall(PetscBTCreate(nroots, rootbt)); 45693defd67SJames Wright for (PetscInt p = 0; p < nroots; p++) { 45793defd67SJames Wright const PetscInt *donor_dof = rootdata + nroots; 45893defd67SJames Wright if (donor_dof[p] == 0) { 45993defd67SJames Wright rootdata[2 * p] = -1; 46093defd67SJames Wright rootdata[2 * p + 1] = -1; 46193defd67SJames Wright continue; 46293defd67SJames Wright } 46393defd67SJames Wright PetscCall(PetscBTSet(*rootbt, p)); 46485de0153SJames 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); 46593defd67SJames Wright PetscInt p_size = point_sizes[p - pStart]; 46693defd67SJames 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); 46793defd67SJames Wright rootdata[2 * p] = root_offset; 46893defd67SJames Wright rootdata[2 * p + 1] = p_size; 46993defd67SJames Wright root_offset += p_size; 47093defd67SJames Wright } 47193defd67SJames Wright PetscCall(PetscSFBcastBegin(face_sf, MPIU_2INT, rootdata, leafdata, MPI_REPLACE)); 47293defd67SJames Wright PetscCall(PetscSFBcastEnd(face_sf, MPIU_2INT, rootdata, leafdata, MPI_REPLACE)); 47393defd67SJames Wright // Count how many leaves we need to communicate the closures 47493defd67SJames Wright PetscInt leaf_offset = 0; 47593defd67SJames Wright for (PetscInt i = 0; i < nleaves; i++) { 47693defd67SJames Wright PetscInt point = filocal[i]; 47793defd67SJames Wright if (leafdata[2 * point + 1] < 0) continue; 47893defd67SJames Wright leaf_offset += leafdata[2 * point + 1]; 47993defd67SJames Wright } 48093defd67SJames Wright 48193defd67SJames Wright PetscSFNode *closure_leaf; 48293defd67SJames Wright PetscCall(PetscMalloc1(leaf_offset, &closure_leaf)); 48393defd67SJames Wright leaf_offset = 0; 48493defd67SJames Wright for (PetscInt i = 0; i < nleaves; i++) { 48593defd67SJames Wright PetscInt point = filocal[i]; 48693defd67SJames Wright PetscInt cl_size = leafdata[2 * point + 1]; 48793defd67SJames Wright if (cl_size < 0) continue; 48893defd67SJames Wright for (PetscInt j = 0; j < cl_size; j++) { 48993defd67SJames Wright closure_leaf[leaf_offset].rank = firemote[i].rank; 49093defd67SJames Wright closure_leaf[leaf_offset].index = leafdata[2 * point] + j; 49193defd67SJames Wright leaf_offset++; 49293defd67SJames Wright } 49393defd67SJames Wright } 49493defd67SJames Wright 49593defd67SJames Wright PetscCall(PetscSFCreate(comm, sf_closure)); 49693defd67SJames Wright PetscCall(PetscSFSetGraph(*sf_closure, root_offset, leaf_offset, NULL, PETSC_USE_POINTER, closure_leaf, PETSC_OWN_POINTER)); 49793defd67SJames Wright *rootbuffersize = root_offset; 49893defd67SJames Wright *leafbuffersize = leaf_offset; 49993defd67SJames Wright PetscCall(PetscFree2(rootdata, leafdata)); 50093defd67SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 50193defd67SJames Wright } 50293defd67SJames Wright 50353b735f8SJames Wright // Determine if `key` is in `array`. `array` does NOT need to be sorted 50453b735f8SJames Wright static inline PetscBool SearchIntArray(PetscInt key, PetscInt array_size, const PetscInt array[]) 50553b735f8SJames Wright { 50653b735f8SJames Wright for (PetscInt i = 0; i < array_size; i++) 50753b735f8SJames Wright if (array[i] == key) return PETSC_TRUE; 50853b735f8SJames Wright return PETSC_FALSE; 50953b735f8SJames Wright } 51053b735f8SJames Wright 51153b735f8SJames Wright // Translate a cone in periodic points to the cone in donor points based on the `periodic2donor` array 51253b735f8SJames Wright static inline PetscErrorCode TranslateConeP2D(const PetscInt periodic_cone[], PetscInt cone_size, const PetscInt periodic2donor[], PetscInt p2d_count, PetscInt p2d_cone[]) 51353b735f8SJames Wright { 51453b735f8SJames Wright PetscFunctionBeginUser; 51553b735f8SJames Wright for (PetscInt p = 0; p < cone_size; p++) { 51653b735f8SJames Wright PetscInt p2d_index = -1; 51753b735f8SJames Wright for (PetscInt p2d = 0; p2d < p2d_count; p2d++) { 51853b735f8SJames Wright if (periodic2donor[p2d * 2] == periodic_cone[p]) p2d_index = p2d; 51953b735f8SJames Wright } 52053b735f8SJames Wright PetscCheck(p2d_index >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not find periodic point in periodic-to-donor array"); 52153b735f8SJames Wright p2d_cone[p] = periodic2donor[2 * p2d_index + 1]; 52253b735f8SJames Wright } 52353b735f8SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 52453b735f8SJames Wright } 52553b735f8SJames Wright 52653b735f8SJames Wright // Corrects the cone order of periodic faces (and their transitive closure's cones) to match their donor face pair. 52753b735f8SJames Wright // 52853b735f8SJames Wright // This is done by: 52953b735f8SJames 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 53053b735f8SJames Wright // - The donor vertices have the isoperiodic transform applied to them such that they should match exactly 53153b735f8SJames Wright // 2. Translating the periodic vertices into the donor vertices point IDs 53253b735f8SJames Wright // 3. Translating the cone of each periodic point into the donor point IDs 53353b735f8SJames Wright // 4. Comparing the periodic-to-donor cone to the donor cone for each point 53453b735f8SJames Wright // 5. Apply the necessary transformation to the periodic cone to make it match the donor cone 53553b735f8SJames Wright static PetscErrorCode DMPlexCorrectOrientationForIsoperiodic(DM dm) 53653b735f8SJames Wright { 53753b735f8SJames Wright MPI_Comm comm; 53853b735f8SJames Wright DM_Plex *plex = (DM_Plex *)dm->data; 53953b735f8SJames Wright PetscInt nroots, nleaves; 540*92a26154SJames Wright PetscInt *local_vec_perm = NULL, local_vec_length = 0, *global_vec_perm = NULL, global_vec_length = 0; 541*92a26154SJames Wright const PetscInt *filocal, coords_field_id = 0; 54253b735f8SJames Wright DM cdm; 543*92a26154SJames Wright PetscSection csection, localSection = NULL; 544*92a26154SJames Wright PetscSF sfNatural_old = NULL; 54553b735f8SJames Wright Vec coordinates; 546*92a26154SJames Wright PetscMPIInt myrank; 54753b735f8SJames Wright PetscBool debug_printing = PETSC_FALSE; 54853b735f8SJames Wright 54953b735f8SJames Wright PetscFunctionBeginUser; 55053b735f8SJames Wright PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 551*92a26154SJames Wright PetscCallMPI(MPI_Comm_rank(comm, &myrank)); 55253b735f8SJames Wright PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 55353b735f8SJames Wright PetscCheck(coordinates, comm, PETSC_ERR_ARG_WRONGSTATE, "DM must have coordinates to setup isoperiodic"); 55453b735f8SJames Wright PetscCall(DMGetCoordinateDM(dm, &cdm)); 55553b735f8SJames Wright PetscCall(DMGetLocalSection(cdm, &csection)); 55653b735f8SJames Wright 557*92a26154SJames Wright PetscCall(DMGetNaturalSF(dm, &sfNatural_old)); 558*92a26154SJames Wright // Prep data for naturalSF correction 559*92a26154SJames Wright if (plex->periodic.num_face_sfs > 0 && sfNatural_old) { 560*92a26154SJames Wright PetscSection globalSection; 561*92a26154SJames Wright PetscSF pointSF, sectionSF; 562*92a26154SJames Wright PetscInt nleaves; 563*92a26154SJames Wright const PetscInt *ilocal; 564*92a26154SJames Wright const PetscSFNode *iremote; 565*92a26154SJames Wright 566*92a26154SJames Wright // Create global section with just pointSF and including constraints 567*92a26154SJames Wright PetscCall(DMGetLocalSection(dm, &localSection)); 568*92a26154SJames Wright PetscCall(DMGetPointSF(dm, &pointSF)); 569*92a26154SJames Wright PetscCall(PetscSectionCreateGlobalSection(localSection, pointSF, PETSC_TRUE, PETSC_TRUE, PETSC_FALSE, &globalSection)); 570*92a26154SJames Wright 571*92a26154SJames Wright // Set local_vec_perm to be negative values when that dof is not owned by the current rank 572*92a26154SJames Wright // Dofs that are owned are set to their corresponding global Vec index 573*92a26154SJames Wright PetscCall(PetscSectionGetStorageSize(globalSection, &global_vec_length)); 574*92a26154SJames Wright PetscCall(PetscSectionGetStorageSize(localSection, &local_vec_length)); 575*92a26154SJames Wright PetscCall(PetscMalloc2(global_vec_length, &global_vec_perm, local_vec_length, &local_vec_perm)); 576*92a26154SJames Wright for (PetscInt i = 0; i < global_vec_length; i++) global_vec_perm[i] = i; 577*92a26154SJames Wright for (PetscInt i = 0; i < local_vec_length; i++) local_vec_perm[i] = -(i + 1); 578*92a26154SJames Wright 579*92a26154SJames Wright PetscCall(PetscSFCreate(comm, §ionSF)); 580*92a26154SJames Wright PetscCall(PetscSFSetGraphSection(sectionSF, localSection, globalSection)); 581*92a26154SJames Wright PetscCall(PetscSFGetGraph(sectionSF, NULL, &nleaves, &ilocal, &iremote)); 582*92a26154SJames Wright for (PetscInt l = 0; l < nleaves; l++) { 583*92a26154SJames Wright if (iremote[l].rank != myrank) continue; 584*92a26154SJames Wright PetscInt local_index = ilocal ? ilocal[l] : l; 585*92a26154SJames Wright local_vec_perm[local_index] = global_vec_perm[iremote[l].index]; 586*92a26154SJames Wright } 587*92a26154SJames Wright PetscCall(PetscSectionDestroy(&globalSection)); 588*92a26154SJames Wright PetscCall(PetscSFDestroy(§ionSF)); 589*92a26154SJames Wright } 590*92a26154SJames Wright 591*92a26154SJames Wright // Create sf_vert_coords and sf_face_cones for communicating donor vertices and cones to periodic faces, respectively 59253b735f8SJames Wright for (PetscInt f = 0; f < plex->periodic.num_face_sfs; f++) { 59353b735f8SJames Wright PetscSF face_sf = plex->periodic.face_sfs[f]; 59453b735f8SJames Wright const PetscScalar(*transform)[4] = (const PetscScalar(*)[4])plex->periodic.transform[f]; 59553b735f8SJames Wright PetscInt *face_vertices_size, *face_cones_size; 59653b735f8SJames Wright PetscInt fStart, fEnd, vStart, vEnd, rootnumvert, leafnumvert, rootconesize, leafconesize, dim; 59753b735f8SJames Wright PetscSF sf_vert_coords, sf_face_cones; 59853b735f8SJames Wright PetscBT rootbt; 59953b735f8SJames Wright 60053b735f8SJames Wright PetscCall(DMGetCoordinateDim(dm, &dim)); 60153b735f8SJames Wright PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 60253b735f8SJames Wright PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd)); 60353b735f8SJames Wright PetscCall(PetscCalloc2(fEnd - fStart, &face_vertices_size, fEnd - fStart, &face_cones_size)); 60453b735f8SJames Wright 60553b735f8SJames Wright // Create SFs to communicate donor vertices and donor cones to periodic faces 60653b735f8SJames Wright for (PetscInt f = fStart, index = 0; f < fEnd; f++, index++) { 60753b735f8SJames Wright PetscInt cl_size, *closure = NULL, num_vertices = 0; 60853b735f8SJames Wright PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &cl_size, &closure)); 60953b735f8SJames Wright for (PetscInt p = 0; p < cl_size; p++) { 61053b735f8SJames Wright PetscInt cl_point = closure[2 * p]; 61153b735f8SJames Wright if (IsPointInsideStratum(cl_point, vStart, vEnd)) num_vertices++; 61253b735f8SJames Wright else { 61353b735f8SJames Wright PetscInt cone_size; 61453b735f8SJames Wright PetscCall(DMPlexGetConeSize(dm, cl_point, &cone_size)); 61553b735f8SJames Wright face_cones_size[index] += cone_size + 2; 61653b735f8SJames Wright } 61753b735f8SJames Wright } 61853b735f8SJames Wright face_vertices_size[index] = num_vertices; 61953b735f8SJames Wright face_cones_size[index] += num_vertices; 62053b735f8SJames Wright PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &cl_size, &closure)); 62153b735f8SJames Wright } 62253b735f8SJames Wright PetscCall(CreateDonorToPeriodicSF(dm, face_sf, fStart, fEnd, face_vertices_size, &rootnumvert, &leafnumvert, &rootbt, &sf_vert_coords)); 62353b735f8SJames Wright PetscCall(PetscBTDestroy(&rootbt)); 62453b735f8SJames Wright PetscCall(CreateDonorToPeriodicSF(dm, face_sf, fStart, fEnd, face_cones_size, &rootconesize, &leafconesize, &rootbt, &sf_face_cones)); 62553b735f8SJames Wright 62653b735f8SJames Wright PetscCall(PetscSFGetGraph(face_sf, &nroots, &nleaves, &filocal, NULL)); 62753b735f8SJames Wright 62853b735f8SJames Wright PetscReal *leaf_donor_coords; 62953b735f8SJames Wright PetscInt *leaf_donor_cones; 63053b735f8SJames Wright 63153b735f8SJames Wright { // Communicate donor coords and cones to the periodic faces 63253b735f8SJames Wright PetscReal *mydonor_vertices; 63353b735f8SJames Wright PetscInt *mydonor_cones; 63453b735f8SJames Wright const PetscScalar *coords_arr; 63553b735f8SJames Wright 63653b735f8SJames Wright PetscCall(PetscCalloc2(rootnumvert * dim, &mydonor_vertices, rootconesize, &mydonor_cones)); 63753b735f8SJames Wright PetscCall(VecGetArrayRead(coordinates, &coords_arr)); 63853b735f8SJames Wright for (PetscInt donor_face = 0, donor_vert_offset = 0, donor_cone_offset = 0; donor_face < nroots; donor_face++) { 63953b735f8SJames Wright if (!PetscBTLookup(rootbt, donor_face)) continue; 64053b735f8SJames Wright PetscInt cl_size, *closure = NULL; 64153b735f8SJames Wright 64253b735f8SJames Wright PetscCall(DMPlexGetTransitiveClosure(dm, donor_face, PETSC_TRUE, &cl_size, &closure)); 64353b735f8SJames Wright // Pack vertex coordinates 64453b735f8SJames Wright for (PetscInt p = 0; p < cl_size; p++) { 64553b735f8SJames Wright PetscInt cl_point = closure[2 * p], dof, offset; 64653b735f8SJames Wright if (!IsPointInsideStratum(cl_point, vStart, vEnd)) continue; 64753b735f8SJames Wright mydonor_cones[donor_cone_offset++] = cl_point; 64853b735f8SJames Wright PetscCall(PetscSectionGetFieldDof(csection, cl_point, coords_field_id, &dof)); 64953b735f8SJames Wright PetscCall(PetscSectionGetFieldOffset(csection, cl_point, coords_field_id, &offset)); 65053b735f8SJames 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); 65153b735f8SJames Wright // Apply isoperiodic transform to donor vertices such that corresponding periodic vertices should match exactly 65253b735f8SJames Wright for (PetscInt d = 0; d < dof; d++) mydonor_vertices[donor_vert_offset * dim + d] = PetscRealPart(coords_arr[offset + d]) + PetscRealPart(transform[d][3]); 65353b735f8SJames Wright donor_vert_offset++; 65453b735f8SJames Wright } 65553b735f8SJames Wright // Pack cones of face points (including face itself) 65653b735f8SJames Wright for (PetscInt p = 0; p < cl_size; p++) { 65753b735f8SJames Wright PetscInt cl_point = closure[2 * p], cone_size, depth; 65853b735f8SJames Wright const PetscInt *cone; 65953b735f8SJames Wright 66053b735f8SJames Wright PetscCall(DMPlexGetConeSize(dm, cl_point, &cone_size)); 66153b735f8SJames Wright PetscCall(DMPlexGetCone(dm, cl_point, &cone)); 66253b735f8SJames Wright PetscCall(DMPlexGetPointDepth(dm, cl_point, &depth)); 66353b735f8SJames Wright if (depth == 0) continue; // don't include vertex depth 66453b735f8SJames Wright mydonor_cones[donor_cone_offset++] = cone_size; 66553b735f8SJames Wright mydonor_cones[donor_cone_offset++] = cl_point; 66653b735f8SJames Wright PetscArraycpy(&mydonor_cones[donor_cone_offset], cone, cone_size); 66753b735f8SJames Wright donor_cone_offset += cone_size; 66853b735f8SJames Wright } 66953b735f8SJames Wright PetscCall(DMPlexRestoreTransitiveClosure(dm, donor_face, PETSC_TRUE, &cl_size, &closure)); 67053b735f8SJames Wright } 67153b735f8SJames Wright PetscCall(VecRestoreArrayRead(coordinates, &coords_arr)); 67253b735f8SJames Wright PetscCall(PetscBTDestroy(&rootbt)); 67353b735f8SJames Wright 67453b735f8SJames Wright MPI_Datatype vertex_unit; 675f370e7f7SJose E. Roman PetscMPIInt n; 676f370e7f7SJose E. Roman PetscCall(PetscMPIIntCast(dim, &n)); 677f370e7f7SJose E. Roman PetscCallMPI(MPI_Type_contiguous(n, MPIU_REAL, &vertex_unit)); 67853b735f8SJames Wright PetscCallMPI(MPI_Type_commit(&vertex_unit)); 67953b735f8SJames Wright PetscCall(PetscMalloc2(leafnumvert * 3, &leaf_donor_coords, leafconesize, &leaf_donor_cones)); 68053b735f8SJames Wright PetscCall(PetscSFBcastBegin(sf_vert_coords, vertex_unit, mydonor_vertices, leaf_donor_coords, MPI_REPLACE)); 68153b735f8SJames Wright PetscCall(PetscSFBcastBegin(sf_face_cones, MPIU_INT, mydonor_cones, leaf_donor_cones, MPI_REPLACE)); 68253b735f8SJames Wright PetscCall(PetscSFBcastEnd(sf_vert_coords, vertex_unit, mydonor_vertices, leaf_donor_coords, MPI_REPLACE)); 68353b735f8SJames Wright PetscCall(PetscSFBcastEnd(sf_face_cones, MPIU_INT, mydonor_cones, leaf_donor_cones, MPI_REPLACE)); 68453b735f8SJames Wright PetscCall(PetscSFDestroy(&sf_vert_coords)); 68553b735f8SJames Wright PetscCall(PetscSFDestroy(&sf_face_cones)); 68653b735f8SJames Wright PetscCallMPI(MPI_Type_free(&vertex_unit)); 68753b735f8SJames Wright PetscCall(PetscFree2(mydonor_vertices, mydonor_cones)); 68853b735f8SJames Wright } 68953b735f8SJames Wright 69053b735f8SJames Wright { // Determine periodic orientation w/rt donor vertices and reorient 69153b735f8SJames Wright PetscReal tol = PetscSqr(PETSC_MACHINE_EPSILON * 1e3); 69253b735f8SJames Wright PetscInt *periodic2donor, dm_depth, maxConeSize; 69353b735f8SJames Wright PetscInt coords_offset = 0, cones_offset = 0; 69453b735f8SJames Wright 69553b735f8SJames Wright PetscCall(DMPlexGetDepth(dm, &dm_depth)); 69653b735f8SJames Wright PetscCall(DMPlexGetMaxSizes(dm, &maxConeSize, NULL)); 69753b735f8SJames Wright PetscCall(DMGetWorkArray(dm, 2 * PetscPowInt(maxConeSize, dm_depth - 1), MPIU_INT, &periodic2donor)); 69853b735f8SJames Wright 69953b735f8SJames Wright // Translate the periodic face vertices into the donor vertices 70053b735f8SJames Wright // Translation stored in periodic2donor 70153b735f8SJames Wright for (PetscInt i = 0; i < nleaves; i++) { 70253b735f8SJames Wright PetscInt periodic_face = filocal[i], cl_size, num_verts = face_vertices_size[periodic_face - fStart]; 70353b735f8SJames Wright PetscInt cones_size = face_cones_size[periodic_face - fStart], p2d_count = 0; 70453b735f8SJames Wright PetscInt *closure = NULL; 70553b735f8SJames Wright 70653b735f8SJames Wright PetscCall(DMPlexGetTransitiveClosure(dm, periodic_face, PETSC_TRUE, &cl_size, &closure)); 70753b735f8SJames Wright for (PetscInt p = 0; p < cl_size; p++) { 70853b735f8SJames Wright PetscInt cl_point = closure[2 * p], coords_size, donor_vertex = -1; 70953b735f8SJames Wright PetscScalar *coords = NULL; 71053b735f8SJames Wright 71153b735f8SJames Wright if (!IsPointInsideStratum(cl_point, vStart, vEnd)) continue; 71253b735f8SJames Wright PetscCall(DMPlexVecGetClosure(dm, csection, coordinates, cl_point, &coords_size, &coords)); 71353b735f8SJames 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); 71453b735f8SJames Wright 71553b735f8SJames Wright for (PetscInt v = 0; v < num_verts; v++) { 71653b735f8SJames Wright PetscReal dist_sqr = 0; 71753b735f8SJames Wright for (PetscInt d = 0; d < coords_size; d++) dist_sqr += PetscSqr(PetscRealPart(coords[d]) - leaf_donor_coords[(v + coords_offset) * dim + d]); 71853b735f8SJames Wright if (dist_sqr < tol) { 71953b735f8SJames Wright donor_vertex = leaf_donor_cones[cones_offset + v]; 72053b735f8SJames Wright break; 72153b735f8SJames Wright } 72253b735f8SJames Wright } 72353b735f8SJames 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); 72453b735f8SJames Wright if (PetscDefined(USE_DEBUG)) { 72553b735f8SJames 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"); 72653b735f8SJames Wright } 72753b735f8SJames Wright 72853b735f8SJames Wright periodic2donor[2 * p2d_count + 0] = cl_point; 72953b735f8SJames Wright periodic2donor[2 * p2d_count + 1] = donor_vertex; 73053b735f8SJames Wright p2d_count++; 73153b735f8SJames Wright PetscCall(DMPlexVecRestoreClosure(dm, csection, coordinates, cl_point, &coords_size, &coords)); 73253b735f8SJames Wright } 73353b735f8SJames Wright coords_offset += num_verts; 73453b735f8SJames Wright PetscCall(DMPlexRestoreTransitiveClosure(dm, periodic_face, PETSC_TRUE, &cl_size, &closure)); 73553b735f8SJames Wright 73653b735f8SJames Wright { // Determine periodic orientation w/rt donor vertices and reorient 73753b735f8SJames Wright PetscInt depth, *p2d_cone, face_is_array[1] = {periodic_face}; 73853b735f8SJames Wright IS *is_arrays, face_is; 73953b735f8SJames Wright PetscSection *section_arrays; 74053b735f8SJames Wright PetscInt *donor_cone_array = &leaf_donor_cones[cones_offset + num_verts]; 74153b735f8SJames Wright 74253b735f8SJames Wright PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 1, face_is_array, PETSC_USE_POINTER, &face_is)); 74353b735f8SJames Wright PetscCall(DMPlexGetConeRecursive(dm, face_is, &depth, &is_arrays, §ion_arrays)); 74453b735f8SJames Wright PetscCall(ISDestroy(&face_is)); 74553b735f8SJames Wright PetscCall(DMGetWorkArray(dm, maxConeSize, MPIU_INT, &p2d_cone)); 74653b735f8SJames Wright for (PetscInt d = 0; d < depth - 1; d++) { 74753b735f8SJames Wright PetscInt pStart, pEnd; 74853b735f8SJames Wright PetscSection section = section_arrays[d]; 74953b735f8SJames Wright const PetscInt *periodic_cone_arrays, *periodic_point_arrays; 75053b735f8SJames Wright 75153b735f8SJames Wright PetscCall(ISGetIndices(is_arrays[d], &periodic_cone_arrays)); 75253b735f8SJames Wright PetscCall(ISGetIndices(is_arrays[d + 1], &periodic_point_arrays)); // Points at d+1 correspond to the cones at d 75353b735f8SJames Wright PetscCall(PetscSectionGetChart(section_arrays[d], &pStart, &pEnd)); 75453b735f8SJames Wright for (PetscInt p = pStart; p < pEnd; p++) { 75553b735f8SJames Wright PetscInt periodic_cone_size, periodic_cone_offset, periodic_point = periodic_point_arrays[p]; 75653b735f8SJames Wright 75753b735f8SJames Wright PetscCall(PetscSectionGetDof(section, p, &periodic_cone_size)); 75853b735f8SJames Wright PetscCall(PetscSectionGetOffset(section, p, &periodic_cone_offset)); 75953b735f8SJames Wright const PetscInt *periodic_cone = &periodic_cone_arrays[periodic_cone_offset]; 76053b735f8SJames Wright PetscCall(TranslateConeP2D(periodic_cone, periodic_cone_size, periodic2donor, p2d_count, p2d_cone)); 76153b735f8SJames Wright 76253b735f8SJames Wright // Find the donor cone that matches the periodic point's cone 76353b735f8SJames Wright PetscInt donor_cone_offset = 0, donor_point = -1, *donor_cone = NULL; 76453b735f8SJames Wright PetscBool cone_matches = PETSC_FALSE; 76553b735f8SJames Wright while (donor_cone_offset < cones_size - num_verts) { 76653b735f8SJames Wright PetscInt donor_cone_size = donor_cone_array[donor_cone_offset]; 76753b735f8SJames Wright donor_point = donor_cone_array[donor_cone_offset + 1]; 76853b735f8SJames Wright donor_cone = &donor_cone_array[donor_cone_offset + 2]; 76953b735f8SJames Wright 77053b735f8SJames Wright if (donor_cone_size != periodic_cone_size) goto next_cone; 77153b735f8SJames Wright for (PetscInt c = 0; c < periodic_cone_size; c++) { 77253b735f8SJames Wright cone_matches = SearchIntArray(donor_cone[c], periodic_cone_size, p2d_cone); 77353b735f8SJames Wright if (!cone_matches) goto next_cone; 77453b735f8SJames Wright } 77553b735f8SJames Wright // Save the found donor cone's point to the translation array. These will be used for higher depth points. 77653b735f8SJames Wright // i.e. we save the edge translations for when we look for face cones 77753b735f8SJames Wright periodic2donor[2 * p2d_count + 0] = periodic_point; 77853b735f8SJames Wright periodic2donor[2 * p2d_count + 1] = donor_point; 77953b735f8SJames Wright p2d_count++; 78053b735f8SJames Wright break; 78153b735f8SJames Wright 78253b735f8SJames Wright next_cone: 78353b735f8SJames Wright donor_cone_offset += donor_cone_size + 2; 78453b735f8SJames Wright } 78553b735f8SJames 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); 78653b735f8SJames Wright 78753b735f8SJames Wright { // Compare the donor cone with the translated periodic cone and reorient 78853b735f8SJames Wright PetscInt ornt; 78953b735f8SJames Wright DMPolytopeType cell_type; 79053b735f8SJames Wright PetscBool found; 79153b735f8SJames Wright PetscCall(DMPlexGetCellType(dm, periodic_point, &cell_type)); 79253b735f8SJames Wright PetscCall(DMPolytopeMatchOrientation(cell_type, donor_cone, p2d_cone, &ornt, &found)); 79353b735f8SJames 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); 79453b735f8SJames Wright if (debug_printing) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Reorienting point %" PetscInt_FMT " by %" PetscInt_FMT "\n", periodic_point, ornt)); 795*92a26154SJames Wright PetscCall(DMPlexOrientPointWithCorrections(dm, periodic_point, ornt, localSection, local_vec_length, local_vec_perm)); 79653b735f8SJames Wright } 79753b735f8SJames Wright } 79853b735f8SJames Wright PetscCall(ISRestoreIndices(is_arrays[d], &periodic_cone_arrays)); 79953b735f8SJames Wright PetscCall(ISRestoreIndices(is_arrays[d + 1], &periodic_point_arrays)); 80053b735f8SJames Wright } 80153b735f8SJames Wright PetscCall(DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &p2d_cone)); 80253b735f8SJames Wright PetscCall(DMPlexRestoreConeRecursive(dm, face_is, &depth, &is_arrays, §ion_arrays)); 80353b735f8SJames Wright } 80453b735f8SJames Wright 80553b735f8SJames Wright PetscCall(DMPlexRestoreTransitiveClosure(dm, periodic_face, PETSC_TRUE, &cl_size, &closure)); 80653b735f8SJames Wright cones_offset += cones_size; 80753b735f8SJames Wright } 80853b735f8SJames Wright PetscCall(DMRestoreWorkArray(dm, 2 * PetscPowInt(maxConeSize, dm_depth - 1), MPIU_INT, &periodic2donor)); 80953b735f8SJames Wright } 810*92a26154SJames Wright // Re-set local coordinates (i.e. destroy global coordinates if they were modified) 8110b6e880eSJames Wright PetscCall(DMSetCoordinatesLocal(dm, coordinates)); 81253b735f8SJames Wright 81353b735f8SJames Wright PetscCall(PetscFree2(leaf_donor_coords, leaf_donor_cones)); 81453b735f8SJames Wright PetscCall(PetscFree2(face_vertices_size, face_cones_size)); 81553b735f8SJames Wright } 816*92a26154SJames Wright 817*92a26154SJames Wright if (sfNatural_old) { // Correct SFNatural based on the permutation of the local vector 818*92a26154SJames Wright PetscSF newglob_to_oldglob_sf, sfNatural_old, sfNatural_new; 819*92a26154SJames Wright PetscSFNode *remote; 820*92a26154SJames Wright 821*92a26154SJames Wright { // Translate permutation of local Vec into permutation of global Vec 822*92a26154SJames Wright PetscInt g = 0; 823*92a26154SJames Wright PetscBT global_vec_check; // Verify that global indices aren't doubled 824*92a26154SJames Wright PetscCall(PetscBTCreate(global_vec_length, &global_vec_check)); 825*92a26154SJames Wright for (PetscInt l = 0; l < local_vec_length; l++) { 826*92a26154SJames Wright PetscInt global_index = local_vec_perm[l]; 827*92a26154SJames Wright if (global_index < 0) continue; 828*92a26154SJames Wright PetscCheck(!PetscBTLookupSet(global_vec_check, global_index), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Found duplicate global index %" PetscInt_FMT " in local_vec_perm", global_index); 829*92a26154SJames Wright global_vec_perm[g++] = global_index; 830*92a26154SJames Wright } 831*92a26154SJames Wright PetscCheck(g == global_vec_length, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Wrong number of non-negative local indices"); 832*92a26154SJames Wright PetscCall(PetscBTDestroy(&global_vec_check)); 833*92a26154SJames Wright } 834*92a26154SJames Wright 835*92a26154SJames Wright PetscCall(PetscMalloc1(global_vec_length, &remote)); 836*92a26154SJames Wright for (PetscInt i = 0; i < global_vec_length; i++) { 837*92a26154SJames Wright remote[i].rank = myrank; 838*92a26154SJames Wright remote[i].index = global_vec_perm[i]; 839*92a26154SJames Wright } 840*92a26154SJames Wright PetscCall(PetscFree2(global_vec_perm, local_vec_perm)); 841*92a26154SJames Wright PetscCall(PetscSFCreate(comm, &newglob_to_oldglob_sf)); 842*92a26154SJames Wright PetscCall(PetscSFSetGraph(newglob_to_oldglob_sf, global_vec_length, global_vec_length, NULL, PETSC_USE_POINTER, remote, PETSC_OWN_POINTER)); 843*92a26154SJames Wright PetscCall(DMGetNaturalSF(dm, &sfNatural_old)); 844*92a26154SJames Wright PetscCall(PetscSFCompose(newglob_to_oldglob_sf, sfNatural_old, &sfNatural_new)); 845*92a26154SJames Wright PetscCall(DMSetNaturalSF(dm, sfNatural_new)); 846*92a26154SJames Wright PetscCall(PetscSFDestroy(&sfNatural_new)); 847*92a26154SJames Wright PetscCall(PetscSFDestroy(&newglob_to_oldglob_sf)); 848*92a26154SJames Wright } 84953b735f8SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 85053b735f8SJames Wright } 85153b735f8SJames Wright 85253b735f8SJames Wright // Start with an SF for a positive depth (e.g., faces) and create a new SF for matched closure. 8536725e60dSJed Brown // 8545f06a3ddSJed Brown // Output Arguments: 8555f06a3ddSJed Brown // 8565f06a3ddSJed Brown // + closure_sf - augmented point SF (see `DMGetPointSF()`) that includes the faces and all points in its closure. This 8575f06a3ddSJed Brown // can be used to create a global section and section SF. 858b83f62b0SJames 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 8595f06a3ddSJed Brown // transformation to periodic dofs; see DMPeriodicCoordinateSetUp_Internal(). 8605f06a3ddSJed Brown // 861b83f62b0SJames Wright static PetscErrorCode DMPlexCreateIsoperiodicPointSF_Private(DM dm, PetscInt num_face_sfs, PetscSF *face_sfs, PetscSF *closure_sf, IS **is_points) 8626725e60dSJed Brown { 8636725e60dSJed Brown MPI_Comm comm; 8646725e60dSJed Brown PetscMPIInt rank; 8656725e60dSJed Brown PetscSF point_sf; 866b83f62b0SJames Wright PetscInt nroots, nleaves; 867b83f62b0SJames Wright const PetscInt *filocal; 868b83f62b0SJames Wright const PetscSFNode *firemote; 8696725e60dSJed Brown 8706725e60dSJed Brown PetscFunctionBegin; 8716725e60dSJed Brown PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 8726725e60dSJed Brown PetscCallMPI(MPI_Comm_rank(comm, &rank)); 8736725e60dSJed Brown PetscCall(DMGetPointSF(dm, &point_sf)); // Point SF has remote points 874b83f62b0SJames Wright PetscCall(PetscMalloc1(num_face_sfs, is_points)); 875b83f62b0SJames Wright 87653b735f8SJames Wright PetscCall(DMPlexCorrectOrientationForIsoperiodic(dm)); 87753b735f8SJames Wright 878b83f62b0SJames Wright for (PetscInt f = 0; f < num_face_sfs; f++) { 879b83f62b0SJames Wright PetscSF face_sf = face_sfs[f]; 88093defd67SJames Wright PetscInt *cl_sizes; 88193defd67SJames Wright PetscInt fStart, fEnd, rootbuffersize, leafbuffersize; 8826725e60dSJed Brown PetscSF sf_closure; 88393defd67SJames Wright PetscBT rootbt; 88493defd67SJames Wright 88593defd67SJames Wright PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd)); 88693defd67SJames Wright PetscCall(PetscMalloc1(fEnd - fStart, &cl_sizes)); 88793defd67SJames Wright for (PetscInt f = fStart, index = 0; f < fEnd; f++, index++) { 88893defd67SJames Wright PetscInt cl_size, *closure = NULL; 88993defd67SJames Wright PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &cl_size, &closure)); 89093defd67SJames Wright cl_sizes[index] = cl_size - 1; 89193defd67SJames Wright PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &cl_size, &closure)); 89293defd67SJames Wright } 89393defd67SJames Wright 89493defd67SJames Wright PetscCall(CreateDonorToPeriodicSF(dm, face_sf, fStart, fEnd, cl_sizes, &rootbuffersize, &leafbuffersize, &rootbt, &sf_closure)); 89593defd67SJames Wright PetscCall(PetscFree(cl_sizes)); 89693defd67SJames Wright PetscCall(PetscSFGetGraph(face_sf, &nroots, &nleaves, &filocal, &firemote)); 8976725e60dSJed Brown 898b83f62b0SJames Wright PetscSFNode *leaf_donor_closure; 899b83f62b0SJames Wright { // Pack root buffer with owner for every point in the root cones 9006725e60dSJed Brown PetscSFNode *donor_closure; 901b83f62b0SJames Wright const PetscInt *pilocal; 902b83f62b0SJames Wright const PetscSFNode *piremote; 903b83f62b0SJames Wright PetscInt npoints; 904b83f62b0SJames Wright 905b83f62b0SJames Wright PetscCall(PetscSFGetGraph(point_sf, NULL, &npoints, &pilocal, &piremote)); 90693defd67SJames Wright PetscCall(PetscCalloc1(rootbuffersize, &donor_closure)); 90793defd67SJames Wright for (PetscInt p = 0, root_offset = 0; p < nroots; p++) { 90893defd67SJames Wright if (!PetscBTLookup(rootbt, p)) continue; 90993defd67SJames Wright PetscInt cl_size, *closure = NULL; 9106725e60dSJed Brown PetscCall(DMPlexGetTransitiveClosure(dm, p, PETSC_TRUE, &cl_size, &closure)); 9116725e60dSJed Brown for (PetscInt j = 1; j < cl_size; j++) { 9126725e60dSJed Brown PetscInt c = closure[2 * j]; 9136725e60dSJed Brown if (pilocal) { 9146725e60dSJed Brown PetscInt found = -1; 9156725e60dSJed Brown if (npoints > 0) PetscCall(PetscFindInt(c, npoints, pilocal, &found)); 9166725e60dSJed Brown if (found >= 0) { 9176725e60dSJed Brown donor_closure[root_offset++] = piremote[found]; 9186725e60dSJed Brown continue; 9196725e60dSJed Brown } 9206725e60dSJed Brown } 9216725e60dSJed Brown // we own c 9226725e60dSJed Brown donor_closure[root_offset].rank = rank; 9236725e60dSJed Brown donor_closure[root_offset].index = c; 9246725e60dSJed Brown root_offset++; 9256725e60dSJed Brown } 9266725e60dSJed Brown PetscCall(DMPlexRestoreTransitiveClosure(dm, p, PETSC_TRUE, &cl_size, &closure)); 9276725e60dSJed Brown } 9286725e60dSJed Brown 92993defd67SJames Wright PetscCall(PetscMalloc1(leafbuffersize, &leaf_donor_closure)); 9306497c311SBarry Smith PetscCall(PetscSFBcastBegin(sf_closure, MPIU_SF_NODE, donor_closure, leaf_donor_closure, MPI_REPLACE)); 9316497c311SBarry Smith PetscCall(PetscSFBcastEnd(sf_closure, MPIU_SF_NODE, donor_closure, leaf_donor_closure, MPI_REPLACE)); 9326725e60dSJed Brown PetscCall(PetscSFDestroy(&sf_closure)); 9336725e60dSJed Brown PetscCall(PetscFree(donor_closure)); 934b83f62b0SJames Wright } 9356725e60dSJed Brown 9366725e60dSJed Brown PetscSFNode *new_iremote; 9376725e60dSJed Brown PetscCall(PetscCalloc1(nroots, &new_iremote)); 9386725e60dSJed Brown for (PetscInt i = 0; i < nroots; i++) new_iremote[i].rank = -1; 9396725e60dSJed Brown // Walk leaves and match vertices 94093defd67SJames Wright for (PetscInt i = 0, leaf_offset = 0; i < nleaves; i++) { 9416725e60dSJed Brown PetscInt point = filocal[i], cl_size; 9426725e60dSJed Brown PetscInt *closure = NULL; 9436725e60dSJed Brown PetscCall(DMPlexGetTransitiveClosure(dm, point, PETSC_TRUE, &cl_size, &closure)); 94453b735f8SJames Wright for (PetscInt j = 1; j < cl_size; j++) { 9456725e60dSJed Brown PetscInt c = closure[2 * j]; 9466725e60dSJed Brown PetscSFNode lc = leaf_donor_closure[leaf_offset]; 9476725e60dSJed Brown // printf("[%d] face %d.%d: %d ?-- (%d,%d)\n", rank, point, j, c, lc.rank, lc.index); 9486725e60dSJed Brown if (new_iremote[c].rank == -1) { 9496725e60dSJed Brown new_iremote[c] = lc; 9505f06a3ddSJed 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"); 9516725e60dSJed Brown leaf_offset++; 9526725e60dSJed Brown } 9536725e60dSJed Brown PetscCall(DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, &cl_size, &closure)); 9546725e60dSJed Brown } 9556725e60dSJed Brown PetscCall(PetscFree(leaf_donor_closure)); 9566725e60dSJed Brown 9576725e60dSJed Brown // Include face points in closure SF 9586725e60dSJed Brown for (PetscInt i = 0; i < nleaves; i++) new_iremote[filocal[i]] = firemote[i]; 9596725e60dSJed Brown // consolidate leaves 96093defd67SJames Wright PetscInt *leafdata; 96193defd67SJames Wright PetscCall(PetscMalloc1(nroots, &leafdata)); 9626725e60dSJed Brown PetscInt num_new_leaves = 0; 9636725e60dSJed Brown for (PetscInt i = 0; i < nroots; i++) { 9646725e60dSJed Brown if (new_iremote[i].rank == -1) continue; 9656725e60dSJed Brown new_iremote[num_new_leaves] = new_iremote[i]; 9666725e60dSJed Brown leafdata[num_new_leaves] = i; 9676725e60dSJed Brown num_new_leaves++; 9686725e60dSJed Brown } 969b83f62b0SJames Wright PetscCall(ISCreateGeneral(PETSC_COMM_SELF, num_new_leaves, leafdata, PETSC_COPY_VALUES, &(*is_points)[f])); 9706725e60dSJed Brown 9716725e60dSJed Brown PetscSF csf; 9726725e60dSJed Brown PetscCall(PetscSFCreate(comm, &csf)); 9736725e60dSJed Brown PetscCall(PetscSFSetGraph(csf, nroots, num_new_leaves, leafdata, PETSC_COPY_VALUES, new_iremote, PETSC_COPY_VALUES)); 9746725e60dSJed Brown PetscCall(PetscFree(new_iremote)); // copy and delete because new_iremote is longer than it needs to be 97593defd67SJames Wright PetscCall(PetscFree(leafdata)); 97693defd67SJames Wright PetscCall(PetscBTDestroy(&rootbt)); 9776725e60dSJed Brown 978b83f62b0SJames Wright PetscInt npoints; 979b83f62b0SJames Wright PetscCall(PetscSFGetGraph(point_sf, NULL, &npoints, NULL, NULL)); 9805f06a3ddSJed Brown if (npoints < 0) { // empty point_sf 9816725e60dSJed Brown *closure_sf = csf; 9825f06a3ddSJed Brown } else { 9835f06a3ddSJed Brown PetscCall(PetscSFMerge(point_sf, csf, closure_sf)); 9845f06a3ddSJed Brown PetscCall(PetscSFDestroy(&csf)); 9855f06a3ddSJed Brown } 986d8b4a066SPierre Jolivet if (f > 0) PetscCall(PetscSFDestroy(&point_sf)); // Only destroy if point_sf is from previous calls to PetscSFMerge 987b83f62b0SJames Wright point_sf = *closure_sf; // Use combined point + isoperiodic SF to define point ownership for further face_sf 988b83f62b0SJames Wright } 9895f06a3ddSJed Brown PetscCall(PetscObjectSetName((PetscObject)*closure_sf, "Composed Periodic Points")); 9903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9916725e60dSJed Brown } 9926725e60dSJed Brown 9935f06a3ddSJed Brown static PetscErrorCode DMGetIsoperiodicPointSF_Plex(DM dm, PetscSF *sf) 9946725e60dSJed Brown { 9956725e60dSJed Brown DM_Plex *plex = (DM_Plex *)dm->data; 9966725e60dSJed Brown 9976725e60dSJed Brown PetscFunctionBegin; 998b83f62b0SJames 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)); 9996725e60dSJed Brown if (sf) *sf = plex->periodic.composed_sf; 10003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10016725e60dSJed Brown } 10026725e60dSJed Brown 10035f06a3ddSJed Brown PetscErrorCode DMPlexMigrateIsoperiodicFaceSF_Internal(DM old_dm, DM dm, PetscSF sf_migration) 10045f06a3ddSJed Brown { 10055f06a3ddSJed Brown DM_Plex *plex = (DM_Plex *)old_dm->data; 1006a45b0079SJames Wright PetscSF sf_point, *new_face_sfs; 10075f06a3ddSJed Brown PetscMPIInt rank; 10085f06a3ddSJed Brown 10095f06a3ddSJed Brown PetscFunctionBegin; 10101fca310dSJames Wright if (!plex->periodic.face_sfs) PetscFunctionReturn(PETSC_SUCCESS); 10115f06a3ddSJed Brown PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 10125f06a3ddSJed Brown PetscCall(DMGetPointSF(dm, &sf_point)); 1013a45b0079SJames Wright PetscCall(PetscMalloc1(plex->periodic.num_face_sfs, &new_face_sfs)); 1014a45b0079SJames Wright 1015a45b0079SJames Wright for (PetscInt f = 0; f < plex->periodic.num_face_sfs; f++) { 10165f06a3ddSJed Brown PetscInt old_npoints, new_npoints, old_nleaf, new_nleaf, point_nleaf; 10175f06a3ddSJed Brown PetscSFNode *new_leafdata, *rootdata, *leafdata; 10185f06a3ddSJed Brown const PetscInt *old_local, *point_local; 10195f06a3ddSJed Brown const PetscSFNode *old_remote, *point_remote; 10206497c311SBarry Smith 1021a45b0079SJames Wright PetscCall(PetscSFGetGraph(plex->periodic.face_sfs[f], &old_npoints, &old_nleaf, &old_local, &old_remote)); 10225f06a3ddSJed Brown PetscCall(PetscSFGetGraph(sf_migration, NULL, &new_nleaf, NULL, NULL)); 10235f06a3ddSJed Brown PetscCall(PetscSFGetGraph(sf_point, &new_npoints, &point_nleaf, &point_local, &point_remote)); 10245f06a3ddSJed Brown PetscAssert(new_nleaf == new_npoints, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected migration leaf space to match new point root space"); 10255f06a3ddSJed Brown PetscCall(PetscMalloc3(old_npoints, &rootdata, old_npoints, &leafdata, new_npoints, &new_leafdata)); 10265f06a3ddSJed Brown 10275f06a3ddSJed Brown // Fill new_leafdata with new owners of all points 10285f06a3ddSJed Brown for (PetscInt i = 0; i < new_npoints; i++) { 10295f06a3ddSJed Brown new_leafdata[i].rank = rank; 10305f06a3ddSJed Brown new_leafdata[i].index = i; 10315f06a3ddSJed Brown } 10325f06a3ddSJed Brown for (PetscInt i = 0; i < point_nleaf; i++) { 10335f06a3ddSJed Brown PetscInt j = point_local[i]; 10345f06a3ddSJed Brown new_leafdata[j] = point_remote[i]; 10355f06a3ddSJed Brown } 10365f06a3ddSJed Brown // REPLACE is okay because every leaf agrees about the new owners 10376497c311SBarry Smith PetscCall(PetscSFReduceBegin(sf_migration, MPIU_SF_NODE, new_leafdata, rootdata, MPI_REPLACE)); 10386497c311SBarry Smith PetscCall(PetscSFReduceEnd(sf_migration, MPIU_SF_NODE, new_leafdata, rootdata, MPI_REPLACE)); 10395f06a3ddSJed Brown // rootdata now contains the new owners 10405f06a3ddSJed Brown 10415f06a3ddSJed Brown // Send to leaves of old space 10425f06a3ddSJed Brown for (PetscInt i = 0; i < old_npoints; i++) { 10435f06a3ddSJed Brown leafdata[i].rank = -1; 10445f06a3ddSJed Brown leafdata[i].index = -1; 10455f06a3ddSJed Brown } 10466497c311SBarry Smith PetscCall(PetscSFBcastBegin(plex->periodic.face_sfs[f], MPIU_SF_NODE, rootdata, leafdata, MPI_REPLACE)); 10476497c311SBarry Smith PetscCall(PetscSFBcastEnd(plex->periodic.face_sfs[f], MPIU_SF_NODE, rootdata, leafdata, MPI_REPLACE)); 10485f06a3ddSJed Brown 10495f06a3ddSJed Brown // Send to new leaf space 10506497c311SBarry Smith PetscCall(PetscSFBcastBegin(sf_migration, MPIU_SF_NODE, leafdata, new_leafdata, MPI_REPLACE)); 10516497c311SBarry Smith PetscCall(PetscSFBcastEnd(sf_migration, MPIU_SF_NODE, leafdata, new_leafdata, MPI_REPLACE)); 10525f06a3ddSJed Brown 10535f06a3ddSJed Brown PetscInt nface = 0, *new_local; 10545f06a3ddSJed Brown PetscSFNode *new_remote; 10555f06a3ddSJed Brown for (PetscInt i = 0; i < new_npoints; i++) nface += (new_leafdata[i].rank >= 0); 10565f06a3ddSJed Brown PetscCall(PetscMalloc1(nface, &new_local)); 10575f06a3ddSJed Brown PetscCall(PetscMalloc1(nface, &new_remote)); 10585f06a3ddSJed Brown nface = 0; 10595f06a3ddSJed Brown for (PetscInt i = 0; i < new_npoints; i++) { 10605f06a3ddSJed Brown if (new_leafdata[i].rank == -1) continue; 10615f06a3ddSJed Brown new_local[nface] = i; 10625f06a3ddSJed Brown new_remote[nface] = new_leafdata[i]; 10635f06a3ddSJed Brown nface++; 10645f06a3ddSJed Brown } 10655f06a3ddSJed Brown PetscCall(PetscFree3(rootdata, leafdata, new_leafdata)); 1066a45b0079SJames Wright PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), &new_face_sfs[f])); 1067a45b0079SJames Wright PetscCall(PetscSFSetGraph(new_face_sfs[f], new_npoints, nface, new_local, PETSC_OWN_POINTER, new_remote, PETSC_OWN_POINTER)); 1068a45b0079SJames Wright { 1069a45b0079SJames Wright char new_face_sf_name[PETSC_MAX_PATH_LEN]; 1070a45b0079SJames Wright PetscCall(PetscSNPrintf(new_face_sf_name, sizeof new_face_sf_name, "Migrated Isoperiodic Faces #%" PetscInt_FMT, f)); 1071a45b0079SJames Wright PetscCall(PetscObjectSetName((PetscObject)new_face_sfs[f], new_face_sf_name)); 1072a45b0079SJames Wright } 1073a45b0079SJames Wright } 1074a45b0079SJames Wright 1075a45b0079SJames Wright PetscCall(DMPlexSetIsoperiodicFaceSF(dm, plex->periodic.num_face_sfs, new_face_sfs)); 1076a45b0079SJames Wright PetscCall(DMPlexSetIsoperiodicFaceTransform(dm, plex->periodic.num_face_sfs, (PetscScalar *)plex->periodic.transform)); 1077a45b0079SJames Wright for (PetscInt f = 0; f < plex->periodic.num_face_sfs; f++) PetscCall(PetscSFDestroy(&new_face_sfs[f])); 1078a45b0079SJames Wright PetscCall(PetscFree(new_face_sfs)); 10793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10805f06a3ddSJed Brown } 10815f06a3ddSJed Brown 10826725e60dSJed Brown PetscErrorCode DMPeriodicCoordinateSetUp_Internal(DM dm) 10836725e60dSJed Brown { 10846725e60dSJed Brown DM_Plex *plex = (DM_Plex *)dm->data; 10856497c311SBarry Smith PetscCount count; 1086ddedc8f6SJames Wright IS isdof; 1087ddedc8f6SJames Wright PetscInt dim; 10884d86920dSPierre Jolivet 10896725e60dSJed Brown PetscFunctionBegin; 10901fca310dSJames Wright if (!plex->periodic.face_sfs) PetscFunctionReturn(PETSC_SUCCESS); 109153b735f8SJames Wright PetscCheck(plex->periodic.periodic_points, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Isoperiodic PointSF must be created before this function is called"); 10926725e60dSJed Brown 109353b735f8SJames Wright PetscCall(DMGetCoordinateDim(dm, &dim)); 1094ddedc8f6SJames Wright dm->periodic.num_affines = plex->periodic.num_face_sfs; 109553b735f8SJames Wright PetscCall(PetscFree2(dm->periodic.affine_to_local, dm->periodic.affine)); 1096ddedc8f6SJames Wright PetscCall(PetscMalloc2(dm->periodic.num_affines, &dm->periodic.affine_to_local, dm->periodic.num_affines, &dm->periodic.affine)); 1097ddedc8f6SJames Wright 1098ddedc8f6SJames Wright for (PetscInt f = 0; f < plex->periodic.num_face_sfs; f++) { 10996497c311SBarry Smith PetscInt npoints, vsize, isize; 11006725e60dSJed Brown const PetscInt *points; 1101ddedc8f6SJames Wright IS is = plex->periodic.periodic_points[f]; 11026725e60dSJed Brown PetscSegBuffer seg; 11036725e60dSJed Brown PetscSection section; 11046497c311SBarry Smith PetscInt *ind; 11056497c311SBarry Smith Vec L, P; 11066497c311SBarry Smith VecType vec_type; 11076497c311SBarry Smith VecScatter scatter; 11086497c311SBarry Smith PetscScalar *x; 11096497c311SBarry Smith 11106725e60dSJed Brown PetscCall(DMGetLocalSection(dm, §ion)); 11116725e60dSJed Brown PetscCall(PetscSegBufferCreate(sizeof(PetscInt), 32, &seg)); 11126725e60dSJed Brown PetscCall(ISGetSize(is, &npoints)); 11136725e60dSJed Brown PetscCall(ISGetIndices(is, &points)); 11146725e60dSJed Brown for (PetscInt i = 0; i < npoints; i++) { 11156725e60dSJed Brown PetscInt point = points[i], off, dof; 11166497c311SBarry Smith 11176725e60dSJed Brown PetscCall(PetscSectionGetOffset(section, point, &off)); 11186725e60dSJed Brown PetscCall(PetscSectionGetDof(section, point, &dof)); 11196725e60dSJed Brown PetscAssert(dof % dim == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unexpected dof %" PetscInt_FMT " not divisible by dimension %" PetscInt_FMT, dof, dim); 11206725e60dSJed Brown for (PetscInt j = 0; j < dof / dim; j++) { 11216725e60dSJed Brown PetscInt *slot; 11226497c311SBarry Smith 11236725e60dSJed Brown PetscCall(PetscSegBufferGetInts(seg, 1, &slot)); 1124729bdd54SJed Brown *slot = off / dim + j; 11256725e60dSJed Brown } 11266725e60dSJed Brown } 11276725e60dSJed Brown PetscCall(PetscSegBufferGetSize(seg, &count)); 11286725e60dSJed Brown PetscCall(PetscSegBufferExtractAlloc(seg, &ind)); 11296725e60dSJed Brown PetscCall(PetscSegBufferDestroy(&seg)); 11306497c311SBarry Smith PetscCall(PetscIntCast(count, &isize)); 11316497c311SBarry Smith PetscCall(ISCreateBlock(PETSC_COMM_SELF, dim, isize, ind, PETSC_OWN_POINTER, &isdof)); 11326497c311SBarry Smith 11336497c311SBarry Smith PetscCall(PetscIntCast(count * dim, &vsize)); 11346725e60dSJed Brown PetscCall(DMGetLocalVector(dm, &L)); 11356725e60dSJed Brown PetscCall(VecCreate(PETSC_COMM_SELF, &P)); 11366497c311SBarry Smith PetscCall(VecSetSizes(P, vsize, vsize)); 11376725e60dSJed Brown PetscCall(VecGetType(L, &vec_type)); 11386725e60dSJed Brown PetscCall(VecSetType(P, vec_type)); 11396725e60dSJed Brown PetscCall(VecScatterCreate(P, NULL, L, isdof, &scatter)); 11406725e60dSJed Brown PetscCall(DMRestoreLocalVector(dm, &L)); 11416725e60dSJed Brown PetscCall(ISDestroy(&isdof)); 11426725e60dSJed Brown 11436725e60dSJed Brown PetscCall(VecGetArrayWrite(P, &x)); 11446497c311SBarry Smith for (PetscCount i = 0; i < count; i++) { 1145ddedc8f6SJames Wright for (PetscInt j = 0; j < dim; j++) x[i * dim + j] = plex->periodic.transform[f][j][3]; 11466725e60dSJed Brown } 11476725e60dSJed Brown PetscCall(VecRestoreArrayWrite(P, &x)); 11486725e60dSJed Brown 1149ddedc8f6SJames Wright dm->periodic.affine_to_local[f] = scatter; 1150ddedc8f6SJames Wright dm->periodic.affine[f] = P; 1151ddedc8f6SJames Wright } 11526725e60dSJed Brown PetscCall(DMGlobalToLocalHookAdd(dm, NULL, DMCoordAddPeriodicOffsets_Private, NULL)); 11533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11546725e60dSJed Brown } 11556725e60dSJed Brown 11563e72e933SJed Brown PetscErrorCode DMPlexCreateBoxMesh_Tensor_SFC_Internal(DM dm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate) 11573e72e933SJed Brown { 11583e72e933SJed Brown PetscInt eextent[3] = {1, 1, 1}, vextent[3] = {1, 1, 1}; 11593e72e933SJed Brown const Ijk closure_1[] = { 11603e72e933SJed Brown {0, 0, 0}, 11613e72e933SJed Brown {1, 0, 0}, 11623e72e933SJed Brown }; 11633e72e933SJed Brown const Ijk closure_2[] = { 11643e72e933SJed Brown {0, 0, 0}, 11653e72e933SJed Brown {1, 0, 0}, 11663e72e933SJed Brown {1, 1, 0}, 11673e72e933SJed Brown {0, 1, 0}, 11683e72e933SJed Brown }; 11693e72e933SJed Brown const Ijk closure_3[] = { 11703e72e933SJed Brown {0, 0, 0}, 11713e72e933SJed Brown {0, 1, 0}, 11723e72e933SJed Brown {1, 1, 0}, 11733e72e933SJed Brown {1, 0, 0}, 11743e72e933SJed Brown {0, 0, 1}, 11753e72e933SJed Brown {1, 0, 1}, 11763e72e933SJed Brown {1, 1, 1}, 11773e72e933SJed Brown {0, 1, 1}, 11783e72e933SJed Brown }; 11793431e603SJed Brown const Ijk *const closure_dim[] = {NULL, closure_1, closure_2, closure_3}; 11803431e603SJed Brown // This must be kept consistent with DMPlexCreateCubeMesh_Internal 11813431e603SJed Brown const PetscInt face_marker_1[] = {1, 2}; 11823431e603SJed Brown const PetscInt face_marker_2[] = {4, 2, 1, 3}; 11833431e603SJed Brown const PetscInt face_marker_3[] = {6, 5, 3, 4, 1, 2}; 11843431e603SJed Brown const PetscInt *const face_marker_dim[] = {NULL, face_marker_1, face_marker_2, face_marker_3}; 11856725e60dSJed Brown // Orient faces so the normal is in the positive axis and the first vertex is the one closest to zero. 11866725e60dSJed Brown // These orientations can be determined by examining cones of a reference quad and hex element. 11876725e60dSJed Brown const PetscInt face_orient_1[] = {0, 0}; 11886725e60dSJed Brown const PetscInt face_orient_2[] = {-1, 0, 0, -1}; 11896725e60dSJed Brown const PetscInt face_orient_3[] = {-2, 0, -2, 1, -2, 0}; 11906725e60dSJed Brown const PetscInt *const face_orient_dim[] = {NULL, face_orient_1, face_orient_2, face_orient_3}; 11913e72e933SJed Brown 11923e72e933SJed Brown PetscFunctionBegin; 1193ea7b3863SJames Wright PetscCall(PetscLogEventBegin(DMPLEX_CreateBoxSFC, dm, 0, 0, 0)); 11944f572ea9SToby Isaac PetscAssertPointer(dm, 1); 11953e72e933SJed Brown PetscValidLogicalCollectiveInt(dm, dim, 2); 11963e72e933SJed Brown PetscCall(DMSetDimension(dm, dim)); 11973e72e933SJed Brown PetscMPIInt rank, size; 11983e72e933SJed Brown PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 11993e72e933SJed Brown PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 12003e72e933SJed Brown for (PetscInt i = 0; i < dim; i++) { 12013e72e933SJed Brown eextent[i] = faces[i]; 12023e72e933SJed Brown vextent[i] = faces[i] + 1; 12033e72e933SJed Brown } 1204c77877e3SJed Brown ZLayout layout; 1205c77877e3SJed Brown PetscCall(ZLayoutCreate(size, eextent, vextent, &layout)); 12063e72e933SJed Brown PetscZSet vset; // set of all vertices in the closure of the owned elements 12073e72e933SJed Brown PetscCall(PetscZSetCreate(&vset)); 12083e72e933SJed Brown PetscInt local_elems = 0; 12093e72e933SJed Brown for (ZCode z = layout.zstarts[rank]; z < layout.zstarts[rank + 1]; z++) { 12103e72e933SJed Brown Ijk loc = ZCodeSplit(z); 12113ba16761SJacob Faibussowitsch if (IjkActive(layout.vextent, loc)) PetscCall(PetscZSetAdd(vset, z)); 12126547ddbdSJed Brown else { 12136547ddbdSJed Brown z += ZStepOct(z); 12146547ddbdSJed Brown continue; 12156547ddbdSJed Brown } 12163e72e933SJed Brown if (IjkActive(layout.eextent, loc)) { 12173e72e933SJed Brown local_elems++; 12183e72e933SJed Brown // Add all neighboring vertices to set 12193e72e933SJed Brown for (PetscInt n = 0; n < PetscPowInt(2, dim); n++) { 12203e72e933SJed Brown Ijk inc = closure_dim[dim][n]; 12213e72e933SJed Brown Ijk nloc = {loc.i + inc.i, loc.j + inc.j, loc.k + inc.k}; 12223e72e933SJed Brown ZCode v = ZEncode(nloc); 12233ba16761SJacob Faibussowitsch PetscCall(PetscZSetAdd(vset, v)); 12243e72e933SJed Brown } 12253e72e933SJed Brown } 12263e72e933SJed Brown } 12273e72e933SJed Brown PetscInt local_verts, off = 0; 12283e72e933SJed Brown ZCode *vert_z; 12293e72e933SJed Brown PetscCall(PetscZSetGetSize(vset, &local_verts)); 12303e72e933SJed Brown PetscCall(PetscMalloc1(local_verts, &vert_z)); 12313e72e933SJed Brown PetscCall(PetscZSetGetElems(vset, &off, vert_z)); 12323e72e933SJed Brown PetscCall(PetscZSetDestroy(&vset)); 12333e72e933SJed Brown // ZCode is unsigned for bitwise convenience, but highest bit should never be set, so can interpret as signed 12343e72e933SJed Brown PetscCall(PetscSortInt64(local_verts, (PetscInt64 *)vert_z)); 12353e72e933SJed Brown 12363e72e933SJed Brown PetscCall(DMPlexSetChart(dm, 0, local_elems + local_verts)); 12373e72e933SJed Brown for (PetscInt e = 0; e < local_elems; e++) PetscCall(DMPlexSetConeSize(dm, e, PetscPowInt(2, dim))); 12383e72e933SJed Brown PetscCall(DMSetUp(dm)); 12393e72e933SJed Brown { 12403e72e933SJed Brown PetscInt e = 0; 12413e72e933SJed Brown for (ZCode z = layout.zstarts[rank]; z < layout.zstarts[rank + 1]; z++) { 12423e72e933SJed Brown Ijk loc = ZCodeSplit(z); 12436547ddbdSJed Brown if (!IjkActive(layout.eextent, loc)) { 12446547ddbdSJed Brown z += ZStepOct(z); 12456547ddbdSJed Brown continue; 12466547ddbdSJed Brown } 12473e72e933SJed Brown PetscInt cone[8], orient[8] = {0}; 12483e72e933SJed Brown for (PetscInt n = 0; n < PetscPowInt(2, dim); n++) { 12493e72e933SJed Brown Ijk inc = closure_dim[dim][n]; 12503e72e933SJed Brown Ijk nloc = {loc.i + inc.i, loc.j + inc.j, loc.k + inc.k}; 12513e72e933SJed Brown ZCode v = ZEncode(nloc); 12523e72e933SJed Brown PetscInt ci = ZCodeFind(v, local_verts, vert_z); 12533e72e933SJed Brown PetscAssert(ci >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Did not find neighbor vertex in set"); 12543e72e933SJed Brown cone[n] = local_elems + ci; 12553e72e933SJed Brown } 12563e72e933SJed Brown PetscCall(DMPlexSetCone(dm, e, cone)); 12573e72e933SJed Brown PetscCall(DMPlexSetConeOrientation(dm, e, orient)); 12583e72e933SJed Brown e++; 12593e72e933SJed Brown } 12603e72e933SJed Brown } 12613e72e933SJed Brown 12623e72e933SJed Brown PetscCall(DMPlexSymmetrize(dm)); 12633e72e933SJed Brown PetscCall(DMPlexStratify(dm)); 12646725e60dSJed Brown 12653e72e933SJed Brown { // Create point SF 12663e72e933SJed Brown PetscSF sf; 12673ba16761SJacob Faibussowitsch PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), &sf)); 12683e72e933SJed Brown PetscInt owned_verts = ZCodeFind(layout.zstarts[rank + 1], local_verts, vert_z); 12693e72e933SJed Brown if (owned_verts < 0) owned_verts = -(owned_verts + 1); // We don't care whether the key was found 12703e72e933SJed Brown PetscInt num_ghosts = local_verts - owned_verts; // Due to sorting, owned vertices always come first 12713e72e933SJed Brown PetscInt *local_ghosts; 12723e72e933SJed Brown PetscSFNode *ghosts; 12733e72e933SJed Brown PetscCall(PetscMalloc1(num_ghosts, &local_ghosts)); 12743e72e933SJed Brown PetscCall(PetscMalloc1(num_ghosts, &ghosts)); 12753e72e933SJed Brown for (PetscInt i = 0; i < num_ghosts;) { 12763e72e933SJed Brown ZCode z = vert_z[owned_verts + i]; 1277835f2295SStefano Zampini PetscMPIInt remote_rank, remote_count = 0; 1278835f2295SStefano Zampini 1279835f2295SStefano Zampini PetscCall(PetscMPIIntCast(ZCodeFind(z, size + 1, layout.zstarts), &remote_rank)); 12803e72e933SJed Brown if (remote_rank < 0) remote_rank = -(remote_rank + 1) - 1; 12813e72e933SJed Brown // We have a new remote rank; find all the ghost indices (which are contiguous in vert_z) 12823e72e933SJed Brown 12833e72e933SJed Brown // Count the elements on remote_rank 12844e2e9504SJed Brown PetscInt remote_elem = ZLayoutElementsOnRank(&layout, remote_rank); 12853e72e933SJed Brown 12863e72e933SJed Brown // Traverse vertices and make ghost links 12873e72e933SJed Brown for (ZCode rz = layout.zstarts[remote_rank]; rz < layout.zstarts[remote_rank + 1]; rz++) { 12883e72e933SJed Brown Ijk loc = ZCodeSplit(rz); 12893e72e933SJed Brown if (rz == z) { 12903e72e933SJed Brown local_ghosts[i] = local_elems + owned_verts + i; 12913e72e933SJed Brown ghosts[i].rank = remote_rank; 12923e72e933SJed Brown ghosts[i].index = remote_elem + remote_count; 12933e72e933SJed Brown i++; 12943e72e933SJed Brown if (i == num_ghosts) break; 12953e72e933SJed Brown z = vert_z[owned_verts + i]; 12963e72e933SJed Brown } 12973e72e933SJed Brown if (IjkActive(layout.vextent, loc)) remote_count++; 12986547ddbdSJed Brown else rz += ZStepOct(rz); 12993e72e933SJed Brown } 13003e72e933SJed Brown } 13013e72e933SJed Brown PetscCall(PetscSFSetGraph(sf, local_elems + local_verts, num_ghosts, local_ghosts, PETSC_OWN_POINTER, ghosts, PETSC_OWN_POINTER)); 13023e72e933SJed Brown PetscCall(PetscObjectSetName((PetscObject)sf, "SFC Point SF")); 13033e72e933SJed Brown PetscCall(DMSetPointSF(dm, sf)); 13043e72e933SJed Brown PetscCall(PetscSFDestroy(&sf)); 13053e72e933SJed Brown } 13063e72e933SJed Brown { 13073e72e933SJed Brown Vec coordinates; 13083e72e933SJed Brown PetscScalar *coords; 13093e72e933SJed Brown PetscSection coord_section; 13103e72e933SJed Brown PetscInt coord_size; 13113e72e933SJed Brown PetscCall(DMGetCoordinateSection(dm, &coord_section)); 13123e72e933SJed Brown PetscCall(PetscSectionSetNumFields(coord_section, 1)); 13133e72e933SJed Brown PetscCall(PetscSectionSetFieldComponents(coord_section, 0, dim)); 13143e72e933SJed Brown PetscCall(PetscSectionSetChart(coord_section, local_elems, local_elems + local_verts)); 13153e72e933SJed Brown for (PetscInt v = 0; v < local_verts; v++) { 13163e72e933SJed Brown PetscInt point = local_elems + v; 13173e72e933SJed Brown PetscCall(PetscSectionSetDof(coord_section, point, dim)); 13183e72e933SJed Brown PetscCall(PetscSectionSetFieldDof(coord_section, point, 0, dim)); 13193e72e933SJed Brown } 13203e72e933SJed Brown PetscCall(PetscSectionSetUp(coord_section)); 13213e72e933SJed Brown PetscCall(PetscSectionGetStorageSize(coord_section, &coord_size)); 13223e72e933SJed Brown PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 13233e72e933SJed Brown PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 13243e72e933SJed Brown PetscCall(VecSetSizes(coordinates, coord_size, PETSC_DETERMINE)); 13253e72e933SJed Brown PetscCall(VecSetBlockSize(coordinates, dim)); 13263e72e933SJed Brown PetscCall(VecSetType(coordinates, VECSTANDARD)); 13273e72e933SJed Brown PetscCall(VecGetArray(coordinates, &coords)); 13283e72e933SJed Brown for (PetscInt v = 0; v < local_verts; v++) { 13293e72e933SJed Brown Ijk loc = ZCodeSplit(vert_z[v]); 13303e72e933SJed Brown coords[v * dim + 0] = lower[0] + loc.i * (upper[0] - lower[0]) / layout.eextent.i; 13313e72e933SJed Brown if (dim > 1) coords[v * dim + 1] = lower[1] + loc.j * (upper[1] - lower[1]) / layout.eextent.j; 13323e72e933SJed Brown if (dim > 2) coords[v * dim + 2] = lower[2] + loc.k * (upper[2] - lower[2]) / layout.eextent.k; 13333e72e933SJed Brown } 13343e72e933SJed Brown PetscCall(VecRestoreArray(coordinates, &coords)); 13353e72e933SJed Brown PetscCall(DMSetCoordinatesLocal(dm, coordinates)); 13363e72e933SJed Brown PetscCall(VecDestroy(&coordinates)); 13373e72e933SJed Brown } 13383e72e933SJed Brown if (interpolate) { 13393431e603SJed Brown PetscCall(DMPlexInterpolateInPlace_Internal(dm)); 13403431e603SJed Brown 13413431e603SJed Brown DMLabel label; 13423431e603SJed Brown PetscCall(DMCreateLabel(dm, "Face Sets")); 13433431e603SJed Brown PetscCall(DMGetLabel(dm, "Face Sets", &label)); 13449ae3492bSJames Wright PetscSegBuffer per_faces[3], donor_face_closure[3], my_donor_faces[3]; 13459ae3492bSJames Wright for (PetscInt i = 0; i < 3; i++) { 13469ae3492bSJames Wright PetscCall(PetscSegBufferCreate(sizeof(PetscInt), 64, &per_faces[i])); 13479ae3492bSJames Wright PetscCall(PetscSegBufferCreate(sizeof(PetscInt), 64, &my_donor_faces[i])); 13489ae3492bSJames Wright PetscCall(PetscSegBufferCreate(sizeof(ZCode), 64 * PetscPowInt(2, dim), &donor_face_closure[i])); 13499ae3492bSJames Wright } 13503431e603SJed Brown PetscInt fStart, fEnd, vStart, vEnd; 13513431e603SJed Brown PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd)); 13523431e603SJed Brown PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 13533431e603SJed Brown for (PetscInt f = fStart; f < fEnd; f++) { 13544e2e9504SJed Brown PetscInt npoints, *points = NULL, num_fverts = 0, fverts[8]; 13553431e603SJed Brown PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &npoints, &points)); 13563431e603SJed Brown PetscInt bc_count[6] = {0}; 13573431e603SJed Brown for (PetscInt i = 0; i < npoints; i++) { 13583431e603SJed Brown PetscInt p = points[2 * i]; 135985de0153SJames Wright if (!IsPointInsideStratum(p, vStart, vEnd)) continue; 13604e2e9504SJed Brown fverts[num_fverts++] = p; 13613431e603SJed Brown Ijk loc = ZCodeSplit(vert_z[p - vStart]); 13623431e603SJed Brown // Convention here matches DMPlexCreateCubeMesh_Internal 13633431e603SJed Brown bc_count[0] += loc.i == 0; 13643431e603SJed Brown bc_count[1] += loc.i == layout.vextent.i - 1; 13653431e603SJed Brown bc_count[2] += loc.j == 0; 13663431e603SJed Brown bc_count[3] += loc.j == layout.vextent.j - 1; 13673431e603SJed Brown bc_count[4] += loc.k == 0; 13683431e603SJed Brown bc_count[5] += loc.k == layout.vextent.k - 1; 13693e72e933SJed Brown } 13703431e603SJed Brown PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &npoints, &points)); 13713431e603SJed Brown for (PetscInt bc = 0, bc_match = 0; bc < 2 * dim; bc++) { 13723431e603SJed Brown if (bc_count[bc] == PetscPowInt(2, dim - 1)) { 13736725e60dSJed Brown PetscCall(DMPlexOrientPoint(dm, f, face_orient_dim[dim][bc])); 13744e2e9504SJed Brown if (periodicity[bc / 2] == DM_BOUNDARY_PERIODIC) { 13754e2e9504SJed Brown PetscInt *put; 13764e2e9504SJed Brown if (bc % 2 == 0) { // donor face; no label 13779ae3492bSJames Wright PetscCall(PetscSegBufferGet(my_donor_faces[bc / 2], 1, &put)); 13784e2e9504SJed Brown *put = f; 13794e2e9504SJed Brown } else { // periodic face 13809ae3492bSJames Wright PetscCall(PetscSegBufferGet(per_faces[bc / 2], 1, &put)); 13814e2e9504SJed Brown *put = f; 13824e2e9504SJed Brown ZCode *zput; 13839ae3492bSJames Wright PetscCall(PetscSegBufferGet(donor_face_closure[bc / 2], num_fverts, &zput)); 13844e2e9504SJed Brown for (PetscInt i = 0; i < num_fverts; i++) { 13854e2e9504SJed Brown Ijk loc = ZCodeSplit(vert_z[fverts[i] - vStart]); 13864e2e9504SJed Brown switch (bc / 2) { 13874e2e9504SJed Brown case 0: 13884e2e9504SJed Brown loc.i = 0; 13894e2e9504SJed Brown break; 13904e2e9504SJed Brown case 1: 13914e2e9504SJed Brown loc.j = 0; 13924e2e9504SJed Brown break; 13934e2e9504SJed Brown case 2: 13944e2e9504SJed Brown loc.k = 0; 13954e2e9504SJed Brown break; 13964e2e9504SJed Brown } 13974e2e9504SJed Brown *zput++ = ZEncode(loc); 13984e2e9504SJed Brown } 13994e2e9504SJed Brown } 14004e2e9504SJed Brown continue; 14014e2e9504SJed Brown } 14023431e603SJed Brown PetscAssert(bc_match == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face matches multiple face sets"); 14033431e603SJed Brown PetscCall(DMLabelSetValue(label, f, face_marker_dim[dim][bc])); 14043431e603SJed Brown bc_match++; 14053431e603SJed Brown } 14063431e603SJed Brown } 14073431e603SJed Brown } 14083431e603SJed Brown // Ensure that the Coordinate DM has our new boundary labels 14093431e603SJed Brown DM cdm; 14103431e603SJed Brown PetscCall(DMGetCoordinateDM(dm, &cdm)); 14113431e603SJed Brown PetscCall(DMCopyLabels(dm, cdm, PETSC_COPY_VALUES, PETSC_FALSE, DM_COPY_LABELS_FAIL)); 14126725e60dSJed Brown if (periodicity[0] == DM_BOUNDARY_PERIODIC || (dim > 1 && periodicity[1] == DM_BOUNDARY_PERIODIC) || (dim > 2 && periodicity[2] == DM_BOUNDARY_PERIODIC)) { 14135f06a3ddSJed Brown PetscCall(DMPlexCreateBoxMesh_Tensor_SFC_Periodicity_Private(dm, &layout, vert_z, per_faces, lower, upper, periodicity, donor_face_closure, my_donor_faces)); 14146725e60dSJed Brown } 14159ae3492bSJames Wright for (PetscInt i = 0; i < 3; i++) { 14169ae3492bSJames Wright PetscCall(PetscSegBufferDestroy(&per_faces[i])); 14179ae3492bSJames Wright PetscCall(PetscSegBufferDestroy(&donor_face_closure[i])); 14189ae3492bSJames Wright PetscCall(PetscSegBufferDestroy(&my_donor_faces[i])); 14199ae3492bSJames Wright } 14203431e603SJed Brown } 14214e2e9504SJed Brown PetscCall(PetscFree(layout.zstarts)); 14223431e603SJed Brown PetscCall(PetscFree(vert_z)); 1423ea7b3863SJames Wright PetscCall(PetscLogEventEnd(DMPLEX_CreateBoxSFC, dm, 0, 0, 0)); 14243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14253e72e933SJed Brown } 14264e2e9504SJed Brown 14274e2e9504SJed Brown /*@ 14285f06a3ddSJed Brown DMPlexSetIsoperiodicFaceSF - Express periodicity from an existing mesh 14294e2e9504SJed Brown 143020f4b53cSBarry Smith Logically Collective 14314e2e9504SJed Brown 14324e2e9504SJed Brown Input Parameters: 14334e2e9504SJed Brown + dm - The `DMPLEX` on which to set periodicity 14341fca310dSJames Wright . num_face_sfs - Number of `PetscSF`s in `face_sfs` 14351fca310dSJames 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. 14364e2e9504SJed Brown 14374e2e9504SJed Brown Level: advanced 14384e2e9504SJed Brown 143953c0d4aeSBarry Smith Note: 14405dca41c3SJed Brown One can use `-dm_plex_shape zbox` to use this mode of periodicity, wherein the periodic points are distinct both globally 14414e2e9504SJed Brown and locally, but are paired when creating a global dof space. 14424e2e9504SJed Brown 14431cc06b55SBarry Smith .seealso: [](ch_unstructured), `DMPLEX`, `DMGetGlobalSection()`, `DMPlexGetIsoperiodicFaceSF()` 14444e2e9504SJed Brown @*/ 14451fca310dSJames Wright PetscErrorCode DMPlexSetIsoperiodicFaceSF(DM dm, PetscInt num_face_sfs, PetscSF *face_sfs) 14464e2e9504SJed Brown { 14474e2e9504SJed Brown DM_Plex *plex = (DM_Plex *)dm->data; 14484d86920dSPierre Jolivet 14494e2e9504SJed Brown PetscFunctionBegin; 14504e2e9504SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1451545efb26SJames Wright if (num_face_sfs) { 1452545efb26SJames Wright PetscAssertPointer(face_sfs, 3); 1453545efb26SJames Wright PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMGetIsoperiodicPointSF_C", DMGetIsoperiodicPointSF_Plex)); 1454545efb26SJames Wright } else PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMGetIsoperiodicPointSF_C", NULL)); 14552b4f33d9SJames Wright if (num_face_sfs == plex->periodic.num_face_sfs && (num_face_sfs == 0 || face_sfs == plex->periodic.face_sfs)) PetscFunctionReturn(PETSC_SUCCESS); 14562b4f33d9SJames Wright PetscCall(DMSetGlobalSection(dm, NULL)); 14571fca310dSJames Wright 14581fca310dSJames Wright for (PetscInt i = 0; i < num_face_sfs; i++) PetscCall(PetscObjectReference((PetscObject)face_sfs[i])); 14591fca310dSJames Wright if (plex->periodic.num_face_sfs > 0) { 14601fca310dSJames Wright for (PetscInt i = 0; i < plex->periodic.num_face_sfs; i++) PetscCall(PetscSFDestroy(&plex->periodic.face_sfs[i])); 14611fca310dSJames Wright PetscCall(PetscFree(plex->periodic.face_sfs)); 14621fca310dSJames Wright } 14631fca310dSJames Wright 14641fca310dSJames Wright plex->periodic.num_face_sfs = num_face_sfs; 14651fca310dSJames Wright PetscCall(PetscCalloc1(num_face_sfs, &plex->periodic.face_sfs)); 14661fca310dSJames Wright for (PetscInt i = 0; i < num_face_sfs; i++) plex->periodic.face_sfs[i] = face_sfs[i]; 14675f06a3ddSJed Brown 14685f06a3ddSJed Brown DM cdm = dm->coordinates[0].dm; // Can't DMGetCoordinateDM because it automatically creates one 14695f06a3ddSJed Brown if (cdm) { 14701fca310dSJames Wright PetscCall(DMPlexSetIsoperiodicFaceSF(cdm, num_face_sfs, face_sfs)); 14711fca310dSJames Wright if (face_sfs) cdm->periodic.setup = DMPeriodicCoordinateSetUp_Internal; 14725f06a3ddSJed Brown } 14733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14744e2e9504SJed Brown } 14754e2e9504SJed Brown 14761fca310dSJames Wright /*@C 14775f06a3ddSJed Brown DMPlexGetIsoperiodicFaceSF - Obtain periodicity for a mesh 14784e2e9504SJed Brown 147920f4b53cSBarry Smith Logically Collective 14804e2e9504SJed Brown 148120f4b53cSBarry Smith Input Parameter: 14824e2e9504SJed Brown . dm - The `DMPLEX` for which to obtain periodic relation 14834e2e9504SJed Brown 14849cde84edSJose E. Roman Output Parameters: 14851fca310dSJames Wright + num_face_sfs - Number of `PetscSF`s in the array 14861fca310dSJames 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. 14874e2e9504SJed Brown 14884e2e9504SJed Brown Level: advanced 14894e2e9504SJed Brown 14901cc06b55SBarry Smith .seealso: [](ch_unstructured), `DMPLEX`, `DMGetGlobalSection()`, `DMPlexSetIsoperiodicFaceSF()` 14914e2e9504SJed Brown @*/ 14921fca310dSJames Wright PetscErrorCode DMPlexGetIsoperiodicFaceSF(DM dm, PetscInt *num_face_sfs, const PetscSF **face_sfs) 14934e2e9504SJed Brown { 149481316d9bSJames Wright PetscBool isPlex; 14954d86920dSPierre Jolivet 14964e2e9504SJed Brown PetscFunctionBegin; 14974e2e9504SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 149881316d9bSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isPlex)); 149981316d9bSJames Wright if (isPlex) { 150081316d9bSJames Wright DM_Plex *plex = (DM_Plex *)dm->data; 150181316d9bSJames Wright if (face_sfs) *face_sfs = plex->periodic.face_sfs; 150281316d9bSJames Wright if (num_face_sfs) *num_face_sfs = plex->periodic.num_face_sfs; 150381316d9bSJames Wright } else { 150481316d9bSJames Wright if (face_sfs) *face_sfs = NULL; 150581316d9bSJames Wright if (num_face_sfs) *num_face_sfs = 0; 150681316d9bSJames Wright } 15073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15084e2e9504SJed Brown } 15096725e60dSJed Brown 15105f06a3ddSJed Brown /*@C 15115f06a3ddSJed Brown DMPlexSetIsoperiodicFaceTransform - set geometric transform from donor to periodic points 15126725e60dSJed Brown 15136725e60dSJed Brown Logically Collective 15146725e60dSJed Brown 151560225df5SJacob Faibussowitsch Input Parameters: 151653c0d4aeSBarry Smith + dm - `DMPLEX` that has been configured with `DMPlexSetIsoperiodicFaceSF()` 15171fca310dSJames Wright . n - Number of transforms in array 15181fca310dSJames Wright - t - Array of 4x4 affine transformation basis. 15196725e60dSJed Brown 152053c0d4aeSBarry Smith Level: advanced 152153c0d4aeSBarry Smith 15225f06a3ddSJed Brown Notes: 15235f06a3ddSJed Brown Affine transforms are 4x4 matrices in which the leading 3x3 block expresses a rotation (or identity for no rotation), 15245f06a3ddSJed Brown the last column contains a translation vector, and the bottom row is all zero except the last entry, which must always 15255f06a3ddSJed Brown be 1. This representation is common in geometric modeling and allows affine transformations to be composed using 1526aaa8cc7dSPierre Jolivet simple matrix multiplication. 15275f06a3ddSJed Brown 15285f06a3ddSJed Brown Although the interface accepts a general affine transform, only affine translation is supported at present. 15295f06a3ddSJed Brown 153060225df5SJacob Faibussowitsch Developer Notes: 15315f06a3ddSJed Brown This interface should be replaced by making BasisTransform public, expanding it to support affine representations, and 15325f06a3ddSJed Brown adding GPU implementations to apply the G2L/L2G transforms. 153353c0d4aeSBarry Smith 15341cc06b55SBarry Smith .seealso: [](ch_unstructured), `DMPLEX`, `DMGetGlobalSection()`, `DMPlexSetIsoperiodicFaceSF()` 15356725e60dSJed Brown @*/ 1536cc4c1da9SBarry Smith PetscErrorCode DMPlexSetIsoperiodicFaceTransform(DM dm, PetscInt n, const PetscScalar t[]) 15376725e60dSJed Brown { 15386725e60dSJed Brown DM_Plex *plex = (DM_Plex *)dm->data; 15394d86920dSPierre Jolivet 15406725e60dSJed Brown PetscFunctionBegin; 15416725e60dSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 15421fca310dSJames 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); 15431fca310dSJames Wright 1544d7d2d1d2SJames Wright PetscCall(PetscFree(plex->periodic.transform)); 15451fca310dSJames Wright PetscCall(PetscMalloc1(n, &plex->periodic.transform)); 15461fca310dSJames Wright for (PetscInt i = 0; i < n; i++) { 15476725e60dSJed Brown for (PetscInt j = 0; j < 4; j++) { 15481fca310dSJames Wright for (PetscInt k = 0; k < 4; k++) { 15491fca310dSJames Wright PetscCheck(j != k || t[i * 16 + j * 4 + k] == 1., PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Rotated transforms not supported"); 15501fca310dSJames Wright plex->periodic.transform[i][j][k] = t[i * 16 + j * 4 + k]; 15511fca310dSJames Wright } 15526725e60dSJed Brown } 15536725e60dSJed Brown } 15543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15556725e60dSJed Brown } 1556