Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
f01bbae76dae1f2d86c4ffb9fb051f31c282cd56
[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 = 1;
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
20   rank = smpi_comm_rank(comm);
21   num_procs = smpi_comm_size(comm);
22
23   MPI_Request *request_array;
24   MPI_Status *status_array;
25
26   request_array = (MPI_Request *) xbt_malloc(pipe_length * sizeof(MPI_Request));
27   status_array = (MPI_Status *) xbt_malloc(pipe_length * sizeof(MPI_Status));
28
29   if (rank != root) {
30     for (i = 0; i < pipe_length; i++) {
31       request_array[i] = smpi_mpi_irecv((char *)buff + (i * increment), segment, data_type, root, tag, comm);
32     }
33     smpi_mpi_waitall(pipe_length, request_array, status_array);
34   }
35
36   else {
37     // Root sends data to all others
38     for (j = 0; j < num_procs; j++) {
39       if (j == rank)
40         continue;
41       else {
42         for (i = 0; i < pipe_length; i++) {
43           smpi_mpi_send((char *)buff + (i * increment), segment, data_type, j, tag, comm);
44         }
45       }
46     }
47
48   }
49
50   free(request_array);
51   free(status_array);
52   return MPI_SUCCESS;
53 }