1 #include <../src/vec/is/sf/impls/basic/sfbasic.h> 2 #include <../src/vec/is/sf/impls/basic/sfpack.h> 3 #include <petsc/private/viewerimpl.h> 4 5 /*===================================================================================*/ 6 /* SF public interface implementations */ 7 /*===================================================================================*/ 8 PETSC_INTERN PetscErrorCode PetscSFSetUp_Basic(PetscSF sf) { 9 PetscSF_Basic *bas = (PetscSF_Basic *)sf->data; 10 PetscInt *rlengths, *ilengths, i, nRemoteRootRanks, nRemoteLeafRanks; 11 PetscMPIInt rank, niranks, *iranks, tag; 12 MPI_Comm comm; 13 MPI_Group group; 14 MPI_Request *rootreqs, *leafreqs; 15 16 PetscFunctionBegin; 17 PetscCallMPI(MPI_Comm_group(PETSC_COMM_SELF, &group)); 18 PetscCall(PetscSFSetUpRanks(sf, group)); 19 PetscCallMPI(MPI_Group_free(&group)); 20 PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 21 PetscCall(PetscObjectGetNewTag((PetscObject)sf, &tag)); 22 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 23 /* 24 * Inform roots about how many leaves and from which ranks 25 */ 26 PetscCall(PetscMalloc1(sf->nranks, &rlengths)); 27 /* Determine number, sending ranks and length of incoming */ 28 for (i = 0; i < sf->nranks; i++) { rlengths[i] = sf->roffset[i + 1] - sf->roffset[i]; /* Number of roots referenced by my leaves; for rank sf->ranks[i] */ } 29 nRemoteRootRanks = sf->nranks - sf->ndranks; 30 PetscCall(PetscCommBuildTwoSided(comm, 1, MPIU_INT, nRemoteRootRanks, sf->ranks + sf->ndranks, rlengths + sf->ndranks, &niranks, &iranks, (void **)&ilengths)); 31 32 /* Sort iranks. See use of VecScatterGetRemoteOrdered_Private() in MatGetBrowsOfAoCols_MPIAIJ() on why. 33 We could sort ranks there at the price of allocating extra working arrays. Presumably, niranks is 34 small and the sorting is cheap. 35 */ 36 PetscCall(PetscSortMPIIntWithIntArray(niranks, iranks, ilengths)); 37 38 /* Partition into distinguished and non-distinguished incoming ranks */ 39 bas->ndiranks = sf->ndranks; 40 bas->niranks = bas->ndiranks + niranks; 41 PetscCall(PetscMalloc2(bas->niranks, &bas->iranks, bas->niranks + 1, &bas->ioffset)); 42 bas->ioffset[0] = 0; 43 for (i = 0; i < bas->ndiranks; i++) { 44 bas->iranks[i] = sf->ranks[i]; 45 bas->ioffset[i + 1] = bas->ioffset[i] + rlengths[i]; 46 } 47 PetscCheck(bas->ndiranks <= 1 && (bas->ndiranks != 1 || bas->iranks[0] == rank), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Broken setup for shared ranks"); 48 for (; i < bas->niranks; i++) { 49 bas->iranks[i] = iranks[i - bas->ndiranks]; 50 bas->ioffset[i + 1] = bas->ioffset[i] + ilengths[i - bas->ndiranks]; 51 } 52 bas->itotal = bas->ioffset[i]; 53 PetscCall(PetscFree(rlengths)); 54 PetscCall(PetscFree(iranks)); 55 PetscCall(PetscFree(ilengths)); 56 57 /* Send leaf identities to roots */ 58 nRemoteLeafRanks = bas->niranks - bas->ndiranks; 59 PetscCall(PetscMalloc1(bas->itotal, &bas->irootloc)); 60 PetscCall(PetscMalloc2(nRemoteLeafRanks, &rootreqs, nRemoteRootRanks, &leafreqs)); 61 for (i = bas->ndiranks; i < bas->niranks; i++) PetscCallMPI(MPIU_Irecv(bas->irootloc + bas->ioffset[i], bas->ioffset[i + 1] - bas->ioffset[i], MPIU_INT, bas->iranks[i], tag, comm, &rootreqs[i - bas->ndiranks])); 62 for (i = 0; i < sf->nranks; i++) { 63 PetscInt npoints = sf->roffset[i + 1] - sf->roffset[i]; 64 if (i < sf->ndranks) { 65 PetscCheck(sf->ranks[i] == rank, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot interpret distinguished leaf rank"); 66 PetscCheck(bas->iranks[0] == rank, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot interpret distinguished root rank"); 67 PetscCheck(npoints == bas->ioffset[1] - bas->ioffset[0], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Distinguished rank exchange has mismatched lengths"); 68 PetscCall(PetscArraycpy(bas->irootloc + bas->ioffset[0], sf->rremote + sf->roffset[i], npoints)); 69 continue; 70 } 71 PetscCallMPI(MPIU_Isend(sf->rremote + sf->roffset[i], npoints, MPIU_INT, sf->ranks[i], tag, comm, &leafreqs[i - sf->ndranks])); 72 } 73 PetscCallMPI(MPI_Waitall(nRemoteLeafRanks, rootreqs, MPI_STATUSES_IGNORE)); 74 PetscCallMPI(MPI_Waitall(nRemoteRootRanks, leafreqs, MPI_STATUSES_IGNORE)); 75 76 sf->nleafreqs = nRemoteRootRanks; 77 bas->nrootreqs = nRemoteLeafRanks; 78 sf->persistent = PETSC_TRUE; 79 80 /* Setup fields related to packing, such as rootbuflen[] */ 81 PetscCall(PetscSFSetUpPackFields(sf)); 82 PetscCall(PetscFree2(rootreqs, leafreqs)); 83 PetscFunctionReturn(0); 84 } 85 86 PETSC_INTERN PetscErrorCode PetscSFReset_Basic(PetscSF sf) { 87 PetscSF_Basic *bas = (PetscSF_Basic *)sf->data; 88 PetscSFLink link = bas->avail, next; 89 90 PetscFunctionBegin; 91 PetscCheck(!bas->inuse, PetscObjectComm((PetscObject)sf), PETSC_ERR_ARG_WRONGSTATE, "Outstanding operation has not been completed"); 92 PetscCall(PetscFree2(bas->iranks, bas->ioffset)); 93 PetscCall(PetscFree(bas->irootloc)); 94 95 #if defined(PETSC_HAVE_DEVICE) 96 for (PetscInt i = 0; i < 2; i++) PetscCall(PetscSFFree(sf, PETSC_MEMTYPE_DEVICE, bas->irootloc_d[i])); 97 #endif 98 99 #if defined(PETSC_HAVE_NVSHMEM) 100 PetscCall(PetscSFReset_Basic_NVSHMEM(sf)); 101 #endif 102 103 for (; link; link = next) { 104 next = link->next; 105 PetscCall(PetscSFLinkDestroy(sf, link)); 106 } 107 bas->avail = NULL; 108 PetscCall(PetscSFResetPackFields(sf)); 109 PetscFunctionReturn(0); 110 } 111 112 PETSC_INTERN PetscErrorCode PetscSFDestroy_Basic(PetscSF sf) { 113 PetscFunctionBegin; 114 PetscCall(PetscSFReset_Basic(sf)); 115 PetscCall(PetscFree(sf->data)); 116 PetscFunctionReturn(0); 117 } 118 119 #if defined(PETSC_USE_SINGLE_LIBRARY) 120 #include <petscmat.h> 121 122 PETSC_INTERN PetscErrorCode PetscSFView_Basic_PatternAndSizes(PetscSF sf, PetscViewer viewer) { 123 PetscSF_Basic *bas = (PetscSF_Basic *)sf->data; 124 PetscInt i, nrootranks, ndrootranks; 125 const PetscInt *rootoffset; 126 PetscMPIInt rank, size; 127 const PetscMPIInt *rootranks; 128 MPI_Comm comm = PetscObjectComm((PetscObject)sf); 129 PetscScalar unitbytes; 130 Mat A; 131 132 PetscFunctionBegin; 133 PetscCallMPI(MPI_Comm_size(comm, &size)); 134 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 135 /* PetscSFView is most useful for the SF used in VecScatterBegin/End in MatMult etc, where we do 136 PetscSFBcast, i.e., roots send data to leaves. We dump the communication pattern into a matrix 137 in senders' view point: how many bytes I will send to my neighbors. 138 139 Looking at a column of the matrix, one can also know how many bytes the rank will receive from others. 140 141 If PetscSFLink bas->inuse is available, we can use that to get tree vertex size. But that would give 142 different interpretations for the same SF for different data types. Since we most care about VecScatter, 143 we uniformly treat each vertex as a PetscScalar. 144 */ 145 unitbytes = (PetscScalar)sizeof(PetscScalar); 146 147 PetscCall(PetscSFGetRootInfo_Basic(sf, &nrootranks, &ndrootranks, &rootranks, &rootoffset, NULL)); 148 PetscCall(MatCreateAIJ(comm, 1, 1, size, size, 1, NULL, nrootranks - ndrootranks, NULL, &A)); 149 PetscCall(MatSetOptionsPrefix(A, "__petsc_internal__")); /* To prevent the internal A from taking any command line options */ 150 for (i = 0; i < nrootranks; i++) PetscCall(MatSetValue(A, (PetscInt)rank, bas->iranks[i], (rootoffset[i + 1] - rootoffset[i]) * unitbytes, INSERT_VALUES)); 151 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 152 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 153 PetscCall(MatView(A, viewer)); 154 PetscCall(MatDestroy(&A)); 155 PetscFunctionReturn(0); 156 } 157 #endif 158 159 PETSC_INTERN PetscErrorCode PetscSFView_Basic(PetscSF sf, PetscViewer viewer) { 160 PetscBool isascii; 161 162 PetscFunctionBegin; 163 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 164 if (isascii && viewer->format != PETSC_VIEWER_ASCII_MATLAB) PetscCall(PetscViewerASCIIPrintf(viewer, " MultiSF sort=%s\n", sf->rankorder ? "rank-order" : "unordered")); 165 #if defined(PETSC_USE_SINGLE_LIBRARY) 166 else { 167 PetscBool isdraw, isbinary; 168 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 169 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 170 if ((isascii && viewer->format == PETSC_VIEWER_ASCII_MATLAB) || isdraw || isbinary) PetscCall(PetscSFView_Basic_PatternAndSizes(sf, viewer)); 171 } 172 #endif 173 PetscFunctionReturn(0); 174 } 175 176 static PetscErrorCode PetscSFBcastBegin_Basic(PetscSF sf, MPI_Datatype unit, PetscMemType rootmtype, const void *rootdata, PetscMemType leafmtype, void *leafdata, MPI_Op op) { 177 PetscSFLink link = NULL; 178 179 PetscFunctionBegin; 180 /* Create a communication link, which provides buffers, MPI requests etc (if MPI is used) */ 181 PetscCall(PetscSFLinkCreate(sf, unit, rootmtype, rootdata, leafmtype, leafdata, op, PETSCSF_BCAST, &link)); 182 /* Pack rootdata to rootbuf for remote communication */ 183 PetscCall(PetscSFLinkPackRootData(sf, link, PETSCSF_REMOTE, rootdata)); 184 /* Start communcation, e.g., post MPI_Isend */ 185 PetscCall(PetscSFLinkStartCommunication(sf, link, PETSCSF_ROOT2LEAF)); 186 /* Do local scatter (i.e., self to self communication), which overlaps with the remote communication above */ 187 PetscCall(PetscSFLinkScatterLocal(sf, link, PETSCSF_ROOT2LEAF, (void *)rootdata, leafdata, op)); 188 PetscFunctionReturn(0); 189 } 190 191 PETSC_INTERN PetscErrorCode PetscSFBcastEnd_Basic(PetscSF sf, MPI_Datatype unit, const void *rootdata, void *leafdata, MPI_Op op) { 192 PetscSFLink link = NULL; 193 194 PetscFunctionBegin; 195 /* Retrieve the link used in XxxBegin() with root/leafdata as key */ 196 PetscCall(PetscSFLinkGetInUse(sf, unit, rootdata, leafdata, PETSC_OWN_POINTER, &link)); 197 /* Finish remote communication, e.g., post MPI_Waitall */ 198 PetscCall(PetscSFLinkFinishCommunication(sf, link, PETSCSF_ROOT2LEAF)); 199 /* Unpack data in leafbuf to leafdata for remote communication */ 200 PetscCall(PetscSFLinkUnpackLeafData(sf, link, PETSCSF_REMOTE, leafdata, op)); 201 /* Recycle the link */ 202 PetscCall(PetscSFLinkReclaim(sf, &link)); 203 PetscFunctionReturn(0); 204 } 205 206 /* Shared by ReduceBegin and FetchAndOpBegin */ 207 static inline PetscErrorCode PetscSFLeafToRootBegin_Basic(PetscSF sf, MPI_Datatype unit, PetscMemType leafmtype, const void *leafdata, PetscMemType rootmtype, void *rootdata, MPI_Op op, PetscSFOperation sfop, PetscSFLink *out) { 208 PetscSFLink link = NULL; 209 210 PetscFunctionBegin; 211 PetscCall(PetscSFLinkCreate(sf, unit, rootmtype, rootdata, leafmtype, leafdata, op, sfop, &link)); 212 PetscCall(PetscSFLinkPackLeafData(sf, link, PETSCSF_REMOTE, leafdata)); 213 PetscCall(PetscSFLinkStartCommunication(sf, link, PETSCSF_LEAF2ROOT)); 214 *out = link; 215 PetscFunctionReturn(0); 216 } 217 218 /* leaf -> root with reduction */ 219 static PetscErrorCode PetscSFReduceBegin_Basic(PetscSF sf, MPI_Datatype unit, PetscMemType leafmtype, const void *leafdata, PetscMemType rootmtype, void *rootdata, MPI_Op op) { 220 PetscSFLink link = NULL; 221 222 PetscFunctionBegin; 223 PetscCall(PetscSFLeafToRootBegin_Basic(sf, unit, leafmtype, leafdata, rootmtype, rootdata, op, PETSCSF_REDUCE, &link)); 224 PetscCall(PetscSFLinkScatterLocal(sf, link, PETSCSF_LEAF2ROOT, rootdata, (void *)leafdata, op)); 225 PetscFunctionReturn(0); 226 } 227 228 PETSC_INTERN PetscErrorCode PetscSFReduceEnd_Basic(PetscSF sf, MPI_Datatype unit, const void *leafdata, void *rootdata, MPI_Op op) { 229 PetscSFLink link = NULL; 230 231 PetscFunctionBegin; 232 PetscCall(PetscSFLinkGetInUse(sf, unit, rootdata, leafdata, PETSC_OWN_POINTER, &link)); 233 PetscCall(PetscSFLinkFinishCommunication(sf, link, PETSCSF_LEAF2ROOT)); 234 PetscCall(PetscSFLinkUnpackRootData(sf, link, PETSCSF_REMOTE, rootdata, op)); 235 PetscCall(PetscSFLinkReclaim(sf, &link)); 236 PetscFunctionReturn(0); 237 } 238 239 PETSC_INTERN PetscErrorCode PetscSFFetchAndOpBegin_Basic(PetscSF sf, MPI_Datatype unit, PetscMemType rootmtype, void *rootdata, PetscMemType leafmtype, const void *leafdata, void *leafupdate, MPI_Op op) { 240 PetscSFLink link = NULL; 241 242 PetscFunctionBegin; 243 PetscCall(PetscSFLeafToRootBegin_Basic(sf, unit, leafmtype, leafdata, rootmtype, rootdata, op, PETSCSF_FETCH, &link)); 244 PetscCall(PetscSFLinkFetchAndOpLocal(sf, link, rootdata, leafdata, leafupdate, op)); 245 PetscFunctionReturn(0); 246 } 247 248 static PetscErrorCode PetscSFFetchAndOpEnd_Basic(PetscSF sf, MPI_Datatype unit, void *rootdata, const void *leafdata, void *leafupdate, MPI_Op op) { 249 PetscSFLink link = NULL; 250 251 PetscFunctionBegin; 252 PetscCall(PetscSFLinkGetInUse(sf, unit, rootdata, leafdata, PETSC_OWN_POINTER, &link)); 253 /* This implementation could be changed to unpack as receives arrive, at the cost of non-determinism */ 254 PetscCall(PetscSFLinkFinishCommunication(sf, link, PETSCSF_LEAF2ROOT)); 255 /* Do fetch-and-op, the (remote) update results are in rootbuf */ 256 PetscCall(PetscSFLinkFetchAndOpRemote(sf, link, rootdata, op)); 257 /* Bcast rootbuf to leafupdate */ 258 PetscCall(PetscSFLinkStartCommunication(sf, link, PETSCSF_ROOT2LEAF)); 259 PetscCall(PetscSFLinkFinishCommunication(sf, link, PETSCSF_ROOT2LEAF)); 260 /* Unpack and insert fetched data into leaves */ 261 PetscCall(PetscSFLinkUnpackLeafData(sf, link, PETSCSF_REMOTE, leafupdate, MPI_REPLACE)); 262 PetscCall(PetscSFLinkReclaim(sf, &link)); 263 PetscFunctionReturn(0); 264 } 265 266 PETSC_INTERN PetscErrorCode PetscSFGetLeafRanks_Basic(PetscSF sf, PetscInt *niranks, const PetscMPIInt **iranks, const PetscInt **ioffset, const PetscInt **irootloc) { 267 PetscSF_Basic *bas = (PetscSF_Basic *)sf->data; 268 269 PetscFunctionBegin; 270 if (niranks) *niranks = bas->niranks; 271 if (iranks) *iranks = bas->iranks; 272 if (ioffset) *ioffset = bas->ioffset; 273 if (irootloc) *irootloc = bas->irootloc; 274 PetscFunctionReturn(0); 275 } 276 277 /* An optimized PetscSFCreateEmbeddedRootSF. We aggresively make use of the established communication on sf. 278 We need one bcast on sf, and no communication anymore to build the embedded sf. Note that selected[] 279 was sorted before calling the routine. 280 */ 281 PETSC_INTERN PetscErrorCode PetscSFCreateEmbeddedRootSF_Basic(PetscSF sf, PetscInt nselected, const PetscInt *selected, PetscSF *newsf) { 282 PetscSF esf; 283 PetscInt esf_nranks, esf_ndranks, *esf_roffset, *esf_rmine, *esf_rremote; 284 PetscInt i, j, p, q, nroots, esf_nleaves, *new_ilocal, nranks, ndranks, niranks, ndiranks, minleaf, maxleaf, maxlocal; 285 char *rootdata, *leafdata, *leafmem; /* Only stores 0 or 1, so we can save memory with char */ 286 PetscMPIInt *esf_ranks; 287 const PetscMPIInt *ranks, *iranks; 288 const PetscInt *roffset, *rmine, *rremote, *ioffset, *irootloc; 289 PetscBool connected; 290 PetscSFNode *new_iremote; 291 PetscSF_Basic *bas; 292 293 PetscFunctionBegin; 294 PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)sf), &esf)); 295 PetscCall(PetscSFSetFromOptions(esf)); 296 PetscCall(PetscSFSetType(esf, PETSCSFBASIC)); /* This optimized routine can only create a basic sf */ 297 298 /* Find out which leaves are still connected to roots in the embedded sf by doing a Bcast */ 299 PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL)); 300 PetscCall(PetscSFGetLeafRange(sf, &minleaf, &maxleaf)); 301 maxlocal = maxleaf - minleaf + 1; 302 PetscCall(PetscCalloc2(nroots, &rootdata, maxlocal, &leafmem)); 303 leafdata = leafmem - minleaf; 304 /* Tag selected roots */ 305 for (i = 0; i < nselected; ++i) rootdata[selected[i]] = 1; 306 307 PetscCall(PetscSFBcastBegin(sf, MPI_CHAR, rootdata, leafdata, MPI_REPLACE)); 308 PetscCall(PetscSFBcastEnd(sf, MPI_CHAR, rootdata, leafdata, MPI_REPLACE)); 309 PetscCall(PetscSFGetLeafInfo_Basic(sf, &nranks, &ndranks, &ranks, &roffset, &rmine, &rremote)); /* Get send info */ 310 esf_nranks = esf_ndranks = esf_nleaves = 0; 311 for (i = 0; i < nranks; i++) { 312 connected = PETSC_FALSE; /* Is this process still connected to this remote root rank? */ 313 for (j = roffset[i]; j < roffset[i + 1]; j++) { 314 if (leafdata[rmine[j]]) { 315 esf_nleaves++; 316 connected = PETSC_TRUE; 317 } 318 } 319 if (connected) { 320 esf_nranks++; 321 if (i < ndranks) esf_ndranks++; 322 } 323 } 324 325 /* Set graph of esf and also set up its outgoing communication (i.e., send info), which is usually done by PetscSFSetUpRanks */ 326 PetscCall(PetscMalloc1(esf_nleaves, &new_ilocal)); 327 PetscCall(PetscMalloc1(esf_nleaves, &new_iremote)); 328 PetscCall(PetscMalloc4(esf_nranks, &esf_ranks, esf_nranks + 1, &esf_roffset, esf_nleaves, &esf_rmine, esf_nleaves, &esf_rremote)); 329 p = 0; /* Counter for connected root ranks */ 330 q = 0; /* Counter for connected leaves */ 331 esf_roffset[0] = 0; 332 for (i = 0; i < nranks; i++) { /* Scan leaf data again to fill esf arrays */ 333 connected = PETSC_FALSE; 334 for (j = roffset[i]; j < roffset[i + 1]; j++) { 335 if (leafdata[rmine[j]]) { 336 esf_rmine[q] = new_ilocal[q] = rmine[j]; 337 esf_rremote[q] = rremote[j]; 338 new_iremote[q].index = rremote[j]; 339 new_iremote[q].rank = ranks[i]; 340 connected = PETSC_TRUE; 341 q++; 342 } 343 } 344 if (connected) { 345 esf_ranks[p] = ranks[i]; 346 esf_roffset[p + 1] = q; 347 p++; 348 } 349 } 350 351 /* SetGraph internally resets the SF, so we only set its fields after the call */ 352 PetscCall(PetscSFSetGraph(esf, nroots, esf_nleaves, new_ilocal, PETSC_OWN_POINTER, new_iremote, PETSC_OWN_POINTER)); 353 esf->nranks = esf_nranks; 354 esf->ndranks = esf_ndranks; 355 esf->ranks = esf_ranks; 356 esf->roffset = esf_roffset; 357 esf->rmine = esf_rmine; 358 esf->rremote = esf_rremote; 359 esf->nleafreqs = esf_nranks - esf_ndranks; 360 361 /* Set up the incoming communication (i.e., recv info) stored in esf->data, which is usually done by PetscSFSetUp_Basic */ 362 bas = (PetscSF_Basic *)esf->data; 363 PetscCall(PetscSFGetRootInfo_Basic(sf, &niranks, &ndiranks, &iranks, &ioffset, &irootloc)); /* Get recv info */ 364 /* Embedded sf always has simpler communication than the original one. We might allocate longer arrays than needed here. But we 365 we do not care since these arrays are usually short. The benefit is we can fill these arrays by just parsing irootloc once. 366 */ 367 PetscCall(PetscMalloc2(niranks, &bas->iranks, niranks + 1, &bas->ioffset)); 368 PetscCall(PetscMalloc1(ioffset[niranks], &bas->irootloc)); 369 bas->niranks = bas->ndiranks = bas->ioffset[0] = 0; 370 p = 0; /* Counter for connected leaf ranks */ 371 q = 0; /* Counter for connected roots */ 372 for (i = 0; i < niranks; i++) { 373 connected = PETSC_FALSE; /* Is the current process still connected to this remote leaf rank? */ 374 for (j = ioffset[i]; j < ioffset[i + 1]; j++) { 375 if (rootdata[irootloc[j]]) { 376 bas->irootloc[q++] = irootloc[j]; 377 connected = PETSC_TRUE; 378 } 379 } 380 if (connected) { 381 bas->niranks++; 382 if (i < ndiranks) bas->ndiranks++; /* Note that order of ranks (including distinguished ranks) is kept */ 383 bas->iranks[p] = iranks[i]; 384 bas->ioffset[p + 1] = q; 385 p++; 386 } 387 } 388 bas->itotal = q; 389 bas->nrootreqs = bas->niranks - bas->ndiranks; 390 esf->persistent = PETSC_TRUE; 391 /* Setup packing related fields */ 392 PetscCall(PetscSFSetUpPackFields(esf)); 393 394 /* Copy from PetscSFSetUp(), since this method wants to skip PetscSFSetUp(). */ 395 #if defined(PETSC_HAVE_CUDA) 396 if (esf->backend == PETSCSF_BACKEND_CUDA) { 397 esf->ops->Malloc = PetscSFMalloc_CUDA; 398 esf->ops->Free = PetscSFFree_CUDA; 399 } 400 #endif 401 402 #if defined(PETSC_HAVE_HIP) 403 /* TODO: Needs debugging */ 404 if (esf->backend == PETSCSF_BACKEND_HIP) { 405 esf->ops->Malloc = PetscSFMalloc_HIP; 406 esf->ops->Free = PetscSFFree_HIP; 407 } 408 #endif 409 410 #if defined(PETSC_HAVE_KOKKOS) 411 if (esf->backend == PETSCSF_BACKEND_KOKKOS) { 412 esf->ops->Malloc = PetscSFMalloc_Kokkos; 413 esf->ops->Free = PetscSFFree_Kokkos; 414 } 415 #endif 416 esf->setupcalled = PETSC_TRUE; /* We have done setup ourselves! */ 417 PetscCall(PetscFree2(rootdata, leafmem)); 418 *newsf = esf; 419 PetscFunctionReturn(0); 420 } 421 422 PETSC_EXTERN PetscErrorCode PetscSFCreate_Basic(PetscSF sf) { 423 PetscSF_Basic *dat; 424 425 PetscFunctionBegin; 426 sf->ops->SetUp = PetscSFSetUp_Basic; 427 sf->ops->Reset = PetscSFReset_Basic; 428 sf->ops->Destroy = PetscSFDestroy_Basic; 429 sf->ops->View = PetscSFView_Basic; 430 sf->ops->BcastBegin = PetscSFBcastBegin_Basic; 431 sf->ops->BcastEnd = PetscSFBcastEnd_Basic; 432 sf->ops->ReduceBegin = PetscSFReduceBegin_Basic; 433 sf->ops->ReduceEnd = PetscSFReduceEnd_Basic; 434 sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Basic; 435 sf->ops->FetchAndOpEnd = PetscSFFetchAndOpEnd_Basic; 436 sf->ops->GetLeafRanks = PetscSFGetLeafRanks_Basic; 437 sf->ops->CreateEmbeddedRootSF = PetscSFCreateEmbeddedRootSF_Basic; 438 439 PetscCall(PetscNewLog(sf, &dat)); 440 sf->data = (void *)dat; 441 PetscFunctionReturn(0); 442 } 443