Eigen  3.2.7
 All Classes Namespaces Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
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 <[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 
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 DenseIndex Index;
41  typedef Matrix<Scalar,AmbientDimAtCompileTime,1,Options> VectorType;
42 
44  inline ParametrizedLine() {}
45 
46  template<int OtherOptions>
48  : m_origin(other.origin()), m_direction(other.direction())
49  {}
50 
53  inline explicit ParametrizedLine(Index _dim) : m_origin(_dim), m_direction(_dim) {}
54 
58  ParametrizedLine(const VectorType& origin, const VectorType& direction)
59  : m_origin(origin), m_direction(direction) {}
60 
61  template <int OtherOptions>
63 
65  static inline ParametrizedLine Through(const VectorType& p0, const VectorType& p1)
66  { return ParametrizedLine(p0, (p1-p0).normalized()); }
67 
68  ~ParametrizedLine() {}
69 
71  inline Index dim() const { return m_direction.size(); }
72 
73  const VectorType& origin() const { return m_origin; }
74  VectorType& origin() { return m_origin; }
75 
76  const VectorType& direction() const { return m_direction; }
77  VectorType& direction() { return m_direction; }
78 
82  RealScalar squaredDistance(const VectorType& p) const
83  {
84  VectorType diff = p - origin();
85  return (diff - direction().dot(diff) * direction()).squaredNorm();
86  }
90  RealScalar distance(const VectorType& p) const { using std::sqrt; return sqrt(squaredDistance(p)); }
91 
94  { return origin() + direction().dot(p-origin()) * direction(); }
95 
96  VectorType pointAt(const Scalar& t) const;
97 
98  template <int OtherOptions>
99  Scalar intersectionParameter(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const;
100 
101  template <int OtherOptions>
102  Scalar intersection(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const;
103 
104  template <int OtherOptions>
105  VectorType intersectionPoint(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const;
106 
112  template<typename NewScalarType>
113  inline typename internal::cast_return_type<ParametrizedLine,
115  {
116  return typename internal::cast_return_type<ParametrizedLine,
118  }
119 
121  template<typename OtherScalarType,int OtherOptions>
123  {
124  m_origin = other.origin().template cast<Scalar>();
125  m_direction = other.direction().template cast<Scalar>();
126  }
127 
133  { return m_origin.isApprox(other.m_origin, prec) && m_direction.isApprox(other.m_direction, prec); }
134 
135 protected:
136 
137  VectorType m_origin, m_direction;
138 };
139 
144 template <typename _Scalar, int _AmbientDim, int _Options>
145 template <int OtherOptions>
147 {
148  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 2)
149  direction() = hyperplane.normal().unitOrthogonal();
150  origin() = -hyperplane.normal()*hyperplane.offset();
151 }
152 
155 template <typename _Scalar, int _AmbientDim, int _Options>
158 {
159  return origin() + (direction()*t);
160 }
161 
164 template <typename _Scalar, int _AmbientDim, int _Options>
165 template <int OtherOptions>
167 {
168  return -(hyperplane.offset()+hyperplane.normal().dot(origin()))
169  / hyperplane.normal().dot(direction());
170 }
171 
172 
176 template <typename _Scalar, int _AmbientDim, int _Options>
177 template <int OtherOptions>
179 {
180  return intersectionParameter(hyperplane);
181 }
182 
185 template <typename _Scalar, int _AmbientDim, int _Options>
186 template <int OtherOptions>
189 {
190  return pointAt(intersectionParameter(hyperplane));
191 }
192 
193 } // end namespace Eigen
194 
195 #endif // EIGEN_PARAMETRIZEDLINE_H
bool isApprox(const ParametrizedLine &other, typename NumTraits< Scalar >::Real prec=NumTraits< Scalar >::dummy_precision()) const
Definition: ParametrizedLine.h:132
VectorType projection(const VectorType &p) const
Definition: ParametrizedLine.h:93
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:88
ParametrizedLine(const VectorType &origin, const VectorType &direction)
Definition: ParametrizedLine.h:58
static ParametrizedLine Through(const VectorType &p0, const VectorType &p1)
Definition: ParametrizedLine.h:65
RealScalar squaredDistance(const VectorType &p) const
Definition: ParametrizedLine.h:82
ParametrizedLine()
Definition: ParametrizedLine.h:44
A hyperplane.
Definition: ForwardDeclarations.h:266
ParametrizedLine(Index _dim)
Definition: ParametrizedLine.h:53
RealScalar distance(const VectorType &p) const
Definition: ParametrizedLine.h:90
ParametrizedLine(const ParametrizedLine< OtherScalarType, AmbientDimAtCompileTime, OtherOptions > &other)
Definition: ParametrizedLine.h:122
VectorType pointAt(const Scalar &t) const
Definition: ParametrizedLine.h:157
ConstNormalReturnType normal() const
Definition: Hyperplane.h:157
const Scalar & offset() const
Definition: Hyperplane.h:167
A parametrized line.
Definition: ForwardDeclarations.h:265
Index dim() const
Definition: ParametrizedLine.h:71
internal::cast_return_type< ParametrizedLine, ParametrizedLine< NewScalarType, AmbientDimAtCompileTime, Options > >::type cast() const
Definition: ParametrizedLine.h:114