Reference documentation for deal.II version 9.0.0
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 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 DoFHandlerImplementation
40  {
41  struct Implementation;
42  }
43  }
44  namespace DoFCellAccessorImplementation
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 
127  static
128  bool
130 
137  static
140 
151  std::vector<active_fe_index_type> active_fe_indices;
152 
165  std::vector<offset_type> dof_offsets;
166 
172  std::vector<types::global_dof_index> dof_indices;
173 
177  std::vector<offset_type> cell_cache_offsets;
178 
184  std::vector<types::global_dof_index> cell_dof_indices_cache;
185 
186  public:
187 
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 memory_consumption () const;
263 
268  template <class Archive>
269  void serialize(Archive &ar,
270  const unsigned int version);
271 
272  private:
281  template <int dim, int spacedim>
282  void compress_data (const ::hp::FECollection<dim,spacedim> &fe_collection);
283 
292  template <int dim, int spacedim>
293  void uncompress_data (const ::hp::FECollection<dim,spacedim> &fe_collection);
294 
295 
314 
315 
320  template <int, int> friend class ::hp::DoFHandler;
321  friend struct ::internal::hp::DoFHandlerImplementation::Implementation;
322  friend struct ::internal::DoFCellAccessorImplementation::Implementation;
323  };
324 
325 
326  // -------------------- template functions --------------------------------
327 
328 
329  inline
331  {
333  }
334 
335 
336 
337  inline
340  {
341  // convert the active_fe_index into a signed type, flip all
342  // bits, and get the unsigned representation back
344  }
345 
346 
347 
348  inline
350  DoFLevel::
351  get_dof_index (const unsigned int obj_index,
352  const unsigned int fe_index,
353  const unsigned int local_index) const
354  {
355  (void)fe_index;
356  Assert (obj_index < dof_offsets.size(),
357  ExcIndexRange (obj_index, 0, dof_offsets.size()));
358 
359  // make sure we are on an object for which DoFs have been
360  // allocated at all
361  Assert (dof_offsets[obj_index] != (offset_type)(-1),
362  ExcMessage ("You are trying to access degree of freedom "
363  "information for an object on which no such "
364  "information is available"));
365 
366  Assert (fe_index == (is_compressed_entry(active_fe_indices[obj_index]) == false ?
367  active_fe_indices[obj_index] :
369  ExcMessage ("FE index does not match that of the present cell"));
370 
371  // see if the dof_indices array has been compressed for this
372  // particular cell
373  if (is_compressed_entry(active_fe_indices[obj_index]) == false)
374  return dof_indices[dof_offsets[obj_index]+local_index];
375  else
376  return dof_indices[dof_offsets[obj_index]]+local_index;
377  }
378 
379 
380 
381  inline
382  void
383  DoFLevel::
384  set_dof_index (const unsigned int obj_index,
385  const unsigned int fe_index,
386  const unsigned int local_index,
387  const types::global_dof_index global_index)
388  {
389  (void)fe_index;
390  Assert (obj_index < dof_offsets.size(),
391  ExcIndexRange (obj_index, 0, dof_offsets.size()));
392 
393  // make sure we are on an
394  // object for which DoFs have
395  // been allocated at all
396  Assert (dof_offsets[obj_index] != (offset_type)(-1),
397  ExcMessage ("You are trying to access degree of freedom "
398  "information for an object on which no such "
399  "information is available"));
400  Assert (is_compressed_entry(active_fe_indices[obj_index]) == false,
401  ExcMessage ("This function can no longer be called after compressing the dof_indices array"));
402  Assert (fe_index == active_fe_indices[obj_index],
403  ExcMessage ("FE index does not match that of the present cell"));
404  dof_indices[dof_offsets[obj_index]+local_index] = global_index;
405  }
406 
407 
408 
409  inline
410  unsigned int
411  DoFLevel::
412  active_fe_index (const unsigned int obj_index) const
413  {
414  Assert (obj_index < active_fe_indices.size(),
415  ExcIndexRange (obj_index, 0, active_fe_indices.size()));
416 
417  if (is_compressed_entry(active_fe_indices[obj_index]) == false)
418  return active_fe_indices[obj_index];
419  else
421  }
422 
423 
424 
425  inline
426  bool
427  DoFLevel::
428  fe_index_is_active (const unsigned int obj_index,
429  const unsigned int fe_index) const
430  {
431  return (fe_index == active_fe_index(obj_index));
432  }
433 
434 
435 
436  inline
437  void
438  DoFLevel::
439  set_active_fe_index (const unsigned int obj_index,
440  const unsigned int fe_index)
441  {
442  Assert (obj_index < active_fe_indices.size(),
443  ExcIndexRange (obj_index, 0, active_fe_indices.size()));
444 
445  // check whether the given fe_index is within the range of
446  // values that we interpret as "not compressed". if not, then
447  // the index is so large that we cannot accept it. (but this
448  // will not likely happen because it requires someone using an
449  // FECollection that has more than 32k entries.)
450  Assert (is_compressed_entry (fe_index) == false,
451  ExcMessage ("You are using an active_fe_index that is larger than an "
452  "internal limitation for these objects. Try to work with "
453  "hp::FECollection objects that have a more modest size."));
454 
455  active_fe_indices[obj_index] = fe_index;
456  }
457 
458 
459 
460  inline
462  DoFLevel::get_cell_cache_start (const unsigned int obj_index,
463  const unsigned int dofs_per_cell) const
464  {
465  (void)dofs_per_cell;
466  Assert ((obj_index < cell_cache_offsets.size())
467  &&
468  (cell_cache_offsets[obj_index]+dofs_per_cell
469  <=
470  cell_dof_indices_cache.size()),
471  ExcMessage("You are trying to access an element of the cache that stores "
472  "the indices of all degrees of freedom that live on one cell. "
473  "However, this element does not exist. Did you forget to call "
474  "DoFHandler::distribute_dofs(), or did you forget to call it "
475  "again after changing the active_fe_index of one of the cells?"));
476 
477  return &cell_dof_indices_cache[cell_cache_offsets[obj_index]];
478  }
479 
480 
481 
482  template <class Archive>
483  inline
484  void
485  DoFLevel::serialize(Archive &ar,
486  const unsigned int)
487  {
488  ar &this->active_fe_indices;
489  ar &this->cell_cache_offsets;
490  ar &this->cell_dof_indices_cache;
491  ar &this->dof_indices;
492  ar &this->dof_offsets;
493  }
494  } // namespace hp
495 
496 } // namespace internal
497 
498 DEAL_II_NAMESPACE_CLOSE
499 
500 #endif
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
void set_active_fe_index(const unsigned int obj_index, const unsigned int fe_index)
Definition: dof_level.h:439
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
const types::global_dof_index * get_cell_cache_start(const unsigned int obj_index, const unsigned int dofs_per_cell) const
Definition: dof_level.h:462
static ::ExceptionBase & ExcMessage(std::string arg1)
unsigned short int active_fe_index_type
Definition: dof_level.h:114
void set_dof_index(const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index, const types::global_dof_index global_index)
Definition: dof_level.h:384
unsigned int global_dof_index
Definition: types.h:88
#define Assert(cond, exc)
Definition: exceptions.h:1142
Definition: hp.h:102
std::vector< offset_type > dof_offsets
Definition: dof_level.h:165
void serialize(Archive &ar, const unsigned int version)
Definition: dof_level.h:485
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
types::global_dof_index get_dof_index(const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index) const
Definition: dof_level.h:351
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
bool fe_index_is_active(const unsigned int obj_index, const unsigned int fe_index) const
Definition: dof_level.h:428
unsigned int active_fe_index(const unsigned int obj_index) const
Definition: dof_level.h:412
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:151