Reference documentation for deal.II version 9.3.3
\(\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\}}\)
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
42template <int spacedim, typename Number>
43void
45 const BoundingBox<spacedim, Number> &other_bbox)
46{
47 for (unsigned int i = 0; i < spacedim; ++i)
48 {
49 this->boundary_points.first[i] =
50 std::min(this->boundary_points.first[i],
51 other_bbox.boundary_points.first[i]);
52 this->boundary_points.second[i] =
53 std::max(this->boundary_points.second[i],
54 other_bbox.boundary_points.second[i]);
55 }
56}
57
58template <int spacedim, typename Number>
61 const BoundingBox<spacedim, Number> &other_bbox) const
62{
63 if (spacedim == 1)
64 {
65 // In dimension 1 if the two bounding box are neighbors
66 // we can merge them
67 if (this->point_inside(other_bbox.boundary_points.first) ||
68 this->point_inside(other_bbox.boundary_points.second))
71 }
72 else
73 {
74 const std::array<Point<spacedim, Number>, 2> bbox1 = {
75 {this->get_boundary_points().first,
76 this->get_boundary_points().second}};
77 const std::array<Point<spacedim, Number>, 2> bbox2 = {
78 {other_bbox.get_boundary_points().first,
79 other_bbox.get_boundary_points().second}};
80
81 // Step 1: testing if the boxes are close enough to intersect
82 for (unsigned int d = 0; d < spacedim; ++d)
83 if (bbox1[0][d] * (1 - std::numeric_limits<Number>::epsilon()) >
84 bbox2[1][d] ||
85 bbox2[0][d] * (1 - std::numeric_limits<Number>::epsilon()) >
86 bbox1[1][d])
88
89 // The boxes intersect: we need to understand now how they intersect.
90 // We begin by computing the intersection:
91 std::array<double, spacedim> intersect_bbox_min;
92 std::array<double, spacedim> intersect_bbox_max;
93 for (unsigned int d = 0; d < spacedim; ++d)
94 {
95 intersect_bbox_min[d] = std::max(bbox1[0][d], bbox2[0][d]);
96 intersect_bbox_max[d] = std::min(bbox1[1][d], bbox2[1][d]);
97 }
98
99 // Finding the intersection's dimension
100 int intersect_dim = spacedim;
101 for (unsigned int d = 0; d < spacedim; ++d)
102 if (std::abs(intersect_bbox_min[d] - intersect_bbox_max[d]) <=
104 (std::abs(intersect_bbox_min[d]) +
105 std::abs(intersect_bbox_max[d])))
106 --intersect_dim;
107
108 if (intersect_dim == 0 || intersect_dim == spacedim - 2)
110
111 // Checking the two mergeable cases: first if the boxes are aligned so
112 // that they can be merged
113 unsigned int not_align_1 = 0, not_align_2 = 0;
114 bool same_direction = true;
115 for (unsigned int d = 0; d < spacedim; ++d)
116 {
117 if (std::abs(bbox2[0][d] - bbox1[0][d]) >
119 (std::abs(bbox2[0][d]) + std::abs(bbox1[0][d])))
120 ++not_align_1;
121 if (std::abs(bbox1[1][d] - bbox2[1][d]) >
123 (std::abs(bbox1[1][d]) + std::abs(bbox2[1][d])))
124 ++not_align_2;
125 if (not_align_1 != not_align_2)
126 {
127 same_direction = false;
128 break;
129 }
130 }
131
132 if (not_align_1 <= 1 && not_align_2 <= 1 && same_direction)
134
135 // Second: one box is contained/equal to the other
136 if ((this->point_inside(bbox2[0]) && this->point_inside(bbox2[1])) ||
137 (other_bbox.point_inside(bbox1[0]) &&
138 other_bbox.point_inside(bbox1[1])))
140
141 // Degenerate and mergeable cases have been found, it remains:
143 }
144}
145
146
147
148template <int spacedim, typename Number>
149double
151{
152 double vol = 1.0;
153 for (unsigned int i = 0; i < spacedim; ++i)
154 vol *= (this->boundary_points.second[i] - this->boundary_points.first[i]);
155 return vol;
156}
157
158
159
160template <int spacedim, typename Number>
161Number
162BoundingBox<spacedim, Number>::lower_bound(const unsigned int direction) const
163{
164 AssertIndexRange(direction, spacedim);
165
166 return boundary_points.first[direction];
167}
168
169
170
171template <int spacedim, typename Number>
172Number
173BoundingBox<spacedim, Number>::upper_bound(const unsigned int direction) const
174{
175 AssertIndexRange(direction, spacedim);
176
177 return boundary_points.second[direction];
178}
180
181
182template <int spacedim, typename Number>
185{
187 for (unsigned int i = 0; i < spacedim; ++i)
188 point[i] = .5 * (boundary_points.first[i] + boundary_points.second[i]);
189
190 return point;
191}
192
193
194
195template <int spacedim, typename Number>
197BoundingBox<spacedim, Number>::bounds(const unsigned int direction) const
198{
199 AssertIndexRange(direction, spacedim);
200
201 std::pair<Point<1, Number>, Point<1, Number>> lower_upper_bounds;
202 lower_upper_bounds.first[0] = lower_bound(direction);
203 lower_upper_bounds.second[0] = upper_bound(direction);
204
205 return BoundingBox<1, Number>(lower_upper_bounds);
206}
207
208
209
210template <int spacedim, typename Number>
211Number
212BoundingBox<spacedim, Number>::side_length(const unsigned int direction) const
213{
214 AssertIndexRange(direction, spacedim);
215
216 return boundary_points.second[direction] - boundary_points.first[direction];
218
219
220
221template <int spacedim, typename Number>
223BoundingBox<spacedim, Number>::vertex(const unsigned int index) const
224{
226
227 const Point<spacedim> unit_cell_vertex =
231 for (unsigned int i = 0; i < spacedim; ++i)
232 point[i] = boundary_points.first[i] + side_length(i) * unit_cell_vertex[i];
233
234 return point;
236
237
238
239template <int spacedim, typename Number>
241BoundingBox<spacedim, Number>::child(const unsigned int index) const
242{
244
245 // Vertex closest to child.
246 const Point<spacedim, Number> parent_vertex = vertex(index);
247 const Point<spacedim, Number> parent_center = center();
248
249 const Point<spacedim> upper_corner_unit_cell =
252
253 const Point<spacedim> lower_corner_unit_cell =
255
256 std::pair<Point<spacedim, Number>, Point<spacedim, Number>>
257 child_lower_upper_corner;
258 for (unsigned int i = 0; i < spacedim; ++i)
259 {
260 const double child_side_length = side_length(i) / 2;
262 const double child_center = (parent_center[i] + parent_vertex[i]) / 2;
263
264 child_lower_upper_corner.first[i] =
265 child_center + child_side_length * (lower_corner_unit_cell[i] - .5);
266 child_lower_upper_corner.second[i] =
267 child_center + child_side_length * (upper_corner_unit_cell[i] - .5);
268 }
269
270 return BoundingBox<spacedim, Number>(child_lower_upper_corner);
272
273
274
275template <int spacedim, typename Number>
276BoundingBox<spacedim - 1, Number>
277BoundingBox<spacedim, Number>::cross_section(const unsigned int direction) const
278{
279 AssertIndexRange(direction, spacedim);
280
281 std::pair<Point<spacedim - 1, Number>, Point<spacedim - 1, Number>>
282 cross_section_lower_upper_corner;
283 for (unsigned int d = 0; d < spacedim - 1; ++d)
284 {
285 const int index_to_write_from =
286 internal::coordinate_to_one_dim_higher<spacedim - 1>(direction, d);
287
288 cross_section_lower_upper_corner.first[d] =
289 boundary_points.first[index_to_write_from];
290
291 cross_section_lower_upper_corner.second[d] =
292 boundary_points.second[index_to_write_from];
294
295 return BoundingBox<spacedim - 1, Number>(cross_section_lower_upper_corner);
296}
297
298
299
300template <int spacedim, typename Number>
303 const Point<spacedim, Number> &point) const
304{
305 auto unit = point;
306 const auto diag = boundary_points.second - boundary_points.first;
307 unit -= boundary_points.first;
308 for (unsigned int d = 0; d < spacedim; ++d)
309 unit[d] /= diag[d];
310 return unit;
311}
312
313
314
315template <int spacedim, typename Number>
318 const Point<spacedim, Number> &point) const
319{
320 auto real = boundary_points.first;
321 const auto diag = boundary_points.second - boundary_points.first;
322 for (unsigned int d = 0; d < spacedim; ++d)
323 real[d] += diag[d] * point[d];
324 return real;
325}
326
327
328
329template <int dim, typename Number>
332{
333 std::pair<Point<dim, Number>, Point<dim, Number>> lower_upper_corner;
334 for (unsigned int i = 0; i < dim; ++i)
335 {
336 lower_upper_corner.second[i] = 1;
337 }
338 return BoundingBox<dim, Number>(lower_upper_corner);
339}
340
341
342#include "bounding_box.inst"
NeighborType
Definition: bounding_box.h:33
std::pair< Point< spacedim, Number >, Point< spacedim, Number > > boundary_points
Definition: bounding_box.h:305
BoundingBox< 1, Number > bounds(const unsigned int direction) const
BoundingBox< dim, Number > create_unit_bounding_box()
NeighborType get_neighbor_type(const BoundingBox< spacedim, Number > &other_bbox) const
Definition: bounding_box.cc:60
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:44
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:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
Point< 3 > center
#define AssertIndexRange(index, range)
Definition: exceptions.h:1690
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition: utilities.cc:188
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)
Iterator lower_bound(Iterator first, Iterator last, const T &val)
Definition: utilities.h:1111
::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)