Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Move all smpi colls to cpp.
[simgrid.git] / src / smpi / colls / bcast-flattree-pipeline.cpp
1 /* Copyright (c) 2013-2014. The SimGrid Team.
2  * All rights reserved.                                                     */
3
4 /* This program is free software; you can redistribute it and/or modify it
5  * under the terms of the license (GNU LGPL) which comes with this package. */
6
7 #include "colls_private.h"
8
9 int flattree_segment_in_byte = 8192;
10
11 int
12 smpi_coll_tuned_bcast_flattree_pipeline(void *buff, int count,
13                                         MPI_Datatype data_type, int root,
14                                         MPI_Comm comm)
15 {
16   int i, j, rank, num_procs;
17   int tag = COLL_TAG_BCAST;
18
19   MPI_Aint extent;
20   extent = smpi_datatype_get_extent(data_type);
21
22   int segment = flattree_segment_in_byte / extent;
23   segment =  segment == 0 ? 1 :segment; 
24   int pipe_length = count / segment;
25   int increment = segment * extent;
26   if (pipe_length==0) {
27     XBT_WARN("MPI_bcast_flattree_pipeline use default MPI_bcast_flattree.");
28     return smpi_coll_tuned_bcast_flattree(buff, count, data_type, root, comm);
29   }
30   rank = smpi_comm_rank(comm);
31   num_procs = smpi_comm_size(comm);
32
33   MPI_Request *request_array;
34   MPI_Status *status_array;
35
36   request_array = (MPI_Request *) xbt_malloc(pipe_length * sizeof(MPI_Request));
37   status_array = (MPI_Status *) xbt_malloc(pipe_length * sizeof(MPI_Status));
38
39   if (rank != root) {
40     for (i = 0; i < pipe_length; i++) {
41       request_array[i] = smpi_mpi_irecv((char *)buff + (i * increment), segment, data_type, root, tag, comm);
42     }
43     smpi_mpi_waitall(pipe_length, request_array, status_array);
44   }
45
46   else {
47     // Root sends data to all others
48     for (j = 0; j < num_procs; j++) {
49       if (j == rank)
50         continue;
51       else {
52         for (i = 0; i < pipe_length; i++) {
53           smpi_mpi_send((char *)buff + (i * increment), segment, data_type, j, tag, comm);
54         }
55       }
56     }
57
58   }
59
60   free(request_array);
61   free(status_array);
62   return MPI_SUCCESS;
63 }