Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
5d5f94537d55bd7f6dec16b1be1bd660cfe48b44
[simgrid.git] / src / smpi / colls / allgatherv / allgatherv-mpich-rdb.cpp
1 /* Copyright (c) 2013-2014. The SimGrid Team.
2  * All rights reserved.                                                     */
3
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. */
6
7         /* Short or medium size message and power-of-two no. of processes. Use
8          * recursive doubling algorithm */
9 #include "../colls_private.h"
10
11 namespace simgrid{
12 namespace smpi{
13
14 int Coll_allgatherv_mpich_rdb::allgatherv (
15   void *sendbuf,
16   int sendcount,
17   MPI_Datatype sendtype,
18   void *recvbuf,
19   int *recvcounts,
20   int *displs,
21   MPI_Datatype recvtype,
22   MPI_Comm comm)
23 {
24   unsigned int  j, i;
25   MPI_Status status;
26   MPI_Aint  recvtype_extent, recvtype_true_extent, recvtype_true_lb;
27   unsigned int curr_cnt, dst, total_count;
28   void *tmp_buf, *tmp_buf_rl;
29   unsigned int mask, dst_tree_root, my_tree_root, position,
30     send_offset, recv_offset, last_recv_cnt=0, nprocs_completed, k,
31     offset, tmp_mask, tree_root;
32
33   unsigned int comm_size = comm->size();
34   unsigned int rank = comm->rank();
35
36   total_count = 0;
37   for (i=0; i<comm_size; i++)
38     total_count += recvcounts[i];
39
40   if (total_count == 0)
41     return MPI_ERR_COUNT;
42
43   recvtype_extent=recvtype->get_extent();
44
45   /* need to receive contiguously into tmp_buf because
46      displs could make the recvbuf noncontiguous */
47
48   recvtype->extent(&recvtype_true_lb, &recvtype_true_extent);
49
50   tmp_buf_rl= (void*)smpi_get_tmp_sendbuffer(total_count*(MAX(recvtype_true_extent,recvtype_extent)));
51
52   /* adjust for potential negative lower bound in datatype */
53   tmp_buf = (void *)((char*)tmp_buf_rl - recvtype_true_lb);
54
55   /* copy local data into right location in tmp_buf */
56   position = 0;
57   for (i=0; i<rank; i++)
58     position += recvcounts[i];
59   if (sendbuf != MPI_IN_PLACE)
60   {
61     Datatype::copy(sendbuf, sendcount, sendtype,
62                        ((char *)tmp_buf + position*
63                         recvtype_extent),
64                        recvcounts[rank], recvtype);
65   }
66   else
67   {
68     /* if in_place specified, local data is found in recvbuf */
69     Datatype::copy(((char *)recvbuf +
70                         displs[rank]*recvtype_extent),
71                        recvcounts[rank], recvtype,
72                        ((char *)tmp_buf + position*
73                         recvtype_extent),
74                        recvcounts[rank], recvtype);
75   }
76   curr_cnt = recvcounts[rank];
77
78   mask = 0x1;
79   i = 0;
80   while (mask < comm_size) {
81     dst = rank ^ mask;
82
83     /* find offset into send and recv buffers. zero out
84        the least significant "i" bits of rank and dst to
85        find root of src and dst subtrees. Use ranks of
86        roots as index to send from and recv into buffer */
87
88     dst_tree_root = dst >> i;
89     dst_tree_root <<= i;
90
91     my_tree_root = rank >> i;
92     my_tree_root <<= i;
93
94     if (dst < comm_size) {
95       send_offset = 0;
96       for (j=0; j<my_tree_root; j++)
97         send_offset += recvcounts[j];
98
99       recv_offset = 0;
100       for (j=0; j<dst_tree_root; j++)
101         recv_offset += recvcounts[j];
102
103       Request::sendrecv(((char *)tmp_buf + send_offset * recvtype_extent),
104                         curr_cnt, recvtype, dst,
105                         COLL_TAG_ALLGATHERV,
106                         ((char *)tmp_buf + recv_offset * recvtype_extent),
107                         total_count - recv_offset, recvtype, dst,
108                         COLL_TAG_ALLGATHERV,
109                         comm, &status);
110       /* for convenience, recv is posted for a bigger amount
111          than will be sent */
112       last_recv_cnt=smpi_mpi_get_count(&status, recvtype);
113       curr_cnt += last_recv_cnt;
114     }
115
116     /* if some processes in this process's subtree in this step
117        did not have any destination process to communicate with
118        because of non-power-of-two, we need to send them the
119        data that they would normally have received from those
120        processes. That is, the haves in this subtree must send to
121        the havenots. We use a logarithmic
122        recursive-halfing algorithm for this. */
123
124     /* This part of the code will not currently be
125        executed because we are not using recursive
126        doubling for non power of two. Mark it as experimental
127        so that it doesn't show up as red in the coverage
128        tests. */
129
130     /* --BEGIN EXPERIMENTAL-- */
131     if (dst_tree_root + mask > comm_size) {
132       nprocs_completed = comm_size - my_tree_root - mask;
133       /* nprocs_completed is the number of processes in this
134          subtree that have all the data. Send data to others
135          in a tree fashion. First find root of current tree
136          that is being divided into two. k is the number of
137          least-significant bits in this process's rank that
138          must be zeroed out to find the rank of the root */
139       j = mask;
140       k = 0;
141       while (j) {
142         j >>= 1;
143         k++;
144       }
145       k--;
146
147       tmp_mask = mask >> 1;
148
149       while (tmp_mask) {
150         dst = rank ^ tmp_mask;
151
152         tree_root = rank >> k;
153         tree_root <<= k;
154
155         /* send only if this proc has data and destination
156            doesn't have data. at any step, multiple processes
157            can send if they have the data */
158         if ((dst > rank) &&
159             (rank < tree_root + nprocs_completed)
160             && (dst >= tree_root + nprocs_completed)) {
161
162           offset = 0;
163           for (j=0; j<(my_tree_root+mask); j++)
164             offset += recvcounts[j];
165           offset *= recvtype_extent;
166
167           Request::send(((char *)tmp_buf + offset),
168                         last_recv_cnt,
169                         recvtype, dst,
170                         COLL_TAG_ALLGATHERV, comm);
171           /* last_recv_cnt was set in the previous
172              receive. that's the amount of data to be
173              sent now. */
174         }
175         /* recv only if this proc. doesn't have data and sender
176            has data */
177         else if ((dst < rank) &&
178                  (dst < tree_root + nprocs_completed) &&
179                  (rank >= tree_root + nprocs_completed)) {
180
181           offset = 0;
182           for (j=0; j<(my_tree_root+mask); j++)
183             offset += recvcounts[j];
184
185           Request::recv(((char *)tmp_buf + offset * recvtype_extent),
186                         total_count - offset, recvtype,
187                         dst, COLL_TAG_ALLGATHERV,
188                         comm, &status);
189           /* for convenience, recv is posted for a
190              bigger amount than will be sent */
191           last_recv_cnt=smpi_mpi_get_count(&status, recvtype);
192           curr_cnt += last_recv_cnt;
193         }
194         tmp_mask >>= 1;
195         k--;
196       }
197     }
198     /* --END EXPERIMENTAL-- */
199
200     mask <<= 1;
201     i++;
202   }
203
204   /* copy data from tmp_buf to recvbuf */
205   position = 0;
206   for (j=0; j<comm_size; j++) {
207     if ((sendbuf != MPI_IN_PLACE) || (j != rank)) {
208       /* not necessary to copy if in_place and
209          j==rank. otherwise copy. */
210       Datatype::copy(((char *)tmp_buf + position*recvtype_extent),
211                          recvcounts[j], recvtype,
212                          ((char *)recvbuf + displs[j]*recvtype_extent),
213                          recvcounts[j], recvtype);
214     }
215     position += recvcounts[j];
216   }
217
218   smpi_free_tmp_buffer(tmp_buf_rl);
219   return MPI_SUCCESS;
220 }
221
222 }
223 }