Logo AND Algorithmique Numérique Distribuée

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