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 /*****************************************************************************
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 * Auther: MPICH / modified by Ahmad Faraj
26 ****************************************************************************/
28 smpi_coll_tuned_alltoall_bruck(void *send_buff, int send_count,
29 MPI_Datatype send_type, void *recv_buff,
30 int recv_count, MPI_Datatype recv_type,
35 MPI_Datatype new_type;
37 int *blocks_length, *disps;
38 int i, src, dst, rank, num_procs, count, remainder, block, position;
39 int pack_size, tag = COLL_TAG_ALLTOALL, pof2 = 1;
43 char *send_ptr = (char *) send_buff;
44 char *recv_ptr = (char *) recv_buff;
46 num_procs = smpi_comm_size(comm);
47 rank = smpi_comm_rank(comm);
49 extent = smpi_datatype_get_extent(recv_type);
51 tmp_buff = (char *) smpi_get_tmp_sendbuffer(num_procs * recv_count * extent);
52 disps = (int *) xbt_malloc(sizeof(int) * num_procs);
53 blocks_length = (int *) xbt_malloc(sizeof(int) * num_procs);
55 smpi_mpi_sendrecv(send_ptr + rank * send_count * extent,
56 (num_procs - rank) * send_count, send_type, rank, tag,
57 recv_ptr, (num_procs - rank) * recv_count, recv_type, rank,
60 smpi_mpi_sendrecv(send_ptr, rank * send_count, send_type, rank, tag,
61 recv_ptr + (num_procs - rank) * recv_count * extent,
62 rank * recv_count, recv_type, rank, tag, comm, &status);
66 MPI_Pack_size(send_count * num_procs, send_type, comm, &pack_size);
68 while (pof2 < num_procs) {
69 dst = (rank + pof2) % num_procs;
70 src = (rank - pof2 + num_procs) % num_procs;
74 for (block = 1; block < num_procs; block++)
76 blocks_length[count] = send_count;
77 disps[count] = block * send_count;
81 MPI_Type_indexed(count, blocks_length, disps, recv_type, &new_type);
82 smpi_datatype_commit(&new_type);
85 MPI_Pack(recv_buff, 1, new_type, tmp_buff, pack_size, &position, comm);
87 smpi_mpi_sendrecv(tmp_buff, position, MPI_PACKED, dst, tag, recv_buff, 1,
88 new_type, src, tag, comm, &status);
89 smpi_datatype_free(&new_type);
97 smpi_mpi_sendrecv(recv_ptr + (rank + 1) * recv_count * extent,
98 (num_procs - rank - 1) * recv_count, send_type,
99 rank, tag, tmp_buff, (num_procs - rank - 1) * recv_count,
100 recv_type, rank, tag, comm, &status);
102 smpi_mpi_sendrecv(recv_ptr, (rank + 1) * recv_count, send_type, rank, tag,
103 tmp_buff + (num_procs - rank - 1) * recv_count * extent,
104 (rank + 1) * recv_count, recv_type, rank, tag, comm, &status);
107 for (i = 0; i < num_procs; i++)
108 smpi_mpi_sendrecv(tmp_buff + i * recv_count * extent, recv_count, send_type,
110 recv_ptr + (num_procs - i - 1) * recv_count * extent,
111 recv_count, recv_type, rank, tag, comm, &status);
113 smpi_free_tmp_buffer(tmp_buff);