Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Add/update copyright notices.
[simgrid.git] / src / smpi / colls / alltoall-2dmesh.c
1 /* Copyright (c) 2013-2014. The SimGrid Team.
2  * All rights reserved.                                                     */
3
4 /* This program is free software; you can redistribute it and/or modify it
5  * under the terms of the license (GNU LGPL) which comes with this package. */
6
7 #include "colls_private.h"
8 #include <math.h>
9
10 /*****************************************************************************
11
12  * Function: alltoall_2dmesh_shoot
13
14  * Return: int
15
16  * Inputs:
17     send_buff: send input buffer
18     send_count: number of elements to send
19     send_type: data type of elements being sent
20     recv_buff: receive output buffer
21     recv_count: number of elements to received
22     recv_type: data type of elements being received
23     comm: communicator
24
25  * Descrp: Function realizes the alltoall operation using the 2dmesh
26            algorithm. It actually performs allgather operation in x dimension
27            then in the y dimension. Each node then extracts the needed data.
28            The communication in each dimension follows "simple."
29  
30  * Auther: Ahmad Faraj
31
32 ****************************************************************************/
33 static int alltoall_check_is_2dmesh(int num, int *i, int *j)
34 {
35   int x, max = num / 2;
36   x = sqrt(num);
37
38   while (x <= max) {
39     if ((num % x) == 0) {
40       *i = x;
41       *j = num / x;
42
43       if (*i > *j) {
44         x = *i;
45         *i = *j;
46         *j = x;
47       }
48
49       return 1;
50     }
51     x++;
52   }
53   return 0;
54 }
55
56 int smpi_coll_tuned_alltoall_2dmesh(void *send_buff, int send_count,
57                                     MPI_Datatype send_type,
58                                     void *recv_buff, int recv_count,
59                                     MPI_Datatype recv_type, MPI_Comm comm)
60 {
61   MPI_Status *statuses, s;
62   MPI_Request *reqs, *req_ptr;;
63   MPI_Aint extent;
64
65   char *tmp_buff1, *tmp_buff2;
66   int i, j, src, dst, rank, num_procs, count, num_reqs;
67   int X, Y, send_offset, recv_offset;
68   int my_row_base, my_col_base, src_row_base, block_size;
69   int tag = COLL_TAG_ALLTOALL;
70
71   rank = smpi_comm_rank(comm);
72   num_procs = smpi_comm_size(comm);
73   extent = smpi_datatype_get_extent(send_type);
74
75   if (!alltoall_check_is_2dmesh(num_procs, &X, &Y))
76     return MPI_ERR_OTHER;
77
78   my_row_base = (rank / Y) * Y;
79   my_col_base = rank % Y;
80
81   block_size = extent * send_count;
82
83   tmp_buff1 = (char *) xbt_malloc(block_size * num_procs * Y);
84   tmp_buff2 = (char *) xbt_malloc(block_size * Y);
85
86   num_reqs = X;
87   if (Y > X)
88     num_reqs = Y;
89
90   statuses = (MPI_Status *) xbt_malloc(num_reqs * sizeof(MPI_Status));
91   reqs = (MPI_Request *) xbt_malloc(num_reqs * sizeof(MPI_Request));
92
93   req_ptr = reqs;
94
95   send_offset = recv_offset = (rank % Y) * block_size * num_procs;
96
97   count = send_count * num_procs;
98
99   for (i = 0; i < Y; i++) {
100     src = i + my_row_base;
101     if (src == rank)
102       continue;
103
104     recv_offset = (src % Y) * block_size * num_procs;
105     *(req_ptr++) = smpi_mpi_irecv(tmp_buff1 + recv_offset, count, recv_type, src, tag, comm);
106   }
107
108   for (i = 0; i < Y; i++) {
109     dst = i + my_row_base;
110     if (dst == rank)
111       continue;
112     smpi_mpi_send(send_buff, count, send_type, dst, tag, comm);
113   }
114
115   smpi_mpi_waitall(Y - 1, reqs, statuses);
116   req_ptr = reqs;
117
118   for (i = 0; i < Y; i++) {
119     send_offset = (rank * block_size) + (i * block_size * num_procs);
120     recv_offset = (my_row_base * block_size) + (i * block_size);
121
122     if (i + my_row_base == rank)
123       smpi_mpi_sendrecv((char *) send_buff + recv_offset, send_count, send_type,
124                    rank, tag,
125                    (char *) recv_buff + recv_offset, recv_count, recv_type,
126                    rank, tag, comm, &s);
127
128     else
129       smpi_mpi_sendrecv(tmp_buff1 + send_offset, send_count, send_type,
130                    rank, tag,
131                    (char *) recv_buff + recv_offset, recv_count, recv_type,
132                    rank, tag, comm, &s);
133   }
134
135
136   for (i = 0; i < X; i++) {
137     src = (i * Y + my_col_base);
138     if (src == rank)
139       continue;
140     src_row_base = (src / Y) * Y;
141
142     *(req_ptr++) = smpi_mpi_irecv((char *) recv_buff + src_row_base * block_size, recv_count * Y,
143               recv_type, src, tag, comm);
144   }
145
146   for (i = 0; i < X; i++) {
147     dst = (i * Y + my_col_base);
148     if (dst == rank)
149       continue;
150
151     recv_offset = 0;
152     for (j = 0; j < Y; j++) {
153       send_offset = (dst + j * num_procs) * block_size;
154
155       if (j + my_row_base == rank)
156         smpi_mpi_sendrecv((char *) send_buff + dst * block_size, send_count,
157                      send_type, rank, tag, tmp_buff2 + recv_offset, recv_count,
158                      recv_type, rank, tag, comm, &s);
159       else
160         smpi_mpi_sendrecv(tmp_buff1 + send_offset, send_count, send_type,
161                      rank, tag,
162                      tmp_buff2 + recv_offset, recv_count, recv_type,
163                      rank, tag, comm, &s);
164
165       recv_offset += block_size;
166     }
167
168     smpi_mpi_send(tmp_buff2, send_count * Y, send_type, dst, tag, comm);
169   }
170   smpi_mpi_waitall(X - 1, reqs, statuses);
171   free(reqs);
172   free(statuses);
173   free(tmp_buff1);
174   free(tmp_buff2);
175   return MPI_SUCCESS;
176 }