13e72e933SJed Brown #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 23e72e933SJed Brown #include <petscsf.h> 33e72e933SJed Brown #include <petsc/private/hashset.h> 43e72e933SJed Brown 53e72e933SJed Brown typedef uint64_t ZCode; 63e72e933SJed Brown 73e72e933SJed Brown PETSC_HASH_SET(ZSet, ZCode, PetscHash_UInt64, PetscHashEqual) 83e72e933SJed Brown 93e72e933SJed Brown typedef struct { 103e72e933SJed Brown PetscInt i, j, k; 113e72e933SJed Brown } Ijk; 123e72e933SJed Brown 133e72e933SJed Brown typedef struct { 143e72e933SJed Brown Ijk eextent; 153e72e933SJed Brown Ijk vextent; 163e72e933SJed Brown PetscMPIInt comm_size; 173e72e933SJed Brown ZCode *zstarts; 183e72e933SJed Brown } ZLayout; 193e72e933SJed Brown 203e72e933SJed Brown static unsigned ZCodeSplit1(ZCode z) 213e72e933SJed Brown { 223e72e933SJed Brown z = ((z & 01001001001001001) | ((z >> 2) & 02002002002002002) | ((z >> 4) & 04004004004004004)); 233e72e933SJed Brown z = (z | (z >> 6) | (z >> 12)) & 0000000777000000777; 243e72e933SJed Brown z = (z | (z >> 18)) & 0777777; 253e72e933SJed Brown return (unsigned)z; 263e72e933SJed Brown } 273e72e933SJed Brown 283e72e933SJed Brown static ZCode ZEncode1(unsigned t) 293e72e933SJed Brown { 303e72e933SJed Brown ZCode z = t; 313e72e933SJed Brown z = (z | (z << 18)) & 0777000000777; 323e72e933SJed Brown z = (z | (z << 6) | (z << 12)) & 07007007007007007; 333e72e933SJed Brown z = (z | (z << 2) | (z << 4)) & 0111111111111111111; 343e72e933SJed Brown return z; 353e72e933SJed Brown } 363e72e933SJed Brown 373e72e933SJed Brown static Ijk ZCodeSplit(ZCode z) 383e72e933SJed Brown { 393e72e933SJed Brown Ijk c; 403e72e933SJed Brown c.i = ZCodeSplit1(z >> 2); 413e72e933SJed Brown c.j = ZCodeSplit1(z >> 1); 423e72e933SJed Brown c.k = ZCodeSplit1(z >> 0); 433e72e933SJed Brown return c; 443e72e933SJed Brown } 453e72e933SJed Brown 463e72e933SJed Brown static ZCode ZEncode(Ijk c) 473e72e933SJed Brown { 483e72e933SJed Brown ZCode z = (ZEncode1(c.i) << 2) | (ZEncode1(c.j) << 1) | ZEncode1(c.k); 493e72e933SJed Brown return z; 503e72e933SJed Brown } 513e72e933SJed Brown 523e72e933SJed Brown static PetscBool IjkActive(Ijk extent, Ijk l) 533e72e933SJed Brown { 543e72e933SJed Brown if (l.i < extent.i && l.j < extent.j && l.k < extent.k) return PETSC_TRUE; 553e72e933SJed Brown return PETSC_FALSE; 563e72e933SJed Brown } 573e72e933SJed Brown 583e72e933SJed Brown // Since element/vertex box extents are typically not equal powers of 2, Z codes that lie within the domain are not contiguous. 593e72e933SJed Brown static ZLayout ZLayoutCreate(PetscMPIInt size, const PetscInt eextent[3], const PetscInt vextent[3]) 603e72e933SJed Brown { 613e72e933SJed Brown ZLayout layout; 623e72e933SJed Brown layout.eextent.i = eextent[0]; 633e72e933SJed Brown layout.eextent.j = eextent[1]; 643e72e933SJed Brown layout.eextent.k = eextent[2]; 653e72e933SJed Brown layout.vextent.i = vextent[0]; 663e72e933SJed Brown layout.vextent.j = vextent[1]; 673e72e933SJed Brown layout.vextent.k = vextent[2]; 683e72e933SJed Brown layout.comm_size = size; 693e72e933SJed Brown PetscMalloc1(size + 1, &layout.zstarts); 703e72e933SJed Brown 713e72e933SJed Brown PetscInt total_elems = eextent[0] * eextent[1] * eextent[2]; 723e72e933SJed Brown ZCode z = 0; 733e72e933SJed Brown layout.zstarts[0] = 0; 743e72e933SJed Brown for (PetscMPIInt r = 0; r < size; r++) { 753e72e933SJed Brown PetscInt elems_needed = (total_elems / size) + (total_elems % size > r), count; 763e72e933SJed Brown for (count = 0; count < elems_needed; z++) { 773e72e933SJed Brown Ijk loc = ZCodeSplit(z); 783e72e933SJed Brown if (IjkActive(layout.eextent, loc)) count++; 793e72e933SJed Brown } 803e72e933SJed Brown // Pick up any extra vertices in the Z ordering before the next rank's first owned element. 813e72e933SJed Brown for (; z <= ZEncode(layout.vextent); z++) { 823e72e933SJed Brown Ijk loc = ZCodeSplit(z); 833e72e933SJed Brown if (IjkActive(layout.eextent, loc)) break; 843e72e933SJed Brown } 853e72e933SJed Brown layout.zstarts[r + 1] = z; 863e72e933SJed Brown } 873e72e933SJed Brown return layout; 883e72e933SJed Brown } 893e72e933SJed Brown 903e72e933SJed Brown PetscInt ZCodeFind(ZCode key, PetscInt n, const ZCode X[]) 913e72e933SJed Brown { 923e72e933SJed Brown PetscInt lo = 0, hi = n; 933e72e933SJed Brown 943e72e933SJed Brown if (n == 0) return -1; 953e72e933SJed Brown while (hi - lo > 1) { 963e72e933SJed Brown PetscInt mid = lo + (hi - lo) / 2; 973e72e933SJed Brown if (key < X[mid]) hi = mid; 983e72e933SJed Brown else lo = mid; 993e72e933SJed Brown } 1003e72e933SJed Brown return key == X[lo] ? lo : -(lo + (key > X[lo]) + 1); 1013e72e933SJed Brown } 1023e72e933SJed Brown 1033e72e933SJed Brown PetscErrorCode DMPlexCreateBoxMesh_Tensor_SFC_Internal(DM dm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate) 1043e72e933SJed Brown { 1053e72e933SJed Brown PetscInt eextent[3] = {1, 1, 1}, vextent[3] = {1, 1, 1}; 1063e72e933SJed Brown const Ijk closure_1[] = { 1073e72e933SJed Brown {0, 0, 0}, 1083e72e933SJed Brown {1, 0, 0}, 1093e72e933SJed Brown }; 1103e72e933SJed Brown const Ijk closure_2[] = { 1113e72e933SJed Brown {0, 0, 0}, 1123e72e933SJed Brown {1, 0, 0}, 1133e72e933SJed Brown {1, 1, 0}, 1143e72e933SJed Brown {0, 1, 0}, 1153e72e933SJed Brown }; 1163e72e933SJed Brown const Ijk closure_3[] = { 1173e72e933SJed Brown {0, 0, 0}, 1183e72e933SJed Brown {0, 1, 0}, 1193e72e933SJed Brown {1, 1, 0}, 1203e72e933SJed Brown {1, 0, 0}, 1213e72e933SJed Brown {0, 0, 1}, 1223e72e933SJed Brown {1, 0, 1}, 1233e72e933SJed Brown {1, 1, 1}, 1243e72e933SJed Brown {0, 1, 1}, 1253e72e933SJed Brown }; 126*3431e603SJed Brown const Ijk *const closure_dim[] = {NULL, closure_1, closure_2, closure_3}; 127*3431e603SJed Brown // This must be kept consistent with DMPlexCreateCubeMesh_Internal 128*3431e603SJed Brown const PetscInt face_marker_1[] = {1, 2}; 129*3431e603SJed Brown const PetscInt face_marker_2[] = {4, 2, 1, 3}; 130*3431e603SJed Brown const PetscInt face_marker_3[] = {6, 5, 3, 4, 1, 2}; 131*3431e603SJed Brown const PetscInt *const face_marker_dim[] = {NULL, face_marker_1, face_marker_2, face_marker_3}; 1323e72e933SJed Brown 1333e72e933SJed Brown PetscFunctionBegin; 1343e72e933SJed Brown PetscValidPointer(dm, 1); 1353e72e933SJed Brown PetscValidLogicalCollectiveInt(dm, dim, 2); 1363e72e933SJed Brown PetscCall(DMSetDimension(dm, dim)); 1373e72e933SJed Brown PetscMPIInt rank, size; 1383e72e933SJed Brown PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 1393e72e933SJed Brown PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 1403e72e933SJed Brown for (PetscInt i = 0; i < dim; i++) { 1413e72e933SJed Brown eextent[i] = faces[i]; 1423e72e933SJed Brown vextent[i] = faces[i] + 1; 1433e72e933SJed Brown } 1443e72e933SJed Brown ZLayout layout = ZLayoutCreate(size, eextent, vextent); 1453e72e933SJed Brown PetscZSet vset; // set of all vertices in the closure of the owned elements 1463e72e933SJed Brown PetscCall(PetscZSetCreate(&vset)); 1473e72e933SJed Brown PetscInt local_elems = 0; 1483e72e933SJed Brown for (ZCode z = layout.zstarts[rank]; z < layout.zstarts[rank + 1]; z++) { 1493e72e933SJed Brown Ijk loc = ZCodeSplit(z); 1503e72e933SJed Brown if (IjkActive(layout.vextent, loc)) PetscZSetAdd(vset, z); 1513e72e933SJed Brown if (IjkActive(layout.eextent, loc)) { 1523e72e933SJed Brown local_elems++; 1533e72e933SJed Brown // Add all neighboring vertices to set 1543e72e933SJed Brown for (PetscInt n = 0; n < PetscPowInt(2, dim); n++) { 1553e72e933SJed Brown Ijk inc = closure_dim[dim][n]; 1563e72e933SJed Brown Ijk nloc = {loc.i + inc.i, loc.j + inc.j, loc.k + inc.k}; 1573e72e933SJed Brown ZCode v = ZEncode(nloc); 1583e72e933SJed Brown PetscZSetAdd(vset, v); 1593e72e933SJed Brown } 1603e72e933SJed Brown } 1613e72e933SJed Brown } 1623e72e933SJed Brown PetscInt local_verts, off = 0; 1633e72e933SJed Brown ZCode *vert_z; 1643e72e933SJed Brown PetscCall(PetscZSetGetSize(vset, &local_verts)); 1653e72e933SJed Brown PetscCall(PetscMalloc1(local_verts, &vert_z)); 1663e72e933SJed Brown PetscCall(PetscZSetGetElems(vset, &off, vert_z)); 1673e72e933SJed Brown PetscCall(PetscZSetDestroy(&vset)); 1683e72e933SJed Brown // ZCode is unsigned for bitwise convenience, but highest bit should never be set, so can interpret as signed 1693e72e933SJed Brown PetscCall(PetscSortInt64(local_verts, (PetscInt64 *)vert_z)); 1703e72e933SJed Brown 1713e72e933SJed Brown PetscCall(DMPlexSetChart(dm, 0, local_elems + local_verts)); 1723e72e933SJed Brown for (PetscInt e = 0; e < local_elems; e++) PetscCall(DMPlexSetConeSize(dm, e, PetscPowInt(2, dim))); 1733e72e933SJed Brown PetscCall(DMSetUp(dm)); 1743e72e933SJed Brown { 1753e72e933SJed Brown PetscInt e = 0; 1763e72e933SJed Brown for (ZCode z = layout.zstarts[rank]; z < layout.zstarts[rank + 1]; z++) { 1773e72e933SJed Brown Ijk loc = ZCodeSplit(z); 1783e72e933SJed Brown if (!IjkActive(layout.eextent, loc)) continue; 1793e72e933SJed Brown PetscInt cone[8], orient[8] = {0}; 1803e72e933SJed Brown for (PetscInt n = 0; n < PetscPowInt(2, dim); n++) { 1813e72e933SJed Brown Ijk inc = closure_dim[dim][n]; 1823e72e933SJed Brown Ijk nloc = {loc.i + inc.i, loc.j + inc.j, loc.k + inc.k}; 1833e72e933SJed Brown ZCode v = ZEncode(nloc); 1843e72e933SJed Brown PetscInt ci = ZCodeFind(v, local_verts, vert_z); 1853e72e933SJed Brown PetscAssert(ci >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Did not find neighbor vertex in set"); 1863e72e933SJed Brown cone[n] = local_elems + ci; 1873e72e933SJed Brown } 1883e72e933SJed Brown PetscCall(DMPlexSetCone(dm, e, cone)); 1893e72e933SJed Brown PetscCall(DMPlexSetConeOrientation(dm, e, orient)); 1903e72e933SJed Brown e++; 1913e72e933SJed Brown } 1923e72e933SJed Brown } 1933e72e933SJed Brown 1943e72e933SJed Brown if (0) { 1953e72e933SJed Brown DMLabel depth; 1963e72e933SJed Brown PetscCall(DMCreateLabel(dm, "depth")); 1973e72e933SJed Brown PetscCall(DMPlexGetDepthLabel(dm, &depth)); 1983e72e933SJed Brown PetscCall(DMLabelSetStratumBounds(depth, 0, local_elems, local_elems + local_verts)); 1993e72e933SJed Brown PetscCall(DMLabelSetStratumBounds(depth, 1, 0, local_elems)); 2003e72e933SJed Brown } else { 2013e72e933SJed Brown PetscCall(DMPlexSymmetrize(dm)); 2023e72e933SJed Brown PetscCall(DMPlexStratify(dm)); 2033e72e933SJed Brown } 2043e72e933SJed Brown { // Create point SF 2053e72e933SJed Brown PetscSF sf; 2063e72e933SJed Brown PetscSFCreate(PetscObjectComm((PetscObject)dm), &sf); 2073e72e933SJed Brown PetscInt owned_verts = ZCodeFind(layout.zstarts[rank + 1], local_verts, vert_z); 2083e72e933SJed Brown if (owned_verts < 0) owned_verts = -(owned_verts + 1); // We don't care whether the key was found 2093e72e933SJed Brown PetscInt num_ghosts = local_verts - owned_verts; // Due to sorting, owned vertices always come first 2103e72e933SJed Brown PetscInt *local_ghosts; 2113e72e933SJed Brown PetscSFNode *ghosts; 2123e72e933SJed Brown PetscCall(PetscMalloc1(num_ghosts, &local_ghosts)); 2133e72e933SJed Brown PetscCall(PetscMalloc1(num_ghosts, &ghosts)); 2143e72e933SJed Brown for (PetscInt i = 0; i < num_ghosts;) { 2153e72e933SJed Brown ZCode z = vert_z[owned_verts + i]; 2163e72e933SJed Brown PetscInt remote_rank = ZCodeFind(z, size + 1, layout.zstarts), remote_count = 0; 2173e72e933SJed Brown if (remote_rank < 0) remote_rank = -(remote_rank + 1) - 1; 2183e72e933SJed Brown // We have a new remote rank; find all the ghost indices (which are contiguous in vert_z) 2193e72e933SJed Brown 2203e72e933SJed Brown // Count the elements on remote_rank 2213e72e933SJed Brown PetscInt remote_elem = 0; 2223e72e933SJed Brown for (ZCode rz = layout.zstarts[remote_rank]; rz < layout.zstarts[remote_rank + 1]; rz++) { 2233e72e933SJed Brown Ijk loc = ZCodeSplit(rz); 2243e72e933SJed Brown if (IjkActive(layout.eextent, loc)) remote_elem++; 2253e72e933SJed Brown } 2263e72e933SJed Brown 2273e72e933SJed Brown // Traverse vertices and make ghost links 2283e72e933SJed Brown for (ZCode rz = layout.zstarts[remote_rank]; rz < layout.zstarts[remote_rank + 1]; rz++) { 2293e72e933SJed Brown Ijk loc = ZCodeSplit(rz); 2303e72e933SJed Brown if (rz == z) { 2313e72e933SJed Brown local_ghosts[i] = local_elems + owned_verts + i; 2323e72e933SJed Brown ghosts[i].rank = remote_rank; 2333e72e933SJed Brown ghosts[i].index = remote_elem + remote_count; 2343e72e933SJed Brown i++; 2353e72e933SJed Brown if (i == num_ghosts) break; 2363e72e933SJed Brown z = vert_z[owned_verts + i]; 2373e72e933SJed Brown } 2383e72e933SJed Brown if (IjkActive(layout.vextent, loc)) remote_count++; 2393e72e933SJed Brown } 2403e72e933SJed Brown } 2413e72e933SJed Brown PetscCall(PetscSFSetGraph(sf, local_elems + local_verts, num_ghosts, local_ghosts, PETSC_OWN_POINTER, ghosts, PETSC_OWN_POINTER)); 2423e72e933SJed Brown PetscCall(PetscObjectSetName((PetscObject)sf, "SFC Point SF")); 2433e72e933SJed Brown PetscCall(DMSetPointSF(dm, sf)); 2443e72e933SJed Brown PetscCall(PetscSFDestroy(&sf)); 2453e72e933SJed Brown } 2463e72e933SJed Brown { 2473e72e933SJed Brown Vec coordinates; 2483e72e933SJed Brown PetscScalar *coords; 2493e72e933SJed Brown PetscSection coord_section; 2503e72e933SJed Brown PetscInt coord_size; 2513e72e933SJed Brown PetscCall(DMGetCoordinateSection(dm, &coord_section)); 2523e72e933SJed Brown PetscCall(PetscSectionSetNumFields(coord_section, 1)); 2533e72e933SJed Brown PetscCall(PetscSectionSetFieldComponents(coord_section, 0, dim)); 2543e72e933SJed Brown PetscCall(PetscSectionSetChart(coord_section, local_elems, local_elems + local_verts)); 2553e72e933SJed Brown for (PetscInt v = 0; v < local_verts; v++) { 2563e72e933SJed Brown PetscInt point = local_elems + v; 2573e72e933SJed Brown PetscCall(PetscSectionSetDof(coord_section, point, dim)); 2583e72e933SJed Brown PetscCall(PetscSectionSetFieldDof(coord_section, point, 0, dim)); 2593e72e933SJed Brown } 2603e72e933SJed Brown PetscCall(PetscSectionSetUp(coord_section)); 2613e72e933SJed Brown PetscCall(PetscSectionGetStorageSize(coord_section, &coord_size)); 2623e72e933SJed Brown PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 2633e72e933SJed Brown PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 2643e72e933SJed Brown PetscCall(VecSetSizes(coordinates, coord_size, PETSC_DETERMINE)); 2653e72e933SJed Brown PetscCall(VecSetBlockSize(coordinates, dim)); 2663e72e933SJed Brown PetscCall(VecSetType(coordinates, VECSTANDARD)); 2673e72e933SJed Brown PetscCall(VecGetArray(coordinates, &coords)); 2683e72e933SJed Brown for (PetscInt v = 0; v < local_verts; v++) { 2693e72e933SJed Brown Ijk loc = ZCodeSplit(vert_z[v]); 2703e72e933SJed Brown coords[v * dim + 0] = lower[0] + loc.i * (upper[0] - lower[0]) / layout.eextent.i; 2713e72e933SJed Brown if (dim > 1) coords[v * dim + 1] = lower[1] + loc.j * (upper[1] - lower[1]) / layout.eextent.j; 2723e72e933SJed Brown if (dim > 2) coords[v * dim + 2] = lower[2] + loc.k * (upper[2] - lower[2]) / layout.eextent.k; 2733e72e933SJed Brown } 2743e72e933SJed Brown PetscCall(VecRestoreArray(coordinates, &coords)); 2753e72e933SJed Brown PetscCall(DMSetCoordinatesLocal(dm, coordinates)); 2763e72e933SJed Brown PetscCall(VecDestroy(&coordinates)); 2773e72e933SJed Brown PetscCall(PetscFree(layout.zstarts)); 2783e72e933SJed Brown } 2793e72e933SJed Brown if (interpolate) { 280*3431e603SJed Brown PetscCall(DMPlexInterpolateInPlace_Internal(dm)); 281*3431e603SJed Brown 282*3431e603SJed Brown DMLabel label; 283*3431e603SJed Brown PetscCall(DMCreateLabel(dm, "Face Sets")); 284*3431e603SJed Brown PetscCall(DMGetLabel(dm, "Face Sets", &label)); 285*3431e603SJed Brown PetscInt fStart, fEnd, vStart, vEnd; 286*3431e603SJed Brown PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd)); 287*3431e603SJed Brown PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 288*3431e603SJed Brown for (PetscInt f = fStart; f < fEnd; f++) { 289*3431e603SJed Brown PetscInt npoints; 290*3431e603SJed Brown PetscInt *points = NULL; 291*3431e603SJed Brown PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &npoints, &points)); 292*3431e603SJed Brown PetscInt bc_count[6] = {0}; 293*3431e603SJed Brown for (PetscInt i = 0; i < npoints; i++) { 294*3431e603SJed Brown PetscInt p = points[2 * i]; 295*3431e603SJed Brown if (p < vStart || vEnd <= p) continue; 296*3431e603SJed Brown Ijk loc = ZCodeSplit(vert_z[p - vStart]); 297*3431e603SJed Brown // Convention here matches DMPlexCreateCubeMesh_Internal 298*3431e603SJed Brown bc_count[0] += loc.i == 0; 299*3431e603SJed Brown bc_count[1] += loc.i == layout.vextent.i - 1; 300*3431e603SJed Brown bc_count[2] += loc.j == 0; 301*3431e603SJed Brown bc_count[3] += loc.j == layout.vextent.j - 1; 302*3431e603SJed Brown bc_count[4] += loc.k == 0; 303*3431e603SJed Brown bc_count[5] += loc.k == layout.vextent.k - 1; 3043e72e933SJed Brown } 305*3431e603SJed Brown PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &npoints, &points)); 306*3431e603SJed Brown for (PetscInt bc = 0, bc_match = 0; bc < 2 * dim; bc++) { 307*3431e603SJed Brown if (bc_count[bc] == PetscPowInt(2, dim - 1)) { 308*3431e603SJed Brown PetscAssert(bc_match == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face matches multiple face sets"); 309*3431e603SJed Brown PetscCall(DMLabelSetValue(label, f, face_marker_dim[dim][bc])); 310*3431e603SJed Brown bc_match++; 311*3431e603SJed Brown } 312*3431e603SJed Brown } 313*3431e603SJed Brown } 314*3431e603SJed Brown // Ensure that the Coordinate DM has our new boundary labels 315*3431e603SJed Brown DM cdm; 316*3431e603SJed Brown PetscCall(DMGetCoordinateDM(dm, &cdm)); 317*3431e603SJed Brown PetscCall(DMCopyLabels(dm, cdm, PETSC_COPY_VALUES, PETSC_FALSE, DM_COPY_LABELS_FAIL)); 318*3431e603SJed Brown } 319*3431e603SJed Brown PetscCall(PetscFree(vert_z)); 3203e72e933SJed Brown PetscFunctionReturn(0); 3213e72e933SJed Brown } 322