--- /dev/null
+New: Utilities::compress() and Utilities::decompress() allow to compress and decompress a string using gzip.
+<br>
+(Luca Heltai, Nicola Giuliani, 2020/01/17)
pack_integers(const std::array<std::uint64_t, dim> &index,
const int bits_per_dim);
+ /**
+ * If the library is configured with ZLIB, then this function compresses the
+ * input string and returns a non-zero terminated string containing the
+ * compressed input.
+ *
+ * If the library was not configured with ZLIB enabled, the returned string
+ * is identical to the input string.
+ *
+ * @param[in] input The string to compress
+ * @param[in] compression_level The compression level at which we compress
+ *
+ * @return A compressed version of the input string
+ *
+ * @authors Luca Heltai, Nicola Giuliani, 2020
+ */
+ std::string
+ compress(
+ const std::string &input,
+ const int compression_level = boost::iostreams::gzip::default_compression);
+
+ /**
+ * If the library is configured with ZLIB, then this function assumes that the
+ * input string has been compressed using the compress() function, and returns
+ * the original decompresses string.
+ *
+ * If the library was not configured with ZLIB enabled, the returned string
+ * is identical to the input string.
+ *
+ * @param[in] compressed_input A compressed string, as returned by the
+ * function compress()
+ *
+ * @return The original uncompressed string.
+ *
+ * @authors Luca Heltai, Nicola Giuliani, 2020
+ */
+ std::string
+ decompress(const std::string &compressed_input);
+
/**
* Convert a number @p value to a string, with as many digits as given to
* fill with leading zeros.
#include <deal.II/base/thread_local_storage.h>
#include <deal.II/base/utilities.h>
+#include <boost/iostreams/copy.hpp>
#include <boost/lexical_cast.hpp>
#include <boost/random.hpp>
+ std::string
+ compress(const std::string &input, const int compression_level)
+ {
+#ifdef DEAL_II_WITH_ZLIB
+ namespace bio = boost::iostreams;
+
+ std::stringstream compressed;
+ std::stringstream origin(input);
+
+ bio::filtering_streambuf<bio::input> out;
+ out.push(bio::gzip_compressor(bio::gzip_params(compression_level)));
+ out.push(origin);
+ bio::copy(out, compressed);
+
+ return compressed.str();
+#else
+ return input;
+#endif
+ }
+
+
+
+ std::string
+ decompress(const std::string &compressed_input)
+ {
+#ifdef DEAL_II_WITH_ZLIB
+ namespace bio = boost::iostreams;
+
+ std::stringstream compressed(compressed_input);
+ std::stringstream decompressed;
+
+ bio::filtering_streambuf<bio::input> out;
+ out.push(bio::gzip_decompressor());
+ out.push(compressed);
+ bio::copy(out, decompressed);
+
+ return decompressed.str();
+#else
+ return compressed_input;
+#endif
+ }
+
+
+
std::string
int_to_string(const unsigned int value, const unsigned int digits)
{
pack_integers<2>(const std::array<std::uint64_t, 2> &, const int);
template std::uint64_t
pack_integers<3>(const std::array<std::uint64_t, 3> &, const int);
+
+
} // namespace Utilities
DEAL_II_NAMESPACE_CLOSE
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2020 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::compress/decompress on a string
+
+#include <deal.II/base/utilities.h>
+
+#include "../tests.h"
+
+int
+main()
+{
+ initlog();
+ std::string input = "deal.II Rocks!";
+ auto compressed = Utilities::compress(input);
+ auto decompressed = Utilities::decompress(compressed);
+ deallog << (decompressed == input ? "OK" : "NOT OK") << std::endl;
+}
--- /dev/null
+
+DEAL::OK