Skip to content

Commit

Permalink
PR review
Browse files Browse the repository at this point in the history
  • Loading branch information
fdrmrc committed Aug 22, 2023
1 parent dda8f18 commit e39c619
Show file tree
Hide file tree
Showing 2 changed files with 40 additions and 49 deletions.
39 changes: 18 additions & 21 deletions include/deal.II/multigrid/mg_transfer_global_coarsening.h
Original file line number Diff line number Diff line change
Expand Up @@ -707,30 +707,27 @@ class MGTwoLevelTransferNonNested<dim,

public:
/**
* AdditionalData structure with the arguments needed by
* RemotePointEvaluation. Default values are the same as the ones described in
* the documentation of RemotePointEvaluation. The last boolean parameter, @p enf_all_points_found is true by defaults and
* checks if RemotePointEvaluation::all_points_found() evaluates to true, i.e.
* all submitted points have been found inside the domain.
* AdditionalData structure that can be used to tweak parameters
* @p tolerance and @p rtree_level related to the search procedure (used
* internally by RemotePointEvaluation) or, in the future, transfer operators
* needed by the non-nested multigrid algorithm. Default values are the same
* as the ones described in the documentation of RemotePointEvaluation.
* The last Boolean parameter, @p enforce_all_points_found is true by default
* and checks if RemotePointEvaluation::all_points_found() evaluates to true,
* i.e. all submitted points have been found inside the domain.
*/
struct AdditionalData
{
AdditionalData(const double tol = 1e-6,
const bool enf_unique_mapping = false,
const unsigned int rtree_l = 0,
const std::function<std::vector<bool>()> &marked_verts = {},
const bool enf_all_points_found = true)
: tolerance(tol)
, enforce_unique_mapping(enf_unique_mapping)
, rtree_level(rtree_l)
, marked_vertices(marked_verts)
, enforce_all_points_found(enf_all_points_found)
AdditionalData(const double tolerance = 1e-6,
const unsigned int rtree_level = 0,
const bool enforce_all_points_found = true)
: tolerance(tolerance)
, rtree_level(rtree_level)
, enforce_all_points_found(enforce_all_points_found)
{}
double tolerance;
bool enforce_unique_mapping;
unsigned int rtree_level;
std::function<std::vector<bool>()> marked_vertices;
bool enforce_all_points_found;
double tolerance;
unsigned int rtree_level;
bool enforce_all_points_found;
};

MGTwoLevelTransferNonNested(const AdditionalData &data = AdditionalData());
Expand Down Expand Up @@ -818,7 +815,7 @@ class MGTwoLevelTransferNonNested<dim,
* Object to evaluate shape functions on one mesh on visited support points of
* the other mesh.
*/
std::shared_ptr<Utilities::MPI::RemotePointEvaluation<dim>> rpe;
Utilities::MPI::RemotePointEvaluation<dim> rpe;

/**
* MappingInfo object needed as Mapping argument by FEPointEvaluation.
Expand Down
50 changes: 22 additions & 28 deletions include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h
Original file line number Diff line number Diff line change
Expand Up @@ -4105,13 +4105,8 @@ template <int dim, typename Number>
MGTwoLevelTransferNonNested<dim, LinearAlgebra::distributed::Vector<Number>>::
MGTwoLevelTransferNonNested(const AdditionalData &data)
: additional_data(data)
{
rpe = std::make_shared<Utilities::MPI::RemotePointEvaluation<dim>>(
data.tolerance,
data.enforce_unique_mapping,
data.rtree_level,
data.marked_vertices);
}
, rpe(data.tolerance, false, data.rtree_level, {})
{}

template <int dim, typename Number>
void
Expand Down Expand Up @@ -4183,19 +4178,18 @@ MGTwoLevelTransferNonNested<dim, LinearAlgebra::distributed::Vector<Number>>::
this->partitioner_fine->global_to_local(global_dof_indices[i]);

// hand points over to RPE
rpe->reinit(points, dof_handler_coarse.get_triangulation(), mapping_coarse);
rpe.reinit(points, dof_handler_coarse.get_triangulation(), mapping_coarse);

if (additional_data.enforce_all_points_found)
AssertThrow(
rpe->all_points_found(),
ExcMessage(
"You requested that all points should be found, but this didn't happen. You can change this option through the AdditionaData struct in the constructor."));
AssertThrow(
!additional_data.enforce_all_points_found || rpe.all_points_found(),
ExcMessage(
"You requested that all points should be found, but this didn't happen. You can change this option through the AdditionaData struct in the constructor."));

// set up MappingInfo for easier data access
mapping_info = internal::fill_mapping_info<dim, Number>(*rpe);
mapping_info = internal::fill_mapping_info<dim, Number>(rpe);

// set up constraints
const auto &cell_data = rpe->get_cell_data();
const auto &cell_data = rpe.get_cell_data();

constraint_info.reinit(dof_handler_coarse,
cell_data.cells.size(),
Expand All @@ -4204,7 +4198,7 @@ MGTwoLevelTransferNonNested<dim, LinearAlgebra::distributed::Vector<Number>>::
for (unsigned int i = 0; i < cell_data.cells.size(); ++i)
{
typename DoFHandler<dim>::active_cell_iterator cell(
&rpe->get_triangulation(),
&rpe.get_triangulation(),
cell_data.cells[i].first,
cell_data.cells[i].second,
&dof_handler_coarse);
Expand Down Expand Up @@ -4314,18 +4308,18 @@ MGTwoLevelTransferNonNested<dim, LinearAlgebra::distributed::Vector<Number>>::
}
};

rpe->template evaluate_and_process<value_type>(evaluation_point_results,
buffer,
evaluation_function);
rpe.template evaluate_and_process<value_type>(evaluation_point_results,
buffer,
evaluation_function);

// Weight operator in case some points are owned by multiple cells.
if (rpe->is_map_unique() == false)
if (rpe.is_map_unique() == false)
{
const auto evaluation_point_results_temp = evaluation_point_results;
evaluation_point_results.clear();
evaluation_point_results.reserve(rpe->get_point_ptrs().size() - 1);
evaluation_point_results.reserve(rpe.get_point_ptrs().size() - 1);

const auto &ptr = rpe->get_point_ptrs();
const auto &ptr = rpe.get_point_ptrs();

for (unsigned int i = 0; i < ptr.size() - 1; ++i)
{
Expand Down Expand Up @@ -4407,7 +4401,7 @@ MGTwoLevelTransferNonNested<dim, LinearAlgebra::distributed::Vector<Number>>::
std::vector<value_type> evaluation_point_results;
std::vector<value_type> buffer;

evaluation_point_results.resize(rpe->get_point_ptrs().size() - 1);
evaluation_point_results.resize(rpe.get_point_ptrs().size() - 1);

for (unsigned int j = 0; j < evaluation_point_results.size(); ++j)
{
Expand Down Expand Up @@ -4442,9 +4436,9 @@ MGTwoLevelTransferNonNested<dim, LinearAlgebra::distributed::Vector<Number>>::
}

// Weight operator in case some points are owned by multiple cells.
if (rpe->is_map_unique() == false)
if (rpe.is_map_unique() == false)
{
const auto &ptr = rpe->get_point_ptrs();
const auto &ptr = rpe.get_point_ptrs();

for (unsigned int i = 0; i < ptr.size() - 1; ++i)
{
Expand Down Expand Up @@ -4489,9 +4483,9 @@ MGTwoLevelTransferNonNested<dim, LinearAlgebra::distributed::Vector<Number>>::
}
};

rpe->template process_and_evaluate<value_type>(evaluation_point_results,
buffer,
evaluation_function);
rpe.template process_and_evaluate<value_type>(evaluation_point_results,
buffer,
evaluation_function);
}


Expand Down

0 comments on commit e39c619

Please sign in to comment.