Doxygen 1.9.8
Toolkit for Adaptive Stochastic Modeling and Non-Intrusive ApproximatioN: Tasmanian v8.2
 
Loading...
Searching...
No Matches
tsgMPIScatterGrid.hpp
Go to the documentation of this file.
1/*
2 * Copyright (c) 2017, Miroslav Stoyanov
3 *
4 * This file is part of
5 * Toolkit for Adaptive Stochastic Modeling And Non-Intrusive ApproximatioN: TASMANIAN
6 *
7 * Redistribution and use in source and binary forms, with or without modification, are permitted provided that the following conditions are met:
8 *
9 * 1. Redistributions of source code must retain the above copyright notice, this list of conditions and the following disclaimer.
10 *
11 * 2. Redistributions in binary form must reproduce the above copyright notice, this list of conditions
12 * and the following disclaimer in the documentation and/or other materials provided with the distribution.
13 *
14 * 3. Neither the name of the copyright holder nor the names of its contributors may be used to endorse
15 * or promote products derived from this software without specific prior written permission.
16 *
17 * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES,
18 * INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED.
19 * IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY,
20 * OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA,
21 * OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY,
22 * OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
23 *
24 * UT-BATTELLE, LLC AND THE UNITED STATES GOVERNMENT MAKE NO REPRESENTATIONS AND DISCLAIM ALL WARRANTIES, BOTH EXPRESSED AND IMPLIED.
25 * THERE ARE NO EXPRESS OR IMPLIED WARRANTIES OF MERCHANTABILITY OR FITNESS FOR A PARTICULAR PURPOSE, OR THAT THE USE OF THE SOFTWARE WILL NOT INFRINGE ANY PATENT,
26 * COPYRIGHT, TRADEMARK, OR OTHER PROPRIETARY RIGHTS, OR THAT THE SOFTWARE WILL ACCOMPLISH THE INTENDED RESULTS OR THAT THE SOFTWARE OR ITS USE WILL NOT RESULT IN INJURY OR DAMAGE.
27 * THE USER ASSUMES RESPONSIBILITY FOR ALL LIABILITIES, PENALTIES, FINES, CLAIMS, CAUSES OF ACTION, AND COSTS AND EXPENSES, CAUSED BY, RESULTING FROM OR ARISING OUT OF,
28 * IN WHOLE OR IN PART THE USE, STORAGE OR DISPOSAL OF THE SOFTWARE.
29 */
30
31#ifndef __TASMANIAN_ADDONS_MPIGRIDSCATTER_HPP
32#define __TASMANIAN_ADDONS_MPIGRIDSCATTER_HPP
33
45#include "tsgAddonsCommon.hpp"
46
56#ifdef Tasmanian_ENABLE_MPI
57
58namespace TasGrid{
59
72class VectorToStreamBuffer : public std::basic_streambuf<char, std::char_traits<char>>{
73public:
75 VectorToStreamBuffer(std::vector<char> &data){
76 setg(&*data.begin(), data.data(), &*data.end());
77 }
78};
79
87inline int getMPIRank(MPI_Comm comm){ int rank; MPI_Comm_rank(comm, &rank); return rank; }
88
126template<bool binary = TasGrid::mode_binary>
127int MPIGridSend(TasmanianSparseGrid const &grid, int destination, int tag, MPI_Comm comm){
128 std::stringstream ss;
129 grid.write(ss, binary);
130 while(ss.str().size() % 16 != 0) ss << " ";
131 return MPI_Send(ss.str().c_str(), (int) (ss.str().size() / 16), MPI_LONG_DOUBLE, destination, tag, comm);
132}
133
154template<bool binary = TasGrid::mode_binary>
155int MPIGridRecv(TasmanianSparseGrid &grid, int source, int tag, MPI_Comm comm, MPI_Status *status = MPI_STATUS_IGNORE){
156 MPI_Status internal_status;
157 if (status == MPI_STATUS_IGNORE) status = &internal_status;
158
159 int short_data_size;
160 MPI_Probe(source, tag, comm, status);
161 MPI_Get_count(status, MPI_LONG_DOUBLE, &short_data_size);
162
163 size_t data_size = Utils::size_mult(short_data_size, 16);
164
165 std::vector<char> buff(data_size);
166 auto result = MPI_Recv(buff.data(), (int) (data_size / 16), MPI_LONG_DOUBLE, source, tag, comm, status);
167
168 VectorToStreamBuffer data_buffer(buff); // do not modify buff after this point
169 std::istream is(&data_buffer);
170 grid.read(is, binary);
171 return result;
172}
173
210template<bool binary = TasGrid::mode_binary>
211int MPIGridBcast(TasmanianSparseGrid &grid, int root, MPI_Comm comm){
212 int me = getMPIRank(comm); // my rank within the comm
213 if (me == root){ // sends the grid
214 std::stringstream ss;
215 grid.write(ss, binary);
216
217 while(ss.str().size() % 16 != 0) ss << " "; // pad with empty chars to align to 16 bytes, i.e., long double
218
219 unsigned long long data_size = (unsigned long long) ss.str().size();
220 auto result = MPI_Bcast(&data_size, 1, MPI_UNSIGNED_LONG_LONG, me, comm);
221 if (result != MPI_SUCCESS) return result;
222
223 return MPI_Bcast(const_cast<char*>(ss.str().c_str()), (int) (data_size / 16), MPI_LONG_DOUBLE, me, comm); // Bcast root does not modify the buffer, this is const-correct
224 }else{ // receives the grid
225 unsigned long long data_size;
226
227 auto result = MPI_Bcast(&data_size, 1, MPI_UNSIGNED_LONG_LONG, root, comm);
228 if (result != MPI_SUCCESS) return result;
229
230 std::vector<char> buff((size_t) data_size);
231 result = MPI_Bcast(buff.data(), (int) (buff.size() / 16), MPI_LONG_DOUBLE, root, comm);
232
233 VectorToStreamBuffer data_buffer(buff); // do not modify buff after this point
234 std::istream is(&data_buffer);
235 grid.read(is, binary);
236 return result;
237 }
238}
239
275template<bool binary = TasGrid::mode_binary>
276int MPIGridScatterOutputs(TasmanianSparseGrid const &source, TasmanianSparseGrid &destination, int root, int tag, MPI_Comm comm){
277 int me = getMPIRank(comm); // my rank within the comm
278
279 if (me == root){ // splitting and sending the grid
280 int num_ranks; MPI_Comm_size(comm, &num_ranks);
281 int num_effective_ranks = std::min(num_ranks, source.getNumOutputs());
282
283 int stride = source.getNumOutputs() / num_effective_ranks;
284 int extras = source.getNumOutputs() % num_effective_ranks;
285
286 // return the starting offset for the given rank
287 auto offset = [&](int rank)->int{ return rank * stride + std::min(rank, extras); };
288
289 for(int rank=0; rank<num_effective_ranks; rank++){
290 if (rank == root){ // this is me, take my own copy of the grid
291 destination = copyGrid(source, offset(rank), offset(rank+1));
292 }else{ // send the grid out
293 auto result = MPIGridSend<binary>(copyGrid(source, offset(rank), offset(rank+1)) , rank, tag, comm);
294 if (result != MPI_SUCCESS) return result;
295 }
296 }
297 for(int rank=num_effective_ranks; rank<num_ranks; rank++){ // if there are any grid remaining, set those to empty
298 if (rank == root){ // this is me, take my own copy of the grid
299 destination = makeEmpty();
300 }else{ // send the grid out
301 auto result = MPIGridSend<binary>(makeEmpty() , rank, tag, comm);
302 if (result != MPI_SUCCESS) return result;
303 }
304 }
305 return MPI_SUCCESS; // if we got here, all was successful
306 }else{ // receiving a grid
307 return MPIGridRecv<binary>(destination, root, tag, comm);
308 }
309}
310
311}
312
313#endif // Tasmanian_ENABLE_MPI
314
315#endif
The master-class that represents an instance of a Tasmanian sparse grid.
Definition TasmanianSparseGrid.hpp:293
int getNumOutputs() const
Return the outputs of the grid, i.e., number of model outputs.
Definition TasmanianSparseGrid.hpp:644
void write(const char *filename, bool binary=mode_binary) const
Write the grid to the given filename using either binary or ASCII format.
void read(const char *filename)
Read the grid from the given filename, automatically detect the format.
Coverts a vector to basic stream-buffer.
Definition tsgMPIScatterGrid.hpp:72
VectorToStreamBuffer(std::vector< char > &data)
Make a stream-buffer from the data vector.
Definition tsgMPIScatterGrid.hpp:75
int MPIGridScatterOutputs(TasmanianSparseGrid const &source, TasmanianSparseGrid &destination, int root, int tag, MPI_Comm comm)
Split the grid across the comm where each rank receives an equal portion of the total outputs.
Definition tsgMPIScatterGrid.hpp:276
int MPIGridRecv(TasmanianSparseGrid &grid, int source, int tag, MPI_Comm comm, MPI_Status *status=MPI_STATUS_IGNORE)
Receive a grid from another process in the MPI comm.
Definition tsgMPIScatterGrid.hpp:155
int MPIGridSend(TasmanianSparseGrid const &grid, int destination, int tag, MPI_Comm comm)
Send a grid to another process in the MPI comm.
Definition tsgMPIScatterGrid.hpp:127
int getMPIRank(MPI_Comm comm)
Utility to return the rank within the given comm.
Definition tsgMPIScatterGrid.hpp:87
int MPIGridBcast(TasmanianSparseGrid &grid, int root, MPI_Comm comm)
Broadcast a grid to all processes in an MPI comm.
Definition tsgMPIScatterGrid.hpp:211
TasmanianSparseGrid makeEmpty()
Returns an empty sparse grid.
Definition TasmanianSparseGrid.hpp:2270
TasmanianSparseGrid copyGrid(TasmanianSparseGrid const &source, int outputs_begin=0, int outputs_end=-1)
Returns a grid that is a copy of the source.
Definition TasmanianSparseGrid.hpp:2382
size_t size_mult(IntA a, IntB b)
Converts two integer-like variables to size_t and returns the product..
Definition tsgUtils.hpp:82
Encapsulates the Tasmanian Sparse Grid module.
Definition TasmanianSparseGrid.hpp:68
Common includes and methods for all addons.