1 /* Copyright (c) 2013-2014. The SimGrid Team.
2 * All rights reserved. */
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. */
7 #include "colls_private.h"
9 /*****************************************************************************
11 Copyright (c) 2006, Ahmad Faraj & Xin Yuan,
14 Redistribution and use in source and binary forms, with or without
15 modification, are permitted provided that the following conditions are met:
17 * Redistributions of source code must retain the above copyright notice,
18 this list of conditions and the following disclaimer.
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.
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.
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.
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: *
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, *
48 *************************************************************************
50 *****************************************************************************/
52 /*****************************************************************************
55 * num: the number of processors in a communicator
59 * descp: takes a number and tries to find a factoring of x*y*z mesh out of it
60 ****************************************************************************/
63 static int is_3dmesh(int num, int *i, int *j, int *k)
69 if ((num % (x * x)) == 0) {
79 /*****************************************************************************
80 * Function: allgather_3dmesh_shoot
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
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
94 ****************************************************************************/
95 int smpi_coll_tuned_allgather_3dmesh(void *send_buff, int send_count,
96 MPI_Datatype send_type, void *recv_buff,
97 int recv_count, MPI_Datatype recv_type,
100 MPI_Request *req, *req_ptr;
103 int i, src, dst, rank, num_procs, block_size, my_z_base;
104 int my_z, X, Y, Z, send_offset, recv_offset;
105 int two_dsize, my_row_base, my_col_base, src_row_base, src_z_base, num_reqs;
106 int tag = COLL_TAG_ALLGATHER;
108 rank = smpi_comm_rank(comm);
109 num_procs = smpi_comm_size(comm);
110 extent = smpi_datatype_get_extent(send_type);
112 if (!is_3dmesh(num_procs, &X, &Y, &Z))
113 THROWF(arg_error,0, "allgather_3dmesh algorithm can't be used with this number of processes! ");
124 my_z = rank / two_dsize;
126 my_row_base = (rank / X) * X;
127 my_col_base = (rank % Y) + (my_z * two_dsize);
128 my_z_base = my_z * two_dsize;
130 block_size = extent * send_count;
132 req = (MPI_Request *) xbt_malloc(num_reqs * sizeof(MPI_Request));
136 // do local allgather/local copy
137 recv_offset = rank * block_size;
138 smpi_datatype_copy(send_buff, send_count, send_type, (char *)recv_buff + recv_offset,
139 recv_count, recv_type);
142 for (i = 0; i < Y; i++) {
143 src = i + my_row_base;
146 recv_offset = src * block_size;
147 *(req_ptr++) = smpi_mpi_irecv((char *)recv_buff + recv_offset, send_count, recv_type, src, tag,
151 for (i = 0; i < Y; i++) {
152 dst = i + my_row_base;
155 smpi_mpi_send(send_buff, send_count, send_type, dst, tag, comm);
158 smpi_mpi_waitall(Y - 1, req, MPI_STATUSES_IGNORE);
161 // do colwise comm, it does not matter here if i*X or i *Y since X == Y
163 for (i = 0; i < X; i++) {
164 src = (i * Y + my_col_base);
168 src_row_base = (src / X) * X;
169 recv_offset = src_row_base * block_size;
170 *(req_ptr++) = smpi_mpi_irecv((char *)recv_buff + recv_offset, recv_count * Y, recv_type, src, tag,
174 send_offset = my_row_base * block_size;
176 for (i = 0; i < X; i++) {
177 dst = (i * Y + my_col_base);
180 smpi_mpi_send((char *)recv_buff + send_offset, send_count * Y, send_type, dst, tag,
184 smpi_mpi_waitall(X - 1, req, MPI_STATUSES_IGNORE);
187 for (i = 1; i < Z; i++) {
188 src = (rank + i * two_dsize) % num_procs;
189 src_z_base = (src / two_dsize) * two_dsize;
191 recv_offset = (src_z_base * block_size);
193 *(req_ptr++) = smpi_mpi_irecv((char *)recv_buff + recv_offset, recv_count * two_dsize, recv_type,
197 for (i = 1; i < Z; i++) {
198 dst = (rank + i * two_dsize) % num_procs;
199 send_offset = my_z_base * block_size;
200 smpi_mpi_send((char *)recv_buff + send_offset, send_count * two_dsize, send_type,
203 smpi_mpi_waitall(Z - 1, req, MPI_STATUSES_IGNORE);