Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
5745f2f19c797b1bc2a34162b0322a00ee92fb6b
[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, MPI_Datatype send_type,
23                void * recv_buff, int recv_count, MPI_Datatype recv_type,
24                MPI_Comm comm)
25 {
26   MPI_Status status;
27   MPI_Aint extent;
28   MPI_Datatype new_type;
29  
30   int * blocks_length, * disps;
31   int i, src, dst, rank, num_procs, count, remainder, block, position;
32   int pack_size, tag = 1, pof2 = 1, success = 1, failure = 0;
33   
34
35   char * tmp_buff;  
36   char * send_ptr = (char *) send_buff;
37   char * recv_ptr = (char *) recv_buff;
38
39   MPI_Comm_size(comm, &num_procs);
40   MPI_Comm_rank(comm, &rank);
41
42   MPI_Type_extent(recv_type, &extent);
43
44   tmp_buff = (char *) malloc (num_procs * recv_count * extent); 
45   if (!tmp_buff)
46     {
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     {
55       printf("alltoall-bruck:61: cannot allocate memory\n");
56       MPI_Finalize();
57       exit(failure);
58     }
59   
60   blocks_length = (int *) malloc(sizeof(int) * num_procs);  
61   if (!blocks_length)
62     {
63       printf("alltoall-bruck:69: cannot allocate memory\n");
64       MPI_Finalize();
65       exit(failure);
66     }
67
68   
69   MPI_Sendrecv(send_ptr + rank * send_count * extent,
70                 (num_procs - rank) * send_count, send_type, rank, tag,
71                 recv_ptr, (num_procs - rank) * recv_count, recv_type, rank,
72                 tag, comm, &status);
73
74   MPI_Sendrecv(send_ptr, rank * send_count, send_type, rank, tag,
75                 recv_ptr + (num_procs - rank) * recv_count * extent,
76                 rank * recv_count, recv_type, rank, tag, comm, &status);  
77
78
79
80   MPI_Pack_size(send_count * num_procs, send_type, comm, &pack_size); 
81
82   while (pof2 < num_procs)
83     {
84       dst = (rank + pof2) % num_procs;
85       src = (rank - pof2 + num_procs) % num_procs;
86
87               
88       count = 0;
89       for (block = 1; block < num_procs; block++)
90         if (block & pof2) 
91           {
92             blocks_length[count] = send_count;
93             disps[count] = block * send_count;
94             count++;
95           }
96           
97       MPI_Type_indexed(count, blocks_length, disps, recv_type, &new_type);
98       MPI_Type_commit(&new_type);
99       
100       position = 0;
101       MPI_Pack(recv_buff, 1, new_type, tmp_buff, pack_size, &position, comm);
102       
103       MPI_Sendrecv(tmp_buff, position, MPI_PACKED, dst, tag, recv_buff, 1,
104                     new_type, src, tag, comm, &status);
105       MPI_Type_free(&new_type);
106       
107       pof2 *= 2;
108     }
109
110   free(disps);
111   free(blocks_length);
112   
113   MPI_Sendrecv(recv_ptr + (rank + 1) * recv_count * extent,
114                (num_procs - rank - 1) * recv_count, send_type,
115                 rank, tag, tmp_buff, (num_procs - rank - 1) * recv_count,
116                 recv_type, rank, tag, comm, &status);
117   
118   MPI_Sendrecv(recv_ptr, (rank + 1) * recv_count, send_type, rank, tag,
119                 tmp_buff + (num_procs - rank - 1) * recv_count * extent,
120                 (rank + 1) * recv_count, recv_type, rank, tag, comm, &status);
121
122       
123   for (i = 0; i < num_procs; i++) 
124     MPI_Sendrecv(tmp_buff + i * recv_count * extent, recv_count, send_type,
125                   rank, tag,
126                   recv_ptr + (num_procs - i - 1) * recv_count * extent,
127                   recv_count, recv_type, rank, tag, comm, &status);
128
129   free(tmp_buff);
130   return success;
131 }