38 template <
int dim,
int spacedim>
40 const double tolerance,
41 const bool enforce_unique_mapping,
42 const unsigned int rtree_level,
43 const std::function<std::vector<bool>()> &marked_vertices)
44 : tolerance(tolerance)
45 , enforce_unique_mapping(enforce_unique_mapping)
46 , rtree_level(rtree_level)
47 , marked_vertices(marked_vertices)
53 template <
int dim,
int spacedim>
56 if (tria_signal.connected())
57 tria_signal.disconnect();
62 template <
int dim,
int spacedim>
69#ifndef DEAL_II_WITH_MPI
75 if (tria_signal.connected())
76 tria_signal.disconnect();
79 tria.
signals.any_change.connect([&]() { this->ready_flag =
false; });
82 this->mapping = &mapping;
84 std::vector<BoundingBox<spacedim>> local_boxes;
85 for (
const auto &cell :
90 const auto local_tree =
pack_rtree(local_boxes);
93 const auto local_reduced_box =
97 const auto global_bboxes =
107 marked_vertices ? marked_vertices() : std::vector<bool>(),
110 enforce_unique_mapping);
112 this->recv_ranks = data.recv_ranks;
113 this->recv_ptrs = data.recv_ptrs;
115 this->send_ranks = data.send_ranks;
116 this->send_ptrs = data.send_ptrs;
118 this->recv_permutation = {};
119 this->recv_permutation.resize(data.recv_components.size());
120 this->point_ptrs.assign(points.size() + 1, 0);
121 for (
unsigned int i = 0; i < data.recv_components.size(); ++i)
124 this->recv_permutation.size());
125 this->recv_permutation[std::get<2>(data.recv_components[i])] = i;
128 this->point_ptrs.size());
129 this->point_ptrs[std::get<1>(data.recv_components[i]) + 1]++;
132 std::tuple<unsigned int, unsigned int> n_owning_processes_default{
134 std::tuple<unsigned int, unsigned int> n_owning_processes_local =
135 n_owning_processes_default;
137 for (
unsigned int i = 0; i < points.size(); ++i)
139 std::get<0>(n_owning_processes_local) =
140 std::min(std::get<0>(n_owning_processes_local),
141 this->point_ptrs[i + 1]);
142 std::get<1>(n_owning_processes_local) =
143 std::max(std::get<1>(n_owning_processes_local),
144 this->point_ptrs[i + 1]);
146 this->point_ptrs[i + 1] += this->point_ptrs[i];
149 const auto n_owning_processes_global =
150 Utilities::MPI::all_reduce<std::tuple<unsigned int, unsigned int>>(
151 n_owning_processes_local,
154 const auto &b) -> std::tuple<unsigned int, unsigned int> {
155 if (a == n_owning_processes_default)
158 if (b == n_owning_processes_default)
161 return std::tuple<unsigned int, unsigned int>{
162 std::min(std::get<0>(a), std::get<0>(b)),
163 std::max(std::get<1>(a), std::get<1>(b))};
166 if (n_owning_processes_global == n_owning_processes_default)
168 unique_mapping =
true;
169 all_points_found_flag =
true;
173 unique_mapping = (std::get<0>(n_owning_processes_global) == 1) &&
174 (std::get<1>(n_owning_processes_global) == 1);
175 all_points_found_flag = std::get<0>(n_owning_processes_global) > 0;
178 Assert(enforce_unique_mapping ==
false || unique_mapping,
182 send_permutation = {};
184 std::pair<int, int> dummy{-1, -1};
185 for (
const auto &i : data.send_components)
187 if (dummy != std::get<0>(i))
189 dummy = std::get<0>(i);
190 cell_data.cells.emplace_back(dummy);
191 cell_data.reference_point_ptrs.emplace_back(
192 cell_data.reference_point_values.size());
195 cell_data.reference_point_values.emplace_back(std::get<3>(i));
196 send_permutation.emplace_back(std::get<5>(i));
199 cell_data.reference_point_ptrs.emplace_back(
200 cell_data.reference_point_values.size());
202 this->ready_flag =
true;
207 template <
int dim,
int spacedim>
208 const std::vector<unsigned int> &
216 template <
int dim,
int spacedim>
220 return unique_mapping;
225 template <
int dim,
int spacedim>
229 return all_points_found_flag;
234 template <
int dim,
int spacedim>
237 const unsigned int i)
const
241 if (all_points_found_flag)
244 return (point_ptrs[i + 1] - point_ptrs[i]) > 0;
249 template <
int dim,
int spacedim>
258 template <
int dim,
int spacedim>
267 template <
int dim,
int spacedim>
277#include "mpi_remote_point_evaluation.inst"
Abstract base class for mapping classes.
virtual BoundingBox< spacedim > get_bounding_box(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
virtual MPI_Comm get_communicator() const
RemotePointEvaluation(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={})
bool all_points_found() const
const Triangulation< dim, spacedim > & get_triangulation() const
const std::vector< unsigned int > & get_point_ptrs() const
bool is_map_unique() const
const Mapping< dim, spacedim > & get_mapping() const
void reinit(const std::vector< Point< spacedim > > &points, const Triangulation< dim, spacedim > &tria, const Mapping< dim, spacedim > &mapping)
bool point_found(const unsigned int i) const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
IteratorRange< active_cell_iterator > active_cell_iterators() const
static ::ExceptionBase & ExcNeedsMPI()
#define Assert(cond, exc)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
std::vector< T > all_gather(const MPI_Comm &comm, const T &object_to_send)
static const unsigned int invalid_unsigned_int
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
std::vector< BoundingBox< boost::geometry::dimension< typename Rtree::indexable_type >::value > > extract_rtree_level(const Rtree &tree, const unsigned int level)
RTree< typename LeafTypeIterator::value_type, IndexType, IndexableGetter > pack_rtree(const LeafTypeIterator &begin, const LeafTypeIterator &end)
const ::Triangulation< dim, spacedim > & tria