Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
use tuned algo here
[simgrid.git] / src / smpi / colls / alltoall-pair-mpi-barrier.c
1 #include "colls_private.h"
2 /*****************************************************************************
3
4  * Function: alltoall_pair_mpi_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 when P is power of two. In each phase of P - 1
18            phases, nodes in pair communicate their data. MPI barriers are
19            inserted between each two phases.
20
21  * Auther: Ahmad Faraj
22
23  ****************************************************************************/
24 int
25 smpi_coll_tuned_alltoall_pair_mpi_barrier(void *send_buff, int send_count,
26                                           MPI_Datatype send_type,
27                                           void *recv_buff, int recv_count,
28                                           MPI_Datatype recv_type, MPI_Comm comm)
29 {
30   MPI_Status s;
31   MPI_Aint send_chunk, recv_chunk;
32   int i, src, dst, rank, num_procs;
33   int tag = COLL_TAG_ALLTOALL;
34   char *send_ptr = (char *) send_buff;
35   char *recv_ptr = (char *) recv_buff;
36
37   rank = smpi_comm_rank(comm);
38   num_procs = smpi_comm_size(comm);
39   send_chunk = smpi_datatype_get_extent(send_type);
40   recv_chunk = smpi_datatype_get_extent(recv_type);
41
42   send_chunk *= send_count;
43   recv_chunk *= recv_count;
44
45   for (i = 0; i < num_procs; i++) {
46     src = dst = rank ^ i;
47     smpi_mpi_barrier(comm);
48     smpi_mpi_sendrecv(send_ptr + dst * send_chunk, send_count, send_type, dst,
49                  tag, recv_ptr + src * recv_chunk, recv_count, recv_type,
50                  src, tag, comm, &s);
51   }
52   return MPI_SUCCESS;
53 }