Reference documentation for deal.II version 8.4.2
fe_poly_tensor.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2005 - 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_tensor_h
17 #define dealii__fe_poly_tensor_h
18 
19 
20 #include <deal.II/lac/full_matrix.h>
21 #include <deal.II/fe/fe.h>
22 #include <deal.II/base/derivative_form.h>
23 #include <deal.II/base/quadrature.h>
24 
25 DEAL_II_NAMESPACE_OPEN
26 
115 template <class PolynomialType, int dim, int spacedim=dim>
116 class FE_PolyTensor : public FiniteElement<dim,spacedim>
117 {
118 public:
125  FE_PolyTensor (const unsigned int degree,
126  const FiniteElementData<dim> &fe_data,
127  const std::vector<bool> &restriction_is_additive_flags,
128  const std::vector<ComponentMask> &nonzero_components);
129 
130  // for documentation, see the FiniteElement base class
131  virtual
133  requires_update_flags (const UpdateFlags update_flags) const;
134 
141  virtual double shape_value (const unsigned int i,
142  const Point<dim> &p) const;
143 
144  // documentation inherited from the base class
145  virtual double shape_value_component (const unsigned int i,
146  const Point<dim> &p,
147  const unsigned int component) const;
148 
155  virtual Tensor<1,dim> shape_grad (const unsigned int i,
156  const Point<dim> &p) const;
157 
158  // documentation inherited from the base class
159  virtual Tensor<1,dim> shape_grad_component (const unsigned int i,
160  const Point<dim> &p,
161  const unsigned int component) const;
162 
169  virtual Tensor<2,dim> shape_grad_grad (const unsigned int i,
170  const Point<dim> &p) const;
171 
172  // documentation inherited from the base class
173  virtual Tensor<2,dim> shape_grad_grad_component (const unsigned int i,
174  const Point<dim> &p,
175  const unsigned int component) const;
176 
177 protected:
183 
184 
185  /* NOTE: The following function has its definition inlined into the class declaration
186  because we otherwise run into a compiler error with MS Visual Studio. */
187  virtual
189  get_data(const UpdateFlags update_flags,
190  const Mapping<dim,spacedim> &/*mapping*/,
191  const Quadrature<dim> &quadrature,
193  {
194  // generate a new data object and
195  // initialize some fields
196  InternalData *data = new InternalData;
197  data->update_each = requires_update_flags(update_flags);
198 
199  const unsigned int n_q_points = quadrature.size();
200 
201  // some scratch arrays
202  std::vector<Tensor<1,dim> > values(0);
203  std::vector<Tensor<2,dim> > grads(0);
204  std::vector<Tensor<3,dim> > grad_grads(0);
205  std::vector<Tensor<4,dim> > third_derivatives(0);
206  std::vector<Tensor<5,dim> > fourth_derivatives(0);
207 
208  if (update_flags & (update_values | update_gradients | update_hessians) )
209  data->sign_change.resize (this->dofs_per_cell);
210 
211  // initialize fields only if really
212  // necessary. otherwise, don't
213  // allocate memory
214  if (update_flags & update_values)
215  {
216  values.resize (this->dofs_per_cell);
217  data->shape_values.reinit (this->dofs_per_cell, n_q_points);
218  if (mapping_type != mapping_none)
219  data->transformed_shape_values.resize (n_q_points);
220  }
221 
222  if (update_flags & update_gradients)
223  {
224  grads.resize (this->dofs_per_cell);
225  data->shape_grads.reinit (this->dofs_per_cell, n_q_points);
226  data->transformed_shape_grads.resize (n_q_points);
227 
228  if ( (mapping_type == mapping_raviart_thomas)
229  ||
230  (mapping_type == mapping_piola)
231  ||
232  (mapping_type == mapping_nedelec)
233  ||
234  (mapping_type == mapping_contravariant))
235  data->untransformed_shape_grads.resize(n_q_points);
236  }
237 
238  if (update_flags & update_hessians)
239  {
240  grad_grads.resize (this->dofs_per_cell);
241  data->shape_grad_grads.reinit (this->dofs_per_cell, n_q_points);
242  data->transformed_shape_hessians.resize (n_q_points);
243  if ( mapping_type != mapping_none )
244  data->untransformed_shape_hessian_tensors.resize(n_q_points);
245  }
246 
247  // Compute shape function values
248  // and derivatives and hessians on
249  // the reference cell.
250  // Make sure, that for the
251  // node values N_i holds
252  // N_i(v_j)=\delta_ij for all basis
253  // functions v_j
254  if (update_flags & (update_values | update_gradients))
255  for (unsigned int k=0; k<n_q_points; ++k)
256  {
257  poly_space.compute(quadrature.point(k),
258  values, grads, grad_grads,
259  third_derivatives,
260  fourth_derivatives);
261 
262  if (update_flags & update_values)
263  {
264  if (inverse_node_matrix.n_cols() == 0)
265  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
266  data->shape_values[i][k] = values[i];
267  else
268  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
269  {
270  Tensor<1,dim> add_values;
271  for (unsigned int j=0; j<this->dofs_per_cell; ++j)
272  add_values += inverse_node_matrix(j,i) * values[j];
273  data->shape_values[i][k] = add_values;
274  }
275  }
276 
277  if (update_flags & update_gradients)
278  {
279  if (inverse_node_matrix.n_cols() == 0)
280  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
281  data->shape_grads[i][k] = grads[i];
282  else
283  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
284  {
285  Tensor<2,dim> add_grads;
286  for (unsigned int j=0; j<this->dofs_per_cell; ++j)
287  add_grads += inverse_node_matrix(j,i) * grads[j];
288  data->shape_grads[i][k] = add_grads;
289  }
290  }
291 
292  if (update_flags & update_hessians)
293  {
294  if (inverse_node_matrix.n_cols() == 0)
295  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
296  data->shape_grad_grads[i][k] = grad_grads[i];
297  else
298  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
299  {
300  Tensor<3,dim> add_grad_grads;
301  for (unsigned int j=0; j<this->dofs_per_cell; ++j)
302  add_grad_grads += inverse_node_matrix(j,i) * grad_grads[j];
303  data->shape_grad_grads[i][k] = add_grad_grads;
304  }
305  }
306 
307  }
308  return data;
309  }
310 
311  virtual
312  void
313  fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
314  const CellSimilarity::Similarity cell_similarity,
315  const Quadrature<dim> &quadrature,
316  const Mapping<dim,spacedim> &mapping,
317  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
318  const ::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
319  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
321 
322  virtual
323  void
324  fill_fe_face_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
325  const unsigned int face_no,
326  const Quadrature<dim-1> &quadrature,
327  const Mapping<dim,spacedim> &mapping,
328  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
329  const ::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
330  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
332 
333  virtual
334  void
335  fill_fe_subface_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
336  const unsigned int face_no,
337  const unsigned int sub_no,
338  const Quadrature<dim-1> &quadrature,
339  const Mapping<dim,spacedim> &mapping,
340  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
341  const ::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
342  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
344 
354  class InternalData : public FiniteElement<dim,spacedim>::InternalDataBase
355  {
356  public:
362 
369 
376 
380  mutable std::vector<double> sign_change;
381  mutable std::vector<Tensor<1, spacedim> > transformed_shape_values;
382  // for shape_gradient computations
383  mutable std::vector<Tensor<2, spacedim > > transformed_shape_grads;
384  mutable std::vector<Tensor<2, dim > > untransformed_shape_grads;
385  // for shape_hessian computations
386  mutable std::vector<Tensor<3, spacedim > > transformed_shape_hessians;
387  mutable std::vector<Tensor<3, dim > > untransformed_shape_hessian_tensors;
388  };
389 
390 
391 
396  PolynomialType poly_space;
397 
410 
417 
421  mutable std::vector<Tensor<1,dim> > cached_values;
422 
426  mutable std::vector<Tensor<2,dim> > cached_grads;
427 
432  mutable std::vector<Tensor<3,dim> > cached_grad_grads;
433 };
434 
435 DEAL_II_NAMESPACE_CLOSE
436 
437 #endif
Shape function values.
std::vector< Tensor< 2, dim > > cached_grads
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:50
MappingType mapping_type
const unsigned int degree
Definition: fe_base.h:299
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 FiniteElement< dim, spacedim >::InternalDataBase * get_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &, const Quadrature< dim > &quadrature, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &) 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:52
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:283
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:2276
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
Definition: table.h:33
FullMatrix< double > inverse_node_matrix
Table< 2, DerivativeForm< 1, dim, spacedim > > shape_grads
const std::vector< bool > restriction_is_additive_flags
Definition: fe.h:2267
UpdateFlags update_each
Definition: fe.h:644