xref: /petsc/src/vec/is/sf/interface/sf.c (revision 5b0d146a4c27a34a26fabb00e248697b7b651ba6)
1af0996ceSBarry Smith #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/
2c4e6a40aSLawrence Mitchell #include <petsc/private/hashseti.h>
395fce210SBarry Smith #include <petscctable.h>
495fce210SBarry Smith 
595fce210SBarry Smith #if defined(PETSC_USE_DEBUG)
695fce210SBarry Smith #  define PetscSFCheckGraphSet(sf,arg) do {                          \
795fce210SBarry Smith     if (PetscUnlikely(!(sf)->graphset))                              \
8dd5b3ca6SJunchao Zhang       SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetGraph() or PetscSFSetGraphWithPattern() on argument %D \"%s\" before %s()",(arg),#sf,PETSC_FUNCTION_NAME); \
995fce210SBarry Smith   } while (0)
1095fce210SBarry Smith #else
1195fce210SBarry Smith #  define PetscSFCheckGraphSet(sf,arg) do {} while (0)
1295fce210SBarry Smith #endif
1395fce210SBarry Smith 
1495fce210SBarry Smith const char *const PetscSFDuplicateOptions[] = {"CONFONLY","RANKS","GRAPH","PetscSFDuplicateOption","PETSCSF_DUPLICATE_",0};
1595fce210SBarry Smith 
168af6ec1cSBarry Smith /*@
1795fce210SBarry Smith    PetscSFCreate - create a star forest communication context
1895fce210SBarry Smith 
19d083f849SBarry Smith    Collective
2095fce210SBarry Smith 
2195fce210SBarry Smith    Input Arguments:
2295fce210SBarry Smith .  comm - communicator on which the star forest will operate
2395fce210SBarry Smith 
2495fce210SBarry Smith    Output Arguments:
2595fce210SBarry Smith .  sf - new star forest context
2695fce210SBarry Smith 
27dd5b3ca6SJunchao Zhang    Options Database Keys:
28dd5b3ca6SJunchao Zhang +  -sf_type basic     -Use MPI persistent Isend/Irecv for communication (Default)
29dd5b3ca6SJunchao Zhang .  -sf_type window    -Use MPI-3 one-sided window for communication
30dd5b3ca6SJunchao Zhang -  -sf_type neighbor  -Use MPI-3 neighborhood collectives for communication
31dd5b3ca6SJunchao Zhang 
3295fce210SBarry Smith    Level: intermediate
3395fce210SBarry Smith 
34dd5b3ca6SJunchao Zhang    Notes:
35dd5b3ca6SJunchao Zhang    When one knows the communication graph is one of the predefined graph, such as MPI_Alltoall, MPI_Allgatherv,
36dd5b3ca6SJunchao Zhang    MPI_Gatherv, one can create a PetscSF and then set its graph with PetscSFSetGraphWithPattern(). These special
37dd5b3ca6SJunchao Zhang    SFs are optimized and they have better performance than general SFs.
38dd5b3ca6SJunchao Zhang 
39dd5b3ca6SJunchao Zhang .seealso: PetscSFSetGraph(), PetscSFSetGraphWithPattern(), PetscSFDestroy()
4095fce210SBarry Smith @*/
4195fce210SBarry Smith PetscErrorCode PetscSFCreate(MPI_Comm comm,PetscSF *sf)
4295fce210SBarry Smith {
4395fce210SBarry Smith   PetscErrorCode ierr;
4495fce210SBarry Smith   PetscSF        b;
4595fce210SBarry Smith 
4695fce210SBarry Smith   PetscFunctionBegin;
4795fce210SBarry Smith   PetscValidPointer(sf,2);
48607a6623SBarry Smith   ierr = PetscSFInitializePackage();CHKERRQ(ierr);
4995fce210SBarry Smith 
5073107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(b,PETSCSF_CLASSID,"PetscSF","Star Forest","PetscSF",comm,PetscSFDestroy,PetscSFView);CHKERRQ(ierr);
5195fce210SBarry Smith 
5295fce210SBarry Smith   b->nroots    = -1;
5395fce210SBarry Smith   b->nleaves   = -1;
5429046d53SLisandro Dalcin   b->minleaf   = PETSC_MAX_INT;
5529046d53SLisandro Dalcin   b->maxleaf   = PETSC_MIN_INT;
5695fce210SBarry Smith   b->nranks    = -1;
5795fce210SBarry Smith   b->rankorder = PETSC_TRUE;
5895fce210SBarry Smith   b->ingroup   = MPI_GROUP_NULL;
5995fce210SBarry Smith   b->outgroup  = MPI_GROUP_NULL;
6095fce210SBarry Smith   b->graphset  = PETSC_FALSE;
6195fce210SBarry Smith 
6295fce210SBarry Smith   *sf = b;
6395fce210SBarry Smith   PetscFunctionReturn(0);
6495fce210SBarry Smith }
6595fce210SBarry Smith 
6629046d53SLisandro Dalcin /*@
6795fce210SBarry Smith    PetscSFReset - Reset a star forest so that different sizes or neighbors can be used
6895fce210SBarry Smith 
6995fce210SBarry Smith    Collective
7095fce210SBarry Smith 
7195fce210SBarry Smith    Input Arguments:
7295fce210SBarry Smith .  sf - star forest
7395fce210SBarry Smith 
7495fce210SBarry Smith    Level: advanced
7595fce210SBarry Smith 
7695fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFDestroy()
7795fce210SBarry Smith @*/
7895fce210SBarry Smith PetscErrorCode PetscSFReset(PetscSF sf)
7995fce210SBarry Smith {
8095fce210SBarry Smith   PetscErrorCode ierr;
8195fce210SBarry Smith 
8295fce210SBarry Smith   PetscFunctionBegin;
8395fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
8479715d56SJed Brown   if (sf->ops->Reset) {ierr = (*sf->ops->Reset)(sf);CHKERRQ(ierr);}
8529046d53SLisandro Dalcin   sf->nroots   = -1;
8629046d53SLisandro Dalcin   sf->nleaves  = -1;
8729046d53SLisandro Dalcin   sf->minleaf  = PETSC_MAX_INT;
8829046d53SLisandro Dalcin   sf->maxleaf  = PETSC_MIN_INT;
8995fce210SBarry Smith   sf->mine     = NULL;
9095fce210SBarry Smith   sf->remote   = NULL;
9129046d53SLisandro Dalcin   sf->graphset = PETSC_FALSE;
9229046d53SLisandro Dalcin   ierr = PetscFree(sf->mine_alloc);CHKERRQ(ierr);
9395fce210SBarry Smith   ierr = PetscFree(sf->remote_alloc);CHKERRQ(ierr);
9421c688dcSJed Brown   sf->nranks = -1;
9529046d53SLisandro Dalcin   ierr = PetscFree4(sf->ranks,sf->roffset,sf->rmine,sf->rremote);CHKERRQ(ierr);
96eb02082bSJunchao Zhang #if defined(PETSC_HAVE_CUDA)
97eb02082bSJunchao Zhang   if (sf->rmine_d) {cudaError_t err = cudaFree(sf->rmine_d);CHKERRCUDA(err);sf->rmine_d=NULL;}
98eb02082bSJunchao Zhang #endif
9929046d53SLisandro Dalcin   sf->degreeknown = PETSC_FALSE;
10095fce210SBarry Smith   ierr = PetscFree(sf->degree);CHKERRQ(ierr);
10195fce210SBarry Smith   if (sf->ingroup  != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->ingroup);CHKERRQ(ierr);}
10295fce210SBarry Smith   if (sf->outgroup != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->outgroup);CHKERRQ(ierr);}
10395fce210SBarry Smith   ierr = PetscSFDestroy(&sf->multi);CHKERRQ(ierr);
104dd5b3ca6SJunchao Zhang   ierr = PetscLayoutDestroy(&sf->map);CHKERRQ(ierr);
10595fce210SBarry Smith   sf->setupcalled = PETSC_FALSE;
10695fce210SBarry Smith   PetscFunctionReturn(0);
10795fce210SBarry Smith }
10895fce210SBarry Smith 
10995fce210SBarry Smith /*@C
11029046d53SLisandro Dalcin    PetscSFSetType - Set the PetscSF communication implementation
11195fce210SBarry Smith 
11295fce210SBarry Smith    Collective on PetscSF
11395fce210SBarry Smith 
11495fce210SBarry Smith    Input Parameters:
11595fce210SBarry Smith +  sf - the PetscSF context
11695fce210SBarry Smith -  type - a known method
11795fce210SBarry Smith 
11895fce210SBarry Smith    Options Database Key:
11995fce210SBarry Smith .  -sf_type <type> - Sets the method; use -help for a list
12070616304SStefano Zampini    of available methods (for instance, window, basic, neighbor)
12195fce210SBarry Smith 
12295fce210SBarry Smith    Notes:
12395fce210SBarry Smith    See "include/petscsf.h" for available methods (for instance)
12495fce210SBarry Smith +    PETSCSFWINDOW - MPI-2/3 one-sided
12595fce210SBarry Smith -    PETSCSFBASIC - basic implementation using MPI-1 two-sided
12695fce210SBarry Smith 
12795fce210SBarry Smith   Level: intermediate
12895fce210SBarry Smith 
12995fce210SBarry Smith .seealso: PetscSFType, PetscSFCreate()
13095fce210SBarry Smith @*/
13195fce210SBarry Smith PetscErrorCode PetscSFSetType(PetscSF sf,PetscSFType type)
13295fce210SBarry Smith {
13395fce210SBarry Smith   PetscErrorCode ierr,(*r)(PetscSF);
13495fce210SBarry Smith   PetscBool      match;
13595fce210SBarry Smith 
13695fce210SBarry Smith   PetscFunctionBegin;
13795fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
13895fce210SBarry Smith   PetscValidCharPointer(type,2);
13995fce210SBarry Smith 
14095fce210SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)sf,type,&match);CHKERRQ(ierr);
14195fce210SBarry Smith   if (match) PetscFunctionReturn(0);
14295fce210SBarry Smith 
143adc40e5bSBarry Smith   ierr = PetscFunctionListFind(PetscSFList,type,&r);CHKERRQ(ierr);
14495fce210SBarry Smith   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested PetscSF type %s",type);
14529046d53SLisandro Dalcin   /* Destroy the previous PetscSF implementation context */
14629046d53SLisandro Dalcin   if (sf->ops->Destroy) {ierr = (*(sf)->ops->Destroy)(sf);CHKERRQ(ierr);}
14795fce210SBarry Smith   ierr = PetscMemzero(sf->ops,sizeof(*sf->ops));CHKERRQ(ierr);
14895fce210SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)sf,type);CHKERRQ(ierr);
14995fce210SBarry Smith   ierr = (*r)(sf);CHKERRQ(ierr);
15095fce210SBarry Smith   PetscFunctionReturn(0);
15195fce210SBarry Smith }
15295fce210SBarry Smith 
15329046d53SLisandro Dalcin /*@C
15429046d53SLisandro Dalcin   PetscSFGetType - Get the PetscSF communication implementation
15529046d53SLisandro Dalcin 
15629046d53SLisandro Dalcin   Not Collective
15729046d53SLisandro Dalcin 
15829046d53SLisandro Dalcin   Input Parameter:
15929046d53SLisandro Dalcin . sf  - the PetscSF context
16029046d53SLisandro Dalcin 
16129046d53SLisandro Dalcin   Output Parameter:
16229046d53SLisandro Dalcin . type - the PetscSF type name
16329046d53SLisandro Dalcin 
16429046d53SLisandro Dalcin   Level: intermediate
16529046d53SLisandro Dalcin 
16629046d53SLisandro Dalcin .seealso: PetscSFSetType(), PetscSFCreate()
16729046d53SLisandro Dalcin @*/
16829046d53SLisandro Dalcin PetscErrorCode PetscSFGetType(PetscSF sf, PetscSFType *type)
16929046d53SLisandro Dalcin {
17029046d53SLisandro Dalcin   PetscFunctionBegin;
17129046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID,1);
17229046d53SLisandro Dalcin   PetscValidPointer(type,2);
17329046d53SLisandro Dalcin   *type = ((PetscObject)sf)->type_name;
17429046d53SLisandro Dalcin   PetscFunctionReturn(0);
17529046d53SLisandro Dalcin }
17629046d53SLisandro Dalcin 
177d36d33e4SMatthew G. Knepley /*@
17895fce210SBarry Smith    PetscSFDestroy - destroy star forest
17995fce210SBarry Smith 
18095fce210SBarry Smith    Collective
18195fce210SBarry Smith 
18295fce210SBarry Smith    Input Arguments:
18395fce210SBarry Smith .  sf - address of star forest
18495fce210SBarry Smith 
18595fce210SBarry Smith    Level: intermediate
18695fce210SBarry Smith 
18795fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFReset()
18895fce210SBarry Smith @*/
18995fce210SBarry Smith PetscErrorCode PetscSFDestroy(PetscSF *sf)
19095fce210SBarry Smith {
19195fce210SBarry Smith   PetscErrorCode ierr;
19295fce210SBarry Smith 
19395fce210SBarry Smith   PetscFunctionBegin;
19495fce210SBarry Smith   if (!*sf) PetscFunctionReturn(0);
19595fce210SBarry Smith   PetscValidHeaderSpecific((*sf),PETSCSF_CLASSID,1);
19629046d53SLisandro Dalcin   if (--((PetscObject)(*sf))->refct > 0) {*sf = NULL; PetscFunctionReturn(0);}
19795fce210SBarry Smith   ierr = PetscSFReset(*sf);CHKERRQ(ierr);
19895fce210SBarry Smith   if ((*sf)->ops->Destroy) {ierr = (*(*sf)->ops->Destroy)(*sf);CHKERRQ(ierr);}
19995fce210SBarry Smith   ierr = PetscHeaderDestroy(sf);CHKERRQ(ierr);
20095fce210SBarry Smith   PetscFunctionReturn(0);
20195fce210SBarry Smith }
20295fce210SBarry Smith 
203c4e6a40aSLawrence Mitchell static PetscErrorCode PetscSFCheckGraphValid_Private(PetscSF sf)
204c4e6a40aSLawrence Mitchell {
205c4e6a40aSLawrence Mitchell #if defined(PETSC_USE_DEBUG)
206c4e6a40aSLawrence Mitchell   PetscInt           i, nleaves;
207c4e6a40aSLawrence Mitchell   PetscMPIInt        size;
208c4e6a40aSLawrence Mitchell   const PetscInt    *ilocal;
209c4e6a40aSLawrence Mitchell   const PetscSFNode *iremote;
210c4e6a40aSLawrence Mitchell   PetscErrorCode     ierr;
211c4e6a40aSLawrence Mitchell 
212c4e6a40aSLawrence Mitchell   PetscFunctionBegin;
213c4e6a40aSLawrence Mitchell   if (!sf->graphset) PetscFunctionReturn(0);
214c4e6a40aSLawrence Mitchell   ierr = PetscSFGetGraph(sf,NULL,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
215c4e6a40aSLawrence Mitchell   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr);
216c4e6a40aSLawrence Mitchell   for (i = 0; i < nleaves; i++) {
217c4e6a40aSLawrence Mitchell     const PetscInt rank = iremote[i].rank;
218c4e6a40aSLawrence Mitchell     const PetscInt remote = iremote[i].index;
219c4e6a40aSLawrence Mitchell     const PetscInt leaf = ilocal ? ilocal[i] : i;
220c4e6a40aSLawrence Mitchell     if (rank < 0 || rank >= size) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Provided rank (%D) for remote %D is invalid, should be in [0, %d)",rank,i,size);
221c4e6a40aSLawrence Mitchell     if (remote < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Provided index (%D) for remote %D is invalid, should be >= 0",remote,i);
222c4e6a40aSLawrence Mitchell     if (leaf < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Provided location (%D) for leaf %D is invalid, should be >= 0",leaf,i);
223c4e6a40aSLawrence Mitchell   }
224c4e6a40aSLawrence Mitchell   PetscFunctionReturn(0);
225c4e6a40aSLawrence Mitchell #else
226c4e6a40aSLawrence Mitchell   PetscFunctionBegin;
227c4e6a40aSLawrence Mitchell   PetscFunctionReturn(0);
228c4e6a40aSLawrence Mitchell #endif
229c4e6a40aSLawrence Mitchell }
230c4e6a40aSLawrence Mitchell 
23195fce210SBarry Smith /*@
23295fce210SBarry Smith    PetscSFSetUp - set up communication structures
23395fce210SBarry Smith 
23495fce210SBarry Smith    Collective
23595fce210SBarry Smith 
23695fce210SBarry Smith    Input Arguments:
23795fce210SBarry Smith .  sf - star forest communication object
23895fce210SBarry Smith 
23995fce210SBarry Smith    Level: beginner
24095fce210SBarry Smith 
24195fce210SBarry Smith .seealso: PetscSFSetFromOptions(), PetscSFSetType()
24295fce210SBarry Smith @*/
24395fce210SBarry Smith PetscErrorCode PetscSFSetUp(PetscSF sf)
24495fce210SBarry Smith {
24595fce210SBarry Smith   PetscErrorCode ierr;
24695fce210SBarry Smith 
24795fce210SBarry Smith   PetscFunctionBegin;
24829046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
24929046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
25095fce210SBarry Smith   if (sf->setupcalled) PetscFunctionReturn(0);
251c4e6a40aSLawrence Mitchell   ierr = PetscSFCheckGraphValid_Private(sf);CHKERRQ(ierr);
25295fce210SBarry Smith   if (!((PetscObject)sf)->type_name) {ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr);}
25329046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_SetUp,sf,0,0,0);CHKERRQ(ierr);
25495fce210SBarry Smith   if (sf->ops->SetUp) {ierr = (*sf->ops->SetUp)(sf);CHKERRQ(ierr);}
25529046d53SLisandro Dalcin   ierr = PetscLogEventEnd(PETSCSF_SetUp,sf,0,0,0);CHKERRQ(ierr);
25695fce210SBarry Smith   sf->setupcalled = PETSC_TRUE;
25795fce210SBarry Smith   PetscFunctionReturn(0);
25895fce210SBarry Smith }
25995fce210SBarry Smith 
2608af6ec1cSBarry Smith /*@
26195fce210SBarry Smith    PetscSFSetFromOptions - set PetscSF options using the options database
26295fce210SBarry Smith 
26395fce210SBarry Smith    Logically Collective
26495fce210SBarry Smith 
26595fce210SBarry Smith    Input Arguments:
26695fce210SBarry Smith .  sf - star forest
26795fce210SBarry Smith 
26895fce210SBarry Smith    Options Database Keys:
26960263706SJed Brown +  -sf_type              - implementation type, see PetscSFSetType()
27051ccb202SJunchao Zhang .  -sf_rank_order        - sort composite points for gathers and scatters in rank order, gathers are non-deterministic otherwise
27151ccb202SJunchao Zhang -  -sf_use_pinned_buffer - use pinned (nonpagable) memory for send/recv buffers on host when communicating GPU data but GPU-aware MPI is not used.
27251ccb202SJunchao Zhang                            Only available for SF types of basic and neighbor.
27395fce210SBarry Smith 
27495fce210SBarry Smith    Level: intermediate
27595fce210SBarry Smith @*/
27695fce210SBarry Smith PetscErrorCode PetscSFSetFromOptions(PetscSF sf)
27795fce210SBarry Smith {
27895fce210SBarry Smith   PetscSFType    deft;
27995fce210SBarry Smith   char           type[256];
28095fce210SBarry Smith   PetscErrorCode ierr;
28195fce210SBarry Smith   PetscBool      flg;
28295fce210SBarry Smith 
28395fce210SBarry Smith   PetscFunctionBegin;
28495fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
28595fce210SBarry Smith   ierr = PetscObjectOptionsBegin((PetscObject)sf);CHKERRQ(ierr);
28695fce210SBarry Smith   deft = ((PetscObject)sf)->type_name ? ((PetscObject)sf)->type_name : PETSCSFBASIC;
28729046d53SLisandro Dalcin   ierr = PetscOptionsFList("-sf_type","PetscSF implementation type","PetscSFSetType",PetscSFList,deft,type,sizeof(type),&flg);CHKERRQ(ierr);
28895fce210SBarry Smith   ierr = PetscSFSetType(sf,flg ? type : deft);CHKERRQ(ierr);
28995fce210SBarry 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);
290e55864a3SBarry Smith   if (sf->ops->SetFromOptions) {ierr = (*sf->ops->SetFromOptions)(PetscOptionsObject,sf);CHKERRQ(ierr);}
29195fce210SBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
29295fce210SBarry Smith   PetscFunctionReturn(0);
29395fce210SBarry Smith }
29495fce210SBarry Smith 
29529046d53SLisandro Dalcin /*@
29695fce210SBarry Smith    PetscSFSetRankOrder - sort multi-points for gathers and scatters by rank order
29795fce210SBarry Smith 
29895fce210SBarry Smith    Logically Collective
29995fce210SBarry Smith 
30095fce210SBarry Smith    Input Arguments:
30195fce210SBarry Smith +  sf - star forest
30295fce210SBarry Smith -  flg - PETSC_TRUE to sort, PETSC_FALSE to skip sorting (lower setup cost, but non-deterministic)
30395fce210SBarry Smith 
30495fce210SBarry Smith    Level: advanced
30595fce210SBarry Smith 
30695fce210SBarry Smith .seealso: PetscSFGatherBegin(), PetscSFScatterBegin()
30795fce210SBarry Smith @*/
30895fce210SBarry Smith PetscErrorCode PetscSFSetRankOrder(PetscSF sf,PetscBool flg)
30995fce210SBarry Smith {
31095fce210SBarry Smith   PetscFunctionBegin;
31195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
31295fce210SBarry Smith   PetscValidLogicalCollectiveBool(sf,flg,2);
31395fce210SBarry Smith   if (sf->multi) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Rank ordering must be set before first call to PetscSFGatherBegin() or PetscSFScatterBegin()");
31495fce210SBarry Smith   sf->rankorder = flg;
31595fce210SBarry Smith   PetscFunctionReturn(0);
31695fce210SBarry Smith }
31795fce210SBarry Smith 
3188af6ec1cSBarry Smith /*@
31995fce210SBarry Smith    PetscSFSetGraph - Set a parallel star forest
32095fce210SBarry Smith 
32195fce210SBarry Smith    Collective
32295fce210SBarry Smith 
32395fce210SBarry Smith    Input Arguments:
32495fce210SBarry Smith +  sf - star forest
32595fce210SBarry Smith .  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
32695fce210SBarry Smith .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
327c4e6a40aSLawrence Mitchell .  ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage (locations must be >= 0, enforced
328c4e6a40aSLawrence Mitchell during setup in debug mode)
32995fce210SBarry Smith .  localmode - copy mode for ilocal
330c4e6a40aSLawrence Mitchell .  iremote - remote locations of root vertices for each leaf on the current process (locations must be >= 0, enforced
331c4e6a40aSLawrence Mitchell during setup in debug mode)
33295fce210SBarry Smith -  remotemode - copy mode for iremote
33395fce210SBarry Smith 
33495fce210SBarry Smith    Level: intermediate
33595fce210SBarry Smith 
33695452b02SPatrick Sanan    Notes:
33795452b02SPatrick Sanan     In Fortran you must use PETSC_COPY_VALUES for localmode and remotemode
33838ab3f8aSBarry Smith 
3392ad1e87fSLisandro Dalcin    Developers Note: Local indices which are the identity permutation in the range [0,nleaves) are discarded as they
3402ad1e87fSLisandro Dalcin    encode contiguous storage. In such case, if localmode is PETSC_OWN_POINTER, the memory is deallocated as it is not
3412ad1e87fSLisandro Dalcin    needed
3422ad1e87fSLisandro Dalcin 
343c4e6a40aSLawrence Mitchell    Developers Note: This object does not necessarily encode a true star forest in the graph theoretic sense, since leaf
344c4e6a40aSLawrence Mitchell    indices are not required to be unique. Some functions, however, rely on unique leaf indices (checked in debug mode).
345c4e6a40aSLawrence Mitchell 
34695fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFGetGraph()
34795fce210SBarry Smith @*/
34895fce210SBarry Smith PetscErrorCode PetscSFSetGraph(PetscSF sf,PetscInt nroots,PetscInt nleaves,const PetscInt *ilocal,PetscCopyMode localmode,const PetscSFNode *iremote,PetscCopyMode remotemode)
34995fce210SBarry Smith {
35095fce210SBarry Smith   PetscErrorCode ierr;
35195fce210SBarry Smith 
35295fce210SBarry Smith   PetscFunctionBegin;
35395fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
35429046d53SLisandro Dalcin   if (nleaves > 0 && ilocal) PetscValidIntPointer(ilocal,4);
35529046d53SLisandro Dalcin   if (nleaves > 0) PetscValidPointer(iremote,6);
35629046d53SLisandro Dalcin   if (nroots  < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nroots %D, cannot be negative",nroots);
35795fce210SBarry Smith   if (nleaves < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nleaves %D, cannot be negative",nleaves);
35829046d53SLisandro Dalcin 
3592a67d2daSStefano Zampini   if (sf->nroots >= 0) { /* Reset only if graph already set */
36095fce210SBarry Smith     ierr = PetscSFReset(sf);CHKERRQ(ierr);
3612a67d2daSStefano Zampini   }
3622a67d2daSStefano Zampini 
36329046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr);
36429046d53SLisandro Dalcin 
36595fce210SBarry Smith   sf->nroots  = nroots;
36695fce210SBarry Smith   sf->nleaves = nleaves;
36729046d53SLisandro Dalcin 
36829046d53SLisandro Dalcin   if (nleaves && ilocal) {
36921c688dcSJed Brown     PetscInt i;
37029046d53SLisandro Dalcin     PetscInt minleaf = PETSC_MAX_INT;
37129046d53SLisandro Dalcin     PetscInt maxleaf = PETSC_MIN_INT;
3722ad1e87fSLisandro Dalcin     int      contiguous = 1;
37329046d53SLisandro Dalcin     for (i=0; i<nleaves; i++) {
37429046d53SLisandro Dalcin       minleaf = PetscMin(minleaf,ilocal[i]);
37529046d53SLisandro Dalcin       maxleaf = PetscMax(maxleaf,ilocal[i]);
3762ad1e87fSLisandro Dalcin       contiguous &= (ilocal[i] == i);
37729046d53SLisandro Dalcin     }
37829046d53SLisandro Dalcin     sf->minleaf = minleaf;
37929046d53SLisandro Dalcin     sf->maxleaf = maxleaf;
3802ad1e87fSLisandro Dalcin     if (contiguous) {
3812ad1e87fSLisandro Dalcin       if (localmode == PETSC_OWN_POINTER) {
3822ad1e87fSLisandro Dalcin         ierr = PetscFree(ilocal);CHKERRQ(ierr);
3832ad1e87fSLisandro Dalcin       }
3842ad1e87fSLisandro Dalcin       ilocal = NULL;
3852ad1e87fSLisandro Dalcin     }
38629046d53SLisandro Dalcin   } else {
38729046d53SLisandro Dalcin     sf->minleaf = 0;
38829046d53SLisandro Dalcin     sf->maxleaf = nleaves - 1;
38929046d53SLisandro Dalcin   }
39029046d53SLisandro Dalcin 
39129046d53SLisandro Dalcin   if (ilocal) {
39295fce210SBarry Smith     switch (localmode) {
39395fce210SBarry Smith     case PETSC_COPY_VALUES:
394785e854fSJed Brown       ierr = PetscMalloc1(nleaves,&sf->mine_alloc);CHKERRQ(ierr);
395580bdb30SBarry Smith       ierr = PetscArraycpy(sf->mine_alloc,ilocal,nleaves);CHKERRQ(ierr);
39695fce210SBarry Smith       sf->mine = sf->mine_alloc;
39795fce210SBarry Smith       break;
39895fce210SBarry Smith     case PETSC_OWN_POINTER:
39995fce210SBarry Smith       sf->mine_alloc = (PetscInt*)ilocal;
40095fce210SBarry Smith       sf->mine       = sf->mine_alloc;
40195fce210SBarry Smith       break;
40295fce210SBarry Smith     case PETSC_USE_POINTER:
40329046d53SLisandro Dalcin       sf->mine_alloc = NULL;
40495fce210SBarry Smith       sf->mine       = (PetscInt*)ilocal;
40595fce210SBarry Smith       break;
40695fce210SBarry Smith     default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown localmode");
40795fce210SBarry Smith     }
40895fce210SBarry Smith   }
40929046d53SLisandro Dalcin 
41095fce210SBarry Smith   switch (remotemode) {
41195fce210SBarry Smith   case PETSC_COPY_VALUES:
412785e854fSJed Brown     ierr = PetscMalloc1(nleaves,&sf->remote_alloc);CHKERRQ(ierr);
413580bdb30SBarry Smith     ierr = PetscArraycpy(sf->remote_alloc,iremote,nleaves);CHKERRQ(ierr);
41495fce210SBarry Smith     sf->remote = sf->remote_alloc;
41595fce210SBarry Smith     break;
41695fce210SBarry Smith   case PETSC_OWN_POINTER:
41795fce210SBarry Smith     sf->remote_alloc = (PetscSFNode*)iremote;
41895fce210SBarry Smith     sf->remote       = sf->remote_alloc;
41995fce210SBarry Smith     break;
42095fce210SBarry Smith   case PETSC_USE_POINTER:
42129046d53SLisandro Dalcin     sf->remote_alloc = NULL;
42295fce210SBarry Smith     sf->remote       = (PetscSFNode*)iremote;
42395fce210SBarry Smith     break;
42495fce210SBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown remotemode");
42595fce210SBarry Smith   }
42695fce210SBarry Smith 
427563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr);
42829046d53SLisandro Dalcin   sf->graphset = PETSC_TRUE;
42995fce210SBarry Smith   PetscFunctionReturn(0);
43095fce210SBarry Smith }
43195fce210SBarry Smith 
43229046d53SLisandro Dalcin /*@
433dd5b3ca6SJunchao Zhang   PetscSFSetGraphWithPattern - Sets the graph of an SF with a specific pattern
434dd5b3ca6SJunchao Zhang 
435dd5b3ca6SJunchao Zhang   Collective
436dd5b3ca6SJunchao Zhang 
437dd5b3ca6SJunchao Zhang   Input Parameters:
438dd5b3ca6SJunchao Zhang + sf      - The PetscSF
439dd5b3ca6SJunchao Zhang . map     - Layout of roots over all processes (insignificant when pattern is PETSCSF_PATTERN_ALLTOALL)
440dd5b3ca6SJunchao Zhang - pattern - One of PETSCSF_PATTERN_ALLGATHER, PETSCSF_PATTERN_GATHER, PETSCSF_PATTERN_ALLTOALL
441dd5b3ca6SJunchao Zhang 
442dd5b3ca6SJunchao Zhang   Notes:
443dd5b3ca6SJunchao Zhang   It is easier to explain PetscSFPattern using vectors. Suppose we have an MPI vector x and its layout is map.
444dd5b3ca6SJunchao Zhang   n and N are local and global sizes of x respectively.
445dd5b3ca6SJunchao Zhang 
446dd5b3ca6SJunchao Zhang   With PETSCSF_PATTERN_ALLGATHER, the routine creates a graph that if one does Bcast on it, it will copy x to
447dd5b3ca6SJunchao Zhang   sequential vectors y on all ranks.
448dd5b3ca6SJunchao Zhang 
449dd5b3ca6SJunchao Zhang   With PETSCSF_PATTERN_GATHER, the routine creates a graph that if one does Bcast on it, it will copy x to a
450dd5b3ca6SJunchao Zhang   sequential vector y on rank 0.
451dd5b3ca6SJunchao Zhang 
452dd5b3ca6SJunchao Zhang   In above cases, entries of x are roots and entries of y are leaves.
453dd5b3ca6SJunchao Zhang 
454dd5b3ca6SJunchao Zhang   With PETSCSF_PATTERN_ALLTOALL, map is insignificant. Suppose NP is size of sf's communicator. The routine
455dd5b3ca6SJunchao Zhang   creates a graph that every rank has NP leaves and NP roots. On rank i, its leaf j is connected to root i
456dd5b3ca6SJunchao Zhang   of rank j. Here 0 <=i,j<NP. It is a kind of MPI_Alltoall with sendcount/recvcount being 1. Note that it does
457dd5b3ca6SJunchao Zhang   not mean one can not send multiple items. One just needs to create a new MPI datatype for the mulptiple data
458dd5b3ca6SJunchao Zhang   items with MPI_Type_contiguous() and use that as the <unit> argument in SF routines.
459dd5b3ca6SJunchao Zhang 
460dd5b3ca6SJunchao Zhang   In this case, roots and leaves are symmetric.
461dd5b3ca6SJunchao Zhang 
462dd5b3ca6SJunchao Zhang   Level: intermediate
463dd5b3ca6SJunchao Zhang  @*/
464dd5b3ca6SJunchao Zhang PetscErrorCode PetscSFSetGraphWithPattern(PetscSF sf,PetscLayout map,PetscSFPattern pattern)
465dd5b3ca6SJunchao Zhang {
466dd5b3ca6SJunchao Zhang   MPI_Comm       comm;
467dd5b3ca6SJunchao Zhang   PetscInt       n,N,res[2];
468dd5b3ca6SJunchao Zhang   PetscMPIInt    rank,size;
469dd5b3ca6SJunchao Zhang   PetscSFType    type;
470dd5b3ca6SJunchao Zhang   PetscErrorCode ierr;
471dd5b3ca6SJunchao Zhang 
472dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
473dd5b3ca6SJunchao Zhang   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
474dd5b3ca6SJunchao Zhang   if (pattern < PETSCSF_PATTERN_ALLGATHER || pattern > PETSCSF_PATTERN_ALLTOALL) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Unsupported PetscSFPattern %D\n",pattern);
475dd5b3ca6SJunchao Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
476dd5b3ca6SJunchao Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
477dd5b3ca6SJunchao Zhang 
478dd5b3ca6SJunchao Zhang   if (pattern == PETSCSF_PATTERN_ALLTOALL) {
479dd5b3ca6SJunchao Zhang     type = PETSCSFALLTOALL;
480dd5b3ca6SJunchao Zhang     ierr = PetscLayoutCreate(comm,&sf->map);CHKERRQ(ierr);
481dd5b3ca6SJunchao Zhang     ierr = PetscLayoutSetLocalSize(sf->map,size);CHKERRQ(ierr);
482dd5b3ca6SJunchao Zhang     ierr = PetscLayoutSetSize(sf->map,((PetscInt)size)*size);CHKERRQ(ierr);
483dd5b3ca6SJunchao Zhang     ierr = PetscLayoutSetUp(sf->map);CHKERRQ(ierr);
484dd5b3ca6SJunchao Zhang   } else {
485dd5b3ca6SJunchao Zhang     ierr   = PetscLayoutGetLocalSize(map,&n);CHKERRQ(ierr);
486dd5b3ca6SJunchao Zhang     ierr   = PetscLayoutGetSize(map,&N);CHKERRQ(ierr);
487dd5b3ca6SJunchao Zhang     res[0] = n;
488dd5b3ca6SJunchao Zhang     res[1] = -n;
489dd5b3ca6SJunchao Zhang     /* Check if n are same over all ranks so that we can optimize it */
490dd5b3ca6SJunchao Zhang     ierr   = MPIU_Allreduce(MPI_IN_PLACE,res,2,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr);
491dd5b3ca6SJunchao Zhang     if (res[0] == -res[1]) { /* same n */
492dd5b3ca6SJunchao Zhang       type = (pattern == PETSCSF_PATTERN_ALLGATHER) ? PETSCSFALLGATHER  : PETSCSFGATHER;
493dd5b3ca6SJunchao Zhang     } else {
494dd5b3ca6SJunchao Zhang       type = (pattern == PETSCSF_PATTERN_ALLGATHER) ? PETSCSFALLGATHERV : PETSCSFGATHERV;
495dd5b3ca6SJunchao Zhang     }
496dd5b3ca6SJunchao Zhang     ierr = PetscLayoutReference(map,&sf->map);CHKERRQ(ierr);
497dd5b3ca6SJunchao Zhang   }
498dd5b3ca6SJunchao Zhang   ierr = PetscSFSetType(sf,type);CHKERRQ(ierr);
499dd5b3ca6SJunchao Zhang 
500dd5b3ca6SJunchao Zhang   sf->pattern = pattern;
501dd5b3ca6SJunchao Zhang   sf->mine    = NULL; /* Contiguous */
502dd5b3ca6SJunchao Zhang 
503dd5b3ca6SJunchao Zhang   /* Set nleaves, nroots here in case user calls PetscSFGetGraph, which is legal to call even before PetscSFSetUp is called.
504dd5b3ca6SJunchao Zhang      Also set other easy stuff.
505dd5b3ca6SJunchao Zhang    */
506dd5b3ca6SJunchao Zhang   if (pattern == PETSCSF_PATTERN_ALLGATHER) {
507dd5b3ca6SJunchao Zhang     sf->nleaves      = N;
508dd5b3ca6SJunchao Zhang     sf->nroots       = n;
509dd5b3ca6SJunchao Zhang     sf->nranks       = size;
510dd5b3ca6SJunchao Zhang     sf->minleaf      = 0;
511dd5b3ca6SJunchao Zhang     sf->maxleaf      = N - 1;
512dd5b3ca6SJunchao Zhang   } else if (pattern == PETSCSF_PATTERN_GATHER) {
513dd5b3ca6SJunchao Zhang     sf->nleaves      = rank ? 0 : N;
514dd5b3ca6SJunchao Zhang     sf->nroots       = n;
515dd5b3ca6SJunchao Zhang     sf->nranks       = rank ? 0 : size;
516dd5b3ca6SJunchao Zhang     sf->minleaf      = 0;
517dd5b3ca6SJunchao Zhang     sf->maxleaf      = rank ? -1 : N - 1;
518dd5b3ca6SJunchao Zhang   } else if (pattern == PETSCSF_PATTERN_ALLTOALL) {
519dd5b3ca6SJunchao Zhang     sf->nleaves      = size;
520dd5b3ca6SJunchao Zhang     sf->nroots       = size;
521dd5b3ca6SJunchao Zhang     sf->nranks       = size;
522dd5b3ca6SJunchao Zhang     sf->minleaf      = 0;
523dd5b3ca6SJunchao Zhang     sf->maxleaf      = size - 1;
524dd5b3ca6SJunchao Zhang   }
525dd5b3ca6SJunchao Zhang   sf->ndranks  = 0; /* We do not need to separate out distinguished ranks for patterned graphs to improve communication performance */
526dd5b3ca6SJunchao Zhang   sf->graphset = PETSC_TRUE;
527dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
528dd5b3ca6SJunchao Zhang }
529dd5b3ca6SJunchao Zhang 
530dd5b3ca6SJunchao Zhang /*@
53195fce210SBarry Smith    PetscSFCreateInverseSF - given a PetscSF in which all vertices have degree 1, creates the inverse map
53295fce210SBarry Smith 
53395fce210SBarry Smith    Collective
53495fce210SBarry Smith 
53595fce210SBarry Smith    Input Arguments:
536dd5b3ca6SJunchao Zhang 
53795fce210SBarry Smith .  sf - star forest to invert
53895fce210SBarry Smith 
53995fce210SBarry Smith    Output Arguments:
54095fce210SBarry Smith .  isf - inverse of sf
54195fce210SBarry Smith    Level: advanced
54295fce210SBarry Smith 
54395fce210SBarry Smith    Notes:
54495fce210SBarry Smith    All roots must have degree 1.
54595fce210SBarry Smith 
54695fce210SBarry Smith    The local space may be a permutation, but cannot be sparse.
54795fce210SBarry Smith 
54895fce210SBarry Smith .seealso: PetscSFSetGraph()
54995fce210SBarry Smith @*/
55095fce210SBarry Smith PetscErrorCode PetscSFCreateInverseSF(PetscSF sf,PetscSF *isf)
55195fce210SBarry Smith {
55295fce210SBarry Smith   PetscErrorCode ierr;
55395fce210SBarry Smith   PetscMPIInt    rank;
55495fce210SBarry Smith   PetscInt       i,nroots,nleaves,maxlocal,count,*newilocal;
55595fce210SBarry Smith   const PetscInt *ilocal;
55695fce210SBarry Smith   PetscSFNode    *roots,*leaves;
55795fce210SBarry Smith 
55895fce210SBarry Smith   PetscFunctionBegin;
55929046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
56029046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
56129046d53SLisandro Dalcin   PetscValidPointer(isf,2);
56229046d53SLisandro Dalcin 
56395fce210SBarry Smith   ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);CHKERRQ(ierr);
56429046d53SLisandro Dalcin   maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */
56529046d53SLisandro Dalcin 
56629046d53SLisandro Dalcin   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
567ae9aee6dSMatthew G. Knepley   ierr = PetscMalloc2(nroots,&roots,maxlocal,&leaves);CHKERRQ(ierr);
568ae9aee6dSMatthew G. Knepley   for (i=0; i<maxlocal; i++) {
56995fce210SBarry Smith     leaves[i].rank  = rank;
57095fce210SBarry Smith     leaves[i].index = i;
57195fce210SBarry Smith   }
57295fce210SBarry Smith   for (i=0; i <nroots; i++) {
57395fce210SBarry Smith     roots[i].rank  = -1;
57495fce210SBarry Smith     roots[i].index = -1;
57595fce210SBarry Smith   }
5768bfbc91cSJed Brown   ierr = PetscSFReduceBegin(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr);
5778bfbc91cSJed Brown   ierr = PetscSFReduceEnd(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr);
57895fce210SBarry Smith 
57995fce210SBarry Smith   /* Check whether our leaves are sparse */
58095fce210SBarry Smith   for (i=0,count=0; i<nroots; i++) if (roots[i].rank >= 0) count++;
58195fce210SBarry Smith   if (count == nroots) newilocal = NULL;
58295fce210SBarry Smith   else {                        /* Index for sparse leaves and compact "roots" array (which is to become our leaves). */
583785e854fSJed Brown     ierr = PetscMalloc1(count,&newilocal);CHKERRQ(ierr);
58495fce210SBarry Smith     for (i=0,count=0; i<nroots; i++) {
58595fce210SBarry Smith       if (roots[i].rank >= 0) {
58695fce210SBarry Smith         newilocal[count]   = i;
58795fce210SBarry Smith         roots[count].rank  = roots[i].rank;
58895fce210SBarry Smith         roots[count].index = roots[i].index;
58995fce210SBarry Smith         count++;
59095fce210SBarry Smith       }
59195fce210SBarry Smith     }
59295fce210SBarry Smith   }
59395fce210SBarry Smith 
59495fce210SBarry Smith   ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,isf);CHKERRQ(ierr);
59595fce210SBarry Smith   ierr = PetscSFSetGraph(*isf,maxlocal,count,newilocal,PETSC_OWN_POINTER,roots,PETSC_COPY_VALUES);CHKERRQ(ierr);
59695fce210SBarry Smith   ierr = PetscFree2(roots,leaves);CHKERRQ(ierr);
59795fce210SBarry Smith   PetscFunctionReturn(0);
59895fce210SBarry Smith }
59995fce210SBarry Smith 
60095fce210SBarry Smith /*@
60195fce210SBarry Smith    PetscSFDuplicate - duplicate a PetscSF, optionally preserving rank connectivity and graph
60295fce210SBarry Smith 
60395fce210SBarry Smith    Collective
60495fce210SBarry Smith 
60595fce210SBarry Smith    Input Arguments:
60695fce210SBarry Smith +  sf - communication object to duplicate
60795fce210SBarry Smith -  opt - PETSCSF_DUPLICATE_CONFONLY, PETSCSF_DUPLICATE_RANKS, or PETSCSF_DUPLICATE_GRAPH (see PetscSFDuplicateOption)
60895fce210SBarry Smith 
60995fce210SBarry Smith    Output Arguments:
61095fce210SBarry Smith .  newsf - new communication object
61195fce210SBarry Smith 
61295fce210SBarry Smith    Level: beginner
61395fce210SBarry Smith 
61495fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetType(), PetscSFSetGraph()
61595fce210SBarry Smith @*/
61695fce210SBarry Smith PetscErrorCode PetscSFDuplicate(PetscSF sf,PetscSFDuplicateOption opt,PetscSF *newsf)
61795fce210SBarry Smith {
61829046d53SLisandro Dalcin   PetscSFType    type;
61995fce210SBarry Smith   PetscErrorCode ierr;
62095fce210SBarry Smith 
62195fce210SBarry Smith   PetscFunctionBegin;
62229046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
62329046d53SLisandro Dalcin   PetscValidLogicalCollectiveEnum(sf,opt,2);
62429046d53SLisandro Dalcin   PetscValidPointer(newsf,3);
62595fce210SBarry Smith   ierr = PetscSFCreate(PetscObjectComm((PetscObject)sf),newsf);CHKERRQ(ierr);
62629046d53SLisandro Dalcin   ierr = PetscSFGetType(sf,&type);CHKERRQ(ierr);
62729046d53SLisandro Dalcin   if (type) {ierr = PetscSFSetType(*newsf,type);CHKERRQ(ierr);}
62895fce210SBarry Smith   if (opt == PETSCSF_DUPLICATE_GRAPH) {
629dd5b3ca6SJunchao Zhang     PetscSFCheckGraphSet(sf,1);
630dd5b3ca6SJunchao Zhang     if (sf->pattern == PETSCSF_PATTERN_GENERAL) {
63195fce210SBarry Smith       PetscInt          nroots,nleaves;
63295fce210SBarry Smith       const PetscInt    *ilocal;
63395fce210SBarry Smith       const PetscSFNode *iremote;
63495fce210SBarry Smith       ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
63595fce210SBarry Smith       ierr = PetscSFSetGraph(*newsf,nroots,nleaves,ilocal,PETSC_COPY_VALUES,iremote,PETSC_COPY_VALUES);CHKERRQ(ierr);
636dd5b3ca6SJunchao Zhang     } else {
637dd5b3ca6SJunchao Zhang       ierr = PetscSFSetGraphWithPattern(*newsf,sf->map,sf->pattern);CHKERRQ(ierr);
638dd5b3ca6SJunchao Zhang     }
63995fce210SBarry Smith   }
64029046d53SLisandro Dalcin   if (sf->ops->Duplicate) {ierr = (*sf->ops->Duplicate)(sf,opt,*newsf);CHKERRQ(ierr);}
64195fce210SBarry Smith   PetscFunctionReturn(0);
64295fce210SBarry Smith }
64395fce210SBarry Smith 
64495fce210SBarry Smith /*@C
64595fce210SBarry Smith    PetscSFGetGraph - Get the graph specifying a parallel star forest
64695fce210SBarry Smith 
64795fce210SBarry Smith    Not Collective
64895fce210SBarry Smith 
64995fce210SBarry Smith    Input Arguments:
65095fce210SBarry Smith .  sf - star forest
65195fce210SBarry Smith 
65295fce210SBarry Smith    Output Arguments:
65395fce210SBarry Smith +  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
65495fce210SBarry Smith .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
65595fce210SBarry Smith .  ilocal - locations of leaves in leafdata buffers
65695fce210SBarry Smith -  iremote - remote locations of root vertices for each leaf on the current process
65795fce210SBarry Smith 
658373e0d91SLisandro Dalcin    Notes:
659373e0d91SLisandro Dalcin    We are not currently requiring that the graph is set, thus returning nroots=-1 if it has not been set yet
660373e0d91SLisandro Dalcin 
661ca797d7aSLawrence Mitchell    When called from Fortran, the returned iremote array is a copy and must deallocated after use. Consequently, if you
662ca797d7aSLawrence Mitchell    want to update the graph, you must call PetscSFSetGraph after modifying the iremote array.
663ca797d7aSLawrence Mitchell 
66495fce210SBarry Smith    Level: intermediate
66595fce210SBarry Smith 
66695fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph()
66795fce210SBarry Smith @*/
66895fce210SBarry Smith PetscErrorCode PetscSFGetGraph(PetscSF sf,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote)
66995fce210SBarry Smith {
670b8dee149SJunchao Zhang   PetscErrorCode ierr;
67195fce210SBarry Smith 
67295fce210SBarry Smith   PetscFunctionBegin;
67395fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
674b8dee149SJunchao Zhang   if (sf->ops->GetGraph) {
675b8dee149SJunchao Zhang     ierr = (sf->ops->GetGraph)(sf,nroots,nleaves,ilocal,iremote);CHKERRQ(ierr);
676b8dee149SJunchao Zhang   } else {
67795fce210SBarry Smith     if (nroots) *nroots = sf->nroots;
67895fce210SBarry Smith     if (nleaves) *nleaves = sf->nleaves;
67995fce210SBarry Smith     if (ilocal) *ilocal = sf->mine;
68095fce210SBarry Smith     if (iremote) *iremote = sf->remote;
681b8dee149SJunchao Zhang   }
68295fce210SBarry Smith   PetscFunctionReturn(0);
68395fce210SBarry Smith }
68495fce210SBarry Smith 
68529046d53SLisandro Dalcin /*@
68695fce210SBarry Smith    PetscSFGetLeafRange - Get the active leaf ranges
68795fce210SBarry Smith 
68895fce210SBarry Smith    Not Collective
68995fce210SBarry Smith 
69095fce210SBarry Smith    Input Arguments:
69195fce210SBarry Smith .  sf - star forest
69295fce210SBarry Smith 
69395fce210SBarry Smith    Output Arguments:
694dd5b3ca6SJunchao Zhang +  minleaf - minimum active leaf on this process. Return 0 if there are no leaves.
695dd5b3ca6SJunchao Zhang -  maxleaf - maximum active leaf on this process. Return -1 if there are no leaves.
69695fce210SBarry Smith 
69795fce210SBarry Smith    Level: developer
69895fce210SBarry Smith 
69995fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph(), PetscSFGetGraph()
70095fce210SBarry Smith @*/
70195fce210SBarry Smith PetscErrorCode PetscSFGetLeafRange(PetscSF sf,PetscInt *minleaf,PetscInt *maxleaf)
70295fce210SBarry Smith {
70395fce210SBarry Smith 
70495fce210SBarry Smith   PetscFunctionBegin;
70595fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
70629046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
70795fce210SBarry Smith   if (minleaf) *minleaf = sf->minleaf;
70895fce210SBarry Smith   if (maxleaf) *maxleaf = sf->maxleaf;
70995fce210SBarry Smith   PetscFunctionReturn(0);
71095fce210SBarry Smith }
71195fce210SBarry Smith 
71295fce210SBarry Smith /*@C
713fe2efc57SMark    PetscSFViewFromOptions - View from Options
714fe2efc57SMark 
715fe2efc57SMark    Collective on PetscSF
716fe2efc57SMark 
717fe2efc57SMark    Input Parameters:
718fe2efc57SMark +  A - the star forest
719736c3998SJose E. Roman .  obj - Optional object
720736c3998SJose E. Roman -  name - command line option
721fe2efc57SMark 
722fe2efc57SMark    Level: intermediate
723fe2efc57SMark .seealso:  PetscSF, PetscSFView, PetscObjectViewFromOptions(), PetscSFCreate()
724fe2efc57SMark @*/
725fe2efc57SMark PetscErrorCode  PetscSFViewFromOptions(PetscSF A,PetscObject obj,const char name[])
726fe2efc57SMark {
727fe2efc57SMark   PetscErrorCode ierr;
728fe2efc57SMark 
729fe2efc57SMark   PetscFunctionBegin;
730fe2efc57SMark   PetscValidHeaderSpecific(A,PETSCSF_CLASSID,1);
731fe2efc57SMark   ierr = PetscObjectViewFromOptions((PetscObject)A,obj,name);CHKERRQ(ierr);
732fe2efc57SMark   PetscFunctionReturn(0);
733fe2efc57SMark }
734fe2efc57SMark 
735fe2efc57SMark /*@C
73695fce210SBarry Smith    PetscSFView - view a star forest
73795fce210SBarry Smith 
73895fce210SBarry Smith    Collective
73995fce210SBarry Smith 
74095fce210SBarry Smith    Input Arguments:
74195fce210SBarry Smith +  sf - star forest
74295fce210SBarry Smith -  viewer - viewer to display graph, for example PETSC_VIEWER_STDOUT_WORLD
74395fce210SBarry Smith 
74495fce210SBarry Smith    Level: beginner
74595fce210SBarry Smith 
74695fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph()
74795fce210SBarry Smith @*/
74895fce210SBarry Smith PetscErrorCode PetscSFView(PetscSF sf,PetscViewer viewer)
74995fce210SBarry Smith {
75095fce210SBarry Smith   PetscErrorCode    ierr;
75195fce210SBarry Smith   PetscBool         iascii;
75295fce210SBarry Smith   PetscViewerFormat format;
75395fce210SBarry Smith 
75495fce210SBarry Smith   PetscFunctionBegin;
75595fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
75695fce210SBarry Smith   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sf),&viewer);CHKERRQ(ierr);}
75795fce210SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
75895fce210SBarry Smith   PetscCheckSameComm(sf,1,viewer,2);
75980153354SVaclav Hapla   if (sf->graphset) {ierr = PetscSFSetUp(sf);CHKERRQ(ierr);}
76095fce210SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
76195fce210SBarry Smith   if (iascii) {
76295fce210SBarry Smith     PetscMPIInt rank;
76381bfa7aaSJed Brown     PetscInt    ii,i,j;
76495fce210SBarry Smith 
765dae58748SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)sf,viewer);CHKERRQ(ierr);
76695fce210SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
76795fce210SBarry Smith     if (sf->ops->View) {ierr = (*sf->ops->View)(sf,viewer);CHKERRQ(ierr);}
768dd5b3ca6SJunchao Zhang     if (sf->pattern == PETSCSF_PATTERN_GENERAL) {
76980153354SVaclav Hapla       if (!sf->graphset) {
77080153354SVaclav Hapla         ierr = PetscViewerASCIIPrintf(viewer,"PetscSFSetGraph() has not been called yet\n");CHKERRQ(ierr);
77180153354SVaclav Hapla         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
77280153354SVaclav Hapla         PetscFunctionReturn(0);
77380153354SVaclav Hapla       }
77495fce210SBarry Smith       ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
7751575c14dSBarry Smith       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
77695fce210SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of roots=%D, leaves=%D, remote ranks=%D\n",rank,sf->nroots,sf->nleaves,sf->nranks);CHKERRQ(ierr);
77795fce210SBarry Smith       for (i=0; i<sf->nleaves; i++) {
77895fce210SBarry 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);
77995fce210SBarry Smith       }
78095fce210SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
78195fce210SBarry Smith       ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
78295fce210SBarry Smith       if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
78381bfa7aaSJed Brown         PetscMPIInt *tmpranks,*perm;
78481bfa7aaSJed Brown         ierr = PetscMalloc2(sf->nranks,&tmpranks,sf->nranks,&perm);CHKERRQ(ierr);
785580bdb30SBarry Smith         ierr = PetscArraycpy(tmpranks,sf->ranks,sf->nranks);CHKERRQ(ierr);
78681bfa7aaSJed Brown         for (i=0; i<sf->nranks; i++) perm[i] = i;
78781bfa7aaSJed Brown         ierr = PetscSortMPIIntWithArray(sf->nranks,tmpranks,perm);CHKERRQ(ierr);
78895fce210SBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Roots referenced by my leaves, by rank\n",rank);CHKERRQ(ierr);
78981bfa7aaSJed Brown         for (ii=0; ii<sf->nranks; ii++) {
79081bfa7aaSJed Brown           i = perm[ii];
7917904a332SBarry Smith           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %d: %D edges\n",rank,sf->ranks[i],sf->roffset[i+1]-sf->roffset[i]);CHKERRQ(ierr);
79295fce210SBarry Smith           for (j=sf->roffset[i]; j<sf->roffset[i+1]; j++) {
79395fce210SBarry Smith             ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d]    %D <- %D\n",rank,sf->rmine[j],sf->rremote[j]);CHKERRQ(ierr);
79495fce210SBarry Smith           }
79595fce210SBarry Smith         }
79681bfa7aaSJed Brown         ierr = PetscFree2(tmpranks,perm);CHKERRQ(ierr);
79795fce210SBarry Smith       }
79895fce210SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
7991575c14dSBarry Smith       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
800dd5b3ca6SJunchao Zhang     }
80195fce210SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
80295fce210SBarry Smith   }
80395fce210SBarry Smith   PetscFunctionReturn(0);
80495fce210SBarry Smith }
80595fce210SBarry Smith 
80695fce210SBarry Smith /*@C
807dec1416fSJunchao Zhang    PetscSFGetRootRanks - Get root ranks and number of vertices referenced by leaves on this process
80895fce210SBarry Smith 
80995fce210SBarry Smith    Not Collective
81095fce210SBarry Smith 
81195fce210SBarry Smith    Input Arguments:
81295fce210SBarry Smith .  sf - star forest
81395fce210SBarry Smith 
81495fce210SBarry Smith    Output Arguments:
81595fce210SBarry Smith +  nranks - number of ranks referenced by local part
81695fce210SBarry Smith .  ranks - array of ranks
81795fce210SBarry Smith .  roffset - offset in rmine/rremote for each rank (length nranks+1)
81895fce210SBarry Smith .  rmine - concatenated array holding local indices referencing each remote rank
81995fce210SBarry Smith -  rremote - concatenated array holding remote indices referenced for each remote rank
82095fce210SBarry Smith 
82195fce210SBarry Smith    Level: developer
82295fce210SBarry Smith 
823dec1416fSJunchao Zhang .seealso: PetscSFGetLeafRanks()
82495fce210SBarry Smith @*/
825dec1416fSJunchao Zhang PetscErrorCode PetscSFGetRootRanks(PetscSF sf,PetscInt *nranks,const PetscMPIInt **ranks,const PetscInt **roffset,const PetscInt **rmine,const PetscInt **rremote)
82695fce210SBarry Smith {
827dec1416fSJunchao Zhang   PetscErrorCode ierr;
82895fce210SBarry Smith 
82995fce210SBarry Smith   PetscFunctionBegin;
83095fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
83129046d53SLisandro Dalcin   if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining ranks");
832dec1416fSJunchao Zhang   if (sf->ops->GetRootRanks) {
833dec1416fSJunchao Zhang     ierr = (sf->ops->GetRootRanks)(sf,nranks,ranks,roffset,rmine,rremote);CHKERRQ(ierr);
834dec1416fSJunchao Zhang   } else {
835dec1416fSJunchao Zhang     /* The generic implementation */
83695fce210SBarry Smith     if (nranks)  *nranks  = sf->nranks;
83795fce210SBarry Smith     if (ranks)   *ranks   = sf->ranks;
83895fce210SBarry Smith     if (roffset) *roffset = sf->roffset;
83995fce210SBarry Smith     if (rmine)   *rmine   = sf->rmine;
84095fce210SBarry Smith     if (rremote) *rremote = sf->rremote;
841dec1416fSJunchao Zhang   }
84295fce210SBarry Smith   PetscFunctionReturn(0);
84395fce210SBarry Smith }
84495fce210SBarry Smith 
8458750ddebSJunchao Zhang /*@C
8468750ddebSJunchao Zhang    PetscSFGetLeafRanks - Get leaf ranks referencing roots on this process
8478750ddebSJunchao Zhang 
8488750ddebSJunchao Zhang    Not Collective
8498750ddebSJunchao Zhang 
8508750ddebSJunchao Zhang    Input Arguments:
8518750ddebSJunchao Zhang .  sf - star forest
8528750ddebSJunchao Zhang 
8538750ddebSJunchao Zhang    Output Arguments:
8548750ddebSJunchao Zhang +  niranks - number of leaf ranks referencing roots on this process
8558750ddebSJunchao Zhang .  iranks - array of ranks
8568750ddebSJunchao Zhang .  ioffset - offset in irootloc for each rank (length niranks+1)
8578750ddebSJunchao Zhang -  irootloc - concatenated array holding local indices of roots referenced by each leaf rank
8588750ddebSJunchao Zhang 
8598750ddebSJunchao Zhang    Level: developer
8608750ddebSJunchao Zhang 
861dec1416fSJunchao Zhang .seealso: PetscSFGetRootRanks()
8628750ddebSJunchao Zhang @*/
8638750ddebSJunchao Zhang PetscErrorCode PetscSFGetLeafRanks(PetscSF sf,PetscInt *niranks,const PetscMPIInt **iranks,const PetscInt **ioffset,const PetscInt **irootloc)
8648750ddebSJunchao Zhang {
8658750ddebSJunchao Zhang   PetscErrorCode ierr;
8668750ddebSJunchao Zhang 
8678750ddebSJunchao Zhang   PetscFunctionBegin;
8688750ddebSJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
8698750ddebSJunchao Zhang   if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining ranks");
8708750ddebSJunchao Zhang   if (sf->ops->GetLeafRanks) {
8718750ddebSJunchao Zhang     ierr = (sf->ops->GetLeafRanks)(sf,niranks,iranks,ioffset,irootloc);CHKERRQ(ierr);
8728750ddebSJunchao Zhang   } else {
8738750ddebSJunchao Zhang     PetscSFType type;
8748750ddebSJunchao Zhang     ierr = PetscSFGetType(sf,&type);CHKERRQ(ierr);
8758750ddebSJunchao Zhang     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSFGetLeafRanks() is not supported on this StarForest type: %s", type);
8768750ddebSJunchao Zhang   }
8778750ddebSJunchao Zhang   PetscFunctionReturn(0);
8788750ddebSJunchao Zhang }
8798750ddebSJunchao Zhang 
880b5a8e515SJed Brown static PetscBool InList(PetscMPIInt needle,PetscMPIInt n,const PetscMPIInt *list) {
881b5a8e515SJed Brown   PetscInt i;
882b5a8e515SJed Brown   for (i=0; i<n; i++) {
883b5a8e515SJed Brown     if (needle == list[i]) return PETSC_TRUE;
884b5a8e515SJed Brown   }
885b5a8e515SJed Brown   return PETSC_FALSE;
886b5a8e515SJed Brown }
887b5a8e515SJed Brown 
88895fce210SBarry Smith /*@C
88921c688dcSJed Brown    PetscSFSetUpRanks - Set up data structures associated with ranks; this is for internal use by PetscSF implementations.
89021c688dcSJed Brown 
89121c688dcSJed Brown    Collective
89221c688dcSJed Brown 
89321c688dcSJed Brown    Input Arguments:
894b5a8e515SJed Brown +  sf - PetscSF to set up; PetscSFSetGraph() must have been called
895b5a8e515SJed Brown -  dgroup - MPI_Group of ranks to be distinguished (e.g., for self or shared memory exchange)
89621c688dcSJed Brown 
89721c688dcSJed Brown    Level: developer
89821c688dcSJed Brown 
899dec1416fSJunchao Zhang .seealso: PetscSFGetRootRanks()
90021c688dcSJed Brown @*/
901b5a8e515SJed Brown PetscErrorCode PetscSFSetUpRanks(PetscSF sf,MPI_Group dgroup)
90221c688dcSJed Brown {
90321c688dcSJed Brown   PetscErrorCode     ierr;
90421c688dcSJed Brown   PetscTable         table;
90521c688dcSJed Brown   PetscTablePosition pos;
906b5a8e515SJed Brown   PetscMPIInt        size,groupsize,*groupranks;
907247e8311SStefano Zampini   PetscInt           *rcount,*ranks;
908247e8311SStefano Zampini   PetscInt           i, irank = -1,orank = -1;
90921c688dcSJed Brown 
91021c688dcSJed Brown   PetscFunctionBegin;
91121c688dcSJed Brown   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
91229046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
91321c688dcSJed Brown   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr);
91421c688dcSJed Brown   ierr = PetscTableCreate(10,size,&table);CHKERRQ(ierr);
91521c688dcSJed Brown   for (i=0; i<sf->nleaves; i++) {
91621c688dcSJed Brown     /* Log 1-based rank */
91721c688dcSJed Brown     ierr = PetscTableAdd(table,sf->remote[i].rank+1,1,ADD_VALUES);CHKERRQ(ierr);
91821c688dcSJed Brown   }
91921c688dcSJed Brown   ierr = PetscTableGetCount(table,&sf->nranks);CHKERRQ(ierr);
92021c688dcSJed Brown   ierr = PetscMalloc4(sf->nranks,&sf->ranks,sf->nranks+1,&sf->roffset,sf->nleaves,&sf->rmine,sf->nleaves,&sf->rremote);CHKERRQ(ierr);
92121c688dcSJed Brown   ierr = PetscMalloc2(sf->nranks,&rcount,sf->nranks,&ranks);CHKERRQ(ierr);
92221c688dcSJed Brown   ierr = PetscTableGetHeadPosition(table,&pos);CHKERRQ(ierr);
92321c688dcSJed Brown   for (i=0; i<sf->nranks; i++) {
92421c688dcSJed Brown     ierr = PetscTableGetNext(table,&pos,&ranks[i],&rcount[i]);CHKERRQ(ierr);
92521c688dcSJed Brown     ranks[i]--;             /* Convert back to 0-based */
92621c688dcSJed Brown   }
92721c688dcSJed Brown   ierr = PetscTableDestroy(&table);CHKERRQ(ierr);
928b5a8e515SJed Brown 
929b5a8e515SJed Brown   /* We expect that dgroup is reliably "small" while nranks could be large */
930b5a8e515SJed Brown   {
9317fb8a5e4SKarl Rupp     MPI_Group group = MPI_GROUP_NULL;
932b5a8e515SJed Brown     PetscMPIInt *dgroupranks;
933b5a8e515SJed Brown     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
934b5a8e515SJed Brown     ierr = MPI_Group_size(dgroup,&groupsize);CHKERRQ(ierr);
935b5a8e515SJed Brown     ierr = PetscMalloc1(groupsize,&dgroupranks);CHKERRQ(ierr);
936b5a8e515SJed Brown     ierr = PetscMalloc1(groupsize,&groupranks);CHKERRQ(ierr);
937b5a8e515SJed Brown     for (i=0; i<groupsize; i++) dgroupranks[i] = i;
938f3fc9a17SJed Brown     if (groupsize) {ierr = MPI_Group_translate_ranks(dgroup,groupsize,dgroupranks,group,groupranks);CHKERRQ(ierr);}
939b5a8e515SJed Brown     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
940b5a8e515SJed Brown     ierr = PetscFree(dgroupranks);CHKERRQ(ierr);
941b5a8e515SJed Brown   }
942b5a8e515SJed Brown 
943b5a8e515SJed Brown   /* Partition ranks[] into distinguished (first sf->ndranks) followed by non-distinguished */
944b5a8e515SJed Brown   for (sf->ndranks=0,i=sf->nranks; sf->ndranks<i; ) {
945b5a8e515SJed Brown     for (i--; sf->ndranks<i; i--) { /* Scan i backward looking for distinguished rank */
946b5a8e515SJed Brown       if (InList(ranks[i],groupsize,groupranks)) break;
947b5a8e515SJed Brown     }
948b5a8e515SJed Brown     for ( ; sf->ndranks<=i; sf->ndranks++) { /* Scan sf->ndranks forward looking for non-distinguished rank */
949b5a8e515SJed Brown       if (!InList(ranks[sf->ndranks],groupsize,groupranks)) break;
950b5a8e515SJed Brown     }
951b5a8e515SJed Brown     if (sf->ndranks < i) {                         /* Swap ranks[sf->ndranks] with ranks[i] */
952b5a8e515SJed Brown       PetscInt tmprank,tmpcount;
953247e8311SStefano Zampini 
954b5a8e515SJed Brown       tmprank             = ranks[i];
955b5a8e515SJed Brown       tmpcount            = rcount[i];
956b5a8e515SJed Brown       ranks[i]            = ranks[sf->ndranks];
957b5a8e515SJed Brown       rcount[i]           = rcount[sf->ndranks];
958b5a8e515SJed Brown       ranks[sf->ndranks]  = tmprank;
959b5a8e515SJed Brown       rcount[sf->ndranks] = tmpcount;
960b5a8e515SJed Brown       sf->ndranks++;
961b5a8e515SJed Brown     }
962b5a8e515SJed Brown   }
963b5a8e515SJed Brown   ierr = PetscFree(groupranks);CHKERRQ(ierr);
964b5a8e515SJed Brown   ierr = PetscSortIntWithArray(sf->ndranks,ranks,rcount);CHKERRQ(ierr);
965b5a8e515SJed Brown   ierr = PetscSortIntWithArray(sf->nranks-sf->ndranks,ranks+sf->ndranks,rcount+sf->ndranks);CHKERRQ(ierr);
96621c688dcSJed Brown   sf->roffset[0] = 0;
96721c688dcSJed Brown   for (i=0; i<sf->nranks; i++) {
96821c688dcSJed Brown     ierr = PetscMPIIntCast(ranks[i],sf->ranks+i);CHKERRQ(ierr);
96921c688dcSJed Brown     sf->roffset[i+1] = sf->roffset[i] + rcount[i];
97021c688dcSJed Brown     rcount[i]        = 0;
97121c688dcSJed Brown   }
972247e8311SStefano Zampini   for (i=0, irank = -1, orank = -1; i<sf->nleaves; i++) {
973247e8311SStefano Zampini     /* short circuit */
974247e8311SStefano Zampini     if (orank != sf->remote[i].rank) {
97521c688dcSJed Brown       /* Search for index of iremote[i].rank in sf->ranks */
976b5a8e515SJed Brown       ierr = PetscFindMPIInt(sf->remote[i].rank,sf->ndranks,sf->ranks,&irank);CHKERRQ(ierr);
977b5a8e515SJed Brown       if (irank < 0) {
978b5a8e515SJed Brown         ierr = PetscFindMPIInt(sf->remote[i].rank,sf->nranks-sf->ndranks,sf->ranks+sf->ndranks,&irank);CHKERRQ(ierr);
979b5a8e515SJed Brown         if (irank >= 0) irank += sf->ndranks;
98021c688dcSJed Brown       }
981247e8311SStefano Zampini       orank = sf->remote[i].rank;
982247e8311SStefano Zampini     }
983b5a8e515SJed Brown     if (irank < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Could not find rank %D in array",sf->remote[i].rank);
98421c688dcSJed Brown     sf->rmine[sf->roffset[irank] + rcount[irank]]   = sf->mine ? sf->mine[i] : i;
98521c688dcSJed Brown     sf->rremote[sf->roffset[irank] + rcount[irank]] = sf->remote[i].index;
98621c688dcSJed Brown     rcount[irank]++;
98721c688dcSJed Brown   }
98821c688dcSJed Brown   ierr = PetscFree2(rcount,ranks);CHKERRQ(ierr);
98921c688dcSJed Brown   PetscFunctionReturn(0);
99021c688dcSJed Brown }
99121c688dcSJed Brown 
99221c688dcSJed Brown /*@C
99395fce210SBarry Smith    PetscSFGetGroups - gets incoming and outgoing process groups
99495fce210SBarry Smith 
99595fce210SBarry Smith    Collective
99695fce210SBarry Smith 
99795fce210SBarry Smith    Input Argument:
99895fce210SBarry Smith .  sf - star forest
99995fce210SBarry Smith 
100095fce210SBarry Smith    Output Arguments:
100195fce210SBarry Smith +  incoming - group of origin processes for incoming edges (leaves that reference my roots)
100295fce210SBarry Smith -  outgoing - group of destination processes for outgoing edges (roots that I reference)
100395fce210SBarry Smith 
100495fce210SBarry Smith    Level: developer
100595fce210SBarry Smith 
100695fce210SBarry Smith .seealso: PetscSFGetWindow(), PetscSFRestoreWindow()
100795fce210SBarry Smith @*/
100895fce210SBarry Smith PetscErrorCode PetscSFGetGroups(PetscSF sf,MPI_Group *incoming,MPI_Group *outgoing)
100995fce210SBarry Smith {
101095fce210SBarry Smith   PetscErrorCode ierr;
10117fb8a5e4SKarl Rupp   MPI_Group      group = MPI_GROUP_NULL;
101295fce210SBarry Smith 
101395fce210SBarry Smith   PetscFunctionBegin;
101444ee17edSStefano Zampini   if (sf->nranks < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUpRanks() before obtaining groups");
101595fce210SBarry Smith   if (sf->ingroup == MPI_GROUP_NULL) {
101695fce210SBarry Smith     PetscInt       i;
101795fce210SBarry Smith     const PetscInt *indegree;
101895fce210SBarry Smith     PetscMPIInt    rank,*outranks,*inranks;
101995fce210SBarry Smith     PetscSFNode    *remote;
102095fce210SBarry Smith     PetscSF        bgcount;
102195fce210SBarry Smith 
102295fce210SBarry Smith     /* Compute the number of incoming ranks */
1023785e854fSJed Brown     ierr = PetscMalloc1(sf->nranks,&remote);CHKERRQ(ierr);
102495fce210SBarry Smith     for (i=0; i<sf->nranks; i++) {
102595fce210SBarry Smith       remote[i].rank  = sf->ranks[i];
102695fce210SBarry Smith       remote[i].index = 0;
102795fce210SBarry Smith     }
102895fce210SBarry Smith     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,&bgcount);CHKERRQ(ierr);
102995fce210SBarry Smith     ierr = PetscSFSetGraph(bgcount,1,sf->nranks,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
103095fce210SBarry Smith     ierr = PetscSFComputeDegreeBegin(bgcount,&indegree);CHKERRQ(ierr);
103195fce210SBarry Smith     ierr = PetscSFComputeDegreeEnd(bgcount,&indegree);CHKERRQ(ierr);
103295fce210SBarry Smith     /* Enumerate the incoming ranks */
1033dcca6d9dSJed Brown     ierr = PetscMalloc2(indegree[0],&inranks,sf->nranks,&outranks);CHKERRQ(ierr);
103495fce210SBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
103595fce210SBarry Smith     for (i=0; i<sf->nranks; i++) outranks[i] = rank;
103695fce210SBarry Smith     ierr = PetscSFGatherBegin(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr);
103795fce210SBarry Smith     ierr = PetscSFGatherEnd(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr);
103895fce210SBarry Smith     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
103995fce210SBarry Smith     ierr = MPI_Group_incl(group,indegree[0],inranks,&sf->ingroup);CHKERRQ(ierr);
104095fce210SBarry Smith     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
104195fce210SBarry Smith     ierr = PetscFree2(inranks,outranks);CHKERRQ(ierr);
104295fce210SBarry Smith     ierr = PetscSFDestroy(&bgcount);CHKERRQ(ierr);
104395fce210SBarry Smith   }
104495fce210SBarry Smith   *incoming = sf->ingroup;
104595fce210SBarry Smith 
104695fce210SBarry Smith   if (sf->outgroup == MPI_GROUP_NULL) {
104795fce210SBarry Smith     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
104895fce210SBarry Smith     ierr = MPI_Group_incl(group,sf->nranks,sf->ranks,&sf->outgroup);CHKERRQ(ierr);
104995fce210SBarry Smith     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
105095fce210SBarry Smith   }
105195fce210SBarry Smith   *outgoing = sf->outgroup;
105295fce210SBarry Smith   PetscFunctionReturn(0);
105395fce210SBarry Smith }
105495fce210SBarry Smith 
105529046d53SLisandro Dalcin /*@
105695fce210SBarry Smith    PetscSFGetMultiSF - gets the inner SF implemeting gathers and scatters
105795fce210SBarry Smith 
105895fce210SBarry Smith    Collective
105995fce210SBarry Smith 
106095fce210SBarry Smith    Input Argument:
106195fce210SBarry Smith .  sf - star forest that may contain roots with 0 or with more than 1 vertex
106295fce210SBarry Smith 
106395fce210SBarry Smith    Output Arguments:
106495fce210SBarry Smith .  multi - star forest with split roots, such that each root has degree exactly 1
106595fce210SBarry Smith 
106695fce210SBarry Smith    Level: developer
106795fce210SBarry Smith 
106895fce210SBarry Smith    Notes:
106995fce210SBarry Smith 
107095fce210SBarry Smith    In most cases, users should use PetscSFGatherBegin() and PetscSFScatterBegin() instead of manipulating multi
107195fce210SBarry Smith    directly. Since multi satisfies the stronger condition that each entry in the global space has exactly one incoming
107295fce210SBarry Smith    edge, it is a candidate for future optimization that might involve its removal.
107395fce210SBarry Smith 
1074673100f5SVaclav Hapla .seealso: PetscSFSetGraph(), PetscSFGatherBegin(), PetscSFScatterBegin(), PetscSFComputeMultiRootOriginalNumbering()
107595fce210SBarry Smith @*/
107695fce210SBarry Smith PetscErrorCode PetscSFGetMultiSF(PetscSF sf,PetscSF *multi)
107795fce210SBarry Smith {
107895fce210SBarry Smith   PetscErrorCode ierr;
107995fce210SBarry Smith 
108095fce210SBarry Smith   PetscFunctionBegin;
108195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
108295fce210SBarry Smith   PetscValidPointer(multi,2);
108395fce210SBarry Smith   if (sf->nroots < 0) {         /* Graph has not been set yet; why do we need this? */
108495fce210SBarry Smith     ierr   = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr);
108595fce210SBarry Smith     *multi = sf->multi;
108695fce210SBarry Smith     PetscFunctionReturn(0);
108795fce210SBarry Smith   }
108895fce210SBarry Smith   if (!sf->multi) {
108995fce210SBarry Smith     const PetscInt *indegree;
10909837ea96SMatthew G. Knepley     PetscInt       i,*inoffset,*outones,*outoffset,maxlocal;
109195fce210SBarry Smith     PetscSFNode    *remote;
109229046d53SLisandro Dalcin     maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */
109395fce210SBarry Smith     ierr = PetscSFComputeDegreeBegin(sf,&indegree);CHKERRQ(ierr);
109495fce210SBarry Smith     ierr = PetscSFComputeDegreeEnd(sf,&indegree);CHKERRQ(ierr);
10959837ea96SMatthew G. Knepley     ierr = PetscMalloc3(sf->nroots+1,&inoffset,maxlocal,&outones,maxlocal,&outoffset);CHKERRQ(ierr);
109695fce210SBarry Smith     inoffset[0] = 0;
109795fce210SBarry Smith     for (i=0; i<sf->nroots; i++) inoffset[i+1] = inoffset[i] + indegree[i];
10989837ea96SMatthew G. Knepley     for (i=0; i<maxlocal; i++) outones[i] = 1;
1099dbd2ff41SBarry Smith     ierr = PetscSFFetchAndOpBegin(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr);
1100dbd2ff41SBarry Smith     ierr = PetscSFFetchAndOpEnd(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr);
110195fce210SBarry Smith     for (i=0; i<sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */
110295fce210SBarry Smith #if defined(PETSC_USE_DEBUG)                                 /* Check that the expected number of increments occurred */
110395fce210SBarry Smith     for (i=0; i<sf->nroots; i++) {
110495fce210SBarry Smith       if (inoffset[i] + indegree[i] != inoffset[i+1]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect result after PetscSFFetchAndOp");
110595fce210SBarry Smith     }
110695fce210SBarry Smith #endif
1107785e854fSJed Brown     ierr = PetscMalloc1(sf->nleaves,&remote);CHKERRQ(ierr);
110895fce210SBarry Smith     for (i=0; i<sf->nleaves; i++) {
110995fce210SBarry Smith       remote[i].rank  = sf->remote[i].rank;
111038e7336fSToby Isaac       remote[i].index = outoffset[sf->mine ? sf->mine[i] : i];
111195fce210SBarry Smith     }
111295fce210SBarry Smith     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr);
111301365b40SToby Isaac     ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
111495fce210SBarry Smith     if (sf->rankorder) {        /* Sort the ranks */
111595fce210SBarry Smith       PetscMPIInt rank;
111695fce210SBarry Smith       PetscInt    *inranks,*newoffset,*outranks,*newoutoffset,*tmpoffset,maxdegree;
111795fce210SBarry Smith       PetscSFNode *newremote;
111895fce210SBarry Smith       ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
111995fce210SBarry Smith       for (i=0,maxdegree=0; i<sf->nroots; i++) maxdegree = PetscMax(maxdegree,indegree[i]);
11209837ea96SMatthew G. Knepley       ierr = PetscMalloc5(sf->multi->nroots,&inranks,sf->multi->nroots,&newoffset,maxlocal,&outranks,maxlocal,&newoutoffset,maxdegree,&tmpoffset);CHKERRQ(ierr);
11219837ea96SMatthew G. Knepley       for (i=0; i<maxlocal; i++) outranks[i] = rank;
11228bfbc91cSJed Brown       ierr = PetscSFReduceBegin(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr);
11238bfbc91cSJed Brown       ierr = PetscSFReduceEnd(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr);
112495fce210SBarry Smith       /* Sort the incoming ranks at each vertex, build the inverse map */
112595fce210SBarry Smith       for (i=0; i<sf->nroots; i++) {
112695fce210SBarry Smith         PetscInt j;
112795fce210SBarry Smith         for (j=0; j<indegree[i]; j++) tmpoffset[j] = j;
112895fce210SBarry Smith         ierr = PetscSortIntWithArray(indegree[i],inranks+inoffset[i],tmpoffset);CHKERRQ(ierr);
112995fce210SBarry Smith         for (j=0; j<indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j;
113095fce210SBarry Smith       }
113195fce210SBarry Smith       ierr = PetscSFBcastBegin(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr);
113295fce210SBarry Smith       ierr = PetscSFBcastEnd(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr);
1133785e854fSJed Brown       ierr = PetscMalloc1(sf->nleaves,&newremote);CHKERRQ(ierr);
113495fce210SBarry Smith       for (i=0; i<sf->nleaves; i++) {
113595fce210SBarry Smith         newremote[i].rank  = sf->remote[i].rank;
113601365b40SToby Isaac         newremote[i].index = newoutoffset[sf->mine ? sf->mine[i] : i];
113795fce210SBarry Smith       }
113801365b40SToby Isaac       ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,newremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
113995fce210SBarry Smith       ierr = PetscFree5(inranks,newoffset,outranks,newoutoffset,tmpoffset);CHKERRQ(ierr);
114095fce210SBarry Smith     }
114195fce210SBarry Smith     ierr = PetscFree3(inoffset,outones,outoffset);CHKERRQ(ierr);
114295fce210SBarry Smith   }
114395fce210SBarry Smith   *multi = sf->multi;
114495fce210SBarry Smith   PetscFunctionReturn(0);
114595fce210SBarry Smith }
114695fce210SBarry Smith 
114795fce210SBarry Smith /*@C
114895fce210SBarry Smith    PetscSFCreateEmbeddedSF - removes edges from all but the selected roots, does not remap indices
114995fce210SBarry Smith 
115095fce210SBarry Smith    Collective
115195fce210SBarry Smith 
115295fce210SBarry Smith    Input Arguments:
115395fce210SBarry Smith +  sf - original star forest
1154ba2a7774SJunchao Zhang .  nselected  - number of selected roots on this process
1155ba2a7774SJunchao Zhang -  selected   - indices of the selected roots on this process
115695fce210SBarry Smith 
115795fce210SBarry Smith    Output Arguments:
115895fce210SBarry Smith .  newsf - new star forest
115995fce210SBarry Smith 
116095fce210SBarry Smith    Level: advanced
116195fce210SBarry Smith 
116295fce210SBarry Smith    Note:
116395fce210SBarry Smith    To use the new PetscSF, it may be necessary to know the indices of the leaves that are still participating. This can
116495fce210SBarry Smith    be done by calling PetscSFGetGraph().
116595fce210SBarry Smith 
116695fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFGetGraph()
116795fce210SBarry Smith @*/
1168ba2a7774SJunchao Zhang PetscErrorCode PetscSFCreateEmbeddedSF(PetscSF sf,PetscInt nselected,const PetscInt *selected,PetscSF *newsf)
116995fce210SBarry Smith {
1170f659e5c7SJunchao Zhang   PetscInt          i,n,*roots,*rootdata,*leafdata,nroots,nleaves,connected_leaves,*new_ilocal;
1171ba2a7774SJunchao Zhang   const PetscSFNode *iremote;
1172f659e5c7SJunchao Zhang   PetscSFNode       *new_iremote;
1173ba2a7774SJunchao Zhang   PetscSF           tmpsf;
1174f659e5c7SJunchao Zhang   MPI_Comm          comm;
11750511a646SMatthew G. Knepley   PetscErrorCode    ierr;
117695fce210SBarry Smith 
117795fce210SBarry Smith   PetscFunctionBegin;
117895fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
117929046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
1180ba2a7774SJunchao Zhang   if (nselected) PetscValidPointer(selected,3);
118195fce210SBarry Smith   PetscValidPointer(newsf,4);
11820511a646SMatthew G. Knepley 
1183f659e5c7SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1184140a1472SStefano Zampini   ierr = PetscLogEventBegin(PETSCSF_EmbedSF,sf,0,0,0);CHKERRQ(ierr);
1185f659e5c7SJunchao Zhang   /* Uniq selected[] and put results in roots[] */
1186f659e5c7SJunchao Zhang   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
1187f659e5c7SJunchao Zhang   ierr = PetscMalloc1(nselected,&roots);CHKERRQ(ierr);
1188dd5b3ca6SJunchao Zhang   ierr = PetscArraycpy(roots,selected,nselected);CHKERRQ(ierr);
1189f659e5c7SJunchao Zhang   ierr = PetscSortedRemoveDupsInt(&nselected,roots);CHKERRQ(ierr);
1190f659e5c7SJunchao Zhang   if (nselected && (roots[0] < 0 || roots[nselected-1] >= sf->nroots)) SETERRQ3(comm,PETSC_ERR_ARG_OUTOFRANGE,"Min/Max root indices %D/%D are not in [0,%D)",roots[0],roots[nselected-1],sf->nroots);
1191f659e5c7SJunchao Zhang 
1192f659e5c7SJunchao Zhang   if (sf->ops->CreateEmbeddedSF) {
1193f659e5c7SJunchao Zhang     ierr = (*sf->ops->CreateEmbeddedSF)(sf,nselected,roots,newsf);CHKERRQ(ierr);
1194f659e5c7SJunchao Zhang   } else {
1195f659e5c7SJunchao Zhang     /* A generic version of creating embedded sf. Note that we called PetscSFSetGraph() twice, which is certainly expensive */
1196ba2a7774SJunchao Zhang     /* Find out which leaves (not leaf data items) are still connected to roots in the embedded sf */
1197f659e5c7SJunchao Zhang     ierr = PetscSFGetGraph(sf,&nroots,&nleaves,NULL,&iremote);CHKERRQ(ierr);
1198ba2a7774SJunchao Zhang     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&tmpsf);CHKERRQ(ierr);
1199ba2a7774SJunchao Zhang     ierr = PetscSFSetGraph(tmpsf,nroots,nleaves,NULL/*contiguous*/,PETSC_USE_POINTER,iremote,PETSC_USE_POINTER);CHKERRQ(ierr);
1200ba2a7774SJunchao Zhang     ierr = PetscCalloc2(nroots,&rootdata,nleaves,&leafdata);CHKERRQ(ierr);
1201f659e5c7SJunchao Zhang     for (i=0; i<nselected; ++i) rootdata[roots[i]] = 1;
1202ba2a7774SJunchao Zhang     ierr = PetscSFBcastBegin(tmpsf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr);
1203ba2a7774SJunchao Zhang     ierr = PetscSFBcastEnd(tmpsf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr);
1204ba2a7774SJunchao Zhang     ierr = PetscSFDestroy(&tmpsf);CHKERRQ(ierr);
1205ba2a7774SJunchao Zhang 
1206ba2a7774SJunchao Zhang     /* Build newsf with leaves that are still connected */
1207f659e5c7SJunchao Zhang     connected_leaves = 0;
1208f659e5c7SJunchao Zhang     for (i=0; i<nleaves; ++i) connected_leaves += leafdata[i];
1209f659e5c7SJunchao Zhang     ierr = PetscMalloc1(connected_leaves,&new_ilocal);CHKERRQ(ierr);
1210f659e5c7SJunchao Zhang     ierr = PetscMalloc1(connected_leaves,&new_iremote);CHKERRQ(ierr);
1211ba2a7774SJunchao Zhang     for (i=0, n=0; i<nleaves; ++i) {
1212ba2a7774SJunchao Zhang       if (leafdata[i]) {
1213ba2a7774SJunchao Zhang         new_ilocal[n]        = sf->mine ? sf->mine[i] : i;
1214ba2a7774SJunchao Zhang         new_iremote[n].rank  = sf->remote[i].rank;
1215ba2a7774SJunchao Zhang         new_iremote[n].index = sf->remote[i].index;
1216fc1ede2bSMatthew G. Knepley         ++n;
121795fce210SBarry Smith       }
121895fce210SBarry Smith     }
1219f659e5c7SJunchao Zhang 
1220f659e5c7SJunchao Zhang     if (n != connected_leaves) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"There is a size mismatch in the SF embedding, %d != %d",n,connected_leaves);
12214cc19a0cSStefano Zampini     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,newsf);CHKERRQ(ierr);
1222f659e5c7SJunchao Zhang     ierr = PetscSFSetGraph(*newsf,nroots,connected_leaves,new_ilocal,PETSC_OWN_POINTER,new_iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
122395fce210SBarry Smith     ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr);
1224f659e5c7SJunchao Zhang   }
1225f659e5c7SJunchao Zhang   ierr = PetscFree(roots);CHKERRQ(ierr);
1226140a1472SStefano Zampini   ierr = PetscSFSetUp(*newsf);CHKERRQ(ierr);
1227140a1472SStefano Zampini   ierr = PetscLogEventEnd(PETSCSF_EmbedSF,sf,0,0,0);CHKERRQ(ierr);
122895fce210SBarry Smith   PetscFunctionReturn(0);
122995fce210SBarry Smith }
123095fce210SBarry Smith 
12312f5fb4c2SMatthew G. Knepley /*@C
12322f5fb4c2SMatthew G. Knepley   PetscSFCreateEmbeddedLeafSF - removes edges from all but the selected leaves, does not remap indices
12332f5fb4c2SMatthew G. Knepley 
12342f5fb4c2SMatthew G. Knepley   Collective
12352f5fb4c2SMatthew G. Knepley 
12362f5fb4c2SMatthew G. Knepley   Input Arguments:
12372f5fb4c2SMatthew G. Knepley + sf - original star forest
1238f659e5c7SJunchao Zhang . nselected  - number of selected leaves on this process
1239f659e5c7SJunchao Zhang - selected   - indices of the selected leaves on this process
12402f5fb4c2SMatthew G. Knepley 
12412f5fb4c2SMatthew G. Knepley   Output Arguments:
12422f5fb4c2SMatthew G. Knepley .  newsf - new star forest
12432f5fb4c2SMatthew G. Knepley 
12442f5fb4c2SMatthew G. Knepley   Level: advanced
12452f5fb4c2SMatthew G. Knepley 
12462f5fb4c2SMatthew G. Knepley .seealso: PetscSFCreateEmbeddedSF(), PetscSFSetGraph(), PetscSFGetGraph()
12472f5fb4c2SMatthew G. Knepley @*/
1248f659e5c7SJunchao Zhang PetscErrorCode PetscSFCreateEmbeddedLeafSF(PetscSF sf,PetscInt nselected,const PetscInt *selected,PetscSF *newsf)
12492f5fb4c2SMatthew G. Knepley {
1250f659e5c7SJunchao Zhang   const PetscSFNode *iremote;
1251f659e5c7SJunchao Zhang   PetscSFNode       *new_iremote;
1252f659e5c7SJunchao Zhang   const PetscInt    *ilocal;
1253f659e5c7SJunchao Zhang   PetscInt          i,nroots,*leaves,*new_ilocal;
1254f659e5c7SJunchao Zhang   MPI_Comm          comm;
12552f5fb4c2SMatthew G. Knepley   PetscErrorCode    ierr;
12562f5fb4c2SMatthew G. Knepley 
12572f5fb4c2SMatthew G. Knepley   PetscFunctionBegin;
12582f5fb4c2SMatthew G. Knepley   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
125929046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
1260f659e5c7SJunchao Zhang   if (nselected) PetscValidPointer(selected,3);
12612f5fb4c2SMatthew G. Knepley   PetscValidPointer(newsf,4);
12622f5fb4c2SMatthew G. Knepley 
1263f659e5c7SJunchao Zhang   /* Uniq selected[] and put results in leaves[] */
1264f659e5c7SJunchao Zhang   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
1265f659e5c7SJunchao Zhang   ierr = PetscMalloc1(nselected,&leaves);CHKERRQ(ierr);
1266dd5b3ca6SJunchao Zhang   ierr = PetscArraycpy(leaves,selected,nselected);CHKERRQ(ierr);
1267f659e5c7SJunchao Zhang   ierr = PetscSortedRemoveDupsInt(&nselected,leaves);CHKERRQ(ierr);
1268f659e5c7SJunchao Zhang   if (nselected && (leaves[0] < 0 || leaves[nselected-1] >= sf->nleaves)) SETERRQ3(comm,PETSC_ERR_ARG_OUTOFRANGE,"Min/Max leaf indices %D/%D are not in [0,%D)",leaves[0],leaves[nselected-1],sf->nleaves);
1269f659e5c7SJunchao Zhang 
1270f659e5c7SJunchao Zhang   /* Optimize the routine only when sf is setup and hence we can reuse sf's communication pattern */
1271f659e5c7SJunchao Zhang   if (sf->setupcalled && sf->ops->CreateEmbeddedLeafSF) {
1272f659e5c7SJunchao Zhang     ierr = (*sf->ops->CreateEmbeddedLeafSF)(sf,nselected,leaves,newsf);CHKERRQ(ierr);
1273f659e5c7SJunchao Zhang   } else {
1274f659e5c7SJunchao Zhang     ierr = PetscSFGetGraph(sf,&nroots,NULL,&ilocal,&iremote);CHKERRQ(ierr);
1275f659e5c7SJunchao Zhang     ierr = PetscMalloc1(nselected,&new_ilocal);CHKERRQ(ierr);
1276f659e5c7SJunchao Zhang     ierr = PetscMalloc1(nselected,&new_iremote);CHKERRQ(ierr);
1277f659e5c7SJunchao Zhang     for (i=0; i<nselected; ++i) {
1278f659e5c7SJunchao Zhang       const PetscInt l     = leaves[i];
1279f659e5c7SJunchao Zhang       new_ilocal[i]        = ilocal ? ilocal[l] : l;
1280f659e5c7SJunchao Zhang       new_iremote[i].rank  = iremote[l].rank;
1281f659e5c7SJunchao Zhang       new_iremote[i].index = iremote[l].index;
12822f5fb4c2SMatthew G. Knepley     }
12834cc19a0cSStefano Zampini     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,newsf);CHKERRQ(ierr);
1284f659e5c7SJunchao Zhang     ierr = PetscSFSetGraph(*newsf,nroots,nselected,new_ilocal,PETSC_OWN_POINTER,new_iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1285f659e5c7SJunchao Zhang   }
1286f659e5c7SJunchao Zhang   ierr = PetscFree(leaves);CHKERRQ(ierr);
12872f5fb4c2SMatthew G. Knepley   PetscFunctionReturn(0);
12882f5fb4c2SMatthew G. Knepley }
12892f5fb4c2SMatthew G. Knepley 
129095fce210SBarry Smith /*@C
12913482bfa8SJunchao Zhang    PetscSFBcastAndOpBegin - begin pointwise broadcast with root value being reduced to leaf value, to be concluded with call to PetscSFBcastAndOpEnd()
12923482bfa8SJunchao Zhang 
12933482bfa8SJunchao Zhang    Collective on PetscSF
12943482bfa8SJunchao Zhang 
12953482bfa8SJunchao Zhang    Input Arguments:
12963482bfa8SJunchao Zhang +  sf - star forest on which to communicate
12973482bfa8SJunchao Zhang .  unit - data type associated with each node
12983482bfa8SJunchao Zhang .  rootdata - buffer to broadcast
12993482bfa8SJunchao Zhang -  op - operation to use for reduction
13003482bfa8SJunchao Zhang 
13013482bfa8SJunchao Zhang    Output Arguments:
13023482bfa8SJunchao Zhang .  leafdata - buffer to be reduced with values from each leaf's respective root
13033482bfa8SJunchao Zhang 
13043482bfa8SJunchao Zhang    Level: intermediate
13053482bfa8SJunchao Zhang 
13063482bfa8SJunchao Zhang .seealso: PetscSFBcastAndOpEnd(), PetscSFBcastBegin(), PetscSFBcastEnd()
13073482bfa8SJunchao Zhang @*/
13083482bfa8SJunchao Zhang PetscErrorCode PetscSFBcastAndOpBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
13093482bfa8SJunchao Zhang {
13103482bfa8SJunchao Zhang   PetscErrorCode ierr;
1311eb02082bSJunchao Zhang   PetscMemType   rootmtype,leafmtype;
13123482bfa8SJunchao Zhang 
13133482bfa8SJunchao Zhang   PetscFunctionBegin;
13143482bfa8SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
13153482bfa8SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
13163482bfa8SJunchao Zhang   ierr = PetscLogEventBegin(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
1317eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
1318eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
1319eb02082bSJunchao Zhang #if defined(PETSC_HAVE_CUDA)
1320eb02082bSJunchao Zhang   /*  Shall we assume rootdata, leafdata are ready to use (instead of being computed by some asynchronous kernels)?
1321eb02082bSJunchao Zhang     To be similar to MPI, I'd like to have this assumption, since MPI does not have a concept of stream.
1322eb02082bSJunchao Zhang     But currently this assumption is not enforecd in Petsc. To be safe, I do synchronization here. Otherwise, if
1323eb02082bSJunchao Zhang     we do not sync now and call the Pack kernel directly on the default NULL stream (assume petsc objects are also
1324eb02082bSJunchao Zhang     computed on it), we have to sync the NULL stream before calling MPI routines. So, it looks a cudaDeviceSynchronize
1325eb02082bSJunchao Zhang     is inevitable. We do it now and put pack/unpack kernels to non-NULL streams.
1326eb02082bSJunchao Zhang    */
1327eb02082bSJunchao Zhang   if (rootmtype == PETSC_MEMTYPE_DEVICE || leafmtype == PETSC_MEMTYPE_DEVICE) {cudaError_t err = cudaDeviceSynchronize();CHKERRCUDA(err);}
1328eb02082bSJunchao Zhang #endif
1329eb02082bSJunchao Zhang   ierr = (*sf->ops->BcastAndOpBegin)(sf,unit,rootmtype,rootdata,leafmtype,leafdata,op);CHKERRQ(ierr);
13303482bfa8SJunchao Zhang   ierr = PetscLogEventEnd(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
13313482bfa8SJunchao Zhang   PetscFunctionReturn(0);
13323482bfa8SJunchao Zhang }
13333482bfa8SJunchao Zhang 
13343482bfa8SJunchao Zhang /*@C
13353482bfa8SJunchao Zhang    PetscSFBcastAndOpEnd - end a broadcast & reduce operation started with PetscSFBcastAndOpBegin()
13363482bfa8SJunchao Zhang 
13373482bfa8SJunchao Zhang    Collective
13383482bfa8SJunchao Zhang 
13393482bfa8SJunchao Zhang    Input Arguments:
13403482bfa8SJunchao Zhang +  sf - star forest
13413482bfa8SJunchao Zhang .  unit - data type
13423482bfa8SJunchao Zhang .  rootdata - buffer to broadcast
13433482bfa8SJunchao Zhang -  op - operation to use for reduction
13443482bfa8SJunchao Zhang 
13453482bfa8SJunchao Zhang    Output Arguments:
13463482bfa8SJunchao Zhang .  leafdata - buffer to be reduced with values from each leaf's respective root
13473482bfa8SJunchao Zhang 
13483482bfa8SJunchao Zhang    Level: intermediate
13493482bfa8SJunchao Zhang 
13503482bfa8SJunchao Zhang .seealso: PetscSFSetGraph(), PetscSFReduceEnd()
13513482bfa8SJunchao Zhang @*/
13523482bfa8SJunchao Zhang PetscErrorCode PetscSFBcastAndOpEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
13533482bfa8SJunchao Zhang {
13543482bfa8SJunchao Zhang   PetscErrorCode ierr;
1355eb02082bSJunchao Zhang   PetscMemType   rootmtype,leafmtype;
13563482bfa8SJunchao Zhang 
13573482bfa8SJunchao Zhang   PetscFunctionBegin;
13583482bfa8SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
13593482bfa8SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
13603482bfa8SJunchao Zhang   ierr = PetscLogEventBegin(PETSCSF_BcastAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
1361eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
1362eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
1363eb02082bSJunchao Zhang   ierr = (*sf->ops->BcastAndOpEnd)(sf,unit,rootmtype,rootdata,leafmtype,leafdata,op);CHKERRQ(ierr);
13643482bfa8SJunchao Zhang   ierr = PetscLogEventEnd(PETSCSF_BcastAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
13653482bfa8SJunchao Zhang   PetscFunctionReturn(0);
13663482bfa8SJunchao Zhang }
13673482bfa8SJunchao Zhang 
13683482bfa8SJunchao Zhang /*@C
136995fce210SBarry Smith    PetscSFReduceBegin - begin reduction of leafdata into rootdata, to be completed with call to PetscSFReduceEnd()
137095fce210SBarry Smith 
137195fce210SBarry Smith    Collective
137295fce210SBarry Smith 
137395fce210SBarry Smith    Input Arguments:
137495fce210SBarry Smith +  sf - star forest
137595fce210SBarry Smith .  unit - data type
137695fce210SBarry Smith .  leafdata - values to reduce
137795fce210SBarry Smith -  op - reduction operation
137895fce210SBarry Smith 
137995fce210SBarry Smith    Output Arguments:
138095fce210SBarry Smith .  rootdata - result of reduction of values from all leaves of each root
138195fce210SBarry Smith 
138295fce210SBarry Smith    Level: intermediate
138395fce210SBarry Smith 
138495fce210SBarry Smith .seealso: PetscSFBcastBegin()
138595fce210SBarry Smith @*/
138695fce210SBarry Smith PetscErrorCode PetscSFReduceBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
138795fce210SBarry Smith {
138895fce210SBarry Smith   PetscErrorCode ierr;
1389eb02082bSJunchao Zhang   PetscMemType   rootmtype,leafmtype;
139095fce210SBarry Smith 
139195fce210SBarry Smith   PetscFunctionBegin;
139295fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
139395fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
139429046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);
1395eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
1396eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
1397eb02082bSJunchao Zhang #if defined(PETSC_HAVE_CUDA)
1398eb02082bSJunchao Zhang   if (rootmtype == PETSC_MEMTYPE_DEVICE || leafmtype == PETSC_MEMTYPE_DEVICE) {cudaError_t err = cudaDeviceSynchronize();CHKERRCUDA(err);}
1399eb02082bSJunchao Zhang #endif
1400eb02082bSJunchao Zhang   ierr = (sf->ops->ReduceBegin)(sf,unit,leafmtype,leafdata,rootmtype,rootdata,op);CHKERRQ(ierr);
1401563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);
140295fce210SBarry Smith   PetscFunctionReturn(0);
140395fce210SBarry Smith }
140495fce210SBarry Smith 
140595fce210SBarry Smith /*@C
140695fce210SBarry Smith    PetscSFReduceEnd - end a reduction operation started with PetscSFReduceBegin()
140795fce210SBarry Smith 
140895fce210SBarry Smith    Collective
140995fce210SBarry Smith 
141095fce210SBarry Smith    Input Arguments:
141195fce210SBarry Smith +  sf - star forest
141295fce210SBarry Smith .  unit - data type
141395fce210SBarry Smith .  leafdata - values to reduce
141495fce210SBarry Smith -  op - reduction operation
141595fce210SBarry Smith 
141695fce210SBarry Smith    Output Arguments:
141795fce210SBarry Smith .  rootdata - result of reduction of values from all leaves of each root
141895fce210SBarry Smith 
141995fce210SBarry Smith    Level: intermediate
142095fce210SBarry Smith 
142195fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFBcastEnd()
142295fce210SBarry Smith @*/
142395fce210SBarry Smith PetscErrorCode PetscSFReduceEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
142495fce210SBarry Smith {
142595fce210SBarry Smith   PetscErrorCode ierr;
1426eb02082bSJunchao Zhang   PetscMemType   rootmtype,leafmtype;
142795fce210SBarry Smith 
142895fce210SBarry Smith   PetscFunctionBegin;
142995fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
143095fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
143129046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);
1432eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
1433eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
1434eb02082bSJunchao Zhang   ierr = (*sf->ops->ReduceEnd)(sf,unit,leafmtype,leafdata,rootmtype,rootdata,op);CHKERRQ(ierr);
1435563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);
143695fce210SBarry Smith   PetscFunctionReturn(0);
143795fce210SBarry Smith }
143895fce210SBarry Smith 
143995fce210SBarry Smith /*@C
1440a1729e3fSJunchao Zhang    PetscSFFetchAndOpBegin - begin operation that fetches values from root and updates atomically by applying operation using my leaf value, to be completed with PetscSFFetchAndOpEnd()
1441a1729e3fSJunchao Zhang 
1442a1729e3fSJunchao Zhang    Collective
1443a1729e3fSJunchao Zhang 
1444a1729e3fSJunchao Zhang    Input Arguments:
1445a1729e3fSJunchao Zhang +  sf - star forest
1446a1729e3fSJunchao Zhang .  unit - data type
1447a1729e3fSJunchao Zhang .  leafdata - leaf values to use in reduction
1448a1729e3fSJunchao Zhang -  op - operation to use for reduction
1449a1729e3fSJunchao Zhang 
1450a1729e3fSJunchao Zhang    Output Arguments:
1451a1729e3fSJunchao Zhang +  rootdata - root values to be updated, input state is seen by first process to perform an update
1452a1729e3fSJunchao Zhang -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
1453a1729e3fSJunchao Zhang 
1454a1729e3fSJunchao Zhang    Level: advanced
1455a1729e3fSJunchao Zhang 
1456a1729e3fSJunchao Zhang    Note:
1457a1729e3fSJunchao Zhang    The update is only atomic at the granularity provided by the hardware. Different roots referenced by the same process
1458a1729e3fSJunchao Zhang    might be updated in a different order. Furthermore, if a composite type is used for the unit datatype, atomicity is
1459a1729e3fSJunchao Zhang    not guaranteed across the whole vertex. Therefore, this function is mostly only used with primitive types such as
1460a1729e3fSJunchao Zhang    integers.
1461a1729e3fSJunchao Zhang 
1462a1729e3fSJunchao Zhang .seealso: PetscSFComputeDegreeBegin(), PetscSFReduceBegin(), PetscSFSetGraph()
1463a1729e3fSJunchao Zhang @*/
1464a1729e3fSJunchao Zhang PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1465a1729e3fSJunchao Zhang {
1466a1729e3fSJunchao Zhang   PetscErrorCode ierr;
1467eb02082bSJunchao Zhang   PetscMemType   rootmtype,leafmtype,leafupdatemtype;
1468a1729e3fSJunchao Zhang 
1469a1729e3fSJunchao Zhang   PetscFunctionBegin;
1470a1729e3fSJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1471a1729e3fSJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1472a1729e3fSJunchao Zhang   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
1473eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
1474eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
1475eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafupdate,&leafupdatemtype);CHKERRQ(ierr);
1476eb02082bSJunchao Zhang   if (leafmtype != leafupdatemtype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for leafdata and leafupdate in different memory types");
1477eb02082bSJunchao Zhang #if defined(PETSC_HAVE_CUDA)
1478eb02082bSJunchao Zhang   if (rootmtype == PETSC_MEMTYPE_DEVICE || leafmtype == PETSC_MEMTYPE_DEVICE) {cudaError_t err = cudaDeviceSynchronize();CHKERRCUDA(err);}
1479eb02082bSJunchao Zhang #endif
1480eb02082bSJunchao Zhang   ierr = (*sf->ops->FetchAndOpBegin)(sf,unit,rootmtype,rootdata,leafmtype,leafdata,leafupdate,op);CHKERRQ(ierr);
1481a1729e3fSJunchao Zhang   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
1482a1729e3fSJunchao Zhang   PetscFunctionReturn(0);
1483a1729e3fSJunchao Zhang }
1484a1729e3fSJunchao Zhang 
1485a1729e3fSJunchao Zhang /*@C
1486a1729e3fSJunchao Zhang    PetscSFFetchAndOpEnd - end operation started in matching call to PetscSFFetchAndOpBegin() to fetch values from roots and update atomically by applying operation using my leaf value
1487a1729e3fSJunchao Zhang 
1488a1729e3fSJunchao Zhang    Collective
1489a1729e3fSJunchao Zhang 
1490a1729e3fSJunchao Zhang    Input Arguments:
1491a1729e3fSJunchao Zhang +  sf - star forest
1492a1729e3fSJunchao Zhang .  unit - data type
1493a1729e3fSJunchao Zhang .  leafdata - leaf values to use in reduction
1494a1729e3fSJunchao Zhang -  op - operation to use for reduction
1495a1729e3fSJunchao Zhang 
1496a1729e3fSJunchao Zhang    Output Arguments:
1497a1729e3fSJunchao Zhang +  rootdata - root values to be updated, input state is seen by first process to perform an update
1498a1729e3fSJunchao Zhang -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
1499a1729e3fSJunchao Zhang 
1500a1729e3fSJunchao Zhang    Level: advanced
1501a1729e3fSJunchao Zhang 
1502a1729e3fSJunchao Zhang .seealso: PetscSFComputeDegreeEnd(), PetscSFReduceEnd(), PetscSFSetGraph()
1503a1729e3fSJunchao Zhang @*/
1504a1729e3fSJunchao Zhang PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1505a1729e3fSJunchao Zhang {
1506a1729e3fSJunchao Zhang   PetscErrorCode ierr;
1507eb02082bSJunchao Zhang   PetscMemType   rootmtype,leafmtype,leafupdatemtype;
1508a1729e3fSJunchao Zhang 
1509a1729e3fSJunchao Zhang   PetscFunctionBegin;
1510a1729e3fSJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1511a1729e3fSJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1512a1729e3fSJunchao Zhang   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
1513eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
1514eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
1515eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafupdate,&leafupdatemtype);CHKERRQ(ierr);
1516eb02082bSJunchao Zhang   if (leafmtype != leafupdatemtype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for leafdata and leafupdate in different memory types");
1517eb02082bSJunchao Zhang   ierr = (*sf->ops->FetchAndOpEnd)(sf,unit,rootmtype,rootdata,leafmtype,leafdata,leafupdate,op);CHKERRQ(ierr);
1518a1729e3fSJunchao Zhang   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
1519a1729e3fSJunchao Zhang   PetscFunctionReturn(0);
1520a1729e3fSJunchao Zhang }
1521a1729e3fSJunchao Zhang 
1522a1729e3fSJunchao Zhang /*@C
152395fce210SBarry Smith    PetscSFComputeDegreeBegin - begin computation of degree for each root vertex, to be completed with PetscSFComputeDegreeEnd()
152495fce210SBarry Smith 
152595fce210SBarry Smith    Collective
152695fce210SBarry Smith 
152795fce210SBarry Smith    Input Arguments:
152895fce210SBarry Smith .  sf - star forest
152995fce210SBarry Smith 
153095fce210SBarry Smith    Output Arguments:
153195fce210SBarry Smith .  degree - degree of each root vertex
153295fce210SBarry Smith 
153395fce210SBarry Smith    Level: advanced
153495fce210SBarry Smith 
1535ffe67aa5SVáclav Hapla    Notes:
1536ffe67aa5SVáclav Hapla    The returned array is owned by PetscSF and automatically freed by PetscSFDestroy(). Hence no need to call PetscFree() on it.
1537ffe67aa5SVáclav Hapla 
153895fce210SBarry Smith .seealso: PetscSFGatherBegin()
153995fce210SBarry Smith @*/
154095fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf,const PetscInt **degree)
154195fce210SBarry Smith {
154295fce210SBarry Smith   PetscErrorCode ierr;
154395fce210SBarry Smith 
154495fce210SBarry Smith   PetscFunctionBegin;
154595fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
154695fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
154795fce210SBarry Smith   PetscValidPointer(degree,2);
1548803bd9e8SMatthew G. Knepley   if (!sf->degreeknown) {
1549*5b0d146aSStefano Zampini     PetscInt i, nroots = sf->nroots, maxlocal;
1550803bd9e8SMatthew G. Knepley     if (sf->degree) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Calls to PetscSFComputeDegreeBegin() cannot be nested.");
1551*5b0d146aSStefano Zampini     maxlocal = sf->maxleaf-sf->minleaf+1;
155229046d53SLisandro Dalcin     ierr = PetscMalloc1(nroots,&sf->degree);CHKERRQ(ierr);
155329046d53SLisandro Dalcin     ierr = PetscMalloc1(PetscMax(maxlocal,1),&sf->degreetmp);CHKERRQ(ierr); /* allocate at least one entry, see check in PetscSFComputeDegreeEnd() */
155429046d53SLisandro Dalcin     for (i=0; i<nroots; i++) sf->degree[i] = 0;
15559837ea96SMatthew G. Knepley     for (i=0; i<maxlocal; i++) sf->degreetmp[i] = 1;
1556*5b0d146aSStefano Zampini     ierr = PetscSFReduceBegin(sf,MPIU_INT,sf->degreetmp-sf->minleaf,sf->degree,MPI_SUM);CHKERRQ(ierr);
155795fce210SBarry Smith   }
155895fce210SBarry Smith   *degree = NULL;
155995fce210SBarry Smith   PetscFunctionReturn(0);
156095fce210SBarry Smith }
156195fce210SBarry Smith 
156295fce210SBarry Smith /*@C
156395fce210SBarry Smith    PetscSFComputeDegreeEnd - complete computation of degree for each root vertex, started with PetscSFComputeDegreeBegin()
156495fce210SBarry Smith 
156595fce210SBarry Smith    Collective
156695fce210SBarry Smith 
156795fce210SBarry Smith    Input Arguments:
156895fce210SBarry Smith .  sf - star forest
156995fce210SBarry Smith 
157095fce210SBarry Smith    Output Arguments:
157195fce210SBarry Smith .  degree - degree of each root vertex
157295fce210SBarry Smith 
157395fce210SBarry Smith    Level: developer
157495fce210SBarry Smith 
1575ffe67aa5SVáclav Hapla    Notes:
1576ffe67aa5SVáclav Hapla    The returned array is owned by PetscSF and automatically freed by PetscSFDestroy(). Hence no need to call PetscFree() on it.
1577ffe67aa5SVáclav Hapla 
157895fce210SBarry Smith .seealso:
157995fce210SBarry Smith @*/
158095fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf,const PetscInt **degree)
158195fce210SBarry Smith {
158295fce210SBarry Smith   PetscErrorCode ierr;
158395fce210SBarry Smith 
158495fce210SBarry Smith   PetscFunctionBegin;
158595fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
158695fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
158729046d53SLisandro Dalcin   PetscValidPointer(degree,2);
158895fce210SBarry Smith   if (!sf->degreeknown) {
158929046d53SLisandro Dalcin     if (!sf->degreetmp) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFComputeDegreeBegin() before PetscSFComputeDegreeEnd()");
1590*5b0d146aSStefano Zampini     ierr = PetscSFReduceEnd(sf,MPIU_INT,sf->degreetmp-sf->minleaf,sf->degree,MPI_SUM);CHKERRQ(ierr);
159195fce210SBarry Smith     ierr = PetscFree(sf->degreetmp);CHKERRQ(ierr);
159295fce210SBarry Smith     sf->degreeknown = PETSC_TRUE;
159395fce210SBarry Smith   }
159495fce210SBarry Smith   *degree = sf->degree;
159595fce210SBarry Smith   PetscFunctionReturn(0);
159695fce210SBarry Smith }
159795fce210SBarry Smith 
1598673100f5SVaclav Hapla 
1599673100f5SVaclav Hapla /*@C
160066dfcd1aSVaclav Hapla    PetscSFComputeMultiRootOriginalNumbering - Returns original numbering of multi-roots (roots of multi-SF returned by PetscSFGetMultiSF()).
160166dfcd1aSVaclav Hapla    Each multi-root is assigned index of the corresponding original root.
1602673100f5SVaclav Hapla 
1603673100f5SVaclav Hapla    Collective
1604673100f5SVaclav Hapla 
1605673100f5SVaclav Hapla    Input Arguments:
1606673100f5SVaclav Hapla +  sf - star forest
1607673100f5SVaclav Hapla -  degree - degree of each root vertex, computed with PetscSFComputeDegreeBegin()/PetscSFComputeDegreeEnd()
1608673100f5SVaclav Hapla 
1609673100f5SVaclav Hapla    Output Arguments:
161066dfcd1aSVaclav Hapla +  nMultiRoots - (optional) number of multi-roots (roots of multi-SF)
161166dfcd1aSVaclav Hapla -  multiRootsOrigNumbering - original indices of multi-roots; length of this array is nMultiRoots
1612673100f5SVaclav Hapla 
1613673100f5SVaclav Hapla    Level: developer
1614673100f5SVaclav Hapla 
1615ffe67aa5SVáclav Hapla    Notes:
1616ffe67aa5SVáclav Hapla    The returned array multiRootsOrigNumbering is newly allocated and should be destroyed with PetscFree() when no longer needed.
1617ffe67aa5SVáclav Hapla 
1618673100f5SVaclav Hapla .seealso: PetscSFComputeDegreeBegin(), PetscSFComputeDegreeEnd(), PetscSFGetMultiSF()
1619673100f5SVaclav Hapla @*/
162066dfcd1aSVaclav Hapla PetscErrorCode PetscSFComputeMultiRootOriginalNumbering(PetscSF sf, const PetscInt degree[], PetscInt *nMultiRoots, PetscInt *multiRootsOrigNumbering[])
1621673100f5SVaclav Hapla {
1622673100f5SVaclav Hapla   PetscSF             msf;
1623673100f5SVaclav Hapla   PetscInt            i, j, k, nroots, nmroots;
1624673100f5SVaclav Hapla   PetscErrorCode      ierr;
1625673100f5SVaclav Hapla 
1626673100f5SVaclav Hapla   PetscFunctionBegin;
1627673100f5SVaclav Hapla   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1628673100f5SVaclav Hapla   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
162929328920SVaclav Hapla   if (nroots) PetscValidIntPointer(degree,2);
163066dfcd1aSVaclav Hapla   if (nMultiRoots) PetscValidIntPointer(nMultiRoots,3);
163166dfcd1aSVaclav Hapla   PetscValidPointer(multiRootsOrigNumbering,4);
1632673100f5SVaclav Hapla   ierr = PetscSFGetMultiSF(sf,&msf);CHKERRQ(ierr);
1633673100f5SVaclav Hapla   ierr = PetscSFGetGraph(msf, &nmroots, NULL, NULL, NULL);CHKERRQ(ierr);
163466dfcd1aSVaclav Hapla   ierr = PetscMalloc1(nmroots, multiRootsOrigNumbering);CHKERRQ(ierr);
1635673100f5SVaclav Hapla   for (i=0,j=0,k=0; i<nroots; i++) {
1636673100f5SVaclav Hapla     if (!degree[i]) continue;
1637673100f5SVaclav Hapla     for (j=0; j<degree[i]; j++,k++) {
163866dfcd1aSVaclav Hapla       (*multiRootsOrigNumbering)[k] = i;
1639673100f5SVaclav Hapla     }
1640673100f5SVaclav Hapla   }
1641673100f5SVaclav Hapla #if defined(PETSC_USE_DEBUG)
1642673100f5SVaclav Hapla   if (PetscUnlikely(k != nmroots)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sanity check fail");
1643673100f5SVaclav Hapla #endif
164466dfcd1aSVaclav Hapla   if (nMultiRoots) *nMultiRoots = nmroots;
1645673100f5SVaclav Hapla   PetscFunctionReturn(0);
1646673100f5SVaclav Hapla }
1647673100f5SVaclav Hapla 
164895fce210SBarry Smith /*@C
164995fce210SBarry Smith    PetscSFGatherBegin - begin pointwise gather of all leaves into multi-roots, to be completed with PetscSFGatherEnd()
165095fce210SBarry Smith 
165195fce210SBarry Smith    Collective
165295fce210SBarry Smith 
165395fce210SBarry Smith    Input Arguments:
165495fce210SBarry Smith +  sf - star forest
165595fce210SBarry Smith .  unit - data type
165695fce210SBarry Smith -  leafdata - leaf data to gather to roots
165795fce210SBarry Smith 
165895fce210SBarry Smith    Output Argument:
165995fce210SBarry Smith .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
166095fce210SBarry Smith 
166195fce210SBarry Smith    Level: intermediate
166295fce210SBarry Smith 
166395fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
166495fce210SBarry Smith @*/
166595fce210SBarry Smith PetscErrorCode PetscSFGatherBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
166695fce210SBarry Smith {
166795fce210SBarry Smith   PetscErrorCode ierr;
166895fce210SBarry Smith   PetscSF        multi;
166995fce210SBarry Smith 
167095fce210SBarry Smith   PetscFunctionBegin;
167195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
167229046d53SLisandro Dalcin   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
167395fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
16748bfbc91cSJed Brown   ierr = PetscSFReduceBegin(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr);
167595fce210SBarry Smith   PetscFunctionReturn(0);
167695fce210SBarry Smith }
167795fce210SBarry Smith 
167895fce210SBarry Smith /*@C
167995fce210SBarry Smith    PetscSFGatherEnd - ends pointwise gather operation that was started with PetscSFGatherBegin()
168095fce210SBarry Smith 
168195fce210SBarry Smith    Collective
168295fce210SBarry Smith 
168395fce210SBarry Smith    Input Arguments:
168495fce210SBarry Smith +  sf - star forest
168595fce210SBarry Smith .  unit - data type
168695fce210SBarry Smith -  leafdata - leaf data to gather to roots
168795fce210SBarry Smith 
168895fce210SBarry Smith    Output Argument:
168995fce210SBarry Smith .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
169095fce210SBarry Smith 
169195fce210SBarry Smith    Level: intermediate
169295fce210SBarry Smith 
169395fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
169495fce210SBarry Smith @*/
169595fce210SBarry Smith PetscErrorCode PetscSFGatherEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
169695fce210SBarry Smith {
169795fce210SBarry Smith   PetscErrorCode ierr;
169895fce210SBarry Smith   PetscSF        multi;
169995fce210SBarry Smith 
170095fce210SBarry Smith   PetscFunctionBegin;
170195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
170295fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
170395fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
17048bfbc91cSJed Brown   ierr = PetscSFReduceEnd(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr);
170595fce210SBarry Smith   PetscFunctionReturn(0);
170695fce210SBarry Smith }
170795fce210SBarry Smith 
170895fce210SBarry Smith /*@C
170995fce210SBarry Smith    PetscSFScatterBegin - begin pointwise scatter operation from multi-roots to leaves, to be completed with PetscSFScatterEnd()
171095fce210SBarry Smith 
171195fce210SBarry Smith    Collective
171295fce210SBarry Smith 
171395fce210SBarry Smith    Input Arguments:
171495fce210SBarry Smith +  sf - star forest
171595fce210SBarry Smith .  unit - data type
171695fce210SBarry Smith -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
171795fce210SBarry Smith 
171895fce210SBarry Smith    Output Argument:
171995fce210SBarry Smith .  leafdata - leaf data to be update with personal data from each respective root
172095fce210SBarry Smith 
172195fce210SBarry Smith    Level: intermediate
172295fce210SBarry Smith 
172395fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
172495fce210SBarry Smith @*/
172595fce210SBarry Smith PetscErrorCode PetscSFScatterBegin(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
172695fce210SBarry Smith {
172795fce210SBarry Smith   PetscErrorCode ierr;
172895fce210SBarry Smith   PetscSF        multi;
172995fce210SBarry Smith 
173095fce210SBarry Smith   PetscFunctionBegin;
173195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
173295fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
173395fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
173495fce210SBarry Smith   ierr = PetscSFBcastBegin(multi,unit,multirootdata,leafdata);CHKERRQ(ierr);
173595fce210SBarry Smith   PetscFunctionReturn(0);
173695fce210SBarry Smith }
173795fce210SBarry Smith 
173895fce210SBarry Smith /*@C
173995fce210SBarry Smith    PetscSFScatterEnd - ends pointwise scatter operation that was started with PetscSFScatterBegin()
174095fce210SBarry Smith 
174195fce210SBarry Smith    Collective
174295fce210SBarry Smith 
174395fce210SBarry Smith    Input Arguments:
174495fce210SBarry Smith +  sf - star forest
174595fce210SBarry Smith .  unit - data type
174695fce210SBarry Smith -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
174795fce210SBarry Smith 
174895fce210SBarry Smith    Output Argument:
174995fce210SBarry Smith .  leafdata - leaf data to be update with personal data from each respective root
175095fce210SBarry Smith 
175195fce210SBarry Smith    Level: intermediate
175295fce210SBarry Smith 
175395fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
175495fce210SBarry Smith @*/
175595fce210SBarry Smith PetscErrorCode PetscSFScatterEnd(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
175695fce210SBarry Smith {
175795fce210SBarry Smith   PetscErrorCode ierr;
175895fce210SBarry Smith   PetscSF        multi;
175995fce210SBarry Smith 
176095fce210SBarry Smith   PetscFunctionBegin;
176195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
176295fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
176395fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
176495fce210SBarry Smith   ierr = PetscSFBcastEnd(multi,unit,multirootdata,leafdata);CHKERRQ(ierr);
176595fce210SBarry Smith   PetscFunctionReturn(0);
176695fce210SBarry Smith }
1767a7b3aa13SAta Mesgarnejad 
1768a072220fSLawrence Mitchell static PetscErrorCode PetscSFCheckLeavesUnique_Private(PetscSF sf)
1769a072220fSLawrence Mitchell {
1770a072220fSLawrence Mitchell #if defined(PETSC_USE_DEBUG)
1771a072220fSLawrence Mitchell   PetscInt        i, n, nleaves;
1772a072220fSLawrence Mitchell   const PetscInt *ilocal = NULL;
1773a072220fSLawrence Mitchell   PetscHSetI      seen;
1774a072220fSLawrence Mitchell   PetscErrorCode  ierr;
1775a072220fSLawrence Mitchell 
1776a072220fSLawrence Mitchell   PetscFunctionBegin;
1777a072220fSLawrence Mitchell   ierr = PetscSFGetGraph(sf,NULL,&nleaves,&ilocal,NULL);CHKERRQ(ierr);
1778a072220fSLawrence Mitchell   ierr = PetscHSetICreate(&seen);CHKERRQ(ierr);
1779a072220fSLawrence Mitchell   for (i = 0; i < nleaves; i++) {
1780a072220fSLawrence Mitchell     const PetscInt leaf = ilocal ? ilocal[i] : i;
1781a072220fSLawrence Mitchell     ierr = PetscHSetIAdd(seen,leaf);CHKERRQ(ierr);
1782a072220fSLawrence Mitchell   }
1783a072220fSLawrence Mitchell   ierr = PetscHSetIGetSize(seen,&n);CHKERRQ(ierr);
1784a072220fSLawrence Mitchell   if (n != nleaves) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Provided leaves have repeated values: all leaves must be unique");
1785a072220fSLawrence Mitchell   ierr = PetscHSetIDestroy(&seen);CHKERRQ(ierr);
1786a072220fSLawrence Mitchell   PetscFunctionReturn(0);
1787a072220fSLawrence Mitchell #else
1788a072220fSLawrence Mitchell   PetscFunctionBegin;
1789a072220fSLawrence Mitchell   PetscFunctionReturn(0);
1790a072220fSLawrence Mitchell #endif
1791a072220fSLawrence Mitchell }
179254729392SStefano Zampini 
1793a7b3aa13SAta Mesgarnejad /*@
179404c0ada0SJunchao Zhang   PetscSFCompose - Compose a new PetscSF by putting the second SF under the first one in a top (roots) down (leaves) view
1795a7b3aa13SAta Mesgarnejad 
1796a7b3aa13SAta Mesgarnejad   Input Parameters:
179754729392SStefano Zampini + sfA - The first PetscSF
179854729392SStefano Zampini - sfB - The second PetscSF
1799a7b3aa13SAta Mesgarnejad 
1800a7b3aa13SAta Mesgarnejad   Output Parameters:
180154729392SStefano Zampini . sfBA - The composite SF
1802a7b3aa13SAta Mesgarnejad 
1803a7b3aa13SAta Mesgarnejad   Level: developer
1804a7b3aa13SAta Mesgarnejad 
1805a072220fSLawrence Mitchell   Notes:
180654729392SStefano Zampini   Currently, the two SFs must be defined on congruent communicators and they must be true star
180754729392SStefano Zampini   forests, i.e. the same leaf is not connected with different roots.
180854729392SStefano Zampini 
180954729392SStefano Zampini   sfA's leaf space and sfB's root space might be partially overlapped. The composition builds
181054729392SStefano Zampini   a graph with sfA's roots and sfB's leaves only when there is a path between them. Unconnected
181154729392SStefano Zampini   nodes (roots or leaves) are not in sfBA. Doing a Bcast on the new SF is equivalent to doing a
181254729392SStefano Zampini   Bcast on sfA, then a Bcast on sfB, on connected nodes.
1813a072220fSLawrence Mitchell 
181404c0ada0SJunchao Zhang .seealso: PetscSF, PetscSFComposeInverse(), PetscSFGetGraph(), PetscSFSetGraph()
1815a7b3aa13SAta Mesgarnejad @*/
1816a7b3aa13SAta Mesgarnejad PetscErrorCode PetscSFCompose(PetscSF sfA,PetscSF sfB,PetscSF *sfBA)
1817a7b3aa13SAta Mesgarnejad {
181804c0ada0SJunchao Zhang   PetscErrorCode    ierr;
1819a7b3aa13SAta Mesgarnejad   const PetscSFNode *remotePointsA,*remotePointsB;
1820d41018fbSJunchao Zhang   PetscSFNode       *remotePointsBA=NULL,*reorderedRemotePointsA = NULL,*leafdataB;
182154729392SStefano Zampini   const PetscInt    *localPointsA,*localPointsB;
182254729392SStefano Zampini   PetscInt          *localPointsBA;
182354729392SStefano Zampini   PetscInt          i,numRootsA,numLeavesA,numRootsB,numLeavesB,minleaf,maxleaf,numLeavesBA;
182454729392SStefano Zampini   PetscBool         denseB;
1825a7b3aa13SAta Mesgarnejad 
1826a7b3aa13SAta Mesgarnejad   PetscFunctionBegin;
1827a7b3aa13SAta Mesgarnejad   PetscValidHeaderSpecific(sfA,PETSCSF_CLASSID,1);
182829046d53SLisandro Dalcin   PetscSFCheckGraphSet(sfA,1);
182929046d53SLisandro Dalcin   PetscValidHeaderSpecific(sfB,PETSCSF_CLASSID,2);
183029046d53SLisandro Dalcin   PetscSFCheckGraphSet(sfB,2);
183154729392SStefano Zampini   PetscCheckSameComm(sfA,1,sfB,2);
183229046d53SLisandro Dalcin   PetscValidPointer(sfBA,3);
183354729392SStefano Zampini   ierr = PetscSFCheckLeavesUnique_Private(sfA);CHKERRQ(ierr);
183454729392SStefano Zampini   ierr = PetscSFCheckLeavesUnique_Private(sfB);CHKERRQ(ierr);
183554729392SStefano Zampini 
1836a7b3aa13SAta Mesgarnejad   ierr = PetscSFGetGraph(sfA,&numRootsA,&numLeavesA,&localPointsA,&remotePointsA);CHKERRQ(ierr);
1837a7b3aa13SAta Mesgarnejad   ierr = PetscSFGetGraph(sfB,&numRootsB,&numLeavesB,&localPointsB,&remotePointsB);CHKERRQ(ierr);
1838d41018fbSJunchao Zhang   if (localPointsA) {
183954729392SStefano Zampini     ierr = PetscMalloc1(numRootsB,&reorderedRemotePointsA);CHKERRQ(ierr);
184054729392SStefano Zampini     for (i=0; i<numRootsB; i++) {
184154729392SStefano Zampini       reorderedRemotePointsA[i].rank = -1;
184254729392SStefano Zampini       reorderedRemotePointsA[i].index = -1;
184354729392SStefano Zampini     }
184454729392SStefano Zampini     for (i=0; i<numLeavesA; i++) {
184554729392SStefano Zampini       if (localPointsA[i] >= numRootsB) continue;
184654729392SStefano Zampini       reorderedRemotePointsA[localPointsA[i]] = remotePointsA[i];
184754729392SStefano Zampini     }
1848d41018fbSJunchao Zhang     remotePointsA = reorderedRemotePointsA;
1849d41018fbSJunchao Zhang   }
1850d41018fbSJunchao Zhang   ierr = PetscSFGetLeafRange(sfB,&minleaf,&maxleaf);CHKERRQ(ierr);
1851d41018fbSJunchao Zhang   ierr = PetscMalloc1(maxleaf-minleaf+1,&leafdataB);CHKERRQ(ierr);
1852d41018fbSJunchao Zhang   ierr = PetscSFBcastBegin(sfB,MPIU_2INT,remotePointsA,leafdataB-minleaf);CHKERRQ(ierr);
1853d41018fbSJunchao Zhang   ierr = PetscSFBcastEnd(sfB,MPIU_2INT,remotePointsA,leafdataB-minleaf);CHKERRQ(ierr);
1854d41018fbSJunchao Zhang   ierr = PetscFree(reorderedRemotePointsA);CHKERRQ(ierr);
1855d41018fbSJunchao Zhang 
185654729392SStefano Zampini   denseB = (PetscBool)!localPointsB;
185754729392SStefano Zampini   for (i=0,numLeavesBA=0; i<numLeavesB; i++) {
185854729392SStefano Zampini     if (leafdataB[localPointsB ? localPointsB[i]-minleaf : i].rank == -1) denseB = PETSC_FALSE;
185954729392SStefano Zampini     else numLeavesBA++;
186054729392SStefano Zampini   }
186154729392SStefano Zampini   if (denseB) {
1862d41018fbSJunchao Zhang     localPointsBA  = NULL;
1863d41018fbSJunchao Zhang     remotePointsBA = leafdataB;
1864d41018fbSJunchao Zhang   } else {
186554729392SStefano Zampini     ierr = PetscMalloc1(numLeavesBA,&localPointsBA);CHKERRQ(ierr);
186654729392SStefano Zampini     ierr = PetscMalloc1(numLeavesBA,&remotePointsBA);CHKERRQ(ierr);
186754729392SStefano Zampini     for (i=0,numLeavesBA=0; i<numLeavesB; i++) {
186854729392SStefano Zampini       const PetscInt l = localPointsB ? localPointsB[i] : i;
186954729392SStefano Zampini 
187054729392SStefano Zampini       if (leafdataB[l-minleaf].rank == -1) continue;
187154729392SStefano Zampini       remotePointsBA[numLeavesBA] = leafdataB[l-minleaf];
187254729392SStefano Zampini       localPointsBA[numLeavesBA] = l;
187354729392SStefano Zampini       numLeavesBA++;
187454729392SStefano Zampini     }
1875d41018fbSJunchao Zhang     ierr = PetscFree(leafdataB);CHKERRQ(ierr);
1876d41018fbSJunchao Zhang   }
187754729392SStefano Zampini   ierr = PetscSFCreate(PetscObjectComm((PetscObject)sfA),sfBA);CHKERRQ(ierr);
187854729392SStefano Zampini   ierr = PetscSFSetGraph(*sfBA,numRootsA,numLeavesBA,localPointsBA,PETSC_OWN_POINTER,remotePointsBA,PETSC_OWN_POINTER);CHKERRQ(ierr);
1879a7b3aa13SAta Mesgarnejad   PetscFunctionReturn(0);
1880a7b3aa13SAta Mesgarnejad }
18811c6ba672SJunchao Zhang 
188204c0ada0SJunchao Zhang /*@
188354729392SStefano Zampini   PetscSFComposeInverse - Compose a new PetscSF by putting the inverse of the second SF under the first one
188404c0ada0SJunchao Zhang 
188504c0ada0SJunchao Zhang   Input Parameters:
188654729392SStefano Zampini + sfA - The first PetscSF
188754729392SStefano Zampini - sfB - The second PetscSF
188804c0ada0SJunchao Zhang 
188904c0ada0SJunchao Zhang   Output Parameters:
189054729392SStefano Zampini . sfBA - The composite SF.
189104c0ada0SJunchao Zhang 
189204c0ada0SJunchao Zhang   Level: developer
189304c0ada0SJunchao Zhang 
189454729392SStefano Zampini   Notes:
189554729392SStefano Zampini   Currently, the two SFs must be defined on congruent communicators and they must be true star
189654729392SStefano Zampini   forests, i.e. the same leaf is not connected with different roots. Even more, all roots of the
189754729392SStefano Zampini   second SF must have a degree of 1, i.e., no roots have more than one leaf connected.
189854729392SStefano Zampini 
189954729392SStefano Zampini   sfA's leaf space and sfB's leaf space might be partially overlapped. The composition builds
190054729392SStefano Zampini   a graph with sfA's roots and sfB's roots only when there is a path between them. Unconnected
190154729392SStefano Zampini   roots are not in sfBA. Doing a Bcast on the new SF is equivalent to doing a Bcast on sfA, then
190254729392SStefano Zampini   a Reduce on sfB, on connected roots.
190354729392SStefano Zampini 
190454729392SStefano Zampini .seealso: PetscSF, PetscSFCompose(), PetscSFGetGraph(), PetscSFSetGraph(), PetscSFCreateInverseSF()
190504c0ada0SJunchao Zhang @*/
190604c0ada0SJunchao Zhang PetscErrorCode PetscSFComposeInverse(PetscSF sfA,PetscSF sfB,PetscSF *sfBA)
190704c0ada0SJunchao Zhang {
190804c0ada0SJunchao Zhang   PetscErrorCode    ierr;
190904c0ada0SJunchao Zhang   const PetscSFNode *remotePointsA,*remotePointsB;
191004c0ada0SJunchao Zhang   PetscSFNode       *remotePointsBA;
191104c0ada0SJunchao Zhang   const PetscInt    *localPointsA,*localPointsB;
191254729392SStefano Zampini   PetscSFNode       *reorderedRemotePointsA = NULL;
191354729392SStefano Zampini   PetscInt          i,numRootsA,numLeavesA,numLeavesBA,numRootsB,numLeavesB,minleaf,maxleaf,*localPointsBA;
1914*5b0d146aSStefano Zampini   MPI_Op            op;
1915*5b0d146aSStefano Zampini #if defined(PETSC_USE_64BIT_INDICES)
1916*5b0d146aSStefano Zampini   PetscBool         iswin;
1917*5b0d146aSStefano Zampini #endif
191804c0ada0SJunchao Zhang 
191904c0ada0SJunchao Zhang   PetscFunctionBegin;
192004c0ada0SJunchao Zhang   PetscValidHeaderSpecific(sfA,PETSCSF_CLASSID,1);
192104c0ada0SJunchao Zhang   PetscSFCheckGraphSet(sfA,1);
192204c0ada0SJunchao Zhang   PetscValidHeaderSpecific(sfB,PETSCSF_CLASSID,2);
192304c0ada0SJunchao Zhang   PetscSFCheckGraphSet(sfB,2);
192454729392SStefano Zampini   PetscCheckSameComm(sfA,1,sfB,2);
192504c0ada0SJunchao Zhang   PetscValidPointer(sfBA,3);
192654729392SStefano Zampini   ierr = PetscSFCheckLeavesUnique_Private(sfA);CHKERRQ(ierr);
192754729392SStefano Zampini   ierr = PetscSFCheckLeavesUnique_Private(sfB);CHKERRQ(ierr);
192854729392SStefano Zampini 
192904c0ada0SJunchao Zhang   ierr = PetscSFGetGraph(sfA, &numRootsA, &numLeavesA, &localPointsA, &remotePointsA);CHKERRQ(ierr);
193004c0ada0SJunchao Zhang   ierr = PetscSFGetGraph(sfB, &numRootsB, &numLeavesB, &localPointsB, &remotePointsB);CHKERRQ(ierr);
1931*5b0d146aSStefano Zampini 
1932*5b0d146aSStefano Zampini   /* TODO: Check roots of sfB have degree of 1 */
1933*5b0d146aSStefano Zampini   /* Once we implement it, we can replace the MPI_MAXLOC
1934*5b0d146aSStefano Zampini      with MPIU_REPLACE. In that case, MPI_MAXLOC and MPIU_REPLACE have the same effect.
1935*5b0d146aSStefano Zampini      We use MPI_MAXLOC only to have a deterministic output from this routine if
1936*5b0d146aSStefano Zampini      the root condition is not meet.
1937*5b0d146aSStefano Zampini    */
1938*5b0d146aSStefano Zampini   op = MPI_MAXLOC;
1939*5b0d146aSStefano Zampini #if defined(PETSC_USE_64BIT_INDICES)
1940*5b0d146aSStefano Zampini   /* we accept a non-deterministic output (if any) with PETSCSFWINDOW, since MPI_MAXLOC cannot operate on MPIU_2INT with MPI_Accumulate */
1941*5b0d146aSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)sfB,PETSCSFWINDOW,&iswin);CHKERRQ(ierr);
1942*5b0d146aSStefano Zampini   if (iswin) op = MPIU_REPLACE;
1943*5b0d146aSStefano Zampini #endif
1944*5b0d146aSStefano Zampini 
194554729392SStefano Zampini   ierr = PetscSFGetLeafRange(sfB, &minleaf, &maxleaf);CHKERRQ(ierr);
194654729392SStefano Zampini   ierr = PetscMalloc1(maxleaf - minleaf + 1,&reorderedRemotePointsA);CHKERRQ(ierr);
194754729392SStefano Zampini   for (i=0; i<maxleaf - minleaf + 1; i++) {
194854729392SStefano Zampini     reorderedRemotePointsA[i].rank = -1;
194954729392SStefano Zampini     reorderedRemotePointsA[i].index = -1;
195054729392SStefano Zampini   }
195154729392SStefano Zampini   if (localPointsA) {
195254729392SStefano Zampini     for (i=0; i<numLeavesA; i++) {
195354729392SStefano Zampini       if (localPointsA[i] > maxleaf || localPointsA[i] < minleaf) continue;
195454729392SStefano Zampini       reorderedRemotePointsA[localPointsA[i] - minleaf] = remotePointsA[i];
195554729392SStefano Zampini     }
195654729392SStefano Zampini   } else {
195754729392SStefano Zampini     for (i=0; i<numLeavesA; i++) {
195854729392SStefano Zampini       if (i > maxleaf || i < minleaf) continue;
195954729392SStefano Zampini       reorderedRemotePointsA[i - minleaf] = remotePointsA[i];
196054729392SStefano Zampini     }
196154729392SStefano Zampini   }
196254729392SStefano Zampini 
196354729392SStefano Zampini   ierr = PetscMalloc1(numRootsB,&localPointsBA);CHKERRQ(ierr);
196404c0ada0SJunchao Zhang   ierr = PetscMalloc1(numRootsB,&remotePointsBA);CHKERRQ(ierr);
196554729392SStefano Zampini   for (i=0; i<numRootsB; i++) {
196654729392SStefano Zampini     remotePointsBA[i].rank = -1;
196754729392SStefano Zampini     remotePointsBA[i].index = -1;
196854729392SStefano Zampini   }
196954729392SStefano Zampini 
1970*5b0d146aSStefano Zampini   ierr = PetscSFReduceBegin(sfB,MPIU_2INT,reorderedRemotePointsA-minleaf,remotePointsBA,op);CHKERRQ(ierr);
1971*5b0d146aSStefano Zampini   ierr = PetscSFReduceEnd(sfB,MPIU_2INT,reorderedRemotePointsA-minleaf,remotePointsBA,op);CHKERRQ(ierr);
197254729392SStefano Zampini   ierr = PetscFree(reorderedRemotePointsA);CHKERRQ(ierr);
197354729392SStefano Zampini   for (i=0,numLeavesBA=0; i<numRootsB; i++) {
197454729392SStefano Zampini     if (remotePointsBA[i].rank == -1) continue;
197554729392SStefano Zampini     remotePointsBA[numLeavesBA].rank = remotePointsBA[i].rank;
197654729392SStefano Zampini     remotePointsBA[numLeavesBA].index = remotePointsBA[i].index;
197754729392SStefano Zampini     localPointsBA[numLeavesBA] = i;
197854729392SStefano Zampini     numLeavesBA++;
197954729392SStefano Zampini   }
198054729392SStefano Zampini   ierr = PetscSFCreate(PetscObjectComm((PetscObject)sfA),sfBA);CHKERRQ(ierr);
198154729392SStefano Zampini   ierr = PetscSFSetGraph(*sfBA,numRootsA,numLeavesBA,localPointsBA,PETSC_OWN_POINTER,remotePointsBA,PETSC_OWN_POINTER);CHKERRQ(ierr);
198204c0ada0SJunchao Zhang   PetscFunctionReturn(0);
198304c0ada0SJunchao Zhang }
198404c0ada0SJunchao Zhang 
19851c6ba672SJunchao Zhang /*
19861c6ba672SJunchao Zhang   PetscSFCreateLocalSF_Private - Creates a local PetscSF that only has intra-process edges of the global PetscSF
19871c6ba672SJunchao Zhang 
19881c6ba672SJunchao Zhang   Input Parameters:
19891c6ba672SJunchao Zhang . sf - The global PetscSF
19901c6ba672SJunchao Zhang 
19911c6ba672SJunchao Zhang   Output Parameters:
19921c6ba672SJunchao Zhang . out - The local PetscSF
19931c6ba672SJunchao Zhang  */
19941c6ba672SJunchao Zhang PetscErrorCode PetscSFCreateLocalSF_Private(PetscSF sf,PetscSF *out)
19951c6ba672SJunchao Zhang {
19961c6ba672SJunchao Zhang   MPI_Comm           comm;
19971c6ba672SJunchao Zhang   PetscMPIInt        myrank;
19981c6ba672SJunchao Zhang   const PetscInt     *ilocal;
19991c6ba672SJunchao Zhang   const PetscSFNode  *iremote;
20001c6ba672SJunchao Zhang   PetscInt           i,j,nroots,nleaves,lnleaves,*lilocal;
20011c6ba672SJunchao Zhang   PetscSFNode        *liremote;
20021c6ba672SJunchao Zhang   PetscSF            lsf;
20031c6ba672SJunchao Zhang   PetscErrorCode     ierr;
20041c6ba672SJunchao Zhang 
20051c6ba672SJunchao Zhang   PetscFunctionBegin;
20061c6ba672SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
20071c6ba672SJunchao Zhang   if (sf->ops->CreateLocalSF) {
20081c6ba672SJunchao Zhang     ierr = (*sf->ops->CreateLocalSF)(sf,out);CHKERRQ(ierr);
20091c6ba672SJunchao Zhang   } else {
20101c6ba672SJunchao Zhang     /* Could use PetscSFCreateEmbeddedLeafSF, but since we know the comm is PETSC_COMM_SELF, we can make it fast */
20111c6ba672SJunchao Zhang     ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
20121c6ba672SJunchao Zhang     ierr = MPI_Comm_rank(comm,&myrank);CHKERRQ(ierr);
20131c6ba672SJunchao Zhang 
20141c6ba672SJunchao Zhang     /* Find out local edges and build a local SF */
20151c6ba672SJunchao Zhang     ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
20161c6ba672SJunchao Zhang     for (i=lnleaves=0; i<nleaves; i++) {if (iremote[i].rank == (PetscInt)myrank) lnleaves++;}
20171c6ba672SJunchao Zhang     ierr = PetscMalloc1(lnleaves,&lilocal);CHKERRQ(ierr);
20181c6ba672SJunchao Zhang     ierr = PetscMalloc1(lnleaves,&liremote);CHKERRQ(ierr);
20191c6ba672SJunchao Zhang 
20201c6ba672SJunchao Zhang     for (i=j=0; i<nleaves; i++) {
20211c6ba672SJunchao Zhang       if (iremote[i].rank == (PetscInt)myrank) {
20221c6ba672SJunchao Zhang         lilocal[j]        = ilocal? ilocal[i] : i; /* ilocal=NULL for contiguous storage */
20231c6ba672SJunchao Zhang         liremote[j].rank  = 0; /* rank in PETSC_COMM_SELF */
20241c6ba672SJunchao Zhang         liremote[j].index = iremote[i].index;
20251c6ba672SJunchao Zhang         j++;
20261c6ba672SJunchao Zhang       }
20271c6ba672SJunchao Zhang     }
20281c6ba672SJunchao Zhang     ierr = PetscSFCreate(PETSC_COMM_SELF,&lsf);CHKERRQ(ierr);
20291c6ba672SJunchao Zhang     ierr = PetscSFSetGraph(lsf,nroots,lnleaves,lilocal,PETSC_OWN_POINTER,liremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
20301c6ba672SJunchao Zhang     ierr = PetscSFSetUp(lsf);CHKERRQ(ierr);
20311c6ba672SJunchao Zhang     *out = lsf;
20321c6ba672SJunchao Zhang   }
20331c6ba672SJunchao Zhang   PetscFunctionReturn(0);
20341c6ba672SJunchao Zhang }
2035dd5b3ca6SJunchao Zhang 
2036dd5b3ca6SJunchao Zhang /* Similar to PetscSFBcast, but only Bcast to leaves on rank 0 */
2037dd5b3ca6SJunchao Zhang PetscErrorCode PetscSFBcastToZero_Private(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
2038dd5b3ca6SJunchao Zhang {
2039dd5b3ca6SJunchao Zhang   PetscErrorCode     ierr;
2040eb02082bSJunchao Zhang   PetscMemType       rootmtype,leafmtype;
2041dd5b3ca6SJunchao Zhang 
2042dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
2043dd5b3ca6SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
2044dd5b3ca6SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
2045dd5b3ca6SJunchao Zhang   ierr = PetscLogEventBegin(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
2046eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
2047eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
2048dd5b3ca6SJunchao Zhang   if (sf->ops->BcastToZero) {
2049eb02082bSJunchao Zhang     ierr = (*sf->ops->BcastToZero)(sf,unit,rootmtype,rootdata,leafmtype,leafdata);CHKERRQ(ierr);
2050dd5b3ca6SJunchao Zhang   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSFBcastToZero_Private is not supported on this SF type");
2051dd5b3ca6SJunchao Zhang   ierr = PetscLogEventEnd(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
2052dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
2053dd5b3ca6SJunchao Zhang }
2054dd5b3ca6SJunchao Zhang 
2055