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 Notes: 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 Notes: 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: `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(PetscLayoutGetLocalSize(layout, &nroots)); 39 PetscCall(PetscLayoutGetRanges(layout, &range)); 40 PetscCall(PetscMalloc1(nleaves, &remote)); 41 if (nleaves) ls = gremote[0] + 1; 42 for (i = 0; i < nleaves; i++) { 43 const PetscInt idx = gremote[i] - ls; 44 if (idx < 0 || idx >= ln) { /* short-circuit the search */ 45 PetscCall(PetscLayoutFindOwnerIndex(layout, gremote[i], &lr, &remote[i].index)); 46 remote[i].rank = lr; 47 ls = range[lr]; 48 ln = range[lr + 1] - ls; 49 } else { 50 remote[i].rank = lr; 51 remote[i].index = idx; 52 } 53 } 54 PetscCall(PetscSFSetGraph(sf, nroots, nleaves, ilocal, localmode, remote, PETSC_OWN_POINTER)); 55 PetscFunctionReturn(0); 56 } 57 58 /*@C 59 PetscSFGetGraphLayout - Get the global indices and PetscLayout that describe this star forest 60 61 Collective 62 63 Input Parameter: 64 . sf - star forest 65 66 Output Parameters: 67 + layout - PetscLayout defining the global space for roots 68 . nleaves - number of leaf vertices on the current process, each of these references a root on any process 69 . ilocal - locations of leaves in leafdata buffers, or NULL for contiguous storage 70 - gremote - root vertices in global numbering corresponding to leaves in ilocal 71 72 Level: intermediate 73 74 Notes: 75 The outputs are such that passing them as inputs to `PetscSFSetGraphLayout()` would lead to the same star forest. 76 The outputs layout and gremote are freshly created each time this function is called, 77 so they need to be freed by user and cannot be qualified as const. 78 79 80 .seealso: `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 PetscSections 109 describing the data layout. 110 111 Input Parameters: 112 + sf - The SF 113 . localSection - PetscSection describing the local data layout 114 - globalSection - PetscSection describing the global data layout 115 116 Level: developer 117 118 .seealso: `PetscSFSetGraph()`, `PetscSFSetGraphLayout()` 119 @*/ 120 PetscErrorCode PetscSFSetGraphSection(PetscSF sf, PetscSection localSection, PetscSection globalSection) 121 { 122 MPI_Comm comm; 123 PetscLayout layout; 124 const PetscInt *ranges; 125 PetscInt *local; 126 PetscSFNode *remote; 127 PetscInt pStart, pEnd, p, nroots, nleaves = 0, l; 128 PetscMPIInt size, rank; 129 130 PetscFunctionBegin; 131 PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 132 PetscValidHeaderSpecific(localSection, PETSC_SECTION_CLASSID, 2); 133 PetscValidHeaderSpecific(globalSection, PETSC_SECTION_CLASSID, 3); 134 135 PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 136 PetscCallMPI(MPI_Comm_size(comm, &size)); 137 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 138 PetscCall(PetscSectionGetChart(globalSection, &pStart, &pEnd)); 139 PetscCall(PetscSectionGetConstrainedStorageSize(globalSection, &nroots)); 140 PetscCall(PetscLayoutCreate(comm, &layout)); 141 PetscCall(PetscLayoutSetBlockSize(layout, 1)); 142 PetscCall(PetscLayoutSetLocalSize(layout, nroots)); 143 PetscCall(PetscLayoutSetUp(layout)); 144 PetscCall(PetscLayoutGetRanges(layout, &ranges)); 145 for (p = pStart; p < pEnd; ++p) { 146 PetscInt gdof, gcdof; 147 148 PetscCall(PetscSectionGetDof(globalSection, p, &gdof)); 149 PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof)); 150 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)); 151 nleaves += gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof; 152 } 153 PetscCall(PetscMalloc1(nleaves, &local)); 154 PetscCall(PetscMalloc1(nleaves, &remote)); 155 for (p = pStart, l = 0; p < pEnd; ++p) { 156 const PetscInt *cind; 157 PetscInt dof, cdof, off, gdof, gcdof, goff, gsize, d, c; 158 159 PetscCall(PetscSectionGetDof(localSection, p, &dof)); 160 PetscCall(PetscSectionGetOffset(localSection, p, &off)); 161 PetscCall(PetscSectionGetConstraintDof(localSection, p, &cdof)); 162 PetscCall(PetscSectionGetConstraintIndices(localSection, p, &cind)); 163 PetscCall(PetscSectionGetDof(globalSection, p, &gdof)); 164 PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof)); 165 PetscCall(PetscSectionGetOffset(globalSection, p, &goff)); 166 if (!gdof) continue; /* Censored point */ 167 gsize = gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof; 168 if (gsize != dof - cdof) { 169 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); 170 cdof = 0; /* Ignore constraints */ 171 } 172 for (d = 0, c = 0; d < dof; ++d) { 173 if ((c < cdof) && (cind[c] == d)) { 174 ++c; 175 continue; 176 } 177 local[l + d - c] = off + d; 178 } 179 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); 180 if (gdof < 0) { 181 for (d = 0; d < gsize; ++d, ++l) { 182 PetscInt offset = -(goff + 1) + d, r; 183 184 PetscCall(PetscFindInt(offset, size + 1, ranges, &r)); 185 if (r < 0) r = -(r + 2); 186 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); 187 remote[l].rank = r; 188 remote[l].index = offset - ranges[r]; 189 } 190 } else { 191 for (d = 0; d < gsize; ++d, ++l) { 192 remote[l].rank = rank; 193 remote[l].index = goff + d - ranges[rank]; 194 } 195 } 196 } 197 PetscCheck(l == nleaves, comm, PETSC_ERR_PLIB, "Iteration error, l %" PetscInt_FMT " != nleaves %" PetscInt_FMT, l, nleaves); 198 PetscCall(PetscLayoutDestroy(&layout)); 199 PetscCall(PetscSFSetGraph(sf, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER)); 200 PetscFunctionReturn(0); 201 } 202 203 /*@C 204 PetscSFDistributeSection - Create a new PetscSection reorganized, moving from the root to the leaves of the SF 205 206 Collective on sf 207 208 Input Parameters: 209 + sf - The SF 210 - rootSection - Section defined on root space 211 212 Output Parameters: 213 + remoteOffsets - root offsets in leaf storage, or NULL 214 - leafSection - Section defined on the leaf space 215 216 Level: advanced 217 218 .seealso: `PetscSFCreate()` 219 @*/ 220 PetscErrorCode PetscSFDistributeSection(PetscSF sf, PetscSection rootSection, PetscInt **remoteOffsets, PetscSection leafSection) 221 { 222 PetscSF embedSF; 223 const PetscInt *indices; 224 IS selected; 225 PetscInt numFields, nroots, rpStart, rpEnd, lpStart = PETSC_MAX_INT, lpEnd = -1, f, c; 226 PetscBool *sub, hasc; 227 228 PetscFunctionBegin; 229 PetscCall(PetscLogEventBegin(PETSCSF_DistSect, sf, 0, 0, 0)); 230 PetscCall(PetscSectionGetNumFields(rootSection, &numFields)); 231 if (numFields) { 232 IS perm; 233 234 /* PetscSectionSetNumFields() calls PetscSectionReset(), which destroys 235 leafSection->perm. To keep this permutation set by the user, we grab 236 the reference before calling PetscSectionSetNumFields() and set it 237 back after. */ 238 PetscCall(PetscSectionGetPermutation(leafSection, &perm)); 239 PetscCall(PetscObjectReference((PetscObject)perm)); 240 PetscCall(PetscSectionSetNumFields(leafSection, numFields)); 241 PetscCall(PetscSectionSetPermutation(leafSection, perm)); 242 PetscCall(ISDestroy(&perm)); 243 } 244 PetscCall(PetscMalloc1(numFields + 2, &sub)); 245 sub[1] = rootSection->bc ? PETSC_TRUE : PETSC_FALSE; 246 for (f = 0; f < numFields; ++f) { 247 PetscSectionSym sym, dsym = NULL; 248 const char *name = NULL; 249 PetscInt numComp = 0; 250 251 sub[2 + f] = rootSection->field[f]->bc ? PETSC_TRUE : PETSC_FALSE; 252 PetscCall(PetscSectionGetFieldComponents(rootSection, f, &numComp)); 253 PetscCall(PetscSectionGetFieldName(rootSection, f, &name)); 254 PetscCall(PetscSectionGetFieldSym(rootSection, f, &sym)); 255 if (sym) PetscCall(PetscSectionSymDistribute(sym, sf, &dsym)); 256 PetscCall(PetscSectionSetFieldComponents(leafSection, f, numComp)); 257 PetscCall(PetscSectionSetFieldName(leafSection, f, name)); 258 PetscCall(PetscSectionSetFieldSym(leafSection, f, dsym)); 259 PetscCall(PetscSectionSymDestroy(&dsym)); 260 for (c = 0; c < rootSection->numFieldComponents[f]; ++c) { 261 PetscCall(PetscSectionGetComponentName(rootSection, f, c, &name)); 262 PetscCall(PetscSectionSetComponentName(leafSection, f, c, name)); 263 } 264 } 265 PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd)); 266 PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL)); 267 rpEnd = PetscMin(rpEnd, nroots); 268 rpEnd = PetscMax(rpStart, rpEnd); 269 /* see if we can avoid creating the embedded SF, since it can cost more than an allreduce */ 270 sub[0] = (PetscBool)(nroots != rpEnd - rpStart); 271 PetscCall(MPIU_Allreduce(MPI_IN_PLACE, sub, 2 + numFields, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)sf))); 272 if (sub[0]) { 273 PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected)); 274 PetscCall(ISGetIndices(selected, &indices)); 275 PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF)); 276 PetscCall(ISRestoreIndices(selected, &indices)); 277 PetscCall(ISDestroy(&selected)); 278 } else { 279 PetscCall(PetscObjectReference((PetscObject)sf)); 280 embedSF = sf; 281 } 282 PetscCall(PetscSFGetLeafRange(embedSF, &lpStart, &lpEnd)); 283 lpEnd++; 284 285 PetscCall(PetscSectionSetChart(leafSection, lpStart, lpEnd)); 286 287 /* Constrained dof section */ 288 hasc = sub[1]; 289 for (f = 0; f < numFields; ++f) hasc = (PetscBool)(hasc || sub[2 + f]); 290 291 /* Could fuse these at the cost of copies and extra allocation */ 292 PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart], MPI_REPLACE)); 293 PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart], MPI_REPLACE)); 294 if (sub[1]) { 295 PetscCall(PetscSectionCheckConstraints_Private(rootSection)); 296 PetscCall(PetscSectionCheckConstraints_Private(leafSection)); 297 PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE)); 298 PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE)); 299 } 300 for (f = 0; f < numFields; ++f) { 301 PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart], MPI_REPLACE)); 302 PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart], MPI_REPLACE)); 303 if (sub[2 + f]) { 304 PetscCall(PetscSectionCheckConstraints_Private(rootSection->field[f])); 305 PetscCall(PetscSectionCheckConstraints_Private(leafSection->field[f])); 306 PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE)); 307 PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE)); 308 } 309 } 310 if (remoteOffsets) { 311 PetscCall(PetscMalloc1(lpEnd - lpStart, remoteOffsets)); 312 PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE)); 313 PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE)); 314 } 315 PetscCall(PetscSectionInvalidateMaxDof_Internal(leafSection)); 316 PetscCall(PetscSectionSetUp(leafSection)); 317 if (hasc) { /* need to communicate bcIndices */ 318 PetscSF bcSF; 319 PetscInt *rOffBc; 320 321 PetscCall(PetscMalloc1(lpEnd - lpStart, &rOffBc)); 322 if (sub[1]) { 323 PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 324 PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 325 PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->bc, rOffBc, leafSection->bc, &bcSF)); 326 PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE)); 327 PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE)); 328 PetscCall(PetscSFDestroy(&bcSF)); 329 } 330 for (f = 0; f < numFields; ++f) { 331 if (sub[2 + f]) { 332 PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 333 PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 334 PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->field[f]->bc, rOffBc, leafSection->field[f]->bc, &bcSF)); 335 PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE)); 336 PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE)); 337 PetscCall(PetscSFDestroy(&bcSF)); 338 } 339 } 340 PetscCall(PetscFree(rOffBc)); 341 } 342 PetscCall(PetscSFDestroy(&embedSF)); 343 PetscCall(PetscFree(sub)); 344 PetscCall(PetscLogEventEnd(PETSCSF_DistSect, sf, 0, 0, 0)); 345 PetscFunctionReturn(0); 346 } 347 348 /*@C 349 PetscSFCreateRemoteOffsets - Create offsets for point data on remote processes 350 351 Collective on sf 352 353 Input Parameters: 354 + sf - The SF 355 . rootSection - Data layout of remote points for outgoing data (this is layout for SF roots) 356 - leafSection - Data layout of local points for incoming data (this is layout for SF leaves) 357 358 Output Parameter: 359 . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL 360 361 Level: developer 362 363 .seealso: `PetscSFCreate()` 364 @*/ 365 PetscErrorCode PetscSFCreateRemoteOffsets(PetscSF sf, PetscSection rootSection, PetscSection leafSection, PetscInt **remoteOffsets) 366 { 367 PetscSF embedSF; 368 const PetscInt *indices; 369 IS selected; 370 PetscInt numRoots, rpStart = 0, rpEnd = 0, lpStart = 0, lpEnd = 0; 371 372 PetscFunctionBegin; 373 *remoteOffsets = NULL; 374 PetscCall(PetscSFGetGraph(sf, &numRoots, NULL, NULL, NULL)); 375 if (numRoots < 0) PetscFunctionReturn(0); 376 PetscCall(PetscLogEventBegin(PETSCSF_RemoteOff, sf, 0, 0, 0)); 377 PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd)); 378 PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd)); 379 PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected)); 380 PetscCall(ISGetIndices(selected, &indices)); 381 PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF)); 382 PetscCall(ISRestoreIndices(selected, &indices)); 383 PetscCall(ISDestroy(&selected)); 384 PetscCall(PetscCalloc1(lpEnd - lpStart, remoteOffsets)); 385 PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE)); 386 PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE)); 387 PetscCall(PetscSFDestroy(&embedSF)); 388 PetscCall(PetscLogEventEnd(PETSCSF_RemoteOff, sf, 0, 0, 0)); 389 PetscFunctionReturn(0); 390 } 391 392 /*@C 393 PetscSFCreateSectionSF - Create an expanded SF of dofs, assuming the input SF relates points 394 395 Collective on sf 396 397 Input Parameters: 398 + sf - The SF 399 . rootSection - Data layout of remote points for outgoing data (this is usually the serial section) 400 . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL 401 - leafSection - Data layout of local points for incoming data (this is the distributed section) 402 403 Output Parameters: 404 - sectionSF - The new SF 405 406 Note: Either rootSection or remoteOffsets can be specified 407 408 Level: advanced 409 410 .seealso: `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: `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 SF 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 Example 1: 601 $ 602 $ rank : 0 1 2 603 $ rootIndices : [1 0 2] [3] [3] 604 $ rootLocalOffset : 100 200 300 605 $ layout : [0 1] [2] [3] 606 $ leafIndices : [0] [2] [0 3] 607 $ leafLocalOffset : 400 500 600 608 $ 609 would build the following SF 610 $ 611 $ [0] 400 <- (0,101) 612 $ [1] 500 <- (0,102) 613 $ [2] 600 <- (0,101) 614 $ [2] 601 <- (2,300) 615 $ 616 Example 2: 617 $ 618 $ rank : 0 1 2 619 $ rootIndices : [1 0 2] [3] [3] 620 $ rootLocalOffset : 100 200 300 621 $ layout : [0 1] [2] [3] 622 $ leafIndices : rootIndices rootIndices rootIndices 623 $ leafLocalOffset : rootLocalOffset rootLocalOffset rootLocalOffset 624 $ 625 would build the following SF 626 $ 627 $ [1] 200 <- (2,300) 628 $ 629 Example 3: 630 $ 631 $ No process requests ownership of global index 1, but no process needs it. 632 $ 633 $ rank : 0 1 2 634 $ numRootIndices : 2 1 1 635 $ rootIndices : [0 2] [3] [3] 636 $ rootLocalOffset : 100 200 300 637 $ layout : [0 1] [2] [3] 638 $ numLeafIndices : 1 1 2 639 $ leafIndices : [0] [2] [0 3] 640 $ leafLocalOffset : 400 500 600 641 $ 642 would build the following SF 643 $ 644 $ [0] 400 <- (0,100) 645 $ [1] 500 <- (0,101) 646 $ [2] 600 <- (0,100) 647 $ [2] 601 <- (2,300) 648 $ 649 650 Notes: 651 The layout parameter represents any partitioning of [0, N), where N is the total number of global indices, and its 652 local size can be set to PETSC_DECIDE. 653 If a global index x lies in the partition owned by process i, each process whose rootIndices contains x requests 654 ownership of x and sends its own rank and the local index of x to process i. 655 If multiple processes request ownership of x, the one with the highest rank is to own x. 656 Process i then broadcasts the ownership information, so that each process whose leafIndices contains x knows the 657 ownership information of x. 658 The output sf is constructed by associating each leaf point to a root point in this way. 659 660 Suppose there is point data ordered according to the global indices and partitioned according to the given layout. 661 The optional output PetscSF object sfA can be used to push such data to leaf points. 662 663 All indices in rootIndices and leafIndices must lie in the layout range. The union (over all processes) of rootIndices 664 must cover that of leafIndices, but need not cover the entire layout. 665 666 If (leafIndices, leafLocalIndices, leafLocalOffset) == (rootIndices, rootLocalIndices, rootLocalOffset), the output 667 star forest is almost identity, so will only include non-trivial part of the map. 668 669 Developer Notes: 670 Current approach of a process of the highest rank gaining the ownership may cause load imbalance; consider using 671 hash(rank, root_local_index) as the bid for the ownership determination. 672 673 Level: advanced 674 675 .seealso: `PetscSFCreate()` 676 @*/ 677 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) 678 { 679 MPI_Comm comm = layout->comm; 680 PetscMPIInt size, rank; 681 PetscSF sf1; 682 PetscSFNode *owners, *buffer, *iremote; 683 PetscInt *ilocal, nleaves, N, n, i; 684 #if defined(PETSC_USE_DEBUG) 685 PetscInt N1; 686 #endif 687 PetscBool flag; 688 689 PetscFunctionBegin; 690 if (rootIndices) PetscValidIntPointer(rootIndices, 3); 691 if (rootLocalIndices) PetscValidIntPointer(rootLocalIndices, 4); 692 if (leafIndices) PetscValidIntPointer(leafIndices, 7); 693 if (leafLocalIndices) PetscValidIntPointer(leafLocalIndices, 8); 694 if (sfA) PetscValidPointer(sfA, 10); 695 PetscValidPointer(sf, 11); 696 PetscCheck(numRootIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numRootIndices (%" PetscInt_FMT ") must be non-negative", numRootIndices); 697 PetscCheck(numLeafIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numLeafIndices (%" PetscInt_FMT ") must be non-negative", numLeafIndices); 698 PetscCallMPI(MPI_Comm_size(comm, &size)); 699 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 700 PetscCall(PetscLayoutSetUp(layout)); 701 PetscCall(PetscLayoutGetSize(layout, &N)); 702 PetscCall(PetscLayoutGetLocalSize(layout, &n)); 703 flag = (PetscBool)(leafIndices == rootIndices); 704 PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &flag, 1, MPIU_BOOL, MPI_LAND, comm)); 705 PetscCheck(!flag || numLeafIndices == numRootIndices, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "leafIndices == rootIndices, but numLeafIndices (%" PetscInt_FMT ") != numRootIndices(%" PetscInt_FMT ")", numLeafIndices, numRootIndices); 706 #if defined(PETSC_USE_DEBUG) 707 N1 = PETSC_MIN_INT; 708 for (i = 0; i < numRootIndices; i++) 709 if (rootIndices[i] > N1) N1 = rootIndices[i]; 710 PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm)); 711 PetscCheck(N1 < N, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Max. root index (%" PetscInt_FMT ") out of layout range [0,%" PetscInt_FMT ")", N1, N); 712 if (!flag) { 713 N1 = PETSC_MIN_INT; 714 for (i = 0; i < numLeafIndices; i++) 715 if (leafIndices[i] > N1) N1 = leafIndices[i]; 716 PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm)); 717 PetscCheck(N1 < N, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Max. leaf index (%" PetscInt_FMT ") out of layout range [0,%" PetscInt_FMT ")", N1, N); 718 } 719 #endif 720 /* Reduce: owners -> buffer */ 721 PetscCall(PetscMalloc1(n, &buffer)); 722 PetscCall(PetscSFCreate(comm, &sf1)); 723 PetscCall(PetscSFSetFromOptions(sf1)); 724 PetscCall(PetscSFSetGraphLayout(sf1, layout, numRootIndices, NULL, PETSC_OWN_POINTER, rootIndices)); 725 PetscCall(PetscMalloc1(numRootIndices, &owners)); 726 for (i = 0; i < numRootIndices; ++i) { 727 owners[i].rank = rank; 728 owners[i].index = rootLocalOffset + (rootLocalIndices ? rootLocalIndices[i] : i); 729 } 730 for (i = 0; i < n; ++i) { 731 buffer[i].index = -1; 732 buffer[i].rank = -1; 733 } 734 PetscCall(PetscSFReduceBegin(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC)); 735 PetscCall(PetscSFReduceEnd(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC)); 736 /* Bcast: buffer -> owners */ 737 if (!flag) { 738 /* leafIndices is different from rootIndices */ 739 PetscCall(PetscFree(owners)); 740 PetscCall(PetscSFSetGraphLayout(sf1, layout, numLeafIndices, NULL, PETSC_OWN_POINTER, leafIndices)); 741 PetscCall(PetscMalloc1(numLeafIndices, &owners)); 742 } 743 PetscCall(PetscSFBcastBegin(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE)); 744 PetscCall(PetscSFBcastEnd(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE)); 745 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]); 746 PetscCall(PetscFree(buffer)); 747 if (sfA) { 748 *sfA = sf1; 749 } else PetscCall(PetscSFDestroy(&sf1)); 750 /* Create sf */ 751 if (flag && rootLocalIndices == leafLocalIndices && leafLocalOffset == rootLocalOffset) { 752 /* leaf space == root space */ 753 for (i = 0, nleaves = 0; i < numLeafIndices; ++i) 754 if (owners[i].rank != rank) ++nleaves; 755 PetscCall(PetscMalloc1(nleaves, &ilocal)); 756 PetscCall(PetscMalloc1(nleaves, &iremote)); 757 for (i = 0, nleaves = 0; i < numLeafIndices; ++i) { 758 if (owners[i].rank != rank) { 759 ilocal[nleaves] = leafLocalOffset + i; 760 iremote[nleaves].rank = owners[i].rank; 761 iremote[nleaves].index = owners[i].index; 762 ++nleaves; 763 } 764 } 765 PetscCall(PetscFree(owners)); 766 } else { 767 nleaves = numLeafIndices; 768 PetscCall(PetscMalloc1(nleaves, &ilocal)); 769 for (i = 0; i < nleaves; ++i) ilocal[i] = leafLocalOffset + (leafLocalIndices ? leafLocalIndices[i] : i); 770 iremote = owners; 771 } 772 PetscCall(PetscSFCreate(comm, sf)); 773 PetscCall(PetscSFSetFromOptions(*sf)); 774 PetscCall(PetscSFSetGraph(*sf, rootLocalOffset + numRootIndices, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER)); 775 PetscFunctionReturn(0); 776 } 777