Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Update copyright headers.
[simgrid.git] / src / smpi / colls / bcast / bcast-scatter-LR-allgather.cpp
1 /* Copyright (c) 2013-2018. The SimGrid Team. All rights reserved.          */
2
3 /* This program is free software; you can redistribute it and/or modify it
4  * under the terms of the license (GNU LGPL) which comes with this package. */
5
6 #include "../colls_private.hpp"
7 #include "smpi_status.hpp"
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 namespace simgrid{
71 namespace smpi{
72 int
73 Coll_bcast_scatter_LR_allgather::bcast(void *buff, int count,
74                                            MPI_Datatype data_type, int root,
75                                            MPI_Comm comm)
76 {
77   MPI_Aint extent;
78   MPI_Status status;
79   int i, src, dst, rank, num_procs;
80   int mask, relative_rank, curr_size, recv_size, send_size, nbytes;
81   int scatter_size, left, right, next_src, *recv_counts, *disps;
82   int tag = COLL_TAG_BCAST;
83
84   rank = comm->rank();
85   num_procs = comm->size();
86   extent = data_type->get_extent();
87
88
89   nbytes = extent * count;
90   scatter_size = (nbytes + num_procs - 1) / num_procs;  // ceiling division
91   curr_size = (rank == root) ? nbytes : 0;      // root starts with all the data
92   relative_rank = (rank >= root) ? rank - root : rank - root + num_procs;
93
94   mask = 0x1;
95   while (mask < num_procs) {
96     if (relative_rank & mask) {
97       src = rank - mask;
98       if (src < 0)
99         src += num_procs;
100       recv_size = nbytes - relative_rank * scatter_size;
101       //  recv_size is larger than what might actually be sent by the
102       //  sender. We don't need compute the exact value because MPI
103       //  allows you to post a larger recv.
104       if (recv_size <= 0)
105         curr_size = 0;          // this process doesn't receive any data
106       // because of uneven division
107       else {
108         Request::recv((char *) buff + relative_rank * scatter_size, recv_size,
109                  MPI_BYTE, src, tag, comm, &status);
110         curr_size = Status::get_count(&status, MPI_BYTE);
111       }
112       break;
113     }
114     mask <<= 1;
115   }
116
117   // This process is responsible for all processes that have bits
118   // set from the LSB upto (but not including) mask.  Because of
119   // the "not including", we start by shifting mask back down
120   // one.
121
122   mask >>= 1;
123   while (mask > 0) {
124     if (relative_rank + mask < num_procs) {
125       send_size = curr_size - scatter_size * mask;
126       // mask is also the size of this process's subtree
127
128       if (send_size > 0) {
129         dst = rank + mask;
130         if (dst >= num_procs)
131           dst -= num_procs;
132         Request::send((char *) buff + scatter_size * (relative_rank + mask),
133                  send_size, MPI_BYTE, dst, tag, comm);
134
135         curr_size -= send_size;
136       }
137     }
138     mask >>= 1;
139   }
140
141   // done scatter now do allgather
142   recv_counts = (int *) xbt_malloc(sizeof(int) * num_procs);
143   disps = (int *) xbt_malloc(sizeof(int) * num_procs);
144
145   for (i = 0; i < num_procs; i++) {
146     recv_counts[i] = nbytes - i * scatter_size;
147     if (recv_counts[i] > scatter_size)
148       recv_counts[i] = scatter_size;
149     if (recv_counts[i] < 0)
150       recv_counts[i] = 0;
151   }
152
153   disps[0] = 0;
154   for (i = 1; i < num_procs; i++)
155     disps[i] = disps[i - 1] + recv_counts[i - 1];
156
157   left = (num_procs + rank - 1) % num_procs;
158   right = (rank + 1) % num_procs;
159
160   src = rank;
161   next_src = left;
162
163   for (i = 1; i < num_procs; i++) {
164     Request::sendrecv((char *) buff + disps[(src - root + num_procs) % num_procs],
165                  recv_counts[(src - root + num_procs) % num_procs],
166                  MPI_BYTE, right, tag,
167                  (char *) buff +
168                  disps[(next_src - root + num_procs) % num_procs],
169                  recv_counts[(next_src - root + num_procs) % num_procs],
170                  MPI_BYTE, left, tag, comm, &status);
171     src = next_src;
172     next_src = (num_procs + next_src - 1) % num_procs;
173   }
174
175
176   free(recv_counts);
177   free(disps);
178
179   return MPI_SUCCESS;
180 }
181
182 }
183 }