AlignedBox.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 //
00006 // This Source Code Form is subject to the terms of the Mozilla
00007 // Public License v. 2.0. If a copy of the MPL was not distributed
00008 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
00009 
00010 #ifndef EIGEN_ALIGNEDBOX_H
00011 #define EIGEN_ALIGNEDBOX_H
00012 
00013 namespace Eigen { 
00014 
00029 template <typename _Scalar, int _AmbientDim>
00030 class AlignedBox
00031 {
00032 public:
00033 EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(_Scalar,_AmbientDim)
00034   enum { AmbientDimAtCompileTime = _AmbientDim };
00035   typedef _Scalar                                   Scalar;
00036   typedef NumTraits<Scalar>                         ScalarTraits;
00037   typedef DenseIndex                                Index;
00038   typedef typename ScalarTraits::Real               RealScalar;
00039   typedef typename ScalarTraits::NonInteger      NonInteger;
00040   typedef Matrix<Scalar,AmbientDimAtCompileTime,1>  VectorType;
00041 
00043   enum CornerType
00044   {
00046     Min=0, Max=1,
00050     BottomLeft=0, BottomRight=1,
00051     TopLeft=2, TopRight=3,
00055     BottomLeftFloor=0, BottomRightFloor=1,
00056     TopLeftFloor=2, TopRightFloor=3,
00057     BottomLeftCeil=4, BottomRightCeil=5,
00058     TopLeftCeil=6, TopRightCeil=7
00060   };
00061 
00062 
00064   inline AlignedBox()
00065   { if (AmbientDimAtCompileTime!=Dynamic) setEmpty(); }
00066 
00068   inline explicit AlignedBox(Index _dim) : m_min(_dim), m_max(_dim)
00069   { setEmpty(); }
00070 
00073   template<typename OtherVectorType1, typename OtherVectorType2>
00074   inline AlignedBox(const OtherVectorType1& _min, const OtherVectorType2& _max) : m_min(_min), m_max(_max) {}
00075 
00077   template<typename Derived>
00078   inline explicit AlignedBox(const MatrixBase<Derived>& p) : m_min(p), m_max(m_min)
00079   { }
00080 
00081   ~AlignedBox() {}
00082 
00084   inline Index dim() const { return AmbientDimAtCompileTime==Dynamic ? m_min.size() : Index(AmbientDimAtCompileTime); }
00085 
00087   inline bool isNull() const { return isEmpty(); }
00088 
00090   inline void setNull() { setEmpty(); }
00091 
00094   inline bool isEmpty() const { return (m_min.array() > m_max.array()).any(); }
00095 
00098   inline void setEmpty()
00099   {
00100     m_min.setConstant( ScalarTraits::highest() );
00101     m_max.setConstant( ScalarTraits::lowest() );
00102   }
00103 
00105   inline const VectorType& (min)() const { return m_min; }
00107   inline VectorType& (min)() { return m_min; }
00109   inline const VectorType& (max)() const { return m_max; }
00111   inline VectorType& (max)() { return m_max; }
00112 
00114   inline const CwiseUnaryOp<internal::scalar_quotient1_op<Scalar>,
00115                             const CwiseBinaryOp<internal::scalar_sum_op<Scalar>, const VectorType, const VectorType> >
00116   center() const
00117   { return (m_min+m_max)/2; }
00118 
00123   inline const CwiseBinaryOp< internal::scalar_difference_op<Scalar>, const VectorType, const VectorType> sizes() const
00124   { return m_max - m_min; }
00125 
00127   inline Scalar volume() const
00128   { return sizes().prod(); }
00129 
00134   inline CwiseBinaryOp< internal::scalar_difference_op<Scalar>, const VectorType, const VectorType> diagonal() const
00135   { return sizes(); }
00136 
00146   inline VectorType corner(CornerType corner) const
00147   {
00148     EIGEN_STATIC_ASSERT(_AmbientDim <= 3, THIS_METHOD_IS_ONLY_FOR_VECTORS_OF_A_SPECIFIC_SIZE);
00149 
00150     VectorType res;
00151 
00152     Index mult = 1;
00153     for(Index d=0; d<dim(); ++d)
00154     {
00155       if( mult & corner ) res[d] = m_max[d];
00156       else                res[d] = m_min[d];
00157       mult *= 2;
00158     }
00159     return res;
00160   }
00161 
00164   inline VectorType sample() const
00165   {
00166     VectorType r;
00167     for(Index d=0; d<dim(); ++d)
00168     {
00169       if(!ScalarTraits::IsInteger)
00170       {
00171         r[d] = m_min[d] + (m_max[d]-m_min[d])
00172              * internal::random<Scalar>(Scalar(0), Scalar(1));
00173       }
00174       else
00175         r[d] = internal::random(m_min[d], m_max[d]);
00176     }
00177     return r;
00178   }
00179 
00181   template<typename Derived>
00182   inline bool contains(const MatrixBase<Derived>& p) const
00183   {
00184     typename internal::nested<Derived,2>::type p_n(p.derived());
00185     return (m_min.array()<=p_n.array()).all() && (p_n.array()<=m_max.array()).all();
00186   }
00187 
00189   inline bool contains(const AlignedBox& b) const
00190   { return (m_min.array()<=(b.min)().array()).all() && ((b.max)().array()<=m_max.array()).all(); }
00191 
00194   inline bool intersects(const AlignedBox& b) const
00195   { return (m_min.array()<=(b.max)().array()).all() && ((b.min)().array()<=m_max.array()).all(); }
00196 
00199   template<typename Derived>
00200   inline AlignedBox& extend(const MatrixBase<Derived>& p)
00201   {
00202     typename internal::nested<Derived,2>::type p_n(p.derived());
00203     m_min = m_min.cwiseMin(p_n);
00204     m_max = m_max.cwiseMax(p_n);
00205     return *this;
00206   }
00207 
00210   inline AlignedBox& extend(const AlignedBox& b)
00211   {
00212     m_min = m_min.cwiseMin(b.m_min);
00213     m_max = m_max.cwiseMax(b.m_max);
00214     return *this;
00215   }
00216 
00220   inline AlignedBox& clamp(const AlignedBox& b)
00221   {
00222     m_min = m_min.cwiseMax(b.m_min);
00223     m_max = m_max.cwiseMin(b.m_max);
00224     return *this;
00225   }
00226 
00230   inline AlignedBox intersection(const AlignedBox& b) const
00231   {return AlignedBox(m_min.cwiseMax(b.m_min), m_max.cwiseMin(b.m_max)); }
00232 
00236   inline AlignedBox merged(const AlignedBox& b) const
00237   { return AlignedBox(m_min.cwiseMin(b.m_min), m_max.cwiseMax(b.m_max)); }
00238 
00240   template<typename Derived>
00241   inline AlignedBox& translate(const MatrixBase<Derived>& a_t)
00242   {
00243     const typename internal::nested<Derived,2>::type t(a_t.derived());
00244     m_min += t;
00245     m_max += t;
00246     return *this;
00247   }
00248 
00253   template<typename Derived>
00254   inline Scalar squaredExteriorDistance(const MatrixBase<Derived>& p) const;
00255 
00260   inline Scalar squaredExteriorDistance(const AlignedBox& b) const;
00261 
00266   template<typename Derived>
00267   inline NonInteger exteriorDistance(const MatrixBase<Derived>& p) const
00268   { using std::sqrt; return sqrt(NonInteger(squaredExteriorDistance(p))); }
00269 
00274   inline NonInteger exteriorDistance(const AlignedBox& b) const
00275   { using std::sqrt; return sqrt(NonInteger(squaredExteriorDistance(b))); }
00276 
00282   template<typename NewScalarType>
00283   inline typename internal::cast_return_type<AlignedBox,
00284            AlignedBox<NewScalarType,AmbientDimAtCompileTime> >::type cast() const
00285   {
00286     return typename internal::cast_return_type<AlignedBox,
00287                     AlignedBox<NewScalarType,AmbientDimAtCompileTime> >::type(*this);
00288   }
00289 
00291   template<typename OtherScalarType>
00292   inline explicit AlignedBox(const AlignedBox<OtherScalarType,AmbientDimAtCompileTime>& other)
00293   {
00294     m_min = (other.min)().template cast<Scalar>();
00295     m_max = (other.max)().template cast<Scalar>();
00296   }
00297 
00302   bool isApprox(const AlignedBox& other, const RealScalar& prec = ScalarTraits::dummy_precision()) const
00303   { return m_min.isApprox(other.m_min, prec) && m_max.isApprox(other.m_max, prec); }
00304 
00305 protected:
00306 
00307   VectorType m_min, m_max;
00308 };
00309 
00310 
00311 
00312 template<typename Scalar,int AmbientDim>
00313 template<typename Derived>
00314 inline Scalar AlignedBox<Scalar,AmbientDim>::squaredExteriorDistance(const MatrixBase<Derived>& a_p) const
00315 {
00316   typename internal::nested<Derived,2*AmbientDim>::type p(a_p.derived());
00317   Scalar dist2(0);
00318   Scalar aux;
00319   for (Index k=0; k<dim(); ++k)
00320   {
00321     if( m_min[k] > p[k] )
00322     {
00323       aux = m_min[k] - p[k];
00324       dist2 += aux*aux;
00325     }
00326     else if( p[k] > m_max[k] )
00327     {
00328       aux = p[k] - m_max[k];
00329       dist2 += aux*aux;
00330     }
00331   }
00332   return dist2;
00333 }
00334 
00335 template<typename Scalar,int AmbientDim>
00336 inline Scalar AlignedBox<Scalar,AmbientDim>::squaredExteriorDistance(const AlignedBox& b) const
00337 {
00338   Scalar dist2(0);
00339   Scalar aux;
00340   for (Index k=0; k<dim(); ++k)
00341   {
00342     if( m_min[k] > b.m_max[k] )
00343     {
00344       aux = m_min[k] - b.m_max[k];
00345       dist2 += aux*aux;
00346     }
00347     else if( b.m_min[k] > m_max[k] )
00348     {
00349       aux = b.m_min[k] - m_max[k];
00350       dist2 += aux*aux;
00351     }
00352   }
00353   return dist2;
00354 }
00355 
00372 #define EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, Size, SizeSuffix)    \
00373                                  \
00374 typedef AlignedBox<Type, Size>   AlignedBox##SizeSuffix##TypeSuffix;
00375 
00376 #define EIGEN_MAKE_TYPEDEFS_ALL_SIZES(Type, TypeSuffix) \
00377 EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, 1, 1) \
00378 EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, 2, 2) \
00379 EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, 3, 3) \
00380 EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, 4, 4) \
00381 EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, Dynamic, X)
00382 
00383 EIGEN_MAKE_TYPEDEFS_ALL_SIZES(int,                  i)
00384 EIGEN_MAKE_TYPEDEFS_ALL_SIZES(float,                f)
00385 EIGEN_MAKE_TYPEDEFS_ALL_SIZES(double,               d)
00386 
00387 #undef EIGEN_MAKE_TYPEDEFS_ALL_SIZES
00388 #undef EIGEN_MAKE_TYPEDEFS
00389 
00390 } // end namespace Eigen
00391 
00392 #endif // EIGEN_ALIGNEDBOX_H


turtlebot_exploration_3d
Author(s): Bona , Shawn
autogenerated on Thu Jun 6 2019 20:57:43