1 #include "colls_private.h"
2 /*****************************************************************************
4 * Function: alltoall_ring_light_barrier
9 send_buff: send input buffer
10 send_count: number of elements to send
11 send_type: data type of elements being sent
12 recv_buff: receive output buffer
13 recv_count: number of elements to received
14 recv_type: data type of elements being received
17 * Descrp: Function works in P - 1 steps. In step i, node j - i -> j -> j + i.
18 Light barriers are inserted between communications in different
23 ****************************************************************************/
25 smpi_coll_tuned_alltoall_ring_light_barrier(void *send_buff, int send_count,
26 MPI_Datatype send_type,
27 void *recv_buff, int recv_count,
28 MPI_Datatype recv_type,
31 MPI_Aint send_chunk, recv_chunk;
33 int i, src, dst, rank, num_procs, next_dst, next_src;
34 int tag = COLL_TAG_ALLTOALL;
36 char send_sync = 'a', recv_sync = 'b';
37 char *send_ptr = (char *) send_buff;
38 char *recv_ptr = (char *) recv_buff;
40 rank = smpi_comm_rank(comm);
41 num_procs = smpi_comm_size(comm);
42 send_chunk = smpi_datatype_get_extent(send_type);
43 recv_chunk = smpi_datatype_get_extent(recv_type);
45 send_chunk *= send_count;
46 recv_chunk *= recv_count;
48 smpi_mpi_sendrecv(send_ptr + rank * send_chunk, send_count, send_type, rank, tag,
49 recv_ptr + rank * recv_chunk, recv_count, recv_type, rank, tag,
52 for (i = 1; i < num_procs; i++) {
53 src = (rank - i + num_procs) % num_procs;
54 dst = (rank + i) % num_procs;
56 smpi_mpi_sendrecv(send_ptr + dst * send_chunk, send_count, send_type,
57 dst, tag, recv_ptr + src * recv_chunk, recv_count,
58 recv_type, src, tag, comm, &s);
60 if ((i + 1) < num_procs) {
61 next_src = (rank - (i + 1) + num_procs) % num_procs;
62 next_dst = (rank + (i + 1) + num_procs) % num_procs;
63 smpi_mpi_sendrecv(&send_sync, 1, MPI_CHAR, next_src, tag,
64 &recv_sync, 1, MPI_CHAR, next_dst, tag, comm, &s);