Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Reindent files before changes.
[simgrid.git] / src / smpi / colls / alltoall-3dmesh.c
1 #include "colls.h"
2 #include <math.h>
3
4 /*****************************************************************************
5
6  * Function: alltoall_3dmesh_shoot
7
8  * Return: int
9
10  * Inputs:
11     send_buff: send input buffer
12     send_count: number of elements to send
13     send_type: data type of elements being sent
14     recv_buff: receive output buffer
15     recv_count: number of elements to received
16     recv_type: data type of elements being received
17     comm: communicator
18
19  * Descrp: Function realizes the alltoall operation using the 3dmesh
20            algorithm. It actually performs allgather operation in x dimension,
21            y dimension, then in z dimension. Each node then extracts the
22            needed data. The communication in all dimension is simple.
23
24  * Auther: Ahmad Faraj
25 ****************************************************************************/
26
27 int alltoall_check_is_3dmesh(int num, int *i, int *j, int *k)
28 {
29   int x, max = num / 3;
30   x = cbrt(num);
31   *i = *j = *k = 0;
32   while (x <= max) {
33     if ((num % (x * x)) == 0) {
34       *i = *j = x;
35       *k = num / (x * x);
36       return 1;
37     }
38     x++;
39   }
40   return 0;
41 }
42
43 int smpi_coll_tuned_alltoall_3dmesh(void *send_buff, int send_count,
44                                     MPI_Datatype send_type,
45                                     void *recv_buff, int recv_count,
46                                     MPI_Datatype recv_type,
47                                     MPI_Comm comm)
48 {
49   MPI_Request *reqs, *req_ptr;
50   MPI_Aint extent;
51   MPI_Status status, *statuses;
52   int i, j, src, dst, rank, num_procs, num_reqs, X, Y, Z, block_size, count;
53   int my_z, two_dsize, my_row_base, my_col_base, my_z_base, src_row_base;
54   int src_z_base, send_offset, recv_offset, tag = 1, failure = 0, success = 1;
55
56   char *tmp_buff1, *tmp_buff2;
57
58   MPI_Comm_rank(comm, &rank);
59   MPI_Comm_size(comm, &num_procs);
60   MPI_Type_extent(send_type, &extent);
61
62   if (!alltoall_check_is_3dmesh(num_procs, &X, &Y, &Z))
63     return failure;
64
65   num_reqs = X;
66   if (Y > X)
67     num_reqs = Y;
68   if (Z > Y)
69     num_reqs = Z;
70
71   two_dsize = X * Y;
72   my_z = rank / two_dsize;
73
74   my_row_base = (rank / X) * X;
75   my_col_base = (rank % Y) + (my_z * two_dsize);
76   my_z_base = my_z * two_dsize;
77
78   block_size = extent * send_count;
79
80   tmp_buff1 = (char *) malloc(block_size * num_procs * two_dsize);
81   if (!tmp_buff1) {
82     printf("alltoall-3Dmesh:97: cannot allocate memory\n");
83     MPI_Finalize();
84     exit(failure);
85   }
86
87   tmp_buff2 = (char *) malloc(block_size * two_dsize);
88   if (!tmp_buff2) {
89     printf("alltoall-3Dmesh:105: cannot allocate memory\n");
90     MPI_Finalize();
91     exit(failure);
92   }
93
94   statuses = (MPI_Status *) malloc(num_reqs * sizeof(MPI_Status));
95   reqs = (MPI_Request *) malloc(num_reqs * sizeof(MPI_Request));
96   if (!reqs) {
97     printf("alltoall-3Dmesh:113: cannot allocate memory\n");
98     MPI_Finalize();
99     exit(failure);
100   }
101
102   req_ptr = reqs;
103
104   send_offset = recv_offset = (rank % two_dsize) * block_size * num_procs;
105
106   MPI_Sendrecv(send_buff, send_count * num_procs, send_type, rank, tag,
107                tmp_buff1 + recv_offset, num_procs * recv_count,
108                recv_type, rank, tag, comm, &status);
109
110   count = send_count * num_procs;
111
112   for (i = 0; i < Y; i++) {
113     src = i + my_row_base;
114     if (src == rank)
115       continue;
116     recv_offset = (src % two_dsize) * block_size * num_procs;
117     MPI_Irecv(tmp_buff1 + recv_offset, count, recv_type, src, tag, comm,
118               req_ptr++);
119   }
120
121   for (i = 0; i < Y; i++) {
122     dst = i + my_row_base;
123     if (dst == rank)
124       continue;
125     MPI_Send(send_buff, count, send_type, dst, tag, comm);
126   }
127
128   MPI_Waitall(Y - 1, reqs, statuses);
129   req_ptr = reqs;
130
131
132   for (i = 0; i < X; i++) {
133     src = (i * Y + my_col_base);
134     if (src == rank)
135       continue;
136
137     src_row_base = (src / X) * X;
138
139     recv_offset = (src_row_base % two_dsize) * block_size * num_procs;
140     MPI_Irecv(tmp_buff1 + recv_offset, recv_count * num_procs * Y,
141               recv_type, src, tag, comm, req_ptr++);
142   }
143
144   send_offset = (my_row_base % two_dsize) * block_size * num_procs;
145   for (i = 0; i < X; i++) {
146     dst = (i * Y + my_col_base);
147     if (dst == rank)
148       continue;
149     MPI_Send(tmp_buff1 + send_offset, send_count * num_procs * Y, send_type,
150              dst, tag, comm);
151   }
152
153   MPI_Waitall(X - 1, reqs, statuses);
154   req_ptr = reqs;
155
156   for (i = 0; i < two_dsize; i++) {
157     send_offset = (rank * block_size) + (i * block_size * num_procs);
158     recv_offset = (my_z_base * block_size) + (i * block_size);
159     MPI_Sendrecv(tmp_buff1 + send_offset, send_count, send_type, rank, tag,
160                  recv_buff + recv_offset, recv_count, recv_type, rank, tag,
161                  comm, &status);
162   }
163
164   for (i = 1; i < Z; i++) {
165     src = (rank + i * two_dsize) % num_procs;
166     src_z_base = (src / two_dsize) * two_dsize;
167
168     recv_offset = (src_z_base * block_size);
169
170     MPI_Irecv(recv_buff + recv_offset, recv_count * two_dsize, recv_type,
171               src, tag, comm, req_ptr++);
172   }
173
174   for (i = 1; i < Z; i++) {
175     dst = (rank + i * two_dsize) % num_procs;
176
177     recv_offset = 0;
178     for (j = 0; j < two_dsize; j++) {
179       send_offset = (dst + j * num_procs) * block_size;
180       MPI_Sendrecv(tmp_buff1 + send_offset, send_count, send_type,
181                    rank, tag, tmp_buff2 + recv_offset, recv_count,
182                    recv_type, rank, tag, comm, &status);
183
184       recv_offset += block_size;
185     }
186
187     MPI_Send(tmp_buff2, send_count * two_dsize, send_type, dst, tag, comm);
188
189   }
190
191   MPI_Waitall(Z - 1, reqs, statuses);
192
193   free(reqs);
194   free(statuses);
195   free(tmp_buff1);
196   free(tmp_buff2);
197   return success;
198 }