Reference documentation for deal.II version 8.5.1
tria_objects.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2006 - 2016 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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/base/memory_consumption.h>
17 #include <deal.II/base/std_cxx11/bind.h>
18 #include <deal.II/grid/tria_objects.h>
19 #include <deal.II/grid/tria.h>
20 #include <deal.II/grid/tria_iterator.h>
21 #include <deal.II/grid/tria_accessor.h>
22 
23 #include <algorithm>
24 #include <functional>
25 
26 
27 
28 DEAL_II_NAMESPACE_OPEN
29 
30 namespace internal
31 {
32  namespace Triangulation
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(),
91  new_size-cells.size(),
92  G ());
93 
94  used.reserve (new_size);
95  used.insert (used.end(),
96  new_size-used.size(),
97  false);
98 
99  user_flags.reserve (new_size);
100  user_flags.insert (user_flags.end(),
101  new_size-user_flags.size(),
102  false);
103 
104  const unsigned int factor = GeometryInfo<G::dimension>::max_children_per_cell / 2;
105  children.reserve (factor*new_size);
106  children.insert (children.end(),
107  factor*new_size-children.size(),
108  -1);
109 
110  if (G::dimension > 1)
111  {
112  refinement_cases.reserve (new_size);
113  refinement_cases.insert (refinement_cases.end(),
114  new_size - refinement_cases.size(),
116  }
117 
118  // first reserve, then resize. Otherwise the std library can decide to allocate
119  // more entries.
120  boundary_or_material_id.reserve (new_size);
121  boundary_or_material_id.resize (new_size);
122 
123  user_data.reserve (new_size);
124  user_data.resize (new_size);
125 
126  manifold_id.reserve (new_size);
127  manifold_id.insert (manifold_id.end(),
128  new_size-manifold_id.size(),
130 
131  }
132 
133  if (n_unused_singles==0)
134  {
135  next_free_single=new_size-1;
136  reverse_order_next_free_single=true;
137  }
138  }
139 
140 
141  template <>
142  template <int dim, int spacedim>
143  typename ::Triangulation<dim,spacedim>::raw_hex_iterator
144  TriaObjects<TriaObject<3> >::next_free_hex (const ::Triangulation<dim,spacedim> &tria,
145  const unsigned int level)
146  {
147  // TODO: Think of a way to ensure that we are using the correct triangulation, i.e. the one containing *this.
148 
149  int pos=next_free_pair,
150  last=used.size()-1;
151  for (; pos<last; ++pos)
152  if (!used[pos])
153  {
154  // this should be a pair slot
155  Assert(!used[pos+1], ExcInternalError());
156  break;
157  }
158  if (pos>=last)
159  // no free slot
160  return tria.end_hex();
161  else
162  next_free_pair=pos+2;
163 
164  return typename ::Triangulation<dim,spacedim>::raw_hex_iterator(&tria,level,pos);
165  }
166 
167 
168  void
169  TriaObjectsHex::reserve_space (const unsigned int new_hexes)
170  {
171  const unsigned int new_size = new_hexes +
172  std::count_if (used.begin(),
173  used.end(),
174  std_cxx11::bind (std::equal_to<bool>(), std_cxx11::_1, true));
175 
176  // see above...
177  if (new_size>cells.size())
178  {
179  cells.reserve (new_size);
180  cells.insert (cells.end(),
181  new_size-cells.size(),
182  TriaObject<3> ());
183 
184  used.reserve (new_size);
185  used.insert (used.end(),
186  new_size-used.size(),
187  false);
188 
189  user_flags.reserve (new_size);
190  user_flags.insert (user_flags.end(),
191  new_size-user_flags.size(),
192  false);
193 
194  children.reserve (4*new_size);
195  children.insert (children.end(),
196  4*new_size-children.size(),
197  -1);
198 
199  // for the following fields, we know exactly how many elements
200  // we need, so first reserve then resize (resize itself, at least
201  // with some compiler libraries, appears to round up the size it
202  // actually reserves)
203  boundary_or_material_id.reserve (new_size);
204  boundary_or_material_id.resize (new_size);
205 
206  manifold_id.reserve (new_size);
207  manifold_id.insert (manifold_id.end(),
208  new_size-manifold_id.size(),
210 
211  user_data.reserve (new_size);
212  user_data.resize (new_size);
213 
215  face_orientations.insert (face_orientations.end(),
217  - face_orientations.size(),
218  true);
219 
220  refinement_cases.reserve (new_size);
221  refinement_cases.insert (refinement_cases.end(),
222  new_size-refinement_cases.size(),
224 
225  face_flips.reserve (new_size * GeometryInfo<3>::faces_per_cell);
226  face_flips.insert (face_flips.end(),
228  - face_flips.size(),
229  false);
230  face_rotations.reserve (new_size * GeometryInfo<3>::faces_per_cell);
231  face_rotations.insert (face_rotations.end(),
233  - face_rotations.size(),
234  false);
235  }
237  }
238 
239 
240  void
241  TriaObjectsQuad3D::reserve_space (const unsigned int new_quads_in_pairs,
242  const unsigned int new_quads_single)
243  {
244  Assert(new_quads_in_pairs%2==0, ExcInternalError());
245 
247  next_free_pair=0;
249 
250  // count the number of objects, of unused single objects and of
251  // unused pairs of objects
252  unsigned int n_quads=0;
253  unsigned int n_unused_pairs=0;
254  unsigned int n_unused_singles=0;
255  for (unsigned int i=0; i<used.size(); ++i)
256  {
257  if (used[i])
258  ++n_quads;
259  else if (i+1<used.size())
260  {
261  if (used[i+1])
262  {
263  ++n_unused_singles;
264  if (next_free_single==0)
266  }
267  else
268  {
269  ++n_unused_pairs;
270  if (next_free_pair==0)
271  next_free_pair=i;
272  ++i;
273  }
274  }
275  else
276  ++n_unused_singles;
277  }
278  Assert(n_quads+2*n_unused_pairs+n_unused_singles==used.size(),
279  ExcInternalError());
280 
281  // how many single quads are needed in addition to n_unused_quads?
282  const int additional_single_quads=
283  new_quads_single-n_unused_singles;
284 
285  unsigned int new_size=
286  used.size() + new_quads_in_pairs - 2*n_unused_pairs;
287  if (additional_single_quads>0)
288  new_size+=additional_single_quads;
289 
290  // see above...
291  if (new_size>cells.size())
292  {
293  // reseve space for the base class
294  TriaObjects<TriaObject<2> >::reserve_space(new_quads_in_pairs,new_quads_single);
295  // reserve the field of the derived class
297  line_orientations.insert (line_orientations.end(),
299  - line_orientations.size(),
300  true);
301  }
302 
303  if (n_unused_singles==0)
304  {
305  next_free_single=new_size-1;
307  }
308  }
309 
310 
311  template<>
312  void
313  TriaObjects<TriaObject<1> >::monitor_memory (const unsigned int) const
314  {
315  Assert (cells.size() == used.size(),
316  ExcMemoryInexact (cells.size(), used.size()));
317  Assert (cells.size() == user_flags.size(),
318  ExcMemoryInexact (cells.size(), user_flags.size()));
319  Assert (cells.size() == children.size(),
320  ExcMemoryInexact (cells.size(), children.size()));
321  Assert (cells.size() == boundary_or_material_id.size(),
322  ExcMemoryInexact (cells.size(), boundary_or_material_id.size()));
323  Assert (cells.size() == manifold_id.size(),
324  ExcMemoryInexact (cells.size(), manifold_id.size()));
325  Assert (cells.size() == user_data.size(),
326  ExcMemoryInexact (cells.size(), user_data.size()));
327  }
328 
329 
330  template<>
331  void
332  TriaObjects<TriaObject<2> >::monitor_memory (const unsigned int) const
333  {
334  Assert (cells.size() == used.size(),
335  ExcMemoryInexact (cells.size(), used.size()));
336  Assert (cells.size() == user_flags.size(),
337  ExcMemoryInexact (cells.size(), user_flags.size()));
338  Assert (2*cells.size() == children.size(),
339  ExcMemoryInexact (cells.size(), children.size()));
340  Assert (cells.size() == refinement_cases.size(),
341  ExcMemoryInexact (cells.size(), refinement_cases.size()));
342  Assert (cells.size() == boundary_or_material_id.size(),
343  ExcMemoryInexact (cells.size(), boundary_or_material_id.size()));
344  Assert (cells.size() == manifold_id.size(),
345  ExcMemoryInexact (cells.size(), manifold_id.size()));
346  Assert (cells.size() == user_data.size(),
347  ExcMemoryInexact (cells.size(), user_data.size()));
348  }
349 
350 
351  void
352  TriaObjectsHex::monitor_memory (const unsigned int) const
353  {
354  Assert (cells.size() == used.size(),
355  ExcMemoryInexact (cells.size(), used.size()));
356  Assert (cells.size() == user_flags.size(),
357  ExcMemoryInexact (cells.size(), user_flags.size()));
358  Assert (4*cells.size() == children.size(),
359  ExcMemoryInexact (cells.size(), children.size()));
360  Assert (cells.size() == boundary_or_material_id.size(),
362  Assert (cells.size() == manifold_id.size(),
363  ExcMemoryInexact (cells.size(), manifold_id.size()));
364  Assert (cells.size() == user_data.size(),
365  ExcMemoryInexact (cells.size(), user_data.size()));
367  == face_orientations.size(),
369  face_orientations.size()));
371  == face_flips.size(),
373  face_flips.size()));
375  == face_rotations.size(),
377  face_rotations.size()));
378  }
379 
380 
381  void
382  TriaObjectsQuad3D::monitor_memory (const unsigned int) const
383  {
384  // check that we have not allocated too much memory. note that bool
385  // vectors allocate their memory in chunks of whole integers, so they
386  // may over-allocate by up to as many elements as an integer has bits
388  == line_orientations.size(),
390  line_orientations.size()));
392 
393  }
394 
395 
396  template <typename G>
397  void
399  {
400  cells.clear();
401  children.clear();
402  refinement_cases.clear();
403  used.clear();
404  user_flags.clear();
405  boundary_or_material_id.clear();
406  manifold_id.clear();
407  user_data.clear();
408  user_data_type = data_unknown;
409  }
410 
411 
412  void
414  {
416  face_orientations.clear();
417  face_flips.clear();
418  face_rotations.clear();
419  }
420 
421 
422  void
424  {
426  line_orientations.clear();
427  }
428 
429 
430  template<typename G>
431  std::size_t
433  {
434  return (MemoryConsumption::memory_consumption (cells) +
438  MemoryConsumption::memory_consumption (boundary_or_material_id) +
440  MemoryConsumption::memory_consumption (refinement_cases) +
441  user_data.capacity() * sizeof(UserData) + sizeof(user_data));
442  }
443 
444 
445  std::size_t
447  {
452  }
453 
454 
455  std::size_t
457  {
460  }
461 
462 
463 
464 // explicit instantiations
465  template class TriaObjects<TriaObject<1> >;
466  template class TriaObjects<TriaObject<2> >;
467 
468 #include "tria_objects.inst"
469  }
470 }
471 
472 DEAL_II_NAMESPACE_CLOSE
473 
const types::manifold_id flat_manifold_id
Definition: types.h:234
static ::ExceptionBase & ExcMemoryInexact(int arg1, int arg2)
void reserve_space(const unsigned int new_quads_in_pairs, const unsigned int new_quads_single=0)
#define Assert(cond, exc)
Definition: exceptions.h:313
std::vector< RefinementCase< TriaObject< 3 > ::dimension > > refinement_cases
Definition: tria_objects.h:94
void reserve_space(const unsigned int new_objs)
void monitor_memory(const unsigned int true_dimension) const
std_cxx11::enable_if< std_cxx11::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
std::vector< BoundaryOrMaterialId > boundary_or_material_id
Definition: tria_objects.h:165
void monitor_memory(const unsigned int true_dimension) const
std::size_t memory_consumption() const
void reserve_space(const unsigned int new_objs_in_pairs, const unsigned int new_objs_single=0)
Definition: tria_objects.cc:36
static ::ExceptionBase & ExcInternalError()