Reference documentation for deal.II version 8.4.2
fe_raviart_thomas.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2003 - 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__fe_raviart_thomas_h
17 #define dealii__fe_raviart_thomas_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/polynomial.h>
23 #include <deal.II/base/tensor_product_polynomials.h>
24 #include <deal.II/base/geometry_info.h>
25 #include <deal.II/fe/fe.h>
26 #include <deal.II/fe/fe_poly_tensor.h>
27 
28 #include <vector>
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 template <int dim, int spacedim> class MappingQ;
33 
34 
37 
105 template <int dim>
107  :
108  public FE_PolyTensor<PolynomialsRaviartThomas<dim>, dim>
109 {
110 public:
114  FE_RaviartThomas (const unsigned int p);
115 
121  virtual std::string get_name () const;
122 
123 
131  virtual bool has_support_on_face (const unsigned int shape_index,
132  const unsigned int face_index) const;
133 
134  virtual void interpolate(std::vector<double> &local_dofs,
135  const std::vector<double> &values) const;
136  virtual void interpolate(std::vector<double> &local_dofs,
137  const std::vector<Vector<double> > &values,
138  unsigned int offset = 0) const;
139  virtual void interpolate(
140  std::vector<double> &local_dofs,
141  const VectorSlice<const std::vector<std::vector<double> > > &values) const;
142 
147  virtual std::pair<Table<2,bool>, std::vector<unsigned int> >
148  get_constant_modes () const;
149 
150  virtual std::size_t memory_consumption () const;
151  virtual FiniteElement<dim> *clone() const;
152 
153 private:
160  static std::vector<unsigned int>
161  get_dpo_vector (const unsigned int degree);
162 
168  void initialize_support_points (const unsigned int rt_degree);
169 
176  void initialize_restriction ();
177 
189 
197 
201  template <int dim1> friend class FE_RaviartThomas;
202 };
203 
204 
205 
245 template <int dim>
247  :
248  public FE_PolyTensor<PolynomialsRaviartThomas<dim>, dim>
249 {
250 public:
254  FE_RaviartThomasNodal (const unsigned int p);
255 
261  virtual std::string get_name () const;
262 
263  virtual FiniteElement<dim> *clone () const;
264 
265  virtual void interpolate(std::vector<double> &local_dofs,
266  const std::vector<double> &values) const;
267  virtual void interpolate(std::vector<double> &local_dofs,
268  const std::vector<Vector<double> > &values,
269  unsigned int offset = 0) const;
270  virtual void interpolate(
271  std::vector<double> &local_dofs,
272  const VectorSlice<const std::vector<std::vector<double> > > &values) const;
273 
274 
275  virtual void get_face_interpolation_matrix (const FiniteElement<dim> &source,
276  FullMatrix<double> &matrix) const;
277 
278  virtual void get_subface_interpolation_matrix (const FiniteElement<dim> &source,
279  const unsigned int subface,
280  FullMatrix<double> &matrix) const;
281  virtual bool hp_constraints_are_implemented () const;
282 
283  virtual std::vector<std::pair<unsigned int, unsigned int> >
284  hp_vertex_dof_identities (const FiniteElement<dim> &fe_other) const;
285 
286  virtual std::vector<std::pair<unsigned int, unsigned int> >
287  hp_line_dof_identities (const FiniteElement<dim> &fe_other) const;
288 
289  virtual std::vector<std::pair<unsigned int, unsigned int> >
290  hp_quad_dof_identities (const FiniteElement<dim> &fe_other) const;
291 
293  compare_for_face_domination (const FiniteElement<dim> &fe_other) const;
294 
295 private:
302  static std::vector<unsigned int>
303  get_dpo_vector (const unsigned int degree);
304 
309  static std::vector<bool>
310  get_ria_vector (const unsigned int degree);
311 
319  virtual bool has_support_on_face (const unsigned int shape_index,
320  const unsigned int face_index) const;
330  void initialize_support_points (const unsigned int rt_degree);
331 };
332 
333 
336 /* -------------- declaration of explicit specializations ------------- */
337 
338 #ifndef DOXYGEN
339 
340 template <>
341 void
343 
344 #endif // DOXYGEN
345 
346 DEAL_II_NAMESPACE_CLOSE
347 
348 #endif
virtual std::string get_name() const
virtual FiniteElement< dim > * clone() const
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe.cc:868
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe.cc:924
const unsigned int degree
Definition: fe_base.h:299
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
friend class FE_RaviartThomas
void initialize_support_points(const unsigned int rt_degree)
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
Definition: fe.cc:885
virtual std::size_t memory_consumption() const
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe.cc:902
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe.cc:913
virtual bool hp_constraints_are_implemented() const
Definition: fe.cc:789
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
virtual void interpolate(std::vector< double > &local_dofs, const std::vector< double > &values) const
Table< 3, double > interior_weights
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe.cc:935
Table< 2, double > boundary_weights