Reference documentation for deal.II version 9.4.1
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
bounding_box.cc
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2017 - 2021 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
18
20
21template <int spacedim, typename Number>
22bool
24 const double tolerance) const
25{
26 for (unsigned int i = 0; i < spacedim; ++i)
27 {
28 // Bottom left-top right convention: the point is outside if it's smaller
29 // than the first or bigger than the second boundary point The bounding
30 // box is defined as a closed set
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])))
37 return false;
38 }
39 return true;
40}
41
42
43
44template <int spacedim, typename Number>
45void
47 const BoundingBox<spacedim, Number> &other_bbox)
48{
49 for (unsigned int i = 0; i < spacedim; ++i)
50 {
51 this->boundary_points.first[i] =
52 std::min(this->boundary_points.first[i],
53 other_bbox.boundary_points.first[i]);
54 this->boundary_points.second[i] =
55 std::max(this->boundary_points.second[i],
56 other_bbox.boundary_points.second[i]);
57 }
58}
59
60
61
62template <int spacedim, typename Number>
65 const BoundingBox<spacedim, Number> &other_bbox) const
66{
67 if (spacedim == 1)
68 {
69 // In dimension 1 if the two bounding box are neighbors
70 // we can merge them
71 if (this->point_inside(other_bbox.boundary_points.first) ||
72 this->point_inside(other_bbox.boundary_points.second))
75 }
76 else
77 {
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 = {
82 {other_bbox.get_boundary_points().first,
83 other_bbox.get_boundary_points().second}};
84
85 // Step 1: testing if the boxes are close enough to intersect
86 for (unsigned int d = 0; d < spacedim; ++d)
87 if (bbox1[0][d] * (1 - std::numeric_limits<Number>::epsilon()) >
88 bbox2[1][d] ||
89 bbox2[0][d] * (1 - std::numeric_limits<Number>::epsilon()) >
90 bbox1[1][d])
92
93 // The boxes intersect: we need to understand now how they intersect.
94 // We begin by computing the intersection:
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)
98 {
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]);
101 }
102
103 // Finding the intersection's dimension
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() *
108 (std::abs(intersect_bbox_min[d]) +
109 std::abs(intersect_bbox_max[d])))
110 --intersect_dim;
111
112 if (intersect_dim == 0 || intersect_dim == spacedim - 2)
114
115 // Checking the two mergeable cases: first if the boxes are aligned so
116 // that they can be merged
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)
120 {
121 if (std::abs(bbox2[0][d] - bbox1[0][d]) >
122 std::numeric_limits<double>::epsilon() *
123 (std::abs(bbox2[0][d]) + std::abs(bbox1[0][d])))
124 ++not_align_1;
125 if (std::abs(bbox1[1][d] - bbox2[1][d]) >
126 std::numeric_limits<double>::epsilon() *
127 (std::abs(bbox1[1][d]) + std::abs(bbox2[1][d])))
128 ++not_align_2;
129 if (not_align_1 != not_align_2)
130 {
131 same_direction = false;
132 break;
133 }
134 }
135
136 if (not_align_1 <= 1 && not_align_2 <= 1 && same_direction)
138
139 // Second: one box is contained/equal to the other
140 if ((this->point_inside(bbox2[0]) && this->point_inside(bbox2[1])) ||
141 (other_bbox.point_inside(bbox1[0]) &&
142 other_bbox.point_inside(bbox1[1])))
144
145 // Degenerate and mergeable cases have been found, it remains:
147 }
148}
149
150
151
152template <int spacedim, typename Number>
153double
155{
156 double vol = 1.0;
157 for (unsigned int i = 0; i < spacedim; ++i)
158 vol *= (this->boundary_points.second[i] - this->boundary_points.first[i]);
159 return vol;
160}
161
162
163
164template <int spacedim, typename Number>
165Number
166BoundingBox<spacedim, Number>::lower_bound(const unsigned int direction) const
167{
168 AssertIndexRange(direction, spacedim);
169
170 return boundary_points.first[direction];
171}
172
173
174
175template <int spacedim, typename Number>
176Number
177BoundingBox<spacedim, Number>::upper_bound(const unsigned int direction) const
178{
179 AssertIndexRange(direction, spacedim);
180
181 return boundary_points.second[direction];
182}
183
184
185
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]);
193
194 return point;
195}
196
198
199template <int spacedim, typename Number>
201BoundingBox<spacedim, Number>::bounds(const unsigned int direction) const
202{
203 AssertIndexRange(direction, spacedim);
204
205 std::pair<Point<1, Number>, Point<1, Number>> lower_upper_bounds;
206 lower_upper_bounds.first[0] = lower_bound(direction);
207 lower_upper_bounds.second[0] = upper_bound(direction);
208
209 return BoundingBox<1, Number>(lower_upper_bounds);
210}
211
212
213
214template <int spacedim, typename Number>
215Number
216BoundingBox<spacedim, Number>::side_length(const unsigned int direction) const
217{
218 AssertIndexRange(direction, spacedim);
219
220 return boundary_points.second[direction] - boundary_points.first[direction];
221}
222
223
224
225template <int spacedim, typename Number>
227BoundingBox<spacedim, Number>::vertex(const unsigned int index) const
228{
230
231 const Point<spacedim> unit_cell_vertex =
235 for (unsigned int i = 0; i < spacedim; ++i)
236 point[i] = boundary_points.first[i] + side_length(i) * unit_cell_vertex[i];
237
238 return point;
240
241
242
243template <int spacedim, typename Number>
245BoundingBox<spacedim, Number>::child(const unsigned int index) const
246{
248
249 // Vertex closest to child.
250 const Point<spacedim, Number> parent_vertex = vertex(index);
251 const Point<spacedim, Number> parent_center = center();
252
253 const Point<spacedim> upper_corner_unit_cell =
256
257 const Point<spacedim> lower_corner_unit_cell =
259
260 std::pair<Point<spacedim, Number>, Point<spacedim, Number>>
261 child_lower_upper_corner;
262 for (unsigned int i = 0; i < spacedim; ++i)
263 {
264 const double child_side_length = side_length(i) / 2;
265
266 const double child_center = (parent_center[i] + parent_vertex[i]) / 2;
267
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);
272 }
273
274 return BoundingBox<spacedim, Number>(child_lower_upper_corner);
275}
276
277
278
279template <int spacedim, typename Number>
280BoundingBox<spacedim - 1, Number>
281BoundingBox<spacedim, Number>::cross_section(const unsigned int direction) const
282{
283 AssertIndexRange(direction, spacedim);
284
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)
288 {
289 const int index_to_write_from =
290 internal::coordinate_to_one_dim_higher<spacedim - 1>(direction, d);
291
292 cross_section_lower_upper_corner.first[d] =
293 boundary_points.first[index_to_write_from];
294
295 cross_section_lower_upper_corner.second[d] =
296 boundary_points.second[index_to_write_from];
297 }
298
299 return BoundingBox<spacedim - 1, Number>(cross_section_lower_upper_corner);
300}
301
302
304template <int spacedim, typename Number>
307 const Point<spacedim, Number> &point) const
308{
309 auto unit = point;
310 const auto diag = boundary_points.second - boundary_points.first;
311 unit -= boundary_points.first;
312 for (unsigned int d = 0; d < spacedim; ++d)
313 unit[d] /= diag[d];
314 return unit;
315}
316
317
318
319template <int spacedim, typename Number>
322 const Point<spacedim, Number> &point) const
323{
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];
328 return real;
329}
330
331
332
333template <int dim, typename Number>
336{
337 std::pair<Point<dim, Number>, Point<dim, Number>> lower_upper_corner;
338 for (unsigned int i = 0; i < dim; ++i)
339 {
340 lower_upper_corner.second[i] = 1;
341 }
342 return BoundingBox<dim, Number>(lower_upper_corner);
343}
344
345
346#include "bounding_box.inst"
NeighborType
Definition: bounding_box.h:33
BoundingBox< dim, Number > create_unit_bounding_box()
std::pair< Point< spacedim, Number >, Point< spacedim, Number > > boundary_points
Definition: bounding_box.h:315
BoundingBox< 1, Number > bounds(const unsigned int direction) const
NeighborType get_neighbor_type(const BoundingBox< spacedim, Number > &other_bbox) const
Definition: bounding_box.cc:64
Point< spacedim, Number > center() const
Number lower_bound(const unsigned int direction) const
void merge_with(const BoundingBox< spacedim, Number > &other_bbox)
Definition: bounding_box.cc:46
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
Definition: bounding_box.cc:23
double volume() 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
Definition: point.h:111
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
Point< 3 > center
#define AssertIndexRange(index, range)
Definition: exceptions.h:1732
::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)