Reference documentation for deal.II version 9.3.3
\(\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 - 2021 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
31#include <array>
32#include <cmath>
33#include <memory>
34
36
37template <typename ElementType, typename MemorySpaceType>
38class ArrayView;
39template <int dim>
40class Quadrature;
41template <int dim, int spacedim>
42class FEValues;
43template <int dim, int spacedim>
44class FEValuesBase;
45template <int dim, int spacedim>
46class FEValues;
47template <int dim, int spacedim>
48class FEFaceValues;
49template <int dim, int spacedim>
50class FESubfaceValues;
51
52
65{
70 mapping_none = 0x0000,
71
76
81
87
93
101
107
116
121
126
137
143
150
151
302template <int dim, int spacedim = dim>
303class Mapping : public Subscriptor
304{
305public:
309 virtual ~Mapping() override = default;
310
320 virtual std::unique_ptr<Mapping<dim, spacedim>>
321 clone() const = 0;
322
337 virtual boost::container::small_vector<Point<spacedim>,
340 const typename Triangulation<dim, spacedim>::cell_iterator &cell) const;
341
365 virtual Point<spacedim>
367 const bool map_center_of_reference_cell = true) const;
368
388 const typename Triangulation<dim, spacedim>::cell_iterator &cell) const;
389
401 virtual bool
403
408 virtual bool
410
427 virtual Point<spacedim>
430 const Point<dim> & p) const = 0;
431
460 virtual Point<dim>
463 const Point<spacedim> & p) const = 0;
464
476 virtual void
479 const ArrayView<const Point<spacedim>> & real_points,
480 const ArrayView<Point<dim>> &unit_points) const;
481
492 Point<dim - 1>
495 const unsigned int face_no,
496 const Point<spacedim> & p) const;
497
512
513
523 "Computing the mapping between a real space point and a point in reference "
524 "space failed, typically because the given point lies outside the cell "
525 "where the inverse mapping is not unique.");
526
536 double,
537 int,
538 << "The image of the mapping applied to cell with center ["
539 << arg1 << "] is distorted. The cell geometry or the "
540 << "mapping are invalid, giving a non-positive volume "
541 << "fraction of " << arg2 << " in quadrature point " << arg3
542 << ".");
543
553public:
620 {
621 public:
627
632
636 virtual ~InternalDataBase() = default;
637
653
657 virtual std::size_t
659 };
660
661
662protected:
686 virtual UpdateFlags
687 requires_update_flags(const UpdateFlags update_flags) const = 0;
688
737 virtual std::unique_ptr<InternalDataBase>
738 get_data(const UpdateFlags update_flags,
739 const Quadrature<dim> &quadrature) const = 0;
740
768 virtual std::unique_ptr<InternalDataBase>
769 get_face_data(const UpdateFlags update_flags,
770 const hp::QCollection<dim - 1> &quadrature) const;
771
775 virtual std::unique_ptr<InternalDataBase>
776 get_face_data(const UpdateFlags update_flags,
777 const Quadrature<dim - 1> &quadrature) const;
778
807 virtual std::unique_ptr<InternalDataBase>
808 get_subface_data(const UpdateFlags update_flags,
809 const Quadrature<dim - 1> &quadrature) const = 0;
810
897 const CellSimilarity::Similarity cell_similarity,
898 const Quadrature<dim> & quadrature,
899 const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
901 &output_data) const = 0;
902
927 virtual void
930 const unsigned int face_no,
931 const hp::QCollection<dim - 1> & quadrature,
932 const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
934 &output_data) const;
935
939 virtual void
942 const unsigned int face_no,
943 const Quadrature<dim - 1> & quadrature,
944 const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
946 &output_data) const;
947
974 virtual void
977 const unsigned int face_no,
978 const unsigned int subface_no,
979 const Quadrature<dim - 1> & quadrature,
980 const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
982 &output_data) const = 0;
983
988public:
1057 virtual void
1058 transform(const ArrayView<const Tensor<1, dim>> & input,
1059 const MappingKind kind,
1061 const ArrayView<Tensor<1, spacedim>> &output) const = 0;
1062
1110 virtual void
1112 const MappingKind kind,
1114 const ArrayView<Tensor<2, spacedim>> &output) const = 0;
1115
1168 virtual void
1169 transform(const ArrayView<const Tensor<2, dim>> & input,
1170 const MappingKind kind,
1172 const ArrayView<Tensor<2, spacedim>> &output) const = 0;
1173
1215 virtual void
1217 const MappingKind kind,
1219 const ArrayView<Tensor<3, spacedim>> &output) const = 0;
1220
1268 virtual void
1269 transform(const ArrayView<const Tensor<3, dim>> & input,
1270 const MappingKind kind,
1272 const ArrayView<Tensor<3, spacedim>> &output) const = 0;
1273
1279 // Give class @p FEValues access to the private <tt>get_...data</tt> and
1280 // <tt>fill_fe_...values</tt> functions.
1281 friend class FEValuesBase<dim, spacedim>;
1282 friend class FEValues<dim, spacedim>;
1283 friend class FEFaceValues<dim, spacedim>;
1284 friend class FESubfaceValues<dim, spacedim>;
1285};
1286
1287
1295template <int dim, int spacedim>
1298
1299
1301
1302#endif
UpdateFlags update_each
Definition: mapping.h:652
InternalDataBase(const InternalDataBase &)=delete
virtual ~InternalDataBase()=default
virtual std::size_t memory_consumption() const
Abstract base class for mapping classes.
Definition: mapping.h:304
virtual ~Mapping() override=default
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 is_compatible_with(const ReferenceCell &reference_cell) 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_data(const UpdateFlags update_flags, const Quadrature< dim > &quadrature) const =0
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const =0
virtual BoundingBox< spacedim > get_bounding_box(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const =0
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) 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 Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const =0
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 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 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 std::unique_ptr< InternalDataBase > get_face_data(const UpdateFlags update_flags, const Quadrature< dim - 1 > &quadrature) 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 bool preserves_vertex_locations() const =0
virtual Point< spacedim > get_center(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const bool map_center_of_reference_cell=true) 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 boost::container::small_vector< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
virtual std::unique_ptr< InternalDataBase > get_face_data(const UpdateFlags update_flags, const hp::QCollection< dim - 1 > &quadrature) const
virtual std::unique_ptr< InternalDataBase > get_subface_data(const UpdateFlags update_flags, const Quadrature< dim - 1 > &quadrature) const =0
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
UpdateFlags
#define DeclException0(Exception0)
Definition: exceptions.h:470
static ::ExceptionBase & ExcDistortedMappedCell(Point< spacedim > arg1, double arg2, int arg3)
static ::ExceptionBase & ExcTransformationFailed()
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:493
static ::ExceptionBase & ExcInvalidData()
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:561
MappingKind
Definition: mapping.h:65
@ mapping_piola
Definition: mapping.h:100
@ mapping_covariant_gradient
Definition: mapping.h:86
@ mapping_covariant
Definition: mapping.h:75
@ mapping_nedelec
Definition: mapping.h:115
@ mapping_contravariant
Definition: mapping.h:80
@ mapping_raviart_thomas
Definition: mapping.h:120
@ mapping_contravariant_hessian
Definition: mapping.h:142
@ mapping_covariant_hessian
Definition: mapping.h:136
@ mapping_contravariant_gradient
Definition: mapping.h:92
@ mapping_none
Definition: mapping.h:70
@ mapping_bdm
Definition: mapping.h:125
@ mapping_piola_gradient
Definition: mapping.h:106
@ mapping_piola_hessian
Definition: mapping.h:148
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
Definition: mapping.cc:240
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation