PolyDEAL
 
Loading...
Searching...
No Matches
fe_agglodgp.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
17#include <deal.II/fe/fe_tools.h>
18
19#include <fe_agglodgp.h>
20
21#include <memory>
22#include <sstream>
23
24
26
27template <int dim, int spacedim>
29 : FE_Poly<dim, spacedim>(
30 PolynomialSpace<dim>(
31 Polynomials::Legendre::generate_complete_basis(degree)),
32 FiniteElementData<dim>(get_dpo_vector(degree),
33 1,
34 degree,
35 FiniteElementData<dim>::L2),
36 std::vector<bool>(
37 FiniteElementData<dim>(get_dpo_vector(degree), 1, degree)
38 .n_dofs_per_cell(),
39 true),
40 std::vector<ComponentMask>(
41 FiniteElementData<dim>(get_dpo_vector(degree), 1, degree)
42 .n_dofs_per_cell(),
43 ComponentMask(std::vector<bool>(1, true))))
44{
45 // Reinit the vectors of restriction and prolongation matrices to the right
46 // sizes
48 // Fill prolongation matrices with embedding operators
49 if (dim == spacedim)
50 {
52 // Fill restriction matrices with L2-projection
54 }
55}
56
57
58template <int dim, int spacedim>
59std::string
61{
62 // note that the FETools::get_fe_by_name function depends on the
63 // particular format of the string this function returns, so they have to be
64 // kept in sync
65
66 std::ostringstream namebuf;
67 namebuf << "FE_AggloDGP<" << Utilities::dim_string(dim, spacedim) << ">("
68 << this->degree << ")";
69
70 return namebuf.str();
71}
72
73
74
75template <int dim, int spacedim>
76std::unique_ptr<FiniteElement<dim, spacedim>>
78{
79 return std::make_unique<FE_AggloDGP<dim, spacedim>>(*this);
80}
81
82
83
84//---------------------------------------------------------------------------
85// Auxiliary functions
86//---------------------------------------------------------------------------
87
88
89template <int dim, int spacedim>
90std::vector<unsigned int>
92{
93 std::vector<unsigned int> dpo(dim + 1, 0U);
94 dpo[dim] = deg + 1;
95 for (unsigned int i = 1; i < dim; ++i)
96 {
97 dpo[dim] *= deg + 1 + i;
98 dpo[dim] /= i + 1;
99 }
100 return dpo;
101}
102
103
104
105template <int dim, int spacedim>
106void
108 const FiniteElement<dim, spacedim> &x_source_fe,
109 FullMatrix<double> &interpolation_matrix,
110 const unsigned int) const
111{
112 // this is only implemented, if the source FE is also a DGP element. in that
113 // case, both elements have no dofs on their faces and the face
114 // interpolation matrix is necessarily empty -- i.e. there isn't much we
115 // need to do here.
116 (void)interpolation_matrix;
118 using FEDGP = FE_AggloDGP<dim, spacedim>;
119 AssertThrow((x_source_fe.get_name().find("FE_AggloDGP<") == 0) ||
120 (dynamic_cast<const FEDGP *>(&x_source_fe) != nullptr),
121 typename FE::ExcInterpolationNotImplemented());
122
123 Assert(interpolation_matrix.m() == 0,
124 ExcDimensionMismatch(interpolation_matrix.m(), 0));
125 Assert(interpolation_matrix.n() == 0,
126 ExcDimensionMismatch(interpolation_matrix.n(), 0));
127}
128
129
130
131template <int dim, int spacedim>
132void
134 const FiniteElement<dim, spacedim> &x_source_fe,
135 const unsigned int,
136 FullMatrix<double> &interpolation_matrix,
137 const unsigned int) const
138{
139 // this is only implemented, if the source FE is also a DGP element. in that
140 // case, both elements have no dofs on their faces and the face
141 // interpolation matrix is necessarily empty -- i.e. there isn't much we
142 // need to do here.
143 (void)interpolation_matrix;
145 using FEDGP = FE_AggloDGP<dim, spacedim>;
146 AssertThrow((x_source_fe.get_name().find("FE_AggloDGP<") == 0) ||
147 (dynamic_cast<const FEDGP *>(&x_source_fe) != nullptr),
148 typename FE::ExcInterpolationNotImplemented());
149
150 Assert(interpolation_matrix.m() == 0,
151 ExcDimensionMismatch(interpolation_matrix.m(), 0));
152 Assert(interpolation_matrix.n() == 0,
153 ExcDimensionMismatch(interpolation_matrix.n(), 0));
154}
155
156
157
158template <int dim, int spacedim>
159bool
164
165
166
167template <int dim, int spacedim>
168std::vector<std::pair<unsigned int, unsigned int>>
170 const FiniteElement<dim, spacedim> &fe_other) const
171{
172 (void)fe_other;
173 return std::vector<std::pair<unsigned int, unsigned int>>();
174}
175
176
177
178template <int dim, int spacedim>
179std::vector<std::pair<unsigned int, unsigned int>>
181 const FiniteElement<dim, spacedim> &fe_other) const
182{
183 // there are no such constraints for DGP elements at all
184 if (dynamic_cast<const FE_AggloDGP<dim, spacedim> *>(&fe_other) != nullptr)
185 return std::vector<std::pair<unsigned int, unsigned int>>();
186 else
187 {
189 return std::vector<std::pair<unsigned int, unsigned int>>();
190 }
191}
192
193
194
195template <int dim, int spacedim>
196std::vector<std::pair<unsigned int, unsigned int>>
198 const FiniteElement<dim, spacedim> &fe_other,
199 const unsigned int) const
200{
201 // there are no such constraints for DGP elements at all
202 if (dynamic_cast<const FE_AggloDGP<dim, spacedim> *>(&fe_other) != nullptr)
203 return std::vector<std::pair<unsigned int, unsigned int>>();
204 else
205 {
207 return std::vector<std::pair<unsigned int, unsigned int>>();
208 }
209}
210
211
212
213template <int dim, int spacedim>
216 const FiniteElement<dim, spacedim> &fe_other,
217 const unsigned int codim) const
218{
219 Assert(codim <= dim, ExcImpossibleInDim(dim));
220
221 // vertex/line/face domination
222 // ---------------------------
223 if (codim > 0)
224 // this is a discontinuous element, so by definition there will
225 // be no constraints wherever this element comes together with
226 // any other kind of element
228
229 // cell domination
230 // ---------------
231 if (const FE_AggloDGP<dim, spacedim> *fe_dgp_other =
232 dynamic_cast<const FE_AggloDGP<dim, spacedim> *>(&fe_other))
233 {
234 if (this->degree < fe_dgp_other->degree)
236 else if (this->degree == fe_dgp_other->degree)
238 else
240 }
241 else if (const FE_Nothing<dim> *fe_nothing =
242 dynamic_cast<const FE_Nothing<dim> *>(&fe_other))
243 {
244 if (fe_nothing->is_dominating())
246 else
247 // the FE_Nothing has no degrees of freedom and it is typically used
248 // in a context where we don't require any continuity along the
249 // interface
251 }
252
255}
256
257
258
259template <int dim, int spacedim>
260bool
262 const unsigned int) const
263{
264 // all shape functions have support on all faces
265 return true;
266}
267
268
269
270template <int dim, int spacedim>
271std::pair<Table<2, bool>, std::vector<unsigned int>>
273{
274 Table<2, bool> constant_modes(1, this->n_dofs_per_cell());
275 constant_modes(0, 0) = true;
276 return std::pair<Table<2, bool>, std::vector<unsigned int>>(
277 constant_modes, std::vector<unsigned int>(1, 0));
278}
279
280
281
282template <int dim, int spacedim>
283std::size_t
289
290
291
292// explicit instantiations
293template class FE_AggloDGP<1>;
294template class FE_AggloDGP<2>;
295template class FE_AggloDGP<3>;
296
297
virtual std::string get_name() const=0
std::vector< std::vector< FullMatrix< double > > > restriction
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
size_type n() const
size_type m() const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
virtual std::string get_name() const override
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
virtual bool hp_constraints_are_implemented() const override
virtual std::size_t memory_consumption() const override
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
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
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
FE_AggloDGP(const unsigned int p)
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
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const override final
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#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)