Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
git checkout hypervisor
[simgrid.git] / src / smpi / colls / bcast-NTSL-Isend.c
1 #include "colls.h"
2
3 static int bcast_NTSL_segment_size_in_byte = 8192;
4
5 /* Non-topology-specific pipelined linear-bcast function 
6    0->1, 1->2 ,2->3, ....., ->last node : in a pipeline fashion
7 */
8 int smpi_coll_tuned_bcast_NTSL_Isend(void *buf, int count, MPI_Datatype datatype,
9                                int root, MPI_Comm comm)
10 {
11   int tag = 50;
12   MPI_Status status;
13   MPI_Request request;
14   MPI_Request *send_request_array;
15   MPI_Request *recv_request_array;
16   MPI_Status *send_status_array;
17   MPI_Status *recv_status_array;
18   int rank, size;
19   int i;
20   MPI_Aint extent;
21   MPI_Type_extent(datatype, &extent);
22
23   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
24   MPI_Comm_size(MPI_COMM_WORLD, &size);
25
26   /* source node and destination nodes (same through out the functions) */
27   int to = (rank + 1) % size;
28   int from = (rank + size - 1) % size;
29
30   /* segment is segment size in number of elements (not bytes) */
31   int segment = bcast_NTSL_segment_size_in_byte / extent;
32
33   /* pipeline length */
34   int pipe_length = count / segment;
35
36   /* use for buffer offset for sending and receiving data = segment size in byte */
37   int increment = segment * extent;
38
39   /* if the input size is not divisible by segment size => 
40      the small remainder will be done with native implementation */
41   int remainder = count % segment;
42
43   /* if root is not zero send to rank zero first
44      this can be modified to make it faster by using logical src, dst.
45    */
46   if (root != 0) {
47     if (rank == root) {
48       MPI_Send(buf, count, datatype, 0, tag, comm);
49     } else if (rank == 0) {
50       MPI_Recv(buf, count, datatype, root, tag, comm, &status);
51     }
52   }
53
54   /* when a message is smaller than a block size => no pipeline */
55   if (count <= segment) {
56     if (rank == 0) {
57       MPI_Send(buf, count, datatype, to, tag, comm);
58     } else if (rank == (size - 1)) {
59       MPI_Irecv(buf, count, datatype, from, tag, comm, &request);
60       MPI_Wait(&request, &status);
61     } else {
62       MPI_Irecv(buf, count, datatype, from, tag, comm, &request);
63       MPI_Wait(&request, &status);
64       MPI_Send(buf, count, datatype, to, tag, comm);
65     }
66     return MPI_SUCCESS;
67   }
68
69   /* pipeline bcast */
70   else {
71     send_request_array =
72         (MPI_Request *) malloc((size + pipe_length) * sizeof(MPI_Request));
73     recv_request_array =
74         (MPI_Request *) malloc((size + pipe_length) * sizeof(MPI_Request));
75     send_status_array =
76         (MPI_Status *) malloc((size + pipe_length) * sizeof(MPI_Status));
77     recv_status_array =
78         (MPI_Status *) malloc((size + pipe_length) * sizeof(MPI_Status));
79
80     /* root send data */
81     if (rank == 0) {
82       for (i = 0; i < pipe_length; i++) {
83         MPI_Isend((char *) buf + (i * increment), segment, datatype, to,
84                   (tag + i), comm, &send_request_array[i]);
85       }
86       MPI_Waitall((pipe_length), send_request_array, send_status_array);
87     }
88
89     /* last node only receive data */
90     else if (rank == (size - 1)) {
91       for (i = 0; i < pipe_length; i++) {
92         MPI_Irecv((char *) buf + (i * increment), segment, datatype, from,
93                   (tag + i), comm, &recv_request_array[i]);
94       }
95       MPI_Waitall((pipe_length), recv_request_array, recv_status_array);
96     }
97
98     /* intermediate nodes relay (receive, then send) data */
99     else {
100       for (i = 0; i < pipe_length; i++) {
101         MPI_Irecv((char *) buf + (i * increment), segment, datatype, from,
102                   (tag + i), comm, &recv_request_array[i]);
103       }
104       for (i = 0; i < pipe_length; i++) {
105         MPI_Wait(&recv_request_array[i], &status);
106         MPI_Isend((char *) buf + (i * increment), segment, datatype, to,
107                   (tag + i), comm, &send_request_array[i]);
108       }
109       MPI_Waitall((pipe_length), send_request_array, send_status_array);
110     }
111
112     free(send_request_array);
113     free(recv_request_array);
114     free(send_status_array);
115     free(recv_status_array);
116   }                             /* end pipeline */
117
118   /* when count is not divisible by block size, use default BCAST for the remainder */
119   if ((remainder != 0) && (count > segment)) {
120     MPI_Bcast((char *) buf + (pipe_length * increment), remainder, datatype,
121               root, comm);
122   }
123
124   return MPI_SUCCESS;
125 }