Skip to content

Commit

Permalink
p:d:GridRefinement: accept other types of triangulations
Browse files Browse the repository at this point in the history
  • Loading branch information
peterrum committed Jun 30, 2022
1 parent 685249d commit 3450abc
Show file tree
Hide file tree
Showing 3 changed files with 34 additions and 22 deletions.
25 changes: 14 additions & 11 deletions include/deal.II/distributed/grid_refinement.h
Original file line number Diff line number Diff line change
Expand Up @@ -98,16 +98,18 @@ namespace parallel
* cycle of the adaptive finite element loop.
*
* In contrast to the functions in namespace dealii::GridRefinement,
* the functions in the current namespace are intended for distributed
* meshes, i.e., objects of type parallel::distributed::Triangulation.
* the functions in the current namespace are intended for parallel
* computations, i.e., computations, e.g., on objects of type
* parallel::distributed::Triangulation.
*
* @ingroup grid
*/
namespace GridRefinement
{
/**
* Like dealii::GridRefinement::refine_and_coarsen_fixed_number, but for
* parallel distributed triangulations.
* Like dealii::GridRefinement::refine_and_coarsen_fixed_number, but
* designed for parallel computations, where each process has only
* information about locally owned cells.
*
* The vector of criteria needs to be a vector of refinement criteria
* for all cells active on the current triangulation, i.e.,
Expand Down Expand Up @@ -152,16 +154,17 @@ namespace parallel
template <int dim, typename Number, int spacedim>
void
refine_and_coarsen_fixed_number(
parallel::distributed::Triangulation<dim, spacedim> &tria,
const dealii::Vector<Number> & criteria,
Triangulation<dim, spacedim> & tria,
const dealii::Vector<Number> & criteria,
const double top_fraction_of_cells,
const double bottom_fraction_of_cells,
const types::global_cell_index max_n_cells =
std::numeric_limits<types::global_cell_index>::max());

/**
* Like dealii::GridRefinement::refine_and_coarsen_fixed_fraction, but
* for parallel distributed triangulations.
* designed for parallel computations, where each process only has
* information about locally owned cells.
*
* The vector of criteria needs to be a vector of refinement criteria
* for all cells active on the current triangulation, i.e.,
Expand Down Expand Up @@ -205,10 +208,10 @@ namespace parallel
template <int dim, typename Number, int spacedim>
void
refine_and_coarsen_fixed_fraction(
parallel::distributed::Triangulation<dim, spacedim> &tria,
const dealii::Vector<Number> & criteria,
const double top_fraction_of_error,
const double bottom_fraction_of_error,
Triangulation<dim, spacedim> &tria,
const dealii::Vector<Number> &criteria,
const double top_fraction_of_error,
const double bottom_fraction_of_error,
const VectorTools::NormType norm_type = VectorTools::NormType::L1_norm);
} // namespace GridRefinement
} // namespace distributed
Expand Down
21 changes: 16 additions & 5 deletions source/distributed/grid_refinement.cc
Original file line number Diff line number Diff line change
Expand Up @@ -42,6 +42,18 @@ DEAL_II_NAMESPACE_OPEN

namespace
{
template <int dim, int spacedim>
unsigned int
n_locally_owned_active_cells(const Triangulation<dim, spacedim> &tria)
{
if (const auto parallel_tria =
dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
&tria))
return parallel_tria->n_locally_owned_active_cells();
else
return tria.n_active_cells();
}

template <typename number>
inline number
max_element(const dealii::Vector<number> &criteria)
Expand Down Expand Up @@ -107,7 +119,7 @@ namespace
Vector<Number> &locally_owned_indicators)
{
Assert(locally_owned_indicators.size() ==
tria.n_locally_owned_active_cells(),
n_locally_owned_active_cells(tria),
ExcInternalError());

unsigned int owned_index = 0;
Expand All @@ -118,7 +130,7 @@ namespace
criteria(cell->active_cell_index());
++owned_index;
}
Assert(owned_index == tria.n_locally_owned_active_cells(),
Assert(owned_index == n_locally_owned_active_cells(tria),
ExcInternalError());
}

Expand Down Expand Up @@ -207,8 +219,7 @@ namespace
{
// first extract from the vector of indicators the ones that correspond
// to cells that we locally own
Vector<Number> locally_owned_indicators(
tria.n_locally_owned_active_cells());
Vector<Number> locally_owned_indicators(n_locally_owned_active_cells(tria));
get_locally_owned_indicators(tria, criteria, locally_owned_indicators);

MPI_Comm mpi_communicator = tria.get_communicator();
Expand Down Expand Up @@ -520,7 +531,7 @@ namespace parallel
// first extract from the vector of indicators the ones that correspond
// to cells that we locally own
Vector<Number> locally_owned_indicators(
tria.n_locally_owned_active_cells());
n_locally_owned_active_cells(tria));
get_locally_owned_indicators(tria, criteria, locally_owned_indicators);

MPI_Comm mpi_communicator = tria.get_communicator();
Expand Down
10 changes: 4 additions & 6 deletions source/distributed/grid_refinement.inst.in
Original file line number Diff line number Diff line change
Expand Up @@ -65,7 +65,7 @@ for (S : REAL_SCALARS; deal_II_dimension : DIMENSIONS)
refine_and_coarsen_fixed_number<deal_II_dimension,
S,
deal_II_dimension>(
parallel::distributed::Triangulation<deal_II_dimension> &,
Triangulation<deal_II_dimension> &,
const dealii::Vector<S> &,
const double,
const double,
Expand All @@ -75,7 +75,7 @@ for (S : REAL_SCALARS; deal_II_dimension : DIMENSIONS)
refine_and_coarsen_fixed_fraction<deal_II_dimension,
S,
deal_II_dimension>(
parallel::distributed::Triangulation<deal_II_dimension> &,
Triangulation<deal_II_dimension> &,
const dealii::Vector<S> &,
const double,
const double,
Expand All @@ -97,8 +97,7 @@ for (S : REAL_SCALARS; deal_II_dimension : DIMENSIONS)
refine_and_coarsen_fixed_number<deal_II_dimension - 1,
S,
deal_II_dimension>(
parallel::distributed::Triangulation<deal_II_dimension - 1,
deal_II_dimension> &,
Triangulation<deal_II_dimension - 1, deal_II_dimension> &,
const dealii::Vector<S> &,
const double,
const double,
Expand All @@ -108,8 +107,7 @@ for (S : REAL_SCALARS; deal_II_dimension : DIMENSIONS)
refine_and_coarsen_fixed_fraction<deal_II_dimension - 1,
S,
deal_II_dimension>(
parallel::distributed::Triangulation<deal_II_dimension - 1,
deal_II_dimension> &,
Triangulation<deal_II_dimension - 1, deal_II_dimension> &,
const dealii::Vector<S> &,
const double,
const double,
Expand Down

0 comments on commit 3450abc

Please sign in to comment.