Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
tria_objects.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2006 - 2018 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 
16 #include <deal.II/base/memory_consumption.h>
17 
18 #include <deal.II/grid/tria.h>
19 #include <deal.II/grid/tria_accessor.h>
20 #include <deal.II/grid/tria_iterator.h>
21 #include <deal.II/grid/tria_objects.h>
22 
23 #include <algorithm>
24 #include <functional>
25 
26 
27 
28 DEAL_II_NAMESPACE_OPEN
29 
30 namespace internal
31 {
32  namespace TriangulationImplementation
33  {
34  template <class G>
35  void
36  TriaObjects<G>::reserve_space(const unsigned int new_objects_in_pairs,
37  const unsigned int new_objects_single)
38  {
39  Assert(new_objects_in_pairs % 2 == 0, ExcInternalError());
40 
41  next_free_single = 0;
42  next_free_pair = 0;
43  reverse_order_next_free_single = false;
44 
45  // count the number of objects, of unused single objects and of
46  // unused pairs of objects
47  unsigned int n_objects = 0;
48  unsigned int n_unused_pairs = 0;
49  unsigned int n_unused_singles = 0;
50  for (unsigned int i = 0; i < used.size(); ++i)
51  {
52  if (used[i])
53  ++n_objects;
54  else if (i + 1 < used.size())
55  {
56  if (used[i + 1])
57  {
58  ++n_unused_singles;
59  if (next_free_single == 0)
60  next_free_single = i;
61  }
62  else
63  {
64  ++n_unused_pairs;
65  if (next_free_pair == 0)
66  next_free_pair = i;
67  ++i;
68  }
69  }
70  else
71  ++n_unused_singles;
72  }
73  Assert(n_objects + 2 * n_unused_pairs + n_unused_singles == used.size(),
75 
76  // how many single objects are needed in addition to
77  // n_unused_objects?
78  const int additional_single_objects =
79  new_objects_single - n_unused_singles;
80 
81  unsigned int new_size =
82  used.size() + new_objects_in_pairs - 2 * n_unused_pairs;
83  if (additional_single_objects > 0)
84  new_size += additional_single_objects;
85 
86  // only allocate space if necessary
87  if (new_size > cells.size())
88  {
89  cells.reserve(new_size);
90  cells.insert(cells.end(), new_size - cells.size(), G());
91 
92  used.reserve(new_size);
93  used.insert(used.end(), new_size - used.size(), false);
94 
95  user_flags.reserve(new_size);
96  user_flags.insert(user_flags.end(),
97  new_size - user_flags.size(),
98  false);
99 
100  const unsigned int factor =
102  children.reserve(factor * new_size);
103  children.insert(children.end(),
104  factor * new_size - children.size(),
105  -1);
106 
107  if (G::dimension > 1)
108  {
109  refinement_cases.reserve(new_size);
110  refinement_cases.insert(
111  refinement_cases.end(),
112  new_size - refinement_cases.size(),
114  }
115 
116  // first reserve, then resize. Otherwise the std library can decide to
117  // allocate more entries.
118  boundary_or_material_id.reserve(new_size);
119  boundary_or_material_id.resize(new_size);
120 
121  user_data.reserve(new_size);
122  user_data.resize(new_size);
123 
124  manifold_id.reserve(new_size);
125  manifold_id.insert(manifold_id.end(),
126  new_size - manifold_id.size(),
128  }
129 
130  if (n_unused_singles == 0)
131  {
132  next_free_single = new_size - 1;
133  reverse_order_next_free_single = true;
134  }
135  }
136 
137 
138  template <>
139  template <int dim, int spacedim>
140  typename ::Triangulation<dim, spacedim>::raw_hex_iterator
141  TriaObjects<TriaObject<3>>::next_free_hex(
142  const ::Triangulation<dim, spacedim> &tria,
143  const unsigned int level)
144  {
145  // TODO: Think of a way to ensure that we are using the correct
146  // triangulation, i.e. the one containing *this.
147 
148  int pos = next_free_pair, last = used.size() - 1;
149  for (; pos < last; ++pos)
150  if (!used[pos])
151  {
152  // this should be a pair slot
153  Assert(!used[pos + 1], ExcInternalError());
154  break;
155  }
156  if (pos >= last)
157  // no free slot
158  return tria.end_hex();
159  else
160  next_free_pair = pos + 2;
161 
162  return
163  typename ::Triangulation<dim, spacedim>::raw_hex_iterator(&tria,
164  level,
165  pos);
166  }
167 
168 
169  void
170  TriaObjectsHex::reserve_space(const unsigned int new_hexes)
171  {
172  const unsigned int new_size =
173  new_hexes + std::count_if(used.begin(),
174  used.end(),
175  std::bind(std::equal_to<bool>(),
176  std::placeholders::_1,
177  true));
178 
179  // see above...
180  if (new_size > cells.size())
181  {
182  cells.reserve(new_size);
183  cells.insert(cells.end(), new_size - cells.size(), TriaObject<3>());
184 
185  used.reserve(new_size);
186  used.insert(used.end(), new_size - used.size(), false);
187 
188  user_flags.reserve(new_size);
189  user_flags.insert(user_flags.end(),
190  new_size - user_flags.size(),
191  false);
192 
193  children.reserve(4 * new_size);
194  children.insert(children.end(), 4 * new_size - children.size(), -1);
195 
196  // for the following fields, we know exactly how many elements
197  // we need, so first reserve then resize (resize itself, at least
198  // with some compiler libraries, appears to round up the size it
199  // actually reserves)
200  boundary_or_material_id.reserve(new_size);
201  boundary_or_material_id.resize(new_size);
202 
203  manifold_id.reserve(new_size);
204  manifold_id.insert(manifold_id.end(),
205  new_size - manifold_id.size(),
207 
208  user_data.reserve(new_size);
209  user_data.resize(new_size);
210 
214  face_orientations.size(),
215  true);
216 
217  refinement_cases.reserve(new_size);
218  refinement_cases.insert(refinement_cases.end(),
219  new_size - refinement_cases.size(),
221 
222  face_flips.reserve(new_size * GeometryInfo<3>::faces_per_cell);
223  face_flips.insert(face_flips.end(),
225  face_flips.size(),
226  false);
228  face_rotations.insert(face_rotations.end(),
230  face_rotations.size(),
231  false);
232  }
234  }
235 
236 
237  void
238  TriaObjectsQuad3D::reserve_space(const unsigned int new_quads_in_pairs,
239  const unsigned int new_quads_single)
240  {
241  Assert(new_quads_in_pairs % 2 == 0, ExcInternalError());
242 
243  next_free_single = 0;
244  next_free_pair = 0;
246 
247  // count the number of objects, of unused single objects and of
248  // unused pairs of objects
249  unsigned int n_quads = 0;
250  unsigned int n_unused_pairs = 0;
251  unsigned int n_unused_singles = 0;
252  for (unsigned int i = 0; i < used.size(); ++i)
253  {
254  if (used[i])
255  ++n_quads;
256  else if (i + 1 < used.size())
257  {
258  if (used[i + 1])
259  {
260  ++n_unused_singles;
261  if (next_free_single == 0)
262  next_free_single = i;
263  }
264  else
265  {
266  ++n_unused_pairs;
267  if (next_free_pair == 0)
268  next_free_pair = i;
269  ++i;
270  }
271  }
272  else
273  ++n_unused_singles;
274  }
275  Assert(n_quads + 2 * n_unused_pairs + n_unused_singles == used.size(),
276  ExcInternalError());
277 
278  // how many single quads are needed in addition to n_unused_quads?
279  const int additional_single_quads = new_quads_single - n_unused_singles;
280 
281  unsigned int new_size =
282  used.size() + new_quads_in_pairs - 2 * n_unused_pairs;
283  if (additional_single_quads > 0)
284  new_size += additional_single_quads;
285 
286  // see above...
287  if (new_size > cells.size())
288  {
289  // reseve space for the base class
290  TriaObjects<TriaObject<2>>::reserve_space(new_quads_in_pairs,
291  new_quads_single);
292  // reserve the field of the derived class
296  line_orientations.size(),
297  true);
298  }
299 
300  if (n_unused_singles == 0)
301  {
302  next_free_single = new_size - 1;
304  }
305  }
306 
307 
308  template <>
309  void
310  TriaObjects<TriaObject<1>>::monitor_memory(const unsigned int) const
311  {
312  Assert(cells.size() == used.size(),
313  ExcMemoryInexact(cells.size(), used.size()));
314  Assert(cells.size() == user_flags.size(),
315  ExcMemoryInexact(cells.size(), user_flags.size()));
316  Assert(cells.size() == children.size(),
317  ExcMemoryInexact(cells.size(), children.size()));
318  Assert(cells.size() == boundary_or_material_id.size(),
319  ExcMemoryInexact(cells.size(), boundary_or_material_id.size()));
320  Assert(cells.size() == manifold_id.size(),
321  ExcMemoryInexact(cells.size(), manifold_id.size()));
322  Assert(cells.size() == user_data.size(),
323  ExcMemoryInexact(cells.size(), user_data.size()));
324  }
325 
326 
327  template <>
328  void
329  TriaObjects<TriaObject<2>>::monitor_memory(const unsigned int) const
330  {
331  Assert(cells.size() == used.size(),
332  ExcMemoryInexact(cells.size(), used.size()));
333  Assert(cells.size() == user_flags.size(),
334  ExcMemoryInexact(cells.size(), user_flags.size()));
335  Assert(2 * cells.size() == children.size(),
336  ExcMemoryInexact(cells.size(), children.size()));
337  Assert(cells.size() == refinement_cases.size(),
338  ExcMemoryInexact(cells.size(), refinement_cases.size()));
339  Assert(cells.size() == boundary_or_material_id.size(),
340  ExcMemoryInexact(cells.size(), boundary_or_material_id.size()));
341  Assert(cells.size() == manifold_id.size(),
342  ExcMemoryInexact(cells.size(), manifold_id.size()));
343  Assert(cells.size() == user_data.size(),
344  ExcMemoryInexact(cells.size(), user_data.size()));
345  }
346 
347 
348  void
349  TriaObjectsHex::monitor_memory(const unsigned int) const
350  {
351  Assert(cells.size() == used.size(),
352  ExcMemoryInexact(cells.size(), used.size()));
353  Assert(cells.size() == user_flags.size(),
354  ExcMemoryInexact(cells.size(), user_flags.size()));
355  Assert(4 * cells.size() == children.size(),
356  ExcMemoryInexact(cells.size(), children.size()));
357  Assert(cells.size() == boundary_or_material_id.size(),
359  Assert(cells.size() == manifold_id.size(),
360  ExcMemoryInexact(cells.size(), manifold_id.size()));
361  Assert(cells.size() == user_data.size(),
362  ExcMemoryInexact(cells.size(), user_data.size()));
364  face_orientations.size(),
366  face_orientations.size()));
368  face_flips.size(),
370  face_flips.size()));
372  face_rotations.size(),
374  face_rotations.size()));
375  }
376 
377 
378  void
379  TriaObjectsQuad3D::monitor_memory(const unsigned int) const
380  {
381  // check that we have not allocated too much memory. note that bool
382  // vectors allocate their memory in chunks of whole integers, so they
383  // may over-allocate by up to as many elements as an integer has bits
385  line_orientations.size(),
387  line_orientations.size()));
389  }
390 
391 
392  template <typename G>
393  void
395  {
396  cells.clear();
397  children.clear();
398  refinement_cases.clear();
399  used.clear();
400  user_flags.clear();
401  boundary_or_material_id.clear();
402  manifold_id.clear();
403  user_data.clear();
404  user_data_type = data_unknown;
405  }
406 
407 
408  void
410  {
412  face_orientations.clear();
413  face_flips.clear();
414  face_rotations.clear();
415  }
416 
417 
418  void
420  {
422  line_orientations.clear();
423  }
424 
425 
426  template <typename G>
427  std::size_t
429  {
434  MemoryConsumption::memory_consumption(boundary_or_material_id) +
436  MemoryConsumption::memory_consumption(refinement_cases) +
437  user_data.capacity() * sizeof(UserData) + sizeof(user_data));
438  }
439 
440 
441  std::size_t
443  {
448  }
449 
450 
451  std::size_t
453  {
456  }
457 
458 
459 
460  // explicit instantiations
461  template class TriaObjects<TriaObject<1>>;
462  template class TriaObjects<TriaObject<2>>;
463 
464 #include "tria_objects.inst"
465  } // namespace TriangulationImplementation
466 } // namespace internal
467 
468 DEAL_II_NAMESPACE_CLOSE
void reserve_space(const unsigned int new_objs)
const types::manifold_id flat_manifold_id
Definition: types.h:246
static ::ExceptionBase & ExcMemoryInexact(int arg1, int arg2)
void monitor_memory(const unsigned int true_dimension) const
void monitor_memory(const unsigned int true_dimension) const
#define Assert(cond, exc)
Definition: exceptions.h:1407
void reserve_space(const unsigned int new_objs_in_pairs, const unsigned int new_objs_single=0)
Definition: tria_objects.cc:36
std::vector< RefinementCase< TriaObject< 3 > ::dimension > > refinement_cases
Definition: tria_objects.h:98
void reserve_space(const unsigned int new_quads_in_pairs, const unsigned int new_quads_single=0)
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
static ::ExceptionBase & ExcInternalError()