1 #include <petsc-private/sfimpl.h> /*I "petscsf.h" I*/ 2 3 typedef struct _n_PetscSFDataLink *PetscSFDataLink; 4 typedef struct _n_PetscSFWinLink *PetscSFWinLink; 5 6 typedef struct { 7 PetscSFWindowSyncType sync; /* FENCE, LOCK, or ACTIVE synchronization */ 8 PetscSFDataLink link; /* List of MPI data types and windows, lazily constructed for each data type */ 9 PetscSFWinLink wins; /* List of active windows */ 10 } PetscSF_Window; 11 12 struct _n_PetscSFDataLink { 13 MPI_Datatype unit; 14 MPI_Datatype *mine; 15 MPI_Datatype *remote; 16 PetscSFDataLink next; 17 }; 18 19 struct _n_PetscSFWinLink { 20 PetscBool inuse; 21 size_t bytes; 22 void *addr; 23 MPI_Win win; 24 PetscBool epoch; 25 PetscSFWinLink next; 26 }; 27 28 const char *const PetscSFWindowSyncTypes[] = {"FENCE","LOCK","ACTIVE","PetscSFWindowSyncType","PETSCSF_WINDOW_SYNC_",0}; 29 30 #undef __FUNCT__ 31 #define __FUNCT__ "PetscSFWindowOpTranslate" 32 /* Built-in MPI_Ops act elementwise inside MPI_Accumulate, but cannot be used with composite types inside collectives (MPI_Allreduce) */ 33 static PetscErrorCode PetscSFWindowOpTranslate(MPI_Op *op) 34 { 35 36 PetscFunctionBegin; 37 if (*op == MPIU_SUM) *op = MPI_SUM; 38 else if (*op == MPIU_MAX) *op = MPI_MAX; 39 else if (*op == MPIU_MIN) *op = MPI_MIN; 40 PetscFunctionReturn(0); 41 } 42 43 #undef __FUNCT__ 44 #define __FUNCT__ "PetscSFWindowGetDataTypes" 45 /*@C 46 PetscSFWindowGetDataTypes - gets composite local and remote data types for each rank 47 48 Not Collective 49 50 Input Arguments: 51 + sf - star forest 52 - unit - data type for each node 53 54 Output Arguments: 55 + localtypes - types describing part of local leaf buffer referencing each remote rank 56 - remotetypes - types describing part of remote root buffer referenced for each remote rank 57 58 Level: developer 59 60 .seealso: PetscSFSetGraph(), PetscSFView() 61 @*/ 62 static PetscErrorCode PetscSFWindowGetDataTypes(PetscSF sf,MPI_Datatype unit,const MPI_Datatype **localtypes,const MPI_Datatype **remotetypes) 63 { 64 PetscSF_Window *w = (PetscSF_Window*)sf->data; 65 PetscErrorCode ierr; 66 PetscSFDataLink link; 67 PetscInt i,nranks; 68 const PetscInt *roffset,*rmine,*rremote; 69 const PetscMPIInt *ranks; 70 71 PetscFunctionBegin; 72 /* Look for types in cache */ 73 for (link=w->link; link; link=link->next) { 74 PetscBool match; 75 ierr = MPIPetsc_Type_compare(unit,link->unit,&match);CHKERRQ(ierr); 76 if (match) { 77 *localtypes = link->mine; 78 *remotetypes = link->remote; 79 PetscFunctionReturn(0); 80 } 81 } 82 83 /* Create new composite types for each send rank */ 84 ierr = PetscSFGetRanks(sf,&nranks,&ranks,&roffset,&rmine,&rremote);CHKERRQ(ierr); 85 ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr); 86 ierr = MPI_Type_dup(unit,&link->unit);CHKERRQ(ierr); 87 ierr = PetscMalloc2(nranks,&link->mine,nranks,&link->remote);CHKERRQ(ierr); 88 for (i=0; i<nranks; i++) { 89 PETSC_UNUSED PetscInt rcount = roffset[i+1] - roffset[i]; 90 PetscMPIInt *rmine,*rremote; 91 #if !defined(PETSC_USE_64BIT_INDICES) 92 rmine = sf->rmine + sf->roffset[i]; 93 rremote = sf->rremote + sf->roffset[i]; 94 #else 95 PetscInt j; 96 ierr = PetscMalloc2(rcount,&rmine,rcount,&rremote);CHKERRQ(ierr); 97 for (j=0; j<rcount; j++) { 98 ierr = PetscMPIIntCast(sf->rmine[sf->roffset[i]+j],rmine+j);CHKERRQ(ierr); 99 ierr = PetscMPIIntCast(sf->rremote[sf->roffset[i]+j],rremote+j);CHKERRQ(ierr); 100 } 101 #endif 102 ierr = MPI_Type_create_indexed_block(rcount,1,rmine,link->unit,&link->mine[i]);CHKERRQ(ierr); 103 ierr = MPI_Type_create_indexed_block(rcount,1,rremote,link->unit,&link->remote[i]);CHKERRQ(ierr); 104 #if defined(PETSC_USE_64BIT_INDICES) 105 ierr = PetscFree2(rmine,rremote);CHKERRQ(ierr); 106 #endif 107 ierr = MPI_Type_commit(&link->mine[i]);CHKERRQ(ierr); 108 ierr = MPI_Type_commit(&link->remote[i]);CHKERRQ(ierr); 109 } 110 link->next = w->link; 111 w->link = link; 112 113 *localtypes = link->mine; 114 *remotetypes = link->remote; 115 PetscFunctionReturn(0); 116 } 117 118 #undef __FUNCT__ 119 #define __FUNCT__ "PetscSFWindowSetSyncType" 120 /*@C 121 PetscSFWindowSetSyncType - set synchrozitaion type for PetscSF communication 122 123 Logically Collective 124 125 Input Arguments: 126 + sf - star forest for communication 127 - sync - synchronization type 128 129 Options Database Key: 130 . -sf_window_sync <sync> - sets the synchronization type FENCE, LOCK, or ACTIVE (see PetscSFWindowSyncType) 131 132 Level: advanced 133 134 .seealso: PetscSFSetFromOptions(), PetscSFWindowGetSyncType() 135 @*/ 136 PetscErrorCode PetscSFWindowSetSyncType(PetscSF sf,PetscSFWindowSyncType sync) 137 { 138 PetscErrorCode ierr; 139 140 PetscFunctionBegin; 141 PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 142 PetscValidLogicalCollectiveEnum(sf,sync,2); 143 ierr = PetscUseMethod(sf,"PetscSFWindowSetSyncType_C",(PetscSF,PetscSFWindowSyncType),(sf,sync));CHKERRQ(ierr); 144 PetscFunctionReturn(0); 145 } 146 147 #undef __FUNCT__ 148 #define __FUNCT__ "PetscSFWindowSetSyncType_Window" 149 static PetscErrorCode PetscSFWindowSetSyncType_Window(PetscSF sf,PetscSFWindowSyncType sync) 150 { 151 PetscSF_Window *w = (PetscSF_Window*)sf->data; 152 153 PetscFunctionBegin; 154 w->sync = sync; 155 PetscFunctionReturn(0); 156 } 157 158 #undef __FUNCT__ 159 #define __FUNCT__ "PetscSFWindowGetSyncType" 160 /*@C 161 PetscSFWindowGetSyncType - get synchrozitaion type for PetscSF communication 162 163 Logically Collective 164 165 Input Argument: 166 . sf - star forest for communication 167 168 Output Argument: 169 . sync - synchronization type 170 171 Level: advanced 172 173 .seealso: PetscSFGetFromOptions(), PetscSFWindowSetSyncType() 174 @*/ 175 PetscErrorCode PetscSFWindowGetSyncType(PetscSF sf,PetscSFWindowSyncType *sync) 176 { 177 PetscErrorCode ierr; 178 179 PetscFunctionBegin; 180 PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 181 PetscValidPointer(sync,2); 182 ierr = PetscTryMethod(sf,"PetscSFWindowGetSyncType_C",(PetscSF,PetscSFWindowSyncType*),(sf,sync));CHKERRQ(ierr); 183 PetscFunctionReturn(0); 184 } 185 186 #undef __FUNCT__ 187 #define __FUNCT__ "PetscSFWindowGetSyncType_Window" 188 static PetscErrorCode PetscSFWindowGetSyncType_Window(PetscSF sf,PetscSFWindowSyncType *sync) 189 { 190 PetscSF_Window *w = (PetscSF_Window*)sf->data; 191 192 PetscFunctionBegin; 193 *sync = w->sync; 194 PetscFunctionReturn(0); 195 } 196 197 #undef __FUNCT__ 198 #define __FUNCT__ "PetscSFGetWindow" 199 /*@C 200 PetscSFGetWindow - Get a window for use with a given data type 201 202 Collective on PetscSF 203 204 Input Arguments: 205 + sf - star forest 206 . unit - data type 207 . array - array to be sent 208 . epoch - PETSC_TRUE to acquire the window and start an epoch, PETSC_FALSE to just acquire the window 209 . fenceassert - assert parameter for call to MPI_Win_fence(), if PETSCSF_WINDOW_SYNC_FENCE 210 . postassert - assert parameter for call to MPI_Win_post(), if PETSCSF_WINDOW_SYNC_ACTIVE 211 - startassert - assert parameter for call to MPI_Win_start(), if PETSCSF_WINDOW_SYNC_ACTIVE 212 213 Output Arguments: 214 . win - window 215 216 Level: developer 217 218 Developer Notes: 219 This currently always creates a new window. This is more synchronous than necessary. An alternative is to try to 220 reuse an existing window created with the same array. Another alternative is to maintain a cache of windows and reuse 221 whichever one is available, by copying the array into it if necessary. 222 223 .seealso: PetscSFGetRanks(), PetscSFWindowGetDataTypes() 224 @*/ 225 static PetscErrorCode PetscSFGetWindow(PetscSF sf,MPI_Datatype unit,void *array,PetscBool epoch,PetscMPIInt fenceassert,PetscMPIInt postassert,PetscMPIInt startassert,MPI_Win *win) 226 { 227 PetscSF_Window *w = (PetscSF_Window*)sf->data; 228 PetscErrorCode ierr; 229 MPI_Aint lb,lb_true,bytes,bytes_true; 230 PetscSFWinLink link; 231 232 PetscFunctionBegin; 233 ierr = MPI_Type_get_extent(unit,&lb,&bytes);CHKERRQ(ierr); 234 ierr = MPI_Type_get_true_extent(unit,&lb_true,&bytes_true);CHKERRQ(ierr); 235 if (lb != 0 || lb_true != 0) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_SUP,"No support for unit type with nonzero lower bound, write petsc-maint@mcs.anl.gov if you want this feature"); 236 if (bytes != bytes_true) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_SUP,"No support for unit type with modified extent, write petsc-maint@mcs.anl.gov if you want this feature"); 237 ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr); 238 239 link->bytes = bytes; 240 link->addr = array; 241 242 ierr = MPI_Win_create(array,(MPI_Aint)bytes*sf->nroots,(PetscMPIInt)bytes,MPI_INFO_NULL,PetscObjectComm((PetscObject)sf),&link->win);CHKERRQ(ierr); 243 244 link->epoch = epoch; 245 link->next = w->wins; 246 link->inuse = PETSC_TRUE; 247 w->wins = link; 248 *win = link->win; 249 250 if (epoch) { 251 switch (w->sync) { 252 case PETSCSF_WINDOW_SYNC_FENCE: 253 ierr = MPI_Win_fence(fenceassert,*win);CHKERRQ(ierr); 254 break; 255 case PETSCSF_WINDOW_SYNC_LOCK: /* Handled outside */ 256 break; 257 case PETSCSF_WINDOW_SYNC_ACTIVE: { 258 MPI_Group ingroup,outgroup; 259 ierr = PetscSFGetGroups(sf,&ingroup,&outgroup);CHKERRQ(ierr); 260 ierr = MPI_Win_post(ingroup,postassert,*win);CHKERRQ(ierr); 261 ierr = MPI_Win_start(outgroup,startassert,*win);CHKERRQ(ierr); 262 } break; 263 default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_PLIB,"Unknown synchronization type"); 264 } 265 } 266 PetscFunctionReturn(0); 267 } 268 269 #undef __FUNCT__ 270 #define __FUNCT__ "PetscSFFindWindow" 271 /*@C 272 PetscSFFindWindow - Finds a window that is already in use 273 274 Not Collective 275 276 Input Arguments: 277 + sf - star forest 278 . unit - data type 279 - array - array with which the window is associated 280 281 Output Arguments: 282 . win - window 283 284 Level: developer 285 286 .seealso: PetscSFGetWindow(), PetscSFRestoreWindow() 287 @*/ 288 static PetscErrorCode PetscSFFindWindow(PetscSF sf,MPI_Datatype unit,const void *array,MPI_Win *win) 289 { 290 PetscSF_Window *w = (PetscSF_Window*)sf->data; 291 PetscSFWinLink link; 292 293 PetscFunctionBegin; 294 *win = MPI_WIN_NULL; 295 for (link=w->wins; link; link=link->next) { 296 if (array == link->addr) { 297 *win = link->win; 298 PetscFunctionReturn(0); 299 } 300 } 301 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Requested window not in use"); 302 PetscFunctionReturn(0); 303 } 304 305 #undef __FUNCT__ 306 #define __FUNCT__ "PetscSFRestoreWindow" 307 /*@C 308 PetscSFRestoreWindow - Restores a window obtained with PetscSFGetWindow() 309 310 Collective 311 312 Input Arguments: 313 + sf - star forest 314 . unit - data type 315 . array - array associated with window 316 . epoch - close an epoch, must match argument to PetscSFGetWindow() 317 - win - window 318 319 Level: developer 320 321 .seealso: PetscSFFindWindow() 322 @*/ 323 static PetscErrorCode PetscSFRestoreWindow(PetscSF sf,MPI_Datatype unit,const void *array,PetscBool epoch,PetscMPIInt fenceassert,MPI_Win *win) 324 { 325 PetscSF_Window *w = (PetscSF_Window*)sf->data; 326 PetscErrorCode ierr; 327 PetscSFWinLink *p,link; 328 329 PetscFunctionBegin; 330 for (p=&w->wins; *p; p=&(*p)->next) { 331 link = *p; 332 if (*win == link->win) { 333 if (array != link->addr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Matched window, but not array"); 334 if (epoch != link->epoch) { 335 if (epoch) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"No epoch to end"); 336 else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Restoring window without ending epoch"); 337 } 338 *p = link->next; 339 goto found; 340 } 341 } 342 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Requested window not in use"); 343 344 found: 345 if (epoch) { 346 switch (w->sync) { 347 case PETSCSF_WINDOW_SYNC_FENCE: 348 ierr = MPI_Win_fence(fenceassert,*win);CHKERRQ(ierr); 349 break; 350 case PETSCSF_WINDOW_SYNC_LOCK: 351 break; /* handled outside */ 352 case PETSCSF_WINDOW_SYNC_ACTIVE: { 353 ierr = MPI_Win_complete(*win);CHKERRQ(ierr); 354 ierr = MPI_Win_wait(*win);CHKERRQ(ierr); 355 } break; 356 default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_PLIB,"Unknown synchronization type"); 357 } 358 } 359 360 ierr = MPI_Win_free(&link->win);CHKERRQ(ierr); 361 ierr = PetscFree(link);CHKERRQ(ierr); 362 *win = MPI_WIN_NULL; 363 PetscFunctionReturn(0); 364 } 365 366 #undef __FUNCT__ 367 #define __FUNCT__ "PetscSFSetUp_Window" 368 static PetscErrorCode PetscSFSetUp_Window(PetscSF sf) 369 { 370 PetscSF_Window *w = (PetscSF_Window*)sf->data; 371 PetscErrorCode ierr; 372 MPI_Group ingroup,outgroup; 373 374 PetscFunctionBegin; 375 switch (w->sync) { 376 case PETSCSF_WINDOW_SYNC_ACTIVE: 377 ierr = PetscSFGetGroups(sf,&ingroup,&outgroup);CHKERRQ(ierr); 378 default: 379 break; 380 } 381 PetscFunctionReturn(0); 382 } 383 384 #undef __FUNCT__ 385 #define __FUNCT__ "PetscSFSetFromOptions_Window" 386 static PetscErrorCode PetscSFSetFromOptions_Window(PetscSF sf) 387 { 388 PetscSF_Window *w = (PetscSF_Window*)sf->data; 389 PetscErrorCode ierr; 390 391 PetscFunctionBegin; 392 ierr = PetscOptionsHead("PetscSF Window options");CHKERRQ(ierr); 393 ierr = PetscOptionsEnum("-sf_window_sync","synchronization type to use for PetscSF Window communication","PetscSFWindowSetSyncType",PetscSFWindowSyncTypes,(PetscEnum)w->sync,(PetscEnum*)&w->sync,NULL);CHKERRQ(ierr); 394 ierr = PetscOptionsTail();CHKERRQ(ierr); 395 PetscFunctionReturn(0); 396 } 397 398 #undef __FUNCT__ 399 #define __FUNCT__ "PetscSFReset_Window" 400 static PetscErrorCode PetscSFReset_Window(PetscSF sf) 401 { 402 PetscSF_Window *w = (PetscSF_Window*)sf->data; 403 PetscErrorCode ierr; 404 PetscSFDataLink link,next; 405 PetscSFWinLink wlink,wnext; 406 PetscInt i; 407 408 PetscFunctionBegin; 409 for (link=w->link; link; link=next) { 410 next = link->next; 411 ierr = MPI_Type_free(&link->unit);CHKERRQ(ierr); 412 for (i=0; i<sf->nranks; i++) { 413 ierr = MPI_Type_free(&link->mine[i]);CHKERRQ(ierr); 414 ierr = MPI_Type_free(&link->remote[i]);CHKERRQ(ierr); 415 } 416 ierr = PetscFree2(link->mine,link->remote);CHKERRQ(ierr); 417 ierr = PetscFree(link);CHKERRQ(ierr); 418 } 419 w->link = NULL; 420 for (wlink=w->wins; wlink; wlink=wnext) { 421 wnext = wlink->next; 422 if (wlink->inuse) SETERRQ1(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Window still in use with address %p",(void*)wlink->addr); 423 ierr = MPI_Win_free(&wlink->win);CHKERRQ(ierr); 424 ierr = PetscFree(wlink);CHKERRQ(ierr); 425 } 426 w->wins = NULL; 427 PetscFunctionReturn(0); 428 } 429 430 #undef __FUNCT__ 431 #define __FUNCT__ "PetscSFDestroy_Window" 432 static PetscErrorCode PetscSFDestroy_Window(PetscSF sf) 433 { 434 PetscErrorCode ierr; 435 436 PetscFunctionBegin; 437 ierr = PetscSFReset_Window(sf);CHKERRQ(ierr); 438 ierr = PetscFree(sf->data);CHKERRQ(ierr); 439 ierr = PetscObjectComposeFunction((PetscObject)sf,"PetscSFWindowSetSyncType_C",NULL);CHKERRQ(ierr); 440 ierr = PetscObjectComposeFunction((PetscObject)sf,"PetscSFWindowGetSyncType_C",NULL);CHKERRQ(ierr); 441 PetscFunctionReturn(0); 442 } 443 444 #undef __FUNCT__ 445 #define __FUNCT__ "PetscSFView_Window" 446 static PetscErrorCode PetscSFView_Window(PetscSF sf,PetscViewer viewer) 447 { 448 PetscSF_Window *w = (PetscSF_Window*)sf->data; 449 PetscErrorCode ierr; 450 PetscBool iascii; 451 452 PetscFunctionBegin; 453 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 454 if (iascii) { 455 ierr = PetscViewerASCIIPrintf(viewer," synchronization=%s sort=%s\n",PetscSFWindowSyncTypes[w->sync],sf->rankorder ? "rank-order" : "unordered");CHKERRQ(ierr); 456 } 457 PetscFunctionReturn(0); 458 } 459 460 #undef __FUNCT__ 461 #define __FUNCT__ "PetscSFDuplicate_Window" 462 static PetscErrorCode PetscSFDuplicate_Window(PetscSF sf,PetscSFDuplicateOption opt,PetscSF newsf) 463 { 464 PetscSF_Window *w = (PetscSF_Window*)sf->data; 465 PetscErrorCode ierr; 466 PetscSFWindowSyncType synctype; 467 468 PetscFunctionBegin; 469 synctype = w->sync; 470 if (!sf->setupcalled) { 471 /* HACK: Must use FENCE or LOCK when called from PetscSFGetGroups() because ACTIVE here would cause recursion. */ 472 synctype = PETSCSF_WINDOW_SYNC_LOCK; 473 } 474 ierr = PetscSFWindowSetSyncType(newsf,synctype);CHKERRQ(ierr); 475 PetscFunctionReturn(0); 476 } 477 478 #undef __FUNCT__ 479 #define __FUNCT__ "PetscSFBcastBegin_Window" 480 static PetscErrorCode PetscSFBcastBegin_Window(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata) 481 { 482 PetscSF_Window *w = (PetscSF_Window*)sf->data; 483 PetscErrorCode ierr; 484 PetscInt i,nranks; 485 const PetscMPIInt *ranks; 486 const MPI_Datatype *mine,*remote; 487 MPI_Win win; 488 489 PetscFunctionBegin; 490 ierr = PetscSFGetRanks(sf,&nranks,&ranks,NULL,NULL,NULL);CHKERRQ(ierr); 491 ierr = PetscSFWindowGetDataTypes(sf,unit,&mine,&remote);CHKERRQ(ierr); 492 ierr = PetscSFGetWindow(sf,unit,(void*)rootdata,PETSC_TRUE,MPI_MODE_NOPUT|MPI_MODE_NOPRECEDE,MPI_MODE_NOPUT,0,&win);CHKERRQ(ierr); 493 for (i=0; i<nranks; i++) { 494 if (w->sync == PETSCSF_WINDOW_SYNC_LOCK) {ierr = MPI_Win_lock(MPI_LOCK_SHARED,ranks[i],MPI_MODE_NOCHECK,win);CHKERRQ(ierr);} 495 ierr = MPI_Get(leafdata,1,mine[i],ranks[i],0,1,remote[i],win);CHKERRQ(ierr); 496 if (w->sync == PETSCSF_WINDOW_SYNC_LOCK) {ierr = MPI_Win_unlock(ranks[i],win);CHKERRQ(ierr);} 497 } 498 PetscFunctionReturn(0); 499 } 500 501 #undef __FUNCT__ 502 #define __FUNCT__ "PetscSFBcastEnd_Window" 503 PetscErrorCode PetscSFBcastEnd_Window(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata) 504 { 505 PetscErrorCode ierr; 506 MPI_Win win; 507 508 PetscFunctionBegin; 509 ierr = PetscSFFindWindow(sf,unit,rootdata,&win);CHKERRQ(ierr); 510 ierr = PetscSFRestoreWindow(sf,unit,rootdata,PETSC_TRUE,MPI_MODE_NOSTORE|MPI_MODE_NOSUCCEED,&win);CHKERRQ(ierr); 511 PetscFunctionReturn(0); 512 } 513 514 #undef __FUNCT__ 515 #define __FUNCT__ "PetscSFReduceBegin_Window" 516 PetscErrorCode PetscSFReduceBegin_Window(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op) 517 { 518 PetscSF_Window *w = (PetscSF_Window*)sf->data; 519 PetscErrorCode ierr; 520 PetscInt i,nranks; 521 const PetscMPIInt *ranks; 522 const MPI_Datatype *mine,*remote; 523 MPI_Win win; 524 525 PetscFunctionBegin; 526 ierr = PetscSFGetRanks(sf,&nranks,&ranks,NULL,NULL,NULL);CHKERRQ(ierr); 527 ierr = PetscSFWindowGetDataTypes(sf,unit,&mine,&remote);CHKERRQ(ierr); 528 ierr = PetscSFWindowOpTranslate(&op);CHKERRQ(ierr); 529 ierr = PetscSFGetWindow(sf,unit,rootdata,PETSC_TRUE,MPI_MODE_NOPRECEDE,0,0,&win);CHKERRQ(ierr); 530 for (i=0; i<nranks; i++) { 531 if (w->sync == PETSCSF_WINDOW_SYNC_LOCK) {ierr = MPI_Win_lock(MPI_LOCK_SHARED,ranks[i],MPI_MODE_NOCHECK,win);CHKERRQ(ierr);} 532 ierr = MPI_Accumulate((void*)leafdata,1,mine[i],ranks[i],0,1,remote[i],op,win);CHKERRQ(ierr); 533 if (w->sync == PETSCSF_WINDOW_SYNC_LOCK) {ierr = MPI_Win_unlock(ranks[i],win);CHKERRQ(ierr);} 534 } 535 PetscFunctionReturn(0); 536 } 537 538 #undef __FUNCT__ 539 #define __FUNCT__ "PetscSFReduceEnd_Window" 540 static PetscErrorCode PetscSFReduceEnd_Window(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op) 541 { 542 PetscSF_Window *w = (PetscSF_Window*)sf->data; 543 PetscErrorCode ierr; 544 MPI_Win win; 545 546 PetscFunctionBegin; 547 if (!w->wins) PetscFunctionReturn(0); 548 ierr = PetscSFFindWindow(sf,unit,rootdata,&win);CHKERRQ(ierr); 549 ierr = MPI_Win_fence(MPI_MODE_NOSUCCEED,win);CHKERRQ(ierr); 550 ierr = PetscSFRestoreWindow(sf,unit,rootdata,PETSC_TRUE,MPI_MODE_NOSUCCEED,&win);CHKERRQ(ierr); 551 PetscFunctionReturn(0); 552 } 553 #undef __FUNCT__ 554 #define __FUNCT__ "PetscSFFetchAndOpBegin_Window" 555 static PetscErrorCode PetscSFFetchAndOpBegin_Window(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op) 556 { 557 PetscErrorCode ierr; 558 PetscInt i,nranks; 559 const PetscMPIInt *ranks; 560 const MPI_Datatype *mine,*remote; 561 MPI_Win win; 562 563 PetscFunctionBegin; 564 ierr = PetscSFGetRanks(sf,&nranks,&ranks,NULL,NULL,NULL);CHKERRQ(ierr); 565 ierr = PetscSFWindowGetDataTypes(sf,unit,&mine,&remote);CHKERRQ(ierr); 566 ierr = PetscSFWindowOpTranslate(&op);CHKERRQ(ierr); 567 ierr = PetscSFGetWindow(sf,unit,rootdata,PETSC_FALSE,0,0,0,&win);CHKERRQ(ierr); 568 for (i=0; i<sf->nranks; i++) { 569 ierr = MPI_Win_lock(MPI_LOCK_EXCLUSIVE,sf->ranks[i],0,win);CHKERRQ(ierr); 570 ierr = MPI_Get(leafupdate,1,mine[i],ranks[i],0,1,remote[i],win);CHKERRQ(ierr); 571 ierr = MPI_Accumulate((void*)leafdata,1,mine[i],ranks[i],0,1,remote[i],op,win);CHKERRQ(ierr); 572 ierr = MPI_Win_unlock(ranks[i],win);CHKERRQ(ierr); 573 } 574 PetscFunctionReturn(0); 575 } 576 577 #undef __FUNCT__ 578 #define __FUNCT__ "PetscSFFetchAndOpEnd_Window" 579 static PetscErrorCode PetscSFFetchAndOpEnd_Window(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op) 580 { 581 PetscErrorCode ierr; 582 MPI_Win win; 583 584 PetscFunctionBegin; 585 ierr = PetscSFFindWindow(sf,unit,rootdata,&win);CHKERRQ(ierr); 586 /* Nothing to do currently because MPI_LOCK_EXCLUSIVE is used in PetscSFFetchAndOpBegin(), rendering this implementation synchronous. */ 587 ierr = PetscSFRestoreWindow(sf,unit,rootdata,PETSC_FALSE,0,&win);CHKERRQ(ierr); 588 PetscFunctionReturn(0); 589 } 590 591 #undef __FUNCT__ 592 #define __FUNCT__ "PetscSFCreate_Window" 593 PETSC_EXTERN PetscErrorCode PetscSFCreate_Window(PetscSF sf) 594 { 595 PetscSF_Window *w = (PetscSF_Window*)sf->data; 596 PetscErrorCode ierr; 597 598 PetscFunctionBegin; 599 sf->ops->SetUp = PetscSFSetUp_Window; 600 sf->ops->SetFromOptions = PetscSFSetFromOptions_Window; 601 sf->ops->Reset = PetscSFReset_Window; 602 sf->ops->Destroy = PetscSFDestroy_Window; 603 sf->ops->View = PetscSFView_Window; 604 sf->ops->Duplicate = PetscSFDuplicate_Window; 605 sf->ops->BcastBegin = PetscSFBcastBegin_Window; 606 sf->ops->BcastEnd = PetscSFBcastEnd_Window; 607 sf->ops->ReduceBegin = PetscSFReduceBegin_Window; 608 sf->ops->ReduceEnd = PetscSFReduceEnd_Window; 609 sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Window; 610 sf->ops->FetchAndOpEnd = PetscSFFetchAndOpEnd_Window; 611 612 ierr = PetscNewLog(sf,&w);CHKERRQ(ierr); 613 sf->data = (void*)w; 614 w->sync = PETSCSF_WINDOW_SYNC_FENCE; 615 616 ierr = PetscObjectComposeFunction((PetscObject)sf,"PetscSFWindowSetSyncType_C",PetscSFWindowSetSyncType_Window);CHKERRQ(ierr); 617 ierr = PetscObjectComposeFunction((PetscObject)sf,"PetscSFWindowGetSyncType_C",PetscSFWindowGetSyncType_Window);CHKERRQ(ierr); 618 619 #if defined(OMPI_MAJOR_VERSION) && (OMPI_MAJOR_VERSION < 1 || (OMPI_MAJOR_VERSION == 1 && OMPI_MINOR_VERSION <= 6)) 620 { 621 PetscBool ackbug = PETSC_FALSE; 622 ierr = PetscOptionsGetBool(NULL,"-acknowledge_ompi_onesided_bug",&ackbug,NULL);CHKERRQ(ierr); 623 if (ackbug) { 624 ierr = PetscInfo(sf,"Acknowledged Open MPI bug, proceeding anyway. Expect memory corruption.\n");CHKERRQ(ierr); 625 } else SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_LIB,"Open MPI is known to be buggy (https://svn.open-mpi.org/trac/ompi/ticket/1905 and 2656), use -acknowledge_ompi_onesided_bug to proceed"); 626 } 627 #endif 628 PetscFunctionReturn(0); 629 } 630