Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
adapt two collectives of starmpi to avoid timing issues, by using only smpi calls...
[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, 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, tmp_mask, tree_root, num_procs_completed;
36   int tag = 1, mask = 1, success = 1, failure = 0, i = 0;
37
38   char *tmp_buff;
39   char *send_ptr = (char *) send_buff;
40   char *recv_ptr = (char *) recv_buff;
41
42   MPI_Comm_size(comm, &num_procs);
43   MPI_Comm_rank(comm, &rank);
44   MPI_Type_extent(send_type, &send_increment);
45   MPI_Type_extent(recv_type, &recv_increment);
46   MPI_Type_extent(recv_type, &extent);
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 *) malloc(max_size);
54   if (!tmp_buff) {
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     dst = rank ^ mask;
68     dst_tree_root = dst >> i;
69     dst_tree_root <<= i;
70     rank_tree_root = rank >> i;
71     rank_tree_root <<= i;
72     send_offset = rank_tree_root * send_increment;
73     recv_offset = dst_tree_root * recv_increment;
74
75     if (dst < num_procs) {
76       MPI_Sendrecv(tmp_buff + send_offset, curr_size, send_type, dst, tag,
77                    tmp_buff + recv_offset, mask * recv_count * num_procs,
78                    recv_type, dst, tag, comm, &status);
79
80       MPI_Get_count(&status, recv_type, &last_recv_count);
81       curr_size += last_recv_count;
82     }
83
84
85     if (dst_tree_root + mask > num_procs) {
86
87       num_procs_completed = num_procs - rank_tree_root - mask;
88       /* num_procs_completed is the number of processes in this
89          subtree that have all the data. Send data to others
90          in a tree fashion. First find root of current tree
91          that is being divided into two. k is the number of
92          least-significant bits in this process's rank that
93          must be zeroed out to find the rank of the root */
94
95       j = mask;
96       k = 0;
97       while (j) {
98         j >>= 1;
99         k++;
100       }
101       k--;
102
103       tmp_mask = mask >> 1;
104
105       while (tmp_mask) {
106         dst = rank ^ tmp_mask;
107
108         tree_root = rank >> k;
109         tree_root <<= k;
110
111         /* send only if this proc has data and destination
112            doesn't have data. at any step, multiple processes
113            can send if they have the data */
114
115         if ((dst > rank)
116             && (rank < tree_root + num_procs_completed)
117             && (dst >= tree_root + num_procs_completed)) {
118           MPI_Send(tmp_buff + dst_tree_root * send_increment,
119                    last_recv_count, send_type, dst, tag, comm);
120
121         }
122
123         /* recv only if this proc. doesn't have data and sender
124            has data */
125
126         else if ((dst < rank)
127                  && (dst < tree_root + num_procs_completed)
128                  && (rank >= tree_root + num_procs_completed)) {
129           MPI_Recv(tmp_buff + dst_tree_root * send_increment,
130                    mask * num_procs * send_count, send_type, dst,
131                    tag, comm, &status);
132
133           MPI_Get_count(&status, send_type, &last_recv_count);
134           curr_size += last_recv_count;
135         }
136
137         tmp_mask >>= 1;
138         k--;
139       }
140     }
141
142     mask <<= 1;
143     i++;
144   }
145
146   for (i = 0; i < num_procs; i++)
147     MPI_Sendrecv(tmp_buff + (rank + i * num_procs) * send_count * extent,
148                  send_count, send_type, rank, tag,
149                  recv_ptr + (i * recv_count * extent),
150                  recv_count, recv_type, rank, tag, comm, &status);
151   free(tmp_buff);
152   return success;
153 }