Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
5abbbc9e8061c699cfe426fe149570239df09e13
[simgrid.git] / src / smpi / colls / alltoall-rdb.c
1 #include "colls.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,
28                                  MPI_Comm comm)
29 {
30   /* MPI variables */
31   MPI_Status status;
32   MPI_Aint send_increment, recv_increment, extent;
33
34   int dst_tree_root, rank_tree_root, send_offset, recv_offset;
35   int rank, num_procs, j, k, dst, curr_size, max_size;
36   int last_recv_count, tmp_mask, tree_root, num_procs_completed;
37   int tag = 1, mask = 1, success = 1, failure = 0, i = 0;
38
39   char *tmp_buff;
40   char *send_ptr = (char *) send_buff;
41   char *recv_ptr = (char *) recv_buff;
42
43   MPI_Comm_size(comm, &num_procs);
44   MPI_Comm_rank(comm, &rank);
45   MPI_Type_extent(send_type, &send_increment);
46   MPI_Type_extent(recv_type, &recv_increment);
47   MPI_Type_extent(recv_type, &extent);
48
49   send_increment *= (send_count * num_procs);
50   recv_increment *= (recv_count * num_procs);
51
52   max_size = num_procs * recv_increment;
53
54   tmp_buff = (char *) malloc(max_size);
55   if (!tmp_buff) {
56     printf("alltoall-rdb:56: cannot allocate memory\n");
57     MPI_Finalize();
58     exit(failure);
59   }
60
61   curr_size = send_count * num_procs;
62
63   MPI_Sendrecv(send_ptr, curr_size, send_type, rank, tag,
64                tmp_buff + (rank * recv_increment),
65                curr_size, recv_type, rank, tag, comm, &status);
66
67   while (mask < num_procs) {
68     dst = rank ^ mask;
69     dst_tree_root = dst >> i;
70     dst_tree_root <<= i;
71     rank_tree_root = rank >> i;
72     rank_tree_root <<= i;
73     send_offset = rank_tree_root * send_increment;
74     recv_offset = dst_tree_root * recv_increment;
75
76     if (dst < num_procs) {
77       MPI_Sendrecv(tmp_buff + send_offset, curr_size, send_type, dst, tag,
78                    tmp_buff + recv_offset, mask * recv_count * num_procs,
79                    recv_type, dst, tag, comm, &status);
80
81       MPI_Get_count(&status, recv_type, &last_recv_count);
82       curr_size += last_recv_count;
83     }
84
85
86     if (dst_tree_root + mask > num_procs) {
87
88       num_procs_completed = num_procs - rank_tree_root - mask;
89       /* num_procs_completed is the number of processes in this
90          subtree that have all the data. Send data to others
91          in a tree fashion. First find root of current tree
92          that is being divided into two. k is the number of
93          least-significant bits in this process's rank that
94          must be zeroed out to find the rank of the root */
95
96       j = mask;
97       k = 0;
98       while (j) {
99         j >>= 1;
100         k++;
101       }
102       k--;
103
104       tmp_mask = mask >> 1;
105
106       while (tmp_mask) {
107         dst = rank ^ tmp_mask;
108
109         tree_root = rank >> k;
110         tree_root <<= k;
111
112         /* send only if this proc has data and destination
113            doesn't have data. at any step, multiple processes
114            can send if they have the data */
115
116         if ((dst > rank)
117             && (rank < tree_root + num_procs_completed)
118             && (dst >= tree_root + num_procs_completed)) {
119           MPI_Send(tmp_buff + dst_tree_root * send_increment,
120                    last_recv_count, send_type, dst, tag, comm);
121
122         }
123
124         /* recv only if this proc. doesn't have data and sender
125            has data */
126
127         else if ((dst < rank)
128                  && (dst < tree_root + num_procs_completed)
129                  && (rank >= tree_root + num_procs_completed)) {
130           MPI_Recv(tmp_buff + dst_tree_root * send_increment,
131                    mask * num_procs * send_count, send_type, dst,
132                    tag, comm, &status);
133
134           MPI_Get_count(&status, send_type, &last_recv_count);
135           curr_size += last_recv_count;
136         }
137
138         tmp_mask >>= 1;
139         k--;
140       }
141     }
142
143     mask <<= 1;
144     i++;
145   }
146
147   for (i = 0; i < num_procs; i++)
148     MPI_Sendrecv(tmp_buff + (rank + i * num_procs) * send_count * extent,
149                  send_count, send_type, rank, tag,
150                  recv_ptr + (i * recv_count * extent),
151                  recv_count, recv_type, rank, tag, comm, &status);
152   free(tmp_buff);
153   return success;
154 }