Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
unify collective tags
[simgrid.git] / src / smpi / colls / bcast-flattree-pipeline.c
1 #include "colls_private.h"
2
3 int flattree_segment_in_byte = 8192;
4
5 int
6 smpi_coll_tuned_bcast_flattree_pipeline(void *buff, int count,
7                                         MPI_Datatype data_type, int root,
8                                         MPI_Comm comm)
9 {
10   int i, j, rank, num_procs;
11   int tag = COLL_TAG_BCAST;
12
13   MPI_Aint extent;
14   extent = smpi_datatype_get_extent(data_type);
15
16   int segment = flattree_segment_in_byte / extent;
17   int pipe_length = count / segment;
18   int increment = segment * extent;
19   if (pipe_length==0) {
20     XBT_WARN("MPI_bcast_flattree_pipeline use default MPI_bcast_flattree.");
21     return smpi_coll_tuned_bcast_flattree(buff, count, data_type, root, comm);
22   }
23   rank = smpi_comm_rank(comm);
24   num_procs = smpi_comm_size(comm);
25
26   MPI_Request *request_array;
27   MPI_Status *status_array;
28
29   request_array = (MPI_Request *) xbt_malloc(pipe_length * sizeof(MPI_Request));
30   status_array = (MPI_Status *) xbt_malloc(pipe_length * sizeof(MPI_Status));
31
32   if (rank != root) {
33     for (i = 0; i < pipe_length; i++) {
34       request_array[i] = smpi_mpi_irecv((char *)buff + (i * increment), segment, data_type, root, tag, comm);
35     }
36     smpi_mpi_waitall(pipe_length, request_array, status_array);
37   }
38
39   else {
40     // Root sends data to all others
41     for (j = 0; j < num_procs; j++) {
42       if (j == rank)
43         continue;
44       else {
45         for (i = 0; i < pipe_length; i++) {
46           smpi_mpi_send((char *)buff + (i * increment), segment, data_type, j, tag, comm);
47         }
48       }
49     }
50
51   }
52
53   free(request_array);
54   free(status_array);
55   return MPI_SUCCESS;
56 }