Hyperplane.h
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 
14 namespace Eigen {
15 
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 
53  inline explicit Hyperplane() {}
54 
55  template<int OtherOptions>
57  : m_coeffs(other.coeffs())
58  {}
59 
62  inline explicit Hyperplane(Index _dim) : m_coeffs(_dim+1) {}
63 
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 
78  inline Hyperplane(const VectorType& n, Scalar d)
79  : m_coeffs(n.size()+1)
80  {
81  normal() = n;
82  offset() = d;
83  }
84 
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 
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 
112  // FIXME to be consitent with the rest this could be implemented as a static Through function ??
114  {
115  normal() = parametrized.direction().unitOrthogonal();
116  offset() = -parametrized.origin().dot(normal());
117  }
118 
119  ~Hyperplane() {}
120 
122  inline Index dim() const { return AmbientDimAtCompileTime==Dynamic ? m_coeffs.size()-1 : Index(AmbientDimAtCompileTime); }
123 
125  void normalize(void)
126  {
127  m_coeffs /= normal().norm();
128  }
129 
133  inline Scalar signedDistance(const VectorType& p) const { return normal().dot(p) + offset(); }
134 
138  inline Scalar absDistance(const VectorType& p) const { return internal::abs(signedDistance(p)); }
139 
142  inline VectorType projection(const VectorType& p) const { return p - signedDistance(p) * normal(); }
143 
147  inline ConstNormalReturnType normal() const { return ConstNormalReturnType(m_coeffs,0,0,dim(),1); }
148 
152  inline NormalReturnType normal() { return NormalReturnType(m_coeffs,0,0,dim(),1); }
153 
157  inline const Scalar& offset() const { return m_coeffs.coeff(dim()); }
158 
161  inline Scalar& offset() { return m_coeffs(dim()); }
162 
166  inline const Coefficients& coeffs() const { return m_coeffs; }
167 
171  inline Coefficients& coeffs() { return m_coeffs; }
172 
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 
206  template<typename XprType>
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 
227  template<int TrOptions>
229  TransformTraits traits = Affine)
230  {
231  transform(t.linear(), traits);
232  offset() -= normal().dot(t.translation());
233  return *this;
234  }
235 
241  template<typename NewScalarType>
242  inline typename internal::cast_return_type<Hyperplane,
244  {
245  return typename internal::cast_return_type<Hyperplane,
247  }
248 
250  template<typename OtherScalarType,int OtherOptions>
252  { m_coeffs = other.coeffs().template cast<Scalar>(); }
253 
258  template<int OtherOptions>
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