Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Remove unused variables.
[simgrid.git] / src / smpi / colls / bcast-SMP-linear.c
1 #include "colls_private.h"
2 #ifndef NUM_CORE
3 #define NUM_CORE 8
4 #endif
5
6 int bcast_SMP_linear_segment_byte = 8192;
7
8 int smpi_coll_tuned_bcast_SMP_linear(void *buf, int count,
9                                      MPI_Datatype datatype, int root,
10                                      MPI_Comm comm)
11 {
12   int tag = 5000;
13   MPI_Status status;
14   MPI_Request request;
15   MPI_Request *request_array;
16   MPI_Status *status_array;
17   int rank, size;
18   int i;
19   MPI_Aint extent;
20   extent = smpi_datatype_get_extent(datatype);
21
22   rank = smpi_comm_rank(comm);
23   size = smpi_comm_size(comm);
24
25   int segment = bcast_SMP_linear_segment_byte / extent;
26   int pipe_length = count / segment;
27   int remainder = count % segment;
28   int increment = segment * extent;
29
30
31   /* leader of each SMP do inter-communication
32      and act as a root for intra-communication */
33   int to_inter = (rank + NUM_CORE) % size;
34   int to_intra = (rank + 1) % size;
35   int from_inter = (rank - NUM_CORE + size) % size;
36   int from_intra = (rank + size - 1) % size;
37
38   // call native when MPI communication size is too small
39   if (size <= NUM_CORE) {
40     XBT_WARN("MPI_bcast_SMP_linear use default MPI_bcast.");              
41     smpi_mpi_bcast(buf, count, datatype, root, comm);
42     return MPI_SUCCESS;            
43   }
44   // if root is not zero send to rank zero first
45   if (root != 0) {
46     if (rank == root)
47       smpi_mpi_send(buf, count, datatype, 0, tag, comm);
48     else if (rank == 0)
49       smpi_mpi_recv(buf, count, datatype, root, tag, comm, &status);
50   }
51   // when a message is smaller than a block size => no pipeline 
52   if (count <= segment) {
53     // case ROOT
54     if (rank == 0) {
55       smpi_mpi_send(buf, count, datatype, to_inter, tag, comm);
56       smpi_mpi_send(buf, count, datatype, to_intra, tag, comm);
57     }
58     // case last ROOT of each SMP
59     else if (rank == (((size - 1) / NUM_CORE) * NUM_CORE)) {
60       request = smpi_mpi_irecv(buf, count, datatype, from_inter, tag, comm);
61       smpi_mpi_wait(&request, &status);
62       smpi_mpi_send(buf, count, datatype, to_intra, tag, comm);
63     }
64     // case intermediate ROOT of each SMP
65     else if (rank % NUM_CORE == 0) {
66       request = smpi_mpi_irecv(buf, count, datatype, from_inter, tag, comm);
67       smpi_mpi_wait(&request, &status);
68       smpi_mpi_send(buf, count, datatype, to_inter, tag, comm);
69       smpi_mpi_send(buf, count, datatype, to_intra, tag, comm);
70     }
71     // case last non-ROOT of each SMP
72     else if (((rank + 1) % NUM_CORE == 0) || (rank == (size - 1))) {
73       request = smpi_mpi_irecv(buf, count, datatype, from_intra, tag, comm);
74       smpi_mpi_wait(&request, &status);
75     }
76     // case intermediate non-ROOT of each SMP
77     else {
78       request = smpi_mpi_irecv(buf, count, datatype, from_intra, tag, comm);
79       smpi_mpi_wait(&request, &status);
80       smpi_mpi_send(buf, count, datatype, to_intra, tag, comm);
81     }
82     return MPI_SUCCESS;
83   }
84   // pipeline bcast
85   else {
86     request_array =
87         (MPI_Request *) xbt_malloc((size + pipe_length) * sizeof(MPI_Request));
88     status_array =
89         (MPI_Status *) xbt_malloc((size + pipe_length) * sizeof(MPI_Status));
90
91     // case ROOT of each SMP
92     if (rank % NUM_CORE == 0) {
93       // case real root
94       if (rank == 0) {
95         for (i = 0; i < pipe_length; i++) {
96           smpi_mpi_send((char *) buf + (i * increment), segment, datatype, to_inter,
97                    (tag + i), comm);
98           smpi_mpi_send((char *) buf + (i * increment), segment, datatype, to_intra,
99                    (tag + i), comm);
100         }
101       }
102       // case last ROOT of each SMP
103       else if (rank == (((size - 1) / NUM_CORE) * NUM_CORE)) {
104         for (i = 0; i < pipe_length; i++) {
105           request_array[i] = smpi_mpi_irecv((char *) buf + (i * increment), segment, datatype,
106                     from_inter, (tag + i), comm);
107         }
108         for (i = 0; i < pipe_length; i++) {
109           smpi_mpi_wait(&request_array[i], &status);
110           smpi_mpi_send((char *) buf + (i * increment), segment, datatype, to_intra,
111                    (tag + i), comm);
112         }
113       }
114       // case intermediate ROOT of each SMP
115       else {
116         for (i = 0; i < pipe_length; i++) {
117           request_array[i] = smpi_mpi_irecv((char *) buf + (i * increment), segment, datatype,
118                     from_inter, (tag + i), comm);
119         }
120         for (i = 0; i < pipe_length; i++) {
121           smpi_mpi_wait(&request_array[i], &status);
122           smpi_mpi_send((char *) buf + (i * increment), segment, datatype, to_inter,
123                    (tag + i), comm);
124           smpi_mpi_send((char *) buf + (i * increment), segment, datatype, to_intra,
125                    (tag + i), comm);
126         }
127       }
128     } else {                    // case last non-ROOT of each SMP
129       if (((rank + 1) % NUM_CORE == 0) || (rank == (size - 1))) {
130         for (i = 0; i < pipe_length; i++) {
131           request_array[i] = smpi_mpi_irecv((char *) buf + (i * increment), segment, datatype,
132                     from_intra, (tag + i), comm);
133         }
134         for (i = 0; i < pipe_length; i++) {
135           smpi_mpi_wait(&request_array[i], &status);
136         }
137       }
138       // case intermediate non-ROOT of each SMP
139       else {
140         for (i = 0; i < pipe_length; i++) {
141           request_array[i] = smpi_mpi_irecv((char *) buf + (i * increment), segment, datatype,
142                     from_intra, (tag + i), comm);
143         }
144         for (i = 0; i < pipe_length; i++) {
145           smpi_mpi_wait(&request_array[i], &status);
146           smpi_mpi_send((char *) buf + (i * increment), segment, datatype, to_intra,
147                    (tag + i), comm);
148         }
149       }
150     }
151     free(request_array);
152     free(status_array);
153   }
154
155   // when count is not divisible by block size, use default BCAST for the remainder
156   if ((remainder != 0) && (count > segment)) {
157     XBT_WARN("MPI_bcast_SMP_linear use default MPI_bcast.");                     
158     smpi_mpi_bcast((char *) buf + (pipe_length * increment), remainder, datatype,
159               root, comm);
160   }
161
162   return MPI_SUCCESS;
163 }