const std::vector<types::global_dof_index> &indices_want,
const MPI_Comm & communicator)
{
+ // step 0) clean vectors from numbers::invalid_dof_index (indicating
+ // padding)
+ std::vector<types::global_dof_index> 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<types::global_dof_index> 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);
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];
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];
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// 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 <deal.II/base/mpi.h>
+#include <deal.II/base/mpi_noncontiguous_partitioner.h>
+
+#include "../tests.h"
+
+using namespace dealii;
+
+void
+test(const MPI_Comm comm,
+ std::vector<types::global_dof_index> index_set_has,
+ std::vector<types::global_dof_index> index_set_want)
+{
+ Utilities::MPI::NoncontiguousPartitioner<double> vector;
+ vector.reinit(index_set_has, index_set_want, comm);
+
+ AlignedVector<double> src(index_set_has.size(), 0);
+ AlignedVector<double> 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<int>(src[i]) << " ";
+ deallog << std::endl;
+ for (size_t i = 0; i < dst.size(); i++)
+ deallog << static_cast<int>(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();
+ }
+}