xref: /petsc/src/vec/is/sf/impls/basic/sfbasic.c (revision 2a74e5f4c5bc84d3530a2446235a7edad941f3f8)
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 {
6*2a74e5f4SJed Brown   void (*Pack)(PetscInt,PetscInt,const PetscInt*,const void*,void*);
7*2a74e5f4SJed Brown   void (*UnpackInsert)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
8*2a74e5f4SJed Brown   void (*UnpackAdd)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
9*2a74e5f4SJed Brown   void (*UnpackMin)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
10*2a74e5f4SJed Brown   void (*UnpackMax)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
11*2a74e5f4SJed Brown   void (*UnpackMinloc)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
12*2a74e5f4SJed Brown   void (*UnpackMaxloc)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
13*2a74e5f4SJed Brown   void (*UnpackMult)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
14*2a74e5f4SJed Brown   void (*UnpackLAND)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
15*2a74e5f4SJed Brown   void (*UnpackBAND)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
16*2a74e5f4SJed Brown   void (*UnpackLOR)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
17*2a74e5f4SJed Brown   void (*UnpackBOR)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
18*2a74e5f4SJed Brown   void (*UnpackLXOR)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
19*2a74e5f4SJed Brown   void (*UnpackBXOR)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
20*2a74e5f4SJed Brown   void (*FetchAndInsert)(PetscInt,PetscInt,const PetscInt*,void*,void*);
21*2a74e5f4SJed Brown   void (*FetchAndAdd)(PetscInt,PetscInt,const PetscInt*,void*,void*);
22*2a74e5f4SJed Brown   void (*FetchAndMin)(PetscInt,PetscInt,const PetscInt*,void*,void*);
23*2a74e5f4SJed Brown   void (*FetchAndMax)(PetscInt,PetscInt,const PetscInt*,void*,void*);
24*2a74e5f4SJed Brown   void (*FetchAndMinloc)(PetscInt,PetscInt,const PetscInt*,void*,void*);
25*2a74e5f4SJed Brown   void (*FetchAndMaxloc)(PetscInt,PetscInt,const PetscInt*,void*,void*);
26*2a74e5f4SJed Brown   void (*FetchAndMult)(PetscInt,PetscInt,const PetscInt*,void*,void*);
27*2a74e5f4SJed Brown   void (*FetchAndLAND)(PetscInt,PetscInt,const PetscInt*,void*,void*);
28*2a74e5f4SJed Brown   void (*FetchAndBAND)(PetscInt,PetscInt,const PetscInt*,void*,void*);
29*2a74e5f4SJed Brown   void (*FetchAndLOR)(PetscInt,PetscInt,const PetscInt*,void*,void*);
30*2a74e5f4SJed Brown   void (*FetchAndBOR)(PetscInt,PetscInt,const PetscInt*,void*,void*);
31*2a74e5f4SJed Brown   void (*FetchAndLXOR)(PetscInt,PetscInt,const PetscInt*,void*,void*);
32*2a74e5f4SJed Brown   void (*FetchAndBXOR)(PetscInt,PetscInt,const PetscInt*,void*,void*);
3395fce210SBarry Smith 
3495fce210SBarry Smith   MPI_Datatype     unit;
3595fce210SBarry Smith   size_t           unitbytes;   /* Number of bytes in a unit */
36*2a74e5f4SJed Brown   PetscInt         bs;          /* Number of basic units in a unit */
3795fce210SBarry Smith   const void       *key;        /* Array used as key for operation */
3895fce210SBarry Smith   char             *root;       /* Packed root data, contiguous by leaf rank */
3995fce210SBarry Smith   char             *leaf;       /* Packed leaf data, contiguous by root rank */
4095fce210SBarry Smith   MPI_Request      *requests;   /* Array of root requests followed by leaf requests */
4195fce210SBarry Smith   PetscSFBasicPack next;
4295fce210SBarry Smith };
4395fce210SBarry Smith 
4495fce210SBarry Smith typedef struct {
4595fce210SBarry Smith   PetscMPIInt      tag;
4695fce210SBarry Smith   PetscInt         niranks;     /* Number of incoming ranks (ranks accessing my roots) */
4795fce210SBarry Smith   PetscMPIInt      *iranks;     /* Array of ranks that reference my roots */
4895fce210SBarry Smith   PetscInt         itotal;      /* Total number of graph edges referencing my roots */
4995fce210SBarry Smith   PetscInt         *ioffset;    /* Array of length niranks+1 holding offset in irootloc[] for each rank */
5095fce210SBarry Smith   PetscInt         *irootloc;   /* Incoming roots referenced by ranks starting at ioffset[rank] */
5195fce210SBarry Smith   PetscSFBasicPack avail;       /* One or more entries per MPI Datatype, lazily constructed */
5295fce210SBarry Smith   PetscSFBasicPack inuse;       /* Buffers being used for transactions that have not yet completed */
5395fce210SBarry Smith } PetscSF_Basic;
5495fce210SBarry Smith 
5595fce210SBarry Smith #if !defined(PETSC_HAVE_MPI_TYPE_DUP) /* Danger: type is not reference counted; subject to ABA problem */
5695fce210SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MPI_Type_dup(MPI_Datatype datatype,MPI_Datatype *newtype)
5795fce210SBarry Smith {
5895fce210SBarry Smith   *newtype = datatype;
5995fce210SBarry Smith   return 0;
6095fce210SBarry Smith }
6195fce210SBarry Smith #endif
6295fce210SBarry Smith 
6395fce210SBarry Smith /*
6495fce210SBarry Smith  * MPI_Reduce_local is not really useful because it can't handle sparse data and it vectorizes "in the wrong direction",
6595fce210SBarry Smith  * therefore we pack data types manually. This section defines packing routines for the standard data types.
6695fce210SBarry Smith  */
6795fce210SBarry Smith 
6895fce210SBarry Smith #define CPPJoin2_exp(a,b) a ## b
6995fce210SBarry Smith #define CPPJoin2(a,b) CPPJoin2_exp(a,b)
7095fce210SBarry Smith #define CPPJoin3_exp_(a,b,c) a ## b ## _ ## c
7195fce210SBarry Smith #define CPPJoin3_(a,b,c) CPPJoin3_exp_(a,b,c)
7295fce210SBarry Smith 
7395fce210SBarry Smith /* Basic types without addition */
74*2a74e5f4SJed Brown #define DEF_PackNoInit(type,BS)                                         \
75*2a74e5f4SJed Brown   static void CPPJoin3_(Pack_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,const void *unpacked,void *packed) { \
7695fce210SBarry Smith     const type *u = (const type*)unpacked;                              \
7795fce210SBarry Smith     type *p = (type*)packed;                                            \
78*2a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
79*2a74e5f4SJed Brown     for (i=0; i<n; i++)                                                 \
80*2a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
81*2a74e5f4SJed Brown         for (k=j; k<j+BS; k++)                                          \
82*2a74e5f4SJed Brown           p[i*bs+k] = u[idx[i]*bs+k];                                   \
8395fce210SBarry Smith   }                                                                     \
84*2a74e5f4SJed Brown   static void CPPJoin3_(UnpackInsert_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
8595fce210SBarry Smith     type *u = (type*)unpacked;                                          \
8695fce210SBarry Smith     const type *p = (const type*)packed;                                \
87*2a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
88*2a74e5f4SJed Brown     for (i=0; i<n; i++)                                                 \
89*2a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
90*2a74e5f4SJed Brown         for (k=j; k<j+BS; k++)                                          \
91*2a74e5f4SJed Brown           u[idx[i]*bs+k] = p[i*bs+k];                                   \
9295fce210SBarry Smith   }                                                                     \
93*2a74e5f4SJed Brown   static void CPPJoin3_(FetchAndInsert_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
9495fce210SBarry Smith     type *u = (type*)unpacked;                                          \
9595fce210SBarry Smith     type *p = (type*)packed;                                            \
96*2a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
9795fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
98*2a74e5f4SJed Brown       PetscInt ii = idx[i];                                             \
99*2a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
100*2a74e5f4SJed Brown         for (k=j; k<j+BS; k++) {                                        \
101*2a74e5f4SJed Brown           type t = u[ii*bs+k];                                          \
102*2a74e5f4SJed Brown           u[ii*bs+k] = p[i*bs+k];                                       \
103*2a74e5f4SJed Brown           p[i*bs+k] = t;                                                \
104*2a74e5f4SJed Brown         }                                                               \
10595fce210SBarry Smith     }                                                                   \
10695fce210SBarry Smith   }
10795fce210SBarry Smith 
10895fce210SBarry Smith /* Basic types defining addition */
109*2a74e5f4SJed Brown #define DEF_PackAddNoInit(type,BS)                                      \
110*2a74e5f4SJed Brown   DEF_PackNoInit(type,BS)                                               \
111*2a74e5f4SJed Brown   static void CPPJoin3_(UnpackAdd_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
11295fce210SBarry Smith     type *u = (type*)unpacked;                                          \
11395fce210SBarry Smith     const type *p = (const type*)packed;                                \
114*2a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
115*2a74e5f4SJed Brown     for (i=0; i<n; i++)                                                 \
116*2a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
117*2a74e5f4SJed Brown         for (k=j; k<j+BS; k++)                                          \
118*2a74e5f4SJed Brown           u[idx[i]*bs+k] += p[i*bs+k];                                  \
11995fce210SBarry Smith   }                                                                     \
120*2a74e5f4SJed Brown   static void CPPJoin3_(FetchAndAdd_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
12195fce210SBarry Smith     type *u = (type*)unpacked;                                          \
12295fce210SBarry Smith     type *p = (type*)packed;                                            \
123*2a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
12495fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
125*2a74e5f4SJed Brown       PetscInt ii = idx[i];                                             \
126*2a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
127*2a74e5f4SJed Brown         for (k=j; k<j+BS; k++) {                                        \
128*2a74e5f4SJed Brown           type t = u[ii*bs+k];                                          \
129*2a74e5f4SJed Brown           u[ii*bs+k] = t + p[i*bs+k];                                   \
130*2a74e5f4SJed Brown           p[i*bs+k] = t;                                                \
13195fce210SBarry Smith         }                                                               \
13253deab39SPeter Brune     }                                                                   \
133*2a74e5f4SJed Brown   }                                                                     \
134*2a74e5f4SJed Brown   static void CPPJoin3_(UnpackMult_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
13553deab39SPeter Brune     type *u = (type*)unpacked;                                          \
13653deab39SPeter Brune     const type *p = (const type*)packed;                                \
137*2a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
138*2a74e5f4SJed Brown     for (i=0; i<n; i++)                                                 \
139*2a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
140*2a74e5f4SJed Brown         for (k=j; k<j+BS; k++)                                          \
141*2a74e5f4SJed Brown           u[idx[i]*bs+k] *= p[i*bs+k];                                  \
14253deab39SPeter Brune   }                                                                     \
143*2a74e5f4SJed Brown   static void CPPJoin3_(FetchAndMult_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
14453deab39SPeter Brune     type *u = (type*)unpacked;                                          \
14553deab39SPeter Brune     type *p = (type*)packed;                                            \
146*2a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
14753deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
148*2a74e5f4SJed Brown       PetscInt ii = idx[i];                                             \
149*2a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
150*2a74e5f4SJed Brown         for (k=j; k<j+BS; k++) {                                        \
151*2a74e5f4SJed Brown           type t = u[ii*bs+k];                                          \
152*2a74e5f4SJed Brown           u[ii*bs+k] = t * p[i*bs+k];                                   \
153*2a74e5f4SJed Brown           p[i*bs+k] = t;                                                \
154*2a74e5f4SJed Brown         }                                                               \
15553deab39SPeter Brune     }                                                                   \
15695fce210SBarry Smith   }
157*2a74e5f4SJed Brown #define DEF_Pack(type,BS)                                               \
158*2a74e5f4SJed Brown   DEF_PackAddNoInit(type,BS)                                            \
159*2a74e5f4SJed Brown   static void CPPJoin3_(PackInit_,type,BS)(PetscSFBasicPack link) {     \
160*2a74e5f4SJed Brown     link->Pack = CPPJoin3_(Pack_,type,BS);                              \
161*2a74e5f4SJed Brown     link->UnpackInsert = CPPJoin3_(UnpackInsert_,type,BS);              \
162*2a74e5f4SJed Brown     link->UnpackAdd = CPPJoin3_(UnpackAdd_,type,BS);                    \
163*2a74e5f4SJed Brown     link->UnpackMult = CPPJoin3_(UnpackMult_,type,BS);                  \
164*2a74e5f4SJed Brown     link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,type,BS);          \
165*2a74e5f4SJed Brown     link->FetchAndAdd = CPPJoin3_(FetchAndAdd_,type,BS);                \
166*2a74e5f4SJed Brown     link->FetchAndMult = CPPJoin3_(FetchAndMult_,type,BS);              \
16795fce210SBarry Smith     link->unitbytes = sizeof(type);                                     \
16895fce210SBarry Smith   }
16995fce210SBarry Smith /* Comparable types */
17095fce210SBarry Smith #define DEF_PackCmp(type)                                               \
171*2a74e5f4SJed Brown   DEF_PackAddNoInit(type,1)                                             \
172*2a74e5f4SJed Brown   static void CPPJoin2(UnpackMax_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
17395fce210SBarry Smith     type *u = (type*)unpacked;                                          \
17495fce210SBarry Smith     const type *p = (const type*)packed;                                \
17595fce210SBarry Smith     PetscInt i;                                                         \
17695fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
17795fce210SBarry Smith       type v = u[idx[i]];                                               \
17895fce210SBarry Smith       u[idx[i]] = PetscMax(v,p[i]);                                     \
17995fce210SBarry Smith     }                                                                   \
18095fce210SBarry Smith   }                                                                     \
181*2a74e5f4SJed Brown   static void CPPJoin2(UnpackMin_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
18295fce210SBarry Smith     type *u = (type*)unpacked;                                          \
18395fce210SBarry Smith     const type *p = (const type*)packed;                                \
18495fce210SBarry Smith     PetscInt i;                                                         \
18595fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
18695fce210SBarry Smith       type v = u[idx[i]];                                               \
18795fce210SBarry Smith       u[idx[i]] = PetscMin(v,p[i]);                                     \
18895fce210SBarry Smith     }                                                                   \
18995fce210SBarry Smith   }                                                                     \
190*2a74e5f4SJed Brown   static void CPPJoin2(FetchAndMax_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
19195fce210SBarry Smith     type *u = (type*)unpacked;                                          \
19295fce210SBarry Smith     type *p = (type*)packed;                                            \
19395fce210SBarry Smith     PetscInt i;                                                         \
19495fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
19595fce210SBarry Smith       PetscInt j = idx[i];                                              \
19695fce210SBarry Smith       type v = u[j];                                                    \
19795fce210SBarry Smith       u[j] = PetscMax(v,p[i]);                                          \
19895fce210SBarry Smith       p[i] = v;                                                         \
19995fce210SBarry Smith     }                                                                   \
20095fce210SBarry Smith   }                                                                     \
201*2a74e5f4SJed Brown   static void CPPJoin2(FetchAndMin_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
20295fce210SBarry Smith     type *u = (type*)unpacked;                                          \
20395fce210SBarry Smith     type *p = (type*)packed;                                            \
20495fce210SBarry Smith     PetscInt i;                                                         \
20595fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
20695fce210SBarry Smith       PetscInt j = idx[i];                                              \
20795fce210SBarry Smith       type v = u[j];                                                    \
20895fce210SBarry Smith       u[j] = PetscMin(v,p[i]);                                          \
20995fce210SBarry Smith       p[i] = v;                                                         \
21095fce210SBarry Smith     }                                                                   \
21195fce210SBarry Smith   }                                                                     \
21295fce210SBarry Smith   static void CPPJoin2(PackInit_,type)(PetscSFBasicPack link) {         \
213*2a74e5f4SJed Brown     link->Pack = CPPJoin3_(Pack_,type,1);                               \
214*2a74e5f4SJed Brown     link->UnpackInsert = CPPJoin3_(UnpackInsert_,type,1);               \
215*2a74e5f4SJed Brown     link->UnpackAdd  = CPPJoin3_(UnpackAdd_,type,1);                    \
21695fce210SBarry Smith     link->UnpackMax  = CPPJoin2(UnpackMax_,type);                       \
21795fce210SBarry Smith     link->UnpackMin  = CPPJoin2(UnpackMin_,type);                       \
218*2a74e5f4SJed Brown     link->UnpackMult = CPPJoin3_(UnpackMult_,type,1);                   \
219*2a74e5f4SJed Brown     link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,type,1);           \
220*2a74e5f4SJed Brown     link->FetchAndAdd = CPPJoin3_(FetchAndAdd_ ,type,1);                \
22195fce210SBarry Smith     link->FetchAndMax = CPPJoin2(FetchAndMax_ ,type);                   \
22295fce210SBarry Smith     link->FetchAndMin = CPPJoin2(FetchAndMin_ ,type);                   \
223*2a74e5f4SJed Brown     link->FetchAndMult = CPPJoin3_(FetchAndMult_,type,1);               \
22495fce210SBarry Smith     link->unitbytes = sizeof(type);                                     \
22595fce210SBarry Smith   }
22695fce210SBarry Smith 
22753deab39SPeter Brune /* Logical Types */
22853deab39SPeter Brune #define DEF_PackLog(type)                                               \
229*2a74e5f4SJed Brown   static void CPPJoin2(UnpackLAND_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
23053deab39SPeter Brune     type *u = (type*)unpacked;                                          \
23153deab39SPeter Brune     const type *p = (const type*)packed;                                \
23253deab39SPeter Brune     PetscInt i;                                                         \
23353deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
23453deab39SPeter Brune       type v = u[idx[i]];                                               \
23553deab39SPeter Brune       u[idx[i]] = v && p[i];                                            \
23653deab39SPeter Brune     }                                                                   \
23753deab39SPeter Brune   }                                                                     \
238*2a74e5f4SJed Brown   static void CPPJoin2(UnpackLOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
23953deab39SPeter Brune     type *u = (type*)unpacked;                                          \
24053deab39SPeter Brune     const type *p = (const type*)packed;                                \
24153deab39SPeter Brune     PetscInt i;                                                         \
24253deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
24353deab39SPeter Brune       type v = u[idx[i]];                                               \
24453deab39SPeter Brune       u[idx[i]] = v || p[i];                                            \
24553deab39SPeter Brune     }                                                                   \
24653deab39SPeter Brune   }                                                                     \
247*2a74e5f4SJed Brown   static void CPPJoin2(UnpackLXOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
24853deab39SPeter Brune     type *u = (type*)unpacked;                                          \
24953deab39SPeter Brune     const type *p = (const type*)packed;                                \
25053deab39SPeter Brune     PetscInt i;                                                         \
25153deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
25253deab39SPeter Brune       type v = u[idx[i]];                                               \
25353deab39SPeter Brune       u[idx[i]] = (!v)!=(!p[i]);                                        \
25453deab39SPeter Brune     }                                                                   \
25553deab39SPeter Brune   }                                                                     \
256*2a74e5f4SJed Brown   static void CPPJoin2(FetchAndLAND_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
25753deab39SPeter Brune     type *u = (type*)unpacked;                                          \
25853deab39SPeter Brune     type *p = (type*)packed;                                            \
25953deab39SPeter Brune     PetscInt i;                                                         \
26053deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
26153deab39SPeter Brune       PetscInt j = idx[i];                                              \
26253deab39SPeter Brune       type v = u[j];                                                    \
26353deab39SPeter Brune       u[j] = v && p[i];                                                 \
26453deab39SPeter Brune       p[i] = v;                                                         \
26553deab39SPeter Brune     }                                                                   \
26653deab39SPeter Brune   }                                                                     \
267*2a74e5f4SJed Brown   static void CPPJoin2(FetchAndLOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
26853deab39SPeter Brune     type *u = (type*)unpacked;                                          \
26953deab39SPeter Brune     type *p = (type*)packed;                                            \
27053deab39SPeter Brune     PetscInt i;                                                         \
27153deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
27253deab39SPeter Brune       PetscInt j = idx[i];                                              \
27353deab39SPeter Brune       type v = u[j];                                                    \
27453deab39SPeter Brune       u[j] = v || p[i];                                                 \
27553deab39SPeter Brune       p[i] = v;                                                         \
27653deab39SPeter Brune     }                                                                   \
27753deab39SPeter Brune   }                                                                     \
278*2a74e5f4SJed Brown   static void CPPJoin2(FetchAndLXOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
27953deab39SPeter Brune     type *u = (type*)unpacked;                                          \
28053deab39SPeter Brune     type *p = (type*)packed;                                            \
28153deab39SPeter Brune     PetscInt i;                                                         \
28253deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
28353deab39SPeter Brune       PetscInt j = idx[i];                                              \
28453deab39SPeter Brune       type v = u[j];                                                    \
28553deab39SPeter Brune       u[j] = (!v)!=(!p[i]);                                             \
28653deab39SPeter Brune       p[i] = v;                                                         \
28753deab39SPeter Brune     }                                                                   \
28853deab39SPeter Brune   }                                                                     \
28953deab39SPeter Brune   static void CPPJoin2(PackInit_Logical_,type)(PetscSFBasicPack link) { \
29053deab39SPeter Brune     link->UnpackLAND = CPPJoin2(UnpackLAND_,type);                      \
29153deab39SPeter Brune     link->UnpackLOR  = CPPJoin2(UnpackLOR_,type);                       \
29253deab39SPeter Brune     link->UnpackLXOR = CPPJoin2(UnpackLXOR_,type);                      \
29353deab39SPeter Brune     link->FetchAndLAND = CPPJoin2(FetchAndLAND_,type);                  \
29453deab39SPeter Brune     link->FetchAndLOR  = CPPJoin2(FetchAndLOR_,type);                   \
29553deab39SPeter Brune     link->FetchAndLXOR = CPPJoin2(FetchAndLXOR_,type);                  \
29653deab39SPeter Brune   }
29753deab39SPeter Brune 
29853deab39SPeter Brune 
29953deab39SPeter Brune /* Bitwise Types */
30053deab39SPeter Brune #define DEF_PackBit(type)                                               \
301*2a74e5f4SJed Brown   static void CPPJoin2(UnpackBAND_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
30253deab39SPeter Brune     type *u = (type*)unpacked;                                          \
30353deab39SPeter Brune     const type *p = (const type*)packed;                                \
30453deab39SPeter Brune     PetscInt i;                                                         \
30553deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
30653deab39SPeter Brune       type v = u[idx[i]];                                               \
30753deab39SPeter Brune       u[idx[i]] = v & p[i];                                             \
30853deab39SPeter Brune     }                                                                   \
30953deab39SPeter Brune   }                                                                     \
310*2a74e5f4SJed Brown   static void CPPJoin2(UnpackBOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
31153deab39SPeter Brune     type *u = (type*)unpacked;                                          \
31253deab39SPeter Brune     const type *p = (const type*)packed;                                \
31353deab39SPeter Brune     PetscInt i;                                                         \
31453deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
31553deab39SPeter Brune       type v = u[idx[i]];                                               \
31653deab39SPeter Brune       u[idx[i]] = v | p[i];                                             \
31753deab39SPeter Brune     }                                                                   \
31853deab39SPeter Brune   }                                                                     \
319*2a74e5f4SJed Brown   static void CPPJoin2(UnpackBXOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
32053deab39SPeter Brune     type *u = (type*)unpacked;                                          \
32153deab39SPeter Brune     const type *p = (const type*)packed;                                \
32253deab39SPeter Brune     PetscInt i;                                                         \
32353deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
32453deab39SPeter Brune       type v = u[idx[i]];                                               \
32553deab39SPeter Brune       u[idx[i]] = v^p[i];                                               \
32653deab39SPeter Brune     }                                                                   \
32753deab39SPeter Brune   }                                                                     \
328*2a74e5f4SJed Brown   static void CPPJoin2(FetchAndBAND_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
32953deab39SPeter Brune     type *u = (type*)unpacked;                                          \
33053deab39SPeter Brune     type *p = (type*)packed;                                            \
33153deab39SPeter Brune     PetscInt i;                                                         \
33253deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
33353deab39SPeter Brune       PetscInt j = idx[i];                                              \
33453deab39SPeter Brune       type v = u[j];                                                    \
33553deab39SPeter Brune       u[j] = v & p[i];                                                  \
33653deab39SPeter Brune       p[i] = v;                                                         \
33753deab39SPeter Brune     }                                                                   \
33853deab39SPeter Brune   }                                                                     \
339*2a74e5f4SJed Brown   static void CPPJoin2(FetchAndBOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
34053deab39SPeter Brune     type *u = (type*)unpacked;                                          \
34153deab39SPeter Brune     type *p = (type*)packed;                                            \
34253deab39SPeter Brune     PetscInt i;                                                         \
34353deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
34453deab39SPeter Brune       PetscInt j = idx[i];                                              \
34553deab39SPeter Brune       type v = u[j];                                                    \
34653deab39SPeter Brune       u[j] = v | p[i];                                                  \
34753deab39SPeter Brune       p[i] = v;                                                         \
34853deab39SPeter Brune     }                                                                   \
34953deab39SPeter Brune   }                                                                     \
350*2a74e5f4SJed Brown   static void CPPJoin2(FetchAndBXOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
35153deab39SPeter Brune     type *u = (type*)unpacked;                                          \
35253deab39SPeter Brune     type *p = (type*)packed;                                            \
35353deab39SPeter Brune     PetscInt i;                                                         \
35453deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
35553deab39SPeter Brune       PetscInt j = idx[i];                                              \
35653deab39SPeter Brune       type v = u[j];                                                    \
35753deab39SPeter Brune       u[j] = v^p[i];                                                    \
35853deab39SPeter Brune       p[i] = v;                                                         \
35953deab39SPeter Brune     }                                                                   \
36053deab39SPeter Brune   }                                                                     \
36153deab39SPeter Brune   static void CPPJoin2(PackInit_Bitwise_,type)(PetscSFBasicPack link) { \
36253deab39SPeter Brune     link->UnpackBAND = CPPJoin2(UnpackBAND_,type);                      \
36353deab39SPeter Brune     link->UnpackBOR  = CPPJoin2(UnpackBOR_,type);                       \
36453deab39SPeter Brune     link->UnpackBXOR = CPPJoin2(UnpackBXOR_,type);                      \
36553deab39SPeter Brune     link->FetchAndBAND = CPPJoin2(FetchAndBAND_,type);                  \
36653deab39SPeter Brune     link->FetchAndBOR  = CPPJoin2(FetchAndBOR_,type);                   \
36753deab39SPeter Brune     link->FetchAndBXOR = CPPJoin2(FetchAndBXOR_,type);                  \
36853deab39SPeter Brune   }
36953deab39SPeter Brune 
37095fce210SBarry Smith /* Pair types */
37195fce210SBarry Smith #define CPPJoinloc_exp(base,op,t1,t2) base ## op ## loc_ ## t1 ## _ ## t2
37295fce210SBarry Smith #define CPPJoinloc(base,op,t1,t2) CPPJoinloc_exp(base,op,t1,t2)
37395fce210SBarry Smith #define PairType(type1,type2) CPPJoin3_(_pairtype_,type1,type2)
37495fce210SBarry Smith #define DEF_UnpackXloc(type1,type2,locname,op)                              \
375*2a74e5f4SJed Brown   static void CPPJoinloc(Unpack,locname,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
37695fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;        \
37795fce210SBarry Smith     const PairType(type1,type2) *p = (const PairType(type1,type2)*)packed; \
37895fce210SBarry Smith     PetscInt i;                                                         \
37995fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
38095fce210SBarry Smith       PetscInt j = idx[i];                                              \
38195fce210SBarry Smith       if (p[i].a op u[j].a) {                                           \
38295fce210SBarry Smith         u[j].a = p[i].a;                                                \
38395fce210SBarry Smith         u[j].b = p[i].b;                                                \
38495fce210SBarry Smith       } else if (u[j].a == p[i].a) {                                    \
38595fce210SBarry Smith         u[j].b = PetscMin(u[j].b,p[i].b);                               \
38695fce210SBarry Smith       }                                                                 \
38795fce210SBarry Smith     }                                                                   \
38895fce210SBarry Smith   }                                                                     \
389*2a74e5f4SJed Brown   static void CPPJoinloc(FetchAnd,locname,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
39095fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;        \
39195fce210SBarry Smith     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;          \
39295fce210SBarry Smith     PetscInt i;                                                         \
39395fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
39495fce210SBarry Smith       PetscInt j = idx[i];                                              \
39595fce210SBarry Smith       PairType(type1,type2) v;                                          \
39695fce210SBarry Smith       v.a = u[j].a;                                                     \
39795fce210SBarry Smith       v.b = u[j].b;                                                     \
39895fce210SBarry Smith       if (p[i].a op u[j].a) {                                           \
39995fce210SBarry Smith         u[j].a = p[i].a;                                                \
40095fce210SBarry Smith         u[j].b = p[i].b;                                                \
40195fce210SBarry Smith       } else if (u[j].a == p[i].a) {                                    \
40295fce210SBarry Smith         u[j].b = PetscMin(u[j].b,p[i].b);                               \
40395fce210SBarry Smith       }                                                                 \
40495fce210SBarry Smith       p[i].a = v.a;                                                     \
40595fce210SBarry Smith       p[i].b = v.b;                                                     \
40695fce210SBarry Smith     }                                                                   \
40795fce210SBarry Smith   }
40895fce210SBarry Smith #define DEF_PackPair(type1,type2)                                       \
40995fce210SBarry Smith   typedef struct {type1 a; type2 b;} PairType(type1,type2);             \
410*2a74e5f4SJed Brown   static void CPPJoin3_(Pack_,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,const void *unpacked,void *packed) { \
41195fce210SBarry Smith     const PairType(type1,type2) *u = (const PairType(type1,type2)*)unpacked; \
41295fce210SBarry Smith     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;          \
41395fce210SBarry Smith     PetscInt i;                                                         \
41495fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
41595fce210SBarry Smith       p[i].a = u[idx[i]].a;                                             \
41695fce210SBarry Smith       p[i].b = u[idx[i]].b;                                             \
41795fce210SBarry Smith     }                                                                   \
41895fce210SBarry Smith   }                                                                     \
419*2a74e5f4SJed Brown   static void CPPJoin3_(UnpackInsert_,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
42095fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;       \
42195fce210SBarry Smith     const PairType(type1,type2) *p = (const PairType(type1,type2)*)packed; \
42295fce210SBarry Smith     PetscInt i;                                                         \
42395fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
42495fce210SBarry Smith       u[idx[i]].a = p[i].a;                                             \
42595fce210SBarry Smith       u[idx[i]].b = p[i].b;                                             \
42695fce210SBarry Smith     }                                                                   \
42795fce210SBarry Smith   }                                                                     \
428*2a74e5f4SJed Brown   static void CPPJoin3_(UnpackAdd_,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
42995fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;       \
43095fce210SBarry Smith     const PairType(type1,type2) *p = (const PairType(type1,type2)*)packed; \
43195fce210SBarry Smith     PetscInt i;                                                         \
43295fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
43395fce210SBarry Smith       u[idx[i]].a += p[i].a;                                            \
43495fce210SBarry Smith       u[idx[i]].b += p[i].b;                                            \
43595fce210SBarry Smith     }                                                                   \
43695fce210SBarry Smith   }                                                                     \
437*2a74e5f4SJed Brown   static void CPPJoin3_(FetchAndInsert_,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
43895fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;        \
43995fce210SBarry Smith     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;          \
44095fce210SBarry Smith     PetscInt i;                                                         \
44195fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
44295fce210SBarry Smith       PetscInt j = idx[i];                                              \
44395fce210SBarry Smith       PairType(type1,type2) v;                                          \
44495fce210SBarry Smith       v.a = u[j].a;                                                     \
44595fce210SBarry Smith       v.b = u[j].b;                                                     \
44695fce210SBarry Smith       u[j].a = p[i].a;                                                  \
44795fce210SBarry Smith       u[j].b = p[i].b;                                                  \
44895fce210SBarry Smith       p[i].a = v.a;                                                     \
44995fce210SBarry Smith       p[i].b = v.b;                                                     \
45095fce210SBarry Smith     }                                                                   \
45195fce210SBarry Smith   }                                                                     \
452*2a74e5f4SJed Brown   static void FetchAndAdd_ ## type1 ## _ ## type2(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
45395fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;       \
45495fce210SBarry Smith     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;         \
45595fce210SBarry Smith     PetscInt i;                                                         \
45695fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
45795fce210SBarry Smith       PetscInt j = idx[i];                                              \
45895fce210SBarry Smith       PairType(type1,type2) v;                                          \
45995fce210SBarry Smith       v.a = u[j].a;                                                     \
46095fce210SBarry Smith       v.b = u[j].b;                                                     \
46195fce210SBarry Smith       u[j].a = v.a + p[i].a;                                            \
46295fce210SBarry Smith       u[j].b = v.b + p[i].b;                                            \
46395fce210SBarry Smith       p[i].a = v.a;                                                     \
46495fce210SBarry Smith       p[i].b = v.b;                                                     \
46595fce210SBarry Smith     }                                                                   \
46695fce210SBarry Smith   }                                                                     \
46795fce210SBarry Smith   DEF_UnpackXloc(type1,type2,Max,>)                                     \
46895fce210SBarry Smith   DEF_UnpackXloc(type1,type2,Min,<)                                     \
46995fce210SBarry Smith   static void CPPJoin3_(PackInit_,type1,type2)(PetscSFBasicPack link) { \
47095fce210SBarry Smith     link->Pack = CPPJoin3_(Pack_,type1,type2);                          \
47195fce210SBarry Smith     link->UnpackInsert = CPPJoin3_(UnpackInsert_,type1,type2);          \
47295fce210SBarry Smith     link->UnpackAdd = CPPJoin3_(UnpackAdd_,type1,type2);                \
47395fce210SBarry Smith     link->UnpackMaxloc = CPPJoin3_(UnpackMaxloc_,type1,type2);          \
47495fce210SBarry Smith     link->UnpackMinloc = CPPJoin3_(UnpackMinloc_,type1,type2);          \
47595fce210SBarry Smith     link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,type1,type2);      \
47695fce210SBarry Smith     link->FetchAndAdd = CPPJoin3_(FetchAndAdd_,type1,type2);            \
47795fce210SBarry Smith     link->FetchAndMaxloc = CPPJoin3_(FetchAndMaxloc_,type1,type2);      \
47895fce210SBarry Smith     link->FetchAndMinloc = CPPJoin3_(FetchAndMinloc_,type1,type2);      \
47995fce210SBarry Smith     link->unitbytes = sizeof(PairType(type1,type2));                    \
48095fce210SBarry Smith   }
48195fce210SBarry Smith 
48295fce210SBarry Smith /* Currently only dumb blocks of data */
48395fce210SBarry Smith #define BlockType(unit,count) CPPJoin3_(_blocktype_,unit,count)
48495fce210SBarry Smith #define DEF_Block(unit,count)                                           \
48595fce210SBarry Smith   typedef struct {unit v[count];} BlockType(unit,count);                \
486*2a74e5f4SJed Brown   DEF_PackNoInit(BlockType(unit,count),1)                               \
48795fce210SBarry Smith   static void CPPJoin3_(PackInit_block_,unit,count)(PetscSFBasicPack link) { \
488*2a74e5f4SJed Brown     link->Pack = CPPJoin3_(Pack_,BlockType(unit,count),1);               \
489*2a74e5f4SJed Brown     link->UnpackInsert = CPPJoin3_(UnpackInsert_,BlockType(unit,count),1); \
490*2a74e5f4SJed Brown     link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,BlockType(unit,count),1); \
49195fce210SBarry Smith     link->unitbytes = sizeof(BlockType(unit,count));                    \
49295fce210SBarry Smith   }
49395fce210SBarry Smith 
49495fce210SBarry Smith DEF_PackCmp(int)
49553deab39SPeter Brune DEF_PackBit(int)
49653deab39SPeter Brune DEF_PackLog(int)
49795fce210SBarry Smith DEF_PackCmp(PetscInt)
49853deab39SPeter Brune DEF_PackBit(PetscInt)
49953deab39SPeter Brune DEF_PackLog(PetscInt)
500*2a74e5f4SJed Brown DEF_Pack(PetscInt,2)
501*2a74e5f4SJed Brown DEF_Pack(PetscInt,3)
502*2a74e5f4SJed Brown DEF_Pack(PetscInt,4)
503*2a74e5f4SJed Brown DEF_Pack(PetscInt,5)
504*2a74e5f4SJed Brown DEF_Pack(PetscInt,7)
50595fce210SBarry Smith DEF_PackCmp(PetscReal)
50653deab39SPeter Brune DEF_PackLog(PetscReal)
507*2a74e5f4SJed Brown DEF_Pack(PetscReal,2)
508*2a74e5f4SJed Brown DEF_Pack(PetscReal,3)
509*2a74e5f4SJed Brown DEF_Pack(PetscReal,4)
510*2a74e5f4SJed Brown DEF_Pack(PetscReal,5)
511*2a74e5f4SJed Brown DEF_Pack(PetscReal,7)
51295fce210SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
513*2a74e5f4SJed Brown DEF_Pack(PetscComplex,1)
514*2a74e5f4SJed Brown DEF_Pack(PetscComplex,2)
515*2a74e5f4SJed Brown DEF_Pack(PetscComplex,3)
516*2a74e5f4SJed Brown DEF_Pack(PetscComplex,4)
517*2a74e5f4SJed Brown DEF_Pack(PetscComplex,5)
518*2a74e5f4SJed Brown DEF_Pack(PetscComplex,7)
51995fce210SBarry Smith #endif
52095fce210SBarry Smith DEF_PackPair(int,int)
52195fce210SBarry Smith DEF_PackPair(PetscInt,PetscInt)
52295fce210SBarry Smith DEF_Block(int,1)
52395fce210SBarry Smith DEF_Block(int,2)
52495fce210SBarry Smith DEF_Block(int,3)
52595fce210SBarry Smith DEF_Block(int,4)
52695fce210SBarry Smith DEF_Block(int,5)
52795fce210SBarry Smith DEF_Block(int,6)
52895fce210SBarry Smith DEF_Block(int,7)
52995fce210SBarry Smith DEF_Block(int,8)
53095fce210SBarry Smith 
53195fce210SBarry Smith #undef __FUNCT__
53295fce210SBarry Smith #define __FUNCT__ "PetscSFSetUp_Basic"
53395fce210SBarry Smith static PetscErrorCode PetscSFSetUp_Basic(PetscSF sf)
53495fce210SBarry Smith {
53595fce210SBarry Smith   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
53695fce210SBarry Smith   PetscErrorCode ierr;
53795fce210SBarry Smith   PetscInt *rlengths,*ilengths,i;
53895fce210SBarry Smith   MPI_Comm comm;
53995fce210SBarry Smith   MPI_Request *rootreqs,*leafreqs;
54095fce210SBarry Smith 
54195fce210SBarry Smith   PetscFunctionBegin;
54295fce210SBarry Smith   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
54395fce210SBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)sf,&bas->tag);CHKERRQ(ierr);
54495fce210SBarry Smith   /*
54595fce210SBarry Smith    * Inform roots about how many leaves and from which ranks
54695fce210SBarry Smith    */
547785e854fSJed Brown   ierr = PetscMalloc1(sf->nranks,&rlengths);CHKERRQ(ierr);
54895fce210SBarry Smith   /* Determine number, sending ranks, and length of incoming  */
54995fce210SBarry Smith   for (i=0; i<sf->nranks; i++) {
55095fce210SBarry Smith     rlengths[i] = sf->roffset[i+1] - sf->roffset[i]; /* Number of roots referenced by my leaves; for rank sf->ranks[i] */
55195fce210SBarry Smith   }
55295fce210SBarry Smith   ierr = PetscCommBuildTwoSided(comm,1,MPIU_INT,sf->nranks,sf->ranks,rlengths,&bas->niranks,&bas->iranks,(void**)&ilengths);CHKERRQ(ierr);
55395fce210SBarry Smith   ierr = PetscFree(rlengths);CHKERRQ(ierr);
55495fce210SBarry Smith 
55595fce210SBarry Smith   /* Send leaf identities to roots */
55695fce210SBarry Smith   for (i=0,bas->itotal=0; i<bas->niranks; i++) bas->itotal += ilengths[i];
557dcca6d9dSJed Brown   ierr = PetscMalloc2(bas->niranks+1,&bas->ioffset,bas->itotal,&bas->irootloc);CHKERRQ(ierr);
558dcca6d9dSJed Brown   ierr = PetscMalloc2(bas->niranks,&rootreqs,sf->nranks,&leafreqs);CHKERRQ(ierr);
55995fce210SBarry Smith   bas->ioffset[0] = 0;
56095fce210SBarry Smith   for (i=0; i<bas->niranks; i++) {
56195fce210SBarry Smith     bas->ioffset[i+1] = bas->ioffset[i] + ilengths[i];
56295fce210SBarry Smith     ierr = MPI_Irecv(bas->irootloc+bas->ioffset[i],ilengths[i],MPIU_INT,bas->iranks[i],bas->tag,comm,&rootreqs[i]);CHKERRQ(ierr);
56395fce210SBarry Smith   }
56495fce210SBarry Smith   for (i=0; i<sf->nranks; i++) {
56595fce210SBarry Smith     PetscMPIInt npoints;
56695fce210SBarry Smith     ierr = PetscMPIIntCast(sf->roffset[i+1] - sf->roffset[i],&npoints);CHKERRQ(ierr);
56795fce210SBarry Smith     ierr = MPI_Isend(sf->rremote+sf->roffset[i],npoints,MPIU_INT,sf->ranks[i],bas->tag,comm,&leafreqs[i]);CHKERRQ(ierr);
56895fce210SBarry Smith   }
569eb9baa12SBarry Smith   ierr = MPI_Waitall(bas->niranks,rootreqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
570eb9baa12SBarry Smith   ierr = MPI_Waitall(sf->nranks,leafreqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
57195fce210SBarry Smith   ierr = PetscFree(ilengths);CHKERRQ(ierr);
572eb9baa12SBarry Smith   ierr = PetscFree2(rootreqs,leafreqs);CHKERRQ(ierr);
57395fce210SBarry Smith   PetscFunctionReturn(0);
57495fce210SBarry Smith }
57595fce210SBarry Smith 
57695fce210SBarry Smith #undef __FUNCT__
57795fce210SBarry Smith #define __FUNCT__ "PetscSFBasicPackTypeSetup"
57895fce210SBarry Smith static PetscErrorCode PetscSFBasicPackTypeSetup(PetscSFBasicPack link,MPI_Datatype unit)
57995fce210SBarry Smith {
58095fce210SBarry Smith   PetscErrorCode ierr;
58195fce210SBarry Smith   PetscBool      isInt,isPetscInt,isPetscReal,is2Int,is2PetscInt;
58295fce210SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
58395fce210SBarry Smith   PetscBool isPetscComplex;
58495fce210SBarry Smith #endif
585*2a74e5f4SJed Brown   PetscInt       nPetscIntContig,nPetscRealContig,nPetscComplexContig;
58695fce210SBarry Smith 
58795fce210SBarry Smith   PetscFunctionBegin;
58895fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPI_INT,&isInt);CHKERRQ(ierr);
58995fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPIU_INT,&isPetscInt);CHKERRQ(ierr);
590*2a74e5f4SJed Brown   ierr = MPIPetsc_Type_compare_contig(unit,MPIU_INT,&nPetscIntContig);CHKERRQ(ierr);
59195fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPIU_REAL,&isPetscReal);CHKERRQ(ierr);
592*2a74e5f4SJed Brown   ierr = MPIPetsc_Type_compare_contig(unit,MPIU_REAL,&nPetscRealContig);CHKERRQ(ierr);
59395fce210SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
59495fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPIU_COMPLEX,&isPetscComplex);CHKERRQ(ierr);
595*2a74e5f4SJed Brown   ierr = MPIPetsc_Type_compare_contig(unit,MPIU_COMPLEX,&nPetscComplexContig);CHKERRQ(ierr);
59695fce210SBarry Smith #endif
59795fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPI_2INT,&is2Int);CHKERRQ(ierr);
59895fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPIU_2INT,&is2PetscInt);CHKERRQ(ierr);
599*2a74e5f4SJed Brown   link->bs = 1;
60053deab39SPeter Brune   if (isInt) {PackInit_int(link); PackInit_Logical_int(link); PackInit_Bitwise_int(link);}
60153deab39SPeter Brune   else if (isPetscInt) {PackInit_PetscInt(link); PackInit_Logical_PetscInt(link); PackInit_Bitwise_PetscInt(link);}
60253deab39SPeter Brune   else if (isPetscReal) {PackInit_PetscReal(link); PackInit_Logical_PetscReal(link);}
60395fce210SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
604*2a74e5f4SJed Brown   else if (isPetscComplex) PackInit_PetscComplex_1(link);
60595fce210SBarry Smith #endif
60695fce210SBarry Smith   else if (is2Int) PackInit_int_int(link);
60795fce210SBarry Smith   else if (is2PetscInt) PackInit_PetscInt_PetscInt(link);
608*2a74e5f4SJed Brown   else if (nPetscIntContig) {
609*2a74e5f4SJed Brown     if (nPetscIntContig%7 == 0) PackInit_PetscInt_7(link);
610*2a74e5f4SJed Brown     else if (nPetscIntContig%5 == 0) PackInit_PetscInt_5(link);
611*2a74e5f4SJed Brown     else if (nPetscIntContig%4 == 0) PackInit_PetscInt_4(link);
612*2a74e5f4SJed Brown     else if (nPetscIntContig%3 == 0) PackInit_PetscInt_3(link);
613*2a74e5f4SJed Brown     else if (nPetscIntContig%2 == 0) PackInit_PetscInt_2(link);
614*2a74e5f4SJed Brown     else PackInit_PetscInt(link);
615*2a74e5f4SJed Brown     link->bs = nPetscIntContig;
616*2a74e5f4SJed Brown     link->unitbytes *= nPetscIntContig;
617*2a74e5f4SJed Brown   } else if (nPetscRealContig) {
618*2a74e5f4SJed Brown     if (nPetscRealContig%7 == 0) PackInit_PetscReal_7(link);
619*2a74e5f4SJed Brown     else if (nPetscRealContig%5 == 0) PackInit_PetscReal_5(link);
620*2a74e5f4SJed Brown     else if (nPetscRealContig%4 == 0) PackInit_PetscReal_4(link);
621*2a74e5f4SJed Brown     else if (nPetscRealContig%3 == 0) PackInit_PetscReal_3(link);
622*2a74e5f4SJed Brown     else if (nPetscRealContig%2 == 0) PackInit_PetscReal_2(link);
623*2a74e5f4SJed Brown     else PackInit_PetscReal(link);
624*2a74e5f4SJed Brown     link->bs = nPetscRealContig;
625*2a74e5f4SJed Brown     link->unitbytes *= nPetscRealContig;
626*2a74e5f4SJed Brown #if defined(PETSC_HAVE_COMPLEX)
627*2a74e5f4SJed Brown   } else if (nPetscComplexContig) {
628*2a74e5f4SJed Brown     if (nPetscComplexContig%7 == 0) PackInit_PetscComplex_7(link);
629*2a74e5f4SJed Brown     else if (nPetscComplexContig%5 == 0) PackInit_PetscComplex_5(link);
630*2a74e5f4SJed Brown     else if (nPetscComplexContig%4 == 0) PackInit_PetscComplex_4(link);
631*2a74e5f4SJed Brown     else if (nPetscComplexContig%3 == 0) PackInit_PetscComplex_3(link);
632*2a74e5f4SJed Brown     else if (nPetscComplexContig%2 == 0) PackInit_PetscComplex_2(link);
633*2a74e5f4SJed Brown     else PackInit_PetscComplex_1(link);
634*2a74e5f4SJed Brown     link->bs = nPetscComplexContig;
635*2a74e5f4SJed Brown     link->unitbytes *= nPetscComplexContig;
636*2a74e5f4SJed Brown #endif
637*2a74e5f4SJed Brown   } else {
63895fce210SBarry Smith     PetscMPIInt bytes;
63995fce210SBarry Smith     ierr = MPI_Type_size(unit,&bytes);CHKERRQ(ierr);
64095fce210SBarry Smith     if (bytes % sizeof(int)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for type size not divisible by %D",sizeof(int));
64195fce210SBarry Smith     switch (bytes / sizeof(int)) {
64295fce210SBarry Smith     case 1: PackInit_block_int_1(link); break;
64395fce210SBarry Smith     case 2: PackInit_block_int_2(link); break;
64495fce210SBarry Smith     case 3: PackInit_block_int_3(link); break;
64595fce210SBarry Smith     case 4: PackInit_block_int_4(link); break;
64695fce210SBarry Smith     case 5: PackInit_block_int_5(link); break;
64795fce210SBarry Smith     case 6: PackInit_block_int_6(link); break;
64895fce210SBarry Smith     case 7: PackInit_block_int_7(link); break;
64995fce210SBarry Smith     case 8: PackInit_block_int_8(link); break;
65095fce210SBarry Smith     default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for arbitrary block sizes");
65195fce210SBarry Smith     }
65295fce210SBarry Smith   }
65395fce210SBarry Smith   ierr = MPI_Type_dup(unit,&link->unit);CHKERRQ(ierr);
65495fce210SBarry Smith   PetscFunctionReturn(0);
65595fce210SBarry Smith }
65695fce210SBarry Smith 
65795fce210SBarry Smith #undef __FUNCT__
65895fce210SBarry Smith #define __FUNCT__ "PetscSFBasicPackGetUnpackOp"
659*2a74e5f4SJed Brown static PetscErrorCode PetscSFBasicPackGetUnpackOp(PetscSF sf,PetscSFBasicPack link,MPI_Op op,void (**UnpackOp)(PetscInt,PetscInt,const PetscInt*,void*,const void*))
66095fce210SBarry Smith {
66195fce210SBarry Smith   PetscFunctionBegin;
66295fce210SBarry Smith   *UnpackOp = NULL;
6638bfbc91cSJed Brown   if (op == MPIU_REPLACE) *UnpackOp = link->UnpackInsert;
66495fce210SBarry Smith   else if (op == MPI_SUM || op == MPIU_SUM) *UnpackOp = link->UnpackAdd;
66553deab39SPeter Brune   else if (op == MPI_PROD) *UnpackOp = link->UnpackMult;
66695fce210SBarry Smith   else if (op == MPI_MAX || op == MPIU_MAX) *UnpackOp = link->UnpackMax;
66795fce210SBarry Smith   else if (op == MPI_MIN || op == MPIU_MIN) *UnpackOp = link->UnpackMin;
66853deab39SPeter Brune   else if (op == MPI_LAND) *UnpackOp = link->UnpackLAND;
66953deab39SPeter Brune   else if (op == MPI_BAND) *UnpackOp = link->UnpackBAND;
67053deab39SPeter Brune   else if (op == MPI_LOR) *UnpackOp = link->UnpackLOR;
67153deab39SPeter Brune   else if (op == MPI_BOR) *UnpackOp = link->UnpackBOR;
67253deab39SPeter Brune   else if (op == MPI_LXOR) *UnpackOp = link->UnpackLXOR;
67353deab39SPeter Brune   else if (op == MPI_BXOR) *UnpackOp = link->UnpackBXOR;
67495fce210SBarry Smith   else if (op == MPI_MAXLOC) *UnpackOp = link->UnpackMaxloc;
67595fce210SBarry Smith   else if (op == MPI_MINLOC) *UnpackOp = link->UnpackMinloc;
67695fce210SBarry Smith   else SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_SUP,"No support for MPI_Op");
67795fce210SBarry Smith   PetscFunctionReturn(0);
67895fce210SBarry Smith }
67995fce210SBarry Smith #undef __FUNCT__
68095fce210SBarry Smith #define __FUNCT__ "PetscSFBasicPackGetFetchAndOp"
681*2a74e5f4SJed Brown static PetscErrorCode PetscSFBasicPackGetFetchAndOp(PetscSF sf,PetscSFBasicPack link,MPI_Op op,void (**FetchAndOp)(PetscInt,PetscInt,const PetscInt*,void*,void*))
68295fce210SBarry Smith {
68395fce210SBarry Smith   PetscFunctionBegin;
68495fce210SBarry Smith   *FetchAndOp = NULL;
6858bfbc91cSJed Brown   if (op == MPIU_REPLACE) *FetchAndOp = link->FetchAndInsert;
68695fce210SBarry Smith   else if (op == MPI_SUM || op == MPIU_SUM) *FetchAndOp = link->FetchAndAdd;
68795fce210SBarry Smith   else if (op == MPI_MAX || op == MPIU_MAX) *FetchAndOp = link->FetchAndMax;
68895fce210SBarry Smith   else if (op == MPI_MIN || op == MPIU_MIN) *FetchAndOp = link->FetchAndMin;
68995fce210SBarry Smith   else if (op == MPI_MAXLOC) *FetchAndOp = link->FetchAndMaxloc;
69095fce210SBarry Smith   else if (op == MPI_MINLOC) *FetchAndOp = link->FetchAndMinloc;
69153deab39SPeter Brune   else if (op == MPI_PROD)   *FetchAndOp = link->FetchAndMult;
69253deab39SPeter Brune   else if (op == MPI_LAND)   *FetchAndOp = link->FetchAndLAND;
69353deab39SPeter Brune   else if (op == MPI_BAND)   *FetchAndOp = link->FetchAndBAND;
69453deab39SPeter Brune   else if (op == MPI_LOR)    *FetchAndOp = link->FetchAndLOR;
69553deab39SPeter Brune   else if (op == MPI_BOR)    *FetchAndOp = link->FetchAndBOR;
69653deab39SPeter Brune   else if (op == MPI_LXOR)   *FetchAndOp = link->FetchAndLXOR;
69753deab39SPeter Brune   else if (op == MPI_BXOR)   *FetchAndOp = link->FetchAndBXOR;
69895fce210SBarry Smith   else SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_SUP,"No support for MPI_Op");
69995fce210SBarry Smith   PetscFunctionReturn(0);
70095fce210SBarry Smith }
70195fce210SBarry Smith 
70295fce210SBarry Smith #undef __FUNCT__
70395fce210SBarry Smith #define __FUNCT__ "PetscSFBasicPackGetReqs"
70495fce210SBarry Smith static PetscErrorCode PetscSFBasicPackGetReqs(PetscSF sf,PetscSFBasicPack link,MPI_Request **rootreqs,MPI_Request **leafreqs)
70595fce210SBarry Smith {
70695fce210SBarry Smith   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
70795fce210SBarry Smith 
70895fce210SBarry Smith   PetscFunctionBegin;
70995fce210SBarry Smith   if (rootreqs) *rootreqs = link->requests;
71095fce210SBarry Smith   if (leafreqs) *leafreqs = link->requests + bas->niranks;
71195fce210SBarry Smith   PetscFunctionReturn(0);
71295fce210SBarry Smith }
71395fce210SBarry Smith 
71495fce210SBarry Smith #undef __FUNCT__
71595fce210SBarry Smith #define __FUNCT__ "PetscSFBasicPackWaitall"
71695fce210SBarry Smith static PetscErrorCode PetscSFBasicPackWaitall(PetscSF sf,PetscSFBasicPack link)
71795fce210SBarry Smith {
71895fce210SBarry Smith   PetscSF_Basic  *bas = (PetscSF_Basic*)sf->data;
71995fce210SBarry Smith   PetscErrorCode ierr;
72095fce210SBarry Smith 
72195fce210SBarry Smith   PetscFunctionBegin;
72295fce210SBarry Smith   ierr = MPI_Waitall(bas->niranks+sf->nranks,link->requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
72395fce210SBarry Smith   PetscFunctionReturn(0);
72495fce210SBarry Smith }
72595fce210SBarry Smith 
72695fce210SBarry Smith #undef __FUNCT__
72795fce210SBarry Smith #define __FUNCT__ "PetscSFBasicGetRootInfo"
72895fce210SBarry Smith static PetscErrorCode PetscSFBasicGetRootInfo(PetscSF sf,PetscInt *nrootranks,const PetscMPIInt **rootranks,const PetscInt **rootoffset,const PetscInt **rootloc)
72995fce210SBarry Smith {
73095fce210SBarry Smith   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
73195fce210SBarry Smith 
73295fce210SBarry Smith   PetscFunctionBegin;
73395fce210SBarry Smith   if (nrootranks) *nrootranks = bas->niranks;
73495fce210SBarry Smith   if (rootranks)  *rootranks  = bas->iranks;
73595fce210SBarry Smith   if (rootoffset) *rootoffset = bas->ioffset;
73695fce210SBarry Smith   if (rootloc)    *rootloc    = bas->irootloc;
73795fce210SBarry Smith   PetscFunctionReturn(0);
73895fce210SBarry Smith }
73995fce210SBarry Smith 
74095fce210SBarry Smith #undef __FUNCT__
74195fce210SBarry Smith #define __FUNCT__ "PetscSFBasicGetLeafInfo"
74295fce210SBarry Smith static PetscErrorCode PetscSFBasicGetLeafInfo(PetscSF sf,PetscInt *nleafranks,const PetscMPIInt **leafranks,const PetscInt **leafoffset,const PetscInt **leafloc)
74395fce210SBarry Smith {
74495fce210SBarry Smith   PetscFunctionBegin;
74595fce210SBarry Smith   if (nleafranks) *nleafranks = sf->nranks;
74695fce210SBarry Smith   if (leafranks)  *leafranks  = sf->ranks;
74795fce210SBarry Smith   if (leafoffset) *leafoffset = sf->roffset;
74895fce210SBarry Smith   if (leafloc)    *leafloc    = sf->rmine;
74995fce210SBarry Smith   PetscFunctionReturn(0);
75095fce210SBarry Smith }
75195fce210SBarry Smith 
75295fce210SBarry Smith #undef __FUNCT__
75395fce210SBarry Smith #define __FUNCT__ "PetscSFBasicGetPack"
75495fce210SBarry Smith static PetscErrorCode PetscSFBasicGetPack(PetscSF sf,MPI_Datatype unit,const void *key,PetscSFBasicPack *mylink)
75595fce210SBarry Smith {
75695fce210SBarry Smith   PetscSF_Basic    *bas = (PetscSF_Basic*)sf->data;
75795fce210SBarry Smith   PetscErrorCode   ierr;
75895fce210SBarry Smith   PetscSFBasicPack link,*p;
75995fce210SBarry Smith   PetscInt         nrootranks,nleafranks;
76095fce210SBarry Smith   const PetscInt   *rootoffset,*leafoffset;
76195fce210SBarry Smith 
76295fce210SBarry Smith   PetscFunctionBegin;
76395fce210SBarry Smith   /* Look for types in cache */
76495fce210SBarry Smith   for (p=&bas->avail; (link=*p); p=&link->next) {
76595fce210SBarry Smith     PetscBool match;
76695fce210SBarry Smith     ierr = MPIPetsc_Type_compare(unit,link->unit,&match);CHKERRQ(ierr);
76795fce210SBarry Smith     if (match) {
76895fce210SBarry Smith       *p = link->next;          /* Remove from available list */
76995fce210SBarry Smith       goto found;
77095fce210SBarry Smith     }
77195fce210SBarry Smith   }
77295fce210SBarry Smith 
77395fce210SBarry Smith   /* Create new composite types for each send rank */
77495fce210SBarry Smith   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,NULL,&rootoffset,NULL);CHKERRQ(ierr);
77595fce210SBarry Smith   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,NULL,&leafoffset,NULL);CHKERRQ(ierr);
776b00a9115SJed Brown   ierr = PetscNew(&link);CHKERRQ(ierr);
77795fce210SBarry Smith   ierr = PetscSFBasicPackTypeSetup(link,unit);CHKERRQ(ierr);
778dcca6d9dSJed Brown   ierr = PetscMalloc2(rootoffset[nrootranks]*link->unitbytes,&link->root,leafoffset[nleafranks]*link->unitbytes,&link->leaf);CHKERRQ(ierr);
779785e854fSJed Brown   ierr = PetscMalloc1((nrootranks+nleafranks),&link->requests);CHKERRQ(ierr);
78095fce210SBarry Smith 
78195fce210SBarry Smith found:
78295fce210SBarry Smith   link->key  = key;
78395fce210SBarry Smith   link->next = bas->inuse;
78495fce210SBarry Smith   bas->inuse = link;
78595fce210SBarry Smith 
78695fce210SBarry Smith   *mylink = link;
78795fce210SBarry Smith   PetscFunctionReturn(0);
78895fce210SBarry Smith }
78995fce210SBarry Smith 
79095fce210SBarry Smith #undef __FUNCT__
79195fce210SBarry Smith #define __FUNCT__ "PetscSFBasicGetPackInUse"
79295fce210SBarry Smith static PetscErrorCode PetscSFBasicGetPackInUse(PetscSF sf,MPI_Datatype unit,const void *key,PetscCopyMode cmode,PetscSFBasicPack *mylink)
79395fce210SBarry Smith {
79495fce210SBarry Smith   PetscSF_Basic    *bas = (PetscSF_Basic*)sf->data;
79595fce210SBarry Smith   PetscErrorCode   ierr;
79695fce210SBarry Smith   PetscSFBasicPack link,*p;
79795fce210SBarry Smith 
79895fce210SBarry Smith   PetscFunctionBegin;
79995fce210SBarry Smith   /* Look for types in cache */
80095fce210SBarry Smith   for (p=&bas->inuse; (link=*p); p=&link->next) {
80195fce210SBarry Smith     PetscBool match;
80295fce210SBarry Smith     ierr = MPIPetsc_Type_compare(unit,link->unit,&match);CHKERRQ(ierr);
80395fce210SBarry Smith     if (match) {
80495fce210SBarry Smith       switch (cmode) {
80595fce210SBarry Smith       case PETSC_OWN_POINTER: *p = link->next; break; /* Remove from inuse list */
80695fce210SBarry Smith       case PETSC_USE_POINTER: break;
80795fce210SBarry Smith       default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"invalid cmode");
80895fce210SBarry Smith       }
80995fce210SBarry Smith       *mylink = link;
81095fce210SBarry Smith       PetscFunctionReturn(0);
81195fce210SBarry Smith     }
81295fce210SBarry Smith   }
81395fce210SBarry Smith   SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Could not find pack");
81495fce210SBarry Smith   PetscFunctionReturn(0);
81595fce210SBarry Smith }
81695fce210SBarry Smith 
81795fce210SBarry Smith #undef __FUNCT__
81895fce210SBarry Smith #define __FUNCT__ "PetscSFBasicReclaimPack"
81995fce210SBarry Smith static PetscErrorCode PetscSFBasicReclaimPack(PetscSF sf,PetscSFBasicPack *link)
82095fce210SBarry Smith {
82195fce210SBarry Smith   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
82295fce210SBarry Smith 
82395fce210SBarry Smith   PetscFunctionBegin;
82495fce210SBarry Smith   (*link)->key  = NULL;
82595fce210SBarry Smith   (*link)->next = bas->avail;
82695fce210SBarry Smith   bas->avail    = *link;
82795fce210SBarry Smith   *link         = NULL;
82895fce210SBarry Smith   PetscFunctionReturn(0);
82995fce210SBarry Smith }
83095fce210SBarry Smith 
83195fce210SBarry Smith #undef __FUNCT__
83295fce210SBarry Smith #define __FUNCT__ "PetscSFSetFromOptions_Basic"
83395fce210SBarry Smith static PetscErrorCode PetscSFSetFromOptions_Basic(PetscSF sf)
83495fce210SBarry Smith {
83595fce210SBarry Smith   PetscErrorCode ierr;
83695fce210SBarry Smith 
83795fce210SBarry Smith   PetscFunctionBegin;
83895fce210SBarry Smith   ierr = PetscOptionsHead("PetscSF Basic options");CHKERRQ(ierr);
83995fce210SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
84095fce210SBarry Smith   PetscFunctionReturn(0);
84195fce210SBarry Smith }
84295fce210SBarry Smith 
84395fce210SBarry Smith #undef __FUNCT__
84495fce210SBarry Smith #define __FUNCT__ "PetscSFReset_Basic"
84595fce210SBarry Smith static PetscErrorCode PetscSFReset_Basic(PetscSF sf)
84695fce210SBarry Smith {
84795fce210SBarry Smith   PetscSF_Basic    *bas = (PetscSF_Basic*)sf->data;
84895fce210SBarry Smith   PetscErrorCode   ierr;
84995fce210SBarry Smith   PetscSFBasicPack link,next;
85095fce210SBarry Smith 
85195fce210SBarry Smith   PetscFunctionBegin;
85295fce210SBarry Smith   ierr = PetscFree(bas->iranks);CHKERRQ(ierr);
85395fce210SBarry Smith   ierr = PetscFree2(bas->ioffset,bas->irootloc);CHKERRQ(ierr);
85495fce210SBarry Smith   if (bas->inuse) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Outstanding operation has not been completed");
85595fce210SBarry Smith   for (link=bas->avail; link; link=next) {
85695fce210SBarry Smith     next = link->next;
857bb86044dSJed Brown #if defined(PETSC_HAVE_MPI_TYPE_DUP)
85895fce210SBarry Smith     ierr = MPI_Type_free(&link->unit);CHKERRQ(ierr);
859bb86044dSJed Brown #endif
86095fce210SBarry Smith     ierr = PetscFree2(link->root,link->leaf);CHKERRQ(ierr);
86195fce210SBarry Smith     ierr = PetscFree(link->requests);CHKERRQ(ierr);
86295fce210SBarry Smith     ierr = PetscFree(link);CHKERRQ(ierr);
86395fce210SBarry Smith   }
86495fce210SBarry Smith   bas->avail = NULL;
86595fce210SBarry Smith   PetscFunctionReturn(0);
86695fce210SBarry Smith }
86795fce210SBarry Smith 
86895fce210SBarry Smith #undef __FUNCT__
86995fce210SBarry Smith #define __FUNCT__ "PetscSFDestroy_Basic"
87095fce210SBarry Smith static PetscErrorCode PetscSFDestroy_Basic(PetscSF sf)
87195fce210SBarry Smith {
87295fce210SBarry Smith   PetscErrorCode ierr;
87395fce210SBarry Smith 
87495fce210SBarry Smith   PetscFunctionBegin;
87595fce210SBarry Smith   ierr = PetscSFReset_Basic(sf);CHKERRQ(ierr);
87695fce210SBarry Smith   ierr = PetscFree(sf->data);CHKERRQ(ierr);
87795fce210SBarry Smith   PetscFunctionReturn(0);
87895fce210SBarry Smith }
87995fce210SBarry Smith 
88095fce210SBarry Smith #undef __FUNCT__
88195fce210SBarry Smith #define __FUNCT__ "PetscSFView_Basic"
88295fce210SBarry Smith static PetscErrorCode PetscSFView_Basic(PetscSF sf,PetscViewer viewer)
88395fce210SBarry Smith {
88495fce210SBarry Smith   /* PetscSF_Basic *bas = (PetscSF_Basic*)sf->data; */
88595fce210SBarry Smith   PetscErrorCode ierr;
88695fce210SBarry Smith   PetscBool      iascii;
88795fce210SBarry Smith 
88895fce210SBarry Smith   PetscFunctionBegin;
88995fce210SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
89095fce210SBarry Smith   if (iascii) {
89195fce210SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  sort=%s\n",sf->rankorder ? "rank-order" : "unordered");CHKERRQ(ierr);
89295fce210SBarry Smith   }
89395fce210SBarry Smith   PetscFunctionReturn(0);
89495fce210SBarry Smith }
89595fce210SBarry Smith 
89695fce210SBarry Smith #undef __FUNCT__
89795fce210SBarry Smith #define __FUNCT__ "PetscSFBcastBegin_Basic"
89895fce210SBarry Smith /* Send from roots to leaves */
89995fce210SBarry Smith static PetscErrorCode PetscSFBcastBegin_Basic(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
90095fce210SBarry Smith {
90195fce210SBarry Smith   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;
90295fce210SBarry Smith   PetscErrorCode    ierr;
90395fce210SBarry Smith   PetscSFBasicPack  link;
90495fce210SBarry Smith   PetscInt          i,nrootranks,nleafranks;
90595fce210SBarry Smith   const PetscInt    *rootoffset,*leafoffset,*rootloc,*leafloc;
90695fce210SBarry Smith   const PetscMPIInt *rootranks,*leafranks;
90795fce210SBarry Smith   MPI_Request       *rootreqs,*leafreqs;
90895fce210SBarry Smith   size_t            unitbytes;
90995fce210SBarry Smith 
91095fce210SBarry Smith   PetscFunctionBegin;
91195fce210SBarry Smith   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,&rootranks,&rootoffset,&rootloc);CHKERRQ(ierr);
91295fce210SBarry Smith   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,&leafranks,&leafoffset,&leafloc);CHKERRQ(ierr);
91395fce210SBarry Smith   ierr = PetscSFBasicGetPack(sf,unit,rootdata,&link);CHKERRQ(ierr);
91495fce210SBarry Smith 
91595fce210SBarry Smith   unitbytes = link->unitbytes;
91695fce210SBarry Smith 
91795fce210SBarry Smith   ierr = PetscSFBasicPackGetReqs(sf,link,&rootreqs,&leafreqs);CHKERRQ(ierr);
91895fce210SBarry Smith   /* Eagerly post leaf receives */
91995fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
92095fce210SBarry Smith     PetscMPIInt n = leafoffset[i+1] - leafoffset[i];
92195fce210SBarry Smith     ierr = MPI_Irecv(link->leaf+leafoffset[i]*unitbytes,n,unit,leafranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&leafreqs[i]);CHKERRQ(ierr);
92295fce210SBarry Smith   }
92395fce210SBarry Smith   /* Pack and send root data */
92495fce210SBarry Smith   for (i=0; i<nrootranks; i++) {
92595fce210SBarry Smith     PetscMPIInt n          = rootoffset[i+1] - rootoffset[i];
92695fce210SBarry Smith     void        *packstart = link->root+rootoffset[i]*unitbytes;
927*2a74e5f4SJed Brown     (*link->Pack)(n,link->bs,rootloc+rootoffset[i],rootdata,packstart);
92895fce210SBarry Smith     ierr = MPI_Isend(packstart,n,unit,rootranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&rootreqs[i]);CHKERRQ(ierr);
92995fce210SBarry Smith   }
93095fce210SBarry Smith   PetscFunctionReturn(0);
93195fce210SBarry Smith }
93295fce210SBarry Smith 
93395fce210SBarry Smith #undef __FUNCT__
93495fce210SBarry Smith #define __FUNCT__ "PetscSFBcastEnd_Basic"
93595fce210SBarry Smith PetscErrorCode PetscSFBcastEnd_Basic(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
93695fce210SBarry Smith {
93795fce210SBarry Smith   PetscErrorCode   ierr;
93895fce210SBarry Smith   PetscSFBasicPack link;
93995fce210SBarry Smith   PetscInt         i,nleafranks;
94095fce210SBarry Smith   const PetscInt   *leafoffset,*leafloc;
94195fce210SBarry Smith 
94295fce210SBarry Smith   PetscFunctionBegin;
94395fce210SBarry Smith   ierr = PetscSFBasicGetPackInUse(sf,unit,rootdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr);
94495fce210SBarry Smith   ierr = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
94595fce210SBarry Smith   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,NULL,&leafoffset,&leafloc);CHKERRQ(ierr);
94695fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
94795fce210SBarry Smith     PetscMPIInt n          = leafoffset[i+1] - leafoffset[i];
94895fce210SBarry Smith     const void  *packstart = link->leaf+leafoffset[i]*link->unitbytes;
949*2a74e5f4SJed Brown     (*link->UnpackInsert)(n,link->bs,leafloc+leafoffset[i],leafdata,packstart);
95095fce210SBarry Smith   }
95195fce210SBarry Smith   ierr = PetscSFBasicReclaimPack(sf,&link);CHKERRQ(ierr);
95295fce210SBarry Smith   PetscFunctionReturn(0);
95395fce210SBarry Smith }
95495fce210SBarry Smith 
95595fce210SBarry Smith #undef __FUNCT__
95695fce210SBarry Smith #define __FUNCT__ "PetscSFReduceBegin_Basic"
95795fce210SBarry Smith /* leaf -> root with reduction */
95895fce210SBarry Smith PetscErrorCode PetscSFReduceBegin_Basic(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
95995fce210SBarry Smith {
96095fce210SBarry Smith   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;
96195fce210SBarry Smith   PetscSFBasicPack  link;
96295fce210SBarry Smith   PetscErrorCode    ierr;
96395fce210SBarry Smith   PetscInt          i,nrootranks,nleafranks;
96495fce210SBarry Smith   const PetscInt    *rootoffset,*leafoffset,*rootloc,*leafloc;
96595fce210SBarry Smith   const PetscMPIInt *rootranks,*leafranks;
96695fce210SBarry Smith   MPI_Request       *rootreqs,*leafreqs;
96795fce210SBarry Smith   size_t            unitbytes;
96895fce210SBarry Smith 
96995fce210SBarry Smith   PetscFunctionBegin;
97095fce210SBarry Smith   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,&rootranks,&rootoffset,&rootloc);CHKERRQ(ierr);
97195fce210SBarry Smith   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,&leafranks,&leafoffset,&leafloc);CHKERRQ(ierr);
97295fce210SBarry Smith   ierr = PetscSFBasicGetPack(sf,unit,rootdata,&link);CHKERRQ(ierr);
97395fce210SBarry Smith 
97495fce210SBarry Smith   unitbytes = link->unitbytes;
97595fce210SBarry Smith 
97695fce210SBarry Smith   ierr = PetscSFBasicPackGetReqs(sf,link,&rootreqs,&leafreqs);CHKERRQ(ierr);
97795fce210SBarry Smith   /* Eagerly post root receives */
97895fce210SBarry Smith   for (i=0; i<nrootranks; i++) {
97995fce210SBarry Smith     PetscMPIInt n = rootoffset[i+1] - rootoffset[i];
98095fce210SBarry Smith     ierr = MPI_Irecv(link->root+rootoffset[i]*unitbytes,n,unit,rootranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&rootreqs[i]);CHKERRQ(ierr);
98195fce210SBarry Smith   }
98295fce210SBarry Smith   /* Pack and send leaf data */
98395fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
98495fce210SBarry Smith     PetscMPIInt n          = leafoffset[i+1] - leafoffset[i];
98595fce210SBarry Smith     void        *packstart = link->leaf+leafoffset[i]*unitbytes;
986*2a74e5f4SJed Brown     (*link->Pack)(n,link->bs,leafloc+leafoffset[i],leafdata,packstart);
98795fce210SBarry Smith     ierr = MPI_Isend(packstart,n,unit,leafranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&leafreqs[i]);CHKERRQ(ierr);
98895fce210SBarry Smith   }
98995fce210SBarry Smith   PetscFunctionReturn(0);
99095fce210SBarry Smith }
99195fce210SBarry Smith 
99295fce210SBarry Smith #undef __FUNCT__
99395fce210SBarry Smith #define __FUNCT__ "PetscSFReduceEnd_Basic"
99495fce210SBarry Smith static PetscErrorCode PetscSFReduceEnd_Basic(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
99595fce210SBarry Smith {
996*2a74e5f4SJed Brown   void             (*UnpackOp)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
99795fce210SBarry Smith   PetscErrorCode   ierr;
99895fce210SBarry Smith   PetscSFBasicPack link;
99995fce210SBarry Smith   PetscInt         i,nrootranks;
100095fce210SBarry Smith   const PetscInt   *rootoffset,*rootloc;
100195fce210SBarry Smith 
100295fce210SBarry Smith   PetscFunctionBegin;
100395fce210SBarry Smith   ierr = PetscSFBasicGetPackInUse(sf,unit,rootdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr);
100495fce210SBarry Smith   /* This implementation could be changed to unpack as receives arrive, at the cost of non-determinism */
100595fce210SBarry Smith   ierr = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
100695fce210SBarry Smith   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,NULL,&rootoffset,&rootloc);CHKERRQ(ierr);
100795fce210SBarry Smith   ierr = PetscSFBasicPackGetUnpackOp(sf,link,op,&UnpackOp);CHKERRQ(ierr);
100895fce210SBarry Smith   for (i=0; i<nrootranks; i++) {
100995fce210SBarry Smith     PetscMPIInt n          = rootoffset[i+1] - rootoffset[i];
101095fce210SBarry Smith     const void  *packstart = link->root+rootoffset[i]*link->unitbytes;
101195fce210SBarry Smith 
1012*2a74e5f4SJed Brown     (*UnpackOp)(n,link->bs,rootloc+rootoffset[i],rootdata,packstart);
101395fce210SBarry Smith   }
101495fce210SBarry Smith   ierr = PetscSFBasicReclaimPack(sf,&link);CHKERRQ(ierr);
101595fce210SBarry Smith   PetscFunctionReturn(0);
101695fce210SBarry Smith }
101795fce210SBarry Smith 
101895fce210SBarry Smith #undef __FUNCT__
101995fce210SBarry Smith #define __FUNCT__ "PetscSFFetchAndOpBegin_Basic"
102095fce210SBarry Smith static PetscErrorCode PetscSFFetchAndOpBegin_Basic(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
102195fce210SBarry Smith {
102295fce210SBarry Smith   PetscErrorCode ierr;
102395fce210SBarry Smith 
102495fce210SBarry Smith   PetscFunctionBegin;
102595fce210SBarry Smith   ierr = PetscSFReduceBegin_Basic(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
102695fce210SBarry Smith   PetscFunctionReturn(0);
102795fce210SBarry Smith }
102895fce210SBarry Smith 
102995fce210SBarry Smith #undef __FUNCT__
103095fce210SBarry Smith #define __FUNCT__ "PetscSFFetchAndOpEnd_Basic"
103195fce210SBarry Smith static PetscErrorCode PetscSFFetchAndOpEnd_Basic(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
103295fce210SBarry Smith {
103395fce210SBarry Smith   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;
1034*2a74e5f4SJed Brown   void              (*FetchAndOp)(PetscInt,PetscInt,const PetscInt*,void*,void*);
103595fce210SBarry Smith   PetscErrorCode    ierr;
103695fce210SBarry Smith   PetscSFBasicPack  link;
103795fce210SBarry Smith   PetscInt          i,nrootranks,nleafranks;
103895fce210SBarry Smith   const PetscInt    *rootoffset,*leafoffset,*rootloc,*leafloc;
103995fce210SBarry Smith   const PetscMPIInt *rootranks,*leafranks;
104095fce210SBarry Smith   MPI_Request       *rootreqs,*leafreqs;
104195fce210SBarry Smith   size_t            unitbytes;
104295fce210SBarry Smith 
104395fce210SBarry Smith   PetscFunctionBegin;
104495fce210SBarry Smith   ierr = PetscSFBasicGetPackInUse(sf,unit,rootdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr);
104595fce210SBarry Smith   /* This implementation could be changed to unpack as receives arrive, at the cost of non-determinism */
104695fce210SBarry Smith   ierr      = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
104795fce210SBarry Smith   unitbytes = link->unitbytes;
104895fce210SBarry Smith   ierr      = PetscSFBasicGetRootInfo(sf,&nrootranks,&rootranks,&rootoffset,&rootloc);CHKERRQ(ierr);
104995fce210SBarry Smith   ierr      = PetscSFBasicGetLeafInfo(sf,&nleafranks,&leafranks,&leafoffset,&leafloc);CHKERRQ(ierr);
105095fce210SBarry Smith   ierr      = PetscSFBasicPackGetReqs(sf,link,&rootreqs,&leafreqs);CHKERRQ(ierr);
105195fce210SBarry Smith   /* Post leaf receives */
105295fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
105395fce210SBarry Smith     PetscMPIInt n = leafoffset[i+1] - leafoffset[i];
105495fce210SBarry Smith     ierr = MPI_Irecv(link->leaf+leafoffset[i]*unitbytes,n,unit,leafranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&leafreqs[i]);CHKERRQ(ierr);
105595fce210SBarry Smith   }
105695fce210SBarry Smith   /* Process local fetch-and-op, post root sends */
105795fce210SBarry Smith   ierr = PetscSFBasicPackGetFetchAndOp(sf,link,op,&FetchAndOp);CHKERRQ(ierr);
105895fce210SBarry Smith   for (i=0; i<nrootranks; i++) {
105995fce210SBarry Smith     PetscMPIInt n          = rootoffset[i+1] - rootoffset[i];
106095fce210SBarry Smith     void        *packstart = link->root+rootoffset[i]*unitbytes;
106195fce210SBarry Smith 
1062*2a74e5f4SJed Brown     (*FetchAndOp)(n,link->bs,rootloc+rootoffset[i],rootdata,packstart);
106395fce210SBarry Smith     ierr = MPI_Isend(packstart,n,unit,rootranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&rootreqs[i]);CHKERRQ(ierr);
106495fce210SBarry Smith   }
106595fce210SBarry Smith   ierr = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
106695fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
106795fce210SBarry Smith     PetscMPIInt n          = leafoffset[i+1] - leafoffset[i];
106895fce210SBarry Smith     const void  *packstart = link->leaf+leafoffset[i]*unitbytes;
1069*2a74e5f4SJed Brown     (*link->UnpackInsert)(n,link->bs,leafloc+leafoffset[i],leafupdate,packstart);
107095fce210SBarry Smith   }
107195fce210SBarry Smith   ierr = PetscSFBasicReclaimPack(sf,&link);CHKERRQ(ierr);
107295fce210SBarry Smith   PetscFunctionReturn(0);
107395fce210SBarry Smith }
107495fce210SBarry Smith 
107595fce210SBarry Smith #undef __FUNCT__
107695fce210SBarry Smith #define __FUNCT__ "PetscSFCreate_Basic"
10778cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscSFCreate_Basic(PetscSF sf)
107895fce210SBarry Smith {
107995fce210SBarry Smith   PetscSF_Basic  *bas = (PetscSF_Basic*)sf->data;
108095fce210SBarry Smith   PetscErrorCode ierr;
108195fce210SBarry Smith 
108295fce210SBarry Smith   PetscFunctionBegin;
108395fce210SBarry Smith   sf->ops->SetUp           = PetscSFSetUp_Basic;
108495fce210SBarry Smith   sf->ops->SetFromOptions  = PetscSFSetFromOptions_Basic;
108595fce210SBarry Smith   sf->ops->Reset           = PetscSFReset_Basic;
108695fce210SBarry Smith   sf->ops->Destroy         = PetscSFDestroy_Basic;
108795fce210SBarry Smith   sf->ops->View            = PetscSFView_Basic;
108895fce210SBarry Smith   sf->ops->BcastBegin      = PetscSFBcastBegin_Basic;
108995fce210SBarry Smith   sf->ops->BcastEnd        = PetscSFBcastEnd_Basic;
109095fce210SBarry Smith   sf->ops->ReduceBegin     = PetscSFReduceBegin_Basic;
109195fce210SBarry Smith   sf->ops->ReduceEnd       = PetscSFReduceEnd_Basic;
109295fce210SBarry Smith   sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Basic;
109395fce210SBarry Smith   sf->ops->FetchAndOpEnd   = PetscSFFetchAndOpEnd_Basic;
109495fce210SBarry Smith 
1095b00a9115SJed Brown   ierr     = PetscNewLog(sf,&bas);CHKERRQ(ierr);
109695fce210SBarry Smith   sf->data = (void*)bas;
109795fce210SBarry Smith   PetscFunctionReturn(0);
109895fce210SBarry Smith }
1099