1af0996ceSBarry Smith #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/ 295fce210SBarry Smith #include <petscctable.h> 395fce210SBarry Smith 495fce210SBarry Smith #if defined(PETSC_USE_DEBUG) 595fce210SBarry Smith # define PetscSFCheckGraphSet(sf,arg) do { \ 695fce210SBarry Smith if (PetscUnlikely(!(sf)->graphset)) \ 795fce210SBarry Smith SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetGraph() on argument %D \"%s\" before %s()",(arg),#sf,PETSC_FUNCTION_NAME); \ 895fce210SBarry Smith } while (0) 995fce210SBarry Smith #else 1095fce210SBarry Smith # define PetscSFCheckGraphSet(sf,arg) do {} while (0) 1195fce210SBarry Smith #endif 1295fce210SBarry Smith 1395fce210SBarry Smith const char *const PetscSFDuplicateOptions[] = {"CONFONLY","RANKS","GRAPH","PetscSFDuplicateOption","PETSCSF_DUPLICATE_",0}; 1495fce210SBarry Smith 158af6ec1cSBarry Smith /*@ 1695fce210SBarry Smith PetscSFCreate - create a star forest communication context 1795fce210SBarry Smith 1895fce210SBarry Smith Not Collective 1995fce210SBarry Smith 2095fce210SBarry Smith Input Arguments: 2195fce210SBarry Smith . comm - communicator on which the star forest will operate 2295fce210SBarry Smith 2395fce210SBarry Smith Output Arguments: 2495fce210SBarry Smith . sf - new star forest context 2595fce210SBarry Smith 2695fce210SBarry Smith Level: intermediate 2795fce210SBarry Smith 2895fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFDestroy() 2995fce210SBarry Smith @*/ 3095fce210SBarry Smith PetscErrorCode PetscSFCreate(MPI_Comm comm,PetscSF *sf) 3195fce210SBarry Smith { 3295fce210SBarry Smith PetscErrorCode ierr; 3395fce210SBarry Smith PetscSF b; 3495fce210SBarry Smith 3595fce210SBarry Smith PetscFunctionBegin; 3695fce210SBarry Smith PetscValidPointer(sf,2); 37607a6623SBarry Smith ierr = PetscSFInitializePackage();CHKERRQ(ierr); 3895fce210SBarry Smith 3973107ff1SLisandro Dalcin ierr = PetscHeaderCreate(b,PETSCSF_CLASSID,"PetscSF","Star Forest","PetscSF",comm,PetscSFDestroy,PetscSFView);CHKERRQ(ierr); 4095fce210SBarry Smith 4195fce210SBarry Smith b->nroots = -1; 4295fce210SBarry Smith b->nleaves = -1; 43*29046d53SLisandro Dalcin b->minleaf = PETSC_MAX_INT; 44*29046d53SLisandro Dalcin b->maxleaf = PETSC_MIN_INT; 4595fce210SBarry Smith b->nranks = -1; 4695fce210SBarry Smith b->rankorder = PETSC_TRUE; 4795fce210SBarry Smith b->ingroup = MPI_GROUP_NULL; 4895fce210SBarry Smith b->outgroup = MPI_GROUP_NULL; 4995fce210SBarry Smith b->graphset = PETSC_FALSE; 5095fce210SBarry Smith 5195fce210SBarry Smith *sf = b; 5295fce210SBarry Smith PetscFunctionReturn(0); 5395fce210SBarry Smith } 5495fce210SBarry Smith 55*29046d53SLisandro Dalcin /*@ 5695fce210SBarry Smith PetscSFReset - Reset a star forest so that different sizes or neighbors can be used 5795fce210SBarry Smith 5895fce210SBarry Smith Collective 5995fce210SBarry Smith 6095fce210SBarry Smith Input Arguments: 6195fce210SBarry Smith . sf - star forest 6295fce210SBarry Smith 6395fce210SBarry Smith Level: advanced 6495fce210SBarry Smith 6595fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFDestroy() 6695fce210SBarry Smith @*/ 6795fce210SBarry Smith PetscErrorCode PetscSFReset(PetscSF sf) 6895fce210SBarry Smith { 6995fce210SBarry Smith PetscErrorCode ierr; 7095fce210SBarry Smith 7195fce210SBarry Smith PetscFunctionBegin; 7295fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 7379715d56SJed Brown if (sf->ops->Reset) {ierr = (*sf->ops->Reset)(sf);CHKERRQ(ierr);} 74*29046d53SLisandro Dalcin sf->nroots = -1; 75*29046d53SLisandro Dalcin sf->nleaves = -1; 76*29046d53SLisandro Dalcin sf->minleaf = PETSC_MAX_INT; 77*29046d53SLisandro Dalcin sf->maxleaf = PETSC_MIN_INT; 7895fce210SBarry Smith sf->mine = NULL; 7995fce210SBarry Smith sf->remote = NULL; 80*29046d53SLisandro Dalcin sf->graphset = PETSC_FALSE; 81*29046d53SLisandro Dalcin ierr = PetscFree(sf->mine_alloc);CHKERRQ(ierr); 8295fce210SBarry Smith ierr = PetscFree(sf->remote_alloc);CHKERRQ(ierr); 8321c688dcSJed Brown sf->nranks = -1; 84*29046d53SLisandro Dalcin ierr = PetscFree4(sf->ranks,sf->roffset,sf->rmine,sf->rremote);CHKERRQ(ierr); 85*29046d53SLisandro Dalcin sf->degreeknown = PETSC_FALSE; 8695fce210SBarry Smith ierr = PetscFree(sf->degree);CHKERRQ(ierr); 8795fce210SBarry Smith if (sf->ingroup != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->ingroup);CHKERRQ(ierr);} 8895fce210SBarry Smith if (sf->outgroup != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->outgroup);CHKERRQ(ierr);} 8995fce210SBarry Smith ierr = PetscSFDestroy(&sf->multi);CHKERRQ(ierr); 9095fce210SBarry Smith sf->setupcalled = PETSC_FALSE; 9195fce210SBarry Smith PetscFunctionReturn(0); 9295fce210SBarry Smith } 9395fce210SBarry Smith 9495fce210SBarry Smith /*@C 95*29046d53SLisandro Dalcin PetscSFSetType - Set the PetscSF communication implementation 9695fce210SBarry Smith 9795fce210SBarry Smith Collective on PetscSF 9895fce210SBarry Smith 9995fce210SBarry Smith Input Parameters: 10095fce210SBarry Smith + sf - the PetscSF context 10195fce210SBarry Smith - type - a known method 10295fce210SBarry Smith 10395fce210SBarry Smith Options Database Key: 10495fce210SBarry Smith . -sf_type <type> - Sets the method; use -help for a list 10595fce210SBarry Smith of available methods (for instance, window, pt2pt, neighbor) 10695fce210SBarry Smith 10795fce210SBarry Smith Notes: 10895fce210SBarry Smith See "include/petscsf.h" for available methods (for instance) 10995fce210SBarry Smith + PETSCSFWINDOW - MPI-2/3 one-sided 11095fce210SBarry Smith - PETSCSFBASIC - basic implementation using MPI-1 two-sided 11195fce210SBarry Smith 11295fce210SBarry Smith Level: intermediate 11395fce210SBarry Smith 11495fce210SBarry Smith .keywords: PetscSF, set, type 11595fce210SBarry Smith 11695fce210SBarry Smith .seealso: PetscSFType, PetscSFCreate() 11795fce210SBarry Smith @*/ 11895fce210SBarry Smith PetscErrorCode PetscSFSetType(PetscSF sf,PetscSFType type) 11995fce210SBarry Smith { 12095fce210SBarry Smith PetscErrorCode ierr,(*r)(PetscSF); 12195fce210SBarry Smith PetscBool match; 12295fce210SBarry Smith 12395fce210SBarry Smith PetscFunctionBegin; 12495fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 12595fce210SBarry Smith PetscValidCharPointer(type,2); 12695fce210SBarry Smith 12795fce210SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)sf,type,&match);CHKERRQ(ierr); 12895fce210SBarry Smith if (match) PetscFunctionReturn(0); 12995fce210SBarry Smith 130adc40e5bSBarry Smith ierr = PetscFunctionListFind(PetscSFList,type,&r);CHKERRQ(ierr); 13195fce210SBarry Smith if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested PetscSF type %s",type); 132*29046d53SLisandro Dalcin /* Destroy the previous PetscSF implementation context */ 133*29046d53SLisandro Dalcin if (sf->ops->Destroy) {ierr = (*(sf)->ops->Destroy)(sf);CHKERRQ(ierr);} 13495fce210SBarry Smith ierr = PetscMemzero(sf->ops,sizeof(*sf->ops));CHKERRQ(ierr); 13595fce210SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)sf,type);CHKERRQ(ierr); 13695fce210SBarry Smith ierr = (*r)(sf);CHKERRQ(ierr); 13795fce210SBarry Smith PetscFunctionReturn(0); 13895fce210SBarry Smith } 13995fce210SBarry Smith 140*29046d53SLisandro Dalcin /*@C 141*29046d53SLisandro Dalcin PetscSFGetType - Get the PetscSF communication implementation 142*29046d53SLisandro Dalcin 143*29046d53SLisandro Dalcin Not Collective 144*29046d53SLisandro Dalcin 145*29046d53SLisandro Dalcin Input Parameter: 146*29046d53SLisandro Dalcin . sf - the PetscSF context 147*29046d53SLisandro Dalcin 148*29046d53SLisandro Dalcin Output Parameter: 149*29046d53SLisandro Dalcin . type - the PetscSF type name 150*29046d53SLisandro Dalcin 151*29046d53SLisandro Dalcin Level: intermediate 152*29046d53SLisandro Dalcin 153*29046d53SLisandro Dalcin .keywords: PetscSF, get, type 154*29046d53SLisandro Dalcin .seealso: PetscSFSetType(), PetscSFCreate() 155*29046d53SLisandro Dalcin @*/ 156*29046d53SLisandro Dalcin PetscErrorCode PetscSFGetType(PetscSF sf, PetscSFType *type) 157*29046d53SLisandro Dalcin { 158*29046d53SLisandro Dalcin PetscFunctionBegin; 159*29046d53SLisandro Dalcin PetscValidHeaderSpecific(sf, PETSCSF_CLASSID,1); 160*29046d53SLisandro Dalcin PetscValidPointer(type,2); 161*29046d53SLisandro Dalcin *type = ((PetscObject)sf)->type_name; 162*29046d53SLisandro Dalcin PetscFunctionReturn(0); 163*29046d53SLisandro Dalcin } 164*29046d53SLisandro Dalcin 165d36d33e4SMatthew G. Knepley /*@ 16695fce210SBarry Smith PetscSFDestroy - destroy star forest 16795fce210SBarry Smith 16895fce210SBarry Smith Collective 16995fce210SBarry Smith 17095fce210SBarry Smith Input Arguments: 17195fce210SBarry Smith . sf - address of star forest 17295fce210SBarry Smith 17395fce210SBarry Smith Level: intermediate 17495fce210SBarry Smith 17595fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFReset() 17695fce210SBarry Smith @*/ 17795fce210SBarry Smith PetscErrorCode PetscSFDestroy(PetscSF *sf) 17895fce210SBarry Smith { 17995fce210SBarry Smith PetscErrorCode ierr; 18095fce210SBarry Smith 18195fce210SBarry Smith PetscFunctionBegin; 18295fce210SBarry Smith if (!*sf) PetscFunctionReturn(0); 18395fce210SBarry Smith PetscValidHeaderSpecific((*sf),PETSCSF_CLASSID,1); 184*29046d53SLisandro Dalcin if (--((PetscObject)(*sf))->refct > 0) {*sf = NULL; PetscFunctionReturn(0);} 18595fce210SBarry Smith ierr = PetscSFReset(*sf);CHKERRQ(ierr); 18695fce210SBarry Smith if ((*sf)->ops->Destroy) {ierr = (*(*sf)->ops->Destroy)(*sf);CHKERRQ(ierr);} 18795fce210SBarry Smith ierr = PetscHeaderDestroy(sf);CHKERRQ(ierr); 18895fce210SBarry Smith PetscFunctionReturn(0); 18995fce210SBarry Smith } 19095fce210SBarry Smith 19195fce210SBarry Smith /*@ 19295fce210SBarry Smith PetscSFSetUp - set up communication structures 19395fce210SBarry Smith 19495fce210SBarry Smith Collective 19595fce210SBarry Smith 19695fce210SBarry Smith Input Arguments: 19795fce210SBarry Smith . sf - star forest communication object 19895fce210SBarry Smith 19995fce210SBarry Smith Level: beginner 20095fce210SBarry Smith 20195fce210SBarry Smith .seealso: PetscSFSetFromOptions(), PetscSFSetType() 20295fce210SBarry Smith @*/ 20395fce210SBarry Smith PetscErrorCode PetscSFSetUp(PetscSF sf) 20495fce210SBarry Smith { 20595fce210SBarry Smith PetscErrorCode ierr; 20695fce210SBarry Smith 20795fce210SBarry Smith PetscFunctionBegin; 208*29046d53SLisandro Dalcin PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 209*29046d53SLisandro Dalcin PetscSFCheckGraphSet(sf,1); 21095fce210SBarry Smith if (sf->setupcalled) PetscFunctionReturn(0); 21195fce210SBarry Smith if (!((PetscObject)sf)->type_name) {ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr);} 212*29046d53SLisandro Dalcin ierr = PetscLogEventBegin(PETSCSF_SetUp,sf,0,0,0);CHKERRQ(ierr); 21395fce210SBarry Smith if (sf->ops->SetUp) {ierr = (*sf->ops->SetUp)(sf);CHKERRQ(ierr);} 214*29046d53SLisandro Dalcin ierr = PetscLogEventEnd(PETSCSF_SetUp,sf,0,0,0);CHKERRQ(ierr); 21595fce210SBarry Smith sf->setupcalled = PETSC_TRUE; 21695fce210SBarry Smith PetscFunctionReturn(0); 21795fce210SBarry Smith } 21895fce210SBarry Smith 2198af6ec1cSBarry Smith /*@ 22095fce210SBarry Smith PetscSFSetFromOptions - set PetscSF options using the options database 22195fce210SBarry Smith 22295fce210SBarry Smith Logically Collective 22395fce210SBarry Smith 22495fce210SBarry Smith Input Arguments: 22595fce210SBarry Smith . sf - star forest 22695fce210SBarry Smith 22795fce210SBarry Smith Options Database Keys: 22860263706SJed Brown + -sf_type - implementation type, see PetscSFSetType() 22960263706SJed Brown - -sf_rank_order - sort composite points for gathers and scatters in rank order, gathers are non-deterministic otherwise 23095fce210SBarry Smith 23195fce210SBarry Smith Level: intermediate 23295fce210SBarry Smith 23395fce210SBarry Smith .keywords: KSP, set, from, options, database 23495fce210SBarry Smith 23595fce210SBarry Smith .seealso: PetscSFWindowSetSyncType() 23695fce210SBarry Smith @*/ 23795fce210SBarry Smith PetscErrorCode PetscSFSetFromOptions(PetscSF sf) 23895fce210SBarry Smith { 23995fce210SBarry Smith PetscSFType deft; 24095fce210SBarry Smith char type[256]; 24195fce210SBarry Smith PetscErrorCode ierr; 24295fce210SBarry Smith PetscBool flg; 24395fce210SBarry Smith 24495fce210SBarry Smith PetscFunctionBegin; 24595fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 24695fce210SBarry Smith ierr = PetscObjectOptionsBegin((PetscObject)sf);CHKERRQ(ierr); 24795fce210SBarry Smith deft = ((PetscObject)sf)->type_name ? ((PetscObject)sf)->type_name : PETSCSFBASIC; 248*29046d53SLisandro Dalcin ierr = PetscOptionsFList("-sf_type","PetscSF implementation type","PetscSFSetType",PetscSFList,deft,type,sizeof(type),&flg);CHKERRQ(ierr); 24995fce210SBarry Smith ierr = PetscSFSetType(sf,flg ? type : deft);CHKERRQ(ierr); 25095fce210SBarry Smith ierr = PetscOptionsBool("-sf_rank_order","sort composite points for gathers and scatters in rank order, gathers are non-deterministic otherwise","PetscSFSetRankOrder",sf->rankorder,&sf->rankorder,NULL);CHKERRQ(ierr); 251e55864a3SBarry Smith if (sf->ops->SetFromOptions) {ierr = (*sf->ops->SetFromOptions)(PetscOptionsObject,sf);CHKERRQ(ierr);} 25295fce210SBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 25395fce210SBarry Smith PetscFunctionReturn(0); 25495fce210SBarry Smith } 25595fce210SBarry Smith 256*29046d53SLisandro Dalcin /*@ 25795fce210SBarry Smith PetscSFSetRankOrder - sort multi-points for gathers and scatters by rank order 25895fce210SBarry Smith 25995fce210SBarry Smith Logically Collective 26095fce210SBarry Smith 26195fce210SBarry Smith Input Arguments: 26295fce210SBarry Smith + sf - star forest 26395fce210SBarry Smith - flg - PETSC_TRUE to sort, PETSC_FALSE to skip sorting (lower setup cost, but non-deterministic) 26495fce210SBarry Smith 26595fce210SBarry Smith Level: advanced 26695fce210SBarry Smith 26795fce210SBarry Smith .seealso: PetscSFGatherBegin(), PetscSFScatterBegin() 26895fce210SBarry Smith @*/ 26995fce210SBarry Smith PetscErrorCode PetscSFSetRankOrder(PetscSF sf,PetscBool flg) 27095fce210SBarry Smith { 27195fce210SBarry Smith 27295fce210SBarry Smith PetscFunctionBegin; 27395fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 27495fce210SBarry Smith PetscValidLogicalCollectiveBool(sf,flg,2); 27595fce210SBarry Smith if (sf->multi) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Rank ordering must be set before first call to PetscSFGatherBegin() or PetscSFScatterBegin()"); 27695fce210SBarry Smith sf->rankorder = flg; 27795fce210SBarry Smith PetscFunctionReturn(0); 27895fce210SBarry Smith } 27995fce210SBarry Smith 2808af6ec1cSBarry Smith /*@ 28195fce210SBarry Smith PetscSFSetGraph - Set a parallel star forest 28295fce210SBarry Smith 28395fce210SBarry Smith Collective 28495fce210SBarry Smith 28595fce210SBarry Smith Input Arguments: 28695fce210SBarry Smith + sf - star forest 28795fce210SBarry Smith . nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves) 28895fce210SBarry Smith . nleaves - number of leaf vertices on the current process, each of these references a root on any process 28995fce210SBarry Smith . ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage 29095fce210SBarry Smith . localmode - copy mode for ilocal 29195fce210SBarry Smith . iremote - remote locations of root vertices for each leaf on the current process 29295fce210SBarry Smith - remotemode - copy mode for iremote 29395fce210SBarry Smith 29495fce210SBarry Smith Level: intermediate 29595fce210SBarry Smith 29638ab3f8aSBarry Smith Notes: In Fortran you must use PETSC_COPY_VALUES for localmode and remotemode 29738ab3f8aSBarry Smith 29895fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFGetGraph() 29995fce210SBarry Smith @*/ 30095fce210SBarry Smith PetscErrorCode PetscSFSetGraph(PetscSF sf,PetscInt nroots,PetscInt nleaves,const PetscInt *ilocal,PetscCopyMode localmode,const PetscSFNode *iremote,PetscCopyMode remotemode) 30195fce210SBarry Smith { 30295fce210SBarry Smith PetscErrorCode ierr; 30395fce210SBarry Smith 30495fce210SBarry Smith PetscFunctionBegin; 30595fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 306*29046d53SLisandro Dalcin if (nleaves > 0 && ilocal) PetscValidIntPointer(ilocal,4); 307*29046d53SLisandro Dalcin if (nleaves > 0) PetscValidPointer(iremote,6); 308*29046d53SLisandro Dalcin if (nroots < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nroots %D, cannot be negative",nroots); 30995fce210SBarry Smith if (nleaves < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nleaves %D, cannot be negative",nleaves); 310*29046d53SLisandro Dalcin 31195fce210SBarry Smith ierr = PetscSFReset(sf);CHKERRQ(ierr); 312*29046d53SLisandro Dalcin ierr = PetscLogEventBegin(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr); 313*29046d53SLisandro Dalcin 31495fce210SBarry Smith sf->nroots = nroots; 31595fce210SBarry Smith sf->nleaves = nleaves; 316*29046d53SLisandro Dalcin 317*29046d53SLisandro Dalcin if (nleaves && ilocal) { 31821c688dcSJed Brown PetscInt i; 319*29046d53SLisandro Dalcin PetscInt minleaf = PETSC_MAX_INT; 320*29046d53SLisandro Dalcin PetscInt maxleaf = PETSC_MIN_INT; 321*29046d53SLisandro Dalcin for (i=0; i<nleaves; i++) { 322*29046d53SLisandro Dalcin minleaf = PetscMin(minleaf,ilocal[i]); 323*29046d53SLisandro Dalcin maxleaf = PetscMax(maxleaf,ilocal[i]); 324*29046d53SLisandro Dalcin } 325*29046d53SLisandro Dalcin sf->minleaf = minleaf; 326*29046d53SLisandro Dalcin sf->maxleaf = maxleaf; 327*29046d53SLisandro Dalcin } else { 328*29046d53SLisandro Dalcin sf->minleaf = 0; 329*29046d53SLisandro Dalcin sf->maxleaf = nleaves - 1; 330*29046d53SLisandro Dalcin } 331*29046d53SLisandro Dalcin 332*29046d53SLisandro Dalcin if (ilocal) { 33395fce210SBarry Smith switch (localmode) { 33495fce210SBarry Smith case PETSC_COPY_VALUES: 335785e854fSJed Brown ierr = PetscMalloc1(nleaves,&sf->mine_alloc);CHKERRQ(ierr); 336*29046d53SLisandro Dalcin ierr = PetscMemcpy(sf->mine_alloc,ilocal,nleaves*sizeof(*ilocal));CHKERRQ(ierr); 33795fce210SBarry Smith sf->mine = sf->mine_alloc; 33895fce210SBarry Smith break; 33995fce210SBarry Smith case PETSC_OWN_POINTER: 34095fce210SBarry Smith sf->mine_alloc = (PetscInt*)ilocal; 34195fce210SBarry Smith sf->mine = sf->mine_alloc; 34295fce210SBarry Smith break; 34395fce210SBarry Smith case PETSC_USE_POINTER: 344*29046d53SLisandro Dalcin sf->mine_alloc = NULL; 34595fce210SBarry Smith sf->mine = (PetscInt*)ilocal; 34695fce210SBarry Smith break; 34795fce210SBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown localmode"); 34895fce210SBarry Smith } 34995fce210SBarry Smith } 350*29046d53SLisandro Dalcin 35195fce210SBarry Smith switch (remotemode) { 35295fce210SBarry Smith case PETSC_COPY_VALUES: 353785e854fSJed Brown ierr = PetscMalloc1(nleaves,&sf->remote_alloc);CHKERRQ(ierr); 354*29046d53SLisandro Dalcin ierr = PetscMemcpy(sf->remote_alloc,iremote,nleaves*sizeof(*iremote));CHKERRQ(ierr); 35595fce210SBarry Smith sf->remote = sf->remote_alloc; 35695fce210SBarry Smith break; 35795fce210SBarry Smith case PETSC_OWN_POINTER: 35895fce210SBarry Smith sf->remote_alloc = (PetscSFNode*)iremote; 35995fce210SBarry Smith sf->remote = sf->remote_alloc; 36095fce210SBarry Smith break; 36195fce210SBarry Smith case PETSC_USE_POINTER: 362*29046d53SLisandro Dalcin sf->remote_alloc = NULL; 36395fce210SBarry Smith sf->remote = (PetscSFNode*)iremote; 36495fce210SBarry Smith break; 36595fce210SBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown remotemode"); 36695fce210SBarry Smith } 36795fce210SBarry Smith 368563ffbabSMatthew G. Knepley ierr = PetscLogEventEnd(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr); 369*29046d53SLisandro Dalcin sf->graphset = PETSC_TRUE; 37095fce210SBarry Smith PetscFunctionReturn(0); 37195fce210SBarry Smith } 37295fce210SBarry Smith 373*29046d53SLisandro Dalcin /*@ 37495fce210SBarry Smith PetscSFCreateInverseSF - given a PetscSF in which all vertices have degree 1, creates the inverse map 37595fce210SBarry Smith 37695fce210SBarry Smith Collective 37795fce210SBarry Smith 37895fce210SBarry Smith Input Arguments: 37995fce210SBarry Smith . sf - star forest to invert 38095fce210SBarry Smith 38195fce210SBarry Smith Output Arguments: 38295fce210SBarry Smith . isf - inverse of sf 38395fce210SBarry Smith 38495fce210SBarry Smith Level: advanced 38595fce210SBarry Smith 38695fce210SBarry Smith Notes: 38795fce210SBarry Smith All roots must have degree 1. 38895fce210SBarry Smith 38995fce210SBarry Smith The local space may be a permutation, but cannot be sparse. 39095fce210SBarry Smith 39195fce210SBarry Smith .seealso: PetscSFSetGraph() 39295fce210SBarry Smith @*/ 39395fce210SBarry Smith PetscErrorCode PetscSFCreateInverseSF(PetscSF sf,PetscSF *isf) 39495fce210SBarry Smith { 39595fce210SBarry Smith PetscErrorCode ierr; 39695fce210SBarry Smith PetscMPIInt rank; 39795fce210SBarry Smith PetscInt i,nroots,nleaves,maxlocal,count,*newilocal; 39895fce210SBarry Smith const PetscInt *ilocal; 39995fce210SBarry Smith PetscSFNode *roots,*leaves; 40095fce210SBarry Smith 40195fce210SBarry Smith PetscFunctionBegin; 402*29046d53SLisandro Dalcin PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 403*29046d53SLisandro Dalcin PetscSFCheckGraphSet(sf,1); 404*29046d53SLisandro Dalcin PetscValidPointer(isf,2); 405*29046d53SLisandro Dalcin 40695fce210SBarry Smith ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);CHKERRQ(ierr); 407*29046d53SLisandro Dalcin maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */ 408*29046d53SLisandro Dalcin 409*29046d53SLisandro Dalcin ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); 410ae9aee6dSMatthew G. Knepley ierr = PetscMalloc2(nroots,&roots,maxlocal,&leaves);CHKERRQ(ierr); 411ae9aee6dSMatthew G. Knepley for (i=0; i<maxlocal; i++) { 41295fce210SBarry Smith leaves[i].rank = rank; 41395fce210SBarry Smith leaves[i].index = i; 41495fce210SBarry Smith } 41595fce210SBarry Smith for (i=0; i <nroots; i++) { 41695fce210SBarry Smith roots[i].rank = -1; 41795fce210SBarry Smith roots[i].index = -1; 41895fce210SBarry Smith } 4198bfbc91cSJed Brown ierr = PetscSFReduceBegin(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr); 4208bfbc91cSJed Brown ierr = PetscSFReduceEnd(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr); 42195fce210SBarry Smith 42295fce210SBarry Smith /* Check whether our leaves are sparse */ 42395fce210SBarry Smith for (i=0,count=0; i<nroots; i++) if (roots[i].rank >= 0) count++; 42495fce210SBarry Smith if (count == nroots) newilocal = NULL; 42595fce210SBarry Smith else { /* Index for sparse leaves and compact "roots" array (which is to become our leaves). */ 426785e854fSJed Brown ierr = PetscMalloc1(count,&newilocal);CHKERRQ(ierr); 42795fce210SBarry Smith for (i=0,count=0; i<nroots; i++) { 42895fce210SBarry Smith if (roots[i].rank >= 0) { 42995fce210SBarry Smith newilocal[count] = i; 43095fce210SBarry Smith roots[count].rank = roots[i].rank; 43195fce210SBarry Smith roots[count].index = roots[i].index; 43295fce210SBarry Smith count++; 43395fce210SBarry Smith } 43495fce210SBarry Smith } 43595fce210SBarry Smith } 43695fce210SBarry Smith 43795fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,isf);CHKERRQ(ierr); 43895fce210SBarry Smith ierr = PetscSFSetGraph(*isf,maxlocal,count,newilocal,PETSC_OWN_POINTER,roots,PETSC_COPY_VALUES);CHKERRQ(ierr); 43995fce210SBarry Smith ierr = PetscFree2(roots,leaves);CHKERRQ(ierr); 44095fce210SBarry Smith PetscFunctionReturn(0); 44195fce210SBarry Smith } 44295fce210SBarry Smith 44395fce210SBarry Smith /*@ 44495fce210SBarry Smith PetscSFDuplicate - duplicate a PetscSF, optionally preserving rank connectivity and graph 44595fce210SBarry Smith 44695fce210SBarry Smith Collective 44795fce210SBarry Smith 44895fce210SBarry Smith Input Arguments: 44995fce210SBarry Smith + sf - communication object to duplicate 45095fce210SBarry Smith - opt - PETSCSF_DUPLICATE_CONFONLY, PETSCSF_DUPLICATE_RANKS, or PETSCSF_DUPLICATE_GRAPH (see PetscSFDuplicateOption) 45195fce210SBarry Smith 45295fce210SBarry Smith Output Arguments: 45395fce210SBarry Smith . newsf - new communication object 45495fce210SBarry Smith 45595fce210SBarry Smith Level: beginner 45695fce210SBarry Smith 45795fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetType(), PetscSFSetGraph() 45895fce210SBarry Smith @*/ 45995fce210SBarry Smith PetscErrorCode PetscSFDuplicate(PetscSF sf,PetscSFDuplicateOption opt,PetscSF *newsf) 46095fce210SBarry Smith { 461*29046d53SLisandro Dalcin PetscSFType type; 46295fce210SBarry Smith PetscErrorCode ierr; 46395fce210SBarry Smith 46495fce210SBarry Smith PetscFunctionBegin; 465*29046d53SLisandro Dalcin PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 466*29046d53SLisandro Dalcin PetscValidLogicalCollectiveEnum(sf,opt,2); 467*29046d53SLisandro Dalcin PetscValidPointer(newsf,3); 46895fce210SBarry Smith ierr = PetscSFCreate(PetscObjectComm((PetscObject)sf),newsf);CHKERRQ(ierr); 469*29046d53SLisandro Dalcin ierr = PetscSFGetType(sf,&type);CHKERRQ(ierr); 470*29046d53SLisandro Dalcin if (type) {ierr = PetscSFSetType(*newsf,type);CHKERRQ(ierr);} 47195fce210SBarry Smith if (opt == PETSCSF_DUPLICATE_GRAPH) { 47295fce210SBarry Smith PetscInt nroots,nleaves; 47395fce210SBarry Smith const PetscInt *ilocal; 47495fce210SBarry Smith const PetscSFNode *iremote; 475*29046d53SLisandro Dalcin PetscSFCheckGraphSet(sf,1); 47695fce210SBarry Smith ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr); 47795fce210SBarry Smith ierr = PetscSFSetGraph(*newsf,nroots,nleaves,ilocal,PETSC_COPY_VALUES,iremote,PETSC_COPY_VALUES);CHKERRQ(ierr); 47895fce210SBarry Smith } 479*29046d53SLisandro Dalcin if (sf->ops->Duplicate) {ierr = (*sf->ops->Duplicate)(sf,opt,*newsf);CHKERRQ(ierr);} 48095fce210SBarry Smith PetscFunctionReturn(0); 48195fce210SBarry Smith } 48295fce210SBarry Smith 48395fce210SBarry Smith /*@C 48495fce210SBarry Smith PetscSFGetGraph - Get the graph specifying a parallel star forest 48595fce210SBarry Smith 48695fce210SBarry Smith Not Collective 48795fce210SBarry Smith 48895fce210SBarry Smith Input Arguments: 48995fce210SBarry Smith . sf - star forest 49095fce210SBarry Smith 49195fce210SBarry Smith Output Arguments: 49295fce210SBarry Smith + nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves) 49395fce210SBarry Smith . nleaves - number of leaf vertices on the current process, each of these references a root on any process 49495fce210SBarry Smith . ilocal - locations of leaves in leafdata buffers 49595fce210SBarry Smith - iremote - remote locations of root vertices for each leaf on the current process 49695fce210SBarry Smith 49795fce210SBarry Smith Level: intermediate 49895fce210SBarry Smith 49995fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph() 50095fce210SBarry Smith @*/ 50195fce210SBarry Smith PetscErrorCode PetscSFGetGraph(PetscSF sf,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote) 50295fce210SBarry Smith { 50395fce210SBarry Smith 50495fce210SBarry Smith PetscFunctionBegin; 50595fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 50695fce210SBarry Smith /* We are not currently requiring that the graph is set, thus returning nroots=-1 if it has not been set */ 507*29046d53SLisandro Dalcin /* PetscSFCheckGraphSet(sf,1); */ 50895fce210SBarry Smith if (nroots) *nroots = sf->nroots; 50995fce210SBarry Smith if (nleaves) *nleaves = sf->nleaves; 51095fce210SBarry Smith if (ilocal) *ilocal = sf->mine; 51195fce210SBarry Smith if (iremote) *iremote = sf->remote; 51295fce210SBarry Smith PetscFunctionReturn(0); 51395fce210SBarry Smith } 51495fce210SBarry Smith 515*29046d53SLisandro Dalcin /*@ 51695fce210SBarry Smith PetscSFGetLeafRange - Get the active leaf ranges 51795fce210SBarry Smith 51895fce210SBarry Smith Not Collective 51995fce210SBarry Smith 52095fce210SBarry Smith Input Arguments: 52195fce210SBarry Smith . sf - star forest 52295fce210SBarry Smith 52395fce210SBarry Smith Output Arguments: 52495fce210SBarry Smith + minleaf - minimum active leaf on this process 52595fce210SBarry Smith - maxleaf - maximum active leaf on this process 52695fce210SBarry Smith 52795fce210SBarry Smith Level: developer 52895fce210SBarry Smith 52995fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph(), PetscSFGetGraph() 53095fce210SBarry Smith @*/ 53195fce210SBarry Smith PetscErrorCode PetscSFGetLeafRange(PetscSF sf,PetscInt *minleaf,PetscInt *maxleaf) 53295fce210SBarry Smith { 53395fce210SBarry Smith 53495fce210SBarry Smith PetscFunctionBegin; 53595fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 536*29046d53SLisandro Dalcin PetscSFCheckGraphSet(sf,1); 53795fce210SBarry Smith if (minleaf) *minleaf = sf->minleaf; 53895fce210SBarry Smith if (maxleaf) *maxleaf = sf->maxleaf; 53995fce210SBarry Smith PetscFunctionReturn(0); 54095fce210SBarry Smith } 54195fce210SBarry Smith 54295fce210SBarry Smith /*@C 54395fce210SBarry Smith PetscSFView - view a star forest 54495fce210SBarry Smith 54595fce210SBarry Smith Collective 54695fce210SBarry Smith 54795fce210SBarry Smith Input Arguments: 54895fce210SBarry Smith + sf - star forest 54995fce210SBarry Smith - viewer - viewer to display graph, for example PETSC_VIEWER_STDOUT_WORLD 55095fce210SBarry Smith 55195fce210SBarry Smith Level: beginner 55295fce210SBarry Smith 55395fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph() 55495fce210SBarry Smith @*/ 55595fce210SBarry Smith PetscErrorCode PetscSFView(PetscSF sf,PetscViewer viewer) 55695fce210SBarry Smith { 55795fce210SBarry Smith PetscErrorCode ierr; 55895fce210SBarry Smith PetscBool iascii; 55995fce210SBarry Smith PetscViewerFormat format; 56095fce210SBarry Smith 56195fce210SBarry Smith PetscFunctionBegin; 56295fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 56395fce210SBarry Smith if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sf),&viewer);CHKERRQ(ierr);} 56495fce210SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 56595fce210SBarry Smith PetscCheckSameComm(sf,1,viewer,2); 56680153354SVaclav Hapla if (sf->graphset) {ierr = PetscSFSetUp(sf);CHKERRQ(ierr);} 56795fce210SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 56895fce210SBarry Smith if (iascii) { 56995fce210SBarry Smith PetscMPIInt rank; 57081bfa7aaSJed Brown PetscInt ii,i,j; 57195fce210SBarry Smith 572dae58748SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)sf,viewer);CHKERRQ(ierr); 57395fce210SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 57495fce210SBarry Smith if (sf->ops->View) {ierr = (*sf->ops->View)(sf,viewer);CHKERRQ(ierr);} 57580153354SVaclav Hapla if (!sf->graphset) { 57680153354SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"PetscSFSetGraph() has not been called yet\n");CHKERRQ(ierr); 57780153354SVaclav Hapla ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 57880153354SVaclav Hapla PetscFunctionReturn(0); 57980153354SVaclav Hapla } 58095fce210SBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); 5811575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 58295fce210SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of roots=%D, leaves=%D, remote ranks=%D\n",rank,sf->nroots,sf->nleaves,sf->nranks);CHKERRQ(ierr); 58395fce210SBarry Smith for (i=0; i<sf->nleaves; i++) { 58495fce210SBarry 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); 58595fce210SBarry Smith } 58695fce210SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 58795fce210SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 58895fce210SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 58981bfa7aaSJed Brown PetscMPIInt *tmpranks,*perm; 59081bfa7aaSJed Brown ierr = PetscMalloc2(sf->nranks,&tmpranks,sf->nranks,&perm);CHKERRQ(ierr); 59181bfa7aaSJed Brown ierr = PetscMemcpy(tmpranks,sf->ranks,sf->nranks*sizeof(tmpranks[0]));CHKERRQ(ierr); 59281bfa7aaSJed Brown for (i=0; i<sf->nranks; i++) perm[i] = i; 59381bfa7aaSJed Brown ierr = PetscSortMPIIntWithArray(sf->nranks,tmpranks,perm);CHKERRQ(ierr); 59495fce210SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Roots referenced by my leaves, by rank\n",rank);CHKERRQ(ierr); 59581bfa7aaSJed Brown for (ii=0; ii<sf->nranks; ii++) { 59681bfa7aaSJed Brown i = perm[ii]; 5977904a332SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %d: %D edges\n",rank,sf->ranks[i],sf->roffset[i+1]-sf->roffset[i]);CHKERRQ(ierr); 59895fce210SBarry Smith for (j=sf->roffset[i]; j<sf->roffset[i+1]; j++) { 59995fce210SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D <- %D\n",rank,sf->rmine[j],sf->rremote[j]);CHKERRQ(ierr); 60095fce210SBarry Smith } 60195fce210SBarry Smith } 60281bfa7aaSJed Brown ierr = PetscFree2(tmpranks,perm);CHKERRQ(ierr); 60395fce210SBarry Smith } 60495fce210SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 6051575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 60695fce210SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 60795fce210SBarry Smith } 60895fce210SBarry Smith PetscFunctionReturn(0); 60995fce210SBarry Smith } 61095fce210SBarry Smith 61195fce210SBarry Smith /*@C 61295fce210SBarry Smith PetscSFGetRanks - Get ranks and number of vertices referenced by leaves on this process 61395fce210SBarry Smith 61495fce210SBarry Smith Not Collective 61595fce210SBarry Smith 61695fce210SBarry Smith Input Arguments: 61795fce210SBarry Smith . sf - star forest 61895fce210SBarry Smith 61995fce210SBarry Smith Output Arguments: 62095fce210SBarry Smith + nranks - number of ranks referenced by local part 62195fce210SBarry Smith . ranks - array of ranks 62295fce210SBarry Smith . roffset - offset in rmine/rremote for each rank (length nranks+1) 62395fce210SBarry Smith . rmine - concatenated array holding local indices referencing each remote rank 62495fce210SBarry Smith - rremote - concatenated array holding remote indices referenced for each remote rank 62595fce210SBarry Smith 62695fce210SBarry Smith Level: developer 62795fce210SBarry Smith 62895fce210SBarry Smith .seealso: PetscSFSetGraph() 62995fce210SBarry Smith @*/ 63095fce210SBarry Smith PetscErrorCode PetscSFGetRanks(PetscSF sf,PetscInt *nranks,const PetscMPIInt **ranks,const PetscInt **roffset,const PetscInt **rmine,const PetscInt **rremote) 63195fce210SBarry Smith { 63295fce210SBarry Smith 63395fce210SBarry Smith PetscFunctionBegin; 63495fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 635*29046d53SLisandro Dalcin if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining ranks"); 63695fce210SBarry Smith if (nranks) *nranks = sf->nranks; 63795fce210SBarry Smith if (ranks) *ranks = sf->ranks; 63895fce210SBarry Smith if (roffset) *roffset = sf->roffset; 63995fce210SBarry Smith if (rmine) *rmine = sf->rmine; 64095fce210SBarry Smith if (rremote) *rremote = sf->rremote; 64195fce210SBarry Smith PetscFunctionReturn(0); 64295fce210SBarry Smith } 64395fce210SBarry Smith 644b5a8e515SJed Brown static PetscBool InList(PetscMPIInt needle,PetscMPIInt n,const PetscMPIInt *list) { 645b5a8e515SJed Brown PetscInt i; 646b5a8e515SJed Brown for (i=0; i<n; i++) { 647b5a8e515SJed Brown if (needle == list[i]) return PETSC_TRUE; 648b5a8e515SJed Brown } 649b5a8e515SJed Brown return PETSC_FALSE; 650b5a8e515SJed Brown } 651b5a8e515SJed Brown 65295fce210SBarry Smith /*@C 65321c688dcSJed Brown PetscSFSetUpRanks - Set up data structures associated with ranks; this is for internal use by PetscSF implementations. 65421c688dcSJed Brown 65521c688dcSJed Brown Collective 65621c688dcSJed Brown 65721c688dcSJed Brown Input Arguments: 658b5a8e515SJed Brown + sf - PetscSF to set up; PetscSFSetGraph() must have been called 659b5a8e515SJed Brown - dgroup - MPI_Group of ranks to be distinguished (e.g., for self or shared memory exchange) 66021c688dcSJed Brown 66121c688dcSJed Brown Level: developer 66221c688dcSJed Brown 66321c688dcSJed Brown .seealso: PetscSFGetRanks() 66421c688dcSJed Brown @*/ 665b5a8e515SJed Brown PetscErrorCode PetscSFSetUpRanks(PetscSF sf,MPI_Group dgroup) 66621c688dcSJed Brown { 66721c688dcSJed Brown PetscErrorCode ierr; 66821c688dcSJed Brown PetscTable table; 66921c688dcSJed Brown PetscTablePosition pos; 670b5a8e515SJed Brown PetscMPIInt size,groupsize,*groupranks; 67121c688dcSJed Brown PetscInt i,*rcount,*ranks; 67221c688dcSJed Brown 67321c688dcSJed Brown PetscFunctionBegin; 67421c688dcSJed Brown PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 675*29046d53SLisandro Dalcin PetscSFCheckGraphSet(sf,1); 67621c688dcSJed Brown ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr); 67721c688dcSJed Brown ierr = PetscTableCreate(10,size,&table);CHKERRQ(ierr); 67821c688dcSJed Brown for (i=0; i<sf->nleaves; i++) { 67921c688dcSJed Brown /* Log 1-based rank */ 68021c688dcSJed Brown ierr = PetscTableAdd(table,sf->remote[i].rank+1,1,ADD_VALUES);CHKERRQ(ierr); 68121c688dcSJed Brown } 68221c688dcSJed Brown ierr = PetscTableGetCount(table,&sf->nranks);CHKERRQ(ierr); 68321c688dcSJed Brown ierr = PetscMalloc4(sf->nranks,&sf->ranks,sf->nranks+1,&sf->roffset,sf->nleaves,&sf->rmine,sf->nleaves,&sf->rremote);CHKERRQ(ierr); 68421c688dcSJed Brown ierr = PetscMalloc2(sf->nranks,&rcount,sf->nranks,&ranks);CHKERRQ(ierr); 68521c688dcSJed Brown ierr = PetscTableGetHeadPosition(table,&pos);CHKERRQ(ierr); 68621c688dcSJed Brown for (i=0; i<sf->nranks; i++) { 68721c688dcSJed Brown ierr = PetscTableGetNext(table,&pos,&ranks[i],&rcount[i]);CHKERRQ(ierr); 68821c688dcSJed Brown ranks[i]--; /* Convert back to 0-based */ 68921c688dcSJed Brown } 69021c688dcSJed Brown ierr = PetscTableDestroy(&table);CHKERRQ(ierr); 691b5a8e515SJed Brown 692b5a8e515SJed Brown /* We expect that dgroup is reliably "small" while nranks could be large */ 693b5a8e515SJed Brown { 6947fb8a5e4SKarl Rupp MPI_Group group = MPI_GROUP_NULL; 695b5a8e515SJed Brown PetscMPIInt *dgroupranks; 696b5a8e515SJed Brown ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr); 697b5a8e515SJed Brown ierr = MPI_Group_size(dgroup,&groupsize);CHKERRQ(ierr); 698b5a8e515SJed Brown ierr = PetscMalloc1(groupsize,&dgroupranks);CHKERRQ(ierr); 699b5a8e515SJed Brown ierr = PetscMalloc1(groupsize,&groupranks);CHKERRQ(ierr); 700b5a8e515SJed Brown for (i=0; i<groupsize; i++) dgroupranks[i] = i; 701f3fc9a17SJed Brown if (groupsize) {ierr = MPI_Group_translate_ranks(dgroup,groupsize,dgroupranks,group,groupranks);CHKERRQ(ierr);} 702b5a8e515SJed Brown ierr = MPI_Group_free(&group);CHKERRQ(ierr); 703b5a8e515SJed Brown ierr = PetscFree(dgroupranks);CHKERRQ(ierr); 704b5a8e515SJed Brown } 705b5a8e515SJed Brown 706b5a8e515SJed Brown /* Partition ranks[] into distinguished (first sf->ndranks) followed by non-distinguished */ 707b5a8e515SJed Brown for (sf->ndranks=0,i=sf->nranks; sf->ndranks<i; ) { 708b5a8e515SJed Brown for (i--; sf->ndranks<i; i--) { /* Scan i backward looking for distinguished rank */ 709b5a8e515SJed Brown if (InList(ranks[i],groupsize,groupranks)) break; 710b5a8e515SJed Brown } 711b5a8e515SJed Brown for ( ; sf->ndranks<=i; sf->ndranks++) { /* Scan sf->ndranks forward looking for non-distinguished rank */ 712b5a8e515SJed Brown if (!InList(ranks[sf->ndranks],groupsize,groupranks)) break; 713b5a8e515SJed Brown } 714b5a8e515SJed Brown if (sf->ndranks < i) { /* Swap ranks[sf->ndranks] with ranks[i] */ 715b5a8e515SJed Brown PetscInt tmprank,tmpcount; 716b5a8e515SJed Brown tmprank = ranks[i]; 717b5a8e515SJed Brown tmpcount = rcount[i]; 718b5a8e515SJed Brown ranks[i] = ranks[sf->ndranks]; 719b5a8e515SJed Brown rcount[i] = rcount[sf->ndranks]; 720b5a8e515SJed Brown ranks[sf->ndranks] = tmprank; 721b5a8e515SJed Brown rcount[sf->ndranks] = tmpcount; 722b5a8e515SJed Brown sf->ndranks++; 723b5a8e515SJed Brown } 724b5a8e515SJed Brown } 725b5a8e515SJed Brown ierr = PetscFree(groupranks);CHKERRQ(ierr); 726b5a8e515SJed Brown ierr = PetscSortIntWithArray(sf->ndranks,ranks,rcount);CHKERRQ(ierr); 727b5a8e515SJed Brown ierr = PetscSortIntWithArray(sf->nranks-sf->ndranks,ranks+sf->ndranks,rcount+sf->ndranks);CHKERRQ(ierr); 72821c688dcSJed Brown sf->roffset[0] = 0; 72921c688dcSJed Brown for (i=0; i<sf->nranks; i++) { 73021c688dcSJed Brown ierr = PetscMPIIntCast(ranks[i],sf->ranks+i);CHKERRQ(ierr); 73121c688dcSJed Brown sf->roffset[i+1] = sf->roffset[i] + rcount[i]; 73221c688dcSJed Brown rcount[i] = 0; 73321c688dcSJed Brown } 73421c688dcSJed Brown for (i=0; i<sf->nleaves; i++) { 735b5a8e515SJed Brown PetscInt irank; 73621c688dcSJed Brown /* Search for index of iremote[i].rank in sf->ranks */ 737b5a8e515SJed Brown ierr = PetscFindMPIInt(sf->remote[i].rank,sf->ndranks,sf->ranks,&irank);CHKERRQ(ierr); 738b5a8e515SJed Brown if (irank < 0) { 739b5a8e515SJed Brown ierr = PetscFindMPIInt(sf->remote[i].rank,sf->nranks-sf->ndranks,sf->ranks+sf->ndranks,&irank);CHKERRQ(ierr); 740b5a8e515SJed Brown if (irank >= 0) irank += sf->ndranks; 74121c688dcSJed Brown } 742b5a8e515SJed Brown if (irank < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Could not find rank %D in array",sf->remote[i].rank); 74321c688dcSJed Brown sf->rmine[sf->roffset[irank] + rcount[irank]] = sf->mine ? sf->mine[i] : i; 74421c688dcSJed Brown sf->rremote[sf->roffset[irank] + rcount[irank]] = sf->remote[i].index; 74521c688dcSJed Brown rcount[irank]++; 74621c688dcSJed Brown } 74721c688dcSJed Brown ierr = PetscFree2(rcount,ranks);CHKERRQ(ierr); 74821c688dcSJed Brown PetscFunctionReturn(0); 74921c688dcSJed Brown } 75021c688dcSJed Brown 75121c688dcSJed Brown /*@C 75295fce210SBarry Smith PetscSFGetGroups - gets incoming and outgoing process groups 75395fce210SBarry Smith 75495fce210SBarry Smith Collective 75595fce210SBarry Smith 75695fce210SBarry Smith Input Argument: 75795fce210SBarry Smith . sf - star forest 75895fce210SBarry Smith 75995fce210SBarry Smith Output Arguments: 76095fce210SBarry Smith + incoming - group of origin processes for incoming edges (leaves that reference my roots) 76195fce210SBarry Smith - outgoing - group of destination processes for outgoing edges (roots that I reference) 76295fce210SBarry Smith 76395fce210SBarry Smith Level: developer 76495fce210SBarry Smith 76595fce210SBarry Smith .seealso: PetscSFGetWindow(), PetscSFRestoreWindow() 76695fce210SBarry Smith @*/ 76795fce210SBarry Smith PetscErrorCode PetscSFGetGroups(PetscSF sf,MPI_Group *incoming,MPI_Group *outgoing) 76895fce210SBarry Smith { 76995fce210SBarry Smith PetscErrorCode ierr; 7707fb8a5e4SKarl Rupp MPI_Group group = MPI_GROUP_NULL; 77195fce210SBarry Smith 77295fce210SBarry Smith PetscFunctionBegin; 773*29046d53SLisandro Dalcin if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining groups"); 77495fce210SBarry Smith if (sf->ingroup == MPI_GROUP_NULL) { 77595fce210SBarry Smith PetscInt i; 77695fce210SBarry Smith const PetscInt *indegree; 77795fce210SBarry Smith PetscMPIInt rank,*outranks,*inranks; 77895fce210SBarry Smith PetscSFNode *remote; 77995fce210SBarry Smith PetscSF bgcount; 78095fce210SBarry Smith 78195fce210SBarry Smith /* Compute the number of incoming ranks */ 782785e854fSJed Brown ierr = PetscMalloc1(sf->nranks,&remote);CHKERRQ(ierr); 78395fce210SBarry Smith for (i=0; i<sf->nranks; i++) { 78495fce210SBarry Smith remote[i].rank = sf->ranks[i]; 78595fce210SBarry Smith remote[i].index = 0; 78695fce210SBarry Smith } 78795fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,&bgcount);CHKERRQ(ierr); 78895fce210SBarry Smith ierr = PetscSFSetGraph(bgcount,1,sf->nranks,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr); 78995fce210SBarry Smith ierr = PetscSFComputeDegreeBegin(bgcount,&indegree);CHKERRQ(ierr); 79095fce210SBarry Smith ierr = PetscSFComputeDegreeEnd(bgcount,&indegree);CHKERRQ(ierr); 79195fce210SBarry Smith 79295fce210SBarry Smith /* Enumerate the incoming ranks */ 793dcca6d9dSJed Brown ierr = PetscMalloc2(indegree[0],&inranks,sf->nranks,&outranks);CHKERRQ(ierr); 79495fce210SBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); 79595fce210SBarry Smith for (i=0; i<sf->nranks; i++) outranks[i] = rank; 79695fce210SBarry Smith ierr = PetscSFGatherBegin(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr); 79795fce210SBarry Smith ierr = PetscSFGatherEnd(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr); 79895fce210SBarry Smith ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr); 79995fce210SBarry Smith ierr = MPI_Group_incl(group,indegree[0],inranks,&sf->ingroup);CHKERRQ(ierr); 80095fce210SBarry Smith ierr = MPI_Group_free(&group);CHKERRQ(ierr); 80195fce210SBarry Smith ierr = PetscFree2(inranks,outranks);CHKERRQ(ierr); 80295fce210SBarry Smith ierr = PetscSFDestroy(&bgcount);CHKERRQ(ierr); 80395fce210SBarry Smith } 80495fce210SBarry Smith *incoming = sf->ingroup; 80595fce210SBarry Smith 80695fce210SBarry Smith if (sf->outgroup == MPI_GROUP_NULL) { 80795fce210SBarry Smith ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr); 80895fce210SBarry Smith ierr = MPI_Group_incl(group,sf->nranks,sf->ranks,&sf->outgroup);CHKERRQ(ierr); 80995fce210SBarry Smith ierr = MPI_Group_free(&group);CHKERRQ(ierr); 81095fce210SBarry Smith } 81195fce210SBarry Smith *outgoing = sf->outgroup; 81295fce210SBarry Smith PetscFunctionReturn(0); 81395fce210SBarry Smith } 81495fce210SBarry Smith 815*29046d53SLisandro Dalcin /*@ 81695fce210SBarry Smith PetscSFGetMultiSF - gets the inner SF implemeting gathers and scatters 81795fce210SBarry Smith 81895fce210SBarry Smith Collective 81995fce210SBarry Smith 82095fce210SBarry Smith Input Argument: 82195fce210SBarry Smith . sf - star forest that may contain roots with 0 or with more than 1 vertex 82295fce210SBarry Smith 82395fce210SBarry Smith Output Arguments: 82495fce210SBarry Smith . multi - star forest with split roots, such that each root has degree exactly 1 82595fce210SBarry Smith 82695fce210SBarry Smith Level: developer 82795fce210SBarry Smith 82895fce210SBarry Smith Notes: 82995fce210SBarry Smith 83095fce210SBarry Smith In most cases, users should use PetscSFGatherBegin() and PetscSFScatterBegin() instead of manipulating multi 83195fce210SBarry Smith directly. Since multi satisfies the stronger condition that each entry in the global space has exactly one incoming 83295fce210SBarry Smith edge, it is a candidate for future optimization that might involve its removal. 83395fce210SBarry Smith 83495fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFGatherBegin(), PetscSFScatterBegin() 83595fce210SBarry Smith @*/ 83695fce210SBarry Smith PetscErrorCode PetscSFGetMultiSF(PetscSF sf,PetscSF *multi) 83795fce210SBarry Smith { 83895fce210SBarry Smith PetscErrorCode ierr; 83995fce210SBarry Smith 84095fce210SBarry Smith PetscFunctionBegin; 84195fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 84295fce210SBarry Smith PetscValidPointer(multi,2); 84395fce210SBarry Smith if (sf->nroots < 0) { /* Graph has not been set yet; why do we need this? */ 84495fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr); 84595fce210SBarry Smith *multi = sf->multi; 84695fce210SBarry Smith PetscFunctionReturn(0); 84795fce210SBarry Smith } 84895fce210SBarry Smith if (!sf->multi) { 84995fce210SBarry Smith const PetscInt *indegree; 8509837ea96SMatthew G. Knepley PetscInt i,*inoffset,*outones,*outoffset,maxlocal; 85195fce210SBarry Smith PetscSFNode *remote; 852*29046d53SLisandro Dalcin maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */ 85395fce210SBarry Smith ierr = PetscSFComputeDegreeBegin(sf,&indegree);CHKERRQ(ierr); 85495fce210SBarry Smith ierr = PetscSFComputeDegreeEnd(sf,&indegree);CHKERRQ(ierr); 8559837ea96SMatthew G. Knepley ierr = PetscMalloc3(sf->nroots+1,&inoffset,maxlocal,&outones,maxlocal,&outoffset);CHKERRQ(ierr); 85695fce210SBarry Smith inoffset[0] = 0; 85795fce210SBarry Smith for (i=0; i<sf->nroots; i++) inoffset[i+1] = inoffset[i] + indegree[i]; 8589837ea96SMatthew G. Knepley for (i=0; i<maxlocal; i++) outones[i] = 1; 859dbd2ff41SBarry Smith ierr = PetscSFFetchAndOpBegin(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr); 860dbd2ff41SBarry Smith ierr = PetscSFFetchAndOpEnd(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr); 86195fce210SBarry Smith for (i=0; i<sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */ 86295fce210SBarry Smith #if 0 86395fce210SBarry Smith #if defined(PETSC_USE_DEBUG) /* Check that the expected number of increments occurred */ 86495fce210SBarry Smith for (i=0; i<sf->nroots; i++) { 86595fce210SBarry Smith if (inoffset[i] + indegree[i] != inoffset[i+1]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect result after PetscSFFetchAndOp"); 86695fce210SBarry Smith } 86795fce210SBarry Smith #endif 86895fce210SBarry Smith #endif 869785e854fSJed Brown ierr = PetscMalloc1(sf->nleaves,&remote);CHKERRQ(ierr); 87095fce210SBarry Smith for (i=0; i<sf->nleaves; i++) { 87195fce210SBarry Smith remote[i].rank = sf->remote[i].rank; 87238e7336fSToby Isaac remote[i].index = outoffset[sf->mine ? sf->mine[i] : i]; 87395fce210SBarry Smith } 87495fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr); 87501365b40SToby Isaac ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr); 87695fce210SBarry Smith if (sf->rankorder) { /* Sort the ranks */ 87795fce210SBarry Smith PetscMPIInt rank; 87895fce210SBarry Smith PetscInt *inranks,*newoffset,*outranks,*newoutoffset,*tmpoffset,maxdegree; 87995fce210SBarry Smith PetscSFNode *newremote; 88095fce210SBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); 88195fce210SBarry Smith for (i=0,maxdegree=0; i<sf->nroots; i++) maxdegree = PetscMax(maxdegree,indegree[i]); 8829837ea96SMatthew G. Knepley ierr = PetscMalloc5(sf->multi->nroots,&inranks,sf->multi->nroots,&newoffset,maxlocal,&outranks,maxlocal,&newoutoffset,maxdegree,&tmpoffset);CHKERRQ(ierr); 8839837ea96SMatthew G. Knepley for (i=0; i<maxlocal; i++) outranks[i] = rank; 8848bfbc91cSJed Brown ierr = PetscSFReduceBegin(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr); 8858bfbc91cSJed Brown ierr = PetscSFReduceEnd(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr); 88695fce210SBarry Smith /* Sort the incoming ranks at each vertex, build the inverse map */ 88795fce210SBarry Smith for (i=0; i<sf->nroots; i++) { 88895fce210SBarry Smith PetscInt j; 88995fce210SBarry Smith for (j=0; j<indegree[i]; j++) tmpoffset[j] = j; 89095fce210SBarry Smith ierr = PetscSortIntWithArray(indegree[i],inranks+inoffset[i],tmpoffset);CHKERRQ(ierr); 89195fce210SBarry Smith for (j=0; j<indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j; 89295fce210SBarry Smith } 89395fce210SBarry Smith ierr = PetscSFBcastBegin(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr); 89495fce210SBarry Smith ierr = PetscSFBcastEnd(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr); 895785e854fSJed Brown ierr = PetscMalloc1(sf->nleaves,&newremote);CHKERRQ(ierr); 89695fce210SBarry Smith for (i=0; i<sf->nleaves; i++) { 89795fce210SBarry Smith newremote[i].rank = sf->remote[i].rank; 89801365b40SToby Isaac newremote[i].index = newoutoffset[sf->mine ? sf->mine[i] : i]; 89995fce210SBarry Smith } 90001365b40SToby Isaac ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,newremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 90195fce210SBarry Smith ierr = PetscFree5(inranks,newoffset,outranks,newoutoffset,tmpoffset);CHKERRQ(ierr); 90295fce210SBarry Smith } 90395fce210SBarry Smith ierr = PetscFree3(inoffset,outones,outoffset);CHKERRQ(ierr); 90495fce210SBarry Smith } 90595fce210SBarry Smith *multi = sf->multi; 90695fce210SBarry Smith PetscFunctionReturn(0); 90795fce210SBarry Smith } 90895fce210SBarry Smith 90995fce210SBarry Smith /*@C 91095fce210SBarry Smith PetscSFCreateEmbeddedSF - removes edges from all but the selected roots, does not remap indices 91195fce210SBarry Smith 91295fce210SBarry Smith Collective 91395fce210SBarry Smith 91495fce210SBarry Smith Input Arguments: 91595fce210SBarry Smith + sf - original star forest 91695fce210SBarry Smith . nroots - number of roots to select on this process 91795fce210SBarry Smith - selected - selected roots on this process 91895fce210SBarry Smith 91995fce210SBarry Smith Output Arguments: 92095fce210SBarry Smith . newsf - new star forest 92195fce210SBarry Smith 92295fce210SBarry Smith Level: advanced 92395fce210SBarry Smith 92495fce210SBarry Smith Note: 92595fce210SBarry Smith To use the new PetscSF, it may be necessary to know the indices of the leaves that are still participating. This can 92695fce210SBarry Smith be done by calling PetscSFGetGraph(). 92795fce210SBarry Smith 92895fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFGetGraph() 92995fce210SBarry Smith @*/ 93095fce210SBarry Smith PetscErrorCode PetscSFCreateEmbeddedSF(PetscSF sf,PetscInt nroots,const PetscInt *selected,PetscSF *newsf) 93195fce210SBarry Smith { 9320511a646SMatthew G. Knepley PetscInt *rootdata, *leafdata, *ilocal; 93395fce210SBarry Smith PetscSFNode *iremote; 934fc1ede2bSMatthew G. Knepley PetscInt leafsize = 0, nleaves = 0, n, i; 9350511a646SMatthew G. Knepley PetscErrorCode ierr; 93695fce210SBarry Smith 93795fce210SBarry Smith PetscFunctionBegin; 93895fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 939*29046d53SLisandro Dalcin PetscSFCheckGraphSet(sf,1); 94095fce210SBarry Smith if (nroots) PetscValidPointer(selected,3); 94195fce210SBarry Smith PetscValidPointer(newsf,4); 9420511a646SMatthew G. Knepley if (sf->mine) for (i = 0; i < sf->nleaves; ++i) {leafsize = PetscMax(leafsize, sf->mine[i]+1);} 9430511a646SMatthew G. Knepley else leafsize = sf->nleaves; 9441795a4d1SJed Brown ierr = PetscCalloc2(sf->nroots,&rootdata,leafsize,&leafdata);CHKERRQ(ierr); 9450511a646SMatthew G. Knepley for (i=0; i<nroots; ++i) rootdata[selected[i]] = 1; 94695fce210SBarry Smith ierr = PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr); 94795fce210SBarry Smith ierr = PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr); 9480511a646SMatthew G. Knepley for (i = 0; i < leafsize; ++i) nleaves += leafdata[i]; 949785e854fSJed Brown ierr = PetscMalloc1(nleaves,&ilocal);CHKERRQ(ierr); 950785e854fSJed Brown ierr = PetscMalloc1(nleaves,&iremote);CHKERRQ(ierr); 951fc1ede2bSMatthew G. Knepley for (i = 0, n = 0; i < sf->nleaves; ++i) { 9520511a646SMatthew G. Knepley const PetscInt lidx = sf->mine ? sf->mine[i] : i; 9530511a646SMatthew G. Knepley 9540511a646SMatthew G. Knepley if (leafdata[lidx]) { 955fc1ede2bSMatthew G. Knepley ilocal[n] = lidx; 956fc1ede2bSMatthew G. Knepley iremote[n].rank = sf->remote[i].rank; 957fc1ede2bSMatthew G. Knepley iremote[n].index = sf->remote[i].index; 958fc1ede2bSMatthew G. Knepley ++n; 95995fce210SBarry Smith } 96095fce210SBarry Smith } 961fc1ede2bSMatthew G. Knepley if (n != nleaves) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "There is a size mismatch in the SF embedding, %d != %d", n, nleaves); 96295fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,newsf);CHKERRQ(ierr); 96395fce210SBarry Smith ierr = PetscSFSetGraph(*newsf,sf->nroots,nleaves,ilocal,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 96495fce210SBarry Smith ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr); 96595fce210SBarry Smith PetscFunctionReturn(0); 96695fce210SBarry Smith } 96795fce210SBarry Smith 9682f5fb4c2SMatthew G. Knepley /*@C 9692f5fb4c2SMatthew G. Knepley PetscSFCreateEmbeddedLeafSF - removes edges from all but the selected leaves, does not remap indices 9702f5fb4c2SMatthew G. Knepley 9712f5fb4c2SMatthew G. Knepley Collective 9722f5fb4c2SMatthew G. Knepley 9732f5fb4c2SMatthew G. Knepley Input Arguments: 9742f5fb4c2SMatthew G. Knepley + sf - original star forest 9752f5fb4c2SMatthew G. Knepley . nleaves - number of leaves to select on this process 9762f5fb4c2SMatthew G. Knepley - selected - selected leaves on this process 9772f5fb4c2SMatthew G. Knepley 9782f5fb4c2SMatthew G. Knepley Output Arguments: 9792f5fb4c2SMatthew G. Knepley . newsf - new star forest 9802f5fb4c2SMatthew G. Knepley 9812f5fb4c2SMatthew G. Knepley Level: advanced 9822f5fb4c2SMatthew G. Knepley 9832f5fb4c2SMatthew G. Knepley .seealso: PetscSFCreateEmbeddedSF(), PetscSFSetGraph(), PetscSFGetGraph() 9842f5fb4c2SMatthew G. Knepley @*/ 9852f5fb4c2SMatthew G. Knepley PetscErrorCode PetscSFCreateEmbeddedLeafSF(PetscSF sf, PetscInt nleaves, const PetscInt *selected, PetscSF *newsf) 9862f5fb4c2SMatthew G. Knepley { 9872f5fb4c2SMatthew G. Knepley PetscSFNode *iremote; 9882f5fb4c2SMatthew G. Knepley PetscInt *ilocal; 9892f5fb4c2SMatthew G. Knepley PetscInt i; 9902f5fb4c2SMatthew G. Knepley PetscErrorCode ierr; 9912f5fb4c2SMatthew G. Knepley 9922f5fb4c2SMatthew G. Knepley PetscFunctionBegin; 9932f5fb4c2SMatthew G. Knepley PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 994*29046d53SLisandro Dalcin PetscSFCheckGraphSet(sf, 1); 9952f5fb4c2SMatthew G. Knepley if (nleaves) PetscValidPointer(selected, 3); 9962f5fb4c2SMatthew G. Knepley PetscValidPointer(newsf, 4); 9972f5fb4c2SMatthew G. Knepley ierr = PetscMalloc1(nleaves, &ilocal);CHKERRQ(ierr); 9982f5fb4c2SMatthew G. Knepley ierr = PetscMalloc1(nleaves, &iremote);CHKERRQ(ierr); 9992f5fb4c2SMatthew G. Knepley for (i = 0; i < nleaves; ++i) { 10002f5fb4c2SMatthew G. Knepley const PetscInt l = selected[i]; 10012f5fb4c2SMatthew G. Knepley 10022f5fb4c2SMatthew G. Knepley ilocal[i] = sf->mine ? sf->mine[l] : l; 10032f5fb4c2SMatthew G. Knepley iremote[i].rank = sf->remote[l].rank; 10042f5fb4c2SMatthew G. Knepley iremote[i].index = sf->remote[l].index; 10052f5fb4c2SMatthew G. Knepley } 10062f5fb4c2SMatthew G. Knepley ierr = PetscSFDuplicate(sf, PETSCSF_DUPLICATE_RANKS, newsf);CHKERRQ(ierr); 10072f5fb4c2SMatthew G. Knepley ierr = PetscSFSetGraph(*newsf, sf->nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr); 10082f5fb4c2SMatthew G. Knepley PetscFunctionReturn(0); 10092f5fb4c2SMatthew G. Knepley } 10102f5fb4c2SMatthew G. Knepley 101195fce210SBarry Smith /*@C 101295fce210SBarry Smith PetscSFBcastBegin - begin pointwise broadcast to be concluded with call to PetscSFBcastEnd() 101395fce210SBarry Smith 101495fce210SBarry Smith Collective on PetscSF 101595fce210SBarry Smith 101695fce210SBarry Smith Input Arguments: 101795fce210SBarry Smith + sf - star forest on which to communicate 101895fce210SBarry Smith . unit - data type associated with each node 101995fce210SBarry Smith - rootdata - buffer to broadcast 102095fce210SBarry Smith 102195fce210SBarry Smith Output Arguments: 102295fce210SBarry Smith . leafdata - buffer to update with values from each leaf's respective root 102395fce210SBarry Smith 102495fce210SBarry Smith Level: intermediate 102595fce210SBarry Smith 102695fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFView(), PetscSFBcastEnd(), PetscSFReduceBegin() 102795fce210SBarry Smith @*/ 102895fce210SBarry Smith PetscErrorCode PetscSFBcastBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata) 102995fce210SBarry Smith { 103095fce210SBarry Smith PetscErrorCode ierr; 103195fce210SBarry Smith 103295fce210SBarry Smith PetscFunctionBegin; 103395fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 103495fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 1035*29046d53SLisandro Dalcin ierr = PetscLogEventBegin(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr); 103695fce210SBarry Smith ierr = (*sf->ops->BcastBegin)(sf,unit,rootdata,leafdata);CHKERRQ(ierr); 1037563ffbabSMatthew G. Knepley ierr = PetscLogEventEnd(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr); 103895fce210SBarry Smith PetscFunctionReturn(0); 103995fce210SBarry Smith } 104095fce210SBarry Smith 104195fce210SBarry Smith /*@C 104295fce210SBarry Smith PetscSFBcastEnd - end a broadcast operation started with PetscSFBcastBegin() 104395fce210SBarry Smith 104495fce210SBarry Smith Collective 104595fce210SBarry Smith 104695fce210SBarry Smith Input Arguments: 104795fce210SBarry Smith + sf - star forest 104895fce210SBarry Smith . unit - data type 104995fce210SBarry Smith - rootdata - buffer to broadcast 105095fce210SBarry Smith 105195fce210SBarry Smith Output Arguments: 105295fce210SBarry Smith . leafdata - buffer to update with values from each leaf's respective root 105395fce210SBarry Smith 105495fce210SBarry Smith Level: intermediate 105595fce210SBarry Smith 105695fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFReduceEnd() 105795fce210SBarry Smith @*/ 105895fce210SBarry Smith PetscErrorCode PetscSFBcastEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata) 105995fce210SBarry Smith { 106095fce210SBarry Smith PetscErrorCode ierr; 106195fce210SBarry Smith 106295fce210SBarry Smith PetscFunctionBegin; 106395fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 106495fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 1065*29046d53SLisandro Dalcin ierr = PetscLogEventBegin(PETSCSF_BcastEnd,sf,0,0,0);CHKERRQ(ierr); 106695fce210SBarry Smith ierr = (*sf->ops->BcastEnd)(sf,unit,rootdata,leafdata);CHKERRQ(ierr); 1067563ffbabSMatthew G. Knepley ierr = PetscLogEventEnd(PETSCSF_BcastEnd,sf,0,0,0);CHKERRQ(ierr); 106895fce210SBarry Smith PetscFunctionReturn(0); 106995fce210SBarry Smith } 107095fce210SBarry Smith 107195fce210SBarry Smith /*@C 107295fce210SBarry Smith PetscSFReduceBegin - begin reduction of leafdata into rootdata, to be completed with call to PetscSFReduceEnd() 107395fce210SBarry Smith 107495fce210SBarry Smith Collective 107595fce210SBarry Smith 107695fce210SBarry Smith Input Arguments: 107795fce210SBarry Smith + sf - star forest 107895fce210SBarry Smith . unit - data type 107995fce210SBarry Smith . leafdata - values to reduce 108095fce210SBarry Smith - op - reduction operation 108195fce210SBarry Smith 108295fce210SBarry Smith Output Arguments: 108395fce210SBarry Smith . rootdata - result of reduction of values from all leaves of each root 108495fce210SBarry Smith 108595fce210SBarry Smith Level: intermediate 108695fce210SBarry Smith 108795fce210SBarry Smith .seealso: PetscSFBcastBegin() 108895fce210SBarry Smith @*/ 108995fce210SBarry Smith PetscErrorCode PetscSFReduceBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op) 109095fce210SBarry Smith { 109195fce210SBarry Smith PetscErrorCode ierr; 109295fce210SBarry Smith 109395fce210SBarry Smith PetscFunctionBegin; 109495fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 109595fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 1096*29046d53SLisandro Dalcin ierr = PetscLogEventBegin(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr); 109795fce210SBarry Smith ierr = (sf->ops->ReduceBegin)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr); 1098563ffbabSMatthew G. Knepley ierr = PetscLogEventEnd(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr); 109995fce210SBarry Smith PetscFunctionReturn(0); 110095fce210SBarry Smith } 110195fce210SBarry Smith 110295fce210SBarry Smith /*@C 110395fce210SBarry Smith PetscSFReduceEnd - end a reduction operation started with PetscSFReduceBegin() 110495fce210SBarry Smith 110595fce210SBarry Smith Collective 110695fce210SBarry Smith 110795fce210SBarry Smith Input Arguments: 110895fce210SBarry Smith + sf - star forest 110995fce210SBarry Smith . unit - data type 111095fce210SBarry Smith . leafdata - values to reduce 111195fce210SBarry Smith - op - reduction operation 111295fce210SBarry Smith 111395fce210SBarry Smith Output Arguments: 111495fce210SBarry Smith . rootdata - result of reduction of values from all leaves of each root 111595fce210SBarry Smith 111695fce210SBarry Smith Level: intermediate 111795fce210SBarry Smith 111895fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFBcastEnd() 111995fce210SBarry Smith @*/ 112095fce210SBarry Smith PetscErrorCode PetscSFReduceEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op) 112195fce210SBarry Smith { 112295fce210SBarry Smith PetscErrorCode ierr; 112395fce210SBarry Smith 112495fce210SBarry Smith PetscFunctionBegin; 112595fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 112695fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 1127*29046d53SLisandro Dalcin ierr = PetscLogEventBegin(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr); 112895fce210SBarry Smith ierr = (*sf->ops->ReduceEnd)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr); 1129563ffbabSMatthew G. Knepley ierr = PetscLogEventEnd(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr); 113095fce210SBarry Smith PetscFunctionReturn(0); 113195fce210SBarry Smith } 113295fce210SBarry Smith 113395fce210SBarry Smith /*@C 113495fce210SBarry Smith PetscSFComputeDegreeBegin - begin computation of degree for each root vertex, to be completed with PetscSFComputeDegreeEnd() 113595fce210SBarry Smith 113695fce210SBarry Smith Collective 113795fce210SBarry Smith 113895fce210SBarry Smith Input Arguments: 113995fce210SBarry Smith . sf - star forest 114095fce210SBarry Smith 114195fce210SBarry Smith Output Arguments: 114295fce210SBarry Smith . degree - degree of each root vertex 114395fce210SBarry Smith 114495fce210SBarry Smith Level: advanced 114595fce210SBarry Smith 114695fce210SBarry Smith .seealso: PetscSFGatherBegin() 114795fce210SBarry Smith @*/ 114895fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf,const PetscInt **degree) 114995fce210SBarry Smith { 115095fce210SBarry Smith PetscErrorCode ierr; 115195fce210SBarry Smith 115295fce210SBarry Smith PetscFunctionBegin; 115395fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 115495fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 115595fce210SBarry Smith PetscValidPointer(degree,2); 1156803bd9e8SMatthew G. Knepley if (!sf->degreeknown) { 1157*29046d53SLisandro Dalcin PetscInt i, nroots = sf->nroots, maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */ 1158803bd9e8SMatthew G. Knepley if (sf->degree) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Calls to PetscSFComputeDegreeBegin() cannot be nested."); 1159*29046d53SLisandro Dalcin ierr = PetscMalloc1(nroots,&sf->degree);CHKERRQ(ierr); 1160*29046d53SLisandro Dalcin ierr = PetscMalloc1(PetscMax(maxlocal,1),&sf->degreetmp);CHKERRQ(ierr); /* allocate at least one entry, see check in PetscSFComputeDegreeEnd() */ 1161*29046d53SLisandro Dalcin for (i=0; i<nroots; i++) sf->degree[i] = 0; 11629837ea96SMatthew G. Knepley for (i=0; i<maxlocal; i++) sf->degreetmp[i] = 1; 1163dbd2ff41SBarry Smith ierr = PetscSFReduceBegin(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);CHKERRQ(ierr); 116495fce210SBarry Smith } 116595fce210SBarry Smith *degree = NULL; 116695fce210SBarry Smith PetscFunctionReturn(0); 116795fce210SBarry Smith } 116895fce210SBarry Smith 116995fce210SBarry Smith /*@C 117095fce210SBarry Smith PetscSFComputeDegreeEnd - complete computation of degree for each root vertex, started with PetscSFComputeDegreeBegin() 117195fce210SBarry Smith 117295fce210SBarry Smith Collective 117395fce210SBarry Smith 117495fce210SBarry Smith Input Arguments: 117595fce210SBarry Smith . sf - star forest 117695fce210SBarry Smith 117795fce210SBarry Smith Output Arguments: 117895fce210SBarry Smith . degree - degree of each root vertex 117995fce210SBarry Smith 118095fce210SBarry Smith Level: developer 118195fce210SBarry Smith 118295fce210SBarry Smith .seealso: 118395fce210SBarry Smith @*/ 118495fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf,const PetscInt **degree) 118595fce210SBarry Smith { 118695fce210SBarry Smith PetscErrorCode ierr; 118795fce210SBarry Smith 118895fce210SBarry Smith PetscFunctionBegin; 118995fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 119095fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 1191*29046d53SLisandro Dalcin PetscValidPointer(degree,2); 119295fce210SBarry Smith if (!sf->degreeknown) { 1193*29046d53SLisandro Dalcin if (!sf->degreetmp) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFComputeDegreeBegin() before PetscSFComputeDegreeEnd()"); 1194dbd2ff41SBarry Smith ierr = PetscSFReduceEnd(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);CHKERRQ(ierr); 119595fce210SBarry Smith ierr = PetscFree(sf->degreetmp);CHKERRQ(ierr); 119695fce210SBarry Smith sf->degreeknown = PETSC_TRUE; 119795fce210SBarry Smith } 119895fce210SBarry Smith *degree = sf->degree; 119995fce210SBarry Smith PetscFunctionReturn(0); 120095fce210SBarry Smith } 120195fce210SBarry Smith 120295fce210SBarry Smith /*@C 120395fce210SBarry Smith PetscSFFetchAndOpBegin - begin operation that fetches values from root and updates atomically by applying operation using my leaf value, to be completed with PetscSFFetchAndOpEnd() 120495fce210SBarry Smith 120595fce210SBarry Smith Collective 120695fce210SBarry Smith 120795fce210SBarry Smith Input Arguments: 120895fce210SBarry Smith + sf - star forest 120995fce210SBarry Smith . unit - data type 121095fce210SBarry Smith . leafdata - leaf values to use in reduction 121195fce210SBarry Smith - op - operation to use for reduction 121295fce210SBarry Smith 121395fce210SBarry Smith Output Arguments: 121495fce210SBarry Smith + rootdata - root values to be updated, input state is seen by first process to perform an update 121595fce210SBarry Smith - leafupdate - state at each leaf's respective root immediately prior to my atomic update 121695fce210SBarry Smith 121795fce210SBarry Smith Level: advanced 121895fce210SBarry Smith 121995fce210SBarry Smith Note: 122095fce210SBarry Smith The update is only atomic at the granularity provided by the hardware. Different roots referenced by the same process 122195fce210SBarry Smith might be updated in a different order. Furthermore, if a composite type is used for the unit datatype, atomicity is 122295fce210SBarry Smith not guaranteed across the whole vertex. Therefore, this function is mostly only used with primitive types such as 122395fce210SBarry Smith integers. 122495fce210SBarry Smith 122595fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFReduceBegin(), PetscSFSetGraph() 122695fce210SBarry Smith @*/ 122795fce210SBarry Smith PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op) 122895fce210SBarry Smith { 122995fce210SBarry Smith PetscErrorCode ierr; 123095fce210SBarry Smith 123195fce210SBarry Smith PetscFunctionBegin; 123295fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 123395fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 1234*29046d53SLisandro Dalcin ierr = PetscLogEventBegin(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr); 123595fce210SBarry Smith ierr = (*sf->ops->FetchAndOpBegin)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr); 1236563ffbabSMatthew G. Knepley ierr = PetscLogEventEnd(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr); 123795fce210SBarry Smith PetscFunctionReturn(0); 123895fce210SBarry Smith } 123995fce210SBarry Smith 124095fce210SBarry Smith /*@C 124195fce210SBarry Smith PetscSFFetchAndOpEnd - end operation started in matching call to PetscSFFetchAndOpBegin() to fetch values from roots and update atomically by applying operation using my leaf value 124295fce210SBarry Smith 124395fce210SBarry Smith Collective 124495fce210SBarry Smith 124595fce210SBarry Smith Input Arguments: 124695fce210SBarry Smith + sf - star forest 124795fce210SBarry Smith . unit - data type 124895fce210SBarry Smith . leafdata - leaf values to use in reduction 124995fce210SBarry Smith - op - operation to use for reduction 125095fce210SBarry Smith 125195fce210SBarry Smith Output Arguments: 125295fce210SBarry Smith + rootdata - root values to be updated, input state is seen by first process to perform an update 125395fce210SBarry Smith - leafupdate - state at each leaf's respective root immediately prior to my atomic update 125495fce210SBarry Smith 125595fce210SBarry Smith Level: advanced 125695fce210SBarry Smith 125795fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFReduceEnd(), PetscSFSetGraph() 125895fce210SBarry Smith @*/ 125995fce210SBarry Smith PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op) 126095fce210SBarry Smith { 126195fce210SBarry Smith PetscErrorCode ierr; 126295fce210SBarry Smith 126395fce210SBarry Smith PetscFunctionBegin; 126495fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 126595fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 1266*29046d53SLisandro Dalcin ierr = PetscLogEventBegin(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr); 126795fce210SBarry Smith ierr = (*sf->ops->FetchAndOpEnd)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr); 1268563ffbabSMatthew G. Knepley ierr = PetscLogEventEnd(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr); 126995fce210SBarry Smith PetscFunctionReturn(0); 127095fce210SBarry Smith } 127195fce210SBarry Smith 127295fce210SBarry Smith /*@C 127395fce210SBarry Smith PetscSFGatherBegin - begin pointwise gather of all leaves into multi-roots, to be completed with PetscSFGatherEnd() 127495fce210SBarry Smith 127595fce210SBarry Smith Collective 127695fce210SBarry Smith 127795fce210SBarry Smith Input Arguments: 127895fce210SBarry Smith + sf - star forest 127995fce210SBarry Smith . unit - data type 128095fce210SBarry Smith - leafdata - leaf data to gather to roots 128195fce210SBarry Smith 128295fce210SBarry Smith Output Argument: 128395fce210SBarry Smith . multirootdata - root buffer to gather into, amount of space per root is equal to its degree 128495fce210SBarry Smith 128595fce210SBarry Smith Level: intermediate 128695fce210SBarry Smith 128795fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin() 128895fce210SBarry Smith @*/ 128995fce210SBarry Smith PetscErrorCode PetscSFGatherBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata) 129095fce210SBarry Smith { 129195fce210SBarry Smith PetscErrorCode ierr; 129295fce210SBarry Smith PetscSF multi; 129395fce210SBarry Smith 129495fce210SBarry Smith PetscFunctionBegin; 129595fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 1296*29046d53SLisandro Dalcin ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 129795fce210SBarry Smith ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr); 12988bfbc91cSJed Brown ierr = PetscSFReduceBegin(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr); 129995fce210SBarry Smith PetscFunctionReturn(0); 130095fce210SBarry Smith } 130195fce210SBarry Smith 130295fce210SBarry Smith /*@C 130395fce210SBarry Smith PetscSFGatherEnd - ends pointwise gather operation that was started with PetscSFGatherBegin() 130495fce210SBarry Smith 130595fce210SBarry Smith Collective 130695fce210SBarry Smith 130795fce210SBarry Smith Input Arguments: 130895fce210SBarry Smith + sf - star forest 130995fce210SBarry Smith . unit - data type 131095fce210SBarry Smith - leafdata - leaf data to gather to roots 131195fce210SBarry Smith 131295fce210SBarry Smith Output Argument: 131395fce210SBarry Smith . multirootdata - root buffer to gather into, amount of space per root is equal to its degree 131495fce210SBarry Smith 131595fce210SBarry Smith Level: intermediate 131695fce210SBarry Smith 131795fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd() 131895fce210SBarry Smith @*/ 131995fce210SBarry Smith PetscErrorCode PetscSFGatherEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata) 132095fce210SBarry Smith { 132195fce210SBarry Smith PetscErrorCode ierr; 132295fce210SBarry Smith PetscSF multi; 132395fce210SBarry Smith 132495fce210SBarry Smith PetscFunctionBegin; 132595fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 132695fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 132795fce210SBarry Smith ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr); 13288bfbc91cSJed Brown ierr = PetscSFReduceEnd(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr); 132995fce210SBarry Smith PetscFunctionReturn(0); 133095fce210SBarry Smith } 133195fce210SBarry Smith 133295fce210SBarry Smith /*@C 133395fce210SBarry Smith PetscSFScatterBegin - begin pointwise scatter operation from multi-roots to leaves, to be completed with PetscSFScatterEnd() 133495fce210SBarry Smith 133595fce210SBarry Smith Collective 133695fce210SBarry Smith 133795fce210SBarry Smith Input Arguments: 133895fce210SBarry Smith + sf - star forest 133995fce210SBarry Smith . unit - data type 134095fce210SBarry Smith - multirootdata - root buffer to send to each leaf, one unit of data per leaf 134195fce210SBarry Smith 134295fce210SBarry Smith Output Argument: 134395fce210SBarry Smith . leafdata - leaf data to be update with personal data from each respective root 134495fce210SBarry Smith 134595fce210SBarry Smith Level: intermediate 134695fce210SBarry Smith 134795fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin() 134895fce210SBarry Smith @*/ 134995fce210SBarry Smith PetscErrorCode PetscSFScatterBegin(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata) 135095fce210SBarry Smith { 135195fce210SBarry Smith PetscErrorCode ierr; 135295fce210SBarry Smith PetscSF multi; 135395fce210SBarry Smith 135495fce210SBarry Smith PetscFunctionBegin; 135595fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 135695fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 135795fce210SBarry Smith ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr); 135895fce210SBarry Smith ierr = PetscSFBcastBegin(multi,unit,multirootdata,leafdata);CHKERRQ(ierr); 135995fce210SBarry Smith PetscFunctionReturn(0); 136095fce210SBarry Smith } 136195fce210SBarry Smith 136295fce210SBarry Smith /*@C 136395fce210SBarry Smith PetscSFScatterEnd - ends pointwise scatter operation that was started with PetscSFScatterBegin() 136495fce210SBarry Smith 136595fce210SBarry Smith Collective 136695fce210SBarry Smith 136795fce210SBarry Smith Input Arguments: 136895fce210SBarry Smith + sf - star forest 136995fce210SBarry Smith . unit - data type 137095fce210SBarry Smith - multirootdata - root buffer to send to each leaf, one unit of data per leaf 137195fce210SBarry Smith 137295fce210SBarry Smith Output Argument: 137395fce210SBarry Smith . leafdata - leaf data to be update with personal data from each respective root 137495fce210SBarry Smith 137595fce210SBarry Smith Level: intermediate 137695fce210SBarry Smith 137795fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd() 137895fce210SBarry Smith @*/ 137995fce210SBarry Smith PetscErrorCode PetscSFScatterEnd(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata) 138095fce210SBarry Smith { 138195fce210SBarry Smith PetscErrorCode ierr; 138295fce210SBarry Smith PetscSF multi; 138395fce210SBarry Smith 138495fce210SBarry Smith PetscFunctionBegin; 138595fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 138695fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 138795fce210SBarry Smith ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr); 138895fce210SBarry Smith ierr = PetscSFBcastEnd(multi,unit,multirootdata,leafdata);CHKERRQ(ierr); 138995fce210SBarry Smith PetscFunctionReturn(0); 139095fce210SBarry Smith } 1391a7b3aa13SAta Mesgarnejad 1392a7b3aa13SAta Mesgarnejad /*@ 1393a7b3aa13SAta Mesgarnejad PetscSFCompose - Compose a new PetscSF equivalent to action to PetscSFs 1394a7b3aa13SAta Mesgarnejad 1395a7b3aa13SAta Mesgarnejad Input Parameters: 1396a7b3aa13SAta Mesgarnejad + sfA - The first PetscSF 1397a7b3aa13SAta Mesgarnejad - sfB - The second PetscSF 1398a7b3aa13SAta Mesgarnejad 1399a7b3aa13SAta Mesgarnejad Output Parameters: 1400a7b3aa13SAta Mesgarnejad . sfBA - equvalent PetscSF for applying A then B 1401a7b3aa13SAta Mesgarnejad 1402a7b3aa13SAta Mesgarnejad Level: developer 1403a7b3aa13SAta Mesgarnejad 1404a7b3aa13SAta Mesgarnejad .seealso: PetscSF, PetscSFGetGraph(), PetscSFSetGraph() 1405a7b3aa13SAta Mesgarnejad @*/ 1406a7b3aa13SAta Mesgarnejad PetscErrorCode PetscSFCompose(PetscSF sfA, PetscSF sfB, PetscSF *sfBA) 1407a7b3aa13SAta Mesgarnejad { 1408a7b3aa13SAta Mesgarnejad MPI_Comm comm; 1409a7b3aa13SAta Mesgarnejad const PetscSFNode *remotePointsA, *remotePointsB; 1410a7b3aa13SAta Mesgarnejad PetscSFNode *remotePointsBA; 1411a7b3aa13SAta Mesgarnejad const PetscInt *localPointsA, *localPointsB; 1412a7b3aa13SAta Mesgarnejad PetscInt numRootsA, numLeavesA, numRootsB, numLeavesB; 1413a7b3aa13SAta Mesgarnejad PetscErrorCode ierr; 1414a7b3aa13SAta Mesgarnejad 1415a7b3aa13SAta Mesgarnejad PetscFunctionBegin; 1416a7b3aa13SAta Mesgarnejad PetscValidHeaderSpecific(sfA, PETSCSF_CLASSID, 1); 1417*29046d53SLisandro Dalcin PetscSFCheckGraphSet(sfA, 1); 1418*29046d53SLisandro Dalcin PetscValidHeaderSpecific(sfB, PETSCSF_CLASSID, 2); 1419*29046d53SLisandro Dalcin PetscSFCheckGraphSet(sfB, 2); 1420*29046d53SLisandro Dalcin PetscValidPointer(sfBA, 3); 1421a7b3aa13SAta Mesgarnejad ierr = PetscObjectGetComm((PetscObject) sfA, &comm);CHKERRQ(ierr); 1422a7b3aa13SAta Mesgarnejad ierr = PetscSFGetGraph(sfA, &numRootsA, &numLeavesA, &localPointsA, &remotePointsA);CHKERRQ(ierr); 1423a7b3aa13SAta Mesgarnejad ierr = PetscSFGetGraph(sfB, &numRootsB, &numLeavesB, &localPointsB, &remotePointsB);CHKERRQ(ierr); 1424a7b3aa13SAta Mesgarnejad ierr = PetscMalloc1(numLeavesB, &remotePointsBA);CHKERRQ(ierr); 1425a7b3aa13SAta Mesgarnejad ierr = PetscSFBcastBegin(sfB, MPIU_2INT, remotePointsA, remotePointsBA);CHKERRQ(ierr); 1426a7b3aa13SAta Mesgarnejad ierr = PetscSFBcastEnd(sfB, MPIU_2INT, remotePointsA, remotePointsBA);CHKERRQ(ierr); 1427a7b3aa13SAta Mesgarnejad ierr = PetscSFCreate(comm, sfBA);CHKERRQ(ierr); 1428a7b3aa13SAta Mesgarnejad ierr = PetscSFSetGraph(*sfBA, numRootsA, numLeavesB, localPointsB, PETSC_COPY_VALUES, remotePointsBA, PETSC_OWN_POINTER);CHKERRQ(ierr); 1429a7b3aa13SAta Mesgarnejad PetscFunctionReturn(0); 1430a7b3aa13SAta Mesgarnejad } 1431