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 23cab54364SBarry Smith Developer Note: 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; 38b114984aSVaclav Hapla PetscCall(PetscLayoutSetUp(layout)); 399566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(layout, &nroots)); 409566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRanges(layout, &range)); 419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &remote)); 42eb58282aSVaclav Hapla if (nleaves) ls = gremote[0] + 1; 43b0c7db22SLisandro Dalcin for (i = 0; i < nleaves; i++) { 44eb58282aSVaclav Hapla const PetscInt idx = gremote[i] - ls; 4538a25198SStefano Zampini if (idx < 0 || idx >= ln) { /* short-circuit the search */ 46eb58282aSVaclav Hapla PetscCall(PetscLayoutFindOwnerIndex(layout, gremote[i], &lr, &remote[i].index)); 4738a25198SStefano Zampini remote[i].rank = lr; 4838a25198SStefano Zampini ls = range[lr]; 4938a25198SStefano Zampini ln = range[lr + 1] - ls; 5038a25198SStefano Zampini } else { 5138a25198SStefano Zampini remote[i].rank = lr; 5238a25198SStefano Zampini remote[i].index = idx; 5338a25198SStefano Zampini } 54b0c7db22SLisandro Dalcin } 559566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, nroots, nleaves, ilocal, localmode, remote, PETSC_OWN_POINTER)); 56b0c7db22SLisandro Dalcin PetscFunctionReturn(0); 57b0c7db22SLisandro Dalcin } 58b0c7db22SLisandro Dalcin 59eb58282aSVaclav Hapla /*@C 60cab54364SBarry Smith PetscSFGetGraphLayout - Get the global indices and `PetscLayout` that describe this star forest 61eb58282aSVaclav Hapla 62eb58282aSVaclav Hapla Collective 63eb58282aSVaclav Hapla 64eb58282aSVaclav Hapla Input Parameter: 65eb58282aSVaclav Hapla . sf - star forest 66eb58282aSVaclav Hapla 67eb58282aSVaclav Hapla Output Parameters: 68cab54364SBarry Smith + layout - `PetscLayout` defining the global space for roots 69eb58282aSVaclav Hapla . nleaves - number of leaf vertices on the current process, each of these references a root on any process 70eb58282aSVaclav Hapla . ilocal - locations of leaves in leafdata buffers, or NULL for contiguous storage 71eb58282aSVaclav Hapla - gremote - root vertices in global numbering corresponding to leaves in ilocal 72eb58282aSVaclav Hapla 73eb58282aSVaclav Hapla Level: intermediate 74eb58282aSVaclav Hapla 75eb58282aSVaclav Hapla Notes: 76eb58282aSVaclav Hapla The outputs are such that passing them as inputs to `PetscSFSetGraphLayout()` would lead to the same star forest. 77eb58282aSVaclav Hapla The outputs layout and gremote are freshly created each time this function is called, 78eb58282aSVaclav Hapla so they need to be freed by user and cannot be qualified as const. 79eb58282aSVaclav Hapla 80cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFSetGraphLayout()`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()` 81eb58282aSVaclav Hapla @*/ 82d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetGraphLayout(PetscSF sf, PetscLayout *layout, PetscInt *nleaves, const PetscInt *ilocal[], PetscInt *gremote[]) 83d71ae5a4SJacob Faibussowitsch { 84eb58282aSVaclav Hapla PetscInt nr, nl; 85eb58282aSVaclav Hapla const PetscSFNode *ir; 86eb58282aSVaclav Hapla PetscLayout lt; 87eb58282aSVaclav Hapla 88eb58282aSVaclav Hapla PetscFunctionBegin; 89eb58282aSVaclav Hapla PetscCall(PetscSFGetGraph(sf, &nr, &nl, ilocal, &ir)); 90eb58282aSVaclav Hapla PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)sf), nr, PETSC_DECIDE, 1, <)); 91eb58282aSVaclav Hapla if (gremote) { 92eb58282aSVaclav Hapla PetscInt i; 93eb58282aSVaclav Hapla const PetscInt *range; 94eb58282aSVaclav Hapla PetscInt *gr; 95eb58282aSVaclav Hapla 96eb58282aSVaclav Hapla PetscCall(PetscLayoutGetRanges(lt, &range)); 97eb58282aSVaclav Hapla PetscCall(PetscMalloc1(nl, &gr)); 98eb58282aSVaclav Hapla for (i = 0; i < nl; i++) gr[i] = range[ir[i].rank] + ir[i].index; 99eb58282aSVaclav Hapla *gremote = gr; 100eb58282aSVaclav Hapla } 101eb58282aSVaclav Hapla if (nleaves) *nleaves = nl; 102eb58282aSVaclav Hapla if (layout) *layout = lt; 103eb58282aSVaclav Hapla else PetscCall(PetscLayoutDestroy(<)); 104eb58282aSVaclav Hapla PetscFunctionReturn(0); 105eb58282aSVaclav Hapla } 106eb58282aSVaclav Hapla 107b0c7db22SLisandro Dalcin /*@ 108cab54364SBarry Smith PetscSFSetGraphSection - Sets the `PetscSF` graph encoding the parallel dof overlap based upon the `PetscSection` describing the data layout. 109b0c7db22SLisandro Dalcin 110b0c7db22SLisandro Dalcin Input Parameters: 111cab54364SBarry Smith + sf - The `PetscSF` 112cab54364SBarry Smith . localSection - `PetscSection` describing the local data layout 113cab54364SBarry Smith - globalSection - `PetscSection` describing the global data layout 114b0c7db22SLisandro Dalcin 115b0c7db22SLisandro Dalcin Level: developer 116b0c7db22SLisandro Dalcin 117cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFSetGraph()`, `PetscSFSetGraphLayout()` 118b0c7db22SLisandro Dalcin @*/ 119d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetGraphSection(PetscSF sf, PetscSection localSection, PetscSection globalSection) 120d71ae5a4SJacob Faibussowitsch { 121b0c7db22SLisandro Dalcin MPI_Comm comm; 122b0c7db22SLisandro Dalcin PetscLayout layout; 123b0c7db22SLisandro Dalcin const PetscInt *ranges; 124b0c7db22SLisandro Dalcin PetscInt *local; 125b0c7db22SLisandro Dalcin PetscSFNode *remote; 126b0c7db22SLisandro Dalcin PetscInt pStart, pEnd, p, nroots, nleaves = 0, l; 127b0c7db22SLisandro Dalcin PetscMPIInt size, rank; 128b0c7db22SLisandro Dalcin 129b0c7db22SLisandro Dalcin PetscFunctionBegin; 130b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 131b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(localSection, PETSC_SECTION_CLASSID, 2); 132b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(globalSection, PETSC_SECTION_CLASSID, 3); 133b0c7db22SLisandro Dalcin 1349566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 1359566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 1369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1379566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(globalSection, &pStart, &pEnd)); 1389566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstrainedStorageSize(globalSection, &nroots)); 1399566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(comm, &layout)); 1409566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(layout, 1)); 1419566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetLocalSize(layout, nroots)); 1429566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(layout)); 1439566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRanges(layout, &ranges)); 144b0c7db22SLisandro Dalcin for (p = pStart; p < pEnd; ++p) { 145b0c7db22SLisandro Dalcin PetscInt gdof, gcdof; 146b0c7db22SLisandro Dalcin 1479566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(globalSection, p, &gdof)); 1489566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof)); 149c9cc58a2SBarry 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)); 150b0c7db22SLisandro Dalcin nleaves += gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof; 151b0c7db22SLisandro Dalcin } 1529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &local)); 1539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &remote)); 154b0c7db22SLisandro Dalcin for (p = pStart, l = 0; p < pEnd; ++p) { 155b0c7db22SLisandro Dalcin const PetscInt *cind; 156b0c7db22SLisandro Dalcin PetscInt dof, cdof, off, gdof, gcdof, goff, gsize, d, c; 157b0c7db22SLisandro Dalcin 1589566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(localSection, p, &dof)); 1599566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(localSection, p, &off)); 1609566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(localSection, p, &cdof)); 1619566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintIndices(localSection, p, &cind)); 1629566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(globalSection, p, &gdof)); 1639566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof)); 1649566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(globalSection, p, &goff)); 165b0c7db22SLisandro Dalcin if (!gdof) continue; /* Censored point */ 166b0c7db22SLisandro Dalcin gsize = gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof; 167b0c7db22SLisandro Dalcin if (gsize != dof - cdof) { 16808401ef6SPierre 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); 169b0c7db22SLisandro Dalcin cdof = 0; /* Ignore constraints */ 170b0c7db22SLisandro Dalcin } 171b0c7db22SLisandro Dalcin for (d = 0, c = 0; d < dof; ++d) { 1729371c9d4SSatish Balay if ((c < cdof) && (cind[c] == d)) { 1739371c9d4SSatish Balay ++c; 1749371c9d4SSatish Balay continue; 1759371c9d4SSatish Balay } 176b0c7db22SLisandro Dalcin local[l + d - c] = off + d; 177b0c7db22SLisandro Dalcin } 17808401ef6SPierre 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); 179b0c7db22SLisandro Dalcin if (gdof < 0) { 180b0c7db22SLisandro Dalcin for (d = 0; d < gsize; ++d, ++l) { 181b0c7db22SLisandro Dalcin PetscInt offset = -(goff + 1) + d, r; 182b0c7db22SLisandro Dalcin 1839566063dSJacob Faibussowitsch PetscCall(PetscFindInt(offset, size + 1, ranges, &r)); 184b0c7db22SLisandro Dalcin if (r < 0) r = -(r + 2); 18508401ef6SPierre 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); 186b0c7db22SLisandro Dalcin remote[l].rank = r; 187b0c7db22SLisandro Dalcin remote[l].index = offset - ranges[r]; 188b0c7db22SLisandro Dalcin } 189b0c7db22SLisandro Dalcin } else { 190b0c7db22SLisandro Dalcin for (d = 0; d < gsize; ++d, ++l) { 191b0c7db22SLisandro Dalcin remote[l].rank = rank; 192b0c7db22SLisandro Dalcin remote[l].index = goff + d - ranges[rank]; 193b0c7db22SLisandro Dalcin } 194b0c7db22SLisandro Dalcin } 195b0c7db22SLisandro Dalcin } 19608401ef6SPierre Jolivet PetscCheck(l == nleaves, comm, PETSC_ERR_PLIB, "Iteration error, l %" PetscInt_FMT " != nleaves %" PetscInt_FMT, l, nleaves); 1979566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&layout)); 1989566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER)); 199b0c7db22SLisandro Dalcin PetscFunctionReturn(0); 200b0c7db22SLisandro Dalcin } 201b0c7db22SLisandro Dalcin 202b0c7db22SLisandro Dalcin /*@C 203cab54364SBarry Smith PetscSFDistributeSection - Create a new `PetscSection` reorganized, moving from the root to the leaves of the `PetscSF` 204b0c7db22SLisandro Dalcin 205c3339decSBarry Smith Collective 206b0c7db22SLisandro Dalcin 207b0c7db22SLisandro Dalcin Input Parameters: 208cab54364SBarry Smith + sf - The `PetscSF` 209b0c7db22SLisandro Dalcin - rootSection - Section defined on root space 210b0c7db22SLisandro Dalcin 211b0c7db22SLisandro Dalcin Output Parameters: 212b0c7db22SLisandro Dalcin + remoteOffsets - root offsets in leaf storage, or NULL 213b0c7db22SLisandro Dalcin - leafSection - Section defined on the leaf space 214b0c7db22SLisandro Dalcin 215b0c7db22SLisandro Dalcin Level: advanced 216b0c7db22SLisandro Dalcin 217*23e9620eSJunchao Zhang Fortran Notes: 218*23e9620eSJunchao Zhang In Fortran, use PetscSFDistributeSectionF90() 219*23e9620eSJunchao Zhang 220cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()` 221b0c7db22SLisandro Dalcin @*/ 222d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFDistributeSection(PetscSF sf, PetscSection rootSection, PetscInt **remoteOffsets, PetscSection leafSection) 223d71ae5a4SJacob Faibussowitsch { 224b0c7db22SLisandro Dalcin PetscSF embedSF; 225b0c7db22SLisandro Dalcin const PetscInt *indices; 226b0c7db22SLisandro Dalcin IS selected; 227b0c7db22SLisandro Dalcin PetscInt numFields, nroots, rpStart, rpEnd, lpStart = PETSC_MAX_INT, lpEnd = -1, f, c; 228b0c7db22SLisandro Dalcin PetscBool *sub, hasc; 229b0c7db22SLisandro Dalcin 230b0c7db22SLisandro Dalcin PetscFunctionBegin; 2319566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_DistSect, sf, 0, 0, 0)); 2329566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(rootSection, &numFields)); 233029e8818Sksagiyam if (numFields) { 234029e8818Sksagiyam IS perm; 235029e8818Sksagiyam 236029e8818Sksagiyam /* PetscSectionSetNumFields() calls PetscSectionReset(), which destroys 237029e8818Sksagiyam leafSection->perm. To keep this permutation set by the user, we grab 238029e8818Sksagiyam the reference before calling PetscSectionSetNumFields() and set it 239029e8818Sksagiyam back after. */ 2409566063dSJacob Faibussowitsch PetscCall(PetscSectionGetPermutation(leafSection, &perm)); 2419566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)perm)); 2429566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(leafSection, numFields)); 2439566063dSJacob Faibussowitsch PetscCall(PetscSectionSetPermutation(leafSection, perm)); 2449566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 245029e8818Sksagiyam } 2469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numFields + 2, &sub)); 247b0c7db22SLisandro Dalcin sub[1] = rootSection->bc ? PETSC_TRUE : PETSC_FALSE; 248b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) { 2492ee82556SMatthew G. Knepley PetscSectionSym sym, dsym = NULL; 250b0c7db22SLisandro Dalcin const char *name = NULL; 251b0c7db22SLisandro Dalcin PetscInt numComp = 0; 252b0c7db22SLisandro Dalcin 253b0c7db22SLisandro Dalcin sub[2 + f] = rootSection->field[f]->bc ? PETSC_TRUE : PETSC_FALSE; 2549566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldComponents(rootSection, f, &numComp)); 2559566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldName(rootSection, f, &name)); 2569566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldSym(rootSection, f, &sym)); 2579566063dSJacob Faibussowitsch if (sym) PetscCall(PetscSectionSymDistribute(sym, sf, &dsym)); 2589566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(leafSection, f, numComp)); 2599566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldName(leafSection, f, name)); 2609566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldSym(leafSection, f, dsym)); 2619566063dSJacob Faibussowitsch PetscCall(PetscSectionSymDestroy(&dsym)); 262b0c7db22SLisandro Dalcin for (c = 0; c < rootSection->numFieldComponents[f]; ++c) { 2639566063dSJacob Faibussowitsch PetscCall(PetscSectionGetComponentName(rootSection, f, c, &name)); 2649566063dSJacob Faibussowitsch PetscCall(PetscSectionSetComponentName(leafSection, f, c, name)); 265b0c7db22SLisandro Dalcin } 266b0c7db22SLisandro Dalcin } 2679566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd)); 2689566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL)); 269b0c7db22SLisandro Dalcin rpEnd = PetscMin(rpEnd, nroots); 270b0c7db22SLisandro Dalcin rpEnd = PetscMax(rpStart, rpEnd); 271b0c7db22SLisandro Dalcin /* see if we can avoid creating the embedded SF, since it can cost more than an allreduce */ 272b0c7db22SLisandro Dalcin sub[0] = (PetscBool)(nroots != rpEnd - rpStart); 2731c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(MPI_IN_PLACE, sub, 2 + numFields, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)sf))); 274b0c7db22SLisandro Dalcin if (sub[0]) { 2759566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected)); 2769566063dSJacob Faibussowitsch PetscCall(ISGetIndices(selected, &indices)); 2779566063dSJacob Faibussowitsch PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF)); 2789566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(selected, &indices)); 2799566063dSJacob Faibussowitsch PetscCall(ISDestroy(&selected)); 280b0c7db22SLisandro Dalcin } else { 2819566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)sf)); 282b0c7db22SLisandro Dalcin embedSF = sf; 283b0c7db22SLisandro Dalcin } 2849566063dSJacob Faibussowitsch PetscCall(PetscSFGetLeafRange(embedSF, &lpStart, &lpEnd)); 285b0c7db22SLisandro Dalcin lpEnd++; 286b0c7db22SLisandro Dalcin 2879566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(leafSection, lpStart, lpEnd)); 288b0c7db22SLisandro Dalcin 289b0c7db22SLisandro Dalcin /* Constrained dof section */ 290b0c7db22SLisandro Dalcin hasc = sub[1]; 291b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) hasc = (PetscBool)(hasc || sub[2 + f]); 292b0c7db22SLisandro Dalcin 293b0c7db22SLisandro Dalcin /* Could fuse these at the cost of copies and extra allocation */ 2949566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart], MPI_REPLACE)); 2959566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart], MPI_REPLACE)); 296b0c7db22SLisandro Dalcin if (sub[1]) { 2977c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(rootSection)); 2987c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(leafSection)); 2999566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE)); 3009566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE)); 301b0c7db22SLisandro Dalcin } 302b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) { 3039566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart], MPI_REPLACE)); 3049566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart], MPI_REPLACE)); 305b0c7db22SLisandro Dalcin if (sub[2 + f]) { 3067c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(rootSection->field[f])); 3077c0883d5SVaclav Hapla PetscCall(PetscSectionCheckConstraints_Private(leafSection->field[f])); 3089566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE)); 3099566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE)); 310b0c7db22SLisandro Dalcin } 311b0c7db22SLisandro Dalcin } 312b0c7db22SLisandro Dalcin if (remoteOffsets) { 3139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lpEnd - lpStart, remoteOffsets)); 3149566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE)); 3159566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE)); 316b0c7db22SLisandro Dalcin } 31769c11d05SVaclav Hapla PetscCall(PetscSectionInvalidateMaxDof_Internal(leafSection)); 3189566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(leafSection)); 319b0c7db22SLisandro Dalcin if (hasc) { /* need to communicate bcIndices */ 320b0c7db22SLisandro Dalcin PetscSF bcSF; 321b0c7db22SLisandro Dalcin PetscInt *rOffBc; 322b0c7db22SLisandro Dalcin 3239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lpEnd - lpStart, &rOffBc)); 324b0c7db22SLisandro Dalcin if (sub[1]) { 3259566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 3269566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 3279566063dSJacob Faibussowitsch PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->bc, rOffBc, leafSection->bc, &bcSF)); 3289566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE)); 3299566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE)); 3309566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&bcSF)); 331b0c7db22SLisandro Dalcin } 332b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) { 333b0c7db22SLisandro Dalcin if (sub[2 + f]) { 3349566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 3359566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE)); 3369566063dSJacob Faibussowitsch PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->field[f]->bc, rOffBc, leafSection->field[f]->bc, &bcSF)); 3379566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE)); 3389566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE)); 3399566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&bcSF)); 340b0c7db22SLisandro Dalcin } 341b0c7db22SLisandro Dalcin } 3429566063dSJacob Faibussowitsch PetscCall(PetscFree(rOffBc)); 343b0c7db22SLisandro Dalcin } 3449566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&embedSF)); 3459566063dSJacob Faibussowitsch PetscCall(PetscFree(sub)); 3469566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_DistSect, sf, 0, 0, 0)); 347b0c7db22SLisandro Dalcin PetscFunctionReturn(0); 348b0c7db22SLisandro Dalcin } 349b0c7db22SLisandro Dalcin 350b0c7db22SLisandro Dalcin /*@C 351b0c7db22SLisandro Dalcin PetscSFCreateRemoteOffsets - Create offsets for point data on remote processes 352b0c7db22SLisandro Dalcin 353c3339decSBarry Smith Collective 354b0c7db22SLisandro Dalcin 355b0c7db22SLisandro Dalcin Input Parameters: 356cab54364SBarry Smith + sf - The `PetscSF` 357b0c7db22SLisandro Dalcin . rootSection - Data layout of remote points for outgoing data (this is layout for SF roots) 358b0c7db22SLisandro Dalcin - leafSection - Data layout of local points for incoming data (this is layout for SF leaves) 359b0c7db22SLisandro Dalcin 360b0c7db22SLisandro Dalcin Output Parameter: 361b0c7db22SLisandro Dalcin . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL 362b0c7db22SLisandro Dalcin 363b0c7db22SLisandro Dalcin Level: developer 364b0c7db22SLisandro Dalcin 365*23e9620eSJunchao Zhang Fortran Notes: 366*23e9620eSJunchao Zhang In Fortran, use PetscSFCreateRemoteOffsetsF90() 367*23e9620eSJunchao Zhang 368cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()` 369b0c7db22SLisandro Dalcin @*/ 370d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateRemoteOffsets(PetscSF sf, PetscSection rootSection, PetscSection leafSection, PetscInt **remoteOffsets) 371d71ae5a4SJacob Faibussowitsch { 372b0c7db22SLisandro Dalcin PetscSF embedSF; 373b0c7db22SLisandro Dalcin const PetscInt *indices; 374b0c7db22SLisandro Dalcin IS selected; 375b0c7db22SLisandro Dalcin PetscInt numRoots, rpStart = 0, rpEnd = 0, lpStart = 0, lpEnd = 0; 376b0c7db22SLisandro Dalcin 377b0c7db22SLisandro Dalcin PetscFunctionBegin; 378b0c7db22SLisandro Dalcin *remoteOffsets = NULL; 3799566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &numRoots, NULL, NULL, NULL)); 380b0c7db22SLisandro Dalcin if (numRoots < 0) PetscFunctionReturn(0); 3819566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_RemoteOff, sf, 0, 0, 0)); 3829566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd)); 3839566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd)); 3849566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected)); 3859566063dSJacob Faibussowitsch PetscCall(ISGetIndices(selected, &indices)); 3869566063dSJacob Faibussowitsch PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF)); 3879566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(selected, &indices)); 3889566063dSJacob Faibussowitsch PetscCall(ISDestroy(&selected)); 3899566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(lpEnd - lpStart, remoteOffsets)); 3909566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE)); 3919566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE)); 3929566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&embedSF)); 3939566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_RemoteOff, sf, 0, 0, 0)); 394b0c7db22SLisandro Dalcin PetscFunctionReturn(0); 395b0c7db22SLisandro Dalcin } 396b0c7db22SLisandro Dalcin 397b0c7db22SLisandro Dalcin /*@C 398cab54364SBarry Smith PetscSFCreateSectionSF - Create an expanded `PetscSF` of dofs, assuming the input `PetscSF` relates points 399b0c7db22SLisandro Dalcin 400c3339decSBarry Smith Collective 401b0c7db22SLisandro Dalcin 402b0c7db22SLisandro Dalcin Input Parameters: 403cab54364SBarry Smith + sf - The `PetscSF` 404b0c7db22SLisandro Dalcin . rootSection - Data layout of remote points for outgoing data (this is usually the serial section) 405b0c7db22SLisandro Dalcin . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL 406b0c7db22SLisandro Dalcin - leafSection - Data layout of local points for incoming data (this is the distributed section) 407b0c7db22SLisandro Dalcin 408b0c7db22SLisandro Dalcin Output Parameters: 409cab54364SBarry Smith - sectionSF - The new `PetscSF` 410b0c7db22SLisandro Dalcin 411b0c7db22SLisandro Dalcin Level: advanced 412b0c7db22SLisandro Dalcin 413*23e9620eSJunchao Zhang Notes: 414cab54364SBarry Smith Either rootSection or remoteOffsets can be specified 415cab54364SBarry Smith 416*23e9620eSJunchao Zhang Fortran Notes: 417*23e9620eSJunchao Zhang In Fortran, use PetscSFCreateSectionSFF90() 418*23e9620eSJunchao Zhang 419cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()` 420b0c7db22SLisandro Dalcin @*/ 421d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateSectionSF(PetscSF sf, PetscSection rootSection, PetscInt remoteOffsets[], PetscSection leafSection, PetscSF *sectionSF) 422d71ae5a4SJacob Faibussowitsch { 423b0c7db22SLisandro Dalcin MPI_Comm comm; 424b0c7db22SLisandro Dalcin const PetscInt *localPoints; 425b0c7db22SLisandro Dalcin const PetscSFNode *remotePoints; 426b0c7db22SLisandro Dalcin PetscInt lpStart, lpEnd; 427b0c7db22SLisandro Dalcin PetscInt numRoots, numSectionRoots, numPoints, numIndices = 0; 428b0c7db22SLisandro Dalcin PetscInt *localIndices; 429b0c7db22SLisandro Dalcin PetscSFNode *remoteIndices; 430b0c7db22SLisandro Dalcin PetscInt i, ind; 431b0c7db22SLisandro Dalcin 432b0c7db22SLisandro Dalcin PetscFunctionBegin; 433b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 434b0c7db22SLisandro Dalcin PetscValidPointer(rootSection, 2); 435b0c7db22SLisandro Dalcin /* Cannot check PetscValidIntPointer(remoteOffsets,3) because it can be NULL if sf does not reference any points in leafSection */ 436b0c7db22SLisandro Dalcin PetscValidPointer(leafSection, 4); 437b0c7db22SLisandro Dalcin PetscValidPointer(sectionSF, 5); 4389566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 4399566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, sectionSF)); 4409566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd)); 4419566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(rootSection, &numSectionRoots)); 4429566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &numRoots, &numPoints, &localPoints, &remotePoints)); 443b0c7db22SLisandro Dalcin if (numRoots < 0) PetscFunctionReturn(0); 4449566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_SectSF, sf, 0, 0, 0)); 445b0c7db22SLisandro Dalcin for (i = 0; i < numPoints; ++i) { 446b0c7db22SLisandro Dalcin PetscInt localPoint = localPoints ? localPoints[i] : i; 447b0c7db22SLisandro Dalcin PetscInt dof; 448b0c7db22SLisandro Dalcin 449b0c7db22SLisandro Dalcin if ((localPoint >= lpStart) && (localPoint < lpEnd)) { 4509566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof)); 45122a18c9fSMatthew G. Knepley numIndices += dof < 0 ? 0 : dof; 452b0c7db22SLisandro Dalcin } 453b0c7db22SLisandro Dalcin } 4549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numIndices, &localIndices)); 4559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numIndices, &remoteIndices)); 456b0c7db22SLisandro Dalcin /* Create new index graph */ 457b0c7db22SLisandro Dalcin for (i = 0, ind = 0; i < numPoints; ++i) { 458b0c7db22SLisandro Dalcin PetscInt localPoint = localPoints ? localPoints[i] : i; 459b0c7db22SLisandro Dalcin PetscInt rank = remotePoints[i].rank; 460b0c7db22SLisandro Dalcin 461b0c7db22SLisandro Dalcin if ((localPoint >= lpStart) && (localPoint < lpEnd)) { 462b0c7db22SLisandro Dalcin PetscInt remoteOffset = remoteOffsets[localPoint - lpStart]; 463b0c7db22SLisandro Dalcin PetscInt loff, dof, d; 464b0c7db22SLisandro Dalcin 4659566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(leafSection, localPoint, &loff)); 4669566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof)); 467b0c7db22SLisandro Dalcin for (d = 0; d < dof; ++d, ++ind) { 468b0c7db22SLisandro Dalcin localIndices[ind] = loff + d; 469b0c7db22SLisandro Dalcin remoteIndices[ind].rank = rank; 470b0c7db22SLisandro Dalcin remoteIndices[ind].index = remoteOffset + d; 471b0c7db22SLisandro Dalcin } 472b0c7db22SLisandro Dalcin } 473b0c7db22SLisandro Dalcin } 47408401ef6SPierre Jolivet PetscCheck(numIndices == ind, comm, PETSC_ERR_PLIB, "Inconsistency in indices, %" PetscInt_FMT " should be %" PetscInt_FMT, ind, numIndices); 4759566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*sectionSF, numSectionRoots, numIndices, localIndices, PETSC_OWN_POINTER, remoteIndices, PETSC_OWN_POINTER)); 4769566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(*sectionSF)); 4779566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_SectSF, sf, 0, 0, 0)); 478b0c7db22SLisandro Dalcin PetscFunctionReturn(0); 479b0c7db22SLisandro Dalcin } 4808fa9e22eSVaclav Hapla 4818fa9e22eSVaclav Hapla /*@C 482cab54364SBarry Smith PetscSFCreateFromLayouts - Creates a parallel star forest mapping two `PetscLayout` objects 4838fa9e22eSVaclav Hapla 4848fa9e22eSVaclav Hapla Collective 4858fa9e22eSVaclav Hapla 4864165533cSJose E. Roman Input Parameters: 487cab54364SBarry Smith + rmap - `PetscLayout` defining the global root space 488cab54364SBarry Smith - lmap - `PetscLayout` defining the global leaf space 4898fa9e22eSVaclav Hapla 4904165533cSJose E. Roman Output Parameter: 4918fa9e22eSVaclav Hapla . sf - The parallel star forest 4928fa9e22eSVaclav Hapla 4938fa9e22eSVaclav Hapla Level: intermediate 4948fa9e22eSVaclav Hapla 495cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()`, `PetscLayoutCreate()`, `PetscSFSetGraphLayout()` 4968fa9e22eSVaclav Hapla @*/ 497d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateFromLayouts(PetscLayout rmap, PetscLayout lmap, PetscSF *sf) 498d71ae5a4SJacob Faibussowitsch { 4998fa9e22eSVaclav Hapla PetscInt i, nroots, nleaves = 0; 5008fa9e22eSVaclav Hapla PetscInt rN, lst, len; 5018fa9e22eSVaclav Hapla PetscMPIInt owner = -1; 5028fa9e22eSVaclav Hapla PetscSFNode *remote; 5038fa9e22eSVaclav Hapla MPI_Comm rcomm = rmap->comm; 5048fa9e22eSVaclav Hapla MPI_Comm lcomm = lmap->comm; 5058fa9e22eSVaclav Hapla PetscMPIInt flg; 5068fa9e22eSVaclav Hapla 5078fa9e22eSVaclav Hapla PetscFunctionBegin; 5088fa9e22eSVaclav Hapla PetscValidPointer(sf, 3); 50928b400f6SJacob Faibussowitsch PetscCheck(rmap->setupcalled, rcomm, PETSC_ERR_ARG_WRONGSTATE, "Root layout not setup"); 51028b400f6SJacob Faibussowitsch PetscCheck(lmap->setupcalled, lcomm, PETSC_ERR_ARG_WRONGSTATE, "Leaf layout not setup"); 5119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(rcomm, lcomm, &flg)); 512c9cc58a2SBarry Smith PetscCheck(flg == MPI_CONGRUENT || flg == MPI_IDENT, rcomm, PETSC_ERR_SUP, "cannot map two layouts with non-matching communicators"); 5139566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(rcomm, sf)); 5149566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(rmap, &nroots)); 5159566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(rmap, &rN)); 5169566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRange(lmap, &lst, &len)); 5179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len - lst, &remote)); 5188fa9e22eSVaclav Hapla for (i = lst; i < len && i < rN; i++) { 51948a46eb9SPierre Jolivet if (owner < -1 || i >= rmap->range[owner + 1]) PetscCall(PetscLayoutFindOwner(rmap, i, &owner)); 5208fa9e22eSVaclav Hapla remote[nleaves].rank = owner; 5218fa9e22eSVaclav Hapla remote[nleaves].index = i - rmap->range[owner]; 5228fa9e22eSVaclav Hapla nleaves++; 5238fa9e22eSVaclav Hapla } 5249566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*sf, nroots, nleaves, NULL, PETSC_OWN_POINTER, remote, PETSC_COPY_VALUES)); 5259566063dSJacob Faibussowitsch PetscCall(PetscFree(remote)); 5268fa9e22eSVaclav Hapla PetscFunctionReturn(0); 5278fa9e22eSVaclav Hapla } 5288fa9e22eSVaclav Hapla 5298fa9e22eSVaclav Hapla /* TODO: handle nooffprocentries like MatZeroRowsMapLocal_Private, since this code is the same */ 530d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLayoutMapLocal(PetscLayout map, PetscInt N, const PetscInt idxs[], PetscInt *on, PetscInt **oidxs, PetscInt **ogidxs) 531d71ae5a4SJacob Faibussowitsch { 5328fa9e22eSVaclav Hapla PetscInt *owners = map->range; 5338fa9e22eSVaclav Hapla PetscInt n = map->n; 5348fa9e22eSVaclav Hapla PetscSF sf; 5358fa9e22eSVaclav Hapla PetscInt *lidxs, *work = NULL; 5368fa9e22eSVaclav Hapla PetscSFNode *ridxs; 5378fa9e22eSVaclav Hapla PetscMPIInt rank, p = 0; 5388fa9e22eSVaclav Hapla PetscInt r, len = 0; 5398fa9e22eSVaclav Hapla 5408fa9e22eSVaclav Hapla PetscFunctionBegin; 5418fa9e22eSVaclav Hapla if (on) *on = 0; /* squelch -Wmaybe-uninitialized */ 5428fa9e22eSVaclav Hapla /* Create SF where leaves are input idxs and roots are owned idxs */ 5439566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(map->comm, &rank)); 5449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &lidxs)); 5458fa9e22eSVaclav Hapla for (r = 0; r < n; ++r) lidxs[r] = -1; 5469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N, &ridxs)); 5478fa9e22eSVaclav Hapla for (r = 0; r < N; ++r) { 5488fa9e22eSVaclav Hapla const PetscInt idx = idxs[r]; 549c9cc58a2SBarry 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); 5508fa9e22eSVaclav Hapla if (idx < owners[p] || owners[p + 1] <= idx) { /* short-circuit the search if the last p owns this idx too */ 5519566063dSJacob Faibussowitsch PetscCall(PetscLayoutFindOwner(map, idx, &p)); 5528fa9e22eSVaclav Hapla } 5538fa9e22eSVaclav Hapla ridxs[r].rank = p; 5548fa9e22eSVaclav Hapla ridxs[r].index = idxs[r] - owners[p]; 5558fa9e22eSVaclav Hapla } 5569566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(map->comm, &sf)); 5579566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, n, N, NULL, PETSC_OWN_POINTER, ridxs, PETSC_OWN_POINTER)); 5589566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR)); 5599566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR)); 5608fa9e22eSVaclav Hapla if (ogidxs) { /* communicate global idxs */ 5618fa9e22eSVaclav Hapla PetscInt cum = 0, start, *work2; 5628fa9e22eSVaclav Hapla 5639566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &work)); 5649566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(N, &work2)); 5659371c9d4SSatish Balay for (r = 0; r < N; ++r) 5669371c9d4SSatish Balay if (idxs[r] >= 0) cum++; 5679566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&cum, &start, 1, MPIU_INT, MPI_SUM, map->comm)); 5688fa9e22eSVaclav Hapla start -= cum; 5698fa9e22eSVaclav Hapla cum = 0; 5709371c9d4SSatish Balay for (r = 0; r < N; ++r) 5719371c9d4SSatish Balay if (idxs[r] >= 0) work2[r] = start + cum++; 5729566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf, MPIU_INT, work2, work, MPI_REPLACE)); 5739566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf, MPIU_INT, work2, work, MPI_REPLACE)); 5749566063dSJacob Faibussowitsch PetscCall(PetscFree(work2)); 5758fa9e22eSVaclav Hapla } 5769566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf)); 5778fa9e22eSVaclav Hapla /* Compress and put in indices */ 5788fa9e22eSVaclav Hapla for (r = 0; r < n; ++r) 5798fa9e22eSVaclav Hapla if (lidxs[r] >= 0) { 5808fa9e22eSVaclav Hapla if (work) work[len] = work[r]; 5818fa9e22eSVaclav Hapla lidxs[len++] = r; 5828fa9e22eSVaclav Hapla } 5838fa9e22eSVaclav Hapla if (on) *on = len; 5848fa9e22eSVaclav Hapla if (oidxs) *oidxs = lidxs; 5858fa9e22eSVaclav Hapla if (ogidxs) *ogidxs = work; 5868fa9e22eSVaclav Hapla PetscFunctionReturn(0); 5878fa9e22eSVaclav Hapla } 588deffd5ebSksagiyam 589deffd5ebSksagiyam /*@ 590cab54364SBarry Smith PetscSFCreateByMatchingIndices - Create `PetscSF` by matching root and leaf indices 591deffd5ebSksagiyam 592deffd5ebSksagiyam Collective 593deffd5ebSksagiyam 594deffd5ebSksagiyam Input Parameters: 595cab54364SBarry Smith + layout - `PetscLayout` defining the global index space and the rank that brokers each index 596deffd5ebSksagiyam . numRootIndices - size of rootIndices 597cab54364SBarry Smith . rootIndices - `PetscInt` array of global indices of which this process requests ownership 598deffd5ebSksagiyam . rootLocalIndices - root local index permutation (NULL if no permutation) 599deffd5ebSksagiyam . rootLocalOffset - offset to be added to root local indices 600deffd5ebSksagiyam . numLeafIndices - size of leafIndices 601cab54364SBarry Smith . leafIndices - `PetscInt` array of global indices with which this process requires data associated 602deffd5ebSksagiyam . leafLocalIndices - leaf local index permutation (NULL if no permutation) 603deffd5ebSksagiyam - leafLocalOffset - offset to be added to leaf local indices 604deffd5ebSksagiyam 605d8d19677SJose E. Roman Output Parameters: 606deffd5ebSksagiyam + sfA - star forest representing the communication pattern from the layout space to the leaf space (NULL if not needed) 607deffd5ebSksagiyam - sf - star forest representing the communication pattern from the root space to the leaf space 608deffd5ebSksagiyam 609cab54364SBarry Smith Level: advanced 610cab54364SBarry Smith 611deffd5ebSksagiyam Example 1: 612cab54364SBarry Smith .vb 613cab54364SBarry Smith rank : 0 1 2 614cab54364SBarry Smith rootIndices : [1 0 2] [3] [3] 615cab54364SBarry Smith rootLocalOffset : 100 200 300 616cab54364SBarry Smith layout : [0 1] [2] [3] 617cab54364SBarry Smith leafIndices : [0] [2] [0 3] 618cab54364SBarry Smith leafLocalOffset : 400 500 600 619cab54364SBarry Smith 620cab54364SBarry Smith would build the following PetscSF 621cab54364SBarry Smith 622cab54364SBarry Smith [0] 400 <- (0,101) 623cab54364SBarry Smith [1] 500 <- (0,102) 624cab54364SBarry Smith [2] 600 <- (0,101) 625cab54364SBarry Smith [2] 601 <- (2,300) 626cab54364SBarry Smith .ve 627cab54364SBarry Smith 628deffd5ebSksagiyam Example 2: 629cab54364SBarry Smith .vb 630cab54364SBarry Smith rank : 0 1 2 631cab54364SBarry Smith rootIndices : [1 0 2] [3] [3] 632cab54364SBarry Smith rootLocalOffset : 100 200 300 633cab54364SBarry Smith layout : [0 1] [2] [3] 634cab54364SBarry Smith leafIndices : rootIndices rootIndices rootIndices 635cab54364SBarry Smith leafLocalOffset : rootLocalOffset rootLocalOffset rootLocalOffset 636cab54364SBarry Smith 637cab54364SBarry Smith would build the following PetscSF 638cab54364SBarry Smith 639cab54364SBarry Smith [1] 200 <- (2,300) 640cab54364SBarry Smith .ve 641cab54364SBarry Smith 642deffd5ebSksagiyam Example 3: 643cab54364SBarry Smith .vb 644cab54364SBarry Smith No process requests ownership of global index 1, but no process needs it. 645cab54364SBarry Smith 646cab54364SBarry Smith rank : 0 1 2 647cab54364SBarry Smith numRootIndices : 2 1 1 648cab54364SBarry Smith rootIndices : [0 2] [3] [3] 649cab54364SBarry Smith rootLocalOffset : 100 200 300 650cab54364SBarry Smith layout : [0 1] [2] [3] 651cab54364SBarry Smith numLeafIndices : 1 1 2 652cab54364SBarry Smith leafIndices : [0] [2] [0 3] 653cab54364SBarry Smith leafLocalOffset : 400 500 600 654cab54364SBarry Smith 655cab54364SBarry Smith would build the following PetscSF 656cab54364SBarry Smith 657cab54364SBarry Smith [0] 400 <- (0,100) 658cab54364SBarry Smith [1] 500 <- (0,101) 659cab54364SBarry Smith [2] 600 <- (0,100) 660cab54364SBarry Smith [2] 601 <- (2,300) 661cab54364SBarry Smith .ve 662deffd5ebSksagiyam 663deffd5ebSksagiyam Notes: 664deffd5ebSksagiyam The layout parameter represents any partitioning of [0, N), where N is the total number of global indices, and its 665cab54364SBarry Smith local size can be set to `PETSC_DECIDE`. 666cab54364SBarry Smith 667deffd5ebSksagiyam If a global index x lies in the partition owned by process i, each process whose rootIndices contains x requests 668deffd5ebSksagiyam ownership of x and sends its own rank and the local index of x to process i. 669deffd5ebSksagiyam If multiple processes request ownership of x, the one with the highest rank is to own x. 670deffd5ebSksagiyam Process i then broadcasts the ownership information, so that each process whose leafIndices contains x knows the 671deffd5ebSksagiyam ownership information of x. 672deffd5ebSksagiyam The output sf is constructed by associating each leaf point to a root point in this way. 673deffd5ebSksagiyam 674deffd5ebSksagiyam Suppose there is point data ordered according to the global indices and partitioned according to the given layout. 675cab54364SBarry Smith The optional output `PetscSF` object sfA can be used to push such data to leaf points. 676deffd5ebSksagiyam 677deffd5ebSksagiyam All indices in rootIndices and leafIndices must lie in the layout range. The union (over all processes) of rootIndices 678deffd5ebSksagiyam must cover that of leafIndices, but need not cover the entire layout. 679deffd5ebSksagiyam 680deffd5ebSksagiyam If (leafIndices, leafLocalIndices, leafLocalOffset) == (rootIndices, rootLocalIndices, rootLocalOffset), the output 681deffd5ebSksagiyam star forest is almost identity, so will only include non-trivial part of the map. 682deffd5ebSksagiyam 683cab54364SBarry Smith Developer Note: 684deffd5ebSksagiyam Current approach of a process of the highest rank gaining the ownership may cause load imbalance; consider using 685deffd5ebSksagiyam hash(rank, root_local_index) as the bid for the ownership determination. 686deffd5ebSksagiyam 687cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()` 688deffd5ebSksagiyam @*/ 689d71ae5a4SJacob 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) 690d71ae5a4SJacob Faibussowitsch { 691deffd5ebSksagiyam MPI_Comm comm = layout->comm; 692deffd5ebSksagiyam PetscMPIInt size, rank; 693deffd5ebSksagiyam PetscSF sf1; 694deffd5ebSksagiyam PetscSFNode *owners, *buffer, *iremote; 695deffd5ebSksagiyam PetscInt *ilocal, nleaves, N, n, i; 696deffd5ebSksagiyam #if defined(PETSC_USE_DEBUG) 697deffd5ebSksagiyam PetscInt N1; 698deffd5ebSksagiyam #endif 699deffd5ebSksagiyam PetscBool flag; 700deffd5ebSksagiyam 701deffd5ebSksagiyam PetscFunctionBegin; 702deffd5ebSksagiyam if (rootIndices) PetscValidIntPointer(rootIndices, 3); 703deffd5ebSksagiyam if (rootLocalIndices) PetscValidIntPointer(rootLocalIndices, 4); 704deffd5ebSksagiyam if (leafIndices) PetscValidIntPointer(leafIndices, 7); 705deffd5ebSksagiyam if (leafLocalIndices) PetscValidIntPointer(leafLocalIndices, 8); 706deffd5ebSksagiyam if (sfA) PetscValidPointer(sfA, 10); 707deffd5ebSksagiyam PetscValidPointer(sf, 11); 70808401ef6SPierre Jolivet PetscCheck(numRootIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numRootIndices (%" PetscInt_FMT ") must be non-negative", numRootIndices); 70908401ef6SPierre Jolivet PetscCheck(numLeafIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numLeafIndices (%" PetscInt_FMT ") must be non-negative", numLeafIndices); 7109566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 7119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 7129566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(layout)); 7139566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(layout, &N)); 7149566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(layout, &n)); 715deffd5ebSksagiyam flag = (PetscBool)(leafIndices == rootIndices); 7161c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &flag, 1, MPIU_BOOL, MPI_LAND, comm)); 717c9cc58a2SBarry Smith PetscCheck(!flag || numLeafIndices == numRootIndices, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "leafIndices == rootIndices, but numLeafIndices (%" PetscInt_FMT ") != numRootIndices(%" PetscInt_FMT ")", numLeafIndices, numRootIndices); 718deffd5ebSksagiyam #if defined(PETSC_USE_DEBUG) 719deffd5ebSksagiyam N1 = PETSC_MIN_INT; 7209371c9d4SSatish Balay for (i = 0; i < numRootIndices; i++) 7219371c9d4SSatish Balay if (rootIndices[i] > N1) N1 = rootIndices[i]; 7229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm)); 72308401ef6SPierre 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); 724deffd5ebSksagiyam if (!flag) { 725deffd5ebSksagiyam N1 = PETSC_MIN_INT; 7269371c9d4SSatish Balay for (i = 0; i < numLeafIndices; i++) 7279371c9d4SSatish Balay if (leafIndices[i] > N1) N1 = leafIndices[i]; 7289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm)); 72908401ef6SPierre 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); 730deffd5ebSksagiyam } 731deffd5ebSksagiyam #endif 732deffd5ebSksagiyam /* Reduce: owners -> buffer */ 7339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &buffer)); 7349566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, &sf1)); 7359566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sf1)); 7369566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraphLayout(sf1, layout, numRootIndices, NULL, PETSC_OWN_POINTER, rootIndices)); 7379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numRootIndices, &owners)); 738deffd5ebSksagiyam for (i = 0; i < numRootIndices; ++i) { 739deffd5ebSksagiyam owners[i].rank = rank; 740deffd5ebSksagiyam owners[i].index = rootLocalOffset + (rootLocalIndices ? rootLocalIndices[i] : i); 741deffd5ebSksagiyam } 742deffd5ebSksagiyam for (i = 0; i < n; ++i) { 743deffd5ebSksagiyam buffer[i].index = -1; 744deffd5ebSksagiyam buffer[i].rank = -1; 745deffd5ebSksagiyam } 7469566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC)); 7479566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC)); 748deffd5ebSksagiyam /* Bcast: buffer -> owners */ 749deffd5ebSksagiyam if (!flag) { 750deffd5ebSksagiyam /* leafIndices is different from rootIndices */ 7519566063dSJacob Faibussowitsch PetscCall(PetscFree(owners)); 7529566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraphLayout(sf1, layout, numLeafIndices, NULL, PETSC_OWN_POINTER, leafIndices)); 7539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numLeafIndices, &owners)); 754deffd5ebSksagiyam } 7559566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE)); 7569566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE)); 75708401ef6SPierre 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]); 7589566063dSJacob Faibussowitsch PetscCall(PetscFree(buffer)); 7599371c9d4SSatish Balay if (sfA) { 7609371c9d4SSatish Balay *sfA = sf1; 7619371c9d4SSatish Balay } else PetscCall(PetscSFDestroy(&sf1)); 762deffd5ebSksagiyam /* Create sf */ 763deffd5ebSksagiyam if (flag && rootLocalIndices == leafLocalIndices && leafLocalOffset == rootLocalOffset) { 764deffd5ebSksagiyam /* leaf space == root space */ 7659371c9d4SSatish Balay for (i = 0, nleaves = 0; i < numLeafIndices; ++i) 7669371c9d4SSatish Balay if (owners[i].rank != rank) ++nleaves; 7679566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &ilocal)); 7689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &iremote)); 769deffd5ebSksagiyam for (i = 0, nleaves = 0; i < numLeafIndices; ++i) { 770deffd5ebSksagiyam if (owners[i].rank != rank) { 771deffd5ebSksagiyam ilocal[nleaves] = leafLocalOffset + i; 772deffd5ebSksagiyam iremote[nleaves].rank = owners[i].rank; 773deffd5ebSksagiyam iremote[nleaves].index = owners[i].index; 774deffd5ebSksagiyam ++nleaves; 775deffd5ebSksagiyam } 776deffd5ebSksagiyam } 7779566063dSJacob Faibussowitsch PetscCall(PetscFree(owners)); 778deffd5ebSksagiyam } else { 779deffd5ebSksagiyam nleaves = numLeafIndices; 7809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &ilocal)); 781ad540459SPierre Jolivet for (i = 0; i < nleaves; ++i) ilocal[i] = leafLocalOffset + (leafLocalIndices ? leafLocalIndices[i] : i); 782deffd5ebSksagiyam iremote = owners; 783deffd5ebSksagiyam } 7849566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, sf)); 7859566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(*sf)); 7869566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*sf, rootLocalOffset + numRootIndices, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER)); 787deffd5ebSksagiyam PetscFunctionReturn(0); 788deffd5ebSksagiyam } 789fbc7033fSJed Brown 790fbc7033fSJed Brown /*@ 791fbc7033fSJed Brown PetscSFMerge - append/merge indices of sfb into sfa, with preference for sfb 792fbc7033fSJed Brown 793fbc7033fSJed Brown Collective 794fbc7033fSJed Brown 795fbc7033fSJed Brown Input Arguments: 796fbc7033fSJed Brown + sfa - default `PetscSF` 797fbc7033fSJed Brown - sfb - additional edges to add/replace edges in sfa 798fbc7033fSJed Brown 799fbc7033fSJed Brown Output Arguments: 800fbc7033fSJed Brown . merged - new `PetscSF` with combined edges 801fbc7033fSJed Brown 802fbc7033fSJed Brown .seealse: `PetscSFCompose()` 803fbc7033fSJed Brown @*/ 804fbc7033fSJed Brown PetscErrorCode PetscSFMerge(PetscSF sfa, PetscSF sfb, PetscSF *merged) 805fbc7033fSJed Brown { 806fbc7033fSJed Brown PetscInt maxleaf; 807fbc7033fSJed Brown 808fbc7033fSJed Brown PetscFunctionBegin; 809fbc7033fSJed Brown PetscValidHeaderSpecific(sfa, PETSCSF_CLASSID, 1); 810fbc7033fSJed Brown PetscValidHeaderSpecific(sfb, PETSCSF_CLASSID, 2); 811fbc7033fSJed Brown PetscCheckSameComm(sfa, 1, sfb, 2); 812fbc7033fSJed Brown PetscValidPointer(merged, 3); 813fbc7033fSJed Brown { 814fbc7033fSJed Brown PetscInt aleaf, bleaf; 815fbc7033fSJed Brown PetscCall(PetscSFGetLeafRange(sfa, NULL, &aleaf)); 816fbc7033fSJed Brown PetscCall(PetscSFGetLeafRange(sfb, NULL, &bleaf)); 817fbc7033fSJed Brown maxleaf = PetscMax(aleaf, bleaf) + 1; // One more than the last index 818fbc7033fSJed Brown } 819fbc7033fSJed Brown PetscInt *clocal, aroots, aleaves, broots, bleaves; 820fbc7033fSJed Brown PetscSFNode *cremote; 821fbc7033fSJed Brown const PetscInt *alocal, *blocal; 822fbc7033fSJed Brown const PetscSFNode *aremote, *bremote; 823fbc7033fSJed Brown PetscCall(PetscMalloc2(maxleaf, &clocal, maxleaf, &cremote)); 824fbc7033fSJed Brown for (PetscInt i = 0; i < maxleaf; i++) clocal[i] = -1; 825fbc7033fSJed Brown PetscCall(PetscSFGetGraph(sfa, &aroots, &aleaves, &alocal, &aremote)); 826fbc7033fSJed Brown PetscCall(PetscSFGetGraph(sfb, &broots, &bleaves, &blocal, &bremote)); 827fbc7033fSJed Brown PetscCheck(aroots == broots, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Both sfa and sfb must have the same root space"); 828fbc7033fSJed Brown for (PetscInt i = 0; i < aleaves; i++) { 829fbc7033fSJed Brown PetscInt a = alocal ? alocal[i] : i; 830fbc7033fSJed Brown clocal[a] = a; 831fbc7033fSJed Brown cremote[a] = aremote[i]; 832fbc7033fSJed Brown } 833fbc7033fSJed Brown for (PetscInt i = 0; i < bleaves; i++) { 834fbc7033fSJed Brown PetscInt b = blocal ? blocal[i] : i; 835fbc7033fSJed Brown clocal[b] = b; 836fbc7033fSJed Brown cremote[b] = bremote[i]; 837fbc7033fSJed Brown } 838fbc7033fSJed Brown PetscInt nleaves = 0; 839fbc7033fSJed Brown for (PetscInt i = 0; i < maxleaf; i++) { 840fbc7033fSJed Brown if (clocal[i] < 0) continue; 841fbc7033fSJed Brown clocal[nleaves] = clocal[i]; 842fbc7033fSJed Brown cremote[nleaves] = cremote[i]; 843fbc7033fSJed Brown nleaves++; 844fbc7033fSJed Brown } 845fbc7033fSJed Brown PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)sfa), merged)); 846fbc7033fSJed Brown PetscCall(PetscSFSetGraph(*merged, aroots, nleaves, clocal, PETSC_COPY_VALUES, cremote, PETSC_COPY_VALUES)); 847fbc7033fSJed Brown PetscCall(PetscFree2(clocal, cremote)); 848fbc7033fSJed Brown PetscFunctionReturn(0); 849fbc7033fSJed Brown } 850