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