Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
9c2f5bcd466ae0fbd014cbf312fbae5c7d460b62
[simgrid.git] / src / smpi / colls / alltoall / alltoall-bruck.cpp
1 /* Copyright (c) 2013-2019. 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 /*****************************************************************************
8
9  * Function: alltoall_bruck
10
11  * Return: int
12
13  * inputs:
14     send_buff: send input buffer
15     send_count: number of elements to send
16     send_type: data type of elements being sent
17     recv_buff: receive output buffer
18     recv_count: number of elements to received
19     recv_type: data type of elements being received
20     comm: communicator
21
22  * Descrp: Function realizes the alltoall operation using the bruck algorithm.
23
24  * Auther: MPICH / modified by Ahmad Faraj
25
26  ****************************************************************************/
27
28 #include "../colls_private.hpp"
29
30 namespace simgrid{
31 namespace smpi{
32
33
34 int
35 Coll_alltoall_bruck::alltoall(const void *send_buff, int send_count,
36                                MPI_Datatype send_type, void *recv_buff,
37                                int recv_count, MPI_Datatype recv_type,
38                                MPI_Comm comm)
39 {
40   MPI_Status status;
41   MPI_Aint extent;
42   MPI_Datatype new_type;
43
44   int i, src, dst, rank, num_procs, count, block, position;
45   int pack_size, tag = COLL_TAG_ALLTOALL, pof2 = 1;
46
47
48   char *tmp_buff;
49   char *send_ptr = (char *) send_buff;
50   char *recv_ptr = (char *) recv_buff;
51
52   num_procs = comm->size();
53   rank = comm->rank();
54
55   extent = recv_type->get_extent();
56
57   tmp_buff = (char *) smpi_get_tmp_sendbuffer(num_procs * recv_count * extent);
58   int* disps         = new int[num_procs];
59   int* blocks_length = new int[num_procs];
60
61   Request::sendrecv(send_ptr + rank * send_count * extent,
62                (num_procs - rank) * send_count, send_type, rank, tag,
63                recv_ptr, (num_procs - rank) * recv_count, recv_type, rank,
64                tag, comm, &status);
65
66   Request::sendrecv(send_ptr, rank * send_count, send_type, rank, tag,
67                recv_ptr + (num_procs - rank) * recv_count * extent,
68                rank * recv_count, recv_type, rank, tag, comm, &status);
69
70
71
72   MPI_Pack_size(send_count * num_procs, send_type, comm, &pack_size);
73
74   while (pof2 < num_procs) {
75     dst = (rank + pof2) % num_procs;
76     src = (rank - pof2 + num_procs) % num_procs;
77
78
79     count = 0;
80     for (block = 1; block < num_procs; block++)
81       if (block & pof2) {
82         blocks_length[count] = send_count;
83         disps[count] = block * send_count;
84         count++;
85       }
86
87     MPI_Type_indexed(count, blocks_length, disps, recv_type, &new_type);
88     new_type->commit();
89
90     position = 0;
91     MPI_Pack(recv_buff, 1, new_type, tmp_buff, pack_size, &position, comm);
92
93     Request::sendrecv(tmp_buff, position, MPI_PACKED, dst, tag, recv_buff, 1,
94                  new_type, src, tag, comm, &status);
95     Datatype::unref(new_type);
96
97     pof2 *= 2;
98   }
99
100   delete[] disps;
101   delete[] blocks_length;
102
103   Request::sendrecv(recv_ptr + (rank + 1) * recv_count * extent,
104                (num_procs - rank - 1) * recv_count, send_type,
105                rank, tag, tmp_buff, (num_procs - rank - 1) * recv_count,
106                recv_type, rank, tag, comm, &status);
107
108   Request::sendrecv(recv_ptr, (rank + 1) * recv_count, send_type, rank, tag,
109                tmp_buff + (num_procs - rank - 1) * recv_count * extent,
110                (rank + 1) * recv_count, recv_type, rank, tag, comm, &status);
111
112
113   for (i = 0; i < num_procs; i++)
114     Request::sendrecv(tmp_buff + i * recv_count * extent, recv_count, send_type,
115                  rank, tag,
116                  recv_ptr + (num_procs - i - 1) * recv_count * extent,
117                  recv_count, recv_type, rank, tag, comm, &status);
118
119   smpi_free_tmp_buffer(tmp_buff);
120   return MPI_SUCCESS;
121 }
122
123 }
124 }