195fce210SBarry Smith #include <petsc-private/sfimpl.h> /*I "petscsf.h" I*/ 295fce210SBarry Smith #include <petscctable.h> 395fce210SBarry Smith 495fce210SBarry Smith #if defined(PETSC_USE_DEBUG) 595fce210SBarry Smith # define PetscSFCheckGraphSet(sf,arg) do { \ 695fce210SBarry Smith if (PetscUnlikely(!(sf)->graphset)) \ 795fce210SBarry Smith SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetGraph() on argument %D \"%s\" before %s()",(arg),#sf,PETSC_FUNCTION_NAME); \ 895fce210SBarry Smith } while (0) 995fce210SBarry Smith #else 1095fce210SBarry Smith # define PetscSFCheckGraphSet(sf,arg) do {} while (0) 1195fce210SBarry Smith #endif 1295fce210SBarry Smith 1395fce210SBarry Smith const char *const PetscSFDuplicateOptions[] = {"CONFONLY","RANKS","GRAPH","PetscSFDuplicateOption","PETSCSF_DUPLICATE_",0}; 1495fce210SBarry Smith 1595fce210SBarry Smith #undef __FUNCT__ 1695fce210SBarry Smith #define __FUNCT__ "PetscSFCreate" 1795fce210SBarry Smith /*@C 1895fce210SBarry Smith PetscSFCreate - create a star forest communication context 1995fce210SBarry Smith 2095fce210SBarry Smith Not Collective 2195fce210SBarry Smith 2295fce210SBarry Smith Input Arguments: 2395fce210SBarry Smith . comm - communicator on which the star forest will operate 2495fce210SBarry Smith 2595fce210SBarry Smith Output Arguments: 2695fce210SBarry Smith . sf - new star forest context 2795fce210SBarry Smith 2895fce210SBarry Smith Level: intermediate 2995fce210SBarry Smith 3095fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFDestroy() 3195fce210SBarry Smith @*/ 3295fce210SBarry Smith PetscErrorCode PetscSFCreate(MPI_Comm comm,PetscSF *sf) 3395fce210SBarry Smith { 3495fce210SBarry Smith PetscErrorCode ierr; 3595fce210SBarry Smith PetscSF b; 3695fce210SBarry Smith 3795fce210SBarry Smith PetscFunctionBegin; 3895fce210SBarry Smith PetscValidPointer(sf,2); 3995fce210SBarry Smith #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 40607a6623SBarry Smith ierr = PetscSFInitializePackage();CHKERRQ(ierr); 4195fce210SBarry Smith #endif 4295fce210SBarry Smith 4395fce210SBarry Smith ierr = PetscHeaderCreate(b,_p_PetscSF,struct _PetscSFOps,PETSCSF_CLASSID,"PetscSF","Star Forest","PetscSF",comm,PetscSFDestroy,PetscSFView);CHKERRQ(ierr); 4495fce210SBarry Smith 4595fce210SBarry Smith b->nroots = -1; 4695fce210SBarry Smith b->nleaves = -1; 4795fce210SBarry Smith b->nranks = -1; 4895fce210SBarry Smith b->rankorder = PETSC_TRUE; 4995fce210SBarry Smith b->ingroup = MPI_GROUP_NULL; 5095fce210SBarry Smith b->outgroup = MPI_GROUP_NULL; 5195fce210SBarry Smith b->graphset = PETSC_FALSE; 5295fce210SBarry Smith 5395fce210SBarry Smith *sf = b; 5495fce210SBarry Smith PetscFunctionReturn(0); 5595fce210SBarry Smith } 5695fce210SBarry Smith 5795fce210SBarry Smith #undef __FUNCT__ 5895fce210SBarry Smith #define __FUNCT__ "PetscSFReset" 5995fce210SBarry Smith /*@C 6095fce210SBarry Smith PetscSFReset - Reset a star forest so that different sizes or neighbors can be used 6195fce210SBarry Smith 6295fce210SBarry Smith Collective 6395fce210SBarry Smith 6495fce210SBarry Smith Input Arguments: 6595fce210SBarry Smith . sf - star forest 6695fce210SBarry Smith 6795fce210SBarry Smith Level: advanced 6895fce210SBarry Smith 6995fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFDestroy() 7095fce210SBarry Smith @*/ 7195fce210SBarry Smith PetscErrorCode PetscSFReset(PetscSF sf) 7295fce210SBarry Smith { 7395fce210SBarry Smith PetscErrorCode ierr; 7495fce210SBarry Smith 7595fce210SBarry Smith PetscFunctionBegin; 7695fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 7795fce210SBarry Smith sf->mine = NULL; 7895fce210SBarry Smith ierr = PetscFree(sf->mine_alloc);CHKERRQ(ierr); 7995fce210SBarry Smith sf->remote = NULL; 8095fce210SBarry Smith ierr = PetscFree(sf->remote_alloc);CHKERRQ(ierr); 8195fce210SBarry Smith ierr = PetscFree4(sf->ranks,sf->roffset,sf->rmine,sf->rremote);CHKERRQ(ierr); 8295fce210SBarry Smith ierr = PetscFree(sf->degree);CHKERRQ(ierr); 8395fce210SBarry Smith if (sf->ingroup != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->ingroup);CHKERRQ(ierr);} 8495fce210SBarry Smith if (sf->outgroup != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->outgroup);CHKERRQ(ierr);} 8595fce210SBarry Smith ierr = PetscSFDestroy(&sf->multi);CHKERRQ(ierr); 8695fce210SBarry Smith sf->graphset = PETSC_FALSE; 8795fce210SBarry Smith if (sf->ops->Reset) {ierr = (*sf->ops->Reset)(sf);CHKERRQ(ierr);} 8895fce210SBarry Smith sf->setupcalled = PETSC_FALSE; 8995fce210SBarry Smith PetscFunctionReturn(0); 9095fce210SBarry Smith } 9195fce210SBarry Smith 9295fce210SBarry Smith #undef __FUNCT__ 9395fce210SBarry Smith #define __FUNCT__ "PetscSFSetType" 9495fce210SBarry Smith /*@C 9595fce210SBarry Smith PetscSFSetType - set the PetscSF communication implementation 9695fce210SBarry Smith 9795fce210SBarry Smith Collective on PetscSF 9895fce210SBarry Smith 9995fce210SBarry Smith Input Parameters: 10095fce210SBarry Smith + sf - the PetscSF context 10195fce210SBarry Smith - type - a known method 10295fce210SBarry Smith 10395fce210SBarry Smith Options Database Key: 10495fce210SBarry Smith . -sf_type <type> - Sets the method; use -help for a list 10595fce210SBarry Smith of available methods (for instance, window, pt2pt, neighbor) 10695fce210SBarry Smith 10795fce210SBarry Smith Notes: 10895fce210SBarry Smith See "include/petscsf.h" for available methods (for instance) 10995fce210SBarry Smith + PETSCSFWINDOW - MPI-2/3 one-sided 11095fce210SBarry Smith - PETSCSFBASIC - basic implementation using MPI-1 two-sided 11195fce210SBarry Smith 11295fce210SBarry Smith Level: intermediate 11395fce210SBarry Smith 11495fce210SBarry Smith .keywords: PetscSF, set, type 11595fce210SBarry Smith 11695fce210SBarry Smith .seealso: PetscSFType, PetscSFCreate() 11795fce210SBarry Smith @*/ 11895fce210SBarry Smith PetscErrorCode PetscSFSetType(PetscSF sf,PetscSFType type) 11995fce210SBarry Smith { 12095fce210SBarry Smith PetscErrorCode ierr,(*r)(PetscSF); 12195fce210SBarry Smith PetscBool match; 12295fce210SBarry Smith 12395fce210SBarry Smith PetscFunctionBegin; 12495fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 12595fce210SBarry Smith PetscValidCharPointer(type,2); 12695fce210SBarry Smith 12795fce210SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)sf,type,&match);CHKERRQ(ierr); 12895fce210SBarry Smith if (match) PetscFunctionReturn(0); 12995fce210SBarry Smith 130adc40e5bSBarry Smith ierr = PetscFunctionListFind(PetscSFList,type,&r);CHKERRQ(ierr); 13195fce210SBarry Smith if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested PetscSF type %s",type); 13295fce210SBarry Smith /* Destroy the previous private PetscSF context */ 13395fce210SBarry Smith if (sf->ops->Destroy) { 13495fce210SBarry Smith ierr = (*(sf)->ops->Destroy)(sf);CHKERRQ(ierr); 13595fce210SBarry Smith } 13695fce210SBarry Smith ierr = PetscMemzero(sf->ops,sizeof(*sf->ops));CHKERRQ(ierr); 13795fce210SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)sf,type);CHKERRQ(ierr); 13895fce210SBarry Smith ierr = (*r)(sf);CHKERRQ(ierr); 13995fce210SBarry Smith PetscFunctionReturn(0); 14095fce210SBarry Smith } 14195fce210SBarry Smith 14295fce210SBarry Smith #undef __FUNCT__ 14395fce210SBarry Smith #define __FUNCT__ "PetscSFDestroy" 14495fce210SBarry Smith /*@C 14595fce210SBarry Smith PetscSFDestroy - destroy star forest 14695fce210SBarry Smith 14795fce210SBarry Smith Collective 14895fce210SBarry Smith 14995fce210SBarry Smith Input Arguments: 15095fce210SBarry Smith . sf - address of star forest 15195fce210SBarry Smith 15295fce210SBarry Smith Level: intermediate 15395fce210SBarry Smith 15495fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFReset() 15595fce210SBarry Smith @*/ 15695fce210SBarry Smith PetscErrorCode PetscSFDestroy(PetscSF *sf) 15795fce210SBarry Smith { 15895fce210SBarry Smith PetscErrorCode ierr; 15995fce210SBarry Smith 16095fce210SBarry Smith PetscFunctionBegin; 16195fce210SBarry Smith if (!*sf) PetscFunctionReturn(0); 16295fce210SBarry Smith PetscValidHeaderSpecific((*sf),PETSCSF_CLASSID,1); 16395fce210SBarry Smith if (--((PetscObject)(*sf))->refct > 0) {*sf = 0; PetscFunctionReturn(0);} 16495fce210SBarry Smith ierr = PetscSFReset(*sf);CHKERRQ(ierr); 16595fce210SBarry Smith if ((*sf)->ops->Destroy) {ierr = (*(*sf)->ops->Destroy)(*sf);CHKERRQ(ierr);} 16695fce210SBarry Smith ierr = PetscHeaderDestroy(sf);CHKERRQ(ierr); 16795fce210SBarry Smith PetscFunctionReturn(0); 16895fce210SBarry Smith } 16995fce210SBarry Smith 17095fce210SBarry Smith #undef __FUNCT__ 17195fce210SBarry Smith #define __FUNCT__ "PetscSFSetUp" 17295fce210SBarry Smith /*@ 17395fce210SBarry Smith PetscSFSetUp - set up communication structures 17495fce210SBarry Smith 17595fce210SBarry Smith Collective 17695fce210SBarry Smith 17795fce210SBarry Smith Input Arguments: 17895fce210SBarry Smith . sf - star forest communication object 17995fce210SBarry Smith 18095fce210SBarry Smith Level: beginner 18195fce210SBarry Smith 18295fce210SBarry Smith .seealso: PetscSFSetFromOptions(), PetscSFSetType() 18395fce210SBarry Smith @*/ 18495fce210SBarry Smith PetscErrorCode PetscSFSetUp(PetscSF sf) 18595fce210SBarry Smith { 18695fce210SBarry Smith PetscErrorCode ierr; 18795fce210SBarry Smith 18895fce210SBarry Smith PetscFunctionBegin; 18995fce210SBarry Smith if (sf->setupcalled) PetscFunctionReturn(0); 19095fce210SBarry Smith if (!((PetscObject)sf)->type_name) {ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr);} 19195fce210SBarry Smith if (sf->ops->SetUp) {ierr = (*sf->ops->SetUp)(sf);CHKERRQ(ierr);} 19295fce210SBarry Smith sf->setupcalled = PETSC_TRUE; 19395fce210SBarry Smith PetscFunctionReturn(0); 19495fce210SBarry Smith } 19595fce210SBarry Smith 19695fce210SBarry Smith #undef __FUNCT__ 19795fce210SBarry Smith #define __FUNCT__ "PetscSFSetFromOptions" 19895fce210SBarry Smith /*@C 19995fce210SBarry Smith PetscSFSetFromOptions - set PetscSF options using the options database 20095fce210SBarry Smith 20195fce210SBarry Smith Logically Collective 20295fce210SBarry Smith 20395fce210SBarry Smith Input Arguments: 20495fce210SBarry Smith . sf - star forest 20595fce210SBarry Smith 20695fce210SBarry Smith Options Database Keys: 20795fce210SBarry Smith . -sf_synchronization - synchronization type used by PetscSF 20895fce210SBarry Smith 20995fce210SBarry Smith Level: intermediate 21095fce210SBarry Smith 21195fce210SBarry Smith .keywords: KSP, set, from, options, database 21295fce210SBarry Smith 21395fce210SBarry Smith .seealso: PetscSFWindowSetSyncType() 21495fce210SBarry Smith @*/ 21595fce210SBarry Smith PetscErrorCode PetscSFSetFromOptions(PetscSF sf) 21695fce210SBarry Smith { 21795fce210SBarry Smith PetscSFType deft; 21895fce210SBarry Smith char type[256]; 21995fce210SBarry Smith PetscErrorCode ierr; 22095fce210SBarry Smith PetscBool flg; 22195fce210SBarry Smith 22295fce210SBarry Smith PetscFunctionBegin; 22395fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 22495fce210SBarry Smith ierr = PetscObjectOptionsBegin((PetscObject)sf);CHKERRQ(ierr); 22595fce210SBarry Smith deft = ((PetscObject)sf)->type_name ? ((PetscObject)sf)->type_name : PETSCSFBASIC; 226adc40e5bSBarry Smith ierr = PetscOptionsList("-sf_type","PetscSF implementation type","PetscSFSetType",PetscSFList,deft,type,256,&flg);CHKERRQ(ierr); 22795fce210SBarry Smith ierr = PetscSFSetType(sf,flg ? type : deft);CHKERRQ(ierr); 22895fce210SBarry 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); 22995fce210SBarry Smith if (sf->ops->SetFromOptions) {ierr = (*sf->ops->SetFromOptions)(sf);CHKERRQ(ierr);} 23095fce210SBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 23195fce210SBarry Smith PetscFunctionReturn(0); 23295fce210SBarry Smith } 23395fce210SBarry Smith 23495fce210SBarry Smith #undef __FUNCT__ 23595fce210SBarry Smith #define __FUNCT__ "PetscSFSetRankOrder" 23695fce210SBarry Smith /*@C 23795fce210SBarry Smith PetscSFSetRankOrder - sort multi-points for gathers and scatters by rank order 23895fce210SBarry Smith 23995fce210SBarry Smith Logically Collective 24095fce210SBarry Smith 24195fce210SBarry Smith Input Arguments: 24295fce210SBarry Smith + sf - star forest 24395fce210SBarry Smith - flg - PETSC_TRUE to sort, PETSC_FALSE to skip sorting (lower setup cost, but non-deterministic) 24495fce210SBarry Smith 24595fce210SBarry Smith Level: advanced 24695fce210SBarry Smith 24795fce210SBarry Smith .seealso: PetscSFGatherBegin(), PetscSFScatterBegin() 24895fce210SBarry Smith @*/ 24995fce210SBarry Smith PetscErrorCode PetscSFSetRankOrder(PetscSF sf,PetscBool flg) 25095fce210SBarry Smith { 25195fce210SBarry Smith 25295fce210SBarry Smith PetscFunctionBegin; 25395fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 25495fce210SBarry Smith PetscValidLogicalCollectiveBool(sf,flg,2); 25595fce210SBarry Smith if (sf->multi) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Rank ordering must be set before first call to PetscSFGatherBegin() or PetscSFScatterBegin()"); 25695fce210SBarry Smith sf->rankorder = flg; 25795fce210SBarry Smith PetscFunctionReturn(0); 25895fce210SBarry Smith } 25995fce210SBarry Smith 26095fce210SBarry Smith #undef __FUNCT__ 26195fce210SBarry Smith #define __FUNCT__ "PetscSFSetGraph" 26295fce210SBarry Smith /*@C 26395fce210SBarry Smith PetscSFSetGraph - Set a parallel star forest 26495fce210SBarry Smith 26595fce210SBarry Smith Collective 26695fce210SBarry Smith 26795fce210SBarry Smith Input Arguments: 26895fce210SBarry Smith + sf - star forest 26995fce210SBarry Smith . nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves) 27095fce210SBarry Smith . nleaves - number of leaf vertices on the current process, each of these references a root on any process 27195fce210SBarry Smith . ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage 27295fce210SBarry Smith . localmode - copy mode for ilocal 27395fce210SBarry Smith . iremote - remote locations of root vertices for each leaf on the current process 27495fce210SBarry Smith - remotemode - copy mode for iremote 27595fce210SBarry Smith 27695fce210SBarry Smith Level: intermediate 27795fce210SBarry Smith 27895fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFGetGraph() 27995fce210SBarry Smith @*/ 28095fce210SBarry Smith PetscErrorCode PetscSFSetGraph(PetscSF sf,PetscInt nroots,PetscInt nleaves,const PetscInt *ilocal,PetscCopyMode localmode,const PetscSFNode *iremote,PetscCopyMode remotemode) 28195fce210SBarry Smith { 28295fce210SBarry Smith PetscErrorCode ierr; 28395fce210SBarry Smith PetscTable table; 28495fce210SBarry Smith PetscTablePosition pos; 28595fce210SBarry Smith PetscMPIInt size; 28695fce210SBarry Smith PetscInt i,*rcount,*ranks; 28795fce210SBarry Smith 28895fce210SBarry Smith PetscFunctionBegin; 28995fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 29095fce210SBarry Smith if (nleaves && ilocal) PetscValidIntPointer(ilocal,4); 29195fce210SBarry Smith if (nleaves) PetscValidPointer(iremote,6); 29295fce210SBarry Smith if (nroots < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"roots %D, cannot be negative",nroots); 29395fce210SBarry Smith if (nleaves < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nleaves %D, cannot be negative",nleaves); 29495fce210SBarry Smith ierr = PetscSFReset(sf);CHKERRQ(ierr); 29595fce210SBarry Smith sf->nroots = nroots; 29695fce210SBarry Smith sf->nleaves = nleaves; 29795fce210SBarry Smith if (ilocal) { 29895fce210SBarry Smith switch (localmode) { 29995fce210SBarry Smith case PETSC_COPY_VALUES: 30095fce210SBarry Smith ierr = PetscMalloc(nleaves*sizeof(*sf->mine),&sf->mine_alloc);CHKERRQ(ierr); 30195fce210SBarry Smith sf->mine = sf->mine_alloc; 30295fce210SBarry Smith ierr = PetscMemcpy(sf->mine,ilocal,nleaves*sizeof(*sf->mine));CHKERRQ(ierr); 30395fce210SBarry Smith sf->minleaf = PETSC_MAX_INT; 30495fce210SBarry Smith sf->maxleaf = PETSC_MIN_INT; 30595fce210SBarry Smith for (i=0; i<nleaves; i++) { 30695fce210SBarry Smith sf->minleaf = PetscMin(sf->minleaf,ilocal[i]); 30795fce210SBarry Smith sf->maxleaf = PetscMax(sf->maxleaf,ilocal[i]); 30895fce210SBarry Smith } 30995fce210SBarry Smith break; 31095fce210SBarry Smith case PETSC_OWN_POINTER: 31195fce210SBarry Smith sf->mine_alloc = (PetscInt*)ilocal; 31295fce210SBarry Smith sf->mine = sf->mine_alloc; 31395fce210SBarry Smith break; 31495fce210SBarry Smith case PETSC_USE_POINTER: 31595fce210SBarry Smith sf->mine = (PetscInt*)ilocal; 31695fce210SBarry Smith break; 31795fce210SBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown localmode"); 31895fce210SBarry Smith } 31995fce210SBarry Smith } 32095fce210SBarry Smith if (!ilocal || nleaves > 0) { 32195fce210SBarry Smith sf->minleaf = 0; 32295fce210SBarry Smith sf->maxleaf = nleaves - 1; 32395fce210SBarry Smith } 32495fce210SBarry Smith switch (remotemode) { 32595fce210SBarry Smith case PETSC_COPY_VALUES: 32695fce210SBarry Smith ierr = PetscMalloc(nleaves*sizeof(*sf->remote),&sf->remote_alloc);CHKERRQ(ierr); 32795fce210SBarry Smith sf->remote = sf->remote_alloc; 32895fce210SBarry Smith ierr = PetscMemcpy(sf->remote,iremote,nleaves*sizeof(*sf->remote));CHKERRQ(ierr); 32995fce210SBarry Smith break; 33095fce210SBarry Smith case PETSC_OWN_POINTER: 33195fce210SBarry Smith sf->remote_alloc = (PetscSFNode*)iremote; 33295fce210SBarry Smith sf->remote = sf->remote_alloc; 33395fce210SBarry Smith break; 33495fce210SBarry Smith case PETSC_USE_POINTER: 33595fce210SBarry Smith sf->remote = (PetscSFNode*)iremote; 33695fce210SBarry Smith break; 33795fce210SBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown remotemode"); 33895fce210SBarry Smith } 33995fce210SBarry Smith 34095fce210SBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr); 34195fce210SBarry Smith ierr = PetscTableCreate(10,size,&table);CHKERRQ(ierr); 34295fce210SBarry Smith for (i=0; i<nleaves; i++) { 34395fce210SBarry Smith /* Log 1-based rank */ 34495fce210SBarry Smith ierr = PetscTableAdd(table,iremote[i].rank+1,1,ADD_VALUES);CHKERRQ(ierr); 34595fce210SBarry Smith } 34695fce210SBarry Smith ierr = PetscTableGetCount(table,&sf->nranks);CHKERRQ(ierr); 34795fce210SBarry Smith ierr = PetscMalloc4(sf->nranks,PetscInt,&sf->ranks,sf->nranks+1,PetscInt,&sf->roffset,nleaves,PetscInt,&sf->rmine,nleaves,PetscInt,&sf->rremote);CHKERRQ(ierr); 34895fce210SBarry Smith ierr = PetscMalloc2(sf->nranks,PetscInt,&rcount,sf->nranks,PetscInt,&ranks);CHKERRQ(ierr); 34995fce210SBarry Smith ierr = PetscTableGetHeadPosition(table,&pos);CHKERRQ(ierr); 35095fce210SBarry Smith for (i=0; i<sf->nranks; i++) { 35195fce210SBarry Smith ierr = PetscTableGetNext(table,&pos,&ranks[i],&rcount[i]);CHKERRQ(ierr); 35295fce210SBarry Smith ranks[i]--; /* Convert back to 0-based */ 35395fce210SBarry Smith } 35495fce210SBarry Smith ierr = PetscTableDestroy(&table);CHKERRQ(ierr); 35595fce210SBarry Smith ierr = PetscSortIntWithArray(sf->nranks,ranks,rcount);CHKERRQ(ierr); 35695fce210SBarry Smith sf->roffset[0] = 0; 35795fce210SBarry Smith for (i=0; i<sf->nranks; i++) { 35895fce210SBarry Smith ierr = PetscMPIIntCast(ranks[i],sf->ranks+i);CHKERRQ(ierr); 35995fce210SBarry Smith sf->roffset[i+1] = sf->roffset[i] + rcount[i]; 36095fce210SBarry Smith rcount[i] = 0; 36195fce210SBarry Smith } 36295fce210SBarry Smith for (i=0; i<nleaves; i++) { 36395fce210SBarry Smith PetscInt lo,hi,irank; 36495fce210SBarry Smith /* Search for index of iremote[i].rank in sf->ranks */ 36595fce210SBarry Smith lo = 0; hi = sf->nranks; 36695fce210SBarry Smith while (hi - lo > 1) { 36795fce210SBarry Smith PetscInt mid = lo + (hi - lo)/2; 36895fce210SBarry Smith if (iremote[i].rank < sf->ranks[mid]) hi = mid; 36995fce210SBarry Smith else lo = mid; 37095fce210SBarry Smith } 37195fce210SBarry Smith if (hi - lo == 1 && iremote[i].rank == sf->ranks[lo]) irank = lo; 37295fce210SBarry Smith else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Could not find rank %D in array",iremote[i].rank); 37395fce210SBarry Smith sf->rmine[sf->roffset[irank] + rcount[irank]] = ilocal ? ilocal[i] : i; 37495fce210SBarry Smith sf->rremote[sf->roffset[irank] + rcount[irank]] = iremote[i].index; 37595fce210SBarry Smith rcount[irank]++; 37695fce210SBarry Smith } 37795fce210SBarry Smith ierr = PetscFree2(rcount,ranks);CHKERRQ(ierr); 37895fce210SBarry Smith #if !defined(PETSC_USE_64BIT_INDICES) 37995fce210SBarry Smith if (nroots == PETSC_DETERMINE) { 38095fce210SBarry Smith /* Jed, if you have a better way to do this, put it in */ 38195fce210SBarry Smith PetscInt *numRankLeaves, *leafOff, *leafIndices, *numRankRoots, *rootOff, *rootIndices, maxRoots = 0; 38295fce210SBarry Smith 38395fce210SBarry Smith /* All to all to determine number of leaf indices from each (you can do this using Scan and asynch messages) */ 38495fce210SBarry Smith ierr = PetscMalloc4(size,PetscInt,&numRankLeaves,size+1,PetscInt,&leafOff,size,PetscInt,&numRankRoots,size+1,PetscInt,&rootOff);CHKERRQ(ierr); 38595fce210SBarry Smith ierr = PetscMemzero(numRankLeaves, size * sizeof(PetscInt));CHKERRQ(ierr); 38695fce210SBarry Smith for (i = 0; i < nleaves; ++i) ++numRankLeaves[iremote[i].rank]; 38795fce210SBarry Smith ierr = MPI_Alltoall(numRankLeaves, 1, MPIU_INT, numRankRoots, 1, MPIU_INT, PetscObjectComm((PetscObject)sf));CHKERRQ(ierr); 38895fce210SBarry Smith /* Could set nroots to this maximum */ 38995fce210SBarry Smith for (i = 0; i < size; ++i) maxRoots += numRankRoots[i]; 39095fce210SBarry Smith 39195fce210SBarry Smith /* Gather all indices */ 39295fce210SBarry Smith ierr = PetscMalloc2(nleaves,PetscInt,&leafIndices,maxRoots,PetscInt,&rootIndices);CHKERRQ(ierr); 39395fce210SBarry Smith leafOff[0] = 0; 39495fce210SBarry Smith for (i = 0; i < size; ++i) leafOff[i+1] = leafOff[i] + numRankLeaves[i]; 39595fce210SBarry Smith for (i = 0; i < nleaves; ++i) leafIndices[leafOff[iremote[i].rank]++] = iremote[i].index; 39695fce210SBarry Smith leafOff[0] = 0; 39795fce210SBarry Smith for (i = 0; i < size; ++i) leafOff[i+1] = leafOff[i] + numRankLeaves[i]; 39895fce210SBarry Smith rootOff[0] = 0; 39995fce210SBarry Smith for (i = 0; i < size; ++i) rootOff[i+1] = rootOff[i] + numRankRoots[i]; 40095fce210SBarry Smith ierr = MPI_Alltoallv(leafIndices, numRankLeaves, leafOff, MPIU_INT, rootIndices, numRankRoots, rootOff, MPIU_INT, PetscObjectComm((PetscObject)sf));CHKERRQ(ierr); 40195fce210SBarry Smith /* Sort and reduce */ 40295fce210SBarry Smith ierr = PetscSortRemoveDupsInt(&maxRoots, rootIndices);CHKERRQ(ierr); 40395fce210SBarry Smith ierr = PetscFree2(leafIndices,rootIndices);CHKERRQ(ierr); 40495fce210SBarry Smith ierr = PetscFree4(numRankLeaves,leafOff,numRankRoots,rootOff);CHKERRQ(ierr); 40595fce210SBarry Smith sf->nroots = maxRoots; 40695fce210SBarry Smith } 40795fce210SBarry Smith #endif 40895fce210SBarry Smith 40995fce210SBarry Smith sf->graphset = PETSC_TRUE; 41095fce210SBarry Smith PetscFunctionReturn(0); 41195fce210SBarry Smith } 41295fce210SBarry Smith 41395fce210SBarry Smith #undef __FUNCT__ 41495fce210SBarry Smith #define __FUNCT__ "PetscSFCreateInverseSF" 41595fce210SBarry Smith /*@C 41695fce210SBarry Smith PetscSFCreateInverseSF - given a PetscSF in which all vertices have degree 1, creates the inverse map 41795fce210SBarry Smith 41895fce210SBarry Smith Collective 41995fce210SBarry Smith 42095fce210SBarry Smith Input Arguments: 42195fce210SBarry Smith . sf - star forest to invert 42295fce210SBarry Smith 42395fce210SBarry Smith Output Arguments: 42495fce210SBarry Smith . isf - inverse of sf 42595fce210SBarry Smith 42695fce210SBarry Smith Level: advanced 42795fce210SBarry Smith 42895fce210SBarry Smith Notes: 42995fce210SBarry Smith All roots must have degree 1. 43095fce210SBarry Smith 43195fce210SBarry Smith The local space may be a permutation, but cannot be sparse. 43295fce210SBarry Smith 43395fce210SBarry Smith .seealso: PetscSFSetGraph() 43495fce210SBarry Smith @*/ 43595fce210SBarry Smith PetscErrorCode PetscSFCreateInverseSF(PetscSF sf,PetscSF *isf) 43695fce210SBarry Smith { 43795fce210SBarry Smith PetscErrorCode ierr; 43895fce210SBarry Smith PetscMPIInt rank; 43995fce210SBarry Smith PetscInt i,nroots,nleaves,maxlocal,count,*newilocal; 44095fce210SBarry Smith const PetscInt *ilocal; 44195fce210SBarry Smith PetscSFNode *roots,*leaves; 44295fce210SBarry Smith 44395fce210SBarry Smith PetscFunctionBegin; 44495fce210SBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); 44595fce210SBarry Smith ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);CHKERRQ(ierr); 44695fce210SBarry Smith for (i=0,maxlocal=0; i<nleaves; i++) maxlocal = PetscMax(maxlocal,(ilocal ? ilocal[i] : i)+1); 44795fce210SBarry Smith ierr = PetscMalloc2(nroots,PetscSFNode,&roots,nleaves,PetscSFNode,&leaves);CHKERRQ(ierr); 44895fce210SBarry Smith for (i=0; i<nleaves; i++) { 44995fce210SBarry Smith leaves[i].rank = rank; 45095fce210SBarry Smith leaves[i].index = i; 45195fce210SBarry Smith } 45295fce210SBarry Smith for (i=0; i <nroots; i++) { 45395fce210SBarry Smith roots[i].rank = -1; 45495fce210SBarry Smith roots[i].index = -1; 45595fce210SBarry Smith } 4568bfbc91cSJed Brown ierr = PetscSFReduceBegin(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr); 4578bfbc91cSJed Brown ierr = PetscSFReduceEnd(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr); 45895fce210SBarry Smith 45995fce210SBarry Smith /* Check whether our leaves are sparse */ 46095fce210SBarry Smith for (i=0,count=0; i<nroots; i++) if (roots[i].rank >= 0) count++; 46195fce210SBarry Smith if (count == nroots) newilocal = NULL; 46295fce210SBarry Smith else { /* Index for sparse leaves and compact "roots" array (which is to become our leaves). */ 46395fce210SBarry Smith ierr = PetscMalloc(count*sizeof(PetscInt),&newilocal);CHKERRQ(ierr); 46495fce210SBarry Smith for (i=0,count=0; i<nroots; i++) { 46595fce210SBarry Smith if (roots[i].rank >= 0) { 46695fce210SBarry Smith newilocal[count] = i; 46795fce210SBarry Smith roots[count].rank = roots[i].rank; 46895fce210SBarry Smith roots[count].index = roots[i].index; 46995fce210SBarry Smith count++; 47095fce210SBarry Smith } 47195fce210SBarry Smith } 47295fce210SBarry Smith } 47395fce210SBarry Smith 47495fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,isf);CHKERRQ(ierr); 47595fce210SBarry Smith ierr = PetscSFSetGraph(*isf,maxlocal,count,newilocal,PETSC_OWN_POINTER,roots,PETSC_COPY_VALUES);CHKERRQ(ierr); 47695fce210SBarry Smith ierr = PetscFree2(roots,leaves);CHKERRQ(ierr); 47795fce210SBarry Smith PetscFunctionReturn(0); 47895fce210SBarry Smith } 47995fce210SBarry Smith 48095fce210SBarry Smith #undef __FUNCT__ 48195fce210SBarry Smith #define __FUNCT__ "PetscSFDuplicate" 48295fce210SBarry Smith /*@ 48395fce210SBarry Smith PetscSFDuplicate - duplicate a PetscSF, optionally preserving rank connectivity and graph 48495fce210SBarry Smith 48595fce210SBarry Smith Collective 48695fce210SBarry Smith 48795fce210SBarry Smith Input Arguments: 48895fce210SBarry Smith + sf - communication object to duplicate 48995fce210SBarry Smith - opt - PETSCSF_DUPLICATE_CONFONLY, PETSCSF_DUPLICATE_RANKS, or PETSCSF_DUPLICATE_GRAPH (see PetscSFDuplicateOption) 49095fce210SBarry Smith 49195fce210SBarry Smith Output Arguments: 49295fce210SBarry Smith . newsf - new communication object 49395fce210SBarry Smith 49495fce210SBarry Smith Level: beginner 49595fce210SBarry Smith 49695fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetType(), PetscSFSetGraph() 49795fce210SBarry Smith @*/ 49895fce210SBarry Smith PetscErrorCode PetscSFDuplicate(PetscSF sf,PetscSFDuplicateOption opt,PetscSF *newsf) 49995fce210SBarry Smith { 50095fce210SBarry Smith PetscErrorCode ierr; 50195fce210SBarry Smith 50295fce210SBarry Smith PetscFunctionBegin; 50395fce210SBarry Smith ierr = PetscSFCreate(PetscObjectComm((PetscObject)sf),newsf);CHKERRQ(ierr); 50495fce210SBarry Smith ierr = PetscSFSetType(*newsf,((PetscObject)sf)->type_name);CHKERRQ(ierr); 50595fce210SBarry Smith if (sf->ops->Duplicate) {ierr = (*sf->ops->Duplicate)(sf,opt,*newsf);CHKERRQ(ierr);} 50695fce210SBarry Smith if (opt == PETSCSF_DUPLICATE_GRAPH) { 50795fce210SBarry Smith PetscInt nroots,nleaves; 50895fce210SBarry Smith const PetscInt *ilocal; 50995fce210SBarry Smith const PetscSFNode *iremote; 51095fce210SBarry Smith ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr); 51195fce210SBarry Smith ierr = PetscSFSetGraph(*newsf,nroots,nleaves,ilocal,PETSC_COPY_VALUES,iremote,PETSC_COPY_VALUES);CHKERRQ(ierr); 51295fce210SBarry Smith } 51395fce210SBarry Smith PetscFunctionReturn(0); 51495fce210SBarry Smith } 51595fce210SBarry Smith 51695fce210SBarry Smith #undef __FUNCT__ 51795fce210SBarry Smith #define __FUNCT__ "PetscSFGetGraph" 51895fce210SBarry Smith /*@C 51995fce210SBarry Smith PetscSFGetGraph - Get the graph specifying a parallel star forest 52095fce210SBarry Smith 52195fce210SBarry Smith Not Collective 52295fce210SBarry Smith 52395fce210SBarry Smith Input Arguments: 52495fce210SBarry Smith . sf - star forest 52595fce210SBarry Smith 52695fce210SBarry Smith Output Arguments: 52795fce210SBarry Smith + nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves) 52895fce210SBarry Smith . nleaves - number of leaf vertices on the current process, each of these references a root on any process 52995fce210SBarry Smith . ilocal - locations of leaves in leafdata buffers 53095fce210SBarry Smith - iremote - remote locations of root vertices for each leaf on the current process 53195fce210SBarry Smith 53295fce210SBarry Smith Level: intermediate 53395fce210SBarry Smith 53495fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph() 53595fce210SBarry Smith @*/ 53695fce210SBarry Smith PetscErrorCode PetscSFGetGraph(PetscSF sf,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote) 53795fce210SBarry Smith { 53895fce210SBarry Smith 53995fce210SBarry Smith PetscFunctionBegin; 54095fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 54195fce210SBarry Smith /* We are not currently requiring that the graph is set, thus returning nroots=-1 if it has not been set */ 54295fce210SBarry Smith /* if (!sf->graphset) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Graph has not been set, must call PetscSFSetGraph()"); */ 54395fce210SBarry Smith if (nroots) *nroots = sf->nroots; 54495fce210SBarry Smith if (nleaves) *nleaves = sf->nleaves; 54595fce210SBarry Smith if (ilocal) *ilocal = sf->mine; 54695fce210SBarry Smith if (iremote) *iremote = sf->remote; 54795fce210SBarry Smith PetscFunctionReturn(0); 54895fce210SBarry Smith } 54995fce210SBarry Smith 55095fce210SBarry Smith #undef __FUNCT__ 55195fce210SBarry Smith #define __FUNCT__ "PetscSFGetLeafRange" 55295fce210SBarry Smith /*@C 55395fce210SBarry Smith PetscSFGetLeafRange - Get the active leaf ranges 55495fce210SBarry Smith 55595fce210SBarry Smith Not Collective 55695fce210SBarry Smith 55795fce210SBarry Smith Input Arguments: 55895fce210SBarry Smith . sf - star forest 55995fce210SBarry Smith 56095fce210SBarry Smith Output Arguments: 56195fce210SBarry Smith + minleaf - minimum active leaf on this process 56295fce210SBarry Smith - maxleaf - maximum active leaf on this process 56395fce210SBarry Smith 56495fce210SBarry Smith Level: developer 56595fce210SBarry Smith 56695fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph(), PetscSFGetGraph() 56795fce210SBarry Smith @*/ 56895fce210SBarry Smith PetscErrorCode PetscSFGetLeafRange(PetscSF sf,PetscInt *minleaf,PetscInt *maxleaf) 56995fce210SBarry Smith { 57095fce210SBarry Smith 57195fce210SBarry Smith PetscFunctionBegin; 57295fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 57395fce210SBarry Smith if (minleaf) *minleaf = sf->minleaf; 57495fce210SBarry Smith if (maxleaf) *maxleaf = sf->maxleaf; 57595fce210SBarry Smith PetscFunctionReturn(0); 57695fce210SBarry Smith } 57795fce210SBarry Smith 57895fce210SBarry Smith #undef __FUNCT__ 57995fce210SBarry Smith #define __FUNCT__ "PetscSFView" 58095fce210SBarry Smith /*@C 58195fce210SBarry Smith PetscSFView - view a star forest 58295fce210SBarry Smith 58395fce210SBarry Smith Collective 58495fce210SBarry Smith 58595fce210SBarry Smith Input Arguments: 58695fce210SBarry Smith + sf - star forest 58795fce210SBarry Smith - viewer - viewer to display graph, for example PETSC_VIEWER_STDOUT_WORLD 58895fce210SBarry Smith 58995fce210SBarry Smith Level: beginner 59095fce210SBarry Smith 59195fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph() 59295fce210SBarry Smith @*/ 59395fce210SBarry Smith PetscErrorCode PetscSFView(PetscSF sf,PetscViewer viewer) 59495fce210SBarry Smith { 59595fce210SBarry Smith PetscErrorCode ierr; 59695fce210SBarry Smith PetscBool iascii; 59795fce210SBarry Smith PetscViewerFormat format; 59895fce210SBarry Smith 59995fce210SBarry Smith PetscFunctionBegin; 60095fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 60195fce210SBarry Smith if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sf),&viewer);CHKERRQ(ierr);} 60295fce210SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 60395fce210SBarry Smith PetscCheckSameComm(sf,1,viewer,2); 60495fce210SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 60595fce210SBarry Smith if (iascii) { 60695fce210SBarry Smith PetscMPIInt rank; 60795fce210SBarry Smith PetscInt i,j; 60895fce210SBarry Smith 609dae58748SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)sf,viewer);CHKERRQ(ierr); 61095fce210SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 61195fce210SBarry Smith if (sf->ops->View) {ierr = (*sf->ops->View)(sf,viewer);CHKERRQ(ierr);} 61295fce210SBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); 61395fce210SBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 61495fce210SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of roots=%D, leaves=%D, remote ranks=%D\n",rank,sf->nroots,sf->nleaves,sf->nranks);CHKERRQ(ierr); 61595fce210SBarry Smith for (i=0; i<sf->nleaves; i++) { 61695fce210SBarry 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); 61795fce210SBarry Smith } 61895fce210SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 61995fce210SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 62095fce210SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 62195fce210SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Roots referenced by my leaves, by rank\n",rank);CHKERRQ(ierr); 62295fce210SBarry Smith for (i=0; i<sf->nranks; i++) { 62395fce210SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D: %D edges\n",rank,sf->ranks[i],sf->roffset[i+1]-sf->roffset[i]);CHKERRQ(ierr); 62495fce210SBarry Smith for (j=sf->roffset[i]; j<sf->roffset[i+1]; j++) { 62595fce210SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D <- %D\n",rank,sf->rmine[j],sf->rremote[j]);CHKERRQ(ierr); 62695fce210SBarry Smith } 62795fce210SBarry Smith } 62895fce210SBarry Smith } 62995fce210SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 63095fce210SBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 63195fce210SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 63295fce210SBarry Smith } 63395fce210SBarry Smith PetscFunctionReturn(0); 63495fce210SBarry Smith } 63595fce210SBarry Smith 63695fce210SBarry Smith #undef __FUNCT__ 63795fce210SBarry Smith #define __FUNCT__ "PetscSFGetRanks" 63895fce210SBarry Smith /*@C 63995fce210SBarry Smith PetscSFGetRanks - Get ranks and number of vertices referenced by leaves on this process 64095fce210SBarry Smith 64195fce210SBarry Smith Not Collective 64295fce210SBarry Smith 64395fce210SBarry Smith Input Arguments: 64495fce210SBarry Smith . sf - star forest 64595fce210SBarry Smith 64695fce210SBarry Smith Output Arguments: 64795fce210SBarry Smith + nranks - number of ranks referenced by local part 64895fce210SBarry Smith . ranks - array of ranks 64995fce210SBarry Smith . roffset - offset in rmine/rremote for each rank (length nranks+1) 65095fce210SBarry Smith . rmine - concatenated array holding local indices referencing each remote rank 65195fce210SBarry Smith - rremote - concatenated array holding remote indices referenced for each remote rank 65295fce210SBarry Smith 65395fce210SBarry Smith Level: developer 65495fce210SBarry Smith 65595fce210SBarry Smith .seealso: PetscSFSetGraph() 65695fce210SBarry Smith @*/ 65795fce210SBarry Smith PetscErrorCode PetscSFGetRanks(PetscSF sf,PetscInt *nranks,const PetscMPIInt **ranks,const PetscInt **roffset,const PetscInt **rmine,const PetscInt **rremote) 65895fce210SBarry Smith { 65995fce210SBarry Smith 66095fce210SBarry Smith PetscFunctionBegin; 66195fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 66295fce210SBarry Smith if (nranks) *nranks = sf->nranks; 66395fce210SBarry Smith if (ranks) *ranks = sf->ranks; 66495fce210SBarry Smith if (roffset) *roffset = sf->roffset; 66595fce210SBarry Smith if (rmine) *rmine = sf->rmine; 66695fce210SBarry Smith if (rremote) *rremote = sf->rremote; 66795fce210SBarry Smith PetscFunctionReturn(0); 66895fce210SBarry Smith } 66995fce210SBarry Smith 67095fce210SBarry Smith #undef __FUNCT__ 67195fce210SBarry Smith #define __FUNCT__ "PetscSFGetGroups" 67295fce210SBarry Smith /*@C 67395fce210SBarry Smith PetscSFGetGroups - gets incoming and outgoing process groups 67495fce210SBarry Smith 67595fce210SBarry Smith Collective 67695fce210SBarry Smith 67795fce210SBarry Smith Input Argument: 67895fce210SBarry Smith . sf - star forest 67995fce210SBarry Smith 68095fce210SBarry Smith Output Arguments: 68195fce210SBarry Smith + incoming - group of origin processes for incoming edges (leaves that reference my roots) 68295fce210SBarry Smith - outgoing - group of destination processes for outgoing edges (roots that I reference) 68395fce210SBarry Smith 68495fce210SBarry Smith Level: developer 68595fce210SBarry Smith 68695fce210SBarry Smith .seealso: PetscSFGetWindow(), PetscSFRestoreWindow() 68795fce210SBarry Smith @*/ 68895fce210SBarry Smith PetscErrorCode PetscSFGetGroups(PetscSF sf,MPI_Group *incoming,MPI_Group *outgoing) 68995fce210SBarry Smith { 69095fce210SBarry Smith PetscErrorCode ierr; 69195fce210SBarry Smith MPI_Group group; 69295fce210SBarry Smith 69395fce210SBarry Smith PetscFunctionBegin; 69495fce210SBarry Smith if (sf->ingroup == MPI_GROUP_NULL) { 69595fce210SBarry Smith PetscInt i; 69695fce210SBarry Smith const PetscInt *indegree; 69795fce210SBarry Smith PetscMPIInt rank,*outranks,*inranks; 69895fce210SBarry Smith PetscSFNode *remote; 69995fce210SBarry Smith PetscSF bgcount; 70095fce210SBarry Smith 70195fce210SBarry Smith /* Compute the number of incoming ranks */ 70295fce210SBarry Smith ierr = PetscMalloc(sf->nranks*sizeof(PetscSFNode),&remote);CHKERRQ(ierr); 70395fce210SBarry Smith for (i=0; i<sf->nranks; i++) { 70495fce210SBarry Smith remote[i].rank = sf->ranks[i]; 70595fce210SBarry Smith remote[i].index = 0; 70695fce210SBarry Smith } 70795fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,&bgcount);CHKERRQ(ierr); 70895fce210SBarry Smith ierr = PetscSFSetGraph(bgcount,1,sf->nranks,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr); 70995fce210SBarry Smith ierr = PetscSFComputeDegreeBegin(bgcount,&indegree);CHKERRQ(ierr); 71095fce210SBarry Smith ierr = PetscSFComputeDegreeEnd(bgcount,&indegree);CHKERRQ(ierr); 71195fce210SBarry Smith 71295fce210SBarry Smith /* Enumerate the incoming ranks */ 71395fce210SBarry Smith ierr = PetscMalloc2(indegree[0],PetscMPIInt,&inranks,sf->nranks,PetscMPIInt,&outranks);CHKERRQ(ierr); 71495fce210SBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); 71595fce210SBarry Smith for (i=0; i<sf->nranks; i++) outranks[i] = rank; 71695fce210SBarry Smith ierr = PetscSFGatherBegin(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr); 71795fce210SBarry Smith ierr = PetscSFGatherEnd(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr); 71895fce210SBarry Smith ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr); 71995fce210SBarry Smith ierr = MPI_Group_incl(group,indegree[0],inranks,&sf->ingroup);CHKERRQ(ierr); 72095fce210SBarry Smith ierr = MPI_Group_free(&group);CHKERRQ(ierr); 72195fce210SBarry Smith ierr = PetscFree2(inranks,outranks);CHKERRQ(ierr); 72295fce210SBarry Smith ierr = PetscSFDestroy(&bgcount);CHKERRQ(ierr); 72395fce210SBarry Smith } 72495fce210SBarry Smith *incoming = sf->ingroup; 72595fce210SBarry Smith 72695fce210SBarry Smith if (sf->outgroup == MPI_GROUP_NULL) { 72795fce210SBarry Smith ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr); 72895fce210SBarry Smith ierr = MPI_Group_incl(group,sf->nranks,sf->ranks,&sf->outgroup);CHKERRQ(ierr); 72995fce210SBarry Smith ierr = MPI_Group_free(&group);CHKERRQ(ierr); 73095fce210SBarry Smith } 73195fce210SBarry Smith *outgoing = sf->outgroup; 73295fce210SBarry Smith PetscFunctionReturn(0); 73395fce210SBarry Smith } 73495fce210SBarry Smith 73595fce210SBarry Smith #undef __FUNCT__ 73695fce210SBarry Smith #define __FUNCT__ "PetscSFGetMultiSF" 73795fce210SBarry Smith /*@C 73895fce210SBarry Smith PetscSFGetMultiSF - gets the inner SF implemeting gathers and scatters 73995fce210SBarry Smith 74095fce210SBarry Smith Collective 74195fce210SBarry Smith 74295fce210SBarry Smith Input Argument: 74395fce210SBarry Smith . sf - star forest that may contain roots with 0 or with more than 1 vertex 74495fce210SBarry Smith 74595fce210SBarry Smith Output Arguments: 74695fce210SBarry Smith . multi - star forest with split roots, such that each root has degree exactly 1 74795fce210SBarry Smith 74895fce210SBarry Smith Level: developer 74995fce210SBarry Smith 75095fce210SBarry Smith Notes: 75195fce210SBarry Smith 75295fce210SBarry Smith In most cases, users should use PetscSFGatherBegin() and PetscSFScatterBegin() instead of manipulating multi 75395fce210SBarry Smith directly. Since multi satisfies the stronger condition that each entry in the global space has exactly one incoming 75495fce210SBarry Smith edge, it is a candidate for future optimization that might involve its removal. 75595fce210SBarry Smith 75695fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFGatherBegin(), PetscSFScatterBegin() 75795fce210SBarry Smith @*/ 75895fce210SBarry Smith PetscErrorCode PetscSFGetMultiSF(PetscSF sf,PetscSF *multi) 75995fce210SBarry Smith { 76095fce210SBarry Smith PetscErrorCode ierr; 76195fce210SBarry Smith 76295fce210SBarry Smith PetscFunctionBegin; 76395fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 76495fce210SBarry Smith PetscValidPointer(multi,2); 76595fce210SBarry Smith if (sf->nroots < 0) { /* Graph has not been set yet; why do we need this? */ 76695fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr); 76795fce210SBarry Smith *multi = sf->multi; 76895fce210SBarry Smith PetscFunctionReturn(0); 76995fce210SBarry Smith } 77095fce210SBarry Smith if (!sf->multi) { 77195fce210SBarry Smith const PetscInt *indegree; 77295fce210SBarry Smith PetscInt i,*inoffset,*outones,*outoffset; 77395fce210SBarry Smith PetscSFNode *remote; 77495fce210SBarry Smith ierr = PetscSFComputeDegreeBegin(sf,&indegree);CHKERRQ(ierr); 77595fce210SBarry Smith ierr = PetscSFComputeDegreeEnd(sf,&indegree);CHKERRQ(ierr); 77695fce210SBarry Smith ierr = PetscMalloc3(sf->nroots+1,PetscInt,&inoffset,sf->nleaves,PetscInt,&outones,sf->nleaves,PetscInt,&outoffset);CHKERRQ(ierr); 77795fce210SBarry Smith inoffset[0] = 0; 77895fce210SBarry Smith #if 1 77995fce210SBarry Smith for (i=0; i<sf->nroots; i++) inoffset[i+1] = PetscMax(i+1, inoffset[i] + indegree[i]); 78095fce210SBarry Smith #else 78195fce210SBarry Smith for (i=0; i<sf->nroots; i++) inoffset[i+1] = inoffset[i] + indegree[i]; 78295fce210SBarry Smith #endif 78395fce210SBarry Smith for (i=0; i<sf->nleaves; i++) outones[i] = 1; 78495fce210SBarry Smith ierr = PetscSFFetchAndOpBegin(sf,MPIU_INT,inoffset,outones,outoffset,MPIU_SUM);CHKERRQ(ierr); 78595fce210SBarry Smith ierr = PetscSFFetchAndOpEnd(sf,MPIU_INT,inoffset,outones,outoffset,MPIU_SUM);CHKERRQ(ierr); 78695fce210SBarry Smith for (i=0; i<sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */ 78795fce210SBarry Smith #if 0 78895fce210SBarry Smith #if defined(PETSC_USE_DEBUG) /* Check that the expected number of increments occurred */ 78995fce210SBarry Smith for (i=0; i<sf->nroots; i++) { 79095fce210SBarry Smith if (inoffset[i] + indegree[i] != inoffset[i+1]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect result after PetscSFFetchAndOp"); 79195fce210SBarry Smith } 79295fce210SBarry Smith #endif 79395fce210SBarry Smith #endif 79495fce210SBarry Smith ierr = PetscMalloc(sf->nleaves*sizeof(*remote),&remote);CHKERRQ(ierr); 79595fce210SBarry Smith for (i=0; i<sf->nleaves; i++) { 79695fce210SBarry Smith remote[i].rank = sf->remote[i].rank; 79795fce210SBarry Smith remote[i].index = outoffset[i]; 79895fce210SBarry Smith } 79995fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr); 80095fce210SBarry Smith ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr); 80195fce210SBarry Smith if (sf->rankorder) { /* Sort the ranks */ 80295fce210SBarry Smith PetscMPIInt rank; 80395fce210SBarry Smith PetscInt *inranks,*newoffset,*outranks,*newoutoffset,*tmpoffset,maxdegree; 80495fce210SBarry Smith PetscSFNode *newremote; 80595fce210SBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); 80695fce210SBarry Smith for (i=0,maxdegree=0; i<sf->nroots; i++) maxdegree = PetscMax(maxdegree,indegree[i]); 80795fce210SBarry Smith ierr = PetscMalloc5(sf->multi->nroots,PetscInt,&inranks,sf->multi->nroots,PetscInt,&newoffset,sf->nleaves,PetscInt,&outranks,sf->nleaves,PetscInt,&newoutoffset,maxdegree,PetscInt,&tmpoffset);CHKERRQ(ierr); 80895fce210SBarry Smith for (i=0; i<sf->nleaves; i++) outranks[i] = rank; 8098bfbc91cSJed Brown ierr = PetscSFReduceBegin(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr); 8108bfbc91cSJed Brown ierr = PetscSFReduceEnd(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr); 81195fce210SBarry Smith /* Sort the incoming ranks at each vertex, build the inverse map */ 81295fce210SBarry Smith for (i=0; i<sf->nroots; i++) { 81395fce210SBarry Smith PetscInt j; 81495fce210SBarry Smith for (j=0; j<indegree[i]; j++) tmpoffset[j] = j; 81595fce210SBarry Smith ierr = PetscSortIntWithArray(indegree[i],inranks+inoffset[i],tmpoffset);CHKERRQ(ierr); 81695fce210SBarry Smith for (j=0; j<indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j; 81795fce210SBarry Smith } 81895fce210SBarry Smith ierr = PetscSFBcastBegin(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr); 81995fce210SBarry Smith ierr = PetscSFBcastEnd(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr); 82095fce210SBarry Smith ierr = PetscMalloc(sf->nleaves*sizeof(*newremote),&newremote);CHKERRQ(ierr); 82195fce210SBarry Smith for (i=0; i<sf->nleaves; i++) { 82295fce210SBarry Smith newremote[i].rank = sf->remote[i].rank; 82395fce210SBarry Smith newremote[i].index = newoutoffset[i]; 82495fce210SBarry Smith } 82595fce210SBarry Smith ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,NULL,PETSC_COPY_VALUES,newremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 82695fce210SBarry Smith ierr = PetscFree5(inranks,newoffset,outranks,newoutoffset,tmpoffset);CHKERRQ(ierr); 82795fce210SBarry Smith } 82895fce210SBarry Smith ierr = PetscFree3(inoffset,outones,outoffset);CHKERRQ(ierr); 82995fce210SBarry Smith } 83095fce210SBarry Smith *multi = sf->multi; 83195fce210SBarry Smith PetscFunctionReturn(0); 83295fce210SBarry Smith } 83395fce210SBarry Smith 83495fce210SBarry Smith #undef __FUNCT__ 83595fce210SBarry Smith #define __FUNCT__ "PetscSFCreateEmbeddedSF" 83695fce210SBarry Smith /*@C 83795fce210SBarry Smith PetscSFCreateEmbeddedSF - removes edges from all but the selected roots, does not remap indices 83895fce210SBarry Smith 83995fce210SBarry Smith Collective 84095fce210SBarry Smith 84195fce210SBarry Smith Input Arguments: 84295fce210SBarry Smith + sf - original star forest 84395fce210SBarry Smith . nroots - number of roots to select on this process 84495fce210SBarry Smith - selected - selected roots on this process 84595fce210SBarry Smith 84695fce210SBarry Smith Output Arguments: 84795fce210SBarry Smith . newsf - new star forest 84895fce210SBarry Smith 84995fce210SBarry Smith Level: advanced 85095fce210SBarry Smith 85195fce210SBarry Smith Note: 85295fce210SBarry Smith To use the new PetscSF, it may be necessary to know the indices of the leaves that are still participating. This can 85395fce210SBarry Smith be done by calling PetscSFGetGraph(). 85495fce210SBarry Smith 85595fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFGetGraph() 85695fce210SBarry Smith @*/ 85795fce210SBarry Smith PetscErrorCode PetscSFCreateEmbeddedSF(PetscSF sf,PetscInt nroots,const PetscInt *selected,PetscSF *newsf) 85895fce210SBarry Smith { 859*0511a646SMatthew G. Knepley PetscInt *rootdata, *leafdata, *ilocal; 86095fce210SBarry Smith PetscSFNode *iremote; 861*0511a646SMatthew G. Knepley PetscInt leafsize = 0, nleaves = 0, i; 862*0511a646SMatthew G. Knepley PetscErrorCode ierr; 86395fce210SBarry Smith 86495fce210SBarry Smith PetscFunctionBegin; 86595fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 86695fce210SBarry Smith if (nroots) PetscValidPointer(selected,3); 86795fce210SBarry Smith PetscValidPointer(newsf,4); 868*0511a646SMatthew G. Knepley if (sf->mine) for (i = 0; i < sf->nleaves; ++i) {leafsize = PetscMax(leafsize, sf->mine[i]+1);} 869*0511a646SMatthew G. Knepley else leafsize = sf->nleaves; 870*0511a646SMatthew G. Knepley ierr = PetscMalloc2(sf->nroots,PetscInt,&rootdata,leafsize,PetscInt,&leafdata);CHKERRQ(ierr); 87195fce210SBarry Smith ierr = PetscMemzero(rootdata,sf->nroots*sizeof(PetscInt));CHKERRQ(ierr); 872*0511a646SMatthew G. Knepley ierr = PetscMemzero(leafdata,leafsize*sizeof(PetscInt));CHKERRQ(ierr); 873*0511a646SMatthew G. Knepley for (i=0; i<nroots; ++i) rootdata[selected[i]] = 1; 87495fce210SBarry Smith ierr = PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr); 87595fce210SBarry Smith ierr = PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr); 87695fce210SBarry Smith 877*0511a646SMatthew G. Knepley for (i = 0; i < leafsize; ++i) nleaves += leafdata[i]; 87895fce210SBarry Smith ierr = PetscMalloc(nleaves*sizeof(PetscInt),&ilocal);CHKERRQ(ierr); 87995fce210SBarry Smith ierr = PetscMalloc(nleaves*sizeof(PetscSFNode),&iremote);CHKERRQ(ierr); 880*0511a646SMatthew G. Knepley for (i = 0, nleaves = 0; i < sf->nleaves; ++i) { 881*0511a646SMatthew G. Knepley const PetscInt lidx = sf->mine ? sf->mine[i] : i; 882*0511a646SMatthew G. Knepley 883*0511a646SMatthew G. Knepley if (leafdata[lidx]) { 884*0511a646SMatthew G. Knepley ilocal[nleaves] = lidx; 88595fce210SBarry Smith iremote[nleaves].rank = sf->remote[i].rank; 88695fce210SBarry Smith iremote[nleaves].index = sf->remote[i].index; 887*0511a646SMatthew G. Knepley ++nleaves; 88895fce210SBarry Smith } 88995fce210SBarry Smith } 89095fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,newsf);CHKERRQ(ierr); 89195fce210SBarry Smith ierr = PetscSFSetGraph(*newsf,sf->nroots,nleaves,ilocal,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 89295fce210SBarry Smith ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr); 89395fce210SBarry Smith PetscFunctionReturn(0); 89495fce210SBarry Smith } 89595fce210SBarry Smith 89695fce210SBarry Smith #undef __FUNCT__ 89795fce210SBarry Smith #define __FUNCT__ "PetscSFBcastBegin" 89895fce210SBarry Smith /*@C 89995fce210SBarry Smith PetscSFBcastBegin - begin pointwise broadcast to be concluded with call to PetscSFBcastEnd() 90095fce210SBarry Smith 90195fce210SBarry Smith Collective on PetscSF 90295fce210SBarry Smith 90395fce210SBarry Smith Input Arguments: 90495fce210SBarry Smith + sf - star forest on which to communicate 90595fce210SBarry Smith . unit - data type associated with each node 90695fce210SBarry Smith - rootdata - buffer to broadcast 90795fce210SBarry Smith 90895fce210SBarry Smith Output Arguments: 90995fce210SBarry Smith . leafdata - buffer to update with values from each leaf's respective root 91095fce210SBarry Smith 91195fce210SBarry Smith Level: intermediate 91295fce210SBarry Smith 91395fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFView(), PetscSFBcastEnd(), PetscSFReduceBegin() 91495fce210SBarry Smith @*/ 91595fce210SBarry Smith PetscErrorCode PetscSFBcastBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata) 91695fce210SBarry Smith { 91795fce210SBarry Smith PetscErrorCode ierr; 91895fce210SBarry Smith 91995fce210SBarry Smith PetscFunctionBegin; 92095fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 92195fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 92295fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 92395fce210SBarry Smith ierr = (*sf->ops->BcastBegin)(sf,unit,rootdata,leafdata);CHKERRQ(ierr); 92495fce210SBarry Smith PetscFunctionReturn(0); 92595fce210SBarry Smith } 92695fce210SBarry Smith 92795fce210SBarry Smith #undef __FUNCT__ 92895fce210SBarry Smith #define __FUNCT__ "PetscSFBcastEnd" 92995fce210SBarry Smith /*@C 93095fce210SBarry Smith PetscSFBcastEnd - end a broadcast operation started with PetscSFBcastBegin() 93195fce210SBarry Smith 93295fce210SBarry Smith Collective 93395fce210SBarry Smith 93495fce210SBarry Smith Input Arguments: 93595fce210SBarry Smith + sf - star forest 93695fce210SBarry Smith . unit - data type 93795fce210SBarry Smith - rootdata - buffer to broadcast 93895fce210SBarry Smith 93995fce210SBarry Smith Output Arguments: 94095fce210SBarry Smith . leafdata - buffer to update with values from each leaf's respective root 94195fce210SBarry Smith 94295fce210SBarry Smith Level: intermediate 94395fce210SBarry Smith 94495fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFReduceEnd() 94595fce210SBarry Smith @*/ 94695fce210SBarry Smith PetscErrorCode PetscSFBcastEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata) 94795fce210SBarry Smith { 94895fce210SBarry Smith PetscErrorCode ierr; 94995fce210SBarry Smith 95095fce210SBarry Smith PetscFunctionBegin; 95195fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 95295fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 95395fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 95495fce210SBarry Smith ierr = (*sf->ops->BcastEnd)(sf,unit,rootdata,leafdata);CHKERRQ(ierr); 95595fce210SBarry Smith PetscFunctionReturn(0); 95695fce210SBarry Smith } 95795fce210SBarry Smith 95895fce210SBarry Smith #undef __FUNCT__ 95995fce210SBarry Smith #define __FUNCT__ "PetscSFReduceBegin" 96095fce210SBarry Smith /*@C 96195fce210SBarry Smith PetscSFReduceBegin - begin reduction of leafdata into rootdata, to be completed with call to PetscSFReduceEnd() 96295fce210SBarry Smith 96395fce210SBarry Smith Collective 96495fce210SBarry Smith 96595fce210SBarry Smith Input Arguments: 96695fce210SBarry Smith + sf - star forest 96795fce210SBarry Smith . unit - data type 96895fce210SBarry Smith . leafdata - values to reduce 96995fce210SBarry Smith - op - reduction operation 97095fce210SBarry Smith 97195fce210SBarry Smith Output Arguments: 97295fce210SBarry Smith . rootdata - result of reduction of values from all leaves of each root 97395fce210SBarry Smith 97495fce210SBarry Smith Level: intermediate 97595fce210SBarry Smith 97695fce210SBarry Smith .seealso: PetscSFBcastBegin() 97795fce210SBarry Smith @*/ 97895fce210SBarry Smith PetscErrorCode PetscSFReduceBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op) 97995fce210SBarry Smith { 98095fce210SBarry Smith PetscErrorCode ierr; 98195fce210SBarry Smith 98295fce210SBarry Smith PetscFunctionBegin; 98395fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 98495fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 98595fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 98695fce210SBarry Smith ierr = (sf->ops->ReduceBegin)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr); 98795fce210SBarry Smith PetscFunctionReturn(0); 98895fce210SBarry Smith } 98995fce210SBarry Smith 99095fce210SBarry Smith #undef __FUNCT__ 99195fce210SBarry Smith #define __FUNCT__ "PetscSFReduceEnd" 99295fce210SBarry Smith /*@C 99395fce210SBarry Smith PetscSFReduceEnd - end a reduction operation started with PetscSFReduceBegin() 99495fce210SBarry Smith 99595fce210SBarry Smith Collective 99695fce210SBarry Smith 99795fce210SBarry Smith Input Arguments: 99895fce210SBarry Smith + sf - star forest 99995fce210SBarry Smith . unit - data type 100095fce210SBarry Smith . leafdata - values to reduce 100195fce210SBarry Smith - op - reduction operation 100295fce210SBarry Smith 100395fce210SBarry Smith Output Arguments: 100495fce210SBarry Smith . rootdata - result of reduction of values from all leaves of each root 100595fce210SBarry Smith 100695fce210SBarry Smith Level: intermediate 100795fce210SBarry Smith 100895fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFBcastEnd() 100995fce210SBarry Smith @*/ 101095fce210SBarry Smith PetscErrorCode PetscSFReduceEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op) 101195fce210SBarry Smith { 101295fce210SBarry Smith PetscErrorCode ierr; 101395fce210SBarry Smith 101495fce210SBarry Smith PetscFunctionBegin; 101595fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 101695fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 101795fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 101895fce210SBarry Smith ierr = (*sf->ops->ReduceEnd)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr); 101995fce210SBarry Smith PetscFunctionReturn(0); 102095fce210SBarry Smith } 102195fce210SBarry Smith 102295fce210SBarry Smith #undef __FUNCT__ 102395fce210SBarry Smith #define __FUNCT__ "PetscSFComputeDegreeBegin" 102495fce210SBarry Smith /*@C 102595fce210SBarry Smith PetscSFComputeDegreeBegin - begin computation of degree for each root vertex, to be completed with PetscSFComputeDegreeEnd() 102695fce210SBarry Smith 102795fce210SBarry Smith Collective 102895fce210SBarry Smith 102995fce210SBarry Smith Input Arguments: 103095fce210SBarry Smith . sf - star forest 103195fce210SBarry Smith 103295fce210SBarry Smith Output Arguments: 103395fce210SBarry Smith . degree - degree of each root vertex 103495fce210SBarry Smith 103595fce210SBarry Smith Level: advanced 103695fce210SBarry Smith 103795fce210SBarry Smith .seealso: PetscSFGatherBegin() 103895fce210SBarry Smith @*/ 103995fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf,const PetscInt **degree) 104095fce210SBarry Smith { 104195fce210SBarry Smith PetscErrorCode ierr; 104295fce210SBarry Smith 104395fce210SBarry Smith PetscFunctionBegin; 104495fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 104595fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 104695fce210SBarry Smith PetscValidPointer(degree,2); 104795fce210SBarry Smith if (!sf->degree) { 104895fce210SBarry Smith PetscInt i; 104995fce210SBarry Smith ierr = PetscMalloc(sf->nroots*sizeof(PetscInt),&sf->degree);CHKERRQ(ierr); 105095fce210SBarry Smith ierr = PetscMalloc(sf->nleaves*sizeof(PetscInt),&sf->degreetmp);CHKERRQ(ierr); 105195fce210SBarry Smith for (i=0; i<sf->nroots; i++) sf->degree[i] = 0; 105295fce210SBarry Smith for (i=0; i<sf->nleaves; i++) sf->degreetmp[i] = 1; 105395fce210SBarry Smith ierr = PetscSFReduceBegin(sf,MPIU_INT,sf->degreetmp,sf->degree,MPIU_SUM);CHKERRQ(ierr); 105495fce210SBarry Smith } 105595fce210SBarry Smith *degree = NULL; 105695fce210SBarry Smith PetscFunctionReturn(0); 105795fce210SBarry Smith } 105895fce210SBarry Smith 105995fce210SBarry Smith #undef __FUNCT__ 106095fce210SBarry Smith #define __FUNCT__ "PetscSFComputeDegreeEnd" 106195fce210SBarry Smith /*@C 106295fce210SBarry Smith PetscSFComputeDegreeEnd - complete computation of degree for each root vertex, started with PetscSFComputeDegreeBegin() 106395fce210SBarry Smith 106495fce210SBarry Smith Collective 106595fce210SBarry Smith 106695fce210SBarry Smith Input Arguments: 106795fce210SBarry Smith . sf - star forest 106895fce210SBarry Smith 106995fce210SBarry Smith Output Arguments: 107095fce210SBarry Smith . degree - degree of each root vertex 107195fce210SBarry Smith 107295fce210SBarry Smith Level: developer 107395fce210SBarry Smith 107495fce210SBarry Smith .seealso: 107595fce210SBarry Smith @*/ 107695fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf,const PetscInt **degree) 107795fce210SBarry Smith { 107895fce210SBarry Smith PetscErrorCode ierr; 107995fce210SBarry Smith 108095fce210SBarry Smith PetscFunctionBegin; 108195fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 108295fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 108395fce210SBarry Smith if (!sf->degreeknown) { 108495fce210SBarry Smith ierr = PetscSFReduceEnd(sf,MPIU_INT,sf->degreetmp,sf->degree,MPIU_SUM);CHKERRQ(ierr); 108595fce210SBarry Smith ierr = PetscFree(sf->degreetmp);CHKERRQ(ierr); 108695fce210SBarry Smith 108795fce210SBarry Smith sf->degreeknown = PETSC_TRUE; 108895fce210SBarry Smith } 108995fce210SBarry Smith *degree = sf->degree; 109095fce210SBarry Smith PetscFunctionReturn(0); 109195fce210SBarry Smith } 109295fce210SBarry Smith 109395fce210SBarry Smith #undef __FUNCT__ 109495fce210SBarry Smith #define __FUNCT__ "PetscSFFetchAndOpBegin" 109595fce210SBarry Smith /*@C 109695fce210SBarry Smith PetscSFFetchAndOpBegin - begin operation that fetches values from root and updates atomically by applying operation using my leaf value, to be completed with PetscSFFetchAndOpEnd() 109795fce210SBarry Smith 109895fce210SBarry Smith Collective 109995fce210SBarry Smith 110095fce210SBarry Smith Input Arguments: 110195fce210SBarry Smith + sf - star forest 110295fce210SBarry Smith . unit - data type 110395fce210SBarry Smith . leafdata - leaf values to use in reduction 110495fce210SBarry Smith - op - operation to use for reduction 110595fce210SBarry Smith 110695fce210SBarry Smith Output Arguments: 110795fce210SBarry Smith + rootdata - root values to be updated, input state is seen by first process to perform an update 110895fce210SBarry Smith - leafupdate - state at each leaf's respective root immediately prior to my atomic update 110995fce210SBarry Smith 111095fce210SBarry Smith Level: advanced 111195fce210SBarry Smith 111295fce210SBarry Smith Note: 111395fce210SBarry Smith The update is only atomic at the granularity provided by the hardware. Different roots referenced by the same process 111495fce210SBarry Smith might be updated in a different order. Furthermore, if a composite type is used for the unit datatype, atomicity is 111595fce210SBarry Smith not guaranteed across the whole vertex. Therefore, this function is mostly only used with primitive types such as 111695fce210SBarry Smith integers. 111795fce210SBarry Smith 111895fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFReduceBegin(), PetscSFSetGraph() 111995fce210SBarry Smith @*/ 112095fce210SBarry Smith PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op) 112195fce210SBarry Smith { 112295fce210SBarry Smith PetscErrorCode ierr; 112395fce210SBarry Smith 112495fce210SBarry Smith PetscFunctionBegin; 112595fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 112695fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 112795fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 112895fce210SBarry Smith ierr = (*sf->ops->FetchAndOpBegin)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr); 112995fce210SBarry Smith PetscFunctionReturn(0); 113095fce210SBarry Smith } 113195fce210SBarry Smith 113295fce210SBarry Smith #undef __FUNCT__ 113395fce210SBarry Smith #define __FUNCT__ "PetscSFFetchAndOpEnd" 113495fce210SBarry Smith /*@C 113595fce210SBarry 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 113695fce210SBarry Smith 113795fce210SBarry Smith Collective 113895fce210SBarry Smith 113995fce210SBarry Smith Input Arguments: 114095fce210SBarry Smith + sf - star forest 114195fce210SBarry Smith . unit - data type 114295fce210SBarry Smith . leafdata - leaf values to use in reduction 114395fce210SBarry Smith - op - operation to use for reduction 114495fce210SBarry Smith 114595fce210SBarry Smith Output Arguments: 114695fce210SBarry Smith + rootdata - root values to be updated, input state is seen by first process to perform an update 114795fce210SBarry Smith - leafupdate - state at each leaf's respective root immediately prior to my atomic update 114895fce210SBarry Smith 114995fce210SBarry Smith Level: advanced 115095fce210SBarry Smith 115195fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFReduceEnd(), PetscSFSetGraph() 115295fce210SBarry Smith @*/ 115395fce210SBarry Smith PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op) 115495fce210SBarry Smith { 115595fce210SBarry Smith PetscErrorCode ierr; 115695fce210SBarry Smith 115795fce210SBarry Smith PetscFunctionBegin; 115895fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 115995fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 116095fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 116195fce210SBarry Smith ierr = (*sf->ops->FetchAndOpEnd)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr); 116295fce210SBarry Smith PetscFunctionReturn(0); 116395fce210SBarry Smith } 116495fce210SBarry Smith 116595fce210SBarry Smith #undef __FUNCT__ 116695fce210SBarry Smith #define __FUNCT__ "PetscSFGatherBegin" 116795fce210SBarry Smith /*@C 116895fce210SBarry Smith PetscSFGatherBegin - begin pointwise gather of all leaves into multi-roots, to be completed with PetscSFGatherEnd() 116995fce210SBarry Smith 117095fce210SBarry Smith Collective 117195fce210SBarry Smith 117295fce210SBarry Smith Input Arguments: 117395fce210SBarry Smith + sf - star forest 117495fce210SBarry Smith . unit - data type 117595fce210SBarry Smith - leafdata - leaf data to gather to roots 117695fce210SBarry Smith 117795fce210SBarry Smith Output Argument: 117895fce210SBarry Smith . multirootdata - root buffer to gather into, amount of space per root is equal to its degree 117995fce210SBarry Smith 118095fce210SBarry Smith Level: intermediate 118195fce210SBarry Smith 118295fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin() 118395fce210SBarry Smith @*/ 118495fce210SBarry Smith PetscErrorCode PetscSFGatherBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata) 118595fce210SBarry Smith { 118695fce210SBarry Smith PetscErrorCode ierr; 118795fce210SBarry Smith PetscSF multi; 118895fce210SBarry Smith 118995fce210SBarry Smith PetscFunctionBegin; 119095fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 119195fce210SBarry Smith ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr); 11928bfbc91cSJed Brown ierr = PetscSFReduceBegin(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr); 119395fce210SBarry Smith PetscFunctionReturn(0); 119495fce210SBarry Smith } 119595fce210SBarry Smith 119695fce210SBarry Smith #undef __FUNCT__ 119795fce210SBarry Smith #define __FUNCT__ "PetscSFGatherEnd" 119895fce210SBarry Smith /*@C 119995fce210SBarry Smith PetscSFGatherEnd - ends pointwise gather operation that was started with PetscSFGatherBegin() 120095fce210SBarry Smith 120195fce210SBarry Smith Collective 120295fce210SBarry Smith 120395fce210SBarry Smith Input Arguments: 120495fce210SBarry Smith + sf - star forest 120595fce210SBarry Smith . unit - data type 120695fce210SBarry Smith - leafdata - leaf data to gather to roots 120795fce210SBarry Smith 120895fce210SBarry Smith Output Argument: 120995fce210SBarry Smith . multirootdata - root buffer to gather into, amount of space per root is equal to its degree 121095fce210SBarry Smith 121195fce210SBarry Smith Level: intermediate 121295fce210SBarry Smith 121395fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd() 121495fce210SBarry Smith @*/ 121595fce210SBarry Smith PetscErrorCode PetscSFGatherEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata) 121695fce210SBarry Smith { 121795fce210SBarry Smith PetscErrorCode ierr; 121895fce210SBarry Smith PetscSF multi; 121995fce210SBarry Smith 122095fce210SBarry Smith PetscFunctionBegin; 122195fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 122295fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 122395fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 122495fce210SBarry Smith ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr); 12258bfbc91cSJed Brown ierr = PetscSFReduceEnd(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr); 122695fce210SBarry Smith PetscFunctionReturn(0); 122795fce210SBarry Smith } 122895fce210SBarry Smith 122995fce210SBarry Smith #undef __FUNCT__ 123095fce210SBarry Smith #define __FUNCT__ "PetscSFScatterBegin" 123195fce210SBarry Smith /*@C 123295fce210SBarry Smith PetscSFScatterBegin - begin pointwise scatter operation from multi-roots to leaves, to be completed with PetscSFScatterEnd() 123395fce210SBarry Smith 123495fce210SBarry Smith Collective 123595fce210SBarry Smith 123695fce210SBarry Smith Input Arguments: 123795fce210SBarry Smith + sf - star forest 123895fce210SBarry Smith . unit - data type 123995fce210SBarry Smith - multirootdata - root buffer to send to each leaf, one unit of data per leaf 124095fce210SBarry Smith 124195fce210SBarry Smith Output Argument: 124295fce210SBarry Smith . leafdata - leaf data to be update with personal data from each respective root 124395fce210SBarry Smith 124495fce210SBarry Smith Level: intermediate 124595fce210SBarry Smith 124695fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin() 124795fce210SBarry Smith @*/ 124895fce210SBarry Smith PetscErrorCode PetscSFScatterBegin(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata) 124995fce210SBarry Smith { 125095fce210SBarry Smith PetscErrorCode ierr; 125195fce210SBarry Smith PetscSF multi; 125295fce210SBarry Smith 125395fce210SBarry Smith PetscFunctionBegin; 125495fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 125595fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 125695fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 125795fce210SBarry Smith ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr); 125895fce210SBarry Smith ierr = PetscSFBcastBegin(multi,unit,multirootdata,leafdata);CHKERRQ(ierr); 125995fce210SBarry Smith PetscFunctionReturn(0); 126095fce210SBarry Smith } 126195fce210SBarry Smith 126295fce210SBarry Smith #undef __FUNCT__ 126395fce210SBarry Smith #define __FUNCT__ "PetscSFScatterEnd" 126495fce210SBarry Smith /*@C 126595fce210SBarry Smith PetscSFScatterEnd - ends pointwise scatter operation that was started with PetscSFScatterBegin() 126695fce210SBarry Smith 126795fce210SBarry Smith Collective 126895fce210SBarry Smith 126995fce210SBarry Smith Input Arguments: 127095fce210SBarry Smith + sf - star forest 127195fce210SBarry Smith . unit - data type 127295fce210SBarry Smith - multirootdata - root buffer to send to each leaf, one unit of data per leaf 127395fce210SBarry Smith 127495fce210SBarry Smith Output Argument: 127595fce210SBarry Smith . leafdata - leaf data to be update with personal data from each respective root 127695fce210SBarry Smith 127795fce210SBarry Smith Level: intermediate 127895fce210SBarry Smith 127995fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd() 128095fce210SBarry Smith @*/ 128195fce210SBarry Smith PetscErrorCode PetscSFScatterEnd(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata) 128295fce210SBarry Smith { 128395fce210SBarry Smith PetscErrorCode ierr; 128495fce210SBarry Smith PetscSF multi; 128595fce210SBarry Smith 128695fce210SBarry Smith PetscFunctionBegin; 128795fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 128895fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 128995fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 129095fce210SBarry Smith ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr); 129195fce210SBarry Smith ierr = PetscSFBcastEnd(multi,unit,multirootdata,leafdata);CHKERRQ(ierr); 129295fce210SBarry Smith PetscFunctionReturn(0); 129395fce210SBarry Smith } 1294