Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00042 template <typename _Scalar, int _AmbientDim>
00043 class ParametrizedLine
00044 {
00045 public:
00046 EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(_Scalar,_AmbientDim)
00047 enum { AmbientDimAtCompileTime = _AmbientDim };
00048 typedef _Scalar Scalar;
00049 typedef typename NumTraits<Scalar>::Real RealScalar;
00050 typedef Matrix<Scalar,AmbientDimAtCompileTime,1> VectorType;
00051
00053 inline explicit ParametrizedLine() {}
00054
00057 inline explicit ParametrizedLine(int _dim) : m_origin(_dim), m_direction(_dim) {}
00058
00062 ParametrizedLine(const VectorType& origin, const VectorType& direction)
00063 : m_origin(origin), m_direction(direction) {}
00064
00065 explicit ParametrizedLine(const Hyperplane<_Scalar, _AmbientDim>& hyperplane);
00066
00068 static inline ParametrizedLine Through(const VectorType& p0, const VectorType& p1)
00069 { return ParametrizedLine(p0, (p1-p0).normalized()); }
00070
00071 ~ParametrizedLine() {}
00072
00074 inline int dim() const { return m_direction.size(); }
00075
00076 const VectorType& origin() const { return m_origin; }
00077 VectorType& origin() { return m_origin; }
00078
00079 const VectorType& direction() const { return m_direction; }
00080 VectorType& direction() { return m_direction; }
00081
00085 RealScalar squaredDistance(const VectorType& p) const
00086 {
00087 VectorType diff = p-origin();
00088 return (diff - diff.eigen2_dot(direction())* direction()).squaredNorm();
00089 }
00093 RealScalar distance(const VectorType& p) const { return ei_sqrt(squaredDistance(p)); }
00094
00096 VectorType projection(const VectorType& p) const
00097 { return origin() + (p-origin()).eigen2_dot(direction()) * direction(); }
00098
00099 Scalar intersection(const Hyperplane<_Scalar, _AmbientDim>& hyperplane);
00100
00106 template<typename NewScalarType>
00107 inline typename internal::cast_return_type<ParametrizedLine,
00108 ParametrizedLine<NewScalarType,AmbientDimAtCompileTime> >::type cast() const
00109 {
00110 return typename internal::cast_return_type<ParametrizedLine,
00111 ParametrizedLine<NewScalarType,AmbientDimAtCompileTime> >::type(*this);
00112 }
00113
00115 template<typename OtherScalarType>
00116 inline explicit ParametrizedLine(const ParametrizedLine<OtherScalarType,AmbientDimAtCompileTime>& other)
00117 {
00118 m_origin = other.origin().template cast<Scalar>();
00119 m_direction = other.direction().template cast<Scalar>();
00120 }
00121
00126 bool isApprox(const ParametrizedLine& other, typename NumTraits<Scalar>::Real prec = precision<Scalar>()) const
00127 { return m_origin.isApprox(other.m_origin, prec) && m_direction.isApprox(other.m_direction, prec); }
00128
00129 protected:
00130
00131 VectorType m_origin, m_direction;
00132 };
00133
00138 template <typename _Scalar, int _AmbientDim>
00139 inline ParametrizedLine<_Scalar, _AmbientDim>::ParametrizedLine(const Hyperplane<_Scalar, _AmbientDim>& hyperplane)
00140 {
00141 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 2)
00142 direction() = hyperplane.normal().unitOrthogonal();
00143 origin() = -hyperplane.normal()*hyperplane.offset();
00144 }
00145
00148 template <typename _Scalar, int _AmbientDim>
00149 inline _Scalar ParametrizedLine<_Scalar, _AmbientDim>::intersection(const Hyperplane<_Scalar, _AmbientDim>& hyperplane)
00150 {
00151 return -(hyperplane.offset()+origin().eigen2_dot(hyperplane.normal()))
00152 /(direction().eigen2_dot(hyperplane.normal()));
00153 }