Home | History | Annotate | Download | only in Geometry
      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 (at) inria.fr>
      5 // Copyright (C) 2008 Benoit Jacob <jacob.benoit.1 (at) 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 
     14 namespace 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   * \param _Scalar the scalar type, i.e., the type of the coefficients
     26   * \param _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   */
     33 template <typename _Scalar, int _AmbientDim, int _Options>
     34 class Hyperplane
     35 {
     36 public:
     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 DenseIndex Index;
     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   inline explicit Hyperplane() {}
     54 
     55   template<int OtherOptions>
     56   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   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   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   inline Hyperplane(const VectorType& n, 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   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   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     result.normal() = (p2 - p0).cross(p1 - p0).normalized();
    104     result.offset() = -p0.dot(result.normal());
    105     return result;
    106   }
    107 
    108   /** Constructs a hyperplane passing through the parametrized line \a parametrized.
    109     * If the dimension of the ambient space is greater than 2, then there isn't uniqueness,
    110     * so an arbitrary choice is made.
    111     */
    112   // FIXME to be consitent with the rest this could be implemented as a static Through function ??
    113   explicit Hyperplane(const ParametrizedLine<Scalar, AmbientDimAtCompileTime>& parametrized)
    114   {
    115     normal() = parametrized.direction().unitOrthogonal();
    116     offset() = -parametrized.origin().dot(normal());
    117   }
    118 
    119   ~Hyperplane() {}
    120 
    121   /** \returns the dimension in which the plane holds */
    122   inline Index dim() const { return AmbientDimAtCompileTime==Dynamic ? m_coeffs.size()-1 : Index(AmbientDimAtCompileTime); }
    123 
    124   /** normalizes \c *this */
    125   void normalize(void)
    126   {
    127     m_coeffs /= normal().norm();
    128   }
    129 
    130   /** \returns the signed distance between the plane \c *this and a point \a p.
    131     * \sa absDistance()
    132     */
    133   inline Scalar signedDistance(const VectorType& p) const { return normal().dot(p) + offset(); }
    134 
    135   /** \returns the absolute distance between the plane \c *this and a point \a p.
    136     * \sa signedDistance()
    137     */
    138   inline Scalar absDistance(const VectorType& p) const { return internal::abs(signedDistance(p)); }
    139 
    140   /** \returns the projection of a point \a p onto the plane \c *this.
    141     */
    142   inline VectorType projection(const VectorType& p) const { return p - signedDistance(p) * normal(); }
    143 
    144   /** \returns a constant reference to the unit normal vector of the plane, which corresponds
    145     * to the linear part of the implicit equation.
    146     */
    147   inline ConstNormalReturnType normal() const { return ConstNormalReturnType(m_coeffs,0,0,dim(),1); }
    148 
    149   /** \returns a non-constant reference to the unit normal vector of the plane, which corresponds
    150     * to the linear part of the implicit equation.
    151     */
    152   inline NormalReturnType normal() { return NormalReturnType(m_coeffs,0,0,dim(),1); }
    153 
    154   /** \returns the distance to the origin, which is also the "constant term" of the implicit equation
    155     * \warning the vector normal is assumed to be normalized.
    156     */
    157   inline const Scalar& offset() const { return m_coeffs.coeff(dim()); }
    158 
    159   /** \returns a non-constant reference to the distance to the origin, which is also the constant part
    160     * of the implicit equation */
    161   inline Scalar& offset() { return m_coeffs(dim()); }
    162 
    163   /** \returns a constant reference to the coefficients c_i of the plane equation:
    164     * \f$ c_0*x_0 + ... + c_{d-1}*x_{d-1} + c_d = 0 \f$
    165     */
    166   inline const Coefficients& coeffs() const { return m_coeffs; }
    167 
    168   /** \returns a non-constant reference to the coefficients c_i of the plane equation:
    169     * \f$ c_0*x_0 + ... + c_{d-1}*x_{d-1} + c_d = 0 \f$
    170     */
    171   inline Coefficients& coeffs() { return m_coeffs; }
    172 
    173   /** \returns the intersection of *this with \a other.
    174     *
    175     * \warning The ambient space must be a plane, i.e. have dimension 2, so that \c *this and \a other are lines.
    176     *
    177     * \note If \a other is approximately parallel to *this, this method will return any point on *this.
    178     */
    179   VectorType intersection(const Hyperplane& other) const
    180   {
    181     EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 2)
    182     Scalar det = coeffs().coeff(0) * other.coeffs().coeff(1) - coeffs().coeff(1) * other.coeffs().coeff(0);
    183     // since the line equations ax+by=c are normalized with a^2+b^2=1, the following tests
    184     // whether the two lines are approximately parallel.
    185     if(internal::isMuchSmallerThan(det, Scalar(1)))
    186     {   // special case where the two lines are approximately parallel. Pick any point on the first line.
    187         if(internal::abs(coeffs().coeff(1))>internal::abs(coeffs().coeff(0)))
    188             return VectorType(coeffs().coeff(1), -coeffs().coeff(2)/coeffs().coeff(1)-coeffs().coeff(0));
    189         else
    190             return VectorType(-coeffs().coeff(2)/coeffs().coeff(0)-coeffs().coeff(1), coeffs().coeff(0));
    191     }
    192     else
    193     {   // general case
    194         Scalar invdet = Scalar(1) / det;
    195         return VectorType(invdet*(coeffs().coeff(1)*other.coeffs().coeff(2)-other.coeffs().coeff(1)*coeffs().coeff(2)),
    196                           invdet*(other.coeffs().coeff(0)*coeffs().coeff(2)-coeffs().coeff(0)*other.coeffs().coeff(2)));
    197     }
    198   }
    199 
    200   /** Applies the transformation matrix \a mat to \c *this and returns a reference to \c *this.
    201     *
    202     * \param mat the Dim x Dim transformation matrix
    203     * \param traits specifies whether the matrix \a mat represents an #Isometry
    204     *               or a more generic #Affine transformation. The default is #Affine.
    205     */
    206   template<typename XprType>
    207   inline Hyperplane& transform(const MatrixBase<XprType>& mat, TransformTraits traits = Affine)
    208   {
    209     if (traits==Affine)
    210       normal() = mat.inverse().transpose() * normal();
    211     else if (traits==Isometry)
    212       normal() = mat * normal();
    213     else
    214     {
    215       eigen_assert(0 && "invalid traits value in Hyperplane::transform()");
    216     }
    217     return *this;
    218   }
    219 
    220   /** Applies the transformation \a t to \c *this and returns a reference to \c *this.
    221     *
    222     * \param t the transformation of dimension Dim
    223     * \param traits specifies whether the transformation \a t represents an #Isometry
    224     *               or a more generic #Affine transformation. The default is #Affine.
    225     *               Other kind of transformations are not supported.
    226     */
    227   template<int TrOptions>
    228   inline Hyperplane& transform(const Transform<Scalar,AmbientDimAtCompileTime,Affine,TrOptions>& t,
    229                                 TransformTraits traits = Affine)
    230   {
    231     transform(t.linear(), traits);
    232     offset() -= normal().dot(t.translation());
    233     return *this;
    234   }
    235 
    236   /** \returns \c *this with scalar type casted to \a NewScalarType
    237     *
    238     * Note that if \a NewScalarType is equal to the current scalar type of \c *this
    239     * then this function smartly returns a const reference to \c *this.
    240     */
    241   template<typename NewScalarType>
    242   inline typename internal::cast_return_type<Hyperplane,
    243            Hyperplane<NewScalarType,AmbientDimAtCompileTime,Options> >::type cast() const
    244   {
    245     return typename internal::cast_return_type<Hyperplane,
    246                     Hyperplane<NewScalarType,AmbientDimAtCompileTime,Options> >::type(*this);
    247   }
    248 
    249   /** Copy constructor with scalar type conversion */
    250   template<typename OtherScalarType,int OtherOptions>
    251   inline explicit Hyperplane(const Hyperplane<OtherScalarType,AmbientDimAtCompileTime,OtherOptions>& other)
    252   { m_coeffs = other.coeffs().template cast<Scalar>(); }
    253 
    254   /** \returns \c true if \c *this is approximately equal to \a other, within the precision
    255     * determined by \a prec.
    256     *
    257     * \sa MatrixBase::isApprox() */
    258   template<int OtherOptions>
    259   bool isApprox(const Hyperplane<Scalar,AmbientDimAtCompileTime,OtherOptions>& other, typename NumTraits<Scalar>::Real prec = NumTraits<Scalar>::dummy_precision()) const
    260   { return m_coeffs.isApprox(other.m_coeffs, prec); }
    261 
    262 protected:
    263 
    264   Coefficients m_coeffs;
    265 };
    266 
    267 } // end namespace Eigen
    268 
    269 #endif // EIGEN_HYPERPLANE_H
    270