xref: /petsc/src/vec/is/sf/interface/sf.c (revision 1c6ba67218fe9020c97a3bfccb426c200ebde0e5)
1af0996ceSBarry Smith #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/
295fce210SBarry Smith #include <petscctable.h>
395fce210SBarry Smith 
495fce210SBarry Smith #if defined(PETSC_USE_DEBUG)
595fce210SBarry Smith #  define PetscSFCheckGraphSet(sf,arg) do {                          \
695fce210SBarry Smith     if (PetscUnlikely(!(sf)->graphset))                              \
795fce210SBarry Smith       SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetGraph() on argument %D \"%s\" before %s()",(arg),#sf,PETSC_FUNCTION_NAME); \
895fce210SBarry Smith   } while (0)
995fce210SBarry Smith #else
1095fce210SBarry Smith #  define PetscSFCheckGraphSet(sf,arg) do {} while (0)
1195fce210SBarry Smith #endif
1295fce210SBarry Smith 
1395fce210SBarry Smith const char *const PetscSFDuplicateOptions[] = {"CONFONLY","RANKS","GRAPH","PetscSFDuplicateOption","PETSCSF_DUPLICATE_",0};
1495fce210SBarry Smith 
158af6ec1cSBarry Smith /*@
1695fce210SBarry Smith    PetscSFCreate - create a star forest communication context
1795fce210SBarry Smith 
18d083f849SBarry Smith    Collective
1995fce210SBarry Smith 
2095fce210SBarry Smith    Input Arguments:
2195fce210SBarry Smith .  comm - communicator on which the star forest will operate
2295fce210SBarry Smith 
2395fce210SBarry Smith    Output Arguments:
2495fce210SBarry Smith .  sf - new star forest context
2595fce210SBarry Smith 
2695fce210SBarry Smith    Level: intermediate
2795fce210SBarry Smith 
2895fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFDestroy()
2995fce210SBarry Smith @*/
3095fce210SBarry Smith PetscErrorCode PetscSFCreate(MPI_Comm comm,PetscSF *sf)
3195fce210SBarry Smith {
3295fce210SBarry Smith   PetscErrorCode ierr;
3395fce210SBarry Smith   PetscSF        b;
3495fce210SBarry Smith 
3595fce210SBarry Smith   PetscFunctionBegin;
3695fce210SBarry Smith   PetscValidPointer(sf,2);
37607a6623SBarry Smith   ierr = PetscSFInitializePackage();CHKERRQ(ierr);
3895fce210SBarry Smith 
3973107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(b,PETSCSF_CLASSID,"PetscSF","Star Forest","PetscSF",comm,PetscSFDestroy,PetscSFView);CHKERRQ(ierr);
4095fce210SBarry Smith 
4195fce210SBarry Smith   b->nroots    = -1;
4295fce210SBarry Smith   b->nleaves   = -1;
4329046d53SLisandro Dalcin   b->minleaf   = PETSC_MAX_INT;
4429046d53SLisandro Dalcin   b->maxleaf   = PETSC_MIN_INT;
4595fce210SBarry Smith   b->nranks    = -1;
4695fce210SBarry Smith   b->rankorder = PETSC_TRUE;
4795fce210SBarry Smith   b->ingroup   = MPI_GROUP_NULL;
4895fce210SBarry Smith   b->outgroup  = MPI_GROUP_NULL;
4995fce210SBarry Smith   b->graphset  = PETSC_FALSE;
5095fce210SBarry Smith 
5195fce210SBarry Smith   *sf = b;
5295fce210SBarry Smith   PetscFunctionReturn(0);
5395fce210SBarry Smith }
5495fce210SBarry Smith 
5529046d53SLisandro Dalcin /*@
5695fce210SBarry Smith    PetscSFReset - Reset a star forest so that different sizes or neighbors can be used
5795fce210SBarry Smith 
5895fce210SBarry Smith    Collective
5995fce210SBarry Smith 
6095fce210SBarry Smith    Input Arguments:
6195fce210SBarry Smith .  sf - star forest
6295fce210SBarry Smith 
6395fce210SBarry Smith    Level: advanced
6495fce210SBarry Smith 
6595fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFDestroy()
6695fce210SBarry Smith @*/
6795fce210SBarry Smith PetscErrorCode PetscSFReset(PetscSF sf)
6895fce210SBarry Smith {
6995fce210SBarry Smith   PetscErrorCode ierr;
7095fce210SBarry Smith 
7195fce210SBarry Smith   PetscFunctionBegin;
7295fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
7379715d56SJed Brown   if (sf->ops->Reset) {ierr = (*sf->ops->Reset)(sf);CHKERRQ(ierr);}
7429046d53SLisandro Dalcin   sf->nroots   = -1;
7529046d53SLisandro Dalcin   sf->nleaves  = -1;
7629046d53SLisandro Dalcin   sf->minleaf  = PETSC_MAX_INT;
7729046d53SLisandro Dalcin   sf->maxleaf  = PETSC_MIN_INT;
7895fce210SBarry Smith   sf->mine     = NULL;
7995fce210SBarry Smith   sf->remote   = NULL;
8029046d53SLisandro Dalcin   sf->graphset = PETSC_FALSE;
8129046d53SLisandro Dalcin   ierr = PetscFree(sf->mine_alloc);CHKERRQ(ierr);
8295fce210SBarry Smith   ierr = PetscFree(sf->remote_alloc);CHKERRQ(ierr);
8321c688dcSJed Brown   sf->nranks = -1;
8429046d53SLisandro Dalcin   ierr = PetscFree4(sf->ranks,sf->roffset,sf->rmine,sf->rremote);CHKERRQ(ierr);
8529046d53SLisandro Dalcin   sf->degreeknown = PETSC_FALSE;
8695fce210SBarry Smith   ierr = PetscFree(sf->degree);CHKERRQ(ierr);
8795fce210SBarry Smith   if (sf->ingroup  != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->ingroup);CHKERRQ(ierr);}
8895fce210SBarry Smith   if (sf->outgroup != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->outgroup);CHKERRQ(ierr);}
8995fce210SBarry Smith   ierr = PetscSFDestroy(&sf->multi);CHKERRQ(ierr);
9095fce210SBarry Smith   sf->setupcalled = PETSC_FALSE;
9195fce210SBarry Smith   PetscFunctionReturn(0);
9295fce210SBarry Smith }
9395fce210SBarry Smith 
9495fce210SBarry Smith /*@C
9529046d53SLisandro Dalcin    PetscSFSetType - Set the PetscSF communication implementation
9695fce210SBarry Smith 
9795fce210SBarry Smith    Collective on PetscSF
9895fce210SBarry Smith 
9995fce210SBarry Smith    Input Parameters:
10095fce210SBarry Smith +  sf - the PetscSF context
10195fce210SBarry Smith -  type - a known method
10295fce210SBarry Smith 
10395fce210SBarry Smith    Options Database Key:
10495fce210SBarry Smith .  -sf_type <type> - Sets the method; use -help for a list
10595fce210SBarry Smith    of available methods (for instance, window, pt2pt, neighbor)
10695fce210SBarry Smith 
10795fce210SBarry Smith    Notes:
10895fce210SBarry Smith    See "include/petscsf.h" for available methods (for instance)
10995fce210SBarry Smith +    PETSCSFWINDOW - MPI-2/3 one-sided
11095fce210SBarry Smith -    PETSCSFBASIC - basic implementation using MPI-1 two-sided
11195fce210SBarry Smith 
11295fce210SBarry Smith   Level: intermediate
11395fce210SBarry Smith 
11495fce210SBarry Smith .seealso: PetscSFType, PetscSFCreate()
11595fce210SBarry Smith @*/
11695fce210SBarry Smith PetscErrorCode PetscSFSetType(PetscSF sf,PetscSFType type)
11795fce210SBarry Smith {
11895fce210SBarry Smith   PetscErrorCode ierr,(*r)(PetscSF);
11995fce210SBarry Smith   PetscBool      match;
12095fce210SBarry Smith 
12195fce210SBarry Smith   PetscFunctionBegin;
12295fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
12395fce210SBarry Smith   PetscValidCharPointer(type,2);
12495fce210SBarry Smith 
12595fce210SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)sf,type,&match);CHKERRQ(ierr);
12695fce210SBarry Smith   if (match) PetscFunctionReturn(0);
12795fce210SBarry Smith 
128adc40e5bSBarry Smith   ierr = PetscFunctionListFind(PetscSFList,type,&r);CHKERRQ(ierr);
12995fce210SBarry Smith   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested PetscSF type %s",type);
13029046d53SLisandro Dalcin   /* Destroy the previous PetscSF implementation context */
13129046d53SLisandro Dalcin   if (sf->ops->Destroy) {ierr = (*(sf)->ops->Destroy)(sf);CHKERRQ(ierr);}
13295fce210SBarry Smith   ierr = PetscMemzero(sf->ops,sizeof(*sf->ops));CHKERRQ(ierr);
13395fce210SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)sf,type);CHKERRQ(ierr);
13495fce210SBarry Smith   ierr = (*r)(sf);CHKERRQ(ierr);
13595fce210SBarry Smith   PetscFunctionReturn(0);
13695fce210SBarry Smith }
13795fce210SBarry Smith 
13829046d53SLisandro Dalcin /*@C
13929046d53SLisandro Dalcin   PetscSFGetType - Get the PetscSF communication implementation
14029046d53SLisandro Dalcin 
14129046d53SLisandro Dalcin   Not Collective
14229046d53SLisandro Dalcin 
14329046d53SLisandro Dalcin   Input Parameter:
14429046d53SLisandro Dalcin . sf  - the PetscSF context
14529046d53SLisandro Dalcin 
14629046d53SLisandro Dalcin   Output Parameter:
14729046d53SLisandro Dalcin . type - the PetscSF type name
14829046d53SLisandro Dalcin 
14929046d53SLisandro Dalcin   Level: intermediate
15029046d53SLisandro Dalcin 
15129046d53SLisandro Dalcin .seealso: PetscSFSetType(), PetscSFCreate()
15229046d53SLisandro Dalcin @*/
15329046d53SLisandro Dalcin PetscErrorCode PetscSFGetType(PetscSF sf, PetscSFType *type)
15429046d53SLisandro Dalcin {
15529046d53SLisandro Dalcin   PetscFunctionBegin;
15629046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID,1);
15729046d53SLisandro Dalcin   PetscValidPointer(type,2);
15829046d53SLisandro Dalcin   *type = ((PetscObject)sf)->type_name;
15929046d53SLisandro Dalcin   PetscFunctionReturn(0);
16029046d53SLisandro Dalcin }
16129046d53SLisandro Dalcin 
162d36d33e4SMatthew G. Knepley /*@
16395fce210SBarry Smith    PetscSFDestroy - destroy star forest
16495fce210SBarry Smith 
16595fce210SBarry Smith    Collective
16695fce210SBarry Smith 
16795fce210SBarry Smith    Input Arguments:
16895fce210SBarry Smith .  sf - address of star forest
16995fce210SBarry Smith 
17095fce210SBarry Smith    Level: intermediate
17195fce210SBarry Smith 
17295fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFReset()
17395fce210SBarry Smith @*/
17495fce210SBarry Smith PetscErrorCode PetscSFDestroy(PetscSF *sf)
17595fce210SBarry Smith {
17695fce210SBarry Smith   PetscErrorCode ierr;
17795fce210SBarry Smith 
17895fce210SBarry Smith   PetscFunctionBegin;
17995fce210SBarry Smith   if (!*sf) PetscFunctionReturn(0);
18095fce210SBarry Smith   PetscValidHeaderSpecific((*sf),PETSCSF_CLASSID,1);
18129046d53SLisandro Dalcin   if (--((PetscObject)(*sf))->refct > 0) {*sf = NULL; PetscFunctionReturn(0);}
18295fce210SBarry Smith   ierr = PetscSFReset(*sf);CHKERRQ(ierr);
18395fce210SBarry Smith   if ((*sf)->ops->Destroy) {ierr = (*(*sf)->ops->Destroy)(*sf);CHKERRQ(ierr);}
18495fce210SBarry Smith   ierr = PetscHeaderDestroy(sf);CHKERRQ(ierr);
18595fce210SBarry Smith   PetscFunctionReturn(0);
18695fce210SBarry Smith }
18795fce210SBarry Smith 
18895fce210SBarry Smith /*@
18995fce210SBarry Smith    PetscSFSetUp - set up communication structures
19095fce210SBarry Smith 
19195fce210SBarry Smith    Collective
19295fce210SBarry Smith 
19395fce210SBarry Smith    Input Arguments:
19495fce210SBarry Smith .  sf - star forest communication object
19595fce210SBarry Smith 
19695fce210SBarry Smith    Level: beginner
19795fce210SBarry Smith 
19895fce210SBarry Smith .seealso: PetscSFSetFromOptions(), PetscSFSetType()
19995fce210SBarry Smith @*/
20095fce210SBarry Smith PetscErrorCode PetscSFSetUp(PetscSF sf)
20195fce210SBarry Smith {
20295fce210SBarry Smith   PetscErrorCode ierr;
20395fce210SBarry Smith 
20495fce210SBarry Smith   PetscFunctionBegin;
20529046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
20629046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
20795fce210SBarry Smith   if (sf->setupcalled) PetscFunctionReturn(0);
20895fce210SBarry Smith   if (!((PetscObject)sf)->type_name) {ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr);}
20929046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_SetUp,sf,0,0,0);CHKERRQ(ierr);
21095fce210SBarry Smith   if (sf->ops->SetUp) {ierr = (*sf->ops->SetUp)(sf);CHKERRQ(ierr);}
21129046d53SLisandro Dalcin   ierr = PetscLogEventEnd(PETSCSF_SetUp,sf,0,0,0);CHKERRQ(ierr);
21295fce210SBarry Smith   sf->setupcalled = PETSC_TRUE;
21395fce210SBarry Smith   PetscFunctionReturn(0);
21495fce210SBarry Smith }
21595fce210SBarry Smith 
2168af6ec1cSBarry Smith /*@
21795fce210SBarry Smith    PetscSFSetFromOptions - set PetscSF options using the options database
21895fce210SBarry Smith 
21995fce210SBarry Smith    Logically Collective
22095fce210SBarry Smith 
22195fce210SBarry Smith    Input Arguments:
22295fce210SBarry Smith .  sf - star forest
22395fce210SBarry Smith 
22495fce210SBarry Smith    Options Database Keys:
22560263706SJed Brown +  -sf_type - implementation type, see PetscSFSetType()
22660263706SJed Brown -  -sf_rank_order - sort composite points for gathers and scatters in rank order, gathers are non-deterministic otherwise
22795fce210SBarry Smith 
22895fce210SBarry Smith    Level: intermediate
22995fce210SBarry Smith 
23095fce210SBarry Smith .seealso: PetscSFWindowSetSyncType()
23195fce210SBarry Smith @*/
23295fce210SBarry Smith PetscErrorCode PetscSFSetFromOptions(PetscSF sf)
23395fce210SBarry Smith {
23495fce210SBarry Smith   PetscSFType    deft;
23595fce210SBarry Smith   char           type[256];
23695fce210SBarry Smith   PetscErrorCode ierr;
23795fce210SBarry Smith   PetscBool      flg;
23895fce210SBarry Smith 
23995fce210SBarry Smith   PetscFunctionBegin;
24095fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
24195fce210SBarry Smith   ierr = PetscObjectOptionsBegin((PetscObject)sf);CHKERRQ(ierr);
24295fce210SBarry Smith   deft = ((PetscObject)sf)->type_name ? ((PetscObject)sf)->type_name : PETSCSFBASIC;
24329046d53SLisandro Dalcin   ierr = PetscOptionsFList("-sf_type","PetscSF implementation type","PetscSFSetType",PetscSFList,deft,type,sizeof(type),&flg);CHKERRQ(ierr);
24495fce210SBarry Smith   ierr = PetscSFSetType(sf,flg ? type : deft);CHKERRQ(ierr);
24595fce210SBarry 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);
246e55864a3SBarry Smith   if (sf->ops->SetFromOptions) {ierr = (*sf->ops->SetFromOptions)(PetscOptionsObject,sf);CHKERRQ(ierr);}
24795fce210SBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
24895fce210SBarry Smith   PetscFunctionReturn(0);
24995fce210SBarry Smith }
25095fce210SBarry Smith 
25129046d53SLisandro Dalcin /*@
25295fce210SBarry Smith    PetscSFSetRankOrder - sort multi-points for gathers and scatters by rank order
25395fce210SBarry Smith 
25495fce210SBarry Smith    Logically Collective
25595fce210SBarry Smith 
25695fce210SBarry Smith    Input Arguments:
25795fce210SBarry Smith +  sf - star forest
25895fce210SBarry Smith -  flg - PETSC_TRUE to sort, PETSC_FALSE to skip sorting (lower setup cost, but non-deterministic)
25995fce210SBarry Smith 
26095fce210SBarry Smith    Level: advanced
26195fce210SBarry Smith 
26295fce210SBarry Smith .seealso: PetscSFGatherBegin(), PetscSFScatterBegin()
26395fce210SBarry Smith @*/
26495fce210SBarry Smith PetscErrorCode PetscSFSetRankOrder(PetscSF sf,PetscBool flg)
26595fce210SBarry Smith {
26695fce210SBarry Smith 
26795fce210SBarry Smith   PetscFunctionBegin;
26895fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
26995fce210SBarry Smith   PetscValidLogicalCollectiveBool(sf,flg,2);
27095fce210SBarry Smith   if (sf->multi) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Rank ordering must be set before first call to PetscSFGatherBegin() or PetscSFScatterBegin()");
27195fce210SBarry Smith   sf->rankorder = flg;
27295fce210SBarry Smith   PetscFunctionReturn(0);
27395fce210SBarry Smith }
27495fce210SBarry Smith 
2758af6ec1cSBarry Smith /*@
27695fce210SBarry Smith    PetscSFSetGraph - Set a parallel star forest
27795fce210SBarry Smith 
27895fce210SBarry Smith    Collective
27995fce210SBarry Smith 
28095fce210SBarry Smith    Input Arguments:
28195fce210SBarry Smith +  sf - star forest
28295fce210SBarry Smith .  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
28395fce210SBarry Smith .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
28495fce210SBarry Smith .  ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage
28595fce210SBarry Smith .  localmode - copy mode for ilocal
28695fce210SBarry Smith .  iremote - remote locations of root vertices for each leaf on the current process
28795fce210SBarry Smith -  remotemode - copy mode for iremote
28895fce210SBarry Smith 
28995fce210SBarry Smith    Level: intermediate
29095fce210SBarry Smith 
29195452b02SPatrick Sanan    Notes:
29295452b02SPatrick Sanan     In Fortran you must use PETSC_COPY_VALUES for localmode and remotemode
29338ab3f8aSBarry Smith 
2942ad1e87fSLisandro Dalcin    Developers Note: Local indices which are the identity permutation in the range [0,nleaves) are discarded as they
2952ad1e87fSLisandro Dalcin    encode contiguous storage. In such case, if localmode is PETSC_OWN_POINTER, the memory is deallocated as it is not
2962ad1e87fSLisandro Dalcin    needed
2972ad1e87fSLisandro Dalcin 
29895fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFGetGraph()
29995fce210SBarry Smith @*/
30095fce210SBarry Smith PetscErrorCode PetscSFSetGraph(PetscSF sf,PetscInt nroots,PetscInt nleaves,const PetscInt *ilocal,PetscCopyMode localmode,const PetscSFNode *iremote,PetscCopyMode remotemode)
30195fce210SBarry Smith {
30295fce210SBarry Smith   PetscErrorCode ierr;
30395fce210SBarry Smith 
30495fce210SBarry Smith   PetscFunctionBegin;
30595fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
30629046d53SLisandro Dalcin   if (nleaves > 0 && ilocal) PetscValidIntPointer(ilocal,4);
30729046d53SLisandro Dalcin   if (nleaves > 0) PetscValidPointer(iremote,6);
30829046d53SLisandro Dalcin   if (nroots  < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nroots %D, cannot be negative",nroots);
30995fce210SBarry Smith   if (nleaves < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nleaves %D, cannot be negative",nleaves);
31029046d53SLisandro Dalcin 
31195fce210SBarry Smith   ierr = PetscSFReset(sf);CHKERRQ(ierr);
31229046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr);
31329046d53SLisandro Dalcin 
31495fce210SBarry Smith   sf->nroots  = nroots;
31595fce210SBarry Smith   sf->nleaves = nleaves;
31629046d53SLisandro Dalcin 
31729046d53SLisandro Dalcin   if (nleaves && ilocal) {
31821c688dcSJed Brown     PetscInt i;
31929046d53SLisandro Dalcin     PetscInt minleaf = PETSC_MAX_INT;
32029046d53SLisandro Dalcin     PetscInt maxleaf = PETSC_MIN_INT;
3212ad1e87fSLisandro Dalcin     int      contiguous = 1;
32229046d53SLisandro Dalcin     for (i=0; i<nleaves; i++) {
32329046d53SLisandro Dalcin       minleaf = PetscMin(minleaf,ilocal[i]);
32429046d53SLisandro Dalcin       maxleaf = PetscMax(maxleaf,ilocal[i]);
3252ad1e87fSLisandro Dalcin       contiguous &= (ilocal[i] == i);
32629046d53SLisandro Dalcin     }
32729046d53SLisandro Dalcin     sf->minleaf = minleaf;
32829046d53SLisandro Dalcin     sf->maxleaf = maxleaf;
3292ad1e87fSLisandro Dalcin     if (contiguous) {
3302ad1e87fSLisandro Dalcin       if (localmode == PETSC_OWN_POINTER) {
3312ad1e87fSLisandro Dalcin         ierr = PetscFree(ilocal);CHKERRQ(ierr);
3322ad1e87fSLisandro Dalcin       }
3332ad1e87fSLisandro Dalcin       ilocal = NULL;
3342ad1e87fSLisandro Dalcin     }
33529046d53SLisandro Dalcin   } else {
33629046d53SLisandro Dalcin     sf->minleaf = 0;
33729046d53SLisandro Dalcin     sf->maxleaf = nleaves - 1;
33829046d53SLisandro Dalcin   }
33929046d53SLisandro Dalcin 
34029046d53SLisandro Dalcin   if (ilocal) {
34195fce210SBarry Smith     switch (localmode) {
34295fce210SBarry Smith     case PETSC_COPY_VALUES:
343785e854fSJed Brown       ierr = PetscMalloc1(nleaves,&sf->mine_alloc);CHKERRQ(ierr);
344580bdb30SBarry Smith       ierr = PetscArraycpy(sf->mine_alloc,ilocal,nleaves);CHKERRQ(ierr);
34595fce210SBarry Smith       sf->mine = sf->mine_alloc;
34695fce210SBarry Smith       break;
34795fce210SBarry Smith     case PETSC_OWN_POINTER:
34895fce210SBarry Smith       sf->mine_alloc = (PetscInt*)ilocal;
34995fce210SBarry Smith       sf->mine       = sf->mine_alloc;
35095fce210SBarry Smith       break;
35195fce210SBarry Smith     case PETSC_USE_POINTER:
35229046d53SLisandro Dalcin       sf->mine_alloc = NULL;
35395fce210SBarry Smith       sf->mine       = (PetscInt*)ilocal;
35495fce210SBarry Smith       break;
35595fce210SBarry Smith     default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown localmode");
35695fce210SBarry Smith     }
35795fce210SBarry Smith   }
35829046d53SLisandro Dalcin 
35995fce210SBarry Smith   switch (remotemode) {
36095fce210SBarry Smith   case PETSC_COPY_VALUES:
361785e854fSJed Brown     ierr = PetscMalloc1(nleaves,&sf->remote_alloc);CHKERRQ(ierr);
362580bdb30SBarry Smith     ierr = PetscArraycpy(sf->remote_alloc,iremote,nleaves);CHKERRQ(ierr);
36395fce210SBarry Smith     sf->remote = sf->remote_alloc;
36495fce210SBarry Smith     break;
36595fce210SBarry Smith   case PETSC_OWN_POINTER:
36695fce210SBarry Smith     sf->remote_alloc = (PetscSFNode*)iremote;
36795fce210SBarry Smith     sf->remote       = sf->remote_alloc;
36895fce210SBarry Smith     break;
36995fce210SBarry Smith   case PETSC_USE_POINTER:
37029046d53SLisandro Dalcin     sf->remote_alloc = NULL;
37195fce210SBarry Smith     sf->remote       = (PetscSFNode*)iremote;
37295fce210SBarry Smith     break;
37395fce210SBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown remotemode");
37495fce210SBarry Smith   }
37595fce210SBarry Smith 
376563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr);
37729046d53SLisandro Dalcin   sf->graphset = PETSC_TRUE;
37895fce210SBarry Smith   PetscFunctionReturn(0);
37995fce210SBarry Smith }
38095fce210SBarry Smith 
38129046d53SLisandro Dalcin /*@
38295fce210SBarry Smith    PetscSFCreateInverseSF - given a PetscSF in which all vertices have degree 1, creates the inverse map
38395fce210SBarry Smith 
38495fce210SBarry Smith    Collective
38595fce210SBarry Smith 
38695fce210SBarry Smith    Input Arguments:
38795fce210SBarry Smith .  sf - star forest to invert
38895fce210SBarry Smith 
38995fce210SBarry Smith    Output Arguments:
39095fce210SBarry Smith .  isf - inverse of sf
39195fce210SBarry Smith 
39295fce210SBarry Smith    Level: advanced
39395fce210SBarry Smith 
39495fce210SBarry Smith    Notes:
39595fce210SBarry Smith    All roots must have degree 1.
39695fce210SBarry Smith 
39795fce210SBarry Smith    The local space may be a permutation, but cannot be sparse.
39895fce210SBarry Smith 
39995fce210SBarry Smith .seealso: PetscSFSetGraph()
40095fce210SBarry Smith @*/
40195fce210SBarry Smith PetscErrorCode PetscSFCreateInverseSF(PetscSF sf,PetscSF *isf)
40295fce210SBarry Smith {
40395fce210SBarry Smith   PetscErrorCode ierr;
40495fce210SBarry Smith   PetscMPIInt    rank;
40595fce210SBarry Smith   PetscInt       i,nroots,nleaves,maxlocal,count,*newilocal;
40695fce210SBarry Smith   const PetscInt *ilocal;
40795fce210SBarry Smith   PetscSFNode    *roots,*leaves;
40895fce210SBarry Smith 
40995fce210SBarry Smith   PetscFunctionBegin;
41029046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
41129046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
41229046d53SLisandro Dalcin   PetscValidPointer(isf,2);
41329046d53SLisandro Dalcin 
41495fce210SBarry Smith   ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);CHKERRQ(ierr);
41529046d53SLisandro Dalcin   maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */
41629046d53SLisandro Dalcin 
41729046d53SLisandro Dalcin   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
418ae9aee6dSMatthew G. Knepley   ierr = PetscMalloc2(nroots,&roots,maxlocal,&leaves);CHKERRQ(ierr);
419ae9aee6dSMatthew G. Knepley   for (i=0; i<maxlocal; i++) {
42095fce210SBarry Smith     leaves[i].rank  = rank;
42195fce210SBarry Smith     leaves[i].index = i;
42295fce210SBarry Smith   }
42395fce210SBarry Smith   for (i=0; i <nroots; i++) {
42495fce210SBarry Smith     roots[i].rank  = -1;
42595fce210SBarry Smith     roots[i].index = -1;
42695fce210SBarry Smith   }
4278bfbc91cSJed Brown   ierr = PetscSFReduceBegin(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr);
4288bfbc91cSJed Brown   ierr = PetscSFReduceEnd(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr);
42995fce210SBarry Smith 
43095fce210SBarry Smith   /* Check whether our leaves are sparse */
43195fce210SBarry Smith   for (i=0,count=0; i<nroots; i++) if (roots[i].rank >= 0) count++;
43295fce210SBarry Smith   if (count == nroots) newilocal = NULL;
43395fce210SBarry Smith   else {                        /* Index for sparse leaves and compact "roots" array (which is to become our leaves). */
434785e854fSJed Brown     ierr = PetscMalloc1(count,&newilocal);CHKERRQ(ierr);
43595fce210SBarry Smith     for (i=0,count=0; i<nroots; i++) {
43695fce210SBarry Smith       if (roots[i].rank >= 0) {
43795fce210SBarry Smith         newilocal[count]   = i;
43895fce210SBarry Smith         roots[count].rank  = roots[i].rank;
43995fce210SBarry Smith         roots[count].index = roots[i].index;
44095fce210SBarry Smith         count++;
44195fce210SBarry Smith       }
44295fce210SBarry Smith     }
44395fce210SBarry Smith   }
44495fce210SBarry Smith 
44595fce210SBarry Smith   ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,isf);CHKERRQ(ierr);
44695fce210SBarry Smith   ierr = PetscSFSetGraph(*isf,maxlocal,count,newilocal,PETSC_OWN_POINTER,roots,PETSC_COPY_VALUES);CHKERRQ(ierr);
44795fce210SBarry Smith   ierr = PetscFree2(roots,leaves);CHKERRQ(ierr);
44895fce210SBarry Smith   PetscFunctionReturn(0);
44995fce210SBarry Smith }
45095fce210SBarry Smith 
45195fce210SBarry Smith /*@
45295fce210SBarry Smith    PetscSFDuplicate - duplicate a PetscSF, optionally preserving rank connectivity and graph
45395fce210SBarry Smith 
45495fce210SBarry Smith    Collective
45595fce210SBarry Smith 
45695fce210SBarry Smith    Input Arguments:
45795fce210SBarry Smith +  sf - communication object to duplicate
45895fce210SBarry Smith -  opt - PETSCSF_DUPLICATE_CONFONLY, PETSCSF_DUPLICATE_RANKS, or PETSCSF_DUPLICATE_GRAPH (see PetscSFDuplicateOption)
45995fce210SBarry Smith 
46095fce210SBarry Smith    Output Arguments:
46195fce210SBarry Smith .  newsf - new communication object
46295fce210SBarry Smith 
46395fce210SBarry Smith    Level: beginner
46495fce210SBarry Smith 
46595fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetType(), PetscSFSetGraph()
46695fce210SBarry Smith @*/
46795fce210SBarry Smith PetscErrorCode PetscSFDuplicate(PetscSF sf,PetscSFDuplicateOption opt,PetscSF *newsf)
46895fce210SBarry Smith {
46929046d53SLisandro Dalcin   PetscSFType    type;
47095fce210SBarry Smith   PetscErrorCode ierr;
47195fce210SBarry Smith 
47295fce210SBarry Smith   PetscFunctionBegin;
47329046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
47429046d53SLisandro Dalcin   PetscValidLogicalCollectiveEnum(sf,opt,2);
47529046d53SLisandro Dalcin   PetscValidPointer(newsf,3);
47695fce210SBarry Smith   ierr = PetscSFCreate(PetscObjectComm((PetscObject)sf),newsf);CHKERRQ(ierr);
47729046d53SLisandro Dalcin   ierr = PetscSFGetType(sf,&type);CHKERRQ(ierr);
47829046d53SLisandro Dalcin   if (type) {ierr = PetscSFSetType(*newsf,type);CHKERRQ(ierr);}
47995fce210SBarry Smith   if (opt == PETSCSF_DUPLICATE_GRAPH) {
48095fce210SBarry Smith     PetscInt          nroots,nleaves;
48195fce210SBarry Smith     const PetscInt    *ilocal;
48295fce210SBarry Smith     const PetscSFNode *iremote;
48329046d53SLisandro Dalcin     PetscSFCheckGraphSet(sf,1);
48495fce210SBarry Smith     ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
48595fce210SBarry Smith     ierr = PetscSFSetGraph(*newsf,nroots,nleaves,ilocal,PETSC_COPY_VALUES,iremote,PETSC_COPY_VALUES);CHKERRQ(ierr);
48695fce210SBarry Smith   }
48729046d53SLisandro Dalcin   if (sf->ops->Duplicate) {ierr = (*sf->ops->Duplicate)(sf,opt,*newsf);CHKERRQ(ierr);}
48895fce210SBarry Smith   PetscFunctionReturn(0);
48995fce210SBarry Smith }
49095fce210SBarry Smith 
49195fce210SBarry Smith /*@C
49295fce210SBarry Smith    PetscSFGetGraph - Get the graph specifying a parallel star forest
49395fce210SBarry Smith 
49495fce210SBarry Smith    Not Collective
49595fce210SBarry Smith 
49695fce210SBarry Smith    Input Arguments:
49795fce210SBarry Smith .  sf - star forest
49895fce210SBarry Smith 
49995fce210SBarry Smith    Output Arguments:
50095fce210SBarry Smith +  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
50195fce210SBarry Smith .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
50295fce210SBarry Smith .  ilocal - locations of leaves in leafdata buffers
50395fce210SBarry Smith -  iremote - remote locations of root vertices for each leaf on the current process
50495fce210SBarry Smith 
505373e0d91SLisandro Dalcin    Notes:
506373e0d91SLisandro Dalcin    We are not currently requiring that the graph is set, thus returning nroots=-1 if it has not been set yet
507373e0d91SLisandro Dalcin 
50895fce210SBarry Smith    Level: intermediate
50995fce210SBarry Smith 
51095fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph()
51195fce210SBarry Smith @*/
51295fce210SBarry Smith PetscErrorCode PetscSFGetGraph(PetscSF sf,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote)
51395fce210SBarry Smith {
51495fce210SBarry Smith 
51595fce210SBarry Smith   PetscFunctionBegin;
51695fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
51795fce210SBarry Smith   if (nroots) *nroots = sf->nroots;
51895fce210SBarry Smith   if (nleaves) *nleaves = sf->nleaves;
51995fce210SBarry Smith   if (ilocal) *ilocal = sf->mine;
52095fce210SBarry Smith   if (iremote) *iremote = sf->remote;
52195fce210SBarry Smith   PetscFunctionReturn(0);
52295fce210SBarry Smith }
52395fce210SBarry Smith 
52429046d53SLisandro Dalcin /*@
52595fce210SBarry Smith    PetscSFGetLeafRange - Get the active leaf ranges
52695fce210SBarry Smith 
52795fce210SBarry Smith    Not Collective
52895fce210SBarry Smith 
52995fce210SBarry Smith    Input Arguments:
53095fce210SBarry Smith .  sf - star forest
53195fce210SBarry Smith 
53295fce210SBarry Smith    Output Arguments:
53395fce210SBarry Smith +  minleaf - minimum active leaf on this process
53495fce210SBarry Smith -  maxleaf - maximum active leaf on this process
53595fce210SBarry Smith 
53695fce210SBarry Smith    Level: developer
53795fce210SBarry Smith 
53895fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph(), PetscSFGetGraph()
53995fce210SBarry Smith @*/
54095fce210SBarry Smith PetscErrorCode PetscSFGetLeafRange(PetscSF sf,PetscInt *minleaf,PetscInt *maxleaf)
54195fce210SBarry Smith {
54295fce210SBarry Smith 
54395fce210SBarry Smith   PetscFunctionBegin;
54495fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
54529046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
54695fce210SBarry Smith   if (minleaf) *minleaf = sf->minleaf;
54795fce210SBarry Smith   if (maxleaf) *maxleaf = sf->maxleaf;
54895fce210SBarry Smith   PetscFunctionReturn(0);
54995fce210SBarry Smith }
55095fce210SBarry Smith 
55195fce210SBarry Smith /*@C
55295fce210SBarry Smith    PetscSFView - view a star forest
55395fce210SBarry Smith 
55495fce210SBarry Smith    Collective
55595fce210SBarry Smith 
55695fce210SBarry Smith    Input Arguments:
55795fce210SBarry Smith +  sf - star forest
55895fce210SBarry Smith -  viewer - viewer to display graph, for example PETSC_VIEWER_STDOUT_WORLD
55995fce210SBarry Smith 
56095fce210SBarry Smith    Level: beginner
56195fce210SBarry Smith 
56295fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph()
56395fce210SBarry Smith @*/
56495fce210SBarry Smith PetscErrorCode PetscSFView(PetscSF sf,PetscViewer viewer)
56595fce210SBarry Smith {
56695fce210SBarry Smith   PetscErrorCode    ierr;
56795fce210SBarry Smith   PetscBool         iascii;
56895fce210SBarry Smith   PetscViewerFormat format;
56995fce210SBarry Smith 
57095fce210SBarry Smith   PetscFunctionBegin;
57195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
57295fce210SBarry Smith   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sf),&viewer);CHKERRQ(ierr);}
57395fce210SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
57495fce210SBarry Smith   PetscCheckSameComm(sf,1,viewer,2);
57580153354SVaclav Hapla   if (sf->graphset) {ierr = PetscSFSetUp(sf);CHKERRQ(ierr);}
57695fce210SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
57795fce210SBarry Smith   if (iascii) {
57895fce210SBarry Smith     PetscMPIInt rank;
57981bfa7aaSJed Brown     PetscInt    ii,i,j;
58095fce210SBarry Smith 
581dae58748SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)sf,viewer);CHKERRQ(ierr);
58295fce210SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
58395fce210SBarry Smith     if (sf->ops->View) {ierr = (*sf->ops->View)(sf,viewer);CHKERRQ(ierr);}
58480153354SVaclav Hapla     if (!sf->graphset) {
58580153354SVaclav Hapla       ierr = PetscViewerASCIIPrintf(viewer,"PetscSFSetGraph() has not been called yet\n");CHKERRQ(ierr);
58680153354SVaclav Hapla       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
58780153354SVaclav Hapla       PetscFunctionReturn(0);
58880153354SVaclav Hapla     }
58995fce210SBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
5901575c14dSBarry Smith     ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
59195fce210SBarry Smith     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of roots=%D, leaves=%D, remote ranks=%D\n",rank,sf->nroots,sf->nleaves,sf->nranks);CHKERRQ(ierr);
59295fce210SBarry Smith     for (i=0; i<sf->nleaves; i++) {
59395fce210SBarry 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);
59495fce210SBarry Smith     }
59595fce210SBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
59695fce210SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
59795fce210SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
59881bfa7aaSJed Brown       PetscMPIInt *tmpranks,*perm;
59981bfa7aaSJed Brown       ierr = PetscMalloc2(sf->nranks,&tmpranks,sf->nranks,&perm);CHKERRQ(ierr);
600580bdb30SBarry Smith       ierr = PetscArraycpy(tmpranks,sf->ranks,sf->nranks);CHKERRQ(ierr);
60181bfa7aaSJed Brown       for (i=0; i<sf->nranks; i++) perm[i] = i;
60281bfa7aaSJed Brown       ierr = PetscSortMPIIntWithArray(sf->nranks,tmpranks,perm);CHKERRQ(ierr);
60395fce210SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Roots referenced by my leaves, by rank\n",rank);CHKERRQ(ierr);
60481bfa7aaSJed Brown       for (ii=0; ii<sf->nranks; ii++) {
60581bfa7aaSJed Brown         i = perm[ii];
6067904a332SBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %d: %D edges\n",rank,sf->ranks[i],sf->roffset[i+1]-sf->roffset[i]);CHKERRQ(ierr);
60795fce210SBarry Smith         for (j=sf->roffset[i]; j<sf->roffset[i+1]; j++) {
60895fce210SBarry Smith           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d]    %D <- %D\n",rank,sf->rmine[j],sf->rremote[j]);CHKERRQ(ierr);
60995fce210SBarry Smith         }
61095fce210SBarry Smith       }
61181bfa7aaSJed Brown       ierr = PetscFree2(tmpranks,perm);CHKERRQ(ierr);
61295fce210SBarry Smith     }
61395fce210SBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
6141575c14dSBarry Smith     ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
61595fce210SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
61695fce210SBarry Smith   }
61795fce210SBarry Smith   PetscFunctionReturn(0);
61895fce210SBarry Smith }
61995fce210SBarry Smith 
62095fce210SBarry Smith /*@C
621dec1416fSJunchao Zhang    PetscSFGetRootRanks - Get root ranks and number of vertices referenced by leaves on this process
62295fce210SBarry Smith 
62395fce210SBarry Smith    Not Collective
62495fce210SBarry Smith 
62595fce210SBarry Smith    Input Arguments:
62695fce210SBarry Smith .  sf - star forest
62795fce210SBarry Smith 
62895fce210SBarry Smith    Output Arguments:
62995fce210SBarry Smith +  nranks - number of ranks referenced by local part
63095fce210SBarry Smith .  ranks - array of ranks
63195fce210SBarry Smith .  roffset - offset in rmine/rremote for each rank (length nranks+1)
63295fce210SBarry Smith .  rmine - concatenated array holding local indices referencing each remote rank
63395fce210SBarry Smith -  rremote - concatenated array holding remote indices referenced for each remote rank
63495fce210SBarry Smith 
63595fce210SBarry Smith    Level: developer
63695fce210SBarry Smith 
637dec1416fSJunchao Zhang .seealso: PetscSFGetLeafRanks()
63895fce210SBarry Smith @*/
639dec1416fSJunchao Zhang PetscErrorCode PetscSFGetRootRanks(PetscSF sf,PetscInt *nranks,const PetscMPIInt **ranks,const PetscInt **roffset,const PetscInt **rmine,const PetscInt **rremote)
64095fce210SBarry Smith {
641dec1416fSJunchao Zhang   PetscErrorCode ierr;
64295fce210SBarry Smith 
64395fce210SBarry Smith   PetscFunctionBegin;
64495fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
64529046d53SLisandro Dalcin   if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining ranks");
646dec1416fSJunchao Zhang   if (sf->ops->GetRootRanks) {
647dec1416fSJunchao Zhang     ierr = (sf->ops->GetRootRanks)(sf,nranks,ranks,roffset,rmine,rremote);CHKERRQ(ierr);
648dec1416fSJunchao Zhang   } else {
649dec1416fSJunchao Zhang     /* The generic implementation */
65095fce210SBarry Smith     if (nranks)  *nranks  = sf->nranks;
65195fce210SBarry Smith     if (ranks)   *ranks   = sf->ranks;
65295fce210SBarry Smith     if (roffset) *roffset = sf->roffset;
65395fce210SBarry Smith     if (rmine)   *rmine   = sf->rmine;
65495fce210SBarry Smith     if (rremote) *rremote = sf->rremote;
655dec1416fSJunchao Zhang   }
65695fce210SBarry Smith   PetscFunctionReturn(0);
65795fce210SBarry Smith }
65895fce210SBarry Smith 
6598750ddebSJunchao Zhang /*@C
6608750ddebSJunchao Zhang    PetscSFGetLeafRanks - Get leaf ranks referencing roots on this process
6618750ddebSJunchao Zhang 
6628750ddebSJunchao Zhang    Not Collective
6638750ddebSJunchao Zhang 
6648750ddebSJunchao Zhang    Input Arguments:
6658750ddebSJunchao Zhang .  sf - star forest
6668750ddebSJunchao Zhang 
6678750ddebSJunchao Zhang    Output Arguments:
6688750ddebSJunchao Zhang +  niranks - number of leaf ranks referencing roots on this process
6698750ddebSJunchao Zhang .  iranks - array of ranks
6708750ddebSJunchao Zhang .  ioffset - offset in irootloc for each rank (length niranks+1)
6718750ddebSJunchao Zhang -  irootloc - concatenated array holding local indices of roots referenced by each leaf rank
6728750ddebSJunchao Zhang 
6738750ddebSJunchao Zhang    Level: developer
6748750ddebSJunchao Zhang 
675dec1416fSJunchao Zhang .seealso: PetscSFGetRootRanks()
6768750ddebSJunchao Zhang @*/
6778750ddebSJunchao Zhang PetscErrorCode PetscSFGetLeafRanks(PetscSF sf,PetscInt *niranks,const PetscMPIInt **iranks,const PetscInt **ioffset,const PetscInt **irootloc)
6788750ddebSJunchao Zhang {
6798750ddebSJunchao Zhang   PetscErrorCode ierr;
6808750ddebSJunchao Zhang 
6818750ddebSJunchao Zhang   PetscFunctionBegin;
6828750ddebSJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
6838750ddebSJunchao Zhang   if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining ranks");
6848750ddebSJunchao Zhang   if (sf->ops->GetLeafRanks) {
6858750ddebSJunchao Zhang     ierr = (sf->ops->GetLeafRanks)(sf,niranks,iranks,ioffset,irootloc);CHKERRQ(ierr);
6868750ddebSJunchao Zhang   } else {
6878750ddebSJunchao Zhang     PetscSFType type;
6888750ddebSJunchao Zhang     ierr = PetscSFGetType(sf,&type);CHKERRQ(ierr);
6898750ddebSJunchao Zhang     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSFGetLeafRanks() is not supported on this StarForest type: %s", type);
6908750ddebSJunchao Zhang   }
6918750ddebSJunchao Zhang   PetscFunctionReturn(0);
6928750ddebSJunchao Zhang }
6938750ddebSJunchao Zhang 
694b5a8e515SJed Brown static PetscBool InList(PetscMPIInt needle,PetscMPIInt n,const PetscMPIInt *list) {
695b5a8e515SJed Brown   PetscInt i;
696b5a8e515SJed Brown   for (i=0; i<n; i++) {
697b5a8e515SJed Brown     if (needle == list[i]) return PETSC_TRUE;
698b5a8e515SJed Brown   }
699b5a8e515SJed Brown   return PETSC_FALSE;
700b5a8e515SJed Brown }
701b5a8e515SJed Brown 
70295fce210SBarry Smith /*@C
70321c688dcSJed Brown    PetscSFSetUpRanks - Set up data structures associated with ranks; this is for internal use by PetscSF implementations.
70421c688dcSJed Brown 
70521c688dcSJed Brown    Collective
70621c688dcSJed Brown 
70721c688dcSJed Brown    Input Arguments:
708b5a8e515SJed Brown +  sf - PetscSF to set up; PetscSFSetGraph() must have been called
709b5a8e515SJed Brown -  dgroup - MPI_Group of ranks to be distinguished (e.g., for self or shared memory exchange)
71021c688dcSJed Brown 
71121c688dcSJed Brown    Level: developer
71221c688dcSJed Brown 
713dec1416fSJunchao Zhang .seealso: PetscSFGetRootRanks()
71421c688dcSJed Brown @*/
715b5a8e515SJed Brown PetscErrorCode PetscSFSetUpRanks(PetscSF sf,MPI_Group dgroup)
71621c688dcSJed Brown {
71721c688dcSJed Brown   PetscErrorCode     ierr;
71821c688dcSJed Brown   PetscTable         table;
71921c688dcSJed Brown   PetscTablePosition pos;
720b5a8e515SJed Brown   PetscMPIInt        size,groupsize,*groupranks;
721247e8311SStefano Zampini   PetscInt           *rcount,*ranks;
722247e8311SStefano Zampini   PetscInt           i, irank = -1,orank = -1;
72321c688dcSJed Brown 
72421c688dcSJed Brown   PetscFunctionBegin;
72521c688dcSJed Brown   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
72629046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
72721c688dcSJed Brown   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr);
72821c688dcSJed Brown   ierr = PetscTableCreate(10,size,&table);CHKERRQ(ierr);
72921c688dcSJed Brown   for (i=0; i<sf->nleaves; i++) {
73021c688dcSJed Brown     /* Log 1-based rank */
73121c688dcSJed Brown     ierr = PetscTableAdd(table,sf->remote[i].rank+1,1,ADD_VALUES);CHKERRQ(ierr);
73221c688dcSJed Brown   }
73321c688dcSJed Brown   ierr = PetscTableGetCount(table,&sf->nranks);CHKERRQ(ierr);
73421c688dcSJed Brown   ierr = PetscMalloc4(sf->nranks,&sf->ranks,sf->nranks+1,&sf->roffset,sf->nleaves,&sf->rmine,sf->nleaves,&sf->rremote);CHKERRQ(ierr);
73521c688dcSJed Brown   ierr = PetscMalloc2(sf->nranks,&rcount,sf->nranks,&ranks);CHKERRQ(ierr);
73621c688dcSJed Brown   ierr = PetscTableGetHeadPosition(table,&pos);CHKERRQ(ierr);
73721c688dcSJed Brown   for (i=0; i<sf->nranks; i++) {
73821c688dcSJed Brown     ierr = PetscTableGetNext(table,&pos,&ranks[i],&rcount[i]);CHKERRQ(ierr);
73921c688dcSJed Brown     ranks[i]--;             /* Convert back to 0-based */
74021c688dcSJed Brown   }
74121c688dcSJed Brown   ierr = PetscTableDestroy(&table);CHKERRQ(ierr);
742b5a8e515SJed Brown 
743b5a8e515SJed Brown   /* We expect that dgroup is reliably "small" while nranks could be large */
744b5a8e515SJed Brown   {
7457fb8a5e4SKarl Rupp     MPI_Group group = MPI_GROUP_NULL;
746b5a8e515SJed Brown     PetscMPIInt *dgroupranks;
747b5a8e515SJed Brown     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
748b5a8e515SJed Brown     ierr = MPI_Group_size(dgroup,&groupsize);CHKERRQ(ierr);
749b5a8e515SJed Brown     ierr = PetscMalloc1(groupsize,&dgroupranks);CHKERRQ(ierr);
750b5a8e515SJed Brown     ierr = PetscMalloc1(groupsize,&groupranks);CHKERRQ(ierr);
751b5a8e515SJed Brown     for (i=0; i<groupsize; i++) dgroupranks[i] = i;
752f3fc9a17SJed Brown     if (groupsize) {ierr = MPI_Group_translate_ranks(dgroup,groupsize,dgroupranks,group,groupranks);CHKERRQ(ierr);}
753b5a8e515SJed Brown     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
754b5a8e515SJed Brown     ierr = PetscFree(dgroupranks);CHKERRQ(ierr);
755b5a8e515SJed Brown   }
756b5a8e515SJed Brown 
757b5a8e515SJed Brown   /* Partition ranks[] into distinguished (first sf->ndranks) followed by non-distinguished */
758b5a8e515SJed Brown   for (sf->ndranks=0,i=sf->nranks; sf->ndranks<i; ) {
759b5a8e515SJed Brown     for (i--; sf->ndranks<i; i--) { /* Scan i backward looking for distinguished rank */
760b5a8e515SJed Brown       if (InList(ranks[i],groupsize,groupranks)) break;
761b5a8e515SJed Brown     }
762b5a8e515SJed Brown     for ( ; sf->ndranks<=i; sf->ndranks++) { /* Scan sf->ndranks forward looking for non-distinguished rank */
763b5a8e515SJed Brown       if (!InList(ranks[sf->ndranks],groupsize,groupranks)) break;
764b5a8e515SJed Brown     }
765b5a8e515SJed Brown     if (sf->ndranks < i) {                         /* Swap ranks[sf->ndranks] with ranks[i] */
766b5a8e515SJed Brown       PetscInt    tmprank,tmpcount;
767247e8311SStefano Zampini 
768b5a8e515SJed Brown       tmprank             = ranks[i];
769b5a8e515SJed Brown       tmpcount            = rcount[i];
770b5a8e515SJed Brown       ranks[i]            = ranks[sf->ndranks];
771b5a8e515SJed Brown       rcount[i]           = rcount[sf->ndranks];
772b5a8e515SJed Brown       ranks[sf->ndranks]  = tmprank;
773b5a8e515SJed Brown       rcount[sf->ndranks] = tmpcount;
774b5a8e515SJed Brown       sf->ndranks++;
775b5a8e515SJed Brown     }
776b5a8e515SJed Brown   }
777b5a8e515SJed Brown   ierr = PetscFree(groupranks);CHKERRQ(ierr);
778b5a8e515SJed Brown   ierr = PetscSortIntWithArray(sf->ndranks,ranks,rcount);CHKERRQ(ierr);
779b5a8e515SJed Brown   ierr = PetscSortIntWithArray(sf->nranks-sf->ndranks,ranks+sf->ndranks,rcount+sf->ndranks);CHKERRQ(ierr);
78021c688dcSJed Brown   sf->roffset[0] = 0;
78121c688dcSJed Brown   for (i=0; i<sf->nranks; i++) {
78221c688dcSJed Brown     ierr = PetscMPIIntCast(ranks[i],sf->ranks+i);CHKERRQ(ierr);
78321c688dcSJed Brown     sf->roffset[i+1] = sf->roffset[i] + rcount[i];
78421c688dcSJed Brown     rcount[i]        = 0;
78521c688dcSJed Brown   }
786247e8311SStefano Zampini   for (i=0, irank = -1, orank = -1; i<sf->nleaves; i++) {
787247e8311SStefano Zampini     /* short circuit */
788247e8311SStefano Zampini     if (orank != sf->remote[i].rank) {
78921c688dcSJed Brown       /* Search for index of iremote[i].rank in sf->ranks */
790b5a8e515SJed Brown       ierr = PetscFindMPIInt(sf->remote[i].rank,sf->ndranks,sf->ranks,&irank);CHKERRQ(ierr);
791b5a8e515SJed Brown       if (irank < 0) {
792b5a8e515SJed Brown         ierr = PetscFindMPIInt(sf->remote[i].rank,sf->nranks-sf->ndranks,sf->ranks+sf->ndranks,&irank);CHKERRQ(ierr);
793b5a8e515SJed Brown         if (irank >= 0) irank += sf->ndranks;
79421c688dcSJed Brown       }
795247e8311SStefano Zampini       orank = sf->remote[i].rank;
796247e8311SStefano Zampini     }
797b5a8e515SJed Brown     if (irank < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Could not find rank %D in array",sf->remote[i].rank);
79821c688dcSJed Brown     sf->rmine[sf->roffset[irank] + rcount[irank]]   = sf->mine ? sf->mine[i] : i;
79921c688dcSJed Brown     sf->rremote[sf->roffset[irank] + rcount[irank]] = sf->remote[i].index;
80021c688dcSJed Brown     rcount[irank]++;
80121c688dcSJed Brown   }
80221c688dcSJed Brown   ierr = PetscFree2(rcount,ranks);CHKERRQ(ierr);
80321c688dcSJed Brown   PetscFunctionReturn(0);
80421c688dcSJed Brown }
80521c688dcSJed Brown 
80621c688dcSJed Brown /*@C
80795fce210SBarry Smith    PetscSFGetGroups - gets incoming and outgoing process groups
80895fce210SBarry Smith 
80995fce210SBarry Smith    Collective
81095fce210SBarry Smith 
81195fce210SBarry Smith    Input Argument:
81295fce210SBarry Smith .  sf - star forest
81395fce210SBarry Smith 
81495fce210SBarry Smith    Output Arguments:
81595fce210SBarry Smith +  incoming - group of origin processes for incoming edges (leaves that reference my roots)
81695fce210SBarry Smith -  outgoing - group of destination processes for outgoing edges (roots that I reference)
81795fce210SBarry Smith 
81895fce210SBarry Smith    Level: developer
81995fce210SBarry Smith 
82095fce210SBarry Smith .seealso: PetscSFGetWindow(), PetscSFRestoreWindow()
82195fce210SBarry Smith @*/
82295fce210SBarry Smith PetscErrorCode PetscSFGetGroups(PetscSF sf,MPI_Group *incoming,MPI_Group *outgoing)
82395fce210SBarry Smith {
82495fce210SBarry Smith   PetscErrorCode ierr;
8257fb8a5e4SKarl Rupp   MPI_Group      group = MPI_GROUP_NULL;
82695fce210SBarry Smith 
82795fce210SBarry Smith   PetscFunctionBegin;
82829046d53SLisandro Dalcin   if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining groups");
82995fce210SBarry Smith   if (sf->ingroup == MPI_GROUP_NULL) {
83095fce210SBarry Smith     PetscInt       i;
83195fce210SBarry Smith     const PetscInt *indegree;
83295fce210SBarry Smith     PetscMPIInt    rank,*outranks,*inranks;
83395fce210SBarry Smith     PetscSFNode    *remote;
83495fce210SBarry Smith     PetscSF        bgcount;
83595fce210SBarry Smith 
83695fce210SBarry Smith     /* Compute the number of incoming ranks */
837785e854fSJed Brown     ierr = PetscMalloc1(sf->nranks,&remote);CHKERRQ(ierr);
83895fce210SBarry Smith     for (i=0; i<sf->nranks; i++) {
83995fce210SBarry Smith       remote[i].rank  = sf->ranks[i];
84095fce210SBarry Smith       remote[i].index = 0;
84195fce210SBarry Smith     }
84295fce210SBarry Smith     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,&bgcount);CHKERRQ(ierr);
84395fce210SBarry Smith     ierr = PetscSFSetGraph(bgcount,1,sf->nranks,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
84495fce210SBarry Smith     ierr = PetscSFComputeDegreeBegin(bgcount,&indegree);CHKERRQ(ierr);
84595fce210SBarry Smith     ierr = PetscSFComputeDegreeEnd(bgcount,&indegree);CHKERRQ(ierr);
84695fce210SBarry Smith 
84795fce210SBarry Smith     /* Enumerate the incoming ranks */
848dcca6d9dSJed Brown     ierr = PetscMalloc2(indegree[0],&inranks,sf->nranks,&outranks);CHKERRQ(ierr);
84995fce210SBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
85095fce210SBarry Smith     for (i=0; i<sf->nranks; i++) outranks[i] = rank;
85195fce210SBarry Smith     ierr = PetscSFGatherBegin(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr);
85295fce210SBarry Smith     ierr = PetscSFGatherEnd(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr);
85395fce210SBarry Smith     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
85495fce210SBarry Smith     ierr = MPI_Group_incl(group,indegree[0],inranks,&sf->ingroup);CHKERRQ(ierr);
85595fce210SBarry Smith     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
85695fce210SBarry Smith     ierr = PetscFree2(inranks,outranks);CHKERRQ(ierr);
85795fce210SBarry Smith     ierr = PetscSFDestroy(&bgcount);CHKERRQ(ierr);
85895fce210SBarry Smith   }
85995fce210SBarry Smith   *incoming = sf->ingroup;
86095fce210SBarry Smith 
86195fce210SBarry Smith   if (sf->outgroup == MPI_GROUP_NULL) {
86295fce210SBarry Smith     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
86395fce210SBarry Smith     ierr = MPI_Group_incl(group,sf->nranks,sf->ranks,&sf->outgroup);CHKERRQ(ierr);
86495fce210SBarry Smith     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
86595fce210SBarry Smith   }
86695fce210SBarry Smith   *outgoing = sf->outgroup;
86795fce210SBarry Smith   PetscFunctionReturn(0);
86895fce210SBarry Smith }
86995fce210SBarry Smith 
87029046d53SLisandro Dalcin /*@
87195fce210SBarry Smith    PetscSFGetMultiSF - gets the inner SF implemeting gathers and scatters
87295fce210SBarry Smith 
87395fce210SBarry Smith    Collective
87495fce210SBarry Smith 
87595fce210SBarry Smith    Input Argument:
87695fce210SBarry Smith .  sf - star forest that may contain roots with 0 or with more than 1 vertex
87795fce210SBarry Smith 
87895fce210SBarry Smith    Output Arguments:
87995fce210SBarry Smith .  multi - star forest with split roots, such that each root has degree exactly 1
88095fce210SBarry Smith 
88195fce210SBarry Smith    Level: developer
88295fce210SBarry Smith 
88395fce210SBarry Smith    Notes:
88495fce210SBarry Smith 
88595fce210SBarry Smith    In most cases, users should use PetscSFGatherBegin() and PetscSFScatterBegin() instead of manipulating multi
88695fce210SBarry Smith    directly. Since multi satisfies the stronger condition that each entry in the global space has exactly one incoming
88795fce210SBarry Smith    edge, it is a candidate for future optimization that might involve its removal.
88895fce210SBarry Smith 
889673100f5SVaclav Hapla .seealso: PetscSFSetGraph(), PetscSFGatherBegin(), PetscSFScatterBegin(), PetscSFComputeMultiRootOriginalNumbering()
89095fce210SBarry Smith @*/
89195fce210SBarry Smith PetscErrorCode PetscSFGetMultiSF(PetscSF sf,PetscSF *multi)
89295fce210SBarry Smith {
89395fce210SBarry Smith   PetscErrorCode ierr;
89495fce210SBarry Smith 
89595fce210SBarry Smith   PetscFunctionBegin;
89695fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
89795fce210SBarry Smith   PetscValidPointer(multi,2);
89895fce210SBarry Smith   if (sf->nroots < 0) {         /* Graph has not been set yet; why do we need this? */
89995fce210SBarry Smith     ierr   = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr);
90095fce210SBarry Smith     *multi = sf->multi;
90195fce210SBarry Smith     PetscFunctionReturn(0);
90295fce210SBarry Smith   }
90395fce210SBarry Smith   if (!sf->multi) {
90495fce210SBarry Smith     const PetscInt *indegree;
9059837ea96SMatthew G. Knepley     PetscInt       i,*inoffset,*outones,*outoffset,maxlocal;
90695fce210SBarry Smith     PetscSFNode    *remote;
90729046d53SLisandro Dalcin     maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */
90895fce210SBarry Smith     ierr = PetscSFComputeDegreeBegin(sf,&indegree);CHKERRQ(ierr);
90995fce210SBarry Smith     ierr = PetscSFComputeDegreeEnd(sf,&indegree);CHKERRQ(ierr);
9109837ea96SMatthew G. Knepley     ierr = PetscMalloc3(sf->nroots+1,&inoffset,maxlocal,&outones,maxlocal,&outoffset);CHKERRQ(ierr);
91195fce210SBarry Smith     inoffset[0] = 0;
91295fce210SBarry Smith     for (i=0; i<sf->nroots; i++) inoffset[i+1] = inoffset[i] + indegree[i];
9139837ea96SMatthew G. Knepley     for (i=0; i<maxlocal; i++) outones[i] = 1;
914dbd2ff41SBarry Smith     ierr = PetscSFFetchAndOpBegin(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr);
915dbd2ff41SBarry Smith     ierr = PetscSFFetchAndOpEnd(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr);
91695fce210SBarry Smith     for (i=0; i<sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */
91795fce210SBarry Smith #if 0
91895fce210SBarry Smith #if defined(PETSC_USE_DEBUG)                                 /* Check that the expected number of increments occurred */
91995fce210SBarry Smith     for (i=0; i<sf->nroots; i++) {
92095fce210SBarry Smith       if (inoffset[i] + indegree[i] != inoffset[i+1]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect result after PetscSFFetchAndOp");
92195fce210SBarry Smith     }
92295fce210SBarry Smith #endif
92395fce210SBarry Smith #endif
924785e854fSJed Brown     ierr = PetscMalloc1(sf->nleaves,&remote);CHKERRQ(ierr);
92595fce210SBarry Smith     for (i=0; i<sf->nleaves; i++) {
92695fce210SBarry Smith       remote[i].rank  = sf->remote[i].rank;
92738e7336fSToby Isaac       remote[i].index = outoffset[sf->mine ? sf->mine[i] : i];
92895fce210SBarry Smith     }
92995fce210SBarry Smith     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr);
93001365b40SToby Isaac     ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
93195fce210SBarry Smith     if (sf->rankorder) {        /* Sort the ranks */
93295fce210SBarry Smith       PetscMPIInt rank;
93395fce210SBarry Smith       PetscInt    *inranks,*newoffset,*outranks,*newoutoffset,*tmpoffset,maxdegree;
93495fce210SBarry Smith       PetscSFNode *newremote;
93595fce210SBarry Smith       ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
93695fce210SBarry Smith       for (i=0,maxdegree=0; i<sf->nroots; i++) maxdegree = PetscMax(maxdegree,indegree[i]);
9379837ea96SMatthew G. Knepley       ierr = PetscMalloc5(sf->multi->nroots,&inranks,sf->multi->nroots,&newoffset,maxlocal,&outranks,maxlocal,&newoutoffset,maxdegree,&tmpoffset);CHKERRQ(ierr);
9389837ea96SMatthew G. Knepley       for (i=0; i<maxlocal; i++) outranks[i] = rank;
9398bfbc91cSJed Brown       ierr = PetscSFReduceBegin(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr);
9408bfbc91cSJed Brown       ierr = PetscSFReduceEnd(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr);
94195fce210SBarry Smith       /* Sort the incoming ranks at each vertex, build the inverse map */
94295fce210SBarry Smith       for (i=0; i<sf->nroots; i++) {
94395fce210SBarry Smith         PetscInt j;
94495fce210SBarry Smith         for (j=0; j<indegree[i]; j++) tmpoffset[j] = j;
94595fce210SBarry Smith         ierr = PetscSortIntWithArray(indegree[i],inranks+inoffset[i],tmpoffset);CHKERRQ(ierr);
94695fce210SBarry Smith         for (j=0; j<indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j;
94795fce210SBarry Smith       }
94895fce210SBarry Smith       ierr = PetscSFBcastBegin(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr);
94995fce210SBarry Smith       ierr = PetscSFBcastEnd(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr);
950785e854fSJed Brown       ierr = PetscMalloc1(sf->nleaves,&newremote);CHKERRQ(ierr);
95195fce210SBarry Smith       for (i=0; i<sf->nleaves; i++) {
95295fce210SBarry Smith         newremote[i].rank  = sf->remote[i].rank;
95301365b40SToby Isaac         newremote[i].index = newoutoffset[sf->mine ? sf->mine[i] : i];
95495fce210SBarry Smith       }
95501365b40SToby Isaac       ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,newremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
95695fce210SBarry Smith       ierr = PetscFree5(inranks,newoffset,outranks,newoutoffset,tmpoffset);CHKERRQ(ierr);
95795fce210SBarry Smith     }
95895fce210SBarry Smith     ierr = PetscFree3(inoffset,outones,outoffset);CHKERRQ(ierr);
95995fce210SBarry Smith   }
96095fce210SBarry Smith   *multi = sf->multi;
96195fce210SBarry Smith   PetscFunctionReturn(0);
96295fce210SBarry Smith }
96395fce210SBarry Smith 
96495fce210SBarry Smith /*@C
96595fce210SBarry Smith    PetscSFCreateEmbeddedSF - removes edges from all but the selected roots, does not remap indices
96695fce210SBarry Smith 
96795fce210SBarry Smith    Collective
96895fce210SBarry Smith 
96995fce210SBarry Smith    Input Arguments:
97095fce210SBarry Smith +  sf - original star forest
971ba2a7774SJunchao Zhang .  nselected  - number of selected roots on this process
972ba2a7774SJunchao Zhang -  selected   - indices of the selected roots on this process
97395fce210SBarry Smith 
97495fce210SBarry Smith    Output Arguments:
97595fce210SBarry Smith .  newsf - new star forest
97695fce210SBarry Smith 
97795fce210SBarry Smith    Level: advanced
97895fce210SBarry Smith 
97995fce210SBarry Smith    Note:
98095fce210SBarry Smith    To use the new PetscSF, it may be necessary to know the indices of the leaves that are still participating. This can
98195fce210SBarry Smith    be done by calling PetscSFGetGraph().
98295fce210SBarry Smith 
98395fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFGetGraph()
98495fce210SBarry Smith @*/
985ba2a7774SJunchao Zhang PetscErrorCode PetscSFCreateEmbeddedSF(PetscSF sf,PetscInt nselected,const PetscInt *selected,PetscSF *newsf)
98695fce210SBarry Smith {
987f659e5c7SJunchao Zhang   PetscInt          i,n,*roots,*rootdata,*leafdata,nroots,nleaves,connected_leaves,*new_ilocal;
988ba2a7774SJunchao Zhang   const PetscSFNode *iremote;
989f659e5c7SJunchao Zhang   PetscSFNode       *new_iremote;
990ba2a7774SJunchao Zhang   PetscSF           tmpsf;
991f659e5c7SJunchao Zhang   MPI_Comm          comm;
9920511a646SMatthew G. Knepley   PetscErrorCode    ierr;
99395fce210SBarry Smith 
99495fce210SBarry Smith   PetscFunctionBegin;
99595fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
99629046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
997ba2a7774SJunchao Zhang   if (nselected) PetscValidPointer(selected,3);
99895fce210SBarry Smith   PetscValidPointer(newsf,4);
9990511a646SMatthew G. Knepley 
1000f659e5c7SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1001f659e5c7SJunchao Zhang 
1002f659e5c7SJunchao Zhang   /* Uniq selected[] and put results in roots[] */
1003f659e5c7SJunchao Zhang   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
1004f659e5c7SJunchao Zhang   ierr = PetscMalloc1(nselected,&roots);CHKERRQ(ierr);
1005f659e5c7SJunchao Zhang   ierr = PetscMemcpy(roots,selected,sizeof(PetscInt)*nselected);CHKERRQ(ierr);
1006f659e5c7SJunchao Zhang   ierr = PetscSortedRemoveDupsInt(&nselected,roots);CHKERRQ(ierr);
1007f659e5c7SJunchao Zhang   if (nselected && (roots[0] < 0 || roots[nselected-1] >= sf->nroots)) SETERRQ3(comm,PETSC_ERR_ARG_OUTOFRANGE,"Min/Max root indices %D/%D are not in [0,%D)",roots[0],roots[nselected-1],sf->nroots);
1008f659e5c7SJunchao Zhang 
1009f659e5c7SJunchao Zhang   if (sf->ops->CreateEmbeddedSF) {
1010f659e5c7SJunchao Zhang     ierr = (*sf->ops->CreateEmbeddedSF)(sf,nselected,roots,newsf);CHKERRQ(ierr);
1011f659e5c7SJunchao Zhang   } else {
1012f659e5c7SJunchao Zhang     /* A generic version of creating embedded sf. Note that we called PetscSFSetGraph() twice, which is certainly expensive */
1013ba2a7774SJunchao Zhang     /* Find out which leaves (not leaf data items) are still connected to roots in the embedded sf */
1014f659e5c7SJunchao Zhang     ierr = PetscSFGetGraph(sf,&nroots,&nleaves,NULL,&iremote);CHKERRQ(ierr);
1015ba2a7774SJunchao Zhang     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&tmpsf);CHKERRQ(ierr);
1016ba2a7774SJunchao Zhang     ierr = PetscSFSetGraph(tmpsf,nroots,nleaves,NULL/*contiguous*/,PETSC_USE_POINTER,iremote,PETSC_USE_POINTER);CHKERRQ(ierr);
1017ba2a7774SJunchao Zhang     ierr = PetscCalloc2(nroots,&rootdata,nleaves,&leafdata);CHKERRQ(ierr);
1018f659e5c7SJunchao Zhang     for (i=0; i<nselected; ++i) rootdata[roots[i]] = 1;
1019ba2a7774SJunchao Zhang     ierr = PetscSFBcastBegin(tmpsf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr);
1020ba2a7774SJunchao Zhang     ierr = PetscSFBcastEnd(tmpsf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr);
1021ba2a7774SJunchao Zhang     ierr = PetscSFDestroy(&tmpsf);CHKERRQ(ierr);
1022ba2a7774SJunchao Zhang 
1023ba2a7774SJunchao Zhang     /* Build newsf with leaves that are still connected */
1024f659e5c7SJunchao Zhang     connected_leaves = 0;
1025f659e5c7SJunchao Zhang     for (i=0; i<nleaves; ++i) connected_leaves += leafdata[i];
1026f659e5c7SJunchao Zhang     ierr = PetscMalloc1(connected_leaves,&new_ilocal);CHKERRQ(ierr);
1027f659e5c7SJunchao Zhang     ierr = PetscMalloc1(connected_leaves,&new_iremote);CHKERRQ(ierr);
1028ba2a7774SJunchao Zhang     for (i=0, n=0; i<nleaves; ++i) {
1029ba2a7774SJunchao Zhang       if (leafdata[i]) {
1030ba2a7774SJunchao Zhang         new_ilocal[n]        = sf->mine ? sf->mine[i] : i;
1031ba2a7774SJunchao Zhang         new_iremote[n].rank  = sf->remote[i].rank;
1032ba2a7774SJunchao Zhang         new_iremote[n].index = sf->remote[i].index;
1033fc1ede2bSMatthew G. Knepley         ++n;
103495fce210SBarry Smith       }
103595fce210SBarry Smith     }
1036f659e5c7SJunchao Zhang 
1037f659e5c7SJunchao Zhang     if (n != connected_leaves) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"There is a size mismatch in the SF embedding, %d != %d",n,connected_leaves);
1038f659e5c7SJunchao Zhang     ierr = PetscSFCreate(comm,newsf);CHKERRQ(ierr);
1039f659e5c7SJunchao Zhang     ierr = PetscSFSetFromOptions(*newsf);CHKERRQ(ierr);
1040f659e5c7SJunchao Zhang     ierr = PetscSFSetGraph(*newsf,nroots,connected_leaves,new_ilocal,PETSC_OWN_POINTER,new_iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
104195fce210SBarry Smith     ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr);
1042f659e5c7SJunchao Zhang   }
1043f659e5c7SJunchao Zhang   ierr = PetscFree(roots);CHKERRQ(ierr);
104495fce210SBarry Smith   PetscFunctionReturn(0);
104595fce210SBarry Smith }
104695fce210SBarry Smith 
10472f5fb4c2SMatthew G. Knepley /*@C
10482f5fb4c2SMatthew G. Knepley   PetscSFCreateEmbeddedLeafSF - removes edges from all but the selected leaves, does not remap indices
10492f5fb4c2SMatthew G. Knepley 
10502f5fb4c2SMatthew G. Knepley   Collective
10512f5fb4c2SMatthew G. Knepley 
10522f5fb4c2SMatthew G. Knepley   Input Arguments:
10532f5fb4c2SMatthew G. Knepley + sf - original star forest
1054f659e5c7SJunchao Zhang . nselected  - number of selected leaves on this process
1055f659e5c7SJunchao Zhang - selected   - indices of the selected leaves on this process
10562f5fb4c2SMatthew G. Knepley 
10572f5fb4c2SMatthew G. Knepley   Output Arguments:
10582f5fb4c2SMatthew G. Knepley .  newsf - new star forest
10592f5fb4c2SMatthew G. Knepley 
10602f5fb4c2SMatthew G. Knepley   Level: advanced
10612f5fb4c2SMatthew G. Knepley 
10622f5fb4c2SMatthew G. Knepley .seealso: PetscSFCreateEmbeddedSF(), PetscSFSetGraph(), PetscSFGetGraph()
10632f5fb4c2SMatthew G. Knepley @*/
1064f659e5c7SJunchao Zhang PetscErrorCode PetscSFCreateEmbeddedLeafSF(PetscSF sf,PetscInt nselected,const PetscInt *selected,PetscSF *newsf)
10652f5fb4c2SMatthew G. Knepley {
1066f659e5c7SJunchao Zhang   const PetscSFNode *iremote;
1067f659e5c7SJunchao Zhang   PetscSFNode       *new_iremote;
1068f659e5c7SJunchao Zhang   const PetscInt    *ilocal;
1069f659e5c7SJunchao Zhang   PetscInt          i,nroots,*leaves,*new_ilocal;
1070f659e5c7SJunchao Zhang   MPI_Comm          comm;
10712f5fb4c2SMatthew G. Knepley   PetscErrorCode    ierr;
10722f5fb4c2SMatthew G. Knepley 
10732f5fb4c2SMatthew G. Knepley   PetscFunctionBegin;
10742f5fb4c2SMatthew G. Knepley   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
107529046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
1076f659e5c7SJunchao Zhang   if (nselected) PetscValidPointer(selected,3);
10772f5fb4c2SMatthew G. Knepley   PetscValidPointer(newsf,4);
10782f5fb4c2SMatthew G. Knepley 
1079f659e5c7SJunchao Zhang   /* Uniq selected[] and put results in leaves[] */
1080f659e5c7SJunchao Zhang   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
1081f659e5c7SJunchao Zhang   ierr = PetscMalloc1(nselected,&leaves);CHKERRQ(ierr);
1082f659e5c7SJunchao Zhang   ierr = PetscMemcpy(leaves,selected,sizeof(PetscInt)*nselected);CHKERRQ(ierr);
1083f659e5c7SJunchao Zhang   ierr = PetscSortedRemoveDupsInt(&nselected,leaves);CHKERRQ(ierr);
1084f659e5c7SJunchao 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);
1085f659e5c7SJunchao Zhang 
1086f659e5c7SJunchao Zhang   /* Optimize the routine only when sf is setup and hence we can reuse sf's communication pattern */
1087f659e5c7SJunchao Zhang   if (sf->setupcalled && sf->ops->CreateEmbeddedLeafSF) {
1088f659e5c7SJunchao Zhang     ierr = (*sf->ops->CreateEmbeddedLeafSF)(sf,nselected,leaves,newsf);CHKERRQ(ierr);
1089f659e5c7SJunchao Zhang   } else {
1090f659e5c7SJunchao Zhang     ierr = PetscSFGetGraph(sf,&nroots,NULL,&ilocal,&iremote);CHKERRQ(ierr);
1091f659e5c7SJunchao Zhang     ierr = PetscMalloc1(nselected,&new_ilocal);CHKERRQ(ierr);
1092f659e5c7SJunchao Zhang     ierr = PetscMalloc1(nselected,&new_iremote);CHKERRQ(ierr);
1093f659e5c7SJunchao Zhang     for (i=0; i<nselected; ++i) {
1094f659e5c7SJunchao Zhang       const PetscInt l     = leaves[i];
1095f659e5c7SJunchao Zhang       new_ilocal[i]        = ilocal ? ilocal[l] : l;
1096f659e5c7SJunchao Zhang       new_iremote[i].rank  = iremote[l].rank;
1097f659e5c7SJunchao Zhang       new_iremote[i].index = iremote[l].index;
10982f5fb4c2SMatthew G. Knepley     }
1099f659e5c7SJunchao Zhang     ierr = PetscSFCreate(comm,newsf);CHKERRQ(ierr);
1100f659e5c7SJunchao Zhang     ierr = PetscSFSetFromOptions(*newsf);CHKERRQ(ierr);
1101f659e5c7SJunchao Zhang     ierr = PetscSFSetGraph(*newsf,nroots,nselected,new_ilocal,PETSC_OWN_POINTER,new_iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1102f659e5c7SJunchao Zhang   }
1103f659e5c7SJunchao Zhang   ierr = PetscFree(leaves);CHKERRQ(ierr);
11042f5fb4c2SMatthew G. Knepley   PetscFunctionReturn(0);
11052f5fb4c2SMatthew G. Knepley }
11062f5fb4c2SMatthew G. Knepley 
110795fce210SBarry Smith /*@C
11083482bfa8SJunchao Zhang    PetscSFBcastAndOpBegin - begin pointwise broadcast with root value being reduced to leaf value, to be concluded with call to PetscSFBcastAndOpEnd()
11093482bfa8SJunchao Zhang 
11103482bfa8SJunchao Zhang    Collective on PetscSF
11113482bfa8SJunchao Zhang 
11123482bfa8SJunchao Zhang    Input Arguments:
11133482bfa8SJunchao Zhang +  sf - star forest on which to communicate
11143482bfa8SJunchao Zhang .  unit - data type associated with each node
11153482bfa8SJunchao Zhang .  rootdata - buffer to broadcast
11163482bfa8SJunchao Zhang -  op - operation to use for reduction
11173482bfa8SJunchao Zhang 
11183482bfa8SJunchao Zhang    Output Arguments:
11193482bfa8SJunchao Zhang .  leafdata - buffer to be reduced with values from each leaf's respective root
11203482bfa8SJunchao Zhang 
11213482bfa8SJunchao Zhang    Level: intermediate
11223482bfa8SJunchao Zhang 
11233482bfa8SJunchao Zhang .seealso: PetscSFBcastAndOpEnd(), PetscSFBcastBegin(), PetscSFBcastEnd()
11243482bfa8SJunchao Zhang @*/
11253482bfa8SJunchao Zhang PetscErrorCode PetscSFBcastAndOpBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
11263482bfa8SJunchao Zhang {
11273482bfa8SJunchao Zhang   PetscErrorCode ierr;
11283482bfa8SJunchao Zhang 
11293482bfa8SJunchao Zhang   PetscFunctionBegin;
11303482bfa8SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
11313482bfa8SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
11323482bfa8SJunchao Zhang   ierr = PetscLogEventBegin(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
11333482bfa8SJunchao Zhang   ierr = (*sf->ops->BcastAndOpBegin)(sf,unit,rootdata,leafdata,op);CHKERRQ(ierr);
11343482bfa8SJunchao Zhang   ierr = PetscLogEventEnd(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
11353482bfa8SJunchao Zhang   PetscFunctionReturn(0);
11363482bfa8SJunchao Zhang }
11373482bfa8SJunchao Zhang 
11383482bfa8SJunchao Zhang /*@C
11393482bfa8SJunchao Zhang    PetscSFBcastAndOpEnd - end a broadcast & reduce operation started with PetscSFBcastAndOpBegin()
11403482bfa8SJunchao Zhang 
11413482bfa8SJunchao Zhang    Collective
11423482bfa8SJunchao Zhang 
11433482bfa8SJunchao Zhang    Input Arguments:
11443482bfa8SJunchao Zhang +  sf - star forest
11453482bfa8SJunchao Zhang .  unit - data type
11463482bfa8SJunchao Zhang .  rootdata - buffer to broadcast
11473482bfa8SJunchao Zhang -  op - operation to use for reduction
11483482bfa8SJunchao Zhang 
11493482bfa8SJunchao Zhang    Output Arguments:
11503482bfa8SJunchao Zhang .  leafdata - buffer to be reduced with values from each leaf's respective root
11513482bfa8SJunchao Zhang 
11523482bfa8SJunchao Zhang    Level: intermediate
11533482bfa8SJunchao Zhang 
11543482bfa8SJunchao Zhang .seealso: PetscSFSetGraph(), PetscSFReduceEnd()
11553482bfa8SJunchao Zhang @*/
11563482bfa8SJunchao Zhang PetscErrorCode PetscSFBcastAndOpEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
11573482bfa8SJunchao Zhang {
11583482bfa8SJunchao Zhang   PetscErrorCode ierr;
11593482bfa8SJunchao Zhang 
11603482bfa8SJunchao Zhang   PetscFunctionBegin;
11613482bfa8SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
11623482bfa8SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
11633482bfa8SJunchao Zhang   ierr = PetscLogEventBegin(PETSCSF_BcastAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
11643482bfa8SJunchao Zhang   ierr = (*sf->ops->BcastAndOpEnd)(sf,unit,rootdata,leafdata,op);CHKERRQ(ierr);
11653482bfa8SJunchao Zhang   ierr = PetscLogEventEnd(PETSCSF_BcastAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
11663482bfa8SJunchao Zhang   PetscFunctionReturn(0);
11673482bfa8SJunchao Zhang }
11683482bfa8SJunchao Zhang 
11693482bfa8SJunchao Zhang /*@C
117095fce210SBarry Smith    PetscSFReduceBegin - begin reduction of leafdata into rootdata, to be completed with call to PetscSFReduceEnd()
117195fce210SBarry Smith 
117295fce210SBarry Smith    Collective
117395fce210SBarry Smith 
117495fce210SBarry Smith    Input Arguments:
117595fce210SBarry Smith +  sf - star forest
117695fce210SBarry Smith .  unit - data type
117795fce210SBarry Smith .  leafdata - values to reduce
117895fce210SBarry Smith -  op - reduction operation
117995fce210SBarry Smith 
118095fce210SBarry Smith    Output Arguments:
118195fce210SBarry Smith .  rootdata - result of reduction of values from all leaves of each root
118295fce210SBarry Smith 
118395fce210SBarry Smith    Level: intermediate
118495fce210SBarry Smith 
118595fce210SBarry Smith .seealso: PetscSFBcastBegin()
118695fce210SBarry Smith @*/
118795fce210SBarry Smith PetscErrorCode PetscSFReduceBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
118895fce210SBarry Smith {
118995fce210SBarry Smith   PetscErrorCode ierr;
119095fce210SBarry Smith 
119195fce210SBarry Smith   PetscFunctionBegin;
119295fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
119395fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
119429046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);
119595fce210SBarry Smith   ierr = (sf->ops->ReduceBegin)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
1196563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);
119795fce210SBarry Smith   PetscFunctionReturn(0);
119895fce210SBarry Smith }
119995fce210SBarry Smith 
120095fce210SBarry Smith /*@C
120195fce210SBarry Smith    PetscSFReduceEnd - end a reduction operation started with PetscSFReduceBegin()
120295fce210SBarry Smith 
120395fce210SBarry Smith    Collective
120495fce210SBarry Smith 
120595fce210SBarry Smith    Input Arguments:
120695fce210SBarry Smith +  sf - star forest
120795fce210SBarry Smith .  unit - data type
120895fce210SBarry Smith .  leafdata - values to reduce
120995fce210SBarry Smith -  op - reduction operation
121095fce210SBarry Smith 
121195fce210SBarry Smith    Output Arguments:
121295fce210SBarry Smith .  rootdata - result of reduction of values from all leaves of each root
121395fce210SBarry Smith 
121495fce210SBarry Smith    Level: intermediate
121595fce210SBarry Smith 
121695fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFBcastEnd()
121795fce210SBarry Smith @*/
121895fce210SBarry Smith PetscErrorCode PetscSFReduceEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
121995fce210SBarry Smith {
122095fce210SBarry Smith   PetscErrorCode ierr;
122195fce210SBarry Smith 
122295fce210SBarry Smith   PetscFunctionBegin;
122395fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
122495fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
122529046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);
122695fce210SBarry Smith   ierr = (*sf->ops->ReduceEnd)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
1227563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);
122895fce210SBarry Smith   PetscFunctionReturn(0);
122995fce210SBarry Smith }
123095fce210SBarry Smith 
123195fce210SBarry Smith /*@C
123295fce210SBarry Smith    PetscSFComputeDegreeBegin - begin computation of degree for each root vertex, to be completed with PetscSFComputeDegreeEnd()
123395fce210SBarry Smith 
123495fce210SBarry Smith    Collective
123595fce210SBarry Smith 
123695fce210SBarry Smith    Input Arguments:
123795fce210SBarry Smith .  sf - star forest
123895fce210SBarry Smith 
123995fce210SBarry Smith    Output Arguments:
124095fce210SBarry Smith .  degree - degree of each root vertex
124195fce210SBarry Smith 
124295fce210SBarry Smith    Level: advanced
124395fce210SBarry Smith 
1244ffe67aa5SVáclav Hapla    Notes:
1245ffe67aa5SVáclav Hapla    The returned array is owned by PetscSF and automatically freed by PetscSFDestroy(). Hence no need to call PetscFree() on it.
1246ffe67aa5SVáclav Hapla 
124795fce210SBarry Smith .seealso: PetscSFGatherBegin()
124895fce210SBarry Smith @*/
124995fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf,const PetscInt **degree)
125095fce210SBarry Smith {
125195fce210SBarry Smith   PetscErrorCode ierr;
125295fce210SBarry Smith 
125395fce210SBarry Smith   PetscFunctionBegin;
125495fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
125595fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
125695fce210SBarry Smith   PetscValidPointer(degree,2);
1257803bd9e8SMatthew G. Knepley   if (!sf->degreeknown) {
125829046d53SLisandro Dalcin     PetscInt i, nroots = sf->nroots, maxlocal = sf->maxleaf+1;  /* TODO: We should use PetscSFGetLeafRange() */
1259803bd9e8SMatthew G. Knepley     if (sf->degree) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Calls to PetscSFComputeDegreeBegin() cannot be nested.");
126029046d53SLisandro Dalcin     ierr = PetscMalloc1(nroots,&sf->degree);CHKERRQ(ierr);
126129046d53SLisandro Dalcin     ierr = PetscMalloc1(PetscMax(maxlocal,1),&sf->degreetmp);CHKERRQ(ierr); /* allocate at least one entry, see check in PetscSFComputeDegreeEnd() */
126229046d53SLisandro Dalcin     for (i=0; i<nroots; i++) sf->degree[i] = 0;
12639837ea96SMatthew G. Knepley     for (i=0; i<maxlocal; i++) sf->degreetmp[i] = 1;
1264dbd2ff41SBarry Smith     ierr = PetscSFReduceBegin(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);CHKERRQ(ierr);
126595fce210SBarry Smith   }
126695fce210SBarry Smith   *degree = NULL;
126795fce210SBarry Smith   PetscFunctionReturn(0);
126895fce210SBarry Smith }
126995fce210SBarry Smith 
127095fce210SBarry Smith /*@C
127195fce210SBarry Smith    PetscSFComputeDegreeEnd - complete computation of degree for each root vertex, started with PetscSFComputeDegreeBegin()
127295fce210SBarry Smith 
127395fce210SBarry Smith    Collective
127495fce210SBarry Smith 
127595fce210SBarry Smith    Input Arguments:
127695fce210SBarry Smith .  sf - star forest
127795fce210SBarry Smith 
127895fce210SBarry Smith    Output Arguments:
127995fce210SBarry Smith .  degree - degree of each root vertex
128095fce210SBarry Smith 
128195fce210SBarry Smith    Level: developer
128295fce210SBarry Smith 
1283ffe67aa5SVáclav Hapla    Notes:
1284ffe67aa5SVáclav Hapla    The returned array is owned by PetscSF and automatically freed by PetscSFDestroy(). Hence no need to call PetscFree() on it.
1285ffe67aa5SVáclav Hapla 
128695fce210SBarry Smith .seealso:
128795fce210SBarry Smith @*/
128895fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf,const PetscInt **degree)
128995fce210SBarry Smith {
129095fce210SBarry Smith   PetscErrorCode ierr;
129195fce210SBarry Smith 
129295fce210SBarry Smith   PetscFunctionBegin;
129395fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
129495fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
129529046d53SLisandro Dalcin   PetscValidPointer(degree,2);
129695fce210SBarry Smith   if (!sf->degreeknown) {
129729046d53SLisandro Dalcin     if (!sf->degreetmp) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFComputeDegreeBegin() before PetscSFComputeDegreeEnd()");
1298dbd2ff41SBarry Smith     ierr = PetscSFReduceEnd(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);CHKERRQ(ierr);
129995fce210SBarry Smith     ierr = PetscFree(sf->degreetmp);CHKERRQ(ierr);
130095fce210SBarry Smith     sf->degreeknown = PETSC_TRUE;
130195fce210SBarry Smith   }
130295fce210SBarry Smith   *degree = sf->degree;
130395fce210SBarry Smith   PetscFunctionReturn(0);
130495fce210SBarry Smith }
130595fce210SBarry Smith 
1306673100f5SVaclav Hapla 
1307673100f5SVaclav Hapla /*@C
130866dfcd1aSVaclav Hapla    PetscSFComputeMultiRootOriginalNumbering - Returns original numbering of multi-roots (roots of multi-SF returned by PetscSFGetMultiSF()).
130966dfcd1aSVaclav Hapla    Each multi-root is assigned index of the corresponding original root.
1310673100f5SVaclav Hapla 
1311673100f5SVaclav Hapla    Collective
1312673100f5SVaclav Hapla 
1313673100f5SVaclav Hapla    Input Arguments:
1314673100f5SVaclav Hapla +  sf - star forest
1315673100f5SVaclav Hapla -  degree - degree of each root vertex, computed with PetscSFComputeDegreeBegin()/PetscSFComputeDegreeEnd()
1316673100f5SVaclav Hapla 
1317673100f5SVaclav Hapla    Output Arguments:
131866dfcd1aSVaclav Hapla +  nMultiRoots - (optional) number of multi-roots (roots of multi-SF)
131966dfcd1aSVaclav Hapla -  multiRootsOrigNumbering - original indices of multi-roots; length of this array is nMultiRoots
1320673100f5SVaclav Hapla 
1321673100f5SVaclav Hapla    Level: developer
1322673100f5SVaclav Hapla 
1323ffe67aa5SVáclav Hapla    Notes:
1324ffe67aa5SVáclav Hapla    The returned array multiRootsOrigNumbering is newly allocated and should be destroyed with PetscFree() when no longer needed.
1325ffe67aa5SVáclav Hapla 
1326673100f5SVaclav Hapla .seealso: PetscSFComputeDegreeBegin(), PetscSFComputeDegreeEnd(), PetscSFGetMultiSF()
1327673100f5SVaclav Hapla @*/
132866dfcd1aSVaclav Hapla PetscErrorCode PetscSFComputeMultiRootOriginalNumbering(PetscSF sf, const PetscInt degree[], PetscInt *nMultiRoots, PetscInt *multiRootsOrigNumbering[])
1329673100f5SVaclav Hapla {
1330673100f5SVaclav Hapla   PetscSF             msf;
1331673100f5SVaclav Hapla   PetscInt            i, j, k, nroots, nmroots;
1332673100f5SVaclav Hapla   PetscErrorCode      ierr;
1333673100f5SVaclav Hapla 
1334673100f5SVaclav Hapla   PetscFunctionBegin;
1335673100f5SVaclav Hapla   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1336673100f5SVaclav Hapla   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
133729328920SVaclav Hapla   if (nroots) PetscValidIntPointer(degree,2);
133866dfcd1aSVaclav Hapla   if (nMultiRoots) PetscValidIntPointer(nMultiRoots,3);
133966dfcd1aSVaclav Hapla   PetscValidPointer(multiRootsOrigNumbering,4);
1340673100f5SVaclav Hapla   ierr = PetscSFGetMultiSF(sf,&msf);CHKERRQ(ierr);
1341673100f5SVaclav Hapla   ierr = PetscSFGetGraph(msf, &nmroots, NULL, NULL, NULL);CHKERRQ(ierr);
134266dfcd1aSVaclav Hapla   ierr = PetscMalloc1(nmroots, multiRootsOrigNumbering);CHKERRQ(ierr);
1343673100f5SVaclav Hapla   for (i=0,j=0,k=0; i<nroots; i++) {
1344673100f5SVaclav Hapla     if (!degree[i]) continue;
1345673100f5SVaclav Hapla     for (j=0; j<degree[i]; j++,k++) {
134666dfcd1aSVaclav Hapla       (*multiRootsOrigNumbering)[k] = i;
1347673100f5SVaclav Hapla     }
1348673100f5SVaclav Hapla   }
1349673100f5SVaclav Hapla #if defined(PETSC_USE_DEBUG)
1350673100f5SVaclav Hapla   if (PetscUnlikely(k != nmroots)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sanity check fail");
1351673100f5SVaclav Hapla #endif
135266dfcd1aSVaclav Hapla   if (nMultiRoots) *nMultiRoots = nmroots;
1353673100f5SVaclav Hapla   PetscFunctionReturn(0);
1354673100f5SVaclav Hapla }
1355673100f5SVaclav Hapla 
135695fce210SBarry Smith /*@C
135795fce210SBarry Smith    PetscSFFetchAndOpBegin - begin operation that fetches values from root and updates atomically by applying operation using my leaf value, to be completed with PetscSFFetchAndOpEnd()
135895fce210SBarry Smith 
135995fce210SBarry Smith    Collective
136095fce210SBarry Smith 
136195fce210SBarry Smith    Input Arguments:
136295fce210SBarry Smith +  sf - star forest
136395fce210SBarry Smith .  unit - data type
136495fce210SBarry Smith .  leafdata - leaf values to use in reduction
136595fce210SBarry Smith -  op - operation to use for reduction
136695fce210SBarry Smith 
136795fce210SBarry Smith    Output Arguments:
136895fce210SBarry Smith +  rootdata - root values to be updated, input state is seen by first process to perform an update
136995fce210SBarry Smith -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
137095fce210SBarry Smith 
137195fce210SBarry Smith    Level: advanced
137295fce210SBarry Smith 
137395fce210SBarry Smith    Note:
137495fce210SBarry Smith    The update is only atomic at the granularity provided by the hardware. Different roots referenced by the same process
137595fce210SBarry Smith    might be updated in a different order. Furthermore, if a composite type is used for the unit datatype, atomicity is
137695fce210SBarry Smith    not guaranteed across the whole vertex. Therefore, this function is mostly only used with primitive types such as
137795fce210SBarry Smith    integers.
137895fce210SBarry Smith 
137995fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFReduceBegin(), PetscSFSetGraph()
138095fce210SBarry Smith @*/
138195fce210SBarry Smith PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
138295fce210SBarry Smith {
138395fce210SBarry Smith   PetscErrorCode ierr;
138495fce210SBarry Smith 
138595fce210SBarry Smith   PetscFunctionBegin;
138695fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
138795fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
138829046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
138995fce210SBarry Smith   ierr = (*sf->ops->FetchAndOpBegin)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr);
1390563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
139195fce210SBarry Smith   PetscFunctionReturn(0);
139295fce210SBarry Smith }
139395fce210SBarry Smith 
139495fce210SBarry Smith /*@C
139595fce210SBarry Smith    PetscSFFetchAndOpEnd - end operation started in matching call to PetscSFFetchAndOpBegin() to fetch values from roots and update atomically by applying operation using my leaf value
139695fce210SBarry Smith 
139795fce210SBarry Smith    Collective
139895fce210SBarry Smith 
139995fce210SBarry Smith    Input Arguments:
140095fce210SBarry Smith +  sf - star forest
140195fce210SBarry Smith .  unit - data type
140295fce210SBarry Smith .  leafdata - leaf values to use in reduction
140395fce210SBarry Smith -  op - operation to use for reduction
140495fce210SBarry Smith 
140595fce210SBarry Smith    Output Arguments:
140695fce210SBarry Smith +  rootdata - root values to be updated, input state is seen by first process to perform an update
140795fce210SBarry Smith -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
140895fce210SBarry Smith 
140995fce210SBarry Smith    Level: advanced
141095fce210SBarry Smith 
141195fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFReduceEnd(), PetscSFSetGraph()
141295fce210SBarry Smith @*/
141395fce210SBarry Smith PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
141495fce210SBarry Smith {
141595fce210SBarry Smith   PetscErrorCode ierr;
141695fce210SBarry Smith 
141795fce210SBarry Smith   PetscFunctionBegin;
141895fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
141995fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
142029046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
142195fce210SBarry Smith   ierr = (*sf->ops->FetchAndOpEnd)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr);
1422563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
142395fce210SBarry Smith   PetscFunctionReturn(0);
142495fce210SBarry Smith }
142595fce210SBarry Smith 
142695fce210SBarry Smith /*@C
142795fce210SBarry Smith    PetscSFGatherBegin - begin pointwise gather of all leaves into multi-roots, to be completed with PetscSFGatherEnd()
142895fce210SBarry Smith 
142995fce210SBarry Smith    Collective
143095fce210SBarry Smith 
143195fce210SBarry Smith    Input Arguments:
143295fce210SBarry Smith +  sf - star forest
143395fce210SBarry Smith .  unit - data type
143495fce210SBarry Smith -  leafdata - leaf data to gather to roots
143595fce210SBarry Smith 
143695fce210SBarry Smith    Output Argument:
143795fce210SBarry Smith .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
143895fce210SBarry Smith 
143995fce210SBarry Smith    Level: intermediate
144095fce210SBarry Smith 
144195fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
144295fce210SBarry Smith @*/
144395fce210SBarry Smith PetscErrorCode PetscSFGatherBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
144495fce210SBarry Smith {
144595fce210SBarry Smith   PetscErrorCode ierr;
144695fce210SBarry Smith   PetscSF        multi;
144795fce210SBarry Smith 
144895fce210SBarry Smith   PetscFunctionBegin;
144995fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
145029046d53SLisandro Dalcin   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
145195fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
14528bfbc91cSJed Brown   ierr = PetscSFReduceBegin(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr);
145395fce210SBarry Smith   PetscFunctionReturn(0);
145495fce210SBarry Smith }
145595fce210SBarry Smith 
145695fce210SBarry Smith /*@C
145795fce210SBarry Smith    PetscSFGatherEnd - ends pointwise gather operation that was started with PetscSFGatherBegin()
145895fce210SBarry Smith 
145995fce210SBarry Smith    Collective
146095fce210SBarry Smith 
146195fce210SBarry Smith    Input Arguments:
146295fce210SBarry Smith +  sf - star forest
146395fce210SBarry Smith .  unit - data type
146495fce210SBarry Smith -  leafdata - leaf data to gather to roots
146595fce210SBarry Smith 
146695fce210SBarry Smith    Output Argument:
146795fce210SBarry Smith .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
146895fce210SBarry Smith 
146995fce210SBarry Smith    Level: intermediate
147095fce210SBarry Smith 
147195fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
147295fce210SBarry Smith @*/
147395fce210SBarry Smith PetscErrorCode PetscSFGatherEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
147495fce210SBarry Smith {
147595fce210SBarry Smith   PetscErrorCode ierr;
147695fce210SBarry Smith   PetscSF        multi;
147795fce210SBarry Smith 
147895fce210SBarry Smith   PetscFunctionBegin;
147995fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
148095fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
148195fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
14828bfbc91cSJed Brown   ierr = PetscSFReduceEnd(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr);
148395fce210SBarry Smith   PetscFunctionReturn(0);
148495fce210SBarry Smith }
148595fce210SBarry Smith 
148695fce210SBarry Smith /*@C
148795fce210SBarry Smith    PetscSFScatterBegin - begin pointwise scatter operation from multi-roots to leaves, to be completed with PetscSFScatterEnd()
148895fce210SBarry Smith 
148995fce210SBarry Smith    Collective
149095fce210SBarry Smith 
149195fce210SBarry Smith    Input Arguments:
149295fce210SBarry Smith +  sf - star forest
149395fce210SBarry Smith .  unit - data type
149495fce210SBarry Smith -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
149595fce210SBarry Smith 
149695fce210SBarry Smith    Output Argument:
149795fce210SBarry Smith .  leafdata - leaf data to be update with personal data from each respective root
149895fce210SBarry Smith 
149995fce210SBarry Smith    Level: intermediate
150095fce210SBarry Smith 
150195fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
150295fce210SBarry Smith @*/
150395fce210SBarry Smith PetscErrorCode PetscSFScatterBegin(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
150495fce210SBarry Smith {
150595fce210SBarry Smith   PetscErrorCode ierr;
150695fce210SBarry Smith   PetscSF        multi;
150795fce210SBarry Smith 
150895fce210SBarry Smith   PetscFunctionBegin;
150995fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
151095fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
151195fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
151295fce210SBarry Smith   ierr = PetscSFBcastBegin(multi,unit,multirootdata,leafdata);CHKERRQ(ierr);
151395fce210SBarry Smith   PetscFunctionReturn(0);
151495fce210SBarry Smith }
151595fce210SBarry Smith 
151695fce210SBarry Smith /*@C
151795fce210SBarry Smith    PetscSFScatterEnd - ends pointwise scatter operation that was started with PetscSFScatterBegin()
151895fce210SBarry Smith 
151995fce210SBarry Smith    Collective
152095fce210SBarry Smith 
152195fce210SBarry Smith    Input Arguments:
152295fce210SBarry Smith +  sf - star forest
152395fce210SBarry Smith .  unit - data type
152495fce210SBarry Smith -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
152595fce210SBarry Smith 
152695fce210SBarry Smith    Output Argument:
152795fce210SBarry Smith .  leafdata - leaf data to be update with personal data from each respective root
152895fce210SBarry Smith 
152995fce210SBarry Smith    Level: intermediate
153095fce210SBarry Smith 
153195fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
153295fce210SBarry Smith @*/
153395fce210SBarry Smith PetscErrorCode PetscSFScatterEnd(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
153495fce210SBarry Smith {
153595fce210SBarry Smith   PetscErrorCode ierr;
153695fce210SBarry Smith   PetscSF        multi;
153795fce210SBarry Smith 
153895fce210SBarry Smith   PetscFunctionBegin;
153995fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
154095fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
154195fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
154295fce210SBarry Smith   ierr = PetscSFBcastEnd(multi,unit,multirootdata,leafdata);CHKERRQ(ierr);
154395fce210SBarry Smith   PetscFunctionReturn(0);
154495fce210SBarry Smith }
1545a7b3aa13SAta Mesgarnejad 
1546a7b3aa13SAta Mesgarnejad /*@
1547a7b3aa13SAta Mesgarnejad   PetscSFCompose - Compose a new PetscSF equivalent to action to PetscSFs
1548a7b3aa13SAta Mesgarnejad 
1549a7b3aa13SAta Mesgarnejad   Input Parameters:
1550a7b3aa13SAta Mesgarnejad + sfA - The first PetscSF
1551a7b3aa13SAta Mesgarnejad - sfB - The second PetscSF
1552a7b3aa13SAta Mesgarnejad 
1553a7b3aa13SAta Mesgarnejad   Output Parameters:
1554a7b3aa13SAta Mesgarnejad . sfBA - equvalent PetscSF for applying A then B
1555a7b3aa13SAta Mesgarnejad 
1556a7b3aa13SAta Mesgarnejad   Level: developer
1557a7b3aa13SAta Mesgarnejad 
1558a7b3aa13SAta Mesgarnejad .seealso: PetscSF, PetscSFGetGraph(), PetscSFSetGraph()
1559a7b3aa13SAta Mesgarnejad @*/
1560a7b3aa13SAta Mesgarnejad PetscErrorCode PetscSFCompose(PetscSF sfA, PetscSF sfB, PetscSF *sfBA)
1561a7b3aa13SAta Mesgarnejad {
1562a7b3aa13SAta Mesgarnejad   MPI_Comm           comm;
1563a7b3aa13SAta Mesgarnejad   const PetscSFNode *remotePointsA, *remotePointsB;
1564a7b3aa13SAta Mesgarnejad   PetscSFNode       *remotePointsBA;
1565a7b3aa13SAta Mesgarnejad   const PetscInt    *localPointsA, *localPointsB;
1566a7b3aa13SAta Mesgarnejad   PetscInt           numRootsA, numLeavesA, numRootsB, numLeavesB;
1567a7b3aa13SAta Mesgarnejad   PetscErrorCode     ierr;
1568a7b3aa13SAta Mesgarnejad 
1569a7b3aa13SAta Mesgarnejad   PetscFunctionBegin;
1570a7b3aa13SAta Mesgarnejad   PetscValidHeaderSpecific(sfA, PETSCSF_CLASSID, 1);
157129046d53SLisandro Dalcin   PetscSFCheckGraphSet(sfA, 1);
157229046d53SLisandro Dalcin   PetscValidHeaderSpecific(sfB, PETSCSF_CLASSID, 2);
157329046d53SLisandro Dalcin   PetscSFCheckGraphSet(sfB, 2);
157429046d53SLisandro Dalcin   PetscValidPointer(sfBA, 3);
1575a7b3aa13SAta Mesgarnejad   ierr = PetscObjectGetComm((PetscObject) sfA, &comm);CHKERRQ(ierr);
1576a7b3aa13SAta Mesgarnejad   ierr = PetscSFGetGraph(sfA, &numRootsA, &numLeavesA, &localPointsA, &remotePointsA);CHKERRQ(ierr);
1577a7b3aa13SAta Mesgarnejad   ierr = PetscSFGetGraph(sfB, &numRootsB, &numLeavesB, &localPointsB, &remotePointsB);CHKERRQ(ierr);
1578a7b3aa13SAta Mesgarnejad   ierr = PetscMalloc1(numLeavesB, &remotePointsBA);CHKERRQ(ierr);
1579a7b3aa13SAta Mesgarnejad   ierr = PetscSFBcastBegin(sfB, MPIU_2INT, remotePointsA, remotePointsBA);CHKERRQ(ierr);
1580a7b3aa13SAta Mesgarnejad   ierr = PetscSFBcastEnd(sfB, MPIU_2INT, remotePointsA, remotePointsBA);CHKERRQ(ierr);
1581a7b3aa13SAta Mesgarnejad   ierr = PetscSFCreate(comm, sfBA);CHKERRQ(ierr);
1582a7b3aa13SAta Mesgarnejad   ierr = PetscSFSetGraph(*sfBA, numRootsA, numLeavesB, localPointsB, PETSC_COPY_VALUES, remotePointsBA, PETSC_OWN_POINTER);CHKERRQ(ierr);
1583a7b3aa13SAta Mesgarnejad   PetscFunctionReturn(0);
1584a7b3aa13SAta Mesgarnejad }
1585*1c6ba672SJunchao Zhang 
1586*1c6ba672SJunchao Zhang /*
1587*1c6ba672SJunchao Zhang   PetscSFCreateLocalSF_Private - Creates a local PetscSF that only has intra-process edges of the global PetscSF
1588*1c6ba672SJunchao Zhang 
1589*1c6ba672SJunchao Zhang   Input Parameters:
1590*1c6ba672SJunchao Zhang . sf - The global PetscSF
1591*1c6ba672SJunchao Zhang 
1592*1c6ba672SJunchao Zhang   Output Parameters:
1593*1c6ba672SJunchao Zhang . out - The local PetscSF
1594*1c6ba672SJunchao Zhang  */
1595*1c6ba672SJunchao Zhang PetscErrorCode PetscSFCreateLocalSF_Private(PetscSF sf,PetscSF *out)
1596*1c6ba672SJunchao Zhang {
1597*1c6ba672SJunchao Zhang   MPI_Comm           comm;
1598*1c6ba672SJunchao Zhang   PetscMPIInt        myrank;
1599*1c6ba672SJunchao Zhang   const PetscInt     *ilocal;
1600*1c6ba672SJunchao Zhang   const PetscSFNode  *iremote;
1601*1c6ba672SJunchao Zhang   PetscInt           i,j,nroots,nleaves,lnleaves,*lilocal;
1602*1c6ba672SJunchao Zhang   PetscSFNode        *liremote;
1603*1c6ba672SJunchao Zhang   PetscSF            lsf;
1604*1c6ba672SJunchao Zhang   PetscErrorCode     ierr;
1605*1c6ba672SJunchao Zhang 
1606*1c6ba672SJunchao Zhang   PetscFunctionBegin;
1607*1c6ba672SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1608*1c6ba672SJunchao Zhang   if (sf->ops->CreateLocalSF) {
1609*1c6ba672SJunchao Zhang     ierr = (*sf->ops->CreateLocalSF)(sf,out);CHKERRQ(ierr);
1610*1c6ba672SJunchao Zhang   } else {
1611*1c6ba672SJunchao Zhang     /* Could use PetscSFCreateEmbeddedLeafSF, but since we know the comm is PETSC_COMM_SELF, we can make it fast */
1612*1c6ba672SJunchao Zhang     ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
1613*1c6ba672SJunchao Zhang     ierr = MPI_Comm_rank(comm,&myrank);CHKERRQ(ierr);
1614*1c6ba672SJunchao Zhang 
1615*1c6ba672SJunchao Zhang     /* Find out local edges and build a local SF */
1616*1c6ba672SJunchao Zhang     ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
1617*1c6ba672SJunchao Zhang     for (i=lnleaves=0; i<nleaves; i++) {if (iremote[i].rank == (PetscInt)myrank) lnleaves++;}
1618*1c6ba672SJunchao Zhang     ierr = PetscMalloc1(lnleaves,&lilocal);CHKERRQ(ierr);
1619*1c6ba672SJunchao Zhang     ierr = PetscMalloc1(lnleaves,&liremote);CHKERRQ(ierr);
1620*1c6ba672SJunchao Zhang 
1621*1c6ba672SJunchao Zhang     for (i=j=0; i<nleaves; i++) {
1622*1c6ba672SJunchao Zhang       if (iremote[i].rank == (PetscInt)myrank) {
1623*1c6ba672SJunchao Zhang         lilocal[j]        = ilocal? ilocal[i] : i; /* ilocal=NULL for contiguous storage */
1624*1c6ba672SJunchao Zhang         liremote[j].rank  = 0; /* rank in PETSC_COMM_SELF */
1625*1c6ba672SJunchao Zhang         liremote[j].index = iremote[i].index;
1626*1c6ba672SJunchao Zhang         j++;
1627*1c6ba672SJunchao Zhang       }
1628*1c6ba672SJunchao Zhang     }
1629*1c6ba672SJunchao Zhang     ierr = PetscSFCreate(PETSC_COMM_SELF,&lsf);CHKERRQ(ierr);
1630*1c6ba672SJunchao Zhang     ierr = PetscSFSetGraph(lsf,nroots,lnleaves,lilocal,PETSC_OWN_POINTER,liremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1631*1c6ba672SJunchao Zhang     ierr = PetscSFSetUp(lsf);CHKERRQ(ierr);
1632*1c6ba672SJunchao Zhang     *out = lsf;
1633*1c6ba672SJunchao Zhang   }
1634*1c6ba672SJunchao Zhang   PetscFunctionReturn(0);
1635*1c6ba672SJunchao Zhang }
1636