Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
[sonar] Use unsigned char* for smpi buffers.
[simgrid.git] / src / smpi / colls / alltoall / alltoall-rdb.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 #include "../colls_private.hpp"
8 #include "smpi_status.hpp"
9
10 /*****************************************************************************
11
12  * Function: alltoall_rdb
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 realizes the allgather operation using the recursive
26            doubling algorithm.
27
28  * Auther: MPICH / slightly modified by Ahmad Faraj.
29
30  ****************************************************************************/
31 namespace simgrid{
32 namespace smpi{
33 int Coll_alltoall_rdb::alltoall(const void *send_buff, int send_count,
34                                  MPI_Datatype send_type,
35                                  void *recv_buff, int recv_count,
36                                  MPI_Datatype recv_type, MPI_Comm comm)
37 {
38   /* MPI variables */
39   MPI_Status status;
40   MPI_Aint send_increment, recv_increment, extent;
41
42   int dst_tree_root, rank_tree_root, send_offset, recv_offset;
43   int rank, num_procs, j, k, dst, curr_size, max_size;
44   int last_recv_count = 0, tmp_mask, tree_root, num_procs_completed;
45   int tag = COLL_TAG_ALLTOALL, mask = 1, i = 0;
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   send_increment = send_type->get_extent();
53   recv_increment = recv_type->get_extent();
54   extent = recv_type->get_extent();
55
56   send_increment *= (send_count * num_procs);
57   recv_increment *= (recv_count * num_procs);
58
59   max_size = num_procs * recv_increment;
60
61   unsigned char* tmp_buff = smpi_get_tmp_sendbuffer(max_size);
62
63   curr_size = send_count * num_procs;
64
65   Request::sendrecv(send_ptr, curr_size, send_type, rank, tag,
66                tmp_buff + (rank * recv_increment),
67                curr_size, recv_type, rank, tag, comm, &status);
68
69   while (mask < num_procs) {
70     dst = rank ^ mask;
71     dst_tree_root = dst >> i;
72     dst_tree_root <<= i;
73     rank_tree_root = rank >> i;
74     rank_tree_root <<= i;
75     send_offset = rank_tree_root * send_increment;
76     recv_offset = dst_tree_root * recv_increment;
77
78     if (dst < num_procs) {
79       Request::sendrecv(tmp_buff + send_offset, curr_size, send_type, dst, tag,
80                    tmp_buff + recv_offset, mask * recv_count * num_procs,
81                    recv_type, dst, tag, comm, &status);
82
83       last_recv_count = Status::get_count(&status, recv_type);
84       curr_size += last_recv_count;
85     }
86
87
88     if (dst_tree_root + mask > num_procs) {
89
90       num_procs_completed = num_procs - rank_tree_root - mask;
91       /* num_procs_completed is the number of processes in this
92          subtree that have all the data. Send data to others
93          in a tree fashion. First find root of current tree
94          that is being divided into two. k is the number of
95          least-significant bits in this process's rank that
96          must be zeroed out to find the rank of the root */
97
98       j = mask;
99       k = 0;
100       while (j) {
101         j >>= 1;
102         k++;
103       }
104       k--;
105
106       tmp_mask = mask >> 1;
107
108       while (tmp_mask) {
109         dst = rank ^ tmp_mask;
110
111         tree_root = rank >> k;
112         tree_root <<= k;
113
114         /* send only if this proc has data and destination
115            doesn't have data. at any step, multiple processes
116            can send if they have the data */
117
118         if ((dst > rank)
119             && (rank < tree_root + num_procs_completed)
120             && (dst >= tree_root + num_procs_completed)) {
121           Request::send(tmp_buff + dst_tree_root * send_increment,
122                    last_recv_count, send_type, dst, tag, comm);
123
124         }
125
126         /* recv only if this proc. doesn't have data and sender
127            has data */
128
129         else if ((dst < rank)
130                  && (dst < tree_root + num_procs_completed)
131                  && (rank >= tree_root + num_procs_completed)) {
132           Request::recv(tmp_buff + dst_tree_root * send_increment,
133                    mask * num_procs * send_count, send_type, dst,
134                    tag, comm, &status);
135
136           last_recv_count = Status::get_count(&status, send_type);
137           curr_size += last_recv_count;
138         }
139
140         tmp_mask >>= 1;
141         k--;
142       }
143     }
144
145     mask <<= 1;
146     i++;
147   }
148
149   for (i = 0; i < num_procs; i++)
150     Request::sendrecv(tmp_buff + (rank + i * num_procs) * send_count * extent,
151                  send_count, send_type, rank, tag,
152                  recv_ptr + (i * recv_count * extent),
153                  recv_count, recv_type, rank, tag, comm, &status);
154   smpi_free_tmp_buffer(tmp_buff);
155   return MPI_SUCCESS;
156 }
157 }
158 }