Skip to content

Commit

Permalink
Allow throwing if not all points have been found
Browse files Browse the repository at this point in the history
  • Loading branch information
fdrmrc committed Aug 21, 2023
1 parent 4cda1db commit dda8f18
Show file tree
Hide file tree
Showing 4 changed files with 25 additions and 6 deletions.
14 changes: 10 additions & 4 deletions include/deal.II/multigrid/mg_transfer_global_coarsening.h
Original file line number Diff line number Diff line change
Expand Up @@ -707,25 +707,30 @@ class MGTwoLevelTransferNonNested<dim,

public:
/**
* AdditionalData structure for construction arguments needed by
* RemotePointEvaluation. Default values are the same as the ones in
* RemotePointEvaluation.
* 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.
*/
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 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)
{}
double tolerance;
bool enforce_unique_mapping;
unsigned int rtree_level;
std::function<std::vector<bool>()> marked_vertices;
bool enforce_all_points_found;
};

MGTwoLevelTransferNonNested(const AdditionalData &data = AdditionalData());
Expand Down Expand Up @@ -773,6 +778,7 @@ class MGTwoLevelTransferNonNested<dim,
memory_consumption() const override;

protected:
AdditionalData additional_data;
/**
* Perform prolongation.
*/
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -4104,6 +4104,7 @@ namespace internal
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,
Expand Down Expand Up @@ -4184,6 +4185,12 @@ MGTwoLevelTransferNonNested<dim, LinearAlgebra::distributed::Vector<Number>>::
// hand points over to RPE
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."));

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

Expand Down
4 changes: 3 additions & 1 deletion tests/multigrid-global-coarsening/non_nested_multigrid_03.cc
Original file line number Diff line number Diff line change
Expand Up @@ -77,10 +77,12 @@ test(const unsigned int n_refinements, const unsigned int fe_degree_fine)
}

// set up transfer operator
typename MGTwoLevelTransferNonNested<dim, VectorType>::AdditionalData data;
data.enforce_all_points_found = false;
for (unsigned int l = min_level; l < max_level; ++l)
{
transfers[l + 1] =
std::make_shared<MGTwoLevelTransferNonNested<dim, VectorType>>();
std::make_shared<MGTwoLevelTransferNonNested<dim, VectorType>>(data);
transfers[l + 1]->reinit(dof_handlers[l + 1],
dof_handlers[l],
mappings[l + 1],
Expand Down
6 changes: 5 additions & 1 deletion tests/multigrid-global-coarsening/non_nested_transfer_02.cc
Original file line number Diff line number Diff line change
Expand Up @@ -80,8 +80,12 @@ do_test(const FiniteElement<dim> & fe_fine,
constraint_fine.close();

// setup transfer operator
typename MGTwoLevelTransferNonNested<
dim,
LinearAlgebra::distributed::Vector<Number>>::AdditionalData data;
data.enforce_all_points_found = false;
MGTwoLevelTransferNonNested<dim, LinearAlgebra::distributed::Vector<Number>>
transfer;
transfer(data);
MappingQ1<dim> mapping_fine, mapping_coarse;
transfer.reinit(dof_handler_fine,
dof_handler_coarse,
Expand Down

0 comments on commit dda8f18

Please sign in to comment.