24#ifdef DEAL_II_BOOST_HAS_BROKEN_HEADER_DEPRECATIONS
25# define BOOST_ALLOW_DEPRECATED_HEADERS
27#include <boost/geometry.hpp>
28#ifdef DEAL_II_BOOST_HAS_BROKEN_HEADER_DEPRECATIONS
29# undef BOOST_ALLOW_DEPRECATED_HEADERS
36template <
int dim,
int spacedim>
37boost::container::small_vector<Point<spacedim>,
42 boost::container::small_vector<Point<spacedim>,
45 for (
const unsigned int i : cell->vertex_indices())
53template <
int dim,
int spacedim>
57 const bool map_center_of_reference_cell)
const
59 if (map_center_of_reference_cell)
62 for (
unsigned int d = 0;
d < dim; ++
d)
63 reference_center[d] = .5;
64 return transform_unit_to_real_cell(cell, reference_center);
68 const auto vertices = get_vertices(cell);
78template <
int dim,
int spacedim>
83 if (preserves_vertex_locations())
84 return cell->bounding_box();
91template <
int dim,
int spacedim>
105template <
int dim,
int spacedim>
113 for (
unsigned int i = 0; i < real_points.size(); ++i)
122 unit_points[i][0] = std::numeric_limits<double>::infinity();
129template <
int dim,
int spacedim>
133 const unsigned int face_no,
143 const unsigned int unit_normal_direction =
148 if (unit_normal_direction == 0)
149 return Point<dim - 1>{unit_cell_pt(1)};
150 else if (unit_normal_direction == 1)
151 return Point<dim - 1>{unit_cell_pt(0)};
155 if (unit_normal_direction == 0)
156 return Point<dim - 1>{unit_cell_pt(1), unit_cell_pt(2)};
157 else if (unit_normal_direction == 1)
158 return Point<dim - 1>{unit_cell_pt(0), unit_cell_pt(2)};
159 else if (unit_normal_direction == 2)
160 return Point<dim - 1>{unit_cell_pt(0), unit_cell_pt(1)};
171template <
int dim,
int spacedim>
175 const unsigned int face_no,
183 fill_fe_face_values(cell, face_no, quadrature[0], internal_data, output_data);
188template <
int dim,
int spacedim>
192 const unsigned int face_no,
199 ExcMessage(
"Use of a deprecated interface, please implement "
200 "fill_fe_face_values taking a hp::QCollection argument"));
210template <
int dim,
int spacedim>
211std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
217 return get_face_data(update_flags, quadrature[0]);
222template <
int dim,
int spacedim>
223std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
229 ExcMessage(
"Use of a deprecated interface, please implement "
230 "fill_fe_face_values taking a hp::QCollection argument"));
234 return std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>();
241template <
int dim,
int spacedim>
248template <
int dim,
int spacedim>
252 return sizeof(*this);
258template <
int dim,
int spacedim>
262 const auto &reference_cells =
triangulation.get_reference_cells();
263 Assert(reference_cells.size() == 1,
265 "This function can only work for triangulations that "
266 "use only a single cell type -- for example, only triangles "
267 "or only quadrilaterals. For mixed meshes, there is no "
268 "single linear mapping object that can be used for all "
269 "cells of the triangulation. The triangulation you are "
270 "passing to this function uses multiple cell types."));
272 return reference_cells.front()
273 .template get_default_linear_mapping<dim, spacedim>();
279#include "mapping.inst"
virtual std::size_t memory_consumption() const
Abstract base class for mapping classes.
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 BoundingBox< spacedim > get_bounding_box(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
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 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_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 Point< spacedim > get_center(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const bool map_center_of_reference_cell=true) const
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
unsigned int size() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
@ update_default
No update.
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
Point< 1 > transform_real_to_unit_cell(const std::array< Point< spacedim >, GeometryInfo< 1 >::vertices_per_cell > &vertices, const Point< spacedim > &p)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation