3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
intersectspointsimplex.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_GEOMETRY_INTERSECTS_POINT_SIMPLEX_HH
25#define DUMUX_GEOMETRY_INTERSECTS_POINT_SIMPLEX_HH
26
27#include <cmath>
28#include <dune/common/fvector.hh>
29#include <dumux/common/math.hh>
30
31namespace Dumux {
32
37template<class ctype, int dimworld, typename std::enable_if_t<(dimworld == 3), int> = 0>
38bool intersectsPointSimplex(const Dune::FieldVector<ctype, dimworld>& point,
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)
43{
44 // Algorithm from http://www.blackpawn.com/texts/pointinpoly/
45 // See also "Real-Time Collision Detection" by Christer Ericson.
46 using GlobalPosition = Dune::FieldVector<ctype, dimworld>;
47 static constexpr ctype eps_ = 1.0e-7;
48
49 // put the tetrahedron points in an array
50 const GlobalPosition *p[4] = {&p0, &p1, &p2, &p3};
51
52 // iterate over all faces
53 for (int i = 0; i < 4; ++i)
54 {
55 // compute all the vectors from vertex (local index 0) to the other points
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];
60 // compute the normal to the facet (cross product)
61 GlobalPosition n1 = crossProduct(v1, v2);
62 n1 /= n1.two_norm();
63 // find out on which side of the plane v and v3 are
64 const auto t1 = n1.dot(v);
65 const auto t2 = n1.dot(v3);
66 // If the point is not exactly on the plane the
67 // points have to be on the same side
68 const auto eps = eps_ * v1.two_norm();
69 if ((t1 > eps || t1 < -eps) && std::signbit(t1) != std::signbit(t2))
70 return false;
71 }
72 return true;
73}
74
79template<class ctype, int dimworld, typename std::enable_if_t<(dimworld == 3), int> = 0>
80bool intersectsPointSimplex(const Dune::FieldVector<ctype, dimworld>& point,
81 const Dune::FieldVector<ctype, dimworld>& p0,
82 const Dune::FieldVector<ctype, dimworld>& p1,
83 const Dune::FieldVector<ctype, dimworld>& p2)
84{
85 // adapted from the algorithm from from "Real-Time Collision Detection" by Christer Ericson,
86 // published by Morgan Kaufmann Publishers, (c) 2005 Elsevier Inc. (Chapter 5.4.2)
87 constexpr ctype eps_ = 1.0e-7;
88
89 // compute the normal of the triangle
90 const auto v1 = p0 - p2;
91 auto n = crossProduct(v1, p1 - p0);
92 const ctype nnorm = n.two_norm();
93 const ctype eps4 = eps_*nnorm*nnorm; // compute an epsilon for later
94 n /= nnorm; // normalize
95
96 // first check if we are in the plane of the triangle
97 // if not we can return early
98 using std::abs;
99 auto x = p0 - point;
100 x /= x.two_norm(); // normalize
101
102 if (abs(x*n) > eps_)
103 return false;
104
105 // translate the triangle so that 'point' is the origin
106 const auto a = p0 - point;
107 const auto b = p1 - point;
108 const auto c = p2 - point;
109
110 // compute the normal vectors for triangles P->A->B and P->B->C
111 const auto u = crossProduct(b, c);
112 const auto v = crossProduct(c, a);
113
114 // they have to point in the same direction or be orthogonal
115 if (u*v < 0.0 - eps4)
116 return false;
117
118 // compute the normal vector for triangle P->C->A
119 const auto w = crossProduct(a, b);
120
121 // it also has to point in the same direction or be orthogonal
122 if (u*w < 0.0 - eps4)
123 return false;
124
125 // check if point is on the line of one of the edges
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;
130
131 // now the point must be in the triangle (or on the faces)
132 return true;
133}
134
139template<class ctype, int dimworld, typename std::enable_if_t<(dimworld == 2), int> = 0>
140bool intersectsPointSimplex(const Dune::FieldVector<ctype, dimworld>& point,
141 const Dune::FieldVector<ctype, dimworld>& p0,
142 const Dune::FieldVector<ctype, dimworld>& p1,
143 const Dune::FieldVector<ctype, dimworld>& p2)
144{
145 static constexpr ctype eps_ = 1.0e-7;
146
147 // Use barycentric coordinates
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_;
156
157 return (s > -eps
158 && t > -eps
159 && (s + t) < 2*A*sign + eps);
160}
161
167template<class ctype, int dimworld, typename std::enable_if_t<(dimworld == 3 || dimworld == 2), int> = 0>
168bool intersectsPointSimplex(const Dune::FieldVector<ctype, dimworld>& point,
169 const Dune::FieldVector<ctype, dimworld>& p0,
170 const Dune::FieldVector<ctype, dimworld>& p1)
171{
172 using GlobalPosition = Dune::FieldVector<ctype, dimworld>;
173 static constexpr ctype eps_ = 1.0e-7;
174
175 // compute the vectors between p0 and the other points
176 const GlobalPosition v1 = p1 - p0;
177 const GlobalPosition v2 = point - p0;
178
179 const ctype v1norm = v1.two_norm();
180 const ctype v2norm = v2.two_norm();
181
182 // early exit if point and p0 are the same
183 if (v2norm < v1norm*eps_)
184 return true;
185
186 // early exit if the point is outside the segment (no epsilon in the
187 // first statement because we already did the above equality check)
188 if (v1.dot(v2) < 0.0 || v2norm > v1norm*(1.0 + eps_))
189 return false;
190
191 // If the area spanned by the 2 vectors is zero, the points are colinear.
192 // If that is the case, the given point is on the segment.
193 const auto n = crossProduct(v1, v2);
194 const auto eps2 = v1norm*v1norm*eps_;
195 if constexpr (dimworld == 3)
196 return n.two_norm2() < eps2*eps2;
197 else
198 {
199 using std::abs;
200 return abs(n) < eps2;
201 }
202}
203
208template<class ctype, int dimworld, typename std::enable_if_t<(dimworld == 1), int> = 0>
209bool intersectsPointSimplex(const Dune::FieldVector<ctype, dimworld>& point,
210 const Dune::FieldVector<ctype, dimworld>& p0,
211 const Dune::FieldVector<ctype, dimworld>& p1)
212{
213 static constexpr ctype eps_ = 1.0e-7;
214
215 // sort the interval so interval[1] is the end and interval[0] the start
216 const ctype *interval[2] = {&p0[0], &p1[0]};
217 if (*interval[0] > *interval[1])
218 std::swap(interval[0], interval[1]);
219
220 const ctype v1 = point[0] - *interval[0];
221 const ctype v2 = *interval[1] - *interval[0]; // always positive
222
223 // the coordinates are the same
224 using std::abs;
225 if (abs(v1) < v2*eps_)
226 return true;
227
228 // the point doesn't coincide with p0
229 // so if p0 and p1 are equal it's not inside
230 if (v2 < 1.0e-30)
231 return false;
232
233 // the point is inside if the length is
234 // smaller than the interval length and the
235 // sign of v1 & v2 are the same
236 using std::signbit;
237 return (!signbit(v1) && abs(v1) < v2*(1.0 + eps_));
238}
239
240} // end namespace Dumux
241
242#endif
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