Doxygen 1.9.8
Toolkit for Adaptive Stochastic Modeling and Non-Intrusive ApproximatioN: Tasmanian v8.2
 
Loading...
Searching...
No Matches
tsgMPIScatterDream.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_MPIDREAMSCATTER_HPP
32#define __TASMANIAN_ADDONS_MPIDREAMSCATTER_HPP
33
45#include "tsgMPIScatterGrid.hpp"
46
56#ifdef Tasmanian_ENABLE_MPI
57
58namespace TasDREAM{
59
85template<class Likelihood>
86int MPILikelihoodSend(Likelihood const &likely, int destination, int tag, MPI_Comm comm, int outputs_begin = 0, int outputs_end = -1){
87 std::stringstream ss;
88 likely.write(ss, outputs_begin, outputs_end);
89 while(ss.str().size() % 16 != 0) ss << " ";
90 return MPI_Send(ss.str().c_str(), (int) (ss.str().size() / 16), MPI_LONG_DOUBLE, destination, tag, comm);
91}
92
115template<class Likelihood>
116int MPILikelihoodRecv(Likelihood &likely, int source, int tag, MPI_Comm comm, MPI_Status *status = MPI_STATUS_IGNORE){
117 MPI_Status internal_status;
118 if (status == MPI_STATUS_IGNORE) status = &internal_status;
119
120 int short_data_size;
121 MPI_Probe(source, tag, comm, status);
122 MPI_Get_count(status, MPI_LONG_DOUBLE, &short_data_size);
123
124 size_t data_size = TasGrid::Utils::size_mult(short_data_size, 16);
125
126 std::vector<char> buff(data_size);
127 auto result = MPI_Recv(buff.data(), (int) (data_size / 16), MPI_LONG_DOUBLE, source, tag, comm, status);
128
129 TasGrid::VectorToStreamBuffer data_buffer(buff); // do not modify buff after this point
130 std::istream is(&data_buffer);
131 likely.read(is);
132 return result;
133}
134
159template<class Likelihood>
160int MPILikelihoodScatter(Likelihood const &source, Likelihood &destination, int root, int tag, MPI_Comm comm){
161 int me = TasGrid::getMPIRank(comm); // my rank within the comm
162
163 if (me == root){ // splitting and sending the grid
164 int num_ranks; MPI_Comm_size(comm, &num_ranks);
165 int num_effective_ranks = std::min(num_ranks, source.getNumOutputs());
166
167 int stride = source.getNumOutputs() / num_effective_ranks;
168 int extras = source.getNumOutputs() % num_effective_ranks;
169
170 // return the starting offset for the given rank
171 auto offset = [&](int rank)->int{ return rank * stride + std::min(rank, extras); };
172
173 for(int rank=0; rank<num_effective_ranks; rank++){
174 if (rank == root){ // this is me, take my own copy of the grid
175 std::stringstream ss;
176 source.write(ss, offset(rank), offset(rank+1));
177 destination.read(ss);
178 }else{ // send the grid out
179 auto result = MPILikelihoodSend(source, rank, tag, comm, offset(rank), offset(rank+1));
180 if (result != MPI_SUCCESS) return result;
181 }
182 }
183 for(int rank=num_effective_ranks; rank<num_ranks; rank++){ // if there are any grid remaining, set those to empty
184 if (rank == root){ // this is me, take my own copy of the grid
185 destination = Likelihood();
186 }else{ // send the grid out
187 auto result = MPILikelihoodSend(Likelihood() , rank, tag, comm);
188 if (result != MPI_SUCCESS) return result;
189 }
190 }
191 return MPI_SUCCESS; // if we got here, all was successful
192 }else{ // receiving a grid
193 return MPILikelihoodRecv(destination, root, tag, comm);
194 }
195}
196
197}
198
199#endif // Tasmanian_ENABLE_MPI
200
201#endif
Coverts a vector to basic stream-buffer.
Definition tsgMPIScatterGrid.hpp:72
int MPILikelihoodRecv(Likelihood &likely, int source, int tag, MPI_Comm comm, MPI_Status *status=MPI_STATUS_IGNORE)
Receive a likelihood from another process in the MPI comm.
Definition tsgMPIScatterDream.hpp:116
int MPILikelihoodScatter(Likelihood const &source, Likelihood &destination, int root, int tag, MPI_Comm comm)
Split the likelihood across the comm where each rank receives an equal portion of the total outputs.
Definition tsgMPIScatterDream.hpp:160
int MPILikelihoodSend(Likelihood const &likely, int destination, int tag, MPI_Comm comm, int outputs_begin=0, int outputs_end=-1)
Send a likelihood to another process in the MPI comm.
Definition tsgMPIScatterDream.hpp:86
int getMPIRank(MPI_Comm comm)
Utility to return the rank within the given comm.
Definition tsgMPIScatterGrid.hpp:87
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 DREAM module.
Definition TasmanianDREAM.hpp:80
Sparse Grids send/receive through MPI.