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 #include "./InternalHeaderCheck.h"
15 
16 namespace Eigen {
17 
31 template <typename Scalar_, int AmbientDim_, int Options_>
33 {
34 public:
36  enum {
37  AmbientDimAtCompileTime = AmbientDim_,
38  Options = Options_
39  };
40  typedef Scalar_ Scalar;
42  typedef Eigen::Index Index;
44 
47 
48  template<int OtherOptions>
50  : m_origin(other.origin()), m_direction(other.direction())
51  {}
52 
55  EIGEN_DEVICE_FUNC inline explicit ParametrizedLine(Index _dim) : m_origin(_dim), m_direction(_dim) {}
56 
62 
63  template <int OtherOptions>
65 
68  { return ParametrizedLine(p0, (p1-p0).normalized()); }
69 
71 
73  EIGEN_DEVICE_FUNC inline Index dim() const { return m_direction.size(); }
74 
75  EIGEN_DEVICE_FUNC const VectorType& origin() const { return m_origin; }
77 
80 
85  {
86  VectorType diff = p - origin();
87  return (diff - direction().dot(diff) * direction()).squaredNorm();
88  }
93 
96  { return origin() + direction().dot(p-origin()) * direction(); }
97 
99 
100  template <int OtherOptions>
102 
103  template <int OtherOptions>
105 
106  template <int OtherOptions>
108 
115  template<typename XprType>
117  {
118  if (traits==Affine)
119  direction() = (mat * direction()).normalized();
120  else if (traits==Isometry)
121  direction() = mat * direction();
122  else
123  {
124  eigen_assert(0 && "invalid traits value in ParametrizedLine::transform()");
125  }
126  origin() = mat * origin();
127  return *this;
128  }
129 
137  template<int TrOptions>
139  TransformTraits traits = Affine)
140  {
141  transform(t.linear(), traits);
142  origin() += t.translation();
143  return *this;
144  }
145 
151  template<typename NewScalarType>
152  EIGEN_DEVICE_FUNC inline typename internal::cast_return_type<ParametrizedLine,
154  {
155  return typename internal::cast_return_type<ParametrizedLine,
157  }
158 
160  template<typename OtherScalarType,int OtherOptions>
162  {
163  m_origin = other.origin().template cast<Scalar>();
164  m_direction = other.direction().template cast<Scalar>();
165  }
166 
172  { return m_origin.isApprox(other.m_origin, prec) && m_direction.isApprox(other.m_direction, prec); }
173 
174 protected:
175 
177 };
178 
183 template <typename Scalar_, int AmbientDim_, int Options_>
184 template <int OtherOptions>
186 {
188  direction() = hyperplane.normal().unitOrthogonal();
189  origin() = -hyperplane.normal()*hyperplane.offset();
190 }
191 
194 template <typename Scalar_, int AmbientDim_, int Options_>
197 {
198  return origin() + (direction()*t);
199 }
200 
203 template <typename Scalar_, int AmbientDim_, int Options_>
204 template <int OtherOptions>
206 {
207  return -(hyperplane.offset()+hyperplane.normal().dot(origin()))
208  / hyperplane.normal().dot(direction());
209 }
210 
211 
215 template <typename Scalar_, int AmbientDim_, int Options_>
216 template <int OtherOptions>
218 {
219  return intersectionParameter(hyperplane);
220 }
221 
224 template <typename Scalar_, int AmbientDim_, int Options_>
225 template <int OtherOptions>
228 {
229  return pointAt(intersectionParameter(hyperplane));
230 }
231 
232 } // end namespace Eigen
233 
234 #endif // EIGEN_PARAMETRIZEDLINE_H
#define EIGEN_USING_STD(FUNC)
Definition: Macros.h:1080
#define EIGEN_DEVICE_FUNC
Definition: Macros.h:883
#define eigen_assert(x)
Definition: Macros.h:902
Vector3f p0
Vector3f p1
#define EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(Scalar, Size)
Definition: Memory.h:921
#define EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(TYPE, SIZE)
Definition: StaticAssert.h:51
float * p
A hyperplane.
Definition: Hyperplane.h:37
ConstNormalReturnType normal() const
Definition: Hyperplane.h:159
const Scalar & offset() const
Definition: Hyperplane.h:169
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
A parametrized line.
ParametrizedLine & transform(const Transform< Scalar, AmbientDimAtCompileTime, Affine, TrOptions > &t, TransformTraits traits=Affine)
Scalar intersection(const Hyperplane< Scalar_, AmbientDim_, OtherOptions > &hyperplane) const
ParametrizedLine & transform(const MatrixBase< XprType > &mat, TransformTraits traits=Affine)
RealScalar squaredDistance(const VectorType &p) const
const VectorType & direction() const
RealScalar distance(const VectorType &p) const
const VectorType & origin() const
Matrix< Scalar, AmbientDimAtCompileTime, 1, Options > VectorType
VectorType intersectionPoint(const Hyperplane< Scalar_, AmbientDim_, OtherOptions > &hyperplane) const
ParametrizedLine(const ParametrizedLine< Scalar, AmbientDimAtCompileTime, OtherOptions > &other)
Scalar intersectionParameter(const Hyperplane< Scalar_, AmbientDim_, OtherOptions > &hyperplane) const
ParametrizedLine(const VectorType &origin, const VectorType &direction)
bool isApprox(const ParametrizedLine &other, const typename NumTraits< Scalar >::Real &prec=NumTraits< Scalar >::dummy_precision()) const
internal::cast_return_type< ParametrizedLine, ParametrizedLine< NewScalarType, AmbientDimAtCompileTime, Options > >::type cast() const
VectorType pointAt(const Scalar &t) const
static ParametrizedLine Through(const VectorType &p0, const VectorType &p1)
NumTraits< Scalar >::Real RealScalar
ParametrizedLine(const ParametrizedLine< OtherScalarType, AmbientDimAtCompileTime, OtherOptions > &other)
VectorType projection(const VectorType &p) const
Represents an homogeneous transformation in a N dimensional space.
Definition: Transform.h:207
ConstLinearPart linear() const
Definition: Transform.h:398
ConstTranslationPart translation() const
Definition: Transform.h:408
TransformTraits
Definition: Constants.h:459
@ Affine
Definition: Constants.h:464
@ Isometry
Definition: Constants.h:461
: InteropHeaders
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:231