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