xref: /petsc/src/vec/is/sf/interface/sf.c (revision 245d98336dcf84d7e5386ee5cb055ed89e78f1a7)
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
12095fce210SBarry Smith    of available methods (for instance, window, pt2pt, 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()
27060263706SJed Brown -  -sf_rank_order - sort composite points for gathers and scatters in rank order, gathers are non-deterministic otherwise
27195fce210SBarry Smith 
27295fce210SBarry Smith    Level: intermediate
27395fce210SBarry Smith 
27495fce210SBarry Smith .seealso: PetscSFWindowSetSyncType()
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 
31195fce210SBarry Smith   PetscFunctionBegin;
31295fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
31395fce210SBarry Smith   PetscValidLogicalCollectiveBool(sf,flg,2);
31495fce210SBarry Smith   if (sf->multi) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Rank ordering must be set before first call to PetscSFGatherBegin() or PetscSFScatterBegin()");
31595fce210SBarry Smith   sf->rankorder = flg;
31695fce210SBarry Smith   PetscFunctionReturn(0);
31795fce210SBarry Smith }
31895fce210SBarry Smith 
3198af6ec1cSBarry Smith /*@
32095fce210SBarry Smith    PetscSFSetGraph - Set a parallel star forest
32195fce210SBarry Smith 
32295fce210SBarry Smith    Collective
32395fce210SBarry Smith 
32495fce210SBarry Smith    Input Arguments:
32595fce210SBarry Smith +  sf - star forest
32695fce210SBarry Smith .  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
32795fce210SBarry Smith .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
328c4e6a40aSLawrence Mitchell .  ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage (locations must be >= 0, enforced
329c4e6a40aSLawrence Mitchell during setup in debug mode)
33095fce210SBarry Smith .  localmode - copy mode for ilocal
331c4e6a40aSLawrence Mitchell .  iremote - remote locations of root vertices for each leaf on the current process (locations must be >= 0, enforced
332c4e6a40aSLawrence Mitchell during setup in debug mode)
33395fce210SBarry Smith -  remotemode - copy mode for iremote
33495fce210SBarry Smith 
33595fce210SBarry Smith    Level: intermediate
33695fce210SBarry Smith 
33795452b02SPatrick Sanan    Notes:
33895452b02SPatrick Sanan     In Fortran you must use PETSC_COPY_VALUES for localmode and remotemode
33938ab3f8aSBarry Smith 
3402ad1e87fSLisandro Dalcin    Developers Note: Local indices which are the identity permutation in the range [0,nleaves) are discarded as they
3412ad1e87fSLisandro Dalcin    encode contiguous storage. In such case, if localmode is PETSC_OWN_POINTER, the memory is deallocated as it is not
3422ad1e87fSLisandro Dalcin    needed
3432ad1e87fSLisandro Dalcin 
344c4e6a40aSLawrence Mitchell    Developers Note: This object does not necessarily encode a true star forest in the graph theoretic sense, since leaf
345c4e6a40aSLawrence Mitchell    indices are not required to be unique. Some functions, however, rely on unique leaf indices (checked in debug mode).
346c4e6a40aSLawrence Mitchell 
34795fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFGetGraph()
34895fce210SBarry Smith @*/
34995fce210SBarry Smith PetscErrorCode PetscSFSetGraph(PetscSF sf,PetscInt nroots,PetscInt nleaves,const PetscInt *ilocal,PetscCopyMode localmode,const PetscSFNode *iremote,PetscCopyMode remotemode)
35095fce210SBarry Smith {
35195fce210SBarry Smith   PetscErrorCode ierr;
35295fce210SBarry Smith 
35395fce210SBarry Smith   PetscFunctionBegin;
35495fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
35529046d53SLisandro Dalcin   if (nleaves > 0 && ilocal) PetscValidIntPointer(ilocal,4);
35629046d53SLisandro Dalcin   if (nleaves > 0) PetscValidPointer(iremote,6);
35729046d53SLisandro Dalcin   if (nroots  < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nroots %D, cannot be negative",nroots);
35895fce210SBarry Smith   if (nleaves < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nleaves %D, cannot be negative",nleaves);
35929046d53SLisandro Dalcin 
36095fce210SBarry Smith   ierr = PetscSFReset(sf);CHKERRQ(ierr);
36129046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr);
36229046d53SLisandro Dalcin 
36395fce210SBarry Smith   sf->nroots  = nroots;
36495fce210SBarry Smith   sf->nleaves = nleaves;
36529046d53SLisandro Dalcin 
36629046d53SLisandro Dalcin   if (nleaves && ilocal) {
36721c688dcSJed Brown     PetscInt i;
36829046d53SLisandro Dalcin     PetscInt minleaf = PETSC_MAX_INT;
36929046d53SLisandro Dalcin     PetscInt maxleaf = PETSC_MIN_INT;
3702ad1e87fSLisandro Dalcin     int      contiguous = 1;
37129046d53SLisandro Dalcin     for (i=0; i<nleaves; i++) {
37229046d53SLisandro Dalcin       minleaf = PetscMin(minleaf,ilocal[i]);
37329046d53SLisandro Dalcin       maxleaf = PetscMax(maxleaf,ilocal[i]);
3742ad1e87fSLisandro Dalcin       contiguous &= (ilocal[i] == i);
37529046d53SLisandro Dalcin     }
37629046d53SLisandro Dalcin     sf->minleaf = minleaf;
37729046d53SLisandro Dalcin     sf->maxleaf = maxleaf;
3782ad1e87fSLisandro Dalcin     if (contiguous) {
3792ad1e87fSLisandro Dalcin       if (localmode == PETSC_OWN_POINTER) {
3802ad1e87fSLisandro Dalcin         ierr = PetscFree(ilocal);CHKERRQ(ierr);
3812ad1e87fSLisandro Dalcin       }
3822ad1e87fSLisandro Dalcin       ilocal = NULL;
3832ad1e87fSLisandro Dalcin     }
38429046d53SLisandro Dalcin   } else {
38529046d53SLisandro Dalcin     sf->minleaf = 0;
38629046d53SLisandro Dalcin     sf->maxleaf = nleaves - 1;
38729046d53SLisandro Dalcin   }
38829046d53SLisandro Dalcin 
38929046d53SLisandro Dalcin   if (ilocal) {
39095fce210SBarry Smith     switch (localmode) {
39195fce210SBarry Smith     case PETSC_COPY_VALUES:
392785e854fSJed Brown       ierr = PetscMalloc1(nleaves,&sf->mine_alloc);CHKERRQ(ierr);
393580bdb30SBarry Smith       ierr = PetscArraycpy(sf->mine_alloc,ilocal,nleaves);CHKERRQ(ierr);
39495fce210SBarry Smith       sf->mine = sf->mine_alloc;
39595fce210SBarry Smith       break;
39695fce210SBarry Smith     case PETSC_OWN_POINTER:
39795fce210SBarry Smith       sf->mine_alloc = (PetscInt*)ilocal;
39895fce210SBarry Smith       sf->mine       = sf->mine_alloc;
39995fce210SBarry Smith       break;
40095fce210SBarry Smith     case PETSC_USE_POINTER:
40129046d53SLisandro Dalcin       sf->mine_alloc = NULL;
40295fce210SBarry Smith       sf->mine       = (PetscInt*)ilocal;
40395fce210SBarry Smith       break;
40495fce210SBarry Smith     default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown localmode");
40595fce210SBarry Smith     }
40695fce210SBarry Smith   }
40729046d53SLisandro Dalcin 
40895fce210SBarry Smith   switch (remotemode) {
40995fce210SBarry Smith   case PETSC_COPY_VALUES:
410785e854fSJed Brown     ierr = PetscMalloc1(nleaves,&sf->remote_alloc);CHKERRQ(ierr);
411580bdb30SBarry Smith     ierr = PetscArraycpy(sf->remote_alloc,iremote,nleaves);CHKERRQ(ierr);
41295fce210SBarry Smith     sf->remote = sf->remote_alloc;
41395fce210SBarry Smith     break;
41495fce210SBarry Smith   case PETSC_OWN_POINTER:
41595fce210SBarry Smith     sf->remote_alloc = (PetscSFNode*)iremote;
41695fce210SBarry Smith     sf->remote       = sf->remote_alloc;
41795fce210SBarry Smith     break;
41895fce210SBarry Smith   case PETSC_USE_POINTER:
41929046d53SLisandro Dalcin     sf->remote_alloc = NULL;
42095fce210SBarry Smith     sf->remote       = (PetscSFNode*)iremote;
42195fce210SBarry Smith     break;
42295fce210SBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown remotemode");
42395fce210SBarry Smith   }
42495fce210SBarry Smith 
425563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr);
42629046d53SLisandro Dalcin   sf->graphset = PETSC_TRUE;
42795fce210SBarry Smith   PetscFunctionReturn(0);
42895fce210SBarry Smith }
42995fce210SBarry Smith 
43029046d53SLisandro Dalcin /*@
431dd5b3ca6SJunchao Zhang   PetscSFSetGraphWithPattern - Sets the graph of an SF with a specific pattern
432dd5b3ca6SJunchao Zhang 
433dd5b3ca6SJunchao Zhang   Collective
434dd5b3ca6SJunchao Zhang 
435dd5b3ca6SJunchao Zhang   Input Parameters:
436dd5b3ca6SJunchao Zhang + sf      - The PetscSF
437dd5b3ca6SJunchao Zhang . map     - Layout of roots over all processes (insignificant when pattern is PETSCSF_PATTERN_ALLTOALL)
438dd5b3ca6SJunchao Zhang - pattern - One of PETSCSF_PATTERN_ALLGATHER, PETSCSF_PATTERN_GATHER, PETSCSF_PATTERN_ALLTOALL
439dd5b3ca6SJunchao Zhang 
440dd5b3ca6SJunchao Zhang   Notes:
441dd5b3ca6SJunchao Zhang   It is easier to explain PetscSFPattern using vectors. Suppose we have an MPI vector x and its layout is map.
442dd5b3ca6SJunchao Zhang   n and N are local and global sizes of x respectively.
443dd5b3ca6SJunchao Zhang 
444dd5b3ca6SJunchao Zhang   With PETSCSF_PATTERN_ALLGATHER, the routine creates a graph that if one does Bcast on it, it will copy x to
445dd5b3ca6SJunchao Zhang   sequential vectors y on all ranks.
446dd5b3ca6SJunchao Zhang 
447dd5b3ca6SJunchao Zhang   With PETSCSF_PATTERN_GATHER, the routine creates a graph that if one does Bcast on it, it will copy x to a
448dd5b3ca6SJunchao Zhang   sequential vector y on rank 0.
449dd5b3ca6SJunchao Zhang 
450dd5b3ca6SJunchao Zhang   In above cases, entries of x are roots and entries of y are leaves.
451dd5b3ca6SJunchao Zhang 
452dd5b3ca6SJunchao Zhang   With PETSCSF_PATTERN_ALLTOALL, map is insignificant. Suppose NP is size of sf's communicator. The routine
453dd5b3ca6SJunchao Zhang   creates a graph that every rank has NP leaves and NP roots. On rank i, its leaf j is connected to root i
454dd5b3ca6SJunchao 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
455dd5b3ca6SJunchao Zhang   not mean one can not send multiple items. One just needs to create a new MPI datatype for the mulptiple data
456dd5b3ca6SJunchao Zhang   items with MPI_Type_contiguous() and use that as the <unit> argument in SF routines.
457dd5b3ca6SJunchao Zhang 
458dd5b3ca6SJunchao Zhang   In this case, roots and leaves are symmetric.
459dd5b3ca6SJunchao Zhang 
460dd5b3ca6SJunchao Zhang   Level: intermediate
461dd5b3ca6SJunchao Zhang  @*/
462dd5b3ca6SJunchao Zhang PetscErrorCode PetscSFSetGraphWithPattern(PetscSF sf,PetscLayout map,PetscSFPattern pattern)
463dd5b3ca6SJunchao Zhang {
464dd5b3ca6SJunchao Zhang   MPI_Comm       comm;
465dd5b3ca6SJunchao Zhang   PetscInt       n,N,res[2];
466dd5b3ca6SJunchao Zhang   PetscMPIInt    rank,size;
467dd5b3ca6SJunchao Zhang   PetscSFType    type;
468dd5b3ca6SJunchao Zhang   PetscErrorCode ierr;
469dd5b3ca6SJunchao Zhang 
470dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
471dd5b3ca6SJunchao Zhang   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
472dd5b3ca6SJunchao Zhang   if (pattern < PETSCSF_PATTERN_ALLGATHER || pattern > PETSCSF_PATTERN_ALLTOALL) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Unsupported PetscSFPattern %D\n",pattern);
473dd5b3ca6SJunchao Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
474dd5b3ca6SJunchao Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
475dd5b3ca6SJunchao Zhang 
476dd5b3ca6SJunchao Zhang   if (pattern == PETSCSF_PATTERN_ALLTOALL) {
477dd5b3ca6SJunchao Zhang     type = PETSCSFALLTOALL;
478dd5b3ca6SJunchao Zhang     ierr = PetscLayoutCreate(comm,&sf->map);CHKERRQ(ierr);
479dd5b3ca6SJunchao Zhang     ierr = PetscLayoutSetLocalSize(sf->map,size);CHKERRQ(ierr);
480dd5b3ca6SJunchao Zhang     ierr = PetscLayoutSetSize(sf->map,((PetscInt)size)*size);CHKERRQ(ierr);
481dd5b3ca6SJunchao Zhang     ierr = PetscLayoutSetUp(sf->map);CHKERRQ(ierr);
482dd5b3ca6SJunchao Zhang   } else {
483dd5b3ca6SJunchao Zhang     ierr   = PetscLayoutGetLocalSize(map,&n);CHKERRQ(ierr);
484dd5b3ca6SJunchao Zhang     ierr   = PetscLayoutGetSize(map,&N);CHKERRQ(ierr);
485dd5b3ca6SJunchao Zhang     res[0] = n;
486dd5b3ca6SJunchao Zhang     res[1] = -n;
487dd5b3ca6SJunchao Zhang     /* Check if n are same over all ranks so that we can optimize it */
488dd5b3ca6SJunchao Zhang     ierr   = MPIU_Allreduce(MPI_IN_PLACE,res,2,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr);
489dd5b3ca6SJunchao Zhang     if (res[0] == -res[1]) { /* same n */
490dd5b3ca6SJunchao Zhang       type = (pattern == PETSCSF_PATTERN_ALLGATHER) ? PETSCSFALLGATHER  : PETSCSFGATHER;
491dd5b3ca6SJunchao Zhang     } else {
492dd5b3ca6SJunchao Zhang       type = (pattern == PETSCSF_PATTERN_ALLGATHER) ? PETSCSFALLGATHERV : PETSCSFGATHERV;
493dd5b3ca6SJunchao Zhang     }
494dd5b3ca6SJunchao Zhang     ierr = PetscLayoutReference(map,&sf->map);CHKERRQ(ierr);
495dd5b3ca6SJunchao Zhang   }
496dd5b3ca6SJunchao Zhang   ierr = PetscSFSetType(sf,type);CHKERRQ(ierr);
497dd5b3ca6SJunchao Zhang 
498dd5b3ca6SJunchao Zhang   sf->pattern = pattern;
499dd5b3ca6SJunchao Zhang   sf->mine    = NULL; /* Contiguous */
500dd5b3ca6SJunchao Zhang 
501dd5b3ca6SJunchao Zhang   /* Set nleaves, nroots here in case user calls PetscSFGetGraph, which is legal to call even before PetscSFSetUp is called.
502dd5b3ca6SJunchao Zhang      Also set other easy stuff.
503dd5b3ca6SJunchao Zhang    */
504dd5b3ca6SJunchao Zhang   if (pattern == PETSCSF_PATTERN_ALLGATHER) {
505dd5b3ca6SJunchao Zhang     sf->nleaves      = N;
506dd5b3ca6SJunchao Zhang     sf->nroots       = n;
507dd5b3ca6SJunchao Zhang     sf->nranks       = size;
508dd5b3ca6SJunchao Zhang     sf->minleaf      = 0;
509dd5b3ca6SJunchao Zhang     sf->maxleaf      = N - 1;
510dd5b3ca6SJunchao Zhang   } else if (pattern == PETSCSF_PATTERN_GATHER) {
511dd5b3ca6SJunchao Zhang     sf->nleaves      = rank ? 0 : N;
512dd5b3ca6SJunchao Zhang     sf->nroots       = n;
513dd5b3ca6SJunchao Zhang     sf->nranks       = rank ? 0 : size;
514dd5b3ca6SJunchao Zhang     sf->minleaf      = 0;
515dd5b3ca6SJunchao Zhang     sf->maxleaf      = rank ? -1 : N - 1;
516dd5b3ca6SJunchao Zhang   } else if (pattern == PETSCSF_PATTERN_ALLTOALL) {
517dd5b3ca6SJunchao Zhang     sf->nleaves      = size;
518dd5b3ca6SJunchao Zhang     sf->nroots       = size;
519dd5b3ca6SJunchao Zhang     sf->nranks       = size;
520dd5b3ca6SJunchao Zhang     sf->minleaf      = 0;
521dd5b3ca6SJunchao Zhang     sf->maxleaf      = size - 1;
522dd5b3ca6SJunchao Zhang   }
523dd5b3ca6SJunchao Zhang   sf->ndranks  = 0; /* We do not need to separate out distinguished ranks for patterned graphs to improve communication performance */
524dd5b3ca6SJunchao Zhang   sf->graphset = PETSC_TRUE;
525dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
526dd5b3ca6SJunchao Zhang }
527dd5b3ca6SJunchao Zhang 
528dd5b3ca6SJunchao Zhang /*@
52995fce210SBarry Smith    PetscSFCreateInverseSF - given a PetscSF in which all vertices have degree 1, creates the inverse map
53095fce210SBarry Smith 
53195fce210SBarry Smith    Collective
53295fce210SBarry Smith 
53395fce210SBarry Smith    Input Arguments:
534dd5b3ca6SJunchao Zhang 
53595fce210SBarry Smith .  sf - star forest to invert
53695fce210SBarry Smith 
53795fce210SBarry Smith    Output Arguments:
53895fce210SBarry Smith .  isf - inverse of sf
53995fce210SBarry Smith    Level: advanced
54095fce210SBarry Smith 
54195fce210SBarry Smith    Notes:
54295fce210SBarry Smith    All roots must have degree 1.
54395fce210SBarry Smith 
54495fce210SBarry Smith    The local space may be a permutation, but cannot be sparse.
54595fce210SBarry Smith 
54695fce210SBarry Smith .seealso: PetscSFSetGraph()
54795fce210SBarry Smith @*/
54895fce210SBarry Smith PetscErrorCode PetscSFCreateInverseSF(PetscSF sf,PetscSF *isf)
54995fce210SBarry Smith {
55095fce210SBarry Smith   PetscErrorCode ierr;
55195fce210SBarry Smith   PetscMPIInt    rank;
55295fce210SBarry Smith   PetscInt       i,nroots,nleaves,maxlocal,count,*newilocal;
55395fce210SBarry Smith   const PetscInt *ilocal;
55495fce210SBarry Smith   PetscSFNode    *roots,*leaves;
55595fce210SBarry Smith 
55695fce210SBarry Smith   PetscFunctionBegin;
55729046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
55829046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
55929046d53SLisandro Dalcin   PetscValidPointer(isf,2);
56029046d53SLisandro Dalcin 
56195fce210SBarry Smith   ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);CHKERRQ(ierr);
56229046d53SLisandro Dalcin   maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */
56329046d53SLisandro Dalcin 
56429046d53SLisandro Dalcin   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
565ae9aee6dSMatthew G. Knepley   ierr = PetscMalloc2(nroots,&roots,maxlocal,&leaves);CHKERRQ(ierr);
566ae9aee6dSMatthew G. Knepley   for (i=0; i<maxlocal; i++) {
56795fce210SBarry Smith     leaves[i].rank  = rank;
56895fce210SBarry Smith     leaves[i].index = i;
56995fce210SBarry Smith   }
57095fce210SBarry Smith   for (i=0; i <nroots; i++) {
57195fce210SBarry Smith     roots[i].rank  = -1;
57295fce210SBarry Smith     roots[i].index = -1;
57395fce210SBarry Smith   }
5748bfbc91cSJed Brown   ierr = PetscSFReduceBegin(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr);
5758bfbc91cSJed Brown   ierr = PetscSFReduceEnd(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr);
57695fce210SBarry Smith 
57795fce210SBarry Smith   /* Check whether our leaves are sparse */
57895fce210SBarry Smith   for (i=0,count=0; i<nroots; i++) if (roots[i].rank >= 0) count++;
57995fce210SBarry Smith   if (count == nroots) newilocal = NULL;
58095fce210SBarry Smith   else {                        /* Index for sparse leaves and compact "roots" array (which is to become our leaves). */
581785e854fSJed Brown     ierr = PetscMalloc1(count,&newilocal);CHKERRQ(ierr);
58295fce210SBarry Smith     for (i=0,count=0; i<nroots; i++) {
58395fce210SBarry Smith       if (roots[i].rank >= 0) {
58495fce210SBarry Smith         newilocal[count]   = i;
58595fce210SBarry Smith         roots[count].rank  = roots[i].rank;
58695fce210SBarry Smith         roots[count].index = roots[i].index;
58795fce210SBarry Smith         count++;
58895fce210SBarry Smith       }
58995fce210SBarry Smith     }
59095fce210SBarry Smith   }
59195fce210SBarry Smith 
59295fce210SBarry Smith   ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,isf);CHKERRQ(ierr);
59395fce210SBarry Smith   ierr = PetscSFSetGraph(*isf,maxlocal,count,newilocal,PETSC_OWN_POINTER,roots,PETSC_COPY_VALUES);CHKERRQ(ierr);
59495fce210SBarry Smith   ierr = PetscFree2(roots,leaves);CHKERRQ(ierr);
59595fce210SBarry Smith   PetscFunctionReturn(0);
59695fce210SBarry Smith }
59795fce210SBarry Smith 
59895fce210SBarry Smith /*@
59995fce210SBarry Smith    PetscSFDuplicate - duplicate a PetscSF, optionally preserving rank connectivity and graph
60095fce210SBarry Smith 
60195fce210SBarry Smith    Collective
60295fce210SBarry Smith 
60395fce210SBarry Smith    Input Arguments:
60495fce210SBarry Smith +  sf - communication object to duplicate
60595fce210SBarry Smith -  opt - PETSCSF_DUPLICATE_CONFONLY, PETSCSF_DUPLICATE_RANKS, or PETSCSF_DUPLICATE_GRAPH (see PetscSFDuplicateOption)
60695fce210SBarry Smith 
60795fce210SBarry Smith    Output Arguments:
60895fce210SBarry Smith .  newsf - new communication object
60995fce210SBarry Smith 
61095fce210SBarry Smith    Level: beginner
61195fce210SBarry Smith 
61295fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetType(), PetscSFSetGraph()
61395fce210SBarry Smith @*/
61495fce210SBarry Smith PetscErrorCode PetscSFDuplicate(PetscSF sf,PetscSFDuplicateOption opt,PetscSF *newsf)
61595fce210SBarry Smith {
61629046d53SLisandro Dalcin   PetscSFType    type;
61795fce210SBarry Smith   PetscErrorCode ierr;
61895fce210SBarry Smith 
61995fce210SBarry Smith   PetscFunctionBegin;
62029046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
62129046d53SLisandro Dalcin   PetscValidLogicalCollectiveEnum(sf,opt,2);
62229046d53SLisandro Dalcin   PetscValidPointer(newsf,3);
62395fce210SBarry Smith   ierr = PetscSFCreate(PetscObjectComm((PetscObject)sf),newsf);CHKERRQ(ierr);
62429046d53SLisandro Dalcin   ierr = PetscSFGetType(sf,&type);CHKERRQ(ierr);
62529046d53SLisandro Dalcin   if (type) {ierr = PetscSFSetType(*newsf,type);CHKERRQ(ierr);}
62695fce210SBarry Smith   if (opt == PETSCSF_DUPLICATE_GRAPH) {
627dd5b3ca6SJunchao Zhang     PetscSFCheckGraphSet(sf,1);
628dd5b3ca6SJunchao Zhang     if (sf->pattern == PETSCSF_PATTERN_GENERAL) {
62995fce210SBarry Smith       PetscInt          nroots,nleaves;
63095fce210SBarry Smith       const PetscInt    *ilocal;
63195fce210SBarry Smith       const PetscSFNode *iremote;
63295fce210SBarry Smith       ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
63395fce210SBarry Smith       ierr = PetscSFSetGraph(*newsf,nroots,nleaves,ilocal,PETSC_COPY_VALUES,iremote,PETSC_COPY_VALUES);CHKERRQ(ierr);
634dd5b3ca6SJunchao Zhang     } else {
635dd5b3ca6SJunchao Zhang       ierr = PetscSFSetGraphWithPattern(*newsf,sf->map,sf->pattern);CHKERRQ(ierr);
636dd5b3ca6SJunchao Zhang     }
63795fce210SBarry Smith   }
63829046d53SLisandro Dalcin   if (sf->ops->Duplicate) {ierr = (*sf->ops->Duplicate)(sf,opt,*newsf);CHKERRQ(ierr);}
63995fce210SBarry Smith   PetscFunctionReturn(0);
64095fce210SBarry Smith }
64195fce210SBarry Smith 
64295fce210SBarry Smith /*@C
64395fce210SBarry Smith    PetscSFGetGraph - Get the graph specifying a parallel star forest
64495fce210SBarry Smith 
64595fce210SBarry Smith    Not Collective
64695fce210SBarry Smith 
64795fce210SBarry Smith    Input Arguments:
64895fce210SBarry Smith .  sf - star forest
64995fce210SBarry Smith 
65095fce210SBarry Smith    Output Arguments:
65195fce210SBarry Smith +  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
65295fce210SBarry Smith .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
65395fce210SBarry Smith .  ilocal - locations of leaves in leafdata buffers
65495fce210SBarry Smith -  iremote - remote locations of root vertices for each leaf on the current process
65595fce210SBarry Smith 
656373e0d91SLisandro Dalcin    Notes:
657373e0d91SLisandro Dalcin    We are not currently requiring that the graph is set, thus returning nroots=-1 if it has not been set yet
658373e0d91SLisandro Dalcin 
659*245d9833Sprj-    When called from Fortran, the returned iremote array is a copy and must be deallocated after use. Consequently, if you
660ca797d7aSLawrence Mitchell    want to update the graph, you must call PetscSFSetGraph after modifying the iremote array.
661ca797d7aSLawrence Mitchell 
66295fce210SBarry Smith    Level: intermediate
66395fce210SBarry Smith 
66495fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph()
66595fce210SBarry Smith @*/
66695fce210SBarry Smith PetscErrorCode PetscSFGetGraph(PetscSF sf,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote)
66795fce210SBarry Smith {
668b8dee149SJunchao Zhang   PetscErrorCode ierr;
66995fce210SBarry Smith 
67095fce210SBarry Smith   PetscFunctionBegin;
67195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
672b8dee149SJunchao Zhang   if (sf->ops->GetGraph) {
673b8dee149SJunchao Zhang     ierr = (sf->ops->GetGraph)(sf,nroots,nleaves,ilocal,iremote);CHKERRQ(ierr);
674b8dee149SJunchao Zhang   } else {
67595fce210SBarry Smith     if (nroots) *nroots = sf->nroots;
67695fce210SBarry Smith     if (nleaves) *nleaves = sf->nleaves;
67795fce210SBarry Smith     if (ilocal) *ilocal = sf->mine;
67895fce210SBarry Smith     if (iremote) *iremote = sf->remote;
679b8dee149SJunchao Zhang   }
68095fce210SBarry Smith   PetscFunctionReturn(0);
68195fce210SBarry Smith }
68295fce210SBarry Smith 
68329046d53SLisandro Dalcin /*@
68495fce210SBarry Smith    PetscSFGetLeafRange - Get the active leaf ranges
68595fce210SBarry Smith 
68695fce210SBarry Smith    Not Collective
68795fce210SBarry Smith 
68895fce210SBarry Smith    Input Arguments:
68995fce210SBarry Smith .  sf - star forest
69095fce210SBarry Smith 
69195fce210SBarry Smith    Output Arguments:
692dd5b3ca6SJunchao Zhang +  minleaf - minimum active leaf on this process. Return 0 if there are no leaves.
693dd5b3ca6SJunchao Zhang -  maxleaf - maximum active leaf on this process. Return -1 if there are no leaves.
69495fce210SBarry Smith 
69595fce210SBarry Smith    Level: developer
69695fce210SBarry Smith 
69795fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph(), PetscSFGetGraph()
69895fce210SBarry Smith @*/
69995fce210SBarry Smith PetscErrorCode PetscSFGetLeafRange(PetscSF sf,PetscInt *minleaf,PetscInt *maxleaf)
70095fce210SBarry Smith {
70195fce210SBarry Smith 
70295fce210SBarry Smith   PetscFunctionBegin;
70395fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
70429046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
70595fce210SBarry Smith   if (minleaf) *minleaf = sf->minleaf;
70695fce210SBarry Smith   if (maxleaf) *maxleaf = sf->maxleaf;
70795fce210SBarry Smith   PetscFunctionReturn(0);
70895fce210SBarry Smith }
70995fce210SBarry Smith 
71095fce210SBarry Smith /*@C
71195fce210SBarry Smith    PetscSFView - view a star forest
71295fce210SBarry Smith 
71395fce210SBarry Smith    Collective
71495fce210SBarry Smith 
71595fce210SBarry Smith    Input Arguments:
71695fce210SBarry Smith +  sf - star forest
71795fce210SBarry Smith -  viewer - viewer to display graph, for example PETSC_VIEWER_STDOUT_WORLD
71895fce210SBarry Smith 
71995fce210SBarry Smith    Level: beginner
72095fce210SBarry Smith 
72195fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph()
72295fce210SBarry Smith @*/
72395fce210SBarry Smith PetscErrorCode PetscSFView(PetscSF sf,PetscViewer viewer)
72495fce210SBarry Smith {
72595fce210SBarry Smith   PetscErrorCode    ierr;
72695fce210SBarry Smith   PetscBool         iascii;
72795fce210SBarry Smith   PetscViewerFormat format;
72895fce210SBarry Smith 
72995fce210SBarry Smith   PetscFunctionBegin;
73095fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
73195fce210SBarry Smith   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sf),&viewer);CHKERRQ(ierr);}
73295fce210SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
73395fce210SBarry Smith   PetscCheckSameComm(sf,1,viewer,2);
73480153354SVaclav Hapla   if (sf->graphset) {ierr = PetscSFSetUp(sf);CHKERRQ(ierr);}
73595fce210SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
73695fce210SBarry Smith   if (iascii) {
73795fce210SBarry Smith     PetscMPIInt rank;
73881bfa7aaSJed Brown     PetscInt    ii,i,j;
73995fce210SBarry Smith 
740dae58748SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)sf,viewer);CHKERRQ(ierr);
74195fce210SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
74295fce210SBarry Smith     if (sf->ops->View) {ierr = (*sf->ops->View)(sf,viewer);CHKERRQ(ierr);}
743dd5b3ca6SJunchao Zhang     if (sf->pattern == PETSCSF_PATTERN_GENERAL) {
74480153354SVaclav Hapla       if (!sf->graphset) {
74580153354SVaclav Hapla         ierr = PetscViewerASCIIPrintf(viewer,"PetscSFSetGraph() has not been called yet\n");CHKERRQ(ierr);
74680153354SVaclav Hapla         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
74780153354SVaclav Hapla         PetscFunctionReturn(0);
74880153354SVaclav Hapla       }
74995fce210SBarry Smith       ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
7501575c14dSBarry Smith       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
75195fce210SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of roots=%D, leaves=%D, remote ranks=%D\n",rank,sf->nroots,sf->nleaves,sf->nranks);CHKERRQ(ierr);
75295fce210SBarry Smith       for (i=0; i<sf->nleaves; i++) {
75395fce210SBarry 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);
75495fce210SBarry Smith       }
75595fce210SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
75695fce210SBarry Smith       ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
75795fce210SBarry Smith       if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
75881bfa7aaSJed Brown         PetscMPIInt *tmpranks,*perm;
75981bfa7aaSJed Brown         ierr = PetscMalloc2(sf->nranks,&tmpranks,sf->nranks,&perm);CHKERRQ(ierr);
760580bdb30SBarry Smith         ierr = PetscArraycpy(tmpranks,sf->ranks,sf->nranks);CHKERRQ(ierr);
76181bfa7aaSJed Brown         for (i=0; i<sf->nranks; i++) perm[i] = i;
76281bfa7aaSJed Brown         ierr = PetscSortMPIIntWithArray(sf->nranks,tmpranks,perm);CHKERRQ(ierr);
76395fce210SBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Roots referenced by my leaves, by rank\n",rank);CHKERRQ(ierr);
76481bfa7aaSJed Brown         for (ii=0; ii<sf->nranks; ii++) {
76581bfa7aaSJed Brown           i = perm[ii];
7667904a332SBarry Smith           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %d: %D edges\n",rank,sf->ranks[i],sf->roffset[i+1]-sf->roffset[i]);CHKERRQ(ierr);
76795fce210SBarry Smith           for (j=sf->roffset[i]; j<sf->roffset[i+1]; j++) {
76895fce210SBarry Smith             ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d]    %D <- %D\n",rank,sf->rmine[j],sf->rremote[j]);CHKERRQ(ierr);
76995fce210SBarry Smith           }
77095fce210SBarry Smith         }
77181bfa7aaSJed Brown         ierr = PetscFree2(tmpranks,perm);CHKERRQ(ierr);
77295fce210SBarry Smith       }
77395fce210SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
7741575c14dSBarry Smith       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
775dd5b3ca6SJunchao Zhang     }
77695fce210SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
77795fce210SBarry Smith   }
77895fce210SBarry Smith   PetscFunctionReturn(0);
77995fce210SBarry Smith }
78095fce210SBarry Smith 
78195fce210SBarry Smith /*@C
782dec1416fSJunchao Zhang    PetscSFGetRootRanks - Get root ranks and number of vertices referenced by leaves on this process
78395fce210SBarry Smith 
78495fce210SBarry Smith    Not Collective
78595fce210SBarry Smith 
78695fce210SBarry Smith    Input Arguments:
78795fce210SBarry Smith .  sf - star forest
78895fce210SBarry Smith 
78995fce210SBarry Smith    Output Arguments:
79095fce210SBarry Smith +  nranks - number of ranks referenced by local part
79195fce210SBarry Smith .  ranks - array of ranks
79295fce210SBarry Smith .  roffset - offset in rmine/rremote for each rank (length nranks+1)
79395fce210SBarry Smith .  rmine - concatenated array holding local indices referencing each remote rank
79495fce210SBarry Smith -  rremote - concatenated array holding remote indices referenced for each remote rank
79595fce210SBarry Smith 
79695fce210SBarry Smith    Level: developer
79795fce210SBarry Smith 
798dec1416fSJunchao Zhang .seealso: PetscSFGetLeafRanks()
79995fce210SBarry Smith @*/
800dec1416fSJunchao Zhang PetscErrorCode PetscSFGetRootRanks(PetscSF sf,PetscInt *nranks,const PetscMPIInt **ranks,const PetscInt **roffset,const PetscInt **rmine,const PetscInt **rremote)
80195fce210SBarry Smith {
802dec1416fSJunchao Zhang   PetscErrorCode ierr;
80395fce210SBarry Smith 
80495fce210SBarry Smith   PetscFunctionBegin;
80595fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
80629046d53SLisandro Dalcin   if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining ranks");
807dec1416fSJunchao Zhang   if (sf->ops->GetRootRanks) {
808dec1416fSJunchao Zhang     ierr = (sf->ops->GetRootRanks)(sf,nranks,ranks,roffset,rmine,rremote);CHKERRQ(ierr);
809dec1416fSJunchao Zhang   } else {
810dec1416fSJunchao Zhang     /* The generic implementation */
81195fce210SBarry Smith     if (nranks)  *nranks  = sf->nranks;
81295fce210SBarry Smith     if (ranks)   *ranks   = sf->ranks;
81395fce210SBarry Smith     if (roffset) *roffset = sf->roffset;
81495fce210SBarry Smith     if (rmine)   *rmine   = sf->rmine;
81595fce210SBarry Smith     if (rremote) *rremote = sf->rremote;
816dec1416fSJunchao Zhang   }
81795fce210SBarry Smith   PetscFunctionReturn(0);
81895fce210SBarry Smith }
81995fce210SBarry Smith 
8208750ddebSJunchao Zhang /*@C
8218750ddebSJunchao Zhang    PetscSFGetLeafRanks - Get leaf ranks referencing roots on this process
8228750ddebSJunchao Zhang 
8238750ddebSJunchao Zhang    Not Collective
8248750ddebSJunchao Zhang 
8258750ddebSJunchao Zhang    Input Arguments:
8268750ddebSJunchao Zhang .  sf - star forest
8278750ddebSJunchao Zhang 
8288750ddebSJunchao Zhang    Output Arguments:
8298750ddebSJunchao Zhang +  niranks - number of leaf ranks referencing roots on this process
8308750ddebSJunchao Zhang .  iranks - array of ranks
8318750ddebSJunchao Zhang .  ioffset - offset in irootloc for each rank (length niranks+1)
8328750ddebSJunchao Zhang -  irootloc - concatenated array holding local indices of roots referenced by each leaf rank
8338750ddebSJunchao Zhang 
8348750ddebSJunchao Zhang    Level: developer
8358750ddebSJunchao Zhang 
836dec1416fSJunchao Zhang .seealso: PetscSFGetRootRanks()
8378750ddebSJunchao Zhang @*/
8388750ddebSJunchao Zhang PetscErrorCode PetscSFGetLeafRanks(PetscSF sf,PetscInt *niranks,const PetscMPIInt **iranks,const PetscInt **ioffset,const PetscInt **irootloc)
8398750ddebSJunchao Zhang {
8408750ddebSJunchao Zhang   PetscErrorCode ierr;
8418750ddebSJunchao Zhang 
8428750ddebSJunchao Zhang   PetscFunctionBegin;
8438750ddebSJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
8448750ddebSJunchao Zhang   if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining ranks");
8458750ddebSJunchao Zhang   if (sf->ops->GetLeafRanks) {
8468750ddebSJunchao Zhang     ierr = (sf->ops->GetLeafRanks)(sf,niranks,iranks,ioffset,irootloc);CHKERRQ(ierr);
8478750ddebSJunchao Zhang   } else {
8488750ddebSJunchao Zhang     PetscSFType type;
8498750ddebSJunchao Zhang     ierr = PetscSFGetType(sf,&type);CHKERRQ(ierr);
8508750ddebSJunchao Zhang     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSFGetLeafRanks() is not supported on this StarForest type: %s", type);
8518750ddebSJunchao Zhang   }
8528750ddebSJunchao Zhang   PetscFunctionReturn(0);
8538750ddebSJunchao Zhang }
8548750ddebSJunchao Zhang 
855b5a8e515SJed Brown static PetscBool InList(PetscMPIInt needle,PetscMPIInt n,const PetscMPIInt *list) {
856b5a8e515SJed Brown   PetscInt i;
857b5a8e515SJed Brown   for (i=0; i<n; i++) {
858b5a8e515SJed Brown     if (needle == list[i]) return PETSC_TRUE;
859b5a8e515SJed Brown   }
860b5a8e515SJed Brown   return PETSC_FALSE;
861b5a8e515SJed Brown }
862b5a8e515SJed Brown 
86395fce210SBarry Smith /*@C
86421c688dcSJed Brown    PetscSFSetUpRanks - Set up data structures associated with ranks; this is for internal use by PetscSF implementations.
86521c688dcSJed Brown 
86621c688dcSJed Brown    Collective
86721c688dcSJed Brown 
86821c688dcSJed Brown    Input Arguments:
869b5a8e515SJed Brown +  sf - PetscSF to set up; PetscSFSetGraph() must have been called
870b5a8e515SJed Brown -  dgroup - MPI_Group of ranks to be distinguished (e.g., for self or shared memory exchange)
87121c688dcSJed Brown 
87221c688dcSJed Brown    Level: developer
87321c688dcSJed Brown 
874dec1416fSJunchao Zhang .seealso: PetscSFGetRootRanks()
87521c688dcSJed Brown @*/
876b5a8e515SJed Brown PetscErrorCode PetscSFSetUpRanks(PetscSF sf,MPI_Group dgroup)
87721c688dcSJed Brown {
87821c688dcSJed Brown   PetscErrorCode     ierr;
87921c688dcSJed Brown   PetscTable         table;
88021c688dcSJed Brown   PetscTablePosition pos;
881b5a8e515SJed Brown   PetscMPIInt        size,groupsize,*groupranks;
882247e8311SStefano Zampini   PetscInt           *rcount,*ranks;
883247e8311SStefano Zampini   PetscInt           i, irank = -1,orank = -1;
88421c688dcSJed Brown 
88521c688dcSJed Brown   PetscFunctionBegin;
88621c688dcSJed Brown   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
88729046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
88821c688dcSJed Brown   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr);
88921c688dcSJed Brown   ierr = PetscTableCreate(10,size,&table);CHKERRQ(ierr);
89021c688dcSJed Brown   for (i=0; i<sf->nleaves; i++) {
89121c688dcSJed Brown     /* Log 1-based rank */
89221c688dcSJed Brown     ierr = PetscTableAdd(table,sf->remote[i].rank+1,1,ADD_VALUES);CHKERRQ(ierr);
89321c688dcSJed Brown   }
89421c688dcSJed Brown   ierr = PetscTableGetCount(table,&sf->nranks);CHKERRQ(ierr);
89521c688dcSJed Brown   ierr = PetscMalloc4(sf->nranks,&sf->ranks,sf->nranks+1,&sf->roffset,sf->nleaves,&sf->rmine,sf->nleaves,&sf->rremote);CHKERRQ(ierr);
89621c688dcSJed Brown   ierr = PetscMalloc2(sf->nranks,&rcount,sf->nranks,&ranks);CHKERRQ(ierr);
89721c688dcSJed Brown   ierr = PetscTableGetHeadPosition(table,&pos);CHKERRQ(ierr);
89821c688dcSJed Brown   for (i=0; i<sf->nranks; i++) {
89921c688dcSJed Brown     ierr = PetscTableGetNext(table,&pos,&ranks[i],&rcount[i]);CHKERRQ(ierr);
90021c688dcSJed Brown     ranks[i]--;             /* Convert back to 0-based */
90121c688dcSJed Brown   }
90221c688dcSJed Brown   ierr = PetscTableDestroy(&table);CHKERRQ(ierr);
903b5a8e515SJed Brown 
904b5a8e515SJed Brown   /* We expect that dgroup is reliably "small" while nranks could be large */
905b5a8e515SJed Brown   {
9067fb8a5e4SKarl Rupp     MPI_Group group = MPI_GROUP_NULL;
907b5a8e515SJed Brown     PetscMPIInt *dgroupranks;
908b5a8e515SJed Brown     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
909b5a8e515SJed Brown     ierr = MPI_Group_size(dgroup,&groupsize);CHKERRQ(ierr);
910b5a8e515SJed Brown     ierr = PetscMalloc1(groupsize,&dgroupranks);CHKERRQ(ierr);
911b5a8e515SJed Brown     ierr = PetscMalloc1(groupsize,&groupranks);CHKERRQ(ierr);
912b5a8e515SJed Brown     for (i=0; i<groupsize; i++) dgroupranks[i] = i;
913f3fc9a17SJed Brown     if (groupsize) {ierr = MPI_Group_translate_ranks(dgroup,groupsize,dgroupranks,group,groupranks);CHKERRQ(ierr);}
914b5a8e515SJed Brown     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
915b5a8e515SJed Brown     ierr = PetscFree(dgroupranks);CHKERRQ(ierr);
916b5a8e515SJed Brown   }
917b5a8e515SJed Brown 
918b5a8e515SJed Brown   /* Partition ranks[] into distinguished (first sf->ndranks) followed by non-distinguished */
919b5a8e515SJed Brown   for (sf->ndranks=0,i=sf->nranks; sf->ndranks<i; ) {
920b5a8e515SJed Brown     for (i--; sf->ndranks<i; i--) { /* Scan i backward looking for distinguished rank */
921b5a8e515SJed Brown       if (InList(ranks[i],groupsize,groupranks)) break;
922b5a8e515SJed Brown     }
923b5a8e515SJed Brown     for ( ; sf->ndranks<=i; sf->ndranks++) { /* Scan sf->ndranks forward looking for non-distinguished rank */
924b5a8e515SJed Brown       if (!InList(ranks[sf->ndranks],groupsize,groupranks)) break;
925b5a8e515SJed Brown     }
926b5a8e515SJed Brown     if (sf->ndranks < i) {                         /* Swap ranks[sf->ndranks] with ranks[i] */
927b5a8e515SJed Brown       PetscInt    tmprank,tmpcount;
928247e8311SStefano Zampini 
929b5a8e515SJed Brown       tmprank             = ranks[i];
930b5a8e515SJed Brown       tmpcount            = rcount[i];
931b5a8e515SJed Brown       ranks[i]            = ranks[sf->ndranks];
932b5a8e515SJed Brown       rcount[i]           = rcount[sf->ndranks];
933b5a8e515SJed Brown       ranks[sf->ndranks]  = tmprank;
934b5a8e515SJed Brown       rcount[sf->ndranks] = tmpcount;
935b5a8e515SJed Brown       sf->ndranks++;
936b5a8e515SJed Brown     }
937b5a8e515SJed Brown   }
938b5a8e515SJed Brown   ierr = PetscFree(groupranks);CHKERRQ(ierr);
939b5a8e515SJed Brown   ierr = PetscSortIntWithArray(sf->ndranks,ranks,rcount);CHKERRQ(ierr);
940b5a8e515SJed Brown   ierr = PetscSortIntWithArray(sf->nranks-sf->ndranks,ranks+sf->ndranks,rcount+sf->ndranks);CHKERRQ(ierr);
94121c688dcSJed Brown   sf->roffset[0] = 0;
94221c688dcSJed Brown   for (i=0; i<sf->nranks; i++) {
94321c688dcSJed Brown     ierr = PetscMPIIntCast(ranks[i],sf->ranks+i);CHKERRQ(ierr);
94421c688dcSJed Brown     sf->roffset[i+1] = sf->roffset[i] + rcount[i];
94521c688dcSJed Brown     rcount[i]        = 0;
94621c688dcSJed Brown   }
947247e8311SStefano Zampini   for (i=0, irank = -1, orank = -1; i<sf->nleaves; i++) {
948247e8311SStefano Zampini     /* short circuit */
949247e8311SStefano Zampini     if (orank != sf->remote[i].rank) {
95021c688dcSJed Brown       /* Search for index of iremote[i].rank in sf->ranks */
951b5a8e515SJed Brown       ierr = PetscFindMPIInt(sf->remote[i].rank,sf->ndranks,sf->ranks,&irank);CHKERRQ(ierr);
952b5a8e515SJed Brown       if (irank < 0) {
953b5a8e515SJed Brown         ierr = PetscFindMPIInt(sf->remote[i].rank,sf->nranks-sf->ndranks,sf->ranks+sf->ndranks,&irank);CHKERRQ(ierr);
954b5a8e515SJed Brown         if (irank >= 0) irank += sf->ndranks;
95521c688dcSJed Brown       }
956247e8311SStefano Zampini       orank = sf->remote[i].rank;
957247e8311SStefano Zampini     }
958b5a8e515SJed Brown     if (irank < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Could not find rank %D in array",sf->remote[i].rank);
95921c688dcSJed Brown     sf->rmine[sf->roffset[irank] + rcount[irank]]   = sf->mine ? sf->mine[i] : i;
96021c688dcSJed Brown     sf->rremote[sf->roffset[irank] + rcount[irank]] = sf->remote[i].index;
96121c688dcSJed Brown     rcount[irank]++;
96221c688dcSJed Brown   }
96321c688dcSJed Brown   ierr = PetscFree2(rcount,ranks);CHKERRQ(ierr);
96421c688dcSJed Brown   PetscFunctionReturn(0);
96521c688dcSJed Brown }
96621c688dcSJed Brown 
96721c688dcSJed Brown /*@C
96895fce210SBarry Smith    PetscSFGetGroups - gets incoming and outgoing process groups
96995fce210SBarry Smith 
97095fce210SBarry Smith    Collective
97195fce210SBarry Smith 
97295fce210SBarry Smith    Input Argument:
97395fce210SBarry Smith .  sf - star forest
97495fce210SBarry Smith 
97595fce210SBarry Smith    Output Arguments:
97695fce210SBarry Smith +  incoming - group of origin processes for incoming edges (leaves that reference my roots)
97795fce210SBarry Smith -  outgoing - group of destination processes for outgoing edges (roots that I reference)
97895fce210SBarry Smith 
97995fce210SBarry Smith    Level: developer
98095fce210SBarry Smith 
98195fce210SBarry Smith .seealso: PetscSFGetWindow(), PetscSFRestoreWindow()
98295fce210SBarry Smith @*/
98395fce210SBarry Smith PetscErrorCode PetscSFGetGroups(PetscSF sf,MPI_Group *incoming,MPI_Group *outgoing)
98495fce210SBarry Smith {
98595fce210SBarry Smith   PetscErrorCode ierr;
9867fb8a5e4SKarl Rupp   MPI_Group      group = MPI_GROUP_NULL;
98795fce210SBarry Smith 
98895fce210SBarry Smith   PetscFunctionBegin;
98929046d53SLisandro Dalcin   if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining groups");
99095fce210SBarry Smith   if (sf->ingroup == MPI_GROUP_NULL) {
99195fce210SBarry Smith     PetscInt       i;
99295fce210SBarry Smith     const PetscInt *indegree;
99395fce210SBarry Smith     PetscMPIInt    rank,*outranks,*inranks;
99495fce210SBarry Smith     PetscSFNode    *remote;
99595fce210SBarry Smith     PetscSF        bgcount;
99695fce210SBarry Smith 
99795fce210SBarry Smith     /* Compute the number of incoming ranks */
998785e854fSJed Brown     ierr = PetscMalloc1(sf->nranks,&remote);CHKERRQ(ierr);
99995fce210SBarry Smith     for (i=0; i<sf->nranks; i++) {
100095fce210SBarry Smith       remote[i].rank  = sf->ranks[i];
100195fce210SBarry Smith       remote[i].index = 0;
100295fce210SBarry Smith     }
100395fce210SBarry Smith     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,&bgcount);CHKERRQ(ierr);
100495fce210SBarry Smith     ierr = PetscSFSetGraph(bgcount,1,sf->nranks,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
100595fce210SBarry Smith     ierr = PetscSFComputeDegreeBegin(bgcount,&indegree);CHKERRQ(ierr);
100695fce210SBarry Smith     ierr = PetscSFComputeDegreeEnd(bgcount,&indegree);CHKERRQ(ierr);
100795fce210SBarry Smith 
100895fce210SBarry Smith     /* Enumerate the incoming ranks */
1009dcca6d9dSJed Brown     ierr = PetscMalloc2(indegree[0],&inranks,sf->nranks,&outranks);CHKERRQ(ierr);
101095fce210SBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
101195fce210SBarry Smith     for (i=0; i<sf->nranks; i++) outranks[i] = rank;
101295fce210SBarry Smith     ierr = PetscSFGatherBegin(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr);
101395fce210SBarry Smith     ierr = PetscSFGatherEnd(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr);
101495fce210SBarry Smith     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
101595fce210SBarry Smith     ierr = MPI_Group_incl(group,indegree[0],inranks,&sf->ingroup);CHKERRQ(ierr);
101695fce210SBarry Smith     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
101795fce210SBarry Smith     ierr = PetscFree2(inranks,outranks);CHKERRQ(ierr);
101895fce210SBarry Smith     ierr = PetscSFDestroy(&bgcount);CHKERRQ(ierr);
101995fce210SBarry Smith   }
102095fce210SBarry Smith   *incoming = sf->ingroup;
102195fce210SBarry Smith 
102295fce210SBarry Smith   if (sf->outgroup == MPI_GROUP_NULL) {
102395fce210SBarry Smith     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
102495fce210SBarry Smith     ierr = MPI_Group_incl(group,sf->nranks,sf->ranks,&sf->outgroup);CHKERRQ(ierr);
102595fce210SBarry Smith     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
102695fce210SBarry Smith   }
102795fce210SBarry Smith   *outgoing = sf->outgroup;
102895fce210SBarry Smith   PetscFunctionReturn(0);
102995fce210SBarry Smith }
103095fce210SBarry Smith 
103129046d53SLisandro Dalcin /*@
103295fce210SBarry Smith    PetscSFGetMultiSF - gets the inner SF implemeting gathers and scatters
103395fce210SBarry Smith 
103495fce210SBarry Smith    Collective
103595fce210SBarry Smith 
103695fce210SBarry Smith    Input Argument:
103795fce210SBarry Smith .  sf - star forest that may contain roots with 0 or with more than 1 vertex
103895fce210SBarry Smith 
103995fce210SBarry Smith    Output Arguments:
104095fce210SBarry Smith .  multi - star forest with split roots, such that each root has degree exactly 1
104195fce210SBarry Smith 
104295fce210SBarry Smith    Level: developer
104395fce210SBarry Smith 
104495fce210SBarry Smith    Notes:
104595fce210SBarry Smith 
104695fce210SBarry Smith    In most cases, users should use PetscSFGatherBegin() and PetscSFScatterBegin() instead of manipulating multi
104795fce210SBarry Smith    directly. Since multi satisfies the stronger condition that each entry in the global space has exactly one incoming
104895fce210SBarry Smith    edge, it is a candidate for future optimization that might involve its removal.
104995fce210SBarry Smith 
1050673100f5SVaclav Hapla .seealso: PetscSFSetGraph(), PetscSFGatherBegin(), PetscSFScatterBegin(), PetscSFComputeMultiRootOriginalNumbering()
105195fce210SBarry Smith @*/
105295fce210SBarry Smith PetscErrorCode PetscSFGetMultiSF(PetscSF sf,PetscSF *multi)
105395fce210SBarry Smith {
105495fce210SBarry Smith   PetscErrorCode ierr;
105595fce210SBarry Smith 
105695fce210SBarry Smith   PetscFunctionBegin;
105795fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
105895fce210SBarry Smith   PetscValidPointer(multi,2);
105995fce210SBarry Smith   if (sf->nroots < 0) {         /* Graph has not been set yet; why do we need this? */
106095fce210SBarry Smith     ierr   = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr);
106195fce210SBarry Smith     *multi = sf->multi;
106295fce210SBarry Smith     PetscFunctionReturn(0);
106395fce210SBarry Smith   }
106495fce210SBarry Smith   if (!sf->multi) {
106595fce210SBarry Smith     const PetscInt *indegree;
10669837ea96SMatthew G. Knepley     PetscInt       i,*inoffset,*outones,*outoffset,maxlocal;
106795fce210SBarry Smith     PetscSFNode    *remote;
106829046d53SLisandro Dalcin     maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */
106995fce210SBarry Smith     ierr = PetscSFComputeDegreeBegin(sf,&indegree);CHKERRQ(ierr);
107095fce210SBarry Smith     ierr = PetscSFComputeDegreeEnd(sf,&indegree);CHKERRQ(ierr);
10719837ea96SMatthew G. Knepley     ierr = PetscMalloc3(sf->nroots+1,&inoffset,maxlocal,&outones,maxlocal,&outoffset);CHKERRQ(ierr);
107295fce210SBarry Smith     inoffset[0] = 0;
107395fce210SBarry Smith     for (i=0; i<sf->nroots; i++) inoffset[i+1] = inoffset[i] + indegree[i];
10749837ea96SMatthew G. Knepley     for (i=0; i<maxlocal; i++) outones[i] = 1;
1075dbd2ff41SBarry Smith     ierr = PetscSFFetchAndOpBegin(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr);
1076dbd2ff41SBarry Smith     ierr = PetscSFFetchAndOpEnd(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr);
107795fce210SBarry Smith     for (i=0; i<sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */
107895fce210SBarry Smith #if 0
107995fce210SBarry Smith #if defined(PETSC_USE_DEBUG)                                 /* Check that the expected number of increments occurred */
108095fce210SBarry Smith     for (i=0; i<sf->nroots; i++) {
108195fce210SBarry Smith       if (inoffset[i] + indegree[i] != inoffset[i+1]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect result after PetscSFFetchAndOp");
108295fce210SBarry Smith     }
108395fce210SBarry Smith #endif
108495fce210SBarry Smith #endif
1085785e854fSJed Brown     ierr = PetscMalloc1(sf->nleaves,&remote);CHKERRQ(ierr);
108695fce210SBarry Smith     for (i=0; i<sf->nleaves; i++) {
108795fce210SBarry Smith       remote[i].rank  = sf->remote[i].rank;
108838e7336fSToby Isaac       remote[i].index = outoffset[sf->mine ? sf->mine[i] : i];
108995fce210SBarry Smith     }
109095fce210SBarry Smith     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr);
109101365b40SToby Isaac     ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
109295fce210SBarry Smith     if (sf->rankorder) {        /* Sort the ranks */
109395fce210SBarry Smith       PetscMPIInt rank;
109495fce210SBarry Smith       PetscInt    *inranks,*newoffset,*outranks,*newoutoffset,*tmpoffset,maxdegree;
109595fce210SBarry Smith       PetscSFNode *newremote;
109695fce210SBarry Smith       ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
109795fce210SBarry Smith       for (i=0,maxdegree=0; i<sf->nroots; i++) maxdegree = PetscMax(maxdegree,indegree[i]);
10989837ea96SMatthew G. Knepley       ierr = PetscMalloc5(sf->multi->nroots,&inranks,sf->multi->nroots,&newoffset,maxlocal,&outranks,maxlocal,&newoutoffset,maxdegree,&tmpoffset);CHKERRQ(ierr);
10999837ea96SMatthew G. Knepley       for (i=0; i<maxlocal; i++) outranks[i] = rank;
11008bfbc91cSJed Brown       ierr = PetscSFReduceBegin(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr);
11018bfbc91cSJed Brown       ierr = PetscSFReduceEnd(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr);
110295fce210SBarry Smith       /* Sort the incoming ranks at each vertex, build the inverse map */
110395fce210SBarry Smith       for (i=0; i<sf->nroots; i++) {
110495fce210SBarry Smith         PetscInt j;
110595fce210SBarry Smith         for (j=0; j<indegree[i]; j++) tmpoffset[j] = j;
110695fce210SBarry Smith         ierr = PetscSortIntWithArray(indegree[i],inranks+inoffset[i],tmpoffset);CHKERRQ(ierr);
110795fce210SBarry Smith         for (j=0; j<indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j;
110895fce210SBarry Smith       }
110995fce210SBarry Smith       ierr = PetscSFBcastBegin(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr);
111095fce210SBarry Smith       ierr = PetscSFBcastEnd(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr);
1111785e854fSJed Brown       ierr = PetscMalloc1(sf->nleaves,&newremote);CHKERRQ(ierr);
111295fce210SBarry Smith       for (i=0; i<sf->nleaves; i++) {
111395fce210SBarry Smith         newremote[i].rank  = sf->remote[i].rank;
111401365b40SToby Isaac         newremote[i].index = newoutoffset[sf->mine ? sf->mine[i] : i];
111595fce210SBarry Smith       }
111601365b40SToby Isaac       ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,newremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
111795fce210SBarry Smith       ierr = PetscFree5(inranks,newoffset,outranks,newoutoffset,tmpoffset);CHKERRQ(ierr);
111895fce210SBarry Smith     }
111995fce210SBarry Smith     ierr = PetscFree3(inoffset,outones,outoffset);CHKERRQ(ierr);
112095fce210SBarry Smith   }
112195fce210SBarry Smith   *multi = sf->multi;
112295fce210SBarry Smith   PetscFunctionReturn(0);
112395fce210SBarry Smith }
112495fce210SBarry Smith 
112595fce210SBarry Smith /*@C
112695fce210SBarry Smith    PetscSFCreateEmbeddedSF - removes edges from all but the selected roots, does not remap indices
112795fce210SBarry Smith 
112895fce210SBarry Smith    Collective
112995fce210SBarry Smith 
113095fce210SBarry Smith    Input Arguments:
113195fce210SBarry Smith +  sf - original star forest
1132ba2a7774SJunchao Zhang .  nselected  - number of selected roots on this process
1133ba2a7774SJunchao Zhang -  selected   - indices of the selected roots on this process
113495fce210SBarry Smith 
113595fce210SBarry Smith    Output Arguments:
113695fce210SBarry Smith .  newsf - new star forest
113795fce210SBarry Smith 
113895fce210SBarry Smith    Level: advanced
113995fce210SBarry Smith 
114095fce210SBarry Smith    Note:
114195fce210SBarry Smith    To use the new PetscSF, it may be necessary to know the indices of the leaves that are still participating. This can
114295fce210SBarry Smith    be done by calling PetscSFGetGraph().
114395fce210SBarry Smith 
114495fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFGetGraph()
114595fce210SBarry Smith @*/
1146ba2a7774SJunchao Zhang PetscErrorCode PetscSFCreateEmbeddedSF(PetscSF sf,PetscInt nselected,const PetscInt *selected,PetscSF *newsf)
114795fce210SBarry Smith {
1148f659e5c7SJunchao Zhang   PetscInt          i,n,*roots,*rootdata,*leafdata,nroots,nleaves,connected_leaves,*new_ilocal;
1149ba2a7774SJunchao Zhang   const PetscSFNode *iremote;
1150f659e5c7SJunchao Zhang   PetscSFNode       *new_iremote;
1151ba2a7774SJunchao Zhang   PetscSF           tmpsf;
1152f659e5c7SJunchao Zhang   MPI_Comm          comm;
11530511a646SMatthew G. Knepley   PetscErrorCode    ierr;
115495fce210SBarry Smith 
115595fce210SBarry Smith   PetscFunctionBegin;
115695fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
115729046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
1158ba2a7774SJunchao Zhang   if (nselected) PetscValidPointer(selected,3);
115995fce210SBarry Smith   PetscValidPointer(newsf,4);
11600511a646SMatthew G. Knepley 
1161f659e5c7SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1162f659e5c7SJunchao Zhang 
1163f659e5c7SJunchao Zhang   /* Uniq selected[] and put results in roots[] */
1164f659e5c7SJunchao Zhang   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
1165f659e5c7SJunchao Zhang   ierr = PetscMalloc1(nselected,&roots);CHKERRQ(ierr);
1166dd5b3ca6SJunchao Zhang   ierr = PetscArraycpy(roots,selected,nselected);CHKERRQ(ierr);
1167f659e5c7SJunchao Zhang   ierr = PetscSortedRemoveDupsInt(&nselected,roots);CHKERRQ(ierr);
1168f659e5c7SJunchao 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);
1169f659e5c7SJunchao Zhang 
1170f659e5c7SJunchao Zhang   if (sf->ops->CreateEmbeddedSF) {
1171f659e5c7SJunchao Zhang     ierr = (*sf->ops->CreateEmbeddedSF)(sf,nselected,roots,newsf);CHKERRQ(ierr);
1172f659e5c7SJunchao Zhang   } else {
1173f659e5c7SJunchao Zhang     /* A generic version of creating embedded sf. Note that we called PetscSFSetGraph() twice, which is certainly expensive */
1174ba2a7774SJunchao Zhang     /* Find out which leaves (not leaf data items) are still connected to roots in the embedded sf */
1175f659e5c7SJunchao Zhang     ierr = PetscSFGetGraph(sf,&nroots,&nleaves,NULL,&iremote);CHKERRQ(ierr);
1176ba2a7774SJunchao Zhang     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&tmpsf);CHKERRQ(ierr);
1177ba2a7774SJunchao Zhang     ierr = PetscSFSetGraph(tmpsf,nroots,nleaves,NULL/*contiguous*/,PETSC_USE_POINTER,iremote,PETSC_USE_POINTER);CHKERRQ(ierr);
1178ba2a7774SJunchao Zhang     ierr = PetscCalloc2(nroots,&rootdata,nleaves,&leafdata);CHKERRQ(ierr);
1179f659e5c7SJunchao Zhang     for (i=0; i<nselected; ++i) rootdata[roots[i]] = 1;
1180ba2a7774SJunchao Zhang     ierr = PetscSFBcastBegin(tmpsf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr);
1181ba2a7774SJunchao Zhang     ierr = PetscSFBcastEnd(tmpsf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr);
1182ba2a7774SJunchao Zhang     ierr = PetscSFDestroy(&tmpsf);CHKERRQ(ierr);
1183ba2a7774SJunchao Zhang 
1184ba2a7774SJunchao Zhang     /* Build newsf with leaves that are still connected */
1185f659e5c7SJunchao Zhang     connected_leaves = 0;
1186f659e5c7SJunchao Zhang     for (i=0; i<nleaves; ++i) connected_leaves += leafdata[i];
1187f659e5c7SJunchao Zhang     ierr = PetscMalloc1(connected_leaves,&new_ilocal);CHKERRQ(ierr);
1188f659e5c7SJunchao Zhang     ierr = PetscMalloc1(connected_leaves,&new_iremote);CHKERRQ(ierr);
1189ba2a7774SJunchao Zhang     for (i=0, n=0; i<nleaves; ++i) {
1190ba2a7774SJunchao Zhang       if (leafdata[i]) {
1191ba2a7774SJunchao Zhang         new_ilocal[n]        = sf->mine ? sf->mine[i] : i;
1192ba2a7774SJunchao Zhang         new_iremote[n].rank  = sf->remote[i].rank;
1193ba2a7774SJunchao Zhang         new_iremote[n].index = sf->remote[i].index;
1194fc1ede2bSMatthew G. Knepley         ++n;
119595fce210SBarry Smith       }
119695fce210SBarry Smith     }
1197f659e5c7SJunchao Zhang 
1198f659e5c7SJunchao 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);
1199f659e5c7SJunchao Zhang     ierr = PetscSFCreate(comm,newsf);CHKERRQ(ierr);
1200f659e5c7SJunchao Zhang     ierr = PetscSFSetFromOptions(*newsf);CHKERRQ(ierr);
1201f659e5c7SJunchao Zhang     ierr = PetscSFSetGraph(*newsf,nroots,connected_leaves,new_ilocal,PETSC_OWN_POINTER,new_iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
120295fce210SBarry Smith     ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr);
1203f659e5c7SJunchao Zhang   }
1204f659e5c7SJunchao Zhang   ierr = PetscFree(roots);CHKERRQ(ierr);
120595fce210SBarry Smith   PetscFunctionReturn(0);
120695fce210SBarry Smith }
120795fce210SBarry Smith 
12082f5fb4c2SMatthew G. Knepley /*@C
12092f5fb4c2SMatthew G. Knepley   PetscSFCreateEmbeddedLeafSF - removes edges from all but the selected leaves, does not remap indices
12102f5fb4c2SMatthew G. Knepley 
12112f5fb4c2SMatthew G. Knepley   Collective
12122f5fb4c2SMatthew G. Knepley 
12132f5fb4c2SMatthew G. Knepley   Input Arguments:
12142f5fb4c2SMatthew G. Knepley + sf - original star forest
1215f659e5c7SJunchao Zhang . nselected  - number of selected leaves on this process
1216f659e5c7SJunchao Zhang - selected   - indices of the selected leaves on this process
12172f5fb4c2SMatthew G. Knepley 
12182f5fb4c2SMatthew G. Knepley   Output Arguments:
12192f5fb4c2SMatthew G. Knepley .  newsf - new star forest
12202f5fb4c2SMatthew G. Knepley 
12212f5fb4c2SMatthew G. Knepley   Level: advanced
12222f5fb4c2SMatthew G. Knepley 
12232f5fb4c2SMatthew G. Knepley .seealso: PetscSFCreateEmbeddedSF(), PetscSFSetGraph(), PetscSFGetGraph()
12242f5fb4c2SMatthew G. Knepley @*/
1225f659e5c7SJunchao Zhang PetscErrorCode PetscSFCreateEmbeddedLeafSF(PetscSF sf,PetscInt nselected,const PetscInt *selected,PetscSF *newsf)
12262f5fb4c2SMatthew G. Knepley {
1227f659e5c7SJunchao Zhang   const PetscSFNode *iremote;
1228f659e5c7SJunchao Zhang   PetscSFNode       *new_iremote;
1229f659e5c7SJunchao Zhang   const PetscInt    *ilocal;
1230f659e5c7SJunchao Zhang   PetscInt          i,nroots,*leaves,*new_ilocal;
1231f659e5c7SJunchao Zhang   MPI_Comm          comm;
12322f5fb4c2SMatthew G. Knepley   PetscErrorCode    ierr;
12332f5fb4c2SMatthew G. Knepley 
12342f5fb4c2SMatthew G. Knepley   PetscFunctionBegin;
12352f5fb4c2SMatthew G. Knepley   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
123629046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
1237f659e5c7SJunchao Zhang   if (nselected) PetscValidPointer(selected,3);
12382f5fb4c2SMatthew G. Knepley   PetscValidPointer(newsf,4);
12392f5fb4c2SMatthew G. Knepley 
1240f659e5c7SJunchao Zhang   /* Uniq selected[] and put results in leaves[] */
1241f659e5c7SJunchao Zhang   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
1242f659e5c7SJunchao Zhang   ierr = PetscMalloc1(nselected,&leaves);CHKERRQ(ierr);
1243dd5b3ca6SJunchao Zhang   ierr = PetscArraycpy(leaves,selected,nselected);CHKERRQ(ierr);
1244f659e5c7SJunchao Zhang   ierr = PetscSortedRemoveDupsInt(&nselected,leaves);CHKERRQ(ierr);
1245f659e5c7SJunchao 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);
1246f659e5c7SJunchao Zhang 
1247f659e5c7SJunchao Zhang   /* Optimize the routine only when sf is setup and hence we can reuse sf's communication pattern */
1248f659e5c7SJunchao Zhang   if (sf->setupcalled && sf->ops->CreateEmbeddedLeafSF) {
1249f659e5c7SJunchao Zhang     ierr = (*sf->ops->CreateEmbeddedLeafSF)(sf,nselected,leaves,newsf);CHKERRQ(ierr);
1250f659e5c7SJunchao Zhang   } else {
1251f659e5c7SJunchao Zhang     ierr = PetscSFGetGraph(sf,&nroots,NULL,&ilocal,&iremote);CHKERRQ(ierr);
1252f659e5c7SJunchao Zhang     ierr = PetscMalloc1(nselected,&new_ilocal);CHKERRQ(ierr);
1253f659e5c7SJunchao Zhang     ierr = PetscMalloc1(nselected,&new_iremote);CHKERRQ(ierr);
1254f659e5c7SJunchao Zhang     for (i=0; i<nselected; ++i) {
1255f659e5c7SJunchao Zhang       const PetscInt l     = leaves[i];
1256f659e5c7SJunchao Zhang       new_ilocal[i]        = ilocal ? ilocal[l] : l;
1257f659e5c7SJunchao Zhang       new_iremote[i].rank  = iremote[l].rank;
1258f659e5c7SJunchao Zhang       new_iremote[i].index = iremote[l].index;
12592f5fb4c2SMatthew G. Knepley     }
1260f659e5c7SJunchao Zhang     ierr = PetscSFCreate(comm,newsf);CHKERRQ(ierr);
1261f659e5c7SJunchao Zhang     ierr = PetscSFSetFromOptions(*newsf);CHKERRQ(ierr);
1262f659e5c7SJunchao Zhang     ierr = PetscSFSetGraph(*newsf,nroots,nselected,new_ilocal,PETSC_OWN_POINTER,new_iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1263f659e5c7SJunchao Zhang   }
1264f659e5c7SJunchao Zhang   ierr = PetscFree(leaves);CHKERRQ(ierr);
12652f5fb4c2SMatthew G. Knepley   PetscFunctionReturn(0);
12662f5fb4c2SMatthew G. Knepley }
12672f5fb4c2SMatthew G. Knepley 
126895fce210SBarry Smith /*@C
12693482bfa8SJunchao Zhang    PetscSFBcastAndOpBegin - begin pointwise broadcast with root value being reduced to leaf value, to be concluded with call to PetscSFBcastAndOpEnd()
12703482bfa8SJunchao Zhang 
12713482bfa8SJunchao Zhang    Collective on PetscSF
12723482bfa8SJunchao Zhang 
12733482bfa8SJunchao Zhang    Input Arguments:
12743482bfa8SJunchao Zhang +  sf - star forest on which to communicate
12753482bfa8SJunchao Zhang .  unit - data type associated with each node
12763482bfa8SJunchao Zhang .  rootdata - buffer to broadcast
12773482bfa8SJunchao Zhang -  op - operation to use for reduction
12783482bfa8SJunchao Zhang 
12793482bfa8SJunchao Zhang    Output Arguments:
12803482bfa8SJunchao Zhang .  leafdata - buffer to be reduced with values from each leaf's respective root
12813482bfa8SJunchao Zhang 
12823482bfa8SJunchao Zhang    Level: intermediate
12833482bfa8SJunchao Zhang 
12843482bfa8SJunchao Zhang .seealso: PetscSFBcastAndOpEnd(), PetscSFBcastBegin(), PetscSFBcastEnd()
12853482bfa8SJunchao Zhang @*/
12863482bfa8SJunchao Zhang PetscErrorCode PetscSFBcastAndOpBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
12873482bfa8SJunchao Zhang {
12883482bfa8SJunchao Zhang   PetscErrorCode ierr;
1289eb02082bSJunchao Zhang   PetscMemType   rootmtype,leafmtype;
12903482bfa8SJunchao Zhang 
12913482bfa8SJunchao Zhang   PetscFunctionBegin;
12923482bfa8SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
12933482bfa8SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
12943482bfa8SJunchao Zhang   ierr = PetscLogEventBegin(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
1295eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
1296eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
1297eb02082bSJunchao Zhang #if defined(PETSC_HAVE_CUDA)
1298eb02082bSJunchao Zhang   /*  Shall we assume rootdata, leafdata are ready to use (instead of being computed by some asynchronous kernels)?
1299eb02082bSJunchao Zhang     To be similar to MPI, I'd like to have this assumption, since MPI does not have a concept of stream.
1300eb02082bSJunchao Zhang     But currently this assumption is not enforecd in Petsc. To be safe, I do synchronization here. Otherwise, if
1301eb02082bSJunchao Zhang     we do not sync now and call the Pack kernel directly on the default NULL stream (assume petsc objects are also
1302eb02082bSJunchao Zhang     computed on it), we have to sync the NULL stream before calling MPI routines. So, it looks a cudaDeviceSynchronize
1303eb02082bSJunchao Zhang     is inevitable. We do it now and put pack/unpack kernels to non-NULL streams.
1304eb02082bSJunchao Zhang    */
1305eb02082bSJunchao Zhang   if (rootmtype == PETSC_MEMTYPE_DEVICE || leafmtype == PETSC_MEMTYPE_DEVICE) {cudaError_t err = cudaDeviceSynchronize();CHKERRCUDA(err);}
1306eb02082bSJunchao Zhang #endif
1307eb02082bSJunchao Zhang   ierr = (*sf->ops->BcastAndOpBegin)(sf,unit,rootmtype,rootdata,leafmtype,leafdata,op);CHKERRQ(ierr);
13083482bfa8SJunchao Zhang   ierr = PetscLogEventEnd(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
13093482bfa8SJunchao Zhang   PetscFunctionReturn(0);
13103482bfa8SJunchao Zhang }
13113482bfa8SJunchao Zhang 
13123482bfa8SJunchao Zhang /*@C
13133482bfa8SJunchao Zhang    PetscSFBcastAndOpEnd - end a broadcast & reduce operation started with PetscSFBcastAndOpBegin()
13143482bfa8SJunchao Zhang 
13153482bfa8SJunchao Zhang    Collective
13163482bfa8SJunchao Zhang 
13173482bfa8SJunchao Zhang    Input Arguments:
13183482bfa8SJunchao Zhang +  sf - star forest
13193482bfa8SJunchao Zhang .  unit - data type
13203482bfa8SJunchao Zhang .  rootdata - buffer to broadcast
13213482bfa8SJunchao Zhang -  op - operation to use for reduction
13223482bfa8SJunchao Zhang 
13233482bfa8SJunchao Zhang    Output Arguments:
13243482bfa8SJunchao Zhang .  leafdata - buffer to be reduced with values from each leaf's respective root
13253482bfa8SJunchao Zhang 
13263482bfa8SJunchao Zhang    Level: intermediate
13273482bfa8SJunchao Zhang 
13283482bfa8SJunchao Zhang .seealso: PetscSFSetGraph(), PetscSFReduceEnd()
13293482bfa8SJunchao Zhang @*/
13303482bfa8SJunchao Zhang PetscErrorCode PetscSFBcastAndOpEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
13313482bfa8SJunchao Zhang {
13323482bfa8SJunchao Zhang   PetscErrorCode ierr;
1333eb02082bSJunchao Zhang   PetscMemType   rootmtype,leafmtype;
13343482bfa8SJunchao Zhang 
13353482bfa8SJunchao Zhang   PetscFunctionBegin;
13363482bfa8SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
13373482bfa8SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
13383482bfa8SJunchao Zhang   ierr = PetscLogEventBegin(PETSCSF_BcastAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
1339eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
1340eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
1341eb02082bSJunchao Zhang   ierr = (*sf->ops->BcastAndOpEnd)(sf,unit,rootmtype,rootdata,leafmtype,leafdata,op);CHKERRQ(ierr);
13423482bfa8SJunchao Zhang   ierr = PetscLogEventEnd(PETSCSF_BcastAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
13433482bfa8SJunchao Zhang   PetscFunctionReturn(0);
13443482bfa8SJunchao Zhang }
13453482bfa8SJunchao Zhang 
13463482bfa8SJunchao Zhang /*@C
134795fce210SBarry Smith    PetscSFReduceBegin - begin reduction of leafdata into rootdata, to be completed with call to PetscSFReduceEnd()
134895fce210SBarry Smith 
134995fce210SBarry Smith    Collective
135095fce210SBarry Smith 
135195fce210SBarry Smith    Input Arguments:
135295fce210SBarry Smith +  sf - star forest
135395fce210SBarry Smith .  unit - data type
135495fce210SBarry Smith .  leafdata - values to reduce
135595fce210SBarry Smith -  op - reduction operation
135695fce210SBarry Smith 
135795fce210SBarry Smith    Output Arguments:
135895fce210SBarry Smith .  rootdata - result of reduction of values from all leaves of each root
135995fce210SBarry Smith 
136095fce210SBarry Smith    Level: intermediate
136195fce210SBarry Smith 
136295fce210SBarry Smith .seealso: PetscSFBcastBegin()
136395fce210SBarry Smith @*/
136495fce210SBarry Smith PetscErrorCode PetscSFReduceBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
136595fce210SBarry Smith {
136695fce210SBarry Smith   PetscErrorCode ierr;
1367eb02082bSJunchao Zhang   PetscMemType   rootmtype,leafmtype;
136895fce210SBarry Smith 
136995fce210SBarry Smith   PetscFunctionBegin;
137095fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
137195fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
137229046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);
1373eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
1374eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
1375eb02082bSJunchao Zhang #if defined(PETSC_HAVE_CUDA)
1376eb02082bSJunchao Zhang   if (rootmtype == PETSC_MEMTYPE_DEVICE || leafmtype == PETSC_MEMTYPE_DEVICE) {cudaError_t err = cudaDeviceSynchronize();CHKERRCUDA(err);}
1377eb02082bSJunchao Zhang #endif
1378eb02082bSJunchao Zhang   ierr = (sf->ops->ReduceBegin)(sf,unit,leafmtype,leafdata,rootmtype,rootdata,op);CHKERRQ(ierr);
1379563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);
138095fce210SBarry Smith   PetscFunctionReturn(0);
138195fce210SBarry Smith }
138295fce210SBarry Smith 
138395fce210SBarry Smith /*@C
138495fce210SBarry Smith    PetscSFReduceEnd - end a reduction operation started with PetscSFReduceBegin()
138595fce210SBarry Smith 
138695fce210SBarry Smith    Collective
138795fce210SBarry Smith 
138895fce210SBarry Smith    Input Arguments:
138995fce210SBarry Smith +  sf - star forest
139095fce210SBarry Smith .  unit - data type
139195fce210SBarry Smith .  leafdata - values to reduce
139295fce210SBarry Smith -  op - reduction operation
139395fce210SBarry Smith 
139495fce210SBarry Smith    Output Arguments:
139595fce210SBarry Smith .  rootdata - result of reduction of values from all leaves of each root
139695fce210SBarry Smith 
139795fce210SBarry Smith    Level: intermediate
139895fce210SBarry Smith 
139995fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFBcastEnd()
140095fce210SBarry Smith @*/
140195fce210SBarry Smith PetscErrorCode PetscSFReduceEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
140295fce210SBarry Smith {
140395fce210SBarry Smith   PetscErrorCode ierr;
1404eb02082bSJunchao Zhang   PetscMemType   rootmtype,leafmtype;
140595fce210SBarry Smith 
140695fce210SBarry Smith   PetscFunctionBegin;
140795fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
140895fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
140929046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);
1410eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
1411eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
1412eb02082bSJunchao Zhang   ierr = (*sf->ops->ReduceEnd)(sf,unit,leafmtype,leafdata,rootmtype,rootdata,op);CHKERRQ(ierr);
1413563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);
141495fce210SBarry Smith   PetscFunctionReturn(0);
141595fce210SBarry Smith }
141695fce210SBarry Smith 
141795fce210SBarry Smith /*@C
1418a1729e3fSJunchao Zhang    PetscSFFetchAndOpBegin - begin operation that fetches values from root and updates atomically by applying operation using my leaf value, to be completed with PetscSFFetchAndOpEnd()
1419a1729e3fSJunchao Zhang 
1420a1729e3fSJunchao Zhang    Collective
1421a1729e3fSJunchao Zhang 
1422a1729e3fSJunchao Zhang    Input Arguments:
1423a1729e3fSJunchao Zhang +  sf - star forest
1424a1729e3fSJunchao Zhang .  unit - data type
1425a1729e3fSJunchao Zhang .  leafdata - leaf values to use in reduction
1426a1729e3fSJunchao Zhang -  op - operation to use for reduction
1427a1729e3fSJunchao Zhang 
1428a1729e3fSJunchao Zhang    Output Arguments:
1429a1729e3fSJunchao Zhang +  rootdata - root values to be updated, input state is seen by first process to perform an update
1430a1729e3fSJunchao Zhang -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
1431a1729e3fSJunchao Zhang 
1432a1729e3fSJunchao Zhang    Level: advanced
1433a1729e3fSJunchao Zhang 
1434a1729e3fSJunchao Zhang    Note:
1435a1729e3fSJunchao Zhang    The update is only atomic at the granularity provided by the hardware. Different roots referenced by the same process
1436a1729e3fSJunchao Zhang    might be updated in a different order. Furthermore, if a composite type is used for the unit datatype, atomicity is
1437a1729e3fSJunchao Zhang    not guaranteed across the whole vertex. Therefore, this function is mostly only used with primitive types such as
1438a1729e3fSJunchao Zhang    integers.
1439a1729e3fSJunchao Zhang 
1440a1729e3fSJunchao Zhang .seealso: PetscSFComputeDegreeBegin(), PetscSFReduceBegin(), PetscSFSetGraph()
1441a1729e3fSJunchao Zhang @*/
1442a1729e3fSJunchao Zhang PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1443a1729e3fSJunchao Zhang {
1444a1729e3fSJunchao Zhang   PetscErrorCode ierr;
1445eb02082bSJunchao Zhang   PetscMemType   rootmtype,leafmtype,leafupdatemtype;
1446a1729e3fSJunchao Zhang 
1447a1729e3fSJunchao Zhang   PetscFunctionBegin;
1448a1729e3fSJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1449a1729e3fSJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1450a1729e3fSJunchao Zhang   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
1451eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
1452eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
1453eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafupdate,&leafupdatemtype);CHKERRQ(ierr);
1454eb02082bSJunchao Zhang   if (leafmtype != leafupdatemtype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for leafdata and leafupdate in different memory types");
1455eb02082bSJunchao Zhang #if defined(PETSC_HAVE_CUDA)
1456eb02082bSJunchao Zhang   if (rootmtype == PETSC_MEMTYPE_DEVICE || leafmtype == PETSC_MEMTYPE_DEVICE) {cudaError_t err = cudaDeviceSynchronize();CHKERRCUDA(err);}
1457eb02082bSJunchao Zhang #endif
1458eb02082bSJunchao Zhang   ierr = (*sf->ops->FetchAndOpBegin)(sf,unit,rootmtype,rootdata,leafmtype,leafdata,leafupdate,op);CHKERRQ(ierr);
1459a1729e3fSJunchao Zhang   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
1460a1729e3fSJunchao Zhang   PetscFunctionReturn(0);
1461a1729e3fSJunchao Zhang }
1462a1729e3fSJunchao Zhang 
1463a1729e3fSJunchao Zhang /*@C
1464a1729e3fSJunchao 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
1465a1729e3fSJunchao Zhang 
1466a1729e3fSJunchao Zhang    Collective
1467a1729e3fSJunchao Zhang 
1468a1729e3fSJunchao Zhang    Input Arguments:
1469a1729e3fSJunchao Zhang +  sf - star forest
1470a1729e3fSJunchao Zhang .  unit - data type
1471a1729e3fSJunchao Zhang .  leafdata - leaf values to use in reduction
1472a1729e3fSJunchao Zhang -  op - operation to use for reduction
1473a1729e3fSJunchao Zhang 
1474a1729e3fSJunchao Zhang    Output Arguments:
1475a1729e3fSJunchao Zhang +  rootdata - root values to be updated, input state is seen by first process to perform an update
1476a1729e3fSJunchao Zhang -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
1477a1729e3fSJunchao Zhang 
1478a1729e3fSJunchao Zhang    Level: advanced
1479a1729e3fSJunchao Zhang 
1480a1729e3fSJunchao Zhang .seealso: PetscSFComputeDegreeEnd(), PetscSFReduceEnd(), PetscSFSetGraph()
1481a1729e3fSJunchao Zhang @*/
1482a1729e3fSJunchao Zhang PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1483a1729e3fSJunchao Zhang {
1484a1729e3fSJunchao Zhang   PetscErrorCode ierr;
1485eb02082bSJunchao Zhang   PetscMemType   rootmtype,leafmtype,leafupdatemtype;
1486a1729e3fSJunchao Zhang 
1487a1729e3fSJunchao Zhang   PetscFunctionBegin;
1488a1729e3fSJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1489a1729e3fSJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1490a1729e3fSJunchao Zhang   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
1491eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
1492eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
1493eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafupdate,&leafupdatemtype);CHKERRQ(ierr);
1494eb02082bSJunchao Zhang   if (leafmtype != leafupdatemtype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for leafdata and leafupdate in different memory types");
1495eb02082bSJunchao Zhang   ierr = (*sf->ops->FetchAndOpEnd)(sf,unit,rootmtype,rootdata,leafmtype,leafdata,leafupdate,op);CHKERRQ(ierr);
1496a1729e3fSJunchao Zhang   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
1497a1729e3fSJunchao Zhang   PetscFunctionReturn(0);
1498a1729e3fSJunchao Zhang }
1499a1729e3fSJunchao Zhang 
1500a1729e3fSJunchao Zhang /*@C
150195fce210SBarry Smith    PetscSFComputeDegreeBegin - begin computation of degree for each root vertex, to be completed with PetscSFComputeDegreeEnd()
150295fce210SBarry Smith 
150395fce210SBarry Smith    Collective
150495fce210SBarry Smith 
150595fce210SBarry Smith    Input Arguments:
150695fce210SBarry Smith .  sf - star forest
150795fce210SBarry Smith 
150895fce210SBarry Smith    Output Arguments:
150995fce210SBarry Smith .  degree - degree of each root vertex
151095fce210SBarry Smith 
151195fce210SBarry Smith    Level: advanced
151295fce210SBarry Smith 
1513ffe67aa5SVáclav Hapla    Notes:
1514ffe67aa5SVáclav Hapla    The returned array is owned by PetscSF and automatically freed by PetscSFDestroy(). Hence no need to call PetscFree() on it.
1515ffe67aa5SVáclav Hapla 
151695fce210SBarry Smith .seealso: PetscSFGatherBegin()
151795fce210SBarry Smith @*/
151895fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf,const PetscInt **degree)
151995fce210SBarry Smith {
152095fce210SBarry Smith   PetscErrorCode ierr;
152195fce210SBarry Smith 
152295fce210SBarry Smith   PetscFunctionBegin;
152395fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
152495fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
152595fce210SBarry Smith   PetscValidPointer(degree,2);
1526803bd9e8SMatthew G. Knepley   if (!sf->degreeknown) {
152729046d53SLisandro Dalcin     PetscInt i, nroots = sf->nroots, maxlocal = sf->maxleaf+1;  /* TODO: We should use PetscSFGetLeafRange() */
1528803bd9e8SMatthew G. Knepley     if (sf->degree) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Calls to PetscSFComputeDegreeBegin() cannot be nested.");
152929046d53SLisandro Dalcin     ierr = PetscMalloc1(nroots,&sf->degree);CHKERRQ(ierr);
153029046d53SLisandro Dalcin     ierr = PetscMalloc1(PetscMax(maxlocal,1),&sf->degreetmp);CHKERRQ(ierr); /* allocate at least one entry, see check in PetscSFComputeDegreeEnd() */
153129046d53SLisandro Dalcin     for (i=0; i<nroots; i++) sf->degree[i] = 0;
15329837ea96SMatthew G. Knepley     for (i=0; i<maxlocal; i++) sf->degreetmp[i] = 1;
1533dbd2ff41SBarry Smith     ierr = PetscSFReduceBegin(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);CHKERRQ(ierr);
153495fce210SBarry Smith   }
153595fce210SBarry Smith   *degree = NULL;
153695fce210SBarry Smith   PetscFunctionReturn(0);
153795fce210SBarry Smith }
153895fce210SBarry Smith 
153995fce210SBarry Smith /*@C
154095fce210SBarry Smith    PetscSFComputeDegreeEnd - complete computation of degree for each root vertex, started with PetscSFComputeDegreeBegin()
154195fce210SBarry Smith 
154295fce210SBarry Smith    Collective
154395fce210SBarry Smith 
154495fce210SBarry Smith    Input Arguments:
154595fce210SBarry Smith .  sf - star forest
154695fce210SBarry Smith 
154795fce210SBarry Smith    Output Arguments:
154895fce210SBarry Smith .  degree - degree of each root vertex
154995fce210SBarry Smith 
155095fce210SBarry Smith    Level: developer
155195fce210SBarry Smith 
1552ffe67aa5SVáclav Hapla    Notes:
1553ffe67aa5SVáclav Hapla    The returned array is owned by PetscSF and automatically freed by PetscSFDestroy(). Hence no need to call PetscFree() on it.
1554ffe67aa5SVáclav Hapla 
155595fce210SBarry Smith .seealso:
155695fce210SBarry Smith @*/
155795fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf,const PetscInt **degree)
155895fce210SBarry Smith {
155995fce210SBarry Smith   PetscErrorCode ierr;
156095fce210SBarry Smith 
156195fce210SBarry Smith   PetscFunctionBegin;
156295fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
156395fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
156429046d53SLisandro Dalcin   PetscValidPointer(degree,2);
156595fce210SBarry Smith   if (!sf->degreeknown) {
156629046d53SLisandro Dalcin     if (!sf->degreetmp) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFComputeDegreeBegin() before PetscSFComputeDegreeEnd()");
1567dbd2ff41SBarry Smith     ierr = PetscSFReduceEnd(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);CHKERRQ(ierr);
156895fce210SBarry Smith     ierr = PetscFree(sf->degreetmp);CHKERRQ(ierr);
156995fce210SBarry Smith     sf->degreeknown = PETSC_TRUE;
157095fce210SBarry Smith   }
157195fce210SBarry Smith   *degree = sf->degree;
157295fce210SBarry Smith   PetscFunctionReturn(0);
157395fce210SBarry Smith }
157495fce210SBarry Smith 
1575673100f5SVaclav Hapla 
1576673100f5SVaclav Hapla /*@C
157766dfcd1aSVaclav Hapla    PetscSFComputeMultiRootOriginalNumbering - Returns original numbering of multi-roots (roots of multi-SF returned by PetscSFGetMultiSF()).
157866dfcd1aSVaclav Hapla    Each multi-root is assigned index of the corresponding original root.
1579673100f5SVaclav Hapla 
1580673100f5SVaclav Hapla    Collective
1581673100f5SVaclav Hapla 
1582673100f5SVaclav Hapla    Input Arguments:
1583673100f5SVaclav Hapla +  sf - star forest
1584673100f5SVaclav Hapla -  degree - degree of each root vertex, computed with PetscSFComputeDegreeBegin()/PetscSFComputeDegreeEnd()
1585673100f5SVaclav Hapla 
1586673100f5SVaclav Hapla    Output Arguments:
158766dfcd1aSVaclav Hapla +  nMultiRoots - (optional) number of multi-roots (roots of multi-SF)
158866dfcd1aSVaclav Hapla -  multiRootsOrigNumbering - original indices of multi-roots; length of this array is nMultiRoots
1589673100f5SVaclav Hapla 
1590673100f5SVaclav Hapla    Level: developer
1591673100f5SVaclav Hapla 
1592ffe67aa5SVáclav Hapla    Notes:
1593ffe67aa5SVáclav Hapla    The returned array multiRootsOrigNumbering is newly allocated and should be destroyed with PetscFree() when no longer needed.
1594ffe67aa5SVáclav Hapla 
1595673100f5SVaclav Hapla .seealso: PetscSFComputeDegreeBegin(), PetscSFComputeDegreeEnd(), PetscSFGetMultiSF()
1596673100f5SVaclav Hapla @*/
159766dfcd1aSVaclav Hapla PetscErrorCode PetscSFComputeMultiRootOriginalNumbering(PetscSF sf, const PetscInt degree[], PetscInt *nMultiRoots, PetscInt *multiRootsOrigNumbering[])
1598673100f5SVaclav Hapla {
1599673100f5SVaclav Hapla   PetscSF             msf;
1600673100f5SVaclav Hapla   PetscInt            i, j, k, nroots, nmroots;
1601673100f5SVaclav Hapla   PetscErrorCode      ierr;
1602673100f5SVaclav Hapla 
1603673100f5SVaclav Hapla   PetscFunctionBegin;
1604673100f5SVaclav Hapla   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1605673100f5SVaclav Hapla   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
160629328920SVaclav Hapla   if (nroots) PetscValidIntPointer(degree,2);
160766dfcd1aSVaclav Hapla   if (nMultiRoots) PetscValidIntPointer(nMultiRoots,3);
160866dfcd1aSVaclav Hapla   PetscValidPointer(multiRootsOrigNumbering,4);
1609673100f5SVaclav Hapla   ierr = PetscSFGetMultiSF(sf,&msf);CHKERRQ(ierr);
1610673100f5SVaclav Hapla   ierr = PetscSFGetGraph(msf, &nmroots, NULL, NULL, NULL);CHKERRQ(ierr);
161166dfcd1aSVaclav Hapla   ierr = PetscMalloc1(nmroots, multiRootsOrigNumbering);CHKERRQ(ierr);
1612673100f5SVaclav Hapla   for (i=0,j=0,k=0; i<nroots; i++) {
1613673100f5SVaclav Hapla     if (!degree[i]) continue;
1614673100f5SVaclav Hapla     for (j=0; j<degree[i]; j++,k++) {
161566dfcd1aSVaclav Hapla       (*multiRootsOrigNumbering)[k] = i;
1616673100f5SVaclav Hapla     }
1617673100f5SVaclav Hapla   }
1618673100f5SVaclav Hapla #if defined(PETSC_USE_DEBUG)
1619673100f5SVaclav Hapla   if (PetscUnlikely(k != nmroots)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sanity check fail");
1620673100f5SVaclav Hapla #endif
162166dfcd1aSVaclav Hapla   if (nMultiRoots) *nMultiRoots = nmroots;
1622673100f5SVaclav Hapla   PetscFunctionReturn(0);
1623673100f5SVaclav Hapla }
1624673100f5SVaclav Hapla 
162595fce210SBarry Smith /*@C
162695fce210SBarry Smith    PetscSFGatherBegin - begin pointwise gather of all leaves into multi-roots, to be completed with PetscSFGatherEnd()
162795fce210SBarry Smith 
162895fce210SBarry Smith    Collective
162995fce210SBarry Smith 
163095fce210SBarry Smith    Input Arguments:
163195fce210SBarry Smith +  sf - star forest
163295fce210SBarry Smith .  unit - data type
163395fce210SBarry Smith -  leafdata - leaf data to gather to roots
163495fce210SBarry Smith 
163595fce210SBarry Smith    Output Argument:
163695fce210SBarry Smith .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
163795fce210SBarry Smith 
163895fce210SBarry Smith    Level: intermediate
163995fce210SBarry Smith 
164095fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
164195fce210SBarry Smith @*/
164295fce210SBarry Smith PetscErrorCode PetscSFGatherBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
164395fce210SBarry Smith {
164495fce210SBarry Smith   PetscErrorCode ierr;
164595fce210SBarry Smith   PetscSF        multi;
164695fce210SBarry Smith 
164795fce210SBarry Smith   PetscFunctionBegin;
164895fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
164929046d53SLisandro Dalcin   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
165095fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
16518bfbc91cSJed Brown   ierr = PetscSFReduceBegin(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr);
165295fce210SBarry Smith   PetscFunctionReturn(0);
165395fce210SBarry Smith }
165495fce210SBarry Smith 
165595fce210SBarry Smith /*@C
165695fce210SBarry Smith    PetscSFGatherEnd - ends pointwise gather operation that was started with PetscSFGatherBegin()
165795fce210SBarry Smith 
165895fce210SBarry Smith    Collective
165995fce210SBarry Smith 
166095fce210SBarry Smith    Input Arguments:
166195fce210SBarry Smith +  sf - star forest
166295fce210SBarry Smith .  unit - data type
166395fce210SBarry Smith -  leafdata - leaf data to gather to roots
166495fce210SBarry Smith 
166595fce210SBarry Smith    Output Argument:
166695fce210SBarry Smith .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
166795fce210SBarry Smith 
166895fce210SBarry Smith    Level: intermediate
166995fce210SBarry Smith 
167095fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
167195fce210SBarry Smith @*/
167295fce210SBarry Smith PetscErrorCode PetscSFGatherEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
167395fce210SBarry Smith {
167495fce210SBarry Smith   PetscErrorCode ierr;
167595fce210SBarry Smith   PetscSF        multi;
167695fce210SBarry Smith 
167795fce210SBarry Smith   PetscFunctionBegin;
167895fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
167995fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
168095fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
16818bfbc91cSJed Brown   ierr = PetscSFReduceEnd(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr);
168295fce210SBarry Smith   PetscFunctionReturn(0);
168395fce210SBarry Smith }
168495fce210SBarry Smith 
168595fce210SBarry Smith /*@C
168695fce210SBarry Smith    PetscSFScatterBegin - begin pointwise scatter operation from multi-roots to leaves, to be completed with PetscSFScatterEnd()
168795fce210SBarry Smith 
168895fce210SBarry Smith    Collective
168995fce210SBarry Smith 
169095fce210SBarry Smith    Input Arguments:
169195fce210SBarry Smith +  sf - star forest
169295fce210SBarry Smith .  unit - data type
169395fce210SBarry Smith -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
169495fce210SBarry Smith 
169595fce210SBarry Smith    Output Argument:
169695fce210SBarry Smith .  leafdata - leaf data to be update with personal data from each respective root
169795fce210SBarry Smith 
169895fce210SBarry Smith    Level: intermediate
169995fce210SBarry Smith 
170095fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
170195fce210SBarry Smith @*/
170295fce210SBarry Smith PetscErrorCode PetscSFScatterBegin(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
170395fce210SBarry Smith {
170495fce210SBarry Smith   PetscErrorCode ierr;
170595fce210SBarry Smith   PetscSF        multi;
170695fce210SBarry Smith 
170795fce210SBarry Smith   PetscFunctionBegin;
170895fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
170995fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
171095fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
171195fce210SBarry Smith   ierr = PetscSFBcastBegin(multi,unit,multirootdata,leafdata);CHKERRQ(ierr);
171295fce210SBarry Smith   PetscFunctionReturn(0);
171395fce210SBarry Smith }
171495fce210SBarry Smith 
171595fce210SBarry Smith /*@C
171695fce210SBarry Smith    PetscSFScatterEnd - ends pointwise scatter operation that was started with PetscSFScatterBegin()
171795fce210SBarry Smith 
171895fce210SBarry Smith    Collective
171995fce210SBarry Smith 
172095fce210SBarry Smith    Input Arguments:
172195fce210SBarry Smith +  sf - star forest
172295fce210SBarry Smith .  unit - data type
172395fce210SBarry Smith -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
172495fce210SBarry Smith 
172595fce210SBarry Smith    Output Argument:
172695fce210SBarry Smith .  leafdata - leaf data to be update with personal data from each respective root
172795fce210SBarry Smith 
172895fce210SBarry Smith    Level: intermediate
172995fce210SBarry Smith 
173095fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
173195fce210SBarry Smith @*/
173295fce210SBarry Smith PetscErrorCode PetscSFScatterEnd(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
173395fce210SBarry Smith {
173495fce210SBarry Smith   PetscErrorCode ierr;
173595fce210SBarry Smith   PetscSF        multi;
173695fce210SBarry Smith 
173795fce210SBarry Smith   PetscFunctionBegin;
173895fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
173995fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
174095fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
174195fce210SBarry Smith   ierr = PetscSFBcastEnd(multi,unit,multirootdata,leafdata);CHKERRQ(ierr);
174295fce210SBarry Smith   PetscFunctionReturn(0);
174395fce210SBarry Smith }
1744a7b3aa13SAta Mesgarnejad 
1745a072220fSLawrence Mitchell static PetscErrorCode PetscSFCheckLeavesUnique_Private(PetscSF sf)
1746a072220fSLawrence Mitchell {
1747a072220fSLawrence Mitchell #if defined(PETSC_USE_DEBUG)
1748a072220fSLawrence Mitchell   PetscInt        i, n, nleaves;
1749a072220fSLawrence Mitchell   const PetscInt *ilocal = NULL;
1750a072220fSLawrence Mitchell   PetscHSetI      seen;
1751a072220fSLawrence Mitchell   PetscErrorCode  ierr;
1752a072220fSLawrence Mitchell 
1753a072220fSLawrence Mitchell   PetscFunctionBegin;
1754a072220fSLawrence Mitchell   ierr = PetscSFGetGraph(sf,NULL,&nleaves,&ilocal,NULL);CHKERRQ(ierr);
1755a072220fSLawrence Mitchell   ierr = PetscHSetICreate(&seen);CHKERRQ(ierr);
1756a072220fSLawrence Mitchell   for (i = 0; i < nleaves; i++) {
1757a072220fSLawrence Mitchell     const PetscInt leaf = ilocal ? ilocal[i] : i;
1758a072220fSLawrence Mitchell     ierr = PetscHSetIAdd(seen,leaf);CHKERRQ(ierr);
1759a072220fSLawrence Mitchell   }
1760a072220fSLawrence Mitchell   ierr = PetscHSetIGetSize(seen,&n);CHKERRQ(ierr);
1761a072220fSLawrence Mitchell   if (n != nleaves) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Provided leaves have repeated values: all leaves must be unique");
1762a072220fSLawrence Mitchell   ierr = PetscHSetIDestroy(&seen);CHKERRQ(ierr);
1763a072220fSLawrence Mitchell   PetscFunctionReturn(0);
1764a072220fSLawrence Mitchell #else
1765a072220fSLawrence Mitchell   PetscFunctionBegin;
1766a072220fSLawrence Mitchell   PetscFunctionReturn(0);
1767a072220fSLawrence Mitchell #endif
1768a072220fSLawrence Mitchell }
1769a7b3aa13SAta Mesgarnejad /*@
177004c0ada0SJunchao Zhang   PetscSFCompose - Compose a new PetscSF by putting the second SF under the first one in a top (roots) down (leaves) view
1771a7b3aa13SAta Mesgarnejad 
1772a7b3aa13SAta Mesgarnejad   Input Parameters:
177304c0ada0SJunchao Zhang + sfA - The first PetscSF, whose local space may be a permutation, but can not be sparse.
177404c0ada0SJunchao Zhang - sfB - The second PetscSF, whose number of roots must be equal to number of leaves of sfA on each processor
1775a7b3aa13SAta Mesgarnejad 
1776a7b3aa13SAta Mesgarnejad   Output Parameters:
177704c0ada0SJunchao Zhang . sfBA - The composite SF. Doing a Bcast on the new SF is equvalent to doing Bcast on sfA, then Bcast on sfB
1778a7b3aa13SAta Mesgarnejad 
1779a7b3aa13SAta Mesgarnejad   Level: developer
1780a7b3aa13SAta Mesgarnejad 
1781a072220fSLawrence Mitchell   Notes:
1782a072220fSLawrence Mitchell   For the resulting composed SF to be valid, the input SFs must be true star forests: the leaves must be unique. This is
1783a072220fSLawrence Mitchell   checked in debug mode.
1784a072220fSLawrence Mitchell 
178504c0ada0SJunchao Zhang .seealso: PetscSF, PetscSFComposeInverse(), PetscSFGetGraph(), PetscSFSetGraph()
1786a7b3aa13SAta Mesgarnejad @*/
1787a7b3aa13SAta Mesgarnejad PetscErrorCode PetscSFCompose(PetscSF sfA,PetscSF sfB,PetscSF *sfBA)
1788a7b3aa13SAta Mesgarnejad {
178904c0ada0SJunchao Zhang   PetscErrorCode    ierr;
1790a7b3aa13SAta Mesgarnejad   MPI_Comm          comm;
1791a7b3aa13SAta Mesgarnejad   const PetscSFNode *remotePointsA,*remotePointsB;
1792d41018fbSJunchao Zhang   PetscSFNode       *remotePointsBA=NULL,*reorderedRemotePointsA = NULL,*leafdataB;
1793d41018fbSJunchao Zhang   const PetscInt    *localPointsA,*localPointsB,*localPointsBA;
1794d41018fbSJunchao Zhang   PetscInt          i,numRootsA,numLeavesA,numRootsB,numLeavesB,minleaf,maxleaf;
1795a7b3aa13SAta Mesgarnejad 
1796a7b3aa13SAta Mesgarnejad   PetscFunctionBegin;
1797a7b3aa13SAta Mesgarnejad   PetscValidHeaderSpecific(sfA,PETSCSF_CLASSID,1);
179829046d53SLisandro Dalcin   PetscSFCheckGraphSet(sfA,1);
179929046d53SLisandro Dalcin   PetscValidHeaderSpecific(sfB,PETSCSF_CLASSID,2);
180029046d53SLisandro Dalcin   PetscSFCheckGraphSet(sfB,2);
180129046d53SLisandro Dalcin   PetscValidPointer(sfBA,3);
1802a7b3aa13SAta Mesgarnejad   ierr = PetscObjectGetComm((PetscObject)sfA,&comm);CHKERRQ(ierr);
1803a7b3aa13SAta Mesgarnejad   ierr = PetscSFGetGraph(sfA,&numRootsA,&numLeavesA,&localPointsA,&remotePointsA);CHKERRQ(ierr);
1804a7b3aa13SAta Mesgarnejad   ierr = PetscSFGetGraph(sfB,&numRootsB,&numLeavesB,&localPointsB,&remotePointsB);CHKERRQ(ierr);
180504c0ada0SJunchao Zhang   ierr = PetscSFGetLeafRange(sfA,&minleaf,&maxleaf);CHKERRQ(ierr);
180604c0ada0SJunchao Zhang   if (numRootsB != numLeavesA) SETERRQ(comm,PETSC_ERR_ARG_INCOMP,"The second SF's number of roots must be equal to the first SF's number of leaves");
1807a072220fSLawrence Mitchell   if (maxleaf+1 != numLeavesA || minleaf) SETERRQ(comm,PETSC_ERR_ARG_INCOMP,"The first SF can not have sparse local space");
1808a072220fSLawrence Mitchell   /* The above check is fast, but not sufficient, since we cannot guarantee that the SF has unique leaves. So in debug
1809a072220fSLawrence Mitchell    mode, check properly. */
1810a072220fSLawrence Mitchell   ierr = PetscSFCheckLeavesUnique_Private(sfA);CHKERRQ(ierr);
1811d41018fbSJunchao Zhang   if (localPointsA) {
1812d41018fbSJunchao Zhang     /* Local space is dense permutation of identity. Need to rewire order of the remote points */
1813d41018fbSJunchao Zhang     ierr = PetscMalloc1(numLeavesA,&reorderedRemotePointsA);CHKERRQ(ierr);
1814d41018fbSJunchao Zhang     for (i=0; i<numLeavesA; i++) reorderedRemotePointsA[localPointsA[i]-minleaf] = remotePointsA[i];
1815d41018fbSJunchao Zhang     remotePointsA = reorderedRemotePointsA;
1816d41018fbSJunchao Zhang   }
1817d41018fbSJunchao Zhang   ierr = PetscSFGetLeafRange(sfB,&minleaf,&maxleaf);CHKERRQ(ierr);
1818d41018fbSJunchao Zhang   ierr = PetscMalloc1(maxleaf-minleaf+1,&leafdataB);CHKERRQ(ierr);
1819d41018fbSJunchao Zhang   ierr = PetscSFBcastBegin(sfB,MPIU_2INT,remotePointsA,leafdataB-minleaf);CHKERRQ(ierr);
1820d41018fbSJunchao Zhang   ierr = PetscSFBcastEnd(sfB,MPIU_2INT,remotePointsA,leafdataB-minleaf);CHKERRQ(ierr);
1821d41018fbSJunchao Zhang   ierr = PetscFree(reorderedRemotePointsA);CHKERRQ(ierr);
1822d41018fbSJunchao Zhang 
1823a072220fSLawrence Mitchell   /* sfB's leaves must be unique, otherwise BcastAndOp(B o A) != BcastAndOp(B) o BcastAndOp(A) */
1824a072220fSLawrence Mitchell   ierr = PetscSFCheckLeavesUnique_Private(sfB);CHKERRQ(ierr);
1825d41018fbSJunchao Zhang   if (minleaf == 0 && maxleaf + 1 == numLeavesB) { /* Local space of sfB is an identity or permutation */
1826d41018fbSJunchao Zhang     localPointsBA  = NULL;
1827d41018fbSJunchao Zhang     remotePointsBA = leafdataB;
1828d41018fbSJunchao Zhang   } else {
1829d41018fbSJunchao Zhang     localPointsBA  = localPointsB;
1830a7b3aa13SAta Mesgarnejad     ierr = PetscMalloc1(numLeavesB,&remotePointsBA);CHKERRQ(ierr);
1831d41018fbSJunchao Zhang     for (i=0; i<numLeavesB; i++) remotePointsBA[i] = leafdataB[localPointsB[i]-minleaf];
1832d41018fbSJunchao Zhang     ierr = PetscFree(leafdataB);CHKERRQ(ierr);
1833d41018fbSJunchao Zhang   }
1834a7b3aa13SAta Mesgarnejad   ierr = PetscSFCreate(comm, sfBA);CHKERRQ(ierr);
1835d41018fbSJunchao Zhang   ierr = PetscSFSetGraph(*sfBA,numRootsA,numLeavesB,localPointsBA,PETSC_COPY_VALUES,remotePointsBA,PETSC_OWN_POINTER);CHKERRQ(ierr);
1836a7b3aa13SAta Mesgarnejad   PetscFunctionReturn(0);
1837a7b3aa13SAta Mesgarnejad }
18381c6ba672SJunchao Zhang 
183904c0ada0SJunchao Zhang /*@
184004c0ada0SJunchao Zhang   PetscSFComposeInverse - Compose a new PetscSF by putting inverse of the second SF under the first one
184104c0ada0SJunchao Zhang 
184204c0ada0SJunchao Zhang   Input Parameters:
184304c0ada0SJunchao Zhang + sfA - The first PetscSF, whose local space may be a permutation, but can not be sparse.
184404c0ada0SJunchao Zhang - sfB - The second PetscSF, whose number of leaves must be equal to number of leaves of sfA on each processor. All roots must have degree 1.
184504c0ada0SJunchao Zhang 
184604c0ada0SJunchao Zhang   Output Parameters:
184704c0ada0SJunchao Zhang . sfBA - The composite SF. Doing a Bcast on the new SF is equvalent to doing Bcast on sfA, then Bcast on inverse of sfB
184804c0ada0SJunchao Zhang 
184904c0ada0SJunchao Zhang   Level: developer
185004c0ada0SJunchao Zhang 
185104c0ada0SJunchao Zhang .seealso: PetscSF, PetscSFCompose(), PetscSFGetGraph(), PetscSFSetGraph()
185204c0ada0SJunchao Zhang @*/
185304c0ada0SJunchao Zhang PetscErrorCode PetscSFComposeInverse(PetscSF sfA,PetscSF sfB,PetscSF *sfBA)
185404c0ada0SJunchao Zhang {
185504c0ada0SJunchao Zhang   PetscErrorCode    ierr;
185604c0ada0SJunchao Zhang   MPI_Comm          comm;
185704c0ada0SJunchao Zhang   const PetscSFNode *remotePointsA,*remotePointsB;
185804c0ada0SJunchao Zhang   PetscSFNode       *remotePointsBA;
185904c0ada0SJunchao Zhang   const PetscInt    *localPointsA,*localPointsB;
186004c0ada0SJunchao Zhang   PetscInt          numRootsA,numLeavesA,numRootsB,numLeavesB,minleaf,maxleaf;
186104c0ada0SJunchao Zhang 
186204c0ada0SJunchao Zhang   PetscFunctionBegin;
186304c0ada0SJunchao Zhang   PetscValidHeaderSpecific(sfA,PETSCSF_CLASSID,1);
186404c0ada0SJunchao Zhang   PetscSFCheckGraphSet(sfA,1);
186504c0ada0SJunchao Zhang   PetscValidHeaderSpecific(sfB,PETSCSF_CLASSID,2);
186604c0ada0SJunchao Zhang   PetscSFCheckGraphSet(sfB,2);
186704c0ada0SJunchao Zhang   PetscValidPointer(sfBA,3);
186804c0ada0SJunchao Zhang   ierr = PetscObjectGetComm((PetscObject) sfA, &comm);CHKERRQ(ierr);
186904c0ada0SJunchao Zhang   ierr = PetscSFGetGraph(sfA, &numRootsA, &numLeavesA, &localPointsA, &remotePointsA);CHKERRQ(ierr);
187004c0ada0SJunchao Zhang   ierr = PetscSFGetGraph(sfB, &numRootsB, &numLeavesB, &localPointsB, &remotePointsB);CHKERRQ(ierr);
187104c0ada0SJunchao Zhang   ierr = PetscSFGetLeafRange(sfA,&minleaf,&maxleaf);CHKERRQ(ierr);
187204c0ada0SJunchao Zhang   if (maxleaf+1-minleaf != numLeavesA) SETERRQ(comm,PETSC_ERR_ARG_INCOMP,"The first SF can not have sparse local space");
187304c0ada0SJunchao Zhang   if (numLeavesA != numLeavesB) SETERRQ(comm,PETSC_ERR_ARG_INCOMP,"The second SF's number of leaves must be equal to the first SF's number of leaves");
187404c0ada0SJunchao Zhang   ierr = PetscMalloc1(numRootsB,&remotePointsBA);CHKERRQ(ierr);
187504c0ada0SJunchao Zhang   ierr = PetscSFReduceBegin(sfB,MPIU_2INT,remotePointsA,remotePointsBA,MPIU_REPLACE);CHKERRQ(ierr);
187604c0ada0SJunchao Zhang   ierr = PetscSFReduceEnd(sfB,MPIU_2INT,remotePointsA,remotePointsBA,MPIU_REPLACE);CHKERRQ(ierr);
187704c0ada0SJunchao Zhang   ierr = PetscSFCreate(comm,sfBA);CHKERRQ(ierr);
187804c0ada0SJunchao Zhang   ierr = PetscSFSetGraph(*sfBA,numRootsA,numRootsB,NULL,PETSC_COPY_VALUES,remotePointsBA,PETSC_OWN_POINTER);CHKERRQ(ierr);
187904c0ada0SJunchao Zhang   PetscFunctionReturn(0);
188004c0ada0SJunchao Zhang }
188104c0ada0SJunchao Zhang 
18821c6ba672SJunchao Zhang /*
18831c6ba672SJunchao Zhang   PetscSFCreateLocalSF_Private - Creates a local PetscSF that only has intra-process edges of the global PetscSF
18841c6ba672SJunchao Zhang 
18851c6ba672SJunchao Zhang   Input Parameters:
18861c6ba672SJunchao Zhang . sf - The global PetscSF
18871c6ba672SJunchao Zhang 
18881c6ba672SJunchao Zhang   Output Parameters:
18891c6ba672SJunchao Zhang . out - The local PetscSF
18901c6ba672SJunchao Zhang  */
18911c6ba672SJunchao Zhang PetscErrorCode PetscSFCreateLocalSF_Private(PetscSF sf,PetscSF *out)
18921c6ba672SJunchao Zhang {
18931c6ba672SJunchao Zhang   MPI_Comm           comm;
18941c6ba672SJunchao Zhang   PetscMPIInt        myrank;
18951c6ba672SJunchao Zhang   const PetscInt     *ilocal;
18961c6ba672SJunchao Zhang   const PetscSFNode  *iremote;
18971c6ba672SJunchao Zhang   PetscInt           i,j,nroots,nleaves,lnleaves,*lilocal;
18981c6ba672SJunchao Zhang   PetscSFNode        *liremote;
18991c6ba672SJunchao Zhang   PetscSF            lsf;
19001c6ba672SJunchao Zhang   PetscErrorCode     ierr;
19011c6ba672SJunchao Zhang 
19021c6ba672SJunchao Zhang   PetscFunctionBegin;
19031c6ba672SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
19041c6ba672SJunchao Zhang   if (sf->ops->CreateLocalSF) {
19051c6ba672SJunchao Zhang     ierr = (*sf->ops->CreateLocalSF)(sf,out);CHKERRQ(ierr);
19061c6ba672SJunchao Zhang   } else {
19071c6ba672SJunchao Zhang     /* Could use PetscSFCreateEmbeddedLeafSF, but since we know the comm is PETSC_COMM_SELF, we can make it fast */
19081c6ba672SJunchao Zhang     ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
19091c6ba672SJunchao Zhang     ierr = MPI_Comm_rank(comm,&myrank);CHKERRQ(ierr);
19101c6ba672SJunchao Zhang 
19111c6ba672SJunchao Zhang     /* Find out local edges and build a local SF */
19121c6ba672SJunchao Zhang     ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
19131c6ba672SJunchao Zhang     for (i=lnleaves=0; i<nleaves; i++) {if (iremote[i].rank == (PetscInt)myrank) lnleaves++;}
19141c6ba672SJunchao Zhang     ierr = PetscMalloc1(lnleaves,&lilocal);CHKERRQ(ierr);
19151c6ba672SJunchao Zhang     ierr = PetscMalloc1(lnleaves,&liremote);CHKERRQ(ierr);
19161c6ba672SJunchao Zhang 
19171c6ba672SJunchao Zhang     for (i=j=0; i<nleaves; i++) {
19181c6ba672SJunchao Zhang       if (iremote[i].rank == (PetscInt)myrank) {
19191c6ba672SJunchao Zhang         lilocal[j]        = ilocal? ilocal[i] : i; /* ilocal=NULL for contiguous storage */
19201c6ba672SJunchao Zhang         liremote[j].rank  = 0; /* rank in PETSC_COMM_SELF */
19211c6ba672SJunchao Zhang         liremote[j].index = iremote[i].index;
19221c6ba672SJunchao Zhang         j++;
19231c6ba672SJunchao Zhang       }
19241c6ba672SJunchao Zhang     }
19251c6ba672SJunchao Zhang     ierr = PetscSFCreate(PETSC_COMM_SELF,&lsf);CHKERRQ(ierr);
19261c6ba672SJunchao Zhang     ierr = PetscSFSetGraph(lsf,nroots,lnleaves,lilocal,PETSC_OWN_POINTER,liremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
19271c6ba672SJunchao Zhang     ierr = PetscSFSetUp(lsf);CHKERRQ(ierr);
19281c6ba672SJunchao Zhang     *out = lsf;
19291c6ba672SJunchao Zhang   }
19301c6ba672SJunchao Zhang   PetscFunctionReturn(0);
19311c6ba672SJunchao Zhang }
1932dd5b3ca6SJunchao Zhang 
1933dd5b3ca6SJunchao Zhang /* Similar to PetscSFBcast, but only Bcast to leaves on rank 0 */
1934dd5b3ca6SJunchao Zhang PetscErrorCode PetscSFBcastToZero_Private(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
1935dd5b3ca6SJunchao Zhang {
1936dd5b3ca6SJunchao Zhang   PetscErrorCode     ierr;
1937eb02082bSJunchao Zhang   PetscMemType       rootmtype,leafmtype;
1938dd5b3ca6SJunchao Zhang 
1939dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
1940dd5b3ca6SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1941dd5b3ca6SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1942dd5b3ca6SJunchao Zhang   ierr = PetscLogEventBegin(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
1943eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
1944eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
1945dd5b3ca6SJunchao Zhang   if (sf->ops->BcastToZero) {
1946eb02082bSJunchao Zhang     ierr = (*sf->ops->BcastToZero)(sf,unit,rootmtype,rootdata,leafmtype,leafdata);CHKERRQ(ierr);
1947dd5b3ca6SJunchao Zhang   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSFBcastToZero_Private is not supported on this SF type");
1948dd5b3ca6SJunchao Zhang   ierr = PetscLogEventEnd(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
1949dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
1950dd5b3ca6SJunchao Zhang }
1951dd5b3ca6SJunchao Zhang 
1952