--- /dev/null
+New: Add new function Utilities::MPI::compute_union(), which computes the union of given input
+vectors of all processes in a given MPI communicator.
+<br>
+(Peter Munch, 2020/04/04)
const IndexSet &indices_to_look_up,
const MPI_Comm &comm);
+ /**
+ * Compute the union of the input vectors @p vec of all processes in the
+ * MPI communicator @p comm.
+ *
+ * @note This is a collective operation. The result will available on all
+ * processes.
+ */
+ template <typename T>
+ std::vector<T>
+ compute_union(const std::vector<T> &vec, const MPI_Comm &comm);
+
#ifndef DOXYGEN
// declaration for an internal function that lives in mpi.templates.h
namespace internal
#include <deal.II/lac/sparse_matrix.h>
#include <deal.II/lac/vector.h>
+#include <set>
#include <vector>
DEAL_II_NAMESPACE_OPEN
internal::all_reduce(MPI_MIN, values, mpi_communicator, minima);
}
+
+
+ template <typename T>
+ std::vector<T>
+ compute_union(const std::vector<T> &vec, const MPI_Comm &comm)
+ {
+#ifdef DEAL_II_WITH_MPI
+ // 1) collect vector entries and create union
+ std::vector<T> result = vec;
+
+ if (this_mpi_process(comm) == 0)
+ {
+ for (unsigned int i = 1; i < n_mpi_processes(comm); i++)
+ {
+ MPI_Status status;
+ const auto ierr_1 = MPI_Probe(MPI_ANY_SOURCE,
+ internal::Tags::compute_union,
+ comm,
+ &status);
+ AssertThrowMPI(ierr_1);
+
+ int amount;
+ const auto ierr_2 =
+ MPI_Get_count(&status,
+ internal::mpi_type_id(vec.data()),
+ &amount);
+ AssertThrowMPI(ierr_2);
+
+ const unsigned int old_size = result.size();
+ result.resize(old_size + amount);
+
+ const auto ierr_3 = MPI_Recv(result.data() + old_size,
+ amount,
+ MPI_INT,
+ status.MPI_SOURCE,
+ status.MPI_TAG,
+ comm,
+ MPI_STATUS_IGNORE);
+ AssertThrowMPI(ierr_3);
+
+ std::sort(result.begin(), result.end());
+ result.erase(std::unique(result.begin(), result.end()),
+ result.end());
+ }
+ }
+ else
+ {
+ const auto ierr = MPI_Send(vec.data(),
+ vec.size(),
+ internal::mpi_type_id(vec.data()),
+ 0,
+ internal::Tags::compute_union,
+ comm);
+ AssertThrowMPI(ierr);
+ }
+
+ // 2) broadcast result
+ int size = result.size();
+ MPI_Bcast(&size, 1, MPI_INT, 0, comm);
+ result.resize(size);
+ MPI_Bcast(
+ result.data(), size, internal::mpi_type_id(vec.data()), 0, comm);
+
+ return result;
+#else
+ (void)comm;
+ return vec;
+#endif
+ }
+
} // end of namespace MPI
} // end of namespace Utilities
/// NoncontiguousPartitioner::update_values
noncontiguous_partitioner_update_values,
+ // Utilities::MPI::compute_union
+ compute_union,
+
};
} // namespace Tags
} // namespace internal
locked = false;
}
+
+ template std::vector<unsigned int>
+ compute_union(const std::vector<unsigned int> &vec, const MPI_Comm &comm);
+
#include "mpi.inst"
} // end of namespace MPI
} // end of namespace Utilities
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2020 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+// Create a distributed triangulation with multigrid levels and copy it.
+
+#include <deal.II/base/mpi.h>
+
+#include "../tests.h"
+
+using namespace dealii;
+
+void
+test(const MPI_Comm comm)
+{
+ std::vector<unsigned int> vector;
+
+ if (Utilities::MPI::this_mpi_process(comm) == 0)
+ vector = std::vector<unsigned int>{2, 5, 7};
+ else
+ vector = std::vector<unsigned int>{5, 8, 9, 10};
+
+ const auto result = Utilities::MPI::compute_union(vector, comm);
+
+ for (auto i : result)
+ deallog << i << " ";
+ deallog << std::endl;
+}
+
+int
+main(int argc, char *argv[])
+{
+ Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+ MPILogInitAll all;
+
+ const MPI_Comm comm = MPI_COMM_WORLD;
+
+ test(comm);
+}
--- /dev/null
+
+DEAL:0::2 5 7 8 9 10
+
+DEAL:1::2 5 7 8 9 10
+