Reference documentation for deal.II version Git 662679e6de 2020-10-27 13:02:35 -0600
\(\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 
27 #include <memory>
28 
30 
33 
73 template <int dim, int spacedim = dim>
74 class FE_Poly : public FiniteElement<dim, spacedim>
75 {
76 public:
81  const FiniteElementData<dim> & fe_data,
82  const std::vector<bool> & restriction_is_additive_flags,
83  const std::vector<ComponentMask> &nonzero_components);
84 
88  FE_Poly(const FE_Poly &fe);
89 
94  unsigned int
95  get_degree() const;
96 
97  // for documentation, see the FiniteElement base class
98  virtual UpdateFlags
99  requires_update_flags(const UpdateFlags update_flags) const override;
100 
105  get_poly_space() const;
106 
112  std::vector<unsigned int>
113  get_poly_space_numbering() const;
114 
119  std::vector<unsigned int>
121 
127  virtual double
128  shape_value(const unsigned int i, const Point<dim> &p) const override;
129 
140  virtual double
141  shape_value_component(const unsigned int i,
142  const Point<dim> & p,
143  const unsigned int component) const override;
144 
150  virtual Tensor<1, dim>
151  shape_grad(const unsigned int i, const Point<dim> &p) const override;
152 
163  virtual Tensor<1, dim>
164  shape_grad_component(const unsigned int i,
165  const Point<dim> & p,
166  const unsigned int component) const override;
167 
173  virtual Tensor<2, dim>
174  shape_grad_grad(const unsigned int i, const Point<dim> &p) const override;
175 
186  virtual Tensor<2, dim>
187  shape_grad_grad_component(const unsigned int i,
188  const Point<dim> & p,
189  const unsigned int component) const override;
190 
196  virtual Tensor<3, dim>
197  shape_3rd_derivative(const unsigned int i,
198  const Point<dim> & p) const override;
199 
210  virtual Tensor<3, dim>
211  shape_3rd_derivative_component(const unsigned int i,
212  const Point<dim> & p,
213  const unsigned int component) const override;
214 
220  virtual Tensor<4, dim>
221  shape_4th_derivative(const unsigned int i,
222  const Point<dim> & p) const override;
223 
234  virtual Tensor<4, dim>
235  shape_4th_derivative_component(const unsigned int i,
236  const Point<dim> & p,
237  const unsigned int component) const override;
238 
242  virtual std::size_t
243  memory_consumption() const override;
244 
245 protected:
246  /*
247  * NOTE: The following function has its definition inlined into the class
248  * declaration because we otherwise run into a compiler error with MS Visual
249  * Studio.
250  */
251 
252 
253  virtual std::unique_ptr<
256  const UpdateFlags update_flags,
257  const Mapping<dim, spacedim> &mapping,
258  const Quadrature<dim> & quadrature,
260  spacedim>
261  &output_data) const override
262  {
263  (void)mapping;
264 
265  // generate a new data object and
266  // initialize some fields
267  std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
268  data_ptr = std::make_unique<InternalData>();
269  auto &data = dynamic_cast<InternalData &>(*data_ptr);
270  data.update_each = requires_update_flags(update_flags);
271 
272  const unsigned int n_q_points = quadrature.size();
273 
274  // initialize some scratch arrays. we need them for the underlying
275  // polynomial to put the values and derivatives of shape functions
276  // to put there, depending on what the user requested
277  std::vector<double> values(
278  update_flags & update_values ? this->n_dofs_per_cell() : 0);
279  std::vector<Tensor<1, dim>> grads(
280  update_flags & update_gradients ? this->n_dofs_per_cell() : 0);
281  std::vector<Tensor<2, dim>> grad_grads(
282  update_flags & update_hessians ? this->n_dofs_per_cell() : 0);
283  std::vector<Tensor<3, dim>> third_derivatives(
284  update_flags & update_3rd_derivatives ? this->n_dofs_per_cell() : 0);
285  std::vector<Tensor<4, dim>>
286  fourth_derivatives; // won't be needed, so leave empty
287 
288  // now also initialize fields the fields of this class's own
289  // temporary storage, depending on what we need for the given
290  // update flags.
291  //
292  // there is one exception from the rule: if we are dealing with
293  // cells (i.e., if this function is not called via
294  // get_(sub)face_data()), then we can already store things in the
295  // final location where FEValues::reinit() later wants to see
296  // things. we then don't need the intermediate space. we determine
297  // whether we are on a cell by asking whether the number of
298  // elements in the output array equals the number of quadrature
299  // points (yes, it's a cell) or not (because in that case the
300  // number of quadrature points we use here equals the number of
301  // quadrature points summed over *all* faces or subfaces, whereas
302  // the number of output slots equals the number of quadrature
303  // points on only *one* face)
304  if ((update_flags & update_values) &&
305  !((output_data.shape_values.n_rows() > 0) &&
306  (output_data.shape_values.n_cols() == n_q_points)))
307  data.shape_values.reinit(this->n_dofs_per_cell(), n_q_points);
308 
309  if (update_flags & update_gradients)
310  data.shape_gradients.reinit(this->n_dofs_per_cell(), n_q_points);
311 
312  if (update_flags & update_hessians)
313  data.shape_hessians.reinit(this->n_dofs_per_cell(), n_q_points);
314 
315  if (update_flags & update_3rd_derivatives)
316  data.shape_3rd_derivatives.reinit(this->n_dofs_per_cell(), n_q_points);
317 
318  // next already fill those fields of which we have information by
319  // now. note that the shape gradients are only those on the unit
320  // cell, and need to be transformed when visiting an actual cell
321  if (update_flags & (update_values | update_gradients | update_hessians |
322  update_3rd_derivatives))
323  for (unsigned int i = 0; i < n_q_points; ++i)
324  {
325  poly_space->evaluate(quadrature.point(i),
326  values,
327  grads,
328  grad_grads,
329  third_derivatives,
330  fourth_derivatives);
331 
332  // the values of shape functions at quadrature points don't change.
333  // consequently, write these values right into the output array if
334  // we can, i.e., if the output array has the correct size. this is
335  // the case on cells. on faces, we already precompute data on *all*
336  // faces and subfaces, but we later on copy only a portion of it
337  // into the output object; in that case, copy the data from all
338  // faces into the scratch object
339  if (update_flags & update_values)
340  if (output_data.shape_values.n_rows() > 0)
341  {
342  if (output_data.shape_values.n_cols() == n_q_points)
343  for (unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
344  output_data.shape_values[k][i] = values[k];
345  else
346  for (unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
347  data.shape_values[k][i] = values[k];
348  }
349 
350  // for everything else, derivatives need to be transformed,
351  // so we write them into our scratch space and only later
352  // copy stuff into where FEValues wants it
353  if (update_flags & update_gradients)
354  for (unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
355  data.shape_gradients[k][i] = grads[k];
356 
357  if (update_flags & update_hessians)
358  for (unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
359  data.shape_hessians[k][i] = grad_grads[k];
360 
361  if (update_flags & update_3rd_derivatives)
362  for (unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
363  data.shape_3rd_derivatives[k][i] = third_derivatives[k];
364  }
365  return data_ptr;
366  }
367 
368  virtual void
370  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
371  const CellSimilarity::Similarity cell_similarity,
372  const Quadrature<dim> & quadrature,
373  const Mapping<dim, spacedim> & mapping,
374  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
375  const ::internal::FEValuesImplementation::MappingRelatedData<dim,
376  spacedim>
377  & mapping_data,
378  const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
380  spacedim>
381  &output_data) const override;
382 
383  virtual void
385  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
386  const unsigned int face_no,
387  const Quadrature<dim - 1> & quadrature,
388  const Mapping<dim, spacedim> & mapping,
389  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
390  const ::internal::FEValuesImplementation::MappingRelatedData<dim,
391  spacedim>
392  & mapping_data,
393  const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
395  spacedim>
396  &output_data) const override;
397 
398  virtual void
400  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
401  const unsigned int face_no,
402  const unsigned int sub_no,
403  const Quadrature<dim - 1> & quadrature,
404  const Mapping<dim, spacedim> & mapping,
405  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
406  const ::internal::FEValuesImplementation::MappingRelatedData<dim,
407  spacedim>
408  & mapping_data,
409  const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
411  spacedim>
412  &output_data) const override;
413 
420  class InternalData : public FiniteElement<dim, spacedim>::InternalDataBase
421  {
422  public:
433 
444 
455 
466  };
467 
484  void
487  &output_data,
489  & mapping_data,
490  const unsigned int n_q_points) const;
491 
514  void
517  &output_data,
519  & mapping_data,
520  const unsigned int n_q_points) const;
521 
522 
526  const std::unique_ptr<ScalarPolynomialsBase<dim>> poly_space;
527 };
528 
532 
533 #endif
Shape function values.
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)
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
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const override
virtual std::size_t memory_consumption() const override
const Point< dim > & point(const unsigned int i) const
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const override
const ScalarPolynomialsBase< dim > & get_poly_space() const
Table< 2, double > shape_values
Definition: fe_poly.h:432
virtual double shape_value(const unsigned int i, const Point< dim > &p) 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
Third derivatives of shape functions.
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) 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
UpdateFlags
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
Abstract base class for mapping classes.
Definition: mapping.h:301
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:369
virtual Tensor< 3, dim > shape_3rd_derivative(const unsigned int i, const Point< dim > &p) const override
Second derivatives of shape functions.
unsigned int size() const
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
virtual Tensor< 4, dim > shape_4th_derivative(const unsigned int i, const Point< dim > &p) const override
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 UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
unsigned int get_degree() const
unsigned int n_dofs_per_cell() const
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:255
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:368
Shape function gradients.
virtual Tensor< 4, dim > shape_4th_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
const std::vector< ComponentMask > nonzero_components
Definition: fe.h:2578
Table< 2, Tensor< 2, dim > > shape_hessians
Definition: fe_poly.h:454
Table< 2, Tensor< 1, dim > > shape_gradients
Definition: fe_poly.h:443
Table< 2, Tensor< 3, dim > > shape_3rd_derivatives
Definition: fe_poly.h:465
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
std::vector< unsigned int > get_poly_space_numbering() const
const std::unique_ptr< ScalarPolynomialsBase< dim > > poly_space
Definition: fe_poly.h:526
const std::vector< bool > restriction_is_additive_flags
Definition: fe.h:2569