Eigen  3.2.7
 All Classes Namespaces Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
AlignedBox.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 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9 
10 #ifndef EIGEN_ALIGNEDBOX_H
11 #define EIGEN_ALIGNEDBOX_H
12 
13 namespace Eigen {
14 
29 template <typename _Scalar, int _AmbientDim>
31 {
32 public:
33 EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(_Scalar,_AmbientDim)
34  enum { AmbientDimAtCompileTime = _AmbientDim };
35  typedef _Scalar Scalar;
37  typedef DenseIndex Index;
38  typedef typename ScalarTraits::Real RealScalar;
39  typedef typename ScalarTraits::NonInteger NonInteger;
41 
44  {
46  Min=0, Max=1,
58  TopLeftCeil=6, TopRightCeil=7
60  };
61 
62 
64  inline AlignedBox()
65  { if (AmbientDimAtCompileTime!=Dynamic) setEmpty(); }
66 
68  inline explicit AlignedBox(Index _dim) : m_min(_dim), m_max(_dim)
69  { setEmpty(); }
70 
73  template<typename OtherVectorType1, typename OtherVectorType2>
74  inline AlignedBox(const OtherVectorType1& _min, const OtherVectorType2& _max) : m_min(_min), m_max(_max) {}
75 
77  template<typename Derived>
78  inline explicit AlignedBox(const MatrixBase<Derived>& p) : m_min(p), m_max(m_min)
79  { }
80 
81  ~AlignedBox() {}
82 
84  inline Index dim() const { return AmbientDimAtCompileTime==Dynamic ? m_min.size() : Index(AmbientDimAtCompileTime); }
85 
87  inline bool isNull() const { return isEmpty(); }
88 
90  inline void setNull() { setEmpty(); }
91 
94  inline bool isEmpty() const { return (m_min.array() > m_max.array()).any(); }
95 
98  inline void setEmpty()
99  {
100  m_min.setConstant( ScalarTraits::highest() );
101  m_max.setConstant( ScalarTraits::lowest() );
102  }
103 
105  inline const VectorType& (min)() const { return m_min; }
107  inline VectorType& (min)() { return m_min; }
109  inline const VectorType& (max)() const { return m_max; }
111  inline VectorType& (max)() { return m_max; }
112 
115  const CwiseBinaryOp<internal::scalar_sum_op<Scalar>, const VectorType, const VectorType> >
116  center() const
117  { return (m_min+m_max)/2; }
118 
124  { return m_max - m_min; }
125 
127  inline Scalar volume() const
128  { return sizes().prod(); }
129 
135  { return sizes(); }
136 
147  {
148  EIGEN_STATIC_ASSERT(_AmbientDim <= 3, THIS_METHOD_IS_ONLY_FOR_VECTORS_OF_A_SPECIFIC_SIZE);
149 
150  VectorType res;
151 
152  Index mult = 1;
153  for(Index d=0; d<dim(); ++d)
154  {
155  if( mult & corner ) res[d] = m_max[d];
156  else res[d] = m_min[d];
157  mult *= 2;
158  }
159  return res;
160  }
161 
164  inline VectorType sample() const
165  {
166  VectorType r(dim());
167  for(Index d=0; d<dim(); ++d)
168  {
169  if(!ScalarTraits::IsInteger)
170  {
171  r[d] = m_min[d] + (m_max[d]-m_min[d])
172  * internal::random<Scalar>(Scalar(0), Scalar(1));
173  }
174  else
175  r[d] = internal::random(m_min[d], m_max[d]);
176  }
177  return r;
178  }
179 
181  template<typename Derived>
182  inline bool contains(const MatrixBase<Derived>& p) const
183  {
184  typename internal::nested<Derived,2>::type p_n(p.derived());
185  return (m_min.array()<=p_n.array()).all() && (p_n.array()<=m_max.array()).all();
186  }
187 
189  inline bool contains(const AlignedBox& b) const
190  { return (m_min.array()<=(b.min)().array()).all() && ((b.max)().array()<=m_max.array()).all(); }
191 
194  inline bool intersects(const AlignedBox& b) const
195  { return (m_min.array()<=(b.max)().array()).all() && ((b.min)().array()<=m_max.array()).all(); }
196 
199  template<typename Derived>
201  {
202  typename internal::nested<Derived,2>::type p_n(p.derived());
203  m_min = m_min.cwiseMin(p_n);
204  m_max = m_max.cwiseMax(p_n);
205  return *this;
206  }
207 
210  inline AlignedBox& extend(const AlignedBox& b)
211  {
212  m_min = m_min.cwiseMin(b.m_min);
213  m_max = m_max.cwiseMax(b.m_max);
214  return *this;
215  }
216 
220  inline AlignedBox& clamp(const AlignedBox& b)
221  {
222  m_min = m_min.cwiseMax(b.m_min);
223  m_max = m_max.cwiseMin(b.m_max);
224  return *this;
225  }
226 
230  inline AlignedBox intersection(const AlignedBox& b) const
231  {return AlignedBox(m_min.cwiseMax(b.m_min), m_max.cwiseMin(b.m_max)); }
232 
236  inline AlignedBox merged(const AlignedBox& b) const
237  { return AlignedBox(m_min.cwiseMin(b.m_min), m_max.cwiseMax(b.m_max)); }
238 
240  template<typename Derived>
242  {
243  const typename internal::nested<Derived,2>::type t(a_t.derived());
244  m_min += t;
245  m_max += t;
246  return *this;
247  }
248 
253  template<typename Derived>
254  inline Scalar squaredExteriorDistance(const MatrixBase<Derived>& p) const;
255 
260  inline Scalar squaredExteriorDistance(const AlignedBox& b) const;
261 
266  template<typename Derived>
267  inline NonInteger exteriorDistance(const MatrixBase<Derived>& p) const
268  { using std::sqrt; return sqrt(NonInteger(squaredExteriorDistance(p))); }
269 
274  inline NonInteger exteriorDistance(const AlignedBox& b) const
275  { using std::sqrt; return sqrt(NonInteger(squaredExteriorDistance(b))); }
276 
282  template<typename NewScalarType>
283  inline typename internal::cast_return_type<AlignedBox,
285  {
286  return typename internal::cast_return_type<AlignedBox,
288  }
289 
291  template<typename OtherScalarType>
293  {
294  m_min = (other.min)().template cast<Scalar>();
295  m_max = (other.max)().template cast<Scalar>();
296  }
297 
302  bool isApprox(const AlignedBox& other, const RealScalar& prec = ScalarTraits::dummy_precision()) const
303  { return m_min.isApprox(other.m_min, prec) && m_max.isApprox(other.m_max, prec); }
304 
305 protected:
306 
307  VectorType m_min, m_max;
308 };
309 
310 
311 
312 template<typename Scalar,int AmbientDim>
313 template<typename Derived>
315 {
316  typename internal::nested<Derived,2*AmbientDim>::type p(a_p.derived());
317  Scalar dist2(0);
318  Scalar aux;
319  for (Index k=0; k<dim(); ++k)
320  {
321  if( m_min[k] > p[k] )
322  {
323  aux = m_min[k] - p[k];
324  dist2 += aux*aux;
325  }
326  else if( p[k] > m_max[k] )
327  {
328  aux = p[k] - m_max[k];
329  dist2 += aux*aux;
330  }
331  }
332  return dist2;
333 }
334 
335 template<typename Scalar,int AmbientDim>
337 {
338  Scalar dist2(0);
339  Scalar aux;
340  for (Index k=0; k<dim(); ++k)
341  {
342  if( m_min[k] > b.m_max[k] )
343  {
344  aux = m_min[k] - b.m_max[k];
345  dist2 += aux*aux;
346  }
347  else if( b.m_min[k] > m_max[k] )
348  {
349  aux = b.m_min[k] - m_max[k];
350  dist2 += aux*aux;
351  }
352  }
353  return dist2;
354 }
355 
372 #define EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, Size, SizeSuffix) \
373  \
374 typedef AlignedBox<Type, Size> AlignedBox##SizeSuffix##TypeSuffix;
375 
376 #define EIGEN_MAKE_TYPEDEFS_ALL_SIZES(Type, TypeSuffix) \
377 EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, 1, 1) \
378 EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, 2, 2) \
379 EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, 3, 3) \
380 EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, 4, 4) \
381 EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, Dynamic, X)
382 
383 EIGEN_MAKE_TYPEDEFS_ALL_SIZES(int, i)
384 EIGEN_MAKE_TYPEDEFS_ALL_SIZES(float, f)
385 EIGEN_MAKE_TYPEDEFS_ALL_SIZES(double, d)
386 
387 #undef EIGEN_MAKE_TYPEDEFS_ALL_SIZES
388 #undef EIGEN_MAKE_TYPEDEFS
389 
390 } // end namespace Eigen
391 
392 #endif // EIGEN_ALIGNEDBOX_H
Definition: AlignedBox.h:46
AlignedBox(const OtherVectorType1 &_min, const OtherVectorType2 &_max)
Definition: AlignedBox.h:74
Definition: AlignedBox.h:58
Definition: AlignedBox.h:50
CwiseBinaryOp< internal::scalar_difference_op< Scalar >, const VectorType, const VectorType > diagonal() const
Definition: AlignedBox.h:134
CornerType
Definition: AlignedBox.h:43
bool contains(const AlignedBox &b) const
Definition: AlignedBox.h:189
AlignedBox()
Definition: AlignedBox.h:64
bool isNull() const
Definition: AlignedBox.h:87
AlignedBox & clamp(const AlignedBox &b)
Definition: AlignedBox.h:220
AlignedBox(Index _dim)
Definition: AlignedBox.h:68
const VectorType &() min() const
Definition: AlignedBox.h:105
Definition: AlignedBox.h:57
VectorType &() max()
Definition: AlignedBox.h:111
VectorType sample() const
Definition: AlignedBox.h:164
Scalar volume() const
Definition: AlignedBox.h:127
const int Dynamic
Definition: Constants.h:21
AlignedBox & translate(const MatrixBase< Derived > &a_t)
Definition: AlignedBox.h:241
const CwiseBinaryOp< internal::scalar_difference_op< Scalar >, const VectorType, const VectorType > sizes() const
Definition: AlignedBox.h:123
An axis aligned box.
Definition: AlignedBox.h:30
bool isApprox(const AlignedBox &other, const RealScalar &prec=ScalarTraits::dummy_precision()) const
Definition: AlignedBox.h:302
Definition: AlignedBox.h:51
void setEmpty()
Definition: AlignedBox.h:98
Generic expression where a coefficient-wise binary operator is applied to two expressions.
Definition: CwiseBinaryOp.h:107
bool isEmpty() const
Definition: AlignedBox.h:94
Definition: AlignedBox.h:51
AlignedBox(const AlignedBox< OtherScalarType, AmbientDimAtCompileTime > &other)
Definition: AlignedBox.h:292
const VectorType &() max() const
Definition: AlignedBox.h:109
Index dim() const
Definition: AlignedBox.h:84
Definition: AlignedBox.h:50
AlignedBox(const MatrixBase< Derived > &p)
Definition: AlignedBox.h:78
internal::cast_return_type< AlignedBox, AlignedBox< NewScalarType, AmbientDimAtCompileTime > >::type cast() const
Definition: AlignedBox.h:284
Definition: AlignedBox.h:46
const CwiseUnaryOp< internal::scalar_quotient1_op< Scalar >, const CwiseBinaryOp< internal::scalar_sum_op< Scalar >, const VectorType, const VectorType > > center() const
Definition: AlignedBox.h:116
bool contains(const MatrixBase< Derived > &p) const
Definition: AlignedBox.h:182
VectorType &() min()
Definition: AlignedBox.h:107
Scalar squaredExteriorDistance(const MatrixBase< Derived > &p) const
Definition: AlignedBox.h:314
Derived & setConstant(Index size, const Scalar &value)
Definition: CwiseNullaryOp.h:348
bool intersects(const AlignedBox &b) const
Definition: AlignedBox.h:194
AlignedBox & extend(const MatrixBase< Derived > &p)
Definition: AlignedBox.h:200
Definition: AlignedBox.h:57
Definition: AlignedBox.h:56
AlignedBox intersection(const AlignedBox &b) const
Definition: AlignedBox.h:230
VectorType corner(CornerType corner) const
Definition: AlignedBox.h:146
NonInteger exteriorDistance(const MatrixBase< Derived > &p) const
Definition: AlignedBox.h:267
Generic expression where a coefficient-wise unary operator is applied to an expression.
Definition: CwiseUnaryOp.h:59
AlignedBox merged(const AlignedBox &b) const
Definition: AlignedBox.h:236
Definition: AlignedBox.h:56
void setNull()
Definition: AlignedBox.h:90
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
Definition: AlignedBox.h:55
AlignedBox & extend(const AlignedBox &b)
Definition: AlignedBox.h:210
NonInteger exteriorDistance(const AlignedBox &b) const
Definition: AlignedBox.h:274
Definition: AlignedBox.h:55