Reference documentation for deal.II version 9.5.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\}}\)
Loading...
Searching...
No Matches
fe_dgq.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2001 - 2023 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
21#include <deal.II/base/mutex.h>
23
24#include <deal.II/fe/fe_poly.h>
25
27
28// Forward declarations
29#ifndef DOXYGEN
30template <int dim, int spacedim>
31class MappingQ;
32template <int dim>
33class Quadrature;
34#endif
35
111template <int dim, int spacedim = dim>
112class FE_DGQ : public FE_Poly<dim, spacedim>
113{
114public:
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,
158 const unsigned int face_no = 0) const override;
159
171 virtual void
173 const FiniteElement<dim, spacedim> &source,
174 const unsigned int subface,
175 FullMatrix<double> & matrix,
176 const unsigned int face_no = 0) const override;
177
195 virtual const FullMatrix<double> &
197 const unsigned int child,
198 const RefinementCase<dim> &refinement_case =
200
222 virtual const FullMatrix<double> &
224 const unsigned int child,
225 const RefinementCase<dim> &refinement_case =
227
251 virtual std::vector<std::pair<unsigned int, unsigned int>>
253 const FiniteElement<dim, spacedim> &fe_other) const override;
254
262 virtual std::vector<std::pair<unsigned int, unsigned int>>
264 const FiniteElement<dim, spacedim> &fe_other) const override;
265
273 virtual std::vector<std::pair<unsigned int, unsigned int>>
275 const unsigned int face_no = 0) const override;
276
285 virtual bool
286 hp_constraints_are_implemented() const override;
287
293 const unsigned int codim = 0) const override final;
294
303 virtual bool
304 has_support_on_face(const unsigned int shape_index,
305 const unsigned int face_index) const override;
306
311 virtual std::pair<Table<2, bool>, std::vector<unsigned int>>
312 get_constant_modes() const override;
313
321 virtual void
323 const std::vector<Vector<double>> &support_point_values,
324 std::vector<double> & nodal_values) const override;
325
326 virtual std::unique_ptr<FiniteElement<dim, spacedim>>
327 clone() const override;
328
329protected:
338 FE_DGQ(const std::vector<Polynomials::Polynomial<double>> &polynomials);
339
340private:
347 static std::vector<unsigned int>
348 get_dpo_vector(const unsigned int degree);
349
366 void
367 rotate_indices(std::vector<unsigned int> &indices,
368 const char direction) const;
369
370 /*
371 * Mutex for protecting initialization of restriction and embedding matrix.
372 */
374
375 // Allow access from other dimensions.
376 template <int dim1, int spacedim1>
377 friend class FE_DGQ;
378
379 // Allow @p MappingQ class to access to build_renumbering function.
380 template <int dim1, int spacedim1>
381 friend class MappingQ;
382};
383
384
385
400template <int dim, int spacedim = dim>
401class FE_DGQArbitraryNodes : public FE_DGQ<dim, spacedim>
402{
403public:
410 FE_DGQArbitraryNodes(const Quadrature<1> &points);
411
417 virtual std::string
418 get_name() const override;
419
427 virtual void
429 const std::vector<Vector<double>> &support_point_values,
430 std::vector<double> & nodal_values) const override;
431 virtual std::unique_ptr<FiniteElement<dim, spacedim>>
432 clone() const override;
433};
434
435
436
452template <int dim, int spacedim = dim>
453class FE_DGQLegendre : public FE_DGQ<dim, spacedim>
454{
455public:
460 FE_DGQLegendre(const unsigned int degree);
461
467 virtual std::pair<Table<2, bool>, std::vector<unsigned int>>
468 get_constant_modes() const override;
469
476 virtual std::string
477 get_name() const override;
478
479 virtual std::unique_ptr<FiniteElement<dim, spacedim>>
480 clone() const override;
481};
482
483
484
501template <int dim, int spacedim = dim>
502class FE_DGQHermite : public FE_DGQ<dim, spacedim>
503{
504public:
509 FE_DGQHermite(const unsigned int degree);
510
517 virtual std::string
518 get_name() const override;
519
520 virtual std::unique_ptr<FiniteElement<dim, spacedim>>
521 clone() const override;
522};
523
524
528
529#endif
virtual std::string get_name() const override
Definition fe_dgq.cc:851
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:951
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
Definition fe_dgq.cc:972
virtual std::string get_name() const override
Definition fe_dgq.cc:1045
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
Definition fe_dgq.cc:1055
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
Definition fe_dgq.cc:1002
virtual std::string get_name() const override
Definition fe_dgq.cc:1016
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
Definition fe_dgq.cc:1026
friend class FE_DGQ
Definition fe_dgq.h:377
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
Definition fe_dgq.cc:733
virtual bool hp_constraints_are_implemented() const override
Definition fe_dgq.cc:562
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_dgq.cc:599
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
Definition fe_dgq.cc:183
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
Definition fe_dgq.cc:824
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:486
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
Definition fe_dgq.cc:355
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const override final
Definition fe_dgq.cc:614
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:571
virtual std::string get_name() const override
Definition fe_dgq.cc:133
void rotate_indices(std::vector< unsigned int > &indices, const char direction) const
Definition fe_dgq.cc:196
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:585
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:149
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
Definition fe_dgq.cc:170
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:410
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_dgq.cc:382
Threads::Mutex mutex
Definition fe_dgq.h:373
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const override
Definition fe_dgq.cc:273
const unsigned int degree
Definition fe_data.h:453
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:472
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:473