Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Spell check.
[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  * Author: 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   char *send_ptr = (char *) send_buff;
48   char *recv_ptr = (char *) recv_buff;
49
50   num_procs = comm->size();
51   rank = comm->rank();
52
53   extent = recv_type->get_extent();
54
55   unsigned char* tmp_buff = smpi_get_tmp_sendbuffer(num_procs * recv_count * extent);
56   int* disps         = new int[num_procs];
57   int* blocks_length = new int[num_procs];
58
59   Request::sendrecv(send_ptr + rank * send_count * extent,
60                (num_procs - rank) * send_count, send_type, rank, tag,
61                recv_ptr, (num_procs - rank) * recv_count, recv_type, rank,
62                tag, comm, &status);
63
64   Request::sendrecv(send_ptr, rank * send_count, send_type, rank, tag,
65                recv_ptr + (num_procs - rank) * recv_count * extent,
66                rank * recv_count, recv_type, rank, tag, comm, &status);
67
68
69
70   MPI_Pack_size(send_count * num_procs, send_type, comm, &pack_size);
71
72   while (pof2 < num_procs) {
73     dst = (rank + pof2) % num_procs;
74     src = (rank - pof2 + num_procs) % num_procs;
75
76
77     count = 0;
78     for (block = 1; block < num_procs; block++)
79       if (block & pof2) {
80         blocks_length[count] = send_count;
81         disps[count] = block * send_count;
82         count++;
83       }
84
85     MPI_Type_indexed(count, blocks_length, disps, recv_type, &new_type);
86     new_type->commit();
87
88     position = 0;
89     MPI_Pack(recv_buff, 1, new_type, tmp_buff, pack_size, &position, comm);
90
91     Request::sendrecv(tmp_buff, position, MPI_PACKED, dst, tag, recv_buff, 1,
92                  new_type, src, tag, comm, &status);
93     Datatype::unref(new_type);
94
95     pof2 *= 2;
96   }
97
98   delete[] disps;
99   delete[] blocks_length;
100
101   Request::sendrecv(recv_ptr + (rank + 1) * recv_count * extent,
102                (num_procs - rank - 1) * recv_count, send_type,
103                rank, tag, tmp_buff, (num_procs - rank - 1) * recv_count,
104                recv_type, rank, tag, comm, &status);
105
106   Request::sendrecv(recv_ptr, (rank + 1) * recv_count, send_type, rank, tag,
107                tmp_buff + (num_procs - rank - 1) * recv_count * extent,
108                (rank + 1) * recv_count, recv_type, rank, tag, comm, &status);
109
110
111   for (i = 0; i < num_procs; i++)
112     Request::sendrecv(tmp_buff + i * recv_count * extent, recv_count, send_type,
113                  rank, tag,
114                  recv_ptr + (num_procs - i - 1) * recv_count * extent,
115                  recv_count, recv_type, rank, tag, comm, &status);
116
117   smpi_free_tmp_buffer(tmp_buff);
118   return MPI_SUCCESS;
119 }
120
121 }
122 }