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