Reference documentation for deal.II version 9.2.0
\(\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\}}\)
mapping_fe_field.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2001 - 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_mapping_fe_h
17 #define dealii_mapping_fe_h
18 
19 
20 #include <deal.II/base/config.h>
21 
24 
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 
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 
131  MappingFEField(const DoFHandlerType & euler_dof_handler,
132  const std::vector<VectorType> &euler_vector,
133  const ComponentMask & mask = ComponentMask());
134 
142  MappingFEField(const DoFHandlerType & euler_dof_handler,
144  const ComponentMask & mask = ComponentMask());
145 
151 
156  virtual std::unique_ptr<Mapping<dim, spacedim>>
157  clone() const override;
158 
164  virtual bool
165  preserves_vertex_locations() const override;
166 
174  virtual std::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell>
176  const override;
177 
183  // for documentation, see the Mapping base class
184  virtual Point<spacedim>
186  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
187  const Point<dim> &p) const override;
188 
189  // for documentation, see the Mapping base class
190  virtual Point<dim>
192  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
193  const Point<spacedim> &p) const override;
194 
204  // for documentation, see the Mapping base class
205  virtual void
206  transform(const ArrayView<const Tensor<1, dim>> & input,
207  const MappingKind kind,
208  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
209  const ArrayView<Tensor<1, spacedim>> &output) const override;
210 
211  // for documentation, see the Mapping base class
212  virtual void
214  const MappingKind kind,
215  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
216  const ArrayView<Tensor<2, spacedim>> &output) const override;
217 
218  // for documentation, see the Mapping base class
219  virtual void
220  transform(const ArrayView<const Tensor<2, dim>> & input,
221  const MappingKind kind,
222  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
223  const ArrayView<Tensor<2, spacedim>> &output) const override;
224 
225  // for documentation, see the Mapping base class
226  virtual void
228  const MappingKind kind,
229  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
230  const ArrayView<Tensor<3, spacedim>> &output) const override;
231 
232  // for documentation, see the Mapping base class
233  virtual void
234  transform(const ArrayView<const Tensor<3, dim>> & input,
235  const MappingKind kind,
236  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
237  const ArrayView<Tensor<3, spacedim>> &output) const override;
238 
247  unsigned int
248  get_degree() const;
249 
255  get_component_mask() const;
256 
261 
262 private:
268  // documentation can be found in Mapping::requires_update_flags()
269  virtual UpdateFlags
270  requires_update_flags(const UpdateFlags update_flags) const override;
271 
272 public:
284  class InternalData : public Mapping<dim, spacedim>::InternalDataBase
285  {
286  public:
291  const ComponentMask & mask);
292 
297  const double &
298  shape(const unsigned int qpoint, const unsigned int shape_nr) const;
299 
303  double &
304  shape(const unsigned int qpoint, const unsigned int shape_nr);
305 
309  const Tensor<1, dim> &
310  derivative(const unsigned int qpoint, const unsigned int shape_nr) const;
311 
316  derivative(const unsigned int qpoint, const unsigned int shape_nr);
317 
321  const Tensor<2, dim> &
322  second_derivative(const unsigned int qpoint,
323  const unsigned int shape_nr) const;
324 
329  second_derivative(const unsigned int qpoint, const unsigned int shape_nr);
330 
334  const Tensor<3, dim> &
335  third_derivative(const unsigned int qpoint,
336  const unsigned int shape_nr) const;
337 
342  third_derivative(const unsigned int qpoint, const unsigned int shape_nr);
343 
347  const Tensor<4, dim> &
348  fourth_derivative(const unsigned int qpoint,
349  const unsigned int shape_nr) const;
350 
355  fourth_derivative(const unsigned int qpoint, const unsigned int shape_nr);
356 
360  virtual std::size_t
361  memory_consumption() const override;
362 
368  std::vector<double> shape_values;
369 
375  std::vector<Tensor<1, dim>> shape_derivatives;
376 
383  std::vector<Tensor<2, dim>> shape_second_derivatives;
384 
391  std::vector<Tensor<3, dim>> shape_third_derivatives;
392 
399  std::vector<Tensor<4, dim>> shape_fourth_derivatives;
400 
414  std::array<std::vector<Tensor<1, dim>>,
417 
424  unsigned int n_shape_functions;
425 
439 
449  mutable std::vector<DerivativeForm<1, dim, spacedim>> covariant;
450 
458  mutable std::vector<DerivativeForm<1, dim, spacedim>> contravariant;
459 
464  mutable std::vector<double> volume_elements;
465 
469  mutable std::vector<std::vector<Tensor<1, spacedim>>> aux;
470 
474  mutable std::vector<types::global_dof_index> local_dof_indices;
475 
479  mutable std::vector<double> local_dof_values;
480  };
481 
482 private:
483  // documentation can be found in Mapping::get_data()
484  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
485  get_data(const UpdateFlags, const Quadrature<dim> &quadrature) const override;
486 
487  // documentation can be found in Mapping::get_face_data()
488  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
489  get_face_data(const UpdateFlags flags,
490  const Quadrature<dim - 1> &quadrature) const override;
491 
492  // documentation can be found in Mapping::get_subface_data()
493  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
494  get_subface_data(const UpdateFlags flags,
495  const Quadrature<dim - 1> &quadrature) const override;
496 
497  // documentation can be found in Mapping::fill_fe_values()
500  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
501  const CellSimilarity::Similarity cell_similarity,
502  const Quadrature<dim> & quadrature,
503  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
505  &output_data) const override;
506 
507  // documentation can be found in Mapping::fill_fe_face_values()
508  virtual void
510  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
511  const unsigned int face_no,
512  const Quadrature<dim - 1> & quadrature,
513  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
515  &output_data) const override;
516 
517  // documentation can be found in Mapping::fill_fe_subface_values()
518  virtual void
520  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
521  const unsigned int face_no,
522  const unsigned int subface_no,
523  const Quadrature<dim - 1> & quadrature,
524  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
526  &output_data) const override;
527 
536  const bool uses_level_dofs;
537 
541  std::vector<
542  SmartPointer<const VectorType,
545 
549  SmartPointer<const DoFHandlerType,
552 
553 private:
570  do_transform_unit_to_real_cell(const InternalData &mdata) const;
571 
572 
587  Point<dim>
589  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
590  const Point<spacedim> & p,
591  const Point<dim> & initial_p_unit,
592  InternalData & mdata) const;
593 
597  void
599  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
601  InternalData &data) const;
602 
606  virtual void
608  const std::vector<Point<dim>> &unit_points,
610  InternalData &data) const;
611 
612  /*
613  * Which components to use for the mapping.
614  */
616 
626  std::vector<unsigned int> fe_to_real;
627 
633 
637  mutable std::mutex fe_values_mutex;
638 
639  void
640  compute_data(const UpdateFlags update_flags,
641  const Quadrature<dim> &q,
642  const unsigned int n_original_q_points,
643  InternalData & data) const;
644 
645  void
646  compute_face_data(const UpdateFlags update_flags,
647  const Quadrature<dim> &q,
648  const unsigned int n_original_q_points,
649  InternalData & data) const;
650 
651 
652  // Declare other MappingFEField classes friends.
653  template <int, int, class, class>
654  friend class MappingFEField;
655 };
656 
659 /* -------------- declaration of explicit specializations ------------- */
660 
661 #ifndef DOXYGEN
662 
663 #endif // DOXYGEN
664 
666 
667 #endif
MappingFEField::InternalData::shape_fourth_derivatives
std::vector< Tensor< 4, dim > > shape_fourth_derivatives
Definition: mapping_fe_field.h:399
MappingFEField::euler_vector
std::vector< SmartPointer< const VectorType, MappingFEField< dim, spacedim, VectorType, DoFHandlerType > > > euler_vector
Definition: mapping_fe_field.h:544
fe_values.h
MappingFEField::transform
virtual void transform(const ArrayView< const Tensor< 1, dim >> &input, const MappingKind kind, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 1, spacedim >> &output) const override
Definition: mapping_fe_field.cc:1968
MappingKind
MappingKind
Definition: mapping.h:62
MappingFEField
Definition: mapping_fe_field.h:83
MappingFEField::InternalData::InternalData
InternalData(const FiniteElement< dim, spacedim > &fe, const ComponentMask &mask)
Definition: mapping_fe_field.cc:64
MappingFEField::InternalData::unit_tangentials
std::array< std::vector< Tensor< 1, dim > >, GeometryInfo< dim >::faces_per_cell *(dim - 1)> unit_tangentials
Definition: mapping_fe_field.h:416
MappingFEField::fe_values_mutex
std::mutex fe_values_mutex
Definition: mapping_fe_field.h:637
MappingFEField::InternalData::third_derivative
const Tensor< 3, dim > & third_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
Definition: mapping_fe_field.cc:147
internal::FEValuesImplementation::MappingRelatedData
Definition: fe_update_flags.h:418
MappingFEField::fe_mask
const ComponentMask fe_mask
Definition: mapping_fe_field.h:615
MappingFEField::compute_face_data
void compute_face_data(const UpdateFlags update_flags, const Quadrature< dim > &q, const unsigned int n_original_q_points, InternalData &data) const
Definition: mapping_fe_field.cc:567
ArrayView
Definition: array_view.h:77
MappingFEField::requires_update_flags
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
Definition: mapping_fe_field.cc:462
MappingFEField::fe_to_real
std::vector< unsigned int > fe_to_real
Definition: mapping_fe_field.h:626
GeometryInfo
Definition: geometry_info.h:1224
VectorType
MappingFEField::fill_fe_values
virtual CellSimilarity::Similarity fill_fe_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
Definition: mapping_fe_field.cc:1537
MappingFEField::InternalData::fourth_derivative
const Tensor< 4, dim > & fourth_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
Definition: mapping_fe_field.cc:170
MappingFEField::InternalData::shape_derivatives
std::vector< Tensor< 1, dim > > shape_derivatives
Definition: mapping_fe_field.h:375
ComponentMask
Definition: component_mask.h:83
CellSimilarity::Similarity
Similarity
Definition: fe_update_flags.h:378
MappingFEField::fe_values
FEValues< dim, spacedim > fe_values
Definition: mapping_fe_field.h:632
MappingFEField::InternalData::local_dof_indices
std::vector< types::global_dof_index > local_dof_indices
Definition: mapping_fe_field.h:474
MappingFEField::InternalData::volume_elements
std::vector< double > volume_elements
Definition: mapping_fe_field.h:464
DoFHandler
Definition: dof_handler.h:205
FEValues
Definition: fe.h:38
MappingFEField::compute_shapes_virtual
virtual void compute_shapes_virtual(const std::vector< Point< dim >> &unit_points, typename MappingFEField< dim, spacedim, VectorType, DoFHandlerType >::InternalData &data) const
Definition: mapping_fe_field.cc:423
DerivativeForm< 1, dim, spacedim >
thread_management.h
Mapping
Abstract base class for mapping classes.
Definition: mapping.h:302
MappingFEField::InternalData::shape_third_derivatives
std::vector< Tensor< 3, dim > > shape_third_derivatives
Definition: mapping_fe_field.h:391
MappingFEField::InternalData::mask
ComponentMask mask
Definition: mapping_fe_field.h:438
MappingFEField::MappingFEField
friend class MappingFEField
Definition: mapping_fe_field.h:654
MappingFEField::InternalData::memory_consumption
virtual std::size_t memory_consumption() const override
Definition: mapping_fe_field.cc:78
MappingFEField::fill_fe_face_values
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 typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
Definition: mapping_fe_field.cc:1744
MappingFEField::InternalData::shape
const double & shape(const unsigned int qpoint, const unsigned int shape_nr) const
Definition: mapping_fe_field.cc:340
Tensor< 1, dim >
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
MappingFEField::do_transform_unit_to_real_cell
Point< spacedim > do_transform_unit_to_real_cell(const InternalData &mdata) const
Definition: mapping_fe_field.cc:2111
MappingFEField::InternalData::covariant
std::vector< DerivativeForm< 1, dim, spacedim > > covariant
Definition: mapping_fe_field.h:449
fe.h
MappingFEField::get_degree
unsigned int get_degree() const
Definition: mapping_fe_field.cc:2292
MappingFEField::get_face_data
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
Definition: mapping_fe_field.cc:624
mg_level_object.h
MappingFEField::euler_dof_handler
SmartPointer< const DoFHandlerType, MappingFEField< dim, spacedim, VectorType, DoFHandlerType > > euler_dof_handler
Definition: mapping_fe_field.h:551
FiniteElement
Definition: fe.h:648
MappingFEField::get_component_mask
ComponentMask get_component_mask() const
Definition: mapping_fe_field.cc:2300
MappingFEField::InternalData::derivative
const Tensor< 1, dim > & derivative(const unsigned int qpoint, const unsigned int shape_nr) const
Definition: mapping_fe_field.cc:100
UpdateFlags
UpdateFlags
Definition: fe_update_flags.h:66
MappingFEField::InternalData::n_shape_functions
unsigned int n_shape_functions
Definition: mapping_fe_field.h:424
MappingFEField::InternalData::shape_second_derivatives
std::vector< Tensor< 2, dim > > shape_second_derivatives
Definition: mapping_fe_field.h:383
SmartPointer
Definition: smartpointer.h:68
MappingFEField::InternalData::local_dof_values
std::vector< double > local_dof_values
Definition: mapping_fe_field.h:479
mapping.h
Mapping::InternalDataBase
Definition: mapping.h:597
MappingFEField::fill_fe_subface_values
virtual void fill_fe_subface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no, const Quadrature< dim - 1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
Definition: mapping_fe_field.cc:1783
MGLevelObject< VectorType >
MappingFEField::compute_data
void compute_data(const UpdateFlags update_flags, const Quadrature< dim > &q, const unsigned int n_original_q_points, InternalData &data) const
Definition: mapping_fe_field.cc:516
vector.h
dof_handler.h
DeclException0
#define DeclException0(Exception0)
Definition: exceptions.h:473
MappingFEField::InternalData::second_derivative
const Tensor< 2, dim > & second_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
Definition: mapping_fe_field.cc:123
MappingFEField::InternalData::aux
std::vector< std::vector< Tensor< 1, spacedim > > > aux
Definition: mapping_fe_field.h:469
Point< spacedim >
MappingFEField::get_data
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
Definition: mapping_fe_field.cc:608
MappingFEField::InternalData::contravariant
std::vector< DerivativeForm< 1, dim, spacedim > > contravariant
Definition: mapping_fe_field.h:458
MappingFEField::transform_real_to_unit_cell
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
Definition: mapping_fe_field.cc:2132
config.h
MappingFEField::uses_level_dofs
const bool uses_level_dofs
Definition: mapping_fe_field.h:536
MappingFEField::transform_unit_to_real_cell
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
Definition: mapping_fe_field.cc:2089
MappingFEField::preserves_vertex_locations
virtual bool preserves_vertex_locations() const override
Definition: mapping_fe_field.cc:353
MappingFEField::ExcInactiveCell
static ::ExceptionBase & ExcInactiveCell()
MappingFEField::get_vertices
virtual std::array< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
Definition: mapping_fe_field.cc:362
MappingFEField::get_subface_data
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
Definition: mapping_fe_field.cc:640
MappingFEField::update_internal_dofs
void update_internal_dofs(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename MappingFEField< dim, spacedim, VectorType, DoFHandlerType >::InternalData &data) const
Definition: mapping_fe_field.cc:2318
MappingFEField::InternalData::shape_values
std::vector< double > shape_values
Definition: mapping_fe_field.h:368
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
MappingFEField::clone
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
Definition: mapping_fe_field.cc:2309
Quadrature
Definition: quadrature.h:85
Vector< double >
MappingFEField::do_transform_real_to_unit_cell
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
Definition: mapping_fe_field.cc:2181
MappingFEField::InternalData
Definition: mapping_fe_field.h:284
TriaIterator
Definition: tria_iterator.h:578