Reference documentation for deal.II version 9.0.0
tria_levels.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 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/grid/tria_levels.h>
18 
19 DEAL_II_NAMESPACE_OPEN
20 
21 
22 namespace internal
23 {
24  namespace TriangulationImplementation
25  {
26  template <int dim>
27  void
28  TriaLevel<dim>::reserve_space (const unsigned int total_cells,
29  const unsigned int dimension,
30  const unsigned int space_dimension)
31  {
32  // we need space for total_cells cells. Maybe we have more already
33  // with those cells which are unused, so only allocate new space if
34  // needed.
35  //
36  // note that all arrays should have equal sizes (checked by
37  // @p{monitor_memory}
38  if (total_cells > refine_flags.size())
39  {
40  refine_flags.reserve (total_cells);
41  refine_flags.insert (refine_flags.end(),
42  total_cells - refine_flags.size(),
44 
45  coarsen_flags.reserve (total_cells);
46  coarsen_flags.insert (coarsen_flags.end(),
47  total_cells - coarsen_flags.size(),
48  false);
49 
50  active_cell_indices.reserve (total_cells);
51  active_cell_indices.insert (active_cell_indices.end(),
52  total_cells - active_cell_indices.size(),
54 
55  subdomain_ids.reserve (total_cells);
56  subdomain_ids.insert (subdomain_ids.end(),
57  total_cells - subdomain_ids.size(),
58  0);
59 
60  level_subdomain_ids.reserve (total_cells);
61  level_subdomain_ids.insert (level_subdomain_ids.end(),
62  total_cells - level_subdomain_ids.size(),
63  0);
64 
65  if (dimension < space_dimension)
66  {
67  direction_flags.reserve (total_cells);
68  direction_flags.insert (direction_flags.end(),
69  total_cells - direction_flags.size(),
70  true);
71  }
72  else
73  direction_flags.clear ();
74 
75  parents.reserve ((int) (total_cells + 1) / 2);
76  parents.insert (parents.end (),
77  (total_cells + 1) / 2 - parents.size (),
78  -1);
79 
80  neighbors.reserve (total_cells*(2*dimension));
81  neighbors.insert (neighbors.end(),
82  total_cells*(2*dimension) - neighbors.size(),
83  std::make_pair(-1,-1));
84  };
85  }
86 
87 
88  template <int dim>
89  void
90  TriaLevel<dim>::monitor_memory (const unsigned int true_dimension) const
91  {
92  (void)true_dimension;
93  Assert (2*true_dimension*refine_flags.size() == neighbors.size(),
94  ExcMemoryInexact (refine_flags.size(), neighbors.size()));
95  Assert (2*true_dimension*coarsen_flags.size() == neighbors.size(),
96  ExcMemoryInexact (coarsen_flags.size(), neighbors.size()));
97  }
98 
99 
100  template <int dim>
101  std::size_t
103  {
104  return (MemoryConsumption::memory_consumption (refine_flags) +
105  MemoryConsumption::memory_consumption (coarsen_flags) +
106  MemoryConsumption::memory_consumption (active_cell_indices) +
108  MemoryConsumption::memory_consumption (subdomain_ids) +
109  MemoryConsumption::memory_consumption (level_subdomain_ids) +
111  MemoryConsumption::memory_consumption (direction_flags) +
113  }
114 
115 // This specialization should be only temporary, until the TriaObjects
116 // classes are straightened out.
117 
118  void
119  TriaLevel<3>::reserve_space (const unsigned int total_cells,
120  const unsigned int dimension,
121  const unsigned int space_dimension)
122  {
123  // we need space for total_cells
124  // cells. Maybe we have more already
125  // with those cells which are unused,
126  // so only allocate new space if needed.
127  //
128  // note that all arrays should have equal
129  // sizes (checked by @p{monitor_memory}
130  if (total_cells > refine_flags.size())
131  {
132  refine_flags.reserve (total_cells);
133  refine_flags.insert (refine_flags.end(),
134  total_cells - refine_flags.size(),
136 
137  coarsen_flags.reserve (total_cells);
138  coarsen_flags.insert (coarsen_flags.end(),
139  total_cells - coarsen_flags.size(),
140  false);
141 
142  active_cell_indices.reserve (total_cells);
144  total_cells - active_cell_indices.size(),
146 
147  subdomain_ids.reserve (total_cells);
148  subdomain_ids.insert (subdomain_ids.end(),
149  total_cells - subdomain_ids.size(),
150  0);
151 
152  level_subdomain_ids.reserve (total_cells);
154  total_cells - level_subdomain_ids.size(),
155  0);
156 
157  if (dimension < space_dimension)
158  {
159  direction_flags.reserve (total_cells);
160  direction_flags.insert (direction_flags.end(),
161  total_cells - direction_flags.size(),
162  true);
163  }
164  else
165  direction_flags.clear ();
166 
167  parents.reserve ((int) (total_cells + 1) / 2);
168  parents.insert (parents.end (),
169  (total_cells + 1) / 2 - parents.size (),
170  -1);
171 
172  neighbors.reserve (total_cells*(2*dimension));
173  neighbors.insert (neighbors.end(),
174  total_cells*(2*dimension) - neighbors.size(),
175  std::make_pair(-1,-1));
176  };
177  }
178 
179 
180  void
181  TriaLevel<3>::monitor_memory (const unsigned int true_dimension) const
182  {
183  (void)true_dimension;
184  Assert (2*true_dimension*refine_flags.size() == neighbors.size(),
185  ExcMemoryInexact (refine_flags.size(), neighbors.size()));
186  Assert (2*true_dimension*coarsen_flags.size() == neighbors.size(),
187  ExcMemoryInexact (coarsen_flags.size(), neighbors.size()));
188  }
189 
190 
191  std::size_t
193  {
202  }
203  }
204 }
205 
206 
209 
210 DEAL_II_NAMESPACE_CLOSE
211 
TriaObjects< TriaObject< dim > > cells
Definition: tria_levels.h:154
static const unsigned int invalid_unsigned_int
Definition: types.h:173
void monitor_memory(const unsigned int true_dimension) const
Definition: tria_levels.cc:90
std::vector< types::subdomain_id > level_subdomain_ids
Definition: tria_levels.h:129
void reserve_space(const unsigned int total_cells, const unsigned int dimension, const unsigned int space_dimension)
Definition: tria_levels.cc:28
std::vector< std::uint8_t > refine_flags
Definition: tria_levels.h:67
std::vector< unsigned int > active_cell_indices
Definition: tria_levels.h:81
std::vector< std::pair< int, int > > neighbors
Definition: tria_levels.h:116
std::vector< types::subdomain_id > subdomain_ids
Definition: tria_levels.h:124
#define Assert(cond, exc)
Definition: exceptions.h:1142
static ::ExceptionBase & ExcMemoryInexact(int arg1, int arg2)
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)