xref: /aosp_15_r20/external/eigen/Eigen/src/Geometry/ParametrizedLine.h (revision bf2c37156dfe67e5dfebd6d394bad8b2ab5804d4)
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008 Gael Guennebaud <[email protected]>
5 // Copyright (C) 2008 Benoit Jacob <[email protected]>
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 
14 namespace 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   */
29 template <typename _Scalar, int _AmbientDim, int _Options>
30 class ParametrizedLine
31 {
32 public:
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 */
ParametrizedLine()44   EIGEN_DEVICE_FUNC inline ParametrizedLine() {}
45 
46   template<int OtherOptions>
ParametrizedLine(const ParametrizedLine<Scalar,AmbientDimAtCompileTime,OtherOptions> & other)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 */
ParametrizedLine(Index _dim)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     */
ParametrizedLine(const VectorType & origin,const VectorType & direction)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. */
Through(const VectorType & p0,const VectorType & p1)65   EIGEN_DEVICE_FUNC static inline ParametrizedLine Through(const VectorType& p0, const VectorType& p1)
66   { return ParametrizedLine(p0, (p1-p0).normalized()); }
67 
~ParametrizedLine()68   EIGEN_DEVICE_FUNC ~ParametrizedLine() {}
69 
70   /** \returns the dimension in which the line holds */
dim()71   EIGEN_DEVICE_FUNC inline Index dim() const { return m_direction.size(); }
72 
origin()73   EIGEN_DEVICE_FUNC const VectorType& origin() const { return m_origin; }
origin()74   EIGEN_DEVICE_FUNC VectorType& origin() { return m_origin; }
75 
direction()76   EIGEN_DEVICE_FUNC const VectorType& direction() const { return m_direction; }
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     */
squaredDistance(const VectorType & p)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     */
distance(const VectorType & p)90   EIGEN_DEVICE_FUNC RealScalar distance(const VectorType& p) const { EIGEN_USING_STD(sqrt) return sqrt(squaredDistance(p)); }
91 
92   /** \returns the projection of a point \a p onto the line \c *this. */
projection(const VectorType & p)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   /** Applies the transformation matrix \a mat to \c *this and returns a reference to \c *this.
108     *
109     * \param mat the Dim x Dim transformation matrix
110     * \param traits specifies whether the matrix \a mat represents an #Isometry
111     *               or a more generic #Affine transformation. The default is #Affine.
112     */
113   template<typename XprType>
114   EIGEN_DEVICE_FUNC inline ParametrizedLine& transform(const MatrixBase<XprType>& mat, TransformTraits traits = Affine)
115   {
116     if (traits==Affine)
117       direction() = (mat * direction()).normalized();
118     else if (traits==Isometry)
119       direction() = mat * direction();
120     else
121     {
122       eigen_assert(0 && "invalid traits value in ParametrizedLine::transform()");
123     }
124     origin() = mat * origin();
125     return *this;
126   }
127 
128   /** Applies the transformation \a t to \c *this and returns a reference to \c *this.
129     *
130     * \param t the transformation of dimension Dim
131     * \param traits specifies whether the transformation \a t represents an #Isometry
132     *               or a more generic #Affine transformation. The default is #Affine.
133     *               Other kind of transformations are not supported.
134     */
135   template<int TrOptions>
136   EIGEN_DEVICE_FUNC inline ParametrizedLine& transform(const Transform<Scalar,AmbientDimAtCompileTime,Affine,TrOptions>& t,
137                                                        TransformTraits traits = Affine)
138   {
139     transform(t.linear(), traits);
140     origin() += t.translation();
141     return *this;
142   }
143 
144 /** \returns \c *this with scalar type casted to \a NewScalarType
145     *
146     * Note that if \a NewScalarType is equal to the current scalar type of \c *this
147     * then this function smartly returns a const reference to \c *this.
148     */
149   template<typename NewScalarType>
150   EIGEN_DEVICE_FUNC inline typename internal::cast_return_type<ParametrizedLine,
cast()151            ParametrizedLine<NewScalarType,AmbientDimAtCompileTime,Options> >::type cast() const
152   {
153     return typename internal::cast_return_type<ParametrizedLine,
154                     ParametrizedLine<NewScalarType,AmbientDimAtCompileTime,Options> >::type(*this);
155   }
156 
157   /** Copy constructor with scalar type conversion */
158   template<typename OtherScalarType,int OtherOptions>
ParametrizedLine(const ParametrizedLine<OtherScalarType,AmbientDimAtCompileTime,OtherOptions> & other)159   EIGEN_DEVICE_FUNC inline explicit ParametrizedLine(const ParametrizedLine<OtherScalarType,AmbientDimAtCompileTime,OtherOptions>& other)
160   {
161     m_origin = other.origin().template cast<Scalar>();
162     m_direction = other.direction().template cast<Scalar>();
163   }
164 
165   /** \returns \c true if \c *this is approximately equal to \a other, within the precision
166     * determined by \a prec.
167     *
168     * \sa MatrixBase::isApprox() */
169   EIGEN_DEVICE_FUNC bool isApprox(const ParametrizedLine& other, const typename NumTraits<Scalar>::Real& prec = NumTraits<Scalar>::dummy_precision()) const
170   { return m_origin.isApprox(other.m_origin, prec) && m_direction.isApprox(other.m_direction, prec); }
171 
172 protected:
173 
174   VectorType m_origin, m_direction;
175 };
176 
177 /** Constructs a parametrized line from a 2D hyperplane
178   *
179   * \warning the ambient space must have dimension 2 such that the hyperplane actually describes a line
180   */
181 template <typename _Scalar, int _AmbientDim, int _Options>
182 template <int OtherOptions>
ParametrizedLine(const Hyperplane<_Scalar,_AmbientDim,OtherOptions> & hyperplane)183 EIGEN_DEVICE_FUNC inline ParametrizedLine<_Scalar, _AmbientDim,_Options>::ParametrizedLine(const Hyperplane<_Scalar, _AmbientDim,OtherOptions>& hyperplane)
184 {
185   EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 2)
186   direction() = hyperplane.normal().unitOrthogonal();
187   origin() = -hyperplane.normal()*hyperplane.offset();
188 }
189 
190 /** \returns the point at \a t along this line
191   */
192 template <typename _Scalar, int _AmbientDim, int _Options>
193 EIGEN_DEVICE_FUNC inline typename ParametrizedLine<_Scalar, _AmbientDim,_Options>::VectorType
pointAt(const _Scalar & t)194 ParametrizedLine<_Scalar, _AmbientDim,_Options>::pointAt(const _Scalar& t) const
195 {
196   return origin() + (direction()*t);
197 }
198 
199 /** \returns the parameter value of the intersection between \c *this and the given \a hyperplane
200   */
201 template <typename _Scalar, int _AmbientDim, int _Options>
202 template <int OtherOptions>
intersectionParameter(const Hyperplane<_Scalar,_AmbientDim,OtherOptions> & hyperplane)203 EIGEN_DEVICE_FUNC inline _Scalar ParametrizedLine<_Scalar, _AmbientDim,_Options>::intersectionParameter(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const
204 {
205   return -(hyperplane.offset()+hyperplane.normal().dot(origin()))
206           / hyperplane.normal().dot(direction());
207 }
208 
209 
210 /** \deprecated use intersectionParameter()
211   * \returns the parameter value of the intersection between \c *this and the given \a hyperplane
212   */
213 template <typename _Scalar, int _AmbientDim, int _Options>
214 template <int OtherOptions>
intersection(const Hyperplane<_Scalar,_AmbientDim,OtherOptions> & hyperplane)215 EIGEN_DEVICE_FUNC inline _Scalar ParametrizedLine<_Scalar, _AmbientDim,_Options>::intersection(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const
216 {
217   return intersectionParameter(hyperplane);
218 }
219 
220 /** \returns the point of the intersection between \c *this and the given hyperplane
221   */
222 template <typename _Scalar, int _AmbientDim, int _Options>
223 template <int OtherOptions>
224 EIGEN_DEVICE_FUNC inline typename ParametrizedLine<_Scalar, _AmbientDim,_Options>::VectorType
intersectionPoint(const Hyperplane<_Scalar,_AmbientDim,OtherOptions> & hyperplane)225 ParametrizedLine<_Scalar, _AmbientDim,_Options>::intersectionPoint(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const
226 {
227   return pointAt(intersectionParameter(hyperplane));
228 }
229 
230 } // end namespace Eigen
231 
232 #endif // EIGEN_PARAMETRIZEDLINE_H
233