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_faces.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2006 - 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_faces_h
17 #define dealii_hp_dof_faces_h
18 
19 
20 #include <deal.II/base/config.h>
21 
23 
24 #include <vector>
25 
27 
28 namespace internal
29 {
30  namespace hp
31  {
94  template <int structdim>
96  {
97  public:
105  std::vector<unsigned int> dof_offsets;
106 
111  std::vector<types::global_dof_index> dofs;
112 
125  template <int dim, int spacedim>
126  void
127  set_dof_index(const ::hp::DoFHandler<dim, spacedim> &dof_handler,
128  const unsigned int obj_index,
129  const unsigned int fe_index,
130  const unsigned int local_index,
132  const unsigned int obj_level);
133 
145  template <int dim, int spacedim>
147  get_dof_index(const ::hp::DoFHandler<dim, spacedim> &dof_handler,
148  const unsigned int obj_index,
149  const unsigned int fe_index,
150  const unsigned int local_index,
151  const unsigned int obj_level) const;
152 
164  template <int dim, int spacedim>
165  unsigned int
167  const ::hp::DoFHandler<dim, spacedim> &dof_handler,
168  const unsigned int obj_index) const;
169 
173  template <int dim, int spacedim>
176  const ::hp::DoFHandler<dim, spacedim> &dof_handler,
177  const unsigned int obj_level,
178  const unsigned int obj_index,
179  const unsigned int n) const;
180 
185  template <int dim, int spacedim>
186  bool
188  const ::hp::DoFHandler<dim, spacedim> &dof_handler,
189  const unsigned int obj_index,
190  const unsigned int fe_index,
191  const unsigned int obj_level) const;
192 
197  std::size_t
198  memory_consumption() const;
199 
204  template <class Archive>
205  void
206  serialize(Archive &ar, const unsigned int version);
207  };
208 
209 
210 
238  template <int dim>
240 
241 
249  template <>
251  {
252  public:
257  std::size_t
258  memory_consumption() const;
259 
264  template <class Archive>
265  void
266  serialize(Archive &ar, const unsigned int version);
267  };
268 
276  template <>
278  {
279  public:
284 
289  std::size_t
290  memory_consumption() const;
291 
296  template <class Archive>
297  void
298  serialize(Archive &ar, const unsigned int version);
299  };
300 
308  template <>
310  {
311  public:
316 
321 
326  std::size_t
327  memory_consumption() const;
328 
333  template <class Archive>
334  void
335  serialize(Archive &ar, const unsigned int version);
336  };
337 
338 
339  // --------------------- inline and template functions ------------------
340  template <class Archive>
341  void
342  DoFIndicesOnFaces<1>::serialize(Archive &, const unsigned int)
343  {}
344 
345 
346 
347  inline std::size_t
349  {
350  return 0;
351  }
352 
353 
354 
355  template <class Archive>
356  void
357  DoFIndicesOnFaces<2>::serialize(Archive &ar, const unsigned int)
358  {
359  ar &lines;
360  }
361 
362 
363 
364  inline std::size_t
366  {
368  }
369 
370 
371 
372  template <class Archive>
373  void
374  DoFIndicesOnFaces<3>::serialize(Archive &ar, const unsigned int)
375  {
376  ar &lines &quads;
377  }
378 
379 
380 
381  inline std::size_t
383  {
386  }
387 
388  template <int structdim>
389  template <int dim, int spacedim>
392  const ::hp::DoFHandler<dim, spacedim> &dof_handler,
393  const unsigned int obj_index,
394  const unsigned int fe_index,
395  const unsigned int local_index,
396  const unsigned int /*obj_level*/) const
397  {
398  Assert((fe_index !=
400  ExcMessage("You need to specify a FE index when working "
401  "with hp DoFHandlers"));
402  AssertIndexRange(fe_index, dof_handler.get_fe_collection().size());
404  local_index,
405  dof_handler.get_fe(fe_index).template n_dofs_per_object<structdim>());
406  AssertIndexRange(obj_index, dof_offsets.size());
407 
408  // make sure we are on an
409  // object for which DoFs have
410  // been allocated at all
411  Assert(dof_offsets[obj_index] != numbers::invalid_unsigned_int,
412  ExcMessage("You are trying to access degree of freedom "
413  "information for an object on which no such "
414  "information is available"));
415 
416  Assert(structdim < dim,
417  ExcMessage("This object can not be used for cells."));
418 
419  // there may be multiple finite elements associated with
420  // this object. hop along the list of index sets until we
421  // find the one with the correct fe_index, and then poke
422  // into that part. trigger an exception if we can't find a
423  // set for this particular fe_index
424  const types::global_dof_index starting_offset = dof_offsets[obj_index];
425  const types::global_dof_index *pointer = &dofs[starting_offset];
426  while (true)
427  {
429  if (*pointer == fe_index)
430  return *(pointer + 1 + local_index);
431  else
432  pointer += static_cast<types::global_dof_index>(
433  dof_handler.get_fe(*pointer)
434  .template n_dofs_per_object<structdim>() +
435  1);
436  }
437  }
438 
439 
440 
441  template <int structdim>
442  template <int dim, int spacedim>
443  inline void
445  const ::hp::DoFHandler<dim, spacedim> &dof_handler,
446  const unsigned int obj_index,
447  const unsigned int fe_index,
448  const unsigned int local_index,
450  const unsigned int /*obj_level*/)
451  {
452  Assert((fe_index !=
454  ExcMessage("You need to specify a FE index when working "
455  "with hp DoFHandlers"));
456  AssertIndexRange(fe_index, dof_handler.get_fe_collection().size());
458  local_index,
459  dof_handler.get_fe(fe_index).template n_dofs_per_object<structdim>());
460  AssertIndexRange(obj_index, dof_offsets.size());
461 
462  // make sure we are on an
463  // object for which DoFs have
464  // been allocated at all
465  Assert(dof_offsets[obj_index] != numbers::invalid_unsigned_int,
466  ExcMessage("You are trying to access degree of freedom "
467  "information for an object on which no such "
468  "information is available"));
469 
470  Assert(structdim < dim,
471  ExcMessage("This object can not be used for cells."));
472 
473  // there may be multiple finite elements associated with
474  // this object. hop along the list of index sets until we
475  // find the one with the correct fe_index, and then poke
476  // into that part. trigger an exception if we can't find a
477  // set for this particular fe_index
478  const types::global_dof_index starting_offset = dof_offsets[obj_index];
479  types::global_dof_index * pointer = &dofs[starting_offset];
480  while (true)
481  {
483  if (*pointer == fe_index)
484  {
485  *(pointer + 1 + local_index) = global_index;
486  return;
487  }
488  else
489  pointer += dof_handler.get_fe(*pointer)
490  .template n_dofs_per_object<structdim>() +
491  1;
492  }
493  }
494 
495 
496 
497  template <int structdim>
498  template <int dim, int spacedim>
499  inline unsigned int
501  const ::hp::DoFHandler<dim, spacedim> &dof_handler,
502  const unsigned int obj_index) const
503  {
504  AssertIndexRange(obj_index, dof_offsets.size());
505 
506  // make sure we are on an
507  // object for which DoFs have
508  // been allocated at all
509  if (dof_offsets[obj_index] == numbers::invalid_unsigned_int)
510  return 0;
511 
512  Assert(structdim < dim,
513  ExcMessage("This object can not be used for cells."));
514 
515  // there may be multiple finite elements associated with this
516  // object. hop along the list of index sets until we find the
517  // one with the correct fe_index, and then poke into that
518  // part. trigger an exception if we can't find a set for this
519  // particular fe_index
520  const unsigned int starting_offset = dof_offsets[obj_index];
521  const types::global_dof_index *pointer = &dofs[starting_offset];
522  unsigned int counter = 0;
523  while (true)
524  {
525  if (*pointer == numbers::invalid_dof_index)
526  // end of list reached
527  return counter;
528  else
529  {
530  ++counter;
531  pointer += dof_handler.get_fe(*pointer)
532  .template n_dofs_per_object<structdim>() +
533  1;
534  }
535  }
536  }
537 
538 
539 
540  template <int structdim>
541  template <int dim, int spacedim>
544  const ::hp::DoFHandler<dim, spacedim> &dof_handler,
545  const unsigned int /*obj_level*/,
546  const unsigned int obj_index,
547  const unsigned int n) const
548  {
549  AssertIndexRange(obj_index, dof_offsets.size());
550 
551  // make sure we are on an
552  // object for which DoFs have
553  // been allocated at all
554  Assert(dof_offsets[obj_index] != numbers::invalid_unsigned_int,
555  ExcMessage("You are trying to access degree of freedom "
556  "information for an object on which no such "
557  "information is available"));
558 
559  Assert(structdim < dim,
560  ExcMessage("This object can not be used for cells."));
561 
562  AssertIndexRange(n, n_active_fe_indices(dof_handler, obj_index));
563 
564  // there may be multiple finite elements associated with
565  // this object. hop along the list of index sets until we
566  // find the one with the correct fe_index, and then poke
567  // into that part. trigger an exception if we can't find a
568  // set for this particular fe_index
569  const unsigned int starting_offset = dof_offsets[obj_index];
570  const types::global_dof_index *pointer = &dofs[starting_offset];
571  unsigned int counter = 0;
572  while (true)
573  {
575 
576  const unsigned int fe_index = *pointer;
577 
578  Assert(fe_index < dof_handler.get_fe_collection().size(),
579  ExcInternalError());
580 
581  if (counter == n)
582  return fe_index;
583 
584  ++counter;
585  pointer += dof_handler.get_fe(fe_index)
586  .template n_dofs_per_object<structdim>() +
587  1;
588  }
589  }
590 
591 
592 
593  template <int structdim>
594  template <int dim, int spacedim>
595  inline bool
597  const ::hp::DoFHandler<dim, spacedim> &dof_handler,
598  const unsigned int obj_index,
599  const unsigned int fe_index,
600  const unsigned int /*obj_level*/) const
601  {
602  AssertIndexRange(obj_index, dof_offsets.size());
603  Assert((fe_index !=
605  ExcMessage("You need to specify a FE index when working "
606  "with hp DoFHandlers"));
607  AssertIndexRange(fe_index, dof_handler.get_fe_collection().size());
608 
609  // make sure we are on an
610  // object for which DoFs have
611  // been allocated at all
612  Assert(dof_offsets[obj_index] != numbers::invalid_unsigned_int,
613  ExcMessage("You are trying to access degree of freedom "
614  "information for an object on which no such "
615  "information is available"));
616 
617  Assert(structdim < dim,
618  ExcMessage("This object can not be used for cells."));
619 
620  // there may be multiple finite elements associated with
621  // this object. hop along the list of index sets until we
622  // find the one with the correct fe_index, and then poke
623  // into that part. trigger an exception if we can't find a
624  // set for this particular fe_index
625  const types::global_dof_index starting_offset = dof_offsets[obj_index];
626  const types::global_dof_index *pointer = &dofs[starting_offset];
627  while (true)
628  {
629  if (*pointer == numbers::invalid_dof_index)
630  // end of list reached
631  return false;
632  else if (*pointer == fe_index)
633  return true;
634  else
635  pointer += static_cast<types::global_dof_index>(
636  dof_handler.get_fe(*pointer)
637  .template n_dofs_per_object<structdim>() +
638  1);
639  }
640  }
641 
642  template <int structdim>
643  template <class Archive>
644  void
646  const unsigned int)
647  {
648  ar &dofs;
649  ar &dof_offsets;
650  }
651 
652 
653  } // namespace hp
654 
655 } // namespace internal
656 
658 
659 #endif
numbers::invalid_dof_index
const types::global_dof_index invalid_dof_index
Definition: types.h:206
internal::hp::DoFIndicesOnFacesOrEdges::nth_active_fe_index
types::global_dof_index nth_active_fe_index(const ::hp::DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_level, const unsigned int obj_index, const unsigned int n) const
Definition: dof_faces.h:543
internal::hp::DoFIndicesOnFacesOrEdges::memory_consumption
std::size_t memory_consumption() const
Definition: dof_faces.cc:30
internal::hp::DoFIndicesOnFacesOrEdges::n_active_fe_indices
unsigned int n_active_fe_indices(const ::hp::DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_index) const
Definition: dof_faces.h:500
internal::hp::DoFIndicesOnFaces< 3 >::lines
internal::hp::DoFIndicesOnFacesOrEdges< 1 > lines
Definition: dof_faces.h:315
AssertIndexRange
#define AssertIndexRange(index, range)
Definition: exceptions.h:1649
hp::DoFHandler
Definition: dof_handler.h:203
internal::hp::DoFIndicesOnFacesOrEdges::serialize
void serialize(Archive &ar, const unsigned int version)
Definition: dof_faces.h:645
internal::hp::DoFIndicesOnFacesOrEdges::fe_index_is_active
bool fe_index_is_active(const ::hp::DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_index, const unsigned int fe_index, const unsigned int obj_level) const
Definition: dof_faces.h:596
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
fe_collection.h
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
MemoryConsumption::memory_consumption
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
Definition: memory_consumption.h:268
internal::hp::DoFIndicesOnFacesOrEdges
Definition: dof_faces.h:95
unsigned int
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
internal::hp::DoFIndicesOnFaces< 3 >::quads
internal::hp::DoFIndicesOnFacesOrEdges< 2 > quads
Definition: dof_faces.h:320
internal::hp::DoFIndicesOnFacesOrEdges::dof_offsets
std::vector< unsigned int > dof_offsets
Definition: dof_faces.h:105
hp
Definition: hp.h:117
numbers::invalid_unsigned_int
static const unsigned int invalid_unsigned_int
Definition: types.h:191
internal::hp::DoFIndicesOnFacesOrEdges::get_dof_index
types::global_dof_index get_dof_index(const ::hp::DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index, const unsigned int obj_level) const
Definition: dof_faces.h:391
internal::hp::DoFIndicesOnFacesOrEdges::set_dof_index
void set_dof_index(const ::hp::DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index, const types::global_dof_index global_index, const unsigned int obj_level)
Definition: dof_faces.h:444
config.h
internal
Definition: aligned_vector.h:369
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
internal::hp::DoFIndicesOnFacesOrEdges::dofs
std::vector< types::global_dof_index > dofs
Definition: dof_faces.h:111
internal::hp::DoFIndicesOnFaces
Definition: dof_faces.h:239
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
internal::hp::DoFIndicesOnFaces< 2 >::lines
internal::hp::DoFIndicesOnFacesOrEdges< 1 > lines
Definition: dof_faces.h:283