Reference documentation for deal.II version 8.5.1
fe_poly.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2016 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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii__fe_poly_h
17 #define dealii__fe_poly_h
18 
19 
20 #include <deal.II/fe/fe.h>
21 #include <deal.II/base/quadrature.h>
22 
23 DEAL_II_NAMESPACE_OPEN
24 
27 
67 template <class PolynomialType, int dim=PolynomialType::dimension, int spacedim=dim>
68 class FE_Poly : public FiniteElement<dim,spacedim>
69 {
70 public:
74  FE_Poly (const PolynomialType &poly_space,
75  const FiniteElementData<dim> &fe_data,
76  const std::vector<bool> &restriction_is_additive_flags,
77  const std::vector<ComponentMask> &nonzero_components);
78 
83  unsigned int get_degree () const;
84 
85  // for documentation, see the FiniteElement base class
86  virtual
88  requires_update_flags (const UpdateFlags update_flags) const;
89 
95  std::vector<unsigned int> get_poly_space_numbering() const;
96 
101  std::vector<unsigned int> get_poly_space_numbering_inverse() const;
102 
108  virtual double shape_value (const unsigned int i,
109  const Point<dim> &p) const;
110 
121  virtual double shape_value_component (const unsigned int i,
122  const Point<dim> &p,
123  const unsigned int component) const;
124 
130  virtual Tensor<1,dim> shape_grad (const unsigned int i,
131  const Point<dim> &p) const;
132 
143  virtual Tensor<1,dim> shape_grad_component (const unsigned int i,
144  const Point<dim> &p,
145  const unsigned int component) const;
146 
152  virtual Tensor<2,dim> shape_grad_grad (const unsigned int i,
153  const Point<dim> &p) const;
154 
165  virtual Tensor<2,dim> shape_grad_grad_component (const unsigned int i,
166  const Point<dim> &p,
167  const unsigned int component) const;
168 
174  virtual Tensor<3,dim> shape_3rd_derivative (const unsigned int i,
175  const Point<dim> &p) const;
176 
187  virtual Tensor<3,dim> shape_3rd_derivative_component (const unsigned int i,
188  const Point<dim> &p,
189  const unsigned int component) const;
190 
196  virtual Tensor<4,dim> shape_4th_derivative (const unsigned int i,
197  const Point<dim> &p) const;
198 
209  virtual Tensor<4,dim> shape_4th_derivative_component (const unsigned int i,
210  const Point<dim> &p,
211  const unsigned int component) const;
212 
213 protected:
214  /*
215  * NOTE: The following function has its definition inlined into the class declaration
216  * because we otherwise run into a compiler error with MS Visual Studio.
217  */
218 
219 
220  virtual
222  get_data(const UpdateFlags update_flags,
223  const Mapping<dim,spacedim> &/*mapping*/,
224  const Quadrature<dim> &quadrature,
226  {
227  // generate a new data object and
228  // initialize some fields
229  InternalData *data = new InternalData;
230  data->update_each = requires_update_flags(update_flags);
231 
232  const unsigned int n_q_points = quadrature.size();
233 
234  // initialize some scratch arrays. we need them for the underlying
235  // polynomial to put the values and derivatives of shape functions
236  // to put there, depending on what the user requested
237  std::vector<double> values(update_flags & update_values ?
238  this->dofs_per_cell : 0);
239  std::vector<Tensor<1,dim> > grads(update_flags & update_gradients ?
240  this->dofs_per_cell : 0);
241  std::vector<Tensor<2,dim> > grad_grads(update_flags & update_hessians ?
242  this->dofs_per_cell : 0);
243  std::vector<Tensor<3,dim> > third_derivatives(update_flags & update_3rd_derivatives ?
244  this->dofs_per_cell : 0);
245  std::vector<Tensor<4,dim> > fourth_derivatives; // won't be needed, so leave empty
246 
247  // now also initialize fields the fields of this class's own
248  // temporary storage, depending on what we need for the given
249  // update flags.
250  //
251  // there is one exception from the rule: if we are dealing with
252  // cells (i.e., if this function is not called via
253  // get_(sub)face_data()), then we can already store things in the
254  // final location where FEValues::reinit() later wants to see
255  // things. we then don't need the intermediate space. we determine
256  // whether we are on a cell by asking whether the number of
257  // elements in the output array equals the number of quadrature
258  // points (yes, it's a cell) or not (because in that case the
259  // number of quadrature points we use here equals the number of
260  // quadrature points summed over *all* faces or subfaces, whereas
261  // the number of output slots equals the number of quadrature
262  // points on only *one* face)
263  if ((update_flags & update_values)
264  &&
265  !((output_data.shape_values.n_rows() > 0)
266  &&
267  (output_data.shape_values.n_cols() == n_q_points)))
268  data->shape_values.reinit (this->dofs_per_cell, n_q_points);
269 
270  if (update_flags & update_gradients)
271  data->shape_gradients.reinit (this->dofs_per_cell, n_q_points);
272 
273  if (update_flags & update_hessians)
274  data->shape_hessians.reinit (this->dofs_per_cell, n_q_points);
275 
276  if (update_flags & update_3rd_derivatives)
277  data->shape_3rd_derivatives.reinit (this->dofs_per_cell, n_q_points);
278 
279  // next already fill those fields of which we have information by
280  // now. note that the shape gradients are only those on the unit
281  // cell, and need to be transformed when visiting an actual cell
282  if (update_flags & (update_values | update_gradients
284  for (unsigned int i=0; i<n_q_points; ++i)
285  {
286  poly_space.compute(quadrature.point(i),
287  values, grads, grad_grads,
288  third_derivatives,
289  fourth_derivatives);
290 
291  // the values of shape functions at quadrature points don't change.
292  // consequently, write these values right into the output array if
293  // we can, i.e., if the output array has the correct size. this is
294  // the case on cells. on faces, we already precompute data on *all*
295  // faces and subfaces, but we later on copy only a portion of it
296  // into the output object; in that case, copy the data from all
297  // faces into the scratch object
298  if (update_flags & update_values)
299  if (output_data.shape_values.n_rows() > 0)
300  {
301  if (output_data.shape_values.n_cols() == n_q_points)
302  for (unsigned int k=0; k<this->dofs_per_cell; ++k)
303  output_data.shape_values[k][i] = values[k];
304  else
305  for (unsigned int k=0; k<this->dofs_per_cell; ++k)
306  data->shape_values[k][i] = values[k];
307  }
308 
309  // for everything else, derivatives need to be transformed,
310  // so we write them into our scratch space and only later
311  // copy stuff into where FEValues wants it
312  if (update_flags & update_gradients)
313  for (unsigned int k=0; k<this->dofs_per_cell; ++k)
314  data->shape_gradients[k][i] = grads[k];
315 
316  if (update_flags & update_hessians)
317  for (unsigned int k=0; k<this->dofs_per_cell; ++k)
318  data->shape_hessians[k][i] = grad_grads[k];
319 
320  if (update_flags & update_3rd_derivatives)
321  for (unsigned int k=0; k<this->dofs_per_cell; ++k)
322  data->shape_3rd_derivatives[k][i] = third_derivatives[k];
323  }
324  return data;
325  }
326 
327  virtual
328  void
329  fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
330  const CellSimilarity::Similarity cell_similarity,
331  const Quadrature<dim> &quadrature,
332  const Mapping<dim,spacedim> &mapping,
333  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
334  const ::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
335  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
337 
338  virtual
339  void
340  fill_fe_face_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
341  const unsigned int face_no,
342  const Quadrature<dim-1> &quadrature,
343  const Mapping<dim,spacedim> &mapping,
344  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
345  const ::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
346  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
348 
349  virtual
350  void
351  fill_fe_subface_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
352  const unsigned int face_no,
353  const unsigned int sub_no,
354  const Quadrature<dim-1> &quadrature,
355  const Mapping<dim,spacedim> &mapping,
356  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
357  const ::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
358  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
360 
367  class InternalData : public FiniteElement<dim,spacedim>::InternalDataBase
368  {
369  public:
380 
391 
402 
413  };
414 
434  void
437  const unsigned int n_q_points,
438  const unsigned int dof) const;
439 
444  PolynomialType poly_space;
445 };
446 
449 DEAL_II_NAMESPACE_CLOSE
450 
451 #endif
Shape function values.
virtual Tensor< 3, dim > shape_3rd_derivative(const unsigned int i, const Point< dim > &p) const
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
const Point< dim > & point(const unsigned int i) const
Table< 2, Tensor< 2, dim > > shape_hessians
Definition: fe_poly.h:401
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
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)
Table< 2, double > shape_values
Definition: fe_poly.h:379
Third derivatives of shape functions.
UpdateFlags
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
Abstract base class for mapping classes.
Definition: dof_tools.h:46
Table< 2, Tensor< 1, dim > > shape_gradients
Definition: fe_poly.h:390
virtual Tensor< 3, dim > shape_3rd_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
Second derivatives of shape functions.
unsigned int size() const
const unsigned int dofs_per_cell
Definition: fe_base.h:295
unsigned int get_degree() const
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const
virtual double shape_value(const unsigned int i, const Point< dim > &p) const
void correct_third_derivatives(internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data, const internal::FEValues::MappingRelatedData< dim, spacedim > &mapping_data, const unsigned int n_q_points, const unsigned int dof) const
Shape function gradients.
virtual Tensor< 4, dim > shape_4th_derivative(const unsigned int i, const Point< dim > &p) const
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const
const std::vector< ComponentMask > nonzero_components
Definition: fe.h:2355
Table< 2, Tensor< 3, dim > > shape_3rd_derivatives
Definition: fe_poly.h:412
virtual Tensor< 4, dim > shape_4th_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
virtual FiniteElement< dim, spacedim >::InternalDataBase * get_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &, const Quadrature< dim > &quadrature, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const
Definition: fe_poly.h:222
PolynomialType poly_space
Definition: fe_poly.h:444
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:2346
UpdateFlags update_each
Definition: fe.h:638