--- /dev/null
+New:Added Utilities::pack and Utilities::unpack, to serialize and unserialize arbitrary
+objects that support boost::serialization/deserialization.
+<br>
+(Luca Heltai, 2017/11/18)
# endif
#endif
+DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
+#include <boost/archive/binary_oarchive.hpp>
+#include <boost/archive/binary_iarchive.hpp>
+#include <boost/serialization/vector.hpp>
+#include <boost/serialization/array.hpp>
+
+#ifdef DEAL_II_WITH_ZLIB
+# include <boost/iostreams/stream.hpp>
+# include <boost/iostreams/filtering_stream.hpp>
+# include <boost/iostreams/device/back_inserter.hpp>
+# include <boost/iostreams/filter/gzip.hpp>
+#endif
+
DEAL_II_NAMESPACE_OPEN
std::vector<unsigned long long int>
invert_permutation (const std::vector<unsigned long long int> &permutation);
+ /**
+ * Given an arbitrary object of type T, use boost::serialization utilities
+ * to pack the object into a vector of characters. The object can be unpacked
+ * using the Utilities::unpack function below.
+ *
+ * If the library has been compiled with ZLIB enabled, then the output buffer
+ * is compressed.
+ *
+ * @author Timo Heister, Wolfgang Bangerth, 2017.
+ */
+ template<typename T>
+ std::vector<char> pack(const T &object);
+
+ /**
+ * Given a vector of characters, obtained through a call to the function
+ * Utilities::pack, restore its content in an object of type T.
+ *
+ * This function uses boost::serialization utilities to unpack the object
+ * from a vector of characters, and it is the inverse of the function
+ * Utilities::pack
+ *
+ * @author Timo Heister, Wolfgang Bangerth, 2017.
+ */
+ template<typename T>
+ T unpack(const std::vector<char> &buffer);
+
/**
* A namespace for utility functions that probe system properties.
*
len = half;
}
}
+
+
+
+ template<typename T>
+ std::vector<char> pack(const T &object)
+ {
+ // set up a buffer and then use it as the target of a compressing
+ // stream into which we serialize the current object
+ std::vector<char> buffer;
+ {
+#ifdef DEAL_II_WITH_ZLIB
+ boost::iostreams::filtering_ostream out;
+ out.push(boost::iostreams::gzip_compressor
+ (boost::iostreams::gzip_params
+ (boost::iostreams::gzip::best_compression)));
+ out.push(boost::iostreams::back_inserter(buffer));
+
+ boost::archive::binary_oarchive archive(out);
+ archive << object;
+ out.flush();
+#else
+ std::ostringstream out;
+ boost::archive::binary_oarchive archive(out);
+ archive << object;
+ const std::string &s = out.str();
+ buffer.reserve(s.size());
+ buffer.assign(s.begin(), s.end());
+#endif
+ }
+ return buffer;
+ }
+
+
+
+ template<typename T>
+ T unpack(const std::vector<char> &buffer)
+ {
+ std::string decompressed_buffer;
+ T object;
+
+ // first decompress the buffer
+ {
+#ifdef DEAL_II_WITH_ZLIB
+ boost::iostreams::filtering_ostream decompressing_stream;
+ decompressing_stream.push(boost::iostreams::gzip_decompressor());
+ decompressing_stream.push(boost::iostreams::back_inserter(decompressed_buffer));
+ decompressing_stream.write (buffer.data(), buffer.size());
+#else
+ decompressed_buffer.assign (buffer.begin(), buffer.end());
+#endif
+ }
+
+ // then restore the object from the buffer
+ std::istringstream in(decompressed_buffer);
+ boost::archive::binary_iarchive archive(in);
+
+ archive >> object;
+ return object;
+ }
}
const unsigned int version);
BOOST_SERIALIZATION_SPLIT_MEMBER()
-
- /**
- * Pack the data that corresponds to this object into a buffer in
- * the form of a vector of chars and return it.
- */
- std::vector<char> pack_data () const;
-
- /**
- * Given a buffer in the form of an array of chars, unpack it and
- * restore the current object to the state that it was when
- * it was packed into said buffer by the pack_data() function.
- */
- void unpack_data (const std::vector<char> &buffer);
-
};
/**
- template <int dim, typename T>
- std::vector<char>
- CellDataTransferBuffer<dim,T>::pack_data () const
- {
- // set up a buffer and then use it as the target of a compressing
- // stream into which we serialize the current object
- std::vector<char> buffer;
- {
-#ifdef DEAL_II_WITH_ZLIB
- boost::iostreams::filtering_ostream out;
- out.push(boost::iostreams::gzip_compressor
- (boost::iostreams::gzip_params
- (boost::iostreams::gzip::best_compression)));
- out.push(boost::iostreams::back_inserter(buffer));
-
- boost::archive::binary_oarchive archive(out);
- archive << *this;
- out.flush();
-#else
- std::ostringstream out;
- boost::archive::binary_oarchive archive(out);
- archive << *this;
- const std::string &s = out.str();
- buffer.reserve(s.size());
- buffer.assign(s.begin(), s.end());
-#endif
- }
-
- return buffer;
- }
-
-
-
- template <int dim, typename T>
- void
- CellDataTransferBuffer<dim,T>::unpack_data (const std::vector<char> &buffer)
- {
- std::string decompressed_buffer;
-
- // first decompress the buffer
- {
-#ifdef DEAL_II_WITH_ZLIB
- boost::iostreams::filtering_ostream decompressing_stream;
- decompressing_stream.push(boost::iostreams::gzip_decompressor());
- decompressing_stream.push(boost::iostreams::back_inserter(decompressed_buffer));
- decompressing_stream.write (buffer.data(), buffer.size());
-#else
- decompressed_buffer.assign (buffer.begin(), buffer.end());
-#endif
- }
-
- // then restore the object from the buffer
- std::istringstream in(decompressed_buffer);
- boost::archive::binary_iarchive archive(in);
-
- archive >> *this;
- }
-
-
-
template <typename DataType, typename MeshType>
void
exchange_cell_data_to_ghosts (const MeshType &mesh,
// pack all the data into the buffer for this recipient and send it.
// keep data around till we can make sure that the packet has been
// received
- sendbuffers[idx] = data.pack_data ();
+ sendbuffers[idx] = Utilities::pack(data);
const int ierr = MPI_Isend(sendbuffers[idx].data(), sendbuffers[idx].size(),
MPI_BYTE, *it,
786, tria->get_communicator(), &requests[idx]);
tria->get_communicator(), &status);
AssertThrowMPI(ierr);
- CellDataTransferBuffer<dim, DataType> cellinfo;
- cellinfo.unpack_data(receive);
+ auto cellinfo = Utilities::unpack<CellDataTransferBuffer<dim, DataType> >(receive);
DataType *data = cellinfo.data.data();
for (unsigned int c=0; c<cellinfo.cell_ids.size(); ++c, ++data)
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2017 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 at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+// test Utilities::pack/unpack on some types.
+
+#include "../tests.h"
+
+#include <deal.II/base/utilities.h>
+#include <deal.II/base/point.h>
+
+template <int dim>
+void test(const unsigned int &size)
+{
+ std::vector<Point<dim> > points(size);
+
+ for (auto &p : points)
+ for (unsigned int i=0; i<dim; ++i)
+ p[i] = Testing::rand()/double(RAND_MAX);
+
+ auto buffer = Utilities::pack(points);
+
+ auto unpacked = Utilities::unpack<std::vector<Point<dim> > >(buffer);
+
+ unsigned int i=0;
+ bool ok = true;
+ for (const auto &p : points)
+ if (p.distance(unpacked[i++]) > 1e-12)
+ {
+ deallog << "NOT OK: "
+ << p << " != " << unpacked[i-1] << std::endl;
+ ok = false;
+ }
+
+ if (ok)
+ deallog << "OK!" << std::endl;
+}
+
+int main()
+{
+ initlog();
+
+ test<1>(10);
+ test<2>(10);
+ test<3>(10);
+}
--- /dev/null
+
+DEAL::OK!
+DEAL::OK!
+DEAL::OK!