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