Reference documentation for deal.II version 8.4.2
dof_level.h
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 #ifndef dealii__hp_dof_level_h
17 #define dealii__hp_dof_level_h
18 
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/base/exceptions.h>
22 
23 #include <vector>
24 
25 
26 DEAL_II_NAMESPACE_OPEN
27 
28 namespace hp
29 {
30  template <int, int> class DoFHandler;
31  template <int, int> class FECollection;
32 }
33 
34 
35 namespace internal
36 {
37  namespace hp
38  {
39  namespace DoFHandler
40  {
41  struct Implementation;
42  }
43  }
44  namespace DoFCellAccessor
45  {
46  struct Implementation;
47  }
48 }
49 
50 
51 namespace internal
52 {
53  namespace hp
54  {
103  class DoFLevel
104  {
105  private:
109  typedef unsigned int offset_type;
110 
114  typedef unsigned short int active_fe_index_type;
115 
120  typedef signed short int signed_active_fe_index_type;
121 
132  std::vector<active_fe_index_type> active_fe_indices;
133 
146  std::vector<offset_type> dof_offsets;
147 
153  std::vector<types::global_dof_index> dof_indices;
154 
158  std::vector<offset_type> cell_cache_offsets;
159 
165  std::vector<types::global_dof_index> cell_dof_indices_cache;
166 
167  public:
168 
181  void
182  set_dof_index (const unsigned int obj_index,
183  const unsigned int fe_index,
184  const unsigned int local_index,
185  const types::global_dof_index global_index);
186 
199  get_dof_index (const unsigned int obj_index,
200  const unsigned int fe_index,
201  const unsigned int local_index) const;
202 
206  unsigned int
207  active_fe_index (const unsigned int obj_index) const;
208 
213  bool
214  fe_index_is_active (const unsigned int obj_index,
215  const unsigned int fe_index) const;
216 
220  void
221  set_active_fe_index (const unsigned int obj_index,
222  const unsigned int fe_index);
223 
236  get_cell_cache_start (const unsigned int obj_index,
237  const unsigned int dofs_per_cell) const;
238 
243  std::size_t memory_consumption () const;
244 
249  template <class Archive>
250  void serialize(Archive &ar,
251  const unsigned int version);
252 
253  private:
262  template <int dim, int spacedim>
263  void compress_data (const ::hp::FECollection<dim,spacedim> &fe_collection);
264 
273  template <int dim, int spacedim>
274  void uncompress_data (const ::hp::FECollection<dim,spacedim> &fe_collection);
275 
280  template <int, int> friend class ::hp::DoFHandler;
281  friend struct ::internal::hp::DoFHandler::Implementation;
282  friend struct ::internal::DoFCellAccessor::Implementation;
283  };
284 
285 
286  // -------------------- template functions --------------------------------
287 
288  inline
290  DoFLevel::
291  get_dof_index (const unsigned int obj_index,
292  const unsigned int fe_index,
293  const unsigned int local_index) const
294  {
295  (void)fe_index;
296  Assert (obj_index < dof_offsets.size(),
297  ExcIndexRange (obj_index, 0, dof_offsets.size()));
298 
299  // make sure we are on an
300  // object for which DoFs have
301  // been allocated at all
302  Assert (dof_offsets[obj_index] != (offset_type)(-1),
303  ExcMessage ("You are trying to access degree of freedom "
304  "information for an object on which no such "
305  "information is available"));
306 
307  Assert (fe_index == active_fe_indices[obj_index],
308  ExcMessage ("FE index does not match that of the present cell"));
309 
310  // see if the dof_indices array has been compressed for this
311  // particular cell
312  if ((signed_active_fe_index_type)active_fe_indices[obj_index]>=0)
313  return dof_indices[dof_offsets[obj_index]+local_index];
314  else
315  return dof_indices[dof_offsets[obj_index]]+local_index;
316  }
317 
318 
319 
320  inline
321  void
322  DoFLevel::
323  set_dof_index (const unsigned int obj_index,
324  const unsigned int fe_index,
325  const unsigned int local_index,
326  const types::global_dof_index global_index)
327  {
328  (void)fe_index;
329  Assert (obj_index < dof_offsets.size(),
330  ExcIndexRange (obj_index, 0, dof_offsets.size()));
331 
332  // make sure we are on an
333  // object for which DoFs have
334  // been allocated at all
335  Assert (dof_offsets[obj_index] != (offset_type)(-1),
336  ExcMessage ("You are trying to access degree of freedom "
337  "information for an object on which no such "
338  "information is available"));
339  Assert ((signed_active_fe_index_type)active_fe_indices[obj_index]>=0,
340  ExcMessage ("This function can no longer be called after compressing the dof_indices array"));
341  Assert (fe_index == active_fe_indices[obj_index],
342  ExcMessage ("FE index does not match that of the present cell"));
343  dof_indices[dof_offsets[obj_index]+local_index] = global_index;
344  }
345 
346 
347 
348  inline
349  unsigned int
350  DoFLevel::
351  active_fe_index (const unsigned int obj_index) const
352  {
353  Assert (obj_index < active_fe_indices.size(),
354  ExcIndexRange (obj_index, 0, active_fe_indices.size()));
355 
356  if (((signed_active_fe_index_type)active_fe_indices[obj_index]) >= 0)
357  return active_fe_indices[obj_index];
358  else
359  return (active_fe_index_type)~(signed_active_fe_index_type)active_fe_indices[obj_index];
360  }
361 
362 
363 
364  inline
365  bool
366  DoFLevel::
367  fe_index_is_active (const unsigned int obj_index,
368  const unsigned int fe_index) const
369  {
370  return (fe_index == active_fe_index(obj_index));
371  }
372 
373 
374  inline
375  void
376  DoFLevel::
377  set_active_fe_index (const unsigned int obj_index,
378  const unsigned int fe_index)
379  {
380  Assert (obj_index < active_fe_indices.size(),
381  ExcIndexRange (obj_index, 0, active_fe_indices.size()));
382 
383  active_fe_indices[obj_index] = fe_index;
384  }
385 
386 
387 
388  inline
390  DoFLevel::get_cell_cache_start (const unsigned int obj_index,
391  const unsigned int dofs_per_cell) const
392  {
393  (void)dofs_per_cell;
394  Assert (obj_index < cell_cache_offsets.size(),
395  ExcInternalError());
396  Assert (cell_cache_offsets[obj_index]+dofs_per_cell
397  <=
398  cell_dof_indices_cache.size(),
399  ExcInternalError());
400 
401  return &cell_dof_indices_cache[cell_cache_offsets[obj_index]];
402  }
403 
404  template <class Archive>
405  inline
406  void
407  DoFLevel::serialize(Archive &ar,
408  const unsigned int)
409  {
410  ar &this->active_fe_indices;
411  ar &this->cell_cache_offsets;
412  ar &this->cell_dof_indices_cache;
413  ar &this->dof_indices;
414  ar &this->dof_offsets;
415  }
416  } // namespace hp
417 
418 } // namespace internal
419 
420 DEAL_II_NAMESPACE_CLOSE
421 
422 #endif
unsigned int offset_type
Definition: dof_level.h:109
::ExceptionBase & ExcMessage(std::string arg1)
std::vector< offset_type > cell_cache_offsets
Definition: dof_level.h:158
unsigned short int active_fe_index_type
Definition: dof_level.h:114
unsigned int global_dof_index
Definition: types.h:88
#define Assert(cond, exc)
Definition: exceptions.h:294
Definition: hp.h:102
std::vector< offset_type > dof_offsets
Definition: dof_level.h:146
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
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