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