Actual source code: sfbasic.c

petsc-3.11.4 2019-09-28
Report Typos and Errors

  2:  #include <petsc/private/sfimpl.h>

  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*);

 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       *rkey,*lkey; /* rootdata and leafdata 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: };

 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;

 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:   MPI_Type_contiguous(1,datatype,newtype); if (ierr) return ierr;
 62:   MPI_Type_commit(newtype); if (ierr) return ierr;
 63:   return MPI_SUCCESS;
 64: }
 65: #endif

 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:  */

 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)

 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:   }

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:   }

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:   }


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:   }

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:   }

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:   }

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)

535: static PetscErrorCode PetscSFSetUp_Basic(PetscSF sf)
536: {
537:   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
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;

547:   MPI_Comm_group(PETSC_COMM_SELF,&group);
548:   PetscSFSetUpRanks(sf,group);
549:   MPI_Group_free(&group);
550:   PetscObjectGetComm((PetscObject)sf,&comm);
551:   PetscObjectGetNewTag((PetscObject)sf,&bas->tag);
552:   MPI_Comm_rank(comm,&rank);
553:   /*
554:    * Inform roots about how many leaves and from which ranks
555:    */
556:   PetscMalloc1(sf->nranks,&rlengths);
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:   PetscCommBuildTwoSided(comm,1,MPIU_INT,sf->nranks-sf->ndranks,sf->ranks+sf->ndranks,rlengths+sf->ndranks,&niranks,&iranks,&ilengths);

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:   PetscSortMPIIntWithIntArray(niranks,iranks,ilengths);

569:   /* Partition into distinguished and non-distinguished incoming ranks */
570:   bas->ndiranks = sf->ndranks;
571:   bas->niranks = bas->ndiranks + niranks;
572:   PetscMalloc2(bas->niranks,&bas->iranks,bas->niranks+1,&bas->ioffset);
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:   PetscFree(rlengths);
584:   PetscFree(iranks);
585:   PetscFree(ilengths);

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");

592:   /* Send leaf identities to roots */
593:   PetscMalloc1(bas->itotal,&bas->irootloc);
594:   PetscMalloc1(bas->niranks-bas->ndiranks+sf->nranks-sf->ndranks,&reqs);
595:   for (i=sf->ndranks; i<sf->nranks; i++, nreqs++) {
596:     PetscMPIInt npoints;
597:     PetscMPIIntCast(sf->roffset[i+1]-sf->roffset[i],&npoints);
598:     MPI_Isend(sf->rremote+sf->roffset[i],npoints,MPIU_INT,sf->ranks[i],bas->tag,comm,&reqs[nreqs]);
599:   }
600:   for (i=bas->ndiranks; i<bas->niranks; i++, nreqs++) {
601:     PetscMPIInt npoints;
602:     PetscMPIIntCast(bas->ioffset[i+1]-bas->ioffset[i],&npoints);
603:     MPI_Irecv(bas->irootloc+bas->ioffset[i],npoints,MPIU_INT,bas->iranks[i],bas->tag,comm,&reqs[nreqs]);
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:     PetscMemcpy(bas->irootloc+bas->ioffset[i],sf->rremote+sf->roffset[i],npoints*sizeof(bas->irootloc[0]));
609:   }
610:   MPI_Waitall(nreqs,reqs,MPI_STATUSES_IGNORE);
611:   PetscFree(reqs);
612:   return(0);
613: }

615: static PetscErrorCode PetscSFBasicPackTypeSetup(PetscSFBasicPack link,MPI_Datatype unit)
616: {
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

627:   MPIPetsc_Type_compare(unit,MPI_INT,&isInt);
628:   MPIPetsc_Type_compare(unit,MPIU_INT,&isPetscInt);
629:   MPIPetsc_Type_compare_contig(unit,MPIU_INT,&nPetscIntContig);
630:   MPIPetsc_Type_compare(unit,MPIU_REAL,&isPetscReal);
631:   MPIPetsc_Type_compare_contig(unit,MPIU_REAL,&nPetscRealContig);
632: #if defined(PETSC_HAVE_COMPLEX)
633:   MPIPetsc_Type_compare(unit,MPIU_COMPLEX,&isPetscComplex);
634:   MPIPetsc_Type_compare_contig(unit,MPIU_COMPLEX,&nPetscComplexContig);
635: #endif
636:   MPIPetsc_Type_compare(unit,MPI_2INT,&is2Int);
637:   MPIPetsc_Type_compare(unit,MPIU_2INT,&is2PetscInt);
638:   MPI_Type_get_envelope(unit,&ni,&na,&nd,&combiner);
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:     MPI_Type_get_extent(unit,&lb,&bytes);
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 {MPI_Type_dup(unit,&link->unit);}
697:   return(0);
698: }

700: static PetscErrorCode PetscSFBasicPackGetUnpackOp(PetscSF sf,PetscSFBasicPack link,MPI_Op op,void (**UnpackOp)(PetscInt,PetscInt,const PetscInt*,void*,const void*))
701: {
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:   return(0);
719: }
720: static PetscErrorCode PetscSFBasicPackGetFetchAndOp(PetscSF sf,PetscSFBasicPack link,MPI_Op op,void (**FetchAndOp)(PetscInt,PetscInt,const PetscInt*,void*,void*))
721: {
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:   return(0);
739: }

741: typedef enum {PETSC_SF_LEAF2../../../../../.._REDUCE, PETSC_SF_../../../../../..2LEAF_BCAST} PetscSFDirection;

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_LEAF2../../../../../.._REDUCE)? 0 : (sf->nranks + bas->niranks); /* reduce reqs are in the front, bcast reqs are at the end */

749:   if (rootreqs) *rootreqs = link->requests + shift;
750:   if (leafreqs) *leafreqs = link->requests + (bas->niranks - bas->ndiranks) + shift;
751:   return(0);
752: }

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_LEAF2../../../../../.._REDUCE)? 0 : (sf->nranks + bas->niranks);

761:   MPI_Waitall(bas->niranks+sf->nranks-(bas->ndiranks+sf->ndranks),link->requests+shift,MPI_STATUSES_IGNORE);
762:   return(0);
763: }

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;

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:   return(0);
776: }

778: static PetscErrorCode PetscSFBasicGetLeafInfo(PetscSF sf,PetscInt *nleafranks,PetscInt *ndleafranks,const PetscMPIInt **leafranks,const PetscInt **leafoffset,const PetscInt **leafloc)
779: {
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:   return(0);
787: }

789: static PetscErrorCode PetscSFBasicGetPack(PetscSF sf,MPI_Datatype unit,const void *rkey,const void *lkey,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,tag;
798:   MPI_Request      *rootreqs,*leafreqs;
799:   PetscBool        match;

802:   /* Look up links in use and error out if there is a match. When both rootdata and leafdata are NULL, ignore
803:      the potential overlapping since this process does not participate in communication. Overlapping is harmless.
804:   */
805:   if (rkey || lkey) {
806:     for (p=&bas->inuse; (link=*p); p=&link->next) {
807:       MPIPetsc_Type_compare(unit,link->unit,&match);
808:       if (match && (rkey == link->rkey) && (lkey == link->lkey)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for overlapped PetscSF communications with the same SF, rootdata, leafdatadata and data type. You can undo the overlap to avoid the error.");
809:     }
810:   }

812:   /* Look up free links */
813:   for (p=&bas->avail; (link=*p); p=&link->next) {
814:     MPIPetsc_Type_compare(unit,link->unit,&match);
815:     if (match) {
816:       *p = link->next; /* Remove from available list */
817:       goto found;
818:     }
819:   }

821:   /* Create new composite types for each send rank */
822:   PetscSFBasicGetRootInfo(sf,&nrootranks,&ndrootranks,NULL,&rootoffset,NULL);
823:   PetscSFBasicGetLeafInfo(sf,&nleafranks,&ndleafranks,NULL,&leafoffset,NULL);
824:   PetscNew(&link);
825:   PetscSFBasicPackTypeSetup(link,unit);
826:   PetscMalloc2(nrootranks,&link->root,nleafranks,&link->leaf);
827:   /* Double the requests. First half are used for reduce (leaf to root) communication, second half for bcast (root to leaf) communication */
828:   half     = nrootranks + nleafranks;
829:   PetscCalloc1(half*2,&link->requests);
830:   rootreqs = link->requests;
831:   leafreqs = link->requests + bas->niranks - bas->ndiranks;
832:   comm     = PetscObjectComm((PetscObject)sf);
833:   PetscCommGetNewTag(comm,&tag);

835:   /* Allocate buffer and then init the persistent communcation */
836:   for (i=0; i<nrootranks; i++) {
837:     PetscMalloc((rootoffset[i+1]-rootoffset[i])*link->unitbytes,&link->root[i]);
838:     if (i >= ndrootranks) {
839:       PetscMPIIntCast(rootoffset[i+1]-rootoffset[i],&n);
840:       MPI_Recv_init(link->root[i],n,unit,bas->iranks[i],tag,comm,&rootreqs[i-ndrootranks]);      /* reduce */
841:       MPI_Send_init(link->root[i],n,unit,bas->iranks[i],tag,comm,&rootreqs[i-ndrootranks+half]); /* bcast  */
842:     }
843:   }
844:   for (i=0; i<nleafranks; i++) {
845:     if (i < ndleafranks) {      /* Leaf buffers for distinguished ranks are pointers directly into root buffers */
846:       if (ndrootranks != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot match distinguished ranks");
847:       link->leaf[i] = link->root[0];
848:       continue;
849:     }
850:     PetscMalloc((leafoffset[i+1]-leafoffset[i])*link->unitbytes,&link->leaf[i]);
851:     PetscMPIIntCast(leafoffset[i+1]-leafoffset[i],&n);
852:     MPI_Send_init(link->leaf[i],n,unit,sf->ranks[i],tag,comm,&leafreqs[i-ndleafranks]);      /* reduce */
853:     MPI_Recv_init(link->leaf[i],n,unit,sf->ranks[i],tag,comm,&leafreqs[i-ndleafranks+half]); /* bcast  */
854:   }

856: found:
857:   link->rkey = rkey;
858:   link->lkey = lkey;
859:   link->next = bas->inuse;
860:   bas->inuse = link;

862:   *mylink = link;
863:   return(0);
864: }

866: static PetscErrorCode PetscSFBasicGetPackInUse(PetscSF sf,MPI_Datatype unit,const void *rkey,const void *lkey,PetscCopyMode cmode,PetscSFBasicPack *mylink)
867: {
868:   PetscSF_Basic    *bas = (PetscSF_Basic*)sf->data;
869:   PetscErrorCode   ierr;
870:   PetscSFBasicPack link,*p;

873:   /* Look for types in cache */
874:   for (p=&bas->inuse; (link=*p); p=&link->next) {
875:     PetscBool match;
876:     MPIPetsc_Type_compare(unit,link->unit,&match);
877:     if (match && (rkey == link->rkey) && (lkey == link->lkey)) {
878:       switch (cmode) {
879:       case PETSC_OWN_POINTER: *p = link->next; break; /* Remove from inuse list */
880:       case PETSC_USE_POINTER: break;
881:       default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"invalid cmode");
882:       }
883:       *mylink = link;
884:       return(0);
885:     }
886:   }
887:   SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Could not find pack");
888:   return(0);
889: }

891: static PetscErrorCode PetscSFBasicReclaimPack(PetscSF sf,PetscSFBasicPack *link)
892: {
893:   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;

896:   (*link)->rkey = NULL;
897:   (*link)->lkey = NULL;
898:   (*link)->next = bas->avail;
899:   bas->avail    = *link;
900:   *link         = NULL;
901:   return(0);
902: }

904: static PetscErrorCode PetscSFSetFromOptions_Basic(PetscOptionItems *PetscOptionsObject,PetscSF sf)
905: {

909:   PetscOptionsHead(PetscOptionsObject,"PetscSF Basic options");
910:   PetscOptionsTail();
911:   return(0);
912: }

914: static PetscErrorCode PetscSFReset_Basic(PetscSF sf)
915: {
916:   PetscSF_Basic    *bas = (PetscSF_Basic*)sf->data;
917:   PetscErrorCode   ierr;
918:   PetscSFBasicPack link,next;

921:   if (bas->inuse) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Outstanding operation has not been completed");
922:   PetscFree2(bas->iranks,bas->ioffset);
923:   PetscFree(bas->irootloc);
924:   for (link=bas->avail; link; link=next) {
925:     PetscInt i;
926:     next = link->next;
927:     if (!link->isbuiltin) {MPI_Type_free(&link->unit);}
928:     for (i=0; i<bas->niranks; i++) {PetscFree(link->root[i]);}
929:     for (i=sf->ndranks; i<sf->nranks; i++) {PetscFree(link->leaf[i]);} /* Free only non-distinguished leaf buffers */
930:     PetscFree2(link->root,link->leaf);
931:     /* Free persistent requests using MPI_Request_free */
932:     for (i=0; i<sf->nranks+bas->niranks-(sf->ndranks+bas->ndiranks); i++) {
933:       MPI_Request_free(&link->requests[i]); /* used in reduce */
934:       MPI_Request_free(&link->requests[sf->nranks+bas->niranks+i]); /* used in bcast */
935:     }
936:     PetscFree(link->requests);
937:     PetscFree(link);
938:   }
939:   bas->avail = NULL;
940:   return(0);
941: }

943: static PetscErrorCode PetscSFDestroy_Basic(PetscSF sf)
944: {

948:   PetscSFReset_Basic(sf);
949:   PetscFree(sf->data);
950:   return(0);
951: }

953: static PetscErrorCode PetscSFView_Basic(PetscSF sf,PetscViewer viewer)
954: {
955:   /* PetscSF_Basic *bas = (PetscSF_Basic*)sf->data; */
957:   PetscBool      iascii;

960:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
961:   if (iascii) {
962:     PetscViewerASCIIPrintf(viewer,"  sort=%s\n",sf->rankorder ? "rank-order" : "unordered");
963:   }
964:   return(0);
965: }

967: static PetscErrorCode PetscSFBcastAndOpBegin_Basic(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
968: {
969:   PetscErrorCode    ierr;
970:   PetscSFBasicPack  link;
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:   PetscMPIInt       n;

978:   PetscSFBasicGetRootInfo(sf,&nrootranks,&ndrootranks,&rootranks,&rootoffset,&rootloc);
979:   PetscSFBasicGetLeafInfo(sf,&nleafranks,&ndleafranks,&leafranks,&leafoffset,&leafloc);
980:   PetscSFBasicGetPack(sf,unit,rootdata,leafdata,&link);

982:   PetscSFBasicPackGetReqs(sf,link,PETSC_SF_../../../../../..2LEAF_BCAST,&rootreqs,&leafreqs);
983:   /* Eagerly post leaf receives, but only from non-distinguished ranks -- distinguished ranks will receive via shared memory */
984:   PetscMPIIntCast(leafoffset[nleafranks]-leafoffset[ndleafranks],&n);
985:   MPI_Startall_irecv(n,unit,nleafranks-ndleafranks,leafreqs);

987:   /* Pack and send root data */
988:   for (i=0; i<nrootranks; i++) {
989:     void *packstart = link->root[i];
990:     PetscMPIIntCast(rootoffset[i+1]-rootoffset[i],&n);
991:     (*link->Pack)(n,link->bs,rootloc+rootoffset[i],rootdata,packstart);
992:     if (i < ndrootranks) continue; /* shared memory */
993:     MPI_Start_isend(n,unit,&rootreqs[i-ndrootranks]);
994:   }
995:   return(0);
996: }

998: PetscErrorCode PetscSFBcastAndOpEnd_Basic(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
999: {
1000:   PetscErrorCode   ierr;
1001:   PetscSFBasicPack link;
1002:   PetscInt         i,nleafranks,ndleafranks;
1003:   const PetscInt   *leafoffset,*leafloc;
1004:   void             (*UnpackOp)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
1005:   PetscMPIInt      typesize = -1;

1008:   PetscSFBasicGetPackInUse(sf,unit,rootdata,leafdata,PETSC_OWN_POINTER,&link);
1009:   PetscSFBasicPackWaitall(sf,link,PETSC_SF_../../../../../..2LEAF_BCAST);
1010:   PetscSFBasicGetLeafInfo(sf,&nleafranks,&ndleafranks,NULL,&leafoffset,&leafloc);
1011:   PetscSFBasicPackGetUnpackOp(sf,link,op,&UnpackOp);

1013:   if (UnpackOp) { typesize = link->unitbytes; }
1014:   else { MPI_Type_size(unit,&typesize); }

1016:   for (i=0; i<nleafranks; i++) {
1017:     PetscMPIInt n   = leafoffset[i+1] - leafoffset[i];
1018:     char *packstart = (char *) link->leaf[i];
1019:     if (UnpackOp) { (*UnpackOp)(n,link->bs,leafloc+leafoffset[i],leafdata,(const void *)packstart); }
1020: #if defined(PETSC_HAVE_MPI_REDUCE_LOCAL)
1021:     else if (n) { /* the op should be defined to operate on the whole datatype, so we ignore link->bs */
1022:       PetscInt j;
1023:       for (j=0; j<n; j++) { MPI_Reduce_local(packstart+j*typesize,((char *) leafdata)+(leafloc[leafoffset[i]+j])*typesize,1,unit,op); }
1024:     }
1025: #else
1026:     else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No unpacking reduction operation for this MPI_Op");
1027: #endif
1028:   }

1030:   PetscSFBasicReclaimPack(sf,&link);
1031:   return(0);
1032: }

1034: /* Send from roots to leaves */
1035: static PetscErrorCode PetscSFBcastBegin_Basic(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
1036: {
1037:   PetscErrorCode   ierr;

1040:   PetscSFBcastAndOpBegin_Basic(sf,unit,rootdata,leafdata,MPI_REPLACE);
1041:   return(0);
1042: }

1044: PetscErrorCode PetscSFBcastEnd_Basic(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
1045: {
1046:   PetscErrorCode   ierr;

1049:   PetscSFBcastAndOpEnd_Basic(sf,unit,rootdata,leafdata,MPI_REPLACE);
1050:   return(0);
1051: }

1053: /* leaf -> root with reduction */
1054: PetscErrorCode PetscSFReduceBegin_Basic(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
1055: {
1056:   PetscSFBasicPack  link;
1057:   PetscErrorCode    ierr;
1058:   PetscInt          i,nrootranks,ndrootranks,nleafranks,ndleafranks;
1059:   const PetscInt    *rootoffset,*leafoffset,*rootloc,*leafloc;
1060:   const PetscMPIInt *rootranks,*leafranks;
1061:   MPI_Request       *rootreqs,*leafreqs;
1062:   PetscMPIInt       n;

1065:   PetscSFBasicGetRootInfo(sf,&nrootranks,&ndrootranks,&rootranks,&rootoffset,&rootloc);
1066:   PetscSFBasicGetLeafInfo(sf,&nleafranks,&ndleafranks,&leafranks,&leafoffset,&leafloc);
1067:   PetscSFBasicGetPack(sf,unit,rootdata,leafdata,&link);

1069:   PetscSFBasicPackGetReqs(sf,link,PETSC_SF_LEAF2../../../../../.._REDUCE,&rootreqs,&leafreqs);
1070:   /* Eagerly post root receives for non-distinguished ranks */
1071:   PetscMPIIntCast(rootoffset[nrootranks]-rootoffset[ndrootranks],&n);
1072:   MPI_Startall_irecv(n,unit,nrootranks-ndrootranks,rootreqs);

1074:   /* Pack and send leaf data */
1075:   for (i=0; i<nleafranks; i++) {
1076:     void *packstart = link->leaf[i];
1077:     PetscMPIIntCast(leafoffset[i+1]-leafoffset[i],&n);
1078:     (*link->Pack)(n,link->bs,leafloc+leafoffset[i],leafdata,packstart);
1079:     if (i < ndleafranks) continue; /* shared memory */
1080:     MPI_Start_isend(n,unit,&leafreqs[i-ndleafranks]);
1081:   }
1082:   return(0);
1083: }

1085: static PetscErrorCode PetscSFReduceEnd_Basic(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
1086: {
1087:   void             (*UnpackOp)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
1088:   PetscErrorCode   ierr;
1089:   PetscSFBasicPack link;
1090:   PetscInt         i,nrootranks;
1091:   PetscMPIInt      typesize = -1;
1092:   const PetscInt   *rootoffset,*rootloc;

1095:   PetscSFBasicGetPackInUse(sf,unit,rootdata,leafdata,PETSC_OWN_POINTER,&link);
1096:   /* This implementation could be changed to unpack as receives arrive, at the cost of non-determinism */
1097:   PetscSFBasicPackWaitall(sf,link,PETSC_SF_LEAF2../../../../../.._REDUCE);
1098:   PetscSFBasicGetRootInfo(sf,&nrootranks,NULL,NULL,&rootoffset,&rootloc);
1099:   PetscSFBasicPackGetUnpackOp(sf,link,op,&UnpackOp);
1100:   if (UnpackOp) {
1101:     typesize = link->unitbytes;
1102:   }
1103:   else {
1104:     MPI_Type_size(unit,&typesize);
1105:   }
1106:   for (i=0; i<nrootranks; i++) {
1107:     PetscMPIInt n   = rootoffset[i+1] - rootoffset[i];
1108:     char *packstart = (char *) link->root[i];

1110:     if (UnpackOp) {
1111:       (*UnpackOp)(n,link->bs,rootloc+rootoffset[i],rootdata,(const void *)packstart);
1112:     }
1113: #if defined(PETSC_HAVE_MPI_REDUCE_LOCAL)
1114:     else if (n) { /* the op should be defined to operate on the whole datatype, so we ignore link->bs */
1115:       PetscInt j;

1117:       for (j = 0; j < n; j++) {
1118:         MPI_Reduce_local(packstart+j*typesize,((char *) rootdata)+(rootloc[rootoffset[i]+j])*typesize,1,unit,op);
1119:       }
1120:     }
1121: #else
1122:     else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No unpacking reduction operation for this MPI_Op");
1123: #endif
1124:   }
1125:   PetscSFBasicReclaimPack(sf,&link);
1126:   return(0);
1127: }

1129: static PetscErrorCode PetscSFFetchAndOpBegin_Basic(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1130: {

1134:   PetscSFReduceBegin_Basic(sf,unit,leafdata,rootdata,op);
1135:   return(0);
1136: }

1138: static PetscErrorCode PetscSFFetchAndOpEnd_Basic(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1139: {
1140:   void              (*FetchAndOp)(PetscInt,PetscInt,const PetscInt*,void*,void*);
1141:   PetscErrorCode    ierr;
1142:   PetscSFBasicPack  link;
1143:   PetscInt          i,nrootranks,ndrootranks,nleafranks,ndleafranks;
1144:   const PetscInt    *rootoffset,*leafoffset,*rootloc,*leafloc;
1145:   const PetscMPIInt *rootranks,*leafranks;
1146:   MPI_Request       *rootreqs,*leafreqs;
1147:   PetscMPIInt       n;

1150:   PetscSFBasicGetPackInUse(sf,unit,rootdata,leafdata,PETSC_OWN_POINTER,&link);
1151:   /* This implementation could be changed to unpack as receives arrive, at the cost of non-determinism */
1152:   PetscSFBasicPackWaitall(sf,link,PETSC_SF_LEAF2../../../../../.._REDUCE);
1153:   PetscSFBasicGetRootInfo(sf,&nrootranks,&ndrootranks,&rootranks,&rootoffset,&rootloc);
1154:   PetscSFBasicGetLeafInfo(sf,&nleafranks,&ndleafranks,&leafranks,&leafoffset,&leafloc);
1155:   PetscSFBasicPackGetReqs(sf,link,PETSC_SF_../../../../../..2LEAF_BCAST,&rootreqs,&leafreqs);
1156:   /* Post leaf receives */
1157:   PetscMPIIntCast(leafoffset[nleafranks]-leafoffset[ndleafranks],&n);
1158:   MPI_Startall_irecv(n,unit,nleafranks-ndleafranks,leafreqs);

1160:   /* Process local fetch-and-op, post root sends */
1161:   PetscSFBasicPackGetFetchAndOp(sf,link,op,&FetchAndOp);
1162:   for (i=0; i<nrootranks; i++) {
1163:     void *packstart = link->root[i];
1164:     PetscMPIIntCast(rootoffset[i+1]-rootoffset[i],&n);
1165:     (*FetchAndOp)(n,link->bs,rootloc+rootoffset[i],rootdata,packstart);
1166:     if (i < ndrootranks) continue; /* shared memory */
1167:     MPI_Start_isend(n,unit,&rootreqs[i-ndrootranks]);
1168:   }
1169:   PetscSFBasicPackWaitall(sf,link,PETSC_SF_../../../../../..2LEAF_BCAST);
1170:   for (i=0; i<nleafranks; i++) {
1171:     const void  *packstart = link->leaf[i];
1172:     PetscMPIIntCast(leafoffset[i+1]-leafoffset[i],&n);
1173:     (*link->UnpackInsert)(n,link->bs,leafloc+leafoffset[i],leafupdate,packstart);
1174:   }
1175:   PetscSFBasicReclaimPack(sf,&link);
1176:   return(0);
1177: }

1179: static PetscErrorCode PetscSFGetLeafRanks_Basic(PetscSF sf,PetscInt *niranks,const PetscMPIInt **iranks,const PetscInt **ioffset,const PetscInt **irootloc)
1180: {
1181:   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;

1184:   if (niranks)  *niranks  = bas->niranks;
1185:   if (iranks)   *iranks   = bas->iranks;
1186:   if (ioffset)  *ioffset  = bas->ioffset;
1187:   if (irootloc) *irootloc = bas->irootloc;
1188:   return(0);
1189: }

1191: PETSC_EXTERN PetscErrorCode PetscSFCreate_Basic(PetscSF sf)
1192: {
1193:   PetscSF_Basic  *bas = (PetscSF_Basic*)sf->data;

1197:   sf->ops->SetUp           = PetscSFSetUp_Basic;
1198:   sf->ops->SetFromOptions  = PetscSFSetFromOptions_Basic;
1199:   sf->ops->Reset           = PetscSFReset_Basic;
1200:   sf->ops->Destroy         = PetscSFDestroy_Basic;
1201:   sf->ops->View            = PetscSFView_Basic;
1202:   sf->ops->BcastBegin      = PetscSFBcastBegin_Basic;
1203:   sf->ops->BcastEnd        = PetscSFBcastEnd_Basic;
1204:   sf->ops->BcastAndOpBegin = PetscSFBcastAndOpBegin_Basic;
1205:   sf->ops->BcastAndOpEnd   = PetscSFBcastAndOpEnd_Basic;
1206:   sf->ops->ReduceBegin     = PetscSFReduceBegin_Basic;
1207:   sf->ops->ReduceEnd       = PetscSFReduceEnd_Basic;
1208:   sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Basic;
1209:   sf->ops->FetchAndOpEnd   = PetscSFFetchAndOpEnd_Basic;
1210:   sf->ops->GetLeafRanks    = PetscSFGetLeafRanks_Basic;

1212:   PetscNewLog(sf,&bas);
1213:   sf->data = (void*)bas;
1214:   return(0);
1215: }