ParametrizedLine.h
Go to the documentation of this file.
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 namespace Eigen {
15 
29 template <typename _Scalar, int _AmbientDim, int _Options>
30 class ParametrizedLine
31 {
32 public:
34  enum {
35  AmbientDimAtCompileTime = _AmbientDim,
36  Options = _Options
37  };
38  typedef _Scalar Scalar;
40  typedef Eigen::Index Index;
42 
45 
46  template<int OtherOptions>
48  : m_origin(other.origin()), m_direction(other.direction())
49  {}
50 
53  EIGEN_DEVICE_FUNC inline explicit ParametrizedLine(Index _dim) : m_origin(_dim), m_direction(_dim) {}
54 
58  EIGEN_DEVICE_FUNC ParametrizedLine(const VectorType& origin, const VectorType& direction)
59  : m_origin(origin), m_direction(direction) {}
60 
61  template <int OtherOptions>
63 
65  EIGEN_DEVICE_FUNC static inline ParametrizedLine Through(const VectorType& p0, const VectorType& p1)
66  { return ParametrizedLine(p0, (p1-p0).normalized()); }
67 
69 
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 
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  }
90  EIGEN_DEVICE_FUNC RealScalar distance(const VectorType& p) const { EIGEN_USING_STD(sqrt) return sqrt(squaredDistance(p)); }
91 
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>
100 
101  template <int OtherOptions>
103 
104  template <int OtherOptions>
106 
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 
135  template<int TrOptions>
137  TransformTraits traits = Affine)
138  {
139  transform(t.linear(), traits);
140  origin() += t.translation();
141  return *this;
142  }
143 
149  template<typename NewScalarType>
152  {
155  }
156 
158  template<typename OtherScalarType,int OtherOptions>
160  {
161  m_origin = other.origin().template cast<Scalar>();
162  m_direction = other.direction().template cast<Scalar>();
163  }
164 
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 
181 template <typename _Scalar, int _AmbientDim, int _Options>
182 template <int OtherOptions>
184 {
186  direction() = hyperplane.normal().unitOrthogonal();
187  origin() = -hyperplane.normal()*hyperplane.offset();
188 }
189 
192 template <typename _Scalar, int _AmbientDim, int _Options>
195 {
196  return origin() + (direction()*t);
197 }
198 
201 template <typename _Scalar, int _AmbientDim, int _Options>
202 template <int OtherOptions>
204 {
205  return -(hyperplane.offset()+hyperplane.normal().dot(origin()))
206  / hyperplane.normal().dot(direction());
207 }
208 
209 
213 template <typename _Scalar, int _AmbientDim, int _Options>
214 template <int OtherOptions>
216 {
217  return intersectionParameter(hyperplane);
218 }
219 
222 template <typename _Scalar, int _AmbientDim, int _Options>
223 template <int OtherOptions>
226 {
227  return pointAt(intersectionParameter(hyperplane));
228 }
229 
230 } // end namespace Eigen
231 
232 #endif // EIGEN_PARAMETRIZEDLINE_H
EIGEN_DEVICE_FUNC Index dim() const
EIGEN_DEVICE_FUNC ~ParametrizedLine()
Vector3f p1
EIGEN_DEVICE_FUNC const VectorType & origin() const
EIGEN_DEVICE_FUNC ParametrizedLine(const VectorType &origin, const VectorType &direction)
void hyperplane(const HyperplaneType &_plane)
EIGEN_DEVICE_FUNC ParametrizedLine()
Namespace containing all symbols from the Eigen library.
Definition: jet.h:637
Vector3f p0
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:232
EIGEN_DEVICE_FUNC Scalar intersection(const Hyperplane< _Scalar, _AmbientDim, OtherOptions > &hyperplane) const
EIGEN_DEVICE_FUNC const VectorType & direction() const
EIGEN_DEVICE_FUNC VectorType pointAt(const Scalar &t) const
#define EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(Scalar, Size)
Definition: Memory.h:842
EIGEN_DEVICE_FUNC internal::cast_return_type< ParametrizedLine, ParametrizedLine< NewScalarType, AmbientDimAtCompileTime, Options > >::type cast() const
Scalar EIGEN_BLAS_FUNC() dot(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy)
EIGEN_DEVICE_FUNC ConstTranslationPart translation() const
Definition: Transform.h:404
TransformTraits
Definition: Constants.h:455
EIGEN_DEVICE_FUNC VectorType projection(const VectorType &p) const
Matrix< Scalar, AmbientDimAtCompileTime, 1, Options > VectorType
EIGEN_DEVICE_FUNC ConstLinearPart linear() const
Definition: Transform.h:394
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:74
#define eigen_assert(x)
Definition: Macros.h:1037
EIGEN_DEVICE_FUNC ConstNormalReturnType normal() const
Definition: Hyperplane.h:157
EIGEN_DEVICE_FUNC ParametrizedLine(const ParametrizedLine< OtherScalarType, AmbientDimAtCompileTime, OtherOptions > &other)
EIGEN_DEVICE_FUNC RealScalar squaredDistance(const VectorType &p) const
EIGEN_DEVICE_FUNC bool isApprox(const ParametrizedLine &other, const typename NumTraits< Scalar >::Real &prec=NumTraits< Scalar >::dummy_precision()) const
EIGEN_DEVICE_FUNC const Scalar & offset() const
Definition: Hyperplane.h:167
#define EIGEN_DEVICE_FUNC
Definition: Macros.h:976
EIGEN_DEVICE_FUNC ParametrizedLine & transform(const Transform< Scalar, AmbientDimAtCompileTime, Affine, TrOptions > &t, TransformTraits traits=Affine)
EIGEN_DEVICE_FUNC ParametrizedLine(const ParametrizedLine< Scalar, AmbientDimAtCompileTime, OtherOptions > &other)
EIGEN_DEVICE_FUNC RealScalar distance(const VectorType &p) const
#define EIGEN_USING_STD(FUNC)
Definition: Macros.h:1185
EIGEN_DEVICE_FUNC ParametrizedLine(Index _dim)
float * p
EIGEN_DEVICE_FUNC VectorType & direction()
NumTraits< Scalar >::Real RealScalar
Jet< T, N > sqrt(const Jet< T, N > &f)
Definition: jet.h:418
static EIGEN_DEVICE_FUNC ParametrizedLine Through(const VectorType &p0, const VectorType &p1)
EIGEN_DEVICE_FUNC ParametrizedLine & transform(const MatrixBase< XprType > &mat, TransformTraits traits=Affine)
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
EIGEN_DEVICE_FUNC VectorType intersectionPoint(const Hyperplane< _Scalar, _AmbientDim, OtherOptions > &hyperplane) const
A parametrized line.
EIGEN_DEVICE_FUNC Scalar intersectionParameter(const Hyperplane< _Scalar, _AmbientDim, OtherOptions > &hyperplane) const
Represents an homogeneous transformation in a N dimensional space.
Point2 t(10, 10)
EIGEN_DEVICE_FUNC VectorType & origin()
Definition: pytypes.h:1370
#define EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(TYPE, SIZE)
Definition: StaticAssert.h:157


gtsam
Author(s):
autogenerated on Tue Jul 4 2023 02:35:01