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