1 #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/ 2 #include <petsc/private/sectionimpl.h> 3 4 /*@ 5 PetscSFSetGraphLayout - Set a parallel star forest via global indices and a `PetscLayout` 6 7 Collective 8 9 Input Parameters: 10 + sf - star forest 11 . layout - `PetscLayout` defining the global space for roots 12 . nleaves - number of leaf vertices on the current process, each of these references a root on any process 13 . ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage 14 . localmode - copy mode for ilocal 15 - gremote - root vertices in global numbering corresponding to leaves in ilocal 16 17 Level: intermediate 18 19 Note: 20 Global indices must lie in [0, N) where N is the global size of layout. 21 Leaf indices in ilocal get sorted; this means the user-provided array gets sorted if localmode is `PETSC_OWN_POINTER`. 22 23 Developer Note: 24 Local indices which are the identity permutation in the range [0,nleaves) are discarded as they 25 encode contiguous storage. In such case, if localmode is `PETSC_OWN_POINTER`, the memory is deallocated as it is not 26 needed 27 28 .seealso: `PetscSF`, `PetscSFGetGraphLayout()`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()` 29 @*/ 30 PetscErrorCode PetscSFSetGraphLayout(PetscSF sf, PetscLayout layout, PetscInt nleaves, PetscInt *ilocal, PetscCopyMode localmode, const PetscInt *gremote) 31 { 32 const PetscInt *range; 33 PetscInt i, nroots, ls = -1, ln = -1; 34 PetscMPIInt lr = -1; 35 PetscSFNode *remote; 36 37 PetscFunctionBegin; 38 PetscCall(PetscLayoutSetUp(layout)); 39 PetscCall(PetscLayoutGetLocalSize(layout, &nroots)); 40 PetscCall(PetscLayoutGetRanges(layout, &range)); 41 PetscCall(PetscMalloc1(nleaves, &remote)); 42 if (nleaves) ls = gremote[0] + 1; 43 for (i = 0; i < nleaves; i++) { 44 const PetscInt idx = gremote[i] - ls; 45 if (idx < 0 || idx >= ln) { /* short-circuit the search */ 46 PetscCall(PetscLayoutFindOwnerIndex(layout, gremote[i], &lr, &remote[i].index)); 47 remote[i].rank = lr; 48 ls = range[lr]; 49 ln = range[lr + 1] - ls; 50 } else { 51 remote[i].rank = lr; 52 remote[i].index = idx; 53 } 54 } 55 PetscCall(PetscSFSetGraph(sf, nroots, nleaves, ilocal, localmode, remote, PETSC_OWN_POINTER)); 56 PetscFunctionReturn(0); 57 } 58 59 /*@C 60 PetscSFGetGraphLayout - Get the global indices and `PetscLayout` that describe this star forest 61 62 Collective 63 64 Input Parameter: 65 . sf - star forest 66 67 Output Parameters: 68 + layout - `PetscLayout` defining the global space for roots 69 . nleaves - number of leaf vertices on the current process, each of these references a root on any process 70 . ilocal - locations of leaves in leafdata buffers, or NULL for contiguous storage 71 - gremote - root vertices in global numbering corresponding to leaves in ilocal 72 73 Level: intermediate 74 75 Notes: 76 The outputs are such that passing them as inputs to `PetscSFSetGraphLayout()` would lead to the same star forest. 77 The outputs layout and gremote are freshly created each time this function is called, 78 so they need to be freed by user and cannot be qualified as const. 79 80 .seealso: `PetscSF`, `PetscSFSetGraphLayout()`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()` 81 @*/ 82 PetscErrorCode PetscSFGetGraphLayout(PetscSF sf, PetscLayout *layout, PetscInt *nleaves, const PetscInt *ilocal[], PetscInt *gremote[]) 83 { 84 PetscInt nr, nl; 85 const PetscSFNode *ir; 86 PetscLayout lt; 87 88 PetscFunctionBegin; 89 PetscCall(PetscSFGetGraph(sf, &nr, &nl, ilocal, &ir)); 90 PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)sf), nr, PETSC_DECIDE, 1, <)); 91 if (gremote) { 92 PetscInt i; 93 const PetscInt *range; 94 PetscInt *gr; 95 96 PetscCall(PetscLayoutGetRanges(lt, &range)); 97 PetscCall(PetscMalloc1(nl, &gr)); 98 for (i = 0; i < nl; i++) gr[i] = range[ir[i].rank] + ir[i].index; 99 *gremote = gr; 100 } 101 if (nleaves) *nleaves = nl; 102 if (layout) *layout = lt; 103 else PetscCall(PetscLayoutDestroy(<)); 104 PetscFunctionReturn(0); 105 } 106 107 /*@ 108 PetscSFSetGraphSection - Sets the `PetscSF` graph encoding the parallel dof overlap based upon the `PetscSection` describing the data layout. 109 110 Input Parameters: 111 + sf - The `PetscSF` 112 . localSection - `PetscSection` describing the local data layout 113 - globalSection - `PetscSection` describing the global data layout 114 115 Level: developer 116 117 .seealso: `PetscSF`, `PetscSFSetGraph()`, `PetscSFSetGraphLayout()` 118 @*/ 119 PetscErrorCode PetscSFSetGraphSection(PetscSF sf, PetscSection localSection, PetscSection globalSection) 120 { 121 MPI_Comm comm; 122 PetscLayout layout; 123 const PetscInt *ranges; 124 PetscInt *local; 125 PetscSFNode *remote; 126 PetscInt pStart, pEnd, p, nroots, nleaves = 0, l; 127 PetscMPIInt size, rank; 128 129 PetscFunctionBegin; 130 PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 131 PetscValidHeaderSpecific(localSection, PETSC_SECTION_CLASSID, 2); 132 PetscValidHeaderSpecific(globalSection, PETSC_SECTION_CLASSID, 3); 133 134 PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 135 PetscCallMPI(MPI_Comm_size(comm, &size)); 136 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 137 PetscCall(PetscSectionGetChart(globalSection, &pStart, &pEnd)); 138 PetscCall(PetscSectionGetConstrainedStorageSize(globalSection, &nroots)); 139 PetscCall(PetscLayoutCreate(comm, &layout)); 140 PetscCall(PetscLayoutSetBlockSize(layout, 1)); 141 PetscCall(PetscLayoutSetLocalSize(layout, nroots)); 142 PetscCall(PetscLayoutSetUp(layout)); 143 PetscCall(PetscLayoutGetRanges(layout, &ranges)); 144 for (p = pStart; p < pEnd; ++p) { 145 PetscInt gdof, gcdof; 146 147 PetscCall(PetscSectionGetDof(globalSection, p, &gdof)); 148 PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof)); 149 PetscCheck(gcdof <= (gdof < 0 ? -(gdof + 1) : gdof), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " has %" PetscInt_FMT " constraints > %" PetscInt_FMT " dof", p, gcdof, (gdof < 0 ? -(gdof + 1) : gdof)); 150 nleaves += gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof; 151 } 152 PetscCall(PetscMalloc1(nleaves, &local)); 153 PetscCall(PetscMalloc1(nleaves, &remote)); 154 for (p = pStart, l = 0; p < pEnd; ++p) { 155 const PetscInt *cind; 156 PetscInt dof, cdof, off, gdof, gcdof, goff, gsize, d, c; 157 158 PetscCall(PetscSectionGetDof(localSection, p, &dof)); 159 PetscCall(PetscSectionGetOffset(localSection, p, &off)); 160 PetscCall(PetscSectionGetConstraintDof(localSection, p, &cdof)); 161 PetscCall(PetscSectionGetConstraintIndices(localSection, p, &cind)); 162 PetscCall(PetscSectionGetDof(globalSection, p, &gdof)); 163 PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof)); 164 PetscCall(PetscSectionGetOffset(globalSection, p, &goff)); 165 if (!gdof) continue; /* Censored point */ 166 gsize = gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof; 167 if (gsize != dof - cdof) { 168 PetscCheck(gsize == dof, comm, PETSC_ERR_ARG_WRONG, "Global dof %" PetscInt_FMT " for point %" PetscInt_FMT " is neither the constrained size %" PetscInt_FMT ", nor the unconstrained %" PetscInt_FMT, gsize, p, dof - cdof, dof); 169 cdof = 0; /* Ignore constraints */ 170 } 171 for (d = 0, c = 0; d < dof; ++d) { 172 if ((c < cdof) && (cind[c] == d)) { 173 ++c; 174 continue; 175 } 176 local[l + d - c] = off + d; 177 } 178 PetscCheck(d - c == gsize, comm, PETSC_ERR_ARG_WRONG, "Point %" PetscInt_FMT ": Global dof %" PetscInt_FMT " != %" PetscInt_FMT " size - number of constraints", p, gsize, d - c); 179 if (gdof < 0) { 180 for (d = 0; d < gsize; ++d, ++l) { 181 PetscInt offset = -(goff + 1) + d, r; 182 183 PetscCall(PetscFindInt(offset, size + 1, ranges, &r)); 184 if (r < 0) r = -(r + 2); 185 PetscCheck(!(r < 0) && !(r >= size), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " mapped to invalid process %" PetscInt_FMT " (%" PetscInt_FMT ", %" PetscInt_FMT ")", p, r, gdof, goff); 186 remote[l].rank = r; 187 remote[l].index = offset - ranges[r]; 188 } 189 } else { 190 for (d = 0; d < gsize; ++d, ++l) { 191 remote[l].rank = rank; 192 remote[l].index = goff + d - ranges[rank]; 193 } 194 } 195 } 196 PetscCheck(l == nleaves, comm, PETSC_ERR_PLIB, "Iteration error, l %" PetscInt_FMT " != nleaves %" PetscInt_FMT, l, nleaves); 197 PetscCall(PetscLayoutDestroy(&layout)); 198 PetscCall(PetscSFSetGraph(sf, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER)); 199 PetscFunctionReturn(0); 200 } 201 202 /*@C 203 PetscSFDistributeSection - Create a new `PetscSection` reorganized, moving from the root to the leaves of the `PetscSF` 204 205 Collective 206 207 Input Parameters: 208 + sf - The `PetscSF` 209 - rootSection - Section defined on root space 210 211 Output Parameters: 212 + remoteOffsets - root offsets in leaf storage, or NULL 213 - leafSection - Section defined on the leaf space 214 215 Level: advanced 216 217 .seealso: `PetscSF`, `PetscSFCreate()` 218 @*/ 219 PetscErrorCode PetscSFDistributeSection(PetscSF sf, PetscSection rootSection, PetscInt **remoteOffsets, PetscSection leafSection) 220 { 221 PetscSF embedSF; 222 const PetscInt *indices; 223 IS selected; 224 PetscInt numFields, nroots, rpStart, rpEnd, lpStart = PETSC_MAX_INT, lpEnd = -1, f, c; 225 PetscBool *sub, hasc; 226 227 PetscFunctionBegin; 228 PetscCall(PetscLogEventBegin(PETSCSF_DistSect, sf, 0, 0, 0)); 229 PetscCall(PetscSectionGetNumFields(rootSection, &numFields)); 230 if (numFields) { 231 IS perm; 232 233 /* PetscSectionSetNumFields() calls PetscSectionReset(), which destroys 234 leafSection->perm. To keep this permutation set by the user, we grab 235 the reference before calling PetscSectionSetNumFields() and set it 236 back after. */ 237 PetscCall(PetscSectionGetPermutation(leafSection, &perm)); 238 PetscCall(PetscObjectReference((PetscObject)perm)); 239 PetscCall(PetscSectionSetNumFields(leafSection, numFields)); 240 PetscCall(PetscSectionSetPermutation(leafSection, perm)); 241 PetscCall(ISDestroy(&perm)); 242 } 243 PetscCall(PetscMalloc1(numFields + 2, &sub)); 244 sub[1] = rootSection->bc ? PETSC_TRUE : PETSC_FALSE; 245 for (f = 0; f < numFields; ++f) { 246 PetscSectionSym sym, dsym = NULL; 247 const char *name = NULL; 248 PetscInt numComp = 0; 249 250 sub[2 + f] = rootSection->field[f]->bc ? PETSC_TRUE : PETSC_FALSE; 251 PetscCall(PetscSectionGetFieldComponents(rootSection, f, &numComp)); 252 PetscCall(PetscSectionGetFieldName(rootSection, f, &name)); 253 PetscCall(PetscSectionGetFieldSym(rootSection, f, &sym)); 254 if (sym) PetscCall(PetscSectionSymDistribute(sym, sf, &dsym)); 255 PetscCall(PetscSectionSetFieldComponents(leafSection, f, numComp)); 256 PetscCall(PetscSectionSetFieldName(leafSection, f, name)); 257 PetscCall(PetscSectionSetFieldSym(leafSection, f, dsym)); 258 PetscCall(PetscSectionSymDestroy(&dsym)); 259 for (c = 0; c < rootSection->numFieldComponents[f]; ++c) { 260 PetscCall(PetscSectionGetComponentName(rootSection, f, c, &name)); 261 PetscCall(PetscSectionSetComponentName(leafSection, f, c, name)); 262 } 263 } 264 PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd)); 265 PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL)); 266 rpEnd = PetscMin(rpEnd, nroots); 267 rpEnd = PetscMax(rpStart, rpEnd); 268 /* see if we can avoid creating the embedded SF, since it can cost more than an allreduce */ 269 sub[0] = (PetscBool)(nroots != rpEnd - rpStart); 270 PetscCall(MPIU_Allreduce(MPI_IN_PLACE, sub, 2 + numFields, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)sf))); 271 if (sub[0]) { 272 PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected)); 273 PetscCall(ISGetIndices(selected, &indices)); 274 PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF)); 275 PetscCall(ISRestoreIndices(selected, &indices)); 276 PetscCall(ISDestroy(&selected)); 277 } else { 278 PetscCall(PetscObjectReference((PetscObject)sf)); 279 embedSF = sf; 280 } 281 PetscCall(PetscSFGetLeafRange(embedSF, &lpStart, &lpEnd)); 282 lpEnd++; 283 284 PetscCall(PetscSectionSetChart(leafSection, lpStart, lpEnd)); 285 286 /* Constrained dof section */ 287 hasc = sub[1]; 288 for (f = 0; f < numFields; ++f) hasc = (PetscBool)(hasc || sub[2 + f]); 289 290 /* Could fuse these at the cost of copies and extra allocation */ 291 PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart], MPI_REPLACE)); 292 PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart], MPI_REPLACE)); 293 if (sub[1]) { 294 PetscCall(PetscSectionCheckConstraints_Private(rootSection)); 295 PetscCall(PetscSectionCheckConstraints_Private(leafSection)); 296 PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE)); 297 PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE)); 298 } 299 for (f = 0; f < numFields; ++f) { 300 PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart], MPI_REPLACE)); 301 PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart], MPI_REPLACE)); 302 if (sub[2 + f]) { 303 PetscCall(PetscSectionCheckConstraints_Private(rootSection->field[f])); 304 PetscCall(PetscSectionCheckConstraints_Private(leafSection->field[f])); 305 PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE)); 306 PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE)); 307 } 308 } 309 if (remoteOffsets) { 310 PetscCall(PetscMalloc1(lpEnd - lpStart, remoteOffsets)); 311 PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE)); 312 PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE)); 313 } 314 PetscCall(PetscSectionInvalidateMaxDof_Internal(leafSection)); 315 PetscCall(PetscSectionSetUp(leafSection)); 316 if (hasc) { /* need to communicate bcIndices */ 317 PetscSF bcSF; 318 PetscInt *rOffBc; 319 320 PetscCall(PetscMalloc1(lpEnd - lpStart, &rOffBc)); 321 if (sub[1]) { 322 PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 323 PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 324 PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->bc, rOffBc, leafSection->bc, &bcSF)); 325 PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE)); 326 PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE)); 327 PetscCall(PetscSFDestroy(&bcSF)); 328 } 329 for (f = 0; f < numFields; ++f) { 330 if (sub[2 + f]) { 331 PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 332 PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 333 PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->field[f]->bc, rOffBc, leafSection->field[f]->bc, &bcSF)); 334 PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE)); 335 PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE)); 336 PetscCall(PetscSFDestroy(&bcSF)); 337 } 338 } 339 PetscCall(PetscFree(rOffBc)); 340 } 341 PetscCall(PetscSFDestroy(&embedSF)); 342 PetscCall(PetscFree(sub)); 343 PetscCall(PetscLogEventEnd(PETSCSF_DistSect, sf, 0, 0, 0)); 344 PetscFunctionReturn(0); 345 } 346 347 /*@C 348 PetscSFCreateRemoteOffsets - Create offsets for point data on remote processes 349 350 Collective 351 352 Input Parameters: 353 + sf - The `PetscSF` 354 . rootSection - Data layout of remote points for outgoing data (this is layout for SF roots) 355 - leafSection - Data layout of local points for incoming data (this is layout for SF leaves) 356 357 Output Parameter: 358 . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL 359 360 Level: developer 361 362 .seealso: `PetscSF`, `PetscSFCreate()` 363 @*/ 364 PetscErrorCode PetscSFCreateRemoteOffsets(PetscSF sf, PetscSection rootSection, PetscSection leafSection, PetscInt **remoteOffsets) 365 { 366 PetscSF embedSF; 367 const PetscInt *indices; 368 IS selected; 369 PetscInt numRoots, rpStart = 0, rpEnd = 0, lpStart = 0, lpEnd = 0; 370 371 PetscFunctionBegin; 372 *remoteOffsets = NULL; 373 PetscCall(PetscSFGetGraph(sf, &numRoots, NULL, NULL, NULL)); 374 if (numRoots < 0) PetscFunctionReturn(0); 375 PetscCall(PetscLogEventBegin(PETSCSF_RemoteOff, sf, 0, 0, 0)); 376 PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd)); 377 PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd)); 378 PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected)); 379 PetscCall(ISGetIndices(selected, &indices)); 380 PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF)); 381 PetscCall(ISRestoreIndices(selected, &indices)); 382 PetscCall(ISDestroy(&selected)); 383 PetscCall(PetscCalloc1(lpEnd - lpStart, remoteOffsets)); 384 PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE)); 385 PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE)); 386 PetscCall(PetscSFDestroy(&embedSF)); 387 PetscCall(PetscLogEventEnd(PETSCSF_RemoteOff, sf, 0, 0, 0)); 388 PetscFunctionReturn(0); 389 } 390 391 /*@C 392 PetscSFCreateSectionSF - Create an expanded `PetscSF` of dofs, assuming the input `PetscSF` relates points 393 394 Collective 395 396 Input Parameters: 397 + sf - The `PetscSF` 398 . rootSection - Data layout of remote points for outgoing data (this is usually the serial section) 399 . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL 400 - leafSection - Data layout of local points for incoming data (this is the distributed section) 401 402 Output Parameters: 403 - sectionSF - The new `PetscSF` 404 405 Level: advanced 406 407 Note: 408 Either rootSection or remoteOffsets can be specified 409 410 .seealso: `PetscSF`, `PetscSFCreate()` 411 @*/ 412 PetscErrorCode PetscSFCreateSectionSF(PetscSF sf, PetscSection rootSection, PetscInt remoteOffsets[], PetscSection leafSection, PetscSF *sectionSF) 413 { 414 MPI_Comm comm; 415 const PetscInt *localPoints; 416 const PetscSFNode *remotePoints; 417 PetscInt lpStart, lpEnd; 418 PetscInt numRoots, numSectionRoots, numPoints, numIndices = 0; 419 PetscInt *localIndices; 420 PetscSFNode *remoteIndices; 421 PetscInt i, ind; 422 423 PetscFunctionBegin; 424 PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 425 PetscValidPointer(rootSection, 2); 426 /* Cannot check PetscValidIntPointer(remoteOffsets,3) because it can be NULL if sf does not reference any points in leafSection */ 427 PetscValidPointer(leafSection, 4); 428 PetscValidPointer(sectionSF, 5); 429 PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 430 PetscCall(PetscSFCreate(comm, sectionSF)); 431 PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd)); 432 PetscCall(PetscSectionGetStorageSize(rootSection, &numSectionRoots)); 433 PetscCall(PetscSFGetGraph(sf, &numRoots, &numPoints, &localPoints, &remotePoints)); 434 if (numRoots < 0) PetscFunctionReturn(0); 435 PetscCall(PetscLogEventBegin(PETSCSF_SectSF, sf, 0, 0, 0)); 436 for (i = 0; i < numPoints; ++i) { 437 PetscInt localPoint = localPoints ? localPoints[i] : i; 438 PetscInt dof; 439 440 if ((localPoint >= lpStart) && (localPoint < lpEnd)) { 441 PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof)); 442 numIndices += dof < 0 ? 0 : dof; 443 } 444 } 445 PetscCall(PetscMalloc1(numIndices, &localIndices)); 446 PetscCall(PetscMalloc1(numIndices, &remoteIndices)); 447 /* Create new index graph */ 448 for (i = 0, ind = 0; i < numPoints; ++i) { 449 PetscInt localPoint = localPoints ? localPoints[i] : i; 450 PetscInt rank = remotePoints[i].rank; 451 452 if ((localPoint >= lpStart) && (localPoint < lpEnd)) { 453 PetscInt remoteOffset = remoteOffsets[localPoint - lpStart]; 454 PetscInt loff, dof, d; 455 456 PetscCall(PetscSectionGetOffset(leafSection, localPoint, &loff)); 457 PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof)); 458 for (d = 0; d < dof; ++d, ++ind) { 459 localIndices[ind] = loff + d; 460 remoteIndices[ind].rank = rank; 461 remoteIndices[ind].index = remoteOffset + d; 462 } 463 } 464 } 465 PetscCheck(numIndices == ind, comm, PETSC_ERR_PLIB, "Inconsistency in indices, %" PetscInt_FMT " should be %" PetscInt_FMT, ind, numIndices); 466 PetscCall(PetscSFSetGraph(*sectionSF, numSectionRoots, numIndices, localIndices, PETSC_OWN_POINTER, remoteIndices, PETSC_OWN_POINTER)); 467 PetscCall(PetscSFSetUp(*sectionSF)); 468 PetscCall(PetscLogEventEnd(PETSCSF_SectSF, sf, 0, 0, 0)); 469 PetscFunctionReturn(0); 470 } 471 472 /*@C 473 PetscSFCreateFromLayouts - Creates a parallel star forest mapping two `PetscLayout` objects 474 475 Collective 476 477 Input Parameters: 478 + rmap - `PetscLayout` defining the global root space 479 - lmap - `PetscLayout` defining the global leaf space 480 481 Output Parameter: 482 . sf - The parallel star forest 483 484 Level: intermediate 485 486 .seealso: `PetscSF`, `PetscSFCreate()`, `PetscLayoutCreate()`, `PetscSFSetGraphLayout()` 487 @*/ 488 PetscErrorCode PetscSFCreateFromLayouts(PetscLayout rmap, PetscLayout lmap, PetscSF *sf) 489 { 490 PetscInt i, nroots, nleaves = 0; 491 PetscInt rN, lst, len; 492 PetscMPIInt owner = -1; 493 PetscSFNode *remote; 494 MPI_Comm rcomm = rmap->comm; 495 MPI_Comm lcomm = lmap->comm; 496 PetscMPIInt flg; 497 498 PetscFunctionBegin; 499 PetscValidPointer(sf, 3); 500 PetscCheck(rmap->setupcalled, rcomm, PETSC_ERR_ARG_WRONGSTATE, "Root layout not setup"); 501 PetscCheck(lmap->setupcalled, lcomm, PETSC_ERR_ARG_WRONGSTATE, "Leaf layout not setup"); 502 PetscCallMPI(MPI_Comm_compare(rcomm, lcomm, &flg)); 503 PetscCheck(flg == MPI_CONGRUENT || flg == MPI_IDENT, rcomm, PETSC_ERR_SUP, "cannot map two layouts with non-matching communicators"); 504 PetscCall(PetscSFCreate(rcomm, sf)); 505 PetscCall(PetscLayoutGetLocalSize(rmap, &nroots)); 506 PetscCall(PetscLayoutGetSize(rmap, &rN)); 507 PetscCall(PetscLayoutGetRange(lmap, &lst, &len)); 508 PetscCall(PetscMalloc1(len - lst, &remote)); 509 for (i = lst; i < len && i < rN; i++) { 510 if (owner < -1 || i >= rmap->range[owner + 1]) PetscCall(PetscLayoutFindOwner(rmap, i, &owner)); 511 remote[nleaves].rank = owner; 512 remote[nleaves].index = i - rmap->range[owner]; 513 nleaves++; 514 } 515 PetscCall(PetscSFSetGraph(*sf, nroots, nleaves, NULL, PETSC_OWN_POINTER, remote, PETSC_COPY_VALUES)); 516 PetscCall(PetscFree(remote)); 517 PetscFunctionReturn(0); 518 } 519 520 /* TODO: handle nooffprocentries like MatZeroRowsMapLocal_Private, since this code is the same */ 521 PetscErrorCode PetscLayoutMapLocal(PetscLayout map, PetscInt N, const PetscInt idxs[], PetscInt *on, PetscInt **oidxs, PetscInt **ogidxs) 522 { 523 PetscInt *owners = map->range; 524 PetscInt n = map->n; 525 PetscSF sf; 526 PetscInt *lidxs, *work = NULL; 527 PetscSFNode *ridxs; 528 PetscMPIInt rank, p = 0; 529 PetscInt r, len = 0; 530 531 PetscFunctionBegin; 532 if (on) *on = 0; /* squelch -Wmaybe-uninitialized */ 533 /* Create SF where leaves are input idxs and roots are owned idxs */ 534 PetscCallMPI(MPI_Comm_rank(map->comm, &rank)); 535 PetscCall(PetscMalloc1(n, &lidxs)); 536 for (r = 0; r < n; ++r) lidxs[r] = -1; 537 PetscCall(PetscMalloc1(N, &ridxs)); 538 for (r = 0; r < N; ++r) { 539 const PetscInt idx = idxs[r]; 540 PetscCheck(idx >= 0 && idx < map->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Index %" PetscInt_FMT " out of range [0,%" PetscInt_FMT ")", idx, map->N); 541 if (idx < owners[p] || owners[p + 1] <= idx) { /* short-circuit the search if the last p owns this idx too */ 542 PetscCall(PetscLayoutFindOwner(map, idx, &p)); 543 } 544 ridxs[r].rank = p; 545 ridxs[r].index = idxs[r] - owners[p]; 546 } 547 PetscCall(PetscSFCreate(map->comm, &sf)); 548 PetscCall(PetscSFSetGraph(sf, n, N, NULL, PETSC_OWN_POINTER, ridxs, PETSC_OWN_POINTER)); 549 PetscCall(PetscSFReduceBegin(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR)); 550 PetscCall(PetscSFReduceEnd(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR)); 551 if (ogidxs) { /* communicate global idxs */ 552 PetscInt cum = 0, start, *work2; 553 554 PetscCall(PetscMalloc1(n, &work)); 555 PetscCall(PetscCalloc1(N, &work2)); 556 for (r = 0; r < N; ++r) 557 if (idxs[r] >= 0) cum++; 558 PetscCallMPI(MPI_Scan(&cum, &start, 1, MPIU_INT, MPI_SUM, map->comm)); 559 start -= cum; 560 cum = 0; 561 for (r = 0; r < N; ++r) 562 if (idxs[r] >= 0) work2[r] = start + cum++; 563 PetscCall(PetscSFReduceBegin(sf, MPIU_INT, work2, work, MPI_REPLACE)); 564 PetscCall(PetscSFReduceEnd(sf, MPIU_INT, work2, work, MPI_REPLACE)); 565 PetscCall(PetscFree(work2)); 566 } 567 PetscCall(PetscSFDestroy(&sf)); 568 /* Compress and put in indices */ 569 for (r = 0; r < n; ++r) 570 if (lidxs[r] >= 0) { 571 if (work) work[len] = work[r]; 572 lidxs[len++] = r; 573 } 574 if (on) *on = len; 575 if (oidxs) *oidxs = lidxs; 576 if (ogidxs) *ogidxs = work; 577 PetscFunctionReturn(0); 578 } 579 580 /*@ 581 PetscSFCreateByMatchingIndices - Create `PetscSF` by matching root and leaf indices 582 583 Collective 584 585 Input Parameters: 586 + layout - `PetscLayout` defining the global index space and the rank that brokers each index 587 . numRootIndices - size of rootIndices 588 . rootIndices - `PetscInt` array of global indices of which this process requests ownership 589 . rootLocalIndices - root local index permutation (NULL if no permutation) 590 . rootLocalOffset - offset to be added to root local indices 591 . numLeafIndices - size of leafIndices 592 . leafIndices - `PetscInt` array of global indices with which this process requires data associated 593 . leafLocalIndices - leaf local index permutation (NULL if no permutation) 594 - leafLocalOffset - offset to be added to leaf local indices 595 596 Output Parameters: 597 + sfA - star forest representing the communication pattern from the layout space to the leaf space (NULL if not needed) 598 - sf - star forest representing the communication pattern from the root space to the leaf space 599 600 Level: advanced 601 602 Example 1: 603 .vb 604 rank : 0 1 2 605 rootIndices : [1 0 2] [3] [3] 606 rootLocalOffset : 100 200 300 607 layout : [0 1] [2] [3] 608 leafIndices : [0] [2] [0 3] 609 leafLocalOffset : 400 500 600 610 611 would build the following PetscSF 612 613 [0] 400 <- (0,101) 614 [1] 500 <- (0,102) 615 [2] 600 <- (0,101) 616 [2] 601 <- (2,300) 617 .ve 618 619 Example 2: 620 .vb 621 rank : 0 1 2 622 rootIndices : [1 0 2] [3] [3] 623 rootLocalOffset : 100 200 300 624 layout : [0 1] [2] [3] 625 leafIndices : rootIndices rootIndices rootIndices 626 leafLocalOffset : rootLocalOffset rootLocalOffset rootLocalOffset 627 628 would build the following PetscSF 629 630 [1] 200 <- (2,300) 631 .ve 632 633 Example 3: 634 .vb 635 No process requests ownership of global index 1, but no process needs it. 636 637 rank : 0 1 2 638 numRootIndices : 2 1 1 639 rootIndices : [0 2] [3] [3] 640 rootLocalOffset : 100 200 300 641 layout : [0 1] [2] [3] 642 numLeafIndices : 1 1 2 643 leafIndices : [0] [2] [0 3] 644 leafLocalOffset : 400 500 600 645 646 would build the following PetscSF 647 648 [0] 400 <- (0,100) 649 [1] 500 <- (0,101) 650 [2] 600 <- (0,100) 651 [2] 601 <- (2,300) 652 .ve 653 654 Notes: 655 The layout parameter represents any partitioning of [0, N), where N is the total number of global indices, and its 656 local size can be set to `PETSC_DECIDE`. 657 658 If a global index x lies in the partition owned by process i, each process whose rootIndices contains x requests 659 ownership of x and sends its own rank and the local index of x to process i. 660 If multiple processes request ownership of x, the one with the highest rank is to own x. 661 Process i then broadcasts the ownership information, so that each process whose leafIndices contains x knows the 662 ownership information of x. 663 The output sf is constructed by associating each leaf point to a root point in this way. 664 665 Suppose there is point data ordered according to the global indices and partitioned according to the given layout. 666 The optional output `PetscSF` object sfA can be used to push such data to leaf points. 667 668 All indices in rootIndices and leafIndices must lie in the layout range. The union (over all processes) of rootIndices 669 must cover that of leafIndices, but need not cover the entire layout. 670 671 If (leafIndices, leafLocalIndices, leafLocalOffset) == (rootIndices, rootLocalIndices, rootLocalOffset), the output 672 star forest is almost identity, so will only include non-trivial part of the map. 673 674 Developer Note: 675 Current approach of a process of the highest rank gaining the ownership may cause load imbalance; consider using 676 hash(rank, root_local_index) as the bid for the ownership determination. 677 678 .seealso: `PetscSF`, `PetscSFCreate()` 679 @*/ 680 PetscErrorCode PetscSFCreateByMatchingIndices(PetscLayout layout, PetscInt numRootIndices, const PetscInt *rootIndices, const PetscInt *rootLocalIndices, PetscInt rootLocalOffset, PetscInt numLeafIndices, const PetscInt *leafIndices, const PetscInt *leafLocalIndices, PetscInt leafLocalOffset, PetscSF *sfA, PetscSF *sf) 681 { 682 MPI_Comm comm = layout->comm; 683 PetscMPIInt size, rank; 684 PetscSF sf1; 685 PetscSFNode *owners, *buffer, *iremote; 686 PetscInt *ilocal, nleaves, N, n, i; 687 #if defined(PETSC_USE_DEBUG) 688 PetscInt N1; 689 #endif 690 PetscBool flag; 691 692 PetscFunctionBegin; 693 if (rootIndices) PetscValidIntPointer(rootIndices, 3); 694 if (rootLocalIndices) PetscValidIntPointer(rootLocalIndices, 4); 695 if (leafIndices) PetscValidIntPointer(leafIndices, 7); 696 if (leafLocalIndices) PetscValidIntPointer(leafLocalIndices, 8); 697 if (sfA) PetscValidPointer(sfA, 10); 698 PetscValidPointer(sf, 11); 699 PetscCheck(numRootIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numRootIndices (%" PetscInt_FMT ") must be non-negative", numRootIndices); 700 PetscCheck(numLeafIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numLeafIndices (%" PetscInt_FMT ") must be non-negative", numLeafIndices); 701 PetscCallMPI(MPI_Comm_size(comm, &size)); 702 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 703 PetscCall(PetscLayoutSetUp(layout)); 704 PetscCall(PetscLayoutGetSize(layout, &N)); 705 PetscCall(PetscLayoutGetLocalSize(layout, &n)); 706 flag = (PetscBool)(leafIndices == rootIndices); 707 PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &flag, 1, MPIU_BOOL, MPI_LAND, comm)); 708 PetscCheck(!flag || numLeafIndices == numRootIndices, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "leafIndices == rootIndices, but numLeafIndices (%" PetscInt_FMT ") != numRootIndices(%" PetscInt_FMT ")", numLeafIndices, numRootIndices); 709 #if defined(PETSC_USE_DEBUG) 710 N1 = PETSC_MIN_INT; 711 for (i = 0; i < numRootIndices; i++) 712 if (rootIndices[i] > N1) N1 = rootIndices[i]; 713 PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm)); 714 PetscCheck(N1 < N, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Max. root index (%" PetscInt_FMT ") out of layout range [0,%" PetscInt_FMT ")", N1, N); 715 if (!flag) { 716 N1 = PETSC_MIN_INT; 717 for (i = 0; i < numLeafIndices; i++) 718 if (leafIndices[i] > N1) N1 = leafIndices[i]; 719 PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm)); 720 PetscCheck(N1 < N, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Max. leaf index (%" PetscInt_FMT ") out of layout range [0,%" PetscInt_FMT ")", N1, N); 721 } 722 #endif 723 /* Reduce: owners -> buffer */ 724 PetscCall(PetscMalloc1(n, &buffer)); 725 PetscCall(PetscSFCreate(comm, &sf1)); 726 PetscCall(PetscSFSetFromOptions(sf1)); 727 PetscCall(PetscSFSetGraphLayout(sf1, layout, numRootIndices, NULL, PETSC_OWN_POINTER, rootIndices)); 728 PetscCall(PetscMalloc1(numRootIndices, &owners)); 729 for (i = 0; i < numRootIndices; ++i) { 730 owners[i].rank = rank; 731 owners[i].index = rootLocalOffset + (rootLocalIndices ? rootLocalIndices[i] : i); 732 } 733 for (i = 0; i < n; ++i) { 734 buffer[i].index = -1; 735 buffer[i].rank = -1; 736 } 737 PetscCall(PetscSFReduceBegin(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC)); 738 PetscCall(PetscSFReduceEnd(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC)); 739 /* Bcast: buffer -> owners */ 740 if (!flag) { 741 /* leafIndices is different from rootIndices */ 742 PetscCall(PetscFree(owners)); 743 PetscCall(PetscSFSetGraphLayout(sf1, layout, numLeafIndices, NULL, PETSC_OWN_POINTER, leafIndices)); 744 PetscCall(PetscMalloc1(numLeafIndices, &owners)); 745 } 746 PetscCall(PetscSFBcastBegin(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE)); 747 PetscCall(PetscSFBcastEnd(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE)); 748 for (i = 0; i < numLeafIndices; ++i) PetscCheck(owners[i].rank >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Global point %" PetscInt_FMT " was unclaimed", leafIndices[i]); 749 PetscCall(PetscFree(buffer)); 750 if (sfA) { 751 *sfA = sf1; 752 } else PetscCall(PetscSFDestroy(&sf1)); 753 /* Create sf */ 754 if (flag && rootLocalIndices == leafLocalIndices && leafLocalOffset == rootLocalOffset) { 755 /* leaf space == root space */ 756 for (i = 0, nleaves = 0; i < numLeafIndices; ++i) 757 if (owners[i].rank != rank) ++nleaves; 758 PetscCall(PetscMalloc1(nleaves, &ilocal)); 759 PetscCall(PetscMalloc1(nleaves, &iremote)); 760 for (i = 0, nleaves = 0; i < numLeafIndices; ++i) { 761 if (owners[i].rank != rank) { 762 ilocal[nleaves] = leafLocalOffset + i; 763 iremote[nleaves].rank = owners[i].rank; 764 iremote[nleaves].index = owners[i].index; 765 ++nleaves; 766 } 767 } 768 PetscCall(PetscFree(owners)); 769 } else { 770 nleaves = numLeafIndices; 771 PetscCall(PetscMalloc1(nleaves, &ilocal)); 772 for (i = 0; i < nleaves; ++i) ilocal[i] = leafLocalOffset + (leafLocalIndices ? leafLocalIndices[i] : i); 773 iremote = owners; 774 } 775 PetscCall(PetscSFCreate(comm, sf)); 776 PetscCall(PetscSFSetFromOptions(*sf)); 777 PetscCall(PetscSFSetGraph(*sf, rootLocalOffset + numRootIndices, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER)); 778 PetscFunctionReturn(0); 779 } 780 781 /*@ 782 PetscSFMerge - append/merge indices of sfb into sfa, with preference for sfb 783 784 Collective 785 786 Input Arguments: 787 + sfa - default `PetscSF` 788 - sfb - additional edges to add/replace edges in sfa 789 790 Output Arguments: 791 . merged - new `PetscSF` with combined edges 792 793 .seealse: `PetscSFCompose()` 794 @*/ 795 PetscErrorCode PetscSFMerge(PetscSF sfa, PetscSF sfb, PetscSF *merged) 796 { 797 PetscInt maxleaf; 798 799 PetscFunctionBegin; 800 PetscValidHeaderSpecific(sfa, PETSCSF_CLASSID, 1); 801 PetscValidHeaderSpecific(sfb, PETSCSF_CLASSID, 2); 802 PetscCheckSameComm(sfa, 1, sfb, 2); 803 PetscValidPointer(merged, 3); 804 { 805 PetscInt aleaf, bleaf; 806 PetscCall(PetscSFGetLeafRange(sfa, NULL, &aleaf)); 807 PetscCall(PetscSFGetLeafRange(sfb, NULL, &bleaf)); 808 maxleaf = PetscMax(aleaf, bleaf) + 1; // One more than the last index 809 } 810 PetscInt *clocal, aroots, aleaves, broots, bleaves; 811 PetscSFNode *cremote; 812 const PetscInt *alocal, *blocal; 813 const PetscSFNode *aremote, *bremote; 814 PetscCall(PetscMalloc2(maxleaf, &clocal, maxleaf, &cremote)); 815 for (PetscInt i = 0; i < maxleaf; i++) clocal[i] = -1; 816 PetscCall(PetscSFGetGraph(sfa, &aroots, &aleaves, &alocal, &aremote)); 817 PetscCall(PetscSFGetGraph(sfb, &broots, &bleaves, &blocal, &bremote)); 818 PetscCheck(aroots == broots, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Both sfa and sfb must have the same root space"); 819 for (PetscInt i = 0; i < aleaves; i++) { 820 PetscInt a = alocal ? alocal[i] : i; 821 clocal[a] = a; 822 cremote[a] = aremote[i]; 823 } 824 for (PetscInt i = 0; i < bleaves; i++) { 825 PetscInt b = blocal ? blocal[i] : i; 826 clocal[b] = b; 827 cremote[b] = bremote[i]; 828 } 829 PetscInt nleaves = 0; 830 for (PetscInt i = 0; i < maxleaf; i++) { 831 if (clocal[i] < 0) continue; 832 clocal[nleaves] = clocal[i]; 833 cremote[nleaves] = cremote[i]; 834 nleaves++; 835 } 836 PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)sfa), merged)); 837 PetscCall(PetscSFSetGraph(*merged, aroots, nleaves, clocal, PETSC_COPY_VALUES, cremote, PETSC_COPY_VALUES)); 838 PetscCall(PetscFree2(clocal, cremote)); 839 PetscFunctionReturn(0); 840 } 841