From: Peter Munch Date: Tue, 28 Jun 2022 19:46:26 +0000 (+0200) Subject: p:d:GridRefinement: accept other types of triangulations X-Git-Tag: v9.5.0-rc1~1093^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=3450abcf7ae997ad4e5ec2e2f7058d52eabd4b80;p=dealii.git p:d:GridRefinement: accept other types of triangulations --- diff --git a/include/deal.II/distributed/grid_refinement.h b/include/deal.II/distributed/grid_refinement.h index 96e73eaff9..bf2adf4c5f 100644 --- a/include/deal.II/distributed/grid_refinement.h +++ b/include/deal.II/distributed/grid_refinement.h @@ -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., @@ -152,8 +154,8 @@ namespace parallel template void refine_and_coarsen_fixed_number( - parallel::distributed::Triangulation &tria, - const dealii::Vector & criteria, + Triangulation & tria, + const dealii::Vector & criteria, const double top_fraction_of_cells, const double bottom_fraction_of_cells, const types::global_cell_index max_n_cells = @@ -161,7 +163,8 @@ namespace parallel /** * 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., @@ -205,10 +208,10 @@ namespace parallel template void refine_and_coarsen_fixed_fraction( - parallel::distributed::Triangulation &tria, - const dealii::Vector & criteria, - const double top_fraction_of_error, - const double bottom_fraction_of_error, + Triangulation &tria, + const dealii::Vector &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 diff --git a/source/distributed/grid_refinement.cc b/source/distributed/grid_refinement.cc index f5100915a2..6a9a5ceff0 100644 --- a/source/distributed/grid_refinement.cc +++ b/source/distributed/grid_refinement.cc @@ -42,6 +42,18 @@ DEAL_II_NAMESPACE_OPEN namespace { + template + unsigned int + n_locally_owned_active_cells(const Triangulation &tria) + { + if (const auto parallel_tria = + dynamic_cast *>( + &tria)) + return parallel_tria->n_locally_owned_active_cells(); + else + return tria.n_active_cells(); + } + template inline number max_element(const dealii::Vector &criteria) @@ -107,7 +119,7 @@ namespace Vector &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; @@ -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()); } @@ -207,8 +219,7 @@ namespace { // first extract from the vector of indicators the ones that correspond // to cells that we locally own - Vector locally_owned_indicators( - tria.n_locally_owned_active_cells()); + Vector 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(); @@ -520,7 +531,7 @@ namespace parallel // first extract from the vector of indicators the ones that correspond // to cells that we locally own Vector 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(); diff --git a/source/distributed/grid_refinement.inst.in b/source/distributed/grid_refinement.inst.in index 48ab446542..8f3a7e025d 100644 --- a/source/distributed/grid_refinement.inst.in +++ b/source/distributed/grid_refinement.inst.in @@ -65,7 +65,7 @@ for (S : REAL_SCALARS; deal_II_dimension : DIMENSIONS) refine_and_coarsen_fixed_number( - parallel::distributed::Triangulation &, + Triangulation &, const dealii::Vector &, const double, const double, @@ -75,7 +75,7 @@ for (S : REAL_SCALARS; deal_II_dimension : DIMENSIONS) refine_and_coarsen_fixed_fraction( - parallel::distributed::Triangulation &, + Triangulation &, const dealii::Vector &, const double, const double, @@ -97,8 +97,7 @@ for (S : REAL_SCALARS; deal_II_dimension : DIMENSIONS) refine_and_coarsen_fixed_number( - parallel::distributed::Triangulation &, + Triangulation &, const dealii::Vector &, const double, const double, @@ -108,8 +107,7 @@ for (S : REAL_SCALARS; deal_II_dimension : DIMENSIONS) refine_and_coarsen_fixed_fraction( - parallel::distributed::Triangulation &, + Triangulation &, const dealii::Vector &, const double, const double,