29#include <meteoio/MeteoIO.h>
101 std::string
name()
const;
116 void gather(
const int& send_value, std::vector<int>& receive_vector,
const size_t& root = 0);
128 void reduce_max(
double& value,
const bool all=
true);
129 void reduce_min(
double& value,
const bool all=
true);
130 void reduce_sum(
double& value,
const bool all=
true);
131 void reduce_sum(
int& value,
const bool all=
true);
142 template<
class T>
static void deserialize(
const void* in,
const size_t& len, T& obj)
144 std::stringstream obj_stream;
145 obj_stream.write(
reinterpret_cast<const char*
>(in), len);
157 template<
class T>
static size_t serialize(
void** out,
const T& obj,
const bool alloc=
false)
159 std::stringstream obj_stream;
161 const size_t len = obj_stream.str().size();
163 *out =
new char[len];
165 obj_stream.read(
reinterpret_cast<char*
>(*out), len);
179 template <
class T>
static void op_sum_func(
void* in,
void* out,
int* , MPI_Datatype* datatype)
182 MPI_Type_size(*datatype, &tmp_size);
186 deserialize(in,
static_cast<size_t>(tmp_size), in_obj);
187 deserialize(out,
static_cast<size_t>(tmp_size), out_obj);
203 template <
class T>
void reduce_sum(T& obj,
const bool all=
true)
205 if (size_ <= 1)
return;
208 MPI_Op_create(op_sum_func<T>,
true, &op);
210 void* in_obj_char=NULL;
211 void* out_obj_char=NULL;
213 const size_t len =
serialize(&in_obj_char, obj,
true);
214 out_obj_char =
new char[len];
216 MPI_Datatype chars_datatype;
217 MPI_Type_contiguous(
static_cast<int>(len), MPI_CHAR, &chars_datatype);
218 MPI_Type_commit(&chars_datatype);
221 MPI_Allreduce(in_obj_char, out_obj_char, 1, chars_datatype, op, MPI_COMM_WORLD);
223 MPI_Reduce(in_obj_char, out_obj_char, 1, chars_datatype, op,
master_rank(), MPI_COMM_WORLD);
228 delete[] (
char*)out_obj_char;
229 delete[] (
char*)in_obj_char;
232 template <
class T>
void reduce_sum(T& ,
const bool =
true) {}
242 template <
class T>
void broadcast(T& obj,
const size_t& root = 0)
244 if (size_ <= 1)
return;
246 std::string obj_string;
249 std::stringstream obj_stream;
251 obj_string.insert(0, obj_stream.str());
257 std::stringstream obj_stream;
258 obj_stream << obj_string;
263 template <
class T>
void broadcast(T& ,
const size_t& root=0) {(void)root;}
274 template <
class T>
void broadcast(std::vector<T>& vec_obj,
const size_t& root = 0)
276 if (size_ <= 1)
return;
278 std::string obj_string;
282 std::stringstream objs_stream;
283 for (
size_t ii=0; ii<vec_obj.size(); ii++)
284 objs_stream << vec_obj[ii];
285 obj_string.insert(0, objs_stream.str());
286 vec_size = vec_obj.size();
295 std::stringstream objs_stream;
296 objs_stream << obj_string;
299 for (
size_t ii=0; ii<vec_size; ii++) {
301 vec_obj.push_back(tmp);
306 template <
class T>
void broadcast(std::vector<T>& ,
const size_t& root=0) {(void)root;}
320 template <
class T>
void scatter(std::vector<T*>& vec_local,
const size_t& root = 0)
322 if (size_ <= 1)
return;
325 for (
size_t ii=1; ii<size_; ii++) {
326 size_t startx, deltax;
329 send((
int)deltax, ii);
330 for (
size_t jj=startx; jj<(startx+deltax); jj++) {
331 std::string obj_string;
332 std::stringstream objs_stream;
334 objs_stream << *(vec_local[jj]);
335 obj_string.insert(0, objs_stream.str());
336 send(obj_string, ii);
337 delete vec_local[jj];
341 size_t startx, deltax;
343 vec_local.resize(deltax);
345 std::string obj_string;
348 vec_local.resize((
size_t)deltax);
351 for (
size_t ii=0; ii<vec_local.size(); ii++) {
352 std::stringstream objs_stream;
355 objs_stream << obj_string;
358 vec_local[ii] =
new T(obj);
363 template <
class T>
void scatter(std::vector<T*>& ,
const size_t& root=0) {(void)root;}
374 template <
class T>
void send(
const std::vector<T*>& vec_local,
const size_t& destination,
const int& tag=0);
375 template <
class T>
void send(
const std::vector<T>& vec_local,
const size_t& destination,
const int& tag=0);
377 template <
class T>
void send(
const std::vector<T*>& ,
const size_t& ,
const int& tag=0) {(void)tag;}
378 template <
class T>
void send(
const std::vector<T>& ,
const size_t& ,
const int& tag=0) {(void)tag;}
389 template <
class T>
void receive(std::vector<T*>& vec_local,
const size_t& source,
const int& tag=0);
390 template <
class T>
void receive(std::vector<T>& vec_local,
const size_t& source,
const int& tag=0);
392 template <
class T>
void receive(std::vector<T*>& ,
const size_t& ,
const int& tag=0) {(void)tag;}
393 template <
class T>
void receive(std::vector<T>& ,
const size_t& ,
const int& tag=0) {(void)tag;}
407 template <
class T>
void gather(std::vector<T*>& vec_local,
const size_t& root = 0)
409 if (size_ <= 1)
return;
411 std::vector<int> vec_sizes;
412 const size_t sum = vec_local.size();
414 gather(vec_local.size(), vec_sizes);
418 std::string obj_string;
419 size_t local_sum = vec_local.size();
421 vec_local.resize(sum);
424 for (
size_t ii=1; ii<vec_sizes.size(); ii++) {
425 for (
size_t jj=0; jj<vec_sizes[ii]; jj++) {
426 std::stringstream objs_stream;
429 objs_stream << obj_string;
431 vec_local[local_sum+jj] =
new T(obj);
433 local_sum += vec_sizes[ii];
437 for (
size_t ii=0; ii<vec_local.size(); ii++) {
438 std::string obj_string;
439 std::stringstream objs_stream;
441 objs_stream << *(vec_local[ii]);
442 obj_string.insert(0, objs_stream.str());
443 send(obj_string, root);
448 template <
class T>
void gather(std::vector<T*>& ,
const size_t& root=0) {(void)root;}
461 void getArraySliceParamsOptim(
const size_t& dimx,
size_t& startx_sub,
size_t& nx_sub,
const mio::DEMObject& dem,
const mio::Grid2DObject& landuse)
464 void getArraySliceParamsOptim(
const size_t& dimx,
const size_t& idx_wk,
size_t& startx_sub,
size_t& nx_sub,
const mio::DEMObject& dem,
const mio::Grid2DObject& landuse);
476 static void getArraySliceParams(
const size_t& dimx,
const size_t& nbworkers,
const size_t& idx_wk,
size_t& startx_sub,
size_t& nx_sub);
487 static void getArraySliceParams(
const size_t& dimx,
const size_t& nbworkers, std::vector<size_t>& offset, std::vector<size_t>& nx);
495 void broadcast(std::string& message,
const size_t& root);
497 void send(std::string& message,
const size_t& recipient,
const int& tag=0);
498 void receive(std::string& message,
const size_t& source,
const int& tag=0);
499 void send(
const int& value,
const size_t& recipient,
const int& tag=0);
500 void receive(
int& value,
const size_t& source,
const int& tag=0);
502 static void checkSuccess(
const int& ierr);
A singleton class that deals with all aspects of parallelization within Alpine3D (MPI,...
Definition: MPIControl.h:44
std::string name() const
Definition: MPIControl.cc:90
void reduce_min(double &value, const bool all=true)
Definition: MPIControl.cc:217
void scatter(std::vector< T * > &vec_local, const size_t &root=0)
Scatter the objects pointed to by vector<T*> by slices to all preocesses. Internally no MPI_Scatterv ...
Definition: MPIControl.h:320
void broadcast(std::vector< T > &vec_obj, const size_t &root=0)
Broadcast a vector of class T objects via MPI or in case MPI is not activated don't do anything.
Definition: MPIControl.h:274
void getArraySliceParamsOptim(const size_t &dimx, size_t &startx_sub, size_t &nx_sub, const mio::DEMObject &dem, const mio::Grid2DObject &landuse)
Definition: MPIControl.h:461
size_t max_threads() const
Definition: MPIControl.cc:66
size_t thread() const
Definition: MPIControl.cc:57
void barrier() const
Definition: MPIControl.cc:247
void broadcast(T &obj, const size_t &root=0)
Broadcast an object via MPI or in case MPI is not activated don't do anything.
Definition: MPIControl.h:242
void gather(const int &send_value, std::vector< int > &receive_vector, const size_t &root=0)
Definition: MPIControl.cc:195
void send(const std::vector< T * > &vec_local, const size_t &destination, const int &tag=0)
Send the objects pointed to by vector<T*> to process #destination.
Definition: MPIControl.cc:299
bool master() const
Definition: MPIControl.cc:95
void reduce_sum(double &value, const bool all=true)
Definition: MPIControl.cc:227
static void getArraySliceParams(const size_t &dimx, const size_t &nbworkers, std::vector< size_t > &offset, std::vector< size_t > &nx)
Returns the parameters for splitting an array in several, balanced sub-arrays. This is mostly usefull...
void gather(std::vector< T * > &vec_local, const size_t &root=0)
Gathers the objects pointed to by vector<T*> from all processes into a vector<T*> on the root node....
Definition: MPIControl.h:407
size_t size() const
Definition: MPIControl.cc:85
void getArraySliceParams(const size_t &dimx, size_t &startx_sub, size_t &nx_sub) const
Helper to split an array for domain decomposition. Given the overall size of an array calculate a sta...
Definition: MPIControl.h:458
void getArraySliceParams(const size_t &dimx, const size_t &idx_wk, size_t &startx_sub, size_t &nx_sub) const
Definition: MPIControl.h:459
static void op_sum_func(void *in, void *out, int *, MPI_Datatype *datatype)
Definition: MPIControl.h:179
void receive(std::vector< T * > &vec_local, const size_t &source, const int &tag=0)
Receive vector of objects from process #source.
Definition: MPIControl.cc:322
static void deserialize(const void *in, const size_t &len, T &obj)
Definition: MPIControl.h:142
void reduce_sum(T &obj, const bool all=true)
Adds up the values of type T from all processes and, if all==true, distributes the sum back to all pr...
Definition: MPIControl.h:203
size_t rank() const
Definition: MPIControl.cc:75
size_t master_rank() const
Definition: MPIControl.cc:80
static MPIControl & instance()
Definition: MPIControl.cc:417
void reduce_max(double &value, const bool all=true)
Definition: MPIControl.cc:207
static size_t serialize(void **out, const T &obj, const bool alloc=false)
Definition: MPIControl.h:157
bool openmp() const
Definition: MPIControl.cc:48