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); 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 18895fce210SBarry Smith /*@C 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 24995fce210SBarry Smith /*@C 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 26595fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFGetGraph() 26695fce210SBarry Smith @*/ 26795fce210SBarry Smith PetscErrorCode PetscSFSetGraph(PetscSF sf,PetscInt nroots,PetscInt nleaves,const PetscInt *ilocal,PetscCopyMode localmode,const PetscSFNode *iremote,PetscCopyMode remotemode) 26895fce210SBarry Smith { 26995fce210SBarry Smith PetscErrorCode ierr; 27095fce210SBarry Smith 27195fce210SBarry Smith PetscFunctionBegin; 27295fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 273563ffbabSMatthew G. Knepley ierr = PetscLogEventBegin(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr); 27495fce210SBarry Smith if (nleaves && ilocal) PetscValidIntPointer(ilocal,4); 27595fce210SBarry Smith if (nleaves) PetscValidPointer(iremote,6); 27695fce210SBarry Smith if (nroots < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"roots %D, cannot be negative",nroots); 27795fce210SBarry Smith if (nleaves < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nleaves %D, cannot be negative",nleaves); 27895fce210SBarry Smith ierr = PetscSFReset(sf);CHKERRQ(ierr); 27995fce210SBarry Smith sf->nroots = nroots; 28095fce210SBarry Smith sf->nleaves = nleaves; 28195fce210SBarry Smith if (ilocal) { 28221c688dcSJed Brown PetscInt i; 28395fce210SBarry Smith switch (localmode) { 28495fce210SBarry Smith case PETSC_COPY_VALUES: 285785e854fSJed Brown ierr = PetscMalloc1(nleaves,&sf->mine_alloc);CHKERRQ(ierr); 28695fce210SBarry Smith sf->mine = sf->mine_alloc; 28795fce210SBarry Smith ierr = PetscMemcpy(sf->mine,ilocal,nleaves*sizeof(*sf->mine));CHKERRQ(ierr); 28895fce210SBarry Smith sf->minleaf = PETSC_MAX_INT; 28995fce210SBarry Smith sf->maxleaf = PETSC_MIN_INT; 29095fce210SBarry Smith for (i=0; i<nleaves; i++) { 29195fce210SBarry Smith sf->minleaf = PetscMin(sf->minleaf,ilocal[i]); 29295fce210SBarry Smith sf->maxleaf = PetscMax(sf->maxleaf,ilocal[i]); 29395fce210SBarry Smith } 29495fce210SBarry Smith break; 29595fce210SBarry Smith case PETSC_OWN_POINTER: 29695fce210SBarry Smith sf->mine_alloc = (PetscInt*)ilocal; 29795fce210SBarry Smith sf->mine = sf->mine_alloc; 29895fce210SBarry Smith break; 29995fce210SBarry Smith case PETSC_USE_POINTER: 30095fce210SBarry Smith sf->mine = (PetscInt*)ilocal; 30195fce210SBarry Smith break; 30295fce210SBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown localmode"); 30395fce210SBarry Smith } 30495fce210SBarry Smith } 30595fce210SBarry Smith if (!ilocal || nleaves > 0) { 30695fce210SBarry Smith sf->minleaf = 0; 30795fce210SBarry Smith sf->maxleaf = nleaves - 1; 30895fce210SBarry Smith } 30995fce210SBarry Smith switch (remotemode) { 31095fce210SBarry Smith case PETSC_COPY_VALUES: 311785e854fSJed Brown ierr = PetscMalloc1(nleaves,&sf->remote_alloc);CHKERRQ(ierr); 31295fce210SBarry Smith sf->remote = sf->remote_alloc; 31395fce210SBarry Smith ierr = PetscMemcpy(sf->remote,iremote,nleaves*sizeof(*sf->remote));CHKERRQ(ierr); 31495fce210SBarry Smith break; 31595fce210SBarry Smith case PETSC_OWN_POINTER: 31695fce210SBarry Smith sf->remote_alloc = (PetscSFNode*)iremote; 31795fce210SBarry Smith sf->remote = sf->remote_alloc; 31895fce210SBarry Smith break; 31995fce210SBarry Smith case PETSC_USE_POINTER: 32095fce210SBarry Smith sf->remote = (PetscSFNode*)iremote; 32195fce210SBarry Smith break; 32295fce210SBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown remotemode"); 32395fce210SBarry Smith } 32495fce210SBarry Smith 32595fce210SBarry Smith sf->graphset = PETSC_TRUE; 326563ffbabSMatthew G. Knepley ierr = PetscLogEventEnd(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr); 32795fce210SBarry Smith PetscFunctionReturn(0); 32895fce210SBarry Smith } 32995fce210SBarry Smith 33095fce210SBarry Smith /*@C 33195fce210SBarry Smith PetscSFCreateInverseSF - given a PetscSF in which all vertices have degree 1, creates the inverse map 33295fce210SBarry Smith 33395fce210SBarry Smith Collective 33495fce210SBarry Smith 33595fce210SBarry Smith Input Arguments: 33695fce210SBarry Smith . sf - star forest to invert 33795fce210SBarry Smith 33895fce210SBarry Smith Output Arguments: 33995fce210SBarry Smith . isf - inverse of sf 34095fce210SBarry Smith 34195fce210SBarry Smith Level: advanced 34295fce210SBarry Smith 34395fce210SBarry Smith Notes: 34495fce210SBarry Smith All roots must have degree 1. 34595fce210SBarry Smith 34695fce210SBarry Smith The local space may be a permutation, but cannot be sparse. 34795fce210SBarry Smith 34895fce210SBarry Smith .seealso: PetscSFSetGraph() 34995fce210SBarry Smith @*/ 35095fce210SBarry Smith PetscErrorCode PetscSFCreateInverseSF(PetscSF sf,PetscSF *isf) 35195fce210SBarry Smith { 35295fce210SBarry Smith PetscErrorCode ierr; 35395fce210SBarry Smith PetscMPIInt rank; 35495fce210SBarry Smith PetscInt i,nroots,nleaves,maxlocal,count,*newilocal; 35595fce210SBarry Smith const PetscInt *ilocal; 35695fce210SBarry Smith PetscSFNode *roots,*leaves; 35795fce210SBarry Smith 35895fce210SBarry Smith PetscFunctionBegin; 35995fce210SBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); 36095fce210SBarry Smith ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);CHKERRQ(ierr); 36195fce210SBarry Smith for (i=0,maxlocal=0; i<nleaves; i++) maxlocal = PetscMax(maxlocal,(ilocal ? ilocal[i] : i)+1); 362ae9aee6dSMatthew G. Knepley ierr = PetscMalloc2(nroots,&roots,maxlocal,&leaves);CHKERRQ(ierr); 363ae9aee6dSMatthew G. Knepley for (i=0; i<maxlocal; i++) { 36495fce210SBarry Smith leaves[i].rank = rank; 36595fce210SBarry Smith leaves[i].index = i; 36695fce210SBarry Smith } 36795fce210SBarry Smith for (i=0; i <nroots; i++) { 36895fce210SBarry Smith roots[i].rank = -1; 36995fce210SBarry Smith roots[i].index = -1; 37095fce210SBarry Smith } 3718bfbc91cSJed Brown ierr = PetscSFReduceBegin(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr); 3728bfbc91cSJed Brown ierr = PetscSFReduceEnd(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr); 37395fce210SBarry Smith 37495fce210SBarry Smith /* Check whether our leaves are sparse */ 37595fce210SBarry Smith for (i=0,count=0; i<nroots; i++) if (roots[i].rank >= 0) count++; 37695fce210SBarry Smith if (count == nroots) newilocal = NULL; 37795fce210SBarry Smith else { /* Index for sparse leaves and compact "roots" array (which is to become our leaves). */ 378785e854fSJed Brown ierr = PetscMalloc1(count,&newilocal);CHKERRQ(ierr); 37995fce210SBarry Smith for (i=0,count=0; i<nroots; i++) { 38095fce210SBarry Smith if (roots[i].rank >= 0) { 38195fce210SBarry Smith newilocal[count] = i; 38295fce210SBarry Smith roots[count].rank = roots[i].rank; 38395fce210SBarry Smith roots[count].index = roots[i].index; 38495fce210SBarry Smith count++; 38595fce210SBarry Smith } 38695fce210SBarry Smith } 38795fce210SBarry Smith } 38895fce210SBarry Smith 38995fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,isf);CHKERRQ(ierr); 39095fce210SBarry Smith ierr = PetscSFSetGraph(*isf,maxlocal,count,newilocal,PETSC_OWN_POINTER,roots,PETSC_COPY_VALUES);CHKERRQ(ierr); 39195fce210SBarry Smith ierr = PetscFree2(roots,leaves);CHKERRQ(ierr); 39295fce210SBarry Smith PetscFunctionReturn(0); 39395fce210SBarry Smith } 39495fce210SBarry Smith 39595fce210SBarry Smith /*@ 39695fce210SBarry Smith PetscSFDuplicate - duplicate a PetscSF, optionally preserving rank connectivity and graph 39795fce210SBarry Smith 39895fce210SBarry Smith Collective 39995fce210SBarry Smith 40095fce210SBarry Smith Input Arguments: 40195fce210SBarry Smith + sf - communication object to duplicate 40295fce210SBarry Smith - opt - PETSCSF_DUPLICATE_CONFONLY, PETSCSF_DUPLICATE_RANKS, or PETSCSF_DUPLICATE_GRAPH (see PetscSFDuplicateOption) 40395fce210SBarry Smith 40495fce210SBarry Smith Output Arguments: 40595fce210SBarry Smith . newsf - new communication object 40695fce210SBarry Smith 40795fce210SBarry Smith Level: beginner 40895fce210SBarry Smith 40995fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetType(), PetscSFSetGraph() 41095fce210SBarry Smith @*/ 41195fce210SBarry Smith PetscErrorCode PetscSFDuplicate(PetscSF sf,PetscSFDuplicateOption opt,PetscSF *newsf) 41295fce210SBarry Smith { 41395fce210SBarry Smith PetscErrorCode ierr; 41495fce210SBarry Smith 41595fce210SBarry Smith PetscFunctionBegin; 41695fce210SBarry Smith ierr = PetscSFCreate(PetscObjectComm((PetscObject)sf),newsf);CHKERRQ(ierr); 41795fce210SBarry Smith ierr = PetscSFSetType(*newsf,((PetscObject)sf)->type_name);CHKERRQ(ierr); 41895fce210SBarry Smith if (sf->ops->Duplicate) {ierr = (*sf->ops->Duplicate)(sf,opt,*newsf);CHKERRQ(ierr);} 41995fce210SBarry Smith if (opt == PETSCSF_DUPLICATE_GRAPH) { 42095fce210SBarry Smith PetscInt nroots,nleaves; 42195fce210SBarry Smith const PetscInt *ilocal; 42295fce210SBarry Smith const PetscSFNode *iremote; 42395fce210SBarry Smith ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr); 42495fce210SBarry Smith ierr = PetscSFSetGraph(*newsf,nroots,nleaves,ilocal,PETSC_COPY_VALUES,iremote,PETSC_COPY_VALUES);CHKERRQ(ierr); 42595fce210SBarry Smith } 42695fce210SBarry Smith PetscFunctionReturn(0); 42795fce210SBarry Smith } 42895fce210SBarry Smith 42995fce210SBarry Smith /*@C 43095fce210SBarry Smith PetscSFGetGraph - Get the graph specifying a parallel star forest 43195fce210SBarry Smith 43295fce210SBarry Smith Not Collective 43395fce210SBarry Smith 43495fce210SBarry Smith Input Arguments: 43595fce210SBarry Smith . sf - star forest 43695fce210SBarry Smith 43795fce210SBarry Smith Output Arguments: 43895fce210SBarry Smith + nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves) 43995fce210SBarry Smith . nleaves - number of leaf vertices on the current process, each of these references a root on any process 44095fce210SBarry Smith . ilocal - locations of leaves in leafdata buffers 44195fce210SBarry Smith - iremote - remote locations of root vertices for each leaf on the current process 44295fce210SBarry Smith 44395fce210SBarry Smith Level: intermediate 44495fce210SBarry Smith 44595fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph() 44695fce210SBarry Smith @*/ 44795fce210SBarry Smith PetscErrorCode PetscSFGetGraph(PetscSF sf,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote) 44895fce210SBarry Smith { 44995fce210SBarry Smith 45095fce210SBarry Smith PetscFunctionBegin; 45195fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 45295fce210SBarry Smith /* We are not currently requiring that the graph is set, thus returning nroots=-1 if it has not been set */ 45395fce210SBarry Smith /* if (!sf->graphset) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Graph has not been set, must call PetscSFSetGraph()"); */ 45495fce210SBarry Smith if (nroots) *nroots = sf->nroots; 45595fce210SBarry Smith if (nleaves) *nleaves = sf->nleaves; 45695fce210SBarry Smith if (ilocal) *ilocal = sf->mine; 45795fce210SBarry Smith if (iremote) *iremote = sf->remote; 45895fce210SBarry Smith PetscFunctionReturn(0); 45995fce210SBarry Smith } 46095fce210SBarry Smith 46195fce210SBarry Smith /*@C 46295fce210SBarry Smith PetscSFGetLeafRange - Get the active leaf ranges 46395fce210SBarry Smith 46495fce210SBarry Smith Not Collective 46595fce210SBarry Smith 46695fce210SBarry Smith Input Arguments: 46795fce210SBarry Smith . sf - star forest 46895fce210SBarry Smith 46995fce210SBarry Smith Output Arguments: 47095fce210SBarry Smith + minleaf - minimum active leaf on this process 47195fce210SBarry Smith - maxleaf - maximum active leaf on this process 47295fce210SBarry Smith 47395fce210SBarry Smith Level: developer 47495fce210SBarry Smith 47595fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph(), PetscSFGetGraph() 47695fce210SBarry Smith @*/ 47795fce210SBarry Smith PetscErrorCode PetscSFGetLeafRange(PetscSF sf,PetscInt *minleaf,PetscInt *maxleaf) 47895fce210SBarry Smith { 47995fce210SBarry Smith 48095fce210SBarry Smith PetscFunctionBegin; 48195fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 48295fce210SBarry Smith if (minleaf) *minleaf = sf->minleaf; 48395fce210SBarry Smith if (maxleaf) *maxleaf = sf->maxleaf; 48495fce210SBarry Smith PetscFunctionReturn(0); 48595fce210SBarry Smith } 48695fce210SBarry Smith 48795fce210SBarry Smith /*@C 48895fce210SBarry Smith PetscSFView - view a star forest 48995fce210SBarry Smith 49095fce210SBarry Smith Collective 49195fce210SBarry Smith 49295fce210SBarry Smith Input Arguments: 49395fce210SBarry Smith + sf - star forest 49495fce210SBarry Smith - viewer - viewer to display graph, for example PETSC_VIEWER_STDOUT_WORLD 49595fce210SBarry Smith 49695fce210SBarry Smith Level: beginner 49795fce210SBarry Smith 49895fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph() 49995fce210SBarry Smith @*/ 50095fce210SBarry Smith PetscErrorCode PetscSFView(PetscSF sf,PetscViewer viewer) 50195fce210SBarry Smith { 50295fce210SBarry Smith PetscErrorCode ierr; 50395fce210SBarry Smith PetscBool iascii; 50495fce210SBarry Smith PetscViewerFormat format; 50595fce210SBarry Smith 50695fce210SBarry Smith PetscFunctionBegin; 50795fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 50895fce210SBarry Smith if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sf),&viewer);CHKERRQ(ierr);} 50995fce210SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 51095fce210SBarry Smith PetscCheckSameComm(sf,1,viewer,2); 511dcca2fcaSJed Brown ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 51295fce210SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 51395fce210SBarry Smith if (iascii) { 51495fce210SBarry Smith PetscMPIInt rank; 51581bfa7aaSJed Brown PetscInt ii,i,j; 51695fce210SBarry Smith 517dae58748SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)sf,viewer);CHKERRQ(ierr); 51895fce210SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 51995fce210SBarry Smith if (sf->ops->View) {ierr = (*sf->ops->View)(sf,viewer);CHKERRQ(ierr);} 52095fce210SBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); 5211575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 52295fce210SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of roots=%D, leaves=%D, remote ranks=%D\n",rank,sf->nroots,sf->nleaves,sf->nranks);CHKERRQ(ierr); 52395fce210SBarry Smith for (i=0; i<sf->nleaves; i++) { 52495fce210SBarry 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); 52595fce210SBarry Smith } 52695fce210SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 52795fce210SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 52895fce210SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 52981bfa7aaSJed Brown PetscMPIInt *tmpranks,*perm; 53081bfa7aaSJed Brown ierr = PetscMalloc2(sf->nranks,&tmpranks,sf->nranks,&perm);CHKERRQ(ierr); 53181bfa7aaSJed Brown ierr = PetscMemcpy(tmpranks,sf->ranks,sf->nranks*sizeof(tmpranks[0]));CHKERRQ(ierr); 53281bfa7aaSJed Brown for (i=0; i<sf->nranks; i++) perm[i] = i; 53381bfa7aaSJed Brown ierr = PetscSortMPIIntWithArray(sf->nranks,tmpranks,perm);CHKERRQ(ierr); 53495fce210SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Roots referenced by my leaves, by rank\n",rank);CHKERRQ(ierr); 53581bfa7aaSJed Brown for (ii=0; ii<sf->nranks; ii++) { 53681bfa7aaSJed Brown i = perm[ii]; 5377904a332SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %d: %D edges\n",rank,sf->ranks[i],sf->roffset[i+1]-sf->roffset[i]);CHKERRQ(ierr); 53895fce210SBarry Smith for (j=sf->roffset[i]; j<sf->roffset[i+1]; j++) { 53995fce210SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D <- %D\n",rank,sf->rmine[j],sf->rremote[j]);CHKERRQ(ierr); 54095fce210SBarry Smith } 54195fce210SBarry Smith } 54281bfa7aaSJed Brown ierr = PetscFree2(tmpranks,perm);CHKERRQ(ierr); 54395fce210SBarry Smith } 54495fce210SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 5451575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 54695fce210SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 54795fce210SBarry Smith } 54895fce210SBarry Smith PetscFunctionReturn(0); 54995fce210SBarry Smith } 55095fce210SBarry Smith 55195fce210SBarry Smith /*@C 55295fce210SBarry Smith PetscSFGetRanks - Get ranks and number of vertices referenced by leaves on this process 55395fce210SBarry Smith 55495fce210SBarry Smith Not Collective 55595fce210SBarry Smith 55695fce210SBarry Smith Input Arguments: 55795fce210SBarry Smith . sf - star forest 55895fce210SBarry Smith 55995fce210SBarry Smith Output Arguments: 56095fce210SBarry Smith + nranks - number of ranks referenced by local part 56195fce210SBarry Smith . ranks - array of ranks 56295fce210SBarry Smith . roffset - offset in rmine/rremote for each rank (length nranks+1) 56395fce210SBarry Smith . rmine - concatenated array holding local indices referencing each remote rank 56495fce210SBarry Smith - rremote - concatenated array holding remote indices referenced for each remote rank 56595fce210SBarry Smith 56695fce210SBarry Smith Level: developer 56795fce210SBarry Smith 56895fce210SBarry Smith .seealso: PetscSFSetGraph() 56995fce210SBarry Smith @*/ 57095fce210SBarry Smith PetscErrorCode PetscSFGetRanks(PetscSF sf,PetscInt *nranks,const PetscMPIInt **ranks,const PetscInt **roffset,const PetscInt **rmine,const PetscInt **rremote) 57195fce210SBarry Smith { 57295fce210SBarry Smith 57395fce210SBarry Smith PetscFunctionBegin; 57495fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 575dcca2fcaSJed Brown if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp before obtaining ranks"); 57695fce210SBarry Smith if (nranks) *nranks = sf->nranks; 57795fce210SBarry Smith if (ranks) *ranks = sf->ranks; 57895fce210SBarry Smith if (roffset) *roffset = sf->roffset; 57995fce210SBarry Smith if (rmine) *rmine = sf->rmine; 58095fce210SBarry Smith if (rremote) *rremote = sf->rremote; 58195fce210SBarry Smith PetscFunctionReturn(0); 58295fce210SBarry Smith } 58395fce210SBarry Smith 584b5a8e515SJed Brown static PetscBool InList(PetscMPIInt needle,PetscMPIInt n,const PetscMPIInt *list) { 585b5a8e515SJed Brown PetscInt i; 586b5a8e515SJed Brown for (i=0; i<n; i++) { 587b5a8e515SJed Brown if (needle == list[i]) return PETSC_TRUE; 588b5a8e515SJed Brown } 589b5a8e515SJed Brown return PETSC_FALSE; 590b5a8e515SJed Brown } 591b5a8e515SJed Brown 59295fce210SBarry Smith /*@C 59321c688dcSJed Brown PetscSFSetUpRanks - Set up data structures associated with ranks; this is for internal use by PetscSF implementations. 59421c688dcSJed Brown 59521c688dcSJed Brown Collective 59621c688dcSJed Brown 59721c688dcSJed Brown Input Arguments: 598b5a8e515SJed Brown + sf - PetscSF to set up; PetscSFSetGraph() must have been called 599b5a8e515SJed Brown - dgroup - MPI_Group of ranks to be distinguished (e.g., for self or shared memory exchange) 60021c688dcSJed Brown 60121c688dcSJed Brown Level: developer 60221c688dcSJed Brown 60321c688dcSJed Brown .seealso: PetscSFGetRanks() 60421c688dcSJed Brown @*/ 605b5a8e515SJed Brown PetscErrorCode PetscSFSetUpRanks(PetscSF sf,MPI_Group dgroup) 60621c688dcSJed Brown { 60721c688dcSJed Brown PetscErrorCode ierr; 60821c688dcSJed Brown PetscTable table; 60921c688dcSJed Brown PetscTablePosition pos; 610b5a8e515SJed Brown PetscMPIInt size,groupsize,*groupranks; 61121c688dcSJed Brown PetscInt i,*rcount,*ranks; 61221c688dcSJed Brown 61321c688dcSJed Brown PetscFunctionBegin; 61421c688dcSJed Brown PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 61521c688dcSJed Brown if (!sf->graphset) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"PetscSFSetGraph() has not been called yet"); 61621c688dcSJed Brown ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr); 61721c688dcSJed Brown ierr = PetscTableCreate(10,size,&table);CHKERRQ(ierr); 61821c688dcSJed Brown for (i=0; i<sf->nleaves; i++) { 61921c688dcSJed Brown /* Log 1-based rank */ 62021c688dcSJed Brown ierr = PetscTableAdd(table,sf->remote[i].rank+1,1,ADD_VALUES);CHKERRQ(ierr); 62121c688dcSJed Brown } 62221c688dcSJed Brown ierr = PetscTableGetCount(table,&sf->nranks);CHKERRQ(ierr); 62321c688dcSJed Brown ierr = PetscMalloc4(sf->nranks,&sf->ranks,sf->nranks+1,&sf->roffset,sf->nleaves,&sf->rmine,sf->nleaves,&sf->rremote);CHKERRQ(ierr); 62421c688dcSJed Brown ierr = PetscMalloc2(sf->nranks,&rcount,sf->nranks,&ranks);CHKERRQ(ierr); 62521c688dcSJed Brown ierr = PetscTableGetHeadPosition(table,&pos);CHKERRQ(ierr); 62621c688dcSJed Brown for (i=0; i<sf->nranks; i++) { 62721c688dcSJed Brown ierr = PetscTableGetNext(table,&pos,&ranks[i],&rcount[i]);CHKERRQ(ierr); 62821c688dcSJed Brown ranks[i]--; /* Convert back to 0-based */ 62921c688dcSJed Brown } 63021c688dcSJed Brown ierr = PetscTableDestroy(&table);CHKERRQ(ierr); 631b5a8e515SJed Brown 632b5a8e515SJed Brown /* We expect that dgroup is reliably "small" while nranks could be large */ 633b5a8e515SJed Brown { 634*7fb8a5e4SKarl Rupp MPI_Group group = MPI_GROUP_NULL; 635b5a8e515SJed Brown PetscMPIInt *dgroupranks; 636b5a8e515SJed Brown ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr); 637b5a8e515SJed Brown ierr = MPI_Group_size(dgroup,&groupsize);CHKERRQ(ierr); 638b5a8e515SJed Brown ierr = PetscMalloc1(groupsize,&dgroupranks);CHKERRQ(ierr); 639b5a8e515SJed Brown ierr = PetscMalloc1(groupsize,&groupranks);CHKERRQ(ierr); 640b5a8e515SJed Brown for (i=0; i<groupsize; i++) dgroupranks[i] = i; 641f3fc9a17SJed Brown if (groupsize) {ierr = MPI_Group_translate_ranks(dgroup,groupsize,dgroupranks,group,groupranks);CHKERRQ(ierr);} 642b5a8e515SJed Brown ierr = MPI_Group_free(&group);CHKERRQ(ierr); 643b5a8e515SJed Brown ierr = PetscFree(dgroupranks);CHKERRQ(ierr); 644b5a8e515SJed Brown } 645b5a8e515SJed Brown 646b5a8e515SJed Brown /* Partition ranks[] into distinguished (first sf->ndranks) followed by non-distinguished */ 647b5a8e515SJed Brown for (sf->ndranks=0,i=sf->nranks; sf->ndranks<i; ) { 648b5a8e515SJed Brown for (i--; sf->ndranks<i; i--) { /* Scan i backward looking for distinguished rank */ 649b5a8e515SJed Brown if (InList(ranks[i],groupsize,groupranks)) break; 650b5a8e515SJed Brown } 651b5a8e515SJed Brown for ( ; sf->ndranks<=i; sf->ndranks++) { /* Scan sf->ndranks forward looking for non-distinguished rank */ 652b5a8e515SJed Brown if (!InList(ranks[sf->ndranks],groupsize,groupranks)) break; 653b5a8e515SJed Brown } 654b5a8e515SJed Brown if (sf->ndranks < i) { /* Swap ranks[sf->ndranks] with ranks[i] */ 655b5a8e515SJed Brown PetscInt tmprank,tmpcount; 656b5a8e515SJed Brown tmprank = ranks[i]; 657b5a8e515SJed Brown tmpcount = rcount[i]; 658b5a8e515SJed Brown ranks[i] = ranks[sf->ndranks]; 659b5a8e515SJed Brown rcount[i] = rcount[sf->ndranks]; 660b5a8e515SJed Brown ranks[sf->ndranks] = tmprank; 661b5a8e515SJed Brown rcount[sf->ndranks] = tmpcount; 662b5a8e515SJed Brown sf->ndranks++; 663b5a8e515SJed Brown } 664b5a8e515SJed Brown } 665b5a8e515SJed Brown ierr = PetscFree(groupranks);CHKERRQ(ierr); 666b5a8e515SJed Brown ierr = PetscSortIntWithArray(sf->ndranks,ranks,rcount);CHKERRQ(ierr); 667b5a8e515SJed Brown ierr = PetscSortIntWithArray(sf->nranks-sf->ndranks,ranks+sf->ndranks,rcount+sf->ndranks);CHKERRQ(ierr); 66821c688dcSJed Brown sf->roffset[0] = 0; 66921c688dcSJed Brown for (i=0; i<sf->nranks; i++) { 67021c688dcSJed Brown ierr = PetscMPIIntCast(ranks[i],sf->ranks+i);CHKERRQ(ierr); 67121c688dcSJed Brown sf->roffset[i+1] = sf->roffset[i] + rcount[i]; 67221c688dcSJed Brown rcount[i] = 0; 67321c688dcSJed Brown } 67421c688dcSJed Brown for (i=0; i<sf->nleaves; i++) { 675b5a8e515SJed Brown PetscInt irank; 67621c688dcSJed Brown /* Search for index of iremote[i].rank in sf->ranks */ 677b5a8e515SJed Brown ierr = PetscFindMPIInt(sf->remote[i].rank,sf->ndranks,sf->ranks,&irank);CHKERRQ(ierr); 678b5a8e515SJed Brown if (irank < 0) { 679b5a8e515SJed Brown ierr = PetscFindMPIInt(sf->remote[i].rank,sf->nranks-sf->ndranks,sf->ranks+sf->ndranks,&irank);CHKERRQ(ierr); 680b5a8e515SJed Brown if (irank >= 0) irank += sf->ndranks; 68121c688dcSJed Brown } 682b5a8e515SJed Brown if (irank < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Could not find rank %D in array",sf->remote[i].rank); 68321c688dcSJed Brown sf->rmine[sf->roffset[irank] + rcount[irank]] = sf->mine ? sf->mine[i] : i; 68421c688dcSJed Brown sf->rremote[sf->roffset[irank] + rcount[irank]] = sf->remote[i].index; 68521c688dcSJed Brown rcount[irank]++; 68621c688dcSJed Brown } 68721c688dcSJed Brown ierr = PetscFree2(rcount,ranks);CHKERRQ(ierr); 68821c688dcSJed Brown PetscFunctionReturn(0); 68921c688dcSJed Brown } 69021c688dcSJed Brown 69121c688dcSJed Brown /*@C 69295fce210SBarry Smith PetscSFGetGroups - gets incoming and outgoing process groups 69395fce210SBarry Smith 69495fce210SBarry Smith Collective 69595fce210SBarry Smith 69695fce210SBarry Smith Input Argument: 69795fce210SBarry Smith . sf - star forest 69895fce210SBarry Smith 69995fce210SBarry Smith Output Arguments: 70095fce210SBarry Smith + incoming - group of origin processes for incoming edges (leaves that reference my roots) 70195fce210SBarry Smith - outgoing - group of destination processes for outgoing edges (roots that I reference) 70295fce210SBarry Smith 70395fce210SBarry Smith Level: developer 70495fce210SBarry Smith 70595fce210SBarry Smith .seealso: PetscSFGetWindow(), PetscSFRestoreWindow() 70695fce210SBarry Smith @*/ 70795fce210SBarry Smith PetscErrorCode PetscSFGetGroups(PetscSF sf,MPI_Group *incoming,MPI_Group *outgoing) 70895fce210SBarry Smith { 70995fce210SBarry Smith PetscErrorCode ierr; 710*7fb8a5e4SKarl Rupp MPI_Group group = MPI_GROUP_NULL; 71195fce210SBarry Smith 71295fce210SBarry Smith PetscFunctionBegin; 71395fce210SBarry Smith if (sf->ingroup == MPI_GROUP_NULL) { 71495fce210SBarry Smith PetscInt i; 71595fce210SBarry Smith const PetscInt *indegree; 71695fce210SBarry Smith PetscMPIInt rank,*outranks,*inranks; 71795fce210SBarry Smith PetscSFNode *remote; 71895fce210SBarry Smith PetscSF bgcount; 71995fce210SBarry Smith 72095fce210SBarry Smith /* Compute the number of incoming ranks */ 721785e854fSJed Brown ierr = PetscMalloc1(sf->nranks,&remote);CHKERRQ(ierr); 72295fce210SBarry Smith for (i=0; i<sf->nranks; i++) { 72395fce210SBarry Smith remote[i].rank = sf->ranks[i]; 72495fce210SBarry Smith remote[i].index = 0; 72595fce210SBarry Smith } 72695fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,&bgcount);CHKERRQ(ierr); 72795fce210SBarry Smith ierr = PetscSFSetGraph(bgcount,1,sf->nranks,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr); 72895fce210SBarry Smith ierr = PetscSFComputeDegreeBegin(bgcount,&indegree);CHKERRQ(ierr); 72995fce210SBarry Smith ierr = PetscSFComputeDegreeEnd(bgcount,&indegree);CHKERRQ(ierr); 73095fce210SBarry Smith 73195fce210SBarry Smith /* Enumerate the incoming ranks */ 732dcca6d9dSJed Brown ierr = PetscMalloc2(indegree[0],&inranks,sf->nranks,&outranks);CHKERRQ(ierr); 73395fce210SBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); 73495fce210SBarry Smith for (i=0; i<sf->nranks; i++) outranks[i] = rank; 73595fce210SBarry Smith ierr = PetscSFGatherBegin(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr); 73695fce210SBarry Smith ierr = PetscSFGatherEnd(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr); 73795fce210SBarry Smith ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr); 73895fce210SBarry Smith ierr = MPI_Group_incl(group,indegree[0],inranks,&sf->ingroup);CHKERRQ(ierr); 73995fce210SBarry Smith ierr = MPI_Group_free(&group);CHKERRQ(ierr); 74095fce210SBarry Smith ierr = PetscFree2(inranks,outranks);CHKERRQ(ierr); 74195fce210SBarry Smith ierr = PetscSFDestroy(&bgcount);CHKERRQ(ierr); 74295fce210SBarry Smith } 74395fce210SBarry Smith *incoming = sf->ingroup; 74495fce210SBarry Smith 74595fce210SBarry Smith if (sf->outgroup == MPI_GROUP_NULL) { 74695fce210SBarry Smith ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr); 74795fce210SBarry Smith ierr = MPI_Group_incl(group,sf->nranks,sf->ranks,&sf->outgroup);CHKERRQ(ierr); 74895fce210SBarry Smith ierr = MPI_Group_free(&group);CHKERRQ(ierr); 74995fce210SBarry Smith } 75095fce210SBarry Smith *outgoing = sf->outgroup; 75195fce210SBarry Smith PetscFunctionReturn(0); 75295fce210SBarry Smith } 75395fce210SBarry Smith 75495fce210SBarry Smith /*@C 75595fce210SBarry Smith PetscSFGetMultiSF - gets the inner SF implemeting gathers and scatters 75695fce210SBarry Smith 75795fce210SBarry Smith Collective 75895fce210SBarry Smith 75995fce210SBarry Smith Input Argument: 76095fce210SBarry Smith . sf - star forest that may contain roots with 0 or with more than 1 vertex 76195fce210SBarry Smith 76295fce210SBarry Smith Output Arguments: 76395fce210SBarry Smith . multi - star forest with split roots, such that each root has degree exactly 1 76495fce210SBarry Smith 76595fce210SBarry Smith Level: developer 76695fce210SBarry Smith 76795fce210SBarry Smith Notes: 76895fce210SBarry Smith 76995fce210SBarry Smith In most cases, users should use PetscSFGatherBegin() and PetscSFScatterBegin() instead of manipulating multi 77095fce210SBarry Smith directly. Since multi satisfies the stronger condition that each entry in the global space has exactly one incoming 77195fce210SBarry Smith edge, it is a candidate for future optimization that might involve its removal. 77295fce210SBarry Smith 77395fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFGatherBegin(), PetscSFScatterBegin() 77495fce210SBarry Smith @*/ 77595fce210SBarry Smith PetscErrorCode PetscSFGetMultiSF(PetscSF sf,PetscSF *multi) 77695fce210SBarry Smith { 77795fce210SBarry Smith PetscErrorCode ierr; 77895fce210SBarry Smith 77995fce210SBarry Smith PetscFunctionBegin; 78095fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 78195fce210SBarry Smith PetscValidPointer(multi,2); 78295fce210SBarry Smith if (sf->nroots < 0) { /* Graph has not been set yet; why do we need this? */ 78395fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr); 78495fce210SBarry Smith *multi = sf->multi; 78595fce210SBarry Smith PetscFunctionReturn(0); 78695fce210SBarry Smith } 78795fce210SBarry Smith if (!sf->multi) { 78895fce210SBarry Smith const PetscInt *indegree; 7899837ea96SMatthew G. Knepley PetscInt i,*inoffset,*outones,*outoffset,maxlocal; 79095fce210SBarry Smith PetscSFNode *remote; 79195fce210SBarry Smith ierr = PetscSFComputeDegreeBegin(sf,&indegree);CHKERRQ(ierr); 79295fce210SBarry Smith ierr = PetscSFComputeDegreeEnd(sf,&indegree);CHKERRQ(ierr); 7939837ea96SMatthew G. Knepley for (i=0,maxlocal=0; i<sf->nleaves; i++) maxlocal = PetscMax(maxlocal,(sf->mine ? sf->mine[i] : i)+1); 7949837ea96SMatthew G. Knepley ierr = PetscMalloc3(sf->nroots+1,&inoffset,maxlocal,&outones,maxlocal,&outoffset);CHKERRQ(ierr); 79595fce210SBarry Smith inoffset[0] = 0; 79695fce210SBarry Smith for (i=0; i<sf->nroots; i++) inoffset[i+1] = inoffset[i] + indegree[i]; 7979837ea96SMatthew G. Knepley for (i=0; i<maxlocal; i++) outones[i] = 1; 798dbd2ff41SBarry Smith ierr = PetscSFFetchAndOpBegin(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr); 799dbd2ff41SBarry Smith ierr = PetscSFFetchAndOpEnd(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr); 80095fce210SBarry Smith for (i=0; i<sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */ 80195fce210SBarry Smith #if 0 80295fce210SBarry Smith #if defined(PETSC_USE_DEBUG) /* Check that the expected number of increments occurred */ 80395fce210SBarry Smith for (i=0; i<sf->nroots; i++) { 80495fce210SBarry Smith if (inoffset[i] + indegree[i] != inoffset[i+1]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect result after PetscSFFetchAndOp"); 80595fce210SBarry Smith } 80695fce210SBarry Smith #endif 80795fce210SBarry Smith #endif 808785e854fSJed Brown ierr = PetscMalloc1(sf->nleaves,&remote);CHKERRQ(ierr); 80995fce210SBarry Smith for (i=0; i<sf->nleaves; i++) { 81095fce210SBarry Smith remote[i].rank = sf->remote[i].rank; 81138e7336fSToby Isaac remote[i].index = outoffset[sf->mine ? sf->mine[i] : i]; 81295fce210SBarry Smith } 81395fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr); 81401365b40SToby Isaac ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr); 81595fce210SBarry Smith if (sf->rankorder) { /* Sort the ranks */ 81695fce210SBarry Smith PetscMPIInt rank; 81795fce210SBarry Smith PetscInt *inranks,*newoffset,*outranks,*newoutoffset,*tmpoffset,maxdegree; 81895fce210SBarry Smith PetscSFNode *newremote; 81995fce210SBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); 82095fce210SBarry Smith for (i=0,maxdegree=0; i<sf->nroots; i++) maxdegree = PetscMax(maxdegree,indegree[i]); 8219837ea96SMatthew G. Knepley ierr = PetscMalloc5(sf->multi->nroots,&inranks,sf->multi->nroots,&newoffset,maxlocal,&outranks,maxlocal,&newoutoffset,maxdegree,&tmpoffset);CHKERRQ(ierr); 8229837ea96SMatthew G. Knepley for (i=0; i<maxlocal; i++) outranks[i] = rank; 8238bfbc91cSJed Brown ierr = PetscSFReduceBegin(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr); 8248bfbc91cSJed Brown ierr = PetscSFReduceEnd(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr); 82595fce210SBarry Smith /* Sort the incoming ranks at each vertex, build the inverse map */ 82695fce210SBarry Smith for (i=0; i<sf->nroots; i++) { 82795fce210SBarry Smith PetscInt j; 82895fce210SBarry Smith for (j=0; j<indegree[i]; j++) tmpoffset[j] = j; 82995fce210SBarry Smith ierr = PetscSortIntWithArray(indegree[i],inranks+inoffset[i],tmpoffset);CHKERRQ(ierr); 83095fce210SBarry Smith for (j=0; j<indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j; 83195fce210SBarry Smith } 83295fce210SBarry Smith ierr = PetscSFBcastBegin(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr); 83395fce210SBarry Smith ierr = PetscSFBcastEnd(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr); 834785e854fSJed Brown ierr = PetscMalloc1(sf->nleaves,&newremote);CHKERRQ(ierr); 83595fce210SBarry Smith for (i=0; i<sf->nleaves; i++) { 83695fce210SBarry Smith newremote[i].rank = sf->remote[i].rank; 83701365b40SToby Isaac newremote[i].index = newoutoffset[sf->mine ? sf->mine[i] : i]; 83895fce210SBarry Smith } 83901365b40SToby Isaac ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,newremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 84095fce210SBarry Smith ierr = PetscFree5(inranks,newoffset,outranks,newoutoffset,tmpoffset);CHKERRQ(ierr); 84195fce210SBarry Smith } 84295fce210SBarry Smith ierr = PetscFree3(inoffset,outones,outoffset);CHKERRQ(ierr); 84395fce210SBarry Smith } 84495fce210SBarry Smith *multi = sf->multi; 84595fce210SBarry Smith PetscFunctionReturn(0); 84695fce210SBarry Smith } 84795fce210SBarry Smith 84895fce210SBarry Smith /*@C 84995fce210SBarry Smith PetscSFCreateEmbeddedSF - removes edges from all but the selected roots, does not remap indices 85095fce210SBarry Smith 85195fce210SBarry Smith Collective 85295fce210SBarry Smith 85395fce210SBarry Smith Input Arguments: 85495fce210SBarry Smith + sf - original star forest 85595fce210SBarry Smith . nroots - number of roots to select on this process 85695fce210SBarry Smith - selected - selected roots on this process 85795fce210SBarry Smith 85895fce210SBarry Smith Output Arguments: 85995fce210SBarry Smith . newsf - new star forest 86095fce210SBarry Smith 86195fce210SBarry Smith Level: advanced 86295fce210SBarry Smith 86395fce210SBarry Smith Note: 86495fce210SBarry Smith To use the new PetscSF, it may be necessary to know the indices of the leaves that are still participating. This can 86595fce210SBarry Smith be done by calling PetscSFGetGraph(). 86695fce210SBarry Smith 86795fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFGetGraph() 86895fce210SBarry Smith @*/ 86995fce210SBarry Smith PetscErrorCode PetscSFCreateEmbeddedSF(PetscSF sf,PetscInt nroots,const PetscInt *selected,PetscSF *newsf) 87095fce210SBarry Smith { 8710511a646SMatthew G. Knepley PetscInt *rootdata, *leafdata, *ilocal; 87295fce210SBarry Smith PetscSFNode *iremote; 873fc1ede2bSMatthew G. Knepley PetscInt leafsize = 0, nleaves = 0, n, i; 8740511a646SMatthew G. Knepley PetscErrorCode ierr; 87595fce210SBarry Smith 87695fce210SBarry Smith PetscFunctionBegin; 87795fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 87895fce210SBarry Smith if (nroots) PetscValidPointer(selected,3); 87995fce210SBarry Smith PetscValidPointer(newsf,4); 8800511a646SMatthew G. Knepley if (sf->mine) for (i = 0; i < sf->nleaves; ++i) {leafsize = PetscMax(leafsize, sf->mine[i]+1);} 8810511a646SMatthew G. Knepley else leafsize = sf->nleaves; 8821795a4d1SJed Brown ierr = PetscCalloc2(sf->nroots,&rootdata,leafsize,&leafdata);CHKERRQ(ierr); 8830511a646SMatthew G. Knepley for (i=0; i<nroots; ++i) rootdata[selected[i]] = 1; 88495fce210SBarry Smith ierr = PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr); 88595fce210SBarry Smith ierr = PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr); 88695fce210SBarry Smith 8870511a646SMatthew G. Knepley for (i = 0; i < leafsize; ++i) nleaves += leafdata[i]; 888785e854fSJed Brown ierr = PetscMalloc1(nleaves,&ilocal);CHKERRQ(ierr); 889785e854fSJed Brown ierr = PetscMalloc1(nleaves,&iremote);CHKERRQ(ierr); 890fc1ede2bSMatthew G. Knepley for (i = 0, n = 0; i < sf->nleaves; ++i) { 8910511a646SMatthew G. Knepley const PetscInt lidx = sf->mine ? sf->mine[i] : i; 8920511a646SMatthew G. Knepley 8930511a646SMatthew G. Knepley if (leafdata[lidx]) { 894fc1ede2bSMatthew G. Knepley ilocal[n] = lidx; 895fc1ede2bSMatthew G. Knepley iremote[n].rank = sf->remote[i].rank; 896fc1ede2bSMatthew G. Knepley iremote[n].index = sf->remote[i].index; 897fc1ede2bSMatthew G. Knepley ++n; 89895fce210SBarry Smith } 89995fce210SBarry Smith } 900fc1ede2bSMatthew 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); 90195fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,newsf);CHKERRQ(ierr); 90295fce210SBarry Smith ierr = PetscSFSetGraph(*newsf,sf->nroots,nleaves,ilocal,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 90395fce210SBarry Smith ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr); 90495fce210SBarry Smith PetscFunctionReturn(0); 90595fce210SBarry Smith } 90695fce210SBarry Smith 9072f5fb4c2SMatthew G. Knepley /*@C 9082f5fb4c2SMatthew G. Knepley PetscSFCreateEmbeddedLeafSF - removes edges from all but the selected leaves, does not remap indices 9092f5fb4c2SMatthew G. Knepley 9102f5fb4c2SMatthew G. Knepley Collective 9112f5fb4c2SMatthew G. Knepley 9122f5fb4c2SMatthew G. Knepley Input Arguments: 9132f5fb4c2SMatthew G. Knepley + sf - original star forest 9142f5fb4c2SMatthew G. Knepley . nleaves - number of leaves to select on this process 9152f5fb4c2SMatthew G. Knepley - selected - selected leaves on this process 9162f5fb4c2SMatthew G. Knepley 9172f5fb4c2SMatthew G. Knepley Output Arguments: 9182f5fb4c2SMatthew G. Knepley . newsf - new star forest 9192f5fb4c2SMatthew G. Knepley 9202f5fb4c2SMatthew G. Knepley Level: advanced 9212f5fb4c2SMatthew G. Knepley 9222f5fb4c2SMatthew G. Knepley .seealso: PetscSFCreateEmbeddedSF(), PetscSFSetGraph(), PetscSFGetGraph() 9232f5fb4c2SMatthew G. Knepley @*/ 9242f5fb4c2SMatthew G. Knepley PetscErrorCode PetscSFCreateEmbeddedLeafSF(PetscSF sf, PetscInt nleaves, const PetscInt *selected, PetscSF *newsf) 9252f5fb4c2SMatthew G. Knepley { 9262f5fb4c2SMatthew G. Knepley PetscSFNode *iremote; 9272f5fb4c2SMatthew G. Knepley PetscInt *ilocal; 9282f5fb4c2SMatthew G. Knepley PetscInt i; 9292f5fb4c2SMatthew G. Knepley PetscErrorCode ierr; 9302f5fb4c2SMatthew G. Knepley 9312f5fb4c2SMatthew G. Knepley PetscFunctionBegin; 9322f5fb4c2SMatthew G. Knepley PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 9332f5fb4c2SMatthew G. Knepley if (nleaves) PetscValidPointer(selected, 3); 9342f5fb4c2SMatthew G. Knepley PetscValidPointer(newsf, 4); 9352f5fb4c2SMatthew G. Knepley ierr = PetscMalloc1(nleaves, &ilocal);CHKERRQ(ierr); 9362f5fb4c2SMatthew G. Knepley ierr = PetscMalloc1(nleaves, &iremote);CHKERRQ(ierr); 9372f5fb4c2SMatthew G. Knepley for (i = 0; i < nleaves; ++i) { 9382f5fb4c2SMatthew G. Knepley const PetscInt l = selected[i]; 9392f5fb4c2SMatthew G. Knepley 9402f5fb4c2SMatthew G. Knepley ilocal[i] = sf->mine ? sf->mine[l] : l; 9412f5fb4c2SMatthew G. Knepley iremote[i].rank = sf->remote[l].rank; 9422f5fb4c2SMatthew G. Knepley iremote[i].index = sf->remote[l].index; 9432f5fb4c2SMatthew G. Knepley } 9442f5fb4c2SMatthew G. Knepley ierr = PetscSFDuplicate(sf, PETSCSF_DUPLICATE_RANKS, newsf);CHKERRQ(ierr); 9452f5fb4c2SMatthew G. Knepley ierr = PetscSFSetGraph(*newsf, sf->nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr); 9462f5fb4c2SMatthew G. Knepley PetscFunctionReturn(0); 9472f5fb4c2SMatthew G. Knepley } 9482f5fb4c2SMatthew G. Knepley 94995fce210SBarry Smith /*@C 95095fce210SBarry Smith PetscSFBcastBegin - begin pointwise broadcast to be concluded with call to PetscSFBcastEnd() 95195fce210SBarry Smith 95295fce210SBarry Smith Collective on PetscSF 95395fce210SBarry Smith 95495fce210SBarry Smith Input Arguments: 95595fce210SBarry Smith + sf - star forest on which to communicate 95695fce210SBarry Smith . unit - data type associated with each node 95795fce210SBarry Smith - rootdata - buffer to broadcast 95895fce210SBarry Smith 95995fce210SBarry Smith Output Arguments: 96095fce210SBarry Smith . leafdata - buffer to update with values from each leaf's respective root 96195fce210SBarry Smith 96295fce210SBarry Smith Level: intermediate 96395fce210SBarry Smith 96495fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFView(), PetscSFBcastEnd(), PetscSFReduceBegin() 96595fce210SBarry Smith @*/ 96695fce210SBarry Smith PetscErrorCode PetscSFBcastBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata) 96795fce210SBarry Smith { 96895fce210SBarry Smith PetscErrorCode ierr; 96995fce210SBarry Smith 97095fce210SBarry Smith PetscFunctionBegin; 97195fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 97295fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 973563ffbabSMatthew G. Knepley ierr = PetscLogEventBegin(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr); 97495fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 97595fce210SBarry Smith ierr = (*sf->ops->BcastBegin)(sf,unit,rootdata,leafdata);CHKERRQ(ierr); 976563ffbabSMatthew G. Knepley ierr = PetscLogEventEnd(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr); 97795fce210SBarry Smith PetscFunctionReturn(0); 97895fce210SBarry Smith } 97995fce210SBarry Smith 98095fce210SBarry Smith /*@C 98195fce210SBarry Smith PetscSFBcastEnd - end a broadcast operation started with PetscSFBcastBegin() 98295fce210SBarry Smith 98395fce210SBarry Smith Collective 98495fce210SBarry Smith 98595fce210SBarry Smith Input Arguments: 98695fce210SBarry Smith + sf - star forest 98795fce210SBarry Smith . unit - data type 98895fce210SBarry Smith - rootdata - buffer to broadcast 98995fce210SBarry Smith 99095fce210SBarry Smith Output Arguments: 99195fce210SBarry Smith . leafdata - buffer to update with values from each leaf's respective root 99295fce210SBarry Smith 99395fce210SBarry Smith Level: intermediate 99495fce210SBarry Smith 99595fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFReduceEnd() 99695fce210SBarry Smith @*/ 99795fce210SBarry Smith PetscErrorCode PetscSFBcastEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata) 99895fce210SBarry Smith { 99995fce210SBarry Smith PetscErrorCode ierr; 100095fce210SBarry Smith 100195fce210SBarry Smith PetscFunctionBegin; 100295fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 100395fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 1004563ffbabSMatthew G. Knepley ierr = PetscLogEventBegin(PETSCSF_BcastEnd,sf,0,0,0);CHKERRQ(ierr); 100595fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 100695fce210SBarry Smith ierr = (*sf->ops->BcastEnd)(sf,unit,rootdata,leafdata);CHKERRQ(ierr); 1007563ffbabSMatthew G. Knepley ierr = PetscLogEventEnd(PETSCSF_BcastEnd,sf,0,0,0);CHKERRQ(ierr); 100895fce210SBarry Smith PetscFunctionReturn(0); 100995fce210SBarry Smith } 101095fce210SBarry Smith 101195fce210SBarry Smith /*@C 101295fce210SBarry Smith PetscSFReduceBegin - begin reduction of leafdata into rootdata, to be completed with call to PetscSFReduceEnd() 101395fce210SBarry Smith 101495fce210SBarry Smith Collective 101595fce210SBarry Smith 101695fce210SBarry Smith Input Arguments: 101795fce210SBarry Smith + sf - star forest 101895fce210SBarry Smith . unit - data type 101995fce210SBarry Smith . leafdata - values to reduce 102095fce210SBarry Smith - op - reduction operation 102195fce210SBarry Smith 102295fce210SBarry Smith Output Arguments: 102395fce210SBarry Smith . rootdata - result of reduction of values from all leaves of each root 102495fce210SBarry Smith 102595fce210SBarry Smith Level: intermediate 102695fce210SBarry Smith 102795fce210SBarry Smith .seealso: PetscSFBcastBegin() 102895fce210SBarry Smith @*/ 102995fce210SBarry Smith PetscErrorCode PetscSFReduceBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op) 103095fce210SBarry Smith { 103195fce210SBarry Smith PetscErrorCode ierr; 103295fce210SBarry Smith 103395fce210SBarry Smith PetscFunctionBegin; 103495fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 103595fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 1036563ffbabSMatthew G. Knepley ierr = PetscLogEventBegin(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr); 103795fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 103895fce210SBarry Smith ierr = (sf->ops->ReduceBegin)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr); 1039563ffbabSMatthew G. Knepley ierr = PetscLogEventEnd(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr); 104095fce210SBarry Smith PetscFunctionReturn(0); 104195fce210SBarry Smith } 104295fce210SBarry Smith 104395fce210SBarry Smith /*@C 104495fce210SBarry Smith PetscSFReduceEnd - end a reduction operation started with PetscSFReduceBegin() 104595fce210SBarry Smith 104695fce210SBarry Smith Collective 104795fce210SBarry Smith 104895fce210SBarry Smith Input Arguments: 104995fce210SBarry Smith + sf - star forest 105095fce210SBarry Smith . unit - data type 105195fce210SBarry Smith . leafdata - values to reduce 105295fce210SBarry Smith - op - reduction operation 105395fce210SBarry Smith 105495fce210SBarry Smith Output Arguments: 105595fce210SBarry Smith . rootdata - result of reduction of values from all leaves of each root 105695fce210SBarry Smith 105795fce210SBarry Smith Level: intermediate 105895fce210SBarry Smith 105995fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFBcastEnd() 106095fce210SBarry Smith @*/ 106195fce210SBarry Smith PetscErrorCode PetscSFReduceEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op) 106295fce210SBarry Smith { 106395fce210SBarry Smith PetscErrorCode ierr; 106495fce210SBarry Smith 106595fce210SBarry Smith PetscFunctionBegin; 106695fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 106795fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 1068563ffbabSMatthew G. Knepley ierr = PetscLogEventBegin(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr); 106995fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 107095fce210SBarry Smith ierr = (*sf->ops->ReduceEnd)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr); 1071563ffbabSMatthew G. Knepley ierr = PetscLogEventEnd(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr); 107295fce210SBarry Smith PetscFunctionReturn(0); 107395fce210SBarry Smith } 107495fce210SBarry Smith 107595fce210SBarry Smith /*@C 107695fce210SBarry Smith PetscSFComputeDegreeBegin - begin computation of degree for each root vertex, to be completed with PetscSFComputeDegreeEnd() 107795fce210SBarry Smith 107895fce210SBarry Smith Collective 107995fce210SBarry Smith 108095fce210SBarry Smith Input Arguments: 108195fce210SBarry Smith . sf - star forest 108295fce210SBarry Smith 108395fce210SBarry Smith Output Arguments: 108495fce210SBarry Smith . degree - degree of each root vertex 108595fce210SBarry Smith 108695fce210SBarry Smith Level: advanced 108795fce210SBarry Smith 108895fce210SBarry Smith .seealso: PetscSFGatherBegin() 108995fce210SBarry Smith @*/ 109095fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf,const PetscInt **degree) 109195fce210SBarry Smith { 109295fce210SBarry Smith PetscErrorCode ierr; 109395fce210SBarry Smith 109495fce210SBarry Smith PetscFunctionBegin; 109595fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 109695fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 109795fce210SBarry Smith PetscValidPointer(degree,2); 1098803bd9e8SMatthew G. Knepley if (!sf->degreeknown) { 10999837ea96SMatthew G. Knepley PetscInt i,maxlocal; 1100803bd9e8SMatthew G. Knepley if (sf->degree) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Calls to PetscSFComputeDegreeBegin() cannot be nested."); 11019837ea96SMatthew G. Knepley for (i=0,maxlocal=0; i<sf->nleaves; i++) maxlocal = PetscMax(maxlocal,(sf->mine ? sf->mine[i] : i)+1); 1102785e854fSJed Brown ierr = PetscMalloc1(sf->nroots,&sf->degree);CHKERRQ(ierr); 11039837ea96SMatthew G. Knepley ierr = PetscMalloc1(maxlocal,&sf->degreetmp);CHKERRQ(ierr); 110495fce210SBarry Smith for (i=0; i<sf->nroots; i++) sf->degree[i] = 0; 11059837ea96SMatthew G. Knepley for (i=0; i<maxlocal; i++) sf->degreetmp[i] = 1; 1106dbd2ff41SBarry Smith ierr = PetscSFReduceBegin(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);CHKERRQ(ierr); 110795fce210SBarry Smith } 110895fce210SBarry Smith *degree = NULL; 110995fce210SBarry Smith PetscFunctionReturn(0); 111095fce210SBarry Smith } 111195fce210SBarry Smith 111295fce210SBarry Smith /*@C 111395fce210SBarry Smith PetscSFComputeDegreeEnd - complete computation of degree for each root vertex, started with PetscSFComputeDegreeBegin() 111495fce210SBarry Smith 111595fce210SBarry Smith Collective 111695fce210SBarry Smith 111795fce210SBarry Smith Input Arguments: 111895fce210SBarry Smith . sf - star forest 111995fce210SBarry Smith 112095fce210SBarry Smith Output Arguments: 112195fce210SBarry Smith . degree - degree of each root vertex 112295fce210SBarry Smith 112395fce210SBarry Smith Level: developer 112495fce210SBarry Smith 112595fce210SBarry Smith .seealso: 112695fce210SBarry Smith @*/ 112795fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf,const PetscInt **degree) 112895fce210SBarry Smith { 112995fce210SBarry Smith PetscErrorCode ierr; 113095fce210SBarry Smith 113195fce210SBarry Smith PetscFunctionBegin; 113295fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 113395fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 113495fce210SBarry Smith if (!sf->degreeknown) { 1135dbd2ff41SBarry Smith ierr = PetscSFReduceEnd(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);CHKERRQ(ierr); 113695fce210SBarry Smith ierr = PetscFree(sf->degreetmp);CHKERRQ(ierr); 113795fce210SBarry Smith 113895fce210SBarry Smith sf->degreeknown = PETSC_TRUE; 113995fce210SBarry Smith } 114095fce210SBarry Smith *degree = sf->degree; 114195fce210SBarry Smith PetscFunctionReturn(0); 114295fce210SBarry Smith } 114395fce210SBarry Smith 114495fce210SBarry Smith /*@C 114595fce210SBarry Smith PetscSFFetchAndOpBegin - begin operation that fetches values from root and updates atomically by applying operation using my leaf value, to be completed with PetscSFFetchAndOpEnd() 114695fce210SBarry Smith 114795fce210SBarry Smith Collective 114895fce210SBarry Smith 114995fce210SBarry Smith Input Arguments: 115095fce210SBarry Smith + sf - star forest 115195fce210SBarry Smith . unit - data type 115295fce210SBarry Smith . leafdata - leaf values to use in reduction 115395fce210SBarry Smith - op - operation to use for reduction 115495fce210SBarry Smith 115595fce210SBarry Smith Output Arguments: 115695fce210SBarry Smith + rootdata - root values to be updated, input state is seen by first process to perform an update 115795fce210SBarry Smith - leafupdate - state at each leaf's respective root immediately prior to my atomic update 115895fce210SBarry Smith 115995fce210SBarry Smith Level: advanced 116095fce210SBarry Smith 116195fce210SBarry Smith Note: 116295fce210SBarry Smith The update is only atomic at the granularity provided by the hardware. Different roots referenced by the same process 116395fce210SBarry Smith might be updated in a different order. Furthermore, if a composite type is used for the unit datatype, atomicity is 116495fce210SBarry Smith not guaranteed across the whole vertex. Therefore, this function is mostly only used with primitive types such as 116595fce210SBarry Smith integers. 116695fce210SBarry Smith 116795fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFReduceBegin(), PetscSFSetGraph() 116895fce210SBarry Smith @*/ 116995fce210SBarry Smith PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op) 117095fce210SBarry Smith { 117195fce210SBarry Smith PetscErrorCode ierr; 117295fce210SBarry Smith 117395fce210SBarry Smith PetscFunctionBegin; 117495fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 117595fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 1176563ffbabSMatthew G. Knepley ierr = PetscLogEventBegin(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr); 117795fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 117895fce210SBarry Smith ierr = (*sf->ops->FetchAndOpBegin)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr); 1179563ffbabSMatthew G. Knepley ierr = PetscLogEventEnd(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr); 118095fce210SBarry Smith PetscFunctionReturn(0); 118195fce210SBarry Smith } 118295fce210SBarry Smith 118395fce210SBarry Smith /*@C 118495fce210SBarry 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 118595fce210SBarry Smith 118695fce210SBarry Smith Collective 118795fce210SBarry Smith 118895fce210SBarry Smith Input Arguments: 118995fce210SBarry Smith + sf - star forest 119095fce210SBarry Smith . unit - data type 119195fce210SBarry Smith . leafdata - leaf values to use in reduction 119295fce210SBarry Smith - op - operation to use for reduction 119395fce210SBarry Smith 119495fce210SBarry Smith Output Arguments: 119595fce210SBarry Smith + rootdata - root values to be updated, input state is seen by first process to perform an update 119695fce210SBarry Smith - leafupdate - state at each leaf's respective root immediately prior to my atomic update 119795fce210SBarry Smith 119895fce210SBarry Smith Level: advanced 119995fce210SBarry Smith 120095fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFReduceEnd(), PetscSFSetGraph() 120195fce210SBarry Smith @*/ 120295fce210SBarry Smith PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op) 120395fce210SBarry Smith { 120495fce210SBarry Smith PetscErrorCode ierr; 120595fce210SBarry Smith 120695fce210SBarry Smith PetscFunctionBegin; 120795fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 120895fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 1209563ffbabSMatthew G. Knepley ierr = PetscLogEventBegin(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr); 121095fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 121195fce210SBarry Smith ierr = (*sf->ops->FetchAndOpEnd)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr); 1212563ffbabSMatthew G. Knepley ierr = PetscLogEventEnd(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr); 121395fce210SBarry Smith PetscFunctionReturn(0); 121495fce210SBarry Smith } 121595fce210SBarry Smith 121695fce210SBarry Smith /*@C 121795fce210SBarry Smith PetscSFGatherBegin - begin pointwise gather of all leaves into multi-roots, to be completed with PetscSFGatherEnd() 121895fce210SBarry Smith 121995fce210SBarry Smith Collective 122095fce210SBarry Smith 122195fce210SBarry Smith Input Arguments: 122295fce210SBarry Smith + sf - star forest 122395fce210SBarry Smith . unit - data type 122495fce210SBarry Smith - leafdata - leaf data to gather to roots 122595fce210SBarry Smith 122695fce210SBarry Smith Output Argument: 122795fce210SBarry Smith . multirootdata - root buffer to gather into, amount of space per root is equal to its degree 122895fce210SBarry Smith 122995fce210SBarry Smith Level: intermediate 123095fce210SBarry Smith 123195fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin() 123295fce210SBarry Smith @*/ 123395fce210SBarry Smith PetscErrorCode PetscSFGatherBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata) 123495fce210SBarry Smith { 123595fce210SBarry Smith PetscErrorCode ierr; 123695fce210SBarry Smith PetscSF multi; 123795fce210SBarry Smith 123895fce210SBarry Smith PetscFunctionBegin; 123995fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 124095fce210SBarry Smith ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr); 12418bfbc91cSJed Brown ierr = PetscSFReduceBegin(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr); 124295fce210SBarry Smith PetscFunctionReturn(0); 124395fce210SBarry Smith } 124495fce210SBarry Smith 124595fce210SBarry Smith /*@C 124695fce210SBarry Smith PetscSFGatherEnd - ends pointwise gather operation that was started with PetscSFGatherBegin() 124795fce210SBarry Smith 124895fce210SBarry Smith Collective 124995fce210SBarry Smith 125095fce210SBarry Smith Input Arguments: 125195fce210SBarry Smith + sf - star forest 125295fce210SBarry Smith . unit - data type 125395fce210SBarry Smith - leafdata - leaf data to gather to roots 125495fce210SBarry Smith 125595fce210SBarry Smith Output Argument: 125695fce210SBarry Smith . multirootdata - root buffer to gather into, amount of space per root is equal to its degree 125795fce210SBarry Smith 125895fce210SBarry Smith Level: intermediate 125995fce210SBarry Smith 126095fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd() 126195fce210SBarry Smith @*/ 126295fce210SBarry Smith PetscErrorCode PetscSFGatherEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata) 126395fce210SBarry Smith { 126495fce210SBarry Smith PetscErrorCode ierr; 126595fce210SBarry Smith PetscSF multi; 126695fce210SBarry Smith 126795fce210SBarry Smith PetscFunctionBegin; 126895fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 126995fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 127095fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 127195fce210SBarry Smith ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr); 12728bfbc91cSJed Brown ierr = PetscSFReduceEnd(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr); 127395fce210SBarry Smith PetscFunctionReturn(0); 127495fce210SBarry Smith } 127595fce210SBarry Smith 127695fce210SBarry Smith /*@C 127795fce210SBarry Smith PetscSFScatterBegin - begin pointwise scatter operation from multi-roots to leaves, to be completed with PetscSFScatterEnd() 127895fce210SBarry Smith 127995fce210SBarry Smith Collective 128095fce210SBarry Smith 128195fce210SBarry Smith Input Arguments: 128295fce210SBarry Smith + sf - star forest 128395fce210SBarry Smith . unit - data type 128495fce210SBarry Smith - multirootdata - root buffer to send to each leaf, one unit of data per leaf 128595fce210SBarry Smith 128695fce210SBarry Smith Output Argument: 128795fce210SBarry Smith . leafdata - leaf data to be update with personal data from each respective root 128895fce210SBarry Smith 128995fce210SBarry Smith Level: intermediate 129095fce210SBarry Smith 129195fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin() 129295fce210SBarry Smith @*/ 129395fce210SBarry Smith PetscErrorCode PetscSFScatterBegin(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata) 129495fce210SBarry Smith { 129595fce210SBarry Smith PetscErrorCode ierr; 129695fce210SBarry Smith PetscSF multi; 129795fce210SBarry Smith 129895fce210SBarry Smith PetscFunctionBegin; 129995fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 130095fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 130195fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 130295fce210SBarry Smith ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr); 130395fce210SBarry Smith ierr = PetscSFBcastBegin(multi,unit,multirootdata,leafdata);CHKERRQ(ierr); 130495fce210SBarry Smith PetscFunctionReturn(0); 130595fce210SBarry Smith } 130695fce210SBarry Smith 130795fce210SBarry Smith /*@C 130895fce210SBarry Smith PetscSFScatterEnd - ends pointwise scatter operation that was started with PetscSFScatterBegin() 130995fce210SBarry Smith 131095fce210SBarry Smith Collective 131195fce210SBarry Smith 131295fce210SBarry Smith Input Arguments: 131395fce210SBarry Smith + sf - star forest 131495fce210SBarry Smith . unit - data type 131595fce210SBarry Smith - multirootdata - root buffer to send to each leaf, one unit of data per leaf 131695fce210SBarry Smith 131795fce210SBarry Smith Output Argument: 131895fce210SBarry Smith . leafdata - leaf data to be update with personal data from each respective root 131995fce210SBarry Smith 132095fce210SBarry Smith Level: intermediate 132195fce210SBarry Smith 132295fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd() 132395fce210SBarry Smith @*/ 132495fce210SBarry Smith PetscErrorCode PetscSFScatterEnd(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata) 132595fce210SBarry Smith { 132695fce210SBarry Smith PetscErrorCode ierr; 132795fce210SBarry Smith PetscSF multi; 132895fce210SBarry Smith 132995fce210SBarry Smith PetscFunctionBegin; 133095fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 133195fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 133295fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 133395fce210SBarry Smith ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr); 133495fce210SBarry Smith ierr = PetscSFBcastEnd(multi,unit,multirootdata,leafdata);CHKERRQ(ierr); 133595fce210SBarry Smith PetscFunctionReturn(0); 133695fce210SBarry Smith } 1337a7b3aa13SAta Mesgarnejad 1338a7b3aa13SAta Mesgarnejad /*@ 1339a7b3aa13SAta Mesgarnejad PetscSFCompose - Compose a new PetscSF equivalent to action to PetscSFs 1340a7b3aa13SAta Mesgarnejad 1341a7b3aa13SAta Mesgarnejad Input Parameters: 1342a7b3aa13SAta Mesgarnejad + sfA - The first PetscSF 1343a7b3aa13SAta Mesgarnejad - sfB - The second PetscSF 1344a7b3aa13SAta Mesgarnejad 1345a7b3aa13SAta Mesgarnejad Output Parameters: 1346a7b3aa13SAta Mesgarnejad . sfBA - equvalent PetscSF for applying A then B 1347a7b3aa13SAta Mesgarnejad 1348a7b3aa13SAta Mesgarnejad Level: developer 1349a7b3aa13SAta Mesgarnejad 1350a7b3aa13SAta Mesgarnejad .seealso: PetscSF, PetscSFGetGraph(), PetscSFSetGraph() 1351a7b3aa13SAta Mesgarnejad @*/ 1352a7b3aa13SAta Mesgarnejad PetscErrorCode PetscSFCompose(PetscSF sfA, PetscSF sfB, PetscSF *sfBA) 1353a7b3aa13SAta Mesgarnejad { 1354a7b3aa13SAta Mesgarnejad MPI_Comm comm; 1355a7b3aa13SAta Mesgarnejad const PetscSFNode *remotePointsA, *remotePointsB; 1356a7b3aa13SAta Mesgarnejad PetscSFNode *remotePointsBA; 1357a7b3aa13SAta Mesgarnejad const PetscInt *localPointsA, *localPointsB; 1358a7b3aa13SAta Mesgarnejad PetscInt numRootsA, numLeavesA, numRootsB, numLeavesB; 1359a7b3aa13SAta Mesgarnejad PetscErrorCode ierr; 1360a7b3aa13SAta Mesgarnejad 1361a7b3aa13SAta Mesgarnejad PetscFunctionBegin; 1362a7b3aa13SAta Mesgarnejad PetscValidHeaderSpecific(sfA, PETSCSF_CLASSID, 1); 1363a7b3aa13SAta Mesgarnejad PetscValidHeaderSpecific(sfB, PETSCSF_CLASSID, 1); 1364a7b3aa13SAta Mesgarnejad ierr = PetscObjectGetComm((PetscObject) sfA, &comm);CHKERRQ(ierr); 1365a7b3aa13SAta Mesgarnejad ierr = PetscSFGetGraph(sfA, &numRootsA, &numLeavesA, &localPointsA, &remotePointsA);CHKERRQ(ierr); 1366a7b3aa13SAta Mesgarnejad ierr = PetscSFGetGraph(sfB, &numRootsB, &numLeavesB, &localPointsB, &remotePointsB);CHKERRQ(ierr); 1367a7b3aa13SAta Mesgarnejad ierr = PetscMalloc1(numLeavesB, &remotePointsBA);CHKERRQ(ierr); 1368a7b3aa13SAta Mesgarnejad ierr = PetscSFBcastBegin(sfB, MPIU_2INT, remotePointsA, remotePointsBA);CHKERRQ(ierr); 1369a7b3aa13SAta Mesgarnejad ierr = PetscSFBcastEnd(sfB, MPIU_2INT, remotePointsA, remotePointsBA);CHKERRQ(ierr); 1370a7b3aa13SAta Mesgarnejad ierr = PetscSFCreate(comm, sfBA);CHKERRQ(ierr); 1371a7b3aa13SAta Mesgarnejad ierr = PetscSFSetGraph(*sfBA, numRootsA, numLeavesB, localPointsB, PETSC_COPY_VALUES, remotePointsBA, PETSC_OWN_POINTER);CHKERRQ(ierr); 1372a7b3aa13SAta Mesgarnejad PetscFunctionReturn(0); 1373a7b3aa13SAta Mesgarnejad } 1374