1 2 #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/ 3 4 typedef struct _n_PetscSFBasicPack *PetscSFBasicPack; 5 struct _n_PetscSFBasicPack { 6 void (*Pack)(PetscInt,PetscInt,const PetscInt*,const void*,void*); 7 void (*UnpackInsert)(PetscInt,PetscInt,const PetscInt*,void*,const void*); 8 void (*UnpackAdd)(PetscInt,PetscInt,const PetscInt*,void*,const void*); 9 void (*UnpackMin)(PetscInt,PetscInt,const PetscInt*,void*,const void*); 10 void (*UnpackMax)(PetscInt,PetscInt,const PetscInt*,void*,const void*); 11 void (*UnpackMinloc)(PetscInt,PetscInt,const PetscInt*,void*,const void*); 12 void (*UnpackMaxloc)(PetscInt,PetscInt,const PetscInt*,void*,const void*); 13 void (*UnpackMult)(PetscInt,PetscInt,const PetscInt*,void*,const void *); 14 void (*UnpackLAND)(PetscInt,PetscInt,const PetscInt*,void*,const void *); 15 void (*UnpackBAND)(PetscInt,PetscInt,const PetscInt*,void*,const void *); 16 void (*UnpackLOR)(PetscInt,PetscInt,const PetscInt*,void*,const void *); 17 void (*UnpackBOR)(PetscInt,PetscInt,const PetscInt*,void*,const void *); 18 void (*UnpackLXOR)(PetscInt,PetscInt,const PetscInt*,void*,const void *); 19 void (*UnpackBXOR)(PetscInt,PetscInt,const PetscInt*,void*,const void *); 20 void (*FetchAndInsert)(PetscInt,PetscInt,const PetscInt*,void*,void*); 21 void (*FetchAndAdd)(PetscInt,PetscInt,const PetscInt*,void*,void*); 22 void (*FetchAndMin)(PetscInt,PetscInt,const PetscInt*,void*,void*); 23 void (*FetchAndMax)(PetscInt,PetscInt,const PetscInt*,void*,void*); 24 void (*FetchAndMinloc)(PetscInt,PetscInt,const PetscInt*,void*,void*); 25 void (*FetchAndMaxloc)(PetscInt,PetscInt,const PetscInt*,void*,void*); 26 void (*FetchAndMult)(PetscInt,PetscInt,const PetscInt*,void*,void*); 27 void (*FetchAndLAND)(PetscInt,PetscInt,const PetscInt*,void*,void*); 28 void (*FetchAndBAND)(PetscInt,PetscInt,const PetscInt*,void*,void*); 29 void (*FetchAndLOR)(PetscInt,PetscInt,const PetscInt*,void*,void*); 30 void (*FetchAndBOR)(PetscInt,PetscInt,const PetscInt*,void*,void*); 31 void (*FetchAndLXOR)(PetscInt,PetscInt,const PetscInt*,void*,void*); 32 void (*FetchAndBXOR)(PetscInt,PetscInt,const PetscInt*,void*,void*); 33 34 MPI_Datatype unit; 35 PetscBool isbuiltin; /* Is unit an MPI builtin datatype? */ 36 size_t unitbytes; /* Number of bytes in a unit */ 37 PetscInt bs; /* Number of basic units in a unit */ 38 const void *key; /* Array used as key for operation */ 39 char **root; /* Packed root data, indexed by leaf rank */ 40 char **leaf; /* Packed leaf data, indexed by root rank */ 41 MPI_Request *requests; /* Array of root requests followed by leaf requests */ 42 PetscSFBasicPack next; 43 }; 44 45 typedef struct { 46 PetscMPIInt tag; 47 PetscMPIInt niranks; /* Number of incoming ranks (ranks accessing my roots) */ 48 PetscMPIInt ndiranks; /* Number of incoming ranks (ranks accessing my roots) in distinguished set */ 49 PetscMPIInt *iranks; /* Array of ranks that reference my roots */ 50 PetscInt itotal; /* Total number of graph edges referencing my roots */ 51 PetscInt *ioffset; /* Array of length niranks+1 holding offset in irootloc[] for each rank */ 52 PetscInt *irootloc; /* Incoming roots referenced by ranks starting at ioffset[rank] */ 53 PetscSFBasicPack avail; /* One or more entries per MPI Datatype, lazily constructed */ 54 PetscSFBasicPack inuse; /* Buffers being used for transactions that have not yet completed */ 55 } PetscSF_Basic; 56 57 #if !defined(PETSC_HAVE_MPI_TYPE_DUP) 58 PETSC_STATIC_INLINE int MPI_Type_dup(MPI_Datatype datatype,MPI_Datatype *newtype) 59 { 60 int ierr; 61 ierr = MPI_Type_contiguous(1,datatype,newtype); if (ierr) return ierr; 62 ierr = MPI_Type_commit(newtype); if (ierr) return ierr; 63 return MPI_SUCCESS; 64 } 65 #endif 66 67 /* 68 * MPI_Reduce_local is not really useful because it can't handle sparse data and it vectorizes "in the wrong direction", 69 * therefore we pack data types manually. This section defines packing routines for the standard data types. 70 */ 71 72 #define CPPJoin2_exp(a,b) a ## b 73 #define CPPJoin2(a,b) CPPJoin2_exp(a,b) 74 #define CPPJoin3_exp_(a,b,c) a ## b ## _ ## c 75 #define CPPJoin3_(a,b,c) CPPJoin3_exp_(a,b,c) 76 77 /* Basic types without addition */ 78 #define DEF_PackNoInit(type,BS) \ 79 static void CPPJoin3_(Pack_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,const void *unpacked,void *packed) { \ 80 const type *u = (const type*)unpacked; \ 81 type *p = (type*)packed; \ 82 PetscInt i,j,k; \ 83 for (i=0; i<n; i++) \ 84 for (j=0; j<bs; j+=BS) \ 85 for (k=j; k<j+BS; k++) \ 86 p[i*bs+k] = u[idx[i]*bs+k]; \ 87 } \ 88 static void CPPJoin3_(UnpackInsert_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \ 89 type *u = (type*)unpacked; \ 90 const type *p = (const type*)packed; \ 91 PetscInt i,j,k; \ 92 for (i=0; i<n; i++) \ 93 for (j=0; j<bs; j+=BS) \ 94 for (k=j; k<j+BS; k++) \ 95 u[idx[i]*bs+k] = p[i*bs+k]; \ 96 } \ 97 static void CPPJoin3_(FetchAndInsert_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \ 98 type *u = (type*)unpacked; \ 99 type *p = (type*)packed; \ 100 PetscInt i,j,k; \ 101 for (i=0; i<n; i++) { \ 102 PetscInt ii = idx[i]; \ 103 for (j=0; j<bs; j+=BS) \ 104 for (k=j; k<j+BS; k++) { \ 105 type t = u[ii*bs+k]; \ 106 u[ii*bs+k] = p[i*bs+k]; \ 107 p[i*bs+k] = t; \ 108 } \ 109 } \ 110 } 111 112 /* Basic types defining addition */ 113 #define DEF_PackAddNoInit(type,BS) \ 114 DEF_PackNoInit(type,BS) \ 115 static void CPPJoin3_(UnpackAdd_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \ 116 type *u = (type*)unpacked; \ 117 const type *p = (const type*)packed; \ 118 PetscInt i,j,k; \ 119 for (i=0; i<n; i++) \ 120 for (j=0; j<bs; j+=BS) \ 121 for (k=j; k<j+BS; k++) \ 122 u[idx[i]*bs+k] += p[i*bs+k]; \ 123 } \ 124 static void CPPJoin3_(FetchAndAdd_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \ 125 type *u = (type*)unpacked; \ 126 type *p = (type*)packed; \ 127 PetscInt i,j,k; \ 128 for (i=0; i<n; i++) { \ 129 PetscInt ii = idx[i]; \ 130 for (j=0; j<bs; j+=BS) \ 131 for (k=j; k<j+BS; k++) { \ 132 type t = u[ii*bs+k]; \ 133 u[ii*bs+k] = t + p[i*bs+k]; \ 134 p[i*bs+k] = t; \ 135 } \ 136 } \ 137 } \ 138 static void CPPJoin3_(UnpackMult_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \ 139 type *u = (type*)unpacked; \ 140 const type *p = (const type*)packed; \ 141 PetscInt i,j,k; \ 142 for (i=0; i<n; i++) \ 143 for (j=0; j<bs; j+=BS) \ 144 for (k=j; k<j+BS; k++) \ 145 u[idx[i]*bs+k] *= p[i*bs+k]; \ 146 } \ 147 static void CPPJoin3_(FetchAndMult_,type,BS)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \ 148 type *u = (type*)unpacked; \ 149 type *p = (type*)packed; \ 150 PetscInt i,j,k; \ 151 for (i=0; i<n; i++) { \ 152 PetscInt ii = idx[i]; \ 153 for (j=0; j<bs; j+=BS) \ 154 for (k=j; k<j+BS; k++) { \ 155 type t = u[ii*bs+k]; \ 156 u[ii*bs+k] = t * p[i*bs+k]; \ 157 p[i*bs+k] = t; \ 158 } \ 159 } \ 160 } 161 #define DEF_Pack(type,BS) \ 162 DEF_PackAddNoInit(type,BS) \ 163 static void CPPJoin3_(PackInit_,type,BS)(PetscSFBasicPack link) { \ 164 link->Pack = CPPJoin3_(Pack_,type,BS); \ 165 link->UnpackInsert = CPPJoin3_(UnpackInsert_,type,BS); \ 166 link->UnpackAdd = CPPJoin3_(UnpackAdd_,type,BS); \ 167 link->UnpackMult = CPPJoin3_(UnpackMult_,type,BS); \ 168 link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,type,BS); \ 169 link->FetchAndAdd = CPPJoin3_(FetchAndAdd_,type,BS); \ 170 link->FetchAndMult = CPPJoin3_(FetchAndMult_,type,BS); \ 171 link->unitbytes = sizeof(type); \ 172 } 173 /* Comparable types */ 174 #define DEF_PackCmp(type) \ 175 DEF_PackAddNoInit(type,1) \ 176 static void CPPJoin2(UnpackMax_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \ 177 type *u = (type*)unpacked; \ 178 const type *p = (const type*)packed; \ 179 PetscInt i; \ 180 for (i=0; i<n; i++) { \ 181 type v = u[idx[i]]; \ 182 u[idx[i]] = PetscMax(v,p[i]); \ 183 } \ 184 } \ 185 static void CPPJoin2(UnpackMin_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \ 186 type *u = (type*)unpacked; \ 187 const type *p = (const type*)packed; \ 188 PetscInt i; \ 189 for (i=0; i<n; i++) { \ 190 type v = u[idx[i]]; \ 191 u[idx[i]] = PetscMin(v,p[i]); \ 192 } \ 193 } \ 194 static void CPPJoin2(FetchAndMax_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \ 195 type *u = (type*)unpacked; \ 196 type *p = (type*)packed; \ 197 PetscInt i; \ 198 for (i=0; i<n; i++) { \ 199 PetscInt j = idx[i]; \ 200 type v = u[j]; \ 201 u[j] = PetscMax(v,p[i]); \ 202 p[i] = v; \ 203 } \ 204 } \ 205 static void CPPJoin2(FetchAndMin_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \ 206 type *u = (type*)unpacked; \ 207 type *p = (type*)packed; \ 208 PetscInt i; \ 209 for (i=0; i<n; i++) { \ 210 PetscInt j = idx[i]; \ 211 type v = u[j]; \ 212 u[j] = PetscMin(v,p[i]); \ 213 p[i] = v; \ 214 } \ 215 } \ 216 static void CPPJoin2(PackInit_,type)(PetscSFBasicPack link) { \ 217 link->Pack = CPPJoin3_(Pack_,type,1); \ 218 link->UnpackInsert = CPPJoin3_(UnpackInsert_,type,1); \ 219 link->UnpackAdd = CPPJoin3_(UnpackAdd_,type,1); \ 220 link->UnpackMax = CPPJoin2(UnpackMax_,type); \ 221 link->UnpackMin = CPPJoin2(UnpackMin_,type); \ 222 link->UnpackMult = CPPJoin3_(UnpackMult_,type,1); \ 223 link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,type,1); \ 224 link->FetchAndAdd = CPPJoin3_(FetchAndAdd_ ,type,1); \ 225 link->FetchAndMax = CPPJoin2(FetchAndMax_ ,type); \ 226 link->FetchAndMin = CPPJoin2(FetchAndMin_ ,type); \ 227 link->FetchAndMult = CPPJoin3_(FetchAndMult_,type,1); \ 228 link->unitbytes = sizeof(type); \ 229 } 230 231 /* Logical Types */ 232 #define DEF_PackLog(type) \ 233 static void CPPJoin2(UnpackLAND_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \ 234 type *u = (type*)unpacked; \ 235 const type *p = (const type*)packed; \ 236 PetscInt i; \ 237 for (i=0; i<n; i++) { \ 238 type v = u[idx[i]]; \ 239 u[idx[i]] = v && p[i]; \ 240 } \ 241 } \ 242 static void CPPJoin2(UnpackLOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \ 243 type *u = (type*)unpacked; \ 244 const type *p = (const type*)packed; \ 245 PetscInt i; \ 246 for (i=0; i<n; i++) { \ 247 type v = u[idx[i]]; \ 248 u[idx[i]] = v || p[i]; \ 249 } \ 250 } \ 251 static void CPPJoin2(UnpackLXOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \ 252 type *u = (type*)unpacked; \ 253 const type *p = (const type*)packed; \ 254 PetscInt i; \ 255 for (i=0; i<n; i++) { \ 256 type v = u[idx[i]]; \ 257 u[idx[i]] = (!v)!=(!p[i]); \ 258 } \ 259 } \ 260 static void CPPJoin2(FetchAndLAND_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \ 261 type *u = (type*)unpacked; \ 262 type *p = (type*)packed; \ 263 PetscInt i; \ 264 for (i=0; i<n; i++) { \ 265 PetscInt j = idx[i]; \ 266 type v = u[j]; \ 267 u[j] = v && p[i]; \ 268 p[i] = v; \ 269 } \ 270 } \ 271 static void CPPJoin2(FetchAndLOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \ 272 type *u = (type*)unpacked; \ 273 type *p = (type*)packed; \ 274 PetscInt i; \ 275 for (i=0; i<n; i++) { \ 276 PetscInt j = idx[i]; \ 277 type v = u[j]; \ 278 u[j] = v || p[i]; \ 279 p[i] = v; \ 280 } \ 281 } \ 282 static void CPPJoin2(FetchAndLXOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \ 283 type *u = (type*)unpacked; \ 284 type *p = (type*)packed; \ 285 PetscInt i; \ 286 for (i=0; i<n; i++) { \ 287 PetscInt j = idx[i]; \ 288 type v = u[j]; \ 289 u[j] = (!v)!=(!p[i]); \ 290 p[i] = v; \ 291 } \ 292 } \ 293 static void CPPJoin2(PackInit_Logical_,type)(PetscSFBasicPack link) { \ 294 link->UnpackLAND = CPPJoin2(UnpackLAND_,type); \ 295 link->UnpackLOR = CPPJoin2(UnpackLOR_,type); \ 296 link->UnpackLXOR = CPPJoin2(UnpackLXOR_,type); \ 297 link->FetchAndLAND = CPPJoin2(FetchAndLAND_,type); \ 298 link->FetchAndLOR = CPPJoin2(FetchAndLOR_,type); \ 299 link->FetchAndLXOR = CPPJoin2(FetchAndLXOR_,type); \ 300 } 301 302 303 /* Bitwise Types */ 304 #define DEF_PackBit(type) \ 305 static void CPPJoin2(UnpackBAND_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \ 306 type *u = (type*)unpacked; \ 307 const type *p = (const type*)packed; \ 308 PetscInt i; \ 309 for (i=0; i<n; i++) { \ 310 type v = u[idx[i]]; \ 311 u[idx[i]] = v & p[i]; \ 312 } \ 313 } \ 314 static void CPPJoin2(UnpackBOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \ 315 type *u = (type*)unpacked; \ 316 const type *p = (const type*)packed; \ 317 PetscInt i; \ 318 for (i=0; i<n; i++) { \ 319 type v = u[idx[i]]; \ 320 u[idx[i]] = v | p[i]; \ 321 } \ 322 } \ 323 static void CPPJoin2(UnpackBXOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \ 324 type *u = (type*)unpacked; \ 325 const type *p = (const type*)packed; \ 326 PetscInt i; \ 327 for (i=0; i<n; i++) { \ 328 type v = u[idx[i]]; \ 329 u[idx[i]] = v^p[i]; \ 330 } \ 331 } \ 332 static void CPPJoin2(FetchAndBAND_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \ 333 type *u = (type*)unpacked; \ 334 type *p = (type*)packed; \ 335 PetscInt i; \ 336 for (i=0; i<n; i++) { \ 337 PetscInt j = idx[i]; \ 338 type v = u[j]; \ 339 u[j] = v & p[i]; \ 340 p[i] = v; \ 341 } \ 342 } \ 343 static void CPPJoin2(FetchAndBOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \ 344 type *u = (type*)unpacked; \ 345 type *p = (type*)packed; \ 346 PetscInt i; \ 347 for (i=0; i<n; i++) { \ 348 PetscInt j = idx[i]; \ 349 type v = u[j]; \ 350 u[j] = v | p[i]; \ 351 p[i] = v; \ 352 } \ 353 } \ 354 static void CPPJoin2(FetchAndBXOR_,type)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \ 355 type *u = (type*)unpacked; \ 356 type *p = (type*)packed; \ 357 PetscInt i; \ 358 for (i=0; i<n; i++) { \ 359 PetscInt j = idx[i]; \ 360 type v = u[j]; \ 361 u[j] = v^p[i]; \ 362 p[i] = v; \ 363 } \ 364 } \ 365 static void CPPJoin2(PackInit_Bitwise_,type)(PetscSFBasicPack link) { \ 366 link->UnpackBAND = CPPJoin2(UnpackBAND_,type); \ 367 link->UnpackBOR = CPPJoin2(UnpackBOR_,type); \ 368 link->UnpackBXOR = CPPJoin2(UnpackBXOR_,type); \ 369 link->FetchAndBAND = CPPJoin2(FetchAndBAND_,type); \ 370 link->FetchAndBOR = CPPJoin2(FetchAndBOR_,type); \ 371 link->FetchAndBXOR = CPPJoin2(FetchAndBXOR_,type); \ 372 } 373 374 /* Pair types */ 375 #define CPPJoinloc_exp(base,op,t1,t2) base ## op ## loc_ ## t1 ## _ ## t2 376 #define CPPJoinloc(base,op,t1,t2) CPPJoinloc_exp(base,op,t1,t2) 377 #define PairType(type1,type2) CPPJoin3_(_pairtype_,type1,type2) 378 #define DEF_UnpackXloc(type1,type2,locname,op) \ 379 static void CPPJoinloc(Unpack,locname,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \ 380 PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked; \ 381 const PairType(type1,type2) *p = (const PairType(type1,type2)*)packed; \ 382 PetscInt i; \ 383 for (i=0; i<n; i++) { \ 384 PetscInt j = idx[i]; \ 385 if (p[i].a op u[j].a) { \ 386 u[j].a = p[i].a; \ 387 u[j].b = p[i].b; \ 388 } else if (u[j].a == p[i].a) { \ 389 u[j].b = PetscMin(u[j].b,p[i].b); \ 390 } \ 391 } \ 392 } \ 393 static void CPPJoinloc(FetchAnd,locname,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \ 394 PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked; \ 395 PairType(type1,type2) *p = (PairType(type1,type2)*)packed; \ 396 PetscInt i; \ 397 for (i=0; i<n; i++) { \ 398 PetscInt j = idx[i]; \ 399 PairType(type1,type2) v; \ 400 v.a = u[j].a; \ 401 v.b = u[j].b; \ 402 if (p[i].a op u[j].a) { \ 403 u[j].a = p[i].a; \ 404 u[j].b = p[i].b; \ 405 } else if (u[j].a == p[i].a) { \ 406 u[j].b = PetscMin(u[j].b,p[i].b); \ 407 } \ 408 p[i].a = v.a; \ 409 p[i].b = v.b; \ 410 } \ 411 } 412 #define DEF_PackPair(type1,type2) \ 413 typedef struct {type1 a; type2 b;} PairType(type1,type2); \ 414 static void CPPJoin3_(Pack_,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,const void *unpacked,void *packed) { \ 415 const PairType(type1,type2) *u = (const PairType(type1,type2)*)unpacked; \ 416 PairType(type1,type2) *p = (PairType(type1,type2)*)packed; \ 417 PetscInt i; \ 418 for (i=0; i<n; i++) { \ 419 p[i].a = u[idx[i]].a; \ 420 p[i].b = u[idx[i]].b; \ 421 } \ 422 } \ 423 static void CPPJoin3_(UnpackInsert_,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \ 424 PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked; \ 425 const PairType(type1,type2) *p = (const PairType(type1,type2)*)packed; \ 426 PetscInt i; \ 427 for (i=0; i<n; i++) { \ 428 u[idx[i]].a = p[i].a; \ 429 u[idx[i]].b = p[i].b; \ 430 } \ 431 } \ 432 static void CPPJoin3_(UnpackAdd_,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,const void *packed) { \ 433 PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked; \ 434 const PairType(type1,type2) *p = (const PairType(type1,type2)*)packed; \ 435 PetscInt i; \ 436 for (i=0; i<n; i++) { \ 437 u[idx[i]].a += p[i].a; \ 438 u[idx[i]].b += p[i].b; \ 439 } \ 440 } \ 441 static void CPPJoin3_(FetchAndInsert_,type1,type2)(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \ 442 PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked; \ 443 PairType(type1,type2) *p = (PairType(type1,type2)*)packed; \ 444 PetscInt i; \ 445 for (i=0; i<n; i++) { \ 446 PetscInt j = idx[i]; \ 447 PairType(type1,type2) v; \ 448 v.a = u[j].a; \ 449 v.b = u[j].b; \ 450 u[j].a = p[i].a; \ 451 u[j].b = p[i].b; \ 452 p[i].a = v.a; \ 453 p[i].b = v.b; \ 454 } \ 455 } \ 456 static void FetchAndAdd_ ## type1 ## _ ## type2(PetscInt n,PetscInt bs,const PetscInt *idx,void *unpacked,void *packed) { \ 457 PairType(type1,type2) *u = (PairType(type1,type2)*)unpacked; \ 458 PairType(type1,type2) *p = (PairType(type1,type2)*)packed; \ 459 PetscInt i; \ 460 for (i=0; i<n; i++) { \ 461 PetscInt j = idx[i]; \ 462 PairType(type1,type2) v; \ 463 v.a = u[j].a; \ 464 v.b = u[j].b; \ 465 u[j].a = v.a + p[i].a; \ 466 u[j].b = v.b + p[i].b; \ 467 p[i].a = v.a; \ 468 p[i].b = v.b; \ 469 } \ 470 } \ 471 DEF_UnpackXloc(type1,type2,Max,>) \ 472 DEF_UnpackXloc(type1,type2,Min,<) \ 473 static void CPPJoin3_(PackInit_,type1,type2)(PetscSFBasicPack link) { \ 474 link->Pack = CPPJoin3_(Pack_,type1,type2); \ 475 link->UnpackInsert = CPPJoin3_(UnpackInsert_,type1,type2); \ 476 link->UnpackAdd = CPPJoin3_(UnpackAdd_,type1,type2); \ 477 link->UnpackMaxloc = CPPJoin3_(UnpackMaxloc_,type1,type2); \ 478 link->UnpackMinloc = CPPJoin3_(UnpackMinloc_,type1,type2); \ 479 link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,type1,type2); \ 480 link->FetchAndAdd = CPPJoin3_(FetchAndAdd_,type1,type2); \ 481 link->FetchAndMaxloc = CPPJoin3_(FetchAndMaxloc_,type1,type2); \ 482 link->FetchAndMinloc = CPPJoin3_(FetchAndMinloc_,type1,type2); \ 483 link->unitbytes = sizeof(PairType(type1,type2)); \ 484 } 485 486 /* Currently only dumb blocks of data */ 487 #define BlockType(unit,count) CPPJoin3_(_blocktype_,unit,count) 488 #define DEF_Block(unit,count) \ 489 typedef struct {unit v[count];} BlockType(unit,count); \ 490 DEF_PackNoInit(BlockType(unit,count),1) \ 491 static void CPPJoin3_(PackInit_block_,unit,count)(PetscSFBasicPack link) { \ 492 link->Pack = CPPJoin3_(Pack_,BlockType(unit,count),1); \ 493 link->UnpackInsert = CPPJoin3_(UnpackInsert_,BlockType(unit,count),1); \ 494 link->FetchAndInsert = CPPJoin3_(FetchAndInsert_,BlockType(unit,count),1); \ 495 link->unitbytes = sizeof(BlockType(unit,count)); \ 496 } 497 498 DEF_PackCmp(int) 499 DEF_PackBit(int) 500 DEF_PackLog(int) 501 DEF_PackCmp(PetscInt) 502 DEF_PackBit(PetscInt) 503 DEF_PackLog(PetscInt) 504 DEF_Pack(PetscInt,2) 505 DEF_Pack(PetscInt,3) 506 DEF_Pack(PetscInt,4) 507 DEF_Pack(PetscInt,5) 508 DEF_Pack(PetscInt,7) 509 DEF_PackCmp(PetscReal) 510 DEF_PackLog(PetscReal) 511 DEF_Pack(PetscReal,2) 512 DEF_Pack(PetscReal,3) 513 DEF_Pack(PetscReal,4) 514 DEF_Pack(PetscReal,5) 515 DEF_Pack(PetscReal,7) 516 #if defined(PETSC_HAVE_COMPLEX) 517 DEF_Pack(PetscComplex,1) 518 DEF_Pack(PetscComplex,2) 519 DEF_Pack(PetscComplex,3) 520 DEF_Pack(PetscComplex,4) 521 DEF_Pack(PetscComplex,5) 522 DEF_Pack(PetscComplex,7) 523 #endif 524 DEF_PackPair(int,int) 525 DEF_PackPair(PetscInt,PetscInt) 526 DEF_Block(int,1) 527 DEF_Block(int,2) 528 DEF_Block(int,3) 529 DEF_Block(int,4) 530 DEF_Block(int,5) 531 DEF_Block(int,6) 532 DEF_Block(int,7) 533 DEF_Block(int,8) 534 535 static PetscErrorCode PetscSFSetUp_Basic(PetscSF sf) 536 { 537 PetscSF_Basic *bas = (PetscSF_Basic*)sf->data; 538 PetscErrorCode ierr; 539 PetscInt *rlengths,*ilengths,i; 540 PetscMPIInt rank,niranks,*iranks; 541 MPI_Comm comm; 542 MPI_Group group; 543 MPI_Request *rootreqs,*leafreqs; 544 545 PetscFunctionBegin; 546 ierr = MPI_Comm_group(PETSC_COMM_SELF,&group);CHKERRQ(ierr); 547 ierr = PetscSFSetUpRanks(sf,group);CHKERRQ(ierr); 548 ierr = MPI_Group_free(&group);CHKERRQ(ierr); 549 ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr); 550 ierr = PetscObjectGetNewTag((PetscObject)sf,&bas->tag);CHKERRQ(ierr); 551 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 552 /* 553 * Inform roots about how many leaves and from which ranks 554 */ 555 ierr = PetscMalloc1(sf->nranks,&rlengths);CHKERRQ(ierr); 556 /* Determine number, sending ranks, and length of incoming */ 557 for (i=0; i<sf->nranks; i++) { 558 rlengths[i] = sf->roffset[i+1] - sf->roffset[i]; /* Number of roots referenced by my leaves; for rank sf->ranks[i] */ 559 } 560 ierr = PetscCommBuildTwoSided(comm,1,MPIU_INT,sf->nranks-sf->ndranks,sf->ranks+sf->ndranks,rlengths+sf->ndranks,&niranks,&iranks,(void**)&ilengths);CHKERRQ(ierr); 561 562 /* Sort iranks. See use of VecScatterGetRemoteOrdered_Private() in MatGetBrowsOfAoCols_MPIAIJ() on why. 563 We could sort ranks there at the price of allocating extra working arrays. Presumably, niranks is 564 small and the sorting is cheap. 565 */ 566 ierr = PetscSortMPIIntWithIntArray(niranks,iranks,ilengths);CHKERRQ(ierr); 567 568 /* Partition into distinguished and non-distinguished incoming ranks */ 569 bas->ndiranks = sf->ndranks; 570 bas->niranks = bas->ndiranks + niranks; 571 ierr = PetscMalloc2(bas->niranks,&bas->iranks,bas->niranks+1,&bas->ioffset);CHKERRQ(ierr); 572 bas->ioffset[0] = 0; 573 for (i=0; i<bas->ndiranks; i++) { 574 bas->iranks[i] = sf->ranks[i]; 575 bas->ioffset[i+1] = bas->ioffset[i] + rlengths[i]; 576 } 577 if (bas->ndiranks > 1 || (bas->ndiranks == 1 && bas->iranks[0] != rank)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Broken setup for shared ranks"); 578 for ( ; i<bas->niranks; i++) { 579 bas->iranks[i] = iranks[i-bas->ndiranks]; 580 bas->ioffset[i+1] = bas->ioffset[i] + ilengths[i-bas->ndiranks]; 581 } 582 bas->itotal = bas->ioffset[i]; 583 ierr = PetscFree(rlengths);CHKERRQ(ierr); 584 ierr = PetscFree(iranks);CHKERRQ(ierr); 585 ierr = PetscFree(ilengths);CHKERRQ(ierr); 586 587 /* Send leaf identities to roots */ 588 ierr = PetscMalloc1(bas->itotal,&bas->irootloc);CHKERRQ(ierr); 589 ierr = PetscMalloc2(bas->niranks-bas->ndiranks,&rootreqs,sf->nranks-sf->ndranks,&leafreqs);CHKERRQ(ierr); 590 for (i=bas->ndiranks; i<bas->niranks; i++) { 591 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); 592 } 593 for (i=0; i<sf->nranks; i++) { 594 PetscMPIInt npoints; 595 ierr = PetscMPIIntCast(sf->roffset[i+1] - sf->roffset[i],&npoints);CHKERRQ(ierr); 596 if (i < sf->ndranks) { 597 if (sf->ranks[i] != rank) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot interpret distinguished leaf rank"); 598 if (bas->iranks[0] != rank) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot interpret distinguished root rank"); 599 if (npoints != bas->ioffset[1]-bas->ioffset[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Distinguished rank exchange has mismatched lengths"); 600 ierr = PetscMemcpy(bas->irootloc+bas->ioffset[0],sf->rremote+sf->roffset[i],npoints*sizeof(bas->irootloc[0]));CHKERRQ(ierr); 601 continue; 602 } 603 ierr = MPI_Isend(sf->rremote+sf->roffset[i],npoints,MPIU_INT,sf->ranks[i],bas->tag,comm,&leafreqs[i-sf->ndranks]);CHKERRQ(ierr); 604 } 605 ierr = MPI_Waitall(bas->niranks-bas->ndiranks,rootreqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 606 ierr = MPI_Waitall(sf->nranks-sf->ndranks,leafreqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 607 ierr = PetscFree2(rootreqs,leafreqs);CHKERRQ(ierr); 608 PetscFunctionReturn(0); 609 } 610 611 static PetscErrorCode PetscSFBasicPackTypeSetup(PetscSFBasicPack link,MPI_Datatype unit) 612 { 613 PetscErrorCode ierr; 614 PetscBool isInt,isPetscInt,isPetscReal,is2Int,is2PetscInt; 615 PetscInt nPetscIntContig,nPetscRealContig; 616 PetscMPIInt ni,na,nd,combiner; 617 #if defined(PETSC_HAVE_COMPLEX) 618 PetscBool isPetscComplex; 619 PetscInt nPetscComplexContig; 620 #endif 621 622 PetscFunctionBegin; 623 ierr = MPIPetsc_Type_compare(unit,MPI_INT,&isInt);CHKERRQ(ierr); 624 ierr = MPIPetsc_Type_compare(unit,MPIU_INT,&isPetscInt);CHKERRQ(ierr); 625 ierr = MPIPetsc_Type_compare_contig(unit,MPIU_INT,&nPetscIntContig);CHKERRQ(ierr); 626 ierr = MPIPetsc_Type_compare(unit,MPIU_REAL,&isPetscReal);CHKERRQ(ierr); 627 ierr = MPIPetsc_Type_compare_contig(unit,MPIU_REAL,&nPetscRealContig);CHKERRQ(ierr); 628 #if defined(PETSC_HAVE_COMPLEX) 629 ierr = MPIPetsc_Type_compare(unit,MPIU_COMPLEX,&isPetscComplex);CHKERRQ(ierr); 630 ierr = MPIPetsc_Type_compare_contig(unit,MPIU_COMPLEX,&nPetscComplexContig);CHKERRQ(ierr); 631 #endif 632 ierr = MPIPetsc_Type_compare(unit,MPI_2INT,&is2Int);CHKERRQ(ierr); 633 ierr = MPIPetsc_Type_compare(unit,MPIU_2INT,&is2PetscInt);CHKERRQ(ierr); 634 ierr = MPI_Type_get_envelope(unit,&ni,&na,&nd,&combiner);CHKERRQ(ierr); 635 link->isbuiltin = (combiner == MPI_COMBINER_NAMED) ? PETSC_TRUE : PETSC_FALSE; 636 link->bs = 1; 637 if (isInt) {PackInit_int(link); PackInit_Logical_int(link); PackInit_Bitwise_int(link);} 638 else if (isPetscInt) {PackInit_PetscInt(link); PackInit_Logical_PetscInt(link); PackInit_Bitwise_PetscInt(link);} 639 else if (isPetscReal) {PackInit_PetscReal(link); PackInit_Logical_PetscReal(link);} 640 #if defined(PETSC_HAVE_COMPLEX) 641 else if (isPetscComplex) PackInit_PetscComplex_1(link); 642 #endif 643 else if (is2Int) PackInit_int_int(link); 644 else if (is2PetscInt) PackInit_PetscInt_PetscInt(link); 645 else if (nPetscIntContig) { 646 if (nPetscIntContig%7 == 0) PackInit_PetscInt_7(link); 647 else if (nPetscIntContig%5 == 0) PackInit_PetscInt_5(link); 648 else if (nPetscIntContig%4 == 0) PackInit_PetscInt_4(link); 649 else if (nPetscIntContig%3 == 0) PackInit_PetscInt_3(link); 650 else if (nPetscIntContig%2 == 0) PackInit_PetscInt_2(link); 651 else PackInit_PetscInt(link); 652 link->bs = nPetscIntContig; 653 link->unitbytes *= nPetscIntContig; 654 } else if (nPetscRealContig) { 655 if (nPetscRealContig%7 == 0) PackInit_PetscReal_7(link); 656 else if (nPetscRealContig%5 == 0) PackInit_PetscReal_5(link); 657 else if (nPetscRealContig%4 == 0) PackInit_PetscReal_4(link); 658 else if (nPetscRealContig%3 == 0) PackInit_PetscReal_3(link); 659 else if (nPetscRealContig%2 == 0) PackInit_PetscReal_2(link); 660 else PackInit_PetscReal(link); 661 link->bs = nPetscRealContig; 662 link->unitbytes *= nPetscRealContig; 663 #if defined(PETSC_HAVE_COMPLEX) 664 } else if (nPetscComplexContig) { 665 if (nPetscComplexContig%7 == 0) PackInit_PetscComplex_7(link); 666 else if (nPetscComplexContig%5 == 0) PackInit_PetscComplex_5(link); 667 else if (nPetscComplexContig%4 == 0) PackInit_PetscComplex_4(link); 668 else if (nPetscComplexContig%3 == 0) PackInit_PetscComplex_3(link); 669 else if (nPetscComplexContig%2 == 0) PackInit_PetscComplex_2(link); 670 else PackInit_PetscComplex_1(link); 671 link->bs = nPetscComplexContig; 672 link->unitbytes *= nPetscComplexContig; 673 #endif 674 } else { 675 MPI_Aint lb,bytes; 676 ierr = MPI_Type_get_extent(unit,&lb,&bytes);CHKERRQ(ierr); 677 if (lb != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Datatype with nonzero lower bound %ld\n",(long)lb); 678 if (bytes % sizeof(int)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for type size not divisible by %D",sizeof(int)); 679 switch (bytes / sizeof(int)) { 680 case 1: PackInit_block_int_1(link); break; 681 case 2: PackInit_block_int_2(link); break; 682 case 3: PackInit_block_int_3(link); break; 683 case 4: PackInit_block_int_4(link); break; 684 case 5: PackInit_block_int_5(link); break; 685 case 6: PackInit_block_int_6(link); break; 686 case 7: PackInit_block_int_7(link); break; 687 case 8: PackInit_block_int_8(link); break; 688 default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for arbitrary block sizes"); 689 } 690 } 691 if (link->isbuiltin) link->unit = unit; /* builtin datatypes are common. Make it fast */ 692 else {ierr = MPI_Type_dup(unit,&link->unit);CHKERRQ(ierr);} 693 PetscFunctionReturn(0); 694 } 695 696 static PetscErrorCode PetscSFBasicPackGetUnpackOp(PetscSF sf,PetscSFBasicPack link,MPI_Op op,void (**UnpackOp)(PetscInt,PetscInt,const PetscInt*,void*,const void*)) 697 { 698 PetscFunctionBegin; 699 *UnpackOp = NULL; 700 if (op == MPIU_REPLACE) *UnpackOp = link->UnpackInsert; 701 else if (op == MPI_SUM || op == MPIU_SUM) *UnpackOp = link->UnpackAdd; 702 else if (op == MPI_PROD) *UnpackOp = link->UnpackMult; 703 else if (op == MPI_MAX || op == MPIU_MAX) *UnpackOp = link->UnpackMax; 704 else if (op == MPI_MIN || op == MPIU_MIN) *UnpackOp = link->UnpackMin; 705 else if (op == MPI_LAND) *UnpackOp = link->UnpackLAND; 706 else if (op == MPI_BAND) *UnpackOp = link->UnpackBAND; 707 else if (op == MPI_LOR) *UnpackOp = link->UnpackLOR; 708 else if (op == MPI_BOR) *UnpackOp = link->UnpackBOR; 709 else if (op == MPI_LXOR) *UnpackOp = link->UnpackLXOR; 710 else if (op == MPI_BXOR) *UnpackOp = link->UnpackBXOR; 711 else if (op == MPI_MAXLOC) *UnpackOp = link->UnpackMaxloc; 712 else if (op == MPI_MINLOC) *UnpackOp = link->UnpackMinloc; 713 else *UnpackOp = NULL; 714 PetscFunctionReturn(0); 715 } 716 static PetscErrorCode PetscSFBasicPackGetFetchAndOp(PetscSF sf,PetscSFBasicPack link,MPI_Op op,void (**FetchAndOp)(PetscInt,PetscInt,const PetscInt*,void*,void*)) 717 { 718 PetscFunctionBegin; 719 *FetchAndOp = NULL; 720 if (op == MPIU_REPLACE) *FetchAndOp = link->FetchAndInsert; 721 else if (op == MPI_SUM || op == MPIU_SUM) *FetchAndOp = link->FetchAndAdd; 722 else if (op == MPI_MAX || op == MPIU_MAX) *FetchAndOp = link->FetchAndMax; 723 else if (op == MPI_MIN || op == MPIU_MIN) *FetchAndOp = link->FetchAndMin; 724 else if (op == MPI_MAXLOC) *FetchAndOp = link->FetchAndMaxloc; 725 else if (op == MPI_MINLOC) *FetchAndOp = link->FetchAndMinloc; 726 else if (op == MPI_PROD) *FetchAndOp = link->FetchAndMult; 727 else if (op == MPI_LAND) *FetchAndOp = link->FetchAndLAND; 728 else if (op == MPI_BAND) *FetchAndOp = link->FetchAndBAND; 729 else if (op == MPI_LOR) *FetchAndOp = link->FetchAndLOR; 730 else if (op == MPI_BOR) *FetchAndOp = link->FetchAndBOR; 731 else if (op == MPI_LXOR) *FetchAndOp = link->FetchAndLXOR; 732 else if (op == MPI_BXOR) *FetchAndOp = link->FetchAndBXOR; 733 else SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_SUP,"No support for MPI_Op"); 734 PetscFunctionReturn(0); 735 } 736 737 static PetscErrorCode PetscSFBasicPackGetReqs(PetscSF sf,PetscSFBasicPack link,MPI_Request **rootreqs,MPI_Request **leafreqs) 738 { 739 PetscSF_Basic *bas = (PetscSF_Basic*)sf->data; 740 741 PetscFunctionBegin; 742 if (rootreqs) *rootreqs = link->requests; 743 if (leafreqs) *leafreqs = link->requests + (bas->niranks - bas->ndiranks); 744 PetscFunctionReturn(0); 745 } 746 747 static PetscErrorCode PetscSFBasicPackWaitall(PetscSF sf,PetscSFBasicPack link) 748 { 749 PetscSF_Basic *bas = (PetscSF_Basic*)sf->data; 750 PetscErrorCode ierr; 751 752 PetscFunctionBegin; 753 ierr = MPI_Waitall(bas->niranks+sf->nranks-(bas->ndiranks+sf->ndranks),link->requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 754 PetscFunctionReturn(0); 755 } 756 757 static PetscErrorCode PetscSFBasicGetRootInfo(PetscSF sf,PetscInt *nrootranks,PetscInt *ndrootranks,const PetscMPIInt **rootranks,const PetscInt **rootoffset,const PetscInt **rootloc) 758 { 759 PetscSF_Basic *bas = (PetscSF_Basic*)sf->data; 760 761 PetscFunctionBegin; 762 if (nrootranks) *nrootranks = bas->niranks; 763 if (ndrootranks) *ndrootranks = bas->ndiranks; 764 if (rootranks) *rootranks = bas->iranks; 765 if (rootoffset) *rootoffset = bas->ioffset; 766 if (rootloc) *rootloc = bas->irootloc; 767 PetscFunctionReturn(0); 768 } 769 770 static PetscErrorCode PetscSFBasicGetLeafInfo(PetscSF sf,PetscInt *nleafranks,PetscInt *ndleafranks,const PetscMPIInt **leafranks,const PetscInt **leafoffset,const PetscInt **leafloc) 771 { 772 PetscFunctionBegin; 773 if (nleafranks) *nleafranks = sf->nranks; 774 if (ndleafranks) *ndleafranks = sf->ndranks; 775 if (leafranks) *leafranks = sf->ranks; 776 if (leafoffset) *leafoffset = sf->roffset; 777 if (leafloc) *leafloc = sf->rmine; 778 PetscFunctionReturn(0); 779 } 780 781 static PetscErrorCode PetscSFBasicGetPack(PetscSF sf,MPI_Datatype unit,const void *key,PetscSFBasicPack *mylink) 782 { 783 PetscSF_Basic *bas = (PetscSF_Basic*)sf->data; 784 PetscErrorCode ierr; 785 PetscSFBasicPack link,*p; 786 PetscInt nrootranks,ndrootranks,nleafranks,ndleafranks,i; 787 const PetscInt *rootoffset,*leafoffset; 788 789 PetscFunctionBegin; 790 /* Look for types in cache */ 791 for (p=&bas->avail; (link=*p); p=&link->next) { 792 PetscBool match; 793 ierr = MPIPetsc_Type_compare(unit,link->unit,&match);CHKERRQ(ierr); 794 if (match) { 795 *p = link->next; /* Remove from available list */ 796 goto found; 797 } 798 } 799 800 /* Create new composite types for each send rank */ 801 ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,&ndrootranks,NULL,&rootoffset,NULL);CHKERRQ(ierr); 802 ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,&ndleafranks,NULL,&leafoffset,NULL);CHKERRQ(ierr); 803 ierr = PetscNew(&link);CHKERRQ(ierr); 804 ierr = PetscSFBasicPackTypeSetup(link,unit);CHKERRQ(ierr); 805 ierr = PetscMalloc2(nrootranks,&link->root,nleafranks,&link->leaf);CHKERRQ(ierr); 806 for (i=0; i<nrootranks; i++) { 807 ierr = PetscMalloc((rootoffset[i+1]-rootoffset[i])*link->unitbytes,&link->root[i]);CHKERRQ(ierr); 808 } 809 for (i=0; i<nleafranks; i++) { 810 if (i < ndleafranks) { /* Leaf buffers for distinguished ranks are pointers directly into root buffers */ 811 if (ndrootranks != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot match distinguished ranks"); 812 link->leaf[i] = link->root[0]; 813 continue; 814 } 815 ierr = PetscMalloc((leafoffset[i+1]-leafoffset[i])*link->unitbytes,&link->leaf[i]);CHKERRQ(ierr); 816 } 817 ierr = PetscCalloc1(nrootranks+nleafranks,&link->requests);CHKERRQ(ierr); 818 819 found: 820 link->key = key; 821 link->next = bas->inuse; 822 bas->inuse = link; 823 824 *mylink = link; 825 PetscFunctionReturn(0); 826 } 827 828 static PetscErrorCode PetscSFBasicGetPackInUse(PetscSF sf,MPI_Datatype unit,const void *key,PetscCopyMode cmode,PetscSFBasicPack *mylink) 829 { 830 PetscSF_Basic *bas = (PetscSF_Basic*)sf->data; 831 PetscErrorCode ierr; 832 PetscSFBasicPack link,*p; 833 834 PetscFunctionBegin; 835 /* Look for types in cache */ 836 for (p=&bas->inuse; (link=*p); p=&link->next) { 837 PetscBool match; 838 ierr = MPIPetsc_Type_compare(unit,link->unit,&match);CHKERRQ(ierr); 839 if (match && (key == link->key)) { 840 switch (cmode) { 841 case PETSC_OWN_POINTER: *p = link->next; break; /* Remove from inuse list */ 842 case PETSC_USE_POINTER: break; 843 default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"invalid cmode"); 844 } 845 *mylink = link; 846 PetscFunctionReturn(0); 847 } 848 } 849 SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Could not find pack"); 850 PetscFunctionReturn(0); 851 } 852 853 static PetscErrorCode PetscSFBasicReclaimPack(PetscSF sf,PetscSFBasicPack *link) 854 { 855 PetscSF_Basic *bas = (PetscSF_Basic*)sf->data; 856 857 PetscFunctionBegin; 858 (*link)->key = NULL; 859 (*link)->next = bas->avail; 860 bas->avail = *link; 861 *link = NULL; 862 PetscFunctionReturn(0); 863 } 864 865 static PetscErrorCode PetscSFSetFromOptions_Basic(PetscOptionItems *PetscOptionsObject,PetscSF sf) 866 { 867 PetscErrorCode ierr; 868 869 PetscFunctionBegin; 870 ierr = PetscOptionsHead(PetscOptionsObject,"PetscSF Basic options");CHKERRQ(ierr); 871 ierr = PetscOptionsTail();CHKERRQ(ierr); 872 PetscFunctionReturn(0); 873 } 874 875 static PetscErrorCode PetscSFReset_Basic(PetscSF sf) 876 { 877 PetscSF_Basic *bas = (PetscSF_Basic*)sf->data; 878 PetscErrorCode ierr; 879 PetscSFBasicPack link,next; 880 881 PetscFunctionBegin; 882 if (bas->inuse) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Outstanding operation has not been completed"); 883 ierr = PetscFree2(bas->iranks,bas->ioffset);CHKERRQ(ierr); 884 ierr = PetscFree(bas->irootloc);CHKERRQ(ierr); 885 for (link=bas->avail; link; link=next) { 886 PetscInt i; 887 next = link->next; 888 if (!link->isbuiltin) {ierr = MPI_Type_free(&link->unit);CHKERRQ(ierr);} 889 for (i=0; i<bas->niranks; i++) {ierr = PetscFree(link->root[i]);CHKERRQ(ierr);} 890 for (i=sf->ndranks; i<sf->nranks; i++) {ierr = PetscFree(link->leaf[i]);CHKERRQ(ierr);} /* Free only non-distinguished leaf buffers */ 891 ierr = PetscFree2(link->root,link->leaf);CHKERRQ(ierr); 892 ierr = PetscFree(link->requests);CHKERRQ(ierr); 893 ierr = PetscFree(link);CHKERRQ(ierr); 894 } 895 bas->avail = NULL; 896 PetscFunctionReturn(0); 897 } 898 899 static PetscErrorCode PetscSFDestroy_Basic(PetscSF sf) 900 { 901 PetscErrorCode ierr; 902 903 PetscFunctionBegin; 904 ierr = PetscSFReset_Basic(sf);CHKERRQ(ierr); 905 ierr = PetscFree(sf->data);CHKERRQ(ierr); 906 PetscFunctionReturn(0); 907 } 908 909 static PetscErrorCode PetscSFView_Basic(PetscSF sf,PetscViewer viewer) 910 { 911 /* PetscSF_Basic *bas = (PetscSF_Basic*)sf->data; */ 912 PetscErrorCode ierr; 913 PetscBool iascii; 914 915 PetscFunctionBegin; 916 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 917 if (iascii) { 918 ierr = PetscViewerASCIIPrintf(viewer," sort=%s\n",sf->rankorder ? "rank-order" : "unordered");CHKERRQ(ierr); 919 } 920 PetscFunctionReturn(0); 921 } 922 923 static PetscErrorCode PetscSFBcastAndOpBegin_Basic(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op) 924 { 925 PetscSF_Basic *bas = (PetscSF_Basic*)sf->data; 926 PetscErrorCode ierr; 927 PetscSFBasicPack link; 928 PetscInt i,nrootranks,ndrootranks,nleafranks,ndleafranks; 929 const PetscInt *rootoffset,*leafoffset,*rootloc,*leafloc; 930 const PetscMPIInt *rootranks,*leafranks; 931 MPI_Request *rootreqs,*leafreqs; 932 933 PetscFunctionBegin; 934 ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,&ndrootranks,&rootranks,&rootoffset,&rootloc);CHKERRQ(ierr); 935 ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,&ndleafranks,&leafranks,&leafoffset,&leafloc);CHKERRQ(ierr); 936 ierr = PetscSFBasicGetPack(sf,unit,rootdata,&link);CHKERRQ(ierr); 937 938 ierr = PetscSFBasicPackGetReqs(sf,link,&rootreqs,&leafreqs);CHKERRQ(ierr); 939 /* Eagerly post leaf receives, but only from non-distinguished ranks -- distinguished ranks will receive via shared memory */ 940 for (i=ndleafranks; i<nleafranks; i++) { 941 PetscMPIInt n = leafoffset[i+1] - leafoffset[i]; 942 ierr = MPI_Irecv(link->leaf[i],n,unit,leafranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&leafreqs[i-ndleafranks]);CHKERRQ(ierr); 943 } 944 /* Pack and send root data */ 945 for (i=0; i<nrootranks; i++) { 946 PetscMPIInt n = rootoffset[i+1] - rootoffset[i]; 947 void *packstart = link->root[i]; 948 (*link->Pack)(n,link->bs,rootloc+rootoffset[i],rootdata,packstart); 949 if (i < ndrootranks) continue; /* shared memory */ 950 ierr = MPI_Isend(packstart,n,unit,rootranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&rootreqs[i-ndrootranks]);CHKERRQ(ierr); 951 } 952 PetscFunctionReturn(0); 953 } 954 955 PetscErrorCode PetscSFBcastAndOpEnd_Basic(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op) 956 { 957 PetscErrorCode ierr; 958 PetscSFBasicPack link; 959 PetscInt i,nleafranks,ndleafranks; 960 const PetscInt *leafoffset,*leafloc; 961 void (*UnpackOp)(PetscInt,PetscInt,const PetscInt*,void*,const void*); 962 PetscMPIInt typesize = -1; 963 964 PetscFunctionBegin; 965 ierr = PetscSFBasicGetPackInUse(sf,unit,rootdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr); 966 ierr = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr); 967 ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,&ndleafranks,NULL,&leafoffset,&leafloc);CHKERRQ(ierr); 968 ierr = PetscSFBasicPackGetUnpackOp(sf,link,op,&UnpackOp);CHKERRQ(ierr); 969 970 if (UnpackOp) { typesize = link->unitbytes; } 971 else { ierr = MPI_Type_size(unit,&typesize);CHKERRQ(ierr); } 972 973 for (i=0; i<nleafranks; i++) { 974 PetscMPIInt n = leafoffset[i+1] - leafoffset[i]; 975 char *packstart = (char *) link->leaf[i]; 976 if (UnpackOp) { (*UnpackOp)(n,link->bs,leafloc+leafoffset[i],leafdata,(const void *)packstart); } 977 #if defined(PETSC_HAVE_MPI_REDUCE_LOCAL) 978 else if (n) { /* the op should be defined to operate on the whole datatype, so we ignore link->bs */ 979 PetscInt j; 980 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); } 981 } 982 #else 983 else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No unpacking reduction operation for this MPI_Op"); 984 #endif 985 } 986 987 ierr = PetscSFBasicReclaimPack(sf,&link);CHKERRQ(ierr); 988 PetscFunctionReturn(0); 989 } 990 991 /* Send from roots to leaves */ 992 static PetscErrorCode PetscSFBcastBegin_Basic(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata) 993 { 994 PetscErrorCode ierr; 995 996 PetscFunctionBegin; 997 ierr = PetscSFBcastAndOpBegin_Basic(sf,unit,rootdata,leafdata,MPI_REPLACE);CHKERRQ(ierr); 998 PetscFunctionReturn(0); 999 } 1000 1001 PetscErrorCode PetscSFBcastEnd_Basic(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata) 1002 { 1003 PetscErrorCode ierr; 1004 1005 PetscFunctionBegin; 1006 ierr = PetscSFBcastAndOpEnd_Basic(sf,unit,rootdata,leafdata,MPI_REPLACE);CHKERRQ(ierr); 1007 PetscFunctionReturn(0); 1008 } 1009 1010 /* leaf -> root with reduction */ 1011 PetscErrorCode PetscSFReduceBegin_Basic(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op) 1012 { 1013 PetscSF_Basic *bas = (PetscSF_Basic*)sf->data; 1014 PetscSFBasicPack link; 1015 PetscErrorCode ierr; 1016 PetscInt i,nrootranks,ndrootranks,nleafranks,ndleafranks; 1017 const PetscInt *rootoffset,*leafoffset,*rootloc,*leafloc; 1018 const PetscMPIInt *rootranks,*leafranks; 1019 MPI_Request *rootreqs,*leafreqs; 1020 1021 PetscFunctionBegin; 1022 ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,&ndrootranks,&rootranks,&rootoffset,&rootloc);CHKERRQ(ierr); 1023 ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,&ndleafranks,&leafranks,&leafoffset,&leafloc);CHKERRQ(ierr); 1024 ierr = PetscSFBasicGetPack(sf,unit,leafdata,&link);CHKERRQ(ierr); 1025 1026 ierr = PetscSFBasicPackGetReqs(sf,link,&rootreqs,&leafreqs);CHKERRQ(ierr); 1027 /* Eagerly post root receives for non-distinguished ranks */ 1028 for (i=ndrootranks; i<nrootranks; i++) { 1029 PetscMPIInt n = rootoffset[i+1] - rootoffset[i]; 1030 ierr = MPI_Irecv(link->root[i],n,unit,rootranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&rootreqs[i-ndrootranks]);CHKERRQ(ierr); 1031 } 1032 /* Pack and send leaf data */ 1033 for (i=0; i<nleafranks; i++) { 1034 PetscMPIInt n = leafoffset[i+1] - leafoffset[i]; 1035 void *packstart = link->leaf[i]; 1036 (*link->Pack)(n,link->bs,leafloc+leafoffset[i],leafdata,packstart); 1037 if (i < ndleafranks) continue; /* shared memory */ 1038 ierr = MPI_Isend(packstart,n,unit,leafranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&leafreqs[i-ndleafranks]);CHKERRQ(ierr); 1039 } 1040 PetscFunctionReturn(0); 1041 } 1042 1043 static PetscErrorCode PetscSFReduceEnd_Basic(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op) 1044 { 1045 void (*UnpackOp)(PetscInt,PetscInt,const PetscInt*,void*,const void*); 1046 PetscErrorCode ierr; 1047 PetscSFBasicPack link; 1048 PetscInt i,nrootranks; 1049 PetscMPIInt typesize = -1; 1050 const PetscInt *rootoffset,*rootloc; 1051 1052 PetscFunctionBegin; 1053 ierr = PetscSFBasicGetPackInUse(sf,unit,leafdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr); 1054 /* This implementation could be changed to unpack as receives arrive, at the cost of non-determinism */ 1055 ierr = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr); 1056 ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,NULL,NULL,&rootoffset,&rootloc);CHKERRQ(ierr); 1057 ierr = PetscSFBasicPackGetUnpackOp(sf,link,op,&UnpackOp);CHKERRQ(ierr); 1058 if (UnpackOp) { 1059 typesize = link->unitbytes; 1060 } 1061 else { 1062 ierr = MPI_Type_size(unit,&typesize);CHKERRQ(ierr); 1063 } 1064 for (i=0; i<nrootranks; i++) { 1065 PetscMPIInt n = rootoffset[i+1] - rootoffset[i]; 1066 char *packstart = (char *) link->root[i]; 1067 1068 if (UnpackOp) { 1069 (*UnpackOp)(n,link->bs,rootloc+rootoffset[i],rootdata,(const void *)packstart); 1070 } 1071 #if defined(PETSC_HAVE_MPI_REDUCE_LOCAL) 1072 else if (n) { /* the op should be defined to operate on the whole datatype, so we ignore link->bs */ 1073 PetscInt j; 1074 1075 for (j = 0; j < n; j++) { 1076 ierr = MPI_Reduce_local(packstart+j*typesize,((char *) rootdata)+(rootloc[rootoffset[i]+j])*typesize,1,unit,op);CHKERRQ(ierr); 1077 } 1078 } 1079 #else 1080 else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No unpacking reduction operation for this MPI_Op"); 1081 #endif 1082 } 1083 ierr = PetscSFBasicReclaimPack(sf,&link);CHKERRQ(ierr); 1084 PetscFunctionReturn(0); 1085 } 1086 1087 static PetscErrorCode PetscSFFetchAndOpBegin_Basic(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op) 1088 { 1089 PetscErrorCode ierr; 1090 1091 PetscFunctionBegin; 1092 ierr = PetscSFReduceBegin_Basic(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr); 1093 PetscFunctionReturn(0); 1094 } 1095 1096 static PetscErrorCode PetscSFFetchAndOpEnd_Basic(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op) 1097 { 1098 PetscSF_Basic *bas = (PetscSF_Basic*)sf->data; 1099 void (*FetchAndOp)(PetscInt,PetscInt,const PetscInt*,void*,void*); 1100 PetscErrorCode ierr; 1101 PetscSFBasicPack link; 1102 PetscInt i,nrootranks,ndrootranks,nleafranks,ndleafranks; 1103 const PetscInt *rootoffset,*leafoffset,*rootloc,*leafloc; 1104 const PetscMPIInt *rootranks,*leafranks; 1105 MPI_Request *rootreqs,*leafreqs; 1106 1107 PetscFunctionBegin; 1108 ierr = PetscSFBasicGetPackInUse(sf,unit,leafdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr); 1109 /* This implementation could be changed to unpack as receives arrive, at the cost of non-determinism */ 1110 ierr = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr); 1111 ierr = PetscSFBasicGetRootInfo(sf,&nrootranks,&ndrootranks,&rootranks,&rootoffset,&rootloc);CHKERRQ(ierr); 1112 ierr = PetscSFBasicGetLeafInfo(sf,&nleafranks,&ndleafranks,&leafranks,&leafoffset,&leafloc);CHKERRQ(ierr); 1113 ierr = PetscSFBasicPackGetReqs(sf,link,&rootreqs,&leafreqs);CHKERRQ(ierr); 1114 /* Post leaf receives */ 1115 for (i=ndleafranks; i<nleafranks; i++) { 1116 PetscMPIInt n = leafoffset[i+1] - leafoffset[i]; 1117 ierr = MPI_Irecv(link->leaf[i],n,unit,leafranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&leafreqs[i-ndleafranks]);CHKERRQ(ierr); 1118 } 1119 /* Process local fetch-and-op, post root sends */ 1120 ierr = PetscSFBasicPackGetFetchAndOp(sf,link,op,&FetchAndOp);CHKERRQ(ierr); 1121 for (i=0; i<nrootranks; i++) { 1122 PetscMPIInt n = rootoffset[i+1] - rootoffset[i]; 1123 void *packstart = link->root[i]; 1124 1125 (*FetchAndOp)(n,link->bs,rootloc+rootoffset[i],rootdata,packstart); 1126 if (i < ndrootranks) continue; /* shared memory */ 1127 ierr = MPI_Isend(packstart,n,unit,rootranks[i],bas->tag,PetscObjectComm((PetscObject)sf),&rootreqs[i-ndrootranks]);CHKERRQ(ierr); 1128 } 1129 ierr = PetscSFBasicPackWaitall(sf,link);CHKERRQ(ierr); 1130 for (i=0; i<nleafranks; i++) { 1131 PetscMPIInt n = leafoffset[i+1] - leafoffset[i]; 1132 const void *packstart = link->leaf[i]; 1133 (*link->UnpackInsert)(n,link->bs,leafloc+leafoffset[i],leafupdate,packstart); 1134 } 1135 ierr = PetscSFBasicReclaimPack(sf,&link);CHKERRQ(ierr); 1136 PetscFunctionReturn(0); 1137 } 1138 1139 static PetscErrorCode PetscSFGetLeafRanks_Basic(PetscSF sf,PetscInt *niranks,const PetscMPIInt **iranks,const PetscInt **ioffset,const PetscInt **irootloc) 1140 { 1141 PetscSF_Basic *bas = (PetscSF_Basic*)sf->data; 1142 1143 PetscFunctionBegin; 1144 if (niranks) *niranks = bas->niranks; 1145 if (iranks) *iranks = bas->iranks; 1146 if (ioffset) *ioffset = bas->ioffset; 1147 if (irootloc) *irootloc = bas->irootloc; 1148 PetscFunctionReturn(0); 1149 } 1150 1151 PETSC_EXTERN PetscErrorCode PetscSFCreate_Basic(PetscSF sf) 1152 { 1153 PetscSF_Basic *bas = (PetscSF_Basic*)sf->data; 1154 PetscErrorCode ierr; 1155 1156 PetscFunctionBegin; 1157 sf->ops->SetUp = PetscSFSetUp_Basic; 1158 sf->ops->SetFromOptions = PetscSFSetFromOptions_Basic; 1159 sf->ops->Reset = PetscSFReset_Basic; 1160 sf->ops->Destroy = PetscSFDestroy_Basic; 1161 sf->ops->View = PetscSFView_Basic; 1162 sf->ops->BcastBegin = PetscSFBcastBegin_Basic; 1163 sf->ops->BcastEnd = PetscSFBcastEnd_Basic; 1164 sf->ops->BcastAndOpBegin = PetscSFBcastAndOpBegin_Basic; 1165 sf->ops->BcastAndOpEnd = PetscSFBcastAndOpEnd_Basic; 1166 sf->ops->ReduceBegin = PetscSFReduceBegin_Basic; 1167 sf->ops->ReduceEnd = PetscSFReduceEnd_Basic; 1168 sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Basic; 1169 sf->ops->FetchAndOpEnd = PetscSFFetchAndOpEnd_Basic; 1170 sf->ops->GetLeafRanks = PetscSFGetLeafRanks_Basic; 1171 1172 ierr = PetscNewLog(sf,&bas);CHKERRQ(ierr); 1173 sf->data = (void*)bas; 1174 PetscFunctionReturn(0); 1175 } 1176