1af0996ceSBarry Smith #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/ 2c4e6a40aSLawrence Mitchell #include <petsc/private/hashseti.h> 395fce210SBarry Smith #include <petscctable.h> 495fce210SBarry Smith 595fce210SBarry Smith #if defined(PETSC_USE_DEBUG) 695fce210SBarry Smith # define PetscSFCheckGraphSet(sf,arg) do { \ 795fce210SBarry Smith if (PetscUnlikely(!(sf)->graphset)) \ 8dd5b3ca6SJunchao Zhang SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetGraph() or PetscSFSetGraphWithPattern() on argument %D \"%s\" before %s()",(arg),#sf,PETSC_FUNCTION_NAME); \ 995fce210SBarry Smith } while (0) 1095fce210SBarry Smith #else 1195fce210SBarry Smith # define PetscSFCheckGraphSet(sf,arg) do {} while (0) 1295fce210SBarry Smith #endif 1395fce210SBarry Smith 1495fce210SBarry Smith const char *const PetscSFDuplicateOptions[] = {"CONFONLY","RANKS","GRAPH","PetscSFDuplicateOption","PETSCSF_DUPLICATE_",0}; 1595fce210SBarry Smith 168af6ec1cSBarry Smith /*@ 1795fce210SBarry Smith PetscSFCreate - create a star forest communication context 1895fce210SBarry Smith 19d083f849SBarry Smith Collective 2095fce210SBarry Smith 2195fce210SBarry Smith Input Arguments: 2295fce210SBarry Smith . comm - communicator on which the star forest will operate 2395fce210SBarry Smith 2495fce210SBarry Smith Output Arguments: 2595fce210SBarry Smith . sf - new star forest context 2695fce210SBarry Smith 27dd5b3ca6SJunchao Zhang Options Database Keys: 28dd5b3ca6SJunchao Zhang + -sf_type basic -Use MPI persistent Isend/Irecv for communication (Default) 29dd5b3ca6SJunchao Zhang . -sf_type window -Use MPI-3 one-sided window for communication 30dd5b3ca6SJunchao Zhang - -sf_type neighbor -Use MPI-3 neighborhood collectives for communication 31dd5b3ca6SJunchao Zhang 3295fce210SBarry Smith Level: intermediate 3395fce210SBarry Smith 34dd5b3ca6SJunchao Zhang Notes: 35dd5b3ca6SJunchao Zhang When one knows the communication graph is one of the predefined graph, such as MPI_Alltoall, MPI_Allgatherv, 36dd5b3ca6SJunchao Zhang MPI_Gatherv, one can create a PetscSF and then set its graph with PetscSFSetGraphWithPattern(). These special 37dd5b3ca6SJunchao Zhang SFs are optimized and they have better performance than general SFs. 38dd5b3ca6SJunchao Zhang 39dd5b3ca6SJunchao Zhang .seealso: PetscSFSetGraph(), PetscSFSetGraphWithPattern(), PetscSFDestroy() 4095fce210SBarry Smith @*/ 4195fce210SBarry Smith PetscErrorCode PetscSFCreate(MPI_Comm comm,PetscSF *sf) 4295fce210SBarry Smith { 4395fce210SBarry Smith PetscErrorCode ierr; 4495fce210SBarry Smith PetscSF b; 4595fce210SBarry Smith 4695fce210SBarry Smith PetscFunctionBegin; 4795fce210SBarry Smith PetscValidPointer(sf,2); 48607a6623SBarry Smith ierr = PetscSFInitializePackage();CHKERRQ(ierr); 4995fce210SBarry Smith 5073107ff1SLisandro Dalcin ierr = PetscHeaderCreate(b,PETSCSF_CLASSID,"PetscSF","Star Forest","PetscSF",comm,PetscSFDestroy,PetscSFView);CHKERRQ(ierr); 5195fce210SBarry Smith 5295fce210SBarry Smith b->nroots = -1; 5395fce210SBarry Smith b->nleaves = -1; 5429046d53SLisandro Dalcin b->minleaf = PETSC_MAX_INT; 5529046d53SLisandro Dalcin b->maxleaf = PETSC_MIN_INT; 5695fce210SBarry Smith b->nranks = -1; 5795fce210SBarry Smith b->rankorder = PETSC_TRUE; 5895fce210SBarry Smith b->ingroup = MPI_GROUP_NULL; 5995fce210SBarry Smith b->outgroup = MPI_GROUP_NULL; 6095fce210SBarry Smith b->graphset = PETSC_FALSE; 6195fce210SBarry Smith 6295fce210SBarry Smith *sf = b; 6395fce210SBarry Smith PetscFunctionReturn(0); 6495fce210SBarry Smith } 6595fce210SBarry Smith 6629046d53SLisandro Dalcin /*@ 6795fce210SBarry Smith PetscSFReset - Reset a star forest so that different sizes or neighbors can be used 6895fce210SBarry Smith 6995fce210SBarry Smith Collective 7095fce210SBarry Smith 7195fce210SBarry Smith Input Arguments: 7295fce210SBarry Smith . sf - star forest 7395fce210SBarry Smith 7495fce210SBarry Smith Level: advanced 7595fce210SBarry Smith 7695fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFDestroy() 7795fce210SBarry Smith @*/ 7895fce210SBarry Smith PetscErrorCode PetscSFReset(PetscSF sf) 7995fce210SBarry Smith { 8095fce210SBarry Smith PetscErrorCode ierr; 8195fce210SBarry Smith 8295fce210SBarry Smith PetscFunctionBegin; 8395fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 8479715d56SJed Brown if (sf->ops->Reset) {ierr = (*sf->ops->Reset)(sf);CHKERRQ(ierr);} 8529046d53SLisandro Dalcin sf->nroots = -1; 8629046d53SLisandro Dalcin sf->nleaves = -1; 8729046d53SLisandro Dalcin sf->minleaf = PETSC_MAX_INT; 8829046d53SLisandro Dalcin sf->maxleaf = PETSC_MIN_INT; 8995fce210SBarry Smith sf->mine = NULL; 9095fce210SBarry Smith sf->remote = NULL; 9129046d53SLisandro Dalcin sf->graphset = PETSC_FALSE; 9229046d53SLisandro Dalcin ierr = PetscFree(sf->mine_alloc);CHKERRQ(ierr); 9395fce210SBarry Smith ierr = PetscFree(sf->remote_alloc);CHKERRQ(ierr); 9421c688dcSJed Brown sf->nranks = -1; 9529046d53SLisandro Dalcin ierr = PetscFree4(sf->ranks,sf->roffset,sf->rmine,sf->rremote);CHKERRQ(ierr); 96eb02082bSJunchao Zhang #if defined(PETSC_HAVE_CUDA) 97cd620004SJunchao Zhang { 98cd620004SJunchao Zhang PetscInt i; 99cd620004SJunchao Zhang for (i=0; i<2; i++) {if (sf->rmine_d[i]) {cudaError_t err = cudaFree(sf->rmine_d[i]);CHKERRCUDA(err);sf->rmine_d[i]=NULL;}} 100cd620004SJunchao Zhang } 101eb02082bSJunchao Zhang #endif 10229046d53SLisandro Dalcin sf->degreeknown = PETSC_FALSE; 10395fce210SBarry Smith ierr = PetscFree(sf->degree);CHKERRQ(ierr); 10495fce210SBarry Smith if (sf->ingroup != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->ingroup);CHKERRQ(ierr);} 10595fce210SBarry Smith if (sf->outgroup != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->outgroup);CHKERRQ(ierr);} 10695fce210SBarry Smith ierr = PetscSFDestroy(&sf->multi);CHKERRQ(ierr); 107dd5b3ca6SJunchao Zhang ierr = PetscLayoutDestroy(&sf->map);CHKERRQ(ierr); 10895fce210SBarry Smith sf->setupcalled = PETSC_FALSE; 10995fce210SBarry Smith PetscFunctionReturn(0); 11095fce210SBarry Smith } 11195fce210SBarry Smith 11295fce210SBarry Smith /*@C 11329046d53SLisandro Dalcin PetscSFSetType - Set the PetscSF communication implementation 11495fce210SBarry Smith 11595fce210SBarry Smith Collective on PetscSF 11695fce210SBarry Smith 11795fce210SBarry Smith Input Parameters: 11895fce210SBarry Smith + sf - the PetscSF context 11995fce210SBarry Smith - type - a known method 12095fce210SBarry Smith 12195fce210SBarry Smith Options Database Key: 12295fce210SBarry Smith . -sf_type <type> - Sets the method; use -help for a list 12370616304SStefano Zampini of available methods (for instance, window, basic, neighbor) 12495fce210SBarry Smith 12595fce210SBarry Smith Notes: 12695fce210SBarry Smith See "include/petscsf.h" for available methods (for instance) 12795fce210SBarry Smith + PETSCSFWINDOW - MPI-2/3 one-sided 12895fce210SBarry Smith - PETSCSFBASIC - basic implementation using MPI-1 two-sided 12995fce210SBarry Smith 13095fce210SBarry Smith Level: intermediate 13195fce210SBarry Smith 13295fce210SBarry Smith .seealso: PetscSFType, PetscSFCreate() 13395fce210SBarry Smith @*/ 13495fce210SBarry Smith PetscErrorCode PetscSFSetType(PetscSF sf,PetscSFType type) 13595fce210SBarry Smith { 13695fce210SBarry Smith PetscErrorCode ierr,(*r)(PetscSF); 13795fce210SBarry Smith PetscBool match; 13895fce210SBarry Smith 13995fce210SBarry Smith PetscFunctionBegin; 14095fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 14195fce210SBarry Smith PetscValidCharPointer(type,2); 14295fce210SBarry Smith 14395fce210SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)sf,type,&match);CHKERRQ(ierr); 14495fce210SBarry Smith if (match) PetscFunctionReturn(0); 14595fce210SBarry Smith 146adc40e5bSBarry Smith ierr = PetscFunctionListFind(PetscSFList,type,&r);CHKERRQ(ierr); 14795fce210SBarry Smith if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested PetscSF type %s",type); 14829046d53SLisandro Dalcin /* Destroy the previous PetscSF implementation context */ 14929046d53SLisandro Dalcin if (sf->ops->Destroy) {ierr = (*(sf)->ops->Destroy)(sf);CHKERRQ(ierr);} 15095fce210SBarry Smith ierr = PetscMemzero(sf->ops,sizeof(*sf->ops));CHKERRQ(ierr); 15195fce210SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)sf,type);CHKERRQ(ierr); 15295fce210SBarry Smith ierr = (*r)(sf);CHKERRQ(ierr); 15395fce210SBarry Smith PetscFunctionReturn(0); 15495fce210SBarry Smith } 15595fce210SBarry Smith 15629046d53SLisandro Dalcin /*@C 15729046d53SLisandro Dalcin PetscSFGetType - Get the PetscSF communication implementation 15829046d53SLisandro Dalcin 15929046d53SLisandro Dalcin Not Collective 16029046d53SLisandro Dalcin 16129046d53SLisandro Dalcin Input Parameter: 16229046d53SLisandro Dalcin . sf - the PetscSF context 16329046d53SLisandro Dalcin 16429046d53SLisandro Dalcin Output Parameter: 16529046d53SLisandro Dalcin . type - the PetscSF type name 16629046d53SLisandro Dalcin 16729046d53SLisandro Dalcin Level: intermediate 16829046d53SLisandro Dalcin 16929046d53SLisandro Dalcin .seealso: PetscSFSetType(), PetscSFCreate() 17029046d53SLisandro Dalcin @*/ 17129046d53SLisandro Dalcin PetscErrorCode PetscSFGetType(PetscSF sf, PetscSFType *type) 17229046d53SLisandro Dalcin { 17329046d53SLisandro Dalcin PetscFunctionBegin; 17429046d53SLisandro Dalcin PetscValidHeaderSpecific(sf, PETSCSF_CLASSID,1); 17529046d53SLisandro Dalcin PetscValidPointer(type,2); 17629046d53SLisandro Dalcin *type = ((PetscObject)sf)->type_name; 17729046d53SLisandro Dalcin PetscFunctionReturn(0); 17829046d53SLisandro Dalcin } 17929046d53SLisandro Dalcin 180d36d33e4SMatthew G. Knepley /*@ 18195fce210SBarry Smith PetscSFDestroy - destroy star forest 18295fce210SBarry Smith 18395fce210SBarry Smith Collective 18495fce210SBarry Smith 18595fce210SBarry Smith Input Arguments: 18695fce210SBarry Smith . sf - address of star forest 18795fce210SBarry Smith 18895fce210SBarry Smith Level: intermediate 18995fce210SBarry Smith 19095fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFReset() 19195fce210SBarry Smith @*/ 19295fce210SBarry Smith PetscErrorCode PetscSFDestroy(PetscSF *sf) 19395fce210SBarry Smith { 19495fce210SBarry Smith PetscErrorCode ierr; 19595fce210SBarry Smith 19695fce210SBarry Smith PetscFunctionBegin; 19795fce210SBarry Smith if (!*sf) PetscFunctionReturn(0); 19895fce210SBarry Smith PetscValidHeaderSpecific((*sf),PETSCSF_CLASSID,1); 19929046d53SLisandro Dalcin if (--((PetscObject)(*sf))->refct > 0) {*sf = NULL; PetscFunctionReturn(0);} 20095fce210SBarry Smith ierr = PetscSFReset(*sf);CHKERRQ(ierr); 20195fce210SBarry Smith if ((*sf)->ops->Destroy) {ierr = (*(*sf)->ops->Destroy)(*sf);CHKERRQ(ierr);} 20295fce210SBarry Smith ierr = PetscHeaderDestroy(sf);CHKERRQ(ierr); 20395fce210SBarry Smith PetscFunctionReturn(0); 20495fce210SBarry Smith } 20595fce210SBarry Smith 206c4e6a40aSLawrence Mitchell static PetscErrorCode PetscSFCheckGraphValid_Private(PetscSF sf) 207c4e6a40aSLawrence Mitchell { 208c4e6a40aSLawrence Mitchell PetscInt i, nleaves; 209c4e6a40aSLawrence Mitchell PetscMPIInt size; 210c4e6a40aSLawrence Mitchell const PetscInt *ilocal; 211c4e6a40aSLawrence Mitchell const PetscSFNode *iremote; 212c4e6a40aSLawrence Mitchell PetscErrorCode ierr; 213c4e6a40aSLawrence Mitchell 214c4e6a40aSLawrence Mitchell PetscFunctionBegin; 215*76bd3646SJed Brown if (!sf->graphset || !PetscDefined(USE_DEBUG)) PetscFunctionReturn(0); 216c4e6a40aSLawrence Mitchell ierr = PetscSFGetGraph(sf,NULL,&nleaves,&ilocal,&iremote);CHKERRQ(ierr); 217c4e6a40aSLawrence Mitchell ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr); 218c4e6a40aSLawrence Mitchell for (i = 0; i < nleaves; i++) { 219c4e6a40aSLawrence Mitchell const PetscInt rank = iremote[i].rank; 220c4e6a40aSLawrence Mitchell const PetscInt remote = iremote[i].index; 221c4e6a40aSLawrence Mitchell const PetscInt leaf = ilocal ? ilocal[i] : i; 222c4e6a40aSLawrence Mitchell if (rank < 0 || rank >= size) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Provided rank (%D) for remote %D is invalid, should be in [0, %d)",rank,i,size); 223c4e6a40aSLawrence Mitchell if (remote < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Provided index (%D) for remote %D is invalid, should be >= 0",remote,i); 224c4e6a40aSLawrence Mitchell if (leaf < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Provided location (%D) for leaf %D is invalid, should be >= 0",leaf,i); 225c4e6a40aSLawrence Mitchell } 226c4e6a40aSLawrence Mitchell PetscFunctionReturn(0); 227c4e6a40aSLawrence Mitchell } 228c4e6a40aSLawrence Mitchell 22995fce210SBarry Smith /*@ 23095fce210SBarry Smith PetscSFSetUp - set up communication structures 23195fce210SBarry Smith 23295fce210SBarry Smith Collective 23395fce210SBarry Smith 23495fce210SBarry Smith Input Arguments: 23595fce210SBarry Smith . sf - star forest communication object 23695fce210SBarry Smith 23795fce210SBarry Smith Level: beginner 23895fce210SBarry Smith 23995fce210SBarry Smith .seealso: PetscSFSetFromOptions(), PetscSFSetType() 24095fce210SBarry Smith @*/ 24195fce210SBarry Smith PetscErrorCode PetscSFSetUp(PetscSF sf) 24295fce210SBarry Smith { 24395fce210SBarry Smith PetscErrorCode ierr; 24495fce210SBarry Smith 24595fce210SBarry Smith PetscFunctionBegin; 24629046d53SLisandro Dalcin PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 24729046d53SLisandro Dalcin PetscSFCheckGraphSet(sf,1); 24895fce210SBarry Smith if (sf->setupcalled) PetscFunctionReturn(0); 249c4e6a40aSLawrence Mitchell ierr = PetscSFCheckGraphValid_Private(sf);CHKERRQ(ierr); 250c2a741eeSJunchao Zhang sf->use_gpu_aware_mpi = use_gpu_aware_mpi; 25195fce210SBarry Smith if (!((PetscObject)sf)->type_name) {ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr);} 25229046d53SLisandro Dalcin ierr = PetscLogEventBegin(PETSCSF_SetUp,sf,0,0,0);CHKERRQ(ierr); 25395fce210SBarry Smith if (sf->ops->SetUp) {ierr = (*sf->ops->SetUp)(sf);CHKERRQ(ierr);} 25429046d53SLisandro Dalcin ierr = PetscLogEventEnd(PETSCSF_SetUp,sf,0,0,0);CHKERRQ(ierr); 25595fce210SBarry Smith sf->setupcalled = PETSC_TRUE; 25695fce210SBarry Smith PetscFunctionReturn(0); 25795fce210SBarry Smith } 25895fce210SBarry Smith 2598af6ec1cSBarry Smith /*@ 26095fce210SBarry Smith PetscSFSetFromOptions - set PetscSF options using the options database 26195fce210SBarry Smith 26295fce210SBarry Smith Logically Collective 26395fce210SBarry Smith 26495fce210SBarry Smith Input Arguments: 26595fce210SBarry Smith . sf - star forest 26695fce210SBarry Smith 26795fce210SBarry Smith Options Database Keys: 26860263706SJed Brown + -sf_type - implementation type, see PetscSFSetType() 26951ccb202SJunchao Zhang . -sf_rank_order - sort composite points for gathers and scatters in rank order, gathers are non-deterministic otherwise 270b85e67b7SJunchao Zhang . -sf_use_default_stream - Assume callers of SF computed the input root/leafdata with the default cuda stream. SF will also 271c2a741eeSJunchao Zhang use the default stream to process data. Therefore, no stream synchronization is needed between SF and its caller (default: true). 272c2a741eeSJunchao Zhang If true, this option only works with -use_cuda_aware_mpi 1. 273c2a741eeSJunchao Zhang - -sf_use_stream_aware_mpi - Assume the underlying MPI is cuda-stream aware and SF won't sync streams for send/recv buffers passed to MPI (default: false). 274c2a741eeSJunchao Zhang If true, this option only works with -use_cuda_aware_mpi 1. 27595fce210SBarry Smith 27695fce210SBarry Smith Level: intermediate 27795fce210SBarry Smith @*/ 27895fce210SBarry Smith PetscErrorCode PetscSFSetFromOptions(PetscSF sf) 27995fce210SBarry Smith { 28095fce210SBarry Smith PetscSFType deft; 28195fce210SBarry Smith char type[256]; 28295fce210SBarry Smith PetscErrorCode ierr; 28395fce210SBarry Smith PetscBool flg; 28495fce210SBarry Smith 28595fce210SBarry Smith PetscFunctionBegin; 28695fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 28795fce210SBarry Smith ierr = PetscObjectOptionsBegin((PetscObject)sf);CHKERRQ(ierr); 28895fce210SBarry Smith deft = ((PetscObject)sf)->type_name ? ((PetscObject)sf)->type_name : PETSCSFBASIC; 28929046d53SLisandro Dalcin ierr = PetscOptionsFList("-sf_type","PetscSF implementation type","PetscSFSetType",PetscSFList,deft,type,sizeof(type),&flg);CHKERRQ(ierr); 29095fce210SBarry Smith ierr = PetscSFSetType(sf,flg ? type : deft);CHKERRQ(ierr); 29195fce210SBarry 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); 292c2a741eeSJunchao Zhang 293c2a741eeSJunchao Zhang #if defined(PETSC_HAVE_CUDA) 294c2a741eeSJunchao Zhang sf->use_default_stream = PETSC_TRUE; /* The assumption is true for PETSc internal use of SF */ 295cd620004SJunchao Zhang ierr = PetscOptionsBool("-sf_use_default_stream","Assume SF's input and output root/leafdata is computed on the default stream","PetscSFSetFromOptions",sf->use_default_stream,&sf->use_default_stream,NULL);CHKERRQ(ierr); 296b85e67b7SJunchao Zhang sf->use_stream_aware_mpi = PETSC_FALSE; 297b85e67b7SJunchao Zhang ierr = PetscOptionsBool("-sf_use_stream_aware_mpi","Assume the underlying MPI is cuda-stream aware","PetscSFSetFromOptions",sf->use_stream_aware_mpi,&sf->use_stream_aware_mpi,NULL);CHKERRQ(ierr); 298c2a741eeSJunchao Zhang #endif 299e55864a3SBarry Smith if (sf->ops->SetFromOptions) {ierr = (*sf->ops->SetFromOptions)(PetscOptionsObject,sf);CHKERRQ(ierr);} 30095fce210SBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 30195fce210SBarry Smith PetscFunctionReturn(0); 30295fce210SBarry Smith } 30395fce210SBarry Smith 30429046d53SLisandro Dalcin /*@ 30595fce210SBarry Smith PetscSFSetRankOrder - sort multi-points for gathers and scatters by rank order 30695fce210SBarry Smith 30795fce210SBarry Smith Logically Collective 30895fce210SBarry Smith 30995fce210SBarry Smith Input Arguments: 31095fce210SBarry Smith + sf - star forest 31195fce210SBarry Smith - flg - PETSC_TRUE to sort, PETSC_FALSE to skip sorting (lower setup cost, but non-deterministic) 31295fce210SBarry Smith 31395fce210SBarry Smith Level: advanced 31495fce210SBarry Smith 31595fce210SBarry Smith .seealso: PetscSFGatherBegin(), PetscSFScatterBegin() 31695fce210SBarry Smith @*/ 31795fce210SBarry Smith PetscErrorCode PetscSFSetRankOrder(PetscSF sf,PetscBool flg) 31895fce210SBarry Smith { 31995fce210SBarry Smith PetscFunctionBegin; 32095fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 32195fce210SBarry Smith PetscValidLogicalCollectiveBool(sf,flg,2); 32295fce210SBarry Smith if (sf->multi) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Rank ordering must be set before first call to PetscSFGatherBegin() or PetscSFScatterBegin()"); 32395fce210SBarry Smith sf->rankorder = flg; 32495fce210SBarry Smith PetscFunctionReturn(0); 32595fce210SBarry Smith } 32695fce210SBarry Smith 3278af6ec1cSBarry Smith /*@ 32895fce210SBarry Smith PetscSFSetGraph - Set a parallel star forest 32995fce210SBarry Smith 33095fce210SBarry Smith Collective 33195fce210SBarry Smith 33295fce210SBarry Smith Input Arguments: 33395fce210SBarry Smith + sf - star forest 33495fce210SBarry Smith . nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves) 33595fce210SBarry Smith . nleaves - number of leaf vertices on the current process, each of these references a root on any process 336c4e6a40aSLawrence Mitchell . ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage (locations must be >= 0, enforced 337c4e6a40aSLawrence Mitchell during setup in debug mode) 33895fce210SBarry Smith . localmode - copy mode for ilocal 339c4e6a40aSLawrence Mitchell . iremote - remote locations of root vertices for each leaf on the current process (locations must be >= 0, enforced 340c4e6a40aSLawrence Mitchell during setup in debug mode) 34195fce210SBarry Smith - remotemode - copy mode for iremote 34295fce210SBarry Smith 34395fce210SBarry Smith Level: intermediate 34495fce210SBarry Smith 34595452b02SPatrick Sanan Notes: 34695452b02SPatrick Sanan In Fortran you must use PETSC_COPY_VALUES for localmode and remotemode 34738ab3f8aSBarry Smith 3482ad1e87fSLisandro Dalcin Developers Note: Local indices which are the identity permutation in the range [0,nleaves) are discarded as they 3492ad1e87fSLisandro Dalcin encode contiguous storage. In such case, if localmode is PETSC_OWN_POINTER, the memory is deallocated as it is not 3502ad1e87fSLisandro Dalcin needed 3512ad1e87fSLisandro Dalcin 352c4e6a40aSLawrence Mitchell Developers Note: This object does not necessarily encode a true star forest in the graph theoretic sense, since leaf 353c4e6a40aSLawrence Mitchell indices are not required to be unique. Some functions, however, rely on unique leaf indices (checked in debug mode). 354c4e6a40aSLawrence Mitchell 35595fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFGetGraph() 35695fce210SBarry Smith @*/ 35795fce210SBarry Smith PetscErrorCode PetscSFSetGraph(PetscSF sf,PetscInt nroots,PetscInt nleaves,const PetscInt *ilocal,PetscCopyMode localmode,const PetscSFNode *iremote,PetscCopyMode remotemode) 35895fce210SBarry Smith { 35995fce210SBarry Smith PetscErrorCode ierr; 36095fce210SBarry Smith 36195fce210SBarry Smith PetscFunctionBegin; 36295fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 36329046d53SLisandro Dalcin if (nleaves > 0 && ilocal) PetscValidIntPointer(ilocal,4); 36429046d53SLisandro Dalcin if (nleaves > 0) PetscValidPointer(iremote,6); 36529046d53SLisandro Dalcin if (nroots < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nroots %D, cannot be negative",nroots); 36695fce210SBarry Smith if (nleaves < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nleaves %D, cannot be negative",nleaves); 36729046d53SLisandro Dalcin 3682a67d2daSStefano Zampini if (sf->nroots >= 0) { /* Reset only if graph already set */ 36995fce210SBarry Smith ierr = PetscSFReset(sf);CHKERRQ(ierr); 3702a67d2daSStefano Zampini } 3712a67d2daSStefano Zampini 37229046d53SLisandro Dalcin ierr = PetscLogEventBegin(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr); 37329046d53SLisandro Dalcin 37495fce210SBarry Smith sf->nroots = nroots; 37595fce210SBarry Smith sf->nleaves = nleaves; 37629046d53SLisandro Dalcin 37729046d53SLisandro Dalcin if (nleaves && ilocal) { 37821c688dcSJed Brown PetscInt i; 37929046d53SLisandro Dalcin PetscInt minleaf = PETSC_MAX_INT; 38029046d53SLisandro Dalcin PetscInt maxleaf = PETSC_MIN_INT; 3812ad1e87fSLisandro Dalcin int contiguous = 1; 38229046d53SLisandro Dalcin for (i=0; i<nleaves; i++) { 38329046d53SLisandro Dalcin minleaf = PetscMin(minleaf,ilocal[i]); 38429046d53SLisandro Dalcin maxleaf = PetscMax(maxleaf,ilocal[i]); 3852ad1e87fSLisandro Dalcin contiguous &= (ilocal[i] == i); 38629046d53SLisandro Dalcin } 38729046d53SLisandro Dalcin sf->minleaf = minleaf; 38829046d53SLisandro Dalcin sf->maxleaf = maxleaf; 3892ad1e87fSLisandro Dalcin if (contiguous) { 3902ad1e87fSLisandro Dalcin if (localmode == PETSC_OWN_POINTER) { 3912ad1e87fSLisandro Dalcin ierr = PetscFree(ilocal);CHKERRQ(ierr); 3922ad1e87fSLisandro Dalcin } 3932ad1e87fSLisandro Dalcin ilocal = NULL; 3942ad1e87fSLisandro Dalcin } 39529046d53SLisandro Dalcin } else { 39629046d53SLisandro Dalcin sf->minleaf = 0; 39729046d53SLisandro Dalcin sf->maxleaf = nleaves - 1; 39829046d53SLisandro Dalcin } 39929046d53SLisandro Dalcin 40029046d53SLisandro Dalcin if (ilocal) { 40195fce210SBarry Smith switch (localmode) { 40295fce210SBarry Smith case PETSC_COPY_VALUES: 403785e854fSJed Brown ierr = PetscMalloc1(nleaves,&sf->mine_alloc);CHKERRQ(ierr); 404580bdb30SBarry Smith ierr = PetscArraycpy(sf->mine_alloc,ilocal,nleaves);CHKERRQ(ierr); 40595fce210SBarry Smith sf->mine = sf->mine_alloc; 40695fce210SBarry Smith break; 40795fce210SBarry Smith case PETSC_OWN_POINTER: 40895fce210SBarry Smith sf->mine_alloc = (PetscInt*)ilocal; 40995fce210SBarry Smith sf->mine = sf->mine_alloc; 41095fce210SBarry Smith break; 41195fce210SBarry Smith case PETSC_USE_POINTER: 41229046d53SLisandro Dalcin sf->mine_alloc = NULL; 41395fce210SBarry Smith sf->mine = (PetscInt*)ilocal; 41495fce210SBarry Smith break; 41595fce210SBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown localmode"); 41695fce210SBarry Smith } 41795fce210SBarry Smith } 41829046d53SLisandro Dalcin 41995fce210SBarry Smith switch (remotemode) { 42095fce210SBarry Smith case PETSC_COPY_VALUES: 421785e854fSJed Brown ierr = PetscMalloc1(nleaves,&sf->remote_alloc);CHKERRQ(ierr); 422580bdb30SBarry Smith ierr = PetscArraycpy(sf->remote_alloc,iremote,nleaves);CHKERRQ(ierr); 42395fce210SBarry Smith sf->remote = sf->remote_alloc; 42495fce210SBarry Smith break; 42595fce210SBarry Smith case PETSC_OWN_POINTER: 42695fce210SBarry Smith sf->remote_alloc = (PetscSFNode*)iremote; 42795fce210SBarry Smith sf->remote = sf->remote_alloc; 42895fce210SBarry Smith break; 42995fce210SBarry Smith case PETSC_USE_POINTER: 43029046d53SLisandro Dalcin sf->remote_alloc = NULL; 43195fce210SBarry Smith sf->remote = (PetscSFNode*)iremote; 43295fce210SBarry Smith break; 43395fce210SBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown remotemode"); 43495fce210SBarry Smith } 43595fce210SBarry Smith 436563ffbabSMatthew G. Knepley ierr = PetscLogEventEnd(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr); 43729046d53SLisandro Dalcin sf->graphset = PETSC_TRUE; 43895fce210SBarry Smith PetscFunctionReturn(0); 43995fce210SBarry Smith } 44095fce210SBarry Smith 44129046d53SLisandro Dalcin /*@ 442dd5b3ca6SJunchao Zhang PetscSFSetGraphWithPattern - Sets the graph of an SF with a specific pattern 443dd5b3ca6SJunchao Zhang 444dd5b3ca6SJunchao Zhang Collective 445dd5b3ca6SJunchao Zhang 446dd5b3ca6SJunchao Zhang Input Parameters: 447dd5b3ca6SJunchao Zhang + sf - The PetscSF 448dd5b3ca6SJunchao Zhang . map - Layout of roots over all processes (insignificant when pattern is PETSCSF_PATTERN_ALLTOALL) 449dd5b3ca6SJunchao Zhang - pattern - One of PETSCSF_PATTERN_ALLGATHER, PETSCSF_PATTERN_GATHER, PETSCSF_PATTERN_ALLTOALL 450dd5b3ca6SJunchao Zhang 451dd5b3ca6SJunchao Zhang Notes: 452dd5b3ca6SJunchao Zhang It is easier to explain PetscSFPattern using vectors. Suppose we have an MPI vector x and its layout is map. 453dd5b3ca6SJunchao Zhang n and N are local and global sizes of x respectively. 454dd5b3ca6SJunchao Zhang 455dd5b3ca6SJunchao Zhang With PETSCSF_PATTERN_ALLGATHER, the routine creates a graph that if one does Bcast on it, it will copy x to 456dd5b3ca6SJunchao Zhang sequential vectors y on all ranks. 457dd5b3ca6SJunchao Zhang 458dd5b3ca6SJunchao Zhang With PETSCSF_PATTERN_GATHER, the routine creates a graph that if one does Bcast on it, it will copy x to a 459dd5b3ca6SJunchao Zhang sequential vector y on rank 0. 460dd5b3ca6SJunchao Zhang 461dd5b3ca6SJunchao Zhang In above cases, entries of x are roots and entries of y are leaves. 462dd5b3ca6SJunchao Zhang 463dd5b3ca6SJunchao Zhang With PETSCSF_PATTERN_ALLTOALL, map is insignificant. Suppose NP is size of sf's communicator. The routine 464dd5b3ca6SJunchao Zhang creates a graph that every rank has NP leaves and NP roots. On rank i, its leaf j is connected to root i 465dd5b3ca6SJunchao Zhang of rank j. Here 0 <=i,j<NP. It is a kind of MPI_Alltoall with sendcount/recvcount being 1. Note that it does 466dd5b3ca6SJunchao Zhang not mean one can not send multiple items. One just needs to create a new MPI datatype for the mulptiple data 467dd5b3ca6SJunchao Zhang items with MPI_Type_contiguous() and use that as the <unit> argument in SF routines. 468dd5b3ca6SJunchao Zhang 469dd5b3ca6SJunchao Zhang In this case, roots and leaves are symmetric. 470dd5b3ca6SJunchao Zhang 471dd5b3ca6SJunchao Zhang Level: intermediate 472dd5b3ca6SJunchao Zhang @*/ 473dd5b3ca6SJunchao Zhang PetscErrorCode PetscSFSetGraphWithPattern(PetscSF sf,PetscLayout map,PetscSFPattern pattern) 474dd5b3ca6SJunchao Zhang { 475dd5b3ca6SJunchao Zhang MPI_Comm comm; 476dd5b3ca6SJunchao Zhang PetscInt n,N,res[2]; 477dd5b3ca6SJunchao Zhang PetscMPIInt rank,size; 478dd5b3ca6SJunchao Zhang PetscSFType type; 479dd5b3ca6SJunchao Zhang PetscErrorCode ierr; 480dd5b3ca6SJunchao Zhang 481dd5b3ca6SJunchao Zhang PetscFunctionBegin; 482dd5b3ca6SJunchao Zhang ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 483dd5b3ca6SJunchao Zhang if (pattern < PETSCSF_PATTERN_ALLGATHER || pattern > PETSCSF_PATTERN_ALLTOALL) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Unsupported PetscSFPattern %D\n",pattern); 484dd5b3ca6SJunchao Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 485dd5b3ca6SJunchao Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 486dd5b3ca6SJunchao Zhang 487dd5b3ca6SJunchao Zhang if (pattern == PETSCSF_PATTERN_ALLTOALL) { 488dd5b3ca6SJunchao Zhang type = PETSCSFALLTOALL; 489dd5b3ca6SJunchao Zhang ierr = PetscLayoutCreate(comm,&sf->map);CHKERRQ(ierr); 490dd5b3ca6SJunchao Zhang ierr = PetscLayoutSetLocalSize(sf->map,size);CHKERRQ(ierr); 491dd5b3ca6SJunchao Zhang ierr = PetscLayoutSetSize(sf->map,((PetscInt)size)*size);CHKERRQ(ierr); 492dd5b3ca6SJunchao Zhang ierr = PetscLayoutSetUp(sf->map);CHKERRQ(ierr); 493dd5b3ca6SJunchao Zhang } else { 494dd5b3ca6SJunchao Zhang ierr = PetscLayoutGetLocalSize(map,&n);CHKERRQ(ierr); 495dd5b3ca6SJunchao Zhang ierr = PetscLayoutGetSize(map,&N);CHKERRQ(ierr); 496dd5b3ca6SJunchao Zhang res[0] = n; 497dd5b3ca6SJunchao Zhang res[1] = -n; 498dd5b3ca6SJunchao Zhang /* Check if n are same over all ranks so that we can optimize it */ 499dd5b3ca6SJunchao Zhang ierr = MPIU_Allreduce(MPI_IN_PLACE,res,2,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr); 500dd5b3ca6SJunchao Zhang if (res[0] == -res[1]) { /* same n */ 501dd5b3ca6SJunchao Zhang type = (pattern == PETSCSF_PATTERN_ALLGATHER) ? PETSCSFALLGATHER : PETSCSFGATHER; 502dd5b3ca6SJunchao Zhang } else { 503dd5b3ca6SJunchao Zhang type = (pattern == PETSCSF_PATTERN_ALLGATHER) ? PETSCSFALLGATHERV : PETSCSFGATHERV; 504dd5b3ca6SJunchao Zhang } 505dd5b3ca6SJunchao Zhang ierr = PetscLayoutReference(map,&sf->map);CHKERRQ(ierr); 506dd5b3ca6SJunchao Zhang } 507dd5b3ca6SJunchao Zhang ierr = PetscSFSetType(sf,type);CHKERRQ(ierr); 508dd5b3ca6SJunchao Zhang 509dd5b3ca6SJunchao Zhang sf->pattern = pattern; 510dd5b3ca6SJunchao Zhang sf->mine = NULL; /* Contiguous */ 511dd5b3ca6SJunchao Zhang 512dd5b3ca6SJunchao Zhang /* Set nleaves, nroots here in case user calls PetscSFGetGraph, which is legal to call even before PetscSFSetUp is called. 513dd5b3ca6SJunchao Zhang Also set other easy stuff. 514dd5b3ca6SJunchao Zhang */ 515dd5b3ca6SJunchao Zhang if (pattern == PETSCSF_PATTERN_ALLGATHER) { 516dd5b3ca6SJunchao Zhang sf->nleaves = N; 517dd5b3ca6SJunchao Zhang sf->nroots = n; 518dd5b3ca6SJunchao Zhang sf->nranks = size; 519dd5b3ca6SJunchao Zhang sf->minleaf = 0; 520dd5b3ca6SJunchao Zhang sf->maxleaf = N - 1; 521dd5b3ca6SJunchao Zhang } else if (pattern == PETSCSF_PATTERN_GATHER) { 522dd5b3ca6SJunchao Zhang sf->nleaves = rank ? 0 : N; 523dd5b3ca6SJunchao Zhang sf->nroots = n; 524dd5b3ca6SJunchao Zhang sf->nranks = rank ? 0 : size; 525dd5b3ca6SJunchao Zhang sf->minleaf = 0; 526dd5b3ca6SJunchao Zhang sf->maxleaf = rank ? -1 : N - 1; 527dd5b3ca6SJunchao Zhang } else if (pattern == PETSCSF_PATTERN_ALLTOALL) { 528dd5b3ca6SJunchao Zhang sf->nleaves = size; 529dd5b3ca6SJunchao Zhang sf->nroots = size; 530dd5b3ca6SJunchao Zhang sf->nranks = size; 531dd5b3ca6SJunchao Zhang sf->minleaf = 0; 532dd5b3ca6SJunchao Zhang sf->maxleaf = size - 1; 533dd5b3ca6SJunchao Zhang } 534dd5b3ca6SJunchao Zhang sf->ndranks = 0; /* We do not need to separate out distinguished ranks for patterned graphs to improve communication performance */ 535dd5b3ca6SJunchao Zhang sf->graphset = PETSC_TRUE; 536dd5b3ca6SJunchao Zhang PetscFunctionReturn(0); 537dd5b3ca6SJunchao Zhang } 538dd5b3ca6SJunchao Zhang 539dd5b3ca6SJunchao Zhang /*@ 54095fce210SBarry Smith PetscSFCreateInverseSF - given a PetscSF in which all vertices have degree 1, creates the inverse map 54195fce210SBarry Smith 54295fce210SBarry Smith Collective 54395fce210SBarry Smith 54495fce210SBarry Smith Input Arguments: 545dd5b3ca6SJunchao Zhang 54695fce210SBarry Smith . sf - star forest to invert 54795fce210SBarry Smith 54895fce210SBarry Smith Output Arguments: 54995fce210SBarry Smith . isf - inverse of sf 55095fce210SBarry Smith Level: advanced 55195fce210SBarry Smith 55295fce210SBarry Smith Notes: 55395fce210SBarry Smith All roots must have degree 1. 55495fce210SBarry Smith 55595fce210SBarry Smith The local space may be a permutation, but cannot be sparse. 55695fce210SBarry Smith 55795fce210SBarry Smith .seealso: PetscSFSetGraph() 55895fce210SBarry Smith @*/ 55995fce210SBarry Smith PetscErrorCode PetscSFCreateInverseSF(PetscSF sf,PetscSF *isf) 56095fce210SBarry Smith { 56195fce210SBarry Smith PetscErrorCode ierr; 56295fce210SBarry Smith PetscMPIInt rank; 56395fce210SBarry Smith PetscInt i,nroots,nleaves,maxlocal,count,*newilocal; 56495fce210SBarry Smith const PetscInt *ilocal; 56595fce210SBarry Smith PetscSFNode *roots,*leaves; 56695fce210SBarry Smith 56795fce210SBarry Smith PetscFunctionBegin; 56829046d53SLisandro Dalcin PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 56929046d53SLisandro Dalcin PetscSFCheckGraphSet(sf,1); 57029046d53SLisandro Dalcin PetscValidPointer(isf,2); 57129046d53SLisandro Dalcin 57295fce210SBarry Smith ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);CHKERRQ(ierr); 57329046d53SLisandro Dalcin maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */ 57429046d53SLisandro Dalcin 57529046d53SLisandro Dalcin ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); 576ae9aee6dSMatthew G. Knepley ierr = PetscMalloc2(nroots,&roots,maxlocal,&leaves);CHKERRQ(ierr); 577ae9aee6dSMatthew G. Knepley for (i=0; i<maxlocal; i++) { 57895fce210SBarry Smith leaves[i].rank = rank; 57995fce210SBarry Smith leaves[i].index = i; 58095fce210SBarry Smith } 58195fce210SBarry Smith for (i=0; i <nroots; i++) { 58295fce210SBarry Smith roots[i].rank = -1; 58395fce210SBarry Smith roots[i].index = -1; 58495fce210SBarry Smith } 5858bfbc91cSJed Brown ierr = PetscSFReduceBegin(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr); 5868bfbc91cSJed Brown ierr = PetscSFReduceEnd(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr); 58795fce210SBarry Smith 58895fce210SBarry Smith /* Check whether our leaves are sparse */ 58995fce210SBarry Smith for (i=0,count=0; i<nroots; i++) if (roots[i].rank >= 0) count++; 59095fce210SBarry Smith if (count == nroots) newilocal = NULL; 59195fce210SBarry Smith else { /* Index for sparse leaves and compact "roots" array (which is to become our leaves). */ 592785e854fSJed Brown ierr = PetscMalloc1(count,&newilocal);CHKERRQ(ierr); 59395fce210SBarry Smith for (i=0,count=0; i<nroots; i++) { 59495fce210SBarry Smith if (roots[i].rank >= 0) { 59595fce210SBarry Smith newilocal[count] = i; 59695fce210SBarry Smith roots[count].rank = roots[i].rank; 59795fce210SBarry Smith roots[count].index = roots[i].index; 59895fce210SBarry Smith count++; 59995fce210SBarry Smith } 60095fce210SBarry Smith } 60195fce210SBarry Smith } 60295fce210SBarry Smith 60395fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,isf);CHKERRQ(ierr); 60495fce210SBarry Smith ierr = PetscSFSetGraph(*isf,maxlocal,count,newilocal,PETSC_OWN_POINTER,roots,PETSC_COPY_VALUES);CHKERRQ(ierr); 60595fce210SBarry Smith ierr = PetscFree2(roots,leaves);CHKERRQ(ierr); 60695fce210SBarry Smith PetscFunctionReturn(0); 60795fce210SBarry Smith } 60895fce210SBarry Smith 60995fce210SBarry Smith /*@ 61095fce210SBarry Smith PetscSFDuplicate - duplicate a PetscSF, optionally preserving rank connectivity and graph 61195fce210SBarry Smith 61295fce210SBarry Smith Collective 61395fce210SBarry Smith 61495fce210SBarry Smith Input Arguments: 61595fce210SBarry Smith + sf - communication object to duplicate 61695fce210SBarry Smith - opt - PETSCSF_DUPLICATE_CONFONLY, PETSCSF_DUPLICATE_RANKS, or PETSCSF_DUPLICATE_GRAPH (see PetscSFDuplicateOption) 61795fce210SBarry Smith 61895fce210SBarry Smith Output Arguments: 61995fce210SBarry Smith . newsf - new communication object 62095fce210SBarry Smith 62195fce210SBarry Smith Level: beginner 62295fce210SBarry Smith 62395fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetType(), PetscSFSetGraph() 62495fce210SBarry Smith @*/ 62595fce210SBarry Smith PetscErrorCode PetscSFDuplicate(PetscSF sf,PetscSFDuplicateOption opt,PetscSF *newsf) 62695fce210SBarry Smith { 62729046d53SLisandro Dalcin PetscSFType type; 62895fce210SBarry Smith PetscErrorCode ierr; 62995fce210SBarry Smith 63095fce210SBarry Smith PetscFunctionBegin; 63129046d53SLisandro Dalcin PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 63229046d53SLisandro Dalcin PetscValidLogicalCollectiveEnum(sf,opt,2); 63329046d53SLisandro Dalcin PetscValidPointer(newsf,3); 63495fce210SBarry Smith ierr = PetscSFCreate(PetscObjectComm((PetscObject)sf),newsf);CHKERRQ(ierr); 63529046d53SLisandro Dalcin ierr = PetscSFGetType(sf,&type);CHKERRQ(ierr); 63629046d53SLisandro Dalcin if (type) {ierr = PetscSFSetType(*newsf,type);CHKERRQ(ierr);} 63795fce210SBarry Smith if (opt == PETSCSF_DUPLICATE_GRAPH) { 638dd5b3ca6SJunchao Zhang PetscSFCheckGraphSet(sf,1); 639dd5b3ca6SJunchao Zhang if (sf->pattern == PETSCSF_PATTERN_GENERAL) { 64095fce210SBarry Smith PetscInt nroots,nleaves; 64195fce210SBarry Smith const PetscInt *ilocal; 64295fce210SBarry Smith const PetscSFNode *iremote; 64395fce210SBarry Smith ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr); 64495fce210SBarry Smith ierr = PetscSFSetGraph(*newsf,nroots,nleaves,ilocal,PETSC_COPY_VALUES,iremote,PETSC_COPY_VALUES);CHKERRQ(ierr); 645dd5b3ca6SJunchao Zhang } else { 646dd5b3ca6SJunchao Zhang ierr = PetscSFSetGraphWithPattern(*newsf,sf->map,sf->pattern);CHKERRQ(ierr); 647dd5b3ca6SJunchao Zhang } 64895fce210SBarry Smith } 64929046d53SLisandro Dalcin if (sf->ops->Duplicate) {ierr = (*sf->ops->Duplicate)(sf,opt,*newsf);CHKERRQ(ierr);} 65095fce210SBarry Smith PetscFunctionReturn(0); 65195fce210SBarry Smith } 65295fce210SBarry Smith 65395fce210SBarry Smith /*@C 65495fce210SBarry Smith PetscSFGetGraph - Get the graph specifying a parallel star forest 65595fce210SBarry Smith 65695fce210SBarry Smith Not Collective 65795fce210SBarry Smith 65895fce210SBarry Smith Input Arguments: 65995fce210SBarry Smith . sf - star forest 66095fce210SBarry Smith 66195fce210SBarry Smith Output Arguments: 66295fce210SBarry Smith + nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves) 66395fce210SBarry Smith . nleaves - number of leaf vertices on the current process, each of these references a root on any process 66495fce210SBarry Smith . ilocal - locations of leaves in leafdata buffers 66595fce210SBarry Smith - iremote - remote locations of root vertices for each leaf on the current process 66695fce210SBarry Smith 667373e0d91SLisandro Dalcin Notes: 668373e0d91SLisandro Dalcin We are not currently requiring that the graph is set, thus returning nroots=-1 if it has not been set yet 669373e0d91SLisandro Dalcin 670245d9833Sprj- When called from Fortran, the returned iremote array is a copy and must be deallocated after use. Consequently, if you 671ca797d7aSLawrence Mitchell want to update the graph, you must call PetscSFSetGraph after modifying the iremote array. 672ca797d7aSLawrence Mitchell 67395fce210SBarry Smith Level: intermediate 67495fce210SBarry Smith 67595fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph() 67695fce210SBarry Smith @*/ 67795fce210SBarry Smith PetscErrorCode PetscSFGetGraph(PetscSF sf,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote) 67895fce210SBarry Smith { 679b8dee149SJunchao Zhang PetscErrorCode ierr; 68095fce210SBarry Smith 68195fce210SBarry Smith PetscFunctionBegin; 68295fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 683b8dee149SJunchao Zhang if (sf->ops->GetGraph) { 684b8dee149SJunchao Zhang ierr = (sf->ops->GetGraph)(sf,nroots,nleaves,ilocal,iremote);CHKERRQ(ierr); 685b8dee149SJunchao Zhang } else { 68695fce210SBarry Smith if (nroots) *nroots = sf->nroots; 68795fce210SBarry Smith if (nleaves) *nleaves = sf->nleaves; 68895fce210SBarry Smith if (ilocal) *ilocal = sf->mine; 68995fce210SBarry Smith if (iremote) *iremote = sf->remote; 690b8dee149SJunchao Zhang } 69195fce210SBarry Smith PetscFunctionReturn(0); 69295fce210SBarry Smith } 69395fce210SBarry Smith 69429046d53SLisandro Dalcin /*@ 69595fce210SBarry Smith PetscSFGetLeafRange - Get the active leaf ranges 69695fce210SBarry Smith 69795fce210SBarry Smith Not Collective 69895fce210SBarry Smith 69995fce210SBarry Smith Input Arguments: 70095fce210SBarry Smith . sf - star forest 70195fce210SBarry Smith 70295fce210SBarry Smith Output Arguments: 703dd5b3ca6SJunchao Zhang + minleaf - minimum active leaf on this process. Return 0 if there are no leaves. 704dd5b3ca6SJunchao Zhang - maxleaf - maximum active leaf on this process. Return -1 if there are no leaves. 70595fce210SBarry Smith 70695fce210SBarry Smith Level: developer 70795fce210SBarry Smith 70895fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph(), PetscSFGetGraph() 70995fce210SBarry Smith @*/ 71095fce210SBarry Smith PetscErrorCode PetscSFGetLeafRange(PetscSF sf,PetscInt *minleaf,PetscInt *maxleaf) 71195fce210SBarry Smith { 71295fce210SBarry Smith 71395fce210SBarry Smith PetscFunctionBegin; 71495fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 71529046d53SLisandro Dalcin PetscSFCheckGraphSet(sf,1); 71695fce210SBarry Smith if (minleaf) *minleaf = sf->minleaf; 71795fce210SBarry Smith if (maxleaf) *maxleaf = sf->maxleaf; 71895fce210SBarry Smith PetscFunctionReturn(0); 71995fce210SBarry Smith } 72095fce210SBarry Smith 72195fce210SBarry Smith /*@C 722fe2efc57SMark PetscSFViewFromOptions - View from Options 723fe2efc57SMark 724fe2efc57SMark Collective on PetscSF 725fe2efc57SMark 726fe2efc57SMark Input Parameters: 727fe2efc57SMark + A - the star forest 728736c3998SJose E. Roman . obj - Optional object 729736c3998SJose E. Roman - name - command line option 730fe2efc57SMark 731fe2efc57SMark Level: intermediate 732fe2efc57SMark .seealso: PetscSF, PetscSFView, PetscObjectViewFromOptions(), PetscSFCreate() 733fe2efc57SMark @*/ 734fe2efc57SMark PetscErrorCode PetscSFViewFromOptions(PetscSF A,PetscObject obj,const char name[]) 735fe2efc57SMark { 736fe2efc57SMark PetscErrorCode ierr; 737fe2efc57SMark 738fe2efc57SMark PetscFunctionBegin; 739fe2efc57SMark PetscValidHeaderSpecific(A,PETSCSF_CLASSID,1); 740fe2efc57SMark ierr = PetscObjectViewFromOptions((PetscObject)A,obj,name);CHKERRQ(ierr); 741fe2efc57SMark PetscFunctionReturn(0); 742fe2efc57SMark } 743fe2efc57SMark 744fe2efc57SMark /*@C 74595fce210SBarry Smith PetscSFView - view a star forest 74695fce210SBarry Smith 74795fce210SBarry Smith Collective 74895fce210SBarry Smith 74995fce210SBarry Smith Input Arguments: 75095fce210SBarry Smith + sf - star forest 75195fce210SBarry Smith - viewer - viewer to display graph, for example PETSC_VIEWER_STDOUT_WORLD 75295fce210SBarry Smith 75395fce210SBarry Smith Level: beginner 75495fce210SBarry Smith 75595fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph() 75695fce210SBarry Smith @*/ 75795fce210SBarry Smith PetscErrorCode PetscSFView(PetscSF sf,PetscViewer viewer) 75895fce210SBarry Smith { 75995fce210SBarry Smith PetscErrorCode ierr; 76095fce210SBarry Smith PetscBool iascii; 76195fce210SBarry Smith PetscViewerFormat format; 76295fce210SBarry Smith 76395fce210SBarry Smith PetscFunctionBegin; 76495fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 76595fce210SBarry Smith if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sf),&viewer);CHKERRQ(ierr);} 76695fce210SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 76795fce210SBarry Smith PetscCheckSameComm(sf,1,viewer,2); 76880153354SVaclav Hapla if (sf->graphset) {ierr = PetscSFSetUp(sf);CHKERRQ(ierr);} 76995fce210SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 77095fce210SBarry Smith if (iascii) { 77195fce210SBarry Smith PetscMPIInt rank; 77281bfa7aaSJed Brown PetscInt ii,i,j; 77395fce210SBarry Smith 774dae58748SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)sf,viewer);CHKERRQ(ierr); 77595fce210SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 77695fce210SBarry Smith if (sf->ops->View) {ierr = (*sf->ops->View)(sf,viewer);CHKERRQ(ierr);} 777dd5b3ca6SJunchao Zhang if (sf->pattern == PETSCSF_PATTERN_GENERAL) { 77880153354SVaclav Hapla if (!sf->graphset) { 77980153354SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"PetscSFSetGraph() has not been called yet\n");CHKERRQ(ierr); 78080153354SVaclav Hapla ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 78180153354SVaclav Hapla PetscFunctionReturn(0); 78280153354SVaclav Hapla } 78395fce210SBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); 7841575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 78595fce210SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of roots=%D, leaves=%D, remote ranks=%D\n",rank,sf->nroots,sf->nleaves,sf->nranks);CHKERRQ(ierr); 78695fce210SBarry Smith for (i=0; i<sf->nleaves; i++) { 78795fce210SBarry 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); 78895fce210SBarry Smith } 78995fce210SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 79095fce210SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 79195fce210SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 79281bfa7aaSJed Brown PetscMPIInt *tmpranks,*perm; 79381bfa7aaSJed Brown ierr = PetscMalloc2(sf->nranks,&tmpranks,sf->nranks,&perm);CHKERRQ(ierr); 794580bdb30SBarry Smith ierr = PetscArraycpy(tmpranks,sf->ranks,sf->nranks);CHKERRQ(ierr); 79581bfa7aaSJed Brown for (i=0; i<sf->nranks; i++) perm[i] = i; 79681bfa7aaSJed Brown ierr = PetscSortMPIIntWithArray(sf->nranks,tmpranks,perm);CHKERRQ(ierr); 79795fce210SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Roots referenced by my leaves, by rank\n",rank);CHKERRQ(ierr); 79881bfa7aaSJed Brown for (ii=0; ii<sf->nranks; ii++) { 79981bfa7aaSJed Brown i = perm[ii]; 8007904a332SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %d: %D edges\n",rank,sf->ranks[i],sf->roffset[i+1]-sf->roffset[i]);CHKERRQ(ierr); 80195fce210SBarry Smith for (j=sf->roffset[i]; j<sf->roffset[i+1]; j++) { 80295fce210SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D <- %D\n",rank,sf->rmine[j],sf->rremote[j]);CHKERRQ(ierr); 80395fce210SBarry Smith } 80495fce210SBarry Smith } 80581bfa7aaSJed Brown ierr = PetscFree2(tmpranks,perm);CHKERRQ(ierr); 80695fce210SBarry Smith } 80795fce210SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 8081575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 809dd5b3ca6SJunchao Zhang } 81095fce210SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 81195fce210SBarry Smith } 81295fce210SBarry Smith PetscFunctionReturn(0); 81395fce210SBarry Smith } 81495fce210SBarry Smith 81595fce210SBarry Smith /*@C 816dec1416fSJunchao Zhang PetscSFGetRootRanks - Get root ranks and number of vertices referenced by leaves on this process 81795fce210SBarry Smith 81895fce210SBarry Smith Not Collective 81995fce210SBarry Smith 82095fce210SBarry Smith Input Arguments: 82195fce210SBarry Smith . sf - star forest 82295fce210SBarry Smith 82395fce210SBarry Smith Output Arguments: 82495fce210SBarry Smith + nranks - number of ranks referenced by local part 82595fce210SBarry Smith . ranks - array of ranks 82695fce210SBarry Smith . roffset - offset in rmine/rremote for each rank (length nranks+1) 82795fce210SBarry Smith . rmine - concatenated array holding local indices referencing each remote rank 82895fce210SBarry Smith - rremote - concatenated array holding remote indices referenced for each remote rank 82995fce210SBarry Smith 83095fce210SBarry Smith Level: developer 83195fce210SBarry Smith 832dec1416fSJunchao Zhang .seealso: PetscSFGetLeafRanks() 83395fce210SBarry Smith @*/ 834dec1416fSJunchao Zhang PetscErrorCode PetscSFGetRootRanks(PetscSF sf,PetscInt *nranks,const PetscMPIInt **ranks,const PetscInt **roffset,const PetscInt **rmine,const PetscInt **rremote) 83595fce210SBarry Smith { 836dec1416fSJunchao Zhang PetscErrorCode ierr; 83795fce210SBarry Smith 83895fce210SBarry Smith PetscFunctionBegin; 83995fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 84029046d53SLisandro Dalcin if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining ranks"); 841dec1416fSJunchao Zhang if (sf->ops->GetRootRanks) { 842dec1416fSJunchao Zhang ierr = (sf->ops->GetRootRanks)(sf,nranks,ranks,roffset,rmine,rremote);CHKERRQ(ierr); 843dec1416fSJunchao Zhang } else { 844dec1416fSJunchao Zhang /* The generic implementation */ 84595fce210SBarry Smith if (nranks) *nranks = sf->nranks; 84695fce210SBarry Smith if (ranks) *ranks = sf->ranks; 84795fce210SBarry Smith if (roffset) *roffset = sf->roffset; 84895fce210SBarry Smith if (rmine) *rmine = sf->rmine; 84995fce210SBarry Smith if (rremote) *rremote = sf->rremote; 850dec1416fSJunchao Zhang } 85195fce210SBarry Smith PetscFunctionReturn(0); 85295fce210SBarry Smith } 85395fce210SBarry Smith 8548750ddebSJunchao Zhang /*@C 8558750ddebSJunchao Zhang PetscSFGetLeafRanks - Get leaf ranks referencing roots on this process 8568750ddebSJunchao Zhang 8578750ddebSJunchao Zhang Not Collective 8588750ddebSJunchao Zhang 8598750ddebSJunchao Zhang Input Arguments: 8608750ddebSJunchao Zhang . sf - star forest 8618750ddebSJunchao Zhang 8628750ddebSJunchao Zhang Output Arguments: 8638750ddebSJunchao Zhang + niranks - number of leaf ranks referencing roots on this process 8648750ddebSJunchao Zhang . iranks - array of ranks 8658750ddebSJunchao Zhang . ioffset - offset in irootloc for each rank (length niranks+1) 8668750ddebSJunchao Zhang - irootloc - concatenated array holding local indices of roots referenced by each leaf rank 8678750ddebSJunchao Zhang 8688750ddebSJunchao Zhang Level: developer 8698750ddebSJunchao Zhang 870dec1416fSJunchao Zhang .seealso: PetscSFGetRootRanks() 8718750ddebSJunchao Zhang @*/ 8728750ddebSJunchao Zhang PetscErrorCode PetscSFGetLeafRanks(PetscSF sf,PetscInt *niranks,const PetscMPIInt **iranks,const PetscInt **ioffset,const PetscInt **irootloc) 8738750ddebSJunchao Zhang { 8748750ddebSJunchao Zhang PetscErrorCode ierr; 8758750ddebSJunchao Zhang 8768750ddebSJunchao Zhang PetscFunctionBegin; 8778750ddebSJunchao Zhang PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 8788750ddebSJunchao Zhang if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining ranks"); 8798750ddebSJunchao Zhang if (sf->ops->GetLeafRanks) { 8808750ddebSJunchao Zhang ierr = (sf->ops->GetLeafRanks)(sf,niranks,iranks,ioffset,irootloc);CHKERRQ(ierr); 8818750ddebSJunchao Zhang } else { 8828750ddebSJunchao Zhang PetscSFType type; 8838750ddebSJunchao Zhang ierr = PetscSFGetType(sf,&type);CHKERRQ(ierr); 8848750ddebSJunchao Zhang SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSFGetLeafRanks() is not supported on this StarForest type: %s", type); 8858750ddebSJunchao Zhang } 8868750ddebSJunchao Zhang PetscFunctionReturn(0); 8878750ddebSJunchao Zhang } 8888750ddebSJunchao Zhang 889b5a8e515SJed Brown static PetscBool InList(PetscMPIInt needle,PetscMPIInt n,const PetscMPIInt *list) { 890b5a8e515SJed Brown PetscInt i; 891b5a8e515SJed Brown for (i=0; i<n; i++) { 892b5a8e515SJed Brown if (needle == list[i]) return PETSC_TRUE; 893b5a8e515SJed Brown } 894b5a8e515SJed Brown return PETSC_FALSE; 895b5a8e515SJed Brown } 896b5a8e515SJed Brown 89795fce210SBarry Smith /*@C 89821c688dcSJed Brown PetscSFSetUpRanks - Set up data structures associated with ranks; this is for internal use by PetscSF implementations. 89921c688dcSJed Brown 90021c688dcSJed Brown Collective 90121c688dcSJed Brown 90221c688dcSJed Brown Input Arguments: 903b5a8e515SJed Brown + sf - PetscSF to set up; PetscSFSetGraph() must have been called 904b5a8e515SJed Brown - dgroup - MPI_Group of ranks to be distinguished (e.g., for self or shared memory exchange) 90521c688dcSJed Brown 90621c688dcSJed Brown Level: developer 90721c688dcSJed Brown 908dec1416fSJunchao Zhang .seealso: PetscSFGetRootRanks() 90921c688dcSJed Brown @*/ 910b5a8e515SJed Brown PetscErrorCode PetscSFSetUpRanks(PetscSF sf,MPI_Group dgroup) 91121c688dcSJed Brown { 91221c688dcSJed Brown PetscErrorCode ierr; 91321c688dcSJed Brown PetscTable table; 91421c688dcSJed Brown PetscTablePosition pos; 915b5a8e515SJed Brown PetscMPIInt size,groupsize,*groupranks; 916247e8311SStefano Zampini PetscInt *rcount,*ranks; 917247e8311SStefano Zampini PetscInt i, irank = -1,orank = -1; 91821c688dcSJed Brown 91921c688dcSJed Brown PetscFunctionBegin; 92021c688dcSJed Brown PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 92129046d53SLisandro Dalcin PetscSFCheckGraphSet(sf,1); 92221c688dcSJed Brown ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr); 92321c688dcSJed Brown ierr = PetscTableCreate(10,size,&table);CHKERRQ(ierr); 92421c688dcSJed Brown for (i=0; i<sf->nleaves; i++) { 92521c688dcSJed Brown /* Log 1-based rank */ 92621c688dcSJed Brown ierr = PetscTableAdd(table,sf->remote[i].rank+1,1,ADD_VALUES);CHKERRQ(ierr); 92721c688dcSJed Brown } 92821c688dcSJed Brown ierr = PetscTableGetCount(table,&sf->nranks);CHKERRQ(ierr); 92921c688dcSJed Brown ierr = PetscMalloc4(sf->nranks,&sf->ranks,sf->nranks+1,&sf->roffset,sf->nleaves,&sf->rmine,sf->nleaves,&sf->rremote);CHKERRQ(ierr); 93021c688dcSJed Brown ierr = PetscMalloc2(sf->nranks,&rcount,sf->nranks,&ranks);CHKERRQ(ierr); 93121c688dcSJed Brown ierr = PetscTableGetHeadPosition(table,&pos);CHKERRQ(ierr); 93221c688dcSJed Brown for (i=0; i<sf->nranks; i++) { 93321c688dcSJed Brown ierr = PetscTableGetNext(table,&pos,&ranks[i],&rcount[i]);CHKERRQ(ierr); 93421c688dcSJed Brown ranks[i]--; /* Convert back to 0-based */ 93521c688dcSJed Brown } 93621c688dcSJed Brown ierr = PetscTableDestroy(&table);CHKERRQ(ierr); 937b5a8e515SJed Brown 938b5a8e515SJed Brown /* We expect that dgroup is reliably "small" while nranks could be large */ 939b5a8e515SJed Brown { 9407fb8a5e4SKarl Rupp MPI_Group group = MPI_GROUP_NULL; 941b5a8e515SJed Brown PetscMPIInt *dgroupranks; 942b5a8e515SJed Brown ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr); 943b5a8e515SJed Brown ierr = MPI_Group_size(dgroup,&groupsize);CHKERRQ(ierr); 944b5a8e515SJed Brown ierr = PetscMalloc1(groupsize,&dgroupranks);CHKERRQ(ierr); 945b5a8e515SJed Brown ierr = PetscMalloc1(groupsize,&groupranks);CHKERRQ(ierr); 946b5a8e515SJed Brown for (i=0; i<groupsize; i++) dgroupranks[i] = i; 947f3fc9a17SJed Brown if (groupsize) {ierr = MPI_Group_translate_ranks(dgroup,groupsize,dgroupranks,group,groupranks);CHKERRQ(ierr);} 948b5a8e515SJed Brown ierr = MPI_Group_free(&group);CHKERRQ(ierr); 949b5a8e515SJed Brown ierr = PetscFree(dgroupranks);CHKERRQ(ierr); 950b5a8e515SJed Brown } 951b5a8e515SJed Brown 952b5a8e515SJed Brown /* Partition ranks[] into distinguished (first sf->ndranks) followed by non-distinguished */ 953b5a8e515SJed Brown for (sf->ndranks=0,i=sf->nranks; sf->ndranks<i; ) { 954b5a8e515SJed Brown for (i--; sf->ndranks<i; i--) { /* Scan i backward looking for distinguished rank */ 955b5a8e515SJed Brown if (InList(ranks[i],groupsize,groupranks)) break; 956b5a8e515SJed Brown } 957b5a8e515SJed Brown for ( ; sf->ndranks<=i; sf->ndranks++) { /* Scan sf->ndranks forward looking for non-distinguished rank */ 958b5a8e515SJed Brown if (!InList(ranks[sf->ndranks],groupsize,groupranks)) break; 959b5a8e515SJed Brown } 960b5a8e515SJed Brown if (sf->ndranks < i) { /* Swap ranks[sf->ndranks] with ranks[i] */ 961b5a8e515SJed Brown PetscInt tmprank,tmpcount; 962247e8311SStefano Zampini 963b5a8e515SJed Brown tmprank = ranks[i]; 964b5a8e515SJed Brown tmpcount = rcount[i]; 965b5a8e515SJed Brown ranks[i] = ranks[sf->ndranks]; 966b5a8e515SJed Brown rcount[i] = rcount[sf->ndranks]; 967b5a8e515SJed Brown ranks[sf->ndranks] = tmprank; 968b5a8e515SJed Brown rcount[sf->ndranks] = tmpcount; 969b5a8e515SJed Brown sf->ndranks++; 970b5a8e515SJed Brown } 971b5a8e515SJed Brown } 972b5a8e515SJed Brown ierr = PetscFree(groupranks);CHKERRQ(ierr); 973b5a8e515SJed Brown ierr = PetscSortIntWithArray(sf->ndranks,ranks,rcount);CHKERRQ(ierr); 974b5a8e515SJed Brown ierr = PetscSortIntWithArray(sf->nranks-sf->ndranks,ranks+sf->ndranks,rcount+sf->ndranks);CHKERRQ(ierr); 97521c688dcSJed Brown sf->roffset[0] = 0; 97621c688dcSJed Brown for (i=0; i<sf->nranks; i++) { 97721c688dcSJed Brown ierr = PetscMPIIntCast(ranks[i],sf->ranks+i);CHKERRQ(ierr); 97821c688dcSJed Brown sf->roffset[i+1] = sf->roffset[i] + rcount[i]; 97921c688dcSJed Brown rcount[i] = 0; 98021c688dcSJed Brown } 981247e8311SStefano Zampini for (i=0, irank = -1, orank = -1; i<sf->nleaves; i++) { 982247e8311SStefano Zampini /* short circuit */ 983247e8311SStefano Zampini if (orank != sf->remote[i].rank) { 98421c688dcSJed Brown /* Search for index of iremote[i].rank in sf->ranks */ 985b5a8e515SJed Brown ierr = PetscFindMPIInt(sf->remote[i].rank,sf->ndranks,sf->ranks,&irank);CHKERRQ(ierr); 986b5a8e515SJed Brown if (irank < 0) { 987b5a8e515SJed Brown ierr = PetscFindMPIInt(sf->remote[i].rank,sf->nranks-sf->ndranks,sf->ranks+sf->ndranks,&irank);CHKERRQ(ierr); 988b5a8e515SJed Brown if (irank >= 0) irank += sf->ndranks; 98921c688dcSJed Brown } 990247e8311SStefano Zampini orank = sf->remote[i].rank; 991247e8311SStefano Zampini } 992b5a8e515SJed Brown if (irank < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Could not find rank %D in array",sf->remote[i].rank); 99321c688dcSJed Brown sf->rmine[sf->roffset[irank] + rcount[irank]] = sf->mine ? sf->mine[i] : i; 99421c688dcSJed Brown sf->rremote[sf->roffset[irank] + rcount[irank]] = sf->remote[i].index; 99521c688dcSJed Brown rcount[irank]++; 99621c688dcSJed Brown } 99721c688dcSJed Brown ierr = PetscFree2(rcount,ranks);CHKERRQ(ierr); 99821c688dcSJed Brown PetscFunctionReturn(0); 99921c688dcSJed Brown } 100021c688dcSJed Brown 100121c688dcSJed Brown /*@C 100295fce210SBarry Smith PetscSFGetGroups - gets incoming and outgoing process groups 100395fce210SBarry Smith 100495fce210SBarry Smith Collective 100595fce210SBarry Smith 100695fce210SBarry Smith Input Argument: 100795fce210SBarry Smith . sf - star forest 100895fce210SBarry Smith 100995fce210SBarry Smith Output Arguments: 101095fce210SBarry Smith + incoming - group of origin processes for incoming edges (leaves that reference my roots) 101195fce210SBarry Smith - outgoing - group of destination processes for outgoing edges (roots that I reference) 101295fce210SBarry Smith 101395fce210SBarry Smith Level: developer 101495fce210SBarry Smith 101595fce210SBarry Smith .seealso: PetscSFGetWindow(), PetscSFRestoreWindow() 101695fce210SBarry Smith @*/ 101795fce210SBarry Smith PetscErrorCode PetscSFGetGroups(PetscSF sf,MPI_Group *incoming,MPI_Group *outgoing) 101895fce210SBarry Smith { 101995fce210SBarry Smith PetscErrorCode ierr; 10207fb8a5e4SKarl Rupp MPI_Group group = MPI_GROUP_NULL; 102195fce210SBarry Smith 102295fce210SBarry Smith PetscFunctionBegin; 102344ee17edSStefano Zampini if (sf->nranks < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUpRanks() before obtaining groups"); 102495fce210SBarry Smith if (sf->ingroup == MPI_GROUP_NULL) { 102595fce210SBarry Smith PetscInt i; 102695fce210SBarry Smith const PetscInt *indegree; 102795fce210SBarry Smith PetscMPIInt rank,*outranks,*inranks; 102895fce210SBarry Smith PetscSFNode *remote; 102995fce210SBarry Smith PetscSF bgcount; 103095fce210SBarry Smith 103195fce210SBarry Smith /* Compute the number of incoming ranks */ 1032785e854fSJed Brown ierr = PetscMalloc1(sf->nranks,&remote);CHKERRQ(ierr); 103395fce210SBarry Smith for (i=0; i<sf->nranks; i++) { 103495fce210SBarry Smith remote[i].rank = sf->ranks[i]; 103595fce210SBarry Smith remote[i].index = 0; 103695fce210SBarry Smith } 103795fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,&bgcount);CHKERRQ(ierr); 103895fce210SBarry Smith ierr = PetscSFSetGraph(bgcount,1,sf->nranks,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr); 103995fce210SBarry Smith ierr = PetscSFComputeDegreeBegin(bgcount,&indegree);CHKERRQ(ierr); 104095fce210SBarry Smith ierr = PetscSFComputeDegreeEnd(bgcount,&indegree);CHKERRQ(ierr); 104195fce210SBarry Smith /* Enumerate the incoming ranks */ 1042dcca6d9dSJed Brown ierr = PetscMalloc2(indegree[0],&inranks,sf->nranks,&outranks);CHKERRQ(ierr); 104395fce210SBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); 104495fce210SBarry Smith for (i=0; i<sf->nranks; i++) outranks[i] = rank; 104595fce210SBarry Smith ierr = PetscSFGatherBegin(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr); 104695fce210SBarry Smith ierr = PetscSFGatherEnd(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr); 104795fce210SBarry Smith ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr); 104895fce210SBarry Smith ierr = MPI_Group_incl(group,indegree[0],inranks,&sf->ingroup);CHKERRQ(ierr); 104995fce210SBarry Smith ierr = MPI_Group_free(&group);CHKERRQ(ierr); 105095fce210SBarry Smith ierr = PetscFree2(inranks,outranks);CHKERRQ(ierr); 105195fce210SBarry Smith ierr = PetscSFDestroy(&bgcount);CHKERRQ(ierr); 105295fce210SBarry Smith } 105395fce210SBarry Smith *incoming = sf->ingroup; 105495fce210SBarry Smith 105595fce210SBarry Smith if (sf->outgroup == MPI_GROUP_NULL) { 105695fce210SBarry Smith ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr); 105795fce210SBarry Smith ierr = MPI_Group_incl(group,sf->nranks,sf->ranks,&sf->outgroup);CHKERRQ(ierr); 105895fce210SBarry Smith ierr = MPI_Group_free(&group);CHKERRQ(ierr); 105995fce210SBarry Smith } 106095fce210SBarry Smith *outgoing = sf->outgroup; 106195fce210SBarry Smith PetscFunctionReturn(0); 106295fce210SBarry Smith } 106395fce210SBarry Smith 106429046d53SLisandro Dalcin /*@ 106595fce210SBarry Smith PetscSFGetMultiSF - gets the inner SF implemeting gathers and scatters 106695fce210SBarry Smith 106795fce210SBarry Smith Collective 106895fce210SBarry Smith 106995fce210SBarry Smith Input Argument: 107095fce210SBarry Smith . sf - star forest that may contain roots with 0 or with more than 1 vertex 107195fce210SBarry Smith 107295fce210SBarry Smith Output Arguments: 107395fce210SBarry Smith . multi - star forest with split roots, such that each root has degree exactly 1 107495fce210SBarry Smith 107595fce210SBarry Smith Level: developer 107695fce210SBarry Smith 107795fce210SBarry Smith Notes: 107895fce210SBarry Smith 107995fce210SBarry Smith In most cases, users should use PetscSFGatherBegin() and PetscSFScatterBegin() instead of manipulating multi 108095fce210SBarry Smith directly. Since multi satisfies the stronger condition that each entry in the global space has exactly one incoming 108195fce210SBarry Smith edge, it is a candidate for future optimization that might involve its removal. 108295fce210SBarry Smith 1083673100f5SVaclav Hapla .seealso: PetscSFSetGraph(), PetscSFGatherBegin(), PetscSFScatterBegin(), PetscSFComputeMultiRootOriginalNumbering() 108495fce210SBarry Smith @*/ 108595fce210SBarry Smith PetscErrorCode PetscSFGetMultiSF(PetscSF sf,PetscSF *multi) 108695fce210SBarry Smith { 108795fce210SBarry Smith PetscErrorCode ierr; 108895fce210SBarry Smith 108995fce210SBarry Smith PetscFunctionBegin; 109095fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 109195fce210SBarry Smith PetscValidPointer(multi,2); 109295fce210SBarry Smith if (sf->nroots < 0) { /* Graph has not been set yet; why do we need this? */ 109395fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr); 109495fce210SBarry Smith *multi = sf->multi; 109595fce210SBarry Smith PetscFunctionReturn(0); 109695fce210SBarry Smith } 109795fce210SBarry Smith if (!sf->multi) { 109895fce210SBarry Smith const PetscInt *indegree; 10999837ea96SMatthew G. Knepley PetscInt i,*inoffset,*outones,*outoffset,maxlocal; 110095fce210SBarry Smith PetscSFNode *remote; 110129046d53SLisandro Dalcin maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */ 110295fce210SBarry Smith ierr = PetscSFComputeDegreeBegin(sf,&indegree);CHKERRQ(ierr); 110395fce210SBarry Smith ierr = PetscSFComputeDegreeEnd(sf,&indegree);CHKERRQ(ierr); 11049837ea96SMatthew G. Knepley ierr = PetscMalloc3(sf->nroots+1,&inoffset,maxlocal,&outones,maxlocal,&outoffset);CHKERRQ(ierr); 110595fce210SBarry Smith inoffset[0] = 0; 110695fce210SBarry Smith for (i=0; i<sf->nroots; i++) inoffset[i+1] = inoffset[i] + indegree[i]; 11079837ea96SMatthew G. Knepley for (i=0; i<maxlocal; i++) outones[i] = 1; 1108dbd2ff41SBarry Smith ierr = PetscSFFetchAndOpBegin(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr); 1109dbd2ff41SBarry Smith ierr = PetscSFFetchAndOpEnd(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr); 111095fce210SBarry Smith for (i=0; i<sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */ 1111*76bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { /* Check that the expected number of increments occurred */ 111295fce210SBarry Smith for (i=0; i<sf->nroots; i++) { 111395fce210SBarry Smith if (inoffset[i] + indegree[i] != inoffset[i+1]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect result after PetscSFFetchAndOp"); 111495fce210SBarry Smith } 1115*76bd3646SJed Brown } 1116785e854fSJed Brown ierr = PetscMalloc1(sf->nleaves,&remote);CHKERRQ(ierr); 111795fce210SBarry Smith for (i=0; i<sf->nleaves; i++) { 111895fce210SBarry Smith remote[i].rank = sf->remote[i].rank; 111938e7336fSToby Isaac remote[i].index = outoffset[sf->mine ? sf->mine[i] : i]; 112095fce210SBarry Smith } 112195fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr); 112201365b40SToby Isaac ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr); 112395fce210SBarry Smith if (sf->rankorder) { /* Sort the ranks */ 112495fce210SBarry Smith PetscMPIInt rank; 112595fce210SBarry Smith PetscInt *inranks,*newoffset,*outranks,*newoutoffset,*tmpoffset,maxdegree; 112695fce210SBarry Smith PetscSFNode *newremote; 112795fce210SBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); 112895fce210SBarry Smith for (i=0,maxdegree=0; i<sf->nroots; i++) maxdegree = PetscMax(maxdegree,indegree[i]); 11299837ea96SMatthew G. Knepley ierr = PetscMalloc5(sf->multi->nroots,&inranks,sf->multi->nroots,&newoffset,maxlocal,&outranks,maxlocal,&newoutoffset,maxdegree,&tmpoffset);CHKERRQ(ierr); 11309837ea96SMatthew G. Knepley for (i=0; i<maxlocal; i++) outranks[i] = rank; 11318bfbc91cSJed Brown ierr = PetscSFReduceBegin(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr); 11328bfbc91cSJed Brown ierr = PetscSFReduceEnd(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr); 113395fce210SBarry Smith /* Sort the incoming ranks at each vertex, build the inverse map */ 113495fce210SBarry Smith for (i=0; i<sf->nroots; i++) { 113595fce210SBarry Smith PetscInt j; 113695fce210SBarry Smith for (j=0; j<indegree[i]; j++) tmpoffset[j] = j; 113795fce210SBarry Smith ierr = PetscSortIntWithArray(indegree[i],inranks+inoffset[i],tmpoffset);CHKERRQ(ierr); 113895fce210SBarry Smith for (j=0; j<indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j; 113995fce210SBarry Smith } 114095fce210SBarry Smith ierr = PetscSFBcastBegin(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr); 114195fce210SBarry Smith ierr = PetscSFBcastEnd(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr); 1142785e854fSJed Brown ierr = PetscMalloc1(sf->nleaves,&newremote);CHKERRQ(ierr); 114395fce210SBarry Smith for (i=0; i<sf->nleaves; i++) { 114495fce210SBarry Smith newremote[i].rank = sf->remote[i].rank; 114501365b40SToby Isaac newremote[i].index = newoutoffset[sf->mine ? sf->mine[i] : i]; 114695fce210SBarry Smith } 114701365b40SToby Isaac ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,newremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 114895fce210SBarry Smith ierr = PetscFree5(inranks,newoffset,outranks,newoutoffset,tmpoffset);CHKERRQ(ierr); 114995fce210SBarry Smith } 115095fce210SBarry Smith ierr = PetscFree3(inoffset,outones,outoffset);CHKERRQ(ierr); 115195fce210SBarry Smith } 115295fce210SBarry Smith *multi = sf->multi; 115395fce210SBarry Smith PetscFunctionReturn(0); 115495fce210SBarry Smith } 115595fce210SBarry Smith 115695fce210SBarry Smith /*@C 115795fce210SBarry Smith PetscSFCreateEmbeddedSF - removes edges from all but the selected roots, does not remap indices 115895fce210SBarry Smith 115995fce210SBarry Smith Collective 116095fce210SBarry Smith 116195fce210SBarry Smith Input Arguments: 116295fce210SBarry Smith + sf - original star forest 1163ba2a7774SJunchao Zhang . nselected - number of selected roots on this process 1164ba2a7774SJunchao Zhang - selected - indices of the selected roots on this process 116595fce210SBarry Smith 116695fce210SBarry Smith Output Arguments: 1167cd620004SJunchao Zhang . esf - new star forest 116895fce210SBarry Smith 116995fce210SBarry Smith Level: advanced 117095fce210SBarry Smith 117195fce210SBarry Smith Note: 117295fce210SBarry Smith To use the new PetscSF, it may be necessary to know the indices of the leaves that are still participating. This can 117395fce210SBarry Smith be done by calling PetscSFGetGraph(). 117495fce210SBarry Smith 117595fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFGetGraph() 117695fce210SBarry Smith @*/ 1177cd620004SJunchao Zhang PetscErrorCode PetscSFCreateEmbeddedSF(PetscSF sf,PetscInt nselected,const PetscInt *selected,PetscSF *esf) 117895fce210SBarry Smith { 1179cd620004SJunchao Zhang PetscInt i,j,n,nroots,nleaves,esf_nleaves,*new_ilocal,minleaf,maxleaf,maxlocal; 1180cd620004SJunchao Zhang const PetscInt *ilocal; 1181cd620004SJunchao Zhang signed char *rootdata,*leafdata,*leafmem; 1182ba2a7774SJunchao Zhang const PetscSFNode *iremote; 1183f659e5c7SJunchao Zhang PetscSFNode *new_iremote; 1184f659e5c7SJunchao Zhang MPI_Comm comm; 11850511a646SMatthew G. Knepley PetscErrorCode ierr; 118695fce210SBarry Smith 118795fce210SBarry Smith PetscFunctionBegin; 118895fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 118929046d53SLisandro Dalcin PetscSFCheckGraphSet(sf,1); 1190ba2a7774SJunchao Zhang if (nselected) PetscValidPointer(selected,3); 1191cd620004SJunchao Zhang PetscValidPointer(esf,4); 11920511a646SMatthew G. Knepley 1193f659e5c7SJunchao Zhang ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 1194cd620004SJunchao Zhang 1195140a1472SStefano Zampini ierr = PetscLogEventBegin(PETSCSF_EmbedSF,sf,0,0,0);CHKERRQ(ierr); 1196f659e5c7SJunchao Zhang ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr); 1197cd620004SJunchao Zhang ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr); 1198cd620004SJunchao Zhang 1199*76bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { /* Error out if selected[] has dups or out of range indices */ 1200*76bd3646SJed Brown 1201cd620004SJunchao Zhang PetscBool dups; 1202cd620004SJunchao Zhang ierr = PetscCheckDupsInt(nselected,selected,&dups);CHKERRQ(ierr); 1203cd620004SJunchao Zhang if (dups) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"selected[] has dups"); 1204cd620004SJunchao Zhang for (i=0; i<nselected; i++) 1205cd620004SJunchao Zhang if (selected[i] < 0 || selected[i] >= nroots) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"selected root indice %D is out of [0,%D)",selected[i],nroots); 1206cd620004SJunchao Zhang } 1207f659e5c7SJunchao Zhang 1208f659e5c7SJunchao Zhang if (sf->ops->CreateEmbeddedSF) { 1209cd620004SJunchao Zhang ierr = (*sf->ops->CreateEmbeddedSF)(sf,nselected,selected,esf);CHKERRQ(ierr); 1210f659e5c7SJunchao Zhang } else { 1211cd620004SJunchao Zhang /* A generic version of creating embedded sf */ 1212cd620004SJunchao Zhang ierr = PetscSFGetLeafRange(sf,&minleaf,&maxleaf);CHKERRQ(ierr); 1213cd620004SJunchao Zhang maxlocal = maxleaf - minleaf + 1; 1214cd620004SJunchao Zhang ierr = PetscCalloc2(nroots,&rootdata,maxlocal,&leafmem);CHKERRQ(ierr); 1215cd620004SJunchao Zhang leafdata = leafmem - minleaf; 1216cd620004SJunchao Zhang /* Tag selected roots and bcast to leaves */ 1217cd620004SJunchao Zhang for (i=0; i<nselected; i++) rootdata[selected[i]] = 1; 1218cd620004SJunchao Zhang ierr = PetscSFBcastBegin(sf,MPI_SIGNED_CHAR,rootdata,leafdata);CHKERRQ(ierr); 1219cd620004SJunchao Zhang ierr = PetscSFBcastEnd(sf,MPI_SIGNED_CHAR,rootdata,leafdata);CHKERRQ(ierr); 1220ba2a7774SJunchao Zhang 1221cd620004SJunchao Zhang /* Build esf with leaves that are still connected */ 1222cd620004SJunchao Zhang esf_nleaves = 0; 1223cd620004SJunchao Zhang for (i=0; i<nleaves; i++) { 1224cd620004SJunchao Zhang j = ilocal ? ilocal[i] : i; 1225cd620004SJunchao Zhang /* esf_nleaves += leafdata[j] should work in theory, but failed with SFWindow bugs 1226cd620004SJunchao Zhang with PetscSFBcast. See https://gitlab.com/petsc/petsc/issues/555 1227cd620004SJunchao Zhang */ 1228cd620004SJunchao Zhang esf_nleaves += (leafdata[j] ? 1 : 0); 1229cd620004SJunchao Zhang } 1230cd620004SJunchao Zhang ierr = PetscMalloc1(esf_nleaves,&new_ilocal);CHKERRQ(ierr); 1231cd620004SJunchao Zhang ierr = PetscMalloc1(esf_nleaves,&new_iremote);CHKERRQ(ierr); 1232cd620004SJunchao Zhang for (i=n=0; i<nleaves; i++) { 1233cd620004SJunchao Zhang j = ilocal ? ilocal[i] : i; 1234cd620004SJunchao Zhang if (leafdata[j]) { 1235cd620004SJunchao Zhang new_ilocal[n] = j; 1236cd620004SJunchao Zhang new_iremote[n].rank = iremote[i].rank; 1237cd620004SJunchao Zhang new_iremote[n].index = iremote[i].index; 1238fc1ede2bSMatthew G. Knepley ++n; 123995fce210SBarry Smith } 124095fce210SBarry Smith } 1241cd620004SJunchao Zhang ierr = PetscSFCreate(comm,esf);CHKERRQ(ierr); 1242cd620004SJunchao Zhang ierr = PetscSFSetFromOptions(*esf);CHKERRQ(ierr); 1243cd620004SJunchao Zhang ierr = PetscSFSetGraph(*esf,nroots,esf_nleaves,new_ilocal,PETSC_OWN_POINTER,new_iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 1244cd620004SJunchao Zhang ierr = PetscFree2(rootdata,leafmem);CHKERRQ(ierr); 1245f659e5c7SJunchao Zhang } 1246140a1472SStefano Zampini ierr = PetscLogEventEnd(PETSCSF_EmbedSF,sf,0,0,0);CHKERRQ(ierr); 124795fce210SBarry Smith PetscFunctionReturn(0); 124895fce210SBarry Smith } 124995fce210SBarry Smith 12502f5fb4c2SMatthew G. Knepley /*@C 12512f5fb4c2SMatthew G. Knepley PetscSFCreateEmbeddedLeafSF - removes edges from all but the selected leaves, does not remap indices 12522f5fb4c2SMatthew G. Knepley 12532f5fb4c2SMatthew G. Knepley Collective 12542f5fb4c2SMatthew G. Knepley 12552f5fb4c2SMatthew G. Knepley Input Arguments: 12562f5fb4c2SMatthew G. Knepley + sf - original star forest 1257f659e5c7SJunchao Zhang . nselected - number of selected leaves on this process 1258f659e5c7SJunchao Zhang - selected - indices of the selected leaves on this process 12592f5fb4c2SMatthew G. Knepley 12602f5fb4c2SMatthew G. Knepley Output Arguments: 12612f5fb4c2SMatthew G. Knepley . newsf - new star forest 12622f5fb4c2SMatthew G. Knepley 12632f5fb4c2SMatthew G. Knepley Level: advanced 12642f5fb4c2SMatthew G. Knepley 12652f5fb4c2SMatthew G. Knepley .seealso: PetscSFCreateEmbeddedSF(), PetscSFSetGraph(), PetscSFGetGraph() 12662f5fb4c2SMatthew G. Knepley @*/ 1267f659e5c7SJunchao Zhang PetscErrorCode PetscSFCreateEmbeddedLeafSF(PetscSF sf,PetscInt nselected,const PetscInt *selected,PetscSF *newsf) 12682f5fb4c2SMatthew G. Knepley { 1269f659e5c7SJunchao Zhang const PetscSFNode *iremote; 1270f659e5c7SJunchao Zhang PetscSFNode *new_iremote; 1271f659e5c7SJunchao Zhang const PetscInt *ilocal; 1272f659e5c7SJunchao Zhang PetscInt i,nroots,*leaves,*new_ilocal; 1273f659e5c7SJunchao Zhang MPI_Comm comm; 12742f5fb4c2SMatthew G. Knepley PetscErrorCode ierr; 12752f5fb4c2SMatthew G. Knepley 12762f5fb4c2SMatthew G. Knepley PetscFunctionBegin; 12772f5fb4c2SMatthew G. Knepley PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 127829046d53SLisandro Dalcin PetscSFCheckGraphSet(sf,1); 1279f659e5c7SJunchao Zhang if (nselected) PetscValidPointer(selected,3); 12802f5fb4c2SMatthew G. Knepley PetscValidPointer(newsf,4); 12812f5fb4c2SMatthew G. Knepley 1282f659e5c7SJunchao Zhang /* Uniq selected[] and put results in leaves[] */ 1283f659e5c7SJunchao Zhang ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr); 1284f659e5c7SJunchao Zhang ierr = PetscMalloc1(nselected,&leaves);CHKERRQ(ierr); 1285dd5b3ca6SJunchao Zhang ierr = PetscArraycpy(leaves,selected,nselected);CHKERRQ(ierr); 1286f659e5c7SJunchao Zhang ierr = PetscSortedRemoveDupsInt(&nselected,leaves);CHKERRQ(ierr); 1287f659e5c7SJunchao Zhang if (nselected && (leaves[0] < 0 || leaves[nselected-1] >= sf->nleaves)) SETERRQ3(comm,PETSC_ERR_ARG_OUTOFRANGE,"Min/Max leaf indices %D/%D are not in [0,%D)",leaves[0],leaves[nselected-1],sf->nleaves); 1288f659e5c7SJunchao Zhang 1289f659e5c7SJunchao Zhang /* Optimize the routine only when sf is setup and hence we can reuse sf's communication pattern */ 1290f659e5c7SJunchao Zhang if (sf->setupcalled && sf->ops->CreateEmbeddedLeafSF) { 1291f659e5c7SJunchao Zhang ierr = (*sf->ops->CreateEmbeddedLeafSF)(sf,nselected,leaves,newsf);CHKERRQ(ierr); 1292f659e5c7SJunchao Zhang } else { 1293f659e5c7SJunchao Zhang ierr = PetscSFGetGraph(sf,&nroots,NULL,&ilocal,&iremote);CHKERRQ(ierr); 1294f659e5c7SJunchao Zhang ierr = PetscMalloc1(nselected,&new_ilocal);CHKERRQ(ierr); 1295f659e5c7SJunchao Zhang ierr = PetscMalloc1(nselected,&new_iremote);CHKERRQ(ierr); 1296f659e5c7SJunchao Zhang for (i=0; i<nselected; ++i) { 1297f659e5c7SJunchao Zhang const PetscInt l = leaves[i]; 1298f659e5c7SJunchao Zhang new_ilocal[i] = ilocal ? ilocal[l] : l; 1299f659e5c7SJunchao Zhang new_iremote[i].rank = iremote[l].rank; 1300f659e5c7SJunchao Zhang new_iremote[i].index = iremote[l].index; 13012f5fb4c2SMatthew G. Knepley } 13024cc19a0cSStefano Zampini ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,newsf);CHKERRQ(ierr); 1303f659e5c7SJunchao Zhang ierr = PetscSFSetGraph(*newsf,nroots,nselected,new_ilocal,PETSC_OWN_POINTER,new_iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 1304f659e5c7SJunchao Zhang } 1305f659e5c7SJunchao Zhang ierr = PetscFree(leaves);CHKERRQ(ierr); 13062f5fb4c2SMatthew G. Knepley PetscFunctionReturn(0); 13072f5fb4c2SMatthew G. Knepley } 13082f5fb4c2SMatthew G. Knepley 130995fce210SBarry Smith /*@C 13103482bfa8SJunchao Zhang PetscSFBcastAndOpBegin - begin pointwise broadcast with root value being reduced to leaf value, to be concluded with call to PetscSFBcastAndOpEnd() 13113482bfa8SJunchao Zhang 13123482bfa8SJunchao Zhang Collective on PetscSF 13133482bfa8SJunchao Zhang 13143482bfa8SJunchao Zhang Input Arguments: 13153482bfa8SJunchao Zhang + sf - star forest on which to communicate 13163482bfa8SJunchao Zhang . unit - data type associated with each node 13173482bfa8SJunchao Zhang . rootdata - buffer to broadcast 13183482bfa8SJunchao Zhang - op - operation to use for reduction 13193482bfa8SJunchao Zhang 13203482bfa8SJunchao Zhang Output Arguments: 13213482bfa8SJunchao Zhang . leafdata - buffer to be reduced with values from each leaf's respective root 13223482bfa8SJunchao Zhang 13233482bfa8SJunchao Zhang Level: intermediate 13243482bfa8SJunchao Zhang 13253482bfa8SJunchao Zhang .seealso: PetscSFBcastAndOpEnd(), PetscSFBcastBegin(), PetscSFBcastEnd() 13263482bfa8SJunchao Zhang @*/ 13273482bfa8SJunchao Zhang PetscErrorCode PetscSFBcastAndOpBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op) 13283482bfa8SJunchao Zhang { 13293482bfa8SJunchao Zhang PetscErrorCode ierr; 1330eb02082bSJunchao Zhang PetscMemType rootmtype,leafmtype; 13313482bfa8SJunchao Zhang 13323482bfa8SJunchao Zhang PetscFunctionBegin; 13333482bfa8SJunchao Zhang PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 13343482bfa8SJunchao Zhang ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 13353482bfa8SJunchao Zhang ierr = PetscLogEventBegin(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr); 1336eb02082bSJunchao Zhang ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr); 1337eb02082bSJunchao Zhang ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr); 1338eb02082bSJunchao Zhang ierr = (*sf->ops->BcastAndOpBegin)(sf,unit,rootmtype,rootdata,leafmtype,leafdata,op);CHKERRQ(ierr); 13393482bfa8SJunchao Zhang ierr = PetscLogEventEnd(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr); 13403482bfa8SJunchao Zhang PetscFunctionReturn(0); 13413482bfa8SJunchao Zhang } 13423482bfa8SJunchao Zhang 13433482bfa8SJunchao Zhang /*@C 13443482bfa8SJunchao Zhang PetscSFBcastAndOpEnd - end a broadcast & reduce operation started with PetscSFBcastAndOpBegin() 13453482bfa8SJunchao Zhang 13463482bfa8SJunchao Zhang Collective 13473482bfa8SJunchao Zhang 13483482bfa8SJunchao Zhang Input Arguments: 13493482bfa8SJunchao Zhang + sf - star forest 13503482bfa8SJunchao Zhang . unit - data type 13513482bfa8SJunchao Zhang . rootdata - buffer to broadcast 13523482bfa8SJunchao Zhang - op - operation to use for reduction 13533482bfa8SJunchao Zhang 13543482bfa8SJunchao Zhang Output Arguments: 13553482bfa8SJunchao Zhang . leafdata - buffer to be reduced with values from each leaf's respective root 13563482bfa8SJunchao Zhang 13573482bfa8SJunchao Zhang Level: intermediate 13583482bfa8SJunchao Zhang 13593482bfa8SJunchao Zhang .seealso: PetscSFSetGraph(), PetscSFReduceEnd() 13603482bfa8SJunchao Zhang @*/ 13613482bfa8SJunchao Zhang PetscErrorCode PetscSFBcastAndOpEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op) 13623482bfa8SJunchao Zhang { 13633482bfa8SJunchao Zhang PetscErrorCode ierr; 13643482bfa8SJunchao Zhang 13653482bfa8SJunchao Zhang PetscFunctionBegin; 13663482bfa8SJunchao Zhang PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 13673482bfa8SJunchao Zhang ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 13683482bfa8SJunchao Zhang ierr = PetscLogEventBegin(PETSCSF_BcastAndOpEnd,sf,0,0,0);CHKERRQ(ierr); 136900816365SJunchao Zhang ierr = (*sf->ops->BcastAndOpEnd)(sf,unit,rootdata,leafdata,op);CHKERRQ(ierr); 13703482bfa8SJunchao Zhang ierr = PetscLogEventEnd(PETSCSF_BcastAndOpEnd,sf,0,0,0);CHKERRQ(ierr); 13713482bfa8SJunchao Zhang PetscFunctionReturn(0); 13723482bfa8SJunchao Zhang } 13733482bfa8SJunchao Zhang 13743482bfa8SJunchao Zhang /*@C 137595fce210SBarry Smith PetscSFReduceBegin - begin reduction of leafdata into rootdata, to be completed with call to PetscSFReduceEnd() 137695fce210SBarry Smith 137795fce210SBarry Smith Collective 137895fce210SBarry Smith 137995fce210SBarry Smith Input Arguments: 138095fce210SBarry Smith + sf - star forest 138195fce210SBarry Smith . unit - data type 138295fce210SBarry Smith . leafdata - values to reduce 138395fce210SBarry Smith - op - reduction operation 138495fce210SBarry Smith 138595fce210SBarry Smith Output Arguments: 138695fce210SBarry Smith . rootdata - result of reduction of values from all leaves of each root 138795fce210SBarry Smith 138895fce210SBarry Smith Level: intermediate 138995fce210SBarry Smith 139095fce210SBarry Smith .seealso: PetscSFBcastBegin() 139195fce210SBarry Smith @*/ 139295fce210SBarry Smith PetscErrorCode PetscSFReduceBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op) 139395fce210SBarry Smith { 139495fce210SBarry Smith PetscErrorCode ierr; 1395eb02082bSJunchao Zhang PetscMemType rootmtype,leafmtype; 139695fce210SBarry Smith 139795fce210SBarry Smith PetscFunctionBegin; 139895fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 139995fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 140029046d53SLisandro Dalcin ierr = PetscLogEventBegin(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr); 1401eb02082bSJunchao Zhang ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr); 1402eb02082bSJunchao Zhang ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr); 1403eb02082bSJunchao Zhang ierr = (sf->ops->ReduceBegin)(sf,unit,leafmtype,leafdata,rootmtype,rootdata,op);CHKERRQ(ierr); 1404563ffbabSMatthew G. Knepley ierr = PetscLogEventEnd(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr); 140595fce210SBarry Smith PetscFunctionReturn(0); 140695fce210SBarry Smith } 140795fce210SBarry Smith 140895fce210SBarry Smith /*@C 140995fce210SBarry Smith PetscSFReduceEnd - end a reduction operation started with PetscSFReduceBegin() 141095fce210SBarry Smith 141195fce210SBarry Smith Collective 141295fce210SBarry Smith 141395fce210SBarry Smith Input Arguments: 141495fce210SBarry Smith + sf - star forest 141595fce210SBarry Smith . unit - data type 141695fce210SBarry Smith . leafdata - values to reduce 141795fce210SBarry Smith - op - reduction operation 141895fce210SBarry Smith 141995fce210SBarry Smith Output Arguments: 142095fce210SBarry Smith . rootdata - result of reduction of values from all leaves of each root 142195fce210SBarry Smith 142295fce210SBarry Smith Level: intermediate 142395fce210SBarry Smith 142495fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFBcastEnd() 142595fce210SBarry Smith @*/ 142695fce210SBarry Smith PetscErrorCode PetscSFReduceEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op) 142795fce210SBarry Smith { 142895fce210SBarry Smith PetscErrorCode ierr; 142995fce210SBarry Smith 143095fce210SBarry Smith PetscFunctionBegin; 143195fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 143295fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 143329046d53SLisandro Dalcin ierr = PetscLogEventBegin(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr); 143400816365SJunchao Zhang ierr = (*sf->ops->ReduceEnd)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr); 1435563ffbabSMatthew G. Knepley ierr = PetscLogEventEnd(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr); 143695fce210SBarry Smith PetscFunctionReturn(0); 143795fce210SBarry Smith } 143895fce210SBarry Smith 143995fce210SBarry Smith /*@C 1440a1729e3fSJunchao Zhang PetscSFFetchAndOpBegin - begin operation that fetches values from root and updates atomically by applying operation using my leaf value, to be completed with PetscSFFetchAndOpEnd() 1441a1729e3fSJunchao Zhang 1442a1729e3fSJunchao Zhang Collective 1443a1729e3fSJunchao Zhang 1444a1729e3fSJunchao Zhang Input Arguments: 1445a1729e3fSJunchao Zhang + sf - star forest 1446a1729e3fSJunchao Zhang . unit - data type 1447a1729e3fSJunchao Zhang . leafdata - leaf values to use in reduction 1448a1729e3fSJunchao Zhang - op - operation to use for reduction 1449a1729e3fSJunchao Zhang 1450a1729e3fSJunchao Zhang Output Arguments: 1451a1729e3fSJunchao Zhang + rootdata - root values to be updated, input state is seen by first process to perform an update 1452a1729e3fSJunchao Zhang - leafupdate - state at each leaf's respective root immediately prior to my atomic update 1453a1729e3fSJunchao Zhang 1454a1729e3fSJunchao Zhang Level: advanced 1455a1729e3fSJunchao Zhang 1456a1729e3fSJunchao Zhang Note: 1457a1729e3fSJunchao Zhang The update is only atomic at the granularity provided by the hardware. Different roots referenced by the same process 1458a1729e3fSJunchao Zhang might be updated in a different order. Furthermore, if a composite type is used for the unit datatype, atomicity is 1459a1729e3fSJunchao Zhang not guaranteed across the whole vertex. Therefore, this function is mostly only used with primitive types such as 1460a1729e3fSJunchao Zhang integers. 1461a1729e3fSJunchao Zhang 1462a1729e3fSJunchao Zhang .seealso: PetscSFComputeDegreeBegin(), PetscSFReduceBegin(), PetscSFSetGraph() 1463a1729e3fSJunchao Zhang @*/ 1464a1729e3fSJunchao Zhang PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op) 1465a1729e3fSJunchao Zhang { 1466a1729e3fSJunchao Zhang PetscErrorCode ierr; 1467eb02082bSJunchao Zhang PetscMemType rootmtype,leafmtype,leafupdatemtype; 1468a1729e3fSJunchao Zhang 1469a1729e3fSJunchao Zhang PetscFunctionBegin; 1470a1729e3fSJunchao Zhang PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 1471a1729e3fSJunchao Zhang ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 1472a1729e3fSJunchao Zhang ierr = PetscLogEventBegin(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr); 1473eb02082bSJunchao Zhang ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr); 1474eb02082bSJunchao Zhang ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr); 1475eb02082bSJunchao Zhang ierr = PetscGetMemType(leafupdate,&leafupdatemtype);CHKERRQ(ierr); 1476eb02082bSJunchao Zhang if (leafmtype != leafupdatemtype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for leafdata and leafupdate in different memory types"); 1477eb02082bSJunchao Zhang ierr = (*sf->ops->FetchAndOpBegin)(sf,unit,rootmtype,rootdata,leafmtype,leafdata,leafupdate,op);CHKERRQ(ierr); 1478a1729e3fSJunchao Zhang ierr = PetscLogEventEnd(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr); 1479a1729e3fSJunchao Zhang PetscFunctionReturn(0); 1480a1729e3fSJunchao Zhang } 1481a1729e3fSJunchao Zhang 1482a1729e3fSJunchao Zhang /*@C 1483a1729e3fSJunchao Zhang PetscSFFetchAndOpEnd - end operation started in matching call to PetscSFFetchAndOpBegin() to fetch values from roots and update atomically by applying operation using my leaf value 1484a1729e3fSJunchao Zhang 1485a1729e3fSJunchao Zhang Collective 1486a1729e3fSJunchao Zhang 1487a1729e3fSJunchao Zhang Input Arguments: 1488a1729e3fSJunchao Zhang + sf - star forest 1489a1729e3fSJunchao Zhang . unit - data type 1490a1729e3fSJunchao Zhang . leafdata - leaf values to use in reduction 1491a1729e3fSJunchao Zhang - op - operation to use for reduction 1492a1729e3fSJunchao Zhang 1493a1729e3fSJunchao Zhang Output Arguments: 1494a1729e3fSJunchao Zhang + rootdata - root values to be updated, input state is seen by first process to perform an update 1495a1729e3fSJunchao Zhang - leafupdate - state at each leaf's respective root immediately prior to my atomic update 1496a1729e3fSJunchao Zhang 1497a1729e3fSJunchao Zhang Level: advanced 1498a1729e3fSJunchao Zhang 1499a1729e3fSJunchao Zhang .seealso: PetscSFComputeDegreeEnd(), PetscSFReduceEnd(), PetscSFSetGraph() 1500a1729e3fSJunchao Zhang @*/ 1501a1729e3fSJunchao Zhang PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op) 1502a1729e3fSJunchao Zhang { 1503a1729e3fSJunchao Zhang PetscErrorCode ierr; 1504a1729e3fSJunchao Zhang 1505a1729e3fSJunchao Zhang PetscFunctionBegin; 1506a1729e3fSJunchao Zhang PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 1507a1729e3fSJunchao Zhang ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 1508a1729e3fSJunchao Zhang ierr = PetscLogEventBegin(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr); 150900816365SJunchao Zhang ierr = (*sf->ops->FetchAndOpEnd)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr); 1510a1729e3fSJunchao Zhang ierr = PetscLogEventEnd(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr); 1511a1729e3fSJunchao Zhang PetscFunctionReturn(0); 1512a1729e3fSJunchao Zhang } 1513a1729e3fSJunchao Zhang 1514a1729e3fSJunchao Zhang /*@C 151595fce210SBarry Smith PetscSFComputeDegreeBegin - begin computation of degree for each root vertex, to be completed with PetscSFComputeDegreeEnd() 151695fce210SBarry Smith 151795fce210SBarry Smith Collective 151895fce210SBarry Smith 151995fce210SBarry Smith Input Arguments: 152095fce210SBarry Smith . sf - star forest 152195fce210SBarry Smith 152295fce210SBarry Smith Output Arguments: 152395fce210SBarry Smith . degree - degree of each root vertex 152495fce210SBarry Smith 152595fce210SBarry Smith Level: advanced 152695fce210SBarry Smith 1527ffe67aa5SVáclav Hapla Notes: 1528ffe67aa5SVáclav Hapla The returned array is owned by PetscSF and automatically freed by PetscSFDestroy(). Hence no need to call PetscFree() on it. 1529ffe67aa5SVáclav Hapla 153095fce210SBarry Smith .seealso: PetscSFGatherBegin() 153195fce210SBarry Smith @*/ 153295fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf,const PetscInt **degree) 153395fce210SBarry Smith { 153495fce210SBarry Smith PetscErrorCode ierr; 153595fce210SBarry Smith 153695fce210SBarry Smith PetscFunctionBegin; 153795fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 153895fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 153995fce210SBarry Smith PetscValidPointer(degree,2); 1540803bd9e8SMatthew G. Knepley if (!sf->degreeknown) { 15415b0d146aSStefano Zampini PetscInt i, nroots = sf->nroots, maxlocal; 1542803bd9e8SMatthew G. Knepley if (sf->degree) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Calls to PetscSFComputeDegreeBegin() cannot be nested."); 15435b0d146aSStefano Zampini maxlocal = sf->maxleaf-sf->minleaf+1; 154429046d53SLisandro Dalcin ierr = PetscMalloc1(nroots,&sf->degree);CHKERRQ(ierr); 154529046d53SLisandro Dalcin ierr = PetscMalloc1(PetscMax(maxlocal,1),&sf->degreetmp);CHKERRQ(ierr); /* allocate at least one entry, see check in PetscSFComputeDegreeEnd() */ 154629046d53SLisandro Dalcin for (i=0; i<nroots; i++) sf->degree[i] = 0; 15479837ea96SMatthew G. Knepley for (i=0; i<maxlocal; i++) sf->degreetmp[i] = 1; 15485b0d146aSStefano Zampini ierr = PetscSFReduceBegin(sf,MPIU_INT,sf->degreetmp-sf->minleaf,sf->degree,MPI_SUM);CHKERRQ(ierr); 154995fce210SBarry Smith } 155095fce210SBarry Smith *degree = NULL; 155195fce210SBarry Smith PetscFunctionReturn(0); 155295fce210SBarry Smith } 155395fce210SBarry Smith 155495fce210SBarry Smith /*@C 155595fce210SBarry Smith PetscSFComputeDegreeEnd - complete computation of degree for each root vertex, started with PetscSFComputeDegreeBegin() 155695fce210SBarry Smith 155795fce210SBarry Smith Collective 155895fce210SBarry Smith 155995fce210SBarry Smith Input Arguments: 156095fce210SBarry Smith . sf - star forest 156195fce210SBarry Smith 156295fce210SBarry Smith Output Arguments: 156395fce210SBarry Smith . degree - degree of each root vertex 156495fce210SBarry Smith 156595fce210SBarry Smith Level: developer 156695fce210SBarry Smith 1567ffe67aa5SVáclav Hapla Notes: 1568ffe67aa5SVáclav Hapla The returned array is owned by PetscSF and automatically freed by PetscSFDestroy(). Hence no need to call PetscFree() on it. 1569ffe67aa5SVáclav Hapla 157095fce210SBarry Smith .seealso: 157195fce210SBarry Smith @*/ 157295fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf,const PetscInt **degree) 157395fce210SBarry Smith { 157495fce210SBarry Smith PetscErrorCode ierr; 157595fce210SBarry Smith 157695fce210SBarry Smith PetscFunctionBegin; 157795fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 157895fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 157929046d53SLisandro Dalcin PetscValidPointer(degree,2); 158095fce210SBarry Smith if (!sf->degreeknown) { 158129046d53SLisandro Dalcin if (!sf->degreetmp) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFComputeDegreeBegin() before PetscSFComputeDegreeEnd()"); 15825b0d146aSStefano Zampini ierr = PetscSFReduceEnd(sf,MPIU_INT,sf->degreetmp-sf->minleaf,sf->degree,MPI_SUM);CHKERRQ(ierr); 158395fce210SBarry Smith ierr = PetscFree(sf->degreetmp);CHKERRQ(ierr); 158495fce210SBarry Smith sf->degreeknown = PETSC_TRUE; 158595fce210SBarry Smith } 158695fce210SBarry Smith *degree = sf->degree; 158795fce210SBarry Smith PetscFunctionReturn(0); 158895fce210SBarry Smith } 158995fce210SBarry Smith 1590673100f5SVaclav Hapla 1591673100f5SVaclav Hapla /*@C 159266dfcd1aSVaclav Hapla PetscSFComputeMultiRootOriginalNumbering - Returns original numbering of multi-roots (roots of multi-SF returned by PetscSFGetMultiSF()). 159366dfcd1aSVaclav Hapla Each multi-root is assigned index of the corresponding original root. 1594673100f5SVaclav Hapla 1595673100f5SVaclav Hapla Collective 1596673100f5SVaclav Hapla 1597673100f5SVaclav Hapla Input Arguments: 1598673100f5SVaclav Hapla + sf - star forest 1599673100f5SVaclav Hapla - degree - degree of each root vertex, computed with PetscSFComputeDegreeBegin()/PetscSFComputeDegreeEnd() 1600673100f5SVaclav Hapla 1601673100f5SVaclav Hapla Output Arguments: 160266dfcd1aSVaclav Hapla + nMultiRoots - (optional) number of multi-roots (roots of multi-SF) 160366dfcd1aSVaclav Hapla - multiRootsOrigNumbering - original indices of multi-roots; length of this array is nMultiRoots 1604673100f5SVaclav Hapla 1605673100f5SVaclav Hapla Level: developer 1606673100f5SVaclav Hapla 1607ffe67aa5SVáclav Hapla Notes: 1608ffe67aa5SVáclav Hapla The returned array multiRootsOrigNumbering is newly allocated and should be destroyed with PetscFree() when no longer needed. 1609ffe67aa5SVáclav Hapla 1610673100f5SVaclav Hapla .seealso: PetscSFComputeDegreeBegin(), PetscSFComputeDegreeEnd(), PetscSFGetMultiSF() 1611673100f5SVaclav Hapla @*/ 161266dfcd1aSVaclav Hapla PetscErrorCode PetscSFComputeMultiRootOriginalNumbering(PetscSF sf, const PetscInt degree[], PetscInt *nMultiRoots, PetscInt *multiRootsOrigNumbering[]) 1613673100f5SVaclav Hapla { 1614673100f5SVaclav Hapla PetscSF msf; 1615673100f5SVaclav Hapla PetscInt i, j, k, nroots, nmroots; 1616673100f5SVaclav Hapla PetscErrorCode ierr; 1617673100f5SVaclav Hapla 1618673100f5SVaclav Hapla PetscFunctionBegin; 1619673100f5SVaclav Hapla PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 1620673100f5SVaclav Hapla ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 162129328920SVaclav Hapla if (nroots) PetscValidIntPointer(degree,2); 162266dfcd1aSVaclav Hapla if (nMultiRoots) PetscValidIntPointer(nMultiRoots,3); 162366dfcd1aSVaclav Hapla PetscValidPointer(multiRootsOrigNumbering,4); 1624673100f5SVaclav Hapla ierr = PetscSFGetMultiSF(sf,&msf);CHKERRQ(ierr); 1625673100f5SVaclav Hapla ierr = PetscSFGetGraph(msf, &nmroots, NULL, NULL, NULL);CHKERRQ(ierr); 162666dfcd1aSVaclav Hapla ierr = PetscMalloc1(nmroots, multiRootsOrigNumbering);CHKERRQ(ierr); 1627673100f5SVaclav Hapla for (i=0,j=0,k=0; i<nroots; i++) { 1628673100f5SVaclav Hapla if (!degree[i]) continue; 1629673100f5SVaclav Hapla for (j=0; j<degree[i]; j++,k++) { 163066dfcd1aSVaclav Hapla (*multiRootsOrigNumbering)[k] = i; 1631673100f5SVaclav Hapla } 1632673100f5SVaclav Hapla } 1633673100f5SVaclav Hapla if (PetscUnlikely(k != nmroots)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sanity check fail"); 163466dfcd1aSVaclav Hapla if (nMultiRoots) *nMultiRoots = nmroots; 1635673100f5SVaclav Hapla PetscFunctionReturn(0); 1636673100f5SVaclav Hapla } 1637673100f5SVaclav Hapla 163895fce210SBarry Smith /*@C 163995fce210SBarry Smith PetscSFGatherBegin - begin pointwise gather of all leaves into multi-roots, to be completed with PetscSFGatherEnd() 164095fce210SBarry Smith 164195fce210SBarry Smith Collective 164295fce210SBarry Smith 164395fce210SBarry Smith Input Arguments: 164495fce210SBarry Smith + sf - star forest 164595fce210SBarry Smith . unit - data type 164695fce210SBarry Smith - leafdata - leaf data to gather to roots 164795fce210SBarry Smith 164895fce210SBarry Smith Output Argument: 164995fce210SBarry Smith . multirootdata - root buffer to gather into, amount of space per root is equal to its degree 165095fce210SBarry Smith 165195fce210SBarry Smith Level: intermediate 165295fce210SBarry Smith 165395fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin() 165495fce210SBarry Smith @*/ 165595fce210SBarry Smith PetscErrorCode PetscSFGatherBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata) 165695fce210SBarry Smith { 165795fce210SBarry Smith PetscErrorCode ierr; 1658a5526d50SJunchao Zhang PetscSF multi = NULL; 165995fce210SBarry Smith 166095fce210SBarry Smith PetscFunctionBegin; 166195fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 166229046d53SLisandro Dalcin ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 166395fce210SBarry Smith ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr); 16648bfbc91cSJed Brown ierr = PetscSFReduceBegin(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr); 166595fce210SBarry Smith PetscFunctionReturn(0); 166695fce210SBarry Smith } 166795fce210SBarry Smith 166895fce210SBarry Smith /*@C 166995fce210SBarry Smith PetscSFGatherEnd - ends pointwise gather operation that was started with PetscSFGatherBegin() 167095fce210SBarry Smith 167195fce210SBarry Smith Collective 167295fce210SBarry Smith 167395fce210SBarry Smith Input Arguments: 167495fce210SBarry Smith + sf - star forest 167595fce210SBarry Smith . unit - data type 167695fce210SBarry Smith - leafdata - leaf data to gather to roots 167795fce210SBarry Smith 167895fce210SBarry Smith Output Argument: 167995fce210SBarry Smith . multirootdata - root buffer to gather into, amount of space per root is equal to its degree 168095fce210SBarry Smith 168195fce210SBarry Smith Level: intermediate 168295fce210SBarry Smith 168395fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd() 168495fce210SBarry Smith @*/ 168595fce210SBarry Smith PetscErrorCode PetscSFGatherEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata) 168695fce210SBarry Smith { 168795fce210SBarry Smith PetscErrorCode ierr; 1688a5526d50SJunchao Zhang PetscSF multi = NULL; 168995fce210SBarry Smith 169095fce210SBarry Smith PetscFunctionBegin; 169195fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 169295fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 169395fce210SBarry Smith ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr); 16948bfbc91cSJed Brown ierr = PetscSFReduceEnd(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr); 169595fce210SBarry Smith PetscFunctionReturn(0); 169695fce210SBarry Smith } 169795fce210SBarry Smith 169895fce210SBarry Smith /*@C 169995fce210SBarry Smith PetscSFScatterBegin - begin pointwise scatter operation from multi-roots to leaves, to be completed with PetscSFScatterEnd() 170095fce210SBarry Smith 170195fce210SBarry Smith Collective 170295fce210SBarry Smith 170395fce210SBarry Smith Input Arguments: 170495fce210SBarry Smith + sf - star forest 170595fce210SBarry Smith . unit - data type 170695fce210SBarry Smith - multirootdata - root buffer to send to each leaf, one unit of data per leaf 170795fce210SBarry Smith 170895fce210SBarry Smith Output Argument: 170995fce210SBarry Smith . leafdata - leaf data to be update with personal data from each respective root 171095fce210SBarry Smith 171195fce210SBarry Smith Level: intermediate 171295fce210SBarry Smith 171395fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin() 171495fce210SBarry Smith @*/ 171595fce210SBarry Smith PetscErrorCode PetscSFScatterBegin(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata) 171695fce210SBarry Smith { 171795fce210SBarry Smith PetscErrorCode ierr; 1718a5526d50SJunchao Zhang PetscSF multi = NULL; 171995fce210SBarry Smith 172095fce210SBarry Smith PetscFunctionBegin; 172195fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 172295fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 172395fce210SBarry Smith ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr); 172495fce210SBarry Smith ierr = PetscSFBcastBegin(multi,unit,multirootdata,leafdata);CHKERRQ(ierr); 172595fce210SBarry Smith PetscFunctionReturn(0); 172695fce210SBarry Smith } 172795fce210SBarry Smith 172895fce210SBarry Smith /*@C 172995fce210SBarry Smith PetscSFScatterEnd - ends pointwise scatter operation that was started with PetscSFScatterBegin() 173095fce210SBarry Smith 173195fce210SBarry Smith Collective 173295fce210SBarry Smith 173395fce210SBarry Smith Input Arguments: 173495fce210SBarry Smith + sf - star forest 173595fce210SBarry Smith . unit - data type 173695fce210SBarry Smith - multirootdata - root buffer to send to each leaf, one unit of data per leaf 173795fce210SBarry Smith 173895fce210SBarry Smith Output Argument: 173995fce210SBarry Smith . leafdata - leaf data to be update with personal data from each respective root 174095fce210SBarry Smith 174195fce210SBarry Smith Level: intermediate 174295fce210SBarry Smith 174395fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd() 174495fce210SBarry Smith @*/ 174595fce210SBarry Smith PetscErrorCode PetscSFScatterEnd(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata) 174695fce210SBarry Smith { 174795fce210SBarry Smith PetscErrorCode ierr; 1748a5526d50SJunchao Zhang PetscSF multi = NULL; 174995fce210SBarry Smith 175095fce210SBarry Smith PetscFunctionBegin; 175195fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 175295fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 175395fce210SBarry Smith ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr); 175495fce210SBarry Smith ierr = PetscSFBcastEnd(multi,unit,multirootdata,leafdata);CHKERRQ(ierr); 175595fce210SBarry Smith PetscFunctionReturn(0); 175695fce210SBarry Smith } 1757a7b3aa13SAta Mesgarnejad 1758a072220fSLawrence Mitchell static PetscErrorCode PetscSFCheckLeavesUnique_Private(PetscSF sf) 1759a072220fSLawrence Mitchell { 1760a072220fSLawrence Mitchell PetscInt i, n, nleaves; 1761a072220fSLawrence Mitchell const PetscInt *ilocal = NULL; 1762a072220fSLawrence Mitchell PetscHSetI seen; 1763a072220fSLawrence Mitchell PetscErrorCode ierr; 1764a072220fSLawrence Mitchell 1765a072220fSLawrence Mitchell PetscFunctionBegin; 1766*76bd3646SJed Brown if (!PetscDefined(USE_DEBUG)) PetscFunctionReturn(0); 1767a072220fSLawrence Mitchell ierr = PetscSFGetGraph(sf,NULL,&nleaves,&ilocal,NULL);CHKERRQ(ierr); 1768a072220fSLawrence Mitchell ierr = PetscHSetICreate(&seen);CHKERRQ(ierr); 1769a072220fSLawrence Mitchell for (i = 0; i < nleaves; i++) { 1770a072220fSLawrence Mitchell const PetscInt leaf = ilocal ? ilocal[i] : i; 1771a072220fSLawrence Mitchell ierr = PetscHSetIAdd(seen,leaf);CHKERRQ(ierr); 1772a072220fSLawrence Mitchell } 1773a072220fSLawrence Mitchell ierr = PetscHSetIGetSize(seen,&n);CHKERRQ(ierr); 1774a072220fSLawrence Mitchell if (n != nleaves) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Provided leaves have repeated values: all leaves must be unique"); 1775a072220fSLawrence Mitchell ierr = PetscHSetIDestroy(&seen);CHKERRQ(ierr); 1776a072220fSLawrence Mitchell PetscFunctionReturn(0); 1777a072220fSLawrence Mitchell } 177854729392SStefano Zampini 1779a7b3aa13SAta Mesgarnejad /*@ 178004c0ada0SJunchao Zhang PetscSFCompose - Compose a new PetscSF by putting the second SF under the first one in a top (roots) down (leaves) view 1781a7b3aa13SAta Mesgarnejad 1782a7b3aa13SAta Mesgarnejad Input Parameters: 178354729392SStefano Zampini + sfA - The first PetscSF 178454729392SStefano Zampini - sfB - The second PetscSF 1785a7b3aa13SAta Mesgarnejad 1786a7b3aa13SAta Mesgarnejad Output Parameters: 178754729392SStefano Zampini . sfBA - The composite SF 1788a7b3aa13SAta Mesgarnejad 1789a7b3aa13SAta Mesgarnejad Level: developer 1790a7b3aa13SAta Mesgarnejad 1791a072220fSLawrence Mitchell Notes: 179254729392SStefano Zampini Currently, the two SFs must be defined on congruent communicators and they must be true star 179354729392SStefano Zampini forests, i.e. the same leaf is not connected with different roots. 179454729392SStefano Zampini 179554729392SStefano Zampini sfA's leaf space and sfB's root space might be partially overlapped. The composition builds 179654729392SStefano Zampini a graph with sfA's roots and sfB's leaves only when there is a path between them. Unconnected 179754729392SStefano Zampini nodes (roots or leaves) are not in sfBA. Doing a Bcast on the new SF is equivalent to doing a 179854729392SStefano Zampini Bcast on sfA, then a Bcast on sfB, on connected nodes. 1799a072220fSLawrence Mitchell 180004c0ada0SJunchao Zhang .seealso: PetscSF, PetscSFComposeInverse(), PetscSFGetGraph(), PetscSFSetGraph() 1801a7b3aa13SAta Mesgarnejad @*/ 1802a7b3aa13SAta Mesgarnejad PetscErrorCode PetscSFCompose(PetscSF sfA,PetscSF sfB,PetscSF *sfBA) 1803a7b3aa13SAta Mesgarnejad { 180404c0ada0SJunchao Zhang PetscErrorCode ierr; 1805a7b3aa13SAta Mesgarnejad const PetscSFNode *remotePointsA,*remotePointsB; 1806d41018fbSJunchao Zhang PetscSFNode *remotePointsBA=NULL,*reorderedRemotePointsA = NULL,*leafdataB; 180754729392SStefano Zampini const PetscInt *localPointsA,*localPointsB; 180854729392SStefano Zampini PetscInt *localPointsBA; 180954729392SStefano Zampini PetscInt i,numRootsA,numLeavesA,numRootsB,numLeavesB,minleaf,maxleaf,numLeavesBA; 181054729392SStefano Zampini PetscBool denseB; 1811a7b3aa13SAta Mesgarnejad 1812a7b3aa13SAta Mesgarnejad PetscFunctionBegin; 1813a7b3aa13SAta Mesgarnejad PetscValidHeaderSpecific(sfA,PETSCSF_CLASSID,1); 181429046d53SLisandro Dalcin PetscSFCheckGraphSet(sfA,1); 181529046d53SLisandro Dalcin PetscValidHeaderSpecific(sfB,PETSCSF_CLASSID,2); 181629046d53SLisandro Dalcin PetscSFCheckGraphSet(sfB,2); 181754729392SStefano Zampini PetscCheckSameComm(sfA,1,sfB,2); 181829046d53SLisandro Dalcin PetscValidPointer(sfBA,3); 181954729392SStefano Zampini ierr = PetscSFCheckLeavesUnique_Private(sfA);CHKERRQ(ierr); 182054729392SStefano Zampini ierr = PetscSFCheckLeavesUnique_Private(sfB);CHKERRQ(ierr); 182154729392SStefano Zampini 1822a7b3aa13SAta Mesgarnejad ierr = PetscSFGetGraph(sfA,&numRootsA,&numLeavesA,&localPointsA,&remotePointsA);CHKERRQ(ierr); 1823a7b3aa13SAta Mesgarnejad ierr = PetscSFGetGraph(sfB,&numRootsB,&numLeavesB,&localPointsB,&remotePointsB);CHKERRQ(ierr); 1824d41018fbSJunchao Zhang if (localPointsA) { 182554729392SStefano Zampini ierr = PetscMalloc1(numRootsB,&reorderedRemotePointsA);CHKERRQ(ierr); 182654729392SStefano Zampini for (i=0; i<numRootsB; i++) { 182754729392SStefano Zampini reorderedRemotePointsA[i].rank = -1; 182854729392SStefano Zampini reorderedRemotePointsA[i].index = -1; 182954729392SStefano Zampini } 183054729392SStefano Zampini for (i=0; i<numLeavesA; i++) { 183154729392SStefano Zampini if (localPointsA[i] >= numRootsB) continue; 183254729392SStefano Zampini reorderedRemotePointsA[localPointsA[i]] = remotePointsA[i]; 183354729392SStefano Zampini } 1834d41018fbSJunchao Zhang remotePointsA = reorderedRemotePointsA; 1835d41018fbSJunchao Zhang } 1836d41018fbSJunchao Zhang ierr = PetscSFGetLeafRange(sfB,&minleaf,&maxleaf);CHKERRQ(ierr); 1837d41018fbSJunchao Zhang ierr = PetscMalloc1(maxleaf-minleaf+1,&leafdataB);CHKERRQ(ierr); 1838d41018fbSJunchao Zhang ierr = PetscSFBcastBegin(sfB,MPIU_2INT,remotePointsA,leafdataB-minleaf);CHKERRQ(ierr); 1839d41018fbSJunchao Zhang ierr = PetscSFBcastEnd(sfB,MPIU_2INT,remotePointsA,leafdataB-minleaf);CHKERRQ(ierr); 1840d41018fbSJunchao Zhang ierr = PetscFree(reorderedRemotePointsA);CHKERRQ(ierr); 1841d41018fbSJunchao Zhang 184254729392SStefano Zampini denseB = (PetscBool)!localPointsB; 184354729392SStefano Zampini for (i=0,numLeavesBA=0; i<numLeavesB; i++) { 184454729392SStefano Zampini if (leafdataB[localPointsB ? localPointsB[i]-minleaf : i].rank == -1) denseB = PETSC_FALSE; 184554729392SStefano Zampini else numLeavesBA++; 184654729392SStefano Zampini } 184754729392SStefano Zampini if (denseB) { 1848d41018fbSJunchao Zhang localPointsBA = NULL; 1849d41018fbSJunchao Zhang remotePointsBA = leafdataB; 1850d41018fbSJunchao Zhang } else { 185154729392SStefano Zampini ierr = PetscMalloc1(numLeavesBA,&localPointsBA);CHKERRQ(ierr); 185254729392SStefano Zampini ierr = PetscMalloc1(numLeavesBA,&remotePointsBA);CHKERRQ(ierr); 185354729392SStefano Zampini for (i=0,numLeavesBA=0; i<numLeavesB; i++) { 185454729392SStefano Zampini const PetscInt l = localPointsB ? localPointsB[i] : i; 185554729392SStefano Zampini 185654729392SStefano Zampini if (leafdataB[l-minleaf].rank == -1) continue; 185754729392SStefano Zampini remotePointsBA[numLeavesBA] = leafdataB[l-minleaf]; 185854729392SStefano Zampini localPointsBA[numLeavesBA] = l; 185954729392SStefano Zampini numLeavesBA++; 186054729392SStefano Zampini } 1861d41018fbSJunchao Zhang ierr = PetscFree(leafdataB);CHKERRQ(ierr); 1862d41018fbSJunchao Zhang } 186354729392SStefano Zampini ierr = PetscSFCreate(PetscObjectComm((PetscObject)sfA),sfBA);CHKERRQ(ierr); 186454729392SStefano Zampini ierr = PetscSFSetGraph(*sfBA,numRootsA,numLeavesBA,localPointsBA,PETSC_OWN_POINTER,remotePointsBA,PETSC_OWN_POINTER);CHKERRQ(ierr); 1865a7b3aa13SAta Mesgarnejad PetscFunctionReturn(0); 1866a7b3aa13SAta Mesgarnejad } 18671c6ba672SJunchao Zhang 186804c0ada0SJunchao Zhang /*@ 186954729392SStefano Zampini PetscSFComposeInverse - Compose a new PetscSF by putting the inverse of the second SF under the first one 187004c0ada0SJunchao Zhang 187104c0ada0SJunchao Zhang Input Parameters: 187254729392SStefano Zampini + sfA - The first PetscSF 187354729392SStefano Zampini - sfB - The second PetscSF 187404c0ada0SJunchao Zhang 187504c0ada0SJunchao Zhang Output Parameters: 187654729392SStefano Zampini . sfBA - The composite SF. 187704c0ada0SJunchao Zhang 187804c0ada0SJunchao Zhang Level: developer 187904c0ada0SJunchao Zhang 188054729392SStefano Zampini Notes: 188154729392SStefano Zampini Currently, the two SFs must be defined on congruent communicators and they must be true star 188254729392SStefano Zampini forests, i.e. the same leaf is not connected with different roots. Even more, all roots of the 188354729392SStefano Zampini second SF must have a degree of 1, i.e., no roots have more than one leaf connected. 188454729392SStefano Zampini 188554729392SStefano Zampini sfA's leaf space and sfB's leaf space might be partially overlapped. The composition builds 188654729392SStefano Zampini a graph with sfA's roots and sfB's roots only when there is a path between them. Unconnected 188754729392SStefano Zampini roots are not in sfBA. Doing a Bcast on the new SF is equivalent to doing a Bcast on sfA, then 188854729392SStefano Zampini a Reduce on sfB, on connected roots. 188954729392SStefano Zampini 189054729392SStefano Zampini .seealso: PetscSF, PetscSFCompose(), PetscSFGetGraph(), PetscSFSetGraph(), PetscSFCreateInverseSF() 189104c0ada0SJunchao Zhang @*/ 189204c0ada0SJunchao Zhang PetscErrorCode PetscSFComposeInverse(PetscSF sfA,PetscSF sfB,PetscSF *sfBA) 189304c0ada0SJunchao Zhang { 189404c0ada0SJunchao Zhang PetscErrorCode ierr; 189504c0ada0SJunchao Zhang const PetscSFNode *remotePointsA,*remotePointsB; 189604c0ada0SJunchao Zhang PetscSFNode *remotePointsBA; 189704c0ada0SJunchao Zhang const PetscInt *localPointsA,*localPointsB; 189854729392SStefano Zampini PetscSFNode *reorderedRemotePointsA = NULL; 189954729392SStefano Zampini PetscInt i,numRootsA,numLeavesA,numLeavesBA,numRootsB,numLeavesB,minleaf,maxleaf,*localPointsBA; 19005b0d146aSStefano Zampini MPI_Op op; 19015b0d146aSStefano Zampini #if defined(PETSC_USE_64BIT_INDICES) 19025b0d146aSStefano Zampini PetscBool iswin; 19035b0d146aSStefano Zampini #endif 190404c0ada0SJunchao Zhang 190504c0ada0SJunchao Zhang PetscFunctionBegin; 190604c0ada0SJunchao Zhang PetscValidHeaderSpecific(sfA,PETSCSF_CLASSID,1); 190704c0ada0SJunchao Zhang PetscSFCheckGraphSet(sfA,1); 190804c0ada0SJunchao Zhang PetscValidHeaderSpecific(sfB,PETSCSF_CLASSID,2); 190904c0ada0SJunchao Zhang PetscSFCheckGraphSet(sfB,2); 191054729392SStefano Zampini PetscCheckSameComm(sfA,1,sfB,2); 191104c0ada0SJunchao Zhang PetscValidPointer(sfBA,3); 191254729392SStefano Zampini ierr = PetscSFCheckLeavesUnique_Private(sfA);CHKERRQ(ierr); 191354729392SStefano Zampini ierr = PetscSFCheckLeavesUnique_Private(sfB);CHKERRQ(ierr); 191454729392SStefano Zampini 191504c0ada0SJunchao Zhang ierr = PetscSFGetGraph(sfA, &numRootsA, &numLeavesA, &localPointsA, &remotePointsA);CHKERRQ(ierr); 191604c0ada0SJunchao Zhang ierr = PetscSFGetGraph(sfB, &numRootsB, &numLeavesB, &localPointsB, &remotePointsB);CHKERRQ(ierr); 19175b0d146aSStefano Zampini 19185b0d146aSStefano Zampini /* TODO: Check roots of sfB have degree of 1 */ 19195b0d146aSStefano Zampini /* Once we implement it, we can replace the MPI_MAXLOC 19205b0d146aSStefano Zampini with MPIU_REPLACE. In that case, MPI_MAXLOC and MPIU_REPLACE have the same effect. 19215b0d146aSStefano Zampini We use MPI_MAXLOC only to have a deterministic output from this routine if 19225b0d146aSStefano Zampini the root condition is not meet. 19235b0d146aSStefano Zampini */ 19245b0d146aSStefano Zampini op = MPI_MAXLOC; 19255b0d146aSStefano Zampini #if defined(PETSC_USE_64BIT_INDICES) 19265b0d146aSStefano Zampini /* we accept a non-deterministic output (if any) with PETSCSFWINDOW, since MPI_MAXLOC cannot operate on MPIU_2INT with MPI_Accumulate */ 19275b0d146aSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)sfB,PETSCSFWINDOW,&iswin);CHKERRQ(ierr); 19285b0d146aSStefano Zampini if (iswin) op = MPIU_REPLACE; 19295b0d146aSStefano Zampini #endif 19305b0d146aSStefano Zampini 193154729392SStefano Zampini ierr = PetscSFGetLeafRange(sfB, &minleaf, &maxleaf);CHKERRQ(ierr); 193254729392SStefano Zampini ierr = PetscMalloc1(maxleaf - minleaf + 1,&reorderedRemotePointsA);CHKERRQ(ierr); 193354729392SStefano Zampini for (i=0; i<maxleaf - minleaf + 1; i++) { 193454729392SStefano Zampini reorderedRemotePointsA[i].rank = -1; 193554729392SStefano Zampini reorderedRemotePointsA[i].index = -1; 193654729392SStefano Zampini } 193754729392SStefano Zampini if (localPointsA) { 193854729392SStefano Zampini for (i=0; i<numLeavesA; i++) { 193954729392SStefano Zampini if (localPointsA[i] > maxleaf || localPointsA[i] < minleaf) continue; 194054729392SStefano Zampini reorderedRemotePointsA[localPointsA[i] - minleaf] = remotePointsA[i]; 194154729392SStefano Zampini } 194254729392SStefano Zampini } else { 194354729392SStefano Zampini for (i=0; i<numLeavesA; i++) { 194454729392SStefano Zampini if (i > maxleaf || i < minleaf) continue; 194554729392SStefano Zampini reorderedRemotePointsA[i - minleaf] = remotePointsA[i]; 194654729392SStefano Zampini } 194754729392SStefano Zampini } 194854729392SStefano Zampini 194954729392SStefano Zampini ierr = PetscMalloc1(numRootsB,&localPointsBA);CHKERRQ(ierr); 195004c0ada0SJunchao Zhang ierr = PetscMalloc1(numRootsB,&remotePointsBA);CHKERRQ(ierr); 195154729392SStefano Zampini for (i=0; i<numRootsB; i++) { 195254729392SStefano Zampini remotePointsBA[i].rank = -1; 195354729392SStefano Zampini remotePointsBA[i].index = -1; 195454729392SStefano Zampini } 195554729392SStefano Zampini 19565b0d146aSStefano Zampini ierr = PetscSFReduceBegin(sfB,MPIU_2INT,reorderedRemotePointsA-minleaf,remotePointsBA,op);CHKERRQ(ierr); 19575b0d146aSStefano Zampini ierr = PetscSFReduceEnd(sfB,MPIU_2INT,reorderedRemotePointsA-minleaf,remotePointsBA,op);CHKERRQ(ierr); 195854729392SStefano Zampini ierr = PetscFree(reorderedRemotePointsA);CHKERRQ(ierr); 195954729392SStefano Zampini for (i=0,numLeavesBA=0; i<numRootsB; i++) { 196054729392SStefano Zampini if (remotePointsBA[i].rank == -1) continue; 196154729392SStefano Zampini remotePointsBA[numLeavesBA].rank = remotePointsBA[i].rank; 196254729392SStefano Zampini remotePointsBA[numLeavesBA].index = remotePointsBA[i].index; 196354729392SStefano Zampini localPointsBA[numLeavesBA] = i; 196454729392SStefano Zampini numLeavesBA++; 196554729392SStefano Zampini } 196654729392SStefano Zampini ierr = PetscSFCreate(PetscObjectComm((PetscObject)sfA),sfBA);CHKERRQ(ierr); 196754729392SStefano Zampini ierr = PetscSFSetGraph(*sfBA,numRootsA,numLeavesBA,localPointsBA,PETSC_OWN_POINTER,remotePointsBA,PETSC_OWN_POINTER);CHKERRQ(ierr); 196804c0ada0SJunchao Zhang PetscFunctionReturn(0); 196904c0ada0SJunchao Zhang } 197004c0ada0SJunchao Zhang 19711c6ba672SJunchao Zhang /* 19721c6ba672SJunchao Zhang PetscSFCreateLocalSF_Private - Creates a local PetscSF that only has intra-process edges of the global PetscSF 19731c6ba672SJunchao Zhang 19741c6ba672SJunchao Zhang Input Parameters: 19751c6ba672SJunchao Zhang . sf - The global PetscSF 19761c6ba672SJunchao Zhang 19771c6ba672SJunchao Zhang Output Parameters: 19781c6ba672SJunchao Zhang . out - The local PetscSF 19791c6ba672SJunchao Zhang */ 19801c6ba672SJunchao Zhang PetscErrorCode PetscSFCreateLocalSF_Private(PetscSF sf,PetscSF *out) 19811c6ba672SJunchao Zhang { 19821c6ba672SJunchao Zhang MPI_Comm comm; 19831c6ba672SJunchao Zhang PetscMPIInt myrank; 19841c6ba672SJunchao Zhang const PetscInt *ilocal; 19851c6ba672SJunchao Zhang const PetscSFNode *iremote; 19861c6ba672SJunchao Zhang PetscInt i,j,nroots,nleaves,lnleaves,*lilocal; 19871c6ba672SJunchao Zhang PetscSFNode *liremote; 19881c6ba672SJunchao Zhang PetscSF lsf; 19891c6ba672SJunchao Zhang PetscErrorCode ierr; 19901c6ba672SJunchao Zhang 19911c6ba672SJunchao Zhang PetscFunctionBegin; 19921c6ba672SJunchao Zhang PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 19931c6ba672SJunchao Zhang if (sf->ops->CreateLocalSF) { 19941c6ba672SJunchao Zhang ierr = (*sf->ops->CreateLocalSF)(sf,out);CHKERRQ(ierr); 19951c6ba672SJunchao Zhang } else { 19961c6ba672SJunchao Zhang /* Could use PetscSFCreateEmbeddedLeafSF, but since we know the comm is PETSC_COMM_SELF, we can make it fast */ 19971c6ba672SJunchao Zhang ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr); 19981c6ba672SJunchao Zhang ierr = MPI_Comm_rank(comm,&myrank);CHKERRQ(ierr); 19991c6ba672SJunchao Zhang 20001c6ba672SJunchao Zhang /* Find out local edges and build a local SF */ 20011c6ba672SJunchao Zhang ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr); 20021c6ba672SJunchao Zhang for (i=lnleaves=0; i<nleaves; i++) {if (iremote[i].rank == (PetscInt)myrank) lnleaves++;} 20031c6ba672SJunchao Zhang ierr = PetscMalloc1(lnleaves,&lilocal);CHKERRQ(ierr); 20041c6ba672SJunchao Zhang ierr = PetscMalloc1(lnleaves,&liremote);CHKERRQ(ierr); 20051c6ba672SJunchao Zhang 20061c6ba672SJunchao Zhang for (i=j=0; i<nleaves; i++) { 20071c6ba672SJunchao Zhang if (iremote[i].rank == (PetscInt)myrank) { 20081c6ba672SJunchao Zhang lilocal[j] = ilocal? ilocal[i] : i; /* ilocal=NULL for contiguous storage */ 20091c6ba672SJunchao Zhang liremote[j].rank = 0; /* rank in PETSC_COMM_SELF */ 20101c6ba672SJunchao Zhang liremote[j].index = iremote[i].index; 20111c6ba672SJunchao Zhang j++; 20121c6ba672SJunchao Zhang } 20131c6ba672SJunchao Zhang } 20141c6ba672SJunchao Zhang ierr = PetscSFCreate(PETSC_COMM_SELF,&lsf);CHKERRQ(ierr); 20151c6ba672SJunchao Zhang ierr = PetscSFSetGraph(lsf,nroots,lnleaves,lilocal,PETSC_OWN_POINTER,liremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 20161c6ba672SJunchao Zhang ierr = PetscSFSetUp(lsf);CHKERRQ(ierr); 20171c6ba672SJunchao Zhang *out = lsf; 20181c6ba672SJunchao Zhang } 20191c6ba672SJunchao Zhang PetscFunctionReturn(0); 20201c6ba672SJunchao Zhang } 2021dd5b3ca6SJunchao Zhang 2022dd5b3ca6SJunchao Zhang /* Similar to PetscSFBcast, but only Bcast to leaves on rank 0 */ 2023dd5b3ca6SJunchao Zhang PetscErrorCode PetscSFBcastToZero_Private(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata) 2024dd5b3ca6SJunchao Zhang { 2025dd5b3ca6SJunchao Zhang PetscErrorCode ierr; 2026eb02082bSJunchao Zhang PetscMemType rootmtype,leafmtype; 2027dd5b3ca6SJunchao Zhang 2028dd5b3ca6SJunchao Zhang PetscFunctionBegin; 2029dd5b3ca6SJunchao Zhang PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 2030dd5b3ca6SJunchao Zhang ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 2031dd5b3ca6SJunchao Zhang ierr = PetscLogEventBegin(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr); 2032eb02082bSJunchao Zhang ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr); 2033eb02082bSJunchao Zhang ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr); 2034dd5b3ca6SJunchao Zhang if (sf->ops->BcastToZero) { 2035eb02082bSJunchao Zhang ierr = (*sf->ops->BcastToZero)(sf,unit,rootmtype,rootdata,leafmtype,leafdata);CHKERRQ(ierr); 2036dd5b3ca6SJunchao Zhang } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSFBcastToZero_Private is not supported on this SF type"); 2037dd5b3ca6SJunchao Zhang ierr = PetscLogEventEnd(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr); 2038dd5b3ca6SJunchao Zhang PetscFunctionReturn(0); 2039dd5b3ca6SJunchao Zhang } 2040dd5b3ca6SJunchao Zhang 2041