Reference documentation for deal.II version GIT relicensing-399-g79d89019c5 2024-04-16 15:00: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_poly.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2004 - 2023 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_h
16#define dealii_fe_poly_h
17
18
19#include <deal.II/base/config.h>
20
23
24#include <deal.II/fe/fe.h>
25
26#include <memory>
27
29
74template <int dim, int spacedim = dim>
75class FE_Poly : public FiniteElement<dim, spacedim>
76{
77public:
82 const FiniteElementData<dim> &fe_data,
83 const std::vector<bool> &restriction_is_additive_flags,
84 const std::vector<ComponentMask> &nonzero_components);
85
89 FE_Poly(const FE_Poly &fe);
90
95 unsigned int
96 get_degree() const;
97
98 // for documentation, see the FiniteElement base class
99 virtual UpdateFlags
100 requires_update_flags(const UpdateFlags update_flags) const override;
101
107
117 std::vector<unsigned int>
119
126 std::vector<unsigned int>
128
134 virtual double
135 shape_value(const unsigned int i, const Point<dim> &p) const override;
136
147 virtual double
148 shape_value_component(const unsigned int i,
149 const Point<dim> &p,
150 const unsigned int component) const override;
151
157 virtual Tensor<1, dim>
158 shape_grad(const unsigned int i, const Point<dim> &p) const override;
159
170 virtual Tensor<1, dim>
171 shape_grad_component(const unsigned int i,
172 const Point<dim> &p,
173 const unsigned int component) const override;
174
180 virtual Tensor<2, dim>
181 shape_grad_grad(const unsigned int i, const Point<dim> &p) const override;
182
193 virtual Tensor<2, dim>
194 shape_grad_grad_component(const unsigned int i,
195 const Point<dim> &p,
196 const unsigned int component) const override;
197
203 virtual Tensor<3, dim>
204 shape_3rd_derivative(const unsigned int i,
205 const Point<dim> &p) const override;
206
217 virtual Tensor<3, dim>
218 shape_3rd_derivative_component(const unsigned int i,
219 const Point<dim> &p,
220 const unsigned int component) const override;
221
227 virtual Tensor<4, dim>
228 shape_4th_derivative(const unsigned int i,
229 const Point<dim> &p) const override;
230
241 virtual Tensor<4, dim>
242 shape_4th_derivative_component(const unsigned int i,
243 const Point<dim> &p,
244 const unsigned int component) const override;
245
249 virtual std::size_t
250 memory_consumption() const override;
251
252protected:
253 /*
254 * NOTE: The following function has its definition inlined into the class
255 * declaration because we otherwise run into a compiler error with MS Visual
256 * Studio.
257 */
258
259
260 virtual std::unique_ptr<
263 const UpdateFlags update_flags,
264 const Mapping<dim, spacedim> &mapping,
265 const Quadrature<dim> &quadrature,
267 spacedim>
268 &output_data) const override
269 {
270 (void)mapping;
271
272 // generate a new data object and
273 // initialize some fields
274 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
275 data_ptr = std::make_unique<InternalData>();
276 auto &data = dynamic_cast<InternalData &>(*data_ptr);
277 data.update_each = requires_update_flags(update_flags);
278
279 const unsigned int n_q_points = quadrature.size();
280
281 // initialize some scratch arrays. we need them for the underlying
282 // polynomial to put the values and derivatives of shape functions
283 // to put there, depending on what the user requested
284 std::vector<double> values(
285 update_flags & update_values ? this->n_dofs_per_cell() : 0);
286 std::vector<Tensor<1, dim>> grads(
287 update_flags & update_gradients ? this->n_dofs_per_cell() : 0);
288 std::vector<Tensor<2, dim>> grad_grads(
289 update_flags & update_hessians ? this->n_dofs_per_cell() : 0);
290 std::vector<Tensor<3, dim>> third_derivatives(
291 update_flags & update_3rd_derivatives ? this->n_dofs_per_cell() : 0);
292 std::vector<Tensor<4, dim>>
293 fourth_derivatives; // won't be needed, so leave empty
294
295 // now also initialize fields the fields of this class's own
296 // temporary storage, depending on what we need for the given
297 // update flags.
298 //
299 // there is one exception from the rule: if we are dealing with
300 // cells (i.e., if this function is not called via
301 // get_(sub)face_data()), then we can already store things in the
302 // final location where FEValues::reinit() later wants to see
303 // things. we then don't need the intermediate space. we determine
304 // whether we are on a cell by asking whether the number of
305 // elements in the output array equals the number of quadrature
306 // points (yes, it's a cell) or not (because in that case the
307 // number of quadrature points we use here equals the number of
308 // quadrature points summed over *all* faces or subfaces, whereas
309 // the number of output slots equals the number of quadrature
310 // points on only *one* face)
311 if ((update_flags & update_values) &&
312 !((output_data.shape_values.n_rows() > 0) &&
313 (output_data.shape_values.n_cols() == n_q_points)))
314 data.shape_values.reinit(this->n_dofs_per_cell(), n_q_points);
315
316 if (update_flags & update_gradients)
317 data.shape_gradients.reinit(this->n_dofs_per_cell(), n_q_points);
318
319 if (update_flags & update_hessians)
320 data.shape_hessians.reinit(this->n_dofs_per_cell(), n_q_points);
321
322 if (update_flags & update_3rd_derivatives)
323 data.shape_3rd_derivatives.reinit(this->n_dofs_per_cell(), n_q_points);
324
325 // next already fill those fields of which we have information by
326 // now. note that the shape gradients are only those on the unit
327 // cell, and need to be transformed when visiting an actual cell
328 if (update_flags & (update_values | update_gradients | update_hessians |
330 for (unsigned int i = 0; i < n_q_points; ++i)
331 {
332 poly_space->evaluate(quadrature.point(i),
333 values,
334 grads,
335 grad_grads,
336 third_derivatives,
337 fourth_derivatives);
338
339 // the values of shape functions at quadrature points don't change.
340 // consequently, write these values right into the output array if
341 // we can, i.e., if the output array has the correct size. this is
342 // the case on cells. on faces, we already precompute data on *all*
343 // faces and subfaces, but we later on copy only a portion of it
344 // into the output object; in that case, copy the data from all
345 // faces into the scratch object
346 if (update_flags & update_values)
347 if (output_data.shape_values.n_rows() > 0)
348 {
349 if (output_data.shape_values.n_cols() == n_q_points)
350 for (unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
351 output_data.shape_values[k][i] = values[k];
352 else
353 for (unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
354 data.shape_values[k][i] = values[k];
355 }
356
357 // for everything else, derivatives need to be transformed,
358 // so we write them into our scratch space and only later
359 // copy stuff into where FEValues wants it
360 if (update_flags & update_gradients)
361 for (unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
362 data.shape_gradients[k][i] = grads[k];
363
364 if (update_flags & update_hessians)
365 for (unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
366 data.shape_hessians[k][i] = grad_grads[k];
367
368 if (update_flags & update_3rd_derivatives)
369 for (unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
370 data.shape_3rd_derivatives[k][i] = third_derivatives[k];
371 }
372 return data_ptr;
373 }
374
375 virtual void
378 const CellSimilarity::Similarity cell_similarity,
379 const Quadrature<dim> &quadrature,
380 const Mapping<dim, spacedim> &mapping,
381 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
383 &mapping_data,
384 const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
386 spacedim>
387 &output_data) const override;
388
389 using FiniteElement<dim, spacedim>::fill_fe_face_values;
390
391 virtual void
394 const unsigned int face_no,
395 const hp::QCollection<dim - 1> &quadrature,
396 const Mapping<dim, spacedim> &mapping,
397 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
399 &mapping_data,
400 const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
402 spacedim>
403 &output_data) const override;
404
405 virtual void
408 const unsigned int face_no,
409 const unsigned int sub_no,
410 const Quadrature<dim - 1> &quadrature,
411 const Mapping<dim, spacedim> &mapping,
412 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
414 &mapping_data,
415 const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
417 spacedim>
418 &output_data) const override;
419
473
490 void
493 &output_data,
495 &mapping_data,
496 const unsigned int n_q_points) const;
497
520 void
523 &output_data,
525 &mapping_data,
526 const unsigned int n_q_points) const;
527
528
532 const std::unique_ptr<ScalarPolynomialsBase<dim>> poly_space;
533};
534
538
539#endif
Table< 2, double > shape_values
Definition fe_poly.h:438
Table< 2, Tensor< 3, dim > > shape_3rd_derivatives
Definition fe_poly.h:471
Table< 2, Tensor< 2, dim > > shape_hessians
Definition fe_poly.h:460
Table< 2, Tensor< 1, dim > > shape_gradients
Definition fe_poly.h:449
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const override
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
FE_Poly(const FE_Poly &fe)
virtual Tensor< 3, dim > shape_3rd_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
std::vector< unsigned int > get_poly_space_numbering_inverse() const
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) 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
Definition fe_poly.h:262
FE_Poly(const ScalarPolynomialsBase< dim > &poly_space, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components)
void correct_hessians(internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const unsigned int n_q_points) const
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
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const override
unsigned int get_degree() const
const ScalarPolynomialsBase< dim > & get_poly_space() const
virtual Tensor< 4, dim > shape_4th_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
void correct_third_derivatives(internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const unsigned int n_q_points) const
virtual double shape_value(const unsigned int i, const Point< dim > &p) const override
const std::unique_ptr< ScalarPolynomialsBase< dim > > poly_space
Definition fe_poly.h:532
virtual Tensor< 3, dim > shape_3rd_derivative(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 UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
virtual Tensor< 4, dim > shape_4th_derivative(const unsigned int i, const Point< dim > &p) const override
std::vector< unsigned int > get_poly_space_numbering() const
virtual std::size_t memory_consumption() 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 Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
unsigned int n_dofs_per_cell() const
const std::vector< bool > restriction_is_additive_flags
Definition fe.h:2564
const std::vector< ComponentMask > nonzero_components
Definition fe.h:2573
Abstract base class for mapping classes.
Definition mapping.h:318
Definition point.h:111
const Point< dim > & point(const unsigned int i) const
unsigned int size() const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
UpdateFlags
@ update_hessians
Second derivatives of shape functions.
@ update_values
Shape function values.
@ update_3rd_derivatives
Third derivatives of shape functions.
@ update_gradients
Shape function gradients.