16#ifndef dealii_base_bounding_box_h
17#define dealii_base_bounding_box_h
120template <
int spacedim,
typename Number =
double>
145 template <
class Container>
247 bounds(
const unsigned int direction)
const;
254 vertex(
const unsigned int index)
const;
261 child(
const unsigned int index)
const;
300 template <
class Archive>
313template <
typename Number>
330 template <
class Container>
340template <
int dim,
typename Number =
double>
382 const int coordinate_in_dim)
386 return (locked_coordinate + coordinate_in_dim + 1) % (dim + 1);
396template <
int spacedim,
typename Number>
402 for (
unsigned int i = 0; i < spacedim; ++i)
403 Assert(boundary_points.first[i] <= boundary_points.second[i],
404 ExcMessage(
"Bounding Box can't be created: the points' "
405 "order should be bottom left, top right!"));
407 this->boundary_points = boundary_points;
412template <
int spacedim,
typename Number>
413template <
class Container>
418 if (points.size() > 0)
420 auto &
min = boundary_points.first;
421 auto &
max = boundary_points.second;
422 std::fill(
min.begin_raw(),
424 std::numeric_limits<Number>::infinity());
425 std::fill(
max.begin_raw(),
427 -std::numeric_limits<Number>::infinity());
430 for (
unsigned int d = 0;
d < spacedim; ++
d)
440template <
int spacedim,
typename Number>
444 return this->boundary_points;
449template <
int spacedim,
typename Number>
453 return this->boundary_points;
458template <
int spacedim,
typename Number>
468template <
int spacedim,
typename Number>
478template <
int spacedim,
typename Number>
482 for (
unsigned int d = 0;
d < spacedim; ++
d)
484 boundary_points.first[
d] -= amount;
485 boundary_points.second[
d] += amount;
486 Assert(boundary_points.first[
d] <= boundary_points.second[
d],
487 ExcMessage(
"Bounding Box can't be shrunk this much: the points' "
488 "order should remain bottom left, top right."));
493template <
int spacedim,
typename Number>
494template <
class Archive>
504template <
typename Number>
512template <
typename Number>
521template <
typename Number>
522template <
class Container>
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)
BoundingBox< dim, Number > create_unit_bounding_box()
NeighborType get_neighbor_type(const BoundingBox< spacedim, Number > &other_bbox) const
Point< spacedim, Number > center() const
Number lower_bound(const unsigned int direction) const
void extend(const Number &amount)
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
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
VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &x, const ::VectorizedArray< Number, width > &y)
VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &x, const ::VectorizedArray< Number, width > &y)
#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)
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > epsilon(const Tensor< 2, dim, Number > &Grad_u)
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)