1b0c7db22SLisandro Dalcin #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/ 2b0c7db22SLisandro Dalcin #include <petsc/private/sectionimpl.h> 3b0c7db22SLisandro Dalcin 4eb58282aSVaclav Hapla /*@ 5cab54364SBarry Smith PetscSFSetGraphLayout - Set a parallel star forest via global indices and a `PetscLayout` 6b0c7db22SLisandro Dalcin 7b0c7db22SLisandro Dalcin Collective 8b0c7db22SLisandro Dalcin 94165533cSJose E. Roman Input Parameters: 10b0c7db22SLisandro Dalcin + sf - star forest 11cab54364SBarry Smith . layout - `PetscLayout` defining the global space for roots 12b0c7db22SLisandro Dalcin . nleaves - number of leaf vertices on the current process, each of these references a root on any process 13b0c7db22SLisandro Dalcin . ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage 14b0c7db22SLisandro Dalcin . localmode - copy mode for ilocal 15eb58282aSVaclav Hapla - gremote - root vertices in global numbering corresponding to leaves in ilocal 16b0c7db22SLisandro Dalcin 17b0c7db22SLisandro Dalcin Level: intermediate 18b0c7db22SLisandro Dalcin 19cab54364SBarry Smith Note: 2086c56f52SVaclav Hapla Global indices must lie in [0, N) where N is the global size of layout. 21cab54364SBarry 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: 2486c56f52SVaclav Hapla 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 28cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFGetGraphLayout()`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()` 29b0c7db22SLisandro Dalcin @*/ 30d71ae5a4SJacob Faibussowitsch 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); 39b114984aSVaclav Hapla PetscCall(PetscLayoutSetUp(layout)); 409566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(layout, &nroots)); 419566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRanges(layout, &range)); 429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &remote)); 43eb58282aSVaclav Hapla if (nleaves) ls = gremote[0] + 1; 44b0c7db22SLisandro Dalcin for (i = 0; i < nleaves; i++) { 45eb58282aSVaclav Hapla const PetscInt idx = gremote[i] - ls; 4638a25198SStefano Zampini if (idx < 0 || idx >= ln) { /* short-circuit the search */ 47eb58282aSVaclav Hapla PetscCall(PetscLayoutFindOwnerIndex(layout, gremote[i], &lr, &remote[i].index)); 4838a25198SStefano Zampini remote[i].rank = lr; 4938a25198SStefano Zampini ls = range[lr]; 5038a25198SStefano Zampini ln = range[lr + 1] - ls; 5138a25198SStefano Zampini } else { 5238a25198SStefano Zampini remote[i].rank = lr; 5338a25198SStefano Zampini remote[i].index = idx; 5438a25198SStefano Zampini } 55b0c7db22SLisandro Dalcin } 569566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, nroots, nleaves, ilocal, localmode, remote, PETSC_OWN_POINTER)); 573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 58b0c7db22SLisandro Dalcin } 59b0c7db22SLisandro Dalcin 60eb58282aSVaclav Hapla /*@C 61cab54364SBarry Smith PetscSFGetGraphLayout - Get the global indices and `PetscLayout` that describe this star forest 62eb58282aSVaclav Hapla 63eb58282aSVaclav Hapla Collective 64eb58282aSVaclav Hapla 65eb58282aSVaclav Hapla Input Parameter: 66eb58282aSVaclav Hapla . sf - star forest 67eb58282aSVaclav Hapla 68eb58282aSVaclav Hapla Output Parameters: 69cab54364SBarry Smith + layout - `PetscLayout` defining the global space for roots 70eb58282aSVaclav Hapla . nleaves - number of leaf vertices on the current process, each of these references a root on any process 71f13dfd9eSBarry Smith . ilocal - locations of leaves in leafdata buffers, or `NULL` for contiguous storage 72eb58282aSVaclav Hapla - gremote - root vertices in global numbering corresponding to leaves in ilocal 73eb58282aSVaclav Hapla 74eb58282aSVaclav Hapla Level: intermediate 75eb58282aSVaclav Hapla 76eb58282aSVaclav Hapla Notes: 77eb58282aSVaclav Hapla The outputs are such that passing them as inputs to `PetscSFSetGraphLayout()` would lead to the same star forest. 78f13dfd9eSBarry Smith The outputs `layout` and `gremote` are freshly created each time this function is called, 79f13dfd9eSBarry Smith so they need to be freed (with `PetscLayoutDestroy()` and `PetscFree()`) by the user. 80eb58282aSVaclav Hapla 81cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFSetGraphLayout()`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()` 82eb58282aSVaclav Hapla @*/ 83d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetGraphLayout(PetscSF sf, PetscLayout *layout, PetscInt *nleaves, const PetscInt *ilocal[], PetscInt *gremote[]) 84d71ae5a4SJacob Faibussowitsch { 85eb58282aSVaclav Hapla PetscInt nr, nl; 86eb58282aSVaclav Hapla const PetscSFNode *ir; 87eb58282aSVaclav Hapla PetscLayout lt; 88eb58282aSVaclav Hapla 89eb58282aSVaclav Hapla PetscFunctionBegin; 90eb58282aSVaclav Hapla PetscCall(PetscSFGetGraph(sf, &nr, &nl, ilocal, &ir)); 91eb58282aSVaclav Hapla PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)sf), nr, PETSC_DECIDE, 1, <)); 92eb58282aSVaclav Hapla if (gremote) { 93eb58282aSVaclav Hapla PetscInt i; 94eb58282aSVaclav Hapla const PetscInt *range; 95eb58282aSVaclav Hapla PetscInt *gr; 96eb58282aSVaclav Hapla 97eb58282aSVaclav Hapla PetscCall(PetscLayoutGetRanges(lt, &range)); 98eb58282aSVaclav Hapla PetscCall(PetscMalloc1(nl, &gr)); 99eb58282aSVaclav Hapla for (i = 0; i < nl; i++) gr[i] = range[ir[i].rank] + ir[i].index; 100eb58282aSVaclav Hapla *gremote = gr; 101eb58282aSVaclav Hapla } 102eb58282aSVaclav Hapla if (nleaves) *nleaves = nl; 103eb58282aSVaclav Hapla if (layout) *layout = lt; 104eb58282aSVaclav Hapla else PetscCall(PetscLayoutDestroy(<)); 1053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 106eb58282aSVaclav Hapla } 107eb58282aSVaclav Hapla 108b0c7db22SLisandro Dalcin /*@ 109cab54364SBarry Smith PetscSFSetGraphSection - Sets the `PetscSF` graph encoding the parallel dof overlap based upon the `PetscSection` describing the data layout. 110b0c7db22SLisandro Dalcin 111b0c7db22SLisandro Dalcin Input Parameters: 112cab54364SBarry Smith + sf - The `PetscSF` 113cab54364SBarry Smith . localSection - `PetscSection` describing the local data layout 114cab54364SBarry Smith - globalSection - `PetscSection` describing the global data layout 115b0c7db22SLisandro Dalcin 116b0c7db22SLisandro Dalcin Level: developer 117b0c7db22SLisandro Dalcin 118cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFSetGraph()`, `PetscSFSetGraphLayout()` 119b0c7db22SLisandro Dalcin @*/ 120d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetGraphSection(PetscSF sf, PetscSection localSection, PetscSection globalSection) 121d71ae5a4SJacob Faibussowitsch { 122b0c7db22SLisandro Dalcin MPI_Comm comm; 123b0c7db22SLisandro Dalcin PetscLayout layout; 124b0c7db22SLisandro Dalcin const PetscInt *ranges; 125b0c7db22SLisandro Dalcin PetscInt *local; 126b0c7db22SLisandro Dalcin PetscSFNode *remote; 127b0c7db22SLisandro Dalcin PetscInt pStart, pEnd, p, nroots, nleaves = 0, l; 128b0c7db22SLisandro Dalcin PetscMPIInt size, rank; 129b0c7db22SLisandro Dalcin 130b0c7db22SLisandro Dalcin PetscFunctionBegin; 131b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 132b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(localSection, PETSC_SECTION_CLASSID, 2); 133b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(globalSection, PETSC_SECTION_CLASSID, 3); 134b0c7db22SLisandro Dalcin 1359566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 1369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 1379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1389566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(globalSection, &pStart, &pEnd)); 1399566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstrainedStorageSize(globalSection, &nroots)); 1409566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(comm, &layout)); 1419566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(layout, 1)); 1429566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetLocalSize(layout, nroots)); 1439566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(layout)); 1449566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRanges(layout, &ranges)); 145b0c7db22SLisandro Dalcin for (p = pStart; p < pEnd; ++p) { 146b0c7db22SLisandro Dalcin PetscInt gdof, gcdof; 147b0c7db22SLisandro Dalcin 1489566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(globalSection, p, &gdof)); 1499566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof)); 15057508eceSPierre 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); 151b0c7db22SLisandro Dalcin nleaves += gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof; 152b0c7db22SLisandro Dalcin } 1539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &local)); 1549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &remote)); 155b0c7db22SLisandro Dalcin for (p = pStart, l = 0; p < pEnd; ++p) { 156b0c7db22SLisandro Dalcin const PetscInt *cind; 157b0c7db22SLisandro Dalcin PetscInt dof, cdof, off, gdof, gcdof, goff, gsize, d, c; 158b0c7db22SLisandro Dalcin 1599566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(localSection, p, &dof)); 1609566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(localSection, p, &off)); 1619566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(localSection, p, &cdof)); 1629566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintIndices(localSection, p, &cind)); 1639566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(globalSection, p, &gdof)); 1649566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof)); 1659566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(globalSection, p, &goff)); 166b0c7db22SLisandro Dalcin if (!gdof) continue; /* Censored point */ 167b0c7db22SLisandro Dalcin gsize = gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof; 168b0c7db22SLisandro Dalcin if (gsize != dof - cdof) { 16908401ef6SPierre 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); 170b0c7db22SLisandro Dalcin cdof = 0; /* Ignore constraints */ 171b0c7db22SLisandro Dalcin } 172b0c7db22SLisandro Dalcin for (d = 0, c = 0; d < dof; ++d) { 1739371c9d4SSatish Balay if ((c < cdof) && (cind[c] == d)) { 1749371c9d4SSatish Balay ++c; 1759371c9d4SSatish Balay continue; 1769371c9d4SSatish Balay } 177b0c7db22SLisandro Dalcin local[l + d - c] = off + d; 178b0c7db22SLisandro Dalcin } 17908401ef6SPierre 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); 180b0c7db22SLisandro Dalcin if (gdof < 0) { 181b0c7db22SLisandro Dalcin for (d = 0; d < gsize; ++d, ++l) { 1826497c311SBarry Smith PetscInt offset = -(goff + 1) + d, ir; 1836497c311SBarry Smith PetscMPIInt r; 184b0c7db22SLisandro Dalcin 1856497c311SBarry Smith PetscCall(PetscFindInt(offset, size + 1, ranges, &ir)); 1866497c311SBarry Smith PetscCall(PetscMPIIntCast(ir, &r)); 187b0c7db22SLisandro Dalcin if (r < 0) r = -(r + 2); 1886497c311SBarry 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); 189b0c7db22SLisandro Dalcin remote[l].rank = r; 190b0c7db22SLisandro Dalcin remote[l].index = offset - ranges[r]; 191b0c7db22SLisandro Dalcin } 192b0c7db22SLisandro Dalcin } else { 193b0c7db22SLisandro Dalcin for (d = 0; d < gsize; ++d, ++l) { 194b0c7db22SLisandro Dalcin remote[l].rank = rank; 195b0c7db22SLisandro Dalcin remote[l].index = goff + d - ranges[rank]; 196b0c7db22SLisandro Dalcin } 197b0c7db22SLisandro Dalcin } 198b0c7db22SLisandro Dalcin } 19908401ef6SPierre Jolivet PetscCheck(l == nleaves, comm, PETSC_ERR_PLIB, "Iteration error, l %" PetscInt_FMT " != nleaves %" PetscInt_FMT, l, nleaves); 2009566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&layout)); 2019566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER)); 2023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 203b0c7db22SLisandro Dalcin } 204b0c7db22SLisandro Dalcin 205b0c7db22SLisandro Dalcin /*@C 206cab54364SBarry Smith PetscSFDistributeSection - Create a new `PetscSection` reorganized, moving from the root to the leaves of the `PetscSF` 207b0c7db22SLisandro Dalcin 208c3339decSBarry Smith Collective 209b0c7db22SLisandro Dalcin 210b0c7db22SLisandro Dalcin Input Parameters: 211cab54364SBarry Smith + sf - The `PetscSF` 212b0c7db22SLisandro Dalcin - rootSection - Section defined on root space 213b0c7db22SLisandro Dalcin 214b0c7db22SLisandro Dalcin Output Parameters: 2156497c311SBarry Smith + remoteOffsets - root offsets in leaf storage, or `NULL` 216b0c7db22SLisandro Dalcin - leafSection - Section defined on the leaf space 217b0c7db22SLisandro Dalcin 218b0c7db22SLisandro Dalcin Level: advanced 219b0c7db22SLisandro Dalcin 22023e9620eSJunchao Zhang Fortran Notes: 22123e9620eSJunchao Zhang In Fortran, use PetscSFDistributeSectionF90() 22223e9620eSJunchao Zhang 223cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()` 224b0c7db22SLisandro Dalcin @*/ 225d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFDistributeSection(PetscSF sf, PetscSection rootSection, PetscInt **remoteOffsets, PetscSection leafSection) 226d71ae5a4SJacob Faibussowitsch { 227b0c7db22SLisandro Dalcin PetscSF embedSF; 228b0c7db22SLisandro Dalcin const PetscInt *indices; 229b0c7db22SLisandro Dalcin IS selected; 2301690c2aeSBarry Smith PetscInt numFields, nroots, rpStart, rpEnd, lpStart = PETSC_INT_MAX, lpEnd = -1, f, c; 231b0c7db22SLisandro Dalcin PetscBool *sub, hasc; 232b0c7db22SLisandro Dalcin 233b0c7db22SLisandro Dalcin PetscFunctionBegin; 2349566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_DistSect, sf, 0, 0, 0)); 2359566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(rootSection, &numFields)); 236029e8818Sksagiyam if (numFields) { 237029e8818Sksagiyam IS perm; 238029e8818Sksagiyam 239029e8818Sksagiyam /* PetscSectionSetNumFields() calls PetscSectionReset(), which destroys 240029e8818Sksagiyam leafSection->perm. To keep this permutation set by the user, we grab 241029e8818Sksagiyam the reference before calling PetscSectionSetNumFields() and set it 242029e8818Sksagiyam back after. */ 2439566063dSJacob Faibussowitsch PetscCall(PetscSectionGetPermutation(leafSection, &perm)); 2449566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)perm)); 2459566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(leafSection, numFields)); 2469566063dSJacob Faibussowitsch PetscCall(PetscSectionSetPermutation(leafSection, perm)); 2479566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 248029e8818Sksagiyam } 2499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numFields + 2, &sub)); 250b0c7db22SLisandro Dalcin sub[1] = rootSection->bc ? PETSC_TRUE : PETSC_FALSE; 251b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) { 2522ee82556SMatthew G. Knepley PetscSectionSym sym, dsym = NULL; 253b0c7db22SLisandro Dalcin const char *name = NULL; 254b0c7db22SLisandro Dalcin PetscInt numComp = 0; 255b0c7db22SLisandro Dalcin 256b0c7db22SLisandro Dalcin sub[2 + f] = rootSection->field[f]->bc ? PETSC_TRUE : PETSC_FALSE; 2579566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldComponents(rootSection, f, &numComp)); 2589566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldName(rootSection, f, &name)); 2599566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldSym(rootSection, f, &sym)); 2609566063dSJacob Faibussowitsch if (sym) PetscCall(PetscSectionSymDistribute(sym, sf, &dsym)); 2619566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(leafSection, f, numComp)); 2629566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldName(leafSection, f, name)); 2639566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldSym(leafSection, f, dsym)); 2649566063dSJacob Faibussowitsch PetscCall(PetscSectionSymDestroy(&dsym)); 265b0c7db22SLisandro Dalcin for (c = 0; c < rootSection->numFieldComponents[f]; ++c) { 2669566063dSJacob Faibussowitsch PetscCall(PetscSectionGetComponentName(rootSection, f, c, &name)); 2679566063dSJacob Faibussowitsch PetscCall(PetscSectionSetComponentName(leafSection, f, c, name)); 268b0c7db22SLisandro Dalcin } 269b0c7db22SLisandro Dalcin } 2709566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd)); 2719566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL)); 272b0c7db22SLisandro Dalcin rpEnd = PetscMin(rpEnd, nroots); 273b0c7db22SLisandro Dalcin rpEnd = PetscMax(rpStart, rpEnd); 274b0c7db22SLisandro Dalcin /* see if we can avoid creating the embedded SF, since it can cost more than an allreduce */ 275b0c7db22SLisandro Dalcin sub[0] = (PetscBool)(nroots != rpEnd - rpStart); 276*e91c04dfSPierre Jolivet PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, sub, 2 + numFields, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)sf))); 277b0c7db22SLisandro Dalcin if (sub[0]) { 2789566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected)); 2799566063dSJacob Faibussowitsch PetscCall(ISGetIndices(selected, &indices)); 2809566063dSJacob Faibussowitsch PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF)); 2819566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(selected, &indices)); 2829566063dSJacob Faibussowitsch PetscCall(ISDestroy(&selected)); 283b0c7db22SLisandro Dalcin } else { 2849566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)sf)); 285b0c7db22SLisandro Dalcin embedSF = sf; 286b0c7db22SLisandro Dalcin } 2879566063dSJacob Faibussowitsch PetscCall(PetscSFGetLeafRange(embedSF, &lpStart, &lpEnd)); 288b0c7db22SLisandro Dalcin lpEnd++; 289b0c7db22SLisandro Dalcin 2909566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(leafSection, lpStart, lpEnd)); 291b0c7db22SLisandro Dalcin 292b0c7db22SLisandro Dalcin /* Constrained dof section */ 293b0c7db22SLisandro Dalcin hasc = sub[1]; 294b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) hasc = (PetscBool)(hasc || sub[2 + f]); 295b0c7db22SLisandro Dalcin 296b0c7db22SLisandro Dalcin /* Could fuse these at the cost of copies and extra allocation */ 29716cd844bSPierre Jolivet PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasDof, -rpStart), PetscSafePointerPlusOffset(leafSection->atlasDof, -lpStart), MPI_REPLACE)); 29816cd844bSPierre Jolivet PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasDof, -rpStart), PetscSafePointerPlusOffset(leafSection->atlasDof, -lpStart), MPI_REPLACE)); 299b0c7db22SLisandro Dalcin if (sub[1]) { 3007c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(rootSection)); 3017c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(leafSection)); 3029566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE)); 3039566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE)); 304b0c7db22SLisandro Dalcin } 305b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) { 30616cd844bSPierre Jolivet PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->field[f]->atlasDof, -rpStart), PetscSafePointerPlusOffset(leafSection->field[f]->atlasDof, -lpStart), MPI_REPLACE)); 30716cd844bSPierre Jolivet PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->field[f]->atlasDof, -rpStart), PetscSafePointerPlusOffset(leafSection->field[f]->atlasDof, -lpStart), MPI_REPLACE)); 308b0c7db22SLisandro Dalcin if (sub[2 + f]) { 3097c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(rootSection->field[f])); 3107c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(leafSection->field[f])); 3119566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE)); 3129566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE)); 313b0c7db22SLisandro Dalcin } 314b0c7db22SLisandro Dalcin } 315b0c7db22SLisandro Dalcin if (remoteOffsets) { 3169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lpEnd - lpStart, remoteOffsets)); 31716cd844bSPierre Jolivet PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(*remoteOffsets, -lpStart), MPI_REPLACE)); 31816cd844bSPierre Jolivet PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(*remoteOffsets, -lpStart), MPI_REPLACE)); 319b0c7db22SLisandro Dalcin } 32069c11d05SVaclav Hapla PetscCall(PetscSectionInvalidateMaxDof_Internal(leafSection)); 3219566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(leafSection)); 322b0c7db22SLisandro Dalcin if (hasc) { /* need to communicate bcIndices */ 323b0c7db22SLisandro Dalcin PetscSF bcSF; 324b0c7db22SLisandro Dalcin PetscInt *rOffBc; 325b0c7db22SLisandro Dalcin 3269566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lpEnd - lpStart, &rOffBc)); 327b0c7db22SLisandro Dalcin if (sub[1]) { 3289566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 3299566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 3309566063dSJacob Faibussowitsch PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->bc, rOffBc, leafSection->bc, &bcSF)); 3319566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE)); 3329566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE)); 3339566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&bcSF)); 334b0c7db22SLisandro Dalcin } 335b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) { 336b0c7db22SLisandro Dalcin if (sub[2 + f]) { 3379566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 3389566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 3399566063dSJacob Faibussowitsch PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->field[f]->bc, rOffBc, leafSection->field[f]->bc, &bcSF)); 3409566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE)); 3419566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE)); 3429566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&bcSF)); 343b0c7db22SLisandro Dalcin } 344b0c7db22SLisandro Dalcin } 3459566063dSJacob Faibussowitsch PetscCall(PetscFree(rOffBc)); 346b0c7db22SLisandro Dalcin } 3479566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&embedSF)); 3489566063dSJacob Faibussowitsch PetscCall(PetscFree(sub)); 3499566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_DistSect, sf, 0, 0, 0)); 3503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 351b0c7db22SLisandro Dalcin } 352b0c7db22SLisandro Dalcin 353b0c7db22SLisandro Dalcin /*@C 354b0c7db22SLisandro Dalcin PetscSFCreateRemoteOffsets - Create offsets for point data on remote processes 355b0c7db22SLisandro Dalcin 356c3339decSBarry Smith Collective 357b0c7db22SLisandro Dalcin 358b0c7db22SLisandro Dalcin Input Parameters: 359cab54364SBarry Smith + sf - The `PetscSF` 360b0c7db22SLisandro Dalcin . rootSection - Data layout of remote points for outgoing data (this is layout for SF roots) 361b0c7db22SLisandro Dalcin - leafSection - Data layout of local points for incoming data (this is layout for SF leaves) 362b0c7db22SLisandro Dalcin 363b0c7db22SLisandro Dalcin Output Parameter: 364b0c7db22SLisandro Dalcin . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL 365b0c7db22SLisandro Dalcin 366b0c7db22SLisandro Dalcin Level: developer 367b0c7db22SLisandro Dalcin 36823e9620eSJunchao Zhang Fortran Notes: 36923e9620eSJunchao Zhang In Fortran, use PetscSFCreateRemoteOffsetsF90() 37023e9620eSJunchao Zhang 371cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()` 372b0c7db22SLisandro Dalcin @*/ 373d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateRemoteOffsets(PetscSF sf, PetscSection rootSection, PetscSection leafSection, PetscInt **remoteOffsets) 374d71ae5a4SJacob Faibussowitsch { 375b0c7db22SLisandro Dalcin PetscSF embedSF; 376b0c7db22SLisandro Dalcin const PetscInt *indices; 377b0c7db22SLisandro Dalcin IS selected; 378b0c7db22SLisandro Dalcin PetscInt numRoots, rpStart = 0, rpEnd = 0, lpStart = 0, lpEnd = 0; 379b0c7db22SLisandro Dalcin 380b0c7db22SLisandro Dalcin PetscFunctionBegin; 381b0c7db22SLisandro Dalcin *remoteOffsets = NULL; 3829566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &numRoots, NULL, NULL, NULL)); 3833ba16761SJacob Faibussowitsch if (numRoots < 0) PetscFunctionReturn(PETSC_SUCCESS); 3849566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_RemoteOff, sf, 0, 0, 0)); 3859566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd)); 3869566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd)); 3879566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected)); 3889566063dSJacob Faibussowitsch PetscCall(ISGetIndices(selected, &indices)); 3899566063dSJacob Faibussowitsch PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF)); 3909566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(selected, &indices)); 3919566063dSJacob Faibussowitsch PetscCall(ISDestroy(&selected)); 3929566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(lpEnd - lpStart, remoteOffsets)); 3938e3a54c0SPierre Jolivet PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(*remoteOffsets, -lpStart), MPI_REPLACE)); 3948e3a54c0SPierre Jolivet PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(*remoteOffsets, -lpStart), MPI_REPLACE)); 3959566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&embedSF)); 3969566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_RemoteOff, sf, 0, 0, 0)); 3973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 398b0c7db22SLisandro Dalcin } 399b0c7db22SLisandro Dalcin 400b0c7db22SLisandro Dalcin /*@C 401cab54364SBarry Smith PetscSFCreateSectionSF - Create an expanded `PetscSF` of dofs, assuming the input `PetscSF` relates points 402b0c7db22SLisandro Dalcin 403c3339decSBarry Smith Collective 404b0c7db22SLisandro Dalcin 405b0c7db22SLisandro Dalcin Input Parameters: 406cab54364SBarry Smith + sf - The `PetscSF` 407b0c7db22SLisandro Dalcin . rootSection - Data layout of remote points for outgoing data (this is usually the serial section) 408b0c7db22SLisandro Dalcin . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL 409b0c7db22SLisandro Dalcin - leafSection - Data layout of local points for incoming data (this is the distributed section) 410b0c7db22SLisandro Dalcin 4112fe279fdSBarry Smith Output Parameter: 4122fe279fdSBarry Smith . sectionSF - The new `PetscSF` 413b0c7db22SLisandro Dalcin 414b0c7db22SLisandro Dalcin Level: advanced 415b0c7db22SLisandro Dalcin 41623e9620eSJunchao Zhang Notes: 417021e28c0SJames Wright `remoteOffsets` can be NULL if `sf` does not reference any points in leafSection 418cab54364SBarry Smith 41923e9620eSJunchao Zhang Fortran Notes: 42023e9620eSJunchao Zhang In Fortran, use PetscSFCreateSectionSFF90() 42123e9620eSJunchao Zhang 422cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()` 423b0c7db22SLisandro Dalcin @*/ 424d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateSectionSF(PetscSF sf, PetscSection rootSection, PetscInt remoteOffsets[], PetscSection leafSection, PetscSF *sectionSF) 425d71ae5a4SJacob Faibussowitsch { 426b0c7db22SLisandro Dalcin MPI_Comm comm; 427b0c7db22SLisandro Dalcin const PetscInt *localPoints; 428b0c7db22SLisandro Dalcin const PetscSFNode *remotePoints; 429b0c7db22SLisandro Dalcin PetscInt lpStart, lpEnd; 430b0c7db22SLisandro Dalcin PetscInt numRoots, numSectionRoots, numPoints, numIndices = 0; 431b0c7db22SLisandro Dalcin PetscInt *localIndices; 432b0c7db22SLisandro Dalcin PetscSFNode *remoteIndices; 433b0c7db22SLisandro Dalcin PetscInt i, ind; 434b0c7db22SLisandro Dalcin 435b0c7db22SLisandro Dalcin PetscFunctionBegin; 436b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 4374f572ea9SToby Isaac PetscAssertPointer(rootSection, 2); 4384f572ea9SToby Isaac /* Cannot check PetscAssertPointer(remoteOffsets,3) because it can be NULL if sf does not reference any points in leafSection */ 4394f572ea9SToby Isaac PetscAssertPointer(leafSection, 4); 4404f572ea9SToby Isaac PetscAssertPointer(sectionSF, 5); 4419566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 4429566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, sectionSF)); 4439566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd)); 4449566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(rootSection, &numSectionRoots)); 4459566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &numRoots, &numPoints, &localPoints, &remotePoints)); 4463ba16761SJacob Faibussowitsch if (numRoots < 0) PetscFunctionReturn(PETSC_SUCCESS); 4479566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_SectSF, sf, 0, 0, 0)); 448b0c7db22SLisandro Dalcin for (i = 0; i < numPoints; ++i) { 449b0c7db22SLisandro Dalcin PetscInt localPoint = localPoints ? localPoints[i] : i; 450b0c7db22SLisandro Dalcin PetscInt dof; 451b0c7db22SLisandro Dalcin 452b0c7db22SLisandro Dalcin if ((localPoint >= lpStart) && (localPoint < lpEnd)) { 4539566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof)); 45422a18c9fSMatthew G. Knepley numIndices += dof < 0 ? 0 : dof; 455b0c7db22SLisandro Dalcin } 456b0c7db22SLisandro Dalcin } 4579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numIndices, &localIndices)); 4589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numIndices, &remoteIndices)); 459b0c7db22SLisandro Dalcin /* Create new index graph */ 460b0c7db22SLisandro Dalcin for (i = 0, ind = 0; i < numPoints; ++i) { 461b0c7db22SLisandro Dalcin PetscInt localPoint = localPoints ? localPoints[i] : i; 462835f2295SStefano Zampini PetscInt rank = remotePoints[i].rank; 463b0c7db22SLisandro Dalcin 464b0c7db22SLisandro Dalcin if ((localPoint >= lpStart) && (localPoint < lpEnd)) { 465b0c7db22SLisandro Dalcin PetscInt remoteOffset = remoteOffsets[localPoint - lpStart]; 466b0c7db22SLisandro Dalcin PetscInt loff, dof, d; 467b0c7db22SLisandro Dalcin 4689566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(leafSection, localPoint, &loff)); 4699566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof)); 470b0c7db22SLisandro Dalcin for (d = 0; d < dof; ++d, ++ind) { 471b0c7db22SLisandro Dalcin localIndices[ind] = loff + d; 472b0c7db22SLisandro Dalcin remoteIndices[ind].rank = rank; 473b0c7db22SLisandro Dalcin remoteIndices[ind].index = remoteOffset + d; 474b0c7db22SLisandro Dalcin } 475b0c7db22SLisandro Dalcin } 476b0c7db22SLisandro Dalcin } 47708401ef6SPierre Jolivet PetscCheck(numIndices == ind, comm, PETSC_ERR_PLIB, "Inconsistency in indices, %" PetscInt_FMT " should be %" PetscInt_FMT, ind, numIndices); 4789566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*sectionSF, numSectionRoots, numIndices, localIndices, PETSC_OWN_POINTER, remoteIndices, PETSC_OWN_POINTER)); 4799566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(*sectionSF)); 4809566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_SectSF, sf, 0, 0, 0)); 4813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 482b0c7db22SLisandro Dalcin } 4838fa9e22eSVaclav Hapla 4848fa9e22eSVaclav Hapla /*@C 485cab54364SBarry Smith PetscSFCreateFromLayouts - Creates a parallel star forest mapping two `PetscLayout` objects 4868fa9e22eSVaclav Hapla 4878fa9e22eSVaclav Hapla Collective 4888fa9e22eSVaclav Hapla 4894165533cSJose E. Roman Input Parameters: 490cab54364SBarry Smith + rmap - `PetscLayout` defining the global root space 491cab54364SBarry Smith - lmap - `PetscLayout` defining the global leaf space 4928fa9e22eSVaclav Hapla 4934165533cSJose E. Roman Output Parameter: 4948fa9e22eSVaclav Hapla . sf - The parallel star forest 4958fa9e22eSVaclav Hapla 4968fa9e22eSVaclav Hapla Level: intermediate 4978fa9e22eSVaclav Hapla 498cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()`, `PetscLayoutCreate()`, `PetscSFSetGraphLayout()` 4998fa9e22eSVaclav Hapla @*/ 500d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateFromLayouts(PetscLayout rmap, PetscLayout lmap, PetscSF *sf) 501d71ae5a4SJacob Faibussowitsch { 5028fa9e22eSVaclav Hapla PetscInt i, nroots, nleaves = 0; 5038fa9e22eSVaclav Hapla PetscInt rN, lst, len; 5048fa9e22eSVaclav Hapla PetscMPIInt owner = -1; 5058fa9e22eSVaclav Hapla PetscSFNode *remote; 5068fa9e22eSVaclav Hapla MPI_Comm rcomm = rmap->comm; 5078fa9e22eSVaclav Hapla MPI_Comm lcomm = lmap->comm; 5088fa9e22eSVaclav Hapla PetscMPIInt flg; 5098fa9e22eSVaclav Hapla 5108fa9e22eSVaclav Hapla PetscFunctionBegin; 5114f572ea9SToby Isaac PetscAssertPointer(sf, 3); 51228b400f6SJacob Faibussowitsch PetscCheck(rmap->setupcalled, rcomm, PETSC_ERR_ARG_WRONGSTATE, "Root layout not setup"); 51328b400f6SJacob Faibussowitsch PetscCheck(lmap->setupcalled, lcomm, PETSC_ERR_ARG_WRONGSTATE, "Leaf layout not setup"); 5149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(rcomm, lcomm, &flg)); 515c9cc58a2SBarry Smith PetscCheck(flg == MPI_CONGRUENT || flg == MPI_IDENT, rcomm, PETSC_ERR_SUP, "cannot map two layouts with non-matching communicators"); 5169566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(rcomm, sf)); 5179566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(rmap, &nroots)); 5189566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(rmap, &rN)); 5199566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRange(lmap, &lst, &len)); 5209566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len - lst, &remote)); 5218fa9e22eSVaclav Hapla for (i = lst; i < len && i < rN; i++) { 52248a46eb9SPierre Jolivet if (owner < -1 || i >= rmap->range[owner + 1]) PetscCall(PetscLayoutFindOwner(rmap, i, &owner)); 5238fa9e22eSVaclav Hapla remote[nleaves].rank = owner; 5248fa9e22eSVaclav Hapla remote[nleaves].index = i - rmap->range[owner]; 5258fa9e22eSVaclav Hapla nleaves++; 5268fa9e22eSVaclav Hapla } 5279566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*sf, nroots, nleaves, NULL, PETSC_OWN_POINTER, remote, PETSC_COPY_VALUES)); 5289566063dSJacob Faibussowitsch PetscCall(PetscFree(remote)); 5293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5308fa9e22eSVaclav Hapla } 5318fa9e22eSVaclav Hapla 5328fa9e22eSVaclav Hapla /* TODO: handle nooffprocentries like MatZeroRowsMapLocal_Private, since this code is the same */ 533d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLayoutMapLocal(PetscLayout map, PetscInt N, const PetscInt idxs[], PetscInt *on, PetscInt **oidxs, PetscInt **ogidxs) 534d71ae5a4SJacob Faibussowitsch { 5358fa9e22eSVaclav Hapla PetscInt *owners = map->range; 5368fa9e22eSVaclav Hapla PetscInt n = map->n; 5378fa9e22eSVaclav Hapla PetscSF sf; 5388fa9e22eSVaclav Hapla PetscInt *lidxs, *work = NULL; 5398fa9e22eSVaclav Hapla PetscSFNode *ridxs; 5408fa9e22eSVaclav Hapla PetscMPIInt rank, p = 0; 5418fa9e22eSVaclav Hapla PetscInt r, len = 0; 5428fa9e22eSVaclav Hapla 5438fa9e22eSVaclav Hapla PetscFunctionBegin; 5448fa9e22eSVaclav Hapla if (on) *on = 0; /* squelch -Wmaybe-uninitialized */ 5458fa9e22eSVaclav Hapla /* Create SF where leaves are input idxs and roots are owned idxs */ 5469566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(map->comm, &rank)); 5479566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &lidxs)); 5488fa9e22eSVaclav Hapla for (r = 0; r < n; ++r) lidxs[r] = -1; 5499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N, &ridxs)); 5508fa9e22eSVaclav Hapla for (r = 0; r < N; ++r) { 5518fa9e22eSVaclav Hapla const PetscInt idx = idxs[r]; 552c9cc58a2SBarry Smith 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); 5538fa9e22eSVaclav Hapla if (idx < owners[p] || owners[p + 1] <= idx) { /* short-circuit the search if the last p owns this idx too */ 5549566063dSJacob Faibussowitsch PetscCall(PetscLayoutFindOwner(map, idx, &p)); 5558fa9e22eSVaclav Hapla } 5568fa9e22eSVaclav Hapla ridxs[r].rank = p; 5578fa9e22eSVaclav Hapla ridxs[r].index = idxs[r] - owners[p]; 5588fa9e22eSVaclav Hapla } 5599566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(map->comm, &sf)); 5609566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, n, N, NULL, PETSC_OWN_POINTER, ridxs, PETSC_OWN_POINTER)); 5619566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR)); 5629566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR)); 5638fa9e22eSVaclav Hapla if (ogidxs) { /* communicate global idxs */ 5648fa9e22eSVaclav Hapla PetscInt cum = 0, start, *work2; 5658fa9e22eSVaclav Hapla 5669566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &work)); 5679566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(N, &work2)); 5689371c9d4SSatish Balay for (r = 0; r < N; ++r) 5699371c9d4SSatish Balay if (idxs[r] >= 0) cum++; 5709566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&cum, &start, 1, MPIU_INT, MPI_SUM, map->comm)); 5718fa9e22eSVaclav Hapla start -= cum; 5728fa9e22eSVaclav Hapla cum = 0; 5739371c9d4SSatish Balay for (r = 0; r < N; ++r) 5749371c9d4SSatish Balay if (idxs[r] >= 0) work2[r] = start + cum++; 5759566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf, MPIU_INT, work2, work, MPI_REPLACE)); 5769566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf, MPIU_INT, work2, work, MPI_REPLACE)); 5779566063dSJacob Faibussowitsch PetscCall(PetscFree(work2)); 5788fa9e22eSVaclav Hapla } 5799566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf)); 5808fa9e22eSVaclav Hapla /* Compress and put in indices */ 5818fa9e22eSVaclav Hapla for (r = 0; r < n; ++r) 5828fa9e22eSVaclav Hapla if (lidxs[r] >= 0) { 5838fa9e22eSVaclav Hapla if (work) work[len] = work[r]; 5848fa9e22eSVaclav Hapla lidxs[len++] = r; 5858fa9e22eSVaclav Hapla } 5868fa9e22eSVaclav Hapla if (on) *on = len; 5878fa9e22eSVaclav Hapla if (oidxs) *oidxs = lidxs; 5888fa9e22eSVaclav Hapla if (ogidxs) *ogidxs = work; 5893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5908fa9e22eSVaclav Hapla } 591deffd5ebSksagiyam 592deffd5ebSksagiyam /*@ 593cab54364SBarry Smith PetscSFCreateByMatchingIndices - Create `PetscSF` by matching root and leaf indices 594deffd5ebSksagiyam 595deffd5ebSksagiyam Collective 596deffd5ebSksagiyam 597deffd5ebSksagiyam Input Parameters: 598cab54364SBarry Smith + layout - `PetscLayout` defining the global index space and the rank that brokers each index 599deffd5ebSksagiyam . numRootIndices - size of rootIndices 600cab54364SBarry Smith . rootIndices - `PetscInt` array of global indices of which this process requests ownership 601deffd5ebSksagiyam . rootLocalIndices - root local index permutation (NULL if no permutation) 602deffd5ebSksagiyam . rootLocalOffset - offset to be added to root local indices 603deffd5ebSksagiyam . numLeafIndices - size of leafIndices 604cab54364SBarry Smith . leafIndices - `PetscInt` array of global indices with which this process requires data associated 605deffd5ebSksagiyam . leafLocalIndices - leaf local index permutation (NULL if no permutation) 606deffd5ebSksagiyam - leafLocalOffset - offset to be added to leaf local indices 607deffd5ebSksagiyam 608d8d19677SJose E. Roman Output Parameters: 609deffd5ebSksagiyam + sfA - star forest representing the communication pattern from the layout space to the leaf space (NULL if not needed) 610deffd5ebSksagiyam - sf - star forest representing the communication pattern from the root space to the leaf space 611deffd5ebSksagiyam 612cab54364SBarry Smith Level: advanced 613cab54364SBarry Smith 614deffd5ebSksagiyam Example 1: 615cab54364SBarry Smith .vb 616cab54364SBarry Smith rank : 0 1 2 617cab54364SBarry Smith rootIndices : [1 0 2] [3] [3] 618cab54364SBarry Smith rootLocalOffset : 100 200 300 619cab54364SBarry Smith layout : [0 1] [2] [3] 620cab54364SBarry Smith leafIndices : [0] [2] [0 3] 621cab54364SBarry Smith leafLocalOffset : 400 500 600 622cab54364SBarry Smith 623cab54364SBarry Smith would build the following PetscSF 624cab54364SBarry Smith 625cab54364SBarry Smith [0] 400 <- (0,101) 626cab54364SBarry Smith [1] 500 <- (0,102) 627cab54364SBarry Smith [2] 600 <- (0,101) 628cab54364SBarry Smith [2] 601 <- (2,300) 629cab54364SBarry Smith .ve 630cab54364SBarry Smith 631deffd5ebSksagiyam Example 2: 632cab54364SBarry Smith .vb 633cab54364SBarry Smith rank : 0 1 2 634cab54364SBarry Smith rootIndices : [1 0 2] [3] [3] 635cab54364SBarry Smith rootLocalOffset : 100 200 300 636cab54364SBarry Smith layout : [0 1] [2] [3] 637cab54364SBarry Smith leafIndices : rootIndices rootIndices rootIndices 638cab54364SBarry Smith leafLocalOffset : rootLocalOffset rootLocalOffset rootLocalOffset 639cab54364SBarry Smith 640cab54364SBarry Smith would build the following PetscSF 641cab54364SBarry Smith 642cab54364SBarry Smith [1] 200 <- (2,300) 643cab54364SBarry Smith .ve 644cab54364SBarry Smith 645deffd5ebSksagiyam Example 3: 646cab54364SBarry Smith .vb 647cab54364SBarry Smith No process requests ownership of global index 1, but no process needs it. 648cab54364SBarry Smith 649cab54364SBarry Smith rank : 0 1 2 650cab54364SBarry Smith numRootIndices : 2 1 1 651cab54364SBarry Smith rootIndices : [0 2] [3] [3] 652cab54364SBarry Smith rootLocalOffset : 100 200 300 653cab54364SBarry Smith layout : [0 1] [2] [3] 654cab54364SBarry Smith numLeafIndices : 1 1 2 655cab54364SBarry Smith leafIndices : [0] [2] [0 3] 656cab54364SBarry Smith leafLocalOffset : 400 500 600 657cab54364SBarry Smith 658cab54364SBarry Smith would build the following PetscSF 659cab54364SBarry Smith 660cab54364SBarry Smith [0] 400 <- (0,100) 661cab54364SBarry Smith [1] 500 <- (0,101) 662cab54364SBarry Smith [2] 600 <- (0,100) 663cab54364SBarry Smith [2] 601 <- (2,300) 664cab54364SBarry Smith .ve 665deffd5ebSksagiyam 666deffd5ebSksagiyam Notes: 667deffd5ebSksagiyam The layout parameter represents any partitioning of [0, N), where N is the total number of global indices, and its 668cab54364SBarry Smith local size can be set to `PETSC_DECIDE`. 669cab54364SBarry Smith 670deffd5ebSksagiyam If a global index x lies in the partition owned by process i, each process whose rootIndices contains x requests 671deffd5ebSksagiyam ownership of x and sends its own rank and the local index of x to process i. 672deffd5ebSksagiyam If multiple processes request ownership of x, the one with the highest rank is to own x. 673deffd5ebSksagiyam Process i then broadcasts the ownership information, so that each process whose leafIndices contains x knows the 674deffd5ebSksagiyam ownership information of x. 675deffd5ebSksagiyam The output sf is constructed by associating each leaf point to a root point in this way. 676deffd5ebSksagiyam 677deffd5ebSksagiyam Suppose there is point data ordered according to the global indices and partitioned according to the given layout. 678cab54364SBarry Smith The optional output `PetscSF` object sfA can be used to push such data to leaf points. 679deffd5ebSksagiyam 680deffd5ebSksagiyam All indices in rootIndices and leafIndices must lie in the layout range. The union (over all processes) of rootIndices 681deffd5ebSksagiyam must cover that of leafIndices, but need not cover the entire layout. 682deffd5ebSksagiyam 683deffd5ebSksagiyam If (leafIndices, leafLocalIndices, leafLocalOffset) == (rootIndices, rootLocalIndices, rootLocalOffset), the output 684deffd5ebSksagiyam star forest is almost identity, so will only include non-trivial part of the map. 685deffd5ebSksagiyam 68638b5cf2dSJacob Faibussowitsch Developer Notes: 687deffd5ebSksagiyam Current approach of a process of the highest rank gaining the ownership may cause load imbalance; consider using 688deffd5ebSksagiyam hash(rank, root_local_index) as the bid for the ownership determination. 689deffd5ebSksagiyam 690cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()` 691deffd5ebSksagiyam @*/ 692d71ae5a4SJacob Faibussowitsch 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) 693d71ae5a4SJacob Faibussowitsch { 694deffd5ebSksagiyam MPI_Comm comm = layout->comm; 695deffd5ebSksagiyam PetscMPIInt size, rank; 696deffd5ebSksagiyam PetscSF sf1; 697deffd5ebSksagiyam PetscSFNode *owners, *buffer, *iremote; 698deffd5ebSksagiyam PetscInt *ilocal, nleaves, N, n, i; 699deffd5ebSksagiyam #if defined(PETSC_USE_DEBUG) 700deffd5ebSksagiyam PetscInt N1; 701deffd5ebSksagiyam #endif 702deffd5ebSksagiyam PetscBool flag; 703deffd5ebSksagiyam 704deffd5ebSksagiyam PetscFunctionBegin; 7054f572ea9SToby Isaac if (rootIndices) PetscAssertPointer(rootIndices, 3); 7064f572ea9SToby Isaac if (rootLocalIndices) PetscAssertPointer(rootLocalIndices, 4); 7074f572ea9SToby Isaac if (leafIndices) PetscAssertPointer(leafIndices, 7); 7084f572ea9SToby Isaac if (leafLocalIndices) PetscAssertPointer(leafLocalIndices, 8); 7094f572ea9SToby Isaac if (sfA) PetscAssertPointer(sfA, 10); 7104f572ea9SToby Isaac PetscAssertPointer(sf, 11); 71108401ef6SPierre Jolivet PetscCheck(numRootIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numRootIndices (%" PetscInt_FMT ") must be non-negative", numRootIndices); 71208401ef6SPierre Jolivet PetscCheck(numLeafIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numLeafIndices (%" PetscInt_FMT ") must be non-negative", numLeafIndices); 7139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 7149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 7159566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(layout)); 7169566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(layout, &N)); 7179566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(layout, &n)); 718deffd5ebSksagiyam flag = (PetscBool)(leafIndices == rootIndices); 719462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &flag, 1, MPIU_BOOL, MPI_LAND, comm)); 720c9cc58a2SBarry Smith PetscCheck(!flag || numLeafIndices == numRootIndices, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "leafIndices == rootIndices, but numLeafIndices (%" PetscInt_FMT ") != numRootIndices(%" PetscInt_FMT ")", numLeafIndices, numRootIndices); 721deffd5ebSksagiyam #if defined(PETSC_USE_DEBUG) 7221690c2aeSBarry Smith N1 = PETSC_INT_MIN; 7239371c9d4SSatish Balay for (i = 0; i < numRootIndices; i++) 7249371c9d4SSatish Balay if (rootIndices[i] > N1) N1 = rootIndices[i]; 725462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm)); 72608401ef6SPierre 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); 727deffd5ebSksagiyam if (!flag) { 7281690c2aeSBarry Smith N1 = PETSC_INT_MIN; 7299371c9d4SSatish Balay for (i = 0; i < numLeafIndices; i++) 7309371c9d4SSatish Balay if (leafIndices[i] > N1) N1 = leafIndices[i]; 731462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm)); 73208401ef6SPierre 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); 733deffd5ebSksagiyam } 734deffd5ebSksagiyam #endif 735deffd5ebSksagiyam /* Reduce: owners -> buffer */ 7369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &buffer)); 7379566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, &sf1)); 7389566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sf1)); 7399566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraphLayout(sf1, layout, numRootIndices, NULL, PETSC_OWN_POINTER, rootIndices)); 7409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numRootIndices, &owners)); 741deffd5ebSksagiyam for (i = 0; i < numRootIndices; ++i) { 742deffd5ebSksagiyam owners[i].rank = rank; 743deffd5ebSksagiyam owners[i].index = rootLocalOffset + (rootLocalIndices ? rootLocalIndices[i] : i); 744deffd5ebSksagiyam } 745deffd5ebSksagiyam for (i = 0; i < n; ++i) { 746deffd5ebSksagiyam buffer[i].index = -1; 747deffd5ebSksagiyam buffer[i].rank = -1; 748deffd5ebSksagiyam } 7496497c311SBarry Smith PetscCall(PetscSFReduceBegin(sf1, MPIU_SF_NODE, owners, buffer, MPI_MAXLOC)); 7506497c311SBarry Smith PetscCall(PetscSFReduceEnd(sf1, MPIU_SF_NODE, owners, buffer, MPI_MAXLOC)); 751deffd5ebSksagiyam /* Bcast: buffer -> owners */ 752deffd5ebSksagiyam if (!flag) { 753deffd5ebSksagiyam /* leafIndices is different from rootIndices */ 7549566063dSJacob Faibussowitsch PetscCall(PetscFree(owners)); 7559566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraphLayout(sf1, layout, numLeafIndices, NULL, PETSC_OWN_POINTER, leafIndices)); 7569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numLeafIndices, &owners)); 757deffd5ebSksagiyam } 7586497c311SBarry Smith PetscCall(PetscSFBcastBegin(sf1, MPIU_SF_NODE, buffer, owners, MPI_REPLACE)); 7596497c311SBarry Smith PetscCall(PetscSFBcastEnd(sf1, MPIU_SF_NODE, buffer, owners, MPI_REPLACE)); 76008401ef6SPierre 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]); 7619566063dSJacob Faibussowitsch PetscCall(PetscFree(buffer)); 7629371c9d4SSatish Balay if (sfA) { 7639371c9d4SSatish Balay *sfA = sf1; 7649371c9d4SSatish Balay } else PetscCall(PetscSFDestroy(&sf1)); 765deffd5ebSksagiyam /* Create sf */ 766deffd5ebSksagiyam if (flag && rootLocalIndices == leafLocalIndices && leafLocalOffset == rootLocalOffset) { 767deffd5ebSksagiyam /* leaf space == root space */ 7689371c9d4SSatish Balay for (i = 0, nleaves = 0; i < numLeafIndices; ++i) 7699371c9d4SSatish Balay if (owners[i].rank != rank) ++nleaves; 7709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &ilocal)); 7719566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &iremote)); 772deffd5ebSksagiyam for (i = 0, nleaves = 0; i < numLeafIndices; ++i) { 773deffd5ebSksagiyam if (owners[i].rank != rank) { 774deffd5ebSksagiyam ilocal[nleaves] = leafLocalOffset + i; 775deffd5ebSksagiyam iremote[nleaves].rank = owners[i].rank; 776deffd5ebSksagiyam iremote[nleaves].index = owners[i].index; 777deffd5ebSksagiyam ++nleaves; 778deffd5ebSksagiyam } 779deffd5ebSksagiyam } 7809566063dSJacob Faibussowitsch PetscCall(PetscFree(owners)); 781deffd5ebSksagiyam } else { 782deffd5ebSksagiyam nleaves = numLeafIndices; 7839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &ilocal)); 784ad540459SPierre Jolivet for (i = 0; i < nleaves; ++i) ilocal[i] = leafLocalOffset + (leafLocalIndices ? leafLocalIndices[i] : i); 785deffd5ebSksagiyam iremote = owners; 786deffd5ebSksagiyam } 7879566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, sf)); 7889566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(*sf)); 7899566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*sf, rootLocalOffset + numRootIndices, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER)); 7903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 791deffd5ebSksagiyam } 792fbc7033fSJed Brown 793fbc7033fSJed Brown /*@ 79453c0d4aeSBarry Smith PetscSFMerge - append/merge indices of `sfb` into `sfa`, with preference for `sfb` 795fbc7033fSJed Brown 796fbc7033fSJed Brown Collective 797fbc7033fSJed Brown 79838b5cf2dSJacob Faibussowitsch Input Parameters: 799fbc7033fSJed Brown + sfa - default `PetscSF` 800fbc7033fSJed Brown - sfb - additional edges to add/replace edges in sfa 801fbc7033fSJed Brown 80238b5cf2dSJacob Faibussowitsch Output Parameter: 803fbc7033fSJed Brown . merged - new `PetscSF` with combined edges 804fbc7033fSJed Brown 80553c0d4aeSBarry Smith Level: intermediate 80653c0d4aeSBarry Smith 80738b5cf2dSJacob Faibussowitsch .seealso: `PetscSFCompose()` 808fbc7033fSJed Brown @*/ 809fbc7033fSJed Brown PetscErrorCode PetscSFMerge(PetscSF sfa, PetscSF sfb, PetscSF *merged) 810fbc7033fSJed Brown { 811fbc7033fSJed Brown PetscInt maxleaf; 812fbc7033fSJed Brown 813fbc7033fSJed Brown PetscFunctionBegin; 814fbc7033fSJed Brown PetscValidHeaderSpecific(sfa, PETSCSF_CLASSID, 1); 815fbc7033fSJed Brown PetscValidHeaderSpecific(sfb, PETSCSF_CLASSID, 2); 816fbc7033fSJed Brown PetscCheckSameComm(sfa, 1, sfb, 2); 8174f572ea9SToby Isaac PetscAssertPointer(merged, 3); 818fbc7033fSJed Brown { 819fbc7033fSJed Brown PetscInt aleaf, bleaf; 820fbc7033fSJed Brown PetscCall(PetscSFGetLeafRange(sfa, NULL, &aleaf)); 821fbc7033fSJed Brown PetscCall(PetscSFGetLeafRange(sfb, NULL, &bleaf)); 822fbc7033fSJed Brown maxleaf = PetscMax(aleaf, bleaf) + 1; // One more than the last index 823fbc7033fSJed Brown } 824fbc7033fSJed Brown PetscInt *clocal, aroots, aleaves, broots, bleaves; 825fbc7033fSJed Brown PetscSFNode *cremote; 826fbc7033fSJed Brown const PetscInt *alocal, *blocal; 827fbc7033fSJed Brown const PetscSFNode *aremote, *bremote; 828fbc7033fSJed Brown PetscCall(PetscMalloc2(maxleaf, &clocal, maxleaf, &cremote)); 829fbc7033fSJed Brown for (PetscInt i = 0; i < maxleaf; i++) clocal[i] = -1; 830fbc7033fSJed Brown PetscCall(PetscSFGetGraph(sfa, &aroots, &aleaves, &alocal, &aremote)); 831fbc7033fSJed Brown PetscCall(PetscSFGetGraph(sfb, &broots, &bleaves, &blocal, &bremote)); 832fbc7033fSJed Brown PetscCheck(aroots == broots, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Both sfa and sfb must have the same root space"); 833fbc7033fSJed Brown for (PetscInt i = 0; i < aleaves; i++) { 834fbc7033fSJed Brown PetscInt a = alocal ? alocal[i] : i; 835fbc7033fSJed Brown clocal[a] = a; 836fbc7033fSJed Brown cremote[a] = aremote[i]; 837fbc7033fSJed Brown } 838fbc7033fSJed Brown for (PetscInt i = 0; i < bleaves; i++) { 839fbc7033fSJed Brown PetscInt b = blocal ? blocal[i] : i; 840fbc7033fSJed Brown clocal[b] = b; 841fbc7033fSJed Brown cremote[b] = bremote[i]; 842fbc7033fSJed Brown } 843fbc7033fSJed Brown PetscInt nleaves = 0; 844fbc7033fSJed Brown for (PetscInt i = 0; i < maxleaf; i++) { 845fbc7033fSJed Brown if (clocal[i] < 0) continue; 846fbc7033fSJed Brown clocal[nleaves] = clocal[i]; 847fbc7033fSJed Brown cremote[nleaves] = cremote[i]; 848fbc7033fSJed Brown nleaves++; 849fbc7033fSJed Brown } 850fbc7033fSJed Brown PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)sfa), merged)); 851fbc7033fSJed Brown PetscCall(PetscSFSetGraph(*merged, aroots, nleaves, clocal, PETSC_COPY_VALUES, cremote, PETSC_COPY_VALUES)); 852fbc7033fSJed Brown PetscCall(PetscFree2(clocal, cremote)); 8533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 854fbc7033fSJed Brown } 8550dd791a8SStefano Zampini 8560dd791a8SStefano Zampini /*@ 8570dd791a8SStefano Zampini PetscSFCreateStridedSF - Create an `PetscSF` to communicate interleaved blocks of data 8580dd791a8SStefano Zampini 8590dd791a8SStefano Zampini Collective 8600dd791a8SStefano Zampini 8610dd791a8SStefano Zampini Input Parameters: 8620dd791a8SStefano Zampini + sf - star forest 8630dd791a8SStefano Zampini . bs - stride 8640dd791a8SStefano Zampini . ldr - leading dimension of root space 8650dd791a8SStefano Zampini - ldl - leading dimension of leaf space 8660dd791a8SStefano Zampini 8670dd791a8SStefano Zampini Output Parameter: 8680dd791a8SStefano Zampini . vsf - the new `PetscSF` 8690dd791a8SStefano Zampini 8700dd791a8SStefano Zampini Level: intermediate 8710dd791a8SStefano Zampini 8720dd791a8SStefano Zampini Notes: 8730dd791a8SStefano Zampini This can be useful to perform communications on blocks of right-hand sides. For example, the calling sequence 8740dd791a8SStefano Zampini .vb 8750dd791a8SStefano Zampini c_datatype *roots, *leaves; 8760dd791a8SStefano Zampini for i in [0,bs) do 8770dd791a8SStefano Zampini PetscSFBcastBegin(sf, mpi_datatype, roots + i*ldr, leaves + i*ldl, op) 8780dd791a8SStefano Zampini PetscSFBcastEnd(sf, mpi_datatype, roots + i*ldr, leaves + i*ldl, op) 8790dd791a8SStefano Zampini .ve 8800dd791a8SStefano Zampini is equivalent to 8810dd791a8SStefano Zampini .vb 8820dd791a8SStefano Zampini c_datatype *roots, *leaves; 8830dd791a8SStefano Zampini PetscSFCreateStridedSF(sf, bs, ldr, ldl, &vsf) 8840dd791a8SStefano Zampini PetscSFBcastBegin(vsf, mpi_datatype, roots, leaves, op) 8850dd791a8SStefano Zampini PetscSFBcastEnd(vsf, mpi_datatype, roots, leaves, op) 8860dd791a8SStefano Zampini .ve 8870dd791a8SStefano Zampini 8880dd791a8SStefano Zampini Developer Notes: 8890dd791a8SStefano Zampini Should this functionality be handled with a new API instead of creating a new object? 8900dd791a8SStefano Zampini 8910dd791a8SStefano Zampini .seealso: `PetscSF`, `PetscSFCreate()`, `PetscSFSetGraph()` 8920dd791a8SStefano Zampini @*/ 8930dd791a8SStefano Zampini PetscErrorCode PetscSFCreateStridedSF(PetscSF sf, PetscInt bs, PetscInt ldr, PetscInt ldl, PetscSF *vsf) 8940dd791a8SStefano Zampini { 8950dd791a8SStefano Zampini PetscSF rankssf; 8960dd791a8SStefano Zampini const PetscSFNode *iremote, *sfrremote; 8970dd791a8SStefano Zampini PetscSFNode *viremote; 8980dd791a8SStefano Zampini const PetscInt *ilocal; 8990dd791a8SStefano Zampini PetscInt *vilocal = NULL, *ldrs; 9000dd791a8SStefano Zampini PetscInt nranks, nr, nl, vnr, vnl, maxl; 9010dd791a8SStefano Zampini PetscMPIInt rank; 9020dd791a8SStefano Zampini MPI_Comm comm; 9030dd791a8SStefano Zampini PetscSFType sftype; 9040dd791a8SStefano Zampini 9050dd791a8SStefano Zampini PetscFunctionBegin; 9060dd791a8SStefano Zampini PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 9070dd791a8SStefano Zampini PetscValidLogicalCollectiveInt(sf, bs, 2); 9080dd791a8SStefano Zampini PetscAssertPointer(vsf, 5); 9090dd791a8SStefano Zampini if (bs == 1) { 9100dd791a8SStefano Zampini PetscCall(PetscObjectReference((PetscObject)sf)); 9110dd791a8SStefano Zampini *vsf = sf; 9120dd791a8SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 9130dd791a8SStefano Zampini } 9140dd791a8SStefano Zampini PetscCall(PetscSFSetUp(sf)); 9150dd791a8SStefano Zampini PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 9160dd791a8SStefano Zampini PetscCallMPI(MPI_Comm_rank(comm, &rank)); 9170dd791a8SStefano Zampini PetscCall(PetscSFGetGraph(sf, &nr, &nl, &ilocal, &iremote)); 9180dd791a8SStefano Zampini PetscCall(PetscSFGetLeafRange(sf, NULL, &maxl)); 9190dd791a8SStefano Zampini maxl += 1; 9200dd791a8SStefano Zampini if (ldl == PETSC_DECIDE) ldl = maxl; 9210dd791a8SStefano Zampini if (ldr == PETSC_DECIDE) ldr = nr; 9220dd791a8SStefano 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); 9230dd791a8SStefano 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); 9240dd791a8SStefano Zampini vnr = nr * bs; 9250dd791a8SStefano Zampini vnl = nl * bs; 9260dd791a8SStefano Zampini PetscCall(PetscMalloc1(vnl, &viremote)); 9270dd791a8SStefano Zampini PetscCall(PetscMalloc1(vnl, &vilocal)); 9280dd791a8SStefano Zampini 9290dd791a8SStefano Zampini /* Communicate root leading dimensions to leaf ranks */ 9300dd791a8SStefano Zampini PetscCall(PetscSFGetRanksSF(sf, &rankssf)); 9310dd791a8SStefano Zampini PetscCall(PetscSFGetGraph(rankssf, NULL, &nranks, NULL, &sfrremote)); 9320dd791a8SStefano Zampini PetscCall(PetscMalloc1(nranks, &ldrs)); 9330dd791a8SStefano Zampini PetscCall(PetscSFBcastBegin(rankssf, MPIU_INT, &ldr, ldrs, MPI_REPLACE)); 9340dd791a8SStefano Zampini PetscCall(PetscSFBcastEnd(rankssf, MPIU_INT, &ldr, ldrs, MPI_REPLACE)); 9350dd791a8SStefano Zampini 9360dd791a8SStefano Zampini for (PetscInt i = 0, rold = -1, lda = -1; i < nl; i++) { 937835f2295SStefano Zampini const PetscInt r = iremote[i].rank; 9380dd791a8SStefano Zampini const PetscInt ii = iremote[i].index; 9390dd791a8SStefano Zampini 9400dd791a8SStefano Zampini if (r == rank) lda = ldr; 9410dd791a8SStefano Zampini else if (rold != r) { 9420dd791a8SStefano Zampini PetscInt j; 9430dd791a8SStefano Zampini 9440dd791a8SStefano Zampini for (j = 0; j < nranks; j++) 9450dd791a8SStefano Zampini if (sfrremote[j].rank == r) break; 946835f2295SStefano Zampini PetscCheck(j < nranks, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unable to locate neighbor rank %" PetscInt_FMT, r); 9470dd791a8SStefano Zampini lda = ldrs[j]; 9480dd791a8SStefano Zampini } 9490dd791a8SStefano Zampini rold = r; 9500dd791a8SStefano Zampini for (PetscInt v = 0; v < bs; v++) { 9510dd791a8SStefano Zampini viremote[v * nl + i].rank = r; 9520dd791a8SStefano Zampini viremote[v * nl + i].index = v * lda + ii; 9530dd791a8SStefano Zampini vilocal[v * nl + i] = v * ldl + (ilocal ? ilocal[i] : i); 9540dd791a8SStefano Zampini } 9550dd791a8SStefano Zampini } 9560dd791a8SStefano Zampini PetscCall(PetscFree(ldrs)); 9570dd791a8SStefano Zampini PetscCall(PetscSFCreate(comm, vsf)); 9580dd791a8SStefano Zampini PetscCall(PetscSFGetType(sf, &sftype)); 9590dd791a8SStefano Zampini PetscCall(PetscSFSetType(*vsf, sftype)); 9600dd791a8SStefano Zampini PetscCall(PetscSFSetGraph(*vsf, vnr, vnl, vilocal, PETSC_OWN_POINTER, viremote, PETSC_OWN_POINTER)); 9610dd791a8SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 9620dd791a8SStefano Zampini } 963