xref: /petsc/src/vec/is/sf/interface/sf.c (revision d083f849a86f1f43e18d534ee43954e2786cb29a)
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 
18*d083f849SBarry 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);
34429046d53SLisandro Dalcin       ierr = PetscMemcpy(sf->mine_alloc,ilocal,nleaves*sizeof(*ilocal));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);
36229046d53SLisandro Dalcin     ierr = PetscMemcpy(sf->remote_alloc,iremote,nleaves*sizeof(*iremote));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);
60081bfa7aaSJed Brown       ierr = PetscMemcpy(tmpranks,sf->ranks,sf->nranks*sizeof(tmpranks[0]));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
62195fce210SBarry Smith    PetscSFGetRanks - Get 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 
63795fce210SBarry Smith .seealso: PetscSFSetGraph()
63895fce210SBarry Smith @*/
63995fce210SBarry Smith PetscErrorCode PetscSFGetRanks(PetscSF sf,PetscInt *nranks,const PetscMPIInt **ranks,const PetscInt **roffset,const PetscInt **rmine,const PetscInt **rremote)
64095fce210SBarry Smith {
64195fce210SBarry Smith 
64295fce210SBarry Smith   PetscFunctionBegin;
64395fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
64429046d53SLisandro Dalcin   if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining ranks");
64595fce210SBarry Smith   if (nranks)  *nranks  = sf->nranks;
64695fce210SBarry Smith   if (ranks)   *ranks   = sf->ranks;
64795fce210SBarry Smith   if (roffset) *roffset = sf->roffset;
64895fce210SBarry Smith   if (rmine)   *rmine   = sf->rmine;
64995fce210SBarry Smith   if (rremote) *rremote = sf->rremote;
65095fce210SBarry Smith   PetscFunctionReturn(0);
65195fce210SBarry Smith }
65295fce210SBarry Smith 
6538750ddebSJunchao Zhang /*@C
6548750ddebSJunchao Zhang    PetscSFGetLeafRanks - Get leaf ranks referencing roots on this process
6558750ddebSJunchao Zhang 
6568750ddebSJunchao Zhang    Not Collective
6578750ddebSJunchao Zhang 
6588750ddebSJunchao Zhang    Input Arguments:
6598750ddebSJunchao Zhang .  sf - star forest
6608750ddebSJunchao Zhang 
6618750ddebSJunchao Zhang    Output Arguments:
6628750ddebSJunchao Zhang +  niranks - number of leaf ranks referencing roots on this process
6638750ddebSJunchao Zhang .  iranks - array of ranks
6648750ddebSJunchao Zhang .  ioffset - offset in irootloc for each rank (length niranks+1)
6658750ddebSJunchao Zhang -  irootloc - concatenated array holding local indices of roots referenced by each leaf rank
6668750ddebSJunchao Zhang 
6678750ddebSJunchao Zhang    Level: developer
6688750ddebSJunchao Zhang 
6698750ddebSJunchao Zhang .seealso: PetscSFGetRanks()
6708750ddebSJunchao Zhang @*/
6718750ddebSJunchao Zhang PetscErrorCode PetscSFGetLeafRanks(PetscSF sf,PetscInt *niranks,const PetscMPIInt **iranks,const PetscInt **ioffset,const PetscInt **irootloc)
6728750ddebSJunchao Zhang {
6738750ddebSJunchao Zhang   PetscErrorCode ierr;
6748750ddebSJunchao Zhang 
6758750ddebSJunchao Zhang   PetscFunctionBegin;
6768750ddebSJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
6778750ddebSJunchao Zhang   if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining ranks");
6788750ddebSJunchao Zhang   if (sf->ops->GetLeafRanks) {
6798750ddebSJunchao Zhang     ierr = (sf->ops->GetLeafRanks)(sf,niranks,iranks,ioffset,irootloc);CHKERRQ(ierr);
6808750ddebSJunchao Zhang   } else {
6818750ddebSJunchao Zhang     PetscSFType type;
6828750ddebSJunchao Zhang     ierr = PetscSFGetType(sf,&type);CHKERRQ(ierr);
6838750ddebSJunchao Zhang     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSFGetLeafRanks() is not supported on this StarForest type: %s", type);
6848750ddebSJunchao Zhang   }
6858750ddebSJunchao Zhang   PetscFunctionReturn(0);
6868750ddebSJunchao Zhang }
6878750ddebSJunchao Zhang 
688b5a8e515SJed Brown static PetscBool InList(PetscMPIInt needle,PetscMPIInt n,const PetscMPIInt *list) {
689b5a8e515SJed Brown   PetscInt i;
690b5a8e515SJed Brown   for (i=0; i<n; i++) {
691b5a8e515SJed Brown     if (needle == list[i]) return PETSC_TRUE;
692b5a8e515SJed Brown   }
693b5a8e515SJed Brown   return PETSC_FALSE;
694b5a8e515SJed Brown }
695b5a8e515SJed Brown 
69695fce210SBarry Smith /*@C
69721c688dcSJed Brown    PetscSFSetUpRanks - Set up data structures associated with ranks; this is for internal use by PetscSF implementations.
69821c688dcSJed Brown 
69921c688dcSJed Brown    Collective
70021c688dcSJed Brown 
70121c688dcSJed Brown    Input Arguments:
702b5a8e515SJed Brown +  sf - PetscSF to set up; PetscSFSetGraph() must have been called
703b5a8e515SJed Brown -  dgroup - MPI_Group of ranks to be distinguished (e.g., for self or shared memory exchange)
70421c688dcSJed Brown 
70521c688dcSJed Brown    Level: developer
70621c688dcSJed Brown 
70721c688dcSJed Brown .seealso: PetscSFGetRanks()
70821c688dcSJed Brown @*/
709b5a8e515SJed Brown PetscErrorCode PetscSFSetUpRanks(PetscSF sf,MPI_Group dgroup)
71021c688dcSJed Brown {
71121c688dcSJed Brown   PetscErrorCode     ierr;
71221c688dcSJed Brown   PetscTable         table;
71321c688dcSJed Brown   PetscTablePosition pos;
714b5a8e515SJed Brown   PetscMPIInt        size,groupsize,*groupranks;
715247e8311SStefano Zampini   PetscInt           *rcount,*ranks;
716247e8311SStefano Zampini   PetscInt           i, irank = -1,orank = -1;
71721c688dcSJed Brown 
71821c688dcSJed Brown   PetscFunctionBegin;
71921c688dcSJed Brown   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
72029046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
72121c688dcSJed Brown   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr);
72221c688dcSJed Brown   ierr = PetscTableCreate(10,size,&table);CHKERRQ(ierr);
72321c688dcSJed Brown   for (i=0; i<sf->nleaves; i++) {
72421c688dcSJed Brown     /* Log 1-based rank */
72521c688dcSJed Brown     ierr = PetscTableAdd(table,sf->remote[i].rank+1,1,ADD_VALUES);CHKERRQ(ierr);
72621c688dcSJed Brown   }
72721c688dcSJed Brown   ierr = PetscTableGetCount(table,&sf->nranks);CHKERRQ(ierr);
72821c688dcSJed Brown   ierr = PetscMalloc4(sf->nranks,&sf->ranks,sf->nranks+1,&sf->roffset,sf->nleaves,&sf->rmine,sf->nleaves,&sf->rremote);CHKERRQ(ierr);
72921c688dcSJed Brown   ierr = PetscMalloc2(sf->nranks,&rcount,sf->nranks,&ranks);CHKERRQ(ierr);
73021c688dcSJed Brown   ierr = PetscTableGetHeadPosition(table,&pos);CHKERRQ(ierr);
73121c688dcSJed Brown   for (i=0; i<sf->nranks; i++) {
73221c688dcSJed Brown     ierr = PetscTableGetNext(table,&pos,&ranks[i],&rcount[i]);CHKERRQ(ierr);
73321c688dcSJed Brown     ranks[i]--;             /* Convert back to 0-based */
73421c688dcSJed Brown   }
73521c688dcSJed Brown   ierr = PetscTableDestroy(&table);CHKERRQ(ierr);
736b5a8e515SJed Brown 
737b5a8e515SJed Brown   /* We expect that dgroup is reliably "small" while nranks could be large */
738b5a8e515SJed Brown   {
7397fb8a5e4SKarl Rupp     MPI_Group group = MPI_GROUP_NULL;
740b5a8e515SJed Brown     PetscMPIInt *dgroupranks;
741b5a8e515SJed Brown     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
742b5a8e515SJed Brown     ierr = MPI_Group_size(dgroup,&groupsize);CHKERRQ(ierr);
743b5a8e515SJed Brown     ierr = PetscMalloc1(groupsize,&dgroupranks);CHKERRQ(ierr);
744b5a8e515SJed Brown     ierr = PetscMalloc1(groupsize,&groupranks);CHKERRQ(ierr);
745b5a8e515SJed Brown     for (i=0; i<groupsize; i++) dgroupranks[i] = i;
746f3fc9a17SJed Brown     if (groupsize) {ierr = MPI_Group_translate_ranks(dgroup,groupsize,dgroupranks,group,groupranks);CHKERRQ(ierr);}
747b5a8e515SJed Brown     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
748b5a8e515SJed Brown     ierr = PetscFree(dgroupranks);CHKERRQ(ierr);
749b5a8e515SJed Brown   }
750b5a8e515SJed Brown 
751b5a8e515SJed Brown   /* Partition ranks[] into distinguished (first sf->ndranks) followed by non-distinguished */
752b5a8e515SJed Brown   for (sf->ndranks=0,i=sf->nranks; sf->ndranks<i; ) {
753b5a8e515SJed Brown     for (i--; sf->ndranks<i; i--) { /* Scan i backward looking for distinguished rank */
754b5a8e515SJed Brown       if (InList(ranks[i],groupsize,groupranks)) break;
755b5a8e515SJed Brown     }
756b5a8e515SJed Brown     for ( ; sf->ndranks<=i; sf->ndranks++) { /* Scan sf->ndranks forward looking for non-distinguished rank */
757b5a8e515SJed Brown       if (!InList(ranks[sf->ndranks],groupsize,groupranks)) break;
758b5a8e515SJed Brown     }
759b5a8e515SJed Brown     if (sf->ndranks < i) {                         /* Swap ranks[sf->ndranks] with ranks[i] */
760b5a8e515SJed Brown       PetscInt    tmprank,tmpcount;
761247e8311SStefano Zampini 
762b5a8e515SJed Brown       tmprank             = ranks[i];
763b5a8e515SJed Brown       tmpcount            = rcount[i];
764b5a8e515SJed Brown       ranks[i]            = ranks[sf->ndranks];
765b5a8e515SJed Brown       rcount[i]           = rcount[sf->ndranks];
766b5a8e515SJed Brown       ranks[sf->ndranks]  = tmprank;
767b5a8e515SJed Brown       rcount[sf->ndranks] = tmpcount;
768b5a8e515SJed Brown       sf->ndranks++;
769b5a8e515SJed Brown     }
770b5a8e515SJed Brown   }
771b5a8e515SJed Brown   ierr = PetscFree(groupranks);CHKERRQ(ierr);
772b5a8e515SJed Brown   ierr = PetscSortIntWithArray(sf->ndranks,ranks,rcount);CHKERRQ(ierr);
773b5a8e515SJed Brown   ierr = PetscSortIntWithArray(sf->nranks-sf->ndranks,ranks+sf->ndranks,rcount+sf->ndranks);CHKERRQ(ierr);
77421c688dcSJed Brown   sf->roffset[0] = 0;
77521c688dcSJed Brown   for (i=0; i<sf->nranks; i++) {
77621c688dcSJed Brown     ierr = PetscMPIIntCast(ranks[i],sf->ranks+i);CHKERRQ(ierr);
77721c688dcSJed Brown     sf->roffset[i+1] = sf->roffset[i] + rcount[i];
77821c688dcSJed Brown     rcount[i]        = 0;
77921c688dcSJed Brown   }
780247e8311SStefano Zampini   for (i=0, irank = -1, orank = -1; i<sf->nleaves; i++) {
781247e8311SStefano Zampini     /* short circuit */
782247e8311SStefano Zampini     if (orank != sf->remote[i].rank) {
78321c688dcSJed Brown       /* Search for index of iremote[i].rank in sf->ranks */
784b5a8e515SJed Brown       ierr = PetscFindMPIInt(sf->remote[i].rank,sf->ndranks,sf->ranks,&irank);CHKERRQ(ierr);
785b5a8e515SJed Brown       if (irank < 0) {
786b5a8e515SJed Brown         ierr = PetscFindMPIInt(sf->remote[i].rank,sf->nranks-sf->ndranks,sf->ranks+sf->ndranks,&irank);CHKERRQ(ierr);
787b5a8e515SJed Brown         if (irank >= 0) irank += sf->ndranks;
78821c688dcSJed Brown       }
789247e8311SStefano Zampini       orank = sf->remote[i].rank;
790247e8311SStefano Zampini     }
791b5a8e515SJed Brown     if (irank < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Could not find rank %D in array",sf->remote[i].rank);
79221c688dcSJed Brown     sf->rmine[sf->roffset[irank] + rcount[irank]]   = sf->mine ? sf->mine[i] : i;
79321c688dcSJed Brown     sf->rremote[sf->roffset[irank] + rcount[irank]] = sf->remote[i].index;
79421c688dcSJed Brown     rcount[irank]++;
79521c688dcSJed Brown   }
79621c688dcSJed Brown   ierr = PetscFree2(rcount,ranks);CHKERRQ(ierr);
79721c688dcSJed Brown   PetscFunctionReturn(0);
79821c688dcSJed Brown }
79921c688dcSJed Brown 
80021c688dcSJed Brown /*@C
80195fce210SBarry Smith    PetscSFGetGroups - gets incoming and outgoing process groups
80295fce210SBarry Smith 
80395fce210SBarry Smith    Collective
80495fce210SBarry Smith 
80595fce210SBarry Smith    Input Argument:
80695fce210SBarry Smith .  sf - star forest
80795fce210SBarry Smith 
80895fce210SBarry Smith    Output Arguments:
80995fce210SBarry Smith +  incoming - group of origin processes for incoming edges (leaves that reference my roots)
81095fce210SBarry Smith -  outgoing - group of destination processes for outgoing edges (roots that I reference)
81195fce210SBarry Smith 
81295fce210SBarry Smith    Level: developer
81395fce210SBarry Smith 
81495fce210SBarry Smith .seealso: PetscSFGetWindow(), PetscSFRestoreWindow()
81595fce210SBarry Smith @*/
81695fce210SBarry Smith PetscErrorCode PetscSFGetGroups(PetscSF sf,MPI_Group *incoming,MPI_Group *outgoing)
81795fce210SBarry Smith {
81895fce210SBarry Smith   PetscErrorCode ierr;
8197fb8a5e4SKarl Rupp   MPI_Group      group = MPI_GROUP_NULL;
82095fce210SBarry Smith 
82195fce210SBarry Smith   PetscFunctionBegin;
82229046d53SLisandro Dalcin   if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining groups");
82395fce210SBarry Smith   if (sf->ingroup == MPI_GROUP_NULL) {
82495fce210SBarry Smith     PetscInt       i;
82595fce210SBarry Smith     const PetscInt *indegree;
82695fce210SBarry Smith     PetscMPIInt    rank,*outranks,*inranks;
82795fce210SBarry Smith     PetscSFNode    *remote;
82895fce210SBarry Smith     PetscSF        bgcount;
82995fce210SBarry Smith 
83095fce210SBarry Smith     /* Compute the number of incoming ranks */
831785e854fSJed Brown     ierr = PetscMalloc1(sf->nranks,&remote);CHKERRQ(ierr);
83295fce210SBarry Smith     for (i=0; i<sf->nranks; i++) {
83395fce210SBarry Smith       remote[i].rank  = sf->ranks[i];
83495fce210SBarry Smith       remote[i].index = 0;
83595fce210SBarry Smith     }
83695fce210SBarry Smith     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,&bgcount);CHKERRQ(ierr);
83795fce210SBarry Smith     ierr = PetscSFSetGraph(bgcount,1,sf->nranks,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
83895fce210SBarry Smith     ierr = PetscSFComputeDegreeBegin(bgcount,&indegree);CHKERRQ(ierr);
83995fce210SBarry Smith     ierr = PetscSFComputeDegreeEnd(bgcount,&indegree);CHKERRQ(ierr);
84095fce210SBarry Smith 
84195fce210SBarry Smith     /* Enumerate the incoming ranks */
842dcca6d9dSJed Brown     ierr = PetscMalloc2(indegree[0],&inranks,sf->nranks,&outranks);CHKERRQ(ierr);
84395fce210SBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
84495fce210SBarry Smith     for (i=0; i<sf->nranks; i++) outranks[i] = rank;
84595fce210SBarry Smith     ierr = PetscSFGatherBegin(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr);
84695fce210SBarry Smith     ierr = PetscSFGatherEnd(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr);
84795fce210SBarry Smith     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
84895fce210SBarry Smith     ierr = MPI_Group_incl(group,indegree[0],inranks,&sf->ingroup);CHKERRQ(ierr);
84995fce210SBarry Smith     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
85095fce210SBarry Smith     ierr = PetscFree2(inranks,outranks);CHKERRQ(ierr);
85195fce210SBarry Smith     ierr = PetscSFDestroy(&bgcount);CHKERRQ(ierr);
85295fce210SBarry Smith   }
85395fce210SBarry Smith   *incoming = sf->ingroup;
85495fce210SBarry Smith 
85595fce210SBarry Smith   if (sf->outgroup == MPI_GROUP_NULL) {
85695fce210SBarry Smith     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
85795fce210SBarry Smith     ierr = MPI_Group_incl(group,sf->nranks,sf->ranks,&sf->outgroup);CHKERRQ(ierr);
85895fce210SBarry Smith     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
85995fce210SBarry Smith   }
86095fce210SBarry Smith   *outgoing = sf->outgroup;
86195fce210SBarry Smith   PetscFunctionReturn(0);
86295fce210SBarry Smith }
86395fce210SBarry Smith 
86429046d53SLisandro Dalcin /*@
86595fce210SBarry Smith    PetscSFGetMultiSF - gets the inner SF implemeting gathers and scatters
86695fce210SBarry Smith 
86795fce210SBarry Smith    Collective
86895fce210SBarry Smith 
86995fce210SBarry Smith    Input Argument:
87095fce210SBarry Smith .  sf - star forest that may contain roots with 0 or with more than 1 vertex
87195fce210SBarry Smith 
87295fce210SBarry Smith    Output Arguments:
87395fce210SBarry Smith .  multi - star forest with split roots, such that each root has degree exactly 1
87495fce210SBarry Smith 
87595fce210SBarry Smith    Level: developer
87695fce210SBarry Smith 
87795fce210SBarry Smith    Notes:
87895fce210SBarry Smith 
87995fce210SBarry Smith    In most cases, users should use PetscSFGatherBegin() and PetscSFScatterBegin() instead of manipulating multi
88095fce210SBarry Smith    directly. Since multi satisfies the stronger condition that each entry in the global space has exactly one incoming
88195fce210SBarry Smith    edge, it is a candidate for future optimization that might involve its removal.
88295fce210SBarry Smith 
883673100f5SVaclav Hapla .seealso: PetscSFSetGraph(), PetscSFGatherBegin(), PetscSFScatterBegin(), PetscSFComputeMultiRootOriginalNumbering()
88495fce210SBarry Smith @*/
88595fce210SBarry Smith PetscErrorCode PetscSFGetMultiSF(PetscSF sf,PetscSF *multi)
88695fce210SBarry Smith {
88795fce210SBarry Smith   PetscErrorCode ierr;
88895fce210SBarry Smith 
88995fce210SBarry Smith   PetscFunctionBegin;
89095fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
89195fce210SBarry Smith   PetscValidPointer(multi,2);
89295fce210SBarry Smith   if (sf->nroots < 0) {         /* Graph has not been set yet; why do we need this? */
89395fce210SBarry Smith     ierr   = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr);
89495fce210SBarry Smith     *multi = sf->multi;
89595fce210SBarry Smith     PetscFunctionReturn(0);
89695fce210SBarry Smith   }
89795fce210SBarry Smith   if (!sf->multi) {
89895fce210SBarry Smith     const PetscInt *indegree;
8999837ea96SMatthew G. Knepley     PetscInt       i,*inoffset,*outones,*outoffset,maxlocal;
90095fce210SBarry Smith     PetscSFNode    *remote;
90129046d53SLisandro Dalcin     maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */
90295fce210SBarry Smith     ierr = PetscSFComputeDegreeBegin(sf,&indegree);CHKERRQ(ierr);
90395fce210SBarry Smith     ierr = PetscSFComputeDegreeEnd(sf,&indegree);CHKERRQ(ierr);
9049837ea96SMatthew G. Knepley     ierr = PetscMalloc3(sf->nroots+1,&inoffset,maxlocal,&outones,maxlocal,&outoffset);CHKERRQ(ierr);
90595fce210SBarry Smith     inoffset[0] = 0;
90695fce210SBarry Smith     for (i=0; i<sf->nroots; i++) inoffset[i+1] = inoffset[i] + indegree[i];
9079837ea96SMatthew G. Knepley     for (i=0; i<maxlocal; i++) outones[i] = 1;
908dbd2ff41SBarry Smith     ierr = PetscSFFetchAndOpBegin(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr);
909dbd2ff41SBarry Smith     ierr = PetscSFFetchAndOpEnd(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr);
91095fce210SBarry Smith     for (i=0; i<sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */
91195fce210SBarry Smith #if 0
91295fce210SBarry Smith #if defined(PETSC_USE_DEBUG)                                 /* Check that the expected number of increments occurred */
91395fce210SBarry Smith     for (i=0; i<sf->nroots; i++) {
91495fce210SBarry Smith       if (inoffset[i] + indegree[i] != inoffset[i+1]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect result after PetscSFFetchAndOp");
91595fce210SBarry Smith     }
91695fce210SBarry Smith #endif
91795fce210SBarry Smith #endif
918785e854fSJed Brown     ierr = PetscMalloc1(sf->nleaves,&remote);CHKERRQ(ierr);
91995fce210SBarry Smith     for (i=0; i<sf->nleaves; i++) {
92095fce210SBarry Smith       remote[i].rank  = sf->remote[i].rank;
92138e7336fSToby Isaac       remote[i].index = outoffset[sf->mine ? sf->mine[i] : i];
92295fce210SBarry Smith     }
92395fce210SBarry Smith     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr);
92401365b40SToby Isaac     ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
92595fce210SBarry Smith     if (sf->rankorder) {        /* Sort the ranks */
92695fce210SBarry Smith       PetscMPIInt rank;
92795fce210SBarry Smith       PetscInt    *inranks,*newoffset,*outranks,*newoutoffset,*tmpoffset,maxdegree;
92895fce210SBarry Smith       PetscSFNode *newremote;
92995fce210SBarry Smith       ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
93095fce210SBarry Smith       for (i=0,maxdegree=0; i<sf->nroots; i++) maxdegree = PetscMax(maxdegree,indegree[i]);
9319837ea96SMatthew G. Knepley       ierr = PetscMalloc5(sf->multi->nroots,&inranks,sf->multi->nroots,&newoffset,maxlocal,&outranks,maxlocal,&newoutoffset,maxdegree,&tmpoffset);CHKERRQ(ierr);
9329837ea96SMatthew G. Knepley       for (i=0; i<maxlocal; i++) outranks[i] = rank;
9338bfbc91cSJed Brown       ierr = PetscSFReduceBegin(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr);
9348bfbc91cSJed Brown       ierr = PetscSFReduceEnd(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr);
93595fce210SBarry Smith       /* Sort the incoming ranks at each vertex, build the inverse map */
93695fce210SBarry Smith       for (i=0; i<sf->nroots; i++) {
93795fce210SBarry Smith         PetscInt j;
93895fce210SBarry Smith         for (j=0; j<indegree[i]; j++) tmpoffset[j] = j;
93995fce210SBarry Smith         ierr = PetscSortIntWithArray(indegree[i],inranks+inoffset[i],tmpoffset);CHKERRQ(ierr);
94095fce210SBarry Smith         for (j=0; j<indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j;
94195fce210SBarry Smith       }
94295fce210SBarry Smith       ierr = PetscSFBcastBegin(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr);
94395fce210SBarry Smith       ierr = PetscSFBcastEnd(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr);
944785e854fSJed Brown       ierr = PetscMalloc1(sf->nleaves,&newremote);CHKERRQ(ierr);
94595fce210SBarry Smith       for (i=0; i<sf->nleaves; i++) {
94695fce210SBarry Smith         newremote[i].rank  = sf->remote[i].rank;
94701365b40SToby Isaac         newremote[i].index = newoutoffset[sf->mine ? sf->mine[i] : i];
94895fce210SBarry Smith       }
94901365b40SToby Isaac       ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,newremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
95095fce210SBarry Smith       ierr = PetscFree5(inranks,newoffset,outranks,newoutoffset,tmpoffset);CHKERRQ(ierr);
95195fce210SBarry Smith     }
95295fce210SBarry Smith     ierr = PetscFree3(inoffset,outones,outoffset);CHKERRQ(ierr);
95395fce210SBarry Smith   }
95495fce210SBarry Smith   *multi = sf->multi;
95595fce210SBarry Smith   PetscFunctionReturn(0);
95695fce210SBarry Smith }
95795fce210SBarry Smith 
95895fce210SBarry Smith /*@C
95995fce210SBarry Smith    PetscSFCreateEmbeddedSF - removes edges from all but the selected roots, does not remap indices
96095fce210SBarry Smith 
96195fce210SBarry Smith    Collective
96295fce210SBarry Smith 
96395fce210SBarry Smith    Input Arguments:
96495fce210SBarry Smith +  sf - original star forest
965ba2a7774SJunchao Zhang .  nselected  - number of selected roots on this process
966ba2a7774SJunchao Zhang -  selected   - indices of the selected roots on this process
96795fce210SBarry Smith 
96895fce210SBarry Smith    Output Arguments:
96995fce210SBarry Smith .  newsf - new star forest
97095fce210SBarry Smith 
97195fce210SBarry Smith    Level: advanced
97295fce210SBarry Smith 
97395fce210SBarry Smith    Note:
97495fce210SBarry Smith    To use the new PetscSF, it may be necessary to know the indices of the leaves that are still participating. This can
97595fce210SBarry Smith    be done by calling PetscSFGetGraph().
97695fce210SBarry Smith 
97795fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFGetGraph()
97895fce210SBarry Smith @*/
979ba2a7774SJunchao Zhang PetscErrorCode PetscSFCreateEmbeddedSF(PetscSF sf,PetscInt nselected,const PetscInt *selected,PetscSF *newsf)
98095fce210SBarry Smith {
981ba2a7774SJunchao Zhang   PetscInt          *rootdata,*leafdata,*new_ilocal;
982ba2a7774SJunchao Zhang   PetscSFNode       *new_iremote;
983ba2a7774SJunchao Zhang   const PetscInt    *ilocal;
984ba2a7774SJunchao Zhang   const PetscSFNode *iremote;
985ba2a7774SJunchao Zhang   PetscInt          nleaves,nroots,n,i,new_nleaves = 0;
986ba2a7774SJunchao Zhang   PetscSF           tmpsf;
9870511a646SMatthew G. Knepley   PetscErrorCode    ierr;
98895fce210SBarry Smith 
98995fce210SBarry Smith   PetscFunctionBegin;
99095fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
99129046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
992ba2a7774SJunchao Zhang   if (nselected) PetscValidPointer(selected,3);
99395fce210SBarry Smith   PetscValidPointer(newsf,4);
994524e35f8SStefano Zampini   ierr = PetscLogEventBegin(PETSCSF_EmbedSF,sf,0,0,0);CHKERRQ(ierr);
9950511a646SMatthew G. Knepley 
996ba2a7774SJunchao Zhang   /* Find out which leaves (not leaf data items) are still connected to roots in the embedded sf */
997ba2a7774SJunchao Zhang   ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
998ba2a7774SJunchao Zhang   ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&tmpsf);CHKERRQ(ierr);
999ba2a7774SJunchao Zhang   ierr = PetscSFSetGraph(tmpsf,nroots,nleaves,NULL/*contiguous*/,PETSC_USE_POINTER,iremote,PETSC_USE_POINTER);CHKERRQ(ierr);
1000ba2a7774SJunchao Zhang   ierr = PetscCalloc2(nroots,&rootdata,nleaves,&leafdata);CHKERRQ(ierr);
1001ba2a7774SJunchao Zhang   for (i=0; i<nselected; ++i) {
100246fbb50bSJunchao Zhang     if (selected[i] < 0 || selected[i] >= nroots) SETERRQ2(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Root index %D is not in [0,%D)",selected[i],nroots);
1003ba2a7774SJunchao Zhang     rootdata[selected[i]] = 1;
1004ba2a7774SJunchao Zhang   }
1005ba2a7774SJunchao Zhang 
1006ba2a7774SJunchao Zhang   ierr = PetscSFBcastBegin(tmpsf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr);
1007ba2a7774SJunchao Zhang   ierr = PetscSFBcastEnd(tmpsf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr);
1008ba2a7774SJunchao Zhang   ierr = PetscSFDestroy(&tmpsf);CHKERRQ(ierr);
1009ba2a7774SJunchao Zhang 
1010ba2a7774SJunchao Zhang   /* Build newsf with leaves that are still connected */
1011ba2a7774SJunchao Zhang   for (i = 0; i < nleaves; ++i) new_nleaves += leafdata[i];
1012ba2a7774SJunchao Zhang   ierr = PetscMalloc1(new_nleaves,&new_ilocal);CHKERRQ(ierr);
1013ba2a7774SJunchao Zhang   ierr = PetscMalloc1(new_nleaves,&new_iremote);CHKERRQ(ierr);
1014ba2a7774SJunchao Zhang   for (i = 0, n = 0; i < nleaves; ++i) {
1015ba2a7774SJunchao Zhang     if (leafdata[i]) {
1016ba2a7774SJunchao Zhang       new_ilocal[n]        = sf->mine ? sf->mine[i] : i;
1017ba2a7774SJunchao Zhang       new_iremote[n].rank  = sf->remote[i].rank;
1018ba2a7774SJunchao Zhang       new_iremote[n].index = sf->remote[i].index;
1019fc1ede2bSMatthew G. Knepley       ++n;
102095fce210SBarry Smith     }
102195fce210SBarry Smith   }
102241bb73faSSatish Balay   if (n != new_nleaves) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"There is a size mismatch in the SF embedding, %D != %D",n,new_nleaves);
1023ba2a7774SJunchao Zhang   ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,newsf);CHKERRQ(ierr);
1024ba2a7774SJunchao Zhang   ierr = PetscSFSetGraph(*newsf,nroots,new_nleaves,new_ilocal,PETSC_OWN_POINTER,new_iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
102595fce210SBarry Smith   ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr);
1026524e35f8SStefano Zampini   ierr = PetscSFSetUp(*newsf);CHKERRQ(ierr);
1027524e35f8SStefano Zampini   ierr = PetscLogEventEnd(PETSCSF_EmbedSF,sf,0,0,0);CHKERRQ(ierr);
102895fce210SBarry Smith   PetscFunctionReturn(0);
102995fce210SBarry Smith }
103095fce210SBarry Smith 
10312f5fb4c2SMatthew G. Knepley /*@C
10322f5fb4c2SMatthew G. Knepley   PetscSFCreateEmbeddedLeafSF - removes edges from all but the selected leaves, does not remap indices
10332f5fb4c2SMatthew G. Knepley 
10342f5fb4c2SMatthew G. Knepley   Collective
10352f5fb4c2SMatthew G. Knepley 
10362f5fb4c2SMatthew G. Knepley   Input Arguments:
10372f5fb4c2SMatthew G. Knepley + sf - original star forest
10382f5fb4c2SMatthew G. Knepley . nleaves - number of leaves to select on this process
10392f5fb4c2SMatthew G. Knepley - selected - selected leaves on this process
10402f5fb4c2SMatthew G. Knepley 
10412f5fb4c2SMatthew G. Knepley   Output Arguments:
10422f5fb4c2SMatthew G. Knepley .  newsf - new star forest
10432f5fb4c2SMatthew G. Knepley 
10442f5fb4c2SMatthew G. Knepley   Level: advanced
10452f5fb4c2SMatthew G. Knepley 
10462f5fb4c2SMatthew G. Knepley .seealso: PetscSFCreateEmbeddedSF(), PetscSFSetGraph(), PetscSFGetGraph()
10472f5fb4c2SMatthew G. Knepley @*/
10482f5fb4c2SMatthew G. Knepley PetscErrorCode PetscSFCreateEmbeddedLeafSF(PetscSF sf, PetscInt nleaves, const PetscInt *selected, PetscSF *newsf)
10492f5fb4c2SMatthew G. Knepley {
10502f5fb4c2SMatthew G. Knepley   PetscSFNode   *iremote;
10512f5fb4c2SMatthew G. Knepley   PetscInt      *ilocal;
10522f5fb4c2SMatthew G. Knepley   PetscInt       i;
10532f5fb4c2SMatthew G. Knepley   PetscErrorCode ierr;
10542f5fb4c2SMatthew G. Knepley 
10552f5fb4c2SMatthew G. Knepley   PetscFunctionBegin;
10562f5fb4c2SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
105729046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf, 1);
10582f5fb4c2SMatthew G. Knepley   if (nleaves) PetscValidPointer(selected, 3);
10592f5fb4c2SMatthew G. Knepley   PetscValidPointer(newsf, 4);
10602f5fb4c2SMatthew G. Knepley   ierr = PetscMalloc1(nleaves, &ilocal);CHKERRQ(ierr);
10612f5fb4c2SMatthew G. Knepley   ierr = PetscMalloc1(nleaves, &iremote);CHKERRQ(ierr);
10622f5fb4c2SMatthew G. Knepley   for (i = 0; i < nleaves; ++i) {
10632f5fb4c2SMatthew G. Knepley     const PetscInt l = selected[i];
10642f5fb4c2SMatthew G. Knepley 
10652f5fb4c2SMatthew G. Knepley     ilocal[i]        = sf->mine ? sf->mine[l] : l;
10662f5fb4c2SMatthew G. Knepley     iremote[i].rank  = sf->remote[l].rank;
10672f5fb4c2SMatthew G. Knepley     iremote[i].index = sf->remote[l].index;
10682f5fb4c2SMatthew G. Knepley   }
10692f5fb4c2SMatthew G. Knepley   ierr = PetscSFDuplicate(sf, PETSCSF_DUPLICATE_RANKS, newsf);CHKERRQ(ierr);
10702f5fb4c2SMatthew G. Knepley   ierr = PetscSFSetGraph(*newsf, sf->nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr);
10712f5fb4c2SMatthew G. Knepley   PetscFunctionReturn(0);
10722f5fb4c2SMatthew G. Knepley }
10732f5fb4c2SMatthew G. Knepley 
107495fce210SBarry Smith /*@C
107595fce210SBarry Smith    PetscSFBcastBegin - begin pointwise broadcast to be concluded with call to PetscSFBcastEnd()
107695fce210SBarry Smith 
107795fce210SBarry Smith    Collective on PetscSF
107895fce210SBarry Smith 
107995fce210SBarry Smith    Input Arguments:
108095fce210SBarry Smith +  sf - star forest on which to communicate
108195fce210SBarry Smith .  unit - data type associated with each node
108295fce210SBarry Smith -  rootdata - buffer to broadcast
108395fce210SBarry Smith 
108495fce210SBarry Smith    Output Arguments:
108595fce210SBarry Smith .  leafdata - buffer to update with values from each leaf's respective root
108695fce210SBarry Smith 
108795fce210SBarry Smith    Level: intermediate
108895fce210SBarry Smith 
108995fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFView(), PetscSFBcastEnd(), PetscSFReduceBegin()
109095fce210SBarry Smith @*/
109195fce210SBarry Smith PetscErrorCode PetscSFBcastBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
109295fce210SBarry Smith {
109395fce210SBarry Smith   PetscErrorCode ierr;
109495fce210SBarry Smith 
109595fce210SBarry Smith   PetscFunctionBegin;
109695fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
109795fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
109829046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr);
109995fce210SBarry Smith   ierr = (*sf->ops->BcastBegin)(sf,unit,rootdata,leafdata);CHKERRQ(ierr);
1100563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr);
110195fce210SBarry Smith   PetscFunctionReturn(0);
110295fce210SBarry Smith }
110395fce210SBarry Smith 
110495fce210SBarry Smith /*@C
110595fce210SBarry Smith    PetscSFBcastEnd - end a broadcast operation started with PetscSFBcastBegin()
110695fce210SBarry Smith 
110795fce210SBarry Smith    Collective
110895fce210SBarry Smith 
110995fce210SBarry Smith    Input Arguments:
111095fce210SBarry Smith +  sf - star forest
111195fce210SBarry Smith .  unit - data type
111295fce210SBarry Smith -  rootdata - buffer to broadcast
111395fce210SBarry Smith 
111495fce210SBarry Smith    Output Arguments:
111595fce210SBarry Smith .  leafdata - buffer to update with values from each leaf's respective root
111695fce210SBarry Smith 
111795fce210SBarry Smith    Level: intermediate
111895fce210SBarry Smith 
111995fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFReduceEnd()
112095fce210SBarry Smith @*/
112195fce210SBarry Smith PetscErrorCode PetscSFBcastEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
112295fce210SBarry Smith {
112395fce210SBarry Smith   PetscErrorCode ierr;
112495fce210SBarry Smith 
112595fce210SBarry Smith   PetscFunctionBegin;
112695fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
112795fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
112829046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_BcastEnd,sf,0,0,0);CHKERRQ(ierr);
112995fce210SBarry Smith   ierr = (*sf->ops->BcastEnd)(sf,unit,rootdata,leafdata);CHKERRQ(ierr);
1130563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_BcastEnd,sf,0,0,0);CHKERRQ(ierr);
113195fce210SBarry Smith   PetscFunctionReturn(0);
113295fce210SBarry Smith }
113395fce210SBarry Smith 
113495fce210SBarry Smith /*@C
11353482bfa8SJunchao Zhang    PetscSFBcastAndOpBegin - begin pointwise broadcast with root value being reduced to leaf value, to be concluded with call to PetscSFBcastAndOpEnd()
11363482bfa8SJunchao Zhang 
11373482bfa8SJunchao Zhang    Collective on PetscSF
11383482bfa8SJunchao Zhang 
11393482bfa8SJunchao Zhang    Input Arguments:
11403482bfa8SJunchao Zhang +  sf - star forest on which to communicate
11413482bfa8SJunchao Zhang .  unit - data type associated with each node
11423482bfa8SJunchao Zhang .  rootdata - buffer to broadcast
11433482bfa8SJunchao Zhang -  op - operation to use for reduction
11443482bfa8SJunchao Zhang 
11453482bfa8SJunchao Zhang    Output Arguments:
11463482bfa8SJunchao Zhang .  leafdata - buffer to be reduced with values from each leaf's respective root
11473482bfa8SJunchao Zhang 
11483482bfa8SJunchao Zhang    Level: intermediate
11493482bfa8SJunchao Zhang 
11503482bfa8SJunchao Zhang .seealso: PetscSFBcastAndOpEnd(), PetscSFBcastBegin(), PetscSFBcastEnd()
11513482bfa8SJunchao Zhang @*/
11523482bfa8SJunchao Zhang PetscErrorCode PetscSFBcastAndOpBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
11533482bfa8SJunchao Zhang {
11543482bfa8SJunchao Zhang   PetscErrorCode ierr;
11553482bfa8SJunchao Zhang 
11563482bfa8SJunchao Zhang   PetscFunctionBegin;
11573482bfa8SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
11583482bfa8SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
11593482bfa8SJunchao Zhang   ierr = PetscLogEventBegin(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
11603482bfa8SJunchao Zhang   ierr = (*sf->ops->BcastAndOpBegin)(sf,unit,rootdata,leafdata,op);CHKERRQ(ierr);
11613482bfa8SJunchao Zhang   ierr = PetscLogEventEnd(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
11623482bfa8SJunchao Zhang   PetscFunctionReturn(0);
11633482bfa8SJunchao Zhang }
11643482bfa8SJunchao Zhang 
11653482bfa8SJunchao Zhang /*@C
11663482bfa8SJunchao Zhang    PetscSFBcastAndOpEnd - end a broadcast & reduce operation started with PetscSFBcastAndOpBegin()
11673482bfa8SJunchao Zhang 
11683482bfa8SJunchao Zhang    Collective
11693482bfa8SJunchao Zhang 
11703482bfa8SJunchao Zhang    Input Arguments:
11713482bfa8SJunchao Zhang +  sf - star forest
11723482bfa8SJunchao Zhang .  unit - data type
11733482bfa8SJunchao Zhang .  rootdata - buffer to broadcast
11743482bfa8SJunchao Zhang -  op - operation to use for reduction
11753482bfa8SJunchao Zhang 
11763482bfa8SJunchao Zhang    Output Arguments:
11773482bfa8SJunchao Zhang .  leafdata - buffer to be reduced with values from each leaf's respective root
11783482bfa8SJunchao Zhang 
11793482bfa8SJunchao Zhang    Level: intermediate
11803482bfa8SJunchao Zhang 
11813482bfa8SJunchao Zhang .seealso: PetscSFSetGraph(), PetscSFReduceEnd()
11823482bfa8SJunchao Zhang @*/
11833482bfa8SJunchao Zhang PetscErrorCode PetscSFBcastAndOpEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
11843482bfa8SJunchao Zhang {
11853482bfa8SJunchao Zhang   PetscErrorCode ierr;
11863482bfa8SJunchao Zhang 
11873482bfa8SJunchao Zhang   PetscFunctionBegin;
11883482bfa8SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
11893482bfa8SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
11903482bfa8SJunchao Zhang   ierr = PetscLogEventBegin(PETSCSF_BcastAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
11913482bfa8SJunchao Zhang   ierr = (*sf->ops->BcastAndOpEnd)(sf,unit,rootdata,leafdata,op);CHKERRQ(ierr);
11923482bfa8SJunchao Zhang   ierr = PetscLogEventEnd(PETSCSF_BcastAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
11933482bfa8SJunchao Zhang   PetscFunctionReturn(0);
11943482bfa8SJunchao Zhang }
11953482bfa8SJunchao Zhang 
11963482bfa8SJunchao Zhang /*@C
119795fce210SBarry Smith    PetscSFReduceBegin - begin reduction of leafdata into rootdata, to be completed with call to PetscSFReduceEnd()
119895fce210SBarry Smith 
119995fce210SBarry Smith    Collective
120095fce210SBarry Smith 
120195fce210SBarry Smith    Input Arguments:
120295fce210SBarry Smith +  sf - star forest
120395fce210SBarry Smith .  unit - data type
120495fce210SBarry Smith .  leafdata - values to reduce
120595fce210SBarry Smith -  op - reduction operation
120695fce210SBarry Smith 
120795fce210SBarry Smith    Output Arguments:
120895fce210SBarry Smith .  rootdata - result of reduction of values from all leaves of each root
120995fce210SBarry Smith 
121095fce210SBarry Smith    Level: intermediate
121195fce210SBarry Smith 
121295fce210SBarry Smith .seealso: PetscSFBcastBegin()
121395fce210SBarry Smith @*/
121495fce210SBarry Smith PetscErrorCode PetscSFReduceBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
121595fce210SBarry Smith {
121695fce210SBarry Smith   PetscErrorCode ierr;
121795fce210SBarry Smith 
121895fce210SBarry Smith   PetscFunctionBegin;
121995fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
122095fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
122129046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);
122295fce210SBarry Smith   ierr = (sf->ops->ReduceBegin)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
1223563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);
122495fce210SBarry Smith   PetscFunctionReturn(0);
122595fce210SBarry Smith }
122695fce210SBarry Smith 
122795fce210SBarry Smith /*@C
122895fce210SBarry Smith    PetscSFReduceEnd - end a reduction operation started with PetscSFReduceBegin()
122995fce210SBarry Smith 
123095fce210SBarry Smith    Collective
123195fce210SBarry Smith 
123295fce210SBarry Smith    Input Arguments:
123395fce210SBarry Smith +  sf - star forest
123495fce210SBarry Smith .  unit - data type
123595fce210SBarry Smith .  leafdata - values to reduce
123695fce210SBarry Smith -  op - reduction operation
123795fce210SBarry Smith 
123895fce210SBarry Smith    Output Arguments:
123995fce210SBarry Smith .  rootdata - result of reduction of values from all leaves of each root
124095fce210SBarry Smith 
124195fce210SBarry Smith    Level: intermediate
124295fce210SBarry Smith 
124395fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFBcastEnd()
124495fce210SBarry Smith @*/
124595fce210SBarry Smith PetscErrorCode PetscSFReduceEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
124695fce210SBarry Smith {
124795fce210SBarry Smith   PetscErrorCode ierr;
124895fce210SBarry Smith 
124995fce210SBarry Smith   PetscFunctionBegin;
125095fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
125195fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
125229046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);
125395fce210SBarry Smith   ierr = (*sf->ops->ReduceEnd)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
1254563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);
125595fce210SBarry Smith   PetscFunctionReturn(0);
125695fce210SBarry Smith }
125795fce210SBarry Smith 
125895fce210SBarry Smith /*@C
125995fce210SBarry Smith    PetscSFComputeDegreeBegin - begin computation of degree for each root vertex, to be completed with PetscSFComputeDegreeEnd()
126095fce210SBarry Smith 
126195fce210SBarry Smith    Collective
126295fce210SBarry Smith 
126395fce210SBarry Smith    Input Arguments:
126495fce210SBarry Smith .  sf - star forest
126595fce210SBarry Smith 
126695fce210SBarry Smith    Output Arguments:
126795fce210SBarry Smith .  degree - degree of each root vertex
126895fce210SBarry Smith 
126995fce210SBarry Smith    Level: advanced
127095fce210SBarry Smith 
1271ffe67aa5SVáclav Hapla    Notes:
1272ffe67aa5SVáclav Hapla    The returned array is owned by PetscSF and automatically freed by PetscSFDestroy(). Hence no need to call PetscFree() on it.
1273ffe67aa5SVáclav Hapla 
127495fce210SBarry Smith .seealso: PetscSFGatherBegin()
127595fce210SBarry Smith @*/
127695fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf,const PetscInt **degree)
127795fce210SBarry Smith {
127895fce210SBarry Smith   PetscErrorCode ierr;
127995fce210SBarry Smith 
128095fce210SBarry Smith   PetscFunctionBegin;
128195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
128295fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
128395fce210SBarry Smith   PetscValidPointer(degree,2);
1284803bd9e8SMatthew G. Knepley   if (!sf->degreeknown) {
128529046d53SLisandro Dalcin     PetscInt i, nroots = sf->nroots, maxlocal = sf->maxleaf+1;  /* TODO: We should use PetscSFGetLeafRange() */
1286803bd9e8SMatthew G. Knepley     if (sf->degree) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Calls to PetscSFComputeDegreeBegin() cannot be nested.");
128729046d53SLisandro Dalcin     ierr = PetscMalloc1(nroots,&sf->degree);CHKERRQ(ierr);
128829046d53SLisandro Dalcin     ierr = PetscMalloc1(PetscMax(maxlocal,1),&sf->degreetmp);CHKERRQ(ierr); /* allocate at least one entry, see check in PetscSFComputeDegreeEnd() */
128929046d53SLisandro Dalcin     for (i=0; i<nroots; i++) sf->degree[i] = 0;
12909837ea96SMatthew G. Knepley     for (i=0; i<maxlocal; i++) sf->degreetmp[i] = 1;
1291dbd2ff41SBarry Smith     ierr = PetscSFReduceBegin(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);CHKERRQ(ierr);
129295fce210SBarry Smith   }
129395fce210SBarry Smith   *degree = NULL;
129495fce210SBarry Smith   PetscFunctionReturn(0);
129595fce210SBarry Smith }
129695fce210SBarry Smith 
129795fce210SBarry Smith /*@C
129895fce210SBarry Smith    PetscSFComputeDegreeEnd - complete computation of degree for each root vertex, started with PetscSFComputeDegreeBegin()
129995fce210SBarry Smith 
130095fce210SBarry Smith    Collective
130195fce210SBarry Smith 
130295fce210SBarry Smith    Input Arguments:
130395fce210SBarry Smith .  sf - star forest
130495fce210SBarry Smith 
130595fce210SBarry Smith    Output Arguments:
130695fce210SBarry Smith .  degree - degree of each root vertex
130795fce210SBarry Smith 
130895fce210SBarry Smith    Level: developer
130995fce210SBarry Smith 
1310ffe67aa5SVáclav Hapla    Notes:
1311ffe67aa5SVáclav Hapla    The returned array is owned by PetscSF and automatically freed by PetscSFDestroy(). Hence no need to call PetscFree() on it.
1312ffe67aa5SVáclav Hapla 
131395fce210SBarry Smith .seealso:
131495fce210SBarry Smith @*/
131595fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf,const PetscInt **degree)
131695fce210SBarry Smith {
131795fce210SBarry Smith   PetscErrorCode ierr;
131895fce210SBarry Smith 
131995fce210SBarry Smith   PetscFunctionBegin;
132095fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
132195fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
132229046d53SLisandro Dalcin   PetscValidPointer(degree,2);
132395fce210SBarry Smith   if (!sf->degreeknown) {
132429046d53SLisandro Dalcin     if (!sf->degreetmp) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFComputeDegreeBegin() before PetscSFComputeDegreeEnd()");
1325dbd2ff41SBarry Smith     ierr = PetscSFReduceEnd(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);CHKERRQ(ierr);
132695fce210SBarry Smith     ierr = PetscFree(sf->degreetmp);CHKERRQ(ierr);
132795fce210SBarry Smith     sf->degreeknown = PETSC_TRUE;
132895fce210SBarry Smith   }
132995fce210SBarry Smith   *degree = sf->degree;
133095fce210SBarry Smith   PetscFunctionReturn(0);
133195fce210SBarry Smith }
133295fce210SBarry Smith 
1333673100f5SVaclav Hapla 
1334673100f5SVaclav Hapla /*@C
133566dfcd1aSVaclav Hapla    PetscSFComputeMultiRootOriginalNumbering - Returns original numbering of multi-roots (roots of multi-SF returned by PetscSFGetMultiSF()).
133666dfcd1aSVaclav Hapla    Each multi-root is assigned index of the corresponding original root.
1337673100f5SVaclav Hapla 
1338673100f5SVaclav Hapla    Collective
1339673100f5SVaclav Hapla 
1340673100f5SVaclav Hapla    Input Arguments:
1341673100f5SVaclav Hapla +  sf - star forest
1342673100f5SVaclav Hapla -  degree - degree of each root vertex, computed with PetscSFComputeDegreeBegin()/PetscSFComputeDegreeEnd()
1343673100f5SVaclav Hapla 
1344673100f5SVaclav Hapla    Output Arguments:
134566dfcd1aSVaclav Hapla +  nMultiRoots - (optional) number of multi-roots (roots of multi-SF)
134666dfcd1aSVaclav Hapla -  multiRootsOrigNumbering - original indices of multi-roots; length of this array is nMultiRoots
1347673100f5SVaclav Hapla 
1348673100f5SVaclav Hapla    Level: developer
1349673100f5SVaclav Hapla 
1350ffe67aa5SVáclav Hapla    Notes:
1351ffe67aa5SVáclav Hapla    The returned array multiRootsOrigNumbering is newly allocated and should be destroyed with PetscFree() when no longer needed.
1352ffe67aa5SVáclav Hapla 
1353673100f5SVaclav Hapla .seealso: PetscSFComputeDegreeBegin(), PetscSFComputeDegreeEnd(), PetscSFGetMultiSF()
1354673100f5SVaclav Hapla @*/
135566dfcd1aSVaclav Hapla PetscErrorCode PetscSFComputeMultiRootOriginalNumbering(PetscSF sf, const PetscInt degree[], PetscInt *nMultiRoots, PetscInt *multiRootsOrigNumbering[])
1356673100f5SVaclav Hapla {
1357673100f5SVaclav Hapla   PetscSF             msf;
1358673100f5SVaclav Hapla   PetscInt            i, j, k, nroots, nmroots;
1359673100f5SVaclav Hapla   PetscErrorCode      ierr;
1360673100f5SVaclav Hapla 
1361673100f5SVaclav Hapla   PetscFunctionBegin;
1362673100f5SVaclav Hapla   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1363673100f5SVaclav Hapla   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
136429328920SVaclav Hapla   if (nroots) PetscValidIntPointer(degree,2);
136566dfcd1aSVaclav Hapla   if (nMultiRoots) PetscValidIntPointer(nMultiRoots,3);
136666dfcd1aSVaclav Hapla   PetscValidPointer(multiRootsOrigNumbering,4);
1367673100f5SVaclav Hapla   ierr = PetscSFGetMultiSF(sf,&msf);CHKERRQ(ierr);
1368673100f5SVaclav Hapla   ierr = PetscSFGetGraph(msf, &nmroots, NULL, NULL, NULL);CHKERRQ(ierr);
136966dfcd1aSVaclav Hapla   ierr = PetscMalloc1(nmroots, multiRootsOrigNumbering);CHKERRQ(ierr);
1370673100f5SVaclav Hapla   for (i=0,j=0,k=0; i<nroots; i++) {
1371673100f5SVaclav Hapla     if (!degree[i]) continue;
1372673100f5SVaclav Hapla     for (j=0; j<degree[i]; j++,k++) {
137366dfcd1aSVaclav Hapla       (*multiRootsOrigNumbering)[k] = i;
1374673100f5SVaclav Hapla     }
1375673100f5SVaclav Hapla   }
1376673100f5SVaclav Hapla #if defined(PETSC_USE_DEBUG)
1377673100f5SVaclav Hapla   if (PetscUnlikely(k != nmroots)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sanity check fail");
1378673100f5SVaclav Hapla #endif
137966dfcd1aSVaclav Hapla   if (nMultiRoots) *nMultiRoots = nmroots;
1380673100f5SVaclav Hapla   PetscFunctionReturn(0);
1381673100f5SVaclav Hapla }
1382673100f5SVaclav Hapla 
138395fce210SBarry Smith /*@C
138495fce210SBarry Smith    PetscSFFetchAndOpBegin - begin operation that fetches values from root and updates atomically by applying operation using my leaf value, to be completed with PetscSFFetchAndOpEnd()
138595fce210SBarry Smith 
138695fce210SBarry Smith    Collective
138795fce210SBarry Smith 
138895fce210SBarry Smith    Input Arguments:
138995fce210SBarry Smith +  sf - star forest
139095fce210SBarry Smith .  unit - data type
139195fce210SBarry Smith .  leafdata - leaf values to use in reduction
139295fce210SBarry Smith -  op - operation to use for reduction
139395fce210SBarry Smith 
139495fce210SBarry Smith    Output Arguments:
139595fce210SBarry Smith +  rootdata - root values to be updated, input state is seen by first process to perform an update
139695fce210SBarry Smith -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
139795fce210SBarry Smith 
139895fce210SBarry Smith    Level: advanced
139995fce210SBarry Smith 
140095fce210SBarry Smith    Note:
140195fce210SBarry Smith    The update is only atomic at the granularity provided by the hardware. Different roots referenced by the same process
140295fce210SBarry Smith    might be updated in a different order. Furthermore, if a composite type is used for the unit datatype, atomicity is
140395fce210SBarry Smith    not guaranteed across the whole vertex. Therefore, this function is mostly only used with primitive types such as
140495fce210SBarry Smith    integers.
140595fce210SBarry Smith 
140695fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFReduceBegin(), PetscSFSetGraph()
140795fce210SBarry Smith @*/
140895fce210SBarry Smith PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
140995fce210SBarry Smith {
141095fce210SBarry Smith   PetscErrorCode ierr;
141195fce210SBarry Smith 
141295fce210SBarry Smith   PetscFunctionBegin;
141395fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
141495fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
141529046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
141695fce210SBarry Smith   ierr = (*sf->ops->FetchAndOpBegin)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr);
1417563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
141895fce210SBarry Smith   PetscFunctionReturn(0);
141995fce210SBarry Smith }
142095fce210SBarry Smith 
142195fce210SBarry Smith /*@C
142295fce210SBarry 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
142395fce210SBarry Smith 
142495fce210SBarry Smith    Collective
142595fce210SBarry Smith 
142695fce210SBarry Smith    Input Arguments:
142795fce210SBarry Smith +  sf - star forest
142895fce210SBarry Smith .  unit - data type
142995fce210SBarry Smith .  leafdata - leaf values to use in reduction
143095fce210SBarry Smith -  op - operation to use for reduction
143195fce210SBarry Smith 
143295fce210SBarry Smith    Output Arguments:
143395fce210SBarry Smith +  rootdata - root values to be updated, input state is seen by first process to perform an update
143495fce210SBarry Smith -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
143595fce210SBarry Smith 
143695fce210SBarry Smith    Level: advanced
143795fce210SBarry Smith 
143895fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFReduceEnd(), PetscSFSetGraph()
143995fce210SBarry Smith @*/
144095fce210SBarry Smith PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
144195fce210SBarry Smith {
144295fce210SBarry Smith   PetscErrorCode ierr;
144395fce210SBarry Smith 
144495fce210SBarry Smith   PetscFunctionBegin;
144595fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
144695fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
144729046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
144895fce210SBarry Smith   ierr = (*sf->ops->FetchAndOpEnd)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr);
1449563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
145095fce210SBarry Smith   PetscFunctionReturn(0);
145195fce210SBarry Smith }
145295fce210SBarry Smith 
145395fce210SBarry Smith /*@C
145495fce210SBarry Smith    PetscSFGatherBegin - begin pointwise gather of all leaves into multi-roots, to be completed with PetscSFGatherEnd()
145595fce210SBarry Smith 
145695fce210SBarry Smith    Collective
145795fce210SBarry Smith 
145895fce210SBarry Smith    Input Arguments:
145995fce210SBarry Smith +  sf - star forest
146095fce210SBarry Smith .  unit - data type
146195fce210SBarry Smith -  leafdata - leaf data to gather to roots
146295fce210SBarry Smith 
146395fce210SBarry Smith    Output Argument:
146495fce210SBarry Smith .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
146595fce210SBarry Smith 
146695fce210SBarry Smith    Level: intermediate
146795fce210SBarry Smith 
146895fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
146995fce210SBarry Smith @*/
147095fce210SBarry Smith PetscErrorCode PetscSFGatherBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
147195fce210SBarry Smith {
147295fce210SBarry Smith   PetscErrorCode ierr;
147395fce210SBarry Smith   PetscSF        multi;
147495fce210SBarry Smith 
147595fce210SBarry Smith   PetscFunctionBegin;
147695fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
147729046d53SLisandro Dalcin   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
147895fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
14798bfbc91cSJed Brown   ierr = PetscSFReduceBegin(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr);
148095fce210SBarry Smith   PetscFunctionReturn(0);
148195fce210SBarry Smith }
148295fce210SBarry Smith 
148395fce210SBarry Smith /*@C
148495fce210SBarry Smith    PetscSFGatherEnd - ends pointwise gather operation that was started with PetscSFGatherBegin()
148595fce210SBarry Smith 
148695fce210SBarry Smith    Collective
148795fce210SBarry Smith 
148895fce210SBarry Smith    Input Arguments:
148995fce210SBarry Smith +  sf - star forest
149095fce210SBarry Smith .  unit - data type
149195fce210SBarry Smith -  leafdata - leaf data to gather to roots
149295fce210SBarry Smith 
149395fce210SBarry Smith    Output Argument:
149495fce210SBarry Smith .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
149595fce210SBarry Smith 
149695fce210SBarry Smith    Level: intermediate
149795fce210SBarry Smith 
149895fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
149995fce210SBarry Smith @*/
150095fce210SBarry Smith PetscErrorCode PetscSFGatherEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
150195fce210SBarry Smith {
150295fce210SBarry Smith   PetscErrorCode ierr;
150395fce210SBarry Smith   PetscSF        multi;
150495fce210SBarry Smith 
150595fce210SBarry Smith   PetscFunctionBegin;
150695fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
150795fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
150895fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
15098bfbc91cSJed Brown   ierr = PetscSFReduceEnd(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr);
151095fce210SBarry Smith   PetscFunctionReturn(0);
151195fce210SBarry Smith }
151295fce210SBarry Smith 
151395fce210SBarry Smith /*@C
151495fce210SBarry Smith    PetscSFScatterBegin - begin pointwise scatter operation from multi-roots to leaves, to be completed with PetscSFScatterEnd()
151595fce210SBarry Smith 
151695fce210SBarry Smith    Collective
151795fce210SBarry Smith 
151895fce210SBarry Smith    Input Arguments:
151995fce210SBarry Smith +  sf - star forest
152095fce210SBarry Smith .  unit - data type
152195fce210SBarry Smith -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
152295fce210SBarry Smith 
152395fce210SBarry Smith    Output Argument:
152495fce210SBarry Smith .  leafdata - leaf data to be update with personal data from each respective root
152595fce210SBarry Smith 
152695fce210SBarry Smith    Level: intermediate
152795fce210SBarry Smith 
152895fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
152995fce210SBarry Smith @*/
153095fce210SBarry Smith PetscErrorCode PetscSFScatterBegin(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
153195fce210SBarry Smith {
153295fce210SBarry Smith   PetscErrorCode ierr;
153395fce210SBarry Smith   PetscSF        multi;
153495fce210SBarry Smith 
153595fce210SBarry Smith   PetscFunctionBegin;
153695fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
153795fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
153895fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
153995fce210SBarry Smith   ierr = PetscSFBcastBegin(multi,unit,multirootdata,leafdata);CHKERRQ(ierr);
154095fce210SBarry Smith   PetscFunctionReturn(0);
154195fce210SBarry Smith }
154295fce210SBarry Smith 
154395fce210SBarry Smith /*@C
154495fce210SBarry Smith    PetscSFScatterEnd - ends pointwise scatter operation that was started with PetscSFScatterBegin()
154595fce210SBarry Smith 
154695fce210SBarry Smith    Collective
154795fce210SBarry Smith 
154895fce210SBarry Smith    Input Arguments:
154995fce210SBarry Smith +  sf - star forest
155095fce210SBarry Smith .  unit - data type
155195fce210SBarry Smith -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
155295fce210SBarry Smith 
155395fce210SBarry Smith    Output Argument:
155495fce210SBarry Smith .  leafdata - leaf data to be update with personal data from each respective root
155595fce210SBarry Smith 
155695fce210SBarry Smith    Level: intermediate
155795fce210SBarry Smith 
155895fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
155995fce210SBarry Smith @*/
156095fce210SBarry Smith PetscErrorCode PetscSFScatterEnd(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
156195fce210SBarry Smith {
156295fce210SBarry Smith   PetscErrorCode ierr;
156395fce210SBarry Smith   PetscSF        multi;
156495fce210SBarry Smith 
156595fce210SBarry Smith   PetscFunctionBegin;
156695fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
156795fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
156895fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
156995fce210SBarry Smith   ierr = PetscSFBcastEnd(multi,unit,multirootdata,leafdata);CHKERRQ(ierr);
157095fce210SBarry Smith   PetscFunctionReturn(0);
157195fce210SBarry Smith }
1572a7b3aa13SAta Mesgarnejad 
1573a7b3aa13SAta Mesgarnejad /*@
1574a7b3aa13SAta Mesgarnejad   PetscSFCompose - Compose a new PetscSF equivalent to action to PetscSFs
1575a7b3aa13SAta Mesgarnejad 
1576a7b3aa13SAta Mesgarnejad   Input Parameters:
1577a7b3aa13SAta Mesgarnejad + sfA - The first PetscSF
1578a7b3aa13SAta Mesgarnejad - sfB - The second PetscSF
1579a7b3aa13SAta Mesgarnejad 
1580a7b3aa13SAta Mesgarnejad   Output Parameters:
1581a7b3aa13SAta Mesgarnejad . sfBA - equvalent PetscSF for applying A then B
1582a7b3aa13SAta Mesgarnejad 
1583a7b3aa13SAta Mesgarnejad   Level: developer
1584a7b3aa13SAta Mesgarnejad 
1585a7b3aa13SAta Mesgarnejad .seealso: PetscSF, PetscSFGetGraph(), PetscSFSetGraph()
1586a7b3aa13SAta Mesgarnejad @*/
1587a7b3aa13SAta Mesgarnejad PetscErrorCode PetscSFCompose(PetscSF sfA, PetscSF sfB, PetscSF *sfBA)
1588a7b3aa13SAta Mesgarnejad {
1589a7b3aa13SAta Mesgarnejad   MPI_Comm           comm;
1590a7b3aa13SAta Mesgarnejad   const PetscSFNode *remotePointsA, *remotePointsB;
1591a7b3aa13SAta Mesgarnejad   PetscSFNode       *remotePointsBA;
1592a7b3aa13SAta Mesgarnejad   const PetscInt    *localPointsA, *localPointsB;
1593a7b3aa13SAta Mesgarnejad   PetscInt           numRootsA, numLeavesA, numRootsB, numLeavesB;
1594a7b3aa13SAta Mesgarnejad   PetscErrorCode     ierr;
1595a7b3aa13SAta Mesgarnejad 
1596a7b3aa13SAta Mesgarnejad   PetscFunctionBegin;
1597a7b3aa13SAta Mesgarnejad   PetscValidHeaderSpecific(sfA, PETSCSF_CLASSID, 1);
159829046d53SLisandro Dalcin   PetscSFCheckGraphSet(sfA, 1);
159929046d53SLisandro Dalcin   PetscValidHeaderSpecific(sfB, PETSCSF_CLASSID, 2);
160029046d53SLisandro Dalcin   PetscSFCheckGraphSet(sfB, 2);
160129046d53SLisandro Dalcin   PetscValidPointer(sfBA, 3);
1602a7b3aa13SAta Mesgarnejad   ierr = PetscObjectGetComm((PetscObject) sfA, &comm);CHKERRQ(ierr);
1603a7b3aa13SAta Mesgarnejad   ierr = PetscSFGetGraph(sfA, &numRootsA, &numLeavesA, &localPointsA, &remotePointsA);CHKERRQ(ierr);
1604a7b3aa13SAta Mesgarnejad   ierr = PetscSFGetGraph(sfB, &numRootsB, &numLeavesB, &localPointsB, &remotePointsB);CHKERRQ(ierr);
1605a7b3aa13SAta Mesgarnejad   ierr = PetscMalloc1(numLeavesB, &remotePointsBA);CHKERRQ(ierr);
1606a7b3aa13SAta Mesgarnejad   ierr = PetscSFBcastBegin(sfB, MPIU_2INT, remotePointsA, remotePointsBA);CHKERRQ(ierr);
1607a7b3aa13SAta Mesgarnejad   ierr = PetscSFBcastEnd(sfB, MPIU_2INT, remotePointsA, remotePointsBA);CHKERRQ(ierr);
1608a7b3aa13SAta Mesgarnejad   ierr = PetscSFCreate(comm, sfBA);CHKERRQ(ierr);
1609a7b3aa13SAta Mesgarnejad   ierr = PetscSFSetGraph(*sfBA, numRootsA, numLeavesB, localPointsB, PETSC_COPY_VALUES, remotePointsBA, PETSC_OWN_POINTER);CHKERRQ(ierr);
1610a7b3aa13SAta Mesgarnejad   PetscFunctionReturn(0);
1611a7b3aa13SAta Mesgarnejad }
1612