deal.II version GIT relicensing-2289-g1e5549a87a 2024-12-21 21:30:00+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_dgp.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2002 - 2024 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
16#include <deal.II/fe/fe_dgp.h>
18#include <deal.II/fe/fe_tools.h>
19
20#include <memory>
21#include <sstream>
22
23
25
26template <int dim, int spacedim>
27FE_DGP<dim, spacedim>::FE_DGP(const unsigned int degree)
28 : FE_Poly<dim, spacedim>(
29 PolynomialSpace<dim>(
30 Polynomials::Legendre::generate_complete_basis(degree)),
31 FiniteElementData<dim>(get_dpo_vector(degree),
32 1,
33 degree,
34 FiniteElementData<dim>::L2),
35 std::vector<bool>(
36 FiniteElementData<dim>(get_dpo_vector(degree), 1, degree)
37 .n_dofs_per_cell(),
38 true),
39 std::vector<ComponentMask>(
40 FiniteElementData<dim>(get_dpo_vector(degree), 1, degree)
41 .n_dofs_per_cell(),
42 ComponentMask(std::vector<bool>(1, true))))
43{
44 // Reinit the vectors of restriction and prolongation matrices to the right
45 // sizes
47 // Fill prolongation matrices with embedding operators
48 if (dim == spacedim)
49 {
51 // Fill restriction matrices with L2-projection
53 }
54}
55
56
57template <int dim, int spacedim>
58std::string
60{
61 // note that the FETools::get_fe_by_name function depends on the
62 // particular format of the string this function returns, so they have to be
63 // kept in sync
64
65 std::ostringstream namebuf;
66 namebuf << "FE_DGP<" << Utilities::dim_string(dim, spacedim) << ">("
67 << this->degree << ")";
68
69 return namebuf.str();
70}
71
72
73
74template <int dim, int spacedim>
75std::unique_ptr<FiniteElement<dim, spacedim>>
77{
78 return std::make_unique<FE_DGP<dim, spacedim>>(*this);
79}
80
81
82
83//---------------------------------------------------------------------------
84// Auxiliary functions
85//---------------------------------------------------------------------------
86
87
88template <int dim, int spacedim>
89std::vector<unsigned int>
91{
92 std::vector<unsigned int> dpo(dim + 1, 0U);
93 dpo[dim] = deg + 1;
94 for (unsigned int i = 1; i < dim; ++i)
95 {
96 dpo[dim] *= deg + 1 + i;
97 dpo[dim] /= i + 1;
98 }
99 return dpo;
100}
101
102
103
104template <int dim, int spacedim>
105void
107 const FiniteElement<dim, spacedim> &x_source_fe,
108 FullMatrix<double> &interpolation_matrix,
109 const unsigned int) const
110{
111 // this is only implemented, if the source FE is also a DGP element. in that
112 // case, both elements have no dofs on their faces and the face
113 // interpolation matrix is necessarily empty -- i.e. there isn't much we
114 // need to do here.
115 (void)interpolation_matrix;
117 using FEDGP = FE_DGP<dim, spacedim>;
118 AssertThrow((x_source_fe.get_name().find("FE_DGP<") == 0) ||
119 (dynamic_cast<const FEDGP *>(&x_source_fe) != nullptr),
120 typename FE::ExcInterpolationNotImplemented());
121
122 Assert(interpolation_matrix.m() == 0,
123 ExcDimensionMismatch(interpolation_matrix.m(), 0));
124 Assert(interpolation_matrix.n() == 0,
125 ExcDimensionMismatch(interpolation_matrix.n(), 0));
126}
127
128
129
130template <int dim, int spacedim>
131void
133 const FiniteElement<dim, spacedim> &x_source_fe,
134 const unsigned int,
135 FullMatrix<double> &interpolation_matrix,
136 const unsigned int) const
137{
138 // this is only implemented, if the source FE is also a DGP element. in that
139 // case, both elements have no dofs on their faces and the face
140 // interpolation matrix is necessarily empty -- i.e. there isn't much we
141 // need to do here.
142 (void)interpolation_matrix;
144 using FEDGP = FE_DGP<dim, spacedim>;
145 AssertThrow((x_source_fe.get_name().find("FE_DGP<") == 0) ||
146 (dynamic_cast<const FEDGP *>(&x_source_fe) != nullptr),
147 typename FE::ExcInterpolationNotImplemented());
148
149 Assert(interpolation_matrix.m() == 0,
150 ExcDimensionMismatch(interpolation_matrix.m(), 0));
151 Assert(interpolation_matrix.n() == 0,
152 ExcDimensionMismatch(interpolation_matrix.n(), 0));
153}
154
155
156
157template <int dim, int spacedim>
158bool
163
164
165
166template <int dim, int spacedim>
167std::vector<std::pair<unsigned int, unsigned int>>
169 const FiniteElement<dim, spacedim> & /*fe_other*/) const
170{
171 // this element is discontinuous, so by definition there can
172 // be no identities between its dofs and those of any neighbor
173 // (of whichever type the neighbor may be -- after all, we have
174 // no face dofs on this side to begin with)
175 return std::vector<std::pair<unsigned int, unsigned int>>();
176}
177
178
179
180template <int dim, int spacedim>
181std::vector<std::pair<unsigned int, unsigned int>>
183 const FiniteElement<dim, spacedim> &fe_other) const
184{
185 // there are no such constraints for DGP elements at all
186 if (dynamic_cast<const FE_DGP<dim, spacedim> *>(&fe_other) != nullptr)
187 return std::vector<std::pair<unsigned int, unsigned int>>();
188 else
189 {
191 return std::vector<std::pair<unsigned int, unsigned int>>();
192 }
193}
194
195
196
197template <int dim, int spacedim>
198std::vector<std::pair<unsigned int, unsigned int>>
200 const FiniteElement<dim, spacedim> &fe_other,
201 const unsigned int) const
202{
203 // there are no such constraints for DGP elements at all
204 if (dynamic_cast<const FE_DGP<dim, spacedim> *>(&fe_other) != nullptr)
205 return std::vector<std::pair<unsigned int, unsigned int>>();
206 else
207 {
209 return std::vector<std::pair<unsigned int, unsigned int>>();
210 }
211}
212
213
214
215template <int dim, int spacedim>
218 const FiniteElement<dim, spacedim> &fe_other,
219 const unsigned int codim) const
220{
221 Assert(codim <= dim, ExcImpossibleInDim(dim));
222
223 // vertex/line/face domination
224 // ---------------------------
225 if (codim > 0)
226 // this is a discontinuous element, so by definition there will
227 // be no constraints wherever this element comes together with
228 // any other kind of element
230
231 // cell domination
232 // ---------------
233 if (const FE_DGP<dim, spacedim> *fe_dgp_other =
234 dynamic_cast<const FE_DGP<dim, spacedim> *>(&fe_other))
235 {
236 if (this->degree < fe_dgp_other->degree)
238 else if (this->degree == fe_dgp_other->degree)
240 else
242 }
243 else if (const FE_Nothing<dim> *fe_nothing =
244 dynamic_cast<const FE_Nothing<dim> *>(&fe_other))
245 {
246 if (fe_nothing->is_dominating())
248 else
249 // the FE_Nothing has no degrees of freedom and it is typically used
250 // in a context where we don't require any continuity along the
251 // interface
253 }
254
257}
258
259
260
261template <int dim, int spacedim>
262bool
264 const unsigned int) const
265{
266 // all shape functions have support on all faces
267 return true;
268}
269
270
271
272template <int dim, int spacedim>
273std::pair<Table<2, bool>, std::vector<unsigned int>>
275{
276 Table<2, bool> constant_modes(1, this->n_dofs_per_cell());
277 constant_modes(0, 0) = true;
278 return std::pair<Table<2, bool>, std::vector<unsigned int>>(
279 constant_modes, std::vector<unsigned int>(1, 0));
280}
281
282
283
284template <int dim, int spacedim>
285std::size_t
291
292
293
294// explicit instantiations
295#include "fe_dgp.inst"
296
297
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
Definition fe_dgp.cc:182
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
Definition fe_dgp.cc:263
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
Definition fe_dgp.cc:274
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_dgp.cc:132
virtual std::string get_name() const override
Definition fe_dgp.cc:59
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
Definition fe_dgp.cc:90
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
Definition fe_dgp.cc:168
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
Definition fe_dgp.cc:106
FE_DGP(const unsigned int p)
Definition fe_dgp.cc:27
virtual std::size_t memory_consumption() const override
Definition fe_dgp.cc:286
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_dgp.cc:199
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
Definition fe_dgp.cc:76
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const override final
Definition fe_dgp.cc:217
virtual bool hp_constraints_are_implemented() const override
Definition fe_dgp.cc:159
virtual std::string get_name() const =0
std::vector< std::vector< FullMatrix< double > > > restriction
Definition fe.h:2450
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
std::vector< std::vector< FullMatrix< double > > > prolongation
Definition fe.h:2464
size_type n() const
size_type m() const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
#define AssertThrow(cond, exc)
#define DEAL_II_NOT_IMPLEMENTED()
void compute_embedding_matrices(const FiniteElement< dim, spacedim > &fe, std::vector< std::vector< FullMatrix< number > > > &matrices, const bool isotropic_only=false, const double threshold=1.e-12)
void compute_projection_matrices(const FiniteElement< dim, spacedim > &fe, std::vector< std::vector< FullMatrix< number > > > &matrices, const bool isotropic_only=false)
std::string dim_string(const int dim, const int spacedim)
Definition utilities.cc:555
STL namespace.