1b0c7db22SLisandro Dalcin #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/ 2b0c7db22SLisandro Dalcin #include <petsc/private/sectionimpl.h> 3b0c7db22SLisandro Dalcin 4*eb58282aSVaclav Hapla /*@ 5b0c7db22SLisandro Dalcin 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 11ddea5d60SJunchao Zhang . 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 15*eb58282aSVaclav Hapla - gremote - root vertices in global numbering corresponding to leaves in ilocal 16b0c7db22SLisandro Dalcin 17b0c7db22SLisandro Dalcin Level: intermediate 18b0c7db22SLisandro Dalcin 1986c56f52SVaclav Hapla Notes: 2086c56f52SVaclav Hapla Global indices must lie in [0, N) where N is the global size of layout. 21db2b9530SVaclav Hapla Leaf indices in ilocal get sorted; this means the user-provided array gets sorted if localmode is PETSC_OWN_POINTER. 2286c56f52SVaclav Hapla 23db2b9530SVaclav Hapla Developer Notes: 2486c56f52SVaclav Hapla Local indices which are the identity permutation in the range [0,nleaves) are discarded as they 25b0c7db22SLisandro Dalcin 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 28*eb58282aSVaclav Hapla .seealso: `PetscSFGetGraphLayout()`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()` 29b0c7db22SLisandro Dalcin @*/ 30*eb58282aSVaclav Hapla PetscErrorCode PetscSFSetGraphLayout(PetscSF sf, PetscLayout layout, PetscInt nleaves, PetscInt *ilocal, PetscCopyMode localmode, const PetscInt *gremote) { 3138a25198SStefano Zampini const PetscInt *range; 3238a25198SStefano Zampini PetscInt i, nroots, ls = -1, ln = -1; 3338a25198SStefano Zampini PetscMPIInt lr = -1; 34b0c7db22SLisandro Dalcin PetscSFNode *remote; 35b0c7db22SLisandro Dalcin 36b0c7db22SLisandro Dalcin PetscFunctionBegin; 379566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(layout, &nroots)); 389566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRanges(layout, &range)); 399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &remote)); 40*eb58282aSVaclav Hapla if (nleaves) ls = gremote[0] + 1; 41b0c7db22SLisandro Dalcin for (i = 0; i < nleaves; i++) { 42*eb58282aSVaclav Hapla const PetscInt idx = gremote[i] - ls; 4338a25198SStefano Zampini if (idx < 0 || idx >= ln) { /* short-circuit the search */ 44*eb58282aSVaclav Hapla PetscCall(PetscLayoutFindOwnerIndex(layout, gremote[i], &lr, &remote[i].index)); 4538a25198SStefano Zampini remote[i].rank = lr; 4638a25198SStefano Zampini ls = range[lr]; 4738a25198SStefano Zampini ln = range[lr + 1] - ls; 4838a25198SStefano Zampini } else { 4938a25198SStefano Zampini remote[i].rank = lr; 5038a25198SStefano Zampini remote[i].index = idx; 5138a25198SStefano Zampini } 52b0c7db22SLisandro Dalcin } 539566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, nroots, nleaves, ilocal, localmode, remote, PETSC_OWN_POINTER)); 54b0c7db22SLisandro Dalcin PetscFunctionReturn(0); 55b0c7db22SLisandro Dalcin } 56b0c7db22SLisandro Dalcin 57*eb58282aSVaclav Hapla /*@C 58*eb58282aSVaclav Hapla PetscSFGetGraphLayout - Get the global indices and PetscLayout that describe this star forest 59*eb58282aSVaclav Hapla 60*eb58282aSVaclav Hapla Collective 61*eb58282aSVaclav Hapla 62*eb58282aSVaclav Hapla Input Parameter: 63*eb58282aSVaclav Hapla . sf - star forest 64*eb58282aSVaclav Hapla 65*eb58282aSVaclav Hapla Output Parameters: 66*eb58282aSVaclav Hapla + layout - PetscLayout defining the global space for roots 67*eb58282aSVaclav Hapla . nleaves - number of leaf vertices on the current process, each of these references a root on any process 68*eb58282aSVaclav Hapla . ilocal - locations of leaves in leafdata buffers, or NULL for contiguous storage 69*eb58282aSVaclav Hapla - gremote - root vertices in global numbering corresponding to leaves in ilocal 70*eb58282aSVaclav Hapla 71*eb58282aSVaclav Hapla Level: intermediate 72*eb58282aSVaclav Hapla 73*eb58282aSVaclav Hapla Notes: 74*eb58282aSVaclav Hapla The outputs are such that passing them as inputs to `PetscSFSetGraphLayout()` would lead to the same star forest. 75*eb58282aSVaclav Hapla The outputs layout and gremote are freshly created each time this function is called, 76*eb58282aSVaclav Hapla so they need to be freed by user and cannot be qualified as const. 77*eb58282aSVaclav Hapla 78*eb58282aSVaclav Hapla 79*eb58282aSVaclav Hapla .seealso: `PetscSFSetGraphLayout()`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()` 80*eb58282aSVaclav Hapla @*/ 81*eb58282aSVaclav Hapla PetscErrorCode PetscSFGetGraphLayout(PetscSF sf, PetscLayout *layout, PetscInt *nleaves, const PetscInt *ilocal[], PetscInt *gremote[]) { 82*eb58282aSVaclav Hapla PetscInt nr, nl; 83*eb58282aSVaclav Hapla const PetscSFNode *ir; 84*eb58282aSVaclav Hapla PetscLayout lt; 85*eb58282aSVaclav Hapla 86*eb58282aSVaclav Hapla PetscFunctionBegin; 87*eb58282aSVaclav Hapla PetscCall(PetscSFGetGraph(sf, &nr, &nl, ilocal, &ir)); 88*eb58282aSVaclav Hapla PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)sf), nr, PETSC_DECIDE, 1, <)); 89*eb58282aSVaclav Hapla if (gremote) { 90*eb58282aSVaclav Hapla PetscInt i; 91*eb58282aSVaclav Hapla const PetscInt *range; 92*eb58282aSVaclav Hapla PetscInt *gr; 93*eb58282aSVaclav Hapla 94*eb58282aSVaclav Hapla PetscCall(PetscLayoutGetRanges(lt, &range)); 95*eb58282aSVaclav Hapla PetscCall(PetscMalloc1(nl, &gr)); 96*eb58282aSVaclav Hapla for (i = 0; i < nl; i++) gr[i] = range[ir[i].rank] + ir[i].index; 97*eb58282aSVaclav Hapla *gremote = gr; 98*eb58282aSVaclav Hapla } 99*eb58282aSVaclav Hapla if (nleaves) *nleaves = nl; 100*eb58282aSVaclav Hapla if (layout) *layout = lt; 101*eb58282aSVaclav Hapla else PetscCall(PetscLayoutDestroy(<)); 102*eb58282aSVaclav Hapla PetscFunctionReturn(0); 103*eb58282aSVaclav Hapla } 104*eb58282aSVaclav Hapla 105b0c7db22SLisandro Dalcin /*@ 106b0c7db22SLisandro Dalcin PetscSFSetGraphSection - Sets the PetscSF graph encoding the parallel dof overlap based upon the PetscSections 107b0c7db22SLisandro Dalcin describing the data layout. 108b0c7db22SLisandro Dalcin 109b0c7db22SLisandro Dalcin Input Parameters: 110b0c7db22SLisandro Dalcin + sf - The SF 111b0c7db22SLisandro Dalcin . localSection - PetscSection describing the local data layout 112b0c7db22SLisandro Dalcin - globalSection - PetscSection describing the global data layout 113b0c7db22SLisandro Dalcin 114b0c7db22SLisandro Dalcin Level: developer 115b0c7db22SLisandro Dalcin 116db781477SPatrick Sanan .seealso: `PetscSFSetGraph()`, `PetscSFSetGraphLayout()` 117b0c7db22SLisandro Dalcin @*/ 1189371c9d4SSatish Balay PetscErrorCode PetscSFSetGraphSection(PetscSF sf, PetscSection localSection, PetscSection globalSection) { 119b0c7db22SLisandro Dalcin MPI_Comm comm; 120b0c7db22SLisandro Dalcin PetscLayout layout; 121b0c7db22SLisandro Dalcin const PetscInt *ranges; 122b0c7db22SLisandro Dalcin PetscInt *local; 123b0c7db22SLisandro Dalcin PetscSFNode *remote; 124b0c7db22SLisandro Dalcin PetscInt pStart, pEnd, p, nroots, nleaves = 0, l; 125b0c7db22SLisandro Dalcin PetscMPIInt size, rank; 126b0c7db22SLisandro Dalcin 127b0c7db22SLisandro Dalcin PetscFunctionBegin; 128b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 129b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(localSection, PETSC_SECTION_CLASSID, 2); 130b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(globalSection, PETSC_SECTION_CLASSID, 3); 131b0c7db22SLisandro Dalcin 1329566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 1339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 1349566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1359566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(globalSection, &pStart, &pEnd)); 1369566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstrainedStorageSize(globalSection, &nroots)); 1379566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(comm, &layout)); 1389566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(layout, 1)); 1399566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetLocalSize(layout, nroots)); 1409566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(layout)); 1419566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRanges(layout, &ranges)); 142b0c7db22SLisandro Dalcin for (p = pStart; p < pEnd; ++p) { 143b0c7db22SLisandro Dalcin PetscInt gdof, gcdof; 144b0c7db22SLisandro Dalcin 1459566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(globalSection, p, &gdof)); 1469566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof)); 147c9cc58a2SBarry Smith 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)); 148b0c7db22SLisandro Dalcin nleaves += gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof; 149b0c7db22SLisandro Dalcin } 1509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &local)); 1519566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &remote)); 152b0c7db22SLisandro Dalcin for (p = pStart, l = 0; p < pEnd; ++p) { 153b0c7db22SLisandro Dalcin const PetscInt *cind; 154b0c7db22SLisandro Dalcin PetscInt dof, cdof, off, gdof, gcdof, goff, gsize, d, c; 155b0c7db22SLisandro Dalcin 1569566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(localSection, p, &dof)); 1579566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(localSection, p, &off)); 1589566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(localSection, p, &cdof)); 1599566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintIndices(localSection, p, &cind)); 1609566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(globalSection, p, &gdof)); 1619566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof)); 1629566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(globalSection, p, &goff)); 163b0c7db22SLisandro Dalcin if (!gdof) continue; /* Censored point */ 164b0c7db22SLisandro Dalcin gsize = gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof; 165b0c7db22SLisandro Dalcin if (gsize != dof - cdof) { 16608401ef6SPierre 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); 167b0c7db22SLisandro Dalcin cdof = 0; /* Ignore constraints */ 168b0c7db22SLisandro Dalcin } 169b0c7db22SLisandro Dalcin for (d = 0, c = 0; d < dof; ++d) { 1709371c9d4SSatish Balay if ((c < cdof) && (cind[c] == d)) { 1719371c9d4SSatish Balay ++c; 1729371c9d4SSatish Balay continue; 1739371c9d4SSatish Balay } 174b0c7db22SLisandro Dalcin local[l + d - c] = off + d; 175b0c7db22SLisandro Dalcin } 17608401ef6SPierre 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); 177b0c7db22SLisandro Dalcin if (gdof < 0) { 178b0c7db22SLisandro Dalcin for (d = 0; d < gsize; ++d, ++l) { 179b0c7db22SLisandro Dalcin PetscInt offset = -(goff + 1) + d, r; 180b0c7db22SLisandro Dalcin 1819566063dSJacob Faibussowitsch PetscCall(PetscFindInt(offset, size + 1, ranges, &r)); 182b0c7db22SLisandro Dalcin if (r < 0) r = -(r + 2); 18308401ef6SPierre Jolivet PetscCheck(!(r < 0) && !(r >= size), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " mapped to invalid process %" PetscInt_FMT " (%" PetscInt_FMT ", %" PetscInt_FMT ")", p, r, gdof, goff); 184b0c7db22SLisandro Dalcin remote[l].rank = r; 185b0c7db22SLisandro Dalcin remote[l].index = offset - ranges[r]; 186b0c7db22SLisandro Dalcin } 187b0c7db22SLisandro Dalcin } else { 188b0c7db22SLisandro Dalcin for (d = 0; d < gsize; ++d, ++l) { 189b0c7db22SLisandro Dalcin remote[l].rank = rank; 190b0c7db22SLisandro Dalcin remote[l].index = goff + d - ranges[rank]; 191b0c7db22SLisandro Dalcin } 192b0c7db22SLisandro Dalcin } 193b0c7db22SLisandro Dalcin } 19408401ef6SPierre Jolivet PetscCheck(l == nleaves, comm, PETSC_ERR_PLIB, "Iteration error, l %" PetscInt_FMT " != nleaves %" PetscInt_FMT, l, nleaves); 1959566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&layout)); 1969566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER)); 197b0c7db22SLisandro Dalcin PetscFunctionReturn(0); 198b0c7db22SLisandro Dalcin } 199b0c7db22SLisandro Dalcin 200b0c7db22SLisandro Dalcin /*@C 201b0c7db22SLisandro Dalcin PetscSFDistributeSection - Create a new PetscSection reorganized, moving from the root to the leaves of the SF 202b0c7db22SLisandro Dalcin 203b0c7db22SLisandro Dalcin Collective on sf 204b0c7db22SLisandro Dalcin 205b0c7db22SLisandro Dalcin Input Parameters: 206b0c7db22SLisandro Dalcin + sf - The SF 207b0c7db22SLisandro Dalcin - rootSection - Section defined on root space 208b0c7db22SLisandro Dalcin 209b0c7db22SLisandro Dalcin Output Parameters: 210b0c7db22SLisandro Dalcin + remoteOffsets - root offsets in leaf storage, or NULL 211b0c7db22SLisandro Dalcin - leafSection - Section defined on the leaf space 212b0c7db22SLisandro Dalcin 213b0c7db22SLisandro Dalcin Level: advanced 214b0c7db22SLisandro Dalcin 215db781477SPatrick Sanan .seealso: `PetscSFCreate()` 216b0c7db22SLisandro Dalcin @*/ 2179371c9d4SSatish Balay PetscErrorCode PetscSFDistributeSection(PetscSF sf, PetscSection rootSection, PetscInt **remoteOffsets, PetscSection leafSection) { 218b0c7db22SLisandro Dalcin PetscSF embedSF; 219b0c7db22SLisandro Dalcin const PetscInt *indices; 220b0c7db22SLisandro Dalcin IS selected; 221b0c7db22SLisandro Dalcin PetscInt numFields, nroots, rpStart, rpEnd, lpStart = PETSC_MAX_INT, lpEnd = -1, f, c; 222b0c7db22SLisandro Dalcin PetscBool *sub, hasc; 223b0c7db22SLisandro Dalcin 224b0c7db22SLisandro Dalcin PetscFunctionBegin; 2259566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_DistSect, sf, 0, 0, 0)); 2269566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(rootSection, &numFields)); 227029e8818Sksagiyam if (numFields) { 228029e8818Sksagiyam IS perm; 229029e8818Sksagiyam 230029e8818Sksagiyam /* PetscSectionSetNumFields() calls PetscSectionReset(), which destroys 231029e8818Sksagiyam leafSection->perm. To keep this permutation set by the user, we grab 232029e8818Sksagiyam the reference before calling PetscSectionSetNumFields() and set it 233029e8818Sksagiyam back after. */ 2349566063dSJacob Faibussowitsch PetscCall(PetscSectionGetPermutation(leafSection, &perm)); 2359566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)perm)); 2369566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(leafSection, numFields)); 2379566063dSJacob Faibussowitsch PetscCall(PetscSectionSetPermutation(leafSection, perm)); 2389566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 239029e8818Sksagiyam } 2409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numFields + 2, &sub)); 241b0c7db22SLisandro Dalcin sub[1] = rootSection->bc ? PETSC_TRUE : PETSC_FALSE; 242b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) { 2432ee82556SMatthew G. Knepley PetscSectionSym sym, dsym = NULL; 244b0c7db22SLisandro Dalcin const char *name = NULL; 245b0c7db22SLisandro Dalcin PetscInt numComp = 0; 246b0c7db22SLisandro Dalcin 247b0c7db22SLisandro Dalcin sub[2 + f] = rootSection->field[f]->bc ? PETSC_TRUE : PETSC_FALSE; 2489566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldComponents(rootSection, f, &numComp)); 2499566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldName(rootSection, f, &name)); 2509566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldSym(rootSection, f, &sym)); 2519566063dSJacob Faibussowitsch if (sym) PetscCall(PetscSectionSymDistribute(sym, sf, &dsym)); 2529566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(leafSection, f, numComp)); 2539566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldName(leafSection, f, name)); 2549566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldSym(leafSection, f, dsym)); 2559566063dSJacob Faibussowitsch PetscCall(PetscSectionSymDestroy(&dsym)); 256b0c7db22SLisandro Dalcin for (c = 0; c < rootSection->numFieldComponents[f]; ++c) { 2579566063dSJacob Faibussowitsch PetscCall(PetscSectionGetComponentName(rootSection, f, c, &name)); 2589566063dSJacob Faibussowitsch PetscCall(PetscSectionSetComponentName(leafSection, f, c, name)); 259b0c7db22SLisandro Dalcin } 260b0c7db22SLisandro Dalcin } 2619566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd)); 2629566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL)); 263b0c7db22SLisandro Dalcin rpEnd = PetscMin(rpEnd, nroots); 264b0c7db22SLisandro Dalcin rpEnd = PetscMax(rpStart, rpEnd); 265b0c7db22SLisandro Dalcin /* see if we can avoid creating the embedded SF, since it can cost more than an allreduce */ 266b0c7db22SLisandro Dalcin sub[0] = (PetscBool)(nroots != rpEnd - rpStart); 2671c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(MPI_IN_PLACE, sub, 2 + numFields, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)sf))); 268b0c7db22SLisandro Dalcin if (sub[0]) { 2699566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected)); 2709566063dSJacob Faibussowitsch PetscCall(ISGetIndices(selected, &indices)); 2719566063dSJacob Faibussowitsch PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF)); 2729566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(selected, &indices)); 2739566063dSJacob Faibussowitsch PetscCall(ISDestroy(&selected)); 274b0c7db22SLisandro Dalcin } else { 2759566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)sf)); 276b0c7db22SLisandro Dalcin embedSF = sf; 277b0c7db22SLisandro Dalcin } 2789566063dSJacob Faibussowitsch PetscCall(PetscSFGetLeafRange(embedSF, &lpStart, &lpEnd)); 279b0c7db22SLisandro Dalcin lpEnd++; 280b0c7db22SLisandro Dalcin 2819566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(leafSection, lpStart, lpEnd)); 282b0c7db22SLisandro Dalcin 283b0c7db22SLisandro Dalcin /* Constrained dof section */ 284b0c7db22SLisandro Dalcin hasc = sub[1]; 285b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) hasc = (PetscBool)(hasc || sub[2 + f]); 286b0c7db22SLisandro Dalcin 287b0c7db22SLisandro Dalcin /* Could fuse these at the cost of copies and extra allocation */ 2889566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart], MPI_REPLACE)); 2899566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart], MPI_REPLACE)); 290b0c7db22SLisandro Dalcin if (sub[1]) { 2917c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(rootSection)); 2927c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(leafSection)); 2939566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE)); 2949566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE)); 295b0c7db22SLisandro Dalcin } 296b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) { 2979566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart], MPI_REPLACE)); 2989566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart], MPI_REPLACE)); 299b0c7db22SLisandro Dalcin if (sub[2 + f]) { 3007c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(rootSection->field[f])); 3017c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(leafSection->field[f])); 3029566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE)); 3039566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE)); 304b0c7db22SLisandro Dalcin } 305b0c7db22SLisandro Dalcin } 306b0c7db22SLisandro Dalcin if (remoteOffsets) { 3079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lpEnd - lpStart, remoteOffsets)); 3089566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE)); 3099566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE)); 310b0c7db22SLisandro Dalcin } 31169c11d05SVaclav Hapla PetscCall(PetscSectionInvalidateMaxDof_Internal(leafSection)); 3129566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(leafSection)); 313b0c7db22SLisandro Dalcin if (hasc) { /* need to communicate bcIndices */ 314b0c7db22SLisandro Dalcin PetscSF bcSF; 315b0c7db22SLisandro Dalcin PetscInt *rOffBc; 316b0c7db22SLisandro Dalcin 3179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lpEnd - lpStart, &rOffBc)); 318b0c7db22SLisandro Dalcin if (sub[1]) { 3199566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 3209566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 3219566063dSJacob Faibussowitsch PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->bc, rOffBc, leafSection->bc, &bcSF)); 3229566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE)); 3239566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE)); 3249566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&bcSF)); 325b0c7db22SLisandro Dalcin } 326b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) { 327b0c7db22SLisandro Dalcin if (sub[2 + f]) { 3289566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 3299566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 3309566063dSJacob Faibussowitsch PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->field[f]->bc, rOffBc, leafSection->field[f]->bc, &bcSF)); 3319566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE)); 3329566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE)); 3339566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&bcSF)); 334b0c7db22SLisandro Dalcin } 335b0c7db22SLisandro Dalcin } 3369566063dSJacob Faibussowitsch PetscCall(PetscFree(rOffBc)); 337b0c7db22SLisandro Dalcin } 3389566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&embedSF)); 3399566063dSJacob Faibussowitsch PetscCall(PetscFree(sub)); 3409566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_DistSect, sf, 0, 0, 0)); 341b0c7db22SLisandro Dalcin PetscFunctionReturn(0); 342b0c7db22SLisandro Dalcin } 343b0c7db22SLisandro Dalcin 344b0c7db22SLisandro Dalcin /*@C 345b0c7db22SLisandro Dalcin PetscSFCreateRemoteOffsets - Create offsets for point data on remote processes 346b0c7db22SLisandro Dalcin 347b0c7db22SLisandro Dalcin Collective on sf 348b0c7db22SLisandro Dalcin 349b0c7db22SLisandro Dalcin Input Parameters: 350b0c7db22SLisandro Dalcin + sf - The SF 351b0c7db22SLisandro Dalcin . rootSection - Data layout of remote points for outgoing data (this is layout for SF roots) 352b0c7db22SLisandro Dalcin - leafSection - Data layout of local points for incoming data (this is layout for SF leaves) 353b0c7db22SLisandro Dalcin 354b0c7db22SLisandro Dalcin Output Parameter: 355b0c7db22SLisandro Dalcin . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL 356b0c7db22SLisandro Dalcin 357b0c7db22SLisandro Dalcin Level: developer 358b0c7db22SLisandro Dalcin 359db781477SPatrick Sanan .seealso: `PetscSFCreate()` 360b0c7db22SLisandro Dalcin @*/ 3619371c9d4SSatish Balay PetscErrorCode PetscSFCreateRemoteOffsets(PetscSF sf, PetscSection rootSection, PetscSection leafSection, PetscInt **remoteOffsets) { 362b0c7db22SLisandro Dalcin PetscSF embedSF; 363b0c7db22SLisandro Dalcin const PetscInt *indices; 364b0c7db22SLisandro Dalcin IS selected; 365b0c7db22SLisandro Dalcin PetscInt numRoots, rpStart = 0, rpEnd = 0, lpStart = 0, lpEnd = 0; 366b0c7db22SLisandro Dalcin 367b0c7db22SLisandro Dalcin PetscFunctionBegin; 368b0c7db22SLisandro Dalcin *remoteOffsets = NULL; 3699566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &numRoots, NULL, NULL, NULL)); 370b0c7db22SLisandro Dalcin if (numRoots < 0) PetscFunctionReturn(0); 3719566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_RemoteOff, sf, 0, 0, 0)); 3729566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd)); 3739566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd)); 3749566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected)); 3759566063dSJacob Faibussowitsch PetscCall(ISGetIndices(selected, &indices)); 3769566063dSJacob Faibussowitsch PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF)); 3779566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(selected, &indices)); 3789566063dSJacob Faibussowitsch PetscCall(ISDestroy(&selected)); 3799566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(lpEnd - lpStart, remoteOffsets)); 3809566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE)); 3819566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE)); 3829566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&embedSF)); 3839566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_RemoteOff, sf, 0, 0, 0)); 384b0c7db22SLisandro Dalcin PetscFunctionReturn(0); 385b0c7db22SLisandro Dalcin } 386b0c7db22SLisandro Dalcin 387b0c7db22SLisandro Dalcin /*@C 388b0c7db22SLisandro Dalcin PetscSFCreateSectionSF - Create an expanded SF of dofs, assuming the input SF relates points 389b0c7db22SLisandro Dalcin 390b0c7db22SLisandro Dalcin Collective on sf 391b0c7db22SLisandro Dalcin 392b0c7db22SLisandro Dalcin Input Parameters: 393b0c7db22SLisandro Dalcin + sf - The SF 394b0c7db22SLisandro Dalcin . rootSection - Data layout of remote points for outgoing data (this is usually the serial section) 395b0c7db22SLisandro Dalcin . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL 396b0c7db22SLisandro Dalcin - leafSection - Data layout of local points for incoming data (this is the distributed section) 397b0c7db22SLisandro Dalcin 398b0c7db22SLisandro Dalcin Output Parameters: 399b0c7db22SLisandro Dalcin - sectionSF - The new SF 400b0c7db22SLisandro Dalcin 401b0c7db22SLisandro Dalcin Note: Either rootSection or remoteOffsets can be specified 402b0c7db22SLisandro Dalcin 403b0c7db22SLisandro Dalcin Level: advanced 404b0c7db22SLisandro Dalcin 405db781477SPatrick Sanan .seealso: `PetscSFCreate()` 406b0c7db22SLisandro Dalcin @*/ 4079371c9d4SSatish Balay PetscErrorCode PetscSFCreateSectionSF(PetscSF sf, PetscSection rootSection, PetscInt remoteOffsets[], PetscSection leafSection, PetscSF *sectionSF) { 408b0c7db22SLisandro Dalcin MPI_Comm comm; 409b0c7db22SLisandro Dalcin const PetscInt *localPoints; 410b0c7db22SLisandro Dalcin const PetscSFNode *remotePoints; 411b0c7db22SLisandro Dalcin PetscInt lpStart, lpEnd; 412b0c7db22SLisandro Dalcin PetscInt numRoots, numSectionRoots, numPoints, numIndices = 0; 413b0c7db22SLisandro Dalcin PetscInt *localIndices; 414b0c7db22SLisandro Dalcin PetscSFNode *remoteIndices; 415b0c7db22SLisandro Dalcin PetscInt i, ind; 416b0c7db22SLisandro Dalcin 417b0c7db22SLisandro Dalcin PetscFunctionBegin; 418b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 419b0c7db22SLisandro Dalcin PetscValidPointer(rootSection, 2); 420b0c7db22SLisandro Dalcin /* Cannot check PetscValidIntPointer(remoteOffsets,3) because it can be NULL if sf does not reference any points in leafSection */ 421b0c7db22SLisandro Dalcin PetscValidPointer(leafSection, 4); 422b0c7db22SLisandro Dalcin PetscValidPointer(sectionSF, 5); 4239566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 4249566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, sectionSF)); 4259566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd)); 4269566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(rootSection, &numSectionRoots)); 4279566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &numRoots, &numPoints, &localPoints, &remotePoints)); 428b0c7db22SLisandro Dalcin if (numRoots < 0) PetscFunctionReturn(0); 4299566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_SectSF, sf, 0, 0, 0)); 430b0c7db22SLisandro Dalcin for (i = 0; i < numPoints; ++i) { 431b0c7db22SLisandro Dalcin PetscInt localPoint = localPoints ? localPoints[i] : i; 432b0c7db22SLisandro Dalcin PetscInt dof; 433b0c7db22SLisandro Dalcin 434b0c7db22SLisandro Dalcin if ((localPoint >= lpStart) && (localPoint < lpEnd)) { 4359566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof)); 43622a18c9fSMatthew G. Knepley numIndices += dof < 0 ? 0 : dof; 437b0c7db22SLisandro Dalcin } 438b0c7db22SLisandro Dalcin } 4399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numIndices, &localIndices)); 4409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numIndices, &remoteIndices)); 441b0c7db22SLisandro Dalcin /* Create new index graph */ 442b0c7db22SLisandro Dalcin for (i = 0, ind = 0; i < numPoints; ++i) { 443b0c7db22SLisandro Dalcin PetscInt localPoint = localPoints ? localPoints[i] : i; 444b0c7db22SLisandro Dalcin PetscInt rank = remotePoints[i].rank; 445b0c7db22SLisandro Dalcin 446b0c7db22SLisandro Dalcin if ((localPoint >= lpStart) && (localPoint < lpEnd)) { 447b0c7db22SLisandro Dalcin PetscInt remoteOffset = remoteOffsets[localPoint - lpStart]; 448b0c7db22SLisandro Dalcin PetscInt loff, dof, d; 449b0c7db22SLisandro Dalcin 4509566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(leafSection, localPoint, &loff)); 4519566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof)); 452b0c7db22SLisandro Dalcin for (d = 0; d < dof; ++d, ++ind) { 453b0c7db22SLisandro Dalcin localIndices[ind] = loff + d; 454b0c7db22SLisandro Dalcin remoteIndices[ind].rank = rank; 455b0c7db22SLisandro Dalcin remoteIndices[ind].index = remoteOffset + d; 456b0c7db22SLisandro Dalcin } 457b0c7db22SLisandro Dalcin } 458b0c7db22SLisandro Dalcin } 45908401ef6SPierre Jolivet PetscCheck(numIndices == ind, comm, PETSC_ERR_PLIB, "Inconsistency in indices, %" PetscInt_FMT " should be %" PetscInt_FMT, ind, numIndices); 4609566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*sectionSF, numSectionRoots, numIndices, localIndices, PETSC_OWN_POINTER, remoteIndices, PETSC_OWN_POINTER)); 4619566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(*sectionSF)); 4629566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_SectSF, sf, 0, 0, 0)); 463b0c7db22SLisandro Dalcin PetscFunctionReturn(0); 464b0c7db22SLisandro Dalcin } 4658fa9e22eSVaclav Hapla 4668fa9e22eSVaclav Hapla /*@C 4678fa9e22eSVaclav Hapla PetscSFCreateFromLayouts - Creates a parallel star forest mapping two PetscLayout objects 4688fa9e22eSVaclav Hapla 4698fa9e22eSVaclav Hapla Collective 4708fa9e22eSVaclav Hapla 4714165533cSJose E. Roman Input Parameters: 4728fa9e22eSVaclav Hapla + rmap - PetscLayout defining the global root space 4738fa9e22eSVaclav Hapla - lmap - PetscLayout defining the global leaf space 4748fa9e22eSVaclav Hapla 4754165533cSJose E. Roman Output Parameter: 4768fa9e22eSVaclav Hapla . sf - The parallel star forest 4778fa9e22eSVaclav Hapla 4788fa9e22eSVaclav Hapla Level: intermediate 4798fa9e22eSVaclav Hapla 480db781477SPatrick Sanan .seealso: `PetscSFCreate()`, `PetscLayoutCreate()`, `PetscSFSetGraphLayout()` 4818fa9e22eSVaclav Hapla @*/ 4829371c9d4SSatish Balay PetscErrorCode PetscSFCreateFromLayouts(PetscLayout rmap, PetscLayout lmap, PetscSF *sf) { 4838fa9e22eSVaclav Hapla PetscInt i, nroots, nleaves = 0; 4848fa9e22eSVaclav Hapla PetscInt rN, lst, len; 4858fa9e22eSVaclav Hapla PetscMPIInt owner = -1; 4868fa9e22eSVaclav Hapla PetscSFNode *remote; 4878fa9e22eSVaclav Hapla MPI_Comm rcomm = rmap->comm; 4888fa9e22eSVaclav Hapla MPI_Comm lcomm = lmap->comm; 4898fa9e22eSVaclav Hapla PetscMPIInt flg; 4908fa9e22eSVaclav Hapla 4918fa9e22eSVaclav Hapla PetscFunctionBegin; 4928fa9e22eSVaclav Hapla PetscValidPointer(sf, 3); 49328b400f6SJacob Faibussowitsch PetscCheck(rmap->setupcalled, rcomm, PETSC_ERR_ARG_WRONGSTATE, "Root layout not setup"); 49428b400f6SJacob Faibussowitsch PetscCheck(lmap->setupcalled, lcomm, PETSC_ERR_ARG_WRONGSTATE, "Leaf layout not setup"); 4959566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(rcomm, lcomm, &flg)); 496c9cc58a2SBarry Smith PetscCheck(flg == MPI_CONGRUENT || flg == MPI_IDENT, rcomm, PETSC_ERR_SUP, "cannot map two layouts with non-matching communicators"); 4979566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(rcomm, sf)); 4989566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(rmap, &nroots)); 4999566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(rmap, &rN)); 5009566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRange(lmap, &lst, &len)); 5019566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len - lst, &remote)); 5028fa9e22eSVaclav Hapla for (i = lst; i < len && i < rN; i++) { 50348a46eb9SPierre Jolivet if (owner < -1 || i >= rmap->range[owner + 1]) PetscCall(PetscLayoutFindOwner(rmap, i, &owner)); 5048fa9e22eSVaclav Hapla remote[nleaves].rank = owner; 5058fa9e22eSVaclav Hapla remote[nleaves].index = i - rmap->range[owner]; 5068fa9e22eSVaclav Hapla nleaves++; 5078fa9e22eSVaclav Hapla } 5089566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*sf, nroots, nleaves, NULL, PETSC_OWN_POINTER, remote, PETSC_COPY_VALUES)); 5099566063dSJacob Faibussowitsch PetscCall(PetscFree(remote)); 5108fa9e22eSVaclav Hapla PetscFunctionReturn(0); 5118fa9e22eSVaclav Hapla } 5128fa9e22eSVaclav Hapla 5138fa9e22eSVaclav Hapla /* TODO: handle nooffprocentries like MatZeroRowsMapLocal_Private, since this code is the same */ 5149371c9d4SSatish Balay PetscErrorCode PetscLayoutMapLocal(PetscLayout map, PetscInt N, const PetscInt idxs[], PetscInt *on, PetscInt **oidxs, PetscInt **ogidxs) { 5158fa9e22eSVaclav Hapla PetscInt *owners = map->range; 5168fa9e22eSVaclav Hapla PetscInt n = map->n; 5178fa9e22eSVaclav Hapla PetscSF sf; 5188fa9e22eSVaclav Hapla PetscInt *lidxs, *work = NULL; 5198fa9e22eSVaclav Hapla PetscSFNode *ridxs; 5208fa9e22eSVaclav Hapla PetscMPIInt rank, p = 0; 5218fa9e22eSVaclav Hapla PetscInt r, len = 0; 5228fa9e22eSVaclav Hapla 5238fa9e22eSVaclav Hapla PetscFunctionBegin; 5248fa9e22eSVaclav Hapla if (on) *on = 0; /* squelch -Wmaybe-uninitialized */ 5258fa9e22eSVaclav Hapla /* Create SF where leaves are input idxs and roots are owned idxs */ 5269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(map->comm, &rank)); 5279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &lidxs)); 5288fa9e22eSVaclav Hapla for (r = 0; r < n; ++r) lidxs[r] = -1; 5299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N, &ridxs)); 5308fa9e22eSVaclav Hapla for (r = 0; r < N; ++r) { 5318fa9e22eSVaclav Hapla const PetscInt idx = idxs[r]; 532c9cc58a2SBarry 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); 5338fa9e22eSVaclav Hapla if (idx < owners[p] || owners[p + 1] <= idx) { /* short-circuit the search if the last p owns this idx too */ 5349566063dSJacob Faibussowitsch PetscCall(PetscLayoutFindOwner(map, idx, &p)); 5358fa9e22eSVaclav Hapla } 5368fa9e22eSVaclav Hapla ridxs[r].rank = p; 5378fa9e22eSVaclav Hapla ridxs[r].index = idxs[r] - owners[p]; 5388fa9e22eSVaclav Hapla } 5399566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(map->comm, &sf)); 5409566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, n, N, NULL, PETSC_OWN_POINTER, ridxs, PETSC_OWN_POINTER)); 5419566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR)); 5429566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR)); 5438fa9e22eSVaclav Hapla if (ogidxs) { /* communicate global idxs */ 5448fa9e22eSVaclav Hapla PetscInt cum = 0, start, *work2; 5458fa9e22eSVaclav Hapla 5469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &work)); 5479566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(N, &work2)); 5489371c9d4SSatish Balay for (r = 0; r < N; ++r) 5499371c9d4SSatish Balay if (idxs[r] >= 0) cum++; 5509566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&cum, &start, 1, MPIU_INT, MPI_SUM, map->comm)); 5518fa9e22eSVaclav Hapla start -= cum; 5528fa9e22eSVaclav Hapla cum = 0; 5539371c9d4SSatish Balay for (r = 0; r < N; ++r) 5549371c9d4SSatish Balay if (idxs[r] >= 0) work2[r] = start + cum++; 5559566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf, MPIU_INT, work2, work, MPI_REPLACE)); 5569566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf, MPIU_INT, work2, work, MPI_REPLACE)); 5579566063dSJacob Faibussowitsch PetscCall(PetscFree(work2)); 5588fa9e22eSVaclav Hapla } 5599566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf)); 5608fa9e22eSVaclav Hapla /* Compress and put in indices */ 5618fa9e22eSVaclav Hapla for (r = 0; r < n; ++r) 5628fa9e22eSVaclav Hapla if (lidxs[r] >= 0) { 5638fa9e22eSVaclav Hapla if (work) work[len] = work[r]; 5648fa9e22eSVaclav Hapla lidxs[len++] = r; 5658fa9e22eSVaclav Hapla } 5668fa9e22eSVaclav Hapla if (on) *on = len; 5678fa9e22eSVaclav Hapla if (oidxs) *oidxs = lidxs; 5688fa9e22eSVaclav Hapla if (ogidxs) *ogidxs = work; 5698fa9e22eSVaclav Hapla PetscFunctionReturn(0); 5708fa9e22eSVaclav Hapla } 571deffd5ebSksagiyam 572deffd5ebSksagiyam /*@ 573deffd5ebSksagiyam PetscSFCreateByMatchingIndices - Create SF by matching root and leaf indices 574deffd5ebSksagiyam 575deffd5ebSksagiyam Collective 576deffd5ebSksagiyam 577deffd5ebSksagiyam Input Parameters: 578deffd5ebSksagiyam + layout - PetscLayout defining the global index space and the rank that brokers each index 579deffd5ebSksagiyam . numRootIndices - size of rootIndices 580deffd5ebSksagiyam . rootIndices - PetscInt array of global indices of which this process requests ownership 581deffd5ebSksagiyam . rootLocalIndices - root local index permutation (NULL if no permutation) 582deffd5ebSksagiyam . rootLocalOffset - offset to be added to root local indices 583deffd5ebSksagiyam . numLeafIndices - size of leafIndices 584deffd5ebSksagiyam . leafIndices - PetscInt array of global indices with which this process requires data associated 585deffd5ebSksagiyam . leafLocalIndices - leaf local index permutation (NULL if no permutation) 586deffd5ebSksagiyam - leafLocalOffset - offset to be added to leaf local indices 587deffd5ebSksagiyam 588d8d19677SJose E. Roman Output Parameters: 589deffd5ebSksagiyam + sfA - star forest representing the communication pattern from the layout space to the leaf space (NULL if not needed) 590deffd5ebSksagiyam - sf - star forest representing the communication pattern from the root space to the leaf space 591deffd5ebSksagiyam 592deffd5ebSksagiyam Example 1: 593deffd5ebSksagiyam $ 594deffd5ebSksagiyam $ rank : 0 1 2 595deffd5ebSksagiyam $ rootIndices : [1 0 2] [3] [3] 596deffd5ebSksagiyam $ rootLocalOffset : 100 200 300 597deffd5ebSksagiyam $ layout : [0 1] [2] [3] 598deffd5ebSksagiyam $ leafIndices : [0] [2] [0 3] 599deffd5ebSksagiyam $ leafLocalOffset : 400 500 600 600deffd5ebSksagiyam $ 6017ef5781fSVaclav Hapla would build the following SF 602deffd5ebSksagiyam $ 603deffd5ebSksagiyam $ [0] 400 <- (0,101) 604deffd5ebSksagiyam $ [1] 500 <- (0,102) 605deffd5ebSksagiyam $ [2] 600 <- (0,101) 606deffd5ebSksagiyam $ [2] 601 <- (2,300) 607deffd5ebSksagiyam $ 608deffd5ebSksagiyam Example 2: 609deffd5ebSksagiyam $ 610deffd5ebSksagiyam $ rank : 0 1 2 611deffd5ebSksagiyam $ rootIndices : [1 0 2] [3] [3] 612deffd5ebSksagiyam $ rootLocalOffset : 100 200 300 613deffd5ebSksagiyam $ layout : [0 1] [2] [3] 614deffd5ebSksagiyam $ leafIndices : rootIndices rootIndices rootIndices 615deffd5ebSksagiyam $ leafLocalOffset : rootLocalOffset rootLocalOffset rootLocalOffset 616deffd5ebSksagiyam $ 6177ef5781fSVaclav Hapla would build the following SF 618deffd5ebSksagiyam $ 619deffd5ebSksagiyam $ [1] 200 <- (2,300) 620deffd5ebSksagiyam $ 621deffd5ebSksagiyam Example 3: 622deffd5ebSksagiyam $ 623deffd5ebSksagiyam $ No process requests ownership of global index 1, but no process needs it. 624deffd5ebSksagiyam $ 625deffd5ebSksagiyam $ rank : 0 1 2 626deffd5ebSksagiyam $ numRootIndices : 2 1 1 627deffd5ebSksagiyam $ rootIndices : [0 2] [3] [3] 628deffd5ebSksagiyam $ rootLocalOffset : 100 200 300 629deffd5ebSksagiyam $ layout : [0 1] [2] [3] 630deffd5ebSksagiyam $ numLeafIndices : 1 1 2 631deffd5ebSksagiyam $ leafIndices : [0] [2] [0 3] 632deffd5ebSksagiyam $ leafLocalOffset : 400 500 600 633deffd5ebSksagiyam $ 6347ef5781fSVaclav Hapla would build the following SF 635deffd5ebSksagiyam $ 636deffd5ebSksagiyam $ [0] 400 <- (0,100) 637deffd5ebSksagiyam $ [1] 500 <- (0,101) 638deffd5ebSksagiyam $ [2] 600 <- (0,100) 639deffd5ebSksagiyam $ [2] 601 <- (2,300) 640deffd5ebSksagiyam $ 641deffd5ebSksagiyam 642deffd5ebSksagiyam Notes: 643deffd5ebSksagiyam The layout parameter represents any partitioning of [0, N), where N is the total number of global indices, and its 644deffd5ebSksagiyam local size can be set to PETSC_DECIDE. 645deffd5ebSksagiyam If a global index x lies in the partition owned by process i, each process whose rootIndices contains x requests 646deffd5ebSksagiyam ownership of x and sends its own rank and the local index of x to process i. 647deffd5ebSksagiyam If multiple processes request ownership of x, the one with the highest rank is to own x. 648deffd5ebSksagiyam Process i then broadcasts the ownership information, so that each process whose leafIndices contains x knows the 649deffd5ebSksagiyam ownership information of x. 650deffd5ebSksagiyam The output sf is constructed by associating each leaf point to a root point in this way. 651deffd5ebSksagiyam 652deffd5ebSksagiyam Suppose there is point data ordered according to the global indices and partitioned according to the given layout. 653deffd5ebSksagiyam The optional output PetscSF object sfA can be used to push such data to leaf points. 654deffd5ebSksagiyam 655deffd5ebSksagiyam All indices in rootIndices and leafIndices must lie in the layout range. The union (over all processes) of rootIndices 656deffd5ebSksagiyam must cover that of leafIndices, but need not cover the entire layout. 657deffd5ebSksagiyam 658deffd5ebSksagiyam If (leafIndices, leafLocalIndices, leafLocalOffset) == (rootIndices, rootLocalIndices, rootLocalOffset), the output 659deffd5ebSksagiyam star forest is almost identity, so will only include non-trivial part of the map. 660deffd5ebSksagiyam 661deffd5ebSksagiyam Developer Notes: 662deffd5ebSksagiyam Current approach of a process of the highest rank gaining the ownership may cause load imbalance; consider using 663deffd5ebSksagiyam hash(rank, root_local_index) as the bid for the ownership determination. 664deffd5ebSksagiyam 665deffd5ebSksagiyam Level: advanced 666deffd5ebSksagiyam 667db781477SPatrick Sanan .seealso: `PetscSFCreate()` 668deffd5ebSksagiyam @*/ 6699371c9d4SSatish Balay 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) { 670deffd5ebSksagiyam MPI_Comm comm = layout->comm; 671deffd5ebSksagiyam PetscMPIInt size, rank; 672deffd5ebSksagiyam PetscSF sf1; 673deffd5ebSksagiyam PetscSFNode *owners, *buffer, *iremote; 674deffd5ebSksagiyam PetscInt *ilocal, nleaves, N, n, i; 675deffd5ebSksagiyam #if defined(PETSC_USE_DEBUG) 676deffd5ebSksagiyam PetscInt N1; 677deffd5ebSksagiyam #endif 678deffd5ebSksagiyam PetscBool flag; 679deffd5ebSksagiyam 680deffd5ebSksagiyam PetscFunctionBegin; 681deffd5ebSksagiyam if (rootIndices) PetscValidIntPointer(rootIndices, 3); 682deffd5ebSksagiyam if (rootLocalIndices) PetscValidIntPointer(rootLocalIndices, 4); 683deffd5ebSksagiyam if (leafIndices) PetscValidIntPointer(leafIndices, 7); 684deffd5ebSksagiyam if (leafLocalIndices) PetscValidIntPointer(leafLocalIndices, 8); 685deffd5ebSksagiyam if (sfA) PetscValidPointer(sfA, 10); 686deffd5ebSksagiyam PetscValidPointer(sf, 11); 68708401ef6SPierre Jolivet PetscCheck(numRootIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numRootIndices (%" PetscInt_FMT ") must be non-negative", numRootIndices); 68808401ef6SPierre Jolivet PetscCheck(numLeafIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numLeafIndices (%" PetscInt_FMT ") must be non-negative", numLeafIndices); 6899566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 6909566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 6919566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(layout)); 6929566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(layout, &N)); 6939566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(layout, &n)); 694deffd5ebSksagiyam flag = (PetscBool)(leafIndices == rootIndices); 6951c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &flag, 1, MPIU_BOOL, MPI_LAND, comm)); 696c9cc58a2SBarry Smith PetscCheck(!flag || numLeafIndices == numRootIndices, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "leafIndices == rootIndices, but numLeafIndices (%" PetscInt_FMT ") != numRootIndices(%" PetscInt_FMT ")", numLeafIndices, numRootIndices); 697deffd5ebSksagiyam #if defined(PETSC_USE_DEBUG) 698deffd5ebSksagiyam N1 = PETSC_MIN_INT; 6999371c9d4SSatish Balay for (i = 0; i < numRootIndices; i++) 7009371c9d4SSatish Balay if (rootIndices[i] > N1) N1 = rootIndices[i]; 7019566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm)); 70208401ef6SPierre 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); 703deffd5ebSksagiyam if (!flag) { 704deffd5ebSksagiyam N1 = PETSC_MIN_INT; 7059371c9d4SSatish Balay for (i = 0; i < numLeafIndices; i++) 7069371c9d4SSatish Balay if (leafIndices[i] > N1) N1 = leafIndices[i]; 7079566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm)); 70808401ef6SPierre 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); 709deffd5ebSksagiyam } 710deffd5ebSksagiyam #endif 711deffd5ebSksagiyam /* Reduce: owners -> buffer */ 7129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &buffer)); 7139566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, &sf1)); 7149566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sf1)); 7159566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraphLayout(sf1, layout, numRootIndices, NULL, PETSC_OWN_POINTER, rootIndices)); 7169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numRootIndices, &owners)); 717deffd5ebSksagiyam for (i = 0; i < numRootIndices; ++i) { 718deffd5ebSksagiyam owners[i].rank = rank; 719deffd5ebSksagiyam owners[i].index = rootLocalOffset + (rootLocalIndices ? rootLocalIndices[i] : i); 720deffd5ebSksagiyam } 721deffd5ebSksagiyam for (i = 0; i < n; ++i) { 722deffd5ebSksagiyam buffer[i].index = -1; 723deffd5ebSksagiyam buffer[i].rank = -1; 724deffd5ebSksagiyam } 7259566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC)); 7269566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC)); 727deffd5ebSksagiyam /* Bcast: buffer -> owners */ 728deffd5ebSksagiyam if (!flag) { 729deffd5ebSksagiyam /* leafIndices is different from rootIndices */ 7309566063dSJacob Faibussowitsch PetscCall(PetscFree(owners)); 7319566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraphLayout(sf1, layout, numLeafIndices, NULL, PETSC_OWN_POINTER, leafIndices)); 7329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numLeafIndices, &owners)); 733deffd5ebSksagiyam } 7349566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE)); 7359566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE)); 73608401ef6SPierre 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]); 7379566063dSJacob Faibussowitsch PetscCall(PetscFree(buffer)); 7389371c9d4SSatish Balay if (sfA) { 7399371c9d4SSatish Balay *sfA = sf1; 7409371c9d4SSatish Balay } else PetscCall(PetscSFDestroy(&sf1)); 741deffd5ebSksagiyam /* Create sf */ 742deffd5ebSksagiyam if (flag && rootLocalIndices == leafLocalIndices && leafLocalOffset == rootLocalOffset) { 743deffd5ebSksagiyam /* leaf space == root space */ 7449371c9d4SSatish Balay for (i = 0, nleaves = 0; i < numLeafIndices; ++i) 7459371c9d4SSatish Balay if (owners[i].rank != rank) ++nleaves; 7469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &ilocal)); 7479566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &iremote)); 748deffd5ebSksagiyam for (i = 0, nleaves = 0; i < numLeafIndices; ++i) { 749deffd5ebSksagiyam if (owners[i].rank != rank) { 750deffd5ebSksagiyam ilocal[nleaves] = leafLocalOffset + i; 751deffd5ebSksagiyam iremote[nleaves].rank = owners[i].rank; 752deffd5ebSksagiyam iremote[nleaves].index = owners[i].index; 753deffd5ebSksagiyam ++nleaves; 754deffd5ebSksagiyam } 755deffd5ebSksagiyam } 7569566063dSJacob Faibussowitsch PetscCall(PetscFree(owners)); 757deffd5ebSksagiyam } else { 758deffd5ebSksagiyam nleaves = numLeafIndices; 7599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &ilocal)); 760ad540459SPierre Jolivet for (i = 0; i < nleaves; ++i) ilocal[i] = leafLocalOffset + (leafLocalIndices ? leafLocalIndices[i] : i); 761deffd5ebSksagiyam iremote = owners; 762deffd5ebSksagiyam } 7639566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, sf)); 7649566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(*sf)); 7659566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*sf, rootLocalOffset + numRootIndices, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER)); 766deffd5ebSksagiyam PetscFunctionReturn(0); 767deffd5ebSksagiyam } 768