]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Rename parallel::distributed::SolutionTransfer::prepare_serialization()
authorDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Wed, 27 Mar 2019 02:10:05 +0000 (22:10 -0400)
committerDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Wed, 27 Mar 2019 02:23:49 +0000 (22:23 -0400)
doc/news/changes/minor/20190326Arndt [new file with mode: 0644]
include/deal.II/distributed/solution_transfer.h
source/distributed/solution_transfer.cc
tests/mpi/p4est_save_02.cc
tests/mpi/p4est_save_03.cc
tests/mpi/p4est_save_04.cc
tests/mpi/p4est_save_06.cc

diff --git a/doc/news/changes/minor/20190326Arndt b/doc/news/changes/minor/20190326Arndt
new file mode 100644 (file)
index 0000000..8e9c7ef
--- /dev/null
@@ -0,0 +1,5 @@
+Changed: parallel::distributed::SolutionTransfer::prepare_serialization() has
+been deprecated in favor of
+parallel::distributed::SolutionTransfer::prepare_for_serialization().
+<br<
+(Daniel Arndt, 2019/03/26)
index 9d61878b2770cdbe0f354059b8ca47e8ceee0087..ebc007c3d2162d6ddb5fa93541fe07dc5a76c9dc 100644 (file)
@@ -43,9 +43,9 @@ namespace parallel
      * <h3>Note on ghost elements</h3> In a parallel computation PETSc or
      * Trilinos vector may contain ghost elements or not. For reading in
      * information with prepare_for_coarsening_and_refinement() or
-     * prepare_serialization() you need to supply vectors with ghost elements,
-     * so that all locally_active elements can be read. On the other hand,
-     * ghosted vectors are generally not writable, so for calls to
+     * prepare_for_serialization() you need to supply vectors with ghost
+     * elements, so that all locally_active elements can be read. On the other
+     * hand, ghosted vectors are generally not writable, so for calls to
      * interpolate() or deserialize() you need to supply distributed vectors
      * without ghost elements. More precisely, during interpolation the
      * current algorithm writes into all locally active degrees of freedom.
@@ -133,7 +133,7 @@ namespace parallel
      * @code
      * parallel::distributed::SolutionTransfer<dim,VectorType>
      *   sol_trans(dof_handler);
-     * sol_trans.prepare_serialization (vector);
+     * sol_trans.prepare_for_serialization (vector);
      *
      * triangulation.save(filename);
      * @endcode
@@ -186,7 +186,7 @@ namespace parallel
      *     sol_trans(hp_dof_handler);
      *
      * hp_dof_handler.prepare_for_serialization_of_active_fe_indices();
-     * sol_trans.prepare_serialization(vector);
+     * sol_trans.prepare_for_serialization(vector);
      *
      * triangulation.save(filename);
      * @endcode
@@ -292,24 +292,42 @@ namespace parallel
       void
       interpolate(VectorType &out);
 
+      /**
+       * Prepare the serialization of the given vector. The serialization is
+       * done by Triangulation::save(). The given vector needs all information
+       * on the locally active DoFs (it must be ghosted). See documentation of
+       * this class for more information.
+       */
+      void
+      prepare_for_serialization(const VectorType &in);
+
+      /**
+       * Same as the function above, only for a list of vectors.
+       */
+      void
+      prepare_for_serialization(const std::vector<const VectorType *> &all_in);
 
       /**
        * Prepare the serialization of the given vector. The serialization is
        * done by Triangulation::save(). The given vector needs all information
        * on the locally active DoFs (it must be ghosted). See documentation of
        * this class for more information.
+       *
+       * @deprecated Use parallel::distributed::SolutionTransfer::prepare_for_serialization() instead.
        */
+      DEAL_II_DEPRECATED
       void
       prepare_serialization(const VectorType &in);
 
-
       /**
        * Same as the function above, only for a list of vectors.
+       *
+       * @deprecated Use parallel::distributed::SolutionTransfer::prepare_for_serialization() instead.
        */
+      DEAL_II_DEPRECATED
       void
       prepare_serialization(const std::vector<const VectorType *> &all_in);
 
-
       /**
        * Execute the deserialization of the given vector. This needs to be
        * done after calling Triangulation::load(). The given vector must be a
index 4c4f5d9f26f18cad37d4dff60df7d4dac01df989..0919d7c798048a4180e169ff5f7e5b2c492d5401 100644 (file)
@@ -109,13 +109,34 @@ namespace parallel
 
 
 
+    template <int dim, typename VectorType, typename DoFHandlerType>
+    void
+    SolutionTransfer<dim, VectorType, DoFHandlerType>::
+      prepare_for_serialization(const VectorType &in)
+    {
+      std::vector<const VectorType *> all_in(1, &in);
+      prepare_for_serialization(all_in);
+    }
+
+
+
+    template <int dim, typename VectorType, typename DoFHandlerType>
+    void
+    SolutionTransfer<dim, VectorType, DoFHandlerType>::
+      prepare_for_serialization(const std::vector<const VectorType *> &all_in)
+    {
+      prepare_for_coarsening_and_refinement(all_in);
+    }
+
+
+
     template <int dim, typename VectorType, typename DoFHandlerType>
     void
     SolutionTransfer<dim, VectorType, DoFHandlerType>::prepare_serialization(
       const VectorType &in)
     {
       std::vector<const VectorType *> all_in(1, &in);
-      prepare_serialization(all_in);
+      prepare_for_serialization(all_in);
     }
 
 
index 1ab383a4d77a067137f62c957ff4be48f33e2ad2..be0738f39f493531d26c4a148b59cbcf5f4473fa 100644 (file)
@@ -97,7 +97,7 @@ test()
     x.compress(VectorOperation::insert);
     solution = x;
 
-    soltrans.prepare_serialization(solution);
+    soltrans.prepare_for_serialization(solution);
 
     tr.save(filename.c_str());
 
index 9780e2dfaeac54edcdf8fcb025c757309708ec4a..4619f2ad7536cb51a594c2351e53bc6b187fc6e2 100644 (file)
@@ -108,8 +108,8 @@ test()
     solution2 = x;
 
 
-    soltrans.prepare_serialization(solution);
-    soltrans2.prepare_serialization(solution2);
+    soltrans.prepare_for_serialization(solution);
+    soltrans2.prepare_for_serialization(solution2);
 
     tr.save(filename.c_str());
 
index c908edf673410c2d513d8ecb5af3de632e988b05..0526ce29c2a62d8cc079a8bd6c9c07575e6da2d0 100644 (file)
@@ -100,7 +100,7 @@ test()
       x.compress(VectorOperation::insert);
       rel_x = x;
 
-      soltrans.prepare_serialization(rel_x);
+      soltrans.prepare_for_serialization(rel_x);
 
       tr.save("file");
       //    tr.write_mesh_vtk("before");
index 453ed583896d3190f28726cb1f185edea8230d2a..15025baeaefe5b8295a97875347711b84c8ecc77 100644 (file)
@@ -107,7 +107,7 @@ test()
       rel_x = x;
 
       dh.prepare_for_serialization_of_active_fe_indices();
-      soltrans.prepare_serialization(rel_x);
+      soltrans.prepare_for_serialization(rel_x);
 
       tr.save("file");
       deallog << "#cells: " << tr.n_global_active_cells() << std::endl

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.