const std::vector<std::array<std::uint64_t, dim>> &points,
const int bits_per_dim = 64);
+ /**
+ * Interleave first @p bits_per_dim bits from each element of @p index
+ * (starting from last) into a single unsigned integer. The function is useful
+ * in debugging and visualization of indices returned by
+ * inverse_Hilbert_space_filling_curve(). Clearly, the following should hold
+ * <code>bits\_per\_dim * dim <= 64</code>.
+ *
+ * @note There is no need to use this function in order to compare indices
+ * returned by inverse_Hilbert_space_filling_curve(), that can be easily done
+ * via <code>std::lexicographical_compare()</code>.
+ */
+ template <int dim>
+ std::uint64_t
+ interleave(const std::array<std::uint64_t, dim> &index,
+ const int bits_per_dim);
+
/**
* Convert a number @p value to a string, with as many digits as given to
* fill with leading zeros.
+ template <int dim>
+ std::uint64_t
+ interleave(const std::array<std::uint64_t, dim> &index,
+ const int bits_per_dim)
+ {
+ using Integer = std::uint64_t;
+
+ AssertIndexRange(bits_per_dim * dim, 65);
+ Assert(bits_per_dim > 0, ExcMessage("bits_per_dim should be positive"));
+
+ const Integer mask = (1 << bits_per_dim) - 1;
+
+ Integer res = 0;
+ for (unsigned int i = 0; i < dim; ++i)
+ {
+ // take @bits_per_dim from each integer and shift them
+ const Integer v = (mask & index[dim - 1 - i]) << (bits_per_dim * i);
+ res |= v;
+ }
+ return res;
+ }
+
+
+
std::string
int_to_string(const unsigned int value, const unsigned int digits)
{
const std::vector<std::array<std::uint64_t, 3>> &,
const int);
+ template std::uint64_t
+ interleave<1>(const std::array<std::uint64_t, 1> &, const int);
+ template std::uint64_t
+ interleave<2>(const std::array<std::uint64_t, 2> &, const int);
+ template std::uint64_t
+ interleave<3>(const std::array<std::uint64_t, 3> &, const int);
} // namespace Utilities
DEAL_II_NAMESPACE_CLOSE
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2018 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.
+//
+// ---------------------------------------------------------------------
+
+
+// test Utilities::interleave() for indices from
+// Utilities::inverse_Hilbert_space_filling_curve in 2D on the following points
+// with bits=2:
+
+// 3 | 5---6 9---10
+// | | | | |
+// 2 | 4 7---8 11
+// | | |
+// 1 | 3---2 13--12
+// | | |
+// 0 | 0---1 14--15
+// |
+// ------------------
+// 0 1 2 3
+
+#include <deal.II/base/utilities.h>
+
+#include "../tests.h"
+
+
+void
+test()
+{
+ const std::vector<Point<2>> points = {Point<2>(0, 0),
+ Point<2>(1, 0),
+ Point<2>(1, 1),
+ Point<2>(0, 1),
+ Point<2>(0, 2),
+ Point<2>(0, 3),
+ Point<2>(1, 3),
+ Point<2>(1, 2),
+ Point<2>(2, 2),
+ Point<2>(2, 3),
+ Point<2>(3, 3),
+ Point<2>(3, 2),
+ Point<2>(3, 1),
+ Point<2>(2, 1),
+ Point<2>(2, 0),
+ Point<2>(3, 0)};
+
+ const int bit_depth = 2;
+ const auto res =
+ Utilities::inverse_Hilbert_space_filling_curve(points, bit_depth);
+
+ for (const auto &p : res)
+ deallog << p[0] << " " << p[1] << " "
+ << Utilities::interleave<2>(p, bit_depth) << std::endl;
+}
+
+int
+main()
+{
+ initlog();
+
+ test();
+}
--- /dev/null
+
+DEAL::0 0 0
+DEAL::0 1 1
+DEAL::0 2 2
+DEAL::0 3 3
+DEAL::1 0 4
+DEAL::1 1 5
+DEAL::1 2 6
+DEAL::1 3 7
+DEAL::2 0 8
+DEAL::2 1 9
+DEAL::2 2 10
+DEAL::2 3 11
+DEAL::3 0 12
+DEAL::3 1 13
+DEAL::3 2 14
+DEAL::3 3 15