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