11 #ifndef EIGEN_PARAMETRIZEDLINE_H 12 #define EIGEN_PARAMETRIZEDLINE_H 29 template <
typename _Scalar,
int _AmbientDim,
int _Options>
30 class ParametrizedLine
33 EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(_Scalar,_AmbientDim)
35 AmbientDimAtCompileTime = _AmbientDim,
38 typedef _Scalar Scalar;
39 typedef typename NumTraits<Scalar>::Real RealScalar;
46 template<
int OtherOptions>
48 : m_origin(other.origin()), m_direction(other.direction())
53 EIGEN_DEVICE_FUNC
inline explicit ParametrizedLine(Index _dim) : m_origin(_dim), m_direction(_dim) {}
58 EIGEN_DEVICE_FUNC
ParametrizedLine(
const VectorType& origin,
const VectorType& direction)
59 : m_origin(origin), m_direction(direction) {}
61 template <
int OtherOptions>
65 EIGEN_DEVICE_FUNC
static inline ParametrizedLine
Through(
const VectorType& p0,
const VectorType& p1)
68 EIGEN_DEVICE_FUNC ~ParametrizedLine() {}
71 EIGEN_DEVICE_FUNC
inline Index
dim()
const {
return m_direction.size(); }
73 EIGEN_DEVICE_FUNC
const VectorType& origin()
const {
return m_origin; }
74 EIGEN_DEVICE_FUNC VectorType& origin() {
return m_origin; }
76 EIGEN_DEVICE_FUNC
const VectorType& direction()
const {
return m_direction; }
77 EIGEN_DEVICE_FUNC VectorType& direction() {
return m_direction; }
84 VectorType diff = p - origin();
85 return (diff - direction().dot(diff) * direction()).squaredNorm();
90 EIGEN_DEVICE_FUNC RealScalar
distance(
const VectorType& p)
const { EIGEN_USING_STD_MATH(sqrt)
return sqrt(
squaredDistance(p)); }
93 EIGEN_DEVICE_FUNC VectorType
projection(
const VectorType& p)
const 94 {
return origin() + direction().dot(p-origin()) * direction(); }
96 EIGEN_DEVICE_FUNC VectorType
pointAt(
const Scalar& t)
const;
98 template <
int OtherOptions>
101 template <
int OtherOptions>
104 template <
int OtherOptions>
112 template<
typename NewScalarType>
121 template<
typename OtherScalarType,
int OtherOptions>
124 m_origin = other.origin().template cast<Scalar>();
125 m_direction = other.direction().template cast<Scalar>();
133 {
return m_origin.isApprox(other.m_origin, prec) && m_direction.isApprox(other.m_direction, prec); }
137 VectorType m_origin, m_direction;
144 template <
typename _Scalar,
int _AmbientDim,
int _Options>
145 template <
int OtherOptions>
148 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 2)
149 direction() = hyperplane.
normal().unitOrthogonal();
155 template <
typename _Scalar,
int _AmbientDim,
int _Options>
159 return origin() + (direction()*t);
164 template <
typename _Scalar,
int _AmbientDim,
int _Options>
165 template <
int OtherOptions>
168 return -(hyperplane.
offset()+hyperplane.
normal().dot(origin()))
169 / hyperplane.
normal().dot(direction());
176 template <
typename _Scalar,
int _AmbientDim,
int _Options>
177 template <
int OtherOptions>
180 return intersectionParameter(hyperplane);
185 template <
typename _Scalar,
int _AmbientDim,
int _Options>
186 template <
int OtherOptions>
190 return pointAt(intersectionParameter(hyperplane));
195 #endif // EIGEN_PARAMETRIZEDLINE_H Definition: XprHelper.h:489
EIGEN_DEVICE_FUNC Index dim() const
Definition: ParametrizedLine.h:71
EIGEN_DEVICE_FUNC ParametrizedLine(const VectorType &origin, const VectorType &direction)
Initializes a parametrized line of direction direction and origin origin.
Definition: ParametrizedLine.h:58
EIGEN_DEVICE_FUNC ParametrizedLine()
Default constructor without initialization.
Definition: ParametrizedLine.h:44
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:85
Holds information about the various numeric (i.e.
Definition: NumTraits.h:150
EIGEN_DEVICE_FUNC VectorType pointAt(const Scalar &t) const
Definition: ParametrizedLine.h:157
EIGEN_DEVICE_FUNC internal::cast_return_type< ParametrizedLine, ParametrizedLine< NewScalarType, AmbientDimAtCompileTime, Options > >::type cast() const
Definition: ParametrizedLine.h:114
Definition: ForwardDeclarations.h:276
EIGEN_DEVICE_FUNC VectorType projection(const VectorType &p) const
Definition: ParametrizedLine.h:93
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
Eigen::Index Index
Definition: ParametrizedLine.h:40
EIGEN_DEVICE_FUNC ConstNormalReturnType normal() const
Definition: Hyperplane.h:157
EIGEN_DEVICE_FUNC ParametrizedLine(const ParametrizedLine< OtherScalarType, AmbientDimAtCompileTime, OtherOptions > &other)
Copy constructor with scalar type conversion.
Definition: ParametrizedLine.h:122
EIGEN_DEVICE_FUNC RealScalar squaredDistance(const VectorType &p) const
Definition: ParametrizedLine.h:82
EIGEN_DEVICE_FUNC bool isApprox(const ParametrizedLine &other, const typename NumTraits< Scalar >::Real &prec=NumTraits< Scalar >::dummy_precision()) const
Definition: ParametrizedLine.h:132
EIGEN_DEVICE_FUNC const Scalar & offset() const
Definition: Hyperplane.h:167
EIGEN_DEVICE_FUNC RealScalar distance(const VectorType &p) const
Definition: ParametrizedLine.h:90
EIGEN_DEVICE_FUNC ParametrizedLine(Index _dim)
Constructs a dynamic-size line with _dim the dimension of the ambient space.
Definition: ParametrizedLine.h:53
static EIGEN_DEVICE_FUNC ParametrizedLine Through(const VectorType &p0, const VectorType &p1)
Constructs a parametrized line going from p0 to p1.
Definition: ParametrizedLine.h:65
EIGEN_DEVICE_FUNC VectorType intersectionPoint(const Hyperplane< _Scalar, _AmbientDim, OtherOptions > &hyperplane) const
Definition: ParametrizedLine.h:188
Definition: ForwardDeclarations.h:275