$darkmode
Eigen  5.0.1-dev
ParametrizedLine.h
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 
14 // IWYU pragma: private
15 #include "./InternalHeaderCheck.h"
16 
17 namespace Eigen {
18 
32 template <typename Scalar_, int AmbientDim_, int Options_>
33 class ParametrizedLine {
34  public:
35  EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(Scalar_, AmbientDim_)
36  enum { AmbientDimAtCompileTime = AmbientDim_, Options = Options_ };
37  typedef Scalar_ Scalar;
38  typedef typename NumTraits<Scalar>::Real RealScalar;
39  typedef Eigen::Index Index;
41 
43  EIGEN_DEVICE_FUNC inline ParametrizedLine() {}
44 
45  template <int OtherOptions>
47  : m_origin(other.origin()), m_direction(other.direction()) {}
48 
51  EIGEN_DEVICE_FUNC inline explicit ParametrizedLine(Index _dim) : m_origin(_dim), m_direction(_dim) {}
52 
56  EIGEN_DEVICE_FUNC ParametrizedLine(const VectorType& origin, const VectorType& direction)
57  : m_origin(origin), m_direction(direction) {}
58 
59  template <int OtherOptions>
60  EIGEN_DEVICE_FUNC explicit ParametrizedLine(const Hyperplane<Scalar_, AmbientDim_, OtherOptions>& hyperplane);
61 
63  EIGEN_DEVICE_FUNC static inline ParametrizedLine Through(const VectorType& p0, const VectorType& p1) {
64  return ParametrizedLine(p0, (p1 - p0).normalized());
65  }
66 
67  EIGEN_DEVICE_FUNC ~ParametrizedLine() {}
68 
70  EIGEN_DEVICE_FUNC inline Index dim() const { return m_direction.size(); }
71 
72  EIGEN_DEVICE_FUNC const VectorType& origin() const { return m_origin; }
73  EIGEN_DEVICE_FUNC VectorType& origin() { return m_origin; }
74 
75  EIGEN_DEVICE_FUNC const VectorType& direction() const { return m_direction; }
76  EIGEN_DEVICE_FUNC VectorType& direction() { return m_direction; }
77 
81  EIGEN_DEVICE_FUNC RealScalar squaredDistance(const VectorType& p) const {
82  VectorType diff = p - origin();
83  return (diff - direction().dot(diff) * direction()).squaredNorm();
84  }
88  EIGEN_DEVICE_FUNC RealScalar distance(const VectorType& p) const {
89  EIGEN_USING_STD(sqrt) return sqrt(squaredDistance(p));
90  }
91 
93  EIGEN_DEVICE_FUNC VectorType projection(const VectorType& p) const {
94  return origin() + direction().dot(p - origin()) * direction();
95  }
96 
97  EIGEN_DEVICE_FUNC VectorType pointAt(const Scalar& t) const;
98 
99  template <int OtherOptions>
100  EIGEN_DEVICE_FUNC Scalar
101  intersectionParameter(const Hyperplane<Scalar_, AmbientDim_, OtherOptions>& hyperplane) const;
102 
103  template <int OtherOptions>
104  EIGEN_DEVICE_FUNC Scalar intersection(const Hyperplane<Scalar_, AmbientDim_, OtherOptions>& hyperplane) const;
105 
106  template <int OtherOptions>
107  EIGEN_DEVICE_FUNC VectorType
109 
116  template <typename XprType>
117  EIGEN_DEVICE_FUNC inline ParametrizedLine& transform(const MatrixBase<XprType>& mat,
118  TransformTraits traits = Affine) {
119  if (traits == Affine)
120  direction() = (mat * direction()).normalized();
121  else if (traits == Isometry)
122  direction() = mat * direction();
123  else {
124  eigen_assert(0 && "invalid traits value in ParametrizedLine::transform()");
125  }
126  origin() = mat * origin();
127  return *this;
128  }
129 
137  template <int TrOptions>
138  EIGEN_DEVICE_FUNC inline ParametrizedLine& transform(
140  transform(t.linear(), traits);
141  origin() += t.translation();
142  return *this;
143  }
144 
150  template <typename NewScalarType>
151  EIGEN_DEVICE_FUNC inline
152  typename internal::cast_return_type<ParametrizedLine,
154  cast() const {
155  return typename internal::cast_return_type<
157  }
158 
160  template <typename OtherScalarType, int OtherOptions>
161  EIGEN_DEVICE_FUNC inline explicit ParametrizedLine(
163  m_origin = other.origin().template cast<Scalar>();
164  m_direction = other.direction().template cast<Scalar>();
165  }
166 
171  EIGEN_DEVICE_FUNC bool isApprox(const ParametrizedLine& other, const typename NumTraits<Scalar>::Real& prec =
173  return m_origin.isApprox(other.m_origin, prec) && m_direction.isApprox(other.m_direction, prec);
174  }
175 
176  protected:
177  VectorType m_origin, m_direction;
178 };
179 
184 template <typename Scalar_, int AmbientDim_, int Options_>
185 template <int OtherOptions>
188  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 2)
189  direction() = hyperplane.normal().unitOrthogonal();
190  origin() = -hyperplane.normal() * hyperplane.offset();
191 }
192 
195 template <typename Scalar_, int AmbientDim_, int Options_>
196 EIGEN_DEVICE_FUNC inline typename ParametrizedLine<Scalar_, AmbientDim_, Options_>::VectorType
198  return origin() + (direction() * t);
199 }
200 
203 template <typename Scalar_, int AmbientDim_, int Options_>
204 template <int OtherOptions>
206  const Hyperplane<Scalar_, AmbientDim_, OtherOptions>& hyperplane) const {
207  return -(hyperplane.offset() + hyperplane.normal().dot(origin())) / hyperplane.normal().dot(direction());
208 }
209 
213 template <typename Scalar_, int AmbientDim_, int Options_>
214 template <int OtherOptions>
216  const Hyperplane<Scalar_, AmbientDim_, OtherOptions>& hyperplane) const {
217  return intersectionParameter(hyperplane);
218 }
219 
222 template <typename Scalar_, int AmbientDim_, int Options_>
223 template <int OtherOptions>
224 EIGEN_DEVICE_FUNC inline typename ParametrizedLine<Scalar_, AmbientDim_, Options_>::VectorType
226  const Hyperplane<Scalar_, AmbientDim_, OtherOptions>& hyperplane) const {
227  return pointAt(intersectionParameter(hyperplane));
228 }
229 
230 } // end namespace Eigen
231 
232 #endif // EIGEN_PARAMETRIZEDLINE_H
ConstLinearPart linear() const
Definition: Transform.h:374
RealScalar distance(const VectorType &p) const
Definition: ParametrizedLine.h:88
ParametrizedLine(const ParametrizedLine< OtherScalarType, AmbientDimAtCompileTime, OtherOptions > &other)
Definition: ParametrizedLine.h:161
VectorType projection(const VectorType &p) const
Definition: ParametrizedLine.h:93
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
ConstNormalReturnType normal() const
Definition: Hyperplane.h:147
Namespace containing all symbols from the Eigen library.
Definition: B01_Experimental.dox:1
Definition: Constants.h:455
Eigen::Index Index
Definition: ParametrizedLine.h:39
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:232
internal::cast_return_type< ParametrizedLine, ParametrizedLine< NewScalarType, AmbientDimAtCompileTime, Options > >::type cast() const
Definition: ParametrizedLine.h:154
ParametrizedLine()
Definition: ParametrizedLine.h:43
static ParametrizedLine Through(const VectorType &p0, const VectorType &p1)
Definition: ParametrizedLine.h:63
const Scalar & offset() const
Definition: Hyperplane.h:159
VectorType intersectionPoint(const Hyperplane< Scalar_, AmbientDim_, OtherOptions > &hyperplane) const
Definition: ParametrizedLine.h:225
Index dim() const
Definition: ParametrizedLine.h:70
ParametrizedLine(const VectorType &origin, const VectorType &direction)
Definition: ParametrizedLine.h:56
ParametrizedLine & transform(const MatrixBase< XprType > &mat, TransformTraits traits=Affine)
Definition: ParametrizedLine.h:117
A hyperplane.
Definition: ForwardDeclarations.h:473
TransformTraits
Definition: Constants.h:453
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
ConstTranslationPart translation() const
Definition: Transform.h:384
VectorType pointAt(const Scalar &t) const
Definition: ParametrizedLine.h:197
bool isApprox(const ParametrizedLine &other, const typename NumTraits< Scalar >::Real &prec=NumTraits< Scalar >::dummy_precision()) const
Definition: ParametrizedLine.h:171
RealScalar squaredDistance(const VectorType &p) const
Definition: ParametrizedLine.h:81
ParametrizedLine & transform(const Transform< Scalar, AmbientDimAtCompileTime, Affine, TrOptions > &t, TransformTraits traits=Affine)
Definition: ParametrizedLine.h:138
Definition: Constants.h:458
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
A parametrized line.
Definition: ForwardDeclarations.h:471
ParametrizedLine(Index _dim)
Definition: ParametrizedLine.h:51
Represents an homogeneous transformation in a N dimensional space.
Definition: ForwardDeclarations.h:469