1b0c7db22SLisandro Dalcin #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/ 2b0c7db22SLisandro Dalcin #include <petsc/private/sectionimpl.h> 3b0c7db22SLisandro Dalcin 4eb58282aSVaclav Hapla /*@ 52f04c522SBarry Smith PetscSFSetGraphLayout - Set a `PetscSF` communication pattern using global indices and a `PetscLayout` 6b0c7db22SLisandro Dalcin 7b0c7db22SLisandro Dalcin Collective 8b0c7db22SLisandro Dalcin 94165533cSJose E. Roman Input Parameters: 10b0c7db22SLisandro Dalcin + sf - star forest 112f04c522SBarry Smith . layout - `PetscLayout` defining the global space for roots, i.e. which roots are owned by each MPI process 122f04c522SBarry Smith . nleaves - number of leaf vertices on the current process, each of these references a root on any MPI process 132f04c522SBarry Smith . ilocal - locations of leaves in leafdata buffers, pass `NULL` for contiguous storage, that is the locations are in [0,`nleaves`) 142f04c522SBarry Smith . localmode - copy mode for `ilocal` 152f04c522SBarry Smith - gremote - root vertices in global numbering corresponding to the leaves 16b0c7db22SLisandro Dalcin 17b0c7db22SLisandro Dalcin Level: intermediate 18b0c7db22SLisandro Dalcin 19cab54364SBarry Smith Note: 202f04c522SBarry Smith Global indices must lie in [0, N) where N is the global size of `layout`. 212f04c522SBarry Smith Leaf indices in `ilocal` get sorted; this means the user-provided array gets sorted if localmode is `PETSC_OWN_POINTER`. 2286c56f52SVaclav Hapla 2338b5cf2dSJacob Faibussowitsch Developer Notes: 242f04c522SBarry Smith Local indices which are the identity permutation in the range [0,`nleaves`) are discarded as they 25cab54364SBarry Smith encode contiguous storage. In such case, if localmode is `PETSC_OWN_POINTER`, the memory is deallocated as it is not 26b0c7db22SLisandro Dalcin needed 27b0c7db22SLisandro Dalcin 282f04c522SBarry Smith .seealso: [](sec_petscsf), `PetscSF`, `PetscSFGetGraphLayout()`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()` 29b0c7db22SLisandro Dalcin @*/ 30cf84bf9eSBarry Smith PetscErrorCode PetscSFSetGraphLayout(PetscSF sf, PetscLayout layout, PetscInt nleaves, PetscInt ilocal[], PetscCopyMode localmode, const PetscInt gremote[]) 31d71ae5a4SJacob Faibussowitsch { 3238a25198SStefano Zampini const PetscInt *range; 3338a25198SStefano Zampini PetscInt i, nroots, ls = -1, ln = -1; 3438a25198SStefano Zampini PetscMPIInt lr = -1; 35b0c7db22SLisandro Dalcin PetscSFNode *remote; 36b0c7db22SLisandro Dalcin 37b0c7db22SLisandro Dalcin PetscFunctionBegin; 38da0802e2SStefano Zampini PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 39*18aa8208SJames Wright PetscAssertPointer(layout, 2); 40*18aa8208SJames Wright if (nleaves > 0 && ilocal) PetscAssertPointer(ilocal, 4); 41*18aa8208SJames Wright if (nleaves > 0) PetscAssertPointer(gremote, 6); 42b114984aSVaclav Hapla PetscCall(PetscLayoutSetUp(layout)); 439566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(layout, &nroots)); 449566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRanges(layout, &range)); 459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &remote)); 46eb58282aSVaclav Hapla if (nleaves) ls = gremote[0] + 1; 47b0c7db22SLisandro Dalcin for (i = 0; i < nleaves; i++) { 48eb58282aSVaclav Hapla const PetscInt idx = gremote[i] - ls; 4938a25198SStefano Zampini if (idx < 0 || idx >= ln) { /* short-circuit the search */ 50eb58282aSVaclav Hapla PetscCall(PetscLayoutFindOwnerIndex(layout, gremote[i], &lr, &remote[i].index)); 5138a25198SStefano Zampini remote[i].rank = lr; 5238a25198SStefano Zampini ls = range[lr]; 5338a25198SStefano Zampini ln = range[lr + 1] - ls; 5438a25198SStefano Zampini } else { 5538a25198SStefano Zampini remote[i].rank = lr; 5638a25198SStefano Zampini remote[i].index = idx; 5738a25198SStefano Zampini } 58b0c7db22SLisandro Dalcin } 599566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, nroots, nleaves, ilocal, localmode, remote, PETSC_OWN_POINTER)); 603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 61b0c7db22SLisandro Dalcin } 62b0c7db22SLisandro Dalcin 63eb58282aSVaclav Hapla /*@C 642f04c522SBarry Smith PetscSFGetGraphLayout - Get the global indices and `PetscLayout` that describe a `PetscSF` 65eb58282aSVaclav Hapla 66eb58282aSVaclav Hapla Collective 67eb58282aSVaclav Hapla 68eb58282aSVaclav Hapla Input Parameter: 69eb58282aSVaclav Hapla . sf - star forest 70eb58282aSVaclav Hapla 71eb58282aSVaclav Hapla Output Parameters: 72cab54364SBarry Smith + layout - `PetscLayout` defining the global space for roots 73eb58282aSVaclav Hapla . nleaves - number of leaf vertices on the current process, each of these references a root on any process 74f13dfd9eSBarry Smith . ilocal - locations of leaves in leafdata buffers, or `NULL` for contiguous storage 752f04c522SBarry Smith - gremote - root vertices in global numbering corresponding to the leaves 76eb58282aSVaclav Hapla 77eb58282aSVaclav Hapla Level: intermediate 78eb58282aSVaclav Hapla 79eb58282aSVaclav Hapla Notes: 80eb58282aSVaclav Hapla The outputs are such that passing them as inputs to `PetscSFSetGraphLayout()` would lead to the same star forest. 81f13dfd9eSBarry Smith The outputs `layout` and `gremote` are freshly created each time this function is called, 82f13dfd9eSBarry Smith so they need to be freed (with `PetscLayoutDestroy()` and `PetscFree()`) by the user. 83eb58282aSVaclav Hapla 842f04c522SBarry Smith .seealso: [](sec_petscsf), `PetscSF`, `PetscSFSetGraphLayout()`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()` 85eb58282aSVaclav Hapla @*/ 86d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetGraphLayout(PetscSF sf, PetscLayout *layout, PetscInt *nleaves, const PetscInt *ilocal[], PetscInt *gremote[]) 87d71ae5a4SJacob Faibussowitsch { 88eb58282aSVaclav Hapla PetscInt nr, nl; 89eb58282aSVaclav Hapla const PetscSFNode *ir; 90eb58282aSVaclav Hapla PetscLayout lt; 91eb58282aSVaclav Hapla 92eb58282aSVaclav Hapla PetscFunctionBegin; 93*18aa8208SJames Wright PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 94*18aa8208SJames Wright if (layout) PetscAssertPointer(layout, 2); 95*18aa8208SJames Wright if (nleaves) PetscAssertPointer(nleaves, 3); 96*18aa8208SJames Wright if (ilocal) PetscAssertPointer(ilocal, 4); 97*18aa8208SJames Wright if (gremote) PetscAssertPointer(gremote, 5); 98eb58282aSVaclav Hapla PetscCall(PetscSFGetGraph(sf, &nr, &nl, ilocal, &ir)); 99eb58282aSVaclav Hapla PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)sf), nr, PETSC_DECIDE, 1, <)); 100eb58282aSVaclav Hapla if (gremote) { 101eb58282aSVaclav Hapla PetscInt i; 102eb58282aSVaclav Hapla const PetscInt *range; 103eb58282aSVaclav Hapla PetscInt *gr; 104eb58282aSVaclav Hapla 105eb58282aSVaclav Hapla PetscCall(PetscLayoutGetRanges(lt, &range)); 106eb58282aSVaclav Hapla PetscCall(PetscMalloc1(nl, &gr)); 107eb58282aSVaclav Hapla for (i = 0; i < nl; i++) gr[i] = range[ir[i].rank] + ir[i].index; 108eb58282aSVaclav Hapla *gremote = gr; 109eb58282aSVaclav Hapla } 110eb58282aSVaclav Hapla if (nleaves) *nleaves = nl; 111eb58282aSVaclav Hapla if (layout) *layout = lt; 112eb58282aSVaclav Hapla else PetscCall(PetscLayoutDestroy(<)); 1133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 114eb58282aSVaclav Hapla } 115eb58282aSVaclav Hapla 116b0c7db22SLisandro Dalcin /*@ 1172f04c522SBarry Smith PetscSFSetGraphSection - Sets the `PetscSF` graph (communication pattern) encoding the parallel dof overlap based upon the `PetscSection` describing the data layout. 118b0c7db22SLisandro Dalcin 119b0c7db22SLisandro Dalcin Input Parameters: 120cab54364SBarry Smith + sf - The `PetscSF` 121cab54364SBarry Smith . localSection - `PetscSection` describing the local data layout 122cab54364SBarry Smith - globalSection - `PetscSection` describing the global data layout 123b0c7db22SLisandro Dalcin 124b0c7db22SLisandro Dalcin Level: developer 125b0c7db22SLisandro Dalcin 1262f04c522SBarry Smith .seealso: [](sec_petscsf), `PetscSF`, `PetscSFSetGraph()`, `PetscSFSetGraphLayout()` 127b0c7db22SLisandro Dalcin @*/ 128d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetGraphSection(PetscSF sf, PetscSection localSection, PetscSection globalSection) 129d71ae5a4SJacob Faibussowitsch { 130b0c7db22SLisandro Dalcin MPI_Comm comm; 131b0c7db22SLisandro Dalcin PetscLayout layout; 132b0c7db22SLisandro Dalcin const PetscInt *ranges; 133b0c7db22SLisandro Dalcin PetscInt *local; 134b0c7db22SLisandro Dalcin PetscSFNode *remote; 135b0c7db22SLisandro Dalcin PetscInt pStart, pEnd, p, nroots, nleaves = 0, l; 136b0c7db22SLisandro Dalcin PetscMPIInt size, rank; 137b0c7db22SLisandro Dalcin 138b0c7db22SLisandro Dalcin PetscFunctionBegin; 139b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 140b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(localSection, PETSC_SECTION_CLASSID, 2); 141b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(globalSection, PETSC_SECTION_CLASSID, 3); 142b0c7db22SLisandro Dalcin 1439566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 1449566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 1459566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1469566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(globalSection, &pStart, &pEnd)); 1479566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstrainedStorageSize(globalSection, &nroots)); 1489566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(comm, &layout)); 1499566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(layout, 1)); 1509566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetLocalSize(layout, nroots)); 1519566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(layout)); 1529566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRanges(layout, &ranges)); 153b0c7db22SLisandro Dalcin for (p = pStart; p < pEnd; ++p) { 154b0c7db22SLisandro Dalcin PetscInt gdof, gcdof; 155b0c7db22SLisandro Dalcin 1569566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(globalSection, p, &gdof)); 1579566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof)); 15857508eceSPierre Jolivet 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); 159b0c7db22SLisandro Dalcin nleaves += gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof; 160b0c7db22SLisandro Dalcin } 1619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &local)); 1629566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &remote)); 163b0c7db22SLisandro Dalcin for (p = pStart, l = 0; p < pEnd; ++p) { 164b0c7db22SLisandro Dalcin const PetscInt *cind; 165b0c7db22SLisandro Dalcin PetscInt dof, cdof, off, gdof, gcdof, goff, gsize, d, c; 166b0c7db22SLisandro Dalcin 1679566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(localSection, p, &dof)); 1689566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(localSection, p, &off)); 1699566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(localSection, p, &cdof)); 1709566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintIndices(localSection, p, &cind)); 1719566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(globalSection, p, &gdof)); 1729566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof)); 1739566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(globalSection, p, &goff)); 174b0c7db22SLisandro Dalcin if (!gdof) continue; /* Censored point */ 175b0c7db22SLisandro Dalcin gsize = gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof; 176b0c7db22SLisandro Dalcin if (gsize != dof - cdof) { 17708401ef6SPierre Jolivet 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); 178b0c7db22SLisandro Dalcin cdof = 0; /* Ignore constraints */ 179b0c7db22SLisandro Dalcin } 180b0c7db22SLisandro Dalcin for (d = 0, c = 0; d < dof; ++d) { 1819371c9d4SSatish Balay if ((c < cdof) && (cind[c] == d)) { 1829371c9d4SSatish Balay ++c; 1839371c9d4SSatish Balay continue; 1849371c9d4SSatish Balay } 185b0c7db22SLisandro Dalcin local[l + d - c] = off + d; 186b0c7db22SLisandro Dalcin } 18708401ef6SPierre Jolivet 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); 188b0c7db22SLisandro Dalcin if (gdof < 0) { 189b0c7db22SLisandro Dalcin for (d = 0; d < gsize; ++d, ++l) { 1906497c311SBarry Smith PetscInt offset = -(goff + 1) + d, ir; 1916497c311SBarry Smith PetscMPIInt r; 192b0c7db22SLisandro Dalcin 1936497c311SBarry Smith PetscCall(PetscFindInt(offset, size + 1, ranges, &ir)); 1946497c311SBarry Smith PetscCall(PetscMPIIntCast(ir, &r)); 195b0c7db22SLisandro Dalcin if (r < 0) r = -(r + 2); 1966497c311SBarry Smith PetscCheck(!(r < 0) && !(r >= size), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " mapped to invalid process %d (%" PetscInt_FMT ", %" PetscInt_FMT ")", p, r, gdof, goff); 197b0c7db22SLisandro Dalcin remote[l].rank = r; 198b0c7db22SLisandro Dalcin remote[l].index = offset - ranges[r]; 199b0c7db22SLisandro Dalcin } 200b0c7db22SLisandro Dalcin } else { 201b0c7db22SLisandro Dalcin for (d = 0; d < gsize; ++d, ++l) { 202b0c7db22SLisandro Dalcin remote[l].rank = rank; 203b0c7db22SLisandro Dalcin remote[l].index = goff + d - ranges[rank]; 204b0c7db22SLisandro Dalcin } 205b0c7db22SLisandro Dalcin } 206b0c7db22SLisandro Dalcin } 20708401ef6SPierre Jolivet PetscCheck(l == nleaves, comm, PETSC_ERR_PLIB, "Iteration error, l %" PetscInt_FMT " != nleaves %" PetscInt_FMT, l, nleaves); 2089566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&layout)); 2099566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER)); 2103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 211b0c7db22SLisandro Dalcin } 212b0c7db22SLisandro Dalcin 213b0c7db22SLisandro Dalcin /*@C 214cab54364SBarry Smith PetscSFDistributeSection - Create a new `PetscSection` reorganized, moving from the root to the leaves of the `PetscSF` 215b0c7db22SLisandro Dalcin 216c3339decSBarry Smith Collective 217b0c7db22SLisandro Dalcin 218b0c7db22SLisandro Dalcin Input Parameters: 219cab54364SBarry Smith + sf - The `PetscSF` 220b0c7db22SLisandro Dalcin - rootSection - Section defined on root space 221b0c7db22SLisandro Dalcin 222b0c7db22SLisandro Dalcin Output Parameters: 223ce78bad3SBarry Smith + remoteOffsets - root offsets in leaf storage, or `NULL`, its length will be the size of the chart of `leafSection` 224b0c7db22SLisandro Dalcin - leafSection - Section defined on the leaf space 225b0c7db22SLisandro Dalcin 226b0c7db22SLisandro Dalcin Level: advanced 227b0c7db22SLisandro Dalcin 228ce78bad3SBarry Smith Note: 229ce78bad3SBarry Smith Caller must `PetscFree()` `remoteOffsets` if it was requested 230ce78bad3SBarry Smith 231ce78bad3SBarry Smith Fortran Note: 232ce78bad3SBarry Smith Use `PetscSFDestroyRemoteOffsets()` when `remoteOffsets` is no longer needed. 23323e9620eSJunchao Zhang 2342f04c522SBarry Smith .seealso: [](sec_petscsf), `PetscSF`, `PetscSFCreate()` 235b0c7db22SLisandro Dalcin @*/ 236ce78bad3SBarry Smith PetscErrorCode PetscSFDistributeSection(PetscSF sf, PetscSection rootSection, PetscInt *remoteOffsets[], PetscSection leafSection) 237d71ae5a4SJacob Faibussowitsch { 238b0c7db22SLisandro Dalcin PetscSF embedSF; 239b0c7db22SLisandro Dalcin const PetscInt *indices; 240b0c7db22SLisandro Dalcin IS selected; 2411690c2aeSBarry Smith PetscInt numFields, nroots, rpStart, rpEnd, lpStart = PETSC_INT_MAX, lpEnd = -1, f, c; 242b0c7db22SLisandro Dalcin PetscBool *sub, hasc; 243b0c7db22SLisandro Dalcin 244b0c7db22SLisandro Dalcin PetscFunctionBegin; 245*18aa8208SJames Wright PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 246*18aa8208SJames Wright PetscValidHeaderSpecific(rootSection, PETSC_SECTION_CLASSID, 2); 247*18aa8208SJames Wright if (remoteOffsets) PetscAssertPointer(remoteOffsets, 3); 248*18aa8208SJames Wright PetscValidHeaderSpecific(leafSection, PETSC_SECTION_CLASSID, 4); 2499566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_DistSect, sf, 0, 0, 0)); 2509566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(rootSection, &numFields)); 251029e8818Sksagiyam if (numFields) { 252029e8818Sksagiyam IS perm; 253029e8818Sksagiyam 254029e8818Sksagiyam /* PetscSectionSetNumFields() calls PetscSectionReset(), which destroys 255029e8818Sksagiyam leafSection->perm. To keep this permutation set by the user, we grab 256029e8818Sksagiyam the reference before calling PetscSectionSetNumFields() and set it 257029e8818Sksagiyam back after. */ 2589566063dSJacob Faibussowitsch PetscCall(PetscSectionGetPermutation(leafSection, &perm)); 2599566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)perm)); 2609566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(leafSection, numFields)); 2619566063dSJacob Faibussowitsch PetscCall(PetscSectionSetPermutation(leafSection, perm)); 2629566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 263029e8818Sksagiyam } 2649566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numFields + 2, &sub)); 265b0c7db22SLisandro Dalcin sub[1] = rootSection->bc ? PETSC_TRUE : PETSC_FALSE; 266b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) { 2672ee82556SMatthew G. Knepley PetscSectionSym sym, dsym = NULL; 268b0c7db22SLisandro Dalcin const char *name = NULL; 269b0c7db22SLisandro Dalcin PetscInt numComp = 0; 270b0c7db22SLisandro Dalcin 271b0c7db22SLisandro Dalcin sub[2 + f] = rootSection->field[f]->bc ? PETSC_TRUE : PETSC_FALSE; 2729566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldComponents(rootSection, f, &numComp)); 2739566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldName(rootSection, f, &name)); 2749566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldSym(rootSection, f, &sym)); 2759566063dSJacob Faibussowitsch if (sym) PetscCall(PetscSectionSymDistribute(sym, sf, &dsym)); 2769566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(leafSection, f, numComp)); 2779566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldName(leafSection, f, name)); 2789566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldSym(leafSection, f, dsym)); 2799566063dSJacob Faibussowitsch PetscCall(PetscSectionSymDestroy(&dsym)); 280b0c7db22SLisandro Dalcin for (c = 0; c < rootSection->numFieldComponents[f]; ++c) { 2819566063dSJacob Faibussowitsch PetscCall(PetscSectionGetComponentName(rootSection, f, c, &name)); 2829566063dSJacob Faibussowitsch PetscCall(PetscSectionSetComponentName(leafSection, f, c, name)); 283b0c7db22SLisandro Dalcin } 284b0c7db22SLisandro Dalcin } 2859566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd)); 2869566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL)); 287b0c7db22SLisandro Dalcin rpEnd = PetscMin(rpEnd, nroots); 288b0c7db22SLisandro Dalcin rpEnd = PetscMax(rpStart, rpEnd); 289b0c7db22SLisandro Dalcin /* see if we can avoid creating the embedded SF, since it can cost more than an allreduce */ 290b0c7db22SLisandro Dalcin sub[0] = (PetscBool)(nroots != rpEnd - rpStart); 2915440e5dcSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, sub, 2 + numFields, MPI_C_BOOL, MPI_LOR, PetscObjectComm((PetscObject)sf))); 292b0c7db22SLisandro Dalcin if (sub[0]) { 2939566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected)); 2949566063dSJacob Faibussowitsch PetscCall(ISGetIndices(selected, &indices)); 2959566063dSJacob Faibussowitsch PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF)); 2969566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(selected, &indices)); 2979566063dSJacob Faibussowitsch PetscCall(ISDestroy(&selected)); 298b0c7db22SLisandro Dalcin } else { 2999566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)sf)); 300b0c7db22SLisandro Dalcin embedSF = sf; 301b0c7db22SLisandro Dalcin } 3029566063dSJacob Faibussowitsch PetscCall(PetscSFGetLeafRange(embedSF, &lpStart, &lpEnd)); 303b0c7db22SLisandro Dalcin lpEnd++; 304b0c7db22SLisandro Dalcin 3059566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(leafSection, lpStart, lpEnd)); 306b0c7db22SLisandro Dalcin 307b0c7db22SLisandro Dalcin /* Constrained dof section */ 308b0c7db22SLisandro Dalcin hasc = sub[1]; 309b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) hasc = (PetscBool)(hasc || sub[2 + f]); 310b0c7db22SLisandro Dalcin 311b0c7db22SLisandro Dalcin /* Could fuse these at the cost of copies and extra allocation */ 31216cd844bSPierre Jolivet PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasDof, -rpStart), PetscSafePointerPlusOffset(leafSection->atlasDof, -lpStart), MPI_REPLACE)); 31316cd844bSPierre Jolivet PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasDof, -rpStart), PetscSafePointerPlusOffset(leafSection->atlasDof, -lpStart), MPI_REPLACE)); 314b0c7db22SLisandro Dalcin if (sub[1]) { 3157c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(rootSection)); 3167c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(leafSection)); 3179566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE)); 3189566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE)); 319b0c7db22SLisandro Dalcin } 320b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) { 32116cd844bSPierre Jolivet PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->field[f]->atlasDof, -rpStart), PetscSafePointerPlusOffset(leafSection->field[f]->atlasDof, -lpStart), MPI_REPLACE)); 32216cd844bSPierre Jolivet PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->field[f]->atlasDof, -rpStart), PetscSafePointerPlusOffset(leafSection->field[f]->atlasDof, -lpStart), MPI_REPLACE)); 323b0c7db22SLisandro Dalcin if (sub[2 + f]) { 3247c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(rootSection->field[f])); 3257c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(leafSection->field[f])); 3269566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE)); 3279566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE)); 328b0c7db22SLisandro Dalcin } 329b0c7db22SLisandro Dalcin } 330b0c7db22SLisandro Dalcin if (remoteOffsets) { 3319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lpEnd - lpStart, remoteOffsets)); 33216cd844bSPierre Jolivet PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(*remoteOffsets, -lpStart), MPI_REPLACE)); 33316cd844bSPierre Jolivet PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(*remoteOffsets, -lpStart), MPI_REPLACE)); 334b0c7db22SLisandro Dalcin } 33569c11d05SVaclav Hapla PetscCall(PetscSectionInvalidateMaxDof_Internal(leafSection)); 3369566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(leafSection)); 337b0c7db22SLisandro Dalcin if (hasc) { /* need to communicate bcIndices */ 338b0c7db22SLisandro Dalcin PetscSF bcSF; 339b0c7db22SLisandro Dalcin PetscInt *rOffBc; 340b0c7db22SLisandro Dalcin 3419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lpEnd - lpStart, &rOffBc)); 342b0c7db22SLisandro Dalcin if (sub[1]) { 3439566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 3449566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 3459566063dSJacob Faibussowitsch PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->bc, rOffBc, leafSection->bc, &bcSF)); 3469566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE)); 3479566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE)); 3489566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&bcSF)); 349b0c7db22SLisandro Dalcin } 350b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) { 351b0c7db22SLisandro Dalcin if (sub[2 + f]) { 3529566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 3539566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 3549566063dSJacob Faibussowitsch PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->field[f]->bc, rOffBc, leafSection->field[f]->bc, &bcSF)); 3559566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE)); 3569566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE)); 3579566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&bcSF)); 358b0c7db22SLisandro Dalcin } 359b0c7db22SLisandro Dalcin } 3609566063dSJacob Faibussowitsch PetscCall(PetscFree(rOffBc)); 361b0c7db22SLisandro Dalcin } 3629566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&embedSF)); 3639566063dSJacob Faibussowitsch PetscCall(PetscFree(sub)); 3649566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_DistSect, sf, 0, 0, 0)); 3653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 366b0c7db22SLisandro Dalcin } 367b0c7db22SLisandro Dalcin 368b0c7db22SLisandro Dalcin /*@C 369b0c7db22SLisandro Dalcin PetscSFCreateRemoteOffsets - Create offsets for point data on remote processes 370b0c7db22SLisandro Dalcin 371c3339decSBarry Smith Collective 372b0c7db22SLisandro Dalcin 373b0c7db22SLisandro Dalcin Input Parameters: 374cab54364SBarry Smith + sf - The `PetscSF` 375ce78bad3SBarry Smith . rootSection - Data layout of remote points for outgoing data (this is layout for roots) 376ce78bad3SBarry Smith - leafSection - Data layout of local points for incoming data (this is layout for leaves) 377b0c7db22SLisandro Dalcin 378b0c7db22SLisandro Dalcin Output Parameter: 379ce78bad3SBarry Smith . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or `NULL` 380b0c7db22SLisandro Dalcin 381b0c7db22SLisandro Dalcin Level: developer 382b0c7db22SLisandro Dalcin 383ce78bad3SBarry Smith Note: 384ce78bad3SBarry Smith Caller must `PetscFree()` `remoteOffsets` if it was requested 385ce78bad3SBarry Smith 386ce78bad3SBarry Smith Fortran Note: 387ce78bad3SBarry Smith Use `PetscSFDestroyRemoteOffsets()` when `remoteOffsets` is no longer needed. 38823e9620eSJunchao Zhang 3892f04c522SBarry Smith .seealso: [](sec_petscsf), `PetscSF`, `PetscSFCreate()` 390b0c7db22SLisandro Dalcin @*/ 391ce78bad3SBarry Smith PetscErrorCode PetscSFCreateRemoteOffsets(PetscSF sf, PetscSection rootSection, PetscSection leafSection, PetscInt *remoteOffsets[]) 392d71ae5a4SJacob Faibussowitsch { 393b0c7db22SLisandro Dalcin PetscSF embedSF; 394b0c7db22SLisandro Dalcin const PetscInt *indices; 395b0c7db22SLisandro Dalcin IS selected; 396b0c7db22SLisandro Dalcin PetscInt numRoots, rpStart = 0, rpEnd = 0, lpStart = 0, lpEnd = 0; 397b0c7db22SLisandro Dalcin 398b0c7db22SLisandro Dalcin PetscFunctionBegin; 399*18aa8208SJames Wright PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 400*18aa8208SJames Wright PetscValidHeaderSpecific(rootSection, PETSC_SECTION_CLASSID, 2); 401*18aa8208SJames Wright PetscValidHeaderSpecific(leafSection, PETSC_SECTION_CLASSID, 3); 402*18aa8208SJames Wright PetscAssertPointer(remoteOffsets, 4); 403b0c7db22SLisandro Dalcin *remoteOffsets = NULL; 4049566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &numRoots, NULL, NULL, NULL)); 4053ba16761SJacob Faibussowitsch if (numRoots < 0) PetscFunctionReturn(PETSC_SUCCESS); 4069566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_RemoteOff, sf, 0, 0, 0)); 4079566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd)); 4089566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd)); 4099566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected)); 4109566063dSJacob Faibussowitsch PetscCall(ISGetIndices(selected, &indices)); 4119566063dSJacob Faibussowitsch PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF)); 4129566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(selected, &indices)); 4139566063dSJacob Faibussowitsch PetscCall(ISDestroy(&selected)); 4149566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(lpEnd - lpStart, remoteOffsets)); 4158e3a54c0SPierre Jolivet PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(*remoteOffsets, -lpStart), MPI_REPLACE)); 4168e3a54c0SPierre Jolivet PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(*remoteOffsets, -lpStart), MPI_REPLACE)); 4179566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&embedSF)); 4189566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_RemoteOff, sf, 0, 0, 0)); 4193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 420b0c7db22SLisandro Dalcin } 421b0c7db22SLisandro Dalcin 422ce78bad3SBarry Smith /*@ 423cab54364SBarry Smith PetscSFCreateSectionSF - Create an expanded `PetscSF` of dofs, assuming the input `PetscSF` relates points 424b0c7db22SLisandro Dalcin 425c3339decSBarry Smith Collective 426b0c7db22SLisandro Dalcin 427b0c7db22SLisandro Dalcin Input Parameters: 428cab54364SBarry Smith + sf - The `PetscSF` 429b0c7db22SLisandro Dalcin . rootSection - Data layout of remote points for outgoing data (this is usually the serial section) 4302f04c522SBarry Smith . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or `NULL` 431b0c7db22SLisandro Dalcin - leafSection - Data layout of local points for incoming data (this is the distributed section) 432b0c7db22SLisandro Dalcin 4332fe279fdSBarry Smith Output Parameter: 4342fe279fdSBarry Smith . sectionSF - The new `PetscSF` 435b0c7db22SLisandro Dalcin 436b0c7db22SLisandro Dalcin Level: advanced 437b0c7db22SLisandro Dalcin 43823e9620eSJunchao Zhang Notes: 439ce78bad3SBarry Smith `remoteOffsets` can be `NULL` if `sf` does not reference any points in leafSection 44023e9620eSJunchao Zhang 4412f04c522SBarry Smith .seealso: [](sec_petscsf), `PetscSF`, `PetscSFCreate()` 442b0c7db22SLisandro Dalcin @*/ 443d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateSectionSF(PetscSF sf, PetscSection rootSection, PetscInt remoteOffsets[], PetscSection leafSection, PetscSF *sectionSF) 444d71ae5a4SJacob Faibussowitsch { 445b0c7db22SLisandro Dalcin MPI_Comm comm; 446b0c7db22SLisandro Dalcin const PetscInt *localPoints; 447b0c7db22SLisandro Dalcin const PetscSFNode *remotePoints; 448b0c7db22SLisandro Dalcin PetscInt lpStart, lpEnd; 449b0c7db22SLisandro Dalcin PetscInt numRoots, numSectionRoots, numPoints, numIndices = 0; 450b0c7db22SLisandro Dalcin PetscInt *localIndices; 451b0c7db22SLisandro Dalcin PetscSFNode *remoteIndices; 452b0c7db22SLisandro Dalcin PetscInt i, ind; 453b0c7db22SLisandro Dalcin 454b0c7db22SLisandro Dalcin PetscFunctionBegin; 455b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 4564f572ea9SToby Isaac PetscAssertPointer(rootSection, 2); 4574f572ea9SToby Isaac /* Cannot check PetscAssertPointer(remoteOffsets,3) because it can be NULL if sf does not reference any points in leafSection */ 4584f572ea9SToby Isaac PetscAssertPointer(leafSection, 4); 4594f572ea9SToby Isaac PetscAssertPointer(sectionSF, 5); 4609566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 4619566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, sectionSF)); 4629566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd)); 4639566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(rootSection, &numSectionRoots)); 4649566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &numRoots, &numPoints, &localPoints, &remotePoints)); 4653ba16761SJacob Faibussowitsch if (numRoots < 0) PetscFunctionReturn(PETSC_SUCCESS); 4669566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_SectSF, sf, 0, 0, 0)); 467b0c7db22SLisandro Dalcin for (i = 0; i < numPoints; ++i) { 468b0c7db22SLisandro Dalcin PetscInt localPoint = localPoints ? localPoints[i] : i; 469b0c7db22SLisandro Dalcin PetscInt dof; 470b0c7db22SLisandro Dalcin 471b0c7db22SLisandro Dalcin if ((localPoint >= lpStart) && (localPoint < lpEnd)) { 4729566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof)); 47322a18c9fSMatthew G. Knepley numIndices += dof < 0 ? 0 : dof; 474b0c7db22SLisandro Dalcin } 475b0c7db22SLisandro Dalcin } 4769566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numIndices, &localIndices)); 4779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numIndices, &remoteIndices)); 478b0c7db22SLisandro Dalcin /* Create new index graph */ 479b0c7db22SLisandro Dalcin for (i = 0, ind = 0; i < numPoints; ++i) { 480b0c7db22SLisandro Dalcin PetscInt localPoint = localPoints ? localPoints[i] : i; 481835f2295SStefano Zampini PetscInt rank = remotePoints[i].rank; 482b0c7db22SLisandro Dalcin 483b0c7db22SLisandro Dalcin if ((localPoint >= lpStart) && (localPoint < lpEnd)) { 484b0c7db22SLisandro Dalcin PetscInt remoteOffset = remoteOffsets[localPoint - lpStart]; 485b0c7db22SLisandro Dalcin PetscInt loff, dof, d; 486b0c7db22SLisandro Dalcin 4879566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(leafSection, localPoint, &loff)); 4889566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof)); 489b0c7db22SLisandro Dalcin for (d = 0; d < dof; ++d, ++ind) { 490b0c7db22SLisandro Dalcin localIndices[ind] = loff + d; 491b0c7db22SLisandro Dalcin remoteIndices[ind].rank = rank; 492b0c7db22SLisandro Dalcin remoteIndices[ind].index = remoteOffset + d; 493b0c7db22SLisandro Dalcin } 494b0c7db22SLisandro Dalcin } 495b0c7db22SLisandro Dalcin } 49608401ef6SPierre Jolivet PetscCheck(numIndices == ind, comm, PETSC_ERR_PLIB, "Inconsistency in indices, %" PetscInt_FMT " should be %" PetscInt_FMT, ind, numIndices); 4979566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*sectionSF, numSectionRoots, numIndices, localIndices, PETSC_OWN_POINTER, remoteIndices, PETSC_OWN_POINTER)); 4989566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(*sectionSF)); 4999566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_SectSF, sf, 0, 0, 0)); 5003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 501b0c7db22SLisandro Dalcin } 5028fa9e22eSVaclav Hapla 5038fa9e22eSVaclav Hapla /*@C 5042f04c522SBarry Smith PetscSFCreateFromLayouts - Creates a parallel star forest mapping between two `PetscLayout` objects 5058fa9e22eSVaclav Hapla 5068fa9e22eSVaclav Hapla Collective 5078fa9e22eSVaclav Hapla 5084165533cSJose E. Roman Input Parameters: 509cab54364SBarry Smith + rmap - `PetscLayout` defining the global root space 510cab54364SBarry Smith - lmap - `PetscLayout` defining the global leaf space 5118fa9e22eSVaclav Hapla 5124165533cSJose E. Roman Output Parameter: 5138fa9e22eSVaclav Hapla . sf - The parallel star forest 5148fa9e22eSVaclav Hapla 5158fa9e22eSVaclav Hapla Level: intermediate 5168fa9e22eSVaclav Hapla 5172f04c522SBarry Smith Notes: 5182f04c522SBarry Smith If the global length of `lmap` differs from the global length of `rmap` then the excess entries are ignored. 5192f04c522SBarry Smith 5202f04c522SBarry Smith The resulting `sf` used with `PetscSFBcastBegin()` and `PetscSFBcastEnd()` merely copies the array entries of `rootdata` to 5212f04c522SBarry Smith `leafdata`; moving them between MPI processes if needed. For example, 5222f04c522SBarry Smith if rmap is [0, 3, 5) and lmap is [0, 2, 6) and `rootdata` is (1, 2, 3) on MPI rank 0 and (4, 5) on MPI rank 1 then the 5232f04c522SBarry Smith `leafdata` would become (1, 2) on MPI rank 0 and (3, 4, 5, x) on MPI rank 1. 5242f04c522SBarry Smith 5252f04c522SBarry Smith .seealso: [](sec_petscsf), `PetscSF`, `PetscLayout`, `PetscSFCreate()`, `PetscSFSetGraph()`, `PetscLayoutCreate()`, `PetscSFSetGraphLayout()` 5268fa9e22eSVaclav Hapla @*/ 527d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateFromLayouts(PetscLayout rmap, PetscLayout lmap, PetscSF *sf) 528d71ae5a4SJacob Faibussowitsch { 5298fa9e22eSVaclav Hapla PetscInt i, nroots, nleaves = 0; 5308fa9e22eSVaclav Hapla PetscInt rN, lst, len; 5318fa9e22eSVaclav Hapla PetscMPIInt owner = -1; 5328fa9e22eSVaclav Hapla PetscSFNode *remote; 5338fa9e22eSVaclav Hapla MPI_Comm rcomm = rmap->comm; 5348fa9e22eSVaclav Hapla MPI_Comm lcomm = lmap->comm; 5358fa9e22eSVaclav Hapla PetscMPIInt flg; 5368fa9e22eSVaclav Hapla 5378fa9e22eSVaclav Hapla PetscFunctionBegin; 538*18aa8208SJames Wright PetscAssertPointer(rmap, 1); 539*18aa8208SJames Wright PetscAssertPointer(lmap, 2); 5404f572ea9SToby Isaac PetscAssertPointer(sf, 3); 54128b400f6SJacob Faibussowitsch PetscCheck(rmap->setupcalled, rcomm, PETSC_ERR_ARG_WRONGSTATE, "Root layout not setup"); 54228b400f6SJacob Faibussowitsch PetscCheck(lmap->setupcalled, lcomm, PETSC_ERR_ARG_WRONGSTATE, "Leaf layout not setup"); 5439566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(rcomm, lcomm, &flg)); 544c9cc58a2SBarry Smith PetscCheck(flg == MPI_CONGRUENT || flg == MPI_IDENT, rcomm, PETSC_ERR_SUP, "cannot map two layouts with non-matching communicators"); 5459566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(rcomm, sf)); 5469566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(rmap, &nroots)); 5479566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(rmap, &rN)); 5489566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRange(lmap, &lst, &len)); 5499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len - lst, &remote)); 5508fa9e22eSVaclav Hapla for (i = lst; i < len && i < rN; i++) { 55148a46eb9SPierre Jolivet if (owner < -1 || i >= rmap->range[owner + 1]) PetscCall(PetscLayoutFindOwner(rmap, i, &owner)); 5528fa9e22eSVaclav Hapla remote[nleaves].rank = owner; 5538fa9e22eSVaclav Hapla remote[nleaves].index = i - rmap->range[owner]; 5548fa9e22eSVaclav Hapla nleaves++; 5558fa9e22eSVaclav Hapla } 5569566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*sf, nroots, nleaves, NULL, PETSC_OWN_POINTER, remote, PETSC_COPY_VALUES)); 5579566063dSJacob Faibussowitsch PetscCall(PetscFree(remote)); 5583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5598fa9e22eSVaclav Hapla } 5608fa9e22eSVaclav Hapla 5618fa9e22eSVaclav Hapla /* TODO: handle nooffprocentries like MatZeroRowsMapLocal_Private, since this code is the same */ 562ce78bad3SBarry Smith PetscErrorCode PetscLayoutMapLocal(PetscLayout map, PetscInt N, const PetscInt idxs[], PetscInt *on, PetscInt *oidxs[], PetscInt *ogidxs[]) 563d71ae5a4SJacob Faibussowitsch { 5648fa9e22eSVaclav Hapla PetscInt *owners = map->range; 5658fa9e22eSVaclav Hapla PetscInt n = map->n; 5668fa9e22eSVaclav Hapla PetscSF sf; 567d61846d3SStefano Zampini PetscInt *lidxs, *work = NULL, *ilocal; 5688fa9e22eSVaclav Hapla PetscSFNode *ridxs; 5698fa9e22eSVaclav Hapla PetscMPIInt rank, p = 0; 570d61846d3SStefano Zampini PetscInt r, len = 0, nleaves = 0; 5718fa9e22eSVaclav Hapla 5728fa9e22eSVaclav Hapla PetscFunctionBegin; 5738fa9e22eSVaclav Hapla if (on) *on = 0; /* squelch -Wmaybe-uninitialized */ 5748fa9e22eSVaclav Hapla /* Create SF where leaves are input idxs and roots are owned idxs */ 5759566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(map->comm, &rank)); 5769566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &lidxs)); 5778fa9e22eSVaclav Hapla for (r = 0; r < n; ++r) lidxs[r] = -1; 5789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N, &ridxs)); 579d61846d3SStefano Zampini PetscCall(PetscMalloc1(N, &ilocal)); 5808fa9e22eSVaclav Hapla for (r = 0; r < N; ++r) { 5818fa9e22eSVaclav Hapla const PetscInt idx = idxs[r]; 582d61846d3SStefano Zampini 583d61846d3SStefano Zampini if (idx < 0) continue; 584d61846d3SStefano Zampini PetscCheck(idx < map->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Index %" PetscInt_FMT " out of range [0,%" PetscInt_FMT ")", idx, map->N); 5858fa9e22eSVaclav Hapla if (idx < owners[p] || owners[p + 1] <= idx) { /* short-circuit the search if the last p owns this idx too */ 5869566063dSJacob Faibussowitsch PetscCall(PetscLayoutFindOwner(map, idx, &p)); 5878fa9e22eSVaclav Hapla } 588d61846d3SStefano Zampini ridxs[nleaves].rank = p; 589d61846d3SStefano Zampini ridxs[nleaves].index = idxs[r] - owners[p]; 590d61846d3SStefano Zampini ilocal[nleaves] = r; 591d61846d3SStefano Zampini nleaves++; 5928fa9e22eSVaclav Hapla } 5939566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(map->comm, &sf)); 594d61846d3SStefano Zampini PetscCall(PetscSFSetGraph(sf, n, nleaves, ilocal, PETSC_OWN_POINTER, ridxs, PETSC_OWN_POINTER)); 5959566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR)); 5969566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR)); 5978fa9e22eSVaclav Hapla if (ogidxs) { /* communicate global idxs */ 5988fa9e22eSVaclav Hapla PetscInt cum = 0, start, *work2; 5998fa9e22eSVaclav Hapla 6009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &work)); 6019566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(N, &work2)); 6029371c9d4SSatish Balay for (r = 0; r < N; ++r) 6039371c9d4SSatish Balay if (idxs[r] >= 0) cum++; 6049566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&cum, &start, 1, MPIU_INT, MPI_SUM, map->comm)); 6058fa9e22eSVaclav Hapla start -= cum; 6068fa9e22eSVaclav Hapla cum = 0; 6079371c9d4SSatish Balay for (r = 0; r < N; ++r) 6089371c9d4SSatish Balay if (idxs[r] >= 0) work2[r] = start + cum++; 6099566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf, MPIU_INT, work2, work, MPI_REPLACE)); 6109566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf, MPIU_INT, work2, work, MPI_REPLACE)); 6119566063dSJacob Faibussowitsch PetscCall(PetscFree(work2)); 6128fa9e22eSVaclav Hapla } 6139566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf)); 6148fa9e22eSVaclav Hapla /* Compress and put in indices */ 6158fa9e22eSVaclav Hapla for (r = 0; r < n; ++r) 6168fa9e22eSVaclav Hapla if (lidxs[r] >= 0) { 6178fa9e22eSVaclav Hapla if (work) work[len] = work[r]; 6188fa9e22eSVaclav Hapla lidxs[len++] = r; 6198fa9e22eSVaclav Hapla } 6208fa9e22eSVaclav Hapla if (on) *on = len; 6218fa9e22eSVaclav Hapla if (oidxs) *oidxs = lidxs; 6228fa9e22eSVaclav Hapla if (ogidxs) *ogidxs = work; 6233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6248fa9e22eSVaclav Hapla } 625deffd5ebSksagiyam 626deffd5ebSksagiyam /*@ 627cab54364SBarry Smith PetscSFCreateByMatchingIndices - Create `PetscSF` by matching root and leaf indices 628deffd5ebSksagiyam 629deffd5ebSksagiyam Collective 630deffd5ebSksagiyam 631deffd5ebSksagiyam Input Parameters: 632cf84bf9eSBarry Smith + layout - `PetscLayout` defining the global index space and the MPI rank that brokers each index 633cf84bf9eSBarry Smith . numRootIndices - size of `rootIndices` 634cf84bf9eSBarry Smith . rootIndices - array of global indices of which this process requests ownership 635cf84bf9eSBarry Smith . rootLocalIndices - root local index permutation (`NULL` if no permutation) 636cf84bf9eSBarry Smith . rootLocalOffset - offset to be added to `rootLocalIndices` 637cf84bf9eSBarry Smith . numLeafIndices - size of `leafIndices` 638cf84bf9eSBarry Smith . leafIndices - array of global indices with which this process requires data associated 639cf84bf9eSBarry Smith . leafLocalIndices - leaf local index permutation (`NULL` if no permutation) 640cf84bf9eSBarry Smith - leafLocalOffset - offset to be added to `leafLocalIndices` 641deffd5ebSksagiyam 642d8d19677SJose E. Roman Output Parameters: 643cf84bf9eSBarry Smith + sfA - star forest representing the communication pattern from the layout space to the leaf space (`NULL` if not needed) 644deffd5ebSksagiyam - sf - star forest representing the communication pattern from the root space to the leaf space 645deffd5ebSksagiyam 646cab54364SBarry Smith Level: advanced 647cab54364SBarry Smith 648deffd5ebSksagiyam Example 1: 649cab54364SBarry Smith .vb 650cab54364SBarry Smith rank : 0 1 2 651cab54364SBarry Smith rootIndices : [1 0 2] [3] [3] 652cab54364SBarry Smith rootLocalOffset : 100 200 300 653cab54364SBarry Smith layout : [0 1] [2] [3] 654cab54364SBarry Smith leafIndices : [0] [2] [0 3] 655cab54364SBarry Smith leafLocalOffset : 400 500 600 656cab54364SBarry Smith 657cab54364SBarry Smith would build the following PetscSF 658cab54364SBarry Smith 659cab54364SBarry Smith [0] 400 <- (0,101) 660cab54364SBarry Smith [1] 500 <- (0,102) 661cab54364SBarry Smith [2] 600 <- (0,101) 662cab54364SBarry Smith [2] 601 <- (2,300) 663cab54364SBarry Smith .ve 664cab54364SBarry Smith 665deffd5ebSksagiyam Example 2: 666cab54364SBarry Smith .vb 667cab54364SBarry Smith rank : 0 1 2 668cab54364SBarry Smith rootIndices : [1 0 2] [3] [3] 669cab54364SBarry Smith rootLocalOffset : 100 200 300 670cab54364SBarry Smith layout : [0 1] [2] [3] 671cab54364SBarry Smith leafIndices : rootIndices rootIndices rootIndices 672cab54364SBarry Smith leafLocalOffset : rootLocalOffset rootLocalOffset rootLocalOffset 673cab54364SBarry Smith 674cab54364SBarry Smith would build the following PetscSF 675cab54364SBarry Smith 676cab54364SBarry Smith [1] 200 <- (2,300) 677cab54364SBarry Smith .ve 678cab54364SBarry Smith 679deffd5ebSksagiyam Example 3: 680cab54364SBarry Smith .vb 681cab54364SBarry Smith No process requests ownership of global index 1, but no process needs it. 682cab54364SBarry Smith 683cab54364SBarry Smith rank : 0 1 2 684cab54364SBarry Smith numRootIndices : 2 1 1 685cab54364SBarry Smith rootIndices : [0 2] [3] [3] 686cab54364SBarry Smith rootLocalOffset : 100 200 300 687cab54364SBarry Smith layout : [0 1] [2] [3] 688cab54364SBarry Smith numLeafIndices : 1 1 2 689cab54364SBarry Smith leafIndices : [0] [2] [0 3] 690cab54364SBarry Smith leafLocalOffset : 400 500 600 691cab54364SBarry Smith 692cab54364SBarry Smith would build the following PetscSF 693cab54364SBarry Smith 694cab54364SBarry Smith [0] 400 <- (0,100) 695cab54364SBarry Smith [1] 500 <- (0,101) 696cab54364SBarry Smith [2] 600 <- (0,100) 697cab54364SBarry Smith [2] 601 <- (2,300) 698cab54364SBarry Smith .ve 699deffd5ebSksagiyam 700deffd5ebSksagiyam Notes: 7012f04c522SBarry Smith `layout` represents any partitioning of [0, N), where N is the total number of global indices, and its 702cab54364SBarry Smith local size can be set to `PETSC_DECIDE`. 703cab54364SBarry Smith 7042f04c522SBarry Smith If a global index x lies in the partition owned by process i, each process whose `rootIndices` contains x requests 705deffd5ebSksagiyam ownership of x and sends its own rank and the local index of x to process i. 706deffd5ebSksagiyam If multiple processes request ownership of x, the one with the highest rank is to own x. 7072f04c522SBarry Smith Process i then broadcasts the ownership information, so that each process whose `leafIndices` contains x knows the 708deffd5ebSksagiyam ownership information of x. 7092f04c522SBarry Smith The output `sf` is constructed by associating each leaf point to a root point in this way. 710deffd5ebSksagiyam 711deffd5ebSksagiyam Suppose there is point data ordered according to the global indices and partitioned according to the given layout. 7122f04c522SBarry Smith The optional output `sfA` can be used to push such data to leaf points. 713deffd5ebSksagiyam 7142f04c522SBarry Smith All indices in `rootIndices` and `leafIndices` must lie in the layout range. The union (over all processes) of `rootIndices` 7152f04c522SBarry Smith must cover that of `leafIndices`, but need not cover the entire layout. 716deffd5ebSksagiyam 717deffd5ebSksagiyam If (leafIndices, leafLocalIndices, leafLocalOffset) == (rootIndices, rootLocalIndices, rootLocalOffset), the output 718deffd5ebSksagiyam star forest is almost identity, so will only include non-trivial part of the map. 719deffd5ebSksagiyam 72038b5cf2dSJacob Faibussowitsch Developer Notes: 721deffd5ebSksagiyam Current approach of a process of the highest rank gaining the ownership may cause load imbalance; consider using 722deffd5ebSksagiyam hash(rank, root_local_index) as the bid for the ownership determination. 723deffd5ebSksagiyam 7242f04c522SBarry Smith .seealso: [](sec_petscsf), `PetscSF`, `PetscSFCreate()` 725deffd5ebSksagiyam @*/ 726cf84bf9eSBarry Smith 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) 727d71ae5a4SJacob Faibussowitsch { 728deffd5ebSksagiyam MPI_Comm comm = layout->comm; 729deffd5ebSksagiyam PetscMPIInt size, rank; 730deffd5ebSksagiyam PetscSF sf1; 731deffd5ebSksagiyam PetscSFNode *owners, *buffer, *iremote; 732deffd5ebSksagiyam PetscInt *ilocal, nleaves, N, n, i; 733deffd5ebSksagiyam #if defined(PETSC_USE_DEBUG) 734deffd5ebSksagiyam PetscInt N1; 735deffd5ebSksagiyam #endif 736deffd5ebSksagiyam PetscBool flag; 737deffd5ebSksagiyam 738deffd5ebSksagiyam PetscFunctionBegin; 739*18aa8208SJames Wright PetscAssertPointer(layout, 1); 7404f572ea9SToby Isaac if (rootIndices) PetscAssertPointer(rootIndices, 3); 7414f572ea9SToby Isaac if (rootLocalIndices) PetscAssertPointer(rootLocalIndices, 4); 7424f572ea9SToby Isaac if (leafIndices) PetscAssertPointer(leafIndices, 7); 7434f572ea9SToby Isaac if (leafLocalIndices) PetscAssertPointer(leafLocalIndices, 8); 7444f572ea9SToby Isaac if (sfA) PetscAssertPointer(sfA, 10); 7454f572ea9SToby Isaac PetscAssertPointer(sf, 11); 74608401ef6SPierre Jolivet PetscCheck(numRootIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numRootIndices (%" PetscInt_FMT ") must be non-negative", numRootIndices); 74708401ef6SPierre Jolivet PetscCheck(numLeafIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numLeafIndices (%" PetscInt_FMT ") must be non-negative", numLeafIndices); 7489566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 7499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 7509566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(layout)); 7519566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(layout, &N)); 7529566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(layout, &n)); 753deffd5ebSksagiyam flag = (PetscBool)(leafIndices == rootIndices); 7545440e5dcSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &flag, 1, MPI_C_BOOL, MPI_LAND, comm)); 755c9cc58a2SBarry Smith PetscCheck(!flag || numLeafIndices == numRootIndices, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "leafIndices == rootIndices, but numLeafIndices (%" PetscInt_FMT ") != numRootIndices(%" PetscInt_FMT ")", numLeafIndices, numRootIndices); 756deffd5ebSksagiyam #if defined(PETSC_USE_DEBUG) 7571690c2aeSBarry Smith N1 = PETSC_INT_MIN; 7589371c9d4SSatish Balay for (i = 0; i < numRootIndices; i++) 7599371c9d4SSatish Balay if (rootIndices[i] > N1) N1 = rootIndices[i]; 760462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm)); 76108401ef6SPierre Jolivet PetscCheck(N1 < N, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Max. root index (%" PetscInt_FMT ") out of layout range [0,%" PetscInt_FMT ")", N1, N); 762deffd5ebSksagiyam if (!flag) { 7631690c2aeSBarry Smith N1 = PETSC_INT_MIN; 7649371c9d4SSatish Balay for (i = 0; i < numLeafIndices; i++) 7659371c9d4SSatish Balay if (leafIndices[i] > N1) N1 = leafIndices[i]; 766462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm)); 76708401ef6SPierre Jolivet PetscCheck(N1 < N, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Max. leaf index (%" PetscInt_FMT ") out of layout range [0,%" PetscInt_FMT ")", N1, N); 768deffd5ebSksagiyam } 769deffd5ebSksagiyam #endif 770deffd5ebSksagiyam /* Reduce: owners -> buffer */ 7719566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &buffer)); 7729566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, &sf1)); 7739566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sf1)); 7749566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraphLayout(sf1, layout, numRootIndices, NULL, PETSC_OWN_POINTER, rootIndices)); 7759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numRootIndices, &owners)); 776deffd5ebSksagiyam for (i = 0; i < numRootIndices; ++i) { 777deffd5ebSksagiyam owners[i].rank = rank; 778deffd5ebSksagiyam owners[i].index = rootLocalOffset + (rootLocalIndices ? rootLocalIndices[i] : i); 779deffd5ebSksagiyam } 780deffd5ebSksagiyam for (i = 0; i < n; ++i) { 781deffd5ebSksagiyam buffer[i].index = -1; 782deffd5ebSksagiyam buffer[i].rank = -1; 783deffd5ebSksagiyam } 7846497c311SBarry Smith PetscCall(PetscSFReduceBegin(sf1, MPIU_SF_NODE, owners, buffer, MPI_MAXLOC)); 7856497c311SBarry Smith PetscCall(PetscSFReduceEnd(sf1, MPIU_SF_NODE, owners, buffer, MPI_MAXLOC)); 786deffd5ebSksagiyam /* Bcast: buffer -> owners */ 787deffd5ebSksagiyam if (!flag) { 788deffd5ebSksagiyam /* leafIndices is different from rootIndices */ 7899566063dSJacob Faibussowitsch PetscCall(PetscFree(owners)); 7909566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraphLayout(sf1, layout, numLeafIndices, NULL, PETSC_OWN_POINTER, leafIndices)); 7919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numLeafIndices, &owners)); 792deffd5ebSksagiyam } 7936497c311SBarry Smith PetscCall(PetscSFBcastBegin(sf1, MPIU_SF_NODE, buffer, owners, MPI_REPLACE)); 7946497c311SBarry Smith PetscCall(PetscSFBcastEnd(sf1, MPIU_SF_NODE, buffer, owners, MPI_REPLACE)); 79508401ef6SPierre Jolivet 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]); 7969566063dSJacob Faibussowitsch PetscCall(PetscFree(buffer)); 7979371c9d4SSatish Balay if (sfA) { 7989371c9d4SSatish Balay *sfA = sf1; 7999371c9d4SSatish Balay } else PetscCall(PetscSFDestroy(&sf1)); 800deffd5ebSksagiyam /* Create sf */ 801deffd5ebSksagiyam if (flag && rootLocalIndices == leafLocalIndices && leafLocalOffset == rootLocalOffset) { 802deffd5ebSksagiyam /* leaf space == root space */ 8039371c9d4SSatish Balay for (i = 0, nleaves = 0; i < numLeafIndices; ++i) 8049371c9d4SSatish Balay if (owners[i].rank != rank) ++nleaves; 8059566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &ilocal)); 8069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &iremote)); 807deffd5ebSksagiyam for (i = 0, nleaves = 0; i < numLeafIndices; ++i) { 808deffd5ebSksagiyam if (owners[i].rank != rank) { 809deffd5ebSksagiyam ilocal[nleaves] = leafLocalOffset + i; 810deffd5ebSksagiyam iremote[nleaves].rank = owners[i].rank; 811deffd5ebSksagiyam iremote[nleaves].index = owners[i].index; 812deffd5ebSksagiyam ++nleaves; 813deffd5ebSksagiyam } 814deffd5ebSksagiyam } 8159566063dSJacob Faibussowitsch PetscCall(PetscFree(owners)); 816deffd5ebSksagiyam } else { 817deffd5ebSksagiyam nleaves = numLeafIndices; 8189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &ilocal)); 819ad540459SPierre Jolivet for (i = 0; i < nleaves; ++i) ilocal[i] = leafLocalOffset + (leafLocalIndices ? leafLocalIndices[i] : i); 820deffd5ebSksagiyam iremote = owners; 821deffd5ebSksagiyam } 8229566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, sf)); 8239566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(*sf)); 8249566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*sf, rootLocalOffset + numRootIndices, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER)); 8253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 826deffd5ebSksagiyam } 827fbc7033fSJed Brown 828fbc7033fSJed Brown /*@ 82953c0d4aeSBarry Smith PetscSFMerge - append/merge indices of `sfb` into `sfa`, with preference for `sfb` 830fbc7033fSJed Brown 831fbc7033fSJed Brown Collective 832fbc7033fSJed Brown 83338b5cf2dSJacob Faibussowitsch Input Parameters: 834fbc7033fSJed Brown + sfa - default `PetscSF` 8352f04c522SBarry Smith - sfb - additional edges to add/replace edges in `sfa` 836fbc7033fSJed Brown 83738b5cf2dSJacob Faibussowitsch Output Parameter: 838fbc7033fSJed Brown . merged - new `PetscSF` with combined edges 839fbc7033fSJed Brown 84053c0d4aeSBarry Smith Level: intermediate 84153c0d4aeSBarry Smith 8422f04c522SBarry Smith .seealso: [](sec_petscsf), `PetscSF`, `PetscSFCompose()` 843fbc7033fSJed Brown @*/ 844fbc7033fSJed Brown PetscErrorCode PetscSFMerge(PetscSF sfa, PetscSF sfb, PetscSF *merged) 845fbc7033fSJed Brown { 846fbc7033fSJed Brown PetscInt maxleaf; 847fbc7033fSJed Brown 848fbc7033fSJed Brown PetscFunctionBegin; 849fbc7033fSJed Brown PetscValidHeaderSpecific(sfa, PETSCSF_CLASSID, 1); 850fbc7033fSJed Brown PetscValidHeaderSpecific(sfb, PETSCSF_CLASSID, 2); 851fbc7033fSJed Brown PetscCheckSameComm(sfa, 1, sfb, 2); 8524f572ea9SToby Isaac PetscAssertPointer(merged, 3); 853fbc7033fSJed Brown { 854fbc7033fSJed Brown PetscInt aleaf, bleaf; 855fbc7033fSJed Brown PetscCall(PetscSFGetLeafRange(sfa, NULL, &aleaf)); 856fbc7033fSJed Brown PetscCall(PetscSFGetLeafRange(sfb, NULL, &bleaf)); 857fbc7033fSJed Brown maxleaf = PetscMax(aleaf, bleaf) + 1; // One more than the last index 858fbc7033fSJed Brown } 859fbc7033fSJed Brown PetscInt *clocal, aroots, aleaves, broots, bleaves; 860fbc7033fSJed Brown PetscSFNode *cremote; 861fbc7033fSJed Brown const PetscInt *alocal, *blocal; 862fbc7033fSJed Brown const PetscSFNode *aremote, *bremote; 863fbc7033fSJed Brown PetscCall(PetscMalloc2(maxleaf, &clocal, maxleaf, &cremote)); 864fbc7033fSJed Brown for (PetscInt i = 0; i < maxleaf; i++) clocal[i] = -1; 865fbc7033fSJed Brown PetscCall(PetscSFGetGraph(sfa, &aroots, &aleaves, &alocal, &aremote)); 866fbc7033fSJed Brown PetscCall(PetscSFGetGraph(sfb, &broots, &bleaves, &blocal, &bremote)); 867fbc7033fSJed Brown PetscCheck(aroots == broots, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Both sfa and sfb must have the same root space"); 868fbc7033fSJed Brown for (PetscInt i = 0; i < aleaves; i++) { 869fbc7033fSJed Brown PetscInt a = alocal ? alocal[i] : i; 870fbc7033fSJed Brown clocal[a] = a; 871fbc7033fSJed Brown cremote[a] = aremote[i]; 872fbc7033fSJed Brown } 873fbc7033fSJed Brown for (PetscInt i = 0; i < bleaves; i++) { 874fbc7033fSJed Brown PetscInt b = blocal ? blocal[i] : i; 875fbc7033fSJed Brown clocal[b] = b; 876fbc7033fSJed Brown cremote[b] = bremote[i]; 877fbc7033fSJed Brown } 878fbc7033fSJed Brown PetscInt nleaves = 0; 879fbc7033fSJed Brown for (PetscInt i = 0; i < maxleaf; i++) { 880fbc7033fSJed Brown if (clocal[i] < 0) continue; 881fbc7033fSJed Brown clocal[nleaves] = clocal[i]; 882fbc7033fSJed Brown cremote[nleaves] = cremote[i]; 883fbc7033fSJed Brown nleaves++; 884fbc7033fSJed Brown } 885fbc7033fSJed Brown PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)sfa), merged)); 886fbc7033fSJed Brown PetscCall(PetscSFSetGraph(*merged, aroots, nleaves, clocal, PETSC_COPY_VALUES, cremote, PETSC_COPY_VALUES)); 887fbc7033fSJed Brown PetscCall(PetscFree2(clocal, cremote)); 8883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 889fbc7033fSJed Brown } 8900dd791a8SStefano Zampini 8910dd791a8SStefano Zampini /*@ 8920dd791a8SStefano Zampini PetscSFCreateStridedSF - Create an `PetscSF` to communicate interleaved blocks of data 8930dd791a8SStefano Zampini 8940dd791a8SStefano Zampini Collective 8950dd791a8SStefano Zampini 8960dd791a8SStefano Zampini Input Parameters: 8970dd791a8SStefano Zampini + sf - star forest 8980dd791a8SStefano Zampini . bs - stride 8990dd791a8SStefano Zampini . ldr - leading dimension of root space 9000dd791a8SStefano Zampini - ldl - leading dimension of leaf space 9010dd791a8SStefano Zampini 9020dd791a8SStefano Zampini Output Parameter: 9030dd791a8SStefano Zampini . vsf - the new `PetscSF` 9040dd791a8SStefano Zampini 9050dd791a8SStefano Zampini Level: intermediate 9060dd791a8SStefano Zampini 9070dd791a8SStefano Zampini Notes: 9082f04c522SBarry Smith This can be useful to perform communications on multiple right-hand sides stored in a Fortran-style two dimensional array. 9092f04c522SBarry Smith For example, the calling sequence 9100dd791a8SStefano Zampini .vb 9110dd791a8SStefano Zampini c_datatype *roots, *leaves; 9120dd791a8SStefano Zampini for i in [0,bs) do 9130dd791a8SStefano Zampini PetscSFBcastBegin(sf, mpi_datatype, roots + i*ldr, leaves + i*ldl, op) 9140dd791a8SStefano Zampini PetscSFBcastEnd(sf, mpi_datatype, roots + i*ldr, leaves + i*ldl, op) 9150dd791a8SStefano Zampini .ve 9160dd791a8SStefano Zampini is equivalent to 9170dd791a8SStefano Zampini .vb 9180dd791a8SStefano Zampini c_datatype *roots, *leaves; 9190dd791a8SStefano Zampini PetscSFCreateStridedSF(sf, bs, ldr, ldl, &vsf) 9200dd791a8SStefano Zampini PetscSFBcastBegin(vsf, mpi_datatype, roots, leaves, op) 9210dd791a8SStefano Zampini PetscSFBcastEnd(vsf, mpi_datatype, roots, leaves, op) 9220dd791a8SStefano Zampini .ve 9230dd791a8SStefano Zampini 9240dd791a8SStefano Zampini Developer Notes: 9250dd791a8SStefano Zampini Should this functionality be handled with a new API instead of creating a new object? 9260dd791a8SStefano Zampini 9272f04c522SBarry Smith .seealso: [](sec_petscsf), `PetscSF`, `PetscSFCreate()`, `PetscSFSetGraph()` 9280dd791a8SStefano Zampini @*/ 9290dd791a8SStefano Zampini PetscErrorCode PetscSFCreateStridedSF(PetscSF sf, PetscInt bs, PetscInt ldr, PetscInt ldl, PetscSF *vsf) 9300dd791a8SStefano Zampini { 9310dd791a8SStefano Zampini PetscSF rankssf; 9320dd791a8SStefano Zampini const PetscSFNode *iremote, *sfrremote; 9330dd791a8SStefano Zampini PetscSFNode *viremote; 9340dd791a8SStefano Zampini const PetscInt *ilocal; 9350dd791a8SStefano Zampini PetscInt *vilocal = NULL, *ldrs; 9360dd791a8SStefano Zampini PetscInt nranks, nr, nl, vnr, vnl, maxl; 9370dd791a8SStefano Zampini PetscMPIInt rank; 9380dd791a8SStefano Zampini MPI_Comm comm; 9390dd791a8SStefano Zampini PetscSFType sftype; 9400dd791a8SStefano Zampini 9410dd791a8SStefano Zampini PetscFunctionBegin; 9420dd791a8SStefano Zampini PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 9430dd791a8SStefano Zampini PetscValidLogicalCollectiveInt(sf, bs, 2); 9440dd791a8SStefano Zampini PetscAssertPointer(vsf, 5); 9450dd791a8SStefano Zampini if (bs == 1) { 9460dd791a8SStefano Zampini PetscCall(PetscObjectReference((PetscObject)sf)); 9470dd791a8SStefano Zampini *vsf = sf; 9480dd791a8SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 9490dd791a8SStefano Zampini } 9500dd791a8SStefano Zampini PetscCall(PetscSFSetUp(sf)); 9510dd791a8SStefano Zampini PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 9520dd791a8SStefano Zampini PetscCallMPI(MPI_Comm_rank(comm, &rank)); 9530dd791a8SStefano Zampini PetscCall(PetscSFGetGraph(sf, &nr, &nl, &ilocal, &iremote)); 9540dd791a8SStefano Zampini PetscCall(PetscSFGetLeafRange(sf, NULL, &maxl)); 9550dd791a8SStefano Zampini maxl += 1; 9560dd791a8SStefano Zampini if (ldl == PETSC_DECIDE) ldl = maxl; 9570dd791a8SStefano Zampini if (ldr == PETSC_DECIDE) ldr = nr; 9580dd791a8SStefano Zampini PetscCheck(ldr >= nr, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid leading dimension %" PetscInt_FMT " must be smaller than number of roots %" PetscInt_FMT, ldr, nr); 9590dd791a8SStefano Zampini PetscCheck(ldl >= maxl, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid leading dimension %" PetscInt_FMT " must be larger than leaf range %" PetscInt_FMT, ldl, maxl - 1); 9600dd791a8SStefano Zampini vnr = nr * bs; 9610dd791a8SStefano Zampini vnl = nl * bs; 9620dd791a8SStefano Zampini PetscCall(PetscMalloc1(vnl, &viremote)); 9630dd791a8SStefano Zampini PetscCall(PetscMalloc1(vnl, &vilocal)); 9640dd791a8SStefano Zampini 9650dd791a8SStefano Zampini /* Communicate root leading dimensions to leaf ranks */ 9660dd791a8SStefano Zampini PetscCall(PetscSFGetRanksSF(sf, &rankssf)); 9670dd791a8SStefano Zampini PetscCall(PetscSFGetGraph(rankssf, NULL, &nranks, NULL, &sfrremote)); 9680dd791a8SStefano Zampini PetscCall(PetscMalloc1(nranks, &ldrs)); 9690dd791a8SStefano Zampini PetscCall(PetscSFBcastBegin(rankssf, MPIU_INT, &ldr, ldrs, MPI_REPLACE)); 9700dd791a8SStefano Zampini PetscCall(PetscSFBcastEnd(rankssf, MPIU_INT, &ldr, ldrs, MPI_REPLACE)); 9710dd791a8SStefano Zampini 9720dd791a8SStefano Zampini for (PetscInt i = 0, rold = -1, lda = -1; i < nl; i++) { 973835f2295SStefano Zampini const PetscInt r = iremote[i].rank; 9740dd791a8SStefano Zampini const PetscInt ii = iremote[i].index; 9750dd791a8SStefano Zampini 9760dd791a8SStefano Zampini if (r == rank) lda = ldr; 9770dd791a8SStefano Zampini else if (rold != r) { 9780dd791a8SStefano Zampini PetscInt j; 9790dd791a8SStefano Zampini 9800dd791a8SStefano Zampini for (j = 0; j < nranks; j++) 9810dd791a8SStefano Zampini if (sfrremote[j].rank == r) break; 982835f2295SStefano Zampini PetscCheck(j < nranks, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unable to locate neighbor rank %" PetscInt_FMT, r); 9830dd791a8SStefano Zampini lda = ldrs[j]; 9840dd791a8SStefano Zampini } 9850dd791a8SStefano Zampini rold = r; 9860dd791a8SStefano Zampini for (PetscInt v = 0; v < bs; v++) { 9870dd791a8SStefano Zampini viremote[v * nl + i].rank = r; 9880dd791a8SStefano Zampini viremote[v * nl + i].index = v * lda + ii; 9890dd791a8SStefano Zampini vilocal[v * nl + i] = v * ldl + (ilocal ? ilocal[i] : i); 9900dd791a8SStefano Zampini } 9910dd791a8SStefano Zampini } 9920dd791a8SStefano Zampini PetscCall(PetscFree(ldrs)); 9930dd791a8SStefano Zampini PetscCall(PetscSFCreate(comm, vsf)); 9940dd791a8SStefano Zampini PetscCall(PetscSFGetType(sf, &sftype)); 9950dd791a8SStefano Zampini PetscCall(PetscSFSetType(*vsf, sftype)); 9960dd791a8SStefano Zampini PetscCall(PetscSFSetGraph(*vsf, vnr, vnl, vilocal, PETSC_OWN_POINTER, viremote, PETSC_OWN_POINTER)); 9970dd791a8SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 9980dd791a8SStefano Zampini } 999