1
// This file is part of Eigen, a lightweight C++ template library
4
// Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
5
// Copyright (C) 2008 Benoit Jacob <jacob.benoit.1@gmail.com>
7
// Eigen is free software; you can redistribute it and/or
8
// modify it under the terms of the GNU Lesser General Public
9
// License as published by the Free Software Foundation; either
10
// version 3 of the License, or (at your option) any later version.
12
// Alternatively, you can redistribute it and/or
13
// modify it under the terms of the GNU General Public License as
14
// published by the Free Software Foundation; either version 2 of
15
// the License, or (at your option) any later version.
17
// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
18
// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
19
// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
20
// GNU General Public License for more details.
22
// You should have received a copy of the GNU Lesser General Public
23
// License and a copy of the GNU General Public License along with
24
// Eigen. If not, see <http://www.gnu.org/licenses/>.
26
#ifndef EIGEN_HYPERPLANE_H
27
#define EIGEN_HYPERPLANE_H
29
/** \geometry_module \ingroup Geometry_Module
35
* A hyperplane is an affine subspace of dimension n-1 in a space of dimension n.
36
* For example, a hyperplane in a plane is a line; a hyperplane in 3-space is a plane.
38
* \param _Scalar the scalar type, i.e., the type of the coefficients
39
* \param _AmbientDim the dimension of the ambient space, can be a compile time value or Dynamic.
40
* Notice that the dimension of the hyperplane is _AmbientDim-1.
42
* This class represents an hyperplane as the zero set of the implicit equation
43
* \f$ n \cdot x + d = 0 \f$ where \f$ n \f$ is a unit normal vector of the plane (linear part)
44
* and \f$ d \f$ is the distance (offset) to the origin.
46
template <typename _Scalar, int _AmbientDim, int _Options>
50
EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(_Scalar,_AmbientDim==Dynamic ? Dynamic : _AmbientDim+1)
52
AmbientDimAtCompileTime = _AmbientDim,
55
typedef _Scalar Scalar;
56
typedef typename NumTraits<Scalar>::Real RealScalar;
57
typedef DenseIndex Index;
58
typedef Matrix<Scalar,AmbientDimAtCompileTime,1> VectorType;
59
typedef Matrix<Scalar,Index(AmbientDimAtCompileTime)==Dynamic
61
: Index(AmbientDimAtCompileTime)+1,1,Options> Coefficients;
62
typedef Block<Coefficients,AmbientDimAtCompileTime,1> NormalReturnType;
63
typedef const Block<const Coefficients,AmbientDimAtCompileTime,1> ConstNormalReturnType;
65
/** Default constructor without initialization */
66
inline explicit Hyperplane() {}
68
template<int OtherOptions>
69
Hyperplane(const Hyperplane<Scalar,AmbientDimAtCompileTime,OtherOptions>& other)
70
: m_coeffs(other.coeffs())
73
/** Constructs a dynamic-size hyperplane with \a _dim the dimension
74
* of the ambient space */
75
inline explicit Hyperplane(Index _dim) : m_coeffs(_dim+1) {}
77
/** Construct a plane from its normal \a n and a point \a e onto the plane.
78
* \warning the vector normal is assumed to be normalized.
80
inline Hyperplane(const VectorType& n, const VectorType& e)
81
: m_coeffs(n.size()+1)
87
/** Constructs a plane from its normal \a n and distance to the origin \a d
88
* such that the algebraic equation of the plane is \f$ n \cdot x + d = 0 \f$.
89
* \warning the vector normal is assumed to be normalized.
91
inline Hyperplane(const VectorType& n, Scalar d)
92
: m_coeffs(n.size()+1)
98
/** Constructs a hyperplane passing through the two points. If the dimension of the ambient space
99
* is greater than 2, then there isn't uniqueness, so an arbitrary choice is made.
101
static inline Hyperplane Through(const VectorType& p0, const VectorType& p1)
103
Hyperplane result(p0.size());
104
result.normal() = (p1 - p0).unitOrthogonal();
105
result.offset() = -p0.dot(result.normal());
109
/** Constructs a hyperplane passing through the three points. The dimension of the ambient space
110
* is required to be exactly 3.
112
static inline Hyperplane Through(const VectorType& p0, const VectorType& p1, const VectorType& p2)
114
EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 3)
115
Hyperplane result(p0.size());
116
result.normal() = (p2 - p0).cross(p1 - p0).normalized();
117
result.offset() = -p0.dot(result.normal());
121
/** Constructs a hyperplane passing through the parametrized line \a parametrized.
122
* If the dimension of the ambient space is greater than 2, then there isn't uniqueness,
123
* so an arbitrary choice is made.
125
// FIXME to be consitent with the rest this could be implemented as a static Through function ??
126
explicit Hyperplane(const ParametrizedLine<Scalar, AmbientDimAtCompileTime>& parametrized)
128
normal() = parametrized.direction().unitOrthogonal();
129
offset() = -parametrized.origin().dot(normal());
134
/** \returns the dimension in which the plane holds */
135
inline Index dim() const { return AmbientDimAtCompileTime==Dynamic ? m_coeffs.size()-1 : Index(AmbientDimAtCompileTime); }
137
/** normalizes \c *this */
140
m_coeffs /= normal().norm();
143
/** \returns the signed distance between the plane \c *this and a point \a p.
146
inline Scalar signedDistance(const VectorType& p) const { return normal().dot(p) + offset(); }
148
/** \returns the absolute distance between the plane \c *this and a point \a p.
149
* \sa signedDistance()
151
inline Scalar absDistance(const VectorType& p) const { return internal::abs(signedDistance(p)); }
153
/** \returns the projection of a point \a p onto the plane \c *this.
155
inline VectorType projection(const VectorType& p) const { return p - signedDistance(p) * normal(); }
157
/** \returns a constant reference to the unit normal vector of the plane, which corresponds
158
* to the linear part of the implicit equation.
160
inline ConstNormalReturnType normal() const { return ConstNormalReturnType(m_coeffs,0,0,dim(),1); }
162
/** \returns a non-constant reference to the unit normal vector of the plane, which corresponds
163
* to the linear part of the implicit equation.
165
inline NormalReturnType normal() { return NormalReturnType(m_coeffs,0,0,dim(),1); }
167
/** \returns the distance to the origin, which is also the "constant term" of the implicit equation
168
* \warning the vector normal is assumed to be normalized.
170
inline const Scalar& offset() const { return m_coeffs.coeff(dim()); }
172
/** \returns a non-constant reference to the distance to the origin, which is also the constant part
173
* of the implicit equation */
174
inline Scalar& offset() { return m_coeffs(dim()); }
176
/** \returns a constant reference to the coefficients c_i of the plane equation:
177
* \f$ c_0*x_0 + ... + c_{d-1}*x_{d-1} + c_d = 0 \f$
179
inline const Coefficients& coeffs() const { return m_coeffs; }
181
/** \returns a non-constant reference to the coefficients c_i of the plane equation:
182
* \f$ c_0*x_0 + ... + c_{d-1}*x_{d-1} + c_d = 0 \f$
184
inline Coefficients& coeffs() { return m_coeffs; }
186
/** \returns the intersection of *this with \a other.
188
* \warning The ambient space must be a plane, i.e. have dimension 2, so that \c *this and \a other are lines.
190
* \note If \a other is approximately parallel to *this, this method will return any point on *this.
192
VectorType intersection(const Hyperplane& other) const
194
EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 2)
195
Scalar det = coeffs().coeff(0) * other.coeffs().coeff(1) - coeffs().coeff(1) * other.coeffs().coeff(0);
196
// since the line equations ax+by=c are normalized with a^2+b^2=1, the following tests
197
// whether the two lines are approximately parallel.
198
if(internal::isMuchSmallerThan(det, Scalar(1)))
199
{ // special case where the two lines are approximately parallel. Pick any point on the first line.
200
if(internal::abs(coeffs().coeff(1))>internal::abs(coeffs().coeff(0)))
201
return VectorType(coeffs().coeff(1), -coeffs().coeff(2)/coeffs().coeff(1)-coeffs().coeff(0));
203
return VectorType(-coeffs().coeff(2)/coeffs().coeff(0)-coeffs().coeff(1), coeffs().coeff(0));
207
Scalar invdet = Scalar(1) / det;
208
return VectorType(invdet*(coeffs().coeff(1)*other.coeffs().coeff(2)-other.coeffs().coeff(1)*coeffs().coeff(2)),
209
invdet*(other.coeffs().coeff(0)*coeffs().coeff(2)-coeffs().coeff(0)*other.coeffs().coeff(2)));
213
/** Applies the transformation matrix \a mat to \c *this and returns a reference to \c *this.
215
* \param mat the Dim x Dim transformation matrix
216
* \param traits specifies whether the matrix \a mat represents an #Isometry
217
* or a more generic #Affine transformation. The default is #Affine.
219
template<typename XprType>
220
inline Hyperplane& transform(const MatrixBase<XprType>& mat, TransformTraits traits = Affine)
223
normal() = mat.inverse().transpose() * normal();
224
else if (traits==Isometry)
225
normal() = mat * normal();
228
eigen_assert("invalid traits value in Hyperplane::transform()");
233
/** Applies the transformation \a t to \c *this and returns a reference to \c *this.
235
* \param t the transformation of dimension Dim
236
* \param traits specifies whether the transformation \a t represents an #Isometry
237
* or a more generic #Affine transformation. The default is #Affine.
238
* Other kind of transformations are not supported.
240
template<int TrOptions>
241
inline Hyperplane& transform(const Transform<Scalar,AmbientDimAtCompileTime,Affine,TrOptions>& t,
242
TransformTraits traits = Affine)
244
transform(t.linear(), traits);
245
offset() -= normal().dot(t.translation());
249
/** \returns \c *this with scalar type casted to \a NewScalarType
251
* Note that if \a NewScalarType is equal to the current scalar type of \c *this
252
* then this function smartly returns a const reference to \c *this.
254
template<typename NewScalarType>
255
inline typename internal::cast_return_type<Hyperplane,
256
Hyperplane<NewScalarType,AmbientDimAtCompileTime,Options> >::type cast() const
258
return typename internal::cast_return_type<Hyperplane,
259
Hyperplane<NewScalarType,AmbientDimAtCompileTime,Options> >::type(*this);
262
/** Copy constructor with scalar type conversion */
263
template<typename OtherScalarType,int OtherOptions>
264
inline explicit Hyperplane(const Hyperplane<OtherScalarType,AmbientDimAtCompileTime,OtherOptions>& other)
265
{ m_coeffs = other.coeffs().template cast<Scalar>(); }
267
/** \returns \c true if \c *this is approximately equal to \a other, within the precision
268
* determined by \a prec.
270
* \sa MatrixBase::isApprox() */
271
template<int OtherOptions>
272
bool isApprox(const Hyperplane<Scalar,AmbientDimAtCompileTime,OtherOptions>& other, typename NumTraits<Scalar>::Real prec = NumTraits<Scalar>::dummy_precision()) const
273
{ return m_coeffs.isApprox(other.m_coeffs, prec); }
277
Coefficients m_coeffs;
280
#endif // EIGEN_HYPERPLANE_H