]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Unpack functions take buffers as call-by-reference.
authorMarc Fehling <marc.fehling@gmx.net>
Fri, 6 Jul 2018 18:15:07 +0000 (12:15 -0600)
committerMarc Fehling <m.fehling@fz-juelich.de>
Sat, 7 Jul 2018 18:40:00 +0000 (12:40 -0600)
Pack functions return the packed data.

include/deal.II/base/quadrature_point_data.h
include/deal.II/distributed/solution_transfer.h
include/deal.II/distributed/tria.h
include/deal.II/particles/particle_handler.h
source/distributed/solution_transfer.cc
source/distributed/tria.cc
source/particles/particle_handler.cc
tests/mpi/attach_data_01.cc
tests/mpi/repartition_02.cc

index 3e07cd955bdadc74df668560adef551f24676265..de84d8d957450fbc0adb91818d751f102f75e642 100644 (file)
@@ -409,13 +409,12 @@ namespace parallel
        * objects that can later be retrieved after refinement, coarsening and
        * repartitioning.
        */
-      void
+      std::vector<char>
       pack_function(
         const typename parallel::distributed::Triangulation<dim>::cell_iterator
           &cell,
         const typename parallel::distributed::Triangulation<dim>::CellStatus
-                           status,
-        std::vector<char> &data);
+          status);
 
       /**
        * A callback function used to unpack the data on the current mesh that
@@ -429,7 +428,7 @@ namespace parallel
         const typename parallel::distributed::Triangulation<dim>::CellStatus
           status,
         const boost::iterator_range<std::vector<char>::const_iterator>
-          data_range);
+          &data_range);
 
       /**
        * FiniteElement used to project data from and to quadrature points.
@@ -798,8 +797,7 @@ namespace parallel
         &ContinuousQuadratureDataTransfer<dim, DataType>::pack_function,
         this,
         std::placeholders::_1,
-        std::placeholders::_2,
-        std::placeholders::_3));
+        std::placeholders::_2));
     }
 
 
@@ -825,13 +823,12 @@ namespace parallel
 
 
     template <int dim, typename DataType>
-    void
+    std::vector<char>
     ContinuousQuadratureDataTransfer<dim, DataType>::pack_function(
       const typename parallel::distributed::Triangulation<dim>::cell_iterator
         &cell,
       const typename parallel::distributed::Triangulation<
-        dim>::CellStatus /*status*/,
-      std::vector<char> &data)
+        dim>::CellStatus /*status*/)
     {
       pack_cell_data(cell, data_storage, matrix_quadrature);
 
@@ -840,7 +837,7 @@ namespace parallel
 
       // to get consistent data sizes on each cell for the fixed size transfer,
       // we won't allow compression
-      Utilities::pack(matrix_dofs, data, /*allow_compression=*/false);
+      return Utilities::pack(matrix_dofs, /*allow_compression=*/false);
     }
 
 
@@ -851,8 +848,9 @@ namespace parallel
       const typename parallel::distributed::Triangulation<dim>::cell_iterator
         &cell,
       const typename parallel::distributed::Triangulation<dim>::CellStatus
-                                                                     status,
-      const boost::iterator_range<std::vector<char>::const_iterator> data_range)
+        status,
+      const boost::iterator_range<std::vector<char>::const_iterator>
+        &data_range)
     {
       Assert((status !=
               parallel::distributed::Triangulation<dim, dim>::CELL_COARSEN),
index 72891ce1cb504636eb5a751f547fdd38513ac8da..756ec39b6383298ab1eb52fffadd9401980bbd7e 100644 (file)
@@ -255,13 +255,12 @@ namespace parallel
        * objects that can later be retrieved after refinement, coarsening and
        * repartitioning.
        */
-      void
+      std::vector<char>
       pack_callback(
         const typename Triangulation<dim, DoFHandlerType::space_dimension>::
           cell_iterator &cell,
         const typename Triangulation<dim, DoFHandlerType::space_dimension>::
-          CellStatus       status,
-        std::vector<char> &data);
+          CellStatus status);
 
       /**
        * A callback function used to unpack the data on the current mesh that
@@ -275,7 +274,7 @@ namespace parallel
         const typename Triangulation<dim, DoFHandlerType::space_dimension>::
           CellStatus status,
         const boost::iterator_range<std::vector<char>::const_iterator>
-                                   data_range,
+          &                        data_range,
         std::vector<VectorType *> &all_out);
 
 
index e3ffca184f2136543e1686115130991c3d797692..4c0cd6986b63e7e3d6a606f4721e48457ebbc655 100644 (file)
@@ -645,7 +645,7 @@ namespace parallel
        * of classes that do this is parallel::distributed::SolutionTransfer
        * where each parallel::distributed::SolutionTransfer object that works
        * on the current Triangulation object then needs to register its intent.
-       * Each of these parties registers a call-back function (the second
+       * Each of these parties registers a callback function (the second
        * argument here, @p pack_callback) that will be called whenever the
        * triangulation's execute_coarsening_and_refinement() or save()
        * functions are called.
@@ -715,8 +715,9 @@ namespace parallel
        *   using save() and load(), then the cell status argument with which
        *   the callback is called will always be `CELL_PERSIST`.
        *
-       * The third argument given to the callback is a reference to the
-       * buffer on which the packed data will be appended at the back.
+       * The callback function is expected to return a memory chunk of the
+       * format `std::vector<char>`, representing the packed data on a
+       * certain cell.
        *
        * @note The purpose of this function is to register intent to
        *   attach data for a single, subsequent call to
@@ -729,9 +730,9 @@ namespace parallel
        */
       unsigned int
       register_data_attach(
-        const std::function<void(const cell_iterator &,
-                                 const CellStatus,
-                                 std::vector<char> &)> &pack_callback);
+        const std::function<std::vector<char>(const cell_iterator &,
+                                              const CellStatus)>
+          &pack_callback);
 
       /**
        * This function is the opposite of register_data_attach(). It is called
@@ -784,10 +785,10 @@ namespace parallel
       void
       notify_ready_to_unpack(
         const unsigned int handle,
-        const std::function<
-          void(const cell_iterator &,
-               const CellStatus,
-               const boost::iterator_range<std::vector<char>::const_iterator>)>
+        const std::function<void(
+          const cell_iterator &,
+          const CellStatus,
+          const boost::iterator_range<std::vector<char>::const_iterator> &)>
           &unpack_callback);
 
       /**
@@ -888,10 +889,9 @@ namespace parallel
          */
         unsigned int n_attached_deserialize;
 
-        using pack_callback_t = std::function<
-          void(typename Triangulation<dim, spacedim>::cell_iterator,
-               CellStatus,
-               std::vector<char> &)>;
+        using pack_callback_t = std::function<std::vector<char>(
+          typename Triangulation<dim, spacedim>::cell_iterator,
+          CellStatus)>;
 
         /**
          * These callback functions will be stored in the order on how they have
@@ -1010,7 +1010,7 @@ namespace parallel
             const typename dealii::Triangulation<dim, spacedim>::cell_iterator
               &,
             const typename dealii::Triangulation<dim, spacedim>::CellStatus &,
-            const boost::iterator_range<std::vector<char>::const_iterator>)>
+            const boost::iterator_range<std::vector<char>::const_iterator> &)>
             &unpack_callback) const;
 
         /**
@@ -1289,10 +1289,10 @@ namespace parallel
        */
       unsigned int
       register_data_attach(
-        const std::function<void(
+        const std::function<std::vector<char>(
           const typename dealii::Triangulation<1, spacedim>::cell_iterator &,
-          const typename dealii::Triangulation<1, spacedim>::CellStatus,
-          std::vector<char> &)> &pack_callback);
+          const typename dealii::Triangulation<1, spacedim>::CellStatus)>
+          &pack_callback);
 
       /**
        * This function is not implemented, but needs to be present for the
@@ -1304,7 +1304,7 @@ namespace parallel
         const std::function<void(
           const typename dealii::Triangulation<1, spacedim>::cell_iterator &,
           const typename dealii::Triangulation<1, spacedim>::CellStatus,
-          const boost::iterator_range<std::vector<char>::const_iterator>)>
+          const boost::iterator_range<std::vector<char>::const_iterator> &)>
           &unpack_callback);
 
       /**
index 0b6f8963fc427f42b434b1b723975380364c9a3f..ade19553b58323fc70436af345eb2aec41d0324e 100644 (file)
@@ -507,11 +507,10 @@ namespace Particles
      * before a refinement step. All particles have to be attached to their
      * cell to be sent around to the new processes.
      */
-    void
+    std::vector<char>
     store_particles(
       const typename Triangulation<dim, spacedim>::cell_iterator &cell,
-      const typename Triangulation<dim, spacedim>::CellStatus     status,
-      std::vector<char> &                                         data) const;
+      const typename Triangulation<dim, spacedim>::CellStatus     status) const;
 
     /**
      * Called by listener functions after a refinement step. The local map
@@ -522,7 +521,7 @@ namespace Particles
       const typename Triangulation<dim, spacedim>::cell_iterator &cell,
       const typename Triangulation<dim, spacedim>::CellStatus     status,
       const boost::iterator_range<std::vector<char>::const_iterator>
-        data_range);
+        &data_range);
   };
 
   /* ---------------------- inline and template functions ------------------ */
index add27f1cdd9bd5b82cf8b398a287fd59caf80c39..349ffbae7a645ff42107038b3c69e895c756c966 100644 (file)
@@ -89,8 +89,7 @@ namespace parallel
         &SolutionTransfer<dim, VectorType, DoFHandlerType>::pack_callback,
         this,
         std::placeholders::_1,
-        std::placeholders::_2,
-        std::placeholders::_3));
+        std::placeholders::_2));
     }
 
 
@@ -202,13 +201,12 @@ namespace parallel
 
 
     template <int dim, typename VectorType, typename DoFHandlerType>
-    void
+    std::vector<char>
     SolutionTransfer<dim, VectorType, DoFHandlerType>::pack_callback(
       const typename Triangulation<dim, DoFHandlerType::space_dimension>::
         cell_iterator &cell_,
       const typename Triangulation<dim, DoFHandlerType::space_dimension>::
-        CellStatus /*status*/,
-      std::vector<char> &data)
+        CellStatus /*status*/)
     {
       typename DoFHandlerType::cell_iterator cell(*cell_, dof_handler);
 
@@ -228,7 +226,7 @@ namespace parallel
 
       // to get consistent data sizes on each cell for the fixed size transfer,
       // we won't allow compression
-      Utilities::pack(dofvalues, data, /*allow_compression=*/false);
+      return Utilities::pack(dofvalues, /*allow_compression=*/false);
     }
 
 
@@ -240,8 +238,9 @@ namespace parallel
         cell_iterator &cell_,
       const typename Triangulation<dim, DoFHandlerType::space_dimension>::
         CellStatus /*status*/,
-      const boost::iterator_range<std::vector<char>::const_iterator> data_range,
-      std::vector<VectorType *> &                                    all_out)
+      const boost::iterator_range<std::vector<char>::const_iterator>
+        &                        data_range,
+      std::vector<VectorType *> &all_out)
     {
       typename DoFHandlerType::cell_iterator cell(*cell_, dof_handler);
 
index c4b959893be0ff5630ba35a62a492da6902d9fb4..7fa7c55198c612b144fc9d33e720d8664e57693e 100644 (file)
@@ -1243,7 +1243,7 @@ namespace parallel
                      callback_it != pack_callbacks.cend();
                      ++callback_it, ++data_it)
                   {
-                    (*callback_it)(dealii_cell, cell_status, *data_it);
+                    *data_it = (*callback_it)(dealii_cell, cell_status);
                   }
               }
           } // loop over quad_cell_relations
@@ -1420,7 +1420,7 @@ namespace parallel
       const std::function<void(
         const typename dealii::Triangulation<dim, spacedim>::cell_iterator &,
         const typename dealii::Triangulation<dim, spacedim>::CellStatus &,
-        const boost::iterator_range<std::vector<char>::const_iterator>)>
+        const boost::iterator_range<std::vector<char>::const_iterator> &)>
         &unpack_callback) const
     {
       Assert(cumulative_sizes_fixed.size() > 0,
@@ -3848,9 +3848,8 @@ namespace parallel
     template <int dim, int spacedim>
     unsigned int
     Triangulation<dim, spacedim>::register_data_attach(
-      const std::function<void(const cell_iterator &,
-                               const CellStatus,
-                               std::vector<char> &)> &pack_callback)
+      const std::function<std::vector<char>(const cell_iterator &,
+                                            const CellStatus)> &pack_callback)
     {
       // add new callback_set
       cell_attached_data.pack_callbacks.push_back(pack_callback);
@@ -3869,7 +3868,7 @@ namespace parallel
       const std::function<
         void(const cell_iterator &,
              const CellStatus,
-             const boost::iterator_range<std::vector<char>::const_iterator>)>
+             const boost::iterator_range<std::vector<char>::const_iterator> &)>
         &unpack_callback)
     {
 #  ifdef DEBUG
@@ -4458,10 +4457,10 @@ namespace parallel
     template <int spacedim>
     unsigned int
     Triangulation<1, spacedim>::register_data_attach(
-      const std::function<
-        void(const typename dealii::Triangulation<1, spacedim>::cell_iterator &,
-             const typename dealii::Triangulation<1, spacedim>::CellStatus,
-             std::vector<char> &)> & /*pack_callback*/)
+      const std::function<std::vector<char>(
+        const typename dealii::Triangulation<1, spacedim>::cell_iterator &,
+        const typename dealii::Triangulation<1, spacedim>::CellStatus)>
+        & /*pack_callback*/)
     {
       Assert(false, ExcNotImplemented());
       return 0;
@@ -4476,7 +4475,7 @@ namespace parallel
       const std::function<
         void(const typename dealii::Triangulation<1, spacedim>::cell_iterator &,
              const typename dealii::Triangulation<1, spacedim>::CellStatus,
-             const boost::iterator_range<std::vector<char>::const_iterator>)>
+             const boost::iterator_range<std::vector<char>::const_iterator> &)>
         & /*unpack_callback*/)
     {
       Assert(false, ExcNotImplemented());
index 0941499e2fccb6d3d8dd5dc340ec75be362dd167..2b4131ef04e9f479e18862e1f57e20a64062257c 100644 (file)
@@ -1037,16 +1037,14 @@ namespace Particles
 
     if (global_max_particles_per_cell > 0)
       {
-        const std::function<
-          void(const typename Triangulation<dim, spacedim>::cell_iterator &,
-               const typename Triangulation<dim, spacedim>::CellStatus,
-               std::vector<char> &)>
+        const std::function<std::vector<char>(
+          const typename Triangulation<dim, spacedim>::cell_iterator &,
+          const typename Triangulation<dim, spacedim>::CellStatus)>
           callback_function =
             std::bind(&ParticleHandler<dim, spacedim>::store_particles,
                       std::cref(*this),
                       std::placeholders::_1,
-                      std::placeholders::_2,
-                      std::placeholders::_3);
+                      std::placeholders::_2);
 
         handle =
           non_const_triangulation->register_data_attach(callback_function);
@@ -1075,16 +1073,14 @@ namespace Particles
     // data correctly. Only do this if something was actually stored.
     if (serialization && (global_max_particles_per_cell > 0))
       {
-        const std::function<
-          void(const typename Triangulation<dim, spacedim>::cell_iterator &,
-               const typename Triangulation<dim, spacedim>::CellStatus,
-               std::vector<char> &)>
+        const std::function<std::vector<char>(
+          const typename Triangulation<dim, spacedim>::cell_iterator &,
+          const typename Triangulation<dim, spacedim>::CellStatus)>
           callback_function =
             std::bind(&ParticleHandler<dim, spacedim>::store_particles,
                       std::cref(*this),
                       std::placeholders::_1,
-                      std::placeholders::_2,
-                      std::placeholders::_3);
+                      std::placeholders::_2);
 
         handle =
           non_const_triangulation->register_data_attach(callback_function);
@@ -1093,10 +1089,10 @@ namespace Particles
     // Check if something was stored and load it
     if (handle != numbers::invalid_unsigned_int)
       {
-        const std::function<
-          void(const typename Triangulation<dim, spacedim>::cell_iterator &,
-               const typename Triangulation<dim, spacedim>::CellStatus,
-               const boost::iterator_range<std::vector<char>::const_iterator>)>
+        const std::function<void(
+          const typename Triangulation<dim, spacedim>::cell_iterator &,
+          const typename Triangulation<dim, spacedim>::CellStatus,
+          const boost::iterator_range<std::vector<char>::const_iterator> &)>
           callback_function =
             std::bind(&ParticleHandler<dim, spacedim>::load_particles,
                       std::ref(*this),
@@ -1117,11 +1113,10 @@ namespace Particles
 
 
   template <int dim, int spacedim>
-  void
+  std::vector<char>
   ParticleHandler<dim, spacedim>::store_particles(
     const typename Triangulation<dim, spacedim>::cell_iterator &cell,
-    const typename Triangulation<dim, spacedim>::CellStatus     status,
-    std::vector<char> &data_vector) const
+    const typename Triangulation<dim, spacedim>::CellStatus     status) const
   {
     // -------------------
     // HOTFIX: resize memory and static_cast std::vector to void*
@@ -1149,9 +1144,8 @@ namespace Particles
       (size_per_particle * global_max_particles_per_cell) *
         (serialization ? 1 : std::pow(2, dim));
 
-    const size_t previous_size = data_vector.size();
-    data_vector.resize(previous_size + transfer_size_per_cell);
-    void *data = static_cast<void *>(data_vector.data() + previous_size);
+    std::vector<char> data_vector(transfer_size_per_cell);
+    void *            data = static_cast<void *>(data_vector.data());
     // --------------------
 
     unsigned int n_particles(0);
@@ -1217,14 +1211,16 @@ namespace Particles
       }
     else
       Assert(false, ExcInternalError());
+
+    return data_vector;
   }
 
   template <int dim, int spacedim>
   void
   ParticleHandler<dim, spacedim>::load_particles(
-    const typename Triangulation<dim, spacedim>::cell_iterator &   cell,
-    const typename Triangulation<dim, spacedim>::CellStatus        status,
-    const boost::iterator_range<std::vector<char>::const_iterator> data_range)
+    const typename Triangulation<dim, spacedim>::cell_iterator &    cell,
+    const typename Triangulation<dim, spacedim>::CellStatus         status,
+    const boost::iterator_range<std::vector<char>::const_iterator> &data_range)
   {
     // -------------------
     // HOTFIX: static_cast std::vector to void*
index 06fdafdca42735eadf9f238d3717eb23a583c8a0..9c556459bd060927fab5adc27993518e207ffcac 100644 (file)
 
 
 template <int dim>
-void
+std::vector<char>
 pack_function(
   const typename parallel::distributed::Triangulation<dim, dim>::cell_iterator
     &cell,
   const typename parallel::distributed::Triangulation<dim, dim>::CellStatus
-                     status,
-  std::vector<char> &data)
+    status)
 {
   static int some_number = 0;
   deallog << "packing cell " << cell->id() << " with data=" << some_number
@@ -63,9 +62,7 @@ pack_function(
       Assert(!cell->has_children(), ExcInternalError());
     }
 
-  Utilities::pack(some_number, data, /*allow_compression=*/false);
-
-  ++some_number;
+  return Utilities::pack(some_number++, /*allow_compression=*/false);
 }
 
 template <int dim>
@@ -74,8 +71,8 @@ unpack_function(
   const typename parallel::distributed::Triangulation<dim, dim>::cell_iterator
     &cell,
   const typename parallel::distributed::Triangulation<dim, dim>::CellStatus
-                                                                 status,
-  const boost::iterator_range<std::vector<char>::const_iterator> data_range)
+                                                                  status,
+  const boost::iterator_range<std::vector<char>::const_iterator> &data_range)
 {
   const int intdata = Utilities::unpack<int>(data_range.begin(),
                                              data_range.end(),
index e01a3fd2843c9975389cd99650c1f90bd3e7f7a5..ce953cc88eae06bec6e4239488ca57a392f564d4 100644 (file)
 
 
 template <int dim>
-void
+std::vector<char>
 pack_function(
   const typename parallel::distributed::Triangulation<dim, dim>::cell_iterator
     &cell,
   const typename parallel::distributed::Triangulation<dim, dim>::CellStatus
-                     status,
-  std::vector<char> &data)
+    status)
 {
   static int some_number = cell->index();
   deallog << "packing cell " << cell->id() << " with data=" << some_number
@@ -63,9 +62,7 @@ pack_function(
       Assert(!cell->has_children(), ExcInternalError());
     }
 
-  Utilities::pack(some_number, data, /*allow_compression=*/false);
-
-  ++some_number;
+  return Utilities::pack(some_number++, /*allow_compression=*/false);
 }
 
 template <int dim>
@@ -74,8 +71,8 @@ unpack_function(
   const typename parallel::distributed::Triangulation<dim, dim>::cell_iterator
     &cell,
   const typename parallel::distributed::Triangulation<dim, dim>::CellStatus
-                                                                 status,
-  const boost::iterator_range<std::vector<char>::const_iterator> data_range)
+                                                                  status,
+  const boost::iterator_range<std::vector<char>::const_iterator> &data_range)
 {
   const int number = Utilities::unpack<int>(data_range.begin(),
                                             data_range.end(),

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.