11 #ifndef EIGEN_HYPERPLANE_H 12 #define EIGEN_HYPERPLANE_H 33 template <
typename _Scalar,
int _AmbientDim,
int _Options>
37 EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(_Scalar,_AmbientDim==
Dynamic ?
Dynamic : _AmbientDim+1)
39 AmbientDimAtCompileTime = _AmbientDim,
42 typedef _Scalar Scalar;
43 typedef typename NumTraits<Scalar>::Real RealScalar;
55 template<
int OtherOptions>
57 : m_coeffs(other.coeffs())
62 EIGEN_DEVICE_FUNC
inline explicit Hyperplane(Index _dim) : m_coeffs(_dim+1) {}
67 EIGEN_DEVICE_FUNC
inline Hyperplane(
const VectorType& n,
const VectorType& e)
68 : m_coeffs(n.size()+1)
78 EIGEN_DEVICE_FUNC
inline Hyperplane(
const VectorType& n,
const Scalar& d)
79 : m_coeffs(n.size()+1)
88 EIGEN_DEVICE_FUNC
static inline Hyperplane
Through(
const VectorType& p0,
const VectorType& p1)
90 Hyperplane result(p0.size());
91 result.
normal() = (p1 - p0).unitOrthogonal();
92 result.offset() = -p0.dot(result.normal());
99 EIGEN_DEVICE_FUNC
static inline Hyperplane
Through(
const VectorType& p0,
const VectorType& p1,
const VectorType& p2)
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();
110 result.normal() = svd.
matrixV().col(2);
113 result.normal() /= norm;
114 result.offset() = -p0.dot(result.normal());
125 normal() = parametrized.direction().unitOrthogonal();
129 EIGEN_DEVICE_FUNC ~Hyperplane() {}
132 EIGEN_DEVICE_FUNC
inline Index
dim()
const {
return AmbientDimAtCompileTime==
Dynamic ? m_coeffs.size()-1 :
Index(AmbientDimAtCompileTime); }
137 m_coeffs /=
normal().norm();
157 EIGEN_DEVICE_FUNC
inline ConstNormalReturnType
normal()
const {
return ConstNormalReturnType(m_coeffs,0,0,
dim(),1); }
162 EIGEN_DEVICE_FUNC
inline NormalReturnType
normal() {
return NormalReturnType(m_coeffs,0,0,
dim(),1); }
167 EIGEN_DEVICE_FUNC
inline const Scalar&
offset()
const {
return m_coeffs.
coeff(
dim()); }
171 EIGEN_DEVICE_FUNC
inline Scalar&
offset() {
return m_coeffs(
dim()); }
176 EIGEN_DEVICE_FUNC
inline const Coefficients&
coeffs()
const {
return m_coeffs; }
181 EIGEN_DEVICE_FUNC
inline Coefficients&
coeffs() {
return m_coeffs; }
189 EIGEN_DEVICE_FUNC VectorType
intersection(
const Hyperplane& other)
const 191 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 2)
195 if(internal::isMuchSmallerThan(det, Scalar(1)))
197 if(numext::abs(
coeffs().coeff(1))>numext::abs(
coeffs().coeff(0)))
204 Scalar invdet = Scalar(1) / det;
205 return VectorType(invdet*(
coeffs().coeff(1)*other.coeffs().coeff(2)-other.coeffs().coeff(1)*
coeffs().
coeff(2)),
216 template<
typename XprType>
225 eigen_assert(0 &&
"invalid traits value in Hyperplane::transform()");
237 template<
int TrOptions>
251 template<
typename NewScalarType>
260 template<
typename OtherScalarType,
int OtherOptions>
262 { m_coeffs = other.coeffs().template cast<Scalar>(); }
268 template<
int OtherOptions>
270 {
return m_coeffs.isApprox(other.m_coeffs, prec); }
274 Coefficients m_coeffs;
279 #endif // EIGEN_HYPERPLANE_H EIGEN_DEVICE_FUNC VectorType intersection(const Hyperplane &other) const
Definition: Hyperplane.h:189
Definition: XprHelper.h:489
EIGEN_DEVICE_FUNC Hyperplane & transform(const MatrixBase< XprType > &mat, TransformTraits traits=Affine)
Applies the transformation matrix mat to *this and returns a reference to *this.
Definition: Hyperplane.h:217
static EIGEN_DEVICE_FUNC Hyperplane Through(const VectorType &p0, const VectorType &p1)
Constructs a hyperplane passing through the two points.
Definition: Hyperplane.h:88
EIGEN_DEVICE_FUNC Index dim() const
Definition: Hyperplane.h:132
EIGEN_DEVICE_FUNC VectorType projection(const VectorType &p) const
Definition: Hyperplane.h:152
EIGEN_DEVICE_FUNC Hyperplane(const Hyperplane< OtherScalarType, AmbientDimAtCompileTime, OtherOptions > &other)
Copy constructor with scalar type conversion.
Definition: Hyperplane.h:261
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:85
Transformation is an isometry.
Definition: Constants.h:447
Holds information about the various numeric (i.e.
Definition: NumTraits.h:150
EIGEN_DEVICE_FUNC Hyperplane(const ParametrizedLine< Scalar, AmbientDimAtCompileTime > ¶metrized)
Constructs a hyperplane passing through the parametrized line parametrized.
Definition: Hyperplane.h:123
EIGEN_DEVICE_FUNC const Coefficients & coeffs() const
Definition: Hyperplane.h:176
EIGEN_DEVICE_FUNC Hyperplane(Index _dim)
Constructs a dynamic-size hyperplane with _dim the dimension of the ambient space.
Definition: Hyperplane.h:62
EIGEN_DEVICE_FUNC Scalar absDistance(const VectorType &p) const
Definition: Hyperplane.h:148
EIGEN_DEVICE_FUNC Hyperplane()
Default constructor without initialization.
Definition: Hyperplane.h:53
Definition: ForwardDeclarations.h:276
TransformTraits
Enum used to specify how a particular transformation is stored in a matrix.
Definition: Constants.h:445
Eigen::Index Index
Definition: Hyperplane.h:44
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Scalar & coeff(Index rowId, Index colId) const
This is an overloaded version of DenseCoeffsBase<Derived,ReadOnlyAccessors>::coeff(Index,Index) const provided to by-pass the creation of an evaluator of the expression, thus saving compilation efforts.
Definition: PlainObjectBase.h:154
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
EIGEN_DEVICE_FUNC Scalar signedDistance(const VectorType &p) const
Definition: Hyperplane.h:143
EIGEN_DEVICE_FUNC ConstNormalReturnType normal() const
Definition: Hyperplane.h:157
EIGEN_DEVICE_FUNC NormalReturnType normal()
Definition: Hyperplane.h:162
EIGEN_DEVICE_FUNC Hyperplane & transform(const Transform< Scalar, AmbientDimAtCompileTime, Affine, TrOptions > &t, TransformTraits traits=Affine)
Applies the transformation t to *this and returns a reference to *this.
Definition: Hyperplane.h:238
EIGEN_DEVICE_FUNC const Scalar & offset() const
Definition: Hyperplane.h:167
EIGEN_DEVICE_FUNC Coefficients & coeffs()
Definition: Hyperplane.h:181
static EIGEN_DEVICE_FUNC Hyperplane Through(const VectorType &p0, const VectorType &p1, const VectorType &p2)
Constructs a hyperplane passing through the three points.
Definition: Hyperplane.h:99
const Inverse< Derived > inverse() const
Definition: InverseImpl.h:335
Expression of a fixed-size or dynamic-size block.
Definition: Block.h:103
EIGEN_DEVICE_FUNC void normalize(void)
normalizes *this
Definition: Hyperplane.h:135
Two-sided Jacobi SVD decomposition of a rectangular matrix.
Definition: ForwardDeclarations.h:258
EIGEN_DEVICE_FUNC Hyperplane(const VectorType &n, const VectorType &e)
Construct a plane from its normal n and a point e onto the plane.
Definition: Hyperplane.h:67
EIGEN_DEVICE_FUNC Hyperplane(const VectorType &n, const Scalar &d)
Constructs a plane from its normal n and distance to the origin d such that the algebraic equation of...
Definition: Hyperplane.h:78
const MatrixVType & matrixV() const
Definition: SVDBase.h:99
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:269
const int Dynamic
This value means that a positive quantity (e.g., a size) is not known at compile-time, and that instead the value is stored in some runtime variable.
Definition: Constants.h:21
Used in JacobiSVD to indicate that the square matrix V is to be computed.
Definition: Constants.h:387
Transformation is an affine transformation stored as a (Dim+1)^2 matrix whose last row is assumed to ...
Definition: Constants.h:450
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
Definition: ForwardDeclarations.h:275
EIGEN_DEVICE_FUNC internal::cast_return_type< Hyperplane, Hyperplane< NewScalarType, AmbientDimAtCompileTime, Options > >::type cast() const
Definition: Hyperplane.h:253
EIGEN_DEVICE_FUNC Scalar & offset()
Definition: Hyperplane.h:171