Add BBox_d

This commit is contained in:
Andreas Fabri 2025-01-10 11:18:38 +00:00
parent 0ab8ccf606
commit 2a1ece6969
8 changed files with 164 additions and 38 deletions

View File

@ -25,7 +25,7 @@
#include <CGAL/Origin.h> #include <CGAL/Origin.h>
#include <CGAL/Bbox_2.h> #include <CGAL/Bbox_2.h>
#include <CGAL/Bbox_3.h> #include <CGAL/Bbox_3.h>
#include <CGAL/Bbox.h> #include <CGAL/Bbox_d.h>
#include <CGAL/Default.h> #include <CGAL/Default.h>
#include <CGAL/tss.h> #include <CGAL/tss.h>
#include <CGAL/type_traits/is_iterator.h> #include <CGAL/type_traits/is_iterator.h>
@ -98,9 +98,9 @@ CGAL_LAZY_FORWARD(Bbox_2)
CGAL_LAZY_FORWARD(Bbox_3) CGAL_LAZY_FORWARD(Bbox_3)
#undef CGAL_LAZY_FORWARD #undef CGAL_LAZY_FORWARD
template<class A, class B> Bbox<A,B> const& approx(Bbox<A,B> const& d) { return d; } template<class A> Bbox_d<A> const& approx(Bbox_d<A> const& d) { return d; }
template<class A, class B> Bbox<A,B> const& exact (Bbox<A,B> const& d) { return d; } template<class A> Bbox_d<A> const& exact (Bbox_d<A> const& d) { return d; }
template<class A, class B> int depth (Bbox<A,B> const& ) { return 0; } template<class A> int depth (Bbox_d<A> const& ) { return 0; }
template<class T>inline std::enable_if_t<std::is_arithmetic<T>::value||std::is_enum<T>::value, T> approx(T d){return d;} template<class T>inline std::enable_if_t<std::is_arithmetic<T>::value||std::is_enum<T>::value, T> approx(T d){return d;}
template<class T>inline std::enable_if_t<std::is_arithmetic<T>::value||std::is_enum<T>::value, T> exact (T d){return d;} template<class T>inline std::enable_if_t<std::is_arithmetic<T>::value||std::is_enum<T>::value, T> exact (T d){return d;}

View File

@ -17,7 +17,7 @@
#include <CGAL/license/Frechet_distance.h> #include <CGAL/license/Frechet_distance.h>
#include <CGAL/Bbox.h> #include <CGAL/Bbox_d.h>
#include <array> #include <array>
namespace CGAL { namespace CGAL {
@ -36,9 +36,9 @@ public:
struct Construct_bbox_d struct Construct_bbox_d
{ {
Bbox<Dimension, double> operator()(const Point_d& p) const Bbox_d<Dimension> operator()(const Point_d& p) const
{ {
Bbox<Dimension, double> bb; Bbox_d<Dimension> bb;
for (int i=0;i<dimension; ++i) for (int i=0;i<dimension; ++i)
{ {
(bb.min)(i)=to_interval(p[i]).first; (bb.min)(i)=to_interval(p[i]).first;

View File

@ -22,7 +22,7 @@
#include <CGAL/Search_traits_d.h> #include <CGAL/Search_traits_d.h>
#include <CGAL/Search_traits.h> #include <CGAL/Search_traits.h>
#include <CGAL/number_utils.h> #include <CGAL/number_utils.h>
#include <CGAL/Bbox.h> #include <CGAL/Bbox_d.h>
#include <array> #include <array>
#include <iterator> #include <iterator>
#include <tuple> #include <tuple>
@ -53,7 +53,7 @@ class FrechetKdTree
struct Point_d { struct Point_d {
using BB = Bbox<Dimension_tag<dim>,FT>; using BB = Bbox_d<Dimension_tag<dim>>;
Point ends[2]; Point ends[2];
BB bbox; BB bbox;
using PP = Concatenate_iterator<typename Point::Cartesian_const_iterator, typename Point::Cartesian_const_iterator>; using PP = Concatenate_iterator<typename Point::Cartesian_const_iterator, typename Point::Cartesian_const_iterator>;
@ -232,7 +232,7 @@ auto FrechetKdTree<Traits>::to_kd_tree_point(const Polyline& curve) -> Point_d
res.ends[0] = curve.front(); res.ends[0] = curve.front();
res.ends[1] = curve.back(); res.ends[1] = curve.back();
for (auto const& point : curve) { for (auto const& point : curve) {
Bbox<Dimension_tag<dim>,FT> bb = bbox(point); Bbox_d<Dimension_tag<dim>> bb = bbox(point);
res.bbox += bb; res.bbox += bb;
} }
return res; return res;

View File

@ -18,7 +18,7 @@
#include <CGAL/Frechet_distance/internal/id.h> #include <CGAL/Frechet_distance/internal/id.h>
#include <CGAL/Kernel/Type_mapper.h> #include <CGAL/Kernel/Type_mapper.h>
#include <CGAL/Bbox.h> #include <CGAL/Bbox_d.h>
#include <CGAL/Bbox_2.h> #include <CGAL/Bbox_2.h>
#include <CGAL/Bbox_3.h> #include <CGAL/Bbox_3.h>
#include <CGAL/Cartesian_converter.h> #include <CGAL/Cartesian_converter.h>
@ -46,8 +46,8 @@ double length_of_diagonal(const Bbox_3& bb)
return sqrt(d).sup(); return sqrt(d).sup();
} }
template<int N, typename T> template<int N>
double length_of_diagonal(const Bbox<Dimension_tag<N>,T>& bb) double length_of_diagonal(const Bbox_d<Dimension_tag<N>>& bb)
{ {
using I = Interval_nt<true>; using I = Interval_nt<true>;
I d = square(I((bb.max)(0)) - I((bb.min)(0))); I d = square(I((bb.max)(0)) - I((bb.min)(0)));
@ -227,7 +227,7 @@ public:
static constexpr int dimension = Traits::Dimension::value; static constexpr int dimension = Traits::Dimension::value;
using Bbox = std::conditional_t<dimension==2, using Bbox = std::conditional_t<dimension==2,
Bbox_2, std::conditional_t<dimension==3, Bbox_2, std::conditional_t<dimension==3,
Bbox_3, ::CGAL::Bbox<typename Traits::Dimension, double>>>; Bbox_3, ::CGAL::Bbox_d<typename Traits::Dimension>>>;
using FT = typename Traits::FT; using FT = typename Traits::FT;
using IFT = std::conditional_t<std::is_floating_point_v<FT>, FT, CGAL::Interval_nt<false>>; using IFT = std::conditional_t<std::is_floating_point_v<FT>, FT, CGAL::Interval_nt<false>>;

View File

@ -0,0 +1,125 @@
namespace CGAL {
/*!
\ingroup PkgKernelDKernelObjs
An object `b` of the class `Bbox_d` is a bounding
box in the d-dimensional Euclidean plane \f$ \E^d\f$. This class is templated with a dimension tag.
\cgalModels{Hashable}
\sa `CGAL::Bbox_2`
\sa `CGAL::Bbox_3`
*/
template <typename DimensionTag>
class Bbox_d {
public:
/// \name Creation
/// @{
/*!
introduces an \em empty bounding box with lower left
corner coordinates at \f$ \infty \f$
and with upper right corner coordinates at
\f$ -\infty \f$, \f$ \infty \f$ being
`std::numeric_limits<double>::%infinity()`.
*/
Bbox_d();
/*!
introduces a d-dimensional bounding box from a 2d bounding box.
\pre the dimension must be 2D
*/
Bbox_d(const Bbox_2& b);
/*!
introduces a d-dimensional bounding box from a range.
\pre the range must have the size of the dimension.
\tparam I an iterator model of `InputIterator` with value type double
*/
template <typename I>
Bbox_d(int d, I b, I e);
/*!
introduces a d-dimensional bounding box from a 3d bounding box.
\pre the dimension must be 3D
*/
Bbox_d(const Bbox_3& b);
/// @}
/// \name Operations
/// @{
/*!
ests for equality.
*/
bool operator==(const Bbox_d &c) const;
/*!
tests for inequality.
*/
bool operator!=(const Bbox_d &q) const;
/*!
returns the dimension.
*/
int dimension() const;
/*!
returns an iterator for the %Cartesian coordinates of the lower left and the upper right corner.
*/
Cartesian_const_iterator cartesian_begin() const;
/*!
returns the past-the-end iterator for the %Cartesian coordinates of the lower left and the upper right corner.
*/
Cartesian_const_iterator cartesian_begin() const;
/*!
returns a bounding box of `b` and `c`.
*/
Bbox_d operator+(const Bbox_d &c) const;
/*!
updates `b` to be the bounding box of `b` and `c` and returns itself.
*/
Bbox_d& operator+=(const Bbox_d &c);
/*!
dilates the bounding box by a specified number of ULP.
*/
void dilate(int dist);
/*!
scales the bounding box by `factor`, while keeping its center fixed.
\pre `factor > 0`
*/
void scale(double factor);
/// @}
}; /* end Bbox_d */
/// \ingroup do_overlap_grp
/// @{
/*!
returns `true`, iff `bb1` and `bb2` overlap, i.e., iff their
intersection is non-empty.
\relates Bbox_d
*/
template <typename DimensionTag>
bool do_overlap(const Bbox_d &bb1, const Bbox_d &bb2);
/// @}
} /* end namespace CGAL */

View File

@ -62,6 +62,7 @@
- `CGAL::Sphere_d<Kernel>` - `CGAL::Sphere_d<Kernel>`
- `CGAL::Iso_box_d<Kernel>` - `CGAL::Iso_box_d<Kernel>`
- `CGAL::Aff_transformation_d<Kernel>` - `CGAL::Aff_transformation_d<Kernel>`
- `CGAL::Bbox_d<Dimension>`
\cgalCRPSection{Global %Kernel Functions} \cgalCRPSection{Global %Kernel Functions}
- `CGAL::affinely_independent()` - `CGAL::affinely_independent()`

View File

@ -9,8 +9,8 @@
// //
// Author(s) : Mathieu Brédif // Author(s) : Mathieu Brédif
#ifndef CGAL_BBOX_H #ifndef CGAL_BBOX_D_H
#define CGAL_BBOX_H #define CGAL_BBOX_D_H
#include <boost/config.hpp> // defines BOOST_PREVENT_MACRO_SUBSTITUTION #include <boost/config.hpp> // defines BOOST_PREVENT_MACRO_SUBSTITUTION
#include <stddef.h> #include <stddef.h>
@ -175,26 +175,26 @@ protected:
} }
template <typename Di, typename T> template <typename Di>
class Bbox; class Bbox_d;
// A fixed D-dimensional axis aligned box // A fixed D-dimensional axis aligned box
template<int N, typename T> template<int N>
class Bbox<Dimension_tag<N>, T> : public Impl::Bbox<std::array<T, N>, Bbox<Dimension_tag<N>,T>> class Bbox_d<Dimension_tag<N>> : public Impl::Bbox<std::array<double, N>, Bbox_d<Dimension_tag<N>>>
{ {
enum { D = N }; enum { D = N };
using array_const_iterator = typename std::array<T, N>::const_iterator; using array_const_iterator = typename std::array<double, N>::const_iterator;
public: public:
using Cartesian_const_iterator = Concatenate_iterator<array_const_iterator,array_const_iterator>; using Cartesian_const_iterator = Concatenate_iterator<array_const_iterator,array_const_iterator>;
inline constexpr int dimension() const { return D; } inline constexpr int dimension() const { return D; }
Bbox(int d = 0 ) { CGAL_assertion(d==N || d==0); this->init(d ); } Bbox_d(int d = 0 ) { CGAL_assertion(d==N || d==0); this->init(d ); }
Bbox(int d, const T& range) { CGAL_assertion(d==N || d==0); this->init(d, range); } Bbox_d(int d, double range) { CGAL_assertion(d==N || d==0); this->init(d, range); }
template <typename I> template <typename I>
Bbox(int d, I b, I e) { CGAL_assertion(d==N || d==0); this->init(d, b, e); } Bbox_d(int d, I b, I e) { CGAL_assertion(d==N || d==0); this->init(d, b, e); }
Bbox(const Bbox_2& bb2){ this->init(bb2);} Bbox_d(const Bbox_2& bb2){ this->init(bb2);}
Bbox(const Bbox_3& bb3){ this->init(bb3);} Bbox_d(const Bbox_3& bb3){ this->init(bb3);}
Cartesian_const_iterator cartesian_begin() const Cartesian_const_iterator cartesian_begin() const
{ {
@ -208,15 +208,15 @@ public:
}; };
// A dynamic D-dimensional axis aligned box // A dynamic D-dimensional axis aligned box
template<typename T> template<>
class Bbox<Dynamic_dimension_tag,T> : public Impl::Bbox<std::vector<T>, Bbox<Dynamic_dimension_tag,T>> class Bbox_d<Dynamic_dimension_tag> : public Impl::Bbox<std::vector<double>, Bbox_d<Dynamic_dimension_tag>>
{ {
public: public:
inline int dimension() const { return this->min_values.size(); } inline int dimension() const { return this->min_values.size(); }
Bbox(int d = 0 ) { init_values(d); this->init(d ); } Bbox_d(int d = 0 ) { init_values(d); this->init(d ); }
Bbox(int d, const T& range) { init_values(d); this->init(d, range); } Bbox_d(int d, double range) { init_values(d); this->init(d, range); }
template <typename I> template <typename I>
Bbox(int d, I b, I e) { init_values(d); this->init(d, b, e); } Bbox_d(int d, I b, I e) { init_values(d); this->init(d, b, e); }
protected: protected:
void init_values(int d) { void init_values(int d) {
@ -243,17 +243,17 @@ std::istream& operator>>(std::istream& in, Impl::Bbox<Container, Derived>& bbox)
return in; return in;
} }
template<int N, typename T> template<int N>
Bbox<Dimension_tag<N>, T> operator+(Bbox<Dimension_tag<N>, T> bbox, const Bbox<Dimension_tag<N>, T>& other) Bbox_d<Dimension_tag<N>> operator+(Bbox_d<Dimension_tag<N>> bbox, const Bbox_d<Dimension_tag<N>>& other)
{ {
bbox += other; bbox += other;
return bbox; return bbox;
} }
template<typename Di, typename T> template<typename Di>
inline inline
bool bool
do_overlap(const Bbox<Di,T>& bb1, const Bbox<Di,T>& bb2) do_overlap(const Bbox_d<Di>& bb1, const Bbox_d<Di>& bb2)
{ {
// check for emptiness ?? // check for emptiness ??
int d = bb1.dimension(); int d = bb1.dimension();
@ -268,4 +268,4 @@ do_overlap(const Bbox<Di,T>& bb1, const Bbox<Di,T>& bb2)
} // namespace CGAL } // namespace CGAL
#endif // CGAL_DDT_BBOX_H #endif // CGAL_BBOX_D_H

View File

@ -14,7 +14,7 @@
#include <CGAL/NewKernel_d/utils.h> #include <CGAL/NewKernel_d/utils.h>
#include <CGAL/Dimension.h> #include <CGAL/Dimension.h>
#include <CGAL/Bbox.h> #include <CGAL/Bbox_d.h>
#include <CGAL/Uncertain.h> #include <CGAL/Uncertain.h>
#include <CGAL/NewKernel_d/store_kernel.h> #include <CGAL/NewKernel_d/store_kernel.h>
#include <CGAL/type_traits/is_iterator.h> #include <CGAL/type_traits/is_iterator.h>
@ -1016,7 +1016,7 @@ template<class R_> struct Construct_bbox : private Store_kernel<R_> {
typedef typename Get_type<R, Point_tag>::type Point; typedef typename Get_type<R, Point_tag>::type Point;
typedef typename Get_functor<R, Construct_ttag<Point_cartesian_const_iterator_tag> >::type CI; typedef typename Get_functor<R, Construct_ttag<Point_cartesian_const_iterator_tag> >::type CI;
typedef Bbox<Dimension,double> result_type; typedef Bbox_d<Dimension> result_type;
typedef Point argument_type; typedef Point argument_type;
result_type operator()(Point const&a)const{ result_type operator()(Point const&a)const{
CI ci(this->kernel()); CI ci(this->kernel());