24#ifndef DUMUX_GEOMETRY_INTERSECTS_POINT_SIMPLEX_HH
25#define DUMUX_GEOMETRY_INTERSECTS_POINT_SIMPLEX_HH
28#include <dune/common/fvector.hh>
37template<
class ctype,
int dimworld,
typename std::enable_if_t<(dimworld == 3),
int> = 0>
39 const Dune::FieldVector<ctype, dimworld>& p0,
40 const Dune::FieldVector<ctype, dimworld>& p1,
41 const Dune::FieldVector<ctype, dimworld>& p2,
42 const Dune::FieldVector<ctype, dimworld>& p3)
46 using GlobalPosition = Dune::FieldVector<ctype, dimworld>;
47 static constexpr ctype eps_ = 1.0e-7;
50 const GlobalPosition *p[4] = {&p0, &p1, &p2, &p3};
53 for (
int i = 0; i < 4; ++i)
56 const GlobalPosition v1 = *p[(i + 1)%4] - *p[i];
57 const GlobalPosition v2 = *p[(i + 2)%4] - *p[i];
58 const GlobalPosition v3 = *p[(i + 3)%4] - *p[i];
59 const GlobalPosition v = point - *p[i];
64 const auto t1 = n1.dot(v);
65 const auto t2 = n1.dot(v3);
68 const auto eps = eps_ * v1.two_norm();
69 if ((t1 > eps || t1 < -eps) && std::signbit(t1) != std::signbit(t2))
79template<
class ctype,
int dimworld,
typename std::enable_if_t<(dimworld == 3),
int> = 0>
81 const Dune::FieldVector<ctype, dimworld>& p0,
82 const Dune::FieldVector<ctype, dimworld>& p1,
83 const Dune::FieldVector<ctype, dimworld>& p2)
87 constexpr ctype eps_ = 1.0e-7;
90 const auto v1 = p0 - p2;
92 const ctype nnorm = n.two_norm();
93 const ctype eps4 = eps_*nnorm*nnorm;
106 const auto a = p0 - point;
107 const auto b = p1 - point;
108 const auto c = p2 - point;
115 if (u*v < 0.0 - eps4)
122 if (u*w < 0.0 - eps4)
126 if (u.two_norm2() < eps4)
127 return b*c < 0.0 + eps_*nnorm;
128 if (v.two_norm2() < eps4)
129 return a*c < 0.0 + eps_*nnorm;
139template<
class ctype,
int dimworld,
typename std::enable_if_t<(dimworld == 2),
int> = 0>
141 const Dune::FieldVector<ctype, dimworld>& p0,
142 const Dune::FieldVector<ctype, dimworld>& p1,
143 const Dune::FieldVector<ctype, dimworld>& p2)
145 static constexpr ctype eps_ = 1.0e-7;
148 const ctype A = 0.5*(-p1[1]*p2[0] + p0[1]*(p2[0] - p1[0])
149 +p1[0]*p2[1] + p0[0]*(p1[1] - p2[1]));
150 const ctype
sign = std::copysign(1.0, A);
151 const ctype s =
sign*(p0[1]*p2[0] + point[0]*(p2[1]-p0[1])
152 -p0[0]*p2[1] + point[1]*(p0[0]-p2[0]));
153 const ctype t =
sign*(p0[0]*p1[1] + point[0]*(p0[1]-p1[1])
154 -p0[1]*p1[0] + point[1]*(p1[0]-p0[0]));
155 const ctype eps =
sign*A*eps_;
159 && (s + t) < 2*A*
sign + eps);
167template<
class ctype,
int dimworld,
typename std::enable_if_t<(dimworld == 3 || dimworld == 2),
int> = 0>
169 const Dune::FieldVector<ctype, dimworld>& p0,
170 const Dune::FieldVector<ctype, dimworld>& p1)
172 using GlobalPosition = Dune::FieldVector<ctype, dimworld>;
173 static constexpr ctype eps_ = 1.0e-7;
176 const GlobalPosition v1 = p1 - p0;
177 const GlobalPosition v2 = point - p0;
179 const ctype v1norm = v1.two_norm();
180 const ctype v2norm = v2.two_norm();
183 if (v2norm < v1norm*eps_)
188 if (v1.dot(v2) < 0.0 || v2norm > v1norm*(1.0 + eps_))
194 const auto eps2 = v1norm*v1norm*eps_;
195 if constexpr (dimworld == 3)
196 return n.two_norm2() < eps2*eps2;
200 return abs(n) < eps2;
208template<
class ctype,
int dimworld,
typename std::enable_if_t<(dimworld == 1),
int> = 0>
210 const Dune::FieldVector<ctype, dimworld>& p0,
211 const Dune::FieldVector<ctype, dimworld>& p1)
213 static constexpr ctype eps_ = 1.0e-7;
216 const ctype *interval[2] = {&p0[0], &p1[0]};
217 if (*interval[0] > *interval[1])
218 std::swap(interval[0], interval[1]);
220 const ctype v1 = point[0] - *interval[0];
221 const ctype v2 = *interval[1] - *interval[0];
225 if (abs(v1) < v2*eps_)
237 return (!signbit(v1) && abs(v1) < v2*(1.0 + eps_));
Define some often used mathematical functions.
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 the tetrahedron (p0, p1, p2, p3) (dimworld is 3)
Definition: intersectspointsimplex.hh:38
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:654
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition: math.hh:641
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29