Reference documentation for deal.II version 8.5.1
fe_raviart_thomas.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2003 - 2017 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 
34 
102 template <int dim>
104  :
105  public FE_PolyTensor<PolynomialsRaviartThomas<dim>, dim>
106 {
107 public:
111  FE_RaviartThomas (const unsigned int p);
112 
118  virtual std::string get_name () const;
119 
120 
128  virtual bool has_support_on_face (const unsigned int shape_index,
129  const unsigned int face_index) const;
130 
131  // documentation inherited from the base class
132  virtual
133  void
134  convert_generalized_support_point_values_to_nodal_values (const std::vector<Vector<double> > &support_point_values,
135  std::vector<double> &nodal_values) const;
136 
137  virtual void interpolate(std::vector<double> &local_dofs,
138  const std::vector<double> &values) const DEAL_II_DEPRECATED;
139 
140  virtual void interpolate(std::vector<double> &local_dofs,
141  const std::vector<Vector<double> > &values,
142  const unsigned int offset = 0) const DEAL_II_DEPRECATED;
143 
144  virtual void interpolate(std::vector<double> &local_dofs,
145  const VectorSlice<const std::vector<std::vector<double> > > &values) const DEAL_II_DEPRECATED;
146 
151  virtual std::pair<Table<2,bool>, std::vector<unsigned int> >
152  get_constant_modes () const;
153 
154  virtual std::size_t memory_consumption () const;
155  virtual FiniteElement<dim> *clone() const;
156 
157 private:
164  static std::vector<unsigned int>
165  get_dpo_vector (const unsigned int degree);
166 
172  void initialize_support_points (const unsigned int rt_degree);
173 
180  void initialize_restriction ();
181 
193 
201 
205  template <int dim1> friend class FE_RaviartThomas;
206 };
207 
208 
209 
249 template <int dim>
251  :
252  public FE_PolyTensor<PolynomialsRaviartThomas<dim>, dim>
253 {
254 public:
258  FE_RaviartThomasNodal (const unsigned int p);
259 
265  virtual std::string get_name () const;
266 
267  virtual FiniteElement<dim> *clone () const;
268 
269  // documentation inherited from the base class
270  virtual
271  void
272  convert_generalized_support_point_values_to_nodal_values (const std::vector<Vector<double> > &support_point_values,
273  std::vector<double> &nodal_values) const;
274  virtual void interpolate(std::vector<double> &local_dofs,
275  const std::vector<double> &values) const DEAL_II_DEPRECATED;
276 
277  virtual void interpolate(std::vector<double> &local_dofs,
278  const std::vector<Vector<double> > &values,
279  const unsigned int offset = 0) const DEAL_II_DEPRECATED;
280 
281  virtual void interpolate(std::vector<double> &local_dofs,
282  const VectorSlice<const std::vector<std::vector<double> > > &values) const DEAL_II_DEPRECATED;
283 
284 
285  virtual void get_face_interpolation_matrix (const FiniteElement<dim> &source,
286  FullMatrix<double> &matrix) const;
287 
288  virtual void get_subface_interpolation_matrix (const FiniteElement<dim> &source,
289  const unsigned int subface,
290  FullMatrix<double> &matrix) const;
291  virtual bool hp_constraints_are_implemented () const;
292 
293  virtual std::vector<std::pair<unsigned int, unsigned int> >
294  hp_vertex_dof_identities (const FiniteElement<dim> &fe_other) const;
295 
296  virtual std::vector<std::pair<unsigned int, unsigned int> >
297  hp_line_dof_identities (const FiniteElement<dim> &fe_other) const;
298 
299  virtual std::vector<std::pair<unsigned int, unsigned int> >
300  hp_quad_dof_identities (const FiniteElement<dim> &fe_other) const;
301 
303  compare_for_face_domination (const FiniteElement<dim> &fe_other) const;
304 
305 private:
312  static std::vector<unsigned int>
313  get_dpo_vector (const unsigned int degree);
314 
319  static std::vector<bool>
320  get_ria_vector (const unsigned int degree);
321 
329  virtual bool has_support_on_face (const unsigned int shape_index,
330  const unsigned int face_index) const;
340  void initialize_support_points (const unsigned int rt_degree);
341 };
342 
343 
346 /* -------------- declaration of explicit specializations ------------- */
347 
348 #ifndef DOXYGEN
349 
350 template <>
351 void
353 
354 #endif // DOXYGEN
355 
356 DEAL_II_NAMESPACE_CLOSE
357 
358 #endif
virtual std::string get_name() const
virtual FiniteElement< dim > * clone() const
virtual void convert_generalized_support_point_values_to_nodal_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const
virtual void interpolate(std::vector< double > &local_dofs, const std::vector< double > &values) const 1
FE_RaviartThomasNodal(const unsigned int p)
virtual void interpolate(std::vector< double > &local_dofs, const std::vector< double > &values) const 1
const unsigned int degree
Definition: fe_base.h:311
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
friend class FE_RaviartThomas
virtual void convert_generalized_support_point_values_to_nodal_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const
virtual std::string get_name() const
void initialize_support_points(const unsigned int rt_degree)
void initialize_support_points(const unsigned int rt_degree)
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 bool hp_constraints_are_implemented() const
virtual FiniteElement< dim > * clone() const
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
Table< 3, double > interior_weights
Table< 2, double > boundary_weights
static std::vector< bool > get_ria_vector(const unsigned int degree)