15#ifndef dealii_base_bounding_box_h
16#define dealii_base_bounding_box_h
134template <
int spacedim,
typename Number =
double>
175 template <
class Container>
209 const double tolerance = std::numeric_limits<Number>::epsilon())
const;
217 const double tolerance = std::numeric_limits<Number>::epsilon())
const;
236 const double tolerance = std::numeric_limits<Number>::epsilon())
const;
308 bounds(
const unsigned int direction)
const;
315 vertex(
const unsigned int index)
const;
322 child(
const unsigned int index)
const;
364 const unsigned int direction)
const;
379 template <
class Archive>
392template <
typename Number>
409 template <
class Container>
419template <
int dim,
typename Number =
double>
461 const int coordinate_in_dim)
465 return (locked_coordinate + coordinate_in_dim + 1) % (dim + 1);
475template <
int spacedim,
typename Number>
483template <
int spacedim,
typename Number>
489 for (
unsigned int i = 0; i < spacedim; ++i)
490 Assert(boundary_points.first[i] <= boundary_points.second[i],
491 ExcMessage(
"Bounding Box can't be created: the points' "
492 "order should be bottom left, top right!"));
494 this->boundary_points = boundary_points;
499template <
int spacedim,
typename Number>
500template <
class Container>
505 if (points.size() > 0)
507 auto &
min = boundary_points.first;
508 auto &
max = boundary_points.second;
509 for (
unsigned int d = 0;
d < spacedim; ++
d)
511 min[
d] = std::numeric_limits<Number>::infinity();
512 max[
d] = -std::numeric_limits<Number>::infinity();
516 for (unsigned
int d = 0;
d < spacedim; ++
d)
526template <
int spacedim,
typename Number>
530 return this->boundary_points;
535template <
int spacedim,
typename Number>
539 return this->boundary_points;
544template <
int spacedim,
typename Number>
554template <
int spacedim,
typename Number>
564template <
int spacedim,
typename Number>
568 for (
unsigned int d = 0;
d < spacedim; ++
d)
570 boundary_points.first[
d] -= amount;
571 boundary_points.second[
d] += amount;
572 Assert(boundary_points.first[d] <= boundary_points.second[d],
573 ExcMessage(
"Bounding Box can't be shrunk this much: the points' "
574 "order should remain bottom left, top right."));
580template <
int spacedim,
typename Number>
593template <
int spacedim,
typename Number>
596 const Number relative_amount)
const
601 for (
unsigned int d = 0;
d < spacedim; ++
d)
604 bb.boundary_points.second[
d] += relative_amount * side_length(d);
605 Assert(bb.boundary_points.first[d] <= bb.boundary_points.second[d],
606 ExcMessage(
"Bounding Box can't be shrunk this much: the points' "
607 "order should remain bottom left, top right."));
615template <
int spacedim,
typename Number>
616template <
class Archive>
626template <
typename Number>
634template <
typename Number>
643template <
typename Number>
644template <
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
bool has_overlap_with(const BoundingBox< spacedim, Number > &other_bbox, const double tolerance=std::numeric_limits< Number >::epsilon()) const
BoundingBox(const std::pair< Point< spacedim, Number >, Point< spacedim, Number > > &boundary_points)
Point< spacedim, Number > center() const
Number lower_bound(const unsigned int direction) const
Number signed_distance(const Point< spacedim, Number > &point, 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< spacedim, Number > create_extended_relative(const Number relative_amount) const
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
BoundingBox< spacedim, Number > & operator=(const BoundingBox< spacedim, Number > &t)=default
NeighborType get_neighbor_type(const BoundingBox< spacedim, Number > &other_bbox, const double tolerance=std::numeric_limits< Number >::epsilon()) const
BoundingBox< spacedim, Number > create_extended(const Number amount) const
Point< spacedim, Number > vertex(const unsigned int index) const
BoundingBox(const Point< spacedim, Number > &point)
BoundingBox< spacedim - 1, Number > cross_section(const unsigned int direction) const
BoundingBox(const BoundingBox< spacedim, Number > &box)=default
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 > &)