root/ompi/mca/coll/tuned/coll_tuned_decision_dynamic.c

/* [<][>][^][v][top][bottom][index][help] */

DEFINITIONS

This source file includes following definitions.
  1. ompi_coll_tuned_allreduce_intra_dec_dynamic
  2. ompi_coll_tuned_alltoall_intra_dec_dynamic
  3. ompi_coll_tuned_alltoallv_intra_dec_dynamic
  4. ompi_coll_tuned_barrier_intra_dec_dynamic
  5. ompi_coll_tuned_bcast_intra_dec_dynamic
  6. ompi_coll_tuned_reduce_intra_dec_dynamic
  7. ompi_coll_tuned_reduce_scatter_intra_dec_dynamic
  8. ompi_coll_tuned_reduce_scatter_block_intra_dec_dynamic
  9. ompi_coll_tuned_allgather_intra_dec_dynamic
  10. ompi_coll_tuned_allgatherv_intra_dec_dynamic
  11. ompi_coll_tuned_gather_intra_dec_dynamic
  12. ompi_coll_tuned_scatter_intra_dec_dynamic
  13. ompi_coll_tuned_exscan_intra_dec_dynamic
  14. ompi_coll_tuned_scan_intra_dec_dynamic

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

/* [<][>][^][v][top][bottom][index][help] */