xref: /petsc/src/vec/is/sf/impls/basic/sfbasic.c (revision b5a8e515b429ceb8098d7b9fc29a8331d9ac04b7)
18cd53115SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/
395fce210SBarry Smith 
495fce210SBarry Smith typedef struct _n_PetscSFBasicPack *PetscSFBasicPack;
595fce210SBarry Smith struct _n_PetscSFBasicPack {
62a74e5f4SJed Brown   void (*Pack)(PetscInt,PetscInt,const PetscInt*,const void*,void*);
72a74e5f4SJed Brown   void (*UnpackInsert)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
82a74e5f4SJed Brown   void (*UnpackAdd)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
92a74e5f4SJed Brown   void (*UnpackMin)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
102a74e5f4SJed Brown   void (*UnpackMax)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
112a74e5f4SJed Brown   void (*UnpackMinloc)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
122a74e5f4SJed Brown   void (*UnpackMaxloc)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
132a74e5f4SJed Brown   void (*UnpackMult)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
142a74e5f4SJed Brown   void (*UnpackLAND)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
152a74e5f4SJed Brown   void (*UnpackBAND)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
162a74e5f4SJed Brown   void (*UnpackLOR)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
172a74e5f4SJed Brown   void (*UnpackBOR)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
182a74e5f4SJed Brown   void (*UnpackLXOR)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
192a74e5f4SJed Brown   void (*UnpackBXOR)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
202a74e5f4SJed Brown   void (*FetchAndInsert)(PetscInt,PetscInt,const PetscInt*,void*,void*);
212a74e5f4SJed Brown   void (*FetchAndAdd)(PetscInt,PetscInt,const PetscInt*,void*,void*);
222a74e5f4SJed Brown   void (*FetchAndMin)(PetscInt,PetscInt,const PetscInt*,void*,void*);
232a74e5f4SJed Brown   void (*FetchAndMax)(PetscInt,PetscInt,const PetscInt*,void*,void*);
242a74e5f4SJed Brown   void (*FetchAndMinloc)(PetscInt,PetscInt,const PetscInt*,void*,void*);
252a74e5f4SJed Brown   void (*FetchAndMaxloc)(PetscInt,PetscInt,const PetscInt*,void*,void*);
262a74e5f4SJed Brown   void (*FetchAndMult)(PetscInt,PetscInt,const PetscInt*,void*,void*);
272a74e5f4SJed Brown   void (*FetchAndLAND)(PetscInt,PetscInt,const PetscInt*,void*,void*);
282a74e5f4SJed Brown   void (*FetchAndBAND)(PetscInt,PetscInt,const PetscInt*,void*,void*);
292a74e5f4SJed Brown   void (*FetchAndLOR)(PetscInt,PetscInt,const PetscInt*,void*,void*);
302a74e5f4SJed Brown   void (*FetchAndBOR)(PetscInt,PetscInt,const PetscInt*,void*,void*);
312a74e5f4SJed Brown   void (*FetchAndLXOR)(PetscInt,PetscInt,const PetscInt*,void*,void*);
322a74e5f4SJed Brown   void (*FetchAndBXOR)(PetscInt,PetscInt,const PetscInt*,void*,void*);
3395fce210SBarry Smith 
3495fce210SBarry Smith   MPI_Datatype     unit;
3595fce210SBarry Smith   size_t           unitbytes;   /* Number of bytes in a unit */
362a74e5f4SJed Brown   PetscInt         bs;          /* Number of basic units in a unit */
3795fce210SBarry Smith   const void       *key;        /* Array used as key for operation */
38bf39f1bfSJed Brown   char             **root;      /* Packed root data, indexed by leaf rank */
39bf39f1bfSJed Brown   char             **leaf;      /* Packed leaf data, indexed by root rank */
4095fce210SBarry Smith   MPI_Request      *requests;   /* Array of root requests followed by leaf requests */
4195fce210SBarry Smith   PetscSFBasicPack next;
4295fce210SBarry Smith };
4395fce210SBarry Smith 
4495fce210SBarry Smith typedef struct {
4595fce210SBarry Smith   PetscMPIInt      tag;
46cf4b5b4fSJed Brown   PetscMPIInt      niranks;     /* Number of incoming ranks (ranks accessing my roots) */
4795fce210SBarry Smith   PetscMPIInt      *iranks;     /* Array of ranks that reference my roots */
4895fce210SBarry Smith   PetscInt         itotal;      /* Total number of graph edges referencing my roots */
4995fce210SBarry Smith   PetscInt         *ioffset;    /* Array of length niranks+1 holding offset in irootloc[] for each rank */
5095fce210SBarry Smith   PetscInt         *irootloc;   /* Incoming roots referenced by ranks starting at ioffset[rank] */
5195fce210SBarry Smith   PetscSFBasicPack avail;       /* One or more entries per MPI Datatype, lazily constructed */
5295fce210SBarry Smith   PetscSFBasicPack inuse;       /* Buffers being used for transactions that have not yet completed */
5395fce210SBarry Smith } PetscSF_Basic;
5495fce210SBarry Smith 
5595fce210SBarry Smith #if !defined(PETSC_HAVE_MPI_TYPE_DUP) /* Danger: type is not reference counted; subject to ABA problem */
5695fce210SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MPI_Type_dup(MPI_Datatype datatype,MPI_Datatype *newtype)
5795fce210SBarry Smith {
5895fce210SBarry Smith   *newtype = datatype;
5995fce210SBarry Smith   return 0;
6095fce210SBarry Smith }
6195fce210SBarry Smith #endif
6295fce210SBarry Smith 
6395fce210SBarry Smith /*
6495fce210SBarry Smith  * MPI_Reduce_local is not really useful because it can't handle sparse data and it vectorizes "in the wrong direction",
6595fce210SBarry Smith  * therefore we pack data types manually. This section defines packing routines for the standard data types.
6695fce210SBarry Smith  */
6795fce210SBarry Smith 
6895fce210SBarry Smith #define CPPJoin2_exp(a,b) a ## b
6995fce210SBarry Smith #define CPPJoin2(a,b) CPPJoin2_exp(a,b)
7095fce210SBarry Smith #define CPPJoin3_exp_(a,b,c) a ## b ## _ ## c
7195fce210SBarry Smith #define CPPJoin3_(a,b,c) CPPJoin3_exp_(a,b,c)
7295fce210SBarry Smith 
7395fce210SBarry Smith /* Basic types without addition */
742a74e5f4SJed Brown #define DEF_PackNoInit(type,BS)                                         \
752a74e5f4SJed Brown   static void CPPJoin3_(Pack_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,const void *unpacked,void *packed) { \
7695fce210SBarry Smith     const type *u = (const type*)unpacked;                              \
7795fce210SBarry Smith     type *p = (type*)packed;                                            \
782a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
792a74e5f4SJed Brown     for (i=0; i<n; i++)                                                 \
802a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
812a74e5f4SJed Brown         for (k=j; k<j+BS; k++)                                          \
822a74e5f4SJed Brown           p[i*bs+k] = u[idx[i]*bs+k];                                   \
8395fce210SBarry Smith   }                                                                     \
842a74e5f4SJed Brown   static void CPPJoin3_(UnpackInsert_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
8595fce210SBarry Smith     type *u = (type*)unpacked;                                          \
8695fce210SBarry Smith     const type *p = (const type*)packed;                                \
872a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
882a74e5f4SJed Brown     for (i=0; i<n; i++)                                                 \
892a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
902a74e5f4SJed Brown         for (k=j; k<j+BS; k++)                                          \
912a74e5f4SJed Brown           u[idx[i]*bs+k] = p[i*bs+k];                                   \
9295fce210SBarry Smith   }                                                                     \
932a74e5f4SJed Brown   static void CPPJoin3_(FetchAndInsert_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
9495fce210SBarry Smith     type *u = (type*)unpacked;                                          \
9595fce210SBarry Smith     type *p = (type*)packed;                                            \
962a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
9795fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
982a74e5f4SJed Brown       PetscInt ii = idx[i];                                             \
992a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
1002a74e5f4SJed Brown         for (k=j; k<j+BS; k++) {                                        \
1012a74e5f4SJed Brown           type t = u[ii*bs+k];                                          \
1022a74e5f4SJed Brown           u[ii*bs+k] = p[i*bs+k];                                       \
1032a74e5f4SJed Brown           p[i*bs+k] = t;                                                \
1042a74e5f4SJed Brown         }                                                               \
10595fce210SBarry Smith     }                                                                   \
10695fce210SBarry Smith   }
10795fce210SBarry Smith 
10895fce210SBarry Smith /* Basic types defining addition */
1092a74e5f4SJed Brown #define DEF_PackAddNoInit(type,BS)                                      \
1102a74e5f4SJed Brown   DEF_PackNoInit(type,BS)                                               \
1112a74e5f4SJed Brown   static void CPPJoin3_(UnpackAdd_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
11295fce210SBarry Smith     type *u = (type*)unpacked;                                          \
11395fce210SBarry Smith     const type *p = (const type*)packed;                                \
1142a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
1152a74e5f4SJed Brown     for (i=0; i<n; i++)                                                 \
1162a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
1172a74e5f4SJed Brown         for (k=j; k<j+BS; k++)                                          \
1182a74e5f4SJed Brown           u[idx[i]*bs+k] += p[i*bs+k];                                  \
11995fce210SBarry Smith   }                                                                     \
1202a74e5f4SJed Brown   static void CPPJoin3_(FetchAndAdd_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
12195fce210SBarry Smith     type *u = (type*)unpacked;                                          \
12295fce210SBarry Smith     type *p = (type*)packed;                                            \
1232a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
12495fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
1252a74e5f4SJed Brown       PetscInt ii = idx[i];                                             \
1262a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
1272a74e5f4SJed Brown         for (k=j; k<j+BS; k++) {                                        \
1282a74e5f4SJed Brown           type t = u[ii*bs+k];                                          \
1292a74e5f4SJed Brown           u[ii*bs+k] = t + p[i*bs+k];                                   \
1302a74e5f4SJed Brown           p[i*bs+k] = t;                                                \
13195fce210SBarry Smith         }                                                               \
13253deab39SPeter Brune     }                                                                   \
1332a74e5f4SJed Brown   }                                                                     \
1342a74e5f4SJed Brown   static void CPPJoin3_(UnpackMult_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
13553deab39SPeter Brune     type *u = (type*)unpacked;                                          \
13653deab39SPeter Brune     const type *p = (const type*)packed;                                \
1372a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
1382a74e5f4SJed Brown     for (i=0; i<n; i++)                                                 \
1392a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
1402a74e5f4SJed Brown         for (k=j; k<j+BS; k++)                                          \
1412a74e5f4SJed Brown           u[idx[i]*bs+k] *= p[i*bs+k];                                  \
14253deab39SPeter Brune   }                                                                     \
1432a74e5f4SJed Brown   static void CPPJoin3_(FetchAndMult_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
14453deab39SPeter Brune     type *u = (type*)unpacked;                                          \
14553deab39SPeter Brune     type *p = (type*)packed;                                            \
1462a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
14753deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
1482a74e5f4SJed Brown       PetscInt ii = idx[i];                                             \
1492a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
1502a74e5f4SJed Brown         for (k=j; k<j+BS; k++) {                                        \
1512a74e5f4SJed Brown           type t = u[ii*bs+k];                                          \
1522a74e5f4SJed Brown           u[ii*bs+k] = t * p[i*bs+k];                                   \
1532a74e5f4SJed Brown           p[i*bs+k] = t;                                                \
1542a74e5f4SJed Brown         }                                                               \
15553deab39SPeter Brune     }                                                                   \
15695fce210SBarry Smith   }
1572a74e5f4SJed Brown #define DEF_Pack(type,BS)                                               \
1582a74e5f4SJed Brown   DEF_PackAddNoInit(type,BS)                                            \
1592a74e5f4SJed Brown   static void CPPJoin3_(PackInit_,type,BS)(PetscSFBasicPack link) {     \
1602a74e5f4SJed Brown     link->Pack = CPPJoin3_(Pack_,type,BS);                              \
1612a74e5f4SJed Brown     link->UnpackInsert = CPPJoin3_(UnpackInsert_,type,BS);              \
1622a74e5f4SJed Brown     link->UnpackAdd = CPPJoin3_(UnpackAdd_,type,BS);                    \
1632a74e5f4SJed Brown     link->UnpackMult = CPPJoin3_(UnpackMult_,type,BS);                  \
1642a74e5f4SJed Brown     link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,type,BS);          \
1652a74e5f4SJed Brown     link->FetchAndAdd = CPPJoin3_(FetchAndAdd_,type,BS);                \
1662a74e5f4SJed Brown     link->FetchAndMult = CPPJoin3_(FetchAndMult_,type,BS);              \
16795fce210SBarry Smith     link->unitbytes = sizeof(type);                                     \
16895fce210SBarry Smith   }
16995fce210SBarry Smith /* Comparable types */
17095fce210SBarry Smith #define DEF_PackCmp(type)                                               \
1712a74e5f4SJed Brown   DEF_PackAddNoInit(type,1)                                             \
1722a74e5f4SJed Brown   static void CPPJoin2(UnpackMax_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
17395fce210SBarry Smith     type *u = (type*)unpacked;                                          \
17495fce210SBarry Smith     const type *p = (const type*)packed;                                \
17595fce210SBarry Smith     PetscInt i;                                                         \
17695fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
17795fce210SBarry Smith       type v = u[idx[i]];                                               \
17895fce210SBarry Smith       u[idx[i]] = PetscMax(v,p[i]);                                     \
17995fce210SBarry Smith     }                                                                   \
18095fce210SBarry Smith   }                                                                     \
1812a74e5f4SJed Brown   static void CPPJoin2(UnpackMin_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
18295fce210SBarry Smith     type *u = (type*)unpacked;                                          \
18395fce210SBarry Smith     const type *p = (const type*)packed;                                \
18495fce210SBarry Smith     PetscInt i;                                                         \
18595fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
18695fce210SBarry Smith       type v = u[idx[i]];                                               \
18795fce210SBarry Smith       u[idx[i]] = PetscMin(v,p[i]);                                     \
18895fce210SBarry Smith     }                                                                   \
18995fce210SBarry Smith   }                                                                     \
1902a74e5f4SJed Brown   static void CPPJoin2(FetchAndMax_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
19195fce210SBarry Smith     type *u = (type*)unpacked;                                          \
19295fce210SBarry Smith     type *p = (type*)packed;                                            \
19395fce210SBarry Smith     PetscInt i;                                                         \
19495fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
19595fce210SBarry Smith       PetscInt j = idx[i];                                              \
19695fce210SBarry Smith       type v = u[j];                                                    \
19795fce210SBarry Smith       u[j] = PetscMax(v,p[i]);                                          \
19895fce210SBarry Smith       p[i] = v;                                                         \
19995fce210SBarry Smith     }                                                                   \
20095fce210SBarry Smith   }                                                                     \
2012a74e5f4SJed Brown   static void CPPJoin2(FetchAndMin_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
20295fce210SBarry Smith     type *u = (type*)unpacked;                                          \
20395fce210SBarry Smith     type *p = (type*)packed;                                            \
20495fce210SBarry Smith     PetscInt i;                                                         \
20595fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
20695fce210SBarry Smith       PetscInt j = idx[i];                                              \
20795fce210SBarry Smith       type v = u[j];                                                    \
20895fce210SBarry Smith       u[j] = PetscMin(v,p[i]);                                          \
20995fce210SBarry Smith       p[i] = v;                                                         \
21095fce210SBarry Smith     }                                                                   \
21195fce210SBarry Smith   }                                                                     \
21295fce210SBarry Smith   static void CPPJoin2(PackInit_,type)(PetscSFBasicPack link) {         \
2132a74e5f4SJed Brown     link->Pack = CPPJoin3_(Pack_,type,1);                               \
2142a74e5f4SJed Brown     link->UnpackInsert = CPPJoin3_(UnpackInsert_,type,1);               \
2152a74e5f4SJed Brown     link->UnpackAdd  = CPPJoin3_(UnpackAdd_,type,1);                    \
21695fce210SBarry Smith     link->UnpackMax  = CPPJoin2(UnpackMax_,type);                       \
21795fce210SBarry Smith     link->UnpackMin  = CPPJoin2(UnpackMin_,type);                       \
2182a74e5f4SJed Brown     link->UnpackMult = CPPJoin3_(UnpackMult_,type,1);                   \
2192a74e5f4SJed Brown     link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,type,1);           \
2202a74e5f4SJed Brown     link->FetchAndAdd = CPPJoin3_(FetchAndAdd_ ,type,1);                \
22195fce210SBarry Smith     link->FetchAndMax = CPPJoin2(FetchAndMax_ ,type);                   \
22295fce210SBarry Smith     link->FetchAndMin = CPPJoin2(FetchAndMin_ ,type);                   \
2232a74e5f4SJed Brown     link->FetchAndMult = CPPJoin3_(FetchAndMult_,type,1);               \
22495fce210SBarry Smith     link->unitbytes = sizeof(type);                                     \
22595fce210SBarry Smith   }
22695fce210SBarry Smith 
22753deab39SPeter Brune /* Logical Types */
22853deab39SPeter Brune #define DEF_PackLog(type)                                               \
2292a74e5f4SJed Brown   static void CPPJoin2(UnpackLAND_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
23053deab39SPeter Brune     type *u = (type*)unpacked;                                          \
23153deab39SPeter Brune     const type *p = (const type*)packed;                                \
23253deab39SPeter Brune     PetscInt i;                                                         \
23353deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
23453deab39SPeter Brune       type v = u[idx[i]];                                               \
23553deab39SPeter Brune       u[idx[i]] = v && p[i];                                            \
23653deab39SPeter Brune     }                                                                   \
23753deab39SPeter Brune   }                                                                     \
2382a74e5f4SJed Brown   static void CPPJoin2(UnpackLOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
23953deab39SPeter Brune     type *u = (type*)unpacked;                                          \
24053deab39SPeter Brune     const type *p = (const type*)packed;                                \
24153deab39SPeter Brune     PetscInt i;                                                         \
24253deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
24353deab39SPeter Brune       type v = u[idx[i]];                                               \
24453deab39SPeter Brune       u[idx[i]] = v || p[i];                                            \
24553deab39SPeter Brune     }                                                                   \
24653deab39SPeter Brune   }                                                                     \
2472a74e5f4SJed Brown   static void CPPJoin2(UnpackLXOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
24853deab39SPeter Brune     type *u = (type*)unpacked;                                          \
24953deab39SPeter Brune     const type *p = (const type*)packed;                                \
25053deab39SPeter Brune     PetscInt i;                                                         \
25153deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
25253deab39SPeter Brune       type v = u[idx[i]];                                               \
25353deab39SPeter Brune       u[idx[i]] = (!v)!=(!p[i]);                                        \
25453deab39SPeter Brune     }                                                                   \
25553deab39SPeter Brune   }                                                                     \
2562a74e5f4SJed Brown   static void CPPJoin2(FetchAndLAND_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
25753deab39SPeter Brune     type *u = (type*)unpacked;                                          \
25853deab39SPeter Brune     type *p = (type*)packed;                                            \
25953deab39SPeter Brune     PetscInt i;                                                         \
26053deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
26153deab39SPeter Brune       PetscInt j = idx[i];                                              \
26253deab39SPeter Brune       type v = u[j];                                                    \
26353deab39SPeter Brune       u[j] = v && p[i];                                                 \
26453deab39SPeter Brune       p[i] = v;                                                         \
26553deab39SPeter Brune     }                                                                   \
26653deab39SPeter Brune   }                                                                     \
2672a74e5f4SJed Brown   static void CPPJoin2(FetchAndLOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
26853deab39SPeter Brune     type *u = (type*)unpacked;                                          \
26953deab39SPeter Brune     type *p = (type*)packed;                                            \
27053deab39SPeter Brune     PetscInt i;                                                         \
27153deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
27253deab39SPeter Brune       PetscInt j = idx[i];                                              \
27353deab39SPeter Brune       type v = u[j];                                                    \
27453deab39SPeter Brune       u[j] = v || p[i];                                                 \
27553deab39SPeter Brune       p[i] = v;                                                         \
27653deab39SPeter Brune     }                                                                   \
27753deab39SPeter Brune   }                                                                     \
2782a74e5f4SJed Brown   static void CPPJoin2(FetchAndLXOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
27953deab39SPeter Brune     type *u = (type*)unpacked;                                          \
28053deab39SPeter Brune     type *p = (type*)packed;                                            \
28153deab39SPeter Brune     PetscInt i;                                                         \
28253deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
28353deab39SPeter Brune       PetscInt j = idx[i];                                              \
28453deab39SPeter Brune       type v = u[j];                                                    \
28553deab39SPeter Brune       u[j] = (!v)!=(!p[i]);                                             \
28653deab39SPeter Brune       p[i] = v;                                                         \
28753deab39SPeter Brune     }                                                                   \
28853deab39SPeter Brune   }                                                                     \
28953deab39SPeter Brune   static void CPPJoin2(PackInit_Logical_,type)(PetscSFBasicPack link) { \
29053deab39SPeter Brune     link->UnpackLAND = CPPJoin2(UnpackLAND_,type);                      \
29153deab39SPeter Brune     link->UnpackLOR  = CPPJoin2(UnpackLOR_,type);                       \
29253deab39SPeter Brune     link->UnpackLXOR = CPPJoin2(UnpackLXOR_,type);                      \
29353deab39SPeter Brune     link->FetchAndLAND = CPPJoin2(FetchAndLAND_,type);                  \
29453deab39SPeter Brune     link->FetchAndLOR  = CPPJoin2(FetchAndLOR_,type);                   \
29553deab39SPeter Brune     link->FetchAndLXOR = CPPJoin2(FetchAndLXOR_,type);                  \
29653deab39SPeter Brune   }
29753deab39SPeter Brune 
29853deab39SPeter Brune 
29953deab39SPeter Brune /* Bitwise Types */
30053deab39SPeter Brune #define DEF_PackBit(type)                                               \
3012a74e5f4SJed Brown   static void CPPJoin2(UnpackBAND_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
30253deab39SPeter Brune     type *u = (type*)unpacked;                                          \
30353deab39SPeter Brune     const type *p = (const type*)packed;                                \
30453deab39SPeter Brune     PetscInt i;                                                         \
30553deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
30653deab39SPeter Brune       type v = u[idx[i]];                                               \
30753deab39SPeter Brune       u[idx[i]] = v & p[i];                                             \
30853deab39SPeter Brune     }                                                                   \
30953deab39SPeter Brune   }                                                                     \
3102a74e5f4SJed Brown   static void CPPJoin2(UnpackBOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
31153deab39SPeter Brune     type *u = (type*)unpacked;                                          \
31253deab39SPeter Brune     const type *p = (const type*)packed;                                \
31353deab39SPeter Brune     PetscInt i;                                                         \
31453deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
31553deab39SPeter Brune       type v = u[idx[i]];                                               \
31653deab39SPeter Brune       u[idx[i]] = v | p[i];                                             \
31753deab39SPeter Brune     }                                                                   \
31853deab39SPeter Brune   }                                                                     \
3192a74e5f4SJed Brown   static void CPPJoin2(UnpackBXOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
32053deab39SPeter Brune     type *u = (type*)unpacked;                                          \
32153deab39SPeter Brune     const type *p = (const type*)packed;                                \
32253deab39SPeter Brune     PetscInt i;                                                         \
32353deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
32453deab39SPeter Brune       type v = u[idx[i]];                                               \
32553deab39SPeter Brune       u[idx[i]] = v^p[i];                                               \
32653deab39SPeter Brune     }                                                                   \
32753deab39SPeter Brune   }                                                                     \
3282a74e5f4SJed Brown   static void CPPJoin2(FetchAndBAND_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
32953deab39SPeter Brune     type *u = (type*)unpacked;                                          \
33053deab39SPeter Brune     type *p = (type*)packed;                                            \
33153deab39SPeter Brune     PetscInt i;                                                         \
33253deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
33353deab39SPeter Brune       PetscInt j = idx[i];                                              \
33453deab39SPeter Brune       type v = u[j];                                                    \
33553deab39SPeter Brune       u[j] = v & p[i];                                                  \
33653deab39SPeter Brune       p[i] = v;                                                         \
33753deab39SPeter Brune     }                                                                   \
33853deab39SPeter Brune   }                                                                     \
3392a74e5f4SJed Brown   static void CPPJoin2(FetchAndBOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
34053deab39SPeter Brune     type *u = (type*)unpacked;                                          \
34153deab39SPeter Brune     type *p = (type*)packed;                                            \
34253deab39SPeter Brune     PetscInt i;                                                         \
34353deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
34453deab39SPeter Brune       PetscInt j = idx[i];                                              \
34553deab39SPeter Brune       type v = u[j];                                                    \
34653deab39SPeter Brune       u[j] = v | p[i];                                                  \
34753deab39SPeter Brune       p[i] = v;                                                         \
34853deab39SPeter Brune     }                                                                   \
34953deab39SPeter Brune   }                                                                     \
3502a74e5f4SJed Brown   static void CPPJoin2(FetchAndBXOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
35153deab39SPeter Brune     type *u = (type*)unpacked;                                          \
35253deab39SPeter Brune     type *p = (type*)packed;                                            \
35353deab39SPeter Brune     PetscInt i;                                                         \
35453deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
35553deab39SPeter Brune       PetscInt j = idx[i];                                              \
35653deab39SPeter Brune       type v = u[j];                                                    \
35753deab39SPeter Brune       u[j] = v^p[i];                                                    \
35853deab39SPeter Brune       p[i] = v;                                                         \
35953deab39SPeter Brune     }                                                                   \
36053deab39SPeter Brune   }                                                                     \
36153deab39SPeter Brune   static void CPPJoin2(PackInit_Bitwise_,type)(PetscSFBasicPack link) { \
36253deab39SPeter Brune     link->UnpackBAND = CPPJoin2(UnpackBAND_,type);                      \
36353deab39SPeter Brune     link->UnpackBOR  = CPPJoin2(UnpackBOR_,type);                       \
36453deab39SPeter Brune     link->UnpackBXOR = CPPJoin2(UnpackBXOR_,type);                      \
36553deab39SPeter Brune     link->FetchAndBAND = CPPJoin2(FetchAndBAND_,type);                  \
36653deab39SPeter Brune     link->FetchAndBOR  = CPPJoin2(FetchAndBOR_,type);                   \
36753deab39SPeter Brune     link->FetchAndBXOR = CPPJoin2(FetchAndBXOR_,type);                  \
36853deab39SPeter Brune   }
36953deab39SPeter Brune 
37095fce210SBarry Smith /* Pair types */
37195fce210SBarry Smith #define CPPJoinloc_exp(base,op,t1,t2) base ## op ## loc_ ## t1 ## _ ## t2
37295fce210SBarry Smith #define CPPJoinloc(base,op,t1,t2) CPPJoinloc_exp(base,op,t1,t2)
37395fce210SBarry Smith #define PairType(type1,type2) CPPJoin3_(_pairtype_,type1,type2)
37495fce210SBarry Smith #define DEF_UnpackXloc(type1,type2,locname,op)                              \
3752a74e5f4SJed Brown   static void CPPJoinloc(Unpack,locname,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
37695fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;        \
37795fce210SBarry Smith     const PairType(type1,type2) *p = (const PairType(type1,type2)*)packed; \
37895fce210SBarry Smith     PetscInt i;                                                         \
37995fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
38095fce210SBarry Smith       PetscInt j = idx[i];                                              \
38195fce210SBarry Smith       if (p[i].a op u[j].a) {                                           \
38295fce210SBarry Smith         u[j].a = p[i].a;                                                \
38395fce210SBarry Smith         u[j].b = p[i].b;                                                \
38495fce210SBarry Smith       } else if (u[j].a == p[i].a) {                                    \
38595fce210SBarry Smith         u[j].b = PetscMin(u[j].b,p[i].b);                               \
38695fce210SBarry Smith       }                                                                 \
38795fce210SBarry Smith     }                                                                   \
38895fce210SBarry Smith   }                                                                     \
3892a74e5f4SJed Brown   static void CPPJoinloc(FetchAnd,locname,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
39095fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;        \
39195fce210SBarry Smith     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;          \
39295fce210SBarry Smith     PetscInt i;                                                         \
39395fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
39495fce210SBarry Smith       PetscInt j = idx[i];                                              \
39595fce210SBarry Smith       PairType(type1,type2) v;                                          \
39695fce210SBarry Smith       v.a = u[j].a;                                                     \
39795fce210SBarry Smith       v.b = u[j].b;                                                     \
39895fce210SBarry Smith       if (p[i].a op u[j].a) {                                           \
39995fce210SBarry Smith         u[j].a = p[i].a;                                                \
40095fce210SBarry Smith         u[j].b = p[i].b;                                                \
40195fce210SBarry Smith       } else if (u[j].a == p[i].a) {                                    \
40295fce210SBarry Smith         u[j].b = PetscMin(u[j].b,p[i].b);                               \
40395fce210SBarry Smith       }                                                                 \
40495fce210SBarry Smith       p[i].a = v.a;                                                     \
40595fce210SBarry Smith       p[i].b = v.b;                                                     \
40695fce210SBarry Smith     }                                                                   \
40795fce210SBarry Smith   }
40895fce210SBarry Smith #define DEF_PackPair(type1,type2)                                       \
40995fce210SBarry Smith   typedef struct {type1 a; type2 b;} PairType(type1,type2);             \
4102a74e5f4SJed Brown   static void CPPJoin3_(Pack_,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,const void *unpacked,void *packed) { \
41195fce210SBarry Smith     const PairType(type1,type2) *u = (const PairType(type1,type2)*)unpacked; \
41295fce210SBarry Smith     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;          \
41395fce210SBarry Smith     PetscInt i;                                                         \
41495fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
41595fce210SBarry Smith       p[i].a = u[idx[i]].a;                                             \
41695fce210SBarry Smith       p[i].b = u[idx[i]].b;                                             \
41795fce210SBarry Smith     }                                                                   \
41895fce210SBarry Smith   }                                                                     \
4192a74e5f4SJed Brown   static void CPPJoin3_(UnpackInsert_,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
42095fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;       \
42195fce210SBarry Smith     const PairType(type1,type2) *p = (const PairType(type1,type2)*)packed; \
42295fce210SBarry Smith     PetscInt i;                                                         \
42395fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
42495fce210SBarry Smith       u[idx[i]].a = p[i].a;                                             \
42595fce210SBarry Smith       u[idx[i]].b = p[i].b;                                             \
42695fce210SBarry Smith     }                                                                   \
42795fce210SBarry Smith   }                                                                     \
4282a74e5f4SJed Brown   static void CPPJoin3_(UnpackAdd_,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
42995fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;       \
43095fce210SBarry Smith     const PairType(type1,type2) *p = (const PairType(type1,type2)*)packed; \
43195fce210SBarry Smith     PetscInt i;                                                         \
43295fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
43395fce210SBarry Smith       u[idx[i]].a += p[i].a;                                            \
43495fce210SBarry Smith       u[idx[i]].b += p[i].b;                                            \
43595fce210SBarry Smith     }                                                                   \
43695fce210SBarry Smith   }                                                                     \
4372a74e5f4SJed Brown   static void CPPJoin3_(FetchAndInsert_,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
43895fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;        \
43995fce210SBarry Smith     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;          \
44095fce210SBarry Smith     PetscInt i;                                                         \
44195fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
44295fce210SBarry Smith       PetscInt j = idx[i];                                              \
44395fce210SBarry Smith       PairType(type1,type2) v;                                          \
44495fce210SBarry Smith       v.a = u[j].a;                                                     \
44595fce210SBarry Smith       v.b = u[j].b;                                                     \
44695fce210SBarry Smith       u[j].a = p[i].a;                                                  \
44795fce210SBarry Smith       u[j].b = p[i].b;                                                  \
44895fce210SBarry Smith       p[i].a = v.a;                                                     \
44995fce210SBarry Smith       p[i].b = v.b;                                                     \
45095fce210SBarry Smith     }                                                                   \
45195fce210SBarry Smith   }                                                                     \
4522a74e5f4SJed Brown   static void FetchAndAdd_ ## type1 ## _ ## type2(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
45395fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;       \
45495fce210SBarry Smith     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;         \
45595fce210SBarry Smith     PetscInt i;                                                         \
45695fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
45795fce210SBarry Smith       PetscInt j = idx[i];                                              \
45895fce210SBarry Smith       PairType(type1,type2) v;                                          \
45995fce210SBarry Smith       v.a = u[j].a;                                                     \
46095fce210SBarry Smith       v.b = u[j].b;                                                     \
46195fce210SBarry Smith       u[j].a = v.a + p[i].a;                                            \
46295fce210SBarry Smith       u[j].b = v.b + p[i].b;                                            \
46395fce210SBarry Smith       p[i].a = v.a;                                                     \
46495fce210SBarry Smith       p[i].b = v.b;                                                     \
46595fce210SBarry Smith     }                                                                   \
46695fce210SBarry Smith   }                                                                     \
46795fce210SBarry Smith   DEF_UnpackXloc(type1,type2,Max,>)                                     \
46895fce210SBarry Smith   DEF_UnpackXloc(type1,type2,Min,<)                                     \
46995fce210SBarry Smith   static void CPPJoin3_(PackInit_,type1,type2)(PetscSFBasicPack link) { \
47095fce210SBarry Smith     link->Pack = CPPJoin3_(Pack_,type1,type2);                          \
47195fce210SBarry Smith     link->UnpackInsert = CPPJoin3_(UnpackInsert_,type1,type2);          \
47295fce210SBarry Smith     link->UnpackAdd = CPPJoin3_(UnpackAdd_,type1,type2);                \
47395fce210SBarry Smith     link->UnpackMaxloc = CPPJoin3_(UnpackMaxloc_,type1,type2);          \
47495fce210SBarry Smith     link->UnpackMinloc = CPPJoin3_(UnpackMinloc_,type1,type2);          \
47595fce210SBarry Smith     link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,type1,type2);      \
47695fce210SBarry Smith     link->FetchAndAdd = CPPJoin3_(FetchAndAdd_,type1,type2);            \
47795fce210SBarry Smith     link->FetchAndMaxloc = CPPJoin3_(FetchAndMaxloc_,type1,type2);      \
47895fce210SBarry Smith     link->FetchAndMinloc = CPPJoin3_(FetchAndMinloc_,type1,type2);      \
47995fce210SBarry Smith     link->unitbytes = sizeof(PairType(type1,type2));                    \
48095fce210SBarry Smith   }
48195fce210SBarry Smith 
48295fce210SBarry Smith /* Currently only dumb blocks of data */
48395fce210SBarry Smith #define BlockType(unit,count) CPPJoin3_(_blocktype_,unit,count)
48495fce210SBarry Smith #define DEF_Block(unit,count)                                           \
48595fce210SBarry Smith   typedef struct {unit v[count];} BlockType(unit,count);                \
4862a74e5f4SJed Brown   DEF_PackNoInit(BlockType(unit,count),1)                               \
48795fce210SBarry Smith   static void CPPJoin3_(PackInit_block_,unit,count)(PetscSFBasicPack link) { \
4882a74e5f4SJed Brown     link->Pack = CPPJoin3_(Pack_,BlockType(unit,count),1);               \
4892a74e5f4SJed Brown     link->UnpackInsert = CPPJoin3_(UnpackInsert_,BlockType(unit,count),1); \
4902a74e5f4SJed Brown     link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,BlockType(unit,count),1); \
49195fce210SBarry Smith     link->unitbytes = sizeof(BlockType(unit,count));                    \
49295fce210SBarry Smith   }
49395fce210SBarry Smith 
49495fce210SBarry Smith DEF_PackCmp(int)
49553deab39SPeter Brune DEF_PackBit(int)
49653deab39SPeter Brune DEF_PackLog(int)
49795fce210SBarry Smith DEF_PackCmp(PetscInt)
49853deab39SPeter Brune DEF_PackBit(PetscInt)
49953deab39SPeter Brune DEF_PackLog(PetscInt)
5002a74e5f4SJed Brown DEF_Pack(PetscInt,2)
5012a74e5f4SJed Brown DEF_Pack(PetscInt,3)
5022a74e5f4SJed Brown DEF_Pack(PetscInt,4)
5032a74e5f4SJed Brown DEF_Pack(PetscInt,5)
5042a74e5f4SJed Brown DEF_Pack(PetscInt,7)
50595fce210SBarry Smith DEF_PackCmp(PetscReal)
50653deab39SPeter Brune DEF_PackLog(PetscReal)
5072a74e5f4SJed Brown DEF_Pack(PetscReal,2)
5082a74e5f4SJed Brown DEF_Pack(PetscReal,3)
5092a74e5f4SJed Brown DEF_Pack(PetscReal,4)
5102a74e5f4SJed Brown DEF_Pack(PetscReal,5)
5112a74e5f4SJed Brown DEF_Pack(PetscReal,7)
51295fce210SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
5132a74e5f4SJed Brown DEF_Pack(PetscComplex,1)
5142a74e5f4SJed Brown DEF_Pack(PetscComplex,2)
5152a74e5f4SJed Brown DEF_Pack(PetscComplex,3)
5162a74e5f4SJed Brown DEF_Pack(PetscComplex,4)
5172a74e5f4SJed Brown DEF_Pack(PetscComplex,5)
5182a74e5f4SJed Brown DEF_Pack(PetscComplex,7)
51995fce210SBarry Smith #endif
52095fce210SBarry Smith DEF_PackPair(int,int)
52195fce210SBarry Smith DEF_PackPair(PetscInt,PetscInt)
52295fce210SBarry Smith DEF_Block(int,1)
52395fce210SBarry Smith DEF_Block(int,2)
52495fce210SBarry Smith DEF_Block(int,3)
52595fce210SBarry Smith DEF_Block(int,4)
52695fce210SBarry Smith DEF_Block(int,5)
52795fce210SBarry Smith DEF_Block(int,6)
52895fce210SBarry Smith DEF_Block(int,7)
52995fce210SBarry Smith DEF_Block(int,8)
53095fce210SBarry Smith 
53195fce210SBarry Smith static PetscErrorCode PetscSFSetUp_Basic(PetscSF sf)
53295fce210SBarry Smith {
53395fce210SBarry Smith   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
53495fce210SBarry Smith   PetscErrorCode ierr;
53595fce210SBarry Smith   PetscInt *rlengths,*ilengths,i;
53695fce210SBarry Smith   MPI_Comm comm;
537*b5a8e515SJed Brown   MPI_Group group;
53895fce210SBarry Smith   MPI_Request *rootreqs,*leafreqs;
53995fce210SBarry Smith 
54095fce210SBarry Smith   PetscFunctionBegin;
541*b5a8e515SJed Brown   ierr = MPI_Comm_group(PETSC_COMM_SELF,&group);CHKERRQ(ierr);
542*b5a8e515SJed Brown   ierr = PetscSFSetUpRanks(sf,group);CHKERRQ(ierr);
543*b5a8e515SJed Brown   ierr = MPI_Group_free(&group);CHKERRQ(ierr);
54495fce210SBarry Smith   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
54595fce210SBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)sf,&bas->tag);CHKERRQ(ierr);
54695fce210SBarry Smith   /*
54795fce210SBarry Smith    * Inform roots about how many leaves and from which ranks
54895fce210SBarry Smith    */
549785e854fSJed Brown   ierr = PetscMalloc1(sf->nranks,&rlengths);CHKERRQ(ierr);
55095fce210SBarry Smith   /* Determine number, sending ranks, and length of incoming  */
55195fce210SBarry Smith   for (i=0; i<sf->nranks; i++) {
55295fce210SBarry Smith     rlengths[i] = sf->roffset[i+1] - sf->roffset[i]; /* Number of roots referenced by my leaves; for rank sf->ranks[i] */
55395fce210SBarry Smith   }
55495fce210SBarry Smith   ierr = PetscCommBuildTwoSided(comm,1,MPIU_INT,sf->nranks,sf->ranks,rlengths,&bas->niranks,&bas->iranks,(void**)&ilengths);CHKERRQ(ierr);
55595fce210SBarry Smith   ierr = PetscFree(rlengths);CHKERRQ(ierr);
55695fce210SBarry Smith 
55795fce210SBarry Smith   /* Send leaf identities to roots */
55895fce210SBarry Smith   for (i=0,bas->itotal=0; i<bas->niranks; i++) bas->itotal += ilengths[i];
559dcca6d9dSJed Brown   ierr = PetscMalloc2(bas->niranks+1,&bas->ioffset,bas->itotal,&bas->irootloc);CHKERRQ(ierr);
560dcca6d9dSJed Brown   ierr = PetscMalloc2(bas->niranks,&rootreqs,sf->nranks,&leafreqs);CHKERRQ(ierr);
56195fce210SBarry Smith   bas->ioffset[0] = 0;
56295fce210SBarry Smith   for (i=0; i<bas->niranks; i++) {
56395fce210SBarry Smith     bas->ioffset[i+1] = bas->ioffset[i] + ilengths[i];
56495fce210SBarry Smith     ierr = MPI_Irecv(bas->irootloc+bas->ioffset[i],ilengths[i],MPIU_INT,bas->iranks[i],bas->tag,comm,&rootreqs[i]);CHKERRQ(ierr);
56595fce210SBarry Smith   }
56695fce210SBarry Smith   for (i=0; i<sf->nranks; i++) {
56795fce210SBarry Smith     PetscMPIInt npoints;
56895fce210SBarry Smith     ierr = PetscMPIIntCast(sf->roffset[i+1] - sf->roffset[i],&npoints);CHKERRQ(ierr);
56995fce210SBarry Smith     ierr = MPI_Isend(sf->rremote+sf->roffset[i],npoints,MPIU_INT,sf->ranks[i],bas->tag,comm,&leafreqs[i]);CHKERRQ(ierr);
57095fce210SBarry Smith   }
571eb9baa12SBarry Smith   ierr = MPI_Waitall(bas->niranks,rootreqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
572eb9baa12SBarry Smith   ierr = MPI_Waitall(sf->nranks,leafreqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
57395fce210SBarry Smith   ierr = PetscFree(ilengths);CHKERRQ(ierr);
574eb9baa12SBarry Smith   ierr = PetscFree2(rootreqs,leafreqs);CHKERRQ(ierr);
57595fce210SBarry Smith   PetscFunctionReturn(0);
57695fce210SBarry Smith }
57795fce210SBarry Smith 
57895fce210SBarry Smith static PetscErrorCode PetscSFBasicPackTypeSetup(PetscSFBasicPack link,MPI_Datatype unit)
57995fce210SBarry Smith {
58095fce210SBarry Smith   PetscErrorCode ierr;
58195fce210SBarry Smith   PetscBool      isInt,isPetscInt,isPetscReal,is2Int,is2PetscInt;
582e37f7bbbSSatish Balay   PetscInt       nPetscIntContig,nPetscRealContig;
58395fce210SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
58495fce210SBarry Smith   PetscBool isPetscComplex;
585e37f7bbbSSatish Balay   PetscInt nPetscComplexContig;
58695fce210SBarry Smith #endif
58795fce210SBarry Smith 
58895fce210SBarry Smith   PetscFunctionBegin;
58995fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPI_INT,&isInt);CHKERRQ(ierr);
59095fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPIU_INT,&isPetscInt);CHKERRQ(ierr);
5912a74e5f4SJed Brown   ierr = MPIPetsc_Type_compare_contig(unit,MPIU_INT,&nPetscIntContig);CHKERRQ(ierr);
59295fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPIU_REAL,&isPetscReal);CHKERRQ(ierr);
5932a74e5f4SJed Brown   ierr = MPIPetsc_Type_compare_contig(unit,MPIU_REAL,&nPetscRealContig);CHKERRQ(ierr);
59495fce210SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
59595fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPIU_COMPLEX,&isPetscComplex);CHKERRQ(ierr);
5962a74e5f4SJed Brown   ierr = MPIPetsc_Type_compare_contig(unit,MPIU_COMPLEX,&nPetscComplexContig);CHKERRQ(ierr);
59795fce210SBarry Smith #endif
59895fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPI_2INT,&is2Int);CHKERRQ(ierr);
59995fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPIU_2INT,&is2PetscInt);CHKERRQ(ierr);
6002a74e5f4SJed Brown   link->bs = 1;
60153deab39SPeter Brune   if (isInt) {PackInit_int(link); PackInit_Logical_int(link); PackInit_Bitwise_int(link);}
60253deab39SPeter Brune   else if (isPetscInt) {PackInit_PetscInt(link); PackInit_Logical_PetscInt(link); PackInit_Bitwise_PetscInt(link);}
60353deab39SPeter Brune   else if (isPetscReal) {PackInit_PetscReal(link); PackInit_Logical_PetscReal(link);}
60495fce210SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
6052a74e5f4SJed Brown   else if (isPetscComplex) PackInit_PetscComplex_1(link);
60695fce210SBarry Smith #endif
60795fce210SBarry Smith   else if (is2Int) PackInit_int_int(link);
60895fce210SBarry Smith   else if (is2PetscInt) PackInit_PetscInt_PetscInt(link);
6092a74e5f4SJed Brown   else if (nPetscIntContig) {
6102a74e5f4SJed Brown     if (nPetscIntContig%7 == 0) PackInit_PetscInt_7(link);
6112a74e5f4SJed Brown     else if (nPetscIntContig%5 == 0) PackInit_PetscInt_5(link);
6122a74e5f4SJed Brown     else if (nPetscIntContig%4 == 0) PackInit_PetscInt_4(link);
6132a74e5f4SJed Brown     else if (nPetscIntContig%3 == 0) PackInit_PetscInt_3(link);
6142a74e5f4SJed Brown     else if (nPetscIntContig%2 == 0) PackInit_PetscInt_2(link);
6152a74e5f4SJed Brown     else PackInit_PetscInt(link);
6162a74e5f4SJed Brown     link->bs = nPetscIntContig;
6172a74e5f4SJed Brown     link->unitbytes *= nPetscIntContig;
6182a74e5f4SJed Brown   } else if (nPetscRealContig) {
6192a74e5f4SJed Brown     if (nPetscRealContig%7 == 0) PackInit_PetscReal_7(link);
6202a74e5f4SJed Brown     else if (nPetscRealContig%5 == 0) PackInit_PetscReal_5(link);
6212a74e5f4SJed Brown     else if (nPetscRealContig%4 == 0) PackInit_PetscReal_4(link);
6222a74e5f4SJed Brown     else if (nPetscRealContig%3 == 0) PackInit_PetscReal_3(link);
6232a74e5f4SJed Brown     else if (nPetscRealContig%2 == 0) PackInit_PetscReal_2(link);
6242a74e5f4SJed Brown     else PackInit_PetscReal(link);
6252a74e5f4SJed Brown     link->bs = nPetscRealContig;
6262a74e5f4SJed Brown     link->unitbytes *= nPetscRealContig;
6272a74e5f4SJed Brown #if defined(PETSC_HAVE_COMPLEX)
6282a74e5f4SJed Brown   } else if (nPetscComplexContig) {
6292a74e5f4SJed Brown     if (nPetscComplexContig%7 == 0) PackInit_PetscComplex_7(link);
6302a74e5f4SJed Brown     else if (nPetscComplexContig%5 == 0) PackInit_PetscComplex_5(link);
6312a74e5f4SJed Brown     else if (nPetscComplexContig%4 == 0) PackInit_PetscComplex_4(link);
6322a74e5f4SJed Brown     else if (nPetscComplexContig%3 == 0) PackInit_PetscComplex_3(link);
6332a74e5f4SJed Brown     else if (nPetscComplexContig%2 == 0) PackInit_PetscComplex_2(link);
6342a74e5f4SJed Brown     else PackInit_PetscComplex_1(link);
6352a74e5f4SJed Brown     link->bs = nPetscComplexContig;
6362a74e5f4SJed Brown     link->unitbytes *= nPetscComplexContig;
6372a74e5f4SJed Brown #endif
6382a74e5f4SJed Brown   } else {
639c3a59b84SJed Brown     MPI_Aint lb,bytes;
640c3a59b84SJed Brown     ierr = MPI_Type_get_extent(unit,&lb,&bytes);CHKERRQ(ierr);
641c3a59b84SJed Brown     if (lb != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Datatype with nonzero lower bound %ld\n",(long)lb);
64295fce210SBarry Smith     if (bytes % sizeof(int)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for type size not divisible by %D",sizeof(int));
64395fce210SBarry Smith     switch (bytes / sizeof(int)) {
64495fce210SBarry Smith     case 1: PackInit_block_int_1(link); break;
64595fce210SBarry Smith     case 2: PackInit_block_int_2(link); break;
64695fce210SBarry Smith     case 3: PackInit_block_int_3(link); break;
64795fce210SBarry Smith     case 4: PackInit_block_int_4(link); break;
64895fce210SBarry Smith     case 5: PackInit_block_int_5(link); break;
64995fce210SBarry Smith     case 6: PackInit_block_int_6(link); break;
65095fce210SBarry Smith     case 7: PackInit_block_int_7(link); break;
65195fce210SBarry Smith     case 8: PackInit_block_int_8(link); break;
65295fce210SBarry Smith     default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for arbitrary block sizes");
65395fce210SBarry Smith     }
65495fce210SBarry Smith   }
65595fce210SBarry Smith   ierr = MPI_Type_dup(unit,&link->unit);CHKERRQ(ierr);
65695fce210SBarry Smith   PetscFunctionReturn(0);
65795fce210SBarry Smith }
65895fce210SBarry Smith 
6592a74e5f4SJed Brown static PetscErrorCode PetscSFBasicPackGetUnpackOp(PetscSF sf,PetscSFBasicPack link,MPI_Op op,void (**UnpackOp)(PetscInt,PetscInt,const PetscInt*,void*,const void*))
66095fce210SBarry Smith {
66195fce210SBarry Smith   PetscFunctionBegin;
66295fce210SBarry Smith   *UnpackOp = NULL;
6638bfbc91cSJed Brown   if (op == MPIU_REPLACE) *UnpackOp = link->UnpackInsert;
66495fce210SBarry Smith   else if (op == MPI_SUM || op == MPIU_SUM) *UnpackOp = link->UnpackAdd;
66553deab39SPeter Brune   else if (op == MPI_PROD) *UnpackOp = link->UnpackMult;
66695fce210SBarry Smith   else if (op == MPI_MAX || op == MPIU_MAX) *UnpackOp = link->UnpackMax;
66795fce210SBarry Smith   else if (op == MPI_MIN || op == MPIU_MIN) *UnpackOp = link->UnpackMin;
66853deab39SPeter Brune   else if (op == MPI_LAND) *UnpackOp = link->UnpackLAND;
66953deab39SPeter Brune   else if (op == MPI_BAND) *UnpackOp = link->UnpackBAND;
67053deab39SPeter Brune   else if (op == MPI_LOR) *UnpackOp = link->UnpackLOR;
67153deab39SPeter Brune   else if (op == MPI_BOR) *UnpackOp = link->UnpackBOR;
67253deab39SPeter Brune   else if (op == MPI_LXOR) *UnpackOp = link->UnpackLXOR;
67353deab39SPeter Brune   else if (op == MPI_BXOR) *UnpackOp = link->UnpackBXOR;
67495fce210SBarry Smith   else if (op == MPI_MAXLOC) *UnpackOp = link->UnpackMaxloc;
67595fce210SBarry Smith   else if (op == MPI_MINLOC) *UnpackOp = link->UnpackMinloc;
6761620da3bSToby Isaac   else *UnpackOp = NULL;
67795fce210SBarry Smith   PetscFunctionReturn(0);
67895fce210SBarry Smith }
6792a74e5f4SJed Brown static PetscErrorCode PetscSFBasicPackGetFetchAndOp(PetscSF sf,PetscSFBasicPack link,MPI_Op op,void (**FetchAndOp)(PetscInt,PetscInt,const PetscInt*,void*,void*))
68095fce210SBarry Smith {
68195fce210SBarry Smith   PetscFunctionBegin;
68295fce210SBarry Smith   *FetchAndOp = NULL;
6838bfbc91cSJed Brown   if (op == MPIU_REPLACE) *FetchAndOp = link->FetchAndInsert;
68495fce210SBarry Smith   else if (op == MPI_SUM || op == MPIU_SUM) *FetchAndOp = link->FetchAndAdd;
68595fce210SBarry Smith   else if (op == MPI_MAX || op == MPIU_MAX) *FetchAndOp = link->FetchAndMax;
68695fce210SBarry Smith   else if (op == MPI_MIN || op == MPIU_MIN) *FetchAndOp = link->FetchAndMin;
68795fce210SBarry Smith   else if (op == MPI_MAXLOC) *FetchAndOp = link->FetchAndMaxloc;
68895fce210SBarry Smith   else if (op == MPI_MINLOC) *FetchAndOp = link->FetchAndMinloc;
68953deab39SPeter Brune   else if (op == MPI_PROD)   *FetchAndOp = link->FetchAndMult;
69053deab39SPeter Brune   else if (op == MPI_LAND)   *FetchAndOp = link->FetchAndLAND;
69153deab39SPeter Brune   else if (op == MPI_BAND)   *FetchAndOp = link->FetchAndBAND;
69253deab39SPeter Brune   else if (op == MPI_LOR)    *FetchAndOp = link->FetchAndLOR;
69353deab39SPeter Brune   else if (op == MPI_BOR)    *FetchAndOp = link->FetchAndBOR;
69453deab39SPeter Brune   else if (op == MPI_LXOR)   *FetchAndOp = link->FetchAndLXOR;
69553deab39SPeter Brune   else if (op == MPI_BXOR)   *FetchAndOp = link->FetchAndBXOR;
69695fce210SBarry Smith   else SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_SUP,"No support for MPI_Op");
69795fce210SBarry Smith   PetscFunctionReturn(0);
69895fce210SBarry Smith }
69995fce210SBarry Smith 
70095fce210SBarry Smith static PetscErrorCode PetscSFBasicPackGetReqs(PetscSF sf,PetscSFBasicPack link,MPI_Request **rootreqs,MPI_Request **leafreqs)
70195fce210SBarry Smith {
70295fce210SBarry Smith   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
70395fce210SBarry Smith 
70495fce210SBarry Smith   PetscFunctionBegin;
70595fce210SBarry Smith   if (rootreqs) *rootreqs = link->requests;
70695fce210SBarry Smith   if (leafreqs) *leafreqs = link->requests + bas->niranks;
70795fce210SBarry Smith   PetscFunctionReturn(0);
70895fce210SBarry Smith }
70995fce210SBarry Smith 
71095fce210SBarry Smith static PetscErrorCode PetscSFBasicPackWaitall(PetscSF sf,PetscSFBasicPack link)
71195fce210SBarry Smith {
71295fce210SBarry Smith   PetscSF_Basic  *bas = (PetscSF_Basic*)sf->data;
71395fce210SBarry Smith   PetscErrorCode ierr;
71495fce210SBarry Smith 
71595fce210SBarry Smith   PetscFunctionBegin;
71695fce210SBarry Smith   ierr = MPI_Waitall(bas->niranks+sf->nranks,link->requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
71795fce210SBarry Smith   PetscFunctionReturn(0);
71895fce210SBarry Smith }
71995fce210SBarry Smith 
72095fce210SBarry Smith static PetscErrorCode PetscSFBasicGetRootInfo(PetscSF sf,PetscInt *nrootranks,const PetscMPIInt **rootranks,const PetscInt **rootoffset,const PetscInt **rootloc)
72195fce210SBarry Smith {
72295fce210SBarry Smith   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
72395fce210SBarry Smith 
72495fce210SBarry Smith   PetscFunctionBegin;
72595fce210SBarry Smith   if (nrootranks) *nrootranks = bas->niranks;
72695fce210SBarry Smith   if (rootranks)  *rootranks  = bas->iranks;
72795fce210SBarry Smith   if (rootoffset) *rootoffset = bas->ioffset;
72895fce210SBarry Smith   if (rootloc)    *rootloc    = bas->irootloc;
72995fce210SBarry Smith   PetscFunctionReturn(0);
73095fce210SBarry Smith }
73195fce210SBarry Smith 
73295fce210SBarry Smith static PetscErrorCode PetscSFBasicGetLeafInfo(PetscSF sf,PetscInt *nleafranks,const PetscMPIInt **leafranks,const PetscInt **leafoffset,const PetscInt **leafloc)
73395fce210SBarry Smith {
73495fce210SBarry Smith   PetscFunctionBegin;
73595fce210SBarry Smith   if (nleafranks) *nleafranks = sf->nranks;
73695fce210SBarry Smith   if (leafranks)  *leafranks  = sf->ranks;
73795fce210SBarry Smith   if (leafoffset) *leafoffset = sf->roffset;
73895fce210SBarry Smith   if (leafloc)    *leafloc    = sf->rmine;
73995fce210SBarry Smith   PetscFunctionReturn(0);
74095fce210SBarry Smith }
74195fce210SBarry Smith 
74295fce210SBarry Smith static PetscErrorCode PetscSFBasicGetPack(PetscSF sf,MPI_Datatype unit,const void *key,PetscSFBasicPack *mylink)
74395fce210SBarry Smith {
74495fce210SBarry Smith   PetscSF_Basic    *bas = (PetscSF_Basic*)sf->data;
74595fce210SBarry Smith   PetscErrorCode   ierr;
74695fce210SBarry Smith   PetscSFBasicPack link,*p;
747bf39f1bfSJed Brown   PetscInt         nrootranks,nleafranks,i;
74895fce210SBarry Smith   const PetscInt   *rootoffset,*leafoffset;
74995fce210SBarry Smith 
75095fce210SBarry Smith   PetscFunctionBegin;
75195fce210SBarry Smith   /* Look for types in cache */
75295fce210SBarry Smith   for (p=&bas->avail; (link=*p); p=&link->next) {
75395fce210SBarry Smith     PetscBool match;
75495fce210SBarry Smith     ierr = MPIPetsc_Type_compare(unit,link->unit,&match);CHKERRQ(ierr);
75595fce210SBarry Smith     if (match) {
75695fce210SBarry Smith       *p = link->next;          /* Remove from available list */
75795fce210SBarry Smith       goto found;
75895fce210SBarry Smith     }
75995fce210SBarry Smith   }
76095fce210SBarry Smith 
76195fce210SBarry Smith   /* Create new composite types for each send rank */
76295fce210SBarry Smith   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,NULL,&rootoffset,NULL);CHKERRQ(ierr);
76395fce210SBarry Smith   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,NULL,&leafoffset,NULL);CHKERRQ(ierr);
764b00a9115SJed Brown   ierr = PetscNew(&link);CHKERRQ(ierr);
76595fce210SBarry Smith   ierr = PetscSFBasicPackTypeSetup(link,unit);CHKERRQ(ierr);
766bf39f1bfSJed Brown   ierr = PetscMalloc2(nrootranks,&link->root,nleafranks,&link->leaf);CHKERRQ(ierr);
767bf39f1bfSJed Brown   for (i=0; i<nrootranks; i++) {
768bf39f1bfSJed Brown     ierr = PetscMalloc((rootoffset[i+1]-rootoffset[i])*link->unitbytes,&link->root[i]);CHKERRQ(ierr);
769bf39f1bfSJed Brown   }
770bf39f1bfSJed Brown   for (i=0; i<nleafranks; i++) {
771bf39f1bfSJed Brown     ierr = PetscMalloc((leafoffset[i+1]-leafoffset[i])*link->unitbytes,&link->leaf[i]);CHKERRQ(ierr);
772bf39f1bfSJed Brown   }
7730f880758SHong Zhang   ierr = PetscCalloc1(nrootranks+nleafranks,&link->requests);CHKERRQ(ierr);
77495fce210SBarry Smith 
77595fce210SBarry Smith found:
77695fce210SBarry Smith   link->key  = key;
77795fce210SBarry Smith   link->next = bas->inuse;
77895fce210SBarry Smith   bas->inuse = link;
77995fce210SBarry Smith 
78095fce210SBarry Smith   *mylink = link;
78195fce210SBarry Smith   PetscFunctionReturn(0);
78295fce210SBarry Smith }
78395fce210SBarry Smith 
78495fce210SBarry Smith static PetscErrorCode PetscSFBasicGetPackInUse(PetscSF sf,MPI_Datatype unit,const void *key,PetscCopyMode cmode,PetscSFBasicPack *mylink)
78595fce210SBarry Smith {
78695fce210SBarry Smith   PetscSF_Basic    *bas = (PetscSF_Basic*)sf->data;
78795fce210SBarry Smith   PetscErrorCode   ierr;
78895fce210SBarry Smith   PetscSFBasicPack link,*p;
78995fce210SBarry Smith 
79095fce210SBarry Smith   PetscFunctionBegin;
79195fce210SBarry Smith   /* Look for types in cache */
79295fce210SBarry Smith   for (p=&bas->inuse; (link=*p); p=&link->next) {
79395fce210SBarry Smith     PetscBool match;
79495fce210SBarry Smith     ierr = MPIPetsc_Type_compare(unit,link->unit,&match);CHKERRQ(ierr);
79545038af3SLawrence Mitchell     if (match && (key == link->key)) {
79695fce210SBarry Smith       switch (cmode) {
79795fce210SBarry Smith       case PETSC_OWN_POINTER: *p = link->next; break; /* Remove from inuse list */
79895fce210SBarry Smith       case PETSC_USE_POINTER: break;
79995fce210SBarry Smith       default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"invalid cmode");
80095fce210SBarry Smith       }
80195fce210SBarry Smith       *mylink = link;
80295fce210SBarry Smith       PetscFunctionReturn(0);
80395fce210SBarry Smith     }
80495fce210SBarry Smith   }
80595fce210SBarry Smith   SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Could not find pack");
80695fce210SBarry Smith   PetscFunctionReturn(0);
80795fce210SBarry Smith }
80895fce210SBarry Smith 
80995fce210SBarry Smith static PetscErrorCode PetscSFBasicReclaimPack(PetscSF sf,PetscSFBasicPack *link)
81095fce210SBarry Smith {
81195fce210SBarry Smith   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
81295fce210SBarry Smith 
81395fce210SBarry Smith   PetscFunctionBegin;
81495fce210SBarry Smith   (*link)->key  = NULL;
81595fce210SBarry Smith   (*link)->next = bas->avail;
81695fce210SBarry Smith   bas->avail    = *link;
81795fce210SBarry Smith   *link         = NULL;
81895fce210SBarry Smith   PetscFunctionReturn(0);
81995fce210SBarry Smith }
82095fce210SBarry Smith 
8214416b707SBarry Smith static PetscErrorCode PetscSFSetFromOptions_Basic(PetscOptionItems *PetscOptionsObject,PetscSF sf)
82295fce210SBarry Smith {
82395fce210SBarry Smith   PetscErrorCode ierr;
82495fce210SBarry Smith 
82595fce210SBarry Smith   PetscFunctionBegin;
826e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"PetscSF Basic options");CHKERRQ(ierr);
82795fce210SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
82895fce210SBarry Smith   PetscFunctionReturn(0);
82995fce210SBarry Smith }
83095fce210SBarry Smith 
83195fce210SBarry Smith static PetscErrorCode PetscSFReset_Basic(PetscSF sf)
83295fce210SBarry Smith {
83395fce210SBarry Smith   PetscSF_Basic    *bas = (PetscSF_Basic*)sf->data;
83495fce210SBarry Smith   PetscErrorCode   ierr;
83595fce210SBarry Smith   PetscSFBasicPack link,next;
83695fce210SBarry Smith 
83795fce210SBarry Smith   PetscFunctionBegin;
83895fce210SBarry Smith   ierr = PetscFree(bas->iranks);CHKERRQ(ierr);
83995fce210SBarry Smith   ierr = PetscFree2(bas->ioffset,bas->irootloc);CHKERRQ(ierr);
84095fce210SBarry Smith   if (bas->inuse) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Outstanding operation has not been completed");
84195fce210SBarry Smith   for (link=bas->avail; link; link=next) {
842bf39f1bfSJed Brown     PetscInt i;
84395fce210SBarry Smith     next = link->next;
844bb86044dSJed Brown #if defined(PETSC_HAVE_MPI_TYPE_DUP)
84595fce210SBarry Smith     ierr = MPI_Type_free(&link->unit);CHKERRQ(ierr);
846bb86044dSJed Brown #endif
847bf39f1bfSJed Brown     for (i=0; i<bas->niranks; i++) {ierr = PetscFree(link->root[i]);CHKERRQ(ierr);}
848bf39f1bfSJed Brown     for (i=0; i<sf->nranks; i++) {ierr = PetscFree(link->leaf[i]);CHKERRQ(ierr);}
84995fce210SBarry Smith     ierr = PetscFree2(link->root,link->leaf);CHKERRQ(ierr);
85095fce210SBarry Smith     ierr = PetscFree(link->requests);CHKERRQ(ierr);
85195fce210SBarry Smith     ierr = PetscFree(link);CHKERRQ(ierr);
85295fce210SBarry Smith   }
85395fce210SBarry Smith   bas->avail = NULL;
85495fce210SBarry Smith   PetscFunctionReturn(0);
85595fce210SBarry Smith }
85695fce210SBarry Smith 
85795fce210SBarry Smith static PetscErrorCode PetscSFDestroy_Basic(PetscSF sf)
85895fce210SBarry Smith {
85995fce210SBarry Smith   PetscErrorCode ierr;
86095fce210SBarry Smith 
86195fce210SBarry Smith   PetscFunctionBegin;
86295fce210SBarry Smith   ierr = PetscFree(sf->data);CHKERRQ(ierr);
86395fce210SBarry Smith   PetscFunctionReturn(0);
86495fce210SBarry Smith }
86595fce210SBarry Smith 
86695fce210SBarry Smith static PetscErrorCode PetscSFView_Basic(PetscSF sf,PetscViewer viewer)
86795fce210SBarry Smith {
86895fce210SBarry Smith   /* PetscSF_Basic *bas = (PetscSF_Basic*)sf->data; */
86995fce210SBarry Smith   PetscErrorCode ierr;
87095fce210SBarry Smith   PetscBool      iascii;
87195fce210SBarry Smith 
87295fce210SBarry Smith   PetscFunctionBegin;
87395fce210SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
87495fce210SBarry Smith   if (iascii) {
87595fce210SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  sort=%s\n",sf->rankorder ? "rank-order" : "unordered");CHKERRQ(ierr);
87695fce210SBarry Smith   }
87795fce210SBarry Smith   PetscFunctionReturn(0);
87895fce210SBarry Smith }
87995fce210SBarry Smith 
88095fce210SBarry Smith /* Send from roots to leaves */
88195fce210SBarry Smith static PetscErrorCode PetscSFBcastBegin_Basic(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
88295fce210SBarry Smith {
88395fce210SBarry Smith   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;
88495fce210SBarry Smith   PetscErrorCode    ierr;
88595fce210SBarry Smith   PetscSFBasicPack  link;
88695fce210SBarry Smith   PetscInt          i,nrootranks,nleafranks;
88795fce210SBarry Smith   const PetscInt    *rootoffset,*leafoffset,*rootloc,*leafloc;
88895fce210SBarry Smith   const PetscMPIInt *rootranks,*leafranks;
88995fce210SBarry Smith   MPI_Request       *rootreqs,*leafreqs;
89095fce210SBarry Smith   size_t            unitbytes;
89195fce210SBarry Smith 
89295fce210SBarry Smith   PetscFunctionBegin;
89395fce210SBarry Smith   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,&rootranks,&rootoffset,&rootloc);CHKERRQ(ierr);
89495fce210SBarry Smith   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,&leafranks,&leafoffset,&leafloc);CHKERRQ(ierr);
89595fce210SBarry Smith   ierr = PetscSFBasicGetPack(sf,unit,rootdata,&link);CHKERRQ(ierr);
89695fce210SBarry Smith 
89795fce210SBarry Smith   unitbytes = link->unitbytes;
89895fce210SBarry Smith 
89995fce210SBarry Smith   ierr = PetscSFBasicPackGetReqs(sf,link,&rootreqs,&leafreqs);CHKERRQ(ierr);
90095fce210SBarry Smith   /* Eagerly post leaf receives */
90195fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
90295fce210SBarry Smith     PetscMPIInt n = leafoffset[i+1] - leafoffset[i];
903bf39f1bfSJed Brown     ierr = MPI_Irecv(link->leaf[i],n,unit,leafranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&leafreqs[i]);CHKERRQ(ierr);
90495fce210SBarry Smith   }
90595fce210SBarry Smith   /* Pack and send root data */
90695fce210SBarry Smith   for (i=0; i<nrootranks; i++) {
90795fce210SBarry Smith     PetscMPIInt n          = rootoffset[i+1] - rootoffset[i];
908bf39f1bfSJed Brown     void        *packstart = link->root[i];
9092a74e5f4SJed Brown     (*link->Pack)(n,link->bs,rootloc+rootoffset[i],rootdata,packstart);
91095fce210SBarry Smith     ierr = MPI_Isend(packstart,n,unit,rootranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&rootreqs[i]);CHKERRQ(ierr);
91195fce210SBarry Smith   }
91295fce210SBarry Smith   PetscFunctionReturn(0);
91395fce210SBarry Smith }
91495fce210SBarry Smith 
91595fce210SBarry Smith PetscErrorCode PetscSFBcastEnd_Basic(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
91695fce210SBarry Smith {
91795fce210SBarry Smith   PetscErrorCode   ierr;
91895fce210SBarry Smith   PetscSFBasicPack link;
91995fce210SBarry Smith   PetscInt         i,nleafranks;
92095fce210SBarry Smith   const PetscInt   *leafoffset,*leafloc;
92195fce210SBarry Smith 
92295fce210SBarry Smith   PetscFunctionBegin;
92395fce210SBarry Smith   ierr = PetscSFBasicGetPackInUse(sf,unit,rootdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr);
92495fce210SBarry Smith   ierr = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
92595fce210SBarry Smith   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,NULL,&leafoffset,&leafloc);CHKERRQ(ierr);
92695fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
92795fce210SBarry Smith     PetscMPIInt n          = leafoffset[i+1] - leafoffset[i];
928bf39f1bfSJed Brown     const void  *packstart = link->leaf[i];
9292a74e5f4SJed Brown     (*link->UnpackInsert)(n,link->bs,leafloc+leafoffset[i],leafdata,packstart);
93095fce210SBarry Smith   }
93195fce210SBarry Smith   ierr = PetscSFBasicReclaimPack(sf,&link);CHKERRQ(ierr);
93295fce210SBarry Smith   PetscFunctionReturn(0);
93395fce210SBarry Smith }
93495fce210SBarry Smith 
93595fce210SBarry Smith /* leaf -> root with reduction */
93695fce210SBarry Smith PetscErrorCode PetscSFReduceBegin_Basic(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
93795fce210SBarry Smith {
93895fce210SBarry Smith   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;
93995fce210SBarry Smith   PetscSFBasicPack  link;
94095fce210SBarry Smith   PetscErrorCode    ierr;
94195fce210SBarry Smith   PetscInt          i,nrootranks,nleafranks;
94295fce210SBarry Smith   const PetscInt    *rootoffset,*leafoffset,*rootloc,*leafloc;
94395fce210SBarry Smith   const PetscMPIInt *rootranks,*leafranks;
94495fce210SBarry Smith   MPI_Request       *rootreqs,*leafreqs;
94595fce210SBarry Smith   size_t            unitbytes;
94695fce210SBarry Smith 
94795fce210SBarry Smith   PetscFunctionBegin;
94895fce210SBarry Smith   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,&rootranks,&rootoffset,&rootloc);CHKERRQ(ierr);
94995fce210SBarry Smith   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,&leafranks,&leafoffset,&leafloc);CHKERRQ(ierr);
95095fce210SBarry Smith   ierr = PetscSFBasicGetPack(sf,unit,rootdata,&link);CHKERRQ(ierr);
95195fce210SBarry Smith 
95295fce210SBarry Smith   unitbytes = link->unitbytes;
95395fce210SBarry Smith 
95495fce210SBarry Smith   ierr = PetscSFBasicPackGetReqs(sf,link,&rootreqs,&leafreqs);CHKERRQ(ierr);
95595fce210SBarry Smith   /* Eagerly post root receives */
95695fce210SBarry Smith   for (i=0; i<nrootranks; i++) {
95795fce210SBarry Smith     PetscMPIInt n = rootoffset[i+1] - rootoffset[i];
958bf39f1bfSJed Brown     ierr = MPI_Irecv(link->root[i],n,unit,rootranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&rootreqs[i]);CHKERRQ(ierr);
95995fce210SBarry Smith   }
96095fce210SBarry Smith   /* Pack and send leaf data */
96195fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
96295fce210SBarry Smith     PetscMPIInt n          = leafoffset[i+1] - leafoffset[i];
963bf39f1bfSJed Brown     void        *packstart = link->leaf[i];
9642a74e5f4SJed Brown     (*link->Pack)(n,link->bs,leafloc+leafoffset[i],leafdata,packstart);
96595fce210SBarry Smith     ierr = MPI_Isend(packstart,n,unit,leafranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&leafreqs[i]);CHKERRQ(ierr);
96695fce210SBarry Smith   }
96795fce210SBarry Smith   PetscFunctionReturn(0);
96895fce210SBarry Smith }
96995fce210SBarry Smith 
97095fce210SBarry Smith static PetscErrorCode PetscSFReduceEnd_Basic(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
97195fce210SBarry Smith {
9722a74e5f4SJed Brown   void             (*UnpackOp)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
97395fce210SBarry Smith   PetscErrorCode   ierr;
97495fce210SBarry Smith   PetscSFBasicPack link;
97595fce210SBarry Smith   PetscInt         i,nrootranks;
9761620da3bSToby Isaac   PetscMPIInt      typesize = -1;
97795fce210SBarry Smith   const PetscInt   *rootoffset,*rootloc;
97895fce210SBarry Smith 
97995fce210SBarry Smith   PetscFunctionBegin;
98095fce210SBarry Smith   ierr = PetscSFBasicGetPackInUse(sf,unit,rootdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr);
98195fce210SBarry Smith   /* This implementation could be changed to unpack as receives arrive, at the cost of non-determinism */
98295fce210SBarry Smith   ierr = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
98395fce210SBarry Smith   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,NULL,&rootoffset,&rootloc);CHKERRQ(ierr);
98495fce210SBarry Smith   ierr = PetscSFBasicPackGetUnpackOp(sf,link,op,&UnpackOp);CHKERRQ(ierr);
9851620da3bSToby Isaac   if (UnpackOp) {
9861620da3bSToby Isaac     typesize = link->unitbytes;
9871620da3bSToby Isaac   }
9881620da3bSToby Isaac   else {
9891620da3bSToby Isaac     ierr = MPI_Type_size(unit,&typesize);CHKERRQ(ierr);
9901620da3bSToby Isaac   }
99195fce210SBarry Smith   for (i=0; i<nrootranks; i++) {
99295fce210SBarry Smith     PetscMPIInt n   = rootoffset[i+1] - rootoffset[i];
993bf39f1bfSJed Brown     char *packstart = (char *) link->root[i];
99495fce210SBarry Smith 
9951620da3bSToby Isaac     if (UnpackOp) {
99688ac5ce8SToby Isaac       (*UnpackOp)(n,link->bs,rootloc+rootoffset[i],rootdata,(const void *)packstart);
99795fce210SBarry Smith     }
998a9ec7e9eSToby Isaac #if PETSC_HAVE_MPI_REDUCE_LOCAL
99978ec4774SToby Isaac     else if (n) { /* the op should be defined to operate on the whole datatype, so we ignore link->bs */
100078ec4774SToby Isaac       PetscInt j;
10011620da3bSToby Isaac 
10021620da3bSToby Isaac       for (j = 0; j < n; j++) {
100378ec4774SToby Isaac         ierr = MPI_Reduce_local(packstart+j*typesize,((char *) rootdata)+(rootloc[rootoffset[i]+j])*typesize,1,unit,op);CHKERRQ(ierr);
10041620da3bSToby Isaac       }
10051620da3bSToby Isaac     }
1006a9ec7e9eSToby Isaac #else
1007a9ec7e9eSToby Isaac     else {
1008a9ec7e9eSToby Isaac       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No unpacking reduction operation for this MPI_Op");
1009a9ec7e9eSToby Isaac     }
1010a9ec7e9eSToby Isaac #endif
10111620da3bSToby Isaac   }
101295fce210SBarry Smith   ierr = PetscSFBasicReclaimPack(sf,&link);CHKERRQ(ierr);
101395fce210SBarry Smith   PetscFunctionReturn(0);
101495fce210SBarry Smith }
101595fce210SBarry Smith 
101695fce210SBarry Smith static PetscErrorCode PetscSFFetchAndOpBegin_Basic(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
101795fce210SBarry Smith {
101895fce210SBarry Smith   PetscErrorCode ierr;
101995fce210SBarry Smith 
102095fce210SBarry Smith   PetscFunctionBegin;
102195fce210SBarry Smith   ierr = PetscSFReduceBegin_Basic(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
102295fce210SBarry Smith   PetscFunctionReturn(0);
102395fce210SBarry Smith }
102495fce210SBarry Smith 
102595fce210SBarry Smith static PetscErrorCode PetscSFFetchAndOpEnd_Basic(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
102695fce210SBarry Smith {
102795fce210SBarry Smith   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;
10282a74e5f4SJed Brown   void              (*FetchAndOp)(PetscInt,PetscInt,const PetscInt*,void*,void*);
102995fce210SBarry Smith   PetscErrorCode    ierr;
103095fce210SBarry Smith   PetscSFBasicPack  link;
103195fce210SBarry Smith   PetscInt          i,nrootranks,nleafranks;
103295fce210SBarry Smith   const PetscInt    *rootoffset,*leafoffset,*rootloc,*leafloc;
103395fce210SBarry Smith   const PetscMPIInt *rootranks,*leafranks;
103495fce210SBarry Smith   MPI_Request       *rootreqs,*leafreqs;
103595fce210SBarry Smith   size_t            unitbytes;
103695fce210SBarry Smith 
103795fce210SBarry Smith   PetscFunctionBegin;
103895fce210SBarry Smith   ierr = PetscSFBasicGetPackInUse(sf,unit,rootdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr);
103995fce210SBarry Smith   /* This implementation could be changed to unpack as receives arrive, at the cost of non-determinism */
104095fce210SBarry Smith   ierr      = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
104195fce210SBarry Smith   unitbytes = link->unitbytes;
104295fce210SBarry Smith   ierr      = PetscSFBasicGetRootInfo(sf,&nrootranks,&rootranks,&rootoffset,&rootloc);CHKERRQ(ierr);
104395fce210SBarry Smith   ierr      = PetscSFBasicGetLeafInfo(sf,&nleafranks,&leafranks,&leafoffset,&leafloc);CHKERRQ(ierr);
104495fce210SBarry Smith   ierr      = PetscSFBasicPackGetReqs(sf,link,&rootreqs,&leafreqs);CHKERRQ(ierr);
104595fce210SBarry Smith   /* Post leaf receives */
104695fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
104795fce210SBarry Smith     PetscMPIInt n = leafoffset[i+1] - leafoffset[i];
1048bf39f1bfSJed Brown     ierr = MPI_Irecv(link->leaf[i],n,unit,leafranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&leafreqs[i]);CHKERRQ(ierr);
104995fce210SBarry Smith   }
105095fce210SBarry Smith   /* Process local fetch-and-op, post root sends */
105195fce210SBarry Smith   ierr = PetscSFBasicPackGetFetchAndOp(sf,link,op,&FetchAndOp);CHKERRQ(ierr);
105295fce210SBarry Smith   for (i=0; i<nrootranks; i++) {
105395fce210SBarry Smith     PetscMPIInt n          = rootoffset[i+1] - rootoffset[i];
1054bf39f1bfSJed Brown     void        *packstart = link->root[i];
105595fce210SBarry Smith 
10562a74e5f4SJed Brown     (*FetchAndOp)(n,link->bs,rootloc+rootoffset[i],rootdata,packstart);
105795fce210SBarry Smith     ierr = MPI_Isend(packstart,n,unit,rootranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&rootreqs[i]);CHKERRQ(ierr);
105895fce210SBarry Smith   }
105995fce210SBarry Smith   ierr = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
106095fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
106195fce210SBarry Smith     PetscMPIInt n          = leafoffset[i+1] - leafoffset[i];
1062bf39f1bfSJed Brown     const void  *packstart = link->leaf[i];
10632a74e5f4SJed Brown     (*link->UnpackInsert)(n,link->bs,leafloc+leafoffset[i],leafupdate,packstart);
106495fce210SBarry Smith   }
106595fce210SBarry Smith   ierr = PetscSFBasicReclaimPack(sf,&link);CHKERRQ(ierr);
106695fce210SBarry Smith   PetscFunctionReturn(0);
106795fce210SBarry Smith }
106895fce210SBarry Smith 
10698cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscSFCreate_Basic(PetscSF sf)
107095fce210SBarry Smith {
107195fce210SBarry Smith   PetscSF_Basic  *bas = (PetscSF_Basic*)sf->data;
107295fce210SBarry Smith   PetscErrorCode ierr;
107395fce210SBarry Smith 
107495fce210SBarry Smith   PetscFunctionBegin;
107595fce210SBarry Smith   sf->ops->SetUp           = PetscSFSetUp_Basic;
107695fce210SBarry Smith   sf->ops->SetFromOptions  = PetscSFSetFromOptions_Basic;
107795fce210SBarry Smith   sf->ops->Reset           = PetscSFReset_Basic;
107895fce210SBarry Smith   sf->ops->Destroy         = PetscSFDestroy_Basic;
107995fce210SBarry Smith   sf->ops->View            = PetscSFView_Basic;
108095fce210SBarry Smith   sf->ops->BcastBegin      = PetscSFBcastBegin_Basic;
108195fce210SBarry Smith   sf->ops->BcastEnd        = PetscSFBcastEnd_Basic;
108295fce210SBarry Smith   sf->ops->ReduceBegin     = PetscSFReduceBegin_Basic;
108395fce210SBarry Smith   sf->ops->ReduceEnd       = PetscSFReduceEnd_Basic;
108495fce210SBarry Smith   sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Basic;
108595fce210SBarry Smith   sf->ops->FetchAndOpEnd   = PetscSFFetchAndOpEnd_Basic;
108695fce210SBarry Smith 
1087b00a9115SJed Brown   ierr     = PetscNewLog(sf,&bas);CHKERRQ(ierr);
108895fce210SBarry Smith   sf->data = (void*)bas;
108995fce210SBarry Smith   PetscFunctionReturn(0);
109095fce210SBarry Smith }
1091