xref: /petsc/src/vec/is/sf/impls/basic/sfbasic.c (revision 77b7e48d5dad8fd77343dff3c4c7fc9e92f9c671)
18cd53115SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/
395fce210SBarry Smith 
495fce210SBarry Smith typedef struct _n_PetscSFBasicPack *PetscSFBasicPack;
595fce210SBarry Smith struct _n_PetscSFBasicPack {
62a74e5f4SJed Brown   void (*Pack)(PetscInt,PetscInt,const PetscInt*,const void*,void*);
72a74e5f4SJed Brown   void (*UnpackInsert)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
82a74e5f4SJed Brown   void (*UnpackAdd)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
92a74e5f4SJed Brown   void (*UnpackMin)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
102a74e5f4SJed Brown   void (*UnpackMax)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
112a74e5f4SJed Brown   void (*UnpackMinloc)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
122a74e5f4SJed Brown   void (*UnpackMaxloc)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
132a74e5f4SJed Brown   void (*UnpackMult)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
142a74e5f4SJed Brown   void (*UnpackLAND)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
152a74e5f4SJed Brown   void (*UnpackBAND)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
162a74e5f4SJed Brown   void (*UnpackLOR)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
172a74e5f4SJed Brown   void (*UnpackBOR)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
182a74e5f4SJed Brown   void (*UnpackLXOR)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
192a74e5f4SJed Brown   void (*UnpackBXOR)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
202a74e5f4SJed Brown   void (*FetchAndInsert)(PetscInt,PetscInt,const PetscInt*,void*,void*);
212a74e5f4SJed Brown   void (*FetchAndAdd)(PetscInt,PetscInt,const PetscInt*,void*,void*);
222a74e5f4SJed Brown   void (*FetchAndMin)(PetscInt,PetscInt,const PetscInt*,void*,void*);
232a74e5f4SJed Brown   void (*FetchAndMax)(PetscInt,PetscInt,const PetscInt*,void*,void*);
242a74e5f4SJed Brown   void (*FetchAndMinloc)(PetscInt,PetscInt,const PetscInt*,void*,void*);
252a74e5f4SJed Brown   void (*FetchAndMaxloc)(PetscInt,PetscInt,const PetscInt*,void*,void*);
262a74e5f4SJed Brown   void (*FetchAndMult)(PetscInt,PetscInt,const PetscInt*,void*,void*);
272a74e5f4SJed Brown   void (*FetchAndLAND)(PetscInt,PetscInt,const PetscInt*,void*,void*);
282a74e5f4SJed Brown   void (*FetchAndBAND)(PetscInt,PetscInt,const PetscInt*,void*,void*);
292a74e5f4SJed Brown   void (*FetchAndLOR)(PetscInt,PetscInt,const PetscInt*,void*,void*);
302a74e5f4SJed Brown   void (*FetchAndBOR)(PetscInt,PetscInt,const PetscInt*,void*,void*);
312a74e5f4SJed Brown   void (*FetchAndLXOR)(PetscInt,PetscInt,const PetscInt*,void*,void*);
322a74e5f4SJed Brown   void (*FetchAndBXOR)(PetscInt,PetscInt,const PetscInt*,void*,void*);
3395fce210SBarry Smith 
3495fce210SBarry Smith   MPI_Datatype     unit;
35*77b7e48dSJunchao Zhang   PetscBool        isbuiltin;   /* Is unit an MPI builtin datatype? */
3695fce210SBarry Smith   size_t           unitbytes;   /* Number of bytes in a unit */
372a74e5f4SJed Brown   PetscInt         bs;          /* Number of basic units in a unit */
3895fce210SBarry Smith   const void       *key;        /* Array used as key for operation */
39bf39f1bfSJed Brown   char             **root;      /* Packed root data, indexed by leaf rank */
40bf39f1bfSJed Brown   char             **leaf;      /* Packed leaf data, indexed by root rank */
4195fce210SBarry Smith   MPI_Request      *requests;   /* Array of root requests followed by leaf requests */
4295fce210SBarry Smith   PetscSFBasicPack next;
4395fce210SBarry Smith };
4495fce210SBarry Smith 
4595fce210SBarry Smith typedef struct {
4695fce210SBarry Smith   PetscMPIInt      tag;
47cf4b5b4fSJed Brown   PetscMPIInt      niranks;     /* Number of incoming ranks (ranks accessing my roots) */
48c943f53fSJed Brown   PetscMPIInt      ndiranks;    /* Number of incoming ranks (ranks accessing my roots) in distinguished set */
4995fce210SBarry Smith   PetscMPIInt      *iranks;     /* Array of ranks that reference my roots */
5095fce210SBarry Smith   PetscInt         itotal;      /* Total number of graph edges referencing my roots */
5195fce210SBarry Smith   PetscInt         *ioffset;    /* Array of length niranks+1 holding offset in irootloc[] for each rank */
5295fce210SBarry Smith   PetscInt         *irootloc;   /* Incoming roots referenced by ranks starting at ioffset[rank] */
5395fce210SBarry Smith   PetscSFBasicPack avail;       /* One or more entries per MPI Datatype, lazily constructed */
5495fce210SBarry Smith   PetscSFBasicPack inuse;       /* Buffers being used for transactions that have not yet completed */
5595fce210SBarry Smith } PetscSF_Basic;
5695fce210SBarry Smith 
575cf5d571SLisandro Dalcin #if !defined(PETSC_HAVE_MPI_TYPE_DUP)
585cf5d571SLisandro Dalcin PETSC_STATIC_INLINE int MPI_Type_dup(MPI_Datatype datatype,MPI_Datatype *newtype)
5995fce210SBarry Smith {
605cf5d571SLisandro Dalcin   int ierr;
615cf5d571SLisandro Dalcin   ierr = MPI_Type_contiguous(1,datatype,newtype); if (ierr) return ierr;
625cf5d571SLisandro Dalcin   ierr = MPI_Type_commit(newtype); if (ierr) return ierr;
635cf5d571SLisandro Dalcin   return MPI_SUCCESS;
6495fce210SBarry Smith }
6595fce210SBarry Smith #endif
6695fce210SBarry Smith 
6795fce210SBarry Smith /*
6895fce210SBarry Smith  * MPI_Reduce_local is not really useful because it can't handle sparse data and it vectorizes "in the wrong direction",
6995fce210SBarry Smith  * therefore we pack data types manually. This section defines packing routines for the standard data types.
7095fce210SBarry Smith  */
7195fce210SBarry Smith 
7295fce210SBarry Smith #define CPPJoin2_exp(a,b) a ## b
7395fce210SBarry Smith #define CPPJoin2(a,b) CPPJoin2_exp(a,b)
7495fce210SBarry Smith #define CPPJoin3_exp_(a,b,c) a ## b ## _ ## c
7595fce210SBarry Smith #define CPPJoin3_(a,b,c) CPPJoin3_exp_(a,b,c)
7695fce210SBarry Smith 
7795fce210SBarry Smith /* Basic types without addition */
782a74e5f4SJed Brown #define DEF_PackNoInit(type,BS)                                         \
792a74e5f4SJed Brown   static void CPPJoin3_(Pack_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,const void *unpacked,void *packed) { \
8095fce210SBarry Smith     const type *u = (const type*)unpacked;                              \
8195fce210SBarry Smith     type *p = (type*)packed;                                            \
822a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
832a74e5f4SJed Brown     for (i=0; i<n; i++)                                                 \
842a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
852a74e5f4SJed Brown         for (k=j; k<j+BS; k++)                                          \
862a74e5f4SJed Brown           p[i*bs+k] = u[idx[i]*bs+k];                                   \
8795fce210SBarry Smith   }                                                                     \
882a74e5f4SJed Brown   static void CPPJoin3_(UnpackInsert_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
8995fce210SBarry Smith     type *u = (type*)unpacked;                                          \
9095fce210SBarry Smith     const type *p = (const type*)packed;                                \
912a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
922a74e5f4SJed Brown     for (i=0; i<n; i++)                                                 \
932a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
942a74e5f4SJed Brown         for (k=j; k<j+BS; k++)                                          \
952a74e5f4SJed Brown           u[idx[i]*bs+k] = p[i*bs+k];                                   \
9695fce210SBarry Smith   }                                                                     \
972a74e5f4SJed Brown   static void CPPJoin3_(FetchAndInsert_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
9895fce210SBarry Smith     type *u = (type*)unpacked;                                          \
9995fce210SBarry Smith     type *p = (type*)packed;                                            \
1002a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
10195fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
1022a74e5f4SJed Brown       PetscInt ii = idx[i];                                             \
1032a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
1042a74e5f4SJed Brown         for (k=j; k<j+BS; k++) {                                        \
1052a74e5f4SJed Brown           type t = u[ii*bs+k];                                          \
1062a74e5f4SJed Brown           u[ii*bs+k] = p[i*bs+k];                                       \
1072a74e5f4SJed Brown           p[i*bs+k] = t;                                                \
1082a74e5f4SJed Brown         }                                                               \
10995fce210SBarry Smith     }                                                                   \
11095fce210SBarry Smith   }
11195fce210SBarry Smith 
11295fce210SBarry Smith /* Basic types defining addition */
1132a74e5f4SJed Brown #define DEF_PackAddNoInit(type,BS)                                      \
1142a74e5f4SJed Brown   DEF_PackNoInit(type,BS)                                               \
1152a74e5f4SJed Brown   static void CPPJoin3_(UnpackAdd_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
11695fce210SBarry Smith     type *u = (type*)unpacked;                                          \
11795fce210SBarry Smith     const type *p = (const type*)packed;                                \
1182a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
1192a74e5f4SJed Brown     for (i=0; i<n; i++)                                                 \
1202a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
1212a74e5f4SJed Brown         for (k=j; k<j+BS; k++)                                          \
1222a74e5f4SJed Brown           u[idx[i]*bs+k] += p[i*bs+k];                                  \
12395fce210SBarry Smith   }                                                                     \
1242a74e5f4SJed Brown   static void CPPJoin3_(FetchAndAdd_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
12595fce210SBarry Smith     type *u = (type*)unpacked;                                          \
12695fce210SBarry Smith     type *p = (type*)packed;                                            \
1272a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
12895fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
1292a74e5f4SJed Brown       PetscInt ii = idx[i];                                             \
1302a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
1312a74e5f4SJed Brown         for (k=j; k<j+BS; k++) {                                        \
1322a74e5f4SJed Brown           type t = u[ii*bs+k];                                          \
1332a74e5f4SJed Brown           u[ii*bs+k] = t + p[i*bs+k];                                   \
1342a74e5f4SJed Brown           p[i*bs+k] = t;                                                \
13595fce210SBarry Smith         }                                                               \
13653deab39SPeter Brune     }                                                                   \
1372a74e5f4SJed Brown   }                                                                     \
1382a74e5f4SJed Brown   static void CPPJoin3_(UnpackMult_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
13953deab39SPeter Brune     type *u = (type*)unpacked;                                          \
14053deab39SPeter Brune     const type *p = (const type*)packed;                                \
1412a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
1422a74e5f4SJed Brown     for (i=0; i<n; i++)                                                 \
1432a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
1442a74e5f4SJed Brown         for (k=j; k<j+BS; k++)                                          \
1452a74e5f4SJed Brown           u[idx[i]*bs+k] *= p[i*bs+k];                                  \
14653deab39SPeter Brune   }                                                                     \
1472a74e5f4SJed Brown   static void CPPJoin3_(FetchAndMult_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
14853deab39SPeter Brune     type *u = (type*)unpacked;                                          \
14953deab39SPeter Brune     type *p = (type*)packed;                                            \
1502a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
15153deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
1522a74e5f4SJed Brown       PetscInt ii = idx[i];                                             \
1532a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
1542a74e5f4SJed Brown         for (k=j; k<j+BS; k++) {                                        \
1552a74e5f4SJed Brown           type t = u[ii*bs+k];                                          \
1562a74e5f4SJed Brown           u[ii*bs+k] = t * p[i*bs+k];                                   \
1572a74e5f4SJed Brown           p[i*bs+k] = t;                                                \
1582a74e5f4SJed Brown         }                                                               \
15953deab39SPeter Brune     }                                                                   \
16095fce210SBarry Smith   }
1612a74e5f4SJed Brown #define DEF_Pack(type,BS)                                               \
1622a74e5f4SJed Brown   DEF_PackAddNoInit(type,BS)                                            \
1632a74e5f4SJed Brown   static void CPPJoin3_(PackInit_,type,BS)(PetscSFBasicPack link) {     \
1642a74e5f4SJed Brown     link->Pack = CPPJoin3_(Pack_,type,BS);                              \
1652a74e5f4SJed Brown     link->UnpackInsert = CPPJoin3_(UnpackInsert_,type,BS);              \
1662a74e5f4SJed Brown     link->UnpackAdd = CPPJoin3_(UnpackAdd_,type,BS);                    \
1672a74e5f4SJed Brown     link->UnpackMult = CPPJoin3_(UnpackMult_,type,BS);                  \
1682a74e5f4SJed Brown     link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,type,BS);          \
1692a74e5f4SJed Brown     link->FetchAndAdd = CPPJoin3_(FetchAndAdd_,type,BS);                \
1702a74e5f4SJed Brown     link->FetchAndMult = CPPJoin3_(FetchAndMult_,type,BS);              \
17195fce210SBarry Smith     link->unitbytes = sizeof(type);                                     \
17295fce210SBarry Smith   }
17395fce210SBarry Smith /* Comparable types */
17495fce210SBarry Smith #define DEF_PackCmp(type)                                               \
1752a74e5f4SJed Brown   DEF_PackAddNoInit(type,1)                                             \
1762a74e5f4SJed Brown   static void CPPJoin2(UnpackMax_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
17795fce210SBarry Smith     type *u = (type*)unpacked;                                          \
17895fce210SBarry Smith     const type *p = (const type*)packed;                                \
17995fce210SBarry Smith     PetscInt i;                                                         \
18095fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
18195fce210SBarry Smith       type v = u[idx[i]];                                               \
18295fce210SBarry Smith       u[idx[i]] = PetscMax(v,p[i]);                                     \
18395fce210SBarry Smith     }                                                                   \
18495fce210SBarry Smith   }                                                                     \
1852a74e5f4SJed Brown   static void CPPJoin2(UnpackMin_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
18695fce210SBarry Smith     type *u = (type*)unpacked;                                          \
18795fce210SBarry Smith     const type *p = (const type*)packed;                                \
18895fce210SBarry Smith     PetscInt i;                                                         \
18995fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
19095fce210SBarry Smith       type v = u[idx[i]];                                               \
19195fce210SBarry Smith       u[idx[i]] = PetscMin(v,p[i]);                                     \
19295fce210SBarry Smith     }                                                                   \
19395fce210SBarry Smith   }                                                                     \
1942a74e5f4SJed Brown   static void CPPJoin2(FetchAndMax_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
19595fce210SBarry Smith     type *u = (type*)unpacked;                                          \
19695fce210SBarry Smith     type *p = (type*)packed;                                            \
19795fce210SBarry Smith     PetscInt i;                                                         \
19895fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
19995fce210SBarry Smith       PetscInt j = idx[i];                                              \
20095fce210SBarry Smith       type v = u[j];                                                    \
20195fce210SBarry Smith       u[j] = PetscMax(v,p[i]);                                          \
20295fce210SBarry Smith       p[i] = v;                                                         \
20395fce210SBarry Smith     }                                                                   \
20495fce210SBarry Smith   }                                                                     \
2052a74e5f4SJed Brown   static void CPPJoin2(FetchAndMin_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
20695fce210SBarry Smith     type *u = (type*)unpacked;                                          \
20795fce210SBarry Smith     type *p = (type*)packed;                                            \
20895fce210SBarry Smith     PetscInt i;                                                         \
20995fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
21095fce210SBarry Smith       PetscInt j = idx[i];                                              \
21195fce210SBarry Smith       type v = u[j];                                                    \
21295fce210SBarry Smith       u[j] = PetscMin(v,p[i]);                                          \
21395fce210SBarry Smith       p[i] = v;                                                         \
21495fce210SBarry Smith     }                                                                   \
21595fce210SBarry Smith   }                                                                     \
21695fce210SBarry Smith   static void CPPJoin2(PackInit_,type)(PetscSFBasicPack link) {         \
2172a74e5f4SJed Brown     link->Pack = CPPJoin3_(Pack_,type,1);                               \
2182a74e5f4SJed Brown     link->UnpackInsert = CPPJoin3_(UnpackInsert_,type,1);               \
2192a74e5f4SJed Brown     link->UnpackAdd  = CPPJoin3_(UnpackAdd_,type,1);                    \
22095fce210SBarry Smith     link->UnpackMax  = CPPJoin2(UnpackMax_,type);                       \
22195fce210SBarry Smith     link->UnpackMin  = CPPJoin2(UnpackMin_,type);                       \
2222a74e5f4SJed Brown     link->UnpackMult = CPPJoin3_(UnpackMult_,type,1);                   \
2232a74e5f4SJed Brown     link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,type,1);           \
2242a74e5f4SJed Brown     link->FetchAndAdd = CPPJoin3_(FetchAndAdd_ ,type,1);                \
22595fce210SBarry Smith     link->FetchAndMax = CPPJoin2(FetchAndMax_ ,type);                   \
22695fce210SBarry Smith     link->FetchAndMin = CPPJoin2(FetchAndMin_ ,type);                   \
2272a74e5f4SJed Brown     link->FetchAndMult = CPPJoin3_(FetchAndMult_,type,1);               \
22895fce210SBarry Smith     link->unitbytes = sizeof(type);                                     \
22995fce210SBarry Smith   }
23095fce210SBarry Smith 
23153deab39SPeter Brune /* Logical Types */
23253deab39SPeter Brune #define DEF_PackLog(type)                                               \
2332a74e5f4SJed Brown   static void CPPJoin2(UnpackLAND_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
23453deab39SPeter Brune     type *u = (type*)unpacked;                                          \
23553deab39SPeter Brune     const type *p = (const type*)packed;                                \
23653deab39SPeter Brune     PetscInt i;                                                         \
23753deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
23853deab39SPeter Brune       type v = u[idx[i]];                                               \
23953deab39SPeter Brune       u[idx[i]] = v && p[i];                                            \
24053deab39SPeter Brune     }                                                                   \
24153deab39SPeter Brune   }                                                                     \
2422a74e5f4SJed Brown   static void CPPJoin2(UnpackLOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
24353deab39SPeter Brune     type *u = (type*)unpacked;                                          \
24453deab39SPeter Brune     const type *p = (const type*)packed;                                \
24553deab39SPeter Brune     PetscInt i;                                                         \
24653deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
24753deab39SPeter Brune       type v = u[idx[i]];                                               \
24853deab39SPeter Brune       u[idx[i]] = v || p[i];                                            \
24953deab39SPeter Brune     }                                                                   \
25053deab39SPeter Brune   }                                                                     \
2512a74e5f4SJed Brown   static void CPPJoin2(UnpackLXOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
25253deab39SPeter Brune     type *u = (type*)unpacked;                                          \
25353deab39SPeter Brune     const type *p = (const type*)packed;                                \
25453deab39SPeter Brune     PetscInt i;                                                         \
25553deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
25653deab39SPeter Brune       type v = u[idx[i]];                                               \
25753deab39SPeter Brune       u[idx[i]] = (!v)!=(!p[i]);                                        \
25853deab39SPeter Brune     }                                                                   \
25953deab39SPeter Brune   }                                                                     \
2602a74e5f4SJed Brown   static void CPPJoin2(FetchAndLAND_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
26153deab39SPeter Brune     type *u = (type*)unpacked;                                          \
26253deab39SPeter Brune     type *p = (type*)packed;                                            \
26353deab39SPeter Brune     PetscInt i;                                                         \
26453deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
26553deab39SPeter Brune       PetscInt j = idx[i];                                              \
26653deab39SPeter Brune       type v = u[j];                                                    \
26753deab39SPeter Brune       u[j] = v && p[i];                                                 \
26853deab39SPeter Brune       p[i] = v;                                                         \
26953deab39SPeter Brune     }                                                                   \
27053deab39SPeter Brune   }                                                                     \
2712a74e5f4SJed Brown   static void CPPJoin2(FetchAndLOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
27253deab39SPeter Brune     type *u = (type*)unpacked;                                          \
27353deab39SPeter Brune     type *p = (type*)packed;                                            \
27453deab39SPeter Brune     PetscInt i;                                                         \
27553deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
27653deab39SPeter Brune       PetscInt j = idx[i];                                              \
27753deab39SPeter Brune       type v = u[j];                                                    \
27853deab39SPeter Brune       u[j] = v || p[i];                                                 \
27953deab39SPeter Brune       p[i] = v;                                                         \
28053deab39SPeter Brune     }                                                                   \
28153deab39SPeter Brune   }                                                                     \
2822a74e5f4SJed Brown   static void CPPJoin2(FetchAndLXOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
28353deab39SPeter Brune     type *u = (type*)unpacked;                                          \
28453deab39SPeter Brune     type *p = (type*)packed;                                            \
28553deab39SPeter Brune     PetscInt i;                                                         \
28653deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
28753deab39SPeter Brune       PetscInt j = idx[i];                                              \
28853deab39SPeter Brune       type v = u[j];                                                    \
28953deab39SPeter Brune       u[j] = (!v)!=(!p[i]);                                             \
29053deab39SPeter Brune       p[i] = v;                                                         \
29153deab39SPeter Brune     }                                                                   \
29253deab39SPeter Brune   }                                                                     \
29353deab39SPeter Brune   static void CPPJoin2(PackInit_Logical_,type)(PetscSFBasicPack link) { \
29453deab39SPeter Brune     link->UnpackLAND = CPPJoin2(UnpackLAND_,type);                      \
29553deab39SPeter Brune     link->UnpackLOR  = CPPJoin2(UnpackLOR_,type);                       \
29653deab39SPeter Brune     link->UnpackLXOR = CPPJoin2(UnpackLXOR_,type);                      \
29753deab39SPeter Brune     link->FetchAndLAND = CPPJoin2(FetchAndLAND_,type);                  \
29853deab39SPeter Brune     link->FetchAndLOR  = CPPJoin2(FetchAndLOR_,type);                   \
29953deab39SPeter Brune     link->FetchAndLXOR = CPPJoin2(FetchAndLXOR_,type);                  \
30053deab39SPeter Brune   }
30153deab39SPeter Brune 
30253deab39SPeter Brune 
30353deab39SPeter Brune /* Bitwise Types */
30453deab39SPeter Brune #define DEF_PackBit(type)                                               \
3052a74e5f4SJed Brown   static void CPPJoin2(UnpackBAND_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
30653deab39SPeter Brune     type *u = (type*)unpacked;                                          \
30753deab39SPeter Brune     const type *p = (const type*)packed;                                \
30853deab39SPeter Brune     PetscInt i;                                                         \
30953deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
31053deab39SPeter Brune       type v = u[idx[i]];                                               \
31153deab39SPeter Brune       u[idx[i]] = v & p[i];                                             \
31253deab39SPeter Brune     }                                                                   \
31353deab39SPeter Brune   }                                                                     \
3142a74e5f4SJed Brown   static void CPPJoin2(UnpackBOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
31553deab39SPeter Brune     type *u = (type*)unpacked;                                          \
31653deab39SPeter Brune     const type *p = (const type*)packed;                                \
31753deab39SPeter Brune     PetscInt i;                                                         \
31853deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
31953deab39SPeter Brune       type v = u[idx[i]];                                               \
32053deab39SPeter Brune       u[idx[i]] = v | p[i];                                             \
32153deab39SPeter Brune     }                                                                   \
32253deab39SPeter Brune   }                                                                     \
3232a74e5f4SJed Brown   static void CPPJoin2(UnpackBXOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
32453deab39SPeter Brune     type *u = (type*)unpacked;                                          \
32553deab39SPeter Brune     const type *p = (const type*)packed;                                \
32653deab39SPeter Brune     PetscInt i;                                                         \
32753deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
32853deab39SPeter Brune       type v = u[idx[i]];                                               \
32953deab39SPeter Brune       u[idx[i]] = v^p[i];                                               \
33053deab39SPeter Brune     }                                                                   \
33153deab39SPeter Brune   }                                                                     \
3322a74e5f4SJed Brown   static void CPPJoin2(FetchAndBAND_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
33353deab39SPeter Brune     type *u = (type*)unpacked;                                          \
33453deab39SPeter Brune     type *p = (type*)packed;                                            \
33553deab39SPeter Brune     PetscInt i;                                                         \
33653deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
33753deab39SPeter Brune       PetscInt j = idx[i];                                              \
33853deab39SPeter Brune       type v = u[j];                                                    \
33953deab39SPeter Brune       u[j] = v & p[i];                                                  \
34053deab39SPeter Brune       p[i] = v;                                                         \
34153deab39SPeter Brune     }                                                                   \
34253deab39SPeter Brune   }                                                                     \
3432a74e5f4SJed Brown   static void CPPJoin2(FetchAndBOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
34453deab39SPeter Brune     type *u = (type*)unpacked;                                          \
34553deab39SPeter Brune     type *p = (type*)packed;                                            \
34653deab39SPeter Brune     PetscInt i;                                                         \
34753deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
34853deab39SPeter Brune       PetscInt j = idx[i];                                              \
34953deab39SPeter Brune       type v = u[j];                                                    \
35053deab39SPeter Brune       u[j] = v | p[i];                                                  \
35153deab39SPeter Brune       p[i] = v;                                                         \
35253deab39SPeter Brune     }                                                                   \
35353deab39SPeter Brune   }                                                                     \
3542a74e5f4SJed Brown   static void CPPJoin2(FetchAndBXOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
35553deab39SPeter Brune     type *u = (type*)unpacked;                                          \
35653deab39SPeter Brune     type *p = (type*)packed;                                            \
35753deab39SPeter Brune     PetscInt i;                                                         \
35853deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
35953deab39SPeter Brune       PetscInt j = idx[i];                                              \
36053deab39SPeter Brune       type v = u[j];                                                    \
36153deab39SPeter Brune       u[j] = v^p[i];                                                    \
36253deab39SPeter Brune       p[i] = v;                                                         \
36353deab39SPeter Brune     }                                                                   \
36453deab39SPeter Brune   }                                                                     \
36553deab39SPeter Brune   static void CPPJoin2(PackInit_Bitwise_,type)(PetscSFBasicPack link) { \
36653deab39SPeter Brune     link->UnpackBAND = CPPJoin2(UnpackBAND_,type);                      \
36753deab39SPeter Brune     link->UnpackBOR  = CPPJoin2(UnpackBOR_,type);                       \
36853deab39SPeter Brune     link->UnpackBXOR = CPPJoin2(UnpackBXOR_,type);                      \
36953deab39SPeter Brune     link->FetchAndBAND = CPPJoin2(FetchAndBAND_,type);                  \
37053deab39SPeter Brune     link->FetchAndBOR  = CPPJoin2(FetchAndBOR_,type);                   \
37153deab39SPeter Brune     link->FetchAndBXOR = CPPJoin2(FetchAndBXOR_,type);                  \
37253deab39SPeter Brune   }
37353deab39SPeter Brune 
37495fce210SBarry Smith /* Pair types */
37595fce210SBarry Smith #define CPPJoinloc_exp(base,op,t1,t2) base ## op ## loc_ ## t1 ## _ ## t2
37695fce210SBarry Smith #define CPPJoinloc(base,op,t1,t2) CPPJoinloc_exp(base,op,t1,t2)
37795fce210SBarry Smith #define PairType(type1,type2) CPPJoin3_(_pairtype_,type1,type2)
37895fce210SBarry Smith #define DEF_UnpackXloc(type1,type2,locname,op)                              \
3792a74e5f4SJed Brown   static void CPPJoinloc(Unpack,locname,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
38095fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;        \
38195fce210SBarry Smith     const PairType(type1,type2) *p = (const PairType(type1,type2)*)packed; \
38295fce210SBarry Smith     PetscInt i;                                                         \
38395fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
38495fce210SBarry Smith       PetscInt j = idx[i];                                              \
38595fce210SBarry Smith       if (p[i].a op u[j].a) {                                           \
38695fce210SBarry Smith         u[j].a = p[i].a;                                                \
38795fce210SBarry Smith         u[j].b = p[i].b;                                                \
38895fce210SBarry Smith       } else if (u[j].a == p[i].a) {                                    \
38995fce210SBarry Smith         u[j].b = PetscMin(u[j].b,p[i].b);                               \
39095fce210SBarry Smith       }                                                                 \
39195fce210SBarry Smith     }                                                                   \
39295fce210SBarry Smith   }                                                                     \
3932a74e5f4SJed Brown   static void CPPJoinloc(FetchAnd,locname,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
39495fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;        \
39595fce210SBarry Smith     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;          \
39695fce210SBarry Smith     PetscInt i;                                                         \
39795fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
39895fce210SBarry Smith       PetscInt j = idx[i];                                              \
39995fce210SBarry Smith       PairType(type1,type2) v;                                          \
40095fce210SBarry Smith       v.a = u[j].a;                                                     \
40195fce210SBarry Smith       v.b = u[j].b;                                                     \
40295fce210SBarry Smith       if (p[i].a op u[j].a) {                                           \
40395fce210SBarry Smith         u[j].a = p[i].a;                                                \
40495fce210SBarry Smith         u[j].b = p[i].b;                                                \
40595fce210SBarry Smith       } else if (u[j].a == p[i].a) {                                    \
40695fce210SBarry Smith         u[j].b = PetscMin(u[j].b,p[i].b);                               \
40795fce210SBarry Smith       }                                                                 \
40895fce210SBarry Smith       p[i].a = v.a;                                                     \
40995fce210SBarry Smith       p[i].b = v.b;                                                     \
41095fce210SBarry Smith     }                                                                   \
41195fce210SBarry Smith   }
41295fce210SBarry Smith #define DEF_PackPair(type1,type2)                                       \
41395fce210SBarry Smith   typedef struct {type1 a; type2 b;} PairType(type1,type2);             \
4142a74e5f4SJed Brown   static void CPPJoin3_(Pack_,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,const void *unpacked,void *packed) { \
41595fce210SBarry Smith     const PairType(type1,type2) *u = (const PairType(type1,type2)*)unpacked; \
41695fce210SBarry Smith     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;          \
41795fce210SBarry Smith     PetscInt i;                                                         \
41895fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
41995fce210SBarry Smith       p[i].a = u[idx[i]].a;                                             \
42095fce210SBarry Smith       p[i].b = u[idx[i]].b;                                             \
42195fce210SBarry Smith     }                                                                   \
42295fce210SBarry Smith   }                                                                     \
4232a74e5f4SJed Brown   static void CPPJoin3_(UnpackInsert_,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
42495fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;       \
42595fce210SBarry Smith     const PairType(type1,type2) *p = (const PairType(type1,type2)*)packed; \
42695fce210SBarry Smith     PetscInt i;                                                         \
42795fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
42895fce210SBarry Smith       u[idx[i]].a = p[i].a;                                             \
42995fce210SBarry Smith       u[idx[i]].b = p[i].b;                                             \
43095fce210SBarry Smith     }                                                                   \
43195fce210SBarry Smith   }                                                                     \
4322a74e5f4SJed Brown   static void CPPJoin3_(UnpackAdd_,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
43395fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;       \
43495fce210SBarry Smith     const PairType(type1,type2) *p = (const PairType(type1,type2)*)packed; \
43595fce210SBarry Smith     PetscInt i;                                                         \
43695fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
43795fce210SBarry Smith       u[idx[i]].a += p[i].a;                                            \
43895fce210SBarry Smith       u[idx[i]].b += p[i].b;                                            \
43995fce210SBarry Smith     }                                                                   \
44095fce210SBarry Smith   }                                                                     \
4412a74e5f4SJed Brown   static void CPPJoin3_(FetchAndInsert_,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
44295fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;        \
44395fce210SBarry Smith     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;          \
44495fce210SBarry Smith     PetscInt i;                                                         \
44595fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
44695fce210SBarry Smith       PetscInt j = idx[i];                                              \
44795fce210SBarry Smith       PairType(type1,type2) v;                                          \
44895fce210SBarry Smith       v.a = u[j].a;                                                     \
44995fce210SBarry Smith       v.b = u[j].b;                                                     \
45095fce210SBarry Smith       u[j].a = p[i].a;                                                  \
45195fce210SBarry Smith       u[j].b = p[i].b;                                                  \
45295fce210SBarry Smith       p[i].a = v.a;                                                     \
45395fce210SBarry Smith       p[i].b = v.b;                                                     \
45495fce210SBarry Smith     }                                                                   \
45595fce210SBarry Smith   }                                                                     \
4562a74e5f4SJed Brown   static void FetchAndAdd_ ## type1 ## _ ## type2(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
45795fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;       \
45895fce210SBarry Smith     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;         \
45995fce210SBarry Smith     PetscInt i;                                                         \
46095fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
46195fce210SBarry Smith       PetscInt j = idx[i];                                              \
46295fce210SBarry Smith       PairType(type1,type2) v;                                          \
46395fce210SBarry Smith       v.a = u[j].a;                                                     \
46495fce210SBarry Smith       v.b = u[j].b;                                                     \
46595fce210SBarry Smith       u[j].a = v.a + p[i].a;                                            \
46695fce210SBarry Smith       u[j].b = v.b + p[i].b;                                            \
46795fce210SBarry Smith       p[i].a = v.a;                                                     \
46895fce210SBarry Smith       p[i].b = v.b;                                                     \
46995fce210SBarry Smith     }                                                                   \
47095fce210SBarry Smith   }                                                                     \
47195fce210SBarry Smith   DEF_UnpackXloc(type1,type2,Max,>)                                     \
47295fce210SBarry Smith   DEF_UnpackXloc(type1,type2,Min,<)                                     \
47395fce210SBarry Smith   static void CPPJoin3_(PackInit_,type1,type2)(PetscSFBasicPack link) { \
47495fce210SBarry Smith     link->Pack = CPPJoin3_(Pack_,type1,type2);                          \
47595fce210SBarry Smith     link->UnpackInsert = CPPJoin3_(UnpackInsert_,type1,type2);          \
47695fce210SBarry Smith     link->UnpackAdd = CPPJoin3_(UnpackAdd_,type1,type2);                \
47795fce210SBarry Smith     link->UnpackMaxloc = CPPJoin3_(UnpackMaxloc_,type1,type2);          \
47895fce210SBarry Smith     link->UnpackMinloc = CPPJoin3_(UnpackMinloc_,type1,type2);          \
47995fce210SBarry Smith     link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,type1,type2);      \
48095fce210SBarry Smith     link->FetchAndAdd = CPPJoin3_(FetchAndAdd_,type1,type2);            \
48195fce210SBarry Smith     link->FetchAndMaxloc = CPPJoin3_(FetchAndMaxloc_,type1,type2);      \
48295fce210SBarry Smith     link->FetchAndMinloc = CPPJoin3_(FetchAndMinloc_,type1,type2);      \
48395fce210SBarry Smith     link->unitbytes = sizeof(PairType(type1,type2));                    \
48495fce210SBarry Smith   }
48595fce210SBarry Smith 
48695fce210SBarry Smith /* Currently only dumb blocks of data */
48795fce210SBarry Smith #define BlockType(unit,count) CPPJoin3_(_blocktype_,unit,count)
48895fce210SBarry Smith #define DEF_Block(unit,count)                                           \
48995fce210SBarry Smith   typedef struct {unit v[count];} BlockType(unit,count);                \
4902a74e5f4SJed Brown   DEF_PackNoInit(BlockType(unit,count),1)                               \
49195fce210SBarry Smith   static void CPPJoin3_(PackInit_block_,unit,count)(PetscSFBasicPack link) { \
4922a74e5f4SJed Brown     link->Pack = CPPJoin3_(Pack_,BlockType(unit,count),1);               \
4932a74e5f4SJed Brown     link->UnpackInsert = CPPJoin3_(UnpackInsert_,BlockType(unit,count),1); \
4942a74e5f4SJed Brown     link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,BlockType(unit,count),1); \
49595fce210SBarry Smith     link->unitbytes = sizeof(BlockType(unit,count));                    \
49695fce210SBarry Smith   }
49795fce210SBarry Smith 
49895fce210SBarry Smith DEF_PackCmp(int)
49953deab39SPeter Brune DEF_PackBit(int)
50053deab39SPeter Brune DEF_PackLog(int)
50195fce210SBarry Smith DEF_PackCmp(PetscInt)
50253deab39SPeter Brune DEF_PackBit(PetscInt)
50353deab39SPeter Brune DEF_PackLog(PetscInt)
5042a74e5f4SJed Brown DEF_Pack(PetscInt,2)
5052a74e5f4SJed Brown DEF_Pack(PetscInt,3)
5062a74e5f4SJed Brown DEF_Pack(PetscInt,4)
5072a74e5f4SJed Brown DEF_Pack(PetscInt,5)
5082a74e5f4SJed Brown DEF_Pack(PetscInt,7)
50995fce210SBarry Smith DEF_PackCmp(PetscReal)
51053deab39SPeter Brune DEF_PackLog(PetscReal)
5112a74e5f4SJed Brown DEF_Pack(PetscReal,2)
5122a74e5f4SJed Brown DEF_Pack(PetscReal,3)
5132a74e5f4SJed Brown DEF_Pack(PetscReal,4)
5142a74e5f4SJed Brown DEF_Pack(PetscReal,5)
5152a74e5f4SJed Brown DEF_Pack(PetscReal,7)
51695fce210SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
5172a74e5f4SJed Brown DEF_Pack(PetscComplex,1)
5182a74e5f4SJed Brown DEF_Pack(PetscComplex,2)
5192a74e5f4SJed Brown DEF_Pack(PetscComplex,3)
5202a74e5f4SJed Brown DEF_Pack(PetscComplex,4)
5212a74e5f4SJed Brown DEF_Pack(PetscComplex,5)
5222a74e5f4SJed Brown DEF_Pack(PetscComplex,7)
52395fce210SBarry Smith #endif
52495fce210SBarry Smith DEF_PackPair(int,int)
52595fce210SBarry Smith DEF_PackPair(PetscInt,PetscInt)
52695fce210SBarry Smith DEF_Block(int,1)
52795fce210SBarry Smith DEF_Block(int,2)
52895fce210SBarry Smith DEF_Block(int,3)
52995fce210SBarry Smith DEF_Block(int,4)
53095fce210SBarry Smith DEF_Block(int,5)
53195fce210SBarry Smith DEF_Block(int,6)
53295fce210SBarry Smith DEF_Block(int,7)
53395fce210SBarry Smith DEF_Block(int,8)
53495fce210SBarry Smith 
53595fce210SBarry Smith static PetscErrorCode PetscSFSetUp_Basic(PetscSF sf)
53695fce210SBarry Smith {
53795fce210SBarry Smith   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
53895fce210SBarry Smith   PetscErrorCode ierr;
53995fce210SBarry Smith   PetscInt *rlengths,*ilengths,i;
540c943f53fSJed Brown   PetscMPIInt rank,niranks,*iranks;
54195fce210SBarry Smith   MPI_Comm comm;
542b5a8e515SJed Brown   MPI_Group group;
54395fce210SBarry Smith   MPI_Request *rootreqs,*leafreqs;
54495fce210SBarry Smith 
54595fce210SBarry Smith   PetscFunctionBegin;
546b5a8e515SJed Brown   ierr = MPI_Comm_group(PETSC_COMM_SELF,&group);CHKERRQ(ierr);
547b5a8e515SJed Brown   ierr = PetscSFSetUpRanks(sf,group);CHKERRQ(ierr);
548b5a8e515SJed Brown   ierr = MPI_Group_free(&group);CHKERRQ(ierr);
54995fce210SBarry Smith   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
55095fce210SBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)sf,&bas->tag);CHKERRQ(ierr);
551c943f53fSJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
55295fce210SBarry Smith   /*
55395fce210SBarry Smith    * Inform roots about how many leaves and from which ranks
55495fce210SBarry Smith    */
555785e854fSJed Brown   ierr = PetscMalloc1(sf->nranks,&rlengths);CHKERRQ(ierr);
55695fce210SBarry Smith   /* Determine number, sending ranks, and length of incoming */
55795fce210SBarry Smith   for (i=0; i<sf->nranks; i++) {
55895fce210SBarry Smith     rlengths[i] = sf->roffset[i+1] - sf->roffset[i]; /* Number of roots referenced by my leaves; for rank sf->ranks[i] */
55995fce210SBarry Smith   }
560c943f53fSJed Brown   ierr = PetscCommBuildTwoSided(comm,1,MPIU_INT,sf->nranks-sf->ndranks,sf->ranks+sf->ndranks,rlengths+sf->ndranks,&niranks,&iranks,(void**)&ilengths);CHKERRQ(ierr);
561c943f53fSJed Brown 
5620b899082SJunchao Zhang   /* Sort iranks. See use of VecScatterGetRemoteOrdered_Private() in MatGetBrowsOfAoCols_MPIAIJ() on why.
5630b899082SJunchao Zhang      We could sort ranks there at the price of allocating extra working arrays. Presumably, niranks is
5640b899082SJunchao Zhang      small and the sorting is cheap.
5650b899082SJunchao Zhang    */
5660b899082SJunchao Zhang   ierr = PetscSortMPIIntWithIntArray(niranks,iranks,ilengths);CHKERRQ(ierr);
5670b899082SJunchao Zhang 
568c943f53fSJed Brown   /* Partition into distinguished and non-distinguished incoming ranks */
569c943f53fSJed Brown   bas->ndiranks = sf->ndranks;
570c943f53fSJed Brown   bas->niranks = bas->ndiranks + niranks;
571c943f53fSJed Brown   ierr = PetscMalloc2(bas->niranks,&bas->iranks,bas->niranks+1,&bas->ioffset);CHKERRQ(ierr);
572c943f53fSJed Brown   bas->ioffset[0] = 0;
573c943f53fSJed Brown   for (i=0; i<bas->ndiranks; i++) {
574c943f53fSJed Brown     bas->iranks[i] = sf->ranks[i];
575c943f53fSJed Brown     bas->ioffset[i+1] = bas->ioffset[i] + rlengths[i];
576c943f53fSJed Brown   }
577c943f53fSJed Brown   if (bas->ndiranks > 1 || (bas->ndiranks == 1 && bas->iranks[0] != rank)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Broken setup for shared ranks");
578c943f53fSJed Brown   for ( ; i<bas->niranks; i++) {
579c943f53fSJed Brown     bas->iranks[i] = iranks[i-bas->ndiranks];
580c943f53fSJed Brown     bas->ioffset[i+1] = bas->ioffset[i] + ilengths[i-bas->ndiranks];
581c943f53fSJed Brown   }
582c943f53fSJed Brown   bas->itotal = bas->ioffset[i];
58395fce210SBarry Smith   ierr = PetscFree(rlengths);CHKERRQ(ierr);
584c943f53fSJed Brown   ierr = PetscFree(iranks);CHKERRQ(ierr);
585c943f53fSJed Brown   ierr = PetscFree(ilengths);CHKERRQ(ierr);
58695fce210SBarry Smith 
58795fce210SBarry Smith   /* Send leaf identities to roots */
588c943f53fSJed Brown   ierr = PetscMalloc1(bas->itotal,&bas->irootloc);CHKERRQ(ierr);
589c943f53fSJed Brown   ierr = PetscMalloc2(bas->niranks-bas->ndiranks,&rootreqs,sf->nranks-sf->ndranks,&leafreqs);CHKERRQ(ierr);
590c943f53fSJed Brown   for (i=bas->ndiranks; i<bas->niranks; i++) {
591c943f53fSJed Brown     ierr = MPI_Irecv(bas->irootloc+bas->ioffset[i],bas->ioffset[i+1]-bas->ioffset[i],MPIU_INT,bas->iranks[i],bas->tag,comm,&rootreqs[i-bas->ndiranks]);CHKERRQ(ierr);
59295fce210SBarry Smith   }
59395fce210SBarry Smith   for (i=0; i<sf->nranks; i++) {
59495fce210SBarry Smith     PetscMPIInt npoints;
59595fce210SBarry Smith     ierr = PetscMPIIntCast(sf->roffset[i+1] - sf->roffset[i],&npoints);CHKERRQ(ierr);
596c943f53fSJed Brown     if (i < sf->ndranks) {
597c943f53fSJed Brown       if (sf->ranks[i] != rank) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot interpret distinguished leaf rank");
598c943f53fSJed Brown       if (bas->iranks[0] != rank) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot interpret distinguished root rank");
599c943f53fSJed Brown       if (npoints != bas->ioffset[1]-bas->ioffset[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Distinguished rank exchange has mismatched lengths");
600c943f53fSJed Brown       ierr = PetscMemcpy(bas->irootloc+bas->ioffset[0],sf->rremote+sf->roffset[i],npoints*sizeof(bas->irootloc[0]));CHKERRQ(ierr);
601c943f53fSJed Brown       continue;
60295fce210SBarry Smith     }
603c943f53fSJed Brown     ierr = MPI_Isend(sf->rremote+sf->roffset[i],npoints,MPIU_INT,sf->ranks[i],bas->tag,comm,&leafreqs[i-sf->ndranks]);CHKERRQ(ierr);
604c943f53fSJed Brown   }
605c943f53fSJed Brown   ierr = MPI_Waitall(bas->niranks-bas->ndiranks,rootreqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
606c943f53fSJed Brown   ierr = MPI_Waitall(sf->nranks-sf->ndranks,leafreqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
607eb9baa12SBarry Smith   ierr = PetscFree2(rootreqs,leafreqs);CHKERRQ(ierr);
60895fce210SBarry Smith   PetscFunctionReturn(0);
60995fce210SBarry Smith }
61095fce210SBarry Smith 
61195fce210SBarry Smith static PetscErrorCode PetscSFBasicPackTypeSetup(PetscSFBasicPack link,MPI_Datatype unit)
61295fce210SBarry Smith {
61395fce210SBarry Smith   PetscErrorCode ierr;
61495fce210SBarry Smith   PetscBool      isInt,isPetscInt,isPetscReal,is2Int,is2PetscInt;
615e37f7bbbSSatish Balay   PetscInt       nPetscIntContig,nPetscRealContig;
616*77b7e48dSJunchao Zhang   PetscMPIInt    ni,na,nd,combiner;
61795fce210SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
61895fce210SBarry Smith   PetscBool isPetscComplex;
619e37f7bbbSSatish Balay   PetscInt nPetscComplexContig;
62095fce210SBarry Smith #endif
62195fce210SBarry Smith 
62295fce210SBarry Smith   PetscFunctionBegin;
62395fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPI_INT,&isInt);CHKERRQ(ierr);
62495fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPIU_INT,&isPetscInt);CHKERRQ(ierr);
6252a74e5f4SJed Brown   ierr = MPIPetsc_Type_compare_contig(unit,MPIU_INT,&nPetscIntContig);CHKERRQ(ierr);
62695fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPIU_REAL,&isPetscReal);CHKERRQ(ierr);
6272a74e5f4SJed Brown   ierr = MPIPetsc_Type_compare_contig(unit,MPIU_REAL,&nPetscRealContig);CHKERRQ(ierr);
62895fce210SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
62995fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPIU_COMPLEX,&isPetscComplex);CHKERRQ(ierr);
6302a74e5f4SJed Brown   ierr = MPIPetsc_Type_compare_contig(unit,MPIU_COMPLEX,&nPetscComplexContig);CHKERRQ(ierr);
63195fce210SBarry Smith #endif
63295fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPI_2INT,&is2Int);CHKERRQ(ierr);
63395fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPIU_2INT,&is2PetscInt);CHKERRQ(ierr);
634*77b7e48dSJunchao Zhang   ierr = MPI_Type_get_envelope(unit,&ni,&na,&nd,&combiner);CHKERRQ(ierr);
635*77b7e48dSJunchao Zhang   link->isbuiltin = (combiner == MPI_COMBINER_NAMED) ? PETSC_TRUE : PETSC_FALSE;
6362a74e5f4SJed Brown   link->bs = 1;
63753deab39SPeter Brune   if (isInt) {PackInit_int(link); PackInit_Logical_int(link); PackInit_Bitwise_int(link);}
63853deab39SPeter Brune   else if (isPetscInt) {PackInit_PetscInt(link); PackInit_Logical_PetscInt(link); PackInit_Bitwise_PetscInt(link);}
63953deab39SPeter Brune   else if (isPetscReal) {PackInit_PetscReal(link); PackInit_Logical_PetscReal(link);}
64095fce210SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
6412a74e5f4SJed Brown   else if (isPetscComplex) PackInit_PetscComplex_1(link);
64295fce210SBarry Smith #endif
64395fce210SBarry Smith   else if (is2Int) PackInit_int_int(link);
64495fce210SBarry Smith   else if (is2PetscInt) PackInit_PetscInt_PetscInt(link);
6452a74e5f4SJed Brown   else if (nPetscIntContig) {
6462a74e5f4SJed Brown     if (nPetscIntContig%7 == 0) PackInit_PetscInt_7(link);
6472a74e5f4SJed Brown     else if (nPetscIntContig%5 == 0) PackInit_PetscInt_5(link);
6482a74e5f4SJed Brown     else if (nPetscIntContig%4 == 0) PackInit_PetscInt_4(link);
6492a74e5f4SJed Brown     else if (nPetscIntContig%3 == 0) PackInit_PetscInt_3(link);
6502a74e5f4SJed Brown     else if (nPetscIntContig%2 == 0) PackInit_PetscInt_2(link);
6512a74e5f4SJed Brown     else PackInit_PetscInt(link);
6522a74e5f4SJed Brown     link->bs = nPetscIntContig;
6532a74e5f4SJed Brown     link->unitbytes *= nPetscIntContig;
6542a74e5f4SJed Brown   } else if (nPetscRealContig) {
6552a74e5f4SJed Brown     if (nPetscRealContig%7 == 0) PackInit_PetscReal_7(link);
6562a74e5f4SJed Brown     else if (nPetscRealContig%5 == 0) PackInit_PetscReal_5(link);
6572a74e5f4SJed Brown     else if (nPetscRealContig%4 == 0) PackInit_PetscReal_4(link);
6582a74e5f4SJed Brown     else if (nPetscRealContig%3 == 0) PackInit_PetscReal_3(link);
6592a74e5f4SJed Brown     else if (nPetscRealContig%2 == 0) PackInit_PetscReal_2(link);
6602a74e5f4SJed Brown     else PackInit_PetscReal(link);
6612a74e5f4SJed Brown     link->bs = nPetscRealContig;
6622a74e5f4SJed Brown     link->unitbytes *= nPetscRealContig;
6632a74e5f4SJed Brown #if defined(PETSC_HAVE_COMPLEX)
6642a74e5f4SJed Brown   } else if (nPetscComplexContig) {
6652a74e5f4SJed Brown     if (nPetscComplexContig%7 == 0) PackInit_PetscComplex_7(link);
6662a74e5f4SJed Brown     else if (nPetscComplexContig%5 == 0) PackInit_PetscComplex_5(link);
6672a74e5f4SJed Brown     else if (nPetscComplexContig%4 == 0) PackInit_PetscComplex_4(link);
6682a74e5f4SJed Brown     else if (nPetscComplexContig%3 == 0) PackInit_PetscComplex_3(link);
6692a74e5f4SJed Brown     else if (nPetscComplexContig%2 == 0) PackInit_PetscComplex_2(link);
6702a74e5f4SJed Brown     else PackInit_PetscComplex_1(link);
6712a74e5f4SJed Brown     link->bs = nPetscComplexContig;
6722a74e5f4SJed Brown     link->unitbytes *= nPetscComplexContig;
6732a74e5f4SJed Brown #endif
6742a74e5f4SJed Brown   } else {
675c3a59b84SJed Brown     MPI_Aint lb,bytes;
676c3a59b84SJed Brown     ierr = MPI_Type_get_extent(unit,&lb,&bytes);CHKERRQ(ierr);
677c3a59b84SJed Brown     if (lb != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Datatype with nonzero lower bound %ld\n",(long)lb);
67895fce210SBarry Smith     if (bytes % sizeof(int)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for type size not divisible by %D",sizeof(int));
67995fce210SBarry Smith     switch (bytes / sizeof(int)) {
68095fce210SBarry Smith     case 1: PackInit_block_int_1(link); break;
68195fce210SBarry Smith     case 2: PackInit_block_int_2(link); break;
68295fce210SBarry Smith     case 3: PackInit_block_int_3(link); break;
68395fce210SBarry Smith     case 4: PackInit_block_int_4(link); break;
68495fce210SBarry Smith     case 5: PackInit_block_int_5(link); break;
68595fce210SBarry Smith     case 6: PackInit_block_int_6(link); break;
68695fce210SBarry Smith     case 7: PackInit_block_int_7(link); break;
68795fce210SBarry Smith     case 8: PackInit_block_int_8(link); break;
68895fce210SBarry Smith     default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for arbitrary block sizes");
68995fce210SBarry Smith     }
69095fce210SBarry Smith   }
691*77b7e48dSJunchao Zhang   if (link->isbuiltin) link->unit = unit; /* builtin datatypes are common. Make it fast */
692*77b7e48dSJunchao Zhang   else {ierr = MPI_Type_dup(unit,&link->unit);CHKERRQ(ierr);}
69395fce210SBarry Smith   PetscFunctionReturn(0);
69495fce210SBarry Smith }
69595fce210SBarry Smith 
6962a74e5f4SJed Brown static PetscErrorCode PetscSFBasicPackGetUnpackOp(PetscSF sf,PetscSFBasicPack link,MPI_Op op,void (**UnpackOp)(PetscInt,PetscInt,const PetscInt*,void*,const void*))
69795fce210SBarry Smith {
69895fce210SBarry Smith   PetscFunctionBegin;
69995fce210SBarry Smith   *UnpackOp = NULL;
7008bfbc91cSJed Brown   if (op == MPIU_REPLACE) *UnpackOp = link->UnpackInsert;
70195fce210SBarry Smith   else if (op == MPI_SUM || op == MPIU_SUM) *UnpackOp = link->UnpackAdd;
70253deab39SPeter Brune   else if (op == MPI_PROD) *UnpackOp = link->UnpackMult;
70395fce210SBarry Smith   else if (op == MPI_MAX || op == MPIU_MAX) *UnpackOp = link->UnpackMax;
70495fce210SBarry Smith   else if (op == MPI_MIN || op == MPIU_MIN) *UnpackOp = link->UnpackMin;
70553deab39SPeter Brune   else if (op == MPI_LAND) *UnpackOp = link->UnpackLAND;
70653deab39SPeter Brune   else if (op == MPI_BAND) *UnpackOp = link->UnpackBAND;
70753deab39SPeter Brune   else if (op == MPI_LOR) *UnpackOp = link->UnpackLOR;
70853deab39SPeter Brune   else if (op == MPI_BOR) *UnpackOp = link->UnpackBOR;
70953deab39SPeter Brune   else if (op == MPI_LXOR) *UnpackOp = link->UnpackLXOR;
71053deab39SPeter Brune   else if (op == MPI_BXOR) *UnpackOp = link->UnpackBXOR;
71195fce210SBarry Smith   else if (op == MPI_MAXLOC) *UnpackOp = link->UnpackMaxloc;
71295fce210SBarry Smith   else if (op == MPI_MINLOC) *UnpackOp = link->UnpackMinloc;
7131620da3bSToby Isaac   else *UnpackOp = NULL;
71495fce210SBarry Smith   PetscFunctionReturn(0);
71595fce210SBarry Smith }
7162a74e5f4SJed Brown static PetscErrorCode PetscSFBasicPackGetFetchAndOp(PetscSF sf,PetscSFBasicPack link,MPI_Op op,void (**FetchAndOp)(PetscInt,PetscInt,const PetscInt*,void*,void*))
71795fce210SBarry Smith {
71895fce210SBarry Smith   PetscFunctionBegin;
71995fce210SBarry Smith   *FetchAndOp = NULL;
7208bfbc91cSJed Brown   if (op == MPIU_REPLACE) *FetchAndOp = link->FetchAndInsert;
72195fce210SBarry Smith   else if (op == MPI_SUM || op == MPIU_SUM) *FetchAndOp = link->FetchAndAdd;
72295fce210SBarry Smith   else if (op == MPI_MAX || op == MPIU_MAX) *FetchAndOp = link->FetchAndMax;
72395fce210SBarry Smith   else if (op == MPI_MIN || op == MPIU_MIN) *FetchAndOp = link->FetchAndMin;
72495fce210SBarry Smith   else if (op == MPI_MAXLOC) *FetchAndOp = link->FetchAndMaxloc;
72595fce210SBarry Smith   else if (op == MPI_MINLOC) *FetchAndOp = link->FetchAndMinloc;
72653deab39SPeter Brune   else if (op == MPI_PROD)   *FetchAndOp = link->FetchAndMult;
72753deab39SPeter Brune   else if (op == MPI_LAND)   *FetchAndOp = link->FetchAndLAND;
72853deab39SPeter Brune   else if (op == MPI_BAND)   *FetchAndOp = link->FetchAndBAND;
72953deab39SPeter Brune   else if (op == MPI_LOR)    *FetchAndOp = link->FetchAndLOR;
73053deab39SPeter Brune   else if (op == MPI_BOR)    *FetchAndOp = link->FetchAndBOR;
73153deab39SPeter Brune   else if (op == MPI_LXOR)   *FetchAndOp = link->FetchAndLXOR;
73253deab39SPeter Brune   else if (op == MPI_BXOR)   *FetchAndOp = link->FetchAndBXOR;
73395fce210SBarry Smith   else SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_SUP,"No support for MPI_Op");
73495fce210SBarry Smith   PetscFunctionReturn(0);
73595fce210SBarry Smith }
73695fce210SBarry Smith 
73795fce210SBarry Smith static PetscErrorCode PetscSFBasicPackGetReqs(PetscSF sf,PetscSFBasicPack link,MPI_Request **rootreqs,MPI_Request **leafreqs)
73895fce210SBarry Smith {
73995fce210SBarry Smith   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
74095fce210SBarry Smith 
74195fce210SBarry Smith   PetscFunctionBegin;
74295fce210SBarry Smith   if (rootreqs) *rootreqs = link->requests;
743c943f53fSJed Brown   if (leafreqs) *leafreqs = link->requests + (bas->niranks - bas->ndiranks);
74495fce210SBarry Smith   PetscFunctionReturn(0);
74595fce210SBarry Smith }
74695fce210SBarry Smith 
74795fce210SBarry Smith static PetscErrorCode PetscSFBasicPackWaitall(PetscSF sf,PetscSFBasicPack link)
74895fce210SBarry Smith {
74995fce210SBarry Smith   PetscSF_Basic  *bas = (PetscSF_Basic*)sf->data;
75095fce210SBarry Smith   PetscErrorCode ierr;
75195fce210SBarry Smith 
75295fce210SBarry Smith   PetscFunctionBegin;
753c943f53fSJed Brown   ierr = MPI_Waitall(bas->niranks+sf->nranks-(bas->ndiranks+sf->ndranks),link->requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
75495fce210SBarry Smith   PetscFunctionReturn(0);
75595fce210SBarry Smith }
75695fce210SBarry Smith 
757c943f53fSJed Brown static PetscErrorCode PetscSFBasicGetRootInfo(PetscSF sf,PetscInt *nrootranks,PetscInt *ndrootranks,const PetscMPIInt **rootranks,const PetscInt **rootoffset,const PetscInt **rootloc)
75895fce210SBarry Smith {
75995fce210SBarry Smith   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
76095fce210SBarry Smith 
76195fce210SBarry Smith   PetscFunctionBegin;
76295fce210SBarry Smith   if (nrootranks)  *nrootranks  = bas->niranks;
763c943f53fSJed Brown   if (ndrootranks) *ndrootranks = bas->ndiranks;
76495fce210SBarry Smith   if (rootranks)   *rootranks   = bas->iranks;
76595fce210SBarry Smith   if (rootoffset)  *rootoffset  = bas->ioffset;
76695fce210SBarry Smith   if (rootloc)     *rootloc     = bas->irootloc;
76795fce210SBarry Smith   PetscFunctionReturn(0);
76895fce210SBarry Smith }
76995fce210SBarry Smith 
770c943f53fSJed Brown static PetscErrorCode PetscSFBasicGetLeafInfo(PetscSF sf,PetscInt *nleafranks,PetscInt *ndleafranks,const PetscMPIInt **leafranks,const PetscInt **leafoffset,const PetscInt **leafloc)
77195fce210SBarry Smith {
77295fce210SBarry Smith   PetscFunctionBegin;
77395fce210SBarry Smith   if (nleafranks)  *nleafranks  = sf->nranks;
774c943f53fSJed Brown   if (ndleafranks) *ndleafranks = sf->ndranks;
77595fce210SBarry Smith   if (leafranks)   *leafranks   = sf->ranks;
77695fce210SBarry Smith   if (leafoffset)  *leafoffset  = sf->roffset;
77795fce210SBarry Smith   if (leafloc)     *leafloc     = sf->rmine;
77895fce210SBarry Smith   PetscFunctionReturn(0);
77995fce210SBarry Smith }
78095fce210SBarry Smith 
78195fce210SBarry Smith static PetscErrorCode PetscSFBasicGetPack(PetscSF sf,MPI_Datatype unit,const void *key,PetscSFBasicPack *mylink)
78295fce210SBarry Smith {
78395fce210SBarry Smith   PetscSF_Basic    *bas = (PetscSF_Basic*)sf->data;
78495fce210SBarry Smith   PetscErrorCode   ierr;
78595fce210SBarry Smith   PetscSFBasicPack link,*p;
786c943f53fSJed Brown   PetscInt         nrootranks,ndrootranks,nleafranks,ndleafranks,i;
78795fce210SBarry Smith   const PetscInt   *rootoffset,*leafoffset;
78895fce210SBarry Smith 
78995fce210SBarry Smith   PetscFunctionBegin;
79095fce210SBarry Smith   /* Look for types in cache */
79195fce210SBarry Smith   for (p=&bas->avail; (link=*p); p=&link->next) {
79295fce210SBarry Smith     PetscBool match;
79395fce210SBarry Smith     ierr = MPIPetsc_Type_compare(unit,link->unit,&match);CHKERRQ(ierr);
79495fce210SBarry Smith     if (match) {
79595fce210SBarry Smith       *p = link->next;          /* Remove from available list */
79695fce210SBarry Smith       goto found;
79795fce210SBarry Smith     }
79895fce210SBarry Smith   }
79995fce210SBarry Smith 
80095fce210SBarry Smith   /* Create new composite types for each send rank */
801c943f53fSJed Brown   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,&ndrootranks,NULL,&rootoffset,NULL);CHKERRQ(ierr);
802c943f53fSJed Brown   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,&ndleafranks,NULL,&leafoffset,NULL);CHKERRQ(ierr);
803b00a9115SJed Brown   ierr = PetscNew(&link);CHKERRQ(ierr);
80495fce210SBarry Smith   ierr = PetscSFBasicPackTypeSetup(link,unit);CHKERRQ(ierr);
805bf39f1bfSJed Brown   ierr = PetscMalloc2(nrootranks,&link->root,nleafranks,&link->leaf);CHKERRQ(ierr);
806bf39f1bfSJed Brown   for (i=0; i<nrootranks; i++) {
807bf39f1bfSJed Brown     ierr = PetscMalloc((rootoffset[i+1]-rootoffset[i])*link->unitbytes,&link->root[i]);CHKERRQ(ierr);
808bf39f1bfSJed Brown   }
809bf39f1bfSJed Brown   for (i=0; i<nleafranks; i++) {
810c943f53fSJed Brown     if (i < ndleafranks) {      /* Leaf buffers for distinguished ranks are pointers directly into root buffers */
811c943f53fSJed Brown       if (ndrootranks != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot match distinguished ranks");
812c943f53fSJed Brown       link->leaf[i] = link->root[0];
813c943f53fSJed Brown       continue;
814c943f53fSJed Brown     }
815bf39f1bfSJed Brown     ierr = PetscMalloc((leafoffset[i+1]-leafoffset[i])*link->unitbytes,&link->leaf[i]);CHKERRQ(ierr);
816bf39f1bfSJed Brown   }
8170f880758SHong Zhang   ierr = PetscCalloc1(nrootranks+nleafranks,&link->requests);CHKERRQ(ierr);
81895fce210SBarry Smith 
81995fce210SBarry Smith found:
82095fce210SBarry Smith   link->key  = key;
82195fce210SBarry Smith   link->next = bas->inuse;
82295fce210SBarry Smith   bas->inuse = link;
82395fce210SBarry Smith 
82495fce210SBarry Smith   *mylink = link;
82595fce210SBarry Smith   PetscFunctionReturn(0);
82695fce210SBarry Smith }
82795fce210SBarry Smith 
82895fce210SBarry Smith static PetscErrorCode PetscSFBasicGetPackInUse(PetscSF sf,MPI_Datatype unit,const void *key,PetscCopyMode cmode,PetscSFBasicPack *mylink)
82995fce210SBarry Smith {
83095fce210SBarry Smith   PetscSF_Basic    *bas = (PetscSF_Basic*)sf->data;
83195fce210SBarry Smith   PetscErrorCode   ierr;
83295fce210SBarry Smith   PetscSFBasicPack link,*p;
83395fce210SBarry Smith 
83495fce210SBarry Smith   PetscFunctionBegin;
83595fce210SBarry Smith   /* Look for types in cache */
83695fce210SBarry Smith   for (p=&bas->inuse; (link=*p); p=&link->next) {
83795fce210SBarry Smith     PetscBool match;
83895fce210SBarry Smith     ierr = MPIPetsc_Type_compare(unit,link->unit,&match);CHKERRQ(ierr);
83945038af3SLawrence Mitchell     if (match && (key == link->key)) {
84095fce210SBarry Smith       switch (cmode) {
84195fce210SBarry Smith       case PETSC_OWN_POINTER: *p = link->next; break; /* Remove from inuse list */
84295fce210SBarry Smith       case PETSC_USE_POINTER: break;
84395fce210SBarry Smith       default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"invalid cmode");
84495fce210SBarry Smith       }
84595fce210SBarry Smith       *mylink = link;
84695fce210SBarry Smith       PetscFunctionReturn(0);
84795fce210SBarry Smith     }
84895fce210SBarry Smith   }
84995fce210SBarry Smith   SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Could not find pack");
85095fce210SBarry Smith   PetscFunctionReturn(0);
85195fce210SBarry Smith }
85295fce210SBarry Smith 
85395fce210SBarry Smith static PetscErrorCode PetscSFBasicReclaimPack(PetscSF sf,PetscSFBasicPack *link)
85495fce210SBarry Smith {
85595fce210SBarry Smith   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
85695fce210SBarry Smith 
85795fce210SBarry Smith   PetscFunctionBegin;
85895fce210SBarry Smith   (*link)->key  = NULL;
85995fce210SBarry Smith   (*link)->next = bas->avail;
86095fce210SBarry Smith   bas->avail    = *link;
86195fce210SBarry Smith   *link         = NULL;
86295fce210SBarry Smith   PetscFunctionReturn(0);
86395fce210SBarry Smith }
86495fce210SBarry Smith 
8654416b707SBarry Smith static PetscErrorCode PetscSFSetFromOptions_Basic(PetscOptionItems *PetscOptionsObject,PetscSF sf)
86695fce210SBarry Smith {
86795fce210SBarry Smith   PetscErrorCode ierr;
86895fce210SBarry Smith 
86995fce210SBarry Smith   PetscFunctionBegin;
870e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"PetscSF Basic options");CHKERRQ(ierr);
87195fce210SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
87295fce210SBarry Smith   PetscFunctionReturn(0);
87395fce210SBarry Smith }
87495fce210SBarry Smith 
87595fce210SBarry Smith static PetscErrorCode PetscSFReset_Basic(PetscSF sf)
87695fce210SBarry Smith {
87795fce210SBarry Smith   PetscSF_Basic    *bas = (PetscSF_Basic*)sf->data;
87895fce210SBarry Smith   PetscErrorCode   ierr;
87995fce210SBarry Smith   PetscSFBasicPack link,next;
88095fce210SBarry Smith 
88195fce210SBarry Smith   PetscFunctionBegin;
88229046d53SLisandro Dalcin   if (bas->inuse) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Outstanding operation has not been completed");
883c943f53fSJed Brown   ierr = PetscFree2(bas->iranks,bas->ioffset);CHKERRQ(ierr);
884c943f53fSJed Brown   ierr = PetscFree(bas->irootloc);CHKERRQ(ierr);
88595fce210SBarry Smith   for (link=bas->avail; link; link=next) {
886bf39f1bfSJed Brown     PetscInt i;
88795fce210SBarry Smith     next = link->next;
888*77b7e48dSJunchao Zhang     if (!link->isbuiltin) {ierr = MPI_Type_free(&link->unit);CHKERRQ(ierr);}
889bf39f1bfSJed Brown     for (i=0; i<bas->niranks; i++) {ierr = PetscFree(link->root[i]);CHKERRQ(ierr);}
890c943f53fSJed Brown     for (i=sf->ndranks; i<sf->nranks; i++) {ierr = PetscFree(link->leaf[i]);CHKERRQ(ierr);} /* Free only non-distinguished leaf buffers */
89195fce210SBarry Smith     ierr = PetscFree2(link->root,link->leaf);CHKERRQ(ierr);
89295fce210SBarry Smith     ierr = PetscFree(link->requests);CHKERRQ(ierr);
89395fce210SBarry Smith     ierr = PetscFree(link);CHKERRQ(ierr);
89495fce210SBarry Smith   }
89595fce210SBarry Smith   bas->avail = NULL;
89695fce210SBarry Smith   PetscFunctionReturn(0);
89795fce210SBarry Smith }
89895fce210SBarry Smith 
89995fce210SBarry Smith static PetscErrorCode PetscSFDestroy_Basic(PetscSF sf)
90095fce210SBarry Smith {
90195fce210SBarry Smith   PetscErrorCode ierr;
90295fce210SBarry Smith 
90395fce210SBarry Smith   PetscFunctionBegin;
90429046d53SLisandro Dalcin   ierr = PetscSFReset_Basic(sf);CHKERRQ(ierr);
90595fce210SBarry Smith   ierr = PetscFree(sf->data);CHKERRQ(ierr);
90695fce210SBarry Smith   PetscFunctionReturn(0);
90795fce210SBarry Smith }
90895fce210SBarry Smith 
90995fce210SBarry Smith static PetscErrorCode PetscSFView_Basic(PetscSF sf,PetscViewer viewer)
91095fce210SBarry Smith {
91195fce210SBarry Smith   /* PetscSF_Basic *bas = (PetscSF_Basic*)sf->data; */
91295fce210SBarry Smith   PetscErrorCode ierr;
91395fce210SBarry Smith   PetscBool      iascii;
91495fce210SBarry Smith 
91595fce210SBarry Smith   PetscFunctionBegin;
91695fce210SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
91795fce210SBarry Smith   if (iascii) {
91895fce210SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  sort=%s\n",sf->rankorder ? "rank-order" : "unordered");CHKERRQ(ierr);
91995fce210SBarry Smith   }
92095fce210SBarry Smith   PetscFunctionReturn(0);
92195fce210SBarry Smith }
92295fce210SBarry Smith 
9233482bfa8SJunchao Zhang static PetscErrorCode PetscSFBcastAndOpBegin_Basic(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
92495fce210SBarry Smith {
92595fce210SBarry Smith   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;
92695fce210SBarry Smith   PetscErrorCode    ierr;
92795fce210SBarry Smith   PetscSFBasicPack  link;
928c943f53fSJed Brown   PetscInt          i,nrootranks,ndrootranks,nleafranks,ndleafranks;
92995fce210SBarry Smith   const PetscInt    *rootoffset,*leafoffset,*rootloc,*leafloc;
93095fce210SBarry Smith   const PetscMPIInt *rootranks,*leafranks;
93195fce210SBarry Smith   MPI_Request       *rootreqs,*leafreqs;
93295fce210SBarry Smith 
93395fce210SBarry Smith   PetscFunctionBegin;
934c943f53fSJed Brown   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,&ndrootranks,&rootranks,&rootoffset,&rootloc);CHKERRQ(ierr);
935c943f53fSJed Brown   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,&ndleafranks,&leafranks,&leafoffset,&leafloc);CHKERRQ(ierr);
93695fce210SBarry Smith   ierr = PetscSFBasicGetPack(sf,unit,rootdata,&link);CHKERRQ(ierr);
93795fce210SBarry Smith 
93895fce210SBarry Smith   ierr = PetscSFBasicPackGetReqs(sf,link,&rootreqs,&leafreqs);CHKERRQ(ierr);
939c943f53fSJed Brown   /* Eagerly post leaf receives, but only from non-distinguished ranks -- distinguished ranks will receive via shared memory */
940c943f53fSJed Brown   for (i=ndleafranks; i<nleafranks; i++) {
94195fce210SBarry Smith     PetscMPIInt n = leafoffset[i+1] - leafoffset[i];
942c943f53fSJed Brown     ierr = MPI_Irecv(link->leaf[i],n,unit,leafranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&leafreqs[i-ndleafranks]);CHKERRQ(ierr);
94395fce210SBarry Smith   }
94495fce210SBarry Smith   /* Pack and send root data */
94595fce210SBarry Smith   for (i=0; i<nrootranks; i++) {
94695fce210SBarry Smith     PetscMPIInt n          = rootoffset[i+1] - rootoffset[i];
947bf39f1bfSJed Brown     void        *packstart = link->root[i];
9482a74e5f4SJed Brown     (*link->Pack)(n,link->bs,rootloc+rootoffset[i],rootdata,packstart);
949c943f53fSJed Brown     if (i < ndrootranks) continue; /* shared memory */
950c943f53fSJed Brown     ierr = MPI_Isend(packstart,n,unit,rootranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&rootreqs[i-ndrootranks]);CHKERRQ(ierr);
95195fce210SBarry Smith   }
95295fce210SBarry Smith   PetscFunctionReturn(0);
95395fce210SBarry Smith }
95495fce210SBarry Smith 
9553482bfa8SJunchao Zhang PetscErrorCode PetscSFBcastAndOpEnd_Basic(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
95695fce210SBarry Smith {
95795fce210SBarry Smith   PetscErrorCode   ierr;
95895fce210SBarry Smith   PetscSFBasicPack link;
959c943f53fSJed Brown   PetscInt         i,nleafranks,ndleafranks;
96095fce210SBarry Smith   const PetscInt   *leafoffset,*leafloc;
9613482bfa8SJunchao Zhang   void             (*UnpackOp)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
9623482bfa8SJunchao Zhang   PetscMPIInt      typesize = -1;
96395fce210SBarry Smith 
96495fce210SBarry Smith   PetscFunctionBegin;
96595fce210SBarry Smith   ierr = PetscSFBasicGetPackInUse(sf,unit,rootdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr);
96695fce210SBarry Smith   ierr = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
967c943f53fSJed Brown   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,&ndleafranks,NULL,&leafoffset,&leafloc);CHKERRQ(ierr);
9683482bfa8SJunchao Zhang   ierr = PetscSFBasicPackGetUnpackOp(sf,link,op,&UnpackOp);CHKERRQ(ierr);
9693482bfa8SJunchao Zhang 
9703482bfa8SJunchao Zhang   if (UnpackOp) { typesize = link->unitbytes; }
9713482bfa8SJunchao Zhang   else { ierr = MPI_Type_size(unit,&typesize);CHKERRQ(ierr); }
9723482bfa8SJunchao Zhang 
97395fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
97495fce210SBarry Smith     PetscMPIInt n   = leafoffset[i+1] - leafoffset[i];
9753482bfa8SJunchao Zhang     char *packstart = (char *) link->leaf[i];
9763482bfa8SJunchao Zhang     if (UnpackOp) { (*UnpackOp)(n,link->bs,leafloc+leafoffset[i],leafdata,(const void *)packstart); }
9773482bfa8SJunchao Zhang #if defined(PETSC_HAVE_MPI_REDUCE_LOCAL)
9783482bfa8SJunchao Zhang     else if (n) { /* the op should be defined to operate on the whole datatype, so we ignore link->bs */
9793482bfa8SJunchao Zhang       PetscInt j;
9803482bfa8SJunchao Zhang       for (j=0; j<n; j++) { ierr = MPI_Reduce_local(packstart+j*typesize,((char *) leafdata)+(leafloc[leafoffset[i]+j])*typesize,1,unit,op);CHKERRQ(ierr); }
98195fce210SBarry Smith     }
9823482bfa8SJunchao Zhang #else
9833482bfa8SJunchao Zhang     else { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No unpacking reduction operation for this MPI_Op"); }
9843482bfa8SJunchao Zhang #endif
9853482bfa8SJunchao Zhang   }
9863482bfa8SJunchao Zhang 
98795fce210SBarry Smith   ierr = PetscSFBasicReclaimPack(sf,&link);CHKERRQ(ierr);
98895fce210SBarry Smith   PetscFunctionReturn(0);
98995fce210SBarry Smith }
99095fce210SBarry Smith 
9913482bfa8SJunchao Zhang /* Send from roots to leaves */
9923482bfa8SJunchao Zhang static PetscErrorCode PetscSFBcastBegin_Basic(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
9933482bfa8SJunchao Zhang {
9943482bfa8SJunchao Zhang   PetscErrorCode   ierr;
9953482bfa8SJunchao Zhang 
9963482bfa8SJunchao Zhang   PetscFunctionBegin;
9973482bfa8SJunchao Zhang   ierr = PetscSFBcastAndOpBegin_Basic(sf,unit,rootdata,leafdata,MPI_REPLACE);CHKERRQ(ierr);
9983482bfa8SJunchao Zhang   PetscFunctionReturn(0);
9993482bfa8SJunchao Zhang }
10003482bfa8SJunchao Zhang 
10013482bfa8SJunchao Zhang PetscErrorCode PetscSFBcastEnd_Basic(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
10023482bfa8SJunchao Zhang {
10033482bfa8SJunchao Zhang   PetscErrorCode   ierr;
10043482bfa8SJunchao Zhang 
10053482bfa8SJunchao Zhang   PetscFunctionBegin;
10063482bfa8SJunchao Zhang   ierr = PetscSFBcastAndOpEnd_Basic(sf,unit,rootdata,leafdata,MPI_REPLACE);CHKERRQ(ierr);
10073482bfa8SJunchao Zhang   PetscFunctionReturn(0);
10083482bfa8SJunchao Zhang }
10093482bfa8SJunchao Zhang 
101095fce210SBarry Smith /* leaf -> root with reduction */
101195fce210SBarry Smith PetscErrorCode PetscSFReduceBegin_Basic(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
101295fce210SBarry Smith {
101395fce210SBarry Smith   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;
101495fce210SBarry Smith   PetscSFBasicPack  link;
101595fce210SBarry Smith   PetscErrorCode    ierr;
1016c943f53fSJed Brown   PetscInt          i,nrootranks,ndrootranks,nleafranks,ndleafranks;
101795fce210SBarry Smith   const PetscInt    *rootoffset,*leafoffset,*rootloc,*leafloc;
101895fce210SBarry Smith   const PetscMPIInt *rootranks,*leafranks;
101995fce210SBarry Smith   MPI_Request       *rootreqs,*leafreqs;
102095fce210SBarry Smith 
102195fce210SBarry Smith   PetscFunctionBegin;
1022c943f53fSJed Brown   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,&ndrootranks,&rootranks,&rootoffset,&rootloc);CHKERRQ(ierr);
1023c943f53fSJed Brown   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,&ndleafranks,&leafranks,&leafoffset,&leafloc);CHKERRQ(ierr);
102495fce210SBarry Smith   ierr = PetscSFBasicGetPack(sf,unit,rootdata,&link);CHKERRQ(ierr);
102595fce210SBarry Smith 
102695fce210SBarry Smith   ierr = PetscSFBasicPackGetReqs(sf,link,&rootreqs,&leafreqs);CHKERRQ(ierr);
1027c943f53fSJed Brown   /* Eagerly post root receives for non-distinguished ranks */
1028c943f53fSJed Brown   for (i=ndrootranks; i<nrootranks; i++) {
102995fce210SBarry Smith     PetscMPIInt n = rootoffset[i+1] - rootoffset[i];
1030c943f53fSJed Brown     ierr = MPI_Irecv(link->root[i],n,unit,rootranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&rootreqs[i-ndrootranks]);CHKERRQ(ierr);
103195fce210SBarry Smith   }
103295fce210SBarry Smith   /* Pack and send leaf data */
103395fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
103495fce210SBarry Smith     PetscMPIInt n          = leafoffset[i+1] - leafoffset[i];
1035bf39f1bfSJed Brown     void        *packstart = link->leaf[i];
10362a74e5f4SJed Brown     (*link->Pack)(n,link->bs,leafloc+leafoffset[i],leafdata,packstart);
1037c943f53fSJed Brown     if (i < ndleafranks) continue; /* shared memory */
1038c943f53fSJed Brown     ierr = MPI_Isend(packstart,n,unit,leafranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&leafreqs[i-ndleafranks]);CHKERRQ(ierr);
103995fce210SBarry Smith   }
104095fce210SBarry Smith   PetscFunctionReturn(0);
104195fce210SBarry Smith }
104295fce210SBarry Smith 
104395fce210SBarry Smith static PetscErrorCode PetscSFReduceEnd_Basic(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
104495fce210SBarry Smith {
10452a74e5f4SJed Brown   void             (*UnpackOp)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
104695fce210SBarry Smith   PetscErrorCode   ierr;
104795fce210SBarry Smith   PetscSFBasicPack link;
104895fce210SBarry Smith   PetscInt         i,nrootranks;
10491620da3bSToby Isaac   PetscMPIInt      typesize = -1;
105095fce210SBarry Smith   const PetscInt   *rootoffset,*rootloc;
105195fce210SBarry Smith 
105295fce210SBarry Smith   PetscFunctionBegin;
105395fce210SBarry Smith   ierr = PetscSFBasicGetPackInUse(sf,unit,rootdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr);
105495fce210SBarry Smith   /* This implementation could be changed to unpack as receives arrive, at the cost of non-determinism */
105595fce210SBarry Smith   ierr = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
1056c943f53fSJed Brown   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,NULL,NULL,&rootoffset,&rootloc);CHKERRQ(ierr);
105795fce210SBarry Smith   ierr = PetscSFBasicPackGetUnpackOp(sf,link,op,&UnpackOp);CHKERRQ(ierr);
10581620da3bSToby Isaac   if (UnpackOp) {
10591620da3bSToby Isaac     typesize = link->unitbytes;
10601620da3bSToby Isaac   }
10611620da3bSToby Isaac   else {
10621620da3bSToby Isaac     ierr = MPI_Type_size(unit,&typesize);CHKERRQ(ierr);
10631620da3bSToby Isaac   }
106495fce210SBarry Smith   for (i=0; i<nrootranks; i++) {
106595fce210SBarry Smith     PetscMPIInt n   = rootoffset[i+1] - rootoffset[i];
1066bf39f1bfSJed Brown     char *packstart = (char *) link->root[i];
106795fce210SBarry Smith 
10681620da3bSToby Isaac     if (UnpackOp) {
106988ac5ce8SToby Isaac       (*UnpackOp)(n,link->bs,rootloc+rootoffset[i],rootdata,(const void *)packstart);
107095fce210SBarry Smith     }
107162795ce3SStefano Zampini #if defined(PETSC_HAVE_MPI_REDUCE_LOCAL)
107278ec4774SToby Isaac     else if (n) { /* the op should be defined to operate on the whole datatype, so we ignore link->bs */
107378ec4774SToby Isaac       PetscInt j;
10741620da3bSToby Isaac 
10751620da3bSToby Isaac       for (j = 0; j < n; j++) {
107678ec4774SToby Isaac         ierr = MPI_Reduce_local(packstart+j*typesize,((char *) rootdata)+(rootloc[rootoffset[i]+j])*typesize,1,unit,op);CHKERRQ(ierr);
10771620da3bSToby Isaac       }
10781620da3bSToby Isaac     }
1079a9ec7e9eSToby Isaac #else
1080a9ec7e9eSToby Isaac     else {
1081a9ec7e9eSToby Isaac       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No unpacking reduction operation for this MPI_Op");
1082a9ec7e9eSToby Isaac     }
1083a9ec7e9eSToby Isaac #endif
10841620da3bSToby Isaac   }
108595fce210SBarry Smith   ierr = PetscSFBasicReclaimPack(sf,&link);CHKERRQ(ierr);
108695fce210SBarry Smith   PetscFunctionReturn(0);
108795fce210SBarry Smith }
108895fce210SBarry Smith 
108995fce210SBarry Smith static PetscErrorCode PetscSFFetchAndOpBegin_Basic(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
109095fce210SBarry Smith {
109195fce210SBarry Smith   PetscErrorCode ierr;
109295fce210SBarry Smith 
109395fce210SBarry Smith   PetscFunctionBegin;
109495fce210SBarry Smith   ierr = PetscSFReduceBegin_Basic(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
109595fce210SBarry Smith   PetscFunctionReturn(0);
109695fce210SBarry Smith }
109795fce210SBarry Smith 
109895fce210SBarry Smith static PetscErrorCode PetscSFFetchAndOpEnd_Basic(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
109995fce210SBarry Smith {
110095fce210SBarry Smith   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;
11012a74e5f4SJed Brown   void              (*FetchAndOp)(PetscInt,PetscInt,const PetscInt*,void*,void*);
110295fce210SBarry Smith   PetscErrorCode    ierr;
110395fce210SBarry Smith   PetscSFBasicPack  link;
1104c943f53fSJed Brown   PetscInt          i,nrootranks,ndrootranks,nleafranks,ndleafranks;
110595fce210SBarry Smith   const PetscInt    *rootoffset,*leafoffset,*rootloc,*leafloc;
110695fce210SBarry Smith   const PetscMPIInt *rootranks,*leafranks;
110795fce210SBarry Smith   MPI_Request       *rootreqs,*leafreqs;
110895fce210SBarry Smith 
110995fce210SBarry Smith   PetscFunctionBegin;
111095fce210SBarry Smith   ierr = PetscSFBasicGetPackInUse(sf,unit,rootdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr);
111195fce210SBarry Smith   /* This implementation could be changed to unpack as receives arrive, at the cost of non-determinism */
111295fce210SBarry Smith   ierr      = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
1113c943f53fSJed Brown   ierr      = PetscSFBasicGetRootInfo(sf,&nrootranks,&ndrootranks,&rootranks,&rootoffset,&rootloc);CHKERRQ(ierr);
1114c943f53fSJed Brown   ierr      = PetscSFBasicGetLeafInfo(sf,&nleafranks,&ndleafranks,&leafranks,&leafoffset,&leafloc);CHKERRQ(ierr);
111595fce210SBarry Smith   ierr      = PetscSFBasicPackGetReqs(sf,link,&rootreqs,&leafreqs);CHKERRQ(ierr);
111695fce210SBarry Smith   /* Post leaf receives */
1117c943f53fSJed Brown   for (i=ndleafranks; i<nleafranks; i++) {
111895fce210SBarry Smith     PetscMPIInt n = leafoffset[i+1] - leafoffset[i];
1119c943f53fSJed Brown     ierr = MPI_Irecv(link->leaf[i],n,unit,leafranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&leafreqs[i-ndleafranks]);CHKERRQ(ierr);
112095fce210SBarry Smith   }
112195fce210SBarry Smith   /* Process local fetch-and-op, post root sends */
112295fce210SBarry Smith   ierr = PetscSFBasicPackGetFetchAndOp(sf,link,op,&FetchAndOp);CHKERRQ(ierr);
112395fce210SBarry Smith   for (i=0; i<nrootranks; i++) {
112495fce210SBarry Smith     PetscMPIInt n          = rootoffset[i+1] - rootoffset[i];
1125bf39f1bfSJed Brown     void        *packstart = link->root[i];
112695fce210SBarry Smith 
11272a74e5f4SJed Brown     (*FetchAndOp)(n,link->bs,rootloc+rootoffset[i],rootdata,packstart);
1128c943f53fSJed Brown     if (i < ndrootranks) continue; /* shared memory */
1129c943f53fSJed Brown     ierr = MPI_Isend(packstart,n,unit,rootranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&rootreqs[i-ndrootranks]);CHKERRQ(ierr);
113095fce210SBarry Smith   }
113195fce210SBarry Smith   ierr = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
113295fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
113395fce210SBarry Smith     PetscMPIInt n          = leafoffset[i+1] - leafoffset[i];
1134bf39f1bfSJed Brown     const void  *packstart = link->leaf[i];
11352a74e5f4SJed Brown     (*link->UnpackInsert)(n,link->bs,leafloc+leafoffset[i],leafupdate,packstart);
113695fce210SBarry Smith   }
113795fce210SBarry Smith   ierr = PetscSFBasicReclaimPack(sf,&link);CHKERRQ(ierr);
113895fce210SBarry Smith   PetscFunctionReturn(0);
113995fce210SBarry Smith }
114095fce210SBarry Smith 
11418750ddebSJunchao Zhang static PetscErrorCode PetscSFGetLeafRanks_Basic(PetscSF sf,PetscInt *niranks,const PetscMPIInt **iranks,const PetscInt **ioffset,const PetscInt **irootloc)
11428750ddebSJunchao Zhang {
11438750ddebSJunchao Zhang   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
11448750ddebSJunchao Zhang 
11458750ddebSJunchao Zhang   PetscFunctionBegin;
11468750ddebSJunchao Zhang   if (niranks)  *niranks  = bas->niranks;
11478750ddebSJunchao Zhang   if (iranks)   *iranks   = bas->iranks;
11488750ddebSJunchao Zhang   if (ioffset)  *ioffset  = bas->ioffset;
11498750ddebSJunchao Zhang   if (irootloc) *irootloc = bas->irootloc;
11508750ddebSJunchao Zhang   PetscFunctionReturn(0);
11518750ddebSJunchao Zhang }
11528750ddebSJunchao Zhang 
11538cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscSFCreate_Basic(PetscSF sf)
115495fce210SBarry Smith {
115595fce210SBarry Smith   PetscSF_Basic  *bas = (PetscSF_Basic*)sf->data;
115695fce210SBarry Smith   PetscErrorCode ierr;
115795fce210SBarry Smith 
115895fce210SBarry Smith   PetscFunctionBegin;
115995fce210SBarry Smith   sf->ops->SetUp           = PetscSFSetUp_Basic;
116095fce210SBarry Smith   sf->ops->SetFromOptions  = PetscSFSetFromOptions_Basic;
116195fce210SBarry Smith   sf->ops->Reset           = PetscSFReset_Basic;
116295fce210SBarry Smith   sf->ops->Destroy         = PetscSFDestroy_Basic;
116395fce210SBarry Smith   sf->ops->View            = PetscSFView_Basic;
116495fce210SBarry Smith   sf->ops->BcastBegin      = PetscSFBcastBegin_Basic;
116595fce210SBarry Smith   sf->ops->BcastEnd        = PetscSFBcastEnd_Basic;
11663482bfa8SJunchao Zhang   sf->ops->BcastAndOpBegin = PetscSFBcastAndOpBegin_Basic;
11673482bfa8SJunchao Zhang   sf->ops->BcastAndOpEnd   = PetscSFBcastAndOpEnd_Basic;
116895fce210SBarry Smith   sf->ops->ReduceBegin     = PetscSFReduceBegin_Basic;
116995fce210SBarry Smith   sf->ops->ReduceEnd       = PetscSFReduceEnd_Basic;
117095fce210SBarry Smith   sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Basic;
117195fce210SBarry Smith   sf->ops->FetchAndOpEnd   = PetscSFFetchAndOpEnd_Basic;
11728750ddebSJunchao Zhang   sf->ops->GetLeafRanks    = PetscSFGetLeafRanks_Basic;
117395fce210SBarry Smith 
1174b00a9115SJed Brown   ierr = PetscNewLog(sf,&bas);CHKERRQ(ierr);
117595fce210SBarry Smith   sf->data = (void*)bas;
117695fce210SBarry Smith   PetscFunctionReturn(0);
117795fce210SBarry Smith }
1178