22#ifndef DUMUX_INTERSECTS_POINT_SIMPLEX_HH
23#define DUMUX_INTERSECTS_POINT_SIMPLEX_HH
26#include <dune/common/fvector.hh>
35template<
class ctype,
int dimworld,
typename std::enable_if_t<(dimworld == 3),
int> = 0>
37 const Dune::FieldVector<ctype, dimworld>& p0,
38 const Dune::FieldVector<ctype, dimworld>& p1,
39 const Dune::FieldVector<ctype, dimworld>& p2,
40 const Dune::FieldVector<ctype, dimworld>& p3)
44 using GlobalPosition = Dune::FieldVector<ctype, dimworld>;
45 static constexpr ctype eps_ = 1.0e-7;
48 const GlobalPosition *p[4] = {&p0, &p1, &p2, &p3};
51 for (
int i = 0; i < 4; ++i)
54 const GlobalPosition v1 = *p[(i + 1)%4] - *p[i];
55 const GlobalPosition v2 = *p[(i + 2)%4] - *p[i];
56 const GlobalPosition v3 = *p[(i + 3)%4] - *p[i];
57 const GlobalPosition v = point - *p[i];
62 const auto t1 = n1.dot(v);
63 const auto t2 = n1.dot(v3);
66 const auto eps = eps_ * v1.two_norm();
67 if ((t1 > eps || t1 < -eps) && std::signbit(t1) != std::signbit(t2))
77template<
class ctype,
int dimworld,
typename std::enable_if_t<(dimworld == 3),
int> = 0>
79 const Dune::FieldVector<ctype, dimworld>& p0,
80 const Dune::FieldVector<ctype, dimworld>& p1,
81 const Dune::FieldVector<ctype, dimworld>& p2)
85 constexpr ctype eps_ = 1.0e-7;
88 const auto v1 = p0 - p2;
90 const ctype nnorm = n.two_norm();
91 const ctype eps4 = eps_*nnorm*nnorm;
104 const auto a = p0 - point;
105 const auto b = p1 - point;
106 const auto c = p2 - point;
113 if (u*v < 0.0 - eps4)
120 if (u*w < 0.0 - eps4)
124 if (u.two_norm2() < eps4)
125 return b*c < 0.0 + eps_*nnorm;
126 if (v.two_norm2() < eps4)
127 return a*c < 0.0 + eps_*nnorm;
137template<
class ctype,
int dimworld,
typename std::enable_if_t<(dimworld == 2),
int> = 0>
139 const Dune::FieldVector<ctype, dimworld>& p0,
140 const Dune::FieldVector<ctype, dimworld>& p1,
141 const Dune::FieldVector<ctype, dimworld>& p2)
143 static constexpr ctype eps_ = 1.0e-7;
146 const ctype A = 0.5*(-p1[1]*p2[0] + p0[1]*(p2[0] - p1[0])
147 +p1[0]*p2[1] + p0[0]*(p1[1] - p2[1]));
148 const ctype
sign = std::copysign(1.0, A);
149 const ctype s =
sign*(p0[1]*p2[0] + point[0]*(p2[1]-p0[1])
150 -p0[0]*p2[1] + point[1]*(p0[0]-p2[0]));
151 const ctype t =
sign*(p0[0]*p1[1] + point[0]*(p0[1]-p1[1])
152 -p0[1]*p1[0] + point[1]*(p1[0]-p0[0]));
153 const ctype eps =
sign*A*eps_;
157 && (s + t) < 2*A*
sign + eps);
165template<
class ctype,
int dimworld,
typename std::enable_if_t<(dimworld == 3 || dimworld == 2),
int> = 0>
167 const Dune::FieldVector<ctype, dimworld>& p0,
168 const Dune::FieldVector<ctype, dimworld>& p1)
170 using GlobalPosition = Dune::FieldVector<ctype, dimworld>;
171 static constexpr ctype eps_ = 1.0e-7;
174 const GlobalPosition v1 = p1 - p0;
175 const GlobalPosition v2 = point - p0;
177 const ctype v1norm = v1.two_norm();
178 const ctype v2norm = v2.two_norm();
181 if (v2norm < v1norm*eps_)
184 return (v1.dot(v2) > v1norm*v2norm*(1.0 - eps_) && v2norm < v1norm*(1.0 + eps_));
191template<
class ctype,
int dimworld,
typename std::enable_if_t<(dimworld == 1),
int> = 0>
193 const Dune::FieldVector<ctype, dimworld>& p0,
194 const Dune::FieldVector<ctype, dimworld>& p1)
196 static constexpr ctype eps_ = 1.0e-7;
199 const ctype *interval[2] = {&p0[0], &p1[0]};
200 if (*interval[0] > *interval[1])
201 std::swap(interval[0], interval[1]);
203 const ctype v1 = point[0] - *interval[0];
204 const ctype v2 = *interval[1] - *interval[0];
208 if (abs(v1) < v2*eps_)
220 return (!signbit(v1) && abs(v1) < v2*(1.0 + eps_));
Define some often used mathematical functions.
Dune::FieldVector< Scalar, 3 > crossProduct(const Dune::FieldVector< Scalar, 3 > &vec1, const Dune::FieldVector< Scalar, 3 > &vec2)
Cross product of two vectors in three-dimensional Euclidean space.
Definition: math.hh:631
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition: math.hh:618
bool intersectsPointSimplex(const Dune::FieldVector< ctype, dimworld > &point, const Dune::FieldVector< ctype, dimworld > &p0, const Dune::FieldVector< ctype, dimworld > &p1, const Dune::FieldVector< ctype, dimworld > &p2, const Dune::FieldVector< ctype, dimworld > &p3)
Find out whether a point is inside a tetrahedron (p0, p1, p2, p3) (dimworld is 3)
Definition: intersectspointsimplex.hh:36