Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
kill all trailling whitespaces
[simgrid.git] / src / smpi / colls / alltoall / alltoall-rdb.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 "src/smpi/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(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 *tmp_buff;
48   char *send_ptr = (char *) send_buff;
49   char *recv_ptr = (char *) recv_buff;
50
51   num_procs = comm->size();
52   rank = comm->rank();
53   send_increment = send_type->get_extent();
54   recv_increment = recv_type->get_extent();
55   extent = recv_type->get_extent();
56
57   send_increment *= (send_count * num_procs);
58   recv_increment *= (recv_count * num_procs);
59
60   max_size = num_procs * recv_increment;
61
62   tmp_buff = (char *) smpi_get_tmp_sendbuffer(max_size);
63
64   curr_size = send_count * num_procs;
65
66   Request::sendrecv(send_ptr, curr_size, send_type, rank, tag,
67                tmp_buff + (rank * recv_increment),
68                curr_size, recv_type, rank, tag, comm, &status);
69
70   while (mask < num_procs) {
71     dst = rank ^ mask;
72     dst_tree_root = dst >> i;
73     dst_tree_root <<= i;
74     rank_tree_root = rank >> i;
75     rank_tree_root <<= i;
76     send_offset = rank_tree_root * send_increment;
77     recv_offset = dst_tree_root * recv_increment;
78
79     if (dst < num_procs) {
80       Request::sendrecv(tmp_buff + send_offset, curr_size, send_type, dst, tag,
81                    tmp_buff + recv_offset, mask * recv_count * num_procs,
82                    recv_type, dst, tag, comm, &status);
83
84       last_recv_count = Status::get_count(&status, recv_type);
85       curr_size += last_recv_count;
86     }
87
88
89     if (dst_tree_root + mask > num_procs) {
90
91       num_procs_completed = num_procs - rank_tree_root - mask;
92       /* num_procs_completed is the number of processes in this
93          subtree that have all the data. Send data to others
94          in a tree fashion. First find root of current tree
95          that is being divided into two. k is the number of
96          least-significant bits in this process's rank that
97          must be zeroed out to find the rank of the root */
98
99       j = mask;
100       k = 0;
101       while (j) {
102         j >>= 1;
103         k++;
104       }
105       k--;
106
107       tmp_mask = mask >> 1;
108
109       while (tmp_mask) {
110         dst = rank ^ tmp_mask;
111
112         tree_root = rank >> k;
113         tree_root <<= k;
114
115         /* send only if this proc has data and destination
116            doesn't have data. at any step, multiple processes
117            can send if they have the data */
118
119         if ((dst > rank)
120             && (rank < tree_root + num_procs_completed)
121             && (dst >= tree_root + num_procs_completed)) {
122           Request::send(tmp_buff + dst_tree_root * send_increment,
123                    last_recv_count, send_type, dst, tag, comm);
124
125         }
126
127         /* recv only if this proc. doesn't have data and sender
128            has data */
129
130         else if ((dst < rank)
131                  && (dst < tree_root + num_procs_completed)
132                  && (rank >= tree_root + num_procs_completed)) {
133           Request::recv(tmp_buff + dst_tree_root * send_increment,
134                    mask * num_procs * send_count, send_type, dst,
135                    tag, comm, &status);
136
137           last_recv_count = Status::get_count(&status, send_type);
138           curr_size += last_recv_count;
139         }
140
141         tmp_mask >>= 1;
142         k--;
143       }
144     }
145
146     mask <<= 1;
147     i++;
148   }
149
150   for (i = 0; i < num_procs; i++)
151     Request::sendrecv(tmp_buff + (rank + i * num_procs) * send_count * extent,
152                  send_count, send_type, rank, tag,
153                  recv_ptr + (i * recv_count * extent),
154                  recv_count, recv_type, rank, tag, comm, &status);
155   smpi_free_tmp_buffer(tmp_buff);
156   return MPI_SUCCESS;
157 }
158 }
159 }