1b0c7db22SLisandro Dalcin #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/ 2b0c7db22SLisandro Dalcin #include <petsc/private/sectionimpl.h> 3b0c7db22SLisandro Dalcin 4eb58282aSVaclav Hapla /*@ 5cab54364SBarry Smith PetscSFSetGraphLayout - Set a parallel star forest via global indices and a `PetscLayout` 6b0c7db22SLisandro Dalcin 7b0c7db22SLisandro Dalcin Collective 8b0c7db22SLisandro Dalcin 94165533cSJose E. Roman Input Parameters: 10b0c7db22SLisandro Dalcin + sf - star forest 11cab54364SBarry Smith . layout - `PetscLayout` defining the global space for roots 12b0c7db22SLisandro Dalcin . nleaves - number of leaf vertices on the current process, each of these references a root on any process 13b0c7db22SLisandro Dalcin . ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage 14b0c7db22SLisandro Dalcin . localmode - copy mode for ilocal 15eb58282aSVaclav Hapla - gremote - root vertices in global numbering corresponding to leaves in ilocal 16b0c7db22SLisandro Dalcin 17b0c7db22SLisandro Dalcin Level: intermediate 18b0c7db22SLisandro Dalcin 19cab54364SBarry Smith Note: 2086c56f52SVaclav Hapla Global indices must lie in [0, N) where N is the global size of layout. 21cab54364SBarry Smith Leaf indices in ilocal get sorted; this means the user-provided array gets sorted if localmode is `PETSC_OWN_POINTER`. 2286c56f52SVaclav Hapla 2338b5cf2dSJacob Faibussowitsch Developer Notes: 2486c56f52SVaclav Hapla Local indices which are the identity permutation in the range [0,nleaves) are discarded as they 25cab54364SBarry Smith encode contiguous storage. In such case, if localmode is `PETSC_OWN_POINTER`, the memory is deallocated as it is not 26b0c7db22SLisandro Dalcin needed 27b0c7db22SLisandro Dalcin 28cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFGetGraphLayout()`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()` 29b0c7db22SLisandro Dalcin @*/ 30d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetGraphLayout(PetscSF sf, PetscLayout layout, PetscInt nleaves, PetscInt *ilocal, PetscCopyMode localmode, const PetscInt *gremote) 31d71ae5a4SJacob Faibussowitsch { 3238a25198SStefano Zampini const PetscInt *range; 3338a25198SStefano Zampini PetscInt i, nroots, ls = -1, ln = -1; 3438a25198SStefano Zampini PetscMPIInt lr = -1; 35b0c7db22SLisandro Dalcin PetscSFNode *remote; 36b0c7db22SLisandro Dalcin 37b0c7db22SLisandro Dalcin PetscFunctionBegin; 38da0802e2SStefano Zampini PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 39b114984aSVaclav Hapla PetscCall(PetscLayoutSetUp(layout)); 409566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(layout, &nroots)); 419566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRanges(layout, &range)); 429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &remote)); 43eb58282aSVaclav Hapla if (nleaves) ls = gremote[0] + 1; 44b0c7db22SLisandro Dalcin for (i = 0; i < nleaves; i++) { 45eb58282aSVaclav Hapla const PetscInt idx = gremote[i] - ls; 4638a25198SStefano Zampini if (idx < 0 || idx >= ln) { /* short-circuit the search */ 47eb58282aSVaclav Hapla PetscCall(PetscLayoutFindOwnerIndex(layout, gremote[i], &lr, &remote[i].index)); 4838a25198SStefano Zampini remote[i].rank = lr; 4938a25198SStefano Zampini ls = range[lr]; 5038a25198SStefano Zampini ln = range[lr + 1] - ls; 5138a25198SStefano Zampini } else { 5238a25198SStefano Zampini remote[i].rank = lr; 5338a25198SStefano Zampini remote[i].index = idx; 5438a25198SStefano Zampini } 55b0c7db22SLisandro Dalcin } 569566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, nroots, nleaves, ilocal, localmode, remote, PETSC_OWN_POINTER)); 573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 58b0c7db22SLisandro Dalcin } 59b0c7db22SLisandro Dalcin 60eb58282aSVaclav Hapla /*@C 61cab54364SBarry Smith PetscSFGetGraphLayout - Get the global indices and `PetscLayout` that describe this star forest 62eb58282aSVaclav Hapla 63eb58282aSVaclav Hapla Collective 64eb58282aSVaclav Hapla 65eb58282aSVaclav Hapla Input Parameter: 66eb58282aSVaclav Hapla . sf - star forest 67eb58282aSVaclav Hapla 68eb58282aSVaclav Hapla Output Parameters: 69cab54364SBarry Smith + layout - `PetscLayout` defining the global space for roots 70eb58282aSVaclav Hapla . nleaves - number of leaf vertices on the current process, each of these references a root on any process 71f13dfd9eSBarry Smith . ilocal - locations of leaves in leafdata buffers, or `NULL` for contiguous storage 72eb58282aSVaclav Hapla - gremote - root vertices in global numbering corresponding to leaves in ilocal 73eb58282aSVaclav Hapla 74eb58282aSVaclav Hapla Level: intermediate 75eb58282aSVaclav Hapla 76eb58282aSVaclav Hapla Notes: 77eb58282aSVaclav Hapla The outputs are such that passing them as inputs to `PetscSFSetGraphLayout()` would lead to the same star forest. 78f13dfd9eSBarry Smith The outputs `layout` and `gremote` are freshly created each time this function is called, 79f13dfd9eSBarry Smith so they need to be freed (with `PetscLayoutDestroy()` and `PetscFree()`) by the user. 80eb58282aSVaclav Hapla 81cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFSetGraphLayout()`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()` 82eb58282aSVaclav Hapla @*/ 83d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetGraphLayout(PetscSF sf, PetscLayout *layout, PetscInt *nleaves, const PetscInt *ilocal[], PetscInt *gremote[]) 84d71ae5a4SJacob Faibussowitsch { 85eb58282aSVaclav Hapla PetscInt nr, nl; 86eb58282aSVaclav Hapla const PetscSFNode *ir; 87eb58282aSVaclav Hapla PetscLayout lt; 88eb58282aSVaclav Hapla 89eb58282aSVaclav Hapla PetscFunctionBegin; 90eb58282aSVaclav Hapla PetscCall(PetscSFGetGraph(sf, &nr, &nl, ilocal, &ir)); 91eb58282aSVaclav Hapla PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)sf), nr, PETSC_DECIDE, 1, <)); 92eb58282aSVaclav Hapla if (gremote) { 93eb58282aSVaclav Hapla PetscInt i; 94eb58282aSVaclav Hapla const PetscInt *range; 95eb58282aSVaclav Hapla PetscInt *gr; 96eb58282aSVaclav Hapla 97eb58282aSVaclav Hapla PetscCall(PetscLayoutGetRanges(lt, &range)); 98eb58282aSVaclav Hapla PetscCall(PetscMalloc1(nl, &gr)); 99eb58282aSVaclav Hapla for (i = 0; i < nl; i++) gr[i] = range[ir[i].rank] + ir[i].index; 100eb58282aSVaclav Hapla *gremote = gr; 101eb58282aSVaclav Hapla } 102eb58282aSVaclav Hapla if (nleaves) *nleaves = nl; 103eb58282aSVaclav Hapla if (layout) *layout = lt; 104eb58282aSVaclav Hapla else PetscCall(PetscLayoutDestroy(<)); 1053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 106eb58282aSVaclav Hapla } 107eb58282aSVaclav Hapla 108b0c7db22SLisandro Dalcin /*@ 109cab54364SBarry Smith PetscSFSetGraphSection - Sets the `PetscSF` graph encoding the parallel dof overlap based upon the `PetscSection` describing the data layout. 110b0c7db22SLisandro Dalcin 111b0c7db22SLisandro Dalcin Input Parameters: 112cab54364SBarry Smith + sf - The `PetscSF` 113cab54364SBarry Smith . localSection - `PetscSection` describing the local data layout 114cab54364SBarry Smith - globalSection - `PetscSection` describing the global data layout 115b0c7db22SLisandro Dalcin 116b0c7db22SLisandro Dalcin Level: developer 117b0c7db22SLisandro Dalcin 118cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFSetGraph()`, `PetscSFSetGraphLayout()` 119b0c7db22SLisandro Dalcin @*/ 120d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetGraphSection(PetscSF sf, PetscSection localSection, PetscSection globalSection) 121d71ae5a4SJacob Faibussowitsch { 122b0c7db22SLisandro Dalcin MPI_Comm comm; 123b0c7db22SLisandro Dalcin PetscLayout layout; 124b0c7db22SLisandro Dalcin const PetscInt *ranges; 125b0c7db22SLisandro Dalcin PetscInt *local; 126b0c7db22SLisandro Dalcin PetscSFNode *remote; 127b0c7db22SLisandro Dalcin PetscInt pStart, pEnd, p, nroots, nleaves = 0, l; 128b0c7db22SLisandro Dalcin PetscMPIInt size, rank; 129b0c7db22SLisandro Dalcin 130b0c7db22SLisandro Dalcin PetscFunctionBegin; 131b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 132b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(localSection, PETSC_SECTION_CLASSID, 2); 133b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(globalSection, PETSC_SECTION_CLASSID, 3); 134b0c7db22SLisandro Dalcin 1359566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 1369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 1379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1389566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(globalSection, &pStart, &pEnd)); 1399566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstrainedStorageSize(globalSection, &nroots)); 1409566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(comm, &layout)); 1419566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(layout, 1)); 1429566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetLocalSize(layout, nroots)); 1439566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(layout)); 1449566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRanges(layout, &ranges)); 145b0c7db22SLisandro Dalcin for (p = pStart; p < pEnd; ++p) { 146b0c7db22SLisandro Dalcin PetscInt gdof, gcdof; 147b0c7db22SLisandro Dalcin 1489566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(globalSection, p, &gdof)); 1499566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof)); 15057508eceSPierre Jolivet PetscCheck(gcdof <= (gdof < 0 ? -(gdof + 1) : gdof), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " has %" PetscInt_FMT " constraints > %" PetscInt_FMT " dof", p, gcdof, gdof < 0 ? -(gdof + 1) : gdof); 151b0c7db22SLisandro Dalcin nleaves += gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof; 152b0c7db22SLisandro Dalcin } 1539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &local)); 1549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &remote)); 155b0c7db22SLisandro Dalcin for (p = pStart, l = 0; p < pEnd; ++p) { 156b0c7db22SLisandro Dalcin const PetscInt *cind; 157b0c7db22SLisandro Dalcin PetscInt dof, cdof, off, gdof, gcdof, goff, gsize, d, c; 158b0c7db22SLisandro Dalcin 1599566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(localSection, p, &dof)); 1609566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(localSection, p, &off)); 1619566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(localSection, p, &cdof)); 1629566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintIndices(localSection, p, &cind)); 1639566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(globalSection, p, &gdof)); 1649566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof)); 1659566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(globalSection, p, &goff)); 166b0c7db22SLisandro Dalcin if (!gdof) continue; /* Censored point */ 167b0c7db22SLisandro Dalcin gsize = gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof; 168b0c7db22SLisandro Dalcin if (gsize != dof - cdof) { 16908401ef6SPierre Jolivet PetscCheck(gsize == dof, comm, PETSC_ERR_ARG_WRONG, "Global dof %" PetscInt_FMT " for point %" PetscInt_FMT " is neither the constrained size %" PetscInt_FMT ", nor the unconstrained %" PetscInt_FMT, gsize, p, dof - cdof, dof); 170b0c7db22SLisandro Dalcin cdof = 0; /* Ignore constraints */ 171b0c7db22SLisandro Dalcin } 172b0c7db22SLisandro Dalcin for (d = 0, c = 0; d < dof; ++d) { 1739371c9d4SSatish Balay if ((c < cdof) && (cind[c] == d)) { 1749371c9d4SSatish Balay ++c; 1759371c9d4SSatish Balay continue; 1769371c9d4SSatish Balay } 177b0c7db22SLisandro Dalcin local[l + d - c] = off + d; 178b0c7db22SLisandro Dalcin } 17908401ef6SPierre Jolivet PetscCheck(d - c == gsize, comm, PETSC_ERR_ARG_WRONG, "Point %" PetscInt_FMT ": Global dof %" PetscInt_FMT " != %" PetscInt_FMT " size - number of constraints", p, gsize, d - c); 180b0c7db22SLisandro Dalcin if (gdof < 0) { 181b0c7db22SLisandro Dalcin for (d = 0; d < gsize; ++d, ++l) { 1826497c311SBarry Smith PetscInt offset = -(goff + 1) + d, ir; 1836497c311SBarry Smith PetscMPIInt r; 184b0c7db22SLisandro Dalcin 1856497c311SBarry Smith PetscCall(PetscFindInt(offset, size + 1, ranges, &ir)); 1866497c311SBarry Smith PetscCall(PetscMPIIntCast(ir, &r)); 187b0c7db22SLisandro Dalcin if (r < 0) r = -(r + 2); 1886497c311SBarry Smith PetscCheck(!(r < 0) && !(r >= size), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " mapped to invalid process %d (%" PetscInt_FMT ", %" PetscInt_FMT ")", p, r, gdof, goff); 189b0c7db22SLisandro Dalcin remote[l].rank = r; 190b0c7db22SLisandro Dalcin remote[l].index = offset - ranges[r]; 191b0c7db22SLisandro Dalcin } 192b0c7db22SLisandro Dalcin } else { 193b0c7db22SLisandro Dalcin for (d = 0; d < gsize; ++d, ++l) { 194b0c7db22SLisandro Dalcin remote[l].rank = rank; 195b0c7db22SLisandro Dalcin remote[l].index = goff + d - ranges[rank]; 196b0c7db22SLisandro Dalcin } 197b0c7db22SLisandro Dalcin } 198b0c7db22SLisandro Dalcin } 19908401ef6SPierre Jolivet PetscCheck(l == nleaves, comm, PETSC_ERR_PLIB, "Iteration error, l %" PetscInt_FMT " != nleaves %" PetscInt_FMT, l, nleaves); 2009566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&layout)); 2019566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER)); 2023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 203b0c7db22SLisandro Dalcin } 204b0c7db22SLisandro Dalcin 205b0c7db22SLisandro Dalcin /*@C 206cab54364SBarry Smith PetscSFDistributeSection - Create a new `PetscSection` reorganized, moving from the root to the leaves of the `PetscSF` 207b0c7db22SLisandro Dalcin 208c3339decSBarry Smith Collective 209b0c7db22SLisandro Dalcin 210b0c7db22SLisandro Dalcin Input Parameters: 211cab54364SBarry Smith + sf - The `PetscSF` 212b0c7db22SLisandro Dalcin - rootSection - Section defined on root space 213b0c7db22SLisandro Dalcin 214b0c7db22SLisandro Dalcin Output Parameters: 2156497c311SBarry Smith + remoteOffsets - root offsets in leaf storage, or `NULL` 216b0c7db22SLisandro Dalcin - leafSection - Section defined on the leaf space 217b0c7db22SLisandro Dalcin 218b0c7db22SLisandro Dalcin Level: advanced 219b0c7db22SLisandro Dalcin 22023e9620eSJunchao Zhang Fortran Notes: 22123e9620eSJunchao Zhang In Fortran, use PetscSFDistributeSectionF90() 22223e9620eSJunchao Zhang 223cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()` 224b0c7db22SLisandro Dalcin @*/ 225d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFDistributeSection(PetscSF sf, PetscSection rootSection, PetscInt **remoteOffsets, PetscSection leafSection) 226d71ae5a4SJacob Faibussowitsch { 227b0c7db22SLisandro Dalcin PetscSF embedSF; 228b0c7db22SLisandro Dalcin const PetscInt *indices; 229b0c7db22SLisandro Dalcin IS selected; 2301690c2aeSBarry Smith PetscInt numFields, nroots, rpStart, rpEnd, lpStart = PETSC_INT_MAX, lpEnd = -1, f, c; 231b0c7db22SLisandro Dalcin PetscBool *sub, hasc; 232*835f2295SStefano Zampini PetscMPIInt msize; 233b0c7db22SLisandro Dalcin 234b0c7db22SLisandro Dalcin PetscFunctionBegin; 2359566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_DistSect, sf, 0, 0, 0)); 2369566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(rootSection, &numFields)); 237029e8818Sksagiyam if (numFields) { 238029e8818Sksagiyam IS perm; 239029e8818Sksagiyam 240029e8818Sksagiyam /* PetscSectionSetNumFields() calls PetscSectionReset(), which destroys 241029e8818Sksagiyam leafSection->perm. To keep this permutation set by the user, we grab 242029e8818Sksagiyam the reference before calling PetscSectionSetNumFields() and set it 243029e8818Sksagiyam back after. */ 2449566063dSJacob Faibussowitsch PetscCall(PetscSectionGetPermutation(leafSection, &perm)); 2459566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)perm)); 2469566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(leafSection, numFields)); 2479566063dSJacob Faibussowitsch PetscCall(PetscSectionSetPermutation(leafSection, perm)); 2489566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 249029e8818Sksagiyam } 2509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numFields + 2, &sub)); 251b0c7db22SLisandro Dalcin sub[1] = rootSection->bc ? PETSC_TRUE : PETSC_FALSE; 252b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) { 2532ee82556SMatthew G. Knepley PetscSectionSym sym, dsym = NULL; 254b0c7db22SLisandro Dalcin const char *name = NULL; 255b0c7db22SLisandro Dalcin PetscInt numComp = 0; 256b0c7db22SLisandro Dalcin 257b0c7db22SLisandro Dalcin sub[2 + f] = rootSection->field[f]->bc ? PETSC_TRUE : PETSC_FALSE; 2589566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldComponents(rootSection, f, &numComp)); 2599566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldName(rootSection, f, &name)); 2609566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldSym(rootSection, f, &sym)); 2619566063dSJacob Faibussowitsch if (sym) PetscCall(PetscSectionSymDistribute(sym, sf, &dsym)); 2629566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(leafSection, f, numComp)); 2639566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldName(leafSection, f, name)); 2649566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldSym(leafSection, f, dsym)); 2659566063dSJacob Faibussowitsch PetscCall(PetscSectionSymDestroy(&dsym)); 266b0c7db22SLisandro Dalcin for (c = 0; c < rootSection->numFieldComponents[f]; ++c) { 2679566063dSJacob Faibussowitsch PetscCall(PetscSectionGetComponentName(rootSection, f, c, &name)); 2689566063dSJacob Faibussowitsch PetscCall(PetscSectionSetComponentName(leafSection, f, c, name)); 269b0c7db22SLisandro Dalcin } 270b0c7db22SLisandro Dalcin } 2719566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd)); 2729566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL)); 273b0c7db22SLisandro Dalcin rpEnd = PetscMin(rpEnd, nroots); 274b0c7db22SLisandro Dalcin rpEnd = PetscMax(rpStart, rpEnd); 275b0c7db22SLisandro Dalcin /* see if we can avoid creating the embedded SF, since it can cost more than an allreduce */ 276b0c7db22SLisandro Dalcin sub[0] = (PetscBool)(nroots != rpEnd - rpStart); 277*835f2295SStefano Zampini PetscCall(PetscMPIIntCast(2 + numFields, &msize)); 278*835f2295SStefano Zampini PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, sub, msize, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)sf))); 279b0c7db22SLisandro Dalcin if (sub[0]) { 2809566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected)); 2819566063dSJacob Faibussowitsch PetscCall(ISGetIndices(selected, &indices)); 2829566063dSJacob Faibussowitsch PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF)); 2839566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(selected, &indices)); 2849566063dSJacob Faibussowitsch PetscCall(ISDestroy(&selected)); 285b0c7db22SLisandro Dalcin } else { 2869566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)sf)); 287b0c7db22SLisandro Dalcin embedSF = sf; 288b0c7db22SLisandro Dalcin } 2899566063dSJacob Faibussowitsch PetscCall(PetscSFGetLeafRange(embedSF, &lpStart, &lpEnd)); 290b0c7db22SLisandro Dalcin lpEnd++; 291b0c7db22SLisandro Dalcin 2929566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(leafSection, lpStart, lpEnd)); 293b0c7db22SLisandro Dalcin 294b0c7db22SLisandro Dalcin /* Constrained dof section */ 295b0c7db22SLisandro Dalcin hasc = sub[1]; 296b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) hasc = (PetscBool)(hasc || sub[2 + f]); 297b0c7db22SLisandro Dalcin 298b0c7db22SLisandro Dalcin /* Could fuse these at the cost of copies and extra allocation */ 29916cd844bSPierre Jolivet PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasDof, -rpStart), PetscSafePointerPlusOffset(leafSection->atlasDof, -lpStart), MPI_REPLACE)); 30016cd844bSPierre Jolivet PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasDof, -rpStart), PetscSafePointerPlusOffset(leafSection->atlasDof, -lpStart), MPI_REPLACE)); 301b0c7db22SLisandro Dalcin if (sub[1]) { 3027c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(rootSection)); 3037c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(leafSection)); 3049566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE)); 3059566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE)); 306b0c7db22SLisandro Dalcin } 307b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) { 30816cd844bSPierre Jolivet PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->field[f]->atlasDof, -rpStart), PetscSafePointerPlusOffset(leafSection->field[f]->atlasDof, -lpStart), MPI_REPLACE)); 30916cd844bSPierre Jolivet PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->field[f]->atlasDof, -rpStart), PetscSafePointerPlusOffset(leafSection->field[f]->atlasDof, -lpStart), MPI_REPLACE)); 310b0c7db22SLisandro Dalcin if (sub[2 + f]) { 3117c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(rootSection->field[f])); 3127c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(leafSection->field[f])); 3139566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE)); 3149566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE)); 315b0c7db22SLisandro Dalcin } 316b0c7db22SLisandro Dalcin } 317b0c7db22SLisandro Dalcin if (remoteOffsets) { 3189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lpEnd - lpStart, remoteOffsets)); 31916cd844bSPierre Jolivet PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(*remoteOffsets, -lpStart), MPI_REPLACE)); 32016cd844bSPierre Jolivet PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(*remoteOffsets, -lpStart), MPI_REPLACE)); 321b0c7db22SLisandro Dalcin } 32269c11d05SVaclav Hapla PetscCall(PetscSectionInvalidateMaxDof_Internal(leafSection)); 3239566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(leafSection)); 324b0c7db22SLisandro Dalcin if (hasc) { /* need to communicate bcIndices */ 325b0c7db22SLisandro Dalcin PetscSF bcSF; 326b0c7db22SLisandro Dalcin PetscInt *rOffBc; 327b0c7db22SLisandro Dalcin 3289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lpEnd - lpStart, &rOffBc)); 329b0c7db22SLisandro Dalcin if (sub[1]) { 3309566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 3319566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 3329566063dSJacob Faibussowitsch PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->bc, rOffBc, leafSection->bc, &bcSF)); 3339566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE)); 3349566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE)); 3359566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&bcSF)); 336b0c7db22SLisandro Dalcin } 337b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) { 338b0c7db22SLisandro Dalcin if (sub[2 + f]) { 3399566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 3409566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 3419566063dSJacob Faibussowitsch PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->field[f]->bc, rOffBc, leafSection->field[f]->bc, &bcSF)); 3429566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE)); 3439566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE)); 3449566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&bcSF)); 345b0c7db22SLisandro Dalcin } 346b0c7db22SLisandro Dalcin } 3479566063dSJacob Faibussowitsch PetscCall(PetscFree(rOffBc)); 348b0c7db22SLisandro Dalcin } 3499566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&embedSF)); 3509566063dSJacob Faibussowitsch PetscCall(PetscFree(sub)); 3519566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_DistSect, sf, 0, 0, 0)); 3523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 353b0c7db22SLisandro Dalcin } 354b0c7db22SLisandro Dalcin 355b0c7db22SLisandro Dalcin /*@C 356b0c7db22SLisandro Dalcin PetscSFCreateRemoteOffsets - Create offsets for point data on remote processes 357b0c7db22SLisandro Dalcin 358c3339decSBarry Smith Collective 359b0c7db22SLisandro Dalcin 360b0c7db22SLisandro Dalcin Input Parameters: 361cab54364SBarry Smith + sf - The `PetscSF` 362b0c7db22SLisandro Dalcin . rootSection - Data layout of remote points for outgoing data (this is layout for SF roots) 363b0c7db22SLisandro Dalcin - leafSection - Data layout of local points for incoming data (this is layout for SF leaves) 364b0c7db22SLisandro Dalcin 365b0c7db22SLisandro Dalcin Output Parameter: 366b0c7db22SLisandro Dalcin . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL 367b0c7db22SLisandro Dalcin 368b0c7db22SLisandro Dalcin Level: developer 369b0c7db22SLisandro Dalcin 37023e9620eSJunchao Zhang Fortran Notes: 37123e9620eSJunchao Zhang In Fortran, use PetscSFCreateRemoteOffsetsF90() 37223e9620eSJunchao Zhang 373cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()` 374b0c7db22SLisandro Dalcin @*/ 375d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateRemoteOffsets(PetscSF sf, PetscSection rootSection, PetscSection leafSection, PetscInt **remoteOffsets) 376d71ae5a4SJacob Faibussowitsch { 377b0c7db22SLisandro Dalcin PetscSF embedSF; 378b0c7db22SLisandro Dalcin const PetscInt *indices; 379b0c7db22SLisandro Dalcin IS selected; 380b0c7db22SLisandro Dalcin PetscInt numRoots, rpStart = 0, rpEnd = 0, lpStart = 0, lpEnd = 0; 381b0c7db22SLisandro Dalcin 382b0c7db22SLisandro Dalcin PetscFunctionBegin; 383b0c7db22SLisandro Dalcin *remoteOffsets = NULL; 3849566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &numRoots, NULL, NULL, NULL)); 3853ba16761SJacob Faibussowitsch if (numRoots < 0) PetscFunctionReturn(PETSC_SUCCESS); 3869566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_RemoteOff, sf, 0, 0, 0)); 3879566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd)); 3889566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd)); 3899566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected)); 3909566063dSJacob Faibussowitsch PetscCall(ISGetIndices(selected, &indices)); 3919566063dSJacob Faibussowitsch PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF)); 3929566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(selected, &indices)); 3939566063dSJacob Faibussowitsch PetscCall(ISDestroy(&selected)); 3949566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(lpEnd - lpStart, remoteOffsets)); 3958e3a54c0SPierre Jolivet PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(*remoteOffsets, -lpStart), MPI_REPLACE)); 3968e3a54c0SPierre Jolivet PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(*remoteOffsets, -lpStart), MPI_REPLACE)); 3979566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&embedSF)); 3989566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_RemoteOff, sf, 0, 0, 0)); 3993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 400b0c7db22SLisandro Dalcin } 401b0c7db22SLisandro Dalcin 402b0c7db22SLisandro Dalcin /*@C 403cab54364SBarry Smith PetscSFCreateSectionSF - Create an expanded `PetscSF` of dofs, assuming the input `PetscSF` relates points 404b0c7db22SLisandro Dalcin 405c3339decSBarry Smith Collective 406b0c7db22SLisandro Dalcin 407b0c7db22SLisandro Dalcin Input Parameters: 408cab54364SBarry Smith + sf - The `PetscSF` 409b0c7db22SLisandro Dalcin . rootSection - Data layout of remote points for outgoing data (this is usually the serial section) 410b0c7db22SLisandro Dalcin . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL 411b0c7db22SLisandro Dalcin - leafSection - Data layout of local points for incoming data (this is the distributed section) 412b0c7db22SLisandro Dalcin 4132fe279fdSBarry Smith Output Parameter: 4142fe279fdSBarry Smith . sectionSF - The new `PetscSF` 415b0c7db22SLisandro Dalcin 416b0c7db22SLisandro Dalcin Level: advanced 417b0c7db22SLisandro Dalcin 41823e9620eSJunchao Zhang Notes: 419021e28c0SJames Wright `remoteOffsets` can be NULL if `sf` does not reference any points in leafSection 420cab54364SBarry Smith 42123e9620eSJunchao Zhang Fortran Notes: 42223e9620eSJunchao Zhang In Fortran, use PetscSFCreateSectionSFF90() 42323e9620eSJunchao Zhang 424cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()` 425b0c7db22SLisandro Dalcin @*/ 426d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateSectionSF(PetscSF sf, PetscSection rootSection, PetscInt remoteOffsets[], PetscSection leafSection, PetscSF *sectionSF) 427d71ae5a4SJacob Faibussowitsch { 428b0c7db22SLisandro Dalcin MPI_Comm comm; 429b0c7db22SLisandro Dalcin const PetscInt *localPoints; 430b0c7db22SLisandro Dalcin const PetscSFNode *remotePoints; 431b0c7db22SLisandro Dalcin PetscInt lpStart, lpEnd; 432b0c7db22SLisandro Dalcin PetscInt numRoots, numSectionRoots, numPoints, numIndices = 0; 433b0c7db22SLisandro Dalcin PetscInt *localIndices; 434b0c7db22SLisandro Dalcin PetscSFNode *remoteIndices; 435b0c7db22SLisandro Dalcin PetscInt i, ind; 436b0c7db22SLisandro Dalcin 437b0c7db22SLisandro Dalcin PetscFunctionBegin; 438b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 4394f572ea9SToby Isaac PetscAssertPointer(rootSection, 2); 4404f572ea9SToby Isaac /* Cannot check PetscAssertPointer(remoteOffsets,3) because it can be NULL if sf does not reference any points in leafSection */ 4414f572ea9SToby Isaac PetscAssertPointer(leafSection, 4); 4424f572ea9SToby Isaac PetscAssertPointer(sectionSF, 5); 4439566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 4449566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, sectionSF)); 4459566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd)); 4469566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(rootSection, &numSectionRoots)); 4479566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &numRoots, &numPoints, &localPoints, &remotePoints)); 4483ba16761SJacob Faibussowitsch if (numRoots < 0) PetscFunctionReturn(PETSC_SUCCESS); 4499566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_SectSF, sf, 0, 0, 0)); 450b0c7db22SLisandro Dalcin for (i = 0; i < numPoints; ++i) { 451b0c7db22SLisandro Dalcin PetscInt localPoint = localPoints ? localPoints[i] : i; 452b0c7db22SLisandro Dalcin PetscInt dof; 453b0c7db22SLisandro Dalcin 454b0c7db22SLisandro Dalcin if ((localPoint >= lpStart) && (localPoint < lpEnd)) { 4559566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof)); 45622a18c9fSMatthew G. Knepley numIndices += dof < 0 ? 0 : dof; 457b0c7db22SLisandro Dalcin } 458b0c7db22SLisandro Dalcin } 4599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numIndices, &localIndices)); 4609566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numIndices, &remoteIndices)); 461b0c7db22SLisandro Dalcin /* Create new index graph */ 462b0c7db22SLisandro Dalcin for (i = 0, ind = 0; i < numPoints; ++i) { 463b0c7db22SLisandro Dalcin PetscInt localPoint = localPoints ? localPoints[i] : i; 464*835f2295SStefano Zampini PetscInt rank = remotePoints[i].rank; 465b0c7db22SLisandro Dalcin 466b0c7db22SLisandro Dalcin if ((localPoint >= lpStart) && (localPoint < lpEnd)) { 467b0c7db22SLisandro Dalcin PetscInt remoteOffset = remoteOffsets[localPoint - lpStart]; 468b0c7db22SLisandro Dalcin PetscInt loff, dof, d; 469b0c7db22SLisandro Dalcin 4709566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(leafSection, localPoint, &loff)); 4719566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof)); 472b0c7db22SLisandro Dalcin for (d = 0; d < dof; ++d, ++ind) { 473b0c7db22SLisandro Dalcin localIndices[ind] = loff + d; 474b0c7db22SLisandro Dalcin remoteIndices[ind].rank = rank; 475b0c7db22SLisandro Dalcin remoteIndices[ind].index = remoteOffset + d; 476b0c7db22SLisandro Dalcin } 477b0c7db22SLisandro Dalcin } 478b0c7db22SLisandro Dalcin } 47908401ef6SPierre Jolivet PetscCheck(numIndices == ind, comm, PETSC_ERR_PLIB, "Inconsistency in indices, %" PetscInt_FMT " should be %" PetscInt_FMT, ind, numIndices); 4809566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*sectionSF, numSectionRoots, numIndices, localIndices, PETSC_OWN_POINTER, remoteIndices, PETSC_OWN_POINTER)); 4819566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(*sectionSF)); 4829566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_SectSF, sf, 0, 0, 0)); 4833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 484b0c7db22SLisandro Dalcin } 4858fa9e22eSVaclav Hapla 4868fa9e22eSVaclav Hapla /*@C 487cab54364SBarry Smith PetscSFCreateFromLayouts - Creates a parallel star forest mapping two `PetscLayout` objects 4888fa9e22eSVaclav Hapla 4898fa9e22eSVaclav Hapla Collective 4908fa9e22eSVaclav Hapla 4914165533cSJose E. Roman Input Parameters: 492cab54364SBarry Smith + rmap - `PetscLayout` defining the global root space 493cab54364SBarry Smith - lmap - `PetscLayout` defining the global leaf space 4948fa9e22eSVaclav Hapla 4954165533cSJose E. Roman Output Parameter: 4968fa9e22eSVaclav Hapla . sf - The parallel star forest 4978fa9e22eSVaclav Hapla 4988fa9e22eSVaclav Hapla Level: intermediate 4998fa9e22eSVaclav Hapla 500cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()`, `PetscLayoutCreate()`, `PetscSFSetGraphLayout()` 5018fa9e22eSVaclav Hapla @*/ 502d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateFromLayouts(PetscLayout rmap, PetscLayout lmap, PetscSF *sf) 503d71ae5a4SJacob Faibussowitsch { 5048fa9e22eSVaclav Hapla PetscInt i, nroots, nleaves = 0; 5058fa9e22eSVaclav Hapla PetscInt rN, lst, len; 5068fa9e22eSVaclav Hapla PetscMPIInt owner = -1; 5078fa9e22eSVaclav Hapla PetscSFNode *remote; 5088fa9e22eSVaclav Hapla MPI_Comm rcomm = rmap->comm; 5098fa9e22eSVaclav Hapla MPI_Comm lcomm = lmap->comm; 5108fa9e22eSVaclav Hapla PetscMPIInt flg; 5118fa9e22eSVaclav Hapla 5128fa9e22eSVaclav Hapla PetscFunctionBegin; 5134f572ea9SToby Isaac PetscAssertPointer(sf, 3); 51428b400f6SJacob Faibussowitsch PetscCheck(rmap->setupcalled, rcomm, PETSC_ERR_ARG_WRONGSTATE, "Root layout not setup"); 51528b400f6SJacob Faibussowitsch PetscCheck(lmap->setupcalled, lcomm, PETSC_ERR_ARG_WRONGSTATE, "Leaf layout not setup"); 5169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(rcomm, lcomm, &flg)); 517c9cc58a2SBarry Smith PetscCheck(flg == MPI_CONGRUENT || flg == MPI_IDENT, rcomm, PETSC_ERR_SUP, "cannot map two layouts with non-matching communicators"); 5189566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(rcomm, sf)); 5199566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(rmap, &nroots)); 5209566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(rmap, &rN)); 5219566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRange(lmap, &lst, &len)); 5229566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len - lst, &remote)); 5238fa9e22eSVaclav Hapla for (i = lst; i < len && i < rN; i++) { 52448a46eb9SPierre Jolivet if (owner < -1 || i >= rmap->range[owner + 1]) PetscCall(PetscLayoutFindOwner(rmap, i, &owner)); 5258fa9e22eSVaclav Hapla remote[nleaves].rank = owner; 5268fa9e22eSVaclav Hapla remote[nleaves].index = i - rmap->range[owner]; 5278fa9e22eSVaclav Hapla nleaves++; 5288fa9e22eSVaclav Hapla } 5299566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*sf, nroots, nleaves, NULL, PETSC_OWN_POINTER, remote, PETSC_COPY_VALUES)); 5309566063dSJacob Faibussowitsch PetscCall(PetscFree(remote)); 5313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5328fa9e22eSVaclav Hapla } 5338fa9e22eSVaclav Hapla 5348fa9e22eSVaclav Hapla /* TODO: handle nooffprocentries like MatZeroRowsMapLocal_Private, since this code is the same */ 535d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLayoutMapLocal(PetscLayout map, PetscInt N, const PetscInt idxs[], PetscInt *on, PetscInt **oidxs, PetscInt **ogidxs) 536d71ae5a4SJacob Faibussowitsch { 5378fa9e22eSVaclav Hapla PetscInt *owners = map->range; 5388fa9e22eSVaclav Hapla PetscInt n = map->n; 5398fa9e22eSVaclav Hapla PetscSF sf; 5408fa9e22eSVaclav Hapla PetscInt *lidxs, *work = NULL; 5418fa9e22eSVaclav Hapla PetscSFNode *ridxs; 5428fa9e22eSVaclav Hapla PetscMPIInt rank, p = 0; 5438fa9e22eSVaclav Hapla PetscInt r, len = 0; 5448fa9e22eSVaclav Hapla 5458fa9e22eSVaclav Hapla PetscFunctionBegin; 5468fa9e22eSVaclav Hapla if (on) *on = 0; /* squelch -Wmaybe-uninitialized */ 5478fa9e22eSVaclav Hapla /* Create SF where leaves are input idxs and roots are owned idxs */ 5489566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(map->comm, &rank)); 5499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &lidxs)); 5508fa9e22eSVaclav Hapla for (r = 0; r < n; ++r) lidxs[r] = -1; 5519566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N, &ridxs)); 5528fa9e22eSVaclav Hapla for (r = 0; r < N; ++r) { 5538fa9e22eSVaclav Hapla const PetscInt idx = idxs[r]; 554c9cc58a2SBarry 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); 5558fa9e22eSVaclav Hapla if (idx < owners[p] || owners[p + 1] <= idx) { /* short-circuit the search if the last p owns this idx too */ 5569566063dSJacob Faibussowitsch PetscCall(PetscLayoutFindOwner(map, idx, &p)); 5578fa9e22eSVaclav Hapla } 5588fa9e22eSVaclav Hapla ridxs[r].rank = p; 5598fa9e22eSVaclav Hapla ridxs[r].index = idxs[r] - owners[p]; 5608fa9e22eSVaclav Hapla } 5619566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(map->comm, &sf)); 5629566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, n, N, NULL, PETSC_OWN_POINTER, ridxs, PETSC_OWN_POINTER)); 5639566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR)); 5649566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR)); 5658fa9e22eSVaclav Hapla if (ogidxs) { /* communicate global idxs */ 5668fa9e22eSVaclav Hapla PetscInt cum = 0, start, *work2; 5678fa9e22eSVaclav Hapla 5689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &work)); 5699566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(N, &work2)); 5709371c9d4SSatish Balay for (r = 0; r < N; ++r) 5719371c9d4SSatish Balay if (idxs[r] >= 0) cum++; 5729566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&cum, &start, 1, MPIU_INT, MPI_SUM, map->comm)); 5738fa9e22eSVaclav Hapla start -= cum; 5748fa9e22eSVaclav Hapla cum = 0; 5759371c9d4SSatish Balay for (r = 0; r < N; ++r) 5769371c9d4SSatish Balay if (idxs[r] >= 0) work2[r] = start + cum++; 5779566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf, MPIU_INT, work2, work, MPI_REPLACE)); 5789566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf, MPIU_INT, work2, work, MPI_REPLACE)); 5799566063dSJacob Faibussowitsch PetscCall(PetscFree(work2)); 5808fa9e22eSVaclav Hapla } 5819566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf)); 5828fa9e22eSVaclav Hapla /* Compress and put in indices */ 5838fa9e22eSVaclav Hapla for (r = 0; r < n; ++r) 5848fa9e22eSVaclav Hapla if (lidxs[r] >= 0) { 5858fa9e22eSVaclav Hapla if (work) work[len] = work[r]; 5868fa9e22eSVaclav Hapla lidxs[len++] = r; 5878fa9e22eSVaclav Hapla } 5888fa9e22eSVaclav Hapla if (on) *on = len; 5898fa9e22eSVaclav Hapla if (oidxs) *oidxs = lidxs; 5908fa9e22eSVaclav Hapla if (ogidxs) *ogidxs = work; 5913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5928fa9e22eSVaclav Hapla } 593deffd5ebSksagiyam 594deffd5ebSksagiyam /*@ 595cab54364SBarry Smith PetscSFCreateByMatchingIndices - Create `PetscSF` by matching root and leaf indices 596deffd5ebSksagiyam 597deffd5ebSksagiyam Collective 598deffd5ebSksagiyam 599deffd5ebSksagiyam Input Parameters: 600cab54364SBarry Smith + layout - `PetscLayout` defining the global index space and the rank that brokers each index 601deffd5ebSksagiyam . numRootIndices - size of rootIndices 602cab54364SBarry Smith . rootIndices - `PetscInt` array of global indices of which this process requests ownership 603deffd5ebSksagiyam . rootLocalIndices - root local index permutation (NULL if no permutation) 604deffd5ebSksagiyam . rootLocalOffset - offset to be added to root local indices 605deffd5ebSksagiyam . numLeafIndices - size of leafIndices 606cab54364SBarry Smith . leafIndices - `PetscInt` array of global indices with which this process requires data associated 607deffd5ebSksagiyam . leafLocalIndices - leaf local index permutation (NULL if no permutation) 608deffd5ebSksagiyam - leafLocalOffset - offset to be added to leaf local indices 609deffd5ebSksagiyam 610d8d19677SJose E. Roman Output Parameters: 611deffd5ebSksagiyam + sfA - star forest representing the communication pattern from the layout space to the leaf space (NULL if not needed) 612deffd5ebSksagiyam - sf - star forest representing the communication pattern from the root space to the leaf space 613deffd5ebSksagiyam 614cab54364SBarry Smith Level: advanced 615cab54364SBarry Smith 616deffd5ebSksagiyam Example 1: 617cab54364SBarry Smith .vb 618cab54364SBarry Smith rank : 0 1 2 619cab54364SBarry Smith rootIndices : [1 0 2] [3] [3] 620cab54364SBarry Smith rootLocalOffset : 100 200 300 621cab54364SBarry Smith layout : [0 1] [2] [3] 622cab54364SBarry Smith leafIndices : [0] [2] [0 3] 623cab54364SBarry Smith leafLocalOffset : 400 500 600 624cab54364SBarry Smith 625cab54364SBarry Smith would build the following PetscSF 626cab54364SBarry Smith 627cab54364SBarry Smith [0] 400 <- (0,101) 628cab54364SBarry Smith [1] 500 <- (0,102) 629cab54364SBarry Smith [2] 600 <- (0,101) 630cab54364SBarry Smith [2] 601 <- (2,300) 631cab54364SBarry Smith .ve 632cab54364SBarry Smith 633deffd5ebSksagiyam Example 2: 634cab54364SBarry Smith .vb 635cab54364SBarry Smith rank : 0 1 2 636cab54364SBarry Smith rootIndices : [1 0 2] [3] [3] 637cab54364SBarry Smith rootLocalOffset : 100 200 300 638cab54364SBarry Smith layout : [0 1] [2] [3] 639cab54364SBarry Smith leafIndices : rootIndices rootIndices rootIndices 640cab54364SBarry Smith leafLocalOffset : rootLocalOffset rootLocalOffset rootLocalOffset 641cab54364SBarry Smith 642cab54364SBarry Smith would build the following PetscSF 643cab54364SBarry Smith 644cab54364SBarry Smith [1] 200 <- (2,300) 645cab54364SBarry Smith .ve 646cab54364SBarry Smith 647deffd5ebSksagiyam Example 3: 648cab54364SBarry Smith .vb 649cab54364SBarry Smith No process requests ownership of global index 1, but no process needs it. 650cab54364SBarry Smith 651cab54364SBarry Smith rank : 0 1 2 652cab54364SBarry Smith numRootIndices : 2 1 1 653cab54364SBarry Smith rootIndices : [0 2] [3] [3] 654cab54364SBarry Smith rootLocalOffset : 100 200 300 655cab54364SBarry Smith layout : [0 1] [2] [3] 656cab54364SBarry Smith numLeafIndices : 1 1 2 657cab54364SBarry Smith leafIndices : [0] [2] [0 3] 658cab54364SBarry Smith leafLocalOffset : 400 500 600 659cab54364SBarry Smith 660cab54364SBarry Smith would build the following PetscSF 661cab54364SBarry Smith 662cab54364SBarry Smith [0] 400 <- (0,100) 663cab54364SBarry Smith [1] 500 <- (0,101) 664cab54364SBarry Smith [2] 600 <- (0,100) 665cab54364SBarry Smith [2] 601 <- (2,300) 666cab54364SBarry Smith .ve 667deffd5ebSksagiyam 668deffd5ebSksagiyam Notes: 669deffd5ebSksagiyam The layout parameter represents any partitioning of [0, N), where N is the total number of global indices, and its 670cab54364SBarry Smith local size can be set to `PETSC_DECIDE`. 671cab54364SBarry Smith 672deffd5ebSksagiyam If a global index x lies in the partition owned by process i, each process whose rootIndices contains x requests 673deffd5ebSksagiyam ownership of x and sends its own rank and the local index of x to process i. 674deffd5ebSksagiyam If multiple processes request ownership of x, the one with the highest rank is to own x. 675deffd5ebSksagiyam Process i then broadcasts the ownership information, so that each process whose leafIndices contains x knows the 676deffd5ebSksagiyam ownership information of x. 677deffd5ebSksagiyam The output sf is constructed by associating each leaf point to a root point in this way. 678deffd5ebSksagiyam 679deffd5ebSksagiyam Suppose there is point data ordered according to the global indices and partitioned according to the given layout. 680cab54364SBarry Smith The optional output `PetscSF` object sfA can be used to push such data to leaf points. 681deffd5ebSksagiyam 682deffd5ebSksagiyam All indices in rootIndices and leafIndices must lie in the layout range. The union (over all processes) of rootIndices 683deffd5ebSksagiyam must cover that of leafIndices, but need not cover the entire layout. 684deffd5ebSksagiyam 685deffd5ebSksagiyam If (leafIndices, leafLocalIndices, leafLocalOffset) == (rootIndices, rootLocalIndices, rootLocalOffset), the output 686deffd5ebSksagiyam star forest is almost identity, so will only include non-trivial part of the map. 687deffd5ebSksagiyam 68838b5cf2dSJacob Faibussowitsch Developer Notes: 689deffd5ebSksagiyam Current approach of a process of the highest rank gaining the ownership may cause load imbalance; consider using 690deffd5ebSksagiyam hash(rank, root_local_index) as the bid for the ownership determination. 691deffd5ebSksagiyam 692cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()` 693deffd5ebSksagiyam @*/ 694d71ae5a4SJacob 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) 695d71ae5a4SJacob Faibussowitsch { 696deffd5ebSksagiyam MPI_Comm comm = layout->comm; 697deffd5ebSksagiyam PetscMPIInt size, rank; 698deffd5ebSksagiyam PetscSF sf1; 699deffd5ebSksagiyam PetscSFNode *owners, *buffer, *iremote; 700deffd5ebSksagiyam PetscInt *ilocal, nleaves, N, n, i; 701deffd5ebSksagiyam #if defined(PETSC_USE_DEBUG) 702deffd5ebSksagiyam PetscInt N1; 703deffd5ebSksagiyam #endif 704deffd5ebSksagiyam PetscBool flag; 705deffd5ebSksagiyam 706deffd5ebSksagiyam PetscFunctionBegin; 7074f572ea9SToby Isaac if (rootIndices) PetscAssertPointer(rootIndices, 3); 7084f572ea9SToby Isaac if (rootLocalIndices) PetscAssertPointer(rootLocalIndices, 4); 7094f572ea9SToby Isaac if (leafIndices) PetscAssertPointer(leafIndices, 7); 7104f572ea9SToby Isaac if (leafLocalIndices) PetscAssertPointer(leafLocalIndices, 8); 7114f572ea9SToby Isaac if (sfA) PetscAssertPointer(sfA, 10); 7124f572ea9SToby Isaac PetscAssertPointer(sf, 11); 71308401ef6SPierre Jolivet PetscCheck(numRootIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numRootIndices (%" PetscInt_FMT ") must be non-negative", numRootIndices); 71408401ef6SPierre Jolivet PetscCheck(numLeafIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numLeafIndices (%" PetscInt_FMT ") must be non-negative", numLeafIndices); 7159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 7169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 7179566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(layout)); 7189566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(layout, &N)); 7199566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(layout, &n)); 720deffd5ebSksagiyam flag = (PetscBool)(leafIndices == rootIndices); 721462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &flag, 1, MPIU_BOOL, MPI_LAND, comm)); 722c9cc58a2SBarry Smith PetscCheck(!flag || numLeafIndices == numRootIndices, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "leafIndices == rootIndices, but numLeafIndices (%" PetscInt_FMT ") != numRootIndices(%" PetscInt_FMT ")", numLeafIndices, numRootIndices); 723deffd5ebSksagiyam #if defined(PETSC_USE_DEBUG) 7241690c2aeSBarry Smith N1 = PETSC_INT_MIN; 7259371c9d4SSatish Balay for (i = 0; i < numRootIndices; i++) 7269371c9d4SSatish Balay if (rootIndices[i] > N1) N1 = rootIndices[i]; 727462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm)); 72808401ef6SPierre 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); 729deffd5ebSksagiyam if (!flag) { 7301690c2aeSBarry Smith N1 = PETSC_INT_MIN; 7319371c9d4SSatish Balay for (i = 0; i < numLeafIndices; i++) 7329371c9d4SSatish Balay if (leafIndices[i] > N1) N1 = leafIndices[i]; 733462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm)); 73408401ef6SPierre 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); 735deffd5ebSksagiyam } 736deffd5ebSksagiyam #endif 737deffd5ebSksagiyam /* Reduce: owners -> buffer */ 7389566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &buffer)); 7399566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, &sf1)); 7409566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sf1)); 7419566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraphLayout(sf1, layout, numRootIndices, NULL, PETSC_OWN_POINTER, rootIndices)); 7429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numRootIndices, &owners)); 743deffd5ebSksagiyam for (i = 0; i < numRootIndices; ++i) { 744deffd5ebSksagiyam owners[i].rank = rank; 745deffd5ebSksagiyam owners[i].index = rootLocalOffset + (rootLocalIndices ? rootLocalIndices[i] : i); 746deffd5ebSksagiyam } 747deffd5ebSksagiyam for (i = 0; i < n; ++i) { 748deffd5ebSksagiyam buffer[i].index = -1; 749deffd5ebSksagiyam buffer[i].rank = -1; 750deffd5ebSksagiyam } 7516497c311SBarry Smith PetscCall(PetscSFReduceBegin(sf1, MPIU_SF_NODE, owners, buffer, MPI_MAXLOC)); 7526497c311SBarry Smith PetscCall(PetscSFReduceEnd(sf1, MPIU_SF_NODE, owners, buffer, MPI_MAXLOC)); 753deffd5ebSksagiyam /* Bcast: buffer -> owners */ 754deffd5ebSksagiyam if (!flag) { 755deffd5ebSksagiyam /* leafIndices is different from rootIndices */ 7569566063dSJacob Faibussowitsch PetscCall(PetscFree(owners)); 7579566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraphLayout(sf1, layout, numLeafIndices, NULL, PETSC_OWN_POINTER, leafIndices)); 7589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numLeafIndices, &owners)); 759deffd5ebSksagiyam } 7606497c311SBarry Smith PetscCall(PetscSFBcastBegin(sf1, MPIU_SF_NODE, buffer, owners, MPI_REPLACE)); 7616497c311SBarry Smith PetscCall(PetscSFBcastEnd(sf1, MPIU_SF_NODE, buffer, owners, MPI_REPLACE)); 76208401ef6SPierre 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]); 7639566063dSJacob Faibussowitsch PetscCall(PetscFree(buffer)); 7649371c9d4SSatish Balay if (sfA) { 7659371c9d4SSatish Balay *sfA = sf1; 7669371c9d4SSatish Balay } else PetscCall(PetscSFDestroy(&sf1)); 767deffd5ebSksagiyam /* Create sf */ 768deffd5ebSksagiyam if (flag && rootLocalIndices == leafLocalIndices && leafLocalOffset == rootLocalOffset) { 769deffd5ebSksagiyam /* leaf space == root space */ 7709371c9d4SSatish Balay for (i = 0, nleaves = 0; i < numLeafIndices; ++i) 7719371c9d4SSatish Balay if (owners[i].rank != rank) ++nleaves; 7729566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &ilocal)); 7739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &iremote)); 774deffd5ebSksagiyam for (i = 0, nleaves = 0; i < numLeafIndices; ++i) { 775deffd5ebSksagiyam if (owners[i].rank != rank) { 776deffd5ebSksagiyam ilocal[nleaves] = leafLocalOffset + i; 777deffd5ebSksagiyam iremote[nleaves].rank = owners[i].rank; 778deffd5ebSksagiyam iremote[nleaves].index = owners[i].index; 779deffd5ebSksagiyam ++nleaves; 780deffd5ebSksagiyam } 781deffd5ebSksagiyam } 7829566063dSJacob Faibussowitsch PetscCall(PetscFree(owners)); 783deffd5ebSksagiyam } else { 784deffd5ebSksagiyam nleaves = numLeafIndices; 7859566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &ilocal)); 786ad540459SPierre Jolivet for (i = 0; i < nleaves; ++i) ilocal[i] = leafLocalOffset + (leafLocalIndices ? leafLocalIndices[i] : i); 787deffd5ebSksagiyam iremote = owners; 788deffd5ebSksagiyam } 7899566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, sf)); 7909566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(*sf)); 7919566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*sf, rootLocalOffset + numRootIndices, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER)); 7923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 793deffd5ebSksagiyam } 794fbc7033fSJed Brown 795fbc7033fSJed Brown /*@ 79653c0d4aeSBarry Smith PetscSFMerge - append/merge indices of `sfb` into `sfa`, with preference for `sfb` 797fbc7033fSJed Brown 798fbc7033fSJed Brown Collective 799fbc7033fSJed Brown 80038b5cf2dSJacob Faibussowitsch Input Parameters: 801fbc7033fSJed Brown + sfa - default `PetscSF` 802fbc7033fSJed Brown - sfb - additional edges to add/replace edges in sfa 803fbc7033fSJed Brown 80438b5cf2dSJacob Faibussowitsch Output Parameter: 805fbc7033fSJed Brown . merged - new `PetscSF` with combined edges 806fbc7033fSJed Brown 80753c0d4aeSBarry Smith Level: intermediate 80853c0d4aeSBarry Smith 80938b5cf2dSJacob Faibussowitsch .seealso: `PetscSFCompose()` 810fbc7033fSJed Brown @*/ 811fbc7033fSJed Brown PetscErrorCode PetscSFMerge(PetscSF sfa, PetscSF sfb, PetscSF *merged) 812fbc7033fSJed Brown { 813fbc7033fSJed Brown PetscInt maxleaf; 814fbc7033fSJed Brown 815fbc7033fSJed Brown PetscFunctionBegin; 816fbc7033fSJed Brown PetscValidHeaderSpecific(sfa, PETSCSF_CLASSID, 1); 817fbc7033fSJed Brown PetscValidHeaderSpecific(sfb, PETSCSF_CLASSID, 2); 818fbc7033fSJed Brown PetscCheckSameComm(sfa, 1, sfb, 2); 8194f572ea9SToby Isaac PetscAssertPointer(merged, 3); 820fbc7033fSJed Brown { 821fbc7033fSJed Brown PetscInt aleaf, bleaf; 822fbc7033fSJed Brown PetscCall(PetscSFGetLeafRange(sfa, NULL, &aleaf)); 823fbc7033fSJed Brown PetscCall(PetscSFGetLeafRange(sfb, NULL, &bleaf)); 824fbc7033fSJed Brown maxleaf = PetscMax(aleaf, bleaf) + 1; // One more than the last index 825fbc7033fSJed Brown } 826fbc7033fSJed Brown PetscInt *clocal, aroots, aleaves, broots, bleaves; 827fbc7033fSJed Brown PetscSFNode *cremote; 828fbc7033fSJed Brown const PetscInt *alocal, *blocal; 829fbc7033fSJed Brown const PetscSFNode *aremote, *bremote; 830fbc7033fSJed Brown PetscCall(PetscMalloc2(maxleaf, &clocal, maxleaf, &cremote)); 831fbc7033fSJed Brown for (PetscInt i = 0; i < maxleaf; i++) clocal[i] = -1; 832fbc7033fSJed Brown PetscCall(PetscSFGetGraph(sfa, &aroots, &aleaves, &alocal, &aremote)); 833fbc7033fSJed Brown PetscCall(PetscSFGetGraph(sfb, &broots, &bleaves, &blocal, &bremote)); 834fbc7033fSJed Brown PetscCheck(aroots == broots, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Both sfa and sfb must have the same root space"); 835fbc7033fSJed Brown for (PetscInt i = 0; i < aleaves; i++) { 836fbc7033fSJed Brown PetscInt a = alocal ? alocal[i] : i; 837fbc7033fSJed Brown clocal[a] = a; 838fbc7033fSJed Brown cremote[a] = aremote[i]; 839fbc7033fSJed Brown } 840fbc7033fSJed Brown for (PetscInt i = 0; i < bleaves; i++) { 841fbc7033fSJed Brown PetscInt b = blocal ? blocal[i] : i; 842fbc7033fSJed Brown clocal[b] = b; 843fbc7033fSJed Brown cremote[b] = bremote[i]; 844fbc7033fSJed Brown } 845fbc7033fSJed Brown PetscInt nleaves = 0; 846fbc7033fSJed Brown for (PetscInt i = 0; i < maxleaf; i++) { 847fbc7033fSJed Brown if (clocal[i] < 0) continue; 848fbc7033fSJed Brown clocal[nleaves] = clocal[i]; 849fbc7033fSJed Brown cremote[nleaves] = cremote[i]; 850fbc7033fSJed Brown nleaves++; 851fbc7033fSJed Brown } 852fbc7033fSJed Brown PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)sfa), merged)); 853fbc7033fSJed Brown PetscCall(PetscSFSetGraph(*merged, aroots, nleaves, clocal, PETSC_COPY_VALUES, cremote, PETSC_COPY_VALUES)); 854fbc7033fSJed Brown PetscCall(PetscFree2(clocal, cremote)); 8553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 856fbc7033fSJed Brown } 8570dd791a8SStefano Zampini 8580dd791a8SStefano Zampini /*@ 8590dd791a8SStefano Zampini PetscSFCreateStridedSF - Create an `PetscSF` to communicate interleaved blocks of data 8600dd791a8SStefano Zampini 8610dd791a8SStefano Zampini Collective 8620dd791a8SStefano Zampini 8630dd791a8SStefano Zampini Input Parameters: 8640dd791a8SStefano Zampini + sf - star forest 8650dd791a8SStefano Zampini . bs - stride 8660dd791a8SStefano Zampini . ldr - leading dimension of root space 8670dd791a8SStefano Zampini - ldl - leading dimension of leaf space 8680dd791a8SStefano Zampini 8690dd791a8SStefano Zampini Output Parameter: 8700dd791a8SStefano Zampini . vsf - the new `PetscSF` 8710dd791a8SStefano Zampini 8720dd791a8SStefano Zampini Level: intermediate 8730dd791a8SStefano Zampini 8740dd791a8SStefano Zampini Notes: 8750dd791a8SStefano Zampini This can be useful to perform communications on blocks of right-hand sides. For example, the calling sequence 8760dd791a8SStefano Zampini .vb 8770dd791a8SStefano Zampini c_datatype *roots, *leaves; 8780dd791a8SStefano Zampini for i in [0,bs) do 8790dd791a8SStefano Zampini PetscSFBcastBegin(sf, mpi_datatype, roots + i*ldr, leaves + i*ldl, op) 8800dd791a8SStefano Zampini PetscSFBcastEnd(sf, mpi_datatype, roots + i*ldr, leaves + i*ldl, op) 8810dd791a8SStefano Zampini .ve 8820dd791a8SStefano Zampini is equivalent to 8830dd791a8SStefano Zampini .vb 8840dd791a8SStefano Zampini c_datatype *roots, *leaves; 8850dd791a8SStefano Zampini PetscSFCreateStridedSF(sf, bs, ldr, ldl, &vsf) 8860dd791a8SStefano Zampini PetscSFBcastBegin(vsf, mpi_datatype, roots, leaves, op) 8870dd791a8SStefano Zampini PetscSFBcastEnd(vsf, mpi_datatype, roots, leaves, op) 8880dd791a8SStefano Zampini .ve 8890dd791a8SStefano Zampini 8900dd791a8SStefano Zampini Developer Notes: 8910dd791a8SStefano Zampini Should this functionality be handled with a new API instead of creating a new object? 8920dd791a8SStefano Zampini 8930dd791a8SStefano Zampini .seealso: `PetscSF`, `PetscSFCreate()`, `PetscSFSetGraph()` 8940dd791a8SStefano Zampini @*/ 8950dd791a8SStefano Zampini PetscErrorCode PetscSFCreateStridedSF(PetscSF sf, PetscInt bs, PetscInt ldr, PetscInt ldl, PetscSF *vsf) 8960dd791a8SStefano Zampini { 8970dd791a8SStefano Zampini PetscSF rankssf; 8980dd791a8SStefano Zampini const PetscSFNode *iremote, *sfrremote; 8990dd791a8SStefano Zampini PetscSFNode *viremote; 9000dd791a8SStefano Zampini const PetscInt *ilocal; 9010dd791a8SStefano Zampini PetscInt *vilocal = NULL, *ldrs; 9020dd791a8SStefano Zampini PetscInt nranks, nr, nl, vnr, vnl, maxl; 9030dd791a8SStefano Zampini PetscMPIInt rank; 9040dd791a8SStefano Zampini MPI_Comm comm; 9050dd791a8SStefano Zampini PetscSFType sftype; 9060dd791a8SStefano Zampini 9070dd791a8SStefano Zampini PetscFunctionBegin; 9080dd791a8SStefano Zampini PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 9090dd791a8SStefano Zampini PetscValidLogicalCollectiveInt(sf, bs, 2); 9100dd791a8SStefano Zampini PetscAssertPointer(vsf, 5); 9110dd791a8SStefano Zampini if (bs == 1) { 9120dd791a8SStefano Zampini PetscCall(PetscObjectReference((PetscObject)sf)); 9130dd791a8SStefano Zampini *vsf = sf; 9140dd791a8SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 9150dd791a8SStefano Zampini } 9160dd791a8SStefano Zampini PetscCall(PetscSFSetUp(sf)); 9170dd791a8SStefano Zampini PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 9180dd791a8SStefano Zampini PetscCallMPI(MPI_Comm_rank(comm, &rank)); 9190dd791a8SStefano Zampini PetscCall(PetscSFGetGraph(sf, &nr, &nl, &ilocal, &iremote)); 9200dd791a8SStefano Zampini PetscCall(PetscSFGetLeafRange(sf, NULL, &maxl)); 9210dd791a8SStefano Zampini maxl += 1; 9220dd791a8SStefano Zampini if (ldl == PETSC_DECIDE) ldl = maxl; 9230dd791a8SStefano Zampini if (ldr == PETSC_DECIDE) ldr = nr; 9240dd791a8SStefano Zampini PetscCheck(ldr >= nr, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid leading dimension %" PetscInt_FMT " must be smaller than number of roots %" PetscInt_FMT, ldr, nr); 9250dd791a8SStefano Zampini PetscCheck(ldl >= maxl, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid leading dimension %" PetscInt_FMT " must be larger than leaf range %" PetscInt_FMT, ldl, maxl - 1); 9260dd791a8SStefano Zampini vnr = nr * bs; 9270dd791a8SStefano Zampini vnl = nl * bs; 9280dd791a8SStefano Zampini PetscCall(PetscMalloc1(vnl, &viremote)); 9290dd791a8SStefano Zampini PetscCall(PetscMalloc1(vnl, &vilocal)); 9300dd791a8SStefano Zampini 9310dd791a8SStefano Zampini /* Communicate root leading dimensions to leaf ranks */ 9320dd791a8SStefano Zampini PetscCall(PetscSFGetRanksSF(sf, &rankssf)); 9330dd791a8SStefano Zampini PetscCall(PetscSFGetGraph(rankssf, NULL, &nranks, NULL, &sfrremote)); 9340dd791a8SStefano Zampini PetscCall(PetscMalloc1(nranks, &ldrs)); 9350dd791a8SStefano Zampini PetscCall(PetscSFBcastBegin(rankssf, MPIU_INT, &ldr, ldrs, MPI_REPLACE)); 9360dd791a8SStefano Zampini PetscCall(PetscSFBcastEnd(rankssf, MPIU_INT, &ldr, ldrs, MPI_REPLACE)); 9370dd791a8SStefano Zampini 9380dd791a8SStefano Zampini for (PetscInt i = 0, rold = -1, lda = -1; i < nl; i++) { 939*835f2295SStefano Zampini const PetscInt r = iremote[i].rank; 9400dd791a8SStefano Zampini const PetscInt ii = iremote[i].index; 9410dd791a8SStefano Zampini 9420dd791a8SStefano Zampini if (r == rank) lda = ldr; 9430dd791a8SStefano Zampini else if (rold != r) { 9440dd791a8SStefano Zampini PetscInt j; 9450dd791a8SStefano Zampini 9460dd791a8SStefano Zampini for (j = 0; j < nranks; j++) 9470dd791a8SStefano Zampini if (sfrremote[j].rank == r) break; 948*835f2295SStefano Zampini PetscCheck(j < nranks, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unable to locate neighbor rank %" PetscInt_FMT, r); 9490dd791a8SStefano Zampini lda = ldrs[j]; 9500dd791a8SStefano Zampini } 9510dd791a8SStefano Zampini rold = r; 9520dd791a8SStefano Zampini for (PetscInt v = 0; v < bs; v++) { 9530dd791a8SStefano Zampini viremote[v * nl + i].rank = r; 9540dd791a8SStefano Zampini viremote[v * nl + i].index = v * lda + ii; 9550dd791a8SStefano Zampini vilocal[v * nl + i] = v * ldl + (ilocal ? ilocal[i] : i); 9560dd791a8SStefano Zampini } 9570dd791a8SStefano Zampini } 9580dd791a8SStefano Zampini PetscCall(PetscFree(ldrs)); 9590dd791a8SStefano Zampini PetscCall(PetscSFCreate(comm, vsf)); 9600dd791a8SStefano Zampini PetscCall(PetscSFGetType(sf, &sftype)); 9610dd791a8SStefano Zampini PetscCall(PetscSFSetType(*vsf, sftype)); 9620dd791a8SStefano Zampini PetscCall(PetscSFSetGraph(*vsf, vnr, vnl, vilocal, PETSC_OWN_POINTER, viremote, PETSC_OWN_POINTER)); 9630dd791a8SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 9640dd791a8SStefano Zampini } 965