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