Reference documentation for deal.II version 8.5.1
fe_q_base.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 2016 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_q_base_h
17 #define dealii__fe_q_base_h
18 
19 #include <deal.II/base/config.h>
20 #include <deal.II/fe/fe_poly.h>
21 #include <deal.II/base/thread_management.h>
22 
23 DEAL_II_NAMESPACE_OPEN
24 
25 
28 
39 template <class PolynomialType, int dim=PolynomialType::dimension, int spacedim=dim>
40 class FE_Q_Base : public FE_Poly<PolynomialType,dim,spacedim>
41 {
42 public:
46  FE_Q_Base (const PolynomialType &poly_space,
47  const FiniteElementData<dim> &fe_data,
48  const std::vector<bool> &restriction_is_additive_flags);
49 
59  virtual void
61  FullMatrix<double> &matrix) const;
62 
63 
73  virtual void
75  FullMatrix<double> &matrix) const;
76 
86  virtual void
88  const unsigned int subface,
89  FullMatrix<double> &matrix) const;
90 
95  virtual bool has_support_on_face (const unsigned int shape_index,
96  const unsigned int face_index) const;
97 
120  virtual const FullMatrix<double> &
121  get_restriction_matrix (const unsigned int child,
123 
150  virtual const FullMatrix<double> &
151  get_prolongation_matrix (const unsigned int child,
153 
192  virtual
193  unsigned int face_to_cell_index (const unsigned int face_dof_index,
194  const unsigned int face,
195  const bool face_orientation = true,
196  const bool face_flip = false,
197  const bool face_rotation = false) const;
198 
203  virtual std::pair<Table<2,bool>, std::vector<unsigned int> >
204  get_constant_modes () const;
205 
219  virtual bool hp_constraints_are_implemented () const;
220 
236  virtual
237  std::vector<std::pair<unsigned int, unsigned int> >
239 
244  virtual
245  std::vector<std::pair<unsigned int, unsigned int> >
246  hp_line_dof_identities (const FiniteElement<dim,spacedim> &fe_other) const;
247 
252  virtual
253  std::vector<std::pair<unsigned int, unsigned int> >
254  hp_quad_dof_identities (const FiniteElement<dim,spacedim> &fe_other) const;
255 
265  virtual
269 
276  "FE_Q can only be used for polynomial degrees "
277  "greater than zero. If you want an element of polynomial "
278  "degree zero, then it cannot be continuous and you "
279  "will want to use FE_DGQ<dim>(0).");
280 
281 protected:
288  static std::vector<unsigned int> get_dpo_vector(const unsigned int degree);
289 
295  void initialize (const std::vector<Point<1> > &support_points_1d);
296 
301  void initialize_constraints (const std::vector<Point<1> > &points);
302 
307  void initialize_unit_support_points (const std::vector<Point<1> > &points);
308 
313  void initialize_unit_face_support_points (const std::vector<Point<1> > &points);
314 
320 
327  struct Implementation;
328 
329  /*
330  * Declare implementation friend.
331  */
332  friend struct FE_Q_Base<PolynomialType,dim,spacedim>::Implementation;
333 
334 private:
335  /*
336  * Mutex for protecting initialization of restriction and embedding matrix.
337  */
338  mutable Threads::Mutex mutex;
339 
340  /*
341  * The highest polynomial degree of the underlying tensor product space
342  * without any enrichment. For FE_Q*(p) this is p. Note that enrichments
343  * may lead to a difference to degree.
344  */
345  const unsigned int q_degree;
346 };
347 
348 
351 DEAL_II_NAMESPACE_CLOSE
352 
353 #endif
FE_Q_Base(const PolynomialType &poly_space, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags)
Definition: fe_q_base.cc:415
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
Definition: fe_q_base.cc:1132
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_q_base.cc:684
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe_q_base.cc:568
const unsigned int degree
Definition: fe_base.h:311
static ::ExceptionBase & ExcFEQCannotHaveDegree0()
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_q_base.cc:788
Definition: point.h:89
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
Definition: fe_q_base.cc:1157
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
Definition: fe_q_base.cc:1351
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:553
virtual bool hp_constraints_are_implemented() const
Definition: fe_q_base.cc:673
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
Definition: fe_q_base.cc:1012
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_q_base.cc:725
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
Definition: fe_q_base.cc:1490
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_q_base.cc:855
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
Definition: fe_q_base.cc:1590
void initialize(const std::vector< Point< 1 > > &support_points_1d)
Definition: fe_q_base.cc:430
void initialize_unit_face_support_points(const std::vector< Point< 1 > > &points)
Definition: fe_q_base.cc:920
void initialize_constraints(const std::vector< Point< 1 > > &points)
Definition: fe_q_base.cc:1147
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
Definition: fe_q_base.cc:581
void initialize_quad_dof_index_permutation()
Definition: fe_q_base.cc:944
void initialize_unit_support_points(const std::vector< Point< 1 > > &points)
Definition: fe_q_base.cc:903
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe_q_base.cc:474
PolynomialType poly_space
Definition: fe_poly.h:444
const std::vector< bool > restriction_is_additive_flags
Definition: fe.h:2346