ParametrizedLine.h 8.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195
  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. #ifndef EIGEN_PARAMETRIZEDLINE_H
  11. #define EIGEN_PARAMETRIZEDLINE_H
  12. namespace Eigen {
  13. /** \geometry_module \ingroup Geometry_Module
  14. *
  15. * \class ParametrizedLine
  16. *
  17. * \brief A parametrized line
  18. *
  19. * A parametrized line is defined by an origin point \f$ \mathbf{o} \f$ and a unit
  20. * direction vector \f$ \mathbf{d} \f$ such that the line corresponds to
  21. * the set \f$ l(t) = \mathbf{o} + t \mathbf{d} \f$, \f$ t \in \mathbf{R} \f$.
  22. *
  23. * \tparam _Scalar the scalar type, i.e., the type of the coefficients
  24. * \tparam _AmbientDim the dimension of the ambient space, can be a compile time value or Dynamic.
  25. */
  26. template <typename _Scalar, int _AmbientDim, int _Options>
  27. class ParametrizedLine
  28. {
  29. public:
  30. EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(_Scalar,_AmbientDim)
  31. enum {
  32. AmbientDimAtCompileTime = _AmbientDim,
  33. Options = _Options
  34. };
  35. typedef _Scalar Scalar;
  36. typedef typename NumTraits<Scalar>::Real RealScalar;
  37. typedef Eigen::Index Index; ///< \deprecated since Eigen 3.3
  38. typedef Matrix<Scalar,AmbientDimAtCompileTime,1,Options> VectorType;
  39. /** Default constructor without initialization */
  40. EIGEN_DEVICE_FUNC inline ParametrizedLine() {}
  41. template<int OtherOptions>
  42. EIGEN_DEVICE_FUNC ParametrizedLine(const ParametrizedLine<Scalar,AmbientDimAtCompileTime,OtherOptions>& other)
  43. : m_origin(other.origin()), m_direction(other.direction())
  44. {}
  45. /** Constructs a dynamic-size line with \a _dim the dimension
  46. * of the ambient space */
  47. EIGEN_DEVICE_FUNC inline explicit ParametrizedLine(Index _dim) : m_origin(_dim), m_direction(_dim) {}
  48. /** Initializes a parametrized line of direction \a direction and origin \a origin.
  49. * \warning the vector direction is assumed to be normalized.
  50. */
  51. EIGEN_DEVICE_FUNC ParametrizedLine(const VectorType& origin, const VectorType& direction)
  52. : m_origin(origin), m_direction(direction) {}
  53. template <int OtherOptions>
  54. EIGEN_DEVICE_FUNC explicit ParametrizedLine(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane);
  55. /** Constructs a parametrized line going from \a p0 to \a p1. */
  56. EIGEN_DEVICE_FUNC static inline ParametrizedLine Through(const VectorType& p0, const VectorType& p1)
  57. { return ParametrizedLine(p0, (p1-p0).normalized()); }
  58. EIGEN_DEVICE_FUNC ~ParametrizedLine() {}
  59. /** \returns the dimension in which the line holds */
  60. EIGEN_DEVICE_FUNC inline Index dim() const { return m_direction.size(); }
  61. EIGEN_DEVICE_FUNC const VectorType& origin() const { return m_origin; }
  62. EIGEN_DEVICE_FUNC VectorType& origin() { return m_origin; }
  63. EIGEN_DEVICE_FUNC const VectorType& direction() const { return m_direction; }
  64. EIGEN_DEVICE_FUNC VectorType& direction() { return m_direction; }
  65. /** \returns the squared distance of a point \a p to its projection onto the line \c *this.
  66. * \sa distance()
  67. */
  68. EIGEN_DEVICE_FUNC RealScalar squaredDistance(const VectorType& p) const
  69. {
  70. VectorType diff = p - origin();
  71. return (diff - direction().dot(diff) * direction()).squaredNorm();
  72. }
  73. /** \returns the distance of a point \a p to its projection onto the line \c *this.
  74. * \sa squaredDistance()
  75. */
  76. EIGEN_DEVICE_FUNC RealScalar distance(const VectorType& p) const { EIGEN_USING_STD_MATH(sqrt) return sqrt(squaredDistance(p)); }
  77. /** \returns the projection of a point \a p onto the line \c *this. */
  78. EIGEN_DEVICE_FUNC VectorType projection(const VectorType& p) const
  79. { return origin() + direction().dot(p-origin()) * direction(); }
  80. EIGEN_DEVICE_FUNC VectorType pointAt(const Scalar& t) const;
  81. template <int OtherOptions>
  82. EIGEN_DEVICE_FUNC Scalar intersectionParameter(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const;
  83. template <int OtherOptions>
  84. EIGEN_DEVICE_FUNC Scalar intersection(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const;
  85. template <int OtherOptions>
  86. EIGEN_DEVICE_FUNC VectorType intersectionPoint(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const;
  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. EIGEN_DEVICE_FUNC inline typename internal::cast_return_type<ParametrizedLine,
  94. ParametrizedLine<NewScalarType,AmbientDimAtCompileTime,Options> >::type cast() const
  95. {
  96. return typename internal::cast_return_type<ParametrizedLine,
  97. ParametrizedLine<NewScalarType,AmbientDimAtCompileTime,Options> >::type(*this);
  98. }
  99. /** Copy constructor with scalar type conversion */
  100. template<typename OtherScalarType,int OtherOptions>
  101. EIGEN_DEVICE_FUNC inline explicit ParametrizedLine(const ParametrizedLine<OtherScalarType,AmbientDimAtCompileTime,OtherOptions>& other)
  102. {
  103. m_origin = other.origin().template cast<Scalar>();
  104. m_direction = other.direction().template cast<Scalar>();
  105. }
  106. /** \returns \c true if \c *this is approximately equal to \a other, within the precision
  107. * determined by \a prec.
  108. *
  109. * \sa MatrixBase::isApprox() */
  110. EIGEN_DEVICE_FUNC bool isApprox(const ParametrizedLine& other, const typename NumTraits<Scalar>::Real& prec = NumTraits<Scalar>::dummy_precision()) const
  111. { return m_origin.isApprox(other.m_origin, prec) && m_direction.isApprox(other.m_direction, prec); }
  112. protected:
  113. VectorType m_origin, m_direction;
  114. };
  115. /** Constructs a parametrized line from a 2D hyperplane
  116. *
  117. * \warning the ambient space must have dimension 2 such that the hyperplane actually describes a line
  118. */
  119. template <typename _Scalar, int _AmbientDim, int _Options>
  120. template <int OtherOptions>
  121. EIGEN_DEVICE_FUNC inline ParametrizedLine<_Scalar, _AmbientDim,_Options>::ParametrizedLine(const Hyperplane<_Scalar, _AmbientDim,OtherOptions>& hyperplane)
  122. {
  123. EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 2)
  124. direction() = hyperplane.normal().unitOrthogonal();
  125. origin() = -hyperplane.normal()*hyperplane.offset();
  126. }
  127. /** \returns the point at \a t along this line
  128. */
  129. template <typename _Scalar, int _AmbientDim, int _Options>
  130. EIGEN_DEVICE_FUNC inline typename ParametrizedLine<_Scalar, _AmbientDim,_Options>::VectorType
  131. ParametrizedLine<_Scalar, _AmbientDim,_Options>::pointAt(const _Scalar& t) const
  132. {
  133. return origin() + (direction()*t);
  134. }
  135. /** \returns the parameter value of the intersection between \c *this and the given \a hyperplane
  136. */
  137. template <typename _Scalar, int _AmbientDim, int _Options>
  138. template <int OtherOptions>
  139. EIGEN_DEVICE_FUNC inline _Scalar ParametrizedLine<_Scalar, _AmbientDim,_Options>::intersectionParameter(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const
  140. {
  141. return -(hyperplane.offset()+hyperplane.normal().dot(origin()))
  142. / hyperplane.normal().dot(direction());
  143. }
  144. /** \deprecated use intersectionParameter()
  145. * \returns the parameter value of the intersection between \c *this and the given \a hyperplane
  146. */
  147. template <typename _Scalar, int _AmbientDim, int _Options>
  148. template <int OtherOptions>
  149. EIGEN_DEVICE_FUNC inline _Scalar ParametrizedLine<_Scalar, _AmbientDim,_Options>::intersection(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const
  150. {
  151. return intersectionParameter(hyperplane);
  152. }
  153. /** \returns the point of the intersection between \c *this and the given hyperplane
  154. */
  155. template <typename _Scalar, int _AmbientDim, int _Options>
  156. template <int OtherOptions>
  157. EIGEN_DEVICE_FUNC inline typename ParametrizedLine<_Scalar, _AmbientDim,_Options>::VectorType
  158. ParametrizedLine<_Scalar, _AmbientDim,_Options>::intersectionPoint(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const
  159. {
  160. return pointAt(intersectionParameter(hyperplane));
  161. }
  162. } // end namespace Eigen
  163. #endif // EIGEN_PARAMETRIZEDLINE_H