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 // IWYU pragma: private
15 #include "./InternalHeaderCheck.h"
16 
17 namespace Eigen {
18 
36 template <typename Scalar_, int AmbientDim_, int Options_>
37 class Hyperplane {
38  public:
40  AmbientDim_ == Dynamic ? Dynamic : AmbientDim_ + 1)
41  enum { AmbientDimAtCompileTime = AmbientDim_, Options = Options_ };
42  typedef Scalar_ Scalar;
44  typedef Eigen::Index Index;
46  typedef Matrix<Scalar, Index(AmbientDimAtCompileTime) == Dynamic ? Dynamic : Index(AmbientDimAtCompileTime) + 1, 1,
47  Options>
51 
54 
55  template <int OtherOptions>
57  : m_coeffs(other.coeffs()) {}
58 
61  EIGEN_DEVICE_FUNC inline explicit Hyperplane(Index _dim) : m_coeffs(_dim + 1) {}
62 
66  EIGEN_DEVICE_FUNC inline Hyperplane(const VectorType& n, const VectorType& e) : m_coeffs(n.size() + 1) {
67  normal() = n;
68  offset() = -n.dot(e);
69  }
70 
75  EIGEN_DEVICE_FUNC inline Hyperplane(const VectorType& n, const Scalar& d) : m_coeffs(n.size() + 1) {
76  normal() = n;
77  offset() = d;
78  }
79 
83  EIGEN_DEVICE_FUNC static inline Hyperplane Through(const VectorType& p0, const VectorType& p1) {
84  Hyperplane result(p0.size());
85  result.normal() = (p1 - p0).unitOrthogonal();
86  result.offset() = -p0.dot(result.normal());
87  return result;
88  }
89 
93  EIGEN_DEVICE_FUNC static inline Hyperplane Through(const VectorType& p0, const VectorType& p1, const VectorType& p2) {
95  Hyperplane result(p0.size());
96  VectorType v0(p2 - p0), v1(p1 - p0);
97  result.normal() = v0.cross(v1);
98  RealScalar norm = result.normal().norm();
99  if (norm <= v0.norm() * v1.norm() * NumTraits<RealScalar>::epsilon()) {
101  m << v0.transpose(), v1.transpose();
103  result.normal() = svd.matrixV().col(2);
104  } else
105  result.normal() /= norm;
106  result.offset() = -p0.dot(result.normal());
107  return result;
108  }
109 
114  // FIXME to be consistent with the rest this could be implemented as a static Through function ??
116  normal() = parametrized.direction().unitOrthogonal();
117  offset() = -parametrized.origin().dot(normal());
118  }
119 
121 
123  EIGEN_DEVICE_FUNC inline Index dim() const {
124  return AmbientDimAtCompileTime == Dynamic ? m_coeffs.size() - 1 : Index(AmbientDimAtCompileTime);
125  }
126 
128  EIGEN_DEVICE_FUNC void normalize(void) { m_coeffs /= normal().norm(); }
129 
133  EIGEN_DEVICE_FUNC inline Scalar signedDistance(const VectorType& p) const { return normal().dot(p) + offset(); }
134 
139 
142  EIGEN_DEVICE_FUNC inline VectorType projection(const VectorType& p) const { return p - signedDistance(p) * normal(); }
143 
148  return ConstNormalReturnType(m_coeffs, 0, 0, dim(), 1);
149  }
150 
155 
159  EIGEN_DEVICE_FUNC inline const Scalar& offset() const { return m_coeffs.coeff(dim()); }
160 
163  EIGEN_DEVICE_FUNC inline Scalar& offset() { return m_coeffs(dim()); }
164 
168  EIGEN_DEVICE_FUNC inline const Coefficients& coeffs() const { return m_coeffs; }
169 
174 
183  Scalar det = coeffs().coeff(0) * other.coeffs().coeff(1) - coeffs().coeff(1) * other.coeffs().coeff(0);
184  // since the line equations ax+by=c are normalized with a^2+b^2=1, the following tests
185  // whether the two lines are approximately parallel.
186  if (internal::isMuchSmallerThan(det, Scalar(1))) { // special case where the two lines are approximately parallel.
187  // Pick any point on the first line.
188  if (numext::abs(coeffs().coeff(1)) > numext::abs(coeffs().coeff(0)))
189  return VectorType(coeffs().coeff(1), -coeffs().coeff(2) / coeffs().coeff(1) - coeffs().coeff(0));
190  else
191  return VectorType(-coeffs().coeff(2) / coeffs().coeff(0) - coeffs().coeff(1), coeffs().coeff(0));
192  } else { // general case
193  Scalar invdet = Scalar(1) / det;
194  return VectorType(
195  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  if (traits == Affine) {
209  normal() = mat.inverse().transpose() * normal();
210  m_coeffs /= normal().norm();
211  } else if (traits == Isometry)
212  normal() = mat * normal();
213  else {
214  eigen_assert(0 && "invalid traits value in Hyperplane::transform()");
215  }
216  return *this;
217  }
218 
226  template <int TrOptions>
228  TransformTraits traits = Affine) {
229  transform(t.linear(), traits);
230  offset() -= normal().dot(t.translation());
231  return *this;
232  }
233 
239  template <typename NewScalarType>
240  EIGEN_DEVICE_FUNC inline
243  cast() const {
244  return
247  }
248 
250  template <typename OtherScalarType, int OtherOptions>
251  EIGEN_DEVICE_FUNC inline explicit Hyperplane(
253  m_coeffs = other.coeffs().template cast<Scalar>();
254  }
255 
260  template <int OtherOptions>
263  const typename NumTraits<Scalar>::Real& prec = NumTraits<Scalar>::dummy_precision()) const {
264  return m_coeffs.isApprox(other.m_coeffs, prec);
265  }
266 
267  protected:
269 };
270 
271 } // end namespace Eigen
272 
273 #endif // EIGEN_HYPERPLANE_H
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
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:892
#define eigen_assert(x)
Definition: Macros.h:910
Vector3f p0
Definition: MatrixBase_all.cpp:2
Vector3f p1
Definition: MatrixBase_all.cpp:2
#define EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(TYPE, SIZE)
Definition: StaticAssert.h:50
float * p
Definition: Tutorial_Map_using.cpp:9
M1<< 1, 2, 3, 4, 5, 6, 7, 8, 9;Map< RowVectorXf > v1(M1.data(), M1.size())
Scalar Scalar int size
Definition: benchVecAdd.cpp:17
Expression of a fixed-size or dynamic-size block.
Definition: Block.h:110
A hyperplane.
Definition: Hyperplane.h:37
EIGEN_DEVICE_FUNC Hyperplane(const VectorType &n, const VectorType &e)
Definition: Hyperplane.h:66
EIGEN_DEVICE_FUNC VectorType projection(const VectorType &p) const
Definition: Hyperplane.h:142
EIGEN_DEVICE_FUNC NormalReturnType normal()
Definition: Hyperplane.h:154
static EIGEN_DEVICE_FUNC Hyperplane Through(const VectorType &p0, const VectorType &p1)
Definition: Hyperplane.h:83
EIGEN_DEVICE_FUNC Hyperplane & transform(const Transform< Scalar, AmbientDimAtCompileTime, Affine, TrOptions > &t, TransformTraits traits=Affine)
Definition: Hyperplane.h:227
EIGEN_DEVICE_FUNC void normalize(void)
Definition: Hyperplane.h:128
EIGEN_DEVICE_FUNC ~Hyperplane()
Definition: Hyperplane.h:120
EIGEN_DEVICE_FUNC const Coefficients & coeffs() const
Definition: Hyperplane.h:168
EIGEN_DEVICE_FUNC internal::cast_return_type< Hyperplane, Hyperplane< NewScalarType, AmbientDimAtCompileTime, Options > >::type cast() const
Definition: Hyperplane.h:243
Matrix< Scalar, Index(AmbientDimAtCompileTime)==Dynamic ? Dynamic :Index(AmbientDimAtCompileTime)+1, 1, Options > Coefficients
Definition: Hyperplane.h:48
EIGEN_DEVICE_FUNC const Scalar & offset() const
Definition: Hyperplane.h:159
EIGEN_DEVICE_FUNC Index dim() const
Definition: Hyperplane.h:123
EIGEN_DEVICE_FUNC ConstNormalReturnType normal() const
Definition: Hyperplane.h:147
Scalar_ Scalar
Definition: Hyperplane.h:41
EIGEN_DEVICE_FUNC VectorType intersection(const Hyperplane &other) const
Definition: Hyperplane.h:181
EIGEN_DEVICE_FUNC Scalar signedDistance(const VectorType &p) const
Definition: Hyperplane.h:133
EIGEN_DEVICE_FUNC Hyperplane(const Hyperplane< Scalar, AmbientDimAtCompileTime, OtherOptions > &other)
Definition: Hyperplane.h:56
NumTraits< Scalar >::Real RealScalar
Definition: Hyperplane.h:43
Block< Coefficients, AmbientDimAtCompileTime, 1 > NormalReturnType
Definition: Hyperplane.h:49
EIGEN_DEVICE_FUNC Hyperplane(const Hyperplane< OtherScalarType, AmbientDimAtCompileTime, OtherOptions > &other)
Definition: Hyperplane.h:251
EIGEN_DEVICE_FUNC Hyperplane(Index _dim)
Definition: Hyperplane.h:61
EIGEN_DEVICE_FUNC Hyperplane(const ParametrizedLine< Scalar, AmbientDimAtCompileTime > &parametrized)
Definition: Hyperplane.h:115
EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(Scalar_, AmbientDim_==Dynamic ? Dynamic :AmbientDim_+1) enum
Definition: Hyperplane.h:39
static EIGEN_DEVICE_FUNC Hyperplane Through(const VectorType &p0, const VectorType &p1, const VectorType &p2)
Definition: Hyperplane.h:93
EIGEN_DEVICE_FUNC Scalar & offset()
Definition: Hyperplane.h:163
EIGEN_DEVICE_FUNC Hyperplane()
Definition: Hyperplane.h:53
EIGEN_DEVICE_FUNC Hyperplane(const VectorType &n, const Scalar &d)
Definition: Hyperplane.h:75
Matrix< Scalar, AmbientDimAtCompileTime, 1 > VectorType
Definition: Hyperplane.h:45
EIGEN_DEVICE_FUNC bool isApprox(const Hyperplane< Scalar, AmbientDimAtCompileTime, OtherOptions > &other, const typename NumTraits< Scalar >::Real &prec=NumTraits< Scalar >::dummy_precision()) const
Definition: Hyperplane.h:261
EIGEN_DEVICE_FUNC Hyperplane & transform(const MatrixBase< XprType > &mat, TransformTraits traits=Affine)
Definition: Hyperplane.h:207
EIGEN_DEVICE_FUNC Coefficients & coeffs()
Definition: Hyperplane.h:173
const Block< const Coefficients, AmbientDimAtCompileTime, 1 > ConstNormalReturnType
Definition: Hyperplane.h:50
Eigen::Index Index
Definition: Hyperplane.h:44
Coefficients m_coeffs
Definition: Hyperplane.h:268
EIGEN_DEVICE_FUNC Scalar absDistance(const VectorType &p) const
Definition: Hyperplane.h:138
Two-sided Jacobi SVD decomposition of a rectangular matrix.
Definition: JacobiSVD.h:500
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
A parametrized line.
Definition: ParametrizedLine.h:33
EIGEN_DEVICE_FUNC const VectorType & direction() const
Definition: ParametrizedLine.h:75
EIGEN_DEVICE_FUNC const VectorType & origin() const
Definition: ParametrizedLine.h:72
EIGEN_DEVICE_FUNC constexpr EIGEN_STRONG_INLINE const Scalar & coeff(Index rowId, Index colId) const
Definition: PlainObjectBase.h:198
TransposeReturnType transpose()
Definition: SparseMatrixBase.h:358
Represents an homogeneous transformation in a N dimensional space.
Definition: Transform.h:192
TransformTraits
Definition: Constants.h:453
@ ComputeFullV
Definition: Constants.h:393
@ Affine
Definition: Constants.h:458
@ Isometry
Definition: Constants.h:455
int * m
Definition: level2_cplx_impl.h:294
EIGEN_DEVICE_FUNC bool isMuchSmallerThan(const Scalar &x, const OtherScalar &y, const typename NumTraits< Scalar >::Real &precision=NumTraits< Scalar >::dummy_precision())
Definition: MathFunctions.h:1916
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE std::enable_if_t< NumTraits< T >::IsSigned||NumTraits< T >::IsComplex, typename NumTraits< T >::Real > abs(const T &x)
Definition: MathFunctions.h:1355
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
const int Dynamic
Definition: Constants.h:25
type
Definition: compute_granudrum_aor.py:141
t
Definition: plotPSD.py:36
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:217
Definition: XprHelper.h:583