xref: /petsc/src/vec/is/sf/impls/basic/sfbasic.c (revision dcca6d9d80ebd869fe6029bd05a3aa9faafef49e)
195fce210SBarry Smith #define PETSC_DESIRE_COMPLEX
295fce210SBarry Smith #include <petsc-private/sfimpl.h> /*I "petscsf.h" I*/
395fce210SBarry Smith 
495fce210SBarry Smith typedef struct _n_PetscSFBasicPack *PetscSFBasicPack;
595fce210SBarry Smith struct _n_PetscSFBasicPack {
695fce210SBarry Smith   void (*Pack)(PetscInt,const PetscInt*,const void*,void*);
795fce210SBarry Smith   void (*UnpackInsert)(PetscInt,const PetscInt*,void*,const void*);
895fce210SBarry Smith   void (*UnpackAdd)(PetscInt,const PetscInt*,void*,const void*);
995fce210SBarry Smith   void (*UnpackMin)(PetscInt,const PetscInt*,void*,const void*);
1095fce210SBarry Smith   void (*UnpackMax)(PetscInt,const PetscInt*,void*,const void*);
1195fce210SBarry Smith   void (*UnpackMinloc)(PetscInt,const PetscInt*,void*,const void*);
1295fce210SBarry Smith   void (*UnpackMaxloc)(PetscInt,const PetscInt*,void*,const void*);
1353deab39SPeter Brune   void (*UnpackMult)(PetscInt,const PetscInt*,void*,const void *);
1453deab39SPeter Brune   void (*UnpackLAND)(PetscInt,const PetscInt*,void*,const void *);
1553deab39SPeter Brune   void (*UnpackBAND)(PetscInt,const PetscInt*,void*,const void *);
1653deab39SPeter Brune   void (*UnpackLOR)(PetscInt,const PetscInt*,void*,const void *);
1753deab39SPeter Brune   void (*UnpackBOR)(PetscInt,const PetscInt*,void*,const void *);
1853deab39SPeter Brune   void (*UnpackLXOR)(PetscInt,const PetscInt*,void*,const void *);
1953deab39SPeter Brune   void (*UnpackBXOR)(PetscInt,const PetscInt*,void*,const void *);
2095fce210SBarry Smith   void (*FetchAndInsert)(PetscInt,const PetscInt*,void*,void*);
2195fce210SBarry Smith   void (*FetchAndAdd)(PetscInt,const PetscInt*,void*,void*);
2295fce210SBarry Smith   void (*FetchAndMin)(PetscInt,const PetscInt*,void*,void*);
2395fce210SBarry Smith   void (*FetchAndMax)(PetscInt,const PetscInt*,void*,void*);
2495fce210SBarry Smith   void (*FetchAndMinloc)(PetscInt,const PetscInt*,void*,void*);
2595fce210SBarry Smith   void (*FetchAndMaxloc)(PetscInt,const PetscInt*,void*,void*);
2653deab39SPeter Brune   void (*FetchAndMult)(PetscInt,const PetscInt*,void*,void*);
2753deab39SPeter Brune   void (*FetchAndLAND)(PetscInt,const PetscInt*,void*,void*);
2853deab39SPeter Brune   void (*FetchAndBAND)(PetscInt,const PetscInt*,void*,void*);
2953deab39SPeter Brune   void (*FetchAndLOR)(PetscInt,const PetscInt*,void*,void*);
3053deab39SPeter Brune   void (*FetchAndBOR)(PetscInt,const PetscInt*,void*,void*);
3153deab39SPeter Brune   void (*FetchAndLXOR)(PetscInt,const PetscInt*,void*,void*);
3253deab39SPeter Brune   void (*FetchAndBXOR)(PetscInt,const PetscInt*,void*,void*);
3395fce210SBarry Smith 
3495fce210SBarry Smith   MPI_Datatype     unit;
3595fce210SBarry Smith   size_t           unitbytes;   /* Number of bytes in a unit */
3695fce210SBarry Smith   const void       *key;        /* Array used as key for operation */
3795fce210SBarry Smith   char             *root;       /* Packed root data, contiguous by leaf rank */
3895fce210SBarry Smith   char             *leaf;       /* Packed leaf data, contiguous by root rank */
3995fce210SBarry Smith   MPI_Request      *requests;   /* Array of root requests followed by leaf requests */
4095fce210SBarry Smith   PetscSFBasicPack next;
4195fce210SBarry Smith };
4295fce210SBarry Smith 
4395fce210SBarry Smith typedef struct {
4495fce210SBarry Smith   PetscMPIInt      tag;
4595fce210SBarry Smith   PetscInt         niranks;     /* Number of incoming ranks (ranks accessing my roots) */
4695fce210SBarry Smith   PetscMPIInt      *iranks;     /* Array of ranks that reference my roots */
4795fce210SBarry Smith   PetscInt         itotal;      /* Total number of graph edges referencing my roots */
4895fce210SBarry Smith   PetscInt         *ioffset;    /* Array of length niranks+1 holding offset in irootloc[] for each rank */
4995fce210SBarry Smith   PetscInt         *irootloc;   /* Incoming roots referenced by ranks starting at ioffset[rank] */
5095fce210SBarry Smith   PetscSFBasicPack avail;       /* One or more entries per MPI Datatype, lazily constructed */
5195fce210SBarry Smith   PetscSFBasicPack inuse;       /* Buffers being used for transactions that have not yet completed */
5295fce210SBarry Smith } PetscSF_Basic;
5395fce210SBarry Smith 
5495fce210SBarry Smith #if !defined(PETSC_HAVE_MPI_TYPE_DUP) /* Danger: type is not reference counted; subject to ABA problem */
5595fce210SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MPI_Type_dup(MPI_Datatype datatype,MPI_Datatype *newtype)
5695fce210SBarry Smith {
5795fce210SBarry Smith   *newtype = datatype;
5895fce210SBarry Smith   return 0;
5995fce210SBarry Smith }
6095fce210SBarry Smith #endif
6195fce210SBarry Smith 
6295fce210SBarry Smith /*
6395fce210SBarry Smith  * MPI_Reduce_local is not really useful because it can't handle sparse data and it vectorizes "in the wrong direction",
6495fce210SBarry Smith  * therefore we pack data types manually. This section defines packing routines for the standard data types.
6595fce210SBarry Smith  */
6695fce210SBarry Smith 
6795fce210SBarry Smith #define CPPJoin2_exp(a,b) a ## b
6895fce210SBarry Smith #define CPPJoin2(a,b) CPPJoin2_exp(a,b)
6995fce210SBarry Smith #define CPPJoin3_exp_(a,b,c) a ## b ## _ ## c
7095fce210SBarry Smith #define CPPJoin3_(a,b,c) CPPJoin3_exp_(a,b,c)
7195fce210SBarry Smith 
7295fce210SBarry Smith /* Basic types without addition */
7395fce210SBarry Smith #define DEF_PackNoInit(type)                                            \
7495fce210SBarry Smith   static void CPPJoin2(Pack_,type)(PetscInt n,const PetscInt *idx,const void *unpacked,void *packed) { \
7595fce210SBarry Smith     const type *u = (const type*)unpacked;                              \
7695fce210SBarry Smith     type *p = (type*)packed;                                            \
7795fce210SBarry Smith     PetscInt i;                                                         \
7895fce210SBarry Smith     for (i=0; i<n; i++) p[i] = u[idx[i]];                               \
7995fce210SBarry Smith   }                                                                     \
8095fce210SBarry Smith   static void CPPJoin2(UnpackInsert_,type)(PetscInt n,const PetscInt *idx,void *unpacked,const void *packed) { \
8195fce210SBarry Smith     type *u = (type*)unpacked;                                          \
8295fce210SBarry Smith     const type *p = (const type*)packed;                                \
8395fce210SBarry Smith     PetscInt i;                                                         \
8495fce210SBarry Smith     for (i=0; i<n; i++) u[idx[i]] = p[i];                               \
8595fce210SBarry Smith   }                                                                     \
8695fce210SBarry Smith   static void CPPJoin2(FetchAndInsert_,type)(PetscInt n,const PetscInt *idx,void *unpacked,void *packed) { \
8795fce210SBarry Smith     type *u = (type*)unpacked;                                          \
8895fce210SBarry Smith     type *p = (type*)packed;                                            \
8995fce210SBarry Smith     PetscInt i;                                                         \
9095fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
9195fce210SBarry Smith       PetscInt j = idx[i];                                              \
9295fce210SBarry Smith       type t = u[j];                                                    \
9395fce210SBarry Smith       u[j] = p[i];                                                      \
9495fce210SBarry Smith       p[i] = t;                                                         \
9595fce210SBarry Smith     }                                                                   \
9695fce210SBarry Smith   }
9795fce210SBarry Smith 
9895fce210SBarry Smith /* Basic types defining addition */
9995fce210SBarry Smith #define DEF_PackAddNoInit(type)                                         \
10095fce210SBarry Smith   DEF_PackNoInit(type)                                                  \
10195fce210SBarry Smith   static void CPPJoin2(UnpackAdd_,type)(PetscInt n,const PetscInt *idx,void *unpacked,const void *packed) { \
10295fce210SBarry Smith     type *u = (type*)unpacked;                                          \
10395fce210SBarry Smith     const type *p = (const type*)packed;                                \
10495fce210SBarry Smith     PetscInt i;                                                         \
10595fce210SBarry Smith     for (i=0; i<n; i++) u[idx[i]] += p[i];                              \
10695fce210SBarry Smith   }                                                                     \
10795fce210SBarry Smith   static void CPPJoin2(FetchAndAdd_,type)(PetscInt n,const PetscInt *idx,void *unpacked,void *packed) { \
10895fce210SBarry Smith     type *u = (type*)unpacked;                                          \
10995fce210SBarry Smith     type *p = (type*)packed;                                            \
11095fce210SBarry Smith     PetscInt i;                                                         \
11195fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
11295fce210SBarry Smith       PetscInt j = idx[i];                                              \
11395fce210SBarry Smith       type t = u[j];                                                    \
11495fce210SBarry Smith       u[j] = t + p[i];                                                  \
11595fce210SBarry Smith       p[i] = t;                                                         \
11695fce210SBarry Smith     }                                                                   \
11753deab39SPeter Brune   }                                                                     \
11853deab39SPeter Brune   static void CPPJoin2(UnpackMult_,type)(PetscInt n,const PetscInt *idx,void *unpacked,const void *packed) { \
11953deab39SPeter Brune     type *u = (type*)unpacked;                                          \
12053deab39SPeter Brune     const type *p = (const type*)packed;                                \
12153deab39SPeter Brune     PetscInt i;                                                         \
12253deab39SPeter Brune     for (i=0; i<n; i++) u[idx[i]] *= p[i];                              \
12353deab39SPeter Brune   }                                                                     \
12453deab39SPeter Brune   static void CPPJoin2(FetchAndMult_,type)(PetscInt n,const PetscInt *idx,void *unpacked,void *packed) { \
12553deab39SPeter Brune     type *u = (type*)unpacked;                                          \
12653deab39SPeter Brune     type *p = (type*)packed;                                            \
12753deab39SPeter Brune     PetscInt i;                                                         \
12853deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
12953deab39SPeter Brune       PetscInt j = idx[i];                                              \
13053deab39SPeter Brune       type t = u[j];                                                    \
13153deab39SPeter Brune       u[j] = t * p[i];                                                  \
13253deab39SPeter Brune       p[i] = t;                                                         \
13353deab39SPeter Brune     }                                                                   \
13495fce210SBarry Smith   }
13595fce210SBarry Smith #define DEF_Pack(type) \
13695fce210SBarry Smith   DEF_PackAddNoInit(type)                                               \
13795fce210SBarry Smith   static void CPPJoin2(PackInit_,type)(PetscSFBasicPack link) {         \
13895fce210SBarry Smith     link->Pack = CPPJoin2(Pack_,type);                                  \
13995fce210SBarry Smith     link->UnpackInsert = CPPJoin2(UnpackInsert_,type);                  \
14095fce210SBarry Smith     link->UnpackAdd = CPPJoin2(UnpackAdd_,type);                        \
14153deab39SPeter Brune     link->UnpackMult = CPPJoin2(UnpackMult_,type);                      \
14295fce210SBarry Smith     link->FetchAndInsert = CPPJoin2(FetchAndInsert_,type);              \
14395fce210SBarry Smith     link->FetchAndAdd = CPPJoin2(FetchAndAdd_,type);                    \
14453deab39SPeter Brune     link->FetchAndMult = CPPJoin2(FetchAndMult_,type);                  \
14595fce210SBarry Smith     link->unitbytes = sizeof(type);                                     \
14695fce210SBarry Smith   }
14795fce210SBarry Smith /* Comparable types */
14895fce210SBarry Smith #define DEF_PackCmp(type)                                               \
14995fce210SBarry Smith   DEF_PackAddNoInit(type)                                               \
15095fce210SBarry Smith   static void CPPJoin2(UnpackMax_,type)(PetscInt n,const PetscInt *idx,void *unpacked,const void *packed) { \
15195fce210SBarry Smith     type *u = (type*)unpacked;                                          \
15295fce210SBarry Smith     const type *p = (const type*)packed;                                \
15395fce210SBarry Smith     PetscInt i;                                                         \
15495fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
15595fce210SBarry Smith       type v = u[idx[i]];                                               \
15695fce210SBarry Smith       u[idx[i]] = PetscMax(v,p[i]);                                     \
15795fce210SBarry Smith     }                                                                   \
15895fce210SBarry Smith   }                                                                     \
15995fce210SBarry Smith   static void CPPJoin2(UnpackMin_,type)(PetscInt n,const PetscInt *idx,void *unpacked,const void *packed) { \
16095fce210SBarry Smith     type *u = (type*)unpacked;                                          \
16195fce210SBarry Smith     const type *p = (const type*)packed;                                \
16295fce210SBarry Smith     PetscInt i;                                                         \
16395fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
16495fce210SBarry Smith       type v = u[idx[i]];                                               \
16595fce210SBarry Smith       u[idx[i]] = PetscMin(v,p[i]);                                     \
16695fce210SBarry Smith     }                                                                   \
16795fce210SBarry Smith   }                                                                     \
16895fce210SBarry Smith   static void CPPJoin2(FetchAndMax_,type)(PetscInt n,const PetscInt *idx,void *unpacked,void *packed) { \
16995fce210SBarry Smith     type *u = (type*)unpacked;                                          \
17095fce210SBarry Smith     type *p = (type*)packed;                                            \
17195fce210SBarry Smith     PetscInt i;                                                         \
17295fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
17395fce210SBarry Smith       PetscInt j = idx[i];                                              \
17495fce210SBarry Smith       type v = u[j];                                                    \
17595fce210SBarry Smith       u[j] = PetscMax(v,p[i]);                                          \
17695fce210SBarry Smith       p[i] = v;                                                         \
17795fce210SBarry Smith     }                                                                   \
17895fce210SBarry Smith   }                                                                     \
17995fce210SBarry Smith   static void CPPJoin2(FetchAndMin_,type)(PetscInt n,const PetscInt *idx,void *unpacked,void *packed) { \
18095fce210SBarry Smith     type *u = (type*)unpacked;                                          \
18195fce210SBarry Smith     type *p = (type*)packed;                                            \
18295fce210SBarry Smith     PetscInt i;                                                         \
18395fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
18495fce210SBarry Smith       PetscInt j = idx[i];                                              \
18595fce210SBarry Smith       type v = u[j];                                                    \
18695fce210SBarry Smith       u[j] = PetscMin(v,p[i]);                                          \
18795fce210SBarry Smith       p[i] = v;                                                         \
18895fce210SBarry Smith     }                                                                   \
18995fce210SBarry Smith   }                                                                     \
19095fce210SBarry Smith   static void CPPJoin2(PackInit_,type)(PetscSFBasicPack link) {         \
19195fce210SBarry Smith     link->Pack = CPPJoin2(Pack_,type);                                  \
19295fce210SBarry Smith     link->UnpackInsert = CPPJoin2(UnpackInsert_,type);                  \
19395fce210SBarry Smith     link->UnpackAdd  = CPPJoin2(UnpackAdd_,type);                       \
19495fce210SBarry Smith     link->UnpackMax  = CPPJoin2(UnpackMax_,type);                       \
19595fce210SBarry Smith     link->UnpackMin  = CPPJoin2(UnpackMin_,type);                       \
19653deab39SPeter Brune     link->UnpackMult = CPPJoin2(UnpackMult_,type);                      \
19795fce210SBarry Smith     link->FetchAndInsert = CPPJoin2(FetchAndInsert_,type);              \
19895fce210SBarry Smith     link->FetchAndAdd = CPPJoin2(FetchAndAdd_ ,type);                   \
19995fce210SBarry Smith     link->FetchAndMax = CPPJoin2(FetchAndMax_ ,type);                   \
20095fce210SBarry Smith     link->FetchAndMin = CPPJoin2(FetchAndMin_ ,type);                   \
20153deab39SPeter Brune     link->FetchAndMult = CPPJoin2(FetchAndMult_,type);                  \
20295fce210SBarry Smith     link->unitbytes = sizeof(type);                                     \
20395fce210SBarry Smith   }
20495fce210SBarry Smith 
20553deab39SPeter Brune /* Logical Types */
20653deab39SPeter Brune #define DEF_PackLog(type)                                               \
20753deab39SPeter Brune   static void CPPJoin2(UnpackLAND_,type)(PetscInt n,const PetscInt *idx,void *unpacked,const void *packed) { \
20853deab39SPeter Brune     type *u = (type*)unpacked;                                          \
20953deab39SPeter Brune     const type *p = (const type*)packed;                                \
21053deab39SPeter Brune     PetscInt i;                                                         \
21153deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
21253deab39SPeter Brune       type v = u[idx[i]];                                               \
21353deab39SPeter Brune       u[idx[i]] = v && p[i];                                            \
21453deab39SPeter Brune     }                                                                   \
21553deab39SPeter Brune   }                                                                     \
21653deab39SPeter Brune   static void CPPJoin2(UnpackLOR_,type)(PetscInt n,const PetscInt *idx,void *unpacked,const void *packed) { \
21753deab39SPeter Brune     type *u = (type*)unpacked;                                          \
21853deab39SPeter Brune     const type *p = (const type*)packed;                                \
21953deab39SPeter Brune     PetscInt i;                                                         \
22053deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
22153deab39SPeter Brune       type v = u[idx[i]];                                               \
22253deab39SPeter Brune       u[idx[i]] = v || p[i];                                            \
22353deab39SPeter Brune     }                                                                   \
22453deab39SPeter Brune   }                                                                     \
22553deab39SPeter Brune   static void CPPJoin2(UnpackLXOR_,type)(PetscInt n,const PetscInt *idx,void *unpacked,const void *packed) { \
22653deab39SPeter Brune     type *u = (type*)unpacked;                                          \
22753deab39SPeter Brune     const type *p = (const type*)packed;                                \
22853deab39SPeter Brune     PetscInt i;                                                         \
22953deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
23053deab39SPeter Brune       type v = u[idx[i]];                                               \
23153deab39SPeter Brune       u[idx[i]] = (!v)!=(!p[i]);                                        \
23253deab39SPeter Brune     }                                                                   \
23353deab39SPeter Brune   }                                                                     \
23453deab39SPeter Brune   static void CPPJoin2(FetchAndLAND_,type)(PetscInt n,const PetscInt *idx,void *unpacked,void *packed) { \
23553deab39SPeter Brune     type *u = (type*)unpacked;                                          \
23653deab39SPeter Brune     type *p = (type*)packed;                                            \
23753deab39SPeter Brune     PetscInt i;                                                         \
23853deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
23953deab39SPeter Brune       PetscInt j = idx[i];                                              \
24053deab39SPeter Brune       type v = u[j];                                                    \
24153deab39SPeter Brune       u[j] = v && p[i];                                                 \
24253deab39SPeter Brune       p[i] = v;                                                         \
24353deab39SPeter Brune     }                                                                   \
24453deab39SPeter Brune   }                                                                     \
24553deab39SPeter Brune   static void CPPJoin2(FetchAndLOR_,type)(PetscInt n,const PetscInt *idx,void *unpacked,void *packed) { \
24653deab39SPeter Brune     type *u = (type*)unpacked;                                          \
24753deab39SPeter Brune     type *p = (type*)packed;                                            \
24853deab39SPeter Brune     PetscInt i;                                                         \
24953deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
25053deab39SPeter Brune       PetscInt j = idx[i];                                              \
25153deab39SPeter Brune       type v = u[j];                                                    \
25253deab39SPeter Brune       u[j] = v || p[i];                                                 \
25353deab39SPeter Brune       p[i] = v;                                                         \
25453deab39SPeter Brune     }                                                                   \
25553deab39SPeter Brune   }                                                                     \
25653deab39SPeter Brune   static void CPPJoin2(FetchAndLXOR_,type)(PetscInt n,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   }                                                                     \
26753deab39SPeter Brune   static void CPPJoin2(PackInit_Logical_,type)(PetscSFBasicPack link) { \
26853deab39SPeter Brune     link->UnpackLAND = CPPJoin2(UnpackLAND_,type);                      \
26953deab39SPeter Brune     link->UnpackLOR  = CPPJoin2(UnpackLOR_,type);                       \
27053deab39SPeter Brune     link->UnpackLXOR = CPPJoin2(UnpackLXOR_,type);                      \
27153deab39SPeter Brune     link->FetchAndLAND = CPPJoin2(FetchAndLAND_,type);                  \
27253deab39SPeter Brune     link->FetchAndLOR  = CPPJoin2(FetchAndLOR_,type);                   \
27353deab39SPeter Brune     link->FetchAndLXOR = CPPJoin2(FetchAndLXOR_,type);                  \
27453deab39SPeter Brune   }
27553deab39SPeter Brune 
27653deab39SPeter Brune 
27753deab39SPeter Brune /* Bitwise Types */
27853deab39SPeter Brune #define DEF_PackBit(type)                                               \
27953deab39SPeter Brune   static void CPPJoin2(UnpackBAND_,type)(PetscInt n,const PetscInt *idx,void *unpacked,const void *packed) { \
28053deab39SPeter Brune     type *u = (type*)unpacked;                                          \
28153deab39SPeter Brune     const type *p = (const type*)packed;                                \
28253deab39SPeter Brune     PetscInt i;                                                         \
28353deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
28453deab39SPeter Brune       type v = u[idx[i]];                                               \
28553deab39SPeter Brune       u[idx[i]] = v & p[i];                                             \
28653deab39SPeter Brune     }                                                                   \
28753deab39SPeter Brune   }                                                                     \
28853deab39SPeter Brune   static void CPPJoin2(UnpackBOR_,type)(PetscInt n,const PetscInt *idx,void *unpacked,const void *packed) { \
28953deab39SPeter Brune     type *u = (type*)unpacked;                                          \
29053deab39SPeter Brune     const type *p = (const type*)packed;                                \
29153deab39SPeter Brune     PetscInt i;                                                         \
29253deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
29353deab39SPeter Brune       type v = u[idx[i]];                                               \
29453deab39SPeter Brune       u[idx[i]] = v | p[i];                                             \
29553deab39SPeter Brune     }                                                                   \
29653deab39SPeter Brune   }                                                                     \
29753deab39SPeter Brune   static void CPPJoin2(UnpackBXOR_,type)(PetscInt n,const PetscInt *idx,void *unpacked,const void *packed) { \
29853deab39SPeter Brune     type *u = (type*)unpacked;                                          \
29953deab39SPeter Brune     const type *p = (const type*)packed;                                \
30053deab39SPeter Brune     PetscInt i;                                                         \
30153deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
30253deab39SPeter Brune       type v = u[idx[i]];                                               \
30353deab39SPeter Brune       u[idx[i]] = v^p[i];                                               \
30453deab39SPeter Brune     }                                                                   \
30553deab39SPeter Brune   }                                                                     \
30653deab39SPeter Brune   static void CPPJoin2(FetchAndBAND_,type)(PetscInt n,const PetscInt *idx,void *unpacked,void *packed) { \
30753deab39SPeter Brune     type *u = (type*)unpacked;                                          \
30853deab39SPeter Brune     type *p = (type*)packed;                                            \
30953deab39SPeter Brune     PetscInt i;                                                         \
31053deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
31153deab39SPeter Brune       PetscInt j = idx[i];                                              \
31253deab39SPeter Brune       type v = u[j];                                                    \
31353deab39SPeter Brune       u[j] = v & p[i];                                                  \
31453deab39SPeter Brune       p[i] = v;                                                         \
31553deab39SPeter Brune     }                                                                   \
31653deab39SPeter Brune   }                                                                     \
31753deab39SPeter Brune   static void CPPJoin2(FetchAndBOR_,type)(PetscInt n,const PetscInt *idx,void *unpacked,void *packed) { \
31853deab39SPeter Brune     type *u = (type*)unpacked;                                          \
31953deab39SPeter Brune     type *p = (type*)packed;                                            \
32053deab39SPeter Brune     PetscInt i;                                                         \
32153deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
32253deab39SPeter Brune       PetscInt j = idx[i];                                              \
32353deab39SPeter Brune       type v = u[j];                                                    \
32453deab39SPeter Brune       u[j] = v | p[i];                                                  \
32553deab39SPeter Brune       p[i] = v;                                                         \
32653deab39SPeter Brune     }                                                                   \
32753deab39SPeter Brune   }                                                                     \
32853deab39SPeter Brune   static void CPPJoin2(FetchAndBXOR_,type)(PetscInt n,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   }                                                                     \
33953deab39SPeter Brune   static void CPPJoin2(PackInit_Bitwise_,type)(PetscSFBasicPack link) { \
34053deab39SPeter Brune     link->UnpackBAND = CPPJoin2(UnpackBAND_,type);                      \
34153deab39SPeter Brune     link->UnpackBOR  = CPPJoin2(UnpackBOR_,type);                       \
34253deab39SPeter Brune     link->UnpackBXOR = CPPJoin2(UnpackBXOR_,type);                      \
34353deab39SPeter Brune     link->FetchAndBAND = CPPJoin2(FetchAndBAND_,type);                  \
34453deab39SPeter Brune     link->FetchAndBOR  = CPPJoin2(FetchAndBOR_,type);                   \
34553deab39SPeter Brune     link->FetchAndBXOR = CPPJoin2(FetchAndBXOR_,type);                  \
34653deab39SPeter Brune   }
34753deab39SPeter Brune 
34895fce210SBarry Smith /* Pair types */
34995fce210SBarry Smith #define CPPJoinloc_exp(base,op,t1,t2) base ## op ## loc_ ## t1 ## _ ## t2
35095fce210SBarry Smith #define CPPJoinloc(base,op,t1,t2) CPPJoinloc_exp(base,op,t1,t2)
35195fce210SBarry Smith #define PairType(type1,type2) CPPJoin3_(_pairtype_,type1,type2)
35295fce210SBarry Smith #define DEF_UnpackXloc(type1,type2,locname,op)                              \
35395fce210SBarry Smith   static void CPPJoinloc(Unpack,locname,type1,type2)(PetscInt n,const PetscInt *idx,void *unpacked,const void *packed) { \
35495fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;        \
35595fce210SBarry Smith     const PairType(type1,type2) *p = (const PairType(type1,type2)*)packed; \
35695fce210SBarry Smith     PetscInt i;                                                         \
35795fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
35895fce210SBarry Smith       PetscInt j = idx[i];                                              \
35995fce210SBarry Smith       if (p[i].a op u[j].a) {                                           \
36095fce210SBarry Smith         u[j].a = p[i].a;                                                \
36195fce210SBarry Smith         u[j].b = p[i].b;                                                \
36295fce210SBarry Smith       } else if (u[j].a == p[i].a) {                                    \
36395fce210SBarry Smith         u[j].b = PetscMin(u[j].b,p[i].b);                               \
36495fce210SBarry Smith       }                                                                 \
36595fce210SBarry Smith     }                                                                   \
36695fce210SBarry Smith   }                                                                     \
36795fce210SBarry Smith   static void CPPJoinloc(FetchAnd,locname,type1,type2)(PetscInt n,const PetscInt *idx,void *unpacked,void *packed) { \
36895fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;        \
36995fce210SBarry Smith     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;          \
37095fce210SBarry Smith     PetscInt i;                                                         \
37195fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
37295fce210SBarry Smith       PetscInt j = idx[i];                                              \
37395fce210SBarry Smith       PairType(type1,type2) v;                                          \
37495fce210SBarry Smith       v.a = u[j].a;                                                     \
37595fce210SBarry Smith       v.b = u[j].b;                                                     \
37695fce210SBarry Smith       if (p[i].a op u[j].a) {                                           \
37795fce210SBarry Smith         u[j].a = p[i].a;                                                \
37895fce210SBarry Smith         u[j].b = p[i].b;                                                \
37995fce210SBarry Smith       } else if (u[j].a == p[i].a) {                                    \
38095fce210SBarry Smith         u[j].b = PetscMin(u[j].b,p[i].b);                               \
38195fce210SBarry Smith       }                                                                 \
38295fce210SBarry Smith       p[i].a = v.a;                                                     \
38395fce210SBarry Smith       p[i].b = v.b;                                                     \
38495fce210SBarry Smith     }                                                                   \
38595fce210SBarry Smith   }
38695fce210SBarry Smith #define DEF_PackPair(type1,type2)                                       \
38795fce210SBarry Smith   typedef struct {type1 a; type2 b;} PairType(type1,type2);             \
38895fce210SBarry Smith   static void CPPJoin3_(Pack_,type1,type2)(PetscInt n,const PetscInt *idx,const void *unpacked,void *packed) { \
38995fce210SBarry Smith     const PairType(type1,type2) *u = (const PairType(type1,type2)*)unpacked; \
39095fce210SBarry Smith     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;          \
39195fce210SBarry Smith     PetscInt i;                                                         \
39295fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
39395fce210SBarry Smith       p[i].a = u[idx[i]].a;                                             \
39495fce210SBarry Smith       p[i].b = u[idx[i]].b;                                             \
39595fce210SBarry Smith     }                                                                   \
39695fce210SBarry Smith   }                                                                     \
39795fce210SBarry Smith   static void CPPJoin3_(UnpackInsert_,type1,type2)(PetscInt n,const PetscInt *idx,void *unpacked,const void *packed) { \
39895fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;       \
39995fce210SBarry Smith     const PairType(type1,type2) *p = (const PairType(type1,type2)*)packed; \
40095fce210SBarry Smith     PetscInt i;                                                         \
40195fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
40295fce210SBarry Smith       u[idx[i]].a = p[i].a;                                             \
40395fce210SBarry Smith       u[idx[i]].b = p[i].b;                                             \
40495fce210SBarry Smith     }                                                                   \
40595fce210SBarry Smith   }                                                                     \
40695fce210SBarry Smith   static void CPPJoin3_(UnpackAdd_,type1,type2)(PetscInt n,const PetscInt *idx,void *unpacked,const void *packed) { \
40795fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;       \
40895fce210SBarry Smith     const PairType(type1,type2) *p = (const PairType(type1,type2)*)packed; \
40995fce210SBarry Smith     PetscInt i;                                                         \
41095fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
41195fce210SBarry Smith       u[idx[i]].a += p[i].a;                                            \
41295fce210SBarry Smith       u[idx[i]].b += p[i].b;                                            \
41395fce210SBarry Smith     }                                                                   \
41495fce210SBarry Smith   }                                                                     \
41595fce210SBarry Smith   static void CPPJoin3_(FetchAndInsert_,type1,type2)(PetscInt n,const PetscInt *idx,void *unpacked,void *packed) { \
41695fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;        \
41795fce210SBarry Smith     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;          \
41895fce210SBarry Smith     PetscInt i;                                                         \
41995fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
42095fce210SBarry Smith       PetscInt j = idx[i];                                              \
42195fce210SBarry Smith       PairType(type1,type2) v;                                          \
42295fce210SBarry Smith       v.a = u[j].a;                                                     \
42395fce210SBarry Smith       v.b = u[j].b;                                                     \
42495fce210SBarry Smith       u[j].a = p[i].a;                                                  \
42595fce210SBarry Smith       u[j].b = p[i].b;                                                  \
42695fce210SBarry Smith       p[i].a = v.a;                                                     \
42795fce210SBarry Smith       p[i].b = v.b;                                                     \
42895fce210SBarry Smith     }                                                                   \
42995fce210SBarry Smith   }                                                                     \
43095fce210SBarry Smith   static void FetchAndAdd_ ## type1 ## _ ## type2(PetscInt n,const PetscInt *idx,void *unpacked,void *packed) { \
43195fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;       \
43295fce210SBarry Smith     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;         \
43395fce210SBarry Smith     PetscInt i;                                                         \
43495fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
43595fce210SBarry Smith       PetscInt j = idx[i];                                              \
43695fce210SBarry Smith       PairType(type1,type2) v;                                          \
43795fce210SBarry Smith       v.a = u[j].a;                                                     \
43895fce210SBarry Smith       v.b = u[j].b;                                                     \
43995fce210SBarry Smith       u[j].a = v.a + p[i].a;                                            \
44095fce210SBarry Smith       u[j].b = v.b + p[i].b;                                            \
44195fce210SBarry Smith       p[i].a = v.a;                                                     \
44295fce210SBarry Smith       p[i].b = v.b;                                                     \
44395fce210SBarry Smith     }                                                                   \
44495fce210SBarry Smith   }                                                                     \
44595fce210SBarry Smith   DEF_UnpackXloc(type1,type2,Max,>)                                     \
44695fce210SBarry Smith   DEF_UnpackXloc(type1,type2,Min,<)                                     \
44795fce210SBarry Smith   static void CPPJoin3_(PackInit_,type1,type2)(PetscSFBasicPack link) { \
44895fce210SBarry Smith     link->Pack = CPPJoin3_(Pack_,type1,type2);                          \
44995fce210SBarry Smith     link->UnpackInsert = CPPJoin3_(UnpackInsert_,type1,type2);          \
45095fce210SBarry Smith     link->UnpackAdd = CPPJoin3_(UnpackAdd_,type1,type2);                \
45195fce210SBarry Smith     link->UnpackMaxloc = CPPJoin3_(UnpackMaxloc_,type1,type2);          \
45295fce210SBarry Smith     link->UnpackMinloc = CPPJoin3_(UnpackMinloc_,type1,type2);          \
45395fce210SBarry Smith     link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,type1,type2);      \
45495fce210SBarry Smith     link->FetchAndAdd = CPPJoin3_(FetchAndAdd_,type1,type2);            \
45595fce210SBarry Smith     link->FetchAndMaxloc = CPPJoin3_(FetchAndMaxloc_,type1,type2);      \
45695fce210SBarry Smith     link->FetchAndMinloc = CPPJoin3_(FetchAndMinloc_,type1,type2);      \
45795fce210SBarry Smith     link->unitbytes = sizeof(PairType(type1,type2));                    \
45895fce210SBarry Smith   }
45995fce210SBarry Smith 
46095fce210SBarry Smith /* Currently only dumb blocks of data */
46195fce210SBarry Smith #define BlockType(unit,count) CPPJoin3_(_blocktype_,unit,count)
46295fce210SBarry Smith #define DEF_Block(unit,count)                                           \
46395fce210SBarry Smith   typedef struct {unit v[count];} BlockType(unit,count);                \
46495fce210SBarry Smith   DEF_PackNoInit(BlockType(unit,count))                                 \
46595fce210SBarry Smith   static void CPPJoin3_(PackInit_block_,unit,count)(PetscSFBasicPack link) { \
46695fce210SBarry Smith     link->Pack = CPPJoin2(Pack_,BlockType(unit,count));                 \
46795fce210SBarry Smith     link->UnpackInsert = CPPJoin2(UnpackInsert_,BlockType(unit,count)); \
46895fce210SBarry Smith     link->FetchAndInsert = CPPJoin2(FetchAndInsert_,BlockType(unit,count)); \
46995fce210SBarry Smith     link->unitbytes = sizeof(BlockType(unit,count));                    \
47095fce210SBarry Smith   }
47195fce210SBarry Smith 
47295fce210SBarry Smith DEF_PackCmp(int)
47353deab39SPeter Brune DEF_PackBit(int)
47453deab39SPeter Brune DEF_PackLog(int)
47595fce210SBarry Smith DEF_PackCmp(PetscInt)
47653deab39SPeter Brune DEF_PackBit(PetscInt)
47753deab39SPeter Brune DEF_PackLog(PetscInt)
47895fce210SBarry Smith DEF_PackCmp(PetscReal)
47953deab39SPeter Brune DEF_PackLog(PetscReal)
48095fce210SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
48195fce210SBarry Smith DEF_Pack(PetscComplex)
48295fce210SBarry Smith #endif
48395fce210SBarry Smith DEF_PackPair(int,int)
48495fce210SBarry Smith DEF_PackPair(PetscInt,PetscInt)
48595fce210SBarry Smith DEF_Block(int,1)
48695fce210SBarry Smith DEF_Block(int,2)
48795fce210SBarry Smith DEF_Block(int,3)
48895fce210SBarry Smith DEF_Block(int,4)
48995fce210SBarry Smith DEF_Block(int,5)
49095fce210SBarry Smith DEF_Block(int,6)
49195fce210SBarry Smith DEF_Block(int,7)
49295fce210SBarry Smith DEF_Block(int,8)
49395fce210SBarry Smith 
49495fce210SBarry Smith #undef __FUNCT__
49595fce210SBarry Smith #define __FUNCT__ "PetscSFSetUp_Basic"
49695fce210SBarry Smith static PetscErrorCode PetscSFSetUp_Basic(PetscSF sf)
49795fce210SBarry Smith {
49895fce210SBarry Smith   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
49995fce210SBarry Smith   PetscErrorCode ierr;
50095fce210SBarry Smith   PetscInt *rlengths,*ilengths,i;
50195fce210SBarry Smith   MPI_Comm comm;
50295fce210SBarry Smith   MPI_Request *rootreqs,*leafreqs;
50395fce210SBarry Smith 
50495fce210SBarry Smith   PetscFunctionBegin;
50595fce210SBarry Smith   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
50695fce210SBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)sf,&bas->tag);CHKERRQ(ierr);
50795fce210SBarry Smith   /*
50895fce210SBarry Smith    * Inform roots about how many leaves and from which ranks
50995fce210SBarry Smith    */
51095fce210SBarry Smith   ierr = PetscMalloc(sf->nranks*sizeof(PetscInt),&rlengths);CHKERRQ(ierr);
51195fce210SBarry Smith   /* Determine number, sending ranks, and length of incoming  */
51295fce210SBarry Smith   for (i=0; i<sf->nranks; i++) {
51395fce210SBarry Smith     rlengths[i] = sf->roffset[i+1] - sf->roffset[i]; /* Number of roots referenced by my leaves; for rank sf->ranks[i] */
51495fce210SBarry Smith   }
51595fce210SBarry Smith   ierr = PetscCommBuildTwoSided(comm,1,MPIU_INT,sf->nranks,sf->ranks,rlengths,&bas->niranks,&bas->iranks,(void**)&ilengths);CHKERRQ(ierr);
51695fce210SBarry Smith   ierr = PetscFree(rlengths);CHKERRQ(ierr);
51795fce210SBarry Smith 
51895fce210SBarry Smith   /* Send leaf identities to roots */
51995fce210SBarry Smith   for (i=0,bas->itotal=0; i<bas->niranks; i++) bas->itotal += ilengths[i];
520*dcca6d9dSJed Brown   ierr = PetscMalloc2(bas->niranks+1,&bas->ioffset,bas->itotal,&bas->irootloc);CHKERRQ(ierr);
521*dcca6d9dSJed Brown   ierr = PetscMalloc2(bas->niranks,&rootreqs,sf->nranks,&leafreqs);CHKERRQ(ierr);
52295fce210SBarry Smith   bas->ioffset[0] = 0;
52395fce210SBarry Smith   for (i=0; i<bas->niranks; i++) {
52495fce210SBarry Smith     bas->ioffset[i+1] = bas->ioffset[i] + ilengths[i];
52595fce210SBarry Smith     ierr = MPI_Irecv(bas->irootloc+bas->ioffset[i],ilengths[i],MPIU_INT,bas->iranks[i],bas->tag,comm,&rootreqs[i]);CHKERRQ(ierr);
52695fce210SBarry Smith   }
52795fce210SBarry Smith   for (i=0; i<sf->nranks; i++) {
52895fce210SBarry Smith     PetscMPIInt npoints;
52995fce210SBarry Smith     ierr = PetscMPIIntCast(sf->roffset[i+1] - sf->roffset[i],&npoints);CHKERRQ(ierr);
53095fce210SBarry Smith     ierr = MPI_Isend(sf->rremote+sf->roffset[i],npoints,MPIU_INT,sf->ranks[i],bas->tag,comm,&leafreqs[i]);CHKERRQ(ierr);
53195fce210SBarry Smith   }
532eb9baa12SBarry Smith   ierr = MPI_Waitall(bas->niranks,rootreqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
533eb9baa12SBarry Smith   ierr = MPI_Waitall(sf->nranks,leafreqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
53495fce210SBarry Smith   ierr = PetscFree(ilengths);CHKERRQ(ierr);
535eb9baa12SBarry Smith   ierr = PetscFree2(rootreqs,leafreqs);CHKERRQ(ierr);
53695fce210SBarry Smith   PetscFunctionReturn(0);
53795fce210SBarry Smith }
53895fce210SBarry Smith 
53995fce210SBarry Smith #undef __FUNCT__
54095fce210SBarry Smith #define __FUNCT__ "PetscSFBasicPackTypeSetup"
54195fce210SBarry Smith static PetscErrorCode PetscSFBasicPackTypeSetup(PetscSFBasicPack link,MPI_Datatype unit)
54295fce210SBarry Smith {
54395fce210SBarry Smith   PetscErrorCode ierr;
54495fce210SBarry Smith   PetscBool      isInt,isPetscInt,isPetscReal,is2Int,is2PetscInt;
54595fce210SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
54695fce210SBarry Smith   PetscBool isPetscComplex;
54795fce210SBarry Smith #endif
54895fce210SBarry Smith 
54995fce210SBarry Smith   PetscFunctionBegin;
55095fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPI_INT,&isInt);CHKERRQ(ierr);
55195fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPIU_INT,&isPetscInt);CHKERRQ(ierr);
55295fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPIU_REAL,&isPetscReal);CHKERRQ(ierr);
55395fce210SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
55495fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPIU_COMPLEX,&isPetscComplex);CHKERRQ(ierr);
55595fce210SBarry Smith #endif
55695fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPI_2INT,&is2Int);CHKERRQ(ierr);
55795fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPIU_2INT,&is2PetscInt);CHKERRQ(ierr);
55853deab39SPeter Brune   if (isInt) {PackInit_int(link); PackInit_Logical_int(link); PackInit_Bitwise_int(link);}
55953deab39SPeter Brune   else if (isPetscInt) {PackInit_PetscInt(link); PackInit_Logical_PetscInt(link); PackInit_Bitwise_PetscInt(link);}
56053deab39SPeter Brune   else if (isPetscReal) {PackInit_PetscReal(link); PackInit_Logical_PetscReal(link);}
56195fce210SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
56295fce210SBarry Smith   else if (isPetscComplex) PackInit_PetscComplex(link);
56395fce210SBarry Smith #endif
56495fce210SBarry Smith   else if (is2Int) PackInit_int_int(link);
56595fce210SBarry Smith   else if (is2PetscInt) PackInit_PetscInt_PetscInt(link);
56695fce210SBarry Smith   else {
56795fce210SBarry Smith     PetscMPIInt bytes;
56895fce210SBarry Smith     ierr = MPI_Type_size(unit,&bytes);CHKERRQ(ierr);
56995fce210SBarry Smith     if (bytes % sizeof(int)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for type size not divisible by %D",sizeof(int));
57095fce210SBarry Smith     switch (bytes / sizeof(int)) {
57195fce210SBarry Smith     case 1: PackInit_block_int_1(link); break;
57295fce210SBarry Smith     case 2: PackInit_block_int_2(link); break;
57395fce210SBarry Smith     case 3: PackInit_block_int_3(link); break;
57495fce210SBarry Smith     case 4: PackInit_block_int_4(link); break;
57595fce210SBarry Smith     case 5: PackInit_block_int_5(link); break;
57695fce210SBarry Smith     case 6: PackInit_block_int_6(link); break;
57795fce210SBarry Smith     case 7: PackInit_block_int_7(link); break;
57895fce210SBarry Smith     case 8: PackInit_block_int_8(link); break;
57995fce210SBarry Smith     default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for arbitrary block sizes");
58095fce210SBarry Smith     }
58195fce210SBarry Smith   }
58295fce210SBarry Smith   ierr = MPI_Type_dup(unit,&link->unit);CHKERRQ(ierr);
58395fce210SBarry Smith   PetscFunctionReturn(0);
58495fce210SBarry Smith }
58595fce210SBarry Smith 
58695fce210SBarry Smith #undef __FUNCT__
58795fce210SBarry Smith #define __FUNCT__ "PetscSFBasicPackGetUnpackOp"
58895fce210SBarry Smith static PetscErrorCode PetscSFBasicPackGetUnpackOp(PetscSF sf,PetscSFBasicPack link,MPI_Op op,void (**UnpackOp)(PetscInt,const PetscInt*,void*,const void*))
58995fce210SBarry Smith {
59095fce210SBarry Smith   PetscFunctionBegin;
59195fce210SBarry Smith   *UnpackOp = NULL;
5928bfbc91cSJed Brown   if (op == MPIU_REPLACE) *UnpackOp = link->UnpackInsert;
59395fce210SBarry Smith   else if (op == MPI_SUM || op == MPIU_SUM) *UnpackOp = link->UnpackAdd;
59453deab39SPeter Brune   else if (op == MPI_PROD) *UnpackOp = link->UnpackMult;
59595fce210SBarry Smith   else if (op == MPI_MAX || op == MPIU_MAX) *UnpackOp = link->UnpackMax;
59695fce210SBarry Smith   else if (op == MPI_MIN || op == MPIU_MIN) *UnpackOp = link->UnpackMin;
59753deab39SPeter Brune   else if (op == MPI_LAND) *UnpackOp = link->UnpackLAND;
59853deab39SPeter Brune   else if (op == MPI_BAND) *UnpackOp = link->UnpackBAND;
59953deab39SPeter Brune   else if (op == MPI_LOR) *UnpackOp = link->UnpackLOR;
60053deab39SPeter Brune   else if (op == MPI_BOR) *UnpackOp = link->UnpackBOR;
60153deab39SPeter Brune   else if (op == MPI_LXOR) *UnpackOp = link->UnpackLXOR;
60253deab39SPeter Brune   else if (op == MPI_BXOR) *UnpackOp = link->UnpackBXOR;
60395fce210SBarry Smith   else if (op == MPI_MAXLOC) *UnpackOp = link->UnpackMaxloc;
60495fce210SBarry Smith   else if (op == MPI_MINLOC) *UnpackOp = link->UnpackMinloc;
60595fce210SBarry Smith   else SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_SUP,"No support for MPI_Op");
60695fce210SBarry Smith   PetscFunctionReturn(0);
60795fce210SBarry Smith }
60895fce210SBarry Smith #undef __FUNCT__
60995fce210SBarry Smith #define __FUNCT__ "PetscSFBasicPackGetFetchAndOp"
61095fce210SBarry Smith static PetscErrorCode PetscSFBasicPackGetFetchAndOp(PetscSF sf,PetscSFBasicPack link,MPI_Op op,void (**FetchAndOp)(PetscInt,const PetscInt*,void*,void*))
61195fce210SBarry Smith {
61295fce210SBarry Smith   PetscFunctionBegin;
61395fce210SBarry Smith   *FetchAndOp = NULL;
6148bfbc91cSJed Brown   if (op == MPIU_REPLACE) *FetchAndOp = link->FetchAndInsert;
61595fce210SBarry Smith   else if (op == MPI_SUM || op == MPIU_SUM) *FetchAndOp = link->FetchAndAdd;
61695fce210SBarry Smith   else if (op == MPI_MAX || op == MPIU_MAX) *FetchAndOp = link->FetchAndMax;
61795fce210SBarry Smith   else if (op == MPI_MIN || op == MPIU_MIN) *FetchAndOp = link->FetchAndMin;
61895fce210SBarry Smith   else if (op == MPI_MAXLOC) *FetchAndOp = link->FetchAndMaxloc;
61995fce210SBarry Smith   else if (op == MPI_MINLOC) *FetchAndOp = link->FetchAndMinloc;
62053deab39SPeter Brune   else if (op == MPI_PROD)   *FetchAndOp = link->FetchAndMult;
62153deab39SPeter Brune   else if (op == MPI_LAND)   *FetchAndOp = link->FetchAndLAND;
62253deab39SPeter Brune   else if (op == MPI_BAND)   *FetchAndOp = link->FetchAndBAND;
62353deab39SPeter Brune   else if (op == MPI_LOR)    *FetchAndOp = link->FetchAndLOR;
62453deab39SPeter Brune   else if (op == MPI_BOR)    *FetchAndOp = link->FetchAndBOR;
62553deab39SPeter Brune   else if (op == MPI_LXOR)   *FetchAndOp = link->FetchAndLXOR;
62653deab39SPeter Brune   else if (op == MPI_BXOR)   *FetchAndOp = link->FetchAndBXOR;
62795fce210SBarry Smith   else SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_SUP,"No support for MPI_Op");
62895fce210SBarry Smith   PetscFunctionReturn(0);
62995fce210SBarry Smith }
63095fce210SBarry Smith 
63195fce210SBarry Smith #undef __FUNCT__
63295fce210SBarry Smith #define __FUNCT__ "PetscSFBasicPackGetReqs"
63395fce210SBarry Smith static PetscErrorCode PetscSFBasicPackGetReqs(PetscSF sf,PetscSFBasicPack link,MPI_Request **rootreqs,MPI_Request **leafreqs)
63495fce210SBarry Smith {
63595fce210SBarry Smith   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
63695fce210SBarry Smith 
63795fce210SBarry Smith   PetscFunctionBegin;
63895fce210SBarry Smith   if (rootreqs) *rootreqs = link->requests;
63995fce210SBarry Smith   if (leafreqs) *leafreqs = link->requests + bas->niranks;
64095fce210SBarry Smith   PetscFunctionReturn(0);
64195fce210SBarry Smith }
64295fce210SBarry Smith 
64395fce210SBarry Smith #undef __FUNCT__
64495fce210SBarry Smith #define __FUNCT__ "PetscSFBasicPackWaitall"
64595fce210SBarry Smith static PetscErrorCode PetscSFBasicPackWaitall(PetscSF sf,PetscSFBasicPack link)
64695fce210SBarry Smith {
64795fce210SBarry Smith   PetscSF_Basic  *bas = (PetscSF_Basic*)sf->data;
64895fce210SBarry Smith   PetscErrorCode ierr;
64995fce210SBarry Smith 
65095fce210SBarry Smith   PetscFunctionBegin;
65195fce210SBarry Smith   ierr = MPI_Waitall(bas->niranks+sf->nranks,link->requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
65295fce210SBarry Smith   PetscFunctionReturn(0);
65395fce210SBarry Smith }
65495fce210SBarry Smith 
65595fce210SBarry Smith #undef __FUNCT__
65695fce210SBarry Smith #define __FUNCT__ "PetscSFBasicGetRootInfo"
65795fce210SBarry Smith static PetscErrorCode PetscSFBasicGetRootInfo(PetscSF sf,PetscInt *nrootranks,const PetscMPIInt **rootranks,const PetscInt **rootoffset,const PetscInt **rootloc)
65895fce210SBarry Smith {
65995fce210SBarry Smith   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
66095fce210SBarry Smith 
66195fce210SBarry Smith   PetscFunctionBegin;
66295fce210SBarry Smith   if (nrootranks) *nrootranks = bas->niranks;
66395fce210SBarry Smith   if (rootranks)  *rootranks  = bas->iranks;
66495fce210SBarry Smith   if (rootoffset) *rootoffset = bas->ioffset;
66595fce210SBarry Smith   if (rootloc)    *rootloc    = bas->irootloc;
66695fce210SBarry Smith   PetscFunctionReturn(0);
66795fce210SBarry Smith }
66895fce210SBarry Smith 
66995fce210SBarry Smith #undef __FUNCT__
67095fce210SBarry Smith #define __FUNCT__ "PetscSFBasicGetLeafInfo"
67195fce210SBarry Smith static PetscErrorCode PetscSFBasicGetLeafInfo(PetscSF sf,PetscInt *nleafranks,const PetscMPIInt **leafranks,const PetscInt **leafoffset,const PetscInt **leafloc)
67295fce210SBarry Smith {
67395fce210SBarry Smith   PetscFunctionBegin;
67495fce210SBarry Smith   if (nleafranks) *nleafranks = sf->nranks;
67595fce210SBarry Smith   if (leafranks)  *leafranks  = sf->ranks;
67695fce210SBarry Smith   if (leafoffset) *leafoffset = sf->roffset;
67795fce210SBarry Smith   if (leafloc)    *leafloc    = sf->rmine;
67895fce210SBarry Smith   PetscFunctionReturn(0);
67995fce210SBarry Smith }
68095fce210SBarry Smith 
68195fce210SBarry Smith #undef __FUNCT__
68295fce210SBarry Smith #define __FUNCT__ "PetscSFBasicGetPack"
68395fce210SBarry Smith static PetscErrorCode PetscSFBasicGetPack(PetscSF sf,MPI_Datatype unit,const void *key,PetscSFBasicPack *mylink)
68495fce210SBarry Smith {
68595fce210SBarry Smith   PetscSF_Basic    *bas = (PetscSF_Basic*)sf->data;
68695fce210SBarry Smith   PetscErrorCode   ierr;
68795fce210SBarry Smith   PetscSFBasicPack link,*p;
68895fce210SBarry Smith   PetscInt         nrootranks,nleafranks;
68995fce210SBarry Smith   const PetscInt   *rootoffset,*leafoffset;
69095fce210SBarry Smith 
69195fce210SBarry Smith   PetscFunctionBegin;
69295fce210SBarry Smith   /* Look for types in cache */
69395fce210SBarry Smith   for (p=&bas->avail; (link=*p); p=&link->next) {
69495fce210SBarry Smith     PetscBool match;
69595fce210SBarry Smith     ierr = MPIPetsc_Type_compare(unit,link->unit,&match);CHKERRQ(ierr);
69695fce210SBarry Smith     if (match) {
69795fce210SBarry Smith       *p = link->next;          /* Remove from available list */
69895fce210SBarry Smith       goto found;
69995fce210SBarry Smith     }
70095fce210SBarry Smith   }
70195fce210SBarry Smith 
70295fce210SBarry Smith   /* Create new composite types for each send rank */
70395fce210SBarry Smith   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,NULL,&rootoffset,NULL);CHKERRQ(ierr);
70495fce210SBarry Smith   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,NULL,&leafoffset,NULL);CHKERRQ(ierr);
70595fce210SBarry Smith   ierr = PetscNew(struct _n_PetscSFBasicPack,&link);CHKERRQ(ierr);
70695fce210SBarry Smith   ierr = PetscSFBasicPackTypeSetup(link,unit);CHKERRQ(ierr);
707*dcca6d9dSJed Brown   ierr = PetscMalloc2(rootoffset[nrootranks]*link->unitbytes,&link->root,leafoffset[nleafranks]*link->unitbytes,&link->leaf);CHKERRQ(ierr);
70895fce210SBarry Smith   ierr = PetscMalloc((nrootranks+nleafranks)*sizeof(MPI_Request),&link->requests);CHKERRQ(ierr);
70995fce210SBarry Smith 
71095fce210SBarry Smith found:
71195fce210SBarry Smith   link->key  = key;
71295fce210SBarry Smith   link->next = bas->inuse;
71395fce210SBarry Smith   bas->inuse = link;
71495fce210SBarry Smith 
71595fce210SBarry Smith   *mylink = link;
71695fce210SBarry Smith   PetscFunctionReturn(0);
71795fce210SBarry Smith }
71895fce210SBarry Smith 
71995fce210SBarry Smith #undef __FUNCT__
72095fce210SBarry Smith #define __FUNCT__ "PetscSFBasicGetPackInUse"
72195fce210SBarry Smith static PetscErrorCode PetscSFBasicGetPackInUse(PetscSF sf,MPI_Datatype unit,const void *key,PetscCopyMode cmode,PetscSFBasicPack *mylink)
72295fce210SBarry Smith {
72395fce210SBarry Smith   PetscSF_Basic    *bas = (PetscSF_Basic*)sf->data;
72495fce210SBarry Smith   PetscErrorCode   ierr;
72595fce210SBarry Smith   PetscSFBasicPack link,*p;
72695fce210SBarry Smith 
72795fce210SBarry Smith   PetscFunctionBegin;
72895fce210SBarry Smith   /* Look for types in cache */
72995fce210SBarry Smith   for (p=&bas->inuse; (link=*p); p=&link->next) {
73095fce210SBarry Smith     PetscBool match;
73195fce210SBarry Smith     ierr = MPIPetsc_Type_compare(unit,link->unit,&match);CHKERRQ(ierr);
73295fce210SBarry Smith     if (match) {
73395fce210SBarry Smith       switch (cmode) {
73495fce210SBarry Smith       case PETSC_OWN_POINTER: *p = link->next; break; /* Remove from inuse list */
73595fce210SBarry Smith       case PETSC_USE_POINTER: break;
73695fce210SBarry Smith       default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"invalid cmode");
73795fce210SBarry Smith       }
73895fce210SBarry Smith       *mylink = link;
73995fce210SBarry Smith       PetscFunctionReturn(0);
74095fce210SBarry Smith     }
74195fce210SBarry Smith   }
74295fce210SBarry Smith   SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Could not find pack");
74395fce210SBarry Smith   PetscFunctionReturn(0);
74495fce210SBarry Smith }
74595fce210SBarry Smith 
74695fce210SBarry Smith #undef __FUNCT__
74795fce210SBarry Smith #define __FUNCT__ "PetscSFBasicReclaimPack"
74895fce210SBarry Smith static PetscErrorCode PetscSFBasicReclaimPack(PetscSF sf,PetscSFBasicPack *link)
74995fce210SBarry Smith {
75095fce210SBarry Smith   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
75195fce210SBarry Smith 
75295fce210SBarry Smith   PetscFunctionBegin;
75395fce210SBarry Smith   (*link)->key  = NULL;
75495fce210SBarry Smith   (*link)->next = bas->avail;
75595fce210SBarry Smith   bas->avail    = *link;
75695fce210SBarry Smith   *link         = NULL;
75795fce210SBarry Smith   PetscFunctionReturn(0);
75895fce210SBarry Smith }
75995fce210SBarry Smith 
76095fce210SBarry Smith #undef __FUNCT__
76195fce210SBarry Smith #define __FUNCT__ "PetscSFSetFromOptions_Basic"
76295fce210SBarry Smith static PetscErrorCode PetscSFSetFromOptions_Basic(PetscSF sf)
76395fce210SBarry Smith {
76495fce210SBarry Smith   PetscErrorCode ierr;
76595fce210SBarry Smith 
76695fce210SBarry Smith   PetscFunctionBegin;
76795fce210SBarry Smith   ierr = PetscOptionsHead("PetscSF Basic options");CHKERRQ(ierr);
76895fce210SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
76995fce210SBarry Smith   PetscFunctionReturn(0);
77095fce210SBarry Smith }
77195fce210SBarry Smith 
77295fce210SBarry Smith #undef __FUNCT__
77395fce210SBarry Smith #define __FUNCT__ "PetscSFReset_Basic"
77495fce210SBarry Smith static PetscErrorCode PetscSFReset_Basic(PetscSF sf)
77595fce210SBarry Smith {
77695fce210SBarry Smith   PetscSF_Basic    *bas = (PetscSF_Basic*)sf->data;
77795fce210SBarry Smith   PetscErrorCode   ierr;
77895fce210SBarry Smith   PetscSFBasicPack link,next;
77995fce210SBarry Smith 
78095fce210SBarry Smith   PetscFunctionBegin;
78195fce210SBarry Smith   ierr = PetscFree(bas->iranks);CHKERRQ(ierr);
78295fce210SBarry Smith   ierr = PetscFree2(bas->ioffset,bas->irootloc);CHKERRQ(ierr);
78395fce210SBarry Smith   if (bas->inuse) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Outstanding operation has not been completed");
78495fce210SBarry Smith   for (link=bas->avail; link; link=next) {
78595fce210SBarry Smith     next = link->next;
786bb86044dSJed Brown #if defined(PETSC_HAVE_MPI_TYPE_DUP)
78795fce210SBarry Smith     ierr = MPI_Type_free(&link->unit);CHKERRQ(ierr);
788bb86044dSJed Brown #endif
78995fce210SBarry Smith     ierr = PetscFree2(link->root,link->leaf);CHKERRQ(ierr);
79095fce210SBarry Smith     ierr = PetscFree(link->requests);CHKERRQ(ierr);
79195fce210SBarry Smith     ierr = PetscFree(link);CHKERRQ(ierr);
79295fce210SBarry Smith   }
79395fce210SBarry Smith   bas->avail = NULL;
79495fce210SBarry Smith   PetscFunctionReturn(0);
79595fce210SBarry Smith }
79695fce210SBarry Smith 
79795fce210SBarry Smith #undef __FUNCT__
79895fce210SBarry Smith #define __FUNCT__ "PetscSFDestroy_Basic"
79995fce210SBarry Smith static PetscErrorCode PetscSFDestroy_Basic(PetscSF sf)
80095fce210SBarry Smith {
80195fce210SBarry Smith   PetscErrorCode ierr;
80295fce210SBarry Smith 
80395fce210SBarry Smith   PetscFunctionBegin;
80495fce210SBarry Smith   ierr = PetscSFReset_Basic(sf);CHKERRQ(ierr);
80595fce210SBarry Smith   ierr = PetscFree(sf->data);CHKERRQ(ierr);
80695fce210SBarry Smith   PetscFunctionReturn(0);
80795fce210SBarry Smith }
80895fce210SBarry Smith 
80995fce210SBarry Smith #undef __FUNCT__
81095fce210SBarry Smith #define __FUNCT__ "PetscSFView_Basic"
81195fce210SBarry Smith static PetscErrorCode PetscSFView_Basic(PetscSF sf,PetscViewer viewer)
81295fce210SBarry Smith {
81395fce210SBarry Smith   /* PetscSF_Basic *bas = (PetscSF_Basic*)sf->data; */
81495fce210SBarry Smith   PetscErrorCode ierr;
81595fce210SBarry Smith   PetscBool      iascii;
81695fce210SBarry Smith 
81795fce210SBarry Smith   PetscFunctionBegin;
81895fce210SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
81995fce210SBarry Smith   if (iascii) {
82095fce210SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  sort=%s\n",sf->rankorder ? "rank-order" : "unordered");CHKERRQ(ierr);
82195fce210SBarry Smith   }
82295fce210SBarry Smith   PetscFunctionReturn(0);
82395fce210SBarry Smith }
82495fce210SBarry Smith 
82595fce210SBarry Smith #undef __FUNCT__
82695fce210SBarry Smith #define __FUNCT__ "PetscSFBcastBegin_Basic"
82795fce210SBarry Smith /* Send from roots to leaves */
82895fce210SBarry Smith static PetscErrorCode PetscSFBcastBegin_Basic(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
82995fce210SBarry Smith {
83095fce210SBarry Smith   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;
83195fce210SBarry Smith   PetscErrorCode    ierr;
83295fce210SBarry Smith   PetscSFBasicPack  link;
83395fce210SBarry Smith   PetscInt          i,nrootranks,nleafranks;
83495fce210SBarry Smith   const PetscInt    *rootoffset,*leafoffset,*rootloc,*leafloc;
83595fce210SBarry Smith   const PetscMPIInt *rootranks,*leafranks;
83695fce210SBarry Smith   MPI_Request       *rootreqs,*leafreqs;
83795fce210SBarry Smith   size_t            unitbytes;
83895fce210SBarry Smith 
83995fce210SBarry Smith   PetscFunctionBegin;
84095fce210SBarry Smith   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,&rootranks,&rootoffset,&rootloc);CHKERRQ(ierr);
84195fce210SBarry Smith   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,&leafranks,&leafoffset,&leafloc);CHKERRQ(ierr);
84295fce210SBarry Smith   ierr = PetscSFBasicGetPack(sf,unit,rootdata,&link);CHKERRQ(ierr);
84395fce210SBarry Smith 
84495fce210SBarry Smith   unitbytes = link->unitbytes;
84595fce210SBarry Smith 
84695fce210SBarry Smith   ierr = PetscSFBasicPackGetReqs(sf,link,&rootreqs,&leafreqs);CHKERRQ(ierr);
84795fce210SBarry Smith   /* Eagerly post leaf receives */
84895fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
84995fce210SBarry Smith     PetscMPIInt n = leafoffset[i+1] - leafoffset[i];
85095fce210SBarry Smith     ierr = MPI_Irecv(link->leaf+leafoffset[i]*unitbytes,n,unit,leafranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&leafreqs[i]);CHKERRQ(ierr);
85195fce210SBarry Smith   }
85295fce210SBarry Smith   /* Pack and send root data */
85395fce210SBarry Smith   for (i=0; i<nrootranks; i++) {
85495fce210SBarry Smith     PetscMPIInt n          = rootoffset[i+1] - rootoffset[i];
85595fce210SBarry Smith     void        *packstart = link->root+rootoffset[i]*unitbytes;
85695fce210SBarry Smith     (*link->Pack)(n,rootloc+rootoffset[i],rootdata,packstart);
85795fce210SBarry Smith     ierr = MPI_Isend(packstart,n,unit,rootranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&rootreqs[i]);CHKERRQ(ierr);
85895fce210SBarry Smith   }
85995fce210SBarry Smith   PetscFunctionReturn(0);
86095fce210SBarry Smith }
86195fce210SBarry Smith 
86295fce210SBarry Smith #undef __FUNCT__
86395fce210SBarry Smith #define __FUNCT__ "PetscSFBcastEnd_Basic"
86495fce210SBarry Smith PetscErrorCode PetscSFBcastEnd_Basic(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
86595fce210SBarry Smith {
86695fce210SBarry Smith   PetscErrorCode   ierr;
86795fce210SBarry Smith   PetscSFBasicPack link;
86895fce210SBarry Smith   PetscInt         i,nleafranks;
86995fce210SBarry Smith   const PetscInt   *leafoffset,*leafloc;
87095fce210SBarry Smith 
87195fce210SBarry Smith   PetscFunctionBegin;
87295fce210SBarry Smith   ierr = PetscSFBasicGetPackInUse(sf,unit,rootdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr);
87395fce210SBarry Smith   ierr = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
87495fce210SBarry Smith   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,NULL,&leafoffset,&leafloc);CHKERRQ(ierr);
87595fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
87695fce210SBarry Smith     PetscMPIInt n          = leafoffset[i+1] - leafoffset[i];
87795fce210SBarry Smith     const void  *packstart = link->leaf+leafoffset[i]*link->unitbytes;
87895fce210SBarry Smith     (*link->UnpackInsert)(n,leafloc+leafoffset[i],leafdata,packstart);
87995fce210SBarry Smith   }
88095fce210SBarry Smith   ierr = PetscSFBasicReclaimPack(sf,&link);CHKERRQ(ierr);
88195fce210SBarry Smith   PetscFunctionReturn(0);
88295fce210SBarry Smith }
88395fce210SBarry Smith 
88495fce210SBarry Smith #undef __FUNCT__
88595fce210SBarry Smith #define __FUNCT__ "PetscSFReduceBegin_Basic"
88695fce210SBarry Smith /* leaf -> root with reduction */
88795fce210SBarry Smith PetscErrorCode PetscSFReduceBegin_Basic(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
88895fce210SBarry Smith {
88995fce210SBarry Smith   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;
89095fce210SBarry Smith   PetscSFBasicPack  link;
89195fce210SBarry Smith   PetscErrorCode    ierr;
89295fce210SBarry Smith   PetscInt          i,nrootranks,nleafranks;
89395fce210SBarry Smith   const PetscInt    *rootoffset,*leafoffset,*rootloc,*leafloc;
89495fce210SBarry Smith   const PetscMPIInt *rootranks,*leafranks;
89595fce210SBarry Smith   MPI_Request       *rootreqs,*leafreqs;
89695fce210SBarry Smith   size_t            unitbytes;
89795fce210SBarry Smith 
89895fce210SBarry Smith   PetscFunctionBegin;
89995fce210SBarry Smith   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,&rootranks,&rootoffset,&rootloc);CHKERRQ(ierr);
90095fce210SBarry Smith   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,&leafranks,&leafoffset,&leafloc);CHKERRQ(ierr);
90195fce210SBarry Smith   ierr = PetscSFBasicGetPack(sf,unit,rootdata,&link);CHKERRQ(ierr);
90295fce210SBarry Smith 
90395fce210SBarry Smith   unitbytes = link->unitbytes;
90495fce210SBarry Smith 
90595fce210SBarry Smith   ierr = PetscSFBasicPackGetReqs(sf,link,&rootreqs,&leafreqs);CHKERRQ(ierr);
90695fce210SBarry Smith   /* Eagerly post root receives */
90795fce210SBarry Smith   for (i=0; i<nrootranks; i++) {
90895fce210SBarry Smith     PetscMPIInt n = rootoffset[i+1] - rootoffset[i];
90995fce210SBarry Smith     ierr = MPI_Irecv(link->root+rootoffset[i]*unitbytes,n,unit,rootranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&rootreqs[i]);CHKERRQ(ierr);
91095fce210SBarry Smith   }
91195fce210SBarry Smith   /* Pack and send leaf data */
91295fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
91395fce210SBarry Smith     PetscMPIInt n          = leafoffset[i+1] - leafoffset[i];
91495fce210SBarry Smith     void        *packstart = link->leaf+leafoffset[i]*unitbytes;
91595fce210SBarry Smith     (*link->Pack)(n,leafloc+leafoffset[i],leafdata,packstart);
91695fce210SBarry Smith     ierr = MPI_Isend(packstart,n,unit,leafranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&leafreqs[i]);CHKERRQ(ierr);
91795fce210SBarry Smith   }
91895fce210SBarry Smith   PetscFunctionReturn(0);
91995fce210SBarry Smith }
92095fce210SBarry Smith 
92195fce210SBarry Smith #undef __FUNCT__
92295fce210SBarry Smith #define __FUNCT__ "PetscSFReduceEnd_Basic"
92395fce210SBarry Smith static PetscErrorCode PetscSFReduceEnd_Basic(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
92495fce210SBarry Smith {
92595fce210SBarry Smith   void             (*UnpackOp)(PetscInt,const PetscInt*,void*,const void*);
92695fce210SBarry Smith   PetscErrorCode   ierr;
92795fce210SBarry Smith   PetscSFBasicPack link;
92895fce210SBarry Smith   PetscInt         i,nrootranks;
92995fce210SBarry Smith   const PetscInt   *rootoffset,*rootloc;
93095fce210SBarry Smith 
93195fce210SBarry Smith   PetscFunctionBegin;
93295fce210SBarry Smith   ierr = PetscSFBasicGetPackInUse(sf,unit,rootdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr);
93395fce210SBarry Smith   /* This implementation could be changed to unpack as receives arrive, at the cost of non-determinism */
93495fce210SBarry Smith   ierr = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
93595fce210SBarry Smith   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,NULL,&rootoffset,&rootloc);CHKERRQ(ierr);
93695fce210SBarry Smith   ierr = PetscSFBasicPackGetUnpackOp(sf,link,op,&UnpackOp);CHKERRQ(ierr);
93795fce210SBarry Smith   for (i=0; i<nrootranks; i++) {
93895fce210SBarry Smith     PetscMPIInt n          = rootoffset[i+1] - rootoffset[i];
93995fce210SBarry Smith     const void  *packstart = link->root+rootoffset[i]*link->unitbytes;
94095fce210SBarry Smith 
94195fce210SBarry Smith     (*UnpackOp)(n,rootloc+rootoffset[i],rootdata,packstart);
94295fce210SBarry Smith   }
94395fce210SBarry Smith   ierr = PetscSFBasicReclaimPack(sf,&link);CHKERRQ(ierr);
94495fce210SBarry Smith   PetscFunctionReturn(0);
94595fce210SBarry Smith }
94695fce210SBarry Smith 
94795fce210SBarry Smith #undef __FUNCT__
94895fce210SBarry Smith #define __FUNCT__ "PetscSFFetchAndOpBegin_Basic"
94995fce210SBarry Smith static PetscErrorCode PetscSFFetchAndOpBegin_Basic(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
95095fce210SBarry Smith {
95195fce210SBarry Smith   PetscErrorCode ierr;
95295fce210SBarry Smith 
95395fce210SBarry Smith   PetscFunctionBegin;
95495fce210SBarry Smith   ierr = PetscSFReduceBegin_Basic(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
95595fce210SBarry Smith   PetscFunctionReturn(0);
95695fce210SBarry Smith }
95795fce210SBarry Smith 
95895fce210SBarry Smith #undef __FUNCT__
95995fce210SBarry Smith #define __FUNCT__ "PetscSFFetchAndOpEnd_Basic"
96095fce210SBarry Smith static PetscErrorCode PetscSFFetchAndOpEnd_Basic(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
96195fce210SBarry Smith {
96295fce210SBarry Smith   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;
96395fce210SBarry Smith   void              (*FetchAndOp)(PetscInt,const PetscInt*,void*,void*);
96495fce210SBarry Smith   PetscErrorCode    ierr;
96595fce210SBarry Smith   PetscSFBasicPack  link;
96695fce210SBarry Smith   PetscInt          i,nrootranks,nleafranks;
96795fce210SBarry Smith   const PetscInt    *rootoffset,*leafoffset,*rootloc,*leafloc;
96895fce210SBarry Smith   const PetscMPIInt *rootranks,*leafranks;
96995fce210SBarry Smith   MPI_Request       *rootreqs,*leafreqs;
97095fce210SBarry Smith   size_t            unitbytes;
97195fce210SBarry Smith 
97295fce210SBarry Smith   PetscFunctionBegin;
97395fce210SBarry Smith   ierr = PetscSFBasicGetPackInUse(sf,unit,rootdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr);
97495fce210SBarry Smith   /* This implementation could be changed to unpack as receives arrive, at the cost of non-determinism */
97595fce210SBarry Smith   ierr      = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
97695fce210SBarry Smith   unitbytes = link->unitbytes;
97795fce210SBarry Smith   ierr      = PetscSFBasicGetRootInfo(sf,&nrootranks,&rootranks,&rootoffset,&rootloc);CHKERRQ(ierr);
97895fce210SBarry Smith   ierr      = PetscSFBasicGetLeafInfo(sf,&nleafranks,&leafranks,&leafoffset,&leafloc);CHKERRQ(ierr);
97995fce210SBarry Smith   ierr      = PetscSFBasicPackGetReqs(sf,link,&rootreqs,&leafreqs);CHKERRQ(ierr);
98095fce210SBarry Smith   /* Post leaf receives */
98195fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
98295fce210SBarry Smith     PetscMPIInt n = leafoffset[i+1] - leafoffset[i];
98395fce210SBarry Smith     ierr = MPI_Irecv(link->leaf+leafoffset[i]*unitbytes,n,unit,leafranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&leafreqs[i]);CHKERRQ(ierr);
98495fce210SBarry Smith   }
98595fce210SBarry Smith   /* Process local fetch-and-op, post root sends */
98695fce210SBarry Smith   ierr = PetscSFBasicPackGetFetchAndOp(sf,link,op,&FetchAndOp);CHKERRQ(ierr);
98795fce210SBarry Smith   for (i=0; i<nrootranks; i++) {
98895fce210SBarry Smith     PetscMPIInt n          = rootoffset[i+1] - rootoffset[i];
98995fce210SBarry Smith     void        *packstart = link->root+rootoffset[i]*unitbytes;
99095fce210SBarry Smith 
99195fce210SBarry Smith     (*FetchAndOp)(n,rootloc+rootoffset[i],rootdata,packstart);
99295fce210SBarry Smith     ierr = MPI_Isend(packstart,n,unit,rootranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&rootreqs[i]);CHKERRQ(ierr);
99395fce210SBarry Smith   }
99495fce210SBarry Smith   ierr = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
99595fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
99695fce210SBarry Smith     PetscMPIInt n          = leafoffset[i+1] - leafoffset[i];
99795fce210SBarry Smith     const void  *packstart = link->leaf+leafoffset[i]*unitbytes;
99895fce210SBarry Smith     (*link->UnpackInsert)(n,leafloc+leafoffset[i],leafupdate,packstart);
99995fce210SBarry Smith   }
100095fce210SBarry Smith   ierr = PetscSFBasicReclaimPack(sf,&link);CHKERRQ(ierr);
100195fce210SBarry Smith   PetscFunctionReturn(0);
100295fce210SBarry Smith }
100395fce210SBarry Smith 
100495fce210SBarry Smith #undef __FUNCT__
100595fce210SBarry Smith #define __FUNCT__ "PetscSFCreate_Basic"
10068cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscSFCreate_Basic(PetscSF sf)
100795fce210SBarry Smith {
100895fce210SBarry Smith   PetscSF_Basic  *bas = (PetscSF_Basic*)sf->data;
100995fce210SBarry Smith   PetscErrorCode ierr;
101095fce210SBarry Smith 
101195fce210SBarry Smith   PetscFunctionBegin;
101295fce210SBarry Smith   sf->ops->SetUp           = PetscSFSetUp_Basic;
101395fce210SBarry Smith   sf->ops->SetFromOptions  = PetscSFSetFromOptions_Basic;
101495fce210SBarry Smith   sf->ops->Reset           = PetscSFReset_Basic;
101595fce210SBarry Smith   sf->ops->Destroy         = PetscSFDestroy_Basic;
101695fce210SBarry Smith   sf->ops->View            = PetscSFView_Basic;
101795fce210SBarry Smith   sf->ops->BcastBegin      = PetscSFBcastBegin_Basic;
101895fce210SBarry Smith   sf->ops->BcastEnd        = PetscSFBcastEnd_Basic;
101995fce210SBarry Smith   sf->ops->ReduceBegin     = PetscSFReduceBegin_Basic;
102095fce210SBarry Smith   sf->ops->ReduceEnd       = PetscSFReduceEnd_Basic;
102195fce210SBarry Smith   sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Basic;
102295fce210SBarry Smith   sf->ops->FetchAndOpEnd   = PetscSFFetchAndOpEnd_Basic;
102395fce210SBarry Smith 
102495fce210SBarry Smith   ierr     = PetscNewLog(sf,PetscSF_Basic,&bas);CHKERRQ(ierr);
102595fce210SBarry Smith   sf->data = (void*)bas;
102695fce210SBarry Smith   PetscFunctionReturn(0);
102795fce210SBarry Smith }
1028