Logo AND Algorithmique Numérique Distribuée

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