xref: /aosp_15_r20/external/eigen/Eigen/src/Geometry/ParametrizedLine.h (revision bf2c37156dfe67e5dfebd6d394bad8b2ab5804d4)
1*bf2c3715SXin Li // This file is part of Eigen, a lightweight C++ template library
2*bf2c3715SXin Li // for linear algebra.
3*bf2c3715SXin Li //
4*bf2c3715SXin Li // Copyright (C) 2008 Gael Guennebaud <[email protected]>
5*bf2c3715SXin Li // Copyright (C) 2008 Benoit Jacob <[email protected]>
6*bf2c3715SXin Li //
7*bf2c3715SXin Li // This Source Code Form is subject to the terms of the Mozilla
8*bf2c3715SXin Li // Public License v. 2.0. If a copy of the MPL was not distributed
9*bf2c3715SXin Li // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10*bf2c3715SXin Li 
11*bf2c3715SXin Li #ifndef EIGEN_PARAMETRIZEDLINE_H
12*bf2c3715SXin Li #define EIGEN_PARAMETRIZEDLINE_H
13*bf2c3715SXin Li 
14*bf2c3715SXin Li namespace Eigen {
15*bf2c3715SXin Li 
16*bf2c3715SXin Li /** \geometry_module \ingroup Geometry_Module
17*bf2c3715SXin Li   *
18*bf2c3715SXin Li   * \class ParametrizedLine
19*bf2c3715SXin Li   *
20*bf2c3715SXin Li   * \brief A parametrized line
21*bf2c3715SXin Li   *
22*bf2c3715SXin Li   * A parametrized line is defined by an origin point \f$ \mathbf{o} \f$ and a unit
23*bf2c3715SXin Li   * direction vector \f$ \mathbf{d} \f$ such that the line corresponds to
24*bf2c3715SXin Li   * the set \f$ l(t) = \mathbf{o} + t \mathbf{d} \f$, \f$ t \in \mathbf{R} \f$.
25*bf2c3715SXin Li   *
26*bf2c3715SXin Li   * \tparam _Scalar the scalar type, i.e., the type of the coefficients
27*bf2c3715SXin Li   * \tparam _AmbientDim the dimension of the ambient space, can be a compile time value or Dynamic.
28*bf2c3715SXin Li   */
29*bf2c3715SXin Li template <typename _Scalar, int _AmbientDim, int _Options>
30*bf2c3715SXin Li class ParametrizedLine
31*bf2c3715SXin Li {
32*bf2c3715SXin Li public:
33*bf2c3715SXin Li   EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(_Scalar,_AmbientDim)
34*bf2c3715SXin Li   enum {
35*bf2c3715SXin Li     AmbientDimAtCompileTime = _AmbientDim,
36*bf2c3715SXin Li     Options = _Options
37*bf2c3715SXin Li   };
38*bf2c3715SXin Li   typedef _Scalar Scalar;
39*bf2c3715SXin Li   typedef typename NumTraits<Scalar>::Real RealScalar;
40*bf2c3715SXin Li   typedef Eigen::Index Index; ///< \deprecated since Eigen 3.3
41*bf2c3715SXin Li   typedef Matrix<Scalar,AmbientDimAtCompileTime,1,Options> VectorType;
42*bf2c3715SXin Li 
43*bf2c3715SXin Li   /** Default constructor without initialization */
ParametrizedLine()44*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline ParametrizedLine() {}
45*bf2c3715SXin Li 
46*bf2c3715SXin Li   template<int OtherOptions>
ParametrizedLine(const ParametrizedLine<Scalar,AmbientDimAtCompileTime,OtherOptions> & other)47*bf2c3715SXin Li   EIGEN_DEVICE_FUNC ParametrizedLine(const ParametrizedLine<Scalar,AmbientDimAtCompileTime,OtherOptions>& other)
48*bf2c3715SXin Li    : m_origin(other.origin()), m_direction(other.direction())
49*bf2c3715SXin Li   {}
50*bf2c3715SXin Li 
51*bf2c3715SXin Li   /** Constructs a dynamic-size line with \a _dim the dimension
52*bf2c3715SXin Li     * of the ambient space */
ParametrizedLine(Index _dim)53*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline explicit ParametrizedLine(Index _dim) : m_origin(_dim), m_direction(_dim) {}
54*bf2c3715SXin Li 
55*bf2c3715SXin Li   /** Initializes a parametrized line of direction \a direction and origin \a origin.
56*bf2c3715SXin Li     * \warning the vector direction is assumed to be normalized.
57*bf2c3715SXin Li     */
ParametrizedLine(const VectorType & origin,const VectorType & direction)58*bf2c3715SXin Li   EIGEN_DEVICE_FUNC ParametrizedLine(const VectorType& origin, const VectorType& direction)
59*bf2c3715SXin Li     : m_origin(origin), m_direction(direction) {}
60*bf2c3715SXin Li 
61*bf2c3715SXin Li   template <int OtherOptions>
62*bf2c3715SXin Li   EIGEN_DEVICE_FUNC explicit ParametrizedLine(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane);
63*bf2c3715SXin Li 
64*bf2c3715SXin Li   /** Constructs a parametrized line going from \a p0 to \a p1. */
Through(const VectorType & p0,const VectorType & p1)65*bf2c3715SXin Li   EIGEN_DEVICE_FUNC static inline ParametrizedLine Through(const VectorType& p0, const VectorType& p1)
66*bf2c3715SXin Li   { return ParametrizedLine(p0, (p1-p0).normalized()); }
67*bf2c3715SXin Li 
~ParametrizedLine()68*bf2c3715SXin Li   EIGEN_DEVICE_FUNC ~ParametrizedLine() {}
69*bf2c3715SXin Li 
70*bf2c3715SXin Li   /** \returns the dimension in which the line holds */
dim()71*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline Index dim() const { return m_direction.size(); }
72*bf2c3715SXin Li 
origin()73*bf2c3715SXin Li   EIGEN_DEVICE_FUNC const VectorType& origin() const { return m_origin; }
origin()74*bf2c3715SXin Li   EIGEN_DEVICE_FUNC VectorType& origin() { return m_origin; }
75*bf2c3715SXin Li 
direction()76*bf2c3715SXin Li   EIGEN_DEVICE_FUNC const VectorType& direction() const { return m_direction; }
direction()77*bf2c3715SXin Li   EIGEN_DEVICE_FUNC VectorType& direction() { return m_direction; }
78*bf2c3715SXin Li 
79*bf2c3715SXin Li   /** \returns the squared distance of a point \a p to its projection onto the line \c *this.
80*bf2c3715SXin Li     * \sa distance()
81*bf2c3715SXin Li     */
squaredDistance(const VectorType & p)82*bf2c3715SXin Li   EIGEN_DEVICE_FUNC RealScalar squaredDistance(const VectorType& p) const
83*bf2c3715SXin Li   {
84*bf2c3715SXin Li     VectorType diff = p - origin();
85*bf2c3715SXin Li     return (diff - direction().dot(diff) * direction()).squaredNorm();
86*bf2c3715SXin Li   }
87*bf2c3715SXin Li   /** \returns the distance of a point \a p to its projection onto the line \c *this.
88*bf2c3715SXin Li     * \sa squaredDistance()
89*bf2c3715SXin Li     */
distance(const VectorType & p)90*bf2c3715SXin Li   EIGEN_DEVICE_FUNC RealScalar distance(const VectorType& p) const { EIGEN_USING_STD(sqrt) return sqrt(squaredDistance(p)); }
91*bf2c3715SXin Li 
92*bf2c3715SXin Li   /** \returns the projection of a point \a p onto the line \c *this. */
projection(const VectorType & p)93*bf2c3715SXin Li   EIGEN_DEVICE_FUNC VectorType projection(const VectorType& p) const
94*bf2c3715SXin Li   { return origin() + direction().dot(p-origin()) * direction(); }
95*bf2c3715SXin Li 
96*bf2c3715SXin Li   EIGEN_DEVICE_FUNC VectorType pointAt(const Scalar& t) const;
97*bf2c3715SXin Li 
98*bf2c3715SXin Li   template <int OtherOptions>
99*bf2c3715SXin Li   EIGEN_DEVICE_FUNC Scalar intersectionParameter(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const;
100*bf2c3715SXin Li 
101*bf2c3715SXin Li   template <int OtherOptions>
102*bf2c3715SXin Li   EIGEN_DEVICE_FUNC Scalar intersection(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const;
103*bf2c3715SXin Li 
104*bf2c3715SXin Li   template <int OtherOptions>
105*bf2c3715SXin Li   EIGEN_DEVICE_FUNC VectorType intersectionPoint(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const;
106*bf2c3715SXin Li 
107*bf2c3715SXin Li   /** Applies the transformation matrix \a mat to \c *this and returns a reference to \c *this.
108*bf2c3715SXin Li     *
109*bf2c3715SXin Li     * \param mat the Dim x Dim transformation matrix
110*bf2c3715SXin Li     * \param traits specifies whether the matrix \a mat represents an #Isometry
111*bf2c3715SXin Li     *               or a more generic #Affine transformation. The default is #Affine.
112*bf2c3715SXin Li     */
113*bf2c3715SXin Li   template<typename XprType>
114*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline ParametrizedLine& transform(const MatrixBase<XprType>& mat, TransformTraits traits = Affine)
115*bf2c3715SXin Li   {
116*bf2c3715SXin Li     if (traits==Affine)
117*bf2c3715SXin Li       direction() = (mat * direction()).normalized();
118*bf2c3715SXin Li     else if (traits==Isometry)
119*bf2c3715SXin Li       direction() = mat * direction();
120*bf2c3715SXin Li     else
121*bf2c3715SXin Li     {
122*bf2c3715SXin Li       eigen_assert(0 && "invalid traits value in ParametrizedLine::transform()");
123*bf2c3715SXin Li     }
124*bf2c3715SXin Li     origin() = mat * origin();
125*bf2c3715SXin Li     return *this;
126*bf2c3715SXin Li   }
127*bf2c3715SXin Li 
128*bf2c3715SXin Li   /** Applies the transformation \a t to \c *this and returns a reference to \c *this.
129*bf2c3715SXin Li     *
130*bf2c3715SXin Li     * \param t the transformation of dimension Dim
131*bf2c3715SXin Li     * \param traits specifies whether the transformation \a t represents an #Isometry
132*bf2c3715SXin Li     *               or a more generic #Affine transformation. The default is #Affine.
133*bf2c3715SXin Li     *               Other kind of transformations are not supported.
134*bf2c3715SXin Li     */
135*bf2c3715SXin Li   template<int TrOptions>
136*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline ParametrizedLine& transform(const Transform<Scalar,AmbientDimAtCompileTime,Affine,TrOptions>& t,
137*bf2c3715SXin Li                                                        TransformTraits traits = Affine)
138*bf2c3715SXin Li   {
139*bf2c3715SXin Li     transform(t.linear(), traits);
140*bf2c3715SXin Li     origin() += t.translation();
141*bf2c3715SXin Li     return *this;
142*bf2c3715SXin Li   }
143*bf2c3715SXin Li 
144*bf2c3715SXin Li /** \returns \c *this with scalar type casted to \a NewScalarType
145*bf2c3715SXin Li     *
146*bf2c3715SXin Li     * Note that if \a NewScalarType is equal to the current scalar type of \c *this
147*bf2c3715SXin Li     * then this function smartly returns a const reference to \c *this.
148*bf2c3715SXin Li     */
149*bf2c3715SXin Li   template<typename NewScalarType>
150*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline typename internal::cast_return_type<ParametrizedLine,
cast()151*bf2c3715SXin Li            ParametrizedLine<NewScalarType,AmbientDimAtCompileTime,Options> >::type cast() const
152*bf2c3715SXin Li   {
153*bf2c3715SXin Li     return typename internal::cast_return_type<ParametrizedLine,
154*bf2c3715SXin Li                     ParametrizedLine<NewScalarType,AmbientDimAtCompileTime,Options> >::type(*this);
155*bf2c3715SXin Li   }
156*bf2c3715SXin Li 
157*bf2c3715SXin Li   /** Copy constructor with scalar type conversion */
158*bf2c3715SXin Li   template<typename OtherScalarType,int OtherOptions>
ParametrizedLine(const ParametrizedLine<OtherScalarType,AmbientDimAtCompileTime,OtherOptions> & other)159*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline explicit ParametrizedLine(const ParametrizedLine<OtherScalarType,AmbientDimAtCompileTime,OtherOptions>& other)
160*bf2c3715SXin Li   {
161*bf2c3715SXin Li     m_origin = other.origin().template cast<Scalar>();
162*bf2c3715SXin Li     m_direction = other.direction().template cast<Scalar>();
163*bf2c3715SXin Li   }
164*bf2c3715SXin Li 
165*bf2c3715SXin Li   /** \returns \c true if \c *this is approximately equal to \a other, within the precision
166*bf2c3715SXin Li     * determined by \a prec.
167*bf2c3715SXin Li     *
168*bf2c3715SXin Li     * \sa MatrixBase::isApprox() */
169*bf2c3715SXin Li   EIGEN_DEVICE_FUNC bool isApprox(const ParametrizedLine& other, const typename NumTraits<Scalar>::Real& prec = NumTraits<Scalar>::dummy_precision()) const
170*bf2c3715SXin Li   { return m_origin.isApprox(other.m_origin, prec) && m_direction.isApprox(other.m_direction, prec); }
171*bf2c3715SXin Li 
172*bf2c3715SXin Li protected:
173*bf2c3715SXin Li 
174*bf2c3715SXin Li   VectorType m_origin, m_direction;
175*bf2c3715SXin Li };
176*bf2c3715SXin Li 
177*bf2c3715SXin Li /** Constructs a parametrized line from a 2D hyperplane
178*bf2c3715SXin Li   *
179*bf2c3715SXin Li   * \warning the ambient space must have dimension 2 such that the hyperplane actually describes a line
180*bf2c3715SXin Li   */
181*bf2c3715SXin Li template <typename _Scalar, int _AmbientDim, int _Options>
182*bf2c3715SXin Li template <int OtherOptions>
ParametrizedLine(const Hyperplane<_Scalar,_AmbientDim,OtherOptions> & hyperplane)183*bf2c3715SXin Li EIGEN_DEVICE_FUNC inline ParametrizedLine<_Scalar, _AmbientDim,_Options>::ParametrizedLine(const Hyperplane<_Scalar, _AmbientDim,OtherOptions>& hyperplane)
184*bf2c3715SXin Li {
185*bf2c3715SXin Li   EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 2)
186*bf2c3715SXin Li   direction() = hyperplane.normal().unitOrthogonal();
187*bf2c3715SXin Li   origin() = -hyperplane.normal()*hyperplane.offset();
188*bf2c3715SXin Li }
189*bf2c3715SXin Li 
190*bf2c3715SXin Li /** \returns the point at \a t along this line
191*bf2c3715SXin Li   */
192*bf2c3715SXin Li template <typename _Scalar, int _AmbientDim, int _Options>
193*bf2c3715SXin Li EIGEN_DEVICE_FUNC inline typename ParametrizedLine<_Scalar, _AmbientDim,_Options>::VectorType
pointAt(const _Scalar & t)194*bf2c3715SXin Li ParametrizedLine<_Scalar, _AmbientDim,_Options>::pointAt(const _Scalar& t) const
195*bf2c3715SXin Li {
196*bf2c3715SXin Li   return origin() + (direction()*t);
197*bf2c3715SXin Li }
198*bf2c3715SXin Li 
199*bf2c3715SXin Li /** \returns the parameter value of the intersection between \c *this and the given \a hyperplane
200*bf2c3715SXin Li   */
201*bf2c3715SXin Li template <typename _Scalar, int _AmbientDim, int _Options>
202*bf2c3715SXin Li template <int OtherOptions>
intersectionParameter(const Hyperplane<_Scalar,_AmbientDim,OtherOptions> & hyperplane)203*bf2c3715SXin Li EIGEN_DEVICE_FUNC inline _Scalar ParametrizedLine<_Scalar, _AmbientDim,_Options>::intersectionParameter(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const
204*bf2c3715SXin Li {
205*bf2c3715SXin Li   return -(hyperplane.offset()+hyperplane.normal().dot(origin()))
206*bf2c3715SXin Li           / hyperplane.normal().dot(direction());
207*bf2c3715SXin Li }
208*bf2c3715SXin Li 
209*bf2c3715SXin Li 
210*bf2c3715SXin Li /** \deprecated use intersectionParameter()
211*bf2c3715SXin Li   * \returns the parameter value of the intersection between \c *this and the given \a hyperplane
212*bf2c3715SXin Li   */
213*bf2c3715SXin Li template <typename _Scalar, int _AmbientDim, int _Options>
214*bf2c3715SXin Li template <int OtherOptions>
intersection(const Hyperplane<_Scalar,_AmbientDim,OtherOptions> & hyperplane)215*bf2c3715SXin Li EIGEN_DEVICE_FUNC inline _Scalar ParametrizedLine<_Scalar, _AmbientDim,_Options>::intersection(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const
216*bf2c3715SXin Li {
217*bf2c3715SXin Li   return intersectionParameter(hyperplane);
218*bf2c3715SXin Li }
219*bf2c3715SXin Li 
220*bf2c3715SXin Li /** \returns the point of the intersection between \c *this and the given hyperplane
221*bf2c3715SXin Li   */
222*bf2c3715SXin Li template <typename _Scalar, int _AmbientDim, int _Options>
223*bf2c3715SXin Li template <int OtherOptions>
224*bf2c3715SXin Li EIGEN_DEVICE_FUNC inline typename ParametrizedLine<_Scalar, _AmbientDim,_Options>::VectorType
intersectionPoint(const Hyperplane<_Scalar,_AmbientDim,OtherOptions> & hyperplane)225*bf2c3715SXin Li ParametrizedLine<_Scalar, _AmbientDim,_Options>::intersectionPoint(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const
226*bf2c3715SXin Li {
227*bf2c3715SXin Li   return pointAt(intersectionParameter(hyperplane));
228*bf2c3715SXin Li }
229*bf2c3715SXin Li 
230*bf2c3715SXin Li } // end namespace Eigen
231*bf2c3715SXin Li 
232*bf2c3715SXin Li #endif // EIGEN_PARAMETRIZEDLINE_H
233