Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
57069ef493cdf5104b10dac86b3c7595a8e6fe40
[simgrid.git] / src / smpi / colls / alltoall-pair.cpp
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 smpi_coll_tuned_alltoall_pair_rma(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_Win win;
38   int assert = 0;
39   int i, dst, rank, num_procs;
40
41   char *send_ptr = (char *) send_buff;
42
43   rank = comm->rank();
44   num_procs = comm->size();
45   send_chunk = smpi_datatype_get_extent(send_type);
46   recv_chunk = smpi_datatype_get_extent(recv_type);
47
48   win=new simgrid::smpi::Win(recv_buff, num_procs * recv_chunk * send_count, recv_chunk, 0,
49                  comm);
50   send_chunk *= send_count;
51   recv_chunk *= recv_count;
52
53   win->fence(assert);
54   for (i = 0; i < num_procs; i++) {
55     dst = rank ^ i;
56     win->put(send_ptr + dst * send_chunk, send_count, send_type, dst,
57             rank /* send_chunk*/, send_count, send_type);
58   }
59   win->fence(assert);
60   delete win;
61   return 0;
62 }
63
64
65 int smpi_coll_tuned_alltoall_pair(void *send_buff, int send_count,
66                                   MPI_Datatype send_type,
67                                   void *recv_buff, int recv_count,
68                                   MPI_Datatype recv_type, MPI_Comm comm)
69 {
70
71   MPI_Aint send_chunk, recv_chunk;
72   MPI_Status s;
73   int i, src, dst, rank, num_procs;
74   int tag = COLL_TAG_ALLTOALL;
75   char *send_ptr = (char *) send_buff;
76   char *recv_ptr = (char *) recv_buff;
77
78   rank = comm->rank();
79   num_procs = comm->size();
80
81   if((num_procs&(num_procs-1)))
82     THROWF(arg_error,0, "alltoall pair algorithm can't be used with non power of two number of processes ! ");
83
84   send_chunk = smpi_datatype_get_extent(send_type);
85   recv_chunk = smpi_datatype_get_extent(recv_type);
86
87   send_chunk *= send_count;
88   recv_chunk *= recv_count;
89
90   for (i = 0; i < num_procs; i++) {
91     src = dst = rank ^ i;
92     Request::sendrecv(send_ptr + dst * send_chunk, send_count, send_type, dst, tag,
93                  recv_ptr + src * recv_chunk, recv_count, recv_type, src, tag,
94                  comm, &s);
95   }
96
97   return MPI_SUCCESS;
98 }