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_poly_face.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2009 - 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#ifndef dealii_fe_poly_face_h
16#define dealii_fe_poly_face_h
17
18
19#include <deal.II/base/config.h>
20
22
23#include <deal.II/fe/fe.h>
24
25#include <memory>
26
28
57template <typename PolynomialType,
58 int dim = PolynomialType::dimension + 1,
59 int spacedim = dim>
60class FE_PolyFace : public FiniteElement<dim, spacedim>
61{
62public:
66 FE_PolyFace(const PolynomialType &poly_space,
67 const FiniteElementData<dim> &fe_data,
68 const std::vector<bool> &restriction_is_additive_flags);
69
74 unsigned int
75 get_degree() const;
76
77 // for documentation, see the FiniteElement base class
78 virtual UpdateFlags
79 requires_update_flags(const UpdateFlags update_flags) const override;
80
81protected:
82 /*
83 * NOTE: The following functions have their definitions inlined into the class
84 * declaration because we otherwise run into a compiler error with MS Visual
85 * Studio.
86 */
87
88
89 virtual std::unique_ptr<
92 const UpdateFlags update_flags,
93 const Mapping<dim, spacedim> &mapping,
94 const Quadrature<dim> &quadrature,
96 spacedim>
97 &output_data) const override
98 {
99 (void)update_flags;
100 (void)mapping;
101 (void)quadrature;
102 (void)output_data;
103 return std::make_unique<InternalData>();
104 }
105
106 using FiniteElement<dim, spacedim>::get_face_data;
107
108 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
110 const UpdateFlags update_flags,
111 const Mapping<dim, spacedim> &mapping,
112 const hp::QCollection<dim - 1> &quadrature,
114 spacedim>
115 &output_data) const override
116 {
117 (void)mapping;
118 (void)output_data;
119 AssertDimension(quadrature.size(), 1);
120
121 // generate a new data object and
122 // initialize some fields
123 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
124 data_ptr = std::make_unique<InternalData>();
125 auto &data = dynamic_cast<InternalData &>(*data_ptr);
126 data.update_each = requires_update_flags(update_flags);
127
128 const unsigned int n_q_points = quadrature[0].size();
129
130 // some scratch arrays
131 std::vector<double> values(0);
132 std::vector<Tensor<1, dim - 1>> grads(0);
133 std::vector<Tensor<2, dim - 1>> grad_grads(0);
134 std::vector<Tensor<3, dim - 1>> empty_vector_of_3rd_order_tensors;
135 std::vector<Tensor<4, dim - 1>> empty_vector_of_4th_order_tensors;
136
137 // initialize fields only if really
138 // necessary. otherwise, don't
139 // allocate memory
140 if (data.update_each & update_values)
141 {
142 values.resize(poly_space.n());
143 data.shape_values.resize(poly_space.n(),
144 std::vector<double>(n_q_points));
145 for (unsigned int i = 0; i < n_q_points; ++i)
146 {
147 poly_space.evaluate(quadrature[0].point(i),
148 values,
149 grads,
150 grad_grads,
151 empty_vector_of_3rd_order_tensors,
152 empty_vector_of_4th_order_tensors);
153
154 for (unsigned int k = 0; k < poly_space.n(); ++k)
155 data.shape_values[k][i] = values[k];
156 }
157 }
158 // No derivatives of this element
159 // are implemented.
160 if (data.update_each & update_gradients ||
161 data.update_each & update_hessians)
162 {
164 }
165
166 return data_ptr;
167 }
168
169 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
171 const UpdateFlags update_flags,
172 const Mapping<dim, spacedim> &mapping,
173 const Quadrature<dim - 1> &quadrature,
175 spacedim>
176 &output_data) const override
177 {
178 return get_face_data(
179 update_flags,
180 mapping,
183 output_data);
184 }
185
186 virtual void
189 const CellSimilarity::Similarity cell_similarity,
190 const Quadrature<dim> &quadrature,
191 const Mapping<dim, spacedim> &mapping,
192 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
194 &mapping_data,
195 const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
197 spacedim>
198 &output_data) const override;
199
200 using FiniteElement<dim, spacedim>::fill_fe_face_values;
201
202 virtual void
205 const unsigned int face_no,
206 const hp::QCollection<dim - 1> &quadrature,
207 const Mapping<dim, spacedim> &mapping,
208 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
210 &mapping_data,
211 const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
213 spacedim>
214 &output_data) const override;
215
216 virtual void
219 const unsigned int face_no,
220 const unsigned int sub_no,
221 const Quadrature<dim - 1> &quadrature,
222 const Mapping<dim, spacedim> &mapping,
223 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
225 &mapping_data,
226 const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
228 spacedim>
229 &output_data) const override;
230
237 class InternalData : public FiniteElement<dim, spacedim>::InternalDataBase
238 {
239 public:
253 std::vector<std::vector<double>> shape_values;
254 };
255
260 PolynomialType poly_space;
261};
262
266
267#endif
std::vector< std::vector< double > > shape_values
PolynomialType poly_space
std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim - 1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) 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
std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const hp::QCollection< dim - 1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
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
FE_PolyFace(const PolynomialType &poly_space, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags)
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
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 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
unsigned int get_degree() const
const std::vector< bool > restriction_is_additive_flags
Definition fe.h:2579
friend class InternalDataBase
Definition fe.h:3075
Abstract base class for mapping classes.
Definition mapping.h:318
unsigned int size() const
Definition collection.h:308
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
#define AssertDimension(dim1, dim2)
UpdateFlags
@ update_hessians
Second derivatives of shape functions.
@ update_values
Shape function values.
@ update_gradients
Shape function gradients.
#define DEAL_II_NOT_IMPLEMENTED()
constexpr const ReferenceCell & get_hypercube()