]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Added serialization benchmark. 11364/head
authorMarc Fehling <marc.fehling@gmx.net>
Tue, 15 Dec 2020 04:00:24 +0000 (21:00 -0700)
committerMarc Fehling <mafehling.git@gmail.com>
Sat, 19 Dec 2020 00:55:16 +0000 (17:55 -0700)
tests/grid/cell_id_08.cc [new file with mode: 0644]
tests/grid/cell_id_08.output [new file with mode: 0644]

diff --git a/tests/grid/cell_id_08.cc b/tests/grid/cell_id_08.cc
new file mode 100644 (file)
index 0000000..47031eb
--- /dev/null
@@ -0,0 +1,189 @@
+// ---------------------------------------------------------------------
+//
+// 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.
+//
+// ---------------------------------------------------------------------
+
+
+// benchmark for creation and serialization of CellId objects
+
+#include <deal.II/base/timer.h>
+#include <deal.II/base/utilities.h>
+
+#include <deal.II/grid/cell_id.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria.h>
+
+#include <boost/archive/binary_iarchive.hpp>
+#include <boost/archive/binary_oarchive.hpp>
+#include <boost/iostreams/device/array.hpp>
+#include <boost/iostreams/device/back_inserter.hpp>
+#include <boost/iostreams/filtering_streambuf.hpp>
+
+#include <iostream>
+#include <vector>
+
+#include "../tests.h"
+
+
+template <typename Container>
+std::size_t
+size_in_bytes(const Container &container)
+{
+  return sizeof(Container) +
+         (sizeof(typename Container::value_type) * container.size());
+}
+
+
+
+template <int dim>
+void
+benchmark(const unsigned int n_global_refinements)
+{
+  // be generous with the estimate
+  const std::size_t size_estimate_per_cellid_in_bytes = 96;
+
+  TimerOutput computing_timer(std::cout,
+                              TimerOutput::summary,
+                              TimerOutput::wall_times);
+
+  Triangulation<dim> tria;
+  GridGenerator::hyper_cube(tria);
+  tria.refine_global(n_global_refinements);
+
+  const unsigned int n_cells = tria.n_active_cells();
+
+
+  //
+  // create
+  //
+  std::vector<CellId> cell_ids;
+  cell_ids.reserve(n_cells);
+  {
+    TimerOutput::Scope t(computing_timer, "cell->id");
+
+    for (const auto &cell : tria.active_cell_iterators())
+      cell_ids.push_back(cell->id());
+  }
+  deallog << "Creation successful." << std::endl
+          << "  Number of CellId objects: " << cell_ids.size() << std::endl;
+
+
+  //
+  // boost binary archive
+  //
+  {
+    std::string buffer;
+    buffer.reserve(n_cells * size_estimate_per_cellid_in_bytes);
+    std::vector<CellId> deserialized;
+    deserialized.reserve(n_cells);
+    {
+      TimerOutput::Scope t(computing_timer, "boost binary_archive save");
+
+      boost::iostreams::filtering_ostreambuf fosb;
+      fosb.push(boost::iostreams::back_inserter(buffer));
+      boost::archive::binary_oarchive oa(fosb);
+      oa &                            cell_ids;
+    }
+    {
+      TimerOutput::Scope t(computing_timer, "boost binary_archive load");
+
+      boost::iostreams::filtering_istreambuf fisb;
+      fisb.push(boost::iostreams::array_source(buffer.data(), buffer.size()));
+      boost::archive::binary_iarchive ia(fisb);
+      ia &                            deserialized;
+    }
+    AssertThrow(cell_ids == deserialized,
+                ExcMessage("Serialization with boost binary archives failed."));
+    deallog << "Serialization with boost binary archives successful."
+            << std::endl
+            << "  Buffer size in bytes: " << size_in_bytes(buffer) << std::endl;
+  }
+
+
+  //
+  // Utilities uncompressed
+  //
+  {
+    std::vector<char> buffer;
+    buffer.reserve(n_cells * size_estimate_per_cellid_in_bytes);
+    std::vector<std::size_t> sizes;
+    sizes.reserve(n_cells);
+    std::vector<CellId> deserialized;
+    deserialized.reserve(n_cells);
+    {
+      TimerOutput::Scope t(computing_timer, "Utilities::pack uncompressed");
+
+      for (const auto &cell_id : cell_ids)
+        sizes.push_back(
+          Utilities::pack(cell_id, buffer, /*allow_compression*/ false));
+    }
+    {
+      TimerOutput::Scope t(computing_timer, "Utilities::unpack uncompressed");
+
+      std::vector<char>::const_iterator cbegin = buffer.cbegin(), cend;
+      for (const auto &size : sizes)
+        {
+          cend = cbegin + size;
+          deserialized.push_back(Utilities::unpack<CellId>(
+            cbegin, cend, /*allow_compression*/ false));
+          cbegin = cend;
+        }
+      Assert(cend == buffer.cend(), ExcInternalError());
+    }
+    AssertThrow(cell_ids == deserialized,
+                ExcMessage("Serialization with Utilities failed."));
+    deallog << "Serialization with Utilities successful." << std::endl
+            << "  Buffer size in bytes: " << size_in_bytes(buffer) << std::endl;
+  }
+
+
+  //
+  // binary representation
+  //
+  {
+    std::vector<CellId::binary_type> buffer;
+    buffer.reserve(n_cells);
+    std::vector<CellId> deserialized;
+    deserialized.reserve(n_cells);
+    {
+      TimerOutput::Scope t(computing_timer, "CellId::to_binary");
+
+      for (const auto &cell_id : cell_ids)
+        buffer.push_back(cell_id.to_binary<dim>());
+    }
+    {
+      TimerOutput::Scope t(computing_timer, "CellId from binary");
+
+      for (const auto &binary : buffer)
+        deserialized.emplace_back(binary);
+    }
+    AssertThrow(cell_ids == deserialized,
+                ExcMessage("Serialization with binary representation failed."));
+    deallog << "Serialization with binary representation successful."
+            << std::endl
+            << "  Buffer size in bytes: " << size_in_bytes(buffer) << std::endl;
+  }
+
+  deallog << "OK" << std::endl;
+}
+
+
+
+int
+main()
+{
+  initlog();
+
+  benchmark<2>(4);
+  // benchmark<3>(8);
+}
diff --git a/tests/grid/cell_id_08.output b/tests/grid/cell_id_08.output
new file mode 100644 (file)
index 0000000..22ff31e
--- /dev/null
@@ -0,0 +1,10 @@
+
+DEAL::Creation successful.
+DEAL::  Number of CellId objects: 256
+DEAL::Serialization with boost binary archives successful.
+DEAL::  Buffer size in bytes: 11875
+DEAL::Serialization with Utilities successful.
+DEAL::  Buffer size in bytes: 10264
+DEAL::Serialization with binary representation successful.
+DEAL::  Buffer size in bytes: 4120
+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.