• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
1 // Copyright (C) 2011 Júlio Hoffimann.
2 
3 // Use, modification and distribution is subject to the Boost Software
4 // License, Version 1.0. (See accompanying file LICENSE_1_0.txt or copy at
5 // http://www.boost.org/LICENSE_1_0.txt)
6 
7 // Message Passing Interface 1.1 -- Section 4.6. Scatterv
8 #ifndef BOOST_MPI_SCATTERV_HPP
9 #define BOOST_MPI_SCATTERV_HPP
10 
11 #include <boost/scoped_array.hpp>
12 #include <boost/mpi/collectives/scatter.hpp>
13 #include <boost/mpi/detail/offsets.hpp>
14 #include <boost/mpi/detail/antiques.hpp>
15 
16 namespace boost { namespace mpi {
17 
18 namespace detail {
19 
20 //////////////////////////////////////////////
21 /// Implementation for MPI primitive types ///
22 //////////////////////////////////////////////
23 
24 // We're scattering from the root for a type that has an associated MPI
25 // datatype, so we'll use MPI_Scatterv to do all of the work.
26 template<typename T>
27 void
scatterv_impl(const communicator & comm,const T * in_values,T * out_values,int out_size,const int * sizes,const int * displs,int root,mpl::true_)28 scatterv_impl(const communicator& comm, const T* in_values, T* out_values, int out_size,
29               const int* sizes, const int* displs, int root, mpl::true_)
30 {
31   assert(!sizes || out_size == sizes[comm.rank()]);
32   assert(bool(sizes) == bool(in_values));
33 
34   scoped_array<int> new_offsets_mem(make_offsets(comm, sizes, displs, root));
35   if (new_offsets_mem) displs = new_offsets_mem.get();
36   MPI_Datatype type = get_mpi_datatype<T>(*in_values);
37   BOOST_MPI_CHECK_RESULT(MPI_Scatterv,
38                          (const_cast<T*>(in_values), const_cast<int*>(sizes),
39                           const_cast<int*>(displs), type,
40                           out_values, out_size, type, root, comm));
41 }
42 
43 // We're scattering from a non-root for a type that has an associated MPI
44 // datatype, so we'll use MPI_Scatterv to do all of the work.
45 template<typename T>
46 void
scatterv_impl(const communicator & comm,T * out_values,int out_size,int root,mpl::true_ is_mpi_type)47 scatterv_impl(const communicator& comm, T* out_values, int out_size, int root,
48               mpl::true_ is_mpi_type)
49 {
50   scatterv_impl(comm, (T const*)0, out_values, out_size,
51                 (const int*)0, (const int*)0, root, is_mpi_type);
52 }
53 
54 //////////////////////////////////////////////////
55 /// Implementation for non MPI primitive types ///
56 //////////////////////////////////////////////////
57 
58 // We're scattering from the root for a type that does not have an
59 // associated MPI datatype, so we'll need to serialize it.
60 template<typename T>
61 void
scatterv_impl(const communicator & comm,const T * in_values,T * out_values,int out_size,int const * sizes,int const * displs,int root,mpl::false_)62 scatterv_impl(const communicator& comm, const T* in_values, T* out_values, int out_size,
63               int const* sizes, int const* displs, int root, mpl::false_)
64 {
65   packed_oarchive::buffer_type sendbuf;
66   bool is_root = comm.rank() == root;
67   int nproc = comm.size();
68   std::vector<int> archsizes;
69   if (is_root) {
70     assert(out_size == sizes[comm.rank()]);
71     archsizes.resize(nproc);
72     std::vector<int> skipped;
73     if (displs) {
74       skipped.resize(nproc);
75       offsets2skipped(sizes, displs, c_data(skipped), nproc);
76       displs = c_data(skipped);
77     }
78     fill_scatter_sendbuf(comm, in_values, sizes, (int const*)0, sendbuf, archsizes);
79   }
80   dispatch_scatter_sendbuf(comm, sendbuf, archsizes, (T const*)0, out_values, out_size, root);
81 }
82 
83 // We're scattering to a non-root for a type that does not have an
84 // associated MPI datatype. input data not needed.
85 // it.
86 template<typename T>
87 void
scatterv_impl(const communicator & comm,T * out_values,int n,int root,mpl::false_ isnt_mpi_type)88 scatterv_impl(const communicator& comm, T* out_values, int n, int root,
89               mpl::false_ isnt_mpi_type)
90 {
91   assert(root != comm.rank());
92   scatterv_impl(comm, (T const*)0, out_values, n, (int const*)0, (int const*)0, root, isnt_mpi_type);
93 }
94 
95 } // end namespace detail
96 
97 template<typename T>
98 void
scatterv(const communicator & comm,const T * in_values,const std::vector<int> & sizes,const std::vector<int> & displs,T * out_values,int out_size,int root)99 scatterv(const communicator& comm, const T* in_values,
100          const std::vector<int>& sizes, const std::vector<int>& displs,
101          T* out_values, int out_size, int root)
102 {
103   using detail::c_data;
104   detail::scatterv_impl(comm, in_values, out_values, out_size, c_data(sizes), c_data(displs),
105                 root, is_mpi_datatype<T>());
106 }
107 
108 template<typename T>
109 void
scatterv(const communicator & comm,const std::vector<T> & in_values,const std::vector<int> & sizes,const std::vector<int> & displs,T * out_values,int out_size,int root)110 scatterv(const communicator& comm, const std::vector<T>& in_values,
111          const std::vector<int>& sizes, const std::vector<int>& displs,
112          T* out_values, int out_size, int root)
113 {
114   using detail::c_data;
115   ::boost::mpi::scatterv(comm, c_data(in_values), sizes, displs,
116                          out_values, out_size, root);
117 }
118 
119 template<typename T>
scatterv(const communicator & comm,T * out_values,int out_size,int root)120 void scatterv(const communicator& comm, T* out_values, int out_size, int root)
121 {
122   BOOST_ASSERT(comm.rank() != root);
123   detail::scatterv_impl(comm, out_values, out_size, root, is_mpi_datatype<T>());
124 }
125 
126 ///////////////////////
127 // common use versions
128 ///////////////////////
129 template<typename T>
130 void
scatterv(const communicator & comm,const T * in_values,const std::vector<int> & sizes,T * out_values,int root)131 scatterv(const communicator& comm, const T* in_values,
132          const std::vector<int>& sizes, T* out_values, int root)
133 {
134   using detail::c_data;
135   detail::scatterv_impl(comm, in_values, out_values, sizes[comm.rank()],
136                         c_data(sizes), (int const*)0,
137                         root, is_mpi_datatype<T>());
138 }
139 
140 template<typename T>
141 void
scatterv(const communicator & comm,const std::vector<T> & in_values,const std::vector<int> & sizes,T * out_values,int root)142 scatterv(const communicator& comm, const std::vector<T>& in_values,
143          const std::vector<int>& sizes, T* out_values, int root)
144 {
145   ::boost::mpi::scatterv(comm, detail::c_data(in_values), sizes, out_values, root);
146 }
147 
148 template<typename T>
149 void
scatterv(const communicator & comm,const T * in_values,T * out_values,int n,int root)150 scatterv(const communicator& comm, const T* in_values,
151          T* out_values, int n, int root)
152 {
153   detail::scatterv_impl(comm, in_values, out_values, n, (int const*)0, (int const*)0,
154                 root, is_mpi_datatype<T>());
155 }
156 
157 template<typename T>
158 void
scatterv(const communicator & comm,const std::vector<T> & in_values,T * out_values,int out_size,int root)159 scatterv(const communicator& comm, const std::vector<T>& in_values,
160          T* out_values, int out_size, int root)
161 {
162   ::boost::mpi::scatterv(comm, detail::c_data(in_values), out_values, out_size, root);
163 }
164 
165 } } // end namespace boost::mpi
166 
167 #endif // BOOST_MPI_SCATTERV_HPP
168