Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
dof_faces.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2006 - 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_faces_h
17 #define dealii_hp_dof_faces_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #include <deal.II/hp/fe_collection.h>
23 
24 #include <vector>
25 
26 DEAL_II_NAMESPACE_OPEN
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,
131  const types::global_dof_index global_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  template <class Archive>
347  void
348  DoFIndicesOnFaces<2>::serialize(Archive &ar, const unsigned int)
349  {
350  ar &lines;
351  }
352 
353 
354  template <class Archive>
355  void
356  DoFIndicesOnFaces<3>::serialize(Archive &ar, const unsigned int)
357  {
358  ar &lines &quads;
359  }
360 
361  template <int structdim>
362  template <int dim, int spacedim>
365  const ::hp::DoFHandler<dim, spacedim> &dof_handler,
366  const unsigned int obj_index,
367  const unsigned int fe_index,
368  const unsigned int local_index,
369  const unsigned int /*obj_level*/) const
370  {
371  Assert((fe_index !=
373  ExcMessage("You need to specify a FE index when working "
374  "with hp DoFHandlers"));
375  Assert(fe_index < dof_handler.get_fe_collection().size(),
376  ExcIndexRange(fe_index,
377  0,
378  dof_handler.get_fe_collection().size()));
379  Assert(local_index < dof_handler.get_fe(fe_index)
380  .template n_dofs_per_object<structdim>(),
381  ExcIndexRange(local_index,
382  0,
383  dof_handler.get_fe(fe_index)
384  .template n_dofs_per_object<structdim>()));
385  Assert(obj_index < dof_offsets.size(),
386  ExcIndexRange(obj_index, 0, dof_offsets.size()));
387 
388  // make sure we are on an
389  // object for which DoFs have
390  // been allocated at all
391  Assert(dof_offsets[obj_index] != numbers::invalid_unsigned_int,
392  ExcMessage("You are trying to access degree of freedom "
393  "information for an object on which no such "
394  "information is available"));
395 
396  Assert(structdim < dim,
397  ExcMessage("This object can not be used for cells."));
398 
399  // there may be multiple finite elements associated with
400  // this object. hop along the list of index sets until we
401  // find the one with the correct fe_index, and then poke
402  // into that part. trigger an exception if we can't find a
403  // set for this particular fe_index
404  const types::global_dof_index starting_offset = dof_offsets[obj_index];
405  const types::global_dof_index *pointer = &dofs[starting_offset];
406  while (true)
407  {
409  if (*pointer == fe_index)
410  return *(pointer + 1 + local_index);
411  else
412  pointer += static_cast<types::global_dof_index>(
413  dof_handler.get_fe(*pointer)
414  .template n_dofs_per_object<structdim>() +
415  1);
416  }
417  }
418 
419 
420 
421  template <int structdim>
422  template <int dim, int spacedim>
423  inline void
425  const ::hp::DoFHandler<dim, spacedim> &dof_handler,
426  const unsigned int obj_index,
427  const unsigned int fe_index,
428  const unsigned int local_index,
429  const types::global_dof_index global_index,
430  const unsigned int /*obj_level*/)
431  {
432  Assert((fe_index !=
434  ExcMessage("You need to specify a FE index when working "
435  "with hp DoFHandlers"));
436  Assert(fe_index < dof_handler.get_fe_collection().size(),
437  ExcIndexRange(fe_index,
438  0,
439  dof_handler.get_fe_collection().size()));
440  Assert(local_index < dof_handler.get_fe(fe_index)
441  .template n_dofs_per_object<structdim>(),
442  ExcIndexRange(local_index,
443  0,
444  dof_handler.get_fe(fe_index)
445  .template n_dofs_per_object<structdim>()));
446  Assert(obj_index < dof_offsets.size(),
447  ExcIndexRange(obj_index, 0, dof_offsets.size()));
448 
449  // make sure we are on an
450  // object for which DoFs have
451  // been allocated at all
452  Assert(dof_offsets[obj_index] != numbers::invalid_unsigned_int,
453  ExcMessage("You are trying to access degree of freedom "
454  "information for an object on which no such "
455  "information is available"));
456 
457  Assert(structdim < dim,
458  ExcMessage("This object can not be used for cells."));
459 
460  // there may be multiple finite elements associated with
461  // this object. hop along the list of index sets until we
462  // find the one with the correct fe_index, and then poke
463  // into that part. trigger an exception if we can't find a
464  // set for this particular fe_index
465  const types::global_dof_index starting_offset = dof_offsets[obj_index];
466  types::global_dof_index * pointer = &dofs[starting_offset];
467  while (true)
468  {
470  if (*pointer == fe_index)
471  {
472  *(pointer + 1 + local_index) = global_index;
473  return;
474  }
475  else
476  pointer += dof_handler.get_fe(*pointer)
477  .template n_dofs_per_object<structdim>() +
478  1;
479  }
480  }
481 
482 
483 
484  template <int structdim>
485  template <int dim, int spacedim>
486  inline unsigned int
488  const ::hp::DoFHandler<dim, spacedim> &dof_handler,
489  const unsigned int obj_index) const
490  {
491  Assert(obj_index < dof_offsets.size(),
492  ExcIndexRange(obj_index, 0, dof_offsets.size()));
493 
494  // make sure we are on an
495  // object for which DoFs have
496  // been allocated at all
497  if (dof_offsets[obj_index] == numbers::invalid_unsigned_int)
498  return 0;
499 
500  Assert(structdim < dim,
501  ExcMessage("This object can not be used for cells."));
502 
503  // there may be multiple finite elements associated with this
504  // object. hop along the list of index sets until we find the
505  // one with the correct fe_index, and then poke into that
506  // part. trigger an exception if we can't find a set for this
507  // particular fe_index
508  const unsigned int starting_offset = dof_offsets[obj_index];
509  const types::global_dof_index *pointer = &dofs[starting_offset];
510  unsigned int counter = 0;
511  while (true)
512  {
513  if (*pointer == numbers::invalid_dof_index)
514  // end of list reached
515  return counter;
516  else
517  {
518  ++counter;
519  pointer += dof_handler.get_fe(*pointer)
520  .template n_dofs_per_object<structdim>() +
521  1;
522  }
523  }
524  }
525 
526 
527 
528  template <int structdim>
529  template <int dim, int spacedim>
532  const ::hp::DoFHandler<dim, spacedim> &dof_handler,
533  const unsigned int /*obj_level*/,
534  const unsigned int obj_index,
535  const unsigned int n) const
536  {
537  Assert(obj_index < dof_offsets.size(),
538  ExcIndexRange(obj_index, 0, dof_offsets.size()));
539 
540  // make sure we are on an
541  // object for which DoFs have
542  // been allocated at all
543  Assert(dof_offsets[obj_index] != numbers::invalid_unsigned_int,
544  ExcMessage("You are trying to access degree of freedom "
545  "information for an object on which no such "
546  "information is available"));
547 
548  Assert(structdim < dim,
549  ExcMessage("This object can not be used for cells."));
550 
551  Assert(n < n_active_fe_indices(dof_handler, obj_index),
552  ExcIndexRange(n, 0, n_active_fe_indices(dof_handler, obj_index)));
553 
554  // there may be multiple finite elements associated with
555  // this object. hop along the list of index sets until we
556  // find the one with the correct fe_index, and then poke
557  // into that part. trigger an exception if we can't find a
558  // set for this particular fe_index
559  const unsigned int starting_offset = dof_offsets[obj_index];
560  const types::global_dof_index *pointer = &dofs[starting_offset];
561  unsigned int counter = 0;
562  while (true)
563  {
565 
566  const unsigned int fe_index = *pointer;
567 
568  Assert(fe_index < dof_handler.get_fe_collection().size(),
569  ExcInternalError());
570 
571  if (counter == n)
572  return fe_index;
573 
574  ++counter;
575  pointer += dof_handler.get_fe(fe_index)
576  .template n_dofs_per_object<structdim>() +
577  1;
578  }
579  }
580 
581 
582 
583  template <int structdim>
584  template <int dim, int spacedim>
585  inline bool
587  const ::hp::DoFHandler<dim, spacedim> &dof_handler,
588  const unsigned int obj_index,
589  const unsigned int fe_index,
590  const unsigned int /*obj_level*/) const
591  {
592  Assert(obj_index < dof_offsets.size(),
593  ExcIndexRange(obj_index,
594  0,
595  static_cast<unsigned int>(dof_offsets.size())));
596  Assert((fe_index !=
598  ExcMessage("You need to specify a FE index when working "
599  "with hp DoFHandlers"));
600  Assert(fe_index < dof_handler.get_fe_collection().size(),
601  ExcIndexRange(fe_index,
602  0,
603  dof_handler.get_fe_collection().size()));
604 
605  // make sure we are on an
606  // object for which DoFs have
607  // been allocated at all
608  Assert(dof_offsets[obj_index] != numbers::invalid_unsigned_int,
609  ExcMessage("You are trying to access degree of freedom "
610  "information for an object on which no such "
611  "information is available"));
612 
613  Assert(structdim < dim,
614  ExcMessage("This object can not be used for cells."));
615 
616  // there may be multiple finite elements associated with
617  // this object. hop along the list of index sets until we
618  // find the one with the correct fe_index, and then poke
619  // into that part. trigger an exception if we can't find a
620  // set for this particular fe_index
621  const types::global_dof_index starting_offset = dof_offsets[obj_index];
622  const types::global_dof_index *pointer = &dofs[starting_offset];
623  while (true)
624  {
625  if (*pointer == numbers::invalid_dof_index)
626  // end of list reached
627  return false;
628  else if (*pointer == fe_index)
629  return true;
630  else
631  pointer += static_cast<types::global_dof_index>(
632  dof_handler.get_fe(*pointer)
633  .template n_dofs_per_object<structdim>() +
634  1);
635  }
636  }
637 
638  template <int structdim>
639  template <class Archive>
640  void
642  const unsigned int)
643  {
644  ar &dofs;
645  ar &dof_offsets;
646  }
647 
648 
649  } // namespace hp
650 
651 } // namespace internal
652 
653 DEAL_II_NAMESPACE_CLOSE
654 
655 #endif
static const unsigned int invalid_unsigned_int
Definition: types.h:173
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:424
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:364
unsigned int n_active_fe_indices(const ::hp::DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_index) const
Definition: dof_faces.h:487
internal::hp::DoFIndicesOnFacesOrEdges< 1 > lines
Definition: dof_faces.h:283
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
Definition: exceptions.h:1407
void serialize(Archive &ar, const unsigned int version)
Definition: dof_faces.h:641
std::size_t memory_consumption() const
Definition: dof_faces.cc:30
std::vector< types::global_dof_index > dofs
Definition: dof_faces.h:111
std::vector< unsigned int > dof_offsets
Definition: dof_faces.h:105
Definition: hp.h:117
internal::hp::DoFIndicesOnFacesOrEdges< 2 > quads
Definition: dof_faces.h:320
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:586
unsigned int global_dof_index
Definition: types.h:89
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:531
const types::global_dof_index invalid_dof_index
Definition: types.h:188
internal::hp::DoFIndicesOnFacesOrEdges< 1 > lines
Definition: dof_faces.h:315
static ::ExceptionBase & ExcInternalError()