1b0c7db22SLisandro Dalcin #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/ 2b0c7db22SLisandro Dalcin #include <petsc/private/sectionimpl.h> 3b0c7db22SLisandro Dalcin 4eb58282aSVaclav 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 15eb58282aSVaclav 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 28eb58282aSVaclav Hapla .seealso: `PetscSFGetGraphLayout()`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()` 29b0c7db22SLisandro Dalcin @*/ 30*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetGraphLayout(PetscSF sf, PetscLayout layout, PetscInt nleaves, PetscInt *ilocal, PetscCopyMode localmode, const PetscInt *gremote) 31*d71ae5a4SJacob 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; 389566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(layout, &nroots)); 399566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRanges(layout, &range)); 409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &remote)); 41eb58282aSVaclav Hapla if (nleaves) ls = gremote[0] + 1; 42b0c7db22SLisandro Dalcin for (i = 0; i < nleaves; i++) { 43eb58282aSVaclav Hapla const PetscInt idx = gremote[i] - ls; 4438a25198SStefano Zampini if (idx < 0 || idx >= ln) { /* short-circuit the search */ 45eb58282aSVaclav Hapla PetscCall(PetscLayoutFindOwnerIndex(layout, gremote[i], &lr, &remote[i].index)); 4638a25198SStefano Zampini remote[i].rank = lr; 4738a25198SStefano Zampini ls = range[lr]; 4838a25198SStefano Zampini ln = range[lr + 1] - ls; 4938a25198SStefano Zampini } else { 5038a25198SStefano Zampini remote[i].rank = lr; 5138a25198SStefano Zampini remote[i].index = idx; 5238a25198SStefano Zampini } 53b0c7db22SLisandro Dalcin } 549566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, nroots, nleaves, ilocal, localmode, remote, PETSC_OWN_POINTER)); 55b0c7db22SLisandro Dalcin PetscFunctionReturn(0); 56b0c7db22SLisandro Dalcin } 57b0c7db22SLisandro Dalcin 58eb58282aSVaclav Hapla /*@C 59eb58282aSVaclav Hapla PetscSFGetGraphLayout - Get the global indices and PetscLayout that describe this star forest 60eb58282aSVaclav Hapla 61eb58282aSVaclav Hapla Collective 62eb58282aSVaclav Hapla 63eb58282aSVaclav Hapla Input Parameter: 64eb58282aSVaclav Hapla . sf - star forest 65eb58282aSVaclav Hapla 66eb58282aSVaclav Hapla Output Parameters: 67eb58282aSVaclav Hapla + layout - PetscLayout defining the global space for roots 68eb58282aSVaclav Hapla . nleaves - number of leaf vertices on the current process, each of these references a root on any process 69eb58282aSVaclav Hapla . ilocal - locations of leaves in leafdata buffers, or NULL for contiguous storage 70eb58282aSVaclav Hapla - gremote - root vertices in global numbering corresponding to leaves in ilocal 71eb58282aSVaclav Hapla 72eb58282aSVaclav Hapla Level: intermediate 73eb58282aSVaclav Hapla 74eb58282aSVaclav Hapla Notes: 75eb58282aSVaclav Hapla The outputs are such that passing them as inputs to `PetscSFSetGraphLayout()` would lead to the same star forest. 76eb58282aSVaclav Hapla The outputs layout and gremote are freshly created each time this function is called, 77eb58282aSVaclav Hapla so they need to be freed by user and cannot be qualified as const. 78eb58282aSVaclav Hapla 79eb58282aSVaclav Hapla 80eb58282aSVaclav Hapla .seealso: `PetscSFSetGraphLayout()`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()` 81eb58282aSVaclav Hapla @*/ 82*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetGraphLayout(PetscSF sf, PetscLayout *layout, PetscInt *nleaves, const PetscInt *ilocal[], PetscInt *gremote[]) 83*d71ae5a4SJacob Faibussowitsch { 84eb58282aSVaclav Hapla PetscInt nr, nl; 85eb58282aSVaclav Hapla const PetscSFNode *ir; 86eb58282aSVaclav Hapla PetscLayout lt; 87eb58282aSVaclav Hapla 88eb58282aSVaclav Hapla PetscFunctionBegin; 89eb58282aSVaclav Hapla PetscCall(PetscSFGetGraph(sf, &nr, &nl, ilocal, &ir)); 90eb58282aSVaclav Hapla PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)sf), nr, PETSC_DECIDE, 1, <)); 91eb58282aSVaclav Hapla if (gremote) { 92eb58282aSVaclav Hapla PetscInt i; 93eb58282aSVaclav Hapla const PetscInt *range; 94eb58282aSVaclav Hapla PetscInt *gr; 95eb58282aSVaclav Hapla 96eb58282aSVaclav Hapla PetscCall(PetscLayoutGetRanges(lt, &range)); 97eb58282aSVaclav Hapla PetscCall(PetscMalloc1(nl, &gr)); 98eb58282aSVaclav Hapla for (i = 0; i < nl; i++) gr[i] = range[ir[i].rank] + ir[i].index; 99eb58282aSVaclav Hapla *gremote = gr; 100eb58282aSVaclav Hapla } 101eb58282aSVaclav Hapla if (nleaves) *nleaves = nl; 102eb58282aSVaclav Hapla if (layout) *layout = lt; 103eb58282aSVaclav Hapla else PetscCall(PetscLayoutDestroy(<)); 104eb58282aSVaclav Hapla PetscFunctionReturn(0); 105eb58282aSVaclav Hapla } 106eb58282aSVaclav Hapla 107b0c7db22SLisandro Dalcin /*@ 108b0c7db22SLisandro Dalcin PetscSFSetGraphSection - Sets the PetscSF graph encoding the parallel dof overlap based upon the PetscSections 109b0c7db22SLisandro Dalcin describing the data layout. 110b0c7db22SLisandro Dalcin 111b0c7db22SLisandro Dalcin Input Parameters: 112b0c7db22SLisandro Dalcin + sf - The SF 113b0c7db22SLisandro Dalcin . localSection - PetscSection describing the local data layout 114b0c7db22SLisandro Dalcin - globalSection - PetscSection describing the global data layout 115b0c7db22SLisandro Dalcin 116b0c7db22SLisandro Dalcin Level: developer 117b0c7db22SLisandro Dalcin 118db781477SPatrick Sanan .seealso: `PetscSFSetGraph()`, `PetscSFSetGraphLayout()` 119b0c7db22SLisandro Dalcin @*/ 120*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetGraphSection(PetscSF sf, PetscSection localSection, PetscSection globalSection) 121*d71ae5a4SJacob 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)); 150c9cc58a2SBarry 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)); 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) { 182b0c7db22SLisandro Dalcin PetscInt offset = -(goff + 1) + d, r; 183b0c7db22SLisandro Dalcin 1849566063dSJacob Faibussowitsch PetscCall(PetscFindInt(offset, size + 1, ranges, &r)); 185b0c7db22SLisandro Dalcin if (r < 0) r = -(r + 2); 18608401ef6SPierre 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); 187b0c7db22SLisandro Dalcin remote[l].rank = r; 188b0c7db22SLisandro Dalcin remote[l].index = offset - ranges[r]; 189b0c7db22SLisandro Dalcin } 190b0c7db22SLisandro Dalcin } else { 191b0c7db22SLisandro Dalcin for (d = 0; d < gsize; ++d, ++l) { 192b0c7db22SLisandro Dalcin remote[l].rank = rank; 193b0c7db22SLisandro Dalcin remote[l].index = goff + d - ranges[rank]; 194b0c7db22SLisandro Dalcin } 195b0c7db22SLisandro Dalcin } 196b0c7db22SLisandro Dalcin } 19708401ef6SPierre Jolivet PetscCheck(l == nleaves, comm, PETSC_ERR_PLIB, "Iteration error, l %" PetscInt_FMT " != nleaves %" PetscInt_FMT, l, nleaves); 1989566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&layout)); 1999566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER)); 200b0c7db22SLisandro Dalcin PetscFunctionReturn(0); 201b0c7db22SLisandro Dalcin } 202b0c7db22SLisandro Dalcin 203b0c7db22SLisandro Dalcin /*@C 204b0c7db22SLisandro Dalcin PetscSFDistributeSection - Create a new PetscSection reorganized, moving from the root to the leaves of the SF 205b0c7db22SLisandro Dalcin 206b0c7db22SLisandro Dalcin Collective on sf 207b0c7db22SLisandro Dalcin 208b0c7db22SLisandro Dalcin Input Parameters: 209b0c7db22SLisandro Dalcin + sf - The SF 210b0c7db22SLisandro Dalcin - rootSection - Section defined on root space 211b0c7db22SLisandro Dalcin 212b0c7db22SLisandro Dalcin Output Parameters: 213b0c7db22SLisandro Dalcin + remoteOffsets - root offsets in leaf storage, or NULL 214b0c7db22SLisandro Dalcin - leafSection - Section defined on the leaf space 215b0c7db22SLisandro Dalcin 216b0c7db22SLisandro Dalcin Level: advanced 217b0c7db22SLisandro Dalcin 218db781477SPatrick Sanan .seealso: `PetscSFCreate()` 219b0c7db22SLisandro Dalcin @*/ 220*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFDistributeSection(PetscSF sf, PetscSection rootSection, PetscInt **remoteOffsets, PetscSection leafSection) 221*d71ae5a4SJacob Faibussowitsch { 222b0c7db22SLisandro Dalcin PetscSF embedSF; 223b0c7db22SLisandro Dalcin const PetscInt *indices; 224b0c7db22SLisandro Dalcin IS selected; 225b0c7db22SLisandro Dalcin PetscInt numFields, nroots, rpStart, rpEnd, lpStart = PETSC_MAX_INT, lpEnd = -1, f, c; 226b0c7db22SLisandro Dalcin PetscBool *sub, hasc; 227b0c7db22SLisandro Dalcin 228b0c7db22SLisandro Dalcin PetscFunctionBegin; 2299566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_DistSect, sf, 0, 0, 0)); 2309566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(rootSection, &numFields)); 231029e8818Sksagiyam if (numFields) { 232029e8818Sksagiyam IS perm; 233029e8818Sksagiyam 234029e8818Sksagiyam /* PetscSectionSetNumFields() calls PetscSectionReset(), which destroys 235029e8818Sksagiyam leafSection->perm. To keep this permutation set by the user, we grab 236029e8818Sksagiyam the reference before calling PetscSectionSetNumFields() and set it 237029e8818Sksagiyam back after. */ 2389566063dSJacob Faibussowitsch PetscCall(PetscSectionGetPermutation(leafSection, &perm)); 2399566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)perm)); 2409566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(leafSection, numFields)); 2419566063dSJacob Faibussowitsch PetscCall(PetscSectionSetPermutation(leafSection, perm)); 2429566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 243029e8818Sksagiyam } 2449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numFields + 2, &sub)); 245b0c7db22SLisandro Dalcin sub[1] = rootSection->bc ? PETSC_TRUE : PETSC_FALSE; 246b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) { 2472ee82556SMatthew G. Knepley PetscSectionSym sym, dsym = NULL; 248b0c7db22SLisandro Dalcin const char *name = NULL; 249b0c7db22SLisandro Dalcin PetscInt numComp = 0; 250b0c7db22SLisandro Dalcin 251b0c7db22SLisandro Dalcin sub[2 + f] = rootSection->field[f]->bc ? PETSC_TRUE : PETSC_FALSE; 2529566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldComponents(rootSection, f, &numComp)); 2539566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldName(rootSection, f, &name)); 2549566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldSym(rootSection, f, &sym)); 2559566063dSJacob Faibussowitsch if (sym) PetscCall(PetscSectionSymDistribute(sym, sf, &dsym)); 2569566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(leafSection, f, numComp)); 2579566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldName(leafSection, f, name)); 2589566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldSym(leafSection, f, dsym)); 2599566063dSJacob Faibussowitsch PetscCall(PetscSectionSymDestroy(&dsym)); 260b0c7db22SLisandro Dalcin for (c = 0; c < rootSection->numFieldComponents[f]; ++c) { 2619566063dSJacob Faibussowitsch PetscCall(PetscSectionGetComponentName(rootSection, f, c, &name)); 2629566063dSJacob Faibussowitsch PetscCall(PetscSectionSetComponentName(leafSection, f, c, name)); 263b0c7db22SLisandro Dalcin } 264b0c7db22SLisandro Dalcin } 2659566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd)); 2669566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL)); 267b0c7db22SLisandro Dalcin rpEnd = PetscMin(rpEnd, nroots); 268b0c7db22SLisandro Dalcin rpEnd = PetscMax(rpStart, rpEnd); 269b0c7db22SLisandro Dalcin /* see if we can avoid creating the embedded SF, since it can cost more than an allreduce */ 270b0c7db22SLisandro Dalcin sub[0] = (PetscBool)(nroots != rpEnd - rpStart); 2711c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(MPI_IN_PLACE, sub, 2 + numFields, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)sf))); 272b0c7db22SLisandro Dalcin if (sub[0]) { 2739566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected)); 2749566063dSJacob Faibussowitsch PetscCall(ISGetIndices(selected, &indices)); 2759566063dSJacob Faibussowitsch PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF)); 2769566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(selected, &indices)); 2779566063dSJacob Faibussowitsch PetscCall(ISDestroy(&selected)); 278b0c7db22SLisandro Dalcin } else { 2799566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)sf)); 280b0c7db22SLisandro Dalcin embedSF = sf; 281b0c7db22SLisandro Dalcin } 2829566063dSJacob Faibussowitsch PetscCall(PetscSFGetLeafRange(embedSF, &lpStart, &lpEnd)); 283b0c7db22SLisandro Dalcin lpEnd++; 284b0c7db22SLisandro Dalcin 2859566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(leafSection, lpStart, lpEnd)); 286b0c7db22SLisandro Dalcin 287b0c7db22SLisandro Dalcin /* Constrained dof section */ 288b0c7db22SLisandro Dalcin hasc = sub[1]; 289b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) hasc = (PetscBool)(hasc || sub[2 + f]); 290b0c7db22SLisandro Dalcin 291b0c7db22SLisandro Dalcin /* Could fuse these at the cost of copies and extra allocation */ 2929566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart], MPI_REPLACE)); 2939566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart], MPI_REPLACE)); 294b0c7db22SLisandro Dalcin if (sub[1]) { 2957c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(rootSection)); 2967c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(leafSection)); 2979566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE)); 2989566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE)); 299b0c7db22SLisandro Dalcin } 300b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) { 3019566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart], MPI_REPLACE)); 3029566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart], MPI_REPLACE)); 303b0c7db22SLisandro Dalcin if (sub[2 + f]) { 3047c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(rootSection->field[f])); 3057c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(leafSection->field[f])); 3069566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE)); 3079566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE)); 308b0c7db22SLisandro Dalcin } 309b0c7db22SLisandro Dalcin } 310b0c7db22SLisandro Dalcin if (remoteOffsets) { 3119566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lpEnd - lpStart, remoteOffsets)); 3129566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE)); 3139566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE)); 314b0c7db22SLisandro Dalcin } 31569c11d05SVaclav Hapla PetscCall(PetscSectionInvalidateMaxDof_Internal(leafSection)); 3169566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(leafSection)); 317b0c7db22SLisandro Dalcin if (hasc) { /* need to communicate bcIndices */ 318b0c7db22SLisandro Dalcin PetscSF bcSF; 319b0c7db22SLisandro Dalcin PetscInt *rOffBc; 320b0c7db22SLisandro Dalcin 3219566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lpEnd - lpStart, &rOffBc)); 322b0c7db22SLisandro Dalcin if (sub[1]) { 3239566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 3249566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 3259566063dSJacob Faibussowitsch PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->bc, rOffBc, leafSection->bc, &bcSF)); 3269566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE)); 3279566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE)); 3289566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&bcSF)); 329b0c7db22SLisandro Dalcin } 330b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) { 331b0c7db22SLisandro Dalcin if (sub[2 + f]) { 3329566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 3339566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 3349566063dSJacob Faibussowitsch PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->field[f]->bc, rOffBc, leafSection->field[f]->bc, &bcSF)); 3359566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE)); 3369566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE)); 3379566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&bcSF)); 338b0c7db22SLisandro Dalcin } 339b0c7db22SLisandro Dalcin } 3409566063dSJacob Faibussowitsch PetscCall(PetscFree(rOffBc)); 341b0c7db22SLisandro Dalcin } 3429566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&embedSF)); 3439566063dSJacob Faibussowitsch PetscCall(PetscFree(sub)); 3449566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_DistSect, sf, 0, 0, 0)); 345b0c7db22SLisandro Dalcin PetscFunctionReturn(0); 346b0c7db22SLisandro Dalcin } 347b0c7db22SLisandro Dalcin 348b0c7db22SLisandro Dalcin /*@C 349b0c7db22SLisandro Dalcin PetscSFCreateRemoteOffsets - Create offsets for point data on remote processes 350b0c7db22SLisandro Dalcin 351b0c7db22SLisandro Dalcin Collective on sf 352b0c7db22SLisandro Dalcin 353b0c7db22SLisandro Dalcin Input Parameters: 354b0c7db22SLisandro Dalcin + sf - The SF 355b0c7db22SLisandro Dalcin . rootSection - Data layout of remote points for outgoing data (this is layout for SF roots) 356b0c7db22SLisandro Dalcin - leafSection - Data layout of local points for incoming data (this is layout for SF leaves) 357b0c7db22SLisandro Dalcin 358b0c7db22SLisandro Dalcin Output Parameter: 359b0c7db22SLisandro Dalcin . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL 360b0c7db22SLisandro Dalcin 361b0c7db22SLisandro Dalcin Level: developer 362b0c7db22SLisandro Dalcin 363db781477SPatrick Sanan .seealso: `PetscSFCreate()` 364b0c7db22SLisandro Dalcin @*/ 365*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateRemoteOffsets(PetscSF sf, PetscSection rootSection, PetscSection leafSection, PetscInt **remoteOffsets) 366*d71ae5a4SJacob Faibussowitsch { 367b0c7db22SLisandro Dalcin PetscSF embedSF; 368b0c7db22SLisandro Dalcin const PetscInt *indices; 369b0c7db22SLisandro Dalcin IS selected; 370b0c7db22SLisandro Dalcin PetscInt numRoots, rpStart = 0, rpEnd = 0, lpStart = 0, lpEnd = 0; 371b0c7db22SLisandro Dalcin 372b0c7db22SLisandro Dalcin PetscFunctionBegin; 373b0c7db22SLisandro Dalcin *remoteOffsets = NULL; 3749566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &numRoots, NULL, NULL, NULL)); 375b0c7db22SLisandro Dalcin if (numRoots < 0) PetscFunctionReturn(0); 3769566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_RemoteOff, sf, 0, 0, 0)); 3779566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd)); 3789566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd)); 3799566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected)); 3809566063dSJacob Faibussowitsch PetscCall(ISGetIndices(selected, &indices)); 3819566063dSJacob Faibussowitsch PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF)); 3829566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(selected, &indices)); 3839566063dSJacob Faibussowitsch PetscCall(ISDestroy(&selected)); 3849566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(lpEnd - lpStart, remoteOffsets)); 3859566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE)); 3869566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE)); 3879566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&embedSF)); 3889566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_RemoteOff, sf, 0, 0, 0)); 389b0c7db22SLisandro Dalcin PetscFunctionReturn(0); 390b0c7db22SLisandro Dalcin } 391b0c7db22SLisandro Dalcin 392b0c7db22SLisandro Dalcin /*@C 393b0c7db22SLisandro Dalcin PetscSFCreateSectionSF - Create an expanded SF of dofs, assuming the input SF relates points 394b0c7db22SLisandro Dalcin 395b0c7db22SLisandro Dalcin Collective on sf 396b0c7db22SLisandro Dalcin 397b0c7db22SLisandro Dalcin Input Parameters: 398b0c7db22SLisandro Dalcin + sf - The SF 399b0c7db22SLisandro Dalcin . rootSection - Data layout of remote points for outgoing data (this is usually the serial section) 400b0c7db22SLisandro Dalcin . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL 401b0c7db22SLisandro Dalcin - leafSection - Data layout of local points for incoming data (this is the distributed section) 402b0c7db22SLisandro Dalcin 403b0c7db22SLisandro Dalcin Output Parameters: 404b0c7db22SLisandro Dalcin - sectionSF - The new SF 405b0c7db22SLisandro Dalcin 406b0c7db22SLisandro Dalcin Note: Either rootSection or remoteOffsets can be specified 407b0c7db22SLisandro Dalcin 408b0c7db22SLisandro Dalcin Level: advanced 409b0c7db22SLisandro Dalcin 410db781477SPatrick Sanan .seealso: `PetscSFCreate()` 411b0c7db22SLisandro Dalcin @*/ 412*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateSectionSF(PetscSF sf, PetscSection rootSection, PetscInt remoteOffsets[], PetscSection leafSection, PetscSF *sectionSF) 413*d71ae5a4SJacob Faibussowitsch { 414b0c7db22SLisandro Dalcin MPI_Comm comm; 415b0c7db22SLisandro Dalcin const PetscInt *localPoints; 416b0c7db22SLisandro Dalcin const PetscSFNode *remotePoints; 417b0c7db22SLisandro Dalcin PetscInt lpStart, lpEnd; 418b0c7db22SLisandro Dalcin PetscInt numRoots, numSectionRoots, numPoints, numIndices = 0; 419b0c7db22SLisandro Dalcin PetscInt *localIndices; 420b0c7db22SLisandro Dalcin PetscSFNode *remoteIndices; 421b0c7db22SLisandro Dalcin PetscInt i, ind; 422b0c7db22SLisandro Dalcin 423b0c7db22SLisandro Dalcin PetscFunctionBegin; 424b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 425b0c7db22SLisandro Dalcin PetscValidPointer(rootSection, 2); 426b0c7db22SLisandro Dalcin /* Cannot check PetscValidIntPointer(remoteOffsets,3) because it can be NULL if sf does not reference any points in leafSection */ 427b0c7db22SLisandro Dalcin PetscValidPointer(leafSection, 4); 428b0c7db22SLisandro Dalcin PetscValidPointer(sectionSF, 5); 4299566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 4309566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, sectionSF)); 4319566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd)); 4329566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(rootSection, &numSectionRoots)); 4339566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &numRoots, &numPoints, &localPoints, &remotePoints)); 434b0c7db22SLisandro Dalcin if (numRoots < 0) PetscFunctionReturn(0); 4359566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_SectSF, sf, 0, 0, 0)); 436b0c7db22SLisandro Dalcin for (i = 0; i < numPoints; ++i) { 437b0c7db22SLisandro Dalcin PetscInt localPoint = localPoints ? localPoints[i] : i; 438b0c7db22SLisandro Dalcin PetscInt dof; 439b0c7db22SLisandro Dalcin 440b0c7db22SLisandro Dalcin if ((localPoint >= lpStart) && (localPoint < lpEnd)) { 4419566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof)); 44222a18c9fSMatthew G. Knepley numIndices += dof < 0 ? 0 : dof; 443b0c7db22SLisandro Dalcin } 444b0c7db22SLisandro Dalcin } 4459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numIndices, &localIndices)); 4469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numIndices, &remoteIndices)); 447b0c7db22SLisandro Dalcin /* Create new index graph */ 448b0c7db22SLisandro Dalcin for (i = 0, ind = 0; i < numPoints; ++i) { 449b0c7db22SLisandro Dalcin PetscInt localPoint = localPoints ? localPoints[i] : i; 450b0c7db22SLisandro Dalcin PetscInt rank = remotePoints[i].rank; 451b0c7db22SLisandro Dalcin 452b0c7db22SLisandro Dalcin if ((localPoint >= lpStart) && (localPoint < lpEnd)) { 453b0c7db22SLisandro Dalcin PetscInt remoteOffset = remoteOffsets[localPoint - lpStart]; 454b0c7db22SLisandro Dalcin PetscInt loff, dof, d; 455b0c7db22SLisandro Dalcin 4569566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(leafSection, localPoint, &loff)); 4579566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof)); 458b0c7db22SLisandro Dalcin for (d = 0; d < dof; ++d, ++ind) { 459b0c7db22SLisandro Dalcin localIndices[ind] = loff + d; 460b0c7db22SLisandro Dalcin remoteIndices[ind].rank = rank; 461b0c7db22SLisandro Dalcin remoteIndices[ind].index = remoteOffset + d; 462b0c7db22SLisandro Dalcin } 463b0c7db22SLisandro Dalcin } 464b0c7db22SLisandro Dalcin } 46508401ef6SPierre Jolivet PetscCheck(numIndices == ind, comm, PETSC_ERR_PLIB, "Inconsistency in indices, %" PetscInt_FMT " should be %" PetscInt_FMT, ind, numIndices); 4669566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*sectionSF, numSectionRoots, numIndices, localIndices, PETSC_OWN_POINTER, remoteIndices, PETSC_OWN_POINTER)); 4679566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(*sectionSF)); 4689566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_SectSF, sf, 0, 0, 0)); 469b0c7db22SLisandro Dalcin PetscFunctionReturn(0); 470b0c7db22SLisandro Dalcin } 4718fa9e22eSVaclav Hapla 4728fa9e22eSVaclav Hapla /*@C 4738fa9e22eSVaclav Hapla PetscSFCreateFromLayouts - Creates a parallel star forest mapping two PetscLayout objects 4748fa9e22eSVaclav Hapla 4758fa9e22eSVaclav Hapla Collective 4768fa9e22eSVaclav Hapla 4774165533cSJose E. Roman Input Parameters: 4788fa9e22eSVaclav Hapla + rmap - PetscLayout defining the global root space 4798fa9e22eSVaclav Hapla - lmap - PetscLayout defining the global leaf space 4808fa9e22eSVaclav Hapla 4814165533cSJose E. Roman Output Parameter: 4828fa9e22eSVaclav Hapla . sf - The parallel star forest 4838fa9e22eSVaclav Hapla 4848fa9e22eSVaclav Hapla Level: intermediate 4858fa9e22eSVaclav Hapla 486db781477SPatrick Sanan .seealso: `PetscSFCreate()`, `PetscLayoutCreate()`, `PetscSFSetGraphLayout()` 4878fa9e22eSVaclav Hapla @*/ 488*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateFromLayouts(PetscLayout rmap, PetscLayout lmap, PetscSF *sf) 489*d71ae5a4SJacob Faibussowitsch { 4908fa9e22eSVaclav Hapla PetscInt i, nroots, nleaves = 0; 4918fa9e22eSVaclav Hapla PetscInt rN, lst, len; 4928fa9e22eSVaclav Hapla PetscMPIInt owner = -1; 4938fa9e22eSVaclav Hapla PetscSFNode *remote; 4948fa9e22eSVaclav Hapla MPI_Comm rcomm = rmap->comm; 4958fa9e22eSVaclav Hapla MPI_Comm lcomm = lmap->comm; 4968fa9e22eSVaclav Hapla PetscMPIInt flg; 4978fa9e22eSVaclav Hapla 4988fa9e22eSVaclav Hapla PetscFunctionBegin; 4998fa9e22eSVaclav Hapla PetscValidPointer(sf, 3); 50028b400f6SJacob Faibussowitsch PetscCheck(rmap->setupcalled, rcomm, PETSC_ERR_ARG_WRONGSTATE, "Root layout not setup"); 50128b400f6SJacob Faibussowitsch PetscCheck(lmap->setupcalled, lcomm, PETSC_ERR_ARG_WRONGSTATE, "Leaf layout not setup"); 5029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(rcomm, lcomm, &flg)); 503c9cc58a2SBarry Smith PetscCheck(flg == MPI_CONGRUENT || flg == MPI_IDENT, rcomm, PETSC_ERR_SUP, "cannot map two layouts with non-matching communicators"); 5049566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(rcomm, sf)); 5059566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(rmap, &nroots)); 5069566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(rmap, &rN)); 5079566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRange(lmap, &lst, &len)); 5089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len - lst, &remote)); 5098fa9e22eSVaclav Hapla for (i = lst; i < len && i < rN; i++) { 51048a46eb9SPierre Jolivet if (owner < -1 || i >= rmap->range[owner + 1]) PetscCall(PetscLayoutFindOwner(rmap, i, &owner)); 5118fa9e22eSVaclav Hapla remote[nleaves].rank = owner; 5128fa9e22eSVaclav Hapla remote[nleaves].index = i - rmap->range[owner]; 5138fa9e22eSVaclav Hapla nleaves++; 5148fa9e22eSVaclav Hapla } 5159566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*sf, nroots, nleaves, NULL, PETSC_OWN_POINTER, remote, PETSC_COPY_VALUES)); 5169566063dSJacob Faibussowitsch PetscCall(PetscFree(remote)); 5178fa9e22eSVaclav Hapla PetscFunctionReturn(0); 5188fa9e22eSVaclav Hapla } 5198fa9e22eSVaclav Hapla 5208fa9e22eSVaclav Hapla /* TODO: handle nooffprocentries like MatZeroRowsMapLocal_Private, since this code is the same */ 521*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLayoutMapLocal(PetscLayout map, PetscInt N, const PetscInt idxs[], PetscInt *on, PetscInt **oidxs, PetscInt **ogidxs) 522*d71ae5a4SJacob Faibussowitsch { 5238fa9e22eSVaclav Hapla PetscInt *owners = map->range; 5248fa9e22eSVaclav Hapla PetscInt n = map->n; 5258fa9e22eSVaclav Hapla PetscSF sf; 5268fa9e22eSVaclav Hapla PetscInt *lidxs, *work = NULL; 5278fa9e22eSVaclav Hapla PetscSFNode *ridxs; 5288fa9e22eSVaclav Hapla PetscMPIInt rank, p = 0; 5298fa9e22eSVaclav Hapla PetscInt r, len = 0; 5308fa9e22eSVaclav Hapla 5318fa9e22eSVaclav Hapla PetscFunctionBegin; 5328fa9e22eSVaclav Hapla if (on) *on = 0; /* squelch -Wmaybe-uninitialized */ 5338fa9e22eSVaclav Hapla /* Create SF where leaves are input idxs and roots are owned idxs */ 5349566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(map->comm, &rank)); 5359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &lidxs)); 5368fa9e22eSVaclav Hapla for (r = 0; r < n; ++r) lidxs[r] = -1; 5379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N, &ridxs)); 5388fa9e22eSVaclav Hapla for (r = 0; r < N; ++r) { 5398fa9e22eSVaclav Hapla const PetscInt idx = idxs[r]; 540c9cc58a2SBarry 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); 5418fa9e22eSVaclav Hapla if (idx < owners[p] || owners[p + 1] <= idx) { /* short-circuit the search if the last p owns this idx too */ 5429566063dSJacob Faibussowitsch PetscCall(PetscLayoutFindOwner(map, idx, &p)); 5438fa9e22eSVaclav Hapla } 5448fa9e22eSVaclav Hapla ridxs[r].rank = p; 5458fa9e22eSVaclav Hapla ridxs[r].index = idxs[r] - owners[p]; 5468fa9e22eSVaclav Hapla } 5479566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(map->comm, &sf)); 5489566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, n, N, NULL, PETSC_OWN_POINTER, ridxs, PETSC_OWN_POINTER)); 5499566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR)); 5509566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR)); 5518fa9e22eSVaclav Hapla if (ogidxs) { /* communicate global idxs */ 5528fa9e22eSVaclav Hapla PetscInt cum = 0, start, *work2; 5538fa9e22eSVaclav Hapla 5549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &work)); 5559566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(N, &work2)); 5569371c9d4SSatish Balay for (r = 0; r < N; ++r) 5579371c9d4SSatish Balay if (idxs[r] >= 0) cum++; 5589566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&cum, &start, 1, MPIU_INT, MPI_SUM, map->comm)); 5598fa9e22eSVaclav Hapla start -= cum; 5608fa9e22eSVaclav Hapla cum = 0; 5619371c9d4SSatish Balay for (r = 0; r < N; ++r) 5629371c9d4SSatish Balay if (idxs[r] >= 0) work2[r] = start + cum++; 5639566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf, MPIU_INT, work2, work, MPI_REPLACE)); 5649566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf, MPIU_INT, work2, work, MPI_REPLACE)); 5659566063dSJacob Faibussowitsch PetscCall(PetscFree(work2)); 5668fa9e22eSVaclav Hapla } 5679566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf)); 5688fa9e22eSVaclav Hapla /* Compress and put in indices */ 5698fa9e22eSVaclav Hapla for (r = 0; r < n; ++r) 5708fa9e22eSVaclav Hapla if (lidxs[r] >= 0) { 5718fa9e22eSVaclav Hapla if (work) work[len] = work[r]; 5728fa9e22eSVaclav Hapla lidxs[len++] = r; 5738fa9e22eSVaclav Hapla } 5748fa9e22eSVaclav Hapla if (on) *on = len; 5758fa9e22eSVaclav Hapla if (oidxs) *oidxs = lidxs; 5768fa9e22eSVaclav Hapla if (ogidxs) *ogidxs = work; 5778fa9e22eSVaclav Hapla PetscFunctionReturn(0); 5788fa9e22eSVaclav Hapla } 579deffd5ebSksagiyam 580deffd5ebSksagiyam /*@ 581deffd5ebSksagiyam PetscSFCreateByMatchingIndices - Create SF by matching root and leaf indices 582deffd5ebSksagiyam 583deffd5ebSksagiyam Collective 584deffd5ebSksagiyam 585deffd5ebSksagiyam Input Parameters: 586deffd5ebSksagiyam + layout - PetscLayout defining the global index space and the rank that brokers each index 587deffd5ebSksagiyam . numRootIndices - size of rootIndices 588deffd5ebSksagiyam . rootIndices - PetscInt array of global indices of which this process requests ownership 589deffd5ebSksagiyam . rootLocalIndices - root local index permutation (NULL if no permutation) 590deffd5ebSksagiyam . rootLocalOffset - offset to be added to root local indices 591deffd5ebSksagiyam . numLeafIndices - size of leafIndices 592deffd5ebSksagiyam . leafIndices - PetscInt array of global indices with which this process requires data associated 593deffd5ebSksagiyam . leafLocalIndices - leaf local index permutation (NULL if no permutation) 594deffd5ebSksagiyam - leafLocalOffset - offset to be added to leaf local indices 595deffd5ebSksagiyam 596d8d19677SJose E. Roman Output Parameters: 597deffd5ebSksagiyam + sfA - star forest representing the communication pattern from the layout space to the leaf space (NULL if not needed) 598deffd5ebSksagiyam - sf - star forest representing the communication pattern from the root space to the leaf space 599deffd5ebSksagiyam 600deffd5ebSksagiyam Example 1: 601deffd5ebSksagiyam $ 602deffd5ebSksagiyam $ rank : 0 1 2 603deffd5ebSksagiyam $ rootIndices : [1 0 2] [3] [3] 604deffd5ebSksagiyam $ rootLocalOffset : 100 200 300 605deffd5ebSksagiyam $ layout : [0 1] [2] [3] 606deffd5ebSksagiyam $ leafIndices : [0] [2] [0 3] 607deffd5ebSksagiyam $ leafLocalOffset : 400 500 600 608deffd5ebSksagiyam $ 6097ef5781fSVaclav Hapla would build the following SF 610deffd5ebSksagiyam $ 611deffd5ebSksagiyam $ [0] 400 <- (0,101) 612deffd5ebSksagiyam $ [1] 500 <- (0,102) 613deffd5ebSksagiyam $ [2] 600 <- (0,101) 614deffd5ebSksagiyam $ [2] 601 <- (2,300) 615deffd5ebSksagiyam $ 616deffd5ebSksagiyam Example 2: 617deffd5ebSksagiyam $ 618deffd5ebSksagiyam $ rank : 0 1 2 619deffd5ebSksagiyam $ rootIndices : [1 0 2] [3] [3] 620deffd5ebSksagiyam $ rootLocalOffset : 100 200 300 621deffd5ebSksagiyam $ layout : [0 1] [2] [3] 622deffd5ebSksagiyam $ leafIndices : rootIndices rootIndices rootIndices 623deffd5ebSksagiyam $ leafLocalOffset : rootLocalOffset rootLocalOffset rootLocalOffset 624deffd5ebSksagiyam $ 6257ef5781fSVaclav Hapla would build the following SF 626deffd5ebSksagiyam $ 627deffd5ebSksagiyam $ [1] 200 <- (2,300) 628deffd5ebSksagiyam $ 629deffd5ebSksagiyam Example 3: 630deffd5ebSksagiyam $ 631deffd5ebSksagiyam $ No process requests ownership of global index 1, but no process needs it. 632deffd5ebSksagiyam $ 633deffd5ebSksagiyam $ rank : 0 1 2 634deffd5ebSksagiyam $ numRootIndices : 2 1 1 635deffd5ebSksagiyam $ rootIndices : [0 2] [3] [3] 636deffd5ebSksagiyam $ rootLocalOffset : 100 200 300 637deffd5ebSksagiyam $ layout : [0 1] [2] [3] 638deffd5ebSksagiyam $ numLeafIndices : 1 1 2 639deffd5ebSksagiyam $ leafIndices : [0] [2] [0 3] 640deffd5ebSksagiyam $ leafLocalOffset : 400 500 600 641deffd5ebSksagiyam $ 6427ef5781fSVaclav Hapla would build the following SF 643deffd5ebSksagiyam $ 644deffd5ebSksagiyam $ [0] 400 <- (0,100) 645deffd5ebSksagiyam $ [1] 500 <- (0,101) 646deffd5ebSksagiyam $ [2] 600 <- (0,100) 647deffd5ebSksagiyam $ [2] 601 <- (2,300) 648deffd5ebSksagiyam $ 649deffd5ebSksagiyam 650deffd5ebSksagiyam Notes: 651deffd5ebSksagiyam The layout parameter represents any partitioning of [0, N), where N is the total number of global indices, and its 652deffd5ebSksagiyam local size can be set to PETSC_DECIDE. 653deffd5ebSksagiyam If a global index x lies in the partition owned by process i, each process whose rootIndices contains x requests 654deffd5ebSksagiyam ownership of x and sends its own rank and the local index of x to process i. 655deffd5ebSksagiyam If multiple processes request ownership of x, the one with the highest rank is to own x. 656deffd5ebSksagiyam Process i then broadcasts the ownership information, so that each process whose leafIndices contains x knows the 657deffd5ebSksagiyam ownership information of x. 658deffd5ebSksagiyam The output sf is constructed by associating each leaf point to a root point in this way. 659deffd5ebSksagiyam 660deffd5ebSksagiyam Suppose there is point data ordered according to the global indices and partitioned according to the given layout. 661deffd5ebSksagiyam The optional output PetscSF object sfA can be used to push such data to leaf points. 662deffd5ebSksagiyam 663deffd5ebSksagiyam All indices in rootIndices and leafIndices must lie in the layout range. The union (over all processes) of rootIndices 664deffd5ebSksagiyam must cover that of leafIndices, but need not cover the entire layout. 665deffd5ebSksagiyam 666deffd5ebSksagiyam If (leafIndices, leafLocalIndices, leafLocalOffset) == (rootIndices, rootLocalIndices, rootLocalOffset), the output 667deffd5ebSksagiyam star forest is almost identity, so will only include non-trivial part of the map. 668deffd5ebSksagiyam 669deffd5ebSksagiyam Developer Notes: 670deffd5ebSksagiyam Current approach of a process of the highest rank gaining the ownership may cause load imbalance; consider using 671deffd5ebSksagiyam hash(rank, root_local_index) as the bid for the ownership determination. 672deffd5ebSksagiyam 673deffd5ebSksagiyam Level: advanced 674deffd5ebSksagiyam 675db781477SPatrick Sanan .seealso: `PetscSFCreate()` 676deffd5ebSksagiyam @*/ 677*d71ae5a4SJacob 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) 678*d71ae5a4SJacob Faibussowitsch { 679deffd5ebSksagiyam MPI_Comm comm = layout->comm; 680deffd5ebSksagiyam PetscMPIInt size, rank; 681deffd5ebSksagiyam PetscSF sf1; 682deffd5ebSksagiyam PetscSFNode *owners, *buffer, *iremote; 683deffd5ebSksagiyam PetscInt *ilocal, nleaves, N, n, i; 684deffd5ebSksagiyam #if defined(PETSC_USE_DEBUG) 685deffd5ebSksagiyam PetscInt N1; 686deffd5ebSksagiyam #endif 687deffd5ebSksagiyam PetscBool flag; 688deffd5ebSksagiyam 689deffd5ebSksagiyam PetscFunctionBegin; 690deffd5ebSksagiyam if (rootIndices) PetscValidIntPointer(rootIndices, 3); 691deffd5ebSksagiyam if (rootLocalIndices) PetscValidIntPointer(rootLocalIndices, 4); 692deffd5ebSksagiyam if (leafIndices) PetscValidIntPointer(leafIndices, 7); 693deffd5ebSksagiyam if (leafLocalIndices) PetscValidIntPointer(leafLocalIndices, 8); 694deffd5ebSksagiyam if (sfA) PetscValidPointer(sfA, 10); 695deffd5ebSksagiyam PetscValidPointer(sf, 11); 69608401ef6SPierre Jolivet PetscCheck(numRootIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numRootIndices (%" PetscInt_FMT ") must be non-negative", numRootIndices); 69708401ef6SPierre Jolivet PetscCheck(numLeafIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numLeafIndices (%" PetscInt_FMT ") must be non-negative", numLeafIndices); 6989566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 6999566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 7009566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(layout)); 7019566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(layout, &N)); 7029566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(layout, &n)); 703deffd5ebSksagiyam flag = (PetscBool)(leafIndices == rootIndices); 7041c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &flag, 1, MPIU_BOOL, MPI_LAND, comm)); 705c9cc58a2SBarry Smith PetscCheck(!flag || numLeafIndices == numRootIndices, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "leafIndices == rootIndices, but numLeafIndices (%" PetscInt_FMT ") != numRootIndices(%" PetscInt_FMT ")", numLeafIndices, numRootIndices); 706deffd5ebSksagiyam #if defined(PETSC_USE_DEBUG) 707deffd5ebSksagiyam N1 = PETSC_MIN_INT; 7089371c9d4SSatish Balay for (i = 0; i < numRootIndices; i++) 7099371c9d4SSatish Balay if (rootIndices[i] > N1) N1 = rootIndices[i]; 7109566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm)); 71108401ef6SPierre 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); 712deffd5ebSksagiyam if (!flag) { 713deffd5ebSksagiyam N1 = PETSC_MIN_INT; 7149371c9d4SSatish Balay for (i = 0; i < numLeafIndices; i++) 7159371c9d4SSatish Balay if (leafIndices[i] > N1) N1 = leafIndices[i]; 7169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm)); 71708401ef6SPierre 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); 718deffd5ebSksagiyam } 719deffd5ebSksagiyam #endif 720deffd5ebSksagiyam /* Reduce: owners -> buffer */ 7219566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &buffer)); 7229566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, &sf1)); 7239566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sf1)); 7249566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraphLayout(sf1, layout, numRootIndices, NULL, PETSC_OWN_POINTER, rootIndices)); 7259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numRootIndices, &owners)); 726deffd5ebSksagiyam for (i = 0; i < numRootIndices; ++i) { 727deffd5ebSksagiyam owners[i].rank = rank; 728deffd5ebSksagiyam owners[i].index = rootLocalOffset + (rootLocalIndices ? rootLocalIndices[i] : i); 729deffd5ebSksagiyam } 730deffd5ebSksagiyam for (i = 0; i < n; ++i) { 731deffd5ebSksagiyam buffer[i].index = -1; 732deffd5ebSksagiyam buffer[i].rank = -1; 733deffd5ebSksagiyam } 7349566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC)); 7359566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC)); 736deffd5ebSksagiyam /* Bcast: buffer -> owners */ 737deffd5ebSksagiyam if (!flag) { 738deffd5ebSksagiyam /* leafIndices is different from rootIndices */ 7399566063dSJacob Faibussowitsch PetscCall(PetscFree(owners)); 7409566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraphLayout(sf1, layout, numLeafIndices, NULL, PETSC_OWN_POINTER, leafIndices)); 7419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numLeafIndices, &owners)); 742deffd5ebSksagiyam } 7439566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE)); 7449566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE)); 74508401ef6SPierre 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]); 7469566063dSJacob Faibussowitsch PetscCall(PetscFree(buffer)); 7479371c9d4SSatish Balay if (sfA) { 7489371c9d4SSatish Balay *sfA = sf1; 7499371c9d4SSatish Balay } else PetscCall(PetscSFDestroy(&sf1)); 750deffd5ebSksagiyam /* Create sf */ 751deffd5ebSksagiyam if (flag && rootLocalIndices == leafLocalIndices && leafLocalOffset == rootLocalOffset) { 752deffd5ebSksagiyam /* leaf space == root space */ 7539371c9d4SSatish Balay for (i = 0, nleaves = 0; i < numLeafIndices; ++i) 7549371c9d4SSatish Balay if (owners[i].rank != rank) ++nleaves; 7559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &ilocal)); 7569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &iremote)); 757deffd5ebSksagiyam for (i = 0, nleaves = 0; i < numLeafIndices; ++i) { 758deffd5ebSksagiyam if (owners[i].rank != rank) { 759deffd5ebSksagiyam ilocal[nleaves] = leafLocalOffset + i; 760deffd5ebSksagiyam iremote[nleaves].rank = owners[i].rank; 761deffd5ebSksagiyam iremote[nleaves].index = owners[i].index; 762deffd5ebSksagiyam ++nleaves; 763deffd5ebSksagiyam } 764deffd5ebSksagiyam } 7659566063dSJacob Faibussowitsch PetscCall(PetscFree(owners)); 766deffd5ebSksagiyam } else { 767deffd5ebSksagiyam nleaves = numLeafIndices; 7689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &ilocal)); 769ad540459SPierre Jolivet for (i = 0; i < nleaves; ++i) ilocal[i] = leafLocalOffset + (leafLocalIndices ? leafLocalIndices[i] : i); 770deffd5ebSksagiyam iremote = owners; 771deffd5ebSksagiyam } 7729566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, sf)); 7739566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(*sf)); 7749566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*sf, rootLocalOffset + numRootIndices, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER)); 775deffd5ebSksagiyam PetscFunctionReturn(0); 776deffd5ebSksagiyam } 777