1*95fce210SBarry Smith #include <petsc-private/sfimpl.h> /*I "petscsf.h" I*/ 2*95fce210SBarry Smith #include <petscctable.h> 3*95fce210SBarry Smith 4*95fce210SBarry Smith #if defined(PETSC_USE_DEBUG) 5*95fce210SBarry Smith # define PetscSFCheckGraphSet(sf,arg) do { \ 6*95fce210SBarry Smith if (PetscUnlikely(!(sf)->graphset)) \ 7*95fce210SBarry Smith SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetGraph() on argument %D \"%s\" before %s()",(arg),#sf,PETSC_FUNCTION_NAME); \ 8*95fce210SBarry Smith } while (0) 9*95fce210SBarry Smith #else 10*95fce210SBarry Smith # define PetscSFCheckGraphSet(sf,arg) do {} while (0) 11*95fce210SBarry Smith #endif 12*95fce210SBarry Smith 13*95fce210SBarry Smith const char *const PetscSFDuplicateOptions[] = {"CONFONLY","RANKS","GRAPH","PetscSFDuplicateOption","PETSCSF_DUPLICATE_",0}; 14*95fce210SBarry Smith 15*95fce210SBarry Smith #undef __FUNCT__ 16*95fce210SBarry Smith #define __FUNCT__ "PetscSFCreate" 17*95fce210SBarry Smith /*@C 18*95fce210SBarry Smith PetscSFCreate - create a star forest communication context 19*95fce210SBarry Smith 20*95fce210SBarry Smith Not Collective 21*95fce210SBarry Smith 22*95fce210SBarry Smith Input Arguments: 23*95fce210SBarry Smith . comm - communicator on which the star forest will operate 24*95fce210SBarry Smith 25*95fce210SBarry Smith Output Arguments: 26*95fce210SBarry Smith . sf - new star forest context 27*95fce210SBarry Smith 28*95fce210SBarry Smith Level: intermediate 29*95fce210SBarry Smith 30*95fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFDestroy() 31*95fce210SBarry Smith @*/ 32*95fce210SBarry Smith PetscErrorCode PetscSFCreate(MPI_Comm comm,PetscSF *sf) 33*95fce210SBarry Smith { 34*95fce210SBarry Smith PetscErrorCode ierr; 35*95fce210SBarry Smith PetscSF b; 36*95fce210SBarry Smith 37*95fce210SBarry Smith PetscFunctionBegin; 38*95fce210SBarry Smith PetscValidPointer(sf,2); 39*95fce210SBarry Smith #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 40*95fce210SBarry Smith ierr = PetscSFInitializePackage(NULL);CHKERRQ(ierr); 41*95fce210SBarry Smith #endif 42*95fce210SBarry Smith 43*95fce210SBarry Smith ierr = PetscHeaderCreate(b,_p_PetscSF,struct _PetscSFOps,PETSCSF_CLASSID,"PetscSF","Star Forest","PetscSF",comm,PetscSFDestroy,PetscSFView);CHKERRQ(ierr); 44*95fce210SBarry Smith 45*95fce210SBarry Smith b->nroots = -1; 46*95fce210SBarry Smith b->nleaves = -1; 47*95fce210SBarry Smith b->nranks = -1; 48*95fce210SBarry Smith b->rankorder = PETSC_TRUE; 49*95fce210SBarry Smith b->ingroup = MPI_GROUP_NULL; 50*95fce210SBarry Smith b->outgroup = MPI_GROUP_NULL; 51*95fce210SBarry Smith b->graphset = PETSC_FALSE; 52*95fce210SBarry Smith 53*95fce210SBarry Smith *sf = b; 54*95fce210SBarry Smith PetscFunctionReturn(0); 55*95fce210SBarry Smith } 56*95fce210SBarry Smith 57*95fce210SBarry Smith #undef __FUNCT__ 58*95fce210SBarry Smith #define __FUNCT__ "PetscSFReset" 59*95fce210SBarry Smith /*@C 60*95fce210SBarry Smith PetscSFReset - Reset a star forest so that different sizes or neighbors can be used 61*95fce210SBarry Smith 62*95fce210SBarry Smith Collective 63*95fce210SBarry Smith 64*95fce210SBarry Smith Input Arguments: 65*95fce210SBarry Smith . sf - star forest 66*95fce210SBarry Smith 67*95fce210SBarry Smith Level: advanced 68*95fce210SBarry Smith 69*95fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFDestroy() 70*95fce210SBarry Smith @*/ 71*95fce210SBarry Smith PetscErrorCode PetscSFReset(PetscSF sf) 72*95fce210SBarry Smith { 73*95fce210SBarry Smith PetscErrorCode ierr; 74*95fce210SBarry Smith 75*95fce210SBarry Smith PetscFunctionBegin; 76*95fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 77*95fce210SBarry Smith sf->mine = NULL; 78*95fce210SBarry Smith ierr = PetscFree(sf->mine_alloc);CHKERRQ(ierr); 79*95fce210SBarry Smith sf->remote = NULL; 80*95fce210SBarry Smith ierr = PetscFree(sf->remote_alloc);CHKERRQ(ierr); 81*95fce210SBarry Smith ierr = PetscFree4(sf->ranks,sf->roffset,sf->rmine,sf->rremote);CHKERRQ(ierr); 82*95fce210SBarry Smith ierr = PetscFree(sf->degree);CHKERRQ(ierr); 83*95fce210SBarry Smith if (sf->ingroup != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->ingroup);CHKERRQ(ierr);} 84*95fce210SBarry Smith if (sf->outgroup != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->outgroup);CHKERRQ(ierr);} 85*95fce210SBarry Smith ierr = PetscSFDestroy(&sf->multi);CHKERRQ(ierr); 86*95fce210SBarry Smith sf->graphset = PETSC_FALSE; 87*95fce210SBarry Smith if (sf->ops->Reset) {ierr = (*sf->ops->Reset)(sf);CHKERRQ(ierr);} 88*95fce210SBarry Smith sf->setupcalled = PETSC_FALSE; 89*95fce210SBarry Smith PetscFunctionReturn(0); 90*95fce210SBarry Smith } 91*95fce210SBarry Smith 92*95fce210SBarry Smith #undef __FUNCT__ 93*95fce210SBarry Smith #define __FUNCT__ "PetscSFSetType" 94*95fce210SBarry Smith /*@C 95*95fce210SBarry Smith PetscSFSetType - set the PetscSF communication implementation 96*95fce210SBarry Smith 97*95fce210SBarry Smith Collective on PetscSF 98*95fce210SBarry Smith 99*95fce210SBarry Smith Input Parameters: 100*95fce210SBarry Smith + sf - the PetscSF context 101*95fce210SBarry Smith - type - a known method 102*95fce210SBarry Smith 103*95fce210SBarry Smith Options Database Key: 104*95fce210SBarry Smith . -sf_type <type> - Sets the method; use -help for a list 105*95fce210SBarry Smith of available methods (for instance, window, pt2pt, neighbor) 106*95fce210SBarry Smith 107*95fce210SBarry Smith Notes: 108*95fce210SBarry Smith See "include/petscsf.h" for available methods (for instance) 109*95fce210SBarry Smith + PETSCSFWINDOW - MPI-2/3 one-sided 110*95fce210SBarry Smith - PETSCSFBASIC - basic implementation using MPI-1 two-sided 111*95fce210SBarry Smith 112*95fce210SBarry Smith Level: intermediate 113*95fce210SBarry Smith 114*95fce210SBarry Smith .keywords: PetscSF, set, type 115*95fce210SBarry Smith 116*95fce210SBarry Smith .seealso: PetscSFType, PetscSFCreate() 117*95fce210SBarry Smith @*/ 118*95fce210SBarry Smith PetscErrorCode PetscSFSetType(PetscSF sf,PetscSFType type) 119*95fce210SBarry Smith { 120*95fce210SBarry Smith PetscErrorCode ierr,(*r)(PetscSF); 121*95fce210SBarry Smith PetscBool match; 122*95fce210SBarry Smith 123*95fce210SBarry Smith PetscFunctionBegin; 124*95fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 125*95fce210SBarry Smith PetscValidCharPointer(type,2); 126*95fce210SBarry Smith 127*95fce210SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)sf,type,&match);CHKERRQ(ierr); 128*95fce210SBarry Smith if (match) PetscFunctionReturn(0); 129*95fce210SBarry Smith 130*95fce210SBarry Smith ierr = PetscFunctionListFind(PetscObjectComm((PetscObject)sf),PetscSFunctionList,type,PETSC_TRUE,(void (**)(void))&r);CHKERRQ(ierr); 131*95fce210SBarry Smith if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested PetscSF type %s",type); 132*95fce210SBarry Smith /* Destroy the previous private PetscSF context */ 133*95fce210SBarry Smith if (sf->ops->Destroy) { 134*95fce210SBarry Smith ierr = (*(sf)->ops->Destroy)(sf);CHKERRQ(ierr); 135*95fce210SBarry Smith } 136*95fce210SBarry Smith ierr = PetscMemzero(sf->ops,sizeof(*sf->ops));CHKERRQ(ierr); 137*95fce210SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)sf,type);CHKERRQ(ierr); 138*95fce210SBarry Smith ierr = (*r)(sf);CHKERRQ(ierr); 139*95fce210SBarry Smith PetscFunctionReturn(0); 140*95fce210SBarry Smith } 141*95fce210SBarry Smith 142*95fce210SBarry Smith #undef __FUNCT__ 143*95fce210SBarry Smith #define __FUNCT__ "PetscSFDestroy" 144*95fce210SBarry Smith /*@C 145*95fce210SBarry Smith PetscSFDestroy - destroy star forest 146*95fce210SBarry Smith 147*95fce210SBarry Smith Collective 148*95fce210SBarry Smith 149*95fce210SBarry Smith Input Arguments: 150*95fce210SBarry Smith . sf - address of star forest 151*95fce210SBarry Smith 152*95fce210SBarry Smith Level: intermediate 153*95fce210SBarry Smith 154*95fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFReset() 155*95fce210SBarry Smith @*/ 156*95fce210SBarry Smith PetscErrorCode PetscSFDestroy(PetscSF *sf) 157*95fce210SBarry Smith { 158*95fce210SBarry Smith PetscErrorCode ierr; 159*95fce210SBarry Smith 160*95fce210SBarry Smith PetscFunctionBegin; 161*95fce210SBarry Smith if (!*sf) PetscFunctionReturn(0); 162*95fce210SBarry Smith PetscValidHeaderSpecific((*sf),PETSCSF_CLASSID,1); 163*95fce210SBarry Smith if (--((PetscObject)(*sf))->refct > 0) {*sf = 0; PetscFunctionReturn(0);} 164*95fce210SBarry Smith ierr = PetscSFReset(*sf);CHKERRQ(ierr); 165*95fce210SBarry Smith if ((*sf)->ops->Destroy) {ierr = (*(*sf)->ops->Destroy)(*sf);CHKERRQ(ierr);} 166*95fce210SBarry Smith ierr = PetscHeaderDestroy(sf);CHKERRQ(ierr); 167*95fce210SBarry Smith PetscFunctionReturn(0); 168*95fce210SBarry Smith } 169*95fce210SBarry Smith 170*95fce210SBarry Smith #undef __FUNCT__ 171*95fce210SBarry Smith #define __FUNCT__ "PetscSFSetUp" 172*95fce210SBarry Smith /*@ 173*95fce210SBarry Smith PetscSFSetUp - set up communication structures 174*95fce210SBarry Smith 175*95fce210SBarry Smith Collective 176*95fce210SBarry Smith 177*95fce210SBarry Smith Input Arguments: 178*95fce210SBarry Smith . sf - star forest communication object 179*95fce210SBarry Smith 180*95fce210SBarry Smith Level: beginner 181*95fce210SBarry Smith 182*95fce210SBarry Smith .seealso: PetscSFSetFromOptions(), PetscSFSetType() 183*95fce210SBarry Smith @*/ 184*95fce210SBarry Smith PetscErrorCode PetscSFSetUp(PetscSF sf) 185*95fce210SBarry Smith { 186*95fce210SBarry Smith PetscErrorCode ierr; 187*95fce210SBarry Smith 188*95fce210SBarry Smith PetscFunctionBegin; 189*95fce210SBarry Smith if (sf->setupcalled) PetscFunctionReturn(0); 190*95fce210SBarry Smith if (!((PetscObject)sf)->type_name) {ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr);} 191*95fce210SBarry Smith if (sf->ops->SetUp) {ierr = (*sf->ops->SetUp)(sf);CHKERRQ(ierr);} 192*95fce210SBarry Smith sf->setupcalled = PETSC_TRUE; 193*95fce210SBarry Smith PetscFunctionReturn(0); 194*95fce210SBarry Smith } 195*95fce210SBarry Smith 196*95fce210SBarry Smith #undef __FUNCT__ 197*95fce210SBarry Smith #define __FUNCT__ "PetscSFSetFromOptions" 198*95fce210SBarry Smith /*@C 199*95fce210SBarry Smith PetscSFSetFromOptions - set PetscSF options using the options database 200*95fce210SBarry Smith 201*95fce210SBarry Smith Logically Collective 202*95fce210SBarry Smith 203*95fce210SBarry Smith Input Arguments: 204*95fce210SBarry Smith . sf - star forest 205*95fce210SBarry Smith 206*95fce210SBarry Smith Options Database Keys: 207*95fce210SBarry Smith . -sf_synchronization - synchronization type used by PetscSF 208*95fce210SBarry Smith 209*95fce210SBarry Smith Level: intermediate 210*95fce210SBarry Smith 211*95fce210SBarry Smith .keywords: KSP, set, from, options, database 212*95fce210SBarry Smith 213*95fce210SBarry Smith .seealso: PetscSFWindowSetSyncType() 214*95fce210SBarry Smith @*/ 215*95fce210SBarry Smith PetscErrorCode PetscSFSetFromOptions(PetscSF sf) 216*95fce210SBarry Smith { 217*95fce210SBarry Smith PetscSFType deft; 218*95fce210SBarry Smith char type[256]; 219*95fce210SBarry Smith PetscErrorCode ierr; 220*95fce210SBarry Smith PetscBool flg; 221*95fce210SBarry Smith 222*95fce210SBarry Smith PetscFunctionBegin; 223*95fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 224*95fce210SBarry Smith ierr = PetscObjectOptionsBegin((PetscObject)sf);CHKERRQ(ierr); 225*95fce210SBarry Smith deft = ((PetscObject)sf)->type_name ? ((PetscObject)sf)->type_name : PETSCSFBASIC; 226*95fce210SBarry Smith ierr = PetscOptionsList("-sf_type","PetscSF implementation type","PetscSFSetType",PetscSFunctionList,deft,type,256,&flg);CHKERRQ(ierr); 227*95fce210SBarry Smith ierr = PetscSFSetType(sf,flg ? type : deft);CHKERRQ(ierr); 228*95fce210SBarry 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); 229*95fce210SBarry Smith if (sf->ops->SetFromOptions) {ierr = (*sf->ops->SetFromOptions)(sf);CHKERRQ(ierr);} 230*95fce210SBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 231*95fce210SBarry Smith PetscFunctionReturn(0); 232*95fce210SBarry Smith } 233*95fce210SBarry Smith 234*95fce210SBarry Smith #undef __FUNCT__ 235*95fce210SBarry Smith #define __FUNCT__ "PetscSFSetRankOrder" 236*95fce210SBarry Smith /*@C 237*95fce210SBarry Smith PetscSFSetRankOrder - sort multi-points for gathers and scatters by rank order 238*95fce210SBarry Smith 239*95fce210SBarry Smith Logically Collective 240*95fce210SBarry Smith 241*95fce210SBarry Smith Input Arguments: 242*95fce210SBarry Smith + sf - star forest 243*95fce210SBarry Smith - flg - PETSC_TRUE to sort, PETSC_FALSE to skip sorting (lower setup cost, but non-deterministic) 244*95fce210SBarry Smith 245*95fce210SBarry Smith Level: advanced 246*95fce210SBarry Smith 247*95fce210SBarry Smith .seealso: PetscSFGatherBegin(), PetscSFScatterBegin() 248*95fce210SBarry Smith @*/ 249*95fce210SBarry Smith PetscErrorCode PetscSFSetRankOrder(PetscSF sf,PetscBool flg) 250*95fce210SBarry Smith { 251*95fce210SBarry Smith 252*95fce210SBarry Smith PetscFunctionBegin; 253*95fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 254*95fce210SBarry Smith PetscValidLogicalCollectiveBool(sf,flg,2); 255*95fce210SBarry Smith if (sf->multi) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Rank ordering must be set before first call to PetscSFGatherBegin() or PetscSFScatterBegin()"); 256*95fce210SBarry Smith sf->rankorder = flg; 257*95fce210SBarry Smith PetscFunctionReturn(0); 258*95fce210SBarry Smith } 259*95fce210SBarry Smith 260*95fce210SBarry Smith #undef __FUNCT__ 261*95fce210SBarry Smith #define __FUNCT__ "PetscSFSetGraph" 262*95fce210SBarry Smith /*@C 263*95fce210SBarry Smith PetscSFSetGraph - Set a parallel star forest 264*95fce210SBarry Smith 265*95fce210SBarry Smith Collective 266*95fce210SBarry Smith 267*95fce210SBarry Smith Input Arguments: 268*95fce210SBarry Smith + sf - star forest 269*95fce210SBarry Smith . nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves) 270*95fce210SBarry Smith . nleaves - number of leaf vertices on the current process, each of these references a root on any process 271*95fce210SBarry Smith . ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage 272*95fce210SBarry Smith . localmode - copy mode for ilocal 273*95fce210SBarry Smith . iremote - remote locations of root vertices for each leaf on the current process 274*95fce210SBarry Smith - remotemode - copy mode for iremote 275*95fce210SBarry Smith 276*95fce210SBarry Smith Level: intermediate 277*95fce210SBarry Smith 278*95fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFGetGraph() 279*95fce210SBarry Smith @*/ 280*95fce210SBarry Smith PetscErrorCode PetscSFSetGraph(PetscSF sf,PetscInt nroots,PetscInt nleaves,const PetscInt *ilocal,PetscCopyMode localmode,const PetscSFNode *iremote,PetscCopyMode remotemode) 281*95fce210SBarry Smith { 282*95fce210SBarry Smith PetscErrorCode ierr; 283*95fce210SBarry Smith PetscTable table; 284*95fce210SBarry Smith PetscTablePosition pos; 285*95fce210SBarry Smith PetscMPIInt size; 286*95fce210SBarry Smith PetscInt i,*rcount,*ranks; 287*95fce210SBarry Smith 288*95fce210SBarry Smith PetscFunctionBegin; 289*95fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 290*95fce210SBarry Smith if (nleaves && ilocal) PetscValidIntPointer(ilocal,4); 291*95fce210SBarry Smith if (nleaves) PetscValidPointer(iremote,6); 292*95fce210SBarry Smith if (nroots < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"roots %D, cannot be negative",nroots); 293*95fce210SBarry Smith if (nleaves < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nleaves %D, cannot be negative",nleaves); 294*95fce210SBarry Smith ierr = PetscSFReset(sf);CHKERRQ(ierr); 295*95fce210SBarry Smith sf->nroots = nroots; 296*95fce210SBarry Smith sf->nleaves = nleaves; 297*95fce210SBarry Smith if (ilocal) { 298*95fce210SBarry Smith switch (localmode) { 299*95fce210SBarry Smith case PETSC_COPY_VALUES: 300*95fce210SBarry Smith ierr = PetscMalloc(nleaves*sizeof(*sf->mine),&sf->mine_alloc);CHKERRQ(ierr); 301*95fce210SBarry Smith sf->mine = sf->mine_alloc; 302*95fce210SBarry Smith ierr = PetscMemcpy(sf->mine,ilocal,nleaves*sizeof(*sf->mine));CHKERRQ(ierr); 303*95fce210SBarry Smith sf->minleaf = PETSC_MAX_INT; 304*95fce210SBarry Smith sf->maxleaf = PETSC_MIN_INT; 305*95fce210SBarry Smith for (i=0; i<nleaves; i++) { 306*95fce210SBarry Smith sf->minleaf = PetscMin(sf->minleaf,ilocal[i]); 307*95fce210SBarry Smith sf->maxleaf = PetscMax(sf->maxleaf,ilocal[i]); 308*95fce210SBarry Smith } 309*95fce210SBarry Smith break; 310*95fce210SBarry Smith case PETSC_OWN_POINTER: 311*95fce210SBarry Smith sf->mine_alloc = (PetscInt*)ilocal; 312*95fce210SBarry Smith sf->mine = sf->mine_alloc; 313*95fce210SBarry Smith break; 314*95fce210SBarry Smith case PETSC_USE_POINTER: 315*95fce210SBarry Smith sf->mine = (PetscInt*)ilocal; 316*95fce210SBarry Smith break; 317*95fce210SBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown localmode"); 318*95fce210SBarry Smith } 319*95fce210SBarry Smith } 320*95fce210SBarry Smith if (!ilocal || nleaves > 0) { 321*95fce210SBarry Smith sf->minleaf = 0; 322*95fce210SBarry Smith sf->maxleaf = nleaves - 1; 323*95fce210SBarry Smith } 324*95fce210SBarry Smith switch (remotemode) { 325*95fce210SBarry Smith case PETSC_COPY_VALUES: 326*95fce210SBarry Smith ierr = PetscMalloc(nleaves*sizeof(*sf->remote),&sf->remote_alloc);CHKERRQ(ierr); 327*95fce210SBarry Smith sf->remote = sf->remote_alloc; 328*95fce210SBarry Smith ierr = PetscMemcpy(sf->remote,iremote,nleaves*sizeof(*sf->remote));CHKERRQ(ierr); 329*95fce210SBarry Smith break; 330*95fce210SBarry Smith case PETSC_OWN_POINTER: 331*95fce210SBarry Smith sf->remote_alloc = (PetscSFNode*)iremote; 332*95fce210SBarry Smith sf->remote = sf->remote_alloc; 333*95fce210SBarry Smith break; 334*95fce210SBarry Smith case PETSC_USE_POINTER: 335*95fce210SBarry Smith sf->remote = (PetscSFNode*)iremote; 336*95fce210SBarry Smith break; 337*95fce210SBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown remotemode"); 338*95fce210SBarry Smith } 339*95fce210SBarry Smith 340*95fce210SBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr); 341*95fce210SBarry Smith ierr = PetscTableCreate(10,size,&table);CHKERRQ(ierr); 342*95fce210SBarry Smith for (i=0; i<nleaves; i++) { 343*95fce210SBarry Smith /* Log 1-based rank */ 344*95fce210SBarry Smith ierr = PetscTableAdd(table,iremote[i].rank+1,1,ADD_VALUES);CHKERRQ(ierr); 345*95fce210SBarry Smith } 346*95fce210SBarry Smith ierr = PetscTableGetCount(table,&sf->nranks);CHKERRQ(ierr); 347*95fce210SBarry Smith ierr = PetscMalloc4(sf->nranks,PetscInt,&sf->ranks,sf->nranks+1,PetscInt,&sf->roffset,nleaves,PetscInt,&sf->rmine,nleaves,PetscInt,&sf->rremote);CHKERRQ(ierr); 348*95fce210SBarry Smith ierr = PetscMalloc2(sf->nranks,PetscInt,&rcount,sf->nranks,PetscInt,&ranks);CHKERRQ(ierr); 349*95fce210SBarry Smith ierr = PetscTableGetHeadPosition(table,&pos);CHKERRQ(ierr); 350*95fce210SBarry Smith for (i=0; i<sf->nranks; i++) { 351*95fce210SBarry Smith ierr = PetscTableGetNext(table,&pos,&ranks[i],&rcount[i]);CHKERRQ(ierr); 352*95fce210SBarry Smith ranks[i]--; /* Convert back to 0-based */ 353*95fce210SBarry Smith } 354*95fce210SBarry Smith ierr = PetscTableDestroy(&table);CHKERRQ(ierr); 355*95fce210SBarry Smith ierr = PetscSortIntWithArray(sf->nranks,ranks,rcount);CHKERRQ(ierr); 356*95fce210SBarry Smith sf->roffset[0] = 0; 357*95fce210SBarry Smith for (i=0; i<sf->nranks; i++) { 358*95fce210SBarry Smith ierr = PetscMPIIntCast(ranks[i],sf->ranks+i);CHKERRQ(ierr); 359*95fce210SBarry Smith sf->roffset[i+1] = sf->roffset[i] + rcount[i]; 360*95fce210SBarry Smith rcount[i] = 0; 361*95fce210SBarry Smith } 362*95fce210SBarry Smith for (i=0; i<nleaves; i++) { 363*95fce210SBarry Smith PetscInt lo,hi,irank; 364*95fce210SBarry Smith /* Search for index of iremote[i].rank in sf->ranks */ 365*95fce210SBarry Smith lo = 0; hi = sf->nranks; 366*95fce210SBarry Smith while (hi - lo > 1) { 367*95fce210SBarry Smith PetscInt mid = lo + (hi - lo)/2; 368*95fce210SBarry Smith if (iremote[i].rank < sf->ranks[mid]) hi = mid; 369*95fce210SBarry Smith else lo = mid; 370*95fce210SBarry Smith } 371*95fce210SBarry Smith if (hi - lo == 1 && iremote[i].rank == sf->ranks[lo]) irank = lo; 372*95fce210SBarry Smith else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Could not find rank %D in array",iremote[i].rank); 373*95fce210SBarry Smith sf->rmine[sf->roffset[irank] + rcount[irank]] = ilocal ? ilocal[i] : i; 374*95fce210SBarry Smith sf->rremote[sf->roffset[irank] + rcount[irank]] = iremote[i].index; 375*95fce210SBarry Smith rcount[irank]++; 376*95fce210SBarry Smith } 377*95fce210SBarry Smith ierr = PetscFree2(rcount,ranks);CHKERRQ(ierr); 378*95fce210SBarry Smith #if !defined(PETSC_USE_64BIT_INDICES) 379*95fce210SBarry Smith if (nroots == PETSC_DETERMINE) { 380*95fce210SBarry Smith /* Jed, if you have a better way to do this, put it in */ 381*95fce210SBarry Smith PetscInt *numRankLeaves, *leafOff, *leafIndices, *numRankRoots, *rootOff, *rootIndices, maxRoots = 0; 382*95fce210SBarry Smith 383*95fce210SBarry Smith /* All to all to determine number of leaf indices from each (you can do this using Scan and asynch messages) */ 384*95fce210SBarry Smith ierr = PetscMalloc4(size,PetscInt,&numRankLeaves,size+1,PetscInt,&leafOff,size,PetscInt,&numRankRoots,size+1,PetscInt,&rootOff);CHKERRQ(ierr); 385*95fce210SBarry Smith ierr = PetscMemzero(numRankLeaves, size * sizeof(PetscInt));CHKERRQ(ierr); 386*95fce210SBarry Smith for (i = 0; i < nleaves; ++i) ++numRankLeaves[iremote[i].rank]; 387*95fce210SBarry Smith ierr = MPI_Alltoall(numRankLeaves, 1, MPIU_INT, numRankRoots, 1, MPIU_INT, PetscObjectComm((PetscObject)sf));CHKERRQ(ierr); 388*95fce210SBarry Smith /* Could set nroots to this maximum */ 389*95fce210SBarry Smith for (i = 0; i < size; ++i) maxRoots += numRankRoots[i]; 390*95fce210SBarry Smith 391*95fce210SBarry Smith /* Gather all indices */ 392*95fce210SBarry Smith ierr = PetscMalloc2(nleaves,PetscInt,&leafIndices,maxRoots,PetscInt,&rootIndices);CHKERRQ(ierr); 393*95fce210SBarry Smith leafOff[0] = 0; 394*95fce210SBarry Smith for (i = 0; i < size; ++i) leafOff[i+1] = leafOff[i] + numRankLeaves[i]; 395*95fce210SBarry Smith for (i = 0; i < nleaves; ++i) leafIndices[leafOff[iremote[i].rank]++] = iremote[i].index; 396*95fce210SBarry Smith leafOff[0] = 0; 397*95fce210SBarry Smith for (i = 0; i < size; ++i) leafOff[i+1] = leafOff[i] + numRankLeaves[i]; 398*95fce210SBarry Smith rootOff[0] = 0; 399*95fce210SBarry Smith for (i = 0; i < size; ++i) rootOff[i+1] = rootOff[i] + numRankRoots[i]; 400*95fce210SBarry Smith ierr = MPI_Alltoallv(leafIndices, numRankLeaves, leafOff, MPIU_INT, rootIndices, numRankRoots, rootOff, MPIU_INT, PetscObjectComm((PetscObject)sf));CHKERRQ(ierr); 401*95fce210SBarry Smith /* Sort and reduce */ 402*95fce210SBarry Smith ierr = PetscSortRemoveDupsInt(&maxRoots, rootIndices);CHKERRQ(ierr); 403*95fce210SBarry Smith ierr = PetscFree2(leafIndices,rootIndices);CHKERRQ(ierr); 404*95fce210SBarry Smith ierr = PetscFree4(numRankLeaves,leafOff,numRankRoots,rootOff);CHKERRQ(ierr); 405*95fce210SBarry Smith sf->nroots = maxRoots; 406*95fce210SBarry Smith } 407*95fce210SBarry Smith #endif 408*95fce210SBarry Smith 409*95fce210SBarry Smith sf->graphset = PETSC_TRUE; 410*95fce210SBarry Smith PetscFunctionReturn(0); 411*95fce210SBarry Smith } 412*95fce210SBarry Smith 413*95fce210SBarry Smith #undef __FUNCT__ 414*95fce210SBarry Smith #define __FUNCT__ "PetscSFCreateInverseSF" 415*95fce210SBarry Smith /*@C 416*95fce210SBarry Smith PetscSFCreateInverseSF - given a PetscSF in which all vertices have degree 1, creates the inverse map 417*95fce210SBarry Smith 418*95fce210SBarry Smith Collective 419*95fce210SBarry Smith 420*95fce210SBarry Smith Input Arguments: 421*95fce210SBarry Smith . sf - star forest to invert 422*95fce210SBarry Smith 423*95fce210SBarry Smith Output Arguments: 424*95fce210SBarry Smith . isf - inverse of sf 425*95fce210SBarry Smith 426*95fce210SBarry Smith Level: advanced 427*95fce210SBarry Smith 428*95fce210SBarry Smith Notes: 429*95fce210SBarry Smith All roots must have degree 1. 430*95fce210SBarry Smith 431*95fce210SBarry Smith The local space may be a permutation, but cannot be sparse. 432*95fce210SBarry Smith 433*95fce210SBarry Smith .seealso: PetscSFSetGraph() 434*95fce210SBarry Smith @*/ 435*95fce210SBarry Smith PetscErrorCode PetscSFCreateInverseSF(PetscSF sf,PetscSF *isf) 436*95fce210SBarry Smith { 437*95fce210SBarry Smith PetscErrorCode ierr; 438*95fce210SBarry Smith PetscMPIInt rank; 439*95fce210SBarry Smith PetscInt i,nroots,nleaves,maxlocal,count,*newilocal; 440*95fce210SBarry Smith const PetscInt *ilocal; 441*95fce210SBarry Smith PetscSFNode *roots,*leaves; 442*95fce210SBarry Smith 443*95fce210SBarry Smith PetscFunctionBegin; 444*95fce210SBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); 445*95fce210SBarry Smith ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);CHKERRQ(ierr); 446*95fce210SBarry Smith for (i=0,maxlocal=0; i<nleaves; i++) maxlocal = PetscMax(maxlocal,(ilocal ? ilocal[i] : i)+1); 447*95fce210SBarry Smith ierr = PetscMalloc2(nroots,PetscSFNode,&roots,nleaves,PetscSFNode,&leaves);CHKERRQ(ierr); 448*95fce210SBarry Smith for (i=0; i<nleaves; i++) { 449*95fce210SBarry Smith leaves[i].rank = rank; 450*95fce210SBarry Smith leaves[i].index = i; 451*95fce210SBarry Smith } 452*95fce210SBarry Smith for (i=0; i <nroots; i++) { 453*95fce210SBarry Smith roots[i].rank = -1; 454*95fce210SBarry Smith roots[i].index = -1; 455*95fce210SBarry Smith } 456*95fce210SBarry Smith ierr = PetscSFReduceBegin(sf,MPIU_2INT,leaves,roots,MPI_REPLACE);CHKERRQ(ierr); 457*95fce210SBarry Smith ierr = PetscSFReduceEnd(sf,MPIU_2INT,leaves,roots,MPI_REPLACE);CHKERRQ(ierr); 458*95fce210SBarry Smith 459*95fce210SBarry Smith /* Check whether our leaves are sparse */ 460*95fce210SBarry Smith for (i=0,count=0; i<nroots; i++) if (roots[i].rank >= 0) count++; 461*95fce210SBarry Smith if (count == nroots) newilocal = NULL; 462*95fce210SBarry Smith else { /* Index for sparse leaves and compact "roots" array (which is to become our leaves). */ 463*95fce210SBarry Smith ierr = PetscMalloc(count*sizeof(PetscInt),&newilocal);CHKERRQ(ierr); 464*95fce210SBarry Smith for (i=0,count=0; i<nroots; i++) { 465*95fce210SBarry Smith if (roots[i].rank >= 0) { 466*95fce210SBarry Smith newilocal[count] = i; 467*95fce210SBarry Smith roots[count].rank = roots[i].rank; 468*95fce210SBarry Smith roots[count].index = roots[i].index; 469*95fce210SBarry Smith count++; 470*95fce210SBarry Smith } 471*95fce210SBarry Smith } 472*95fce210SBarry Smith } 473*95fce210SBarry Smith 474*95fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,isf);CHKERRQ(ierr); 475*95fce210SBarry Smith ierr = PetscSFSetGraph(*isf,maxlocal,count,newilocal,PETSC_OWN_POINTER,roots,PETSC_COPY_VALUES);CHKERRQ(ierr); 476*95fce210SBarry Smith ierr = PetscFree2(roots,leaves);CHKERRQ(ierr); 477*95fce210SBarry Smith PetscFunctionReturn(0); 478*95fce210SBarry Smith } 479*95fce210SBarry Smith 480*95fce210SBarry Smith #undef __FUNCT__ 481*95fce210SBarry Smith #define __FUNCT__ "PetscSFDuplicate" 482*95fce210SBarry Smith /*@ 483*95fce210SBarry Smith PetscSFDuplicate - duplicate a PetscSF, optionally preserving rank connectivity and graph 484*95fce210SBarry Smith 485*95fce210SBarry Smith Collective 486*95fce210SBarry Smith 487*95fce210SBarry Smith Input Arguments: 488*95fce210SBarry Smith + sf - communication object to duplicate 489*95fce210SBarry Smith - opt - PETSCSF_DUPLICATE_CONFONLY, PETSCSF_DUPLICATE_RANKS, or PETSCSF_DUPLICATE_GRAPH (see PetscSFDuplicateOption) 490*95fce210SBarry Smith 491*95fce210SBarry Smith Output Arguments: 492*95fce210SBarry Smith . newsf - new communication object 493*95fce210SBarry Smith 494*95fce210SBarry Smith Level: beginner 495*95fce210SBarry Smith 496*95fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetType(), PetscSFSetGraph() 497*95fce210SBarry Smith @*/ 498*95fce210SBarry Smith PetscErrorCode PetscSFDuplicate(PetscSF sf,PetscSFDuplicateOption opt,PetscSF *newsf) 499*95fce210SBarry Smith { 500*95fce210SBarry Smith PetscErrorCode ierr; 501*95fce210SBarry Smith 502*95fce210SBarry Smith PetscFunctionBegin; 503*95fce210SBarry Smith ierr = PetscSFCreate(PetscObjectComm((PetscObject)sf),newsf);CHKERRQ(ierr); 504*95fce210SBarry Smith ierr = PetscSFSetType(*newsf,((PetscObject)sf)->type_name);CHKERRQ(ierr); 505*95fce210SBarry Smith if (sf->ops->Duplicate) {ierr = (*sf->ops->Duplicate)(sf,opt,*newsf);CHKERRQ(ierr);} 506*95fce210SBarry Smith if (opt == PETSCSF_DUPLICATE_GRAPH) { 507*95fce210SBarry Smith PetscInt nroots,nleaves; 508*95fce210SBarry Smith const PetscInt *ilocal; 509*95fce210SBarry Smith const PetscSFNode *iremote; 510*95fce210SBarry Smith ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr); 511*95fce210SBarry Smith ierr = PetscSFSetGraph(*newsf,nroots,nleaves,ilocal,PETSC_COPY_VALUES,iremote,PETSC_COPY_VALUES);CHKERRQ(ierr); 512*95fce210SBarry Smith } 513*95fce210SBarry Smith PetscFunctionReturn(0); 514*95fce210SBarry Smith } 515*95fce210SBarry Smith 516*95fce210SBarry Smith #undef __FUNCT__ 517*95fce210SBarry Smith #define __FUNCT__ "PetscSFGetGraph" 518*95fce210SBarry Smith /*@C 519*95fce210SBarry Smith PetscSFGetGraph - Get the graph specifying a parallel star forest 520*95fce210SBarry Smith 521*95fce210SBarry Smith Not Collective 522*95fce210SBarry Smith 523*95fce210SBarry Smith Input Arguments: 524*95fce210SBarry Smith . sf - star forest 525*95fce210SBarry Smith 526*95fce210SBarry Smith Output Arguments: 527*95fce210SBarry Smith + nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves) 528*95fce210SBarry Smith . nleaves - number of leaf vertices on the current process, each of these references a root on any process 529*95fce210SBarry Smith . ilocal - locations of leaves in leafdata buffers 530*95fce210SBarry Smith - iremote - remote locations of root vertices for each leaf on the current process 531*95fce210SBarry Smith 532*95fce210SBarry Smith Level: intermediate 533*95fce210SBarry Smith 534*95fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph() 535*95fce210SBarry Smith @*/ 536*95fce210SBarry Smith PetscErrorCode PetscSFGetGraph(PetscSF sf,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote) 537*95fce210SBarry Smith { 538*95fce210SBarry Smith 539*95fce210SBarry Smith PetscFunctionBegin; 540*95fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 541*95fce210SBarry Smith /* We are not currently requiring that the graph is set, thus returning nroots=-1 if it has not been set */ 542*95fce210SBarry Smith /* if (!sf->graphset) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Graph has not been set, must call PetscSFSetGraph()"); */ 543*95fce210SBarry Smith if (nroots) *nroots = sf->nroots; 544*95fce210SBarry Smith if (nleaves) *nleaves = sf->nleaves; 545*95fce210SBarry Smith if (ilocal) *ilocal = sf->mine; 546*95fce210SBarry Smith if (iremote) *iremote = sf->remote; 547*95fce210SBarry Smith PetscFunctionReturn(0); 548*95fce210SBarry Smith } 549*95fce210SBarry Smith 550*95fce210SBarry Smith #undef __FUNCT__ 551*95fce210SBarry Smith #define __FUNCT__ "PetscSFGetLeafRange" 552*95fce210SBarry Smith /*@C 553*95fce210SBarry Smith PetscSFGetLeafRange - Get the active leaf ranges 554*95fce210SBarry Smith 555*95fce210SBarry Smith Not Collective 556*95fce210SBarry Smith 557*95fce210SBarry Smith Input Arguments: 558*95fce210SBarry Smith . sf - star forest 559*95fce210SBarry Smith 560*95fce210SBarry Smith Output Arguments: 561*95fce210SBarry Smith + minleaf - minimum active leaf on this process 562*95fce210SBarry Smith - maxleaf - maximum active leaf on this process 563*95fce210SBarry Smith 564*95fce210SBarry Smith Level: developer 565*95fce210SBarry Smith 566*95fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph(), PetscSFGetGraph() 567*95fce210SBarry Smith @*/ 568*95fce210SBarry Smith PetscErrorCode PetscSFGetLeafRange(PetscSF sf,PetscInt *minleaf,PetscInt *maxleaf) 569*95fce210SBarry Smith { 570*95fce210SBarry Smith 571*95fce210SBarry Smith PetscFunctionBegin; 572*95fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 573*95fce210SBarry Smith if (minleaf) *minleaf = sf->minleaf; 574*95fce210SBarry Smith if (maxleaf) *maxleaf = sf->maxleaf; 575*95fce210SBarry Smith PetscFunctionReturn(0); 576*95fce210SBarry Smith } 577*95fce210SBarry Smith 578*95fce210SBarry Smith #undef __FUNCT__ 579*95fce210SBarry Smith #define __FUNCT__ "PetscSFView" 580*95fce210SBarry Smith /*@C 581*95fce210SBarry Smith PetscSFView - view a star forest 582*95fce210SBarry Smith 583*95fce210SBarry Smith Collective 584*95fce210SBarry Smith 585*95fce210SBarry Smith Input Arguments: 586*95fce210SBarry Smith + sf - star forest 587*95fce210SBarry Smith - viewer - viewer to display graph, for example PETSC_VIEWER_STDOUT_WORLD 588*95fce210SBarry Smith 589*95fce210SBarry Smith Level: beginner 590*95fce210SBarry Smith 591*95fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph() 592*95fce210SBarry Smith @*/ 593*95fce210SBarry Smith PetscErrorCode PetscSFView(PetscSF sf,PetscViewer viewer) 594*95fce210SBarry Smith { 595*95fce210SBarry Smith PetscErrorCode ierr; 596*95fce210SBarry Smith PetscBool iascii; 597*95fce210SBarry Smith PetscViewerFormat format; 598*95fce210SBarry Smith 599*95fce210SBarry Smith PetscFunctionBegin; 600*95fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 601*95fce210SBarry Smith if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sf),&viewer);CHKERRQ(ierr);} 602*95fce210SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 603*95fce210SBarry Smith PetscCheckSameComm(sf,1,viewer,2); 604*95fce210SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 605*95fce210SBarry Smith if (iascii) { 606*95fce210SBarry Smith PetscMPIInt rank; 607*95fce210SBarry Smith PetscInt i,j; 608*95fce210SBarry Smith 609*95fce210SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)sf,viewer,"Star Forest Object");CHKERRQ(ierr); 610*95fce210SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 611*95fce210SBarry Smith if (sf->ops->View) {ierr = (*sf->ops->View)(sf,viewer);CHKERRQ(ierr);} 612*95fce210SBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); 613*95fce210SBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 614*95fce210SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of roots=%D, leaves=%D, remote ranks=%D\n",rank,sf->nroots,sf->nleaves,sf->nranks);CHKERRQ(ierr); 615*95fce210SBarry Smith for (i=0; i<sf->nleaves; i++) { 616*95fce210SBarry 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); 617*95fce210SBarry Smith } 618*95fce210SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 619*95fce210SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 620*95fce210SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 621*95fce210SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Roots referenced by my leaves, by rank\n",rank);CHKERRQ(ierr); 622*95fce210SBarry Smith for (i=0; i<sf->nranks; i++) { 623*95fce210SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D: %D edges\n",rank,sf->ranks[i],sf->roffset[i+1]-sf->roffset[i]);CHKERRQ(ierr); 624*95fce210SBarry Smith for (j=sf->roffset[i]; j<sf->roffset[i+1]; j++) { 625*95fce210SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D <- %D\n",rank,sf->rmine[j],sf->rremote[j]);CHKERRQ(ierr); 626*95fce210SBarry Smith } 627*95fce210SBarry Smith } 628*95fce210SBarry Smith } 629*95fce210SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 630*95fce210SBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 631*95fce210SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 632*95fce210SBarry Smith } 633*95fce210SBarry Smith PetscFunctionReturn(0); 634*95fce210SBarry Smith } 635*95fce210SBarry Smith 636*95fce210SBarry Smith #undef __FUNCT__ 637*95fce210SBarry Smith #define __FUNCT__ "PetscSFGetRanks" 638*95fce210SBarry Smith /*@C 639*95fce210SBarry Smith PetscSFGetRanks - Get ranks and number of vertices referenced by leaves on this process 640*95fce210SBarry Smith 641*95fce210SBarry Smith Not Collective 642*95fce210SBarry Smith 643*95fce210SBarry Smith Input Arguments: 644*95fce210SBarry Smith . sf - star forest 645*95fce210SBarry Smith 646*95fce210SBarry Smith Output Arguments: 647*95fce210SBarry Smith + nranks - number of ranks referenced by local part 648*95fce210SBarry Smith . ranks - array of ranks 649*95fce210SBarry Smith . roffset - offset in rmine/rremote for each rank (length nranks+1) 650*95fce210SBarry Smith . rmine - concatenated array holding local indices referencing each remote rank 651*95fce210SBarry Smith - rremote - concatenated array holding remote indices referenced for each remote rank 652*95fce210SBarry Smith 653*95fce210SBarry Smith Level: developer 654*95fce210SBarry Smith 655*95fce210SBarry Smith .seealso: PetscSFSetGraph() 656*95fce210SBarry Smith @*/ 657*95fce210SBarry Smith PetscErrorCode PetscSFGetRanks(PetscSF sf,PetscInt *nranks,const PetscMPIInt **ranks,const PetscInt **roffset,const PetscInt **rmine,const PetscInt **rremote) 658*95fce210SBarry Smith { 659*95fce210SBarry Smith 660*95fce210SBarry Smith PetscFunctionBegin; 661*95fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 662*95fce210SBarry Smith if (nranks) *nranks = sf->nranks; 663*95fce210SBarry Smith if (ranks) *ranks = sf->ranks; 664*95fce210SBarry Smith if (roffset) *roffset = sf->roffset; 665*95fce210SBarry Smith if (rmine) *rmine = sf->rmine; 666*95fce210SBarry Smith if (rremote) *rremote = sf->rremote; 667*95fce210SBarry Smith PetscFunctionReturn(0); 668*95fce210SBarry Smith } 669*95fce210SBarry Smith 670*95fce210SBarry Smith #undef __FUNCT__ 671*95fce210SBarry Smith #define __FUNCT__ "PetscSFGetGroups" 672*95fce210SBarry Smith /*@C 673*95fce210SBarry Smith PetscSFGetGroups - gets incoming and outgoing process groups 674*95fce210SBarry Smith 675*95fce210SBarry Smith Collective 676*95fce210SBarry Smith 677*95fce210SBarry Smith Input Argument: 678*95fce210SBarry Smith . sf - star forest 679*95fce210SBarry Smith 680*95fce210SBarry Smith Output Arguments: 681*95fce210SBarry Smith + incoming - group of origin processes for incoming edges (leaves that reference my roots) 682*95fce210SBarry Smith - outgoing - group of destination processes for outgoing edges (roots that I reference) 683*95fce210SBarry Smith 684*95fce210SBarry Smith Level: developer 685*95fce210SBarry Smith 686*95fce210SBarry Smith .seealso: PetscSFGetWindow(), PetscSFRestoreWindow() 687*95fce210SBarry Smith @*/ 688*95fce210SBarry Smith PetscErrorCode PetscSFGetGroups(PetscSF sf,MPI_Group *incoming,MPI_Group *outgoing) 689*95fce210SBarry Smith { 690*95fce210SBarry Smith PetscErrorCode ierr; 691*95fce210SBarry Smith MPI_Group group; 692*95fce210SBarry Smith 693*95fce210SBarry Smith PetscFunctionBegin; 694*95fce210SBarry Smith if (sf->ingroup == MPI_GROUP_NULL) { 695*95fce210SBarry Smith PetscInt i; 696*95fce210SBarry Smith const PetscInt *indegree; 697*95fce210SBarry Smith PetscMPIInt rank,*outranks,*inranks; 698*95fce210SBarry Smith PetscSFNode *remote; 699*95fce210SBarry Smith PetscSF bgcount; 700*95fce210SBarry Smith 701*95fce210SBarry Smith /* Compute the number of incoming ranks */ 702*95fce210SBarry Smith ierr = PetscMalloc(sf->nranks*sizeof(PetscSFNode),&remote);CHKERRQ(ierr); 703*95fce210SBarry Smith for (i=0; i<sf->nranks; i++) { 704*95fce210SBarry Smith remote[i].rank = sf->ranks[i]; 705*95fce210SBarry Smith remote[i].index = 0; 706*95fce210SBarry Smith } 707*95fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,&bgcount);CHKERRQ(ierr); 708*95fce210SBarry Smith ierr = PetscSFSetGraph(bgcount,1,sf->nranks,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr); 709*95fce210SBarry Smith ierr = PetscSFComputeDegreeBegin(bgcount,&indegree);CHKERRQ(ierr); 710*95fce210SBarry Smith ierr = PetscSFComputeDegreeEnd(bgcount,&indegree);CHKERRQ(ierr); 711*95fce210SBarry Smith 712*95fce210SBarry Smith /* Enumerate the incoming ranks */ 713*95fce210SBarry Smith ierr = PetscMalloc2(indegree[0],PetscMPIInt,&inranks,sf->nranks,PetscMPIInt,&outranks);CHKERRQ(ierr); 714*95fce210SBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); 715*95fce210SBarry Smith for (i=0; i<sf->nranks; i++) outranks[i] = rank; 716*95fce210SBarry Smith ierr = PetscSFGatherBegin(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr); 717*95fce210SBarry Smith ierr = PetscSFGatherEnd(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr); 718*95fce210SBarry Smith ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr); 719*95fce210SBarry Smith ierr = MPI_Group_incl(group,indegree[0],inranks,&sf->ingroup);CHKERRQ(ierr); 720*95fce210SBarry Smith ierr = MPI_Group_free(&group);CHKERRQ(ierr); 721*95fce210SBarry Smith ierr = PetscFree2(inranks,outranks);CHKERRQ(ierr); 722*95fce210SBarry Smith ierr = PetscSFDestroy(&bgcount);CHKERRQ(ierr); 723*95fce210SBarry Smith } 724*95fce210SBarry Smith *incoming = sf->ingroup; 725*95fce210SBarry Smith 726*95fce210SBarry Smith if (sf->outgroup == MPI_GROUP_NULL) { 727*95fce210SBarry Smith ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr); 728*95fce210SBarry Smith ierr = MPI_Group_incl(group,sf->nranks,sf->ranks,&sf->outgroup);CHKERRQ(ierr); 729*95fce210SBarry Smith ierr = MPI_Group_free(&group);CHKERRQ(ierr); 730*95fce210SBarry Smith } 731*95fce210SBarry Smith *outgoing = sf->outgroup; 732*95fce210SBarry Smith PetscFunctionReturn(0); 733*95fce210SBarry Smith } 734*95fce210SBarry Smith 735*95fce210SBarry Smith #undef __FUNCT__ 736*95fce210SBarry Smith #define __FUNCT__ "PetscSFGetMultiSF" 737*95fce210SBarry Smith /*@C 738*95fce210SBarry Smith PetscSFGetMultiSF - gets the inner SF implemeting gathers and scatters 739*95fce210SBarry Smith 740*95fce210SBarry Smith Collective 741*95fce210SBarry Smith 742*95fce210SBarry Smith Input Argument: 743*95fce210SBarry Smith . sf - star forest that may contain roots with 0 or with more than 1 vertex 744*95fce210SBarry Smith 745*95fce210SBarry Smith Output Arguments: 746*95fce210SBarry Smith . multi - star forest with split roots, such that each root has degree exactly 1 747*95fce210SBarry Smith 748*95fce210SBarry Smith Level: developer 749*95fce210SBarry Smith 750*95fce210SBarry Smith Notes: 751*95fce210SBarry Smith 752*95fce210SBarry Smith In most cases, users should use PetscSFGatherBegin() and PetscSFScatterBegin() instead of manipulating multi 753*95fce210SBarry Smith directly. Since multi satisfies the stronger condition that each entry in the global space has exactly one incoming 754*95fce210SBarry Smith edge, it is a candidate for future optimization that might involve its removal. 755*95fce210SBarry Smith 756*95fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFGatherBegin(), PetscSFScatterBegin() 757*95fce210SBarry Smith @*/ 758*95fce210SBarry Smith PetscErrorCode PetscSFGetMultiSF(PetscSF sf,PetscSF *multi) 759*95fce210SBarry Smith { 760*95fce210SBarry Smith PetscErrorCode ierr; 761*95fce210SBarry Smith 762*95fce210SBarry Smith PetscFunctionBegin; 763*95fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 764*95fce210SBarry Smith PetscValidPointer(multi,2); 765*95fce210SBarry Smith if (sf->nroots < 0) { /* Graph has not been set yet; why do we need this? */ 766*95fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr); 767*95fce210SBarry Smith *multi = sf->multi; 768*95fce210SBarry Smith PetscFunctionReturn(0); 769*95fce210SBarry Smith } 770*95fce210SBarry Smith if (!sf->multi) { 771*95fce210SBarry Smith const PetscInt *indegree; 772*95fce210SBarry Smith PetscInt i,*inoffset,*outones,*outoffset; 773*95fce210SBarry Smith PetscSFNode *remote; 774*95fce210SBarry Smith ierr = PetscSFComputeDegreeBegin(sf,&indegree);CHKERRQ(ierr); 775*95fce210SBarry Smith ierr = PetscSFComputeDegreeEnd(sf,&indegree);CHKERRQ(ierr); 776*95fce210SBarry Smith ierr = PetscMalloc3(sf->nroots+1,PetscInt,&inoffset,sf->nleaves,PetscInt,&outones,sf->nleaves,PetscInt,&outoffset);CHKERRQ(ierr); 777*95fce210SBarry Smith inoffset[0] = 0; 778*95fce210SBarry Smith #if 1 779*95fce210SBarry Smith for (i=0; i<sf->nroots; i++) inoffset[i+1] = PetscMax(i+1, inoffset[i] + indegree[i]); 780*95fce210SBarry Smith #else 781*95fce210SBarry Smith for (i=0; i<sf->nroots; i++) inoffset[i+1] = inoffset[i] + indegree[i]; 782*95fce210SBarry Smith #endif 783*95fce210SBarry Smith for (i=0; i<sf->nleaves; i++) outones[i] = 1; 784*95fce210SBarry Smith ierr = PetscSFFetchAndOpBegin(sf,MPIU_INT,inoffset,outones,outoffset,MPIU_SUM);CHKERRQ(ierr); 785*95fce210SBarry Smith ierr = PetscSFFetchAndOpEnd(sf,MPIU_INT,inoffset,outones,outoffset,MPIU_SUM);CHKERRQ(ierr); 786*95fce210SBarry Smith for (i=0; i<sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */ 787*95fce210SBarry Smith #if 0 788*95fce210SBarry Smith #if defined(PETSC_USE_DEBUG) /* Check that the expected number of increments occurred */ 789*95fce210SBarry Smith for (i=0; i<sf->nroots; i++) { 790*95fce210SBarry Smith if (inoffset[i] + indegree[i] != inoffset[i+1]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect result after PetscSFFetchAndOp"); 791*95fce210SBarry Smith } 792*95fce210SBarry Smith #endif 793*95fce210SBarry Smith #endif 794*95fce210SBarry Smith ierr = PetscMalloc(sf->nleaves*sizeof(*remote),&remote);CHKERRQ(ierr); 795*95fce210SBarry Smith for (i=0; i<sf->nleaves; i++) { 796*95fce210SBarry Smith remote[i].rank = sf->remote[i].rank; 797*95fce210SBarry Smith remote[i].index = outoffset[i]; 798*95fce210SBarry Smith } 799*95fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr); 800*95fce210SBarry Smith ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr); 801*95fce210SBarry Smith if (sf->rankorder) { /* Sort the ranks */ 802*95fce210SBarry Smith PetscMPIInt rank; 803*95fce210SBarry Smith PetscInt *inranks,*newoffset,*outranks,*newoutoffset,*tmpoffset,maxdegree; 804*95fce210SBarry Smith PetscSFNode *newremote; 805*95fce210SBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); 806*95fce210SBarry Smith for (i=0,maxdegree=0; i<sf->nroots; i++) maxdegree = PetscMax(maxdegree,indegree[i]); 807*95fce210SBarry Smith ierr = PetscMalloc5(sf->multi->nroots,PetscInt,&inranks,sf->multi->nroots,PetscInt,&newoffset,sf->nleaves,PetscInt,&outranks,sf->nleaves,PetscInt,&newoutoffset,maxdegree,PetscInt,&tmpoffset);CHKERRQ(ierr); 808*95fce210SBarry Smith for (i=0; i<sf->nleaves; i++) outranks[i] = rank; 809*95fce210SBarry Smith ierr = PetscSFReduceBegin(sf->multi,MPIU_INT,outranks,inranks,MPI_REPLACE);CHKERRQ(ierr); 810*95fce210SBarry Smith ierr = PetscSFReduceEnd(sf->multi,MPIU_INT,outranks,inranks,MPI_REPLACE);CHKERRQ(ierr); 811*95fce210SBarry Smith /* Sort the incoming ranks at each vertex, build the inverse map */ 812*95fce210SBarry Smith for (i=0; i<sf->nroots; i++) { 813*95fce210SBarry Smith PetscInt j; 814*95fce210SBarry Smith for (j=0; j<indegree[i]; j++) tmpoffset[j] = j; 815*95fce210SBarry Smith ierr = PetscSortIntWithArray(indegree[i],inranks+inoffset[i],tmpoffset);CHKERRQ(ierr); 816*95fce210SBarry Smith for (j=0; j<indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j; 817*95fce210SBarry Smith } 818*95fce210SBarry Smith ierr = PetscSFBcastBegin(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr); 819*95fce210SBarry Smith ierr = PetscSFBcastEnd(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr); 820*95fce210SBarry Smith ierr = PetscMalloc(sf->nleaves*sizeof(*newremote),&newremote);CHKERRQ(ierr); 821*95fce210SBarry Smith for (i=0; i<sf->nleaves; i++) { 822*95fce210SBarry Smith newremote[i].rank = sf->remote[i].rank; 823*95fce210SBarry Smith newremote[i].index = newoutoffset[i]; 824*95fce210SBarry Smith } 825*95fce210SBarry Smith ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,NULL,PETSC_COPY_VALUES,newremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 826*95fce210SBarry Smith ierr = PetscFree5(inranks,newoffset,outranks,newoutoffset,tmpoffset);CHKERRQ(ierr); 827*95fce210SBarry Smith } 828*95fce210SBarry Smith ierr = PetscFree3(inoffset,outones,outoffset);CHKERRQ(ierr); 829*95fce210SBarry Smith } 830*95fce210SBarry Smith *multi = sf->multi; 831*95fce210SBarry Smith PetscFunctionReturn(0); 832*95fce210SBarry Smith } 833*95fce210SBarry Smith 834*95fce210SBarry Smith #undef __FUNCT__ 835*95fce210SBarry Smith #define __FUNCT__ "PetscSFCreateEmbeddedSF" 836*95fce210SBarry Smith /*@C 837*95fce210SBarry Smith PetscSFCreateEmbeddedSF - removes edges from all but the selected roots, does not remap indices 838*95fce210SBarry Smith 839*95fce210SBarry Smith Collective 840*95fce210SBarry Smith 841*95fce210SBarry Smith Input Arguments: 842*95fce210SBarry Smith + sf - original star forest 843*95fce210SBarry Smith . nroots - number of roots to select on this process 844*95fce210SBarry Smith - selected - selected roots on this process 845*95fce210SBarry Smith 846*95fce210SBarry Smith Output Arguments: 847*95fce210SBarry Smith . newsf - new star forest 848*95fce210SBarry Smith 849*95fce210SBarry Smith Level: advanced 850*95fce210SBarry Smith 851*95fce210SBarry Smith Note: 852*95fce210SBarry Smith To use the new PetscSF, it may be necessary to know the indices of the leaves that are still participating. This can 853*95fce210SBarry Smith be done by calling PetscSFGetGraph(). 854*95fce210SBarry Smith 855*95fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFGetGraph() 856*95fce210SBarry Smith @*/ 857*95fce210SBarry Smith PetscErrorCode PetscSFCreateEmbeddedSF(PetscSF sf,PetscInt nroots,const PetscInt *selected,PetscSF *newsf) 858*95fce210SBarry Smith { 859*95fce210SBarry Smith PetscErrorCode ierr; 860*95fce210SBarry Smith PetscInt i,nleaves,*ilocal,*rootdata,*leafdata; 861*95fce210SBarry Smith PetscSFNode *iremote; 862*95fce210SBarry Smith 863*95fce210SBarry Smith PetscFunctionBegin; 864*95fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 865*95fce210SBarry Smith if (nroots) PetscValidPointer(selected,3); 866*95fce210SBarry Smith PetscValidPointer(newsf,4); 867*95fce210SBarry Smith ierr = PetscMalloc2(sf->nroots,PetscInt,&rootdata,sf->nleaves,PetscInt,&leafdata);CHKERRQ(ierr); 868*95fce210SBarry Smith ierr = PetscMemzero(rootdata,sf->nroots*sizeof(PetscInt));CHKERRQ(ierr); 869*95fce210SBarry Smith ierr = PetscMemzero(leafdata,sf->nleaves*sizeof(PetscInt));CHKERRQ(ierr); 870*95fce210SBarry Smith for (i=0; i<nroots; i++) rootdata[selected[i]] = 1; 871*95fce210SBarry Smith ierr = PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr); 872*95fce210SBarry Smith ierr = PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr); 873*95fce210SBarry Smith 874*95fce210SBarry Smith for (i=0,nleaves=0; i<sf->nleaves; i++) nleaves += leafdata[i]; 875*95fce210SBarry Smith ierr = PetscMalloc(nleaves*sizeof(PetscInt),&ilocal);CHKERRQ(ierr); 876*95fce210SBarry Smith ierr = PetscMalloc(nleaves*sizeof(PetscSFNode),&iremote);CHKERRQ(ierr); 877*95fce210SBarry Smith for (i=0,nleaves=0; i<sf->nleaves; i++) { 878*95fce210SBarry Smith if (leafdata[i]) { 879*95fce210SBarry Smith ilocal[nleaves] = sf->mine ? sf->mine[i] : i; 880*95fce210SBarry Smith iremote[nleaves].rank = sf->remote[i].rank; 881*95fce210SBarry Smith iremote[nleaves].index = sf->remote[i].index; 882*95fce210SBarry Smith nleaves++; 883*95fce210SBarry Smith } 884*95fce210SBarry Smith } 885*95fce210SBarry Smith ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,newsf);CHKERRQ(ierr); 886*95fce210SBarry Smith ierr = PetscSFSetGraph(*newsf,sf->nroots,nleaves,ilocal,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 887*95fce210SBarry Smith ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr); 888*95fce210SBarry Smith PetscFunctionReturn(0); 889*95fce210SBarry Smith } 890*95fce210SBarry Smith 891*95fce210SBarry Smith #undef __FUNCT__ 892*95fce210SBarry Smith #define __FUNCT__ "PetscSFBcastBegin" 893*95fce210SBarry Smith /*@C 894*95fce210SBarry Smith PetscSFBcastBegin - begin pointwise broadcast to be concluded with call to PetscSFBcastEnd() 895*95fce210SBarry Smith 896*95fce210SBarry Smith Collective on PetscSF 897*95fce210SBarry Smith 898*95fce210SBarry Smith Input Arguments: 899*95fce210SBarry Smith + sf - star forest on which to communicate 900*95fce210SBarry Smith . unit - data type associated with each node 901*95fce210SBarry Smith - rootdata - buffer to broadcast 902*95fce210SBarry Smith 903*95fce210SBarry Smith Output Arguments: 904*95fce210SBarry Smith . leafdata - buffer to update with values from each leaf's respective root 905*95fce210SBarry Smith 906*95fce210SBarry Smith Level: intermediate 907*95fce210SBarry Smith 908*95fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFView(), PetscSFBcastEnd(), PetscSFReduceBegin() 909*95fce210SBarry Smith @*/ 910*95fce210SBarry Smith PetscErrorCode PetscSFBcastBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata) 911*95fce210SBarry Smith { 912*95fce210SBarry Smith PetscErrorCode ierr; 913*95fce210SBarry Smith 914*95fce210SBarry Smith PetscFunctionBegin; 915*95fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 916*95fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 917*95fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 918*95fce210SBarry Smith ierr = (*sf->ops->BcastBegin)(sf,unit,rootdata,leafdata);CHKERRQ(ierr); 919*95fce210SBarry Smith PetscFunctionReturn(0); 920*95fce210SBarry Smith } 921*95fce210SBarry Smith 922*95fce210SBarry Smith #undef __FUNCT__ 923*95fce210SBarry Smith #define __FUNCT__ "PetscSFBcastEnd" 924*95fce210SBarry Smith /*@C 925*95fce210SBarry Smith PetscSFBcastEnd - end a broadcast operation started with PetscSFBcastBegin() 926*95fce210SBarry Smith 927*95fce210SBarry Smith Collective 928*95fce210SBarry Smith 929*95fce210SBarry Smith Input Arguments: 930*95fce210SBarry Smith + sf - star forest 931*95fce210SBarry Smith . unit - data type 932*95fce210SBarry Smith - rootdata - buffer to broadcast 933*95fce210SBarry Smith 934*95fce210SBarry Smith Output Arguments: 935*95fce210SBarry Smith . leafdata - buffer to update with values from each leaf's respective root 936*95fce210SBarry Smith 937*95fce210SBarry Smith Level: intermediate 938*95fce210SBarry Smith 939*95fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFReduceEnd() 940*95fce210SBarry Smith @*/ 941*95fce210SBarry Smith PetscErrorCode PetscSFBcastEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata) 942*95fce210SBarry Smith { 943*95fce210SBarry Smith PetscErrorCode ierr; 944*95fce210SBarry Smith 945*95fce210SBarry Smith PetscFunctionBegin; 946*95fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 947*95fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 948*95fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 949*95fce210SBarry Smith ierr = (*sf->ops->BcastEnd)(sf,unit,rootdata,leafdata);CHKERRQ(ierr); 950*95fce210SBarry Smith PetscFunctionReturn(0); 951*95fce210SBarry Smith } 952*95fce210SBarry Smith 953*95fce210SBarry Smith #undef __FUNCT__ 954*95fce210SBarry Smith #define __FUNCT__ "PetscSFReduceBegin" 955*95fce210SBarry Smith /*@C 956*95fce210SBarry Smith PetscSFReduceBegin - begin reduction of leafdata into rootdata, to be completed with call to PetscSFReduceEnd() 957*95fce210SBarry Smith 958*95fce210SBarry Smith Collective 959*95fce210SBarry Smith 960*95fce210SBarry Smith Input Arguments: 961*95fce210SBarry Smith + sf - star forest 962*95fce210SBarry Smith . unit - data type 963*95fce210SBarry Smith . leafdata - values to reduce 964*95fce210SBarry Smith - op - reduction operation 965*95fce210SBarry Smith 966*95fce210SBarry Smith Output Arguments: 967*95fce210SBarry Smith . rootdata - result of reduction of values from all leaves of each root 968*95fce210SBarry Smith 969*95fce210SBarry Smith Level: intermediate 970*95fce210SBarry Smith 971*95fce210SBarry Smith .seealso: PetscSFBcastBegin() 972*95fce210SBarry Smith @*/ 973*95fce210SBarry Smith PetscErrorCode PetscSFReduceBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op) 974*95fce210SBarry Smith { 975*95fce210SBarry Smith PetscErrorCode ierr; 976*95fce210SBarry Smith 977*95fce210SBarry Smith PetscFunctionBegin; 978*95fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 979*95fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 980*95fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 981*95fce210SBarry Smith ierr = (sf->ops->ReduceBegin)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr); 982*95fce210SBarry Smith PetscFunctionReturn(0); 983*95fce210SBarry Smith } 984*95fce210SBarry Smith 985*95fce210SBarry Smith #undef __FUNCT__ 986*95fce210SBarry Smith #define __FUNCT__ "PetscSFReduceEnd" 987*95fce210SBarry Smith /*@C 988*95fce210SBarry Smith PetscSFReduceEnd - end a reduction operation started with PetscSFReduceBegin() 989*95fce210SBarry Smith 990*95fce210SBarry Smith Collective 991*95fce210SBarry Smith 992*95fce210SBarry Smith Input Arguments: 993*95fce210SBarry Smith + sf - star forest 994*95fce210SBarry Smith . unit - data type 995*95fce210SBarry Smith . leafdata - values to reduce 996*95fce210SBarry Smith - op - reduction operation 997*95fce210SBarry Smith 998*95fce210SBarry Smith Output Arguments: 999*95fce210SBarry Smith . rootdata - result of reduction of values from all leaves of each root 1000*95fce210SBarry Smith 1001*95fce210SBarry Smith Level: intermediate 1002*95fce210SBarry Smith 1003*95fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFBcastEnd() 1004*95fce210SBarry Smith @*/ 1005*95fce210SBarry Smith PetscErrorCode PetscSFReduceEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op) 1006*95fce210SBarry Smith { 1007*95fce210SBarry Smith PetscErrorCode ierr; 1008*95fce210SBarry Smith 1009*95fce210SBarry Smith PetscFunctionBegin; 1010*95fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 1011*95fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 1012*95fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 1013*95fce210SBarry Smith ierr = (*sf->ops->ReduceEnd)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr); 1014*95fce210SBarry Smith PetscFunctionReturn(0); 1015*95fce210SBarry Smith } 1016*95fce210SBarry Smith 1017*95fce210SBarry Smith #undef __FUNCT__ 1018*95fce210SBarry Smith #define __FUNCT__ "PetscSFComputeDegreeBegin" 1019*95fce210SBarry Smith /*@C 1020*95fce210SBarry Smith PetscSFComputeDegreeBegin - begin computation of degree for each root vertex, to be completed with PetscSFComputeDegreeEnd() 1021*95fce210SBarry Smith 1022*95fce210SBarry Smith Collective 1023*95fce210SBarry Smith 1024*95fce210SBarry Smith Input Arguments: 1025*95fce210SBarry Smith . sf - star forest 1026*95fce210SBarry Smith 1027*95fce210SBarry Smith Output Arguments: 1028*95fce210SBarry Smith . degree - degree of each root vertex 1029*95fce210SBarry Smith 1030*95fce210SBarry Smith Level: advanced 1031*95fce210SBarry Smith 1032*95fce210SBarry Smith .seealso: PetscSFGatherBegin() 1033*95fce210SBarry Smith @*/ 1034*95fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf,const PetscInt **degree) 1035*95fce210SBarry Smith { 1036*95fce210SBarry Smith PetscErrorCode ierr; 1037*95fce210SBarry Smith 1038*95fce210SBarry Smith PetscFunctionBegin; 1039*95fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 1040*95fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 1041*95fce210SBarry Smith PetscValidPointer(degree,2); 1042*95fce210SBarry Smith if (!sf->degree) { 1043*95fce210SBarry Smith PetscInt i; 1044*95fce210SBarry Smith ierr = PetscMalloc(sf->nroots*sizeof(PetscInt),&sf->degree);CHKERRQ(ierr); 1045*95fce210SBarry Smith ierr = PetscMalloc(sf->nleaves*sizeof(PetscInt),&sf->degreetmp);CHKERRQ(ierr); 1046*95fce210SBarry Smith for (i=0; i<sf->nroots; i++) sf->degree[i] = 0; 1047*95fce210SBarry Smith for (i=0; i<sf->nleaves; i++) sf->degreetmp[i] = 1; 1048*95fce210SBarry Smith ierr = PetscSFReduceBegin(sf,MPIU_INT,sf->degreetmp,sf->degree,MPIU_SUM);CHKERRQ(ierr); 1049*95fce210SBarry Smith } 1050*95fce210SBarry Smith *degree = NULL; 1051*95fce210SBarry Smith PetscFunctionReturn(0); 1052*95fce210SBarry Smith } 1053*95fce210SBarry Smith 1054*95fce210SBarry Smith #undef __FUNCT__ 1055*95fce210SBarry Smith #define __FUNCT__ "PetscSFComputeDegreeEnd" 1056*95fce210SBarry Smith /*@C 1057*95fce210SBarry Smith PetscSFComputeDegreeEnd - complete computation of degree for each root vertex, started with PetscSFComputeDegreeBegin() 1058*95fce210SBarry Smith 1059*95fce210SBarry Smith Collective 1060*95fce210SBarry Smith 1061*95fce210SBarry Smith Input Arguments: 1062*95fce210SBarry Smith . sf - star forest 1063*95fce210SBarry Smith 1064*95fce210SBarry Smith Output Arguments: 1065*95fce210SBarry Smith . degree - degree of each root vertex 1066*95fce210SBarry Smith 1067*95fce210SBarry Smith Level: developer 1068*95fce210SBarry Smith 1069*95fce210SBarry Smith .seealso: 1070*95fce210SBarry Smith @*/ 1071*95fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf,const PetscInt **degree) 1072*95fce210SBarry Smith { 1073*95fce210SBarry Smith PetscErrorCode ierr; 1074*95fce210SBarry Smith 1075*95fce210SBarry Smith PetscFunctionBegin; 1076*95fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 1077*95fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 1078*95fce210SBarry Smith if (!sf->degreeknown) { 1079*95fce210SBarry Smith ierr = PetscSFReduceEnd(sf,MPIU_INT,sf->degreetmp,sf->degree,MPIU_SUM);CHKERRQ(ierr); 1080*95fce210SBarry Smith ierr = PetscFree(sf->degreetmp);CHKERRQ(ierr); 1081*95fce210SBarry Smith 1082*95fce210SBarry Smith sf->degreeknown = PETSC_TRUE; 1083*95fce210SBarry Smith } 1084*95fce210SBarry Smith *degree = sf->degree; 1085*95fce210SBarry Smith PetscFunctionReturn(0); 1086*95fce210SBarry Smith } 1087*95fce210SBarry Smith 1088*95fce210SBarry Smith #undef __FUNCT__ 1089*95fce210SBarry Smith #define __FUNCT__ "PetscSFFetchAndOpBegin" 1090*95fce210SBarry Smith /*@C 1091*95fce210SBarry Smith PetscSFFetchAndOpBegin - begin operation that fetches values from root and updates atomically by applying operation using my leaf value, to be completed with PetscSFFetchAndOpEnd() 1092*95fce210SBarry Smith 1093*95fce210SBarry Smith Collective 1094*95fce210SBarry Smith 1095*95fce210SBarry Smith Input Arguments: 1096*95fce210SBarry Smith + sf - star forest 1097*95fce210SBarry Smith . unit - data type 1098*95fce210SBarry Smith . leafdata - leaf values to use in reduction 1099*95fce210SBarry Smith - op - operation to use for reduction 1100*95fce210SBarry Smith 1101*95fce210SBarry Smith Output Arguments: 1102*95fce210SBarry Smith + rootdata - root values to be updated, input state is seen by first process to perform an update 1103*95fce210SBarry Smith - leafupdate - state at each leaf's respective root immediately prior to my atomic update 1104*95fce210SBarry Smith 1105*95fce210SBarry Smith Level: advanced 1106*95fce210SBarry Smith 1107*95fce210SBarry Smith Note: 1108*95fce210SBarry Smith The update is only atomic at the granularity provided by the hardware. Different roots referenced by the same process 1109*95fce210SBarry Smith might be updated in a different order. Furthermore, if a composite type is used for the unit datatype, atomicity is 1110*95fce210SBarry Smith not guaranteed across the whole vertex. Therefore, this function is mostly only used with primitive types such as 1111*95fce210SBarry Smith integers. 1112*95fce210SBarry Smith 1113*95fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFReduceBegin(), PetscSFSetGraph() 1114*95fce210SBarry Smith @*/ 1115*95fce210SBarry Smith PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op) 1116*95fce210SBarry Smith { 1117*95fce210SBarry Smith PetscErrorCode ierr; 1118*95fce210SBarry Smith 1119*95fce210SBarry Smith PetscFunctionBegin; 1120*95fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 1121*95fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 1122*95fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 1123*95fce210SBarry Smith ierr = (*sf->ops->FetchAndOpBegin)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr); 1124*95fce210SBarry Smith PetscFunctionReturn(0); 1125*95fce210SBarry Smith } 1126*95fce210SBarry Smith 1127*95fce210SBarry Smith #undef __FUNCT__ 1128*95fce210SBarry Smith #define __FUNCT__ "PetscSFFetchAndOpEnd" 1129*95fce210SBarry Smith /*@C 1130*95fce210SBarry 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 1131*95fce210SBarry Smith 1132*95fce210SBarry Smith Collective 1133*95fce210SBarry Smith 1134*95fce210SBarry Smith Input Arguments: 1135*95fce210SBarry Smith + sf - star forest 1136*95fce210SBarry Smith . unit - data type 1137*95fce210SBarry Smith . leafdata - leaf values to use in reduction 1138*95fce210SBarry Smith - op - operation to use for reduction 1139*95fce210SBarry Smith 1140*95fce210SBarry Smith Output Arguments: 1141*95fce210SBarry Smith + rootdata - root values to be updated, input state is seen by first process to perform an update 1142*95fce210SBarry Smith - leafupdate - state at each leaf's respective root immediately prior to my atomic update 1143*95fce210SBarry Smith 1144*95fce210SBarry Smith Level: advanced 1145*95fce210SBarry Smith 1146*95fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFReduceEnd(), PetscSFSetGraph() 1147*95fce210SBarry Smith @*/ 1148*95fce210SBarry Smith PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op) 1149*95fce210SBarry Smith { 1150*95fce210SBarry Smith PetscErrorCode ierr; 1151*95fce210SBarry Smith 1152*95fce210SBarry Smith PetscFunctionBegin; 1153*95fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 1154*95fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 1155*95fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 1156*95fce210SBarry Smith ierr = (*sf->ops->FetchAndOpEnd)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr); 1157*95fce210SBarry Smith PetscFunctionReturn(0); 1158*95fce210SBarry Smith } 1159*95fce210SBarry Smith 1160*95fce210SBarry Smith #undef __FUNCT__ 1161*95fce210SBarry Smith #define __FUNCT__ "PetscSFGatherBegin" 1162*95fce210SBarry Smith /*@C 1163*95fce210SBarry Smith PetscSFGatherBegin - begin pointwise gather of all leaves into multi-roots, to be completed with PetscSFGatherEnd() 1164*95fce210SBarry Smith 1165*95fce210SBarry Smith Collective 1166*95fce210SBarry Smith 1167*95fce210SBarry Smith Input Arguments: 1168*95fce210SBarry Smith + sf - star forest 1169*95fce210SBarry Smith . unit - data type 1170*95fce210SBarry Smith - leafdata - leaf data to gather to roots 1171*95fce210SBarry Smith 1172*95fce210SBarry Smith Output Argument: 1173*95fce210SBarry Smith . multirootdata - root buffer to gather into, amount of space per root is equal to its degree 1174*95fce210SBarry Smith 1175*95fce210SBarry Smith Level: intermediate 1176*95fce210SBarry Smith 1177*95fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin() 1178*95fce210SBarry Smith @*/ 1179*95fce210SBarry Smith PetscErrorCode PetscSFGatherBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata) 1180*95fce210SBarry Smith { 1181*95fce210SBarry Smith PetscErrorCode ierr; 1182*95fce210SBarry Smith PetscSF multi; 1183*95fce210SBarry Smith 1184*95fce210SBarry Smith PetscFunctionBegin; 1185*95fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 1186*95fce210SBarry Smith ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr); 1187*95fce210SBarry Smith ierr = PetscSFReduceBegin(multi,unit,leafdata,multirootdata,MPI_REPLACE);CHKERRQ(ierr); 1188*95fce210SBarry Smith PetscFunctionReturn(0); 1189*95fce210SBarry Smith } 1190*95fce210SBarry Smith 1191*95fce210SBarry Smith #undef __FUNCT__ 1192*95fce210SBarry Smith #define __FUNCT__ "PetscSFGatherEnd" 1193*95fce210SBarry Smith /*@C 1194*95fce210SBarry Smith PetscSFGatherEnd - ends pointwise gather operation that was started with PetscSFGatherBegin() 1195*95fce210SBarry Smith 1196*95fce210SBarry Smith Collective 1197*95fce210SBarry Smith 1198*95fce210SBarry Smith Input Arguments: 1199*95fce210SBarry Smith + sf - star forest 1200*95fce210SBarry Smith . unit - data type 1201*95fce210SBarry Smith - leafdata - leaf data to gather to roots 1202*95fce210SBarry Smith 1203*95fce210SBarry Smith Output Argument: 1204*95fce210SBarry Smith . multirootdata - root buffer to gather into, amount of space per root is equal to its degree 1205*95fce210SBarry Smith 1206*95fce210SBarry Smith Level: intermediate 1207*95fce210SBarry Smith 1208*95fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd() 1209*95fce210SBarry Smith @*/ 1210*95fce210SBarry Smith PetscErrorCode PetscSFGatherEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata) 1211*95fce210SBarry Smith { 1212*95fce210SBarry Smith PetscErrorCode ierr; 1213*95fce210SBarry Smith PetscSF multi; 1214*95fce210SBarry Smith 1215*95fce210SBarry Smith PetscFunctionBegin; 1216*95fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 1217*95fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 1218*95fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 1219*95fce210SBarry Smith ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr); 1220*95fce210SBarry Smith ierr = PetscSFReduceEnd(multi,unit,leafdata,multirootdata,MPI_REPLACE);CHKERRQ(ierr); 1221*95fce210SBarry Smith PetscFunctionReturn(0); 1222*95fce210SBarry Smith } 1223*95fce210SBarry Smith 1224*95fce210SBarry Smith #undef __FUNCT__ 1225*95fce210SBarry Smith #define __FUNCT__ "PetscSFScatterBegin" 1226*95fce210SBarry Smith /*@C 1227*95fce210SBarry Smith PetscSFScatterBegin - begin pointwise scatter operation from multi-roots to leaves, to be completed with PetscSFScatterEnd() 1228*95fce210SBarry Smith 1229*95fce210SBarry Smith Collective 1230*95fce210SBarry Smith 1231*95fce210SBarry Smith Input Arguments: 1232*95fce210SBarry Smith + sf - star forest 1233*95fce210SBarry Smith . unit - data type 1234*95fce210SBarry Smith - multirootdata - root buffer to send to each leaf, one unit of data per leaf 1235*95fce210SBarry Smith 1236*95fce210SBarry Smith Output Argument: 1237*95fce210SBarry Smith . leafdata - leaf data to be update with personal data from each respective root 1238*95fce210SBarry Smith 1239*95fce210SBarry Smith Level: intermediate 1240*95fce210SBarry Smith 1241*95fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin() 1242*95fce210SBarry Smith @*/ 1243*95fce210SBarry Smith PetscErrorCode PetscSFScatterBegin(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata) 1244*95fce210SBarry Smith { 1245*95fce210SBarry Smith PetscErrorCode ierr; 1246*95fce210SBarry Smith PetscSF multi; 1247*95fce210SBarry Smith 1248*95fce210SBarry Smith PetscFunctionBegin; 1249*95fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 1250*95fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 1251*95fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 1252*95fce210SBarry Smith ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr); 1253*95fce210SBarry Smith ierr = PetscSFBcastBegin(multi,unit,multirootdata,leafdata);CHKERRQ(ierr); 1254*95fce210SBarry Smith PetscFunctionReturn(0); 1255*95fce210SBarry Smith } 1256*95fce210SBarry Smith 1257*95fce210SBarry Smith #undef __FUNCT__ 1258*95fce210SBarry Smith #define __FUNCT__ "PetscSFScatterEnd" 1259*95fce210SBarry Smith /*@C 1260*95fce210SBarry Smith PetscSFScatterEnd - ends pointwise scatter operation that was started with PetscSFScatterBegin() 1261*95fce210SBarry Smith 1262*95fce210SBarry Smith Collective 1263*95fce210SBarry Smith 1264*95fce210SBarry Smith Input Arguments: 1265*95fce210SBarry Smith + sf - star forest 1266*95fce210SBarry Smith . unit - data type 1267*95fce210SBarry Smith - multirootdata - root buffer to send to each leaf, one unit of data per leaf 1268*95fce210SBarry Smith 1269*95fce210SBarry Smith Output Argument: 1270*95fce210SBarry Smith . leafdata - leaf data to be update with personal data from each respective root 1271*95fce210SBarry Smith 1272*95fce210SBarry Smith Level: intermediate 1273*95fce210SBarry Smith 1274*95fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd() 1275*95fce210SBarry Smith @*/ 1276*95fce210SBarry Smith PetscErrorCode PetscSFScatterEnd(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata) 1277*95fce210SBarry Smith { 1278*95fce210SBarry Smith PetscErrorCode ierr; 1279*95fce210SBarry Smith PetscSF multi; 1280*95fce210SBarry Smith 1281*95fce210SBarry Smith PetscFunctionBegin; 1282*95fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 1283*95fce210SBarry Smith PetscSFCheckGraphSet(sf,1); 1284*95fce210SBarry Smith ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 1285*95fce210SBarry Smith ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr); 1286*95fce210SBarry Smith ierr = PetscSFBcastEnd(multi,unit,multirootdata,leafdata);CHKERRQ(ierr); 1287*95fce210SBarry Smith PetscFunctionReturn(0); 1288*95fce210SBarry Smith } 1289