Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Update copyright notices
[simgrid.git] / src / smpi / colls / bcast-NTSB.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 bcast_NTSB_segment_size_in_byte = 8192;
10
11 int smpi_coll_tuned_bcast_NTSB(void *buf, int count, MPI_Datatype datatype,
12                                int root, MPI_Comm comm)
13 {
14   int tag = COLL_TAG_BCAST;
15   MPI_Status status;
16   int rank, size;
17   int i;
18
19   MPI_Request *send_request_array;
20   MPI_Request *recv_request_array;
21   MPI_Status *send_status_array;
22   MPI_Status *recv_status_array;
23
24   MPI_Aint extent;
25   extent = smpi_datatype_get_extent(datatype);
26
27   rank = smpi_comm_rank(comm);
28   size = smpi_comm_size(comm);
29
30   /* source node and destination nodes (same through out the functions) */
31   int from = (rank - 1) / 2;
32   int to_left = rank * 2 + 1;
33   int to_right = rank * 2 + 2;
34   if (to_left >= size)
35     to_left = -1;
36   if (to_right >= size)
37     to_right = -1;
38
39   /* segment is segment size in number of elements (not bytes) */
40   int segment = bcast_NTSB_segment_size_in_byte / extent;
41   segment =  segment == 0 ? 1 :segment; 
42   /* pipeline length */
43   int pipe_length = count / segment;
44
45   /* use for buffer offset for sending and receiving data = segment size in byte */
46   int increment = segment * extent;
47
48   /* if the input size is not divisible by segment size => 
49      the small remainder will be done with native implementation */
50   int remainder = count % segment;
51
52   /* if root is not zero send to rank zero first */
53   if (root != 0) {
54     if (rank == root) {
55       smpi_mpi_send(buf, count, datatype, 0, tag, comm);
56     } else if (rank == 0) {
57       smpi_mpi_recv(buf, count, datatype, root, tag, comm, &status);
58     }
59   }
60
61   /* when a message is smaller than a block size => no pipeline */
62   if (count <= segment) {
63
64     /* case: root */
65     if (rank == 0) {
66       /* case root has only a left child */
67       if (to_right == -1) {
68         smpi_mpi_send(buf, count, datatype, to_left, tag, comm);
69       }
70       /* case root has both left and right children */
71       else {
72         smpi_mpi_send(buf, count, datatype, to_left, tag, comm);
73         smpi_mpi_send(buf, count, datatype, to_right, tag, comm);
74       }
75     }
76
77     /* case: leaf ==> receive only */
78     else if (to_left == -1) {
79       smpi_mpi_recv(buf, count, datatype, from, tag, comm, &status);
80     }
81
82     /* case: intermidiate node with only left child ==> relay message */
83     else if (to_right == -1) {
84       smpi_mpi_recv(buf, count, datatype, from, tag, comm, &status);
85       smpi_mpi_send(buf, count, datatype, to_left, tag, comm);
86     }
87
88     /* case: intermidiate node with both left and right children ==> relay message */
89     else {
90       smpi_mpi_recv(buf, count, datatype, from, tag, comm, &status);
91       smpi_mpi_send(buf, count, datatype, to_left, tag, comm);
92       smpi_mpi_send(buf, count, datatype, to_right, tag, comm);
93     }
94     return MPI_SUCCESS;
95   }
96   // pipelining
97   else {
98
99     send_request_array =
100         (MPI_Request *) xbt_malloc(2 * (size + pipe_length) * sizeof(MPI_Request));
101     recv_request_array =
102         (MPI_Request *) xbt_malloc((size + pipe_length) * sizeof(MPI_Request));
103     send_status_array =
104         (MPI_Status *) xbt_malloc(2 * (size + pipe_length) * sizeof(MPI_Status));
105     recv_status_array =
106         (MPI_Status *) xbt_malloc((size + pipe_length) * sizeof(MPI_Status));
107
108
109
110     /* case: root */
111     if (rank == 0) {
112       /* case root has only a left child */
113       if (to_right == -1) {
114         for (i = 0; i < pipe_length; i++) {
115           send_request_array[i] = smpi_mpi_isend((char *) buf + (i * increment), segment, datatype, to_left,
116                     tag + i, comm);
117         }
118         smpi_mpi_waitall((pipe_length), send_request_array, send_status_array);
119       }
120       /* case root has both left and right children */
121       else {
122         for (i = 0; i < pipe_length; i++) {
123           send_request_array[i] = smpi_mpi_isend((char *) buf + (i * increment), segment, datatype, to_left,
124                     tag + i, comm);
125           send_request_array[i + pipe_length] = smpi_mpi_isend((char *) buf + (i * increment), segment, datatype, to_right,
126                     tag + i, comm);
127         }
128         smpi_mpi_waitall((2 * pipe_length), send_request_array, send_status_array);
129       }
130     }
131
132     /* case: leaf ==> receive only */
133     else if (to_left == -1) {
134       for (i = 0; i < pipe_length; i++) {
135         recv_request_array[i] = smpi_mpi_irecv((char *) buf + (i * increment), segment, datatype, from,
136                   tag + i, comm);
137       }
138       smpi_mpi_waitall((pipe_length), recv_request_array, recv_status_array);
139     }
140
141     /* case: intermidiate node with only left child ==> relay message */
142     else if (to_right == -1) {
143       for (i = 0; i < pipe_length; i++) {
144         recv_request_array[i] = smpi_mpi_irecv((char *) buf + (i * increment), segment, datatype, from,
145                   tag + i, comm);
146       }
147       for (i = 0; i < pipe_length; i++) {
148         smpi_mpi_wait(&recv_request_array[i], &status);
149         send_request_array[i] = smpi_mpi_isend((char *) buf + (i * increment), segment, datatype, to_left,
150                   tag + i, comm);
151       }
152       smpi_mpi_waitall(pipe_length, send_request_array, send_status_array);
153
154     }
155     /* case: intermidiate node with both left and right children ==> relay message */
156     else {
157       for (i = 0; i < pipe_length; i++) {
158         recv_request_array[i] = smpi_mpi_irecv((char *) buf + (i * increment), segment, datatype, from,
159                   tag + i, comm);
160       }
161       for (i = 0; i < pipe_length; i++) {
162         smpi_mpi_wait(&recv_request_array[i], &status);
163         send_request_array[i] = smpi_mpi_isend((char *) buf + (i * increment), segment, datatype, to_left,
164                   tag + i, comm);
165         send_request_array[i + pipe_length] = smpi_mpi_isend((char *) buf + (i * increment), segment, datatype, to_right,
166                   tag + i, comm);
167       }
168       smpi_mpi_waitall((2 * pipe_length), send_request_array, send_status_array);
169     }
170
171     free(send_request_array);
172     free(recv_request_array);
173     free(send_status_array);
174     free(recv_status_array);
175   }                             /* end pipeline */
176
177   /* when count is not divisible by block size, use default BCAST for the remainder */
178   if ((remainder != 0) && (count > segment)) {
179     XBT_WARN("MPI_bcast_NTSB use default MPI_bcast.");            
180     smpi_mpi_bcast((char *) buf + (pipe_length * increment), remainder, datatype,
181               root, comm);
182   }
183
184   return MPI_SUCCESS;
185 }