ALPSCore reference
mpi_bcast.hpp
Go to the documentation of this file.
1 /*
2  * Copyright (C) 1998-2018 ALPS Collaboration. See COPYRIGHT.TXT
3  * All rights reserved. Use is subject to license terms. See LICENSE.TXT
4  * For use in publications, see ACKNOWLEDGE.TXT
5  */
6 
11 #ifndef ALPS_GF_MPI_BCAST_HPP_c030bec39d4b43b9a24a16b5805f542d
12 #define ALPS_GF_MPI_BCAST_HPP_c030bec39d4b43b9a24a16b5805f542d
13 
14 #include <alps/utilities/mpi.hpp>
15 #include <iostream>
16 namespace alps {
17  namespace gf {
18  namespace detail {
20 
21  template <typename T>
22  void broadcast(const alps::mpi::communicator& comm, std::vector<T>& data, int root) {
23  typedef std::vector<T> data_type;
24  typedef typename data_type::size_type size_type;
25  size_type root_sz=data.size();
26  alps::mpi::broadcast(comm, root_sz, root);
27  data.resize(root_sz);
28  alps::mpi::broadcast(comm, &data[0], root_sz, root);
29  }
30 
32 
39  template <typename T, size_t N>
40  void broadcast(const alps::mpi::communicator& comm, boost::multi_array<T,N>& data, int root)
41  {
42  typedef boost::multi_array<T,N> data_type;
43  typedef typename data_type::index index_type;
44  typedef typename data_type::size_type size_type;
45 
46  int rank=comm.rank();
47  const bool is_root=(rank==root);
48 
49  try {
50  // NOTE: questionable; boost::multi_array does not document comparison of storage orders
51  if (! (data.storage_order()==boost::c_storage_order()) )
52  throw std::logic_error("Unsupported storage order in multi_array broadcast at rank #"+
53  std::to_string(rank));
54 
55  // Compare dimensions with root. Normally should not be needed,
56  // and incurs extra communication cost ==> enabled only in debug mode.
57 #ifndef BOOST_DISABLE_ASSERTS
58  {
59  size_type ndim=N;
60  alps::mpi::broadcast(comm, ndim, root);
61  if (ndim!=N) {
62  throw std::logic_error("Different multi_array dimensions in broadcast:\n"
63  "root (rank #"+
64  std::to_string(root)+
65  ") expects N="+
66  std::to_string(ndim)+
67  ", rank #"+
68  std::to_string(rank)+
69  " has N="+std::to_string(N));
70  }
71  }
72 #endif
73  } catch (const std::exception& exc) {
74  // Abort here and now -- otherwise we end up with unmatched broadcast.
75  std::cerr << "Exception in broadcast of boost::multi_array data on rank " << rank << std::endl
76  << exc.what()
77  << "\nAborting." << std::endl;
78  MPI_Abort(MPI_COMM_WORLD,1);
79 
80  // FIXME: In debug mode we could use MPI_Allreduce() and notify other processes
81  // of the error. Has to be done carefully to avoid unmatched broadcasts.
82  }
83 
84  // Broadcast the bases
85  boost::array<index_type,N> bases;
86  const index_type* root_bases=data.index_bases();
87  if (is_root) {
88  std::copy(root_bases, root_bases+N, bases.begin()); // FIXME: this copy is not needed if done carefully
89  }
90  alps::mpi::broadcast(comm, bases.data(), N, root);
91 
92  // Broadcast the array shape
93  boost::array<size_type,N> shape;
94  if (is_root) {
95  const size_type* root_shape=data.shape();
96  std::copy(root_shape, root_shape+N, shape.begin()); // FIXME: this copy is not needed if done carefully
97  }
98  alps::mpi::broadcast(comm, shape.data(), N, root);
99 
100  if (! is_root) {
101  data.resize(shape);
102  data.reindex(bases);
103  }
104 
105  size_t nelements=data.num_elements();
106 
107  // This is an additional broadcast, but we need to ensure MPI broadcast correctness,
108  // otherwise it would be a hell to debug down the road figuring out what exactly went wrong.
109  unsigned long nelements_root=nelements;
110  alps::mpi::broadcast(comm, nelements_root, root);
111  if (nelements_root!=nelements) {
112  // Should never happen, but if it has,
113  // our best course is to abort here and now.
114  std::cerr << "Broadcast of incompatible boost::multi_array data detected on rank " << rank
115  << ".\nRoot sends " << nelements_root << " elements,"
116  << " this process expects " << nelements << " elements."
117  << "\nAborting..."
118  << std::endl;
119  MPI_Abort(MPI_COMM_WORLD, 1);
120  }
121  alps::mpi::broadcast(comm, &data(bases), nelements, root);
122  }
123 
124  } // detail::
125  } // gf::
126 } // alps::
127 
128 #endif /* ALPS_GF_MPI_BCAST_HPP_c030bec39d4b43b9a24a16b5805f542d */
index_mesh::index_type index
Definition: mesh.hpp:1247
void broadcast(C const &c, P &p, int r=0)
Definition: api.hpp:56
Header for object-oriented interface to MPI (similar to boost::mpi)
Encapsulation of an MPI communicator and some communicator-related operations.
Definition: mpi.hpp:111
int rank() const
Returns process rank in this communicator.
Definition: mpi.hpp:156
void broadcast(const communicator &comm, T *vals, std::size_t count, int root)
Broadcasts array vals of a primitive type T, length count on communicator comm with root root ...
Definition: mpi.hpp:270