1 // This file is part of Eigen, a lightweight C++ template library 2 // for linear algebra. Eigen itself is part of the KDE project. 3 // 4 // Copyright (C) 2008 Gael Guennebaud <g.gael (at) free.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 // 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 ParametrizedLine 18 * 19 * \brief A parametrized line 20 * 21 * A parametrized line is defined by an origin point \f$ \mathbf{o} \f$ and a unit 22 * direction vector \f$ \mathbf{d} \f$ such that the line corresponds to 23 * the set \f$ l(t) = \mathbf{o} + t \mathbf{d} \f$, \f$ l \in \mathbf{R} \f$. 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 */ 28 template <typename _Scalar, int _AmbientDim> 29 class ParametrizedLine 30 { 31 public: 32 EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(_Scalar,_AmbientDim) 33 enum { AmbientDimAtCompileTime = _AmbientDim }; 34 typedef _Scalar Scalar; 35 typedef typename NumTraits<Scalar>::Real RealScalar; 36 typedef Matrix<Scalar,AmbientDimAtCompileTime,1> VectorType; 37 38 /** Default constructor without initialization */ 39 inline explicit ParametrizedLine() {} 40 41 /** Constructs a dynamic-size line with \a _dim the dimension 42 * of the ambient space */ 43 inline explicit ParametrizedLine(int _dim) : m_origin(_dim), m_direction(_dim) {} 44 45 /** Initializes a parametrized line of direction \a direction and origin \a origin. 46 * \warning the vector direction is assumed to be normalized. 47 */ 48 ParametrizedLine(const VectorType& origin, const VectorType& direction) 49 : m_origin(origin), m_direction(direction) {} 50 51 explicit ParametrizedLine(const Hyperplane<_Scalar, _AmbientDim>& hyperplane); 52 53 /** Constructs a parametrized line going from \a p0 to \a p1. */ 54 static inline ParametrizedLine Through(const VectorType& p0, const VectorType& p1) 55 { return ParametrizedLine(p0, (p1-p0).normalized()); } 56 57 ~ParametrizedLine() {} 58 59 /** \returns the dimension in which the line holds */ 60 inline int dim() const { return m_direction.size(); } 61 62 const VectorType& origin() const { return m_origin; } 63 VectorType& origin() { return m_origin; } 64 65 const VectorType& direction() const { return m_direction; } 66 VectorType& direction() { return m_direction; } 67 68 /** \returns the squared distance of a point \a p to its projection onto the line \c *this. 69 * \sa distance() 70 */ 71 RealScalar squaredDistance(const VectorType& p) const 72 { 73 VectorType diff = p-origin(); 74 return (diff - diff.eigen2_dot(direction())* direction()).squaredNorm(); 75 } 76 /** \returns the distance of a point \a p to its projection onto the line \c *this. 77 * \sa squaredDistance() 78 */ 79 RealScalar distance(const VectorType& p) const { return ei_sqrt(squaredDistance(p)); } 80 81 /** \returns the projection of a point \a p onto the line \c *this. */ 82 VectorType projection(const VectorType& p) const 83 { return origin() + (p-origin()).eigen2_dot(direction()) * direction(); } 84 85 Scalar intersection(const Hyperplane<_Scalar, _AmbientDim>& hyperplane); 86 87 /** \returns \c *this with scalar type casted to \a NewScalarType 88 * 89 * Note that if \a NewScalarType is equal to the current scalar type of \c *this 90 * then this function smartly returns a const reference to \c *this. 91 */ 92 template<typename NewScalarType> 93 inline typename internal::cast_return_type<ParametrizedLine, 94 ParametrizedLine<NewScalarType,AmbientDimAtCompileTime> >::type cast() const 95 { 96 return typename internal::cast_return_type<ParametrizedLine, 97 ParametrizedLine<NewScalarType,AmbientDimAtCompileTime> >::type(*this); 98 } 99 100 /** Copy constructor with scalar type conversion */ 101 template<typename OtherScalarType> 102 inline explicit ParametrizedLine(const ParametrizedLine<OtherScalarType,AmbientDimAtCompileTime>& other) 103 { 104 m_origin = other.origin().template cast<Scalar>(); 105 m_direction = other.direction().template cast<Scalar>(); 106 } 107 108 /** \returns \c true if \c *this is approximately equal to \a other, within the precision 109 * determined by \a prec. 110 * 111 * \sa MatrixBase::isApprox() */ 112 bool isApprox(const ParametrizedLine& other, typename NumTraits<Scalar>::Real prec = precision<Scalar>()) const 113 { return m_origin.isApprox(other.m_origin, prec) && m_direction.isApprox(other.m_direction, prec); } 114 115 protected: 116 117 VectorType m_origin, m_direction; 118 }; 119 120 /** Constructs a parametrized line from a 2D hyperplane 121 * 122 * \warning the ambient space must have dimension 2 such that the hyperplane actually describes a line 123 */ 124 template <typename _Scalar, int _AmbientDim> 125 inline ParametrizedLine<_Scalar, _AmbientDim>::ParametrizedLine(const Hyperplane<_Scalar, _AmbientDim>& hyperplane) 126 { 127 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 2) 128 direction() = hyperplane.normal().unitOrthogonal(); 129 origin() = -hyperplane.normal()*hyperplane.offset(); 130 } 131 132 /** \returns the parameter value of the intersection between \c *this and the given hyperplane 133 */ 134 template <typename _Scalar, int _AmbientDim> 135 inline _Scalar ParametrizedLine<_Scalar, _AmbientDim>::intersection(const Hyperplane<_Scalar, _AmbientDim>& hyperplane) 136 { 137 return -(hyperplane.offset()+origin().eigen2_dot(hyperplane.normal())) 138 /(direction().eigen2_dot(hyperplane.normal())); 139 } 140 141 } // end namespace Eigen 142