16#ifndef dealii_base_bounding_box_h
17#define dealii_base_bounding_box_h
130template <
int spacedim,
typename Number =
double>
155 template <
class Container>
208 const double tolerance = std::numeric_limits<Number>::epsilon())
const;
257 bounds(
const unsigned int direction)
const;
264 vertex(
const unsigned int index)
const;
271 child(
const unsigned int index)
const;
310 template <
class Archive>
323template <
typename Number>
340 template <
class Container>
350template <
int dim,
typename Number =
double>
392 const int coordinate_in_dim)
396 return (locked_coordinate + coordinate_in_dim + 1) % (dim + 1);
406template <
int spacedim,
typename Number>
412 for (
unsigned int i = 0; i < spacedim; ++i)
413 Assert(boundary_points.first[i] <= boundary_points.second[i],
414 ExcMessage(
"Bounding Box can't be created: the points' "
415 "order should be bottom left, top right!"));
417 this->boundary_points = boundary_points;
422template <
int spacedim,
typename Number>
423template <
class Container>
428 if (points.size() > 0)
430 auto &
min = boundary_points.first;
431 auto &
max = boundary_points.second;
432 std::fill(
min.begin_raw(),
434 std::numeric_limits<Number>::infinity());
435 std::fill(
max.begin_raw(),
437 -std::numeric_limits<Number>::infinity());
440 for (
unsigned int d = 0;
d < spacedim; ++
d)
450template <
int spacedim,
typename Number>
454 return this->boundary_points;
459template <
int spacedim,
typename Number>
463 return this->boundary_points;
468template <
int spacedim,
typename Number>
478template <
int spacedim,
typename Number>
488template <
int spacedim,
typename Number>
492 for (
unsigned int d = 0;
d < spacedim; ++
d)
494 boundary_points.first[
d] -= amount;
495 boundary_points.second[
d] += amount;
496 Assert(boundary_points.first[d] <= boundary_points.second[d],
497 ExcMessage(
"Bounding Box can't be shrunk this much: the points' "
498 "order should remain bottom left, top right."));
503template <
int spacedim,
typename Number>
504template <
class Archive>
514template <
typename Number>
522template <
typename Number>
531template <
typename Number>
532template <
class Container>
BoundingBox< dim, Number > create_unit_bounding_box()
BoundingBox(const std::pair< Point< 0, Number >, Point< 0, Number > > &)
BoundingBox(const Container &)
int coordinate_to_one_dim_higher(const int locked_coordinate, const int coordinate_in_dim)
std::pair< Point< spacedim, Number >, Point< spacedim, Number > > boundary_points
BoundingBox< 1, Number > bounds(const unsigned int direction) const
BoundingBox(const std::pair< Point< spacedim, Number >, Point< spacedim, Number > > &boundary_points)
NeighborType get_neighbor_type(const BoundingBox< spacedim, Number > &other_bbox) const
Point< spacedim, Number > center() const
Number lower_bound(const unsigned int direction) const
bool operator==(const BoundingBox< spacedim, Number > &box) const
void serialize(Archive &ar, const unsigned int version)
void merge_with(const BoundingBox< spacedim, Number > &other_bbox)
std::pair< Point< spacedim, Number >, Point< spacedim, Number > > & get_boundary_points()
bool point_inside(const Point< spacedim, Number > &p, const double tolerance=std::numeric_limits< Number >::epsilon()) const
void extend(const Number amount)
BoundingBox(const Container &points)
Point< spacedim, Number > real_to_unit(const Point< spacedim, Number > &point) const
const std::pair< Point< spacedim, Number >, Point< spacedim, Number > > & get_boundary_points() const
Number side_length(const unsigned int direction) const
BoundingBox< spacedim, Number > child(const unsigned int index) const
bool operator!=(const BoundingBox< spacedim, Number > &box) const
Point< spacedim, Number > vertex(const unsigned int index) const
BoundingBox< spacedim - 1, Number > cross_section(const unsigned int direction) const
Number upper_bound(const unsigned int direction) const
Point< spacedim, Number > unit_to_real(const Point< spacedim, Number > &point) const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)