xref: /petsc/src/vec/is/sf/impls/basic/sfbasic.c (revision 992568525dc78450d29be7cf7c589d143c34945b) !
1 
2 #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/
3 
4 typedef struct _n_PetscSFBasicPack *PetscSFBasicPack;
5 struct _n_PetscSFBasicPack {
6   void (*Pack)(PetscInt,PetscInt,const PetscInt*,const void*,void*);
7   void (*UnpackInsert)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
8   void (*UnpackAdd)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
9   void (*UnpackMin)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
10   void (*UnpackMax)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
11   void (*UnpackMinloc)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
12   void (*UnpackMaxloc)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
13   void (*UnpackMult)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
14   void (*UnpackLAND)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
15   void (*UnpackBAND)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
16   void (*UnpackLOR)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
17   void (*UnpackBOR)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
18   void (*UnpackLXOR)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
19   void (*UnpackBXOR)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
20   void (*FetchAndInsert)(PetscInt,PetscInt,const PetscInt*,void*,void*);
21   void (*FetchAndAdd)(PetscInt,PetscInt,const PetscInt*,void*,void*);
22   void (*FetchAndMin)(PetscInt,PetscInt,const PetscInt*,void*,void*);
23   void (*FetchAndMax)(PetscInt,PetscInt,const PetscInt*,void*,void*);
24   void (*FetchAndMinloc)(PetscInt,PetscInt,const PetscInt*,void*,void*);
25   void (*FetchAndMaxloc)(PetscInt,PetscInt,const PetscInt*,void*,void*);
26   void (*FetchAndMult)(PetscInt,PetscInt,const PetscInt*,void*,void*);
27   void (*FetchAndLAND)(PetscInt,PetscInt,const PetscInt*,void*,void*);
28   void (*FetchAndBAND)(PetscInt,PetscInt,const PetscInt*,void*,void*);
29   void (*FetchAndLOR)(PetscInt,PetscInt,const PetscInt*,void*,void*);
30   void (*FetchAndBOR)(PetscInt,PetscInt,const PetscInt*,void*,void*);
31   void (*FetchAndLXOR)(PetscInt,PetscInt,const PetscInt*,void*,void*);
32   void (*FetchAndBXOR)(PetscInt,PetscInt,const PetscInt*,void*,void*);
33 
34   MPI_Datatype     unit;
35   size_t           unitbytes;   /* Number of bytes in a unit */
36   PetscInt         bs;          /* Number of basic units in a unit */
37   const void       *key;        /* Array used as key for operation */
38   char             **root;      /* Packed root data, indexed by leaf rank */
39   char             **leaf;      /* Packed leaf data, indexed by root rank */
40   MPI_Request      *requests;   /* Array of root requests followed by leaf requests */
41   PetscSFBasicPack next;
42 };
43 
44 typedef struct {
45   PetscMPIInt      tag;
46   PetscMPIInt      niranks;     /* Number of incoming ranks (ranks accessing my roots) */
47   PetscMPIInt      ndiranks;    /* Number of incoming ranks (ranks accessing my roots) in distinguished set */
48   PetscMPIInt      *iranks;     /* Array of ranks that reference my roots */
49   PetscInt         itotal;      /* Total number of graph edges referencing my roots */
50   PetscInt         *ioffset;    /* Array of length niranks+1 holding offset in irootloc[] for each rank */
51   PetscInt         *irootloc;   /* Incoming roots referenced by ranks starting at ioffset[rank] */
52   PetscSFBasicPack avail;       /* One or more entries per MPI Datatype, lazily constructed */
53   PetscSFBasicPack inuse;       /* Buffers being used for transactions that have not yet completed */
54 } PetscSF_Basic;
55 
56 #if !defined(PETSC_HAVE_MPI_TYPE_DUP) && !defined(PETSC_HAVE_MPIUNI)
57 /* Danger: type is not reference counted; subject to ABA problem */
58 PETSC_STATIC_INLINE PetscErrorCode MPI_Type_dup(MPI_Datatype datatype,MPI_Datatype *newtype)
59 {
60   *newtype = datatype;
61   return 0;
62 }
63 #endif
64 
65 /*
66  * MPI_Reduce_local is not really useful because it can't handle sparse data and it vectorizes "in the wrong direction",
67  * therefore we pack data types manually. This section defines packing routines for the standard data types.
68  */
69 
70 #define CPPJoin2_exp(a,b) a ## b
71 #define CPPJoin2(a,b) CPPJoin2_exp(a,b)
72 #define CPPJoin3_exp_(a,b,c) a ## b ## _ ## c
73 #define CPPJoin3_(a,b,c) CPPJoin3_exp_(a,b,c)
74 
75 /* Basic types without addition */
76 #define DEF_PackNoInit(type,BS)                                         \
77   static void CPPJoin3_(Pack_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,const void *unpacked,void *packed) { \
78     const type *u = (const type*)unpacked;                              \
79     type *p = (type*)packed;                                            \
80     PetscInt i,j,k;                                                     \
81     for (i=0; i<n; i++)                                                 \
82       for (j=0; j<bs; j+=BS)                                            \
83         for (k=j; k<j+BS; k++)                                          \
84           p[i*bs+k] = u[idx[i]*bs+k];                                   \
85   }                                                                     \
86   static void CPPJoin3_(UnpackInsert_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
87     type *u = (type*)unpacked;                                          \
88     const type *p = (const type*)packed;                                \
89     PetscInt i,j,k;                                                     \
90     for (i=0; i<n; i++)                                                 \
91       for (j=0; j<bs; j+=BS)                                            \
92         for (k=j; k<j+BS; k++)                                          \
93           u[idx[i]*bs+k] = p[i*bs+k];                                   \
94   }                                                                     \
95   static void CPPJoin3_(FetchAndInsert_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
96     type *u = (type*)unpacked;                                          \
97     type *p = (type*)packed;                                            \
98     PetscInt i,j,k;                                                     \
99     for (i=0; i<n; i++) {                                               \
100       PetscInt ii = idx[i];                                             \
101       for (j=0; j<bs; j+=BS)                                            \
102         for (k=j; k<j+BS; k++) {                                        \
103           type t = u[ii*bs+k];                                          \
104           u[ii*bs+k] = p[i*bs+k];                                       \
105           p[i*bs+k] = t;                                                \
106         }                                                               \
107     }                                                                   \
108   }
109 
110 /* Basic types defining addition */
111 #define DEF_PackAddNoInit(type,BS)                                      \
112   DEF_PackNoInit(type,BS)                                               \
113   static void CPPJoin3_(UnpackAdd_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
114     type *u = (type*)unpacked;                                          \
115     const type *p = (const type*)packed;                                \
116     PetscInt i,j,k;                                                     \
117     for (i=0; i<n; i++)                                                 \
118       for (j=0; j<bs; j+=BS)                                            \
119         for (k=j; k<j+BS; k++)                                          \
120           u[idx[i]*bs+k] += p[i*bs+k];                                  \
121   }                                                                     \
122   static void CPPJoin3_(FetchAndAdd_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
123     type *u = (type*)unpacked;                                          \
124     type *p = (type*)packed;                                            \
125     PetscInt i,j,k;                                                     \
126     for (i=0; i<n; i++) {                                               \
127       PetscInt ii = idx[i];                                             \
128       for (j=0; j<bs; j+=BS)                                            \
129         for (k=j; k<j+BS; k++) {                                        \
130           type t = u[ii*bs+k];                                          \
131           u[ii*bs+k] = t + p[i*bs+k];                                   \
132           p[i*bs+k] = t;                                                \
133         }                                                               \
134     }                                                                   \
135   }                                                                     \
136   static void CPPJoin3_(UnpackMult_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
137     type *u = (type*)unpacked;                                          \
138     const type *p = (const type*)packed;                                \
139     PetscInt i,j,k;                                                     \
140     for (i=0; i<n; i++)                                                 \
141       for (j=0; j<bs; j+=BS)                                            \
142         for (k=j; k<j+BS; k++)                                          \
143           u[idx[i]*bs+k] *= p[i*bs+k];                                  \
144   }                                                                     \
145   static void CPPJoin3_(FetchAndMult_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
146     type *u = (type*)unpacked;                                          \
147     type *p = (type*)packed;                                            \
148     PetscInt i,j,k;                                                     \
149     for (i=0; i<n; i++) {                                               \
150       PetscInt ii = idx[i];                                             \
151       for (j=0; j<bs; j+=BS)                                            \
152         for (k=j; k<j+BS; k++) {                                        \
153           type t = u[ii*bs+k];                                          \
154           u[ii*bs+k] = t * p[i*bs+k];                                   \
155           p[i*bs+k] = t;                                                \
156         }                                                               \
157     }                                                                   \
158   }
159 #define DEF_Pack(type,BS)                                               \
160   DEF_PackAddNoInit(type,BS)                                            \
161   static void CPPJoin3_(PackInit_,type,BS)(PetscSFBasicPack link) {     \
162     link->Pack = CPPJoin3_(Pack_,type,BS);                              \
163     link->UnpackInsert = CPPJoin3_(UnpackInsert_,type,BS);              \
164     link->UnpackAdd = CPPJoin3_(UnpackAdd_,type,BS);                    \
165     link->UnpackMult = CPPJoin3_(UnpackMult_,type,BS);                  \
166     link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,type,BS);          \
167     link->FetchAndAdd = CPPJoin3_(FetchAndAdd_,type,BS);                \
168     link->FetchAndMult = CPPJoin3_(FetchAndMult_,type,BS);              \
169     link->unitbytes = sizeof(type);                                     \
170   }
171 /* Comparable types */
172 #define DEF_PackCmp(type)                                               \
173   DEF_PackAddNoInit(type,1)                                             \
174   static void CPPJoin2(UnpackMax_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
175     type *u = (type*)unpacked;                                          \
176     const type *p = (const type*)packed;                                \
177     PetscInt i;                                                         \
178     for (i=0; i<n; i++) {                                               \
179       type v = u[idx[i]];                                               \
180       u[idx[i]] = PetscMax(v,p[i]);                                     \
181     }                                                                   \
182   }                                                                     \
183   static void CPPJoin2(UnpackMin_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
184     type *u = (type*)unpacked;                                          \
185     const type *p = (const type*)packed;                                \
186     PetscInt i;                                                         \
187     for (i=0; i<n; i++) {                                               \
188       type v = u[idx[i]];                                               \
189       u[idx[i]] = PetscMin(v,p[i]);                                     \
190     }                                                                   \
191   }                                                                     \
192   static void CPPJoin2(FetchAndMax_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
193     type *u = (type*)unpacked;                                          \
194     type *p = (type*)packed;                                            \
195     PetscInt i;                                                         \
196     for (i=0; i<n; i++) {                                               \
197       PetscInt j = idx[i];                                              \
198       type v = u[j];                                                    \
199       u[j] = PetscMax(v,p[i]);                                          \
200       p[i] = v;                                                         \
201     }                                                                   \
202   }                                                                     \
203   static void CPPJoin2(FetchAndMin_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
204     type *u = (type*)unpacked;                                          \
205     type *p = (type*)packed;                                            \
206     PetscInt i;                                                         \
207     for (i=0; i<n; i++) {                                               \
208       PetscInt j = idx[i];                                              \
209       type v = u[j];                                                    \
210       u[j] = PetscMin(v,p[i]);                                          \
211       p[i] = v;                                                         \
212     }                                                                   \
213   }                                                                     \
214   static void CPPJoin2(PackInit_,type)(PetscSFBasicPack link) {         \
215     link->Pack = CPPJoin3_(Pack_,type,1);                               \
216     link->UnpackInsert = CPPJoin3_(UnpackInsert_,type,1);               \
217     link->UnpackAdd  = CPPJoin3_(UnpackAdd_,type,1);                    \
218     link->UnpackMax  = CPPJoin2(UnpackMax_,type);                       \
219     link->UnpackMin  = CPPJoin2(UnpackMin_,type);                       \
220     link->UnpackMult = CPPJoin3_(UnpackMult_,type,1);                   \
221     link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,type,1);           \
222     link->FetchAndAdd = CPPJoin3_(FetchAndAdd_ ,type,1);                \
223     link->FetchAndMax = CPPJoin2(FetchAndMax_ ,type);                   \
224     link->FetchAndMin = CPPJoin2(FetchAndMin_ ,type);                   \
225     link->FetchAndMult = CPPJoin3_(FetchAndMult_,type,1);               \
226     link->unitbytes = sizeof(type);                                     \
227   }
228 
229 /* Logical Types */
230 #define DEF_PackLog(type)                                               \
231   static void CPPJoin2(UnpackLAND_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
232     type *u = (type*)unpacked;                                          \
233     const type *p = (const type*)packed;                                \
234     PetscInt i;                                                         \
235     for (i=0; i<n; i++) {                                               \
236       type v = u[idx[i]];                                               \
237       u[idx[i]] = v && p[i];                                            \
238     }                                                                   \
239   }                                                                     \
240   static void CPPJoin2(UnpackLOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
241     type *u = (type*)unpacked;                                          \
242     const type *p = (const type*)packed;                                \
243     PetscInt i;                                                         \
244     for (i=0; i<n; i++) {                                               \
245       type v = u[idx[i]];                                               \
246       u[idx[i]] = v || p[i];                                            \
247     }                                                                   \
248   }                                                                     \
249   static void CPPJoin2(UnpackLXOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
250     type *u = (type*)unpacked;                                          \
251     const type *p = (const type*)packed;                                \
252     PetscInt i;                                                         \
253     for (i=0; i<n; i++) {                                               \
254       type v = u[idx[i]];                                               \
255       u[idx[i]] = (!v)!=(!p[i]);                                        \
256     }                                                                   \
257   }                                                                     \
258   static void CPPJoin2(FetchAndLAND_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
259     type *u = (type*)unpacked;                                          \
260     type *p = (type*)packed;                                            \
261     PetscInt i;                                                         \
262     for (i=0; i<n; i++) {                                               \
263       PetscInt j = idx[i];                                              \
264       type v = u[j];                                                    \
265       u[j] = v && p[i];                                                 \
266       p[i] = v;                                                         \
267     }                                                                   \
268   }                                                                     \
269   static void CPPJoin2(FetchAndLOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
270     type *u = (type*)unpacked;                                          \
271     type *p = (type*)packed;                                            \
272     PetscInt i;                                                         \
273     for (i=0; i<n; i++) {                                               \
274       PetscInt j = idx[i];                                              \
275       type v = u[j];                                                    \
276       u[j] = v || p[i];                                                 \
277       p[i] = v;                                                         \
278     }                                                                   \
279   }                                                                     \
280   static void CPPJoin2(FetchAndLXOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
281     type *u = (type*)unpacked;                                          \
282     type *p = (type*)packed;                                            \
283     PetscInt i;                                                         \
284     for (i=0; i<n; i++) {                                               \
285       PetscInt j = idx[i];                                              \
286       type v = u[j];                                                    \
287       u[j] = (!v)!=(!p[i]);                                             \
288       p[i] = v;                                                         \
289     }                                                                   \
290   }                                                                     \
291   static void CPPJoin2(PackInit_Logical_,type)(PetscSFBasicPack link) { \
292     link->UnpackLAND = CPPJoin2(UnpackLAND_,type);                      \
293     link->UnpackLOR  = CPPJoin2(UnpackLOR_,type);                       \
294     link->UnpackLXOR = CPPJoin2(UnpackLXOR_,type);                      \
295     link->FetchAndLAND = CPPJoin2(FetchAndLAND_,type);                  \
296     link->FetchAndLOR  = CPPJoin2(FetchAndLOR_,type);                   \
297     link->FetchAndLXOR = CPPJoin2(FetchAndLXOR_,type);                  \
298   }
299 
300 
301 /* Bitwise Types */
302 #define DEF_PackBit(type)                                               \
303   static void CPPJoin2(UnpackBAND_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
304     type *u = (type*)unpacked;                                          \
305     const type *p = (const type*)packed;                                \
306     PetscInt i;                                                         \
307     for (i=0; i<n; i++) {                                               \
308       type v = u[idx[i]];                                               \
309       u[idx[i]] = v & p[i];                                             \
310     }                                                                   \
311   }                                                                     \
312   static void CPPJoin2(UnpackBOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
313     type *u = (type*)unpacked;                                          \
314     const type *p = (const type*)packed;                                \
315     PetscInt i;                                                         \
316     for (i=0; i<n; i++) {                                               \
317       type v = u[idx[i]];                                               \
318       u[idx[i]] = v | p[i];                                             \
319     }                                                                   \
320   }                                                                     \
321   static void CPPJoin2(UnpackBXOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
322     type *u = (type*)unpacked;                                          \
323     const type *p = (const type*)packed;                                \
324     PetscInt i;                                                         \
325     for (i=0; i<n; i++) {                                               \
326       type v = u[idx[i]];                                               \
327       u[idx[i]] = v^p[i];                                               \
328     }                                                                   \
329   }                                                                     \
330   static void CPPJoin2(FetchAndBAND_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
331     type *u = (type*)unpacked;                                          \
332     type *p = (type*)packed;                                            \
333     PetscInt i;                                                         \
334     for (i=0; i<n; i++) {                                               \
335       PetscInt j = idx[i];                                              \
336       type v = u[j];                                                    \
337       u[j] = v & p[i];                                                  \
338       p[i] = v;                                                         \
339     }                                                                   \
340   }                                                                     \
341   static void CPPJoin2(FetchAndBOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
342     type *u = (type*)unpacked;                                          \
343     type *p = (type*)packed;                                            \
344     PetscInt i;                                                         \
345     for (i=0; i<n; i++) {                                               \
346       PetscInt j = idx[i];                                              \
347       type v = u[j];                                                    \
348       u[j] = v | p[i];                                                  \
349       p[i] = v;                                                         \
350     }                                                                   \
351   }                                                                     \
352   static void CPPJoin2(FetchAndBXOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
353     type *u = (type*)unpacked;                                          \
354     type *p = (type*)packed;                                            \
355     PetscInt i;                                                         \
356     for (i=0; i<n; i++) {                                               \
357       PetscInt j = idx[i];                                              \
358       type v = u[j];                                                    \
359       u[j] = v^p[i];                                                    \
360       p[i] = v;                                                         \
361     }                                                                   \
362   }                                                                     \
363   static void CPPJoin2(PackInit_Bitwise_,type)(PetscSFBasicPack link) { \
364     link->UnpackBAND = CPPJoin2(UnpackBAND_,type);                      \
365     link->UnpackBOR  = CPPJoin2(UnpackBOR_,type);                       \
366     link->UnpackBXOR = CPPJoin2(UnpackBXOR_,type);                      \
367     link->FetchAndBAND = CPPJoin2(FetchAndBAND_,type);                  \
368     link->FetchAndBOR  = CPPJoin2(FetchAndBOR_,type);                   \
369     link->FetchAndBXOR = CPPJoin2(FetchAndBXOR_,type);                  \
370   }
371 
372 /* Pair types */
373 #define CPPJoinloc_exp(base,op,t1,t2) base ## op ## loc_ ## t1 ## _ ## t2
374 #define CPPJoinloc(base,op,t1,t2) CPPJoinloc_exp(base,op,t1,t2)
375 #define PairType(type1,type2) CPPJoin3_(_pairtype_,type1,type2)
376 #define DEF_UnpackXloc(type1,type2,locname,op)                              \
377   static void CPPJoinloc(Unpack,locname,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
378     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;        \
379     const PairType(type1,type2) *p = (const PairType(type1,type2)*)packed; \
380     PetscInt i;                                                         \
381     for (i=0; i<n; i++) {                                               \
382       PetscInt j = idx[i];                                              \
383       if (p[i].a op u[j].a) {                                           \
384         u[j].a = p[i].a;                                                \
385         u[j].b = p[i].b;                                                \
386       } else if (u[j].a == p[i].a) {                                    \
387         u[j].b = PetscMin(u[j].b,p[i].b);                               \
388       }                                                                 \
389     }                                                                   \
390   }                                                                     \
391   static void CPPJoinloc(FetchAnd,locname,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
392     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;        \
393     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;          \
394     PetscInt i;                                                         \
395     for (i=0; i<n; i++) {                                               \
396       PetscInt j = idx[i];                                              \
397       PairType(type1,type2) v;                                          \
398       v.a = u[j].a;                                                     \
399       v.b = u[j].b;                                                     \
400       if (p[i].a op u[j].a) {                                           \
401         u[j].a = p[i].a;                                                \
402         u[j].b = p[i].b;                                                \
403       } else if (u[j].a == p[i].a) {                                    \
404         u[j].b = PetscMin(u[j].b,p[i].b);                               \
405       }                                                                 \
406       p[i].a = v.a;                                                     \
407       p[i].b = v.b;                                                     \
408     }                                                                   \
409   }
410 #define DEF_PackPair(type1,type2)                                       \
411   typedef struct {type1 a; type2 b;} PairType(type1,type2);             \
412   static void CPPJoin3_(Pack_,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,const void *unpacked,void *packed) { \
413     const PairType(type1,type2) *u = (const PairType(type1,type2)*)unpacked; \
414     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;          \
415     PetscInt i;                                                         \
416     for (i=0; i<n; i++) {                                               \
417       p[i].a = u[idx[i]].a;                                             \
418       p[i].b = u[idx[i]].b;                                             \
419     }                                                                   \
420   }                                                                     \
421   static void CPPJoin3_(UnpackInsert_,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
422     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;       \
423     const PairType(type1,type2) *p = (const PairType(type1,type2)*)packed; \
424     PetscInt i;                                                         \
425     for (i=0; i<n; i++) {                                               \
426       u[idx[i]].a = p[i].a;                                             \
427       u[idx[i]].b = p[i].b;                                             \
428     }                                                                   \
429   }                                                                     \
430   static void CPPJoin3_(UnpackAdd_,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
431     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;       \
432     const PairType(type1,type2) *p = (const PairType(type1,type2)*)packed; \
433     PetscInt i;                                                         \
434     for (i=0; i<n; i++) {                                               \
435       u[idx[i]].a += p[i].a;                                            \
436       u[idx[i]].b += p[i].b;                                            \
437     }                                                                   \
438   }                                                                     \
439   static void CPPJoin3_(FetchAndInsert_,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
440     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;        \
441     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;          \
442     PetscInt i;                                                         \
443     for (i=0; i<n; i++) {                                               \
444       PetscInt j = idx[i];                                              \
445       PairType(type1,type2) v;                                          \
446       v.a = u[j].a;                                                     \
447       v.b = u[j].b;                                                     \
448       u[j].a = p[i].a;                                                  \
449       u[j].b = p[i].b;                                                  \
450       p[i].a = v.a;                                                     \
451       p[i].b = v.b;                                                     \
452     }                                                                   \
453   }                                                                     \
454   static void FetchAndAdd_ ## type1 ## _ ## type2(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
455     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;       \
456     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;         \
457     PetscInt i;                                                         \
458     for (i=0; i<n; i++) {                                               \
459       PetscInt j = idx[i];                                              \
460       PairType(type1,type2) v;                                          \
461       v.a = u[j].a;                                                     \
462       v.b = u[j].b;                                                     \
463       u[j].a = v.a + p[i].a;                                            \
464       u[j].b = v.b + p[i].b;                                            \
465       p[i].a = v.a;                                                     \
466       p[i].b = v.b;                                                     \
467     }                                                                   \
468   }                                                                     \
469   DEF_UnpackXloc(type1,type2,Max,>)                                     \
470   DEF_UnpackXloc(type1,type2,Min,<)                                     \
471   static void CPPJoin3_(PackInit_,type1,type2)(PetscSFBasicPack link) { \
472     link->Pack = CPPJoin3_(Pack_,type1,type2);                          \
473     link->UnpackInsert = CPPJoin3_(UnpackInsert_,type1,type2);          \
474     link->UnpackAdd = CPPJoin3_(UnpackAdd_,type1,type2);                \
475     link->UnpackMaxloc = CPPJoin3_(UnpackMaxloc_,type1,type2);          \
476     link->UnpackMinloc = CPPJoin3_(UnpackMinloc_,type1,type2);          \
477     link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,type1,type2);      \
478     link->FetchAndAdd = CPPJoin3_(FetchAndAdd_,type1,type2);            \
479     link->FetchAndMaxloc = CPPJoin3_(FetchAndMaxloc_,type1,type2);      \
480     link->FetchAndMinloc = CPPJoin3_(FetchAndMinloc_,type1,type2);      \
481     link->unitbytes = sizeof(PairType(type1,type2));                    \
482   }
483 
484 /* Currently only dumb blocks of data */
485 #define BlockType(unit,count) CPPJoin3_(_blocktype_,unit,count)
486 #define DEF_Block(unit,count)                                           \
487   typedef struct {unit v[count];} BlockType(unit,count);                \
488   DEF_PackNoInit(BlockType(unit,count),1)                               \
489   static void CPPJoin3_(PackInit_block_,unit,count)(PetscSFBasicPack link) { \
490     link->Pack = CPPJoin3_(Pack_,BlockType(unit,count),1);               \
491     link->UnpackInsert = CPPJoin3_(UnpackInsert_,BlockType(unit,count),1); \
492     link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,BlockType(unit,count),1); \
493     link->unitbytes = sizeof(BlockType(unit,count));                    \
494   }
495 
496 DEF_PackCmp(int)
497 DEF_PackBit(int)
498 DEF_PackLog(int)
499 DEF_PackCmp(PetscInt)
500 DEF_PackBit(PetscInt)
501 DEF_PackLog(PetscInt)
502 DEF_Pack(PetscInt,2)
503 DEF_Pack(PetscInt,3)
504 DEF_Pack(PetscInt,4)
505 DEF_Pack(PetscInt,5)
506 DEF_Pack(PetscInt,7)
507 DEF_PackCmp(PetscReal)
508 DEF_PackLog(PetscReal)
509 DEF_Pack(PetscReal,2)
510 DEF_Pack(PetscReal,3)
511 DEF_Pack(PetscReal,4)
512 DEF_Pack(PetscReal,5)
513 DEF_Pack(PetscReal,7)
514 #if defined(PETSC_HAVE_COMPLEX)
515 DEF_Pack(PetscComplex,1)
516 DEF_Pack(PetscComplex,2)
517 DEF_Pack(PetscComplex,3)
518 DEF_Pack(PetscComplex,4)
519 DEF_Pack(PetscComplex,5)
520 DEF_Pack(PetscComplex,7)
521 #endif
522 DEF_PackPair(int,int)
523 DEF_PackPair(PetscInt,PetscInt)
524 DEF_Block(int,1)
525 DEF_Block(int,2)
526 DEF_Block(int,3)
527 DEF_Block(int,4)
528 DEF_Block(int,5)
529 DEF_Block(int,6)
530 DEF_Block(int,7)
531 DEF_Block(int,8)
532 
533 static PetscErrorCode PetscSFSetUp_Basic(PetscSF sf)
534 {
535   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
536   PetscErrorCode ierr;
537   PetscInt *rlengths,*ilengths,i;
538   PetscMPIInt rank,niranks,*iranks;
539   MPI_Comm comm;
540   MPI_Group group;
541   MPI_Request *rootreqs,*leafreqs;
542 
543   PetscFunctionBegin;
544   ierr = MPI_Comm_group(PETSC_COMM_SELF,&group);CHKERRQ(ierr);
545   ierr = PetscSFSetUpRanks(sf,group);CHKERRQ(ierr);
546   ierr = MPI_Group_free(&group);CHKERRQ(ierr);
547   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
548   ierr = PetscObjectGetNewTag((PetscObject)sf,&bas->tag);CHKERRQ(ierr);
549   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
550   /*
551    * Inform roots about how many leaves and from which ranks
552    */
553   ierr = PetscMalloc1(sf->nranks,&rlengths);CHKERRQ(ierr);
554   /* Determine number, sending ranks, and length of incoming */
555   for (i=0; i<sf->nranks; i++) {
556     rlengths[i] = sf->roffset[i+1] - sf->roffset[i]; /* Number of roots referenced by my leaves; for rank sf->ranks[i] */
557   }
558   ierr = PetscCommBuildTwoSided(comm,1,MPIU_INT,sf->nranks-sf->ndranks,sf->ranks+sf->ndranks,rlengths+sf->ndranks,&niranks,&iranks,(void**)&ilengths);CHKERRQ(ierr);
559 
560   /* Partition into distinguished and non-distinguished incoming ranks */
561   bas->ndiranks = sf->ndranks;
562   bas->niranks = bas->ndiranks + niranks;
563   ierr = PetscMalloc2(bas->niranks,&bas->iranks,bas->niranks+1,&bas->ioffset);CHKERRQ(ierr);
564   bas->ioffset[0] = 0;
565   for (i=0; i<bas->ndiranks; i++) {
566     bas->iranks[i] = sf->ranks[i];
567     bas->ioffset[i+1] = bas->ioffset[i] + rlengths[i];
568   }
569   if (bas->ndiranks > 1 || (bas->ndiranks == 1 && bas->iranks[0] != rank)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Broken setup for shared ranks");
570   for ( ; i<bas->niranks; i++) {
571     bas->iranks[i] = iranks[i-bas->ndiranks];
572     bas->ioffset[i+1] = bas->ioffset[i] + ilengths[i-bas->ndiranks];
573   }
574   bas->itotal = bas->ioffset[i];
575   ierr = PetscFree(rlengths);CHKERRQ(ierr);
576   ierr = PetscFree(iranks);CHKERRQ(ierr);
577   ierr = PetscFree(ilengths);CHKERRQ(ierr);
578 
579   /* Send leaf identities to roots */
580   ierr = PetscMalloc1(bas->itotal,&bas->irootloc);CHKERRQ(ierr);
581   ierr = PetscMalloc2(bas->niranks-bas->ndiranks,&rootreqs,sf->nranks-sf->ndranks,&leafreqs);CHKERRQ(ierr);
582   for (i=bas->ndiranks; i<bas->niranks; i++) {
583     ierr = MPI_Irecv(bas->irootloc+bas->ioffset[i],bas->ioffset[i+1]-bas->ioffset[i],MPIU_INT,bas->iranks[i],bas->tag,comm,&rootreqs[i-bas->ndiranks]);CHKERRQ(ierr);
584   }
585   for (i=0; i<sf->nranks; i++) {
586     PetscMPIInt npoints;
587     ierr = PetscMPIIntCast(sf->roffset[i+1] - sf->roffset[i],&npoints);CHKERRQ(ierr);
588     if (i < sf->ndranks) {
589       if (sf->ranks[i] != rank) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot interpret distinguished leaf rank");
590       if (bas->iranks[0] != rank) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot interpret distinguished root rank");
591       if (npoints != bas->ioffset[1]-bas->ioffset[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Distinguished rank exchange has mismatched lengths");
592       ierr = PetscMemcpy(bas->irootloc+bas->ioffset[0],sf->rremote+sf->roffset[i],npoints*sizeof(bas->irootloc[0]));CHKERRQ(ierr);
593       continue;
594     }
595     ierr = MPI_Isend(sf->rremote+sf->roffset[i],npoints,MPIU_INT,sf->ranks[i],bas->tag,comm,&leafreqs[i-sf->ndranks]);CHKERRQ(ierr);
596   }
597   ierr = MPI_Waitall(bas->niranks-bas->ndiranks,rootreqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
598   ierr = MPI_Waitall(sf->nranks-sf->ndranks,leafreqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
599   ierr = PetscFree2(rootreqs,leafreqs);CHKERRQ(ierr);
600   PetscFunctionReturn(0);
601 }
602 
603 static PetscErrorCode PetscSFBasicPackTypeSetup(PetscSFBasicPack link,MPI_Datatype unit)
604 {
605   PetscErrorCode ierr;
606   PetscBool      isInt,isPetscInt,isPetscReal,is2Int,is2PetscInt;
607   PetscInt       nPetscIntContig,nPetscRealContig;
608 #if defined(PETSC_HAVE_COMPLEX)
609   PetscBool isPetscComplex;
610   PetscInt nPetscComplexContig;
611 #endif
612 
613   PetscFunctionBegin;
614   ierr = MPIPetsc_Type_compare(unit,MPI_INT,&isInt);CHKERRQ(ierr);
615   ierr = MPIPetsc_Type_compare(unit,MPIU_INT,&isPetscInt);CHKERRQ(ierr);
616   ierr = MPIPetsc_Type_compare_contig(unit,MPIU_INT,&nPetscIntContig);CHKERRQ(ierr);
617   ierr = MPIPetsc_Type_compare(unit,MPIU_REAL,&isPetscReal);CHKERRQ(ierr);
618   ierr = MPIPetsc_Type_compare_contig(unit,MPIU_REAL,&nPetscRealContig);CHKERRQ(ierr);
619 #if defined(PETSC_HAVE_COMPLEX)
620   ierr = MPIPetsc_Type_compare(unit,MPIU_COMPLEX,&isPetscComplex);CHKERRQ(ierr);
621   ierr = MPIPetsc_Type_compare_contig(unit,MPIU_COMPLEX,&nPetscComplexContig);CHKERRQ(ierr);
622 #endif
623   ierr = MPIPetsc_Type_compare(unit,MPI_2INT,&is2Int);CHKERRQ(ierr);
624   ierr = MPIPetsc_Type_compare(unit,MPIU_2INT,&is2PetscInt);CHKERRQ(ierr);
625   link->bs = 1;
626   if (isInt) {PackInit_int(link); PackInit_Logical_int(link); PackInit_Bitwise_int(link);}
627   else if (isPetscInt) {PackInit_PetscInt(link); PackInit_Logical_PetscInt(link); PackInit_Bitwise_PetscInt(link);}
628   else if (isPetscReal) {PackInit_PetscReal(link); PackInit_Logical_PetscReal(link);}
629 #if defined(PETSC_HAVE_COMPLEX)
630   else if (isPetscComplex) PackInit_PetscComplex_1(link);
631 #endif
632   else if (is2Int) PackInit_int_int(link);
633   else if (is2PetscInt) PackInit_PetscInt_PetscInt(link);
634   else if (nPetscIntContig) {
635     if (nPetscIntContig%7 == 0) PackInit_PetscInt_7(link);
636     else if (nPetscIntContig%5 == 0) PackInit_PetscInt_5(link);
637     else if (nPetscIntContig%4 == 0) PackInit_PetscInt_4(link);
638     else if (nPetscIntContig%3 == 0) PackInit_PetscInt_3(link);
639     else if (nPetscIntContig%2 == 0) PackInit_PetscInt_2(link);
640     else PackInit_PetscInt(link);
641     link->bs = nPetscIntContig;
642     link->unitbytes *= nPetscIntContig;
643   } else if (nPetscRealContig) {
644     if (nPetscRealContig%7 == 0) PackInit_PetscReal_7(link);
645     else if (nPetscRealContig%5 == 0) PackInit_PetscReal_5(link);
646     else if (nPetscRealContig%4 == 0) PackInit_PetscReal_4(link);
647     else if (nPetscRealContig%3 == 0) PackInit_PetscReal_3(link);
648     else if (nPetscRealContig%2 == 0) PackInit_PetscReal_2(link);
649     else PackInit_PetscReal(link);
650     link->bs = nPetscRealContig;
651     link->unitbytes *= nPetscRealContig;
652 #if defined(PETSC_HAVE_COMPLEX)
653   } else if (nPetscComplexContig) {
654     if (nPetscComplexContig%7 == 0) PackInit_PetscComplex_7(link);
655     else if (nPetscComplexContig%5 == 0) PackInit_PetscComplex_5(link);
656     else if (nPetscComplexContig%4 == 0) PackInit_PetscComplex_4(link);
657     else if (nPetscComplexContig%3 == 0) PackInit_PetscComplex_3(link);
658     else if (nPetscComplexContig%2 == 0) PackInit_PetscComplex_2(link);
659     else PackInit_PetscComplex_1(link);
660     link->bs = nPetscComplexContig;
661     link->unitbytes *= nPetscComplexContig;
662 #endif
663   } else {
664     MPI_Aint lb,bytes;
665     ierr = MPI_Type_get_extent(unit,&lb,&bytes);CHKERRQ(ierr);
666     if (lb != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Datatype with nonzero lower bound %ld\n",(long)lb);
667     if (bytes % sizeof(int)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for type size not divisible by %D",sizeof(int));
668     switch (bytes / sizeof(int)) {
669     case 1: PackInit_block_int_1(link); break;
670     case 2: PackInit_block_int_2(link); break;
671     case 3: PackInit_block_int_3(link); break;
672     case 4: PackInit_block_int_4(link); break;
673     case 5: PackInit_block_int_5(link); break;
674     case 6: PackInit_block_int_6(link); break;
675     case 7: PackInit_block_int_7(link); break;
676     case 8: PackInit_block_int_8(link); break;
677     default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for arbitrary block sizes");
678     }
679   }
680   ierr = MPI_Type_dup(unit,&link->unit);CHKERRQ(ierr);
681   PetscFunctionReturn(0);
682 }
683 
684 static PetscErrorCode PetscSFBasicPackGetUnpackOp(PetscSF sf,PetscSFBasicPack link,MPI_Op op,void (**UnpackOp)(PetscInt,PetscInt,const PetscInt*,void*,const void*))
685 {
686   PetscFunctionBegin;
687   *UnpackOp = NULL;
688   if (op == MPIU_REPLACE) *UnpackOp = link->UnpackInsert;
689   else if (op == MPI_SUM || op == MPIU_SUM) *UnpackOp = link->UnpackAdd;
690   else if (op == MPI_PROD) *UnpackOp = link->UnpackMult;
691   else if (op == MPI_MAX || op == MPIU_MAX) *UnpackOp = link->UnpackMax;
692   else if (op == MPI_MIN || op == MPIU_MIN) *UnpackOp = link->UnpackMin;
693   else if (op == MPI_LAND) *UnpackOp = link->UnpackLAND;
694   else if (op == MPI_BAND) *UnpackOp = link->UnpackBAND;
695   else if (op == MPI_LOR) *UnpackOp = link->UnpackLOR;
696   else if (op == MPI_BOR) *UnpackOp = link->UnpackBOR;
697   else if (op == MPI_LXOR) *UnpackOp = link->UnpackLXOR;
698   else if (op == MPI_BXOR) *UnpackOp = link->UnpackBXOR;
699   else if (op == MPI_MAXLOC) *UnpackOp = link->UnpackMaxloc;
700   else if (op == MPI_MINLOC) *UnpackOp = link->UnpackMinloc;
701   else *UnpackOp = NULL;
702   PetscFunctionReturn(0);
703 }
704 static PetscErrorCode PetscSFBasicPackGetFetchAndOp(PetscSF sf,PetscSFBasicPack link,MPI_Op op,void (**FetchAndOp)(PetscInt,PetscInt,const PetscInt*,void*,void*))
705 {
706   PetscFunctionBegin;
707   *FetchAndOp = NULL;
708   if (op == MPIU_REPLACE) *FetchAndOp = link->FetchAndInsert;
709   else if (op == MPI_SUM || op == MPIU_SUM) *FetchAndOp = link->FetchAndAdd;
710   else if (op == MPI_MAX || op == MPIU_MAX) *FetchAndOp = link->FetchAndMax;
711   else if (op == MPI_MIN || op == MPIU_MIN) *FetchAndOp = link->FetchAndMin;
712   else if (op == MPI_MAXLOC) *FetchAndOp = link->FetchAndMaxloc;
713   else if (op == MPI_MINLOC) *FetchAndOp = link->FetchAndMinloc;
714   else if (op == MPI_PROD)   *FetchAndOp = link->FetchAndMult;
715   else if (op == MPI_LAND)   *FetchAndOp = link->FetchAndLAND;
716   else if (op == MPI_BAND)   *FetchAndOp = link->FetchAndBAND;
717   else if (op == MPI_LOR)    *FetchAndOp = link->FetchAndLOR;
718   else if (op == MPI_BOR)    *FetchAndOp = link->FetchAndBOR;
719   else if (op == MPI_LXOR)   *FetchAndOp = link->FetchAndLXOR;
720   else if (op == MPI_BXOR)   *FetchAndOp = link->FetchAndBXOR;
721   else SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_SUP,"No support for MPI_Op");
722   PetscFunctionReturn(0);
723 }
724 
725 static PetscErrorCode PetscSFBasicPackGetReqs(PetscSF sf,PetscSFBasicPack link,MPI_Request **rootreqs,MPI_Request **leafreqs)
726 {
727   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
728 
729   PetscFunctionBegin;
730   if (rootreqs) *rootreqs = link->requests;
731   if (leafreqs) *leafreqs = link->requests + (bas->niranks - bas->ndiranks);
732   PetscFunctionReturn(0);
733 }
734 
735 static PetscErrorCode PetscSFBasicPackWaitall(PetscSF sf,PetscSFBasicPack link)
736 {
737   PetscSF_Basic  *bas = (PetscSF_Basic*)sf->data;
738   PetscErrorCode ierr;
739 
740   PetscFunctionBegin;
741   ierr = MPI_Waitall(bas->niranks+sf->nranks-(bas->ndiranks+sf->ndranks),link->requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
742   PetscFunctionReturn(0);
743 }
744 
745 static PetscErrorCode PetscSFBasicGetRootInfo(PetscSF sf,PetscInt *nrootranks,PetscInt *ndrootranks,const PetscMPIInt **rootranks,const PetscInt **rootoffset,const PetscInt **rootloc)
746 {
747   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
748 
749   PetscFunctionBegin;
750   if (nrootranks)  *nrootranks  = bas->niranks;
751   if (ndrootranks) *ndrootranks = bas->ndiranks;
752   if (rootranks)   *rootranks   = bas->iranks;
753   if (rootoffset)  *rootoffset  = bas->ioffset;
754   if (rootloc)     *rootloc     = bas->irootloc;
755   PetscFunctionReturn(0);
756 }
757 
758 static PetscErrorCode PetscSFBasicGetLeafInfo(PetscSF sf,PetscInt *nleafranks,PetscInt *ndleafranks,const PetscMPIInt **leafranks,const PetscInt **leafoffset,const PetscInt **leafloc)
759 {
760   PetscFunctionBegin;
761   if (nleafranks)  *nleafranks  = sf->nranks;
762   if (ndleafranks) *ndleafranks = sf->ndranks;
763   if (leafranks)   *leafranks   = sf->ranks;
764   if (leafoffset)  *leafoffset  = sf->roffset;
765   if (leafloc)     *leafloc     = sf->rmine;
766   PetscFunctionReturn(0);
767 }
768 
769 static PetscErrorCode PetscSFBasicGetPack(PetscSF sf,MPI_Datatype unit,const void *key,PetscSFBasicPack *mylink)
770 {
771   PetscSF_Basic    *bas = (PetscSF_Basic*)sf->data;
772   PetscErrorCode   ierr;
773   PetscSFBasicPack link,*p;
774   PetscInt         nrootranks,ndrootranks,nleafranks,ndleafranks,i;
775   const PetscInt   *rootoffset,*leafoffset;
776 
777   PetscFunctionBegin;
778   /* Look for types in cache */
779   for (p=&bas->avail; (link=*p); p=&link->next) {
780     PetscBool match;
781     ierr = MPIPetsc_Type_compare(unit,link->unit,&match);CHKERRQ(ierr);
782     if (match) {
783       *p = link->next;          /* Remove from available list */
784       goto found;
785     }
786   }
787 
788   /* Create new composite types for each send rank */
789   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,&ndrootranks,NULL,&rootoffset,NULL);CHKERRQ(ierr);
790   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,&ndleafranks,NULL,&leafoffset,NULL);CHKERRQ(ierr);
791   ierr = PetscNew(&link);CHKERRQ(ierr);
792   ierr = PetscSFBasicPackTypeSetup(link,unit);CHKERRQ(ierr);
793   ierr = PetscMalloc2(nrootranks,&link->root,nleafranks,&link->leaf);CHKERRQ(ierr);
794   for (i=0; i<nrootranks; i++) {
795     ierr = PetscMalloc((rootoffset[i+1]-rootoffset[i])*link->unitbytes,&link->root[i]);CHKERRQ(ierr);
796   }
797   for (i=0; i<nleafranks; i++) {
798     if (i < ndleafranks) {      /* Leaf buffers for distinguished ranks are pointers directly into root buffers */
799       if (ndrootranks != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot match distinguished ranks");
800       link->leaf[i] = link->root[0];
801       continue;
802     }
803     ierr = PetscMalloc((leafoffset[i+1]-leafoffset[i])*link->unitbytes,&link->leaf[i]);CHKERRQ(ierr);
804   }
805   ierr = PetscCalloc1(nrootranks+nleafranks,&link->requests);CHKERRQ(ierr);
806 
807 found:
808   link->key  = key;
809   link->next = bas->inuse;
810   bas->inuse = link;
811 
812   *mylink = link;
813   PetscFunctionReturn(0);
814 }
815 
816 static PetscErrorCode PetscSFBasicGetPackInUse(PetscSF sf,MPI_Datatype unit,const void *key,PetscCopyMode cmode,PetscSFBasicPack *mylink)
817 {
818   PetscSF_Basic    *bas = (PetscSF_Basic*)sf->data;
819   PetscErrorCode   ierr;
820   PetscSFBasicPack link,*p;
821 
822   PetscFunctionBegin;
823   /* Look for types in cache */
824   for (p=&bas->inuse; (link=*p); p=&link->next) {
825     PetscBool match;
826     ierr = MPIPetsc_Type_compare(unit,link->unit,&match);CHKERRQ(ierr);
827     if (match && (key == link->key)) {
828       switch (cmode) {
829       case PETSC_OWN_POINTER: *p = link->next; break; /* Remove from inuse list */
830       case PETSC_USE_POINTER: break;
831       default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"invalid cmode");
832       }
833       *mylink = link;
834       PetscFunctionReturn(0);
835     }
836   }
837   SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Could not find pack");
838   PetscFunctionReturn(0);
839 }
840 
841 static PetscErrorCode PetscSFBasicReclaimPack(PetscSF sf,PetscSFBasicPack *link)
842 {
843   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
844 
845   PetscFunctionBegin;
846   (*link)->key  = NULL;
847   (*link)->next = bas->avail;
848   bas->avail    = *link;
849   *link         = NULL;
850   PetscFunctionReturn(0);
851 }
852 
853 static PetscErrorCode PetscSFSetFromOptions_Basic(PetscOptionItems *PetscOptionsObject,PetscSF sf)
854 {
855   PetscErrorCode ierr;
856 
857   PetscFunctionBegin;
858   ierr = PetscOptionsHead(PetscOptionsObject,"PetscSF Basic options");CHKERRQ(ierr);
859   ierr = PetscOptionsTail();CHKERRQ(ierr);
860   PetscFunctionReturn(0);
861 }
862 
863 static PetscErrorCode PetscSFReset_Basic(PetscSF sf)
864 {
865   PetscSF_Basic    *bas = (PetscSF_Basic*)sf->data;
866   PetscErrorCode   ierr;
867   PetscSFBasicPack link,next;
868 
869   PetscFunctionBegin;
870   ierr = PetscFree2(bas->iranks,bas->ioffset);CHKERRQ(ierr);
871   ierr = PetscFree(bas->irootloc);CHKERRQ(ierr);
872   if (bas->inuse) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Outstanding operation has not been completed");
873   for (link=bas->avail; link; link=next) {
874     PetscInt i;
875     next = link->next;
876 #if defined(PETSC_HAVE_MPI_TYPE_DUP) || defined(PETSC_HAVE_MPIUNI)
877     ierr = MPI_Type_free(&link->unit);CHKERRQ(ierr);
878 #endif
879     for (i=0; i<bas->niranks; i++) {ierr = PetscFree(link->root[i]);CHKERRQ(ierr);}
880     for (i=sf->ndranks; i<sf->nranks; i++) {ierr = PetscFree(link->leaf[i]);CHKERRQ(ierr);} /* Free only non-distinguished leaf buffers */
881     ierr = PetscFree2(link->root,link->leaf);CHKERRQ(ierr);
882     ierr = PetscFree(link->requests);CHKERRQ(ierr);
883     ierr = PetscFree(link);CHKERRQ(ierr);
884   }
885   bas->avail = NULL;
886   PetscFunctionReturn(0);
887 }
888 
889 static PetscErrorCode PetscSFDestroy_Basic(PetscSF sf)
890 {
891   PetscErrorCode ierr;
892 
893   PetscFunctionBegin;
894   ierr = PetscFree(sf->data);CHKERRQ(ierr);
895   PetscFunctionReturn(0);
896 }
897 
898 static PetscErrorCode PetscSFView_Basic(PetscSF sf,PetscViewer viewer)
899 {
900   /* PetscSF_Basic *bas = (PetscSF_Basic*)sf->data; */
901   PetscErrorCode ierr;
902   PetscBool      iascii;
903 
904   PetscFunctionBegin;
905   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
906   if (iascii) {
907     ierr = PetscViewerASCIIPrintf(viewer,"  sort=%s\n",sf->rankorder ? "rank-order" : "unordered");CHKERRQ(ierr);
908   }
909   PetscFunctionReturn(0);
910 }
911 
912 /* Send from roots to leaves */
913 static PetscErrorCode PetscSFBcastBegin_Basic(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
914 {
915   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;
916   PetscErrorCode    ierr;
917   PetscSFBasicPack  link;
918   PetscInt          i,nrootranks,ndrootranks,nleafranks,ndleafranks;
919   const PetscInt    *rootoffset,*leafoffset,*rootloc,*leafloc;
920   const PetscMPIInt *rootranks,*leafranks;
921   MPI_Request       *rootreqs,*leafreqs;
922 
923   PetscFunctionBegin;
924   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,&ndrootranks,&rootranks,&rootoffset,&rootloc);CHKERRQ(ierr);
925   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,&ndleafranks,&leafranks,&leafoffset,&leafloc);CHKERRQ(ierr);
926   ierr = PetscSFBasicGetPack(sf,unit,rootdata,&link);CHKERRQ(ierr);
927 
928   ierr = PetscSFBasicPackGetReqs(sf,link,&rootreqs,&leafreqs);CHKERRQ(ierr);
929   /* Eagerly post leaf receives, but only from non-distinguished ranks -- distinguished ranks will receive via shared memory */
930   for (i=ndleafranks; i<nleafranks; i++) {
931     PetscMPIInt n = leafoffset[i+1] - leafoffset[i];
932     ierr = MPI_Irecv(link->leaf[i],n,unit,leafranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&leafreqs[i-ndleafranks]);CHKERRQ(ierr);
933   }
934   /* Pack and send root data */
935   for (i=0; i<nrootranks; i++) {
936     PetscMPIInt n          = rootoffset[i+1] - rootoffset[i];
937     void        *packstart = link->root[i];
938     (*link->Pack)(n,link->bs,rootloc+rootoffset[i],rootdata,packstart);
939     if (i < ndrootranks) continue; /* shared memory */
940     ierr = MPI_Isend(packstart,n,unit,rootranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&rootreqs[i-ndrootranks]);CHKERRQ(ierr);
941   }
942   PetscFunctionReturn(0);
943 }
944 
945 PetscErrorCode PetscSFBcastEnd_Basic(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
946 {
947   PetscErrorCode   ierr;
948   PetscSFBasicPack link;
949   PetscInt         i,nleafranks,ndleafranks;
950   const PetscInt   *leafoffset,*leafloc;
951 
952   PetscFunctionBegin;
953   ierr = PetscSFBasicGetPackInUse(sf,unit,rootdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr);
954   ierr = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
955   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,&ndleafranks,NULL,&leafoffset,&leafloc);CHKERRQ(ierr);
956   for (i=0; i<nleafranks; i++) {
957     PetscMPIInt n          = leafoffset[i+1] - leafoffset[i];
958     const void  *packstart = link->leaf[i];
959     (*link->UnpackInsert)(n,link->bs,leafloc+leafoffset[i],leafdata,packstart);
960   }
961   ierr = PetscSFBasicReclaimPack(sf,&link);CHKERRQ(ierr);
962   PetscFunctionReturn(0);
963 }
964 
965 /* leaf -> root with reduction */
966 PetscErrorCode PetscSFReduceBegin_Basic(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
967 {
968   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;
969   PetscSFBasicPack  link;
970   PetscErrorCode    ierr;
971   PetscInt          i,nrootranks,ndrootranks,nleafranks,ndleafranks;
972   const PetscInt    *rootoffset,*leafoffset,*rootloc,*leafloc;
973   const PetscMPIInt *rootranks,*leafranks;
974   MPI_Request       *rootreqs,*leafreqs;
975 
976   PetscFunctionBegin;
977   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,&ndrootranks,&rootranks,&rootoffset,&rootloc);CHKERRQ(ierr);
978   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,&ndleafranks,&leafranks,&leafoffset,&leafloc);CHKERRQ(ierr);
979   ierr = PetscSFBasicGetPack(sf,unit,rootdata,&link);CHKERRQ(ierr);
980 
981   ierr = PetscSFBasicPackGetReqs(sf,link,&rootreqs,&leafreqs);CHKERRQ(ierr);
982   /* Eagerly post root receives for non-distinguished ranks */
983   for (i=ndrootranks; i<nrootranks; i++) {
984     PetscMPIInt n = rootoffset[i+1] - rootoffset[i];
985     ierr = MPI_Irecv(link->root[i],n,unit,rootranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&rootreqs[i-ndrootranks]);CHKERRQ(ierr);
986   }
987   /* Pack and send leaf data */
988   for (i=0; i<nleafranks; i++) {
989     PetscMPIInt n          = leafoffset[i+1] - leafoffset[i];
990     void        *packstart = link->leaf[i];
991     (*link->Pack)(n,link->bs,leafloc+leafoffset[i],leafdata,packstart);
992     if (i < ndleafranks) continue; /* shared memory */
993     ierr = MPI_Isend(packstart,n,unit,leafranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&leafreqs[i-ndleafranks]);CHKERRQ(ierr);
994   }
995   PetscFunctionReturn(0);
996 }
997 
998 static PetscErrorCode PetscSFReduceEnd_Basic(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
999 {
1000   void             (*UnpackOp)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
1001   PetscErrorCode   ierr;
1002   PetscSFBasicPack link;
1003   PetscInt         i,nrootranks;
1004   PetscMPIInt      typesize = -1;
1005   const PetscInt   *rootoffset,*rootloc;
1006 
1007   PetscFunctionBegin;
1008   ierr = PetscSFBasicGetPackInUse(sf,unit,rootdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr);
1009   /* This implementation could be changed to unpack as receives arrive, at the cost of non-determinism */
1010   ierr = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
1011   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,NULL,NULL,&rootoffset,&rootloc);CHKERRQ(ierr);
1012   ierr = PetscSFBasicPackGetUnpackOp(sf,link,op,&UnpackOp);CHKERRQ(ierr);
1013   if (UnpackOp) {
1014     typesize = link->unitbytes;
1015   }
1016   else {
1017     ierr = MPI_Type_size(unit,&typesize);CHKERRQ(ierr);
1018   }
1019   for (i=0; i<nrootranks; i++) {
1020     PetscMPIInt n   = rootoffset[i+1] - rootoffset[i];
1021     char *packstart = (char *) link->root[i];
1022 
1023     if (UnpackOp) {
1024       (*UnpackOp)(n,link->bs,rootloc+rootoffset[i],rootdata,(const void *)packstart);
1025     }
1026 #if PETSC_HAVE_MPI_REDUCE_LOCAL
1027     else if (n) { /* the op should be defined to operate on the whole datatype, so we ignore link->bs */
1028       PetscInt j;
1029 
1030       for (j = 0; j < n; j++) {
1031         ierr = MPI_Reduce_local(packstart+j*typesize,((char *) rootdata)+(rootloc[rootoffset[i]+j])*typesize,1,unit,op);CHKERRQ(ierr);
1032       }
1033     }
1034 #else
1035     else {
1036       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No unpacking reduction operation for this MPI_Op");
1037     }
1038 #endif
1039   }
1040   ierr = PetscSFBasicReclaimPack(sf,&link);CHKERRQ(ierr);
1041   PetscFunctionReturn(0);
1042 }
1043 
1044 static PetscErrorCode PetscSFFetchAndOpBegin_Basic(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1045 {
1046   PetscErrorCode ierr;
1047 
1048   PetscFunctionBegin;
1049   ierr = PetscSFReduceBegin_Basic(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
1050   PetscFunctionReturn(0);
1051 }
1052 
1053 static PetscErrorCode PetscSFFetchAndOpEnd_Basic(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1054 {
1055   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;
1056   void              (*FetchAndOp)(PetscInt,PetscInt,const PetscInt*,void*,void*);
1057   PetscErrorCode    ierr;
1058   PetscSFBasicPack  link;
1059   PetscInt          i,nrootranks,ndrootranks,nleafranks,ndleafranks;
1060   const PetscInt    *rootoffset,*leafoffset,*rootloc,*leafloc;
1061   const PetscMPIInt *rootranks,*leafranks;
1062   MPI_Request       *rootreqs,*leafreqs;
1063 
1064   PetscFunctionBegin;
1065   ierr = PetscSFBasicGetPackInUse(sf,unit,rootdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr);
1066   /* This implementation could be changed to unpack as receives arrive, at the cost of non-determinism */
1067   ierr      = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
1068   ierr      = PetscSFBasicGetRootInfo(sf,&nrootranks,&ndrootranks,&rootranks,&rootoffset,&rootloc);CHKERRQ(ierr);
1069   ierr      = PetscSFBasicGetLeafInfo(sf,&nleafranks,&ndleafranks,&leafranks,&leafoffset,&leafloc);CHKERRQ(ierr);
1070   ierr      = PetscSFBasicPackGetReqs(sf,link,&rootreqs,&leafreqs);CHKERRQ(ierr);
1071   /* Post leaf receives */
1072   for (i=ndleafranks; i<nleafranks; i++) {
1073     PetscMPIInt n = leafoffset[i+1] - leafoffset[i];
1074     ierr = MPI_Irecv(link->leaf[i],n,unit,leafranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&leafreqs[i-ndleafranks]);CHKERRQ(ierr);
1075   }
1076   /* Process local fetch-and-op, post root sends */
1077   ierr = PetscSFBasicPackGetFetchAndOp(sf,link,op,&FetchAndOp);CHKERRQ(ierr);
1078   for (i=0; i<nrootranks; i++) {
1079     PetscMPIInt n          = rootoffset[i+1] - rootoffset[i];
1080     void        *packstart = link->root[i];
1081 
1082     (*FetchAndOp)(n,link->bs,rootloc+rootoffset[i],rootdata,packstart);
1083     if (i < ndrootranks) continue; /* shared memory */
1084     ierr = MPI_Isend(packstart,n,unit,rootranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&rootreqs[i-ndrootranks]);CHKERRQ(ierr);
1085   }
1086   ierr = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
1087   for (i=0; i<nleafranks; i++) {
1088     PetscMPIInt n          = leafoffset[i+1] - leafoffset[i];
1089     const void  *packstart = link->leaf[i];
1090     (*link->UnpackInsert)(n,link->bs,leafloc+leafoffset[i],leafupdate,packstart);
1091   }
1092   ierr = PetscSFBasicReclaimPack(sf,&link);CHKERRQ(ierr);
1093   PetscFunctionReturn(0);
1094 }
1095 
1096 PETSC_EXTERN PetscErrorCode PetscSFCreate_Basic(PetscSF sf)
1097 {
1098   PetscSF_Basic  *bas = (PetscSF_Basic*)sf->data;
1099   PetscErrorCode ierr;
1100 
1101   PetscFunctionBegin;
1102   sf->ops->SetUp           = PetscSFSetUp_Basic;
1103   sf->ops->SetFromOptions  = PetscSFSetFromOptions_Basic;
1104   sf->ops->Reset           = PetscSFReset_Basic;
1105   sf->ops->Destroy         = PetscSFDestroy_Basic;
1106   sf->ops->View            = PetscSFView_Basic;
1107   sf->ops->BcastBegin      = PetscSFBcastBegin_Basic;
1108   sf->ops->BcastEnd        = PetscSFBcastEnd_Basic;
1109   sf->ops->ReduceBegin     = PetscSFReduceBegin_Basic;
1110   sf->ops->ReduceEnd       = PetscSFReduceEnd_Basic;
1111   sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Basic;
1112   sf->ops->FetchAndOpEnd   = PetscSFFetchAndOpEnd_Basic;
1113 
1114   ierr     = PetscNewLog(sf,&bas);CHKERRQ(ierr);
1115   sf->data = (void*)bas;
1116   PetscFunctionReturn(0);
1117 }
1118