1*3e72e933SJed Brown #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2*3e72e933SJed Brown #include <petscsf.h> 3*3e72e933SJed Brown #include <petsc/private/hashset.h> 4*3e72e933SJed Brown 5*3e72e933SJed Brown typedef uint64_t ZCode; 6*3e72e933SJed Brown 7*3e72e933SJed Brown PETSC_HASH_SET(ZSet, ZCode, PetscHash_UInt64, PetscHashEqual) 8*3e72e933SJed Brown 9*3e72e933SJed Brown typedef struct { 10*3e72e933SJed Brown PetscInt i, j, k; 11*3e72e933SJed Brown } Ijk; 12*3e72e933SJed Brown 13*3e72e933SJed Brown typedef struct { 14*3e72e933SJed Brown Ijk eextent; 15*3e72e933SJed Brown Ijk vextent; 16*3e72e933SJed Brown PetscMPIInt comm_size; 17*3e72e933SJed Brown ZCode *zstarts; 18*3e72e933SJed Brown } ZLayout; 19*3e72e933SJed Brown 20*3e72e933SJed Brown static unsigned ZCodeSplit1(ZCode z) 21*3e72e933SJed Brown { 22*3e72e933SJed Brown z = ((z & 01001001001001001) | ((z >> 2) & 02002002002002002) | ((z >> 4) & 04004004004004004)); 23*3e72e933SJed Brown z = (z | (z >> 6) | (z >> 12)) & 0000000777000000777; 24*3e72e933SJed Brown z = (z | (z >> 18)) & 0777777; 25*3e72e933SJed Brown return (unsigned)z; 26*3e72e933SJed Brown } 27*3e72e933SJed Brown 28*3e72e933SJed Brown static ZCode ZEncode1(unsigned t) 29*3e72e933SJed Brown { 30*3e72e933SJed Brown ZCode z = t; 31*3e72e933SJed Brown z = (z | (z << 18)) & 0777000000777; 32*3e72e933SJed Brown z = (z | (z << 6) | (z << 12)) & 07007007007007007; 33*3e72e933SJed Brown z = (z | (z << 2) | (z << 4)) & 0111111111111111111; 34*3e72e933SJed Brown return z; 35*3e72e933SJed Brown } 36*3e72e933SJed Brown 37*3e72e933SJed Brown static Ijk ZCodeSplit(ZCode z) 38*3e72e933SJed Brown { 39*3e72e933SJed Brown Ijk c; 40*3e72e933SJed Brown c.i = ZCodeSplit1(z >> 2); 41*3e72e933SJed Brown c.j = ZCodeSplit1(z >> 1); 42*3e72e933SJed Brown c.k = ZCodeSplit1(z >> 0); 43*3e72e933SJed Brown return c; 44*3e72e933SJed Brown } 45*3e72e933SJed Brown 46*3e72e933SJed Brown static ZCode ZEncode(Ijk c) 47*3e72e933SJed Brown { 48*3e72e933SJed Brown ZCode z = (ZEncode1(c.i) << 2) | (ZEncode1(c.j) << 1) | ZEncode1(c.k); 49*3e72e933SJed Brown return z; 50*3e72e933SJed Brown } 51*3e72e933SJed Brown 52*3e72e933SJed Brown static PetscBool IjkActive(Ijk extent, Ijk l) 53*3e72e933SJed Brown { 54*3e72e933SJed Brown if (l.i < extent.i && l.j < extent.j && l.k < extent.k) return PETSC_TRUE; 55*3e72e933SJed Brown return PETSC_FALSE; 56*3e72e933SJed Brown } 57*3e72e933SJed Brown 58*3e72e933SJed Brown // Since element/vertex box extents are typically not equal powers of 2, Z codes that lie within the domain are not contiguous. 59*3e72e933SJed Brown static ZLayout ZLayoutCreate(PetscMPIInt size, const PetscInt eextent[3], const PetscInt vextent[3]) 60*3e72e933SJed Brown { 61*3e72e933SJed Brown ZLayout layout; 62*3e72e933SJed Brown layout.eextent.i = eextent[0]; 63*3e72e933SJed Brown layout.eextent.j = eextent[1]; 64*3e72e933SJed Brown layout.eextent.k = eextent[2]; 65*3e72e933SJed Brown layout.vextent.i = vextent[0]; 66*3e72e933SJed Brown layout.vextent.j = vextent[1]; 67*3e72e933SJed Brown layout.vextent.k = vextent[2]; 68*3e72e933SJed Brown layout.comm_size = size; 69*3e72e933SJed Brown PetscMalloc1(size + 1, &layout.zstarts); 70*3e72e933SJed Brown 71*3e72e933SJed Brown PetscInt total_elems = eextent[0] * eextent[1] * eextent[2]; 72*3e72e933SJed Brown ZCode z = 0; 73*3e72e933SJed Brown layout.zstarts[0] = 0; 74*3e72e933SJed Brown for (PetscMPIInt r = 0; r < size; r++) { 75*3e72e933SJed Brown PetscInt elems_needed = (total_elems / size) + (total_elems % size > r), count; 76*3e72e933SJed Brown for (count = 0; count < elems_needed; z++) { 77*3e72e933SJed Brown Ijk loc = ZCodeSplit(z); 78*3e72e933SJed Brown if (IjkActive(layout.eextent, loc)) count++; 79*3e72e933SJed Brown } 80*3e72e933SJed Brown // Pick up any extra vertices in the Z ordering before the next rank's first owned element. 81*3e72e933SJed Brown for (; z <= ZEncode(layout.vextent); z++) { 82*3e72e933SJed Brown Ijk loc = ZCodeSplit(z); 83*3e72e933SJed Brown if (IjkActive(layout.eextent, loc)) break; 84*3e72e933SJed Brown } 85*3e72e933SJed Brown layout.zstarts[r + 1] = z; 86*3e72e933SJed Brown } 87*3e72e933SJed Brown return layout; 88*3e72e933SJed Brown } 89*3e72e933SJed Brown 90*3e72e933SJed Brown PetscInt ZCodeFind(ZCode key, PetscInt n, const ZCode X[]) 91*3e72e933SJed Brown { 92*3e72e933SJed Brown PetscInt lo = 0, hi = n; 93*3e72e933SJed Brown 94*3e72e933SJed Brown if (n == 0) return -1; 95*3e72e933SJed Brown while (hi - lo > 1) { 96*3e72e933SJed Brown PetscInt mid = lo + (hi - lo) / 2; 97*3e72e933SJed Brown if (key < X[mid]) hi = mid; 98*3e72e933SJed Brown else lo = mid; 99*3e72e933SJed Brown } 100*3e72e933SJed Brown return key == X[lo] ? lo : -(lo + (key > X[lo]) + 1); 101*3e72e933SJed Brown } 102*3e72e933SJed Brown 103*3e72e933SJed Brown PetscErrorCode DMPlexCreateBoxMesh_Tensor_SFC_Internal(DM dm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate) 104*3e72e933SJed Brown { 105*3e72e933SJed Brown PetscInt eextent[3] = {1, 1, 1}, vextent[3] = {1, 1, 1}; 106*3e72e933SJed Brown const Ijk closure_1[] = { 107*3e72e933SJed Brown {0, 0, 0}, 108*3e72e933SJed Brown {1, 0, 0}, 109*3e72e933SJed Brown }; 110*3e72e933SJed Brown const Ijk closure_2[] = { 111*3e72e933SJed Brown {0, 0, 0}, 112*3e72e933SJed Brown {1, 0, 0}, 113*3e72e933SJed Brown {1, 1, 0}, 114*3e72e933SJed Brown {0, 1, 0}, 115*3e72e933SJed Brown }; 116*3e72e933SJed Brown const Ijk closure_3[] = { 117*3e72e933SJed Brown {0, 0, 0}, 118*3e72e933SJed Brown {0, 1, 0}, 119*3e72e933SJed Brown {1, 1, 0}, 120*3e72e933SJed Brown {1, 0, 0}, 121*3e72e933SJed Brown {0, 0, 1}, 122*3e72e933SJed Brown {1, 0, 1}, 123*3e72e933SJed Brown {1, 1, 1}, 124*3e72e933SJed Brown {0, 1, 1}, 125*3e72e933SJed Brown }; 126*3e72e933SJed Brown const Ijk *closure_dim[] = {NULL, closure_1, closure_2, closure_3}; 127*3e72e933SJed Brown 128*3e72e933SJed Brown PetscFunctionBegin; 129*3e72e933SJed Brown PetscValidPointer(dm, 1); 130*3e72e933SJed Brown PetscValidLogicalCollectiveInt(dm, dim, 2); 131*3e72e933SJed Brown PetscCall(DMSetDimension(dm, dim)); 132*3e72e933SJed Brown PetscMPIInt rank, size; 133*3e72e933SJed Brown PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 134*3e72e933SJed Brown PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 135*3e72e933SJed Brown for (PetscInt i = 0; i < dim; i++) { 136*3e72e933SJed Brown eextent[i] = faces[i]; 137*3e72e933SJed Brown vextent[i] = faces[i] + 1; 138*3e72e933SJed Brown } 139*3e72e933SJed Brown ZLayout layout = ZLayoutCreate(size, eextent, vextent); 140*3e72e933SJed Brown PetscZSet vset; // set of all vertices in the closure of the owned elements 141*3e72e933SJed Brown PetscCall(PetscZSetCreate(&vset)); 142*3e72e933SJed Brown PetscInt local_elems = 0; 143*3e72e933SJed Brown for (ZCode z = layout.zstarts[rank]; z < layout.zstarts[rank + 1]; z++) { 144*3e72e933SJed Brown Ijk loc = ZCodeSplit(z); 145*3e72e933SJed Brown if (IjkActive(layout.vextent, loc)) PetscZSetAdd(vset, z); 146*3e72e933SJed Brown if (IjkActive(layout.eextent, loc)) { 147*3e72e933SJed Brown local_elems++; 148*3e72e933SJed Brown // Add all neighboring vertices to set 149*3e72e933SJed Brown for (PetscInt n = 0; n < PetscPowInt(2, dim); n++) { 150*3e72e933SJed Brown Ijk inc = closure_dim[dim][n]; 151*3e72e933SJed Brown Ijk nloc = {loc.i + inc.i, loc.j + inc.j, loc.k + inc.k}; 152*3e72e933SJed Brown ZCode v = ZEncode(nloc); 153*3e72e933SJed Brown PetscZSetAdd(vset, v); 154*3e72e933SJed Brown } 155*3e72e933SJed Brown } 156*3e72e933SJed Brown } 157*3e72e933SJed Brown PetscInt local_verts, off = 0; 158*3e72e933SJed Brown ZCode *vert_z; 159*3e72e933SJed Brown PetscCall(PetscZSetGetSize(vset, &local_verts)); 160*3e72e933SJed Brown PetscCall(PetscMalloc1(local_verts, &vert_z)); 161*3e72e933SJed Brown PetscCall(PetscZSetGetElems(vset, &off, vert_z)); 162*3e72e933SJed Brown PetscCall(PetscZSetDestroy(&vset)); 163*3e72e933SJed Brown // ZCode is unsigned for bitwise convenience, but highest bit should never be set, so can interpret as signed 164*3e72e933SJed Brown PetscCall(PetscSortInt64(local_verts, (PetscInt64 *)vert_z)); 165*3e72e933SJed Brown 166*3e72e933SJed Brown PetscCall(DMPlexSetChart(dm, 0, local_elems + local_verts)); 167*3e72e933SJed Brown for (PetscInt e = 0; e < local_elems; e++) PetscCall(DMPlexSetConeSize(dm, e, PetscPowInt(2, dim))); 168*3e72e933SJed Brown PetscCall(DMSetUp(dm)); 169*3e72e933SJed Brown { 170*3e72e933SJed Brown PetscInt e = 0; 171*3e72e933SJed Brown for (ZCode z = layout.zstarts[rank]; z < layout.zstarts[rank + 1]; z++) { 172*3e72e933SJed Brown Ijk loc = ZCodeSplit(z); 173*3e72e933SJed Brown if (!IjkActive(layout.eextent, loc)) continue; 174*3e72e933SJed Brown PetscInt cone[8], orient[8] = {0}; 175*3e72e933SJed Brown for (PetscInt n = 0; n < PetscPowInt(2, dim); n++) { 176*3e72e933SJed Brown Ijk inc = closure_dim[dim][n]; 177*3e72e933SJed Brown Ijk nloc = {loc.i + inc.i, loc.j + inc.j, loc.k + inc.k}; 178*3e72e933SJed Brown ZCode v = ZEncode(nloc); 179*3e72e933SJed Brown PetscInt ci = ZCodeFind(v, local_verts, vert_z); 180*3e72e933SJed Brown PetscAssert(ci >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Did not find neighbor vertex in set"); 181*3e72e933SJed Brown cone[n] = local_elems + ci; 182*3e72e933SJed Brown } 183*3e72e933SJed Brown PetscCall(DMPlexSetCone(dm, e, cone)); 184*3e72e933SJed Brown PetscCall(DMPlexSetConeOrientation(dm, e, orient)); 185*3e72e933SJed Brown e++; 186*3e72e933SJed Brown } 187*3e72e933SJed Brown } 188*3e72e933SJed Brown 189*3e72e933SJed Brown if (0) { 190*3e72e933SJed Brown DMLabel depth; 191*3e72e933SJed Brown PetscCall(DMCreateLabel(dm, "depth")); 192*3e72e933SJed Brown PetscCall(DMPlexGetDepthLabel(dm, &depth)); 193*3e72e933SJed Brown PetscCall(DMLabelSetStratumBounds(depth, 0, local_elems, local_elems + local_verts)); 194*3e72e933SJed Brown PetscCall(DMLabelSetStratumBounds(depth, 1, 0, local_elems)); 195*3e72e933SJed Brown } else { 196*3e72e933SJed Brown PetscCall(DMPlexSymmetrize(dm)); 197*3e72e933SJed Brown PetscCall(DMPlexStratify(dm)); 198*3e72e933SJed Brown } 199*3e72e933SJed Brown { // Create point SF 200*3e72e933SJed Brown PetscSF sf; 201*3e72e933SJed Brown PetscSFCreate(PetscObjectComm((PetscObject)dm), &sf); 202*3e72e933SJed Brown PetscInt owned_verts = ZCodeFind(layout.zstarts[rank + 1], local_verts, vert_z); 203*3e72e933SJed Brown if (owned_verts < 0) owned_verts = -(owned_verts + 1); // We don't care whether the key was found 204*3e72e933SJed Brown PetscInt num_ghosts = local_verts - owned_verts; // Due to sorting, owned vertices always come first 205*3e72e933SJed Brown PetscInt *local_ghosts; 206*3e72e933SJed Brown PetscSFNode *ghosts; 207*3e72e933SJed Brown PetscCall(PetscMalloc1(num_ghosts, &local_ghosts)); 208*3e72e933SJed Brown PetscCall(PetscMalloc1(num_ghosts, &ghosts)); 209*3e72e933SJed Brown for (PetscInt i = 0; i < num_ghosts;) { 210*3e72e933SJed Brown ZCode z = vert_z[owned_verts + i]; 211*3e72e933SJed Brown PetscInt remote_rank = ZCodeFind(z, size + 1, layout.zstarts), remote_count = 0; 212*3e72e933SJed Brown if (remote_rank < 0) remote_rank = -(remote_rank + 1) - 1; 213*3e72e933SJed Brown // We have a new remote rank; find all the ghost indices (which are contiguous in vert_z) 214*3e72e933SJed Brown 215*3e72e933SJed Brown // Count the elements on remote_rank 216*3e72e933SJed Brown PetscInt remote_elem = 0; 217*3e72e933SJed Brown for (ZCode rz = layout.zstarts[remote_rank]; rz < layout.zstarts[remote_rank + 1]; rz++) { 218*3e72e933SJed Brown Ijk loc = ZCodeSplit(rz); 219*3e72e933SJed Brown if (IjkActive(layout.eextent, loc)) remote_elem++; 220*3e72e933SJed Brown } 221*3e72e933SJed Brown 222*3e72e933SJed Brown // Traverse vertices and make ghost links 223*3e72e933SJed Brown for (ZCode rz = layout.zstarts[remote_rank]; rz < layout.zstarts[remote_rank + 1]; rz++) { 224*3e72e933SJed Brown Ijk loc = ZCodeSplit(rz); 225*3e72e933SJed Brown if (rz == z) { 226*3e72e933SJed Brown local_ghosts[i] = local_elems + owned_verts + i; 227*3e72e933SJed Brown ghosts[i].rank = remote_rank; 228*3e72e933SJed Brown ghosts[i].index = remote_elem + remote_count; 229*3e72e933SJed Brown i++; 230*3e72e933SJed Brown if (i == num_ghosts) break; 231*3e72e933SJed Brown z = vert_z[owned_verts + i]; 232*3e72e933SJed Brown } 233*3e72e933SJed Brown if (IjkActive(layout.vextent, loc)) remote_count++; 234*3e72e933SJed Brown } 235*3e72e933SJed Brown } 236*3e72e933SJed Brown PetscCall(PetscSFSetGraph(sf, local_elems + local_verts, num_ghosts, local_ghosts, PETSC_OWN_POINTER, ghosts, PETSC_OWN_POINTER)); 237*3e72e933SJed Brown PetscCall(PetscObjectSetName((PetscObject)sf, "SFC Point SF")); 238*3e72e933SJed Brown PetscCall(DMSetPointSF(dm, sf)); 239*3e72e933SJed Brown PetscCall(PetscSFDestroy(&sf)); 240*3e72e933SJed Brown } 241*3e72e933SJed Brown 242*3e72e933SJed Brown { 243*3e72e933SJed Brown Vec coordinates; 244*3e72e933SJed Brown PetscScalar *coords; 245*3e72e933SJed Brown PetscSection coord_section; 246*3e72e933SJed Brown PetscInt coord_size; 247*3e72e933SJed Brown DM cdm; 248*3e72e933SJed Brown PetscCall(DMGetCoordinateSection(dm, &coord_section)); 249*3e72e933SJed Brown PetscCall(PetscSectionSetNumFields(coord_section, 1)); 250*3e72e933SJed Brown PetscCall(PetscSectionSetFieldComponents(coord_section, 0, dim)); 251*3e72e933SJed Brown PetscCall(PetscSectionSetChart(coord_section, local_elems, local_elems + local_verts)); 252*3e72e933SJed Brown for (PetscInt v = 0; v < local_verts; v++) { 253*3e72e933SJed Brown PetscInt point = local_elems + v; 254*3e72e933SJed Brown PetscCall(PetscSectionSetDof(coord_section, point, dim)); 255*3e72e933SJed Brown PetscCall(PetscSectionSetFieldDof(coord_section, point, 0, dim)); 256*3e72e933SJed Brown } 257*3e72e933SJed Brown PetscCall(PetscSectionSetUp(coord_section)); 258*3e72e933SJed Brown PetscCall(PetscSectionGetStorageSize(coord_section, &coord_size)); 259*3e72e933SJed Brown PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 260*3e72e933SJed Brown PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 261*3e72e933SJed Brown PetscCall(VecSetSizes(coordinates, coord_size, PETSC_DETERMINE)); 262*3e72e933SJed Brown PetscCall(VecSetBlockSize(coordinates, dim)); 263*3e72e933SJed Brown PetscCall(VecSetType(coordinates, VECSTANDARD)); 264*3e72e933SJed Brown PetscCall(VecGetArray(coordinates, &coords)); 265*3e72e933SJed Brown for (PetscInt v = 0; v < local_verts; v++) { 266*3e72e933SJed Brown Ijk loc = ZCodeSplit(vert_z[v]); 267*3e72e933SJed Brown coords[v * dim + 0] = lower[0] + loc.i * (upper[0] - lower[0]) / layout.eextent.i; 268*3e72e933SJed Brown if (dim > 1) coords[v * dim + 1] = lower[1] + loc.j * (upper[1] - lower[1]) / layout.eextent.j; 269*3e72e933SJed Brown if (dim > 2) coords[v * dim + 2] = lower[2] + loc.k * (upper[2] - lower[2]) / layout.eextent.k; 270*3e72e933SJed Brown } 271*3e72e933SJed Brown PetscCall(VecRestoreArray(coordinates, &coords)); 272*3e72e933SJed Brown PetscCall(DMSetCoordinatesLocal(dm, coordinates)); 273*3e72e933SJed Brown PetscCall(VecDestroy(&coordinates)); 274*3e72e933SJed Brown PetscCall(PetscFree(vert_z)); 275*3e72e933SJed Brown PetscCall(PetscFree(layout.zstarts)); 276*3e72e933SJed Brown 277*3e72e933SJed Brown PetscFE fe; 278*3e72e933SJed Brown PetscCall(DMGetCoordinateDM(dm, &cdm)); 279*3e72e933SJed Brown PetscCall(PetscFECreateLagrange(PetscObjectComm((PetscObject)dm), dim, dim, PETSC_FALSE, 1, PETSC_DETERMINE, &fe)); 280*3e72e933SJed Brown PetscCall(DMSetField(cdm, 0, NULL, (PetscObject)fe)); 281*3e72e933SJed Brown PetscCall(PetscFEDestroy(&fe)); 282*3e72e933SJed Brown PetscCall(DMCreateDS(cdm)); 283*3e72e933SJed Brown } 284*3e72e933SJed Brown if (interpolate) { 285*3e72e933SJed Brown DM dmint; 286*3e72e933SJed Brown PetscCall(DMPlexInterpolate(dm, &dmint)); 287*3e72e933SJed Brown PetscCall(DMPlexReplace_Internal(dm, &dmint)); 288*3e72e933SJed Brown } 289*3e72e933SJed Brown PetscFunctionReturn(0); 290*3e72e933SJed Brown } 291