1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773
|
/*
* Copyright (c) 2004-2005 The Trustees of Indiana University and Indiana
* University Research and Technology
* Corporation. All rights reserved.
* Copyright (c) 2004-2020 The University of Tennessee and The University
* of Tennessee Research Foundation. All rights
* reserved.
* Copyright (c) 2004-2005 High Performance Computing Center Stuttgart,
* University of Stuttgart. All rights reserved.
* Copyright (c) 2004-2005 The Regents of the University of California.
* All rights reserved.
* Copyright (c) 2008 Sun Microsystems, Inc. All rights reserved.
* Copyright (c) 2015-2018 Research Organization for Information Science
* and Technology (RIST). All rights reserved.
* Copyright (c) 2020 Amazon.com, Inc. or its affiliates.
* All Rights reserved.
* $COPYRIGHT$
*
* Additional copyrights may follow
*
* $HEADER$
*/
#include "ompi_config.h"
#include "mpi.h"
#include "ompi/constants.h"
#include "ompi/datatype/ompi_datatype.h"
#include "ompi/communicator/communicator.h"
#include "ompi/mca/coll/base/base.h"
#include "ompi/mca/coll/coll.h"
#include "ompi/mca/coll/base/coll_tags.h"
#include "coll_tuned.h"
/*
* Notes on evaluation rules and ordering
*
* The order is:
* use file based rules if presented (-coll_tuned_dynamic_rules_filename = rules)
* Else
* use forced rules (-coll_tuned_dynamic_ALG_intra_algorithm = algorithm-number)
* Else
* use fixed (compiled) rule set (or nested ifs)
*
*/
/*
* allreduce_intra
*
* Function: - allreduce using other MPI collectives
* Accepts: - same as MPI_Allreduce()
* Returns: - MPI_SUCCESS or error code
*/
int
ompi_coll_tuned_allreduce_intra_dec_dynamic (const void *sbuf, void *rbuf, int count,
struct ompi_datatype_t *dtype,
struct ompi_op_t *op,
struct ompi_communicator_t *comm,
mca_coll_base_module_t *module)
{
mca_coll_tuned_module_t *tuned_module = (mca_coll_tuned_module_t*) module;
OPAL_OUTPUT((ompi_coll_tuned_stream, "ompi_coll_tuned_allreduce_intra_dec_dynamic"));
/* Check first if an algorithm is set explicitly for this collective */
if (tuned_module->user_forced[ALLREDUCE].algorithm) {
return ompi_coll_tuned_allreduce_intra_do_this(sbuf, rbuf, count, dtype, op, comm, module,
tuned_module->user_forced[ALLREDUCE].algorithm,
tuned_module->user_forced[ALLREDUCE].tree_fanout,
tuned_module->user_forced[ALLREDUCE].segsize);
}
/* check to see if we have some filebased rules */
if (tuned_module->com_rules[ALLREDUCE]) {
/* we do, so calc the message size or what ever we need and use this for the evaluation */
int alg, faninout, segsize, ignoreme;
size_t dsize;
ompi_datatype_type_size (dtype, &dsize);
dsize *= count;
alg = ompi_coll_tuned_get_target_method_params (tuned_module->com_rules[ALLREDUCE],
dsize, &faninout, &segsize, &ignoreme);
if (alg) {
/* we have found a valid choice from the file based rules for this message size */
return ompi_coll_tuned_allreduce_intra_do_this (sbuf, rbuf, count, dtype, op,
comm, module,
alg, faninout, segsize);
} /* found a method */
} /*end if any com rules to check */
return ompi_coll_tuned_allreduce_intra_dec_fixed (sbuf, rbuf, count, dtype, op,
comm, module);
}
/*
* alltoall_intra_dec
*
* Function: - seletects alltoall algorithm to use
* Accepts: - same arguments as MPI_Alltoall()
* Returns: - MPI_SUCCESS or error code (passed from the alltoall implementation)
*/
int ompi_coll_tuned_alltoall_intra_dec_dynamic(const void *sbuf, int scount,
struct ompi_datatype_t *sdtype,
void* rbuf, int rcount,
struct ompi_datatype_t *rdtype,
struct ompi_communicator_t *comm,
mca_coll_base_module_t *module)
{
mca_coll_tuned_module_t *tuned_module = (mca_coll_tuned_module_t*) module;
OPAL_OUTPUT((ompi_coll_tuned_stream, "ompi_coll_tuned_alltoall_intra_dec_dynamic"));
/* Check first if an algorithm is set explicitly for this collective */
if (tuned_module->user_forced[ALLTOALL].algorithm) {
return ompi_coll_tuned_alltoall_intra_do_this(sbuf, scount, sdtype,
rbuf, rcount, rdtype,
comm, module,
tuned_module->user_forced[ALLTOALL].algorithm,
tuned_module->user_forced[ALLTOALL].tree_fanout,
tuned_module->user_forced[ALLTOALL].segsize,
tuned_module->user_forced[ALLTOALL].max_requests);
}
/* check to see if we have some filebased rules */
if (tuned_module->com_rules[ALLTOALL]) {
/* we do, so calc the message size or what ever we need and use this for the evaluation */
int comsize;
int alg, faninout, segsize, max_requests;
size_t dsize;
ompi_datatype_type_size (sdtype, &dsize);
comsize = ompi_comm_size(comm);
dsize *= (ptrdiff_t)comsize * (ptrdiff_t)scount;
alg = ompi_coll_tuned_get_target_method_params (tuned_module->com_rules[ALLTOALL],
dsize, &faninout, &segsize, &max_requests);
if (alg) {
/* we have found a valid choice from the file based rules for this message size */
return ompi_coll_tuned_alltoall_intra_do_this (sbuf, scount, sdtype,
rbuf, rcount, rdtype,
comm, module,
alg, faninout, segsize, max_requests);
} /* found a method */
} /*end if any com rules to check */
return ompi_coll_tuned_alltoall_intra_dec_fixed (sbuf, scount, sdtype,
rbuf, rcount, rdtype,
comm, module);
}
/*
* Function: - selects alltoallv algorithm to use
* Accepts: - same arguments as MPI_Alltoallv()
* Returns: - MPI_SUCCESS or error code
*/
int ompi_coll_tuned_alltoallv_intra_dec_dynamic(const void *sbuf, const int *scounts, const int *sdisps,
struct ompi_datatype_t *sdtype,
void* rbuf, const int *rcounts, const int *rdisps,
struct ompi_datatype_t *rdtype,
struct ompi_communicator_t *comm,
mca_coll_base_module_t *module)
{
mca_coll_tuned_module_t *tuned_module = (mca_coll_tuned_module_t*) module;
OPAL_OUTPUT((ompi_coll_tuned_stream, "ompi_coll_tuned_alltoallv_intra_dec_dynamic"));
/* Check first if an algorithm is set explicitly for this collective */
if (tuned_module->user_forced[ALLTOALLV].algorithm) {
return ompi_coll_tuned_alltoallv_intra_do_this(sbuf, scounts, sdisps, sdtype,
rbuf, rcounts, rdisps, rdtype,
comm, module,
tuned_module->user_forced[ALLTOALLV].algorithm);
}
/**
* check to see if we have some filebased rules. As we don't have global
* knowledge about the total amount of data, use the first available rule.
* This allow the users to specify the alltoallv algorithm to be used only
* based on the communicator size.
*/
if (tuned_module->com_rules[ALLTOALLV]) {
int alg, faninout, segsize, max_requests;
alg = ompi_coll_tuned_get_target_method_params (tuned_module->com_rules[ALLTOALLV],
0, &faninout, &segsize, &max_requests);
if (alg) {
/* we have found a valid choice from the file based rules for this message size */
return ompi_coll_tuned_alltoallv_intra_do_this (sbuf, scounts, sdisps, sdtype,
rbuf, rcounts, rdisps, rdtype,
comm, module,
alg);
} /* found a method */
} /*end if any com rules to check */
return ompi_coll_tuned_alltoallv_intra_dec_fixed(sbuf, scounts, sdisps, sdtype,
rbuf, rcounts, rdisps, rdtype,
comm, module);
}
/*
* barrier_intra_dec
*
* Function: - seletects barrier algorithm to use
* Accepts: - same arguments as MPI_Barrier()
* Returns: - MPI_SUCCESS or error code (passed from the barrier implementation)
*/
int ompi_coll_tuned_barrier_intra_dec_dynamic(struct ompi_communicator_t *comm,
mca_coll_base_module_t *module)
{
mca_coll_tuned_module_t *tuned_module = (mca_coll_tuned_module_t*) module;
OPAL_OUTPUT((ompi_coll_tuned_stream,"ompi_coll_tuned_barrier_intra_dec_dynamic"));
/* Check first if an algorithm is set explicitly for this collective */
if (tuned_module->user_forced[BARRIER].algorithm) {
return ompi_coll_tuned_barrier_intra_do_this(comm, module,
tuned_module->user_forced[BARRIER].algorithm,
tuned_module->user_forced[BARRIER].tree_fanout,
tuned_module->user_forced[BARRIER].segsize);
}
/* check to see if we have some filebased rules */
if (tuned_module->com_rules[BARRIER]) {
/* we do, so calc the message size or what ever we need and use this for the evaluation */
int alg, faninout, segsize, ignoreme;
alg = ompi_coll_tuned_get_target_method_params (tuned_module->com_rules[BARRIER],
0, &faninout, &segsize, &ignoreme);
if (alg) {
/* we have found a valid choice from the file based rules for this message size */
return ompi_coll_tuned_barrier_intra_do_this (comm, module,
alg, faninout, segsize);
} /* found a method */
} /*end if any com rules to check */
return ompi_coll_tuned_barrier_intra_dec_fixed (comm, module);
}
/*
* bcast_intra_dec
*
* Function: - selects broadcast algorithm to use
* Accepts: - same arguments as MPI_Bcast()
* Returns: - MPI_SUCCESS or error code (passed from the bcast implementation)
*/
int ompi_coll_tuned_bcast_intra_dec_dynamic(void *buf, int count,
struct ompi_datatype_t *dtype, int root,
struct ompi_communicator_t *comm,
mca_coll_base_module_t *module)
{
mca_coll_tuned_module_t *tuned_module = (mca_coll_tuned_module_t*) module;
OPAL_OUTPUT((ompi_coll_tuned_stream, "coll:tuned:bcast_intra_dec_dynamic"));
/* Check first if an algorithm is set explicitly for this collective */
if (tuned_module->user_forced[BCAST].algorithm) {
return ompi_coll_tuned_bcast_intra_do_this(buf, count, dtype,
root, comm, module,
tuned_module->user_forced[BCAST].algorithm,
tuned_module->user_forced[BCAST].chain_fanout,
tuned_module->user_forced[BCAST].segsize);
}
/* check to see if we have some filebased rules */
if (tuned_module->com_rules[BCAST]) {
/* we do, so calc the message size or what ever we need and use this for the evaluation */
int alg, faninout, segsize, ignoreme;
size_t dsize;
ompi_datatype_type_size (dtype, &dsize);
dsize *= count;
alg = ompi_coll_tuned_get_target_method_params (tuned_module->com_rules[BCAST],
dsize, &faninout, &segsize, &ignoreme);
if (alg) {
/* we have found a valid choice from the file based rules for this message size */
return ompi_coll_tuned_bcast_intra_do_this (buf, count, dtype, root,
comm, module,
alg, faninout, segsize);
} /* found a method */
} /*end if any com rules to check */
return ompi_coll_tuned_bcast_intra_dec_fixed (buf, count, dtype, root,
comm, module);
}
/*
* reduce_intra_dec
*
* Function: - seletects reduce algorithm to use
* Accepts: - same arguments as MPI_reduce()
* Returns: - MPI_SUCCESS or error code (passed from the reduce implementation)
*
*/
int ompi_coll_tuned_reduce_intra_dec_dynamic( const void *sbuf, void *rbuf,
int count, struct ompi_datatype_t* dtype,
struct ompi_op_t* op, int root,
struct ompi_communicator_t* comm,
mca_coll_base_module_t *module)
{
mca_coll_tuned_module_t *tuned_module = (mca_coll_tuned_module_t*) module;
OPAL_OUTPUT((ompi_coll_tuned_stream, "coll:tuned:reduce_intra_dec_dynamic"));
/* Check first if an algorithm is set explicitly for this collective */
if (tuned_module->user_forced[REDUCE].algorithm) {
return ompi_coll_tuned_reduce_intra_do_this(sbuf, rbuf, count, dtype,
op, root, comm, module,
tuned_module->user_forced[REDUCE].algorithm,
tuned_module->user_forced[REDUCE].chain_fanout,
tuned_module->user_forced[REDUCE].segsize,
tuned_module->user_forced[REDUCE].max_requests);
}
/* check to see if we have some filebased rules */
if (tuned_module->com_rules[REDUCE]) {
/* we do, so calc the message size or what ever we need and use this for the evaluation */
int alg, faninout, segsize, max_requests;
size_t dsize;
ompi_datatype_type_size(dtype, &dsize);
dsize *= count;
alg = ompi_coll_tuned_get_target_method_params (tuned_module->com_rules[REDUCE],
dsize, &faninout, &segsize, &max_requests);
if (alg) {
/* we have found a valid choice from the file based rules for this message size */
return ompi_coll_tuned_reduce_intra_do_this (sbuf, rbuf, count, dtype,
op, root, comm, module,
alg, faninout,
segsize, max_requests);
} /* found a method */
} /*end if any com rules to check */
return ompi_coll_tuned_reduce_intra_dec_fixed (sbuf, rbuf, count, dtype,
op, root, comm, module);
}
/*
* reduce_scatter_intra_dec
*
* Function: - seletects reduce_scatter algorithm to use
* Accepts: - same arguments as MPI_Reduce_scatter()
* Returns: - MPI_SUCCESS or error code (passed from
* the reduce_scatter implementation)
*
*/
int ompi_coll_tuned_reduce_scatter_intra_dec_dynamic(const void *sbuf, void *rbuf,
const int *rcounts,
struct ompi_datatype_t *dtype,
struct ompi_op_t *op,
struct ompi_communicator_t *comm,
mca_coll_base_module_t *module)
{
mca_coll_tuned_module_t *tuned_module = (mca_coll_tuned_module_t*) module;
OPAL_OUTPUT((ompi_coll_tuned_stream, "coll:tuned:reduce_scatter_intra_dec_dynamic"));
/* Check first if an algorithm is set explicitly for this collective */
if (tuned_module->user_forced[REDUCESCATTER].algorithm) {
return ompi_coll_tuned_reduce_scatter_intra_do_this(sbuf, rbuf, rcounts, dtype,
op, comm, module,
tuned_module->user_forced[REDUCESCATTER].algorithm,
tuned_module->user_forced[REDUCESCATTER].chain_fanout,
tuned_module->user_forced[REDUCESCATTER].segsize);
}
/* check to see if we have some filebased rules */
if (tuned_module->com_rules[REDUCESCATTER]) {
/* we do, so calc the message size or what ever we need and use
this for the evaluation */
int alg, faninout, segsize, ignoreme, i, count, size;
size_t dsize;
size = ompi_comm_size(comm);
for (i = 0, count = 0; i < size; i++) { count += rcounts[i];}
ompi_datatype_type_size (dtype, &dsize);
dsize *= count;
alg = ompi_coll_tuned_get_target_method_params (tuned_module->com_rules[REDUCESCATTER],
dsize, &faninout,
&segsize, &ignoreme);
if (alg) {
/* we have found a valid choice from the file based rules for this message size */
return ompi_coll_tuned_reduce_scatter_intra_do_this (sbuf, rbuf, rcounts, dtype,
op, comm, module,
alg, faninout, segsize);
} /* found a method */
} /*end if any com rules to check */
return ompi_coll_tuned_reduce_scatter_intra_dec_fixed (sbuf, rbuf, rcounts,
dtype, op, comm, module);
}
/*
* reduce_scatter_block_intra_dec
*
* Function: - seletects reduce_scatter_block algorithm to use
* Accepts: - same arguments as MPI_Reduce_scatter_block()
* Returns: - MPI_SUCCESS or error code (passed from
* the reduce_scatter implementation)
*
*/
int ompi_coll_tuned_reduce_scatter_block_intra_dec_dynamic(const void *sbuf, void *rbuf,
int rcount,
struct ompi_datatype_t *dtype,
struct ompi_op_t *op,
struct ompi_communicator_t *comm,
mca_coll_base_module_t *module)
{
mca_coll_tuned_module_t *tuned_module = (mca_coll_tuned_module_t*) module;
OPAL_OUTPUT((ompi_coll_tuned_stream, "coll:tuned:reduce_scatter_block_intra_dec_dynamic"));
/* Check first if an algorithm is set explicitly for this collective */
if (tuned_module->user_forced[REDUCESCATTERBLOCK].algorithm) {
return ompi_coll_tuned_reduce_scatter_block_intra_do_this(sbuf, rbuf, rcount, dtype,
op, comm, module,
tuned_module->user_forced[REDUCESCATTERBLOCK].algorithm,
tuned_module->user_forced[REDUCESCATTERBLOCK].chain_fanout,
tuned_module->user_forced[REDUCESCATTERBLOCK].segsize);
}
/* check to see if we have some filebased rules */
if (tuned_module->com_rules[REDUCESCATTERBLOCK]) {
/* we do, so calc the message size or what ever we need and use
this for the evaluation */
int alg, faninout, segsize, ignoreme, size;
size_t dsize;
size = ompi_comm_size(comm);
ompi_datatype_type_size (dtype, &dsize);
dsize *= rcount * size;
alg = ompi_coll_tuned_get_target_method_params(tuned_module->com_rules[REDUCESCATTERBLOCK],
dsize, &faninout,
&segsize, &ignoreme);
if (alg) {
/* we have found a valid choice from the file based rules for this message size */
return ompi_coll_tuned_reduce_scatter_block_intra_do_this (sbuf, rbuf, rcount, dtype,
op, comm, module,
alg, faninout, segsize);
} /* found a method */
} /* end if any com rules to check */
return ompi_coll_tuned_reduce_scatter_block_intra_dec_fixed (sbuf, rbuf, rcount,
dtype, op, comm, module);
}
/*
* allgather_intra_dec
*
* Function: - seletects allgather algorithm to use
* Accepts: - same arguments as MPI_Allgather()
* Returns: - MPI_SUCCESS or error code (passed from the selected
* allgather function).
*/
int ompi_coll_tuned_allgather_intra_dec_dynamic(const void *sbuf, int scount,
struct ompi_datatype_t *sdtype,
void* rbuf, int rcount,
struct ompi_datatype_t *rdtype,
struct ompi_communicator_t *comm,
mca_coll_base_module_t *module)
{
mca_coll_tuned_module_t *tuned_module = (mca_coll_tuned_module_t*) module;
OPAL_OUTPUT((ompi_coll_tuned_stream,
"ompi_coll_tuned_allgather_intra_dec_dynamic"));
/* Check first if an algorithm is set explicitly for this collective */
if (tuned_module->user_forced[ALLGATHER].algorithm) {
/* User-forced algorithm */
return ompi_coll_tuned_allgather_intra_do_this(sbuf, scount, sdtype,
rbuf, rcount, rdtype,
comm, module,
tuned_module->user_forced[ALLGATHER].algorithm,
tuned_module->user_forced[ALLGATHER].tree_fanout,
tuned_module->user_forced[ALLGATHER].segsize);
}
if (tuned_module->com_rules[ALLGATHER]) {
/* We have file based rules:
- calculate message size and other necessary information */
int comsize;
int alg, faninout, segsize, ignoreme;
size_t dsize;
ompi_datatype_type_size (sdtype, &dsize);
comsize = ompi_comm_size(comm);
dsize *= (ptrdiff_t)comsize * (ptrdiff_t)scount;
alg = ompi_coll_tuned_get_target_method_params (tuned_module->com_rules[ALLGATHER],
dsize, &faninout, &segsize, &ignoreme);
if (alg) {
/* we have found a valid choice from the file based rules for
this message size */
return ompi_coll_tuned_allgather_intra_do_this (sbuf, scount, sdtype,
rbuf, rcount, rdtype,
comm, module,
alg, faninout, segsize);
}
}
/* Use default decision */
return ompi_coll_tuned_allgather_intra_dec_fixed (sbuf, scount, sdtype,
rbuf, rcount, rdtype,
comm, module);
}
/*
* allgatherv_intra_dec
*
* Function: - seletects allgatherv algorithm to use
* Accepts: - same arguments as MPI_Allgatherv()
* Returns: - MPI_SUCCESS or error code (passed from the selected
* allgatherv function).
*/
int ompi_coll_tuned_allgatherv_intra_dec_dynamic(const void *sbuf, int scount,
struct ompi_datatype_t *sdtype,
void* rbuf, const int *rcounts,
const int *rdispls,
struct ompi_datatype_t *rdtype,
struct ompi_communicator_t *comm,
mca_coll_base_module_t *module)
{
mca_coll_tuned_module_t *tuned_module = (mca_coll_tuned_module_t*) module;
OPAL_OUTPUT((ompi_coll_tuned_stream,
"ompi_coll_tuned_allgatherv_intra_dec_dynamic"));
/* Check first if an algorithm is set explicitly for this collective */
if (tuned_module->user_forced[ALLGATHERV].algorithm) {
/* User-forced algorithm */
return ompi_coll_tuned_allgatherv_intra_do_this(sbuf, scount, sdtype,
rbuf, rcounts, rdispls, rdtype,
comm, module,
tuned_module->user_forced[ALLGATHERV].algorithm,
tuned_module->user_forced[ALLGATHERV].tree_fanout,
tuned_module->user_forced[ALLGATHERV].segsize);
}
if (tuned_module->com_rules[ALLGATHERV]) {
/* We have file based rules:
- calculate message size and other necessary information */
int comsize, i;
int alg, faninout, segsize, ignoreme;
size_t dsize, total_size, per_rank_size;
comsize = ompi_comm_size(comm);
ompi_datatype_type_size (sdtype, &dsize);
total_size = 0;
for (i = 0; i < comsize; i++) { total_size += dsize * rcounts[i]; }
per_rank_size = total_size / comsize;
alg = ompi_coll_tuned_get_target_method_params (tuned_module->com_rules[ALLGATHERV],
per_rank_size, &faninout, &segsize, &ignoreme);
if (alg) {
/* we have found a valid choice from the file based rules for
this message size */
return ompi_coll_tuned_allgatherv_intra_do_this (sbuf, scount, sdtype,
rbuf, rcounts,
rdispls, rdtype,
comm, module,
alg, faninout, segsize);
}
}
/* Use default decision */
return ompi_coll_tuned_allgatherv_intra_dec_fixed (sbuf, scount, sdtype,
rbuf, rcounts,
rdispls, rdtype,
comm, module);
}
int ompi_coll_tuned_gather_intra_dec_dynamic(const void *sbuf, int scount,
struct ompi_datatype_t *sdtype,
void* rbuf, int rcount,
struct ompi_datatype_t *rdtype,
int root,
struct ompi_communicator_t *comm,
mca_coll_base_module_t *module)
{
mca_coll_tuned_module_t *tuned_module = (mca_coll_tuned_module_t*) module;
OPAL_OUTPUT((ompi_coll_tuned_stream,
"ompi_coll_tuned_gather_intra_dec_dynamic"));
/* Check first if an algorithm is set explicitly for this collective */
if (tuned_module->user_forced[GATHER].algorithm) {
return ompi_coll_tuned_gather_intra_do_this(sbuf, scount, sdtype,
rbuf, rcount, rdtype,
root, comm, module,
tuned_module->user_forced[GATHER].algorithm,
tuned_module->user_forced[GATHER].tree_fanout,
tuned_module->user_forced[GATHER].segsize);
}
/**
* check to see if we have some filebased rules.
*/
if (tuned_module->com_rules[GATHER]) {
int comsize, alg, faninout, segsize, max_requests;
size_t dsize;
comsize = ompi_comm_size(comm);
ompi_datatype_type_size (sdtype, &dsize);
dsize *= scount * comsize;
alg = ompi_coll_tuned_get_target_method_params (tuned_module->com_rules[GATHER],
dsize, &faninout, &segsize, &max_requests);
if (alg) {
/* we have found a valid choice from the file based rules for this message size */
return ompi_coll_tuned_gather_intra_do_this (sbuf, scount, sdtype,
rbuf, rcount, rdtype,
root, comm, module,
alg, faninout, segsize);
} /* found a method */
} /*end if any com rules to check */
return ompi_coll_tuned_gather_intra_dec_fixed (sbuf, scount, sdtype,
rbuf, rcount, rdtype,
root, comm, module);
}
int ompi_coll_tuned_scatter_intra_dec_dynamic(const void *sbuf, int scount,
struct ompi_datatype_t *sdtype,
void* rbuf, int rcount,
struct ompi_datatype_t *rdtype,
int root, struct ompi_communicator_t *comm,
mca_coll_base_module_t *module)
{
mca_coll_tuned_module_t *tuned_module = (mca_coll_tuned_module_t*) module;
OPAL_OUTPUT((ompi_coll_tuned_stream,
"ompi_coll_tuned_scatter_intra_dec_dynamic"));
/* Check first if an algorithm is set explicitly for this collective */
if (tuned_module->user_forced[SCATTER].algorithm) {
return ompi_coll_tuned_scatter_intra_do_this(sbuf, scount, sdtype,
rbuf, rcount, rdtype,
root, comm, module,
tuned_module->user_forced[SCATTER].algorithm,
tuned_module->user_forced[SCATTER].chain_fanout,
tuned_module->user_forced[SCATTER].segsize);
}
/**
* check to see if we have some filebased rules.
*/
if (tuned_module->com_rules[SCATTER]) {
int comsize, alg, faninout, segsize, max_requests;
size_t dsize;
comsize = ompi_comm_size(comm);
ompi_datatype_type_size (sdtype, &dsize);
dsize *= scount * comsize;
alg = ompi_coll_tuned_get_target_method_params (tuned_module->com_rules[SCATTER],
dsize, &faninout, &segsize, &max_requests);
if (alg) {
/* we have found a valid choice from the file based rules for this message size */
return ompi_coll_tuned_scatter_intra_do_this (sbuf, scount, sdtype,
rbuf, rcount, rdtype,
root, comm, module,
alg, faninout, segsize);
} /* found a method */
} /*end if any com rules to check */
return ompi_coll_tuned_scatter_intra_dec_fixed (sbuf, scount, sdtype,
rbuf, rcount, rdtype,
root, comm, module);
}
int ompi_coll_tuned_exscan_intra_dec_dynamic(const void *sbuf, void* rbuf, int count,
struct ompi_datatype_t *dtype,
struct ompi_op_t *op,
struct ompi_communicator_t *comm,
mca_coll_base_module_t *module)
{
mca_coll_tuned_module_t *tuned_module = (mca_coll_tuned_module_t*) module;
OPAL_OUTPUT((ompi_coll_tuned_stream,
"ompi_coll_tuned_exscan_intra_dec_dynamic"));
/* Check first if an algorithm is set explicitly for this collective */
if (tuned_module->user_forced[EXSCAN].algorithm) {
return ompi_coll_tuned_exscan_intra_do_this(sbuf, rbuf, count, dtype,
op, comm, module,
tuned_module->user_forced[EXSCAN].algorithm);
}
/**
* check to see if we have some filebased rules.
*/
if (tuned_module->com_rules[EXSCAN]) {
int comsize, alg, faninout, segsize, max_requests;
size_t dsize;
comsize = ompi_comm_size(comm);
ompi_datatype_type_size (dtype, &dsize);
dsize *= comsize;
alg = ompi_coll_tuned_get_target_method_params (tuned_module->com_rules[EXSCAN],
dsize, &faninout, &segsize, &max_requests);
if (alg) {
/* we have found a valid choice from the file based rules for this message size */
return ompi_coll_tuned_exscan_intra_do_this (sbuf, rbuf, count, dtype,
op, comm, module,
alg);
} /* found a method */
} /*end if any com rules to check */
return ompi_coll_base_exscan_intra_linear(sbuf, rbuf, count, dtype,
op, comm, module);
}
int ompi_coll_tuned_scan_intra_dec_dynamic(const void *sbuf, void* rbuf, int count,
struct ompi_datatype_t *dtype,
struct ompi_op_t *op,
struct ompi_communicator_t *comm,
mca_coll_base_module_t *module)
{
mca_coll_tuned_module_t *tuned_module = (mca_coll_tuned_module_t*) module;
OPAL_OUTPUT((ompi_coll_tuned_stream,
"ompi_coll_tuned_scan_intra_dec_dynamic"));
/* Check first if an algorithm is set explicitly for this collective */
if (tuned_module->user_forced[SCAN].algorithm) {
return ompi_coll_tuned_scan_intra_do_this(sbuf, rbuf, count, dtype,
op, comm, module,
tuned_module->user_forced[SCAN].algorithm);
}
/**
* check to see if we have some filebased rules.
*/
if (tuned_module->com_rules[SCAN]) {
int comsize, alg, faninout, segsize, max_requests;
size_t dsize;
comsize = ompi_comm_size(comm);
ompi_datatype_type_size (dtype, &dsize);
dsize *= comsize;
alg = ompi_coll_tuned_get_target_method_params (tuned_module->com_rules[SCAN],
dsize, &faninout, &segsize, &max_requests);
if (alg) {
/* we have found a valid choice from the file based rules for this message size */
return ompi_coll_tuned_scan_intra_do_this (sbuf, rbuf, count, dtype,
op, comm, module,
alg);
} /* found a method */
} /*end if any com rules to check */
return ompi_coll_base_scan_intra_linear(sbuf, rbuf, count, dtype,
op, comm, module);
}
|