Reference documentation for deal.II version 9.2.0
\(\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_dgq.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2001 - 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_dgq_h
17 #define dealii_fe_dgq_h
18 
19 #include <deal.II/base/config.h>
20 
23 
24 #include <deal.II/fe/fe_poly.h>
25 
27 
28 // Forward declarations
29 #ifndef DOXYGEN
30 template <int dim, int spacedim>
31 class MappingQ;
32 template <int dim>
33 class Quadrature;
34 #endif
35 
38 
111 template <int dim, int spacedim = dim>
112 class FE_DGQ : public FE_Poly<TensorProductPolynomials<dim>, dim, spacedim>
113 {
114 public:
121  FE_DGQ(const unsigned int p);
122 
128  virtual std::string
129  get_name() const override;
130 
140  virtual void
142  FullMatrix<double> &matrix) const override;
143 
155  virtual void
157  FullMatrix<double> &matrix) const override;
158 
170  virtual void
172  const unsigned int subface,
173  FullMatrix<double> &matrix) const override;
174 
192  virtual const FullMatrix<double> &
194  const unsigned int child,
195  const RefinementCase<dim> &refinement_case =
197 
219  virtual const FullMatrix<double> &
221  const unsigned int child,
222  const RefinementCase<dim> &refinement_case =
224 
248  virtual std::vector<std::pair<unsigned int, unsigned int>>
250  const FiniteElement<dim, spacedim> &fe_other) const override;
251 
259  virtual std::vector<std::pair<unsigned int, unsigned int>>
261  const FiniteElement<dim, spacedim> &fe_other) const override;
262 
270  virtual std::vector<std::pair<unsigned int, unsigned int>>
272  const FiniteElement<dim, spacedim> &fe_other) const override;
273 
282  virtual bool
283  hp_constraints_are_implemented() const override;
284 
290  const unsigned int codim = 0) const override final;
291 
300  virtual bool
301  has_support_on_face(const unsigned int shape_index,
302  const unsigned int face_index) const override;
303 
308  virtual std::pair<Table<2, bool>, std::vector<unsigned int>>
309  get_constant_modes() const override;
310 
318  virtual void
320  const std::vector<Vector<double>> &support_point_values,
321  std::vector<double> & nodal_values) const override;
322 
323  virtual std::unique_ptr<FiniteElement<dim, spacedim>>
324  clone() const override;
325 
326 protected:
335  FE_DGQ(const std::vector<Polynomials::Polynomial<double>> &polynomials);
336 
337 private:
344  static std::vector<unsigned int>
345  get_dpo_vector(const unsigned int degree);
346 
363  void
364  rotate_indices(std::vector<unsigned int> &indices,
365  const char direction) const;
366 
367  /*
368  * Mutex for protecting initialization of restriction and embedding matrix.
369  */
371 
372  // Allow access from other dimensions.
373  template <int dim1, int spacedim1>
374  friend class FE_DGQ;
375 
376  // Allow @p MappingQ class to access to build_renumbering function.
377  template <int dim1, int spacedim1>
378  friend class MappingQ;
379 };
380 
381 
382 
399 template <int dim, int spacedim = dim>
400 class FE_DGQArbitraryNodes : public FE_DGQ<dim, spacedim>
401 {
402 public:
409  FE_DGQArbitraryNodes(const Quadrature<1> &points);
410 
416  virtual std::string
417  get_name() const override;
418 
426  virtual void
428  const std::vector<Vector<double>> &support_point_values,
429  std::vector<double> & nodal_values) const override;
430  virtual std::unique_ptr<FiniteElement<dim, spacedim>>
431  clone() const override;
432 };
433 
434 
435 
446 template <int dim, int spacedim = dim>
447 class FE_DGQLegendre : public FE_DGQ<dim, spacedim>
448 {
449 public:
454  FE_DGQLegendre(const unsigned int degree);
455 
461  virtual std::pair<Table<2, bool>, std::vector<unsigned int>>
462  get_constant_modes() const override;
463 
470  virtual std::string
471  get_name() const override;
472 
473  virtual std::unique_ptr<FiniteElement<dim, spacedim>>
474  clone() const override;
475 };
476 
477 
478 
497 template <int dim, int spacedim = dim>
498 class FE_DGQHermite : public FE_DGQ<dim, spacedim>
499 {
500 public:
505  FE_DGQHermite(const unsigned int degree);
506 
513  virtual std::string
514  get_name() const override;
515 
516  virtual std::unique_ptr<FiniteElement<dim, spacedim>>
517  clone() const override;
518 };
519 
520 
524 
525 #endif
tensor_product_polynomials.h
MappingQ
Definition: mapping_manifold.h:33
FE_DGQ
Definition: fe_dgq.h:112
Threads::Mutex
Definition: thread_management.h:91
FE_DGQLegendre::get_name
virtual std::string get_name() const override
Definition: fe_dgq.cc:987
FE_DGQ::has_support_on_face
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
Definition: fe_dgq.cc:709
FE_DGQ::get_restriction_matrix
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const override
Definition: fe_dgq.cc:474
fe_poly.h
FE_DGQHermite
Definition: fe_dgq.h:498
FE_DGQ::rotate_indices
void rotate_indices(std::vector< unsigned int > &indices, const char direction) const
Definition: fe_dgq.cc:190
FE_DGQ::clone
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
Definition: fe_dgq.cc:164
FE_DGQLegendre
Definition: fe_dgq.h:447
FE_DGQArbitraryNodes::get_name
virtual std::string get_name() const override
Definition: fe_dgq.cc:827
FiniteElementData::degree
const unsigned int degree
Definition: fe_base.h:298
thread_management.h
FE_DGQ::hp_line_dof_identities
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
Definition: fe_dgq.cc:572
FE_DGQArbitraryNodes::FE_DGQArbitraryNodes
FE_DGQArbitraryNodes(const Quadrature< 1 > &points)
Definition: fe_dgq.cc:813
FiniteElementDomination::Domination
Domination
Definition: fe_base.h:83
FE_DGQLegendre::clone
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
Definition: fe_dgq.cc:997
FE_DGQArbitraryNodes::clone
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
Definition: fe_dgq.cc:945
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
FE_DGQ::get_name
virtual std::string get_name() const override
Definition: fe_dgq.cc:127
RefinementCase
Definition: geometry_info.h:795
LAPACKSupport::matrix
@ matrix
Contents is actually a matrix.
Definition: lapack_support.h:60
FE_DGQ::get_subface_interpolation_matrix
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const override
Definition: fe_dgq.cc:372
FE_DGQ::hp_vertex_dof_identities
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
Definition: fe_dgq.cc:558
Polynomials::Polynomial< double >
FiniteElement
Definition: fe.h:648
FE_DGQ::compare_for_domination
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const override final
Definition: fe_dgq.cc:600
FE_DGQArbitraryNodes
Definition: fe_dgq.h:400
FE_DGQLegendre::get_constant_modes
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
Definition: fe_dgq.cc:973
FE_DGQ::get_prolongation_matrix
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const override
Definition: fe_dgq.cc:399
FE_DGQ::convert_generalized_support_point_values_to_dof_values
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
Definition: fe_dgq.cc:143
FE_DGQ::get_face_interpolation_matrix
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const override
Definition: fe_dgq.cc:346
FE_Poly
Definition: fe_poly.h:76
FE_DGQHermite::get_name
virtual std::string get_name() const override
Definition: fe_dgq.cc:1016
FE_DGQ::get_dpo_vector
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
Definition: fe_dgq.cc:177
config.h
FE_DGQHermite::FE_DGQHermite
FE_DGQHermite(const unsigned int degree)
Definition: fe_dgq.cc:1007
FullMatrix< double >
FE_DGQ::get_interpolation_matrix
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const override
Definition: fe_dgq.cc:267
FE_DGQLegendre::FE_DGQLegendre
FE_DGQLegendre(const unsigned int degree)
Definition: fe_dgq.cc:964
FE_DGQHermite::clone
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
Definition: fe_dgq.cc:1026
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
FE_DGQ::hp_constraints_are_implemented
virtual bool hp_constraints_are_implemented() const override
Definition: fe_dgq.cc:549
FE_DGQArbitraryNodes::convert_generalized_support_point_values_to_dof_values
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
Definition: fe_dgq.cc:924
Quadrature
Definition: quadrature.h:85
Vector< double >
FE_DGQ::mutex
Threads::Mutex mutex
Definition: fe_dgq.h:370
FE_DGQ::get_constant_modes
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
Definition: fe_dgq.cc:800
FE_DGQ::FE_DGQ
friend class FE_DGQ
Definition: fe_dgq.h:374
FE_DGQ::hp_quad_dof_identities
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
Definition: fe_dgq.cc:586