Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Add/update copyright notices.
[simgrid.git] / src / smpi / colls / alltoall-pair-one-barrier.c
1 /* Copyright (c) 2013-2014. The SimGrid Team.
2  * All rights reserved.                                                     */
3
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. */
6
7 #include "colls_private.h"
8 /*****************************************************************************
9
10  * Function: alltoall_pair
11
12  * Return: int
13
14  * Inputs:
15     send_buff: send input buffer
16     send_count: number of elements to send
17     send_type: data type of elements being sent
18     recv_buff: receive output buffer
19     recv_count: number of elements to received
20     recv_type: data type of elements being received
21     comm: communicator
22
23  * Descrp: Function works when P is power of two. In each phase of P - 1
24            phases, nodes in pair communicate their data.
25
26  * Auther: Ahmad Faraj
27
28  ****************************************************************************/
29 int
30 smpi_coll_tuned_alltoall_pair_one_barrier(void *send_buff, int send_count,
31                                           MPI_Datatype send_type,
32                                           void *recv_buff, int recv_count,
33                                           MPI_Datatype recv_type, MPI_Comm comm)
34 {
35
36   MPI_Aint send_chunk, recv_chunk;
37   MPI_Status s;
38   int i, src, dst, rank, num_procs;
39   int tag = COLL_TAG_ALLTOALL;
40
41   char *send_ptr = (char *) send_buff;
42   char *recv_ptr = (char *) recv_buff;
43
44   rank = smpi_comm_rank(comm);
45   num_procs = smpi_comm_size(comm);
46
47   if((num_procs&(num_procs-1)))
48     THROWF(arg_error,0, "alltoall pair algorithm can't be used with non power of two number of processes ! ");
49
50   send_chunk = smpi_datatype_get_extent(send_type);
51   recv_chunk = smpi_datatype_get_extent(recv_type);
52
53   send_chunk *= send_count;
54   recv_chunk *= recv_count;
55
56   mpi_coll_barrier_fun(comm);
57   for (i = 0; i < num_procs; i++) {
58     src = dst = rank ^ i;
59     smpi_mpi_sendrecv(send_ptr + dst * send_chunk, send_count, send_type, dst,
60                  tag, recv_ptr + src * recv_chunk, recv_count, recv_type,
61                  src, tag, comm, &s);
62   }
63
64   return MPI_SUCCESS;
65 }