1b0c7db22SLisandro Dalcin #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/ 2b0c7db22SLisandro Dalcin #include <petsc/private/sectionimpl.h> 3b0c7db22SLisandro Dalcin 4b0c7db22SLisandro Dalcin /*@C 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 1586c56f52SVaclav Hapla - iremote - 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 28db781477SPatrick Sanan .seealso: `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()` 29b0c7db22SLisandro Dalcin @*/ 309371c9d4SSatish Balay PetscErrorCode PetscSFSetGraphLayout(PetscSF sf, PetscLayout layout, PetscInt nleaves, PetscInt *ilocal, PetscCopyMode localmode, const PetscInt *iremote) { 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)); 4038a25198SStefano Zampini if (nleaves) { ls = iremote[0] + 1; } 41b0c7db22SLisandro Dalcin for (i = 0; i < nleaves; i++) { 4238a25198SStefano Zampini const PetscInt idx = iremote[i] - ls; 4338a25198SStefano Zampini if (idx < 0 || idx >= ln) { /* short-circuit the search */ 449566063dSJacob Faibussowitsch PetscCall(PetscLayoutFindOwnerIndex(layout, iremote[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 57b0c7db22SLisandro Dalcin /*@ 58b0c7db22SLisandro Dalcin PetscSFSetGraphSection - Sets the PetscSF graph encoding the parallel dof overlap based upon the PetscSections 59b0c7db22SLisandro Dalcin describing the data layout. 60b0c7db22SLisandro Dalcin 61b0c7db22SLisandro Dalcin Input Parameters: 62b0c7db22SLisandro Dalcin + sf - The SF 63b0c7db22SLisandro Dalcin . localSection - PetscSection describing the local data layout 64b0c7db22SLisandro Dalcin - globalSection - PetscSection describing the global data layout 65b0c7db22SLisandro Dalcin 66b0c7db22SLisandro Dalcin Level: developer 67b0c7db22SLisandro Dalcin 68db781477SPatrick Sanan .seealso: `PetscSFSetGraph()`, `PetscSFSetGraphLayout()` 69b0c7db22SLisandro Dalcin @*/ 709371c9d4SSatish Balay PetscErrorCode PetscSFSetGraphSection(PetscSF sf, PetscSection localSection, PetscSection globalSection) { 71b0c7db22SLisandro Dalcin MPI_Comm comm; 72b0c7db22SLisandro Dalcin PetscLayout layout; 73b0c7db22SLisandro Dalcin const PetscInt *ranges; 74b0c7db22SLisandro Dalcin PetscInt *local; 75b0c7db22SLisandro Dalcin PetscSFNode *remote; 76b0c7db22SLisandro Dalcin PetscInt pStart, pEnd, p, nroots, nleaves = 0, l; 77b0c7db22SLisandro Dalcin PetscMPIInt size, rank; 78b0c7db22SLisandro Dalcin 79b0c7db22SLisandro Dalcin PetscFunctionBegin; 80b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 81b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(localSection, PETSC_SECTION_CLASSID, 2); 82b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(globalSection, PETSC_SECTION_CLASSID, 3); 83b0c7db22SLisandro Dalcin 849566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 859566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 869566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 879566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(globalSection, &pStart, &pEnd)); 889566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstrainedStorageSize(globalSection, &nroots)); 899566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(comm, &layout)); 909566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(layout, 1)); 919566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetLocalSize(layout, nroots)); 929566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(layout)); 939566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRanges(layout, &ranges)); 94b0c7db22SLisandro Dalcin for (p = pStart; p < pEnd; ++p) { 95b0c7db22SLisandro Dalcin PetscInt gdof, gcdof; 96b0c7db22SLisandro Dalcin 979566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(globalSection, p, &gdof)); 989566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof)); 99c9cc58a2SBarry 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)); 100b0c7db22SLisandro Dalcin nleaves += gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof; 101b0c7db22SLisandro Dalcin } 1029566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &local)); 1039566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &remote)); 104b0c7db22SLisandro Dalcin for (p = pStart, l = 0; p < pEnd; ++p) { 105b0c7db22SLisandro Dalcin const PetscInt *cind; 106b0c7db22SLisandro Dalcin PetscInt dof, cdof, off, gdof, gcdof, goff, gsize, d, c; 107b0c7db22SLisandro Dalcin 1089566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(localSection, p, &dof)); 1099566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(localSection, p, &off)); 1109566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(localSection, p, &cdof)); 1119566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintIndices(localSection, p, &cind)); 1129566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(globalSection, p, &gdof)); 1139566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof)); 1149566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(globalSection, p, &goff)); 115b0c7db22SLisandro Dalcin if (!gdof) continue; /* Censored point */ 116b0c7db22SLisandro Dalcin gsize = gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof; 117b0c7db22SLisandro Dalcin if (gsize != dof - cdof) { 11808401ef6SPierre 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); 119b0c7db22SLisandro Dalcin cdof = 0; /* Ignore constraints */ 120b0c7db22SLisandro Dalcin } 121b0c7db22SLisandro Dalcin for (d = 0, c = 0; d < dof; ++d) { 1229371c9d4SSatish Balay if ((c < cdof) && (cind[c] == d)) { 1239371c9d4SSatish Balay ++c; 1249371c9d4SSatish Balay continue; 1259371c9d4SSatish Balay } 126b0c7db22SLisandro Dalcin local[l + d - c] = off + d; 127b0c7db22SLisandro Dalcin } 12808401ef6SPierre 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); 129b0c7db22SLisandro Dalcin if (gdof < 0) { 130b0c7db22SLisandro Dalcin for (d = 0; d < gsize; ++d, ++l) { 131b0c7db22SLisandro Dalcin PetscInt offset = -(goff + 1) + d, r; 132b0c7db22SLisandro Dalcin 1339566063dSJacob Faibussowitsch PetscCall(PetscFindInt(offset, size + 1, ranges, &r)); 134b0c7db22SLisandro Dalcin if (r < 0) r = -(r + 2); 13508401ef6SPierre 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); 136b0c7db22SLisandro Dalcin remote[l].rank = r; 137b0c7db22SLisandro Dalcin remote[l].index = offset - ranges[r]; 138b0c7db22SLisandro Dalcin } 139b0c7db22SLisandro Dalcin } else { 140b0c7db22SLisandro Dalcin for (d = 0; d < gsize; ++d, ++l) { 141b0c7db22SLisandro Dalcin remote[l].rank = rank; 142b0c7db22SLisandro Dalcin remote[l].index = goff + d - ranges[rank]; 143b0c7db22SLisandro Dalcin } 144b0c7db22SLisandro Dalcin } 145b0c7db22SLisandro Dalcin } 14608401ef6SPierre Jolivet PetscCheck(l == nleaves, comm, PETSC_ERR_PLIB, "Iteration error, l %" PetscInt_FMT " != nleaves %" PetscInt_FMT, l, nleaves); 1479566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&layout)); 1489566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER)); 149b0c7db22SLisandro Dalcin PetscFunctionReturn(0); 150b0c7db22SLisandro Dalcin } 151b0c7db22SLisandro Dalcin 152b0c7db22SLisandro Dalcin /*@C 153b0c7db22SLisandro Dalcin PetscSFDistributeSection - Create a new PetscSection reorganized, moving from the root to the leaves of the SF 154b0c7db22SLisandro Dalcin 155b0c7db22SLisandro Dalcin Collective on sf 156b0c7db22SLisandro Dalcin 157b0c7db22SLisandro Dalcin Input Parameters: 158b0c7db22SLisandro Dalcin + sf - The SF 159b0c7db22SLisandro Dalcin - rootSection - Section defined on root space 160b0c7db22SLisandro Dalcin 161b0c7db22SLisandro Dalcin Output Parameters: 162b0c7db22SLisandro Dalcin + remoteOffsets - root offsets in leaf storage, or NULL 163b0c7db22SLisandro Dalcin - leafSection - Section defined on the leaf space 164b0c7db22SLisandro Dalcin 165b0c7db22SLisandro Dalcin Level: advanced 166b0c7db22SLisandro Dalcin 167db781477SPatrick Sanan .seealso: `PetscSFCreate()` 168b0c7db22SLisandro Dalcin @*/ 1699371c9d4SSatish Balay PetscErrorCode PetscSFDistributeSection(PetscSF sf, PetscSection rootSection, PetscInt **remoteOffsets, PetscSection leafSection) { 170b0c7db22SLisandro Dalcin PetscSF embedSF; 171b0c7db22SLisandro Dalcin const PetscInt *indices; 172b0c7db22SLisandro Dalcin IS selected; 173b0c7db22SLisandro Dalcin PetscInt numFields, nroots, rpStart, rpEnd, lpStart = PETSC_MAX_INT, lpEnd = -1, f, c; 174b0c7db22SLisandro Dalcin PetscBool *sub, hasc; 175b0c7db22SLisandro Dalcin 176b0c7db22SLisandro Dalcin PetscFunctionBegin; 1779566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_DistSect, sf, 0, 0, 0)); 1789566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(rootSection, &numFields)); 179029e8818Sksagiyam if (numFields) { 180029e8818Sksagiyam IS perm; 181029e8818Sksagiyam 182029e8818Sksagiyam /* PetscSectionSetNumFields() calls PetscSectionReset(), which destroys 183029e8818Sksagiyam leafSection->perm. To keep this permutation set by the user, we grab 184029e8818Sksagiyam the reference before calling PetscSectionSetNumFields() and set it 185029e8818Sksagiyam back after. */ 1869566063dSJacob Faibussowitsch PetscCall(PetscSectionGetPermutation(leafSection, &perm)); 1879566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)perm)); 1889566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(leafSection, numFields)); 1899566063dSJacob Faibussowitsch PetscCall(PetscSectionSetPermutation(leafSection, perm)); 1909566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 191029e8818Sksagiyam } 1929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numFields + 2, &sub)); 193b0c7db22SLisandro Dalcin sub[1] = rootSection->bc ? PETSC_TRUE : PETSC_FALSE; 194b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) { 1952ee82556SMatthew G. Knepley PetscSectionSym sym, dsym = NULL; 196b0c7db22SLisandro Dalcin const char *name = NULL; 197b0c7db22SLisandro Dalcin PetscInt numComp = 0; 198b0c7db22SLisandro Dalcin 199b0c7db22SLisandro Dalcin sub[2 + f] = rootSection->field[f]->bc ? PETSC_TRUE : PETSC_FALSE; 2009566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldComponents(rootSection, f, &numComp)); 2019566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldName(rootSection, f, &name)); 2029566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldSym(rootSection, f, &sym)); 2039566063dSJacob Faibussowitsch if (sym) PetscCall(PetscSectionSymDistribute(sym, sf, &dsym)); 2049566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(leafSection, f, numComp)); 2059566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldName(leafSection, f, name)); 2069566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldSym(leafSection, f, dsym)); 2079566063dSJacob Faibussowitsch PetscCall(PetscSectionSymDestroy(&dsym)); 208b0c7db22SLisandro Dalcin for (c = 0; c < rootSection->numFieldComponents[f]; ++c) { 2099566063dSJacob Faibussowitsch PetscCall(PetscSectionGetComponentName(rootSection, f, c, &name)); 2109566063dSJacob Faibussowitsch PetscCall(PetscSectionSetComponentName(leafSection, f, c, name)); 211b0c7db22SLisandro Dalcin } 212b0c7db22SLisandro Dalcin } 2139566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd)); 2149566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL)); 215b0c7db22SLisandro Dalcin rpEnd = PetscMin(rpEnd, nroots); 216b0c7db22SLisandro Dalcin rpEnd = PetscMax(rpStart, rpEnd); 217b0c7db22SLisandro Dalcin /* see if we can avoid creating the embedded SF, since it can cost more than an allreduce */ 218b0c7db22SLisandro Dalcin sub[0] = (PetscBool)(nroots != rpEnd - rpStart); 2191c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(MPI_IN_PLACE, sub, 2 + numFields, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)sf))); 220b0c7db22SLisandro Dalcin if (sub[0]) { 2219566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected)); 2229566063dSJacob Faibussowitsch PetscCall(ISGetIndices(selected, &indices)); 2239566063dSJacob Faibussowitsch PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF)); 2249566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(selected, &indices)); 2259566063dSJacob Faibussowitsch PetscCall(ISDestroy(&selected)); 226b0c7db22SLisandro Dalcin } else { 2279566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)sf)); 228b0c7db22SLisandro Dalcin embedSF = sf; 229b0c7db22SLisandro Dalcin } 2309566063dSJacob Faibussowitsch PetscCall(PetscSFGetLeafRange(embedSF, &lpStart, &lpEnd)); 231b0c7db22SLisandro Dalcin lpEnd++; 232b0c7db22SLisandro Dalcin 2339566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(leafSection, lpStart, lpEnd)); 234b0c7db22SLisandro Dalcin 235b0c7db22SLisandro Dalcin /* Constrained dof section */ 236b0c7db22SLisandro Dalcin hasc = sub[1]; 237b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) hasc = (PetscBool)(hasc || sub[2 + f]); 238b0c7db22SLisandro Dalcin 239b0c7db22SLisandro Dalcin /* Could fuse these at the cost of copies and extra allocation */ 2409566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart], MPI_REPLACE)); 2419566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart], MPI_REPLACE)); 242b0c7db22SLisandro Dalcin if (sub[1]) { 2437c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(rootSection)); 2447c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(leafSection)); 2459566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE)); 2469566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE)); 247b0c7db22SLisandro Dalcin } 248b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) { 2499566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart], MPI_REPLACE)); 2509566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart], MPI_REPLACE)); 251b0c7db22SLisandro Dalcin if (sub[2 + f]) { 2527c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(rootSection->field[f])); 2537c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(leafSection->field[f])); 2549566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE)); 2559566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE)); 256b0c7db22SLisandro Dalcin } 257b0c7db22SLisandro Dalcin } 258b0c7db22SLisandro Dalcin if (remoteOffsets) { 2599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lpEnd - lpStart, remoteOffsets)); 2609566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE)); 2619566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE)); 262b0c7db22SLisandro Dalcin } 26369c11d05SVaclav Hapla PetscCall(PetscSectionInvalidateMaxDof_Internal(leafSection)); 2649566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(leafSection)); 265b0c7db22SLisandro Dalcin if (hasc) { /* need to communicate bcIndices */ 266b0c7db22SLisandro Dalcin PetscSF bcSF; 267b0c7db22SLisandro Dalcin PetscInt *rOffBc; 268b0c7db22SLisandro Dalcin 2699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lpEnd - lpStart, &rOffBc)); 270b0c7db22SLisandro Dalcin if (sub[1]) { 2719566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 2729566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 2739566063dSJacob Faibussowitsch PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->bc, rOffBc, leafSection->bc, &bcSF)); 2749566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE)); 2759566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE)); 2769566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&bcSF)); 277b0c7db22SLisandro Dalcin } 278b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) { 279b0c7db22SLisandro Dalcin if (sub[2 + f]) { 2809566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 2819566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 2829566063dSJacob Faibussowitsch PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->field[f]->bc, rOffBc, leafSection->field[f]->bc, &bcSF)); 2839566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE)); 2849566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE)); 2859566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&bcSF)); 286b0c7db22SLisandro Dalcin } 287b0c7db22SLisandro Dalcin } 2889566063dSJacob Faibussowitsch PetscCall(PetscFree(rOffBc)); 289b0c7db22SLisandro Dalcin } 2909566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&embedSF)); 2919566063dSJacob Faibussowitsch PetscCall(PetscFree(sub)); 2929566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_DistSect, sf, 0, 0, 0)); 293b0c7db22SLisandro Dalcin PetscFunctionReturn(0); 294b0c7db22SLisandro Dalcin } 295b0c7db22SLisandro Dalcin 296b0c7db22SLisandro Dalcin /*@C 297b0c7db22SLisandro Dalcin PetscSFCreateRemoteOffsets - Create offsets for point data on remote processes 298b0c7db22SLisandro Dalcin 299b0c7db22SLisandro Dalcin Collective on sf 300b0c7db22SLisandro Dalcin 301b0c7db22SLisandro Dalcin Input Parameters: 302b0c7db22SLisandro Dalcin + sf - The SF 303b0c7db22SLisandro Dalcin . rootSection - Data layout of remote points for outgoing data (this is layout for SF roots) 304b0c7db22SLisandro Dalcin - leafSection - Data layout of local points for incoming data (this is layout for SF leaves) 305b0c7db22SLisandro Dalcin 306b0c7db22SLisandro Dalcin Output Parameter: 307b0c7db22SLisandro Dalcin . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL 308b0c7db22SLisandro Dalcin 309b0c7db22SLisandro Dalcin Level: developer 310b0c7db22SLisandro Dalcin 311db781477SPatrick Sanan .seealso: `PetscSFCreate()` 312b0c7db22SLisandro Dalcin @*/ 3139371c9d4SSatish Balay PetscErrorCode PetscSFCreateRemoteOffsets(PetscSF sf, PetscSection rootSection, PetscSection leafSection, PetscInt **remoteOffsets) { 314b0c7db22SLisandro Dalcin PetscSF embedSF; 315b0c7db22SLisandro Dalcin const PetscInt *indices; 316b0c7db22SLisandro Dalcin IS selected; 317b0c7db22SLisandro Dalcin PetscInt numRoots, rpStart = 0, rpEnd = 0, lpStart = 0, lpEnd = 0; 318b0c7db22SLisandro Dalcin 319b0c7db22SLisandro Dalcin PetscFunctionBegin; 320b0c7db22SLisandro Dalcin *remoteOffsets = NULL; 3219566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &numRoots, NULL, NULL, NULL)); 322b0c7db22SLisandro Dalcin if (numRoots < 0) PetscFunctionReturn(0); 3239566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_RemoteOff, sf, 0, 0, 0)); 3249566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd)); 3259566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd)); 3269566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected)); 3279566063dSJacob Faibussowitsch PetscCall(ISGetIndices(selected, &indices)); 3289566063dSJacob Faibussowitsch PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF)); 3299566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(selected, &indices)); 3309566063dSJacob Faibussowitsch PetscCall(ISDestroy(&selected)); 3319566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(lpEnd - lpStart, remoteOffsets)); 3329566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE)); 3339566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE)); 3349566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&embedSF)); 3359566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_RemoteOff, sf, 0, 0, 0)); 336b0c7db22SLisandro Dalcin PetscFunctionReturn(0); 337b0c7db22SLisandro Dalcin } 338b0c7db22SLisandro Dalcin 339b0c7db22SLisandro Dalcin /*@C 340b0c7db22SLisandro Dalcin PetscSFCreateSectionSF - Create an expanded SF of dofs, assuming the input SF relates points 341b0c7db22SLisandro Dalcin 342b0c7db22SLisandro Dalcin Collective on sf 343b0c7db22SLisandro Dalcin 344b0c7db22SLisandro Dalcin Input Parameters: 345b0c7db22SLisandro Dalcin + sf - The SF 346b0c7db22SLisandro Dalcin . rootSection - Data layout of remote points for outgoing data (this is usually the serial section) 347b0c7db22SLisandro Dalcin . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL 348b0c7db22SLisandro Dalcin - leafSection - Data layout of local points for incoming data (this is the distributed section) 349b0c7db22SLisandro Dalcin 350b0c7db22SLisandro Dalcin Output Parameters: 351b0c7db22SLisandro Dalcin - sectionSF - The new SF 352b0c7db22SLisandro Dalcin 353b0c7db22SLisandro Dalcin Note: Either rootSection or remoteOffsets can be specified 354b0c7db22SLisandro Dalcin 355b0c7db22SLisandro Dalcin Level: advanced 356b0c7db22SLisandro Dalcin 357db781477SPatrick Sanan .seealso: `PetscSFCreate()` 358b0c7db22SLisandro Dalcin @*/ 3599371c9d4SSatish Balay PetscErrorCode PetscSFCreateSectionSF(PetscSF sf, PetscSection rootSection, PetscInt remoteOffsets[], PetscSection leafSection, PetscSF *sectionSF) { 360b0c7db22SLisandro Dalcin MPI_Comm comm; 361b0c7db22SLisandro Dalcin const PetscInt *localPoints; 362b0c7db22SLisandro Dalcin const PetscSFNode *remotePoints; 363b0c7db22SLisandro Dalcin PetscInt lpStart, lpEnd; 364b0c7db22SLisandro Dalcin PetscInt numRoots, numSectionRoots, numPoints, numIndices = 0; 365b0c7db22SLisandro Dalcin PetscInt *localIndices; 366b0c7db22SLisandro Dalcin PetscSFNode *remoteIndices; 367b0c7db22SLisandro Dalcin PetscInt i, ind; 368b0c7db22SLisandro Dalcin 369b0c7db22SLisandro Dalcin PetscFunctionBegin; 370b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 371b0c7db22SLisandro Dalcin PetscValidPointer(rootSection, 2); 372b0c7db22SLisandro Dalcin /* Cannot check PetscValidIntPointer(remoteOffsets,3) because it can be NULL if sf does not reference any points in leafSection */ 373b0c7db22SLisandro Dalcin PetscValidPointer(leafSection, 4); 374b0c7db22SLisandro Dalcin PetscValidPointer(sectionSF, 5); 3759566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 3769566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, sectionSF)); 3779566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd)); 3789566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(rootSection, &numSectionRoots)); 3799566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &numRoots, &numPoints, &localPoints, &remotePoints)); 380b0c7db22SLisandro Dalcin if (numRoots < 0) PetscFunctionReturn(0); 3819566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_SectSF, sf, 0, 0, 0)); 382b0c7db22SLisandro Dalcin for (i = 0; i < numPoints; ++i) { 383b0c7db22SLisandro Dalcin PetscInt localPoint = localPoints ? localPoints[i] : i; 384b0c7db22SLisandro Dalcin PetscInt dof; 385b0c7db22SLisandro Dalcin 386b0c7db22SLisandro Dalcin if ((localPoint >= lpStart) && (localPoint < lpEnd)) { 3879566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof)); 38822a18c9fSMatthew G. Knepley numIndices += dof < 0 ? 0 : dof; 389b0c7db22SLisandro Dalcin } 390b0c7db22SLisandro Dalcin } 3919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numIndices, &localIndices)); 3929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numIndices, &remoteIndices)); 393b0c7db22SLisandro Dalcin /* Create new index graph */ 394b0c7db22SLisandro Dalcin for (i = 0, ind = 0; i < numPoints; ++i) { 395b0c7db22SLisandro Dalcin PetscInt localPoint = localPoints ? localPoints[i] : i; 396b0c7db22SLisandro Dalcin PetscInt rank = remotePoints[i].rank; 397b0c7db22SLisandro Dalcin 398b0c7db22SLisandro Dalcin if ((localPoint >= lpStart) && (localPoint < lpEnd)) { 399b0c7db22SLisandro Dalcin PetscInt remoteOffset = remoteOffsets[localPoint - lpStart]; 400b0c7db22SLisandro Dalcin PetscInt loff, dof, d; 401b0c7db22SLisandro Dalcin 4029566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(leafSection, localPoint, &loff)); 4039566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof)); 404b0c7db22SLisandro Dalcin for (d = 0; d < dof; ++d, ++ind) { 405b0c7db22SLisandro Dalcin localIndices[ind] = loff + d; 406b0c7db22SLisandro Dalcin remoteIndices[ind].rank = rank; 407b0c7db22SLisandro Dalcin remoteIndices[ind].index = remoteOffset + d; 408b0c7db22SLisandro Dalcin } 409b0c7db22SLisandro Dalcin } 410b0c7db22SLisandro Dalcin } 41108401ef6SPierre Jolivet PetscCheck(numIndices == ind, comm, PETSC_ERR_PLIB, "Inconsistency in indices, %" PetscInt_FMT " should be %" PetscInt_FMT, ind, numIndices); 4129566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*sectionSF, numSectionRoots, numIndices, localIndices, PETSC_OWN_POINTER, remoteIndices, PETSC_OWN_POINTER)); 4139566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(*sectionSF)); 4149566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_SectSF, sf, 0, 0, 0)); 415b0c7db22SLisandro Dalcin PetscFunctionReturn(0); 416b0c7db22SLisandro Dalcin } 4178fa9e22eSVaclav Hapla 4188fa9e22eSVaclav Hapla /*@C 4198fa9e22eSVaclav Hapla PetscSFCreateFromLayouts - Creates a parallel star forest mapping two PetscLayout objects 4208fa9e22eSVaclav Hapla 4218fa9e22eSVaclav Hapla Collective 4228fa9e22eSVaclav Hapla 4234165533cSJose E. Roman Input Parameters: 4248fa9e22eSVaclav Hapla + rmap - PetscLayout defining the global root space 4258fa9e22eSVaclav Hapla - lmap - PetscLayout defining the global leaf space 4268fa9e22eSVaclav Hapla 4274165533cSJose E. Roman Output Parameter: 4288fa9e22eSVaclav Hapla . sf - The parallel star forest 4298fa9e22eSVaclav Hapla 4308fa9e22eSVaclav Hapla Level: intermediate 4318fa9e22eSVaclav Hapla 432db781477SPatrick Sanan .seealso: `PetscSFCreate()`, `PetscLayoutCreate()`, `PetscSFSetGraphLayout()` 4338fa9e22eSVaclav Hapla @*/ 4349371c9d4SSatish Balay PetscErrorCode PetscSFCreateFromLayouts(PetscLayout rmap, PetscLayout lmap, PetscSF *sf) { 4358fa9e22eSVaclav Hapla PetscInt i, nroots, nleaves = 0; 4368fa9e22eSVaclav Hapla PetscInt rN, lst, len; 4378fa9e22eSVaclav Hapla PetscMPIInt owner = -1; 4388fa9e22eSVaclav Hapla PetscSFNode *remote; 4398fa9e22eSVaclav Hapla MPI_Comm rcomm = rmap->comm; 4408fa9e22eSVaclav Hapla MPI_Comm lcomm = lmap->comm; 4418fa9e22eSVaclav Hapla PetscMPIInt flg; 4428fa9e22eSVaclav Hapla 4438fa9e22eSVaclav Hapla PetscFunctionBegin; 4448fa9e22eSVaclav Hapla PetscValidPointer(sf, 3); 44528b400f6SJacob Faibussowitsch PetscCheck(rmap->setupcalled, rcomm, PETSC_ERR_ARG_WRONGSTATE, "Root layout not setup"); 44628b400f6SJacob Faibussowitsch PetscCheck(lmap->setupcalled, lcomm, PETSC_ERR_ARG_WRONGSTATE, "Leaf layout not setup"); 4479566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(rcomm, lcomm, &flg)); 448c9cc58a2SBarry Smith PetscCheck(flg == MPI_CONGRUENT || flg == MPI_IDENT, rcomm, PETSC_ERR_SUP, "cannot map two layouts with non-matching communicators"); 4499566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(rcomm, sf)); 4509566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(rmap, &nroots)); 4519566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(rmap, &rN)); 4529566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRange(lmap, &lst, &len)); 4539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len - lst, &remote)); 4548fa9e22eSVaclav Hapla for (i = lst; i < len && i < rN; i++) { 455*48a46eb9SPierre Jolivet if (owner < -1 || i >= rmap->range[owner + 1]) PetscCall(PetscLayoutFindOwner(rmap, i, &owner)); 4568fa9e22eSVaclav Hapla remote[nleaves].rank = owner; 4578fa9e22eSVaclav Hapla remote[nleaves].index = i - rmap->range[owner]; 4588fa9e22eSVaclav Hapla nleaves++; 4598fa9e22eSVaclav Hapla } 4609566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*sf, nroots, nleaves, NULL, PETSC_OWN_POINTER, remote, PETSC_COPY_VALUES)); 4619566063dSJacob Faibussowitsch PetscCall(PetscFree(remote)); 4628fa9e22eSVaclav Hapla PetscFunctionReturn(0); 4638fa9e22eSVaclav Hapla } 4648fa9e22eSVaclav Hapla 4658fa9e22eSVaclav Hapla /* TODO: handle nooffprocentries like MatZeroRowsMapLocal_Private, since this code is the same */ 4669371c9d4SSatish Balay PetscErrorCode PetscLayoutMapLocal(PetscLayout map, PetscInt N, const PetscInt idxs[], PetscInt *on, PetscInt **oidxs, PetscInt **ogidxs) { 4678fa9e22eSVaclav Hapla PetscInt *owners = map->range; 4688fa9e22eSVaclav Hapla PetscInt n = map->n; 4698fa9e22eSVaclav Hapla PetscSF sf; 4708fa9e22eSVaclav Hapla PetscInt *lidxs, *work = NULL; 4718fa9e22eSVaclav Hapla PetscSFNode *ridxs; 4728fa9e22eSVaclav Hapla PetscMPIInt rank, p = 0; 4738fa9e22eSVaclav Hapla PetscInt r, len = 0; 4748fa9e22eSVaclav Hapla 4758fa9e22eSVaclav Hapla PetscFunctionBegin; 4768fa9e22eSVaclav Hapla if (on) *on = 0; /* squelch -Wmaybe-uninitialized */ 4778fa9e22eSVaclav Hapla /* Create SF where leaves are input idxs and roots are owned idxs */ 4789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(map->comm, &rank)); 4799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &lidxs)); 4808fa9e22eSVaclav Hapla for (r = 0; r < n; ++r) lidxs[r] = -1; 4819566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N, &ridxs)); 4828fa9e22eSVaclav Hapla for (r = 0; r < N; ++r) { 4838fa9e22eSVaclav Hapla const PetscInt idx = idxs[r]; 484c9cc58a2SBarry 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); 4858fa9e22eSVaclav Hapla if (idx < owners[p] || owners[p + 1] <= idx) { /* short-circuit the search if the last p owns this idx too */ 4869566063dSJacob Faibussowitsch PetscCall(PetscLayoutFindOwner(map, idx, &p)); 4878fa9e22eSVaclav Hapla } 4888fa9e22eSVaclav Hapla ridxs[r].rank = p; 4898fa9e22eSVaclav Hapla ridxs[r].index = idxs[r] - owners[p]; 4908fa9e22eSVaclav Hapla } 4919566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(map->comm, &sf)); 4929566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, n, N, NULL, PETSC_OWN_POINTER, ridxs, PETSC_OWN_POINTER)); 4939566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR)); 4949566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR)); 4958fa9e22eSVaclav Hapla if (ogidxs) { /* communicate global idxs */ 4968fa9e22eSVaclav Hapla PetscInt cum = 0, start, *work2; 4978fa9e22eSVaclav Hapla 4989566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &work)); 4999566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(N, &work2)); 5009371c9d4SSatish Balay for (r = 0; r < N; ++r) 5019371c9d4SSatish Balay if (idxs[r] >= 0) cum++; 5029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&cum, &start, 1, MPIU_INT, MPI_SUM, map->comm)); 5038fa9e22eSVaclav Hapla start -= cum; 5048fa9e22eSVaclav Hapla cum = 0; 5059371c9d4SSatish Balay for (r = 0; r < N; ++r) 5069371c9d4SSatish Balay if (idxs[r] >= 0) work2[r] = start + cum++; 5079566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf, MPIU_INT, work2, work, MPI_REPLACE)); 5089566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf, MPIU_INT, work2, work, MPI_REPLACE)); 5099566063dSJacob Faibussowitsch PetscCall(PetscFree(work2)); 5108fa9e22eSVaclav Hapla } 5119566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf)); 5128fa9e22eSVaclav Hapla /* Compress and put in indices */ 5138fa9e22eSVaclav Hapla for (r = 0; r < n; ++r) 5148fa9e22eSVaclav Hapla if (lidxs[r] >= 0) { 5158fa9e22eSVaclav Hapla if (work) work[len] = work[r]; 5168fa9e22eSVaclav Hapla lidxs[len++] = r; 5178fa9e22eSVaclav Hapla } 5188fa9e22eSVaclav Hapla if (on) *on = len; 5198fa9e22eSVaclav Hapla if (oidxs) *oidxs = lidxs; 5208fa9e22eSVaclav Hapla if (ogidxs) *ogidxs = work; 5218fa9e22eSVaclav Hapla PetscFunctionReturn(0); 5228fa9e22eSVaclav Hapla } 523deffd5ebSksagiyam 524deffd5ebSksagiyam /*@ 525deffd5ebSksagiyam PetscSFCreateByMatchingIndices - Create SF by matching root and leaf indices 526deffd5ebSksagiyam 527deffd5ebSksagiyam Collective 528deffd5ebSksagiyam 529deffd5ebSksagiyam Input Parameters: 530deffd5ebSksagiyam + layout - PetscLayout defining the global index space and the rank that brokers each index 531deffd5ebSksagiyam . numRootIndices - size of rootIndices 532deffd5ebSksagiyam . rootIndices - PetscInt array of global indices of which this process requests ownership 533deffd5ebSksagiyam . rootLocalIndices - root local index permutation (NULL if no permutation) 534deffd5ebSksagiyam . rootLocalOffset - offset to be added to root local indices 535deffd5ebSksagiyam . numLeafIndices - size of leafIndices 536deffd5ebSksagiyam . leafIndices - PetscInt array of global indices with which this process requires data associated 537deffd5ebSksagiyam . leafLocalIndices - leaf local index permutation (NULL if no permutation) 538deffd5ebSksagiyam - leafLocalOffset - offset to be added to leaf local indices 539deffd5ebSksagiyam 540d8d19677SJose E. Roman Output Parameters: 541deffd5ebSksagiyam + sfA - star forest representing the communication pattern from the layout space to the leaf space (NULL if not needed) 542deffd5ebSksagiyam - sf - star forest representing the communication pattern from the root space to the leaf space 543deffd5ebSksagiyam 544deffd5ebSksagiyam Example 1: 545deffd5ebSksagiyam $ 546deffd5ebSksagiyam $ rank : 0 1 2 547deffd5ebSksagiyam $ rootIndices : [1 0 2] [3] [3] 548deffd5ebSksagiyam $ rootLocalOffset : 100 200 300 549deffd5ebSksagiyam $ layout : [0 1] [2] [3] 550deffd5ebSksagiyam $ leafIndices : [0] [2] [0 3] 551deffd5ebSksagiyam $ leafLocalOffset : 400 500 600 552deffd5ebSksagiyam $ 5537ef5781fSVaclav Hapla would build the following SF 554deffd5ebSksagiyam $ 555deffd5ebSksagiyam $ [0] 400 <- (0,101) 556deffd5ebSksagiyam $ [1] 500 <- (0,102) 557deffd5ebSksagiyam $ [2] 600 <- (0,101) 558deffd5ebSksagiyam $ [2] 601 <- (2,300) 559deffd5ebSksagiyam $ 560deffd5ebSksagiyam Example 2: 561deffd5ebSksagiyam $ 562deffd5ebSksagiyam $ rank : 0 1 2 563deffd5ebSksagiyam $ rootIndices : [1 0 2] [3] [3] 564deffd5ebSksagiyam $ rootLocalOffset : 100 200 300 565deffd5ebSksagiyam $ layout : [0 1] [2] [3] 566deffd5ebSksagiyam $ leafIndices : rootIndices rootIndices rootIndices 567deffd5ebSksagiyam $ leafLocalOffset : rootLocalOffset rootLocalOffset rootLocalOffset 568deffd5ebSksagiyam $ 5697ef5781fSVaclav Hapla would build the following SF 570deffd5ebSksagiyam $ 571deffd5ebSksagiyam $ [1] 200 <- (2,300) 572deffd5ebSksagiyam $ 573deffd5ebSksagiyam Example 3: 574deffd5ebSksagiyam $ 575deffd5ebSksagiyam $ No process requests ownership of global index 1, but no process needs it. 576deffd5ebSksagiyam $ 577deffd5ebSksagiyam $ rank : 0 1 2 578deffd5ebSksagiyam $ numRootIndices : 2 1 1 579deffd5ebSksagiyam $ rootIndices : [0 2] [3] [3] 580deffd5ebSksagiyam $ rootLocalOffset : 100 200 300 581deffd5ebSksagiyam $ layout : [0 1] [2] [3] 582deffd5ebSksagiyam $ numLeafIndices : 1 1 2 583deffd5ebSksagiyam $ leafIndices : [0] [2] [0 3] 584deffd5ebSksagiyam $ leafLocalOffset : 400 500 600 585deffd5ebSksagiyam $ 5867ef5781fSVaclav Hapla would build the following SF 587deffd5ebSksagiyam $ 588deffd5ebSksagiyam $ [0] 400 <- (0,100) 589deffd5ebSksagiyam $ [1] 500 <- (0,101) 590deffd5ebSksagiyam $ [2] 600 <- (0,100) 591deffd5ebSksagiyam $ [2] 601 <- (2,300) 592deffd5ebSksagiyam $ 593deffd5ebSksagiyam 594deffd5ebSksagiyam Notes: 595deffd5ebSksagiyam The layout parameter represents any partitioning of [0, N), where N is the total number of global indices, and its 596deffd5ebSksagiyam local size can be set to PETSC_DECIDE. 597deffd5ebSksagiyam If a global index x lies in the partition owned by process i, each process whose rootIndices contains x requests 598deffd5ebSksagiyam ownership of x and sends its own rank and the local index of x to process i. 599deffd5ebSksagiyam If multiple processes request ownership of x, the one with the highest rank is to own x. 600deffd5ebSksagiyam Process i then broadcasts the ownership information, so that each process whose leafIndices contains x knows the 601deffd5ebSksagiyam ownership information of x. 602deffd5ebSksagiyam The output sf is constructed by associating each leaf point to a root point in this way. 603deffd5ebSksagiyam 604deffd5ebSksagiyam Suppose there is point data ordered according to the global indices and partitioned according to the given layout. 605deffd5ebSksagiyam The optional output PetscSF object sfA can be used to push such data to leaf points. 606deffd5ebSksagiyam 607deffd5ebSksagiyam All indices in rootIndices and leafIndices must lie in the layout range. The union (over all processes) of rootIndices 608deffd5ebSksagiyam must cover that of leafIndices, but need not cover the entire layout. 609deffd5ebSksagiyam 610deffd5ebSksagiyam If (leafIndices, leafLocalIndices, leafLocalOffset) == (rootIndices, rootLocalIndices, rootLocalOffset), the output 611deffd5ebSksagiyam star forest is almost identity, so will only include non-trivial part of the map. 612deffd5ebSksagiyam 613deffd5ebSksagiyam Developer Notes: 614deffd5ebSksagiyam Current approach of a process of the highest rank gaining the ownership may cause load imbalance; consider using 615deffd5ebSksagiyam hash(rank, root_local_index) as the bid for the ownership determination. 616deffd5ebSksagiyam 617deffd5ebSksagiyam Level: advanced 618deffd5ebSksagiyam 619db781477SPatrick Sanan .seealso: `PetscSFCreate()` 620deffd5ebSksagiyam @*/ 6219371c9d4SSatish 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) { 622deffd5ebSksagiyam MPI_Comm comm = layout->comm; 623deffd5ebSksagiyam PetscMPIInt size, rank; 624deffd5ebSksagiyam PetscSF sf1; 625deffd5ebSksagiyam PetscSFNode *owners, *buffer, *iremote; 626deffd5ebSksagiyam PetscInt *ilocal, nleaves, N, n, i; 627deffd5ebSksagiyam #if defined(PETSC_USE_DEBUG) 628deffd5ebSksagiyam PetscInt N1; 629deffd5ebSksagiyam #endif 630deffd5ebSksagiyam PetscBool flag; 631deffd5ebSksagiyam 632deffd5ebSksagiyam PetscFunctionBegin; 633deffd5ebSksagiyam if (rootIndices) PetscValidIntPointer(rootIndices, 3); 634deffd5ebSksagiyam if (rootLocalIndices) PetscValidIntPointer(rootLocalIndices, 4); 635deffd5ebSksagiyam if (leafIndices) PetscValidIntPointer(leafIndices, 7); 636deffd5ebSksagiyam if (leafLocalIndices) PetscValidIntPointer(leafLocalIndices, 8); 637deffd5ebSksagiyam if (sfA) PetscValidPointer(sfA, 10); 638deffd5ebSksagiyam PetscValidPointer(sf, 11); 63908401ef6SPierre Jolivet PetscCheck(numRootIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numRootIndices (%" PetscInt_FMT ") must be non-negative", numRootIndices); 64008401ef6SPierre Jolivet PetscCheck(numLeafIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numLeafIndices (%" PetscInt_FMT ") must be non-negative", numLeafIndices); 6419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 6429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 6439566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(layout)); 6449566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(layout, &N)); 6459566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(layout, &n)); 646deffd5ebSksagiyam flag = (PetscBool)(leafIndices == rootIndices); 6471c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &flag, 1, MPIU_BOOL, MPI_LAND, comm)); 648c9cc58a2SBarry Smith PetscCheck(!flag || numLeafIndices == numRootIndices, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "leafIndices == rootIndices, but numLeafIndices (%" PetscInt_FMT ") != numRootIndices(%" PetscInt_FMT ")", numLeafIndices, numRootIndices); 649deffd5ebSksagiyam #if defined(PETSC_USE_DEBUG) 650deffd5ebSksagiyam N1 = PETSC_MIN_INT; 6519371c9d4SSatish Balay for (i = 0; i < numRootIndices; i++) 6529371c9d4SSatish Balay if (rootIndices[i] > N1) N1 = rootIndices[i]; 6539566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm)); 65408401ef6SPierre 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); 655deffd5ebSksagiyam if (!flag) { 656deffd5ebSksagiyam N1 = PETSC_MIN_INT; 6579371c9d4SSatish Balay for (i = 0; i < numLeafIndices; i++) 6589371c9d4SSatish Balay if (leafIndices[i] > N1) N1 = leafIndices[i]; 6599566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm)); 66008401ef6SPierre 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); 661deffd5ebSksagiyam } 662deffd5ebSksagiyam #endif 663deffd5ebSksagiyam /* Reduce: owners -> buffer */ 6649566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &buffer)); 6659566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, &sf1)); 6669566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sf1)); 6679566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraphLayout(sf1, layout, numRootIndices, NULL, PETSC_OWN_POINTER, rootIndices)); 6689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numRootIndices, &owners)); 669deffd5ebSksagiyam for (i = 0; i < numRootIndices; ++i) { 670deffd5ebSksagiyam owners[i].rank = rank; 671deffd5ebSksagiyam owners[i].index = rootLocalOffset + (rootLocalIndices ? rootLocalIndices[i] : i); 672deffd5ebSksagiyam } 673deffd5ebSksagiyam for (i = 0; i < n; ++i) { 674deffd5ebSksagiyam buffer[i].index = -1; 675deffd5ebSksagiyam buffer[i].rank = -1; 676deffd5ebSksagiyam } 6779566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC)); 6789566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC)); 679deffd5ebSksagiyam /* Bcast: buffer -> owners */ 680deffd5ebSksagiyam if (!flag) { 681deffd5ebSksagiyam /* leafIndices is different from rootIndices */ 6829566063dSJacob Faibussowitsch PetscCall(PetscFree(owners)); 6839566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraphLayout(sf1, layout, numLeafIndices, NULL, PETSC_OWN_POINTER, leafIndices)); 6849566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numLeafIndices, &owners)); 685deffd5ebSksagiyam } 6869566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE)); 6879566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE)); 68808401ef6SPierre 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]); 6899566063dSJacob Faibussowitsch PetscCall(PetscFree(buffer)); 6909371c9d4SSatish Balay if (sfA) { 6919371c9d4SSatish Balay *sfA = sf1; 6929371c9d4SSatish Balay } else PetscCall(PetscSFDestroy(&sf1)); 693deffd5ebSksagiyam /* Create sf */ 694deffd5ebSksagiyam if (flag && rootLocalIndices == leafLocalIndices && leafLocalOffset == rootLocalOffset) { 695deffd5ebSksagiyam /* leaf space == root space */ 6969371c9d4SSatish Balay for (i = 0, nleaves = 0; i < numLeafIndices; ++i) 6979371c9d4SSatish Balay if (owners[i].rank != rank) ++nleaves; 6989566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &ilocal)); 6999566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &iremote)); 700deffd5ebSksagiyam for (i = 0, nleaves = 0; i < numLeafIndices; ++i) { 701deffd5ebSksagiyam if (owners[i].rank != rank) { 702deffd5ebSksagiyam ilocal[nleaves] = leafLocalOffset + i; 703deffd5ebSksagiyam iremote[nleaves].rank = owners[i].rank; 704deffd5ebSksagiyam iremote[nleaves].index = owners[i].index; 705deffd5ebSksagiyam ++nleaves; 706deffd5ebSksagiyam } 707deffd5ebSksagiyam } 7089566063dSJacob Faibussowitsch PetscCall(PetscFree(owners)); 709deffd5ebSksagiyam } else { 710deffd5ebSksagiyam nleaves = numLeafIndices; 7119566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &ilocal)); 712deffd5ebSksagiyam for (i = 0; i < nleaves; ++i) { ilocal[i] = leafLocalOffset + (leafLocalIndices ? leafLocalIndices[i] : i); } 713deffd5ebSksagiyam iremote = owners; 714deffd5ebSksagiyam } 7159566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, sf)); 7169566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(*sf)); 7179566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*sf, rootLocalOffset + numRootIndices, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER)); 718deffd5ebSksagiyam PetscFunctionReturn(0); 719deffd5ebSksagiyam } 720