xref: /petsc/src/vec/is/sf/impls/basic/sfbasic.c (revision c943f53f8bfa9306c8e67b2fced8d88cb4128a92)
18cd53115SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/
395fce210SBarry Smith 
495fce210SBarry Smith typedef struct _n_PetscSFBasicPack *PetscSFBasicPack;
595fce210SBarry Smith struct _n_PetscSFBasicPack {
62a74e5f4SJed Brown   void (*Pack)(PetscInt,PetscInt,const PetscInt*,const void*,void*);
72a74e5f4SJed Brown   void (*UnpackInsert)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
82a74e5f4SJed Brown   void (*UnpackAdd)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
92a74e5f4SJed Brown   void (*UnpackMin)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
102a74e5f4SJed Brown   void (*UnpackMax)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
112a74e5f4SJed Brown   void (*UnpackMinloc)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
122a74e5f4SJed Brown   void (*UnpackMaxloc)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
132a74e5f4SJed Brown   void (*UnpackMult)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
142a74e5f4SJed Brown   void (*UnpackLAND)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
152a74e5f4SJed Brown   void (*UnpackBAND)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
162a74e5f4SJed Brown   void (*UnpackLOR)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
172a74e5f4SJed Brown   void (*UnpackBOR)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
182a74e5f4SJed Brown   void (*UnpackLXOR)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
192a74e5f4SJed Brown   void (*UnpackBXOR)(PetscInt,PetscInt,const PetscInt*,void*,const void *);
202a74e5f4SJed Brown   void (*FetchAndInsert)(PetscInt,PetscInt,const PetscInt*,void*,void*);
212a74e5f4SJed Brown   void (*FetchAndAdd)(PetscInt,PetscInt,const PetscInt*,void*,void*);
222a74e5f4SJed Brown   void (*FetchAndMin)(PetscInt,PetscInt,const PetscInt*,void*,void*);
232a74e5f4SJed Brown   void (*FetchAndMax)(PetscInt,PetscInt,const PetscInt*,void*,void*);
242a74e5f4SJed Brown   void (*FetchAndMinloc)(PetscInt,PetscInt,const PetscInt*,void*,void*);
252a74e5f4SJed Brown   void (*FetchAndMaxloc)(PetscInt,PetscInt,const PetscInt*,void*,void*);
262a74e5f4SJed Brown   void (*FetchAndMult)(PetscInt,PetscInt,const PetscInt*,void*,void*);
272a74e5f4SJed Brown   void (*FetchAndLAND)(PetscInt,PetscInt,const PetscInt*,void*,void*);
282a74e5f4SJed Brown   void (*FetchAndBAND)(PetscInt,PetscInt,const PetscInt*,void*,void*);
292a74e5f4SJed Brown   void (*FetchAndLOR)(PetscInt,PetscInt,const PetscInt*,void*,void*);
302a74e5f4SJed Brown   void (*FetchAndBOR)(PetscInt,PetscInt,const PetscInt*,void*,void*);
312a74e5f4SJed Brown   void (*FetchAndLXOR)(PetscInt,PetscInt,const PetscInt*,void*,void*);
322a74e5f4SJed Brown   void (*FetchAndBXOR)(PetscInt,PetscInt,const PetscInt*,void*,void*);
3395fce210SBarry Smith 
3495fce210SBarry Smith   MPI_Datatype     unit;
3595fce210SBarry Smith   size_t           unitbytes;   /* Number of bytes in a unit */
362a74e5f4SJed Brown   PetscInt         bs;          /* Number of basic units in a unit */
3795fce210SBarry Smith   const void       *key;        /* Array used as key for operation */
38bf39f1bfSJed Brown   char             **root;      /* Packed root data, indexed by leaf rank */
39bf39f1bfSJed Brown   char             **leaf;      /* Packed leaf data, indexed by root rank */
4095fce210SBarry Smith   MPI_Request      *requests;   /* Array of root requests followed by leaf requests */
4195fce210SBarry Smith   PetscSFBasicPack next;
4295fce210SBarry Smith };
4395fce210SBarry Smith 
4495fce210SBarry Smith typedef struct {
4595fce210SBarry Smith   PetscMPIInt      tag;
46cf4b5b4fSJed Brown   PetscMPIInt      niranks;     /* Number of incoming ranks (ranks accessing my roots) */
47*c943f53fSJed Brown   PetscMPIInt      ndiranks;    /* Number of incoming ranks (ranks accessing my roots) in distinguished set */
4895fce210SBarry Smith   PetscMPIInt      *iranks;     /* Array of ranks that reference my roots */
4995fce210SBarry Smith   PetscInt         itotal;      /* Total number of graph edges referencing my roots */
5095fce210SBarry Smith   PetscInt         *ioffset;    /* Array of length niranks+1 holding offset in irootloc[] for each rank */
5195fce210SBarry Smith   PetscInt         *irootloc;   /* Incoming roots referenced by ranks starting at ioffset[rank] */
5295fce210SBarry Smith   PetscSFBasicPack avail;       /* One or more entries per MPI Datatype, lazily constructed */
5395fce210SBarry Smith   PetscSFBasicPack inuse;       /* Buffers being used for transactions that have not yet completed */
5495fce210SBarry Smith } PetscSF_Basic;
5595fce210SBarry Smith 
5695fce210SBarry Smith #if !defined(PETSC_HAVE_MPI_TYPE_DUP) /* Danger: type is not reference counted; subject to ABA problem */
5795fce210SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MPI_Type_dup(MPI_Datatype datatype,MPI_Datatype *newtype)
5895fce210SBarry Smith {
5995fce210SBarry Smith   *newtype = datatype;
6095fce210SBarry Smith   return 0;
6195fce210SBarry Smith }
6295fce210SBarry Smith #endif
6395fce210SBarry Smith 
6495fce210SBarry Smith /*
6595fce210SBarry Smith  * MPI_Reduce_local is not really useful because it can't handle sparse data and it vectorizes "in the wrong direction",
6695fce210SBarry Smith  * therefore we pack data types manually. This section defines packing routines for the standard data types.
6795fce210SBarry Smith  */
6895fce210SBarry Smith 
6995fce210SBarry Smith #define CPPJoin2_exp(a,b) a ## b
7095fce210SBarry Smith #define CPPJoin2(a,b) CPPJoin2_exp(a,b)
7195fce210SBarry Smith #define CPPJoin3_exp_(a,b,c) a ## b ## _ ## c
7295fce210SBarry Smith #define CPPJoin3_(a,b,c) CPPJoin3_exp_(a,b,c)
7395fce210SBarry Smith 
7495fce210SBarry Smith /* Basic types without addition */
752a74e5f4SJed Brown #define DEF_PackNoInit(type,BS)                                         \
762a74e5f4SJed Brown   static void CPPJoin3_(Pack_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,const void *unpacked,void *packed) { \
7795fce210SBarry Smith     const type *u = (const type*)unpacked;                              \
7895fce210SBarry Smith     type *p = (type*)packed;                                            \
792a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
802a74e5f4SJed Brown     for (i=0; i<n; i++)                                                 \
812a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
822a74e5f4SJed Brown         for (k=j; k<j+BS; k++)                                          \
832a74e5f4SJed Brown           p[i*bs+k] = u[idx[i]*bs+k];                                   \
8495fce210SBarry Smith   }                                                                     \
852a74e5f4SJed Brown   static void CPPJoin3_(UnpackInsert_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
8695fce210SBarry Smith     type *u = (type*)unpacked;                                          \
8795fce210SBarry Smith     const type *p = (const type*)packed;                                \
882a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
892a74e5f4SJed Brown     for (i=0; i<n; i++)                                                 \
902a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
912a74e5f4SJed Brown         for (k=j; k<j+BS; k++)                                          \
922a74e5f4SJed Brown           u[idx[i]*bs+k] = p[i*bs+k];                                   \
9395fce210SBarry Smith   }                                                                     \
942a74e5f4SJed Brown   static void CPPJoin3_(FetchAndInsert_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
9595fce210SBarry Smith     type *u = (type*)unpacked;                                          \
9695fce210SBarry Smith     type *p = (type*)packed;                                            \
972a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
9895fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
992a74e5f4SJed Brown       PetscInt ii = idx[i];                                             \
1002a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
1012a74e5f4SJed Brown         for (k=j; k<j+BS; k++) {                                        \
1022a74e5f4SJed Brown           type t = u[ii*bs+k];                                          \
1032a74e5f4SJed Brown           u[ii*bs+k] = p[i*bs+k];                                       \
1042a74e5f4SJed Brown           p[i*bs+k] = t;                                                \
1052a74e5f4SJed Brown         }                                                               \
10695fce210SBarry Smith     }                                                                   \
10795fce210SBarry Smith   }
10895fce210SBarry Smith 
10995fce210SBarry Smith /* Basic types defining addition */
1102a74e5f4SJed Brown #define DEF_PackAddNoInit(type,BS)                                      \
1112a74e5f4SJed Brown   DEF_PackNoInit(type,BS)                                               \
1122a74e5f4SJed Brown   static void CPPJoin3_(UnpackAdd_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
11395fce210SBarry Smith     type *u = (type*)unpacked;                                          \
11495fce210SBarry Smith     const type *p = (const type*)packed;                                \
1152a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
1162a74e5f4SJed Brown     for (i=0; i<n; i++)                                                 \
1172a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
1182a74e5f4SJed Brown         for (k=j; k<j+BS; k++)                                          \
1192a74e5f4SJed Brown           u[idx[i]*bs+k] += p[i*bs+k];                                  \
12095fce210SBarry Smith   }                                                                     \
1212a74e5f4SJed Brown   static void CPPJoin3_(FetchAndAdd_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
12295fce210SBarry Smith     type *u = (type*)unpacked;                                          \
12395fce210SBarry Smith     type *p = (type*)packed;                                            \
1242a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
12595fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
1262a74e5f4SJed Brown       PetscInt ii = idx[i];                                             \
1272a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
1282a74e5f4SJed Brown         for (k=j; k<j+BS; k++) {                                        \
1292a74e5f4SJed Brown           type t = u[ii*bs+k];                                          \
1302a74e5f4SJed Brown           u[ii*bs+k] = t + p[i*bs+k];                                   \
1312a74e5f4SJed Brown           p[i*bs+k] = t;                                                \
13295fce210SBarry Smith         }                                                               \
13353deab39SPeter Brune     }                                                                   \
1342a74e5f4SJed Brown   }                                                                     \
1352a74e5f4SJed Brown   static void CPPJoin3_(UnpackMult_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
13653deab39SPeter Brune     type *u = (type*)unpacked;                                          \
13753deab39SPeter Brune     const type *p = (const type*)packed;                                \
1382a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
1392a74e5f4SJed Brown     for (i=0; i<n; i++)                                                 \
1402a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
1412a74e5f4SJed Brown         for (k=j; k<j+BS; k++)                                          \
1422a74e5f4SJed Brown           u[idx[i]*bs+k] *= p[i*bs+k];                                  \
14353deab39SPeter Brune   }                                                                     \
1442a74e5f4SJed Brown   static void CPPJoin3_(FetchAndMult_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
14553deab39SPeter Brune     type *u = (type*)unpacked;                                          \
14653deab39SPeter Brune     type *p = (type*)packed;                                            \
1472a74e5f4SJed Brown     PetscInt i,j,k;                                                     \
14853deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
1492a74e5f4SJed Brown       PetscInt ii = idx[i];                                             \
1502a74e5f4SJed Brown       for (j=0; j<bs; j+=BS)                                            \
1512a74e5f4SJed Brown         for (k=j; k<j+BS; k++) {                                        \
1522a74e5f4SJed Brown           type t = u[ii*bs+k];                                          \
1532a74e5f4SJed Brown           u[ii*bs+k] = t * p[i*bs+k];                                   \
1542a74e5f4SJed Brown           p[i*bs+k] = t;                                                \
1552a74e5f4SJed Brown         }                                                               \
15653deab39SPeter Brune     }                                                                   \
15795fce210SBarry Smith   }
1582a74e5f4SJed Brown #define DEF_Pack(type,BS)                                               \
1592a74e5f4SJed Brown   DEF_PackAddNoInit(type,BS)                                            \
1602a74e5f4SJed Brown   static void CPPJoin3_(PackInit_,type,BS)(PetscSFBasicPack link) {     \
1612a74e5f4SJed Brown     link->Pack = CPPJoin3_(Pack_,type,BS);                              \
1622a74e5f4SJed Brown     link->UnpackInsert = CPPJoin3_(UnpackInsert_,type,BS);              \
1632a74e5f4SJed Brown     link->UnpackAdd = CPPJoin3_(UnpackAdd_,type,BS);                    \
1642a74e5f4SJed Brown     link->UnpackMult = CPPJoin3_(UnpackMult_,type,BS);                  \
1652a74e5f4SJed Brown     link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,type,BS);          \
1662a74e5f4SJed Brown     link->FetchAndAdd = CPPJoin3_(FetchAndAdd_,type,BS);                \
1672a74e5f4SJed Brown     link->FetchAndMult = CPPJoin3_(FetchAndMult_,type,BS);              \
16895fce210SBarry Smith     link->unitbytes = sizeof(type);                                     \
16995fce210SBarry Smith   }
17095fce210SBarry Smith /* Comparable types */
17195fce210SBarry Smith #define DEF_PackCmp(type)                                               \
1722a74e5f4SJed Brown   DEF_PackAddNoInit(type,1)                                             \
1732a74e5f4SJed Brown   static void CPPJoin2(UnpackMax_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
17495fce210SBarry Smith     type *u = (type*)unpacked;                                          \
17595fce210SBarry Smith     const type *p = (const type*)packed;                                \
17695fce210SBarry Smith     PetscInt i;                                                         \
17795fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
17895fce210SBarry Smith       type v = u[idx[i]];                                               \
17995fce210SBarry Smith       u[idx[i]] = PetscMax(v,p[i]);                                     \
18095fce210SBarry Smith     }                                                                   \
18195fce210SBarry Smith   }                                                                     \
1822a74e5f4SJed Brown   static void CPPJoin2(UnpackMin_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
18395fce210SBarry Smith     type *u = (type*)unpacked;                                          \
18495fce210SBarry Smith     const type *p = (const type*)packed;                                \
18595fce210SBarry Smith     PetscInt i;                                                         \
18695fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
18795fce210SBarry Smith       type v = u[idx[i]];                                               \
18895fce210SBarry Smith       u[idx[i]] = PetscMin(v,p[i]);                                     \
18995fce210SBarry Smith     }                                                                   \
19095fce210SBarry Smith   }                                                                     \
1912a74e5f4SJed Brown   static void CPPJoin2(FetchAndMax_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
19295fce210SBarry Smith     type *u = (type*)unpacked;                                          \
19395fce210SBarry Smith     type *p = (type*)packed;                                            \
19495fce210SBarry Smith     PetscInt i;                                                         \
19595fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
19695fce210SBarry Smith       PetscInt j = idx[i];                                              \
19795fce210SBarry Smith       type v = u[j];                                                    \
19895fce210SBarry Smith       u[j] = PetscMax(v,p[i]);                                          \
19995fce210SBarry Smith       p[i] = v;                                                         \
20095fce210SBarry Smith     }                                                                   \
20195fce210SBarry Smith   }                                                                     \
2022a74e5f4SJed Brown   static void CPPJoin2(FetchAndMin_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
20395fce210SBarry Smith     type *u = (type*)unpacked;                                          \
20495fce210SBarry Smith     type *p = (type*)packed;                                            \
20595fce210SBarry Smith     PetscInt i;                                                         \
20695fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
20795fce210SBarry Smith       PetscInt j = idx[i];                                              \
20895fce210SBarry Smith       type v = u[j];                                                    \
20995fce210SBarry Smith       u[j] = PetscMin(v,p[i]);                                          \
21095fce210SBarry Smith       p[i] = v;                                                         \
21195fce210SBarry Smith     }                                                                   \
21295fce210SBarry Smith   }                                                                     \
21395fce210SBarry Smith   static void CPPJoin2(PackInit_,type)(PetscSFBasicPack link) {         \
2142a74e5f4SJed Brown     link->Pack = CPPJoin3_(Pack_,type,1);                               \
2152a74e5f4SJed Brown     link->UnpackInsert = CPPJoin3_(UnpackInsert_,type,1);               \
2162a74e5f4SJed Brown     link->UnpackAdd  = CPPJoin3_(UnpackAdd_,type,1);                    \
21795fce210SBarry Smith     link->UnpackMax  = CPPJoin2(UnpackMax_,type);                       \
21895fce210SBarry Smith     link->UnpackMin  = CPPJoin2(UnpackMin_,type);                       \
2192a74e5f4SJed Brown     link->UnpackMult = CPPJoin3_(UnpackMult_,type,1);                   \
2202a74e5f4SJed Brown     link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,type,1);           \
2212a74e5f4SJed Brown     link->FetchAndAdd = CPPJoin3_(FetchAndAdd_ ,type,1);                \
22295fce210SBarry Smith     link->FetchAndMax = CPPJoin2(FetchAndMax_ ,type);                   \
22395fce210SBarry Smith     link->FetchAndMin = CPPJoin2(FetchAndMin_ ,type);                   \
2242a74e5f4SJed Brown     link->FetchAndMult = CPPJoin3_(FetchAndMult_,type,1);               \
22595fce210SBarry Smith     link->unitbytes = sizeof(type);                                     \
22695fce210SBarry Smith   }
22795fce210SBarry Smith 
22853deab39SPeter Brune /* Logical Types */
22953deab39SPeter Brune #define DEF_PackLog(type)                                               \
2302a74e5f4SJed Brown   static void CPPJoin2(UnpackLAND_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
23153deab39SPeter Brune     type *u = (type*)unpacked;                                          \
23253deab39SPeter Brune     const type *p = (const type*)packed;                                \
23353deab39SPeter Brune     PetscInt i;                                                         \
23453deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
23553deab39SPeter Brune       type v = u[idx[i]];                                               \
23653deab39SPeter Brune       u[idx[i]] = v && p[i];                                            \
23753deab39SPeter Brune     }                                                                   \
23853deab39SPeter Brune   }                                                                     \
2392a74e5f4SJed Brown   static void CPPJoin2(UnpackLOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
24053deab39SPeter Brune     type *u = (type*)unpacked;                                          \
24153deab39SPeter Brune     const type *p = (const type*)packed;                                \
24253deab39SPeter Brune     PetscInt i;                                                         \
24353deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
24453deab39SPeter Brune       type v = u[idx[i]];                                               \
24553deab39SPeter Brune       u[idx[i]] = v || p[i];                                            \
24653deab39SPeter Brune     }                                                                   \
24753deab39SPeter Brune   }                                                                     \
2482a74e5f4SJed Brown   static void CPPJoin2(UnpackLXOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
24953deab39SPeter Brune     type *u = (type*)unpacked;                                          \
25053deab39SPeter Brune     const type *p = (const type*)packed;                                \
25153deab39SPeter Brune     PetscInt i;                                                         \
25253deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
25353deab39SPeter Brune       type v = u[idx[i]];                                               \
25453deab39SPeter Brune       u[idx[i]] = (!v)!=(!p[i]);                                        \
25553deab39SPeter Brune     }                                                                   \
25653deab39SPeter Brune   }                                                                     \
2572a74e5f4SJed Brown   static void CPPJoin2(FetchAndLAND_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
25853deab39SPeter Brune     type *u = (type*)unpacked;                                          \
25953deab39SPeter Brune     type *p = (type*)packed;                                            \
26053deab39SPeter Brune     PetscInt i;                                                         \
26153deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
26253deab39SPeter Brune       PetscInt j = idx[i];                                              \
26353deab39SPeter Brune       type v = u[j];                                                    \
26453deab39SPeter Brune       u[j] = v && p[i];                                                 \
26553deab39SPeter Brune       p[i] = v;                                                         \
26653deab39SPeter Brune     }                                                                   \
26753deab39SPeter Brune   }                                                                     \
2682a74e5f4SJed Brown   static void CPPJoin2(FetchAndLOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
26953deab39SPeter Brune     type *u = (type*)unpacked;                                          \
27053deab39SPeter Brune     type *p = (type*)packed;                                            \
27153deab39SPeter Brune     PetscInt i;                                                         \
27253deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
27353deab39SPeter Brune       PetscInt j = idx[i];                                              \
27453deab39SPeter Brune       type v = u[j];                                                    \
27553deab39SPeter Brune       u[j] = v || p[i];                                                 \
27653deab39SPeter Brune       p[i] = v;                                                         \
27753deab39SPeter Brune     }                                                                   \
27853deab39SPeter Brune   }                                                                     \
2792a74e5f4SJed Brown   static void CPPJoin2(FetchAndLXOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
28053deab39SPeter Brune     type *u = (type*)unpacked;                                          \
28153deab39SPeter Brune     type *p = (type*)packed;                                            \
28253deab39SPeter Brune     PetscInt i;                                                         \
28353deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
28453deab39SPeter Brune       PetscInt j = idx[i];                                              \
28553deab39SPeter Brune       type v = u[j];                                                    \
28653deab39SPeter Brune       u[j] = (!v)!=(!p[i]);                                             \
28753deab39SPeter Brune       p[i] = v;                                                         \
28853deab39SPeter Brune     }                                                                   \
28953deab39SPeter Brune   }                                                                     \
29053deab39SPeter Brune   static void CPPJoin2(PackInit_Logical_,type)(PetscSFBasicPack link) { \
29153deab39SPeter Brune     link->UnpackLAND = CPPJoin2(UnpackLAND_,type);                      \
29253deab39SPeter Brune     link->UnpackLOR  = CPPJoin2(UnpackLOR_,type);                       \
29353deab39SPeter Brune     link->UnpackLXOR = CPPJoin2(UnpackLXOR_,type);                      \
29453deab39SPeter Brune     link->FetchAndLAND = CPPJoin2(FetchAndLAND_,type);                  \
29553deab39SPeter Brune     link->FetchAndLOR  = CPPJoin2(FetchAndLOR_,type);                   \
29653deab39SPeter Brune     link->FetchAndLXOR = CPPJoin2(FetchAndLXOR_,type);                  \
29753deab39SPeter Brune   }
29853deab39SPeter Brune 
29953deab39SPeter Brune 
30053deab39SPeter Brune /* Bitwise Types */
30153deab39SPeter Brune #define DEF_PackBit(type)                                               \
3022a74e5f4SJed Brown   static void CPPJoin2(UnpackBAND_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
30353deab39SPeter Brune     type *u = (type*)unpacked;                                          \
30453deab39SPeter Brune     const type *p = (const type*)packed;                                \
30553deab39SPeter Brune     PetscInt i;                                                         \
30653deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
30753deab39SPeter Brune       type v = u[idx[i]];                                               \
30853deab39SPeter Brune       u[idx[i]] = v & p[i];                                             \
30953deab39SPeter Brune     }                                                                   \
31053deab39SPeter Brune   }                                                                     \
3112a74e5f4SJed Brown   static void CPPJoin2(UnpackBOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
31253deab39SPeter Brune     type *u = (type*)unpacked;                                          \
31353deab39SPeter Brune     const type *p = (const type*)packed;                                \
31453deab39SPeter Brune     PetscInt i;                                                         \
31553deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
31653deab39SPeter Brune       type v = u[idx[i]];                                               \
31753deab39SPeter Brune       u[idx[i]] = v | p[i];                                             \
31853deab39SPeter Brune     }                                                                   \
31953deab39SPeter Brune   }                                                                     \
3202a74e5f4SJed Brown   static void CPPJoin2(UnpackBXOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
32153deab39SPeter Brune     type *u = (type*)unpacked;                                          \
32253deab39SPeter Brune     const type *p = (const type*)packed;                                \
32353deab39SPeter Brune     PetscInt i;                                                         \
32453deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
32553deab39SPeter Brune       type v = u[idx[i]];                                               \
32653deab39SPeter Brune       u[idx[i]] = v^p[i];                                               \
32753deab39SPeter Brune     }                                                                   \
32853deab39SPeter Brune   }                                                                     \
3292a74e5f4SJed Brown   static void CPPJoin2(FetchAndBAND_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
33053deab39SPeter Brune     type *u = (type*)unpacked;                                          \
33153deab39SPeter Brune     type *p = (type*)packed;                                            \
33253deab39SPeter Brune     PetscInt i;                                                         \
33353deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
33453deab39SPeter Brune       PetscInt j = idx[i];                                              \
33553deab39SPeter Brune       type v = u[j];                                                    \
33653deab39SPeter Brune       u[j] = v & p[i];                                                  \
33753deab39SPeter Brune       p[i] = v;                                                         \
33853deab39SPeter Brune     }                                                                   \
33953deab39SPeter Brune   }                                                                     \
3402a74e5f4SJed Brown   static void CPPJoin2(FetchAndBOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
34153deab39SPeter Brune     type *u = (type*)unpacked;                                          \
34253deab39SPeter Brune     type *p = (type*)packed;                                            \
34353deab39SPeter Brune     PetscInt i;                                                         \
34453deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
34553deab39SPeter Brune       PetscInt j = idx[i];                                              \
34653deab39SPeter Brune       type v = u[j];                                                    \
34753deab39SPeter Brune       u[j] = v | p[i];                                                  \
34853deab39SPeter Brune       p[i] = v;                                                         \
34953deab39SPeter Brune     }                                                                   \
35053deab39SPeter Brune   }                                                                     \
3512a74e5f4SJed Brown   static void CPPJoin2(FetchAndBXOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
35253deab39SPeter Brune     type *u = (type*)unpacked;                                          \
35353deab39SPeter Brune     type *p = (type*)packed;                                            \
35453deab39SPeter Brune     PetscInt i;                                                         \
35553deab39SPeter Brune     for (i=0; i<n; i++) {                                               \
35653deab39SPeter Brune       PetscInt j = idx[i];                                              \
35753deab39SPeter Brune       type v = u[j];                                                    \
35853deab39SPeter Brune       u[j] = v^p[i];                                                    \
35953deab39SPeter Brune       p[i] = v;                                                         \
36053deab39SPeter Brune     }                                                                   \
36153deab39SPeter Brune   }                                                                     \
36253deab39SPeter Brune   static void CPPJoin2(PackInit_Bitwise_,type)(PetscSFBasicPack link) { \
36353deab39SPeter Brune     link->UnpackBAND = CPPJoin2(UnpackBAND_,type);                      \
36453deab39SPeter Brune     link->UnpackBOR  = CPPJoin2(UnpackBOR_,type);                       \
36553deab39SPeter Brune     link->UnpackBXOR = CPPJoin2(UnpackBXOR_,type);                      \
36653deab39SPeter Brune     link->FetchAndBAND = CPPJoin2(FetchAndBAND_,type);                  \
36753deab39SPeter Brune     link->FetchAndBOR  = CPPJoin2(FetchAndBOR_,type);                   \
36853deab39SPeter Brune     link->FetchAndBXOR = CPPJoin2(FetchAndBXOR_,type);                  \
36953deab39SPeter Brune   }
37053deab39SPeter Brune 
37195fce210SBarry Smith /* Pair types */
37295fce210SBarry Smith #define CPPJoinloc_exp(base,op,t1,t2) base ## op ## loc_ ## t1 ## _ ## t2
37395fce210SBarry Smith #define CPPJoinloc(base,op,t1,t2) CPPJoinloc_exp(base,op,t1,t2)
37495fce210SBarry Smith #define PairType(type1,type2) CPPJoin3_(_pairtype_,type1,type2)
37595fce210SBarry Smith #define DEF_UnpackXloc(type1,type2,locname,op)                              \
3762a74e5f4SJed Brown   static void CPPJoinloc(Unpack,locname,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
37795fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;        \
37895fce210SBarry Smith     const PairType(type1,type2) *p = (const PairType(type1,type2)*)packed; \
37995fce210SBarry Smith     PetscInt i;                                                         \
38095fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
38195fce210SBarry Smith       PetscInt j = idx[i];                                              \
38295fce210SBarry Smith       if (p[i].a op u[j].a) {                                           \
38395fce210SBarry Smith         u[j].a = p[i].a;                                                \
38495fce210SBarry Smith         u[j].b = p[i].b;                                                \
38595fce210SBarry Smith       } else if (u[j].a == p[i].a) {                                    \
38695fce210SBarry Smith         u[j].b = PetscMin(u[j].b,p[i].b);                               \
38795fce210SBarry Smith       }                                                                 \
38895fce210SBarry Smith     }                                                                   \
38995fce210SBarry Smith   }                                                                     \
3902a74e5f4SJed Brown   static void CPPJoinloc(FetchAnd,locname,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
39195fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;        \
39295fce210SBarry Smith     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;          \
39395fce210SBarry Smith     PetscInt i;                                                         \
39495fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
39595fce210SBarry Smith       PetscInt j = idx[i];                                              \
39695fce210SBarry Smith       PairType(type1,type2) v;                                          \
39795fce210SBarry Smith       v.a = u[j].a;                                                     \
39895fce210SBarry Smith       v.b = u[j].b;                                                     \
39995fce210SBarry Smith       if (p[i].a op u[j].a) {                                           \
40095fce210SBarry Smith         u[j].a = p[i].a;                                                \
40195fce210SBarry Smith         u[j].b = p[i].b;                                                \
40295fce210SBarry Smith       } else if (u[j].a == p[i].a) {                                    \
40395fce210SBarry Smith         u[j].b = PetscMin(u[j].b,p[i].b);                               \
40495fce210SBarry Smith       }                                                                 \
40595fce210SBarry Smith       p[i].a = v.a;                                                     \
40695fce210SBarry Smith       p[i].b = v.b;                                                     \
40795fce210SBarry Smith     }                                                                   \
40895fce210SBarry Smith   }
40995fce210SBarry Smith #define DEF_PackPair(type1,type2)                                       \
41095fce210SBarry Smith   typedef struct {type1 a; type2 b;} PairType(type1,type2);             \
4112a74e5f4SJed Brown   static void CPPJoin3_(Pack_,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,const void *unpacked,void *packed) { \
41295fce210SBarry Smith     const PairType(type1,type2) *u = (const PairType(type1,type2)*)unpacked; \
41395fce210SBarry Smith     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;          \
41495fce210SBarry Smith     PetscInt i;                                                         \
41595fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
41695fce210SBarry Smith       p[i].a = u[idx[i]].a;                                             \
41795fce210SBarry Smith       p[i].b = u[idx[i]].b;                                             \
41895fce210SBarry Smith     }                                                                   \
41995fce210SBarry Smith   }                                                                     \
4202a74e5f4SJed Brown   static void CPPJoin3_(UnpackInsert_,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
42195fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;       \
42295fce210SBarry Smith     const PairType(type1,type2) *p = (const PairType(type1,type2)*)packed; \
42395fce210SBarry Smith     PetscInt i;                                                         \
42495fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
42595fce210SBarry Smith       u[idx[i]].a = p[i].a;                                             \
42695fce210SBarry Smith       u[idx[i]].b = p[i].b;                                             \
42795fce210SBarry Smith     }                                                                   \
42895fce210SBarry Smith   }                                                                     \
4292a74e5f4SJed Brown   static void CPPJoin3_(UnpackAdd_,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \
43095fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;       \
43195fce210SBarry Smith     const PairType(type1,type2) *p = (const PairType(type1,type2)*)packed; \
43295fce210SBarry Smith     PetscInt i;                                                         \
43395fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
43495fce210SBarry Smith       u[idx[i]].a += p[i].a;                                            \
43595fce210SBarry Smith       u[idx[i]].b += p[i].b;                                            \
43695fce210SBarry Smith     }                                                                   \
43795fce210SBarry Smith   }                                                                     \
4382a74e5f4SJed Brown   static void CPPJoin3_(FetchAndInsert_,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
43995fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;        \
44095fce210SBarry Smith     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;          \
44195fce210SBarry Smith     PetscInt i;                                                         \
44295fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
44395fce210SBarry Smith       PetscInt j = idx[i];                                              \
44495fce210SBarry Smith       PairType(type1,type2) v;                                          \
44595fce210SBarry Smith       v.a = u[j].a;                                                     \
44695fce210SBarry Smith       v.b = u[j].b;                                                     \
44795fce210SBarry Smith       u[j].a = p[i].a;                                                  \
44895fce210SBarry Smith       u[j].b = p[i].b;                                                  \
44995fce210SBarry Smith       p[i].a = v.a;                                                     \
45095fce210SBarry Smith       p[i].b = v.b;                                                     \
45195fce210SBarry Smith     }                                                                   \
45295fce210SBarry Smith   }                                                                     \
4532a74e5f4SJed Brown   static void FetchAndAdd_ ## type1 ## _ ## type2(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \
45495fce210SBarry Smith     PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked;       \
45595fce210SBarry Smith     PairType(type1,type2) *p = (PairType(type1,type2)*)packed;         \
45695fce210SBarry Smith     PetscInt i;                                                         \
45795fce210SBarry Smith     for (i=0; i<n; i++) {                                               \
45895fce210SBarry Smith       PetscInt j = idx[i];                                              \
45995fce210SBarry Smith       PairType(type1,type2) v;                                          \
46095fce210SBarry Smith       v.a = u[j].a;                                                     \
46195fce210SBarry Smith       v.b = u[j].b;                                                     \
46295fce210SBarry Smith       u[j].a = v.a + p[i].a;                                            \
46395fce210SBarry Smith       u[j].b = v.b + p[i].b;                                            \
46495fce210SBarry Smith       p[i].a = v.a;                                                     \
46595fce210SBarry Smith       p[i].b = v.b;                                                     \
46695fce210SBarry Smith     }                                                                   \
46795fce210SBarry Smith   }                                                                     \
46895fce210SBarry Smith   DEF_UnpackXloc(type1,type2,Max,>)                                     \
46995fce210SBarry Smith   DEF_UnpackXloc(type1,type2,Min,<)                                     \
47095fce210SBarry Smith   static void CPPJoin3_(PackInit_,type1,type2)(PetscSFBasicPack link) { \
47195fce210SBarry Smith     link->Pack = CPPJoin3_(Pack_,type1,type2);                          \
47295fce210SBarry Smith     link->UnpackInsert = CPPJoin3_(UnpackInsert_,type1,type2);          \
47395fce210SBarry Smith     link->UnpackAdd = CPPJoin3_(UnpackAdd_,type1,type2);                \
47495fce210SBarry Smith     link->UnpackMaxloc = CPPJoin3_(UnpackMaxloc_,type1,type2);          \
47595fce210SBarry Smith     link->UnpackMinloc = CPPJoin3_(UnpackMinloc_,type1,type2);          \
47695fce210SBarry Smith     link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,type1,type2);      \
47795fce210SBarry Smith     link->FetchAndAdd = CPPJoin3_(FetchAndAdd_,type1,type2);            \
47895fce210SBarry Smith     link->FetchAndMaxloc = CPPJoin3_(FetchAndMaxloc_,type1,type2);      \
47995fce210SBarry Smith     link->FetchAndMinloc = CPPJoin3_(FetchAndMinloc_,type1,type2);      \
48095fce210SBarry Smith     link->unitbytes = sizeof(PairType(type1,type2));                    \
48195fce210SBarry Smith   }
48295fce210SBarry Smith 
48395fce210SBarry Smith /* Currently only dumb blocks of data */
48495fce210SBarry Smith #define BlockType(unit,count) CPPJoin3_(_blocktype_,unit,count)
48595fce210SBarry Smith #define DEF_Block(unit,count)                                           \
48695fce210SBarry Smith   typedef struct {unit v[count];} BlockType(unit,count);                \
4872a74e5f4SJed Brown   DEF_PackNoInit(BlockType(unit,count),1)                               \
48895fce210SBarry Smith   static void CPPJoin3_(PackInit_block_,unit,count)(PetscSFBasicPack link) { \
4892a74e5f4SJed Brown     link->Pack = CPPJoin3_(Pack_,BlockType(unit,count),1);               \
4902a74e5f4SJed Brown     link->UnpackInsert = CPPJoin3_(UnpackInsert_,BlockType(unit,count),1); \
4912a74e5f4SJed Brown     link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,BlockType(unit,count),1); \
49295fce210SBarry Smith     link->unitbytes = sizeof(BlockType(unit,count));                    \
49395fce210SBarry Smith   }
49495fce210SBarry Smith 
49595fce210SBarry Smith DEF_PackCmp(int)
49653deab39SPeter Brune DEF_PackBit(int)
49753deab39SPeter Brune DEF_PackLog(int)
49895fce210SBarry Smith DEF_PackCmp(PetscInt)
49953deab39SPeter Brune DEF_PackBit(PetscInt)
50053deab39SPeter Brune DEF_PackLog(PetscInt)
5012a74e5f4SJed Brown DEF_Pack(PetscInt,2)
5022a74e5f4SJed Brown DEF_Pack(PetscInt,3)
5032a74e5f4SJed Brown DEF_Pack(PetscInt,4)
5042a74e5f4SJed Brown DEF_Pack(PetscInt,5)
5052a74e5f4SJed Brown DEF_Pack(PetscInt,7)
50695fce210SBarry Smith DEF_PackCmp(PetscReal)
50753deab39SPeter Brune DEF_PackLog(PetscReal)
5082a74e5f4SJed Brown DEF_Pack(PetscReal,2)
5092a74e5f4SJed Brown DEF_Pack(PetscReal,3)
5102a74e5f4SJed Brown DEF_Pack(PetscReal,4)
5112a74e5f4SJed Brown DEF_Pack(PetscReal,5)
5122a74e5f4SJed Brown DEF_Pack(PetscReal,7)
51395fce210SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
5142a74e5f4SJed Brown DEF_Pack(PetscComplex,1)
5152a74e5f4SJed Brown DEF_Pack(PetscComplex,2)
5162a74e5f4SJed Brown DEF_Pack(PetscComplex,3)
5172a74e5f4SJed Brown DEF_Pack(PetscComplex,4)
5182a74e5f4SJed Brown DEF_Pack(PetscComplex,5)
5192a74e5f4SJed Brown DEF_Pack(PetscComplex,7)
52095fce210SBarry Smith #endif
52195fce210SBarry Smith DEF_PackPair(int,int)
52295fce210SBarry Smith DEF_PackPair(PetscInt,PetscInt)
52395fce210SBarry Smith DEF_Block(int,1)
52495fce210SBarry Smith DEF_Block(int,2)
52595fce210SBarry Smith DEF_Block(int,3)
52695fce210SBarry Smith DEF_Block(int,4)
52795fce210SBarry Smith DEF_Block(int,5)
52895fce210SBarry Smith DEF_Block(int,6)
52995fce210SBarry Smith DEF_Block(int,7)
53095fce210SBarry Smith DEF_Block(int,8)
53195fce210SBarry Smith 
53295fce210SBarry Smith static PetscErrorCode PetscSFSetUp_Basic(PetscSF sf)
53395fce210SBarry Smith {
53495fce210SBarry Smith   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
53595fce210SBarry Smith   PetscErrorCode ierr;
53695fce210SBarry Smith   PetscInt *rlengths,*ilengths,i;
537*c943f53fSJed Brown   PetscMPIInt rank,niranks,*iranks;
53895fce210SBarry Smith   MPI_Comm comm;
539b5a8e515SJed Brown   MPI_Group group;
54095fce210SBarry Smith   MPI_Request *rootreqs,*leafreqs;
54195fce210SBarry Smith 
54295fce210SBarry Smith   PetscFunctionBegin;
543b5a8e515SJed Brown   ierr = MPI_Comm_group(PETSC_COMM_SELF,&group);CHKERRQ(ierr);
544b5a8e515SJed Brown   ierr = PetscSFSetUpRanks(sf,group);CHKERRQ(ierr);
545b5a8e515SJed Brown   ierr = MPI_Group_free(&group);CHKERRQ(ierr);
54695fce210SBarry Smith   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
54795fce210SBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)sf,&bas->tag);CHKERRQ(ierr);
548*c943f53fSJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
54995fce210SBarry Smith   /*
55095fce210SBarry Smith    * Inform roots about how many leaves and from which ranks
55195fce210SBarry Smith    */
552785e854fSJed Brown   ierr = PetscMalloc1(sf->nranks,&rlengths);CHKERRQ(ierr);
55395fce210SBarry Smith   /* Determine number, sending ranks, and length of incoming */
55495fce210SBarry Smith   for (i=0; i<sf->nranks; i++) {
55595fce210SBarry Smith     rlengths[i] = sf->roffset[i+1] - sf->roffset[i]; /* Number of roots referenced by my leaves; for rank sf->ranks[i] */
55695fce210SBarry Smith   }
557*c943f53fSJed Brown   ierr = PetscCommBuildTwoSided(comm,1,MPIU_INT,sf->nranks-sf->ndranks,sf->ranks+sf->ndranks,rlengths+sf->ndranks,&niranks,&iranks,(void**)&ilengths);CHKERRQ(ierr);
558*c943f53fSJed Brown 
559*c943f53fSJed Brown   /* Partition into distinguished and non-distinguished incoming ranks */
560*c943f53fSJed Brown   bas->ndiranks = sf->ndranks;
561*c943f53fSJed Brown   bas->niranks = bas->ndiranks + niranks;
562*c943f53fSJed Brown   ierr = PetscMalloc2(bas->niranks,&bas->iranks,bas->niranks+1,&bas->ioffset);CHKERRQ(ierr);
563*c943f53fSJed Brown   bas->ioffset[0] = 0;
564*c943f53fSJed Brown   for (i=0; i<bas->ndiranks; i++) {
565*c943f53fSJed Brown     bas->iranks[i] = sf->ranks[i];
566*c943f53fSJed Brown     bas->ioffset[i+1] = bas->ioffset[i] + rlengths[i];
567*c943f53fSJed Brown   }
568*c943f53fSJed Brown   if (bas->ndiranks > 1 || (bas->ndiranks == 1 && bas->iranks[0] != rank)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Broken setup for shared ranks");
569*c943f53fSJed Brown   for ( ; i<bas->niranks; i++) {
570*c943f53fSJed Brown     bas->iranks[i] = iranks[i-bas->ndiranks];
571*c943f53fSJed Brown     bas->ioffset[i+1] = bas->ioffset[i] + ilengths[i-bas->ndiranks];
572*c943f53fSJed Brown   }
573*c943f53fSJed Brown   bas->itotal = bas->ioffset[i];
57495fce210SBarry Smith   ierr = PetscFree(rlengths);CHKERRQ(ierr);
575*c943f53fSJed Brown   ierr = PetscFree(iranks);CHKERRQ(ierr);
576*c943f53fSJed Brown   ierr = PetscFree(ilengths);CHKERRQ(ierr);
57795fce210SBarry Smith 
57895fce210SBarry Smith   /* Send leaf identities to roots */
579*c943f53fSJed Brown   ierr = PetscMalloc1(bas->itotal,&bas->irootloc);CHKERRQ(ierr);
580*c943f53fSJed Brown   ierr = PetscMalloc2(bas->niranks-bas->ndiranks,&rootreqs,sf->nranks-sf->ndranks,&leafreqs);CHKERRQ(ierr);
581*c943f53fSJed Brown   for (i=bas->ndiranks; i<bas->niranks; i++) {
582*c943f53fSJed 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);
58395fce210SBarry Smith   }
58495fce210SBarry Smith   for (i=0; i<sf->nranks; i++) {
58595fce210SBarry Smith     PetscMPIInt npoints;
58695fce210SBarry Smith     ierr = PetscMPIIntCast(sf->roffset[i+1] - sf->roffset[i],&npoints);CHKERRQ(ierr);
587*c943f53fSJed Brown     if (i < sf->ndranks) {
588*c943f53fSJed Brown       if (sf->ranks[i] != rank) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot interpret distinguished leaf rank");
589*c943f53fSJed Brown       if (bas->iranks[0] != rank) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot interpret distinguished root rank");
590*c943f53fSJed Brown       if (npoints != bas->ioffset[1]-bas->ioffset[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Distinguished rank exchange has mismatched lengths");
591*c943f53fSJed Brown       ierr = PetscMemcpy(bas->irootloc+bas->ioffset[0],sf->rremote+sf->roffset[i],npoints*sizeof(bas->irootloc[0]));CHKERRQ(ierr);
592*c943f53fSJed Brown       continue;
59395fce210SBarry Smith     }
594*c943f53fSJed Brown     ierr = MPI_Isend(sf->rremote+sf->roffset[i],npoints,MPIU_INT,sf->ranks[i],bas->tag,comm,&leafreqs[i-sf->ndranks]);CHKERRQ(ierr);
595*c943f53fSJed Brown   }
596*c943f53fSJed Brown   ierr = MPI_Waitall(bas->niranks-bas->ndiranks,rootreqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
597*c943f53fSJed Brown   ierr = MPI_Waitall(sf->nranks-sf->ndranks,leafreqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
598eb9baa12SBarry Smith   ierr = PetscFree2(rootreqs,leafreqs);CHKERRQ(ierr);
59995fce210SBarry Smith   PetscFunctionReturn(0);
60095fce210SBarry Smith }
60195fce210SBarry Smith 
60295fce210SBarry Smith static PetscErrorCode PetscSFBasicPackTypeSetup(PetscSFBasicPack link,MPI_Datatype unit)
60395fce210SBarry Smith {
60495fce210SBarry Smith   PetscErrorCode ierr;
60595fce210SBarry Smith   PetscBool      isInt,isPetscInt,isPetscReal,is2Int,is2PetscInt;
606e37f7bbbSSatish Balay   PetscInt       nPetscIntContig,nPetscRealContig;
60795fce210SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
60895fce210SBarry Smith   PetscBool isPetscComplex;
609e37f7bbbSSatish Balay   PetscInt nPetscComplexContig;
61095fce210SBarry Smith #endif
61195fce210SBarry Smith 
61295fce210SBarry Smith   PetscFunctionBegin;
61395fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPI_INT,&isInt);CHKERRQ(ierr);
61495fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPIU_INT,&isPetscInt);CHKERRQ(ierr);
6152a74e5f4SJed Brown   ierr = MPIPetsc_Type_compare_contig(unit,MPIU_INT,&nPetscIntContig);CHKERRQ(ierr);
61695fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPIU_REAL,&isPetscReal);CHKERRQ(ierr);
6172a74e5f4SJed Brown   ierr = MPIPetsc_Type_compare_contig(unit,MPIU_REAL,&nPetscRealContig);CHKERRQ(ierr);
61895fce210SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
61995fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPIU_COMPLEX,&isPetscComplex);CHKERRQ(ierr);
6202a74e5f4SJed Brown   ierr = MPIPetsc_Type_compare_contig(unit,MPIU_COMPLEX,&nPetscComplexContig);CHKERRQ(ierr);
62195fce210SBarry Smith #endif
62295fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPI_2INT,&is2Int);CHKERRQ(ierr);
62395fce210SBarry Smith   ierr = MPIPetsc_Type_compare(unit,MPIU_2INT,&is2PetscInt);CHKERRQ(ierr);
6242a74e5f4SJed Brown   link->bs = 1;
62553deab39SPeter Brune   if (isInt) {PackInit_int(link); PackInit_Logical_int(link); PackInit_Bitwise_int(link);}
62653deab39SPeter Brune   else if (isPetscInt) {PackInit_PetscInt(link); PackInit_Logical_PetscInt(link); PackInit_Bitwise_PetscInt(link);}
62753deab39SPeter Brune   else if (isPetscReal) {PackInit_PetscReal(link); PackInit_Logical_PetscReal(link);}
62895fce210SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
6292a74e5f4SJed Brown   else if (isPetscComplex) PackInit_PetscComplex_1(link);
63095fce210SBarry Smith #endif
63195fce210SBarry Smith   else if (is2Int) PackInit_int_int(link);
63295fce210SBarry Smith   else if (is2PetscInt) PackInit_PetscInt_PetscInt(link);
6332a74e5f4SJed Brown   else if (nPetscIntContig) {
6342a74e5f4SJed Brown     if (nPetscIntContig%7 == 0) PackInit_PetscInt_7(link);
6352a74e5f4SJed Brown     else if (nPetscIntContig%5 == 0) PackInit_PetscInt_5(link);
6362a74e5f4SJed Brown     else if (nPetscIntContig%4 == 0) PackInit_PetscInt_4(link);
6372a74e5f4SJed Brown     else if (nPetscIntContig%3 == 0) PackInit_PetscInt_3(link);
6382a74e5f4SJed Brown     else if (nPetscIntContig%2 == 0) PackInit_PetscInt_2(link);
6392a74e5f4SJed Brown     else PackInit_PetscInt(link);
6402a74e5f4SJed Brown     link->bs = nPetscIntContig;
6412a74e5f4SJed Brown     link->unitbytes *= nPetscIntContig;
6422a74e5f4SJed Brown   } else if (nPetscRealContig) {
6432a74e5f4SJed Brown     if (nPetscRealContig%7 == 0) PackInit_PetscReal_7(link);
6442a74e5f4SJed Brown     else if (nPetscRealContig%5 == 0) PackInit_PetscReal_5(link);
6452a74e5f4SJed Brown     else if (nPetscRealContig%4 == 0) PackInit_PetscReal_4(link);
6462a74e5f4SJed Brown     else if (nPetscRealContig%3 == 0) PackInit_PetscReal_3(link);
6472a74e5f4SJed Brown     else if (nPetscRealContig%2 == 0) PackInit_PetscReal_2(link);
6482a74e5f4SJed Brown     else PackInit_PetscReal(link);
6492a74e5f4SJed Brown     link->bs = nPetscRealContig;
6502a74e5f4SJed Brown     link->unitbytes *= nPetscRealContig;
6512a74e5f4SJed Brown #if defined(PETSC_HAVE_COMPLEX)
6522a74e5f4SJed Brown   } else if (nPetscComplexContig) {
6532a74e5f4SJed Brown     if (nPetscComplexContig%7 == 0) PackInit_PetscComplex_7(link);
6542a74e5f4SJed Brown     else if (nPetscComplexContig%5 == 0) PackInit_PetscComplex_5(link);
6552a74e5f4SJed Brown     else if (nPetscComplexContig%4 == 0) PackInit_PetscComplex_4(link);
6562a74e5f4SJed Brown     else if (nPetscComplexContig%3 == 0) PackInit_PetscComplex_3(link);
6572a74e5f4SJed Brown     else if (nPetscComplexContig%2 == 0) PackInit_PetscComplex_2(link);
6582a74e5f4SJed Brown     else PackInit_PetscComplex_1(link);
6592a74e5f4SJed Brown     link->bs = nPetscComplexContig;
6602a74e5f4SJed Brown     link->unitbytes *= nPetscComplexContig;
6612a74e5f4SJed Brown #endif
6622a74e5f4SJed Brown   } else {
663c3a59b84SJed Brown     MPI_Aint lb,bytes;
664c3a59b84SJed Brown     ierr = MPI_Type_get_extent(unit,&lb,&bytes);CHKERRQ(ierr);
665c3a59b84SJed Brown     if (lb != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Datatype with nonzero lower bound %ld\n",(long)lb);
66695fce210SBarry Smith     if (bytes % sizeof(int)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for type size not divisible by %D",sizeof(int));
66795fce210SBarry Smith     switch (bytes / sizeof(int)) {
66895fce210SBarry Smith     case 1: PackInit_block_int_1(link); break;
66995fce210SBarry Smith     case 2: PackInit_block_int_2(link); break;
67095fce210SBarry Smith     case 3: PackInit_block_int_3(link); break;
67195fce210SBarry Smith     case 4: PackInit_block_int_4(link); break;
67295fce210SBarry Smith     case 5: PackInit_block_int_5(link); break;
67395fce210SBarry Smith     case 6: PackInit_block_int_6(link); break;
67495fce210SBarry Smith     case 7: PackInit_block_int_7(link); break;
67595fce210SBarry Smith     case 8: PackInit_block_int_8(link); break;
67695fce210SBarry Smith     default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for arbitrary block sizes");
67795fce210SBarry Smith     }
67895fce210SBarry Smith   }
67995fce210SBarry Smith   ierr = MPI_Type_dup(unit,&link->unit);CHKERRQ(ierr);
68095fce210SBarry Smith   PetscFunctionReturn(0);
68195fce210SBarry Smith }
68295fce210SBarry Smith 
6832a74e5f4SJed Brown static PetscErrorCode PetscSFBasicPackGetUnpackOp(PetscSF sf,PetscSFBasicPack link,MPI_Op op,void (**UnpackOp)(PetscInt,PetscInt,const PetscInt*,void*,const void*))
68495fce210SBarry Smith {
68595fce210SBarry Smith   PetscFunctionBegin;
68695fce210SBarry Smith   *UnpackOp = NULL;
6878bfbc91cSJed Brown   if (op == MPIU_REPLACE) *UnpackOp = link->UnpackInsert;
68895fce210SBarry Smith   else if (op == MPI_SUM || op == MPIU_SUM) *UnpackOp = link->UnpackAdd;
68953deab39SPeter Brune   else if (op == MPI_PROD) *UnpackOp = link->UnpackMult;
69095fce210SBarry Smith   else if (op == MPI_MAX || op == MPIU_MAX) *UnpackOp = link->UnpackMax;
69195fce210SBarry Smith   else if (op == MPI_MIN || op == MPIU_MIN) *UnpackOp = link->UnpackMin;
69253deab39SPeter Brune   else if (op == MPI_LAND) *UnpackOp = link->UnpackLAND;
69353deab39SPeter Brune   else if (op == MPI_BAND) *UnpackOp = link->UnpackBAND;
69453deab39SPeter Brune   else if (op == MPI_LOR) *UnpackOp = link->UnpackLOR;
69553deab39SPeter Brune   else if (op == MPI_BOR) *UnpackOp = link->UnpackBOR;
69653deab39SPeter Brune   else if (op == MPI_LXOR) *UnpackOp = link->UnpackLXOR;
69753deab39SPeter Brune   else if (op == MPI_BXOR) *UnpackOp = link->UnpackBXOR;
69895fce210SBarry Smith   else if (op == MPI_MAXLOC) *UnpackOp = link->UnpackMaxloc;
69995fce210SBarry Smith   else if (op == MPI_MINLOC) *UnpackOp = link->UnpackMinloc;
7001620da3bSToby Isaac   else *UnpackOp = NULL;
70195fce210SBarry Smith   PetscFunctionReturn(0);
70295fce210SBarry Smith }
7032a74e5f4SJed Brown static PetscErrorCode PetscSFBasicPackGetFetchAndOp(PetscSF sf,PetscSFBasicPack link,MPI_Op op,void (**FetchAndOp)(PetscInt,PetscInt,const PetscInt*,void*,void*))
70495fce210SBarry Smith {
70595fce210SBarry Smith   PetscFunctionBegin;
70695fce210SBarry Smith   *FetchAndOp = NULL;
7078bfbc91cSJed Brown   if (op == MPIU_REPLACE) *FetchAndOp = link->FetchAndInsert;
70895fce210SBarry Smith   else if (op == MPI_SUM || op == MPIU_SUM) *FetchAndOp = link->FetchAndAdd;
70995fce210SBarry Smith   else if (op == MPI_MAX || op == MPIU_MAX) *FetchAndOp = link->FetchAndMax;
71095fce210SBarry Smith   else if (op == MPI_MIN || op == MPIU_MIN) *FetchAndOp = link->FetchAndMin;
71195fce210SBarry Smith   else if (op == MPI_MAXLOC) *FetchAndOp = link->FetchAndMaxloc;
71295fce210SBarry Smith   else if (op == MPI_MINLOC) *FetchAndOp = link->FetchAndMinloc;
71353deab39SPeter Brune   else if (op == MPI_PROD)   *FetchAndOp = link->FetchAndMult;
71453deab39SPeter Brune   else if (op == MPI_LAND)   *FetchAndOp = link->FetchAndLAND;
71553deab39SPeter Brune   else if (op == MPI_BAND)   *FetchAndOp = link->FetchAndBAND;
71653deab39SPeter Brune   else if (op == MPI_LOR)    *FetchAndOp = link->FetchAndLOR;
71753deab39SPeter Brune   else if (op == MPI_BOR)    *FetchAndOp = link->FetchAndBOR;
71853deab39SPeter Brune   else if (op == MPI_LXOR)   *FetchAndOp = link->FetchAndLXOR;
71953deab39SPeter Brune   else if (op == MPI_BXOR)   *FetchAndOp = link->FetchAndBXOR;
72095fce210SBarry Smith   else SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_SUP,"No support for MPI_Op");
72195fce210SBarry Smith   PetscFunctionReturn(0);
72295fce210SBarry Smith }
72395fce210SBarry Smith 
72495fce210SBarry Smith static PetscErrorCode PetscSFBasicPackGetReqs(PetscSF sf,PetscSFBasicPack link,MPI_Request **rootreqs,MPI_Request **leafreqs)
72595fce210SBarry Smith {
72695fce210SBarry Smith   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
72795fce210SBarry Smith 
72895fce210SBarry Smith   PetscFunctionBegin;
72995fce210SBarry Smith   if (rootreqs) *rootreqs = link->requests;
730*c943f53fSJed Brown   if (leafreqs) *leafreqs = link->requests + (bas->niranks - bas->ndiranks);
73195fce210SBarry Smith   PetscFunctionReturn(0);
73295fce210SBarry Smith }
73395fce210SBarry Smith 
73495fce210SBarry Smith static PetscErrorCode PetscSFBasicPackWaitall(PetscSF sf,PetscSFBasicPack link)
73595fce210SBarry Smith {
73695fce210SBarry Smith   PetscSF_Basic  *bas = (PetscSF_Basic*)sf->data;
73795fce210SBarry Smith   PetscErrorCode ierr;
73895fce210SBarry Smith 
73995fce210SBarry Smith   PetscFunctionBegin;
740*c943f53fSJed Brown   ierr = MPI_Waitall(bas->niranks+sf->nranks-(bas->ndiranks+sf->ndranks),link->requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
74195fce210SBarry Smith   PetscFunctionReturn(0);
74295fce210SBarry Smith }
74395fce210SBarry Smith 
744*c943f53fSJed Brown static PetscErrorCode PetscSFBasicGetRootInfo(PetscSF sf,PetscInt *nrootranks,PetscInt *ndrootranks,const PetscMPIInt **rootranks,const PetscInt **rootoffset,const PetscInt **rootloc)
74595fce210SBarry Smith {
74695fce210SBarry Smith   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
74795fce210SBarry Smith 
74895fce210SBarry Smith   PetscFunctionBegin;
74995fce210SBarry Smith   if (nrootranks)  *nrootranks  = bas->niranks;
750*c943f53fSJed Brown   if (ndrootranks) *ndrootranks = bas->ndiranks;
75195fce210SBarry Smith   if (rootranks)   *rootranks   = bas->iranks;
75295fce210SBarry Smith   if (rootoffset)  *rootoffset  = bas->ioffset;
75395fce210SBarry Smith   if (rootloc)     *rootloc     = bas->irootloc;
75495fce210SBarry Smith   PetscFunctionReturn(0);
75595fce210SBarry Smith }
75695fce210SBarry Smith 
757*c943f53fSJed Brown static PetscErrorCode PetscSFBasicGetLeafInfo(PetscSF sf,PetscInt *nleafranks,PetscInt *ndleafranks,const PetscMPIInt **leafranks,const PetscInt **leafoffset,const PetscInt **leafloc)
75895fce210SBarry Smith {
75995fce210SBarry Smith   PetscFunctionBegin;
76095fce210SBarry Smith   if (nleafranks)  *nleafranks  = sf->nranks;
761*c943f53fSJed Brown   if (ndleafranks) *ndleafranks = sf->ndranks;
76295fce210SBarry Smith   if (leafranks)   *leafranks   = sf->ranks;
76395fce210SBarry Smith   if (leafoffset)  *leafoffset  = sf->roffset;
76495fce210SBarry Smith   if (leafloc)     *leafloc     = sf->rmine;
76595fce210SBarry Smith   PetscFunctionReturn(0);
76695fce210SBarry Smith }
76795fce210SBarry Smith 
76895fce210SBarry Smith static PetscErrorCode PetscSFBasicGetPack(PetscSF sf,MPI_Datatype unit,const void *key,PetscSFBasicPack *mylink)
76995fce210SBarry Smith {
77095fce210SBarry Smith   PetscSF_Basic    *bas = (PetscSF_Basic*)sf->data;
77195fce210SBarry Smith   PetscErrorCode   ierr;
77295fce210SBarry Smith   PetscSFBasicPack link,*p;
773*c943f53fSJed Brown   PetscInt         nrootranks,ndrootranks,nleafranks,ndleafranks,i;
77495fce210SBarry Smith   const PetscInt   *rootoffset,*leafoffset;
77595fce210SBarry Smith 
77695fce210SBarry Smith   PetscFunctionBegin;
77795fce210SBarry Smith   /* Look for types in cache */
77895fce210SBarry Smith   for (p=&bas->avail; (link=*p); p=&link->next) {
77995fce210SBarry Smith     PetscBool match;
78095fce210SBarry Smith     ierr = MPIPetsc_Type_compare(unit,link->unit,&match);CHKERRQ(ierr);
78195fce210SBarry Smith     if (match) {
78295fce210SBarry Smith       *p = link->next;          /* Remove from available list */
78395fce210SBarry Smith       goto found;
78495fce210SBarry Smith     }
78595fce210SBarry Smith   }
78695fce210SBarry Smith 
78795fce210SBarry Smith   /* Create new composite types for each send rank */
788*c943f53fSJed Brown   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,&ndrootranks,NULL,&rootoffset,NULL);CHKERRQ(ierr);
789*c943f53fSJed Brown   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,&ndleafranks,NULL,&leafoffset,NULL);CHKERRQ(ierr);
790b00a9115SJed Brown   ierr = PetscNew(&link);CHKERRQ(ierr);
79195fce210SBarry Smith   ierr = PetscSFBasicPackTypeSetup(link,unit);CHKERRQ(ierr);
792bf39f1bfSJed Brown   ierr = PetscMalloc2(nrootranks,&link->root,nleafranks,&link->leaf);CHKERRQ(ierr);
793bf39f1bfSJed Brown   for (i=0; i<nrootranks; i++) {
794bf39f1bfSJed Brown     ierr = PetscMalloc((rootoffset[i+1]-rootoffset[i])*link->unitbytes,&link->root[i]);CHKERRQ(ierr);
795bf39f1bfSJed Brown   }
796bf39f1bfSJed Brown   for (i=0; i<nleafranks; i++) {
797*c943f53fSJed Brown     if (i < ndleafranks) {      /* Leaf buffers for distinguished ranks are pointers directly into root buffers */
798*c943f53fSJed Brown       if (ndrootranks != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot match distinguished ranks");
799*c943f53fSJed Brown       link->leaf[i] = link->root[0];
800*c943f53fSJed Brown       continue;
801*c943f53fSJed Brown     }
802bf39f1bfSJed Brown     ierr = PetscMalloc((leafoffset[i+1]-leafoffset[i])*link->unitbytes,&link->leaf[i]);CHKERRQ(ierr);
803bf39f1bfSJed Brown   }
8040f880758SHong Zhang   ierr = PetscCalloc1(nrootranks+nleafranks,&link->requests);CHKERRQ(ierr);
80595fce210SBarry Smith 
80695fce210SBarry Smith found:
80795fce210SBarry Smith   link->key  = key;
80895fce210SBarry Smith   link->next = bas->inuse;
80995fce210SBarry Smith   bas->inuse = link;
81095fce210SBarry Smith 
81195fce210SBarry Smith   *mylink = link;
81295fce210SBarry Smith   PetscFunctionReturn(0);
81395fce210SBarry Smith }
81495fce210SBarry Smith 
81595fce210SBarry Smith static PetscErrorCode PetscSFBasicGetPackInUse(PetscSF sf,MPI_Datatype unit,const void *key,PetscCopyMode cmode,PetscSFBasicPack *mylink)
81695fce210SBarry Smith {
81795fce210SBarry Smith   PetscSF_Basic    *bas = (PetscSF_Basic*)sf->data;
81895fce210SBarry Smith   PetscErrorCode   ierr;
81995fce210SBarry Smith   PetscSFBasicPack link,*p;
82095fce210SBarry Smith 
82195fce210SBarry Smith   PetscFunctionBegin;
82295fce210SBarry Smith   /* Look for types in cache */
82395fce210SBarry Smith   for (p=&bas->inuse; (link=*p); p=&link->next) {
82495fce210SBarry Smith     PetscBool match;
82595fce210SBarry Smith     ierr = MPIPetsc_Type_compare(unit,link->unit,&match);CHKERRQ(ierr);
82645038af3SLawrence Mitchell     if (match && (key == link->key)) {
82795fce210SBarry Smith       switch (cmode) {
82895fce210SBarry Smith       case PETSC_OWN_POINTER: *p = link->next; break; /* Remove from inuse list */
82995fce210SBarry Smith       case PETSC_USE_POINTER: break;
83095fce210SBarry Smith       default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"invalid cmode");
83195fce210SBarry Smith       }
83295fce210SBarry Smith       *mylink = link;
83395fce210SBarry Smith       PetscFunctionReturn(0);
83495fce210SBarry Smith     }
83595fce210SBarry Smith   }
83695fce210SBarry Smith   SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Could not find pack");
83795fce210SBarry Smith   PetscFunctionReturn(0);
83895fce210SBarry Smith }
83995fce210SBarry Smith 
84095fce210SBarry Smith static PetscErrorCode PetscSFBasicReclaimPack(PetscSF sf,PetscSFBasicPack *link)
84195fce210SBarry Smith {
84295fce210SBarry Smith   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;
84395fce210SBarry Smith 
84495fce210SBarry Smith   PetscFunctionBegin;
84595fce210SBarry Smith   (*link)->key  = NULL;
84695fce210SBarry Smith   (*link)->next = bas->avail;
84795fce210SBarry Smith   bas->avail    = *link;
84895fce210SBarry Smith   *link         = NULL;
84995fce210SBarry Smith   PetscFunctionReturn(0);
85095fce210SBarry Smith }
85195fce210SBarry Smith 
8524416b707SBarry Smith static PetscErrorCode PetscSFSetFromOptions_Basic(PetscOptionItems *PetscOptionsObject,PetscSF sf)
85395fce210SBarry Smith {
85495fce210SBarry Smith   PetscErrorCode ierr;
85595fce210SBarry Smith 
85695fce210SBarry Smith   PetscFunctionBegin;
857e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"PetscSF Basic options");CHKERRQ(ierr);
85895fce210SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
85995fce210SBarry Smith   PetscFunctionReturn(0);
86095fce210SBarry Smith }
86195fce210SBarry Smith 
86295fce210SBarry Smith static PetscErrorCode PetscSFReset_Basic(PetscSF sf)
86395fce210SBarry Smith {
86495fce210SBarry Smith   PetscSF_Basic    *bas = (PetscSF_Basic*)sf->data;
86595fce210SBarry Smith   PetscErrorCode   ierr;
86695fce210SBarry Smith   PetscSFBasicPack link,next;
86795fce210SBarry Smith 
86895fce210SBarry Smith   PetscFunctionBegin;
869*c943f53fSJed Brown   ierr = PetscFree2(bas->iranks,bas->ioffset);CHKERRQ(ierr);
870*c943f53fSJed Brown   ierr = PetscFree(bas->irootloc);CHKERRQ(ierr);
87195fce210SBarry Smith   if (bas->inuse) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Outstanding operation has not been completed");
87295fce210SBarry Smith   for (link=bas->avail; link; link=next) {
873bf39f1bfSJed Brown     PetscInt i;
87495fce210SBarry Smith     next = link->next;
875bb86044dSJed Brown #if defined(PETSC_HAVE_MPI_TYPE_DUP)
87695fce210SBarry Smith     ierr = MPI_Type_free(&link->unit);CHKERRQ(ierr);
877bb86044dSJed Brown #endif
878bf39f1bfSJed Brown     for (i=0; i<bas->niranks; i++) {ierr = PetscFree(link->root[i]);CHKERRQ(ierr);}
879*c943f53fSJed Brown     for (i=sf->ndranks; i<sf->nranks; i++) {ierr = PetscFree(link->leaf[i]);CHKERRQ(ierr);} /* Free only non-distinguished leaf buffers */
88095fce210SBarry Smith     ierr = PetscFree2(link->root,link->leaf);CHKERRQ(ierr);
88195fce210SBarry Smith     ierr = PetscFree(link->requests);CHKERRQ(ierr);
88295fce210SBarry Smith     ierr = PetscFree(link);CHKERRQ(ierr);
88395fce210SBarry Smith   }
88495fce210SBarry Smith   bas->avail = NULL;
88595fce210SBarry Smith   PetscFunctionReturn(0);
88695fce210SBarry Smith }
88795fce210SBarry Smith 
88895fce210SBarry Smith static PetscErrorCode PetscSFDestroy_Basic(PetscSF sf)
88995fce210SBarry Smith {
89095fce210SBarry Smith   PetscErrorCode ierr;
89195fce210SBarry Smith 
89295fce210SBarry Smith   PetscFunctionBegin;
89395fce210SBarry Smith   ierr = PetscFree(sf->data);CHKERRQ(ierr);
89495fce210SBarry Smith   PetscFunctionReturn(0);
89595fce210SBarry Smith }
89695fce210SBarry Smith 
89795fce210SBarry Smith static PetscErrorCode PetscSFView_Basic(PetscSF sf,PetscViewer viewer)
89895fce210SBarry Smith {
89995fce210SBarry Smith   /* PetscSF_Basic *bas = (PetscSF_Basic*)sf->data; */
90095fce210SBarry Smith   PetscErrorCode ierr;
90195fce210SBarry Smith   PetscBool      iascii;
90295fce210SBarry Smith 
90395fce210SBarry Smith   PetscFunctionBegin;
90495fce210SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
90595fce210SBarry Smith   if (iascii) {
90695fce210SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  sort=%s\n",sf->rankorder ? "rank-order" : "unordered");CHKERRQ(ierr);
90795fce210SBarry Smith   }
90895fce210SBarry Smith   PetscFunctionReturn(0);
90995fce210SBarry Smith }
91095fce210SBarry Smith 
91195fce210SBarry Smith /* Send from roots to leaves */
91295fce210SBarry Smith static PetscErrorCode PetscSFBcastBegin_Basic(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
91395fce210SBarry Smith {
91495fce210SBarry Smith   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;
91595fce210SBarry Smith   PetscErrorCode    ierr;
91695fce210SBarry Smith   PetscSFBasicPack  link;
917*c943f53fSJed Brown   PetscInt          i,nrootranks,ndrootranks,nleafranks,ndleafranks;
91895fce210SBarry Smith   const PetscInt    *rootoffset,*leafoffset,*rootloc,*leafloc;
91995fce210SBarry Smith   const PetscMPIInt *rootranks,*leafranks;
92095fce210SBarry Smith   MPI_Request       *rootreqs,*leafreqs;
92195fce210SBarry Smith 
92295fce210SBarry Smith   PetscFunctionBegin;
923*c943f53fSJed Brown   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,&ndrootranks,&rootranks,&rootoffset,&rootloc);CHKERRQ(ierr);
924*c943f53fSJed Brown   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,&ndleafranks,&leafranks,&leafoffset,&leafloc);CHKERRQ(ierr);
92595fce210SBarry Smith   ierr = PetscSFBasicGetPack(sf,unit,rootdata,&link);CHKERRQ(ierr);
92695fce210SBarry Smith 
92795fce210SBarry Smith   ierr = PetscSFBasicPackGetReqs(sf,link,&rootreqs,&leafreqs);CHKERRQ(ierr);
928*c943f53fSJed Brown   /* Eagerly post leaf receives, but only from non-distinguished ranks -- distinguished ranks will receive via shared memory */
929*c943f53fSJed Brown   for (i=ndleafranks; i<nleafranks; i++) {
93095fce210SBarry Smith     PetscMPIInt n = leafoffset[i+1] - leafoffset[i];
931*c943f53fSJed Brown     ierr = MPI_Irecv(link->leaf[i],n,unit,leafranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&leafreqs[i-ndleafranks]);CHKERRQ(ierr);
93295fce210SBarry Smith   }
93395fce210SBarry Smith   /* Pack and send root data */
93495fce210SBarry Smith   for (i=0; i<nrootranks; i++) {
93595fce210SBarry Smith     PetscMPIInt n          = rootoffset[i+1] - rootoffset[i];
936bf39f1bfSJed Brown     void        *packstart = link->root[i];
9372a74e5f4SJed Brown     (*link->Pack)(n,link->bs,rootloc+rootoffset[i],rootdata,packstart);
938*c943f53fSJed Brown     if (i < ndrootranks) continue; /* shared memory */
939*c943f53fSJed Brown     ierr = MPI_Isend(packstart,n,unit,rootranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&rootreqs[i-ndrootranks]);CHKERRQ(ierr);
94095fce210SBarry Smith   }
94195fce210SBarry Smith   PetscFunctionReturn(0);
94295fce210SBarry Smith }
94395fce210SBarry Smith 
94495fce210SBarry Smith PetscErrorCode PetscSFBcastEnd_Basic(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
94595fce210SBarry Smith {
94695fce210SBarry Smith   PetscErrorCode   ierr;
94795fce210SBarry Smith   PetscSFBasicPack link;
948*c943f53fSJed Brown   PetscInt         i,nleafranks,ndleafranks;
94995fce210SBarry Smith   const PetscInt   *leafoffset,*leafloc;
95095fce210SBarry Smith 
95195fce210SBarry Smith   PetscFunctionBegin;
95295fce210SBarry Smith   ierr = PetscSFBasicGetPackInUse(sf,unit,rootdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr);
95395fce210SBarry Smith   ierr = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
954*c943f53fSJed Brown   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,&ndleafranks,NULL,&leafoffset,&leafloc);CHKERRQ(ierr);
95595fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
95695fce210SBarry Smith     PetscMPIInt n          = leafoffset[i+1] - leafoffset[i];
957bf39f1bfSJed Brown     const void  *packstart = link->leaf[i];
9582a74e5f4SJed Brown     (*link->UnpackInsert)(n,link->bs,leafloc+leafoffset[i],leafdata,packstart);
95995fce210SBarry Smith   }
96095fce210SBarry Smith   ierr = PetscSFBasicReclaimPack(sf,&link);CHKERRQ(ierr);
96195fce210SBarry Smith   PetscFunctionReturn(0);
96295fce210SBarry Smith }
96395fce210SBarry Smith 
96495fce210SBarry Smith /* leaf -> root with reduction */
96595fce210SBarry Smith PetscErrorCode PetscSFReduceBegin_Basic(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
96695fce210SBarry Smith {
96795fce210SBarry Smith   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;
96895fce210SBarry Smith   PetscSFBasicPack  link;
96995fce210SBarry Smith   PetscErrorCode    ierr;
970*c943f53fSJed Brown   PetscInt          i,nrootranks,ndrootranks,nleafranks,ndleafranks;
97195fce210SBarry Smith   const PetscInt    *rootoffset,*leafoffset,*rootloc,*leafloc;
97295fce210SBarry Smith   const PetscMPIInt *rootranks,*leafranks;
97395fce210SBarry Smith   MPI_Request       *rootreqs,*leafreqs;
97495fce210SBarry Smith 
97595fce210SBarry Smith   PetscFunctionBegin;
976*c943f53fSJed Brown   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,&ndrootranks,&rootranks,&rootoffset,&rootloc);CHKERRQ(ierr);
977*c943f53fSJed Brown   ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,&ndleafranks,&leafranks,&leafoffset,&leafloc);CHKERRQ(ierr);
97895fce210SBarry Smith   ierr = PetscSFBasicGetPack(sf,unit,rootdata,&link);CHKERRQ(ierr);
97995fce210SBarry Smith 
98095fce210SBarry Smith   ierr = PetscSFBasicPackGetReqs(sf,link,&rootreqs,&leafreqs);CHKERRQ(ierr);
981*c943f53fSJed Brown   /* Eagerly post root receives for non-distinguished ranks */
982*c943f53fSJed Brown   for (i=ndrootranks; i<nrootranks; i++) {
98395fce210SBarry Smith     PetscMPIInt n = rootoffset[i+1] - rootoffset[i];
984*c943f53fSJed Brown     ierr = MPI_Irecv(link->root[i],n,unit,rootranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&rootreqs[i-ndrootranks]);CHKERRQ(ierr);
98595fce210SBarry Smith   }
98695fce210SBarry Smith   /* Pack and send leaf data */
98795fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
98895fce210SBarry Smith     PetscMPIInt n          = leafoffset[i+1] - leafoffset[i];
989bf39f1bfSJed Brown     void        *packstart = link->leaf[i];
9902a74e5f4SJed Brown     (*link->Pack)(n,link->bs,leafloc+leafoffset[i],leafdata,packstart);
991*c943f53fSJed Brown     if (i < ndleafranks) continue; /* shared memory */
992*c943f53fSJed Brown     ierr = MPI_Isend(packstart,n,unit,leafranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&leafreqs[i-ndleafranks]);CHKERRQ(ierr);
99395fce210SBarry Smith   }
99495fce210SBarry Smith   PetscFunctionReturn(0);
99595fce210SBarry Smith }
99695fce210SBarry Smith 
99795fce210SBarry Smith static PetscErrorCode PetscSFReduceEnd_Basic(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
99895fce210SBarry Smith {
9992a74e5f4SJed Brown   void             (*UnpackOp)(PetscInt,PetscInt,const PetscInt*,void*,const void*);
100095fce210SBarry Smith   PetscErrorCode   ierr;
100195fce210SBarry Smith   PetscSFBasicPack link;
100295fce210SBarry Smith   PetscInt         i,nrootranks;
10031620da3bSToby Isaac   PetscMPIInt      typesize = -1;
100495fce210SBarry Smith   const PetscInt   *rootoffset,*rootloc;
100595fce210SBarry Smith 
100695fce210SBarry Smith   PetscFunctionBegin;
100795fce210SBarry Smith   ierr = PetscSFBasicGetPackInUse(sf,unit,rootdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr);
100895fce210SBarry Smith   /* This implementation could be changed to unpack as receives arrive, at the cost of non-determinism */
100995fce210SBarry Smith   ierr = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
1010*c943f53fSJed Brown   ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,NULL,NULL,&rootoffset,&rootloc);CHKERRQ(ierr);
101195fce210SBarry Smith   ierr = PetscSFBasicPackGetUnpackOp(sf,link,op,&UnpackOp);CHKERRQ(ierr);
10121620da3bSToby Isaac   if (UnpackOp) {
10131620da3bSToby Isaac     typesize = link->unitbytes;
10141620da3bSToby Isaac   }
10151620da3bSToby Isaac   else {
10161620da3bSToby Isaac     ierr = MPI_Type_size(unit,&typesize);CHKERRQ(ierr);
10171620da3bSToby Isaac   }
101895fce210SBarry Smith   for (i=0; i<nrootranks; i++) {
101995fce210SBarry Smith     PetscMPIInt n   = rootoffset[i+1] - rootoffset[i];
1020bf39f1bfSJed Brown     char *packstart = (char *) link->root[i];
102195fce210SBarry Smith 
10221620da3bSToby Isaac     if (UnpackOp) {
102388ac5ce8SToby Isaac       (*UnpackOp)(n,link->bs,rootloc+rootoffset[i],rootdata,(const void *)packstart);
102495fce210SBarry Smith     }
1025a9ec7e9eSToby Isaac #if PETSC_HAVE_MPI_REDUCE_LOCAL
102678ec4774SToby Isaac     else if (n) { /* the op should be defined to operate on the whole datatype, so we ignore link->bs */
102778ec4774SToby Isaac       PetscInt j;
10281620da3bSToby Isaac 
10291620da3bSToby Isaac       for (j = 0; j < n; j++) {
103078ec4774SToby Isaac         ierr = MPI_Reduce_local(packstart+j*typesize,((char *) rootdata)+(rootloc[rootoffset[i]+j])*typesize,1,unit,op);CHKERRQ(ierr);
10311620da3bSToby Isaac       }
10321620da3bSToby Isaac     }
1033a9ec7e9eSToby Isaac #else
1034a9ec7e9eSToby Isaac     else {
1035a9ec7e9eSToby Isaac       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No unpacking reduction operation for this MPI_Op");
1036a9ec7e9eSToby Isaac     }
1037a9ec7e9eSToby Isaac #endif
10381620da3bSToby Isaac   }
103995fce210SBarry Smith   ierr = PetscSFBasicReclaimPack(sf,&link);CHKERRQ(ierr);
104095fce210SBarry Smith   PetscFunctionReturn(0);
104195fce210SBarry Smith }
104295fce210SBarry Smith 
104395fce210SBarry Smith static PetscErrorCode PetscSFFetchAndOpBegin_Basic(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
104495fce210SBarry Smith {
104595fce210SBarry Smith   PetscErrorCode ierr;
104695fce210SBarry Smith 
104795fce210SBarry Smith   PetscFunctionBegin;
104895fce210SBarry Smith   ierr = PetscSFReduceBegin_Basic(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
104995fce210SBarry Smith   PetscFunctionReturn(0);
105095fce210SBarry Smith }
105195fce210SBarry Smith 
105295fce210SBarry Smith static PetscErrorCode PetscSFFetchAndOpEnd_Basic(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
105395fce210SBarry Smith {
105495fce210SBarry Smith   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;
10552a74e5f4SJed Brown   void              (*FetchAndOp)(PetscInt,PetscInt,const PetscInt*,void*,void*);
105695fce210SBarry Smith   PetscErrorCode    ierr;
105795fce210SBarry Smith   PetscSFBasicPack  link;
1058*c943f53fSJed Brown   PetscInt          i,nrootranks,ndrootranks,nleafranks,ndleafranks;
105995fce210SBarry Smith   const PetscInt    *rootoffset,*leafoffset,*rootloc,*leafloc;
106095fce210SBarry Smith   const PetscMPIInt *rootranks,*leafranks;
106195fce210SBarry Smith   MPI_Request       *rootreqs,*leafreqs;
106295fce210SBarry Smith 
106395fce210SBarry Smith   PetscFunctionBegin;
106495fce210SBarry Smith   ierr = PetscSFBasicGetPackInUse(sf,unit,rootdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr);
106595fce210SBarry Smith   /* This implementation could be changed to unpack as receives arrive, at the cost of non-determinism */
106695fce210SBarry Smith   ierr      = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
1067*c943f53fSJed Brown   ierr      = PetscSFBasicGetRootInfo(sf,&nrootranks,&ndrootranks,&rootranks,&rootoffset,&rootloc);CHKERRQ(ierr);
1068*c943f53fSJed Brown   ierr      = PetscSFBasicGetLeafInfo(sf,&nleafranks,&ndleafranks,&leafranks,&leafoffset,&leafloc);CHKERRQ(ierr);
106995fce210SBarry Smith   ierr      = PetscSFBasicPackGetReqs(sf,link,&rootreqs,&leafreqs);CHKERRQ(ierr);
107095fce210SBarry Smith   /* Post leaf receives */
1071*c943f53fSJed Brown   for (i=ndleafranks; i<nleafranks; i++) {
107295fce210SBarry Smith     PetscMPIInt n = leafoffset[i+1] - leafoffset[i];
1073*c943f53fSJed Brown     ierr = MPI_Irecv(link->leaf[i],n,unit,leafranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&leafreqs[i-ndleafranks]);CHKERRQ(ierr);
107495fce210SBarry Smith   }
107595fce210SBarry Smith   /* Process local fetch-and-op, post root sends */
107695fce210SBarry Smith   ierr = PetscSFBasicPackGetFetchAndOp(sf,link,op,&FetchAndOp);CHKERRQ(ierr);
107795fce210SBarry Smith   for (i=0; i<nrootranks; i++) {
107895fce210SBarry Smith     PetscMPIInt n          = rootoffset[i+1] - rootoffset[i];
1079bf39f1bfSJed Brown     void        *packstart = link->root[i];
108095fce210SBarry Smith 
10812a74e5f4SJed Brown     (*FetchAndOp)(n,link->bs,rootloc+rootoffset[i],rootdata,packstart);
1082*c943f53fSJed Brown     if (i < ndrootranks) continue; /* shared memory */
1083*c943f53fSJed Brown     ierr = MPI_Isend(packstart,n,unit,rootranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&rootreqs[i-ndrootranks]);CHKERRQ(ierr);
108495fce210SBarry Smith   }
108595fce210SBarry Smith   ierr = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr);
108695fce210SBarry Smith   for (i=0; i<nleafranks; i++) {
108795fce210SBarry Smith     PetscMPIInt n          = leafoffset[i+1] - leafoffset[i];
1088bf39f1bfSJed Brown     const void  *packstart = link->leaf[i];
10892a74e5f4SJed Brown     (*link->UnpackInsert)(n,link->bs,leafloc+leafoffset[i],leafupdate,packstart);
109095fce210SBarry Smith   }
109195fce210SBarry Smith   ierr = PetscSFBasicReclaimPack(sf,&link);CHKERRQ(ierr);
109295fce210SBarry Smith   PetscFunctionReturn(0);
109395fce210SBarry Smith }
109495fce210SBarry Smith 
10958cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscSFCreate_Basic(PetscSF sf)
109695fce210SBarry Smith {
109795fce210SBarry Smith   PetscSF_Basic  *bas = (PetscSF_Basic*)sf->data;
109895fce210SBarry Smith   PetscErrorCode ierr;
109995fce210SBarry Smith 
110095fce210SBarry Smith   PetscFunctionBegin;
110195fce210SBarry Smith   sf->ops->SetUp           = PetscSFSetUp_Basic;
110295fce210SBarry Smith   sf->ops->SetFromOptions  = PetscSFSetFromOptions_Basic;
110395fce210SBarry Smith   sf->ops->Reset           = PetscSFReset_Basic;
110495fce210SBarry Smith   sf->ops->Destroy         = PetscSFDestroy_Basic;
110595fce210SBarry Smith   sf->ops->View            = PetscSFView_Basic;
110695fce210SBarry Smith   sf->ops->BcastBegin      = PetscSFBcastBegin_Basic;
110795fce210SBarry Smith   sf->ops->BcastEnd        = PetscSFBcastEnd_Basic;
110895fce210SBarry Smith   sf->ops->ReduceBegin     = PetscSFReduceBegin_Basic;
110995fce210SBarry Smith   sf->ops->ReduceEnd       = PetscSFReduceEnd_Basic;
111095fce210SBarry Smith   sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Basic;
111195fce210SBarry Smith   sf->ops->FetchAndOpEnd   = PetscSFFetchAndOpEnd_Basic;
111295fce210SBarry Smith 
1113b00a9115SJed Brown   ierr     = PetscNewLog(sf,&bas);CHKERRQ(ierr);
111495fce210SBarry Smith   sf->data = (void*)bas;
111595fce210SBarry Smith   PetscFunctionReturn(0);
111695fce210SBarry Smith }
1117