Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Merge branch 'mc++' into mc-merge
[simgrid.git] / src / smpi / colls / bcast-flattree-pipeline.c
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   int pipe_length = count / segment;
24   int increment = segment * extent;
25   if (pipe_length==0) {
26     XBT_WARN("MPI_bcast_flattree_pipeline use default MPI_bcast_flattree.");
27     return smpi_coll_tuned_bcast_flattree(buff, count, data_type, root, comm);
28   }
29   rank = smpi_comm_rank(comm);
30   num_procs = smpi_comm_size(comm);
31
32   MPI_Request *request_array;
33   MPI_Status *status_array;
34
35   request_array = (MPI_Request *) xbt_malloc(pipe_length * sizeof(MPI_Request));
36   status_array = (MPI_Status *) xbt_malloc(pipe_length * sizeof(MPI_Status));
37
38   if (rank != root) {
39     for (i = 0; i < pipe_length; i++) {
40       request_array[i] = smpi_mpi_irecv((char *)buff + (i * increment), segment, data_type, root, tag, comm);
41     }
42     smpi_mpi_waitall(pipe_length, request_array, status_array);
43   }
44
45   else {
46     // Root sends data to all others
47     for (j = 0; j < num_procs; j++) {
48       if (j == rank)
49         continue;
50       else {
51         for (i = 0; i < pipe_length; i++) {
52           smpi_mpi_send((char *)buff + (i * increment), segment, data_type, j, tag, comm);
53         }
54       }
55     }
56
57   }
58
59   free(request_array);
60   free(status_array);
61   return MPI_SUCCESS;
62 }