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