+ /**
+ * Type of the 9-bit representation of the refinement configuration.
+ */
+ using compressed_constraint_kind = std::uint8_t;
+
+
+
+ /**
+ * Value of compressed_constraint_kind for the unconstrained case.
+ */
+ constexpr compressed_constraint_kind
+ unconstrained_compressed_constraint_kind = 0;
+
+
+
/**
* Check if the combinations of the bits in @p kind_in are valid.
*/
inline bool
- check(const ConstraintKinds &kind_in, const unsigned int dim)
+ check(const ConstraintKinds kind_in, const unsigned int dim)
{
const std::uint16_t kind = static_cast<std::uint16_t>(kind_in);
const std::uint16_t subcell = (kind >> 0) & 7;
+ /**
+ * Compress a 9-bit representation of the refinement configuration
+ * to a 8-bit representation.
+ */
+ inline compressed_constraint_kind
+ compress(const ConstraintKinds kind_in, const unsigned int dim)
+ {
+ Assert(check(kind_in, dim), ExcInternalError());
+
+ if (dim == 2)
+ return static_cast<compressed_constraint_kind>(kind_in);
+
+ if (kind_in == ConstraintKinds::unconstrained)
+ return unconstrained_compressed_constraint_kind;
+
+ const std::uint16_t kind = static_cast<std::uint16_t>(kind_in);
+ const std::uint16_t subcell = (kind >> 0) & 7;
+ const std::uint16_t face = (kind >> 3) & 7;
+ const std::uint16_t edge = (kind >> 6) & 7;
+
+ return subcell + ((face > 0) << 3) + ((edge > 0) << 4) +
+ (std::max(face, edge) << 5);
+ }
+
+
+
+ /**
+ * Decompress a 8-bit representation of the refinement configuration
+ * to a 9-bit representation.
+ */
+ inline ConstraintKinds
+ decompress(const compressed_constraint_kind kind_in, const unsigned int dim)
+ {
+ if (dim == 2)
+ return static_cast<ConstraintKinds>(kind_in);
+
+ if (kind_in == unconstrained_compressed_constraint_kind)
+ return ConstraintKinds::unconstrained;
+
+ const std::uint16_t subcell = (kind_in >> 0) & 7;
+ const std::uint16_t flag_0 = (kind_in >> 3) & 3;
+ const std::uint16_t flag_1 = (kind_in >> 5) & 7;
+
+ const auto result = static_cast<ConstraintKinds>(
+ subcell + (((flag_0 & 0b01) ? flag_1 : 0) << 3) +
+ (((flag_0 & 0b10) ? flag_1 : 0) << 6));
+
+ Assert(check(result, dim), ExcInternalError());
+
+ return result;
+ }
+
+
+
/**
* Return the memory consumption in bytes of this enum class.
*/
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2022 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+
+// Compress/decompress ConstraintKinds.
+
+#include <deal.II/dofs/dof_tools.h>
+
+#include <deal.II/fe/fe_q.h>
+
+#include <deal.II/grid/grid_generator.h>
+
+#include <deal.II/matrix_free/evaluation_kernels_hanging_nodes.h>
+#include <deal.II/matrix_free/fe_evaluation.h>
+#include <deal.II/matrix_free/matrix_free.h>
+
+#include <deal.II/numerics/vector_tools.h>
+
+#include "../tests.h"
+
+int
+main()
+{
+ initlog();
+
+ using namespace dealii::internal::MatrixFreeFunctions;
+
+ const unsigned int dim = 3;
+
+ const auto check = [&](const auto subcell, const auto face, const auto edge) {
+ const auto kind =
+ static_cast<ConstraintKinds>(subcell + (face << 3) + (edge << 6));
+
+ Assert(kind == decompress(compress(kind, dim), dim), ExcInternalError());
+ };
+
+ check(0, 0, 0);
+
+ for (unsigned int subcell = 0; subcell < 8; ++subcell)
+ {
+ for (unsigned int face = 1; face < 8; ++face)
+ check(subcell, face, 0);
+
+ for (unsigned int edge = 1; edge < 8; ++edge)
+ check(subcell, 0, edge);
+
+ for (unsigned int d = 0; d < dim; ++d)
+ check(subcell, 1 << d, 1 << d);
+ }
+
+ deallog << "OK!" << std::endl;
+}