Reference documentation for deal.II version 9.1.0-pre
dof_level.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2003 - 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.md at
12 // the top level directory of deal.II.
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 
22 #include <deal.II/base/exceptions.h>
23 
24 #include <vector>
25 
26 
27 DEAL_II_NAMESPACE_OPEN
28 
29 namespace hp
30 {
31  template <int, int>
32  class DoFHandler;
33  template <int, int>
34  class FECollection;
35 } // namespace hp
36 
37 
38 namespace internal
39 {
40  namespace hp
41  {
42  namespace DoFHandlerImplementation
43  {
44  struct Implementation;
45  }
46  } // namespace hp
47  namespace DoFCellAccessorImplementation
48  {
49  struct Implementation;
50  }
51 } // namespace internal
52 
53 
54 namespace internal
55 {
56  namespace hp
57  {
106  class DoFLevel
107  {
108  private:
112  using offset_type = unsigned int;
113 
117  using active_fe_index_type = unsigned short int;
118 
123  using signed_active_fe_index_type = signed short int;
124 
130  static bool
131  is_compressed_entry(const active_fe_index_type active_fe_index);
132 
139  static active_fe_index_type
140  get_toggled_compression_state(const active_fe_index_type active_fe_index);
141 
152  std::vector<active_fe_index_type> active_fe_indices;
153 
166  std::vector<offset_type> dof_offsets;
167 
173  std::vector<types::global_dof_index> dof_indices;
174 
178  std::vector<offset_type> cell_cache_offsets;
179 
185  std::vector<types::global_dof_index> cell_dof_indices_cache;
186 
187  public:
200  void
201  set_dof_index(const unsigned int obj_index,
202  const unsigned int fe_index,
203  const unsigned int local_index,
204  const types::global_dof_index global_index);
205 
218  get_dof_index(const unsigned int obj_index,
219  const unsigned int fe_index,
220  const unsigned int local_index) const;
221 
225  unsigned int
226  active_fe_index(const unsigned int obj_index) const;
227 
232  bool
233  fe_index_is_active(const unsigned int obj_index,
234  const unsigned int fe_index) const;
235 
239  void
240  set_active_fe_index(const unsigned int obj_index,
241  const unsigned int fe_index);
242 
255  get_cell_cache_start(const unsigned int obj_index,
256  const unsigned int dofs_per_cell) const;
257 
262  std::size_t
263  memory_consumption() const;
264 
269  template <class Archive>
270  void
271  serialize(Archive &ar, const unsigned int version);
272 
273  private:
282  template <int dim, int spacedim>
283  void
284  compress_data(
285  const ::hp::FECollection<dim, spacedim> &fe_collection);
286 
295  template <int dim, int spacedim>
296  void
297  uncompress_data(
298  const ::hp::FECollection<dim, spacedim> &fe_collection);
299 
300 
318  void
319  normalize_active_fe_indices();
320 
321 
326  template <int, int>
327  friend class ::hp::DoFHandler;
328  friend struct ::internal::hp::DoFHandlerImplementation::
329  Implementation;
330  friend struct ::internal::DoFCellAccessorImplementation::
331  Implementation;
332  };
333 
334 
335  // -------------------- template functions --------------------------------
336 
337 
338  inline bool
339  DoFLevel::is_compressed_entry(const active_fe_index_type active_fe_index)
340  {
341  return ((signed_active_fe_index_type)active_fe_index < 0);
342  }
343 
344 
345 
347  DoFLevel::get_toggled_compression_state(
348  const active_fe_index_type active_fe_index)
349  {
350  // convert the active_fe_index into a signed type, flip all
351  // bits, and get the unsigned representation back
352  return (active_fe_index_type) ~(
353  signed_active_fe_index_type)active_fe_index;
354  }
355 
356 
357 
359  DoFLevel::get_dof_index(const unsigned int obj_index,
360  const unsigned int fe_index,
361  const unsigned int local_index) const
362  {
363  (void)fe_index;
364  Assert(obj_index < dof_offsets.size(),
365  ExcIndexRange(obj_index, 0, dof_offsets.size()));
366 
367  // make sure we are on an object for which DoFs have been
368  // allocated at all
369  Assert(dof_offsets[obj_index] != (offset_type)(-1),
370  ExcMessage("You are trying to access degree of freedom "
371  "information for an object on which no such "
372  "information is available"));
373 
374  Assert(fe_index ==
375  (is_compressed_entry(active_fe_indices[obj_index]) == false ?
376  active_fe_indices[obj_index] :
377  get_toggled_compression_state(active_fe_indices[obj_index])),
378  ExcMessage("FE index does not match that of the present cell"));
379 
380  // see if the dof_indices array has been compressed for this
381  // particular cell
382  if (is_compressed_entry(active_fe_indices[obj_index]) == false)
383  return dof_indices[dof_offsets[obj_index] + local_index];
384  else
385  return dof_indices[dof_offsets[obj_index]] + local_index;
386  }
387 
388 
389 
390  inline void
391  DoFLevel::set_dof_index(const unsigned int obj_index,
392  const unsigned int fe_index,
393  const unsigned int local_index,
394  const types::global_dof_index global_index)
395  {
396  (void)fe_index;
397  Assert(obj_index < dof_offsets.size(),
398  ExcIndexRange(obj_index, 0, dof_offsets.size()));
399 
400  // make sure we are on an
401  // object for which DoFs have
402  // been allocated at all
403  Assert(dof_offsets[obj_index] != (offset_type)(-1),
404  ExcMessage("You are trying to access degree of freedom "
405  "information for an object on which no such "
406  "information is available"));
407  Assert(
408  is_compressed_entry(active_fe_indices[obj_index]) == false,
409  ExcMessage(
410  "This function can no longer be called after compressing the dof_indices array"));
411  Assert(fe_index == active_fe_indices[obj_index],
412  ExcMessage("FE index does not match that of the present cell"));
413  dof_indices[dof_offsets[obj_index] + local_index] = global_index;
414  }
415 
416 
417 
418  inline unsigned int
419  DoFLevel::active_fe_index(const unsigned int obj_index) const
420  {
421  Assert(obj_index < active_fe_indices.size(),
422  ExcIndexRange(obj_index, 0, active_fe_indices.size()));
423 
424  if (is_compressed_entry(active_fe_indices[obj_index]) == false)
425  return active_fe_indices[obj_index];
426  else
427  return get_toggled_compression_state(active_fe_indices[obj_index]);
428  }
429 
430 
431 
432  inline bool
433  DoFLevel::fe_index_is_active(const unsigned int obj_index,
434  const unsigned int fe_index) const
435  {
436  return (fe_index == active_fe_index(obj_index));
437  }
438 
439 
440 
441  inline void
442  DoFLevel::set_active_fe_index(const unsigned int obj_index,
443  const unsigned int fe_index)
444  {
445  Assert(obj_index < active_fe_indices.size(),
446  ExcIndexRange(obj_index, 0, active_fe_indices.size()));
447 
448  // check whether the given fe_index is within the range of
449  // values that we interpret as "not compressed". if not, then
450  // the index is so large that we cannot accept it. (but this
451  // will not likely happen because it requires someone using an
452  // FECollection that has more than 32k entries.)
453  Assert(is_compressed_entry(fe_index) == false,
454  ExcMessage(
455  "You are using an active_fe_index that is larger than an "
456  "internal limitation for these objects. Try to work with "
457  "hp::FECollection objects that have a more modest size."));
458 
459  active_fe_indices[obj_index] = fe_index;
460  }
461 
462 
463 
464  inline const types::global_dof_index *
465  DoFLevel::get_cell_cache_start(const unsigned int obj_index,
466  const unsigned int dofs_per_cell) const
467  {
468  (void)dofs_per_cell;
469  Assert(
470  (obj_index < cell_cache_offsets.size()) &&
471  (cell_cache_offsets[obj_index] + dofs_per_cell <=
472  cell_dof_indices_cache.size()),
473  ExcMessage(
474  "You are trying to access an element of the cache that stores "
475  "the indices of all degrees of freedom that live on one cell. "
476  "However, this element does not exist. Did you forget to call "
477  "DoFHandler::distribute_dofs(), or did you forget to call it "
478  "again after changing the active_fe_index of one of the cells?"));
479 
480  return &cell_dof_indices_cache[cell_cache_offsets[obj_index]];
481  }
482 
483 
484 
485  template <class Archive>
486  inline void
487  DoFLevel::serialize(Archive &ar, const unsigned int)
488  {
489  ar & this->active_fe_indices;
490  ar & this->cell_cache_offsets;
491  ar & this->cell_dof_indices_cache;
492  ar & this->dof_indices;
493  ar & this->dof_offsets;
494  }
495  } // namespace hp
496 
497 } // namespace internal
498 
499 DEAL_II_NAMESPACE_CLOSE
500 
501 #endif
unsigned int offset_type
Definition: dof_level.h:112
std::vector< offset_type > cell_cache_offsets
Definition: dof_level.h:178
static::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
unsigned long long int global_dof_index
Definition: types.h:72
unsigned short int active_fe_index_type
Definition: dof_level.h:117
static::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
Definition: exceptions.h:1227
Definition: hp.h:102
std::vector< offset_type > dof_offsets
Definition: dof_level.h:166
std::vector< types::global_dof_index > dof_indices
Definition: dof_level.h:173
signed short int signed_active_fe_index_type
Definition: dof_level.h:123
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