Reference documentation for deal.II version Git 4fbb5374f0 2021-01-22 15:02:13 -0500
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
fe_raviart_thomas.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2003 - 2019 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 
24 #include <deal.II/base/table.h>
26 
27 #include <deal.II/fe/fe.h>
29 
30 #include <vector>
31 
33 
36 
126 template <int dim>
127 class FE_RaviartThomas : public FE_PolyTensor<dim>
128 {
129 public:
133  FE_RaviartThomas(const unsigned int p);
134 
140  virtual std::string
141  get_name() const override;
142 
143  // documentation inherited from the base class
144  virtual std::unique_ptr<FiniteElement<dim, dim>>
145  clone() const override;
146 
154  virtual bool
155  has_support_on_face(const unsigned int shape_index,
156  const unsigned int face_index) const override;
157 
158  // documentation inherited from the base class
159  virtual void
161  const std::vector<Vector<double>> &support_point_values,
162  std::vector<double> & nodal_values) const override;
163 
168  virtual std::pair<Table<2, bool>, std::vector<unsigned int>>
169  get_constant_modes() const override;
170 
171  virtual std::size_t
172  memory_consumption() const override;
173 
174 private:
181  static std::vector<unsigned int>
182  get_dpo_vector(const unsigned int degree);
183 
189  void
190  initialize_support_points(const unsigned int rt_degree);
191 
198  void
200 
212 
220 
221  // Allow access from other dimensions.
222  template <int dim1>
223  friend class FE_RaviartThomas;
224 };
225 
226 
227 
265 template <int dim>
267 {
268 public:
272  FE_RaviartThomasNodal(const unsigned int p);
273 
279  virtual std::string
280  get_name() const override;
281 
282  // documentation inherited from the base class
283  virtual std::unique_ptr<FiniteElement<dim, dim>>
284  clone() const override;
285 
286  virtual void
288  const std::vector<Vector<double>> &support_point_values,
289  std::vector<double> & nodal_values) const override;
290 
291  virtual void
294  const unsigned int face_no = 0) const override;
295 
296  virtual void
298  const FiniteElement<dim> &source,
299  const unsigned int subface,
300  FullMatrix<double> & matrix,
301  const unsigned int face_no = 0) const override;
302  virtual bool
303  hp_constraints_are_implemented() const override;
304 
305  virtual std::vector<std::pair<unsigned int, unsigned int>>
306  hp_vertex_dof_identities(const FiniteElement<dim> &fe_other) const override;
307 
308  virtual std::vector<std::pair<unsigned int, unsigned int>>
309  hp_line_dof_identities(const FiniteElement<dim> &fe_other) const override;
310 
311  virtual std::vector<std::pair<unsigned int, unsigned int>>
313  const unsigned int face_no = 0) const override;
314 
320  const unsigned int codim = 0) const override final;
321 
322 private:
329  static std::vector<unsigned int>
330  get_dpo_vector(const unsigned int degree);
331 
336  static std::vector<bool>
337  get_ria_vector(const unsigned int degree);
338 
346  virtual bool
347  has_support_on_face(const unsigned int shape_index,
348  const unsigned int face_index) const override;
358  void
359  initialize_support_points(const unsigned int rt_degree);
360 };
361 
362 
365 /* -------------- declaration of explicit specializations ------------- */
366 
367 #ifndef DOXYGEN
368 
369 template <>
370 void
372 
373 #endif // DOXYGEN
374 
376 
377 #endif
Contents is actually a matrix.
const unsigned int degree
Definition: fe_base.h:435
virtual std::size_t memory_consumption() const override
virtual std::string get_name() const override
friend class FE_RaviartThomas
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix, const unsigned int face_no=0) const
Definition: fe.cc:947
void initialize_support_points(const unsigned int rt_degree)
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:380
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe.cc:982
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix, const unsigned int face_no=0) const
Definition: fe.cc:964
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const override
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe.cc:993
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 FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const
Definition: fe.cc:1016
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:379
virtual bool hp_constraints_are_implemented() const
Definition: fe.cc:860
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
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other, const unsigned int face_no=0) const
Definition: fe.cc:1004
Table< 3, double > interior_weights
Table< 2, double > boundary_weights