Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Triangulation: implement global_*_cell_index_partitioner() #14270

Merged
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
15 changes: 3 additions & 12 deletions include/deal.II/distributed/tria_base.h
Original file line number Diff line number Diff line change
Expand Up @@ -202,21 +202,12 @@ namespace parallel
const std::set<types::subdomain_id> &
level_ghost_owners() const;

/**
* Return partitioner for the global indices of the cells on the active
* level of the triangulation, which is returned by the function
* CellAccessor::global_active_cell_index().
*/
const std::weak_ptr<const Utilities::MPI::Partitioner>
global_active_cell_index_partitioner() const;
global_active_cell_index_partitioner() const override;
drwells marked this conversation as resolved.
Show resolved Hide resolved

/**
* Return partitioner for the global indices of the cells on the given @p
* level of the triangulation, which is returned by the function
* CellAccessor::global_level_cell_index().
*/
const std::weak_ptr<const Utilities::MPI::Partitioner>
global_level_cell_index_partitioner(const unsigned int level) const;
global_level_cell_index_partitioner(
const unsigned int level) const override;

/**
* @copydoc dealii::Triangulation::get_boundary_ids()
Expand Down
29 changes: 29 additions & 0 deletions include/deal.II/grid/tria.h
Original file line number Diff line number Diff line change
Expand Up @@ -21,6 +21,7 @@

#include <deal.II/base/geometry_info.h>
#include <deal.II/base/iterator_range.h>
#include <deal.II/base/partitioner.h>
#include <deal.II/base/point.h>
#include <deal.II/base/smartpointer.h>
#include <deal.II/base/subscriptor.h>
Expand Down Expand Up @@ -176,6 +177,18 @@ namespace internal
*/
std::vector<unsigned int> n_active_lines_level;

/**
* Partitioner for the global active cell indices.
*/
std::shared_ptr<const Utilities::MPI::Partitioner>
active_cell_index_partitioner;

/**
* Partitioner for the global level cell indices for each level.
*/
std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
level_cell_index_partitioners;

/**
* Constructor. Set values to zero by default.
*/
Expand Down Expand Up @@ -1615,6 +1628,22 @@ class Triangulation : public Subscriptor
virtual MPI_Comm
get_communicator() const;

/**
* Return the partitioner for the global indices of the cells on the active
* level of the triangulation, which is returned by the function
* CellAccessor::global_active_cell_index().
*/
virtual const std::weak_ptr<const Utilities::MPI::Partitioner>
global_active_cell_index_partitioner() const;

/**
* Return the partitioner for the global indices of the cells on the given @p
* level of the triangulation, which is returned by the function
* CellAccessor::global_level_cell_index().
*/
virtual const std::weak_ptr<const Utilities::MPI::Partitioner>
global_level_cell_index_partitioner(const unsigned int level) const;

/**
* Set the mesh smoothing to @p mesh_smoothing. This overrides the
* MeshSmoothing given to the constructor. It is allowed to call this
Expand Down
52 changes: 47 additions & 5 deletions source/grid/tria.cc
Original file line number Diff line number Diff line change
Expand Up @@ -1757,7 +1757,7 @@ namespace internal
*/
template <int dim, int spacedim>
static void
compute_number_cache(
compute_number_cache_dim(
const Triangulation<dim, spacedim> & triangulation,
const unsigned int level_objects,
internal::TriangulationImplementation::NumberCache<1> &number_cache)
Expand Down Expand Up @@ -1845,7 +1845,7 @@ namespace internal
*/
template <int dim, int spacedim>
static void
compute_number_cache(
compute_number_cache_dim(
const Triangulation<dim, spacedim> & triangulation,
const unsigned int level_objects,
internal::TriangulationImplementation::NumberCache<2> &number_cache)
Expand All @@ -1858,7 +1858,7 @@ namespace internal
void (*)(const Triangulation<dim, spacedim> &,
const unsigned int,
internal::TriangulationImplementation::NumberCache<1> &)>(
&compute_number_cache<dim, spacedim>),
&compute_number_cache_dim<dim, spacedim>),
triangulation,
level_objects,
static_cast<internal::TriangulationImplementation::NumberCache<1> &>(
Expand Down Expand Up @@ -1952,7 +1952,7 @@ namespace internal
*/
template <int dim, int spacedim>
static void
compute_number_cache(
compute_number_cache_dim(
const Triangulation<dim, spacedim> & triangulation,
const unsigned int level_objects,
internal::TriangulationImplementation::NumberCache<3> &number_cache)
Expand All @@ -1965,7 +1965,7 @@ namespace internal
void (*)(const Triangulation<dim, spacedim> &,
const unsigned int,
internal::TriangulationImplementation::NumberCache<2> &)>(
&compute_number_cache<dim, spacedim>),
&compute_number_cache_dim<dim, spacedim>),
triangulation,
level_objects,
static_cast<internal::TriangulationImplementation::NumberCache<2> &>(
Expand Down Expand Up @@ -2043,6 +2043,27 @@ namespace internal
}


template <int dim, int spacedim>
static void
compute_number_cache(
const Triangulation<dim, spacedim> & triangulation,
const unsigned int level_objects,
internal::TriangulationImplementation::NumberCache<dim> &number_cache)
{
compute_number_cache_dim(triangulation, level_objects, number_cache);

number_cache.active_cell_index_partitioner =
std::make_shared<const Utilities::MPI::Partitioner>(
triangulation.n_active_cells());

number_cache.level_cell_index_partitioners.resize(
triangulation.n_levels());
for (unsigned int level = 0; level < triangulation.n_levels(); ++level)
number_cache.level_cell_index_partitioners[level] =
std::make_shared<const Utilities::MPI::Partitioner>(
triangulation.n_cells(level));
}


template <int spacedim>
static void
Expand Down Expand Up @@ -11207,6 +11228,27 @@ Triangulation<dim, spacedim>::get_communicator() const



template <int dim, int spacedim>
const std::weak_ptr<const Utilities::MPI::Partitioner>
Triangulation<dim, spacedim>::global_active_cell_index_partitioner() const
{
return number_cache.active_cell_index_partitioner;
}



template <int dim, int spacedim>
const std::weak_ptr<const Utilities::MPI::Partitioner>
Triangulation<dim, spacedim>::global_level_cell_index_partitioner(
const unsigned int level) const
{
AssertIndexRange(level, this->n_levels());

return number_cache.level_cell_index_partitioners[level];
}



template <int dim, int spacedim>
void
Triangulation<dim, spacedim>::set_mesh_smoothing(
Expand Down