Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Merge branch 'hypervisor' of scm.gforge.inria.fr:/gitroot/simgrid/simgrid into hypervisor
[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 = 1, pof2 = 1, success = 1, failure = 0;
34
35
36   char *tmp_buff;
37   char *send_ptr = (char *) send_buff;
38   char *recv_ptr = (char *) recv_buff;
39
40   MPI_Comm_size(comm, &num_procs);
41   MPI_Comm_rank(comm, &rank);
42
43   MPI_Type_extent(recv_type, &extent);
44
45   tmp_buff = (char *) malloc(num_procs * recv_count * extent);
46   if (!tmp_buff) {
47     printf("alltoall-bruck:53: cannot allocate memory\n");
48     MPI_Finalize();
49     exit(failure);
50   }
51
52   disps = (int *) malloc(sizeof(int) * num_procs);
53   if (!disps) {
54     printf("alltoall-bruck:61: cannot allocate memory\n");
55     MPI_Finalize();
56     exit(failure);
57   }
58
59   blocks_length = (int *) malloc(sizeof(int) * num_procs);
60   if (!blocks_length) {
61     printf("alltoall-bruck:69: cannot allocate memory\n");
62     MPI_Finalize();
63     exit(failure);
64   }
65
66
67   MPI_Sendrecv(send_ptr + rank * send_count * extent,
68                (num_procs - rank) * send_count, send_type, rank, tag,
69                recv_ptr, (num_procs - rank) * recv_count, recv_type, rank,
70                tag, comm, &status);
71
72   MPI_Sendrecv(send_ptr, rank * send_count, send_type, rank, tag,
73                recv_ptr + (num_procs - rank) * recv_count * extent,
74                rank * recv_count, recv_type, rank, tag, comm, &status);
75
76
77
78   MPI_Pack_size(send_count * num_procs, send_type, comm, &pack_size);
79
80   while (pof2 < num_procs) {
81     dst = (rank + pof2) % num_procs;
82     src = (rank - pof2 + num_procs) % num_procs;
83
84
85     count = 0;
86     for (block = 1; block < num_procs; block++)
87       if (block & pof2) {
88         blocks_length[count] = send_count;
89         disps[count] = block * send_count;
90         count++;
91       }
92
93     MPI_Type_indexed(count, blocks_length, disps, recv_type, &new_type);
94     MPI_Type_commit(&new_type);
95
96     position = 0;
97     MPI_Pack(recv_buff, 1, new_type, tmp_buff, pack_size, &position, comm);
98
99     MPI_Sendrecv(tmp_buff, position, MPI_PACKED, dst, tag, recv_buff, 1,
100                  new_type, src, tag, comm, &status);
101     MPI_Type_free(&new_type);
102
103     pof2 *= 2;
104   }
105
106   free(disps);
107   free(blocks_length);
108
109   MPI_Sendrecv(recv_ptr + (rank + 1) * recv_count * extent,
110                (num_procs - rank - 1) * recv_count, send_type,
111                rank, tag, tmp_buff, (num_procs - rank - 1) * recv_count,
112                recv_type, rank, tag, comm, &status);
113
114   MPI_Sendrecv(recv_ptr, (rank + 1) * recv_count, send_type, rank, tag,
115                tmp_buff + (num_procs - rank - 1) * recv_count * extent,
116                (rank + 1) * recv_count, recv_type, rank, tag, comm, &status);
117
118
119   for (i = 0; i < num_procs; i++)
120     MPI_Sendrecv(tmp_buff + i * recv_count * extent, recv_count, send_type,
121                  rank, tag,
122                  recv_ptr + (num_procs - i - 1) * recv_count * extent,
123                  recv_count, recv_type, rank, tag, comm, &status);
124
125   free(tmp_buff);
126   return success;
127 }