1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008 Gael Guennebaud <g.gael@free.fr>
5 // Copyright (C) 2008 Benoit Jacob <jacob.benoit.1@gmail.com>
6 //
7 // This Source Code Form is subject to the terms of the Mozilla
8 // Public License v. 2.0. If a copy of the MPL was not distributed
9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10 
11 // no include guard, we'll include this twice from All.h from Eigen2Support, and it's internal anyway
12 
13 namespace Eigen {
14 
15 /** \geometry_module \ingroup Geometry_Module
16   *
17   * \class Hyperplane
18   *
19   * \brief A hyperplane
20   *
21   * A hyperplane is an affine subspace of dimension n-1 in a space of dimension n.
22   * For example, a hyperplane in a plane is a line; a hyperplane in 3-space is a plane.
23   *
24   * \param _Scalar the scalar type, i.e., the type of the coefficients
25   * \param _AmbientDim the dimension of the ambient space, can be a compile time value or Dynamic.
26   *             Notice that the dimension of the hyperplane is _AmbientDim-1.
27   *
28   * This class represents an hyperplane as the zero set of the implicit equation
29   * \f$ n \cdot x + d = 0 \f$ where \f$ n \f$ is a unit normal vector of the plane (linear part)
30   * and \f$ d \f$ is the distance (offset) to the origin.
31   */
32 template <typename _Scalar, int _AmbientDim>
33 class Hyperplane
34 {
35 public:
36   EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(_Scalar,_AmbientDim==Dynamic ? Dynamic : _AmbientDim+1)
37   enum { AmbientDimAtCompileTime = _AmbientDim };
38   typedef _Scalar Scalar;
39   typedef typename NumTraits<Scalar>::Real RealScalar;
40   typedef Matrix<Scalar,AmbientDimAtCompileTime,1> VectorType;
41   typedef Matrix<Scalar,int(AmbientDimAtCompileTime)==Dynamic
42                         ? Dynamic
43                         : int(AmbientDimAtCompileTime)+1,1> Coefficients;
44   typedef Block<Coefficients,AmbientDimAtCompileTime,1> NormalReturnType;
45 
46   /** Default constructor without initialization */
Hyperplane()47   inline Hyperplane() {}
48 
49   /** Constructs a dynamic-size hyperplane with \a _dim the dimension
50     * of the ambient space */
Hyperplane(int _dim)51   inline explicit Hyperplane(int _dim) : m_coeffs(_dim+1) {}
52 
53   /** Construct a plane from its normal \a n and a point \a e onto the plane.
54     * \warning the vector normal is assumed to be normalized.
55     */
Hyperplane(const VectorType & n,const VectorType & e)56   inline Hyperplane(const VectorType& n, const VectorType& e)
57     : m_coeffs(n.size()+1)
58   {
59     normal() = n;
60     offset() = -e.eigen2_dot(n);
61   }
62 
63   /** Constructs a plane from its normal \a n and distance to the origin \a d
64     * such that the algebraic equation of the plane is \f$ n \cdot x + d = 0 \f$.
65     * \warning the vector normal is assumed to be normalized.
66     */
Hyperplane(const VectorType & n,Scalar d)67   inline Hyperplane(const VectorType& n, Scalar d)
68     : m_coeffs(n.size()+1)
69   {
70     normal() = n;
71     offset() = d;
72   }
73 
74   /** Constructs a hyperplane passing through the two points. If the dimension of the ambient space
75     * is greater than 2, then there isn't uniqueness, so an arbitrary choice is made.
76     */
Through(const VectorType & p0,const VectorType & p1)77   static inline Hyperplane Through(const VectorType& p0, const VectorType& p1)
78   {
79     Hyperplane result(p0.size());
80     result.normal() = (p1 - p0).unitOrthogonal();
81     result.offset() = -result.normal().eigen2_dot(p0);
82     return result;
83   }
84 
85   /** Constructs a hyperplane passing through the three points. The dimension of the ambient space
86     * is required to be exactly 3.
87     */
Through(const VectorType & p0,const VectorType & p1,const VectorType & p2)88   static inline Hyperplane Through(const VectorType& p0, const VectorType& p1, const VectorType& p2)
89   {
90     EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 3)
91     Hyperplane result(p0.size());
92     result.normal() = (p2 - p0).cross(p1 - p0).normalized();
93     result.offset() = -result.normal().eigen2_dot(p0);
94     return result;
95   }
96 
97   /** Constructs a hyperplane passing through the parametrized line \a parametrized.
98     * If the dimension of the ambient space is greater than 2, then there isn't uniqueness,
99     * so an arbitrary choice is made.
100     */
101   // FIXME to be consitent with the rest this could be implemented as a static Through function ??
Hyperplane(const ParametrizedLine<Scalar,AmbientDimAtCompileTime> & parametrized)102   explicit Hyperplane(const ParametrizedLine<Scalar, AmbientDimAtCompileTime>& parametrized)
103   {
104     normal() = parametrized.direction().unitOrthogonal();
105     offset() = -normal().eigen2_dot(parametrized.origin());
106   }
107 
~Hyperplane()108   ~Hyperplane() {}
109 
110   /** \returns the dimension in which the plane holds */
dim()111   inline int dim() const { return int(AmbientDimAtCompileTime)==Dynamic ? m_coeffs.size()-1 : int(AmbientDimAtCompileTime); }
112 
113   /** normalizes \c *this */
normalize(void)114   void normalize(void)
115   {
116     m_coeffs /= normal().norm();
117   }
118 
119   /** \returns the signed distance between the plane \c *this and a point \a p.
120     * \sa absDistance()
121     */
signedDistance(const VectorType & p)122   inline Scalar signedDistance(const VectorType& p) const { return p.eigen2_dot(normal()) + offset(); }
123 
124   /** \returns the absolute distance between the plane \c *this and a point \a p.
125     * \sa signedDistance()
126     */
absDistance(const VectorType & p)127   inline Scalar absDistance(const VectorType& p) const { return ei_abs(signedDistance(p)); }
128 
129   /** \returns the projection of a point \a p onto the plane \c *this.
130     */
projection(const VectorType & p)131   inline VectorType projection(const VectorType& p) const { return p - signedDistance(p) * normal(); }
132 
133   /** \returns a constant reference to the unit normal vector of the plane, which corresponds
134     * to the linear part of the implicit equation.
135     */
normal()136   inline const NormalReturnType normal() const { return NormalReturnType(*const_cast<Coefficients*>(&m_coeffs),0,0,dim(),1); }
137 
138   /** \returns a non-constant reference to the unit normal vector of the plane, which corresponds
139     * to the linear part of the implicit equation.
140     */
normal()141   inline NormalReturnType normal() { return NormalReturnType(m_coeffs,0,0,dim(),1); }
142 
143   /** \returns the distance to the origin, which is also the "constant term" of the implicit equation
144     * \warning the vector normal is assumed to be normalized.
145     */
offset()146   inline const Scalar& offset() const { return m_coeffs.coeff(dim()); }
147 
148   /** \returns a non-constant reference to the distance to the origin, which is also the constant part
149     * of the implicit equation */
offset()150   inline Scalar& offset() { return m_coeffs(dim()); }
151 
152   /** \returns a constant reference to the coefficients c_i of the plane equation:
153     * \f$ c_0*x_0 + ... + c_{d-1}*x_{d-1} + c_d = 0 \f$
154     */
coeffs()155   inline const Coefficients& coeffs() const { return m_coeffs; }
156 
157   /** \returns a non-constant reference to the coefficients c_i of the plane equation:
158     * \f$ c_0*x_0 + ... + c_{d-1}*x_{d-1} + c_d = 0 \f$
159     */
coeffs()160   inline Coefficients& coeffs() { return m_coeffs; }
161 
162   /** \returns the intersection of *this with \a other.
163     *
164     * \warning The ambient space must be a plane, i.e. have dimension 2, so that \c *this and \a other are lines.
165     *
166     * \note If \a other is approximately parallel to *this, this method will return any point on *this.
167     */
intersection(const Hyperplane & other)168   VectorType intersection(const Hyperplane& other)
169   {
170     EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 2)
171     Scalar det = coeffs().coeff(0) * other.coeffs().coeff(1) - coeffs().coeff(1) * other.coeffs().coeff(0);
172     // since the line equations ax+by=c are normalized with a^2+b^2=1, the following tests
173     // whether the two lines are approximately parallel.
174     if(ei_isMuchSmallerThan(det, Scalar(1)))
175     {   // special case where the two lines are approximately parallel. Pick any point on the first line.
176         if(ei_abs(coeffs().coeff(1))>ei_abs(coeffs().coeff(0)))
177             return VectorType(coeffs().coeff(1), -coeffs().coeff(2)/coeffs().coeff(1)-coeffs().coeff(0));
178         else
179             return VectorType(-coeffs().coeff(2)/coeffs().coeff(0)-coeffs().coeff(1), coeffs().coeff(0));
180     }
181     else
182     {   // general case
183         Scalar invdet = Scalar(1) / det;
184         return VectorType(invdet*(coeffs().coeff(1)*other.coeffs().coeff(2)-other.coeffs().coeff(1)*coeffs().coeff(2)),
185                           invdet*(other.coeffs().coeff(0)*coeffs().coeff(2)-coeffs().coeff(0)*other.coeffs().coeff(2)));
186     }
187   }
188 
189   /** Applies the transformation matrix \a mat to \c *this and returns a reference to \c *this.
190     *
191     * \param mat the Dim x Dim transformation matrix
192     * \param traits specifies whether the matrix \a mat represents an Isometry
193     *               or a more generic Affine transformation. The default is Affine.
194     */
195   template<typename XprType>
196   inline Hyperplane& transform(const MatrixBase<XprType>& mat, TransformTraits traits = Affine)
197   {
198     if (traits==Affine)
199       normal() = mat.inverse().transpose() * normal();
200     else if (traits==Isometry)
201       normal() = mat * normal();
202     else
203     {
204       ei_assert("invalid traits value in Hyperplane::transform()");
205     }
206     return *this;
207   }
208 
209   /** Applies the transformation \a t to \c *this and returns a reference to \c *this.
210     *
211     * \param t the transformation of dimension Dim
212     * \param traits specifies whether the transformation \a t represents an Isometry
213     *               or a more generic Affine transformation. The default is Affine.
214     *               Other kind of transformations are not supported.
215     */
216   inline Hyperplane& transform(const Transform<Scalar,AmbientDimAtCompileTime>& t,
217                                 TransformTraits traits = Affine)
218   {
219     transform(t.linear(), traits);
220     offset() -= t.translation().eigen2_dot(normal());
221     return *this;
222   }
223 
224   /** \returns \c *this with scalar type casted to \a NewScalarType
225     *
226     * Note that if \a NewScalarType is equal to the current scalar type of \c *this
227     * then this function smartly returns a const reference to \c *this.
228     */
229   template<typename NewScalarType>
230   inline typename internal::cast_return_type<Hyperplane,
cast()231            Hyperplane<NewScalarType,AmbientDimAtCompileTime> >::type cast() const
232   {
233     return typename internal::cast_return_type<Hyperplane,
234                     Hyperplane<NewScalarType,AmbientDimAtCompileTime> >::type(*this);
235   }
236 
237   /** Copy constructor with scalar type conversion */
238   template<typename OtherScalarType>
Hyperplane(const Hyperplane<OtherScalarType,AmbientDimAtCompileTime> & other)239   inline explicit Hyperplane(const Hyperplane<OtherScalarType,AmbientDimAtCompileTime>& other)
240   { m_coeffs = other.coeffs().template cast<Scalar>(); }
241 
242   /** \returns \c true if \c *this is approximately equal to \a other, within the precision
243     * determined by \a prec.
244     *
245     * \sa MatrixBase::isApprox() */
246   bool isApprox(const Hyperplane& other, typename NumTraits<Scalar>::Real prec = precision<Scalar>()) const
247   { return m_coeffs.isApprox(other.m_coeffs, prec); }
248 
249 protected:
250 
251   Coefficients m_coeffs;
252 };
253 
254 } // end namespace Eigen
255