Geometry/AlignedBox.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
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 
27 template <typename _Scalar, int _AmbientDim>
28 class AlignedBox
29 {
30 public:
32  enum { AmbientDimAtCompileTime = _AmbientDim };
33  typedef _Scalar Scalar;
35  typedef DenseIndex Index;
36  typedef typename ScalarTraits::Real RealScalar;
37  typedef typename ScalarTraits::NonInteger NonInteger;
39 
42  {
44  Min=0, Max=1,
45 
49 
55  };
56 
57 
59  inline AlignedBox()
61 
63  inline explicit AlignedBox(Index _dim) : m_min(_dim), m_max(_dim)
64  { setEmpty(); }
65 
67  template<typename OtherVectorType1, typename OtherVectorType2>
68  inline AlignedBox(const OtherVectorType1& _min, const OtherVectorType2& _max) : m_min(_min), m_max(_max) {}
69 
71  template<typename Derived>
72  inline explicit AlignedBox(const MatrixBase<Derived>& a_p)
73  {
74  typename internal::nested<Derived,2>::type p(a_p.derived());
75  m_min = p;
76  m_max = p;
77  }
78 
80 
82  inline Index dim() const { return AmbientDimAtCompileTime==Dynamic ? m_min.size() : Index(AmbientDimAtCompileTime); }
83 
85  inline bool isNull() const { return isEmpty(); }
86 
88  inline void setNull() { setEmpty(); }
89 
91  inline bool isEmpty() const { return (m_min.array() > m_max.array()).any(); }
92 
94  inline void setEmpty()
95  {
96  m_min.setConstant( ScalarTraits::highest() );
97  m_max.setConstant( ScalarTraits::lowest() );
98  }
99 
101  inline const VectorType& (min)() const { return m_min; }
103  inline VectorType& (min)() { return m_min; }
105  inline const VectorType& (max)() const { return m_max; }
107  inline VectorType& (max)() { return m_max; }
108 
111  const CwiseBinaryOp<internal::scalar_sum_op<Scalar>, const VectorType, const VectorType> >
112  center() const
113  { return (m_min+m_max)/2; }
114 
119  inline const CwiseBinaryOp< internal::scalar_difference_op<Scalar>, const VectorType, const VectorType> sizes() const
120  { return m_max - m_min; }
121 
123  inline Scalar volume() const
124  { return sizes().prod(); }
125 
131  { return sizes(); }
132 
142  inline VectorType corner(CornerType corner) const
143  {
144  EIGEN_STATIC_ASSERT(_AmbientDim <= 3, THIS_METHOD_IS_ONLY_FOR_VECTORS_OF_A_SPECIFIC_SIZE);
145 
146  VectorType res;
147 
148  Index mult = 1;
149  for(Index d=0; d<dim(); ++d)
150  {
151  if( mult & corner ) res[d] = m_max[d];
152  else res[d] = m_min[d];
153  mult *= 2;
154  }
155  return res;
156  }
157 
160  inline VectorType sample() const
161  {
162  VectorType r;
163  for(Index d=0; d<dim(); ++d)
164  {
166  {
167  r[d] = m_min[d] + (m_max[d]-m_min[d])
168  * internal::random<Scalar>(Scalar(0), Scalar(1));
169  }
170  else
171  r[d] = internal::random(m_min[d], m_max[d]);
172  }
173  return r;
174  }
175 
177  template<typename Derived>
178  inline bool contains(const MatrixBase<Derived>& a_p) const
179  {
180  typename internal::nested<Derived,2>::type p(a_p.derived());
181  return (m_min.array()<=p.array()).all() && (p.array()<=m_max.array()).all();
182  }
183 
185  inline bool contains(const AlignedBox& b) const
186  { return (m_min.array()<=(b.min)().array()).all() && ((b.max)().array()<=m_max.array()).all(); }
187 
189  template<typename Derived>
191  {
192  typename internal::nested<Derived,2>::type p(a_p.derived());
193  m_min = m_min.cwiseMin(p);
194  m_max = m_max.cwiseMax(p);
195  return *this;
196  }
197 
199  inline AlignedBox& extend(const AlignedBox& b)
200  {
201  m_min = m_min.cwiseMin(b.m_min);
202  m_max = m_max.cwiseMax(b.m_max);
203  return *this;
204  }
205 
207  inline AlignedBox& clamp(const AlignedBox& b)
208  {
209  m_min = m_min.cwiseMax(b.m_min);
210  m_max = m_max.cwiseMin(b.m_max);
211  return *this;
212  }
213 
215  inline AlignedBox intersection(const AlignedBox& b) const
216  {return AlignedBox(m_min.cwiseMax(b.m_min), m_max.cwiseMin(b.m_max)); }
217 
219  inline AlignedBox merged(const AlignedBox& b) const
220  { return AlignedBox(m_min.cwiseMin(b.m_min), m_max.cwiseMax(b.m_max)); }
221 
223  template<typename Derived>
225  {
226  const typename internal::nested<Derived,2>::type t(a_t.derived());
227  m_min += t;
228  m_max += t;
229  return *this;
230  }
231 
236  template<typename Derived>
237  inline Scalar squaredExteriorDistance(const MatrixBase<Derived>& a_p) const;
238 
243  inline Scalar squaredExteriorDistance(const AlignedBox& b) const;
244 
249  template<typename Derived>
250  inline NonInteger exteriorDistance(const MatrixBase<Derived>& p) const
251  { using std::sqrt; return sqrt(NonInteger(squaredExteriorDistance(p))); }
252 
257  inline NonInteger exteriorDistance(const AlignedBox& b) const
258  { using std::sqrt; return sqrt(NonInteger(squaredExteriorDistance(b))); }
259 
265  template<typename NewScalarType>
266  inline typename internal::cast_return_type<AlignedBox,
268  {
269  return typename internal::cast_return_type<AlignedBox,
271  }
272 
274  template<typename OtherScalarType>
276  {
277  m_min = (other.min)().template cast<Scalar>();
278  m_max = (other.max)().template cast<Scalar>();
279  }
280 
285  bool isApprox(const AlignedBox& other, const RealScalar& prec = ScalarTraits::dummy_precision()) const
286  { return m_min.isApprox(other.m_min, prec) && m_max.isApprox(other.m_max, prec); }
287 
288 protected:
289 
290  VectorType m_min, m_max;
291 };
292 
293 
294 
295 template<typename Scalar,int AmbientDim>
296 template<typename Derived>
298 {
299  typename internal::nested<Derived,2*AmbientDim>::type p(a_p.derived());
300  Scalar dist2(0);
301  Scalar aux;
302  for (Index k=0; k<dim(); ++k)
303  {
304  if( m_min[k] > p[k] )
305  {
306  aux = m_min[k] - p[k];
307  dist2 += aux*aux;
308  }
309  else if( p[k] > m_max[k] )
310  {
311  aux = p[k] - m_max[k];
312  dist2 += aux*aux;
313  }
314  }
315  return dist2;
316 }
317 
318 template<typename Scalar,int AmbientDim>
319 inline Scalar AlignedBox<Scalar,AmbientDim>::squaredExteriorDistance(const AlignedBox& b) const
320 {
321  Scalar dist2(0);
322  Scalar aux;
323  for (Index k=0; k<dim(); ++k)
324  {
325  if( m_min[k] > b.m_max[k] )
326  {
327  aux = m_min[k] - b.m_max[k];
328  dist2 += aux*aux;
329  }
330  else if( b.m_min[k] > m_max[k] )
331  {
332  aux = b.m_min[k] - m_max[k];
333  dist2 += aux*aux;
334  }
335  }
336  return dist2;
337 }
338 
355 #define EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, Size, SizeSuffix) \
356  \
357 typedef AlignedBox<Type, Size> AlignedBox##SizeSuffix##TypeSuffix;
358 
359 #define EIGEN_MAKE_TYPEDEFS_ALL_SIZES(Type, TypeSuffix) \
360 EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, 1, 1) \
361 EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, 2, 2) \
362 EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, 3, 3) \
363 EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, 4, 4) \
364 EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, Dynamic, X)
365 
369 
370 #undef EIGEN_MAKE_TYPEDEFS_ALL_SIZES
371 #undef EIGEN_MAKE_TYPEDEFS
372 
373 } // end namespace Eigen
374 
375 #endif // EIGEN_ALIGNEDBOX_H
AlignedBox(const MatrixBase< Derived > &a_p)
ScalarTraits::NonInteger NonInteger
AlignedBox intersection(const AlignedBox &b) const
IntermediateState sqrt(const Expression &arg)
AlignedBox & clamp(const AlignedBox &b)
NonInteger exteriorDistance(const AlignedBox &b) const
ScalarTraits::Real RealScalar
NonInteger exteriorDistance(const MatrixBase< Derived > &p) const
internal::cast_return_type< AlignedBox, AlignedBox< NewScalarType, AmbientDimAtCompileTime > >::type cast() const
iterative scaling algorithm to equilibrate rows and column norms in matrices
Definition: matrix.hpp:471
VectorType sample() const
#define EIGEN_STATIC_ASSERT(CONDITION, MSG)
Definition: StaticAssert.h:111
Scalar squaredExteriorDistance(const VectorType &p) const
AlignedBox merged(const AlignedBox &b) const
const VectorType &() max() const
AlignedBox & extend(const AlignedBox &b)
Generic expression where a coefficient-wise binary operator is applied to two expressions.
const CwiseBinaryOp< internal::scalar_difference_op< Scalar >, const VectorType, const VectorType > sizes() const
const VectorType &() min() const
#define AlignedBox
Definition: All.h:55
VectorType corner(CornerType corner) const
const CwiseUnaryOp< internal::scalar_quotient1_op< Scalar >, const CwiseBinaryOp< internal::scalar_sum_op< Scalar >, const VectorType, const VectorType > > center() const
bool isApprox(const AlignedBox &other, const RealScalar &prec=ScalarTraits::dummy_precision()) const
#define EIGEN_MAKE_TYPEDEFS_ALL_SIZES(Type, TypeSuffix)
Derived & setConstant(Index size, const Scalar &value)
VectorType &() max()
EIGEN_DEFAULT_DENSE_INDEX_TYPE DenseIndex
Definition: XprHelper.h:27
AlignedBox & translate(const MatrixBase< Derived > &a_t)
NumTraits< Scalar > ScalarTraits
VectorType &() min()
CwiseBinaryOp< internal::scalar_difference_op< Scalar >, const VectorType, const VectorType > diagonal() const
Scalar volume() const
AlignedBox(const OtherVectorType1 &_min, const OtherVectorType2 &_max)
Generic expression where a coefficient-wise unary operator is applied to an expression.
Definition: CwiseUnaryOp.h:59
#define EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(Scalar, Size)
AlignedBox & extend(const MatrixBase< Derived > &a_p)
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
bool contains(const MatrixBase< Derived > &a_p) const
Matrix< Scalar, AmbientDimAtCompileTime, 1 > VectorType
AlignedBox(const AlignedBox< OtherScalarType, AmbientDimAtCompileTime > &other)
bool contains(const AlignedBox &b) const


acado
Author(s): Milan Vukov, Rien Quirynen
autogenerated on Mon Jun 10 2019 12:34:27