xref: /petsc/src/vec/is/sf/interface/sf.c (revision b85e67b730cce8d297c053b3c61d3977f008dd9a)
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)
97cd620004SJunchao Zhang   {
98cd620004SJunchao Zhang   PetscInt  i;
99cd620004SJunchao Zhang   for (i=0; i<2; i++) {if (sf->rmine_d[i]) {cudaError_t err = cudaFree(sf->rmine_d[i]);CHKERRCUDA(err);sf->rmine_d[i]=NULL;}}
100cd620004SJunchao Zhang   }
101eb02082bSJunchao Zhang #endif
10229046d53SLisandro Dalcin   sf->degreeknown = PETSC_FALSE;
10395fce210SBarry Smith   ierr = PetscFree(sf->degree);CHKERRQ(ierr);
10495fce210SBarry Smith   if (sf->ingroup  != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->ingroup);CHKERRQ(ierr);}
10595fce210SBarry Smith   if (sf->outgroup != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->outgroup);CHKERRQ(ierr);}
10695fce210SBarry Smith   ierr = PetscSFDestroy(&sf->multi);CHKERRQ(ierr);
107dd5b3ca6SJunchao Zhang   ierr = PetscLayoutDestroy(&sf->map);CHKERRQ(ierr);
10895fce210SBarry Smith   sf->setupcalled = PETSC_FALSE;
10995fce210SBarry Smith   PetscFunctionReturn(0);
11095fce210SBarry Smith }
11195fce210SBarry Smith 
11295fce210SBarry Smith /*@C
11329046d53SLisandro Dalcin    PetscSFSetType - Set the PetscSF communication implementation
11495fce210SBarry Smith 
11595fce210SBarry Smith    Collective on PetscSF
11695fce210SBarry Smith 
11795fce210SBarry Smith    Input Parameters:
11895fce210SBarry Smith +  sf - the PetscSF context
11995fce210SBarry Smith -  type - a known method
12095fce210SBarry Smith 
12195fce210SBarry Smith    Options Database Key:
12295fce210SBarry Smith .  -sf_type <type> - Sets the method; use -help for a list
12370616304SStefano Zampini    of available methods (for instance, window, basic, neighbor)
12495fce210SBarry Smith 
12595fce210SBarry Smith    Notes:
12695fce210SBarry Smith    See "include/petscsf.h" for available methods (for instance)
12795fce210SBarry Smith +    PETSCSFWINDOW - MPI-2/3 one-sided
12895fce210SBarry Smith -    PETSCSFBASIC - basic implementation using MPI-1 two-sided
12995fce210SBarry Smith 
13095fce210SBarry Smith   Level: intermediate
13195fce210SBarry Smith 
13295fce210SBarry Smith .seealso: PetscSFType, PetscSFCreate()
13395fce210SBarry Smith @*/
13495fce210SBarry Smith PetscErrorCode PetscSFSetType(PetscSF sf,PetscSFType type)
13595fce210SBarry Smith {
13695fce210SBarry Smith   PetscErrorCode ierr,(*r)(PetscSF);
13795fce210SBarry Smith   PetscBool      match;
13895fce210SBarry Smith 
13995fce210SBarry Smith   PetscFunctionBegin;
14095fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
14195fce210SBarry Smith   PetscValidCharPointer(type,2);
14295fce210SBarry Smith 
14395fce210SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)sf,type,&match);CHKERRQ(ierr);
14495fce210SBarry Smith   if (match) PetscFunctionReturn(0);
14595fce210SBarry Smith 
146adc40e5bSBarry Smith   ierr = PetscFunctionListFind(PetscSFList,type,&r);CHKERRQ(ierr);
14795fce210SBarry Smith   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested PetscSF type %s",type);
14829046d53SLisandro Dalcin   /* Destroy the previous PetscSF implementation context */
14929046d53SLisandro Dalcin   if (sf->ops->Destroy) {ierr = (*(sf)->ops->Destroy)(sf);CHKERRQ(ierr);}
15095fce210SBarry Smith   ierr = PetscMemzero(sf->ops,sizeof(*sf->ops));CHKERRQ(ierr);
15195fce210SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)sf,type);CHKERRQ(ierr);
15295fce210SBarry Smith   ierr = (*r)(sf);CHKERRQ(ierr);
15395fce210SBarry Smith   PetscFunctionReturn(0);
15495fce210SBarry Smith }
15595fce210SBarry Smith 
15629046d53SLisandro Dalcin /*@C
15729046d53SLisandro Dalcin   PetscSFGetType - Get the PetscSF communication implementation
15829046d53SLisandro Dalcin 
15929046d53SLisandro Dalcin   Not Collective
16029046d53SLisandro Dalcin 
16129046d53SLisandro Dalcin   Input Parameter:
16229046d53SLisandro Dalcin . sf  - the PetscSF context
16329046d53SLisandro Dalcin 
16429046d53SLisandro Dalcin   Output Parameter:
16529046d53SLisandro Dalcin . type - the PetscSF type name
16629046d53SLisandro Dalcin 
16729046d53SLisandro Dalcin   Level: intermediate
16829046d53SLisandro Dalcin 
16929046d53SLisandro Dalcin .seealso: PetscSFSetType(), PetscSFCreate()
17029046d53SLisandro Dalcin @*/
17129046d53SLisandro Dalcin PetscErrorCode PetscSFGetType(PetscSF sf, PetscSFType *type)
17229046d53SLisandro Dalcin {
17329046d53SLisandro Dalcin   PetscFunctionBegin;
17429046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID,1);
17529046d53SLisandro Dalcin   PetscValidPointer(type,2);
17629046d53SLisandro Dalcin   *type = ((PetscObject)sf)->type_name;
17729046d53SLisandro Dalcin   PetscFunctionReturn(0);
17829046d53SLisandro Dalcin }
17929046d53SLisandro Dalcin 
180d36d33e4SMatthew G. Knepley /*@
18195fce210SBarry Smith    PetscSFDestroy - destroy star forest
18295fce210SBarry Smith 
18395fce210SBarry Smith    Collective
18495fce210SBarry Smith 
18595fce210SBarry Smith    Input Arguments:
18695fce210SBarry Smith .  sf - address of star forest
18795fce210SBarry Smith 
18895fce210SBarry Smith    Level: intermediate
18995fce210SBarry Smith 
19095fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFReset()
19195fce210SBarry Smith @*/
19295fce210SBarry Smith PetscErrorCode PetscSFDestroy(PetscSF *sf)
19395fce210SBarry Smith {
19495fce210SBarry Smith   PetscErrorCode ierr;
19595fce210SBarry Smith 
19695fce210SBarry Smith   PetscFunctionBegin;
19795fce210SBarry Smith   if (!*sf) PetscFunctionReturn(0);
19895fce210SBarry Smith   PetscValidHeaderSpecific((*sf),PETSCSF_CLASSID,1);
19929046d53SLisandro Dalcin   if (--((PetscObject)(*sf))->refct > 0) {*sf = NULL; PetscFunctionReturn(0);}
20095fce210SBarry Smith   ierr = PetscSFReset(*sf);CHKERRQ(ierr);
20195fce210SBarry Smith   if ((*sf)->ops->Destroy) {ierr = (*(*sf)->ops->Destroy)(*sf);CHKERRQ(ierr);}
20295fce210SBarry Smith   ierr = PetscHeaderDestroy(sf);CHKERRQ(ierr);
20395fce210SBarry Smith   PetscFunctionReturn(0);
20495fce210SBarry Smith }
20595fce210SBarry Smith 
206c4e6a40aSLawrence Mitchell static PetscErrorCode PetscSFCheckGraphValid_Private(PetscSF sf)
207c4e6a40aSLawrence Mitchell {
208c4e6a40aSLawrence Mitchell #if defined(PETSC_USE_DEBUG)
209c4e6a40aSLawrence Mitchell   PetscInt           i, nleaves;
210c4e6a40aSLawrence Mitchell   PetscMPIInt        size;
211c4e6a40aSLawrence Mitchell   const PetscInt    *ilocal;
212c4e6a40aSLawrence Mitchell   const PetscSFNode *iremote;
213c4e6a40aSLawrence Mitchell   PetscErrorCode     ierr;
214c4e6a40aSLawrence Mitchell 
215c4e6a40aSLawrence Mitchell   PetscFunctionBegin;
216c4e6a40aSLawrence Mitchell   if (!sf->graphset) PetscFunctionReturn(0);
217c4e6a40aSLawrence Mitchell   ierr = PetscSFGetGraph(sf,NULL,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
218c4e6a40aSLawrence Mitchell   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr);
219c4e6a40aSLawrence Mitchell   for (i = 0; i < nleaves; i++) {
220c4e6a40aSLawrence Mitchell     const PetscInt rank = iremote[i].rank;
221c4e6a40aSLawrence Mitchell     const PetscInt remote = iremote[i].index;
222c4e6a40aSLawrence Mitchell     const PetscInt leaf = ilocal ? ilocal[i] : i;
223c4e6a40aSLawrence 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);
224c4e6a40aSLawrence 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);
225c4e6a40aSLawrence 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);
226c4e6a40aSLawrence Mitchell   }
227c4e6a40aSLawrence Mitchell   PetscFunctionReturn(0);
228c4e6a40aSLawrence Mitchell #else
229c4e6a40aSLawrence Mitchell   PetscFunctionBegin;
230c4e6a40aSLawrence Mitchell   PetscFunctionReturn(0);
231c4e6a40aSLawrence Mitchell #endif
232c4e6a40aSLawrence Mitchell }
233c4e6a40aSLawrence Mitchell 
23495fce210SBarry Smith /*@
23595fce210SBarry Smith    PetscSFSetUp - set up communication structures
23695fce210SBarry Smith 
23795fce210SBarry Smith    Collective
23895fce210SBarry Smith 
23995fce210SBarry Smith    Input Arguments:
24095fce210SBarry Smith .  sf - star forest communication object
24195fce210SBarry Smith 
24295fce210SBarry Smith    Level: beginner
24395fce210SBarry Smith 
24495fce210SBarry Smith .seealso: PetscSFSetFromOptions(), PetscSFSetType()
24595fce210SBarry Smith @*/
24695fce210SBarry Smith PetscErrorCode PetscSFSetUp(PetscSF sf)
24795fce210SBarry Smith {
24895fce210SBarry Smith   PetscErrorCode ierr;
24995fce210SBarry Smith 
25095fce210SBarry Smith   PetscFunctionBegin;
25129046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
25229046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
25395fce210SBarry Smith   if (sf->setupcalled) PetscFunctionReturn(0);
254c4e6a40aSLawrence Mitchell   ierr = PetscSFCheckGraphValid_Private(sf);CHKERRQ(ierr);
25595fce210SBarry Smith   if (!((PetscObject)sf)->type_name) {ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr);}
25629046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_SetUp,sf,0,0,0);CHKERRQ(ierr);
25795fce210SBarry Smith   if (sf->ops->SetUp) {ierr = (*sf->ops->SetUp)(sf);CHKERRQ(ierr);}
25829046d53SLisandro Dalcin   ierr = PetscLogEventEnd(PETSCSF_SetUp,sf,0,0,0);CHKERRQ(ierr);
25995fce210SBarry Smith   sf->setupcalled = PETSC_TRUE;
26095fce210SBarry Smith   PetscFunctionReturn(0);
26195fce210SBarry Smith }
26295fce210SBarry Smith 
2638af6ec1cSBarry Smith /*@
26495fce210SBarry Smith    PetscSFSetFromOptions - set PetscSF options using the options database
26595fce210SBarry Smith 
26695fce210SBarry Smith    Logically Collective
26795fce210SBarry Smith 
26895fce210SBarry Smith    Input Arguments:
26995fce210SBarry Smith .  sf - star forest
27095fce210SBarry Smith 
27195fce210SBarry Smith    Options Database Keys:
27260263706SJed Brown +  -sf_type               - implementation type, see PetscSFSetType()
27351ccb202SJunchao Zhang .  -sf_rank_order         - sort composite points for gathers and scatters in rank order, gathers are non-deterministic otherwise
274*b85e67b7SJunchao Zhang .  -sf_use_default_stream - Assume callers of SF computed the input root/leafdata with the default cuda stream. SF will also
275cd620004SJunchao Zhang                             use the default stream to process data. Therefore, no stream synchronization is needed between SF and its caller.
276*b85e67b7SJunchao Zhang -  -sf_use_stream_aware_mpi  - Assumes the underlying MPI is cuda-stream aware and SF won't sync streams for send/recv buffers passed to MPI.
27795fce210SBarry Smith 
27895fce210SBarry Smith    Level: intermediate
27995fce210SBarry Smith @*/
28095fce210SBarry Smith PetscErrorCode PetscSFSetFromOptions(PetscSF sf)
28195fce210SBarry Smith {
28295fce210SBarry Smith   PetscSFType    deft;
28395fce210SBarry Smith   char           type[256];
28495fce210SBarry Smith   PetscErrorCode ierr;
28595fce210SBarry Smith   PetscBool      flg;
28695fce210SBarry Smith 
28795fce210SBarry Smith   PetscFunctionBegin;
28895fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
28995fce210SBarry Smith   ierr = PetscObjectOptionsBegin((PetscObject)sf);CHKERRQ(ierr);
29095fce210SBarry Smith   deft = ((PetscObject)sf)->type_name ? ((PetscObject)sf)->type_name : PETSCSFBASIC;
29129046d53SLisandro Dalcin   ierr = PetscOptionsFList("-sf_type","PetscSF implementation type","PetscSFSetType",PetscSFList,deft,type,sizeof(type),&flg);CHKERRQ(ierr);
29295fce210SBarry Smith   ierr = PetscSFSetType(sf,flg ? type : deft);CHKERRQ(ierr);
29395fce210SBarry 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);
294cd620004SJunchao Zhang   sf->use_default_stream = PETSC_FALSE;
295cd620004SJunchao Zhang   ierr = PetscOptionsBool("-sf_use_default_stream","Assume SF's input and output root/leafdata is computed on the default stream","PetscSFSetFromOptions",sf->use_default_stream,&sf->use_default_stream,NULL);CHKERRQ(ierr);
296*b85e67b7SJunchao Zhang   sf->use_stream_aware_mpi = PETSC_FALSE;
297*b85e67b7SJunchao Zhang   ierr = PetscOptionsBool("-sf_use_stream_aware_mpi","Assume the underlying MPI is cuda-stream aware","PetscSFSetFromOptions",sf->use_stream_aware_mpi,&sf->use_stream_aware_mpi,NULL);CHKERRQ(ierr);
298e55864a3SBarry Smith   if (sf->ops->SetFromOptions) {ierr = (*sf->ops->SetFromOptions)(PetscOptionsObject,sf);CHKERRQ(ierr);}
29995fce210SBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
30095fce210SBarry Smith   PetscFunctionReturn(0);
30195fce210SBarry Smith }
30295fce210SBarry Smith 
30329046d53SLisandro Dalcin /*@
30495fce210SBarry Smith    PetscSFSetRankOrder - sort multi-points for gathers and scatters by rank order
30595fce210SBarry Smith 
30695fce210SBarry Smith    Logically Collective
30795fce210SBarry Smith 
30895fce210SBarry Smith    Input Arguments:
30995fce210SBarry Smith +  sf - star forest
31095fce210SBarry Smith -  flg - PETSC_TRUE to sort, PETSC_FALSE to skip sorting (lower setup cost, but non-deterministic)
31195fce210SBarry Smith 
31295fce210SBarry Smith    Level: advanced
31395fce210SBarry Smith 
31495fce210SBarry Smith .seealso: PetscSFGatherBegin(), PetscSFScatterBegin()
31595fce210SBarry Smith @*/
31695fce210SBarry Smith PetscErrorCode PetscSFSetRankOrder(PetscSF sf,PetscBool flg)
31795fce210SBarry Smith {
31895fce210SBarry Smith   PetscFunctionBegin;
31995fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
32095fce210SBarry Smith   PetscValidLogicalCollectiveBool(sf,flg,2);
32195fce210SBarry Smith   if (sf->multi) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Rank ordering must be set before first call to PetscSFGatherBegin() or PetscSFScatterBegin()");
32295fce210SBarry Smith   sf->rankorder = flg;
32395fce210SBarry Smith   PetscFunctionReturn(0);
32495fce210SBarry Smith }
32595fce210SBarry Smith 
3268af6ec1cSBarry Smith /*@
32795fce210SBarry Smith    PetscSFSetGraph - Set a parallel star forest
32895fce210SBarry Smith 
32995fce210SBarry Smith    Collective
33095fce210SBarry Smith 
33195fce210SBarry Smith    Input Arguments:
33295fce210SBarry Smith +  sf - star forest
33395fce210SBarry Smith .  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
33495fce210SBarry Smith .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
335c4e6a40aSLawrence Mitchell .  ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage (locations must be >= 0, enforced
336c4e6a40aSLawrence Mitchell during setup in debug mode)
33795fce210SBarry Smith .  localmode - copy mode for ilocal
338c4e6a40aSLawrence Mitchell .  iremote - remote locations of root vertices for each leaf on the current process (locations must be >= 0, enforced
339c4e6a40aSLawrence Mitchell during setup in debug mode)
34095fce210SBarry Smith -  remotemode - copy mode for iremote
34195fce210SBarry Smith 
34295fce210SBarry Smith    Level: intermediate
34395fce210SBarry Smith 
34495452b02SPatrick Sanan    Notes:
34595452b02SPatrick Sanan     In Fortran you must use PETSC_COPY_VALUES for localmode and remotemode
34638ab3f8aSBarry Smith 
3472ad1e87fSLisandro Dalcin    Developers Note: Local indices which are the identity permutation in the range [0,nleaves) are discarded as they
3482ad1e87fSLisandro Dalcin    encode contiguous storage. In such case, if localmode is PETSC_OWN_POINTER, the memory is deallocated as it is not
3492ad1e87fSLisandro Dalcin    needed
3502ad1e87fSLisandro Dalcin 
351c4e6a40aSLawrence Mitchell    Developers Note: This object does not necessarily encode a true star forest in the graph theoretic sense, since leaf
352c4e6a40aSLawrence Mitchell    indices are not required to be unique. Some functions, however, rely on unique leaf indices (checked in debug mode).
353c4e6a40aSLawrence Mitchell 
35495fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFGetGraph()
35595fce210SBarry Smith @*/
35695fce210SBarry Smith PetscErrorCode PetscSFSetGraph(PetscSF sf,PetscInt nroots,PetscInt nleaves,const PetscInt *ilocal,PetscCopyMode localmode,const PetscSFNode *iremote,PetscCopyMode remotemode)
35795fce210SBarry Smith {
35895fce210SBarry Smith   PetscErrorCode ierr;
35995fce210SBarry Smith 
36095fce210SBarry Smith   PetscFunctionBegin;
36195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
36229046d53SLisandro Dalcin   if (nleaves > 0 && ilocal) PetscValidIntPointer(ilocal,4);
36329046d53SLisandro Dalcin   if (nleaves > 0) PetscValidPointer(iremote,6);
36429046d53SLisandro Dalcin   if (nroots  < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nroots %D, cannot be negative",nroots);
36595fce210SBarry Smith   if (nleaves < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nleaves %D, cannot be negative",nleaves);
36629046d53SLisandro Dalcin 
3672a67d2daSStefano Zampini   if (sf->nroots >= 0) { /* Reset only if graph already set */
36895fce210SBarry Smith     ierr = PetscSFReset(sf);CHKERRQ(ierr);
3692a67d2daSStefano Zampini   }
3702a67d2daSStefano Zampini 
37129046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr);
37229046d53SLisandro Dalcin 
37395fce210SBarry Smith   sf->nroots  = nroots;
37495fce210SBarry Smith   sf->nleaves = nleaves;
37529046d53SLisandro Dalcin 
37629046d53SLisandro Dalcin   if (nleaves && ilocal) {
37721c688dcSJed Brown     PetscInt i;
37829046d53SLisandro Dalcin     PetscInt minleaf = PETSC_MAX_INT;
37929046d53SLisandro Dalcin     PetscInt maxleaf = PETSC_MIN_INT;
3802ad1e87fSLisandro Dalcin     int      contiguous = 1;
38129046d53SLisandro Dalcin     for (i=0; i<nleaves; i++) {
38229046d53SLisandro Dalcin       minleaf = PetscMin(minleaf,ilocal[i]);
38329046d53SLisandro Dalcin       maxleaf = PetscMax(maxleaf,ilocal[i]);
3842ad1e87fSLisandro Dalcin       contiguous &= (ilocal[i] == i);
38529046d53SLisandro Dalcin     }
38629046d53SLisandro Dalcin     sf->minleaf = minleaf;
38729046d53SLisandro Dalcin     sf->maxleaf = maxleaf;
3882ad1e87fSLisandro Dalcin     if (contiguous) {
3892ad1e87fSLisandro Dalcin       if (localmode == PETSC_OWN_POINTER) {
3902ad1e87fSLisandro Dalcin         ierr = PetscFree(ilocal);CHKERRQ(ierr);
3912ad1e87fSLisandro Dalcin       }
3922ad1e87fSLisandro Dalcin       ilocal = NULL;
3932ad1e87fSLisandro Dalcin     }
39429046d53SLisandro Dalcin   } else {
39529046d53SLisandro Dalcin     sf->minleaf = 0;
39629046d53SLisandro Dalcin     sf->maxleaf = nleaves - 1;
39729046d53SLisandro Dalcin   }
39829046d53SLisandro Dalcin 
39929046d53SLisandro Dalcin   if (ilocal) {
40095fce210SBarry Smith     switch (localmode) {
40195fce210SBarry Smith     case PETSC_COPY_VALUES:
402785e854fSJed Brown       ierr = PetscMalloc1(nleaves,&sf->mine_alloc);CHKERRQ(ierr);
403580bdb30SBarry Smith       ierr = PetscArraycpy(sf->mine_alloc,ilocal,nleaves);CHKERRQ(ierr);
40495fce210SBarry Smith       sf->mine = sf->mine_alloc;
40595fce210SBarry Smith       break;
40695fce210SBarry Smith     case PETSC_OWN_POINTER:
40795fce210SBarry Smith       sf->mine_alloc = (PetscInt*)ilocal;
40895fce210SBarry Smith       sf->mine       = sf->mine_alloc;
40995fce210SBarry Smith       break;
41095fce210SBarry Smith     case PETSC_USE_POINTER:
41129046d53SLisandro Dalcin       sf->mine_alloc = NULL;
41295fce210SBarry Smith       sf->mine       = (PetscInt*)ilocal;
41395fce210SBarry Smith       break;
41495fce210SBarry Smith     default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown localmode");
41595fce210SBarry Smith     }
41695fce210SBarry Smith   }
41729046d53SLisandro Dalcin 
41895fce210SBarry Smith   switch (remotemode) {
41995fce210SBarry Smith   case PETSC_COPY_VALUES:
420785e854fSJed Brown     ierr = PetscMalloc1(nleaves,&sf->remote_alloc);CHKERRQ(ierr);
421580bdb30SBarry Smith     ierr = PetscArraycpy(sf->remote_alloc,iremote,nleaves);CHKERRQ(ierr);
42295fce210SBarry Smith     sf->remote = sf->remote_alloc;
42395fce210SBarry Smith     break;
42495fce210SBarry Smith   case PETSC_OWN_POINTER:
42595fce210SBarry Smith     sf->remote_alloc = (PetscSFNode*)iremote;
42695fce210SBarry Smith     sf->remote       = sf->remote_alloc;
42795fce210SBarry Smith     break;
42895fce210SBarry Smith   case PETSC_USE_POINTER:
42929046d53SLisandro Dalcin     sf->remote_alloc = NULL;
43095fce210SBarry Smith     sf->remote       = (PetscSFNode*)iremote;
43195fce210SBarry Smith     break;
43295fce210SBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown remotemode");
43395fce210SBarry Smith   }
43495fce210SBarry Smith 
435563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr);
43629046d53SLisandro Dalcin   sf->graphset = PETSC_TRUE;
43795fce210SBarry Smith   PetscFunctionReturn(0);
43895fce210SBarry Smith }
43995fce210SBarry Smith 
44029046d53SLisandro Dalcin /*@
441dd5b3ca6SJunchao Zhang   PetscSFSetGraphWithPattern - Sets the graph of an SF with a specific pattern
442dd5b3ca6SJunchao Zhang 
443dd5b3ca6SJunchao Zhang   Collective
444dd5b3ca6SJunchao Zhang 
445dd5b3ca6SJunchao Zhang   Input Parameters:
446dd5b3ca6SJunchao Zhang + sf      - The PetscSF
447dd5b3ca6SJunchao Zhang . map     - Layout of roots over all processes (insignificant when pattern is PETSCSF_PATTERN_ALLTOALL)
448dd5b3ca6SJunchao Zhang - pattern - One of PETSCSF_PATTERN_ALLGATHER, PETSCSF_PATTERN_GATHER, PETSCSF_PATTERN_ALLTOALL
449dd5b3ca6SJunchao Zhang 
450dd5b3ca6SJunchao Zhang   Notes:
451dd5b3ca6SJunchao Zhang   It is easier to explain PetscSFPattern using vectors. Suppose we have an MPI vector x and its layout is map.
452dd5b3ca6SJunchao Zhang   n and N are local and global sizes of x respectively.
453dd5b3ca6SJunchao Zhang 
454dd5b3ca6SJunchao Zhang   With PETSCSF_PATTERN_ALLGATHER, the routine creates a graph that if one does Bcast on it, it will copy x to
455dd5b3ca6SJunchao Zhang   sequential vectors y on all ranks.
456dd5b3ca6SJunchao Zhang 
457dd5b3ca6SJunchao Zhang   With PETSCSF_PATTERN_GATHER, the routine creates a graph that if one does Bcast on it, it will copy x to a
458dd5b3ca6SJunchao Zhang   sequential vector y on rank 0.
459dd5b3ca6SJunchao Zhang 
460dd5b3ca6SJunchao Zhang   In above cases, entries of x are roots and entries of y are leaves.
461dd5b3ca6SJunchao Zhang 
462dd5b3ca6SJunchao Zhang   With PETSCSF_PATTERN_ALLTOALL, map is insignificant. Suppose NP is size of sf's communicator. The routine
463dd5b3ca6SJunchao Zhang   creates a graph that every rank has NP leaves and NP roots. On rank i, its leaf j is connected to root i
464dd5b3ca6SJunchao 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
465dd5b3ca6SJunchao Zhang   not mean one can not send multiple items. One just needs to create a new MPI datatype for the mulptiple data
466dd5b3ca6SJunchao Zhang   items with MPI_Type_contiguous() and use that as the <unit> argument in SF routines.
467dd5b3ca6SJunchao Zhang 
468dd5b3ca6SJunchao Zhang   In this case, roots and leaves are symmetric.
469dd5b3ca6SJunchao Zhang 
470dd5b3ca6SJunchao Zhang   Level: intermediate
471dd5b3ca6SJunchao Zhang  @*/
472dd5b3ca6SJunchao Zhang PetscErrorCode PetscSFSetGraphWithPattern(PetscSF sf,PetscLayout map,PetscSFPattern pattern)
473dd5b3ca6SJunchao Zhang {
474dd5b3ca6SJunchao Zhang   MPI_Comm       comm;
475dd5b3ca6SJunchao Zhang   PetscInt       n,N,res[2];
476dd5b3ca6SJunchao Zhang   PetscMPIInt    rank,size;
477dd5b3ca6SJunchao Zhang   PetscSFType    type;
478dd5b3ca6SJunchao Zhang   PetscErrorCode ierr;
479dd5b3ca6SJunchao Zhang 
480dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
481dd5b3ca6SJunchao Zhang   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
482dd5b3ca6SJunchao Zhang   if (pattern < PETSCSF_PATTERN_ALLGATHER || pattern > PETSCSF_PATTERN_ALLTOALL) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Unsupported PetscSFPattern %D\n",pattern);
483dd5b3ca6SJunchao Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
484dd5b3ca6SJunchao Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
485dd5b3ca6SJunchao Zhang 
486dd5b3ca6SJunchao Zhang   if (pattern == PETSCSF_PATTERN_ALLTOALL) {
487dd5b3ca6SJunchao Zhang     type = PETSCSFALLTOALL;
488dd5b3ca6SJunchao Zhang     ierr = PetscLayoutCreate(comm,&sf->map);CHKERRQ(ierr);
489dd5b3ca6SJunchao Zhang     ierr = PetscLayoutSetLocalSize(sf->map,size);CHKERRQ(ierr);
490dd5b3ca6SJunchao Zhang     ierr = PetscLayoutSetSize(sf->map,((PetscInt)size)*size);CHKERRQ(ierr);
491dd5b3ca6SJunchao Zhang     ierr = PetscLayoutSetUp(sf->map);CHKERRQ(ierr);
492dd5b3ca6SJunchao Zhang   } else {
493dd5b3ca6SJunchao Zhang     ierr   = PetscLayoutGetLocalSize(map,&n);CHKERRQ(ierr);
494dd5b3ca6SJunchao Zhang     ierr   = PetscLayoutGetSize(map,&N);CHKERRQ(ierr);
495dd5b3ca6SJunchao Zhang     res[0] = n;
496dd5b3ca6SJunchao Zhang     res[1] = -n;
497dd5b3ca6SJunchao Zhang     /* Check if n are same over all ranks so that we can optimize it */
498dd5b3ca6SJunchao Zhang     ierr   = MPIU_Allreduce(MPI_IN_PLACE,res,2,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr);
499dd5b3ca6SJunchao Zhang     if (res[0] == -res[1]) { /* same n */
500dd5b3ca6SJunchao Zhang       type = (pattern == PETSCSF_PATTERN_ALLGATHER) ? PETSCSFALLGATHER  : PETSCSFGATHER;
501dd5b3ca6SJunchao Zhang     } else {
502dd5b3ca6SJunchao Zhang       type = (pattern == PETSCSF_PATTERN_ALLGATHER) ? PETSCSFALLGATHERV : PETSCSFGATHERV;
503dd5b3ca6SJunchao Zhang     }
504dd5b3ca6SJunchao Zhang     ierr = PetscLayoutReference(map,&sf->map);CHKERRQ(ierr);
505dd5b3ca6SJunchao Zhang   }
506dd5b3ca6SJunchao Zhang   ierr = PetscSFSetType(sf,type);CHKERRQ(ierr);
507dd5b3ca6SJunchao Zhang 
508dd5b3ca6SJunchao Zhang   sf->pattern = pattern;
509dd5b3ca6SJunchao Zhang   sf->mine    = NULL; /* Contiguous */
510dd5b3ca6SJunchao Zhang 
511dd5b3ca6SJunchao Zhang   /* Set nleaves, nroots here in case user calls PetscSFGetGraph, which is legal to call even before PetscSFSetUp is called.
512dd5b3ca6SJunchao Zhang      Also set other easy stuff.
513dd5b3ca6SJunchao Zhang    */
514dd5b3ca6SJunchao Zhang   if (pattern == PETSCSF_PATTERN_ALLGATHER) {
515dd5b3ca6SJunchao Zhang     sf->nleaves      = N;
516dd5b3ca6SJunchao Zhang     sf->nroots       = n;
517dd5b3ca6SJunchao Zhang     sf->nranks       = size;
518dd5b3ca6SJunchao Zhang     sf->minleaf      = 0;
519dd5b3ca6SJunchao Zhang     sf->maxleaf      = N - 1;
520dd5b3ca6SJunchao Zhang   } else if (pattern == PETSCSF_PATTERN_GATHER) {
521dd5b3ca6SJunchao Zhang     sf->nleaves      = rank ? 0 : N;
522dd5b3ca6SJunchao Zhang     sf->nroots       = n;
523dd5b3ca6SJunchao Zhang     sf->nranks       = rank ? 0 : size;
524dd5b3ca6SJunchao Zhang     sf->minleaf      = 0;
525dd5b3ca6SJunchao Zhang     sf->maxleaf      = rank ? -1 : N - 1;
526dd5b3ca6SJunchao Zhang   } else if (pattern == PETSCSF_PATTERN_ALLTOALL) {
527dd5b3ca6SJunchao Zhang     sf->nleaves      = size;
528dd5b3ca6SJunchao Zhang     sf->nroots       = size;
529dd5b3ca6SJunchao Zhang     sf->nranks       = size;
530dd5b3ca6SJunchao Zhang     sf->minleaf      = 0;
531dd5b3ca6SJunchao Zhang     sf->maxleaf      = size - 1;
532dd5b3ca6SJunchao Zhang   }
533dd5b3ca6SJunchao Zhang   sf->ndranks  = 0; /* We do not need to separate out distinguished ranks for patterned graphs to improve communication performance */
534dd5b3ca6SJunchao Zhang   sf->graphset = PETSC_TRUE;
535dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
536dd5b3ca6SJunchao Zhang }
537dd5b3ca6SJunchao Zhang 
538dd5b3ca6SJunchao Zhang /*@
53995fce210SBarry Smith    PetscSFCreateInverseSF - given a PetscSF in which all vertices have degree 1, creates the inverse map
54095fce210SBarry Smith 
54195fce210SBarry Smith    Collective
54295fce210SBarry Smith 
54395fce210SBarry Smith    Input Arguments:
544dd5b3ca6SJunchao Zhang 
54595fce210SBarry Smith .  sf - star forest to invert
54695fce210SBarry Smith 
54795fce210SBarry Smith    Output Arguments:
54895fce210SBarry Smith .  isf - inverse of sf
54995fce210SBarry Smith    Level: advanced
55095fce210SBarry Smith 
55195fce210SBarry Smith    Notes:
55295fce210SBarry Smith    All roots must have degree 1.
55395fce210SBarry Smith 
55495fce210SBarry Smith    The local space may be a permutation, but cannot be sparse.
55595fce210SBarry Smith 
55695fce210SBarry Smith .seealso: PetscSFSetGraph()
55795fce210SBarry Smith @*/
55895fce210SBarry Smith PetscErrorCode PetscSFCreateInverseSF(PetscSF sf,PetscSF *isf)
55995fce210SBarry Smith {
56095fce210SBarry Smith   PetscErrorCode ierr;
56195fce210SBarry Smith   PetscMPIInt    rank;
56295fce210SBarry Smith   PetscInt       i,nroots,nleaves,maxlocal,count,*newilocal;
56395fce210SBarry Smith   const PetscInt *ilocal;
56495fce210SBarry Smith   PetscSFNode    *roots,*leaves;
56595fce210SBarry Smith 
56695fce210SBarry Smith   PetscFunctionBegin;
56729046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
56829046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
56929046d53SLisandro Dalcin   PetscValidPointer(isf,2);
57029046d53SLisandro Dalcin 
57195fce210SBarry Smith   ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);CHKERRQ(ierr);
57229046d53SLisandro Dalcin   maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */
57329046d53SLisandro Dalcin 
57429046d53SLisandro Dalcin   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
575ae9aee6dSMatthew G. Knepley   ierr = PetscMalloc2(nroots,&roots,maxlocal,&leaves);CHKERRQ(ierr);
576ae9aee6dSMatthew G. Knepley   for (i=0; i<maxlocal; i++) {
57795fce210SBarry Smith     leaves[i].rank  = rank;
57895fce210SBarry Smith     leaves[i].index = i;
57995fce210SBarry Smith   }
58095fce210SBarry Smith   for (i=0; i <nroots; i++) {
58195fce210SBarry Smith     roots[i].rank  = -1;
58295fce210SBarry Smith     roots[i].index = -1;
58395fce210SBarry Smith   }
5848bfbc91cSJed Brown   ierr = PetscSFReduceBegin(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr);
5858bfbc91cSJed Brown   ierr = PetscSFReduceEnd(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr);
58695fce210SBarry Smith 
58795fce210SBarry Smith   /* Check whether our leaves are sparse */
58895fce210SBarry Smith   for (i=0,count=0; i<nroots; i++) if (roots[i].rank >= 0) count++;
58995fce210SBarry Smith   if (count == nroots) newilocal = NULL;
59095fce210SBarry Smith   else {                        /* Index for sparse leaves and compact "roots" array (which is to become our leaves). */
591785e854fSJed Brown     ierr = PetscMalloc1(count,&newilocal);CHKERRQ(ierr);
59295fce210SBarry Smith     for (i=0,count=0; i<nroots; i++) {
59395fce210SBarry Smith       if (roots[i].rank >= 0) {
59495fce210SBarry Smith         newilocal[count]   = i;
59595fce210SBarry Smith         roots[count].rank  = roots[i].rank;
59695fce210SBarry Smith         roots[count].index = roots[i].index;
59795fce210SBarry Smith         count++;
59895fce210SBarry Smith       }
59995fce210SBarry Smith     }
60095fce210SBarry Smith   }
60195fce210SBarry Smith 
60295fce210SBarry Smith   ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,isf);CHKERRQ(ierr);
60395fce210SBarry Smith   ierr = PetscSFSetGraph(*isf,maxlocal,count,newilocal,PETSC_OWN_POINTER,roots,PETSC_COPY_VALUES);CHKERRQ(ierr);
60495fce210SBarry Smith   ierr = PetscFree2(roots,leaves);CHKERRQ(ierr);
60595fce210SBarry Smith   PetscFunctionReturn(0);
60695fce210SBarry Smith }
60795fce210SBarry Smith 
60895fce210SBarry Smith /*@
60995fce210SBarry Smith    PetscSFDuplicate - duplicate a PetscSF, optionally preserving rank connectivity and graph
61095fce210SBarry Smith 
61195fce210SBarry Smith    Collective
61295fce210SBarry Smith 
61395fce210SBarry Smith    Input Arguments:
61495fce210SBarry Smith +  sf - communication object to duplicate
61595fce210SBarry Smith -  opt - PETSCSF_DUPLICATE_CONFONLY, PETSCSF_DUPLICATE_RANKS, or PETSCSF_DUPLICATE_GRAPH (see PetscSFDuplicateOption)
61695fce210SBarry Smith 
61795fce210SBarry Smith    Output Arguments:
61895fce210SBarry Smith .  newsf - new communication object
61995fce210SBarry Smith 
62095fce210SBarry Smith    Level: beginner
62195fce210SBarry Smith 
62295fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetType(), PetscSFSetGraph()
62395fce210SBarry Smith @*/
62495fce210SBarry Smith PetscErrorCode PetscSFDuplicate(PetscSF sf,PetscSFDuplicateOption opt,PetscSF *newsf)
62595fce210SBarry Smith {
62629046d53SLisandro Dalcin   PetscSFType    type;
62795fce210SBarry Smith   PetscErrorCode ierr;
62895fce210SBarry Smith 
62995fce210SBarry Smith   PetscFunctionBegin;
63029046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
63129046d53SLisandro Dalcin   PetscValidLogicalCollectiveEnum(sf,opt,2);
63229046d53SLisandro Dalcin   PetscValidPointer(newsf,3);
63395fce210SBarry Smith   ierr = PetscSFCreate(PetscObjectComm((PetscObject)sf),newsf);CHKERRQ(ierr);
63429046d53SLisandro Dalcin   ierr = PetscSFGetType(sf,&type);CHKERRQ(ierr);
63529046d53SLisandro Dalcin   if (type) {ierr = PetscSFSetType(*newsf,type);CHKERRQ(ierr);}
63695fce210SBarry Smith   if (opt == PETSCSF_DUPLICATE_GRAPH) {
637dd5b3ca6SJunchao Zhang     PetscSFCheckGraphSet(sf,1);
638dd5b3ca6SJunchao Zhang     if (sf->pattern == PETSCSF_PATTERN_GENERAL) {
63995fce210SBarry Smith       PetscInt          nroots,nleaves;
64095fce210SBarry Smith       const PetscInt    *ilocal;
64195fce210SBarry Smith       const PetscSFNode *iremote;
64295fce210SBarry Smith       ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
64395fce210SBarry Smith       ierr = PetscSFSetGraph(*newsf,nroots,nleaves,ilocal,PETSC_COPY_VALUES,iremote,PETSC_COPY_VALUES);CHKERRQ(ierr);
644dd5b3ca6SJunchao Zhang     } else {
645dd5b3ca6SJunchao Zhang       ierr = PetscSFSetGraphWithPattern(*newsf,sf->map,sf->pattern);CHKERRQ(ierr);
646dd5b3ca6SJunchao Zhang     }
64795fce210SBarry Smith   }
64829046d53SLisandro Dalcin   if (sf->ops->Duplicate) {ierr = (*sf->ops->Duplicate)(sf,opt,*newsf);CHKERRQ(ierr);}
64995fce210SBarry Smith   PetscFunctionReturn(0);
65095fce210SBarry Smith }
65195fce210SBarry Smith 
65295fce210SBarry Smith /*@C
65395fce210SBarry Smith    PetscSFGetGraph - Get the graph specifying a parallel star forest
65495fce210SBarry Smith 
65595fce210SBarry Smith    Not Collective
65695fce210SBarry Smith 
65795fce210SBarry Smith    Input Arguments:
65895fce210SBarry Smith .  sf - star forest
65995fce210SBarry Smith 
66095fce210SBarry Smith    Output Arguments:
66195fce210SBarry Smith +  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
66295fce210SBarry Smith .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
66395fce210SBarry Smith .  ilocal - locations of leaves in leafdata buffers
66495fce210SBarry Smith -  iremote - remote locations of root vertices for each leaf on the current process
66595fce210SBarry Smith 
666373e0d91SLisandro Dalcin    Notes:
667373e0d91SLisandro Dalcin    We are not currently requiring that the graph is set, thus returning nroots=-1 if it has not been set yet
668373e0d91SLisandro Dalcin 
669245d9833Sprj-    When called from Fortran, the returned iremote array is a copy and must be deallocated after use. Consequently, if you
670ca797d7aSLawrence Mitchell    want to update the graph, you must call PetscSFSetGraph after modifying the iremote array.
671ca797d7aSLawrence Mitchell 
67295fce210SBarry Smith    Level: intermediate
67395fce210SBarry Smith 
67495fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph()
67595fce210SBarry Smith @*/
67695fce210SBarry Smith PetscErrorCode PetscSFGetGraph(PetscSF sf,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote)
67795fce210SBarry Smith {
678b8dee149SJunchao Zhang   PetscErrorCode ierr;
67995fce210SBarry Smith 
68095fce210SBarry Smith   PetscFunctionBegin;
68195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
682b8dee149SJunchao Zhang   if (sf->ops->GetGraph) {
683b8dee149SJunchao Zhang     ierr = (sf->ops->GetGraph)(sf,nroots,nleaves,ilocal,iremote);CHKERRQ(ierr);
684b8dee149SJunchao Zhang   } else {
68595fce210SBarry Smith     if (nroots) *nroots = sf->nroots;
68695fce210SBarry Smith     if (nleaves) *nleaves = sf->nleaves;
68795fce210SBarry Smith     if (ilocal) *ilocal = sf->mine;
68895fce210SBarry Smith     if (iremote) *iremote = sf->remote;
689b8dee149SJunchao Zhang   }
69095fce210SBarry Smith   PetscFunctionReturn(0);
69195fce210SBarry Smith }
69295fce210SBarry Smith 
69329046d53SLisandro Dalcin /*@
69495fce210SBarry Smith    PetscSFGetLeafRange - Get the active leaf ranges
69595fce210SBarry Smith 
69695fce210SBarry Smith    Not Collective
69795fce210SBarry Smith 
69895fce210SBarry Smith    Input Arguments:
69995fce210SBarry Smith .  sf - star forest
70095fce210SBarry Smith 
70195fce210SBarry Smith    Output Arguments:
702dd5b3ca6SJunchao Zhang +  minleaf - minimum active leaf on this process. Return 0 if there are no leaves.
703dd5b3ca6SJunchao Zhang -  maxleaf - maximum active leaf on this process. Return -1 if there are no leaves.
70495fce210SBarry Smith 
70595fce210SBarry Smith    Level: developer
70695fce210SBarry Smith 
70795fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph(), PetscSFGetGraph()
70895fce210SBarry Smith @*/
70995fce210SBarry Smith PetscErrorCode PetscSFGetLeafRange(PetscSF sf,PetscInt *minleaf,PetscInt *maxleaf)
71095fce210SBarry Smith {
71195fce210SBarry Smith 
71295fce210SBarry Smith   PetscFunctionBegin;
71395fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
71429046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
71595fce210SBarry Smith   if (minleaf) *minleaf = sf->minleaf;
71695fce210SBarry Smith   if (maxleaf) *maxleaf = sf->maxleaf;
71795fce210SBarry Smith   PetscFunctionReturn(0);
71895fce210SBarry Smith }
71995fce210SBarry Smith 
72095fce210SBarry Smith /*@C
721fe2efc57SMark    PetscSFViewFromOptions - View from Options
722fe2efc57SMark 
723fe2efc57SMark    Collective on PetscSF
724fe2efc57SMark 
725fe2efc57SMark    Input Parameters:
726fe2efc57SMark +  A - the star forest
727736c3998SJose E. Roman .  obj - Optional object
728736c3998SJose E. Roman -  name - command line option
729fe2efc57SMark 
730fe2efc57SMark    Level: intermediate
731fe2efc57SMark .seealso:  PetscSF, PetscSFView, PetscObjectViewFromOptions(), PetscSFCreate()
732fe2efc57SMark @*/
733fe2efc57SMark PetscErrorCode  PetscSFViewFromOptions(PetscSF A,PetscObject obj,const char name[])
734fe2efc57SMark {
735fe2efc57SMark   PetscErrorCode ierr;
736fe2efc57SMark 
737fe2efc57SMark   PetscFunctionBegin;
738fe2efc57SMark   PetscValidHeaderSpecific(A,PETSCSF_CLASSID,1);
739fe2efc57SMark   ierr = PetscObjectViewFromOptions((PetscObject)A,obj,name);CHKERRQ(ierr);
740fe2efc57SMark   PetscFunctionReturn(0);
741fe2efc57SMark }
742fe2efc57SMark 
743fe2efc57SMark /*@C
74495fce210SBarry Smith    PetscSFView - view a star forest
74595fce210SBarry Smith 
74695fce210SBarry Smith    Collective
74795fce210SBarry Smith 
74895fce210SBarry Smith    Input Arguments:
74995fce210SBarry Smith +  sf - star forest
75095fce210SBarry Smith -  viewer - viewer to display graph, for example PETSC_VIEWER_STDOUT_WORLD
75195fce210SBarry Smith 
75295fce210SBarry Smith    Level: beginner
75395fce210SBarry Smith 
75495fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph()
75595fce210SBarry Smith @*/
75695fce210SBarry Smith PetscErrorCode PetscSFView(PetscSF sf,PetscViewer viewer)
75795fce210SBarry Smith {
75895fce210SBarry Smith   PetscErrorCode    ierr;
75995fce210SBarry Smith   PetscBool         iascii;
76095fce210SBarry Smith   PetscViewerFormat format;
76195fce210SBarry Smith 
76295fce210SBarry Smith   PetscFunctionBegin;
76395fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
76495fce210SBarry Smith   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sf),&viewer);CHKERRQ(ierr);}
76595fce210SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
76695fce210SBarry Smith   PetscCheckSameComm(sf,1,viewer,2);
76780153354SVaclav Hapla   if (sf->graphset) {ierr = PetscSFSetUp(sf);CHKERRQ(ierr);}
76895fce210SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
76995fce210SBarry Smith   if (iascii) {
77095fce210SBarry Smith     PetscMPIInt rank;
77181bfa7aaSJed Brown     PetscInt    ii,i,j;
77295fce210SBarry Smith 
773dae58748SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)sf,viewer);CHKERRQ(ierr);
77495fce210SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
77595fce210SBarry Smith     if (sf->ops->View) {ierr = (*sf->ops->View)(sf,viewer);CHKERRQ(ierr);}
776dd5b3ca6SJunchao Zhang     if (sf->pattern == PETSCSF_PATTERN_GENERAL) {
77780153354SVaclav Hapla       if (!sf->graphset) {
77880153354SVaclav Hapla         ierr = PetscViewerASCIIPrintf(viewer,"PetscSFSetGraph() has not been called yet\n");CHKERRQ(ierr);
77980153354SVaclav Hapla         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
78080153354SVaclav Hapla         PetscFunctionReturn(0);
78180153354SVaclav Hapla       }
78295fce210SBarry Smith       ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
7831575c14dSBarry Smith       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
78495fce210SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of roots=%D, leaves=%D, remote ranks=%D\n",rank,sf->nroots,sf->nleaves,sf->nranks);CHKERRQ(ierr);
78595fce210SBarry Smith       for (i=0; i<sf->nleaves; i++) {
78695fce210SBarry 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);
78795fce210SBarry Smith       }
78895fce210SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
78995fce210SBarry Smith       ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
79095fce210SBarry Smith       if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
79181bfa7aaSJed Brown         PetscMPIInt *tmpranks,*perm;
79281bfa7aaSJed Brown         ierr = PetscMalloc2(sf->nranks,&tmpranks,sf->nranks,&perm);CHKERRQ(ierr);
793580bdb30SBarry Smith         ierr = PetscArraycpy(tmpranks,sf->ranks,sf->nranks);CHKERRQ(ierr);
79481bfa7aaSJed Brown         for (i=0; i<sf->nranks; i++) perm[i] = i;
79581bfa7aaSJed Brown         ierr = PetscSortMPIIntWithArray(sf->nranks,tmpranks,perm);CHKERRQ(ierr);
79695fce210SBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Roots referenced by my leaves, by rank\n",rank);CHKERRQ(ierr);
79781bfa7aaSJed Brown         for (ii=0; ii<sf->nranks; ii++) {
79881bfa7aaSJed Brown           i = perm[ii];
7997904a332SBarry Smith           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %d: %D edges\n",rank,sf->ranks[i],sf->roffset[i+1]-sf->roffset[i]);CHKERRQ(ierr);
80095fce210SBarry Smith           for (j=sf->roffset[i]; j<sf->roffset[i+1]; j++) {
80195fce210SBarry Smith             ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d]    %D <- %D\n",rank,sf->rmine[j],sf->rremote[j]);CHKERRQ(ierr);
80295fce210SBarry Smith           }
80395fce210SBarry Smith         }
80481bfa7aaSJed Brown         ierr = PetscFree2(tmpranks,perm);CHKERRQ(ierr);
80595fce210SBarry Smith       }
80695fce210SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
8071575c14dSBarry Smith       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
808dd5b3ca6SJunchao Zhang     }
80995fce210SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
81095fce210SBarry Smith   }
81195fce210SBarry Smith   PetscFunctionReturn(0);
81295fce210SBarry Smith }
81395fce210SBarry Smith 
81495fce210SBarry Smith /*@C
815dec1416fSJunchao Zhang    PetscSFGetRootRanks - Get root ranks and number of vertices referenced by leaves on this process
81695fce210SBarry Smith 
81795fce210SBarry Smith    Not Collective
81895fce210SBarry Smith 
81995fce210SBarry Smith    Input Arguments:
82095fce210SBarry Smith .  sf - star forest
82195fce210SBarry Smith 
82295fce210SBarry Smith    Output Arguments:
82395fce210SBarry Smith +  nranks - number of ranks referenced by local part
82495fce210SBarry Smith .  ranks - array of ranks
82595fce210SBarry Smith .  roffset - offset in rmine/rremote for each rank (length nranks+1)
82695fce210SBarry Smith .  rmine - concatenated array holding local indices referencing each remote rank
82795fce210SBarry Smith -  rremote - concatenated array holding remote indices referenced for each remote rank
82895fce210SBarry Smith 
82995fce210SBarry Smith    Level: developer
83095fce210SBarry Smith 
831dec1416fSJunchao Zhang .seealso: PetscSFGetLeafRanks()
83295fce210SBarry Smith @*/
833dec1416fSJunchao Zhang PetscErrorCode PetscSFGetRootRanks(PetscSF sf,PetscInt *nranks,const PetscMPIInt **ranks,const PetscInt **roffset,const PetscInt **rmine,const PetscInt **rremote)
83495fce210SBarry Smith {
835dec1416fSJunchao Zhang   PetscErrorCode ierr;
83695fce210SBarry Smith 
83795fce210SBarry Smith   PetscFunctionBegin;
83895fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
83929046d53SLisandro Dalcin   if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining ranks");
840dec1416fSJunchao Zhang   if (sf->ops->GetRootRanks) {
841dec1416fSJunchao Zhang     ierr = (sf->ops->GetRootRanks)(sf,nranks,ranks,roffset,rmine,rremote);CHKERRQ(ierr);
842dec1416fSJunchao Zhang   } else {
843dec1416fSJunchao Zhang     /* The generic implementation */
84495fce210SBarry Smith     if (nranks)  *nranks  = sf->nranks;
84595fce210SBarry Smith     if (ranks)   *ranks   = sf->ranks;
84695fce210SBarry Smith     if (roffset) *roffset = sf->roffset;
84795fce210SBarry Smith     if (rmine)   *rmine   = sf->rmine;
84895fce210SBarry Smith     if (rremote) *rremote = sf->rremote;
849dec1416fSJunchao Zhang   }
85095fce210SBarry Smith   PetscFunctionReturn(0);
85195fce210SBarry Smith }
85295fce210SBarry Smith 
8538750ddebSJunchao Zhang /*@C
8548750ddebSJunchao Zhang    PetscSFGetLeafRanks - Get leaf ranks referencing roots on this process
8558750ddebSJunchao Zhang 
8568750ddebSJunchao Zhang    Not Collective
8578750ddebSJunchao Zhang 
8588750ddebSJunchao Zhang    Input Arguments:
8598750ddebSJunchao Zhang .  sf - star forest
8608750ddebSJunchao Zhang 
8618750ddebSJunchao Zhang    Output Arguments:
8628750ddebSJunchao Zhang +  niranks - number of leaf ranks referencing roots on this process
8638750ddebSJunchao Zhang .  iranks - array of ranks
8648750ddebSJunchao Zhang .  ioffset - offset in irootloc for each rank (length niranks+1)
8658750ddebSJunchao Zhang -  irootloc - concatenated array holding local indices of roots referenced by each leaf rank
8668750ddebSJunchao Zhang 
8678750ddebSJunchao Zhang    Level: developer
8688750ddebSJunchao Zhang 
869dec1416fSJunchao Zhang .seealso: PetscSFGetRootRanks()
8708750ddebSJunchao Zhang @*/
8718750ddebSJunchao Zhang PetscErrorCode PetscSFGetLeafRanks(PetscSF sf,PetscInt *niranks,const PetscMPIInt **iranks,const PetscInt **ioffset,const PetscInt **irootloc)
8728750ddebSJunchao Zhang {
8738750ddebSJunchao Zhang   PetscErrorCode ierr;
8748750ddebSJunchao Zhang 
8758750ddebSJunchao Zhang   PetscFunctionBegin;
8768750ddebSJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
8778750ddebSJunchao Zhang   if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining ranks");
8788750ddebSJunchao Zhang   if (sf->ops->GetLeafRanks) {
8798750ddebSJunchao Zhang     ierr = (sf->ops->GetLeafRanks)(sf,niranks,iranks,ioffset,irootloc);CHKERRQ(ierr);
8808750ddebSJunchao Zhang   } else {
8818750ddebSJunchao Zhang     PetscSFType type;
8828750ddebSJunchao Zhang     ierr = PetscSFGetType(sf,&type);CHKERRQ(ierr);
8838750ddebSJunchao Zhang     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSFGetLeafRanks() is not supported on this StarForest type: %s", type);
8848750ddebSJunchao Zhang   }
8858750ddebSJunchao Zhang   PetscFunctionReturn(0);
8868750ddebSJunchao Zhang }
8878750ddebSJunchao Zhang 
888b5a8e515SJed Brown static PetscBool InList(PetscMPIInt needle,PetscMPIInt n,const PetscMPIInt *list) {
889b5a8e515SJed Brown   PetscInt i;
890b5a8e515SJed Brown   for (i=0; i<n; i++) {
891b5a8e515SJed Brown     if (needle == list[i]) return PETSC_TRUE;
892b5a8e515SJed Brown   }
893b5a8e515SJed Brown   return PETSC_FALSE;
894b5a8e515SJed Brown }
895b5a8e515SJed Brown 
89695fce210SBarry Smith /*@C
89721c688dcSJed Brown    PetscSFSetUpRanks - Set up data structures associated with ranks; this is for internal use by PetscSF implementations.
89821c688dcSJed Brown 
89921c688dcSJed Brown    Collective
90021c688dcSJed Brown 
90121c688dcSJed Brown    Input Arguments:
902b5a8e515SJed Brown +  sf - PetscSF to set up; PetscSFSetGraph() must have been called
903b5a8e515SJed Brown -  dgroup - MPI_Group of ranks to be distinguished (e.g., for self or shared memory exchange)
90421c688dcSJed Brown 
90521c688dcSJed Brown    Level: developer
90621c688dcSJed Brown 
907dec1416fSJunchao Zhang .seealso: PetscSFGetRootRanks()
90821c688dcSJed Brown @*/
909b5a8e515SJed Brown PetscErrorCode PetscSFSetUpRanks(PetscSF sf,MPI_Group dgroup)
91021c688dcSJed Brown {
91121c688dcSJed Brown   PetscErrorCode     ierr;
91221c688dcSJed Brown   PetscTable         table;
91321c688dcSJed Brown   PetscTablePosition pos;
914b5a8e515SJed Brown   PetscMPIInt        size,groupsize,*groupranks;
915247e8311SStefano Zampini   PetscInt           *rcount,*ranks;
916247e8311SStefano Zampini   PetscInt           i, irank = -1,orank = -1;
91721c688dcSJed Brown 
91821c688dcSJed Brown   PetscFunctionBegin;
91921c688dcSJed Brown   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
92029046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
92121c688dcSJed Brown   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr);
92221c688dcSJed Brown   ierr = PetscTableCreate(10,size,&table);CHKERRQ(ierr);
92321c688dcSJed Brown   for (i=0; i<sf->nleaves; i++) {
92421c688dcSJed Brown     /* Log 1-based rank */
92521c688dcSJed Brown     ierr = PetscTableAdd(table,sf->remote[i].rank+1,1,ADD_VALUES);CHKERRQ(ierr);
92621c688dcSJed Brown   }
92721c688dcSJed Brown   ierr = PetscTableGetCount(table,&sf->nranks);CHKERRQ(ierr);
92821c688dcSJed Brown   ierr = PetscMalloc4(sf->nranks,&sf->ranks,sf->nranks+1,&sf->roffset,sf->nleaves,&sf->rmine,sf->nleaves,&sf->rremote);CHKERRQ(ierr);
92921c688dcSJed Brown   ierr = PetscMalloc2(sf->nranks,&rcount,sf->nranks,&ranks);CHKERRQ(ierr);
93021c688dcSJed Brown   ierr = PetscTableGetHeadPosition(table,&pos);CHKERRQ(ierr);
93121c688dcSJed Brown   for (i=0; i<sf->nranks; i++) {
93221c688dcSJed Brown     ierr = PetscTableGetNext(table,&pos,&ranks[i],&rcount[i]);CHKERRQ(ierr);
93321c688dcSJed Brown     ranks[i]--;             /* Convert back to 0-based */
93421c688dcSJed Brown   }
93521c688dcSJed Brown   ierr = PetscTableDestroy(&table);CHKERRQ(ierr);
936b5a8e515SJed Brown 
937b5a8e515SJed Brown   /* We expect that dgroup is reliably "small" while nranks could be large */
938b5a8e515SJed Brown   {
9397fb8a5e4SKarl Rupp     MPI_Group group = MPI_GROUP_NULL;
940b5a8e515SJed Brown     PetscMPIInt *dgroupranks;
941b5a8e515SJed Brown     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
942b5a8e515SJed Brown     ierr = MPI_Group_size(dgroup,&groupsize);CHKERRQ(ierr);
943b5a8e515SJed Brown     ierr = PetscMalloc1(groupsize,&dgroupranks);CHKERRQ(ierr);
944b5a8e515SJed Brown     ierr = PetscMalloc1(groupsize,&groupranks);CHKERRQ(ierr);
945b5a8e515SJed Brown     for (i=0; i<groupsize; i++) dgroupranks[i] = i;
946f3fc9a17SJed Brown     if (groupsize) {ierr = MPI_Group_translate_ranks(dgroup,groupsize,dgroupranks,group,groupranks);CHKERRQ(ierr);}
947b5a8e515SJed Brown     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
948b5a8e515SJed Brown     ierr = PetscFree(dgroupranks);CHKERRQ(ierr);
949b5a8e515SJed Brown   }
950b5a8e515SJed Brown 
951b5a8e515SJed Brown   /* Partition ranks[] into distinguished (first sf->ndranks) followed by non-distinguished */
952b5a8e515SJed Brown   for (sf->ndranks=0,i=sf->nranks; sf->ndranks<i; ) {
953b5a8e515SJed Brown     for (i--; sf->ndranks<i; i--) { /* Scan i backward looking for distinguished rank */
954b5a8e515SJed Brown       if (InList(ranks[i],groupsize,groupranks)) break;
955b5a8e515SJed Brown     }
956b5a8e515SJed Brown     for ( ; sf->ndranks<=i; sf->ndranks++) { /* Scan sf->ndranks forward looking for non-distinguished rank */
957b5a8e515SJed Brown       if (!InList(ranks[sf->ndranks],groupsize,groupranks)) break;
958b5a8e515SJed Brown     }
959b5a8e515SJed Brown     if (sf->ndranks < i) {                         /* Swap ranks[sf->ndranks] with ranks[i] */
960b5a8e515SJed Brown       PetscInt tmprank,tmpcount;
961247e8311SStefano Zampini 
962b5a8e515SJed Brown       tmprank             = ranks[i];
963b5a8e515SJed Brown       tmpcount            = rcount[i];
964b5a8e515SJed Brown       ranks[i]            = ranks[sf->ndranks];
965b5a8e515SJed Brown       rcount[i]           = rcount[sf->ndranks];
966b5a8e515SJed Brown       ranks[sf->ndranks]  = tmprank;
967b5a8e515SJed Brown       rcount[sf->ndranks] = tmpcount;
968b5a8e515SJed Brown       sf->ndranks++;
969b5a8e515SJed Brown     }
970b5a8e515SJed Brown   }
971b5a8e515SJed Brown   ierr = PetscFree(groupranks);CHKERRQ(ierr);
972b5a8e515SJed Brown   ierr = PetscSortIntWithArray(sf->ndranks,ranks,rcount);CHKERRQ(ierr);
973b5a8e515SJed Brown   ierr = PetscSortIntWithArray(sf->nranks-sf->ndranks,ranks+sf->ndranks,rcount+sf->ndranks);CHKERRQ(ierr);
97421c688dcSJed Brown   sf->roffset[0] = 0;
97521c688dcSJed Brown   for (i=0; i<sf->nranks; i++) {
97621c688dcSJed Brown     ierr = PetscMPIIntCast(ranks[i],sf->ranks+i);CHKERRQ(ierr);
97721c688dcSJed Brown     sf->roffset[i+1] = sf->roffset[i] + rcount[i];
97821c688dcSJed Brown     rcount[i]        = 0;
97921c688dcSJed Brown   }
980247e8311SStefano Zampini   for (i=0, irank = -1, orank = -1; i<sf->nleaves; i++) {
981247e8311SStefano Zampini     /* short circuit */
982247e8311SStefano Zampini     if (orank != sf->remote[i].rank) {
98321c688dcSJed Brown       /* Search for index of iremote[i].rank in sf->ranks */
984b5a8e515SJed Brown       ierr = PetscFindMPIInt(sf->remote[i].rank,sf->ndranks,sf->ranks,&irank);CHKERRQ(ierr);
985b5a8e515SJed Brown       if (irank < 0) {
986b5a8e515SJed Brown         ierr = PetscFindMPIInt(sf->remote[i].rank,sf->nranks-sf->ndranks,sf->ranks+sf->ndranks,&irank);CHKERRQ(ierr);
987b5a8e515SJed Brown         if (irank >= 0) irank += sf->ndranks;
98821c688dcSJed Brown       }
989247e8311SStefano Zampini       orank = sf->remote[i].rank;
990247e8311SStefano Zampini     }
991b5a8e515SJed Brown     if (irank < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Could not find rank %D in array",sf->remote[i].rank);
99221c688dcSJed Brown     sf->rmine[sf->roffset[irank] + rcount[irank]]   = sf->mine ? sf->mine[i] : i;
99321c688dcSJed Brown     sf->rremote[sf->roffset[irank] + rcount[irank]] = sf->remote[i].index;
99421c688dcSJed Brown     rcount[irank]++;
99521c688dcSJed Brown   }
99621c688dcSJed Brown   ierr = PetscFree2(rcount,ranks);CHKERRQ(ierr);
99721c688dcSJed Brown   PetscFunctionReturn(0);
99821c688dcSJed Brown }
99921c688dcSJed Brown 
100021c688dcSJed Brown /*@C
100195fce210SBarry Smith    PetscSFGetGroups - gets incoming and outgoing process groups
100295fce210SBarry Smith 
100395fce210SBarry Smith    Collective
100495fce210SBarry Smith 
100595fce210SBarry Smith    Input Argument:
100695fce210SBarry Smith .  sf - star forest
100795fce210SBarry Smith 
100895fce210SBarry Smith    Output Arguments:
100995fce210SBarry Smith +  incoming - group of origin processes for incoming edges (leaves that reference my roots)
101095fce210SBarry Smith -  outgoing - group of destination processes for outgoing edges (roots that I reference)
101195fce210SBarry Smith 
101295fce210SBarry Smith    Level: developer
101395fce210SBarry Smith 
101495fce210SBarry Smith .seealso: PetscSFGetWindow(), PetscSFRestoreWindow()
101595fce210SBarry Smith @*/
101695fce210SBarry Smith PetscErrorCode PetscSFGetGroups(PetscSF sf,MPI_Group *incoming,MPI_Group *outgoing)
101795fce210SBarry Smith {
101895fce210SBarry Smith   PetscErrorCode ierr;
10197fb8a5e4SKarl Rupp   MPI_Group      group = MPI_GROUP_NULL;
102095fce210SBarry Smith 
102195fce210SBarry Smith   PetscFunctionBegin;
102244ee17edSStefano Zampini   if (sf->nranks < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUpRanks() before obtaining groups");
102395fce210SBarry Smith   if (sf->ingroup == MPI_GROUP_NULL) {
102495fce210SBarry Smith     PetscInt       i;
102595fce210SBarry Smith     const PetscInt *indegree;
102695fce210SBarry Smith     PetscMPIInt    rank,*outranks,*inranks;
102795fce210SBarry Smith     PetscSFNode    *remote;
102895fce210SBarry Smith     PetscSF        bgcount;
102995fce210SBarry Smith 
103095fce210SBarry Smith     /* Compute the number of incoming ranks */
1031785e854fSJed Brown     ierr = PetscMalloc1(sf->nranks,&remote);CHKERRQ(ierr);
103295fce210SBarry Smith     for (i=0; i<sf->nranks; i++) {
103395fce210SBarry Smith       remote[i].rank  = sf->ranks[i];
103495fce210SBarry Smith       remote[i].index = 0;
103595fce210SBarry Smith     }
103695fce210SBarry Smith     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,&bgcount);CHKERRQ(ierr);
103795fce210SBarry Smith     ierr = PetscSFSetGraph(bgcount,1,sf->nranks,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
103895fce210SBarry Smith     ierr = PetscSFComputeDegreeBegin(bgcount,&indegree);CHKERRQ(ierr);
103995fce210SBarry Smith     ierr = PetscSFComputeDegreeEnd(bgcount,&indegree);CHKERRQ(ierr);
104095fce210SBarry Smith     /* Enumerate the incoming ranks */
1041dcca6d9dSJed Brown     ierr = PetscMalloc2(indegree[0],&inranks,sf->nranks,&outranks);CHKERRQ(ierr);
104295fce210SBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
104395fce210SBarry Smith     for (i=0; i<sf->nranks; i++) outranks[i] = rank;
104495fce210SBarry Smith     ierr = PetscSFGatherBegin(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr);
104595fce210SBarry Smith     ierr = PetscSFGatherEnd(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr);
104695fce210SBarry Smith     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
104795fce210SBarry Smith     ierr = MPI_Group_incl(group,indegree[0],inranks,&sf->ingroup);CHKERRQ(ierr);
104895fce210SBarry Smith     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
104995fce210SBarry Smith     ierr = PetscFree2(inranks,outranks);CHKERRQ(ierr);
105095fce210SBarry Smith     ierr = PetscSFDestroy(&bgcount);CHKERRQ(ierr);
105195fce210SBarry Smith   }
105295fce210SBarry Smith   *incoming = sf->ingroup;
105395fce210SBarry Smith 
105495fce210SBarry Smith   if (sf->outgroup == MPI_GROUP_NULL) {
105595fce210SBarry Smith     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
105695fce210SBarry Smith     ierr = MPI_Group_incl(group,sf->nranks,sf->ranks,&sf->outgroup);CHKERRQ(ierr);
105795fce210SBarry Smith     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
105895fce210SBarry Smith   }
105995fce210SBarry Smith   *outgoing = sf->outgroup;
106095fce210SBarry Smith   PetscFunctionReturn(0);
106195fce210SBarry Smith }
106295fce210SBarry Smith 
106329046d53SLisandro Dalcin /*@
106495fce210SBarry Smith    PetscSFGetMultiSF - gets the inner SF implemeting gathers and scatters
106595fce210SBarry Smith 
106695fce210SBarry Smith    Collective
106795fce210SBarry Smith 
106895fce210SBarry Smith    Input Argument:
106995fce210SBarry Smith .  sf - star forest that may contain roots with 0 or with more than 1 vertex
107095fce210SBarry Smith 
107195fce210SBarry Smith    Output Arguments:
107295fce210SBarry Smith .  multi - star forest with split roots, such that each root has degree exactly 1
107395fce210SBarry Smith 
107495fce210SBarry Smith    Level: developer
107595fce210SBarry Smith 
107695fce210SBarry Smith    Notes:
107795fce210SBarry Smith 
107895fce210SBarry Smith    In most cases, users should use PetscSFGatherBegin() and PetscSFScatterBegin() instead of manipulating multi
107995fce210SBarry Smith    directly. Since multi satisfies the stronger condition that each entry in the global space has exactly one incoming
108095fce210SBarry Smith    edge, it is a candidate for future optimization that might involve its removal.
108195fce210SBarry Smith 
1082673100f5SVaclav Hapla .seealso: PetscSFSetGraph(), PetscSFGatherBegin(), PetscSFScatterBegin(), PetscSFComputeMultiRootOriginalNumbering()
108395fce210SBarry Smith @*/
108495fce210SBarry Smith PetscErrorCode PetscSFGetMultiSF(PetscSF sf,PetscSF *multi)
108595fce210SBarry Smith {
108695fce210SBarry Smith   PetscErrorCode ierr;
108795fce210SBarry Smith 
108895fce210SBarry Smith   PetscFunctionBegin;
108995fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
109095fce210SBarry Smith   PetscValidPointer(multi,2);
109195fce210SBarry Smith   if (sf->nroots < 0) {         /* Graph has not been set yet; why do we need this? */
109295fce210SBarry Smith     ierr   = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr);
109395fce210SBarry Smith     *multi = sf->multi;
109495fce210SBarry Smith     PetscFunctionReturn(0);
109595fce210SBarry Smith   }
109695fce210SBarry Smith   if (!sf->multi) {
109795fce210SBarry Smith     const PetscInt *indegree;
10989837ea96SMatthew G. Knepley     PetscInt       i,*inoffset,*outones,*outoffset,maxlocal;
109995fce210SBarry Smith     PetscSFNode    *remote;
110029046d53SLisandro Dalcin     maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */
110195fce210SBarry Smith     ierr = PetscSFComputeDegreeBegin(sf,&indegree);CHKERRQ(ierr);
110295fce210SBarry Smith     ierr = PetscSFComputeDegreeEnd(sf,&indegree);CHKERRQ(ierr);
11039837ea96SMatthew G. Knepley     ierr = PetscMalloc3(sf->nroots+1,&inoffset,maxlocal,&outones,maxlocal,&outoffset);CHKERRQ(ierr);
110495fce210SBarry Smith     inoffset[0] = 0;
110595fce210SBarry Smith     for (i=0; i<sf->nroots; i++) inoffset[i+1] = inoffset[i] + indegree[i];
11069837ea96SMatthew G. Knepley     for (i=0; i<maxlocal; i++) outones[i] = 1;
1107dbd2ff41SBarry Smith     ierr = PetscSFFetchAndOpBegin(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr);
1108dbd2ff41SBarry Smith     ierr = PetscSFFetchAndOpEnd(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr);
110995fce210SBarry Smith     for (i=0; i<sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */
111095fce210SBarry Smith #if defined(PETSC_USE_DEBUG)                                 /* Check that the expected number of increments occurred */
111195fce210SBarry Smith     for (i=0; i<sf->nroots; i++) {
111295fce210SBarry Smith       if (inoffset[i] + indegree[i] != inoffset[i+1]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect result after PetscSFFetchAndOp");
111395fce210SBarry Smith     }
111495fce210SBarry Smith #endif
1115785e854fSJed Brown     ierr = PetscMalloc1(sf->nleaves,&remote);CHKERRQ(ierr);
111695fce210SBarry Smith     for (i=0; i<sf->nleaves; i++) {
111795fce210SBarry Smith       remote[i].rank  = sf->remote[i].rank;
111838e7336fSToby Isaac       remote[i].index = outoffset[sf->mine ? sf->mine[i] : i];
111995fce210SBarry Smith     }
112095fce210SBarry Smith     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr);
112101365b40SToby Isaac     ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
112295fce210SBarry Smith     if (sf->rankorder) {        /* Sort the ranks */
112395fce210SBarry Smith       PetscMPIInt rank;
112495fce210SBarry Smith       PetscInt    *inranks,*newoffset,*outranks,*newoutoffset,*tmpoffset,maxdegree;
112595fce210SBarry Smith       PetscSFNode *newremote;
112695fce210SBarry Smith       ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
112795fce210SBarry Smith       for (i=0,maxdegree=0; i<sf->nroots; i++) maxdegree = PetscMax(maxdegree,indegree[i]);
11289837ea96SMatthew G. Knepley       ierr = PetscMalloc5(sf->multi->nroots,&inranks,sf->multi->nroots,&newoffset,maxlocal,&outranks,maxlocal,&newoutoffset,maxdegree,&tmpoffset);CHKERRQ(ierr);
11299837ea96SMatthew G. Knepley       for (i=0; i<maxlocal; i++) outranks[i] = rank;
11308bfbc91cSJed Brown       ierr = PetscSFReduceBegin(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr);
11318bfbc91cSJed Brown       ierr = PetscSFReduceEnd(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr);
113295fce210SBarry Smith       /* Sort the incoming ranks at each vertex, build the inverse map */
113395fce210SBarry Smith       for (i=0; i<sf->nroots; i++) {
113495fce210SBarry Smith         PetscInt j;
113595fce210SBarry Smith         for (j=0; j<indegree[i]; j++) tmpoffset[j] = j;
113695fce210SBarry Smith         ierr = PetscSortIntWithArray(indegree[i],inranks+inoffset[i],tmpoffset);CHKERRQ(ierr);
113795fce210SBarry Smith         for (j=0; j<indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j;
113895fce210SBarry Smith       }
113995fce210SBarry Smith       ierr = PetscSFBcastBegin(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr);
114095fce210SBarry Smith       ierr = PetscSFBcastEnd(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr);
1141785e854fSJed Brown       ierr = PetscMalloc1(sf->nleaves,&newremote);CHKERRQ(ierr);
114295fce210SBarry Smith       for (i=0; i<sf->nleaves; i++) {
114395fce210SBarry Smith         newremote[i].rank  = sf->remote[i].rank;
114401365b40SToby Isaac         newremote[i].index = newoutoffset[sf->mine ? sf->mine[i] : i];
114595fce210SBarry Smith       }
114601365b40SToby Isaac       ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,newremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
114795fce210SBarry Smith       ierr = PetscFree5(inranks,newoffset,outranks,newoutoffset,tmpoffset);CHKERRQ(ierr);
114895fce210SBarry Smith     }
114995fce210SBarry Smith     ierr = PetscFree3(inoffset,outones,outoffset);CHKERRQ(ierr);
115095fce210SBarry Smith   }
115195fce210SBarry Smith   *multi = sf->multi;
115295fce210SBarry Smith   PetscFunctionReturn(0);
115395fce210SBarry Smith }
115495fce210SBarry Smith 
115595fce210SBarry Smith /*@C
115695fce210SBarry Smith    PetscSFCreateEmbeddedSF - removes edges from all but the selected roots, does not remap indices
115795fce210SBarry Smith 
115895fce210SBarry Smith    Collective
115995fce210SBarry Smith 
116095fce210SBarry Smith    Input Arguments:
116195fce210SBarry Smith +  sf - original star forest
1162ba2a7774SJunchao Zhang .  nselected  - number of selected roots on this process
1163ba2a7774SJunchao Zhang -  selected   - indices of the selected roots on this process
116495fce210SBarry Smith 
116595fce210SBarry Smith    Output Arguments:
1166cd620004SJunchao Zhang .  esf - new star forest
116795fce210SBarry Smith 
116895fce210SBarry Smith    Level: advanced
116995fce210SBarry Smith 
117095fce210SBarry Smith    Note:
117195fce210SBarry Smith    To use the new PetscSF, it may be necessary to know the indices of the leaves that are still participating. This can
117295fce210SBarry Smith    be done by calling PetscSFGetGraph().
117395fce210SBarry Smith 
117495fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFGetGraph()
117595fce210SBarry Smith @*/
1176cd620004SJunchao Zhang PetscErrorCode PetscSFCreateEmbeddedSF(PetscSF sf,PetscInt nselected,const PetscInt *selected,PetscSF *esf)
117795fce210SBarry Smith {
1178cd620004SJunchao Zhang   PetscInt          i,j,n,nroots,nleaves,esf_nleaves,*new_ilocal,minleaf,maxleaf,maxlocal;
1179cd620004SJunchao Zhang   const PetscInt    *ilocal;
1180cd620004SJunchao Zhang   signed char       *rootdata,*leafdata,*leafmem;
1181ba2a7774SJunchao Zhang   const PetscSFNode *iremote;
1182f659e5c7SJunchao Zhang   PetscSFNode       *new_iremote;
1183f659e5c7SJunchao Zhang   MPI_Comm          comm;
11840511a646SMatthew G. Knepley   PetscErrorCode    ierr;
118595fce210SBarry Smith 
118695fce210SBarry Smith   PetscFunctionBegin;
118795fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
118829046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
1189ba2a7774SJunchao Zhang   if (nselected) PetscValidPointer(selected,3);
1190cd620004SJunchao Zhang   PetscValidPointer(esf,4);
11910511a646SMatthew G. Knepley 
1192f659e5c7SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1193cd620004SJunchao Zhang 
1194140a1472SStefano Zampini   ierr = PetscLogEventBegin(PETSCSF_EmbedSF,sf,0,0,0);CHKERRQ(ierr);
1195f659e5c7SJunchao Zhang   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
1196cd620004SJunchao Zhang   ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
1197cd620004SJunchao Zhang 
1198cd620004SJunchao Zhang #if defined(PETSC_USE_DEBUG)
1199cd620004SJunchao Zhang   /* Error out if selected[] has dups or  out of range indices */
1200cd620004SJunchao Zhang   {
1201cd620004SJunchao Zhang     PetscBool dups;
1202cd620004SJunchao Zhang     ierr = PetscCheckDupsInt(nselected,selected,&dups);CHKERRQ(ierr);
1203cd620004SJunchao Zhang     if (dups) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"selected[] has dups");
1204cd620004SJunchao Zhang     for (i=0; i<nselected; i++)
1205cd620004SJunchao Zhang       if (selected[i] < 0 || selected[i] >= nroots) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"selected root indice %D is out of [0,%D)",selected[i],nroots);
1206cd620004SJunchao Zhang   }
1207cd620004SJunchao Zhang #endif
1208f659e5c7SJunchao Zhang 
1209f659e5c7SJunchao Zhang   if (sf->ops->CreateEmbeddedSF) {
1210cd620004SJunchao Zhang     ierr = (*sf->ops->CreateEmbeddedSF)(sf,nselected,selected,esf);CHKERRQ(ierr);
1211f659e5c7SJunchao Zhang   } else {
1212cd620004SJunchao Zhang     /* A generic version of creating embedded sf */
1213cd620004SJunchao Zhang     ierr = PetscSFGetLeafRange(sf,&minleaf,&maxleaf);CHKERRQ(ierr);
1214cd620004SJunchao Zhang     maxlocal = maxleaf - minleaf + 1;
1215cd620004SJunchao Zhang     ierr = PetscCalloc2(nroots,&rootdata,maxlocal,&leafmem);CHKERRQ(ierr);
1216cd620004SJunchao Zhang     leafdata = leafmem - minleaf;
1217cd620004SJunchao Zhang     /* Tag selected roots and bcast to leaves */
1218cd620004SJunchao Zhang     for (i=0; i<nselected; i++) rootdata[selected[i]] = 1;
1219cd620004SJunchao Zhang     ierr = PetscSFBcastBegin(sf,MPI_SIGNED_CHAR,rootdata,leafdata);CHKERRQ(ierr);
1220cd620004SJunchao Zhang     ierr = PetscSFBcastEnd(sf,MPI_SIGNED_CHAR,rootdata,leafdata);CHKERRQ(ierr);
1221ba2a7774SJunchao Zhang 
1222cd620004SJunchao Zhang     /* Build esf with leaves that are still connected */
1223cd620004SJunchao Zhang     esf_nleaves = 0;
1224cd620004SJunchao Zhang     for (i=0; i<nleaves; i++) {
1225cd620004SJunchao Zhang       j = ilocal ? ilocal[i] : i;
1226cd620004SJunchao Zhang       /* esf_nleaves += leafdata[j] should work in theory, but failed with SFWindow bugs
1227cd620004SJunchao Zhang          with PetscSFBcast. See https://gitlab.com/petsc/petsc/issues/555
1228cd620004SJunchao Zhang       */
1229cd620004SJunchao Zhang       esf_nleaves += (leafdata[j] ? 1 : 0);
1230cd620004SJunchao Zhang     }
1231cd620004SJunchao Zhang     ierr = PetscMalloc1(esf_nleaves,&new_ilocal);CHKERRQ(ierr);
1232cd620004SJunchao Zhang     ierr = PetscMalloc1(esf_nleaves,&new_iremote);CHKERRQ(ierr);
1233cd620004SJunchao Zhang     for (i=n=0; i<nleaves; i++) {
1234cd620004SJunchao Zhang       j = ilocal ? ilocal[i] : i;
1235cd620004SJunchao Zhang       if (leafdata[j]) {
1236cd620004SJunchao Zhang         new_ilocal[n]        = j;
1237cd620004SJunchao Zhang         new_iremote[n].rank  = iremote[i].rank;
1238cd620004SJunchao Zhang         new_iremote[n].index = iremote[i].index;
1239fc1ede2bSMatthew G. Knepley         ++n;
124095fce210SBarry Smith       }
124195fce210SBarry Smith     }
1242cd620004SJunchao Zhang     ierr = PetscSFCreate(comm,esf);CHKERRQ(ierr);
1243cd620004SJunchao Zhang     ierr = PetscSFSetFromOptions(*esf);CHKERRQ(ierr);
1244cd620004SJunchao Zhang     ierr = PetscSFSetGraph(*esf,nroots,esf_nleaves,new_ilocal,PETSC_OWN_POINTER,new_iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1245cd620004SJunchao Zhang     ierr = PetscFree2(rootdata,leafmem);CHKERRQ(ierr);
1246f659e5c7SJunchao Zhang   }
1247140a1472SStefano Zampini   ierr = PetscLogEventEnd(PETSCSF_EmbedSF,sf,0,0,0);CHKERRQ(ierr);
124895fce210SBarry Smith   PetscFunctionReturn(0);
124995fce210SBarry Smith }
125095fce210SBarry Smith 
12512f5fb4c2SMatthew G. Knepley /*@C
12522f5fb4c2SMatthew G. Knepley   PetscSFCreateEmbeddedLeafSF - removes edges from all but the selected leaves, does not remap indices
12532f5fb4c2SMatthew G. Knepley 
12542f5fb4c2SMatthew G. Knepley   Collective
12552f5fb4c2SMatthew G. Knepley 
12562f5fb4c2SMatthew G. Knepley   Input Arguments:
12572f5fb4c2SMatthew G. Knepley + sf - original star forest
1258f659e5c7SJunchao Zhang . nselected  - number of selected leaves on this process
1259f659e5c7SJunchao Zhang - selected   - indices of the selected leaves on this process
12602f5fb4c2SMatthew G. Knepley 
12612f5fb4c2SMatthew G. Knepley   Output Arguments:
12622f5fb4c2SMatthew G. Knepley .  newsf - new star forest
12632f5fb4c2SMatthew G. Knepley 
12642f5fb4c2SMatthew G. Knepley   Level: advanced
12652f5fb4c2SMatthew G. Knepley 
12662f5fb4c2SMatthew G. Knepley .seealso: PetscSFCreateEmbeddedSF(), PetscSFSetGraph(), PetscSFGetGraph()
12672f5fb4c2SMatthew G. Knepley @*/
1268f659e5c7SJunchao Zhang PetscErrorCode PetscSFCreateEmbeddedLeafSF(PetscSF sf,PetscInt nselected,const PetscInt *selected,PetscSF *newsf)
12692f5fb4c2SMatthew G. Knepley {
1270f659e5c7SJunchao Zhang   const PetscSFNode *iremote;
1271f659e5c7SJunchao Zhang   PetscSFNode       *new_iremote;
1272f659e5c7SJunchao Zhang   const PetscInt    *ilocal;
1273f659e5c7SJunchao Zhang   PetscInt          i,nroots,*leaves,*new_ilocal;
1274f659e5c7SJunchao Zhang   MPI_Comm          comm;
12752f5fb4c2SMatthew G. Knepley   PetscErrorCode    ierr;
12762f5fb4c2SMatthew G. Knepley 
12772f5fb4c2SMatthew G. Knepley   PetscFunctionBegin;
12782f5fb4c2SMatthew G. Knepley   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
127929046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
1280f659e5c7SJunchao Zhang   if (nselected) PetscValidPointer(selected,3);
12812f5fb4c2SMatthew G. Knepley   PetscValidPointer(newsf,4);
12822f5fb4c2SMatthew G. Knepley 
1283f659e5c7SJunchao Zhang   /* Uniq selected[] and put results in leaves[] */
1284f659e5c7SJunchao Zhang   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
1285f659e5c7SJunchao Zhang   ierr = PetscMalloc1(nselected,&leaves);CHKERRQ(ierr);
1286dd5b3ca6SJunchao Zhang   ierr = PetscArraycpy(leaves,selected,nselected);CHKERRQ(ierr);
1287f659e5c7SJunchao Zhang   ierr = PetscSortedRemoveDupsInt(&nselected,leaves);CHKERRQ(ierr);
1288f659e5c7SJunchao 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);
1289f659e5c7SJunchao Zhang 
1290f659e5c7SJunchao Zhang   /* Optimize the routine only when sf is setup and hence we can reuse sf's communication pattern */
1291f659e5c7SJunchao Zhang   if (sf->setupcalled && sf->ops->CreateEmbeddedLeafSF) {
1292f659e5c7SJunchao Zhang     ierr = (*sf->ops->CreateEmbeddedLeafSF)(sf,nselected,leaves,newsf);CHKERRQ(ierr);
1293f659e5c7SJunchao Zhang   } else {
1294f659e5c7SJunchao Zhang     ierr = PetscSFGetGraph(sf,&nroots,NULL,&ilocal,&iremote);CHKERRQ(ierr);
1295f659e5c7SJunchao Zhang     ierr = PetscMalloc1(nselected,&new_ilocal);CHKERRQ(ierr);
1296f659e5c7SJunchao Zhang     ierr = PetscMalloc1(nselected,&new_iremote);CHKERRQ(ierr);
1297f659e5c7SJunchao Zhang     for (i=0; i<nselected; ++i) {
1298f659e5c7SJunchao Zhang       const PetscInt l     = leaves[i];
1299f659e5c7SJunchao Zhang       new_ilocal[i]        = ilocal ? ilocal[l] : l;
1300f659e5c7SJunchao Zhang       new_iremote[i].rank  = iremote[l].rank;
1301f659e5c7SJunchao Zhang       new_iremote[i].index = iremote[l].index;
13022f5fb4c2SMatthew G. Knepley     }
13034cc19a0cSStefano Zampini     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,newsf);CHKERRQ(ierr);
1304f659e5c7SJunchao Zhang     ierr = PetscSFSetGraph(*newsf,nroots,nselected,new_ilocal,PETSC_OWN_POINTER,new_iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1305f659e5c7SJunchao Zhang   }
1306f659e5c7SJunchao Zhang   ierr = PetscFree(leaves);CHKERRQ(ierr);
13072f5fb4c2SMatthew G. Knepley   PetscFunctionReturn(0);
13082f5fb4c2SMatthew G. Knepley }
13092f5fb4c2SMatthew G. Knepley 
131095fce210SBarry Smith /*@C
13113482bfa8SJunchao Zhang    PetscSFBcastAndOpBegin - begin pointwise broadcast with root value being reduced to leaf value, to be concluded with call to PetscSFBcastAndOpEnd()
13123482bfa8SJunchao Zhang 
13133482bfa8SJunchao Zhang    Collective on PetscSF
13143482bfa8SJunchao Zhang 
13153482bfa8SJunchao Zhang    Input Arguments:
13163482bfa8SJunchao Zhang +  sf - star forest on which to communicate
13173482bfa8SJunchao Zhang .  unit - data type associated with each node
13183482bfa8SJunchao Zhang .  rootdata - buffer to broadcast
13193482bfa8SJunchao Zhang -  op - operation to use for reduction
13203482bfa8SJunchao Zhang 
13213482bfa8SJunchao Zhang    Output Arguments:
13223482bfa8SJunchao Zhang .  leafdata - buffer to be reduced with values from each leaf's respective root
13233482bfa8SJunchao Zhang 
13243482bfa8SJunchao Zhang    Level: intermediate
13253482bfa8SJunchao Zhang 
13263482bfa8SJunchao Zhang .seealso: PetscSFBcastAndOpEnd(), PetscSFBcastBegin(), PetscSFBcastEnd()
13273482bfa8SJunchao Zhang @*/
13283482bfa8SJunchao Zhang PetscErrorCode PetscSFBcastAndOpBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
13293482bfa8SJunchao Zhang {
13303482bfa8SJunchao Zhang   PetscErrorCode ierr;
1331eb02082bSJunchao Zhang   PetscMemType   rootmtype,leafmtype;
13323482bfa8SJunchao Zhang 
13333482bfa8SJunchao Zhang   PetscFunctionBegin;
13343482bfa8SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
13353482bfa8SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
13363482bfa8SJunchao Zhang   ierr = PetscLogEventBegin(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
1337eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
1338eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
1339eb02082bSJunchao Zhang   ierr = (*sf->ops->BcastAndOpBegin)(sf,unit,rootmtype,rootdata,leafmtype,leafdata,op);CHKERRQ(ierr);
13403482bfa8SJunchao Zhang   ierr = PetscLogEventEnd(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
13413482bfa8SJunchao Zhang   PetscFunctionReturn(0);
13423482bfa8SJunchao Zhang }
13433482bfa8SJunchao Zhang 
13443482bfa8SJunchao Zhang /*@C
13453482bfa8SJunchao Zhang    PetscSFBcastAndOpEnd - end a broadcast & reduce operation started with PetscSFBcastAndOpBegin()
13463482bfa8SJunchao Zhang 
13473482bfa8SJunchao Zhang    Collective
13483482bfa8SJunchao Zhang 
13493482bfa8SJunchao Zhang    Input Arguments:
13503482bfa8SJunchao Zhang +  sf - star forest
13513482bfa8SJunchao Zhang .  unit - data type
13523482bfa8SJunchao Zhang .  rootdata - buffer to broadcast
13533482bfa8SJunchao Zhang -  op - operation to use for reduction
13543482bfa8SJunchao Zhang 
13553482bfa8SJunchao Zhang    Output Arguments:
13563482bfa8SJunchao Zhang .  leafdata - buffer to be reduced with values from each leaf's respective root
13573482bfa8SJunchao Zhang 
13583482bfa8SJunchao Zhang    Level: intermediate
13593482bfa8SJunchao Zhang 
13603482bfa8SJunchao Zhang .seealso: PetscSFSetGraph(), PetscSFReduceEnd()
13613482bfa8SJunchao Zhang @*/
13623482bfa8SJunchao Zhang PetscErrorCode PetscSFBcastAndOpEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
13633482bfa8SJunchao Zhang {
13643482bfa8SJunchao Zhang   PetscErrorCode ierr;
13653482bfa8SJunchao Zhang 
13663482bfa8SJunchao Zhang   PetscFunctionBegin;
13673482bfa8SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
13683482bfa8SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
13693482bfa8SJunchao Zhang   ierr = PetscLogEventBegin(PETSCSF_BcastAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
137000816365SJunchao Zhang   ierr = (*sf->ops->BcastAndOpEnd)(sf,unit,rootdata,leafdata,op);CHKERRQ(ierr);
13713482bfa8SJunchao Zhang   ierr = PetscLogEventEnd(PETSCSF_BcastAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
13723482bfa8SJunchao Zhang   PetscFunctionReturn(0);
13733482bfa8SJunchao Zhang }
13743482bfa8SJunchao Zhang 
13753482bfa8SJunchao Zhang /*@C
137695fce210SBarry Smith    PetscSFReduceBegin - begin reduction of leafdata into rootdata, to be completed with call to PetscSFReduceEnd()
137795fce210SBarry Smith 
137895fce210SBarry Smith    Collective
137995fce210SBarry Smith 
138095fce210SBarry Smith    Input Arguments:
138195fce210SBarry Smith +  sf - star forest
138295fce210SBarry Smith .  unit - data type
138395fce210SBarry Smith .  leafdata - values to reduce
138495fce210SBarry Smith -  op - reduction operation
138595fce210SBarry Smith 
138695fce210SBarry Smith    Output Arguments:
138795fce210SBarry Smith .  rootdata - result of reduction of values from all leaves of each root
138895fce210SBarry Smith 
138995fce210SBarry Smith    Level: intermediate
139095fce210SBarry Smith 
139195fce210SBarry Smith .seealso: PetscSFBcastBegin()
139295fce210SBarry Smith @*/
139395fce210SBarry Smith PetscErrorCode PetscSFReduceBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
139495fce210SBarry Smith {
139595fce210SBarry Smith   PetscErrorCode ierr;
1396eb02082bSJunchao Zhang   PetscMemType   rootmtype,leafmtype;
139795fce210SBarry Smith 
139895fce210SBarry Smith   PetscFunctionBegin;
139995fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
140095fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
140129046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);
1402eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
1403eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
1404eb02082bSJunchao Zhang   ierr = (sf->ops->ReduceBegin)(sf,unit,leafmtype,leafdata,rootmtype,rootdata,op);CHKERRQ(ierr);
1405563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);
140695fce210SBarry Smith   PetscFunctionReturn(0);
140795fce210SBarry Smith }
140895fce210SBarry Smith 
140995fce210SBarry Smith /*@C
141095fce210SBarry Smith    PetscSFReduceEnd - end a reduction operation started with PetscSFReduceBegin()
141195fce210SBarry Smith 
141295fce210SBarry Smith    Collective
141395fce210SBarry Smith 
141495fce210SBarry Smith    Input Arguments:
141595fce210SBarry Smith +  sf - star forest
141695fce210SBarry Smith .  unit - data type
141795fce210SBarry Smith .  leafdata - values to reduce
141895fce210SBarry Smith -  op - reduction operation
141995fce210SBarry Smith 
142095fce210SBarry Smith    Output Arguments:
142195fce210SBarry Smith .  rootdata - result of reduction of values from all leaves of each root
142295fce210SBarry Smith 
142395fce210SBarry Smith    Level: intermediate
142495fce210SBarry Smith 
142595fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFBcastEnd()
142695fce210SBarry Smith @*/
142795fce210SBarry Smith PetscErrorCode PetscSFReduceEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
142895fce210SBarry Smith {
142995fce210SBarry Smith   PetscErrorCode ierr;
143095fce210SBarry Smith 
143195fce210SBarry Smith   PetscFunctionBegin;
143295fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
143395fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
143429046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);
143500816365SJunchao Zhang   ierr = (*sf->ops->ReduceEnd)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
1436563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);
143795fce210SBarry Smith   PetscFunctionReturn(0);
143895fce210SBarry Smith }
143995fce210SBarry Smith 
144095fce210SBarry Smith /*@C
1441a1729e3fSJunchao Zhang    PetscSFFetchAndOpBegin - begin operation that fetches values from root and updates atomically by applying operation using my leaf value, to be completed with PetscSFFetchAndOpEnd()
1442a1729e3fSJunchao Zhang 
1443a1729e3fSJunchao Zhang    Collective
1444a1729e3fSJunchao Zhang 
1445a1729e3fSJunchao Zhang    Input Arguments:
1446a1729e3fSJunchao Zhang +  sf - star forest
1447a1729e3fSJunchao Zhang .  unit - data type
1448a1729e3fSJunchao Zhang .  leafdata - leaf values to use in reduction
1449a1729e3fSJunchao Zhang -  op - operation to use for reduction
1450a1729e3fSJunchao Zhang 
1451a1729e3fSJunchao Zhang    Output Arguments:
1452a1729e3fSJunchao Zhang +  rootdata - root values to be updated, input state is seen by first process to perform an update
1453a1729e3fSJunchao Zhang -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
1454a1729e3fSJunchao Zhang 
1455a1729e3fSJunchao Zhang    Level: advanced
1456a1729e3fSJunchao Zhang 
1457a1729e3fSJunchao Zhang    Note:
1458a1729e3fSJunchao Zhang    The update is only atomic at the granularity provided by the hardware. Different roots referenced by the same process
1459a1729e3fSJunchao Zhang    might be updated in a different order. Furthermore, if a composite type is used for the unit datatype, atomicity is
1460a1729e3fSJunchao Zhang    not guaranteed across the whole vertex. Therefore, this function is mostly only used with primitive types such as
1461a1729e3fSJunchao Zhang    integers.
1462a1729e3fSJunchao Zhang 
1463a1729e3fSJunchao Zhang .seealso: PetscSFComputeDegreeBegin(), PetscSFReduceBegin(), PetscSFSetGraph()
1464a1729e3fSJunchao Zhang @*/
1465a1729e3fSJunchao Zhang PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1466a1729e3fSJunchao Zhang {
1467a1729e3fSJunchao Zhang   PetscErrorCode ierr;
1468eb02082bSJunchao Zhang   PetscMemType   rootmtype,leafmtype,leafupdatemtype;
1469a1729e3fSJunchao Zhang 
1470a1729e3fSJunchao Zhang   PetscFunctionBegin;
1471a1729e3fSJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1472a1729e3fSJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1473a1729e3fSJunchao Zhang   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
1474eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
1475eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
1476eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafupdate,&leafupdatemtype);CHKERRQ(ierr);
1477eb02082bSJunchao Zhang   if (leafmtype != leafupdatemtype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for leafdata and leafupdate in different memory types");
1478eb02082bSJunchao Zhang   ierr = (*sf->ops->FetchAndOpBegin)(sf,unit,rootmtype,rootdata,leafmtype,leafdata,leafupdate,op);CHKERRQ(ierr);
1479a1729e3fSJunchao Zhang   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
1480a1729e3fSJunchao Zhang   PetscFunctionReturn(0);
1481a1729e3fSJunchao Zhang }
1482a1729e3fSJunchao Zhang 
1483a1729e3fSJunchao Zhang /*@C
1484a1729e3fSJunchao 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
1485a1729e3fSJunchao Zhang 
1486a1729e3fSJunchao Zhang    Collective
1487a1729e3fSJunchao Zhang 
1488a1729e3fSJunchao Zhang    Input Arguments:
1489a1729e3fSJunchao Zhang +  sf - star forest
1490a1729e3fSJunchao Zhang .  unit - data type
1491a1729e3fSJunchao Zhang .  leafdata - leaf values to use in reduction
1492a1729e3fSJunchao Zhang -  op - operation to use for reduction
1493a1729e3fSJunchao Zhang 
1494a1729e3fSJunchao Zhang    Output Arguments:
1495a1729e3fSJunchao Zhang +  rootdata - root values to be updated, input state is seen by first process to perform an update
1496a1729e3fSJunchao Zhang -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
1497a1729e3fSJunchao Zhang 
1498a1729e3fSJunchao Zhang    Level: advanced
1499a1729e3fSJunchao Zhang 
1500a1729e3fSJunchao Zhang .seealso: PetscSFComputeDegreeEnd(), PetscSFReduceEnd(), PetscSFSetGraph()
1501a1729e3fSJunchao Zhang @*/
1502a1729e3fSJunchao Zhang PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1503a1729e3fSJunchao Zhang {
1504a1729e3fSJunchao Zhang   PetscErrorCode ierr;
1505a1729e3fSJunchao Zhang 
1506a1729e3fSJunchao Zhang   PetscFunctionBegin;
1507a1729e3fSJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1508a1729e3fSJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1509a1729e3fSJunchao Zhang   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
151000816365SJunchao Zhang   ierr = (*sf->ops->FetchAndOpEnd)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr);
1511a1729e3fSJunchao Zhang   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
1512a1729e3fSJunchao Zhang   PetscFunctionReturn(0);
1513a1729e3fSJunchao Zhang }
1514a1729e3fSJunchao Zhang 
1515a1729e3fSJunchao Zhang /*@C
151695fce210SBarry Smith    PetscSFComputeDegreeBegin - begin computation of degree for each root vertex, to be completed with PetscSFComputeDegreeEnd()
151795fce210SBarry Smith 
151895fce210SBarry Smith    Collective
151995fce210SBarry Smith 
152095fce210SBarry Smith    Input Arguments:
152195fce210SBarry Smith .  sf - star forest
152295fce210SBarry Smith 
152395fce210SBarry Smith    Output Arguments:
152495fce210SBarry Smith .  degree - degree of each root vertex
152595fce210SBarry Smith 
152695fce210SBarry Smith    Level: advanced
152795fce210SBarry Smith 
1528ffe67aa5SVáclav Hapla    Notes:
1529ffe67aa5SVáclav Hapla    The returned array is owned by PetscSF and automatically freed by PetscSFDestroy(). Hence no need to call PetscFree() on it.
1530ffe67aa5SVáclav Hapla 
153195fce210SBarry Smith .seealso: PetscSFGatherBegin()
153295fce210SBarry Smith @*/
153395fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf,const PetscInt **degree)
153495fce210SBarry Smith {
153595fce210SBarry Smith   PetscErrorCode ierr;
153695fce210SBarry Smith 
153795fce210SBarry Smith   PetscFunctionBegin;
153895fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
153995fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
154095fce210SBarry Smith   PetscValidPointer(degree,2);
1541803bd9e8SMatthew G. Knepley   if (!sf->degreeknown) {
15425b0d146aSStefano Zampini     PetscInt i, nroots = sf->nroots, maxlocal;
1543803bd9e8SMatthew G. Knepley     if (sf->degree) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Calls to PetscSFComputeDegreeBegin() cannot be nested.");
15445b0d146aSStefano Zampini     maxlocal = sf->maxleaf-sf->minleaf+1;
154529046d53SLisandro Dalcin     ierr = PetscMalloc1(nroots,&sf->degree);CHKERRQ(ierr);
154629046d53SLisandro Dalcin     ierr = PetscMalloc1(PetscMax(maxlocal,1),&sf->degreetmp);CHKERRQ(ierr); /* allocate at least one entry, see check in PetscSFComputeDegreeEnd() */
154729046d53SLisandro Dalcin     for (i=0; i<nroots; i++) sf->degree[i] = 0;
15489837ea96SMatthew G. Knepley     for (i=0; i<maxlocal; i++) sf->degreetmp[i] = 1;
15495b0d146aSStefano Zampini     ierr = PetscSFReduceBegin(sf,MPIU_INT,sf->degreetmp-sf->minleaf,sf->degree,MPI_SUM);CHKERRQ(ierr);
155095fce210SBarry Smith   }
155195fce210SBarry Smith   *degree = NULL;
155295fce210SBarry Smith   PetscFunctionReturn(0);
155395fce210SBarry Smith }
155495fce210SBarry Smith 
155595fce210SBarry Smith /*@C
155695fce210SBarry Smith    PetscSFComputeDegreeEnd - complete computation of degree for each root vertex, started with PetscSFComputeDegreeBegin()
155795fce210SBarry Smith 
155895fce210SBarry Smith    Collective
155995fce210SBarry Smith 
156095fce210SBarry Smith    Input Arguments:
156195fce210SBarry Smith .  sf - star forest
156295fce210SBarry Smith 
156395fce210SBarry Smith    Output Arguments:
156495fce210SBarry Smith .  degree - degree of each root vertex
156595fce210SBarry Smith 
156695fce210SBarry Smith    Level: developer
156795fce210SBarry Smith 
1568ffe67aa5SVáclav Hapla    Notes:
1569ffe67aa5SVáclav Hapla    The returned array is owned by PetscSF and automatically freed by PetscSFDestroy(). Hence no need to call PetscFree() on it.
1570ffe67aa5SVáclav Hapla 
157195fce210SBarry Smith .seealso:
157295fce210SBarry Smith @*/
157395fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf,const PetscInt **degree)
157495fce210SBarry Smith {
157595fce210SBarry Smith   PetscErrorCode ierr;
157695fce210SBarry Smith 
157795fce210SBarry Smith   PetscFunctionBegin;
157895fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
157995fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
158029046d53SLisandro Dalcin   PetscValidPointer(degree,2);
158195fce210SBarry Smith   if (!sf->degreeknown) {
158229046d53SLisandro Dalcin     if (!sf->degreetmp) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFComputeDegreeBegin() before PetscSFComputeDegreeEnd()");
15835b0d146aSStefano Zampini     ierr = PetscSFReduceEnd(sf,MPIU_INT,sf->degreetmp-sf->minleaf,sf->degree,MPI_SUM);CHKERRQ(ierr);
158495fce210SBarry Smith     ierr = PetscFree(sf->degreetmp);CHKERRQ(ierr);
158595fce210SBarry Smith     sf->degreeknown = PETSC_TRUE;
158695fce210SBarry Smith   }
158795fce210SBarry Smith   *degree = sf->degree;
158895fce210SBarry Smith   PetscFunctionReturn(0);
158995fce210SBarry Smith }
159095fce210SBarry Smith 
1591673100f5SVaclav Hapla 
1592673100f5SVaclav Hapla /*@C
159366dfcd1aSVaclav Hapla    PetscSFComputeMultiRootOriginalNumbering - Returns original numbering of multi-roots (roots of multi-SF returned by PetscSFGetMultiSF()).
159466dfcd1aSVaclav Hapla    Each multi-root is assigned index of the corresponding original root.
1595673100f5SVaclav Hapla 
1596673100f5SVaclav Hapla    Collective
1597673100f5SVaclav Hapla 
1598673100f5SVaclav Hapla    Input Arguments:
1599673100f5SVaclav Hapla +  sf - star forest
1600673100f5SVaclav Hapla -  degree - degree of each root vertex, computed with PetscSFComputeDegreeBegin()/PetscSFComputeDegreeEnd()
1601673100f5SVaclav Hapla 
1602673100f5SVaclav Hapla    Output Arguments:
160366dfcd1aSVaclav Hapla +  nMultiRoots - (optional) number of multi-roots (roots of multi-SF)
160466dfcd1aSVaclav Hapla -  multiRootsOrigNumbering - original indices of multi-roots; length of this array is nMultiRoots
1605673100f5SVaclav Hapla 
1606673100f5SVaclav Hapla    Level: developer
1607673100f5SVaclav Hapla 
1608ffe67aa5SVáclav Hapla    Notes:
1609ffe67aa5SVáclav Hapla    The returned array multiRootsOrigNumbering is newly allocated and should be destroyed with PetscFree() when no longer needed.
1610ffe67aa5SVáclav Hapla 
1611673100f5SVaclav Hapla .seealso: PetscSFComputeDegreeBegin(), PetscSFComputeDegreeEnd(), PetscSFGetMultiSF()
1612673100f5SVaclav Hapla @*/
161366dfcd1aSVaclav Hapla PetscErrorCode PetscSFComputeMultiRootOriginalNumbering(PetscSF sf, const PetscInt degree[], PetscInt *nMultiRoots, PetscInt *multiRootsOrigNumbering[])
1614673100f5SVaclav Hapla {
1615673100f5SVaclav Hapla   PetscSF             msf;
1616673100f5SVaclav Hapla   PetscInt            i, j, k, nroots, nmroots;
1617673100f5SVaclav Hapla   PetscErrorCode      ierr;
1618673100f5SVaclav Hapla 
1619673100f5SVaclav Hapla   PetscFunctionBegin;
1620673100f5SVaclav Hapla   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1621673100f5SVaclav Hapla   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
162229328920SVaclav Hapla   if (nroots) PetscValidIntPointer(degree,2);
162366dfcd1aSVaclav Hapla   if (nMultiRoots) PetscValidIntPointer(nMultiRoots,3);
162466dfcd1aSVaclav Hapla   PetscValidPointer(multiRootsOrigNumbering,4);
1625673100f5SVaclav Hapla   ierr = PetscSFGetMultiSF(sf,&msf);CHKERRQ(ierr);
1626673100f5SVaclav Hapla   ierr = PetscSFGetGraph(msf, &nmroots, NULL, NULL, NULL);CHKERRQ(ierr);
162766dfcd1aSVaclav Hapla   ierr = PetscMalloc1(nmroots, multiRootsOrigNumbering);CHKERRQ(ierr);
1628673100f5SVaclav Hapla   for (i=0,j=0,k=0; i<nroots; i++) {
1629673100f5SVaclav Hapla     if (!degree[i]) continue;
1630673100f5SVaclav Hapla     for (j=0; j<degree[i]; j++,k++) {
163166dfcd1aSVaclav Hapla       (*multiRootsOrigNumbering)[k] = i;
1632673100f5SVaclav Hapla     }
1633673100f5SVaclav Hapla   }
1634673100f5SVaclav Hapla #if defined(PETSC_USE_DEBUG)
1635673100f5SVaclav Hapla   if (PetscUnlikely(k != nmroots)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sanity check fail");
1636673100f5SVaclav Hapla #endif
163766dfcd1aSVaclav Hapla   if (nMultiRoots) *nMultiRoots = nmroots;
1638673100f5SVaclav Hapla   PetscFunctionReturn(0);
1639673100f5SVaclav Hapla }
1640673100f5SVaclav Hapla 
164195fce210SBarry Smith /*@C
164295fce210SBarry Smith    PetscSFGatherBegin - begin pointwise gather of all leaves into multi-roots, to be completed with PetscSFGatherEnd()
164395fce210SBarry Smith 
164495fce210SBarry Smith    Collective
164595fce210SBarry Smith 
164695fce210SBarry Smith    Input Arguments:
164795fce210SBarry Smith +  sf - star forest
164895fce210SBarry Smith .  unit - data type
164995fce210SBarry Smith -  leafdata - leaf data to gather to roots
165095fce210SBarry Smith 
165195fce210SBarry Smith    Output Argument:
165295fce210SBarry Smith .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
165395fce210SBarry Smith 
165495fce210SBarry Smith    Level: intermediate
165595fce210SBarry Smith 
165695fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
165795fce210SBarry Smith @*/
165895fce210SBarry Smith PetscErrorCode PetscSFGatherBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
165995fce210SBarry Smith {
166095fce210SBarry Smith   PetscErrorCode ierr;
166195fce210SBarry Smith   PetscSF        multi;
166295fce210SBarry Smith 
166395fce210SBarry Smith   PetscFunctionBegin;
166495fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
166529046d53SLisandro Dalcin   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
166695fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
16678bfbc91cSJed Brown   ierr = PetscSFReduceBegin(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr);
166895fce210SBarry Smith   PetscFunctionReturn(0);
166995fce210SBarry Smith }
167095fce210SBarry Smith 
167195fce210SBarry Smith /*@C
167295fce210SBarry Smith    PetscSFGatherEnd - ends pointwise gather operation that was started with PetscSFGatherBegin()
167395fce210SBarry Smith 
167495fce210SBarry Smith    Collective
167595fce210SBarry Smith 
167695fce210SBarry Smith    Input Arguments:
167795fce210SBarry Smith +  sf - star forest
167895fce210SBarry Smith .  unit - data type
167995fce210SBarry Smith -  leafdata - leaf data to gather to roots
168095fce210SBarry Smith 
168195fce210SBarry Smith    Output Argument:
168295fce210SBarry Smith .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
168395fce210SBarry Smith 
168495fce210SBarry Smith    Level: intermediate
168595fce210SBarry Smith 
168695fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
168795fce210SBarry Smith @*/
168895fce210SBarry Smith PetscErrorCode PetscSFGatherEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
168995fce210SBarry Smith {
169095fce210SBarry Smith   PetscErrorCode ierr;
169195fce210SBarry Smith   PetscSF        multi;
169295fce210SBarry Smith 
169395fce210SBarry Smith   PetscFunctionBegin;
169495fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
169595fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
169695fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
16978bfbc91cSJed Brown   ierr = PetscSFReduceEnd(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr);
169895fce210SBarry Smith   PetscFunctionReturn(0);
169995fce210SBarry Smith }
170095fce210SBarry Smith 
170195fce210SBarry Smith /*@C
170295fce210SBarry Smith    PetscSFScatterBegin - begin pointwise scatter operation from multi-roots to leaves, to be completed with PetscSFScatterEnd()
170395fce210SBarry Smith 
170495fce210SBarry Smith    Collective
170595fce210SBarry Smith 
170695fce210SBarry Smith    Input Arguments:
170795fce210SBarry Smith +  sf - star forest
170895fce210SBarry Smith .  unit - data type
170995fce210SBarry Smith -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
171095fce210SBarry Smith 
171195fce210SBarry Smith    Output Argument:
171295fce210SBarry Smith .  leafdata - leaf data to be update with personal data from each respective root
171395fce210SBarry Smith 
171495fce210SBarry Smith    Level: intermediate
171595fce210SBarry Smith 
171695fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
171795fce210SBarry Smith @*/
171895fce210SBarry Smith PetscErrorCode PetscSFScatterBegin(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
171995fce210SBarry Smith {
172095fce210SBarry Smith   PetscErrorCode ierr;
172195fce210SBarry Smith   PetscSF        multi;
172295fce210SBarry Smith 
172395fce210SBarry Smith   PetscFunctionBegin;
172495fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
172595fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
172695fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
172795fce210SBarry Smith   ierr = PetscSFBcastBegin(multi,unit,multirootdata,leafdata);CHKERRQ(ierr);
172895fce210SBarry Smith   PetscFunctionReturn(0);
172995fce210SBarry Smith }
173095fce210SBarry Smith 
173195fce210SBarry Smith /*@C
173295fce210SBarry Smith    PetscSFScatterEnd - ends pointwise scatter operation that was started with PetscSFScatterBegin()
173395fce210SBarry Smith 
173495fce210SBarry Smith    Collective
173595fce210SBarry Smith 
173695fce210SBarry Smith    Input Arguments:
173795fce210SBarry Smith +  sf - star forest
173895fce210SBarry Smith .  unit - data type
173995fce210SBarry Smith -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
174095fce210SBarry Smith 
174195fce210SBarry Smith    Output Argument:
174295fce210SBarry Smith .  leafdata - leaf data to be update with personal data from each respective root
174395fce210SBarry Smith 
174495fce210SBarry Smith    Level: intermediate
174595fce210SBarry Smith 
174695fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
174795fce210SBarry Smith @*/
174895fce210SBarry Smith PetscErrorCode PetscSFScatterEnd(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
174995fce210SBarry Smith {
175095fce210SBarry Smith   PetscErrorCode ierr;
175195fce210SBarry Smith   PetscSF        multi;
175295fce210SBarry Smith 
175395fce210SBarry Smith   PetscFunctionBegin;
175495fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
175595fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
175695fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
175795fce210SBarry Smith   ierr = PetscSFBcastEnd(multi,unit,multirootdata,leafdata);CHKERRQ(ierr);
175895fce210SBarry Smith   PetscFunctionReturn(0);
175995fce210SBarry Smith }
1760a7b3aa13SAta Mesgarnejad 
1761a072220fSLawrence Mitchell static PetscErrorCode PetscSFCheckLeavesUnique_Private(PetscSF sf)
1762a072220fSLawrence Mitchell {
1763a072220fSLawrence Mitchell #if defined(PETSC_USE_DEBUG)
1764a072220fSLawrence Mitchell   PetscInt        i, n, nleaves;
1765a072220fSLawrence Mitchell   const PetscInt *ilocal = NULL;
1766a072220fSLawrence Mitchell   PetscHSetI      seen;
1767a072220fSLawrence Mitchell   PetscErrorCode  ierr;
1768a072220fSLawrence Mitchell 
1769a072220fSLawrence Mitchell   PetscFunctionBegin;
1770a072220fSLawrence Mitchell   ierr = PetscSFGetGraph(sf,NULL,&nleaves,&ilocal,NULL);CHKERRQ(ierr);
1771a072220fSLawrence Mitchell   ierr = PetscHSetICreate(&seen);CHKERRQ(ierr);
1772a072220fSLawrence Mitchell   for (i = 0; i < nleaves; i++) {
1773a072220fSLawrence Mitchell     const PetscInt leaf = ilocal ? ilocal[i] : i;
1774a072220fSLawrence Mitchell     ierr = PetscHSetIAdd(seen,leaf);CHKERRQ(ierr);
1775a072220fSLawrence Mitchell   }
1776a072220fSLawrence Mitchell   ierr = PetscHSetIGetSize(seen,&n);CHKERRQ(ierr);
1777a072220fSLawrence Mitchell   if (n != nleaves) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Provided leaves have repeated values: all leaves must be unique");
1778a072220fSLawrence Mitchell   ierr = PetscHSetIDestroy(&seen);CHKERRQ(ierr);
1779a072220fSLawrence Mitchell   PetscFunctionReturn(0);
1780a072220fSLawrence Mitchell #else
1781a072220fSLawrence Mitchell   PetscFunctionBegin;
1782a072220fSLawrence Mitchell   PetscFunctionReturn(0);
1783a072220fSLawrence Mitchell #endif
1784a072220fSLawrence Mitchell }
178554729392SStefano Zampini 
1786a7b3aa13SAta Mesgarnejad /*@
178704c0ada0SJunchao Zhang   PetscSFCompose - Compose a new PetscSF by putting the second SF under the first one in a top (roots) down (leaves) view
1788a7b3aa13SAta Mesgarnejad 
1789a7b3aa13SAta Mesgarnejad   Input Parameters:
179054729392SStefano Zampini + sfA - The first PetscSF
179154729392SStefano Zampini - sfB - The second PetscSF
1792a7b3aa13SAta Mesgarnejad 
1793a7b3aa13SAta Mesgarnejad   Output Parameters:
179454729392SStefano Zampini . sfBA - The composite SF
1795a7b3aa13SAta Mesgarnejad 
1796a7b3aa13SAta Mesgarnejad   Level: developer
1797a7b3aa13SAta Mesgarnejad 
1798a072220fSLawrence Mitchell   Notes:
179954729392SStefano Zampini   Currently, the two SFs must be defined on congruent communicators and they must be true star
180054729392SStefano Zampini   forests, i.e. the same leaf is not connected with different roots.
180154729392SStefano Zampini 
180254729392SStefano Zampini   sfA's leaf space and sfB's root space might be partially overlapped. The composition builds
180354729392SStefano Zampini   a graph with sfA's roots and sfB's leaves only when there is a path between them. Unconnected
180454729392SStefano Zampini   nodes (roots or leaves) are not in sfBA. Doing a Bcast on the new SF is equivalent to doing a
180554729392SStefano Zampini   Bcast on sfA, then a Bcast on sfB, on connected nodes.
1806a072220fSLawrence Mitchell 
180704c0ada0SJunchao Zhang .seealso: PetscSF, PetscSFComposeInverse(), PetscSFGetGraph(), PetscSFSetGraph()
1808a7b3aa13SAta Mesgarnejad @*/
1809a7b3aa13SAta Mesgarnejad PetscErrorCode PetscSFCompose(PetscSF sfA,PetscSF sfB,PetscSF *sfBA)
1810a7b3aa13SAta Mesgarnejad {
181104c0ada0SJunchao Zhang   PetscErrorCode    ierr;
1812a7b3aa13SAta Mesgarnejad   const PetscSFNode *remotePointsA,*remotePointsB;
1813d41018fbSJunchao Zhang   PetscSFNode       *remotePointsBA=NULL,*reorderedRemotePointsA = NULL,*leafdataB;
181454729392SStefano Zampini   const PetscInt    *localPointsA,*localPointsB;
181554729392SStefano Zampini   PetscInt          *localPointsBA;
181654729392SStefano Zampini   PetscInt          i,numRootsA,numLeavesA,numRootsB,numLeavesB,minleaf,maxleaf,numLeavesBA;
181754729392SStefano Zampini   PetscBool         denseB;
1818a7b3aa13SAta Mesgarnejad 
1819a7b3aa13SAta Mesgarnejad   PetscFunctionBegin;
1820a7b3aa13SAta Mesgarnejad   PetscValidHeaderSpecific(sfA,PETSCSF_CLASSID,1);
182129046d53SLisandro Dalcin   PetscSFCheckGraphSet(sfA,1);
182229046d53SLisandro Dalcin   PetscValidHeaderSpecific(sfB,PETSCSF_CLASSID,2);
182329046d53SLisandro Dalcin   PetscSFCheckGraphSet(sfB,2);
182454729392SStefano Zampini   PetscCheckSameComm(sfA,1,sfB,2);
182529046d53SLisandro Dalcin   PetscValidPointer(sfBA,3);
182654729392SStefano Zampini   ierr = PetscSFCheckLeavesUnique_Private(sfA);CHKERRQ(ierr);
182754729392SStefano Zampini   ierr = PetscSFCheckLeavesUnique_Private(sfB);CHKERRQ(ierr);
182854729392SStefano Zampini 
1829a7b3aa13SAta Mesgarnejad   ierr = PetscSFGetGraph(sfA,&numRootsA,&numLeavesA,&localPointsA,&remotePointsA);CHKERRQ(ierr);
1830a7b3aa13SAta Mesgarnejad   ierr = PetscSFGetGraph(sfB,&numRootsB,&numLeavesB,&localPointsB,&remotePointsB);CHKERRQ(ierr);
1831d41018fbSJunchao Zhang   if (localPointsA) {
183254729392SStefano Zampini     ierr = PetscMalloc1(numRootsB,&reorderedRemotePointsA);CHKERRQ(ierr);
183354729392SStefano Zampini     for (i=0; i<numRootsB; i++) {
183454729392SStefano Zampini       reorderedRemotePointsA[i].rank = -1;
183554729392SStefano Zampini       reorderedRemotePointsA[i].index = -1;
183654729392SStefano Zampini     }
183754729392SStefano Zampini     for (i=0; i<numLeavesA; i++) {
183854729392SStefano Zampini       if (localPointsA[i] >= numRootsB) continue;
183954729392SStefano Zampini       reorderedRemotePointsA[localPointsA[i]] = remotePointsA[i];
184054729392SStefano Zampini     }
1841d41018fbSJunchao Zhang     remotePointsA = reorderedRemotePointsA;
1842d41018fbSJunchao Zhang   }
1843d41018fbSJunchao Zhang   ierr = PetscSFGetLeafRange(sfB,&minleaf,&maxleaf);CHKERRQ(ierr);
1844d41018fbSJunchao Zhang   ierr = PetscMalloc1(maxleaf-minleaf+1,&leafdataB);CHKERRQ(ierr);
1845d41018fbSJunchao Zhang   ierr = PetscSFBcastBegin(sfB,MPIU_2INT,remotePointsA,leafdataB-minleaf);CHKERRQ(ierr);
1846d41018fbSJunchao Zhang   ierr = PetscSFBcastEnd(sfB,MPIU_2INT,remotePointsA,leafdataB-minleaf);CHKERRQ(ierr);
1847d41018fbSJunchao Zhang   ierr = PetscFree(reorderedRemotePointsA);CHKERRQ(ierr);
1848d41018fbSJunchao Zhang 
184954729392SStefano Zampini   denseB = (PetscBool)!localPointsB;
185054729392SStefano Zampini   for (i=0,numLeavesBA=0; i<numLeavesB; i++) {
185154729392SStefano Zampini     if (leafdataB[localPointsB ? localPointsB[i]-minleaf : i].rank == -1) denseB = PETSC_FALSE;
185254729392SStefano Zampini     else numLeavesBA++;
185354729392SStefano Zampini   }
185454729392SStefano Zampini   if (denseB) {
1855d41018fbSJunchao Zhang     localPointsBA  = NULL;
1856d41018fbSJunchao Zhang     remotePointsBA = leafdataB;
1857d41018fbSJunchao Zhang   } else {
185854729392SStefano Zampini     ierr = PetscMalloc1(numLeavesBA,&localPointsBA);CHKERRQ(ierr);
185954729392SStefano Zampini     ierr = PetscMalloc1(numLeavesBA,&remotePointsBA);CHKERRQ(ierr);
186054729392SStefano Zampini     for (i=0,numLeavesBA=0; i<numLeavesB; i++) {
186154729392SStefano Zampini       const PetscInt l = localPointsB ? localPointsB[i] : i;
186254729392SStefano Zampini 
186354729392SStefano Zampini       if (leafdataB[l-minleaf].rank == -1) continue;
186454729392SStefano Zampini       remotePointsBA[numLeavesBA] = leafdataB[l-minleaf];
186554729392SStefano Zampini       localPointsBA[numLeavesBA] = l;
186654729392SStefano Zampini       numLeavesBA++;
186754729392SStefano Zampini     }
1868d41018fbSJunchao Zhang     ierr = PetscFree(leafdataB);CHKERRQ(ierr);
1869d41018fbSJunchao Zhang   }
187054729392SStefano Zampini   ierr = PetscSFCreate(PetscObjectComm((PetscObject)sfA),sfBA);CHKERRQ(ierr);
187154729392SStefano Zampini   ierr = PetscSFSetGraph(*sfBA,numRootsA,numLeavesBA,localPointsBA,PETSC_OWN_POINTER,remotePointsBA,PETSC_OWN_POINTER);CHKERRQ(ierr);
1872a7b3aa13SAta Mesgarnejad   PetscFunctionReturn(0);
1873a7b3aa13SAta Mesgarnejad }
18741c6ba672SJunchao Zhang 
187504c0ada0SJunchao Zhang /*@
187654729392SStefano Zampini   PetscSFComposeInverse - Compose a new PetscSF by putting the inverse of the second SF under the first one
187704c0ada0SJunchao Zhang 
187804c0ada0SJunchao Zhang   Input Parameters:
187954729392SStefano Zampini + sfA - The first PetscSF
188054729392SStefano Zampini - sfB - The second PetscSF
188104c0ada0SJunchao Zhang 
188204c0ada0SJunchao Zhang   Output Parameters:
188354729392SStefano Zampini . sfBA - The composite SF.
188404c0ada0SJunchao Zhang 
188504c0ada0SJunchao Zhang   Level: developer
188604c0ada0SJunchao Zhang 
188754729392SStefano Zampini   Notes:
188854729392SStefano Zampini   Currently, the two SFs must be defined on congruent communicators and they must be true star
188954729392SStefano Zampini   forests, i.e. the same leaf is not connected with different roots. Even more, all roots of the
189054729392SStefano Zampini   second SF must have a degree of 1, i.e., no roots have more than one leaf connected.
189154729392SStefano Zampini 
189254729392SStefano Zampini   sfA's leaf space and sfB's leaf space might be partially overlapped. The composition builds
189354729392SStefano Zampini   a graph with sfA's roots and sfB's roots only when there is a path between them. Unconnected
189454729392SStefano Zampini   roots are not in sfBA. Doing a Bcast on the new SF is equivalent to doing a Bcast on sfA, then
189554729392SStefano Zampini   a Reduce on sfB, on connected roots.
189654729392SStefano Zampini 
189754729392SStefano Zampini .seealso: PetscSF, PetscSFCompose(), PetscSFGetGraph(), PetscSFSetGraph(), PetscSFCreateInverseSF()
189804c0ada0SJunchao Zhang @*/
189904c0ada0SJunchao Zhang PetscErrorCode PetscSFComposeInverse(PetscSF sfA,PetscSF sfB,PetscSF *sfBA)
190004c0ada0SJunchao Zhang {
190104c0ada0SJunchao Zhang   PetscErrorCode    ierr;
190204c0ada0SJunchao Zhang   const PetscSFNode *remotePointsA,*remotePointsB;
190304c0ada0SJunchao Zhang   PetscSFNode       *remotePointsBA;
190404c0ada0SJunchao Zhang   const PetscInt    *localPointsA,*localPointsB;
190554729392SStefano Zampini   PetscSFNode       *reorderedRemotePointsA = NULL;
190654729392SStefano Zampini   PetscInt          i,numRootsA,numLeavesA,numLeavesBA,numRootsB,numLeavesB,minleaf,maxleaf,*localPointsBA;
19075b0d146aSStefano Zampini   MPI_Op            op;
19085b0d146aSStefano Zampini #if defined(PETSC_USE_64BIT_INDICES)
19095b0d146aSStefano Zampini   PetscBool         iswin;
19105b0d146aSStefano Zampini #endif
191104c0ada0SJunchao Zhang 
191204c0ada0SJunchao Zhang   PetscFunctionBegin;
191304c0ada0SJunchao Zhang   PetscValidHeaderSpecific(sfA,PETSCSF_CLASSID,1);
191404c0ada0SJunchao Zhang   PetscSFCheckGraphSet(sfA,1);
191504c0ada0SJunchao Zhang   PetscValidHeaderSpecific(sfB,PETSCSF_CLASSID,2);
191604c0ada0SJunchao Zhang   PetscSFCheckGraphSet(sfB,2);
191754729392SStefano Zampini   PetscCheckSameComm(sfA,1,sfB,2);
191804c0ada0SJunchao Zhang   PetscValidPointer(sfBA,3);
191954729392SStefano Zampini   ierr = PetscSFCheckLeavesUnique_Private(sfA);CHKERRQ(ierr);
192054729392SStefano Zampini   ierr = PetscSFCheckLeavesUnique_Private(sfB);CHKERRQ(ierr);
192154729392SStefano Zampini 
192204c0ada0SJunchao Zhang   ierr = PetscSFGetGraph(sfA, &numRootsA, &numLeavesA, &localPointsA, &remotePointsA);CHKERRQ(ierr);
192304c0ada0SJunchao Zhang   ierr = PetscSFGetGraph(sfB, &numRootsB, &numLeavesB, &localPointsB, &remotePointsB);CHKERRQ(ierr);
19245b0d146aSStefano Zampini 
19255b0d146aSStefano Zampini   /* TODO: Check roots of sfB have degree of 1 */
19265b0d146aSStefano Zampini   /* Once we implement it, we can replace the MPI_MAXLOC
19275b0d146aSStefano Zampini      with MPIU_REPLACE. In that case, MPI_MAXLOC and MPIU_REPLACE have the same effect.
19285b0d146aSStefano Zampini      We use MPI_MAXLOC only to have a deterministic output from this routine if
19295b0d146aSStefano Zampini      the root condition is not meet.
19305b0d146aSStefano Zampini    */
19315b0d146aSStefano Zampini   op = MPI_MAXLOC;
19325b0d146aSStefano Zampini #if defined(PETSC_USE_64BIT_INDICES)
19335b0d146aSStefano Zampini   /* we accept a non-deterministic output (if any) with PETSCSFWINDOW, since MPI_MAXLOC cannot operate on MPIU_2INT with MPI_Accumulate */
19345b0d146aSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)sfB,PETSCSFWINDOW,&iswin);CHKERRQ(ierr);
19355b0d146aSStefano Zampini   if (iswin) op = MPIU_REPLACE;
19365b0d146aSStefano Zampini #endif
19375b0d146aSStefano Zampini 
193854729392SStefano Zampini   ierr = PetscSFGetLeafRange(sfB, &minleaf, &maxleaf);CHKERRQ(ierr);
193954729392SStefano Zampini   ierr = PetscMalloc1(maxleaf - minleaf + 1,&reorderedRemotePointsA);CHKERRQ(ierr);
194054729392SStefano Zampini   for (i=0; i<maxleaf - minleaf + 1; i++) {
194154729392SStefano Zampini     reorderedRemotePointsA[i].rank = -1;
194254729392SStefano Zampini     reorderedRemotePointsA[i].index = -1;
194354729392SStefano Zampini   }
194454729392SStefano Zampini   if (localPointsA) {
194554729392SStefano Zampini     for (i=0; i<numLeavesA; i++) {
194654729392SStefano Zampini       if (localPointsA[i] > maxleaf || localPointsA[i] < minleaf) continue;
194754729392SStefano Zampini       reorderedRemotePointsA[localPointsA[i] - minleaf] = remotePointsA[i];
194854729392SStefano Zampini     }
194954729392SStefano Zampini   } else {
195054729392SStefano Zampini     for (i=0; i<numLeavesA; i++) {
195154729392SStefano Zampini       if (i > maxleaf || i < minleaf) continue;
195254729392SStefano Zampini       reorderedRemotePointsA[i - minleaf] = remotePointsA[i];
195354729392SStefano Zampini     }
195454729392SStefano Zampini   }
195554729392SStefano Zampini 
195654729392SStefano Zampini   ierr = PetscMalloc1(numRootsB,&localPointsBA);CHKERRQ(ierr);
195704c0ada0SJunchao Zhang   ierr = PetscMalloc1(numRootsB,&remotePointsBA);CHKERRQ(ierr);
195854729392SStefano Zampini   for (i=0; i<numRootsB; i++) {
195954729392SStefano Zampini     remotePointsBA[i].rank = -1;
196054729392SStefano Zampini     remotePointsBA[i].index = -1;
196154729392SStefano Zampini   }
196254729392SStefano Zampini 
19635b0d146aSStefano Zampini   ierr = PetscSFReduceBegin(sfB,MPIU_2INT,reorderedRemotePointsA-minleaf,remotePointsBA,op);CHKERRQ(ierr);
19645b0d146aSStefano Zampini   ierr = PetscSFReduceEnd(sfB,MPIU_2INT,reorderedRemotePointsA-minleaf,remotePointsBA,op);CHKERRQ(ierr);
196554729392SStefano Zampini   ierr = PetscFree(reorderedRemotePointsA);CHKERRQ(ierr);
196654729392SStefano Zampini   for (i=0,numLeavesBA=0; i<numRootsB; i++) {
196754729392SStefano Zampini     if (remotePointsBA[i].rank == -1) continue;
196854729392SStefano Zampini     remotePointsBA[numLeavesBA].rank = remotePointsBA[i].rank;
196954729392SStefano Zampini     remotePointsBA[numLeavesBA].index = remotePointsBA[i].index;
197054729392SStefano Zampini     localPointsBA[numLeavesBA] = i;
197154729392SStefano Zampini     numLeavesBA++;
197254729392SStefano Zampini   }
197354729392SStefano Zampini   ierr = PetscSFCreate(PetscObjectComm((PetscObject)sfA),sfBA);CHKERRQ(ierr);
197454729392SStefano Zampini   ierr = PetscSFSetGraph(*sfBA,numRootsA,numLeavesBA,localPointsBA,PETSC_OWN_POINTER,remotePointsBA,PETSC_OWN_POINTER);CHKERRQ(ierr);
197504c0ada0SJunchao Zhang   PetscFunctionReturn(0);
197604c0ada0SJunchao Zhang }
197704c0ada0SJunchao Zhang 
19781c6ba672SJunchao Zhang /*
19791c6ba672SJunchao Zhang   PetscSFCreateLocalSF_Private - Creates a local PetscSF that only has intra-process edges of the global PetscSF
19801c6ba672SJunchao Zhang 
19811c6ba672SJunchao Zhang   Input Parameters:
19821c6ba672SJunchao Zhang . sf - The global PetscSF
19831c6ba672SJunchao Zhang 
19841c6ba672SJunchao Zhang   Output Parameters:
19851c6ba672SJunchao Zhang . out - The local PetscSF
19861c6ba672SJunchao Zhang  */
19871c6ba672SJunchao Zhang PetscErrorCode PetscSFCreateLocalSF_Private(PetscSF sf,PetscSF *out)
19881c6ba672SJunchao Zhang {
19891c6ba672SJunchao Zhang   MPI_Comm           comm;
19901c6ba672SJunchao Zhang   PetscMPIInt        myrank;
19911c6ba672SJunchao Zhang   const PetscInt     *ilocal;
19921c6ba672SJunchao Zhang   const PetscSFNode  *iremote;
19931c6ba672SJunchao Zhang   PetscInt           i,j,nroots,nleaves,lnleaves,*lilocal;
19941c6ba672SJunchao Zhang   PetscSFNode        *liremote;
19951c6ba672SJunchao Zhang   PetscSF            lsf;
19961c6ba672SJunchao Zhang   PetscErrorCode     ierr;
19971c6ba672SJunchao Zhang 
19981c6ba672SJunchao Zhang   PetscFunctionBegin;
19991c6ba672SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
20001c6ba672SJunchao Zhang   if (sf->ops->CreateLocalSF) {
20011c6ba672SJunchao Zhang     ierr = (*sf->ops->CreateLocalSF)(sf,out);CHKERRQ(ierr);
20021c6ba672SJunchao Zhang   } else {
20031c6ba672SJunchao Zhang     /* Could use PetscSFCreateEmbeddedLeafSF, but since we know the comm is PETSC_COMM_SELF, we can make it fast */
20041c6ba672SJunchao Zhang     ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
20051c6ba672SJunchao Zhang     ierr = MPI_Comm_rank(comm,&myrank);CHKERRQ(ierr);
20061c6ba672SJunchao Zhang 
20071c6ba672SJunchao Zhang     /* Find out local edges and build a local SF */
20081c6ba672SJunchao Zhang     ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
20091c6ba672SJunchao Zhang     for (i=lnleaves=0; i<nleaves; i++) {if (iremote[i].rank == (PetscInt)myrank) lnleaves++;}
20101c6ba672SJunchao Zhang     ierr = PetscMalloc1(lnleaves,&lilocal);CHKERRQ(ierr);
20111c6ba672SJunchao Zhang     ierr = PetscMalloc1(lnleaves,&liremote);CHKERRQ(ierr);
20121c6ba672SJunchao Zhang 
20131c6ba672SJunchao Zhang     for (i=j=0; i<nleaves; i++) {
20141c6ba672SJunchao Zhang       if (iremote[i].rank == (PetscInt)myrank) {
20151c6ba672SJunchao Zhang         lilocal[j]        = ilocal? ilocal[i] : i; /* ilocal=NULL for contiguous storage */
20161c6ba672SJunchao Zhang         liremote[j].rank  = 0; /* rank in PETSC_COMM_SELF */
20171c6ba672SJunchao Zhang         liremote[j].index = iremote[i].index;
20181c6ba672SJunchao Zhang         j++;
20191c6ba672SJunchao Zhang       }
20201c6ba672SJunchao Zhang     }
20211c6ba672SJunchao Zhang     ierr = PetscSFCreate(PETSC_COMM_SELF,&lsf);CHKERRQ(ierr);
20221c6ba672SJunchao Zhang     ierr = PetscSFSetGraph(lsf,nroots,lnleaves,lilocal,PETSC_OWN_POINTER,liremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
20231c6ba672SJunchao Zhang     ierr = PetscSFSetUp(lsf);CHKERRQ(ierr);
20241c6ba672SJunchao Zhang     *out = lsf;
20251c6ba672SJunchao Zhang   }
20261c6ba672SJunchao Zhang   PetscFunctionReturn(0);
20271c6ba672SJunchao Zhang }
2028dd5b3ca6SJunchao Zhang 
2029dd5b3ca6SJunchao Zhang /* Similar to PetscSFBcast, but only Bcast to leaves on rank 0 */
2030dd5b3ca6SJunchao Zhang PetscErrorCode PetscSFBcastToZero_Private(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
2031dd5b3ca6SJunchao Zhang {
2032dd5b3ca6SJunchao Zhang   PetscErrorCode     ierr;
2033eb02082bSJunchao Zhang   PetscMemType       rootmtype,leafmtype;
2034dd5b3ca6SJunchao Zhang 
2035dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
2036dd5b3ca6SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
2037dd5b3ca6SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
2038dd5b3ca6SJunchao Zhang   ierr = PetscLogEventBegin(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
2039eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
2040eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
2041dd5b3ca6SJunchao Zhang   if (sf->ops->BcastToZero) {
2042eb02082bSJunchao Zhang     ierr = (*sf->ops->BcastToZero)(sf,unit,rootmtype,rootdata,leafmtype,leafdata);CHKERRQ(ierr);
2043dd5b3ca6SJunchao Zhang   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSFBcastToZero_Private is not supported on this SF type");
2044dd5b3ca6SJunchao Zhang   ierr = PetscLogEventEnd(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
2045dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
2046dd5b3ca6SJunchao Zhang }
2047dd5b3ca6SJunchao Zhang 
2048