Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
fe_poly.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2019 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/quadrature.h>
21 #include <deal.II/base/std_cxx14/memory.h>
22 
23 #include <deal.II/fe/fe.h>
24 
25 DEAL_II_NAMESPACE_OPEN
26 
29 
69 template <class PolynomialType,
70  int dim = PolynomialType::dimension,
71  int spacedim = dim>
72 class FE_Poly : public FiniteElement<dim, spacedim>
73 {
74 public:
78  FE_Poly(const PolynomialType & poly_space,
79  const FiniteElementData<dim> & fe_data,
80  const std::vector<bool> & restriction_is_additive_flags,
81  const std::vector<ComponentMask> &nonzero_components);
82 
87  unsigned int
88  get_degree() const;
89 
90  // for documentation, see the FiniteElement base class
91  virtual UpdateFlags
92  requires_update_flags(const UpdateFlags update_flags) const override;
93 
99  std::vector<unsigned int>
100  get_poly_space_numbering() const;
101 
106  std::vector<unsigned int>
108 
114  virtual double
115  shape_value(const unsigned int i, const Point<dim> &p) const override;
116 
127  virtual double
128  shape_value_component(const unsigned int i,
129  const Point<dim> & p,
130  const unsigned int component) const override;
131 
137  virtual Tensor<1, dim>
138  shape_grad(const unsigned int i, const Point<dim> &p) const override;
139 
150  virtual Tensor<1, dim>
151  shape_grad_component(const unsigned int i,
152  const Point<dim> & p,
153  const unsigned int component) const override;
154 
160  virtual Tensor<2, dim>
161  shape_grad_grad(const unsigned int i, const Point<dim> &p) const override;
162 
173  virtual Tensor<2, dim>
174  shape_grad_grad_component(const unsigned int i,
175  const Point<dim> & p,
176  const unsigned int component) const override;
177 
183  virtual Tensor<3, dim>
184  shape_3rd_derivative(const unsigned int i,
185  const Point<dim> & p) const override;
186 
197  virtual Tensor<3, dim>
198  shape_3rd_derivative_component(const unsigned int i,
199  const Point<dim> & p,
200  const unsigned int component) const override;
201 
207  virtual Tensor<4, dim>
208  shape_4th_derivative(const unsigned int i,
209  const Point<dim> & p) const override;
210 
221  virtual Tensor<4, dim>
222  shape_4th_derivative_component(const unsigned int i,
223  const Point<dim> & p,
224  const unsigned int component) const override;
225 
226 protected:
227  /*
228  * NOTE: The following function has its definition inlined into the class
229  * declaration because we otherwise run into a compiler error with MS Visual
230  * Studio.
231  */
232 
233 
234  virtual std::unique_ptr<
237  const UpdateFlags update_flags,
238  const Mapping<dim, spacedim> & /*mapping*/,
239  const Quadrature<dim> &quadrature,
241  spacedim>
242  &output_data) const override
243  {
244  // generate a new data object and
245  // initialize some fields
246  std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
247  data_ptr = std_cxx14::make_unique<InternalData>();
248  auto &data = dynamic_cast<InternalData &>(*data_ptr);
249  data.update_each = requires_update_flags(update_flags);
250 
251  const unsigned int n_q_points = quadrature.size();
252 
253  // initialize some scratch arrays. we need them for the underlying
254  // polynomial to put the values and derivatives of shape functions
255  // to put there, depending on what the user requested
256  std::vector<double> values(
257  update_flags & update_values ? this->dofs_per_cell : 0);
258  std::vector<Tensor<1, dim>> grads(
259  update_flags & update_gradients ? this->dofs_per_cell : 0);
260  std::vector<Tensor<2, dim>> grad_grads(
261  update_flags & update_hessians ? this->dofs_per_cell : 0);
262  std::vector<Tensor<3, dim>> third_derivatives(
263  update_flags & update_3rd_derivatives ? this->dofs_per_cell : 0);
264  std::vector<Tensor<4, dim>>
265  fourth_derivatives; // won't be needed, so leave empty
266 
267  // now also initialize fields the fields of this class's own
268  // temporary storage, depending on what we need for the given
269  // update flags.
270  //
271  // there is one exception from the rule: if we are dealing with
272  // cells (i.e., if this function is not called via
273  // get_(sub)face_data()), then we can already store things in the
274  // final location where FEValues::reinit() later wants to see
275  // things. we then don't need the intermediate space. we determine
276  // whether we are on a cell by asking whether the number of
277  // elements in the output array equals the number of quadrature
278  // points (yes, it's a cell) or not (because in that case the
279  // number of quadrature points we use here equals the number of
280  // quadrature points summed over *all* faces or subfaces, whereas
281  // the number of output slots equals the number of quadrature
282  // points on only *one* face)
283  if ((update_flags & update_values) &&
284  !((output_data.shape_values.n_rows() > 0) &&
285  (output_data.shape_values.n_cols() == n_q_points)))
286  data.shape_values.reinit(this->dofs_per_cell, n_q_points);
287 
288  if (update_flags & update_gradients)
289  data.shape_gradients.reinit(this->dofs_per_cell, n_q_points);
290 
291  if (update_flags & update_hessians)
292  data.shape_hessians.reinit(this->dofs_per_cell, n_q_points);
293 
294  if (update_flags & update_3rd_derivatives)
295  data.shape_3rd_derivatives.reinit(this->dofs_per_cell, n_q_points);
296 
297  // next already fill those fields of which we have information by
298  // now. note that the shape gradients are only those on the unit
299  // cell, and need to be transformed when visiting an actual cell
300  if (update_flags & (update_values | update_gradients | update_hessians |
302  for (unsigned int i = 0; i < n_q_points; ++i)
303  {
304  poly_space.compute(quadrature.point(i),
305  values,
306  grads,
307  grad_grads,
308  third_derivatives,
309  fourth_derivatives);
310 
311  // the values of shape functions at quadrature points don't change.
312  // consequently, write these values right into the output array if
313  // we can, i.e., if the output array has the correct size. this is
314  // the case on cells. on faces, we already precompute data on *all*
315  // faces and subfaces, but we later on copy only a portion of it
316  // into the output object; in that case, copy the data from all
317  // faces into the scratch object
318  if (update_flags & update_values)
319  if (output_data.shape_values.n_rows() > 0)
320  {
321  if (output_data.shape_values.n_cols() == n_q_points)
322  for (unsigned int k = 0; k < this->dofs_per_cell; ++k)
323  output_data.shape_values[k][i] = values[k];
324  else
325  for (unsigned int k = 0; k < this->dofs_per_cell; ++k)
326  data.shape_values[k][i] = values[k];
327  }
328 
329  // for everything else, derivatives need to be transformed,
330  // so we write them into our scratch space and only later
331  // copy stuff into where FEValues wants it
332  if (update_flags & update_gradients)
333  for (unsigned int k = 0; k < this->dofs_per_cell; ++k)
334  data.shape_gradients[k][i] = grads[k];
335 
336  if (update_flags & update_hessians)
337  for (unsigned int k = 0; k < this->dofs_per_cell; ++k)
338  data.shape_hessians[k][i] = grad_grads[k];
339 
340  if (update_flags & update_3rd_derivatives)
341  for (unsigned int k = 0; k < this->dofs_per_cell; ++k)
342  data.shape_3rd_derivatives[k][i] = third_derivatives[k];
343  }
344  return data_ptr;
345  }
346 
347  virtual void
348  fill_fe_values(
349  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
350  const CellSimilarity::Similarity cell_similarity,
351  const Quadrature<dim> & quadrature,
352  const Mapping<dim, spacedim> & mapping,
353  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
354  const ::internal::FEValuesImplementation::MappingRelatedData<dim,
355  spacedim>
356  & mapping_data,
357  const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
359  spacedim>
360  &output_data) const override;
361 
362  virtual void
363  fill_fe_face_values(
364  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
365  const unsigned int face_no,
366  const Quadrature<dim - 1> & quadrature,
367  const Mapping<dim, spacedim> & mapping,
368  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
369  const ::internal::FEValuesImplementation::MappingRelatedData<dim,
370  spacedim>
371  & mapping_data,
372  const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
374  spacedim>
375  &output_data) const override;
376 
377  virtual void
378  fill_fe_subface_values(
379  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
380  const unsigned int face_no,
381  const unsigned int sub_no,
382  const Quadrature<dim - 1> & quadrature,
383  const Mapping<dim, spacedim> & mapping,
384  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
385  const ::internal::FEValuesImplementation::MappingRelatedData<dim,
386  spacedim>
387  & mapping_data,
388  const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
390  spacedim>
391  &output_data) const override;
392 
399  class InternalData : public FiniteElement<dim, spacedim>::InternalDataBase
400  {
401  public:
412 
423 
434 
445  };
446 
469  void
472  &output_data,
474  & mapping_data,
475  const unsigned int n_q_points,
476  const unsigned int dof) const;
477 
482  PolynomialType poly_space;
483 };
484 
487 DEAL_II_NAMESPACE_CLOSE
488 
489 #endif
Shape function values.
virtual double shape_value(const unsigned int i, const Point< dim > &p) const override
Table< 2, Tensor< 2, dim > > shape_hessians
Definition: fe_poly.h:433
virtual Tensor< 4, dim > shape_4th_derivative(const unsigned int i, const Point< dim > &p) const override
const Point< dim > & point(const unsigned int i) const
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const override
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
Table< 2, Tensor< 1, dim > > shape_gradients
Definition: fe_poly.h:422
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)
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:236
Third derivatives of shape functions.
UpdateFlags
Abstract base class for mapping classes.
Definition: dof_tools.h:57
Table< 2, Tensor< 3, dim > > shape_3rd_derivatives
Definition: fe_poly.h:444
Second derivatives of shape functions.
unsigned int size() const
const unsigned int dofs_per_cell
Definition: fe_base.h:282
virtual Tensor< 3, dim > shape_3rd_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
Table< 2, double > shape_values
Definition: fe_poly.h:411
unsigned int get_degree() const
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 unsigned int dof) const
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
virtual Tensor< 4, dim > shape_4th_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
Shape function gradients.
virtual Tensor< 3, dim > shape_3rd_derivative(const unsigned int i, const Point< dim > &p) const override
const std::vector< ComponentMask > nonzero_components
Definition: fe.h:2618
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
PolynomialType poly_space
Definition: fe_poly.h:482
std::vector< unsigned int > get_poly_space_numbering_inverse() const
std::vector< unsigned int > get_poly_space_numbering() const
const std::vector< bool > restriction_is_additive_flags
Definition: fe.h:2609
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override