Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
SMPI colls in not really C++. But cleaner than before.
[simgrid.git] / src / smpi / colls / allgather / allgather-3dmesh.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 #include "../colls_private.h"
8
9 /*****************************************************************************
10
11 Copyright (c) 2006, Ahmad Faraj & Xin Yuan,
12 All rights reserved.
13
14 Redistribution and use in source and binary forms, with or without
15 modification, are permitted provided that the following conditions are met:
16
17   * Redistributions of source code must retain the above copyright notice,
18     this list of conditions and the following disclaimer.
19
20   * Redistributions in binary form must reproduce the above copyright notice,
21     this list of conditions and the following disclaimer in the documentation
22     and/or other materials provided with the distribution.
23
24   * Neither the name of the Florida State University nor the names of its
25     contributors may be used to endorse or promote products derived from this
26     software without specific prior written permission.
27
28 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
29 ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
30 WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
31 DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR
32 ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
33 (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
34 LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON
35 ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
36 (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
37 SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
38
39   *************************************************************************
40   *     Any results obtained from executing this software require the     *
41   *     acknowledgment and citation of the software and its owners.       *
42   *     The full citation is given below:                                 *
43   *                                                                       *
44   *     A. Faraj and X. Yuan. "Automatic Generation and Tuning of MPI     *
45   *     Collective Communication Routines." The 19th ACM International    *
46   *     Conference on Supercomputing (ICS), Cambridge, Massachusetts,     *
47   *     June 20-22, 2005.                                                 *
48   *************************************************************************
49
50 *****************************************************************************/
51
52 /*****************************************************************************
53  * Function: is_2dmesh
54  * return: int
55  * num: the number of processors in a communicator
56  * i: x dimension
57  * j: y dimension
58  * k: z dimension
59  * descp: takes a number and tries to find a factoring of x*y*z mesh out of it
60  ****************************************************************************/
61 #ifndef THREED
62 #define THREED
63 static int is_3dmesh(int num, int *i, int *j, int *k)
64 {
65   int x, max = num / 3;
66   x = cbrt(num);
67   *i = *j = *k = 0;
68   while (x <= max) {
69     if ((num % (x * x)) == 0) {
70       *i = *j = x;
71       *k = num / (x * x);
72       return 1;
73     }
74     x++;
75   }
76   return 0;
77 }
78 #endif
79 /*****************************************************************************
80  * Function: allgather_3dmesh_shoot
81  * return: int
82  * send_buff: send input buffer
83  * send_count: number of elements to send
84  * send_type: data type of elements being sent
85  * recv_buff: receive output buffer
86  * recv_count: number of elements to received
87  * recv_type: data type of elements being received
88  * comm: communication
89  * Descrp: Function realizes the allgather operation using the 2dmesh
90  * algorithm. Allgather ommunication occurs first in the x dimension, y
91  * dimension, and then in the z dimension. Communication in each dimension
92  * follows "simple"
93  * Auther: Ahmad Faraj
94 ****************************************************************************/
95 namespace simgrid{
96 namespace smpi{
97
98
99 int Coll_allgather_3dmesh::allgather(void *send_buff, int send_count,
100                                      MPI_Datatype send_type, void *recv_buff,
101                                      int recv_count, MPI_Datatype recv_type,
102                                      MPI_Comm comm)
103 {
104   MPI_Request *req, *req_ptr;
105   MPI_Aint extent;
106
107   int i, src, dst, rank, num_procs, block_size, my_z_base;
108   int my_z, X, Y, Z, send_offset, recv_offset;
109   int two_dsize, my_row_base, my_col_base, src_row_base, src_z_base, num_reqs;
110   int tag = COLL_TAG_ALLGATHER;
111
112   rank = comm->rank();
113   num_procs = comm->size();
114   extent = send_type->get_extent();
115
116   if (!is_3dmesh(num_procs, &X, &Y, &Z))
117     THROWF(arg_error,0, "allgather_3dmesh algorithm can't be used with this number of processes! ");
118
119
120   num_reqs = X;
121
122   if (Y > X)
123     num_reqs = Y;
124   if (Z > Y)
125     num_reqs = Z;
126
127   two_dsize = X * Y;
128   my_z = rank / two_dsize;
129
130   my_row_base = (rank / X) * X;
131   my_col_base = (rank % Y) + (my_z * two_dsize);
132   my_z_base = my_z * two_dsize;
133
134   block_size = extent * send_count;
135
136   req = (MPI_Request *) xbt_malloc(num_reqs * sizeof(MPI_Request));
137
138   req_ptr = req;
139
140   // do local allgather/local copy 
141   recv_offset = rank * block_size;
142   Datatype::copy(send_buff, send_count, send_type, (char *)recv_buff + recv_offset,
143                  recv_count, recv_type);
144
145   // do rowwise comm 
146   for (i = 0; i < Y; i++) {
147     src = i + my_row_base;
148     if (src == rank)
149       continue;
150     recv_offset = src * block_size;
151     *(req_ptr++) = Request::irecv((char *)recv_buff + recv_offset, send_count, recv_type, src, tag,
152                comm);
153   }
154
155   for (i = 0; i < Y; i++) {
156     dst = i + my_row_base;
157     if (dst == rank)
158       continue;
159     Request::send(send_buff, send_count, send_type, dst, tag, comm);
160   }
161
162   Request::waitall(Y - 1, req, MPI_STATUSES_IGNORE);
163   req_ptr = req;
164
165   // do colwise comm, it does not matter here if i*X or i *Y since X == Y
166
167   for (i = 0; i < X; i++) {
168     src = (i * Y + my_col_base);
169     if (src == rank)
170       continue;
171
172     src_row_base = (src / X) * X;
173     recv_offset = src_row_base * block_size;
174     *(req_ptr++) = Request::irecv((char *)recv_buff + recv_offset, recv_count * Y, recv_type, src, tag,
175                comm);
176   }
177
178   send_offset = my_row_base * block_size;
179
180   for (i = 0; i < X; i++) {
181     dst = (i * Y + my_col_base);
182     if (dst == rank)
183       continue;
184     Request::send((char *)recv_buff + send_offset, send_count * Y, send_type, dst, tag,
185               comm);
186   }
187
188   Request::waitall(X - 1, req, MPI_STATUSES_IGNORE);
189   req_ptr = req;
190
191   for (i = 1; i < Z; i++) {
192     src = (rank + i * two_dsize) % num_procs;
193     src_z_base = (src / two_dsize) * two_dsize;
194
195     recv_offset = (src_z_base * block_size);
196
197     *(req_ptr++) = Request::irecv((char *)recv_buff + recv_offset, recv_count * two_dsize, recv_type,
198                src, tag, comm);
199   }
200
201   for (i = 1; i < Z; i++) {
202     dst = (rank + i * two_dsize) % num_procs;
203     send_offset = my_z_base * block_size;
204     Request::send((char *)recv_buff + send_offset, send_count * two_dsize, send_type,
205               dst, tag, comm);
206   }
207   Request::waitall(Z - 1, req, MPI_STATUSES_IGNORE);
208
209   free(req);
210
211   return MPI_SUCCESS;
212 }
213
214
215 }
216 }