--- /dev/null
+New: Add AdditionalData struct for configuring the RemotePointEvaluation class.
+The latter is passed into a new constructor.
+Deprecated: The old constructor of RemotePointEvaluation class taking the
+parameters value-by-value is marked as deprecated.
+<br>
+(Magdalena Schreter-Fleischhacker, Peter Munch, 2023/12/06)
class RemotePointEvaluation
{
public:
+ /**
+ * AdditionalData structure that can be used to tweak parameters
+ * of RemotePointEvaluation.
+ */
+ struct AdditionalData
+ {
+ public:
+ /**
+ * Constructor.
+ */
+ AdditionalData(
+ const double tolerance = 1e-6,
+ const bool enforce_unique_mapping = false,
+ const unsigned int rtree_level = 0,
+ const std::function<std::vector<bool>()> &marked_vertices = {});
+
+ /**
+ * Tolerance in terms of unit cell coordinates for determining all cells
+ * around a point passed to RemotePointEvaluation during reinit().
+ * Depending on the problem, it might be necessary to adjust the
+ * tolerance in order to be able to identify a cell. Floating point
+ * arithmetic implies that a point will, in general, not lie exactly on
+ * a vertex, edge, or face.
+ */
+ double tolerance;
+
+ /**
+ * Enforce unique mapping, i.e., (one-to-one) relation of points and
+ * cells.
+ */
+ bool enforce_unique_mapping;
+
+ /**
+ * RTree level to be used during the construction of the bounding boxes.
+ */
+ unsigned int rtree_level;
+
+ /**
+ * Function that marks relevant vertices to make search of active cells
+ * around point more efficient.
+ */
+ std::function<std::vector<bool>()> marked_vertices;
+ };
+
/**
* Constructor.
*
+ * @param additional_data Configure options for RemotePointEvaluation.
+ */
+ RemotePointEvaluation(
+ const AdditionalData &additional_data = AdditionalData());
+
+ /**
+ * Constructor. This constructor is deprecated. Use the other constructor
+ * taking AdditionalData instead.
+ *
* @param tolerance Tolerance in terms of unit cell coordinates for
* determining all cells around a point passed to the class during
* reinit(). Depending on the problem, it might be necessary to adjust
* @param rtree_level RTree level to be used during the construction of the bounding boxes.
* @param marked_vertices Function that marks relevant vertices to make search
* of active cells around point more efficient.
+ *
+ * @deprecated
*/
+ DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
+ "Use the constructor with AdditionalData struct.")
RemotePointEvaluation(
- const double tolerance = 1e-6,
+ const double tolerance,
const bool enforce_unique_mapping = false,
const unsigned int rtree_level = 0,
const std::function<std::vector<bool>()> &marked_vertices = {});
private:
/**
- * Tolerance to be used while determining the surrounding cells of a
- * point.
- */
- const double tolerance;
-
- /**
- * Enforce unique mapping, i.e., (one-to-one) relation of points and
- * cells.
- */
- const bool enforce_unique_mapping;
-
- /**
- * RTree level to be used during the construction of the bounding boxes.
- */
- const unsigned int rtree_level;
-
- /**
- * Function that marks relevant vertices to make search of active cells
- * around point more efficient.
+ * Additional data with basic settings.
*/
- const std::function<std::vector<bool>()> marked_vertices;
+ const AdditionalData additional_data;
/**
* Storage for the status of the triangulation signal.
MGTwoLevelTransferNonNested<dim, LinearAlgebra::distributed::Vector<Number>>::
MGTwoLevelTransferNonNested(const AdditionalData &data)
: additional_data(data)
- , rpe(data.tolerance, false, data.rtree_level, {})
+ , rpe(typename Utilities::MPI::RemotePointEvaluation<dim>::AdditionalData(
+ data.tolerance,
+ false,
+ data.rtree_level,
+ {}))
{}
template <int dim, typename Number>
namespace MPI
{
template <int dim, int spacedim>
- RemotePointEvaluation<dim, spacedim>::RemotePointEvaluation(
+ RemotePointEvaluation<dim, spacedim>::AdditionalData::AdditionalData(
const double tolerance,
const bool enforce_unique_mapping,
const unsigned int rtree_level,
, enforce_unique_mapping(enforce_unique_mapping)
, rtree_level(rtree_level)
, marked_vertices(marked_vertices)
+ {}
+
+
+
+ template <int dim, int spacedim>
+ RemotePointEvaluation<dim, spacedim>::RemotePointEvaluation(
+ const AdditionalData &additional_data)
+ : additional_data(additional_data)
+ , ready_flag(false)
+ {}
+
+
+
+ template <int dim, int spacedim>
+ RemotePointEvaluation<dim, spacedim>::RemotePointEvaluation(
+ const double tolerance,
+ const bool enforce_unique_mapping,
+ const unsigned int rtree_level,
+ const std::function<std::vector<bool>()> &marked_vertices)
+ : additional_data(tolerance,
+ enforce_unique_mapping,
+ rtree_level,
+ marked_vertices)
, ready_flag(false)
{}
// compress r-tree to a minimal set of bounding boxes
std::vector<std::vector<BoundingBox<spacedim>>> global_bboxes(1);
- global_bboxes[0] = extract_rtree_level(local_tree, rtree_level);
+ global_bboxes[0] =
+ extract_rtree_level(local_tree, additional_data.rtree_level);
const GridTools::Cache<dim, spacedim> cache(tria, mapping);
cache,
points,
global_bboxes,
- marked_vertices ? marked_vertices() : std::vector<bool>(),
- tolerance,
+ additional_data.marked_vertices ? additional_data.marked_vertices() :
+ std::vector<bool>(),
+ additional_data.tolerance,
true,
- enforce_unique_mapping);
+ additional_data.enforce_unique_mapping);
this->reinit(data, tria, mapping);
#endif
all_points_found_flag = std::get<0>(n_owning_processes_global) > 0;
}
- Assert(enforce_unique_mapping == false || unique_mapping,
+ Assert(additional_data.enforce_unique_mapping == false || unique_mapping,
ExcInternalError());
cell_data = std::make_unique<CellData>(tria);
for (unsigned int i = 0; i <= 4; ++i)
evaluation_points.emplace_back(i * 0.25, j * 0.25);
- Utilities::MPI::RemotePointEvaluation<dim> eval(1e-6, enforce_unique_map);
+ typename Utilities::MPI::RemotePointEvaluation<dim>::AdditionalData
+ additional_data;
+ additional_data.enforce_unique_mapping = enforce_unique_map;
+ additional_data.tolerance = 1e-6;
+ Utilities::MPI::RemotePointEvaluation<dim> eval(additional_data);
const auto result_avg =
VectorTools::point_values<1>(mapping,
// should be assigned to rank 0 for unique mapping
evaluation_points.emplace_back(0.5, 0.5 + std::pow<double>(10.0, -i));
- Utilities::MPI::RemotePointEvaluation<dim> eval(1e-6, enforce_unique_map);
+ typename Utilities::MPI::RemotePointEvaluation<dim>::AdditionalData
+ additional_data;
+ additional_data.enforce_unique_mapping = enforce_unique_map;
+ additional_data.tolerance = 1e-6;
+
+ Utilities::MPI::RemotePointEvaluation<dim> eval(additional_data);
const auto result_avg =
VectorTools::point_values<1>(mapping,
MappingQ1<dim> mapping;
- Utilities::MPI::RemotePointEvaluation<dim> eval(1.e-6, enforce_unique_map);
+ typename Utilities::MPI::RemotePointEvaluation<dim>::AdditionalData
+ additional_data;
+ additional_data.enforce_unique_mapping = enforce_unique_map;
+
+ Utilities::MPI::RemotePointEvaluation<dim> eval(additional_data);
eval.reinit(evaluation_points, tria, mapping);
}
// initialize RPE without any marked points
- dealii::Utilities::MPI::RemotePointEvaluation<dim> rpe(tolerance);
+ typename Utilities::MPI::RemotePointEvaluation<dim>::AdditionalData
+ additional_data;
+ additional_data.tolerance = tolerance;
+ dealii::Utilities::MPI::RemotePointEvaluation<dim> rpe(additional_data);
rpe.reinit(points, tria, mapping);
unsigned int n_points_not_found_rpe = 0;
// initialize RPE with all points marked
std::vector<bool> marked_vertices(tria.n_vertices(), true);
- dealii::Utilities::MPI::RemotePointEvaluation<dim> rpe2(
- tolerance, false, 0, [marked_vertices]() { return marked_vertices; });
+
+ typename Utilities::MPI::RemotePointEvaluation<dim>::AdditionalData
+ additional_data2(tolerance, false, 0, [marked_vertices]() {
+ return marked_vertices;
+ });
+
+ dealii::Utilities::MPI::RemotePointEvaluation<dim> rpe2(additional_data2);
rpe2.reinit(points, tria, mapping);