From: Peter Munch Date: Sat, 26 Jun 2021 15:30:16 +0000 (+0200) Subject: Add Triangulation::n_global_coarse_cells() X-Git-Tag: v9.4.0-rc1~1188^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=refs%2Fpull%2F12502%2Fhead;p=dealii.git Add Triangulation::n_global_coarse_cells() --- diff --git a/doc/news/changes/minor/20210626Munch b/doc/news/changes/minor/20210626Munch new file mode 100644 index 0000000000..c57fdc5e46 --- /dev/null +++ b/doc/news/changes/minor/20210626Munch @@ -0,0 +1,6 @@ +New: The new function Triangulation::n_global_coarse_cells() allows to +query the global number of coarse cells. This function is particularly useful +for parallel::fullydistributed::Triangulation, since in this case the coarse +cells might differ on each process. +
+(Peter Munch, 2021/06/26) diff --git a/include/deal.II/distributed/fully_distributed_tria.h b/include/deal.II/distributed/fully_distributed_tria.h index 27414c756e..326dec1120 100644 --- a/include/deal.II/distributed/fully_distributed_tria.h +++ b/include/deal.II/distributed/fully_distributed_tria.h @@ -278,6 +278,9 @@ namespace parallel virtual void update_cell_relations() override; + 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 6e34a49e7a..076a4e6e06 100644 --- a/include/deal.II/distributed/tria_base.h +++ b/include/deal.II/distributed/tria_base.h @@ -289,6 +289,9 @@ namespace parallel communicate_locally_moved_vertices( const std::vector &vertex_locally_moved); + virtual types::coarse_cell_id + n_global_coarse_cells() const override; + protected: /** * MPI communicator to be used for the triangulation. We create a unique @@ -318,22 +321,31 @@ namespace parallel * Number of locally owned active cells of this MPI rank. */ unsigned int n_locally_owned_active_cells; + /** * The total number of active cells (sum of @p * n_locally_owned_active_cells). */ types::global_cell_index n_global_active_cells; + + /** + * Number of global coarse cells. + */ + types::coarse_cell_id number_of_global_coarse_cells; + /** * The global number of levels computed as the maximum number of levels * taken over all MPI ranks, so n_levels()<=n_global_levels = * max(n_levels() on proc i). */ unsigned int n_global_levels; + /** * A set containing the subdomain_id (MPI rank) of the owners of the * ghost cells on this processor. */ std::set ghost_owners; + /** * A set containing the MPI ranks of the owners of the level ghost cells * on this processor (for all levels). diff --git a/include/deal.II/grid/tria.h b/include/deal.II/grid/tria.h index b954f00508..265b70d0ae 100644 --- a/include/deal.II/grid/tria.h +++ b/include/deal.II/grid/tria.h @@ -3082,6 +3082,13 @@ public: unsigned int n_active_cells(const unsigned int level) const; + /** + * Return the total number of coarse cells. If the coarse mesh is replicated + * on each process, this simply returns n_cells(0). + */ + virtual types::coarse_cell_id + n_global_coarse_cells() const; + /** * Return the total number of used faces, active or not. In 2D, the result * equals n_lines(), in 3D it equals n_quads(), while in 1D it equals diff --git a/include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h b/include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h index 58ef8e9192..cbb38d2de2 100644 --- a/include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h +++ b/include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h @@ -551,12 +551,14 @@ namespace internal , mg_level_fine(mg_level_fine) , communicator( mesh_fine.get_communicator() /*TODO: fix for different comms*/) - , cell_id_translator(n_coarse_cells(mesh_fine), - n_global_levels(mesh_fine)) + , cell_id_translator( + mesh_fine.get_triangulation().n_global_coarse_cells(), + mesh_fine.get_triangulation().n_global_levels()) { - AssertDimension(n_coarse_cells(mesh_fine), n_coarse_cells(mesh_coarse)); - AssertIndexRange(n_global_levels(mesh_coarse), - n_global_levels(mesh_fine) + 1); + AssertDimension(mesh_fine.get_triangulation().n_global_coarse_cells(), + mesh_coarse.get_triangulation().n_global_coarse_cells()); + AssertIndexRange(mesh_coarse.get_triangulation().n_global_levels(), + mesh_fine.get_triangulation().n_global_levels() + 1); } void @@ -969,25 +971,6 @@ namespace internal std::map>> map; - - static unsigned int - n_coarse_cells(const DoFHandler &mesh) - { - types::coarse_cell_id n_coarse_cells = 0; - - for (const auto &cell : mesh.get_triangulation().active_cell_iterators()) - if (!cell->is_artificial()) - n_coarse_cells = - std::max(n_coarse_cells, cell->id().get_coarse_cell_id()); - - return Utilities::MPI::max(n_coarse_cells, mesh.get_communicator()) + 1; - } - - static unsigned int - n_global_levels(const DoFHandler &mesh) - { - return mesh.get_triangulation().n_global_levels(); - } }; template diff --git a/source/distributed/fully_distributed_tria.cc b/source/distributed/fully_distributed_tria.cc index 95ed52fcce..22033a0afa 100644 --- a/source/distributed/fully_distributed_tria.cc +++ b/source/distributed/fully_distributed_tria.cc @@ -701,6 +701,32 @@ namespace parallel } + + template + void + Triangulation::update_number_cache() + { + dealii::parallel::TriangulationBase::update_number_cache(); + + // additionally update the number of global coarse cells + types::coarse_cell_id number_of_global_coarse_cells = 0; + + for (const auto &cell : this->active_cell_iterators()) + if (!cell->is_artificial()) + number_of_global_coarse_cells = + std::max(number_of_global_coarse_cells, + cell->id().get_coarse_cell_id()); + + number_of_global_coarse_cells = + Utilities::MPI::max(number_of_global_coarse_cells, + this->mpi_communicator) + + 1; + + this->number_cache.number_of_global_coarse_cells = + number_of_global_coarse_cells; + } + + } // namespace fullydistributed } // namespace parallel diff --git a/source/distributed/repartitioning_policy_tools.cc b/source/distributed/repartitioning_policy_tools.cc index a7e35036f4..31e4739d15 100644 --- a/source/distributed/repartitioning_policy_tools.cc +++ b/source/distributed/repartitioning_policy_tools.cc @@ -28,31 +28,6 @@ namespace RepartitioningPolicyTools { namespace { - template - unsigned int - compute_n_coarse_cells(const MeshType &mesh) - { - types::coarse_cell_id n_coarse_cells = 0; - - for (const auto &cell : - mesh.get_triangulation().cell_iterators_on_level(0)) - n_coarse_cells = - std::max(n_coarse_cells, cell->id().get_coarse_cell_id()); - - return Utilities::MPI::max(n_coarse_cells, mesh.get_communicator()) + 1; - } - - - - template - unsigned int - compute_n_global_levels(const MeshType &mesh) - { - return mesh.get_triangulation().n_global_levels(); - } - - - template void add_indices_recursevly_for_first_child_policy( @@ -83,8 +58,8 @@ namespace RepartitioningPolicyTools template FirstChildPolicy::FirstChildPolicy( const Triangulation &tria_fine) - : n_coarse_cells(compute_n_coarse_cells(tria_fine)) - , n_global_levels(compute_n_global_levels(tria_fine)) + : n_coarse_cells(tria_fine.n_global_coarse_cells()) + , n_global_levels(tria_fine.n_global_levels()) { Assert( tria_fine.all_reference_cells_are_hyper_cube(), diff --git a/source/distributed/tria_base.cc b/source/distributed/tria_base.cc index 864e982435..0035491eb9 100644 --- a/source/distributed/tria_base.cc +++ b/source/distributed/tria_base.cc @@ -111,6 +111,7 @@ namespace parallel template TriangulationBase::NumberCache::NumberCache() : n_locally_owned_active_cells(0) + , number_of_global_coarse_cells(0) , n_global_levels(0) {} @@ -273,6 +274,8 @@ namespace parallel ExcInternalError()); } + this->number_cache.number_of_global_coarse_cells = this->n_cells(0); + // reset global cell ids this->reset_global_cell_indices(); } @@ -649,6 +652,15 @@ namespace parallel + template + types::coarse_cell_id + TriangulationBase::n_global_coarse_cells() const + { + return number_cache.number_of_global_coarse_cells; + } + + + template DistributedTriangulationBase::DistributedTriangulationBase( const MPI_Comm &mpi_communicator, diff --git a/source/grid/tria.cc b/source/grid/tria.cc index b5c3724aca..079b622b6c 100644 --- a/source/grid/tria.cc +++ b/source/grid/tria.cc @@ -13830,7 +13830,12 @@ Triangulation::n_global_active_cells() const return n_active_cells(); } - +template +types::coarse_cell_id +Triangulation::n_global_coarse_cells() const +{ + return n_cells(0); +} template unsigned int