Reference documentation for deal.II version 9.4.1
\(\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 - 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_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
30template <int dim, int spacedim>
31class MappingQ;
32template <int dim>
33class Quadrature;
34#endif
35
38
109template <int dim, int spacedim = dim>
110class FE_DGQ : public FE_Poly<dim, spacedim>
111{
112public:
119 FE_DGQ(const unsigned int p);
120
126 virtual std::string
127 get_name() const override;
128
138 virtual void
140 FullMatrix<double> &matrix) const override;
141
153 virtual void
155 FullMatrix<double> & matrix,
156 const unsigned int face_no = 0) const override;
157
169 virtual void
171 const FiniteElement<dim, spacedim> &source,
172 const unsigned int subface,
173 FullMatrix<double> & matrix,
174 const unsigned int face_no = 0) const override;
175
193 virtual const FullMatrix<double> &
195 const unsigned int child,
196 const RefinementCase<dim> &refinement_case =
198
220 virtual const FullMatrix<double> &
222 const unsigned int child,
223 const RefinementCase<dim> &refinement_case =
225
249 virtual std::vector<std::pair<unsigned int, unsigned int>>
251 const FiniteElement<dim, spacedim> &fe_other) const override;
252
260 virtual std::vector<std::pair<unsigned int, unsigned int>>
262 const FiniteElement<dim, spacedim> &fe_other) const override;
263
271 virtual std::vector<std::pair<unsigned int, unsigned int>>
273 const unsigned int face_no = 0) const override;
274
283 virtual bool
284 hp_constraints_are_implemented() const override;
285
291 const unsigned int codim = 0) const override final;
292
301 virtual bool
302 has_support_on_face(const unsigned int shape_index,
303 const unsigned int face_index) const override;
304
309 virtual std::pair<Table<2, bool>, std::vector<unsigned int>>
310 get_constant_modes() const override;
311
319 virtual void
321 const std::vector<Vector<double>> &support_point_values,
322 std::vector<double> & nodal_values) const override;
323
324 virtual std::unique_ptr<FiniteElement<dim, spacedim>>
325 clone() const override;
326
327protected:
336 FE_DGQ(const std::vector<Polynomials::Polynomial<double>> &polynomials);
337
338private:
345 static std::vector<unsigned int>
346 get_dpo_vector(const unsigned int degree);
347
364 void
365 rotate_indices(std::vector<unsigned int> &indices,
366 const char direction) const;
367
368 /*
369 * Mutex for protecting initialization of restriction and embedding matrix.
370 */
372
373 // Allow access from other dimensions.
374 template <int dim1, int spacedim1>
375 friend class FE_DGQ;
376
377 // Allow @p MappingQ class to access to build_renumbering function.
378 template <int dim1, int spacedim1>
379 friend class MappingQ;
380};
381
382
383
398template <int dim, int spacedim = dim>
399class FE_DGQArbitraryNodes : public FE_DGQ<dim, spacedim>
400{
401public:
408 FE_DGQArbitraryNodes(const Quadrature<1> &points);
409
415 virtual std::string
416 get_name() const override;
417
425 virtual void
427 const std::vector<Vector<double>> &support_point_values,
428 std::vector<double> & nodal_values) const override;
429 virtual std::unique_ptr<FiniteElement<dim, spacedim>>
430 clone() const override;
431};
432
433
434
450template <int dim, int spacedim = dim>
451class FE_DGQLegendre : public FE_DGQ<dim, spacedim>
452{
453public:
458 FE_DGQLegendre(const unsigned int degree);
459
465 virtual std::pair<Table<2, bool>, std::vector<unsigned int>>
466 get_constant_modes() const override;
467
474 virtual std::string
475 get_name() const override;
476
477 virtual std::unique_ptr<FiniteElement<dim, spacedim>>
478 clone() const override;
479};
480
481
482
499template <int dim, int spacedim = dim>
500class FE_DGQHermite : public FE_DGQ<dim, spacedim>
501{
502public:
507 FE_DGQHermite(const unsigned int degree);
508
515 virtual std::string
516 get_name() const override;
517
518 virtual std::unique_ptr<FiniteElement<dim, spacedim>>
519 clone() const override;
520};
521
522
526
527#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
Definition: fe_dgq.h:111
friend class FE_DGQ
Definition: fe_dgq.h:375
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:371
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:449
Definition: vector.h:109
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443