Logo AND Algorithmique Numérique Distribuée

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