]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add boost serialization for ConstructionData 8748/head
authorPeter Munch <peterrmuench@gmail.com>
Fri, 13 Sep 2019 08:04:29 +0000 (10:04 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Thu, 19 Sep 2019 03:49:59 +0000 (05:49 +0200)
doc/news/changes/minor/20190913PeterMunch-1 [new file with mode: 0644]
include/deal.II/distributed/fully_distributed_tria.h
tests/fullydistributed_grids/tests.h
tests/serialization/parallel_fullydistributed_construction_data_1.cc [new file with mode: 0644]
tests/serialization/parallel_fullydistributed_construction_data_1.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20190913PeterMunch-1 b/doc/news/changes/minor/20190913PeterMunch-1
new file mode 100644 (file)
index 0000000..486bf7d
--- /dev/null
@@ -0,0 +1,4 @@
+New: Add Boost serialization function to dealii::fullydistributed::CellData and
+dealii::fullydistributed::ConstructionData.
+<br>
+(Peter Munch, 2019/09/13)
index eef0b55d2d5bc36c838552b15039734f750823c5..b035b298bcbac83cab5264722775301da1946669 100644 (file)
@@ -85,6 +85,13 @@ namespace parallel
        */
       CellData() = default;
 
+      /**
+       * Boost serialization function
+       */
+      template <class Archive>
+      void
+      serialize(Archive &ar, const unsigned int /*version*/);
+
       /**
        * Unique CellID of the cell.
        */
@@ -118,7 +125,8 @@ namespace parallel
        *
        * @note Only used for 3D.
        */
-      std::array<types::material_id, GeometryInfo<dim>::quads_per_cell>
+      std::array<types::material_id,
+                 dim == 1 ? 1 : GeometryInfo<3>::quads_per_cell>
         manifold_quad_ids;
 
       /**
@@ -139,6 +147,13 @@ namespace parallel
     template <int dim, int spacedim>
     struct ConstructionData
     {
+      /**
+       * Boost serialization function
+       */
+      template <class Archive>
+      void
+      serialize(Archive &ar, const unsigned int /*version*/);
+
       /**
        * Cells of the locally-relevant coarse-grid triangulation.
        */
@@ -385,6 +400,42 @@ namespace parallel
         currently_processing_prepare_coarsening_and_refinement_for_internal_usage;
     };
 
+
+#ifndef DOXYGEN
+
+    template <int dim>
+    template <class Archive>
+    void
+    CellData<dim>::serialize(Archive &ar, const unsigned int /*version*/)
+    {
+      ar &id;
+      ar &subdomain_id;
+      ar &level_subdomain_id;
+      ar &manifold_id;
+      if (dim >= 2)
+        ar &manifold_line_ids;
+      if (dim >= 3)
+        ar &manifold_quad_ids;
+      ar &boundary_ids;
+    }
+
+
+    template <int dim, int spacedim>
+    template <class Archive>
+    void
+    ConstructionData<dim, spacedim>::serialize(Archive &ar,
+                                               const unsigned int /*version*/)
+    {
+      ar &coarse_cells;
+      ar &coarse_cell_vertices;
+      ar &coarse_cell_index_to_coarse_cell_id;
+      ar &cell_infos;
+      ar &settings;
+    }
+
+
+#endif
+
   } // namespace fullydistributed
 } // namespace parallel
 
index 67cdbd71bb696520d2e901a8848aa314ad792448..4ba696fc90282845c6de5599e156d90fea45fa5b 100644 (file)
@@ -72,4 +72,55 @@ print_statistics(const DoFHandler<dim, spacedim> &dof_handler,
   deallog << std::endl;
 }
 
+namespace dealii
+{
+  namespace parallel
+  {
+    namespace fullydistributed
+    {
+      template <int dim>
+      bool
+      operator==(const CellData<dim> &t1, const CellData<dim> &t2)
+      {
+        if (t1.id != t2.id)
+          return false;
+        if (t1.subdomain_id != t2.subdomain_id)
+          return false;
+        if (t1.level_subdomain_id != t2.level_subdomain_id)
+          return false;
+        if (t1.manifold_id != t2.manifold_id)
+          return false;
+        if (dim >= 2 && t1.manifold_line_ids != t2.manifold_line_ids)
+          return false;
+        if (dim >= 3 && t1.manifold_quad_ids != t2.manifold_quad_ids)
+          return false;
+        if (t1.boundary_ids != t2.boundary_ids)
+          return false;
+
+        return true;
+      }
+
+      template <int dim, int spacedim>
+      bool
+      operator==(const ConstructionData<dim, spacedim> &t1,
+                 const ConstructionData<dim, spacedim> &t2)
+      {
+        if (t1.coarse_cells != t2.coarse_cells)
+          return false;
+        if (t1.coarse_cell_vertices != t2.coarse_cell_vertices)
+          return false;
+        if (t1.coarse_cell_index_to_coarse_cell_id !=
+            t2.coarse_cell_index_to_coarse_cell_id)
+          return false;
+        if (t1.cell_infos != t2.cell_infos)
+          return false;
+        if (t1.settings != t2.settings)
+          return false;
+
+        return true;
+      }
+    } // namespace fullydistributed
+  }   // namespace parallel
+} // namespace dealii
+
 #endif
diff --git a/tests/serialization/parallel_fullydistributed_construction_data_1.cc b/tests/serialization/parallel_fullydistributed_construction_data_1.cc
new file mode 100644 (file)
index 0000000..67d52e1
--- /dev/null
@@ -0,0 +1,79 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2019 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.
+//
+// ---------------------------------------------------------------------
+
+
+// check serialization for ConstructionData<dim, spacedim>
+
+#include <deal.II/distributed/fully_distributed_tria.h>
+#include <deal.II/distributed/fully_distributed_tria_util.h>
+#include <deal.II/distributed/tria.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_tools.h>
+#include <deal.II/grid/tria.h>
+
+#include <boost/serialization/vector.hpp>
+
+#include "../fullydistributed_grids/tests.h"
+#include "serialization.h"
+
+using namespace dealii;
+
+template <int dim>
+void
+test(MPI_Comm comm)
+{
+  // create serial triangulation
+  Triangulation<dim> basetria;
+  GridGenerator::hyper_cube(basetria);
+  basetria.refine_global(1);
+
+  GridTools::partition_triangulation_zorder(
+    Utilities::MPI::n_mpi_processes(comm), basetria);
+
+  auto t1 = parallel::fullydistributed::Utilities::
+    create_construction_data_from_triangulation(basetria, comm);
+
+  // compare equal ConstructionDatas
+  auto t2 = t1;
+  verify(t1, t2);
+
+  basetria.refine_global(1);
+
+  GridTools::partition_triangulation_zorder(
+    Utilities::MPI::n_mpi_processes(comm), basetria);
+
+  auto t3 = parallel::fullydistributed::Utilities::
+    create_construction_data_from_triangulation(basetria, comm);
+
+  // compare different ConstructionDatas
+  verify(t1, t3);
+}
+
+
+int
+main(int argc, char *argv[])
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+  MPILogInitAll                    all;
+  deallog << std::setprecision(3);
+  const MPI_Comm comm = MPI_COMM_WORLD;
+
+  test<1>(comm);
+  test<2>(comm);
+  test<3>(comm);
+
+  deallog << "OK" << std::endl;
+}
diff --git a/tests/serialization/parallel_fullydistributed_construction_data_1.output b/tests/serialization/parallel_fullydistributed_construction_data_1.output
new file mode 100644 (file)
index 0000000..a636cf5
--- /dev/null
@@ -0,0 +1,14 @@
+
+DEAL:0::0 0 0 0 2 0 0 0 2 0 2 0 0 4294967295 2 2 1 0 0 4294967295 0 0 3 0 0 0 0 0 1 0 0 0.00000000000000000e+00 1 1.00000000000000000e+00 1 5.00000000000000000e-01 2 0 0 1 0 0 1 0 0 0 2 0 0 0 0 0 4 0 1 0 0 0 4294967295 4294967295 0 0 1 0 0 0 0 0 4 1 1 0 0 0 4294967295 4294967295 1 0 1 1 0
+
+DEAL:0::0 0 0 0 2 0 0 0 2 0 2 0 0 4294967295 2 2 1 0 0 4294967295 0 0 3 0 0 0 0 0 1 0 0 0.00000000000000000e+00 1 1.00000000000000000e+00 1 5.00000000000000000e-01 2 0 0 1 0 0 1 0 0 0 2 0 0 0 0 0 4 0 1 0 0 0 4294967295 4294967295 0 0 1 0 0 0 0 0 4 1 1 0 0 0 4294967295 4294967295 1 0 1 1 0
+
+DEAL:0::0 0 0 0 4 0 0 0 4 0 4 5 8 0 0 4294967295 4 4 1 8 6 0 0 4294967295 4 5 8 2 7 0 0 4294967295 4 8 6 7 3 0 0 4294967295 0 0 9 0 0 0 0 0 2 0 0 0.00000000000000000e+00 0.00000000000000000e+00 2 1.00000000000000000e+00 0.00000000000000000e+00 2 0.00000000000000000e+00 1.00000000000000000e+00 2 1.00000000000000000e+00 1.00000000000000000e+00 2 5.00000000000000000e-01 0.00000000000000000e+00 2 0.00000000000000000e+00 5.00000000000000000e-01 2 1.00000000000000000e+00 5.00000000000000000e-01 2 5.00000000000000000e-01 1.00000000000000000e+00 2 5.00000000000000000e-01 5.00000000000000000e-01 4 0 0 1 2 3 0 0 1 0 0 0 4 0 0 0 0 0 4 0 2 0 0 0 4294967295 4294967295 4 4294967295 4294967295 4294967295 4294967295 0 0 2 0 0 0 0 0 2 0 4 1 2 0 0 0 4294967295 4294967295 4 4294967295 4294967295 4294967295 4294967295 2 0 1 0 2 0 4 2 2 0 0 0 4294967295 4294967295 4 4294967295 4294967295 4294967295 4294967295 2 0 0 0 3 0 4 3 2 0 0 0 4294967295 4294967295 4 4294967295 4294967295 4294967295 4294967295 2 0 1 0 3 0 0
+
+DEAL:0::0 0 0 0 4 0 0 0 4 0 4 5 8 0 0 4294967295 4 4 1 8 6 0 0 4294967295 4 5 8 2 7 0 0 4294967295 4 8 6 7 3 0 0 4294967295 0 0 9 0 0 0 0 0 2 0 0 0.00000000000000000e+00 0.00000000000000000e+00 2 1.00000000000000000e+00 0.00000000000000000e+00 2 0.00000000000000000e+00 1.00000000000000000e+00 2 1.00000000000000000e+00 1.00000000000000000e+00 2 5.00000000000000000e-01 0.00000000000000000e+00 2 0.00000000000000000e+00 5.00000000000000000e-01 2 1.00000000000000000e+00 5.00000000000000000e-01 2 5.00000000000000000e-01 1.00000000000000000e+00 2 5.00000000000000000e-01 5.00000000000000000e-01 4 0 0 1 2 3 0 0 1 0 0 0 4 0 0 0 0 0 4 0 2 0 0 0 4294967295 4294967295 4 4294967295 4294967295 4294967295 4294967295 0 0 2 0 0 0 0 0 2 0 4 1 2 0 0 0 4294967295 4294967295 4 4294967295 4294967295 4294967295 4294967295 2 0 1 0 2 0 4 2 2 0 0 0 4294967295 4294967295 4 4294967295 4294967295 4294967295 4294967295 2 0 0 0 3 0 4 3 2 0 0 0 4294967295 4294967295 4 4294967295 4294967295 4294967295 4294967295 2 0 1 0 3 0 0
+
+DEAL:0::0 0 0 0 8 0 0 0 8 0 8 9 21 10 20 22 26 0 0 4294967295 8 8 1 21 11 20 12 26 23 0 0 4294967295 8 9 21 2 13 22 26 14 24 0 0 4294967295 8 21 11 13 3 26 23 24 15 0 0 4294967295 8 10 20 22 26 4 16 17 25 0 0 4294967295 8 20 12 26 23 16 5 25 18 0 0 4294967295 8 22 26 14 24 17 25 6 19 0 0 4294967295 8 26 23 24 15 25 18 19 7 0 0 4294967295 0 0 27 0 0 0 0 0 3 0 0 0.00000000000000000e+00 0.00000000000000000e+00 0.00000000000000000e+00 3 1.00000000000000000e+00 0.00000000000000000e+00 0.00000000000000000e+00 3 0.00000000000000000e+00 1.00000000000000000e+00 0.00000000000000000e+00 3 1.00000000000000000e+00 1.00000000000000000e+00 0.00000000000000000e+00 3 0.00000000000000000e+00 0.00000000000000000e+00 1.00000000000000000e+00 3 1.00000000000000000e+00 0.00000000000000000e+00 1.00000000000000000e+00 3 0.00000000000000000e+00 1.00000000000000000e+00 1.00000000000000000e+00 3 1.00000000000000000e+00 1.00000000000000000e+00 1.00000000000000000e+00 3 5.00000000000000000e-01 0.00000000000000000e+00 0.00000000000000000e+00 3 0.00000000000000000e+00 5.00000000000000000e-01 0.00000000000000000e+00 3 0.00000000000000000e+00 0.00000000000000000e+00 5.00000000000000000e-01 3 1.00000000000000000e+00 5.00000000000000000e-01 0.00000000000000000e+00 3 1.00000000000000000e+00 0.00000000000000000e+00 5.00000000000000000e-01 3 5.00000000000000000e-01 1.00000000000000000e+00 0.00000000000000000e+00 3 0.00000000000000000e+00 1.00000000000000000e+00 5.00000000000000000e-01 3 1.00000000000000000e+00 1.00000000000000000e+00 5.00000000000000000e-01 3 5.00000000000000000e-01 0.00000000000000000e+00 1.00000000000000000e+00 3 0.00000000000000000e+00 5.00000000000000000e-01 1.00000000000000000e+00 3 1.00000000000000000e+00 5.00000000000000000e-01 1.00000000000000000e+00 3 5.00000000000000000e-01 1.00000000000000000e+00 1.00000000000000000e+00 3 5.00000000000000000e-01 0.00000000000000000e+00 5.00000000000000000e-01 3 5.00000000000000000e-01 5.00000000000000000e-01 0.00000000000000000e+00 3 0.00000000000000000e+00 5.00000000000000000e-01 5.00000000000000000e-01 3 1.00000000000000000e+00 5.00000000000000000e-01 5.00000000000000000e-01 3 5.00000000000000000e-01 1.00000000000000000e+00 5.00000000000000000e-01 3 5.00000000000000000e-01 5.00000000000000000e-01 1.00000000000000000e+00 3 5.00000000000000000e-01 5.00000000000000000e-01 5.00000000000000000e-01 8 0 0 1 2 3 4 5 6 7 0 0 1 0 0 0 8 0 0 0 0 0 4 0 3 0 0 0 4294967295 4294967295 0 0 12 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 0 0 6 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 0 0 3 0 0 0 0 0 2 0 4 0 4 1 3 0 0 0 4294967295 4294967295 12 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 6 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 3 0 1 0 2 0 4 0 4 2 3 0 0 0 4294967295 4294967295 12 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 6 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 3 0 0 0 3 0 4 0 4 3 3 0 0 0 4294967295 4294967295 12 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 6 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 3 0 1 0 3 0 4 0 4 4 3 0 0 0 4294967295 4294967295 12 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 6 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 3 0 0 0 2 0 5 0 4 5 3 0 0 0 4294967295 4294967295 12 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 6 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 3 0 1 0 2 0 5 0 4 6 3 0 0 0 4294967295 4294967295 12 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 6 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 3 0 0 0 3 0 5 0 4 7 3 0 0 0 4294967295 4294967295 12 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 6 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 3 0 1 0 3 0 5 0 0
+
+DEAL:0::0 0 0 0 8 0 0 0 8 0 8 9 21 10 20 22 26 0 0 4294967295 8 8 1 21 11 20 12 26 23 0 0 4294967295 8 9 21 2 13 22 26 14 24 0 0 4294967295 8 21 11 13 3 26 23 24 15 0 0 4294967295 8 10 20 22 26 4 16 17 25 0 0 4294967295 8 20 12 26 23 16 5 25 18 0 0 4294967295 8 22 26 14 24 17 25 6 19 0 0 4294967295 8 26 23 24 15 25 18 19 7 0 0 4294967295 0 0 27 0 0 0 0 0 3 0 0 0.00000000000000000e+00 0.00000000000000000e+00 0.00000000000000000e+00 3 1.00000000000000000e+00 0.00000000000000000e+00 0.00000000000000000e+00 3 0.00000000000000000e+00 1.00000000000000000e+00 0.00000000000000000e+00 3 1.00000000000000000e+00 1.00000000000000000e+00 0.00000000000000000e+00 3 0.00000000000000000e+00 0.00000000000000000e+00 1.00000000000000000e+00 3 1.00000000000000000e+00 0.00000000000000000e+00 1.00000000000000000e+00 3 0.00000000000000000e+00 1.00000000000000000e+00 1.00000000000000000e+00 3 1.00000000000000000e+00 1.00000000000000000e+00 1.00000000000000000e+00 3 5.00000000000000000e-01 0.00000000000000000e+00 0.00000000000000000e+00 3 0.00000000000000000e+00 5.00000000000000000e-01 0.00000000000000000e+00 3 0.00000000000000000e+00 0.00000000000000000e+00 5.00000000000000000e-01 3 1.00000000000000000e+00 5.00000000000000000e-01 0.00000000000000000e+00 3 1.00000000000000000e+00 0.00000000000000000e+00 5.00000000000000000e-01 3 5.00000000000000000e-01 1.00000000000000000e+00 0.00000000000000000e+00 3 0.00000000000000000e+00 1.00000000000000000e+00 5.00000000000000000e-01 3 1.00000000000000000e+00 1.00000000000000000e+00 5.00000000000000000e-01 3 5.00000000000000000e-01 0.00000000000000000e+00 1.00000000000000000e+00 3 0.00000000000000000e+00 5.00000000000000000e-01 1.00000000000000000e+00 3 1.00000000000000000e+00 5.00000000000000000e-01 1.00000000000000000e+00 3 5.00000000000000000e-01 1.00000000000000000e+00 1.00000000000000000e+00 3 5.00000000000000000e-01 0.00000000000000000e+00 5.00000000000000000e-01 3 5.00000000000000000e-01 5.00000000000000000e-01 0.00000000000000000e+00 3 0.00000000000000000e+00 5.00000000000000000e-01 5.00000000000000000e-01 3 1.00000000000000000e+00 5.00000000000000000e-01 5.00000000000000000e-01 3 5.00000000000000000e-01 1.00000000000000000e+00 5.00000000000000000e-01 3 5.00000000000000000e-01 5.00000000000000000e-01 1.00000000000000000e+00 3 5.00000000000000000e-01 5.00000000000000000e-01 5.00000000000000000e-01 8 0 0 1 2 3 4 5 6 7 0 0 1 0 0 0 8 0 0 0 0 0 4 0 3 0 0 0 4294967295 4294967295 0 0 12 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 0 0 6 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 0 0 3 0 0 0 0 0 2 0 4 0 4 1 3 0 0 0 4294967295 4294967295 12 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 6 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 3 0 1 0 2 0 4 0 4 2 3 0 0 0 4294967295 4294967295 12 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 6 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 3 0 0 0 3 0 4 0 4 3 3 0 0 0 4294967295 4294967295 12 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 6 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 3 0 1 0 3 0 4 0 4 4 3 0 0 0 4294967295 4294967295 12 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 6 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 3 0 0 0 2 0 5 0 4 5 3 0 0 0 4294967295 4294967295 12 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 6 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 3 0 1 0 2 0 5 0 4 6 3 0 0 0 4294967295 4294967295 12 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 6 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 3 0 0 0 3 0 5 0 4 7 3 0 0 0 4294967295 4294967295 12 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 6 4294967295 4294967295 4294967295 4294967295 4294967295 4294967295 3 0 1 0 3 0 5 0 0
+
+DEAL:0::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.