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_PARAMETRIZEDLINE_H
12#define EIGEN_PARAMETRIZEDLINE_H
13
14namespace Eigen {
15
16/** \geometry_module \ingroup Geometry_Module
17 *
18 * \class ParametrizedLine
19 *
20 * \brief A parametrized line
21 *
22 * A parametrized line is defined by an origin point \f$ \mathbf{o} \f$ and a unit
23 * direction vector \f$ \mathbf{d} \f$ such that the line corresponds to
24 * the set \f$ l(t) = \mathbf{o} + t \mathbf{d} \f$, \f$ t \in \mathbf{R} \f$.
25 *
26 * \tparam _Scalar the scalar type, i.e., the type of the coefficients
27 * \tparam _AmbientDim the dimension of the ambient space, can be a compile time value or Dynamic.
28 */
29template <typename _Scalar, int _AmbientDim, int _Options>
30class ParametrizedLine
31{
32public:
33 EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(_Scalar,_AmbientDim)
34 enum {
35 AmbientDimAtCompileTime = _AmbientDim,
36 Options = _Options
37 };
38 typedef _Scalar Scalar;
39 typedef typename NumTraits<Scalar>::Real RealScalar;
40 typedef Eigen::Index Index; ///< \deprecated since Eigen 3.3
41 typedef Matrix<Scalar,AmbientDimAtCompileTime,1,Options> VectorType;
42
43 /** Default constructor without initialization */
44 EIGEN_DEVICE_FUNC inline ParametrizedLine() {}
45
46 template<int OtherOptions>
47 EIGEN_DEVICE_FUNC ParametrizedLine(const ParametrizedLine<Scalar,AmbientDimAtCompileTime,OtherOptions>& other)
48 : m_origin(other.origin()), m_direction(other.direction())
49 {}
50
51 /** Constructs a dynamic-size line with \a _dim the dimension
52 * of the ambient space */
53 EIGEN_DEVICE_FUNC inline explicit ParametrizedLine(Index _dim) : m_origin(_dim), m_direction(_dim) {}
54
55 /** Initializes a parametrized line of direction \a direction and origin \a origin.
56 * \warning the vector direction is assumed to be normalized.
57 */
58 EIGEN_DEVICE_FUNC ParametrizedLine(const VectorType& origin, const VectorType& direction)
59 : m_origin(origin), m_direction(direction) {}
60
61 template <int OtherOptions>
62 EIGEN_DEVICE_FUNC explicit ParametrizedLine(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane);
63
64 /** Constructs a parametrized line going from \a p0 to \a p1. */
65 EIGEN_DEVICE_FUNC static inline ParametrizedLine Through(const VectorType& p0, const VectorType& p1)
66 { return ParametrizedLine(p0, (p1-p0).normalized()); }
67
68 EIGEN_DEVICE_FUNC ~ParametrizedLine() {}
69
70 /** \returns the dimension in which the line holds */
71 EIGEN_DEVICE_FUNC inline Index dim() const { return m_direction.size(); }
72
73 EIGEN_DEVICE_FUNC const VectorType& origin() const { return m_origin; }
74 EIGEN_DEVICE_FUNC VectorType& origin() { return m_origin; }
75
76 EIGEN_DEVICE_FUNC const VectorType& direction() const { return m_direction; }
77 EIGEN_DEVICE_FUNC VectorType& direction() { return m_direction; }
78
79 /** \returns the squared distance of a point \a p to its projection onto the line \c *this.
80 * \sa distance()
81 */
82 EIGEN_DEVICE_FUNC RealScalar squaredDistance(const VectorType& p) const
83 {
84 VectorType diff = p - origin();
85 return (diff - direction().dot(diff) * direction()).squaredNorm();
86 }
87 /** \returns the distance of a point \a p to its projection onto the line \c *this.
88 * \sa squaredDistance()
89 */
90 EIGEN_DEVICE_FUNC RealScalar distance(const VectorType& p) const { EIGEN_USING_STD_MATH(sqrt) return sqrt(squaredDistance(p)); }
91
92 /** \returns the projection of a point \a p onto the line \c *this. */
93 EIGEN_DEVICE_FUNC VectorType projection(const VectorType& p) const
94 { return origin() + direction().dot(p-origin()) * direction(); }
95
96 EIGEN_DEVICE_FUNC VectorType pointAt(const Scalar& t) const;
97
98 template <int OtherOptions>
99 EIGEN_DEVICE_FUNC Scalar intersectionParameter(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const;
100
101 template <int OtherOptions>
102 EIGEN_DEVICE_FUNC Scalar intersection(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const;
103
104 template <int OtherOptions>
105 EIGEN_DEVICE_FUNC VectorType intersectionPoint(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const;
106
107 /** \returns \c *this with scalar type casted to \a NewScalarType
108 *
109 * Note that if \a NewScalarType is equal to the current scalar type of \c *this
110 * then this function smartly returns a const reference to \c *this.
111 */
112 template<typename NewScalarType>
113 EIGEN_DEVICE_FUNC inline typename internal::cast_return_type<ParametrizedLine,
114 ParametrizedLine<NewScalarType,AmbientDimAtCompileTime,Options> >::type cast() const
115 {
116 return typename internal::cast_return_type<ParametrizedLine,
117 ParametrizedLine<NewScalarType,AmbientDimAtCompileTime,Options> >::type(*this);
118 }
119
120 /** Copy constructor with scalar type conversion */
121 template<typename OtherScalarType,int OtherOptions>
122 EIGEN_DEVICE_FUNC inline explicit ParametrizedLine(const ParametrizedLine<OtherScalarType,AmbientDimAtCompileTime,OtherOptions>& other)
123 {
124 m_origin = other.origin().template cast<Scalar>();
125 m_direction = other.direction().template cast<Scalar>();
126 }
127
128 /** \returns \c true if \c *this is approximately equal to \a other, within the precision
129 * determined by \a prec.
130 *
131 * \sa MatrixBase::isApprox() */
132 EIGEN_DEVICE_FUNC bool isApprox(const ParametrizedLine& other, const typename NumTraits<Scalar>::Real& prec = NumTraits<Scalar>::dummy_precision()) const
133 { return m_origin.isApprox(other.m_origin, prec) && m_direction.isApprox(other.m_direction, prec); }
134
135protected:
136
137 VectorType m_origin, m_direction;
138};
139
140/** Constructs a parametrized line from a 2D hyperplane
141 *
142 * \warning the ambient space must have dimension 2 such that the hyperplane actually describes a line
143 */
144template <typename _Scalar, int _AmbientDim, int _Options>
145template <int OtherOptions>
146EIGEN_DEVICE_FUNC inline ParametrizedLine<_Scalar, _AmbientDim,_Options>::ParametrizedLine(const Hyperplane<_Scalar, _AmbientDim,OtherOptions>& hyperplane)
147{
148 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 2)
149 direction() = hyperplane.normal().unitOrthogonal();
150 origin() = -hyperplane.normal()*hyperplane.offset();
151}
152
153/** \returns the point at \a t along this line
154 */
155template <typename _Scalar, int _AmbientDim, int _Options>
156EIGEN_DEVICE_FUNC inline typename ParametrizedLine<_Scalar, _AmbientDim,_Options>::VectorType
157ParametrizedLine<_Scalar, _AmbientDim,_Options>::pointAt(const _Scalar& t) const
158{
159 return origin() + (direction()*t);
160}
161
162/** \returns the parameter value of the intersection between \c *this and the given \a hyperplane
163 */
164template <typename _Scalar, int _AmbientDim, int _Options>
165template <int OtherOptions>
166EIGEN_DEVICE_FUNC inline _Scalar ParametrizedLine<_Scalar, _AmbientDim,_Options>::intersectionParameter(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const
167{
168 return -(hyperplane.offset()+hyperplane.normal().dot(origin()))
169 / hyperplane.normal().dot(direction());
170}
171
172
173/** \deprecated use intersectionParameter()
174 * \returns the parameter value of the intersection between \c *this and the given \a hyperplane
175 */
176template <typename _Scalar, int _AmbientDim, int _Options>
177template <int OtherOptions>
178EIGEN_DEVICE_FUNC inline _Scalar ParametrizedLine<_Scalar, _AmbientDim,_Options>::intersection(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const
179{
180 return intersectionParameter(hyperplane);
181}
182
183/** \returns the point of the intersection between \c *this and the given hyperplane
184 */
185template <typename _Scalar, int _AmbientDim, int _Options>
186template <int OtherOptions>
187EIGEN_DEVICE_FUNC inline typename ParametrizedLine<_Scalar, _AmbientDim,_Options>::VectorType
188ParametrizedLine<_Scalar, _AmbientDim,_Options>::intersectionPoint(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const
189{
190 return pointAt(intersectionParameter(hyperplane));
191}
192
193} // end namespace Eigen
194
195#endif // EIGEN_PARAMETRIZEDLINE_H
196