Reference documentation for deal.II version 9.1.0-pre
dof_level.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2003 - 2017 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/hp/dof_level.h>
19 #include <deal.II/hp/fe_collection.h>
20 
21 #include <iostream>
22 
23 DEAL_II_NAMESPACE_OPEN
24 
25 namespace internal
26 {
27  namespace hp
28  {
29  template <int dim, int spacedim>
30  void
32  const ::hp::FECollection<dim, spacedim> &fe_collection)
33  {
34  (void)fe_collection;
35 
36  if (dof_offsets.size() == 0 || dof_indices.size() == 0)
37  return;
38 
39  // in a first run through, count how many new slots we need in the
40  // dof_indices array after compression. note that the 'cell'
41  // counter is incremented inside the loop
42  unsigned int new_size = 0;
43  for (unsigned int cell = 0; cell < dof_offsets.size();)
44  // see if this cell is active on the current level
45  if (dof_offsets[cell] != (offset_type)(-1))
46  {
47  // find the next cell active on this level
48  unsigned int next_cell = cell + 1;
49  while ((next_cell < dof_offsets.size()) &&
50  (dof_offsets[next_cell] == (offset_type)(-1)))
51  ++next_cell;
52 
53  const unsigned int next_offset =
54  (next_cell < dof_offsets.size() ? dof_offsets[next_cell] :
55  dof_indices.size());
56 
57  Assert(next_offset - dof_offsets[cell] ==
58  fe_collection[active_fe_indices[cell]]
59  .template n_dofs_per_object<dim>(),
61 
62  // see if the range of dofs for this cell can be compressed and if
63  // so how many slots we have to store for them
64  if (next_offset > dof_offsets[cell])
65  {
66  bool compressible = true;
67  for (unsigned int j = dof_offsets[cell] + 1; j < next_offset;
68  ++j)
69  if (dof_indices[j] != dof_indices[j - 1] + 1)
70  {
71  compressible = false;
72  break;
73  }
74  if (compressible == true)
75  new_size += 1;
76  else
77  new_size += (next_offset - dof_offsets[cell]);
78  }
79 
80  // then move on to the next cell
81  cell = next_cell;
82  }
83  else
84  ++cell;
85 
86  // now allocate the new array and copy into it whatever we need
87  std::vector<types::global_dof_index> new_dof_indices;
88  new_dof_indices.reserve(new_size);
89  std::vector<offset_type> new_dof_offsets(dof_offsets.size(),
90  (offset_type)(-1));
91  for (unsigned int cell = 0; cell < dof_offsets.size();)
92  // see if this cell is active on the current level
93  if (dof_offsets[cell] != (offset_type)(-1))
94  {
95  // find the next cell active on this level
96  unsigned int next_cell = cell + 1;
97  while ((next_cell < dof_offsets.size()) &&
98  (dof_offsets[next_cell] == (offset_type)(-1)))
99  ++next_cell;
100 
101  const unsigned int next_offset =
102  (next_cell < dof_offsets.size() ? dof_offsets[next_cell] :
103  dof_indices.size());
104 
105  Assert(next_offset - dof_offsets[cell] ==
106  fe_collection[active_fe_indices[cell]]
107  .template n_dofs_per_object<dim>(),
108  ExcInternalError());
109 
110  new_dof_offsets[cell] = new_dof_indices.size();
111 
112  // see if the range of dofs for this cell can be compressed and if
113  // so how many slots we have to store for them
114  if (next_offset > dof_offsets[cell])
115  {
116  bool compressible = true;
117  for (unsigned int j = dof_offsets[cell] + 1; j < next_offset;
118  ++j)
119  if (dof_indices[j] != dof_indices[j - 1] + 1)
120  {
121  compressible = false;
122  break;
123  }
124 
125  // if this cell is compressible, then copy the first index and
126  // mark this in the dof_offsets array
127  if (compressible == true)
128  {
129  new_dof_indices.push_back(dof_indices[dof_offsets[cell]]);
130 
131  // make sure that the current active_fe_index indicates
132  // that this entry hasn't been compressed yet
134  false,
135  ExcInternalError());
136 
137  // then mark the compression
138  active_fe_indices[cell] =
140  }
141  else
142  for (unsigned int i = dof_offsets[cell]; i < next_offset; ++i)
143  new_dof_indices.push_back(dof_indices[i]);
144  }
145 
146  // then move on to the next cell
147  cell = next_cell;
148  }
149  else
150  ++cell;
151 
152  // finally swap old and new content
153  Assert(new_dof_indices.size() == new_size, ExcInternalError());
154  dof_indices.swap(new_dof_indices);
155  dof_offsets.swap(new_dof_offsets);
156  }
157 
158 
159 
160  template <int dim, int spacedim>
161  void
163  const ::hp::FECollection<dim, spacedim> &fe_collection)
164  {
165  if (dof_offsets.size() == 0 || dof_indices.size() == 0)
166  return;
167 
168  // in a first run through, count how many new slots we need in the
169  // dof_indices array after uncompression.
170  unsigned int new_size = 0;
171  for (unsigned int cell = 0; cell < dof_offsets.size(); ++cell)
172  if (dof_offsets[cell] != (offset_type)(-1))
173  {
174  // we know now that the slot for this cell is used. extract the
175  // active_fe_index for it and see how many entries we need
176  new_size += fe_collection[active_fe_index(cell)]
177  .template n_dofs_per_object<dim>();
178  }
179 
180  // now allocate the new array and copy into it whatever we need
181  std::vector<types::global_dof_index> new_dof_indices;
182  new_dof_indices.reserve(new_size);
183  std::vector<offset_type> new_dof_offsets(dof_offsets.size(),
184  (offset_type)(-1));
185  for (unsigned int cell = 0; cell < dof_offsets.size();)
186  // see if this cell is active on the current level
187  if (dof_offsets[cell] != (offset_type)(-1))
188  {
189  // find the next cell active on this level
190  unsigned int next_cell = cell + 1;
191  while ((next_cell < dof_offsets.size()) &&
192  (dof_offsets[next_cell] == (offset_type)(-1)))
193  ++next_cell;
194 
195  const unsigned int next_offset =
196  (next_cell < dof_offsets.size() ? dof_offsets[next_cell] :
197  dof_indices.size());
198 
199  // set offset for this cell
200  new_dof_offsets[cell] = new_dof_indices.size();
201 
202  // see if we need to uncompress this set of dofs
203  if (is_compressed_entry(active_fe_indices[cell]) == false)
204  {
205  // apparently not. simply copy them
206  Assert(next_offset - dof_offsets[cell] ==
207  fe_collection[active_fe_indices[cell]]
208  .template n_dofs_per_object<dim>(),
209  ExcInternalError());
210  for (unsigned int i = dof_offsets[cell]; i < next_offset; ++i)
211  new_dof_indices.push_back(dof_indices[i]);
212  }
213  else
214  {
215  // apparently so. uncompress
216  Assert(next_offset - dof_offsets[cell] == 1,
217  ExcInternalError());
218  const unsigned int dofs_per_object =
219  fe_collection[get_toggled_compression_state(
220  active_fe_indices[cell])]
221  .template n_dofs_per_object<dim>();
222  for (unsigned int i = 0; i < dofs_per_object; ++i)
223  new_dof_indices.push_back(dof_indices[dof_offsets[cell]] + i);
224 
225  // then mark the uncompression
226  active_fe_indices[cell] =
228  }
229 
230  // then move on to the next cell
231  cell = next_cell;
232  }
233  else
234  ++cell;
235 
236  // verify correct size, then swap arrays
237  Assert(new_dof_indices.size() == new_size, ExcInternalError());
238  dof_indices.swap(new_dof_indices);
239  dof_offsets.swap(new_dof_offsets);
240  }
241 
242 
243 
244  std::size_t
246  {
252  }
253 
254 
255 
256  void
258  {
259  for (unsigned int i = 0; i < active_fe_indices.size(); ++i)
261  active_fe_indices[i] =
263  }
264 
265 
266 
267  // explicit instantiations
268  template void
269  DoFLevel::compress_data(const ::hp::FECollection<1, 1> &);
270  template void
271  DoFLevel::compress_data(const ::hp::FECollection<1, 2> &);
272  template void
273  DoFLevel::compress_data(const ::hp::FECollection<1, 3> &);
274  template void
275  DoFLevel::compress_data(const ::hp::FECollection<2, 2> &);
276  template void
277  DoFLevel::compress_data(const ::hp::FECollection<2, 3> &);
278  template void
279  DoFLevel::compress_data(const ::hp::FECollection<3, 3> &);
280 
281  template void
282  DoFLevel::uncompress_data(const ::hp::FECollection<1, 1> &);
283  template void
284  DoFLevel::uncompress_data(const ::hp::FECollection<1, 2> &);
285  template void
286  DoFLevel::uncompress_data(const ::hp::FECollection<1, 3> &);
287  template void
288  DoFLevel::uncompress_data(const ::hp::FECollection<2, 2> &);
289  template void
290  DoFLevel::uncompress_data(const ::hp::FECollection<2, 3> &);
291  template void
292  DoFLevel::uncompress_data(const ::hp::FECollection<3, 3> &);
293  } // namespace hp
294 } // namespace internal
295 
296 DEAL_II_NAMESPACE_CLOSE
void uncompress_data(const ::hp::FECollection< dim, spacedim > &fe_collection)
Definition: dof_level.cc:162
static active_fe_index_type get_toggled_compression_state(const active_fe_index_type active_fe_index)
Definition: dof_level.h:347
unsigned int offset_type
Definition: dof_level.h:112
std::vector< offset_type > cell_cache_offsets
Definition: dof_level.h:178
void compress_data(const ::hp::FECollection< dim, spacedim > &fe_collection)
Definition: dof_level.cc:31
#define Assert(cond, exc)
Definition: exceptions.h:1227
Definition: hp.h:102
std::vector< offset_type > dof_offsets
Definition: dof_level.h:166
static bool is_compressed_entry(const active_fe_index_type active_fe_index)
Definition: dof_level.h:339
void normalize_active_fe_indices()
Definition: dof_level.cc:257
std::vector< types::global_dof_index > dof_indices
Definition: dof_level.h:173
std::size_t memory_consumption() const
Definition: dof_level.cc:245
std::vector< types::global_dof_index > cell_dof_indices_cache
Definition: dof_level.h:185
std::vector< active_fe_index_type > active_fe_indices
Definition: dof_level.h:152
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
static::ExceptionBase & ExcInternalError()
unsigned int active_fe_index(const unsigned int obj_index) const
Definition: dof_level.h:419