Reference documentation for deal.II version 9.2.0
\(\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 - 2020 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 // ---------------------------------------------------------------------
17 
19 
20 template <int spacedim, typename Number>
21 bool
23  const Point<spacedim, Number> &p) const
24 {
25  for (unsigned int i = 0; i < spacedim; ++i)
26  {
27  // Bottom left-top right convention: the point is outside if it's smaller
28  // than the first or bigger than the second boundary point The bounding
29  // box is defined as a closed set
31  (std::abs(this->boundary_points.first[i] + p[i])) <
32  this->boundary_points.first[i] - p[i] ||
34  (std::abs(this->boundary_points.second[i] + p[i])) <
35  p[i] - this->boundary_points.second[i])
36  return false;
37  }
38  return true;
39 }
40 
41 template <int spacedim, typename Number>
42 void
44  const BoundingBox<spacedim, Number> &other_bbox)
45 {
46  for (unsigned int i = 0; i < spacedim; ++i)
47  {
48  this->boundary_points.first[i] =
49  std::min(this->boundary_points.first[i],
50  other_bbox.boundary_points.first[i]);
51  this->boundary_points.second[i] =
52  std::max(this->boundary_points.second[i],
53  other_bbox.boundary_points.second[i]);
54  }
55 }
56 
57 template <int spacedim, typename Number>
60  const BoundingBox<spacedim, Number> &other_bbox) const
61 {
62  if (spacedim == 1)
63  {
64  // In dimension 1 if the two bounding box are neighbors
65  // we can merge them
66  if (this->point_inside(other_bbox.boundary_points.first) ||
67  this->point_inside(other_bbox.boundary_points.second))
70  }
71  else
72  {
73  std::vector<Point<spacedim, Number>> bbox1;
74  bbox1.push_back(this->get_boundary_points().first);
75  bbox1.push_back(this->get_boundary_points().second);
76  std::vector<Point<spacedim, Number>> bbox2;
77  bbox2.push_back(other_bbox.get_boundary_points().first);
78  bbox2.push_back(other_bbox.get_boundary_points().second);
79 
80  // Step 1: testing if the boxes are close enough to intersect
81  for (unsigned int d = 0; d < spacedim; ++d)
82  if (bbox1[0][d] * (1 - std::numeric_limits<Number>::epsilon()) >
83  bbox2[1][d] ||
84  bbox2[0][d] * (1 - std::numeric_limits<Number>::epsilon()) >
85  bbox1[1][d])
87 
88  // The boxes intersect: we need to understand now how they intersect.
89  // We begin by computing the intersection:
90  std::vector<double> intersect_bbox_min;
91  std::vector<double> intersect_bbox_max;
92  for (unsigned int d = 0; d < spacedim; ++d)
93  {
94  intersect_bbox_min.push_back(std::max(bbox1[0][d], bbox2[0][d]));
95  intersect_bbox_max.push_back(std::min(bbox1[1][d], bbox2[1][d]));
96  }
97 
98  // Finding the intersection's dimension
99 
100  unsigned 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 + 2 == spacedim)
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 
148 template <int spacedim, typename Number>
149 std::pair<Point<spacedim, Number>, Point<spacedim, Number>> &
151 {
152  return this->boundary_points;
153 }
154 
155 
156 
157 template <int spacedim, typename Number>
158 const std::pair<Point<spacedim, Number>, Point<spacedim, Number>> &
160 {
161  return this->boundary_points;
162 }
163 
164 
165 
166 template <int spacedim, typename Number>
167 double
169 {
170  double vol = 1.0;
171  for (unsigned int i = 0; i < spacedim; ++i)
172  vol *= (this->boundary_points.second[i] - this->boundary_points.first[i]);
173  return vol;
174 }
175 
176 
177 
178 template <int spacedim, typename Number>
179 Number
180 BoundingBox<spacedim, Number>::lower_bound(const unsigned int direction) const
181 {
182  AssertIndexRange(direction, spacedim);
183 
184  return boundary_points.first[direction];
185 }
186 
187 
188 
189 template <int spacedim, typename Number>
190 Number
191 BoundingBox<spacedim, Number>::upper_bound(const unsigned int direction) const
192 {
193  AssertIndexRange(direction, spacedim);
194 
195  return boundary_points.second[direction];
196 }
197 
198 
199 
200 template <int spacedim, typename Number>
203 {
205  for (unsigned int i = 0; i < spacedim; ++i)
206  point[i] = .5 * (boundary_points.first[i] + boundary_points.second[i]);
207 
208  return point;
209 }
210 
211 
212 
213 template <int spacedim, typename Number>
215 BoundingBox<spacedim, Number>::bounds(const unsigned int direction) const
216 {
217  AssertIndexRange(direction, spacedim);
218 
219  std::pair<Point<1, Number>, Point<1, Number>> lower_upper_bounds;
220  lower_upper_bounds.first[0] = lower_bound(direction);
221  lower_upper_bounds.second[0] = upper_bound(direction);
222 
223  return BoundingBox<1, Number>(lower_upper_bounds);
224 }
225 
226 
227 
228 template <int spacedim, typename Number>
229 Number
230 BoundingBox<spacedim, Number>::side_length(const unsigned int direction) const
231 {
232  AssertIndexRange(direction, spacedim);
233 
234  return boundary_points.second[direction] - boundary_points.first[direction];
235 }
236 
237 
238 
239 template <int spacedim, typename Number>
241 BoundingBox<spacedim, Number>::vertex(const unsigned int index) const
242 {
244 
245  const Point<spacedim> unit_cell_vertex =
247 
249  for (unsigned int i = 0; i < spacedim; ++i)
250  point[i] = boundary_points.first[i] + side_length(i) * unit_cell_vertex[i];
251 
252  return point;
253 }
254 
255 
256 
257 template <int spacedim, typename Number>
259 BoundingBox<spacedim, Number>::child(const unsigned int index) const
260 {
262 
263  // Vertex closest to child.
264  const Point<spacedim, Number> parent_vertex = vertex(index);
265  const Point<spacedim, Number> parent_center = center();
266 
267  const Point<spacedim> upper_corner_unit_cell =
270 
271  const Point<spacedim> lower_corner_unit_cell =
273 
274  std::pair<Point<spacedim, Number>, Point<spacedim, Number>>
275  child_lower_upper_corner;
276  for (unsigned int i = 0; i < spacedim; ++i)
277  {
278  const double child_side_length = side_length(i) / 2;
279 
280  const double child_center = (parent_center[i] + parent_vertex[i]) / 2;
281 
282  child_lower_upper_corner.first[i] =
283  child_center + child_side_length * (lower_corner_unit_cell[i] - .5);
284  child_lower_upper_corner.second[i] =
285  child_center + child_side_length * (upper_corner_unit_cell[i] - .5);
286  }
287 
288  return BoundingBox<spacedim, Number>(child_lower_upper_corner);
289 }
290 
291 
292 
293 template <int spacedim, typename Number>
294 BoundingBox<spacedim - 1, Number>
295 BoundingBox<spacedim, Number>::cross_section(const unsigned int direction) const
296 {
297  AssertIndexRange(direction, spacedim);
298 
299  std::pair<Point<spacedim - 1, Number>, Point<spacedim - 1, Number>>
300  cross_section_lower_upper_corner;
301  for (unsigned int d = 0; d < spacedim - 1; ++d)
302  {
303  const int index_to_write_from =
304  internal::coordinate_to_one_dim_higher<spacedim - 1>(direction, d);
305 
306  cross_section_lower_upper_corner.first[d] =
307  boundary_points.first[index_to_write_from];
308 
309  cross_section_lower_upper_corner.second[d] =
310  boundary_points.second[index_to_write_from];
311  }
312 
313  return BoundingBox<spacedim - 1, Number>(cross_section_lower_upper_corner);
314 }
315 
316 
317 
318 template <int dim, typename Number>
321 {
322  std::pair<Point<dim, Number>, Point<dim, Number>> lower_upper_corner;
323  for (unsigned int i = 0; i < dim; ++i)
324  {
325  lower_upper_corner.second[i] = 1;
326  }
327  return BoundingBox<dim, Number>(lower_upper_corner);
328 }
329 
330 
331 #include "bounding_box.inst"
Physics::Elasticity::Kinematics::epsilon
SymmetricTensor< 2, dim, Number > epsilon(const Tensor< 2, dim, Number > &Grad_u)
bounding_box.h
BoundingBox::volume
double volume() const
Definition: bounding_box.cc:168
GeometryInfo::unit_cell_vertex
static Point< dim > unit_cell_vertex(const unsigned int vertex)
BoundingBox::lower_bound
Number lower_bound(const unsigned int direction) const
Definition: bounding_box.cc:180
BoundingBox::child
BoundingBox< spacedim, Number > child(const unsigned int index) const
Definition: bounding_box.cc:259
GeometryInfo
Definition: geometry_info.h:1224
AssertIndexRange
#define AssertIndexRange(index, range)
Definition: exceptions.h:1649
NeighborType::attached_neighbors
@ attached_neighbors
NeighborType
NeighborType
Definition: bounding_box.h:37
NeighborType::mergeable_neighbors
@ mergeable_neighbors
BoundingBox::center
Point< spacedim, Number > center() const
Definition: bounding_box.cc:202
second
Point< 2 > second
Definition: grid_out.cc:4353
Physics::Elasticity::Kinematics::d
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
OpenCASCADE::point
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition: utilities.cc:188
BoundingBox
Definition: bounding_box.h:128
NeighborType::not_neighbors
@ not_neighbors
BoundingBox::merge_with
void merge_with(const BoundingBox< spacedim, Number > &other_bbox)
Definition: bounding_box.cc:43
geometry_info.h
BoundingBox::upper_bound
Number upper_bound(const unsigned int direction) const
Definition: bounding_box.cc:191
Utilities::lower_bound
Iterator lower_bound(Iterator first, Iterator last, const T &val)
Definition: utilities.h:1102
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
BoundingBox::cross_section
BoundingBox< spacedim - 1, Number > cross_section(const unsigned int direction) const
Definition: bounding_box.cc:295
BoundingBox::create_unit_bounding_box
BoundingBox< dim, Number > create_unit_bounding_box()
Definition: bounding_box.cc:320
BoundingBox::point_inside
bool point_inside(const Point< spacedim, Number > &p) const
Definition: bounding_box.cc:22
BoundingBox::boundary_points
std::pair< Point< spacedim, Number >, Point< spacedim, Number > > boundary_points
Definition: bounding_box.h:271
BoundingBox::get_neighbor_type
NeighborType get_neighbor_type(const BoundingBox< spacedim, Number > &other_bbox) const
Definition: bounding_box.cc:59
Utilities::MPI::min
T min(const T &t, const MPI_Comm &mpi_communicator)
Point
Definition: point.h:111
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
BoundingBox::side_length
Number side_length(const unsigned int direction) const
Definition: bounding_box.cc:230
first
Point< 2 > first
Definition: grid_out.cc:4352
BoundingBox::bounds
BoundingBox< 1, Number > bounds(const unsigned int direction) const
Definition: bounding_box.cc:215
BoundingBox::get_boundary_points
std::pair< Point< spacedim, Number >, Point< spacedim, Number > > & get_boundary_points()
Definition: bounding_box.cc:150
BoundingBox::vertex
Point< spacedim, Number > vertex(const unsigned int index) const
Definition: bounding_box.cc:241
center
Point< 3 > center
Definition: data_out_base.cc:169
NeighborType::simple_neighbors
@ simple_neighbors
Utilities::MPI::max
T max(const T &t, const MPI_Comm &mpi_communicator)