Reference documentation for deal.II version 8.4.2
dof_level.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2003 - 2015 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  return;
33 
34  if (dof_offsets.size() == 0 || dof_indices.size()==0)
35  return;
36 
37  // in a first run through, count how many new slots we need in the
38  // dof_indices array after compression. note that the 'cell'
39  // counter is incremented inside the loop
40  unsigned int new_size = 0;
41  for (unsigned int cell=0; cell<dof_offsets.size(); )
42  // see if this cell is active on the current level
43  if (dof_offsets[cell] != (offset_type)(-1))
44  {
45  // find the next cell active on this level
46  unsigned int next_cell = cell+1;
47  while ((next_cell<dof_offsets.size()) &&
48  (dof_offsets[next_cell] == (offset_type)(-1)))
49  ++next_cell;
50 
51  const unsigned int next_offset = (next_cell < dof_offsets.size() ?
52  dof_offsets[next_cell] :
53  dof_indices.size());
54 
55  Assert (next_offset-dof_offsets[cell] == fe_collection[active_fe_indices[cell]].template n_dofs_per_object<dim>(),
56  ExcInternalError());
57 
58  // see if the range of dofs for this cell can be compressed and if so
59  // how many slots we have to store for them
60  if (next_offset > dof_offsets[cell])
61  {
62  bool compressible = true;
63  for (unsigned int j=dof_offsets[cell]+1; j<next_offset; ++j)
64  if (dof_indices[j] != dof_indices[j-1]+1)
65  {
66  compressible = false;
67  break;
68  }
69  if (compressible == true)
70  new_size += 1;
71  else
72  new_size += (next_offset-dof_offsets[cell]);
73  }
74 
75  // then move on to the next cell
76  cell = next_cell;
77  }
78  else
79  ++cell;
80 
81  // now allocate the new array and copy into it whatever we need
82  std::vector<types::global_dof_index> new_dof_indices;
83  new_dof_indices.reserve(new_size);
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>(),
99  ExcInternalError());
100 
101  // see if the range of dofs for this cell can be compressed and if so
102  // how many slots we have to store for them
103  if (next_offset > dof_offsets[cell])
104  {
105  bool compressible = true;
106  for (unsigned int j=dof_offsets[cell]+1; j<next_offset; ++j)
107  if (dof_indices[j] != dof_indices[j-1]+1)
108  {
109  compressible = false;
110  break;
111  }
112 
113  // if this cell is compressible, then copy the first index and mark this
114  // in the dof_offsets array
115  if (compressible == true)
116  {
117  new_dof_indices.push_back (dof_indices[dof_offsets[cell]]);
118 
119  // make sure that the current active_fe_index indicates
120  // that this entry hasn't been compressed yet
121  Assert ((signed_active_fe_index_type)active_fe_indices[cell] >= 0, ExcInternalError());
122 
123  // then mark the compression
125  }
126  else
127  for (unsigned int i=dof_offsets[cell]; i<next_offset; ++i)
128  new_dof_indices.push_back (dof_indices[i]);
129  }
130 
131  // then move on to the next cell
132  cell = next_cell;
133  }
134  else
135  ++cell;
136 
137  // finally swap old and new content
138  Assert (new_dof_indices.size() == new_size, ExcInternalError());
139  dof_indices.swap (new_dof_indices);
140  }
141 
142 
143 
144  template <int dim, int spacedim>
145  void
146  DoFLevel::uncompress_data(const ::hp::FECollection<dim,spacedim> &fe_collection)
147  {
148  return;
149 
150  if (dof_offsets.size() == 0 || dof_indices.size()==0)
151  return;
152 
153  // in a first run through, count how many new slots we need in the
154  // dof_indices array after uncompression.
155  unsigned int new_size = 0;
156  for (unsigned int cell=0; cell<dof_offsets.size(); ++cell)
157  if (dof_offsets[cell] != (offset_type)(-1))
158  {
159  // we know now that the slot for this cell is used. extract the
160  // active_fe_index for it and see how many entries we need
161  new_size += fe_collection[active_fe_index(cell)].template n_dofs_per_object<dim>();
162  }
163 
164  // now allocate the new array and copy into it whatever we need
165  std::vector<types::global_dof_index> new_dof_indices;
166  new_dof_indices.reserve(new_size);
167  std::vector<offset_type> new_dof_offsets (dof_offsets.size(), (offset_type)(-1));
168  for (unsigned int cell=0; cell<dof_offsets.size(); )
169  // see if this cell is active on the current level
170  if (dof_offsets[cell] != (offset_type)(-1))
171  {
172  // find the next cell active on this level
173  unsigned int next_cell = cell+1;
174  while ((next_cell<dof_offsets.size()) &&
175  (dof_offsets[next_cell] == (offset_type)(-1)))
176  ++next_cell;
177 
178  const unsigned int next_offset = (next_cell < dof_offsets.size() ?
179  dof_offsets[next_cell] :
180  dof_indices.size());
181 
182  // set offset for this cell
183  new_dof_offsets[cell] = new_dof_indices.size();
184 
185  // see if we need to uncompress this set of dofs
187  {
188  // apparently not. simply copy them
189  Assert (next_offset-dof_offsets[cell] == fe_collection[active_fe_indices[cell]].template n_dofs_per_object<dim>(),
190  ExcInternalError());
191  for (unsigned int i=dof_offsets[cell]; i<next_offset; ++i)
192  new_dof_indices.push_back (dof_indices[i]);
193  }
194  else
195  {
196  // apparently so. uncompress
197  Assert (next_offset-dof_offsets[cell] == 1,
198  ExcInternalError());
199  for (unsigned int i=0; i<fe_collection[active_fe_indices[cell]].template n_dofs_per_object<dim>(); ++i)
200  new_dof_indices.push_back (dof_indices[dof_offsets[cell]]+i);
201  }
202 
203  // then move on to the next cell
204  cell = next_cell;
205  }
206  else
207  ++cell;
208 
209  // verify correct size, then swap arrays
210  Assert (new_dof_indices.size() == new_size, ExcInternalError());
211  dof_indices.swap (new_dof_indices);
212  dof_offsets.swap (new_dof_offsets);
213  }
214 
215 
216  std::size_t
218  {
224  }
225 
226 
227  // explicit instantiations
228  template void DoFLevel::compress_data(const ::hp::FECollection<1,1> &);
229  template void DoFLevel::compress_data(const ::hp::FECollection<1,2> &);
230  template void DoFLevel::compress_data(const ::hp::FECollection<1,3> &);
231  template void DoFLevel::compress_data(const ::hp::FECollection<2,2> &);
232  template void DoFLevel::compress_data(const ::hp::FECollection<2,3> &);
233  template void DoFLevel::compress_data(const ::hp::FECollection<3,3> &);
234 
235  template void DoFLevel::uncompress_data(const ::hp::FECollection<1,1> &);
236  template void DoFLevel::uncompress_data(const ::hp::FECollection<1,2> &);
237  template void DoFLevel::uncompress_data(const ::hp::FECollection<1,3> &);
238  template void DoFLevel::uncompress_data(const ::hp::FECollection<2,2> &);
239  template void DoFLevel::uncompress_data(const ::hp::FECollection<2,3> &);
240  template void DoFLevel::uncompress_data(const ::hp::FECollection<3,3> &);
241  }
242 }
243 
244 DEAL_II_NAMESPACE_CLOSE
void uncompress_data(const ::hp::FECollection< dim, spacedim > &fe_collection)
Definition: dof_level.cc:146
unsigned int offset_type
Definition: dof_level.h:109
std::size_t memory_consumption() const
Definition: dof_level.cc:217
std::vector< offset_type > cell_cache_offsets
Definition: dof_level.h:158
void compress_data(const ::hp::FECollection< dim, spacedim > &fe_collection)
Definition: dof_level.cc:29
unsigned short int active_fe_index_type
Definition: dof_level.h:114
#define Assert(cond, exc)
Definition: exceptions.h:294
Definition: hp.h:102
std::vector< offset_type > dof_offsets
Definition: dof_level.h:146
std_cxx11::enable_if< std_cxx11::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
std::vector< types::global_dof_index > dof_indices
Definition: dof_level.h:153
std::vector< types::global_dof_index > cell_dof_indices_cache
Definition: dof_level.h:165
unsigned int active_fe_index(const unsigned int obj_index) const
Definition: dof_level.h:351
signed short int signed_active_fe_index_type
Definition: dof_level.h:120
std::vector< active_fe_index_type > active_fe_indices
Definition: dof_level.h:132