xref: /petsc/src/vec/is/sf/interface/sf.c (revision 79715d568eda8fb3c2b41070aa5caa67af195ae0)
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 
1895fce210SBarry Smith /*@C
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);
74*79715d56SJed 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);
8095fce210SBarry Smith   ierr       = PetscFree(sf->degree);CHKERRQ(ierr);
8195fce210SBarry Smith   if (sf->ingroup  != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->ingroup);CHKERRQ(ierr);}
8295fce210SBarry Smith   if (sf->outgroup != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->outgroup);CHKERRQ(ierr);}
8395fce210SBarry Smith   ierr         = PetscSFDestroy(&sf->multi);CHKERRQ(ierr);
8495fce210SBarry Smith   sf->graphset = PETSC_FALSE;
8595fce210SBarry Smith   sf->setupcalled = PETSC_FALSE;
8695fce210SBarry Smith   PetscFunctionReturn(0);
8795fce210SBarry Smith }
8895fce210SBarry Smith 
8995fce210SBarry Smith /*@C
9095fce210SBarry Smith    PetscSFSetType - set the PetscSF communication implementation
9195fce210SBarry Smith 
9295fce210SBarry Smith    Collective on PetscSF
9395fce210SBarry Smith 
9495fce210SBarry Smith    Input Parameters:
9595fce210SBarry Smith +  sf - the PetscSF context
9695fce210SBarry Smith -  type - a known method
9795fce210SBarry Smith 
9895fce210SBarry Smith    Options Database Key:
9995fce210SBarry Smith .  -sf_type <type> - Sets the method; use -help for a list
10095fce210SBarry Smith    of available methods (for instance, window, pt2pt, neighbor)
10195fce210SBarry Smith 
10295fce210SBarry Smith    Notes:
10395fce210SBarry Smith    See "include/petscsf.h" for available methods (for instance)
10495fce210SBarry Smith +    PETSCSFWINDOW - MPI-2/3 one-sided
10595fce210SBarry Smith -    PETSCSFBASIC - basic implementation using MPI-1 two-sided
10695fce210SBarry Smith 
10795fce210SBarry Smith   Level: intermediate
10895fce210SBarry Smith 
10995fce210SBarry Smith .keywords: PetscSF, set, type
11095fce210SBarry Smith 
11195fce210SBarry Smith .seealso: PetscSFType, PetscSFCreate()
11295fce210SBarry Smith @*/
11395fce210SBarry Smith PetscErrorCode PetscSFSetType(PetscSF sf,PetscSFType type)
11495fce210SBarry Smith {
11595fce210SBarry Smith   PetscErrorCode ierr,(*r)(PetscSF);
11695fce210SBarry Smith   PetscBool      match;
11795fce210SBarry Smith 
11895fce210SBarry Smith   PetscFunctionBegin;
11995fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
12095fce210SBarry Smith   PetscValidCharPointer(type,2);
12195fce210SBarry Smith 
12295fce210SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)sf,type,&match);CHKERRQ(ierr);
12395fce210SBarry Smith   if (match) PetscFunctionReturn(0);
12495fce210SBarry Smith 
125adc40e5bSBarry Smith   ierr = PetscFunctionListFind(PetscSFList,type,&r);CHKERRQ(ierr);
12695fce210SBarry Smith   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested PetscSF type %s",type);
12795fce210SBarry Smith   /* Destroy the previous private PetscSF context */
12895fce210SBarry Smith   if (sf->ops->Destroy) {
12995fce210SBarry Smith     ierr = (*(sf)->ops->Destroy)(sf);CHKERRQ(ierr);
13095fce210SBarry Smith   }
13195fce210SBarry Smith   ierr = PetscMemzero(sf->ops,sizeof(*sf->ops));CHKERRQ(ierr);
13295fce210SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)sf,type);CHKERRQ(ierr);
13395fce210SBarry Smith   ierr = (*r)(sf);CHKERRQ(ierr);
13495fce210SBarry Smith   PetscFunctionReturn(0);
13595fce210SBarry Smith }
13695fce210SBarry Smith 
137d36d33e4SMatthew G. Knepley /*@
13895fce210SBarry Smith    PetscSFDestroy - destroy star forest
13995fce210SBarry Smith 
14095fce210SBarry Smith    Collective
14195fce210SBarry Smith 
14295fce210SBarry Smith    Input Arguments:
14395fce210SBarry Smith .  sf - address of star forest
14495fce210SBarry Smith 
14595fce210SBarry Smith    Level: intermediate
14695fce210SBarry Smith 
14795fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFReset()
14895fce210SBarry Smith @*/
14995fce210SBarry Smith PetscErrorCode PetscSFDestroy(PetscSF *sf)
15095fce210SBarry Smith {
15195fce210SBarry Smith   PetscErrorCode ierr;
15295fce210SBarry Smith 
15395fce210SBarry Smith   PetscFunctionBegin;
15495fce210SBarry Smith   if (!*sf) PetscFunctionReturn(0);
15595fce210SBarry Smith   PetscValidHeaderSpecific((*sf),PETSCSF_CLASSID,1);
15695fce210SBarry Smith   if (--((PetscObject)(*sf))->refct > 0) {*sf = 0; PetscFunctionReturn(0);}
15795fce210SBarry Smith   ierr = PetscSFReset(*sf);CHKERRQ(ierr);
15895fce210SBarry Smith   if ((*sf)->ops->Destroy) {ierr = (*(*sf)->ops->Destroy)(*sf);CHKERRQ(ierr);}
15995fce210SBarry Smith   ierr = PetscHeaderDestroy(sf);CHKERRQ(ierr);
16095fce210SBarry Smith   PetscFunctionReturn(0);
16195fce210SBarry Smith }
16295fce210SBarry Smith 
16395fce210SBarry Smith /*@
16495fce210SBarry Smith    PetscSFSetUp - set up communication structures
16595fce210SBarry Smith 
16695fce210SBarry Smith    Collective
16795fce210SBarry Smith 
16895fce210SBarry Smith    Input Arguments:
16995fce210SBarry Smith .  sf - star forest communication object
17095fce210SBarry Smith 
17195fce210SBarry Smith    Level: beginner
17295fce210SBarry Smith 
17395fce210SBarry Smith .seealso: PetscSFSetFromOptions(), PetscSFSetType()
17495fce210SBarry Smith @*/
17595fce210SBarry Smith PetscErrorCode PetscSFSetUp(PetscSF sf)
17695fce210SBarry Smith {
17795fce210SBarry Smith   PetscErrorCode ierr;
17895fce210SBarry Smith 
17995fce210SBarry Smith   PetscFunctionBegin;
18095fce210SBarry Smith   if (sf->setupcalled) PetscFunctionReturn(0);
18195fce210SBarry Smith   if (!((PetscObject)sf)->type_name) {ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr);}
18295fce210SBarry Smith   if (sf->ops->SetUp) {ierr = (*sf->ops->SetUp)(sf);CHKERRQ(ierr);}
18395fce210SBarry Smith   sf->setupcalled = PETSC_TRUE;
18495fce210SBarry Smith   PetscFunctionReturn(0);
18595fce210SBarry Smith }
18695fce210SBarry Smith 
18795fce210SBarry Smith /*@C
18895fce210SBarry Smith    PetscSFSetFromOptions - set PetscSF options using the options database
18995fce210SBarry Smith 
19095fce210SBarry Smith    Logically Collective
19195fce210SBarry Smith 
19295fce210SBarry Smith    Input Arguments:
19395fce210SBarry Smith .  sf - star forest
19495fce210SBarry Smith 
19595fce210SBarry Smith    Options Database Keys:
19660263706SJed Brown +  -sf_type - implementation type, see PetscSFSetType()
19760263706SJed Brown -  -sf_rank_order - sort composite points for gathers and scatters in rank order, gathers are non-deterministic otherwise
19895fce210SBarry Smith 
19995fce210SBarry Smith    Level: intermediate
20095fce210SBarry Smith 
20195fce210SBarry Smith .keywords: KSP, set, from, options, database
20295fce210SBarry Smith 
20395fce210SBarry Smith .seealso: PetscSFWindowSetSyncType()
20495fce210SBarry Smith @*/
20595fce210SBarry Smith PetscErrorCode PetscSFSetFromOptions(PetscSF sf)
20695fce210SBarry Smith {
20795fce210SBarry Smith   PetscSFType    deft;
20895fce210SBarry Smith   char           type[256];
20995fce210SBarry Smith   PetscErrorCode ierr;
21095fce210SBarry Smith   PetscBool      flg;
21195fce210SBarry Smith 
21295fce210SBarry Smith   PetscFunctionBegin;
21395fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
21495fce210SBarry Smith   ierr = PetscObjectOptionsBegin((PetscObject)sf);CHKERRQ(ierr);
21595fce210SBarry Smith   deft = ((PetscObject)sf)->type_name ? ((PetscObject)sf)->type_name : PETSCSFBASIC;
216a264d7a6SBarry Smith   ierr = PetscOptionsFList("-sf_type","PetscSF implementation type","PetscSFSetType",PetscSFList,deft,type,256,&flg);CHKERRQ(ierr);
21795fce210SBarry Smith   ierr = PetscSFSetType(sf,flg ? type : deft);CHKERRQ(ierr);
21895fce210SBarry 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);
219e55864a3SBarry Smith   if (sf->ops->SetFromOptions) {ierr = (*sf->ops->SetFromOptions)(PetscOptionsObject,sf);CHKERRQ(ierr);}
22095fce210SBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
22195fce210SBarry Smith   PetscFunctionReturn(0);
22295fce210SBarry Smith }
22395fce210SBarry Smith 
22495fce210SBarry Smith /*@C
22595fce210SBarry Smith    PetscSFSetRankOrder - sort multi-points for gathers and scatters by rank order
22695fce210SBarry Smith 
22795fce210SBarry Smith    Logically Collective
22895fce210SBarry Smith 
22995fce210SBarry Smith    Input Arguments:
23095fce210SBarry Smith +  sf - star forest
23195fce210SBarry Smith -  flg - PETSC_TRUE to sort, PETSC_FALSE to skip sorting (lower setup cost, but non-deterministic)
23295fce210SBarry Smith 
23395fce210SBarry Smith    Level: advanced
23495fce210SBarry Smith 
23595fce210SBarry Smith .seealso: PetscSFGatherBegin(), PetscSFScatterBegin()
23695fce210SBarry Smith @*/
23795fce210SBarry Smith PetscErrorCode PetscSFSetRankOrder(PetscSF sf,PetscBool flg)
23895fce210SBarry Smith {
23995fce210SBarry Smith 
24095fce210SBarry Smith   PetscFunctionBegin;
24195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
24295fce210SBarry Smith   PetscValidLogicalCollectiveBool(sf,flg,2);
24395fce210SBarry Smith   if (sf->multi) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Rank ordering must be set before first call to PetscSFGatherBegin() or PetscSFScatterBegin()");
24495fce210SBarry Smith   sf->rankorder = flg;
24595fce210SBarry Smith   PetscFunctionReturn(0);
24695fce210SBarry Smith }
24795fce210SBarry Smith 
24895fce210SBarry Smith /*@C
24995fce210SBarry Smith    PetscSFSetGraph - Set a parallel star forest
25095fce210SBarry Smith 
25195fce210SBarry Smith    Collective
25295fce210SBarry Smith 
25395fce210SBarry Smith    Input Arguments:
25495fce210SBarry Smith +  sf - star forest
25595fce210SBarry Smith .  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
25695fce210SBarry Smith .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
25795fce210SBarry Smith .  ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage
25895fce210SBarry Smith .  localmode - copy mode for ilocal
25995fce210SBarry Smith .  iremote - remote locations of root vertices for each leaf on the current process
26095fce210SBarry Smith -  remotemode - copy mode for iremote
26195fce210SBarry Smith 
26295fce210SBarry Smith    Level: intermediate
26395fce210SBarry Smith 
26495fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFGetGraph()
26595fce210SBarry Smith @*/
26695fce210SBarry Smith PetscErrorCode PetscSFSetGraph(PetscSF sf,PetscInt nroots,PetscInt nleaves,const PetscInt *ilocal,PetscCopyMode localmode,const PetscSFNode *iremote,PetscCopyMode remotemode)
26795fce210SBarry Smith {
26895fce210SBarry Smith   PetscErrorCode     ierr;
26995fce210SBarry Smith   PetscTable         table;
27095fce210SBarry Smith   PetscTablePosition pos;
27195fce210SBarry Smith   PetscMPIInt        size;
27295fce210SBarry Smith   PetscInt           i,*rcount,*ranks;
27395fce210SBarry Smith 
27495fce210SBarry Smith   PetscFunctionBegin;
27595fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
276563ffbabSMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr);
27795fce210SBarry Smith   if (nleaves && ilocal) PetscValidIntPointer(ilocal,4);
27895fce210SBarry Smith   if (nleaves) PetscValidPointer(iremote,6);
27995fce210SBarry Smith   if (nroots < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"roots %D, cannot be negative",nroots);
28095fce210SBarry Smith   if (nleaves < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nleaves %D, cannot be negative",nleaves);
28195fce210SBarry Smith   ierr        = PetscSFReset(sf);CHKERRQ(ierr);
28295fce210SBarry Smith   sf->nroots  = nroots;
28395fce210SBarry Smith   sf->nleaves = nleaves;
28495fce210SBarry Smith   if (ilocal) {
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   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr);
32895fce210SBarry Smith   ierr = PetscTableCreate(10,size,&table);CHKERRQ(ierr);
32995fce210SBarry Smith   for (i=0; i<nleaves; i++) {
33095fce210SBarry Smith     /* Log 1-based rank */
33195fce210SBarry Smith     ierr = PetscTableAdd(table,iremote[i].rank+1,1,ADD_VALUES);CHKERRQ(ierr);
33295fce210SBarry Smith   }
33395fce210SBarry Smith   ierr = PetscTableGetCount(table,&sf->nranks);CHKERRQ(ierr);
334dcca6d9dSJed Brown   ierr = PetscMalloc4(sf->nranks,&sf->ranks,sf->nranks+1,&sf->roffset,nleaves,&sf->rmine,nleaves,&sf->rremote);CHKERRQ(ierr);
335dcca6d9dSJed Brown   ierr = PetscMalloc2(sf->nranks,&rcount,sf->nranks,&ranks);CHKERRQ(ierr);
33695fce210SBarry Smith   ierr = PetscTableGetHeadPosition(table,&pos);CHKERRQ(ierr);
33795fce210SBarry Smith   for (i=0; i<sf->nranks; i++) {
33895fce210SBarry Smith     ierr = PetscTableGetNext(table,&pos,&ranks[i],&rcount[i]);CHKERRQ(ierr);
33995fce210SBarry Smith     ranks[i]--;             /* Convert back to 0-based */
34095fce210SBarry Smith   }
34195fce210SBarry Smith   ierr = PetscTableDestroy(&table);CHKERRQ(ierr);
34295fce210SBarry Smith   ierr = PetscSortIntWithArray(sf->nranks,ranks,rcount);CHKERRQ(ierr);
34395fce210SBarry Smith   sf->roffset[0] = 0;
34495fce210SBarry Smith   for (i=0; i<sf->nranks; i++) {
34595fce210SBarry Smith     ierr = PetscMPIIntCast(ranks[i],sf->ranks+i);CHKERRQ(ierr);
34695fce210SBarry Smith     sf->roffset[i+1] = sf->roffset[i] + rcount[i];
34795fce210SBarry Smith     rcount[i]        = 0;
34895fce210SBarry Smith   }
34995fce210SBarry Smith   for (i=0; i<nleaves; i++) {
35095fce210SBarry Smith     PetscInt lo,hi,irank;
35195fce210SBarry Smith     /* Search for index of iremote[i].rank in sf->ranks */
35295fce210SBarry Smith     lo = 0; hi = sf->nranks;
35395fce210SBarry Smith     while (hi - lo > 1) {
35495fce210SBarry Smith       PetscInt mid = lo + (hi - lo)/2;
35595fce210SBarry Smith       if (iremote[i].rank < sf->ranks[mid]) hi = mid;
35695fce210SBarry Smith       else                                  lo = mid;
35795fce210SBarry Smith     }
35895fce210SBarry Smith     if (hi - lo == 1 && iremote[i].rank == sf->ranks[lo]) irank = lo;
35995fce210SBarry Smith     else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Could not find rank %D in array",iremote[i].rank);
36095fce210SBarry Smith     sf->rmine[sf->roffset[irank] + rcount[irank]]   = ilocal ? ilocal[i] : i;
36195fce210SBarry Smith     sf->rremote[sf->roffset[irank] + rcount[irank]] = iremote[i].index;
36295fce210SBarry Smith     rcount[irank]++;
36395fce210SBarry Smith   }
36495fce210SBarry Smith   ierr = PetscFree2(rcount,ranks);CHKERRQ(ierr);
36595fce210SBarry Smith 
36695fce210SBarry Smith   sf->graphset = PETSC_TRUE;
367563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr);
36895fce210SBarry Smith   PetscFunctionReturn(0);
36995fce210SBarry Smith }
37095fce210SBarry Smith 
37195fce210SBarry Smith /*@C
37295fce210SBarry Smith    PetscSFCreateInverseSF - given a PetscSF in which all vertices have degree 1, creates the inverse map
37395fce210SBarry Smith 
37495fce210SBarry Smith    Collective
37595fce210SBarry Smith 
37695fce210SBarry Smith    Input Arguments:
37795fce210SBarry Smith .  sf - star forest to invert
37895fce210SBarry Smith 
37995fce210SBarry Smith    Output Arguments:
38095fce210SBarry Smith .  isf - inverse of sf
38195fce210SBarry Smith 
38295fce210SBarry Smith    Level: advanced
38395fce210SBarry Smith 
38495fce210SBarry Smith    Notes:
38595fce210SBarry Smith    All roots must have degree 1.
38695fce210SBarry Smith 
38795fce210SBarry Smith    The local space may be a permutation, but cannot be sparse.
38895fce210SBarry Smith 
38995fce210SBarry Smith .seealso: PetscSFSetGraph()
39095fce210SBarry Smith @*/
39195fce210SBarry Smith PetscErrorCode PetscSFCreateInverseSF(PetscSF sf,PetscSF *isf)
39295fce210SBarry Smith {
39395fce210SBarry Smith   PetscErrorCode ierr;
39495fce210SBarry Smith   PetscMPIInt    rank;
39595fce210SBarry Smith   PetscInt       i,nroots,nleaves,maxlocal,count,*newilocal;
39695fce210SBarry Smith   const PetscInt *ilocal;
39795fce210SBarry Smith   PetscSFNode    *roots,*leaves;
39895fce210SBarry Smith 
39995fce210SBarry Smith   PetscFunctionBegin;
40095fce210SBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
40195fce210SBarry Smith   ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);CHKERRQ(ierr);
40295fce210SBarry Smith   for (i=0,maxlocal=0; i<nleaves; i++) maxlocal = PetscMax(maxlocal,(ilocal ? ilocal[i] : i)+1);
403ae9aee6dSMatthew G. Knepley   ierr = PetscMalloc2(nroots,&roots,maxlocal,&leaves);CHKERRQ(ierr);
404ae9aee6dSMatthew G. Knepley   for (i=0; i<maxlocal; i++) {
40595fce210SBarry Smith     leaves[i].rank  = rank;
40695fce210SBarry Smith     leaves[i].index = i;
40795fce210SBarry Smith   }
40895fce210SBarry Smith   for (i=0; i <nroots; i++) {
40995fce210SBarry Smith     roots[i].rank  = -1;
41095fce210SBarry Smith     roots[i].index = -1;
41195fce210SBarry Smith   }
4128bfbc91cSJed Brown   ierr = PetscSFReduceBegin(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr);
4138bfbc91cSJed Brown   ierr = PetscSFReduceEnd(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr);
41495fce210SBarry Smith 
41595fce210SBarry Smith   /* Check whether our leaves are sparse */
41695fce210SBarry Smith   for (i=0,count=0; i<nroots; i++) if (roots[i].rank >= 0) count++;
41795fce210SBarry Smith   if (count == nroots) newilocal = NULL;
41895fce210SBarry Smith   else {                        /* Index for sparse leaves and compact "roots" array (which is to become our leaves). */
419785e854fSJed Brown     ierr = PetscMalloc1(count,&newilocal);CHKERRQ(ierr);
42095fce210SBarry Smith     for (i=0,count=0; i<nroots; i++) {
42195fce210SBarry Smith       if (roots[i].rank >= 0) {
42295fce210SBarry Smith         newilocal[count]   = i;
42395fce210SBarry Smith         roots[count].rank  = roots[i].rank;
42495fce210SBarry Smith         roots[count].index = roots[i].index;
42595fce210SBarry Smith         count++;
42695fce210SBarry Smith       }
42795fce210SBarry Smith     }
42895fce210SBarry Smith   }
42995fce210SBarry Smith 
43095fce210SBarry Smith   ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,isf);CHKERRQ(ierr);
43195fce210SBarry Smith   ierr = PetscSFSetGraph(*isf,maxlocal,count,newilocal,PETSC_OWN_POINTER,roots,PETSC_COPY_VALUES);CHKERRQ(ierr);
43295fce210SBarry Smith   ierr = PetscFree2(roots,leaves);CHKERRQ(ierr);
43395fce210SBarry Smith   PetscFunctionReturn(0);
43495fce210SBarry Smith }
43595fce210SBarry Smith 
43695fce210SBarry Smith /*@
43795fce210SBarry Smith    PetscSFDuplicate - duplicate a PetscSF, optionally preserving rank connectivity and graph
43895fce210SBarry Smith 
43995fce210SBarry Smith    Collective
44095fce210SBarry Smith 
44195fce210SBarry Smith    Input Arguments:
44295fce210SBarry Smith +  sf - communication object to duplicate
44395fce210SBarry Smith -  opt - PETSCSF_DUPLICATE_CONFONLY, PETSCSF_DUPLICATE_RANKS, or PETSCSF_DUPLICATE_GRAPH (see PetscSFDuplicateOption)
44495fce210SBarry Smith 
44595fce210SBarry Smith    Output Arguments:
44695fce210SBarry Smith .  newsf - new communication object
44795fce210SBarry Smith 
44895fce210SBarry Smith    Level: beginner
44995fce210SBarry Smith 
45095fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetType(), PetscSFSetGraph()
45195fce210SBarry Smith @*/
45295fce210SBarry Smith PetscErrorCode PetscSFDuplicate(PetscSF sf,PetscSFDuplicateOption opt,PetscSF *newsf)
45395fce210SBarry Smith {
45495fce210SBarry Smith   PetscErrorCode ierr;
45595fce210SBarry Smith 
45695fce210SBarry Smith   PetscFunctionBegin;
45795fce210SBarry Smith   ierr = PetscSFCreate(PetscObjectComm((PetscObject)sf),newsf);CHKERRQ(ierr);
45895fce210SBarry Smith   ierr = PetscSFSetType(*newsf,((PetscObject)sf)->type_name);CHKERRQ(ierr);
45995fce210SBarry Smith   if (sf->ops->Duplicate) {ierr = (*sf->ops->Duplicate)(sf,opt,*newsf);CHKERRQ(ierr);}
46095fce210SBarry Smith   if (opt == PETSCSF_DUPLICATE_GRAPH) {
46195fce210SBarry Smith     PetscInt          nroots,nleaves;
46295fce210SBarry Smith     const PetscInt    *ilocal;
46395fce210SBarry Smith     const PetscSFNode *iremote;
46495fce210SBarry Smith     ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
46595fce210SBarry Smith     ierr = PetscSFSetGraph(*newsf,nroots,nleaves,ilocal,PETSC_COPY_VALUES,iremote,PETSC_COPY_VALUES);CHKERRQ(ierr);
46695fce210SBarry Smith   }
46795fce210SBarry Smith   PetscFunctionReturn(0);
46895fce210SBarry Smith }
46995fce210SBarry Smith 
47095fce210SBarry Smith /*@C
47195fce210SBarry Smith    PetscSFGetGraph - Get the graph specifying a parallel star forest
47295fce210SBarry Smith 
47395fce210SBarry Smith    Not Collective
47495fce210SBarry Smith 
47595fce210SBarry Smith    Input Arguments:
47695fce210SBarry Smith .  sf - star forest
47795fce210SBarry Smith 
47895fce210SBarry Smith    Output Arguments:
47995fce210SBarry Smith +  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
48095fce210SBarry Smith .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
48195fce210SBarry Smith .  ilocal - locations of leaves in leafdata buffers
48295fce210SBarry Smith -  iremote - remote locations of root vertices for each leaf on the current process
48395fce210SBarry Smith 
48495fce210SBarry Smith    Level: intermediate
48595fce210SBarry Smith 
48695fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph()
48795fce210SBarry Smith @*/
48895fce210SBarry Smith PetscErrorCode PetscSFGetGraph(PetscSF sf,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote)
48995fce210SBarry Smith {
49095fce210SBarry Smith 
49195fce210SBarry Smith   PetscFunctionBegin;
49295fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
49395fce210SBarry Smith   /* We are not currently requiring that the graph is set, thus returning nroots=-1 if it has not been set */
49495fce210SBarry Smith   /* if (!sf->graphset) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Graph has not been set, must call PetscSFSetGraph()"); */
49595fce210SBarry Smith   if (nroots) *nroots = sf->nroots;
49695fce210SBarry Smith   if (nleaves) *nleaves = sf->nleaves;
49795fce210SBarry Smith   if (ilocal) *ilocal = sf->mine;
49895fce210SBarry Smith   if (iremote) *iremote = sf->remote;
49995fce210SBarry Smith   PetscFunctionReturn(0);
50095fce210SBarry Smith }
50195fce210SBarry Smith 
50295fce210SBarry Smith /*@C
50395fce210SBarry Smith    PetscSFGetLeafRange - Get the active leaf ranges
50495fce210SBarry Smith 
50595fce210SBarry Smith    Not Collective
50695fce210SBarry Smith 
50795fce210SBarry Smith    Input Arguments:
50895fce210SBarry Smith .  sf - star forest
50995fce210SBarry Smith 
51095fce210SBarry Smith    Output Arguments:
51195fce210SBarry Smith +  minleaf - minimum active leaf on this process
51295fce210SBarry Smith -  maxleaf - maximum active leaf on this process
51395fce210SBarry Smith 
51495fce210SBarry Smith    Level: developer
51595fce210SBarry Smith 
51695fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph(), PetscSFGetGraph()
51795fce210SBarry Smith @*/
51895fce210SBarry Smith PetscErrorCode PetscSFGetLeafRange(PetscSF sf,PetscInt *minleaf,PetscInt *maxleaf)
51995fce210SBarry Smith {
52095fce210SBarry Smith 
52195fce210SBarry Smith   PetscFunctionBegin;
52295fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
52395fce210SBarry Smith   if (minleaf) *minleaf = sf->minleaf;
52495fce210SBarry Smith   if (maxleaf) *maxleaf = sf->maxleaf;
52595fce210SBarry Smith   PetscFunctionReturn(0);
52695fce210SBarry Smith }
52795fce210SBarry Smith 
52895fce210SBarry Smith /*@C
52995fce210SBarry Smith    PetscSFView - view a star forest
53095fce210SBarry Smith 
53195fce210SBarry Smith    Collective
53295fce210SBarry Smith 
53395fce210SBarry Smith    Input Arguments:
53495fce210SBarry Smith +  sf - star forest
53595fce210SBarry Smith -  viewer - viewer to display graph, for example PETSC_VIEWER_STDOUT_WORLD
53695fce210SBarry Smith 
53795fce210SBarry Smith    Level: beginner
53895fce210SBarry Smith 
53995fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph()
54095fce210SBarry Smith @*/
54195fce210SBarry Smith PetscErrorCode PetscSFView(PetscSF sf,PetscViewer viewer)
54295fce210SBarry Smith {
54395fce210SBarry Smith   PetscErrorCode    ierr;
54495fce210SBarry Smith   PetscBool         iascii;
54595fce210SBarry Smith   PetscViewerFormat format;
54695fce210SBarry Smith 
54795fce210SBarry Smith   PetscFunctionBegin;
54895fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
54995fce210SBarry Smith   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sf),&viewer);CHKERRQ(ierr);}
55095fce210SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
55195fce210SBarry Smith   PetscCheckSameComm(sf,1,viewer,2);
55295fce210SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
55395fce210SBarry Smith   if (iascii) {
55495fce210SBarry Smith     PetscMPIInt rank;
55595fce210SBarry Smith     PetscInt    i,j;
55695fce210SBarry Smith 
557dae58748SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)sf,viewer);CHKERRQ(ierr);
55895fce210SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
55995fce210SBarry Smith     if (sf->ops->View) {ierr = (*sf->ops->View)(sf,viewer);CHKERRQ(ierr);}
56095fce210SBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
5611575c14dSBarry Smith     ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
56295fce210SBarry Smith     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of roots=%D, leaves=%D, remote ranks=%D\n",rank,sf->nroots,sf->nleaves,sf->nranks);CHKERRQ(ierr);
56395fce210SBarry Smith     for (i=0; i<sf->nleaves; i++) {
56495fce210SBarry 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);
56595fce210SBarry Smith     }
56695fce210SBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
56795fce210SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
56895fce210SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
56995fce210SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Roots referenced by my leaves, by rank\n",rank);CHKERRQ(ierr);
57095fce210SBarry Smith       for (i=0; i<sf->nranks; i++) {
5717904a332SBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %d: %D edges\n",rank,sf->ranks[i],sf->roffset[i+1]-sf->roffset[i]);CHKERRQ(ierr);
57295fce210SBarry Smith         for (j=sf->roffset[i]; j<sf->roffset[i+1]; j++) {
57395fce210SBarry Smith           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d]    %D <- %D\n",rank,sf->rmine[j],sf->rremote[j]);CHKERRQ(ierr);
57495fce210SBarry Smith         }
57595fce210SBarry Smith       }
57695fce210SBarry Smith     }
57795fce210SBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
5781575c14dSBarry Smith     ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
57995fce210SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
58095fce210SBarry Smith   }
58195fce210SBarry Smith   PetscFunctionReturn(0);
58295fce210SBarry Smith }
58395fce210SBarry Smith 
58495fce210SBarry Smith /*@C
58595fce210SBarry Smith    PetscSFGetRanks - Get ranks and number of vertices referenced by leaves on this process
58695fce210SBarry Smith 
58795fce210SBarry Smith    Not Collective
58895fce210SBarry Smith 
58995fce210SBarry Smith    Input Arguments:
59095fce210SBarry Smith .  sf - star forest
59195fce210SBarry Smith 
59295fce210SBarry Smith    Output Arguments:
59395fce210SBarry Smith +  nranks - number of ranks referenced by local part
59495fce210SBarry Smith .  ranks - array of ranks
59595fce210SBarry Smith .  roffset - offset in rmine/rremote for each rank (length nranks+1)
59695fce210SBarry Smith .  rmine - concatenated array holding local indices referencing each remote rank
59795fce210SBarry Smith -  rremote - concatenated array holding remote indices referenced for each remote rank
59895fce210SBarry Smith 
59995fce210SBarry Smith    Level: developer
60095fce210SBarry Smith 
60195fce210SBarry Smith .seealso: PetscSFSetGraph()
60295fce210SBarry Smith @*/
60395fce210SBarry Smith PetscErrorCode PetscSFGetRanks(PetscSF sf,PetscInt *nranks,const PetscMPIInt **ranks,const PetscInt **roffset,const PetscInt **rmine,const PetscInt **rremote)
60495fce210SBarry Smith {
60595fce210SBarry Smith 
60695fce210SBarry Smith   PetscFunctionBegin;
60795fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
60895fce210SBarry Smith   if (nranks)  *nranks  = sf->nranks;
60995fce210SBarry Smith   if (ranks)   *ranks   = sf->ranks;
61095fce210SBarry Smith   if (roffset) *roffset = sf->roffset;
61195fce210SBarry Smith   if (rmine)   *rmine   = sf->rmine;
61295fce210SBarry Smith   if (rremote) *rremote = sf->rremote;
61395fce210SBarry Smith   PetscFunctionReturn(0);
61495fce210SBarry Smith }
61595fce210SBarry Smith 
61695fce210SBarry Smith /*@C
61795fce210SBarry Smith    PetscSFGetGroups - gets incoming and outgoing process groups
61895fce210SBarry Smith 
61995fce210SBarry Smith    Collective
62095fce210SBarry Smith 
62195fce210SBarry Smith    Input Argument:
62295fce210SBarry Smith .  sf - star forest
62395fce210SBarry Smith 
62495fce210SBarry Smith    Output Arguments:
62595fce210SBarry Smith +  incoming - group of origin processes for incoming edges (leaves that reference my roots)
62695fce210SBarry Smith -  outgoing - group of destination processes for outgoing edges (roots that I reference)
62795fce210SBarry Smith 
62895fce210SBarry Smith    Level: developer
62995fce210SBarry Smith 
63095fce210SBarry Smith .seealso: PetscSFGetWindow(), PetscSFRestoreWindow()
63195fce210SBarry Smith @*/
63295fce210SBarry Smith PetscErrorCode PetscSFGetGroups(PetscSF sf,MPI_Group *incoming,MPI_Group *outgoing)
63395fce210SBarry Smith {
63495fce210SBarry Smith   PetscErrorCode ierr;
63595fce210SBarry Smith   MPI_Group      group;
63695fce210SBarry Smith 
63795fce210SBarry Smith   PetscFunctionBegin;
63895fce210SBarry Smith   if (sf->ingroup == MPI_GROUP_NULL) {
63995fce210SBarry Smith     PetscInt       i;
64095fce210SBarry Smith     const PetscInt *indegree;
64195fce210SBarry Smith     PetscMPIInt    rank,*outranks,*inranks;
64295fce210SBarry Smith     PetscSFNode    *remote;
64395fce210SBarry Smith     PetscSF        bgcount;
64495fce210SBarry Smith 
64595fce210SBarry Smith     /* Compute the number of incoming ranks */
646785e854fSJed Brown     ierr = PetscMalloc1(sf->nranks,&remote);CHKERRQ(ierr);
64795fce210SBarry Smith     for (i=0; i<sf->nranks; i++) {
64895fce210SBarry Smith       remote[i].rank  = sf->ranks[i];
64995fce210SBarry Smith       remote[i].index = 0;
65095fce210SBarry Smith     }
65195fce210SBarry Smith     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,&bgcount);CHKERRQ(ierr);
65295fce210SBarry Smith     ierr = PetscSFSetGraph(bgcount,1,sf->nranks,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
65395fce210SBarry Smith     ierr = PetscSFComputeDegreeBegin(bgcount,&indegree);CHKERRQ(ierr);
65495fce210SBarry Smith     ierr = PetscSFComputeDegreeEnd(bgcount,&indegree);CHKERRQ(ierr);
65595fce210SBarry Smith 
65695fce210SBarry Smith     /* Enumerate the incoming ranks */
657dcca6d9dSJed Brown     ierr = PetscMalloc2(indegree[0],&inranks,sf->nranks,&outranks);CHKERRQ(ierr);
65895fce210SBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
65995fce210SBarry Smith     for (i=0; i<sf->nranks; i++) outranks[i] = rank;
66095fce210SBarry Smith     ierr = PetscSFGatherBegin(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr);
66195fce210SBarry Smith     ierr = PetscSFGatherEnd(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr);
66295fce210SBarry Smith     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
66395fce210SBarry Smith     ierr = MPI_Group_incl(group,indegree[0],inranks,&sf->ingroup);CHKERRQ(ierr);
66495fce210SBarry Smith     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
66595fce210SBarry Smith     ierr = PetscFree2(inranks,outranks);CHKERRQ(ierr);
66695fce210SBarry Smith     ierr = PetscSFDestroy(&bgcount);CHKERRQ(ierr);
66795fce210SBarry Smith   }
66895fce210SBarry Smith   *incoming = sf->ingroup;
66995fce210SBarry Smith 
67095fce210SBarry Smith   if (sf->outgroup == MPI_GROUP_NULL) {
67195fce210SBarry Smith     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
67295fce210SBarry Smith     ierr = MPI_Group_incl(group,sf->nranks,sf->ranks,&sf->outgroup);CHKERRQ(ierr);
67395fce210SBarry Smith     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
67495fce210SBarry Smith   }
67595fce210SBarry Smith   *outgoing = sf->outgroup;
67695fce210SBarry Smith   PetscFunctionReturn(0);
67795fce210SBarry Smith }
67895fce210SBarry Smith 
67995fce210SBarry Smith /*@C
68095fce210SBarry Smith    PetscSFGetMultiSF - gets the inner SF implemeting gathers and scatters
68195fce210SBarry Smith 
68295fce210SBarry Smith    Collective
68395fce210SBarry Smith 
68495fce210SBarry Smith    Input Argument:
68595fce210SBarry Smith .  sf - star forest that may contain roots with 0 or with more than 1 vertex
68695fce210SBarry Smith 
68795fce210SBarry Smith    Output Arguments:
68895fce210SBarry Smith .  multi - star forest with split roots, such that each root has degree exactly 1
68995fce210SBarry Smith 
69095fce210SBarry Smith    Level: developer
69195fce210SBarry Smith 
69295fce210SBarry Smith    Notes:
69395fce210SBarry Smith 
69495fce210SBarry Smith    In most cases, users should use PetscSFGatherBegin() and PetscSFScatterBegin() instead of manipulating multi
69595fce210SBarry Smith    directly. Since multi satisfies the stronger condition that each entry in the global space has exactly one incoming
69695fce210SBarry Smith    edge, it is a candidate for future optimization that might involve its removal.
69795fce210SBarry Smith 
69895fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFGatherBegin(), PetscSFScatterBegin()
69995fce210SBarry Smith @*/
70095fce210SBarry Smith PetscErrorCode PetscSFGetMultiSF(PetscSF sf,PetscSF *multi)
70195fce210SBarry Smith {
70295fce210SBarry Smith   PetscErrorCode ierr;
70395fce210SBarry Smith 
70495fce210SBarry Smith   PetscFunctionBegin;
70595fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
70695fce210SBarry Smith   PetscValidPointer(multi,2);
70795fce210SBarry Smith   if (sf->nroots < 0) {         /* Graph has not been set yet; why do we need this? */
70895fce210SBarry Smith     ierr   = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr);
70995fce210SBarry Smith     *multi = sf->multi;
71095fce210SBarry Smith     PetscFunctionReturn(0);
71195fce210SBarry Smith   }
71295fce210SBarry Smith   if (!sf->multi) {
71395fce210SBarry Smith     const PetscInt *indegree;
7149837ea96SMatthew G. Knepley     PetscInt       i,*inoffset,*outones,*outoffset,maxlocal;
71595fce210SBarry Smith     PetscSFNode    *remote;
71695fce210SBarry Smith     ierr        = PetscSFComputeDegreeBegin(sf,&indegree);CHKERRQ(ierr);
71795fce210SBarry Smith     ierr        = PetscSFComputeDegreeEnd(sf,&indegree);CHKERRQ(ierr);
7189837ea96SMatthew G. Knepley     for (i=0,maxlocal=0; i<sf->nleaves; i++) maxlocal = PetscMax(maxlocal,(sf->mine ? sf->mine[i] : i)+1);
7199837ea96SMatthew G. Knepley     ierr        = PetscMalloc3(sf->nroots+1,&inoffset,maxlocal,&outones,maxlocal,&outoffset);CHKERRQ(ierr);
72095fce210SBarry Smith     inoffset[0] = 0;
72195fce210SBarry Smith     for (i=0; i<sf->nroots; i++) inoffset[i+1] = inoffset[i] + indegree[i];
7229837ea96SMatthew G. Knepley     for (i=0; i<maxlocal; i++) outones[i] = 1;
723dbd2ff41SBarry Smith     ierr = PetscSFFetchAndOpBegin(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr);
724dbd2ff41SBarry Smith     ierr = PetscSFFetchAndOpEnd(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr);
72595fce210SBarry Smith     for (i=0; i<sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */
72695fce210SBarry Smith #if 0
72795fce210SBarry Smith #if defined(PETSC_USE_DEBUG)                                 /* Check that the expected number of increments occurred */
72895fce210SBarry Smith     for (i=0; i<sf->nroots; i++) {
72995fce210SBarry Smith       if (inoffset[i] + indegree[i] != inoffset[i+1]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect result after PetscSFFetchAndOp");
73095fce210SBarry Smith     }
73195fce210SBarry Smith #endif
73295fce210SBarry Smith #endif
733785e854fSJed Brown     ierr = PetscMalloc1(sf->nleaves,&remote);CHKERRQ(ierr);
73495fce210SBarry Smith     for (i=0; i<sf->nleaves; i++) {
73595fce210SBarry Smith       remote[i].rank  = sf->remote[i].rank;
73638e7336fSToby Isaac       remote[i].index = outoffset[sf->mine ? sf->mine[i] : i];
73795fce210SBarry Smith     }
73895fce210SBarry Smith     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr);
73901365b40SToby Isaac     ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
74095fce210SBarry Smith     if (sf->rankorder) {        /* Sort the ranks */
74195fce210SBarry Smith       PetscMPIInt rank;
74295fce210SBarry Smith       PetscInt    *inranks,*newoffset,*outranks,*newoutoffset,*tmpoffset,maxdegree;
74395fce210SBarry Smith       PetscSFNode *newremote;
74495fce210SBarry Smith       ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
74595fce210SBarry Smith       for (i=0,maxdegree=0; i<sf->nroots; i++) maxdegree = PetscMax(maxdegree,indegree[i]);
7469837ea96SMatthew G. Knepley       ierr = PetscMalloc5(sf->multi->nroots,&inranks,sf->multi->nroots,&newoffset,maxlocal,&outranks,maxlocal,&newoutoffset,maxdegree,&tmpoffset);CHKERRQ(ierr);
7479837ea96SMatthew G. Knepley       for (i=0; i<maxlocal; i++) outranks[i] = rank;
7488bfbc91cSJed Brown       ierr = PetscSFReduceBegin(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr);
7498bfbc91cSJed Brown       ierr = PetscSFReduceEnd(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr);
75095fce210SBarry Smith       /* Sort the incoming ranks at each vertex, build the inverse map */
75195fce210SBarry Smith       for (i=0; i<sf->nroots; i++) {
75295fce210SBarry Smith         PetscInt j;
75395fce210SBarry Smith         for (j=0; j<indegree[i]; j++) tmpoffset[j] = j;
75495fce210SBarry Smith         ierr = PetscSortIntWithArray(indegree[i],inranks+inoffset[i],tmpoffset);CHKERRQ(ierr);
75595fce210SBarry Smith         for (j=0; j<indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j;
75695fce210SBarry Smith       }
75795fce210SBarry Smith       ierr = PetscSFBcastBegin(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr);
75895fce210SBarry Smith       ierr = PetscSFBcastEnd(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr);
759785e854fSJed Brown       ierr = PetscMalloc1(sf->nleaves,&newremote);CHKERRQ(ierr);
76095fce210SBarry Smith       for (i=0; i<sf->nleaves; i++) {
76195fce210SBarry Smith         newremote[i].rank  = sf->remote[i].rank;
76201365b40SToby Isaac         newremote[i].index = newoutoffset[sf->mine ? sf->mine[i] : i];
76395fce210SBarry Smith       }
76401365b40SToby Isaac       ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,newremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
76595fce210SBarry Smith       ierr = PetscFree5(inranks,newoffset,outranks,newoutoffset,tmpoffset);CHKERRQ(ierr);
76695fce210SBarry Smith     }
76795fce210SBarry Smith     ierr = PetscFree3(inoffset,outones,outoffset);CHKERRQ(ierr);
76895fce210SBarry Smith   }
76995fce210SBarry Smith   *multi = sf->multi;
77095fce210SBarry Smith   PetscFunctionReturn(0);
77195fce210SBarry Smith }
77295fce210SBarry Smith 
77395fce210SBarry Smith /*@C
77495fce210SBarry Smith    PetscSFCreateEmbeddedSF - removes edges from all but the selected roots, does not remap indices
77595fce210SBarry Smith 
77695fce210SBarry Smith    Collective
77795fce210SBarry Smith 
77895fce210SBarry Smith    Input Arguments:
77995fce210SBarry Smith +  sf - original star forest
78095fce210SBarry Smith .  nroots - number of roots to select on this process
78195fce210SBarry Smith -  selected - selected roots on this process
78295fce210SBarry Smith 
78395fce210SBarry Smith    Output Arguments:
78495fce210SBarry Smith .  newsf - new star forest
78595fce210SBarry Smith 
78695fce210SBarry Smith    Level: advanced
78795fce210SBarry Smith 
78895fce210SBarry Smith    Note:
78995fce210SBarry Smith    To use the new PetscSF, it may be necessary to know the indices of the leaves that are still participating. This can
79095fce210SBarry Smith    be done by calling PetscSFGetGraph().
79195fce210SBarry Smith 
79295fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFGetGraph()
79395fce210SBarry Smith @*/
79495fce210SBarry Smith PetscErrorCode PetscSFCreateEmbeddedSF(PetscSF sf,PetscInt nroots,const PetscInt *selected,PetscSF *newsf)
79595fce210SBarry Smith {
7960511a646SMatthew G. Knepley   PetscInt      *rootdata, *leafdata, *ilocal;
79795fce210SBarry Smith   PetscSFNode   *iremote;
798fc1ede2bSMatthew G. Knepley   PetscInt       leafsize = 0, nleaves = 0, n, i;
7990511a646SMatthew G. Knepley   PetscErrorCode ierr;
80095fce210SBarry Smith 
80195fce210SBarry Smith   PetscFunctionBegin;
80295fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
80395fce210SBarry Smith   if (nroots) PetscValidPointer(selected,3);
80495fce210SBarry Smith   PetscValidPointer(newsf,4);
8050511a646SMatthew G. Knepley   if (sf->mine) for (i = 0; i < sf->nleaves; ++i) {leafsize = PetscMax(leafsize, sf->mine[i]+1);}
8060511a646SMatthew G. Knepley   else leafsize = sf->nleaves;
8071795a4d1SJed Brown   ierr = PetscCalloc2(sf->nroots,&rootdata,leafsize,&leafdata);CHKERRQ(ierr);
8080511a646SMatthew G. Knepley   for (i=0; i<nroots; ++i) rootdata[selected[i]] = 1;
80995fce210SBarry Smith   ierr = PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr);
81095fce210SBarry Smith   ierr = PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr);
81195fce210SBarry Smith 
8120511a646SMatthew G. Knepley   for (i = 0; i < leafsize; ++i) nleaves += leafdata[i];
813785e854fSJed Brown   ierr = PetscMalloc1(nleaves,&ilocal);CHKERRQ(ierr);
814785e854fSJed Brown   ierr = PetscMalloc1(nleaves,&iremote);CHKERRQ(ierr);
815fc1ede2bSMatthew G. Knepley   for (i = 0, n = 0; i < sf->nleaves; ++i) {
8160511a646SMatthew G. Knepley     const PetscInt lidx = sf->mine ? sf->mine[i] : i;
8170511a646SMatthew G. Knepley 
8180511a646SMatthew G. Knepley     if (leafdata[lidx]) {
819fc1ede2bSMatthew G. Knepley       ilocal[n]        = lidx;
820fc1ede2bSMatthew G. Knepley       iremote[n].rank  = sf->remote[i].rank;
821fc1ede2bSMatthew G. Knepley       iremote[n].index = sf->remote[i].index;
822fc1ede2bSMatthew G. Knepley       ++n;
82395fce210SBarry Smith     }
82495fce210SBarry Smith   }
825fc1ede2bSMatthew 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);
82695fce210SBarry Smith   ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,newsf);CHKERRQ(ierr);
82795fce210SBarry Smith   ierr = PetscSFSetGraph(*newsf,sf->nroots,nleaves,ilocal,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
82895fce210SBarry Smith   ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr);
82995fce210SBarry Smith   PetscFunctionReturn(0);
83095fce210SBarry Smith }
83195fce210SBarry Smith 
8322f5fb4c2SMatthew G. Knepley /*@C
8332f5fb4c2SMatthew G. Knepley   PetscSFCreateEmbeddedLeafSF - removes edges from all but the selected leaves, does not remap indices
8342f5fb4c2SMatthew G. Knepley 
8352f5fb4c2SMatthew G. Knepley   Collective
8362f5fb4c2SMatthew G. Knepley 
8372f5fb4c2SMatthew G. Knepley   Input Arguments:
8382f5fb4c2SMatthew G. Knepley + sf - original star forest
8392f5fb4c2SMatthew G. Knepley . nleaves - number of leaves to select on this process
8402f5fb4c2SMatthew G. Knepley - selected - selected leaves on this process
8412f5fb4c2SMatthew G. Knepley 
8422f5fb4c2SMatthew G. Knepley   Output Arguments:
8432f5fb4c2SMatthew G. Knepley .  newsf - new star forest
8442f5fb4c2SMatthew G. Knepley 
8452f5fb4c2SMatthew G. Knepley   Level: advanced
8462f5fb4c2SMatthew G. Knepley 
8472f5fb4c2SMatthew G. Knepley .seealso: PetscSFCreateEmbeddedSF(), PetscSFSetGraph(), PetscSFGetGraph()
8482f5fb4c2SMatthew G. Knepley @*/
8492f5fb4c2SMatthew G. Knepley PetscErrorCode PetscSFCreateEmbeddedLeafSF(PetscSF sf, PetscInt nleaves, const PetscInt *selected, PetscSF *newsf)
8502f5fb4c2SMatthew G. Knepley {
8512f5fb4c2SMatthew G. Knepley   PetscSFNode   *iremote;
8522f5fb4c2SMatthew G. Knepley   PetscInt      *ilocal;
8532f5fb4c2SMatthew G. Knepley   PetscInt       i;
8542f5fb4c2SMatthew G. Knepley   PetscErrorCode ierr;
8552f5fb4c2SMatthew G. Knepley 
8562f5fb4c2SMatthew G. Knepley   PetscFunctionBegin;
8572f5fb4c2SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
8582f5fb4c2SMatthew G. Knepley   if (nleaves) PetscValidPointer(selected, 3);
8592f5fb4c2SMatthew G. Knepley   PetscValidPointer(newsf, 4);
8602f5fb4c2SMatthew G. Knepley   ierr = PetscMalloc1(nleaves, &ilocal);CHKERRQ(ierr);
8612f5fb4c2SMatthew G. Knepley   ierr = PetscMalloc1(nleaves, &iremote);CHKERRQ(ierr);
8622f5fb4c2SMatthew G. Knepley   for (i = 0; i < nleaves; ++i) {
8632f5fb4c2SMatthew G. Knepley     const PetscInt l = selected[i];
8642f5fb4c2SMatthew G. Knepley 
8652f5fb4c2SMatthew G. Knepley     ilocal[i]        = sf->mine ? sf->mine[l] : l;
8662f5fb4c2SMatthew G. Knepley     iremote[i].rank  = sf->remote[l].rank;
8672f5fb4c2SMatthew G. Knepley     iremote[i].index = sf->remote[l].index;
8682f5fb4c2SMatthew G. Knepley   }
8692f5fb4c2SMatthew G. Knepley   ierr = PetscSFDuplicate(sf, PETSCSF_DUPLICATE_RANKS, newsf);CHKERRQ(ierr);
8702f5fb4c2SMatthew G. Knepley   ierr = PetscSFSetGraph(*newsf, sf->nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr);
8712f5fb4c2SMatthew G. Knepley   PetscFunctionReturn(0);
8722f5fb4c2SMatthew G. Knepley }
8732f5fb4c2SMatthew G. Knepley 
87495fce210SBarry Smith /*@C
87595fce210SBarry Smith    PetscSFBcastBegin - begin pointwise broadcast to be concluded with call to PetscSFBcastEnd()
87695fce210SBarry Smith 
87795fce210SBarry Smith    Collective on PetscSF
87895fce210SBarry Smith 
87995fce210SBarry Smith    Input Arguments:
88095fce210SBarry Smith +  sf - star forest on which to communicate
88195fce210SBarry Smith .  unit - data type associated with each node
88295fce210SBarry Smith -  rootdata - buffer to broadcast
88395fce210SBarry Smith 
88495fce210SBarry Smith    Output Arguments:
88595fce210SBarry Smith .  leafdata - buffer to update with values from each leaf's respective root
88695fce210SBarry Smith 
88795fce210SBarry Smith    Level: intermediate
88895fce210SBarry Smith 
88995fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFView(), PetscSFBcastEnd(), PetscSFReduceBegin()
89095fce210SBarry Smith @*/
89195fce210SBarry Smith PetscErrorCode PetscSFBcastBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
89295fce210SBarry Smith {
89395fce210SBarry Smith   PetscErrorCode ierr;
89495fce210SBarry Smith 
89595fce210SBarry Smith   PetscFunctionBegin;
89695fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
89795fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
898563ffbabSMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr);
89995fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
90095fce210SBarry Smith   ierr = (*sf->ops->BcastBegin)(sf,unit,rootdata,leafdata);CHKERRQ(ierr);
901563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr);
90295fce210SBarry Smith   PetscFunctionReturn(0);
90395fce210SBarry Smith }
90495fce210SBarry Smith 
90595fce210SBarry Smith /*@C
90695fce210SBarry Smith    PetscSFBcastEnd - end a broadcast operation started with PetscSFBcastBegin()
90795fce210SBarry Smith 
90895fce210SBarry Smith    Collective
90995fce210SBarry Smith 
91095fce210SBarry Smith    Input Arguments:
91195fce210SBarry Smith +  sf - star forest
91295fce210SBarry Smith .  unit - data type
91395fce210SBarry Smith -  rootdata - buffer to broadcast
91495fce210SBarry Smith 
91595fce210SBarry Smith    Output Arguments:
91695fce210SBarry Smith .  leafdata - buffer to update with values from each leaf's respective root
91795fce210SBarry Smith 
91895fce210SBarry Smith    Level: intermediate
91995fce210SBarry Smith 
92095fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFReduceEnd()
92195fce210SBarry Smith @*/
92295fce210SBarry Smith PetscErrorCode PetscSFBcastEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
92395fce210SBarry Smith {
92495fce210SBarry Smith   PetscErrorCode ierr;
92595fce210SBarry Smith 
92695fce210SBarry Smith   PetscFunctionBegin;
92795fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
92895fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
929563ffbabSMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCSF_BcastEnd,sf,0,0,0);CHKERRQ(ierr);
93095fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
93195fce210SBarry Smith   ierr = (*sf->ops->BcastEnd)(sf,unit,rootdata,leafdata);CHKERRQ(ierr);
932563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_BcastEnd,sf,0,0,0);CHKERRQ(ierr);
93395fce210SBarry Smith   PetscFunctionReturn(0);
93495fce210SBarry Smith }
93595fce210SBarry Smith 
93695fce210SBarry Smith /*@C
93795fce210SBarry Smith    PetscSFReduceBegin - begin reduction of leafdata into rootdata, to be completed with call to PetscSFReduceEnd()
93895fce210SBarry Smith 
93995fce210SBarry Smith    Collective
94095fce210SBarry Smith 
94195fce210SBarry Smith    Input Arguments:
94295fce210SBarry Smith +  sf - star forest
94395fce210SBarry Smith .  unit - data type
94495fce210SBarry Smith .  leafdata - values to reduce
94595fce210SBarry Smith -  op - reduction operation
94695fce210SBarry Smith 
94795fce210SBarry Smith    Output Arguments:
94895fce210SBarry Smith .  rootdata - result of reduction of values from all leaves of each root
94995fce210SBarry Smith 
95095fce210SBarry Smith    Level: intermediate
95195fce210SBarry Smith 
95295fce210SBarry Smith .seealso: PetscSFBcastBegin()
95395fce210SBarry Smith @*/
95495fce210SBarry Smith PetscErrorCode PetscSFReduceBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
95595fce210SBarry Smith {
95695fce210SBarry Smith   PetscErrorCode ierr;
95795fce210SBarry Smith 
95895fce210SBarry Smith   PetscFunctionBegin;
95995fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
96095fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
961563ffbabSMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);
96295fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
96395fce210SBarry Smith   ierr = (sf->ops->ReduceBegin)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
964563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);
96595fce210SBarry Smith   PetscFunctionReturn(0);
96695fce210SBarry Smith }
96795fce210SBarry Smith 
96895fce210SBarry Smith /*@C
96995fce210SBarry Smith    PetscSFReduceEnd - end a reduction operation started with PetscSFReduceBegin()
97095fce210SBarry Smith 
97195fce210SBarry Smith    Collective
97295fce210SBarry Smith 
97395fce210SBarry Smith    Input Arguments:
97495fce210SBarry Smith +  sf - star forest
97595fce210SBarry Smith .  unit - data type
97695fce210SBarry Smith .  leafdata - values to reduce
97795fce210SBarry Smith -  op - reduction operation
97895fce210SBarry Smith 
97995fce210SBarry Smith    Output Arguments:
98095fce210SBarry Smith .  rootdata - result of reduction of values from all leaves of each root
98195fce210SBarry Smith 
98295fce210SBarry Smith    Level: intermediate
98395fce210SBarry Smith 
98495fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFBcastEnd()
98595fce210SBarry Smith @*/
98695fce210SBarry Smith PetscErrorCode PetscSFReduceEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
98795fce210SBarry Smith {
98895fce210SBarry Smith   PetscErrorCode ierr;
98995fce210SBarry Smith 
99095fce210SBarry Smith   PetscFunctionBegin;
99195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
99295fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
993563ffbabSMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);
99495fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
99595fce210SBarry Smith   ierr = (*sf->ops->ReduceEnd)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
996563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);
99795fce210SBarry Smith   PetscFunctionReturn(0);
99895fce210SBarry Smith }
99995fce210SBarry Smith 
100095fce210SBarry Smith /*@C
100195fce210SBarry Smith    PetscSFComputeDegreeBegin - begin computation of degree for each root vertex, to be completed with PetscSFComputeDegreeEnd()
100295fce210SBarry Smith 
100395fce210SBarry Smith    Collective
100495fce210SBarry Smith 
100595fce210SBarry Smith    Input Arguments:
100695fce210SBarry Smith .  sf - star forest
100795fce210SBarry Smith 
100895fce210SBarry Smith    Output Arguments:
100995fce210SBarry Smith .  degree - degree of each root vertex
101095fce210SBarry Smith 
101195fce210SBarry Smith    Level: advanced
101295fce210SBarry Smith 
101395fce210SBarry Smith .seealso: PetscSFGatherBegin()
101495fce210SBarry Smith @*/
101595fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf,const PetscInt **degree)
101695fce210SBarry Smith {
101795fce210SBarry Smith   PetscErrorCode ierr;
101895fce210SBarry Smith 
101995fce210SBarry Smith   PetscFunctionBegin;
102095fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
102195fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
102295fce210SBarry Smith   PetscValidPointer(degree,2);
1023803bd9e8SMatthew G. Knepley   if (!sf->degreeknown) {
10249837ea96SMatthew G. Knepley     PetscInt i,maxlocal;
1025803bd9e8SMatthew G. Knepley     if (sf->degree) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Calls to PetscSFComputeDegreeBegin() cannot be nested.");
10269837ea96SMatthew G. Knepley     for (i=0,maxlocal=0; i<sf->nleaves; i++) maxlocal = PetscMax(maxlocal,(sf->mine ? sf->mine[i] : i)+1);
1027785e854fSJed Brown     ierr = PetscMalloc1(sf->nroots,&sf->degree);CHKERRQ(ierr);
10289837ea96SMatthew G. Knepley     ierr = PetscMalloc1(maxlocal,&sf->degreetmp);CHKERRQ(ierr);
102995fce210SBarry Smith     for (i=0; i<sf->nroots; i++) sf->degree[i] = 0;
10309837ea96SMatthew G. Knepley     for (i=0; i<maxlocal; i++) sf->degreetmp[i] = 1;
1031dbd2ff41SBarry Smith     ierr = PetscSFReduceBegin(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);CHKERRQ(ierr);
103295fce210SBarry Smith   }
103395fce210SBarry Smith   *degree = NULL;
103495fce210SBarry Smith   PetscFunctionReturn(0);
103595fce210SBarry Smith }
103695fce210SBarry Smith 
103795fce210SBarry Smith /*@C
103895fce210SBarry Smith    PetscSFComputeDegreeEnd - complete computation of degree for each root vertex, started with PetscSFComputeDegreeBegin()
103995fce210SBarry Smith 
104095fce210SBarry Smith    Collective
104195fce210SBarry Smith 
104295fce210SBarry Smith    Input Arguments:
104395fce210SBarry Smith .  sf - star forest
104495fce210SBarry Smith 
104595fce210SBarry Smith    Output Arguments:
104695fce210SBarry Smith .  degree - degree of each root vertex
104795fce210SBarry Smith 
104895fce210SBarry Smith    Level: developer
104995fce210SBarry Smith 
105095fce210SBarry Smith .seealso:
105195fce210SBarry Smith @*/
105295fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf,const PetscInt **degree)
105395fce210SBarry Smith {
105495fce210SBarry Smith   PetscErrorCode ierr;
105595fce210SBarry Smith 
105695fce210SBarry Smith   PetscFunctionBegin;
105795fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
105895fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
105995fce210SBarry Smith   if (!sf->degreeknown) {
1060dbd2ff41SBarry Smith     ierr = PetscSFReduceEnd(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);CHKERRQ(ierr);
106195fce210SBarry Smith     ierr = PetscFree(sf->degreetmp);CHKERRQ(ierr);
106295fce210SBarry Smith 
106395fce210SBarry Smith     sf->degreeknown = PETSC_TRUE;
106495fce210SBarry Smith   }
106595fce210SBarry Smith   *degree = sf->degree;
106695fce210SBarry Smith   PetscFunctionReturn(0);
106795fce210SBarry Smith }
106895fce210SBarry Smith 
106995fce210SBarry Smith /*@C
107095fce210SBarry Smith    PetscSFFetchAndOpBegin - begin operation that fetches values from root and updates atomically by applying operation using my leaf value, to be completed with PetscSFFetchAndOpEnd()
107195fce210SBarry Smith 
107295fce210SBarry Smith    Collective
107395fce210SBarry Smith 
107495fce210SBarry Smith    Input Arguments:
107595fce210SBarry Smith +  sf - star forest
107695fce210SBarry Smith .  unit - data type
107795fce210SBarry Smith .  leafdata - leaf values to use in reduction
107895fce210SBarry Smith -  op - operation to use for reduction
107995fce210SBarry Smith 
108095fce210SBarry Smith    Output Arguments:
108195fce210SBarry Smith +  rootdata - root values to be updated, input state is seen by first process to perform an update
108295fce210SBarry Smith -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
108395fce210SBarry Smith 
108495fce210SBarry Smith    Level: advanced
108595fce210SBarry Smith 
108695fce210SBarry Smith    Note:
108795fce210SBarry Smith    The update is only atomic at the granularity provided by the hardware. Different roots referenced by the same process
108895fce210SBarry Smith    might be updated in a different order. Furthermore, if a composite type is used for the unit datatype, atomicity is
108995fce210SBarry Smith    not guaranteed across the whole vertex. Therefore, this function is mostly only used with primitive types such as
109095fce210SBarry Smith    integers.
109195fce210SBarry Smith 
109295fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFReduceBegin(), PetscSFSetGraph()
109395fce210SBarry Smith @*/
109495fce210SBarry Smith PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
109595fce210SBarry Smith {
109695fce210SBarry Smith   PetscErrorCode ierr;
109795fce210SBarry Smith 
109895fce210SBarry Smith   PetscFunctionBegin;
109995fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
110095fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
1101563ffbabSMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
110295fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
110395fce210SBarry Smith   ierr = (*sf->ops->FetchAndOpBegin)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr);
1104563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
110595fce210SBarry Smith   PetscFunctionReturn(0);
110695fce210SBarry Smith }
110795fce210SBarry Smith 
110895fce210SBarry Smith /*@C
110995fce210SBarry 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
111095fce210SBarry Smith 
111195fce210SBarry Smith    Collective
111295fce210SBarry Smith 
111395fce210SBarry Smith    Input Arguments:
111495fce210SBarry Smith +  sf - star forest
111595fce210SBarry Smith .  unit - data type
111695fce210SBarry Smith .  leafdata - leaf values to use in reduction
111795fce210SBarry Smith -  op - operation to use for reduction
111895fce210SBarry Smith 
111995fce210SBarry Smith    Output Arguments:
112095fce210SBarry Smith +  rootdata - root values to be updated, input state is seen by first process to perform an update
112195fce210SBarry Smith -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
112295fce210SBarry Smith 
112395fce210SBarry Smith    Level: advanced
112495fce210SBarry Smith 
112595fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFReduceEnd(), PetscSFSetGraph()
112695fce210SBarry Smith @*/
112795fce210SBarry Smith PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
112895fce210SBarry Smith {
112995fce210SBarry Smith   PetscErrorCode ierr;
113095fce210SBarry Smith 
113195fce210SBarry Smith   PetscFunctionBegin;
113295fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
113395fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
1134563ffbabSMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
113595fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
113695fce210SBarry Smith   ierr = (*sf->ops->FetchAndOpEnd)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr);
1137563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
113895fce210SBarry Smith   PetscFunctionReturn(0);
113995fce210SBarry Smith }
114095fce210SBarry Smith 
114195fce210SBarry Smith /*@C
114295fce210SBarry Smith    PetscSFGatherBegin - begin pointwise gather of all leaves into multi-roots, to be completed with PetscSFGatherEnd()
114395fce210SBarry Smith 
114495fce210SBarry Smith    Collective
114595fce210SBarry Smith 
114695fce210SBarry Smith    Input Arguments:
114795fce210SBarry Smith +  sf - star forest
114895fce210SBarry Smith .  unit - data type
114995fce210SBarry Smith -  leafdata - leaf data to gather to roots
115095fce210SBarry Smith 
115195fce210SBarry Smith    Output Argument:
115295fce210SBarry Smith .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
115395fce210SBarry Smith 
115495fce210SBarry Smith    Level: intermediate
115595fce210SBarry Smith 
115695fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
115795fce210SBarry Smith @*/
115895fce210SBarry Smith PetscErrorCode PetscSFGatherBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
115995fce210SBarry Smith {
116095fce210SBarry Smith   PetscErrorCode ierr;
116195fce210SBarry Smith   PetscSF        multi;
116295fce210SBarry Smith 
116395fce210SBarry Smith   PetscFunctionBegin;
116495fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
116595fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
11668bfbc91cSJed Brown   ierr = PetscSFReduceBegin(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr);
116795fce210SBarry Smith   PetscFunctionReturn(0);
116895fce210SBarry Smith }
116995fce210SBarry Smith 
117095fce210SBarry Smith /*@C
117195fce210SBarry Smith    PetscSFGatherEnd - ends pointwise gather operation that was started with PetscSFGatherBegin()
117295fce210SBarry Smith 
117395fce210SBarry Smith    Collective
117495fce210SBarry Smith 
117595fce210SBarry Smith    Input Arguments:
117695fce210SBarry Smith +  sf - star forest
117795fce210SBarry Smith .  unit - data type
117895fce210SBarry Smith -  leafdata - leaf data to gather to roots
117995fce210SBarry Smith 
118095fce210SBarry Smith    Output Argument:
118195fce210SBarry Smith .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
118295fce210SBarry Smith 
118395fce210SBarry Smith    Level: intermediate
118495fce210SBarry Smith 
118595fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
118695fce210SBarry Smith @*/
118795fce210SBarry Smith PetscErrorCode PetscSFGatherEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
118895fce210SBarry Smith {
118995fce210SBarry Smith   PetscErrorCode ierr;
119095fce210SBarry Smith   PetscSF        multi;
119195fce210SBarry Smith 
119295fce210SBarry Smith   PetscFunctionBegin;
119395fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
119495fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
119595fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
119695fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
11978bfbc91cSJed Brown   ierr = PetscSFReduceEnd(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr);
119895fce210SBarry Smith   PetscFunctionReturn(0);
119995fce210SBarry Smith }
120095fce210SBarry Smith 
120195fce210SBarry Smith /*@C
120295fce210SBarry Smith    PetscSFScatterBegin - begin pointwise scatter operation from multi-roots to leaves, to be completed with PetscSFScatterEnd()
120395fce210SBarry Smith 
120495fce210SBarry Smith    Collective
120595fce210SBarry Smith 
120695fce210SBarry Smith    Input Arguments:
120795fce210SBarry Smith +  sf - star forest
120895fce210SBarry Smith .  unit - data type
120995fce210SBarry Smith -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
121095fce210SBarry Smith 
121195fce210SBarry Smith    Output Argument:
121295fce210SBarry Smith .  leafdata - leaf data to be update with personal data from each respective root
121395fce210SBarry Smith 
121495fce210SBarry Smith    Level: intermediate
121595fce210SBarry Smith 
121695fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
121795fce210SBarry Smith @*/
121895fce210SBarry Smith PetscErrorCode PetscSFScatterBegin(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
121995fce210SBarry Smith {
122095fce210SBarry Smith   PetscErrorCode ierr;
122195fce210SBarry Smith   PetscSF        multi;
122295fce210SBarry Smith 
122395fce210SBarry Smith   PetscFunctionBegin;
122495fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
122595fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
122695fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
122795fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
122895fce210SBarry Smith   ierr = PetscSFBcastBegin(multi,unit,multirootdata,leafdata);CHKERRQ(ierr);
122995fce210SBarry Smith   PetscFunctionReturn(0);
123095fce210SBarry Smith }
123195fce210SBarry Smith 
123295fce210SBarry Smith /*@C
123395fce210SBarry Smith    PetscSFScatterEnd - ends pointwise scatter operation that was started with PetscSFScatterBegin()
123495fce210SBarry Smith 
123595fce210SBarry Smith    Collective
123695fce210SBarry Smith 
123795fce210SBarry Smith    Input Arguments:
123895fce210SBarry Smith +  sf - star forest
123995fce210SBarry Smith .  unit - data type
124095fce210SBarry Smith -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
124195fce210SBarry Smith 
124295fce210SBarry Smith    Output Argument:
124395fce210SBarry Smith .  leafdata - leaf data to be update with personal data from each respective root
124495fce210SBarry Smith 
124595fce210SBarry Smith    Level: intermediate
124695fce210SBarry Smith 
124795fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
124895fce210SBarry Smith @*/
124995fce210SBarry Smith PetscErrorCode PetscSFScatterEnd(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
125095fce210SBarry Smith {
125195fce210SBarry Smith   PetscErrorCode ierr;
125295fce210SBarry Smith   PetscSF        multi;
125395fce210SBarry Smith 
125495fce210SBarry Smith   PetscFunctionBegin;
125595fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
125695fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
125795fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
125895fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
125995fce210SBarry Smith   ierr = PetscSFBcastEnd(multi,unit,multirootdata,leafdata);CHKERRQ(ierr);
126095fce210SBarry Smith   PetscFunctionReturn(0);
126195fce210SBarry Smith }
1262a7b3aa13SAta Mesgarnejad 
1263a7b3aa13SAta Mesgarnejad /*@
1264a7b3aa13SAta Mesgarnejad   PetscSFCompose - Compose a new PetscSF equivalent to action to PetscSFs
1265a7b3aa13SAta Mesgarnejad 
1266a7b3aa13SAta Mesgarnejad   Input Parameters:
1267a7b3aa13SAta Mesgarnejad + sfA - The first PetscSF
1268a7b3aa13SAta Mesgarnejad - sfB - The second PetscSF
1269a7b3aa13SAta Mesgarnejad 
1270a7b3aa13SAta Mesgarnejad   Output Parameters:
1271a7b3aa13SAta Mesgarnejad . sfBA - equvalent PetscSF for applying A then B
1272a7b3aa13SAta Mesgarnejad 
1273a7b3aa13SAta Mesgarnejad   Level: developer
1274a7b3aa13SAta Mesgarnejad 
1275a7b3aa13SAta Mesgarnejad .seealso: PetscSF, PetscSFGetGraph(), PetscSFSetGraph()
1276a7b3aa13SAta Mesgarnejad @*/
1277a7b3aa13SAta Mesgarnejad PetscErrorCode PetscSFCompose(PetscSF sfA, PetscSF sfB, PetscSF *sfBA)
1278a7b3aa13SAta Mesgarnejad {
1279a7b3aa13SAta Mesgarnejad   MPI_Comm           comm;
1280a7b3aa13SAta Mesgarnejad   const PetscSFNode *remotePointsA, *remotePointsB;
1281a7b3aa13SAta Mesgarnejad   PetscSFNode       *remotePointsBA;
1282a7b3aa13SAta Mesgarnejad   const PetscInt    *localPointsA, *localPointsB;
1283a7b3aa13SAta Mesgarnejad   PetscInt           numRootsA, numLeavesA, numRootsB, numLeavesB;
1284a7b3aa13SAta Mesgarnejad   PetscErrorCode     ierr;
1285a7b3aa13SAta Mesgarnejad 
1286a7b3aa13SAta Mesgarnejad   PetscFunctionBegin;
1287a7b3aa13SAta Mesgarnejad   PetscValidHeaderSpecific(sfA, PETSCSF_CLASSID, 1);
1288a7b3aa13SAta Mesgarnejad   PetscValidHeaderSpecific(sfB, PETSCSF_CLASSID, 1);
1289a7b3aa13SAta Mesgarnejad   ierr = PetscObjectGetComm((PetscObject) sfA, &comm);CHKERRQ(ierr);
1290a7b3aa13SAta Mesgarnejad   ierr = PetscSFGetGraph(sfA, &numRootsA, &numLeavesA, &localPointsA, &remotePointsA);CHKERRQ(ierr);
1291a7b3aa13SAta Mesgarnejad   ierr = PetscSFGetGraph(sfB, &numRootsB, &numLeavesB, &localPointsB, &remotePointsB);CHKERRQ(ierr);
1292a7b3aa13SAta Mesgarnejad   ierr = PetscMalloc1(numLeavesB, &remotePointsBA);CHKERRQ(ierr);
1293a7b3aa13SAta Mesgarnejad   ierr = PetscSFBcastBegin(sfB, MPIU_2INT, remotePointsA, remotePointsBA);CHKERRQ(ierr);
1294a7b3aa13SAta Mesgarnejad   ierr = PetscSFBcastEnd(sfB, MPIU_2INT, remotePointsA, remotePointsBA);CHKERRQ(ierr);
1295a7b3aa13SAta Mesgarnejad   ierr = PetscSFCreate(comm, sfBA);CHKERRQ(ierr);
1296a7b3aa13SAta Mesgarnejad   ierr = PetscSFSetGraph(*sfBA, numRootsA, numLeavesB, localPointsB, PETSC_COPY_VALUES, remotePointsBA, PETSC_OWN_POINTER);CHKERRQ(ierr);
1297a7b3aa13SAta Mesgarnejad   PetscFunctionReturn(0);
1298a7b3aa13SAta Mesgarnejad }
1299