1 /* Copyright (c) 2013-2019. 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 /*****************************************************************************
9 * Function: alltoall_bruck
14 send_buff: send input buffer
15 send_count: number of elements to send
16 send_type: data type of elements being sent
17 recv_buff: receive output buffer
18 recv_count: number of elements to received
19 recv_type: data type of elements being received
22 * Descrp: Function realizes the alltoall operation using the bruck algorithm.
24 * Author: MPICH / modified by Ahmad Faraj
26 ****************************************************************************/
28 #include "../colls_private.hpp"
35 alltoall__bruck(const void *send_buff, int send_count,
36 MPI_Datatype send_type, void *recv_buff,
37 int recv_count, MPI_Datatype recv_type,
42 MPI_Datatype new_type;
44 int i, src, dst, rank, num_procs, count, block, position;
45 int pack_size, tag = COLL_TAG_ALLTOALL, pof2 = 1;
47 char *send_ptr = (char *) send_buff;
48 char *recv_ptr = (char *) recv_buff;
50 num_procs = comm->size();
53 extent = recv_type->get_extent();
55 unsigned char* tmp_buff = smpi_get_tmp_sendbuffer(num_procs * recv_count * extent);
56 int* disps = new int[num_procs];
57 int* blocks_length = new int[num_procs];
59 Request::sendrecv(send_ptr + rank * send_count * extent,
60 (num_procs - rank) * send_count, send_type, rank, tag,
61 recv_ptr, (num_procs - rank) * recv_count, recv_type, rank,
64 Request::sendrecv(send_ptr, rank * send_count, send_type, rank, tag,
65 recv_ptr + (num_procs - rank) * recv_count * extent,
66 rank * recv_count, recv_type, rank, tag, comm, &status);
70 MPI_Pack_size(send_count * num_procs, send_type, comm, &pack_size);
72 while (pof2 < num_procs) {
73 dst = (rank + pof2) % num_procs;
74 src = (rank - pof2 + num_procs) % num_procs;
78 for (block = 1; block < num_procs; block++)
80 blocks_length[count] = send_count;
81 disps[count] = block * send_count;
85 MPI_Type_indexed(count, blocks_length, disps, recv_type, &new_type);
89 MPI_Pack(recv_buff, 1, new_type, tmp_buff, pack_size, &position, comm);
91 Request::sendrecv(tmp_buff, position, MPI_PACKED, dst, tag, recv_buff, 1,
92 new_type, src, tag, comm, &status);
93 Datatype::unref(new_type);
99 delete[] blocks_length;
101 Request::sendrecv(recv_ptr + (rank + 1) * recv_count * extent,
102 (num_procs - rank - 1) * recv_count, send_type,
103 rank, tag, tmp_buff, (num_procs - rank - 1) * recv_count,
104 recv_type, rank, tag, comm, &status);
106 Request::sendrecv(recv_ptr, (rank + 1) * recv_count, send_type, rank, tag,
107 tmp_buff + (num_procs - rank - 1) * recv_count * extent,
108 (rank + 1) * recv_count, recv_type, rank, tag, comm, &status);
111 for (i = 0; i < num_procs; i++)
112 Request::sendrecv(tmp_buff + i * recv_count * extent, recv_count, send_type,
114 recv_ptr + (num_procs - i - 1) * recv_count * extent,
115 recv_count, recv_type, rank, tag, comm, &status);
117 smpi_free_tmp_buffer(tmp_buff);