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