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