Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
b5bd8a52fb3ddf471c65c4515170c0e1d9032a34
[simgrid.git] / src / smpi / colls / alltoall-rdb.c
1 #include "colls.h"
2 /*****************************************************************************
3
4  * Function: alltoall_rdb
5
6  * Return: int
7
8  * inputs: 
9     send_buff: send input buffer
10     send_count: number of elements to send
11     send_type: data type of elements being sent
12     recv_buff: receive output buffer
13     recv_count: number of elements to received
14     recv_type: data type of elements being received
15     comm: communicator
16
17  * Descrp: Function realizes the allgather operation using the recursive
18            doubling algorithm.
19
20  * Auther: MPICH / slightly modified by Ahmad Faraj.  
21
22  ****************************************************************************/
23 int
24 smpi_coll_tuned_alltoall_rdb(void * send_buff, int send_count, MPI_Datatype send_type,
25              void * recv_buff, int recv_count, MPI_Datatype recv_type,
26              MPI_Comm comm)
27 {
28   /* MPI variables */
29   MPI_Status status;
30   MPI_Aint send_increment, recv_increment, extent;
31
32   int dst_tree_root, rank_tree_root, send_offset, recv_offset;
33   int rank, num_procs, j, k, dst, curr_size, max_size;
34   int last_recv_count, tmp_mask, tree_root, num_procs_completed;
35   int tag = 1, mask = 1, success = 1, failure = 0, c = 0, i = 0;
36
37   char * tmp_buff;
38   char * send_ptr = (char *) send_buff;
39   char * recv_ptr = (char *) recv_buff;
40
41   MPI_Comm_size(comm, &num_procs);
42   MPI_Comm_rank(comm, &rank);
43   MPI_Type_extent(send_type, &send_increment);
44   MPI_Type_extent(recv_type, &recv_increment);
45   MPI_Type_extent(recv_type, &extent);
46
47   send_increment *= (send_count * num_procs);
48   recv_increment *= (recv_count * num_procs);
49    
50   max_size = num_procs * recv_increment;
51  
52   tmp_buff = (char * ) malloc(max_size);
53   if (!tmp_buff)
54     {
55       printf("alltoall-rdb:56: cannot allocate memory\n");
56       MPI_Finalize();
57       exit(failure);
58     }
59
60   curr_size = send_count * num_procs;
61
62   MPI_Sendrecv(send_ptr, curr_size, send_type, rank, tag,
63                tmp_buff + (rank * recv_increment),
64                curr_size, recv_type, rank, tag, comm, &status);
65
66   while (mask < num_procs)
67     {
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         {
78           MPI_Sendrecv(tmp_buff + send_offset, curr_size, send_type, dst, tag,
79                         tmp_buff + recv_offset, mask * recv_count * num_procs,
80                         recv_type, dst, tag, comm, &status);
81
82           MPI_Get_count(&status, recv_type, &last_recv_count);
83           curr_size += last_recv_count;
84         }
85      
86      
87       if (dst_tree_root + mask > num_procs)
88         {
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             {
102               j >>= 1;
103               k++;
104             }
105           k--;
106
107           tmp_mask = mask >> 1;
108                   
109           while (tmp_mask)
110             {
111               dst = rank ^ tmp_mask;
112                       
113               tree_root = rank >> k;
114               tree_root <<= k;
115                       
116               /* send only if this proc has data and destination
117                  doesn't have data. at any step, multiple processes
118                  can send if they have the data */
119
120               if ((dst > rank)
121                   && (rank < tree_root + num_procs_completed)
122                   && (dst >= tree_root + num_procs_completed))
123                 {
124                   MPI_Send(tmp_buff + dst_tree_root * send_increment,
125                             last_recv_count, send_type, dst, tag, comm);
126
127                 }
128            
129               /* recv only if this proc. doesn't have data and sender
130                  has data */
131                       
132               else if ((dst < rank)
133                        && (dst < tree_root + num_procs_completed)
134                        && (rank >= tree_root + num_procs_completed))
135                 {
136                   MPI_Recv(tmp_buff + dst_tree_root * send_increment,
137                             mask * num_procs * send_count, send_type, dst,
138                             tag, comm, &status);
139
140                   MPI_Get_count(&status, send_type, &last_recv_count);
141                   curr_size += last_recv_count;
142                 }
143
144               tmp_mask >>= 1;
145               k--;
146             }
147         }
148
149       mask <<= 1;
150       i++;
151     }
152
153   for (i = 0; i < num_procs; i++)
154      MPI_Sendrecv(tmp_buff + (rank + i * num_procs) * send_count * extent,
155                    send_count, send_type, rank, tag,
156                    recv_ptr + (i * recv_count * extent),
157                    recv_count, recv_type, rank, tag, comm, &status);
158   free(tmp_buff);
159   return success;
160 }