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