xref: /petsc/src/vec/is/sf/impls/basic/sfbasic.c (revision 8cc058d9cd56c1ccb3be12a47760ddfc446aaffc)
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*);
1395fce210SBarry Smith   void (*FetchAndInsert)(PetscInt,const PetscInt*,void*,void*);
1495fce210SBarry Smith   void (*FetchAndAdd)(PetscInt,const PetscInt*,void*,void*);
1595fce210SBarry Smith   void (*FetchAndMin)(PetscInt,const PetscInt*,void*,void*);
1695fce210SBarry Smith   void (*FetchAndMax)(PetscInt,const PetscInt*,void*,void*);
1795fce210SBarry Smith   void (*FetchAndMinloc)(PetscInt,const PetscInt*,void*,void*);
1895fce210SBarry Smith   void (*FetchAndMaxloc)(PetscInt,const PetscInt*,void*,void*);
1995fce210SBarry Smith 
2095fce210SBarry Smith   MPI_Datatype     unit;
2195fce210SBarry Smith   size_t           unitbytes;   /* Number of bytes in a unit */
2295fce210SBarry Smith   const void       *key;        /* Array used as key for operation */
2395fce210SBarry Smith   char             *root;       /* Packed root data, contiguous by leaf rank */
2495fce210SBarry Smith   char             *leaf;       /* Packed leaf data, contiguous by root rank */
2595fce210SBarry Smith   MPI_Request      *requests;   /* Array of root requests followed by leaf requests */
2695fce210SBarry Smith   PetscSFBasicPack next;
2795fce210SBarry Smith };
2895fce210SBarry Smith 
2995fce210SBarry Smith typedef struct {
3095fce210SBarry Smith   PetscMPIInt      tag;
3195fce210SBarry Smith   PetscInt         niranks;     /* Number of incoming ranks (ranks accessing my roots) */
3295fce210SBarry Smith   PetscMPIInt      *iranks;     /* Array of ranks that reference my roots */
3395fce210SBarry Smith   PetscInt         itotal;      /* Total number of graph edges referencing my roots */
3495fce210SBarry Smith   PetscInt         *ioffset;    /* Array of length niranks+1 holding offset in irootloc[] for each rank */
3595fce210SBarry Smith   PetscInt         *irootloc;   /* Incoming roots referenced by ranks starting at ioffset[rank] */
3695fce210SBarry Smith   PetscSFBasicPack avail;       /* One or more entries per MPI Datatype, lazily constructed */
3795fce210SBarry Smith   PetscSFBasicPack inuse;       /* Buffers being used for transactions that have not yet completed */
3895fce210SBarry Smith } PetscSF_Basic;
3995fce210SBarry Smith 
4095fce210SBarry Smith #if !defined(PETSC_HAVE_MPI_TYPE_DUP) /* Danger: type is not reference counted; subject to ABA problem */
4195fce210SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MPI_Type_dup(MPI_Datatype datatype,MPI_Datatype *newtype)
4295fce210SBarry Smith {
4395fce210SBarry Smith   *newtype = datatype;
4495fce210SBarry Smith   return 0;
4595fce210SBarry Smith }
4695fce210SBarry Smith #endif
4795fce210SBarry Smith 
4895fce210SBarry Smith /*
4995fce210SBarry Smith  * MPI_Reduce_local is not really useful because it can't handle sparse data and it vectorizes "in the wrong direction",
5095fce210SBarry Smith  * therefore we pack data types manually. This section defines packing routines for the standard data types.
5195fce210SBarry Smith  */
5295fce210SBarry Smith 
5395fce210SBarry Smith #define CPPJoin2_exp(a,b) a ## b
5495fce210SBarry Smith #define CPPJoin2(a,b) CPPJoin2_exp(a,b)
5595fce210SBarry Smith #define CPPJoin3_exp_(a,b,c) a ## b ## _ ## c
5695fce210SBarry Smith #define CPPJoin3_(a,b,c) CPPJoin3_exp_(a,b,c)
5795fce210SBarry Smith 
5895fce210SBarry Smith /* Basic types without addition */
5995fce210SBarry Smith #define DEF_PackNoInit(type)                                            \
6095fce210SBarry Smith   static void CPPJoin2(Pack_,type)(PetscInt n,const PetscInt *idx,const void *unpacked,void *packed) { \
6195fce210SBarry Smith     const type *u = (const type*)unpacked;                              \
6295fce210SBarry Smith     type *p = (type*)packed;                                            \
6395fce210SBarry Smith     PetscInt i;                                                         \
6495fce210SBarry Smith     for (i=0; i<n; i++) p[i] = u[idx[i]];                               \
6595fce210SBarry Smith   }                                                                     \
6695fce210SBarry Smith   static void CPPJoin2(UnpackInsert_,type)(PetscInt n,const PetscInt *idx,void *unpacked,const void *packed) { \
6795fce210SBarry Smith     type *u = (type*)unpacked;                                          \
6895fce210SBarry Smith     const type *p = (const type*)packed;                                \
6995fce210SBarry Smith     PetscInt i;                                                         \
7095fce210SBarry Smith     for (i=0; i<n; i++) u[idx[i]] = p[i];                               \
7195fce210SBarry Smith   }                                                                     \
7295fce210SBarry Smith   static void CPPJoin2(FetchAndInsert_,type)(PetscInt n,const PetscInt *idx,void *unpacked,void *packed) { \
7395fce210SBarry Smith     type *u = (type*)unpacked;                                          \
7495fce210SBarry Smith     type *p = (type*)packed;                                            \
7595fce210SBarry Smith     PetscInt i;                                                         \
7695fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
7795fce210SBarry Smith       PetscInt j = idx[i];                                              \
7895fce210SBarry Smith       type t = u[j];                                                    \
7995fce210SBarry Smith       u[j] = p[i];                                                      \
8095fce210SBarry Smith       p[i] = t;                                                         \
8195fce210SBarry Smith     }                                                                   \
8295fce210SBarry Smith   }
8395fce210SBarry Smith 
8495fce210SBarry Smith /* Basic types defining addition */
8595fce210SBarry Smith #define DEF_PackAddNoInit(type)                                         \
8695fce210SBarry Smith   DEF_PackNoInit(type)                                                  \
8795fce210SBarry Smith   static void CPPJoin2(UnpackAdd_,type)(PetscInt n,const PetscInt *idx,void *unpacked,const void *packed) { \
8895fce210SBarry Smith     type *u = (type*)unpacked;                                          \
8995fce210SBarry Smith     const type *p = (const type*)packed;                                \
9095fce210SBarry Smith     PetscInt i;                                                         \
9195fce210SBarry Smith     for (i=0; i<n; i++) u[idx[i]] += p[i];                              \
9295fce210SBarry Smith   }                                                                     \
9395fce210SBarry Smith   static void CPPJoin2(FetchAndAdd_,type)(PetscInt n,const PetscInt *idx,void *unpacked,void *packed) { \
9495fce210SBarry Smith     type *u = (type*)unpacked;                                          \
9595fce210SBarry Smith     type *p = (type*)packed;                                            \
9695fce210SBarry Smith     PetscInt i;                                                         \
9795fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
9895fce210SBarry Smith       PetscInt j = idx[i];                                              \
9995fce210SBarry Smith       type t = u[j];                                                    \
10095fce210SBarry Smith       u[j] = t + p[i];                                                  \
10195fce210SBarry Smith       p[i] = t;                                                         \
10295fce210SBarry Smith     }                                                                   \
10395fce210SBarry Smith   }
10495fce210SBarry Smith #define DEF_Pack(type) \
10595fce210SBarry Smith   DEF_PackAddNoInit(type)                                               \
10695fce210SBarry Smith   static void CPPJoin2(PackInit_,type)(PetscSFBasicPack link) {         \
10795fce210SBarry Smith     link->Pack = CPPJoin2(Pack_,type);                                  \
10895fce210SBarry Smith     link->UnpackInsert = CPPJoin2(UnpackInsert_,type);                  \
10995fce210SBarry Smith     link->UnpackAdd = CPPJoin2(UnpackAdd_,type);                        \
11095fce210SBarry Smith     link->FetchAndInsert = CPPJoin2(FetchAndInsert_,type);              \
11195fce210SBarry Smith     link->FetchAndAdd = CPPJoin2(FetchAndAdd_,type);                    \
11295fce210SBarry Smith     link->unitbytes = sizeof(type);                                     \
11395fce210SBarry Smith   }
11495fce210SBarry Smith /* Comparable types */
11595fce210SBarry Smith #define DEF_PackCmp(type)                                               \
11695fce210SBarry Smith   DEF_PackAddNoInit(type)                                               \
11795fce210SBarry Smith   static void CPPJoin2(UnpackMax_,type)(PetscInt n,const PetscInt *idx,void *unpacked,const void *packed) { \
11895fce210SBarry Smith     type *u = (type*)unpacked;                                          \
11995fce210SBarry Smith     const type *p = (const type*)packed;                                \
12095fce210SBarry Smith     PetscInt i;                                                         \
12195fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
12295fce210SBarry Smith       type v = u[idx[i]];                                               \
12395fce210SBarry Smith       u[idx[i]] = PetscMax(v,p[i]);                                     \
12495fce210SBarry Smith     }                                                                   \
12595fce210SBarry Smith   }                                                                     \
12695fce210SBarry Smith   static void CPPJoin2(UnpackMin_,type)(PetscInt n,const PetscInt *idx,void *unpacked,const void *packed) { \
12795fce210SBarry Smith     type *u = (type*)unpacked;                                          \
12895fce210SBarry Smith     const type *p = (const type*)packed;                                \
12995fce210SBarry Smith     PetscInt i;                                                         \
13095fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
13195fce210SBarry Smith       type v = u[idx[i]];                                               \
13295fce210SBarry Smith       u[idx[i]] = PetscMin(v,p[i]);                                     \
13395fce210SBarry Smith     }                                                                   \
13495fce210SBarry Smith   }                                                                     \
13595fce210SBarry Smith   static void CPPJoin2(FetchAndMax_,type)(PetscInt n,const PetscInt *idx,void *unpacked,void *packed) { \
13695fce210SBarry Smith     type *u = (type*)unpacked;                                          \
13795fce210SBarry Smith     type *p = (type*)packed;                                            \
13895fce210SBarry Smith     PetscInt i;                                                         \
13995fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
14095fce210SBarry Smith       PetscInt j = idx[i];                                              \
14195fce210SBarry Smith       type v = u[j];                                                    \
14295fce210SBarry Smith       u[j] = PetscMax(v,p[i]);                                          \
14395fce210SBarry Smith       p[i] = v;                                                         \
14495fce210SBarry Smith     }                                                                   \
14595fce210SBarry Smith   }                                                                     \
14695fce210SBarry Smith   static void CPPJoin2(FetchAndMin_,type)(PetscInt n,const PetscInt *idx,void *unpacked,void *packed) { \
14795fce210SBarry Smith     type *u = (type*)unpacked;                                          \
14895fce210SBarry Smith     type *p = (type*)packed;                                            \
14995fce210SBarry Smith     PetscInt i;                                                         \
15095fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
15195fce210SBarry Smith       PetscInt j = idx[i];                                              \
15295fce210SBarry Smith       type v = u[j];                                                    \
15395fce210SBarry Smith       u[j] = PetscMin(v,p[i]);                                          \
15495fce210SBarry Smith       p[i] = v;                                                         \
15595fce210SBarry Smith     }                                                                   \
15695fce210SBarry Smith   }                                                                     \
15795fce210SBarry Smith   static void CPPJoin2(PackInit_,type)(PetscSFBasicPack link) {         \
15895fce210SBarry Smith     link->Pack = CPPJoin2(Pack_,type);                                  \
15995fce210SBarry Smith     link->UnpackInsert = CPPJoin2(UnpackInsert_,type);                  \
16095fce210SBarry Smith     link->UnpackAdd = CPPJoin2(UnpackAdd_,type);                        \
16195fce210SBarry Smith     link->UnpackMax = CPPJoin2(UnpackMax_,type);                        \
16295fce210SBarry Smith     link->UnpackMin = CPPJoin2(UnpackMin_,type);                        \
16395fce210SBarry Smith     link->FetchAndInsert = CPPJoin2(FetchAndInsert_,type);              \
16495fce210SBarry Smith     link->FetchAndAdd = CPPJoin2(FetchAndAdd_ ,type);                   \
16595fce210SBarry Smith     link->FetchAndMax = CPPJoin2(FetchAndMax_ ,type);                   \
16695fce210SBarry Smith     link->FetchAndMin = CPPJoin2(FetchAndMin_ ,type);                   \
16795fce210SBarry Smith     link->unitbytes = sizeof(type);                                     \
16895fce210SBarry Smith   }
16995fce210SBarry Smith 
17095fce210SBarry Smith /* Pair types */
17195fce210SBarry Smith #define CPPJoinloc_exp(base,op,t1,t2) base ## op ## loc_ ## t1 ## _ ## t2
17295fce210SBarry Smith #define CPPJoinloc(base,op,t1,t2) CPPJoinloc_exp(base,op,t1,t2)
17395fce210SBarry Smith #define PairType(type1,type2) CPPJoin3_(_pairtype_,type1,type2)
17495fce210SBarry Smith #define DEF_UnpackXloc(type1,type2,locname,op)                              \
17595fce210SBarry Smith   static void CPPJoinloc(Unpack,locname,type1,type2)(PetscInt n,const PetscInt *idx,void *unpacked,const void *packed) { \
17695fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;        \
17795fce210SBarry Smith     const PairType(type1,type2) *p = (const PairType(type1,type2)*)packed; \
17895fce210SBarry Smith     PetscInt i;                                                         \
17995fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
18095fce210SBarry Smith       PetscInt j = idx[i];                                              \
18195fce210SBarry Smith       if (p[i].a op u[j].a) {                                           \
18295fce210SBarry Smith         u[j].a = p[i].a;                                                \
18395fce210SBarry Smith         u[j].b = p[i].b;                                                \
18495fce210SBarry Smith       } else if (u[j].a == p[i].a) {                                    \
18595fce210SBarry Smith         u[j].b = PetscMin(u[j].b,p[i].b);                               \
18695fce210SBarry Smith       }                                                                 \
18795fce210SBarry Smith     }                                                                   \
18895fce210SBarry Smith   }                                                                     \
18995fce210SBarry Smith   static void CPPJoinloc(FetchAnd,locname,type1,type2)(PetscInt n,const PetscInt *idx,void *unpacked,void *packed) { \
19095fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;        \
19195fce210SBarry Smith     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;          \
19295fce210SBarry Smith     PetscInt i;                                                         \
19395fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
19495fce210SBarry Smith       PetscInt j = idx[i];                                              \
19595fce210SBarry Smith       PairType(type1,type2) v;                                          \
19695fce210SBarry Smith       v.a = u[j].a;                                                     \
19795fce210SBarry Smith       v.b = u[j].b;                                                     \
19895fce210SBarry Smith       if (p[i].a op u[j].a) {                                           \
19995fce210SBarry Smith         u[j].a = p[i].a;                                                \
20095fce210SBarry Smith         u[j].b = p[i].b;                                                \
20195fce210SBarry Smith       } else if (u[j].a == p[i].a) {                                    \
20295fce210SBarry Smith         u[j].b = PetscMin(u[j].b,p[i].b);                               \
20395fce210SBarry Smith       }                                                                 \
20495fce210SBarry Smith       p[i].a = v.a;                                                     \
20595fce210SBarry Smith       p[i].b = v.b;                                                     \
20695fce210SBarry Smith     }                                                                   \
20795fce210SBarry Smith   }
20895fce210SBarry Smith #define DEF_PackPair(type1,type2)                                       \
20995fce210SBarry Smith   typedef struct {type1 a; type2 b;} PairType(type1,type2);             \
21095fce210SBarry Smith   static void CPPJoin3_(Pack_,type1,type2)(PetscInt n,const PetscInt *idx,const void *unpacked,void *packed) { \
21195fce210SBarry Smith     const PairType(type1,type2) *u = (const PairType(type1,type2)*)unpacked; \
21295fce210SBarry Smith     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;          \
21395fce210SBarry Smith     PetscInt i;                                                         \
21495fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
21595fce210SBarry Smith       p[i].a = u[idx[i]].a;                                             \
21695fce210SBarry Smith       p[i].b = u[idx[i]].b;                                             \
21795fce210SBarry Smith     }                                                                   \
21895fce210SBarry Smith   }                                                                     \
21995fce210SBarry Smith   static void CPPJoin3_(UnpackInsert_,type1,type2)(PetscInt n,const PetscInt *idx,void *unpacked,const void *packed) { \
22095fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;       \
22195fce210SBarry Smith     const PairType(type1,type2) *p = (const PairType(type1,type2)*)packed; \
22295fce210SBarry Smith     PetscInt i;                                                         \
22395fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
22495fce210SBarry Smith       u[idx[i]].a = p[i].a;                                             \
22595fce210SBarry Smith       u[idx[i]].b = p[i].b;                                             \
22695fce210SBarry Smith     }                                                                   \
22795fce210SBarry Smith   }                                                                     \
22895fce210SBarry Smith   static void CPPJoin3_(UnpackAdd_,type1,type2)(PetscInt n,const PetscInt *idx,void *unpacked,const void *packed) { \
22995fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;       \
23095fce210SBarry Smith     const PairType(type1,type2) *p = (const PairType(type1,type2)*)packed; \
23195fce210SBarry Smith     PetscInt i;                                                         \
23295fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
23395fce210SBarry Smith       u[idx[i]].a += p[i].a;                                            \
23495fce210SBarry Smith       u[idx[i]].b += p[i].b;                                            \
23595fce210SBarry Smith     }                                                                   \
23695fce210SBarry Smith   }                                                                     \
23795fce210SBarry Smith   static void CPPJoin3_(FetchAndInsert_,type1,type2)(PetscInt n,const PetscInt *idx,void *unpacked,void *packed) { \
23895fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;        \
23995fce210SBarry Smith     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;          \
24095fce210SBarry Smith     PetscInt i;                                                         \
24195fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
24295fce210SBarry Smith       PetscInt j = idx[i];                                              \
24395fce210SBarry Smith       PairType(type1,type2) v;                                          \
24495fce210SBarry Smith       v.a = u[j].a;                                                     \
24595fce210SBarry Smith       v.b = u[j].b;                                                     \
24695fce210SBarry Smith       u[j].a = p[i].a;                                                  \
24795fce210SBarry Smith       u[j].b = p[i].b;                                                  \
24895fce210SBarry Smith       p[i].a = v.a;                                                     \
24995fce210SBarry Smith       p[i].b = v.b;                                                     \
25095fce210SBarry Smith     }                                                                   \
25195fce210SBarry Smith   }                                                                     \
25295fce210SBarry Smith   static void FetchAndAdd_ ## type1 ## _ ## type2(PetscInt n,const PetscInt *idx,void *unpacked,void *packed) { \
25395fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;       \
25495fce210SBarry Smith     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;         \
25595fce210SBarry Smith     PetscInt i;                                                         \
25695fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
25795fce210SBarry Smith       PetscInt j = idx[i];                                              \
25895fce210SBarry Smith       PairType(type1,type2) v;                                          \
25995fce210SBarry Smith       v.a = u[j].a;                                                     \
26095fce210SBarry Smith       v.b = u[j].b;                                                     \
26195fce210SBarry Smith       u[j].a = v.a + p[i].a;                                            \
26295fce210SBarry Smith       u[j].b = v.b + p[i].b;                                            \
26395fce210SBarry Smith       p[i].a = v.a;                                                     \
26495fce210SBarry Smith       p[i].b = v.b;                                                     \
26595fce210SBarry Smith     }                                                                   \
26695fce210SBarry Smith   }                                                                     \
26795fce210SBarry Smith   DEF_UnpackXloc(type1,type2,Max,>)                                     \
26895fce210SBarry Smith   DEF_UnpackXloc(type1,type2,Min,<)                                     \
26995fce210SBarry Smith   static void CPPJoin3_(PackInit_,type1,type2)(PetscSFBasicPack link) { \
27095fce210SBarry Smith     link->Pack = CPPJoin3_(Pack_,type1,type2);                          \
27195fce210SBarry Smith     link->UnpackInsert = CPPJoin3_(UnpackInsert_,type1,type2);          \
27295fce210SBarry Smith     link->UnpackAdd = CPPJoin3_(UnpackAdd_,type1,type2);                \
27395fce210SBarry Smith     link->UnpackMaxloc = CPPJoin3_(UnpackMaxloc_,type1,type2);          \
27495fce210SBarry Smith     link->UnpackMinloc = CPPJoin3_(UnpackMinloc_,type1,type2);          \
27595fce210SBarry Smith     link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,type1,type2);      \
27695fce210SBarry Smith     link->FetchAndAdd = CPPJoin3_(FetchAndAdd_,type1,type2);            \
27795fce210SBarry Smith     link->FetchAndMaxloc = CPPJoin3_(FetchAndMaxloc_,type1,type2);      \
27895fce210SBarry Smith     link->FetchAndMinloc = CPPJoin3_(FetchAndMinloc_,type1,type2);      \
27995fce210SBarry Smith     link->unitbytes = sizeof(PairType(type1,type2));                    \
28095fce210SBarry Smith   }
28195fce210SBarry Smith 
28295fce210SBarry Smith /* Currently only dumb blocks of data */
28395fce210SBarry Smith #define BlockType(unit,count) CPPJoin3_(_blocktype_,unit,count)
28495fce210SBarry Smith #define DEF_Block(unit,count)                                           \
28595fce210SBarry Smith   typedef struct {unit v[count];} BlockType(unit,count);                \
28695fce210SBarry Smith   DEF_PackNoInit(BlockType(unit,count))                                 \
28795fce210SBarry Smith   static void CPPJoin3_(PackInit_block_,unit,count)(PetscSFBasicPack link) { \
28895fce210SBarry Smith     link->Pack = CPPJoin2(Pack_,BlockType(unit,count));                 \
28995fce210SBarry Smith     link->UnpackInsert = CPPJoin2(UnpackInsert_,BlockType(unit,count)); \
29095fce210SBarry Smith     link->FetchAndInsert = CPPJoin2(FetchAndInsert_,BlockType(unit,count)); \
29195fce210SBarry Smith     link->unitbytes = sizeof(BlockType(unit,count));                    \
29295fce210SBarry Smith   }
29395fce210SBarry Smith 
29495fce210SBarry Smith DEF_PackCmp(int)
29595fce210SBarry Smith DEF_PackCmp(PetscInt)
29695fce210SBarry Smith DEF_PackCmp(PetscReal)
29795fce210SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
29895fce210SBarry Smith DEF_Pack(PetscComplex)
29995fce210SBarry Smith #endif
30095fce210SBarry Smith DEF_PackPair(int,int)
30195fce210SBarry Smith DEF_PackPair(PetscInt,PetscInt)
30295fce210SBarry Smith DEF_Block(int,1)
30395fce210SBarry Smith DEF_Block(int,2)
30495fce210SBarry Smith DEF_Block(int,3)
30595fce210SBarry Smith DEF_Block(int,4)
30695fce210SBarry Smith DEF_Block(int,5)
30795fce210SBarry Smith DEF_Block(int,6)
30895fce210SBarry Smith DEF_Block(int,7)
30995fce210SBarry Smith DEF_Block(int,8)
31095fce210SBarry Smith 
31195fce210SBarry Smith #undef __FUNCT__
31295fce210SBarry Smith #define __FUNCT__ "PetscSFSetUp_Basic"
31395fce210SBarry Smith static PetscErrorCode PetscSFSetUp_Basic(PetscSF sf)
31495fce210SBarry Smith {
31595fce210SBarry Smith   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
31695fce210SBarry Smith   PetscErrorCode ierr;
31795fce210SBarry Smith   PetscInt *rlengths,*ilengths,i;
31895fce210SBarry Smith   MPI_Comm comm;
31995fce210SBarry Smith   MPI_Request *rootreqs,*leafreqs;
32095fce210SBarry Smith 
32195fce210SBarry Smith   PetscFunctionBegin;
32295fce210SBarry Smith   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
32395fce210SBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)sf,&bas->tag);CHKERRQ(ierr);
32495fce210SBarry Smith   /*
32595fce210SBarry Smith    * Inform roots about how many leaves and from which ranks
32695fce210SBarry Smith    */
32795fce210SBarry Smith   ierr = PetscMalloc(sf->nranks*sizeof(PetscInt),&rlengths);CHKERRQ(ierr);
32895fce210SBarry Smith   /* Determine number, sending ranks, and length of incoming  */
32995fce210SBarry Smith   for (i=0; i<sf->nranks; i++) {
33095fce210SBarry Smith     rlengths[i] = sf->roffset[i+1] - sf->roffset[i]; /* Number of roots referenced by my leaves; for rank sf->ranks[i] */
33195fce210SBarry Smith   }
33295fce210SBarry Smith   ierr = PetscCommBuildTwoSided(comm,1,MPIU_INT,sf->nranks,sf->ranks,rlengths,&bas->niranks,&bas->iranks,(void**)&ilengths);CHKERRQ(ierr);
33395fce210SBarry Smith   ierr = PetscFree(rlengths);CHKERRQ(ierr);
33495fce210SBarry Smith 
33595fce210SBarry Smith   /* Send leaf identities to roots */
33695fce210SBarry Smith   for (i=0,bas->itotal=0; i<bas->niranks; i++) bas->itotal += ilengths[i];
33795fce210SBarry Smith   ierr = PetscMalloc2(bas->niranks+1,PetscInt,&bas->ioffset,bas->itotal,PetscInt,&bas->irootloc);CHKERRQ(ierr);
33895fce210SBarry Smith   ierr = PetscMalloc((bas->niranks+sf->nranks)*sizeof(MPI_Request),&rootreqs);CHKERRQ(ierr);
33995fce210SBarry Smith 
34095fce210SBarry Smith   leafreqs = rootreqs + bas->niranks;
34195fce210SBarry Smith   bas->ioffset[0] = 0;
34295fce210SBarry Smith   for (i=0; i<bas->niranks; i++) {
34395fce210SBarry Smith     bas->ioffset[i+1] = bas->ioffset[i] + ilengths[i];
34495fce210SBarry Smith     ierr = MPI_Irecv(bas->irootloc+bas->ioffset[i],ilengths[i],MPIU_INT,bas->iranks[i],bas->tag,comm,&rootreqs[i]);CHKERRQ(ierr);
34595fce210SBarry Smith   }
34695fce210SBarry Smith   for (i=0; i<sf->nranks; i++) {
34795fce210SBarry Smith     PetscMPIInt npoints;
34895fce210SBarry Smith     ierr = PetscMPIIntCast(sf->roffset[i+1] - sf->roffset[i],&npoints);CHKERRQ(ierr);
34995fce210SBarry Smith     ierr = MPI_Isend(sf->rremote+sf->roffset[i],npoints,MPIU_INT,sf->ranks[i],bas->tag,comm,&leafreqs[i]);CHKERRQ(ierr);
35095fce210SBarry Smith   }
35195fce210SBarry Smith   ierr = MPI_Waitall(sf->nranks+bas->niranks,rootreqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
35295fce210SBarry Smith   ierr = PetscFree(ilengths);CHKERRQ(ierr);
35395fce210SBarry Smith   ierr = PetscFree(rootreqs);CHKERRQ(ierr);
35495fce210SBarry Smith   PetscFunctionReturn(0);
35595fce210SBarry Smith }
35695fce210SBarry Smith 
35795fce210SBarry Smith #undef __FUNCT__
35895fce210SBarry Smith #define __FUNCT__ "PetscSFBasicPackTypeSetup"
35995fce210SBarry Smith static PetscErrorCode PetscSFBasicPackTypeSetup(PetscSFBasicPack link,MPI_Datatype unit)
36095fce210SBarry Smith {
36195fce210SBarry Smith   PetscErrorCode ierr;
36295fce210SBarry Smith   PetscBool      isInt,isPetscInt,isPetscReal,is2Int,is2PetscInt;
36395fce210SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
36495fce210SBarry Smith   PetscBool isPetscComplex;
36595fce210SBarry Smith #endif
36695fce210SBarry Smith 
36795fce210SBarry Smith   PetscFunctionBegin;
36895fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPI_INT,&isInt);CHKERRQ(ierr);
36995fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPIU_INT,&isPetscInt);CHKERRQ(ierr);
37095fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPIU_REAL,&isPetscReal);CHKERRQ(ierr);
37195fce210SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
37295fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPIU_COMPLEX,&isPetscComplex);CHKERRQ(ierr);
37395fce210SBarry Smith #endif
37495fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPI_2INT,&is2Int);CHKERRQ(ierr);
37595fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPIU_2INT,&is2PetscInt);CHKERRQ(ierr);
37695fce210SBarry Smith   if (isInt) PackInit_int(link);
37795fce210SBarry Smith   else if (isPetscInt) PackInit_PetscInt(link);
37895fce210SBarry Smith   else if (isPetscReal) PackInit_PetscReal(link);
37995fce210SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
38095fce210SBarry Smith   else if (isPetscComplex) PackInit_PetscComplex(link);
38195fce210SBarry Smith #endif
38295fce210SBarry Smith   else if (is2Int) PackInit_int_int(link);
38395fce210SBarry Smith   else if (is2PetscInt) PackInit_PetscInt_PetscInt(link);
38495fce210SBarry Smith   else {
38595fce210SBarry Smith     PetscMPIInt bytes;
38695fce210SBarry Smith     ierr = MPI_Type_size(unit,&bytes);CHKERRQ(ierr);
38795fce210SBarry Smith     if (bytes % sizeof(int)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for type size not divisible by %D",sizeof(int));
38895fce210SBarry Smith     switch (bytes / sizeof(int)) {
38995fce210SBarry Smith     case 1: PackInit_block_int_1(link); break;
39095fce210SBarry Smith     case 2: PackInit_block_int_2(link); break;
39195fce210SBarry Smith     case 3: PackInit_block_int_3(link); break;
39295fce210SBarry Smith     case 4: PackInit_block_int_4(link); break;
39395fce210SBarry Smith     case 5: PackInit_block_int_5(link); break;
39495fce210SBarry Smith     case 6: PackInit_block_int_6(link); break;
39595fce210SBarry Smith     case 7: PackInit_block_int_7(link); break;
39695fce210SBarry Smith     case 8: PackInit_block_int_8(link); break;
39795fce210SBarry Smith     default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for arbitrary block sizes");
39895fce210SBarry Smith     }
39995fce210SBarry Smith   }
40095fce210SBarry Smith   ierr = MPI_Type_dup(unit,&link->unit);CHKERRQ(ierr);
40195fce210SBarry Smith   PetscFunctionReturn(0);
40295fce210SBarry Smith }
40395fce210SBarry Smith 
40495fce210SBarry Smith #undef __FUNCT__
40595fce210SBarry Smith #define __FUNCT__ "PetscSFBasicPackGetUnpackOp"
40695fce210SBarry Smith static PetscErrorCode PetscSFBasicPackGetUnpackOp(PetscSF sf,PetscSFBasicPack link,MPI_Op op,void (**UnpackOp)(PetscInt,const PetscInt*,void*,const void*))
40795fce210SBarry Smith {
40895fce210SBarry Smith   PetscFunctionBegin;
40995fce210SBarry Smith   *UnpackOp = NULL;
41095fce210SBarry Smith   if (op == MPI_REPLACE) *UnpackOp = link->UnpackInsert;
41195fce210SBarry Smith   else if (op == MPI_SUM || op == MPIU_SUM) *UnpackOp = link->UnpackAdd;
41295fce210SBarry Smith   else if (op == MPI_MAX || op == MPIU_MAX) *UnpackOp = link->UnpackMax;
41395fce210SBarry Smith   else if (op == MPI_MIN || op == MPIU_MIN) *UnpackOp = link->UnpackMin;
41495fce210SBarry Smith   else if (op == MPI_MAXLOC) *UnpackOp = link->UnpackMaxloc;
41595fce210SBarry Smith   else if (op == MPI_MINLOC) *UnpackOp = link->UnpackMinloc;
41695fce210SBarry Smith   else SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_SUP,"No support for MPI_Op");
41795fce210SBarry Smith   PetscFunctionReturn(0);
41895fce210SBarry Smith }
41995fce210SBarry Smith #undef __FUNCT__
42095fce210SBarry Smith #define __FUNCT__ "PetscSFBasicPackGetFetchAndOp"
42195fce210SBarry Smith static PetscErrorCode PetscSFBasicPackGetFetchAndOp(PetscSF sf,PetscSFBasicPack link,MPI_Op op,void (**FetchAndOp)(PetscInt,const PetscInt*,void*,void*))
42295fce210SBarry Smith {
42395fce210SBarry Smith   PetscFunctionBegin;
42495fce210SBarry Smith   *FetchAndOp = NULL;
42595fce210SBarry Smith   if (op == MPI_REPLACE) *FetchAndOp = link->FetchAndInsert;
42695fce210SBarry Smith   else if (op == MPI_SUM || op == MPIU_SUM) *FetchAndOp = link->FetchAndAdd;
42795fce210SBarry Smith   else if (op == MPI_MAX || op == MPIU_MAX) *FetchAndOp = link->FetchAndMax;
42895fce210SBarry Smith   else if (op == MPI_MIN || op == MPIU_MIN) *FetchAndOp = link->FetchAndMin;
42995fce210SBarry Smith   else if (op == MPI_MAXLOC) *FetchAndOp = link->FetchAndMaxloc;
43095fce210SBarry Smith   else if (op == MPI_MINLOC) *FetchAndOp = link->FetchAndMinloc;
43195fce210SBarry Smith   else SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_SUP,"No support for MPI_Op");
43295fce210SBarry Smith   PetscFunctionReturn(0);
43395fce210SBarry Smith }
43495fce210SBarry Smith 
43595fce210SBarry Smith #undef __FUNCT__
43695fce210SBarry Smith #define __FUNCT__ "PetscSFBasicPackGetReqs"
43795fce210SBarry Smith static PetscErrorCode PetscSFBasicPackGetReqs(PetscSF sf,PetscSFBasicPack link,MPI_Request **rootreqs,MPI_Request **leafreqs)
43895fce210SBarry Smith {
43995fce210SBarry Smith   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
44095fce210SBarry Smith 
44195fce210SBarry Smith   PetscFunctionBegin;
44295fce210SBarry Smith   if (rootreqs) *rootreqs = link->requests;
44395fce210SBarry Smith   if (leafreqs) *leafreqs = link->requests + bas->niranks;
44495fce210SBarry Smith   PetscFunctionReturn(0);
44595fce210SBarry Smith }
44695fce210SBarry Smith 
44795fce210SBarry Smith #undef __FUNCT__
44895fce210SBarry Smith #define __FUNCT__ "PetscSFBasicPackWaitall"
44995fce210SBarry Smith static PetscErrorCode PetscSFBasicPackWaitall(PetscSF sf,PetscSFBasicPack link)
45095fce210SBarry Smith {
45195fce210SBarry Smith   PetscSF_Basic  *bas = (PetscSF_Basic*)sf->data;
45295fce210SBarry Smith   PetscErrorCode ierr;
45395fce210SBarry Smith 
45495fce210SBarry Smith   PetscFunctionBegin;
45595fce210SBarry Smith   ierr = MPI_Waitall(bas->niranks+sf->nranks,link->requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
45695fce210SBarry Smith   PetscFunctionReturn(0);
45795fce210SBarry Smith }
45895fce210SBarry Smith 
45995fce210SBarry Smith #undef __FUNCT__
46095fce210SBarry Smith #define __FUNCT__ "PetscSFBasicGetRootInfo"
46195fce210SBarry Smith static PetscErrorCode PetscSFBasicGetRootInfo(PetscSF sf,PetscInt *nrootranks,const PetscMPIInt **rootranks,const PetscInt **rootoffset,const PetscInt **rootloc)
46295fce210SBarry Smith {
46395fce210SBarry Smith   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
46495fce210SBarry Smith 
46595fce210SBarry Smith   PetscFunctionBegin;
46695fce210SBarry Smith   if (nrootranks) *nrootranks = bas->niranks;
46795fce210SBarry Smith   if (rootranks)  *rootranks  = bas->iranks;
46895fce210SBarry Smith   if (rootoffset) *rootoffset = bas->ioffset;
46995fce210SBarry Smith   if (rootloc)    *rootloc    = bas->irootloc;
47095fce210SBarry Smith   PetscFunctionReturn(0);
47195fce210SBarry Smith }
47295fce210SBarry Smith 
47395fce210SBarry Smith #undef __FUNCT__
47495fce210SBarry Smith #define __FUNCT__ "PetscSFBasicGetLeafInfo"
47595fce210SBarry Smith static PetscErrorCode PetscSFBasicGetLeafInfo(PetscSF sf,PetscInt *nleafranks,const PetscMPIInt **leafranks,const PetscInt **leafoffset,const PetscInt **leafloc)
47695fce210SBarry Smith {
47795fce210SBarry Smith   PetscFunctionBegin;
47895fce210SBarry Smith   if (nleafranks) *nleafranks = sf->nranks;
47995fce210SBarry Smith   if (leafranks)  *leafranks  = sf->ranks;
48095fce210SBarry Smith   if (leafoffset) *leafoffset = sf->roffset;
48195fce210SBarry Smith   if (leafloc)    *leafloc    = sf->rmine;
48295fce210SBarry Smith   PetscFunctionReturn(0);
48395fce210SBarry Smith }
48495fce210SBarry Smith 
48595fce210SBarry Smith #undef __FUNCT__
48695fce210SBarry Smith #define __FUNCT__ "PetscSFBasicGetPack"
48795fce210SBarry Smith static PetscErrorCode PetscSFBasicGetPack(PetscSF sf,MPI_Datatype unit,const void *key,PetscSFBasicPack *mylink)
48895fce210SBarry Smith {
48995fce210SBarry Smith   PetscSF_Basic    *bas = (PetscSF_Basic*)sf->data;
49095fce210SBarry Smith   PetscErrorCode   ierr;
49195fce210SBarry Smith   PetscSFBasicPack link,*p;
49295fce210SBarry Smith   PetscInt         nrootranks,nleafranks;
49395fce210SBarry Smith   const PetscInt   *rootoffset,*leafoffset;
49495fce210SBarry Smith 
49595fce210SBarry Smith   PetscFunctionBegin;
49695fce210SBarry Smith   /* Look for types in cache */
49795fce210SBarry Smith   for (p=&bas->avail; (link=*p); p=&link->next) {
49895fce210SBarry Smith     PetscBool match;
49995fce210SBarry Smith     ierr = MPIPetsc_Type_compare(unit,link->unit,&match);CHKERRQ(ierr);
50095fce210SBarry Smith     if (match) {
50195fce210SBarry Smith       *p = link->next;          /* Remove from available list */
50295fce210SBarry Smith       goto found;
50395fce210SBarry Smith     }
50495fce210SBarry Smith   }
50595fce210SBarry Smith 
50695fce210SBarry Smith   /* Create new composite types for each send rank */
50795fce210SBarry Smith   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,NULL,&rootoffset,NULL);CHKERRQ(ierr);
50895fce210SBarry Smith   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,NULL,&leafoffset,NULL);CHKERRQ(ierr);
50995fce210SBarry Smith   ierr = PetscNew(struct _n_PetscSFBasicPack,&link);CHKERRQ(ierr);
51095fce210SBarry Smith   ierr = PetscSFBasicPackTypeSetup(link,unit);CHKERRQ(ierr);
51195fce210SBarry Smith   ierr = PetscMalloc2(rootoffset[nrootranks]*link->unitbytes,char,&link->root,leafoffset[nleafranks]*link->unitbytes,char,&link->leaf);CHKERRQ(ierr);
51295fce210SBarry Smith   ierr = PetscMalloc((nrootranks+nleafranks)*sizeof(MPI_Request),&link->requests);CHKERRQ(ierr);
51395fce210SBarry Smith 
51495fce210SBarry Smith found:
51595fce210SBarry Smith   link->key  = key;
51695fce210SBarry Smith   link->next = bas->inuse;
51795fce210SBarry Smith   bas->inuse = link;
51895fce210SBarry Smith 
51995fce210SBarry Smith   *mylink = link;
52095fce210SBarry Smith   PetscFunctionReturn(0);
52195fce210SBarry Smith }
52295fce210SBarry Smith 
52395fce210SBarry Smith #undef __FUNCT__
52495fce210SBarry Smith #define __FUNCT__ "PetscSFBasicGetPackInUse"
52595fce210SBarry Smith static PetscErrorCode PetscSFBasicGetPackInUse(PetscSF sf,MPI_Datatype unit,const void *key,PetscCopyMode cmode,PetscSFBasicPack *mylink)
52695fce210SBarry Smith {
52795fce210SBarry Smith   PetscSF_Basic    *bas = (PetscSF_Basic*)sf->data;
52895fce210SBarry Smith   PetscErrorCode   ierr;
52995fce210SBarry Smith   PetscSFBasicPack link,*p;
53095fce210SBarry Smith 
53195fce210SBarry Smith   PetscFunctionBegin;
53295fce210SBarry Smith   /* Look for types in cache */
53395fce210SBarry Smith   for (p=&bas->inuse; (link=*p); p=&link->next) {
53495fce210SBarry Smith     PetscBool match;
53595fce210SBarry Smith     ierr = MPIPetsc_Type_compare(unit,link->unit,&match);CHKERRQ(ierr);
53695fce210SBarry Smith     if (match) {
53795fce210SBarry Smith       switch (cmode) {
53895fce210SBarry Smith       case PETSC_OWN_POINTER: *p = link->next; break; /* Remove from inuse list */
53995fce210SBarry Smith       case PETSC_USE_POINTER: break;
54095fce210SBarry Smith       default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"invalid cmode");
54195fce210SBarry Smith       }
54295fce210SBarry Smith       *mylink = link;
54395fce210SBarry Smith       PetscFunctionReturn(0);
54495fce210SBarry Smith     }
54595fce210SBarry Smith   }
54695fce210SBarry Smith   SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Could not find pack");
54795fce210SBarry Smith   PetscFunctionReturn(0);
54895fce210SBarry Smith }
54995fce210SBarry Smith 
55095fce210SBarry Smith #undef __FUNCT__
55195fce210SBarry Smith #define __FUNCT__ "PetscSFBasicReclaimPack"
55295fce210SBarry Smith static PetscErrorCode PetscSFBasicReclaimPack(PetscSF sf,PetscSFBasicPack *link)
55395fce210SBarry Smith {
55495fce210SBarry Smith   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
55595fce210SBarry Smith 
55695fce210SBarry Smith   PetscFunctionBegin;
55795fce210SBarry Smith   (*link)->key  = NULL;
55895fce210SBarry Smith   (*link)->next = bas->avail;
55995fce210SBarry Smith   bas->avail    = *link;
56095fce210SBarry Smith   *link         = NULL;
56195fce210SBarry Smith   PetscFunctionReturn(0);
56295fce210SBarry Smith }
56395fce210SBarry Smith 
56495fce210SBarry Smith #undef __FUNCT__
56595fce210SBarry Smith #define __FUNCT__ "PetscSFSetFromOptions_Basic"
56695fce210SBarry Smith static PetscErrorCode PetscSFSetFromOptions_Basic(PetscSF sf)
56795fce210SBarry Smith {
56895fce210SBarry Smith   PetscErrorCode ierr;
56995fce210SBarry Smith 
57095fce210SBarry Smith   PetscFunctionBegin;
57195fce210SBarry Smith   ierr = PetscOptionsHead("PetscSF Basic options");CHKERRQ(ierr);
57295fce210SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
57395fce210SBarry Smith   PetscFunctionReturn(0);
57495fce210SBarry Smith }
57595fce210SBarry Smith 
57695fce210SBarry Smith #undef __FUNCT__
57795fce210SBarry Smith #define __FUNCT__ "PetscSFReset_Basic"
57895fce210SBarry Smith static PetscErrorCode PetscSFReset_Basic(PetscSF sf)
57995fce210SBarry Smith {
58095fce210SBarry Smith   PetscSF_Basic    *bas = (PetscSF_Basic*)sf->data;
58195fce210SBarry Smith   PetscErrorCode   ierr;
58295fce210SBarry Smith   PetscSFBasicPack link,next;
58395fce210SBarry Smith 
58495fce210SBarry Smith   PetscFunctionBegin;
58595fce210SBarry Smith   ierr = PetscFree(bas->iranks);CHKERRQ(ierr);
58695fce210SBarry Smith   ierr = PetscFree2(bas->ioffset,bas->irootloc);CHKERRQ(ierr);
58795fce210SBarry Smith   if (bas->inuse) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Outstanding operation has not been completed");
58895fce210SBarry Smith   for (link=bas->avail; link; link=next) {
58995fce210SBarry Smith     next = link->next;
59095fce210SBarry Smith     ierr = MPI_Type_free(&link->unit);CHKERRQ(ierr);
59195fce210SBarry Smith     ierr = PetscFree2(link->root,link->leaf);CHKERRQ(ierr);
59295fce210SBarry Smith     ierr = PetscFree(link->requests);CHKERRQ(ierr);
59395fce210SBarry Smith     ierr = PetscFree(link);CHKERRQ(ierr);
59495fce210SBarry Smith   }
59595fce210SBarry Smith   bas->avail = NULL;
59695fce210SBarry Smith   PetscFunctionReturn(0);
59795fce210SBarry Smith }
59895fce210SBarry Smith 
59995fce210SBarry Smith #undef __FUNCT__
60095fce210SBarry Smith #define __FUNCT__ "PetscSFDestroy_Basic"
60195fce210SBarry Smith static PetscErrorCode PetscSFDestroy_Basic(PetscSF sf)
60295fce210SBarry Smith {
60395fce210SBarry Smith   PetscErrorCode ierr;
60495fce210SBarry Smith 
60595fce210SBarry Smith   PetscFunctionBegin;
60695fce210SBarry Smith   ierr = PetscSFReset_Basic(sf);CHKERRQ(ierr);
60795fce210SBarry Smith   ierr = PetscFree(sf->data);CHKERRQ(ierr);
60895fce210SBarry Smith   PetscFunctionReturn(0);
60995fce210SBarry Smith }
61095fce210SBarry Smith 
61195fce210SBarry Smith #undef __FUNCT__
61295fce210SBarry Smith #define __FUNCT__ "PetscSFView_Basic"
61395fce210SBarry Smith static PetscErrorCode PetscSFView_Basic(PetscSF sf,PetscViewer viewer)
61495fce210SBarry Smith {
61595fce210SBarry Smith   /* PetscSF_Basic *bas = (PetscSF_Basic*)sf->data; */
61695fce210SBarry Smith   PetscErrorCode ierr;
61795fce210SBarry Smith   PetscBool      iascii;
61895fce210SBarry Smith 
61995fce210SBarry Smith   PetscFunctionBegin;
62095fce210SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
62195fce210SBarry Smith   if (iascii) {
62295fce210SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  sort=%s\n",sf->rankorder ? "rank-order" : "unordered");CHKERRQ(ierr);
62395fce210SBarry Smith   }
62495fce210SBarry Smith   PetscFunctionReturn(0);
62595fce210SBarry Smith }
62695fce210SBarry Smith 
62795fce210SBarry Smith #undef __FUNCT__
62895fce210SBarry Smith #define __FUNCT__ "PetscSFBcastBegin_Basic"
62995fce210SBarry Smith /* Send from roots to leaves */
63095fce210SBarry Smith static PetscErrorCode PetscSFBcastBegin_Basic(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
63195fce210SBarry Smith {
63295fce210SBarry Smith   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;
63395fce210SBarry Smith   PetscErrorCode    ierr;
63495fce210SBarry Smith   PetscSFBasicPack  link;
63595fce210SBarry Smith   PetscInt          i,nrootranks,nleafranks;
63695fce210SBarry Smith   const PetscInt    *rootoffset,*leafoffset,*rootloc,*leafloc;
63795fce210SBarry Smith   const PetscMPIInt *rootranks,*leafranks;
63895fce210SBarry Smith   MPI_Request       *rootreqs,*leafreqs;
63995fce210SBarry Smith   size_t            unitbytes;
64095fce210SBarry Smith 
64195fce210SBarry Smith   PetscFunctionBegin;
64295fce210SBarry Smith   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,&rootranks,&rootoffset,&rootloc);CHKERRQ(ierr);
64395fce210SBarry Smith   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,&leafranks,&leafoffset,&leafloc);CHKERRQ(ierr);
64495fce210SBarry Smith   ierr = PetscSFBasicGetPack(sf,unit,rootdata,&link);CHKERRQ(ierr);
64595fce210SBarry Smith 
64695fce210SBarry Smith   unitbytes = link->unitbytes;
64795fce210SBarry Smith 
64895fce210SBarry Smith   ierr = PetscSFBasicPackGetReqs(sf,link,&rootreqs,&leafreqs);CHKERRQ(ierr);
64995fce210SBarry Smith   /* Eagerly post leaf receives */
65095fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
65195fce210SBarry Smith     PetscMPIInt n = leafoffset[i+1] - leafoffset[i];
65295fce210SBarry Smith     ierr = MPI_Irecv(link->leaf+leafoffset[i]*unitbytes,n,unit,leafranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&leafreqs[i]);CHKERRQ(ierr);
65395fce210SBarry Smith   }
65495fce210SBarry Smith   /* Pack and send root data */
65595fce210SBarry Smith   for (i=0; i<nrootranks; i++) {
65695fce210SBarry Smith     PetscMPIInt n          = rootoffset[i+1] - rootoffset[i];
65795fce210SBarry Smith     void        *packstart = link->root+rootoffset[i]*unitbytes;
65895fce210SBarry Smith     (*link->Pack)(n,rootloc+rootoffset[i],rootdata,packstart);
65995fce210SBarry Smith     ierr = MPI_Isend(packstart,n,unit,rootranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&rootreqs[i]);CHKERRQ(ierr);
66095fce210SBarry Smith   }
66195fce210SBarry Smith   PetscFunctionReturn(0);
66295fce210SBarry Smith }
66395fce210SBarry Smith 
66495fce210SBarry Smith #undef __FUNCT__
66595fce210SBarry Smith #define __FUNCT__ "PetscSFBcastEnd_Basic"
66695fce210SBarry Smith PetscErrorCode PetscSFBcastEnd_Basic(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
66795fce210SBarry Smith {
66895fce210SBarry Smith   PetscErrorCode   ierr;
66995fce210SBarry Smith   PetscSFBasicPack link;
67095fce210SBarry Smith   PetscInt         i,nleafranks;
67195fce210SBarry Smith   const PetscInt   *leafoffset,*leafloc;
67295fce210SBarry Smith 
67395fce210SBarry Smith   PetscFunctionBegin;
67495fce210SBarry Smith   ierr = PetscSFBasicGetPackInUse(sf,unit,rootdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr);
67595fce210SBarry Smith   ierr = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
67695fce210SBarry Smith   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,NULL,&leafoffset,&leafloc);CHKERRQ(ierr);
67795fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
67895fce210SBarry Smith     PetscMPIInt n          = leafoffset[i+1] - leafoffset[i];
67995fce210SBarry Smith     const void  *packstart = link->leaf+leafoffset[i]*link->unitbytes;
68095fce210SBarry Smith     (*link->UnpackInsert)(n,leafloc+leafoffset[i],leafdata,packstart);
68195fce210SBarry Smith   }
68295fce210SBarry Smith   ierr = PetscSFBasicReclaimPack(sf,&link);CHKERRQ(ierr);
68395fce210SBarry Smith   PetscFunctionReturn(0);
68495fce210SBarry Smith }
68595fce210SBarry Smith 
68695fce210SBarry Smith #undef __FUNCT__
68795fce210SBarry Smith #define __FUNCT__ "PetscSFReduceBegin_Basic"
68895fce210SBarry Smith /* leaf -> root with reduction */
68995fce210SBarry Smith PetscErrorCode PetscSFReduceBegin_Basic(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
69095fce210SBarry Smith {
69195fce210SBarry Smith   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;
69295fce210SBarry Smith   PetscSFBasicPack  link;
69395fce210SBarry Smith   PetscErrorCode    ierr;
69495fce210SBarry Smith   PetscInt          i,nrootranks,nleafranks;
69595fce210SBarry Smith   const PetscInt    *rootoffset,*leafoffset,*rootloc,*leafloc;
69695fce210SBarry Smith   const PetscMPIInt *rootranks,*leafranks;
69795fce210SBarry Smith   MPI_Request       *rootreqs,*leafreqs;
69895fce210SBarry Smith   size_t            unitbytes;
69995fce210SBarry Smith 
70095fce210SBarry Smith   PetscFunctionBegin;
70195fce210SBarry Smith   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,&rootranks,&rootoffset,&rootloc);CHKERRQ(ierr);
70295fce210SBarry Smith   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,&leafranks,&leafoffset,&leafloc);CHKERRQ(ierr);
70395fce210SBarry Smith   ierr = PetscSFBasicGetPack(sf,unit,rootdata,&link);CHKERRQ(ierr);
70495fce210SBarry Smith 
70595fce210SBarry Smith   unitbytes = link->unitbytes;
70695fce210SBarry Smith 
70795fce210SBarry Smith   ierr = PetscSFBasicPackGetReqs(sf,link,&rootreqs,&leafreqs);CHKERRQ(ierr);
70895fce210SBarry Smith   /* Eagerly post root receives */
70995fce210SBarry Smith   for (i=0; i<nrootranks; i++) {
71095fce210SBarry Smith     PetscMPIInt n = rootoffset[i+1] - rootoffset[i];
71195fce210SBarry Smith     ierr = MPI_Irecv(link->root+rootoffset[i]*unitbytes,n,unit,rootranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&rootreqs[i]);CHKERRQ(ierr);
71295fce210SBarry Smith   }
71395fce210SBarry Smith   /* Pack and send leaf data */
71495fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
71595fce210SBarry Smith     PetscMPIInt n          = leafoffset[i+1] - leafoffset[i];
71695fce210SBarry Smith     void        *packstart = link->leaf+leafoffset[i]*unitbytes;
71795fce210SBarry Smith     (*link->Pack)(n,leafloc+leafoffset[i],leafdata,packstart);
71895fce210SBarry Smith     ierr = MPI_Isend(packstart,n,unit,leafranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&leafreqs[i]);CHKERRQ(ierr);
71995fce210SBarry Smith   }
72095fce210SBarry Smith   PetscFunctionReturn(0);
72195fce210SBarry Smith }
72295fce210SBarry Smith 
72395fce210SBarry Smith #undef __FUNCT__
72495fce210SBarry Smith #define __FUNCT__ "PetscSFReduceEnd_Basic"
72595fce210SBarry Smith static PetscErrorCode PetscSFReduceEnd_Basic(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
72695fce210SBarry Smith {
72795fce210SBarry Smith   void             (*UnpackOp)(PetscInt,const PetscInt*,void*,const void*);
72895fce210SBarry Smith   PetscErrorCode   ierr;
72995fce210SBarry Smith   PetscSFBasicPack link;
73095fce210SBarry Smith   PetscInt         i,nrootranks;
73195fce210SBarry Smith   const PetscInt   *rootoffset,*rootloc;
73295fce210SBarry Smith 
73395fce210SBarry Smith   PetscFunctionBegin;
73495fce210SBarry Smith   ierr = PetscSFBasicGetPackInUse(sf,unit,rootdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr);
73595fce210SBarry Smith   /* This implementation could be changed to unpack as receives arrive, at the cost of non-determinism */
73695fce210SBarry Smith   ierr = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
73795fce210SBarry Smith   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,NULL,&rootoffset,&rootloc);CHKERRQ(ierr);
73895fce210SBarry Smith   ierr = PetscSFBasicPackGetUnpackOp(sf,link,op,&UnpackOp);CHKERRQ(ierr);
73995fce210SBarry Smith   for (i=0; i<nrootranks; i++) {
74095fce210SBarry Smith     PetscMPIInt n          = rootoffset[i+1] - rootoffset[i];
74195fce210SBarry Smith     const void  *packstart = link->root+rootoffset[i]*link->unitbytes;
74295fce210SBarry Smith 
74395fce210SBarry Smith     (*UnpackOp)(n,rootloc+rootoffset[i],rootdata,packstart);
74495fce210SBarry Smith   }
74595fce210SBarry Smith   ierr = PetscSFBasicReclaimPack(sf,&link);CHKERRQ(ierr);
74695fce210SBarry Smith   PetscFunctionReturn(0);
74795fce210SBarry Smith }
74895fce210SBarry Smith 
74995fce210SBarry Smith #undef __FUNCT__
75095fce210SBarry Smith #define __FUNCT__ "PetscSFFetchAndOpBegin_Basic"
75195fce210SBarry Smith static PetscErrorCode PetscSFFetchAndOpBegin_Basic(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
75295fce210SBarry Smith {
75395fce210SBarry Smith   PetscErrorCode ierr;
75495fce210SBarry Smith 
75595fce210SBarry Smith   PetscFunctionBegin;
75695fce210SBarry Smith   ierr = PetscSFReduceBegin_Basic(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
75795fce210SBarry Smith   PetscFunctionReturn(0);
75895fce210SBarry Smith }
75995fce210SBarry Smith 
76095fce210SBarry Smith #undef __FUNCT__
76195fce210SBarry Smith #define __FUNCT__ "PetscSFFetchAndOpEnd_Basic"
76295fce210SBarry Smith static PetscErrorCode PetscSFFetchAndOpEnd_Basic(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
76395fce210SBarry Smith {
76495fce210SBarry Smith   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;
76595fce210SBarry Smith   void              (*FetchAndOp)(PetscInt,const PetscInt*,void*,void*);
76695fce210SBarry Smith   PetscErrorCode    ierr;
76795fce210SBarry Smith   PetscSFBasicPack  link;
76895fce210SBarry Smith   PetscInt          i,nrootranks,nleafranks;
76995fce210SBarry Smith   const PetscInt    *rootoffset,*leafoffset,*rootloc,*leafloc;
77095fce210SBarry Smith   const PetscMPIInt *rootranks,*leafranks;
77195fce210SBarry Smith   MPI_Request       *rootreqs,*leafreqs;
77295fce210SBarry Smith   size_t            unitbytes;
77395fce210SBarry Smith 
77495fce210SBarry Smith   PetscFunctionBegin;
77595fce210SBarry Smith   ierr = PetscSFBasicGetPackInUse(sf,unit,rootdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr);
77695fce210SBarry Smith   /* This implementation could be changed to unpack as receives arrive, at the cost of non-determinism */
77795fce210SBarry Smith   ierr      = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
77895fce210SBarry Smith   unitbytes = link->unitbytes;
77995fce210SBarry Smith   ierr      = PetscSFBasicGetRootInfo(sf,&nrootranks,&rootranks,&rootoffset,&rootloc);CHKERRQ(ierr);
78095fce210SBarry Smith   ierr      = PetscSFBasicGetLeafInfo(sf,&nleafranks,&leafranks,&leafoffset,&leafloc);CHKERRQ(ierr);
78195fce210SBarry Smith   ierr      = PetscSFBasicPackGetReqs(sf,link,&rootreqs,&leafreqs);CHKERRQ(ierr);
78295fce210SBarry Smith   /* Post leaf receives */
78395fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
78495fce210SBarry Smith     PetscMPIInt n = leafoffset[i+1] - leafoffset[i];
78595fce210SBarry Smith     ierr = MPI_Irecv(link->leaf+leafoffset[i]*unitbytes,n,unit,leafranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&leafreqs[i]);CHKERRQ(ierr);
78695fce210SBarry Smith   }
78795fce210SBarry Smith   /* Process local fetch-and-op, post root sends */
78895fce210SBarry Smith   ierr = PetscSFBasicPackGetFetchAndOp(sf,link,op,&FetchAndOp);CHKERRQ(ierr);
78995fce210SBarry Smith   for (i=0; i<nrootranks; i++) {
79095fce210SBarry Smith     PetscMPIInt n          = rootoffset[i+1] - rootoffset[i];
79195fce210SBarry Smith     void        *packstart = link->root+rootoffset[i]*unitbytes;
79295fce210SBarry Smith 
79395fce210SBarry Smith     (*FetchAndOp)(n,rootloc+rootoffset[i],rootdata,packstart);
79495fce210SBarry Smith     ierr = MPI_Isend(packstart,n,unit,rootranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&rootreqs[i]);CHKERRQ(ierr);
79595fce210SBarry Smith   }
79695fce210SBarry Smith   ierr = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
79795fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
79895fce210SBarry Smith     PetscMPIInt n          = leafoffset[i+1] - leafoffset[i];
79995fce210SBarry Smith     const void  *packstart = link->leaf+leafoffset[i]*unitbytes;
80095fce210SBarry Smith     (*link->UnpackInsert)(n,leafloc+leafoffset[i],leafupdate,packstart);
80195fce210SBarry Smith   }
80295fce210SBarry Smith   ierr = PetscSFBasicReclaimPack(sf,&link);CHKERRQ(ierr);
80395fce210SBarry Smith   PetscFunctionReturn(0);
80495fce210SBarry Smith }
80595fce210SBarry Smith 
80695fce210SBarry Smith #undef __FUNCT__
80795fce210SBarry Smith #define __FUNCT__ "PetscSFCreate_Basic"
808*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscSFCreate_Basic(PetscSF sf)
80995fce210SBarry Smith {
81095fce210SBarry Smith   PetscSF_Basic  *bas = (PetscSF_Basic*)sf->data;
81195fce210SBarry Smith   PetscErrorCode ierr;
81295fce210SBarry Smith 
81395fce210SBarry Smith   PetscFunctionBegin;
81495fce210SBarry Smith   sf->ops->SetUp           = PetscSFSetUp_Basic;
81595fce210SBarry Smith   sf->ops->SetFromOptions  = PetscSFSetFromOptions_Basic;
81695fce210SBarry Smith   sf->ops->Reset           = PetscSFReset_Basic;
81795fce210SBarry Smith   sf->ops->Destroy         = PetscSFDestroy_Basic;
81895fce210SBarry Smith   sf->ops->View            = PetscSFView_Basic;
81995fce210SBarry Smith   sf->ops->BcastBegin      = PetscSFBcastBegin_Basic;
82095fce210SBarry Smith   sf->ops->BcastEnd        = PetscSFBcastEnd_Basic;
82195fce210SBarry Smith   sf->ops->ReduceBegin     = PetscSFReduceBegin_Basic;
82295fce210SBarry Smith   sf->ops->ReduceEnd       = PetscSFReduceEnd_Basic;
82395fce210SBarry Smith   sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Basic;
82495fce210SBarry Smith   sf->ops->FetchAndOpEnd   = PetscSFFetchAndOpEnd_Basic;
82595fce210SBarry Smith 
82695fce210SBarry Smith   ierr     = PetscNewLog(sf,PetscSF_Basic,&bas);CHKERRQ(ierr);
82795fce210SBarry Smith   sf->data = (void*)bas;
82895fce210SBarry Smith   PetscFunctionReturn(0);
82995fce210SBarry Smith }
830