Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
fe_values.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2003 - 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_fe_values_h
17 #define dealii_hp_fe_values_h
18 
19 #include <deal.II/base/config.h>
20 
21 #include <deal.II/fe/fe.h>
22 #include <deal.II/fe/fe_values.h>
23 
24 #include <deal.II/hp/fe_collection.h>
25 #include <deal.II/hp/mapping_collection.h>
26 #include <deal.II/hp/q_collection.h>
27 
28 #include <map>
29 #include <memory>
30 
31 DEAL_II_NAMESPACE_OPEN
32 
33 template <int dim, int spacedim>
34 class FiniteElement;
35 
36 
37 
38 namespace internal
39 {
40  namespace hp
41  {
62  template <int dim, int q_dim, class FEValuesType>
64  {
65  public:
71  const ::hp::MappingCollection<dim, FEValuesType::space_dimension>
73  const ::hp::FECollection<dim, FEValuesType::space_dimension>
74  & fe_collection,
75  const ::hp::QCollection<q_dim> &q_collection,
83  const ::hp::FECollection<dim, FEValuesType::space_dimension>
84  & fe_collection,
85  const ::hp::QCollection<q_dim> &q_collection,
87 
92  const ::hp::FECollection<dim, FEValuesType::space_dimension> &
93  get_fe_collection() const;
94 
98  const ::hp::MappingCollection<dim, FEValuesType::space_dimension> &
99  get_mapping_collection() const;
100 
104  const ::hp::QCollection<q_dim> &
106 
111  get_update_flags() const;
112 
119  const FEValuesType &
120  get_present_fe_values() const;
121 
122  protected:
130  FEValuesType &
131  select_fe_values(const unsigned int fe_index,
132  const unsigned int mapping_index,
133  const unsigned int q_index);
134 
135  protected:
139  const SmartPointer<
140  const ::hp::FECollection<dim, FEValuesType::space_dimension>,
143 
147  const SmartPointer<
148  const ::hp::MappingCollection<dim, FEValuesType::space_dimension>,
151 
155  const ::hp::QCollection<q_dim> q_collection;
156 
157  private:
170 
176 
181  };
182 
183  } // namespace hp
184 
185 } // namespace internal
186 
187 
188 namespace hp
189 {
238  template <int dim, int spacedim = dim>
239  class FEValues : public ::internal::hp::
240  FEValuesBase<dim, dim, ::FEValues<dim, spacedim>>
241  {
242  public:
243  static const unsigned int dimension = dim;
244 
245  static const unsigned int space_dimension = spacedim;
246 
256  FEValues(
257  const ::hp::MappingCollection<dim, spacedim> &mapping_collection,
258  const ::hp::FECollection<dim, spacedim> & fe_collection,
259  const ::hp::QCollection<dim> & q_collection,
260  const UpdateFlags update_flags);
261 
262 
274  FEValues(const hp::FECollection<dim, spacedim> &fe_collection,
275  const hp::QCollection<dim> & q_collection,
276  const UpdateFlags update_flags);
277 
278 
324  template <typename DoFHandlerType, bool lda>
325  void
327  const unsigned int q_index = numbers::invalid_unsigned_int,
328  const unsigned int mapping_index = numbers::invalid_unsigned_int,
329  const unsigned int fe_index = numbers::invalid_unsigned_int);
330 
345  void
347  const unsigned int q_index = numbers::invalid_unsigned_int,
348  const unsigned int mapping_index = numbers::invalid_unsigned_int,
349  const unsigned int fe_index = numbers::invalid_unsigned_int);
350  };
351 
352 
353 
378  template <int dim, int spacedim = dim>
380  : public ::internal::hp::
381  FEValuesBase<dim, dim - 1, ::FEFaceValues<dim, spacedim>>
382  {
383  public:
393  FEFaceValues(const hp::MappingCollection<dim, spacedim> &mapping_collection,
394  const hp::FECollection<dim, spacedim> & fe_collection,
395  const hp::QCollection<dim - 1> & q_collection,
396  const UpdateFlags update_flags);
397 
398 
410  FEFaceValues(const hp::FECollection<dim, spacedim> &fe_collection,
411  const hp::QCollection<dim - 1> & q_collection,
412  const UpdateFlags update_flags);
413 
459  template <typename DoFHandlerType, bool lda>
460  void
462  const unsigned int face_no,
463  const unsigned int q_index = numbers::invalid_unsigned_int,
464  const unsigned int mapping_index = numbers::invalid_unsigned_int,
465  const unsigned int fe_index = numbers::invalid_unsigned_int);
466 
481  void
483  const unsigned int face_no,
484  const unsigned int q_index = numbers::invalid_unsigned_int,
485  const unsigned int mapping_index = numbers::invalid_unsigned_int,
486  const unsigned int fe_index = numbers::invalid_unsigned_int);
487  };
488 
489 
490 
498  template <int dim, int spacedim = dim>
500  : public ::internal::hp::
501  FEValuesBase<dim, dim - 1, ::FESubfaceValues<dim, spacedim>>
502  {
503  public:
514  const hp::MappingCollection<dim, spacedim> &mapping_collection,
515  const hp::FECollection<dim, spacedim> & fe_collection,
516  const hp::QCollection<dim - 1> & q_collection,
517  const UpdateFlags update_flags);
518 
519 
532  const hp::QCollection<dim - 1> & q_collection,
533  const UpdateFlags update_flags);
534 
570  template <typename DoFHandlerType, bool lda>
571  void
573  const unsigned int face_no,
574  const unsigned int subface_no,
575  const unsigned int q_index = numbers::invalid_unsigned_int,
576  const unsigned int mapping_index = numbers::invalid_unsigned_int,
577  const unsigned int fe_index = numbers::invalid_unsigned_int);
578 
593  void
595  const unsigned int face_no,
596  const unsigned int subface_no,
597  const unsigned int q_index = numbers::invalid_unsigned_int,
598  const unsigned int mapping_index = numbers::invalid_unsigned_int,
599  const unsigned int fe_index = numbers::invalid_unsigned_int);
600  };
601 
602 } // namespace hp
603 
604 
605 // -------------- inline and template functions --------------
606 
607 namespace internal
608 {
609  namespace hp
610  {
611  template <int dim, int q_dim, class FEValuesType>
612  inline const FEValuesType &
614  {
615  return *fe_values_table(present_fe_values_index);
616  }
617 
618 
619 
620  template <int dim, int q_dim, class FEValuesType>
621  inline const ::hp::FECollection<dim, FEValuesType::space_dimension> &
623  {
624  return *fe_collection;
625  }
626 
627 
628 
629  template <int dim, int q_dim, class FEValuesType>
630  inline const ::hp::MappingCollection<dim,
631  FEValuesType::space_dimension> &
633  {
634  return *mapping_collection;
635  }
636 
637 
638 
639  template <int dim, int q_dim, class FEValuesType>
640  inline const ::hp::QCollection<q_dim> &
642  {
643  return q_collection;
644  }
645 
646 
647 
648  template <int dim, int q_dim, class FEValuesType>
651  {
652  return update_flags;
653  }
654  } // namespace hp
655 
656 } // namespace internal
657 
658 DEAL_II_NAMESPACE_CLOSE
659 
660 #endif
static const unsigned int invalid_unsigned_int
Definition: types.h:173
const ::hp::MappingCollection< dim, FEValuesType::space_dimension > & get_mapping_collection() const
Definition: fe_values.h:632
FEValuesType & select_fe_values(const unsigned int fe_index, const unsigned int mapping_index, const unsigned int q_index)
Definition: fe_values.cc:71
const UpdateFlags update_flags
Definition: fe_values.h:180
TableIndices< 3 > present_fe_values_index
Definition: fe_values.h:175
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType, lda >> cell, const unsigned int face_no, const unsigned int subface_no, const unsigned int q_index=numbers::invalid_unsigned_int, const unsigned int mapping_index=numbers::invalid_unsigned_int, const unsigned int fe_index=numbers::invalid_unsigned_int)
Definition: fe_values.cc:390
const SmartPointer< const ::hp::FECollection< dim, FEValuesType::space_dimension >, FEValuesBase< dim, q_dim, FEValuesType > > fe_collection
Definition: fe_values.h:142
const ::hp::FECollection< dim, FEValuesType::space_dimension > & get_fe_collection() const
Definition: fe_values.h:622
const FEValuesType & get_present_fe_values() const
Definition: fe_values.h:613
FEValuesBase(const ::hp::MappingCollection< dim, FEValuesType::space_dimension > &mapping_collection, const ::hp::FECollection< dim, FEValuesType::space_dimension > &fe_collection, const ::hp::QCollection< q_dim > &q_collection, const ::UpdateFlags update_flags)
const ::hp::QCollection< q_dim > q_collection
Definition: fe_values.h:155
UpdateFlags
FEFaceValues(const hp::MappingCollection< dim, spacedim > &mapping_collection, const hp::FECollection< dim, spacedim > &fe_collection, const hp::QCollection< dim - 1 > &q_collection, const UpdateFlags update_flags)
Definition: fe_values.cc:235
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType, lda >> cell, const unsigned int q_index=numbers::invalid_unsigned_int, const unsigned int mapping_index=numbers::invalid_unsigned_int, const unsigned int fe_index=numbers::invalid_unsigned_int)
Definition: fe_values.cc:142
const ::hp::QCollection< q_dim > & get_quadrature_collection() const
Definition: fe_values.h:641
Definition: hp.h:117
::Table< 3, std::shared_ptr< FEValuesType > > fe_values_table
Definition: fe_values.h:169
FEValues(const ::hp::MappingCollection< dim, spacedim > &mapping_collection, const ::hp::FECollection< dim, spacedim > &fe_collection, const ::hp::QCollection< dim > &q_collection, const UpdateFlags update_flags)
Definition: table.h:37
UpdateFlags get_update_flags() const
Definition: fe_values.h:650
const SmartPointer< const ::hp::MappingCollection< dim, FEValuesType::space_dimension >, FEValuesBase< dim, q_dim, FEValuesType > > mapping_collection
Definition: fe_values.h:150
FESubfaceValues(const hp::MappingCollection< dim, spacedim > &mapping_collection, const hp::FECollection< dim, spacedim > &fe_collection, const hp::QCollection< dim - 1 > &q_collection, const UpdateFlags update_flags)
Definition: fe_values.cc:360
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType, lda >> cell, const unsigned int face_no, const unsigned int q_index=numbers::invalid_unsigned_int, const unsigned int mapping_index=numbers::invalid_unsigned_int, const unsigned int fe_index=numbers::invalid_unsigned_int)
Definition: fe_values.cc:265