Logo AND Algorithmique Numérique Distribuée

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