Reference documentation for deal.II version Git 3f1f337db3 2021-10-23 13:19:02 -0600
\(\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_q_base.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 2021 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_q_base_h
17 #define dealii_fe_q_base_h
18 
19 #include <deal.II/base/config.h>
20 
22 
23 #include <deal.II/fe/fe_poly.h>
24 
26 
27 
30 
37 template <int dim, int spacedim = dim>
38 class FE_Q_Base : public FE_Poly<dim, spacedim>
39 {
40 public:
45  const FiniteElementData<dim> & fe_data,
46  const std::vector<bool> & restriction_is_additive_flags);
47 
57  virtual void
59  FullMatrix<double> &matrix) const override;
60 
61 
71  virtual void
73  FullMatrix<double> & matrix,
74  const unsigned int face_no = 0) const override;
75 
85  virtual void
87  const FiniteElement<dim, spacedim> &source,
88  const unsigned int subface,
89  FullMatrix<double> & matrix,
90  const unsigned int face_no = 0) const override;
91 
96  virtual bool
97  has_support_on_face(const unsigned int shape_index,
98  const unsigned int face_index) const override;
99 
122  virtual const FullMatrix<double> &
124  const unsigned int child,
125  const RefinementCase<dim> &refinement_case =
127 
154  virtual const FullMatrix<double> &
156  const unsigned int child,
157  const RefinementCase<dim> &refinement_case =
159 
198  virtual unsigned int
199  face_to_cell_index(const unsigned int face_dof_index,
200  const unsigned int face,
201  const bool face_orientation = true,
202  const bool face_flip = false,
203  const bool face_rotation = false) const override;
204 
209  virtual std::pair<Table<2, bool>, std::vector<unsigned int>>
210  get_constant_modes() const override;
211 
225  virtual bool
226  hp_constraints_are_implemented() const override;
227 
243  virtual std::vector<std::pair<unsigned int, unsigned int>>
245  const FiniteElement<dim, spacedim> &fe_other) const override;
246 
251  virtual std::vector<std::pair<unsigned int, unsigned int>>
253  const FiniteElement<dim, spacedim> &fe_other) const override;
254 
259  virtual std::vector<std::pair<unsigned int, unsigned int>>
261  const unsigned int face_no = 0) const override;
262 
264 
271  "FE_Q can only be used for polynomial degrees "
272  "greater than zero. If you want an element of polynomial "
273  "degree zero, then it cannot be continuous and you "
274  "will want to use FE_DGQ<dim>(0).");
275 
276 protected:
283  static std::vector<unsigned int>
284  get_dpo_vector(const unsigned int degree);
285 
291  void
292  initialize(const std::vector<Point<1>> &support_points_1d);
293 
298  void
299  initialize_constraints(const std::vector<Point<1>> &points);
300 
305  void
306  initialize_unit_support_points(const std::vector<Point<1>> &points);
307 
312  void
313  initialize_unit_face_support_points(const std::vector<Point<1>> &points);
314 
319  void
321 
329 
330  // Declare implementation friend.
331  friend struct FE_Q_Base<dim, spacedim>::Implementation;
332 
333 private:
338 
344  const unsigned int q_degree;
345 };
346 
347 
351 
352 #endif
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
Definition: fe_q_base.cc:739
void initialize_unit_face_support_points(const std::vector< Point< 1 >> &points)
Definition: fe_q_base.cc:979
void initialize_unit_support_points(const std::vector< Point< 1 >> &points)
Definition: fe_q_base.cc:954
Contents is actually a matrix.
virtual bool hp_constraints_are_implemented() const override
Definition: fe_q_base.cc:730
const unsigned int degree
Definition: fe_base.h:436
void initialize_quad_dof_index_permutation()
Definition: fe_q_base.cc:1018
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const override
Definition: fe_q_base.cc:1441
void initialize(const std::vector< Point< 1 >> &support_points_1d)
Definition: fe_q_base.cc:435
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 override
Definition: fe_q_base.cc:626
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
Definition: fe_q_base.cc:783
FE_Q_Base(const ScalarPolynomialsBase< dim > &poly_space, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags)
Definition: fe_q_base.cc:416
void initialize_constraints(const std::vector< Point< 1 >> &points)
Definition: fe_q_base.cc:1232
virtual unsigned int face_to_cell_index(const unsigned int face_dof_index, const unsigned int face, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false) const override
Definition: fe_q_base.cc:1103
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:487
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const override
Definition: fe_q_base.cc:511
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:401
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
Definition: fe_q_base.cc:1589
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
Definition: fe_q_base.cc:611
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:400
Threads::Mutex mutex
Definition: fe_q_base.h:337
const unsigned int q_degree
Definition: fe_q_base.h:344
static ::ExceptionBase & ExcFEQCannotHaveDegree0()
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
Definition: fe_q_base.cc:1691
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 override
Definition: fe_q_base.cc:878
const std::unique_ptr< ScalarPolynomialsBase< dim > > poly_space
Definition: fe_poly.h:534
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
Definition: fe_q_base.cc:1218
const std::vector< bool > restriction_is_additive_flags
Definition: fe.h:2572
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const override
Definition: fe_q_base.cc:1242