preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
RTreeAdapter.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <Eigen/Core>
4#include <boost/geometry.hpp>
5#include "mesh/Edge.hpp"
6#include "mesh/Mesh.hpp"
8#include "mesh/Vertex.hpp"
9#include "utils/assertion.hpp"
10
11namespace precice::mesh {
12class Triangle;
13} // namespace precice::mesh
14
15namespace pm = precice::mesh;
16namespace bg = boost::geometry;
17
19
21/*
22 * This adapts every VectorXd to a 3d point. For non-existing dimensions, zero is returned.
23 */
24template <>
25struct tag<Eigen::VectorXd> {
26 using type = point_tag;
27};
28template <>
29struct coordinate_type<Eigen::VectorXd> {
30 using type = double;
31};
32template <>
33struct coordinate_system<Eigen::VectorXd> {
34 using type = cs::cartesian;
35};
36template <>
37struct dimension<Eigen::VectorXd> : boost::mpl::int_<3> {
38};
39
40template <size_t Dimension>
41struct access<Eigen::VectorXd, Dimension> {
42 static double get(Eigen::VectorXd const &p)
43 {
44 if (Dimension >= static_cast<size_t>(p.rows())) {
45 return 0;
46 } else {
47 return p(Dimension);
48 }
49 }
50
51 static void set(Eigen::VectorXd &p, double const &value)
52 {
53 if (Dimension >= static_cast<size_t>(p.rows())) {
54 Eigen::VectorXd tmp(Dimension);
55 std::copy_n(p.data(), Dimension - 1, tmp.data());
56 p = std::move(tmp);
57 }
58 p(Dimension) = value;
59 }
60};
61
63template <>
64struct tag<pm::Vertex::RawCoords> {
65 using type = point_tag;
66};
67template <>
68struct coordinate_type<pm::Vertex::RawCoords> {
69 using type = double;
70};
71template <>
72struct coordinate_system<pm::Vertex::RawCoords> {
73 using type = cs::cartesian;
74};
75template <>
76struct dimension<pm::Vertex::RawCoords> : boost::mpl::int_<3> {
77};
78
79template <size_t Dimension>
80struct access<pm::Vertex::RawCoords, Dimension> {
81 static double get(pm::Vertex::RawCoords const &p)
82 {
83 return p[Dimension];
84 }
85
86 static void set(pm::Vertex::RawCoords &p, double const &value)
87 {
88 p[Dimension] = value;
89 }
90};
91
92BOOST_CONCEPT_ASSERT((bg::concepts::Point<pm::Vertex::RawCoords>) );
93
95/*
96 * This adapts every Vertex to a 3d point. For non-existing dimensions, zero is returned.
97 */
98template <>
99struct tag<pm::Vertex> {
100 using type = point_tag;
101};
102template <>
103struct coordinate_type<pm::Vertex> {
104 using type = double;
105};
106template <>
107struct coordinate_system<pm::Vertex> {
108 using type = cs::cartesian;
109};
110template <>
111struct dimension<pm::Vertex> : boost::mpl::int_<3> {
112};
113
114template <size_t Dimension>
115struct access<pm::Vertex, Dimension> {
116 static double get(pm::Vertex const &p)
117 {
118 return p.rawCoords()[Dimension];
119 }
120
121 static void set(pm::Vertex &p, double const &value)
122 {
123 PRECICE_UNREACHABLE("Boost.Geometry is not allowed to write to mesh::Vertex.");
124 }
125};
126
127BOOST_CONCEPT_ASSERT((concepts::Point<pm::Vertex>) );
128
134template <>
135struct tag<pm::Edge> {
136 using type = segment_tag;
137};
138template <>
139struct point_type<pm::Edge> {
140 using type = pm::Vertex::RawCoords;
141};
142
143template <size_t Index, size_t Dimension>
144struct indexed_access<pm::Edge, Index, Dimension> {
145 static_assert((Index <= 1), "Valid Indices are {0, 1}");
146 static_assert((Dimension <= 2), "Valid Dimensions are {0, 1, 2}");
147
148 static double get(pm::Edge const &e)
149 {
150 return access<point_type<pm::Edge>::type, Dimension>::get(e.vertex(Index).rawCoords());
151 }
152
153 static void set(pm::Edge &e, double const &value)
154 {
155 PRECICE_UNREACHABLE("Boost.Geometry is not allowed to write to mesh::Edge.");
156 }
157};
158
164template <>
165struct tag<pm::Triangle> {
166 using type = ring_tag;
167};
168template <>
169struct point_order<pm::Triangle> {
170 static const order_selector value = clockwise;
171};
172template <>
173struct closure<pm::Triangle> {
174 static const closure_selector value = open;
175};
176
177// BOOST_CONCEPT_ASSERT( (bg::concepts::Ring<pm::Triangle>));
178
179} // namespace boost::geometry::traits
180
181namespace precice::query {
182
184using RTreeBox = boost::geometry::model::box<pm::Vertex::RawCoords>;
185
186inline RTreeBox makeBox(const pm::Vertex::RawCoords &min, const pm::Vertex::RawCoords &max)
187{
188 return {min, max};
189}
190
191inline pm::Vertex::RawCoords eigenToRaw(const Eigen::VectorXd &v)
192{
193 const auto size = v.size();
194 PRECICE_ASSERT(size == 2 || size == 3, size);
195 pm::Vertex::RawCoords r{0.0, 0.0, 0.0};
196 std::copy_n(v.data(), size, r.data());
197 return r;
198}
199
200inline Eigen::VectorXd rawToEigen(const pm::Vertex::RawCoords &v)
201{
202 Eigen::VectorXd r(3);
203 std::copy(v.begin(), v.end(), r.data());
204 return r;
205}
206
207inline RTreeBox makeBox(const Eigen::VectorXd &min, const Eigen::VectorXd &max)
208{
209 return {eigenToRaw(min), eigenToRaw(max)};
210}
211
212// Overload for a tetrahedron
214{
215
217 for (int i = 0; i < 4; ++i) {
218 box.expandBy(tetra.vertex(i));
219 }
220
221 // Convert to Boost type
222 return makeBox(box.minCorner(), box.maxCorner());
223}
224
225namespace impl {
226
228using RTreeParameters = boost::geometry::index::rstar<16>;
229
231template <class T>
233
234template <>
238
239template <>
243
244template <>
248
249template <>
250struct PrimitiveTraits<mesh::Tetrahedron> {
252};
253
255template <typename Container>
257 using size_type = typename Container::container::size_type;
258 using cref = const typename Container::value_type &;
259 Container const &container;
260
261public:
263
264 explicit PtrVectorIndexable(Container const &c)
265 : container(c)
266 {
267 }
268
270 {
271 return container[i];
272 }
273};
274
276template <typename Container>
278 using size_type = typename Container::size_type;
279 using cref = const typename Container::value_type &;
280 Container const &container;
281
282public:
284
285 explicit VectorIndexable(Container const &c)
286 : container(c)
287 {
288 }
289
291 {
292 return container[i];
293 }
294};
295
296template <typename Primitive>
298private:
299 template <typename T, typename = typename std::enable_if<
301 typename boost::geometry::traits::tag<T>::type,
302 boost::geometry::point_tag>::value,
304 static std::true_type test(char *);
305 template <typename T, typename = typename std::enable_if<
307 typename boost::geometry::traits::tag<T>::type,
308 boost::geometry::segment_tag>::value,
310 static std::true_type test(int *);
311 template <typename T, typename = typename std::enable_if<
313 typename boost::geometry::traits::tag<T>::type,
314 boost::geometry::box_tag>::value,
316 static std::true_type test(void *);
317
318 template <typename T>
320
321public:
322 using type = decltype(test<Primitive>(nullptr));
323};
324
325template <class Primitive>
327};
328
330template <class Primitive>
333 using MeshContainerIndex = typename MeshContainer::size_type;
334
335 using IndexType = typename std::conditional<
339
343 boost::geometry::index::indexable<IndexType>>::type;
344
345 using RTree = boost::geometry::index::rtree<IndexType, RTreeParameters, IndexGetter>;
347};
348
349} // namespace impl
350} // namespace precice::query
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
#define PRECICE_UNREACHABLE(...)
Definition assertion.hpp:93
An axis-aligned bounding box around a (partition of a) mesh.
Eigen::VectorXd maxCorner() const
the max corner of the bounding box
Eigen::VectorXd minCorner() const
the min corner of the bounding box
void expandBy(const BoundingBox &otherBB)
Expand bounding box using another bounding box.
Linear edge of a mesh, defined by two Vertex objects.
Definition Edge.hpp:15
Vertex & vertex(int i)
Returns the edge's vertex with index 0 or 1.
Definition Edge.hpp:76
std::deque< Triangle > TriangleContainer
Definition Mesh.hpp:42
std::deque< Tetrahedron > TetraContainer
Definition Mesh.hpp:43
std::deque< Edge > EdgeContainer
Definition Mesh.hpp:41
std::deque< Vertex > VertexContainer
Definition Mesh.hpp:40
Tetrahedron of a mesh, defined by 4 vertices.
Vertex & vertex(int i)
Returns tetrahedron vertex with index 0, 1, 2 or 3.
int getDimensions() const
Returns dimensionalty of space the Tetrahedron is embedded in.
Vertex of a mesh.
Definition Vertex.hpp:16
const RawCoords & rawCoords() const
Direct access to the coordinates.
Definition Vertex.hpp:121
static std::true_type test(void *)
static std::true_type test(char *)
decltype(test< Primitive >(nullptr)) type
Makes a utils::PtrVector indexable and thus be usable in boost::geometry::rtree.
typename Container::container::size_type size_type
const typename Container::value_type & cref
result_type operator()(size_type i) const
Makes a std::vector indexable and thus be usable in boost::geometry::rtree.
const typename Container::value_type & cref
result_type operator()(size_type i) const
typename Container::size_type size_type
T copy(T... args)
T copy_n(T... args)
BOOST_CONCEPT_ASSERT((bg::concepts::Point< pm::Vertex::RawCoords >))
provides Mesh, Data and primitives.
boost::geometry::index::rstar< 16 > RTreeParameters
The general rtree parameter type used in precice.
contains geometrical queries.
RTreeBox makeBox(const pm::Vertex::RawCoords &min, const pm::Vertex::RawCoords &max)
pm::Vertex::RawCoords eigenToRaw(const Eigen::VectorXd &v)
Eigen::VectorXd rawToEigen(const pm::Vertex::RawCoords &v)
boost::geometry::model::box< pm::Vertex::RawCoords > RTreeBox
The RTree box type.
static std::unique_ptr< precice::Participant > impl
Definition preciceC.cpp:21
static void set(Eigen::VectorXd &p, double const &value)
static void set(pm::Vertex &p, double const &value)
static void set(pm::Vertex::RawCoords &p, double const &value)
Type trait to extract information based on the type of a Primitive.
The type traits of a rtree based on a Primitive.
typename PrimitiveTraits< Primitive >::MeshContainer MeshContainer
typename std::conditional< IsDirectIndexable< Primitive >::value, MeshContainerIndex, std::pair< RTreeBox, MeshContainerIndex > >::type IndexType
typename MeshContainer::size_type MeshContainerIndex
typename std::conditional< IsDirectIndexable< Primitive >::value, impl::VectorIndexable< MeshContainer >, boost::geometry::index::indexable< IndexType > >::type IndexGetter
boost::geometry::index::rtree< IndexType, RTreeParameters, IndexGetter > RTree