Reference documentation for deal.II version 8.4.2
fe_dg_vector.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2010 - 2016 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__fe_dg_vector_h
17 #define dealii__fe_dg_vector_h
18 
19 #include <deal.II/base/config.h>
20 #include <deal.II/base/table.h>
21 #include <deal.II/base/polynomials_raviart_thomas.h>
22 #include <deal.II/base/polynomials_nedelec.h>
23 #include <deal.II/base/polynomials_bdm.h>
24 #include <deal.II/base/polynomial.h>
25 #include <deal.II/base/tensor_product_polynomials.h>
26 #include <deal.II/base/geometry_info.h>
27 #include <deal.II/fe/fe.h>
28 #include <deal.II/fe/fe_poly_tensor.h>
29 
30 #include <vector>
31 
32 DEAL_II_NAMESPACE_OPEN
33 
34 template <int dim, int spacedim> class MappingQ;
35 
36 
55 template <class PolynomialType, int dim, int spacedim=dim>
57  :
58  public FE_PolyTensor<PolynomialType, dim, spacedim>
59 {
60 public:
64  FE_DGVector (const unsigned int p, MappingType m);
65 public:
66 
68 
74  virtual std::string get_name () const;
75 
76 
83  virtual bool has_support_on_face (const unsigned int shape_index,
84  const unsigned int face_index) const;
85 
86  virtual void interpolate(std::vector<double> &local_dofs,
87  const std::vector<double> &values) const;
88  virtual void interpolate(std::vector<double> &local_dofs,
89  const std::vector<Vector<double> > &values,
90  unsigned int offset = 0) const;
91  virtual void interpolate(
92  std::vector<double> &local_dofs,
93  const VectorSlice<const std::vector<std::vector<double> > > &values) const;
94  virtual std::size_t memory_consumption () const;
95 
96 private:
103  static std::vector<unsigned int>
104  get_dpo_vector (const unsigned int degree);
105 
115  void initialize_support_points (const unsigned int degree);
116 
123  void initialize_restriction ();
124 
132  {
133  public:
145  std::vector<std::vector<Tensor<1,dim> > > shape_values;
146 
156  std::vector<std::vector<Tensor<2,dim> > > shape_gradients;
157  };
158  Table<3, double> interior_weights;
159 };
160 
161 
162 
170 template <int dim, int spacedim=dim>
171 class FE_DGNedelec : public FE_DGVector<PolynomialsNedelec<dim>, dim, spacedim>
172 {
173 public:
178  FE_DGNedelec (const unsigned int p);
179 
185  virtual std::string get_name () const;
186 };
187 
188 
189 
198 template <int dim, int spacedim=dim>
199 class FE_DGRaviartThomas : public FE_DGVector<PolynomialsRaviartThomas<dim>, dim, spacedim>
200 {
201 public:
205  FE_DGRaviartThomas (const unsigned int p);
206 
212  virtual std::string get_name () const;
213 };
214 
215 
216 
224 template <int dim, int spacedim=dim>
225 class FE_DGBDM : public FE_DGVector<PolynomialsBDM<dim>, dim, spacedim>
226 {
227 public:
231  FE_DGBDM (const unsigned int p);
232 
238  virtual std::string get_name () const;
239 };
240 
241 
242 DEAL_II_NAMESPACE_CLOSE
243 
244 #endif
void initialize_restriction()
MappingType
Definition: mapping.h:50
FE_DGVector(const unsigned int p, MappingType m)
const unsigned int degree
Definition: fe_base.h:299
std::vector< std::vector< Tensor< 1, dim > > > shape_values
Definition: fe_dg_vector.h:145
virtual std::size_t memory_consumption() const
virtual std::string get_name() const
virtual void interpolate(std::vector< double > &local_dofs, const std::vector< double > &values) const
FiniteElement< dim, spacedim > * clone() const
void initialize_support_points(const unsigned int degree)
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
std::vector< std::vector< Tensor< 2, dim > > > shape_gradients
Definition: fe_dg_vector.h:156