Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Add collectives for allgather, allreduce, bcast and reduce
[simgrid.git] / src / smpi / colls / allgather-RDB.c
1 #include "colls.h"
2
3 int
4 smpi_coll_tuned_allgather_RDB(void *sbuf, int send_count,
5                                               MPI_Datatype send_type,
6                                               void *rbuf, int recv_count,
7                                               MPI_Datatype recv_type,
8                                               MPI_Comm comm)
9 {
10   // MPI variables
11   MPI_Status status;
12   MPI_Aint send_chunk, recv_chunk;
13
14   // local int variables
15   int i, j, k, dst, rank, num_procs, send_offset, recv_offset, tree_root;
16   int dst_tree_root, rank_tree_root, last_recv_count, num_procs_completed;
17   int offset, tmp_mask;
18   int tag = 1;
19   int mask = 1;
20   int success = 0;
21   int curr_count = recv_count;
22
23   // local string variables
24   char *send_ptr = (char *) sbuf;
25   char *recv_ptr = (char *) rbuf;
26
27   // get size of the communicator, followed by rank 
28   MPI_Comm_size(comm, &num_procs);
29   MPI_Comm_rank(comm, &rank);
30
31   // get size of single element's type for send buffer and recv buffer
32   MPI_Type_extent(send_type, &send_chunk);
33   MPI_Type_extent(recv_type, &recv_chunk);
34
35   // multiply size of each element by number of elements to send or recv
36   send_chunk *= send_count;
37   recv_chunk *= recv_count;
38
39   // perform a local copy
40   MPI_Sendrecv(send_ptr, send_count, send_type, rank, tag,
41                recv_ptr + rank * recv_chunk, recv_count, recv_type, rank, tag,
42                comm, &status);
43
44   i = 0;
45   while (mask < num_procs) {
46     dst = rank ^ mask;
47     dst_tree_root = dst >> i;
48     dst_tree_root <<= i;
49     rank_tree_root = rank >> i;
50     rank_tree_root <<= i;
51     send_offset = rank_tree_root * send_chunk;
52     recv_offset = dst_tree_root * recv_chunk;
53
54     if (dst < num_procs) {
55       MPI_Sendrecv(recv_ptr + send_offset, curr_count, send_type, dst,
56                    tag, recv_ptr + recv_offset, mask * recv_count,
57                    recv_type, dst, tag, comm, &status);
58       MPI_Get_count(&status, recv_type, &last_recv_count);
59       curr_count += last_recv_count;
60     }
61
62     if (dst_tree_root + mask > num_procs) {
63       num_procs_completed = num_procs - rank_tree_root - mask;
64       /* num_procs_completed is the number of processes in this
65          subtree that have all the data. Send data to others
66          in a tree fashion. First find root of current tree
67          that is being divided into two. k is the number of
68          least-significant bits in this process's rank that
69          must be zeroed out to find the rank of the root */
70
71       j = mask;
72       k = 0;
73       while (j) {
74         j >>= 1;
75         k++;
76       }
77       k--;
78
79       offset = recv_chunk * (rank_tree_root + mask);
80       tmp_mask = mask >> 1;
81
82       while (tmp_mask) {
83         dst = rank ^ tmp_mask;
84
85         tree_root = rank >> k;
86         tree_root <<= k;
87
88         /* send only if this proc has data and destination
89            doesn't have data. at any step, multiple processes
90            can send if they have the data */
91         if ((dst > rank)
92             && (rank < tree_root + num_procs_completed)
93             && (dst >= tree_root + num_procs_completed)) {
94           MPI_Send(recv_ptr + offset, last_recv_count, recv_type, dst,
95                    tag, comm);
96
97           /* last_recv_cnt was set in the previous
98              receive. that's the amount of data to be
99              sent now. */
100         }
101         /* recv only if this proc. doesn't have data and sender
102            has data */
103         else if ((dst < rank)
104                  && (dst < tree_root + num_procs_completed)
105                  && (rank >= tree_root + num_procs_completed)) {
106           MPI_Recv(recv_ptr + offset,
107                    recv_count * num_procs_completed,
108                    recv_type, dst, tag, comm, &status);
109           // num_procs_completed is also equal to the no. of processes
110           // whose data we don't have
111           MPI_Get_count(&status, recv_type, &last_recv_count);
112           curr_count += last_recv_count;
113         }
114         tmp_mask >>= 1;
115         k--;
116       }
117     }
118
119     mask <<= 1;
120     i++;
121   }
122
123   return success;
124 }