xref: /petsc/src/vec/is/sf/interface/sf.c (revision a2b725a8db0d6bf6cc2a1c6df7dd8029aadfff6e)
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 .  nselected  - number of selected roots on this process
966 -  selected   - indices of the 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 nselected,const PetscInt *selected,PetscSF *newsf)
980 {
981   PetscInt          *rootdata,*leafdata,*new_ilocal;
982   PetscSFNode       *new_iremote;
983   const PetscInt    *ilocal;
984   const PetscSFNode *iremote;
985   PetscInt          nleaves,nroots,n,i,new_nleaves = 0;
986   PetscSF           tmpsf;
987   PetscErrorCode    ierr;
988 
989   PetscFunctionBegin;
990   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
991   PetscSFCheckGraphSet(sf,1);
992   if (nselected) PetscValidPointer(selected,3);
993   PetscValidPointer(newsf,4);
994 
995   /* Find out which leaves (not leaf data items) are still connected to roots in the embedded sf */
996   ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
997   ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&tmpsf);CHKERRQ(ierr);
998   ierr = PetscSFSetGraph(tmpsf,nroots,nleaves,NULL/*contiguous*/,PETSC_USE_POINTER,iremote,PETSC_USE_POINTER);CHKERRQ(ierr);
999   ierr = PetscCalloc2(nroots,&rootdata,nleaves,&leafdata);CHKERRQ(ierr);
1000   for (i=0; i<nselected; ++i) {
1001     if (selected[i] < 0 || selected[i] >= nroots) SETERRQ2(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Root index %D is not in [0,%D)",selected[i],nroots);
1002     rootdata[selected[i]] = 1;
1003   }
1004 
1005   ierr = PetscSFBcastBegin(tmpsf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr);
1006   ierr = PetscSFBcastEnd(tmpsf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr);
1007   ierr = PetscSFDestroy(&tmpsf);CHKERRQ(ierr);
1008 
1009   /* Build newsf with leaves that are still connected */
1010   for (i = 0; i < nleaves; ++i) new_nleaves += leafdata[i];
1011   ierr = PetscMalloc1(new_nleaves,&new_ilocal);CHKERRQ(ierr);
1012   ierr = PetscMalloc1(new_nleaves,&new_iremote);CHKERRQ(ierr);
1013   for (i = 0, n = 0; i < nleaves; ++i) {
1014     if (leafdata[i]) {
1015       new_ilocal[n]        = sf->mine ? sf->mine[i] : i;
1016       new_iremote[n].rank  = sf->remote[i].rank;
1017       new_iremote[n].index = sf->remote[i].index;
1018       ++n;
1019     }
1020   }
1021   if (n != new_nleaves) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"There is a size mismatch in the SF embedding, %d != %d",n,new_nleaves);
1022   ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,newsf);CHKERRQ(ierr);
1023   ierr = PetscSFSetGraph(*newsf,nroots,new_nleaves,new_ilocal,PETSC_OWN_POINTER,new_iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1024   ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr);
1025   PetscFunctionReturn(0);
1026 }
1027 
1028 /*@C
1029   PetscSFCreateEmbeddedLeafSF - removes edges from all but the selected leaves, does not remap indices
1030 
1031   Collective
1032 
1033   Input Arguments:
1034 + sf - original star forest
1035 . nleaves - number of leaves to select on this process
1036 - selected - selected leaves on this process
1037 
1038   Output Arguments:
1039 .  newsf - new star forest
1040 
1041   Level: advanced
1042 
1043 .seealso: PetscSFCreateEmbeddedSF(), PetscSFSetGraph(), PetscSFGetGraph()
1044 @*/
1045 PetscErrorCode PetscSFCreateEmbeddedLeafSF(PetscSF sf, PetscInt nleaves, const PetscInt *selected, PetscSF *newsf)
1046 {
1047   PetscSFNode   *iremote;
1048   PetscInt      *ilocal;
1049   PetscInt       i;
1050   PetscErrorCode ierr;
1051 
1052   PetscFunctionBegin;
1053   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
1054   PetscSFCheckGraphSet(sf, 1);
1055   if (nleaves) PetscValidPointer(selected, 3);
1056   PetscValidPointer(newsf, 4);
1057   ierr = PetscMalloc1(nleaves, &ilocal);CHKERRQ(ierr);
1058   ierr = PetscMalloc1(nleaves, &iremote);CHKERRQ(ierr);
1059   for (i = 0; i < nleaves; ++i) {
1060     const PetscInt l = selected[i];
1061 
1062     ilocal[i]        = sf->mine ? sf->mine[l] : l;
1063     iremote[i].rank  = sf->remote[l].rank;
1064     iremote[i].index = sf->remote[l].index;
1065   }
1066   ierr = PetscSFDuplicate(sf, PETSCSF_DUPLICATE_RANKS, newsf);CHKERRQ(ierr);
1067   ierr = PetscSFSetGraph(*newsf, sf->nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr);
1068   PetscFunctionReturn(0);
1069 }
1070 
1071 /*@C
1072    PetscSFBcastBegin - begin pointwise broadcast to be concluded with call to PetscSFBcastEnd()
1073 
1074    Collective on PetscSF
1075 
1076    Input Arguments:
1077 +  sf - star forest on which to communicate
1078 .  unit - data type associated with each node
1079 -  rootdata - buffer to broadcast
1080 
1081    Output Arguments:
1082 .  leafdata - buffer to update with values from each leaf's respective root
1083 
1084    Level: intermediate
1085 
1086 .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFView(), PetscSFBcastEnd(), PetscSFReduceBegin()
1087 @*/
1088 PetscErrorCode PetscSFBcastBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
1089 {
1090   PetscErrorCode ierr;
1091 
1092   PetscFunctionBegin;
1093   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1094   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1095   ierr = PetscLogEventBegin(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr);
1096   ierr = (*sf->ops->BcastBegin)(sf,unit,rootdata,leafdata);CHKERRQ(ierr);
1097   ierr = PetscLogEventEnd(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr);
1098   PetscFunctionReturn(0);
1099 }
1100 
1101 /*@C
1102    PetscSFBcastEnd - end a broadcast operation started with PetscSFBcastBegin()
1103 
1104    Collective
1105 
1106    Input Arguments:
1107 +  sf - star forest
1108 .  unit - data type
1109 -  rootdata - buffer to broadcast
1110 
1111    Output Arguments:
1112 .  leafdata - buffer to update with values from each leaf's respective root
1113 
1114    Level: intermediate
1115 
1116 .seealso: PetscSFSetGraph(), PetscSFReduceEnd()
1117 @*/
1118 PetscErrorCode PetscSFBcastEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
1119 {
1120   PetscErrorCode ierr;
1121 
1122   PetscFunctionBegin;
1123   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1124   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1125   ierr = PetscLogEventBegin(PETSCSF_BcastEnd,sf,0,0,0);CHKERRQ(ierr);
1126   ierr = (*sf->ops->BcastEnd)(sf,unit,rootdata,leafdata);CHKERRQ(ierr);
1127   ierr = PetscLogEventEnd(PETSCSF_BcastEnd,sf,0,0,0);CHKERRQ(ierr);
1128   PetscFunctionReturn(0);
1129 }
1130 
1131 /*@C
1132    PetscSFBcastAndOpBegin - begin pointwise broadcast with root value being reduced to leaf value, to be concluded with call to PetscSFBcastAndOpEnd()
1133 
1134    Collective on PetscSF
1135 
1136    Input Arguments:
1137 +  sf - star forest on which to communicate
1138 .  unit - data type associated with each node
1139 .  rootdata - buffer to broadcast
1140 -  op - operation to use for reduction
1141 
1142    Output Arguments:
1143 .  leafdata - buffer to be reduced with values from each leaf's respective root
1144 
1145    Level: intermediate
1146 
1147 .seealso: PetscSFBcastAndOpEnd(), PetscSFBcastBegin(), PetscSFBcastEnd()
1148 @*/
1149 PetscErrorCode PetscSFBcastAndOpBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
1150 {
1151   PetscErrorCode ierr;
1152 
1153   PetscFunctionBegin;
1154   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1155   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1156   ierr = PetscLogEventBegin(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
1157   ierr = (*sf->ops->BcastAndOpBegin)(sf,unit,rootdata,leafdata,op);CHKERRQ(ierr);
1158   ierr = PetscLogEventEnd(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
1159   PetscFunctionReturn(0);
1160 }
1161 
1162 /*@C
1163    PetscSFBcastAndOpEnd - end a broadcast & reduce operation started with PetscSFBcastAndOpBegin()
1164 
1165    Collective
1166 
1167    Input Arguments:
1168 +  sf - star forest
1169 .  unit - data type
1170 .  rootdata - buffer to broadcast
1171 -  op - operation to use for reduction
1172 
1173    Output Arguments:
1174 .  leafdata - buffer to be reduced with values from each leaf's respective root
1175 
1176    Level: intermediate
1177 
1178 .seealso: PetscSFSetGraph(), PetscSFReduceEnd()
1179 @*/
1180 PetscErrorCode PetscSFBcastAndOpEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
1181 {
1182   PetscErrorCode ierr;
1183 
1184   PetscFunctionBegin;
1185   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1186   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1187   ierr = PetscLogEventBegin(PETSCSF_BcastAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
1188   ierr = (*sf->ops->BcastAndOpEnd)(sf,unit,rootdata,leafdata,op);CHKERRQ(ierr);
1189   ierr = PetscLogEventEnd(PETSCSF_BcastAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
1190   PetscFunctionReturn(0);
1191 }
1192 
1193 /*@C
1194    PetscSFReduceBegin - begin reduction of leafdata into rootdata, to be completed with call to PetscSFReduceEnd()
1195 
1196    Collective
1197 
1198    Input Arguments:
1199 +  sf - star forest
1200 .  unit - data type
1201 .  leafdata - values to reduce
1202 -  op - reduction operation
1203 
1204    Output Arguments:
1205 .  rootdata - result of reduction of values from all leaves of each root
1206 
1207    Level: intermediate
1208 
1209 .seealso: PetscSFBcastBegin()
1210 @*/
1211 PetscErrorCode PetscSFReduceBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
1212 {
1213   PetscErrorCode ierr;
1214 
1215   PetscFunctionBegin;
1216   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1217   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1218   ierr = PetscLogEventBegin(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);
1219   ierr = (sf->ops->ReduceBegin)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
1220   ierr = PetscLogEventEnd(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);
1221   PetscFunctionReturn(0);
1222 }
1223 
1224 /*@C
1225    PetscSFReduceEnd - end a reduction operation started with PetscSFReduceBegin()
1226 
1227    Collective
1228 
1229    Input Arguments:
1230 +  sf - star forest
1231 .  unit - data type
1232 .  leafdata - values to reduce
1233 -  op - reduction operation
1234 
1235    Output Arguments:
1236 .  rootdata - result of reduction of values from all leaves of each root
1237 
1238    Level: intermediate
1239 
1240 .seealso: PetscSFSetGraph(), PetscSFBcastEnd()
1241 @*/
1242 PetscErrorCode PetscSFReduceEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
1243 {
1244   PetscErrorCode ierr;
1245 
1246   PetscFunctionBegin;
1247   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1248   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1249   ierr = PetscLogEventBegin(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);
1250   ierr = (*sf->ops->ReduceEnd)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
1251   ierr = PetscLogEventEnd(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);
1252   PetscFunctionReturn(0);
1253 }
1254 
1255 /*@C
1256    PetscSFComputeDegreeBegin - begin computation of degree for each root vertex, to be completed with PetscSFComputeDegreeEnd()
1257 
1258    Collective
1259 
1260    Input Arguments:
1261 .  sf - star forest
1262 
1263    Output Arguments:
1264 .  degree - degree of each root vertex
1265 
1266    Level: advanced
1267 
1268 .seealso: PetscSFGatherBegin()
1269 @*/
1270 PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf,const PetscInt **degree)
1271 {
1272   PetscErrorCode ierr;
1273 
1274   PetscFunctionBegin;
1275   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1276   PetscSFCheckGraphSet(sf,1);
1277   PetscValidPointer(degree,2);
1278   if (!sf->degreeknown) {
1279     PetscInt i, nroots = sf->nroots, maxlocal = sf->maxleaf+1;  /* TODO: We should use PetscSFGetLeafRange() */
1280     if (sf->degree) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Calls to PetscSFComputeDegreeBegin() cannot be nested.");
1281     ierr = PetscMalloc1(nroots,&sf->degree);CHKERRQ(ierr);
1282     ierr = PetscMalloc1(PetscMax(maxlocal,1),&sf->degreetmp);CHKERRQ(ierr); /* allocate at least one entry, see check in PetscSFComputeDegreeEnd() */
1283     for (i=0; i<nroots; i++) sf->degree[i] = 0;
1284     for (i=0; i<maxlocal; i++) sf->degreetmp[i] = 1;
1285     ierr = PetscSFReduceBegin(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);CHKERRQ(ierr);
1286   }
1287   *degree = NULL;
1288   PetscFunctionReturn(0);
1289 }
1290 
1291 /*@C
1292    PetscSFComputeDegreeEnd - complete computation of degree for each root vertex, started with PetscSFComputeDegreeBegin()
1293 
1294    Collective
1295 
1296    Input Arguments:
1297 .  sf - star forest
1298 
1299    Output Arguments:
1300 .  degree - degree of each root vertex
1301 
1302    Level: developer
1303 
1304 .seealso:
1305 @*/
1306 PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf,const PetscInt **degree)
1307 {
1308   PetscErrorCode ierr;
1309 
1310   PetscFunctionBegin;
1311   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1312   PetscSFCheckGraphSet(sf,1);
1313   PetscValidPointer(degree,2);
1314   if (!sf->degreeknown) {
1315     if (!sf->degreetmp) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFComputeDegreeBegin() before PetscSFComputeDegreeEnd()");
1316     ierr = PetscSFReduceEnd(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);CHKERRQ(ierr);
1317     ierr = PetscFree(sf->degreetmp);CHKERRQ(ierr);
1318     sf->degreeknown = PETSC_TRUE;
1319   }
1320   *degree = sf->degree;
1321   PetscFunctionReturn(0);
1322 }
1323 
1324 
1325 /*@C
1326    PetscSFComputeMultiRootOriginalNumbering - Returns original numbering of multi-roots (roots of multi-SF returned by PetscSFGetMultiSF()). Each multi-root is assigned index of its original root.
1327 
1328    Collective
1329 
1330    Input Arguments:
1331 +  sf - star forest
1332 -  degree - degree of each root vertex, computed with PetscSFComputeDegreeBegin()/PetscSFComputeDegreeEnd()
1333 
1334    Output Arguments:
1335 .  mRootsOrigNumbering - original indices of multi-roots; length of the array is equal to the number of multi-roots (roots of multi-SF)
1336 
1337    Level: developer
1338 
1339 .seealso: PetscSFComputeDegreeBegin(), PetscSFComputeDegreeEnd(), PetscSFGetMultiSF()
1340 @*/
1341 PetscErrorCode PetscSFComputeMultiRootOriginalNumbering(PetscSF sf, const PetscInt degree[], PetscInt *mRootsOrigNumbering[])
1342 {
1343   PetscSF             msf;
1344   PetscInt            i, j, k, nroots, nmroots;
1345   PetscErrorCode      ierr;
1346 
1347   PetscFunctionBegin;
1348   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1349   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
1350   if (nroots) PetscValidIntPointer(degree,2);
1351   PetscValidPointer(mRootsOrigNumbering,3);
1352   ierr = PetscSFGetMultiSF(sf,&msf);CHKERRQ(ierr);
1353   ierr = PetscSFGetGraph(msf, &nmroots, NULL, NULL, NULL);CHKERRQ(ierr);
1354   ierr = PetscMalloc1(nmroots, mRootsOrigNumbering);CHKERRQ(ierr);
1355   for (i=0,j=0,k=0; i<nroots; i++) {
1356     if (!degree[i]) continue;
1357     for (j=0; j<degree[i]; j++,k++) {
1358       (*mRootsOrigNumbering)[k] = i;
1359     }
1360   }
1361 #if defined(PETSC_USE_DEBUG)
1362   if (PetscUnlikely(k != nmroots)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sanity check fail");
1363 #endif
1364   PetscFunctionReturn(0);
1365 }
1366 
1367 /*@C
1368    PetscSFFetchAndOpBegin - begin operation that fetches values from root and updates atomically by applying operation using my leaf value, to be completed with PetscSFFetchAndOpEnd()
1369 
1370    Collective
1371 
1372    Input Arguments:
1373 +  sf - star forest
1374 .  unit - data type
1375 .  leafdata - leaf values to use in reduction
1376 -  op - operation to use for reduction
1377 
1378    Output Arguments:
1379 +  rootdata - root values to be updated, input state is seen by first process to perform an update
1380 -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
1381 
1382    Level: advanced
1383 
1384    Note:
1385    The update is only atomic at the granularity provided by the hardware. Different roots referenced by the same process
1386    might be updated in a different order. Furthermore, if a composite type is used for the unit datatype, atomicity is
1387    not guaranteed across the whole vertex. Therefore, this function is mostly only used with primitive types such as
1388    integers.
1389 
1390 .seealso: PetscSFComputeDegreeBegin(), PetscSFReduceBegin(), PetscSFSetGraph()
1391 @*/
1392 PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1393 {
1394   PetscErrorCode ierr;
1395 
1396   PetscFunctionBegin;
1397   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1398   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1399   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
1400   ierr = (*sf->ops->FetchAndOpBegin)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr);
1401   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
1402   PetscFunctionReturn(0);
1403 }
1404 
1405 /*@C
1406    PetscSFFetchAndOpEnd - end operation started in matching call to PetscSFFetchAndOpBegin() to fetch values from roots and update atomically by applying operation using my leaf value
1407 
1408    Collective
1409 
1410    Input Arguments:
1411 +  sf - star forest
1412 .  unit - data type
1413 .  leafdata - leaf values to use in reduction
1414 -  op - operation to use for reduction
1415 
1416    Output Arguments:
1417 +  rootdata - root values to be updated, input state is seen by first process to perform an update
1418 -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
1419 
1420    Level: advanced
1421 
1422 .seealso: PetscSFComputeDegreeEnd(), PetscSFReduceEnd(), PetscSFSetGraph()
1423 @*/
1424 PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1425 {
1426   PetscErrorCode ierr;
1427 
1428   PetscFunctionBegin;
1429   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1430   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1431   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
1432   ierr = (*sf->ops->FetchAndOpEnd)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr);
1433   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
1434   PetscFunctionReturn(0);
1435 }
1436 
1437 /*@C
1438    PetscSFGatherBegin - begin pointwise gather of all leaves into multi-roots, to be completed with PetscSFGatherEnd()
1439 
1440    Collective
1441 
1442    Input Arguments:
1443 +  sf - star forest
1444 .  unit - data type
1445 -  leafdata - leaf data to gather to roots
1446 
1447    Output Argument:
1448 .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
1449 
1450    Level: intermediate
1451 
1452 .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
1453 @*/
1454 PetscErrorCode PetscSFGatherBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
1455 {
1456   PetscErrorCode ierr;
1457   PetscSF        multi;
1458 
1459   PetscFunctionBegin;
1460   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1461   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1462   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
1463   ierr = PetscSFReduceBegin(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr);
1464   PetscFunctionReturn(0);
1465 }
1466 
1467 /*@C
1468    PetscSFGatherEnd - ends pointwise gather operation that was started with PetscSFGatherBegin()
1469 
1470    Collective
1471 
1472    Input Arguments:
1473 +  sf - star forest
1474 .  unit - data type
1475 -  leafdata - leaf data to gather to roots
1476 
1477    Output Argument:
1478 .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
1479 
1480    Level: intermediate
1481 
1482 .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
1483 @*/
1484 PetscErrorCode PetscSFGatherEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
1485 {
1486   PetscErrorCode ierr;
1487   PetscSF        multi;
1488 
1489   PetscFunctionBegin;
1490   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1491   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1492   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
1493   ierr = PetscSFReduceEnd(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr);
1494   PetscFunctionReturn(0);
1495 }
1496 
1497 /*@C
1498    PetscSFScatterBegin - begin pointwise scatter operation from multi-roots to leaves, to be completed with PetscSFScatterEnd()
1499 
1500    Collective
1501 
1502    Input Arguments:
1503 +  sf - star forest
1504 .  unit - data type
1505 -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
1506 
1507    Output Argument:
1508 .  leafdata - leaf data to be update with personal data from each respective root
1509 
1510    Level: intermediate
1511 
1512 .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
1513 @*/
1514 PetscErrorCode PetscSFScatterBegin(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
1515 {
1516   PetscErrorCode ierr;
1517   PetscSF        multi;
1518 
1519   PetscFunctionBegin;
1520   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1521   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1522   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
1523   ierr = PetscSFBcastBegin(multi,unit,multirootdata,leafdata);CHKERRQ(ierr);
1524   PetscFunctionReturn(0);
1525 }
1526 
1527 /*@C
1528    PetscSFScatterEnd - ends pointwise scatter operation that was started with PetscSFScatterBegin()
1529 
1530    Collective
1531 
1532    Input Arguments:
1533 +  sf - star forest
1534 .  unit - data type
1535 -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
1536 
1537    Output Argument:
1538 .  leafdata - leaf data to be update with personal data from each respective root
1539 
1540    Level: intermediate
1541 
1542 .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
1543 @*/
1544 PetscErrorCode PetscSFScatterEnd(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
1545 {
1546   PetscErrorCode ierr;
1547   PetscSF        multi;
1548 
1549   PetscFunctionBegin;
1550   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1551   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1552   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
1553   ierr = PetscSFBcastEnd(multi,unit,multirootdata,leafdata);CHKERRQ(ierr);
1554   PetscFunctionReturn(0);
1555 }
1556 
1557 /*@
1558   PetscSFCompose - Compose a new PetscSF equivalent to action to PetscSFs
1559 
1560   Input Parameters:
1561 + sfA - The first PetscSF
1562 - sfB - The second PetscSF
1563 
1564   Output Parameters:
1565 . sfBA - equvalent PetscSF for applying A then B
1566 
1567   Level: developer
1568 
1569 .seealso: PetscSF, PetscSFGetGraph(), PetscSFSetGraph()
1570 @*/
1571 PetscErrorCode PetscSFCompose(PetscSF sfA, PetscSF sfB, PetscSF *sfBA)
1572 {
1573   MPI_Comm           comm;
1574   const PetscSFNode *remotePointsA, *remotePointsB;
1575   PetscSFNode       *remotePointsBA;
1576   const PetscInt    *localPointsA, *localPointsB;
1577   PetscInt           numRootsA, numLeavesA, numRootsB, numLeavesB;
1578   PetscErrorCode     ierr;
1579 
1580   PetscFunctionBegin;
1581   PetscValidHeaderSpecific(sfA, PETSCSF_CLASSID, 1);
1582   PetscSFCheckGraphSet(sfA, 1);
1583   PetscValidHeaderSpecific(sfB, PETSCSF_CLASSID, 2);
1584   PetscSFCheckGraphSet(sfB, 2);
1585   PetscValidPointer(sfBA, 3);
1586   ierr = PetscObjectGetComm((PetscObject) sfA, &comm);CHKERRQ(ierr);
1587   ierr = PetscSFGetGraph(sfA, &numRootsA, &numLeavesA, &localPointsA, &remotePointsA);CHKERRQ(ierr);
1588   ierr = PetscSFGetGraph(sfB, &numRootsB, &numLeavesB, &localPointsB, &remotePointsB);CHKERRQ(ierr);
1589   ierr = PetscMalloc1(numLeavesB, &remotePointsBA);CHKERRQ(ierr);
1590   ierr = PetscSFBcastBegin(sfB, MPIU_2INT, remotePointsA, remotePointsBA);CHKERRQ(ierr);
1591   ierr = PetscSFBcastEnd(sfB, MPIU_2INT, remotePointsA, remotePointsBA);CHKERRQ(ierr);
1592   ierr = PetscSFCreate(comm, sfBA);CHKERRQ(ierr);
1593   ierr = PetscSFSetGraph(*sfBA, numRootsA, numLeavesB, localPointsB, PETSC_COPY_VALUES, remotePointsBA, PETSC_OWN_POINTER);CHKERRQ(ierr);
1594   PetscFunctionReturn(0);
1595 }
1596