Hyperplane.h
Go to the documentation of this file.
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 #include "./InternalHeaderCheck.h"
15 
16 namespace Eigen {
17 
35 template <typename Scalar_, int AmbientDim_, int Options_>
37 {
38 public:
39  EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(Scalar_,AmbientDim_==Dynamic ? Dynamic : AmbientDim_+1)
40  enum {
41  AmbientDimAtCompileTime = AmbientDim_,
42  Options = Options_
43  };
44  typedef Scalar_ Scalar;
46  typedef Eigen::Index Index;
49  ? Dynamic
53 
56 
57  template<int OtherOptions>
59  : m_coeffs(other.coeffs())
60  {}
61 
64  EIGEN_DEVICE_FUNC inline explicit Hyperplane(Index _dim) : m_coeffs(_dim+1) {}
65 
70  : m_coeffs(n.size()+1)
71  {
72  normal() = n;
73  offset() = -n.dot(e);
74  }
75 
80  EIGEN_DEVICE_FUNC inline Hyperplane(const VectorType& n, const Scalar& d)
81  : m_coeffs(n.size()+1)
82  {
83  normal() = n;
84  offset() = d;
85  }
86 
91  {
92  Hyperplane result(p0.size());
93  result.normal() = (p1 - p0).unitOrthogonal();
94  result.offset() = -p0.dot(result.normal());
95  return result;
96  }
97 
101  EIGEN_DEVICE_FUNC static inline Hyperplane Through(const VectorType& p0, const VectorType& p1, const VectorType& p2)
102  {
104  Hyperplane result(p0.size());
105  VectorType v0(p2 - p0), v1(p1 - p0);
106  result.normal() = v0.cross(v1);
107  RealScalar norm = result.normal().norm();
108  if(norm <= v0.norm() * v1.norm() * NumTraits<RealScalar>::epsilon())
109  {
110  Matrix<Scalar,2,3> m; m << v0.transpose(), v1.transpose();
112  result.normal() = svd.matrixV().col(2);
113  }
114  else
115  result.normal() /= norm;
116  result.offset() = -p0.dot(result.normal());
117  return result;
118  }
119 
124  // FIXME to be consistent with the rest this could be implemented as a static Through function ??
126  {
127  normal() = parametrized.direction().unitOrthogonal();
128  offset() = -parametrized.origin().dot(normal());
129  }
130 
132 
135 
138  {
139  m_coeffs /= normal().norm();
140  }
141 
145  EIGEN_DEVICE_FUNC inline Scalar signedDistance(const VectorType& p) const { return normal().dot(p) + offset(); }
146 
151 
154  EIGEN_DEVICE_FUNC inline VectorType projection(const VectorType& p) const { return p - signedDistance(p) * normal(); }
155 
160 
165 
169  EIGEN_DEVICE_FUNC inline const Scalar& offset() const { return m_coeffs.coeff(dim()); }
170 
173  EIGEN_DEVICE_FUNC inline Scalar& offset() { return m_coeffs(dim()); }
174 
178  EIGEN_DEVICE_FUNC inline const Coefficients& coeffs() const { return m_coeffs; }
179 
184 
192  {
194  Scalar det = coeffs().coeff(0) * other.coeffs().coeff(1) - coeffs().coeff(1) * other.coeffs().coeff(0);
195  // since the line equations ax+by=c are normalized with a^2+b^2=1, the following tests
196  // whether the two lines are approximately parallel.
198  { // special case where the two lines are approximately parallel. Pick any point on the first line.
199  if(numext::abs(coeffs().coeff(1))>numext::abs(coeffs().coeff(0)))
200  return VectorType(coeffs().coeff(1), -coeffs().coeff(2)/coeffs().coeff(1)-coeffs().coeff(0));
201  else
202  return VectorType(-coeffs().coeff(2)/coeffs().coeff(0)-coeffs().coeff(1), coeffs().coeff(0));
203  }
204  else
205  { // general case
206  Scalar invdet = Scalar(1) / det;
207  return VectorType(invdet*(coeffs().coeff(1)*other.coeffs().coeff(2)-other.coeffs().coeff(1)*coeffs().coeff(2)),
208  invdet*(other.coeffs().coeff(0)*coeffs().coeff(2)-coeffs().coeff(0)*other.coeffs().coeff(2)));
209  }
210  }
211 
218  template<typename XprType>
220  {
221  if (traits==Affine)
222  {
223  normal() = mat.inverse().transpose() * normal();
224  m_coeffs /= normal().norm();
225  }
226  else if (traits==Isometry)
227  normal() = mat * normal();
228  else
229  {
230  eigen_assert(0 && "invalid traits value in Hyperplane::transform()");
231  }
232  return *this;
233  }
234 
242  template<int TrOptions>
244  TransformTraits traits = Affine)
245  {
246  transform(t.linear(), traits);
247  offset() -= normal().dot(t.translation());
248  return *this;
249  }
250 
256  template<typename NewScalarType>
257  EIGEN_DEVICE_FUNC inline typename internal::cast_return_type<Hyperplane,
259  {
260  return typename internal::cast_return_type<Hyperplane,
262  }
263 
265  template<typename OtherScalarType,int OtherOptions>
267  { m_coeffs = other.coeffs().template cast<Scalar>(); }
268 
273  template<int OtherOptions>
275  { return m_coeffs.isApprox(other.m_coeffs, prec); }
276 
277 protected:
278 
280 };
281 
282 } // end namespace Eigen
283 
284 #endif // EIGEN_HYPERPLANE_H
Matrix3f m
int n
Array< double, 1, 3 > e(1./3., 0.5, 2.)
cout<< "Here is the matrix m:"<< endl<< m<< endl;JacobiSVD< MatrixXf, ComputeThinU|ComputeThinV > svd(m)
#define EIGEN_DEVICE_FUNC
Definition: Macros.h:883
#define eigen_assert(x)
Definition: Macros.h:902
Vector3f p0
Vector3f p1
#define EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(Scalar, Size)
Definition: Memory.h:921
#define EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(TYPE, SIZE)
Definition: StaticAssert.h:51
float * p
M1<< 1, 2, 3, 4, 5, 6, 7, 8, 9;Map< RowVectorXf > v1(M1.data(), M1.size())
Expression of a fixed-size or dynamic-size block.
Definition: Block.h:107
A hyperplane.
Definition: Hyperplane.h:37
static Hyperplane Through(const VectorType &p0, const VectorType &p1)
Definition: Hyperplane.h:90
VectorType intersection(const Hyperplane &other) const
Definition: Hyperplane.h:191
ConstNormalReturnType normal() const
Definition: Hyperplane.h:159
Hyperplane & transform(const Transform< Scalar, AmbientDimAtCompileTime, Affine, TrOptions > &t, TransformTraits traits=Affine)
Definition: Hyperplane.h:243
Scalar signedDistance(const VectorType &p) const
Definition: Hyperplane.h:145
Block< Coefficients, AmbientDimAtCompileTime, 1 > NormalReturnType
Definition: Hyperplane.h:51
Hyperplane(const Hyperplane< Scalar, AmbientDimAtCompileTime, OtherOptions > &other)
Definition: Hyperplane.h:58
internal::cast_return_type< Hyperplane, Hyperplane< NewScalarType, AmbientDimAtCompileTime, Options > >::type cast() const
Definition: Hyperplane.h:258
Hyperplane(const VectorType &n, const VectorType &e)
Definition: Hyperplane.h:69
NormalReturnType normal()
Definition: Hyperplane.h:164
Scalar & offset()
Definition: Hyperplane.h:173
const Block< const Coefficients, AmbientDimAtCompileTime, 1 > ConstNormalReturnType
Definition: Hyperplane.h:52
const Scalar & offset() const
Definition: Hyperplane.h:169
Index dim() const
Definition: Hyperplane.h:134
Scalar absDistance(const VectorType &p) const
Definition: Hyperplane.h:150
static Hyperplane Through(const VectorType &p0, const VectorType &p1, const VectorType &p2)
Definition: Hyperplane.h:101
Hyperplane(const VectorType &n, const Scalar &d)
Definition: Hyperplane.h:80
bool isApprox(const Hyperplane< Scalar, AmbientDimAtCompileTime, OtherOptions > &other, const typename NumTraits< Scalar >::Real &prec=NumTraits< Scalar >::dummy_precision()) const
Definition: Hyperplane.h:274
Hyperplane(const Hyperplane< OtherScalarType, AmbientDimAtCompileTime, OtherOptions > &other)
Definition: Hyperplane.h:266
Coefficients & coeffs()
Definition: Hyperplane.h:183
NumTraits< Scalar >::Real RealScalar
Definition: Hyperplane.h:45
void normalize(void)
Definition: Hyperplane.h:137
Matrix< Scalar, Index(AmbientDimAtCompileTime)==Dynamic ? Dynamic :Index(AmbientDimAtCompileTime)+1, 1, Options > Coefficients
Definition: Hyperplane.h:50
const Coefficients & coeffs() const
Definition: Hyperplane.h:178
Hyperplane(Index _dim)
Definition: Hyperplane.h:64
Matrix< Scalar, AmbientDimAtCompileTime, 1 > VectorType
Definition: Hyperplane.h:47
Hyperplane(const ParametrizedLine< Scalar, AmbientDimAtCompileTime > &parametrized)
Definition: Hyperplane.h:125
VectorType projection(const VectorType &p) const
Definition: Hyperplane.h:154
Eigen::Index Index
Definition: Hyperplane.h:46
Coefficients m_coeffs
Definition: Hyperplane.h:279
Hyperplane & transform(const MatrixBase< XprType > &mat, TransformTraits traits=Affine)
Definition: Hyperplane.h:219
Two-sided Jacobi SVD decomposition of a rectangular matrix.
Definition: JacobiSVD.h:514
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
A parametrized line.
const VectorType & direction() const
const VectorType & origin() const
constexpr const Scalar & coeff(Index rowId, Index colId) const
Represents an homogeneous transformation in a N dimensional space.
Definition: Transform.h:207
ConstLinearPart linear() const
Definition: Transform.h:398
ConstTranslationPart translation() const
Definition: Transform.h:408
TransformTraits
Definition: Constants.h:459
@ ComputeFullV
Definition: Constants.h:399
@ Affine
Definition: Constants.h:464
@ Isometry
Definition: Constants.h:461
bool isMuchSmallerThan(const Scalar &x, const OtherScalar &y, const typename NumTraits< Scalar >::Real &precision=NumTraits< Scalar >::dummy_precision())
EIGEN_ALWAYS_INLINE std::enable_if_t< NumTraits< T >::IsSigned||NumTraits< T >::IsComplex, typename NumTraits< T >::Real > abs(const T &x)
: InteropHeaders
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
const int Dynamic
Definition: Constants.h:24
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:231