Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
fe_dg_vector.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2010 - 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_fe_dg_vector_h
17 #define dealii_fe_dg_vector_h
18 
19 #include <deal.II/base/config.h>
20 
21 #include <deal.II/base/geometry_info.h>
22 #include <deal.II/base/polynomial.h>
23 #include <deal.II/base/polynomials_bdm.h>
24 #include <deal.II/base/polynomials_nedelec.h>
25 #include <deal.II/base/polynomials_raviart_thomas.h>
26 #include <deal.II/base/table.h>
27 #include <deal.II/base/tensor_product_polynomials.h>
28 
29 #include <deal.II/fe/fe.h>
30 #include <deal.II/fe/fe_poly_tensor.h>
31 
32 #include <vector>
33 
34 DEAL_II_NAMESPACE_OPEN
35 
54 template <class PolynomialType, int dim, int spacedim = dim>
55 class FE_DGVector : public FE_PolyTensor<PolynomialType, dim, spacedim>
56 {
57 public:
61  FE_DGVector(const unsigned int p, MappingType m);
62 
68  virtual std::string
69  get_name() const override;
70 
71  virtual std::unique_ptr<FiniteElement<dim, spacedim>>
72  clone() const override;
73 
80  virtual bool
81  has_support_on_face(const unsigned int shape_index,
82  const unsigned int face_index) const override;
83 
84  virtual std::size_t
85  memory_consumption() const override;
86 
87 private:
94  static std::vector<unsigned int>
95  get_dpo_vector(const unsigned int degree);
96 
104  {
105  public:
117  std::vector<std::vector<Tensor<1, dim>>> shape_values;
118 
128  std::vector<std::vector<Tensor<2, dim>>> shape_gradients;
129  };
130  Table<3, double> interior_weights;
131 };
132 
133 
134 
142 template <int dim, int spacedim = dim>
143 class FE_DGNedelec : public FE_DGVector<PolynomialsNedelec<dim>, dim, spacedim>
144 {
145 public:
150  FE_DGNedelec(const unsigned int p);
151 
157  virtual std::string
158  get_name() const override;
159 };
160 
161 
162 
171 template <int dim, int spacedim = dim>
173  : public FE_DGVector<PolynomialsRaviartThomas<dim>, dim, spacedim>
174 {
175 public:
179  FE_DGRaviartThomas(const unsigned int p);
180 
186  virtual std::string
187  get_name() const override;
188 };
189 
190 
191 
199 template <int dim, int spacedim = dim>
200 class FE_DGBDM : public FE_DGVector<PolynomialsBDM<dim>, dim, spacedim>
201 {
202 public:
206  FE_DGBDM(const unsigned int p);
207 
213  virtual std::string
214  get_name() const override;
215 };
216 
217 
218 DEAL_II_NAMESPACE_CLOSE
219 
220 #endif
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
virtual std::string get_name() const override
FE_DGNedelec(const unsigned int p)
Definition: fe_dg_vector.cc:29
MappingType
Definition: mapping.h:61
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
FE_DGRaviartThomas(const unsigned int p)
Definition: fe_dg_vector.cc:52
FE_DGVector(const unsigned int p, MappingType m)
FE_DGBDM(const unsigned int p)
Definition: fe_dg_vector.cc:79
virtual std::string get_name() const override
Definition: fe_dg_vector.cc:36
const unsigned int degree
Definition: fe_base.h:298
std::vector< std::vector< Tensor< 2, dim > > > shape_gradients
Definition: fe_dg_vector.h:128
virtual std::string get_name() const override
Definition: fe_dg_vector.cc:86
virtual std::size_t memory_consumption() const override
virtual std::string get_name() const override
Definition: fe_dg_vector.cc:61
std::vector< std::vector< Tensor< 1, dim > > > shape_values
Definition: fe_dg_vector.h:117
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)