From 2c05303e11aafbfc1bb44d21669c96922c73c10b Mon Sep 17 00:00:00 2001 From: Peter Munch Date: Sat, 9 Apr 2022 11:44:03 +0200 Subject: [PATCH] RPE: refactor determination of unique mapping and all points found --- .../base/mpi_remote_point_evaluation.h | 5 + source/base/mpi_remote_point_evaluation.cc | 60 ++++-- tests/remote_point_evaluation/mapping_04.cc | 197 ++++++++++++++++++ ...mapping_04.with_p4est=true.mpirun=2.output | 3 + 4 files changed, 252 insertions(+), 13 deletions(-) create mode 100644 tests/remote_point_evaluation/mapping_04.cc create mode 100644 tests/remote_point_evaluation/mapping_04.with_p4est=true.mpirun=2.output diff --git a/include/deal.II/base/mpi_remote_point_evaluation.h b/include/deal.II/base/mpi_remote_point_evaluation.h index cf1429a536..51a0a75778 100644 --- a/include/deal.II/base/mpi_remote_point_evaluation.h +++ b/include/deal.II/base/mpi_remote_point_evaluation.h @@ -242,6 +242,11 @@ namespace Utilities */ bool unique_mapping; + /** + * Cache if all points passed in during reinit() have been found. + */ + bool all_points_found_flag; + /** * Since for each point multiple or no results can be available, the * pointers in this vector indicate the first and last entry associated diff --git a/source/base/mpi_remote_point_evaluation.cc b/source/base/mpi_remote_point_evaluation.cc index dd1907609a..671cd9e467 100644 --- a/source/base/mpi_remote_point_evaluation.cc +++ b/source/base/mpi_remote_point_evaluation.cc @@ -130,14 +130,52 @@ namespace Utilities this->point_ptrs[std::get<1>(data.recv_components[i]) + 1]++; } - unique_mapping = true; + std::tuple n_owning_processes_default{ + numbers::invalid_unsigned_int, 0}; + std::tuple n_owning_processes_local = + n_owning_processes_default; + for (unsigned int i = 0; i < points.size(); ++i) { - if (unique_mapping && this->point_ptrs[i + 1] != 1) - unique_mapping = false; + std::get<0>(n_owning_processes_local) = + std::min(std::get<0>(n_owning_processes_local), + this->point_ptrs[i + 1]); + std::get<1>(n_owning_processes_local) = + std::max(std::get<1>(n_owning_processes_local), + this->point_ptrs[i + 1]); + this->point_ptrs[i + 1] += this->point_ptrs[i]; } + const auto n_owning_processes_global = + Utilities::MPI::all_reduce>( + n_owning_processes_local, + tria.get_communicator(), + [&](const auto &a, + const auto &b) -> std::tuple { + if (a == n_owning_processes_default) + return b; + + if (b == n_owning_processes_default) + return a; + + return std::tuple{ + std::min(std::get<0>(a), std::get<0>(b)), + std::max(std::get<1>(a), std::get<1>(b))}; + }); + + if (n_owning_processes_global == n_owning_processes_default) + { + unique_mapping = true; + all_points_found_flag = true; + } + else + { + unique_mapping = (std::get<0>(n_owning_processes_global) == 1) && + (std::get<1>(n_owning_processes_global) == 1); + all_points_found_flag = std::get<0>(n_owning_processes_global) > 0; + } + Assert(enforce_unique_mapping == false || unique_mapping, ExcInternalError()); @@ -189,15 +227,7 @@ namespace Utilities bool RemotePointEvaluation::all_points_found() const { - if (is_map_unique()) - return true; - - if (point_ptrs.size() > 0) - for (unsigned int i = 0; i < point_ptrs.size() - 1; ++i) - if (point_found(i) == false) - return false; - - return true; + return all_points_found_flag; } @@ -208,7 +238,11 @@ namespace Utilities const unsigned int i) const { AssertIndexRange(i, point_ptrs.size() - 1); - return (point_ptrs[i + 1] - point_ptrs[i]) > 0; + + if (all_points_found_flag) + return true; + else + return (point_ptrs[i + 1] - point_ptrs[i]) > 0; } diff --git a/tests/remote_point_evaluation/mapping_04.cc b/tests/remote_point_evaluation/mapping_04.cc new file mode 100644 index 0000000000..99c9d261dc --- /dev/null +++ b/tests/remote_point_evaluation/mapping_04.cc @@ -0,0 +1,197 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2021 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::RemotePointEvaluation if points have been found +// and if the mapping is unique. + +#include + +#include + +#include +#include + +#include + +#include + +#include + +#include "../tests.h" + +using namespace dealii; + +template +void +do_test(const bool enforce_unique_map, + const std::vector> &evaluation_points, + const std::pair & expected_result) +{ + parallel::distributed::Triangulation tria(MPI_COMM_WORLD); + GridGenerator::subdivided_hyper_rectangle(tria, {2, 1}, {0, 0}, {2, 1}); + + MappingQ1 mapping; + + Utilities::MPI::RemotePointEvaluation eval(1.e-6, enforce_unique_map); + + eval.reinit(evaluation_points, tria, mapping); + + Assert(eval.is_map_unique() == expected_result.first, ExcInternalError()); + Assert(eval.all_points_found() == expected_result.second, ExcInternalError()); +} + + +template +void +test() +{ + const auto my_rank = Utilities::MPI::this_mpi_process(MPI_COMM_WORLD); + + { + // (-, -) + std::vector> evaluation_points; + + do_test(false, evaluation_points, {true, true}); + } + + { + // (not found, -) + std::vector> evaluation_points; + + if (my_rank == 0) + evaluation_points.emplace_back(2.5, 0.5); + + do_test(false, evaluation_points, {false, false}); + } + + { + // (unique, -) + std::vector> evaluation_points; + + if (my_rank == 0) + evaluation_points.emplace_back(1.5, 0.5); + + do_test(false, evaluation_points, {true, true}); + } + + { + // (non-unique, -) + std::vector> evaluation_points; + + if (my_rank == 0) + evaluation_points.emplace_back(1.0, 0.5); + + do_test(false, evaluation_points, {false, true}); + } + + { + // (enforced-unique, -) + std::vector> evaluation_points; + + if (my_rank == 0) + evaluation_points.emplace_back(1.0, 0.5); + + do_test(true, evaluation_points, {true, true}); + } + + { + // (-, not found) + std::vector> evaluation_points; + + if (my_rank == 1) + evaluation_points.emplace_back(2.5, 0.5); + + do_test(false, evaluation_points, {false, false}); + } + + { + // (-, unique) + std::vector> evaluation_points; + + if (my_rank == 1) + evaluation_points.emplace_back(0.5, 0.5); + + do_test(false, evaluation_points, {true, true}); + } + + { + // (-, non-unique) + std::vector> evaluation_points; + + if (my_rank == 1) + evaluation_points.emplace_back(1.0, 0.5); + + do_test(false, evaluation_points, {false, true}); + } + + { + // (-, enforced-unique) + std::vector> evaluation_points; + + if (my_rank == 1) + evaluation_points.emplace_back(1.0, 0.5); + + do_test(true, evaluation_points, {true, true}); + } + + { + // (unique, not found) + std::vector> evaluation_points; + + if (my_rank == 0) + evaluation_points.emplace_back(1.5, 0.5); + if (my_rank == 1) + evaluation_points.emplace_back(2.5, 0.5); + + do_test(false, evaluation_points, {false, false}); + } + + { + // (unique, not unique) + std::vector> evaluation_points; + + if (my_rank == 0) + evaluation_points.emplace_back(1.5, 0.5); + if (my_rank == 1) + evaluation_points.emplace_back(1.0, 0.5); + + do_test(false, evaluation_points, {false, true}); + } + + { + // (unique, enforced unique) + std::vector> evaluation_points; + + if (my_rank == 0) + evaluation_points.emplace_back(1.5, 0.5); + if (my_rank == 1) + evaluation_points.emplace_back(1.0, 0.5); + + do_test(true, evaluation_points, {true, true}); + } +} + + + +int +main(int argc, char **argv) +{ + Utilities::MPI::MPI_InitFinalize mpi(argc, argv, 1); + MPILogInitAll all; + + AssertDimension(Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD), 2); + + test<2>(); +} diff --git a/tests/remote_point_evaluation/mapping_04.with_p4est=true.mpirun=2.output b/tests/remote_point_evaluation/mapping_04.with_p4est=true.mpirun=2.output new file mode 100644 index 0000000000..b28b04f643 --- /dev/null +++ b/tests/remote_point_evaluation/mapping_04.with_p4est=true.mpirun=2.output @@ -0,0 +1,3 @@ + + + -- 2.39.5