#include <deal.II/base/point.h>
#include <deal.II/base/utilities.h>
-DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
-#include <boost/geometry/algorithms/envelope.hpp>
-#include <boost/geometry/geometries/multi_point.hpp>
-#if DEAL_II_BOOST_VERSION_GTE(1, 75, 0)
-# include <boost/geometry/strategies/envelope/cartesian.hpp>
-#endif
-DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
-
DEAL_II_NAMESPACE_OPEN
/**
const std::pair<Point<spacedim, Number>, Point<spacedim, Number>> &
get_boundary_points() const;
+ /**
+ * Test for equality.
+ */
+ bool
+ operator==(const BoundingBox<spacedim, Number> &box) const;
+
+ /**
+ * Test for inequality.
+ */
+ bool
+ operator!=(const BoundingBox<spacedim, Number> &box) const;
+
/**
* Check if the current object and @p other_bbox are neighbors, i.e. if the boxes
* have dimension spacedim, check if their intersection is non empty.
template <class Container>
inline BoundingBox<spacedim, Number>::BoundingBox(const Container &points)
{
- boost::geometry::envelope(
- boost::geometry::model::multi_point<Point<spacedim, Number>>(points.begin(),
- points.end()),
- *this);
+ // Use the default constructor in case points is empty instead of setting
+ // things to +oo and -oo
+ if (points.size() > 0)
+ {
+ auto &min = boundary_points.first;
+ auto &max = boundary_points.second;
+ std::fill(min.begin_raw(),
+ min.end_raw(),
+ std::numeric_limits<Number>::infinity());
+ std::fill(max.begin_raw(),
+ max.end_raw(),
+ -std::numeric_limits<Number>::infinity());
+
+ for (const Point<spacedim, Number> &point : points)
+ for (unsigned int d = 0; d < spacedim; ++d)
+ {
+ min[d] = std::min(min[d], point[d]);
+ max[d] = std::max(max[d], point[d]);
+ }
+ }
+}
+
+
+
+template <int spacedim, typename Number>
+inline bool
+BoundingBox<spacedim, Number>::
+operator==(const BoundingBox<spacedim, Number> &box) const
+{
+ return boundary_points == box.boundary_points;
+}
+
+
+
+template <int spacedim, typename Number>
+inline bool
+BoundingBox<spacedim, Number>::
+operator!=(const BoundingBox<spacedim, Number> &box) const
+{
+ return boundary_points != box.boundary_points;
}
<< " is inside: " << b.point_inside(test_points[i]) << std::endl;
deallog << std::endl;
+
+ // Verify that we get the same box when we use all those points together
+ {
+ test_points.push_back(boundaries.first);
+ test_points.push_back(boundaries.second);
+
+ BoundingBox<spacedim> b2(test_points);
+ deallog << "Boxes should be equal : " << (b2 == b) << std::endl;
+ }
+
test_points.clear();
// To create outside points we take a non-convex combination
for (unsigned int i = 0; i < test_points.size(); ++i)
deallog << test_points[i]
<< " is inside: " << b.point_inside(test_points[i]) << std::endl;
+
+ // Similarly, verify that we get different boxes since some points are
+ // outside:
+ {
+ BoundingBox<spacedim> b2(test_points);
+ deallog << "Boxes should not be equal : " << (b2 != b) << std::endl;
+ }
deallog << std::endl;
}
DEAL::0.320000 is inside: 1
DEAL::0.200000 is inside: 1
DEAL::
+DEAL::Boxes should be equal : 1
DEAL::Points outside:
DEAL::-2.50000 is inside: 0
DEAL::5.00000 is inside: 0
DEAL::-7.50000 is inside: 0
DEAL::10.0000 is inside: 0
DEAL::-12.5000 is inside: 0
+DEAL::Boxes should not be equal : 1
DEAL::
DEAL::
DEAL::Test for dimension 2
DEAL::0.320000 0.320000 is inside: 1
DEAL::0.200000 0.00000 is inside: 1
DEAL::
+DEAL::Boxes should be equal : 1
DEAL::Points outside:
DEAL::-2.50000 -4.00000 is inside: 0
DEAL::5.00000 8.00000 is inside: 0
DEAL::-7.50000 -12.0000 is inside: 0
DEAL::10.0000 16.0000 is inside: 0
DEAL::-12.5000 -20.0000 is inside: 0
+DEAL::Boxes should not be equal : 1
DEAL::
DEAL::
DEAL::Test for dimension 3
DEAL::0.320000 0.320000 0.320000 is inside: 1
DEAL::0.200000 0.00000 -0.200000 is inside: 1
DEAL::
+DEAL::Boxes should be equal : 1
DEAL::Points outside:
DEAL::-2.50000 -4.00000 -5.50000 is inside: 0
DEAL::5.00000 8.00000 11.0000 is inside: 0
DEAL::-7.50000 -12.0000 -16.5000 is inside: 0
DEAL::10.0000 16.0000 22.0000 is inside: 0
DEAL::-12.5000 -20.0000 -27.5000 is inside: 0
+DEAL::Boxes should not be equal : 1
DEAL::
DEAL::
DEAL::Test for dimension 3, unitary box