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