aboutsummaryrefslogtreecommitdiff
path: root/Eigen/src/Geometry/AlignedBox.h
blob: 8e186d57a34ed91738589bc694eb6a8a1c94b0c0 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
//
// This Source Code Form is subject to the terms of the Mozilla
// Public License v. 2.0. If a copy of the MPL was not distributed
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.

#ifndef EIGEN_ALIGNEDBOX_H
#define EIGEN_ALIGNEDBOX_H

namespace Eigen { 

/** \geometry_module \ingroup Geometry_Module
  *
  *
  * \class AlignedBox
  *
  * \brief An axis aligned box
  *
  * \param _Scalar the type of the scalar coefficients
  * \param _AmbientDim the dimension of the ambient space, can be a compile time value or Dynamic.
  *
  * This class represents an axis aligned box as a pair of the minimal and maximal corners.
  */
template <typename _Scalar, int _AmbientDim>
class AlignedBox
{
public:
EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(_Scalar,_AmbientDim)
  enum { AmbientDimAtCompileTime = _AmbientDim };
  typedef _Scalar                                   Scalar;
  typedef NumTraits<Scalar>                         ScalarTraits;
  typedef DenseIndex                                Index;
  typedef typename ScalarTraits::Real               RealScalar;
  typedef typename ScalarTraits::NonInteger      NonInteger;
  typedef Matrix<Scalar,AmbientDimAtCompileTime,1>  VectorType;

  /** Define constants to name the corners of a 1D, 2D or 3D axis aligned bounding box */
  enum CornerType
  {
    /** 1D names */
    Min=0, Max=1,

    /** Added names for 2D */
    BottomLeft=0, BottomRight=1,
    TopLeft=2, TopRight=3,

    /** Added names for 3D */
    BottomLeftFloor=0, BottomRightFloor=1,
    TopLeftFloor=2, TopRightFloor=3,
    BottomLeftCeil=4, BottomRightCeil=5,
    TopLeftCeil=6, TopRightCeil=7
  };


  /** Default constructor initializing a null box. */
  inline AlignedBox()
  { if (AmbientDimAtCompileTime!=Dynamic) setEmpty(); }

  /** Constructs a null box with \a _dim the dimension of the ambient space. */
  inline explicit AlignedBox(Index _dim) : m_min(_dim), m_max(_dim)
  { setEmpty(); }

  /** Constructs a box with extremities \a _min and \a _max. */
  template<typename OtherVectorType1, typename OtherVectorType2>
  inline AlignedBox(const OtherVectorType1& _min, const OtherVectorType2& _max) : m_min(_min), m_max(_max) {}

  /** Constructs a box containing a single point \a p. */
  template<typename Derived>
  inline explicit AlignedBox(const MatrixBase<Derived>& a_p)
  {
    typename internal::nested<Derived,2>::type p(a_p.derived());
    m_min = p;
    m_max = p;
  }

  ~AlignedBox() {}

  /** \returns the dimension in which the box holds */
  inline Index dim() const { return AmbientDimAtCompileTime==Dynamic ? m_min.size() : Index(AmbientDimAtCompileTime); }

  /** \deprecated use isEmpty */
  inline bool isNull() const { return isEmpty(); }

  /** \deprecated use setEmpty */
  inline void setNull() { setEmpty(); }

  /** \returns true if the box is empty. */
  inline bool isEmpty() const { return (m_min.array() > m_max.array()).any(); }

  /** Makes \c *this an empty box. */
  inline void setEmpty()
  {
    m_min.setConstant( ScalarTraits::highest() );
    m_max.setConstant( ScalarTraits::lowest() );
  }

  /** \returns the minimal corner */
  inline const VectorType& (min)() const { return m_min; }
  /** \returns a non const reference to the minimal corner */
  inline VectorType& (min)() { return m_min; }
  /** \returns the maximal corner */
  inline const VectorType& (max)() const { return m_max; }
  /** \returns a non const reference to the maximal corner */
  inline VectorType& (max)() { return m_max; }

  /** \returns the center of the box */
  inline const CwiseUnaryOp<internal::scalar_quotient1_op<Scalar>,
                            const CwiseBinaryOp<internal::scalar_sum_op<Scalar>, const VectorType, const VectorType> >
  center() const
  { return (m_min+m_max)/2; }

  /** \returns the lengths of the sides of the bounding box.
    * Note that this function does not get the same
    * result for integral or floating scalar types: see
    */
  inline const CwiseBinaryOp< internal::scalar_difference_op<Scalar>, const VectorType, const VectorType> sizes() const
  { return m_max - m_min; }

  /** \returns the volume of the bounding box */
  inline Scalar volume() const
  { return sizes().prod(); }

  /** \returns an expression for the bounding box diagonal vector
    * if the length of the diagonal is needed: diagonal().norm()
    * will provide it.
    */
  inline CwiseBinaryOp< internal::scalar_difference_op<Scalar>, const VectorType, const VectorType> diagonal() const
  { return sizes(); }

  /** \returns the vertex of the bounding box at the corner defined by
    * the corner-id corner. It works only for a 1D, 2D or 3D bounding box.
    * For 1D bounding boxes corners are named by 2 enum constants:
    * BottomLeft and BottomRight.
    * For 2D bounding boxes, corners are named by 4 enum constants:
    * BottomLeft, BottomRight, TopLeft, TopRight.
    * For 3D bounding boxes, the following names are added:
    * BottomLeftCeil, BottomRightCeil, TopLeftCeil, TopRightCeil.
    */
  inline VectorType corner(CornerType corner) const
  {
    EIGEN_STATIC_ASSERT(_AmbientDim <= 3, THIS_METHOD_IS_ONLY_FOR_VECTORS_OF_A_SPECIFIC_SIZE);

    VectorType res;

    Index mult = 1;
    for(Index d=0; d<dim(); ++d)
    {
      if( mult & corner ) res[d] = m_max[d];
      else                res[d] = m_min[d];
      mult *= 2;
    }
    return res;
  }

  /** \returns a random point inside the bounding box sampled with
   * a uniform distribution */
  inline VectorType sample() const
  {
    VectorType r;
    for(Index d=0; d<dim(); ++d)
    {
      if(!ScalarTraits::IsInteger)
      {
        r[d] = m_min[d] + (m_max[d]-m_min[d])
             * internal::random<Scalar>(Scalar(0), Scalar(1));
      }
      else
        r[d] = internal::random(m_min[d], m_max[d]);
    }
    return r;
  }

  /** \returns true if the point \a p is inside the box \c *this. */
  template<typename Derived>
  inline bool contains(const MatrixBase<Derived>& a_p) const
  {
    typename internal::nested<Derived,2>::type p(a_p.derived());
    return (m_min.array()<=p.array()).all() && (p.array()<=m_max.array()).all();
  }

  /** \returns true if the box \a b is entirely inside the box \c *this. */
  inline bool contains(const AlignedBox& b) const
  { return (m_min.array()<=(b.min)().array()).all() && ((b.max)().array()<=m_max.array()).all(); }

  /** Extends \c *this such that it contains the point \a p and returns a reference to \c *this. */
  template<typename Derived>
  inline AlignedBox& extend(const MatrixBase<Derived>& a_p)
  {
    typename internal::nested<Derived,2>::type p(a_p.derived());
    m_min = m_min.cwiseMin(p);
    m_max = m_max.cwiseMax(p);
    return *this;
  }

  /** Extends \c *this such that it contains the box \a b and returns a reference to \c *this. */
  inline AlignedBox& extend(const AlignedBox& b)
  {
    m_min = m_min.cwiseMin(b.m_min);
    m_max = m_max.cwiseMax(b.m_max);
    return *this;
  }

  /** Clamps \c *this by the box \a b and returns a reference to \c *this. */
  inline AlignedBox& clamp(const AlignedBox& b)
  {
    m_min = m_min.cwiseMax(b.m_min);
    m_max = m_max.cwiseMin(b.m_max);
    return *this;
  }

  /** Returns an AlignedBox that is the intersection of \a b and \c *this */
  inline AlignedBox intersection(const AlignedBox& b) const
  {return AlignedBox(m_min.cwiseMax(b.m_min), m_max.cwiseMin(b.m_max)); }

  /** Returns an AlignedBox that is the union of \a b and \c *this */
  inline AlignedBox merged(const AlignedBox& b) const
  { return AlignedBox(m_min.cwiseMin(b.m_min), m_max.cwiseMax(b.m_max)); }

  /** Translate \c *this by the vector \a t and returns a reference to \c *this. */
  template<typename Derived>
  inline AlignedBox& translate(const MatrixBase<Derived>& a_t)
  {
    const typename internal::nested<Derived,2>::type t(a_t.derived());
    m_min += t;
    m_max += t;
    return *this;
  }

  /** \returns the squared distance between the point \a p and the box \c *this,
    * and zero if \a p is inside the box.
    * \sa exteriorDistance()
    */
  template<typename Derived>
  inline Scalar squaredExteriorDistance(const MatrixBase<Derived>& a_p) const;

  /** \returns the squared distance between the boxes \a b and \c *this,
    * and zero if the boxes intersect.
    * \sa exteriorDistance()
    */
  inline Scalar squaredExteriorDistance(const AlignedBox& b) const;

  /** \returns the distance between the point \a p and the box \c *this,
    * and zero if \a p is inside the box.
    * \sa squaredExteriorDistance()
    */
  template<typename Derived>
  inline NonInteger exteriorDistance(const MatrixBase<Derived>& p) const
  { using std::sqrt; return sqrt(NonInteger(squaredExteriorDistance(p))); }

  /** \returns the distance between the boxes \a b and \c *this,
    * and zero if the boxes intersect.
    * \sa squaredExteriorDistance()
    */
  inline NonInteger exteriorDistance(const AlignedBox& b) const
  { using std::sqrt; return sqrt(NonInteger(squaredExteriorDistance(b))); }

  /** \returns \c *this with scalar type casted to \a NewScalarType
    *
    * Note that if \a NewScalarType is equal to the current scalar type of \c *this
    * then this function smartly returns a const reference to \c *this.
    */
  template<typename NewScalarType>
  inline typename internal::cast_return_type<AlignedBox,
           AlignedBox<NewScalarType,AmbientDimAtCompileTime> >::type cast() const
  {
    return typename internal::cast_return_type<AlignedBox,
                    AlignedBox<NewScalarType,AmbientDimAtCompileTime> >::type(*this);
  }

  /** Copy constructor with scalar type conversion */
  template<typename OtherScalarType>
  inline explicit AlignedBox(const AlignedBox<OtherScalarType,AmbientDimAtCompileTime>& other)
  {
    m_min = (other.min)().template cast<Scalar>();
    m_max = (other.max)().template cast<Scalar>();
  }

  /** \returns \c true if \c *this is approximately equal to \a other, within the precision
    * determined by \a prec.
    *
    * \sa MatrixBase::isApprox() */
  bool isApprox(const AlignedBox& other, const RealScalar& prec = ScalarTraits::dummy_precision()) const
  { return m_min.isApprox(other.m_min, prec) && m_max.isApprox(other.m_max, prec); }

protected:

  VectorType m_min, m_max;
};



template<typename Scalar,int AmbientDim>
template<typename Derived>
inline Scalar AlignedBox<Scalar,AmbientDim>::squaredExteriorDistance(const MatrixBase<Derived>& a_p) const
{
  typename internal::nested<Derived,2*AmbientDim>::type p(a_p.derived());
  Scalar dist2(0);
  Scalar aux;
  for (Index k=0; k<dim(); ++k)
  {
    if( m_min[k] > p[k] )
    {
      aux = m_min[k] - p[k];
      dist2 += aux*aux;
    }
    else if( p[k] > m_max[k] )
    {
      aux = p[k] - m_max[k];
      dist2 += aux*aux;
    }
  }
  return dist2;
}

template<typename Scalar,int AmbientDim>
inline Scalar AlignedBox<Scalar,AmbientDim>::squaredExteriorDistance(const AlignedBox& b) const
{
  Scalar dist2(0);
  Scalar aux;
  for (Index k=0; k<dim(); ++k)
  {
    if( m_min[k] > b.m_max[k] )
    {
      aux = m_min[k] - b.m_max[k];
      dist2 += aux*aux;
    }
    else if( b.m_min[k] > m_max[k] )
    {
      aux = b.m_min[k] - m_max[k];
      dist2 += aux*aux;
    }
  }
  return dist2;
}

/** \defgroup alignedboxtypedefs Global aligned box typedefs
  *
  * \ingroup Geometry_Module
  *
  * Eigen defines several typedef shortcuts for most common aligned box types.
  *
  * The general patterns are the following:
  *
  * \c AlignedBoxSizeType where \c Size can be \c 1, \c 2,\c 3,\c 4 for fixed size boxes or \c X for dynamic size,
  * and where \c Type can be \c i for integer, \c f for float, \c d for double.
  *
  * For example, \c AlignedBox3d is a fixed-size 3x3 aligned box type of doubles, and \c AlignedBoxXf is a dynamic-size aligned box of floats.
  *
  * \sa class AlignedBox
  */

#define EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, Size, SizeSuffix)    \
/** \ingroup alignedboxtypedefs */                                 \
typedef AlignedBox<Type, Size>   AlignedBox##SizeSuffix##TypeSuffix;

#define EIGEN_MAKE_TYPEDEFS_ALL_SIZES(Type, TypeSuffix) \
EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, 1, 1) \
EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, 2, 2) \
EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, 3, 3) \
EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, 4, 4) \
EIGEN_MAKE_TYPEDEFS(Type, TypeSuffix, Dynamic, X)

EIGEN_MAKE_TYPEDEFS_ALL_SIZES(int,                  i)
EIGEN_MAKE_TYPEDEFS_ALL_SIZES(float,                f)
EIGEN_MAKE_TYPEDEFS_ALL_SIZES(double,               d)

#undef EIGEN_MAKE_TYPEDEFS_ALL_SIZES
#undef EIGEN_MAKE_TYPEDEFS

} // end namespace Eigen

#endif // EIGEN_ALIGNEDBOX_H