permute_by_combined_orientation(const ArrayView<const T> &vertices,
const unsigned char orientation) const;
+ /**
+ * Return the inverse orientation. This is the value such that calling
+ * permute_by_combined_orientation() with <tt>o</tt> and then calling it again
+ * with get_inverse_combined_orientation(o) is the identity operation.
+ */
+ unsigned char
+ get_inverse_combined_orientation(const unsigned char orientation) const;
+
/**
* Return a vector of faces a given @p vertex_index belongs to.
*/
+inline unsigned char
+ReferenceCell::get_inverse_combined_orientation(
+ const unsigned char orientation) const
+{
+ switch (this->kind)
+ {
+ case ReferenceCells::Vertex:
+ // Things are always default-oriented in 1D
+ return orientation;
+
+ case ReferenceCells::Line:
+ // the 1d orientations are the identity and a flip: i.e., the identity
+ // and an involutory mapping
+ return orientation;
+
+ case ReferenceCells::Triangle:
+ {
+ AssertIndexRange(orientation, 6);
+ constexpr std::array<unsigned char, 6> inverses{{0, 1, 2, 5, 4, 3}};
+ return inverses[orientation];
+ }
+ case ReferenceCells::Quadrilateral:
+ {
+ AssertIndexRange(orientation, 8);
+ constexpr std::array<unsigned char, 8> inverses{
+ {0, 1, 2, 7, 4, 5, 6, 3}};
+ return inverses[orientation];
+ }
+ default:
+ DEAL_II_NOT_IMPLEMENTED();
+ }
+
+ return std::numeric_limits<unsigned char>::max();
+}
+
+
+
template <>
unsigned int
ReferenceCell::vtk_lexicographic_to_node_index<0>(
--- /dev/null
+// ------------------------------------------------------------------------
+//
+// SPDX-License-Identifier: LGPL-2.1-or-later
+// Copyright (C) 2024 - 2024 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// Part of the source code is dual licensed under Apache-2.0 WITH
+// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
+// governing the source code and code contributions can be found in
+// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
+//
+// ------------------------------------------------------------------------
+
+
+// Test that we get the correct reversed orientations (i.e., the composition of
+// permutations is the identity operation)
+
+#include <deal.II/grid/reference_cell.h>
+#include <deal.II/grid/tria_orientation.h>
+
+#include "../tests.h"
+
+int
+main()
+{
+ initlog();
+
+ deallog << "lines" << std::endl;
+ {
+ const std::array<unsigned char, 2> inverse_permutations{{0u, 1u}};
+
+ for (unsigned char o = 0; o < 2; ++o)
+ {
+ deallog << "o = " << int(o) << std::endl;
+ std::array<unsigned int, 2> vs{{0u, 1u}};
+
+ const auto result1 =
+ ReferenceCells::Line.permute_by_combined_orientation(
+ make_const_array_view(vs), o);
+
+ ArrayView<const unsigned int> view1(result1.data(), result1.size());
+ const auto result2 =
+ ReferenceCells::Line.permute_by_combined_orientation(
+ view1, inverse_permutations[o]);
+
+ for (const auto &v : result2)
+ deallog << " " << v << std::endl;
+ }
+ }
+
+ deallog << "triangles" << std::endl;
+ {
+ const std::array<unsigned char, 6> inverse_permutations{
+ {0u, 1u, 2u, 5u, 4u, 3u}};
+
+ for (unsigned char o = 0; o < 6; ++o)
+ {
+ deallog << "o = " << int(o) << std::endl;
+ std::array<unsigned int, 3> vs{{0u, 1u, 2u}};
+
+ const auto result1 =
+ ReferenceCells::Triangle.permute_by_combined_orientation(
+ make_const_array_view(vs), o);
+
+ ArrayView<const unsigned int> view1(result1.data(), result1.size());
+ const auto result2 =
+ ReferenceCells::Triangle.permute_by_combined_orientation(
+ view1, inverse_permutations[o]);
+
+ for (const auto &v : result2)
+ deallog << " " << v << std::endl;
+ }
+ }
+
+ deallog << "quadrilaterals" << std::endl;
+ {
+ const std::array<unsigned char, 8> inverse_permutations{
+ {0u, 1u, 2u, 7u, 4u, 5u, 6u, 3u}};
+
+ for (unsigned char o = 0; o < 8; ++o)
+ {
+ deallog << "o = " << int(o) << std::endl;
+ std::array<unsigned int, 4> vs{{0u, 1u, 2u, 3u}};
+
+ const auto result1 =
+ ReferenceCells::Quadrilateral.permute_by_combined_orientation(
+ make_const_array_view(vs), o);
+
+ ArrayView<const unsigned int> view1(result1.data(), result1.size());
+ const auto result2 =
+ ReferenceCells::Quadrilateral.permute_by_combined_orientation(
+ view1, inverse_permutations[o]);
+
+ for (const auto &v : result2)
+ deallog << " " << v << std::endl;
+ }
+ }
+
+ // Verify that the manual version created the same results.
+ deallog << "quadrilaterals (manual)" << std::endl;
+ {
+ for (unsigned char o = 0; o < 8; ++o)
+ {
+ deallog << "o = " << int(o) << std::endl;
+ std::array<unsigned int, 4> vs{{0u, 1u, 2u, 3u}};
+
+ const auto result1 =
+ ReferenceCells::Quadrilateral.permute_by_combined_orientation(
+ make_const_array_view(vs), o);
+
+ ArrayView<const unsigned int> view1(result1.data(), result1.size());
+
+ auto [orientation, rotation, flip] =
+ internal::split_face_orientation(o);
+ flip = orientation ? rotation ^ flip : flip;
+ const auto result2 =
+ ReferenceCells::Quadrilateral.permute_by_combined_orientation(
+ view1,
+ internal::combined_face_orientation(orientation, rotation, flip));
+
+ for (const auto &v : result2)
+ deallog << " " << v << std::endl;
+ }
+ }
+}