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