Reference documentation for deal.II version 9.1.0-pre
mapping_fe_field.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2001 - 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.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_mapping_fe_h
17 #define dealii_mapping_fe_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #include <deal.II/base/table.h>
23 #include <deal.II/base/thread_management.h>
24 
25 #include <deal.II/dofs/dof_handler.h>
26 
27 #include <deal.II/fe/fe.h>
28 #include <deal.II/fe/fe_values.h>
29 #include <deal.II/fe/mapping.h>
30 
31 #include <deal.II/lac/vector.h>
32 
33 #include <array>
34 
35 
36 DEAL_II_NAMESPACE_OPEN
37 
38 
41 
79 template <int dim,
80  int spacedim = dim,
81  typename VectorType = Vector<double>,
82  typename DoFHandlerType = DoFHandler<dim, spacedim>>
83 class MappingFEField : public Mapping<dim, spacedim>
84 {
85 public:
118  MappingFEField(const DoFHandlerType &euler_dof_handler,
119  const VectorType & euler_vector,
120  const ComponentMask & mask = ComponentMask());
121 
127 
132  virtual std::unique_ptr<Mapping<dim, spacedim>>
133  clone() const override;
134 
140  virtual bool
141  preserves_vertex_locations() const override;
142 
150  virtual std::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell>
152  const override;
153 
159  // for documentation, see the Mapping base class
160  virtual Point<spacedim>
162  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
163  const Point<dim> &p) const override;
164 
165  // for documentation, see the Mapping base class
166  virtual Point<dim>
168  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
169  const Point<spacedim> &p) const override;
170 
180  // for documentation, see the Mapping base class
181  virtual void
182  transform(const ArrayView<const Tensor<1, dim>> & input,
183  const MappingType type,
184  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
185  const ArrayView<Tensor<1, spacedim>> &output) const override;
186 
187  // for documentation, see the Mapping base class
188  virtual void
190  const MappingType type,
191  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
192  const ArrayView<Tensor<2, spacedim>> &output) const override;
193 
194  // for documentation, see the Mapping base class
195  virtual void
196  transform(const ArrayView<const Tensor<2, dim>> & input,
197  const MappingType type,
198  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
199  const ArrayView<Tensor<2, spacedim>> &output) const override;
200 
201  // for documentation, see the Mapping base class
202  virtual void
204  const MappingType type,
205  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
206  const ArrayView<Tensor<3, spacedim>> &output) const override;
207 
208  // for documentation, see the Mapping base class
209  virtual void
210  transform(const ArrayView<const Tensor<3, dim>> & input,
211  const MappingType type,
212  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
213  const ArrayView<Tensor<3, spacedim>> &output) const override;
214 
223  unsigned int
224  get_degree() const;
225 
231  get_component_mask() const;
232 
237 
238 private:
244  // documentation can be found in Mapping::requires_update_flags()
245  virtual UpdateFlags
246  requires_update_flags(const UpdateFlags update_flags) const override;
247 
248 public:
260  class InternalData : public Mapping<dim, spacedim>::InternalDataBase
261  {
262  public:
267  const ComponentMask & mask);
268 
273  const double &
274  shape(const unsigned int qpoint, const unsigned int shape_nr) const;
275 
279  double &
280  shape(const unsigned int qpoint, const unsigned int shape_nr);
281 
285  const Tensor<1, dim> &
286  derivative(const unsigned int qpoint, const unsigned int shape_nr) const;
287 
292  derivative(const unsigned int qpoint, const unsigned int shape_nr);
293 
297  const Tensor<2, dim> &
298  second_derivative(const unsigned int qpoint,
299  const unsigned int shape_nr) const;
300 
305  second_derivative(const unsigned int qpoint, const unsigned int shape_nr);
306 
310  const Tensor<3, dim> &
311  third_derivative(const unsigned int qpoint,
312  const unsigned int shape_nr) const;
313 
318  third_derivative(const unsigned int qpoint, const unsigned int shape_nr);
319 
323  const Tensor<4, dim> &
324  fourth_derivative(const unsigned int qpoint,
325  const unsigned int shape_nr) const;
326 
331  fourth_derivative(const unsigned int qpoint, const unsigned int shape_nr);
332 
336  virtual std::size_t
337  memory_consumption() const override;
338 
344  std::vector<double> shape_values;
345 
351  std::vector<Tensor<1, dim>> shape_derivatives;
352 
359  std::vector<Tensor<2, dim>> shape_second_derivatives;
360 
367  std::vector<Tensor<3, dim>> shape_third_derivatives;
368 
375  std::vector<Tensor<4, dim>> shape_fourth_derivatives;
376 
390  std::array<std::vector<Tensor<1, dim>>,
393 
400  unsigned int n_shape_functions;
401 
415 
425  mutable std::vector<DerivativeForm<1, dim, spacedim>> covariant;
426 
434  mutable std::vector<DerivativeForm<1, dim, spacedim>> contravariant;
435 
440  mutable std::vector<double> volume_elements;
441 
445  mutable std::vector<std::vector<Tensor<1, spacedim>>> aux;
446 
450  mutable std::vector<types::global_dof_index> local_dof_indices;
451 
455  mutable std::vector<double> local_dof_values;
456  };
457 
458 private:
459  // documentation can be found in Mapping::get_data()
460  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
461  get_data(const UpdateFlags, const Quadrature<dim> &quadrature) const override;
462 
463  // documentation can be found in Mapping::get_face_data()
464  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
465  get_face_data(const UpdateFlags flags,
466  const Quadrature<dim - 1> &quadrature) const override;
467 
468  // documentation can be found in Mapping::get_subface_data()
469  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
470  get_subface_data(const UpdateFlags flags,
471  const Quadrature<dim - 1> &quadrature) const override;
472 
473  // documentation can be found in Mapping::fill_fe_values()
475  fill_fe_values(
476  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
477  const CellSimilarity::Similarity cell_similarity,
478  const Quadrature<dim> & quadrature,
479  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
481  &output_data) const override;
482 
483  // documentation can be found in Mapping::fill_fe_face_values()
484  virtual void
485  fill_fe_face_values(
486  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
487  const unsigned int face_no,
488  const Quadrature<dim - 1> & quadrature,
489  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
491  &output_data) const override;
492 
493  // documentation can be found in Mapping::fill_fe_subface_values()
494  virtual void
495  fill_fe_subface_values(
496  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
497  const unsigned int face_no,
498  const unsigned int subface_no,
499  const Quadrature<dim - 1> & quadrature,
500  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
502  &output_data) const override;
503 
512  SmartPointer<const VectorType,
515 
519  SmartPointer<const DoFHandlerType,
522 
523 private:
540  do_transform_unit_to_real_cell(const InternalData &mdata) const;
541 
542 
557  Point<dim>
559  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
560  const Point<spacedim> & p,
561  const Point<dim> & initial_p_unit,
562  InternalData & mdata) const;
563 
567  void
569  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
571  InternalData &data) const;
572 
576  virtual void
578  const std::vector<Point<dim>> &unit_points,
580  InternalData &data) const;
581 
582  /*
583  * Which components to use for the mapping.
584  */
585  const ComponentMask fe_mask;
586 
596  std::vector<unsigned int> fe_to_real;
597 
603 
608 
609  void
610  compute_data(const UpdateFlags update_flags,
611  const Quadrature<dim> &q,
612  const unsigned int n_original_q_points,
613  InternalData & data) const;
614 
615  void
616  compute_face_data(const UpdateFlags update_flags,
617  const Quadrature<dim> &q,
618  const unsigned int n_original_q_points,
619  InternalData & data) const;
620 
621 
625  template <int, int, class, class>
626  friend class MappingFEField;
627 };
628 
631 /* -------------- declaration of explicit specializations ------------- */
632 
633 #ifndef DOXYGEN
634 
635 #endif // DOXYGEN
636 
637 DEAL_II_NAMESPACE_CLOSE
638 
639 #endif
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
SmartPointer< const VectorType, MappingFEField< dim, spacedim, VectorType, DoFHandlerType > > euler_vector
MappingType
Definition: mapping.h:61
std::vector< unsigned int > fe_to_real
ComponentMask get_component_mask() const
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const override
Point< dim > do_transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p, const Point< dim > &initial_p_unit, InternalData &mdata) const
FEValues< dim, spacedim > fe_values
std::vector< Tensor< 3, dim > > shape_third_derivatives
std::vector< Tensor< 1, dim > > shape_derivatives
std::vector< Tensor< 2, dim > > shape_second_derivatives
UpdateFlags
Threads::Mutex fe_values_mutex
std::vector< Tensor< 4, dim > > shape_fourth_derivatives
Abstract base class for mapping classes.
Definition: dof_tools.h:57
virtual void transform(const ArrayView< const Tensor< 1, dim >> &input, const MappingType type, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 1, spacedim >> &output) const override
const Tensor< 3, dim > & third_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
#define DeclException0(Exception0)
Definition: exceptions.h:385
std::vector< double > volume_elements
void update_internal_dofs(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename MappingFEField< dim, spacedim, VectorType, DoFHandlerType >::InternalData &data) const
const Tensor< 2, dim > & second_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
std::vector< DerivativeForm< 1, dim, spacedim > > covariant
const double & shape(const unsigned int qpoint, const unsigned int shape_nr) const
std::vector< double > local_dof_values
virtual std::array< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
std::vector< DerivativeForm< 1, dim, spacedim > > contravariant
std::vector< double > shape_values
std::vector< std::vector< Tensor< 1, spacedim > > > aux
friend class MappingFEField
std::vector< types::global_dof_index > local_dof_indices
unsigned int get_degree() const
const Tensor< 1, dim > & derivative(const unsigned int qpoint, const unsigned int shape_nr) const
SmartPointer< const DoFHandlerType, MappingFEField< dim, spacedim, VectorType, DoFHandlerType > > euler_dof_handler
const Tensor< 4, dim > & fourth_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
Definition: fe.h:36
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const override
Point< spacedim > do_transform_unit_to_real_cell(const InternalData &mdata) const
virtual bool preserves_vertex_locations() const override
virtual void compute_shapes_virtual(const std::vector< Point< dim >> &unit_points, typename MappingFEField< dim, spacedim, VectorType, DoFHandlerType >::InternalData &data) const
static::ExceptionBase & ExcInactiveCell()
std::array< std::vector< Tensor< 1, dim > >, GeometryInfo< dim >::faces_per_cell *(dim-1)> unit_tangentials
virtual std::size_t memory_consumption() const override
InternalData(const FiniteElement< dim, spacedim > &fe, const ComponentMask &mask)
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override