1af0996ceSBarry Smith #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/ 295fce210SBarry Smith #include <petscctable.h> 395fce210SBarry Smith 495fce210SBarry Smith #if defined(PETSC_USE_DEBUG) 595fce210SBarry Smith # define PetscSFCheckGraphSet(sf,arg) do { \ 695fce210SBarry Smith if (PetscUnlikely(!(sf)->graphset)) \ 795fce210SBarry Smith SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetGraph() on argument %D \"%s\" before %s()",(arg),#sf,PETSC_FUNCTION_NAME); \ 895fce210SBarry Smith } while (0) 995fce210SBarry Smith #else 1095fce210SBarry Smith # define PetscSFCheckGraphSet(sf,arg) do {} while (0) 1195fce210SBarry Smith #endif 1295fce210SBarry Smith 1395fce210SBarry Smith const char *const PetscSFDuplicateOptions[] = {"CONFONLY","RANKS","GRAPH","PetscSFDuplicateOption","PETSCSF_DUPLICATE_",0}; 1495fce210SBarry Smith 158af6ec1cSBarry Smith /*@ 1695fce210SBarry Smith PetscSFCreate - create a star forest communication context 1795fce210SBarry Smith 18217044c2SLisandro Dalcin Collective on MPI_Comm 1995fce210SBarry Smith 2095fce210SBarry Smith Input Arguments: 2195fce210SBarry Smith . comm - communicator on which the star forest will operate 2295fce210SBarry Smith 2395fce210SBarry Smith Output Arguments: 2495fce210SBarry Smith . sf - new star forest context 2595fce210SBarry Smith 2695fce210SBarry Smith Level: intermediate 2795fce210SBarry Smith 2895fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFDestroy() 2995fce210SBarry Smith @*/ 3095fce210SBarry Smith PetscErrorCode PetscSFCreate(MPI_Comm comm,PetscSF *sf) 3195fce210SBarry Smith { 3295fce210SBarry Smith PetscErrorCode ierr; 3395fce210SBarry Smith PetscSF b; 3495fce210SBarry Smith 3595fce210SBarry Smith PetscFunctionBegin; 3695fce210SBarry Smith PetscValidPointer(sf,2); 37607a6623SBarry Smith ierr = PetscSFInitializePackage();CHKERRQ(ierr); 3895fce210SBarry Smith 3973107ff1SLisandro Dalcin ierr = PetscHeaderCreate(b,PETSCSF_CLASSID,"PetscSF","Star Forest","PetscSF",comm,PetscSFDestroy,PetscSFView);CHKERRQ(ierr); 4095fce210SBarry Smith 4195fce210SBarry Smith b->nroots = -1; 4295fce210SBarry Smith b->nleaves = -1; 4329046d53SLisandro Dalcin b->minleaf = PETSC_MAX_INT; 4429046d53SLisandro Dalcin b->maxleaf = PETSC_MIN_INT; 4595fce210SBarry Smith b->nranks = -1; 4695fce210SBarry Smith b->rankorder = PETSC_TRUE; 4795fce210SBarry Smith b->ingroup = MPI_GROUP_NULL; 4895fce210SBarry Smith b->outgroup = MPI_GROUP_NULL; 4995fce210SBarry Smith b->graphset = PETSC_FALSE; 5095fce210SBarry Smith 5195fce210SBarry Smith *sf = b; 5295fce210SBarry Smith PetscFunctionReturn(0); 5395fce210SBarry Smith } 5495fce210SBarry Smith 5529046d53SLisandro Dalcin /*@ 5695fce210SBarry Smith PetscSFReset - Reset a star forest so that different sizes or neighbors can be used 5795fce210SBarry Smith 5895fce210SBarry Smith Collective 5995fce210SBarry Smith 6095fce210SBarry Smith Input Arguments: 6195fce210SBarry Smith . sf - star forest 6295fce210SBarry Smith 6395fce210SBarry Smith Level: advanced 6495fce210SBarry Smith 6595fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFDestroy() 6695fce210SBarry Smith @*/ 6795fce210SBarry Smith PetscErrorCode PetscSFReset(PetscSF sf) 6895fce210SBarry Smith { 6995fce210SBarry Smith PetscErrorCode ierr; 7095fce210SBarry Smith 7195fce210SBarry Smith PetscFunctionBegin; 7295fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 7379715d56SJed Brown if (sf->ops->Reset) {ierr = (*sf->ops->Reset)(sf);CHKERRQ(ierr);} 7429046d53SLisandro Dalcin sf->nroots = -1; 7529046d53SLisandro Dalcin sf->nleaves = -1; 7629046d53SLisandro Dalcin sf->minleaf = PETSC_MAX_INT; 7729046d53SLisandro Dalcin sf->maxleaf = PETSC_MIN_INT; 7895fce210SBarry Smith sf->mine = NULL; 7995fce210SBarry Smith sf->remote = NULL; 8029046d53SLisandro Dalcin sf->graphset = PETSC_FALSE; 8129046d53SLisandro Dalcin ierr = PetscFree(sf->mine_alloc);CHKERRQ(ierr); 8295fce210SBarry Smith ierr = PetscFree(sf->remote_alloc);CHKERRQ(ierr); 8321c688dcSJed Brown sf->nranks = -1; 8429046d53SLisandro Dalcin ierr = PetscFree4(sf->ranks,sf->roffset,sf->rmine,sf->rremote);CHKERRQ(ierr); 8529046d53SLisandro Dalcin sf->degreeknown = PETSC_FALSE; 8695fce210SBarry Smith ierr = PetscFree(sf->degree);CHKERRQ(ierr); 8795fce210SBarry Smith if (sf->ingroup != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->ingroup);CHKERRQ(ierr);} 8895fce210SBarry Smith if (sf->outgroup != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->outgroup);CHKERRQ(ierr);} 8995fce210SBarry Smith ierr = PetscSFDestroy(&sf->multi);CHKERRQ(ierr); 9095fce210SBarry Smith sf->setupcalled = PETSC_FALSE; 9195fce210SBarry Smith PetscFunctionReturn(0); 9295fce210SBarry Smith } 9395fce210SBarry Smith 9495fce210SBarry Smith /*@C 9529046d53SLisandro Dalcin PetscSFSetType - Set the PetscSF communication implementation 9695fce210SBarry Smith 9795fce210SBarry Smith Collective on PetscSF 9895fce210SBarry Smith 9995fce210SBarry Smith Input Parameters: 10095fce210SBarry Smith + sf - the PetscSF context 10195fce210SBarry Smith - type - a known method 10295fce210SBarry Smith 10395fce210SBarry Smith Options Database Key: 10495fce210SBarry Smith . -sf_type <type> - Sets the method; use -help for a list 10595fce210SBarry Smith of available methods (for instance, window, pt2pt, neighbor) 10695fce210SBarry Smith 10795fce210SBarry Smith Notes: 10895fce210SBarry Smith See "include/petscsf.h" for available methods (for instance) 10995fce210SBarry Smith + PETSCSFWINDOW - MPI-2/3 one-sided 11095fce210SBarry Smith - PETSCSFBASIC - basic implementation using MPI-1 two-sided 11195fce210SBarry Smith 11295fce210SBarry Smith Level: intermediate 11395fce210SBarry Smith 11495fce210SBarry Smith .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); 13229046d53SLisandro Dalcin /* Destroy the previous PetscSF implementation context */ 13329046d53SLisandro Dalcin if (sf->ops->Destroy) {ierr = (*(sf)->ops->Destroy)(sf);CHKERRQ(ierr);} 13495fce210SBarry Smith ierr = PetscMemzero(sf->ops,sizeof(*sf->ops));CHKERRQ(ierr); 13595fce210SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)sf,type);CHKERRQ(ierr); 13695fce210SBarry Smith ierr = (*r)(sf);CHKERRQ(ierr); 13795fce210SBarry Smith PetscFunctionReturn(0); 13895fce210SBarry Smith } 13995fce210SBarry Smith 14029046d53SLisandro Dalcin /*@C 14129046d53SLisandro Dalcin PetscSFGetType - Get the PetscSF communication implementation 14229046d53SLisandro Dalcin 14329046d53SLisandro Dalcin Not Collective 14429046d53SLisandro Dalcin 14529046d53SLisandro Dalcin Input Parameter: 14629046d53SLisandro Dalcin . sf - the PetscSF context 14729046d53SLisandro Dalcin 14829046d53SLisandro Dalcin Output Parameter: 14929046d53SLisandro Dalcin . type - the PetscSF type name 15029046d53SLisandro Dalcin 15129046d53SLisandro Dalcin Level: intermediate 15229046d53SLisandro Dalcin 15329046d53SLisandro Dalcin .keywords: PetscSF, get, type 15429046d53SLisandro Dalcin .seealso: PetscSFSetType(), PetscSFCreate() 15529046d53SLisandro Dalcin @*/ 15629046d53SLisandro Dalcin PetscErrorCode PetscSFGetType(PetscSF sf, PetscSFType *type) 15729046d53SLisandro Dalcin { 15829046d53SLisandro Dalcin PetscFunctionBegin; 15929046d53SLisandro Dalcin PetscValidHeaderSpecific(sf, PETSCSF_CLASSID,1); 16029046d53SLisandro Dalcin PetscValidPointer(type,2); 16129046d53SLisandro Dalcin *type = ((PetscObject)sf)->type_name; 16229046d53SLisandro Dalcin PetscFunctionReturn(0); 16329046d53SLisandro Dalcin } 16429046d53SLisandro Dalcin 165d36d33e4SMatthew G. Knepley /*@ 16695fce210SBarry Smith PetscSFDestroy - destroy star forest 16795fce210SBarry Smith 16895fce210SBarry Smith Collective 16995fce210SBarry Smith 17095fce210SBarry Smith Input Arguments: 17195fce210SBarry Smith . sf - address of star forest 17295fce210SBarry Smith 17395fce210SBarry Smith Level: intermediate 17495fce210SBarry Smith 17595fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFReset() 17695fce210SBarry Smith @*/ 17795fce210SBarry Smith PetscErrorCode PetscSFDestroy(PetscSF *sf) 17895fce210SBarry Smith { 17995fce210SBarry Smith PetscErrorCode ierr; 18095fce210SBarry Smith 18195fce210SBarry Smith PetscFunctionBegin; 18295fce210SBarry Smith if (!*sf) PetscFunctionReturn(0); 18395fce210SBarry Smith PetscValidHeaderSpecific((*sf),PETSCSF_CLASSID,1); 18429046d53SLisandro Dalcin if (--((PetscObject)(*sf))->refct > 0) {*sf = NULL; PetscFunctionReturn(0);} 18595fce210SBarry Smith ierr = PetscSFReset(*sf);CHKERRQ(ierr); 18695fce210SBarry Smith if ((*sf)->ops->Destroy) {ierr = (*(*sf)->ops->Destroy)(*sf);CHKERRQ(ierr);} 18795fce210SBarry Smith ierr = PetscHeaderDestroy(sf);CHKERRQ(ierr); 18895fce210SBarry Smith PetscFunctionReturn(0); 18995fce210SBarry Smith } 19095fce210SBarry Smith 19195fce210SBarry Smith /*@ 19295fce210SBarry Smith PetscSFSetUp - set up communication structures 19395fce210SBarry Smith 19495fce210SBarry Smith Collective 19595fce210SBarry Smith 19695fce210SBarry Smith Input Arguments: 19795fce210SBarry Smith . sf - star forest communication object 19895fce210SBarry Smith 19995fce210SBarry Smith Level: beginner 20095fce210SBarry Smith 20195fce210SBarry Smith .seealso: PetscSFSetFromOptions(), PetscSFSetType() 20295fce210SBarry Smith @*/ 20395fce210SBarry Smith PetscErrorCode PetscSFSetUp(PetscSF sf) 20495fce210SBarry Smith { 20595fce210SBarry Smith PetscErrorCode ierr; 20695fce210SBarry Smith 20795fce210SBarry Smith PetscFunctionBegin; 20829046d53SLisandro Dalcin PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 20929046d53SLisandro Dalcin PetscSFCheckGraphSet(sf,1); 21095fce210SBarry Smith if (sf->setupcalled) PetscFunctionReturn(0); 21195fce210SBarry Smith if (!((PetscObject)sf)->type_name) {ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr);} 21229046d53SLisandro Dalcin ierr = PetscLogEventBegin(PETSCSF_SetUp,sf,0,0,0);CHKERRQ(ierr); 21395fce210SBarry Smith if (sf->ops->SetUp) {ierr = (*sf->ops->SetUp)(sf);CHKERRQ(ierr);} 21429046d53SLisandro Dalcin ierr = PetscLogEventEnd(PETSCSF_SetUp,sf,0,0,0);CHKERRQ(ierr); 21595fce210SBarry Smith sf->setupcalled = PETSC_TRUE; 21695fce210SBarry Smith PetscFunctionReturn(0); 21795fce210SBarry Smith } 21895fce210SBarry Smith 2198af6ec1cSBarry Smith /*@ 22095fce210SBarry Smith PetscSFSetFromOptions - set PetscSF options using the options database 22195fce210SBarry Smith 22295fce210SBarry Smith Logically Collective 22395fce210SBarry Smith 22495fce210SBarry Smith Input Arguments: 22595fce210SBarry Smith . sf - star forest 22695fce210SBarry Smith 22795fce210SBarry Smith Options Database Keys: 22860263706SJed Brown + -sf_type - implementation type, see PetscSFSetType() 22960263706SJed Brown - -sf_rank_order - sort composite points for gathers and scatters in rank order, gathers are non-deterministic otherwise 23095fce210SBarry Smith 23195fce210SBarry Smith Level: intermediate 23295fce210SBarry Smith 23395fce210SBarry Smith .keywords: KSP, set, from, options, database 23495fce210SBarry Smith 23595fce210SBarry Smith .seealso: PetscSFWindowSetSyncType() 23695fce210SBarry Smith @*/ 23795fce210SBarry Smith PetscErrorCode PetscSFSetFromOptions(PetscSF sf) 23895fce210SBarry Smith { 23995fce210SBarry Smith PetscSFType deft; 24095fce210SBarry Smith char type[256]; 24195fce210SBarry Smith PetscErrorCode ierr; 24295fce210SBarry Smith PetscBool flg; 24395fce210SBarry Smith 24495fce210SBarry Smith PetscFunctionBegin; 24595fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 24695fce210SBarry Smith ierr = PetscObjectOptionsBegin((PetscObject)sf);CHKERRQ(ierr); 24795fce210SBarry Smith deft = ((PetscObject)sf)->type_name ? ((PetscObject)sf)->type_name : PETSCSFBASIC; 24829046d53SLisandro Dalcin ierr = PetscOptionsFList("-sf_type","PetscSF implementation type","PetscSFSetType",PetscSFList,deft,type,sizeof(type),&flg);CHKERRQ(ierr); 24995fce210SBarry Smith ierr = PetscSFSetType(sf,flg ? type : deft);CHKERRQ(ierr); 25095fce210SBarry 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); 251e55864a3SBarry Smith if (sf->ops->SetFromOptions) {ierr = (*sf->ops->SetFromOptions)(PetscOptionsObject,sf);CHKERRQ(ierr);} 25295fce210SBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 25395fce210SBarry Smith PetscFunctionReturn(0); 25495fce210SBarry Smith } 25595fce210SBarry Smith 25629046d53SLisandro Dalcin /*@ 25795fce210SBarry Smith PetscSFSetRankOrder - sort multi-points for gathers and scatters by rank order 25895fce210SBarry Smith 25995fce210SBarry Smith Logically Collective 26095fce210SBarry Smith 26195fce210SBarry Smith Input Arguments: 26295fce210SBarry Smith + sf - star forest 26395fce210SBarry Smith - flg - PETSC_TRUE to sort, PETSC_FALSE to skip sorting (lower setup cost, but non-deterministic) 26495fce210SBarry Smith 26595fce210SBarry Smith Level: advanced 26695fce210SBarry Smith 26795fce210SBarry Smith .seealso: PetscSFGatherBegin(), PetscSFScatterBegin() 26895fce210SBarry Smith @*/ 26995fce210SBarry Smith PetscErrorCode PetscSFSetRankOrder(PetscSF sf,PetscBool flg) 27095fce210SBarry Smith { 27195fce210SBarry Smith 27295fce210SBarry Smith PetscFunctionBegin; 27395fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 27495fce210SBarry Smith PetscValidLogicalCollectiveBool(sf,flg,2); 27595fce210SBarry Smith if (sf->multi) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Rank ordering must be set before first call to PetscSFGatherBegin() or PetscSFScatterBegin()"); 27695fce210SBarry Smith sf->rankorder = flg; 27795fce210SBarry Smith PetscFunctionReturn(0); 27895fce210SBarry Smith } 27995fce210SBarry Smith 2808af6ec1cSBarry Smith /*@ 28195fce210SBarry Smith PetscSFSetGraph - Set a parallel star forest 28295fce210SBarry Smith 28395fce210SBarry Smith Collective 28495fce210SBarry Smith 28595fce210SBarry Smith Input Arguments: 28695fce210SBarry Smith + sf - star forest 28795fce210SBarry Smith . nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves) 28895fce210SBarry Smith . nleaves - number of leaf vertices on the current process, each of these references a root on any process 28995fce210SBarry Smith . ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage 29095fce210SBarry Smith . localmode - copy mode for ilocal 29195fce210SBarry Smith . iremote - remote locations of root vertices for each leaf on the current process 29295fce210SBarry Smith - remotemode - copy mode for iremote 29395fce210SBarry Smith 29495fce210SBarry Smith Level: intermediate 29595fce210SBarry Smith 29695452b02SPatrick Sanan Notes: 29795452b02SPatrick Sanan In Fortran you must use PETSC_COPY_VALUES for localmode and remotemode 29838ab3f8aSBarry Smith 2992ad1e87fSLisandro Dalcin Developers Note: Local indices which are the identity permutation in the range [0,nleaves) are discarded as they 3002ad1e87fSLisandro Dalcin encode contiguous storage. In such case, if localmode is PETSC_OWN_POINTER, the memory is deallocated as it is not 3012ad1e87fSLisandro Dalcin needed 3022ad1e87fSLisandro Dalcin 30395fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFGetGraph() 30495fce210SBarry Smith @*/ 30595fce210SBarry Smith PetscErrorCode PetscSFSetGraph(PetscSF sf,PetscInt nroots,PetscInt nleaves,const PetscInt *ilocal,PetscCopyMode localmode,const PetscSFNode *iremote,PetscCopyMode remotemode) 30695fce210SBarry Smith { 30795fce210SBarry Smith PetscErrorCode ierr; 30895fce210SBarry Smith 30995fce210SBarry Smith PetscFunctionBegin; 31095fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 31129046d53SLisandro Dalcin if (nleaves > 0 && ilocal) PetscValidIntPointer(ilocal,4); 31229046d53SLisandro Dalcin if (nleaves > 0) PetscValidPointer(iremote,6); 31329046d53SLisandro Dalcin if (nroots < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nroots %D, cannot be negative",nroots); 31495fce210SBarry Smith if (nleaves < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nleaves %D, cannot be negative",nleaves); 31529046d53SLisandro Dalcin 31695fce210SBarry Smith ierr = PetscSFReset(sf);CHKERRQ(ierr); 31729046d53SLisandro Dalcin ierr = PetscLogEventBegin(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr); 31829046d53SLisandro Dalcin 31995fce210SBarry Smith sf->nroots = nroots; 32095fce210SBarry Smith sf->nleaves = nleaves; 32129046d53SLisandro Dalcin 32229046d53SLisandro Dalcin if (nleaves && ilocal) { 32321c688dcSJed Brown PetscInt i; 32429046d53SLisandro Dalcin PetscInt minleaf = PETSC_MAX_INT; 32529046d53SLisandro Dalcin PetscInt maxleaf = PETSC_MIN_INT; 3262ad1e87fSLisandro Dalcin int contiguous = 1; 32729046d53SLisandro Dalcin for (i=0; i<nleaves; i++) { 32829046d53SLisandro Dalcin minleaf = PetscMin(minleaf,ilocal[i]); 32929046d53SLisandro Dalcin maxleaf = PetscMax(maxleaf,ilocal[i]); 3302ad1e87fSLisandro Dalcin contiguous &= (ilocal[i] == i); 33129046d53SLisandro Dalcin } 33229046d53SLisandro Dalcin sf->minleaf = minleaf; 33329046d53SLisandro Dalcin sf->maxleaf = maxleaf; 3342ad1e87fSLisandro Dalcin if (contiguous) { 3352ad1e87fSLisandro Dalcin if (localmode == PETSC_OWN_POINTER) { 3362ad1e87fSLisandro Dalcin ierr = PetscFree(ilocal);CHKERRQ(ierr); 3372ad1e87fSLisandro Dalcin } 3382ad1e87fSLisandro Dalcin ilocal = NULL; 3392ad1e87fSLisandro Dalcin } 34029046d53SLisandro Dalcin } else { 34129046d53SLisandro Dalcin sf->minleaf = 0; 34229046d53SLisandro Dalcin sf->maxleaf = nleaves - 1; 34329046d53SLisandro Dalcin } 34429046d53SLisandro Dalcin 34529046d53SLisandro Dalcin if (ilocal) { 34695fce210SBarry Smith switch (localmode) { 34795fce210SBarry Smith case PETSC_COPY_VALUES: 348785e854fSJed Brown ierr = PetscMalloc1(nleaves,&sf->mine_alloc);CHKERRQ(ierr); 34929046d53SLisandro Dalcin ierr = PetscMemcpy(sf->mine_alloc,ilocal,nleaves*sizeof(*ilocal));CHKERRQ(ierr); 35095fce210SBarry Smith sf->mine = sf->mine_alloc; 35195fce210SBarry Smith break; 35295fce210SBarry Smith case PETSC_OWN_POINTER: 35395fce210SBarry Smith sf->mine_alloc = (PetscInt*)ilocal; 35495fce210SBarry Smith sf->mine = sf->mine_alloc; 35595fce210SBarry Smith break; 35695fce210SBarry Smith case PETSC_USE_POINTER: 35729046d53SLisandro Dalcin sf->mine_alloc = NULL; 35895fce210SBarry Smith sf->mine = (PetscInt*)ilocal; 35995fce210SBarry Smith break; 36095fce210SBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown localmode"); 36195fce210SBarry Smith } 36295fce210SBarry Smith } 36329046d53SLisandro Dalcin 36495fce210SBarry Smith switch (remotemode) { 36595fce210SBarry Smith case PETSC_COPY_VALUES: 366785e854fSJed Brown ierr = PetscMalloc1(nleaves,&sf->remote_alloc);CHKERRQ(ierr); 36729046d53SLisandro Dalcin ierr = PetscMemcpy(sf->remote_alloc,iremote,nleaves*sizeof(*iremote));CHKERRQ(ierr); 36895fce210SBarry Smith sf->remote = sf->remote_alloc; 36995fce210SBarry Smith break; 37095fce210SBarry Smith case PETSC_OWN_POINTER: 37195fce210SBarry Smith sf->remote_alloc = (PetscSFNode*)iremote; 37295fce210SBarry Smith sf->remote = sf->remote_alloc; 37395fce210SBarry Smith break; 37495fce210SBarry Smith case PETSC_USE_POINTER: 37529046d53SLisandro Dalcin sf->remote_alloc = NULL; 37695fce210SBarry Smith sf->remote = (PetscSFNode*)iremote; 37795fce210SBarry Smith break; 37895fce210SBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown remotemode"); 37995fce210SBarry Smith } 38095fce210SBarry Smith 381563ffbabSMatthew G. Knepley ierr = PetscLogEventEnd(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr); 38229046d53SLisandro Dalcin sf->graphset = PETSC_TRUE; 38395fce210SBarry Smith PetscFunctionReturn(0); 38495fce210SBarry Smith } 38595fce210SBarry Smith 38629046d53SLisandro Dalcin /*@ 38795fce210SBarry Smith PetscSFCreateInverseSF - given a PetscSF in which all vertices have degree 1, creates the inverse map 38895fce210SBarry Smith 38995fce210SBarry Smith Collective 39095fce210SBarry Smith 39195fce210SBarry Smith Input Arguments: 39295fce210SBarry Smith . sf - star forest to invert 39395fce210SBarry Smith 39495fce210SBarry Smith Output Arguments: 39595fce210SBarry Smith . isf - inverse of sf 39695fce210SBarry Smith 39795fce210SBarry Smith Level: advanced 39895fce210SBarry Smith 39995fce210SBarry Smith Notes: 40095fce210SBarry Smith All roots must have degree 1. 40195fce210SBarry Smith 40295fce210SBarry Smith The local space may be a permutation, but cannot be sparse. 40395fce210SBarry Smith 40495fce210SBarry Smith .seealso: PetscSFSetGraph() 40595fce210SBarry Smith @*/ 40695fce210SBarry Smith PetscErrorCode PetscSFCreateInverseSF(PetscSF sf,PetscSF *isf) 40795fce210SBarry Smith { 40895fce210SBarry Smith PetscErrorCode ierr; 40995fce210SBarry Smith PetscMPIInt rank; 41095fce210SBarry Smith PetscInt i,nroots,nleaves,maxlocal,count,*newilocal; 41195fce210SBarry Smith const PetscInt *ilocal; 41295fce210SBarry Smith PetscSFNode *roots,*leaves; 41395fce210SBarry Smith 41495fce210SBarry Smith PetscFunctionBegin; 41529046d53SLisandro Dalcin PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 41629046d53SLisandro Dalcin PetscSFCheckGraphSet(sf,1); 41729046d53SLisandro Dalcin PetscValidPointer(isf,2); 41829046d53SLisandro Dalcin 41995fce210SBarry Smith ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);CHKERRQ(ierr); 42029046d53SLisandro Dalcin maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */ 42129046d53SLisandro Dalcin 42229046d53SLisandro Dalcin ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); 423ae9aee6dSMatthew G. Knepley ierr = PetscMalloc2(nroots,&roots,maxlocal,&leaves);CHKERRQ(ierr); 424ae9aee6dSMatthew G. Knepley for (i=0; i<maxlocal; i++) { 42595fce210SBarry Smith leaves[i].rank = rank; 42695fce210SBarry Smith leaves[i].index = i; 42795fce210SBarry Smith } 42895fce210SBarry Smith for (i=0; i <nroots; i++) { 42995fce210SBarry Smith roots[i].rank = -1; 43095fce210SBarry Smith roots[i].index = -1; 43195fce210SBarry Smith } 4328bfbc91cSJed Brown ierr = PetscSFReduceBegin(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr); 4338bfbc91cSJed Brown ierr = PetscSFReduceEnd(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr); 43495fce210SBarry Smith 43595fce210SBarry Smith /* Check whether our leaves are sparse */ 43695fce210SBarry Smith for (i=0,count=0; i<nroots; i++) if (roots[i].rank >= 0) count++; 43795fce210SBarry Smith if (count == nroots) newilocal = NULL; 43895fce210SBarry Smith else { /* Index for sparse leaves and compact "roots" array (which is to become our leaves). */ 439785e854fSJed Brown ierr = PetscMalloc1(count,&newilocal);CHKERRQ(ierr); 44095fce210SBarry Smith for (i=0,count=0; i<nroots; i++) { 44195fce210SBarry Smith if (roots[i].rank >= 0) { 44295fce210SBarry Smith newilocal[count] = i; 44395fce210SBarry Smith roots[count].rank = roots[i].rank; 44495fce210SBarry Smith roots[count].index = roots[i].index; 44595fce210SBarry Smith count++; 44695fce210SBarry Smith } 44795fce210SBarry Smith } 44895fce210SBarry Smith } 44995fce210SBarry Smith 45095fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,isf);CHKERRQ(ierr); 45195fce210SBarry Smith ierr = PetscSFSetGraph(*isf,maxlocal,count,newilocal,PETSC_OWN_POINTER,roots,PETSC_COPY_VALUES);CHKERRQ(ierr); 45295fce210SBarry Smith ierr = PetscFree2(roots,leaves);CHKERRQ(ierr); 45395fce210SBarry Smith PetscFunctionReturn(0); 45495fce210SBarry Smith } 45595fce210SBarry Smith 45695fce210SBarry Smith /*@ 45795fce210SBarry Smith PetscSFDuplicate - duplicate a PetscSF, optionally preserving rank connectivity and graph 45895fce210SBarry Smith 45995fce210SBarry Smith Collective 46095fce210SBarry Smith 46195fce210SBarry Smith Input Arguments: 46295fce210SBarry Smith + sf - communication object to duplicate 46395fce210SBarry Smith - opt - PETSCSF_DUPLICATE_CONFONLY, PETSCSF_DUPLICATE_RANKS, or PETSCSF_DUPLICATE_GRAPH (see PetscSFDuplicateOption) 46495fce210SBarry Smith 46595fce210SBarry Smith Output Arguments: 46695fce210SBarry Smith . newsf - new communication object 46795fce210SBarry Smith 46895fce210SBarry Smith Level: beginner 46995fce210SBarry Smith 47095fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetType(), PetscSFSetGraph() 47195fce210SBarry Smith @*/ 47295fce210SBarry Smith PetscErrorCode PetscSFDuplicate(PetscSF sf,PetscSFDuplicateOption opt,PetscSF *newsf) 47395fce210SBarry Smith { 47429046d53SLisandro Dalcin PetscSFType type; 47595fce210SBarry Smith PetscErrorCode ierr; 47695fce210SBarry Smith 47795fce210SBarry Smith PetscFunctionBegin; 47829046d53SLisandro Dalcin PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 47929046d53SLisandro Dalcin PetscValidLogicalCollectiveEnum(sf,opt,2); 48029046d53SLisandro Dalcin PetscValidPointer(newsf,3); 48195fce210SBarry Smith ierr = PetscSFCreate(PetscObjectComm((PetscObject)sf),newsf);CHKERRQ(ierr); 48229046d53SLisandro Dalcin ierr = PetscSFGetType(sf,&type);CHKERRQ(ierr); 48329046d53SLisandro Dalcin if (type) {ierr = PetscSFSetType(*newsf,type);CHKERRQ(ierr);} 48495fce210SBarry Smith if (opt == PETSCSF_DUPLICATE_GRAPH) { 48595fce210SBarry Smith PetscInt nroots,nleaves; 48695fce210SBarry Smith const PetscInt *ilocal; 48795fce210SBarry Smith const PetscSFNode *iremote; 48829046d53SLisandro Dalcin PetscSFCheckGraphSet(sf,1); 48995fce210SBarry Smith ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr); 49095fce210SBarry Smith ierr = PetscSFSetGraph(*newsf,nroots,nleaves,ilocal,PETSC_COPY_VALUES,iremote,PETSC_COPY_VALUES);CHKERRQ(ierr); 49195fce210SBarry Smith } 49229046d53SLisandro Dalcin if (sf->ops->Duplicate) {ierr = (*sf->ops->Duplicate)(sf,opt,*newsf);CHKERRQ(ierr);} 49395fce210SBarry Smith PetscFunctionReturn(0); 49495fce210SBarry Smith } 49595fce210SBarry Smith 49695fce210SBarry Smith /*@C 49795fce210SBarry Smith PetscSFGetGraph - Get the graph specifying a parallel star forest 49895fce210SBarry Smith 49995fce210SBarry Smith Not Collective 50095fce210SBarry Smith 50195fce210SBarry Smith Input Arguments: 50295fce210SBarry Smith . sf - star forest 50395fce210SBarry Smith 50495fce210SBarry Smith Output Arguments: 50595fce210SBarry Smith + nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves) 50695fce210SBarry Smith . nleaves - number of leaf vertices on the current process, each of these references a root on any process 50795fce210SBarry Smith . ilocal - locations of leaves in leafdata buffers 50895fce210SBarry Smith - iremote - remote locations of root vertices for each leaf on the current process 50995fce210SBarry Smith 510373e0d91SLisandro Dalcin Notes: 511373e0d91SLisandro Dalcin We are not currently requiring that the graph is set, thus returning nroots=-1 if it has not been set yet 512373e0d91SLisandro Dalcin 51395fce210SBarry Smith Level: intermediate 51495fce210SBarry Smith 51595fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph() 51695fce210SBarry Smith @*/ 51795fce210SBarry Smith PetscErrorCode PetscSFGetGraph(PetscSF sf,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote) 51895fce210SBarry Smith { 51995fce210SBarry Smith 52095fce210SBarry Smith PetscFunctionBegin; 52195fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 52295fce210SBarry Smith if (nroots) *nroots = sf->nroots; 52395fce210SBarry Smith if (nleaves) *nleaves = sf->nleaves; 52495fce210SBarry Smith if (ilocal) *ilocal = sf->mine; 52595fce210SBarry Smith if (iremote) *iremote = sf->remote; 52695fce210SBarry Smith PetscFunctionReturn(0); 52795fce210SBarry Smith } 52895fce210SBarry Smith 52929046d53SLisandro Dalcin /*@ 53095fce210SBarry Smith PetscSFGetLeafRange - Get the active leaf ranges 53195fce210SBarry Smith 53295fce210SBarry Smith Not Collective 53395fce210SBarry Smith 53495fce210SBarry Smith Input Arguments: 53595fce210SBarry Smith . sf - star forest 53695fce210SBarry Smith 53795fce210SBarry Smith Output Arguments: 53895fce210SBarry Smith + minleaf - minimum active leaf on this process 53995fce210SBarry Smith - maxleaf - maximum active leaf on this process 54095fce210SBarry Smith 54195fce210SBarry Smith Level: developer 54295fce210SBarry Smith 54395fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph(), PetscSFGetGraph() 54495fce210SBarry Smith @*/ 54595fce210SBarry Smith PetscErrorCode PetscSFGetLeafRange(PetscSF sf,PetscInt *minleaf,PetscInt *maxleaf) 54695fce210SBarry Smith { 54795fce210SBarry Smith 54895fce210SBarry Smith PetscFunctionBegin; 54995fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 55029046d53SLisandro Dalcin PetscSFCheckGraphSet(sf,1); 55195fce210SBarry Smith if (minleaf) *minleaf = sf->minleaf; 55295fce210SBarry Smith if (maxleaf) *maxleaf = sf->maxleaf; 55395fce210SBarry Smith PetscFunctionReturn(0); 55495fce210SBarry Smith } 55595fce210SBarry Smith 55695fce210SBarry Smith /*@C 55795fce210SBarry Smith PetscSFView - view a star forest 55895fce210SBarry Smith 55995fce210SBarry Smith Collective 56095fce210SBarry Smith 56195fce210SBarry Smith Input Arguments: 56295fce210SBarry Smith + sf - star forest 56395fce210SBarry Smith - viewer - viewer to display graph, for example PETSC_VIEWER_STDOUT_WORLD 56495fce210SBarry Smith 56595fce210SBarry Smith Level: beginner 56695fce210SBarry Smith 56795fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph() 56895fce210SBarry Smith @*/ 56995fce210SBarry Smith PetscErrorCode PetscSFView(PetscSF sf,PetscViewer viewer) 57095fce210SBarry Smith { 57195fce210SBarry Smith PetscErrorCode ierr; 57295fce210SBarry Smith PetscBool iascii; 57395fce210SBarry Smith PetscViewerFormat format; 57495fce210SBarry Smith 57595fce210SBarry Smith PetscFunctionBegin; 57695fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 57795fce210SBarry Smith if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sf),&viewer);CHKERRQ(ierr);} 57895fce210SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 57995fce210SBarry Smith PetscCheckSameComm(sf,1,viewer,2); 58080153354SVaclav Hapla if (sf->graphset) {ierr = PetscSFSetUp(sf);CHKERRQ(ierr);} 58195fce210SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 58295fce210SBarry Smith if (iascii) { 58395fce210SBarry Smith PetscMPIInt rank; 58481bfa7aaSJed Brown PetscInt ii,i,j; 58595fce210SBarry Smith 586dae58748SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)sf,viewer);CHKERRQ(ierr); 58795fce210SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 58895fce210SBarry Smith if (sf->ops->View) {ierr = (*sf->ops->View)(sf,viewer);CHKERRQ(ierr);} 58980153354SVaclav Hapla if (!sf->graphset) { 59080153354SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"PetscSFSetGraph() has not been called yet\n");CHKERRQ(ierr); 59180153354SVaclav Hapla ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 59280153354SVaclav Hapla PetscFunctionReturn(0); 59380153354SVaclav Hapla } 59495fce210SBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); 5951575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 59695fce210SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of roots=%D, leaves=%D, remote ranks=%D\n",rank,sf->nroots,sf->nleaves,sf->nranks);CHKERRQ(ierr); 59795fce210SBarry Smith for (i=0; i<sf->nleaves; i++) { 59895fce210SBarry 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); 59995fce210SBarry Smith } 60095fce210SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 60195fce210SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 60295fce210SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 60381bfa7aaSJed Brown PetscMPIInt *tmpranks,*perm; 60481bfa7aaSJed Brown ierr = PetscMalloc2(sf->nranks,&tmpranks,sf->nranks,&perm);CHKERRQ(ierr); 60581bfa7aaSJed Brown ierr = PetscMemcpy(tmpranks,sf->ranks,sf->nranks*sizeof(tmpranks[0]));CHKERRQ(ierr); 60681bfa7aaSJed Brown for (i=0; i<sf->nranks; i++) perm[i] = i; 60781bfa7aaSJed Brown ierr = PetscSortMPIIntWithArray(sf->nranks,tmpranks,perm);CHKERRQ(ierr); 60895fce210SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Roots referenced by my leaves, by rank\n",rank);CHKERRQ(ierr); 60981bfa7aaSJed Brown for (ii=0; ii<sf->nranks; ii++) { 61081bfa7aaSJed Brown i = perm[ii]; 6117904a332SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %d: %D edges\n",rank,sf->ranks[i],sf->roffset[i+1]-sf->roffset[i]);CHKERRQ(ierr); 61295fce210SBarry Smith for (j=sf->roffset[i]; j<sf->roffset[i+1]; j++) { 61395fce210SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D <- %D\n",rank,sf->rmine[j],sf->rremote[j]);CHKERRQ(ierr); 61495fce210SBarry Smith } 61595fce210SBarry Smith } 61681bfa7aaSJed Brown ierr = PetscFree2(tmpranks,perm);CHKERRQ(ierr); 61795fce210SBarry Smith } 61895fce210SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 6191575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 62095fce210SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 62195fce210SBarry Smith } 62295fce210SBarry Smith PetscFunctionReturn(0); 62395fce210SBarry Smith } 62495fce210SBarry Smith 62595fce210SBarry Smith /*@C 62695fce210SBarry Smith PetscSFGetRanks - Get ranks and number of vertices referenced by leaves on this process 62795fce210SBarry Smith 62895fce210SBarry Smith Not Collective 62995fce210SBarry Smith 63095fce210SBarry Smith Input Arguments: 63195fce210SBarry Smith . sf - star forest 63295fce210SBarry Smith 63395fce210SBarry Smith Output Arguments: 63495fce210SBarry Smith + nranks - number of ranks referenced by local part 63595fce210SBarry Smith . ranks - array of ranks 63695fce210SBarry Smith . roffset - offset in rmine/rremote for each rank (length nranks+1) 63795fce210SBarry Smith . rmine - concatenated array holding local indices referencing each remote rank 63895fce210SBarry Smith - rremote - concatenated array holding remote indices referenced for each remote rank 63995fce210SBarry Smith 64095fce210SBarry Smith Level: developer 64195fce210SBarry Smith 64295fce210SBarry Smith .seealso: PetscSFSetGraph() 64395fce210SBarry Smith @*/ 64495fce210SBarry Smith PetscErrorCode PetscSFGetRanks(PetscSF sf,PetscInt *nranks,const PetscMPIInt **ranks,const PetscInt **roffset,const PetscInt **rmine,const PetscInt **rremote) 64595fce210SBarry Smith { 64695fce210SBarry Smith 64795fce210SBarry Smith PetscFunctionBegin; 64895fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 64929046d53SLisandro Dalcin if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining ranks"); 65095fce210SBarry Smith if (nranks) *nranks = sf->nranks; 65195fce210SBarry Smith if (ranks) *ranks = sf->ranks; 65295fce210SBarry Smith if (roffset) *roffset = sf->roffset; 65395fce210SBarry Smith if (rmine) *rmine = sf->rmine; 65495fce210SBarry Smith if (rremote) *rremote = sf->rremote; 65595fce210SBarry Smith PetscFunctionReturn(0); 65695fce210SBarry Smith } 65795fce210SBarry Smith 658b5a8e515SJed Brown static PetscBool InList(PetscMPIInt needle,PetscMPIInt n,const PetscMPIInt *list) { 659b5a8e515SJed Brown PetscInt i; 660b5a8e515SJed Brown for (i=0; i<n; i++) { 661b5a8e515SJed Brown if (needle == list[i]) return PETSC_TRUE; 662b5a8e515SJed Brown } 663b5a8e515SJed Brown return PETSC_FALSE; 664b5a8e515SJed Brown } 665b5a8e515SJed Brown 66695fce210SBarry Smith /*@C 66721c688dcSJed Brown PetscSFSetUpRanks - Set up data structures associated with ranks; this is for internal use by PetscSF implementations. 66821c688dcSJed Brown 66921c688dcSJed Brown Collective 67021c688dcSJed Brown 67121c688dcSJed Brown Input Arguments: 672b5a8e515SJed Brown + sf - PetscSF to set up; PetscSFSetGraph() must have been called 673b5a8e515SJed Brown - dgroup - MPI_Group of ranks to be distinguished (e.g., for self or shared memory exchange) 67421c688dcSJed Brown 67521c688dcSJed Brown Level: developer 67621c688dcSJed Brown 67721c688dcSJed Brown .seealso: PetscSFGetRanks() 67821c688dcSJed Brown @*/ 679b5a8e515SJed Brown PetscErrorCode PetscSFSetUpRanks(PetscSF sf,MPI_Group dgroup) 68021c688dcSJed Brown { 68121c688dcSJed Brown PetscErrorCode ierr; 68221c688dcSJed Brown PetscTable table; 68321c688dcSJed Brown PetscTablePosition pos; 684b5a8e515SJed Brown PetscMPIInt size,groupsize,*groupranks; 68521c688dcSJed Brown PetscInt i,*rcount,*ranks; 68621c688dcSJed Brown 68721c688dcSJed Brown PetscFunctionBegin; 68821c688dcSJed Brown PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 68929046d53SLisandro Dalcin PetscSFCheckGraphSet(sf,1); 69021c688dcSJed Brown ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr); 69121c688dcSJed Brown ierr = PetscTableCreate(10,size,&table);CHKERRQ(ierr); 69221c688dcSJed Brown for (i=0; i<sf->nleaves; i++) { 69321c688dcSJed Brown /* Log 1-based rank */ 69421c688dcSJed Brown ierr = PetscTableAdd(table,sf->remote[i].rank+1,1,ADD_VALUES);CHKERRQ(ierr); 69521c688dcSJed Brown } 69621c688dcSJed Brown ierr = PetscTableGetCount(table,&sf->nranks);CHKERRQ(ierr); 69721c688dcSJed Brown ierr = PetscMalloc4(sf->nranks,&sf->ranks,sf->nranks+1,&sf->roffset,sf->nleaves,&sf->rmine,sf->nleaves,&sf->rremote);CHKERRQ(ierr); 69821c688dcSJed Brown ierr = PetscMalloc2(sf->nranks,&rcount,sf->nranks,&ranks);CHKERRQ(ierr); 69921c688dcSJed Brown ierr = PetscTableGetHeadPosition(table,&pos);CHKERRQ(ierr); 70021c688dcSJed Brown for (i=0; i<sf->nranks; i++) { 70121c688dcSJed Brown ierr = PetscTableGetNext(table,&pos,&ranks[i],&rcount[i]);CHKERRQ(ierr); 70221c688dcSJed Brown ranks[i]--; /* Convert back to 0-based */ 70321c688dcSJed Brown } 70421c688dcSJed Brown ierr = PetscTableDestroy(&table);CHKERRQ(ierr); 705b5a8e515SJed Brown 706b5a8e515SJed Brown /* We expect that dgroup is reliably "small" while nranks could be large */ 707b5a8e515SJed Brown { 7087fb8a5e4SKarl Rupp MPI_Group group = MPI_GROUP_NULL; 709b5a8e515SJed Brown PetscMPIInt *dgroupranks; 710b5a8e515SJed Brown ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr); 711b5a8e515SJed Brown ierr = MPI_Group_size(dgroup,&groupsize);CHKERRQ(ierr); 712b5a8e515SJed Brown ierr = PetscMalloc1(groupsize,&dgroupranks);CHKERRQ(ierr); 713b5a8e515SJed Brown ierr = PetscMalloc1(groupsize,&groupranks);CHKERRQ(ierr); 714b5a8e515SJed Brown for (i=0; i<groupsize; i++) dgroupranks[i] = i; 715f3fc9a17SJed Brown if (groupsize) {ierr = MPI_Group_translate_ranks(dgroup,groupsize,dgroupranks,group,groupranks);CHKERRQ(ierr);} 716b5a8e515SJed Brown ierr = MPI_Group_free(&group);CHKERRQ(ierr); 717b5a8e515SJed Brown ierr = PetscFree(dgroupranks);CHKERRQ(ierr); 718b5a8e515SJed Brown } 719b5a8e515SJed Brown 720b5a8e515SJed Brown /* Partition ranks[] into distinguished (first sf->ndranks) followed by non-distinguished */ 721b5a8e515SJed Brown for (sf->ndranks=0,i=sf->nranks; sf->ndranks<i; ) { 722b5a8e515SJed Brown for (i--; sf->ndranks<i; i--) { /* Scan i backward looking for distinguished rank */ 723b5a8e515SJed Brown if (InList(ranks[i],groupsize,groupranks)) break; 724b5a8e515SJed Brown } 725b5a8e515SJed Brown for ( ; sf->ndranks<=i; sf->ndranks++) { /* Scan sf->ndranks forward looking for non-distinguished rank */ 726b5a8e515SJed Brown if (!InList(ranks[sf->ndranks],groupsize,groupranks)) break; 727b5a8e515SJed Brown } 728b5a8e515SJed Brown if (sf->ndranks < i) { /* Swap ranks[sf->ndranks] with ranks[i] */ 729b5a8e515SJed Brown PetscInt tmprank,tmpcount; 730b5a8e515SJed Brown tmprank = ranks[i]; 731b5a8e515SJed Brown tmpcount = rcount[i]; 732b5a8e515SJed Brown ranks[i] = ranks[sf->ndranks]; 733b5a8e515SJed Brown rcount[i] = rcount[sf->ndranks]; 734b5a8e515SJed Brown ranks[sf->ndranks] = tmprank; 735b5a8e515SJed Brown rcount[sf->ndranks] = tmpcount; 736b5a8e515SJed Brown sf->ndranks++; 737b5a8e515SJed Brown } 738b5a8e515SJed Brown } 739b5a8e515SJed Brown ierr = PetscFree(groupranks);CHKERRQ(ierr); 740b5a8e515SJed Brown ierr = PetscSortIntWithArray(sf->ndranks,ranks,rcount);CHKERRQ(ierr); 741b5a8e515SJed Brown ierr = PetscSortIntWithArray(sf->nranks-sf->ndranks,ranks+sf->ndranks,rcount+sf->ndranks);CHKERRQ(ierr); 74221c688dcSJed Brown sf->roffset[0] = 0; 74321c688dcSJed Brown for (i=0; i<sf->nranks; i++) { 74421c688dcSJed Brown ierr = PetscMPIIntCast(ranks[i],sf->ranks+i);CHKERRQ(ierr); 74521c688dcSJed Brown sf->roffset[i+1] = sf->roffset[i] + rcount[i]; 74621c688dcSJed Brown rcount[i] = 0; 74721c688dcSJed Brown } 74821c688dcSJed Brown for (i=0; i<sf->nleaves; i++) { 749b5a8e515SJed Brown PetscInt irank; 75021c688dcSJed Brown /* Search for index of iremote[i].rank in sf->ranks */ 751b5a8e515SJed Brown ierr = PetscFindMPIInt(sf->remote[i].rank,sf->ndranks,sf->ranks,&irank);CHKERRQ(ierr); 752b5a8e515SJed Brown if (irank < 0) { 753b5a8e515SJed Brown ierr = PetscFindMPIInt(sf->remote[i].rank,sf->nranks-sf->ndranks,sf->ranks+sf->ndranks,&irank);CHKERRQ(ierr); 754b5a8e515SJed Brown if (irank >= 0) irank += sf->ndranks; 75521c688dcSJed Brown } 756b5a8e515SJed Brown if (irank < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Could not find rank %D in array",sf->remote[i].rank); 75721c688dcSJed Brown sf->rmine[sf->roffset[irank] + rcount[irank]] = sf->mine ? sf->mine[i] : i; 75821c688dcSJed Brown sf->rremote[sf->roffset[irank] + rcount[irank]] = sf->remote[i].index; 75921c688dcSJed Brown rcount[irank]++; 76021c688dcSJed Brown } 76121c688dcSJed Brown ierr = PetscFree2(rcount,ranks);CHKERRQ(ierr); 76221c688dcSJed Brown PetscFunctionReturn(0); 76321c688dcSJed Brown } 76421c688dcSJed Brown 76521c688dcSJed Brown /*@C 76695fce210SBarry Smith PetscSFGetGroups - gets incoming and outgoing process groups 76795fce210SBarry Smith 76895fce210SBarry Smith Collective 76995fce210SBarry Smith 77095fce210SBarry Smith Input Argument: 77195fce210SBarry Smith . sf - star forest 77295fce210SBarry Smith 77395fce210SBarry Smith Output Arguments: 77495fce210SBarry Smith + incoming - group of origin processes for incoming edges (leaves that reference my roots) 77595fce210SBarry Smith - outgoing - group of destination processes for outgoing edges (roots that I reference) 77695fce210SBarry Smith 77795fce210SBarry Smith Level: developer 77895fce210SBarry Smith 77995fce210SBarry Smith .seealso: PetscSFGetWindow(), PetscSFRestoreWindow() 78095fce210SBarry Smith @*/ 78195fce210SBarry Smith PetscErrorCode PetscSFGetGroups(PetscSF sf,MPI_Group *incoming,MPI_Group *outgoing) 78295fce210SBarry Smith { 78395fce210SBarry Smith PetscErrorCode ierr; 7847fb8a5e4SKarl Rupp MPI_Group group = MPI_GROUP_NULL; 78595fce210SBarry Smith 78695fce210SBarry Smith PetscFunctionBegin; 78729046d53SLisandro Dalcin if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining groups"); 78895fce210SBarry Smith if (sf->ingroup == MPI_GROUP_NULL) { 78995fce210SBarry Smith PetscInt i; 79095fce210SBarry Smith const PetscInt *indegree; 79195fce210SBarry Smith PetscMPIInt rank,*outranks,*inranks; 79295fce210SBarry Smith PetscSFNode *remote; 79395fce210SBarry Smith PetscSF bgcount; 79495fce210SBarry Smith 79595fce210SBarry Smith /* Compute the number of incoming ranks */ 796785e854fSJed Brown ierr = PetscMalloc1(sf->nranks,&remote);CHKERRQ(ierr); 79795fce210SBarry Smith for (i=0; i<sf->nranks; i++) { 79895fce210SBarry Smith remote[i].rank = sf->ranks[i]; 79995fce210SBarry Smith remote[i].index = 0; 80095fce210SBarry Smith } 80195fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,&bgcount);CHKERRQ(ierr); 80295fce210SBarry Smith ierr = PetscSFSetGraph(bgcount,1,sf->nranks,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr); 80395fce210SBarry Smith ierr = PetscSFComputeDegreeBegin(bgcount,&indegree);CHKERRQ(ierr); 80495fce210SBarry Smith ierr = PetscSFComputeDegreeEnd(bgcount,&indegree);CHKERRQ(ierr); 80595fce210SBarry Smith 80695fce210SBarry Smith /* Enumerate the incoming ranks */ 807dcca6d9dSJed Brown ierr = PetscMalloc2(indegree[0],&inranks,sf->nranks,&outranks);CHKERRQ(ierr); 80895fce210SBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); 80995fce210SBarry Smith for (i=0; i<sf->nranks; i++) outranks[i] = rank; 81095fce210SBarry Smith ierr = PetscSFGatherBegin(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr); 81195fce210SBarry Smith ierr = PetscSFGatherEnd(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr); 81295fce210SBarry Smith ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr); 81395fce210SBarry Smith ierr = MPI_Group_incl(group,indegree[0],inranks,&sf->ingroup);CHKERRQ(ierr); 81495fce210SBarry Smith ierr = MPI_Group_free(&group);CHKERRQ(ierr); 81595fce210SBarry Smith ierr = PetscFree2(inranks,outranks);CHKERRQ(ierr); 81695fce210SBarry Smith ierr = PetscSFDestroy(&bgcount);CHKERRQ(ierr); 81795fce210SBarry Smith } 81895fce210SBarry Smith *incoming = sf->ingroup; 81995fce210SBarry Smith 82095fce210SBarry Smith if (sf->outgroup == MPI_GROUP_NULL) { 82195fce210SBarry Smith ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr); 82295fce210SBarry Smith ierr = MPI_Group_incl(group,sf->nranks,sf->ranks,&sf->outgroup);CHKERRQ(ierr); 82395fce210SBarry Smith ierr = MPI_Group_free(&group);CHKERRQ(ierr); 82495fce210SBarry Smith } 82595fce210SBarry Smith *outgoing = sf->outgroup; 82695fce210SBarry Smith PetscFunctionReturn(0); 82795fce210SBarry Smith } 82895fce210SBarry Smith 82929046d53SLisandro Dalcin /*@ 83095fce210SBarry Smith PetscSFGetMultiSF - gets the inner SF implemeting gathers and scatters 83195fce210SBarry Smith 83295fce210SBarry Smith Collective 83395fce210SBarry Smith 83495fce210SBarry Smith Input Argument: 83595fce210SBarry Smith . sf - star forest that may contain roots with 0 or with more than 1 vertex 83695fce210SBarry Smith 83795fce210SBarry Smith Output Arguments: 83895fce210SBarry Smith . multi - star forest with split roots, such that each root has degree exactly 1 83995fce210SBarry Smith 84095fce210SBarry Smith Level: developer 84195fce210SBarry Smith 84295fce210SBarry Smith Notes: 84395fce210SBarry Smith 84495fce210SBarry Smith In most cases, users should use PetscSFGatherBegin() and PetscSFScatterBegin() instead of manipulating multi 84595fce210SBarry Smith directly. Since multi satisfies the stronger condition that each entry in the global space has exactly one incoming 84695fce210SBarry Smith edge, it is a candidate for future optimization that might involve its removal. 84795fce210SBarry Smith 848673100f5SVaclav Hapla .seealso: PetscSFSetGraph(), PetscSFGatherBegin(), PetscSFScatterBegin(), PetscSFComputeMultiRootOriginalNumbering() 84995fce210SBarry Smith @*/ 85095fce210SBarry Smith PetscErrorCode PetscSFGetMultiSF(PetscSF sf,PetscSF *multi) 85195fce210SBarry Smith { 85295fce210SBarry Smith PetscErrorCode ierr; 85395fce210SBarry Smith 85495fce210SBarry Smith PetscFunctionBegin; 85595fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 85695fce210SBarry Smith PetscValidPointer(multi,2); 85795fce210SBarry Smith if (sf->nroots < 0) { /* Graph has not been set yet; why do we need this? */ 85895fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr); 85995fce210SBarry Smith *multi = sf->multi; 86095fce210SBarry Smith PetscFunctionReturn(0); 86195fce210SBarry Smith } 86295fce210SBarry Smith if (!sf->multi) { 86395fce210SBarry Smith const PetscInt *indegree; 8649837ea96SMatthew G. Knepley PetscInt i,*inoffset,*outones,*outoffset,maxlocal; 86595fce210SBarry Smith PetscSFNode *remote; 86629046d53SLisandro Dalcin maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */ 86795fce210SBarry Smith ierr = PetscSFComputeDegreeBegin(sf,&indegree);CHKERRQ(ierr); 86895fce210SBarry Smith ierr = PetscSFComputeDegreeEnd(sf,&indegree);CHKERRQ(ierr); 8699837ea96SMatthew G. Knepley ierr = PetscMalloc3(sf->nroots+1,&inoffset,maxlocal,&outones,maxlocal,&outoffset);CHKERRQ(ierr); 87095fce210SBarry Smith inoffset[0] = 0; 87195fce210SBarry Smith for (i=0; i<sf->nroots; i++) inoffset[i+1] = inoffset[i] + indegree[i]; 8729837ea96SMatthew G. Knepley for (i=0; i<maxlocal; i++) outones[i] = 1; 873dbd2ff41SBarry Smith ierr = PetscSFFetchAndOpBegin(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr); 874dbd2ff41SBarry Smith ierr = PetscSFFetchAndOpEnd(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr); 87595fce210SBarry Smith for (i=0; i<sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */ 87695fce210SBarry Smith #if 0 87795fce210SBarry Smith #if defined(PETSC_USE_DEBUG) /* Check that the expected number of increments occurred */ 87895fce210SBarry Smith for (i=0; i<sf->nroots; i++) { 87995fce210SBarry Smith if (inoffset[i] + indegree[i] != inoffset[i+1]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect result after PetscSFFetchAndOp"); 88095fce210SBarry Smith } 88195fce210SBarry Smith #endif 88295fce210SBarry Smith #endif 883785e854fSJed Brown ierr = PetscMalloc1(sf->nleaves,&remote);CHKERRQ(ierr); 88495fce210SBarry Smith for (i=0; i<sf->nleaves; i++) { 88595fce210SBarry Smith remote[i].rank = sf->remote[i].rank; 88638e7336fSToby Isaac remote[i].index = outoffset[sf->mine ? sf->mine[i] : i]; 88795fce210SBarry Smith } 88895fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr); 88901365b40SToby Isaac ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr); 89095fce210SBarry Smith if (sf->rankorder) { /* Sort the ranks */ 89195fce210SBarry Smith PetscMPIInt rank; 89295fce210SBarry Smith PetscInt *inranks,*newoffset,*outranks,*newoutoffset,*tmpoffset,maxdegree; 89395fce210SBarry Smith PetscSFNode *newremote; 89495fce210SBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); 89595fce210SBarry Smith for (i=0,maxdegree=0; i<sf->nroots; i++) maxdegree = PetscMax(maxdegree,indegree[i]); 8969837ea96SMatthew G. Knepley ierr = PetscMalloc5(sf->multi->nroots,&inranks,sf->multi->nroots,&newoffset,maxlocal,&outranks,maxlocal,&newoutoffset,maxdegree,&tmpoffset);CHKERRQ(ierr); 8979837ea96SMatthew G. Knepley for (i=0; i<maxlocal; i++) outranks[i] = rank; 8988bfbc91cSJed Brown ierr = PetscSFReduceBegin(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr); 8998bfbc91cSJed Brown ierr = PetscSFReduceEnd(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr); 90095fce210SBarry Smith /* Sort the incoming ranks at each vertex, build the inverse map */ 90195fce210SBarry Smith for (i=0; i<sf->nroots; i++) { 90295fce210SBarry Smith PetscInt j; 90395fce210SBarry Smith for (j=0; j<indegree[i]; j++) tmpoffset[j] = j; 90495fce210SBarry Smith ierr = PetscSortIntWithArray(indegree[i],inranks+inoffset[i],tmpoffset);CHKERRQ(ierr); 90595fce210SBarry Smith for (j=0; j<indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j; 90695fce210SBarry Smith } 90795fce210SBarry Smith ierr = PetscSFBcastBegin(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr); 90895fce210SBarry Smith ierr = PetscSFBcastEnd(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr); 909785e854fSJed Brown ierr = PetscMalloc1(sf->nleaves,&newremote);CHKERRQ(ierr); 91095fce210SBarry Smith for (i=0; i<sf->nleaves; i++) { 91195fce210SBarry Smith newremote[i].rank = sf->remote[i].rank; 91201365b40SToby Isaac newremote[i].index = newoutoffset[sf->mine ? sf->mine[i] : i]; 91395fce210SBarry Smith } 91401365b40SToby Isaac ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,newremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 91595fce210SBarry Smith ierr = PetscFree5(inranks,newoffset,outranks,newoutoffset,tmpoffset);CHKERRQ(ierr); 91695fce210SBarry Smith } 91795fce210SBarry Smith ierr = PetscFree3(inoffset,outones,outoffset);CHKERRQ(ierr); 91895fce210SBarry Smith } 91995fce210SBarry Smith *multi = sf->multi; 92095fce210SBarry Smith PetscFunctionReturn(0); 92195fce210SBarry Smith } 92295fce210SBarry Smith 92395fce210SBarry Smith /*@C 92495fce210SBarry Smith PetscSFCreateEmbeddedSF - removes edges from all but the selected roots, does not remap indices 92595fce210SBarry Smith 92695fce210SBarry Smith Collective 92795fce210SBarry Smith 92895fce210SBarry Smith Input Arguments: 92995fce210SBarry Smith + sf - original star forest 93095fce210SBarry Smith . nroots - number of roots to select on this process 93195fce210SBarry Smith - selected - selected roots on this process 93295fce210SBarry Smith 93395fce210SBarry Smith Output Arguments: 93495fce210SBarry Smith . newsf - new star forest 93595fce210SBarry Smith 93695fce210SBarry Smith Level: advanced 93795fce210SBarry Smith 93895fce210SBarry Smith Note: 93995fce210SBarry Smith To use the new PetscSF, it may be necessary to know the indices of the leaves that are still participating. This can 94095fce210SBarry Smith be done by calling PetscSFGetGraph(). 94195fce210SBarry Smith 94295fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFGetGraph() 94395fce210SBarry Smith @*/ 94495fce210SBarry Smith PetscErrorCode PetscSFCreateEmbeddedSF(PetscSF sf,PetscInt nroots,const PetscInt *selected,PetscSF *newsf) 94595fce210SBarry Smith { 9460511a646SMatthew G. Knepley PetscInt *rootdata, *leafdata, *ilocal; 94795fce210SBarry Smith PetscSFNode *iremote; 948fc1ede2bSMatthew G. Knepley PetscInt leafsize = 0, nleaves = 0, n, i; 9490511a646SMatthew G. Knepley PetscErrorCode ierr; 95095fce210SBarry Smith 95195fce210SBarry Smith PetscFunctionBegin; 95295fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 95329046d53SLisandro Dalcin PetscSFCheckGraphSet(sf,1); 95495fce210SBarry Smith if (nroots) PetscValidPointer(selected,3); 95595fce210SBarry Smith PetscValidPointer(newsf,4); 9560511a646SMatthew G. Knepley if (sf->mine) for (i = 0; i < sf->nleaves; ++i) {leafsize = PetscMax(leafsize, sf->mine[i]+1);} 9570511a646SMatthew G. Knepley else leafsize = sf->nleaves; 9581795a4d1SJed Brown ierr = PetscCalloc2(sf->nroots,&rootdata,leafsize,&leafdata);CHKERRQ(ierr); 9590511a646SMatthew G. Knepley for (i=0; i<nroots; ++i) rootdata[selected[i]] = 1; 96095fce210SBarry Smith ierr = PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr); 96195fce210SBarry Smith ierr = PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr); 9620511a646SMatthew G. Knepley for (i = 0; i < leafsize; ++i) nleaves += leafdata[i]; 963785e854fSJed Brown ierr = PetscMalloc1(nleaves,&ilocal);CHKERRQ(ierr); 964785e854fSJed Brown ierr = PetscMalloc1(nleaves,&iremote);CHKERRQ(ierr); 965fc1ede2bSMatthew G. Knepley for (i = 0, n = 0; i < sf->nleaves; ++i) { 9660511a646SMatthew G. Knepley const PetscInt lidx = sf->mine ? sf->mine[i] : i; 9670511a646SMatthew G. Knepley 9680511a646SMatthew G. Knepley if (leafdata[lidx]) { 969fc1ede2bSMatthew G. Knepley ilocal[n] = lidx; 970fc1ede2bSMatthew G. Knepley iremote[n].rank = sf->remote[i].rank; 971fc1ede2bSMatthew G. Knepley iremote[n].index = sf->remote[i].index; 972fc1ede2bSMatthew G. Knepley ++n; 97395fce210SBarry Smith } 97495fce210SBarry Smith } 975fc1ede2bSMatthew 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); 97695fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,newsf);CHKERRQ(ierr); 97795fce210SBarry Smith ierr = PetscSFSetGraph(*newsf,sf->nroots,nleaves,ilocal,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 97895fce210SBarry Smith ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr); 97995fce210SBarry Smith PetscFunctionReturn(0); 98095fce210SBarry Smith } 98195fce210SBarry Smith 9822f5fb4c2SMatthew G. Knepley /*@C 9832f5fb4c2SMatthew G. Knepley PetscSFCreateEmbeddedLeafSF - removes edges from all but the selected leaves, does not remap indices 9842f5fb4c2SMatthew G. Knepley 9852f5fb4c2SMatthew G. Knepley Collective 9862f5fb4c2SMatthew G. Knepley 9872f5fb4c2SMatthew G. Knepley Input Arguments: 9882f5fb4c2SMatthew G. Knepley + sf - original star forest 9892f5fb4c2SMatthew G. Knepley . nleaves - number of leaves to select on this process 9902f5fb4c2SMatthew G. Knepley - selected - selected leaves on this process 9912f5fb4c2SMatthew G. Knepley 9922f5fb4c2SMatthew G. Knepley Output Arguments: 9932f5fb4c2SMatthew G. Knepley . newsf - new star forest 9942f5fb4c2SMatthew G. Knepley 9952f5fb4c2SMatthew G. Knepley Level: advanced 9962f5fb4c2SMatthew G. Knepley 9972f5fb4c2SMatthew G. Knepley .seealso: PetscSFCreateEmbeddedSF(), PetscSFSetGraph(), PetscSFGetGraph() 9982f5fb4c2SMatthew G. Knepley @*/ 9992f5fb4c2SMatthew G. Knepley PetscErrorCode PetscSFCreateEmbeddedLeafSF(PetscSF sf, PetscInt nleaves, const PetscInt *selected, PetscSF *newsf) 10002f5fb4c2SMatthew G. Knepley { 10012f5fb4c2SMatthew G. Knepley PetscSFNode *iremote; 10022f5fb4c2SMatthew G. Knepley PetscInt *ilocal; 10032f5fb4c2SMatthew G. Knepley PetscInt i; 10042f5fb4c2SMatthew G. Knepley PetscErrorCode ierr; 10052f5fb4c2SMatthew G. Knepley 10062f5fb4c2SMatthew G. Knepley PetscFunctionBegin; 10072f5fb4c2SMatthew G. Knepley PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 100829046d53SLisandro Dalcin PetscSFCheckGraphSet(sf, 1); 10092f5fb4c2SMatthew G. Knepley if (nleaves) PetscValidPointer(selected, 3); 10102f5fb4c2SMatthew G. Knepley PetscValidPointer(newsf, 4); 10112f5fb4c2SMatthew G. Knepley ierr = PetscMalloc1(nleaves, &ilocal);CHKERRQ(ierr); 10122f5fb4c2SMatthew G. Knepley ierr = PetscMalloc1(nleaves, &iremote);CHKERRQ(ierr); 10132f5fb4c2SMatthew G. Knepley for (i = 0; i < nleaves; ++i) { 10142f5fb4c2SMatthew G. Knepley const PetscInt l = selected[i]; 10152f5fb4c2SMatthew G. Knepley 10162f5fb4c2SMatthew G. Knepley ilocal[i] = sf->mine ? sf->mine[l] : l; 10172f5fb4c2SMatthew G. Knepley iremote[i].rank = sf->remote[l].rank; 10182f5fb4c2SMatthew G. Knepley iremote[i].index = sf->remote[l].index; 10192f5fb4c2SMatthew G. Knepley } 10202f5fb4c2SMatthew G. Knepley ierr = PetscSFDuplicate(sf, PETSCSF_DUPLICATE_RANKS, newsf);CHKERRQ(ierr); 10212f5fb4c2SMatthew G. Knepley ierr = PetscSFSetGraph(*newsf, sf->nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr); 10222f5fb4c2SMatthew G. Knepley PetscFunctionReturn(0); 10232f5fb4c2SMatthew G. Knepley } 10242f5fb4c2SMatthew G. Knepley 102595fce210SBarry Smith /*@C 102695fce210SBarry Smith PetscSFBcastBegin - begin pointwise broadcast to be concluded with call to PetscSFBcastEnd() 102795fce210SBarry Smith 102895fce210SBarry Smith Collective on PetscSF 102995fce210SBarry Smith 103095fce210SBarry Smith Input Arguments: 103195fce210SBarry Smith + sf - star forest on which to communicate 103295fce210SBarry Smith . unit - data type associated with each node 103395fce210SBarry Smith - rootdata - buffer to broadcast 103495fce210SBarry Smith 103595fce210SBarry Smith Output Arguments: 103695fce210SBarry Smith . leafdata - buffer to update with values from each leaf's respective root 103795fce210SBarry Smith 103895fce210SBarry Smith Level: intermediate 103995fce210SBarry Smith 104095fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFView(), PetscSFBcastEnd(), PetscSFReduceBegin() 104195fce210SBarry Smith @*/ 104295fce210SBarry Smith PetscErrorCode PetscSFBcastBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata) 104395fce210SBarry Smith { 104495fce210SBarry Smith PetscErrorCode ierr; 104595fce210SBarry Smith 104695fce210SBarry Smith PetscFunctionBegin; 104795fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 104895fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 104929046d53SLisandro Dalcin ierr = PetscLogEventBegin(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr); 105095fce210SBarry Smith ierr = (*sf->ops->BcastBegin)(sf,unit,rootdata,leafdata);CHKERRQ(ierr); 1051563ffbabSMatthew G. Knepley ierr = PetscLogEventEnd(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr); 105295fce210SBarry Smith PetscFunctionReturn(0); 105395fce210SBarry Smith } 105495fce210SBarry Smith 105595fce210SBarry Smith /*@C 105695fce210SBarry Smith PetscSFBcastEnd - end a broadcast operation started with PetscSFBcastBegin() 105795fce210SBarry Smith 105895fce210SBarry Smith Collective 105995fce210SBarry Smith 106095fce210SBarry Smith Input Arguments: 106195fce210SBarry Smith + sf - star forest 106295fce210SBarry Smith . unit - data type 106395fce210SBarry Smith - rootdata - buffer to broadcast 106495fce210SBarry Smith 106595fce210SBarry Smith Output Arguments: 106695fce210SBarry Smith . leafdata - buffer to update with values from each leaf's respective root 106795fce210SBarry Smith 106895fce210SBarry Smith Level: intermediate 106995fce210SBarry Smith 107095fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFReduceEnd() 107195fce210SBarry Smith @*/ 107295fce210SBarry Smith PetscErrorCode PetscSFBcastEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata) 107395fce210SBarry Smith { 107495fce210SBarry Smith PetscErrorCode ierr; 107595fce210SBarry Smith 107695fce210SBarry Smith PetscFunctionBegin; 107795fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 107895fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 107929046d53SLisandro Dalcin ierr = PetscLogEventBegin(PETSCSF_BcastEnd,sf,0,0,0);CHKERRQ(ierr); 108095fce210SBarry Smith ierr = (*sf->ops->BcastEnd)(sf,unit,rootdata,leafdata);CHKERRQ(ierr); 1081563ffbabSMatthew G. Knepley ierr = PetscLogEventEnd(PETSCSF_BcastEnd,sf,0,0,0);CHKERRQ(ierr); 108295fce210SBarry Smith PetscFunctionReturn(0); 108395fce210SBarry Smith } 108495fce210SBarry Smith 108595fce210SBarry Smith /*@C 1086*3482bfa8SJunchao Zhang PetscSFBcastAndOpBegin - begin pointwise broadcast with root value being reduced to leaf value, to be concluded with call to PetscSFBcastAndOpEnd() 1087*3482bfa8SJunchao Zhang 1088*3482bfa8SJunchao Zhang Collective on PetscSF 1089*3482bfa8SJunchao Zhang 1090*3482bfa8SJunchao Zhang Input Arguments: 1091*3482bfa8SJunchao Zhang + sf - star forest on which to communicate 1092*3482bfa8SJunchao Zhang . unit - data type associated with each node 1093*3482bfa8SJunchao Zhang . rootdata - buffer to broadcast 1094*3482bfa8SJunchao Zhang - op - operation to use for reduction 1095*3482bfa8SJunchao Zhang 1096*3482bfa8SJunchao Zhang Output Arguments: 1097*3482bfa8SJunchao Zhang . leafdata - buffer to be reduced with values from each leaf's respective root 1098*3482bfa8SJunchao Zhang 1099*3482bfa8SJunchao Zhang Level: intermediate 1100*3482bfa8SJunchao Zhang 1101*3482bfa8SJunchao Zhang .seealso: PetscSFBcastAndOpEnd(), PetscSFBcastBegin(), PetscSFBcastEnd() 1102*3482bfa8SJunchao Zhang @*/ 1103*3482bfa8SJunchao Zhang PetscErrorCode PetscSFBcastAndOpBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op) 1104*3482bfa8SJunchao Zhang { 1105*3482bfa8SJunchao Zhang PetscErrorCode ierr; 1106*3482bfa8SJunchao Zhang 1107*3482bfa8SJunchao Zhang PetscFunctionBegin; 1108*3482bfa8SJunchao Zhang PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 1109*3482bfa8SJunchao Zhang ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 1110*3482bfa8SJunchao Zhang ierr = PetscLogEventBegin(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr); 1111*3482bfa8SJunchao Zhang ierr = (*sf->ops->BcastAndOpBegin)(sf,unit,rootdata,leafdata,op);CHKERRQ(ierr); 1112*3482bfa8SJunchao Zhang ierr = PetscLogEventEnd(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr); 1113*3482bfa8SJunchao Zhang PetscFunctionReturn(0); 1114*3482bfa8SJunchao Zhang } 1115*3482bfa8SJunchao Zhang 1116*3482bfa8SJunchao Zhang /*@C 1117*3482bfa8SJunchao Zhang PetscSFBcastAndOpEnd - end a broadcast & reduce operation started with PetscSFBcastAndOpBegin() 1118*3482bfa8SJunchao Zhang 1119*3482bfa8SJunchao Zhang Collective 1120*3482bfa8SJunchao Zhang 1121*3482bfa8SJunchao Zhang Input Arguments: 1122*3482bfa8SJunchao Zhang + sf - star forest 1123*3482bfa8SJunchao Zhang . unit - data type 1124*3482bfa8SJunchao Zhang . rootdata - buffer to broadcast 1125*3482bfa8SJunchao Zhang - op - operation to use for reduction 1126*3482bfa8SJunchao Zhang 1127*3482bfa8SJunchao Zhang Output Arguments: 1128*3482bfa8SJunchao Zhang . leafdata - buffer to be reduced with values from each leaf's respective root 1129*3482bfa8SJunchao Zhang 1130*3482bfa8SJunchao Zhang Level: intermediate 1131*3482bfa8SJunchao Zhang 1132*3482bfa8SJunchao Zhang .seealso: PetscSFSetGraph(), PetscSFReduceEnd() 1133*3482bfa8SJunchao Zhang @*/ 1134*3482bfa8SJunchao Zhang PetscErrorCode PetscSFBcastAndOpEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op) 1135*3482bfa8SJunchao Zhang { 1136*3482bfa8SJunchao Zhang PetscErrorCode ierr; 1137*3482bfa8SJunchao Zhang 1138*3482bfa8SJunchao Zhang PetscFunctionBegin; 1139*3482bfa8SJunchao Zhang PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 1140*3482bfa8SJunchao Zhang ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 1141*3482bfa8SJunchao Zhang ierr = PetscLogEventBegin(PETSCSF_BcastAndOpEnd,sf,0,0,0);CHKERRQ(ierr); 1142*3482bfa8SJunchao Zhang ierr = (*sf->ops->BcastAndOpEnd)(sf,unit,rootdata,leafdata,op);CHKERRQ(ierr); 1143*3482bfa8SJunchao Zhang ierr = PetscLogEventEnd(PETSCSF_BcastAndOpEnd,sf,0,0,0);CHKERRQ(ierr); 1144*3482bfa8SJunchao Zhang PetscFunctionReturn(0); 1145*3482bfa8SJunchao Zhang } 1146*3482bfa8SJunchao Zhang 1147*3482bfa8SJunchao Zhang /*@C 114895fce210SBarry Smith PetscSFReduceBegin - begin reduction of leafdata into rootdata, to be completed with call to PetscSFReduceEnd() 114995fce210SBarry Smith 115095fce210SBarry Smith Collective 115195fce210SBarry Smith 115295fce210SBarry Smith Input Arguments: 115395fce210SBarry Smith + sf - star forest 115495fce210SBarry Smith . unit - data type 115595fce210SBarry Smith . leafdata - values to reduce 115695fce210SBarry Smith - op - reduction operation 115795fce210SBarry Smith 115895fce210SBarry Smith Output Arguments: 115995fce210SBarry Smith . rootdata - result of reduction of values from all leaves of each root 116095fce210SBarry Smith 116195fce210SBarry Smith Level: intermediate 116295fce210SBarry Smith 116395fce210SBarry Smith .seealso: PetscSFBcastBegin() 116495fce210SBarry Smith @*/ 116595fce210SBarry Smith PetscErrorCode PetscSFReduceBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op) 116695fce210SBarry Smith { 116795fce210SBarry Smith PetscErrorCode ierr; 116895fce210SBarry Smith 116995fce210SBarry Smith PetscFunctionBegin; 117095fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 117195fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 117229046d53SLisandro Dalcin ierr = PetscLogEventBegin(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr); 117395fce210SBarry Smith ierr = (sf->ops->ReduceBegin)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr); 1174563ffbabSMatthew G. Knepley ierr = PetscLogEventEnd(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr); 117595fce210SBarry Smith PetscFunctionReturn(0); 117695fce210SBarry Smith } 117795fce210SBarry Smith 117895fce210SBarry Smith /*@C 117995fce210SBarry Smith PetscSFReduceEnd - end a reduction operation started with PetscSFReduceBegin() 118095fce210SBarry Smith 118195fce210SBarry Smith Collective 118295fce210SBarry Smith 118395fce210SBarry Smith Input Arguments: 118495fce210SBarry Smith + sf - star forest 118595fce210SBarry Smith . unit - data type 118695fce210SBarry Smith . leafdata - values to reduce 118795fce210SBarry Smith - op - reduction operation 118895fce210SBarry Smith 118995fce210SBarry Smith Output Arguments: 119095fce210SBarry Smith . rootdata - result of reduction of values from all leaves of each root 119195fce210SBarry Smith 119295fce210SBarry Smith Level: intermediate 119395fce210SBarry Smith 119495fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFBcastEnd() 119595fce210SBarry Smith @*/ 119695fce210SBarry Smith PetscErrorCode PetscSFReduceEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op) 119795fce210SBarry Smith { 119895fce210SBarry Smith PetscErrorCode ierr; 119995fce210SBarry Smith 120095fce210SBarry Smith PetscFunctionBegin; 120195fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 120295fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 120329046d53SLisandro Dalcin ierr = PetscLogEventBegin(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr); 120495fce210SBarry Smith ierr = (*sf->ops->ReduceEnd)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr); 1205563ffbabSMatthew G. Knepley ierr = PetscLogEventEnd(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr); 120695fce210SBarry Smith PetscFunctionReturn(0); 120795fce210SBarry Smith } 120895fce210SBarry Smith 120995fce210SBarry Smith /*@C 121095fce210SBarry Smith PetscSFComputeDegreeBegin - begin computation of degree for each root vertex, to be completed with PetscSFComputeDegreeEnd() 121195fce210SBarry Smith 121295fce210SBarry Smith Collective 121395fce210SBarry Smith 121495fce210SBarry Smith Input Arguments: 121595fce210SBarry Smith . sf - star forest 121695fce210SBarry Smith 121795fce210SBarry Smith Output Arguments: 121895fce210SBarry Smith . degree - degree of each root vertex 121995fce210SBarry Smith 122095fce210SBarry Smith Level: advanced 122195fce210SBarry Smith 122295fce210SBarry Smith .seealso: PetscSFGatherBegin() 122395fce210SBarry Smith @*/ 122495fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf,const PetscInt **degree) 122595fce210SBarry Smith { 122695fce210SBarry Smith PetscErrorCode ierr; 122795fce210SBarry Smith 122895fce210SBarry Smith PetscFunctionBegin; 122995fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 123095fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 123195fce210SBarry Smith PetscValidPointer(degree,2); 1232803bd9e8SMatthew G. Knepley if (!sf->degreeknown) { 123329046d53SLisandro Dalcin PetscInt i, nroots = sf->nroots, maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */ 1234803bd9e8SMatthew G. Knepley if (sf->degree) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Calls to PetscSFComputeDegreeBegin() cannot be nested."); 123529046d53SLisandro Dalcin ierr = PetscMalloc1(nroots,&sf->degree);CHKERRQ(ierr); 123629046d53SLisandro Dalcin ierr = PetscMalloc1(PetscMax(maxlocal,1),&sf->degreetmp);CHKERRQ(ierr); /* allocate at least one entry, see check in PetscSFComputeDegreeEnd() */ 123729046d53SLisandro Dalcin for (i=0; i<nroots; i++) sf->degree[i] = 0; 12389837ea96SMatthew G. Knepley for (i=0; i<maxlocal; i++) sf->degreetmp[i] = 1; 1239dbd2ff41SBarry Smith ierr = PetscSFReduceBegin(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);CHKERRQ(ierr); 124095fce210SBarry Smith } 124195fce210SBarry Smith *degree = NULL; 124295fce210SBarry Smith PetscFunctionReturn(0); 124395fce210SBarry Smith } 124495fce210SBarry Smith 124595fce210SBarry Smith /*@C 124695fce210SBarry Smith PetscSFComputeDegreeEnd - complete computation of degree for each root vertex, started with PetscSFComputeDegreeBegin() 124795fce210SBarry Smith 124895fce210SBarry Smith Collective 124995fce210SBarry Smith 125095fce210SBarry Smith Input Arguments: 125195fce210SBarry Smith . sf - star forest 125295fce210SBarry Smith 125395fce210SBarry Smith Output Arguments: 125495fce210SBarry Smith . degree - degree of each root vertex 125595fce210SBarry Smith 125695fce210SBarry Smith Level: developer 125795fce210SBarry Smith 125895fce210SBarry Smith .seealso: 125995fce210SBarry Smith @*/ 126095fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf,const PetscInt **degree) 126195fce210SBarry Smith { 126295fce210SBarry Smith PetscErrorCode ierr; 126395fce210SBarry Smith 126495fce210SBarry Smith PetscFunctionBegin; 126595fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 126695fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 126729046d53SLisandro Dalcin PetscValidPointer(degree,2); 126895fce210SBarry Smith if (!sf->degreeknown) { 126929046d53SLisandro Dalcin if (!sf->degreetmp) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFComputeDegreeBegin() before PetscSFComputeDegreeEnd()"); 1270dbd2ff41SBarry Smith ierr = PetscSFReduceEnd(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);CHKERRQ(ierr); 127195fce210SBarry Smith ierr = PetscFree(sf->degreetmp);CHKERRQ(ierr); 127295fce210SBarry Smith sf->degreeknown = PETSC_TRUE; 127395fce210SBarry Smith } 127495fce210SBarry Smith *degree = sf->degree; 127595fce210SBarry Smith PetscFunctionReturn(0); 127695fce210SBarry Smith } 127795fce210SBarry Smith 1278673100f5SVaclav Hapla 1279673100f5SVaclav Hapla /*@C 1280673100f5SVaclav Hapla PetscSFComputeMultiRootOriginalNumbering - Returns original numbering of multi-roots (roots of multi-SF returned by PetscSFGetMultiSF()). Each multi-root is assigned index of its original root. 1281673100f5SVaclav Hapla 1282673100f5SVaclav Hapla Collective 1283673100f5SVaclav Hapla 1284673100f5SVaclav Hapla Input Arguments: 1285673100f5SVaclav Hapla + sf - star forest 1286673100f5SVaclav Hapla - degree - degree of each root vertex, computed with PetscSFComputeDegreeBegin()/PetscSFComputeDegreeEnd() 1287673100f5SVaclav Hapla 1288673100f5SVaclav Hapla Output Arguments: 1289673100f5SVaclav Hapla . mRootsOrigNumbering - original indices of multi-roots; length of the array is equal to the number of multi-roots (roots of multi-SF) 1290673100f5SVaclav Hapla 1291673100f5SVaclav Hapla Level: developer 1292673100f5SVaclav Hapla 1293673100f5SVaclav Hapla .seealso: PetscSFComputeDegreeBegin(), PetscSFComputeDegreeEnd(), PetscSFGetMultiSF() 1294673100f5SVaclav Hapla @*/ 1295673100f5SVaclav Hapla PetscErrorCode PetscSFComputeMultiRootOriginalNumbering(PetscSF sf, const PetscInt degree[], PetscInt *mRootsOrigNumbering[]) 1296673100f5SVaclav Hapla { 1297673100f5SVaclav Hapla PetscSF msf; 1298673100f5SVaclav Hapla PetscInt i, j, k, nroots, nmroots; 1299673100f5SVaclav Hapla PetscErrorCode ierr; 1300673100f5SVaclav Hapla 1301673100f5SVaclav Hapla PetscFunctionBegin; 1302673100f5SVaclav Hapla PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 1303673100f5SVaclav Hapla ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 130429328920SVaclav Hapla if (nroots) PetscValidIntPointer(degree,2); 130529328920SVaclav Hapla PetscValidPointer(mRootsOrigNumbering,3); 1306673100f5SVaclav Hapla ierr = PetscSFGetMultiSF(sf,&msf);CHKERRQ(ierr); 1307673100f5SVaclav Hapla ierr = PetscSFGetGraph(msf, &nmroots, NULL, NULL, NULL);CHKERRQ(ierr); 1308673100f5SVaclav Hapla ierr = PetscMalloc1(nmroots, mRootsOrigNumbering);CHKERRQ(ierr); 1309673100f5SVaclav Hapla for (i=0,j=0,k=0; i<nroots; i++) { 1310673100f5SVaclav Hapla if (!degree[i]) continue; 1311673100f5SVaclav Hapla for (j=0; j<degree[i]; j++,k++) { 1312673100f5SVaclav Hapla (*mRootsOrigNumbering)[k] = i; 1313673100f5SVaclav Hapla } 1314673100f5SVaclav Hapla } 1315673100f5SVaclav Hapla #if defined(PETSC_USE_DEBUG) 1316673100f5SVaclav Hapla if (PetscUnlikely(k != nmroots)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sanity check fail"); 1317673100f5SVaclav Hapla #endif 1318673100f5SVaclav Hapla PetscFunctionReturn(0); 1319673100f5SVaclav Hapla } 1320673100f5SVaclav Hapla 132195fce210SBarry Smith /*@C 132295fce210SBarry Smith PetscSFFetchAndOpBegin - begin operation that fetches values from root and updates atomically by applying operation using my leaf value, to be completed with PetscSFFetchAndOpEnd() 132395fce210SBarry Smith 132495fce210SBarry Smith Collective 132595fce210SBarry Smith 132695fce210SBarry Smith Input Arguments: 132795fce210SBarry Smith + sf - star forest 132895fce210SBarry Smith . unit - data type 132995fce210SBarry Smith . leafdata - leaf values to use in reduction 133095fce210SBarry Smith - op - operation to use for reduction 133195fce210SBarry Smith 133295fce210SBarry Smith Output Arguments: 133395fce210SBarry Smith + rootdata - root values to be updated, input state is seen by first process to perform an update 133495fce210SBarry Smith - leafupdate - state at each leaf's respective root immediately prior to my atomic update 133595fce210SBarry Smith 133695fce210SBarry Smith Level: advanced 133795fce210SBarry Smith 133895fce210SBarry Smith Note: 133995fce210SBarry Smith The update is only atomic at the granularity provided by the hardware. Different roots referenced by the same process 134095fce210SBarry Smith might be updated in a different order. Furthermore, if a composite type is used for the unit datatype, atomicity is 134195fce210SBarry Smith not guaranteed across the whole vertex. Therefore, this function is mostly only used with primitive types such as 134295fce210SBarry Smith integers. 134395fce210SBarry Smith 134495fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFReduceBegin(), PetscSFSetGraph() 134595fce210SBarry Smith @*/ 134695fce210SBarry Smith PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op) 134795fce210SBarry Smith { 134895fce210SBarry Smith PetscErrorCode ierr; 134995fce210SBarry Smith 135095fce210SBarry Smith PetscFunctionBegin; 135195fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 135295fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 135329046d53SLisandro Dalcin ierr = PetscLogEventBegin(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr); 135495fce210SBarry Smith ierr = (*sf->ops->FetchAndOpBegin)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr); 1355563ffbabSMatthew G. Knepley ierr = PetscLogEventEnd(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr); 135695fce210SBarry Smith PetscFunctionReturn(0); 135795fce210SBarry Smith } 135895fce210SBarry Smith 135995fce210SBarry Smith /*@C 136095fce210SBarry 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 136195fce210SBarry Smith 136295fce210SBarry Smith Collective 136395fce210SBarry Smith 136495fce210SBarry Smith Input Arguments: 136595fce210SBarry Smith + sf - star forest 136695fce210SBarry Smith . unit - data type 136795fce210SBarry Smith . leafdata - leaf values to use in reduction 136895fce210SBarry Smith - op - operation to use for reduction 136995fce210SBarry Smith 137095fce210SBarry Smith Output Arguments: 137195fce210SBarry Smith + rootdata - root values to be updated, input state is seen by first process to perform an update 137295fce210SBarry Smith - leafupdate - state at each leaf's respective root immediately prior to my atomic update 137395fce210SBarry Smith 137495fce210SBarry Smith Level: advanced 137595fce210SBarry Smith 137695fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFReduceEnd(), PetscSFSetGraph() 137795fce210SBarry Smith @*/ 137895fce210SBarry Smith PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op) 137995fce210SBarry Smith { 138095fce210SBarry Smith PetscErrorCode ierr; 138195fce210SBarry Smith 138295fce210SBarry Smith PetscFunctionBegin; 138395fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 138495fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 138529046d53SLisandro Dalcin ierr = PetscLogEventBegin(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr); 138695fce210SBarry Smith ierr = (*sf->ops->FetchAndOpEnd)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr); 1387563ffbabSMatthew G. Knepley ierr = PetscLogEventEnd(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr); 138895fce210SBarry Smith PetscFunctionReturn(0); 138995fce210SBarry Smith } 139095fce210SBarry Smith 139195fce210SBarry Smith /*@C 139295fce210SBarry Smith PetscSFGatherBegin - begin pointwise gather of all leaves into multi-roots, to be completed with PetscSFGatherEnd() 139395fce210SBarry Smith 139495fce210SBarry Smith Collective 139595fce210SBarry Smith 139695fce210SBarry Smith Input Arguments: 139795fce210SBarry Smith + sf - star forest 139895fce210SBarry Smith . unit - data type 139995fce210SBarry Smith - leafdata - leaf data to gather to roots 140095fce210SBarry Smith 140195fce210SBarry Smith Output Argument: 140295fce210SBarry Smith . multirootdata - root buffer to gather into, amount of space per root is equal to its degree 140395fce210SBarry Smith 140495fce210SBarry Smith Level: intermediate 140595fce210SBarry Smith 140695fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin() 140795fce210SBarry Smith @*/ 140895fce210SBarry Smith PetscErrorCode PetscSFGatherBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata) 140995fce210SBarry Smith { 141095fce210SBarry Smith PetscErrorCode ierr; 141195fce210SBarry Smith PetscSF multi; 141295fce210SBarry Smith 141395fce210SBarry Smith PetscFunctionBegin; 141495fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 141529046d53SLisandro Dalcin ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 141695fce210SBarry Smith ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr); 14178bfbc91cSJed Brown ierr = PetscSFReduceBegin(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr); 141895fce210SBarry Smith PetscFunctionReturn(0); 141995fce210SBarry Smith } 142095fce210SBarry Smith 142195fce210SBarry Smith /*@C 142295fce210SBarry Smith PetscSFGatherEnd - ends pointwise gather operation that was started with PetscSFGatherBegin() 142395fce210SBarry Smith 142495fce210SBarry Smith Collective 142595fce210SBarry Smith 142695fce210SBarry Smith Input Arguments: 142795fce210SBarry Smith + sf - star forest 142895fce210SBarry Smith . unit - data type 142995fce210SBarry Smith - leafdata - leaf data to gather to roots 143095fce210SBarry Smith 143195fce210SBarry Smith Output Argument: 143295fce210SBarry Smith . multirootdata - root buffer to gather into, amount of space per root is equal to its degree 143395fce210SBarry Smith 143495fce210SBarry Smith Level: intermediate 143595fce210SBarry Smith 143695fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd() 143795fce210SBarry Smith @*/ 143895fce210SBarry Smith PetscErrorCode PetscSFGatherEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata) 143995fce210SBarry Smith { 144095fce210SBarry Smith PetscErrorCode ierr; 144195fce210SBarry Smith PetscSF multi; 144295fce210SBarry Smith 144395fce210SBarry Smith PetscFunctionBegin; 144495fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 144595fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 144695fce210SBarry Smith ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr); 14478bfbc91cSJed Brown ierr = PetscSFReduceEnd(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr); 144895fce210SBarry Smith PetscFunctionReturn(0); 144995fce210SBarry Smith } 145095fce210SBarry Smith 145195fce210SBarry Smith /*@C 145295fce210SBarry Smith PetscSFScatterBegin - begin pointwise scatter operation from multi-roots to leaves, to be completed with PetscSFScatterEnd() 145395fce210SBarry Smith 145495fce210SBarry Smith Collective 145595fce210SBarry Smith 145695fce210SBarry Smith Input Arguments: 145795fce210SBarry Smith + sf - star forest 145895fce210SBarry Smith . unit - data type 145995fce210SBarry Smith - multirootdata - root buffer to send to each leaf, one unit of data per leaf 146095fce210SBarry Smith 146195fce210SBarry Smith Output Argument: 146295fce210SBarry Smith . leafdata - leaf data to be update with personal data from each respective root 146395fce210SBarry Smith 146495fce210SBarry Smith Level: intermediate 146595fce210SBarry Smith 146695fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin() 146795fce210SBarry Smith @*/ 146895fce210SBarry Smith PetscErrorCode PetscSFScatterBegin(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata) 146995fce210SBarry Smith { 147095fce210SBarry Smith PetscErrorCode ierr; 147195fce210SBarry Smith PetscSF multi; 147295fce210SBarry Smith 147395fce210SBarry Smith PetscFunctionBegin; 147495fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 147595fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 147695fce210SBarry Smith ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr); 147795fce210SBarry Smith ierr = PetscSFBcastBegin(multi,unit,multirootdata,leafdata);CHKERRQ(ierr); 147895fce210SBarry Smith PetscFunctionReturn(0); 147995fce210SBarry Smith } 148095fce210SBarry Smith 148195fce210SBarry Smith /*@C 148295fce210SBarry Smith PetscSFScatterEnd - ends pointwise scatter operation that was started with PetscSFScatterBegin() 148395fce210SBarry Smith 148495fce210SBarry Smith Collective 148595fce210SBarry Smith 148695fce210SBarry Smith Input Arguments: 148795fce210SBarry Smith + sf - star forest 148895fce210SBarry Smith . unit - data type 148995fce210SBarry Smith - multirootdata - root buffer to send to each leaf, one unit of data per leaf 149095fce210SBarry Smith 149195fce210SBarry Smith Output Argument: 149295fce210SBarry Smith . leafdata - leaf data to be update with personal data from each respective root 149395fce210SBarry Smith 149495fce210SBarry Smith Level: intermediate 149595fce210SBarry Smith 149695fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd() 149795fce210SBarry Smith @*/ 149895fce210SBarry Smith PetscErrorCode PetscSFScatterEnd(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata) 149995fce210SBarry Smith { 150095fce210SBarry Smith PetscErrorCode ierr; 150195fce210SBarry Smith PetscSF multi; 150295fce210SBarry Smith 150395fce210SBarry Smith PetscFunctionBegin; 150495fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 150595fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 150695fce210SBarry Smith ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr); 150795fce210SBarry Smith ierr = PetscSFBcastEnd(multi,unit,multirootdata,leafdata);CHKERRQ(ierr); 150895fce210SBarry Smith PetscFunctionReturn(0); 150995fce210SBarry Smith } 1510a7b3aa13SAta Mesgarnejad 1511a7b3aa13SAta Mesgarnejad /*@ 1512a7b3aa13SAta Mesgarnejad PetscSFCompose - Compose a new PetscSF equivalent to action to PetscSFs 1513a7b3aa13SAta Mesgarnejad 1514a7b3aa13SAta Mesgarnejad Input Parameters: 1515a7b3aa13SAta Mesgarnejad + sfA - The first PetscSF 1516a7b3aa13SAta Mesgarnejad - sfB - The second PetscSF 1517a7b3aa13SAta Mesgarnejad 1518a7b3aa13SAta Mesgarnejad Output Parameters: 1519a7b3aa13SAta Mesgarnejad . sfBA - equvalent PetscSF for applying A then B 1520a7b3aa13SAta Mesgarnejad 1521a7b3aa13SAta Mesgarnejad Level: developer 1522a7b3aa13SAta Mesgarnejad 1523a7b3aa13SAta Mesgarnejad .seealso: PetscSF, PetscSFGetGraph(), PetscSFSetGraph() 1524a7b3aa13SAta Mesgarnejad @*/ 1525a7b3aa13SAta Mesgarnejad PetscErrorCode PetscSFCompose(PetscSF sfA, PetscSF sfB, PetscSF *sfBA) 1526a7b3aa13SAta Mesgarnejad { 1527a7b3aa13SAta Mesgarnejad MPI_Comm comm; 1528a7b3aa13SAta Mesgarnejad const PetscSFNode *remotePointsA, *remotePointsB; 1529a7b3aa13SAta Mesgarnejad PetscSFNode *remotePointsBA; 1530a7b3aa13SAta Mesgarnejad const PetscInt *localPointsA, *localPointsB; 1531a7b3aa13SAta Mesgarnejad PetscInt numRootsA, numLeavesA, numRootsB, numLeavesB; 1532a7b3aa13SAta Mesgarnejad PetscErrorCode ierr; 1533a7b3aa13SAta Mesgarnejad 1534a7b3aa13SAta Mesgarnejad PetscFunctionBegin; 1535a7b3aa13SAta Mesgarnejad PetscValidHeaderSpecific(sfA, PETSCSF_CLASSID, 1); 153629046d53SLisandro Dalcin PetscSFCheckGraphSet(sfA, 1); 153729046d53SLisandro Dalcin PetscValidHeaderSpecific(sfB, PETSCSF_CLASSID, 2); 153829046d53SLisandro Dalcin PetscSFCheckGraphSet(sfB, 2); 153929046d53SLisandro Dalcin PetscValidPointer(sfBA, 3); 1540a7b3aa13SAta Mesgarnejad ierr = PetscObjectGetComm((PetscObject) sfA, &comm);CHKERRQ(ierr); 1541a7b3aa13SAta Mesgarnejad ierr = PetscSFGetGraph(sfA, &numRootsA, &numLeavesA, &localPointsA, &remotePointsA);CHKERRQ(ierr); 1542a7b3aa13SAta Mesgarnejad ierr = PetscSFGetGraph(sfB, &numRootsB, &numLeavesB, &localPointsB, &remotePointsB);CHKERRQ(ierr); 1543a7b3aa13SAta Mesgarnejad ierr = PetscMalloc1(numLeavesB, &remotePointsBA);CHKERRQ(ierr); 1544a7b3aa13SAta Mesgarnejad ierr = PetscSFBcastBegin(sfB, MPIU_2INT, remotePointsA, remotePointsBA);CHKERRQ(ierr); 1545a7b3aa13SAta Mesgarnejad ierr = PetscSFBcastEnd(sfB, MPIU_2INT, remotePointsA, remotePointsBA);CHKERRQ(ierr); 1546a7b3aa13SAta Mesgarnejad ierr = PetscSFCreate(comm, sfBA);CHKERRQ(ierr); 1547a7b3aa13SAta Mesgarnejad ierr = PetscSFSetGraph(*sfBA, numRootsA, numLeavesB, localPointsB, PETSC_COPY_VALUES, remotePointsBA, PETSC_OWN_POINTER);CHKERRQ(ierr); 1548a7b3aa13SAta Mesgarnejad PetscFunctionReturn(0); 1549a7b3aa13SAta Mesgarnejad } 1550