xref: /petsc/src/vec/is/sf/utils/sfutils.c (revision 147403d98689287189d69e47992b7c8152b2c9da)
1 #include <petsc/private/sfimpl.h>       /*I  "petscsf.h"   I*/
2 #include <petsc/private/sectionimpl.h>
3 
4 /*@C
5    PetscSFSetGraphLayout - Set a parallel star forest via global indices and a PetscLayout
6 
7    Collective
8 
9    Input Parameters:
10 +  sf - star forest
11 .  layout - PetscLayout defining the global space for roots
12 .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
13 .  ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage
14 .  localmode - copy mode for ilocal
15 -  iremote - root vertices in global numbering corresponding to leaves in ilocal
16 
17    Level: intermediate
18 
19    Notes:
20    Global indices must lie in [0, N) where N is the global size of layout.
21    Leaf indices in ilocal get sorted; this means the user-provided array gets sorted if localmode is PETSC_OWN_POINTER.
22 
23    Developer Notes:
24    Local indices which are the identity permutation in the range [0,nleaves) are discarded as they
25    encode contiguous storage. In such case, if localmode is PETSC_OWN_POINTER, the memory is deallocated as it is not
26    needed
27 
28 .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph(), PetscSFGetGraph()
29 @*/
30 PetscErrorCode PetscSFSetGraphLayout(PetscSF sf,PetscLayout layout,PetscInt nleaves,PetscInt *ilocal,PetscCopyMode localmode,const PetscInt *iremote)
31 {
32   PetscErrorCode ierr;
33   const PetscInt *range;
34   PetscInt       i, nroots, ls = -1, ln = -1;
35   PetscMPIInt    lr = -1;
36   PetscSFNode    *remote;
37 
38   PetscFunctionBegin;
39   ierr = PetscLayoutGetLocalSize(layout,&nroots);CHKERRQ(ierr);
40   ierr = PetscLayoutGetRanges(layout,&range);CHKERRQ(ierr);
41   ierr = PetscMalloc1(nleaves,&remote);CHKERRQ(ierr);
42   if (nleaves) { ls = iremote[0] + 1; }
43   for (i=0; i<nleaves; i++) {
44     const PetscInt idx = iremote[i] - ls;
45     if (idx < 0 || idx >= ln) { /* short-circuit the search */
46       ierr = PetscLayoutFindOwnerIndex(layout,iremote[i],&lr,&remote[i].index);CHKERRQ(ierr);
47       remote[i].rank = lr;
48       ls = range[lr];
49       ln = range[lr+1] - ls;
50     } else {
51       remote[i].rank  = lr;
52       remote[i].index = idx;
53     }
54   }
55   ierr = PetscSFSetGraph(sf,nroots,nleaves,ilocal,localmode,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
56   PetscFunctionReturn(0);
57 }
58 
59 /*@
60   PetscSFSetGraphSection - Sets the PetscSF graph encoding the parallel dof overlap based upon the PetscSections
61   describing the data layout.
62 
63   Input Parameters:
64 + sf - The SF
65 . localSection - PetscSection describing the local data layout
66 - globalSection - PetscSection describing the global data layout
67 
68   Level: developer
69 
70 .seealso: PetscSFSetGraph(), PetscSFSetGraphLayout()
71 @*/
72 PetscErrorCode PetscSFSetGraphSection(PetscSF sf, PetscSection localSection, PetscSection globalSection)
73 {
74   MPI_Comm       comm;
75   PetscLayout    layout;
76   const PetscInt *ranges;
77   PetscInt       *local;
78   PetscSFNode    *remote;
79   PetscInt       pStart, pEnd, p, nroots, nleaves = 0, l;
80   PetscMPIInt    size, rank;
81   PetscErrorCode ierr;
82 
83   PetscFunctionBegin;
84   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
85   PetscValidHeaderSpecific(localSection, PETSC_SECTION_CLASSID, 2);
86   PetscValidHeaderSpecific(globalSection, PETSC_SECTION_CLASSID, 3);
87 
88   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
89   ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
90   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
91   ierr = PetscSectionGetChart(globalSection, &pStart, &pEnd);CHKERRQ(ierr);
92   ierr = PetscSectionGetConstrainedStorageSize(globalSection, &nroots);CHKERRQ(ierr);
93   ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr);
94   ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr);
95   ierr = PetscLayoutSetLocalSize(layout, nroots);CHKERRQ(ierr);
96   ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr);
97   ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr);
98   for (p = pStart; p < pEnd; ++p) {
99     PetscInt gdof, gcdof;
100 
101     ierr = PetscSectionGetDof(globalSection, p, &gdof);CHKERRQ(ierr);
102     ierr = PetscSectionGetConstraintDof(globalSection, p, &gcdof);CHKERRQ(ierr);
103     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));
104     nleaves += gdof < 0 ? -(gdof+1)-gcdof : gdof-gcdof;
105   }
106   ierr = PetscMalloc1(nleaves, &local);CHKERRQ(ierr);
107   ierr = PetscMalloc1(nleaves, &remote);CHKERRQ(ierr);
108   for (p = pStart, l = 0; p < pEnd; ++p) {
109     const PetscInt *cind;
110     PetscInt       dof, cdof, off, gdof, gcdof, goff, gsize, d, c;
111 
112     ierr = PetscSectionGetDof(localSection, p, &dof);CHKERRQ(ierr);
113     ierr = PetscSectionGetOffset(localSection, p, &off);CHKERRQ(ierr);
114     ierr = PetscSectionGetConstraintDof(localSection, p, &cdof);CHKERRQ(ierr);
115     ierr = PetscSectionGetConstraintIndices(localSection, p, &cind);CHKERRQ(ierr);
116     ierr = PetscSectionGetDof(globalSection, p, &gdof);CHKERRQ(ierr);
117     ierr = PetscSectionGetConstraintDof(globalSection, p, &gcdof);CHKERRQ(ierr);
118     ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr);
119     if (!gdof) continue; /* Censored point */
120     gsize = gdof < 0 ? -(gdof+1)-gcdof : gdof-gcdof;
121     if (gsize != dof-cdof) {
122       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);
123       cdof = 0; /* Ignore constraints */
124     }
125     for (d = 0, c = 0; d < dof; ++d) {
126       if ((c < cdof) && (cind[c] == d)) {++c; continue;}
127       local[l+d-c] = off+d;
128     }
129     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);
130     if (gdof < 0) {
131       for (d = 0; d < gsize; ++d, ++l) {
132         PetscInt offset = -(goff+1) + d, r;
133 
134         ierr = PetscFindInt(offset,size+1,ranges,&r);CHKERRQ(ierr);
135         if (r < 0) r = -(r+2);
136         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);
137         remote[l].rank  = r;
138         remote[l].index = offset - ranges[r];
139       }
140     } else {
141       for (d = 0; d < gsize; ++d, ++l) {
142         remote[l].rank  = rank;
143         remote[l].index = goff+d - ranges[rank];
144       }
145     }
146   }
147   PetscCheckFalse(l != nleaves,comm, PETSC_ERR_PLIB, "Iteration error, l %" PetscInt_FMT " != nleaves %" PetscInt_FMT, l, nleaves);
148   ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr);
149   ierr = PetscSFSetGraph(sf, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER);CHKERRQ(ierr);
150   PetscFunctionReturn(0);
151 }
152 
153 static PetscErrorCode PetscSectionCheckConstraints_Static(PetscSection s)
154 {
155   PetscErrorCode ierr;
156 
157   PetscFunctionBegin;
158   if (!s->bc) {
159     ierr = PetscSectionCreate(PETSC_COMM_SELF, &s->bc);CHKERRQ(ierr);
160     ierr = PetscSectionSetChart(s->bc, s->pStart, s->pEnd);CHKERRQ(ierr);
161   }
162   PetscFunctionReturn(0);
163 }
164 
165 /*@C
166   PetscSFDistributeSection - Create a new PetscSection reorganized, moving from the root to the leaves of the SF
167 
168   Collective on sf
169 
170   Input Parameters:
171 + sf - The SF
172 - rootSection - Section defined on root space
173 
174   Output Parameters:
175 + remoteOffsets - root offsets in leaf storage, or NULL
176 - leafSection - Section defined on the leaf space
177 
178   Level: advanced
179 
180 .seealso: PetscSFCreate()
181 @*/
182 PetscErrorCode PetscSFDistributeSection(PetscSF sf, PetscSection rootSection, PetscInt **remoteOffsets, PetscSection leafSection)
183 {
184   PetscSF        embedSF;
185   const PetscInt *indices;
186   IS             selected;
187   PetscInt       numFields, nroots, rpStart, rpEnd, lpStart = PETSC_MAX_INT, lpEnd = -1, f, c;
188   PetscBool      *sub, hasc;
189   PetscErrorCode ierr;
190 
191   PetscFunctionBegin;
192   ierr = PetscLogEventBegin(PETSCSF_DistSect,sf,0,0,0);CHKERRQ(ierr);
193   ierr = PetscSectionGetNumFields(rootSection, &numFields);CHKERRQ(ierr);
194   if (numFields) {
195     IS perm;
196 
197     /* PetscSectionSetNumFields() calls PetscSectionReset(), which destroys
198        leafSection->perm. To keep this permutation set by the user, we grab
199        the reference before calling PetscSectionSetNumFields() and set it
200        back after. */
201     ierr = PetscSectionGetPermutation(leafSection, &perm);CHKERRQ(ierr);
202     ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
203     ierr = PetscSectionSetNumFields(leafSection, numFields);CHKERRQ(ierr);
204     ierr = PetscSectionSetPermutation(leafSection, perm);CHKERRQ(ierr);
205     ierr = ISDestroy(&perm);CHKERRQ(ierr);
206   }
207   ierr = PetscMalloc1(numFields+2, &sub);CHKERRQ(ierr);
208   sub[1] = rootSection->bc ? PETSC_TRUE : PETSC_FALSE;
209   for (f = 0; f < numFields; ++f) {
210     PetscSectionSym sym, dsym = NULL;
211     const char      *name   = NULL;
212     PetscInt        numComp = 0;
213 
214     sub[2 + f] = rootSection->field[f]->bc ? PETSC_TRUE : PETSC_FALSE;
215     ierr = PetscSectionGetFieldComponents(rootSection, f, &numComp);CHKERRQ(ierr);
216     ierr = PetscSectionGetFieldName(rootSection, f, &name);CHKERRQ(ierr);
217     ierr = PetscSectionGetFieldSym(rootSection, f, &sym);CHKERRQ(ierr);
218     if (sym) {ierr = PetscSectionSymDistribute(sym, sf, &dsym);CHKERRQ(ierr);}
219     ierr = PetscSectionSetFieldComponents(leafSection, f, numComp);CHKERRQ(ierr);
220     ierr = PetscSectionSetFieldName(leafSection, f, name);CHKERRQ(ierr);
221     ierr = PetscSectionSetFieldSym(leafSection, f, dsym);CHKERRQ(ierr);
222     ierr = PetscSectionSymDestroy(&dsym);CHKERRQ(ierr);
223     for (c = 0; c < rootSection->numFieldComponents[f]; ++c) {
224       ierr = PetscSectionGetComponentName(rootSection, f, c, &name);CHKERRQ(ierr);
225       ierr = PetscSectionSetComponentName(leafSection, f, c, name);CHKERRQ(ierr);
226     }
227   }
228   ierr = PetscSectionGetChart(rootSection, &rpStart, &rpEnd);CHKERRQ(ierr);
229   ierr = PetscSFGetGraph(sf,&nroots,NULL,NULL,NULL);CHKERRQ(ierr);
230   rpEnd = PetscMin(rpEnd,nroots);
231   rpEnd = PetscMax(rpStart,rpEnd);
232   /* see if we can avoid creating the embedded SF, since it can cost more than an allreduce */
233   sub[0] = (PetscBool)(nroots != rpEnd - rpStart);
234   ierr = MPIU_Allreduce(MPI_IN_PLACE, sub, 2+numFields, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)sf));CHKERRMPI(ierr);
235   if (sub[0]) {
236     ierr = ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected);CHKERRQ(ierr);
237     ierr = ISGetIndices(selected, &indices);CHKERRQ(ierr);
238     ierr = PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF);CHKERRQ(ierr);
239     ierr = ISRestoreIndices(selected, &indices);CHKERRQ(ierr);
240     ierr = ISDestroy(&selected);CHKERRQ(ierr);
241   } else {
242     ierr = PetscObjectReference((PetscObject)sf);CHKERRQ(ierr);
243     embedSF = sf;
244   }
245   ierr = PetscSFGetLeafRange(embedSF, &lpStart, &lpEnd);CHKERRQ(ierr);
246   lpEnd++;
247 
248   ierr = PetscSectionSetChart(leafSection, lpStart, lpEnd);CHKERRQ(ierr);
249 
250   /* Constrained dof section */
251   hasc = sub[1];
252   for (f = 0; f < numFields; ++f) hasc = (PetscBool)(hasc || sub[2+f]);
253 
254   /* Could fuse these at the cost of copies and extra allocation */
255   ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart],MPI_REPLACE);CHKERRQ(ierr);
256   ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart],MPI_REPLACE);CHKERRQ(ierr);
257   if (sub[1]) {
258     ierr = PetscSectionCheckConstraints_Static(rootSection);CHKERRQ(ierr);
259     ierr = PetscSectionCheckConstraints_Static(leafSection);CHKERRQ(ierr);
260     ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart],MPI_REPLACE);CHKERRQ(ierr);
261     ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart],MPI_REPLACE);CHKERRQ(ierr);
262   }
263   for (f = 0; f < numFields; ++f) {
264     ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart],MPI_REPLACE);CHKERRQ(ierr);
265     ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart],MPI_REPLACE);CHKERRQ(ierr);
266     if (sub[2+f]) {
267       ierr = PetscSectionCheckConstraints_Static(rootSection->field[f]);CHKERRQ(ierr);
268       ierr = PetscSectionCheckConstraints_Static(leafSection->field[f]);CHKERRQ(ierr);
269       ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart],MPI_REPLACE);CHKERRQ(ierr);
270       ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart],MPI_REPLACE);CHKERRQ(ierr);
271     }
272   }
273   if (remoteOffsets) {
274     ierr = PetscMalloc1(lpEnd - lpStart, remoteOffsets);CHKERRQ(ierr);
275     ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart],MPI_REPLACE);CHKERRQ(ierr);
276     ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart],MPI_REPLACE);CHKERRQ(ierr);
277   }
278   ierr = PetscSectionSetUp(leafSection);CHKERRQ(ierr);
279   if (hasc) { /* need to communicate bcIndices */
280     PetscSF  bcSF;
281     PetscInt *rOffBc;
282 
283     ierr = PetscMalloc1(lpEnd - lpStart, &rOffBc);CHKERRQ(ierr);
284     if (sub[1]) {
285       ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart],MPI_REPLACE);CHKERRQ(ierr);
286       ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart],MPI_REPLACE);CHKERRQ(ierr);
287       ierr = PetscSFCreateSectionSF(embedSF, rootSection->bc, rOffBc, leafSection->bc, &bcSF);CHKERRQ(ierr);
288       ierr = PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices,MPI_REPLACE);CHKERRQ(ierr);
289       ierr = PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices,MPI_REPLACE);CHKERRQ(ierr);
290       ierr = PetscSFDestroy(&bcSF);CHKERRQ(ierr);
291     }
292     for (f = 0; f < numFields; ++f) {
293       if (sub[2+f]) {
294         ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart],MPI_REPLACE);CHKERRQ(ierr);
295         ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart],MPI_REPLACE);CHKERRQ(ierr);
296         ierr = PetscSFCreateSectionSF(embedSF, rootSection->field[f]->bc, rOffBc, leafSection->field[f]->bc, &bcSF);CHKERRQ(ierr);
297         ierr = PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices,MPI_REPLACE);CHKERRQ(ierr);
298         ierr = PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices,MPI_REPLACE);CHKERRQ(ierr);
299         ierr = PetscSFDestroy(&bcSF);CHKERRQ(ierr);
300       }
301     }
302     ierr = PetscFree(rOffBc);CHKERRQ(ierr);
303   }
304   ierr = PetscSFDestroy(&embedSF);CHKERRQ(ierr);
305   ierr = PetscFree(sub);CHKERRQ(ierr);
306   ierr = PetscLogEventEnd(PETSCSF_DistSect,sf,0,0,0);CHKERRQ(ierr);
307   PetscFunctionReturn(0);
308 }
309 
310 /*@C
311   PetscSFCreateRemoteOffsets - Create offsets for point data on remote processes
312 
313   Collective on sf
314 
315   Input Parameters:
316 + sf          - The SF
317 . rootSection - Data layout of remote points for outgoing data (this is layout for SF roots)
318 - leafSection - Data layout of local points for incoming data  (this is layout for SF leaves)
319 
320   Output Parameter:
321 . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL
322 
323   Level: developer
324 
325 .seealso: PetscSFCreate()
326 @*/
327 PetscErrorCode PetscSFCreateRemoteOffsets(PetscSF sf, PetscSection rootSection, PetscSection leafSection, PetscInt **remoteOffsets)
328 {
329   PetscSF         embedSF;
330   const PetscInt *indices;
331   IS              selected;
332   PetscInt        numRoots, rpStart = 0, rpEnd = 0, lpStart = 0, lpEnd = 0;
333   PetscErrorCode  ierr;
334 
335   PetscFunctionBegin;
336   *remoteOffsets = NULL;
337   ierr = PetscSFGetGraph(sf, &numRoots, NULL, NULL, NULL);CHKERRQ(ierr);
338   if (numRoots < 0) PetscFunctionReturn(0);
339   ierr = PetscLogEventBegin(PETSCSF_RemoteOff,sf,0,0,0);CHKERRQ(ierr);
340   ierr = PetscSectionGetChart(rootSection, &rpStart, &rpEnd);CHKERRQ(ierr);
341   ierr = PetscSectionGetChart(leafSection, &lpStart, &lpEnd);CHKERRQ(ierr);
342   ierr = ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected);CHKERRQ(ierr);
343   ierr = ISGetIndices(selected, &indices);CHKERRQ(ierr);
344   ierr = PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF);CHKERRQ(ierr);
345   ierr = ISRestoreIndices(selected, &indices);CHKERRQ(ierr);
346   ierr = ISDestroy(&selected);CHKERRQ(ierr);
347   ierr = PetscCalloc1(lpEnd - lpStart, remoteOffsets);CHKERRQ(ierr);
348   ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart],MPI_REPLACE);CHKERRQ(ierr);
349   ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart],MPI_REPLACE);CHKERRQ(ierr);
350   ierr = PetscSFDestroy(&embedSF);CHKERRQ(ierr);
351   ierr = PetscLogEventEnd(PETSCSF_RemoteOff,sf,0,0,0);CHKERRQ(ierr);
352   PetscFunctionReturn(0);
353 }
354 
355 /*@C
356   PetscSFCreateSectionSF - Create an expanded SF of dofs, assuming the input SF relates points
357 
358   Collective on sf
359 
360   Input Parameters:
361 + sf - The SF
362 . rootSection - Data layout of remote points for outgoing data (this is usually the serial section)
363 . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL
364 - leafSection - Data layout of local points for incoming data  (this is the distributed section)
365 
366   Output Parameters:
367 - sectionSF - The new SF
368 
369   Note: Either rootSection or remoteOffsets can be specified
370 
371   Level: advanced
372 
373 .seealso: PetscSFCreate()
374 @*/
375 PetscErrorCode PetscSFCreateSectionSF(PetscSF sf, PetscSection rootSection, PetscInt remoteOffsets[], PetscSection leafSection, PetscSF *sectionSF)
376 {
377   MPI_Comm          comm;
378   const PetscInt    *localPoints;
379   const PetscSFNode *remotePoints;
380   PetscInt          lpStart, lpEnd;
381   PetscInt          numRoots, numSectionRoots, numPoints, numIndices = 0;
382   PetscInt          *localIndices;
383   PetscSFNode       *remoteIndices;
384   PetscInt          i, ind;
385   PetscErrorCode    ierr;
386 
387   PetscFunctionBegin;
388   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
389   PetscValidPointer(rootSection,2);
390   /* Cannot check PetscValidIntPointer(remoteOffsets,3) because it can be NULL if sf does not reference any points in leafSection */
391   PetscValidPointer(leafSection,4);
392   PetscValidPointer(sectionSF,5);
393   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
394   ierr = PetscSFCreate(comm, sectionSF);CHKERRQ(ierr);
395   ierr = PetscSectionGetChart(leafSection, &lpStart, &lpEnd);CHKERRQ(ierr);
396   ierr = PetscSectionGetStorageSize(rootSection, &numSectionRoots);CHKERRQ(ierr);
397   ierr = PetscSFGetGraph(sf, &numRoots, &numPoints, &localPoints, &remotePoints);CHKERRQ(ierr);
398   if (numRoots < 0) PetscFunctionReturn(0);
399   ierr = PetscLogEventBegin(PETSCSF_SectSF,sf,0,0,0);CHKERRQ(ierr);
400   for (i = 0; i < numPoints; ++i) {
401     PetscInt localPoint = localPoints ? localPoints[i] : i;
402     PetscInt dof;
403 
404     if ((localPoint >= lpStart) && (localPoint < lpEnd)) {
405       ierr = PetscSectionGetDof(leafSection, localPoint, &dof);CHKERRQ(ierr);
406       numIndices += dof;
407     }
408   }
409   ierr = PetscMalloc1(numIndices, &localIndices);CHKERRQ(ierr);
410   ierr = PetscMalloc1(numIndices, &remoteIndices);CHKERRQ(ierr);
411   /* Create new index graph */
412   for (i = 0, ind = 0; i < numPoints; ++i) {
413     PetscInt localPoint = localPoints ? localPoints[i] : i;
414     PetscInt rank       = remotePoints[i].rank;
415 
416     if ((localPoint >= lpStart) && (localPoint < lpEnd)) {
417       PetscInt remoteOffset = remoteOffsets[localPoint-lpStart];
418       PetscInt loff, dof, d;
419 
420       ierr = PetscSectionGetOffset(leafSection, localPoint, &loff);CHKERRQ(ierr);
421       ierr = PetscSectionGetDof(leafSection, localPoint, &dof);CHKERRQ(ierr);
422       for (d = 0; d < dof; ++d, ++ind) {
423         localIndices[ind]        = loff+d;
424         remoteIndices[ind].rank  = rank;
425         remoteIndices[ind].index = remoteOffset+d;
426       }
427     }
428   }
429   PetscCheckFalse(numIndices != ind,comm, PETSC_ERR_PLIB, "Inconsistency in indices, %" PetscInt_FMT " should be %" PetscInt_FMT, ind, numIndices);
430   ierr = PetscSFSetGraph(*sectionSF, numSectionRoots, numIndices, localIndices, PETSC_OWN_POINTER, remoteIndices, PETSC_OWN_POINTER);CHKERRQ(ierr);
431   ierr = PetscSFSetUp(*sectionSF);CHKERRQ(ierr);
432   ierr = PetscLogEventEnd(PETSCSF_SectSF,sf,0,0,0);CHKERRQ(ierr);
433   PetscFunctionReturn(0);
434 }
435 
436 /*@C
437    PetscSFCreateFromLayouts - Creates a parallel star forest mapping two PetscLayout objects
438 
439    Collective
440 
441    Input Parameters:
442 +  rmap - PetscLayout defining the global root space
443 -  lmap - PetscLayout defining the global leaf space
444 
445    Output Parameter:
446 .  sf - The parallel star forest
447 
448    Level: intermediate
449 
450 .seealso: PetscSFCreate(), PetscLayoutCreate(), PetscSFSetGraphLayout()
451 @*/
452 PetscErrorCode PetscSFCreateFromLayouts(PetscLayout rmap, PetscLayout lmap, PetscSF* sf)
453 {
454   PetscErrorCode ierr;
455   PetscInt       i,nroots,nleaves = 0;
456   PetscInt       rN, lst, len;
457   PetscMPIInt    owner = -1;
458   PetscSFNode    *remote;
459   MPI_Comm       rcomm = rmap->comm;
460   MPI_Comm       lcomm = lmap->comm;
461   PetscMPIInt    flg;
462 
463   PetscFunctionBegin;
464   PetscValidPointer(sf,3);
465   PetscCheckFalse(!rmap->setupcalled,rcomm,PETSC_ERR_ARG_WRONGSTATE,"Root layout not setup");
466   PetscCheckFalse(!lmap->setupcalled,lcomm,PETSC_ERR_ARG_WRONGSTATE,"Leaf layout not setup");
467   ierr = MPI_Comm_compare(rcomm,lcomm,&flg);CHKERRMPI(ierr);
468   PetscCheckFalse(flg != MPI_CONGRUENT && flg != MPI_IDENT,rcomm,PETSC_ERR_SUP,"cannot map two layouts with non-matching communicators");
469   ierr = PetscSFCreate(rcomm,sf);CHKERRQ(ierr);
470   ierr = PetscLayoutGetLocalSize(rmap,&nroots);CHKERRQ(ierr);
471   ierr = PetscLayoutGetSize(rmap,&rN);CHKERRQ(ierr);
472   ierr = PetscLayoutGetRange(lmap,&lst,&len);CHKERRQ(ierr);
473   ierr = PetscMalloc1(len-lst,&remote);CHKERRQ(ierr);
474   for (i = lst; i < len && i < rN; i++) {
475     if (owner < -1 || i >= rmap->range[owner+1]) {
476       ierr = PetscLayoutFindOwner(rmap,i,&owner);CHKERRQ(ierr);
477     }
478     remote[nleaves].rank  = owner;
479     remote[nleaves].index = i - rmap->range[owner];
480     nleaves++;
481   }
482   ierr = PetscSFSetGraph(*sf,nroots,nleaves,NULL,PETSC_OWN_POINTER,remote,PETSC_COPY_VALUES);CHKERRQ(ierr);
483   ierr = PetscFree(remote);CHKERRQ(ierr);
484   PetscFunctionReturn(0);
485 }
486 
487 /* TODO: handle nooffprocentries like MatZeroRowsMapLocal_Private, since this code is the same */
488 PetscErrorCode PetscLayoutMapLocal(PetscLayout map,PetscInt N,const PetscInt idxs[], PetscInt *on,PetscInt **oidxs,PetscInt **ogidxs)
489 {
490   PetscInt      *owners = map->range;
491   PetscInt       n      = map->n;
492   PetscSF        sf;
493   PetscInt      *lidxs,*work = NULL;
494   PetscSFNode   *ridxs;
495   PetscMPIInt    rank, p = 0;
496   PetscInt       r, len = 0;
497   PetscErrorCode ierr;
498 
499   PetscFunctionBegin;
500   if (on) *on = 0;              /* squelch -Wmaybe-uninitialized */
501   /* Create SF where leaves are input idxs and roots are owned idxs */
502   ierr = MPI_Comm_rank(map->comm,&rank);CHKERRMPI(ierr);
503   ierr = PetscMalloc1(n,&lidxs);CHKERRQ(ierr);
504   for (r = 0; r < n; ++r) lidxs[r] = -1;
505   ierr = PetscMalloc1(N,&ridxs);CHKERRQ(ierr);
506   for (r = 0; r < N; ++r) {
507     const PetscInt idx = idxs[r];
508     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);
509     if (idx < owners[p] || owners[p+1] <= idx) { /* short-circuit the search if the last p owns this idx too */
510       ierr = PetscLayoutFindOwner(map,idx,&p);CHKERRQ(ierr);
511     }
512     ridxs[r].rank = p;
513     ridxs[r].index = idxs[r] - owners[p];
514   }
515   ierr = PetscSFCreate(map->comm,&sf);CHKERRQ(ierr);
516   ierr = PetscSFSetGraph(sf,n,N,NULL,PETSC_OWN_POINTER,ridxs,PETSC_OWN_POINTER);CHKERRQ(ierr);
517   ierr = PetscSFReduceBegin(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr);
518   ierr = PetscSFReduceEnd(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr);
519   if (ogidxs) { /* communicate global idxs */
520     PetscInt cum = 0,start,*work2;
521 
522     ierr = PetscMalloc1(n,&work);CHKERRQ(ierr);
523     ierr = PetscCalloc1(N,&work2);CHKERRQ(ierr);
524     for (r = 0; r < N; ++r) if (idxs[r] >=0) cum++;
525     ierr = MPI_Scan(&cum,&start,1,MPIU_INT,MPI_SUM,map->comm);CHKERRMPI(ierr);
526     start -= cum;
527     cum = 0;
528     for (r = 0; r < N; ++r) if (idxs[r] >=0) work2[r] = start+cum++;
529     ierr = PetscSFReduceBegin(sf,MPIU_INT,work2,work,MPI_REPLACE);CHKERRQ(ierr);
530     ierr = PetscSFReduceEnd(sf,MPIU_INT,work2,work,MPI_REPLACE);CHKERRQ(ierr);
531     ierr = PetscFree(work2);CHKERRQ(ierr);
532   }
533   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
534   /* Compress and put in indices */
535   for (r = 0; r < n; ++r)
536     if (lidxs[r] >= 0) {
537       if (work) work[len] = work[r];
538       lidxs[len++] = r;
539     }
540   if (on) *on = len;
541   if (oidxs) *oidxs = lidxs;
542   if (ogidxs) *ogidxs = work;
543   PetscFunctionReturn(0);
544 }
545 
546 /*@
547   PetscSFCreateByMatchingIndices - Create SF by matching root and leaf indices
548 
549   Collective
550 
551   Input Parameters:
552 + layout - PetscLayout defining the global index space and the rank that brokers each index
553 . numRootIndices - size of rootIndices
554 . rootIndices - PetscInt array of global indices of which this process requests ownership
555 . rootLocalIndices - root local index permutation (NULL if no permutation)
556 . rootLocalOffset - offset to be added to root local indices
557 . numLeafIndices - size of leafIndices
558 . leafIndices - PetscInt array of global indices with which this process requires data associated
559 . leafLocalIndices - leaf local index permutation (NULL if no permutation)
560 - leafLocalOffset - offset to be added to leaf local indices
561 
562   Output Parameters:
563 + sfA - star forest representing the communication pattern from the layout space to the leaf space (NULL if not needed)
564 - sf - star forest representing the communication pattern from the root space to the leaf space
565 
566   Example 1:
567 $
568 $  rank             : 0            1            2
569 $  rootIndices      : [1 0 2]      [3]          [3]
570 $  rootLocalOffset  : 100          200          300
571 $  layout           : [0 1]        [2]          [3]
572 $  leafIndices      : [0]          [2]          [0 3]
573 $  leafLocalOffset  : 400          500          600
574 $
575 would build the following SF
576 $
577 $  [0] 400 <- (0,101)
578 $  [1] 500 <- (0,102)
579 $  [2] 600 <- (0,101)
580 $  [2] 601 <- (2,300)
581 $
582   Example 2:
583 $
584 $  rank             : 0               1               2
585 $  rootIndices      : [1 0 2]         [3]             [3]
586 $  rootLocalOffset  : 100             200             300
587 $  layout           : [0 1]           [2]             [3]
588 $  leafIndices      : rootIndices     rootIndices     rootIndices
589 $  leafLocalOffset  : rootLocalOffset rootLocalOffset rootLocalOffset
590 $
591 would build the following SF
592 $
593 $  [1] 200 <- (2,300)
594 $
595   Example 3:
596 $
597 $  No process requests ownership of global index 1, but no process needs it.
598 $
599 $  rank             : 0            1            2
600 $  numRootIndices   : 2            1            1
601 $  rootIndices      : [0 2]        [3]          [3]
602 $  rootLocalOffset  : 100          200          300
603 $  layout           : [0 1]        [2]          [3]
604 $  numLeafIndices   : 1            1            2
605 $  leafIndices      : [0]          [2]          [0 3]
606 $  leafLocalOffset  : 400          500          600
607 $
608 would build the following SF
609 $
610 $  [0] 400 <- (0,100)
611 $  [1] 500 <- (0,101)
612 $  [2] 600 <- (0,100)
613 $  [2] 601 <- (2,300)
614 $
615 
616   Notes:
617   The layout parameter represents any partitioning of [0, N), where N is the total number of global indices, and its
618   local size can be set to PETSC_DECIDE.
619   If a global index x lies in the partition owned by process i, each process whose rootIndices contains x requests
620   ownership of x and sends its own rank and the local index of x to process i.
621   If multiple processes request ownership of x, the one with the highest rank is to own x.
622   Process i then broadcasts the ownership information, so that each process whose leafIndices contains x knows the
623   ownership information of x.
624   The output sf is constructed by associating each leaf point to a root point in this way.
625 
626   Suppose there is point data ordered according to the global indices and partitioned according to the given layout.
627   The optional output PetscSF object sfA can be used to push such data to leaf points.
628 
629   All indices in rootIndices and leafIndices must lie in the layout range. The union (over all processes) of rootIndices
630   must cover that of leafIndices, but need not cover the entire layout.
631 
632   If (leafIndices, leafLocalIndices, leafLocalOffset) == (rootIndices, rootLocalIndices, rootLocalOffset), the output
633   star forest is almost identity, so will only include non-trivial part of the map.
634 
635   Developer Notes:
636   Current approach of a process of the highest rank gaining the ownership may cause load imbalance; consider using
637   hash(rank, root_local_index) as the bid for the ownership determination.
638 
639   Level: advanced
640 
641 .seealso: PetscSFCreate()
642 @*/
643 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)
644 {
645   MPI_Comm        comm = layout->comm;
646   PetscMPIInt     size, rank;
647   PetscSF         sf1;
648   PetscSFNode    *owners, *buffer, *iremote;
649   PetscInt       *ilocal, nleaves, N, n, i;
650 #if defined(PETSC_USE_DEBUG)
651   PetscInt        N1;
652 #endif
653   PetscBool       flag;
654   PetscErrorCode  ierr;
655 
656   PetscFunctionBegin;
657   if (rootIndices)      PetscValidIntPointer(rootIndices,3);
658   if (rootLocalIndices) PetscValidIntPointer(rootLocalIndices,4);
659   if (leafIndices)      PetscValidIntPointer(leafIndices,7);
660   if (leafLocalIndices) PetscValidIntPointer(leafLocalIndices,8);
661   if (sfA)              PetscValidPointer(sfA,10);
662   PetscValidPointer(sf,11);
663   PetscCheckFalse(numRootIndices < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numRootIndices (%" PetscInt_FMT ") must be non-negative", numRootIndices);
664   PetscCheckFalse(numLeafIndices < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numLeafIndices (%" PetscInt_FMT ") must be non-negative", numLeafIndices);
665   ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
666   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
667   ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr);
668   ierr = PetscLayoutGetSize(layout, &N);CHKERRQ(ierr);
669   ierr = PetscLayoutGetLocalSize(layout, &n);CHKERRQ(ierr);
670   flag = (PetscBool)(leafIndices == rootIndices);
671   ierr = MPIU_Allreduce(MPI_IN_PLACE, &flag, 1, MPIU_BOOL, MPI_LAND, comm);CHKERRMPI(ierr);
672   PetscCheckFalse(flag && numLeafIndices != numRootIndices,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "leafIndices == rootIndices, but numLeafIndices (%" PetscInt_FMT ") != numRootIndices(%" PetscInt_FMT ")", numLeafIndices, numRootIndices);
673 #if defined(PETSC_USE_DEBUG)
674   N1 = PETSC_MIN_INT;
675   for (i = 0; i < numRootIndices; i++) if (rootIndices[i] > N1) N1 = rootIndices[i];
676   ierr = MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm);CHKERRMPI(ierr);
677   PetscCheckFalse(N1 >= N,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Max. root index (%" PetscInt_FMT ") out of layout range [0,%" PetscInt_FMT ")", N1, N);
678   if (!flag) {
679     N1 = PETSC_MIN_INT;
680     for (i = 0; i < numLeafIndices; i++) if (leafIndices[i] > N1) N1 = leafIndices[i];
681     ierr = MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm);CHKERRMPI(ierr);
682     PetscCheckFalse(N1 >= N,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Max. leaf index (%" PetscInt_FMT ") out of layout range [0,%" PetscInt_FMT ")", N1, N);
683   }
684 #endif
685   /* Reduce: owners -> buffer */
686   ierr = PetscMalloc1(n, &buffer);CHKERRQ(ierr);
687   ierr = PetscSFCreate(comm, &sf1);CHKERRQ(ierr);
688   ierr = PetscSFSetFromOptions(sf1);CHKERRQ(ierr);
689   ierr = PetscSFSetGraphLayout(sf1, layout, numRootIndices, NULL, PETSC_OWN_POINTER, rootIndices);CHKERRQ(ierr);
690   ierr = PetscMalloc1(numRootIndices, &owners);CHKERRQ(ierr);
691   for (i = 0; i < numRootIndices; ++i) {
692     owners[i].rank = rank;
693     owners[i].index = rootLocalOffset + (rootLocalIndices ? rootLocalIndices[i] : i);
694   }
695   for (i = 0; i < n; ++i) {
696     buffer[i].index = -1;
697     buffer[i].rank = -1;
698   }
699   ierr = PetscSFReduceBegin(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC);CHKERRQ(ierr);
700   ierr = PetscSFReduceEnd(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC);CHKERRQ(ierr);
701   /* Bcast: buffer -> owners */
702   if (!flag) {
703     /* leafIndices is different from rootIndices */
704     ierr = PetscFree(owners);CHKERRQ(ierr);
705     ierr = PetscSFSetGraphLayout(sf1, layout, numLeafIndices, NULL, PETSC_OWN_POINTER, leafIndices);CHKERRQ(ierr);
706     ierr = PetscMalloc1(numLeafIndices, &owners);CHKERRQ(ierr);
707   }
708   ierr = PetscSFBcastBegin(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE);CHKERRQ(ierr);
709   ierr = PetscSFBcastEnd(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE);CHKERRQ(ierr);
710   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]);
711   ierr = PetscFree(buffer);CHKERRQ(ierr);
712   if (sfA) {*sfA = sf1;}
713   else     {ierr = PetscSFDestroy(&sf1);CHKERRQ(ierr);}
714   /* Create sf */
715   if (flag && rootLocalIndices == leafLocalIndices && leafLocalOffset == rootLocalOffset) {
716     /* leaf space == root space */
717     for (i = 0, nleaves = 0; i < numLeafIndices; ++i) if (owners[i].rank != rank) ++nleaves;
718     ierr = PetscMalloc1(nleaves, &ilocal);CHKERRQ(ierr);
719     ierr = PetscMalloc1(nleaves, &iremote);CHKERRQ(ierr);
720     for (i = 0, nleaves = 0; i < numLeafIndices; ++i) {
721       if (owners[i].rank != rank) {
722         ilocal[nleaves]        = leafLocalOffset + i;
723         iremote[nleaves].rank  = owners[i].rank;
724         iremote[nleaves].index = owners[i].index;
725         ++nleaves;
726       }
727     }
728     ierr = PetscFree(owners);CHKERRQ(ierr);
729   } else {
730     nleaves = numLeafIndices;
731     ierr = PetscMalloc1(nleaves, &ilocal);CHKERRQ(ierr);
732     for (i = 0; i < nleaves; ++i) {ilocal[i] = leafLocalOffset + (leafLocalIndices ? leafLocalIndices[i] : i);}
733     iremote = owners;
734   }
735   ierr = PetscSFCreate(comm, sf);CHKERRQ(ierr);
736   ierr = PetscSFSetFromOptions(*sf);CHKERRQ(ierr);
737   ierr = PetscSFSetGraph(*sf, rootLocalOffset + numRootIndices, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr);
738   PetscFunctionReturn(0);
739 }
740