Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Add mpi alltoallv pair_light_barrier pair_mpi_barrier pair_one_barrier
[simgrid.git] / src / smpi / colls / alltoallv-pair-one-barrier.c
1 #include "colls_private.h"
2 /*****************************************************************************
3
4  * Function: alltoall_pair
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.
19
20  * Auther: Ahmad Faraj
21
22  ****************************************************************************/
23 int
24 smpi_coll_tuned_alltoallv_pair_one_barrier(void *send_buff, int *send_counts, int *send_disps,
25                                           MPI_Datatype send_type,
26                                           void *recv_buff,  int *recv_counts, int *recv_disps,                                                                                  MPI_Datatype recv_type, MPI_Comm comm)
27 {
28
29   MPI_Aint send_chunk, recv_chunk;
30   MPI_Status s;
31   int i, src, dst, rank, num_procs;
32   int tag = 1;
33
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   smpi_mpi_barrier(comm);
43   for (i = 0; i < num_procs; i++) {
44     src = dst = rank ^ i;
45     smpi_mpi_sendrecv(send_ptr + send_disps[dst] * send_chunk, send_counts[dst], send_type, dst,
46                  tag, recv_ptr + recv_disps[src] * recv_chunk, recv_counts[src], recv_type,
47                  src, tag, comm, &s);
48   }
49
50   return MPI_SUCCESS;
51 }