Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Define buffer on the stack here.
[simgrid.git] / src / smpi / colls / alltoall-pair-light-barrier.c
1 #include "colls_private.h"
2 /*****************************************************************************
3
4  * Function: alltoall_pair_light_barrier
5
6  * Return: int
7
8  * Inputs:
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
15     comm: communicator
16
17  * Descrp: Function works in P - 1 steps. In step i, node j exchanges data
18            with node i ^ j. Light barriers are inserted between
19            communications in different phases.
20
21  * Auther: Ahmad Faraj
22
23  ****************************************************************************/
24 int
25 smpi_coll_tuned_alltoall_pair_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,
29                                             MPI_Comm comm)
30 {
31   MPI_Aint send_chunk, recv_chunk;
32   MPI_Status s;
33   int i, src, dst, rank, num_procs, next_partner;
34   int tag = COLL_TAG_ALLTOALL;     /*, failure = 0; */
35
36   char send_sync = 'a', recv_sync = 'b';
37   char *send_ptr = (char *) send_buff;
38   char *recv_ptr = (char *) recv_buff;
39
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);
44
45   send_chunk *= send_count;
46   recv_chunk *= recv_count;
47
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,
50                comm, &s);
51
52   for (i = 1; i < num_procs; i++) {
53     src = dst = rank ^ i;
54
55     smpi_mpi_sendrecv(send_ptr + dst * send_chunk, send_count, send_type,
56                  dst, tag, recv_ptr + src * recv_chunk, recv_count,
57                  recv_type, src, tag, comm, &s);
58
59     if ((i + 1) < num_procs) {
60       next_partner = rank ^ (i + 1);
61       smpi_mpi_sendrecv(&send_sync, 1, MPI_CHAR, next_partner, tag,
62                    &recv_sync, 1, MPI_CHAR, next_partner, tag, comm, &s);
63     }
64   }
65   return MPI_SUCCESS;
66 }