15#ifndef dealii_mapping_h
16#define dealii_mapping_h
39template <
typename ElementType,
typename MemorySpaceType>
43template <
int dim,
int spacedim>
45template <
int dim,
int spacedim>
47template <
int dim,
int spacedim>
49template <
int dim,
int spacedim>
51template <
int dim,
int spacedim>
59 template <
int dim,
int spacedim>
62 template <
int dim,
int spacedim,
typename Number>
316template <
int dim,
int spacedim = dim>
334 virtual std::unique_ptr<Mapping<dim, spacedim>>
351 virtual boost::container::small_vector<Point<spacedim>,
364 boost::container::small_vector<Point<spacedim>,
367 const unsigned int face_no)
const;
394 const bool map_barycenter_of_reference_cell =
true)
const;
522 const unsigned int face_no,
550 "Computing the mapping between a real space point and a point in reference "
551 "space failed, typically because the given point lies outside the cell "
552 "where the inverse mapping is not unique.");
565 <<
"The image of the mapping applied to cell with center ["
566 << arg1 <<
"] is distorted. The cell geometry or the "
567 <<
"mapping are invalid, giving a non-positive volume "
568 <<
"fraction of " << arg2 <<
" in quadrature point " << arg3
775 virtual std::unique_ptr<InternalDataBase>
806 virtual std::unique_ptr<InternalDataBase>
813 virtual std::unique_ptr<InternalDataBase>
845 virtual std::unique_ptr<InternalDataBase>
939 &output_data)
const = 0;
968 const unsigned int face_no,
980 const unsigned int face_no,
1015 const unsigned int face_no,
1016 const unsigned int subface_no,
1020 &output_data)
const = 0;
1034 &output_data)
const;
1334 friend class FEValues<dim, spacedim>;
1339 template <
int,
int,
typename>
1351template <
int dim,
int spacedim>
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.
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 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 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
boost::container::small_vector< Point< spacedim >, GeometryInfo< dim >::vertices_per_face > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no) 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 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 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
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#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)
@ mapping_covariant_gradient
@ mapping_contravariant_hessian
@ mapping_covariant_hessian
@ mapping_contravariant_gradient
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation