xref: /petsc/src/vec/is/sf/interface/sf.c (revision a8c5bd36e19d36ea12c52b76bb2fa9fb5f6fb596)
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           i,*rcount,*ranks;
721 
722   PetscFunctionBegin;
723   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
724   PetscSFCheckGraphSet(sf,1);
725   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr);
726   ierr = PetscTableCreate(10,size,&table);CHKERRQ(ierr);
727   for (i=0; i<sf->nleaves; i++) {
728     /* Log 1-based rank */
729     ierr = PetscTableAdd(table,sf->remote[i].rank+1,1,ADD_VALUES);CHKERRQ(ierr);
730   }
731   ierr = PetscTableGetCount(table,&sf->nranks);CHKERRQ(ierr);
732   ierr = PetscMalloc4(sf->nranks,&sf->ranks,sf->nranks+1,&sf->roffset,sf->nleaves,&sf->rmine,sf->nleaves,&sf->rremote);CHKERRQ(ierr);
733   ierr = PetscMalloc2(sf->nranks,&rcount,sf->nranks,&ranks);CHKERRQ(ierr);
734   ierr = PetscTableGetHeadPosition(table,&pos);CHKERRQ(ierr);
735   for (i=0; i<sf->nranks; i++) {
736     ierr = PetscTableGetNext(table,&pos,&ranks[i],&rcount[i]);CHKERRQ(ierr);
737     ranks[i]--;             /* Convert back to 0-based */
738   }
739   ierr = PetscTableDestroy(&table);CHKERRQ(ierr);
740 
741   /* We expect that dgroup is reliably "small" while nranks could be large */
742   {
743     MPI_Group group = MPI_GROUP_NULL;
744     PetscMPIInt *dgroupranks;
745     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
746     ierr = MPI_Group_size(dgroup,&groupsize);CHKERRQ(ierr);
747     ierr = PetscMalloc1(groupsize,&dgroupranks);CHKERRQ(ierr);
748     ierr = PetscMalloc1(groupsize,&groupranks);CHKERRQ(ierr);
749     for (i=0; i<groupsize; i++) dgroupranks[i] = i;
750     if (groupsize) {ierr = MPI_Group_translate_ranks(dgroup,groupsize,dgroupranks,group,groupranks);CHKERRQ(ierr);}
751     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
752     ierr = PetscFree(dgroupranks);CHKERRQ(ierr);
753   }
754 
755   /* Partition ranks[] into distinguished (first sf->ndranks) followed by non-distinguished */
756   for (sf->ndranks=0,i=sf->nranks; sf->ndranks<i; ) {
757     for (i--; sf->ndranks<i; i--) { /* Scan i backward looking for distinguished rank */
758       if (InList(ranks[i],groupsize,groupranks)) break;
759     }
760     for ( ; sf->ndranks<=i; sf->ndranks++) { /* Scan sf->ndranks forward looking for non-distinguished rank */
761       if (!InList(ranks[sf->ndranks],groupsize,groupranks)) break;
762     }
763     if (sf->ndranks < i) {                         /* Swap ranks[sf->ndranks] with ranks[i] */
764       PetscInt    tmprank,tmpcount;
765       tmprank = ranks[i];
766       tmpcount = rcount[i];
767       ranks[i] = ranks[sf->ndranks];
768       rcount[i] = rcount[sf->ndranks];
769       ranks[sf->ndranks] = tmprank;
770       rcount[sf->ndranks] = tmpcount;
771       sf->ndranks++;
772     }
773   }
774   ierr = PetscFree(groupranks);CHKERRQ(ierr);
775   ierr = PetscSortIntWithArray(sf->ndranks,ranks,rcount);CHKERRQ(ierr);
776   ierr = PetscSortIntWithArray(sf->nranks-sf->ndranks,ranks+sf->ndranks,rcount+sf->ndranks);CHKERRQ(ierr);
777   sf->roffset[0] = 0;
778   for (i=0; i<sf->nranks; i++) {
779     ierr = PetscMPIIntCast(ranks[i],sf->ranks+i);CHKERRQ(ierr);
780     sf->roffset[i+1] = sf->roffset[i] + rcount[i];
781     rcount[i]        = 0;
782   }
783   for (i=0; i<sf->nleaves; i++) {
784     PetscInt irank;
785     /* Search for index of iremote[i].rank in sf->ranks */
786     ierr = PetscFindMPIInt(sf->remote[i].rank,sf->ndranks,sf->ranks,&irank);CHKERRQ(ierr);
787     if (irank < 0) {
788       ierr = PetscFindMPIInt(sf->remote[i].rank,sf->nranks-sf->ndranks,sf->ranks+sf->ndranks,&irank);CHKERRQ(ierr);
789       if (irank >= 0) irank += sf->ndranks;
790     }
791     if (irank < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Could not find rank %D in array",sf->remote[i].rank);
792     sf->rmine[sf->roffset[irank] + rcount[irank]]   = sf->mine ? sf->mine[i] : i;
793     sf->rremote[sf->roffset[irank] + rcount[irank]] = sf->remote[i].index;
794     rcount[irank]++;
795   }
796   ierr = PetscFree2(rcount,ranks);CHKERRQ(ierr);
797   PetscFunctionReturn(0);
798 }
799 
800 /*@C
801    PetscSFGetGroups - gets incoming and outgoing process groups
802 
803    Collective
804 
805    Input Argument:
806 .  sf - star forest
807 
808    Output Arguments:
809 +  incoming - group of origin processes for incoming edges (leaves that reference my roots)
810 -  outgoing - group of destination processes for outgoing edges (roots that I reference)
811 
812    Level: developer
813 
814 .seealso: PetscSFGetWindow(), PetscSFRestoreWindow()
815 @*/
816 PetscErrorCode PetscSFGetGroups(PetscSF sf,MPI_Group *incoming,MPI_Group *outgoing)
817 {
818   PetscErrorCode ierr;
819   MPI_Group      group = MPI_GROUP_NULL;
820 
821   PetscFunctionBegin;
822   if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining groups");
823   if (sf->ingroup == MPI_GROUP_NULL) {
824     PetscInt       i;
825     const PetscInt *indegree;
826     PetscMPIInt    rank,*outranks,*inranks;
827     PetscSFNode    *remote;
828     PetscSF        bgcount;
829 
830     /* Compute the number of incoming ranks */
831     ierr = PetscMalloc1(sf->nranks,&remote);CHKERRQ(ierr);
832     for (i=0; i<sf->nranks; i++) {
833       remote[i].rank  = sf->ranks[i];
834       remote[i].index = 0;
835     }
836     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,&bgcount);CHKERRQ(ierr);
837     ierr = PetscSFSetGraph(bgcount,1,sf->nranks,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
838     ierr = PetscSFComputeDegreeBegin(bgcount,&indegree);CHKERRQ(ierr);
839     ierr = PetscSFComputeDegreeEnd(bgcount,&indegree);CHKERRQ(ierr);
840 
841     /* Enumerate the incoming ranks */
842     ierr = PetscMalloc2(indegree[0],&inranks,sf->nranks,&outranks);CHKERRQ(ierr);
843     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
844     for (i=0; i<sf->nranks; i++) outranks[i] = rank;
845     ierr = PetscSFGatherBegin(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr);
846     ierr = PetscSFGatherEnd(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr);
847     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
848     ierr = MPI_Group_incl(group,indegree[0],inranks,&sf->ingroup);CHKERRQ(ierr);
849     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
850     ierr = PetscFree2(inranks,outranks);CHKERRQ(ierr);
851     ierr = PetscSFDestroy(&bgcount);CHKERRQ(ierr);
852   }
853   *incoming = sf->ingroup;
854 
855   if (sf->outgroup == MPI_GROUP_NULL) {
856     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
857     ierr = MPI_Group_incl(group,sf->nranks,sf->ranks,&sf->outgroup);CHKERRQ(ierr);
858     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
859   }
860   *outgoing = sf->outgroup;
861   PetscFunctionReturn(0);
862 }
863 
864 /*@
865    PetscSFGetMultiSF - gets the inner SF implemeting gathers and scatters
866 
867    Collective
868 
869    Input Argument:
870 .  sf - star forest that may contain roots with 0 or with more than 1 vertex
871 
872    Output Arguments:
873 .  multi - star forest with split roots, such that each root has degree exactly 1
874 
875    Level: developer
876 
877    Notes:
878 
879    In most cases, users should use PetscSFGatherBegin() and PetscSFScatterBegin() instead of manipulating multi
880    directly. Since multi satisfies the stronger condition that each entry in the global space has exactly one incoming
881    edge, it is a candidate for future optimization that might involve its removal.
882 
883 .seealso: PetscSFSetGraph(), PetscSFGatherBegin(), PetscSFScatterBegin(), PetscSFComputeMultiRootOriginalNumbering()
884 @*/
885 PetscErrorCode PetscSFGetMultiSF(PetscSF sf,PetscSF *multi)
886 {
887   PetscErrorCode ierr;
888 
889   PetscFunctionBegin;
890   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
891   PetscValidPointer(multi,2);
892   if (sf->nroots < 0) {         /* Graph has not been set yet; why do we need this? */
893     ierr   = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr);
894     *multi = sf->multi;
895     PetscFunctionReturn(0);
896   }
897   if (!sf->multi) {
898     const PetscInt *indegree;
899     PetscInt       i,*inoffset,*outones,*outoffset,maxlocal;
900     PetscSFNode    *remote;
901     maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */
902     ierr = PetscSFComputeDegreeBegin(sf,&indegree);CHKERRQ(ierr);
903     ierr = PetscSFComputeDegreeEnd(sf,&indegree);CHKERRQ(ierr);
904     ierr = PetscMalloc3(sf->nroots+1,&inoffset,maxlocal,&outones,maxlocal,&outoffset);CHKERRQ(ierr);
905     inoffset[0] = 0;
906     for (i=0; i<sf->nroots; i++) inoffset[i+1] = inoffset[i] + indegree[i];
907     for (i=0; i<maxlocal; i++) outones[i] = 1;
908     ierr = PetscSFFetchAndOpBegin(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr);
909     ierr = PetscSFFetchAndOpEnd(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr);
910     for (i=0; i<sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */
911 #if 0
912 #if defined(PETSC_USE_DEBUG)                                 /* Check that the expected number of increments occurred */
913     for (i=0; i<sf->nroots; i++) {
914       if (inoffset[i] + indegree[i] != inoffset[i+1]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect result after PetscSFFetchAndOp");
915     }
916 #endif
917 #endif
918     ierr = PetscMalloc1(sf->nleaves,&remote);CHKERRQ(ierr);
919     for (i=0; i<sf->nleaves; i++) {
920       remote[i].rank  = sf->remote[i].rank;
921       remote[i].index = outoffset[sf->mine ? sf->mine[i] : i];
922     }
923     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr);
924     ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
925     if (sf->rankorder) {        /* Sort the ranks */
926       PetscMPIInt rank;
927       PetscInt    *inranks,*newoffset,*outranks,*newoutoffset,*tmpoffset,maxdegree;
928       PetscSFNode *newremote;
929       ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
930       for (i=0,maxdegree=0; i<sf->nroots; i++) maxdegree = PetscMax(maxdegree,indegree[i]);
931       ierr = PetscMalloc5(sf->multi->nroots,&inranks,sf->multi->nroots,&newoffset,maxlocal,&outranks,maxlocal,&newoutoffset,maxdegree,&tmpoffset);CHKERRQ(ierr);
932       for (i=0; i<maxlocal; i++) outranks[i] = rank;
933       ierr = PetscSFReduceBegin(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr);
934       ierr = PetscSFReduceEnd(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr);
935       /* Sort the incoming ranks at each vertex, build the inverse map */
936       for (i=0; i<sf->nroots; i++) {
937         PetscInt j;
938         for (j=0; j<indegree[i]; j++) tmpoffset[j] = j;
939         ierr = PetscSortIntWithArray(indegree[i],inranks+inoffset[i],tmpoffset);CHKERRQ(ierr);
940         for (j=0; j<indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j;
941       }
942       ierr = PetscSFBcastBegin(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr);
943       ierr = PetscSFBcastEnd(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr);
944       ierr = PetscMalloc1(sf->nleaves,&newremote);CHKERRQ(ierr);
945       for (i=0; i<sf->nleaves; i++) {
946         newremote[i].rank  = sf->remote[i].rank;
947         newremote[i].index = newoutoffset[sf->mine ? sf->mine[i] : i];
948       }
949       ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,newremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
950       ierr = PetscFree5(inranks,newoffset,outranks,newoutoffset,tmpoffset);CHKERRQ(ierr);
951     }
952     ierr = PetscFree3(inoffset,outones,outoffset);CHKERRQ(ierr);
953   }
954   *multi = sf->multi;
955   PetscFunctionReturn(0);
956 }
957 
958 /*@C
959    PetscSFCreateEmbeddedSF - removes edges from all but the selected roots, does not remap indices
960 
961    Collective
962 
963    Input Arguments:
964 +  sf - original star forest
965 .  nroots - number of roots to select on this process
966 -  selected - selected roots on this process
967 
968    Output Arguments:
969 .  newsf - new star forest
970 
971    Level: advanced
972 
973    Note:
974    To use the new PetscSF, it may be necessary to know the indices of the leaves that are still participating. This can
975    be done by calling PetscSFGetGraph().
976 
977 .seealso: PetscSFSetGraph(), PetscSFGetGraph()
978 @*/
979 PetscErrorCode PetscSFCreateEmbeddedSF(PetscSF sf,PetscInt nroots,const PetscInt *selected,PetscSF *newsf)
980 {
981   PetscInt      *rootdata, *leafdata, *ilocal;
982   PetscSFNode   *iremote;
983   PetscInt       leafsize = 0, nleaves = 0, n, i;
984   PetscErrorCode ierr;
985 
986   PetscFunctionBegin;
987   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
988   PetscSFCheckGraphSet(sf,1);
989   if (nroots) PetscValidPointer(selected,3);
990   PetscValidPointer(newsf,4);
991   ierr = PetscLogEventBegin(PETSCSF_EmbedSF,sf,0,0,0);CHKERRQ(ierr);
992   if (sf->mine) for (i = 0; i < sf->nleaves; ++i) {leafsize = PetscMax(leafsize, sf->mine[i]+1);}
993   else leafsize = sf->nleaves;
994   ierr = PetscCalloc2(sf->nroots,&rootdata,leafsize,&leafdata);CHKERRQ(ierr);
995   for (i=0; i<nroots; ++i) rootdata[selected[i]] = 1;
996   ierr = PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr);
997   ierr = PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr);
998   for (i = 0; i < leafsize; ++i) nleaves += leafdata[i];
999   ierr = PetscMalloc1(nleaves,&ilocal);CHKERRQ(ierr);
1000   ierr = PetscMalloc1(nleaves,&iremote);CHKERRQ(ierr);
1001   for (i = 0, n = 0; i < sf->nleaves; ++i) {
1002     const PetscInt lidx = sf->mine ? sf->mine[i] : i;
1003 
1004     if (leafdata[lidx]) {
1005       ilocal[n]        = lidx;
1006       iremote[n].rank  = sf->remote[i].rank;
1007       iremote[n].index = sf->remote[i].index;
1008       ++n;
1009     }
1010   }
1011   if (n != nleaves) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "There is a size mismatch in the SF embedding, %d != %d", n, nleaves);
1012   ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,newsf);CHKERRQ(ierr);
1013   ierr = PetscSFSetGraph(*newsf,sf->nroots,nleaves,ilocal,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1014   ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr);
1015   ierr = PetscSFSetUp(*newsf);CHKERRQ(ierr);
1016   ierr = PetscLogEventEnd(PETSCSF_EmbedSF,sf,0,0,0);CHKERRQ(ierr);
1017   PetscFunctionReturn(0);
1018 }
1019 
1020 /*@C
1021   PetscSFCreateEmbeddedLeafSF - removes edges from all but the selected leaves, does not remap indices
1022 
1023   Collective
1024 
1025   Input Arguments:
1026 + sf - original star forest
1027 . nleaves - number of leaves to select on this process
1028 - selected - selected leaves on this process
1029 
1030   Output Arguments:
1031 .  newsf - new star forest
1032 
1033   Level: advanced
1034 
1035 .seealso: PetscSFCreateEmbeddedSF(), PetscSFSetGraph(), PetscSFGetGraph()
1036 @*/
1037 PetscErrorCode PetscSFCreateEmbeddedLeafSF(PetscSF sf, PetscInt nleaves, const PetscInt *selected, PetscSF *newsf)
1038 {
1039   PetscSFNode   *iremote;
1040   PetscInt      *ilocal;
1041   PetscInt       i;
1042   PetscErrorCode ierr;
1043 
1044   PetscFunctionBegin;
1045   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
1046   PetscSFCheckGraphSet(sf, 1);
1047   if (nleaves) PetscValidPointer(selected, 3);
1048   PetscValidPointer(newsf, 4);
1049   ierr = PetscMalloc1(nleaves, &ilocal);CHKERRQ(ierr);
1050   ierr = PetscMalloc1(nleaves, &iremote);CHKERRQ(ierr);
1051   for (i = 0; i < nleaves; ++i) {
1052     const PetscInt l = selected[i];
1053 
1054     ilocal[i]        = sf->mine ? sf->mine[l] : l;
1055     iremote[i].rank  = sf->remote[l].rank;
1056     iremote[i].index = sf->remote[l].index;
1057   }
1058   ierr = PetscSFDuplicate(sf, PETSCSF_DUPLICATE_RANKS, newsf);CHKERRQ(ierr);
1059   ierr = PetscSFSetGraph(*newsf, sf->nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr);
1060   PetscFunctionReturn(0);
1061 }
1062 
1063 /*@C
1064    PetscSFBcastBegin - begin pointwise broadcast to be concluded with call to PetscSFBcastEnd()
1065 
1066    Collective on PetscSF
1067 
1068    Input Arguments:
1069 +  sf - star forest on which to communicate
1070 .  unit - data type associated with each node
1071 -  rootdata - buffer to broadcast
1072 
1073    Output Arguments:
1074 .  leafdata - buffer to update with values from each leaf's respective root
1075 
1076    Level: intermediate
1077 
1078 .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFView(), PetscSFBcastEnd(), PetscSFReduceBegin()
1079 @*/
1080 PetscErrorCode PetscSFBcastBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
1081 {
1082   PetscErrorCode ierr;
1083 
1084   PetscFunctionBegin;
1085   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1086   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1087   ierr = PetscLogEventBegin(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr);
1088   ierr = (*sf->ops->BcastBegin)(sf,unit,rootdata,leafdata);CHKERRQ(ierr);
1089   ierr = PetscLogEventEnd(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr);
1090   PetscFunctionReturn(0);
1091 }
1092 
1093 /*@C
1094    PetscSFBcastEnd - end a broadcast operation started with PetscSFBcastBegin()
1095 
1096    Collective
1097 
1098    Input Arguments:
1099 +  sf - star forest
1100 .  unit - data type
1101 -  rootdata - buffer to broadcast
1102 
1103    Output Arguments:
1104 .  leafdata - buffer to update with values from each leaf's respective root
1105 
1106    Level: intermediate
1107 
1108 .seealso: PetscSFSetGraph(), PetscSFReduceEnd()
1109 @*/
1110 PetscErrorCode PetscSFBcastEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
1111 {
1112   PetscErrorCode ierr;
1113 
1114   PetscFunctionBegin;
1115   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1116   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1117   ierr = PetscLogEventBegin(PETSCSF_BcastEnd,sf,0,0,0);CHKERRQ(ierr);
1118   ierr = (*sf->ops->BcastEnd)(sf,unit,rootdata,leafdata);CHKERRQ(ierr);
1119   ierr = PetscLogEventEnd(PETSCSF_BcastEnd,sf,0,0,0);CHKERRQ(ierr);
1120   PetscFunctionReturn(0);
1121 }
1122 
1123 /*@C
1124    PetscSFBcastAndOpBegin - begin pointwise broadcast with root value being reduced to leaf value, to be concluded with call to PetscSFBcastAndOpEnd()
1125 
1126    Collective on PetscSF
1127 
1128    Input Arguments:
1129 +  sf - star forest on which to communicate
1130 .  unit - data type associated with each node
1131 .  rootdata - buffer to broadcast
1132 -  op - operation to use for reduction
1133 
1134    Output Arguments:
1135 .  leafdata - buffer to be reduced with values from each leaf's respective root
1136 
1137    Level: intermediate
1138 
1139 .seealso: PetscSFBcastAndOpEnd(), PetscSFBcastBegin(), PetscSFBcastEnd()
1140 @*/
1141 PetscErrorCode PetscSFBcastAndOpBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
1142 {
1143   PetscErrorCode ierr;
1144 
1145   PetscFunctionBegin;
1146   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1147   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1148   ierr = PetscLogEventBegin(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
1149   ierr = (*sf->ops->BcastAndOpBegin)(sf,unit,rootdata,leafdata,op);CHKERRQ(ierr);
1150   ierr = PetscLogEventEnd(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
1151   PetscFunctionReturn(0);
1152 }
1153 
1154 /*@C
1155    PetscSFBcastAndOpEnd - end a broadcast & reduce operation started with PetscSFBcastAndOpBegin()
1156 
1157    Collective
1158 
1159    Input Arguments:
1160 +  sf - star forest
1161 .  unit - data type
1162 .  rootdata - buffer to broadcast
1163 -  op - operation to use for reduction
1164 
1165    Output Arguments:
1166 .  leafdata - buffer to be reduced with values from each leaf's respective root
1167 
1168    Level: intermediate
1169 
1170 .seealso: PetscSFSetGraph(), PetscSFReduceEnd()
1171 @*/
1172 PetscErrorCode PetscSFBcastAndOpEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
1173 {
1174   PetscErrorCode ierr;
1175 
1176   PetscFunctionBegin;
1177   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1178   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1179   ierr = PetscLogEventBegin(PETSCSF_BcastAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
1180   ierr = (*sf->ops->BcastAndOpEnd)(sf,unit,rootdata,leafdata,op);CHKERRQ(ierr);
1181   ierr = PetscLogEventEnd(PETSCSF_BcastAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
1182   PetscFunctionReturn(0);
1183 }
1184 
1185 /*@C
1186    PetscSFReduceBegin - begin reduction of leafdata into rootdata, to be completed with call to PetscSFReduceEnd()
1187 
1188    Collective
1189 
1190    Input Arguments:
1191 +  sf - star forest
1192 .  unit - data type
1193 .  leafdata - values to reduce
1194 -  op - reduction operation
1195 
1196    Output Arguments:
1197 .  rootdata - result of reduction of values from all leaves of each root
1198 
1199    Level: intermediate
1200 
1201 .seealso: PetscSFBcastBegin()
1202 @*/
1203 PetscErrorCode PetscSFReduceBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
1204 {
1205   PetscErrorCode ierr;
1206 
1207   PetscFunctionBegin;
1208   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1209   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1210   ierr = PetscLogEventBegin(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);
1211   ierr = (sf->ops->ReduceBegin)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
1212   ierr = PetscLogEventEnd(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);
1213   PetscFunctionReturn(0);
1214 }
1215 
1216 /*@C
1217    PetscSFReduceEnd - end a reduction operation started with PetscSFReduceBegin()
1218 
1219    Collective
1220 
1221    Input Arguments:
1222 +  sf - star forest
1223 .  unit - data type
1224 .  leafdata - values to reduce
1225 -  op - reduction operation
1226 
1227    Output Arguments:
1228 .  rootdata - result of reduction of values from all leaves of each root
1229 
1230    Level: intermediate
1231 
1232 .seealso: PetscSFSetGraph(), PetscSFBcastEnd()
1233 @*/
1234 PetscErrorCode PetscSFReduceEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
1235 {
1236   PetscErrorCode ierr;
1237 
1238   PetscFunctionBegin;
1239   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1240   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1241   ierr = PetscLogEventBegin(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);
1242   ierr = (*sf->ops->ReduceEnd)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
1243   ierr = PetscLogEventEnd(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);
1244   PetscFunctionReturn(0);
1245 }
1246 
1247 /*@C
1248    PetscSFComputeDegreeBegin - begin computation of degree for each root vertex, to be completed with PetscSFComputeDegreeEnd()
1249 
1250    Collective
1251 
1252    Input Arguments:
1253 .  sf - star forest
1254 
1255    Output Arguments:
1256 .  degree - degree of each root vertex
1257 
1258    Level: advanced
1259 
1260    Notes:
1261    The returned array is owned by PetscSF and automatically freed by PetscSFDestroy(). Hence no need to call PetscFree() on it.
1262 
1263 .seealso: PetscSFGatherBegin()
1264 @*/
1265 PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf,const PetscInt **degree)
1266 {
1267   PetscErrorCode ierr;
1268 
1269   PetscFunctionBegin;
1270   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1271   PetscSFCheckGraphSet(sf,1);
1272   PetscValidPointer(degree,2);
1273   if (!sf->degreeknown) {
1274     PetscInt i, nroots = sf->nroots, maxlocal = sf->maxleaf+1;  /* TODO: We should use PetscSFGetLeafRange() */
1275     if (sf->degree) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Calls to PetscSFComputeDegreeBegin() cannot be nested.");
1276     ierr = PetscMalloc1(nroots,&sf->degree);CHKERRQ(ierr);
1277     ierr = PetscMalloc1(PetscMax(maxlocal,1),&sf->degreetmp);CHKERRQ(ierr); /* allocate at least one entry, see check in PetscSFComputeDegreeEnd() */
1278     for (i=0; i<nroots; i++) sf->degree[i] = 0;
1279     for (i=0; i<maxlocal; i++) sf->degreetmp[i] = 1;
1280     ierr = PetscSFReduceBegin(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);CHKERRQ(ierr);
1281   }
1282   *degree = NULL;
1283   PetscFunctionReturn(0);
1284 }
1285 
1286 /*@C
1287    PetscSFComputeDegreeEnd - complete computation of degree for each root vertex, started with PetscSFComputeDegreeBegin()
1288 
1289    Collective
1290 
1291    Input Arguments:
1292 .  sf - star forest
1293 
1294    Output Arguments:
1295 .  degree - degree of each root vertex
1296 
1297    Level: developer
1298 
1299    Notes:
1300    The returned array is owned by PetscSF and automatically freed by PetscSFDestroy(). Hence no need to call PetscFree() on it.
1301 
1302 .seealso:
1303 @*/
1304 PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf,const PetscInt **degree)
1305 {
1306   PetscErrorCode ierr;
1307 
1308   PetscFunctionBegin;
1309   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1310   PetscSFCheckGraphSet(sf,1);
1311   PetscValidPointer(degree,2);
1312   if (!sf->degreeknown) {
1313     if (!sf->degreetmp) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFComputeDegreeBegin() before PetscSFComputeDegreeEnd()");
1314     ierr = PetscSFReduceEnd(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);CHKERRQ(ierr);
1315     ierr = PetscFree(sf->degreetmp);CHKERRQ(ierr);
1316     sf->degreeknown = PETSC_TRUE;
1317   }
1318   *degree = sf->degree;
1319   PetscFunctionReturn(0);
1320 }
1321 
1322 
1323 /*@C
1324    PetscSFComputeMultiRootOriginalNumbering - Returns original numbering of multi-roots (roots of multi-SF returned by PetscSFGetMultiSF()).
1325    Each multi-root is assigned index of the corresponding original root.
1326 
1327    Collective
1328 
1329    Input Arguments:
1330 +  sf - star forest
1331 -  degree - degree of each root vertex, computed with PetscSFComputeDegreeBegin()/PetscSFComputeDegreeEnd()
1332 
1333    Output Arguments:
1334 +  nMultiRoots - (optional) number of multi-roots (roots of multi-SF)
1335 -  multiRootsOrigNumbering - original indices of multi-roots; length of this array is nMultiRoots
1336 
1337    Level: developer
1338 
1339    Notes:
1340    The returned array multiRootsOrigNumbering is newly allocated and should be destroyed with PetscFree() when no longer needed.
1341 
1342 .seealso: PetscSFComputeDegreeBegin(), PetscSFComputeDegreeEnd(), PetscSFGetMultiSF()
1343 @*/
1344 PetscErrorCode PetscSFComputeMultiRootOriginalNumbering(PetscSF sf, const PetscInt degree[], PetscInt *nMultiRoots, PetscInt *multiRootsOrigNumbering[])
1345 {
1346   PetscSF             msf;
1347   PetscInt            i, j, k, nroots, nmroots;
1348   PetscErrorCode      ierr;
1349 
1350   PetscFunctionBegin;
1351   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1352   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
1353   if (nroots) PetscValidIntPointer(degree,2);
1354   if (nMultiRoots) PetscValidIntPointer(nMultiRoots,3);
1355   PetscValidPointer(multiRootsOrigNumbering,4);
1356   ierr = PetscSFGetMultiSF(sf,&msf);CHKERRQ(ierr);
1357   ierr = PetscSFGetGraph(msf, &nmroots, NULL, NULL, NULL);CHKERRQ(ierr);
1358   ierr = PetscMalloc1(nmroots, multiRootsOrigNumbering);CHKERRQ(ierr);
1359   for (i=0,j=0,k=0; i<nroots; i++) {
1360     if (!degree[i]) continue;
1361     for (j=0; j<degree[i]; j++,k++) {
1362       (*multiRootsOrigNumbering)[k] = i;
1363     }
1364   }
1365 #if defined(PETSC_USE_DEBUG)
1366   if (PetscUnlikely(k != nmroots)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sanity check fail");
1367 #endif
1368   if (nMultiRoots) *nMultiRoots = nmroots;
1369   PetscFunctionReturn(0);
1370 }
1371 
1372 /*@C
1373    PetscSFFetchAndOpBegin - begin operation that fetches values from root and updates atomically by applying operation using my leaf value, to be completed with PetscSFFetchAndOpEnd()
1374 
1375    Collective
1376 
1377    Input Arguments:
1378 +  sf - star forest
1379 .  unit - data type
1380 .  leafdata - leaf values to use in reduction
1381 -  op - operation to use for reduction
1382 
1383    Output Arguments:
1384 +  rootdata - root values to be updated, input state is seen by first process to perform an update
1385 -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
1386 
1387    Level: advanced
1388 
1389    Note:
1390    The update is only atomic at the granularity provided by the hardware. Different roots referenced by the same process
1391    might be updated in a different order. Furthermore, if a composite type is used for the unit datatype, atomicity is
1392    not guaranteed across the whole vertex. Therefore, this function is mostly only used with primitive types such as
1393    integers.
1394 
1395 .seealso: PetscSFComputeDegreeBegin(), PetscSFReduceBegin(), PetscSFSetGraph()
1396 @*/
1397 PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1398 {
1399   PetscErrorCode ierr;
1400 
1401   PetscFunctionBegin;
1402   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1403   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1404   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
1405   ierr = (*sf->ops->FetchAndOpBegin)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr);
1406   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
1407   PetscFunctionReturn(0);
1408 }
1409 
1410 /*@C
1411    PetscSFFetchAndOpEnd - end operation started in matching call to PetscSFFetchAndOpBegin() to fetch values from roots and update atomically by applying operation using my leaf value
1412 
1413    Collective
1414 
1415    Input Arguments:
1416 +  sf - star forest
1417 .  unit - data type
1418 .  leafdata - leaf values to use in reduction
1419 -  op - operation to use for reduction
1420 
1421    Output Arguments:
1422 +  rootdata - root values to be updated, input state is seen by first process to perform an update
1423 -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
1424 
1425    Level: advanced
1426 
1427 .seealso: PetscSFComputeDegreeEnd(), PetscSFReduceEnd(), PetscSFSetGraph()
1428 @*/
1429 PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1430 {
1431   PetscErrorCode ierr;
1432 
1433   PetscFunctionBegin;
1434   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1435   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1436   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
1437   ierr = (*sf->ops->FetchAndOpEnd)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr);
1438   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
1439   PetscFunctionReturn(0);
1440 }
1441 
1442 /*@C
1443    PetscSFGatherBegin - begin pointwise gather of all leaves into multi-roots, to be completed with PetscSFGatherEnd()
1444 
1445    Collective
1446 
1447    Input Arguments:
1448 +  sf - star forest
1449 .  unit - data type
1450 -  leafdata - leaf data to gather to roots
1451 
1452    Output Argument:
1453 .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
1454 
1455    Level: intermediate
1456 
1457 .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
1458 @*/
1459 PetscErrorCode PetscSFGatherBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
1460 {
1461   PetscErrorCode ierr;
1462   PetscSF        multi;
1463 
1464   PetscFunctionBegin;
1465   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1466   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1467   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
1468   ierr = PetscSFReduceBegin(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr);
1469   PetscFunctionReturn(0);
1470 }
1471 
1472 /*@C
1473    PetscSFGatherEnd - ends pointwise gather operation that was started with PetscSFGatherBegin()
1474 
1475    Collective
1476 
1477    Input Arguments:
1478 +  sf - star forest
1479 .  unit - data type
1480 -  leafdata - leaf data to gather to roots
1481 
1482    Output Argument:
1483 .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
1484 
1485    Level: intermediate
1486 
1487 .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
1488 @*/
1489 PetscErrorCode PetscSFGatherEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
1490 {
1491   PetscErrorCode ierr;
1492   PetscSF        multi;
1493 
1494   PetscFunctionBegin;
1495   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1496   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1497   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
1498   ierr = PetscSFReduceEnd(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr);
1499   PetscFunctionReturn(0);
1500 }
1501 
1502 /*@C
1503    PetscSFScatterBegin - begin pointwise scatter operation from multi-roots to leaves, to be completed with PetscSFScatterEnd()
1504 
1505    Collective
1506 
1507    Input Arguments:
1508 +  sf - star forest
1509 .  unit - data type
1510 -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
1511 
1512    Output Argument:
1513 .  leafdata - leaf data to be update with personal data from each respective root
1514 
1515    Level: intermediate
1516 
1517 .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
1518 @*/
1519 PetscErrorCode PetscSFScatterBegin(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
1520 {
1521   PetscErrorCode ierr;
1522   PetscSF        multi;
1523 
1524   PetscFunctionBegin;
1525   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1526   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1527   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
1528   ierr = PetscSFBcastBegin(multi,unit,multirootdata,leafdata);CHKERRQ(ierr);
1529   PetscFunctionReturn(0);
1530 }
1531 
1532 /*@C
1533    PetscSFScatterEnd - ends pointwise scatter operation that was started with PetscSFScatterBegin()
1534 
1535    Collective
1536 
1537    Input Arguments:
1538 +  sf - star forest
1539 .  unit - data type
1540 -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
1541 
1542    Output Argument:
1543 .  leafdata - leaf data to be update with personal data from each respective root
1544 
1545    Level: intermediate
1546 
1547 .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
1548 @*/
1549 PetscErrorCode PetscSFScatterEnd(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
1550 {
1551   PetscErrorCode ierr;
1552   PetscSF        multi;
1553 
1554   PetscFunctionBegin;
1555   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1556   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1557   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
1558   ierr = PetscSFBcastEnd(multi,unit,multirootdata,leafdata);CHKERRQ(ierr);
1559   PetscFunctionReturn(0);
1560 }
1561 
1562 /*@
1563   PetscSFCompose - Compose a new PetscSF equivalent to action to PetscSFs
1564 
1565   Input Parameters:
1566 + sfA - The first PetscSF
1567 - sfB - The second PetscSF
1568 
1569   Output Parameters:
1570 . sfBA - equvalent PetscSF for applying A then B
1571 
1572   Level: developer
1573 
1574 .seealso: PetscSF, PetscSFGetGraph(), PetscSFSetGraph()
1575 @*/
1576 PetscErrorCode PetscSFCompose(PetscSF sfA, PetscSF sfB, PetscSF *sfBA)
1577 {
1578   MPI_Comm           comm;
1579   const PetscSFNode *remotePointsA, *remotePointsB;
1580   PetscSFNode       *remotePointsBA;
1581   const PetscInt    *localPointsA, *localPointsB;
1582   PetscInt           numRootsA, numLeavesA, numRootsB, numLeavesB;
1583   PetscErrorCode     ierr;
1584 
1585   PetscFunctionBegin;
1586   PetscValidHeaderSpecific(sfA, PETSCSF_CLASSID, 1);
1587   PetscSFCheckGraphSet(sfA, 1);
1588   PetscValidHeaderSpecific(sfB, PETSCSF_CLASSID, 2);
1589   PetscSFCheckGraphSet(sfB, 2);
1590   PetscValidPointer(sfBA, 3);
1591   ierr = PetscObjectGetComm((PetscObject) sfA, &comm);CHKERRQ(ierr);
1592   ierr = PetscSFGetGraph(sfA, &numRootsA, &numLeavesA, &localPointsA, &remotePointsA);CHKERRQ(ierr);
1593   ierr = PetscSFGetGraph(sfB, &numRootsB, &numLeavesB, &localPointsB, &remotePointsB);CHKERRQ(ierr);
1594   ierr = PetscMalloc1(numLeavesB, &remotePointsBA);CHKERRQ(ierr);
1595   ierr = PetscSFBcastBegin(sfB, MPIU_2INT, remotePointsA, remotePointsBA);CHKERRQ(ierr);
1596   ierr = PetscSFBcastEnd(sfB, MPIU_2INT, remotePointsA, remotePointsBA);CHKERRQ(ierr);
1597   ierr = PetscSFCreate(comm, sfBA);CHKERRQ(ierr);
1598   ierr = PetscSFSetGraph(*sfBA, numRootsA, numLeavesB, localPointsB, PETSC_COPY_VALUES, remotePointsBA, PETSC_OWN_POINTER);CHKERRQ(ierr);
1599   PetscFunctionReturn(0);
1600 }
1601