Logo AND Algorithmique Numérique Distribuée

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