1 /* Copyright (c) 2013-2014. The SimGrid Team.
2 * All rights reserved. */
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. */
7 #include "colls_private.h"
9 int bcast_NTSB_segment_size_in_byte = 8192;
11 int smpi_coll_tuned_bcast_NTSB(void *buf, int count, MPI_Datatype datatype,
12 int root, MPI_Comm comm)
14 int tag = COLL_TAG_BCAST;
19 MPI_Request *send_request_array;
20 MPI_Request *recv_request_array;
21 MPI_Status *send_status_array;
22 MPI_Status *recv_status_array;
25 extent = smpi_datatype_get_extent(datatype);
27 rank = smpi_comm_rank(comm);
28 size = smpi_comm_size(comm);
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;
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;
43 int pipe_length = count / segment;
45 /* use for buffer offset for sending and receiving data = segment size in byte */
46 int increment = segment * extent;
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;
52 /* if root is not zero send to rank zero first */
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);
61 /* when a message is smaller than a block size => no pipeline */
62 if (count <= segment) {
66 /* case root has only a left child */
68 smpi_mpi_send(buf, count, datatype, to_left, tag, comm);
70 /* case root has both left and right children */
72 smpi_mpi_send(buf, count, datatype, to_left, tag, comm);
73 smpi_mpi_send(buf, count, datatype, to_right, tag, comm);
77 /* case: leaf ==> receive only */
78 else if (to_left == -1) {
79 smpi_mpi_recv(buf, count, datatype, from, tag, comm, &status);
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);
88 /* case: intermidiate node with both left and right children ==> relay message */
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);
100 (MPI_Request *) xbt_malloc(2 * (size + pipe_length) * sizeof(MPI_Request));
102 (MPI_Request *) xbt_malloc((size + pipe_length) * sizeof(MPI_Request));
104 (MPI_Status *) xbt_malloc(2 * (size + pipe_length) * sizeof(MPI_Status));
106 (MPI_Status *) xbt_malloc((size + pipe_length) * sizeof(MPI_Status));
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,
118 smpi_mpi_waitall((pipe_length), send_request_array, send_status_array);
120 /* case root has both left and right children */
122 for (i = 0; i < pipe_length; i++) {
123 send_request_array[i] = smpi_mpi_isend((char *) buf + (i * increment), segment, datatype, to_left,
125 send_request_array[i + pipe_length] = smpi_mpi_isend((char *) buf + (i * increment), segment, datatype, to_right,
128 smpi_mpi_waitall((2 * pipe_length), send_request_array, send_status_array);
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,
138 smpi_mpi_waitall((pipe_length), recv_request_array, recv_status_array);
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,
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,
152 smpi_mpi_waitall(pipe_length, send_request_array, send_status_array);
155 /* case: intermidiate node with both left and right children ==> relay message */
157 for (i = 0; i < pipe_length; i++) {
158 recv_request_array[i] = smpi_mpi_irecv((char *) buf + (i * increment), segment, datatype, from,
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,
165 send_request_array[i + pipe_length] = smpi_mpi_isend((char *) buf + (i * increment), segment, datatype, to_right,
168 smpi_mpi_waitall((2 * pipe_length), send_request_array, send_status_array);
171 free(send_request_array);
172 free(recv_request_array);
173 free(send_status_array);
174 free(recv_status_array);
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,