ParametrizedLine.h
Go to the documentation of this file.
00001 // This file is part of Eigen, a lightweight C++ template library
00002 // for linear algebra.
00003 //
00004 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
00005 // Copyright (C) 2008 Benoit Jacob <jacob.benoit.1@gmail.com>
00006 //
00007 // Eigen is free software; you can redistribute it and/or
00008 // modify it under the terms of the GNU Lesser General Public
00009 // License as published by the Free Software Foundation; either
00010 // version 3 of the License, or (at your option) any later version.
00011 //
00012 // Alternatively, you can redistribute it and/or
00013 // modify it under the terms of the GNU General Public License as
00014 // published by the Free Software Foundation; either version 2 of
00015 // the License, or (at your option) any later version.
00016 //
00017 // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
00018 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
00019 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
00020 // GNU General Public License for more details.
00021 //
00022 // You should have received a copy of the GNU Lesser General Public
00023 // License and a copy of the GNU General Public License along with
00024 // Eigen. If not, see <http://www.gnu.org/licenses/>.
00025 
00026 #ifndef EIGEN_PARAMETRIZEDLINE_H
00027 #define EIGEN_PARAMETRIZEDLINE_H
00028 
00029 namespace Eigen { 
00030 
00044 template <typename _Scalar, int _AmbientDim, int _Options>
00045 class ParametrizedLine
00046 {
00047 public:
00048   EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(_Scalar,_AmbientDim)
00049   enum {
00050     AmbientDimAtCompileTime = _AmbientDim,
00051     Options = _Options
00052   };
00053   typedef _Scalar Scalar;
00054   typedef typename NumTraits<Scalar>::Real RealScalar;
00055   typedef DenseIndex Index;
00056   typedef Matrix<Scalar,AmbientDimAtCompileTime,1,Options> VectorType;
00057 
00059   inline explicit ParametrizedLine() {}
00060   
00061   template<int OtherOptions>
00062   ParametrizedLine(const ParametrizedLine<Scalar,AmbientDimAtCompileTime,OtherOptions>& other)
00063    : m_origin(other.origin()), m_direction(other.direction())
00064   {}
00065 
00068   inline explicit ParametrizedLine(Index _dim) : m_origin(_dim), m_direction(_dim) {}
00069 
00073   ParametrizedLine(const VectorType& origin, const VectorType& direction)
00074     : m_origin(origin), m_direction(direction) {}
00075 
00076   template <int OtherOptions>
00077   explicit ParametrizedLine(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane);
00078 
00080   static inline ParametrizedLine Through(const VectorType& p0, const VectorType& p1)
00081   { return ParametrizedLine(p0, (p1-p0).normalized()); }
00082 
00083   ~ParametrizedLine() {}
00084 
00086   inline Index dim() const { return m_direction.size(); }
00087 
00088   const VectorType& origin() const { return m_origin; }
00089   VectorType& origin() { return m_origin; }
00090 
00091   const VectorType& direction() const { return m_direction; }
00092   VectorType& direction() { return m_direction; }
00093 
00097   RealScalar squaredDistance(const VectorType& p) const
00098   {
00099     VectorType diff = p - origin();
00100     return (diff - direction().dot(diff) * direction()).squaredNorm();
00101   }
00105   RealScalar distance(const VectorType& p) const { return internal::sqrt(squaredDistance(p)); }
00106 
00108   VectorType projection(const VectorType& p) const
00109   { return origin() + direction().dot(p-origin()) * direction(); }
00110 
00111   VectorType pointAt( Scalar t ) const;
00112   
00113   template <int OtherOptions>
00114   Scalar intersectionParameter(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const;
00115  
00116   template <int OtherOptions>
00117   Scalar intersection(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const;
00118   
00119   template <int OtherOptions>
00120   VectorType intersectionPoint(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const;
00121 
00127   template<typename NewScalarType>
00128   inline typename internal::cast_return_type<ParametrizedLine,
00129            ParametrizedLine<NewScalarType,AmbientDimAtCompileTime,Options> >::type cast() const
00130   {
00131     return typename internal::cast_return_type<ParametrizedLine,
00132                     ParametrizedLine<NewScalarType,AmbientDimAtCompileTime,Options> >::type(*this);
00133   }
00134 
00136   template<typename OtherScalarType,int OtherOptions>
00137   inline explicit ParametrizedLine(const ParametrizedLine<OtherScalarType,AmbientDimAtCompileTime,OtherOptions>& other)
00138   {
00139     m_origin = other.origin().template cast<Scalar>();
00140     m_direction = other.direction().template cast<Scalar>();
00141   }
00142 
00147   bool isApprox(const ParametrizedLine& other, typename NumTraits<Scalar>::Real prec = NumTraits<Scalar>::dummy_precision()) const
00148   { return m_origin.isApprox(other.m_origin, prec) && m_direction.isApprox(other.m_direction, prec); }
00149 
00150 protected:
00151 
00152   VectorType m_origin, m_direction;
00153 };
00154 
00159 template <typename _Scalar, int _AmbientDim, int _Options>
00160 template <int OtherOptions>
00161 inline ParametrizedLine<_Scalar, _AmbientDim,_Options>::ParametrizedLine(const Hyperplane<_Scalar, _AmbientDim,OtherOptions>& hyperplane)
00162 {
00163   EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 2)
00164   direction() = hyperplane.normal().unitOrthogonal();
00165   origin() = -hyperplane.normal()*hyperplane.offset();
00166 }
00167 
00170 template <typename _Scalar, int _AmbientDim, int _Options>
00171 inline typename ParametrizedLine<_Scalar, _AmbientDim,_Options>::VectorType
00172 ParametrizedLine<_Scalar, _AmbientDim,_Options>::pointAt( _Scalar t ) const
00173 {
00174   return origin() + (direction()*t); 
00175 }
00176 
00179 template <typename _Scalar, int _AmbientDim, int _Options>
00180 template <int OtherOptions>
00181 inline _Scalar ParametrizedLine<_Scalar, _AmbientDim,_Options>::intersectionParameter(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const
00182 {
00183   return -(hyperplane.offset()+hyperplane.normal().dot(origin()))
00184           / hyperplane.normal().dot(direction());
00185 }
00186 
00187 
00191 template <typename _Scalar, int _AmbientDim, int _Options>
00192 template <int OtherOptions>
00193 inline _Scalar ParametrizedLine<_Scalar, _AmbientDim,_Options>::intersection(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const
00194 {
00195   return intersectionParameter(hyperplane);
00196 }
00197 
00200 template <typename _Scalar, int _AmbientDim, int _Options>
00201 template <int OtherOptions>
00202 inline typename ParametrizedLine<_Scalar, _AmbientDim,_Options>::VectorType
00203 ParametrizedLine<_Scalar, _AmbientDim,_Options>::intersectionPoint(const Hyperplane<_Scalar, _AmbientDim, OtherOptions>& hyperplane) const
00204 {
00205   return pointAt(intersectionParameter(hyperplane));
00206 }
00207 
00208 } // end namespace Eigen
00209 
00210 #endif // EIGEN_PARAMETRIZEDLINE_H