4#include <boost/geometry.hpp>
29struct tag<Eigen::VectorXd> {
33struct coordinate_type<Eigen::VectorXd> {
37struct coordinate_system<Eigen::VectorXd> {
38 using type = cs::cartesian;
41struct dimension<Eigen::VectorXd> : boost::mpl::int_<3> {
44template <
size_t Dimension>
45struct access<Eigen::VectorXd, Dimension> {
46 static double get(Eigen::VectorXd
const &p)
48 if (Dimension >=
static_cast<size_t>(p.rows())) {
55 static void set(Eigen::VectorXd &p,
double const &value)
57 if (Dimension >=
static_cast<size_t>(p.rows())) {
58 Eigen::VectorXd tmp(Dimension);
68struct tag<
pm::Vertex::RawCoords> {
72struct coordinate_type<
pm::Vertex::RawCoords> {
76struct coordinate_system<
pm::Vertex::RawCoords> {
77 using type = cs::cartesian;
80struct dimension<
pm::Vertex::RawCoords> : boost::mpl::int_<3> {
83template <
size_t Dimension>
84struct access<
pm::Vertex::RawCoords, Dimension> {
85 static double get(pm::Vertex::RawCoords
const &p)
90 static void set(pm::Vertex::RawCoords &p,
double const &value)
103struct tag<
pm::Vertex> {
107struct coordinate_type<
pm::Vertex> {
111struct coordinate_system<
pm::Vertex> {
115struct dimension<
pm::Vertex> : boost::mpl::int_<3> {
118template <
size_t Dimension>
119struct access<
pm::Vertex, Dimension> {
139struct tag<
pm::Edge> {
143struct point_type<
pm::Edge> {
144 using type = pm::Vertex::RawCoords;
147template <
size_t Index,
size_t Dimension>
148struct indexed_access<
pm::Edge, Index, Dimension> {
149 static_assert((Index <= 1),
"Valid Indices are {0, 1}");
150 static_assert((Dimension <= 2),
"Valid Dimensions are {0, 1, 2}");
154 return access<point_type<pm::Edge>::type, Dimension>::get(e.
vertex(Index).
rawCoords());
169struct tag<
pm::Triangle> {
173struct point_order<
pm::Triangle> {
174 static const order_selector value = clockwise;
177struct closure<
pm::Triangle> {
178 static const closure_selector value = open;
191using RTreeBox = boost::geometry::model::box<pm::Vertex::RawCoords>;
193inline RTreeBox makeBox(
const pm::Vertex::RawCoords &min,
const pm::Vertex::RawCoords &max)
198inline pm::Vertex::RawCoords
eigenToRaw(
const Eigen::VectorXd &v)
200 const auto size = v.size();
202 pm::Vertex::RawCoords r{0.0, 0.0, 0.0};
207inline Eigen::VectorXd
rawToEigen(
const pm::Vertex::RawCoords &v)
209 Eigen::VectorXd r(3);
224 for (
int i = 0; i < 4; ++i) {
262template <
typename Container>
264 using size_type =
typename Container::container::size_type;
265 using cref =
const typename Container::value_type &;
283template <
typename Container>
286 using cref =
const typename Container::value_type &;
303template <
typename Primitive>
308 typename boost::geometry::traits::tag<T>::type,
309 boost::geometry::point_tag>::value,
314 typename boost::geometry::traits::tag<T>::type,
315 boost::geometry::segment_tag>::value,
320 typename boost::geometry::traits::tag<T>::type,
321 boost::geometry::box_tag>::value,
325 template <
typename T>
329 using type =
decltype(test<Primitive>(
nullptr));
332template <
class Primitive>
337template <
class Primitive>
350 boost::geometry::index::indexable<IndexType>>::type;
352 using RTree = boost::geometry::index::rtree<IndexType, RTreeParameters, IndexGetter>;
#define PRECICE_ASSERT(...)
#define PRECICE_UNREACHABLE(...)
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.
Vertex & vertex(int i)
Returns the edge's vertex with index 0 or 1.
std::deque< Triangle > TriangleContainer
std::deque< Tetrahedron > TetraContainer
std::deque< Edge > EdgeContainer
std::deque< Vertex > VertexContainer
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.
const RawCoords & rawCoords() const
Direct access to the coordinates.
static std::true_type test(void *)
static std::true_type test(int *)
static std::false_type test(...)
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
PtrVectorIndexable(Container const &c)
Container const & container
result_type operator()(size_type i) const
Makes a std::vector indexable and thus be usable in boost::geometry::rtree.
VectorIndexable(Container const &c)
const typename Container::value_type & cref
result_type operator()(size_type i) const
typename Container::size_type size_type
Container const & container
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.
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.
Main namespace of the precice library.
static std::unique_ptr< precice::Participant > impl
static void set(Eigen::VectorXd &p, double const &value)
static double get(Eigen::VectorXd const &p)
static void set(pm::Vertex &p, double const &value)
static double get(pm::Vertex const &p)
static double get(pm::Vertex::RawCoords const &p)
static void set(pm::Vertex::RawCoords &p, double const &value)
static double get(pm::Edge const &e)
static void set(pm::Edge &e, double const &value)
pm::Vertex::RawCoords type
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