Actual source code: sf.c
petsc-3.5.4 2015-05-23
1: #include <petsc-private/sfimpl.h> /*I "petscsf.h" I*/
2: #include <petscctable.h>
4: /* Logging support */
5: PetscLogEvent PETSCSF_SetGraph, PETSCSF_BcastBegin, PETSCSF_BcastEnd, PETSCSF_ReduceBegin, PETSCSF_ReduceEnd, PETSCSF_FetchAndOpBegin, PETSCSF_FetchAndOpEnd;
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
16: const char *const PetscSFDuplicateOptions[] = {"CONFONLY","RANKS","GRAPH","PetscSFDuplicateOption","PETSCSF_DUPLICATE_",0};
20: /*@C
21: PetscSFCreate - create a star forest communication context
23: Not Collective
25: Input Arguments:
26: . comm - communicator on which the star forest will operate
28: Output Arguments:
29: . sf - new star forest context
31: Level: intermediate
33: .seealso: PetscSFSetGraph(), PetscSFDestroy()
34: @*/
35: PetscErrorCode PetscSFCreate(MPI_Comm comm,PetscSF *sf)
36: {
38: PetscSF b;
42: PetscSFInitializePackage();
44: PetscHeaderCreate(b,_p_PetscSF,struct _PetscSFOps,PETSCSF_CLASSID,"PetscSF","Star Forest","PetscSF",comm,PetscSFDestroy,PetscSFView);
46: b->nroots = -1;
47: b->nleaves = -1;
48: b->nranks = -1;
49: b->rankorder = PETSC_TRUE;
50: b->ingroup = MPI_GROUP_NULL;
51: b->outgroup = MPI_GROUP_NULL;
52: b->graphset = PETSC_FALSE;
54: *sf = b;
55: return(0);
56: }
60: /*@C
61: PetscSFReset - Reset a star forest so that different sizes or neighbors can be used
63: Collective
65: Input Arguments:
66: . sf - star forest
68: Level: advanced
70: .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFDestroy()
71: @*/
72: PetscErrorCode PetscSFReset(PetscSF sf)
73: {
78: sf->mine = NULL;
79: PetscFree(sf->mine_alloc);
80: sf->remote = NULL;
81: PetscFree(sf->remote_alloc);
82: PetscFree4(sf->ranks,sf->roffset,sf->rmine,sf->rremote);
83: PetscFree(sf->degree);
84: if (sf->ingroup != MPI_GROUP_NULL) {MPI_Group_free(&sf->ingroup);}
85: if (sf->outgroup != MPI_GROUP_NULL) {MPI_Group_free(&sf->outgroup);}
86: PetscSFDestroy(&sf->multi);
87: sf->graphset = PETSC_FALSE;
88: if (sf->ops->Reset) {(*sf->ops->Reset)(sf);}
89: sf->setupcalled = PETSC_FALSE;
90: return(0);
91: }
95: /*@C
96: PetscSFSetType - set the PetscSF communication implementation
98: Collective on PetscSF
100: Input Parameters:
101: + sf - the PetscSF context
102: - type - a known method
104: Options Database Key:
105: . -sf_type <type> - Sets the method; use -help for a list
106: of available methods (for instance, window, pt2pt, neighbor)
108: Notes:
109: See "include/petscsf.h" for available methods (for instance)
110: + PETSCSFWINDOW - MPI-2/3 one-sided
111: - PETSCSFBASIC - basic implementation using MPI-1 two-sided
113: Level: intermediate
115: .keywords: PetscSF, set, type
117: .seealso: PetscSFType, PetscSFCreate()
118: @*/
119: PetscErrorCode PetscSFSetType(PetscSF sf,PetscSFType type)
120: {
121: PetscErrorCode ierr,(*r)(PetscSF);
122: PetscBool match;
128: PetscObjectTypeCompare((PetscObject)sf,type,&match);
129: if (match) return(0);
131: PetscFunctionListFind(PetscSFList,type,&r);
132: if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested PetscSF type %s",type);
133: /* Destroy the previous private PetscSF context */
134: if (sf->ops->Destroy) {
135: (*(sf)->ops->Destroy)(sf);
136: }
137: PetscMemzero(sf->ops,sizeof(*sf->ops));
138: PetscObjectChangeTypeName((PetscObject)sf,type);
139: (*r)(sf);
140: return(0);
141: }
145: /*@C
146: PetscSFDestroy - destroy star forest
148: Collective
150: Input Arguments:
151: . sf - address of star forest
153: Level: intermediate
155: .seealso: PetscSFCreate(), PetscSFReset()
156: @*/
157: PetscErrorCode PetscSFDestroy(PetscSF *sf)
158: {
162: if (!*sf) return(0);
164: if (--((PetscObject)(*sf))->refct > 0) {*sf = 0; return(0);}
165: PetscSFReset(*sf);
166: if ((*sf)->ops->Destroy) {(*(*sf)->ops->Destroy)(*sf);}
167: PetscHeaderDestroy(sf);
168: return(0);
169: }
173: /*@
174: PetscSFSetUp - set up communication structures
176: Collective
178: Input Arguments:
179: . sf - star forest communication object
181: Level: beginner
183: .seealso: PetscSFSetFromOptions(), PetscSFSetType()
184: @*/
185: PetscErrorCode PetscSFSetUp(PetscSF sf)
186: {
190: if (sf->setupcalled) return(0);
191: if (!((PetscObject)sf)->type_name) {PetscSFSetType(sf,PETSCSFBASIC);}
192: if (sf->ops->SetUp) {(*sf->ops->SetUp)(sf);}
193: sf->setupcalled = PETSC_TRUE;
194: return(0);
195: }
199: /*@C
200: PetscSFSetFromOptions - set PetscSF options using the options database
202: Logically Collective
204: Input Arguments:
205: . sf - star forest
207: Options Database Keys:
208: + -sf_type - implementation type, see PetscSFSetType()
209: - -sf_rank_order - sort composite points for gathers and scatters in rank order, gathers are non-deterministic otherwise
211: Level: intermediate
213: .keywords: KSP, set, from, options, database
215: .seealso: PetscSFWindowSetSyncType()
216: @*/
217: PetscErrorCode PetscSFSetFromOptions(PetscSF sf)
218: {
219: PetscSFType deft;
220: char type[256];
222: PetscBool flg;
226: PetscObjectOptionsBegin((PetscObject)sf);
227: deft = ((PetscObject)sf)->type_name ? ((PetscObject)sf)->type_name : PETSCSFBASIC;
228: PetscOptionsFList("-sf_type","PetscSF implementation type","PetscSFSetType",PetscSFList,deft,type,256,&flg);
229: PetscSFSetType(sf,flg ? type : deft);
230: 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);
231: if (sf->ops->SetFromOptions) {(*sf->ops->SetFromOptions)(sf);}
232: PetscOptionsEnd();
233: return(0);
234: }
238: /*@C
239: PetscSFSetRankOrder - sort multi-points for gathers and scatters by rank order
241: Logically Collective
243: Input Arguments:
244: + sf - star forest
245: - flg - PETSC_TRUE to sort, PETSC_FALSE to skip sorting (lower setup cost, but non-deterministic)
247: Level: advanced
249: .seealso: PetscSFGatherBegin(), PetscSFScatterBegin()
250: @*/
251: PetscErrorCode PetscSFSetRankOrder(PetscSF sf,PetscBool flg)
252: {
257: if (sf->multi) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Rank ordering must be set before first call to PetscSFGatherBegin() or PetscSFScatterBegin()");
258: sf->rankorder = flg;
259: return(0);
260: }
264: /*@C
265: PetscSFSetGraph - Set a parallel star forest
267: Collective
269: Input Arguments:
270: + sf - star forest
271: . nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
272: . nleaves - number of leaf vertices on the current process, each of these references a root on any process
273: . ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage
274: . localmode - copy mode for ilocal
275: . iremote - remote locations of root vertices for each leaf on the current process
276: - remotemode - copy mode for iremote
278: Level: intermediate
280: .seealso: PetscSFCreate(), PetscSFView(), PetscSFGetGraph()
281: @*/
282: PetscErrorCode PetscSFSetGraph(PetscSF sf,PetscInt nroots,PetscInt nleaves,const PetscInt *ilocal,PetscCopyMode localmode,const PetscSFNode *iremote,PetscCopyMode remotemode)
283: {
284: PetscErrorCode ierr;
285: PetscTable table;
286: PetscTablePosition pos;
287: PetscMPIInt size;
288: PetscInt i,*rcount,*ranks;
292: PetscLogEventBegin(PETSCSF_SetGraph,sf,0,0,0);
295: if (nroots < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"roots %D, cannot be negative",nroots);
296: if (nleaves < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nleaves %D, cannot be negative",nleaves);
297: PetscSFReset(sf);
298: sf->nroots = nroots;
299: sf->nleaves = nleaves;
300: if (ilocal) {
301: switch (localmode) {
302: case PETSC_COPY_VALUES:
303: PetscMalloc1(nleaves,&sf->mine_alloc);
304: sf->mine = sf->mine_alloc;
305: PetscMemcpy(sf->mine,ilocal,nleaves*sizeof(*sf->mine));
306: sf->minleaf = PETSC_MAX_INT;
307: sf->maxleaf = PETSC_MIN_INT;
308: for (i=0; i<nleaves; i++) {
309: sf->minleaf = PetscMin(sf->minleaf,ilocal[i]);
310: sf->maxleaf = PetscMax(sf->maxleaf,ilocal[i]);
311: }
312: break;
313: case PETSC_OWN_POINTER:
314: sf->mine_alloc = (PetscInt*)ilocal;
315: sf->mine = sf->mine_alloc;
316: break;
317: case PETSC_USE_POINTER:
318: sf->mine = (PetscInt*)ilocal;
319: break;
320: default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown localmode");
321: }
322: }
323: if (!ilocal || nleaves > 0) {
324: sf->minleaf = 0;
325: sf->maxleaf = nleaves - 1;
326: }
327: switch (remotemode) {
328: case PETSC_COPY_VALUES:
329: PetscMalloc1(nleaves,&sf->remote_alloc);
330: sf->remote = sf->remote_alloc;
331: PetscMemcpy(sf->remote,iremote,nleaves*sizeof(*sf->remote));
332: break;
333: case PETSC_OWN_POINTER:
334: sf->remote_alloc = (PetscSFNode*)iremote;
335: sf->remote = sf->remote_alloc;
336: break;
337: case PETSC_USE_POINTER:
338: sf->remote = (PetscSFNode*)iremote;
339: break;
340: default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown remotemode");
341: }
343: MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);
344: PetscTableCreate(10,size,&table);
345: for (i=0; i<nleaves; i++) {
346: /* Log 1-based rank */
347: PetscTableAdd(table,iremote[i].rank+1,1,ADD_VALUES);
348: }
349: PetscTableGetCount(table,&sf->nranks);
350: PetscMalloc4(sf->nranks,&sf->ranks,sf->nranks+1,&sf->roffset,nleaves,&sf->rmine,nleaves,&sf->rremote);
351: PetscMalloc2(sf->nranks,&rcount,sf->nranks,&ranks);
352: PetscTableGetHeadPosition(table,&pos);
353: for (i=0; i<sf->nranks; i++) {
354: PetscTableGetNext(table,&pos,&ranks[i],&rcount[i]);
355: ranks[i]--; /* Convert back to 0-based */
356: }
357: PetscTableDestroy(&table);
358: PetscSortIntWithArray(sf->nranks,ranks,rcount);
359: sf->roffset[0] = 0;
360: for (i=0; i<sf->nranks; i++) {
361: PetscMPIIntCast(ranks[i],sf->ranks+i);
362: sf->roffset[i+1] = sf->roffset[i] + rcount[i];
363: rcount[i] = 0;
364: }
365: for (i=0; i<nleaves; i++) {
366: PetscInt lo,hi,irank;
367: /* Search for index of iremote[i].rank in sf->ranks */
368: lo = 0; hi = sf->nranks;
369: while (hi - lo > 1) {
370: PetscInt mid = lo + (hi - lo)/2;
371: if (iremote[i].rank < sf->ranks[mid]) hi = mid;
372: else lo = mid;
373: }
374: if (hi - lo == 1 && iremote[i].rank == sf->ranks[lo]) irank = lo;
375: else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Could not find rank %D in array",iremote[i].rank);
376: sf->rmine[sf->roffset[irank] + rcount[irank]] = ilocal ? ilocal[i] : i;
377: sf->rremote[sf->roffset[irank] + rcount[irank]] = iremote[i].index;
378: rcount[irank]++;
379: }
380: PetscFree2(rcount,ranks);
381: #if !defined(PETSC_USE_64BIT_INDICES)
382: if (nroots == PETSC_DETERMINE) {
383: /* Jed, if you have a better way to do this, put it in */
384: PetscInt *numRankLeaves, *leafOff, *leafIndices, *numRankRoots, *rootOff, *rootIndices, maxRoots = 0;
386: /* All to all to determine number of leaf indices from each (you can do this using Scan and asynch messages) */
387: PetscMalloc4(size,&numRankLeaves,size+1,&leafOff,size,&numRankRoots,size+1,&rootOff);
388: PetscMemzero(numRankLeaves, size * sizeof(PetscInt));
389: for (i = 0; i < nleaves; ++i) ++numRankLeaves[iremote[i].rank];
390: MPI_Alltoall(numRankLeaves, 1, MPIU_INT, numRankRoots, 1, MPIU_INT, PetscObjectComm((PetscObject)sf));
391: /* Could set nroots to this maximum */
392: for (i = 0; i < size; ++i) maxRoots += numRankRoots[i];
394: /* Gather all indices */
395: PetscMalloc2(nleaves,&leafIndices,maxRoots,&rootIndices);
396: leafOff[0] = 0;
397: for (i = 0; i < size; ++i) leafOff[i+1] = leafOff[i] + numRankLeaves[i];
398: for (i = 0; i < nleaves; ++i) leafIndices[leafOff[iremote[i].rank]++] = iremote[i].index;
399: leafOff[0] = 0;
400: for (i = 0; i < size; ++i) leafOff[i+1] = leafOff[i] + numRankLeaves[i];
401: rootOff[0] = 0;
402: for (i = 0; i < size; ++i) rootOff[i+1] = rootOff[i] + numRankRoots[i];
403: MPI_Alltoallv(leafIndices, numRankLeaves, leafOff, MPIU_INT, rootIndices, numRankRoots, rootOff, MPIU_INT, PetscObjectComm((PetscObject)sf));
404: /* Sort and reduce */
405: PetscSortRemoveDupsInt(&maxRoots, rootIndices);
406: PetscFree2(leafIndices,rootIndices);
407: PetscFree4(numRankLeaves,leafOff,numRankRoots,rootOff);
408: sf->nroots = maxRoots;
409: }
410: #endif
412: sf->graphset = PETSC_TRUE;
413: PetscLogEventEnd(PETSCSF_SetGraph,sf,0,0,0);
414: return(0);
415: }
419: /*@C
420: PetscSFCreateInverseSF - given a PetscSF in which all vertices have degree 1, creates the inverse map
422: Collective
424: Input Arguments:
425: . sf - star forest to invert
427: Output Arguments:
428: . isf - inverse of sf
430: Level: advanced
432: Notes:
433: All roots must have degree 1.
435: The local space may be a permutation, but cannot be sparse.
437: .seealso: PetscSFSetGraph()
438: @*/
439: PetscErrorCode PetscSFCreateInverseSF(PetscSF sf,PetscSF *isf)
440: {
442: PetscMPIInt rank;
443: PetscInt i,nroots,nleaves,maxlocal,count,*newilocal;
444: const PetscInt *ilocal;
445: PetscSFNode *roots,*leaves;
448: MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);
449: PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);
450: for (i=0,maxlocal=0; i<nleaves; i++) maxlocal = PetscMax(maxlocal,(ilocal ? ilocal[i] : i)+1);
451: PetscMalloc2(nroots,&roots,nleaves,&leaves);
452: for (i=0; i<nleaves; i++) {
453: leaves[i].rank = rank;
454: leaves[i].index = i;
455: }
456: for (i=0; i <nroots; i++) {
457: roots[i].rank = -1;
458: roots[i].index = -1;
459: }
460: PetscSFReduceBegin(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);
461: PetscSFReduceEnd(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);
463: /* Check whether our leaves are sparse */
464: for (i=0,count=0; i<nroots; i++) if (roots[i].rank >= 0) count++;
465: if (count == nroots) newilocal = NULL;
466: else { /* Index for sparse leaves and compact "roots" array (which is to become our leaves). */
467: PetscMalloc1(count,&newilocal);
468: for (i=0,count=0; i<nroots; i++) {
469: if (roots[i].rank >= 0) {
470: newilocal[count] = i;
471: roots[count].rank = roots[i].rank;
472: roots[count].index = roots[i].index;
473: count++;
474: }
475: }
476: }
478: PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,isf);
479: PetscSFSetGraph(*isf,maxlocal,count,newilocal,PETSC_OWN_POINTER,roots,PETSC_COPY_VALUES);
480: PetscFree2(roots,leaves);
481: return(0);
482: }
486: /*@
487: PetscSFDuplicate - duplicate a PetscSF, optionally preserving rank connectivity and graph
489: Collective
491: Input Arguments:
492: + sf - communication object to duplicate
493: - opt - PETSCSF_DUPLICATE_CONFONLY, PETSCSF_DUPLICATE_RANKS, or PETSCSF_DUPLICATE_GRAPH (see PetscSFDuplicateOption)
495: Output Arguments:
496: . newsf - new communication object
498: Level: beginner
500: .seealso: PetscSFCreate(), PetscSFSetType(), PetscSFSetGraph()
501: @*/
502: PetscErrorCode PetscSFDuplicate(PetscSF sf,PetscSFDuplicateOption opt,PetscSF *newsf)
503: {
507: PetscSFCreate(PetscObjectComm((PetscObject)sf),newsf);
508: PetscSFSetType(*newsf,((PetscObject)sf)->type_name);
509: if (sf->ops->Duplicate) {(*sf->ops->Duplicate)(sf,opt,*newsf);}
510: if (opt == PETSCSF_DUPLICATE_GRAPH) {
511: PetscInt nroots,nleaves;
512: const PetscInt *ilocal;
513: const PetscSFNode *iremote;
514: PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);
515: PetscSFSetGraph(*newsf,nroots,nleaves,ilocal,PETSC_COPY_VALUES,iremote,PETSC_COPY_VALUES);
516: }
517: return(0);
518: }
522: /*@C
523: PetscSFGetGraph - Get the graph specifying a parallel star forest
525: Not Collective
527: Input Arguments:
528: . sf - star forest
530: Output Arguments:
531: + nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
532: . nleaves - number of leaf vertices on the current process, each of these references a root on any process
533: . ilocal - locations of leaves in leafdata buffers
534: - iremote - remote locations of root vertices for each leaf on the current process
536: Level: intermediate
538: .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph()
539: @*/
540: PetscErrorCode PetscSFGetGraph(PetscSF sf,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote)
541: {
545: /* We are not currently requiring that the graph is set, thus returning nroots=-1 if it has not been set */
546: /* if (!sf->graphset) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Graph has not been set, must call PetscSFSetGraph()"); */
547: if (nroots) *nroots = sf->nroots;
548: if (nleaves) *nleaves = sf->nleaves;
549: if (ilocal) *ilocal = sf->mine;
550: if (iremote) *iremote = sf->remote;
551: return(0);
552: }
556: /*@C
557: PetscSFGetLeafRange - Get the active leaf ranges
559: Not Collective
561: Input Arguments:
562: . sf - star forest
564: Output Arguments:
565: + minleaf - minimum active leaf on this process
566: - maxleaf - maximum active leaf on this process
568: Level: developer
570: .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph(), PetscSFGetGraph()
571: @*/
572: PetscErrorCode PetscSFGetLeafRange(PetscSF sf,PetscInt *minleaf,PetscInt *maxleaf)
573: {
577: if (minleaf) *minleaf = sf->minleaf;
578: if (maxleaf) *maxleaf = sf->maxleaf;
579: return(0);
580: }
584: /*@C
585: PetscSFView - view a star forest
587: Collective
589: Input Arguments:
590: + sf - star forest
591: - viewer - viewer to display graph, for example PETSC_VIEWER_STDOUT_WORLD
593: Level: beginner
595: .seealso: PetscSFCreate(), PetscSFSetGraph()
596: @*/
597: PetscErrorCode PetscSFView(PetscSF sf,PetscViewer viewer)
598: {
599: PetscErrorCode ierr;
600: PetscBool iascii;
601: PetscViewerFormat format;
605: if (!viewer) {PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sf),&viewer);}
608: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
609: if (iascii) {
610: PetscMPIInt rank;
611: PetscInt i,j;
613: PetscObjectPrintClassNamePrefixType((PetscObject)sf,viewer);
614: PetscViewerASCIIPushTab(viewer);
615: if (sf->ops->View) {(*sf->ops->View)(sf,viewer);}
616: MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);
617: PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);
618: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of roots=%D, leaves=%D, remote ranks=%D\n",rank,sf->nroots,sf->nleaves,sf->nranks);
619: for (i=0; i<sf->nleaves; i++) {
620: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D <- (%D,%D)\n",rank,sf->mine ? sf->mine[i] : i,sf->remote[i].rank,sf->remote[i].index);
621: }
622: PetscViewerFlush(viewer);
623: PetscViewerGetFormat(viewer,&format);
624: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
625: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Roots referenced by my leaves, by rank\n",rank);
626: for (i=0; i<sf->nranks; i++) {
627: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %d: %D edges\n",rank,sf->ranks[i],sf->roffset[i+1]-sf->roffset[i]);
628: for (j=sf->roffset[i]; j<sf->roffset[i+1]; j++) {
629: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D <- %D\n",rank,sf->rmine[j],sf->rremote[j]);
630: }
631: }
632: }
633: PetscViewerFlush(viewer);
634: PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);
635: PetscViewerASCIIPopTab(viewer);
636: }
637: return(0);
638: }
642: /*@C
643: PetscSFGetRanks - Get ranks and number of vertices referenced by leaves on this process
645: Not Collective
647: Input Arguments:
648: . sf - star forest
650: Output Arguments:
651: + nranks - number of ranks referenced by local part
652: . ranks - array of ranks
653: . roffset - offset in rmine/rremote for each rank (length nranks+1)
654: . rmine - concatenated array holding local indices referencing each remote rank
655: - rremote - concatenated array holding remote indices referenced for each remote rank
657: Level: developer
659: .seealso: PetscSFSetGraph()
660: @*/
661: PetscErrorCode PetscSFGetRanks(PetscSF sf,PetscInt *nranks,const PetscMPIInt **ranks,const PetscInt **roffset,const PetscInt **rmine,const PetscInt **rremote)
662: {
666: if (nranks) *nranks = sf->nranks;
667: if (ranks) *ranks = sf->ranks;
668: if (roffset) *roffset = sf->roffset;
669: if (rmine) *rmine = sf->rmine;
670: if (rremote) *rremote = sf->rremote;
671: return(0);
672: }
676: /*@C
677: PetscSFGetGroups - gets incoming and outgoing process groups
679: Collective
681: Input Argument:
682: . sf - star forest
684: Output Arguments:
685: + incoming - group of origin processes for incoming edges (leaves that reference my roots)
686: - outgoing - group of destination processes for outgoing edges (roots that I reference)
688: Level: developer
690: .seealso: PetscSFGetWindow(), PetscSFRestoreWindow()
691: @*/
692: PetscErrorCode PetscSFGetGroups(PetscSF sf,MPI_Group *incoming,MPI_Group *outgoing)
693: {
695: MPI_Group group;
698: if (sf->ingroup == MPI_GROUP_NULL) {
699: PetscInt i;
700: const PetscInt *indegree;
701: PetscMPIInt rank,*outranks,*inranks;
702: PetscSFNode *remote;
703: PetscSF bgcount;
705: /* Compute the number of incoming ranks */
706: PetscMalloc1(sf->nranks,&remote);
707: for (i=0; i<sf->nranks; i++) {
708: remote[i].rank = sf->ranks[i];
709: remote[i].index = 0;
710: }
711: PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,&bgcount);
712: PetscSFSetGraph(bgcount,1,sf->nranks,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);
713: PetscSFComputeDegreeBegin(bgcount,&indegree);
714: PetscSFComputeDegreeEnd(bgcount,&indegree);
716: /* Enumerate the incoming ranks */
717: PetscMalloc2(indegree[0],&inranks,sf->nranks,&outranks);
718: MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);
719: for (i=0; i<sf->nranks; i++) outranks[i] = rank;
720: PetscSFGatherBegin(bgcount,MPI_INT,outranks,inranks);
721: PetscSFGatherEnd(bgcount,MPI_INT,outranks,inranks);
722: MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);
723: MPI_Group_incl(group,indegree[0],inranks,&sf->ingroup);
724: MPI_Group_free(&group);
725: PetscFree2(inranks,outranks);
726: PetscSFDestroy(&bgcount);
727: }
728: *incoming = sf->ingroup;
730: if (sf->outgroup == MPI_GROUP_NULL) {
731: MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);
732: MPI_Group_incl(group,sf->nranks,sf->ranks,&sf->outgroup);
733: MPI_Group_free(&group);
734: }
735: *outgoing = sf->outgroup;
736: return(0);
737: }
741: /*@C
742: PetscSFGetMultiSF - gets the inner SF implemeting gathers and scatters
744: Collective
746: Input Argument:
747: . sf - star forest that may contain roots with 0 or with more than 1 vertex
749: Output Arguments:
750: . multi - star forest with split roots, such that each root has degree exactly 1
752: Level: developer
754: Notes:
756: In most cases, users should use PetscSFGatherBegin() and PetscSFScatterBegin() instead of manipulating multi
757: directly. Since multi satisfies the stronger condition that each entry in the global space has exactly one incoming
758: edge, it is a candidate for future optimization that might involve its removal.
760: .seealso: PetscSFSetGraph(), PetscSFGatherBegin(), PetscSFScatterBegin()
761: @*/
762: PetscErrorCode PetscSFGetMultiSF(PetscSF sf,PetscSF *multi)
763: {
769: if (sf->nroots < 0) { /* Graph has not been set yet; why do we need this? */
770: PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);
771: *multi = sf->multi;
772: return(0);
773: }
774: if (!sf->multi) {
775: const PetscInt *indegree;
776: PetscInt i,*inoffset,*outones,*outoffset;
777: PetscSFNode *remote;
778: PetscSFComputeDegreeBegin(sf,&indegree);
779: PetscSFComputeDegreeEnd(sf,&indegree);
780: PetscMalloc3(sf->nroots+1,&inoffset,sf->nleaves,&outones,sf->nleaves,&outoffset);
781: inoffset[0] = 0;
782: #if 1
783: for (i=0; i<sf->nroots; i++) inoffset[i+1] = PetscMax(i+1, inoffset[i] + indegree[i]);
784: #else
785: for (i=0; i<sf->nroots; i++) inoffset[i+1] = inoffset[i] + indegree[i];
786: #endif
787: for (i=0; i<sf->nleaves; i++) outones[i] = 1;
788: PetscSFFetchAndOpBegin(sf,MPIU_INT,inoffset,outones,outoffset,MPIU_SUM);
789: PetscSFFetchAndOpEnd(sf,MPIU_INT,inoffset,outones,outoffset,MPIU_SUM);
790: for (i=0; i<sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */
791: #if 0
792: #if defined(PETSC_USE_DEBUG) /* Check that the expected number of increments occurred */
793: for (i=0; i<sf->nroots; i++) {
794: if (inoffset[i] + indegree[i] != inoffset[i+1]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect result after PetscSFFetchAndOp");
795: }
796: #endif
797: #endif
798: PetscMalloc1(sf->nleaves,&remote);
799: for (i=0; i<sf->nleaves; i++) {
800: remote[i].rank = sf->remote[i].rank;
801: remote[i].index = outoffset[i];
802: }
803: PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);
804: PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);
805: if (sf->rankorder) { /* Sort the ranks */
806: PetscMPIInt rank;
807: PetscInt *inranks,*newoffset,*outranks,*newoutoffset,*tmpoffset,maxdegree;
808: PetscSFNode *newremote;
809: MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);
810: for (i=0,maxdegree=0; i<sf->nroots; i++) maxdegree = PetscMax(maxdegree,indegree[i]);
811: PetscMalloc5(sf->multi->nroots,&inranks,sf->multi->nroots,&newoffset,sf->nleaves,&outranks,sf->nleaves,&newoutoffset,maxdegree,&tmpoffset);
812: for (i=0; i<sf->nleaves; i++) outranks[i] = rank;
813: PetscSFReduceBegin(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);
814: PetscSFReduceEnd(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);
815: /* Sort the incoming ranks at each vertex, build the inverse map */
816: for (i=0; i<sf->nroots; i++) {
817: PetscInt j;
818: for (j=0; j<indegree[i]; j++) tmpoffset[j] = j;
819: PetscSortIntWithArray(indegree[i],inranks+inoffset[i],tmpoffset);
820: for (j=0; j<indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j;
821: }
822: PetscSFBcastBegin(sf->multi,MPIU_INT,newoffset,newoutoffset);
823: PetscSFBcastEnd(sf->multi,MPIU_INT,newoffset,newoutoffset);
824: PetscMalloc1(sf->nleaves,&newremote);
825: for (i=0; i<sf->nleaves; i++) {
826: newremote[i].rank = sf->remote[i].rank;
827: newremote[i].index = newoutoffset[i];
828: }
829: PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,NULL,PETSC_COPY_VALUES,newremote,PETSC_OWN_POINTER);
830: PetscFree5(inranks,newoffset,outranks,newoutoffset,tmpoffset);
831: }
832: PetscFree3(inoffset,outones,outoffset);
833: }
834: *multi = sf->multi;
835: return(0);
836: }
840: /*@C
841: PetscSFCreateEmbeddedSF - removes edges from all but the selected roots, does not remap indices
843: Collective
845: Input Arguments:
846: + sf - original star forest
847: . nroots - number of roots to select on this process
848: - selected - selected roots on this process
850: Output Arguments:
851: . newsf - new star forest
853: Level: advanced
855: Note:
856: To use the new PetscSF, it may be necessary to know the indices of the leaves that are still participating. This can
857: be done by calling PetscSFGetGraph().
859: .seealso: PetscSFSetGraph(), PetscSFGetGraph()
860: @*/
861: PetscErrorCode PetscSFCreateEmbeddedSF(PetscSF sf,PetscInt nroots,const PetscInt *selected,PetscSF *newsf)
862: {
863: PetscInt *rootdata, *leafdata, *ilocal;
864: PetscSFNode *iremote;
865: PetscInt leafsize = 0, nleaves = 0, n, i;
872: if (sf->mine) for (i = 0; i < sf->nleaves; ++i) {leafsize = PetscMax(leafsize, sf->mine[i]+1);}
873: else leafsize = sf->nleaves;
874: PetscCalloc2(sf->nroots,&rootdata,leafsize,&leafdata);
875: for (i=0; i<nroots; ++i) rootdata[selected[i]] = 1;
876: PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata);
877: PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata);
879: for (i = 0; i < leafsize; ++i) nleaves += leafdata[i];
880: PetscMalloc1(nleaves,&ilocal);
881: PetscMalloc1(nleaves,&iremote);
882: for (i = 0, n = 0; i < sf->nleaves; ++i) {
883: const PetscInt lidx = sf->mine ? sf->mine[i] : i;
885: if (leafdata[lidx]) {
886: ilocal[n] = lidx;
887: iremote[n].rank = sf->remote[i].rank;
888: iremote[n].index = sf->remote[i].index;
889: ++n;
890: }
891: }
892: if (n != nleaves) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "There is a size mismatch in the SF embedding, %d != %d", n, nleaves);
893: PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,newsf);
894: PetscSFSetGraph(*newsf,sf->nroots,nleaves,ilocal,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);
895: PetscFree2(rootdata,leafdata);
896: return(0);
897: }
901: /*@C
902: PetscSFBcastBegin - begin pointwise broadcast to be concluded with call to PetscSFBcastEnd()
904: Collective on PetscSF
906: Input Arguments:
907: + sf - star forest on which to communicate
908: . unit - data type associated with each node
909: - rootdata - buffer to broadcast
911: Output Arguments:
912: . leafdata - buffer to update with values from each leaf's respective root
914: Level: intermediate
916: .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFView(), PetscSFBcastEnd(), PetscSFReduceBegin()
917: @*/
918: PetscErrorCode PetscSFBcastBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
919: {
924: PetscSFCheckGraphSet(sf,1);
925: PetscLogEventBegin(PETSCSF_BcastBegin,sf,0,0,0);
926: PetscSFSetUp(sf);
927: (*sf->ops->BcastBegin)(sf,unit,rootdata,leafdata);
928: PetscLogEventEnd(PETSCSF_BcastBegin,sf,0,0,0);
929: return(0);
930: }
934: /*@C
935: PetscSFBcastEnd - end a broadcast operation started with PetscSFBcastBegin()
937: Collective
939: Input Arguments:
940: + sf - star forest
941: . unit - data type
942: - rootdata - buffer to broadcast
944: Output Arguments:
945: . leafdata - buffer to update with values from each leaf's respective root
947: Level: intermediate
949: .seealso: PetscSFSetGraph(), PetscSFReduceEnd()
950: @*/
951: PetscErrorCode PetscSFBcastEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
952: {
957: PetscSFCheckGraphSet(sf,1);
958: PetscLogEventBegin(PETSCSF_BcastEnd,sf,0,0,0);
959: PetscSFSetUp(sf);
960: (*sf->ops->BcastEnd)(sf,unit,rootdata,leafdata);
961: PetscLogEventEnd(PETSCSF_BcastEnd,sf,0,0,0);
962: return(0);
963: }
967: /*@C
968: PetscSFReduceBegin - begin reduction of leafdata into rootdata, to be completed with call to PetscSFReduceEnd()
970: Collective
972: Input Arguments:
973: + sf - star forest
974: . unit - data type
975: . leafdata - values to reduce
976: - op - reduction operation
978: Output Arguments:
979: . rootdata - result of reduction of values from all leaves of each root
981: Level: intermediate
983: .seealso: PetscSFBcastBegin()
984: @*/
985: PetscErrorCode PetscSFReduceBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
986: {
991: PetscSFCheckGraphSet(sf,1);
992: PetscLogEventBegin(PETSCSF_ReduceBegin,sf,0,0,0);
993: PetscSFSetUp(sf);
994: (sf->ops->ReduceBegin)(sf,unit,leafdata,rootdata,op);
995: PetscLogEventEnd(PETSCSF_ReduceBegin,sf,0,0,0);
996: return(0);
997: }
1001: /*@C
1002: PetscSFReduceEnd - end a reduction operation started with PetscSFReduceBegin()
1004: Collective
1006: Input Arguments:
1007: + sf - star forest
1008: . unit - data type
1009: . leafdata - values to reduce
1010: - op - reduction operation
1012: Output Arguments:
1013: . rootdata - result of reduction of values from all leaves of each root
1015: Level: intermediate
1017: .seealso: PetscSFSetGraph(), PetscSFBcastEnd()
1018: @*/
1019: PetscErrorCode PetscSFReduceEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
1020: {
1025: PetscSFCheckGraphSet(sf,1);
1026: PetscLogEventBegin(PETSCSF_ReduceEnd,sf,0,0,0);
1027: PetscSFSetUp(sf);
1028: (*sf->ops->ReduceEnd)(sf,unit,leafdata,rootdata,op);
1029: PetscLogEventEnd(PETSCSF_ReduceEnd,sf,0,0,0);
1030: return(0);
1031: }
1035: /*@C
1036: PetscSFComputeDegreeBegin - begin computation of degree for each root vertex, to be completed with PetscSFComputeDegreeEnd()
1038: Collective
1040: Input Arguments:
1041: . sf - star forest
1043: Output Arguments:
1044: . degree - degree of each root vertex
1046: Level: advanced
1048: .seealso: PetscSFGatherBegin()
1049: @*/
1050: PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf,const PetscInt **degree)
1051: {
1056: PetscSFCheckGraphSet(sf,1);
1058: if (!sf->degree) {
1059: PetscInt i;
1060: PetscMalloc1(sf->nroots,&sf->degree);
1061: PetscMalloc1(sf->nleaves,&sf->degreetmp);
1062: for (i=0; i<sf->nroots; i++) sf->degree[i] = 0;
1063: for (i=0; i<sf->nleaves; i++) sf->degreetmp[i] = 1;
1064: PetscSFReduceBegin(sf,MPIU_INT,sf->degreetmp,sf->degree,MPIU_SUM);
1065: }
1066: *degree = NULL;
1067: return(0);
1068: }
1072: /*@C
1073: PetscSFComputeDegreeEnd - complete computation of degree for each root vertex, started with PetscSFComputeDegreeBegin()
1075: Collective
1077: Input Arguments:
1078: . sf - star forest
1080: Output Arguments:
1081: . degree - degree of each root vertex
1083: Level: developer
1085: .seealso:
1086: @*/
1087: PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf,const PetscInt **degree)
1088: {
1093: PetscSFCheckGraphSet(sf,1);
1094: if (!sf->degreeknown) {
1095: PetscSFReduceEnd(sf,MPIU_INT,sf->degreetmp,sf->degree,MPIU_SUM);
1096: PetscFree(sf->degreetmp);
1098: sf->degreeknown = PETSC_TRUE;
1099: }
1100: *degree = sf->degree;
1101: return(0);
1102: }
1106: /*@C
1107: PetscSFFetchAndOpBegin - begin operation that fetches values from root and updates atomically by applying operation using my leaf value, to be completed with PetscSFFetchAndOpEnd()
1109: Collective
1111: Input Arguments:
1112: + sf - star forest
1113: . unit - data type
1114: . leafdata - leaf values to use in reduction
1115: - op - operation to use for reduction
1117: Output Arguments:
1118: + rootdata - root values to be updated, input state is seen by first process to perform an update
1119: - leafupdate - state at each leaf's respective root immediately prior to my atomic update
1121: Level: advanced
1123: Note:
1124: The update is only atomic at the granularity provided by the hardware. Different roots referenced by the same process
1125: might be updated in a different order. Furthermore, if a composite type is used for the unit datatype, atomicity is
1126: not guaranteed across the whole vertex. Therefore, this function is mostly only used with primitive types such as
1127: integers.
1129: .seealso: PetscSFComputeDegreeBegin(), PetscSFReduceBegin(), PetscSFSetGraph()
1130: @*/
1131: PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1132: {
1137: PetscSFCheckGraphSet(sf,1);
1138: PetscLogEventBegin(PETSCSF_FetchAndOpBegin,sf,0,0,0);
1139: PetscSFSetUp(sf);
1140: (*sf->ops->FetchAndOpBegin)(sf,unit,rootdata,leafdata,leafupdate,op);
1141: PetscLogEventEnd(PETSCSF_FetchAndOpBegin,sf,0,0,0);
1142: return(0);
1143: }
1147: /*@C
1148: PetscSFFetchAndOpEnd - end operation started in matching call to PetscSFFetchAndOpBegin() to fetch values from roots and update atomically by applying operation using my leaf value
1150: Collective
1152: Input Arguments:
1153: + sf - star forest
1154: . unit - data type
1155: . leafdata - leaf values to use in reduction
1156: - op - operation to use for reduction
1158: Output Arguments:
1159: + rootdata - root values to be updated, input state is seen by first process to perform an update
1160: - leafupdate - state at each leaf's respective root immediately prior to my atomic update
1162: Level: advanced
1164: .seealso: PetscSFComputeDegreeEnd(), PetscSFReduceEnd(), PetscSFSetGraph()
1165: @*/
1166: PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1167: {
1172: PetscSFCheckGraphSet(sf,1);
1173: PetscLogEventBegin(PETSCSF_FetchAndOpEnd,sf,0,0,0);
1174: PetscSFSetUp(sf);
1175: (*sf->ops->FetchAndOpEnd)(sf,unit,rootdata,leafdata,leafupdate,op);
1176: PetscLogEventEnd(PETSCSF_FetchAndOpEnd,sf,0,0,0);
1177: return(0);
1178: }
1182: /*@C
1183: PetscSFGatherBegin - begin pointwise gather of all leaves into multi-roots, to be completed with PetscSFGatherEnd()
1185: Collective
1187: Input Arguments:
1188: + sf - star forest
1189: . unit - data type
1190: - leafdata - leaf data to gather to roots
1192: Output Argument:
1193: . multirootdata - root buffer to gather into, amount of space per root is equal to its degree
1195: Level: intermediate
1197: .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
1198: @*/
1199: PetscErrorCode PetscSFGatherBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
1200: {
1202: PetscSF multi;
1206: PetscSFGetMultiSF(sf,&multi);
1207: PetscSFReduceBegin(multi,unit,leafdata,multirootdata,MPIU_REPLACE);
1208: return(0);
1209: }
1213: /*@C
1214: PetscSFGatherEnd - ends pointwise gather operation that was started with PetscSFGatherBegin()
1216: Collective
1218: Input Arguments:
1219: + sf - star forest
1220: . unit - data type
1221: - leafdata - leaf data to gather to roots
1223: Output Argument:
1224: . multirootdata - root buffer to gather into, amount of space per root is equal to its degree
1226: Level: intermediate
1228: .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
1229: @*/
1230: PetscErrorCode PetscSFGatherEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
1231: {
1233: PetscSF multi;
1237: PetscSFCheckGraphSet(sf,1);
1238: PetscSFSetUp(sf);
1239: PetscSFGetMultiSF(sf,&multi);
1240: PetscSFReduceEnd(multi,unit,leafdata,multirootdata,MPIU_REPLACE);
1241: return(0);
1242: }
1246: /*@C
1247: PetscSFScatterBegin - begin pointwise scatter operation from multi-roots to leaves, to be completed with PetscSFScatterEnd()
1249: Collective
1251: Input Arguments:
1252: + sf - star forest
1253: . unit - data type
1254: - multirootdata - root buffer to send to each leaf, one unit of data per leaf
1256: Output Argument:
1257: . leafdata - leaf data to be update with personal data from each respective root
1259: Level: intermediate
1261: .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
1262: @*/
1263: PetscErrorCode PetscSFScatterBegin(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
1264: {
1266: PetscSF multi;
1270: PetscSFCheckGraphSet(sf,1);
1271: PetscSFSetUp(sf);
1272: PetscSFGetMultiSF(sf,&multi);
1273: PetscSFBcastBegin(multi,unit,multirootdata,leafdata);
1274: return(0);
1275: }
1279: /*@C
1280: PetscSFScatterEnd - ends pointwise scatter operation that was started with PetscSFScatterBegin()
1282: Collective
1284: Input Arguments:
1285: + sf - star forest
1286: . unit - data type
1287: - multirootdata - root buffer to send to each leaf, one unit of data per leaf
1289: Output Argument:
1290: . leafdata - leaf data to be update with personal data from each respective root
1292: Level: intermediate
1294: .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
1295: @*/
1296: PetscErrorCode PetscSFScatterEnd(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
1297: {
1299: PetscSF multi;
1303: PetscSFCheckGraphSet(sf,1);
1304: PetscSFSetUp(sf);
1305: PetscSFGetMultiSF(sf,&multi);
1306: PetscSFBcastEnd(multi,unit,multirootdata,leafdata);
1307: return(0);
1308: }