xref: /petsc/src/vec/is/sf/interface/sf.c (revision dbd2ff41708ad711a91e64cc194ae80bfcaa9b20)
195fce210SBarry Smith #include <petsc-private/sfimpl.h> /*I "petscsf.h" I*/
295fce210SBarry Smith #include <petscctable.h>
395fce210SBarry Smith 
4563ffbabSMatthew G. Knepley /* Logging support */
5563ffbabSMatthew G. Knepley PetscLogEvent PETSCSF_SetGraph, PETSCSF_BcastBegin, PETSCSF_BcastEnd, PETSCSF_ReduceBegin, PETSCSF_ReduceEnd, PETSCSF_FetchAndOpBegin, PETSCSF_FetchAndOpEnd;
6563ffbabSMatthew G. Knepley 
795fce210SBarry Smith #if defined(PETSC_USE_DEBUG)
895fce210SBarry Smith #  define PetscSFCheckGraphSet(sf,arg) do {                          \
995fce210SBarry Smith     if (PetscUnlikely(!(sf)->graphset))                              \
1095fce210SBarry Smith       SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetGraph() on argument %D \"%s\" before %s()",(arg),#sf,PETSC_FUNCTION_NAME); \
1195fce210SBarry Smith   } while (0)
1295fce210SBarry Smith #else
1395fce210SBarry Smith #  define PetscSFCheckGraphSet(sf,arg) do {} while (0)
1495fce210SBarry Smith #endif
1595fce210SBarry Smith 
1695fce210SBarry Smith const char *const PetscSFDuplicateOptions[] = {"CONFONLY","RANKS","GRAPH","PetscSFDuplicateOption","PETSCSF_DUPLICATE_",0};
1795fce210SBarry Smith 
1895fce210SBarry Smith #undef __FUNCT__
1995fce210SBarry Smith #define __FUNCT__ "PetscSFCreate"
2095fce210SBarry Smith /*@C
2195fce210SBarry Smith    PetscSFCreate - create a star forest communication context
2295fce210SBarry Smith 
2395fce210SBarry Smith    Not Collective
2495fce210SBarry Smith 
2595fce210SBarry Smith    Input Arguments:
2695fce210SBarry Smith .  comm - communicator on which the star forest will operate
2795fce210SBarry Smith 
2895fce210SBarry Smith    Output Arguments:
2995fce210SBarry Smith .  sf - new star forest context
3095fce210SBarry Smith 
3195fce210SBarry Smith    Level: intermediate
3295fce210SBarry Smith 
3395fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFDestroy()
3495fce210SBarry Smith @*/
3595fce210SBarry Smith PetscErrorCode PetscSFCreate(MPI_Comm comm,PetscSF *sf)
3695fce210SBarry Smith {
3795fce210SBarry Smith   PetscErrorCode ierr;
3895fce210SBarry Smith   PetscSF        b;
3995fce210SBarry Smith 
4095fce210SBarry Smith   PetscFunctionBegin;
4195fce210SBarry Smith   PetscValidPointer(sf,2);
42607a6623SBarry Smith   ierr = PetscSFInitializePackage();CHKERRQ(ierr);
4395fce210SBarry Smith 
4495fce210SBarry Smith   ierr = PetscHeaderCreate(b,_p_PetscSF,struct _PetscSFOps,PETSCSF_CLASSID,"PetscSF","Star Forest","PetscSF",comm,PetscSFDestroy,PetscSFView);CHKERRQ(ierr);
4595fce210SBarry Smith 
4695fce210SBarry Smith   b->nroots    = -1;
4795fce210SBarry Smith   b->nleaves   = -1;
4895fce210SBarry Smith   b->nranks    = -1;
4995fce210SBarry Smith   b->rankorder = PETSC_TRUE;
5095fce210SBarry Smith   b->ingroup   = MPI_GROUP_NULL;
5195fce210SBarry Smith   b->outgroup  = MPI_GROUP_NULL;
5295fce210SBarry Smith   b->graphset  = PETSC_FALSE;
5395fce210SBarry Smith 
5495fce210SBarry Smith   *sf = b;
5595fce210SBarry Smith   PetscFunctionReturn(0);
5695fce210SBarry Smith }
5795fce210SBarry Smith 
5895fce210SBarry Smith #undef __FUNCT__
5995fce210SBarry Smith #define __FUNCT__ "PetscSFReset"
6095fce210SBarry Smith /*@C
6195fce210SBarry Smith    PetscSFReset - Reset a star forest so that different sizes or neighbors can be used
6295fce210SBarry Smith 
6395fce210SBarry Smith    Collective
6495fce210SBarry Smith 
6595fce210SBarry Smith    Input Arguments:
6695fce210SBarry Smith .  sf - star forest
6795fce210SBarry Smith 
6895fce210SBarry Smith    Level: advanced
6995fce210SBarry Smith 
7095fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFDestroy()
7195fce210SBarry Smith @*/
7295fce210SBarry Smith PetscErrorCode PetscSFReset(PetscSF sf)
7395fce210SBarry Smith {
7495fce210SBarry Smith   PetscErrorCode ierr;
7595fce210SBarry Smith 
7695fce210SBarry Smith   PetscFunctionBegin;
7795fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
7895fce210SBarry Smith   sf->mine   = NULL;
7995fce210SBarry Smith   ierr       = PetscFree(sf->mine_alloc);CHKERRQ(ierr);
8095fce210SBarry Smith   sf->remote = NULL;
8195fce210SBarry Smith   ierr       = PetscFree(sf->remote_alloc);CHKERRQ(ierr);
8295fce210SBarry Smith   ierr       = PetscFree4(sf->ranks,sf->roffset,sf->rmine,sf->rremote);CHKERRQ(ierr);
8395fce210SBarry Smith   ierr       = PetscFree(sf->degree);CHKERRQ(ierr);
8495fce210SBarry Smith   if (sf->ingroup  != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->ingroup);CHKERRQ(ierr);}
8595fce210SBarry Smith   if (sf->outgroup != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->outgroup);CHKERRQ(ierr);}
8695fce210SBarry Smith   ierr         = PetscSFDestroy(&sf->multi);CHKERRQ(ierr);
8795fce210SBarry Smith   sf->graphset = PETSC_FALSE;
8895fce210SBarry Smith   if (sf->ops->Reset) {ierr = (*sf->ops->Reset)(sf);CHKERRQ(ierr);}
8995fce210SBarry Smith   sf->setupcalled = PETSC_FALSE;
9095fce210SBarry Smith   PetscFunctionReturn(0);
9195fce210SBarry Smith }
9295fce210SBarry Smith 
9395fce210SBarry Smith #undef __FUNCT__
9495fce210SBarry Smith #define __FUNCT__ "PetscSFSetType"
9595fce210SBarry Smith /*@C
9695fce210SBarry Smith    PetscSFSetType - set the PetscSF communication implementation
9795fce210SBarry Smith 
9895fce210SBarry Smith    Collective on PetscSF
9995fce210SBarry Smith 
10095fce210SBarry Smith    Input Parameters:
10195fce210SBarry Smith +  sf - the PetscSF context
10295fce210SBarry Smith -  type - a known method
10395fce210SBarry Smith 
10495fce210SBarry Smith    Options Database Key:
10595fce210SBarry Smith .  -sf_type <type> - Sets the method; use -help for a list
10695fce210SBarry Smith    of available methods (for instance, window, pt2pt, neighbor)
10795fce210SBarry Smith 
10895fce210SBarry Smith    Notes:
10995fce210SBarry Smith    See "include/petscsf.h" for available methods (for instance)
11095fce210SBarry Smith +    PETSCSFWINDOW - MPI-2/3 one-sided
11195fce210SBarry Smith -    PETSCSFBASIC - basic implementation using MPI-1 two-sided
11295fce210SBarry Smith 
11395fce210SBarry Smith   Level: intermediate
11495fce210SBarry Smith 
11595fce210SBarry Smith .keywords: PetscSF, set, type
11695fce210SBarry Smith 
11795fce210SBarry Smith .seealso: PetscSFType, PetscSFCreate()
11895fce210SBarry Smith @*/
11995fce210SBarry Smith PetscErrorCode PetscSFSetType(PetscSF sf,PetscSFType type)
12095fce210SBarry Smith {
12195fce210SBarry Smith   PetscErrorCode ierr,(*r)(PetscSF);
12295fce210SBarry Smith   PetscBool      match;
12395fce210SBarry Smith 
12495fce210SBarry Smith   PetscFunctionBegin;
12595fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
12695fce210SBarry Smith   PetscValidCharPointer(type,2);
12795fce210SBarry Smith 
12895fce210SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)sf,type,&match);CHKERRQ(ierr);
12995fce210SBarry Smith   if (match) PetscFunctionReturn(0);
13095fce210SBarry Smith 
131adc40e5bSBarry Smith   ierr = PetscFunctionListFind(PetscSFList,type,&r);CHKERRQ(ierr);
13295fce210SBarry Smith   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested PetscSF type %s",type);
13395fce210SBarry Smith   /* Destroy the previous private PetscSF context */
13495fce210SBarry Smith   if (sf->ops->Destroy) {
13595fce210SBarry Smith     ierr = (*(sf)->ops->Destroy)(sf);CHKERRQ(ierr);
13695fce210SBarry Smith   }
13795fce210SBarry Smith   ierr = PetscMemzero(sf->ops,sizeof(*sf->ops));CHKERRQ(ierr);
13895fce210SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)sf,type);CHKERRQ(ierr);
13995fce210SBarry Smith   ierr = (*r)(sf);CHKERRQ(ierr);
14095fce210SBarry Smith   PetscFunctionReturn(0);
14195fce210SBarry Smith }
14295fce210SBarry Smith 
14395fce210SBarry Smith #undef __FUNCT__
14495fce210SBarry Smith #define __FUNCT__ "PetscSFDestroy"
145d36d33e4SMatthew G. Knepley /*@
14695fce210SBarry Smith    PetscSFDestroy - destroy star forest
14795fce210SBarry Smith 
14895fce210SBarry Smith    Collective
14995fce210SBarry Smith 
15095fce210SBarry Smith    Input Arguments:
15195fce210SBarry Smith .  sf - address of star forest
15295fce210SBarry Smith 
15395fce210SBarry Smith    Level: intermediate
15495fce210SBarry Smith 
15595fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFReset()
15695fce210SBarry Smith @*/
15795fce210SBarry Smith PetscErrorCode PetscSFDestroy(PetscSF *sf)
15895fce210SBarry Smith {
15995fce210SBarry Smith   PetscErrorCode ierr;
16095fce210SBarry Smith 
16195fce210SBarry Smith   PetscFunctionBegin;
16295fce210SBarry Smith   if (!*sf) PetscFunctionReturn(0);
16395fce210SBarry Smith   PetscValidHeaderSpecific((*sf),PETSCSF_CLASSID,1);
16495fce210SBarry Smith   if (--((PetscObject)(*sf))->refct > 0) {*sf = 0; PetscFunctionReturn(0);}
16595fce210SBarry Smith   ierr = PetscSFReset(*sf);CHKERRQ(ierr);
16695fce210SBarry Smith   if ((*sf)->ops->Destroy) {ierr = (*(*sf)->ops->Destroy)(*sf);CHKERRQ(ierr);}
16795fce210SBarry Smith   ierr = PetscHeaderDestroy(sf);CHKERRQ(ierr);
16895fce210SBarry Smith   PetscFunctionReturn(0);
16995fce210SBarry Smith }
17095fce210SBarry Smith 
17195fce210SBarry Smith #undef __FUNCT__
17295fce210SBarry Smith #define __FUNCT__ "PetscSFSetUp"
17395fce210SBarry Smith /*@
17495fce210SBarry Smith    PetscSFSetUp - set up communication structures
17595fce210SBarry Smith 
17695fce210SBarry Smith    Collective
17795fce210SBarry Smith 
17895fce210SBarry Smith    Input Arguments:
17995fce210SBarry Smith .  sf - star forest communication object
18095fce210SBarry Smith 
18195fce210SBarry Smith    Level: beginner
18295fce210SBarry Smith 
18395fce210SBarry Smith .seealso: PetscSFSetFromOptions(), PetscSFSetType()
18495fce210SBarry Smith @*/
18595fce210SBarry Smith PetscErrorCode PetscSFSetUp(PetscSF sf)
18695fce210SBarry Smith {
18795fce210SBarry Smith   PetscErrorCode ierr;
18895fce210SBarry Smith 
18995fce210SBarry Smith   PetscFunctionBegin;
19095fce210SBarry Smith   if (sf->setupcalled) PetscFunctionReturn(0);
19195fce210SBarry Smith   if (!((PetscObject)sf)->type_name) {ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr);}
19295fce210SBarry Smith   if (sf->ops->SetUp) {ierr = (*sf->ops->SetUp)(sf);CHKERRQ(ierr);}
19395fce210SBarry Smith   sf->setupcalled = PETSC_TRUE;
19495fce210SBarry Smith   PetscFunctionReturn(0);
19595fce210SBarry Smith }
19695fce210SBarry Smith 
19795fce210SBarry Smith #undef __FUNCT__
19895fce210SBarry Smith #define __FUNCT__ "PetscSFSetFromOptions"
19995fce210SBarry Smith /*@C
20095fce210SBarry Smith    PetscSFSetFromOptions - set PetscSF options using the options database
20195fce210SBarry Smith 
20295fce210SBarry Smith    Logically Collective
20395fce210SBarry Smith 
20495fce210SBarry Smith    Input Arguments:
20595fce210SBarry Smith .  sf - star forest
20695fce210SBarry Smith 
20795fce210SBarry Smith    Options Database Keys:
20860263706SJed Brown +  -sf_type - implementation type, see PetscSFSetType()
20960263706SJed Brown -  -sf_rank_order - sort composite points for gathers and scatters in rank order, gathers are non-deterministic otherwise
21095fce210SBarry Smith 
21195fce210SBarry Smith    Level: intermediate
21295fce210SBarry Smith 
21395fce210SBarry Smith .keywords: KSP, set, from, options, database
21495fce210SBarry Smith 
21595fce210SBarry Smith .seealso: PetscSFWindowSetSyncType()
21695fce210SBarry Smith @*/
21795fce210SBarry Smith PetscErrorCode PetscSFSetFromOptions(PetscSF sf)
21895fce210SBarry Smith {
21995fce210SBarry Smith   PetscSFType    deft;
22095fce210SBarry Smith   char           type[256];
22195fce210SBarry Smith   PetscErrorCode ierr;
22295fce210SBarry Smith   PetscBool      flg;
22395fce210SBarry Smith 
22495fce210SBarry Smith   PetscFunctionBegin;
22595fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
22695fce210SBarry Smith   ierr = PetscObjectOptionsBegin((PetscObject)sf);CHKERRQ(ierr);
22795fce210SBarry Smith   deft = ((PetscObject)sf)->type_name ? ((PetscObject)sf)->type_name : PETSCSFBASIC;
228a264d7a6SBarry Smith   ierr = PetscOptionsFList("-sf_type","PetscSF implementation type","PetscSFSetType",PetscSFList,deft,type,256,&flg);CHKERRQ(ierr);
22995fce210SBarry Smith   ierr = PetscSFSetType(sf,flg ? type : deft);CHKERRQ(ierr);
23095fce210SBarry Smith   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);
231e55864a3SBarry Smith   if (sf->ops->SetFromOptions) {ierr = (*sf->ops->SetFromOptions)(PetscOptionsObject,sf);CHKERRQ(ierr);}
23295fce210SBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
23395fce210SBarry Smith   PetscFunctionReturn(0);
23495fce210SBarry Smith }
23595fce210SBarry Smith 
23695fce210SBarry Smith #undef __FUNCT__
23795fce210SBarry Smith #define __FUNCT__ "PetscSFSetRankOrder"
23895fce210SBarry Smith /*@C
23995fce210SBarry Smith    PetscSFSetRankOrder - sort multi-points for gathers and scatters by rank order
24095fce210SBarry Smith 
24195fce210SBarry Smith    Logically Collective
24295fce210SBarry Smith 
24395fce210SBarry Smith    Input Arguments:
24495fce210SBarry Smith +  sf - star forest
24595fce210SBarry Smith -  flg - PETSC_TRUE to sort, PETSC_FALSE to skip sorting (lower setup cost, but non-deterministic)
24695fce210SBarry Smith 
24795fce210SBarry Smith    Level: advanced
24895fce210SBarry Smith 
24995fce210SBarry Smith .seealso: PetscSFGatherBegin(), PetscSFScatterBegin()
25095fce210SBarry Smith @*/
25195fce210SBarry Smith PetscErrorCode PetscSFSetRankOrder(PetscSF sf,PetscBool flg)
25295fce210SBarry Smith {
25395fce210SBarry Smith 
25495fce210SBarry Smith   PetscFunctionBegin;
25595fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
25695fce210SBarry Smith   PetscValidLogicalCollectiveBool(sf,flg,2);
25795fce210SBarry Smith   if (sf->multi) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Rank ordering must be set before first call to PetscSFGatherBegin() or PetscSFScatterBegin()");
25895fce210SBarry Smith   sf->rankorder = flg;
25995fce210SBarry Smith   PetscFunctionReturn(0);
26095fce210SBarry Smith }
26195fce210SBarry Smith 
26295fce210SBarry Smith #undef __FUNCT__
26395fce210SBarry Smith #define __FUNCT__ "PetscSFSetGraph"
26495fce210SBarry Smith /*@C
26595fce210SBarry Smith    PetscSFSetGraph - Set a parallel star forest
26695fce210SBarry Smith 
26795fce210SBarry Smith    Collective
26895fce210SBarry Smith 
26995fce210SBarry Smith    Input Arguments:
27095fce210SBarry Smith +  sf - star forest
27195fce210SBarry Smith .  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
27295fce210SBarry Smith .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
27395fce210SBarry Smith .  ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage
27495fce210SBarry Smith .  localmode - copy mode for ilocal
27595fce210SBarry Smith .  iremote - remote locations of root vertices for each leaf on the current process
27695fce210SBarry Smith -  remotemode - copy mode for iremote
27795fce210SBarry Smith 
27895fce210SBarry Smith    Level: intermediate
27995fce210SBarry Smith 
28095fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFGetGraph()
28195fce210SBarry Smith @*/
28295fce210SBarry Smith PetscErrorCode PetscSFSetGraph(PetscSF sf,PetscInt nroots,PetscInt nleaves,const PetscInt *ilocal,PetscCopyMode localmode,const PetscSFNode *iremote,PetscCopyMode remotemode)
28395fce210SBarry Smith {
28495fce210SBarry Smith   PetscErrorCode     ierr;
28595fce210SBarry Smith   PetscTable         table;
28695fce210SBarry Smith   PetscTablePosition pos;
28795fce210SBarry Smith   PetscMPIInt        size;
28895fce210SBarry Smith   PetscInt           i,*rcount,*ranks;
28995fce210SBarry Smith 
29095fce210SBarry Smith   PetscFunctionBegin;
29195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
292563ffbabSMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr);
29395fce210SBarry Smith   if (nleaves && ilocal) PetscValidIntPointer(ilocal,4);
29495fce210SBarry Smith   if (nleaves) PetscValidPointer(iremote,6);
29595fce210SBarry Smith   if (nroots < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"roots %D, cannot be negative",nroots);
29695fce210SBarry Smith   if (nleaves < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nleaves %D, cannot be negative",nleaves);
29795fce210SBarry Smith   ierr        = PetscSFReset(sf);CHKERRQ(ierr);
29895fce210SBarry Smith   sf->nroots  = nroots;
29995fce210SBarry Smith   sf->nleaves = nleaves;
30095fce210SBarry Smith   if (ilocal) {
30195fce210SBarry Smith     switch (localmode) {
30295fce210SBarry Smith     case PETSC_COPY_VALUES:
303785e854fSJed Brown       ierr        = PetscMalloc1(nleaves,&sf->mine_alloc);CHKERRQ(ierr);
30495fce210SBarry Smith       sf->mine    = sf->mine_alloc;
30595fce210SBarry Smith       ierr        = PetscMemcpy(sf->mine,ilocal,nleaves*sizeof(*sf->mine));CHKERRQ(ierr);
30695fce210SBarry Smith       sf->minleaf = PETSC_MAX_INT;
30795fce210SBarry Smith       sf->maxleaf = PETSC_MIN_INT;
30895fce210SBarry Smith       for (i=0; i<nleaves; i++) {
30995fce210SBarry Smith         sf->minleaf = PetscMin(sf->minleaf,ilocal[i]);
31095fce210SBarry Smith         sf->maxleaf = PetscMax(sf->maxleaf,ilocal[i]);
31195fce210SBarry Smith       }
31295fce210SBarry Smith       break;
31395fce210SBarry Smith     case PETSC_OWN_POINTER:
31495fce210SBarry Smith       sf->mine_alloc = (PetscInt*)ilocal;
31595fce210SBarry Smith       sf->mine       = sf->mine_alloc;
31695fce210SBarry Smith       break;
31795fce210SBarry Smith     case PETSC_USE_POINTER:
31895fce210SBarry Smith       sf->mine = (PetscInt*)ilocal;
31995fce210SBarry Smith       break;
32095fce210SBarry Smith     default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown localmode");
32195fce210SBarry Smith     }
32295fce210SBarry Smith   }
32395fce210SBarry Smith   if (!ilocal || nleaves > 0) {
32495fce210SBarry Smith     sf->minleaf = 0;
32595fce210SBarry Smith     sf->maxleaf = nleaves - 1;
32695fce210SBarry Smith   }
32795fce210SBarry Smith   switch (remotemode) {
32895fce210SBarry Smith   case PETSC_COPY_VALUES:
329785e854fSJed Brown     ierr       = PetscMalloc1(nleaves,&sf->remote_alloc);CHKERRQ(ierr);
33095fce210SBarry Smith     sf->remote = sf->remote_alloc;
33195fce210SBarry Smith     ierr       = PetscMemcpy(sf->remote,iremote,nleaves*sizeof(*sf->remote));CHKERRQ(ierr);
33295fce210SBarry Smith     break;
33395fce210SBarry Smith   case PETSC_OWN_POINTER:
33495fce210SBarry Smith     sf->remote_alloc = (PetscSFNode*)iremote;
33595fce210SBarry Smith     sf->remote       = sf->remote_alloc;
33695fce210SBarry Smith     break;
33795fce210SBarry Smith   case PETSC_USE_POINTER:
33895fce210SBarry Smith     sf->remote = (PetscSFNode*)iremote;
33995fce210SBarry Smith     break;
34095fce210SBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown remotemode");
34195fce210SBarry Smith   }
34295fce210SBarry Smith 
34395fce210SBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr);
34495fce210SBarry Smith   ierr = PetscTableCreate(10,size,&table);CHKERRQ(ierr);
34595fce210SBarry Smith   for (i=0; i<nleaves; i++) {
34695fce210SBarry Smith     /* Log 1-based rank */
34795fce210SBarry Smith     ierr = PetscTableAdd(table,iremote[i].rank+1,1,ADD_VALUES);CHKERRQ(ierr);
34895fce210SBarry Smith   }
34995fce210SBarry Smith   ierr = PetscTableGetCount(table,&sf->nranks);CHKERRQ(ierr);
350dcca6d9dSJed Brown   ierr = PetscMalloc4(sf->nranks,&sf->ranks,sf->nranks+1,&sf->roffset,nleaves,&sf->rmine,nleaves,&sf->rremote);CHKERRQ(ierr);
351dcca6d9dSJed Brown   ierr = PetscMalloc2(sf->nranks,&rcount,sf->nranks,&ranks);CHKERRQ(ierr);
35295fce210SBarry Smith   ierr = PetscTableGetHeadPosition(table,&pos);CHKERRQ(ierr);
35395fce210SBarry Smith   for (i=0; i<sf->nranks; i++) {
35495fce210SBarry Smith     ierr = PetscTableGetNext(table,&pos,&ranks[i],&rcount[i]);CHKERRQ(ierr);
35595fce210SBarry Smith     ranks[i]--;             /* Convert back to 0-based */
35695fce210SBarry Smith   }
35795fce210SBarry Smith   ierr = PetscTableDestroy(&table);CHKERRQ(ierr);
35895fce210SBarry Smith   ierr = PetscSortIntWithArray(sf->nranks,ranks,rcount);CHKERRQ(ierr);
35995fce210SBarry Smith   sf->roffset[0] = 0;
36095fce210SBarry Smith   for (i=0; i<sf->nranks; i++) {
36195fce210SBarry Smith     ierr = PetscMPIIntCast(ranks[i],sf->ranks+i);CHKERRQ(ierr);
36295fce210SBarry Smith     sf->roffset[i+1] = sf->roffset[i] + rcount[i];
36395fce210SBarry Smith     rcount[i]        = 0;
36495fce210SBarry Smith   }
36595fce210SBarry Smith   for (i=0; i<nleaves; i++) {
36695fce210SBarry Smith     PetscInt lo,hi,irank;
36795fce210SBarry Smith     /* Search for index of iremote[i].rank in sf->ranks */
36895fce210SBarry Smith     lo = 0; hi = sf->nranks;
36995fce210SBarry Smith     while (hi - lo > 1) {
37095fce210SBarry Smith       PetscInt mid = lo + (hi - lo)/2;
37195fce210SBarry Smith       if (iremote[i].rank < sf->ranks[mid]) hi = mid;
37295fce210SBarry Smith       else                                  lo = mid;
37395fce210SBarry Smith     }
37495fce210SBarry Smith     if (hi - lo == 1 && iremote[i].rank == sf->ranks[lo]) irank = lo;
37595fce210SBarry Smith     else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Could not find rank %D in array",iremote[i].rank);
37695fce210SBarry Smith     sf->rmine[sf->roffset[irank] + rcount[irank]]   = ilocal ? ilocal[i] : i;
37795fce210SBarry Smith     sf->rremote[sf->roffset[irank] + rcount[irank]] = iremote[i].index;
37895fce210SBarry Smith     rcount[irank]++;
37995fce210SBarry Smith   }
38095fce210SBarry Smith   ierr = PetscFree2(rcount,ranks);CHKERRQ(ierr);
38195fce210SBarry Smith #if !defined(PETSC_USE_64BIT_INDICES)
38295fce210SBarry Smith   if (nroots == PETSC_DETERMINE) {
38395fce210SBarry Smith     /* Jed, if you have a better way to do this, put it in */
38495fce210SBarry Smith     PetscInt *numRankLeaves, *leafOff, *leafIndices, *numRankRoots, *rootOff, *rootIndices, maxRoots = 0;
38595fce210SBarry Smith 
38695fce210SBarry Smith     /* All to all to determine number of leaf indices from each (you can do this using Scan and asynch messages) */
387dcca6d9dSJed Brown     ierr = PetscMalloc4(size,&numRankLeaves,size+1,&leafOff,size,&numRankRoots,size+1,&rootOff);CHKERRQ(ierr);
38895fce210SBarry Smith     ierr = PetscMemzero(numRankLeaves, size * sizeof(PetscInt));CHKERRQ(ierr);
38995fce210SBarry Smith     for (i = 0; i < nleaves; ++i) ++numRankLeaves[iremote[i].rank];
39095fce210SBarry Smith     ierr = MPI_Alltoall(numRankLeaves, 1, MPIU_INT, numRankRoots, 1, MPIU_INT, PetscObjectComm((PetscObject)sf));CHKERRQ(ierr);
39195fce210SBarry Smith     /* Could set nroots to this maximum */
39295fce210SBarry Smith     for (i = 0; i < size; ++i) maxRoots += numRankRoots[i];
39395fce210SBarry Smith 
39495fce210SBarry Smith     /* Gather all indices */
395dcca6d9dSJed Brown     ierr = PetscMalloc2(nleaves,&leafIndices,maxRoots,&rootIndices);CHKERRQ(ierr);
39695fce210SBarry Smith     leafOff[0] = 0;
39795fce210SBarry Smith     for (i = 0; i < size; ++i) leafOff[i+1] = leafOff[i] + numRankLeaves[i];
39895fce210SBarry Smith     for (i = 0; i < nleaves; ++i) leafIndices[leafOff[iremote[i].rank]++] = iremote[i].index;
39995fce210SBarry Smith     leafOff[0] = 0;
40095fce210SBarry Smith     for (i = 0; i < size; ++i) leafOff[i+1] = leafOff[i] + numRankLeaves[i];
40195fce210SBarry Smith     rootOff[0] = 0;
40295fce210SBarry Smith     for (i = 0; i < size; ++i) rootOff[i+1] = rootOff[i] + numRankRoots[i];
40395fce210SBarry Smith     ierr = MPI_Alltoallv(leafIndices, numRankLeaves, leafOff, MPIU_INT, rootIndices, numRankRoots, rootOff, MPIU_INT, PetscObjectComm((PetscObject)sf));CHKERRQ(ierr);
40495fce210SBarry Smith     /* Sort and reduce */
40595fce210SBarry Smith     ierr       = PetscSortRemoveDupsInt(&maxRoots, rootIndices);CHKERRQ(ierr);
40695fce210SBarry Smith     ierr       = PetscFree2(leafIndices,rootIndices);CHKERRQ(ierr);
40795fce210SBarry Smith     ierr       = PetscFree4(numRankLeaves,leafOff,numRankRoots,rootOff);CHKERRQ(ierr);
40895fce210SBarry Smith     sf->nroots = maxRoots;
40995fce210SBarry Smith   }
41095fce210SBarry Smith #endif
41195fce210SBarry Smith 
41295fce210SBarry Smith   sf->graphset = PETSC_TRUE;
413563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr);
41495fce210SBarry Smith   PetscFunctionReturn(0);
41595fce210SBarry Smith }
41695fce210SBarry Smith 
41795fce210SBarry Smith #undef __FUNCT__
41895fce210SBarry Smith #define __FUNCT__ "PetscSFCreateInverseSF"
41995fce210SBarry Smith /*@C
42095fce210SBarry Smith    PetscSFCreateInverseSF - given a PetscSF in which all vertices have degree 1, creates the inverse map
42195fce210SBarry Smith 
42295fce210SBarry Smith    Collective
42395fce210SBarry Smith 
42495fce210SBarry Smith    Input Arguments:
42595fce210SBarry Smith .  sf - star forest to invert
42695fce210SBarry Smith 
42795fce210SBarry Smith    Output Arguments:
42895fce210SBarry Smith .  isf - inverse of sf
42995fce210SBarry Smith 
43095fce210SBarry Smith    Level: advanced
43195fce210SBarry Smith 
43295fce210SBarry Smith    Notes:
43395fce210SBarry Smith    All roots must have degree 1.
43495fce210SBarry Smith 
43595fce210SBarry Smith    The local space may be a permutation, but cannot be sparse.
43695fce210SBarry Smith 
43795fce210SBarry Smith .seealso: PetscSFSetGraph()
43895fce210SBarry Smith @*/
43995fce210SBarry Smith PetscErrorCode PetscSFCreateInverseSF(PetscSF sf,PetscSF *isf)
44095fce210SBarry Smith {
44195fce210SBarry Smith   PetscErrorCode ierr;
44295fce210SBarry Smith   PetscMPIInt    rank;
44395fce210SBarry Smith   PetscInt       i,nroots,nleaves,maxlocal,count,*newilocal;
44495fce210SBarry Smith   const PetscInt *ilocal;
44595fce210SBarry Smith   PetscSFNode    *roots,*leaves;
44695fce210SBarry Smith 
44795fce210SBarry Smith   PetscFunctionBegin;
44895fce210SBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
44995fce210SBarry Smith   ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);CHKERRQ(ierr);
45095fce210SBarry Smith   for (i=0,maxlocal=0; i<nleaves; i++) maxlocal = PetscMax(maxlocal,(ilocal ? ilocal[i] : i)+1);
451dcca6d9dSJed Brown   ierr = PetscMalloc2(nroots,&roots,nleaves,&leaves);CHKERRQ(ierr);
45295fce210SBarry Smith   for (i=0; i<nleaves; i++) {
45395fce210SBarry Smith     leaves[i].rank  = rank;
45495fce210SBarry Smith     leaves[i].index = i;
45595fce210SBarry Smith   }
45695fce210SBarry Smith   for (i=0; i <nroots; i++) {
45795fce210SBarry Smith     roots[i].rank  = -1;
45895fce210SBarry Smith     roots[i].index = -1;
45995fce210SBarry Smith   }
4608bfbc91cSJed Brown   ierr = PetscSFReduceBegin(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr);
4618bfbc91cSJed Brown   ierr = PetscSFReduceEnd(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr);
46295fce210SBarry Smith 
46395fce210SBarry Smith   /* Check whether our leaves are sparse */
46495fce210SBarry Smith   for (i=0,count=0; i<nroots; i++) if (roots[i].rank >= 0) count++;
46595fce210SBarry Smith   if (count == nroots) newilocal = NULL;
46695fce210SBarry Smith   else {                        /* Index for sparse leaves and compact "roots" array (which is to become our leaves). */
467785e854fSJed Brown     ierr = PetscMalloc1(count,&newilocal);CHKERRQ(ierr);
46895fce210SBarry Smith     for (i=0,count=0; i<nroots; i++) {
46995fce210SBarry Smith       if (roots[i].rank >= 0) {
47095fce210SBarry Smith         newilocal[count]   = i;
47195fce210SBarry Smith         roots[count].rank  = roots[i].rank;
47295fce210SBarry Smith         roots[count].index = roots[i].index;
47395fce210SBarry Smith         count++;
47495fce210SBarry Smith       }
47595fce210SBarry Smith     }
47695fce210SBarry Smith   }
47795fce210SBarry Smith 
47895fce210SBarry Smith   ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,isf);CHKERRQ(ierr);
47995fce210SBarry Smith   ierr = PetscSFSetGraph(*isf,maxlocal,count,newilocal,PETSC_OWN_POINTER,roots,PETSC_COPY_VALUES);CHKERRQ(ierr);
48095fce210SBarry Smith   ierr = PetscFree2(roots,leaves);CHKERRQ(ierr);
48195fce210SBarry Smith   PetscFunctionReturn(0);
48295fce210SBarry Smith }
48395fce210SBarry Smith 
48495fce210SBarry Smith #undef __FUNCT__
48595fce210SBarry Smith #define __FUNCT__ "PetscSFDuplicate"
48695fce210SBarry Smith /*@
48795fce210SBarry Smith    PetscSFDuplicate - duplicate a PetscSF, optionally preserving rank connectivity and graph
48895fce210SBarry Smith 
48995fce210SBarry Smith    Collective
49095fce210SBarry Smith 
49195fce210SBarry Smith    Input Arguments:
49295fce210SBarry Smith +  sf - communication object to duplicate
49395fce210SBarry Smith -  opt - PETSCSF_DUPLICATE_CONFONLY, PETSCSF_DUPLICATE_RANKS, or PETSCSF_DUPLICATE_GRAPH (see PetscSFDuplicateOption)
49495fce210SBarry Smith 
49595fce210SBarry Smith    Output Arguments:
49695fce210SBarry Smith .  newsf - new communication object
49795fce210SBarry Smith 
49895fce210SBarry Smith    Level: beginner
49995fce210SBarry Smith 
50095fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetType(), PetscSFSetGraph()
50195fce210SBarry Smith @*/
50295fce210SBarry Smith PetscErrorCode PetscSFDuplicate(PetscSF sf,PetscSFDuplicateOption opt,PetscSF *newsf)
50395fce210SBarry Smith {
50495fce210SBarry Smith   PetscErrorCode ierr;
50595fce210SBarry Smith 
50695fce210SBarry Smith   PetscFunctionBegin;
50795fce210SBarry Smith   ierr = PetscSFCreate(PetscObjectComm((PetscObject)sf),newsf);CHKERRQ(ierr);
50895fce210SBarry Smith   ierr = PetscSFSetType(*newsf,((PetscObject)sf)->type_name);CHKERRQ(ierr);
50995fce210SBarry Smith   if (sf->ops->Duplicate) {ierr = (*sf->ops->Duplicate)(sf,opt,*newsf);CHKERRQ(ierr);}
51095fce210SBarry Smith   if (opt == PETSCSF_DUPLICATE_GRAPH) {
51195fce210SBarry Smith     PetscInt          nroots,nleaves;
51295fce210SBarry Smith     const PetscInt    *ilocal;
51395fce210SBarry Smith     const PetscSFNode *iremote;
51495fce210SBarry Smith     ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
51595fce210SBarry Smith     ierr = PetscSFSetGraph(*newsf,nroots,nleaves,ilocal,PETSC_COPY_VALUES,iremote,PETSC_COPY_VALUES);CHKERRQ(ierr);
51695fce210SBarry Smith   }
51795fce210SBarry Smith   PetscFunctionReturn(0);
51895fce210SBarry Smith }
51995fce210SBarry Smith 
52095fce210SBarry Smith #undef __FUNCT__
52195fce210SBarry Smith #define __FUNCT__ "PetscSFGetGraph"
52295fce210SBarry Smith /*@C
52395fce210SBarry Smith    PetscSFGetGraph - Get the graph specifying a parallel star forest
52495fce210SBarry Smith 
52595fce210SBarry Smith    Not Collective
52695fce210SBarry Smith 
52795fce210SBarry Smith    Input Arguments:
52895fce210SBarry Smith .  sf - star forest
52995fce210SBarry Smith 
53095fce210SBarry Smith    Output Arguments:
53195fce210SBarry Smith +  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
53295fce210SBarry Smith .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
53395fce210SBarry Smith .  ilocal - locations of leaves in leafdata buffers
53495fce210SBarry Smith -  iremote - remote locations of root vertices for each leaf on the current process
53595fce210SBarry Smith 
53695fce210SBarry Smith    Level: intermediate
53795fce210SBarry Smith 
53895fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph()
53995fce210SBarry Smith @*/
54095fce210SBarry Smith PetscErrorCode PetscSFGetGraph(PetscSF sf,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote)
54195fce210SBarry Smith {
54295fce210SBarry Smith 
54395fce210SBarry Smith   PetscFunctionBegin;
54495fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
54595fce210SBarry Smith   /* We are not currently requiring that the graph is set, thus returning nroots=-1 if it has not been set */
54695fce210SBarry Smith   /* if (!sf->graphset) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Graph has not been set, must call PetscSFSetGraph()"); */
54795fce210SBarry Smith   if (nroots) *nroots = sf->nroots;
54895fce210SBarry Smith   if (nleaves) *nleaves = sf->nleaves;
54995fce210SBarry Smith   if (ilocal) *ilocal = sf->mine;
55095fce210SBarry Smith   if (iremote) *iremote = sf->remote;
55195fce210SBarry Smith   PetscFunctionReturn(0);
55295fce210SBarry Smith }
55395fce210SBarry Smith 
55495fce210SBarry Smith #undef __FUNCT__
55595fce210SBarry Smith #define __FUNCT__ "PetscSFGetLeafRange"
55695fce210SBarry Smith /*@C
55795fce210SBarry Smith    PetscSFGetLeafRange - Get the active leaf ranges
55895fce210SBarry Smith 
55995fce210SBarry Smith    Not Collective
56095fce210SBarry Smith 
56195fce210SBarry Smith    Input Arguments:
56295fce210SBarry Smith .  sf - star forest
56395fce210SBarry Smith 
56495fce210SBarry Smith    Output Arguments:
56595fce210SBarry Smith +  minleaf - minimum active leaf on this process
56695fce210SBarry Smith -  maxleaf - maximum active leaf on this process
56795fce210SBarry Smith 
56895fce210SBarry Smith    Level: developer
56995fce210SBarry Smith 
57095fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph(), PetscSFGetGraph()
57195fce210SBarry Smith @*/
57295fce210SBarry Smith PetscErrorCode PetscSFGetLeafRange(PetscSF sf,PetscInt *minleaf,PetscInt *maxleaf)
57395fce210SBarry Smith {
57495fce210SBarry Smith 
57595fce210SBarry Smith   PetscFunctionBegin;
57695fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
57795fce210SBarry Smith   if (minleaf) *minleaf = sf->minleaf;
57895fce210SBarry Smith   if (maxleaf) *maxleaf = sf->maxleaf;
57995fce210SBarry Smith   PetscFunctionReturn(0);
58095fce210SBarry Smith }
58195fce210SBarry Smith 
58295fce210SBarry Smith #undef __FUNCT__
58395fce210SBarry Smith #define __FUNCT__ "PetscSFView"
58495fce210SBarry Smith /*@C
58595fce210SBarry Smith    PetscSFView - view a star forest
58695fce210SBarry Smith 
58795fce210SBarry Smith    Collective
58895fce210SBarry Smith 
58995fce210SBarry Smith    Input Arguments:
59095fce210SBarry Smith +  sf - star forest
59195fce210SBarry Smith -  viewer - viewer to display graph, for example PETSC_VIEWER_STDOUT_WORLD
59295fce210SBarry Smith 
59395fce210SBarry Smith    Level: beginner
59495fce210SBarry Smith 
59595fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph()
59695fce210SBarry Smith @*/
59795fce210SBarry Smith PetscErrorCode PetscSFView(PetscSF sf,PetscViewer viewer)
59895fce210SBarry Smith {
59995fce210SBarry Smith   PetscErrorCode    ierr;
60095fce210SBarry Smith   PetscBool         iascii;
60195fce210SBarry Smith   PetscViewerFormat format;
60295fce210SBarry Smith 
60395fce210SBarry Smith   PetscFunctionBegin;
60495fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
60595fce210SBarry Smith   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sf),&viewer);CHKERRQ(ierr);}
60695fce210SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
60795fce210SBarry Smith   PetscCheckSameComm(sf,1,viewer,2);
60895fce210SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
60995fce210SBarry Smith   if (iascii) {
61095fce210SBarry Smith     PetscMPIInt rank;
61195fce210SBarry Smith     PetscInt    i,j;
61295fce210SBarry Smith 
613dae58748SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)sf,viewer);CHKERRQ(ierr);
61495fce210SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
61595fce210SBarry Smith     if (sf->ops->View) {ierr = (*sf->ops->View)(sf,viewer);CHKERRQ(ierr);}
61695fce210SBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
61795fce210SBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
61895fce210SBarry Smith     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of roots=%D, leaves=%D, remote ranks=%D\n",rank,sf->nroots,sf->nleaves,sf->nranks);CHKERRQ(ierr);
61995fce210SBarry Smith     for (i=0; i<sf->nleaves; i++) {
62095fce210SBarry Smith       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);
62195fce210SBarry Smith     }
62295fce210SBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
62395fce210SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
62495fce210SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
62595fce210SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Roots referenced by my leaves, by rank\n",rank);CHKERRQ(ierr);
62695fce210SBarry Smith       for (i=0; i<sf->nranks; i++) {
6277904a332SBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %d: %D edges\n",rank,sf->ranks[i],sf->roffset[i+1]-sf->roffset[i]);CHKERRQ(ierr);
62895fce210SBarry Smith         for (j=sf->roffset[i]; j<sf->roffset[i+1]; j++) {
62995fce210SBarry Smith           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d]    %D <- %D\n",rank,sf->rmine[j],sf->rremote[j]);CHKERRQ(ierr);
63095fce210SBarry Smith         }
63195fce210SBarry Smith       }
63295fce210SBarry Smith     }
63395fce210SBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
63495fce210SBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
63595fce210SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
63695fce210SBarry Smith   }
63795fce210SBarry Smith   PetscFunctionReturn(0);
63895fce210SBarry Smith }
63995fce210SBarry Smith 
64095fce210SBarry Smith #undef __FUNCT__
64195fce210SBarry Smith #define __FUNCT__ "PetscSFGetRanks"
64295fce210SBarry Smith /*@C
64395fce210SBarry Smith    PetscSFGetRanks - Get ranks and number of vertices referenced by leaves on this process
64495fce210SBarry Smith 
64595fce210SBarry Smith    Not Collective
64695fce210SBarry Smith 
64795fce210SBarry Smith    Input Arguments:
64895fce210SBarry Smith .  sf - star forest
64995fce210SBarry Smith 
65095fce210SBarry Smith    Output Arguments:
65195fce210SBarry Smith +  nranks - number of ranks referenced by local part
65295fce210SBarry Smith .  ranks - array of ranks
65395fce210SBarry Smith .  roffset - offset in rmine/rremote for each rank (length nranks+1)
65495fce210SBarry Smith .  rmine - concatenated array holding local indices referencing each remote rank
65595fce210SBarry Smith -  rremote - concatenated array holding remote indices referenced for each remote rank
65695fce210SBarry Smith 
65795fce210SBarry Smith    Level: developer
65895fce210SBarry Smith 
65995fce210SBarry Smith .seealso: PetscSFSetGraph()
66095fce210SBarry Smith @*/
66195fce210SBarry Smith PetscErrorCode PetscSFGetRanks(PetscSF sf,PetscInt *nranks,const PetscMPIInt **ranks,const PetscInt **roffset,const PetscInt **rmine,const PetscInt **rremote)
66295fce210SBarry Smith {
66395fce210SBarry Smith 
66495fce210SBarry Smith   PetscFunctionBegin;
66595fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
66695fce210SBarry Smith   if (nranks)  *nranks  = sf->nranks;
66795fce210SBarry Smith   if (ranks)   *ranks   = sf->ranks;
66895fce210SBarry Smith   if (roffset) *roffset = sf->roffset;
66995fce210SBarry Smith   if (rmine)   *rmine   = sf->rmine;
67095fce210SBarry Smith   if (rremote) *rremote = sf->rremote;
67195fce210SBarry Smith   PetscFunctionReturn(0);
67295fce210SBarry Smith }
67395fce210SBarry Smith 
67495fce210SBarry Smith #undef __FUNCT__
67595fce210SBarry Smith #define __FUNCT__ "PetscSFGetGroups"
67695fce210SBarry Smith /*@C
67795fce210SBarry Smith    PetscSFGetGroups - gets incoming and outgoing process groups
67895fce210SBarry Smith 
67995fce210SBarry Smith    Collective
68095fce210SBarry Smith 
68195fce210SBarry Smith    Input Argument:
68295fce210SBarry Smith .  sf - star forest
68395fce210SBarry Smith 
68495fce210SBarry Smith    Output Arguments:
68595fce210SBarry Smith +  incoming - group of origin processes for incoming edges (leaves that reference my roots)
68695fce210SBarry Smith -  outgoing - group of destination processes for outgoing edges (roots that I reference)
68795fce210SBarry Smith 
68895fce210SBarry Smith    Level: developer
68995fce210SBarry Smith 
69095fce210SBarry Smith .seealso: PetscSFGetWindow(), PetscSFRestoreWindow()
69195fce210SBarry Smith @*/
69295fce210SBarry Smith PetscErrorCode PetscSFGetGroups(PetscSF sf,MPI_Group *incoming,MPI_Group *outgoing)
69395fce210SBarry Smith {
69495fce210SBarry Smith   PetscErrorCode ierr;
69595fce210SBarry Smith   MPI_Group      group;
69695fce210SBarry Smith 
69795fce210SBarry Smith   PetscFunctionBegin;
69895fce210SBarry Smith   if (sf->ingroup == MPI_GROUP_NULL) {
69995fce210SBarry Smith     PetscInt       i;
70095fce210SBarry Smith     const PetscInt *indegree;
70195fce210SBarry Smith     PetscMPIInt    rank,*outranks,*inranks;
70295fce210SBarry Smith     PetscSFNode    *remote;
70395fce210SBarry Smith     PetscSF        bgcount;
70495fce210SBarry Smith 
70595fce210SBarry Smith     /* Compute the number of incoming ranks */
706785e854fSJed Brown     ierr = PetscMalloc1(sf->nranks,&remote);CHKERRQ(ierr);
70795fce210SBarry Smith     for (i=0; i<sf->nranks; i++) {
70895fce210SBarry Smith       remote[i].rank  = sf->ranks[i];
70995fce210SBarry Smith       remote[i].index = 0;
71095fce210SBarry Smith     }
71195fce210SBarry Smith     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,&bgcount);CHKERRQ(ierr);
71295fce210SBarry Smith     ierr = PetscSFSetGraph(bgcount,1,sf->nranks,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
71395fce210SBarry Smith     ierr = PetscSFComputeDegreeBegin(bgcount,&indegree);CHKERRQ(ierr);
71495fce210SBarry Smith     ierr = PetscSFComputeDegreeEnd(bgcount,&indegree);CHKERRQ(ierr);
71595fce210SBarry Smith 
71695fce210SBarry Smith     /* Enumerate the incoming ranks */
717dcca6d9dSJed Brown     ierr = PetscMalloc2(indegree[0],&inranks,sf->nranks,&outranks);CHKERRQ(ierr);
71895fce210SBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
71995fce210SBarry Smith     for (i=0; i<sf->nranks; i++) outranks[i] = rank;
72095fce210SBarry Smith     ierr = PetscSFGatherBegin(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr);
72195fce210SBarry Smith     ierr = PetscSFGatherEnd(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr);
72295fce210SBarry Smith     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
72395fce210SBarry Smith     ierr = MPI_Group_incl(group,indegree[0],inranks,&sf->ingroup);CHKERRQ(ierr);
72495fce210SBarry Smith     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
72595fce210SBarry Smith     ierr = PetscFree2(inranks,outranks);CHKERRQ(ierr);
72695fce210SBarry Smith     ierr = PetscSFDestroy(&bgcount);CHKERRQ(ierr);
72795fce210SBarry Smith   }
72895fce210SBarry Smith   *incoming = sf->ingroup;
72995fce210SBarry Smith 
73095fce210SBarry Smith   if (sf->outgroup == MPI_GROUP_NULL) {
73195fce210SBarry Smith     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
73295fce210SBarry Smith     ierr = MPI_Group_incl(group,sf->nranks,sf->ranks,&sf->outgroup);CHKERRQ(ierr);
73395fce210SBarry Smith     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
73495fce210SBarry Smith   }
73595fce210SBarry Smith   *outgoing = sf->outgroup;
73695fce210SBarry Smith   PetscFunctionReturn(0);
73795fce210SBarry Smith }
73895fce210SBarry Smith 
73995fce210SBarry Smith #undef __FUNCT__
74095fce210SBarry Smith #define __FUNCT__ "PetscSFGetMultiSF"
74195fce210SBarry Smith /*@C
74295fce210SBarry Smith    PetscSFGetMultiSF - gets the inner SF implemeting gathers and scatters
74395fce210SBarry Smith 
74495fce210SBarry Smith    Collective
74595fce210SBarry Smith 
74695fce210SBarry Smith    Input Argument:
74795fce210SBarry Smith .  sf - star forest that may contain roots with 0 or with more than 1 vertex
74895fce210SBarry Smith 
74995fce210SBarry Smith    Output Arguments:
75095fce210SBarry Smith .  multi - star forest with split roots, such that each root has degree exactly 1
75195fce210SBarry Smith 
75295fce210SBarry Smith    Level: developer
75395fce210SBarry Smith 
75495fce210SBarry Smith    Notes:
75595fce210SBarry Smith 
75695fce210SBarry Smith    In most cases, users should use PetscSFGatherBegin() and PetscSFScatterBegin() instead of manipulating multi
75795fce210SBarry Smith    directly. Since multi satisfies the stronger condition that each entry in the global space has exactly one incoming
75895fce210SBarry Smith    edge, it is a candidate for future optimization that might involve its removal.
75995fce210SBarry Smith 
76095fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFGatherBegin(), PetscSFScatterBegin()
76195fce210SBarry Smith @*/
76295fce210SBarry Smith PetscErrorCode PetscSFGetMultiSF(PetscSF sf,PetscSF *multi)
76395fce210SBarry Smith {
76495fce210SBarry Smith   PetscErrorCode ierr;
76595fce210SBarry Smith 
76695fce210SBarry Smith   PetscFunctionBegin;
76795fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
76895fce210SBarry Smith   PetscValidPointer(multi,2);
76995fce210SBarry Smith   if (sf->nroots < 0) {         /* Graph has not been set yet; why do we need this? */
77095fce210SBarry Smith     ierr   = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr);
77195fce210SBarry Smith     *multi = sf->multi;
77295fce210SBarry Smith     PetscFunctionReturn(0);
77395fce210SBarry Smith   }
77495fce210SBarry Smith   if (!sf->multi) {
77595fce210SBarry Smith     const PetscInt *indegree;
7769837ea96SMatthew G. Knepley     PetscInt       i,*inoffset,*outones,*outoffset,maxlocal;
77795fce210SBarry Smith     PetscSFNode    *remote;
77895fce210SBarry Smith     ierr        = PetscSFComputeDegreeBegin(sf,&indegree);CHKERRQ(ierr);
77995fce210SBarry Smith     ierr        = PetscSFComputeDegreeEnd(sf,&indegree);CHKERRQ(ierr);
7809837ea96SMatthew G. Knepley     for (i=0,maxlocal=0; i<sf->nleaves; i++) maxlocal = PetscMax(maxlocal,(sf->mine ? sf->mine[i] : i)+1);
7819837ea96SMatthew G. Knepley     ierr        = PetscMalloc3(sf->nroots+1,&inoffset,maxlocal,&outones,maxlocal,&outoffset);CHKERRQ(ierr);
78295fce210SBarry Smith     inoffset[0] = 0;
78395fce210SBarry Smith     for (i=0; i<sf->nroots; i++) inoffset[i+1] = inoffset[i] + indegree[i];
7849837ea96SMatthew G. Knepley     for (i=0; i<maxlocal; i++) outones[i] = 1;
785*dbd2ff41SBarry Smith     ierr = PetscSFFetchAndOpBegin(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr);
786*dbd2ff41SBarry Smith     ierr = PetscSFFetchAndOpEnd(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr);
78795fce210SBarry Smith     for (i=0; i<sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */
78895fce210SBarry Smith #if 0
78995fce210SBarry Smith #if defined(PETSC_USE_DEBUG)                                 /* Check that the expected number of increments occurred */
79095fce210SBarry Smith     for (i=0; i<sf->nroots; i++) {
79195fce210SBarry Smith       if (inoffset[i] + indegree[i] != inoffset[i+1]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect result after PetscSFFetchAndOp");
79295fce210SBarry Smith     }
79395fce210SBarry Smith #endif
79495fce210SBarry Smith #endif
795785e854fSJed Brown     ierr = PetscMalloc1(sf->nleaves,&remote);CHKERRQ(ierr);
79695fce210SBarry Smith     for (i=0; i<sf->nleaves; i++) {
79795fce210SBarry Smith       remote[i].rank  = sf->remote[i].rank;
79838e7336fSToby Isaac       remote[i].index = outoffset[sf->mine ? sf->mine[i] : i];
79995fce210SBarry Smith     }
80095fce210SBarry Smith     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr);
80195fce210SBarry Smith     ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
80295fce210SBarry Smith     if (sf->rankorder) {        /* Sort the ranks */
80395fce210SBarry Smith       PetscMPIInt rank;
80495fce210SBarry Smith       PetscInt    *inranks,*newoffset,*outranks,*newoutoffset,*tmpoffset,maxdegree;
80595fce210SBarry Smith       PetscSFNode *newremote;
80695fce210SBarry Smith       ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
80795fce210SBarry Smith       for (i=0,maxdegree=0; i<sf->nroots; i++) maxdegree = PetscMax(maxdegree,indegree[i]);
8089837ea96SMatthew G. Knepley       ierr = PetscMalloc5(sf->multi->nroots,&inranks,sf->multi->nroots,&newoffset,maxlocal,&outranks,maxlocal,&newoutoffset,maxdegree,&tmpoffset);CHKERRQ(ierr);
8099837ea96SMatthew G. Knepley       for (i=0; i<maxlocal; i++) outranks[i] = rank;
8108bfbc91cSJed Brown       ierr = PetscSFReduceBegin(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr);
8118bfbc91cSJed Brown       ierr = PetscSFReduceEnd(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr);
81295fce210SBarry Smith       /* Sort the incoming ranks at each vertex, build the inverse map */
81395fce210SBarry Smith       for (i=0; i<sf->nroots; i++) {
81495fce210SBarry Smith         PetscInt j;
81595fce210SBarry Smith         for (j=0; j<indegree[i]; j++) tmpoffset[j] = j;
81695fce210SBarry Smith         ierr = PetscSortIntWithArray(indegree[i],inranks+inoffset[i],tmpoffset);CHKERRQ(ierr);
81795fce210SBarry Smith         for (j=0; j<indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j;
81895fce210SBarry Smith       }
81995fce210SBarry Smith       ierr = PetscSFBcastBegin(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr);
82095fce210SBarry Smith       ierr = PetscSFBcastEnd(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr);
821785e854fSJed Brown       ierr = PetscMalloc1(sf->nleaves,&newremote);CHKERRQ(ierr);
82295fce210SBarry Smith       for (i=0; i<sf->nleaves; i++) {
82395fce210SBarry Smith         newremote[i].rank  = sf->remote[i].rank;
82495fce210SBarry Smith         newremote[i].index = newoutoffset[i];
82595fce210SBarry Smith       }
82695fce210SBarry Smith       ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,NULL,PETSC_COPY_VALUES,newremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
82795fce210SBarry Smith       ierr = PetscFree5(inranks,newoffset,outranks,newoutoffset,tmpoffset);CHKERRQ(ierr);
82895fce210SBarry Smith     }
82995fce210SBarry Smith     ierr = PetscFree3(inoffset,outones,outoffset);CHKERRQ(ierr);
83095fce210SBarry Smith   }
83195fce210SBarry Smith   *multi = sf->multi;
83295fce210SBarry Smith   PetscFunctionReturn(0);
83395fce210SBarry Smith }
83495fce210SBarry Smith 
83595fce210SBarry Smith #undef __FUNCT__
83695fce210SBarry Smith #define __FUNCT__ "PetscSFCreateEmbeddedSF"
83795fce210SBarry Smith /*@C
83895fce210SBarry Smith    PetscSFCreateEmbeddedSF - removes edges from all but the selected roots, does not remap indices
83995fce210SBarry Smith 
84095fce210SBarry Smith    Collective
84195fce210SBarry Smith 
84295fce210SBarry Smith    Input Arguments:
84395fce210SBarry Smith +  sf - original star forest
84495fce210SBarry Smith .  nroots - number of roots to select on this process
84595fce210SBarry Smith -  selected - selected roots on this process
84695fce210SBarry Smith 
84795fce210SBarry Smith    Output Arguments:
84895fce210SBarry Smith .  newsf - new star forest
84995fce210SBarry Smith 
85095fce210SBarry Smith    Level: advanced
85195fce210SBarry Smith 
85295fce210SBarry Smith    Note:
85395fce210SBarry Smith    To use the new PetscSF, it may be necessary to know the indices of the leaves that are still participating. This can
85495fce210SBarry Smith    be done by calling PetscSFGetGraph().
85595fce210SBarry Smith 
85695fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFGetGraph()
85795fce210SBarry Smith @*/
85895fce210SBarry Smith PetscErrorCode PetscSFCreateEmbeddedSF(PetscSF sf,PetscInt nroots,const PetscInt *selected,PetscSF *newsf)
85995fce210SBarry Smith {
8600511a646SMatthew G. Knepley   PetscInt      *rootdata, *leafdata, *ilocal;
86195fce210SBarry Smith   PetscSFNode   *iremote;
862fc1ede2bSMatthew G. Knepley   PetscInt       leafsize = 0, nleaves = 0, n, i;
8630511a646SMatthew G. Knepley   PetscErrorCode ierr;
86495fce210SBarry Smith 
86595fce210SBarry Smith   PetscFunctionBegin;
86695fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
86795fce210SBarry Smith   if (nroots) PetscValidPointer(selected,3);
86895fce210SBarry Smith   PetscValidPointer(newsf,4);
8690511a646SMatthew G. Knepley   if (sf->mine) for (i = 0; i < sf->nleaves; ++i) {leafsize = PetscMax(leafsize, sf->mine[i]+1);}
8700511a646SMatthew G. Knepley   else leafsize = sf->nleaves;
8711795a4d1SJed Brown   ierr = PetscCalloc2(sf->nroots,&rootdata,leafsize,&leafdata);CHKERRQ(ierr);
8720511a646SMatthew G. Knepley   for (i=0; i<nroots; ++i) rootdata[selected[i]] = 1;
87395fce210SBarry Smith   ierr = PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr);
87495fce210SBarry Smith   ierr = PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr);
87595fce210SBarry Smith 
8760511a646SMatthew G. Knepley   for (i = 0; i < leafsize; ++i) nleaves += leafdata[i];
877785e854fSJed Brown   ierr = PetscMalloc1(nleaves,&ilocal);CHKERRQ(ierr);
878785e854fSJed Brown   ierr = PetscMalloc1(nleaves,&iremote);CHKERRQ(ierr);
879fc1ede2bSMatthew G. Knepley   for (i = 0, n = 0; i < sf->nleaves; ++i) {
8800511a646SMatthew G. Knepley     const PetscInt lidx = sf->mine ? sf->mine[i] : i;
8810511a646SMatthew G. Knepley 
8820511a646SMatthew G. Knepley     if (leafdata[lidx]) {
883fc1ede2bSMatthew G. Knepley       ilocal[n]        = lidx;
884fc1ede2bSMatthew G. Knepley       iremote[n].rank  = sf->remote[i].rank;
885fc1ede2bSMatthew G. Knepley       iremote[n].index = sf->remote[i].index;
886fc1ede2bSMatthew G. Knepley       ++n;
88795fce210SBarry Smith     }
88895fce210SBarry Smith   }
889fc1ede2bSMatthew G. Knepley   if (n != nleaves) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "There is a size mismatch in the SF embedding, %d != %d", n, nleaves);
89095fce210SBarry Smith   ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,newsf);CHKERRQ(ierr);
89195fce210SBarry Smith   ierr = PetscSFSetGraph(*newsf,sf->nroots,nleaves,ilocal,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
89295fce210SBarry Smith   ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr);
89395fce210SBarry Smith   PetscFunctionReturn(0);
89495fce210SBarry Smith }
89595fce210SBarry Smith 
89695fce210SBarry Smith #undef __FUNCT__
89795fce210SBarry Smith #define __FUNCT__ "PetscSFBcastBegin"
89895fce210SBarry Smith /*@C
89995fce210SBarry Smith    PetscSFBcastBegin - begin pointwise broadcast to be concluded with call to PetscSFBcastEnd()
90095fce210SBarry Smith 
90195fce210SBarry Smith    Collective on PetscSF
90295fce210SBarry Smith 
90395fce210SBarry Smith    Input Arguments:
90495fce210SBarry Smith +  sf - star forest on which to communicate
90595fce210SBarry Smith .  unit - data type associated with each node
90695fce210SBarry Smith -  rootdata - buffer to broadcast
90795fce210SBarry Smith 
90895fce210SBarry Smith    Output Arguments:
90995fce210SBarry Smith .  leafdata - buffer to update with values from each leaf's respective root
91095fce210SBarry Smith 
91195fce210SBarry Smith    Level: intermediate
91295fce210SBarry Smith 
91395fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFView(), PetscSFBcastEnd(), PetscSFReduceBegin()
91495fce210SBarry Smith @*/
91595fce210SBarry Smith PetscErrorCode PetscSFBcastBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
91695fce210SBarry Smith {
91795fce210SBarry Smith   PetscErrorCode ierr;
91895fce210SBarry Smith 
91995fce210SBarry Smith   PetscFunctionBegin;
92095fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
92195fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
922563ffbabSMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr);
92395fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
92495fce210SBarry Smith   ierr = (*sf->ops->BcastBegin)(sf,unit,rootdata,leafdata);CHKERRQ(ierr);
925563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr);
92695fce210SBarry Smith   PetscFunctionReturn(0);
92795fce210SBarry Smith }
92895fce210SBarry Smith 
92995fce210SBarry Smith #undef __FUNCT__
93095fce210SBarry Smith #define __FUNCT__ "PetscSFBcastEnd"
93195fce210SBarry Smith /*@C
93295fce210SBarry Smith    PetscSFBcastEnd - end a broadcast operation started with PetscSFBcastBegin()
93395fce210SBarry Smith 
93495fce210SBarry Smith    Collective
93595fce210SBarry Smith 
93695fce210SBarry Smith    Input Arguments:
93795fce210SBarry Smith +  sf - star forest
93895fce210SBarry Smith .  unit - data type
93995fce210SBarry Smith -  rootdata - buffer to broadcast
94095fce210SBarry Smith 
94195fce210SBarry Smith    Output Arguments:
94295fce210SBarry Smith .  leafdata - buffer to update with values from each leaf's respective root
94395fce210SBarry Smith 
94495fce210SBarry Smith    Level: intermediate
94595fce210SBarry Smith 
94695fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFReduceEnd()
94795fce210SBarry Smith @*/
94895fce210SBarry Smith PetscErrorCode PetscSFBcastEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
94995fce210SBarry Smith {
95095fce210SBarry Smith   PetscErrorCode ierr;
95195fce210SBarry Smith 
95295fce210SBarry Smith   PetscFunctionBegin;
95395fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
95495fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
955563ffbabSMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCSF_BcastEnd,sf,0,0,0);CHKERRQ(ierr);
95695fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
95795fce210SBarry Smith   ierr = (*sf->ops->BcastEnd)(sf,unit,rootdata,leafdata);CHKERRQ(ierr);
958563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_BcastEnd,sf,0,0,0);CHKERRQ(ierr);
95995fce210SBarry Smith   PetscFunctionReturn(0);
96095fce210SBarry Smith }
96195fce210SBarry Smith 
96295fce210SBarry Smith #undef __FUNCT__
96395fce210SBarry Smith #define __FUNCT__ "PetscSFReduceBegin"
96495fce210SBarry Smith /*@C
96595fce210SBarry Smith    PetscSFReduceBegin - begin reduction of leafdata into rootdata, to be completed with call to PetscSFReduceEnd()
96695fce210SBarry Smith 
96795fce210SBarry Smith    Collective
96895fce210SBarry Smith 
96995fce210SBarry Smith    Input Arguments:
97095fce210SBarry Smith +  sf - star forest
97195fce210SBarry Smith .  unit - data type
97295fce210SBarry Smith .  leafdata - values to reduce
97395fce210SBarry Smith -  op - reduction operation
97495fce210SBarry Smith 
97595fce210SBarry Smith    Output Arguments:
97695fce210SBarry Smith .  rootdata - result of reduction of values from all leaves of each root
97795fce210SBarry Smith 
97895fce210SBarry Smith    Level: intermediate
97995fce210SBarry Smith 
98095fce210SBarry Smith .seealso: PetscSFBcastBegin()
98195fce210SBarry Smith @*/
98295fce210SBarry Smith PetscErrorCode PetscSFReduceBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
98395fce210SBarry Smith {
98495fce210SBarry Smith   PetscErrorCode ierr;
98595fce210SBarry Smith 
98695fce210SBarry Smith   PetscFunctionBegin;
98795fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
98895fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
989563ffbabSMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);
99095fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
99195fce210SBarry Smith   ierr = (sf->ops->ReduceBegin)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
992563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);
99395fce210SBarry Smith   PetscFunctionReturn(0);
99495fce210SBarry Smith }
99595fce210SBarry Smith 
99695fce210SBarry Smith #undef __FUNCT__
99795fce210SBarry Smith #define __FUNCT__ "PetscSFReduceEnd"
99895fce210SBarry Smith /*@C
99995fce210SBarry Smith    PetscSFReduceEnd - end a reduction operation started with PetscSFReduceBegin()
100095fce210SBarry Smith 
100195fce210SBarry Smith    Collective
100295fce210SBarry Smith 
100395fce210SBarry Smith    Input Arguments:
100495fce210SBarry Smith +  sf - star forest
100595fce210SBarry Smith .  unit - data type
100695fce210SBarry Smith .  leafdata - values to reduce
100795fce210SBarry Smith -  op - reduction operation
100895fce210SBarry Smith 
100995fce210SBarry Smith    Output Arguments:
101095fce210SBarry Smith .  rootdata - result of reduction of values from all leaves of each root
101195fce210SBarry Smith 
101295fce210SBarry Smith    Level: intermediate
101395fce210SBarry Smith 
101495fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFBcastEnd()
101595fce210SBarry Smith @*/
101695fce210SBarry Smith PetscErrorCode PetscSFReduceEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
101795fce210SBarry Smith {
101895fce210SBarry Smith   PetscErrorCode ierr;
101995fce210SBarry Smith 
102095fce210SBarry Smith   PetscFunctionBegin;
102195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
102295fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
1023563ffbabSMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);
102495fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
102595fce210SBarry Smith   ierr = (*sf->ops->ReduceEnd)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
1026563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);
102795fce210SBarry Smith   PetscFunctionReturn(0);
102895fce210SBarry Smith }
102995fce210SBarry Smith 
103095fce210SBarry Smith #undef __FUNCT__
103195fce210SBarry Smith #define __FUNCT__ "PetscSFComputeDegreeBegin"
103295fce210SBarry Smith /*@C
103395fce210SBarry Smith    PetscSFComputeDegreeBegin - begin computation of degree for each root vertex, to be completed with PetscSFComputeDegreeEnd()
103495fce210SBarry Smith 
103595fce210SBarry Smith    Collective
103695fce210SBarry Smith 
103795fce210SBarry Smith    Input Arguments:
103895fce210SBarry Smith .  sf - star forest
103995fce210SBarry Smith 
104095fce210SBarry Smith    Output Arguments:
104195fce210SBarry Smith .  degree - degree of each root vertex
104295fce210SBarry Smith 
104395fce210SBarry Smith    Level: advanced
104495fce210SBarry Smith 
104595fce210SBarry Smith .seealso: PetscSFGatherBegin()
104695fce210SBarry Smith @*/
104795fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf,const PetscInt **degree)
104895fce210SBarry Smith {
104995fce210SBarry Smith   PetscErrorCode ierr;
105095fce210SBarry Smith 
105195fce210SBarry Smith   PetscFunctionBegin;
105295fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
105395fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
105495fce210SBarry Smith   PetscValidPointer(degree,2);
1055803bd9e8SMatthew G. Knepley   if (!sf->degreeknown) {
10569837ea96SMatthew G. Knepley     PetscInt i,maxlocal;
1057803bd9e8SMatthew G. Knepley     if (sf->degree) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Calls to PetscSFComputeDegreeBegin() cannot be nested.");
10589837ea96SMatthew G. Knepley     for (i=0,maxlocal=0; i<sf->nleaves; i++) maxlocal = PetscMax(maxlocal,(sf->mine ? sf->mine[i] : i)+1);
1059785e854fSJed Brown     ierr = PetscMalloc1(sf->nroots,&sf->degree);CHKERRQ(ierr);
10609837ea96SMatthew G. Knepley     ierr = PetscMalloc1(maxlocal,&sf->degreetmp);CHKERRQ(ierr);
106195fce210SBarry Smith     for (i=0; i<sf->nroots; i++) sf->degree[i] = 0;
10629837ea96SMatthew G. Knepley     for (i=0; i<maxlocal; i++) sf->degreetmp[i] = 1;
1063*dbd2ff41SBarry Smith     ierr = PetscSFReduceBegin(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);CHKERRQ(ierr);
106495fce210SBarry Smith   }
106595fce210SBarry Smith   *degree = NULL;
106695fce210SBarry Smith   PetscFunctionReturn(0);
106795fce210SBarry Smith }
106895fce210SBarry Smith 
106995fce210SBarry Smith #undef __FUNCT__
107095fce210SBarry Smith #define __FUNCT__ "PetscSFComputeDegreeEnd"
107195fce210SBarry Smith /*@C
107295fce210SBarry Smith    PetscSFComputeDegreeEnd - complete computation of degree for each root vertex, started with PetscSFComputeDegreeBegin()
107395fce210SBarry Smith 
107495fce210SBarry Smith    Collective
107595fce210SBarry Smith 
107695fce210SBarry Smith    Input Arguments:
107795fce210SBarry Smith .  sf - star forest
107895fce210SBarry Smith 
107995fce210SBarry Smith    Output Arguments:
108095fce210SBarry Smith .  degree - degree of each root vertex
108195fce210SBarry Smith 
108295fce210SBarry Smith    Level: developer
108395fce210SBarry Smith 
108495fce210SBarry Smith .seealso:
108595fce210SBarry Smith @*/
108695fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf,const PetscInt **degree)
108795fce210SBarry Smith {
108895fce210SBarry Smith   PetscErrorCode ierr;
108995fce210SBarry Smith 
109095fce210SBarry Smith   PetscFunctionBegin;
109195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
109295fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
109395fce210SBarry Smith   if (!sf->degreeknown) {
1094*dbd2ff41SBarry Smith     ierr = PetscSFReduceEnd(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);CHKERRQ(ierr);
109595fce210SBarry Smith     ierr = PetscFree(sf->degreetmp);CHKERRQ(ierr);
109695fce210SBarry Smith 
109795fce210SBarry Smith     sf->degreeknown = PETSC_TRUE;
109895fce210SBarry Smith   }
109995fce210SBarry Smith   *degree = sf->degree;
110095fce210SBarry Smith   PetscFunctionReturn(0);
110195fce210SBarry Smith }
110295fce210SBarry Smith 
110395fce210SBarry Smith #undef __FUNCT__
110495fce210SBarry Smith #define __FUNCT__ "PetscSFFetchAndOpBegin"
110595fce210SBarry Smith /*@C
110695fce210SBarry Smith    PetscSFFetchAndOpBegin - begin operation that fetches values from root and updates atomically by applying operation using my leaf value, to be completed with PetscSFFetchAndOpEnd()
110795fce210SBarry Smith 
110895fce210SBarry Smith    Collective
110995fce210SBarry Smith 
111095fce210SBarry Smith    Input Arguments:
111195fce210SBarry Smith +  sf - star forest
111295fce210SBarry Smith .  unit - data type
111395fce210SBarry Smith .  leafdata - leaf values to use in reduction
111495fce210SBarry Smith -  op - operation to use for reduction
111595fce210SBarry Smith 
111695fce210SBarry Smith    Output Arguments:
111795fce210SBarry Smith +  rootdata - root values to be updated, input state is seen by first process to perform an update
111895fce210SBarry Smith -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
111995fce210SBarry Smith 
112095fce210SBarry Smith    Level: advanced
112195fce210SBarry Smith 
112295fce210SBarry Smith    Note:
112395fce210SBarry Smith    The update is only atomic at the granularity provided by the hardware. Different roots referenced by the same process
112495fce210SBarry Smith    might be updated in a different order. Furthermore, if a composite type is used for the unit datatype, atomicity is
112595fce210SBarry Smith    not guaranteed across the whole vertex. Therefore, this function is mostly only used with primitive types such as
112695fce210SBarry Smith    integers.
112795fce210SBarry Smith 
112895fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFReduceBegin(), PetscSFSetGraph()
112995fce210SBarry Smith @*/
113095fce210SBarry Smith PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
113195fce210SBarry Smith {
113295fce210SBarry Smith   PetscErrorCode ierr;
113395fce210SBarry Smith 
113495fce210SBarry Smith   PetscFunctionBegin;
113595fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
113695fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
1137563ffbabSMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
113895fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
113995fce210SBarry Smith   ierr = (*sf->ops->FetchAndOpBegin)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr);
1140563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
114195fce210SBarry Smith   PetscFunctionReturn(0);
114295fce210SBarry Smith }
114395fce210SBarry Smith 
114495fce210SBarry Smith #undef __FUNCT__
114595fce210SBarry Smith #define __FUNCT__ "PetscSFFetchAndOpEnd"
114695fce210SBarry Smith /*@C
114795fce210SBarry Smith    PetscSFFetchAndOpEnd - end operation started in matching call to PetscSFFetchAndOpBegin() to fetch values from roots and update atomically by applying operation using my leaf value
114895fce210SBarry Smith 
114995fce210SBarry Smith    Collective
115095fce210SBarry Smith 
115195fce210SBarry Smith    Input Arguments:
115295fce210SBarry Smith +  sf - star forest
115395fce210SBarry Smith .  unit - data type
115495fce210SBarry Smith .  leafdata - leaf values to use in reduction
115595fce210SBarry Smith -  op - operation to use for reduction
115695fce210SBarry Smith 
115795fce210SBarry Smith    Output Arguments:
115895fce210SBarry Smith +  rootdata - root values to be updated, input state is seen by first process to perform an update
115995fce210SBarry Smith -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
116095fce210SBarry Smith 
116195fce210SBarry Smith    Level: advanced
116295fce210SBarry Smith 
116395fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFReduceEnd(), PetscSFSetGraph()
116495fce210SBarry Smith @*/
116595fce210SBarry Smith PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
116695fce210SBarry Smith {
116795fce210SBarry Smith   PetscErrorCode ierr;
116895fce210SBarry Smith 
116995fce210SBarry Smith   PetscFunctionBegin;
117095fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
117195fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
1172563ffbabSMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
117395fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
117495fce210SBarry Smith   ierr = (*sf->ops->FetchAndOpEnd)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr);
1175563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
117695fce210SBarry Smith   PetscFunctionReturn(0);
117795fce210SBarry Smith }
117895fce210SBarry Smith 
117995fce210SBarry Smith #undef __FUNCT__
118095fce210SBarry Smith #define __FUNCT__ "PetscSFGatherBegin"
118195fce210SBarry Smith /*@C
118295fce210SBarry Smith    PetscSFGatherBegin - begin pointwise gather of all leaves into multi-roots, to be completed with PetscSFGatherEnd()
118395fce210SBarry Smith 
118495fce210SBarry Smith    Collective
118595fce210SBarry Smith 
118695fce210SBarry Smith    Input Arguments:
118795fce210SBarry Smith +  sf - star forest
118895fce210SBarry Smith .  unit - data type
118995fce210SBarry Smith -  leafdata - leaf data to gather to roots
119095fce210SBarry Smith 
119195fce210SBarry Smith    Output Argument:
119295fce210SBarry Smith .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
119395fce210SBarry Smith 
119495fce210SBarry Smith    Level: intermediate
119595fce210SBarry Smith 
119695fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
119795fce210SBarry Smith @*/
119895fce210SBarry Smith PetscErrorCode PetscSFGatherBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
119995fce210SBarry Smith {
120095fce210SBarry Smith   PetscErrorCode ierr;
120195fce210SBarry Smith   PetscSF        multi;
120295fce210SBarry Smith 
120395fce210SBarry Smith   PetscFunctionBegin;
120495fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
120595fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
12068bfbc91cSJed Brown   ierr = PetscSFReduceBegin(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr);
120795fce210SBarry Smith   PetscFunctionReturn(0);
120895fce210SBarry Smith }
120995fce210SBarry Smith 
121095fce210SBarry Smith #undef __FUNCT__
121195fce210SBarry Smith #define __FUNCT__ "PetscSFGatherEnd"
121295fce210SBarry Smith /*@C
121395fce210SBarry Smith    PetscSFGatherEnd - ends pointwise gather operation that was started with PetscSFGatherBegin()
121495fce210SBarry Smith 
121595fce210SBarry Smith    Collective
121695fce210SBarry Smith 
121795fce210SBarry Smith    Input Arguments:
121895fce210SBarry Smith +  sf - star forest
121995fce210SBarry Smith .  unit - data type
122095fce210SBarry Smith -  leafdata - leaf data to gather to roots
122195fce210SBarry Smith 
122295fce210SBarry Smith    Output Argument:
122395fce210SBarry Smith .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
122495fce210SBarry Smith 
122595fce210SBarry Smith    Level: intermediate
122695fce210SBarry Smith 
122795fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
122895fce210SBarry Smith @*/
122995fce210SBarry Smith PetscErrorCode PetscSFGatherEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
123095fce210SBarry Smith {
123195fce210SBarry Smith   PetscErrorCode ierr;
123295fce210SBarry Smith   PetscSF        multi;
123395fce210SBarry Smith 
123495fce210SBarry Smith   PetscFunctionBegin;
123595fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
123695fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
123795fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
123895fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
12398bfbc91cSJed Brown   ierr = PetscSFReduceEnd(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr);
124095fce210SBarry Smith   PetscFunctionReturn(0);
124195fce210SBarry Smith }
124295fce210SBarry Smith 
124395fce210SBarry Smith #undef __FUNCT__
124495fce210SBarry Smith #define __FUNCT__ "PetscSFScatterBegin"
124595fce210SBarry Smith /*@C
124695fce210SBarry Smith    PetscSFScatterBegin - begin pointwise scatter operation from multi-roots to leaves, to be completed with PetscSFScatterEnd()
124795fce210SBarry Smith 
124895fce210SBarry Smith    Collective
124995fce210SBarry Smith 
125095fce210SBarry Smith    Input Arguments:
125195fce210SBarry Smith +  sf - star forest
125295fce210SBarry Smith .  unit - data type
125395fce210SBarry Smith -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
125495fce210SBarry Smith 
125595fce210SBarry Smith    Output Argument:
125695fce210SBarry Smith .  leafdata - leaf data to be update with personal data from each respective root
125795fce210SBarry Smith 
125895fce210SBarry Smith    Level: intermediate
125995fce210SBarry Smith 
126095fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
126195fce210SBarry Smith @*/
126295fce210SBarry Smith PetscErrorCode PetscSFScatterBegin(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
126395fce210SBarry Smith {
126495fce210SBarry Smith   PetscErrorCode ierr;
126595fce210SBarry Smith   PetscSF        multi;
126695fce210SBarry Smith 
126795fce210SBarry Smith   PetscFunctionBegin;
126895fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
126995fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
127095fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
127195fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
127295fce210SBarry Smith   ierr = PetscSFBcastBegin(multi,unit,multirootdata,leafdata);CHKERRQ(ierr);
127395fce210SBarry Smith   PetscFunctionReturn(0);
127495fce210SBarry Smith }
127595fce210SBarry Smith 
127695fce210SBarry Smith #undef __FUNCT__
127795fce210SBarry Smith #define __FUNCT__ "PetscSFScatterEnd"
127895fce210SBarry Smith /*@C
127995fce210SBarry Smith    PetscSFScatterEnd - ends pointwise scatter operation that was started with PetscSFScatterBegin()
128095fce210SBarry Smith 
128195fce210SBarry Smith    Collective
128295fce210SBarry Smith 
128395fce210SBarry Smith    Input Arguments:
128495fce210SBarry Smith +  sf - star forest
128595fce210SBarry Smith .  unit - data type
128695fce210SBarry Smith -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
128795fce210SBarry Smith 
128895fce210SBarry Smith    Output Argument:
128995fce210SBarry Smith .  leafdata - leaf data to be update with personal data from each respective root
129095fce210SBarry Smith 
129195fce210SBarry Smith    Level: intermediate
129295fce210SBarry Smith 
129395fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
129495fce210SBarry Smith @*/
129595fce210SBarry Smith PetscErrorCode PetscSFScatterEnd(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
129695fce210SBarry Smith {
129795fce210SBarry Smith   PetscErrorCode ierr;
129895fce210SBarry Smith   PetscSF        multi;
129995fce210SBarry Smith 
130095fce210SBarry Smith   PetscFunctionBegin;
130195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
130295fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
130395fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
130495fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
130595fce210SBarry Smith   ierr = PetscSFBcastEnd(multi,unit,multirootdata,leafdata);CHKERRQ(ierr);
130695fce210SBarry Smith   PetscFunctionReturn(0);
130795fce210SBarry Smith }
1308