57 VTK_QUADRATIC_EDGE = 21,
58 VTK_QUADRATIC_TRIANGLE = 22,
59 VTK_QUADRATIC_QUAD = 23,
60 VTK_QUADRATIC_TETRA = 24,
61 VTK_QUADRATIC_HEXAHEDRON = 25,
62 VTK_QUADRATIC_WEDGE = 26,
63 VTK_QUADRATIC_PYRAMID = 27,
65 VTK_LAGRANGE_CURVE = 68,
66 VTK_LAGRANGE_TRIANGLE = 69,
67 VTK_LAGRANGE_QUADRILATERAL = 70,
68 VTK_LAGRANGE_TETRAHEDRON = 71,
69 VTK_LAGRANGE_HEXAHEDRON = 72,
70 VTK_LAGRANGE_WEDGE = 73,
71 VTK_LAGRANGE_PYRAMID = 74,
125 const unsigned int subface_no)
const
127 if constexpr (dim == 3)
189 static const unsigned int
213 static const RefinementCase<dim - 1> rotated_refinement_case[4] = {
218 const auto [face_orientation, face_rotation, face_flip] =
222 equivalent_refine_case[subface_case][subface_no];
223 const unsigned int equivalent_subface_no =
224 equivalent_subface_number[subface_case][subface_no];
231 (face_orientation == face_rotation ?
235 const unsigned int final_subface_no =
237 final_refinement_case),
239 equivalent_subface_no,
245 return std::make_pair(final_subface_no, final_refinement_case);
249 (void)combined_face_orientation;
260template <
int dim,
int spacedim>
261std::unique_ptr<Mapping<dim, spacedim>>
267 return std::make_unique<MappingQ<dim, spacedim>>(degree);
269 return std::make_unique<MappingFE<dim, spacedim>>(
272 return std::make_unique<MappingFE<dim, spacedim>>(
275 return std::make_unique<MappingFE<dim, spacedim>>(
282 return std::make_unique<MappingQ<dim, spacedim>>(degree);
287template <
int dim,
int spacedim>
355 const auto create_quadrature = [](
const ReferenceCell &reference_cell) {
356 std::vector<Point<dim>> vertices(reference_cell.n_vertices());
357 for (
const unsigned int v : reference_cell.vertex_indices())
358 vertices[v] = reference_cell.vertex<dim>(v);
404 constexpr std::array<unsigned int, 4> exodus_to_deal{{0, 1, 3, 2}};
405 return exodus_to_deal[vertex_n];
411 constexpr std::array<unsigned int, 8> exodus_to_deal{
412 {0, 1, 3, 2, 4, 5, 7, 6}};
413 return exodus_to_deal[vertex_n];
417 constexpr std::array<unsigned int, 6> exodus_to_deal{
419 return exodus_to_deal[vertex_n];
423 constexpr std::array<unsigned int, 5> exodus_to_deal{{0, 1, 3, 2, 4}};
424 return exodus_to_deal[vertex_n];
449 constexpr std::array<unsigned int, 4> exodus_to_deal{{2, 1, 3, 0}};
450 return exodus_to_deal[face_n];
454 constexpr std::array<unsigned int, 4> exodus_to_deal{{1, 3, 2, 0}};
455 return exodus_to_deal[face_n];
459 constexpr std::array<unsigned int, 6> exodus_to_deal{
461 return exodus_to_deal[face_n];
465 constexpr std::array<unsigned int, 6> exodus_to_deal{{3, 4, 2, 0, 1}};
466 return exodus_to_deal[face_n];
470 constexpr std::array<unsigned int, 5> exodus_to_deal{{3, 2, 4, 1, 0}};
471 return exodus_to_deal[face_n];
501 constexpr std::array<unsigned int, 4> unv_to_deal{{1, 0, 2, 3}};
502 return unv_to_deal[vertex_n];
506 constexpr std::array<unsigned int, 8> unv_to_deal{
507 {6, 7, 5, 4, 2, 3, 1, 0}};
508 return unv_to_deal[vertex_n];
524 return VTKCellType::VTK_VERTEX;
526 return VTKCellType::VTK_LINE;
528 return VTKCellType::VTK_TRIANGLE;
530 return VTKCellType::VTK_QUAD;
532 return VTKCellType::VTK_TETRA;
534 return VTKCellType::VTK_PYRAMID;
536 return VTKCellType::VTK_WEDGE;
538 return VTKCellType::VTK_HEXAHEDRON;
540 return VTKCellType::VTK_INVALID;
545 return VTKCellType::VTK_INVALID;
556 return VTKCellType::VTK_VERTEX;
558 return VTKCellType::VTK_QUADRATIC_EDGE;
560 return VTKCellType::VTK_QUADRATIC_TRIANGLE;
562 return VTKCellType::VTK_QUADRATIC_QUAD;
564 return VTKCellType::VTK_QUADRATIC_TETRA;
566 return VTKCellType::VTK_QUADRATIC_PYRAMID;
568 return VTKCellType::VTK_QUADRATIC_WEDGE;
570 return VTKCellType::VTK_QUADRATIC_HEXAHEDRON;
572 return VTKCellType::VTK_INVALID;
577 return VTKCellType::VTK_INVALID;
588 return VTKCellType::VTK_VERTEX;
590 return VTKCellType::VTK_LAGRANGE_CURVE;
592 return VTKCellType::VTK_LAGRANGE_TRIANGLE;
594 return VTKCellType::VTK_LAGRANGE_QUADRILATERAL;
596 return VTKCellType::VTK_LAGRANGE_TETRAHEDRON;
598 return VTKCellType::VTK_LAGRANGE_PYRAMID;
600 return VTKCellType::VTK_LAGRANGE_WEDGE;
602 return VTKCellType::VTK_LAGRANGE_HEXAHEDRON;
604 return VTKCellType::VTK_INVALID;
609 return VTKCellType::VTK_INVALID;
617 const std::array<unsigned, 0> &,
618 const std::array<unsigned, 0> &,
634 const std::array<unsigned, 1> &node_indices,
635 const std::array<unsigned, 1> &nodes_per_direction,
638 const unsigned int i = node_indices[0];
640 const bool ibdy = (i == 0 || i == nodes_per_direction[0]);
642 const int nbdy = (ibdy ? 1 : 0);
649 const int offset = 2;
650 return (i - 1) + offset;
662 const std::array<unsigned, 2> &node_indices,
663 const std::array<unsigned, 2> &nodes_per_direction,
668 const unsigned int i = node_indices[0];
669 const unsigned int j = node_indices[1];
671 const bool ibdy = (i == 0 || i == nodes_per_direction[0]);
672 const bool jbdy = (j == 0 || j == nodes_per_direction[1]);
674 const int nbdy = (ibdy ? 1 : 0) + (jbdy ? 1 : 0);
678 return (i != 0u ? (j != 0u ? 2 : 1) : (j != 0u ? 3 : 0));
688 nodes_per_direction[0] - 1 + nodes_per_direction[1] - 1 :
696 (i != 0u ? nodes_per_direction[0] - 1 :
697 2 * (nodes_per_direction[0] - 1) +
698 nodes_per_direction[1] - 1) +
703 offset += 2 * (nodes_per_direction[0] - 1 + nodes_per_direction[1] - 1);
705 return offset + (i - 1) + (nodes_per_direction[0] - 1) * ((j - 1));
723 const std::array<unsigned, 3> &node_indices,
724 const std::array<unsigned, 3> &nodes_per_direction,
725 const bool legacy_format)
const
729 const unsigned int i = node_indices[0];
730 const unsigned int j = node_indices[1];
731 const unsigned int k = node_indices[2];
733 const bool ibdy = (i == 0 || i == nodes_per_direction[0]);
734 const bool jbdy = (j == 0 || j == nodes_per_direction[1]);
735 const bool kbdy = (k == 0 || k == nodes_per_direction[2]);
737 const int nbdy = (ibdy ? 1 : 0) + (jbdy ? 1 : 0) + (kbdy ? 1 : 0);
741 return (i != 0u ? (j != 0u ? 2 : 1) : (j != 0u ? 3 : 0)) +
752 nodes_per_direction[0] - 1 + nodes_per_direction[1] - 1 :
754 (k != 0u ? 2 * (nodes_per_direction[0] - 1 +
755 nodes_per_direction[1] - 1) :
762 (i != 0u ? nodes_per_direction[0] - 1 :
763 2 * (nodes_per_direction[0] - 1) +
764 nodes_per_direction[1] - 1) +
765 (k != 0u ? 2 * (nodes_per_direction[0] - 1 +
766 nodes_per_direction[1] - 1) :
772 4 * (nodes_per_direction[0] - 1) + 4 * (nodes_per_direction[1] - 1);
775 (nodes_per_direction[2] - 1) *
776 (i != 0u ? (j != 0u ? 3 : 1) : (j != 0u ? 2 : 0)) +
780 (nodes_per_direction[2] - 1) *
781 (i != 0u ? (j != 0u ? 2 : 1) : (j != 0u ? 3 : 0)) +
785 offset += 4 * (nodes_per_direction[0] - 1 + nodes_per_direction[1] - 1 +
786 nodes_per_direction[2] - 1);
791 return (j - 1) + ((nodes_per_direction[1] - 1) * (k - 1)) +
792 (i != 0u ? (nodes_per_direction[1] - 1) *
793 (nodes_per_direction[2] - 1) :
797 offset += 2 * (nodes_per_direction[1] - 1) * (nodes_per_direction[2] - 1);
800 return (i - 1) + ((nodes_per_direction[0] - 1) * (k - 1)) +
801 (j != 0u ? (nodes_per_direction[2] - 1) *
802 (nodes_per_direction[0] - 1) :
806 offset += 2 * (nodes_per_direction[2] - 1) * (nodes_per_direction[0] - 1);
808 return (i - 1) + ((nodes_per_direction[0] - 1) * (j - 1)) +
810 (nodes_per_direction[0] - 1) * (nodes_per_direction[1] - 1) :
816 offset += 2 * ((nodes_per_direction[1] - 1) * (nodes_per_direction[2] - 1) +
817 (nodes_per_direction[2] - 1) * (nodes_per_direction[0] - 1) +
818 (nodes_per_direction[0] - 1) * (nodes_per_direction[1] - 1));
819 return offset + (i - 1) +
820 (nodes_per_direction[0] - 1) *
821 ((j - 1) + (nodes_per_direction[1] - 1) * ((k - 1)));
847 static constexpr std::array<unsigned int, 4> index_translation_table =
849 return index_translation_table[vertex_index];
855 static constexpr std::array<unsigned int, 5> index_translation_table =
857 return index_translation_table[vertex_index];
863 static constexpr std::array<unsigned int, 8> index_translation_table =
864 {{0, 1, 3, 2, 4, 5, 7, 6}};
865 return index_translation_table[vertex_index];
970 std::pair<Point<dim>,
double>
981 const double t = ((x1 - x0) * (p - x0)) / ((x1 - x0).norm_square());
988 const auto p2 = x0 + t * (x1 - x0);
989 return std::make_pair(p2, p2.distance_square(p));
997 std::pair<Point<dim>,
double>
998 project_to_quad(
const std::array<
Point<dim>, 3> & ,
1004 std::numeric_limits<double>::signaling_NaN());
1025 std::pair<Point<3>,
double>
1026 project_to_quad(
const std::array<
Point<3>, 3> &vertices,
1037 std::array<Point<3>, 3> shifted_vertices = vertices;
1039 for (
Point<3> &shifted_vertex : shifted_vertices)
1040 shifted_vertex +=
shift;
1049 const Point<3> vertex = shifted_vertices[0];
1059 e0 = shifted_vertices[1] - shifted_vertices[0];
1060 e1 = shifted_vertices[2] - shifted_vertices[0];
1065 e0 = shifted_vertices[1] - shifted_vertices[0];
1066 e1 = shifted_vertices[2] - shifted_vertices[0];
1071 const double c0 = e0 * (shifted_p - vertex) / e0.
norm_square();
1072 const double c1 = e1 * (shifted_p - vertex) / e1.norm_square();
1073 const Point<3> projected_shifted_p = vertex + c0 * e0 + c1 * e1;
1075 bool in_quad =
false;
1079 for (
unsigned int i = 0; i < 3; ++i)
1080 shifted_vertex_matrix[i] = shifted_vertices[i];
1083 bool is_convex_combination =
true;
1084 for (
unsigned int i = 0; i < 3; ++i)
1085 is_convex_combination = is_convex_combination &&
1086 (0.0 <= combination_coordinates[i]) &&
1087 (combination_coordinates[i] <= 1.0);
1088 in_quad = is_convex_combination;
1091 in_quad = (0.0 <= c0 && c0 <= 1.0 && 0.0 <= c1 && c1 <= 1.0);
1094 return std::make_pair(projected_shifted_p - shift,
1097 return std::make_pair(
Point<3>(), std::numeric_limits<double>::max());
1120 unsigned int closest_vertex_no = 0;
1121 double closest_vertex_distance_square =
vertex<dim>(0).distance_square(p);
1122 for (
unsigned int i = 1; i <
n_vertices(); ++i)
1124 const double new_vertex_distance_square =
1126 if (new_vertex_distance_square < closest_vertex_distance_square)
1128 closest_vertex_no = i;
1129 closest_vertex_distance_square = new_vertex_distance_square;
1133 double min_distance_square = std::numeric_limits<double>::max();
1136 for (
const unsigned int face_no :
1142 auto pair = project_to_line(v0, v1, p);
1143 if (pair.second < min_distance_square)
1145 result = pair.first;
1146 min_distance_square = pair.second;
1163 const std::array<unsigned int, 5> all_pyramid_faces{{0, 1, 2, 3, 4}};
1167 for (
const unsigned int face_no : faces)
1172 std::array<Point<dim>, 3> vertices;
1173 for (
unsigned int vertex_no = 0; vertex_no < 3; ++vertex_no)
1177 auto pair = project_to_quad(vertices, p, face_cell);
1178 if (pair.second < min_distance_square)
1180 result = pair.first;
1181 min_distance_square = pair.second;
1185 for (
const unsigned int face_no :
1189 for (
const unsigned int face_line_no : face_cell.line_indices())
1191 const auto cell_line_no =
1199 auto pair = project_to_line(v0, v1, p);
1200 if (pair.second < min_distance_square)
1202 result = pair.first;
1203 min_distance_square = pair.second;
1209 Assert(min_distance_square < std::numeric_limits<double>::max(),
1218 constexpr unsigned int x_index = 0;
1219 constexpr unsigned int y_index = (dim >= 2 ? 1 : 0);
1220 constexpr unsigned int z_index = (dim >= 3 ? 2 : 0);
1231 for (
unsigned int d = 0; d < dim; ++d)
1232 result[d] = std::clamp(result[d], 0.0, 1.0);
1236 result[x_index] = std::clamp(result[x_index], 0.0, 1.0);
1238 std::clamp(result[y_index], 0.0, 1.0 - result[x_index]);
1241 result[x_index] = std::clamp(result[x_index], 0.0, 1.0);
1243 std::clamp(result[y_index], 0.0, 1.0 - result[x_index]);
1245 std::clamp(result[z_index],
1247 1.0 - result[x_index] - result[y_index]);
1251 result[x_index] = std::clamp(result[x_index], 0.0, 1.0);
1253 std::clamp(result[y_index], 0.0, 1.0 - result[x_index]);
1254 result[z_index] = std::clamp(result[z_index], 0.0, 1.0);
1258 result[x_index] = std::clamp(result[x_index], -1.0, 1.0);
1259 result[y_index] = std::clamp(result[y_index], -1.0, 1.0);
1262 const auto x_abs =
std::abs(result[x_index]);
1263 const auto y_abs =
std::abs(result[y_index]);
1266 result[z_index] = std::clamp(result[z_index], 0.0, 1.0 - x_abs);
1268 result[z_index] = std::clamp(result[z_index], 0.0, 1.0 - y_abs);
1293 out << static_cast<unsigned int>(reference_cell.kind);
1307 reference_cell.kind =
static_cast<decltype(reference_cell.
kind)
>(
value);
1321 "The reference cell kind just read does not correspond to one of the "
1322 "valid choices. There must be an error."));
1334#include "grid/reference_cell.inst"
Implementation of the classic affine transformation mapping used for simplices.
Abstract base class for mapping classes.
constexpr numbers::NumberTraits< Number >::real_type distance_square(const Point< dim, Number > &p) const
ArrayView< const unsigned int > faces_for_given_vertex(const unsigned int vertex_index) const
unsigned int n_vertices() const
Quadrature< dim > get_gauss_type_quadrature(const unsigned n_points_1d) const
bool is_hyper_cube() const
unsigned int vtk_linear_type() const
std::pair< unsigned int, RefinementCase< dim - 1 > > equivalent_refinement_case(const types::geometric_orientation combined_face_orientation, const internal::SubfaceCase< dim > subface_case, const unsigned int subface_no) const
Point< dim > vertex(const unsigned int v) const
const Quadrature< dim > & get_nodal_type_quadrature() const
unsigned int face_to_cell_lines(const unsigned int face, const unsigned int line, const types::geometric_orientation face_orientation) const
static constexpr ndarray< unsigned int, 2, 2 > line_vertex_permutations
unsigned int vtk_lexicographic_to_node_index(const std::array< unsigned, dim > &node_indices, const std::array< unsigned, dim > &nodes_per_direction, const bool legacy_format) const
bool contains_point(const Point< dim > &p, const double tolerance=0) const
unsigned int gmsh_element_type() const
static constexpr ndarray< unsigned int, 8, 4 > quadrilateral_vertex_permutations
unsigned int n_faces() const
unsigned int line_to_cell_vertices(const unsigned int line, const unsigned int vertex) 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
unsigned int vtk_vertex_to_deal_vertex(const unsigned int vertex_index) const
unsigned int face_to_cell_vertices(const unsigned int face, const unsigned int vertex, const types::geometric_orientation face_orientation) const
std::unique_ptr< Mapping< dim, spacedim > > get_default_mapping(const unsigned int degree) const
unsigned int vtk_quadratic_type() const
unsigned int vtk_lagrange_type() 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
const Mapping< dim, spacedim > & get_default_linear_mapping() const
std::string to_string() const
Point< dim > closest_point(const Point< dim > &p) const
static constexpr ndarray< unsigned int, 6, 3 > triangle_vertex_permutations
constexpr numbers::NumberTraits< Number >::real_type norm_square() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_ASSERT_UNREACHABLE()
#define DEAL_II_NOT_IMPLEMENTED()
static ::ExceptionBase & ExcIO()
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
const bool IsBlockVector< VectorType >::value
MappingQ< dim, spacedim > StaticMappingQ1< dim, spacedim >::mapping
constexpr ReferenceCell Triangle
constexpr ReferenceCell Hexahedron
constexpr ReferenceCell Pyramid
constexpr ReferenceCell Wedge
constexpr ReferenceCell Vertex
constexpr ReferenceCell Invalid
constexpr ReferenceCell Quadrilateral
constexpr ReferenceCell Tetrahedron
constexpr ReferenceCell Line
std::tuple< bool, bool, bool > split_face_orientation(const types::geometric_orientation combined_face_orientation)
constexpr unsigned int invalid_unsigned_int
constexpr types::geometric_orientation default_geometric_orientation
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
unsigned char geometric_orientation
typename internal::ndarray::HelperArray< T, Ns... >::type ndarray
std::istream & operator>>(std::istream &in, ReferenceCell &reference_cell)
std::ostream & operator<<(std::ostream &out, const ReferenceCell &reference_cell)
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 constexpr unsigned int max_children_per_face
static MappingQ< dim, spacedim > mapping
DEAL_II_HOST constexpr SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &)