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