Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Begin to add a MVAPICH2 collectives selector. Alltoall, Allgather and gather done.
[simgrid.git] / src / smpi / colls / alltoall-pair.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
11  * Function: alltoall_pair
12
13  * Return: int
14
15  * Inputs:
16     send_buff: send input buffer
17     send_count: number of elements to send
18     send_type: data type of elements being sent
19     recv_buff: receive output buffer
20     recv_count: number of elements to received
21     recv_type: data type of elements being received
22     comm: communicator
23
24  * Descrp: Function works when P is power of two. In each phase of P - 1
25            phases, nodes in pair communicate their data.
26
27  * Auther: Ahmad Faraj
28
29  ****************************************************************************/
30 /*
31 int alltoall_pair(void *send_buff, int send_count, MPI_Datatype send_type,
32                   void *recv_buff, int recv_count, MPI_Datatype recv_type,
33                   MPI_Comm comm)
34 {
35
36   MPI_Aint send_chunk, recv_chunk;
37   MPI_Status s;
38   MPI_Win win;
39   int assert = 0;
40   int i, src, dst, rank, num_procs;
41   int tag = 1, success = 1, failure = 0, pof2 = 1;
42
43   char *send_ptr = (char *) send_buff;
44   char *recv_ptr = (char *) recv_buff;
45
46   rank = smpi_comm_rank(comm);
47   num_procs = smpi_comm_size(comm);
48   send_chunk = smpi_datatype_get_extent(send_type);
49   recv_chunk = smpi_datatype_get_extent(recv_type);
50
51   MPI_Win_create(recv_buff, num_procs * recv_chunk * send_count, recv_chunk, 0,
52                  comm, &win);
53   send_chunk *= send_count;
54   recv_chunk *= recv_count;
55
56   MPI_Win_fence(assert, win);
57   for (i = 0; i < num_procs; i++) {
58     src = dst = rank ^ i;
59     MPI_Put(send_ptr + dst * send_chunk, send_count, send_type, dst,
60             rank * send_chunk, send_count, send_type, win);
61   }
62   MPI_Win_fence(assert, win);
63   MPI_Win_free(&win);
64   return 0;
65 }
66 */
67
68 int smpi_coll_tuned_alltoall_pair(void *send_buff, int send_count,
69                                   MPI_Datatype send_type,
70                                   void *recv_buff, int recv_count,
71                                   MPI_Datatype recv_type, MPI_Comm comm)
72 {
73
74   MPI_Aint send_chunk, recv_chunk;
75   MPI_Status s;
76   int i, src, dst, rank, num_procs;
77   int tag = COLL_TAG_ALLTOALL;
78   char *send_ptr = (char *) send_buff;
79   char *recv_ptr = (char *) recv_buff;
80
81   rank = smpi_comm_rank(comm);
82   num_procs = smpi_comm_size(comm);
83
84   if((num_procs&(num_procs-1)))
85     THROWF(arg_error,0, "alltoall pair algorithm can't be used with non power of two number of processes ! ");
86
87   send_chunk = smpi_datatype_get_extent(send_type);
88   recv_chunk = smpi_datatype_get_extent(recv_type);
89
90   send_chunk *= send_count;
91   recv_chunk *= recv_count;
92
93   for (i = 0; i < num_procs; i++) {
94     src = dst = rank ^ i;
95     smpi_mpi_sendrecv(send_ptr + dst * send_chunk, send_count, send_type, dst, tag,
96                  recv_ptr + src * recv_chunk, recv_count, recv_type, src, tag,
97                  comm, &s);
98   }
99
100   return MPI_SUCCESS;
101 }