From 6f2f9c4a2f192f1587734c7cbe647a5169da4ccb Mon Sep 17 00:00:00 2001 From: Timo Heister Date: Tue, 22 Oct 2019 22:01:14 +0200 Subject: [PATCH] add DuplicatedCommunicator --- include/deal.II/base/mpi.h | 55 +++++++++++ tests/mpi/communicator_duplicated.cc | 91 +++++++++++++++++++ .../communicator_duplicated.mpirun=2.output | 2 + 3 files changed, 148 insertions(+) create mode 100644 tests/mpi/communicator_duplicated.cc create mode 100644 tests/mpi/communicator_duplicated.mpirun=2.output diff --git a/include/deal.II/base/mpi.h b/include/deal.II/base/mpi.h index 2c03fabdf3..5f749a40e8 100644 --- a/include/deal.II/base/mpi.h +++ b/include/deal.II/base/mpi.h @@ -211,6 +211,61 @@ namespace Utilities void free_communicator(MPI_Comm &mpi_communicator); + /** + * Helper class to automatically duplicate and free an MPI + * @ref GlossMPICommunicator "communicator". + * + * This class duplicates the communicator given in the constructor + * using duplicate_communicator() and frees it automatically when + * this object gets destroyed by calling free_communicator(). You + * can access the wrapped communicator using operator*. + * + * This class exists to easily allow duplicating communicators without + * having to worry when and how to free it after usage. + */ + class DuplicatedCommunicator + { + public: + /** + * Create a duplicate of the given @p communicator. + */ + explicit DuplicatedCommunicator(const MPI_Comm &communicator) + : comm (duplicate_communicator(communicator)) + {} + + /** + * The destructor will free the communicator automatically. + */ + ~DuplicatedCommunicator() + { + free_communicator(comm); + } + + /** + * Access the stored communicator. + */ + const MPI_Comm &operator*() const + { + return comm; + } + + /** + * Do not allow making copies. + */ + DuplicatedCommunicator(const DuplicatedCommunicator &) = delete; + + /** + * Do not allow assignment of this class. + */ + DuplicatedCommunicator& operator=(const DuplicatedCommunicator &) = delete; + + private: + /** + * The communicator of course. + */ + MPI_Comm comm; + }; + /** * If @p comm is an intracommunicator, this function returns a new * communicator @p newcomm with communication group defined by the diff --git a/tests/mpi/communicator_duplicated.cc b/tests/mpi/communicator_duplicated.cc new file mode 100644 index 0000000000..10b408d00f --- /dev/null +++ b/tests/mpi/communicator_duplicated.cc @@ -0,0 +1,91 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2017 - 2018 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. +// +// --------------------------------------------------------------------- + + +// check MPI::duplicate_communicator, free_communicator, and +// DuplicatedCommunicator + + +#include + +#include "../tests.h" + + +void test(const MPI_Comm comm) +{ + int tag = 12345; + + const auto my_rank = Utilities::MPI::this_mpi_process(comm); + + MPI_Comm comm2 = Utilities::MPI::duplicate_communicator(comm); + Utilities::MPI::DuplicatedCommunicator pcomm3(comm); + + if (my_rank==1) + { + int value[3]={1,2,3}; + int dest = 0; + + MPI_Send(&value[0], + 1, + MPI_UNSIGNED, + dest, + tag, + comm); + MPI_Send(&value[1], + 1, + MPI_UNSIGNED, + dest, + tag, + comm2); + MPI_Send(&value[2], + 1, + MPI_UNSIGNED, + dest, + tag, + *pcomm3); + } + + if (my_rank == 0) + { + int value[3]; + int src = 1; + + // receive in reverse order, if duplication of communicators worked, + // these won't be mixed up! + MPI_Recv(&value[2], 1, MPI_UNSIGNED, + src, tag, *pcomm3, MPI_STATUS_IGNORE); + MPI_Recv(&value[1], 1, MPI_UNSIGNED, + src, tag, comm2, MPI_STATUS_IGNORE); + MPI_Recv(&value[0], 1, MPI_UNSIGNED, + src, tag, comm, MPI_STATUS_IGNORE); + + deallog << value[0] << " " << value[1] << " " << value[2] << std::endl; + } + + Utilities::MPI::free_communicator(comm2); +} + + + +int +main(int argc, char **argv) +{ + Utilities::MPI::MPI_InitFinalize mpi_initialization( + argc, argv, testing_max_num_threads()); + + mpi_initlog(); + + test(MPI_COMM_WORLD); +} diff --git a/tests/mpi/communicator_duplicated.mpirun=2.output b/tests/mpi/communicator_duplicated.mpirun=2.output new file mode 100644 index 0000000000..ed62a0ef21 --- /dev/null +++ b/tests/mpi/communicator_duplicated.mpirun=2.output @@ -0,0 +1,2 @@ + +DEAL::1 2 3 -- 2.39.5