using VectorizedArrayType = VectorizedArray<Number, 1>;
public:
+ /**
+ * AdditionalData structure for construction arguments needed by
+ * RemotePointEvaluation. Default values are the same as the ones in
+ * RemotePointEvaluation.
+ */
+ 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 = {})
+ : tolerance(tol)
+ , enforce_unique_mapping(enf_unique_mapping)
+ , rtree_level(rtree_l)
+ , marked_vertices(marked_verts)
+ {}
+ double tolerance;
+ bool enforce_unique_mapping;
+ unsigned int rtree_level;
+ std::function<std::vector<bool>()> marked_vertices;
+ };
+
+ MGTwoLevelTransferNonNested(const AdditionalData &data = AdditionalData());
+
/**
* Set up transfer operator between the given DoFHandler objects (
* @p dof_handler_fine and @p dof_handler_coarse).
* Object to evaluate shape functions on one mesh on visited support points of
* the other mesh.
*/
- Utilities::MPI::RemotePointEvaluation<dim> rpe;
+ std::shared_ptr<Utilities::MPI::RemotePointEvaluation<dim>> rpe;
/**
* MappingInfo object needed as Mapping argument by FEPointEvaluation.
} // namespace internal
+
+template <int dim, typename Number>
+MGTwoLevelTransferNonNested<dim, LinearAlgebra::distributed::Vector<Number>>::
+ MGTwoLevelTransferNonNested(const AdditionalData &data)
+{
+ rpe = std::make_shared<Utilities::MPI::RemotePointEvaluation<dim>>(
+ data.tolerance,
+ data.enforce_unique_mapping,
+ data.rtree_level,
+ data.marked_vertices);
+}
+
template <int dim, typename Number>
void
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);
// 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(),
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);
}
};
- 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)
{
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)
{
}
// 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)
{
}
};
- 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);
}