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