Reference documentation for deal.II version 9.2.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\}}\)
fe_poly.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2020 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 #ifndef dealii_fe_poly_h
17 #define dealii_fe_poly_h
18 
19 
20 #include <deal.II/base/config.h>
21 
24 
25 #include <deal.II/fe/fe.h>
26 
28 
31 
73 template <class PolynomialType,
74  int dim = PolynomialType::dimension,
75  int spacedim = dim>
76 class FE_Poly : public FiniteElement<dim, spacedim>
77 {
78 public:
83  const FiniteElementData<dim> & fe_data,
84  const std::vector<bool> & restriction_is_additive_flags,
85  const std::vector<ComponentMask> &nonzero_components);
86 
91  unsigned int
92  get_degree() const;
93 
94  // for documentation, see the FiniteElement base class
95  virtual UpdateFlags
96  requires_update_flags(const UpdateFlags update_flags) const override;
97 
103  std::vector<unsigned int>
104  get_poly_space_numbering() const;
105 
110  std::vector<unsigned int>
112 
118  virtual double
119  shape_value(const unsigned int i, const Point<dim> &p) const override;
120 
131  virtual double
132  shape_value_component(const unsigned int i,
133  const Point<dim> & p,
134  const unsigned int component) const override;
135 
141  virtual Tensor<1, dim>
142  shape_grad(const unsigned int i, const Point<dim> &p) const override;
143 
154  virtual Tensor<1, dim>
155  shape_grad_component(const unsigned int i,
156  const Point<dim> & p,
157  const unsigned int component) const override;
158 
164  virtual Tensor<2, dim>
165  shape_grad_grad(const unsigned int i, const Point<dim> &p) const override;
166 
177  virtual Tensor<2, dim>
178  shape_grad_grad_component(const unsigned int i,
179  const Point<dim> & p,
180  const unsigned int component) const override;
181 
187  virtual Tensor<3, dim>
188  shape_3rd_derivative(const unsigned int i,
189  const Point<dim> & p) const override;
190 
201  virtual Tensor<3, dim>
202  shape_3rd_derivative_component(const unsigned int i,
203  const Point<dim> & p,
204  const unsigned int component) const override;
205 
211  virtual Tensor<4, dim>
212  shape_4th_derivative(const unsigned int i,
213  const Point<dim> & p) const override;
214 
225  virtual Tensor<4, dim>
226  shape_4th_derivative_component(const unsigned int i,
227  const Point<dim> & p,
228  const unsigned int component) const override;
229 
233  virtual std::size_t
234  memory_consumption() const override;
235 
236 protected:
237  /*
238  * NOTE: The following function has its definition inlined into the class
239  * declaration because we otherwise run into a compiler error with MS Visual
240  * Studio.
241  */
242 
243 
244  virtual std::unique_ptr<
247  const UpdateFlags update_flags,
248  const Mapping<dim, spacedim> & /*mapping*/,
249  const Quadrature<dim> &quadrature,
251  spacedim>
252  &output_data) const override
253  {
254  // generate a new data object and
255  // initialize some fields
256  std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
257  data_ptr = std_cxx14::make_unique<InternalData>();
258  auto &data = dynamic_cast<InternalData &>(*data_ptr);
259  data.update_each = requires_update_flags(update_flags);
260 
261  const unsigned int n_q_points = quadrature.size();
262 
263  // initialize some scratch arrays. we need them for the underlying
264  // polynomial to put the values and derivatives of shape functions
265  // to put there, depending on what the user requested
266  std::vector<double> values(
267  update_flags & update_values ? this->dofs_per_cell : 0);
268  std::vector<Tensor<1, dim>> grads(
269  update_flags & update_gradients ? this->dofs_per_cell : 0);
270  std::vector<Tensor<2, dim>> grad_grads(
271  update_flags & update_hessians ? this->dofs_per_cell : 0);
272  std::vector<Tensor<3, dim>> third_derivatives(
273  update_flags & update_3rd_derivatives ? this->dofs_per_cell : 0);
274  std::vector<Tensor<4, dim>>
275  fourth_derivatives; // won't be needed, so leave empty
276 
277  // now also initialize fields the fields of this class's own
278  // temporary storage, depending on what we need for the given
279  // update flags.
280  //
281  // there is one exception from the rule: if we are dealing with
282  // cells (i.e., if this function is not called via
283  // get_(sub)face_data()), then we can already store things in the
284  // final location where FEValues::reinit() later wants to see
285  // things. we then don't need the intermediate space. we determine
286  // whether we are on a cell by asking whether the number of
287  // elements in the output array equals the number of quadrature
288  // points (yes, it's a cell) or not (because in that case the
289  // number of quadrature points we use here equals the number of
290  // quadrature points summed over *all* faces or subfaces, whereas
291  // the number of output slots equals the number of quadrature
292  // points on only *one* face)
293  if ((update_flags & update_values) &&
294  !((output_data.shape_values.n_rows() > 0) &&
295  (output_data.shape_values.n_cols() == n_q_points)))
296  data.shape_values.reinit(this->dofs_per_cell, n_q_points);
297 
298  if (update_flags & update_gradients)
299  data.shape_gradients.reinit(this->dofs_per_cell, n_q_points);
300 
301  if (update_flags & update_hessians)
302  data.shape_hessians.reinit(this->dofs_per_cell, n_q_points);
303 
304  if (update_flags & update_3rd_derivatives)
305  data.shape_3rd_derivatives.reinit(this->dofs_per_cell, n_q_points);
306 
307  // next already fill those fields of which we have information by
308  // now. note that the shape gradients are only those on the unit
309  // cell, and need to be transformed when visiting an actual cell
310  if (update_flags & (update_values | update_gradients | update_hessians |
312  for (unsigned int i = 0; i < n_q_points; ++i)
313  {
314  poly_space.evaluate(quadrature.point(i),
315  values,
316  grads,
317  grad_grads,
318  third_derivatives,
319  fourth_derivatives);
320 
321  // the values of shape functions at quadrature points don't change.
322  // consequently, write these values right into the output array if
323  // we can, i.e., if the output array has the correct size. this is
324  // the case on cells. on faces, we already precompute data on *all*
325  // faces and subfaces, but we later on copy only a portion of it
326  // into the output object; in that case, copy the data from all
327  // faces into the scratch object
328  if (update_flags & update_values)
329  if (output_data.shape_values.n_rows() > 0)
330  {
331  if (output_data.shape_values.n_cols() == n_q_points)
332  for (unsigned int k = 0; k < this->dofs_per_cell; ++k)
333  output_data.shape_values[k][i] = values[k];
334  else
335  for (unsigned int k = 0; k < this->dofs_per_cell; ++k)
336  data.shape_values[k][i] = values[k];
337  }
338 
339  // for everything else, derivatives need to be transformed,
340  // so we write them into our scratch space and only later
341  // copy stuff into where FEValues wants it
342  if (update_flags & update_gradients)
343  for (unsigned int k = 0; k < this->dofs_per_cell; ++k)
344  data.shape_gradients[k][i] = grads[k];
345 
346  if (update_flags & update_hessians)
347  for (unsigned int k = 0; k < this->dofs_per_cell; ++k)
348  data.shape_hessians[k][i] = grad_grads[k];
349 
350  if (update_flags & update_3rd_derivatives)
351  for (unsigned int k = 0; k < this->dofs_per_cell; ++k)
352  data.shape_3rd_derivatives[k][i] = third_derivatives[k];
353  }
354  return data_ptr;
355  }
356 
357  virtual void
359  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
360  const CellSimilarity::Similarity cell_similarity,
361  const Quadrature<dim> & quadrature,
362  const Mapping<dim, spacedim> & mapping,
363  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
364  const ::internal::FEValuesImplementation::MappingRelatedData<dim,
365  spacedim>
366  & mapping_data,
367  const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
369  spacedim>
370  &output_data) const override;
371 
372  virtual void
374  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
375  const unsigned int face_no,
376  const Quadrature<dim - 1> & quadrature,
377  const Mapping<dim, spacedim> & mapping,
378  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
379  const ::internal::FEValuesImplementation::MappingRelatedData<dim,
380  spacedim>
381  & mapping_data,
382  const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
384  spacedim>
385  &output_data) const override;
386 
387  virtual void
389  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
390  const unsigned int face_no,
391  const unsigned int sub_no,
392  const Quadrature<dim - 1> & quadrature,
393  const Mapping<dim, spacedim> & mapping,
394  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
395  const ::internal::FEValuesImplementation::MappingRelatedData<dim,
396  spacedim>
397  & mapping_data,
398  const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
400  spacedim>
401  &output_data) const override;
402 
409  class InternalData : public FiniteElement<dim, spacedim>::InternalDataBase
410  {
411  public:
422 
433 
444 
455  };
456 
473  void
476  &output_data,
478  & mapping_data,
479  const unsigned int n_q_points) const;
480 
503  void
506  &output_data,
508  & mapping_data,
509  const unsigned int n_q_points) const;
510 
511 
517 };
518 
522 
523 #endif
FE_Poly::correct_hessians
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
FE_Poly::shape_grad
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const override
FE_Poly::shape_value
virtual double shape_value(const unsigned int i, const Point< dim > &p) const override
FE_Poly::poly_space
PolynomialType poly_space
Definition: fe_poly.h:516
FE_Poly::InternalData
Definition: fe_poly.h:409
internal::FEValuesImplementation::FiniteElementRelatedData
Definition: fe_update_flags.h:524
FE_Poly::get_degree
unsigned int get_degree() const
FE_Poly::InternalData::shape_3rd_derivatives
Table< 2, Tensor< 3, dim > > shape_3rd_derivatives
Definition: fe_poly.h:454
internal::FEValuesImplementation::MappingRelatedData
Definition: fe_update_flags.h:418
update_3rd_derivatives
@ update_3rd_derivatives
Third derivatives of shape functions.
Definition: fe_update_flags.h:96
FE_Poly::get_data
virtual std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &, const Quadrature< dim > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
Definition: fe_poly.h:246
FE_Poly::correct_third_derivatives
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
FiniteElementData
Definition: fe_base.h:147
FE_Poly::shape_grad_component
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
FiniteElement< PolynomialType::dimension, PolynomialType::dimension >::InternalDataBase
friend class InternalDataBase
Definition: fe.h:3048
FE_Poly::requires_update_flags
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
FE_Poly::InternalData::shape_values
Table< 2, double > shape_values
Definition: fe_poly.h:421
CellSimilarity::Similarity
Similarity
Definition: fe_update_flags.h:378
update_values
@ update_values
Shape function values.
Definition: fe_update_flags.h:78
FE_Poly::shape_grad_grad_component
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
FE_Poly::shape_4th_derivative_component
virtual Tensor< 4, dim > shape_4th_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
FE_Poly::shape_3rd_derivative
virtual Tensor< 3, dim > shape_3rd_derivative(const unsigned int i, const Point< dim > &p) const override
Table< 2, double >
Quadrature::point
const Point< dim > & point(const unsigned int i) const
Mapping
Abstract base class for mapping classes.
Definition: mapping.h:302
FE_Poly::fill_fe_values
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_Poly::shape_4th_derivative
virtual Tensor< 4, dim > shape_4th_derivative(const unsigned int i, const Point< dim > &p) const override
FE_Poly::shape_value_component
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
FiniteElementData::dofs_per_cell
const unsigned int dofs_per_cell
Definition: fe_base.h:282
Tensor< 1, dim >
FiniteElement< PolynomialType::dimension, PolynomialType::dimension >::nonzero_components
const std::vector< ComponentMask > nonzero_components
Definition: fe.h:2590
FE_Poly::get_poly_space_numbering
std::vector< unsigned int > get_poly_space_numbering() const
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
update_gradients
@ update_gradients
Shape function gradients.
Definition: fe_update_flags.h:84
fe.h
FE_Poly::memory_consumption
virtual std::size_t memory_consumption() const override
FE_Poly::InternalData::shape_hessians
Table< 2, Tensor< 2, dim > > shape_hessians
Definition: fe_poly.h:443
FE_Poly::InternalData::shape_gradients
Table< 2, Tensor< 1, dim > > shape_gradients
Definition: fe_poly.h:432
FE_Poly::shape_3rd_derivative_component
virtual Tensor< 3, dim > shape_3rd_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
FiniteElement
Definition: fe.h:648
UpdateFlags
UpdateFlags
Definition: fe_update_flags.h:66
FE_Poly::FE_Poly
FE_Poly(const PolynomialType &poly_space, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components)
update_hessians
@ update_hessians
Second derivatives of shape functions.
Definition: fe_update_flags.h:90
FE_Poly::fill_fe_subface_values
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
Mapping::InternalDataBase
Definition: mapping.h:597
FE_Poly
Definition: fe_poly.h:76
FiniteElement::InternalDataBase
Definition: fe.h:682
Point< dim >
quadrature.h
FiniteElement< PolynomialType::dimension, PolynomialType::dimension >::restriction_is_additive_flags
const std::vector< bool > restriction_is_additive_flags
Definition: fe.h:2581
config.h
FE_Poly::fill_fe_face_values
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_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
Quadrature::size
unsigned int size() const
memory.h
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
Quadrature
Definition: quadrature.h:85
FE_Poly::shape_grad_grad
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const override
TriaIterator
Definition: tria_iterator.h:578
FE_Poly::get_poly_space_numbering_inverse
std::vector< unsigned int > get_poly_space_numbering_inverse() const
PolynomialType