Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
remove warning
[simgrid.git] / src / smpi / colls / alltoall-rdb.c
1 #include "colls_private.h"
2
3 /*****************************************************************************
4
5  * Function: alltoall_rdb
6
7  * Return: int
8
9  * inputs: 
10     send_buff: send input buffer
11     send_count: number of elements to send
12     send_type: data type of elements being sent
13     recv_buff: receive output buffer
14     recv_count: number of elements to received
15     recv_type: data type of elements being received
16     comm: communicator
17
18  * Descrp: Function realizes the allgather operation using the recursive
19            doubling algorithm.
20
21  * Auther: MPICH / slightly modified by Ahmad Faraj.  
22
23  ****************************************************************************/
24 int smpi_coll_tuned_alltoall_rdb(void *send_buff, int send_count,
25                                  MPI_Datatype send_type,
26                                  void *recv_buff, int recv_count,
27                                  MPI_Datatype recv_type, MPI_Comm comm)
28 {
29   /* MPI variables */
30   MPI_Status status;
31   MPI_Aint send_increment, recv_increment, extent;
32
33   int dst_tree_root, rank_tree_root, send_offset, recv_offset;
34   int rank, num_procs, j, k, dst, curr_size, max_size;
35   int last_recv_count = 0, tmp_mask, tree_root, num_procs_completed;
36   int tag = 1, mask = 1, i = 0;
37
38   char *tmp_buff;
39   char *send_ptr = (char *) send_buff;
40   char *recv_ptr = (char *) recv_buff;
41
42   num_procs = smpi_comm_size(comm);
43   rank = smpi_comm_rank(comm);
44   send_increment = smpi_datatype_get_extent(send_type);
45   recv_increment = smpi_datatype_get_extent(recv_type);
46   extent = smpi_datatype_get_extent(recv_type);
47
48   send_increment *= (send_count * num_procs);
49   recv_increment *= (recv_count * num_procs);
50
51   max_size = num_procs * recv_increment;
52
53   tmp_buff = (char *) xbt_malloc(max_size);
54
55   curr_size = send_count * num_procs;
56
57   smpi_mpi_sendrecv(send_ptr, curr_size, send_type, rank, tag,
58                tmp_buff + (rank * recv_increment),
59                curr_size, recv_type, rank, tag, comm, &status);
60
61   while (mask < num_procs) {
62     dst = rank ^ mask;
63     dst_tree_root = dst >> i;
64     dst_tree_root <<= i;
65     rank_tree_root = rank >> i;
66     rank_tree_root <<= i;
67     send_offset = rank_tree_root * send_increment;
68     recv_offset = dst_tree_root * recv_increment;
69
70     if (dst < num_procs) {
71       smpi_mpi_sendrecv(tmp_buff + send_offset, curr_size, send_type, dst, tag,
72                    tmp_buff + recv_offset, mask * recv_count * num_procs,
73                    recv_type, dst, tag, comm, &status);
74
75       last_recv_count = smpi_mpi_get_count(&status, recv_type);
76       curr_size += last_recv_count;
77     }
78
79
80     if (dst_tree_root + mask > num_procs) {
81
82       num_procs_completed = num_procs - rank_tree_root - mask;
83       /* num_procs_completed is the number of processes in this
84          subtree that have all the data. Send data to others
85          in a tree fashion. First find root of current tree
86          that is being divided into two. k is the number of
87          least-significant bits in this process's rank that
88          must be zeroed out to find the rank of the root */
89
90       j = mask;
91       k = 0;
92       while (j) {
93         j >>= 1;
94         k++;
95       }
96       k--;
97
98       tmp_mask = mask >> 1;
99
100       while (tmp_mask) {
101         dst = rank ^ tmp_mask;
102
103         tree_root = rank >> k;
104         tree_root <<= k;
105
106         /* send only if this proc has data and destination
107            doesn't have data. at any step, multiple processes
108            can send if they have the data */
109
110         if ((dst > rank)
111             && (rank < tree_root + num_procs_completed)
112             && (dst >= tree_root + num_procs_completed)) {
113           smpi_mpi_send(tmp_buff + dst_tree_root * send_increment,
114                    last_recv_count, send_type, dst, tag, comm);
115
116         }
117
118         /* recv only if this proc. doesn't have data and sender
119            has data */
120
121         else if ((dst < rank)
122                  && (dst < tree_root + num_procs_completed)
123                  && (rank >= tree_root + num_procs_completed)) {
124           smpi_mpi_recv(tmp_buff + dst_tree_root * send_increment,
125                    mask * num_procs * send_count, send_type, dst,
126                    tag, comm, &status);
127
128           last_recv_count = smpi_mpi_get_count(&status, send_type);
129           curr_size += last_recv_count;
130         }
131
132         tmp_mask >>= 1;
133         k--;
134       }
135     }
136
137     mask <<= 1;
138     i++;
139   }
140
141   for (i = 0; i < num_procs; i++)
142     smpi_mpi_sendrecv(tmp_buff + (rank + i * num_procs) * send_count * extent,
143                  send_count, send_type, rank, tag,
144                  recv_ptr + (i * recv_count * extent),
145                  recv_count, recv_type, rank, tag, comm, &status);
146   free(tmp_buff);
147   return MPI_SUCCESS;
148 }