Hyperplane.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 // This Source Code Form is subject to the terms of the Mozilla
00008 // Public License v. 2.0. If a copy of the MPL was not distributed
00009 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
00010 
00011 #ifndef EIGEN_HYPERPLANE_H
00012 #define EIGEN_HYPERPLANE_H
00013 
00014 namespace Eigen { 
00015 
00033 template <typename _Scalar, int _AmbientDim, int _Options>
00034 class Hyperplane
00035 {
00036 public:
00037   EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(_Scalar,_AmbientDim==Dynamic ? Dynamic : _AmbientDim+1)
00038   enum {
00039     AmbientDimAtCompileTime = _AmbientDim,
00040     Options = _Options
00041   };
00042   typedef _Scalar Scalar;
00043   typedef typename NumTraits<Scalar>::Real RealScalar;
00044   typedef DenseIndex Index;
00045   typedef Matrix<Scalar,AmbientDimAtCompileTime,1> VectorType;
00046   typedef Matrix<Scalar,Index(AmbientDimAtCompileTime)==Dynamic
00047                         ? Dynamic
00048                         : Index(AmbientDimAtCompileTime)+1,1,Options> Coefficients;
00049   typedef Block<Coefficients,AmbientDimAtCompileTime,1> NormalReturnType;
00050   typedef const Block<const Coefficients,AmbientDimAtCompileTime,1> ConstNormalReturnType;
00051 
00053   inline Hyperplane() {}
00054   
00055   template<int OtherOptions>
00056   Hyperplane(const Hyperplane<Scalar,AmbientDimAtCompileTime,OtherOptions>& other)
00057    : m_coeffs(other.coeffs())
00058   {}
00059 
00062   inline explicit Hyperplane(Index _dim) : m_coeffs(_dim+1) {}
00063 
00067   inline Hyperplane(const VectorType& n, const VectorType& e)
00068     : m_coeffs(n.size()+1)
00069   {
00070     normal() = n;
00071     offset() = -n.dot(e);
00072   }
00073 
00078   inline Hyperplane(const VectorType& n, const Scalar& d)
00079     : m_coeffs(n.size()+1)
00080   {
00081     normal() = n;
00082     offset() = d;
00083   }
00084 
00088   static inline Hyperplane Through(const VectorType& p0, const VectorType& p1)
00089   {
00090     Hyperplane result(p0.size());
00091     result.normal() = (p1 - p0).unitOrthogonal();
00092     result.offset() = -p0.dot(result.normal());
00093     return result;
00094   }
00095 
00099   static inline Hyperplane Through(const VectorType& p0, const VectorType& p1, const VectorType& p2)
00100   {
00101     EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 3)
00102     Hyperplane result(p0.size());
00103     result.normal() = (p2 - p0).cross(p1 - p0).normalized();
00104     result.offset() = -p0.dot(result.normal());
00105     return result;
00106   }
00107 
00112   // FIXME to be consitent with the rest this could be implemented as a static Through function ??
00113   explicit Hyperplane(const ParametrizedLine<Scalar, AmbientDimAtCompileTime>& parametrized)
00114   {
00115     normal() = parametrized.direction().unitOrthogonal();
00116     offset() = -parametrized.origin().dot(normal());
00117   }
00118 
00119   ~Hyperplane() {}
00120 
00122   inline Index dim() const { return AmbientDimAtCompileTime==Dynamic ? m_coeffs.size()-1 : Index(AmbientDimAtCompileTime); }
00123 
00125   void normalize(void)
00126   {
00127     m_coeffs /= normal().norm();
00128   }
00129 
00133   inline Scalar signedDistance(const VectorType& p) const { return normal().dot(p) + offset(); }
00134 
00138   inline Scalar absDistance(const VectorType& p) const { using std::abs; return abs(signedDistance(p)); }
00139 
00142   inline VectorType projection(const VectorType& p) const { return p - signedDistance(p) * normal(); }
00143 
00147   inline ConstNormalReturnType normal() const { return ConstNormalReturnType(m_coeffs,0,0,dim(),1); }
00148 
00152   inline NormalReturnType normal() { return NormalReturnType(m_coeffs,0,0,dim(),1); }
00153 
00157   inline const Scalar& offset() const { return m_coeffs.coeff(dim()); }
00158 
00161   inline Scalar& offset() { return m_coeffs(dim()); }
00162 
00166   inline const Coefficients& coeffs() const { return m_coeffs; }
00167 
00171   inline Coefficients& coeffs() { return m_coeffs; }
00172 
00179   VectorType intersection(const Hyperplane& other) const
00180   {
00181     using std::abs;
00182     EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 2)
00183     Scalar det = coeffs().coeff(0) * other.coeffs().coeff(1) - coeffs().coeff(1) * other.coeffs().coeff(0);
00184     // since the line equations ax+by=c are normalized with a^2+b^2=1, the following tests
00185     // whether the two lines are approximately parallel.
00186     if(internal::isMuchSmallerThan(det, Scalar(1)))
00187     {   // special case where the two lines are approximately parallel. Pick any point on the first line.
00188         if(abs(coeffs().coeff(1))>abs(coeffs().coeff(0)))
00189             return VectorType(coeffs().coeff(1), -coeffs().coeff(2)/coeffs().coeff(1)-coeffs().coeff(0));
00190         else
00191             return VectorType(-coeffs().coeff(2)/coeffs().coeff(0)-coeffs().coeff(1), coeffs().coeff(0));
00192     }
00193     else
00194     {   // general case
00195         Scalar invdet = Scalar(1) / det;
00196         return VectorType(invdet*(coeffs().coeff(1)*other.coeffs().coeff(2)-other.coeffs().coeff(1)*coeffs().coeff(2)),
00197                           invdet*(other.coeffs().coeff(0)*coeffs().coeff(2)-coeffs().coeff(0)*other.coeffs().coeff(2)));
00198     }
00199   }
00200 
00207   template<typename XprType>
00208   inline Hyperplane& transform(const MatrixBase<XprType>& mat, TransformTraits traits = Affine)
00209   {
00210     if (traits==Affine)
00211       normal() = mat.inverse().transpose() * normal();
00212     else if (traits==Isometry)
00213       normal() = mat * normal();
00214     else
00215     {
00216       eigen_assert(0 && "invalid traits value in Hyperplane::transform()");
00217     }
00218     return *this;
00219   }
00220 
00228   template<int TrOptions>
00229   inline Hyperplane& transform(const Transform<Scalar,AmbientDimAtCompileTime,Affine,TrOptions>& t,
00230                                 TransformTraits traits = Affine)
00231   {
00232     transform(t.linear(), traits);
00233     offset() -= normal().dot(t.translation());
00234     return *this;
00235   }
00236 
00242   template<typename NewScalarType>
00243   inline typename internal::cast_return_type<Hyperplane,
00244            Hyperplane<NewScalarType,AmbientDimAtCompileTime,Options> >::type cast() const
00245   {
00246     return typename internal::cast_return_type<Hyperplane,
00247                     Hyperplane<NewScalarType,AmbientDimAtCompileTime,Options> >::type(*this);
00248   }
00249 
00251   template<typename OtherScalarType,int OtherOptions>
00252   inline explicit Hyperplane(const Hyperplane<OtherScalarType,AmbientDimAtCompileTime,OtherOptions>& other)
00253   { m_coeffs = other.coeffs().template cast<Scalar>(); }
00254 
00259   template<int OtherOptions>
00260   bool isApprox(const Hyperplane<Scalar,AmbientDimAtCompileTime,OtherOptions>& other, const typename NumTraits<Scalar>::Real& prec = NumTraits<Scalar>::dummy_precision()) const
00261   { return m_coeffs.isApprox(other.m_coeffs, prec); }
00262 
00263 protected:
00264 
00265   Coefficients m_coeffs;
00266 };
00267 
00268 } // end namespace Eigen
00269 
00270 #endif // EIGEN_HYPERPLANE_H


acado
Author(s): Milan Vukov, Rien Quirynen
autogenerated on Thu Aug 27 2015 11:58:32