From: Peter Munch Date: Sun, 8 Mar 2020 21:56:30 +0000 (+0100) Subject: Enable Utilities::MPI::NoncontiguousPartitioner to handle padding X-Git-Tag: v9.2.0-rc1~407^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=336eab7efbcde9f61197a9672ab3ff997eb1910b;p=dealii.git Enable Utilities::MPI::NoncontiguousPartitioner to handle padding --- diff --git a/include/deal.II/base/mpi_noncontiguous_partitioner.templates.h b/include/deal.II/base/mpi_noncontiguous_partitioner.templates.h index 5487bc823d..085c66772a 100644 --- a/include/deal.II/base/mpi_noncontiguous_partitioner.templates.h +++ b/include/deal.II/base/mpi_noncontiguous_partitioner.templates.h @@ -183,26 +183,49 @@ namespace Utilities const std::vector &indices_want, const MPI_Comm & communicator) { + // step 0) clean vectors from numbers::invalid_dof_index (indicating + // padding) + std::vector indices_has_clean; + indices_has_clean.reserve(indices_has.size()); + + for (const auto &i : indices_has) + if (i != numbers::invalid_dof_index) + indices_has_clean.push_back(i); + + std::vector indices_want_clean; + indices_want_clean.reserve(indices_has.size()); + + for (const auto &i : indices_want) + if (i != numbers::invalid_dof_index) + indices_want_clean.push_back(i); + // step 0) determine "number of degrees of freedom" needed for IndexSet - const types::global_dof_index n_dofs = Utilities::MPI::max( - std::max( - [&indices_has]() { - const auto it = max_element(indices_has.begin(), indices_has.end()); - return it == indices_has.end() ? 0 : (*it + 1); - }(), - [&indices_want]() { - const auto it = - max_element(indices_want.begin(), indices_want.end()); - return it == indices_want.end() ? 0 : (*it + 1); - }()), - communicator); + const types::global_dof_index local_n_dofs_has = + indices_has_clean.empty() ? + 0 : + (*std::max_element(indices_has_clean.begin(), + indices_has_clean.end()) + + 1); + + const types::global_dof_index local_n_dofs_want = + indices_want_clean.empty() ? + 0 : + (*std::max_element(indices_want_clean.begin(), + indices_want_clean.end()) + + 1); + + const types::global_dof_index n_dofs = + Utilities::MPI::max(std::max(local_n_dofs_has, local_n_dofs_want), + communicator); // step 1) convert vectors to indexsets (sorted!) IndexSet index_set_has(n_dofs); - index_set_has.add_indices(indices_has.begin(), indices_has.end()); + index_set_has.add_indices(indices_has_clean.begin(), + indices_has_clean.end()); IndexSet index_set_want(n_dofs); - index_set_want.add_indices(indices_want.begin(), indices_want.end()); + index_set_want.add_indices(indices_want_clean.begin(), + indices_want_clean.end()); // step 2) setup internal data structures with indexset this->reinit(index_set_has, index_set_want, communicator); @@ -214,7 +237,8 @@ namespace Utilities index_set_has.n_elements()); for (types::global_dof_index i = 0; i < indices_has.size(); i++) - temp_map_send[index_set_has.index_within_set(indices_has[i])] = i; + if (indices_has[i] != numbers::invalid_dof_index) + temp_map_send[index_set_has.index_within_set(indices_has[i])] = i; for (auto &i : send_indices) i = temp_map_send[i]; @@ -225,7 +249,8 @@ namespace Utilities index_set_want.n_elements()); for (types::global_dof_index i = 0; i < indices_want.size(); i++) - temp_map_recv[index_set_want.index_within_set(indices_want[i])] = i; + if (indices_want[i] != numbers::invalid_dof_index) + temp_map_recv[index_set_want.index_within_set(indices_want[i])] = i; for (auto &i : recv_indices) i = temp_map_recv[i]; diff --git a/tests/base/mpi_noncontiguous_partitioner_03.cc b/tests/base/mpi_noncontiguous_partitioner_03.cc new file mode 100644 index 0000000000..1ee78adb01 --- /dev/null +++ b/tests/base/mpi_noncontiguous_partitioner_03.cc @@ -0,0 +1,89 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2019 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. +// +// --------------------------------------------------------------------- + + +// Test Utilities::MPI::NoncontiguousPartitioner for padding. + +#include +#include + +#include "../tests.h" + +using namespace dealii; + +void +test(const MPI_Comm comm, + std::vector index_set_has, + std::vector index_set_want) +{ + Utilities::MPI::NoncontiguousPartitioner vector; + vector.reinit(index_set_has, index_set_want, comm); + + AlignedVector src(index_set_has.size(), 0); + AlignedVector dst(index_set_want.size(), 0); + + for (unsigned int i = 0; i < index_set_has.size(); i++) + src[i] = Utilities::MPI::this_mpi_process(comm) * 100 + i; + + vector.update_values(dst, src); + + for (size_t i = 0; i < src.size(); i++) + deallog << static_cast(src[i]) << " "; + deallog << std::endl; + for (size_t i = 0; i < dst.size(); i++) + deallog << static_cast(dst[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; + + const unsigned int rank = Utilities::MPI::this_mpi_process(comm); + + { + deallog.push("padding-non"); + + if (rank == 0) + test(comm, {0, 1, 2, 3}, {4, 5, 6, 7}); + else + test(comm, {4, 5, 6, 7}, {0, 1, 2, 3}); + deallog.pop(); + } + + { + deallog.push("padding-src"); + + if (rank == 0) + test(comm, {0, 1, numbers::invalid_dof_index, 2, 3}, {4, 5, 6, 7}); + else + test(comm, {4, 5, 6, 7}, {0, 1, 2, 3}); + deallog.pop(); + } + + { + deallog.push("padding-dst"); + + if (rank == 0) + test(comm, {0, 1, 2, 3}, {4, 5, numbers::invalid_dof_index, 6, 7}); + else + test(comm, {4, 5, 6, 7}, {0, 1, 2, 3}); + deallog.pop(); + } +} diff --git a/tests/base/mpi_noncontiguous_partitioner_03.mpirun=2.output b/tests/base/mpi_noncontiguous_partitioner_03.mpirun=2.output new file mode 100644 index 0000000000..bfef1b6851 --- /dev/null +++ b/tests/base/mpi_noncontiguous_partitioner_03.mpirun=2.output @@ -0,0 +1,15 @@ + +DEAL:0:padding-non::0 1 2 3 +DEAL:0:padding-non::100 101 102 103 +DEAL:0:padding-src::0 1 2 3 4 +DEAL:0:padding-src::100 101 102 103 +DEAL:0:padding-dst::0 1 2 3 +DEAL:0:padding-dst::100 101 0 102 103 + +DEAL:1:padding-non::100 101 102 103 +DEAL:1:padding-non::0 1 2 3 +DEAL:1:padding-src::100 101 102 103 +DEAL:1:padding-src::0 1 3 4 +DEAL:1:padding-dst::100 101 102 103 +DEAL:1:padding-dst::0 1 2 3 +