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