Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
leaks -- with smp algos
[simgrid.git] / src / smpi / colls / bcast-scatter-LR-allgather.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
9 /*****************************************************************************
10
11 Copyright (c) 2006, Ahmad Faraj & Xin Yuan,
12 All rights reserved.
13
14 Redistribution and use in source and binary forms, with or without
15 modification, are permitted provided that the following conditions are met:
16
17   * Redistributions of source code must retain the above copyright notice,
18     this list of conditions and the following disclaimer.
19
20   * Redistributions in binary form must reproduce the above copyright notice,
21     this list of conditions and the following disclaimer in the documentation
22     and/or other materials provided with the distribution.
23
24   * Neither the name of the Florida State University nor the names of its
25     contributors may be used to endorse or promote products derived from this
26     software without specific prior written permission.
27
28 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
29 ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
30 WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
31 DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR
32 ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
33 (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
34 LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON
35 ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
36 (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
37 SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
38
39   *************************************************************************
40   *     Any results obtained from executing this software require the     *
41   *     acknowledgment and citation of the software and its owners.       *
42   *     The full citation is given below:                                 *
43   *                                                                       *
44   *     A. Faraj and X. Yuan. "Automatic Generation and Tuning of MPI     *
45   *     Collective Communication Routines." The 19th ACM International    *
46   *     Conference on Supercomputing (ICS), Cambridge, Massachusetts,     *
47   *     June 20-22, 2005.                                                 *
48   *************************************************************************
49
50 *****************************************************************************/
51
52 /*****************************************************************************
53
54  * Function: bcast_scatter_LR_allgather
55
56  * Return: int
57
58  * Inputs:
59     buff: send input buffer
60     count: number of elements to send
61     data_type: data type of elements being sent
62     root: source of data
63     comm: communicator
64
65  * Descrp: broadcasts using a scatter followed by LR allgather.
66
67  * Auther: MPIH / modified by Ahmad Faraj
68
69  ****************************************************************************/
70 int
71 smpi_coll_tuned_bcast_scatter_LR_allgather(void *buff, int count,
72                                            MPI_Datatype data_type, int root,
73                                            MPI_Comm comm)
74 {
75   MPI_Aint extent;
76   MPI_Status status;
77   int i, src, dst, rank, num_procs;
78   int mask, relative_rank, curr_size, recv_size, send_size, nbytes;
79   int scatter_size, left, right, next_src, *recv_counts, *disps;
80   int tag = COLL_TAG_BCAST;
81
82   rank = smpi_comm_rank(comm);
83   num_procs = smpi_comm_size(comm);
84   extent = smpi_datatype_get_extent(data_type);
85
86
87   nbytes = extent * count;
88   scatter_size = (nbytes + num_procs - 1) / num_procs;  // ceiling division 
89   curr_size = (rank == root) ? nbytes : 0;      // root starts with all the data
90   relative_rank = (rank >= root) ? rank - root : rank - root + num_procs;
91
92   mask = 0x1;
93   while (mask < num_procs) {
94     if (relative_rank & mask) {
95       src = rank - mask;
96       if (src < 0)
97         src += num_procs;
98       recv_size = nbytes - relative_rank * scatter_size;
99       //  recv_size is larger than what might actually be sent by the
100       //  sender. We don't need compute the exact value because MPI
101       //  allows you to post a larger recv.
102       if (recv_size <= 0)
103         curr_size = 0;          // this process doesn't receive any data
104       // because of uneven division 
105       else {
106         smpi_mpi_recv((char *) buff + relative_rank * scatter_size, recv_size,
107                  MPI_BYTE, src, tag, comm, &status);
108         curr_size = smpi_mpi_get_count(&status, MPI_BYTE);
109       }
110       break;
111     }
112     mask <<= 1;
113   }
114
115   // This process is responsible for all processes that have bits
116   // set from the LSB upto (but not including) mask.  Because of
117   // the "not including", we start by shifting mask back down
118   // one.
119
120   mask >>= 1;
121   while (mask > 0) {
122     if (relative_rank + mask < num_procs) {
123       send_size = curr_size - scatter_size * mask;
124       // mask is also the size of this process's subtree 
125
126       if (send_size > 0) {
127         dst = rank + mask;
128         if (dst >= num_procs)
129           dst -= num_procs;
130         smpi_mpi_send((char *) buff + scatter_size * (relative_rank + mask),
131                  send_size, MPI_BYTE, dst, tag, comm);
132
133         curr_size -= send_size;
134       }
135     }
136     mask >>= 1;
137   }
138
139   // done scatter now do allgather
140   recv_counts = (int *) xbt_malloc(sizeof(int) * num_procs);
141   disps = (int *) xbt_malloc(sizeof(int) * num_procs);
142
143   for (i = 0; i < num_procs; i++) {
144     recv_counts[i] = nbytes - i * scatter_size;
145     if (recv_counts[i] > scatter_size)
146       recv_counts[i] = scatter_size;
147     if (recv_counts[i] < 0)
148       recv_counts[i] = 0;
149   }
150
151   disps[0] = 0;
152   for (i = 1; i < num_procs; i++)
153     disps[i] = disps[i - 1] + recv_counts[i - 1];
154
155   left = (num_procs + rank - 1) % num_procs;
156   right = (rank + 1) % num_procs;
157
158   src = rank;
159   next_src = left;
160
161   for (i = 1; i < num_procs; i++) {
162     smpi_mpi_sendrecv((char *) buff + disps[(src - root + num_procs) % num_procs],
163                  recv_counts[(src - root + num_procs) % num_procs],
164                  MPI_BYTE, right, tag,
165                  (char *) buff +
166                  disps[(next_src - root + num_procs) % num_procs],
167                  recv_counts[(next_src - root + num_procs) % num_procs],
168                  MPI_BYTE, left, tag, comm, &status);
169     src = next_src;
170     next_src = (num_procs + next_src - 1) % num_procs;
171   }
172
173
174   free(recv_counts);
175   free(disps);
176
177   return MPI_SUCCESS;
178 }