Reference documentation for deal.II version 9.6.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_dg_vector.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2010 - 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_dg_vector_h
16#define dealii_fe_dg_vector_h
17
18#include <deal.II/base/config.h>
19
25#include <deal.II/base/table.h>
27
28#include <deal.II/fe/fe.h>
30
31#include <vector>
32
34
51template <typename PolynomialType, int dim, int spacedim = dim>
52class FE_DGVector : public FE_PolyTensor<dim, spacedim>
53{
54public:
58 FE_DGVector(const unsigned int p, MappingKind m);
59
66 virtual std::string
67 get_name() const override;
68
69 virtual std::unique_ptr<FiniteElement<dim, spacedim>>
70 clone() const override;
71
78 virtual bool
79 has_support_on_face(const unsigned int shape_index,
80 const unsigned int face_index) const override;
81
82 virtual std::size_t
83 memory_consumption() const override;
84
85private:
92 static std::vector<unsigned int>
93 get_dpo_vector(const unsigned int degree);
94
102 {
103 public:
115 std::vector<std::vector<Tensor<1, dim>>> shape_values;
116
126 std::vector<std::vector<Tensor<2, dim>>> shape_gradients;
127 };
129};
130
131
132
142template <int dim, int spacedim = dim>
143class FE_DGNedelec : public FE_DGVector<PolynomialsNedelec<dim>, dim, spacedim>
144{
145public:
150 FE_DGNedelec(const unsigned int p);
151
157 virtual std::string
158 get_name() const override;
159};
160
161
162
173template <int dim, int spacedim = dim>
175 : public FE_DGVector<PolynomialsRaviartThomas<dim>, dim, spacedim>
176{
177public:
181 FE_DGRaviartThomas(const unsigned int p);
182
188 virtual std::string
189 get_name() const override;
190};
191
192
193
204template <int dim, int spacedim = dim>
205class FE_DGBDM : public FE_DGVector<PolynomialsBDM<dim>, dim, spacedim>
206{
207public:
211 FE_DGBDM(const unsigned int p);
212
218 virtual std::string
219 get_name() const override;
220};
221
222
224
225#endif
virtual std::string get_name() const override
FE_DGBDM(const unsigned int p)
virtual std::string get_name() const override
FE_DGNedelec(const unsigned int p)
virtual std::string get_name() const override
FE_DGRaviartThomas(const unsigned int p)
std::vector< std::vector< Tensor< 1, dim > > > shape_values
std::vector< std::vector< Tensor< 2, dim > > > shape_gradients
FE_DGVector(const unsigned int p, MappingKind m)
virtual std::string get_name() const override
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
Table< 3, double > interior_weights
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
virtual std::size_t memory_consumption() const override
const unsigned int degree
Definition fe_data.h:452
friend class InternalDataBase
Definition fe.h:3075
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
MappingKind
Definition mapping.h:79