Reference documentation for deal.II version 9.0.0
fe_poly_tensor.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2005 - 2018 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_tensor_h
17 #define dealii_fe_poly_tensor_h
18 
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/lac/full_matrix.h>
22 #include <deal.II/fe/fe.h>
23 #include <deal.II/base/derivative_form.h>
24 #include <deal.II/base/quadrature.h>
25 #include <deal.II/base/std_cxx14/memory.h>
26 #include <deal.II/base/thread_management.h>
27 
28 
29 DEAL_II_NAMESPACE_OPEN
30 
139 template <class PolynomialType, int dim, int spacedim=dim>
140 class FE_PolyTensor : public FiniteElement<dim,spacedim>
141 {
142 public:
149  FE_PolyTensor (const unsigned int degree,
150  const FiniteElementData<dim> &fe_data,
151  const std::vector<bool> &restriction_is_additive_flags,
152  const std::vector<ComponentMask> &nonzero_components);
153 
154  // for documentation, see the FiniteElement base class
155  virtual
157  requires_update_flags (const UpdateFlags update_flags) const;
158 
165  virtual double shape_value (const unsigned int i,
166  const Point<dim> &p) const;
167 
168  // documentation inherited from the base class
169  virtual double shape_value_component (const unsigned int i,
170  const Point<dim> &p,
171  const unsigned int component) const;
172 
179  virtual Tensor<1,dim> shape_grad (const unsigned int i,
180  const Point<dim> &p) const;
181 
182  // documentation inherited from the base class
183  virtual Tensor<1,dim> shape_grad_component (const unsigned int i,
184  const Point<dim> &p,
185  const unsigned int component) const;
186 
193  virtual Tensor<2,dim> shape_grad_grad (const unsigned int i,
194  const Point<dim> &p) const;
195 
196  // documentation inherited from the base class
197  virtual Tensor<2,dim> shape_grad_grad_component (const unsigned int i,
198  const Point<dim> &p,
199  const unsigned int component) const;
200 
201 protected:
207 
208 
209  /* NOTE: The following function has its definition inlined into the class declaration
210  because we otherwise run into a compiler error with MS Visual Studio. */
211  virtual
212  std::unique_ptr<typename FiniteElement<dim,spacedim>::InternalDataBase>
213  get_data(const UpdateFlags update_flags,
214  const Mapping<dim,spacedim> &/*mapping*/,
215  const Quadrature<dim> &quadrature,
217  {
218  // generate a new data object and
219  // initialize some fields
220  auto data = std_cxx14::make_unique<InternalData>();
221  data->update_each = requires_update_flags(update_flags);
222 
223  const unsigned int n_q_points = quadrature.size();
224 
225  // some scratch arrays
226  std::vector<Tensor<1,dim> > values(0);
227  std::vector<Tensor<2,dim> > grads(0);
228  std::vector<Tensor<3,dim> > grad_grads(0);
229  std::vector<Tensor<4,dim> > third_derivatives(0);
230  std::vector<Tensor<5,dim> > fourth_derivatives(0);
231 
232  if (update_flags & (update_values | update_gradients | update_hessians) )
233  data->sign_change.resize (this->dofs_per_cell);
234 
235  // initialize fields only if really
236  // necessary. otherwise, don't
237  // allocate memory
238  if (update_flags & update_values)
239  {
240  values.resize (this->dofs_per_cell);
241  data->shape_values.reinit (this->dofs_per_cell, n_q_points);
242  if (mapping_type != mapping_none)
243  data->transformed_shape_values.resize (n_q_points);
244  }
245 
246  if (update_flags & update_gradients)
247  {
248  grads.resize (this->dofs_per_cell);
249  data->shape_grads.reinit (this->dofs_per_cell, n_q_points);
250  data->transformed_shape_grads.resize (n_q_points);
251 
253  ||
255  ||
257  ||
259  data->untransformed_shape_grads.resize(n_q_points);
260  }
261 
262  if (update_flags & update_hessians)
263  {
264  grad_grads.resize (this->dofs_per_cell);
265  data->shape_grad_grads.reinit (this->dofs_per_cell, n_q_points);
266  data->transformed_shape_hessians.resize (n_q_points);
267  if ( mapping_type != mapping_none )
268  data->untransformed_shape_hessian_tensors.resize(n_q_points);
269  }
270 
271  // Compute shape function values
272  // and derivatives and hessians on
273  // the reference cell.
274  // Make sure, that for the
275  // node values N_i holds
276  // N_i(v_j)=\delta_ij for all basis
277  // functions v_j
278  if (update_flags & (update_values | update_gradients))
279  for (unsigned int k=0; k<n_q_points; ++k)
280  {
281  poly_space.compute(quadrature.point(k),
282  values, grads, grad_grads,
283  third_derivatives,
284  fourth_derivatives);
285 
286  if (update_flags & update_values)
287  {
288  if (inverse_node_matrix.n_cols() == 0)
289  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
290  data->shape_values[i][k] = values[i];
291  else
292  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
293  {
294  Tensor<1,dim> add_values;
295  for (unsigned int j=0; j<this->dofs_per_cell; ++j)
296  add_values += inverse_node_matrix(j,i) * values[j];
297  data->shape_values[i][k] = add_values;
298  }
299  }
300 
301  if (update_flags & update_gradients)
302  {
303  if (inverse_node_matrix.n_cols() == 0)
304  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
305  data->shape_grads[i][k] = grads[i];
306  else
307  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
308  {
309  Tensor<2,dim> add_grads;
310  for (unsigned int j=0; j<this->dofs_per_cell; ++j)
311  add_grads += inverse_node_matrix(j,i) * grads[j];
312  data->shape_grads[i][k] = add_grads;
313  }
314  }
315 
316  if (update_flags & update_hessians)
317  {
318  if (inverse_node_matrix.n_cols() == 0)
319  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
320  data->shape_grad_grads[i][k] = grad_grads[i];
321  else
322  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
323  {
324  Tensor<3,dim> add_grad_grads;
325  for (unsigned int j=0; j<this->dofs_per_cell; ++j)
326  add_grad_grads += inverse_node_matrix(j,i) * grad_grads[j];
327  data->shape_grad_grads[i][k] = add_grad_grads;
328  }
329  }
330 
331  }
332  return std::move(data);
333  }
334 
335  virtual
336  void
337  fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
338  const CellSimilarity::Similarity cell_similarity,
339  const Quadrature<dim> &quadrature,
340  const Mapping<dim,spacedim> &mapping,
341  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
342  const ::internal::FEValuesImplementation::MappingRelatedData<dim, spacedim> &mapping_data,
343  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
345 
346  virtual
347  void
348  fill_fe_face_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
349  const unsigned int face_no,
350  const Quadrature<dim-1> &quadrature,
351  const Mapping<dim,spacedim> &mapping,
352  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
353  const ::internal::FEValuesImplementation::MappingRelatedData<dim, spacedim> &mapping_data,
354  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
356 
357  virtual
358  void
359  fill_fe_subface_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
360  const unsigned int face_no,
361  const unsigned int sub_no,
362  const Quadrature<dim-1> &quadrature,
363  const Mapping<dim,spacedim> &mapping,
364  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
365  const ::internal::FEValuesImplementation::MappingRelatedData<dim, spacedim> &mapping_data,
366  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
368 
378  class InternalData : public FiniteElement<dim,spacedim>::InternalDataBase
379  {
380  public:
386 
393 
400 
404  mutable std::vector<double> sign_change;
405  mutable std::vector<Tensor<1, spacedim> > transformed_shape_values;
406  // for shape_gradient computations
407  mutable std::vector<Tensor<2, spacedim > > transformed_shape_grads;
408  mutable std::vector<Tensor<2, dim > > untransformed_shape_grads;
409  // for shape_hessian computations
410  mutable std::vector<Tensor<3, spacedim > > transformed_shape_hessians;
411  mutable std::vector<Tensor<3, dim > > untransformed_shape_hessian_tensors;
412  };
413 
414 
415 
420  PolynomialType poly_space;
421 
434 
439 
446 
450  mutable std::vector<Tensor<1,dim> > cached_values;
451 
455  mutable std::vector<Tensor<2,dim> > cached_grads;
456 
461  mutable std::vector<Tensor<3,dim> > cached_grad_grads;
462 };
463 
464 DEAL_II_NAMESPACE_CLOSE
465 
466 #endif
Shape function values.
std::vector< Tensor< 2, dim > > cached_grads
Threads::Mutex cache_mutex
Point< dim > cached_point
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const
MappingType
Definition: mapping.h:51
MappingType mapping_type
const unsigned int degree
Definition: fe_base.h:311
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 > &) const
FE_PolyTensor(const unsigned int degree, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components)
const Point< dim > & point(const unsigned int i) const
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
virtual double shape_value(const unsigned int i, const Point< dim > &p) const
PolynomialType poly_space
Table< 2, Tensor< 1, dim > > shape_values
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const
UpdateFlags
Abstract base class for mapping classes.
Definition: dof_tools.h:46
Table< 2, DerivativeForm< 2, dim, spacedim > > shape_grad_grads
std::vector< double > sign_change
Second derivatives of shape functions.
unsigned int size() const
const unsigned int dofs_per_cell
Definition: fe_base.h:295
std::vector< Tensor< 1, dim > > cached_values
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const
std::vector< Tensor< 3, dim > > cached_grad_grads
Shape function gradients.
const std::vector< ComponentMask > nonzero_components
Definition: fe.h:2550
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
Definition: table.h:34
FullMatrix< double > inverse_node_matrix
Table< 2, DerivativeForm< 1, dim, spacedim > > shape_grads
const std::vector< bool > restriction_is_additive_flags
Definition: fe.h:2541