Reference documentation for deal.II version 9.1.0-pre
fe_raviart_thomas.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_fe_raviart_thomas_h
17 #define dealii_fe_raviart_thomas_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_raviart_thomas.h>
24 #include <deal.II/base/table.h>
25 #include <deal.II/base/tensor_product_polynomials.h>
26 
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 
36 
104 template <int dim>
106  : public FE_PolyTensor<PolynomialsRaviartThomas<dim>, dim>
107 {
108 public:
112  FE_RaviartThomas(const unsigned int p);
113 
119  virtual std::string
120  get_name() const override;
121 
122  // documentation inherited from the base class
123  virtual std::unique_ptr<FiniteElement<dim, dim>>
124  clone() const override;
125 
133  virtual bool
134  has_support_on_face(const unsigned int shape_index,
135  const unsigned int face_index) const override;
136 
137  // documentation inherited from the base class
138  virtual void
140  const std::vector<Vector<double>> &support_point_values,
141  std::vector<double> & nodal_values) const override;
142 
147  virtual std::pair<Table<2, bool>, std::vector<unsigned int>>
148  get_constant_modes() const override;
149 
150  virtual std::size_t
151  memory_consumption() const override;
152 
153 private:
160  static std::vector<unsigned int>
161  get_dpo_vector(const unsigned int degree);
162 
168  void
169  initialize_support_points(const unsigned int rt_degree);
170 
177  void
179 
191 
199 
203  template <int dim1>
204  friend class FE_RaviartThomas;
205 };
206 
207 
208 
248 template <int dim>
250  : public FE_PolyTensor<PolynomialsRaviartThomas<dim>, dim>
251 {
252 public:
256  FE_RaviartThomasNodal(const unsigned int p);
257 
263  virtual std::string
264  get_name() const override;
265 
266  // documentation inherited from the base class
267  virtual std::unique_ptr<FiniteElement<dim, dim>>
268  clone() const override;
269 
270  virtual void
272  const std::vector<Vector<double>> &support_point_values,
273  std::vector<double> & nodal_values) const override;
274 
275  virtual void
277  FullMatrix<double> &matrix) const override;
278 
279  virtual void
281  const unsigned int subface,
282  FullMatrix<double> &matrix) const override;
283  virtual bool
284  hp_constraints_are_implemented() const override;
285 
286  virtual std::vector<std::pair<unsigned int, unsigned int>>
287  hp_vertex_dof_identities(const FiniteElement<dim> &fe_other) const override;
288 
289  virtual std::vector<std::pair<unsigned int, unsigned int>>
290  hp_line_dof_identities(const FiniteElement<dim> &fe_other) const override;
291 
292  virtual std::vector<std::pair<unsigned int, unsigned int>>
293  hp_quad_dof_identities(const FiniteElement<dim> &fe_other) const override;
294 
297  const FiniteElement<dim> &fe_other) const override;
298 
299 private:
306  static std::vector<unsigned int>
307  get_dpo_vector(const unsigned int degree);
308 
313  static std::vector<bool>
314  get_ria_vector(const unsigned int degree);
315 
323  virtual bool
324  has_support_on_face(const unsigned int shape_index,
325  const unsigned int face_index) const override;
335  void
336  initialize_support_points(const unsigned int rt_degree);
337 };
338 
339 
342 /* -------------- declaration of explicit specializations ------------- */
343 
344 #ifndef DOXYGEN
345 
346 template <>
347 void
349 
350 #endif // DOXYGEN
351 
352 DEAL_II_NAMESPACE_CLOSE
353 
354 #endif
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
const unsigned int degree
Definition: fe_base.h:313
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
virtual bool hp_constraints_are_implemented() const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
virtual std::size_t memory_consumption() const override
virtual std::string get_name() const override
friend class FE_RaviartThomas
void initialize_support_points(const unsigned int rt_degree)
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const override
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double >> &support_point_values, std::vector< double > &nodal_values) const override
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
Table< 3, double > interior_weights
Table< 2, double > boundary_weights