Hyperplane.h 12 KB

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