21template <
int spacedim,
typename Number>
24 const double tolerance)
const
26 for (
unsigned int i = 0; i < spacedim; ++i)
31 if ((p[i] < this->boundary_points.first[i] -
32 tolerance *
std::abs(this->boundary_points.second[i] -
33 this->boundary_points.first[i])) ||
34 (p[i] > this->boundary_points.second[i] +
35 tolerance *
std::abs(this->boundary_points.second[i] -
36 this->boundary_points.first[i])))
44template <
int spacedim,
typename Number>
49 for (
unsigned int i = 0; i < spacedim; ++i)
51 this->boundary_points.first[i] =
52 std::min(this->boundary_points.first[i],
54 this->boundary_points.second[i] =
55 std::max(this->boundary_points.second[i],
62template <
int spacedim,
typename Number>
78 const std::array<Point<spacedim, Number>, 2> bbox1 = {
79 {this->get_boundary_points().first,
80 this->get_boundary_points().second}};
81 const std::array<Point<spacedim, Number>, 2> bbox2 = {
86 for (
unsigned int d = 0; d < spacedim; ++d)
87 if (bbox1[0][d] * (1 - std::numeric_limits<Number>::epsilon()) >
89 bbox2[0][d] * (1 - std::numeric_limits<Number>::epsilon()) >
95 std::array<double, spacedim> intersect_bbox_min;
96 std::array<double, spacedim> intersect_bbox_max;
97 for (
unsigned int d = 0; d < spacedim; ++d)
99 intersect_bbox_min[d] =
std::max(bbox1[0][d], bbox2[0][d]);
100 intersect_bbox_max[d] =
std::min(bbox1[1][d], bbox2[1][d]);
104 int intersect_dim = spacedim;
105 for (
unsigned int d = 0; d < spacedim; ++d)
106 if (
std::abs(intersect_bbox_min[d] - intersect_bbox_max[d]) <=
107 std::numeric_limits<Number>::epsilon() *
112 if (intersect_dim == 0 || intersect_dim == spacedim - 2)
117 unsigned int not_align_1 = 0, not_align_2 = 0;
118 bool same_direction =
true;
119 for (
unsigned int d = 0; d < spacedim; ++d)
121 if (
std::abs(bbox2[0][d] - bbox1[0][d]) >
122 std::numeric_limits<double>::epsilon() *
125 if (
std::abs(bbox1[1][d] - bbox2[1][d]) >
126 std::numeric_limits<double>::epsilon() *
129 if (not_align_1 != not_align_2)
131 same_direction =
false;
136 if (not_align_1 <= 1 && not_align_2 <= 1 && same_direction)
140 if ((this->point_inside(bbox2[0]) && this->point_inside(bbox2[1])) ||
152template <
int spacedim,
typename Number>
157 for (
unsigned int i = 0; i < spacedim; ++i)
158 vol *= (this->boundary_points.second[i] - this->boundary_points.first[i]);
164template <
int spacedim,
typename Number>
170 return boundary_points.first[direction];
175template <
int spacedim,
typename Number>
181 return boundary_points.second[direction];
186template <
int spacedim,
typename Number>
191 for (
unsigned int i = 0; i < spacedim; ++i)
192 point[i] = .5 * (boundary_points.first[i] + boundary_points.second[i]);
199template <
int spacedim,
typename Number>
206 lower_upper_bounds.first[0] = lower_bound(direction);
207 lower_upper_bounds.second[0] = upper_bound(direction);
214template <
int spacedim,
typename Number>
220 return boundary_points.second[direction] - boundary_points.first[direction];
225template <
int spacedim,
typename Number>
235 for (
unsigned int i = 0; i < spacedim; ++i)
236 point[i] = boundary_points.first[i] + side_length(i) * unit_cell_vertex[i];
243template <
int spacedim,
typename Number>
261 child_lower_upper_corner;
262 for (
unsigned int i = 0; i < spacedim; ++i)
264 const double child_side_length = side_length(i) / 2;
266 const double child_center = (parent_center[i] + parent_vertex[i]) / 2;
268 child_lower_upper_corner.first[i] =
269 child_center + child_side_length * (lower_corner_unit_cell[i] - .5);
270 child_lower_upper_corner.second[i] =
271 child_center + child_side_length * (upper_corner_unit_cell[i] - .5);
279template <
int spacedim,
typename Number>
285 std::pair<
Point<spacedim - 1, Number>,
Point<spacedim - 1, Number>>
286 cross_section_lower_upper_corner;
287 for (
unsigned int d = 0; d < spacedim - 1; ++d)
289 const int index_to_write_from =
290 internal::coordinate_to_one_dim_higher<spacedim - 1>(direction, d);
292 cross_section_lower_upper_corner.first[d] =
293 boundary_points.first[index_to_write_from];
295 cross_section_lower_upper_corner.second[d] =
296 boundary_points.second[index_to_write_from];
299 return BoundingBox<spacedim - 1, Number>(cross_section_lower_upper_corner);
304template <
int spacedim,
typename Number>
310 const auto diag = boundary_points.second - boundary_points.first;
311 unit -= boundary_points.first;
312 for (
unsigned int d = 0; d < spacedim; ++d)
319template <
int spacedim,
typename Number>
324 auto real = boundary_points.first;
325 const auto diag = boundary_points.second - boundary_points.first;
326 for (
unsigned int d = 0; d < spacedim; ++d)
327 real[d] += diag[d] * point[d];
333template <
int dim,
typename Number>
338 for (
unsigned int i = 0; i < dim; ++i)
340 lower_upper_corner.second[i] = 1;
346#include "bounding_box.inst"
BoundingBox< dim, Number > create_unit_bounding_box()
std::pair< Point< spacedim, Number >, Point< spacedim, Number > > boundary_points
BoundingBox< 1, Number > bounds(const unsigned int direction) const
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 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
Point< spacedim, Number > real_to_unit(const Point< spacedim, Number > &point) const
Number side_length(const unsigned int direction) const
BoundingBox< spacedim, Number > child(const unsigned int index) 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 AssertIndexRange(index, range)
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
static Point< dim > unit_cell_vertex(const unsigned int vertex)