From: Peter Munch Date: Sun, 26 Apr 2020 13:48:45 +0000 (+0200) Subject: Simplify parallel::TriangulationBase::update_number_cache X-Git-Tag: v9.2.0-rc1~177^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=refs%2Fpull%2F9962%2Fhead;p=dealii.git Simplify parallel::TriangulationBase::update_number_cache --- diff --git a/include/deal.II/distributed/fully_distributed_tria.h b/include/deal.II/distributed/fully_distributed_tria.h index 8e1678c4cf..6d98768c2c 100644 --- a/include/deal.II/distributed/fully_distributed_tria.h +++ b/include/deal.II/distributed/fully_distributed_tria.h @@ -242,13 +242,6 @@ namespace parallel const unsigned int coarse_cell_index) const override; private: - /** - * Override the function to update the number cache so we can fill data - * like @p level_ghost_owners. - */ - virtual void - update_number_cache() override; - /** * store the Settings. */ diff --git a/include/deal.II/distributed/shared_tria.h b/include/deal.II/distributed/shared_tria.h index cde151c46d..f4604e548a 100644 --- a/include/deal.II/distributed/shared_tria.h +++ b/include/deal.II/distributed/shared_tria.h @@ -342,13 +342,6 @@ namespace parallel with_artificial_cells() const; private: - /** - * Override the function to update the number cache so we can fill data - * like @p level_ghost_owners. - */ - virtual void - update_number_cache() override; - /** * Settings */ diff --git a/include/deal.II/distributed/tria.h b/include/deal.II/distributed/tria.h index d6d219527d..bb133907d1 100644 --- a/include/deal.II/distributed/tria.h +++ b/include/deal.II/distributed/tria.h @@ -861,13 +861,6 @@ namespace parallel private: - /** - * Override the function to update the number cache so we can fill data - * like @p level_ghost_owners. - */ - virtual void - update_number_cache() override; - /** * store the Settings. */ diff --git a/include/deal.II/distributed/tria_base.h b/include/deal.II/distributed/tria_base.h index a4276bfbb6..818c23be10 100644 --- a/include/deal.II/distributed/tria_base.h +++ b/include/deal.II/distributed/tria_base.h @@ -301,12 +301,6 @@ namespace parallel */ virtual void update_number_cache(); - - /** - * Store MPI ranks of level ghost owners of this processor on all levels. - */ - void - fill_level_ghost_owners(); }; /** diff --git a/source/distributed/fully_distributed_tria.cc b/source/distributed/fully_distributed_tria.cc index 7fec8ccd2f..bded33a0e2 100644 --- a/source/distributed/fully_distributed_tria.cc +++ b/source/distributed/fully_distributed_tria.cc @@ -173,7 +173,7 @@ namespace parallel } } - update_number_cache(); + this->update_number_cache(); } @@ -265,19 +265,6 @@ namespace parallel - template - void - Triangulation::update_number_cache() - { - parallel::Triangulation::update_number_cache(); - - if (settings & - TriangulationDescription::Settings::construct_multigrid_hierarchy) - parallel::Triangulation::fill_level_ghost_owners(); - } - - - template void Triangulation::execute_coarsening_and_refinement() diff --git a/source/distributed/shared_tria.cc b/source/distributed/shared_tria.cc index a21102ad37..e2f55d35de 100644 --- a/source/distributed/shared_tria.cc +++ b/source/distributed/shared_tria.cc @@ -419,18 +419,6 @@ namespace parallel partition(); this->update_number_cache(); } - - - - template - void - Triangulation::update_number_cache() - { - parallel::TriangulationBase::update_number_cache(); - - if (settings & construct_multigrid_hierarchy) - parallel::TriangulationBase::fill_level_ghost_owners(); - } } // namespace shared } // namespace parallel diff --git a/source/distributed/tria.cc b/source/distributed/tria.cc index 039bcf1943..b7a4d53140 100644 --- a/source/distributed/tria.cc +++ b/source/distributed/tria.cc @@ -2914,18 +2914,6 @@ namespace parallel - template - void - Triangulation::update_number_cache() - { - parallel::TriangulationBase::update_number_cache(); - - if (settings & construct_multigrid_hierarchy) - parallel::TriangulationBase::fill_level_ghost_owners(); - } - - - template typename dealii::internal::p4est::types::tree * Triangulation::init_tree( diff --git a/source/distributed/tria_base.cc b/source/distributed/tria_base.cc index e12ddedbf9..7282d17092 100644 --- a/source/distributed/tria_base.cc +++ b/source/distributed/tria_base.cc @@ -215,92 +215,90 @@ namespace parallel number_cache.n_global_levels = Utilities::MPI::max(this->n_levels(), this->mpi_communicator); - } - - - template - void - TriangulationBase::fill_level_ghost_owners() - { - number_cache.level_ghost_owners.clear(); - - // if there is nothing to do, then do nothing - if (this->n_levels() == 0) - return; - - // find level ghost owners - for (typename Triangulation::cell_iterator cell = - this->begin(); - cell != this->end(); - ++cell) - if (cell->level_subdomain_id() != numbers::artificial_subdomain_id && - cell->level_subdomain_id() != this->locally_owned_subdomain()) - this->number_cache.level_ghost_owners.insert( - cell->level_subdomain_id()); + // Store MPI ranks of level ghost owners of this processor on all levels. + if (this->is_multilevel_hierarchy_constructed() == true) + { + number_cache.level_ghost_owners.clear(); + + // if there is nothing to do, then do nothing + if (this->n_levels() == 0) + return; + + // find level ghost owners + for (typename Triangulation::cell_iterator cell = + this->begin(); + cell != this->end(); + ++cell) + if (cell->level_subdomain_id() != numbers::artificial_subdomain_id && + cell->level_subdomain_id() != this->locally_owned_subdomain()) + this->number_cache.level_ghost_owners.insert( + cell->level_subdomain_id()); # ifdef DEBUG - // Check that level_ghost_owners is symmetric by sending a message to - // everyone - { - int ierr = MPI_Barrier(this->mpi_communicator); - AssertThrowMPI(ierr); - - const int mpi_tag = Utilities::MPI::internal::Tags:: - triangulation_base_fill_level_ghost_owners; - - // important: preallocate to avoid (re)allocation: - std::vector requests( - this->number_cache.level_ghost_owners.size()); - unsigned int dummy = 0; - unsigned int req_counter = 0; - - for (std::set::iterator it = - this->number_cache.level_ghost_owners.begin(); - it != this->number_cache.level_ghost_owners.end(); - ++it, ++req_counter) - { - ierr = MPI_Isend(&dummy, - 1, - MPI_UNSIGNED, - *it, - mpi_tag, - this->mpi_communicator, - &requests[req_counter]); - AssertThrowMPI(ierr); - } - - for (std::set::iterator it = - this->number_cache.level_ghost_owners.begin(); - it != this->number_cache.level_ghost_owners.end(); - ++it) + // Check that level_ghost_owners is symmetric by sending a message to + // everyone { - unsigned int dummy; - ierr = MPI_Recv(&dummy, - 1, - MPI_UNSIGNED, - *it, - mpi_tag, - this->mpi_communicator, - MPI_STATUS_IGNORE); + int ierr = MPI_Barrier(this->mpi_communicator); AssertThrowMPI(ierr); - } - if (requests.size() > 0) - { - ierr = - MPI_Waitall(requests.size(), requests.data(), MPI_STATUSES_IGNORE); + const int mpi_tag = Utilities::MPI::internal::Tags:: + triangulation_base_fill_level_ghost_owners; + + // important: preallocate to avoid (re)allocation: + std::vector requests( + this->number_cache.level_ghost_owners.size()); + unsigned int dummy = 0; + unsigned int req_counter = 0; + + for (std::set::iterator it = + this->number_cache.level_ghost_owners.begin(); + it != this->number_cache.level_ghost_owners.end(); + ++it, ++req_counter) + { + ierr = MPI_Isend(&dummy, + 1, + MPI_UNSIGNED, + *it, + mpi_tag, + this->mpi_communicator, + &requests[req_counter]); + AssertThrowMPI(ierr); + } + + for (std::set::iterator it = + this->number_cache.level_ghost_owners.begin(); + it != this->number_cache.level_ghost_owners.end(); + ++it) + { + unsigned int dummy; + ierr = MPI_Recv(&dummy, + 1, + MPI_UNSIGNED, + *it, + mpi_tag, + this->mpi_communicator, + MPI_STATUS_IGNORE); + AssertThrowMPI(ierr); + } + + if (requests.size() > 0) + { + ierr = MPI_Waitall(requests.size(), + requests.data(), + MPI_STATUSES_IGNORE); + AssertThrowMPI(ierr); + } + + ierr = MPI_Barrier(this->mpi_communicator); AssertThrowMPI(ierr); } - - ierr = MPI_Barrier(this->mpi_communicator); - AssertThrowMPI(ierr); - } # endif - Assert(this->number_cache.level_ghost_owners.size() < - Utilities::MPI::n_mpi_processes(this->mpi_communicator), - ExcInternalError()); + Assert(this->number_cache.level_ghost_owners.size() < + Utilities::MPI::n_mpi_processes(this->mpi_communicator), + ExcInternalError()); + } } #else @@ -312,13 +310,6 @@ namespace parallel Assert(false, ExcNotImplemented()); } - template - void - TriangulationBase::fill_level_ghost_owners() - { - Assert(false, ExcNotImplemented()); - } - #endif template