xref: /petsc/src/vec/is/sf/impls/window/sfwindow.c (revision 785e854f82a3c614b452fca2cf5ad4f2afe8bdde)
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_synchronization <sync> - sets the synchronization type
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   for (link=w->wins; link; link=link->next) {
295     if (array == link->addr) {
296       *win = link->win;
297       PetscFunctionReturn(0);
298     }
299   }
300   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Requested window not in use");
301   PetscFunctionReturn(0);
302 }
303 
304 #undef __FUNCT__
305 #define __FUNCT__ "PetscSFRestoreWindow"
306 /*@C
307    PetscSFRestoreWindow - Restores a window obtained with PetscSFGetWindow()
308 
309    Collective
310 
311    Input Arguments:
312 +  sf - star forest
313 .  unit - data type
314 .  array - array associated with window
315 .  epoch - close an epoch, must match argument to PetscSFGetWindow()
316 -  win - window
317 
318    Level: developer
319 
320 .seealso: PetscSFFindWindow()
321 @*/
322 static PetscErrorCode PetscSFRestoreWindow(PetscSF sf,MPI_Datatype unit,const void *array,PetscBool epoch,PetscMPIInt fenceassert,MPI_Win *win)
323 {
324   PetscSF_Window *w = (PetscSF_Window*)sf->data;
325   PetscErrorCode ierr;
326   PetscSFWinLink *p,link;
327 
328   PetscFunctionBegin;
329   for (p=&w->wins; *p; p=&(*p)->next) {
330     link = *p;
331     if (*win == link->win) {
332       if (array != link->addr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Matched window, but not array");
333       if (epoch != link->epoch) {
334         if (epoch) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"No epoch to end");
335         else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Restoring window without ending epoch");
336       }
337       *p = link->next;
338       goto found;
339     }
340   }
341   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Requested window not in use");
342 
343 found:
344   if (epoch) {
345     switch (w->sync) {
346     case PETSCSF_WINDOW_SYNC_FENCE:
347       ierr = MPI_Win_fence(fenceassert,*win);CHKERRQ(ierr);
348       break;
349     case PETSCSF_WINDOW_SYNC_LOCK:
350       break;                    /* handled outside */
351     case PETSCSF_WINDOW_SYNC_ACTIVE: {
352       ierr = MPI_Win_complete(*win);CHKERRQ(ierr);
353       ierr = MPI_Win_wait(*win);CHKERRQ(ierr);
354     } break;
355     default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_PLIB,"Unknown synchronization type");
356     }
357   }
358 
359   ierr = MPI_Win_free(&link->win);CHKERRQ(ierr);
360   ierr = PetscFree(link);CHKERRQ(ierr);
361   *win = MPI_WIN_NULL;
362   PetscFunctionReturn(0);
363 }
364 
365 #undef __FUNCT__
366 #define __FUNCT__ "PetscSFSetUp_Window"
367 static PetscErrorCode PetscSFSetUp_Window(PetscSF sf)
368 {
369   PetscSF_Window *w = (PetscSF_Window*)sf->data;
370   PetscErrorCode ierr;
371   MPI_Group      ingroup,outgroup;
372 
373   PetscFunctionBegin;
374   switch (w->sync) {
375   case PETSCSF_WINDOW_SYNC_ACTIVE:
376     ierr = PetscSFGetGroups(sf,&ingroup,&outgroup);CHKERRQ(ierr);
377   default:
378     break;
379   }
380   PetscFunctionReturn(0);
381 }
382 
383 #undef __FUNCT__
384 #define __FUNCT__ "PetscSFSetFromOptions_Window"
385 static PetscErrorCode PetscSFSetFromOptions_Window(PetscSF sf)
386 {
387   PetscSF_Window *w = (PetscSF_Window*)sf->data;
388   PetscErrorCode ierr;
389 
390   PetscFunctionBegin;
391   ierr = PetscOptionsHead("PetscSF Window options");CHKERRQ(ierr);
392   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   ierr = PetscOptionsTail();CHKERRQ(ierr);
394   PetscFunctionReturn(0);
395 }
396 
397 #undef __FUNCT__
398 #define __FUNCT__ "PetscSFReset_Window"
399 static PetscErrorCode PetscSFReset_Window(PetscSF sf)
400 {
401   PetscSF_Window  *w = (PetscSF_Window*)sf->data;
402   PetscErrorCode  ierr;
403   PetscSFDataLink link,next;
404   PetscSFWinLink  wlink,wnext;
405   PetscInt        i;
406 
407   PetscFunctionBegin;
408   for (link=w->link; link; link=next) {
409     next = link->next;
410     ierr = MPI_Type_free(&link->unit);CHKERRQ(ierr);
411     for (i=0; i<sf->nranks; i++) {
412       ierr = MPI_Type_free(&link->mine[i]);CHKERRQ(ierr);
413       ierr = MPI_Type_free(&link->remote[i]);CHKERRQ(ierr);
414     }
415     ierr = PetscFree2(link->mine,link->remote);CHKERRQ(ierr);
416     ierr = PetscFree(link);CHKERRQ(ierr);
417   }
418   w->link = NULL;
419   for (wlink=w->wins; wlink; wlink=wnext) {
420     wnext = wlink->next;
421     if (wlink->inuse) SETERRQ1(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Window still in use with address %p",(void*)wlink->addr);
422     ierr = MPI_Win_free(&wlink->win);CHKERRQ(ierr);
423     ierr = PetscFree(wlink);CHKERRQ(ierr);
424   }
425   w->wins = NULL;
426   PetscFunctionReturn(0);
427 }
428 
429 #undef __FUNCT__
430 #define __FUNCT__ "PetscSFDestroy_Window"
431 static PetscErrorCode PetscSFDestroy_Window(PetscSF sf)
432 {
433   PetscErrorCode ierr;
434 
435   PetscFunctionBegin;
436   ierr = PetscSFReset_Window(sf);CHKERRQ(ierr);
437   ierr = PetscFree(sf->data);CHKERRQ(ierr);
438   ierr = PetscObjectComposeFunction((PetscObject)sf,"PetscSFWindowSetSyncType_C",NULL);CHKERRQ(ierr);
439   ierr = PetscObjectComposeFunction((PetscObject)sf,"PetscSFWindowGetSyncType_C",NULL);CHKERRQ(ierr);
440   PetscFunctionReturn(0);
441 }
442 
443 #undef __FUNCT__
444 #define __FUNCT__ "PetscSFView_Window"
445 static PetscErrorCode PetscSFView_Window(PetscSF sf,PetscViewer viewer)
446 {
447   PetscSF_Window *w = (PetscSF_Window*)sf->data;
448   PetscErrorCode ierr;
449   PetscBool      iascii;
450 
451   PetscFunctionBegin;
452   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
453   if (iascii) {
454     ierr = PetscViewerASCIIPrintf(viewer,"  synchronization=%s sort=%s\n",PetscSFWindowSyncTypes[w->sync],sf->rankorder ? "rank-order" : "unordered");CHKERRQ(ierr);
455   }
456   PetscFunctionReturn(0);
457 }
458 
459 #undef __FUNCT__
460 #define __FUNCT__ "PetscSFDuplicate_Window"
461 static PetscErrorCode PetscSFDuplicate_Window(PetscSF sf,PetscSFDuplicateOption opt,PetscSF newsf)
462 {
463   PetscSF_Window        *w = (PetscSF_Window*)sf->data;
464   PetscErrorCode        ierr;
465   PetscSFWindowSyncType synctype;
466 
467   PetscFunctionBegin;
468   synctype = w->sync;
469   if (!sf->setupcalled) {
470     /* HACK: Must use FENCE or LOCK when called from PetscSFGetGroups() because ACTIVE here would cause recursion. */
471     synctype = PETSCSF_WINDOW_SYNC_LOCK;
472   }
473   ierr = PetscSFWindowSetSyncType(newsf,synctype);CHKERRQ(ierr);
474   PetscFunctionReturn(0);
475 }
476 
477 #undef __FUNCT__
478 #define __FUNCT__ "PetscSFBcastBegin_Window"
479 static PetscErrorCode PetscSFBcastBegin_Window(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
480 {
481   PetscSF_Window     *w = (PetscSF_Window*)sf->data;
482   PetscErrorCode     ierr;
483   PetscInt           i,nranks;
484   const PetscMPIInt  *ranks;
485   const MPI_Datatype *mine,*remote;
486   MPI_Win            win;
487 
488   PetscFunctionBegin;
489   ierr = PetscSFGetRanks(sf,&nranks,&ranks,NULL,NULL,NULL);CHKERRQ(ierr);
490   ierr = PetscSFWindowGetDataTypes(sf,unit,&mine,&remote);CHKERRQ(ierr);
491   ierr = PetscSFGetWindow(sf,unit,(void*)rootdata,PETSC_TRUE,MPI_MODE_NOPUT|MPI_MODE_NOPRECEDE,MPI_MODE_NOPUT,0,&win);CHKERRQ(ierr);
492   for (i=0; i<nranks; i++) {
493     if (w->sync == PETSCSF_WINDOW_SYNC_LOCK) {ierr = MPI_Win_lock(MPI_LOCK_SHARED,ranks[i],MPI_MODE_NOCHECK,win);CHKERRQ(ierr);}
494     ierr = MPI_Get(leafdata,1,mine[i],ranks[i],0,1,remote[i],win);CHKERRQ(ierr);
495     if (w->sync == PETSCSF_WINDOW_SYNC_LOCK) {ierr = MPI_Win_unlock(ranks[i],win);CHKERRQ(ierr);}
496   }
497   PetscFunctionReturn(0);
498 }
499 
500 #undef __FUNCT__
501 #define __FUNCT__ "PetscSFBcastEnd_Window"
502 PetscErrorCode PetscSFBcastEnd_Window(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
503 {
504   PetscErrorCode ierr;
505   MPI_Win        win;
506 
507   PetscFunctionBegin;
508   ierr = PetscSFFindWindow(sf,unit,rootdata,&win);CHKERRQ(ierr);
509   ierr = PetscSFRestoreWindow(sf,unit,rootdata,PETSC_TRUE,MPI_MODE_NOSTORE|MPI_MODE_NOSUCCEED,&win);CHKERRQ(ierr);
510   PetscFunctionReturn(0);
511 }
512 
513 #undef __FUNCT__
514 #define __FUNCT__ "PetscSFReduceBegin_Window"
515 PetscErrorCode PetscSFReduceBegin_Window(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
516 {
517   PetscSF_Window     *w = (PetscSF_Window*)sf->data;
518   PetscErrorCode     ierr;
519   PetscInt           i,nranks;
520   const PetscMPIInt  *ranks;
521   const MPI_Datatype *mine,*remote;
522   MPI_Win            win;
523 
524   PetscFunctionBegin;
525   ierr = PetscSFGetRanks(sf,&nranks,&ranks,NULL,NULL,NULL);CHKERRQ(ierr);
526   ierr = PetscSFWindowGetDataTypes(sf,unit,&mine,&remote);CHKERRQ(ierr);
527   ierr = PetscSFWindowOpTranslate(&op);CHKERRQ(ierr);
528   ierr = PetscSFGetWindow(sf,unit,rootdata,PETSC_TRUE,MPI_MODE_NOPRECEDE,0,0,&win);CHKERRQ(ierr);
529   for (i=0; i<nranks; i++) {
530     if (w->sync == PETSCSF_WINDOW_SYNC_LOCK) {ierr = MPI_Win_lock(MPI_LOCK_SHARED,ranks[i],MPI_MODE_NOCHECK,win);CHKERRQ(ierr);}
531     ierr = MPI_Accumulate((void*)leafdata,1,mine[i],ranks[i],0,1,remote[i],op,win);CHKERRQ(ierr);
532     if (w->sync == PETSCSF_WINDOW_SYNC_LOCK) {ierr = MPI_Win_unlock(ranks[i],win);CHKERRQ(ierr);}
533   }
534   PetscFunctionReturn(0);
535 }
536 
537 #undef __FUNCT__
538 #define __FUNCT__ "PetscSFReduceEnd_Window"
539 static PetscErrorCode PetscSFReduceEnd_Window(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
540 {
541   PetscSF_Window *w = (PetscSF_Window*)sf->data;
542   PetscErrorCode ierr;
543   MPI_Win        win;
544 
545   PetscFunctionBegin;
546   if (!w->wins) PetscFunctionReturn(0);
547   ierr = PetscSFFindWindow(sf,unit,rootdata,&win);CHKERRQ(ierr);
548   ierr = MPI_Win_fence(MPI_MODE_NOSUCCEED,win);CHKERRQ(ierr);
549   ierr = PetscSFRestoreWindow(sf,unit,rootdata,PETSC_TRUE,MPI_MODE_NOSUCCEED,&win);CHKERRQ(ierr);
550   PetscFunctionReturn(0);
551 }
552 #undef __FUNCT__
553 #define __FUNCT__ "PetscSFFetchAndOpBegin_Window"
554 static PetscErrorCode PetscSFFetchAndOpBegin_Window(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
555 {
556   PetscErrorCode     ierr;
557   PetscInt           i,nranks;
558   const PetscMPIInt  *ranks;
559   const MPI_Datatype *mine,*remote;
560   MPI_Win            win;
561 
562   PetscFunctionBegin;
563   ierr = PetscSFGetRanks(sf,&nranks,&ranks,NULL,NULL,NULL);CHKERRQ(ierr);
564   ierr = PetscSFWindowGetDataTypes(sf,unit,&mine,&remote);CHKERRQ(ierr);
565   ierr = PetscSFWindowOpTranslate(&op);CHKERRQ(ierr);
566   ierr = PetscSFGetWindow(sf,unit,rootdata,PETSC_FALSE,0,0,0,&win);CHKERRQ(ierr);
567   for (i=0; i<sf->nranks; i++) {
568     ierr = MPI_Win_lock(MPI_LOCK_EXCLUSIVE,sf->ranks[i],0,win);CHKERRQ(ierr);
569     ierr = MPI_Get(leafupdate,1,mine[i],ranks[i],0,1,remote[i],win);CHKERRQ(ierr);
570     ierr = MPI_Accumulate((void*)leafdata,1,mine[i],ranks[i],0,1,remote[i],op,win);CHKERRQ(ierr);
571     ierr = MPI_Win_unlock(ranks[i],win);CHKERRQ(ierr);
572   }
573   PetscFunctionReturn(0);
574 }
575 
576 #undef __FUNCT__
577 #define __FUNCT__ "PetscSFFetchAndOpEnd_Window"
578 static PetscErrorCode PetscSFFetchAndOpEnd_Window(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
579 {
580   PetscErrorCode ierr;
581   MPI_Win        win;
582 
583   PetscFunctionBegin;
584   ierr = PetscSFFindWindow(sf,unit,rootdata,&win);CHKERRQ(ierr);
585   /* Nothing to do currently because MPI_LOCK_EXCLUSIVE is used in PetscSFFetchAndOpBegin(), rendering this implementation synchronous. */
586   ierr = PetscSFRestoreWindow(sf,unit,rootdata,PETSC_FALSE,0,&win);CHKERRQ(ierr);
587   PetscFunctionReturn(0);
588 }
589 
590 #undef __FUNCT__
591 #define __FUNCT__ "PetscSFCreate_Window"
592 PETSC_EXTERN PetscErrorCode PetscSFCreate_Window(PetscSF sf)
593 {
594   PetscSF_Window *w = (PetscSF_Window*)sf->data;
595   PetscErrorCode ierr;
596 
597   PetscFunctionBegin;
598   sf->ops->SetUp           = PetscSFSetUp_Window;
599   sf->ops->SetFromOptions  = PetscSFSetFromOptions_Window;
600   sf->ops->Reset           = PetscSFReset_Window;
601   sf->ops->Destroy         = PetscSFDestroy_Window;
602   sf->ops->View            = PetscSFView_Window;
603   sf->ops->Duplicate       = PetscSFDuplicate_Window;
604   sf->ops->BcastBegin      = PetscSFBcastBegin_Window;
605   sf->ops->BcastEnd        = PetscSFBcastEnd_Window;
606   sf->ops->ReduceBegin     = PetscSFReduceBegin_Window;
607   sf->ops->ReduceEnd       = PetscSFReduceEnd_Window;
608   sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Window;
609   sf->ops->FetchAndOpEnd   = PetscSFFetchAndOpEnd_Window;
610 
611   ierr     = PetscNewLog(sf,PetscSF_Window,&w);CHKERRQ(ierr);
612   sf->data = (void*)w;
613   w->sync  = PETSCSF_WINDOW_SYNC_FENCE;
614 
615   ierr = PetscObjectComposeFunction((PetscObject)sf,"PetscSFWindowSetSyncType_C",PetscSFWindowSetSyncType_Window);CHKERRQ(ierr);
616   ierr = PetscObjectComposeFunction((PetscObject)sf,"PetscSFWindowGetSyncType_C",PetscSFWindowGetSyncType_Window);CHKERRQ(ierr);
617 
618 #if defined(OMPI_MAJOR_VERSION) && (OMPI_MAJOR_VERSION < 1 || (OMPI_MAJOR_VERSION == 1 && OMPI_MINOR_VERSION <= 6))
619   {
620     PetscBool ackbug = PETSC_FALSE;
621     ierr = PetscOptionsGetBool(NULL,"-acknowledge_ompi_onesided_bug",&ackbug,NULL);CHKERRQ(ierr);
622     if (ackbug) {
623       ierr = PetscInfo(sf,"Acknowledged Open MPI bug, proceeding anyway. Expect memory corruption.");CHKERRQ(ierr);
624     } 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   }
626 #endif
627   PetscFunctionReturn(0);
628 }
629