16#ifndef dealii_tria_reference_cell_h
17#define dealii_tria_reference_cell_h
32template <
int dim,
int spacedim>
59 DEAL_II_CONSTEXPR ::ReferenceCell
161 const unsigned int i)
const;
172 template <
int dim,
int spacedim = dim>
173 std::unique_ptr<Mapping<dim, spacedim>>
187 template <
int dim,
int spacedim = dim>
334 const unsigned int subface_n,
335 const unsigned char face_orientation = 1)
const;
346 std::array<unsigned int, 2>
358 std::array<unsigned int, 2>
366 const unsigned int line,
367 const unsigned char face_orientation)
const;
374 const unsigned int vertex,
375 const unsigned char face_orientation)
const;
382 const unsigned int face,
383 const unsigned char face_orientation)
const;
390 const unsigned int face,
391 const unsigned char face_orientation)
const;
405 const unsigned char face_orientation,
406 const unsigned char line_orientation)
const;
428 const unsigned int i)
const;
441 template <
typename T, std::
size_t N>
444 const std::array<T, N> &vertices_1)
const;
449 template <
typename T, std::
size_t N>
452 const unsigned int orientation)
const;
526 constexpr operator std::uint8_t()
const;
545 template <
class Archive>
547 serialize(Archive &archive,
const unsigned int );
583inline constexpr ReferenceCell::operator std::uint8_t()
const
610 inline DEAL_II_CONSTEXPR ::ReferenceCell
615 Assert((kind ==
static_cast<std::uint8_t
>(-1)) || (kind < 8),
654 static_cast<std::uint8_t
>(-1));
684template <
class Archive>
715 {{{0, 2}}, {{0, 1}}, {{1, 2}}}};
717 return table[vertex];
723 {{{0, 1, 2}}, {{0, 1, 3}}, {{0, 2, 3}}, {{1, 2, 3}}}};
725 return table[vertex];
737 return table[vertex];
749 return {&table[vertex][0], vertex == 4 ? 4u : 3u};
943 const unsigned int face,
944 const unsigned int subface,
945 const unsigned char face_orientation_raw)
const
960 {{{0, 1}}, {{1, 2}}, {{2, 0}}}};
962 return subcells[face][subface];
1011inline std::array<unsigned int, 2>
1013 const unsigned int vertex)
const
1028 {{{0, 0}}, {{0, 1}}, {{1, 1}}}};
1030 return table[vertex];
1039 {{{0, 0}}, {{0, 1}}, {{0, 2}}, {{1, 2}}}};
1041 return table[vertex];
1046 {{{0, 0}}, {{0, 1}}, {{0, 2}}, {{0, 3}}, {{1, 2}}}};
1048 return table[vertex];
1053 {{{0, 1}}, {{0, 0}}, {{0, 2}}, {{1, 0}}, {{1, 1}}, {{1, 2}}}};
1055 return table[vertex];
1068inline std::array<unsigned int, 2>
1070 const unsigned int line)
const
1092 static const std::array<unsigned int, 2> table[6] = {
1093 {{0, 0}}, {{0, 1}}, {{0, 2}}, {{1, 1}}, {{1, 2}}, {{2, 1}}};
1099 static const std::array<unsigned int, 2> table[8] = {{{0, 0}},
1112 static const std::array<unsigned int, 2> table[9] = {{{0, 0}},
1137 const unsigned int line,
1138 const unsigned char face_orientation)
const
1172 {{{0, 1, 2}}, {{0, 3, 4}}, {{2, 5, 3}}, {{1, 4, 5}}}};
1203 const unsigned int vertex,
1204 const unsigned char face_orientation)
const
1225 {{{0, 1}}, {{1, 2}}, {{2, 0}}}};
1227 return table[face][face_orientation ? vertex : (1 - vertex)];
1241 {{{0, 1, 2}}, {{1, 0, 3}}, {{0, 2, 3}}, {{2, 1, 3}}}};
1244 vertex, face, face_orientation)];
1256 vertex, face, face_orientation)];
1268 vertex, face, face_orientation)];
1288 const unsigned int vertex,
1289 const unsigned int face,
1290 const unsigned char face_orientation)
const
1307 return table[face_orientation][vertex];
1323 return table[face_orientation][vertex];
1344 return table[face_orientation][vertex];
1366 return table[face_orientation][vertex];
1386 const unsigned int line,
1387 const unsigned int face,
1388 const unsigned char face_orientation)
const
1418 return table[face_orientation][line];
1439 return table[face_orientation][line];
1461 return table[face_orientation][line];
1560 " does not correspond to a known reference cell type."));
1569 const unsigned int i)
const
1572 if (*
this == ReferenceCells::get_hypercube<dim>())
1625 const double Q14 = 0.25;
1628 const double r = xi[
std::min(0, dim - 1)];
1629 const double s = xi[
std::min(1, dim - 1)];
1630 const double t = xi[
std::min(2, dim - 1)];
1632 if (
fabs(t - 1.0) > 1.0e-14)
1634 ration = (r * s * t) / (1.0 - t);
1642 return Q14 * ((1.0 - r) * (1.0 - s) - t + ration);
1644 return Q14 * ((1.0 + r) * (1.0 - s) - t - ration);
1646 return Q14 * ((1.0 - r) * (1.0 + s) - t - ration);
1648 return Q14 * ((1.0 + r) * (1.0 + s) - t + ration);
1663 const unsigned int i)
const
1666 if (*
this == ReferenceCells::get_hypercube<dim>())
1693 const unsigned int i)
const
1698 if (*
this == ReferenceCells::get_hypercube<dim>())
1706 static const std::array<Tensor<1, dim>, 3> table = {
1711 return table[face_no];
1725 +
std::pow(1.0 / 3.0, 1.0 / 4.0))}}}};
1727 return table[face_no][i];
1740 return table[face_no][i];
1756 return table[face_no][i];
1787 unit_tangential_vectors<dim>(face_no, 1));
1799 const unsigned int line,
1800 const unsigned char face_orientation_raw,
1801 const unsigned char line_orientation)
const
1805 static const bool bool_table[2][2][2][2] = {
1817 {{
true,
false}, {
false,
true}}}};
1823 return (
static_cast<bool>(line_orientation) ==
1824 bool_table[line / 2][face_orientation][face_flip][face_rotation]);
1836 template <
typename T, std::
size_t N>
1866 for (
unsigned int i = 0; i < n_vertices; ++i)
1869 if (i + 1 != n_vertices)
1873 out <<
"] is not a permutation of [";
1875 for (
unsigned int i = 0; i < n_vertices; ++i)
1878 if (i + 1 != n_vertices)
1882 out <<
"]." << std::endl;
1904template <
typename T, std::
size_t N>
1907 const std::array<T, N> &vertices_1)
const
1912 const std::array<T, 2> i{{vertices_0[0], vertices_0[1]}};
1913 const std::array<T, 2> j{{vertices_1[0], vertices_1[1]}};
1916 if (i == std::array<T, 2>{{j[0], j[1]}})
1920 if (i == std::array<T, 2>{{j[1], j[0]}})
1925 const std::array<T, 3> i{{vertices_0[0], vertices_0[1], vertices_0[2]}};
1926 const std::array<T, 3> j{{vertices_1[0], vertices_1[1], vertices_1[2]}};
1929 if (i == std::array<T, 3>{{j[0], j[1], j[2]}})
1933 if (i == std::array<T, 3>{{j[1], j[2], j[0]}})
1937 if (i == std::array<T, 3>{{j[2], j[0], j[1]}})
1941 if (i == std::array<T, 3>{{j[0], j[2], j[1]}})
1945 if (i == std::array<T, 3>{{j[2], j[1], j[0]}})
1949 if (i == std::array<T, 3>{{j[1], j[0], j[2]}})
1954 const std::array<T, 4> i{
1955 {vertices_0[0], vertices_0[1], vertices_0[2], vertices_0[3]}};
1956 const std::array<T, 4> j{
1957 {vertices_1[0], vertices_1[1], vertices_1[2], vertices_1[3]}};
1960 if (i == std::array<T, 4>{{j[0], j[1], j[2], j[3]}})
1964 if (i == std::array<T, 4>{{j[2], j[0], j[3], j[1]}})
1968 if (i == std::array<T, 4>{{j[3], j[2], j[1], j[0]}})
1972 if (i == std::array<T, 4>{{j[1], j[3], j[0], j[2]}})
1976 if (i == std::array<T, 4>{{j[0], j[2], j[1], j[3]}})
1980 if (i == std::array<T, 4>{{j[2], j[3], j[0], j[1]}})
1984 if (i == std::array<T, 4>{{j[3], j[1], j[2], j[0]}})
1988 if (i == std::array<T, 4>{{j[1], j[0], j[3], j[2]}})
1999template <
typename T, std::
size_t N>
2000inline std::array<T, N>
2003 const unsigned int orientation)
const
2005 std::array<T, 4> temp;
2009 switch (orientation)
2023 switch (orientation)
2049 switch (orientation)
2084 std::array<T, N> temp_;
2085 std::copy_n(temp.begin(),
N, temp_.begin());
Abstract base class for mapping classes.
ArrayView< const unsigned int > faces_for_given_vertex(const unsigned int vertex_index) const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > vertex_indices() const
unsigned int n_vertices() const
bool is_hyper_cube() const
unsigned int vtk_linear_type() const
bool standard_vs_true_line_orientation(const unsigned int line, const unsigned char face_orientation, const unsigned char line_orientation) const
void serialize(Archive &archive, const unsigned int)
unsigned int child_cell_on_face(const unsigned int face_n, const unsigned int subface_n, const unsigned char face_orientation=1) const
unsigned int face_to_cell_vertices(const unsigned int face, const unsigned int vertex, const unsigned char face_orientation) const
unsigned char compute_orientation(const std::array< T, N > &vertices_0, const std::array< T, N > &vertices_1) const
std::array< unsigned int, 2 > standard_vertex_to_face_and_vertex_index(const unsigned int vertex) const
Quadrature< dim > get_gauss_type_quadrature(const unsigned n_points_1D) const
unsigned int standard_to_real_face_vertex(const unsigned int vertex, const unsigned int face, const unsigned char face_orientation) const
std::array< unsigned int, 2 > standard_line_to_face_and_line_index(const unsigned int line) const
double d_linear_shape_function(const Point< dim > &xi, const unsigned int i) const
Tensor< 1, dim > unit_tangential_vectors(const unsigned int face_no, const unsigned int i) const
std::array< T, N > permute_according_orientation(const std::array< T, N > &vertices, const unsigned int orientation) const
unsigned int n_faces() const
constexpr ReferenceCell()
unsigned int exodusii_vertex_to_deal_vertex(const unsigned int vertex_n) const
unsigned int unv_vertex_to_deal_vertex(const unsigned int vertex_n) const
std::unique_ptr< Mapping< dim, spacedim > > get_default_mapping(const unsigned int degree) const
unsigned int vtk_quadratic_type() const
unsigned int n_lines() const
unsigned int vtk_lagrange_type() const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > line_indices() const
unsigned int get_dimension() const
ReferenceCell face_reference_cell(const unsigned int face_no) const
unsigned int exodusii_face_to_deal_face(const unsigned int face_n) const
static ReferenceCell n_vertices_to_type(const int dim, const unsigned int n_vertices)
unsigned int standard_to_real_face_line(const unsigned int line, const unsigned int face, const unsigned char face_orientation) const
const Mapping< dim, spacedim > & get_default_linear_mapping() const
std::string to_string() const
unsigned int face_to_cell_lines(const unsigned int face, const unsigned int line, const unsigned char face_orientation) const
constexpr bool operator!=(const ReferenceCell &type) const
constexpr bool operator==(const ReferenceCell &type) const
const Quadrature< dim > & get_nodal_type_quadrature() const
Tensor< 1, dim > unit_normal_vectors(const unsigned int face_no) const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > face_indices() const
Tensor< 1, dim > d_linear_shape_function_gradient(const Point< dim > &xi, const unsigned int i) const
const ::ReferenceCell entity_type
const std::array< T, N > vertices_0
const std::array< T, N > vertices_1
virtual void print_info(std::ostream &out) const override
virtual ~NoPermutation() noexcept override=default
NoPermutation(const ::ReferenceCell &entity_type, const std::array< T, N > &vertices_0, const std::array< T, N > &vertices_1)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_CONSTEXPR
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
std::string to_string(const T &t)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
Expression fabs(const Expression &x)
constexpr const ReferenceCell Tetrahedron
constexpr const ReferenceCell Quadrilateral
constexpr const ReferenceCell Wedge
constexpr const ReferenceCell Pyramid
constexpr const ReferenceCell Invalid
constexpr const ReferenceCell Triangle
constexpr const ReferenceCell Hexahedron
constexpr const ReferenceCell & get_hypercube()
constexpr const ReferenceCell Vertex
constexpr const ReferenceCell Line
constexpr const ReferenceCell & get_simplex()
bool get_bit(const unsigned char number, const unsigned int n)
constexpr ::ReferenceCell make_reference_cell_from_int(const std::uint8_t kind)
static const unsigned int invalid_unsigned_int
boost::integer_range< IncrementableType > iota_view
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
typename internal::ndarray::HelperArray< T, Ns... >::type ndarray
static unsigned int standard_to_real_line_vertex(const unsigned int vertex, const bool line_orientation=true)
static unsigned int child_cell_on_face(const RefinementCase< dim > &ref_case, const unsigned int face, const unsigned int subface, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false, const RefinementCase< dim - 1 > &face_refinement_case=RefinementCase< dim - 1 >::isotropic_refinement)
static unsigned int standard_to_real_face_vertex(const unsigned int vertex, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static unsigned int face_to_cell_vertices(const unsigned int face, const unsigned int vertex, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static std::array< unsigned int, 2 > standard_hex_line_to_quad_line_index(const unsigned int line)
static unsigned int standard_to_real_face_line(const unsigned int line, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static std::array< unsigned int, 2 > standard_hex_vertex_to_quad_vertex_index(const unsigned int vertex)
static std::array< unsigned int, 2 > standard_quad_vertex_to_line_vertex_index(const unsigned int vertex)
static double d_linear_shape_function(const Point< dim > &xi, const unsigned int i)
static unsigned int face_to_cell_lines(const unsigned int face, const unsigned int line, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static Tensor< 1, dim > d_linear_shape_function_gradient(const Point< dim > &xi, const unsigned int i)
constexpr Tensor< 1, dim, Number > cross_product_2d(const Tensor< 1, dim, Number > &src)
constexpr Tensor< 1, dim, typename ProductType< Number1, Number2 >::type > cross_product_3d(const Tensor< 1, dim, Number1 > &src1, const Tensor< 1, dim, Number2 > &src2)