Logo AND Algorithmique Numérique Distribuée

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