Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
fix build without tracing
[simgrid.git] / src / smpi / colls / alltoall-bruck.c
1 /*****************************************************************************
2
3  * Function: alltoall_bruck
4
5  * Return: int
6
7  * inputs:
8     send_buff: send input buffer
9     send_count: number of elements to send
10     send_type: data type of elements being sent
11     recv_buff: receive output buffer
12     recv_count: number of elements to received
13     recv_type: data type of elements being received
14     comm: communicator
15
16  * Descrp: Function realizes the alltoall operation using the bruck algorithm.
17
18  * Auther: MPICH / modified by Ahmad Faraj
19
20  ****************************************************************************/
21 int
22 smpi_coll_tuned_alltoall_bruck(void *send_buff, int send_count,
23                                MPI_Datatype send_type, void *recv_buff,
24                                int recv_count, MPI_Datatype recv_type,
25                                MPI_Comm comm)
26 {
27   MPI_Status status;
28   MPI_Aint extent;
29   MPI_Datatype new_type;
30
31   int *blocks_length, *disps;
32   int i, src, dst, rank, num_procs, count, remainder, block, position;
33   int pack_size, tag = COLL_TAG_ALLTOALL, pof2 = 1;
34
35
36   char *tmp_buff;
37   char *send_ptr = (char *) send_buff;
38   char *recv_ptr = (char *) recv_buff;
39
40   num_procs = smpi_comm_size(comm);
41   rank = smpi_comm_rank(comm);
42
43   extent = smpi_datatype_get_extent(recv_type);
44
45   tmp_buff = (char *) xbt_malloc(num_procs * recv_count * extent);
46   disps = (int *) xbt_malloc(sizeof(int) * num_procs);
47   blocks_length = (int *) xbt_malloc(sizeof(int) * num_procs);
48
49   smpi_mpi_sendrecv(send_ptr + rank * send_count * extent,
50                (num_procs - rank) * send_count, send_type, rank, tag,
51                recv_ptr, (num_procs - rank) * recv_count, recv_type, rank,
52                tag, comm, &status);
53
54   smpi_mpi_sendrecv(send_ptr, rank * send_count, send_type, rank, tag,
55                recv_ptr + (num_procs - rank) * recv_count * extent,
56                rank * recv_count, recv_type, rank, tag, comm, &status);
57
58
59
60   MPI_Pack_size(send_count * num_procs, send_type, comm, &pack_size);
61
62   while (pof2 < num_procs) {
63     dst = (rank + pof2) % num_procs;
64     src = (rank - pof2 + num_procs) % num_procs;
65
66
67     count = 0;
68     for (block = 1; block < num_procs; block++)
69       if (block & pof2) {
70         blocks_length[count] = send_count;
71         disps[count] = block * send_count;
72         count++;
73       }
74
75     MPI_Type_indexed(count, blocks_length, disps, recv_type, &new_type);
76     smpi_datatype_commit(&new_type);
77
78     position = 0;
79     MPI_Pack(recv_buff, 1, new_type, tmp_buff, pack_size, &position, comm);
80
81     smpi_mpi_sendrecv(tmp_buff, position, MPI_PACKED, dst, tag, recv_buff, 1,
82                  new_type, src, tag, comm, &status);
83     smpi_datatype_free(&new_type);
84
85     pof2 *= 2;
86   }
87
88   free(disps);
89   free(blocks_length);
90
91   smpi_mpi_sendrecv(recv_ptr + (rank + 1) * recv_count * extent,
92                (num_procs - rank - 1) * recv_count, send_type,
93                rank, tag, tmp_buff, (num_procs - rank - 1) * recv_count,
94                recv_type, rank, tag, comm, &status);
95
96   smpi_mpi_sendrecv(recv_ptr, (rank + 1) * recv_count, send_type, rank, tag,
97                tmp_buff + (num_procs - rank - 1) * recv_count * extent,
98                (rank + 1) * recv_count, recv_type, rank, tag, comm, &status);
99
100
101   for (i = 0; i < num_procs; i++)
102     smpi_mpi_sendrecv(tmp_buff + i * recv_count * extent, recv_count, send_type,
103                  rank, tag,
104                  recv_ptr + (num_procs - i - 1) * recv_count * extent,
105                  recv_count, recv_type, rank, tag, comm, &status);
106
107   free(tmp_buff);
108   return MPI_SUCCESS;
109 }