Skip to content

Commit

Permalink
Merge pull request #15899 from luca-heltai/rpe_AdditionalData_non_nested
Browse files Browse the repository at this point in the history
Expose RemotePointEvaluation parameters in non-nested multigrid
  • Loading branch information
peterrum committed Aug 24, 2023
2 parents b5564cc + f11181e commit 26c7c44
Show file tree
Hide file tree
Showing 4 changed files with 73 additions and 2 deletions.
51 changes: 51 additions & 0 deletions include/deal.II/multigrid/mg_transfer_global_coarsening.h
Original file line number Diff line number Diff line change
Expand Up @@ -706,6 +706,56 @@ class MGTwoLevelTransferNonNested<dim,
using VectorizedArrayType = VectorizedArray<Number, 1>;

public:
/**
* AdditionalData structure that can be used to tweak parameters
* related to the search procedure (used internally by RemotePointEvaluation)
* or, in the future, transfer operators needed by the non-nested multigrid
* algorithm.
*/
struct AdditionalData
{
/**
* Constructor. By default, the @p tolerance and @p rtree_level parameters
* are set to the default values used in the constructor of
* RemotePointEvaluation, i.e. 1e-6 and 0, respectively. The last Boolean
* parameter @p enforce_all_points_found is true by default and checks
* that all points submitted internally to RemotePointEvaluation::reinit()
* have been 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)
{}

/**
* Tolerance parameter. See the constructor of RemotePointEvaluation for
* more details.
*/
double tolerance;

/**
* RTree level parameter. See the constructor of RemotePointEvaluation for
* more details.
*
*/
unsigned int rtree_level;

/**
* If set to true, it checks if RemotePointEvaluation::all_points_found()
* evaluates to true internally during the each call to reinit() from one
* level to the next one, ensuring that all submitted points have been found
* inside the domain.
*
*/
bool enforce_all_points_found;
};

MGTwoLevelTransferNonNested(const AdditionalData &data = AdditionalData());

/**
* Set up transfer operator between the given DoFHandler objects (
* @p dof_handler_fine and @p dof_handler_coarse).
Expand Down Expand Up @@ -749,6 +799,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 @@ -4738,6 +4738,14 @@ namespace internal
} // namespace internal



template <int dim, typename Number>
MGTwoLevelTransferNonNested<dim, LinearAlgebra::distributed::Vector<Number>>::
MGTwoLevelTransferNonNested(const AdditionalData &data)
: additional_data(data)
, rpe(data.tolerance, false, data.rtree_level, {})
{}

template <int dim, typename Number>
void
MGTwoLevelTransferNonNested<dim, LinearAlgebra::distributed::Vector<Number>>::
Expand Down Expand Up @@ -4810,6 +4818,12 @@ MGTwoLevelTransferNonNested<dim, LinearAlgebra::distributed::Vector<Number>>::
// hand points over to RPE
rpe.reinit(points, dof_handler_coarse.get_triangulation(), mapping_coarse);

AssertThrow(
!additional_data.enforce_all_points_found || rpe.all_points_found(),
ExcMessage(
"You requested that all points should be found, but this didn'thappen."
" 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 26c7c44

Please sign in to comment.