]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Added support for tuple serialization in Utilities. 5500/head
authorLuca Heltai <luca.heltai@sissa.it>
Sat, 18 Nov 2017 19:18:27 +0000 (20:18 +0100)
committerLuca Heltai <luca.heltai@sissa.it>
Sat, 18 Nov 2017 19:18:27 +0000 (20:18 +0100)
include/deal.II/base/utilities.h
tests/base/utilities_pack_unpack_02.cc [new file with mode: 0644]
tests/base/utilities_pack_unpack_02.output [new file with mode: 0644]

index 4c44c34643e75673f053514b3cd81658185cf8da..9aef3ed0feac8c71ddfe0c12d4f9d72f5dadd29b 100644 (file)
@@ -23,6 +23,7 @@
 #include <utility>
 #include <functional>
 #include <string>
+#include <tuple>
 
 #ifdef DEAL_II_WITH_TRILINOS
 #  include <Epetra_Comm.h>
@@ -833,4 +834,43 @@ namespace Utilities
 
 DEAL_II_NAMESPACE_CLOSE
 
+#ifndef DOXYGEN
+namespace boost
+{
+  namespace serialization
+  {
+
+    // Provides boost and c++11 with a way to serialize tuples and pairs automatically
+    template<int N>
+    struct Serialize
+    {
+      template<class Archive, typename... Args>
+      static void serialize(Archive &ar, std::tuple<Args...> &t, const unsigned int version)
+      {
+        ar &std::get<N-1>(t);
+        Serialize<N-1>::serialize(ar, t, version);
+      }
+    };
+
+    template<>
+    struct Serialize<0>
+    {
+      template<class Archive, typename... Args>
+      static void serialize(Archive &ar, std::tuple<Args...> &t, const unsigned int version)
+      {
+        (void) ar;
+        (void) t;
+        (void) version;
+      }
+    };
+
+    template<class Archive, typename... Args>
+    void serialize(Archive &ar, std::tuple<Args...> &t, const unsigned int version)
+    {
+      Serialize<sizeof...(Args)>::serialize(ar, t, version);
+    }
+  }
+}
+#endif
+
 #endif
diff --git a/tests/base/utilities_pack_unpack_02.cc b/tests/base/utilities_pack_unpack_02.cc
new file mode 100644 (file)
index 0000000..194354e
--- /dev/null
@@ -0,0 +1,68 @@
+// ---------------------------------------------------------------------
+//
+// 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>
+#include <boost/serialization/utility.hpp>
+
+#include <tuple>
+
+template <int dim>
+void test(const unsigned int &size)
+{
+  std::vector<Point<dim> > points(size);
+  auto a_pair = std::make_pair(1, 3.14);
+
+  for (auto &p : points)
+    for (unsigned int i=0; i<dim; ++i)
+      p[i] = Testing::rand()/double(RAND_MAX);
+
+  auto a_tuple = std::make_tuple(a_pair, points);
+
+  auto buffer = Utilities::pack(a_tuple);
+
+  auto tuple_unpacked = Utilities::unpack<decltype(a_tuple)>(buffer);
+
+  auto pair_unpacked = std::get<0>(tuple_unpacked);
+  auto points_unpacked = std::get<1>(tuple_unpacked);
+
+  unsigned int i=0;
+  bool ok = (pair_unpacked == a_pair);
+
+  for (const auto &p : points)
+    if (p.distance(points_unpacked[i++]) > 1e-12)
+      {
+        deallog << "NOT OK: "
+                << p << " != " << points_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);
+}
diff --git a/tests/base/utilities_pack_unpack_02.output b/tests/base/utilities_pack_unpack_02.output
new file mode 100644 (file)
index 0000000..7ca52c4
--- /dev/null
@@ -0,0 +1,4 @@
+
+DEAL::OK!
+DEAL::OK!
+DEAL::OK!

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.