Reference documentation for deal.II version GIT relicensing-464-g14f7274e4d 2024-04-22 15:40: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_nothing.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2009 - 2022 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
18#include <memory>
19
21
22
23template <int dim, int spacedim>
25 const unsigned int n_components,
26 const bool dominate)
27 : FiniteElement<dim, spacedim>(
28 FiniteElementData<dim>(std::vector<unsigned>(dim + 1, 0),
29 type,
30 n_components,
31 0,
32 FiniteElementData<dim>::unknown),
33 std::vector<bool>(),
34 std::vector<ComponentMask>())
35 , dominate(dominate)
36{
38 ExcMessage("A finite element needs to have at least one "
39 "vector component."));
40
41 // in most other elements we have to set up all sorts of stuff
42 // here. there isn't much that we have to do here; in particular,
43 // we can simply leave the restriction and prolongation matrices
44 // empty since their proper size is in fact zero given that the
45 // element here has no degrees of freedom
46}
47
48
49
50template <int dim, int spacedim>
51FE_Nothing<dim, spacedim>::FE_Nothing(const unsigned int n_components,
52 const bool dominate)
53 : FE_Nothing<dim, spacedim>(ReferenceCells::get_hypercube<dim>(),
54 n_components,
55 dominate)
56{}
57
58
59
60template <int dim, int spacedim>
61std::unique_ptr<FiniteElement<dim, spacedim>>
63{
64 return std::make_unique<FE_Nothing<dim, spacedim>>(*this);
65}
66
67
68
69template <int dim, int spacedim>
70std::string
72{
73 std::ostringstream namebuf;
74 namebuf << "FE_Nothing<" << Utilities::dim_string(dim, spacedim) << ">(";
75
76 std::vector<std::string> name_components;
77 if (this->reference_cell() != ReferenceCells::get_hypercube<dim>())
78 name_components.push_back(this->reference_cell().to_string());
79 if (this->n_components() > 1)
80 name_components.push_back(std::to_string(this->n_components()));
81 if (dominate)
82 name_components.emplace_back("dominating");
83
84 for (const std::string &comp : name_components)
85 {
86 namebuf << comp;
87 if (comp != name_components.back())
88 namebuf << ", ";
89 }
90 namebuf << ")";
91
92 return namebuf.str();
93}
94
95
96
97template <int dim, int spacedim>
100{
101 return flags;
102}
103
104
105
106template <int dim, int spacedim>
107double
109 const Point<dim> & /*p*/) const
110{
111 Assert(false, ExcMessage("This element has no shape functions."));
112 return 0;
113}
114
115
116
117template <int dim, int spacedim>
118std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
120 const UpdateFlags /*update_flags*/,
121 const Mapping<dim, spacedim> & /*mapping*/,
122 const Quadrature<dim> & /*quadrature*/,
124 spacedim>
125 & /*output_data*/) const
126{
127 // Create a default data object. Normally we would then
128 // need to resize things to hold the appropriate numbers
129 // of dofs, but in this case all data fields are empty.
130 return std::make_unique<
132}
133
134
135
136template <int dim, int spacedim>
137void
141 const Quadrature<dim> &,
147 spacedim>
148 &) const
150 // leave data fields empty
151}
152
153
154
155template <int dim, int spacedim>
156void
159 const unsigned int,
166 spacedim>
167 &) const
168{
169 // leave data fields empty
170}
171
172
173
174template <int dim, int spacedim>
175void
178 const unsigned int,
179 const unsigned int,
180 const Quadrature<dim - 1> &,
186 spacedim>
187 &) const
188{
189 // leave data fields empty
190}
191
193
194template <int dim, int spacedim>
195bool
197{
198 return dominate;
199}
200
201
202
203template <int dim, int spacedim>
207 const unsigned int codim) const
209 Assert(codim <= dim, ExcImpossibleInDim(dim));
210 (void)codim;
211
212 if (!dominate)
213 // if FE_Nothing does not dominate, there are no requirements
215 else if (dynamic_cast<const FE_Nothing<dim> *>(&fe) != nullptr)
216 // if it does and the other is FE_Nothing, either can dominate
218 else
219 // otherwise we dominate whatever FE is provided
221}
223
224
225template <int dim, int spacedim>
226std::vector<std::pair<unsigned int, unsigned int>>
228 const FiniteElement<dim, spacedim> & /*fe_other*/) const
229{
230 // the FE_Nothing has no
231 // degrees of freedom, so there
232 // are no equivalencies to be
233 // recorded
234 return std::vector<std::pair<unsigned int, unsigned int>>();
235}
236
237
238template <int dim, int spacedim>
239std::vector<std::pair<unsigned int, unsigned int>>
241 const FiniteElement<dim, spacedim> & /*fe_other*/) const
242{
243 // the FE_Nothing has no
244 // degrees of freedom, so there
245 // are no equivalencies to be
246 // recorded
247 return std::vector<std::pair<unsigned int, unsigned int>>();
248}
249
250
251template <int dim, int spacedim>
252std::vector<std::pair<unsigned int, unsigned int>>
254 const FiniteElement<dim, spacedim> & /*fe_other*/,
255 const unsigned int) const
256{
257 // the FE_Nothing has no
258 // degrees of freedom, so there
259 // are no equivalencies to be
260 // recorded
261 return std::vector<std::pair<unsigned int, unsigned int>>();
262}
263
264
265template <int dim, int spacedim>
266bool
271
273
274template <int dim, int spacedim>
275void
277 const FiniteElement<dim, spacedim> & /*source_fe*/,
278 FullMatrix<double> &interpolation_matrix) const
279{
280 // Since this element has no dofs,
281 // the interpolation matrix is necessarily empty.
282 (void)interpolation_matrix;
283
284 Assert(interpolation_matrix.m() == 0,
285 ExcDimensionMismatch(interpolation_matrix.m(), 0));
286 Assert(interpolation_matrix.n() == 0,
287 ExcDimensionMismatch(interpolation_matrix.n(), 0));
288}
289
290
291
292template <int dim, int spacedim>
293void
295 const FiniteElement<dim, spacedim> & /*source_fe*/,
296 FullMatrix<double> &interpolation_matrix,
297 const unsigned int) const
298{
299 // since this element has no face dofs, the
300 // interpolation matrix is necessarily empty
301 (void)interpolation_matrix;
302
303 Assert(interpolation_matrix.m() == 0,
304 ExcDimensionMismatch(interpolation_matrix.m(), 0));
305 Assert(interpolation_matrix.n() == 0,
306 ExcDimensionMismatch(interpolation_matrix.m(), 0));
307}
308
309
310
311template <int dim, int spacedim>
312void
314 const FiniteElement<dim, spacedim> & /*source_fe*/,
315 const unsigned int /*index*/,
316 FullMatrix<double> &interpolation_matrix,
317 const unsigned int) const
318{
319 // since this element has no face dofs, the
320 // interpolation matrix is necessarily empty
322 (void)interpolation_matrix;
323 Assert(interpolation_matrix.m() == 0,
324 ExcDimensionMismatch(interpolation_matrix.m(), 0));
325 Assert(interpolation_matrix.n() == 0,
326 ExcDimensionMismatch(interpolation_matrix.m(), 0));
327}
328
329
330
331template <int dim, int spacedim>
332std::pair<Table<2, bool>, std::vector<unsigned int>>
335 // since this element has no dofs, there are no constant modes
336 return {Table<2, bool>{}, std::vector<unsigned int>{}};
337}
338
339
341// explicit instantiations
342#include "fe_nothing.inst"
343
344
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
Definition fe_nothing.cc:71
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source_fe, const unsigned int index, FullMatrix< double > &interpolation_matrix, const unsigned int face_no=0) const override
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
Definition fe_nothing.cc:62
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source_fe, FullMatrix< double > &interpolation_matrix, const unsigned int face_no=0) const override
FE_Nothing(const ReferenceCell &type, const unsigned int n_components=1, const bool dominate=false)
Definition fe_nothing.cc:24
virtual void fill_fe_subface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
Definition fe_nothing.cc:99
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
virtual double shape_value(const unsigned int i, const Point< dim > &p) const override
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const hp::QCollection< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source_fe, FullMatrix< double > &interpolation_matrix) const override
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const override final
virtual void fill_fe_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
bool is_dominating() const
virtual std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
virtual bool hp_constraints_are_implemented() 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_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other, const unsigned int face_no=0) const override
unsigned int n_components() const
size_type n() const
size_type m() const
Abstract base class for mapping classes.
Definition mapping.h:318
Definition point.h:111
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
UpdateFlags
std::string dim_string(const int dim, const int spacedim)
Definition utilities.cc:555
STL namespace.