Eigen  3.2.7
 All Classes Namespaces Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
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 <[email protected]>
5 // Copyright (C) 2008 Benoit Jacob <[email protected]>
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 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, const 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  VectorType v0(p2 - p0), v1(p1 - p0);
104  result.normal() = v0.cross(v1);
105  RealScalar norm = result.normal().norm();
106  if(norm <= v0.norm() * v1.norm() * NumTraits<RealScalar>::epsilon())
107  {
108  Matrix<Scalar,2,3> m; m << v0.transpose(), v1.transpose();
110  result.normal() = svd.matrixV().col(2);
111  }
112  else
113  result.normal() /= norm;
114  result.offset() = -p0.dot(result.normal());
115  return result;
116  }
117 
122  // FIXME to be consitent with the rest this could be implemented as a static Through function ??
124  {
125  normal() = parametrized.direction().unitOrthogonal();
126  offset() = -parametrized.origin().dot(normal());
127  }
128 
129  ~Hyperplane() {}
130 
132  inline Index dim() const { return AmbientDimAtCompileTime==Dynamic ? m_coeffs.size()-1 : Index(AmbientDimAtCompileTime); }
133 
135  void normalize(void)
136  {
137  m_coeffs /= normal().norm();
138  }
139 
143  inline Scalar signedDistance(const VectorType& p) const { return normal().dot(p) + offset(); }
144 
148  inline Scalar absDistance(const VectorType& p) const { using std::abs; return abs(signedDistance(p)); }
149 
152  inline VectorType projection(const VectorType& p) const { return p - signedDistance(p) * normal(); }
153 
157  inline ConstNormalReturnType normal() const { return ConstNormalReturnType(m_coeffs,0,0,dim(),1); }
158 
162  inline NormalReturnType normal() { return NormalReturnType(m_coeffs,0,0,dim(),1); }
163 
167  inline const Scalar& offset() const { return m_coeffs.coeff(dim()); }
168 
171  inline Scalar& offset() { return m_coeffs(dim()); }
172 
176  inline const Coefficients& coeffs() const { return m_coeffs; }
177 
181  inline Coefficients& coeffs() { return m_coeffs; }
182 
189  VectorType intersection(const Hyperplane& other) const
190  {
191  using std::abs;
192  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 2)
193  Scalar det = coeffs().coeff(0) * other.coeffs().coeff(1) - coeffs().coeff(1) * other.coeffs().coeff(0);
194  // since the line equations ax+by=c are normalized with a^2+b^2=1, the following tests
195  // whether the two lines are approximately parallel.
196  if(internal::isMuchSmallerThan(det, Scalar(1)))
197  { // special case where the two lines are approximately parallel. Pick any point on the first line.
198  if(abs(coeffs().coeff(1))>abs(coeffs().coeff(0)))
199  return VectorType(coeffs().coeff(1), -coeffs().coeff(2)/coeffs().coeff(1)-coeffs().coeff(0));
200  else
201  return VectorType(-coeffs().coeff(2)/coeffs().coeff(0)-coeffs().coeff(1), coeffs().coeff(0));
202  }
203  else
204  { // general case
205  Scalar invdet = Scalar(1) / det;
206  return VectorType(invdet*(coeffs().coeff(1)*other.coeffs().coeff(2)-other.coeffs().coeff(1)*coeffs().coeff(2)),
207  invdet*(other.coeffs().coeff(0)*coeffs().coeff(2)-coeffs().coeff(0)*other.coeffs().coeff(2)));
208  }
209  }
210 
217  template<typename XprType>
219  {
220  if (traits==Affine)
221  normal() = mat.inverse().transpose() * normal();
222  else if (traits==Isometry)
223  normal() = mat * normal();
224  else
225  {
226  eigen_assert(0 && "invalid traits value in Hyperplane::transform()");
227  }
228  return *this;
229  }
230 
238  template<int TrOptions>
240  TransformTraits traits = Affine)
241  {
242  transform(t.linear(), traits);
243  offset() -= normal().dot(t.translation());
244  return *this;
245  }
246 
252  template<typename NewScalarType>
253  inline typename internal::cast_return_type<Hyperplane,
255  {
256  return typename internal::cast_return_type<Hyperplane,
258  }
259 
261  template<typename OtherScalarType,int OtherOptions>
263  { m_coeffs = other.coeffs().template cast<Scalar>(); }
264 
269  template<int OtherOptions>
271  { return m_coeffs.isApprox(other.m_coeffs, prec); }
272 
273 protected:
274 
275  Coefficients m_coeffs;
276 };
277 
278 } // end namespace Eigen
279 
280 #endif // EIGEN_HYPERPLANE_H
ConstLinearPart linear() const
Definition: Transform.h:379
NormalReturnType normal()
Definition: Hyperplane.h:162
Hyperplane & transform(const Transform< Scalar, AmbientDimAtCompileTime, Affine, TrOptions > &t, TransformTraits traits=Affine)
Definition: Hyperplane.h:239
Hyperplane(const VectorType &n, const VectorType &e)
Definition: Hyperplane.h:67
Hyperplane(const ParametrizedLine< Scalar, AmbientDimAtCompileTime > &parametrized)
Definition: Hyperplane.h:123
Hyperplane()
Definition: Hyperplane.h:53
internal::cast_return_type< Hyperplane, Hyperplane< NewScalarType, AmbientDimAtCompileTime, Options > >::type cast() const
Definition: Hyperplane.h:254
const Coefficients & coeffs() const
Definition: Hyperplane.h:176
static Hyperplane Through(const VectorType &p0, const VectorType &p1, const VectorType &p2)
Definition: Hyperplane.h:99
Definition: Constants.h:394
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:88
const int Dynamic
Definition: Constants.h:21
Definition: Constants.h:331
Index dim() const
Definition: Hyperplane.h:132
Definition: Constants.h:391
VectorType intersection(const Hyperplane &other) const
Definition: Hyperplane.h:189
Hyperplane(const Hyperplane< OtherScalarType, AmbientDimAtCompileTime, OtherOptions > &other)
Definition: Hyperplane.h:262
Hyperplane(Index _dim)
Definition: Hyperplane.h:62
static Hyperplane Through(const VectorType &p0, const VectorType &p1)
Definition: Hyperplane.h:88
TransformTraits
Definition: Constants.h:389
const internal::inverse_impl< Derived > inverse() const
Definition: Inverse.h:320
Scalar & offset()
Definition: Hyperplane.h:171
A hyperplane.
Definition: ForwardDeclarations.h:266
Coefficients & coeffs()
Definition: Hyperplane.h:181
Scalar signedDistance(const VectorType &p) const
Definition: Hyperplane.h:143
Hyperplane(const VectorType &n, const Scalar &d)
Definition: Hyperplane.h:78
Hyperplane & transform(const MatrixBase< XprType > &mat, TransformTraits traits=Affine)
Definition: Hyperplane.h:218
const MatrixVType & matrixV() const
Definition: JacobiSVD.h:629
Expression of a fixed-size or dynamic-size block.
Definition: Block.h:103
Eigen::Transpose< Derived > transpose()
Definition: Transpose.h:199
Two-sided Jacobi SVD decomposition of a rectangular matrix.
Definition: ForwardDeclarations.h:224
VectorType projection(const VectorType &p) const
Definition: Hyperplane.h:152
ConstNormalReturnType normal() const
Definition: Hyperplane.h:157
const Scalar & offset() const
Definition: Hyperplane.h:167
void normalize(void)
Definition: Hyperplane.h:135
Scalar absDistance(const VectorType &p) const
Definition: Hyperplane.h:148
ConstTranslationPart translation() const
Definition: Transform.h:389
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
A parametrized line.
Definition: ForwardDeclarations.h:265
Represents an homogeneous transformation in a N dimensional space.
Definition: ForwardDeclarations.h:264
bool isApprox(const Hyperplane< Scalar, AmbientDimAtCompileTime, OtherOptions > &other, const typename NumTraits< Scalar >::Real &prec=NumTraits< Scalar >::dummy_precision()) const
Definition: Hyperplane.h:270