Hyperplane.h revision 2b8756b6f1de65d3f8bffab45be6c44ceb7411fc
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.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#ifndef EIGEN_HYPERPLANE_H
12#define EIGEN_HYPERPLANE_H
13
14namespace Eigen {
15
16/** \geometry_module \ingroup Geometry_Module
17  *
18  * \class Hyperplane
19  *
20  * \brief A hyperplane
21  *
22  * A hyperplane is an affine subspace of dimension n-1 in a space of dimension n.
23  * For example, a hyperplane in a plane is a line; a hyperplane in 3-space is a plane.
24  *
25  * \tparam _Scalar the scalar type, i.e., the type of the coefficients
26  * \tparam _AmbientDim the dimension of the ambient space, can be a compile time value or Dynamic.
27  *             Notice that the dimension of the hyperplane is _AmbientDim-1.
28  *
29  * This class represents an hyperplane as the zero set of the implicit equation
30  * \f$ n \cdot x + d = 0 \f$ where \f$ n \f$ is a unit normal vector of the plane (linear part)
31  * and \f$ d \f$ is the distance (offset) to the origin.
32  */
33template <typename _Scalar, int _AmbientDim, int _Options>
34class Hyperplane
35{
36public:
37  EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(_Scalar,_AmbientDim==Dynamic ? Dynamic : _AmbientDim+1)
38  enum {
39    AmbientDimAtCompileTime = _AmbientDim,
40    Options = _Options
41  };
42  typedef _Scalar Scalar;
43  typedef typename NumTraits<Scalar>::Real RealScalar;
44  typedef Eigen::Index Index; ///< \deprecated since Eigen 3.3
45  typedef Matrix<Scalar,AmbientDimAtCompileTime,1> VectorType;
46  typedef Matrix<Scalar,Index(AmbientDimAtCompileTime)==Dynamic
47                        ? Dynamic
48                        : Index(AmbientDimAtCompileTime)+1,1,Options> Coefficients;
49  typedef Block<Coefficients,AmbientDimAtCompileTime,1> NormalReturnType;
50  typedef const Block<const Coefficients,AmbientDimAtCompileTime,1> ConstNormalReturnType;
51
52  /** Default constructor without initialization */
53  EIGEN_DEVICE_FUNC inline Hyperplane() {}
54
55  template<int OtherOptions>
56  EIGEN_DEVICE_FUNC Hyperplane(const Hyperplane<Scalar,AmbientDimAtCompileTime,OtherOptions>& other)
57   : m_coeffs(other.coeffs())
58  {}
59
60  /** Constructs a dynamic-size hyperplane with \a _dim the dimension
61    * of the ambient space */
62  EIGEN_DEVICE_FUNC inline explicit Hyperplane(Index _dim) : m_coeffs(_dim+1) {}
63
64  /** Construct a plane from its normal \a n and a point \a e onto the plane.
65    * \warning the vector normal is assumed to be normalized.
66    */
67  EIGEN_DEVICE_FUNC inline Hyperplane(const VectorType& n, const VectorType& e)
68    : m_coeffs(n.size()+1)
69  {
70    normal() = n;
71    offset() = -n.dot(e);
72  }
73
74  /** Constructs a plane from its normal \a n and distance to the origin \a d
75    * such that the algebraic equation of the plane is \f$ n \cdot x + d = 0 \f$.
76    * \warning the vector normal is assumed to be normalized.
77    */
78  EIGEN_DEVICE_FUNC inline Hyperplane(const VectorType& n, const Scalar& d)
79    : m_coeffs(n.size()+1)
80  {
81    normal() = n;
82    offset() = d;
83  }
84
85  /** Constructs a hyperplane passing through the two points. If the dimension of the ambient space
86    * is greater than 2, then there isn't uniqueness, so an arbitrary choice is made.
87    */
88  EIGEN_DEVICE_FUNC static inline Hyperplane Through(const VectorType& p0, const VectorType& p1)
89  {
90    Hyperplane result(p0.size());
91    result.normal() = (p1 - p0).unitOrthogonal();
92    result.offset() = -p0.dot(result.normal());
93    return result;
94  }
95
96  /** Constructs a hyperplane passing through the three points. The dimension of the ambient space
97    * is required to be exactly 3.
98    */
99  EIGEN_DEVICE_FUNC static inline Hyperplane Through(const VectorType& p0, const VectorType& p1, const VectorType& p2)
100  {
101    EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 3)
102    Hyperplane result(p0.size());
103    VectorType v0(p2 - p0), v1(p1 - p0);
104    result.normal() = v0.cross(v1);
105    RealScalar norm = result.normal().norm();
106    if(norm <= v0.norm() * v1.norm() * NumTraits<RealScalar>::epsilon())
107    {
108      Matrix<Scalar,2,3> m; m << v0.transpose(), v1.transpose();
109      JacobiSVD<Matrix<Scalar,2,3> > svd(m, ComputeFullV);
110      result.normal() = svd.matrixV().col(2);
111    }
112    else
113      result.normal() /= norm;
114    result.offset() = -p0.dot(result.normal());
115    return result;
116  }
117
118  /** Constructs a hyperplane passing through the parametrized line \a parametrized.
119    * If the dimension of the ambient space is greater than 2, then there isn't uniqueness,
120    * so an arbitrary choice is made.
121    */
122  // FIXME to be consitent with the rest this could be implemented as a static Through function ??
123  EIGEN_DEVICE_FUNC explicit Hyperplane(const ParametrizedLine<Scalar, AmbientDimAtCompileTime>& parametrized)
124  {
125    normal() = parametrized.direction().unitOrthogonal();
126    offset() = -parametrized.origin().dot(normal());
127  }
128
129  EIGEN_DEVICE_FUNC ~Hyperplane() {}
130
131  /** \returns the dimension in which the plane holds */
132  EIGEN_DEVICE_FUNC inline Index dim() const { return AmbientDimAtCompileTime==Dynamic ? m_coeffs.size()-1 : Index(AmbientDimAtCompileTime); }
133
134  /** normalizes \c *this */
135  EIGEN_DEVICE_FUNC void normalize(void)
136  {
137    m_coeffs /= normal().norm();
138  }
139
140  /** \returns the signed distance between the plane \c *this and a point \a p.
141    * \sa absDistance()
142    */
143  EIGEN_DEVICE_FUNC inline Scalar signedDistance(const VectorType& p) const { return normal().dot(p) + offset(); }
144
145  /** \returns the absolute distance between the plane \c *this and a point \a p.
146    * \sa signedDistance()
147    */
148  EIGEN_DEVICE_FUNC inline Scalar absDistance(const VectorType& p) const { return numext::abs(signedDistance(p)); }
149
150  /** \returns the projection of a point \a p onto the plane \c *this.
151    */
152  EIGEN_DEVICE_FUNC inline VectorType projection(const VectorType& p) const { return p - signedDistance(p) * normal(); }
153
154  /** \returns a constant reference to the unit normal vector of the plane, which corresponds
155    * to the linear part of the implicit equation.
156    */
157  EIGEN_DEVICE_FUNC inline ConstNormalReturnType normal() const { return ConstNormalReturnType(m_coeffs,0,0,dim(),1); }
158
159  /** \returns a non-constant reference to the unit normal vector of the plane, which corresponds
160    * to the linear part of the implicit equation.
161    */
162  EIGEN_DEVICE_FUNC inline NormalReturnType normal() { return NormalReturnType(m_coeffs,0,0,dim(),1); }
163
164  /** \returns the distance to the origin, which is also the "constant term" of the implicit equation
165    * \warning the vector normal is assumed to be normalized.
166    */
167  EIGEN_DEVICE_FUNC inline const Scalar& offset() const { return m_coeffs.coeff(dim()); }
168
169  /** \returns a non-constant reference to the distance to the origin, which is also the constant part
170    * of the implicit equation */
171  EIGEN_DEVICE_FUNC inline Scalar& offset() { return m_coeffs(dim()); }
172
173  /** \returns a constant reference to the coefficients c_i of the plane equation:
174    * \f$ c_0*x_0 + ... + c_{d-1}*x_{d-1} + c_d = 0 \f$
175    */
176  EIGEN_DEVICE_FUNC inline const Coefficients& coeffs() const { return m_coeffs; }
177
178  /** \returns a non-constant reference to the coefficients c_i of the plane equation:
179    * \f$ c_0*x_0 + ... + c_{d-1}*x_{d-1} + c_d = 0 \f$
180    */
181  EIGEN_DEVICE_FUNC inline Coefficients& coeffs() { return m_coeffs; }
182
183  /** \returns the intersection of *this with \a other.
184    *
185    * \warning The ambient space must be a plane, i.e. have dimension 2, so that \c *this and \a other are lines.
186    *
187    * \note If \a other is approximately parallel to *this, this method will return any point on *this.
188    */
189  EIGEN_DEVICE_FUNC VectorType intersection(const Hyperplane& other) const
190  {
191    EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 2)
192    Scalar det = coeffs().coeff(0) * other.coeffs().coeff(1) - coeffs().coeff(1) * other.coeffs().coeff(0);
193    // since the line equations ax+by=c are normalized with a^2+b^2=1, the following tests
194    // whether the two lines are approximately parallel.
195    if(internal::isMuchSmallerThan(det, Scalar(1)))
196    {   // special case where the two lines are approximately parallel. Pick any point on the first line.
197        if(numext::abs(coeffs().coeff(1))>numext::abs(coeffs().coeff(0)))
198            return VectorType(coeffs().coeff(1), -coeffs().coeff(2)/coeffs().coeff(1)-coeffs().coeff(0));
199        else
200            return VectorType(-coeffs().coeff(2)/coeffs().coeff(0)-coeffs().coeff(1), coeffs().coeff(0));
201    }
202    else
203    {   // general case
204        Scalar invdet = Scalar(1) / det;
205        return VectorType(invdet*(coeffs().coeff(1)*other.coeffs().coeff(2)-other.coeffs().coeff(1)*coeffs().coeff(2)),
206                          invdet*(other.coeffs().coeff(0)*coeffs().coeff(2)-coeffs().coeff(0)*other.coeffs().coeff(2)));
207    }
208  }
209
210  /** Applies the transformation matrix \a mat to \c *this and returns a reference to \c *this.
211    *
212    * \param mat the Dim x Dim transformation matrix
213    * \param traits specifies whether the matrix \a mat represents an #Isometry
214    *               or a more generic #Affine transformation. The default is #Affine.
215    */
216  template<typename XprType>
217  EIGEN_DEVICE_FUNC inline Hyperplane& transform(const MatrixBase<XprType>& mat, TransformTraits traits = Affine)
218  {
219    if (traits==Affine)
220    {
221      normal() = mat.inverse().transpose() * normal();
222      m_coeffs /= normal().norm();
223    }
224    else if (traits==Isometry)
225      normal() = mat * normal();
226    else
227    {
228      eigen_assert(0 && "invalid traits value in Hyperplane::transform()");
229    }
230    return *this;
231  }
232
233  /** Applies the transformation \a t to \c *this and returns a reference to \c *this.
234    *
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.
239    */
240  template<int TrOptions>
241  EIGEN_DEVICE_FUNC inline Hyperplane& transform(const Transform<Scalar,AmbientDimAtCompileTime,Affine,TrOptions>& t,
242                                TransformTraits traits = Affine)
243  {
244    transform(t.linear(), traits);
245    offset() -= normal().dot(t.translation());
246    return *this;
247  }
248
249  /** \returns \c *this with scalar type casted to \a NewScalarType
250    *
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.
253    */
254  template<typename NewScalarType>
255  EIGEN_DEVICE_FUNC inline typename internal::cast_return_type<Hyperplane,
256           Hyperplane<NewScalarType,AmbientDimAtCompileTime,Options> >::type cast() const
257  {
258    return typename internal::cast_return_type<Hyperplane,
259                    Hyperplane<NewScalarType,AmbientDimAtCompileTime,Options> >::type(*this);
260  }
261
262  /** Copy constructor with scalar type conversion */
263  template<typename OtherScalarType,int OtherOptions>
264  EIGEN_DEVICE_FUNC inline explicit Hyperplane(const Hyperplane<OtherScalarType,AmbientDimAtCompileTime,OtherOptions>& other)
265  { m_coeffs = other.coeffs().template cast<Scalar>(); }
266
267  /** \returns \c true if \c *this is approximately equal to \a other, within the precision
268    * determined by \a prec.
269    *
270    * \sa MatrixBase::isApprox() */
271  template<int OtherOptions>
272  EIGEN_DEVICE_FUNC bool isApprox(const Hyperplane<Scalar,AmbientDimAtCompileTime,OtherOptions>& other, const typename NumTraits<Scalar>::Real& prec = NumTraits<Scalar>::dummy_precision()) const
273  { return m_coeffs.isApprox(other.m_coeffs, prec); }
274
275protected:
276
277  Coefficients m_coeffs;
278};
279
280} // end namespace Eigen
281
282#endif // EIGEN_HYPERPLANE_H
283