Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
dof_level.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2003 - 2019 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 
136  static bool
138 
145  static active_fe_index_type
147 
158  std::vector<active_fe_index_type> active_fe_indices;
159 
172  std::vector<active_fe_index_type> future_fe_indices;
173 
186  std::vector<offset_type> dof_offsets;
187 
193  std::vector<types::global_dof_index> dof_indices;
194 
198  std::vector<offset_type> cell_cache_offsets;
199 
205  std::vector<types::global_dof_index> cell_dof_indices_cache;
206 
207  public:
220  void
221  set_dof_index(const unsigned int obj_index,
222  const unsigned int fe_index,
223  const unsigned int local_index,
224  const types::global_dof_index global_index);
225 
238  get_dof_index(const unsigned int obj_index,
239  const unsigned int fe_index,
240  const unsigned int local_index) const;
241 
245  unsigned int
246  active_fe_index(const unsigned int obj_index) const;
247 
252  bool
253  fe_index_is_active(const unsigned int obj_index,
254  const unsigned int fe_index) const;
255 
259  void
260  set_active_fe_index(const unsigned int obj_index,
261  const unsigned int fe_index);
262 
267  unsigned int
268  future_fe_index(const unsigned int obj_index) const;
269 
273  void
274  set_future_fe_index(const unsigned int obj_index,
275  const unsigned int fe_index);
276 
280  bool
281  future_fe_index_set(const unsigned int obj_index) const;
282 
286  void
287  clear_future_fe_index(const unsigned int obj_index);
288 
301  get_cell_cache_start(const unsigned int obj_index,
302  const unsigned int dofs_per_cell) const;
303 
308  std::size_t
309  memory_consumption() const;
310 
315  template <class Archive>
316  void
317  serialize(Archive &ar, const unsigned int version);
318 
319  private:
328  template <int dim, int spacedim>
329  void
331  const ::hp::FECollection<dim, spacedim> &fe_collection);
332 
341  template <int dim, int spacedim>
342  void
344  const ::hp::FECollection<dim, spacedim> &fe_collection);
345 
346 
364  void
366 
367 
372  template <int, int>
373  friend class ::hp::DoFHandler;
374  friend struct ::internal::hp::DoFHandlerImplementation::
375  Implementation;
376  friend struct ::internal::DoFCellAccessorImplementation::
377  Implementation;
378  };
379 
380 
381  // -------------------- template functions --------------------------------
382 
383 
384  inline bool
386  {
387  return (static_cast<signed_active_fe_index_type>(active_fe_index) < 0);
388  }
389 
390 
391 
394  const active_fe_index_type active_fe_index)
395  {
396  // convert the active_fe_index into a signed type, flip all
397  // bits, and get the unsigned representation back
398  return static_cast<active_fe_index_type>(
400  }
401 
402 
403 
405  DoFLevel::get_dof_index(const unsigned int obj_index,
406  const unsigned int fe_index,
407  const unsigned int local_index) const
408  {
409  (void)fe_index;
410  Assert(obj_index < dof_offsets.size(),
411  ExcIndexRange(obj_index, 0, dof_offsets.size()));
412 
413  // make sure we are on an object for which DoFs have been
414  // allocated at all
415  Assert(dof_offsets[obj_index] != static_cast<offset_type>(-1),
416  ExcMessage("You are trying to access degree of freedom "
417  "information for an object on which no such "
418  "information is available"));
419 
420  Assert(fe_index ==
421  (is_compressed_entry(active_fe_indices[obj_index]) == false ?
422  active_fe_indices[obj_index] :
424  ExcMessage("FE index does not match that of the present cell"));
425 
426  // see if the dof_indices array has been compressed for this
427  // particular cell
428  if (is_compressed_entry(active_fe_indices[obj_index]) == false)
429  return dof_indices[dof_offsets[obj_index] + local_index];
430  else
431  return dof_indices[dof_offsets[obj_index]] + local_index;
432  }
433 
434 
435 
436  inline void
437  DoFLevel::set_dof_index(const unsigned int obj_index,
438  const unsigned int fe_index,
439  const unsigned int local_index,
440  const types::global_dof_index global_index)
441  {
442  (void)fe_index;
443  Assert(obj_index < dof_offsets.size(),
444  ExcIndexRange(obj_index, 0, dof_offsets.size()));
445 
446  // make sure we are on an
447  // object for which DoFs have
448  // been allocated at all
449  Assert(dof_offsets[obj_index] != static_cast<offset_type>(-1),
450  ExcMessage("You are trying to access degree of freedom "
451  "information for an object on which no such "
452  "information is available"));
453  Assert(
454  is_compressed_entry(active_fe_indices[obj_index]) == false,
455  ExcMessage(
456  "This function can no longer be called after compressing the dof_indices array"));
457  Assert(fe_index == active_fe_indices[obj_index],
458  ExcMessage("FE index does not match that of the present cell"));
459  dof_indices[dof_offsets[obj_index] + local_index] = global_index;
460  }
461 
462 
463 
464  inline unsigned int
465  DoFLevel::active_fe_index(const unsigned int obj_index) const
466  {
467  Assert(obj_index < active_fe_indices.size(),
468  ExcIndexRange(obj_index, 0, active_fe_indices.size()));
469 
470  if (is_compressed_entry(active_fe_indices[obj_index]) == false)
471  return active_fe_indices[obj_index];
472  else
474  }
475 
476 
477 
478  inline bool
479  DoFLevel::fe_index_is_active(const unsigned int obj_index,
480  const unsigned int fe_index) const
481  {
482  return (fe_index == active_fe_index(obj_index));
483  }
484 
485 
486 
487  inline void
488  DoFLevel::set_active_fe_index(const unsigned int obj_index,
489  const unsigned int fe_index)
490  {
491  Assert(obj_index < active_fe_indices.size(),
492  ExcIndexRange(obj_index, 0, active_fe_indices.size()));
493 
494  // check whether the given fe_index is within the range of
495  // values that we interpret as "not compressed". if not, then
496  // the index is so large that we cannot accept it. (but this
497  // will not likely happen because it requires someone using an
498  // FECollection that has more than 32k entries.)
499  Assert(is_compressed_entry(fe_index) == false,
500  ExcMessage(
501  "You are using an active_fe_index that is larger than an "
502  "internal limitation for these objects. Try to work with "
503  "hp::FECollection objects that have a more modest size."));
504  Assert(fe_index != invalid_active_fe_index,
505  ExcMessage(
506  "You are using an active_fe_index that is reserved for "
507  "internal purposes for these objects. Try to work with "
508  "hp::FECollection objects that have a more modest size."));
509 
510  active_fe_indices[obj_index] = fe_index;
511  }
512 
513 
514 
515  inline unsigned int
516  DoFLevel::future_fe_index(const unsigned int obj_index) const
517  {
518  Assert(obj_index < future_fe_indices.size(),
519  ExcIndexRange(obj_index, 0, future_fe_indices.size()));
520 
521  if (future_fe_index_set(obj_index))
522  return future_fe_indices[obj_index];
523 
524  return active_fe_index(obj_index);
525  }
526 
527 
528 
529  inline void
530  DoFLevel::set_future_fe_index(const unsigned int obj_index,
531  const unsigned int fe_index)
532  {
533  Assert(obj_index < future_fe_indices.size(),
534  ExcIndexRange(obj_index, 0, future_fe_indices.size()));
535 
536  // check whether the given fe_index is within the range of
537  // values that we interpret as "not compressed". if not, then
538  // the index is so large that we cannot accept it. (but this
539  // will not likely happen because it requires someone using an
540  // FECollection that has more than 32k entries.)
541  Assert(is_compressed_entry(fe_index) == false,
542  ExcMessage(
543  "You are using a future_fe_index that is larger than an "
544  "internal limitation for these objects. Try to work with "
545  "hp::FECollection objects that have a more modest size."));
546  Assert(fe_index != invalid_active_fe_index,
547  ExcMessage(
548  "You are using a future_fe_index that is reserved for "
549  "internal purposes for these objects. Try to work with "
550  "hp::FECollection objects that have a more modest size."));
551 
552  future_fe_indices[obj_index] = fe_index;
553  }
554 
555 
556 
557  inline bool
558  DoFLevel::future_fe_index_set(const unsigned int obj_index) const
559  {
560  Assert(obj_index < future_fe_indices.size(),
561  ExcIndexRange(obj_index, 0, future_fe_indices.size()));
562 
563  return (future_fe_indices[obj_index] != invalid_active_fe_index);
564  }
565 
566 
567 
568  inline void
569  DoFLevel::clear_future_fe_index(const unsigned int obj_index)
570  {
571  Assert(obj_index < future_fe_indices.size(),
572  ExcIndexRange(obj_index, 0, future_fe_indices.size()));
573 
575  }
576 
577 
578 
579  inline const types::global_dof_index *
580  DoFLevel::get_cell_cache_start(const unsigned int obj_index,
581  const unsigned int dofs_per_cell) const
582  {
583  (void)dofs_per_cell;
584  Assert(
585  (obj_index < cell_cache_offsets.size()) &&
586  (cell_cache_offsets[obj_index] + dofs_per_cell <=
587  cell_dof_indices_cache.size()),
588  ExcMessage(
589  "You are trying to access an element of the cache that stores "
590  "the indices of all degrees of freedom that live on one cell. "
591  "However, this element does not exist. Did you forget to call "
592  "DoFHandler::distribute_dofs(), or did you forget to call it "
593  "again after changing the active_fe_index of one of the cells?"));
594 
595  return &cell_dof_indices_cache[cell_cache_offsets[obj_index]];
596  }
597 
598 
599 
600  template <class Archive>
601  inline void
602  DoFLevel::serialize(Archive &ar, const unsigned int)
603  {
604  ar & this->active_fe_indices;
605  ar & this->cell_cache_offsets;
606  ar & this->cell_dof_indices_cache;
607  ar & this->dof_indices;
608  ar & this->dof_offsets;
609  ar & this->future_fe_indices;
610  }
611  } // namespace hp
612 
613 } // namespace internal
614 
615 DEAL_II_NAMESPACE_CLOSE
616 
617 #endif
void uncompress_data(const ::hp::FECollection< dim, spacedim > &fe_collection)
Definition: dof_level.cc:167
std::size_t memory_consumption() const
Definition: dof_level.cc:250
static active_fe_index_type get_toggled_compression_state(const active_fe_index_type active_fe_index)
Definition: dof_level.h:393
unsigned int offset_type
Definition: dof_level.h:112
std::vector< offset_type > cell_cache_offsets
Definition: dof_level.h:198
void compress_data(const ::hp::FECollection< dim, spacedim > &fe_collection)
Definition: dof_level.cc:36
void set_active_fe_index(const unsigned int obj_index, const unsigned int fe_index)
Definition: dof_level.h:488
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:580
unsigned short int active_fe_index_type
Definition: dof_level.h:117
bool future_fe_index_set(const unsigned int obj_index) const
Definition: dof_level.h:558
unsigned int future_fe_index(const unsigned int obj_index) const
Definition: dof_level.h:516
void clear_future_fe_index(const unsigned int obj_index)
Definition: dof_level.h:569
static ::ExceptionBase & ExcMessage(std::string arg1)
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:437
#define Assert(cond, exc)
Definition: exceptions.h:1407
Definition: hp.h:117
std::vector< offset_type > dof_offsets
Definition: dof_level.h:186
void serialize(Archive &ar, const unsigned int version)
Definition: dof_level.h:602
static bool is_compressed_entry(const active_fe_index_type active_fe_index)
Definition: dof_level.h:385
void normalize_active_fe_indices()
Definition: dof_level.cc:263
unsigned int global_dof_index
Definition: types.h:89
void set_future_fe_index(const unsigned int obj_index, const unsigned int fe_index)
Definition: dof_level.h:530
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:405
std::vector< types::global_dof_index > dof_indices
Definition: dof_level.h:193
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:205
bool fe_index_is_active(const unsigned int obj_index, const unsigned int fe_index) const
Definition: dof_level.h:479
unsigned int active_fe_index(const unsigned int obj_index) const
Definition: dof_level.h:465
static const active_fe_index_type invalid_active_fe_index
Definition: dof_level.h:129
std::vector< active_fe_index_type > future_fe_indices
Definition: dof_level.h:172
std::vector< active_fe_index_type > active_fe_indices
Definition: dof_level.h:158