Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
dof_level.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2003 - 2020 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 
23 
24 #include <vector>
25 
26 
28 
29 // Forward declarations
30 #ifndef DOXYGEN
31 namespace hp
32 {
33  template <int, int>
34  class DoFHandler;
35  template <int, int>
36  class FECollection;
37 } // namespace hp
38 
39 
40 namespace internal
41 {
42  namespace hp
43  {
44  namespace DoFHandlerImplementation
45  {
46  struct Implementation;
47  }
48  } // namespace hp
49  namespace DoFCellAccessorImplementation
50  {
51  struct Implementation;
52  }
53 } // namespace internal
54 #endif
55 
56 
57 namespace internal
58 {
59  namespace hp
60  {
109  class DoFLevel
110  {
111  private:
115  using offset_type = unsigned int;
116 
120  using active_fe_index_type = unsigned short int;
121 
126  using signed_active_fe_index_type = signed short int;
127 
133 
139  static bool
141 
148  static active_fe_index_type
150 
161  std::vector<active_fe_index_type> active_fe_indices;
162 
175  std::vector<active_fe_index_type> future_fe_indices;
176 
189  std::vector<offset_type> dof_offsets;
190 
196  std::vector<types::global_dof_index> dof_indices;
197 
201  std::vector<offset_type> cell_cache_offsets;
202 
208  std::vector<types::global_dof_index> cell_dof_indices_cache;
209 
210  public:
223  void
224  set_dof_index(const unsigned int obj_index,
225  const unsigned int fe_index,
226  const unsigned int local_index,
228 
241  get_dof_index(const unsigned int obj_index,
242  const unsigned int fe_index,
243  const unsigned int local_index) const;
244 
248  unsigned int
249  active_fe_index(const unsigned int obj_index) const;
250 
255  bool
256  fe_index_is_active(const unsigned int obj_index,
257  const unsigned int fe_index) const;
258 
262  void
263  set_active_fe_index(const unsigned int obj_index,
264  const unsigned int fe_index);
265 
270  unsigned int
271  future_fe_index(const unsigned int obj_index) const;
272 
276  void
277  set_future_fe_index(const unsigned int obj_index,
278  const unsigned int fe_index);
279 
283  bool
284  future_fe_index_set(const unsigned int obj_index) const;
285 
289  void
290  clear_future_fe_index(const unsigned int obj_index);
291 
304  get_cell_cache_start(const unsigned int obj_index,
305  const unsigned int dofs_per_cell) const;
306 
311  std::size_t
312  memory_consumption() const;
313 
318  template <class Archive>
319  void
320  serialize(Archive &ar, const unsigned int version);
321 
322  private:
331  template <int dim, int spacedim>
332  void
334  const ::hp::FECollection<dim, spacedim> &fe_collection);
335 
344  template <int dim, int spacedim>
345  void
347  const ::hp::FECollection<dim, spacedim> &fe_collection);
348 
349 
367  void
369 
370 
371  // Make hp::DoFHandler and its auxiliary class a friend since it is the
372  // class that needs to create these data structures.
373  template <int, int>
374  friend class ::hp::DoFHandler;
375  friend struct ::internal::hp::DoFHandlerImplementation::
376  Implementation;
377  friend struct ::internal::DoFCellAccessorImplementation::
378  Implementation;
379  };
380 
381 
382  // -------------------- template functions --------------------------------
383 
384 
385  inline bool
387  {
388  return (static_cast<signed_active_fe_index_type>(active_fe_index) < 0);
389  }
390 
391 
392 
395  const active_fe_index_type active_fe_index)
396  {
397  // convert the active_fe_index into a signed type, flip all
398  // bits, and get the unsigned representation back
399  return static_cast<active_fe_index_type>(
401  }
402 
403 
404 
406  DoFLevel::get_dof_index(const unsigned int obj_index,
407  const unsigned int fe_index,
408  const unsigned int local_index) const
409  {
410  (void)fe_index;
411  AssertIndexRange(obj_index, 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,
441  {
442  (void)fe_index;
443  AssertIndexRange(obj_index, dof_offsets.size());
444 
445  // make sure we are on an
446  // object for which DoFs have
447  // been allocated at all
448  Assert(dof_offsets[obj_index] != static_cast<offset_type>(-1),
449  ExcMessage("You are trying to access degree of freedom "
450  "information for an object on which no such "
451  "information is available"));
452  Assert(
453  is_compressed_entry(active_fe_indices[obj_index]) == false,
454  ExcMessage(
455  "This function can no longer be called after compressing the dof_indices array"));
456  Assert(fe_index == active_fe_indices[obj_index],
457  ExcMessage("FE index does not match that of the present cell"));
458  dof_indices[dof_offsets[obj_index] + local_index] = global_index;
459  }
460 
461 
462 
463  inline unsigned int
464  DoFLevel::active_fe_index(const unsigned int obj_index) const
465  {
466  AssertIndexRange(obj_index, active_fe_indices.size());
467 
468  if (is_compressed_entry(active_fe_indices[obj_index]) == false)
469  return active_fe_indices[obj_index];
470  else
472  }
473 
474 
475 
476  inline bool
477  DoFLevel::fe_index_is_active(const unsigned int obj_index,
478  const unsigned int fe_index) const
479  {
480  return (fe_index == active_fe_index(obj_index));
481  }
482 
483 
484 
485  inline void
486  DoFLevel::set_active_fe_index(const unsigned int obj_index,
487  const unsigned int fe_index)
488  {
489  AssertIndexRange(obj_index, active_fe_indices.size());
490 
491  // check whether the given fe_index is within the range of
492  // values that we interpret as "not compressed". if not, then
493  // the index is so large that we cannot accept it. (but this
494  // will not likely happen because it requires someone using an
495  // FECollection that has more than 32k entries.)
496  Assert(is_compressed_entry(fe_index) == false,
497  ExcMessage(
498  "You are using an active_fe_index that is larger than an "
499  "internal limitation for these objects. Try to work with "
500  "hp::FECollection objects that have a more modest size."));
501  Assert(fe_index != invalid_active_fe_index,
502  ExcMessage(
503  "You are using an active_fe_index that is reserved for "
504  "internal purposes for these objects. Try to work with "
505  "hp::FECollection objects that have a more modest size."));
506 
507  active_fe_indices[obj_index] = fe_index;
508  }
509 
510 
511 
512  inline unsigned int
513  DoFLevel::future_fe_index(const unsigned int obj_index) const
514  {
515  AssertIndexRange(obj_index, future_fe_indices.size());
516 
517  if (future_fe_index_set(obj_index))
518  return future_fe_indices[obj_index];
519 
520  return active_fe_index(obj_index);
521  }
522 
523 
524 
525  inline void
526  DoFLevel::set_future_fe_index(const unsigned int obj_index,
527  const unsigned int fe_index)
528  {
529  AssertIndexRange(obj_index, future_fe_indices.size());
530 
531  // check whether the given fe_index is within the range of
532  // values that we interpret as "not compressed". if not, then
533  // the index is so large that we cannot accept it. (but this
534  // will not likely happen because it requires someone using an
535  // FECollection that has more than 32k entries.)
536  Assert(is_compressed_entry(fe_index) == false,
537  ExcMessage(
538  "You are using a future_fe_index that is larger than an "
539  "internal limitation for these objects. Try to work with "
540  "hp::FECollection objects that have a more modest size."));
541  Assert(fe_index != invalid_active_fe_index,
542  ExcMessage(
543  "You are using a future_fe_index that is reserved for "
544  "internal purposes for these objects. Try to work with "
545  "hp::FECollection objects that have a more modest size."));
546 
547  future_fe_indices[obj_index] = fe_index;
548  }
549 
550 
551 
552  inline bool
553  DoFLevel::future_fe_index_set(const unsigned int obj_index) const
554  {
555  AssertIndexRange(obj_index, future_fe_indices.size());
556 
557  return (future_fe_indices[obj_index] != invalid_active_fe_index);
558  }
559 
560 
561 
562  inline void
563  DoFLevel::clear_future_fe_index(const unsigned int obj_index)
564  {
565  AssertIndexRange(obj_index, future_fe_indices.size());
566 
568  }
569 
570 
571 
572  inline const types::global_dof_index *
573  DoFLevel::get_cell_cache_start(const unsigned int obj_index,
574  const unsigned int dofs_per_cell) const
575  {
576  (void)dofs_per_cell;
577  Assert(
578  (obj_index < cell_cache_offsets.size()) &&
579  (cell_cache_offsets[obj_index] + dofs_per_cell <=
580  cell_dof_indices_cache.size()),
581  ExcMessage(
582  "You are trying to access an element of the cache that stores "
583  "the indices of all degrees of freedom that live on one cell. "
584  "However, this element does not exist. Did you forget to call "
585  "DoFHandler::distribute_dofs(), or did you forget to call it "
586  "again after changing the active_fe_index of one of the cells?"));
587 
588  return &cell_dof_indices_cache[cell_cache_offsets[obj_index]];
589  }
590 
591 
592 
593  template <class Archive>
594  inline void
595  DoFLevel::serialize(Archive &ar, const unsigned int)
596  {
597  ar & this->active_fe_indices;
598  ar & this->cell_cache_offsets;
599  ar & this->cell_dof_indices_cache;
600  ar & this->dof_indices;
601  ar & this->dof_offsets;
602  ar & this->future_fe_indices;
603  }
604  } // namespace hp
605 
606 } // namespace internal
607 
609 
610 #endif
internal::hp::DoFLevel::dof_indices
std::vector< types::global_dof_index > dof_indices
Definition: dof_level.h:196
internal::hp::DoFLevel::future_fe_index
unsigned int future_fe_index(const unsigned int obj_index) const
Definition: dof_level.h:513
internal::hp::DoFLevel::dof_offsets
std::vector< offset_type > dof_offsets
Definition: dof_level.h:189
internal::hp::DoFLevel::set_dof_index
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
internal::hp::DoFLevel::normalize_active_fe_indices
void normalize_active_fe_indices()
Definition: dof_level.cc:263
internal::hp::DoFLevel::signed_active_fe_index_type
signed short int signed_active_fe_index_type
Definition: dof_level.h:126
internal::hp::DoFLevel::get_dof_index
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:406
internal::hp::DoFLevel::fe_index_is_active
bool fe_index_is_active(const unsigned int obj_index, const unsigned int fe_index) const
Definition: dof_level.h:477
internal::hp::DoFLevel::cell_cache_offsets
std::vector< offset_type > cell_cache_offsets
Definition: dof_level.h:201
internal::hp::DoFLevel::set_future_fe_index
void set_future_fe_index(const unsigned int obj_index, const unsigned int fe_index)
Definition: dof_level.h:526
AssertIndexRange
#define AssertIndexRange(index, range)
Definition: exceptions.h:1649
internal::hp::DoFLevel::invalid_active_fe_index
static const active_fe_index_type invalid_active_fe_index
Definition: dof_level.h:132
internal::hp::DoFLevel::is_compressed_entry
static bool is_compressed_entry(const active_fe_index_type active_fe_index)
Definition: dof_level.h:386
DoFHandler
Definition: dof_handler.h:205
internal::hp::DoFLevel
Definition: dof_level.h:109
internal::hp::DoFLevel::get_cell_cache_start
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:573
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
internal::hp::DoFLevel::active_fe_index
unsigned int active_fe_index(const unsigned int obj_index) const
Definition: dof_level.h:464
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
internal::hp::DoFLevel::serialize
void serialize(Archive &ar, const unsigned int version)
Definition: dof_level.h:595
internal::hp::DoFLevel::future_fe_indices
std::vector< active_fe_index_type > future_fe_indices
Definition: dof_level.h:175
internal::hp::DoFLevel::get_toggled_compression_state
static active_fe_index_type get_toggled_compression_state(const active_fe_index_type active_fe_index)
Definition: dof_level.h:394
internal::hp::DoFLevel::set_active_fe_index
void set_active_fe_index(const unsigned int obj_index, const unsigned int fe_index)
Definition: dof_level.h:486
internal::hp::DoFLevel::memory_consumption
std::size_t memory_consumption() const
Definition: dof_level.cc:250
exceptions.h
internal::hp::DoFLevel::active_fe_indices
std::vector< active_fe_index_type > active_fe_indices
Definition: dof_level.h:161
internal::hp::DoFLevel::cell_dof_indices_cache
std::vector< types::global_dof_index > cell_dof_indices_cache
Definition: dof_level.h:208
unsigned int
internal::hp::DoFLevel::future_fe_index_set
bool future_fe_index_set(const unsigned int obj_index) const
Definition: dof_level.h:553
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
internal::hp::DoFLevel::uncompress_data
void uncompress_data(const ::hp::FECollection< dim, spacedim > &fe_collection)
Definition: dof_level.cc:167
hp
Definition: hp.h:117
internal::hp::DoFLevel::active_fe_index_type
unsigned short int active_fe_index_type
Definition: dof_level.h:120
config.h
internal
Definition: aligned_vector.h:369
internal::hp::DoFLevel::clear_future_fe_index
void clear_future_fe_index(const unsigned int obj_index)
Definition: dof_level.h:563
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
TrilinosWrappers::global_index
TrilinosWrappers::types::int_type global_index(const Epetra_BlockMap &map, const ::types::global_dof_index i)
Definition: trilinos_index_access.h:82
int
internal::hp::DoFLevel::compress_data
void compress_data(const ::hp::FECollection< dim, spacedim > &fe_collection)
Definition: dof_level.cc:36