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