195fce210SBarry Smith #include <petsc-private/sfimpl.h> /*I "petscsf.h" I*/ 295fce210SBarry Smith 395fce210SBarry Smith typedef struct _n_PetscSFDataLink *PetscSFDataLink; 495fce210SBarry Smith typedef struct _n_PetscSFWinLink *PetscSFWinLink; 595fce210SBarry Smith 695fce210SBarry Smith typedef struct { 795fce210SBarry Smith PetscSFWindowSyncType sync; /* FENCE, LOCK, or ACTIVE synchronization */ 895fce210SBarry Smith PetscSFDataLink link; /* List of MPI data types and windows, lazily constructed for each data type */ 995fce210SBarry Smith PetscSFWinLink wins; /* List of active windows */ 1095fce210SBarry Smith } PetscSF_Window; 1195fce210SBarry Smith 1295fce210SBarry Smith struct _n_PetscSFDataLink { 1395fce210SBarry Smith MPI_Datatype unit; 1495fce210SBarry Smith MPI_Datatype *mine; 1595fce210SBarry Smith MPI_Datatype *remote; 1695fce210SBarry Smith PetscSFDataLink next; 1795fce210SBarry Smith }; 1895fce210SBarry Smith 1995fce210SBarry Smith struct _n_PetscSFWinLink { 2095fce210SBarry Smith PetscBool inuse; 2195fce210SBarry Smith size_t bytes; 2295fce210SBarry Smith void *addr; 2395fce210SBarry Smith MPI_Win win; 2495fce210SBarry Smith PetscBool epoch; 2595fce210SBarry Smith PetscSFWinLink next; 2695fce210SBarry Smith }; 2795fce210SBarry Smith 2895fce210SBarry Smith const char *const PetscSFWindowSyncTypes[] = {"FENCE","LOCK","ACTIVE","PetscSFWindowSyncType","PETSCSF_WINDOW_SYNC_",0}; 2995fce210SBarry Smith 3095fce210SBarry Smith #undef __FUNCT__ 3195fce210SBarry Smith #define __FUNCT__ "PetscSFWindowOpTranslate" 3295fce210SBarry Smith /* Built-in MPI_Ops act elementwise inside MPI_Accumulate, but cannot be used with composite types inside collectives (MPI_Allreduce) */ 3395fce210SBarry Smith static PetscErrorCode PetscSFWindowOpTranslate(MPI_Op *op) 3495fce210SBarry Smith { 3595fce210SBarry Smith 3695fce210SBarry Smith PetscFunctionBegin; 3795fce210SBarry Smith if (*op == MPIU_SUM) *op = MPI_SUM; 3895fce210SBarry Smith else if (*op == MPIU_MAX) *op = MPI_MAX; 3995fce210SBarry Smith else if (*op == MPIU_MIN) *op = MPI_MIN; 4095fce210SBarry Smith PetscFunctionReturn(0); 4195fce210SBarry Smith } 4295fce210SBarry Smith 4395fce210SBarry Smith #undef __FUNCT__ 4495fce210SBarry Smith #define __FUNCT__ "PetscSFWindowGetDataTypes" 4595fce210SBarry Smith /*@C 4695fce210SBarry Smith PetscSFWindowGetDataTypes - gets composite local and remote data types for each rank 4795fce210SBarry Smith 4895fce210SBarry Smith Not Collective 4995fce210SBarry Smith 5095fce210SBarry Smith Input Arguments: 5195fce210SBarry Smith + sf - star forest 5295fce210SBarry Smith - unit - data type for each node 5395fce210SBarry Smith 5495fce210SBarry Smith Output Arguments: 5595fce210SBarry Smith + localtypes - types describing part of local leaf buffer referencing each remote rank 5695fce210SBarry Smith - remotetypes - types describing part of remote root buffer referenced for each remote rank 5795fce210SBarry Smith 5895fce210SBarry Smith Level: developer 5995fce210SBarry Smith 6095fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFView() 6195fce210SBarry Smith @*/ 6295fce210SBarry Smith static PetscErrorCode PetscSFWindowGetDataTypes(PetscSF sf,MPI_Datatype unit,const MPI_Datatype **localtypes,const MPI_Datatype **remotetypes) 6395fce210SBarry Smith { 6495fce210SBarry Smith PetscSF_Window *w = (PetscSF_Window*)sf->data; 6595fce210SBarry Smith PetscErrorCode ierr; 6695fce210SBarry Smith PetscSFDataLink link; 6795fce210SBarry Smith PetscInt i,nranks; 6895fce210SBarry Smith const PetscInt *roffset,*rmine,*rremote; 6995fce210SBarry Smith const PetscMPIInt *ranks; 7095fce210SBarry Smith 7195fce210SBarry Smith PetscFunctionBegin; 7295fce210SBarry Smith /* Look for types in cache */ 7395fce210SBarry Smith for (link=w->link; link; link=link->next) { 7495fce210SBarry Smith PetscBool match; 7595fce210SBarry Smith ierr = MPIPetsc_Type_compare(unit,link->unit,&match);CHKERRQ(ierr); 7695fce210SBarry Smith if (match) { 7795fce210SBarry Smith *localtypes = link->mine; 7895fce210SBarry Smith *remotetypes = link->remote; 7995fce210SBarry Smith PetscFunctionReturn(0); 8095fce210SBarry Smith } 8195fce210SBarry Smith } 8295fce210SBarry Smith 8395fce210SBarry Smith /* Create new composite types for each send rank */ 8495fce210SBarry Smith ierr = PetscSFGetRanks(sf,&nranks,&ranks,&roffset,&rmine,&rremote);CHKERRQ(ierr); 8595fce210SBarry Smith ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr); 8695fce210SBarry Smith ierr = MPI_Type_dup(unit,&link->unit);CHKERRQ(ierr); 87dcca6d9dSJed Brown ierr = PetscMalloc2(nranks,&link->mine,nranks,&link->remote);CHKERRQ(ierr); 8895fce210SBarry Smith for (i=0; i<nranks; i++) { 8995fce210SBarry Smith PETSC_UNUSED PetscInt rcount = roffset[i+1] - roffset[i]; 9095fce210SBarry Smith PetscMPIInt *rmine,*rremote; 9195fce210SBarry Smith #if !defined(PETSC_USE_64BIT_INDICES) 9295fce210SBarry Smith rmine = sf->rmine + sf->roffset[i]; 9395fce210SBarry Smith rremote = sf->rremote + sf->roffset[i]; 9495fce210SBarry Smith #else 9595fce210SBarry Smith PetscInt j; 96dcca6d9dSJed Brown ierr = PetscMalloc2(rcount,&rmine,rcount,&rremote);CHKERRQ(ierr); 9795fce210SBarry Smith for (j=0; j<rcount; j++) { 9895fce210SBarry Smith ierr = PetscMPIIntCast(sf->rmine[sf->roffset[i]+j],rmine+j);CHKERRQ(ierr); 9995fce210SBarry Smith ierr = PetscMPIIntCast(sf->rremote[sf->roffset[i]+j],rremote+j);CHKERRQ(ierr); 10095fce210SBarry Smith } 10195fce210SBarry Smith #endif 10295fce210SBarry Smith ierr = MPI_Type_create_indexed_block(rcount,1,rmine,link->unit,&link->mine[i]);CHKERRQ(ierr); 10395fce210SBarry Smith ierr = MPI_Type_create_indexed_block(rcount,1,rremote,link->unit,&link->remote[i]);CHKERRQ(ierr); 10495fce210SBarry Smith #if defined(PETSC_USE_64BIT_INDICES) 10595fce210SBarry Smith ierr = PetscFree2(rmine,rremote);CHKERRQ(ierr); 10695fce210SBarry Smith #endif 10795fce210SBarry Smith ierr = MPI_Type_commit(&link->mine[i]);CHKERRQ(ierr); 10895fce210SBarry Smith ierr = MPI_Type_commit(&link->remote[i]);CHKERRQ(ierr); 10995fce210SBarry Smith } 11095fce210SBarry Smith link->next = w->link; 11195fce210SBarry Smith w->link = link; 11295fce210SBarry Smith 11395fce210SBarry Smith *localtypes = link->mine; 11495fce210SBarry Smith *remotetypes = link->remote; 11595fce210SBarry Smith PetscFunctionReturn(0); 11695fce210SBarry Smith } 11795fce210SBarry Smith 11895fce210SBarry Smith #undef __FUNCT__ 11995fce210SBarry Smith #define __FUNCT__ "PetscSFWindowSetSyncType" 12095fce210SBarry Smith /*@C 12195fce210SBarry Smith PetscSFWindowSetSyncType - set synchrozitaion type for PetscSF communication 12295fce210SBarry Smith 12395fce210SBarry Smith Logically Collective 12495fce210SBarry Smith 12595fce210SBarry Smith Input Arguments: 12695fce210SBarry Smith + sf - star forest for communication 12795fce210SBarry Smith - sync - synchronization type 12895fce210SBarry Smith 12995fce210SBarry Smith Options Database Key: 13095fce210SBarry Smith . -sf_synchronization <sync> - sets the synchronization type 13195fce210SBarry Smith 13295fce210SBarry Smith Level: advanced 13395fce210SBarry Smith 13495fce210SBarry Smith .seealso: PetscSFSetFromOptions(), PetscSFWindowGetSyncType() 13595fce210SBarry Smith @*/ 13695fce210SBarry Smith PetscErrorCode PetscSFWindowSetSyncType(PetscSF sf,PetscSFWindowSyncType sync) 13795fce210SBarry Smith { 13895fce210SBarry Smith PetscErrorCode ierr; 13995fce210SBarry Smith 14095fce210SBarry Smith PetscFunctionBegin; 14195fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 14295fce210SBarry Smith PetscValidLogicalCollectiveEnum(sf,sync,2); 14395fce210SBarry Smith ierr = PetscUseMethod(sf,"PetscSFWindowSetSyncType_C",(PetscSF,PetscSFWindowSyncType),(sf,sync));CHKERRQ(ierr); 14495fce210SBarry Smith PetscFunctionReturn(0); 14595fce210SBarry Smith } 14695fce210SBarry Smith 14795fce210SBarry Smith #undef __FUNCT__ 14895fce210SBarry Smith #define __FUNCT__ "PetscSFWindowSetSyncType_Window" 149f7a08781SBarry Smith static PetscErrorCode PetscSFWindowSetSyncType_Window(PetscSF sf,PetscSFWindowSyncType sync) 15095fce210SBarry Smith { 15195fce210SBarry Smith PetscSF_Window *w = (PetscSF_Window*)sf->data; 15295fce210SBarry Smith 15395fce210SBarry Smith PetscFunctionBegin; 15495fce210SBarry Smith w->sync = sync; 15595fce210SBarry Smith PetscFunctionReturn(0); 15695fce210SBarry Smith } 15795fce210SBarry Smith 15895fce210SBarry Smith #undef __FUNCT__ 15995fce210SBarry Smith #define __FUNCT__ "PetscSFWindowGetSyncType" 16095fce210SBarry Smith /*@C 16195fce210SBarry Smith PetscSFWindowGetSyncType - get synchrozitaion type for PetscSF communication 16295fce210SBarry Smith 16395fce210SBarry Smith Logically Collective 16495fce210SBarry Smith 16595fce210SBarry Smith Input Argument: 16695fce210SBarry Smith . sf - star forest for communication 16795fce210SBarry Smith 16895fce210SBarry Smith Output Argument: 16995fce210SBarry Smith . sync - synchronization type 17095fce210SBarry Smith 17195fce210SBarry Smith Level: advanced 17295fce210SBarry Smith 17395fce210SBarry Smith .seealso: PetscSFGetFromOptions(), PetscSFWindowSetSyncType() 17495fce210SBarry Smith @*/ 17595fce210SBarry Smith PetscErrorCode PetscSFWindowGetSyncType(PetscSF sf,PetscSFWindowSyncType *sync) 17695fce210SBarry Smith { 17795fce210SBarry Smith PetscErrorCode ierr; 17895fce210SBarry Smith 17995fce210SBarry Smith PetscFunctionBegin; 18095fce210SBarry Smith PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 18195fce210SBarry Smith PetscValidPointer(sync,2); 18295fce210SBarry Smith ierr = PetscTryMethod(sf,"PetscSFWindowGetSyncType_C",(PetscSF,PetscSFWindowSyncType*),(sf,sync));CHKERRQ(ierr); 18395fce210SBarry Smith PetscFunctionReturn(0); 18495fce210SBarry Smith } 18595fce210SBarry Smith 18695fce210SBarry Smith #undef __FUNCT__ 18795fce210SBarry Smith #define __FUNCT__ "PetscSFWindowGetSyncType_Window" 188f7a08781SBarry Smith static PetscErrorCode PetscSFWindowGetSyncType_Window(PetscSF sf,PetscSFWindowSyncType *sync) 18995fce210SBarry Smith { 19095fce210SBarry Smith PetscSF_Window *w = (PetscSF_Window*)sf->data; 19195fce210SBarry Smith 19295fce210SBarry Smith PetscFunctionBegin; 19395fce210SBarry Smith *sync = w->sync; 19495fce210SBarry Smith PetscFunctionReturn(0); 19595fce210SBarry Smith } 19695fce210SBarry Smith 19795fce210SBarry Smith #undef __FUNCT__ 19895fce210SBarry Smith #define __FUNCT__ "PetscSFGetWindow" 19995fce210SBarry Smith /*@C 20095fce210SBarry Smith PetscSFGetWindow - Get a window for use with a given data type 20195fce210SBarry Smith 20295fce210SBarry Smith Collective on PetscSF 20395fce210SBarry Smith 20495fce210SBarry Smith Input Arguments: 20595fce210SBarry Smith + sf - star forest 20695fce210SBarry Smith . unit - data type 20795fce210SBarry Smith . array - array to be sent 20895fce210SBarry Smith . epoch - PETSC_TRUE to acquire the window and start an epoch, PETSC_FALSE to just acquire the window 20995fce210SBarry Smith . fenceassert - assert parameter for call to MPI_Win_fence(), if PETSCSF_WINDOW_SYNC_FENCE 21095fce210SBarry Smith . postassert - assert parameter for call to MPI_Win_post(), if PETSCSF_WINDOW_SYNC_ACTIVE 21195fce210SBarry Smith - startassert - assert parameter for call to MPI_Win_start(), if PETSCSF_WINDOW_SYNC_ACTIVE 21295fce210SBarry Smith 21395fce210SBarry Smith Output Arguments: 21495fce210SBarry Smith . win - window 21595fce210SBarry Smith 21695fce210SBarry Smith Level: developer 21795fce210SBarry Smith 21895fce210SBarry Smith Developer Notes: 21995fce210SBarry Smith This currently always creates a new window. This is more synchronous than necessary. An alternative is to try to 22095fce210SBarry Smith reuse an existing window created with the same array. Another alternative is to maintain a cache of windows and reuse 22195fce210SBarry Smith whichever one is available, by copying the array into it if necessary. 22295fce210SBarry Smith 22395fce210SBarry Smith .seealso: PetscSFGetRanks(), PetscSFWindowGetDataTypes() 22495fce210SBarry Smith @*/ 22595fce210SBarry Smith static PetscErrorCode PetscSFGetWindow(PetscSF sf,MPI_Datatype unit,void *array,PetscBool epoch,PetscMPIInt fenceassert,PetscMPIInt postassert,PetscMPIInt startassert,MPI_Win *win) 22695fce210SBarry Smith { 22795fce210SBarry Smith PetscSF_Window *w = (PetscSF_Window*)sf->data; 22895fce210SBarry Smith PetscErrorCode ierr; 22995fce210SBarry Smith MPI_Aint lb,lb_true,bytes,bytes_true; 23095fce210SBarry Smith PetscSFWinLink link; 23195fce210SBarry Smith 23295fce210SBarry Smith PetscFunctionBegin; 23395fce210SBarry Smith ierr = MPI_Type_get_extent(unit,&lb,&bytes);CHKERRQ(ierr); 23495fce210SBarry Smith ierr = MPI_Type_get_true_extent(unit,&lb_true,&bytes_true);CHKERRQ(ierr); 23595fce210SBarry Smith 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"); 23695fce210SBarry Smith 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"); 23795fce210SBarry Smith ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr); 23895fce210SBarry Smith 23995fce210SBarry Smith link->bytes = bytes; 24095fce210SBarry Smith link->addr = array; 24195fce210SBarry Smith 24295fce210SBarry Smith ierr = MPI_Win_create(array,(MPI_Aint)bytes*sf->nroots,(PetscMPIInt)bytes,MPI_INFO_NULL,PetscObjectComm((PetscObject)sf),&link->win);CHKERRQ(ierr); 24395fce210SBarry Smith 24495fce210SBarry Smith link->epoch = epoch; 24595fce210SBarry Smith link->next = w->wins; 24695fce210SBarry Smith link->inuse = PETSC_TRUE; 24795fce210SBarry Smith w->wins = link; 24895fce210SBarry Smith *win = link->win; 24995fce210SBarry Smith 25095fce210SBarry Smith if (epoch) { 25195fce210SBarry Smith switch (w->sync) { 25295fce210SBarry Smith case PETSCSF_WINDOW_SYNC_FENCE: 25395fce210SBarry Smith ierr = MPI_Win_fence(fenceassert,*win);CHKERRQ(ierr); 25495fce210SBarry Smith break; 25595fce210SBarry Smith case PETSCSF_WINDOW_SYNC_LOCK: /* Handled outside */ 25695fce210SBarry Smith break; 25795fce210SBarry Smith case PETSCSF_WINDOW_SYNC_ACTIVE: { 25895fce210SBarry Smith MPI_Group ingroup,outgroup; 25995fce210SBarry Smith ierr = PetscSFGetGroups(sf,&ingroup,&outgroup);CHKERRQ(ierr); 26095fce210SBarry Smith ierr = MPI_Win_post(ingroup,postassert,*win);CHKERRQ(ierr); 26195fce210SBarry Smith ierr = MPI_Win_start(outgroup,startassert,*win);CHKERRQ(ierr); 26295fce210SBarry Smith } break; 26395fce210SBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_PLIB,"Unknown synchronization type"); 26495fce210SBarry Smith } 26595fce210SBarry Smith } 26695fce210SBarry Smith PetscFunctionReturn(0); 26795fce210SBarry Smith } 26895fce210SBarry Smith 26995fce210SBarry Smith #undef __FUNCT__ 27095fce210SBarry Smith #define __FUNCT__ "PetscSFFindWindow" 27195fce210SBarry Smith /*@C 27295fce210SBarry Smith PetscSFFindWindow - Finds a window that is already in use 27395fce210SBarry Smith 27495fce210SBarry Smith Not Collective 27595fce210SBarry Smith 27695fce210SBarry Smith Input Arguments: 27795fce210SBarry Smith + sf - star forest 27895fce210SBarry Smith . unit - data type 27995fce210SBarry Smith - array - array with which the window is associated 28095fce210SBarry Smith 28195fce210SBarry Smith Output Arguments: 28295fce210SBarry Smith . win - window 28395fce210SBarry Smith 28495fce210SBarry Smith Level: developer 28595fce210SBarry Smith 28695fce210SBarry Smith .seealso: PetscSFGetWindow(), PetscSFRestoreWindow() 28795fce210SBarry Smith @*/ 28895fce210SBarry Smith static PetscErrorCode PetscSFFindWindow(PetscSF sf,MPI_Datatype unit,const void *array,MPI_Win *win) 28995fce210SBarry Smith { 29095fce210SBarry Smith PetscSF_Window *w = (PetscSF_Window*)sf->data; 29195fce210SBarry Smith PetscSFWinLink link; 29295fce210SBarry Smith 29395fce210SBarry Smith PetscFunctionBegin; 29495fce210SBarry Smith for (link=w->wins; link; link=link->next) { 29595fce210SBarry Smith if (array == link->addr) { 29695fce210SBarry Smith *win = link->win; 29795fce210SBarry Smith PetscFunctionReturn(0); 29895fce210SBarry Smith } 29995fce210SBarry Smith } 30095fce210SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Requested window not in use"); 30195fce210SBarry Smith PetscFunctionReturn(0); 30295fce210SBarry Smith } 30395fce210SBarry Smith 30495fce210SBarry Smith #undef __FUNCT__ 30595fce210SBarry Smith #define __FUNCT__ "PetscSFRestoreWindow" 30695fce210SBarry Smith /*@C 30795fce210SBarry Smith PetscSFRestoreWindow - Restores a window obtained with PetscSFGetWindow() 30895fce210SBarry Smith 30995fce210SBarry Smith Collective 31095fce210SBarry Smith 31195fce210SBarry Smith Input Arguments: 31295fce210SBarry Smith + sf - star forest 31395fce210SBarry Smith . unit - data type 31495fce210SBarry Smith . array - array associated with window 31595fce210SBarry Smith . epoch - close an epoch, must match argument to PetscSFGetWindow() 31695fce210SBarry Smith - win - window 31795fce210SBarry Smith 31895fce210SBarry Smith Level: developer 31995fce210SBarry Smith 32095fce210SBarry Smith .seealso: PetscSFFindWindow() 32195fce210SBarry Smith @*/ 32295fce210SBarry Smith static PetscErrorCode PetscSFRestoreWindow(PetscSF sf,MPI_Datatype unit,const void *array,PetscBool epoch,PetscMPIInt fenceassert,MPI_Win *win) 32395fce210SBarry Smith { 32495fce210SBarry Smith PetscSF_Window *w = (PetscSF_Window*)sf->data; 32595fce210SBarry Smith PetscErrorCode ierr; 32695fce210SBarry Smith PetscSFWinLink *p,link; 32795fce210SBarry Smith 32895fce210SBarry Smith PetscFunctionBegin; 32995fce210SBarry Smith for (p=&w->wins; *p; p=&(*p)->next) { 33095fce210SBarry Smith link = *p; 33195fce210SBarry Smith if (*win == link->win) { 33295fce210SBarry Smith if (array != link->addr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Matched window, but not array"); 33395fce210SBarry Smith if (epoch != link->epoch) { 33495fce210SBarry Smith if (epoch) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"No epoch to end"); 33595fce210SBarry Smith else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Restoring window without ending epoch"); 33695fce210SBarry Smith } 33795fce210SBarry Smith *p = link->next; 33895fce210SBarry Smith goto found; 33995fce210SBarry Smith } 34095fce210SBarry Smith } 34195fce210SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Requested window not in use"); 34295fce210SBarry Smith 34395fce210SBarry Smith found: 34495fce210SBarry Smith if (epoch) { 34595fce210SBarry Smith switch (w->sync) { 34695fce210SBarry Smith case PETSCSF_WINDOW_SYNC_FENCE: 34795fce210SBarry Smith ierr = MPI_Win_fence(fenceassert,*win);CHKERRQ(ierr); 34895fce210SBarry Smith break; 34995fce210SBarry Smith case PETSCSF_WINDOW_SYNC_LOCK: 35095fce210SBarry Smith break; /* handled outside */ 35195fce210SBarry Smith case PETSCSF_WINDOW_SYNC_ACTIVE: { 35295fce210SBarry Smith ierr = MPI_Win_complete(*win);CHKERRQ(ierr); 35395fce210SBarry Smith ierr = MPI_Win_wait(*win);CHKERRQ(ierr); 35495fce210SBarry Smith } break; 35595fce210SBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_PLIB,"Unknown synchronization type"); 35695fce210SBarry Smith } 35795fce210SBarry Smith } 35895fce210SBarry Smith 35995fce210SBarry Smith ierr = MPI_Win_free(&link->win);CHKERRQ(ierr); 36095fce210SBarry Smith ierr = PetscFree(link);CHKERRQ(ierr); 36195fce210SBarry Smith *win = MPI_WIN_NULL; 36295fce210SBarry Smith PetscFunctionReturn(0); 36395fce210SBarry Smith } 36495fce210SBarry Smith 36595fce210SBarry Smith #undef __FUNCT__ 36695fce210SBarry Smith #define __FUNCT__ "PetscSFSetUp_Window" 36795fce210SBarry Smith static PetscErrorCode PetscSFSetUp_Window(PetscSF sf) 36895fce210SBarry Smith { 36995fce210SBarry Smith PetscSF_Window *w = (PetscSF_Window*)sf->data; 37095fce210SBarry Smith PetscErrorCode ierr; 37195fce210SBarry Smith MPI_Group ingroup,outgroup; 37295fce210SBarry Smith 37395fce210SBarry Smith PetscFunctionBegin; 37495fce210SBarry Smith switch (w->sync) { 37595fce210SBarry Smith case PETSCSF_WINDOW_SYNC_ACTIVE: 37695fce210SBarry Smith ierr = PetscSFGetGroups(sf,&ingroup,&outgroup);CHKERRQ(ierr); 37795fce210SBarry Smith default: 37895fce210SBarry Smith break; 37995fce210SBarry Smith } 38095fce210SBarry Smith PetscFunctionReturn(0); 38195fce210SBarry Smith } 38295fce210SBarry Smith 38395fce210SBarry Smith #undef __FUNCT__ 38495fce210SBarry Smith #define __FUNCT__ "PetscSFSetFromOptions_Window" 38595fce210SBarry Smith static PetscErrorCode PetscSFSetFromOptions_Window(PetscSF sf) 38695fce210SBarry Smith { 38795fce210SBarry Smith PetscSF_Window *w = (PetscSF_Window*)sf->data; 38895fce210SBarry Smith PetscErrorCode ierr; 38995fce210SBarry Smith 39095fce210SBarry Smith PetscFunctionBegin; 39195fce210SBarry Smith ierr = PetscOptionsHead("PetscSF Window options");CHKERRQ(ierr); 39295fce210SBarry Smith ierr = PetscOptionsEnum("-sf_window_sync","synchronization type to use for PetscSF Window communication","PetscSFWindowSetSyncType",PetscSFWindowSyncTypes,(PetscEnum)w->sync,(PetscEnum*)&w->sync,NULL);CHKERRQ(ierr); 39395fce210SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 39495fce210SBarry Smith PetscFunctionReturn(0); 39595fce210SBarry Smith } 39695fce210SBarry Smith 39795fce210SBarry Smith #undef __FUNCT__ 39895fce210SBarry Smith #define __FUNCT__ "PetscSFReset_Window" 39995fce210SBarry Smith static PetscErrorCode PetscSFReset_Window(PetscSF sf) 40095fce210SBarry Smith { 40195fce210SBarry Smith PetscSF_Window *w = (PetscSF_Window*)sf->data; 40295fce210SBarry Smith PetscErrorCode ierr; 40395fce210SBarry Smith PetscSFDataLink link,next; 40495fce210SBarry Smith PetscSFWinLink wlink,wnext; 40595fce210SBarry Smith PetscInt i; 40695fce210SBarry Smith 40795fce210SBarry Smith PetscFunctionBegin; 40895fce210SBarry Smith for (link=w->link; link; link=next) { 40995fce210SBarry Smith next = link->next; 41095fce210SBarry Smith ierr = MPI_Type_free(&link->unit);CHKERRQ(ierr); 41195fce210SBarry Smith for (i=0; i<sf->nranks; i++) { 41295fce210SBarry Smith ierr = MPI_Type_free(&link->mine[i]);CHKERRQ(ierr); 41395fce210SBarry Smith ierr = MPI_Type_free(&link->remote[i]);CHKERRQ(ierr); 41495fce210SBarry Smith } 41595fce210SBarry Smith ierr = PetscFree2(link->mine,link->remote);CHKERRQ(ierr); 41695fce210SBarry Smith ierr = PetscFree(link);CHKERRQ(ierr); 41795fce210SBarry Smith } 41895fce210SBarry Smith w->link = NULL; 41995fce210SBarry Smith for (wlink=w->wins; wlink; wlink=wnext) { 42095fce210SBarry Smith wnext = wlink->next; 42195fce210SBarry Smith if (wlink->inuse) SETERRQ1(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Window still in use with address %p",(void*)wlink->addr); 42295fce210SBarry Smith ierr = MPI_Win_free(&wlink->win);CHKERRQ(ierr); 42395fce210SBarry Smith ierr = PetscFree(wlink);CHKERRQ(ierr); 42495fce210SBarry Smith } 42595fce210SBarry Smith w->wins = NULL; 42695fce210SBarry Smith PetscFunctionReturn(0); 42795fce210SBarry Smith } 42895fce210SBarry Smith 42995fce210SBarry Smith #undef __FUNCT__ 43095fce210SBarry Smith #define __FUNCT__ "PetscSFDestroy_Window" 43195fce210SBarry Smith static PetscErrorCode PetscSFDestroy_Window(PetscSF sf) 43295fce210SBarry Smith { 43395fce210SBarry Smith PetscErrorCode ierr; 43495fce210SBarry Smith 43595fce210SBarry Smith PetscFunctionBegin; 43695fce210SBarry Smith ierr = PetscSFReset_Window(sf);CHKERRQ(ierr); 43795fce210SBarry Smith ierr = PetscFree(sf->data);CHKERRQ(ierr); 438bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)sf,"PetscSFWindowSetSyncType_C",NULL);CHKERRQ(ierr); 439bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)sf,"PetscSFWindowGetSyncType_C",NULL);CHKERRQ(ierr); 44095fce210SBarry Smith PetscFunctionReturn(0); 44195fce210SBarry Smith } 44295fce210SBarry Smith 44395fce210SBarry Smith #undef __FUNCT__ 44495fce210SBarry Smith #define __FUNCT__ "PetscSFView_Window" 44595fce210SBarry Smith static PetscErrorCode PetscSFView_Window(PetscSF sf,PetscViewer viewer) 44695fce210SBarry Smith { 44795fce210SBarry Smith PetscSF_Window *w = (PetscSF_Window*)sf->data; 44895fce210SBarry Smith PetscErrorCode ierr; 44995fce210SBarry Smith PetscBool iascii; 45095fce210SBarry Smith 45195fce210SBarry Smith PetscFunctionBegin; 45295fce210SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 45395fce210SBarry Smith if (iascii) { 45495fce210SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," synchronization=%s sort=%s\n",PetscSFWindowSyncTypes[w->sync],sf->rankorder ? "rank-order" : "unordered");CHKERRQ(ierr); 45595fce210SBarry Smith } 45695fce210SBarry Smith PetscFunctionReturn(0); 45795fce210SBarry Smith } 45895fce210SBarry Smith 45995fce210SBarry Smith #undef __FUNCT__ 46095fce210SBarry Smith #define __FUNCT__ "PetscSFDuplicate_Window" 46195fce210SBarry Smith static PetscErrorCode PetscSFDuplicate_Window(PetscSF sf,PetscSFDuplicateOption opt,PetscSF newsf) 46295fce210SBarry Smith { 46395fce210SBarry Smith PetscSF_Window *w = (PetscSF_Window*)sf->data; 46495fce210SBarry Smith PetscErrorCode ierr; 46595fce210SBarry Smith PetscSFWindowSyncType synctype; 46695fce210SBarry Smith 46795fce210SBarry Smith PetscFunctionBegin; 46895fce210SBarry Smith synctype = w->sync; 46995fce210SBarry Smith if (!sf->setupcalled) { 47095fce210SBarry Smith /* HACK: Must use FENCE or LOCK when called from PetscSFGetGroups() because ACTIVE here would cause recursion. */ 47195fce210SBarry Smith synctype = PETSCSF_WINDOW_SYNC_LOCK; 47295fce210SBarry Smith } 47395fce210SBarry Smith ierr = PetscSFWindowSetSyncType(newsf,synctype);CHKERRQ(ierr); 47495fce210SBarry Smith PetscFunctionReturn(0); 47595fce210SBarry Smith } 47695fce210SBarry Smith 47795fce210SBarry Smith #undef __FUNCT__ 47895fce210SBarry Smith #define __FUNCT__ "PetscSFBcastBegin_Window" 47995fce210SBarry Smith static PetscErrorCode PetscSFBcastBegin_Window(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata) 48095fce210SBarry Smith { 48195fce210SBarry Smith PetscSF_Window *w = (PetscSF_Window*)sf->data; 48295fce210SBarry Smith PetscErrorCode ierr; 48395fce210SBarry Smith PetscInt i,nranks; 48495fce210SBarry Smith const PetscMPIInt *ranks; 48595fce210SBarry Smith const MPI_Datatype *mine,*remote; 48695fce210SBarry Smith MPI_Win win; 48795fce210SBarry Smith 48895fce210SBarry Smith PetscFunctionBegin; 48995fce210SBarry Smith ierr = PetscSFGetRanks(sf,&nranks,&ranks,NULL,NULL,NULL);CHKERRQ(ierr); 49095fce210SBarry Smith ierr = PetscSFWindowGetDataTypes(sf,unit,&mine,&remote);CHKERRQ(ierr); 49195fce210SBarry Smith ierr = PetscSFGetWindow(sf,unit,(void*)rootdata,PETSC_TRUE,MPI_MODE_NOPUT|MPI_MODE_NOPRECEDE,MPI_MODE_NOPUT,0,&win);CHKERRQ(ierr); 49295fce210SBarry Smith for (i=0; i<nranks; i++) { 49395fce210SBarry Smith if (w->sync == PETSCSF_WINDOW_SYNC_LOCK) {ierr = MPI_Win_lock(MPI_LOCK_SHARED,ranks[i],MPI_MODE_NOCHECK,win);CHKERRQ(ierr);} 49495fce210SBarry Smith ierr = MPI_Get(leafdata,1,mine[i],ranks[i],0,1,remote[i],win);CHKERRQ(ierr); 49595fce210SBarry Smith if (w->sync == PETSCSF_WINDOW_SYNC_LOCK) {ierr = MPI_Win_unlock(ranks[i],win);CHKERRQ(ierr);} 49695fce210SBarry Smith } 49795fce210SBarry Smith PetscFunctionReturn(0); 49895fce210SBarry Smith } 49995fce210SBarry Smith 50095fce210SBarry Smith #undef __FUNCT__ 50195fce210SBarry Smith #define __FUNCT__ "PetscSFBcastEnd_Window" 50295fce210SBarry Smith PetscErrorCode PetscSFBcastEnd_Window(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata) 50395fce210SBarry Smith { 50495fce210SBarry Smith PetscErrorCode ierr; 50595fce210SBarry Smith MPI_Win win; 50695fce210SBarry Smith 50795fce210SBarry Smith PetscFunctionBegin; 50895fce210SBarry Smith ierr = PetscSFFindWindow(sf,unit,rootdata,&win);CHKERRQ(ierr); 50995fce210SBarry Smith ierr = PetscSFRestoreWindow(sf,unit,rootdata,PETSC_TRUE,MPI_MODE_NOSTORE|MPI_MODE_NOSUCCEED,&win);CHKERRQ(ierr); 51095fce210SBarry Smith PetscFunctionReturn(0); 51195fce210SBarry Smith } 51295fce210SBarry Smith 51395fce210SBarry Smith #undef __FUNCT__ 51495fce210SBarry Smith #define __FUNCT__ "PetscSFReduceBegin_Window" 51595fce210SBarry Smith PetscErrorCode PetscSFReduceBegin_Window(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op) 51695fce210SBarry Smith { 51795fce210SBarry Smith PetscSF_Window *w = (PetscSF_Window*)sf->data; 51895fce210SBarry Smith PetscErrorCode ierr; 51995fce210SBarry Smith PetscInt i,nranks; 52095fce210SBarry Smith const PetscMPIInt *ranks; 52195fce210SBarry Smith const MPI_Datatype *mine,*remote; 52295fce210SBarry Smith MPI_Win win; 52395fce210SBarry Smith 52495fce210SBarry Smith PetscFunctionBegin; 52595fce210SBarry Smith ierr = PetscSFGetRanks(sf,&nranks,&ranks,NULL,NULL,NULL);CHKERRQ(ierr); 52695fce210SBarry Smith ierr = PetscSFWindowGetDataTypes(sf,unit,&mine,&remote);CHKERRQ(ierr); 52795fce210SBarry Smith ierr = PetscSFWindowOpTranslate(&op);CHKERRQ(ierr); 52895fce210SBarry Smith ierr = PetscSFGetWindow(sf,unit,rootdata,PETSC_TRUE,MPI_MODE_NOPRECEDE,0,0,&win);CHKERRQ(ierr); 52995fce210SBarry Smith for (i=0; i<nranks; i++) { 53095fce210SBarry Smith if (w->sync == PETSCSF_WINDOW_SYNC_LOCK) {ierr = MPI_Win_lock(MPI_LOCK_SHARED,ranks[i],MPI_MODE_NOCHECK,win);CHKERRQ(ierr);} 53195fce210SBarry Smith ierr = MPI_Accumulate((void*)leafdata,1,mine[i],ranks[i],0,1,remote[i],op,win);CHKERRQ(ierr); 53295fce210SBarry Smith if (w->sync == PETSCSF_WINDOW_SYNC_LOCK) {ierr = MPI_Win_unlock(ranks[i],win);CHKERRQ(ierr);} 53395fce210SBarry Smith } 53495fce210SBarry Smith PetscFunctionReturn(0); 53595fce210SBarry Smith } 53695fce210SBarry Smith 53795fce210SBarry Smith #undef __FUNCT__ 53895fce210SBarry Smith #define __FUNCT__ "PetscSFReduceEnd_Window" 53995fce210SBarry Smith static PetscErrorCode PetscSFReduceEnd_Window(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op) 54095fce210SBarry Smith { 54195fce210SBarry Smith PetscSF_Window *w = (PetscSF_Window*)sf->data; 54295fce210SBarry Smith PetscErrorCode ierr; 54395fce210SBarry Smith MPI_Win win; 54495fce210SBarry Smith 54595fce210SBarry Smith PetscFunctionBegin; 54695fce210SBarry Smith if (!w->wins) PetscFunctionReturn(0); 54795fce210SBarry Smith ierr = PetscSFFindWindow(sf,unit,rootdata,&win);CHKERRQ(ierr); 54895fce210SBarry Smith ierr = MPI_Win_fence(MPI_MODE_NOSUCCEED,win);CHKERRQ(ierr); 54995fce210SBarry Smith ierr = PetscSFRestoreWindow(sf,unit,rootdata,PETSC_TRUE,MPI_MODE_NOSUCCEED,&win);CHKERRQ(ierr); 55095fce210SBarry Smith PetscFunctionReturn(0); 55195fce210SBarry Smith } 55295fce210SBarry Smith #undef __FUNCT__ 55395fce210SBarry Smith #define __FUNCT__ "PetscSFFetchAndOpBegin_Window" 55495fce210SBarry Smith static PetscErrorCode PetscSFFetchAndOpBegin_Window(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op) 55595fce210SBarry Smith { 55695fce210SBarry Smith PetscErrorCode ierr; 55795fce210SBarry Smith PetscInt i,nranks; 55895fce210SBarry Smith const PetscMPIInt *ranks; 55995fce210SBarry Smith const MPI_Datatype *mine,*remote; 56095fce210SBarry Smith MPI_Win win; 56195fce210SBarry Smith 56295fce210SBarry Smith PetscFunctionBegin; 56395fce210SBarry Smith ierr = PetscSFGetRanks(sf,&nranks,&ranks,NULL,NULL,NULL);CHKERRQ(ierr); 56495fce210SBarry Smith ierr = PetscSFWindowGetDataTypes(sf,unit,&mine,&remote);CHKERRQ(ierr); 56595fce210SBarry Smith ierr = PetscSFWindowOpTranslate(&op);CHKERRQ(ierr); 56695fce210SBarry Smith ierr = PetscSFGetWindow(sf,unit,rootdata,PETSC_FALSE,0,0,0,&win);CHKERRQ(ierr); 56795fce210SBarry Smith for (i=0; i<sf->nranks; i++) { 56895fce210SBarry Smith ierr = MPI_Win_lock(MPI_LOCK_EXCLUSIVE,sf->ranks[i],0,win);CHKERRQ(ierr); 56995fce210SBarry Smith ierr = MPI_Get(leafupdate,1,mine[i],ranks[i],0,1,remote[i],win);CHKERRQ(ierr); 57095fce210SBarry Smith ierr = MPI_Accumulate((void*)leafdata,1,mine[i],ranks[i],0,1,remote[i],op,win);CHKERRQ(ierr); 57195fce210SBarry Smith ierr = MPI_Win_unlock(ranks[i],win);CHKERRQ(ierr); 57295fce210SBarry Smith } 57395fce210SBarry Smith PetscFunctionReturn(0); 57495fce210SBarry Smith } 57595fce210SBarry Smith 57695fce210SBarry Smith #undef __FUNCT__ 57795fce210SBarry Smith #define __FUNCT__ "PetscSFFetchAndOpEnd_Window" 57895fce210SBarry Smith static PetscErrorCode PetscSFFetchAndOpEnd_Window(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op) 57995fce210SBarry Smith { 58095fce210SBarry Smith PetscErrorCode ierr; 58195fce210SBarry Smith MPI_Win win; 58295fce210SBarry Smith 58395fce210SBarry Smith PetscFunctionBegin; 58495fce210SBarry Smith ierr = PetscSFFindWindow(sf,unit,rootdata,&win);CHKERRQ(ierr); 58595fce210SBarry Smith /* Nothing to do currently because MPI_LOCK_EXCLUSIVE is used in PetscSFFetchAndOpBegin(), rendering this implementation synchronous. */ 58695fce210SBarry Smith ierr = PetscSFRestoreWindow(sf,unit,rootdata,PETSC_FALSE,0,&win);CHKERRQ(ierr); 58795fce210SBarry Smith PetscFunctionReturn(0); 58895fce210SBarry Smith } 58995fce210SBarry Smith 59095fce210SBarry Smith #undef __FUNCT__ 59195fce210SBarry Smith #define __FUNCT__ "PetscSFCreate_Window" 5928cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscSFCreate_Window(PetscSF sf) 59395fce210SBarry Smith { 59495fce210SBarry Smith PetscSF_Window *w = (PetscSF_Window*)sf->data; 59595fce210SBarry Smith PetscErrorCode ierr; 59695fce210SBarry Smith 59795fce210SBarry Smith PetscFunctionBegin; 59895fce210SBarry Smith sf->ops->SetUp = PetscSFSetUp_Window; 59995fce210SBarry Smith sf->ops->SetFromOptions = PetscSFSetFromOptions_Window; 60095fce210SBarry Smith sf->ops->Reset = PetscSFReset_Window; 60195fce210SBarry Smith sf->ops->Destroy = PetscSFDestroy_Window; 60295fce210SBarry Smith sf->ops->View = PetscSFView_Window; 60395fce210SBarry Smith sf->ops->Duplicate = PetscSFDuplicate_Window; 60495fce210SBarry Smith sf->ops->BcastBegin = PetscSFBcastBegin_Window; 60595fce210SBarry Smith sf->ops->BcastEnd = PetscSFBcastEnd_Window; 60695fce210SBarry Smith sf->ops->ReduceBegin = PetscSFReduceBegin_Window; 60795fce210SBarry Smith sf->ops->ReduceEnd = PetscSFReduceEnd_Window; 60895fce210SBarry Smith sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Window; 60995fce210SBarry Smith sf->ops->FetchAndOpEnd = PetscSFFetchAndOpEnd_Window; 61095fce210SBarry Smith 611*b00a9115SJed Brown ierr = PetscNewLog(sf,&w);CHKERRQ(ierr); 61295fce210SBarry Smith sf->data = (void*)w; 61395fce210SBarry Smith w->sync = PETSCSF_WINDOW_SYNC_FENCE; 61495fce210SBarry Smith 615bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)sf,"PetscSFWindowSetSyncType_C",PetscSFWindowSetSyncType_Window);CHKERRQ(ierr); 616bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)sf,"PetscSFWindowGetSyncType_C",PetscSFWindowGetSyncType_Window);CHKERRQ(ierr); 61795fce210SBarry Smith 61895fce210SBarry Smith #if defined(OMPI_MAJOR_VERSION) && (OMPI_MAJOR_VERSION < 1 || (OMPI_MAJOR_VERSION == 1 && OMPI_MINOR_VERSION <= 6)) 61995fce210SBarry Smith { 62095fce210SBarry Smith PetscBool ackbug = PETSC_FALSE; 62195fce210SBarry Smith ierr = PetscOptionsGetBool(NULL,"-acknowledge_ompi_onesided_bug",&ackbug,NULL);CHKERRQ(ierr); 62295fce210SBarry Smith if (ackbug) { 62395fce210SBarry Smith ierr = PetscInfo(sf,"Acknowledged Open MPI bug, proceeding anyway. Expect memory corruption.");CHKERRQ(ierr); 62495fce210SBarry Smith } 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"); 62595fce210SBarry Smith } 62695fce210SBarry Smith #endif 62795fce210SBarry Smith PetscFunctionReturn(0); 62895fce210SBarry Smith } 629