Reference documentation for deal.II version GIT a2efd28e10 2023-02-01 14:45:02+00:00
\(\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.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 2022 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_h
17 #define dealii_mapping_h
18 
19 
20 #include <deal.II/base/config.h>
21 
24 
27 
28 #include <deal.II/grid/tria.h>
29 
31 
33 
34 #include <array>
35 #include <cmath>
36 #include <memory>
37 
39 
40 template <typename ElementType, typename MemorySpaceType>
41 class ArrayView;
42 template <int dim>
43 class Quadrature;
44 template <int dim, int spacedim>
45 class FEValues;
46 template <int dim, int spacedim>
47 class FEValuesBase;
48 template <int dim, int spacedim>
49 class FEValues;
50 template <int dim, int spacedim>
51 class FEFaceValues;
52 template <int dim, int spacedim>
53 class FESubfaceValues;
54 namespace NonMatching
55 {
56  template <int dim>
58 }
59 
60 
73 {
78  mapping_none = 0x0000,
79 
84 
89 
95 
101 
108  mapping_piola = 0x0100,
109 
115 
123  mapping_nedelec = 0x0200,
124 
129 
134 
145 
151 
157 };
158 
159 
310 template <int dim, int spacedim = dim>
311 class Mapping : public Subscriptor
312 {
313 public:
317  virtual ~Mapping() override = default;
318 
328  virtual std::unique_ptr<Mapping<dim, spacedim>>
329  clone() const = 0;
330 
345  virtual boost::container::small_vector<Point<spacedim>,
348  const typename Triangulation<dim, spacedim>::cell_iterator &cell) const;
349 
373  virtual Point<spacedim>
375  const bool map_center_of_reference_cell = true) const;
376 
394  virtual BoundingBox<spacedim>
396  const typename Triangulation<dim, spacedim>::cell_iterator &cell) const;
397 
409  virtual bool
411 
416  virtual bool
418 
435  virtual Point<spacedim>
437  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
438  const Point<dim> & p) const = 0;
439 
468  virtual Point<dim>
470  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
471  const Point<spacedim> & p) const = 0;
472 
484  virtual void
486  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
487  const ArrayView<const Point<spacedim>> & real_points,
488  const ArrayView<Point<dim>> &unit_points) const;
489 
500  Point<dim - 1>
502  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
503  const unsigned int face_no,
504  const Point<spacedim> & p) const;
505 
520 
521 
531  "Computing the mapping between a real space point and a point in reference "
532  "space failed, typically because the given point lies outside the cell "
533  "where the inverse mapping is not unique.");
534 
544  double,
545  int,
546  << "The image of the mapping applied to cell with center ["
547  << arg1 << "] is distorted. The cell geometry or the "
548  << "mapping are invalid, giving a non-positive volume "
549  << "fraction of " << arg2 << " in quadrature point " << arg3
550  << '.');
551 
561 public:
628  {
629  public:
635 
640 
644  virtual ~InternalDataBase() = default;
645 
661 
665  virtual std::size_t
667  };
668 
669 
670 protected:
694  virtual UpdateFlags
695  requires_update_flags(const UpdateFlags update_flags) const = 0;
696 
745  virtual std::unique_ptr<InternalDataBase>
746  get_data(const UpdateFlags update_flags,
747  const Quadrature<dim> &quadrature) const = 0;
748 
776  virtual std::unique_ptr<InternalDataBase>
777  get_face_data(const UpdateFlags update_flags,
778  const hp::QCollection<dim - 1> &quadrature) const;
779 
783  virtual std::unique_ptr<InternalDataBase>
784  get_face_data(const UpdateFlags update_flags,
785  const Quadrature<dim - 1> &quadrature) const;
786 
815  virtual std::unique_ptr<InternalDataBase>
816  get_subface_data(const UpdateFlags update_flags,
817  const Quadrature<dim - 1> &quadrature) const = 0;
818 
904  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
905  const CellSimilarity::Similarity cell_similarity,
906  const Quadrature<dim> & quadrature,
907  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
909  &output_data) const = 0;
910 
935  virtual void
937  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
938  const unsigned int face_no,
939  const hp::QCollection<dim - 1> & quadrature,
940  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
942  &output_data) const;
943 
947  virtual void
949  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
950  const unsigned int face_no,
951  const Quadrature<dim - 1> & quadrature,
952  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
954  &output_data) const;
955 
982  virtual void
984  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
985  const unsigned int face_no,
986  const unsigned int subface_no,
987  const Quadrature<dim - 1> & quadrature,
988  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
990  &output_data) const = 0;
991 
998  virtual void
1000  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
1002  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
1004  &output_data) const;
1005 
1010 public:
1079  virtual void
1080  transform(const ArrayView<const Tensor<1, dim>> & input,
1081  const MappingKind kind,
1082  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
1083  const ArrayView<Tensor<1, spacedim>> &output) const = 0;
1084 
1132  virtual void
1134  const MappingKind kind,
1135  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
1136  const ArrayView<Tensor<2, spacedim>> &output) const = 0;
1137 
1190  virtual void
1191  transform(const ArrayView<const Tensor<2, dim>> & input,
1192  const MappingKind kind,
1193  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
1194  const ArrayView<Tensor<2, spacedim>> &output) const = 0;
1195 
1237  virtual void
1239  const MappingKind kind,
1240  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
1241  const ArrayView<Tensor<3, spacedim>> &output) const = 0;
1242 
1290  virtual void
1291  transform(const ArrayView<const Tensor<3, dim>> & input,
1292  const MappingKind kind,
1293  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
1294  const ArrayView<Tensor<3, spacedim>> &output) const = 0;
1295 
1301  // Give class @p FEValues access to the private <tt>get_...data</tt> and
1302  // <tt>fill_fe_...values</tt> functions.
1303  friend class FEValuesBase<dim, spacedim>;
1304  friend class FEValues<dim, spacedim>;
1305  friend class FEFaceValues<dim, spacedim>;
1306  friend class FESubfaceValues<dim, spacedim>;
1307  friend class NonMatching::FEImmersedSurfaceValues<dim>;
1308 };
1309 
1310 
1318 template <int dim, int spacedim>
1319 const Mapping<dim, spacedim> &
1321 
1322 
1324 
1325 #endif
UpdateFlags update_each
Definition: mapping.h:660
InternalDataBase(const InternalDataBase &)=delete
virtual ~InternalDataBase()=default
virtual std::size_t memory_consumption() const
Abstract base class for mapping classes.
Definition: mapping.h:312
virtual ~Mapping() override=default
virtual bool is_compatible_with(const ReferenceCell &reference_cell) const =0
virtual void transform_points_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const ArrayView< const Point< spacedim >> &real_points, const ArrayView< Point< dim >> &unit_points) const
virtual std::unique_ptr< InternalDataBase > get_data(const UpdateFlags update_flags, const Quadrature< dim > &quadrature) const =0
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const =0
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 =0
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const =0
virtual void fill_fe_immersed_surface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const NonMatching::ImmersedSurfaceQuadrature< dim > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const hp::QCollection< dim - 1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
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 =0
virtual std::unique_ptr< InternalDataBase > get_face_data(const UpdateFlags update_flags, const Quadrature< dim - 1 > &quadrature) const
virtual void transform(const ArrayView< const DerivativeForm< 1, dim, spacedim >> &input, const MappingKind kind, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 2, spacedim >> &output) const =0
virtual std::unique_ptr< InternalDataBase > get_subface_data(const UpdateFlags update_flags, const Quadrature< dim - 1 > &quadrature) const =0
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
virtual void transform(const ArrayView< const Tensor< 2, dim >> &input, const MappingKind kind, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 2, spacedim >> &output) const =0
virtual void transform(const ArrayView< const DerivativeForm< 2, dim, spacedim >> &input, const MappingKind kind, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 3, spacedim >> &output) const =0
virtual boost::container::small_vector< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
Point< dim - 1 > project_real_point_to_unit_point_on_face(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Point< spacedim > &p) const
virtual bool preserves_vertex_locations() const =0
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 =0
virtual BoundingBox< spacedim > get_bounding_box(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const =0
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const =0
virtual std::unique_ptr< InternalDataBase > get_face_data(const UpdateFlags update_flags, const hp::QCollection< dim - 1 > &quadrature) const
virtual Point< spacedim > get_center(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const bool map_center_of_reference_cell=true) const
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:461
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:462
#define DeclException0(Exception0)
Definition: exceptions.h:464
static ::ExceptionBase & ExcTransformationFailed()
static ::ExceptionBase & ExcDistortedMappedCell(Point< spacedim > arg1, double arg2, int arg3)
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:487
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:555
static ::ExceptionBase & ExcInvalidData()
UpdateFlags
MappingKind
Definition: mapping.h:73
@ mapping_piola
Definition: mapping.h:108
@ mapping_covariant_gradient
Definition: mapping.h:94
@ mapping_covariant
Definition: mapping.h:83
@ mapping_nedelec
Definition: mapping.h:123
@ mapping_contravariant
Definition: mapping.h:88
@ mapping_raviart_thomas
Definition: mapping.h:128
@ mapping_contravariant_hessian
Definition: mapping.h:150
@ mapping_covariant_hessian
Definition: mapping.h:144
@ mapping_contravariant_gradient
Definition: mapping.h:100
@ mapping_none
Definition: mapping.h:78
@ mapping_bdm
Definition: mapping.h:133
@ mapping_piola_gradient
Definition: mapping.h:114
@ mapping_piola_hessian
Definition: mapping.h:156
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
Definition: mapping.cc:261
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation