Reference documentation for deal.II version GIT relicensing-249-g48dc7357c7 2024-03-29 12:30:02+00:00
\(\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// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2001 - 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_fe_dgq_h
16#define dealii_fe_dgq_h
17
18#include <deal.II/base/config.h>
19
20#include <deal.II/base/mutex.h>
22
23#include <deal.II/fe/fe_poly.h>
24
26
27// Forward declarations
28#ifndef DOXYGEN
29template <int dim, int spacedim>
30class MappingQ;
31template <int dim>
32class Quadrature;
33#endif
34
110template <int dim, int spacedim = dim>
111class FE_DGQ : public FE_Poly<dim, spacedim>
112{
113public:
120 FE_DGQ(const unsigned int p);
121
127 virtual std::string
128 get_name() const override;
129
139 virtual void
141 FullMatrix<double> &matrix) const override;
142
154 virtual void
156 FullMatrix<double> &matrix,
157 const unsigned int face_no = 0) const override;
158
170 virtual void
172 const FiniteElement<dim, spacedim> &source,
173 const unsigned int subface,
174 FullMatrix<double> &matrix,
175 const unsigned int face_no = 0) const override;
176
194 virtual const FullMatrix<double> &
196 const unsigned int child,
197 const RefinementCase<dim> &refinement_case =
199
221 virtual const FullMatrix<double> &
223 const unsigned int child,
224 const RefinementCase<dim> &refinement_case =
226
250 virtual std::vector<std::pair<unsigned int, unsigned int>>
252 const FiniteElement<dim, spacedim> &fe_other) const override;
253
261 virtual std::vector<std::pair<unsigned int, unsigned int>>
263 const FiniteElement<dim, spacedim> &fe_other) const override;
264
272 virtual std::vector<std::pair<unsigned int, unsigned int>>
274 const unsigned int face_no = 0) const override;
275
284 virtual bool
285 hp_constraints_are_implemented() const override;
286
292 const unsigned int codim = 0) const override final;
293
302 virtual bool
303 has_support_on_face(const unsigned int shape_index,
304 const unsigned int face_index) const override;
305
310 virtual std::pair<Table<2, bool>, std::vector<unsigned int>>
311 get_constant_modes() const override;
312
320 virtual void
322 const std::vector<Vector<double>> &support_point_values,
323 std::vector<double> &nodal_values) const override;
324
325 virtual std::unique_ptr<FiniteElement<dim, spacedim>>
326 clone() const override;
327
328protected:
337 FE_DGQ(const std::vector<Polynomials::Polynomial<double>> &polynomials);
338
339private:
346 static std::vector<unsigned int>
347 get_dpo_vector(const unsigned int degree);
348
365 void
366 rotate_indices(std::vector<unsigned int> &indices,
367 const char direction) const;
368
375
376 // Allow access from other dimensions.
377 template <int dim1, int spacedim1>
378 friend class FE_DGQ;
379
380 // Allow @p MappingQ class to access to build_renumbering function.
381 template <int dim1, int spacedim1>
382 friend class MappingQ;
383};
384
385
386
401template <int dim, int spacedim = dim>
402class FE_DGQArbitraryNodes : public FE_DGQ<dim, spacedim>
403{
404public:
411 FE_DGQArbitraryNodes(const Quadrature<1> &points);
412
418 virtual std::string
419 get_name() const override;
420
428 virtual void
430 const std::vector<Vector<double>> &support_point_values,
431 std::vector<double> &nodal_values) const override;
432 virtual std::unique_ptr<FiniteElement<dim, spacedim>>
433 clone() const override;
434};
435
436
437
453template <int dim, int spacedim = dim>
454class FE_DGQLegendre : public FE_DGQ<dim, spacedim>
455{
456public:
461 FE_DGQLegendre(const unsigned int degree);
462
468 virtual std::pair<Table<2, bool>, std::vector<unsigned int>>
469 get_constant_modes() const override;
470
477 virtual std::string
478 get_name() const override;
479
480 virtual std::unique_ptr<FiniteElement<dim, spacedim>>
481 clone() const override;
482};
483
484
485
502template <int dim, int spacedim = dim>
503class FE_DGQHermite : public FE_DGQ<dim, spacedim>
504{
505public:
510 FE_DGQHermite(const unsigned int degree);
511
518 virtual std::string
519 get_name() const override;
520
521 virtual std::unique_ptr<FiniteElement<dim, spacedim>>
522 clone() const override;
523};
524
525
529
530#endif
virtual std::string get_name() const override
Definition fe_dgq.cc:865
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:965
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
Definition fe_dgq.cc:986
virtual std::string get_name() const override
Definition fe_dgq.cc:1059
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
Definition fe_dgq.cc:1069
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
Definition fe_dgq.cc:1016
virtual std::string get_name() const override
Definition fe_dgq.cc:1030
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
Definition fe_dgq.cc:1040
friend class FE_DGQ
Definition fe_dgq.h:378
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
Definition fe_dgq.cc:747
virtual bool hp_constraints_are_implemented() const override
Definition fe_dgq.cc:576
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:613
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
Definition fe_dgq.cc:182
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
Definition fe_dgq.cc:838
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:500
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:369
Threads::Mutex restriction_matrix_mutex
Definition fe_dgq.h:373
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const override final
Definition fe_dgq.cc:628
Threads::Mutex prolongation_matrix_mutex
Definition fe_dgq.h:374
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:585
virtual std::string get_name() const override
Definition fe_dgq.cc:132
void rotate_indices(std::vector< unsigned int > &indices, const char direction) const
Definition fe_dgq.cc:195
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:599
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:148
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
Definition fe_dgq.cc:169
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:424
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:396
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const override
Definition fe_dgq.cc:272
const unsigned int degree
Definition fe_data.h:452
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503