Reference documentation for deal.II version GIT dad323def1 2022-06-25 19:00: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 
26 
27 #include <deal.II/grid/tria.h>
28 
30 
32 
33 #include <array>
34 #include <cmath>
35 #include <memory>
36 
38 
39 template <typename ElementType, typename MemorySpaceType>
40 class ArrayView;
41 template <int dim>
42 class Quadrature;
43 template <int dim, int spacedim>
44 class FEValues;
45 template <int dim, int spacedim>
46 class FEValuesBase;
47 template <int dim, int spacedim>
48 class FEValues;
49 template <int dim, int spacedim>
50 class FEFaceValues;
51 template <int dim, int spacedim>
52 class FESubfaceValues;
53 namespace NonMatching
54 {
55  template <int dim>
57 }
58 
59 
72 {
77  mapping_none = 0x0000,
78 
83 
88 
94 
100 
107  mapping_piola = 0x0100,
108 
114 
122  mapping_nedelec = 0x0200,
123 
128 
133 
144 
150 
156 };
157 
158 
309 template <int dim, int spacedim = dim>
310 class Mapping : public Subscriptor
311 {
312 public:
316  virtual ~Mapping() override = default;
317 
327  virtual std::unique_ptr<Mapping<dim, spacedim>>
328  clone() const = 0;
329 
344  virtual boost::container::small_vector<Point<spacedim>,
347  const typename Triangulation<dim, spacedim>::cell_iterator &cell) const;
348 
372  virtual Point<spacedim>
374  const bool map_center_of_reference_cell = true) const;
375 
393  virtual BoundingBox<spacedim>
395  const typename Triangulation<dim, spacedim>::cell_iterator &cell) const;
396 
408  virtual bool
410 
415  virtual bool
417 
434  virtual Point<spacedim>
436  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
437  const Point<dim> & p) const = 0;
438 
467  virtual Point<dim>
469  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
470  const Point<spacedim> & p) const = 0;
471 
483  virtual void
485  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
486  const ArrayView<const Point<spacedim>> & real_points,
487  const ArrayView<Point<dim>> &unit_points) const;
488 
499  Point<dim - 1>
501  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
502  const unsigned int face_no,
503  const Point<spacedim> & p) const;
504 
519 
520 
530  "Computing the mapping between a real space point and a point in reference "
531  "space failed, typically because the given point lies outside the cell "
532  "where the inverse mapping is not unique.");
533 
543  double,
544  int,
545  << "The image of the mapping applied to cell with center ["
546  << arg1 << "] is distorted. The cell geometry or the "
547  << "mapping are invalid, giving a non-positive volume "
548  << "fraction of " << arg2 << " in quadrature point " << arg3
549  << '.');
550 
560 public:
627  {
628  public:
634 
639 
643  virtual ~InternalDataBase() = default;
644 
660 
664  virtual std::size_t
666  };
667 
668 
669 protected:
693  virtual UpdateFlags
694  requires_update_flags(const UpdateFlags update_flags) const = 0;
695 
744  virtual std::unique_ptr<InternalDataBase>
745  get_data(const UpdateFlags update_flags,
746  const Quadrature<dim> &quadrature) const = 0;
747 
775  virtual std::unique_ptr<InternalDataBase>
776  get_face_data(const UpdateFlags update_flags,
777  const hp::QCollection<dim - 1> &quadrature) const;
778 
782  virtual std::unique_ptr<InternalDataBase>
783  get_face_data(const UpdateFlags update_flags,
784  const Quadrature<dim - 1> &quadrature) const;
785 
814  virtual std::unique_ptr<InternalDataBase>
815  get_subface_data(const UpdateFlags update_flags,
816  const Quadrature<dim - 1> &quadrature) const = 0;
817 
903  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
904  const CellSimilarity::Similarity cell_similarity,
905  const Quadrature<dim> & quadrature,
906  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
908  &output_data) const = 0;
909 
934  virtual void
936  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
937  const unsigned int face_no,
938  const hp::QCollection<dim - 1> & quadrature,
939  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
941  &output_data) const;
942 
946  virtual void
948  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
949  const unsigned int face_no,
950  const Quadrature<dim - 1> & quadrature,
951  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
953  &output_data) const;
954 
981  virtual void
983  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
984  const unsigned int face_no,
985  const unsigned int subface_no,
986  const Quadrature<dim - 1> & quadrature,
987  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
989  &output_data) const = 0;
990 
997  virtual void
999  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
1001  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
1003  &output_data) const;
1004 
1009 public:
1078  virtual void
1079  transform(const ArrayView<const Tensor<1, dim>> & input,
1080  const MappingKind kind,
1081  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
1082  const ArrayView<Tensor<1, spacedim>> &output) const = 0;
1083 
1131  virtual void
1133  const MappingKind kind,
1134  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
1135  const ArrayView<Tensor<2, spacedim>> &output) const = 0;
1136 
1189  virtual void
1190  transform(const ArrayView<const Tensor<2, dim>> & input,
1191  const MappingKind kind,
1192  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
1193  const ArrayView<Tensor<2, spacedim>> &output) const = 0;
1194 
1236  virtual void
1238  const MappingKind kind,
1239  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
1240  const ArrayView<Tensor<3, spacedim>> &output) const = 0;
1241 
1289  virtual void
1290  transform(const ArrayView<const Tensor<3, dim>> & input,
1291  const MappingKind kind,
1292  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
1293  const ArrayView<Tensor<3, spacedim>> &output) const = 0;
1294 
1300  // Give class @p FEValues access to the private <tt>get_...data</tt> and
1301  // <tt>fill_fe_...values</tt> functions.
1302  friend class FEValuesBase<dim, spacedim>;
1303  friend class FEValues<dim, spacedim>;
1304  friend class FEFaceValues<dim, spacedim>;
1305  friend class FESubfaceValues<dim, spacedim>;
1306  friend class NonMatching::FEImmersedSurfaceValues<dim>;
1307 };
1308 
1309 
1317 template <int dim, int spacedim>
1318 const Mapping<dim, spacedim> &
1320 
1321 
1323 
1324 #endif
UpdateFlags update_each
Definition: mapping.h:659
InternalDataBase(const InternalDataBase &)=delete
virtual ~InternalDataBase()=default
virtual std::size_t memory_consumption() const
Abstract base class for mapping classes.
Definition: mapping.h:311
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 UpdateFlags requires_update_flags(const UpdateFlags update_flags) const =0
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 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 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 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
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 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
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 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:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
UpdateFlags
#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()
MappingKind
Definition: mapping.h:72
@ mapping_piola
Definition: mapping.h:107
@ mapping_covariant_gradient
Definition: mapping.h:93
@ mapping_covariant
Definition: mapping.h:82
@ mapping_nedelec
Definition: mapping.h:122
@ mapping_contravariant
Definition: mapping.h:87
@ mapping_raviart_thomas
Definition: mapping.h:127
@ mapping_contravariant_hessian
Definition: mapping.h:149
@ mapping_covariant_hessian
Definition: mapping.h:143
@ mapping_contravariant_gradient
Definition: mapping.h:99
@ mapping_none
Definition: mapping.h:77
@ mapping_bdm
Definition: mapping.h:132
@ mapping_piola_gradient
Definition: mapping.h:113
@ mapping_piola_hessian
Definition: mapping.h:155
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
Definition: mapping.cc:260
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation