xref: /petsc/src/vec/is/sf/interface/sf.c (revision 80153354bc763d7777ef13196fd03cfc4d2e02f9)
1af0996ceSBarry Smith #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/
295fce210SBarry Smith #include <petscctable.h>
395fce210SBarry Smith 
4563ffbabSMatthew G. Knepley /* Logging support */
5563ffbabSMatthew G. Knepley PetscLogEvent PETSCSF_SetGraph, PETSCSF_BcastBegin, PETSCSF_BcastEnd, PETSCSF_ReduceBegin, PETSCSF_ReduceEnd, PETSCSF_FetchAndOpBegin, PETSCSF_FetchAndOpEnd;
6563ffbabSMatthew G. Knepley 
795fce210SBarry Smith #if defined(PETSC_USE_DEBUG)
895fce210SBarry Smith #  define PetscSFCheckGraphSet(sf,arg) do {                          \
995fce210SBarry Smith     if (PetscUnlikely(!(sf)->graphset))                              \
1095fce210SBarry Smith       SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetGraph() on argument %D \"%s\" before %s()",(arg),#sf,PETSC_FUNCTION_NAME); \
1195fce210SBarry Smith   } while (0)
1295fce210SBarry Smith #else
1395fce210SBarry Smith #  define PetscSFCheckGraphSet(sf,arg) do {} while (0)
1495fce210SBarry Smith #endif
1595fce210SBarry Smith 
1695fce210SBarry Smith const char *const PetscSFDuplicateOptions[] = {"CONFONLY","RANKS","GRAPH","PetscSFDuplicateOption","PETSCSF_DUPLICATE_",0};
1795fce210SBarry Smith 
188af6ec1cSBarry Smith /*@
1995fce210SBarry Smith    PetscSFCreate - create a star forest communication context
2095fce210SBarry Smith 
2195fce210SBarry Smith    Not Collective
2295fce210SBarry Smith 
2395fce210SBarry Smith    Input Arguments:
2495fce210SBarry Smith .  comm - communicator on which the star forest will operate
2595fce210SBarry Smith 
2695fce210SBarry Smith    Output Arguments:
2795fce210SBarry Smith .  sf - new star forest context
2895fce210SBarry Smith 
2995fce210SBarry Smith    Level: intermediate
3095fce210SBarry Smith 
3195fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFDestroy()
3295fce210SBarry Smith @*/
3395fce210SBarry Smith PetscErrorCode PetscSFCreate(MPI_Comm comm,PetscSF *sf)
3495fce210SBarry Smith {
3595fce210SBarry Smith   PetscErrorCode ierr;
3695fce210SBarry Smith   PetscSF        b;
3795fce210SBarry Smith 
3895fce210SBarry Smith   PetscFunctionBegin;
3995fce210SBarry Smith   PetscValidPointer(sf,2);
40607a6623SBarry Smith   ierr = PetscSFInitializePackage();CHKERRQ(ierr);
4195fce210SBarry Smith 
4273107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(b,PETSCSF_CLASSID,"PetscSF","Star Forest","PetscSF",comm,PetscSFDestroy,PetscSFView);CHKERRQ(ierr);
4395fce210SBarry Smith 
4495fce210SBarry Smith   b->nroots    = -1;
4595fce210SBarry Smith   b->nleaves   = -1;
4695fce210SBarry Smith   b->nranks    = -1;
4795fce210SBarry Smith   b->rankorder = PETSC_TRUE;
4895fce210SBarry Smith   b->ingroup   = MPI_GROUP_NULL;
4995fce210SBarry Smith   b->outgroup  = MPI_GROUP_NULL;
5095fce210SBarry Smith   b->graphset  = PETSC_FALSE;
5195fce210SBarry Smith 
5295fce210SBarry Smith   *sf = b;
5395fce210SBarry Smith   PetscFunctionReturn(0);
5495fce210SBarry Smith }
5595fce210SBarry Smith 
5695fce210SBarry Smith /*@C
5795fce210SBarry Smith    PetscSFReset - Reset a star forest so that different sizes or neighbors can be used
5895fce210SBarry Smith 
5995fce210SBarry Smith    Collective
6095fce210SBarry Smith 
6195fce210SBarry Smith    Input Arguments:
6295fce210SBarry Smith .  sf - star forest
6395fce210SBarry Smith 
6495fce210SBarry Smith    Level: advanced
6595fce210SBarry Smith 
6695fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFDestroy()
6795fce210SBarry Smith @*/
6895fce210SBarry Smith PetscErrorCode PetscSFReset(PetscSF sf)
6995fce210SBarry Smith {
7095fce210SBarry Smith   PetscErrorCode ierr;
7195fce210SBarry Smith 
7295fce210SBarry Smith   PetscFunctionBegin;
7395fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
7479715d56SJed Brown   if (sf->ops->Reset) {ierr = (*sf->ops->Reset)(sf);CHKERRQ(ierr);}
7595fce210SBarry Smith   sf->mine   = NULL;
7695fce210SBarry Smith   ierr       = PetscFree(sf->mine_alloc);CHKERRQ(ierr);
7795fce210SBarry Smith   sf->remote = NULL;
7895fce210SBarry Smith   ierr       = PetscFree(sf->remote_alloc);CHKERRQ(ierr);
7995fce210SBarry Smith   ierr       = PetscFree4(sf->ranks,sf->roffset,sf->rmine,sf->rremote);CHKERRQ(ierr);
8021c688dcSJed Brown   sf->nranks = -1;
8195fce210SBarry Smith   ierr       = PetscFree(sf->degree);CHKERRQ(ierr);
8295fce210SBarry Smith   if (sf->ingroup  != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->ingroup);CHKERRQ(ierr);}
8395fce210SBarry Smith   if (sf->outgroup != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->outgroup);CHKERRQ(ierr);}
8495fce210SBarry Smith   ierr         = PetscSFDestroy(&sf->multi);CHKERRQ(ierr);
8595fce210SBarry Smith   sf->graphset = PETSC_FALSE;
8695fce210SBarry Smith   sf->setupcalled = PETSC_FALSE;
8795fce210SBarry Smith   PetscFunctionReturn(0);
8895fce210SBarry Smith }
8995fce210SBarry Smith 
9095fce210SBarry Smith /*@C
9195fce210SBarry Smith    PetscSFSetType - set the PetscSF communication implementation
9295fce210SBarry Smith 
9395fce210SBarry Smith    Collective on PetscSF
9495fce210SBarry Smith 
9595fce210SBarry Smith    Input Parameters:
9695fce210SBarry Smith +  sf - the PetscSF context
9795fce210SBarry Smith -  type - a known method
9895fce210SBarry Smith 
9995fce210SBarry Smith    Options Database Key:
10095fce210SBarry Smith .  -sf_type <type> - Sets the method; use -help for a list
10195fce210SBarry Smith    of available methods (for instance, window, pt2pt, neighbor)
10295fce210SBarry Smith 
10395fce210SBarry Smith    Notes:
10495fce210SBarry Smith    See "include/petscsf.h" for available methods (for instance)
10595fce210SBarry Smith +    PETSCSFWINDOW - MPI-2/3 one-sided
10695fce210SBarry Smith -    PETSCSFBASIC - basic implementation using MPI-1 two-sided
10795fce210SBarry Smith 
10895fce210SBarry Smith   Level: intermediate
10995fce210SBarry Smith 
11095fce210SBarry Smith .keywords: PetscSF, set, type
11195fce210SBarry Smith 
11295fce210SBarry Smith .seealso: PetscSFType, PetscSFCreate()
11395fce210SBarry Smith @*/
11495fce210SBarry Smith PetscErrorCode PetscSFSetType(PetscSF sf,PetscSFType type)
11595fce210SBarry Smith {
11695fce210SBarry Smith   PetscErrorCode ierr,(*r)(PetscSF);
11795fce210SBarry Smith   PetscBool      match;
11895fce210SBarry Smith 
11995fce210SBarry Smith   PetscFunctionBegin;
12095fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
12195fce210SBarry Smith   PetscValidCharPointer(type,2);
12295fce210SBarry Smith 
12395fce210SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)sf,type,&match);CHKERRQ(ierr);
12495fce210SBarry Smith   if (match) PetscFunctionReturn(0);
12595fce210SBarry Smith 
126adc40e5bSBarry Smith   ierr = PetscFunctionListFind(PetscSFList,type,&r);CHKERRQ(ierr);
12795fce210SBarry Smith   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested PetscSF type %s",type);
12895fce210SBarry Smith   /* Destroy the previous private PetscSF context */
12995fce210SBarry Smith   if (sf->ops->Destroy) {
13095fce210SBarry Smith     ierr = (*(sf)->ops->Destroy)(sf);CHKERRQ(ierr);
13195fce210SBarry Smith   }
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 
138d36d33e4SMatthew G. Knepley /*@
13995fce210SBarry Smith    PetscSFDestroy - destroy star forest
14095fce210SBarry Smith 
14195fce210SBarry Smith    Collective
14295fce210SBarry Smith 
14395fce210SBarry Smith    Input Arguments:
14495fce210SBarry Smith .  sf - address of star forest
14595fce210SBarry Smith 
14695fce210SBarry Smith    Level: intermediate
14795fce210SBarry Smith 
14895fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFReset()
14995fce210SBarry Smith @*/
15095fce210SBarry Smith PetscErrorCode PetscSFDestroy(PetscSF *sf)
15195fce210SBarry Smith {
15295fce210SBarry Smith   PetscErrorCode ierr;
15395fce210SBarry Smith 
15495fce210SBarry Smith   PetscFunctionBegin;
15595fce210SBarry Smith   if (!*sf) PetscFunctionReturn(0);
15695fce210SBarry Smith   PetscValidHeaderSpecific((*sf),PETSCSF_CLASSID,1);
15795fce210SBarry Smith   if (--((PetscObject)(*sf))->refct > 0) {*sf = 0; PetscFunctionReturn(0);}
15895fce210SBarry Smith   ierr = PetscSFReset(*sf);CHKERRQ(ierr);
15995fce210SBarry Smith   if ((*sf)->ops->Destroy) {ierr = (*(*sf)->ops->Destroy)(*sf);CHKERRQ(ierr);}
16095fce210SBarry Smith   ierr = PetscHeaderDestroy(sf);CHKERRQ(ierr);
16195fce210SBarry Smith   PetscFunctionReturn(0);
16295fce210SBarry Smith }
16395fce210SBarry Smith 
16495fce210SBarry Smith /*@
16595fce210SBarry Smith    PetscSFSetUp - set up communication structures
16695fce210SBarry Smith 
16795fce210SBarry Smith    Collective
16895fce210SBarry Smith 
16995fce210SBarry Smith    Input Arguments:
17095fce210SBarry Smith .  sf - star forest communication object
17195fce210SBarry Smith 
17295fce210SBarry Smith    Level: beginner
17395fce210SBarry Smith 
17495fce210SBarry Smith .seealso: PetscSFSetFromOptions(), PetscSFSetType()
17595fce210SBarry Smith @*/
17695fce210SBarry Smith PetscErrorCode PetscSFSetUp(PetscSF sf)
17795fce210SBarry Smith {
17895fce210SBarry Smith   PetscErrorCode ierr;
17995fce210SBarry Smith 
18095fce210SBarry Smith   PetscFunctionBegin;
18195fce210SBarry Smith   if (sf->setupcalled) PetscFunctionReturn(0);
18295fce210SBarry Smith   if (!((PetscObject)sf)->type_name) {ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr);}
18395fce210SBarry Smith   if (sf->ops->SetUp) {ierr = (*sf->ops->SetUp)(sf);CHKERRQ(ierr);}
18495fce210SBarry Smith   sf->setupcalled = PETSC_TRUE;
18595fce210SBarry Smith   PetscFunctionReturn(0);
18695fce210SBarry Smith }
18795fce210SBarry Smith 
1888af6ec1cSBarry Smith /*@
18995fce210SBarry Smith    PetscSFSetFromOptions - set PetscSF options using the options database
19095fce210SBarry Smith 
19195fce210SBarry Smith    Logically Collective
19295fce210SBarry Smith 
19395fce210SBarry Smith    Input Arguments:
19495fce210SBarry Smith .  sf - star forest
19595fce210SBarry Smith 
19695fce210SBarry Smith    Options Database Keys:
19760263706SJed Brown +  -sf_type - implementation type, see PetscSFSetType()
19860263706SJed Brown -  -sf_rank_order - sort composite points for gathers and scatters in rank order, gathers are non-deterministic otherwise
19995fce210SBarry Smith 
20095fce210SBarry Smith    Level: intermediate
20195fce210SBarry Smith 
20295fce210SBarry Smith .keywords: KSP, set, from, options, database
20395fce210SBarry Smith 
20495fce210SBarry Smith .seealso: PetscSFWindowSetSyncType()
20595fce210SBarry Smith @*/
20695fce210SBarry Smith PetscErrorCode PetscSFSetFromOptions(PetscSF sf)
20795fce210SBarry Smith {
20895fce210SBarry Smith   PetscSFType    deft;
20995fce210SBarry Smith   char           type[256];
21095fce210SBarry Smith   PetscErrorCode ierr;
21195fce210SBarry Smith   PetscBool      flg;
21295fce210SBarry Smith 
21395fce210SBarry Smith   PetscFunctionBegin;
21495fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
21595fce210SBarry Smith   ierr = PetscObjectOptionsBegin((PetscObject)sf);CHKERRQ(ierr);
21695fce210SBarry Smith   deft = ((PetscObject)sf)->type_name ? ((PetscObject)sf)->type_name : PETSCSFBASIC;
217a264d7a6SBarry Smith   ierr = PetscOptionsFList("-sf_type","PetscSF implementation type","PetscSFSetType",PetscSFList,deft,type,256,&flg);CHKERRQ(ierr);
21895fce210SBarry Smith   ierr = PetscSFSetType(sf,flg ? type : deft);CHKERRQ(ierr);
21995fce210SBarry 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);
220e55864a3SBarry Smith   if (sf->ops->SetFromOptions) {ierr = (*sf->ops->SetFromOptions)(PetscOptionsObject,sf);CHKERRQ(ierr);}
22195fce210SBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
22295fce210SBarry Smith   PetscFunctionReturn(0);
22395fce210SBarry Smith }
22495fce210SBarry Smith 
22595fce210SBarry Smith /*@C
22695fce210SBarry Smith    PetscSFSetRankOrder - sort multi-points for gathers and scatters by rank order
22795fce210SBarry Smith 
22895fce210SBarry Smith    Logically Collective
22995fce210SBarry Smith 
23095fce210SBarry Smith    Input Arguments:
23195fce210SBarry Smith +  sf - star forest
23295fce210SBarry Smith -  flg - PETSC_TRUE to sort, PETSC_FALSE to skip sorting (lower setup cost, but non-deterministic)
23395fce210SBarry Smith 
23495fce210SBarry Smith    Level: advanced
23595fce210SBarry Smith 
23695fce210SBarry Smith .seealso: PetscSFGatherBegin(), PetscSFScatterBegin()
23795fce210SBarry Smith @*/
23895fce210SBarry Smith PetscErrorCode PetscSFSetRankOrder(PetscSF sf,PetscBool flg)
23995fce210SBarry Smith {
24095fce210SBarry Smith 
24195fce210SBarry Smith   PetscFunctionBegin;
24295fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
24395fce210SBarry Smith   PetscValidLogicalCollectiveBool(sf,flg,2);
24495fce210SBarry Smith   if (sf->multi) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Rank ordering must be set before first call to PetscSFGatherBegin() or PetscSFScatterBegin()");
24595fce210SBarry Smith   sf->rankorder = flg;
24695fce210SBarry Smith   PetscFunctionReturn(0);
24795fce210SBarry Smith }
24895fce210SBarry Smith 
2498af6ec1cSBarry Smith /*@
25095fce210SBarry Smith    PetscSFSetGraph - Set a parallel star forest
25195fce210SBarry Smith 
25295fce210SBarry Smith    Collective
25395fce210SBarry Smith 
25495fce210SBarry Smith    Input Arguments:
25595fce210SBarry Smith +  sf - star forest
25695fce210SBarry Smith .  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
25795fce210SBarry Smith .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
25895fce210SBarry Smith .  ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage
25995fce210SBarry Smith .  localmode - copy mode for ilocal
26095fce210SBarry Smith .  iremote - remote locations of root vertices for each leaf on the current process
26195fce210SBarry Smith -  remotemode - copy mode for iremote
26295fce210SBarry Smith 
26395fce210SBarry Smith    Level: intermediate
26495fce210SBarry Smith 
26538ab3f8aSBarry Smith    Notes: In Fortran you must use PETSC_COPY_VALUES for localmode and remotemode
26638ab3f8aSBarry Smith 
26795fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFGetGraph()
26895fce210SBarry Smith @*/
26995fce210SBarry Smith PetscErrorCode PetscSFSetGraph(PetscSF sf,PetscInt nroots,PetscInt nleaves,const PetscInt *ilocal,PetscCopyMode localmode,const PetscSFNode *iremote,PetscCopyMode remotemode)
27095fce210SBarry Smith {
27195fce210SBarry Smith   PetscErrorCode     ierr;
27295fce210SBarry Smith 
27395fce210SBarry Smith   PetscFunctionBegin;
27495fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
275563ffbabSMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr);
27695fce210SBarry Smith   if (nleaves && ilocal) PetscValidIntPointer(ilocal,4);
27795fce210SBarry Smith   if (nleaves) PetscValidPointer(iremote,6);
27895fce210SBarry Smith   if (nroots < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"roots %D, cannot be negative",nroots);
27995fce210SBarry Smith   if (nleaves < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nleaves %D, cannot be negative",nleaves);
28095fce210SBarry Smith   ierr        = PetscSFReset(sf);CHKERRQ(ierr);
28195fce210SBarry Smith   sf->nroots  = nroots;
28295fce210SBarry Smith   sf->nleaves = nleaves;
28395fce210SBarry Smith   if (ilocal) {
28421c688dcSJed Brown     PetscInt i;
28595fce210SBarry Smith     switch (localmode) {
28695fce210SBarry Smith     case PETSC_COPY_VALUES:
287785e854fSJed Brown       ierr        = PetscMalloc1(nleaves,&sf->mine_alloc);CHKERRQ(ierr);
28895fce210SBarry Smith       sf->mine    = sf->mine_alloc;
28995fce210SBarry Smith       ierr        = PetscMemcpy(sf->mine,ilocal,nleaves*sizeof(*sf->mine));CHKERRQ(ierr);
29095fce210SBarry Smith       sf->minleaf = PETSC_MAX_INT;
29195fce210SBarry Smith       sf->maxleaf = PETSC_MIN_INT;
29295fce210SBarry Smith       for (i=0; i<nleaves; i++) {
29395fce210SBarry Smith         sf->minleaf = PetscMin(sf->minleaf,ilocal[i]);
29495fce210SBarry Smith         sf->maxleaf = PetscMax(sf->maxleaf,ilocal[i]);
29595fce210SBarry Smith       }
29695fce210SBarry Smith       break;
29795fce210SBarry Smith     case PETSC_OWN_POINTER:
29895fce210SBarry Smith       sf->mine_alloc = (PetscInt*)ilocal;
29995fce210SBarry Smith       sf->mine       = sf->mine_alloc;
30095fce210SBarry Smith       break;
30195fce210SBarry Smith     case PETSC_USE_POINTER:
30295fce210SBarry Smith       sf->mine = (PetscInt*)ilocal;
30395fce210SBarry Smith       break;
30495fce210SBarry Smith     default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown localmode");
30595fce210SBarry Smith     }
30695fce210SBarry Smith   }
30795fce210SBarry Smith   if (!ilocal || nleaves > 0) {
30895fce210SBarry Smith     sf->minleaf = 0;
30995fce210SBarry Smith     sf->maxleaf = nleaves - 1;
31095fce210SBarry Smith   }
31195fce210SBarry Smith   switch (remotemode) {
31295fce210SBarry Smith   case PETSC_COPY_VALUES:
313785e854fSJed Brown     ierr       = PetscMalloc1(nleaves,&sf->remote_alloc);CHKERRQ(ierr);
31495fce210SBarry Smith     sf->remote = sf->remote_alloc;
31595fce210SBarry Smith     ierr       = PetscMemcpy(sf->remote,iremote,nleaves*sizeof(*sf->remote));CHKERRQ(ierr);
31695fce210SBarry Smith     break;
31795fce210SBarry Smith   case PETSC_OWN_POINTER:
31895fce210SBarry Smith     sf->remote_alloc = (PetscSFNode*)iremote;
31995fce210SBarry Smith     sf->remote       = sf->remote_alloc;
32095fce210SBarry Smith     break;
32195fce210SBarry Smith   case PETSC_USE_POINTER:
32295fce210SBarry Smith     sf->remote = (PetscSFNode*)iremote;
32395fce210SBarry Smith     break;
32495fce210SBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown remotemode");
32595fce210SBarry Smith   }
32695fce210SBarry Smith 
32795fce210SBarry Smith   sf->graphset = PETSC_TRUE;
328563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr);
32995fce210SBarry Smith   PetscFunctionReturn(0);
33095fce210SBarry Smith }
33195fce210SBarry Smith 
33295fce210SBarry Smith /*@C
33395fce210SBarry Smith    PetscSFCreateInverseSF - given a PetscSF in which all vertices have degree 1, creates the inverse map
33495fce210SBarry Smith 
33595fce210SBarry Smith    Collective
33695fce210SBarry Smith 
33795fce210SBarry Smith    Input Arguments:
33895fce210SBarry Smith .  sf - star forest to invert
33995fce210SBarry Smith 
34095fce210SBarry Smith    Output Arguments:
34195fce210SBarry Smith .  isf - inverse of sf
34295fce210SBarry Smith 
34395fce210SBarry Smith    Level: advanced
34495fce210SBarry Smith 
34595fce210SBarry Smith    Notes:
34695fce210SBarry Smith    All roots must have degree 1.
34795fce210SBarry Smith 
34895fce210SBarry Smith    The local space may be a permutation, but cannot be sparse.
34995fce210SBarry Smith 
35095fce210SBarry Smith .seealso: PetscSFSetGraph()
35195fce210SBarry Smith @*/
35295fce210SBarry Smith PetscErrorCode PetscSFCreateInverseSF(PetscSF sf,PetscSF *isf)
35395fce210SBarry Smith {
35495fce210SBarry Smith   PetscErrorCode ierr;
35595fce210SBarry Smith   PetscMPIInt    rank;
35695fce210SBarry Smith   PetscInt       i,nroots,nleaves,maxlocal,count,*newilocal;
35795fce210SBarry Smith   const PetscInt *ilocal;
35895fce210SBarry Smith   PetscSFNode    *roots,*leaves;
35995fce210SBarry Smith 
36095fce210SBarry Smith   PetscFunctionBegin;
36195fce210SBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
36295fce210SBarry Smith   ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);CHKERRQ(ierr);
36395fce210SBarry Smith   for (i=0,maxlocal=0; i<nleaves; i++) maxlocal = PetscMax(maxlocal,(ilocal ? ilocal[i] : i)+1);
364ae9aee6dSMatthew G. Knepley   ierr = PetscMalloc2(nroots,&roots,maxlocal,&leaves);CHKERRQ(ierr);
365ae9aee6dSMatthew G. Knepley   for (i=0; i<maxlocal; i++) {
36695fce210SBarry Smith     leaves[i].rank  = rank;
36795fce210SBarry Smith     leaves[i].index = i;
36895fce210SBarry Smith   }
36995fce210SBarry Smith   for (i=0; i <nroots; i++) {
37095fce210SBarry Smith     roots[i].rank  = -1;
37195fce210SBarry Smith     roots[i].index = -1;
37295fce210SBarry Smith   }
3738bfbc91cSJed Brown   ierr = PetscSFReduceBegin(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr);
3748bfbc91cSJed Brown   ierr = PetscSFReduceEnd(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr);
37595fce210SBarry Smith 
37695fce210SBarry Smith   /* Check whether our leaves are sparse */
37795fce210SBarry Smith   for (i=0,count=0; i<nroots; i++) if (roots[i].rank >= 0) count++;
37895fce210SBarry Smith   if (count == nroots) newilocal = NULL;
37995fce210SBarry Smith   else {                        /* Index for sparse leaves and compact "roots" array (which is to become our leaves). */
380785e854fSJed Brown     ierr = PetscMalloc1(count,&newilocal);CHKERRQ(ierr);
38195fce210SBarry Smith     for (i=0,count=0; i<nroots; i++) {
38295fce210SBarry Smith       if (roots[i].rank >= 0) {
38395fce210SBarry Smith         newilocal[count]   = i;
38495fce210SBarry Smith         roots[count].rank  = roots[i].rank;
38595fce210SBarry Smith         roots[count].index = roots[i].index;
38695fce210SBarry Smith         count++;
38795fce210SBarry Smith       }
38895fce210SBarry Smith     }
38995fce210SBarry Smith   }
39095fce210SBarry Smith 
39195fce210SBarry Smith   ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,isf);CHKERRQ(ierr);
39295fce210SBarry Smith   ierr = PetscSFSetGraph(*isf,maxlocal,count,newilocal,PETSC_OWN_POINTER,roots,PETSC_COPY_VALUES);CHKERRQ(ierr);
39395fce210SBarry Smith   ierr = PetscFree2(roots,leaves);CHKERRQ(ierr);
39495fce210SBarry Smith   PetscFunctionReturn(0);
39595fce210SBarry Smith }
39695fce210SBarry Smith 
39795fce210SBarry Smith /*@
39895fce210SBarry Smith    PetscSFDuplicate - duplicate a PetscSF, optionally preserving rank connectivity and graph
39995fce210SBarry Smith 
40095fce210SBarry Smith    Collective
40195fce210SBarry Smith 
40295fce210SBarry Smith    Input Arguments:
40395fce210SBarry Smith +  sf - communication object to duplicate
40495fce210SBarry Smith -  opt - PETSCSF_DUPLICATE_CONFONLY, PETSCSF_DUPLICATE_RANKS, or PETSCSF_DUPLICATE_GRAPH (see PetscSFDuplicateOption)
40595fce210SBarry Smith 
40695fce210SBarry Smith    Output Arguments:
40795fce210SBarry Smith .  newsf - new communication object
40895fce210SBarry Smith 
40995fce210SBarry Smith    Level: beginner
41095fce210SBarry Smith 
41195fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetType(), PetscSFSetGraph()
41295fce210SBarry Smith @*/
41395fce210SBarry Smith PetscErrorCode PetscSFDuplicate(PetscSF sf,PetscSFDuplicateOption opt,PetscSF *newsf)
41495fce210SBarry Smith {
41595fce210SBarry Smith   PetscErrorCode ierr;
41695fce210SBarry Smith 
41795fce210SBarry Smith   PetscFunctionBegin;
41895fce210SBarry Smith   ierr = PetscSFCreate(PetscObjectComm((PetscObject)sf),newsf);CHKERRQ(ierr);
41995fce210SBarry Smith   ierr = PetscSFSetType(*newsf,((PetscObject)sf)->type_name);CHKERRQ(ierr);
42095fce210SBarry Smith   if (sf->ops->Duplicate) {ierr = (*sf->ops->Duplicate)(sf,opt,*newsf);CHKERRQ(ierr);}
42195fce210SBarry Smith   if (opt == PETSCSF_DUPLICATE_GRAPH) {
42295fce210SBarry Smith     PetscInt          nroots,nleaves;
42395fce210SBarry Smith     const PetscInt    *ilocal;
42495fce210SBarry Smith     const PetscSFNode *iremote;
42595fce210SBarry Smith     ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
42695fce210SBarry Smith     ierr = PetscSFSetGraph(*newsf,nroots,nleaves,ilocal,PETSC_COPY_VALUES,iremote,PETSC_COPY_VALUES);CHKERRQ(ierr);
42795fce210SBarry Smith   }
42895fce210SBarry Smith   PetscFunctionReturn(0);
42995fce210SBarry Smith }
43095fce210SBarry Smith 
43195fce210SBarry Smith /*@C
43295fce210SBarry Smith    PetscSFGetGraph - Get the graph specifying a parallel star forest
43395fce210SBarry Smith 
43495fce210SBarry Smith    Not Collective
43595fce210SBarry Smith 
43695fce210SBarry Smith    Input Arguments:
43795fce210SBarry Smith .  sf - star forest
43895fce210SBarry Smith 
43995fce210SBarry Smith    Output Arguments:
44095fce210SBarry Smith +  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
44195fce210SBarry Smith .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
44295fce210SBarry Smith .  ilocal - locations of leaves in leafdata buffers
44395fce210SBarry Smith -  iremote - remote locations of root vertices for each leaf on the current process
44495fce210SBarry Smith 
44595fce210SBarry Smith    Level: intermediate
44695fce210SBarry Smith 
44795fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph()
44895fce210SBarry Smith @*/
44995fce210SBarry Smith PetscErrorCode PetscSFGetGraph(PetscSF sf,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote)
45095fce210SBarry Smith {
45195fce210SBarry Smith 
45295fce210SBarry Smith   PetscFunctionBegin;
45395fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
45495fce210SBarry Smith   /* We are not currently requiring that the graph is set, thus returning nroots=-1 if it has not been set */
45595fce210SBarry Smith   /* if (!sf->graphset) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Graph has not been set, must call PetscSFSetGraph()"); */
45695fce210SBarry Smith   if (nroots) *nroots = sf->nroots;
45795fce210SBarry Smith   if (nleaves) *nleaves = sf->nleaves;
45895fce210SBarry Smith   if (ilocal) *ilocal = sf->mine;
45995fce210SBarry Smith   if (iremote) *iremote = sf->remote;
46095fce210SBarry Smith   PetscFunctionReturn(0);
46195fce210SBarry Smith }
46295fce210SBarry Smith 
46395fce210SBarry Smith /*@C
46495fce210SBarry Smith    PetscSFGetLeafRange - Get the active leaf ranges
46595fce210SBarry Smith 
46695fce210SBarry Smith    Not Collective
46795fce210SBarry Smith 
46895fce210SBarry Smith    Input Arguments:
46995fce210SBarry Smith .  sf - star forest
47095fce210SBarry Smith 
47195fce210SBarry Smith    Output Arguments:
47295fce210SBarry Smith +  minleaf - minimum active leaf on this process
47395fce210SBarry Smith -  maxleaf - maximum active leaf on this process
47495fce210SBarry Smith 
47595fce210SBarry Smith    Level: developer
47695fce210SBarry Smith 
47795fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph(), PetscSFGetGraph()
47895fce210SBarry Smith @*/
47995fce210SBarry Smith PetscErrorCode PetscSFGetLeafRange(PetscSF sf,PetscInt *minleaf,PetscInt *maxleaf)
48095fce210SBarry Smith {
48195fce210SBarry Smith 
48295fce210SBarry Smith   PetscFunctionBegin;
48395fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
48495fce210SBarry Smith   if (minleaf) *minleaf = sf->minleaf;
48595fce210SBarry Smith   if (maxleaf) *maxleaf = sf->maxleaf;
48695fce210SBarry Smith   PetscFunctionReturn(0);
48795fce210SBarry Smith }
48895fce210SBarry Smith 
48995fce210SBarry Smith /*@C
49095fce210SBarry Smith    PetscSFView - view a star forest
49195fce210SBarry Smith 
49295fce210SBarry Smith    Collective
49395fce210SBarry Smith 
49495fce210SBarry Smith    Input Arguments:
49595fce210SBarry Smith +  sf - star forest
49695fce210SBarry Smith -  viewer - viewer to display graph, for example PETSC_VIEWER_STDOUT_WORLD
49795fce210SBarry Smith 
49895fce210SBarry Smith    Level: beginner
49995fce210SBarry Smith 
50095fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph()
50195fce210SBarry Smith @*/
50295fce210SBarry Smith PetscErrorCode PetscSFView(PetscSF sf,PetscViewer viewer)
50395fce210SBarry Smith {
50495fce210SBarry Smith   PetscErrorCode    ierr;
50595fce210SBarry Smith   PetscBool         iascii;
50695fce210SBarry Smith   PetscViewerFormat format;
50795fce210SBarry Smith 
50895fce210SBarry Smith   PetscFunctionBegin;
50995fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
51095fce210SBarry Smith   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sf),&viewer);CHKERRQ(ierr);}
51195fce210SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
51295fce210SBarry Smith   PetscCheckSameComm(sf,1,viewer,2);
513*80153354SVaclav Hapla   if (sf->graphset) {ierr = PetscSFSetUp(sf);CHKERRQ(ierr);}
51495fce210SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
51595fce210SBarry Smith   if (iascii) {
51695fce210SBarry Smith     PetscMPIInt rank;
51781bfa7aaSJed Brown     PetscInt    ii,i,j;
51895fce210SBarry Smith 
519dae58748SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)sf,viewer);CHKERRQ(ierr);
52095fce210SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
52195fce210SBarry Smith     if (sf->ops->View) {ierr = (*sf->ops->View)(sf,viewer);CHKERRQ(ierr);}
522*80153354SVaclav Hapla     if (!sf->graphset) {
523*80153354SVaclav Hapla       ierr = PetscViewerASCIIPrintf(viewer,"PetscSFSetGraph() has not been called yet\n");CHKERRQ(ierr);
524*80153354SVaclav Hapla       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
525*80153354SVaclav Hapla       PetscFunctionReturn(0);
526*80153354SVaclav Hapla     }
52795fce210SBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
5281575c14dSBarry Smith     ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
52995fce210SBarry Smith     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of roots=%D, leaves=%D, remote ranks=%D\n",rank,sf->nroots,sf->nleaves,sf->nranks);CHKERRQ(ierr);
53095fce210SBarry Smith     for (i=0; i<sf->nleaves; i++) {
53195fce210SBarry 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);
53295fce210SBarry Smith     }
53395fce210SBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
53495fce210SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
53595fce210SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
53681bfa7aaSJed Brown       PetscMPIInt *tmpranks,*perm;
53781bfa7aaSJed Brown       ierr = PetscMalloc2(sf->nranks,&tmpranks,sf->nranks,&perm);CHKERRQ(ierr);
53881bfa7aaSJed Brown       ierr = PetscMemcpy(tmpranks,sf->ranks,sf->nranks*sizeof(tmpranks[0]));CHKERRQ(ierr);
53981bfa7aaSJed Brown       for (i=0; i<sf->nranks; i++) perm[i] = i;
54081bfa7aaSJed Brown       ierr = PetscSortMPIIntWithArray(sf->nranks,tmpranks,perm);CHKERRQ(ierr);
54195fce210SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Roots referenced by my leaves, by rank\n",rank);CHKERRQ(ierr);
54281bfa7aaSJed Brown       for (ii=0; ii<sf->nranks; ii++) {
54381bfa7aaSJed Brown         i = perm[ii];
5447904a332SBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %d: %D edges\n",rank,sf->ranks[i],sf->roffset[i+1]-sf->roffset[i]);CHKERRQ(ierr);
54595fce210SBarry Smith         for (j=sf->roffset[i]; j<sf->roffset[i+1]; j++) {
54695fce210SBarry Smith           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d]    %D <- %D\n",rank,sf->rmine[j],sf->rremote[j]);CHKERRQ(ierr);
54795fce210SBarry Smith         }
54895fce210SBarry Smith       }
54981bfa7aaSJed Brown       ierr = PetscFree2(tmpranks,perm);CHKERRQ(ierr);
55095fce210SBarry Smith     }
55195fce210SBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
5521575c14dSBarry Smith     ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
55395fce210SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
55495fce210SBarry Smith   }
55595fce210SBarry Smith   PetscFunctionReturn(0);
55695fce210SBarry Smith }
55795fce210SBarry Smith 
55895fce210SBarry Smith /*@C
55995fce210SBarry Smith    PetscSFGetRanks - Get ranks and number of vertices referenced by leaves on this process
56095fce210SBarry Smith 
56195fce210SBarry Smith    Not Collective
56295fce210SBarry Smith 
56395fce210SBarry Smith    Input Arguments:
56495fce210SBarry Smith .  sf - star forest
56595fce210SBarry Smith 
56695fce210SBarry Smith    Output Arguments:
56795fce210SBarry Smith +  nranks - number of ranks referenced by local part
56895fce210SBarry Smith .  ranks - array of ranks
56995fce210SBarry Smith .  roffset - offset in rmine/rremote for each rank (length nranks+1)
57095fce210SBarry Smith .  rmine - concatenated array holding local indices referencing each remote rank
57195fce210SBarry Smith -  rremote - concatenated array holding remote indices referenced for each remote rank
57295fce210SBarry Smith 
57395fce210SBarry Smith    Level: developer
57495fce210SBarry Smith 
57595fce210SBarry Smith .seealso: PetscSFSetGraph()
57695fce210SBarry Smith @*/
57795fce210SBarry Smith PetscErrorCode PetscSFGetRanks(PetscSF sf,PetscInt *nranks,const PetscMPIInt **ranks,const PetscInt **roffset,const PetscInt **rmine,const PetscInt **rremote)
57895fce210SBarry Smith {
57995fce210SBarry Smith 
58095fce210SBarry Smith   PetscFunctionBegin;
58195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
582dcca2fcaSJed Brown   if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp before obtaining ranks");
58395fce210SBarry Smith   if (nranks)  *nranks  = sf->nranks;
58495fce210SBarry Smith   if (ranks)   *ranks   = sf->ranks;
58595fce210SBarry Smith   if (roffset) *roffset = sf->roffset;
58695fce210SBarry Smith   if (rmine)   *rmine   = sf->rmine;
58795fce210SBarry Smith   if (rremote) *rremote = sf->rremote;
58895fce210SBarry Smith   PetscFunctionReturn(0);
58995fce210SBarry Smith }
59095fce210SBarry Smith 
591b5a8e515SJed Brown static PetscBool InList(PetscMPIInt needle,PetscMPIInt n,const PetscMPIInt *list) {
592b5a8e515SJed Brown   PetscInt i;
593b5a8e515SJed Brown   for (i=0; i<n; i++) {
594b5a8e515SJed Brown     if (needle == list[i]) return PETSC_TRUE;
595b5a8e515SJed Brown   }
596b5a8e515SJed Brown   return PETSC_FALSE;
597b5a8e515SJed Brown }
598b5a8e515SJed Brown 
59995fce210SBarry Smith /*@C
60021c688dcSJed Brown    PetscSFSetUpRanks - Set up data structures associated with ranks; this is for internal use by PetscSF implementations.
60121c688dcSJed Brown 
60221c688dcSJed Brown    Collective
60321c688dcSJed Brown 
60421c688dcSJed Brown    Input Arguments:
605b5a8e515SJed Brown +  sf - PetscSF to set up; PetscSFSetGraph() must have been called
606b5a8e515SJed Brown -  dgroup - MPI_Group of ranks to be distinguished (e.g., for self or shared memory exchange)
60721c688dcSJed Brown 
60821c688dcSJed Brown    Level: developer
60921c688dcSJed Brown 
61021c688dcSJed Brown .seealso: PetscSFGetRanks()
61121c688dcSJed Brown @*/
612b5a8e515SJed Brown PetscErrorCode PetscSFSetUpRanks(PetscSF sf,MPI_Group dgroup)
61321c688dcSJed Brown {
61421c688dcSJed Brown   PetscErrorCode     ierr;
61521c688dcSJed Brown   PetscTable         table;
61621c688dcSJed Brown   PetscTablePosition pos;
617b5a8e515SJed Brown   PetscMPIInt        size,groupsize,*groupranks;
61821c688dcSJed Brown   PetscInt           i,*rcount,*ranks;
61921c688dcSJed Brown 
62021c688dcSJed Brown   PetscFunctionBegin;
62121c688dcSJed Brown   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
62221c688dcSJed Brown   if (!sf->graphset) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"PetscSFSetGraph() has not been called yet");
62321c688dcSJed Brown   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr);
62421c688dcSJed Brown   ierr = PetscTableCreate(10,size,&table);CHKERRQ(ierr);
62521c688dcSJed Brown   for (i=0; i<sf->nleaves; i++) {
62621c688dcSJed Brown     /* Log 1-based rank */
62721c688dcSJed Brown     ierr = PetscTableAdd(table,sf->remote[i].rank+1,1,ADD_VALUES);CHKERRQ(ierr);
62821c688dcSJed Brown   }
62921c688dcSJed Brown   ierr = PetscTableGetCount(table,&sf->nranks);CHKERRQ(ierr);
63021c688dcSJed Brown   ierr = PetscMalloc4(sf->nranks,&sf->ranks,sf->nranks+1,&sf->roffset,sf->nleaves,&sf->rmine,sf->nleaves,&sf->rremote);CHKERRQ(ierr);
63121c688dcSJed Brown   ierr = PetscMalloc2(sf->nranks,&rcount,sf->nranks,&ranks);CHKERRQ(ierr);
63221c688dcSJed Brown   ierr = PetscTableGetHeadPosition(table,&pos);CHKERRQ(ierr);
63321c688dcSJed Brown   for (i=0; i<sf->nranks; i++) {
63421c688dcSJed Brown     ierr = PetscTableGetNext(table,&pos,&ranks[i],&rcount[i]);CHKERRQ(ierr);
63521c688dcSJed Brown     ranks[i]--;             /* Convert back to 0-based */
63621c688dcSJed Brown   }
63721c688dcSJed Brown   ierr = PetscTableDestroy(&table);CHKERRQ(ierr);
638b5a8e515SJed Brown 
639b5a8e515SJed Brown   /* We expect that dgroup is reliably "small" while nranks could be large */
640b5a8e515SJed Brown   {
6417fb8a5e4SKarl Rupp     MPI_Group group = MPI_GROUP_NULL;
642b5a8e515SJed Brown     PetscMPIInt *dgroupranks;
643b5a8e515SJed Brown     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
644b5a8e515SJed Brown     ierr = MPI_Group_size(dgroup,&groupsize);CHKERRQ(ierr);
645b5a8e515SJed Brown     ierr = PetscMalloc1(groupsize,&dgroupranks);CHKERRQ(ierr);
646b5a8e515SJed Brown     ierr = PetscMalloc1(groupsize,&groupranks);CHKERRQ(ierr);
647b5a8e515SJed Brown     for (i=0; i<groupsize; i++) dgroupranks[i] = i;
648f3fc9a17SJed Brown     if (groupsize) {ierr = MPI_Group_translate_ranks(dgroup,groupsize,dgroupranks,group,groupranks);CHKERRQ(ierr);}
649b5a8e515SJed Brown     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
650b5a8e515SJed Brown     ierr = PetscFree(dgroupranks);CHKERRQ(ierr);
651b5a8e515SJed Brown   }
652b5a8e515SJed Brown 
653b5a8e515SJed Brown   /* Partition ranks[] into distinguished (first sf->ndranks) followed by non-distinguished */
654b5a8e515SJed Brown   for (sf->ndranks=0,i=sf->nranks; sf->ndranks<i; ) {
655b5a8e515SJed Brown     for (i--; sf->ndranks<i; i--) { /* Scan i backward looking for distinguished rank */
656b5a8e515SJed Brown       if (InList(ranks[i],groupsize,groupranks)) break;
657b5a8e515SJed Brown     }
658b5a8e515SJed Brown     for ( ; sf->ndranks<=i; sf->ndranks++) { /* Scan sf->ndranks forward looking for non-distinguished rank */
659b5a8e515SJed Brown       if (!InList(ranks[sf->ndranks],groupsize,groupranks)) break;
660b5a8e515SJed Brown     }
661b5a8e515SJed Brown     if (sf->ndranks < i) {                         /* Swap ranks[sf->ndranks] with ranks[i] */
662b5a8e515SJed Brown       PetscInt    tmprank,tmpcount;
663b5a8e515SJed Brown       tmprank = ranks[i];
664b5a8e515SJed Brown       tmpcount = rcount[i];
665b5a8e515SJed Brown       ranks[i] = ranks[sf->ndranks];
666b5a8e515SJed Brown       rcount[i] = rcount[sf->ndranks];
667b5a8e515SJed Brown       ranks[sf->ndranks] = tmprank;
668b5a8e515SJed Brown       rcount[sf->ndranks] = tmpcount;
669b5a8e515SJed Brown       sf->ndranks++;
670b5a8e515SJed Brown     }
671b5a8e515SJed Brown   }
672b5a8e515SJed Brown   ierr = PetscFree(groupranks);CHKERRQ(ierr);
673b5a8e515SJed Brown   ierr = PetscSortIntWithArray(sf->ndranks,ranks,rcount);CHKERRQ(ierr);
674b5a8e515SJed Brown   ierr = PetscSortIntWithArray(sf->nranks-sf->ndranks,ranks+sf->ndranks,rcount+sf->ndranks);CHKERRQ(ierr);
67521c688dcSJed Brown   sf->roffset[0] = 0;
67621c688dcSJed Brown   for (i=0; i<sf->nranks; i++) {
67721c688dcSJed Brown     ierr = PetscMPIIntCast(ranks[i],sf->ranks+i);CHKERRQ(ierr);
67821c688dcSJed Brown     sf->roffset[i+1] = sf->roffset[i] + rcount[i];
67921c688dcSJed Brown     rcount[i]        = 0;
68021c688dcSJed Brown   }
68121c688dcSJed Brown   for (i=0; i<sf->nleaves; i++) {
682b5a8e515SJed Brown     PetscInt irank;
68321c688dcSJed Brown     /* Search for index of iremote[i].rank in sf->ranks */
684b5a8e515SJed Brown     ierr = PetscFindMPIInt(sf->remote[i].rank,sf->ndranks,sf->ranks,&irank);CHKERRQ(ierr);
685b5a8e515SJed Brown     if (irank < 0) {
686b5a8e515SJed Brown       ierr = PetscFindMPIInt(sf->remote[i].rank,sf->nranks-sf->ndranks,sf->ranks+sf->ndranks,&irank);CHKERRQ(ierr);
687b5a8e515SJed Brown       if (irank >= 0) irank += sf->ndranks;
68821c688dcSJed Brown     }
689b5a8e515SJed Brown     if (irank < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Could not find rank %D in array",sf->remote[i].rank);
69021c688dcSJed Brown     sf->rmine[sf->roffset[irank] + rcount[irank]]   = sf->mine ? sf->mine[i] : i;
69121c688dcSJed Brown     sf->rremote[sf->roffset[irank] + rcount[irank]] = sf->remote[i].index;
69221c688dcSJed Brown     rcount[irank]++;
69321c688dcSJed Brown   }
69421c688dcSJed Brown   ierr = PetscFree2(rcount,ranks);CHKERRQ(ierr);
69521c688dcSJed Brown   PetscFunctionReturn(0);
69621c688dcSJed Brown }
69721c688dcSJed Brown 
69821c688dcSJed Brown /*@C
69995fce210SBarry Smith    PetscSFGetGroups - gets incoming and outgoing process groups
70095fce210SBarry Smith 
70195fce210SBarry Smith    Collective
70295fce210SBarry Smith 
70395fce210SBarry Smith    Input Argument:
70495fce210SBarry Smith .  sf - star forest
70595fce210SBarry Smith 
70695fce210SBarry Smith    Output Arguments:
70795fce210SBarry Smith +  incoming - group of origin processes for incoming edges (leaves that reference my roots)
70895fce210SBarry Smith -  outgoing - group of destination processes for outgoing edges (roots that I reference)
70995fce210SBarry Smith 
71095fce210SBarry Smith    Level: developer
71195fce210SBarry Smith 
71295fce210SBarry Smith .seealso: PetscSFGetWindow(), PetscSFRestoreWindow()
71395fce210SBarry Smith @*/
71495fce210SBarry Smith PetscErrorCode PetscSFGetGroups(PetscSF sf,MPI_Group *incoming,MPI_Group *outgoing)
71595fce210SBarry Smith {
71695fce210SBarry Smith   PetscErrorCode ierr;
7177fb8a5e4SKarl Rupp   MPI_Group      group = MPI_GROUP_NULL;
71895fce210SBarry Smith 
71995fce210SBarry Smith   PetscFunctionBegin;
72095fce210SBarry Smith   if (sf->ingroup == MPI_GROUP_NULL) {
72195fce210SBarry Smith     PetscInt       i;
72295fce210SBarry Smith     const PetscInt *indegree;
72395fce210SBarry Smith     PetscMPIInt    rank,*outranks,*inranks;
72495fce210SBarry Smith     PetscSFNode    *remote;
72595fce210SBarry Smith     PetscSF        bgcount;
72695fce210SBarry Smith 
72795fce210SBarry Smith     /* Compute the number of incoming ranks */
728785e854fSJed Brown     ierr = PetscMalloc1(sf->nranks,&remote);CHKERRQ(ierr);
72995fce210SBarry Smith     for (i=0; i<sf->nranks; i++) {
73095fce210SBarry Smith       remote[i].rank  = sf->ranks[i];
73195fce210SBarry Smith       remote[i].index = 0;
73295fce210SBarry Smith     }
73395fce210SBarry Smith     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,&bgcount);CHKERRQ(ierr);
73495fce210SBarry Smith     ierr = PetscSFSetGraph(bgcount,1,sf->nranks,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
73595fce210SBarry Smith     ierr = PetscSFComputeDegreeBegin(bgcount,&indegree);CHKERRQ(ierr);
73695fce210SBarry Smith     ierr = PetscSFComputeDegreeEnd(bgcount,&indegree);CHKERRQ(ierr);
73795fce210SBarry Smith 
73895fce210SBarry Smith     /* Enumerate the incoming ranks */
739dcca6d9dSJed Brown     ierr = PetscMalloc2(indegree[0],&inranks,sf->nranks,&outranks);CHKERRQ(ierr);
74095fce210SBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
74195fce210SBarry Smith     for (i=0; i<sf->nranks; i++) outranks[i] = rank;
74295fce210SBarry Smith     ierr = PetscSFGatherBegin(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr);
74395fce210SBarry Smith     ierr = PetscSFGatherEnd(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr);
74495fce210SBarry Smith     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
74595fce210SBarry Smith     ierr = MPI_Group_incl(group,indegree[0],inranks,&sf->ingroup);CHKERRQ(ierr);
74695fce210SBarry Smith     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
74795fce210SBarry Smith     ierr = PetscFree2(inranks,outranks);CHKERRQ(ierr);
74895fce210SBarry Smith     ierr = PetscSFDestroy(&bgcount);CHKERRQ(ierr);
74995fce210SBarry Smith   }
75095fce210SBarry Smith   *incoming = sf->ingroup;
75195fce210SBarry Smith 
75295fce210SBarry Smith   if (sf->outgroup == MPI_GROUP_NULL) {
75395fce210SBarry Smith     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
75495fce210SBarry Smith     ierr = MPI_Group_incl(group,sf->nranks,sf->ranks,&sf->outgroup);CHKERRQ(ierr);
75595fce210SBarry Smith     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
75695fce210SBarry Smith   }
75795fce210SBarry Smith   *outgoing = sf->outgroup;
75895fce210SBarry Smith   PetscFunctionReturn(0);
75995fce210SBarry Smith }
76095fce210SBarry Smith 
76195fce210SBarry Smith /*@C
76295fce210SBarry Smith    PetscSFGetMultiSF - gets the inner SF implemeting gathers and scatters
76395fce210SBarry Smith 
76495fce210SBarry Smith    Collective
76595fce210SBarry Smith 
76695fce210SBarry Smith    Input Argument:
76795fce210SBarry Smith .  sf - star forest that may contain roots with 0 or with more than 1 vertex
76895fce210SBarry Smith 
76995fce210SBarry Smith    Output Arguments:
77095fce210SBarry Smith .  multi - star forest with split roots, such that each root has degree exactly 1
77195fce210SBarry Smith 
77295fce210SBarry Smith    Level: developer
77395fce210SBarry Smith 
77495fce210SBarry Smith    Notes:
77595fce210SBarry Smith 
77695fce210SBarry Smith    In most cases, users should use PetscSFGatherBegin() and PetscSFScatterBegin() instead of manipulating multi
77795fce210SBarry Smith    directly. Since multi satisfies the stronger condition that each entry in the global space has exactly one incoming
77895fce210SBarry Smith    edge, it is a candidate for future optimization that might involve its removal.
77995fce210SBarry Smith 
78095fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFGatherBegin(), PetscSFScatterBegin()
78195fce210SBarry Smith @*/
78295fce210SBarry Smith PetscErrorCode PetscSFGetMultiSF(PetscSF sf,PetscSF *multi)
78395fce210SBarry Smith {
78495fce210SBarry Smith   PetscErrorCode ierr;
78595fce210SBarry Smith 
78695fce210SBarry Smith   PetscFunctionBegin;
78795fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
78895fce210SBarry Smith   PetscValidPointer(multi,2);
78995fce210SBarry Smith   if (sf->nroots < 0) {         /* Graph has not been set yet; why do we need this? */
79095fce210SBarry Smith     ierr   = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr);
79195fce210SBarry Smith     *multi = sf->multi;
79295fce210SBarry Smith     PetscFunctionReturn(0);
79395fce210SBarry Smith   }
79495fce210SBarry Smith   if (!sf->multi) {
79595fce210SBarry Smith     const PetscInt *indegree;
7969837ea96SMatthew G. Knepley     PetscInt       i,*inoffset,*outones,*outoffset,maxlocal;
79795fce210SBarry Smith     PetscSFNode    *remote;
79895fce210SBarry Smith     ierr        = PetscSFComputeDegreeBegin(sf,&indegree);CHKERRQ(ierr);
79995fce210SBarry Smith     ierr        = PetscSFComputeDegreeEnd(sf,&indegree);CHKERRQ(ierr);
8009837ea96SMatthew G. Knepley     for (i=0,maxlocal=0; i<sf->nleaves; i++) maxlocal = PetscMax(maxlocal,(sf->mine ? sf->mine[i] : i)+1);
8019837ea96SMatthew G. Knepley     ierr        = PetscMalloc3(sf->nroots+1,&inoffset,maxlocal,&outones,maxlocal,&outoffset);CHKERRQ(ierr);
80295fce210SBarry Smith     inoffset[0] = 0;
80395fce210SBarry Smith     for (i=0; i<sf->nroots; i++) inoffset[i+1] = inoffset[i] + indegree[i];
8049837ea96SMatthew G. Knepley     for (i=0; i<maxlocal; i++) outones[i] = 1;
805dbd2ff41SBarry Smith     ierr = PetscSFFetchAndOpBegin(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr);
806dbd2ff41SBarry Smith     ierr = PetscSFFetchAndOpEnd(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr);
80795fce210SBarry Smith     for (i=0; i<sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */
80895fce210SBarry Smith #if 0
80995fce210SBarry Smith #if defined(PETSC_USE_DEBUG)                                 /* Check that the expected number of increments occurred */
81095fce210SBarry Smith     for (i=0; i<sf->nroots; i++) {
81195fce210SBarry Smith       if (inoffset[i] + indegree[i] != inoffset[i+1]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect result after PetscSFFetchAndOp");
81295fce210SBarry Smith     }
81395fce210SBarry Smith #endif
81495fce210SBarry Smith #endif
815785e854fSJed Brown     ierr = PetscMalloc1(sf->nleaves,&remote);CHKERRQ(ierr);
81695fce210SBarry Smith     for (i=0; i<sf->nleaves; i++) {
81795fce210SBarry Smith       remote[i].rank  = sf->remote[i].rank;
81838e7336fSToby Isaac       remote[i].index = outoffset[sf->mine ? sf->mine[i] : i];
81995fce210SBarry Smith     }
82095fce210SBarry Smith     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr);
82101365b40SToby Isaac     ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
82295fce210SBarry Smith     if (sf->rankorder) {        /* Sort the ranks */
82395fce210SBarry Smith       PetscMPIInt rank;
82495fce210SBarry Smith       PetscInt    *inranks,*newoffset,*outranks,*newoutoffset,*tmpoffset,maxdegree;
82595fce210SBarry Smith       PetscSFNode *newremote;
82695fce210SBarry Smith       ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
82795fce210SBarry Smith       for (i=0,maxdegree=0; i<sf->nroots; i++) maxdegree = PetscMax(maxdegree,indegree[i]);
8289837ea96SMatthew G. Knepley       ierr = PetscMalloc5(sf->multi->nroots,&inranks,sf->multi->nroots,&newoffset,maxlocal,&outranks,maxlocal,&newoutoffset,maxdegree,&tmpoffset);CHKERRQ(ierr);
8299837ea96SMatthew G. Knepley       for (i=0; i<maxlocal; i++) outranks[i] = rank;
8308bfbc91cSJed Brown       ierr = PetscSFReduceBegin(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr);
8318bfbc91cSJed Brown       ierr = PetscSFReduceEnd(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr);
83295fce210SBarry Smith       /* Sort the incoming ranks at each vertex, build the inverse map */
83395fce210SBarry Smith       for (i=0; i<sf->nroots; i++) {
83495fce210SBarry Smith         PetscInt j;
83595fce210SBarry Smith         for (j=0; j<indegree[i]; j++) tmpoffset[j] = j;
83695fce210SBarry Smith         ierr = PetscSortIntWithArray(indegree[i],inranks+inoffset[i],tmpoffset);CHKERRQ(ierr);
83795fce210SBarry Smith         for (j=0; j<indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j;
83895fce210SBarry Smith       }
83995fce210SBarry Smith       ierr = PetscSFBcastBegin(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr);
84095fce210SBarry Smith       ierr = PetscSFBcastEnd(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr);
841785e854fSJed Brown       ierr = PetscMalloc1(sf->nleaves,&newremote);CHKERRQ(ierr);
84295fce210SBarry Smith       for (i=0; i<sf->nleaves; i++) {
84395fce210SBarry Smith         newremote[i].rank  = sf->remote[i].rank;
84401365b40SToby Isaac         newremote[i].index = newoutoffset[sf->mine ? sf->mine[i] : i];
84595fce210SBarry Smith       }
84601365b40SToby Isaac       ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,newremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
84795fce210SBarry Smith       ierr = PetscFree5(inranks,newoffset,outranks,newoutoffset,tmpoffset);CHKERRQ(ierr);
84895fce210SBarry Smith     }
84995fce210SBarry Smith     ierr = PetscFree3(inoffset,outones,outoffset);CHKERRQ(ierr);
85095fce210SBarry Smith   }
85195fce210SBarry Smith   *multi = sf->multi;
85295fce210SBarry Smith   PetscFunctionReturn(0);
85395fce210SBarry Smith }
85495fce210SBarry Smith 
85595fce210SBarry Smith /*@C
85695fce210SBarry Smith    PetscSFCreateEmbeddedSF - removes edges from all but the selected roots, does not remap indices
85795fce210SBarry Smith 
85895fce210SBarry Smith    Collective
85995fce210SBarry Smith 
86095fce210SBarry Smith    Input Arguments:
86195fce210SBarry Smith +  sf - original star forest
86295fce210SBarry Smith .  nroots - number of roots to select on this process
86395fce210SBarry Smith -  selected - selected roots on this process
86495fce210SBarry Smith 
86595fce210SBarry Smith    Output Arguments:
86695fce210SBarry Smith .  newsf - new star forest
86795fce210SBarry Smith 
86895fce210SBarry Smith    Level: advanced
86995fce210SBarry Smith 
87095fce210SBarry Smith    Note:
87195fce210SBarry Smith    To use the new PetscSF, it may be necessary to know the indices of the leaves that are still participating. This can
87295fce210SBarry Smith    be done by calling PetscSFGetGraph().
87395fce210SBarry Smith 
87495fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFGetGraph()
87595fce210SBarry Smith @*/
87695fce210SBarry Smith PetscErrorCode PetscSFCreateEmbeddedSF(PetscSF sf,PetscInt nroots,const PetscInt *selected,PetscSF *newsf)
87795fce210SBarry Smith {
8780511a646SMatthew G. Knepley   PetscInt      *rootdata, *leafdata, *ilocal;
87995fce210SBarry Smith   PetscSFNode   *iremote;
880fc1ede2bSMatthew G. Knepley   PetscInt       leafsize = 0, nleaves = 0, n, i;
8810511a646SMatthew G. Knepley   PetscErrorCode ierr;
88295fce210SBarry Smith 
88395fce210SBarry Smith   PetscFunctionBegin;
88495fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
88595fce210SBarry Smith   if (nroots) PetscValidPointer(selected,3);
88695fce210SBarry Smith   PetscValidPointer(newsf,4);
8870511a646SMatthew G. Knepley   if (sf->mine) for (i = 0; i < sf->nleaves; ++i) {leafsize = PetscMax(leafsize, sf->mine[i]+1);}
8880511a646SMatthew G. Knepley   else leafsize = sf->nleaves;
8891795a4d1SJed Brown   ierr = PetscCalloc2(sf->nroots,&rootdata,leafsize,&leafdata);CHKERRQ(ierr);
8900511a646SMatthew G. Knepley   for (i=0; i<nroots; ++i) rootdata[selected[i]] = 1;
89195fce210SBarry Smith   ierr = PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr);
89295fce210SBarry Smith   ierr = PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr);
89395fce210SBarry Smith 
8940511a646SMatthew G. Knepley   for (i = 0; i < leafsize; ++i) nleaves += leafdata[i];
895785e854fSJed Brown   ierr = PetscMalloc1(nleaves,&ilocal);CHKERRQ(ierr);
896785e854fSJed Brown   ierr = PetscMalloc1(nleaves,&iremote);CHKERRQ(ierr);
897fc1ede2bSMatthew G. Knepley   for (i = 0, n = 0; i < sf->nleaves; ++i) {
8980511a646SMatthew G. Knepley     const PetscInt lidx = sf->mine ? sf->mine[i] : i;
8990511a646SMatthew G. Knepley 
9000511a646SMatthew G. Knepley     if (leafdata[lidx]) {
901fc1ede2bSMatthew G. Knepley       ilocal[n]        = lidx;
902fc1ede2bSMatthew G. Knepley       iremote[n].rank  = sf->remote[i].rank;
903fc1ede2bSMatthew G. Knepley       iremote[n].index = sf->remote[i].index;
904fc1ede2bSMatthew G. Knepley       ++n;
90595fce210SBarry Smith     }
90695fce210SBarry Smith   }
907fc1ede2bSMatthew G. Knepley   if (n != nleaves) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "There is a size mismatch in the SF embedding, %d != %d", n, nleaves);
90895fce210SBarry Smith   ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,newsf);CHKERRQ(ierr);
90995fce210SBarry Smith   ierr = PetscSFSetGraph(*newsf,sf->nroots,nleaves,ilocal,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
91095fce210SBarry Smith   ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr);
91195fce210SBarry Smith   PetscFunctionReturn(0);
91295fce210SBarry Smith }
91395fce210SBarry Smith 
9142f5fb4c2SMatthew G. Knepley /*@C
9152f5fb4c2SMatthew G. Knepley   PetscSFCreateEmbeddedLeafSF - removes edges from all but the selected leaves, does not remap indices
9162f5fb4c2SMatthew G. Knepley 
9172f5fb4c2SMatthew G. Knepley   Collective
9182f5fb4c2SMatthew G. Knepley 
9192f5fb4c2SMatthew G. Knepley   Input Arguments:
9202f5fb4c2SMatthew G. Knepley + sf - original star forest
9212f5fb4c2SMatthew G. Knepley . nleaves - number of leaves to select on this process
9222f5fb4c2SMatthew G. Knepley - selected - selected leaves on this process
9232f5fb4c2SMatthew G. Knepley 
9242f5fb4c2SMatthew G. Knepley   Output Arguments:
9252f5fb4c2SMatthew G. Knepley .  newsf - new star forest
9262f5fb4c2SMatthew G. Knepley 
9272f5fb4c2SMatthew G. Knepley   Level: advanced
9282f5fb4c2SMatthew G. Knepley 
9292f5fb4c2SMatthew G. Knepley .seealso: PetscSFCreateEmbeddedSF(), PetscSFSetGraph(), PetscSFGetGraph()
9302f5fb4c2SMatthew G. Knepley @*/
9312f5fb4c2SMatthew G. Knepley PetscErrorCode PetscSFCreateEmbeddedLeafSF(PetscSF sf, PetscInt nleaves, const PetscInt *selected, PetscSF *newsf)
9322f5fb4c2SMatthew G. Knepley {
9332f5fb4c2SMatthew G. Knepley   PetscSFNode   *iremote;
9342f5fb4c2SMatthew G. Knepley   PetscInt      *ilocal;
9352f5fb4c2SMatthew G. Knepley   PetscInt       i;
9362f5fb4c2SMatthew G. Knepley   PetscErrorCode ierr;
9372f5fb4c2SMatthew G. Knepley 
9382f5fb4c2SMatthew G. Knepley   PetscFunctionBegin;
9392f5fb4c2SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
9402f5fb4c2SMatthew G. Knepley   if (nleaves) PetscValidPointer(selected, 3);
9412f5fb4c2SMatthew G. Knepley   PetscValidPointer(newsf, 4);
9422f5fb4c2SMatthew G. Knepley   ierr = PetscMalloc1(nleaves, &ilocal);CHKERRQ(ierr);
9432f5fb4c2SMatthew G. Knepley   ierr = PetscMalloc1(nleaves, &iremote);CHKERRQ(ierr);
9442f5fb4c2SMatthew G. Knepley   for (i = 0; i < nleaves; ++i) {
9452f5fb4c2SMatthew G. Knepley     const PetscInt l = selected[i];
9462f5fb4c2SMatthew G. Knepley 
9472f5fb4c2SMatthew G. Knepley     ilocal[i]        = sf->mine ? sf->mine[l] : l;
9482f5fb4c2SMatthew G. Knepley     iremote[i].rank  = sf->remote[l].rank;
9492f5fb4c2SMatthew G. Knepley     iremote[i].index = sf->remote[l].index;
9502f5fb4c2SMatthew G. Knepley   }
9512f5fb4c2SMatthew G. Knepley   ierr = PetscSFDuplicate(sf, PETSCSF_DUPLICATE_RANKS, newsf);CHKERRQ(ierr);
9522f5fb4c2SMatthew G. Knepley   ierr = PetscSFSetGraph(*newsf, sf->nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr);
9532f5fb4c2SMatthew G. Knepley   PetscFunctionReturn(0);
9542f5fb4c2SMatthew G. Knepley }
9552f5fb4c2SMatthew G. Knepley 
95695fce210SBarry Smith /*@C
95795fce210SBarry Smith    PetscSFBcastBegin - begin pointwise broadcast to be concluded with call to PetscSFBcastEnd()
95895fce210SBarry Smith 
95995fce210SBarry Smith    Collective on PetscSF
96095fce210SBarry Smith 
96195fce210SBarry Smith    Input Arguments:
96295fce210SBarry Smith +  sf - star forest on which to communicate
96395fce210SBarry Smith .  unit - data type associated with each node
96495fce210SBarry Smith -  rootdata - buffer to broadcast
96595fce210SBarry Smith 
96695fce210SBarry Smith    Output Arguments:
96795fce210SBarry Smith .  leafdata - buffer to update with values from each leaf's respective root
96895fce210SBarry Smith 
96995fce210SBarry Smith    Level: intermediate
97095fce210SBarry Smith 
97195fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFView(), PetscSFBcastEnd(), PetscSFReduceBegin()
97295fce210SBarry Smith @*/
97395fce210SBarry Smith PetscErrorCode PetscSFBcastBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
97495fce210SBarry Smith {
97595fce210SBarry Smith   PetscErrorCode ierr;
97695fce210SBarry Smith 
97795fce210SBarry Smith   PetscFunctionBegin;
97895fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
97995fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
980563ffbabSMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr);
98195fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
98295fce210SBarry Smith   ierr = (*sf->ops->BcastBegin)(sf,unit,rootdata,leafdata);CHKERRQ(ierr);
983563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr);
98495fce210SBarry Smith   PetscFunctionReturn(0);
98595fce210SBarry Smith }
98695fce210SBarry Smith 
98795fce210SBarry Smith /*@C
98895fce210SBarry Smith    PetscSFBcastEnd - end a broadcast operation started with PetscSFBcastBegin()
98995fce210SBarry Smith 
99095fce210SBarry Smith    Collective
99195fce210SBarry Smith 
99295fce210SBarry Smith    Input Arguments:
99395fce210SBarry Smith +  sf - star forest
99495fce210SBarry Smith .  unit - data type
99595fce210SBarry Smith -  rootdata - buffer to broadcast
99695fce210SBarry Smith 
99795fce210SBarry Smith    Output Arguments:
99895fce210SBarry Smith .  leafdata - buffer to update with values from each leaf's respective root
99995fce210SBarry Smith 
100095fce210SBarry Smith    Level: intermediate
100195fce210SBarry Smith 
100295fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFReduceEnd()
100395fce210SBarry Smith @*/
100495fce210SBarry Smith PetscErrorCode PetscSFBcastEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
100595fce210SBarry Smith {
100695fce210SBarry Smith   PetscErrorCode ierr;
100795fce210SBarry Smith 
100895fce210SBarry Smith   PetscFunctionBegin;
100995fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
101095fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
1011563ffbabSMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCSF_BcastEnd,sf,0,0,0);CHKERRQ(ierr);
101295fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
101395fce210SBarry Smith   ierr = (*sf->ops->BcastEnd)(sf,unit,rootdata,leafdata);CHKERRQ(ierr);
1014563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_BcastEnd,sf,0,0,0);CHKERRQ(ierr);
101595fce210SBarry Smith   PetscFunctionReturn(0);
101695fce210SBarry Smith }
101795fce210SBarry Smith 
101895fce210SBarry Smith /*@C
101995fce210SBarry Smith    PetscSFReduceBegin - begin reduction of leafdata into rootdata, to be completed with call to PetscSFReduceEnd()
102095fce210SBarry Smith 
102195fce210SBarry Smith    Collective
102295fce210SBarry Smith 
102395fce210SBarry Smith    Input Arguments:
102495fce210SBarry Smith +  sf - star forest
102595fce210SBarry Smith .  unit - data type
102695fce210SBarry Smith .  leafdata - values to reduce
102795fce210SBarry Smith -  op - reduction operation
102895fce210SBarry Smith 
102995fce210SBarry Smith    Output Arguments:
103095fce210SBarry Smith .  rootdata - result of reduction of values from all leaves of each root
103195fce210SBarry Smith 
103295fce210SBarry Smith    Level: intermediate
103395fce210SBarry Smith 
103495fce210SBarry Smith .seealso: PetscSFBcastBegin()
103595fce210SBarry Smith @*/
103695fce210SBarry Smith PetscErrorCode PetscSFReduceBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
103795fce210SBarry Smith {
103895fce210SBarry Smith   PetscErrorCode ierr;
103995fce210SBarry Smith 
104095fce210SBarry Smith   PetscFunctionBegin;
104195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
104295fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
1043563ffbabSMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);
104495fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
104595fce210SBarry Smith   ierr = (sf->ops->ReduceBegin)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
1046563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);
104795fce210SBarry Smith   PetscFunctionReturn(0);
104895fce210SBarry Smith }
104995fce210SBarry Smith 
105095fce210SBarry Smith /*@C
105195fce210SBarry Smith    PetscSFReduceEnd - end a reduction operation started with PetscSFReduceBegin()
105295fce210SBarry Smith 
105395fce210SBarry Smith    Collective
105495fce210SBarry Smith 
105595fce210SBarry Smith    Input Arguments:
105695fce210SBarry Smith +  sf - star forest
105795fce210SBarry Smith .  unit - data type
105895fce210SBarry Smith .  leafdata - values to reduce
105995fce210SBarry Smith -  op - reduction operation
106095fce210SBarry Smith 
106195fce210SBarry Smith    Output Arguments:
106295fce210SBarry Smith .  rootdata - result of reduction of values from all leaves of each root
106395fce210SBarry Smith 
106495fce210SBarry Smith    Level: intermediate
106595fce210SBarry Smith 
106695fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFBcastEnd()
106795fce210SBarry Smith @*/
106895fce210SBarry Smith PetscErrorCode PetscSFReduceEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
106995fce210SBarry Smith {
107095fce210SBarry Smith   PetscErrorCode ierr;
107195fce210SBarry Smith 
107295fce210SBarry Smith   PetscFunctionBegin;
107395fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
107495fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
1075563ffbabSMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);
107695fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
107795fce210SBarry Smith   ierr = (*sf->ops->ReduceEnd)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
1078563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);
107995fce210SBarry Smith   PetscFunctionReturn(0);
108095fce210SBarry Smith }
108195fce210SBarry Smith 
108295fce210SBarry Smith /*@C
108395fce210SBarry Smith    PetscSFComputeDegreeBegin - begin computation of degree for each root vertex, to be completed with PetscSFComputeDegreeEnd()
108495fce210SBarry Smith 
108595fce210SBarry Smith    Collective
108695fce210SBarry Smith 
108795fce210SBarry Smith    Input Arguments:
108895fce210SBarry Smith .  sf - star forest
108995fce210SBarry Smith 
109095fce210SBarry Smith    Output Arguments:
109195fce210SBarry Smith .  degree - degree of each root vertex
109295fce210SBarry Smith 
109395fce210SBarry Smith    Level: advanced
109495fce210SBarry Smith 
109595fce210SBarry Smith .seealso: PetscSFGatherBegin()
109695fce210SBarry Smith @*/
109795fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf,const PetscInt **degree)
109895fce210SBarry Smith {
109995fce210SBarry Smith   PetscErrorCode ierr;
110095fce210SBarry Smith 
110195fce210SBarry Smith   PetscFunctionBegin;
110295fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
110395fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
110495fce210SBarry Smith   PetscValidPointer(degree,2);
1105803bd9e8SMatthew G. Knepley   if (!sf->degreeknown) {
11069837ea96SMatthew G. Knepley     PetscInt i,maxlocal;
1107803bd9e8SMatthew G. Knepley     if (sf->degree) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Calls to PetscSFComputeDegreeBegin() cannot be nested.");
11089837ea96SMatthew G. Knepley     for (i=0,maxlocal=0; i<sf->nleaves; i++) maxlocal = PetscMax(maxlocal,(sf->mine ? sf->mine[i] : i)+1);
1109785e854fSJed Brown     ierr = PetscMalloc1(sf->nroots,&sf->degree);CHKERRQ(ierr);
11109837ea96SMatthew G. Knepley     ierr = PetscMalloc1(maxlocal,&sf->degreetmp);CHKERRQ(ierr);
111195fce210SBarry Smith     for (i=0; i<sf->nroots; i++) sf->degree[i] = 0;
11129837ea96SMatthew G. Knepley     for (i=0; i<maxlocal; i++) sf->degreetmp[i] = 1;
1113dbd2ff41SBarry Smith     ierr = PetscSFReduceBegin(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);CHKERRQ(ierr);
111495fce210SBarry Smith   }
111595fce210SBarry Smith   *degree = NULL;
111695fce210SBarry Smith   PetscFunctionReturn(0);
111795fce210SBarry Smith }
111895fce210SBarry Smith 
111995fce210SBarry Smith /*@C
112095fce210SBarry Smith    PetscSFComputeDegreeEnd - complete computation of degree for each root vertex, started with PetscSFComputeDegreeBegin()
112195fce210SBarry Smith 
112295fce210SBarry Smith    Collective
112395fce210SBarry Smith 
112495fce210SBarry Smith    Input Arguments:
112595fce210SBarry Smith .  sf - star forest
112695fce210SBarry Smith 
112795fce210SBarry Smith    Output Arguments:
112895fce210SBarry Smith .  degree - degree of each root vertex
112995fce210SBarry Smith 
113095fce210SBarry Smith    Level: developer
113195fce210SBarry Smith 
113295fce210SBarry Smith .seealso:
113395fce210SBarry Smith @*/
113495fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf,const PetscInt **degree)
113595fce210SBarry Smith {
113695fce210SBarry Smith   PetscErrorCode ierr;
113795fce210SBarry Smith 
113895fce210SBarry Smith   PetscFunctionBegin;
113995fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
114095fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
114195fce210SBarry Smith   if (!sf->degreeknown) {
1142dbd2ff41SBarry Smith     ierr = PetscSFReduceEnd(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);CHKERRQ(ierr);
114395fce210SBarry Smith     ierr = PetscFree(sf->degreetmp);CHKERRQ(ierr);
114495fce210SBarry Smith 
114595fce210SBarry Smith     sf->degreeknown = PETSC_TRUE;
114695fce210SBarry Smith   }
114795fce210SBarry Smith   *degree = sf->degree;
114895fce210SBarry Smith   PetscFunctionReturn(0);
114995fce210SBarry Smith }
115095fce210SBarry Smith 
115195fce210SBarry Smith /*@C
115295fce210SBarry Smith    PetscSFFetchAndOpBegin - begin operation that fetches values from root and updates atomically by applying operation using my leaf value, to be completed with PetscSFFetchAndOpEnd()
115395fce210SBarry Smith 
115495fce210SBarry Smith    Collective
115595fce210SBarry Smith 
115695fce210SBarry Smith    Input Arguments:
115795fce210SBarry Smith +  sf - star forest
115895fce210SBarry Smith .  unit - data type
115995fce210SBarry Smith .  leafdata - leaf values to use in reduction
116095fce210SBarry Smith -  op - operation to use for reduction
116195fce210SBarry Smith 
116295fce210SBarry Smith    Output Arguments:
116395fce210SBarry Smith +  rootdata - root values to be updated, input state is seen by first process to perform an update
116495fce210SBarry Smith -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
116595fce210SBarry Smith 
116695fce210SBarry Smith    Level: advanced
116795fce210SBarry Smith 
116895fce210SBarry Smith    Note:
116995fce210SBarry Smith    The update is only atomic at the granularity provided by the hardware. Different roots referenced by the same process
117095fce210SBarry Smith    might be updated in a different order. Furthermore, if a composite type is used for the unit datatype, atomicity is
117195fce210SBarry Smith    not guaranteed across the whole vertex. Therefore, this function is mostly only used with primitive types such as
117295fce210SBarry Smith    integers.
117395fce210SBarry Smith 
117495fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFReduceBegin(), PetscSFSetGraph()
117595fce210SBarry Smith @*/
117695fce210SBarry Smith PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
117795fce210SBarry Smith {
117895fce210SBarry Smith   PetscErrorCode ierr;
117995fce210SBarry Smith 
118095fce210SBarry Smith   PetscFunctionBegin;
118195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
118295fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
1183563ffbabSMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
118495fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
118595fce210SBarry Smith   ierr = (*sf->ops->FetchAndOpBegin)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr);
1186563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
118795fce210SBarry Smith   PetscFunctionReturn(0);
118895fce210SBarry Smith }
118995fce210SBarry Smith 
119095fce210SBarry Smith /*@C
119195fce210SBarry 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
119295fce210SBarry Smith 
119395fce210SBarry Smith    Collective
119495fce210SBarry Smith 
119595fce210SBarry Smith    Input Arguments:
119695fce210SBarry Smith +  sf - star forest
119795fce210SBarry Smith .  unit - data type
119895fce210SBarry Smith .  leafdata - leaf values to use in reduction
119995fce210SBarry Smith -  op - operation to use for reduction
120095fce210SBarry Smith 
120195fce210SBarry Smith    Output Arguments:
120295fce210SBarry Smith +  rootdata - root values to be updated, input state is seen by first process to perform an update
120395fce210SBarry Smith -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
120495fce210SBarry Smith 
120595fce210SBarry Smith    Level: advanced
120695fce210SBarry Smith 
120795fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFReduceEnd(), PetscSFSetGraph()
120895fce210SBarry Smith @*/
120995fce210SBarry Smith PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
121095fce210SBarry Smith {
121195fce210SBarry Smith   PetscErrorCode ierr;
121295fce210SBarry Smith 
121395fce210SBarry Smith   PetscFunctionBegin;
121495fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
121595fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
1216563ffbabSMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
121795fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
121895fce210SBarry Smith   ierr = (*sf->ops->FetchAndOpEnd)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr);
1219563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
122095fce210SBarry Smith   PetscFunctionReturn(0);
122195fce210SBarry Smith }
122295fce210SBarry Smith 
122395fce210SBarry Smith /*@C
122495fce210SBarry Smith    PetscSFGatherBegin - begin pointwise gather of all leaves into multi-roots, to be completed with PetscSFGatherEnd()
122595fce210SBarry Smith 
122695fce210SBarry Smith    Collective
122795fce210SBarry Smith 
122895fce210SBarry Smith    Input Arguments:
122995fce210SBarry Smith +  sf - star forest
123095fce210SBarry Smith .  unit - data type
123195fce210SBarry Smith -  leafdata - leaf data to gather to roots
123295fce210SBarry Smith 
123395fce210SBarry Smith    Output Argument:
123495fce210SBarry Smith .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
123595fce210SBarry Smith 
123695fce210SBarry Smith    Level: intermediate
123795fce210SBarry Smith 
123895fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
123995fce210SBarry Smith @*/
124095fce210SBarry Smith PetscErrorCode PetscSFGatherBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
124195fce210SBarry Smith {
124295fce210SBarry Smith   PetscErrorCode ierr;
124395fce210SBarry Smith   PetscSF        multi;
124495fce210SBarry Smith 
124595fce210SBarry Smith   PetscFunctionBegin;
124695fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
124795fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
12488bfbc91cSJed Brown   ierr = PetscSFReduceBegin(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr);
124995fce210SBarry Smith   PetscFunctionReturn(0);
125095fce210SBarry Smith }
125195fce210SBarry Smith 
125295fce210SBarry Smith /*@C
125395fce210SBarry Smith    PetscSFGatherEnd - ends pointwise gather operation that was started with PetscSFGatherBegin()
125495fce210SBarry Smith 
125595fce210SBarry Smith    Collective
125695fce210SBarry Smith 
125795fce210SBarry Smith    Input Arguments:
125895fce210SBarry Smith +  sf - star forest
125995fce210SBarry Smith .  unit - data type
126095fce210SBarry Smith -  leafdata - leaf data to gather to roots
126195fce210SBarry Smith 
126295fce210SBarry Smith    Output Argument:
126395fce210SBarry Smith .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
126495fce210SBarry Smith 
126595fce210SBarry Smith    Level: intermediate
126695fce210SBarry Smith 
126795fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
126895fce210SBarry Smith @*/
126995fce210SBarry Smith PetscErrorCode PetscSFGatherEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
127095fce210SBarry Smith {
127195fce210SBarry Smith   PetscErrorCode ierr;
127295fce210SBarry Smith   PetscSF        multi;
127395fce210SBarry Smith 
127495fce210SBarry Smith   PetscFunctionBegin;
127595fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
127695fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
127795fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
127895fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
12798bfbc91cSJed Brown   ierr = PetscSFReduceEnd(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr);
128095fce210SBarry Smith   PetscFunctionReturn(0);
128195fce210SBarry Smith }
128295fce210SBarry Smith 
128395fce210SBarry Smith /*@C
128495fce210SBarry Smith    PetscSFScatterBegin - begin pointwise scatter operation from multi-roots to leaves, to be completed with PetscSFScatterEnd()
128595fce210SBarry Smith 
128695fce210SBarry Smith    Collective
128795fce210SBarry Smith 
128895fce210SBarry Smith    Input Arguments:
128995fce210SBarry Smith +  sf - star forest
129095fce210SBarry Smith .  unit - data type
129195fce210SBarry Smith -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
129295fce210SBarry Smith 
129395fce210SBarry Smith    Output Argument:
129495fce210SBarry Smith .  leafdata - leaf data to be update with personal data from each respective root
129595fce210SBarry Smith 
129695fce210SBarry Smith    Level: intermediate
129795fce210SBarry Smith 
129895fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
129995fce210SBarry Smith @*/
130095fce210SBarry Smith PetscErrorCode PetscSFScatterBegin(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
130195fce210SBarry Smith {
130295fce210SBarry Smith   PetscErrorCode ierr;
130395fce210SBarry Smith   PetscSF        multi;
130495fce210SBarry Smith 
130595fce210SBarry Smith   PetscFunctionBegin;
130695fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
130795fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
130895fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
130995fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
131095fce210SBarry Smith   ierr = PetscSFBcastBegin(multi,unit,multirootdata,leafdata);CHKERRQ(ierr);
131195fce210SBarry Smith   PetscFunctionReturn(0);
131295fce210SBarry Smith }
131395fce210SBarry Smith 
131495fce210SBarry Smith /*@C
131595fce210SBarry Smith    PetscSFScatterEnd - ends pointwise scatter operation that was started with PetscSFScatterBegin()
131695fce210SBarry Smith 
131795fce210SBarry Smith    Collective
131895fce210SBarry Smith 
131995fce210SBarry Smith    Input Arguments:
132095fce210SBarry Smith +  sf - star forest
132195fce210SBarry Smith .  unit - data type
132295fce210SBarry Smith -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
132395fce210SBarry Smith 
132495fce210SBarry Smith    Output Argument:
132595fce210SBarry Smith .  leafdata - leaf data to be update with personal data from each respective root
132695fce210SBarry Smith 
132795fce210SBarry Smith    Level: intermediate
132895fce210SBarry Smith 
132995fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
133095fce210SBarry Smith @*/
133195fce210SBarry Smith PetscErrorCode PetscSFScatterEnd(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
133295fce210SBarry Smith {
133395fce210SBarry Smith   PetscErrorCode ierr;
133495fce210SBarry Smith   PetscSF        multi;
133595fce210SBarry Smith 
133695fce210SBarry Smith   PetscFunctionBegin;
133795fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
133895fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
133995fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
134095fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
134195fce210SBarry Smith   ierr = PetscSFBcastEnd(multi,unit,multirootdata,leafdata);CHKERRQ(ierr);
134295fce210SBarry Smith   PetscFunctionReturn(0);
134395fce210SBarry Smith }
1344a7b3aa13SAta Mesgarnejad 
1345a7b3aa13SAta Mesgarnejad /*@
1346a7b3aa13SAta Mesgarnejad   PetscSFCompose - Compose a new PetscSF equivalent to action to PetscSFs
1347a7b3aa13SAta Mesgarnejad 
1348a7b3aa13SAta Mesgarnejad   Input Parameters:
1349a7b3aa13SAta Mesgarnejad + sfA - The first PetscSF
1350a7b3aa13SAta Mesgarnejad - sfB - The second PetscSF
1351a7b3aa13SAta Mesgarnejad 
1352a7b3aa13SAta Mesgarnejad   Output Parameters:
1353a7b3aa13SAta Mesgarnejad . sfBA - equvalent PetscSF for applying A then B
1354a7b3aa13SAta Mesgarnejad 
1355a7b3aa13SAta Mesgarnejad   Level: developer
1356a7b3aa13SAta Mesgarnejad 
1357a7b3aa13SAta Mesgarnejad .seealso: PetscSF, PetscSFGetGraph(), PetscSFSetGraph()
1358a7b3aa13SAta Mesgarnejad @*/
1359a7b3aa13SAta Mesgarnejad PetscErrorCode PetscSFCompose(PetscSF sfA, PetscSF sfB, PetscSF *sfBA)
1360a7b3aa13SAta Mesgarnejad {
1361a7b3aa13SAta Mesgarnejad   MPI_Comm           comm;
1362a7b3aa13SAta Mesgarnejad   const PetscSFNode *remotePointsA, *remotePointsB;
1363a7b3aa13SAta Mesgarnejad   PetscSFNode       *remotePointsBA;
1364a7b3aa13SAta Mesgarnejad   const PetscInt    *localPointsA, *localPointsB;
1365a7b3aa13SAta Mesgarnejad   PetscInt           numRootsA, numLeavesA, numRootsB, numLeavesB;
1366a7b3aa13SAta Mesgarnejad   PetscErrorCode     ierr;
1367a7b3aa13SAta Mesgarnejad 
1368a7b3aa13SAta Mesgarnejad   PetscFunctionBegin;
1369a7b3aa13SAta Mesgarnejad   PetscValidHeaderSpecific(sfA, PETSCSF_CLASSID, 1);
1370a7b3aa13SAta Mesgarnejad   PetscValidHeaderSpecific(sfB, PETSCSF_CLASSID, 1);
1371a7b3aa13SAta Mesgarnejad   ierr = PetscObjectGetComm((PetscObject) sfA, &comm);CHKERRQ(ierr);
1372a7b3aa13SAta Mesgarnejad   ierr = PetscSFGetGraph(sfA, &numRootsA, &numLeavesA, &localPointsA, &remotePointsA);CHKERRQ(ierr);
1373a7b3aa13SAta Mesgarnejad   ierr = PetscSFGetGraph(sfB, &numRootsB, &numLeavesB, &localPointsB, &remotePointsB);CHKERRQ(ierr);
1374a7b3aa13SAta Mesgarnejad   ierr = PetscMalloc1(numLeavesB, &remotePointsBA);CHKERRQ(ierr);
1375a7b3aa13SAta Mesgarnejad   ierr = PetscSFBcastBegin(sfB, MPIU_2INT, remotePointsA, remotePointsBA);CHKERRQ(ierr);
1376a7b3aa13SAta Mesgarnejad   ierr = PetscSFBcastEnd(sfB, MPIU_2INT, remotePointsA, remotePointsBA);CHKERRQ(ierr);
1377a7b3aa13SAta Mesgarnejad   ierr = PetscSFCreate(comm, sfBA);CHKERRQ(ierr);
1378a7b3aa13SAta Mesgarnejad   ierr = PetscSFSetGraph(*sfBA, numRootsA, numLeavesB, localPointsB, PETSC_COPY_VALUES, remotePointsBA, PETSC_OWN_POINTER);CHKERRQ(ierr);
1379a7b3aa13SAta Mesgarnejad   PetscFunctionReturn(0);
1380a7b3aa13SAta Mesgarnejad }
1381