Reference documentation for deal.II version 8.5.1
dof_faces.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2006 - 2015 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_faces_h
17 #define dealii__hp_dof_faces_h
18 
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/hp/fe_collection.h>
22 
23 #include <vector>
24 
25 DEAL_II_NAMESPACE_OPEN
26 
27 namespace hp
28 {
29  template <int dim, int spacedim>
30  class FECollection;
31 }
32 
33 
34 namespace internal
35 {
36  namespace hp
37  {
100  template <int structdim>
102  {
103  public:
111  std::vector<unsigned int> dof_offsets;
112 
117  std::vector<types::global_dof_index> dofs;
118 
131  template <int dim, int spacedim>
132  void
133  set_dof_index (const ::hp::DoFHandler<dim,spacedim> &dof_handler,
134  const unsigned int obj_index,
135  const unsigned int fe_index,
136  const unsigned int local_index,
137  const types::global_dof_index global_index,
138  const unsigned int obj_level);
139 
151  template <int dim, int spacedim>
153  get_dof_index (const ::hp::DoFHandler<dim,spacedim> &dof_handler,
154  const unsigned int obj_index,
155  const unsigned int fe_index,
156  const unsigned int local_index,
157  const unsigned int obj_level) const;
158 
170  template <int dim, int spacedim>
171  unsigned int
172  n_active_fe_indices (const ::hp::DoFHandler<dim,spacedim> &dof_handler,
173  const unsigned int obj_index) const;
174 
178  template <int dim, int spacedim>
180  nth_active_fe_index (const ::hp::DoFHandler<dim,spacedim> &dof_handler,
181  const unsigned int obj_level,
182  const unsigned int obj_index,
183  const unsigned int n) const;
184 
189  template <int dim, int spacedim>
190  bool
191  fe_index_is_active (const ::hp::DoFHandler<dim,spacedim> &dof_handler,
192  const unsigned int obj_index,
193  const unsigned int fe_index,
194  const unsigned int obj_level) const;
195 
200  std::size_t memory_consumption () const;
201 
206  template <class Archive>
207  void serialize(Archive &ar,
208  const unsigned int version);
209  };
210 
211 
212 
240  template<int dim>
242 
243 
251  template<>
253  {
254  public:
259  std::size_t memory_consumption () const;
260 
265  template <class Archive>
266  void serialize(Archive &ar,
267  const unsigned int version);
268  };
269 
277  template<>
279  {
280  public:
285 
290  std::size_t memory_consumption () const;
291 
296  template <class Archive>
297  void serialize(Archive &ar,
298  const unsigned int version);
299  };
300 
308  template<>
310  {
311  public:
316 
321 
326  std::size_t memory_consumption () const;
327 
332  template <class Archive>
333  void serialize(Archive &ar,
334  const unsigned int version);
335  };
336 
337 
338  // --------------------- inline and template functions ------------------
339  template <class Archive>
341  const unsigned int)
342  {}
343 
344 
345  template <class Archive>
347  const unsigned int)
348  {
349  ar &lines;
350  }
351 
352 
353  template <class Archive>
355  const unsigned int)
356  {
357  ar &lines &quads;
358  }
359 
360  template <int structdim>
361  template <int dim, int spacedim>
362  inline
365  get_dof_index (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  {
372  ExcMessage ("You need to specify a FE index when working "
373  "with hp DoFHandlers"));
374  Assert (fe_index < dof_handler.get_fe().size(),
375  ExcIndexRange (fe_index, 0, dof_handler.get_fe().size()));
376  Assert (local_index <
377  dof_handler.get_fe()[fe_index].template n_dofs_per_object<structdim>(),
378  ExcIndexRange(local_index, 0,
379  dof_handler.get_fe()[fe_index]
380  .template n_dofs_per_object<structdim>()));
381  Assert (obj_index < dof_offsets.size(),
382  ExcIndexRange (obj_index, 0, dof_offsets.size()));
383 
384  // make sure we are on an
385  // object for which DoFs have
386  // been allocated at all
387  Assert (dof_offsets[obj_index] != numbers::invalid_unsigned_int,
388  ExcMessage ("You are trying to access degree of freedom "
389  "information for an object on which no such "
390  "information is available"));
391 
392  Assert (structdim<dim, ExcMessage ("This object can not be used for cells."));
393 
394  // there may be multiple finite elements associated with
395  // this object. hop along the list of index sets until we
396  // find the one with the correct fe_index, and then poke
397  // into that part. trigger an exception if we can't find a
398  // set for this particular fe_index
399  const types::global_dof_index starting_offset = dof_offsets[obj_index];
400  const types::global_dof_index *pointer = &dofs[starting_offset];
401  while (true)
402  {
403  Assert (*pointer != numbers::invalid_dof_index,
404  ExcInternalError());
405  if (*pointer == fe_index)
406  return *(pointer + 1 + local_index);
407  else
408  pointer += static_cast<types::global_dof_index>(
409  dof_handler.get_fe()[*pointer]
410  .template n_dofs_per_object<structdim>() + 1);
411  }
412  }
413 
414 
415 
416  template <int structdim>
417  template <int dim, int spacedim>
418  inline
419  void
421  set_dof_index (const ::hp::DoFHandler<dim,spacedim> &dof_handler,
422  const unsigned int obj_index,
423  const unsigned int fe_index,
424  const unsigned int local_index,
425  const types::global_dof_index global_index,
426  const unsigned int /*obj_level*/)
427  {
429  ExcMessage ("You need to specify a FE index when working "
430  "with hp DoFHandlers"));
431  Assert (fe_index < dof_handler.get_fe().size(),
432  ExcIndexRange (fe_index, 0, dof_handler.get_fe().size()));
433  Assert (local_index <
434  dof_handler.get_fe()[fe_index].template n_dofs_per_object<structdim>(),
435  ExcIndexRange(local_index, 0,
436  dof_handler.get_fe()[fe_index]
437  .template n_dofs_per_object<structdim>()));
438  Assert (obj_index < dof_offsets.size(),
439  ExcIndexRange (obj_index, 0, dof_offsets.size()));
440 
441  // make sure we are on an
442  // object for which DoFs have
443  // been allocated at all
444  Assert (dof_offsets[obj_index] != numbers::invalid_unsigned_int,
445  ExcMessage ("You are trying to access degree of freedom "
446  "information for an object on which no such "
447  "information is available"));
448 
449  Assert (structdim<dim, ExcMessage ("This object can not be used for cells."));
450 
451  // there may be multiple finite elements associated with
452  // this object. hop along the list of index sets until we
453  // find the one with the correct fe_index, and then poke
454  // into that part. trigger an exception if we can't find a
455  // set for this particular fe_index
456  const types::global_dof_index starting_offset = dof_offsets[obj_index];
457  types::global_dof_index *pointer = &dofs[starting_offset];
458  while (true)
459  {
460  Assert (*pointer != numbers::invalid_dof_index,
461  ExcInternalError());
462  if (*pointer == fe_index)
463  {
464  *(pointer + 1 + local_index) = global_index;
465  return;
466  }
467  else
468  pointer += dof_handler.get_fe()[*pointer]
469  .template n_dofs_per_object<structdim>() + 1;
470  }
471  }
472 
473 
474 
475  template <int structdim>
476  template <int dim, int spacedim>
477  inline
478  unsigned int
480  n_active_fe_indices (const ::hp::DoFHandler<dim,spacedim> &dof_handler,
481  const unsigned int obj_index) const
482  {
483  Assert (obj_index < dof_offsets.size(),
484  ExcIndexRange (obj_index, 0, dof_offsets.size()));
485 
486  // make sure we are on an
487  // object for which DoFs have
488  // been allocated at all
489  if (dof_offsets[obj_index] == numbers::invalid_unsigned_int)
490  return 0;
491 
492  Assert (structdim<dim, ExcMessage ("This object can not be used for cells."));
493 
494  // there may be multiple finite elements associated with this
495  // object. hop along the list of index sets until we find the
496  // one with the correct fe_index, and then poke into that
497  // part. trigger an exception if we can't find a set for this
498  // particular fe_index
499  const unsigned int starting_offset = dof_offsets[obj_index];
500  const types::global_dof_index *pointer = &dofs[starting_offset];
501  unsigned int counter = 0;
502  while (true)
503  {
504  if (*pointer == numbers::invalid_dof_index)
505  // end of list reached
506  return counter;
507  else
508  {
509  ++counter;
510  pointer += dof_handler.get_fe()[*pointer]
511  .template n_dofs_per_object<structdim>() + 1;
512  }
513  }
514  }
515 
516 
517 
518  template <int structdim>
519  template <int dim, int spacedim>
520  inline
523  nth_active_fe_index (const ::hp::DoFHandler<dim,spacedim> &dof_handler,
524  const unsigned int /*obj_level*/,
525  const unsigned int obj_index,
526  const unsigned int n) const
527  {
528  Assert (obj_index < dof_offsets.size(),
529  ExcIndexRange (obj_index, 0, dof_offsets.size()));
530 
531  // make sure we are on an
532  // object for which DoFs have
533  // been allocated at all
534  Assert (dof_offsets[obj_index] != numbers::invalid_unsigned_int,
535  ExcMessage ("You are trying to access degree of freedom "
536  "information for an object on which no such "
537  "information is available"));
538 
539  Assert (structdim<dim, ExcMessage ("This object can not be used for cells."));
540 
541  Assert (n < n_active_fe_indices(dof_handler, obj_index),
542  ExcIndexRange (n, 0,
543  n_active_fe_indices(dof_handler, obj_index)));
544 
545  // there may be multiple finite elements associated with
546  // this object. hop along the list of index sets until we
547  // find the one with the correct fe_index, and then poke
548  // into that part. trigger an exception if we can't find a
549  // set for this particular fe_index
550  const unsigned int starting_offset = dof_offsets[obj_index];
551  const types::global_dof_index *pointer = &dofs[starting_offset];
552  unsigned int counter = 0;
553  while (true)
554  {
555  Assert (*pointer != numbers::invalid_dof_index,
556  ExcInternalError());
557 
558  const unsigned int fe_index = *pointer;
559 
560  Assert (fe_index < dof_handler.get_fe().size(),
561  ExcInternalError());
562 
563  if (counter == n)
564  return fe_index;
565 
566  ++counter;
567  pointer += dof_handler.get_fe()[fe_index]
568  .template n_dofs_per_object<structdim>() + 1;
569  }
570  }
571 
572 
573 
574  template <int structdim>
575  template <int dim, int spacedim>
576  inline
577  bool
579  fe_index_is_active (const ::hp::DoFHandler<dim,spacedim> &dof_handler,
580  const unsigned int obj_index,
581  const unsigned int fe_index,
582  const unsigned int /*obj_level*/) const
583  {
584  Assert (obj_index < dof_offsets.size(),
585  ExcIndexRange (obj_index, 0, static_cast<unsigned int>(dof_offsets.size())));
587  ExcMessage ("You need to specify a FE index when working "
588  "with hp DoFHandlers"));
589  Assert (fe_index < dof_handler.get_fe().size(),
590  ExcIndexRange (fe_index, 0, dof_handler.get_fe().size()));
591 
592  // make sure we are on an
593  // object for which DoFs have
594  // been allocated at all
595  Assert (dof_offsets[obj_index] != numbers::invalid_unsigned_int,
596  ExcMessage ("You are trying to access degree of freedom "
597  "information for an object on which no such "
598  "information is available"));
599 
600  Assert (structdim<dim, ExcMessage ("This object can not be used for cells."));
601 
602  // there may be multiple finite elements associated with
603  // this object. hop along the list of index sets until we
604  // find the one with the correct fe_index, and then poke
605  // into that part. trigger an exception if we can't find a
606  // set for this particular fe_index
607  const types::global_dof_index starting_offset = dof_offsets[obj_index];
608  const types::global_dof_index *pointer = &dofs[starting_offset];
609  while (true)
610  {
611  if (*pointer == numbers::invalid_dof_index)
612  // end of list reached
613  return false;
614  else if (*pointer == fe_index)
615  return true;
616  else
617  pointer += static_cast<types::global_dof_index>(
618  dof_handler.get_fe()[*pointer]
619  .template n_dofs_per_object<structdim>()+1);
620  }
621  }
622 
623  template <int structdim>
624  template <class Archive>
626  const unsigned int)
627  {
628  ar &dofs;
629  ar &dof_offsets;
630  }
631 
632 
633  } // namespace hp
634 
635 } // namespace internal
636 
637 DEAL_II_NAMESPACE_CLOSE
638 
639 #endif
static const unsigned int invalid_unsigned_int
Definition: types.h:170
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:421
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:365
unsigned int n_active_fe_indices(const ::hp::DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_index) const
Definition: dof_faces.h:480
internal::hp::DoFIndicesOnFacesOrEdges< 1 > lines
Definition: dof_faces.h:284
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
static ::ExceptionBase & ExcMessage(std::string arg1)
unsigned int global_dof_index
Definition: types.h:88
#define Assert(cond, exc)
Definition: exceptions.h:313
void serialize(Archive &ar, const unsigned int version)
Definition: dof_faces.h:625
std::size_t memory_consumption() const
Definition: dof_faces.cc:29
std::vector< types::global_dof_index > dofs
Definition: dof_faces.h:117
std::vector< unsigned int > dof_offsets
Definition: dof_faces.h:111
Definition: hp.h:102
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:579
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:523
const types::global_dof_index invalid_dof_index
Definition: types.h:184
internal::hp::DoFIndicesOnFacesOrEdges< 1 > lines
Definition: dof_faces.h:315
static ::ExceptionBase & ExcInternalError()