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