11 #ifndef EIGEN_PARAMETRIZEDLINE_H
12 #define EIGEN_PARAMETRIZEDLINE_H
14 #include "./InternalHeaderCheck.h"
31 template <
typename Scalar_,
int AmbientDim_,
int Options_>
35 EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(Scalar_,AmbientDim_)
37 AmbientDimAtCompileTime = AmbientDim_,
40 typedef Scalar_ Scalar;
48 template<
int OtherOptions>
50 : m_origin(other.origin()), m_direction(other.direction())
61 : m_origin(origin), m_direction(direction) {}
63 template <
int OtherOptions>
73 EIGEN_DEVICE_FUNC
inline Index dim()
const {
return m_direction.size(); }
75 EIGEN_DEVICE_FUNC
const VectorType& origin()
const {
return m_origin; }
76 EIGEN_DEVICE_FUNC VectorType& origin() {
return m_origin; }
78 EIGEN_DEVICE_FUNC
const VectorType& direction()
const {
return m_direction; }
79 EIGEN_DEVICE_FUNC VectorType& direction() {
return m_direction; }
87 return (diff - direction().dot(diff) * direction()).squaredNorm();
96 {
return origin() + direction().dot(p-origin()) * direction(); }
98 EIGEN_DEVICE_FUNC VectorType
pointAt(
const Scalar& t)
const;
100 template <
int OtherOptions>
103 template <
int OtherOptions>
106 template <
int OtherOptions>
115 template<
typename XprType>
119 direction() = (mat * direction()).normalized();
121 direction() = mat * direction();
124 eigen_assert(0 &&
"invalid traits value in ParametrizedLine::transform()");
126 origin() = mat * origin();
137 template<
int TrOptions>
151 template<
typename NewScalarType>
152 EIGEN_DEVICE_FUNC
inline typename internal::cast_return_type<
ParametrizedLine,
160 template<
typename OtherScalarType,
int OtherOptions>
163 m_origin = other.origin().template cast<Scalar>();
164 m_direction = other.direction().template cast<Scalar>();
172 {
return m_origin.isApprox(other.m_origin, prec) && m_direction.isApprox(other.m_direction, prec); }
176 VectorType m_origin, m_direction;
183 template <
typename Scalar_,
int AmbientDim_,
int Options_>
184 template <
int OtherOptions>
187 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(
VectorType, 2)
188 direction() = hyperplane.
normal().unitOrthogonal();
194 template <
typename Scalar_,
int AmbientDim_,
int Options_>
198 return origin() + (direction()*t);
203 template <
typename Scalar_,
int AmbientDim_,
int Options_>
204 template <
int OtherOptions>
207 return -(hyperplane.
offset()+hyperplane.
normal().dot(origin()))
208 / hyperplane.
normal().dot(direction());
215 template <
typename Scalar_,
int AmbientDim_,
int Options_>
216 template <
int OtherOptions>
219 return intersectionParameter(hyperplane);
224 template <
typename Scalar_,
int AmbientDim_,
int Options_>
225 template <
int OtherOptions>
229 return pointAt(intersectionParameter(hyperplane));
A hyperplane.
Definition: Hyperplane.h:37
ConstNormalReturnType normal() const
Definition: Hyperplane.h:159
const Scalar & offset() const
Definition: Hyperplane.h:169
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
A parametrized line.
Definition: ParametrizedLine.h:33
ParametrizedLine & transform(const Transform< Scalar, AmbientDimAtCompileTime, Affine, TrOptions > &t, TransformTraits traits=Affine)
Definition: ParametrizedLine.h:138
ParametrizedLine & transform(const MatrixBase< XprType > &mat, TransformTraits traits=Affine)
Definition: ParametrizedLine.h:116
RealScalar squaredDistance(const VectorType &p) const
Definition: ParametrizedLine.h:84
RealScalar distance(const VectorType &p) const
Definition: ParametrizedLine.h:92
VectorType intersectionPoint(const Hyperplane< Scalar_, AmbientDim_, OtherOptions > &hyperplane) const
Definition: ParametrizedLine.h:227
ParametrizedLine(const VectorType &origin, const VectorType &direction)
Definition: ParametrizedLine.h:60
ParametrizedLine(Index _dim)
Definition: ParametrizedLine.h:55
bool isApprox(const ParametrizedLine &other, const typename NumTraits< Scalar >::Real &prec=NumTraits< Scalar >::dummy_precision()) const
Definition: ParametrizedLine.h:171
internal::cast_return_type< ParametrizedLine, ParametrizedLine< NewScalarType, AmbientDimAtCompileTime, Options > >::type cast() const
Definition: ParametrizedLine.h:153
VectorType pointAt(const Scalar &t) const
Definition: ParametrizedLine.h:196
ParametrizedLine()
Definition: ParametrizedLine.h:46
static ParametrizedLine Through(const VectorType &p0, const VectorType &p1)
Definition: ParametrizedLine.h:67
Index dim() const
Definition: ParametrizedLine.h:73
Eigen::Index Index
Definition: ParametrizedLine.h:42
ParametrizedLine(const ParametrizedLine< OtherScalarType, AmbientDimAtCompileTime, OtherOptions > &other)
Definition: ParametrizedLine.h:161
VectorType projection(const VectorType &p) const
Definition: ParametrizedLine.h:95
TransformTraits
Definition: Constants.h:457
@ Affine
Definition: Constants.h:462
@ Isometry
Definition: Constants.h:459
Namespace containing all symbols from the Eigen library.
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:59
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:231