xref: /petsc/src/vec/is/sf/interface/sf.c (revision 2f6eced2a19e978d64f27de66fbc6c26cc5c7934)
1 #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/
2 #include <petscctable.h>
3 
4 /* Logging support */
5 PetscLogEvent PETSCSF_SetGraph, PETSCSF_BcastBegin, PETSCSF_BcastEnd, PETSCSF_ReduceBegin, PETSCSF_ReduceEnd, PETSCSF_FetchAndOpBegin, PETSCSF_FetchAndOpEnd;
6 
7 #if defined(PETSC_USE_DEBUG)
8 #  define PetscSFCheckGraphSet(sf,arg) do {                          \
9     if (PetscUnlikely(!(sf)->graphset))                              \
10       SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetGraph() on argument %D \"%s\" before %s()",(arg),#sf,PETSC_FUNCTION_NAME); \
11   } while (0)
12 #else
13 #  define PetscSFCheckGraphSet(sf,arg) do {} while (0)
14 #endif
15 
16 const char *const PetscSFDuplicateOptions[] = {"CONFONLY","RANKS","GRAPH","PetscSFDuplicateOption","PETSCSF_DUPLICATE_",0};
17 
18 /*@C
19    PetscSFCreate - create a star forest communication context
20 
21    Not Collective
22 
23    Input Arguments:
24 .  comm - communicator on which the star forest will operate
25 
26    Output Arguments:
27 .  sf - new star forest context
28 
29    Level: intermediate
30 
31 .seealso: PetscSFSetGraph(), PetscSFDestroy()
32 @*/
33 PetscErrorCode PetscSFCreate(MPI_Comm comm,PetscSF *sf)
34 {
35   PetscErrorCode ierr;
36   PetscSF        b;
37 
38   PetscFunctionBegin;
39   PetscValidPointer(sf,2);
40   ierr = PetscSFInitializePackage();CHKERRQ(ierr);
41 
42   ierr = PetscHeaderCreate(b,PETSCSF_CLASSID,"PetscSF","Star Forest","PetscSF",comm,PetscSFDestroy,PetscSFView);CHKERRQ(ierr);
43 
44   b->nroots    = -1;
45   b->nleaves   = -1;
46   b->nranks    = -1;
47   b->rankorder = PETSC_TRUE;
48   b->ingroup   = MPI_GROUP_NULL;
49   b->outgroup  = MPI_GROUP_NULL;
50   b->graphset  = PETSC_FALSE;
51 
52   *sf = b;
53   PetscFunctionReturn(0);
54 }
55 
56 /*@C
57    PetscSFReset - Reset a star forest so that different sizes or neighbors can be used
58 
59    Collective
60 
61    Input Arguments:
62 .  sf - star forest
63 
64    Level: advanced
65 
66 .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFDestroy()
67 @*/
68 PetscErrorCode PetscSFReset(PetscSF sf)
69 {
70   PetscErrorCode ierr;
71 
72   PetscFunctionBegin;
73   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
74   sf->mine   = NULL;
75   ierr       = PetscFree(sf->mine_alloc);CHKERRQ(ierr);
76   sf->remote = NULL;
77   ierr       = PetscFree(sf->remote_alloc);CHKERRQ(ierr);
78   ierr       = PetscFree4(sf->ranks,sf->roffset,sf->rmine,sf->rremote);CHKERRQ(ierr);
79   ierr       = PetscFree(sf->degree);CHKERRQ(ierr);
80   if (sf->ingroup  != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->ingroup);CHKERRQ(ierr);}
81   if (sf->outgroup != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->outgroup);CHKERRQ(ierr);}
82   ierr         = PetscSFDestroy(&sf->multi);CHKERRQ(ierr);
83   sf->graphset = PETSC_FALSE;
84   if (sf->ops->Reset) {ierr = (*sf->ops->Reset)(sf);CHKERRQ(ierr);}
85   sf->setupcalled = PETSC_FALSE;
86   PetscFunctionReturn(0);
87 }
88 
89 /*@C
90    PetscSFSetType - set the PetscSF communication implementation
91 
92    Collective on PetscSF
93 
94    Input Parameters:
95 +  sf - the PetscSF context
96 -  type - a known method
97 
98    Options Database Key:
99 .  -sf_type <type> - Sets the method; use -help for a list
100    of available methods (for instance, window, pt2pt, neighbor)
101 
102    Notes:
103    See "include/petscsf.h" for available methods (for instance)
104 +    PETSCSFWINDOW - MPI-2/3 one-sided
105 -    PETSCSFBASIC - basic implementation using MPI-1 two-sided
106 
107   Level: intermediate
108 
109 .keywords: PetscSF, set, type
110 
111 .seealso: PetscSFType, PetscSFCreate()
112 @*/
113 PetscErrorCode PetscSFSetType(PetscSF sf,PetscSFType type)
114 {
115   PetscErrorCode ierr,(*r)(PetscSF);
116   PetscBool      match;
117 
118   PetscFunctionBegin;
119   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
120   PetscValidCharPointer(type,2);
121 
122   ierr = PetscObjectTypeCompare((PetscObject)sf,type,&match);CHKERRQ(ierr);
123   if (match) PetscFunctionReturn(0);
124 
125   ierr = PetscFunctionListFind(PetscSFList,type,&r);CHKERRQ(ierr);
126   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested PetscSF type %s",type);
127   /* Destroy the previous private PetscSF context */
128   if (sf->ops->Destroy) {
129     ierr = (*(sf)->ops->Destroy)(sf);CHKERRQ(ierr);
130   }
131   ierr = PetscMemzero(sf->ops,sizeof(*sf->ops));CHKERRQ(ierr);
132   ierr = PetscObjectChangeTypeName((PetscObject)sf,type);CHKERRQ(ierr);
133   ierr = (*r)(sf);CHKERRQ(ierr);
134   PetscFunctionReturn(0);
135 }
136 
137 /*@
138    PetscSFDestroy - destroy star forest
139 
140    Collective
141 
142    Input Arguments:
143 .  sf - address of star forest
144 
145    Level: intermediate
146 
147 .seealso: PetscSFCreate(), PetscSFReset()
148 @*/
149 PetscErrorCode PetscSFDestroy(PetscSF *sf)
150 {
151   PetscErrorCode ierr;
152 
153   PetscFunctionBegin;
154   if (!*sf) PetscFunctionReturn(0);
155   PetscValidHeaderSpecific((*sf),PETSCSF_CLASSID,1);
156   if (--((PetscObject)(*sf))->refct > 0) {*sf = 0; PetscFunctionReturn(0);}
157   ierr = PetscSFReset(*sf);CHKERRQ(ierr);
158   if ((*sf)->ops->Destroy) {ierr = (*(*sf)->ops->Destroy)(*sf);CHKERRQ(ierr);}
159   ierr = PetscHeaderDestroy(sf);CHKERRQ(ierr);
160   PetscFunctionReturn(0);
161 }
162 
163 /*@
164    PetscSFSetUp - set up communication structures
165 
166    Collective
167 
168    Input Arguments:
169 .  sf - star forest communication object
170 
171    Level: beginner
172 
173 .seealso: PetscSFSetFromOptions(), PetscSFSetType()
174 @*/
175 PetscErrorCode PetscSFSetUp(PetscSF sf)
176 {
177   PetscErrorCode ierr;
178 
179   PetscFunctionBegin;
180   if (sf->setupcalled) PetscFunctionReturn(0);
181   if (!((PetscObject)sf)->type_name) {ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr);}
182   if (sf->ops->SetUp) {ierr = (*sf->ops->SetUp)(sf);CHKERRQ(ierr);}
183   sf->setupcalled = PETSC_TRUE;
184   PetscFunctionReturn(0);
185 }
186 
187 /*@C
188    PetscSFSetFromOptions - set PetscSF options using the options database
189 
190    Logically Collective
191 
192    Input Arguments:
193 .  sf - star forest
194 
195    Options Database Keys:
196 +  -sf_type - implementation type, see PetscSFSetType()
197 -  -sf_rank_order - sort composite points for gathers and scatters in rank order, gathers are non-deterministic otherwise
198 
199    Level: intermediate
200 
201 .keywords: KSP, set, from, options, database
202 
203 .seealso: PetscSFWindowSetSyncType()
204 @*/
205 PetscErrorCode PetscSFSetFromOptions(PetscSF sf)
206 {
207   PetscSFType    deft;
208   char           type[256];
209   PetscErrorCode ierr;
210   PetscBool      flg;
211 
212   PetscFunctionBegin;
213   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
214   ierr = PetscObjectOptionsBegin((PetscObject)sf);CHKERRQ(ierr);
215   deft = ((PetscObject)sf)->type_name ? ((PetscObject)sf)->type_name : PETSCSFBASIC;
216   ierr = PetscOptionsFList("-sf_type","PetscSF implementation type","PetscSFSetType",PetscSFList,deft,type,256,&flg);CHKERRQ(ierr);
217   ierr = PetscSFSetType(sf,flg ? type : deft);CHKERRQ(ierr);
218   ierr = PetscOptionsBool("-sf_rank_order","sort composite points for gathers and scatters in rank order, gathers are non-deterministic otherwise","PetscSFSetRankOrder",sf->rankorder,&sf->rankorder,NULL);CHKERRQ(ierr);
219   if (sf->ops->SetFromOptions) {ierr = (*sf->ops->SetFromOptions)(PetscOptionsObject,sf);CHKERRQ(ierr);}
220   ierr = PetscOptionsEnd();CHKERRQ(ierr);
221   PetscFunctionReturn(0);
222 }
223 
224 /*@C
225    PetscSFSetRankOrder - sort multi-points for gathers and scatters by rank order
226 
227    Logically Collective
228 
229    Input Arguments:
230 +  sf - star forest
231 -  flg - PETSC_TRUE to sort, PETSC_FALSE to skip sorting (lower setup cost, but non-deterministic)
232 
233    Level: advanced
234 
235 .seealso: PetscSFGatherBegin(), PetscSFScatterBegin()
236 @*/
237 PetscErrorCode PetscSFSetRankOrder(PetscSF sf,PetscBool flg)
238 {
239 
240   PetscFunctionBegin;
241   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
242   PetscValidLogicalCollectiveBool(sf,flg,2);
243   if (sf->multi) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Rank ordering must be set before first call to PetscSFGatherBegin() or PetscSFScatterBegin()");
244   sf->rankorder = flg;
245   PetscFunctionReturn(0);
246 }
247 
248 /*@C
249    PetscSFSetGraph - Set a parallel star forest
250 
251    Collective
252 
253    Input Arguments:
254 +  sf - star forest
255 .  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
256 .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
257 .  ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage
258 .  localmode - copy mode for ilocal
259 .  iremote - remote locations of root vertices for each leaf on the current process
260 -  remotemode - copy mode for iremote
261 
262    Level: intermediate
263 
264 .seealso: PetscSFCreate(), PetscSFView(), PetscSFGetGraph()
265 @*/
266 PetscErrorCode PetscSFSetGraph(PetscSF sf,PetscInt nroots,PetscInt nleaves,const PetscInt *ilocal,PetscCopyMode localmode,const PetscSFNode *iremote,PetscCopyMode remotemode)
267 {
268   PetscErrorCode     ierr;
269   PetscTable         table;
270   PetscTablePosition pos;
271   PetscMPIInt        size;
272   PetscInt           i,*rcount,*ranks;
273 
274   PetscFunctionBegin;
275   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
276   ierr = PetscLogEventBegin(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr);
277   if (nleaves && ilocal) PetscValidIntPointer(ilocal,4);
278   if (nleaves) PetscValidPointer(iremote,6);
279   if (nroots < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"roots %D, cannot be negative",nroots);
280   if (nleaves < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nleaves %D, cannot be negative",nleaves);
281   ierr        = PetscSFReset(sf);CHKERRQ(ierr);
282   sf->nroots  = nroots;
283   sf->nleaves = nleaves;
284   if (ilocal) {
285     switch (localmode) {
286     case PETSC_COPY_VALUES:
287       ierr        = PetscMalloc1(nleaves,&sf->mine_alloc);CHKERRQ(ierr);
288       sf->mine    = sf->mine_alloc;
289       ierr        = PetscMemcpy(sf->mine,ilocal,nleaves*sizeof(*sf->mine));CHKERRQ(ierr);
290       sf->minleaf = PETSC_MAX_INT;
291       sf->maxleaf = PETSC_MIN_INT;
292       for (i=0; i<nleaves; i++) {
293         sf->minleaf = PetscMin(sf->minleaf,ilocal[i]);
294         sf->maxleaf = PetscMax(sf->maxleaf,ilocal[i]);
295       }
296       break;
297     case PETSC_OWN_POINTER:
298       sf->mine_alloc = (PetscInt*)ilocal;
299       sf->mine       = sf->mine_alloc;
300       break;
301     case PETSC_USE_POINTER:
302       sf->mine = (PetscInt*)ilocal;
303       break;
304     default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown localmode");
305     }
306   }
307   if (!ilocal || nleaves > 0) {
308     sf->minleaf = 0;
309     sf->maxleaf = nleaves - 1;
310   }
311   switch (remotemode) {
312   case PETSC_COPY_VALUES:
313     ierr       = PetscMalloc1(nleaves,&sf->remote_alloc);CHKERRQ(ierr);
314     sf->remote = sf->remote_alloc;
315     ierr       = PetscMemcpy(sf->remote,iremote,nleaves*sizeof(*sf->remote));CHKERRQ(ierr);
316     break;
317   case PETSC_OWN_POINTER:
318     sf->remote_alloc = (PetscSFNode*)iremote;
319     sf->remote       = sf->remote_alloc;
320     break;
321   case PETSC_USE_POINTER:
322     sf->remote = (PetscSFNode*)iremote;
323     break;
324   default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown remotemode");
325   }
326 
327   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr);
328   ierr = PetscTableCreate(10,size,&table);CHKERRQ(ierr);
329   for (i=0; i<nleaves; i++) {
330     /* Log 1-based rank */
331     ierr = PetscTableAdd(table,iremote[i].rank+1,1,ADD_VALUES);CHKERRQ(ierr);
332   }
333   ierr = PetscTableGetCount(table,&sf->nranks);CHKERRQ(ierr);
334   ierr = PetscMalloc4(sf->nranks,&sf->ranks,sf->nranks+1,&sf->roffset,nleaves,&sf->rmine,nleaves,&sf->rremote);CHKERRQ(ierr);
335   ierr = PetscMalloc2(sf->nranks,&rcount,sf->nranks,&ranks);CHKERRQ(ierr);
336   ierr = PetscTableGetHeadPosition(table,&pos);CHKERRQ(ierr);
337   for (i=0; i<sf->nranks; i++) {
338     ierr = PetscTableGetNext(table,&pos,&ranks[i],&rcount[i]);CHKERRQ(ierr);
339     ranks[i]--;             /* Convert back to 0-based */
340   }
341   ierr = PetscTableDestroy(&table);CHKERRQ(ierr);
342   ierr = PetscSortIntWithArray(sf->nranks,ranks,rcount);CHKERRQ(ierr);
343   sf->roffset[0] = 0;
344   for (i=0; i<sf->nranks; i++) {
345     ierr = PetscMPIIntCast(ranks[i],sf->ranks+i);CHKERRQ(ierr);
346     sf->roffset[i+1] = sf->roffset[i] + rcount[i];
347     rcount[i]        = 0;
348   }
349   for (i=0; i<nleaves; i++) {
350     PetscInt lo,hi,irank;
351     /* Search for index of iremote[i].rank in sf->ranks */
352     lo = 0; hi = sf->nranks;
353     while (hi - lo > 1) {
354       PetscInt mid = lo + (hi - lo)/2;
355       if (iremote[i].rank < sf->ranks[mid]) hi = mid;
356       else                                  lo = mid;
357     }
358     if (hi - lo == 1 && iremote[i].rank == sf->ranks[lo]) irank = lo;
359     else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Could not find rank %D in array",iremote[i].rank);
360     sf->rmine[sf->roffset[irank] + rcount[irank]]   = ilocal ? ilocal[i] : i;
361     sf->rremote[sf->roffset[irank] + rcount[irank]] = iremote[i].index;
362     rcount[irank]++;
363   }
364   ierr = PetscFree2(rcount,ranks);CHKERRQ(ierr);
365 #if !defined(PETSC_USE_64BIT_INDICES)
366   if (nroots == PETSC_DETERMINE) {
367     /* Jed, if you have a better way to do this, put it in */
368     PetscInt *numRankLeaves, *leafOff, *leafIndices, *numRankRoots, *rootOff, *rootIndices, maxRoots = 0;
369 
370     /* All to all to determine number of leaf indices from each (you can do this using Scan and asynch messages) */
371     ierr = PetscMalloc4(size,&numRankLeaves,size+1,&leafOff,size,&numRankRoots,size+1,&rootOff);CHKERRQ(ierr);
372     ierr = PetscMemzero(numRankLeaves, size * sizeof(PetscInt));CHKERRQ(ierr);
373     for (i = 0; i < nleaves; ++i) ++numRankLeaves[iremote[i].rank];
374     ierr = MPI_Alltoall(numRankLeaves, 1, MPIU_INT, numRankRoots, 1, MPIU_INT, PetscObjectComm((PetscObject)sf));CHKERRQ(ierr);
375     /* Could set nroots to this maximum */
376     for (i = 0; i < size; ++i) maxRoots += numRankRoots[i];
377 
378     /* Gather all indices */
379     ierr = PetscMalloc2(nleaves,&leafIndices,maxRoots,&rootIndices);CHKERRQ(ierr);
380     leafOff[0] = 0;
381     for (i = 0; i < size; ++i) leafOff[i+1] = leafOff[i] + numRankLeaves[i];
382     for (i = 0; i < nleaves; ++i) leafIndices[leafOff[iremote[i].rank]++] = iremote[i].index;
383     leafOff[0] = 0;
384     for (i = 0; i < size; ++i) leafOff[i+1] = leafOff[i] + numRankLeaves[i];
385     rootOff[0] = 0;
386     for (i = 0; i < size; ++i) rootOff[i+1] = rootOff[i] + numRankRoots[i];
387     ierr = MPI_Alltoallv(leafIndices, numRankLeaves, leafOff, MPIU_INT, rootIndices, numRankRoots, rootOff, MPIU_INT, PetscObjectComm((PetscObject)sf));CHKERRQ(ierr);
388     /* Sort and reduce */
389     ierr       = PetscSortRemoveDupsInt(&maxRoots, rootIndices);CHKERRQ(ierr);
390     ierr       = PetscFree2(leafIndices,rootIndices);CHKERRQ(ierr);
391     ierr       = PetscFree4(numRankLeaves,leafOff,numRankRoots,rootOff);CHKERRQ(ierr);
392     sf->nroots = maxRoots;
393   }
394 #endif
395 
396   sf->graphset = PETSC_TRUE;
397   ierr = PetscLogEventEnd(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr);
398   PetscFunctionReturn(0);
399 }
400 
401 /*@C
402    PetscSFCreateInverseSF - given a PetscSF in which all vertices have degree 1, creates the inverse map
403 
404    Collective
405 
406    Input Arguments:
407 .  sf - star forest to invert
408 
409    Output Arguments:
410 .  isf - inverse of sf
411 
412    Level: advanced
413 
414    Notes:
415    All roots must have degree 1.
416 
417    The local space may be a permutation, but cannot be sparse.
418 
419 .seealso: PetscSFSetGraph()
420 @*/
421 PetscErrorCode PetscSFCreateInverseSF(PetscSF sf,PetscSF *isf)
422 {
423   PetscErrorCode ierr;
424   PetscMPIInt    rank;
425   PetscInt       i,nroots,nleaves,maxlocal,count,*newilocal;
426   const PetscInt *ilocal;
427   PetscSFNode    *roots,*leaves;
428 
429   PetscFunctionBegin;
430   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
431   ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);CHKERRQ(ierr);
432   for (i=0,maxlocal=0; i<nleaves; i++) maxlocal = PetscMax(maxlocal,(ilocal ? ilocal[i] : i)+1);
433   ierr = PetscMalloc2(nroots,&roots,maxlocal,&leaves);CHKERRQ(ierr);
434   for (i=0; i<maxlocal; i++) {
435     leaves[i].rank  = rank;
436     leaves[i].index = i;
437   }
438   for (i=0; i <nroots; i++) {
439     roots[i].rank  = -1;
440     roots[i].index = -1;
441   }
442   ierr = PetscSFReduceBegin(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr);
443   ierr = PetscSFReduceEnd(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr);
444 
445   /* Check whether our leaves are sparse */
446   for (i=0,count=0; i<nroots; i++) if (roots[i].rank >= 0) count++;
447   if (count == nroots) newilocal = NULL;
448   else {                        /* Index for sparse leaves and compact "roots" array (which is to become our leaves). */
449     ierr = PetscMalloc1(count,&newilocal);CHKERRQ(ierr);
450     for (i=0,count=0; i<nroots; i++) {
451       if (roots[i].rank >= 0) {
452         newilocal[count]   = i;
453         roots[count].rank  = roots[i].rank;
454         roots[count].index = roots[i].index;
455         count++;
456       }
457     }
458   }
459 
460   ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,isf);CHKERRQ(ierr);
461   ierr = PetscSFSetGraph(*isf,maxlocal,count,newilocal,PETSC_OWN_POINTER,roots,PETSC_COPY_VALUES);CHKERRQ(ierr);
462   ierr = PetscFree2(roots,leaves);CHKERRQ(ierr);
463   PetscFunctionReturn(0);
464 }
465 
466 /*@
467    PetscSFDuplicate - duplicate a PetscSF, optionally preserving rank connectivity and graph
468 
469    Collective
470 
471    Input Arguments:
472 +  sf - communication object to duplicate
473 -  opt - PETSCSF_DUPLICATE_CONFONLY, PETSCSF_DUPLICATE_RANKS, or PETSCSF_DUPLICATE_GRAPH (see PetscSFDuplicateOption)
474 
475    Output Arguments:
476 .  newsf - new communication object
477 
478    Level: beginner
479 
480 .seealso: PetscSFCreate(), PetscSFSetType(), PetscSFSetGraph()
481 @*/
482 PetscErrorCode PetscSFDuplicate(PetscSF sf,PetscSFDuplicateOption opt,PetscSF *newsf)
483 {
484   PetscErrorCode ierr;
485 
486   PetscFunctionBegin;
487   ierr = PetscSFCreate(PetscObjectComm((PetscObject)sf),newsf);CHKERRQ(ierr);
488   ierr = PetscSFSetType(*newsf,((PetscObject)sf)->type_name);CHKERRQ(ierr);
489   if (sf->ops->Duplicate) {ierr = (*sf->ops->Duplicate)(sf,opt,*newsf);CHKERRQ(ierr);}
490   if (opt == PETSCSF_DUPLICATE_GRAPH) {
491     PetscInt          nroots,nleaves;
492     const PetscInt    *ilocal;
493     const PetscSFNode *iremote;
494     ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
495     ierr = PetscSFSetGraph(*newsf,nroots,nleaves,ilocal,PETSC_COPY_VALUES,iremote,PETSC_COPY_VALUES);CHKERRQ(ierr);
496   }
497   PetscFunctionReturn(0);
498 }
499 
500 /*@C
501    PetscSFGetGraph - Get the graph specifying a parallel star forest
502 
503    Not Collective
504 
505    Input Arguments:
506 .  sf - star forest
507 
508    Output Arguments:
509 +  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
510 .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
511 .  ilocal - locations of leaves in leafdata buffers
512 -  iremote - remote locations of root vertices for each leaf on the current process
513 
514    Level: intermediate
515 
516 .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph()
517 @*/
518 PetscErrorCode PetscSFGetGraph(PetscSF sf,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote)
519 {
520 
521   PetscFunctionBegin;
522   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
523   /* We are not currently requiring that the graph is set, thus returning nroots=-1 if it has not been set */
524   /* if (!sf->graphset) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Graph has not been set, must call PetscSFSetGraph()"); */
525   if (nroots) *nroots = sf->nroots;
526   if (nleaves) *nleaves = sf->nleaves;
527   if (ilocal) *ilocal = sf->mine;
528   if (iremote) *iremote = sf->remote;
529   PetscFunctionReturn(0);
530 }
531 
532 /*@C
533    PetscSFGetLeafRange - Get the active leaf ranges
534 
535    Not Collective
536 
537    Input Arguments:
538 .  sf - star forest
539 
540    Output Arguments:
541 +  minleaf - minimum active leaf on this process
542 -  maxleaf - maximum active leaf on this process
543 
544    Level: developer
545 
546 .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph(), PetscSFGetGraph()
547 @*/
548 PetscErrorCode PetscSFGetLeafRange(PetscSF sf,PetscInt *minleaf,PetscInt *maxleaf)
549 {
550 
551   PetscFunctionBegin;
552   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
553   if (minleaf) *minleaf = sf->minleaf;
554   if (maxleaf) *maxleaf = sf->maxleaf;
555   PetscFunctionReturn(0);
556 }
557 
558 /*@C
559    PetscSFView - view a star forest
560 
561    Collective
562 
563    Input Arguments:
564 +  sf - star forest
565 -  viewer - viewer to display graph, for example PETSC_VIEWER_STDOUT_WORLD
566 
567    Level: beginner
568 
569 .seealso: PetscSFCreate(), PetscSFSetGraph()
570 @*/
571 PetscErrorCode PetscSFView(PetscSF sf,PetscViewer viewer)
572 {
573   PetscErrorCode    ierr;
574   PetscBool         iascii;
575   PetscViewerFormat format;
576 
577   PetscFunctionBegin;
578   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
579   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sf),&viewer);CHKERRQ(ierr);}
580   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
581   PetscCheckSameComm(sf,1,viewer,2);
582   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
583   if (iascii) {
584     PetscMPIInt rank;
585     PetscInt    i,j;
586 
587     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)sf,viewer);CHKERRQ(ierr);
588     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
589     if (sf->ops->View) {ierr = (*sf->ops->View)(sf,viewer);CHKERRQ(ierr);}
590     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
591     ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
592     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of roots=%D, leaves=%D, remote ranks=%D\n",rank,sf->nroots,sf->nleaves,sf->nranks);CHKERRQ(ierr);
593     for (i=0; i<sf->nleaves; i++) {
594       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D <- (%D,%D)\n",rank,sf->mine ? sf->mine[i] : i,sf->remote[i].rank,sf->remote[i].index);CHKERRQ(ierr);
595     }
596     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
597     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
598     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
599       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Roots referenced by my leaves, by rank\n",rank);CHKERRQ(ierr);
600       for (i=0; i<sf->nranks; i++) {
601         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %d: %D edges\n",rank,sf->ranks[i],sf->roffset[i+1]-sf->roffset[i]);CHKERRQ(ierr);
602         for (j=sf->roffset[i]; j<sf->roffset[i+1]; j++) {
603           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d]    %D <- %D\n",rank,sf->rmine[j],sf->rremote[j]);CHKERRQ(ierr);
604         }
605       }
606     }
607     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
608     ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
609     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
610   }
611   PetscFunctionReturn(0);
612 }
613 
614 /*@C
615    PetscSFGetRanks - Get ranks and number of vertices referenced by leaves on this process
616 
617    Not Collective
618 
619    Input Arguments:
620 .  sf - star forest
621 
622    Output Arguments:
623 +  nranks - number of ranks referenced by local part
624 .  ranks - array of ranks
625 .  roffset - offset in rmine/rremote for each rank (length nranks+1)
626 .  rmine - concatenated array holding local indices referencing each remote rank
627 -  rremote - concatenated array holding remote indices referenced for each remote rank
628 
629    Level: developer
630 
631 .seealso: PetscSFSetGraph()
632 @*/
633 PetscErrorCode PetscSFGetRanks(PetscSF sf,PetscInt *nranks,const PetscMPIInt **ranks,const PetscInt **roffset,const PetscInt **rmine,const PetscInt **rremote)
634 {
635 
636   PetscFunctionBegin;
637   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
638   if (nranks)  *nranks  = sf->nranks;
639   if (ranks)   *ranks   = sf->ranks;
640   if (roffset) *roffset = sf->roffset;
641   if (rmine)   *rmine   = sf->rmine;
642   if (rremote) *rremote = sf->rremote;
643   PetscFunctionReturn(0);
644 }
645 
646 /*@C
647    PetscSFGetGroups - gets incoming and outgoing process groups
648 
649    Collective
650 
651    Input Argument:
652 .  sf - star forest
653 
654    Output Arguments:
655 +  incoming - group of origin processes for incoming edges (leaves that reference my roots)
656 -  outgoing - group of destination processes for outgoing edges (roots that I reference)
657 
658    Level: developer
659 
660 .seealso: PetscSFGetWindow(), PetscSFRestoreWindow()
661 @*/
662 PetscErrorCode PetscSFGetGroups(PetscSF sf,MPI_Group *incoming,MPI_Group *outgoing)
663 {
664   PetscErrorCode ierr;
665   MPI_Group      group;
666 
667   PetscFunctionBegin;
668   if (sf->ingroup == MPI_GROUP_NULL) {
669     PetscInt       i;
670     const PetscInt *indegree;
671     PetscMPIInt    rank,*outranks,*inranks;
672     PetscSFNode    *remote;
673     PetscSF        bgcount;
674 
675     /* Compute the number of incoming ranks */
676     ierr = PetscMalloc1(sf->nranks,&remote);CHKERRQ(ierr);
677     for (i=0; i<sf->nranks; i++) {
678       remote[i].rank  = sf->ranks[i];
679       remote[i].index = 0;
680     }
681     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,&bgcount);CHKERRQ(ierr);
682     ierr = PetscSFSetGraph(bgcount,1,sf->nranks,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
683     ierr = PetscSFComputeDegreeBegin(bgcount,&indegree);CHKERRQ(ierr);
684     ierr = PetscSFComputeDegreeEnd(bgcount,&indegree);CHKERRQ(ierr);
685 
686     /* Enumerate the incoming ranks */
687     ierr = PetscMalloc2(indegree[0],&inranks,sf->nranks,&outranks);CHKERRQ(ierr);
688     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
689     for (i=0; i<sf->nranks; i++) outranks[i] = rank;
690     ierr = PetscSFGatherBegin(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr);
691     ierr = PetscSFGatherEnd(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr);
692     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
693     ierr = MPI_Group_incl(group,indegree[0],inranks,&sf->ingroup);CHKERRQ(ierr);
694     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
695     ierr = PetscFree2(inranks,outranks);CHKERRQ(ierr);
696     ierr = PetscSFDestroy(&bgcount);CHKERRQ(ierr);
697   }
698   *incoming = sf->ingroup;
699 
700   if (sf->outgroup == MPI_GROUP_NULL) {
701     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
702     ierr = MPI_Group_incl(group,sf->nranks,sf->ranks,&sf->outgroup);CHKERRQ(ierr);
703     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
704   }
705   *outgoing = sf->outgroup;
706   PetscFunctionReturn(0);
707 }
708 
709 /*@C
710    PetscSFGetMultiSF - gets the inner SF implemeting gathers and scatters
711 
712    Collective
713 
714    Input Argument:
715 .  sf - star forest that may contain roots with 0 or with more than 1 vertex
716 
717    Output Arguments:
718 .  multi - star forest with split roots, such that each root has degree exactly 1
719 
720    Level: developer
721 
722    Notes:
723 
724    In most cases, users should use PetscSFGatherBegin() and PetscSFScatterBegin() instead of manipulating multi
725    directly. Since multi satisfies the stronger condition that each entry in the global space has exactly one incoming
726    edge, it is a candidate for future optimization that might involve its removal.
727 
728 .seealso: PetscSFSetGraph(), PetscSFGatherBegin(), PetscSFScatterBegin()
729 @*/
730 PetscErrorCode PetscSFGetMultiSF(PetscSF sf,PetscSF *multi)
731 {
732   PetscErrorCode ierr;
733 
734   PetscFunctionBegin;
735   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
736   PetscValidPointer(multi,2);
737   if (sf->nroots < 0) {         /* Graph has not been set yet; why do we need this? */
738     ierr   = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr);
739     *multi = sf->multi;
740     PetscFunctionReturn(0);
741   }
742   if (!sf->multi) {
743     const PetscInt *indegree;
744     PetscInt       i,*inoffset,*outones,*outoffset,maxlocal;
745     PetscSFNode    *remote;
746     ierr        = PetscSFComputeDegreeBegin(sf,&indegree);CHKERRQ(ierr);
747     ierr        = PetscSFComputeDegreeEnd(sf,&indegree);CHKERRQ(ierr);
748     for (i=0,maxlocal=0; i<sf->nleaves; i++) maxlocal = PetscMax(maxlocal,(sf->mine ? sf->mine[i] : i)+1);
749     ierr        = PetscMalloc3(sf->nroots+1,&inoffset,maxlocal,&outones,maxlocal,&outoffset);CHKERRQ(ierr);
750     inoffset[0] = 0;
751     for (i=0; i<sf->nroots; i++) inoffset[i+1] = inoffset[i] + indegree[i];
752     for (i=0; i<maxlocal; i++) outones[i] = 1;
753     ierr = PetscSFFetchAndOpBegin(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr);
754     ierr = PetscSFFetchAndOpEnd(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr);
755     for (i=0; i<sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */
756 #if 0
757 #if defined(PETSC_USE_DEBUG)                                 /* Check that the expected number of increments occurred */
758     for (i=0; i<sf->nroots; i++) {
759       if (inoffset[i] + indegree[i] != inoffset[i+1]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect result after PetscSFFetchAndOp");
760     }
761 #endif
762 #endif
763     ierr = PetscMalloc1(sf->nleaves,&remote);CHKERRQ(ierr);
764     for (i=0; i<sf->nleaves; i++) {
765       remote[i].rank  = sf->remote[i].rank;
766       remote[i].index = outoffset[sf->mine ? sf->mine[i] : i];
767     }
768     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr);
769     ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
770     if (sf->rankorder) {        /* Sort the ranks */
771       PetscMPIInt rank;
772       PetscInt    *inranks,*newoffset,*outranks,*newoutoffset,*tmpoffset,maxdegree;
773       PetscSFNode *newremote;
774       ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
775       for (i=0,maxdegree=0; i<sf->nroots; i++) maxdegree = PetscMax(maxdegree,indegree[i]);
776       ierr = PetscMalloc5(sf->multi->nroots,&inranks,sf->multi->nroots,&newoffset,maxlocal,&outranks,maxlocal,&newoutoffset,maxdegree,&tmpoffset);CHKERRQ(ierr);
777       for (i=0; i<maxlocal; i++) outranks[i] = rank;
778       ierr = PetscSFReduceBegin(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr);
779       ierr = PetscSFReduceEnd(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr);
780       /* Sort the incoming ranks at each vertex, build the inverse map */
781       for (i=0; i<sf->nroots; i++) {
782         PetscInt j;
783         for (j=0; j<indegree[i]; j++) tmpoffset[j] = j;
784         ierr = PetscSortIntWithArray(indegree[i],inranks+inoffset[i],tmpoffset);CHKERRQ(ierr);
785         for (j=0; j<indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j;
786       }
787       ierr = PetscSFBcastBegin(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr);
788       ierr = PetscSFBcastEnd(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr);
789       ierr = PetscMalloc1(sf->nleaves,&newremote);CHKERRQ(ierr);
790       for (i=0; i<sf->nleaves; i++) {
791         newremote[i].rank  = sf->remote[i].rank;
792         newremote[i].index = newoutoffset[sf->mine ? sf->mine[i] : i];
793       }
794       ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,newremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
795       ierr = PetscFree5(inranks,newoffset,outranks,newoutoffset,tmpoffset);CHKERRQ(ierr);
796     }
797     ierr = PetscFree3(inoffset,outones,outoffset);CHKERRQ(ierr);
798   }
799   *multi = sf->multi;
800   PetscFunctionReturn(0);
801 }
802 
803 /*@C
804    PetscSFCreateEmbeddedSF - removes edges from all but the selected roots, does not remap indices
805 
806    Collective
807 
808    Input Arguments:
809 +  sf - original star forest
810 .  nroots - number of roots to select on this process
811 -  selected - selected roots on this process
812 
813    Output Arguments:
814 .  newsf - new star forest
815 
816    Level: advanced
817 
818    Note:
819    To use the new PetscSF, it may be necessary to know the indices of the leaves that are still participating. This can
820    be done by calling PetscSFGetGraph().
821 
822 .seealso: PetscSFSetGraph(), PetscSFGetGraph()
823 @*/
824 PetscErrorCode PetscSFCreateEmbeddedSF(PetscSF sf,PetscInt nroots,const PetscInt *selected,PetscSF *newsf)
825 {
826   PetscInt      *rootdata, *leafdata, *ilocal;
827   PetscSFNode   *iremote;
828   PetscInt       leafsize = 0, nleaves = 0, n, i;
829   PetscErrorCode ierr;
830 
831   PetscFunctionBegin;
832   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
833   if (nroots) PetscValidPointer(selected,3);
834   PetscValidPointer(newsf,4);
835   if (sf->mine) for (i = 0; i < sf->nleaves; ++i) {leafsize = PetscMax(leafsize, sf->mine[i]+1);}
836   else leafsize = sf->nleaves;
837   ierr = PetscCalloc2(sf->nroots,&rootdata,leafsize,&leafdata);CHKERRQ(ierr);
838   for (i=0; i<nroots; ++i) rootdata[selected[i]] = 1;
839   ierr = PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr);
840   ierr = PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr);
841 
842   for (i = 0; i < leafsize; ++i) nleaves += leafdata[i];
843   ierr = PetscMalloc1(nleaves,&ilocal);CHKERRQ(ierr);
844   ierr = PetscMalloc1(nleaves,&iremote);CHKERRQ(ierr);
845   for (i = 0, n = 0; i < sf->nleaves; ++i) {
846     const PetscInt lidx = sf->mine ? sf->mine[i] : i;
847 
848     if (leafdata[lidx]) {
849       ilocal[n]        = lidx;
850       iremote[n].rank  = sf->remote[i].rank;
851       iremote[n].index = sf->remote[i].index;
852       ++n;
853     }
854   }
855   if (n != nleaves) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "There is a size mismatch in the SF embedding, %d != %d", n, nleaves);
856   ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,newsf);CHKERRQ(ierr);
857   ierr = PetscSFSetGraph(*newsf,sf->nroots,nleaves,ilocal,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
858   ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr);
859   PetscFunctionReturn(0);
860 }
861 
862 /*@C
863   PetscSFCreateEmbeddedLeafSF - removes edges from all but the selected leaves, does not remap indices
864 
865   Collective
866 
867   Input Arguments:
868 + sf - original star forest
869 . nleaves - number of leaves to select on this process
870 - selected - selected leaves on this process
871 
872   Output Arguments:
873 .  newsf - new star forest
874 
875   Level: advanced
876 
877 .seealso: PetscSFCreateEmbeddedSF(), PetscSFSetGraph(), PetscSFGetGraph()
878 @*/
879 PetscErrorCode PetscSFCreateEmbeddedLeafSF(PetscSF sf, PetscInt nleaves, const PetscInt *selected, PetscSF *newsf)
880 {
881   PetscSFNode   *iremote;
882   PetscInt      *ilocal;
883   PetscInt       i;
884   PetscErrorCode ierr;
885 
886   PetscFunctionBegin;
887   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
888   if (nleaves) PetscValidPointer(selected, 3);
889   PetscValidPointer(newsf, 4);
890   ierr = PetscMalloc1(nleaves, &ilocal);CHKERRQ(ierr);
891   ierr = PetscMalloc1(nleaves, &iremote);CHKERRQ(ierr);
892   for (i = 0; i < nleaves; ++i) {
893     const PetscInt l = selected[i];
894 
895     ilocal[i]        = sf->mine ? sf->mine[l] : l;
896     iremote[i].rank  = sf->remote[l].rank;
897     iremote[i].index = sf->remote[l].index;
898   }
899   ierr = PetscSFDuplicate(sf, PETSCSF_DUPLICATE_RANKS, newsf);CHKERRQ(ierr);
900   ierr = PetscSFSetGraph(*newsf, sf->nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr);
901   PetscFunctionReturn(0);
902 }
903 
904 /*@C
905    PetscSFBcastBegin - begin pointwise broadcast to be concluded with call to PetscSFBcastEnd()
906 
907    Collective on PetscSF
908 
909    Input Arguments:
910 +  sf - star forest on which to communicate
911 .  unit - data type associated with each node
912 -  rootdata - buffer to broadcast
913 
914    Output Arguments:
915 .  leafdata - buffer to update with values from each leaf's respective root
916 
917    Level: intermediate
918 
919 .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFView(), PetscSFBcastEnd(), PetscSFReduceBegin()
920 @*/
921 PetscErrorCode PetscSFBcastBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
922 {
923   PetscErrorCode ierr;
924 
925   PetscFunctionBegin;
926   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
927   PetscSFCheckGraphSet(sf,1);
928   ierr = PetscLogEventBegin(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr);
929   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
930   ierr = (*sf->ops->BcastBegin)(sf,unit,rootdata,leafdata);CHKERRQ(ierr);
931   ierr = PetscLogEventEnd(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr);
932   PetscFunctionReturn(0);
933 }
934 
935 /*@C
936    PetscSFBcastEnd - end a broadcast operation started with PetscSFBcastBegin()
937 
938    Collective
939 
940    Input Arguments:
941 +  sf - star forest
942 .  unit - data type
943 -  rootdata - buffer to broadcast
944 
945    Output Arguments:
946 .  leafdata - buffer to update with values from each leaf's respective root
947 
948    Level: intermediate
949 
950 .seealso: PetscSFSetGraph(), PetscSFReduceEnd()
951 @*/
952 PetscErrorCode PetscSFBcastEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
953 {
954   PetscErrorCode ierr;
955 
956   PetscFunctionBegin;
957   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
958   PetscSFCheckGraphSet(sf,1);
959   ierr = PetscLogEventBegin(PETSCSF_BcastEnd,sf,0,0,0);CHKERRQ(ierr);
960   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
961   ierr = (*sf->ops->BcastEnd)(sf,unit,rootdata,leafdata);CHKERRQ(ierr);
962   ierr = PetscLogEventEnd(PETSCSF_BcastEnd,sf,0,0,0);CHKERRQ(ierr);
963   PetscFunctionReturn(0);
964 }
965 
966 /*@C
967    PetscSFReduceBegin - begin reduction of leafdata into rootdata, to be completed with call to PetscSFReduceEnd()
968 
969    Collective
970 
971    Input Arguments:
972 +  sf - star forest
973 .  unit - data type
974 .  leafdata - values to reduce
975 -  op - reduction operation
976 
977    Output Arguments:
978 .  rootdata - result of reduction of values from all leaves of each root
979 
980    Level: intermediate
981 
982 .seealso: PetscSFBcastBegin()
983 @*/
984 PetscErrorCode PetscSFReduceBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
985 {
986   PetscErrorCode ierr;
987 
988   PetscFunctionBegin;
989   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
990   PetscSFCheckGraphSet(sf,1);
991   ierr = PetscLogEventBegin(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);
992   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
993   ierr = (sf->ops->ReduceBegin)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
994   ierr = PetscLogEventEnd(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);
995   PetscFunctionReturn(0);
996 }
997 
998 /*@C
999    PetscSFReduceEnd - end a reduction operation started with PetscSFReduceBegin()
1000 
1001    Collective
1002 
1003    Input Arguments:
1004 +  sf - star forest
1005 .  unit - data type
1006 .  leafdata - values to reduce
1007 -  op - reduction operation
1008 
1009    Output Arguments:
1010 .  rootdata - result of reduction of values from all leaves of each root
1011 
1012    Level: intermediate
1013 
1014 .seealso: PetscSFSetGraph(), PetscSFBcastEnd()
1015 @*/
1016 PetscErrorCode PetscSFReduceEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
1017 {
1018   PetscErrorCode ierr;
1019 
1020   PetscFunctionBegin;
1021   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1022   PetscSFCheckGraphSet(sf,1);
1023   ierr = PetscLogEventBegin(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);
1024   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1025   ierr = (*sf->ops->ReduceEnd)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
1026   ierr = PetscLogEventEnd(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);
1027   PetscFunctionReturn(0);
1028 }
1029 
1030 /*@C
1031    PetscSFComputeDegreeBegin - begin computation of degree for each root vertex, to be completed with PetscSFComputeDegreeEnd()
1032 
1033    Collective
1034 
1035    Input Arguments:
1036 .  sf - star forest
1037 
1038    Output Arguments:
1039 .  degree - degree of each root vertex
1040 
1041    Level: advanced
1042 
1043 .seealso: PetscSFGatherBegin()
1044 @*/
1045 PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf,const PetscInt **degree)
1046 {
1047   PetscErrorCode ierr;
1048 
1049   PetscFunctionBegin;
1050   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1051   PetscSFCheckGraphSet(sf,1);
1052   PetscValidPointer(degree,2);
1053   if (!sf->degreeknown) {
1054     PetscInt i,maxlocal;
1055     if (sf->degree) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Calls to PetscSFComputeDegreeBegin() cannot be nested.");
1056     for (i=0,maxlocal=0; i<sf->nleaves; i++) maxlocal = PetscMax(maxlocal,(sf->mine ? sf->mine[i] : i)+1);
1057     ierr = PetscMalloc1(sf->nroots,&sf->degree);CHKERRQ(ierr);
1058     ierr = PetscMalloc1(maxlocal,&sf->degreetmp);CHKERRQ(ierr);
1059     for (i=0; i<sf->nroots; i++) sf->degree[i] = 0;
1060     for (i=0; i<maxlocal; i++) sf->degreetmp[i] = 1;
1061     ierr = PetscSFReduceBegin(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);CHKERRQ(ierr);
1062   }
1063   *degree = NULL;
1064   PetscFunctionReturn(0);
1065 }
1066 
1067 /*@C
1068    PetscSFComputeDegreeEnd - complete computation of degree for each root vertex, started with PetscSFComputeDegreeBegin()
1069 
1070    Collective
1071 
1072    Input Arguments:
1073 .  sf - star forest
1074 
1075    Output Arguments:
1076 .  degree - degree of each root vertex
1077 
1078    Level: developer
1079 
1080 .seealso:
1081 @*/
1082 PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf,const PetscInt **degree)
1083 {
1084   PetscErrorCode ierr;
1085 
1086   PetscFunctionBegin;
1087   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1088   PetscSFCheckGraphSet(sf,1);
1089   if (!sf->degreeknown) {
1090     ierr = PetscSFReduceEnd(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);CHKERRQ(ierr);
1091     ierr = PetscFree(sf->degreetmp);CHKERRQ(ierr);
1092 
1093     sf->degreeknown = PETSC_TRUE;
1094   }
1095   *degree = sf->degree;
1096   PetscFunctionReturn(0);
1097 }
1098 
1099 /*@C
1100    PetscSFFetchAndOpBegin - begin operation that fetches values from root and updates atomically by applying operation using my leaf value, to be completed with PetscSFFetchAndOpEnd()
1101 
1102    Collective
1103 
1104    Input Arguments:
1105 +  sf - star forest
1106 .  unit - data type
1107 .  leafdata - leaf values to use in reduction
1108 -  op - operation to use for reduction
1109 
1110    Output Arguments:
1111 +  rootdata - root values to be updated, input state is seen by first process to perform an update
1112 -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
1113 
1114    Level: advanced
1115 
1116    Note:
1117    The update is only atomic at the granularity provided by the hardware. Different roots referenced by the same process
1118    might be updated in a different order. Furthermore, if a composite type is used for the unit datatype, atomicity is
1119    not guaranteed across the whole vertex. Therefore, this function is mostly only used with primitive types such as
1120    integers.
1121 
1122 .seealso: PetscSFComputeDegreeBegin(), PetscSFReduceBegin(), PetscSFSetGraph()
1123 @*/
1124 PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1125 {
1126   PetscErrorCode ierr;
1127 
1128   PetscFunctionBegin;
1129   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1130   PetscSFCheckGraphSet(sf,1);
1131   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
1132   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1133   ierr = (*sf->ops->FetchAndOpBegin)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr);
1134   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
1135   PetscFunctionReturn(0);
1136 }
1137 
1138 /*@C
1139    PetscSFFetchAndOpEnd - end operation started in matching call to PetscSFFetchAndOpBegin() to fetch values from roots and update atomically by applying operation using my leaf value
1140 
1141    Collective
1142 
1143    Input Arguments:
1144 +  sf - star forest
1145 .  unit - data type
1146 .  leafdata - leaf values to use in reduction
1147 -  op - operation to use for reduction
1148 
1149    Output Arguments:
1150 +  rootdata - root values to be updated, input state is seen by first process to perform an update
1151 -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
1152 
1153    Level: advanced
1154 
1155 .seealso: PetscSFComputeDegreeEnd(), PetscSFReduceEnd(), PetscSFSetGraph()
1156 @*/
1157 PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1158 {
1159   PetscErrorCode ierr;
1160 
1161   PetscFunctionBegin;
1162   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1163   PetscSFCheckGraphSet(sf,1);
1164   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
1165   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1166   ierr = (*sf->ops->FetchAndOpEnd)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr);
1167   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
1168   PetscFunctionReturn(0);
1169 }
1170 
1171 /*@C
1172    PetscSFGatherBegin - begin pointwise gather of all leaves into multi-roots, to be completed with PetscSFGatherEnd()
1173 
1174    Collective
1175 
1176    Input Arguments:
1177 +  sf - star forest
1178 .  unit - data type
1179 -  leafdata - leaf data to gather to roots
1180 
1181    Output Argument:
1182 .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
1183 
1184    Level: intermediate
1185 
1186 .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
1187 @*/
1188 PetscErrorCode PetscSFGatherBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
1189 {
1190   PetscErrorCode ierr;
1191   PetscSF        multi;
1192 
1193   PetscFunctionBegin;
1194   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1195   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
1196   ierr = PetscSFReduceBegin(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr);
1197   PetscFunctionReturn(0);
1198 }
1199 
1200 /*@C
1201    PetscSFGatherEnd - ends pointwise gather operation that was started with PetscSFGatherBegin()
1202 
1203    Collective
1204 
1205    Input Arguments:
1206 +  sf - star forest
1207 .  unit - data type
1208 -  leafdata - leaf data to gather to roots
1209 
1210    Output Argument:
1211 .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
1212 
1213    Level: intermediate
1214 
1215 .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
1216 @*/
1217 PetscErrorCode PetscSFGatherEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
1218 {
1219   PetscErrorCode ierr;
1220   PetscSF        multi;
1221 
1222   PetscFunctionBegin;
1223   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1224   PetscSFCheckGraphSet(sf,1);
1225   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1226   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
1227   ierr = PetscSFReduceEnd(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr);
1228   PetscFunctionReturn(0);
1229 }
1230 
1231 /*@C
1232    PetscSFScatterBegin - begin pointwise scatter operation from multi-roots to leaves, to be completed with PetscSFScatterEnd()
1233 
1234    Collective
1235 
1236    Input Arguments:
1237 +  sf - star forest
1238 .  unit - data type
1239 -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
1240 
1241    Output Argument:
1242 .  leafdata - leaf data to be update with personal data from each respective root
1243 
1244    Level: intermediate
1245 
1246 .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
1247 @*/
1248 PetscErrorCode PetscSFScatterBegin(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
1249 {
1250   PetscErrorCode ierr;
1251   PetscSF        multi;
1252 
1253   PetscFunctionBegin;
1254   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1255   PetscSFCheckGraphSet(sf,1);
1256   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1257   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
1258   ierr = PetscSFBcastBegin(multi,unit,multirootdata,leafdata);CHKERRQ(ierr);
1259   PetscFunctionReturn(0);
1260 }
1261 
1262 /*@C
1263    PetscSFScatterEnd - ends pointwise scatter operation that was started with PetscSFScatterBegin()
1264 
1265    Collective
1266 
1267    Input Arguments:
1268 +  sf - star forest
1269 .  unit - data type
1270 -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
1271 
1272    Output Argument:
1273 .  leafdata - leaf data to be update with personal data from each respective root
1274 
1275    Level: intermediate
1276 
1277 .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
1278 @*/
1279 PetscErrorCode PetscSFScatterEnd(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
1280 {
1281   PetscErrorCode ierr;
1282   PetscSF        multi;
1283 
1284   PetscFunctionBegin;
1285   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1286   PetscSFCheckGraphSet(sf,1);
1287   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1288   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
1289   ierr = PetscSFBcastEnd(multi,unit,multirootdata,leafdata);CHKERRQ(ierr);
1290   PetscFunctionReturn(0);
1291 }
1292 
1293 /*@
1294   PetscSFCompose - Compose a new PetscSF equivalent to action to PetscSFs
1295 
1296   Input Parameters:
1297 + sfA - The first PetscSF
1298 - sfB - The second PetscSF
1299 
1300   Output Parameters:
1301 . sfBA - equvalent PetscSF for applying A then B
1302 
1303   Level: developer
1304 
1305 .seealso: PetscSF, PetscSFGetGraph(), PetscSFSetGraph()
1306 @*/
1307 PetscErrorCode PetscSFCompose(PetscSF sfA, PetscSF sfB, PetscSF *sfBA)
1308 {
1309   MPI_Comm           comm;
1310   const PetscSFNode *remotePointsA, *remotePointsB;
1311   PetscSFNode       *remotePointsBA;
1312   const PetscInt    *localPointsA, *localPointsB;
1313   PetscInt           numRootsA, numLeavesA, numRootsB, numLeavesB;
1314   PetscErrorCode     ierr;
1315 
1316   PetscFunctionBegin;
1317   PetscValidHeaderSpecific(sfA, PETSCSF_CLASSID, 1);
1318   PetscValidHeaderSpecific(sfB, PETSCSF_CLASSID, 1);
1319   ierr = PetscObjectGetComm((PetscObject) sfA, &comm);CHKERRQ(ierr);
1320   ierr = PetscSFGetGraph(sfA, &numRootsA, &numLeavesA, &localPointsA, &remotePointsA);CHKERRQ(ierr);
1321   ierr = PetscSFGetGraph(sfB, &numRootsB, &numLeavesB, &localPointsB, &remotePointsB);CHKERRQ(ierr);
1322   ierr = PetscMalloc1(numLeavesB, &remotePointsBA);CHKERRQ(ierr);
1323   ierr = PetscSFBcastBegin(sfB, MPIU_2INT, remotePointsA, remotePointsBA);CHKERRQ(ierr);
1324   ierr = PetscSFBcastEnd(sfB, MPIU_2INT, remotePointsA, remotePointsBA);CHKERRQ(ierr);
1325   ierr = PetscSFCreate(comm, sfBA);CHKERRQ(ierr);
1326   ierr = PetscSFSetGraph(*sfBA, numRootsA, numLeavesB, localPointsB, PETSC_COPY_VALUES, remotePointsBA, PETSC_OWN_POINTER);CHKERRQ(ierr);
1327   PetscFunctionReturn(0);
1328 }
1329