Loading [MathJax]/extensions/TeX/newcommand.js
 deal.II version GIT relicensing-3075-gc235bd4825 2025-04-15 08:10:00+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\}}
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages Concepts
mapping.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2001 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_mapping_h
16#define dealii_mapping_h
17
18
19#include <deal.II/base/config.h>
20
23
26
27#include <deal.II/grid/tria.h>
28
30
32
33#include <boost/container/small_vector.hpp>
34
35#include <array>
36#include <cmath>
37#include <memory>
38
40
41template <typename ElementType, typename MemorySpaceType>
42class ArrayView;
43template <int dim>
44class Quadrature;
45template <int dim, int spacedim>
46class FEValues;
47template <int dim, int spacedim>
48class FEValuesBase;
49template <int dim, int spacedim>
50class FEValues;
51template <int dim, int spacedim>
52class FEFaceValues;
53template <int dim, int spacedim>
54class FESubfaceValues;
55namespace NonMatching
56{
57 template <int dim>
59 namespace internal
60 {
61 template <int dim, int spacedim>
63 }
64 template <int dim, int spacedim, typename Number>
65 class MappingInfo;
66} // namespace NonMatching
67
68
166
167
318template <int dim, int spacedim = dim>
320{
321public:
325 virtual ~Mapping() override = default;
326
336 virtual std::unique_ptr<Mapping<dim, spacedim>>
337 clone() const = 0;
338
353 virtual boost::container::small_vector<Point<spacedim>,
354#ifndef _MSC_VER
355 ReferenceCells::max_n_vertices<dim>()
356#else
358#endif
359 >
361 const typename Triangulation<dim, spacedim>::cell_iterator &cell) const;
362
371 boost::container::small_vector<Point<spacedim>,
372#ifndef _MSC_VER
374#else
376#endif
377 >
379 const unsigned int face_no) const;
380
404 virtual Point<spacedim>
406 const bool map_barycenter_of_reference_cell = true) const;
407
427 const typename Triangulation<dim, spacedim>::cell_iterator &cell) const;
428
440 virtual bool
442
447 virtual bool
448 is_compatible_with(const ReferenceCell &reference_cell) const = 0;
449
466 virtual Point<spacedim>
469 const Point<dim> &p) const = 0;
470
499 virtual Point<dim>
502 const Point<spacedim> &p) const = 0;
503
516 virtual void
519 const ArrayView<const Point<spacedim>> &real_points,
520 const ArrayView<Point<dim>> &unit_points) const;
521
532 Point<dim - 1>
535 const unsigned int face_no,
536 const Point<spacedim> &p) const;
537
552
553
563 "Computing the mapping between a real space point and a point in reference "
564 "space failed, typically because the given point lies outside the cell "
565 "where the inverse mapping is not unique.");
566
576 double,
577 int,
578 << "The image of the mapping applied to cell with center ["
579 << arg1 << "] is distorted. The cell geometry or the "
580 << "mapping are invalid, giving a non-positive volume "
581 << "fraction of " << arg2 << " in quadrature point " << arg3
582 << '.');
583
593public:
660 {
661 public:
667
672
682 virtual void
683 reinit(const UpdateFlags update_flags, const Quadrature<dim> &quadrature);
684
688 virtual ~InternalDataBase() = default;
689
705
709 virtual std::size_t
711 };
712
713protected:
737 virtual UpdateFlags
738 requires_update_flags(const UpdateFlags update_flags) const = 0;
739
788 virtual std::unique_ptr<InternalDataBase>
789 get_data(const UpdateFlags update_flags,
790 const Quadrature<dim> &quadrature) const = 0;
791
819 virtual std::unique_ptr<InternalDataBase>
820 get_face_data(const UpdateFlags update_flags,
821 const hp::QCollection<dim - 1> &quadrature) const;
822
826 virtual std::unique_ptr<InternalDataBase>
827 get_face_data(const UpdateFlags update_flags,
828 const Quadrature<dim - 1> &quadrature) const;
829
858 virtual std::unique_ptr<InternalDataBase>
859 get_subface_data(const UpdateFlags update_flags,
860 const Quadrature<dim - 1> &quadrature) const = 0;
861
948 const CellSimilarity::Similarity cell_similarity,
949 const Quadrature<dim> &quadrature,
950 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
952 &output_data) const = 0;
953
978 virtual void
981 const unsigned int face_no,
982 const hp::QCollection<dim - 1> &quadrature,
983 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
985 &output_data) const;
986
990 virtual void
993 const unsigned int face_no,
994 const Quadrature<dim - 1> &quadrature,
995 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
997 &output_data) const;
998
1025 virtual void
1028 const unsigned int face_no,
1029 const unsigned int subface_no,
1030 const Quadrature<dim - 1> &quadrature,
1031 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
1033 &output_data) const = 0;
1034
1041 virtual void
1045 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
1047 &output_data) const;
1048
1053public:
1122 virtual void
1124 const MappingKind kind,
1126 const ArrayView<Tensor<1, spacedim>> &output) const = 0;
1127
1175 virtual void
1177 const MappingKind kind,
1179 const ArrayView<Tensor<2, spacedim>> &output) const = 0;
1180
1233 virtual void
1235 const MappingKind kind,
1237 const ArrayView<Tensor<2, spacedim>> &output) const = 0;
1238
1280 virtual void
1282 const MappingKind kind,
1284 const ArrayView<Tensor<3, spacedim>> &output) const = 0;
1285
1333 virtual void
1335 const MappingKind kind,
1337 const ArrayView<Tensor<3, spacedim>> &output) const = 0;
1338
1344 // Give class @p FEValues access to the private <tt>get_...data</tt> and
1345 // <tt>fill_fe_...values</tt> functions.
1346 friend class FEValuesBase<dim, spacedim>;
1347 friend class FEValues<dim, spacedim>;
1348 friend class FEFaceValues<dim, spacedim>;
1349 friend class FESubfaceValues<dim, spacedim>;
1350 friend class NonMatching::FEImmersedSurfaceValues<dim>;
1351 friend class NonMatching::internal::ComputeMappingDataHelper<dim, spacedim>;
1352 template <int, int, typename>
1354};
1355
1356
1364template <int dim, int spacedim>
1367
1368
1370
1371#endif
UpdateFlags update_each
Definition mapping.h:704
InternalDataBase(const InternalDataBase &)=delete
virtual ~InternalDataBase()=default
virtual std::size_t memory_consumption() const
virtual void reinit(const UpdateFlags update_flags, const Quadrature< dim > &quadrature)
Abstract base class for mapping classes.
Definition mapping.h:320
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
boost::container::small_vector< Point< spacedim >, ReferenceCells::max_n_vertices< dim - 1 >() > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no) 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_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 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 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 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 >, ReferenceCells::max_n_vertices< dim >() > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
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 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 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 Point< spacedim > get_center(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const bool map_barycenter_of_reference_cell=true) 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
Definition point.h:113
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:35
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:36
#define DeclException0(Exception0)
static ::ExceptionBase & ExcDistortedMappedCell(Point< spacedim > arg1, double arg2, int arg3)
static ::ExceptionBase & ExcTransformationFailed()
#define DeclExceptionMsg(Exception, defaulttext)
static ::ExceptionBase & ExcInvalidData()
#define DeclException3(Exception3, type1, type2, type3, outsequence)
UpdateFlags
MappingKind
Definition mapping.h:81
@ mapping_piola
Definition mapping.h:116
@ mapping_covariant_gradient
Definition mapping.h:102
@ mapping_covariant
Definition mapping.h:91
@ mapping_nedelec
Definition mapping.h:131
@ mapping_contravariant
Definition mapping.h:96
@ mapping_raviart_thomas
Definition mapping.h:136
@ mapping_contravariant_hessian
Definition mapping.h:158
@ mapping_covariant_hessian
Definition mapping.h:152
@ mapping_contravariant_gradient
Definition mapping.h:108
@ mapping_none
Definition mapping.h:86
@ mapping_bdm
Definition mapping.h:141
@ mapping_piola_gradient
Definition mapping.h:122
@ mapping_piola_hessian
Definition mapping.h:164
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
Definition mapping.cc:316
constexpr unsigned int max_n_vertices()
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation