-New: Add new function Utilities::MPI::compute_union(), which computes the union of given input
-vectors of all processes in a given MPI communicator.
+New: Add new function Utilities::MPI::compute_set_union(), which computes the union of given input
+sets and vectors of all processes in a given MPI communicator.
<br>
(Peter Munch, 2020/04/04)
*/
template <typename T>
std::vector<T>
- compute_union(const std::vector<T> &vec, const MPI_Comm &comm);
+ compute_set_union(const std::vector<T> &vec, const MPI_Comm &comm);
+
+ /**
+ * The same as above but for std::set.
+ */
+ template <typename T>
+ std::set<T>
+ compute_set_union(const std::set<T> &set, const MPI_Comm &comm);
#ifndef DOXYGEN
// declaration for an internal function that lives in mpi.templates.h
template <typename T>
std::vector<T>
- compute_union(const std::vector<T> &vec, const MPI_Comm &comm)
+ compute_set_union(const std::vector<T> &vec, const MPI_Comm &comm)
{
#ifdef DEAL_II_WITH_MPI
// 1) collect vector entries and create union
const auto ierr_3 = MPI_Recv(result.data() + old_size,
amount,
- MPI_INT,
+ internal::mpi_type_id(vec.data()),
status.MPI_SOURCE,
status.MPI_TAG,
comm,
#endif
}
+
+
+ template <typename T>
+ std::set<T>
+ compute_set_union(const std::set<T> &set_in, const MPI_Comm &comm)
+ {
+ // convert vector to set
+ std::vector<T> vector_in(set_in.begin(), set_in.end());
+
+ // perform operation to vector
+ const std::vector<T> vector_out = compute_set_union(vector_in, comm);
+
+ // convert vector to set
+ return std::set<T>(vector_out.begin(), vector_out.end());
+ }
+
} // end of namespace MPI
} // end of namespace Utilities
template std::vector<unsigned int>
- compute_union(const std::vector<unsigned int> &vec, const MPI_Comm &comm);
+ compute_set_union(const std::vector<unsigned int> &vec,
+ const MPI_Comm & comm);
+
+
+ template std::set<unsigned int>
+ compute_set_union(const std::set<unsigned int> &set, const MPI_Comm &comm);
#include "mpi.inst"
} // end of namespace MPI
// ---------------------------------------------------------------------
-// Create a distributed triangulation with multigrid levels and copy it.
+// Test the function compute_set_union() for set::vector and std::set.
#include <deal.II/base/mpi.h>
else
vector = std::vector<unsigned int>{5, 8, 9, 10};
- const auto result = Utilities::MPI::compute_union(vector, comm);
+ // test function for vector
+ {
+ const auto result = Utilities::MPI::compute_set_union(vector, comm);
- for (auto i : result)
- deallog << i << " ";
- deallog << std::endl;
+ for (auto i : result)
+ deallog << i << " ";
+ deallog << std::endl;
+ }
+
+ // test function for set
+ {
+ const auto result = Utilities::MPI::compute_set_union(
+ std::set<unsigned int>(vector.begin(), vector.end()), comm);
+
+ for (auto i : result)
+ deallog << i << " ";
+ deallog << std::endl;
+ }
}
int
+DEAL:0::2 5 7 8 9 10
DEAL:0::2 5 7 8 9 10
DEAL:1::2 5 7 8 9 10
+DEAL:1::2 5 7 8 9 10