1b0c7db22SLisandro Dalcin #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/ 2b0c7db22SLisandro Dalcin #include <petsc/private/sectionimpl.h> 3b0c7db22SLisandro Dalcin 4b0c7db22SLisandro Dalcin /*@C 5b0c7db22SLisandro Dalcin PetscSFSetGraphLayout - Set a parallel star forest via global indices and a PetscLayout 6b0c7db22SLisandro Dalcin 7b0c7db22SLisandro Dalcin Collective 8b0c7db22SLisandro Dalcin 94165533cSJose E. Roman Input Parameters: 10b0c7db22SLisandro Dalcin + sf - star forest 11ddea5d60SJunchao Zhang . layout - PetscLayout defining the global space for roots 12b0c7db22SLisandro Dalcin . nleaves - number of leaf vertices on the current process, each of these references a root on any process 13b0c7db22SLisandro Dalcin . ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage 14b0c7db22SLisandro Dalcin . localmode - copy mode for ilocal 1586c56f52SVaclav Hapla - iremote - root vertices in global numbering corresponding to leaves in ilocal 16b0c7db22SLisandro Dalcin 17b0c7db22SLisandro Dalcin Level: intermediate 18b0c7db22SLisandro Dalcin 1986c56f52SVaclav Hapla Notes: 2086c56f52SVaclav Hapla Global indices must lie in [0, N) where N is the global size of layout. 2186c56f52SVaclav Hapla 2286c56f52SVaclav Hapla Developers Note: 2386c56f52SVaclav Hapla Local indices which are the identity permutation in the range [0,nleaves) are discarded as they 24b0c7db22SLisandro Dalcin encode contiguous storage. In such case, if localmode is PETSC_OWN_POINTER, the memory is deallocated as it is not 25b0c7db22SLisandro Dalcin needed 26b0c7db22SLisandro Dalcin 27b0c7db22SLisandro Dalcin .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph(), PetscSFGetGraph() 28b0c7db22SLisandro Dalcin @*/ 29b0c7db22SLisandro Dalcin PetscErrorCode PetscSFSetGraphLayout(PetscSF sf,PetscLayout layout,PetscInt nleaves,const PetscInt *ilocal,PetscCopyMode localmode,const PetscInt *iremote) 30b0c7db22SLisandro Dalcin { 31b0c7db22SLisandro Dalcin PetscErrorCode ierr; 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; 38b0c7db22SLisandro Dalcin ierr = PetscLayoutGetLocalSize(layout,&nroots);CHKERRQ(ierr); 3938a25198SStefano Zampini ierr = PetscLayoutGetRanges(layout,&range);CHKERRQ(ierr); 40b0c7db22SLisandro Dalcin ierr = PetscMalloc1(nleaves,&remote);CHKERRQ(ierr); 4138a25198SStefano Zampini if (nleaves) { ls = iremote[0] + 1; } 42b0c7db22SLisandro Dalcin for (i=0; i<nleaves; i++) { 4338a25198SStefano Zampini const PetscInt idx = iremote[i] - ls; 4438a25198SStefano Zampini if (idx < 0 || idx >= ln) { /* short-circuit the search */ 4538a25198SStefano Zampini ierr = PetscLayoutFindOwnerIndex(layout,iremote[i],&lr,&remote[i].index);CHKERRQ(ierr); 4638a25198SStefano Zampini remote[i].rank = lr; 4738a25198SStefano Zampini ls = range[lr]; 4838a25198SStefano Zampini ln = range[lr+1] - ls; 4938a25198SStefano Zampini } else { 5038a25198SStefano Zampini remote[i].rank = lr; 5138a25198SStefano Zampini remote[i].index = idx; 5238a25198SStefano Zampini } 53b0c7db22SLisandro Dalcin } 54b0c7db22SLisandro Dalcin ierr = PetscSFSetGraph(sf,nroots,nleaves,ilocal,localmode,remote,PETSC_OWN_POINTER);CHKERRQ(ierr); 55b0c7db22SLisandro Dalcin PetscFunctionReturn(0); 56b0c7db22SLisandro Dalcin } 57b0c7db22SLisandro Dalcin 58b0c7db22SLisandro Dalcin /*@ 59b0c7db22SLisandro Dalcin PetscSFSetGraphSection - Sets the PetscSF graph encoding the parallel dof overlap based upon the PetscSections 60b0c7db22SLisandro Dalcin describing the data layout. 61b0c7db22SLisandro Dalcin 62b0c7db22SLisandro Dalcin Input Parameters: 63b0c7db22SLisandro Dalcin + sf - The SF 64b0c7db22SLisandro Dalcin . localSection - PetscSection describing the local data layout 65b0c7db22SLisandro Dalcin - globalSection - PetscSection describing the global data layout 66b0c7db22SLisandro Dalcin 67b0c7db22SLisandro Dalcin Level: developer 68b0c7db22SLisandro Dalcin 69b0c7db22SLisandro Dalcin .seealso: PetscSFSetGraph(), PetscSFSetGraphLayout() 70b0c7db22SLisandro Dalcin @*/ 71b0c7db22SLisandro Dalcin PetscErrorCode PetscSFSetGraphSection(PetscSF sf, PetscSection localSection, PetscSection globalSection) 72b0c7db22SLisandro Dalcin { 73b0c7db22SLisandro Dalcin MPI_Comm comm; 74b0c7db22SLisandro Dalcin PetscLayout layout; 75b0c7db22SLisandro Dalcin const PetscInt *ranges; 76b0c7db22SLisandro Dalcin PetscInt *local; 77b0c7db22SLisandro Dalcin PetscSFNode *remote; 78b0c7db22SLisandro Dalcin PetscInt pStart, pEnd, p, nroots, nleaves = 0, l; 79b0c7db22SLisandro Dalcin PetscMPIInt size, rank; 80b0c7db22SLisandro Dalcin PetscErrorCode ierr; 81b0c7db22SLisandro Dalcin 82b0c7db22SLisandro Dalcin PetscFunctionBegin; 83b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 84b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(localSection, PETSC_SECTION_CLASSID, 2); 85b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(globalSection, PETSC_SECTION_CLASSID, 3); 86b0c7db22SLisandro Dalcin 87b0c7db22SLisandro Dalcin ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr); 8855b25c41SPierre Jolivet ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr); 8955b25c41SPierre Jolivet ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 90b0c7db22SLisandro Dalcin ierr = PetscSectionGetChart(globalSection, &pStart, &pEnd);CHKERRQ(ierr); 91b0c7db22SLisandro Dalcin ierr = PetscSectionGetConstrainedStorageSize(globalSection, &nroots);CHKERRQ(ierr); 92b0c7db22SLisandro Dalcin ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr); 93b0c7db22SLisandro Dalcin ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr); 94b0c7db22SLisandro Dalcin ierr = PetscLayoutSetLocalSize(layout, nroots);CHKERRQ(ierr); 95b0c7db22SLisandro Dalcin ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr); 96b0c7db22SLisandro Dalcin ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr); 97b0c7db22SLisandro Dalcin for (p = pStart; p < pEnd; ++p) { 98b0c7db22SLisandro Dalcin PetscInt gdof, gcdof; 99b0c7db22SLisandro Dalcin 100b0c7db22SLisandro Dalcin ierr = PetscSectionGetDof(globalSection, p, &gdof);CHKERRQ(ierr); 101b0c7db22SLisandro Dalcin ierr = PetscSectionGetConstraintDof(globalSection, p, &gcdof);CHKERRQ(ierr); 1022c71b3e2SJacob Faibussowitsch PetscCheckFalse(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)); 103b0c7db22SLisandro Dalcin nleaves += gdof < 0 ? -(gdof+1)-gcdof : gdof-gcdof; 104b0c7db22SLisandro Dalcin } 105b0c7db22SLisandro Dalcin ierr = PetscMalloc1(nleaves, &local);CHKERRQ(ierr); 106b0c7db22SLisandro Dalcin ierr = PetscMalloc1(nleaves, &remote);CHKERRQ(ierr); 107b0c7db22SLisandro Dalcin for (p = pStart, l = 0; p < pEnd; ++p) { 108b0c7db22SLisandro Dalcin const PetscInt *cind; 109b0c7db22SLisandro Dalcin PetscInt dof, cdof, off, gdof, gcdof, goff, gsize, d, c; 110b0c7db22SLisandro Dalcin 111b0c7db22SLisandro Dalcin ierr = PetscSectionGetDof(localSection, p, &dof);CHKERRQ(ierr); 112b0c7db22SLisandro Dalcin ierr = PetscSectionGetOffset(localSection, p, &off);CHKERRQ(ierr); 113b0c7db22SLisandro Dalcin ierr = PetscSectionGetConstraintDof(localSection, p, &cdof);CHKERRQ(ierr); 114b0c7db22SLisandro Dalcin ierr = PetscSectionGetConstraintIndices(localSection, p, &cind);CHKERRQ(ierr); 115b0c7db22SLisandro Dalcin ierr = PetscSectionGetDof(globalSection, p, &gdof);CHKERRQ(ierr); 116b0c7db22SLisandro Dalcin ierr = PetscSectionGetConstraintDof(globalSection, p, &gcdof);CHKERRQ(ierr); 117b0c7db22SLisandro Dalcin ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr); 118b0c7db22SLisandro Dalcin if (!gdof) continue; /* Censored point */ 119b0c7db22SLisandro Dalcin gsize = gdof < 0 ? -(gdof+1)-gcdof : gdof-gcdof; 120b0c7db22SLisandro Dalcin if (gsize != dof-cdof) { 1212c71b3e2SJacob Faibussowitsch PetscCheckFalse(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); 122b0c7db22SLisandro Dalcin cdof = 0; /* Ignore constraints */ 123b0c7db22SLisandro Dalcin } 124b0c7db22SLisandro Dalcin for (d = 0, c = 0; d < dof; ++d) { 125b0c7db22SLisandro Dalcin if ((c < cdof) && (cind[c] == d)) {++c; continue;} 126b0c7db22SLisandro Dalcin local[l+d-c] = off+d; 127b0c7db22SLisandro Dalcin } 1282c71b3e2SJacob Faibussowitsch PetscCheckFalse(d-c != gsize,comm, PETSC_ERR_ARG_WRONG, "Point %" PetscInt_FMT ": Global dof %" PetscInt_FMT " != %" PetscInt_FMT " size - number of constraints", p, gsize, d-c); 129b0c7db22SLisandro Dalcin if (gdof < 0) { 130b0c7db22SLisandro Dalcin for (d = 0; d < gsize; ++d, ++l) { 131b0c7db22SLisandro Dalcin PetscInt offset = -(goff+1) + d, r; 132b0c7db22SLisandro Dalcin 133b0c7db22SLisandro Dalcin ierr = PetscFindInt(offset,size+1,ranges,&r);CHKERRQ(ierr); 134b0c7db22SLisandro Dalcin if (r < 0) r = -(r+2); 1352c71b3e2SJacob Faibussowitsch PetscCheckFalse((r < 0) || (r >= size),PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " mapped to invalid process %" PetscInt_FMT " (%" PetscInt_FMT ", %" PetscInt_FMT ")", p, r, gdof, goff); 136b0c7db22SLisandro Dalcin remote[l].rank = r; 137b0c7db22SLisandro Dalcin remote[l].index = offset - ranges[r]; 138b0c7db22SLisandro Dalcin } 139b0c7db22SLisandro Dalcin } else { 140b0c7db22SLisandro Dalcin for (d = 0; d < gsize; ++d, ++l) { 141b0c7db22SLisandro Dalcin remote[l].rank = rank; 142b0c7db22SLisandro Dalcin remote[l].index = goff+d - ranges[rank]; 143b0c7db22SLisandro Dalcin } 144b0c7db22SLisandro Dalcin } 145b0c7db22SLisandro Dalcin } 1462c71b3e2SJacob Faibussowitsch PetscCheckFalse(l != nleaves,comm, PETSC_ERR_PLIB, "Iteration error, l %" PetscInt_FMT " != nleaves %" PetscInt_FMT, l, nleaves); 147b0c7db22SLisandro Dalcin ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr); 148b0c7db22SLisandro Dalcin ierr = PetscSFSetGraph(sf, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER);CHKERRQ(ierr); 149b0c7db22SLisandro Dalcin PetscFunctionReturn(0); 150b0c7db22SLisandro Dalcin } 151b0c7db22SLisandro Dalcin 152b0c7db22SLisandro Dalcin static PetscErrorCode PetscSectionCheckConstraints_Static(PetscSection s) 153b0c7db22SLisandro Dalcin { 154b0c7db22SLisandro Dalcin PetscErrorCode ierr; 155b0c7db22SLisandro Dalcin 156b0c7db22SLisandro Dalcin PetscFunctionBegin; 157b0c7db22SLisandro Dalcin if (!s->bc) { 158b0c7db22SLisandro Dalcin ierr = PetscSectionCreate(PETSC_COMM_SELF, &s->bc);CHKERRQ(ierr); 159b0c7db22SLisandro Dalcin ierr = PetscSectionSetChart(s->bc, s->pStart, s->pEnd);CHKERRQ(ierr); 160b0c7db22SLisandro Dalcin } 161b0c7db22SLisandro Dalcin PetscFunctionReturn(0); 162b0c7db22SLisandro Dalcin } 163b0c7db22SLisandro Dalcin 164b0c7db22SLisandro Dalcin /*@C 165b0c7db22SLisandro Dalcin PetscSFDistributeSection - Create a new PetscSection reorganized, moving from the root to the leaves of the SF 166b0c7db22SLisandro Dalcin 167b0c7db22SLisandro Dalcin Collective on sf 168b0c7db22SLisandro Dalcin 169b0c7db22SLisandro Dalcin Input Parameters: 170b0c7db22SLisandro Dalcin + sf - The SF 171b0c7db22SLisandro Dalcin - rootSection - Section defined on root space 172b0c7db22SLisandro Dalcin 173b0c7db22SLisandro Dalcin Output Parameters: 174b0c7db22SLisandro Dalcin + remoteOffsets - root offsets in leaf storage, or NULL 175b0c7db22SLisandro Dalcin - leafSection - Section defined on the leaf space 176b0c7db22SLisandro Dalcin 177b0c7db22SLisandro Dalcin Level: advanced 178b0c7db22SLisandro Dalcin 179b0c7db22SLisandro Dalcin .seealso: PetscSFCreate() 180b0c7db22SLisandro Dalcin @*/ 181b0c7db22SLisandro Dalcin PetscErrorCode PetscSFDistributeSection(PetscSF sf, PetscSection rootSection, PetscInt **remoteOffsets, PetscSection leafSection) 182b0c7db22SLisandro Dalcin { 183b0c7db22SLisandro Dalcin PetscSF embedSF; 184b0c7db22SLisandro Dalcin const PetscInt *indices; 185b0c7db22SLisandro Dalcin IS selected; 186b0c7db22SLisandro Dalcin PetscInt numFields, nroots, rpStart, rpEnd, lpStart = PETSC_MAX_INT, lpEnd = -1, f, c; 187b0c7db22SLisandro Dalcin PetscBool *sub, hasc; 188b0c7db22SLisandro Dalcin PetscErrorCode ierr; 189b0c7db22SLisandro Dalcin 190b0c7db22SLisandro Dalcin PetscFunctionBegin; 191b0c7db22SLisandro Dalcin ierr = PetscLogEventBegin(PETSCSF_DistSect,sf,0,0,0);CHKERRQ(ierr); 192b0c7db22SLisandro Dalcin ierr = PetscSectionGetNumFields(rootSection, &numFields);CHKERRQ(ierr); 193029e8818Sksagiyam if (numFields) { 194029e8818Sksagiyam IS perm; 195029e8818Sksagiyam 196029e8818Sksagiyam /* PetscSectionSetNumFields() calls PetscSectionReset(), which destroys 197029e8818Sksagiyam leafSection->perm. To keep this permutation set by the user, we grab 198029e8818Sksagiyam the reference before calling PetscSectionSetNumFields() and set it 199029e8818Sksagiyam back after. */ 200029e8818Sksagiyam ierr = PetscSectionGetPermutation(leafSection, &perm);CHKERRQ(ierr); 201029e8818Sksagiyam ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 202029e8818Sksagiyam ierr = PetscSectionSetNumFields(leafSection, numFields);CHKERRQ(ierr); 203029e8818Sksagiyam ierr = PetscSectionSetPermutation(leafSection, perm);CHKERRQ(ierr); 204029e8818Sksagiyam ierr = ISDestroy(&perm);CHKERRQ(ierr); 205029e8818Sksagiyam } 206b0c7db22SLisandro Dalcin ierr = PetscMalloc1(numFields+2, &sub);CHKERRQ(ierr); 207b0c7db22SLisandro Dalcin sub[1] = rootSection->bc ? PETSC_TRUE : PETSC_FALSE; 208b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) { 209*2ee82556SMatthew G. Knepley PetscSectionSym sym, dsym = NULL; 210b0c7db22SLisandro Dalcin const char *name = NULL; 211b0c7db22SLisandro Dalcin PetscInt numComp = 0; 212b0c7db22SLisandro Dalcin 213b0c7db22SLisandro Dalcin sub[2 + f] = rootSection->field[f]->bc ? PETSC_TRUE : PETSC_FALSE; 214b0c7db22SLisandro Dalcin ierr = PetscSectionGetFieldComponents(rootSection, f, &numComp);CHKERRQ(ierr); 215b0c7db22SLisandro Dalcin ierr = PetscSectionGetFieldName(rootSection, f, &name);CHKERRQ(ierr); 216b0c7db22SLisandro Dalcin ierr = PetscSectionGetFieldSym(rootSection, f, &sym);CHKERRQ(ierr); 217*2ee82556SMatthew G. Knepley if (sym) {ierr = PetscSectionSymDistribute(sym, sf, &dsym);CHKERRQ(ierr);} 218b0c7db22SLisandro Dalcin ierr = PetscSectionSetFieldComponents(leafSection, f, numComp);CHKERRQ(ierr); 219b0c7db22SLisandro Dalcin ierr = PetscSectionSetFieldName(leafSection, f, name);CHKERRQ(ierr); 220*2ee82556SMatthew G. Knepley ierr = PetscSectionSetFieldSym(leafSection, f, dsym);CHKERRQ(ierr); 221*2ee82556SMatthew G. Knepley ierr = PetscSectionSymDestroy(&dsym);CHKERRQ(ierr); 222b0c7db22SLisandro Dalcin for (c = 0; c < rootSection->numFieldComponents[f]; ++c) { 223b0c7db22SLisandro Dalcin ierr = PetscSectionGetComponentName(rootSection, f, c, &name);CHKERRQ(ierr); 224b0c7db22SLisandro Dalcin ierr = PetscSectionSetComponentName(leafSection, f, c, name);CHKERRQ(ierr); 225b0c7db22SLisandro Dalcin } 226b0c7db22SLisandro Dalcin } 227b0c7db22SLisandro Dalcin ierr = PetscSectionGetChart(rootSection, &rpStart, &rpEnd);CHKERRQ(ierr); 228b0c7db22SLisandro Dalcin ierr = PetscSFGetGraph(sf,&nroots,NULL,NULL,NULL);CHKERRQ(ierr); 229b0c7db22SLisandro Dalcin rpEnd = PetscMin(rpEnd,nroots); 230b0c7db22SLisandro Dalcin rpEnd = PetscMax(rpStart,rpEnd); 231b0c7db22SLisandro Dalcin /* see if we can avoid creating the embedded SF, since it can cost more than an allreduce */ 232b0c7db22SLisandro Dalcin sub[0] = (PetscBool)(nroots != rpEnd - rpStart); 233820f2d46SBarry Smith ierr = MPIU_Allreduce(MPI_IN_PLACE, sub, 2+numFields, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)sf));CHKERRMPI(ierr); 234b0c7db22SLisandro Dalcin if (sub[0]) { 235b0c7db22SLisandro Dalcin ierr = ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected);CHKERRQ(ierr); 236b0c7db22SLisandro Dalcin ierr = ISGetIndices(selected, &indices);CHKERRQ(ierr); 23772502a1fSJunchao Zhang ierr = PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF);CHKERRQ(ierr); 238b0c7db22SLisandro Dalcin ierr = ISRestoreIndices(selected, &indices);CHKERRQ(ierr); 239b0c7db22SLisandro Dalcin ierr = ISDestroy(&selected);CHKERRQ(ierr); 240b0c7db22SLisandro Dalcin } else { 241b0c7db22SLisandro Dalcin ierr = PetscObjectReference((PetscObject)sf);CHKERRQ(ierr); 242b0c7db22SLisandro Dalcin embedSF = sf; 243b0c7db22SLisandro Dalcin } 244b0c7db22SLisandro Dalcin ierr = PetscSFGetLeafRange(embedSF, &lpStart, &lpEnd);CHKERRQ(ierr); 245b0c7db22SLisandro Dalcin lpEnd++; 246b0c7db22SLisandro Dalcin 247b0c7db22SLisandro Dalcin ierr = PetscSectionSetChart(leafSection, lpStart, lpEnd);CHKERRQ(ierr); 248b0c7db22SLisandro Dalcin 249b0c7db22SLisandro Dalcin /* Constrained dof section */ 250b0c7db22SLisandro Dalcin hasc = sub[1]; 251b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) hasc = (PetscBool)(hasc || sub[2+f]); 252b0c7db22SLisandro Dalcin 253b0c7db22SLisandro Dalcin /* Could fuse these at the cost of copies and extra allocation */ 254ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 255ad227feaSJunchao Zhang ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 256b0c7db22SLisandro Dalcin if (sub[1]) { 257b0c7db22SLisandro Dalcin ierr = PetscSectionCheckConstraints_Static(rootSection);CHKERRQ(ierr); 258b0c7db22SLisandro Dalcin ierr = PetscSectionCheckConstraints_Static(leafSection);CHKERRQ(ierr); 259ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 260ad227feaSJunchao Zhang ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 261b0c7db22SLisandro Dalcin } 262b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) { 263ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 264ad227feaSJunchao Zhang ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 265b0c7db22SLisandro Dalcin if (sub[2+f]) { 266b0c7db22SLisandro Dalcin ierr = PetscSectionCheckConstraints_Static(rootSection->field[f]);CHKERRQ(ierr); 267b0c7db22SLisandro Dalcin ierr = PetscSectionCheckConstraints_Static(leafSection->field[f]);CHKERRQ(ierr); 268ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 269ad227feaSJunchao Zhang ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 270b0c7db22SLisandro Dalcin } 271b0c7db22SLisandro Dalcin } 272b0c7db22SLisandro Dalcin if (remoteOffsets) { 273b0c7db22SLisandro Dalcin ierr = PetscMalloc1(lpEnd - lpStart, remoteOffsets);CHKERRQ(ierr); 274ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 275ad227feaSJunchao Zhang ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 276b0c7db22SLisandro Dalcin } 277b0c7db22SLisandro Dalcin ierr = PetscSectionSetUp(leafSection);CHKERRQ(ierr); 278b0c7db22SLisandro Dalcin if (hasc) { /* need to communicate bcIndices */ 279b0c7db22SLisandro Dalcin PetscSF bcSF; 280b0c7db22SLisandro Dalcin PetscInt *rOffBc; 281b0c7db22SLisandro Dalcin 282b0c7db22SLisandro Dalcin ierr = PetscMalloc1(lpEnd - lpStart, &rOffBc);CHKERRQ(ierr); 283b0c7db22SLisandro Dalcin if (sub[1]) { 284ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 285ad227feaSJunchao Zhang ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 286b0c7db22SLisandro Dalcin ierr = PetscSFCreateSectionSF(embedSF, rootSection->bc, rOffBc, leafSection->bc, &bcSF);CHKERRQ(ierr); 287ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices,MPI_REPLACE);CHKERRQ(ierr); 288ad227feaSJunchao Zhang ierr = PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices,MPI_REPLACE);CHKERRQ(ierr); 289b0c7db22SLisandro Dalcin ierr = PetscSFDestroy(&bcSF);CHKERRQ(ierr); 290b0c7db22SLisandro Dalcin } 291b0c7db22SLisandro Dalcin for (f = 0; f < numFields; ++f) { 292b0c7db22SLisandro Dalcin if (sub[2+f]) { 293ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 294ad227feaSJunchao Zhang ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 295b0c7db22SLisandro Dalcin ierr = PetscSFCreateSectionSF(embedSF, rootSection->field[f]->bc, rOffBc, leafSection->field[f]->bc, &bcSF);CHKERRQ(ierr); 296ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices,MPI_REPLACE);CHKERRQ(ierr); 297ad227feaSJunchao Zhang ierr = PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices,MPI_REPLACE);CHKERRQ(ierr); 298b0c7db22SLisandro Dalcin ierr = PetscSFDestroy(&bcSF);CHKERRQ(ierr); 299b0c7db22SLisandro Dalcin } 300b0c7db22SLisandro Dalcin } 301b0c7db22SLisandro Dalcin ierr = PetscFree(rOffBc);CHKERRQ(ierr); 302b0c7db22SLisandro Dalcin } 303b0c7db22SLisandro Dalcin ierr = PetscSFDestroy(&embedSF);CHKERRQ(ierr); 304b0c7db22SLisandro Dalcin ierr = PetscFree(sub);CHKERRQ(ierr); 305b0c7db22SLisandro Dalcin ierr = PetscLogEventEnd(PETSCSF_DistSect,sf,0,0,0);CHKERRQ(ierr); 306b0c7db22SLisandro Dalcin PetscFunctionReturn(0); 307b0c7db22SLisandro Dalcin } 308b0c7db22SLisandro Dalcin 309b0c7db22SLisandro Dalcin /*@C 310b0c7db22SLisandro Dalcin PetscSFCreateRemoteOffsets - Create offsets for point data on remote processes 311b0c7db22SLisandro Dalcin 312b0c7db22SLisandro Dalcin Collective on sf 313b0c7db22SLisandro Dalcin 314b0c7db22SLisandro Dalcin Input Parameters: 315b0c7db22SLisandro Dalcin + sf - The SF 316b0c7db22SLisandro Dalcin . rootSection - Data layout of remote points for outgoing data (this is layout for SF roots) 317b0c7db22SLisandro Dalcin - leafSection - Data layout of local points for incoming data (this is layout for SF leaves) 318b0c7db22SLisandro Dalcin 319b0c7db22SLisandro Dalcin Output Parameter: 320b0c7db22SLisandro Dalcin . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL 321b0c7db22SLisandro Dalcin 322b0c7db22SLisandro Dalcin Level: developer 323b0c7db22SLisandro Dalcin 324b0c7db22SLisandro Dalcin .seealso: PetscSFCreate() 325b0c7db22SLisandro Dalcin @*/ 326b0c7db22SLisandro Dalcin PetscErrorCode PetscSFCreateRemoteOffsets(PetscSF sf, PetscSection rootSection, PetscSection leafSection, PetscInt **remoteOffsets) 327b0c7db22SLisandro Dalcin { 328b0c7db22SLisandro Dalcin PetscSF embedSF; 329b0c7db22SLisandro Dalcin const PetscInt *indices; 330b0c7db22SLisandro Dalcin IS selected; 331b0c7db22SLisandro Dalcin PetscInt numRoots, rpStart = 0, rpEnd = 0, lpStart = 0, lpEnd = 0; 332b0c7db22SLisandro Dalcin PetscErrorCode ierr; 333b0c7db22SLisandro Dalcin 334b0c7db22SLisandro Dalcin PetscFunctionBegin; 335b0c7db22SLisandro Dalcin *remoteOffsets = NULL; 336b0c7db22SLisandro Dalcin ierr = PetscSFGetGraph(sf, &numRoots, NULL, NULL, NULL);CHKERRQ(ierr); 337b0c7db22SLisandro Dalcin if (numRoots < 0) PetscFunctionReturn(0); 338b0c7db22SLisandro Dalcin ierr = PetscLogEventBegin(PETSCSF_RemoteOff,sf,0,0,0);CHKERRQ(ierr); 339b0c7db22SLisandro Dalcin ierr = PetscSectionGetChart(rootSection, &rpStart, &rpEnd);CHKERRQ(ierr); 340b0c7db22SLisandro Dalcin ierr = PetscSectionGetChart(leafSection, &lpStart, &lpEnd);CHKERRQ(ierr); 341b0c7db22SLisandro Dalcin ierr = ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected);CHKERRQ(ierr); 342b0c7db22SLisandro Dalcin ierr = ISGetIndices(selected, &indices);CHKERRQ(ierr); 34372502a1fSJunchao Zhang ierr = PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF);CHKERRQ(ierr); 344b0c7db22SLisandro Dalcin ierr = ISRestoreIndices(selected, &indices);CHKERRQ(ierr); 345b0c7db22SLisandro Dalcin ierr = ISDestroy(&selected);CHKERRQ(ierr); 346b0c7db22SLisandro Dalcin ierr = PetscCalloc1(lpEnd - lpStart, remoteOffsets);CHKERRQ(ierr); 347ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 348ad227feaSJunchao Zhang ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 349b0c7db22SLisandro Dalcin ierr = PetscSFDestroy(&embedSF);CHKERRQ(ierr); 350b0c7db22SLisandro Dalcin ierr = PetscLogEventEnd(PETSCSF_RemoteOff,sf,0,0,0);CHKERRQ(ierr); 351b0c7db22SLisandro Dalcin PetscFunctionReturn(0); 352b0c7db22SLisandro Dalcin } 353b0c7db22SLisandro Dalcin 354b0c7db22SLisandro Dalcin /*@C 355b0c7db22SLisandro Dalcin PetscSFCreateSectionSF - Create an expanded SF of dofs, assuming the input SF relates points 356b0c7db22SLisandro Dalcin 357b0c7db22SLisandro Dalcin Collective on sf 358b0c7db22SLisandro Dalcin 359b0c7db22SLisandro Dalcin Input Parameters: 360b0c7db22SLisandro Dalcin + sf - The SF 361b0c7db22SLisandro Dalcin . rootSection - Data layout of remote points for outgoing data (this is usually the serial section) 362b0c7db22SLisandro Dalcin . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL 363b0c7db22SLisandro Dalcin - leafSection - Data layout of local points for incoming data (this is the distributed section) 364b0c7db22SLisandro Dalcin 365b0c7db22SLisandro Dalcin Output Parameters: 366b0c7db22SLisandro Dalcin - sectionSF - The new SF 367b0c7db22SLisandro Dalcin 368b0c7db22SLisandro Dalcin Note: Either rootSection or remoteOffsets can be specified 369b0c7db22SLisandro Dalcin 370b0c7db22SLisandro Dalcin Level: advanced 371b0c7db22SLisandro Dalcin 372b0c7db22SLisandro Dalcin .seealso: PetscSFCreate() 373b0c7db22SLisandro Dalcin @*/ 374b0c7db22SLisandro Dalcin PetscErrorCode PetscSFCreateSectionSF(PetscSF sf, PetscSection rootSection, PetscInt remoteOffsets[], PetscSection leafSection, PetscSF *sectionSF) 375b0c7db22SLisandro Dalcin { 376b0c7db22SLisandro Dalcin MPI_Comm comm; 377b0c7db22SLisandro Dalcin const PetscInt *localPoints; 378b0c7db22SLisandro Dalcin const PetscSFNode *remotePoints; 379b0c7db22SLisandro Dalcin PetscInt lpStart, lpEnd; 380b0c7db22SLisandro Dalcin PetscInt numRoots, numSectionRoots, numPoints, numIndices = 0; 381b0c7db22SLisandro Dalcin PetscInt *localIndices; 382b0c7db22SLisandro Dalcin PetscSFNode *remoteIndices; 383b0c7db22SLisandro Dalcin PetscInt i, ind; 384b0c7db22SLisandro Dalcin PetscErrorCode ierr; 385b0c7db22SLisandro Dalcin 386b0c7db22SLisandro Dalcin PetscFunctionBegin; 387b0c7db22SLisandro Dalcin PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 388b0c7db22SLisandro Dalcin PetscValidPointer(rootSection,2); 389b0c7db22SLisandro Dalcin /* Cannot check PetscValidIntPointer(remoteOffsets,3) because it can be NULL if sf does not reference any points in leafSection */ 390b0c7db22SLisandro Dalcin PetscValidPointer(leafSection,4); 391b0c7db22SLisandro Dalcin PetscValidPointer(sectionSF,5); 392b0c7db22SLisandro Dalcin ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr); 393b0c7db22SLisandro Dalcin ierr = PetscSFCreate(comm, sectionSF);CHKERRQ(ierr); 394b0c7db22SLisandro Dalcin ierr = PetscSectionGetChart(leafSection, &lpStart, &lpEnd);CHKERRQ(ierr); 395b0c7db22SLisandro Dalcin ierr = PetscSectionGetStorageSize(rootSection, &numSectionRoots);CHKERRQ(ierr); 396b0c7db22SLisandro Dalcin ierr = PetscSFGetGraph(sf, &numRoots, &numPoints, &localPoints, &remotePoints);CHKERRQ(ierr); 397b0c7db22SLisandro Dalcin if (numRoots < 0) PetscFunctionReturn(0); 398b0c7db22SLisandro Dalcin ierr = PetscLogEventBegin(PETSCSF_SectSF,sf,0,0,0);CHKERRQ(ierr); 399b0c7db22SLisandro Dalcin for (i = 0; i < numPoints; ++i) { 400b0c7db22SLisandro Dalcin PetscInt localPoint = localPoints ? localPoints[i] : i; 401b0c7db22SLisandro Dalcin PetscInt dof; 402b0c7db22SLisandro Dalcin 403b0c7db22SLisandro Dalcin if ((localPoint >= lpStart) && (localPoint < lpEnd)) { 404b0c7db22SLisandro Dalcin ierr = PetscSectionGetDof(leafSection, localPoint, &dof);CHKERRQ(ierr); 405b0c7db22SLisandro Dalcin numIndices += dof; 406b0c7db22SLisandro Dalcin } 407b0c7db22SLisandro Dalcin } 408b0c7db22SLisandro Dalcin ierr = PetscMalloc1(numIndices, &localIndices);CHKERRQ(ierr); 409b0c7db22SLisandro Dalcin ierr = PetscMalloc1(numIndices, &remoteIndices);CHKERRQ(ierr); 410b0c7db22SLisandro Dalcin /* Create new index graph */ 411b0c7db22SLisandro Dalcin for (i = 0, ind = 0; i < numPoints; ++i) { 412b0c7db22SLisandro Dalcin PetscInt localPoint = localPoints ? localPoints[i] : i; 413b0c7db22SLisandro Dalcin PetscInt rank = remotePoints[i].rank; 414b0c7db22SLisandro Dalcin 415b0c7db22SLisandro Dalcin if ((localPoint >= lpStart) && (localPoint < lpEnd)) { 416b0c7db22SLisandro Dalcin PetscInt remoteOffset = remoteOffsets[localPoint-lpStart]; 417b0c7db22SLisandro Dalcin PetscInt loff, dof, d; 418b0c7db22SLisandro Dalcin 419b0c7db22SLisandro Dalcin ierr = PetscSectionGetOffset(leafSection, localPoint, &loff);CHKERRQ(ierr); 420b0c7db22SLisandro Dalcin ierr = PetscSectionGetDof(leafSection, localPoint, &dof);CHKERRQ(ierr); 421b0c7db22SLisandro Dalcin for (d = 0; d < dof; ++d, ++ind) { 422b0c7db22SLisandro Dalcin localIndices[ind] = loff+d; 423b0c7db22SLisandro Dalcin remoteIndices[ind].rank = rank; 424b0c7db22SLisandro Dalcin remoteIndices[ind].index = remoteOffset+d; 425b0c7db22SLisandro Dalcin } 426b0c7db22SLisandro Dalcin } 427b0c7db22SLisandro Dalcin } 4282c71b3e2SJacob Faibussowitsch PetscCheckFalse(numIndices != ind,comm, PETSC_ERR_PLIB, "Inconsistency in indices, %" PetscInt_FMT " should be %" PetscInt_FMT, ind, numIndices); 429b0c7db22SLisandro Dalcin ierr = PetscSFSetGraph(*sectionSF, numSectionRoots, numIndices, localIndices, PETSC_OWN_POINTER, remoteIndices, PETSC_OWN_POINTER);CHKERRQ(ierr); 430b0c7db22SLisandro Dalcin ierr = PetscSFSetUp(*sectionSF);CHKERRQ(ierr); 431b0c7db22SLisandro Dalcin ierr = PetscLogEventEnd(PETSCSF_SectSF,sf,0,0,0);CHKERRQ(ierr); 432b0c7db22SLisandro Dalcin PetscFunctionReturn(0); 433b0c7db22SLisandro Dalcin } 4348fa9e22eSVaclav Hapla 4358fa9e22eSVaclav Hapla /*@C 4368fa9e22eSVaclav Hapla PetscSFCreateFromLayouts - Creates a parallel star forest mapping two PetscLayout objects 4378fa9e22eSVaclav Hapla 4388fa9e22eSVaclav Hapla Collective 4398fa9e22eSVaclav Hapla 4404165533cSJose E. Roman Input Parameters: 4418fa9e22eSVaclav Hapla + rmap - PetscLayout defining the global root space 4428fa9e22eSVaclav Hapla - lmap - PetscLayout defining the global leaf space 4438fa9e22eSVaclav Hapla 4444165533cSJose E. Roman Output Parameter: 4458fa9e22eSVaclav Hapla . sf - The parallel star forest 4468fa9e22eSVaclav Hapla 4478fa9e22eSVaclav Hapla Level: intermediate 4488fa9e22eSVaclav Hapla 4498fa9e22eSVaclav Hapla .seealso: PetscSFCreate(), PetscLayoutCreate(), PetscSFSetGraphLayout() 4508fa9e22eSVaclav Hapla @*/ 4518fa9e22eSVaclav Hapla PetscErrorCode PetscSFCreateFromLayouts(PetscLayout rmap, PetscLayout lmap, PetscSF* sf) 4528fa9e22eSVaclav Hapla { 4538fa9e22eSVaclav Hapla PetscErrorCode ierr; 4548fa9e22eSVaclav Hapla PetscInt i,nroots,nleaves = 0; 4558fa9e22eSVaclav Hapla PetscInt rN, lst, len; 4568fa9e22eSVaclav Hapla PetscMPIInt owner = -1; 4578fa9e22eSVaclav Hapla PetscSFNode *remote; 4588fa9e22eSVaclav Hapla MPI_Comm rcomm = rmap->comm; 4598fa9e22eSVaclav Hapla MPI_Comm lcomm = lmap->comm; 4608fa9e22eSVaclav Hapla PetscMPIInt flg; 4618fa9e22eSVaclav Hapla 4628fa9e22eSVaclav Hapla PetscFunctionBegin; 4638fa9e22eSVaclav Hapla PetscValidPointer(sf,3); 4642c71b3e2SJacob Faibussowitsch PetscCheckFalse(!rmap->setupcalled,rcomm,PETSC_ERR_ARG_WRONGSTATE,"Root layout not setup"); 4652c71b3e2SJacob Faibussowitsch PetscCheckFalse(!lmap->setupcalled,lcomm,PETSC_ERR_ARG_WRONGSTATE,"Leaf layout not setup"); 4668fa9e22eSVaclav Hapla ierr = MPI_Comm_compare(rcomm,lcomm,&flg);CHKERRMPI(ierr); 4672c71b3e2SJacob Faibussowitsch PetscCheckFalse(flg != MPI_CONGRUENT && flg != MPI_IDENT,rcomm,PETSC_ERR_SUP,"cannot map two layouts with non-matching communicators"); 4688fa9e22eSVaclav Hapla ierr = PetscSFCreate(rcomm,sf);CHKERRQ(ierr); 4698fa9e22eSVaclav Hapla ierr = PetscLayoutGetLocalSize(rmap,&nroots);CHKERRQ(ierr); 4708fa9e22eSVaclav Hapla ierr = PetscLayoutGetSize(rmap,&rN);CHKERRQ(ierr); 4718fa9e22eSVaclav Hapla ierr = PetscLayoutGetRange(lmap,&lst,&len);CHKERRQ(ierr); 4728fa9e22eSVaclav Hapla ierr = PetscMalloc1(len-lst,&remote);CHKERRQ(ierr); 4738fa9e22eSVaclav Hapla for (i = lst; i < len && i < rN; i++) { 4748fa9e22eSVaclav Hapla if (owner < -1 || i >= rmap->range[owner+1]) { 4758fa9e22eSVaclav Hapla ierr = PetscLayoutFindOwner(rmap,i,&owner);CHKERRQ(ierr); 4768fa9e22eSVaclav Hapla } 4778fa9e22eSVaclav Hapla remote[nleaves].rank = owner; 4788fa9e22eSVaclav Hapla remote[nleaves].index = i - rmap->range[owner]; 4798fa9e22eSVaclav Hapla nleaves++; 4808fa9e22eSVaclav Hapla } 4818fa9e22eSVaclav Hapla ierr = PetscSFSetGraph(*sf,nroots,nleaves,NULL,PETSC_OWN_POINTER,remote,PETSC_COPY_VALUES);CHKERRQ(ierr); 4828fa9e22eSVaclav Hapla ierr = PetscFree(remote);CHKERRQ(ierr); 4838fa9e22eSVaclav Hapla PetscFunctionReturn(0); 4848fa9e22eSVaclav Hapla } 4858fa9e22eSVaclav Hapla 4868fa9e22eSVaclav Hapla /* TODO: handle nooffprocentries like MatZeroRowsMapLocal_Private, since this code is the same */ 4878fa9e22eSVaclav Hapla PetscErrorCode PetscLayoutMapLocal(PetscLayout map,PetscInt N,const PetscInt idxs[], PetscInt *on,PetscInt **oidxs,PetscInt **ogidxs) 4888fa9e22eSVaclav Hapla { 4898fa9e22eSVaclav Hapla PetscInt *owners = map->range; 4908fa9e22eSVaclav Hapla PetscInt n = map->n; 4918fa9e22eSVaclav Hapla PetscSF sf; 4928fa9e22eSVaclav Hapla PetscInt *lidxs,*work = NULL; 4938fa9e22eSVaclav Hapla PetscSFNode *ridxs; 4948fa9e22eSVaclav Hapla PetscMPIInt rank, p = 0; 4958fa9e22eSVaclav Hapla PetscInt r, len = 0; 4968fa9e22eSVaclav Hapla PetscErrorCode ierr; 4978fa9e22eSVaclav Hapla 4988fa9e22eSVaclav Hapla PetscFunctionBegin; 4998fa9e22eSVaclav Hapla if (on) *on = 0; /* squelch -Wmaybe-uninitialized */ 5008fa9e22eSVaclav Hapla /* Create SF where leaves are input idxs and roots are owned idxs */ 5018fa9e22eSVaclav Hapla ierr = MPI_Comm_rank(map->comm,&rank);CHKERRMPI(ierr); 5028fa9e22eSVaclav Hapla ierr = PetscMalloc1(n,&lidxs);CHKERRQ(ierr); 5038fa9e22eSVaclav Hapla for (r = 0; r < n; ++r) lidxs[r] = -1; 5048fa9e22eSVaclav Hapla ierr = PetscMalloc1(N,&ridxs);CHKERRQ(ierr); 5058fa9e22eSVaclav Hapla for (r = 0; r < N; ++r) { 5068fa9e22eSVaclav Hapla const PetscInt idx = idxs[r]; 5072c71b3e2SJacob Faibussowitsch PetscCheckFalse(idx < 0 || map->N <= idx,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index %" PetscInt_FMT " out of range [0,%" PetscInt_FMT ")",idx,map->N); 5088fa9e22eSVaclav Hapla if (idx < owners[p] || owners[p+1] <= idx) { /* short-circuit the search if the last p owns this idx too */ 5098fa9e22eSVaclav Hapla ierr = PetscLayoutFindOwner(map,idx,&p);CHKERRQ(ierr); 5108fa9e22eSVaclav Hapla } 5118fa9e22eSVaclav Hapla ridxs[r].rank = p; 5128fa9e22eSVaclav Hapla ridxs[r].index = idxs[r] - owners[p]; 5138fa9e22eSVaclav Hapla } 5148fa9e22eSVaclav Hapla ierr = PetscSFCreate(map->comm,&sf);CHKERRQ(ierr); 5158fa9e22eSVaclav Hapla ierr = PetscSFSetGraph(sf,n,N,NULL,PETSC_OWN_POINTER,ridxs,PETSC_OWN_POINTER);CHKERRQ(ierr); 5168fa9e22eSVaclav Hapla ierr = PetscSFReduceBegin(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr); 5178fa9e22eSVaclav Hapla ierr = PetscSFReduceEnd(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr); 5188fa9e22eSVaclav Hapla if (ogidxs) { /* communicate global idxs */ 5198fa9e22eSVaclav Hapla PetscInt cum = 0,start,*work2; 5208fa9e22eSVaclav Hapla 5218fa9e22eSVaclav Hapla ierr = PetscMalloc1(n,&work);CHKERRQ(ierr); 5228fa9e22eSVaclav Hapla ierr = PetscCalloc1(N,&work2);CHKERRQ(ierr); 5238fa9e22eSVaclav Hapla for (r = 0; r < N; ++r) if (idxs[r] >=0) cum++; 5248fa9e22eSVaclav Hapla ierr = MPI_Scan(&cum,&start,1,MPIU_INT,MPI_SUM,map->comm);CHKERRMPI(ierr); 5258fa9e22eSVaclav Hapla start -= cum; 5268fa9e22eSVaclav Hapla cum = 0; 5278fa9e22eSVaclav Hapla for (r = 0; r < N; ++r) if (idxs[r] >=0) work2[r] = start+cum++; 52883df288dSJunchao Zhang ierr = PetscSFReduceBegin(sf,MPIU_INT,work2,work,MPI_REPLACE);CHKERRQ(ierr); 52983df288dSJunchao Zhang ierr = PetscSFReduceEnd(sf,MPIU_INT,work2,work,MPI_REPLACE);CHKERRQ(ierr); 5308fa9e22eSVaclav Hapla ierr = PetscFree(work2);CHKERRQ(ierr); 5318fa9e22eSVaclav Hapla } 5328fa9e22eSVaclav Hapla ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 5338fa9e22eSVaclav Hapla /* Compress and put in indices */ 5348fa9e22eSVaclav Hapla for (r = 0; r < n; ++r) 5358fa9e22eSVaclav Hapla if (lidxs[r] >= 0) { 5368fa9e22eSVaclav Hapla if (work) work[len] = work[r]; 5378fa9e22eSVaclav Hapla lidxs[len++] = r; 5388fa9e22eSVaclav Hapla } 5398fa9e22eSVaclav Hapla if (on) *on = len; 5408fa9e22eSVaclav Hapla if (oidxs) *oidxs = lidxs; 5418fa9e22eSVaclav Hapla if (ogidxs) *ogidxs = work; 5428fa9e22eSVaclav Hapla PetscFunctionReturn(0); 5438fa9e22eSVaclav Hapla } 544deffd5ebSksagiyam 545deffd5ebSksagiyam /*@ 546deffd5ebSksagiyam PetscSFCreateByMatchingIndices - Create SF by matching root and leaf indices 547deffd5ebSksagiyam 548deffd5ebSksagiyam Collective 549deffd5ebSksagiyam 550deffd5ebSksagiyam Input Parameters: 551deffd5ebSksagiyam + layout - PetscLayout defining the global index space and the rank that brokers each index 552deffd5ebSksagiyam . numRootIndices - size of rootIndices 553deffd5ebSksagiyam . rootIndices - PetscInt array of global indices of which this process requests ownership 554deffd5ebSksagiyam . rootLocalIndices - root local index permutation (NULL if no permutation) 555deffd5ebSksagiyam . rootLocalOffset - offset to be added to root local indices 556deffd5ebSksagiyam . numLeafIndices - size of leafIndices 557deffd5ebSksagiyam . leafIndices - PetscInt array of global indices with which this process requires data associated 558deffd5ebSksagiyam . leafLocalIndices - leaf local index permutation (NULL if no permutation) 559deffd5ebSksagiyam - leafLocalOffset - offset to be added to leaf local indices 560deffd5ebSksagiyam 561d8d19677SJose E. Roman Output Parameters: 562deffd5ebSksagiyam + sfA - star forest representing the communication pattern from the layout space to the leaf space (NULL if not needed) 563deffd5ebSksagiyam - sf - star forest representing the communication pattern from the root space to the leaf space 564deffd5ebSksagiyam 565deffd5ebSksagiyam Example 1: 566deffd5ebSksagiyam $ 567deffd5ebSksagiyam $ rank : 0 1 2 568deffd5ebSksagiyam $ rootIndices : [1 0 2] [3] [3] 569deffd5ebSksagiyam $ rootLocalOffset : 100 200 300 570deffd5ebSksagiyam $ layout : [0 1] [2] [3] 571deffd5ebSksagiyam $ leafIndices : [0] [2] [0 3] 572deffd5ebSksagiyam $ leafLocalOffset : 400 500 600 573deffd5ebSksagiyam $ 5747ef5781fSVaclav Hapla would build the following SF 575deffd5ebSksagiyam $ 576deffd5ebSksagiyam $ [0] 400 <- (0,101) 577deffd5ebSksagiyam $ [1] 500 <- (0,102) 578deffd5ebSksagiyam $ [2] 600 <- (0,101) 579deffd5ebSksagiyam $ [2] 601 <- (2,300) 580deffd5ebSksagiyam $ 581deffd5ebSksagiyam Example 2: 582deffd5ebSksagiyam $ 583deffd5ebSksagiyam $ rank : 0 1 2 584deffd5ebSksagiyam $ rootIndices : [1 0 2] [3] [3] 585deffd5ebSksagiyam $ rootLocalOffset : 100 200 300 586deffd5ebSksagiyam $ layout : [0 1] [2] [3] 587deffd5ebSksagiyam $ leafIndices : rootIndices rootIndices rootIndices 588deffd5ebSksagiyam $ leafLocalOffset : rootLocalOffset rootLocalOffset rootLocalOffset 589deffd5ebSksagiyam $ 5907ef5781fSVaclav Hapla would build the following SF 591deffd5ebSksagiyam $ 592deffd5ebSksagiyam $ [1] 200 <- (2,300) 593deffd5ebSksagiyam $ 594deffd5ebSksagiyam Example 3: 595deffd5ebSksagiyam $ 596deffd5ebSksagiyam $ No process requests ownership of global index 1, but no process needs it. 597deffd5ebSksagiyam $ 598deffd5ebSksagiyam $ rank : 0 1 2 599deffd5ebSksagiyam $ numRootIndices : 2 1 1 600deffd5ebSksagiyam $ rootIndices : [0 2] [3] [3] 601deffd5ebSksagiyam $ rootLocalOffset : 100 200 300 602deffd5ebSksagiyam $ layout : [0 1] [2] [3] 603deffd5ebSksagiyam $ numLeafIndices : 1 1 2 604deffd5ebSksagiyam $ leafIndices : [0] [2] [0 3] 605deffd5ebSksagiyam $ leafLocalOffset : 400 500 600 606deffd5ebSksagiyam $ 6077ef5781fSVaclav Hapla would build the following SF 608deffd5ebSksagiyam $ 609deffd5ebSksagiyam $ [0] 400 <- (0,100) 610deffd5ebSksagiyam $ [1] 500 <- (0,101) 611deffd5ebSksagiyam $ [2] 600 <- (0,100) 612deffd5ebSksagiyam $ [2] 601 <- (2,300) 613deffd5ebSksagiyam $ 614deffd5ebSksagiyam 615deffd5ebSksagiyam Notes: 616deffd5ebSksagiyam The layout parameter represents any partitioning of [0, N), where N is the total number of global indices, and its 617deffd5ebSksagiyam local size can be set to PETSC_DECIDE. 618deffd5ebSksagiyam If a global index x lies in the partition owned by process i, each process whose rootIndices contains x requests 619deffd5ebSksagiyam ownership of x and sends its own rank and the local index of x to process i. 620deffd5ebSksagiyam If multiple processes request ownership of x, the one with the highest rank is to own x. 621deffd5ebSksagiyam Process i then broadcasts the ownership information, so that each process whose leafIndices contains x knows the 622deffd5ebSksagiyam ownership information of x. 623deffd5ebSksagiyam The output sf is constructed by associating each leaf point to a root point in this way. 624deffd5ebSksagiyam 625deffd5ebSksagiyam Suppose there is point data ordered according to the global indices and partitioned according to the given layout. 626deffd5ebSksagiyam The optional output PetscSF object sfA can be used to push such data to leaf points. 627deffd5ebSksagiyam 628deffd5ebSksagiyam All indices in rootIndices and leafIndices must lie in the layout range. The union (over all processes) of rootIndices 629deffd5ebSksagiyam must cover that of leafIndices, but need not cover the entire layout. 630deffd5ebSksagiyam 631deffd5ebSksagiyam If (leafIndices, leafLocalIndices, leafLocalOffset) == (rootIndices, rootLocalIndices, rootLocalOffset), the output 632deffd5ebSksagiyam star forest is almost identity, so will only include non-trivial part of the map. 633deffd5ebSksagiyam 634deffd5ebSksagiyam Developer Notes: 635deffd5ebSksagiyam Current approach of a process of the highest rank gaining the ownership may cause load imbalance; consider using 636deffd5ebSksagiyam hash(rank, root_local_index) as the bid for the ownership determination. 637deffd5ebSksagiyam 638deffd5ebSksagiyam Level: advanced 639deffd5ebSksagiyam 640deffd5ebSksagiyam .seealso: PetscSFCreate() 641deffd5ebSksagiyam @*/ 642deffd5ebSksagiyam 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) 643deffd5ebSksagiyam { 644deffd5ebSksagiyam MPI_Comm comm = layout->comm; 645deffd5ebSksagiyam PetscMPIInt size, rank; 646deffd5ebSksagiyam PetscSF sf1; 647deffd5ebSksagiyam PetscSFNode *owners, *buffer, *iremote; 648deffd5ebSksagiyam PetscInt *ilocal, nleaves, N, n, i; 649deffd5ebSksagiyam #if defined(PETSC_USE_DEBUG) 650deffd5ebSksagiyam PetscInt N1; 651deffd5ebSksagiyam #endif 652deffd5ebSksagiyam PetscBool flag; 653deffd5ebSksagiyam PetscErrorCode ierr; 654deffd5ebSksagiyam 655deffd5ebSksagiyam PetscFunctionBegin; 656deffd5ebSksagiyam if (rootIndices) PetscValidIntPointer(rootIndices,3); 657deffd5ebSksagiyam if (rootLocalIndices) PetscValidIntPointer(rootLocalIndices,4); 658deffd5ebSksagiyam if (leafIndices) PetscValidIntPointer(leafIndices,7); 659deffd5ebSksagiyam if (leafLocalIndices) PetscValidIntPointer(leafLocalIndices,8); 660deffd5ebSksagiyam if (sfA) PetscValidPointer(sfA,10); 661deffd5ebSksagiyam PetscValidPointer(sf,11); 6622c71b3e2SJacob Faibussowitsch PetscCheckFalse(numRootIndices < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numRootIndices (%" PetscInt_FMT ") must be non-negative", numRootIndices); 6632c71b3e2SJacob Faibussowitsch PetscCheckFalse(numLeafIndices < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numLeafIndices (%" PetscInt_FMT ") must be non-negative", numLeafIndices); 664deffd5ebSksagiyam ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr); 665deffd5ebSksagiyam ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 666deffd5ebSksagiyam ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr); 667deffd5ebSksagiyam ierr = PetscLayoutGetSize(layout, &N);CHKERRQ(ierr); 668deffd5ebSksagiyam ierr = PetscLayoutGetLocalSize(layout, &n);CHKERRQ(ierr); 669deffd5ebSksagiyam flag = (PetscBool)(leafIndices == rootIndices); 670820f2d46SBarry Smith ierr = MPIU_Allreduce(MPI_IN_PLACE, &flag, 1, MPIU_BOOL, MPI_LAND, comm);CHKERRMPI(ierr); 6712c71b3e2SJacob Faibussowitsch PetscCheckFalse(flag && numLeafIndices != numRootIndices,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "leafIndices == rootIndices, but numLeafIndices (%" PetscInt_FMT ") != numRootIndices(%" PetscInt_FMT ")", numLeafIndices, numRootIndices); 672deffd5ebSksagiyam #if defined(PETSC_USE_DEBUG) 673deffd5ebSksagiyam N1 = PETSC_MIN_INT; 674deffd5ebSksagiyam for (i = 0; i < numRootIndices; i++) if (rootIndices[i] > N1) N1 = rootIndices[i]; 675deffd5ebSksagiyam ierr = MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm);CHKERRMPI(ierr); 6762c71b3e2SJacob Faibussowitsch PetscCheckFalse(N1 >= N,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Max. root index (%" PetscInt_FMT ") out of layout range [0,%" PetscInt_FMT ")", N1, N); 677deffd5ebSksagiyam if (!flag) { 678deffd5ebSksagiyam N1 = PETSC_MIN_INT; 679deffd5ebSksagiyam for (i = 0; i < numLeafIndices; i++) if (leafIndices[i] > N1) N1 = leafIndices[i]; 680deffd5ebSksagiyam ierr = MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm);CHKERRMPI(ierr); 6812c71b3e2SJacob Faibussowitsch PetscCheckFalse(N1 >= N,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Max. leaf index (%" PetscInt_FMT ") out of layout range [0,%" PetscInt_FMT ")", N1, N); 682deffd5ebSksagiyam } 683deffd5ebSksagiyam #endif 684deffd5ebSksagiyam /* Reduce: owners -> buffer */ 685deffd5ebSksagiyam ierr = PetscMalloc1(n, &buffer);CHKERRQ(ierr); 686deffd5ebSksagiyam ierr = PetscSFCreate(comm, &sf1);CHKERRQ(ierr); 687deffd5ebSksagiyam ierr = PetscSFSetFromOptions(sf1);CHKERRQ(ierr); 688deffd5ebSksagiyam ierr = PetscSFSetGraphLayout(sf1, layout, numRootIndices, NULL, PETSC_OWN_POINTER, rootIndices);CHKERRQ(ierr); 689deffd5ebSksagiyam ierr = PetscMalloc1(numRootIndices, &owners);CHKERRQ(ierr); 690deffd5ebSksagiyam for (i = 0; i < numRootIndices; ++i) { 691deffd5ebSksagiyam owners[i].rank = rank; 692deffd5ebSksagiyam owners[i].index = rootLocalOffset + (rootLocalIndices ? rootLocalIndices[i] : i); 693deffd5ebSksagiyam } 694deffd5ebSksagiyam for (i = 0; i < n; ++i) { 695deffd5ebSksagiyam buffer[i].index = -1; 696deffd5ebSksagiyam buffer[i].rank = -1; 697deffd5ebSksagiyam } 698deffd5ebSksagiyam ierr = PetscSFReduceBegin(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC);CHKERRQ(ierr); 699deffd5ebSksagiyam ierr = PetscSFReduceEnd(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC);CHKERRQ(ierr); 700deffd5ebSksagiyam /* Bcast: buffer -> owners */ 701deffd5ebSksagiyam if (!flag) { 702deffd5ebSksagiyam /* leafIndices is different from rootIndices */ 703deffd5ebSksagiyam ierr = PetscFree(owners);CHKERRQ(ierr); 704deffd5ebSksagiyam ierr = PetscSFSetGraphLayout(sf1, layout, numLeafIndices, NULL, PETSC_OWN_POINTER, leafIndices);CHKERRQ(ierr); 705deffd5ebSksagiyam ierr = PetscMalloc1(numLeafIndices, &owners);CHKERRQ(ierr); 706deffd5ebSksagiyam } 707deffd5ebSksagiyam ierr = PetscSFBcastBegin(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE);CHKERRQ(ierr); 708deffd5ebSksagiyam ierr = PetscSFBcastEnd(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE);CHKERRQ(ierr); 7092c71b3e2SJacob Faibussowitsch for (i = 0; i < numLeafIndices; ++i) PetscCheckFalse(owners[i].rank < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Global point %" PetscInt_FMT " was unclaimed", leafIndices[i]); 710deffd5ebSksagiyam ierr = PetscFree(buffer);CHKERRQ(ierr); 711deffd5ebSksagiyam if (sfA) {*sfA = sf1;} 712deffd5ebSksagiyam else {ierr = PetscSFDestroy(&sf1);CHKERRQ(ierr);} 713deffd5ebSksagiyam /* Create sf */ 714deffd5ebSksagiyam if (flag && rootLocalIndices == leafLocalIndices && leafLocalOffset == rootLocalOffset) { 715deffd5ebSksagiyam /* leaf space == root space */ 716deffd5ebSksagiyam for (i = 0, nleaves = 0; i < numLeafIndices; ++i) if (owners[i].rank != rank) ++nleaves; 717deffd5ebSksagiyam ierr = PetscMalloc1(nleaves, &ilocal);CHKERRQ(ierr); 718deffd5ebSksagiyam ierr = PetscMalloc1(nleaves, &iremote);CHKERRQ(ierr); 719deffd5ebSksagiyam for (i = 0, nleaves = 0; i < numLeafIndices; ++i) { 720deffd5ebSksagiyam if (owners[i].rank != rank) { 721deffd5ebSksagiyam ilocal[nleaves] = leafLocalOffset + i; 722deffd5ebSksagiyam iremote[nleaves].rank = owners[i].rank; 723deffd5ebSksagiyam iremote[nleaves].index = owners[i].index; 724deffd5ebSksagiyam ++nleaves; 725deffd5ebSksagiyam } 726deffd5ebSksagiyam } 727deffd5ebSksagiyam ierr = PetscFree(owners);CHKERRQ(ierr); 728deffd5ebSksagiyam } else { 729deffd5ebSksagiyam nleaves = numLeafIndices; 730deffd5ebSksagiyam ierr = PetscMalloc1(nleaves, &ilocal);CHKERRQ(ierr); 731deffd5ebSksagiyam for (i = 0; i < nleaves; ++i) {ilocal[i] = leafLocalOffset + (leafLocalIndices ? leafLocalIndices[i] : i);} 732deffd5ebSksagiyam iremote = owners; 733deffd5ebSksagiyam } 734deffd5ebSksagiyam ierr = PetscSFCreate(comm, sf);CHKERRQ(ierr); 735deffd5ebSksagiyam ierr = PetscSFSetFromOptions(*sf);CHKERRQ(ierr); 736deffd5ebSksagiyam ierr = PetscSFSetGraph(*sf, rootLocalOffset + numRootIndices, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr); 737deffd5ebSksagiyam PetscFunctionReturn(0); 738deffd5ebSksagiyam } 739