]> https://gitweb.dealii.org/ - dealii.git/commitdiff
cell_weight: deprecate old signal `cell_weight` in favor of `weight`.
authorMarc Fehling <mafehling.git@gmail.com>
Fri, 11 Mar 2022 20:28:34 +0000 (13:28 -0700)
committerMarc Fehling <mafehling.git@gmail.com>
Mon, 14 Mar 2022 03:05:39 +0000 (21:05 -0600)
17 files changed:
doc/news/changes/incompatibilities/20220225Fehling
examples/step-68/step-68.cc
include/deal.II/distributed/cell_weights.h
include/deal.II/distributed/tria.h
include/deal.II/grid/grid_tools.h
include/deal.II/grid/tria.h
source/distributed/cell_weights.cc
source/distributed/tria.cc
source/grid/grid_tools.cc
tests/mpi/cell_weights_01_back_and_forth_01.cc
tests/mpi/cell_weights_01_back_and_forth_02.cc
tests/mpi/cell_weights_02.cc
tests/mpi/cell_weights_03.cc
tests/mpi/cell_weights_04.cc
tests/mpi/cell_weights_05.cc
tests/mpi/cell_weights_06.cc
tests/particles/step-68.cc

index 4ef1b4f57dc66afc0fbe88777143d3cca12c66b9..8ce72bc5d80d32238791bb3a1aa9b10604b719b2 100644 (file)
@@ -1,10 +1,13 @@
 Changed: For weighted load balancing with parallel::distributed::Triangulation
 objects, an intial weight of `1000` will no longer be assigned to each cell.
 <br>
-You can invoke the old behavior by connecting a function to the signal
-that returns the base weight like this:
+The old signal Triangulation::Signals::cell_weight has been deprecated.
+Please use the new signal Triangulation::Signals::weight instead.
+<br>
+You can invoke the old behavior by connecting a function to the new
+signal that returns the base weight like this:
 @code{.cc}
-triangulation.signals.cell_weight.connect(
+triangulation.signals.weight.connect(
   [](const typename parallel::distributed::Triangulation<dim>::cell_iterator &,
      const typename parallel::distributed::Triangulation<dim>::CellStatus)
     -> unsigned int { return 1000; });
index bde1d7d5a47b91762b9bcb2bc00b303a72e96b24..eb06cd335d782350ec5949fb0c9a5b0dd7a758cd 100644 (file)
@@ -322,7 +322,7 @@ namespace Step68
   // return value of this function (representing "work for this cell") is
   // calculated based on the number of particles in the current cell.
   // The function is
-  // connected to the cell_weight() signal inside the triangulation, and will be
+  // connected to the `weight` signal inside the triangulation, and will be
   // called once per cell, whenever the triangulation repartitions the domain
   // between ranks (the connection is created inside the
   // generate_particles() function of this class).
@@ -398,7 +398,7 @@ namespace Step68
     // be created once, so we might as well have set it up in the constructor
     // of this class, but for the purpose of this example we want to group the
     // particle related instructions.
-    background_triangulation.signals.cell_weight.connect(
+    background_triangulation.signals.weight.connect(
       [&](
         const typename parallel::distributed::Triangulation<dim>::cell_iterator
           &cell,
index 8fce0349f3cf18ce75cb4862035d97caa40b7494..79342bc6da43d3c223f2c4bdabdebed5d90d02fb 100644 (file)
@@ -65,7 +65,7 @@ namespace parallel
    * this case, an analogous code example looks as follows.
    * @code
    * boost::signals2::connection connection =
-   *   hp_dof_handler.get_triangulation().signals.cell_weight.connect(
+   *   hp_dof_handler.get_triangulation().signals.weight.connect(
    *     parallel::CellWeights<dim, spacedim>::make_weighting_callback(
    *       hp_dof_handler,
    *       parallel::CellWeights<dim, spacedim>::ndofs_weighting({1, 1}));
@@ -73,7 +73,7 @@ namespace parallel
    *
    * The use of this class is demonstrated in step-75.
    *
-   * @note See Triangulation::Signals::cell_weight for more information on
+   * @note See Triangulation::Signals::weight for more information on
    * weighting and load balancing.
    *
    * @note Be aware that this class connects the weight function to the
@@ -85,7 +85,7 @@ namespace parallel
    *
    * @note A hp::FECollection needs to be attached to your DoFHandler object via
    * DoFHandler::distribute_dofs() <em>once before</em> the
-   * Triangulation::Signals::cell_weight signal will be triggered. Otherwise,
+   * Triangulation::Signals::weight signal will be triggered. Otherwise,
    * your DoFHandler does not know many degrees of freedom your cells have. In
    * other words, you need to call DoFHandler::distribute_dofs() once before you
    * call
@@ -197,13 +197,13 @@ namespace parallel
 
   private:
     /**
-     * A connection to the corresponding cell_weight signal of the Triangulation
+     * A connection to the corresponding weight signal of the Triangulation
      * which is attached to the DoFHandler.
      */
     boost::signals2::connection connection;
 
     /**
-     * A callback function that will be connected to the cell_weight signal of
+     * A callback function that will be connected to the weight signal of
      * the @p triangulation, to which the @p dof_handler is attached. Ultimately
      * returns the weight for each cell, determined by the @p weighting_function
      * provided as a parameter. Returns zero if @p dof_handler has not been
index e0f62e886528f1e1a2202bb64dcf367a6097c703..a3df6166560c17ce47b66e6f641ac2642a9b720a 100644 (file)
@@ -499,7 +499,7 @@ namespace parallel
        * @note This function by default partitions the mesh in such a way that
        * the number of cells on all processors is roughly equal. If you want
        * to set weights for partitioning, e.g. because some cells are more
-       * expensive to compute than others, you can use the signal cell_weight
+       * expensive to compute than others, you can use the signal weight
        * as documented in the dealii::Triangulation class. This function will
        * check whether a function is connected to the signal and if so use it.
        * If you prefer to repartition the mesh yourself at user-defined
@@ -508,7 +508,7 @@ namespace parallel
        * flag to the constructor, which ensures that calling the current
        * function only refines and coarsens the triangulation, but doesn't
        * partition it. You can then call the repartition() function manually.
-       * The usage of the cell_weights signal is identical in both cases, if a
+       * The usage of the weight signal is identical in both cases, if a
        * function is connected to the signal it will be used to balance the
        * calculated weights, otherwise the number of cells is balanced.
        */
@@ -542,7 +542,7 @@ namespace parallel
        * the same way as execute_coarsening_and_refinement() with respect to
        * dealing with data movement (SolutionTransfer, etc.).
        *
-       * @note If no function is connected to the cell_weight signal described
+       * @note If no function is connected to the weight signal described
        * in the dealii::Triangulation class, this function will balance the
        * number of cells on each processor. If one or more functions are
        * connected, it will calculate the sum of the weights and balance the
index d02bd22642d1a1474f4243188c4ffb53108ac077..858090e6dc29a91ca1614b1423211f38dcb31b5d 100644 (file)
@@ -2163,7 +2163,7 @@ namespace GridTools
    * single-partition case without packages installed, and only requires them
    * installed when multiple partitions are required.
    *
-   * @note If the @p cell_weight signal has been attached to the @p triangulation,
+   * @note If the @p weight signal has been attached to the @p triangulation,
    * then this will be used and passed to the partitioner.
    */
   template <int dim, int spacedim>
@@ -2233,7 +2233,7 @@ namespace GridTools
    * case like this, partitioning algorithm may sometimes make bad decisions and
    * you may want to build your own connectivity graph.
    *
-   * @note If the @p cell_weight signal has been attached to the @p triangulation,
+   * @note If the @p weight signal has been attached to the @p triangulation,
    * then this will be used and passed to the partitioner.
    */
   template <int dim, int spacedim>
index 7badd095d0723f6d67932feb1b6deb828e322e48..669e2657ec0b9c28175288e9fff4b7ff7925d4e9 100644 (file)
@@ -2044,7 +2044,7 @@ public:
   };
 
   /**
-   * A structure used to accumulate the results of the cell_weights slot
+   * A structure used to accumulate the results of the weight signal slot
    * functions below. It takes an iterator range and returns the sum of
    * values.
    */
@@ -2203,7 +2203,126 @@ public:
     boost::signals2::signal<unsigned int(const cell_iterator &,
                                          const CellStatus),
                             CellWeightSum<unsigned int>>
-      cell_weight;
+      weight;
+
+#ifndef DOXYGEN
+    /**
+     * Legacy constructor that connects deprecated signals to their new ones.
+     */
+    Signals()
+      : cell_weight(weight)
+    {}
+
+    /**
+     * Legacy signal emulation to deprecate the old signal.
+     */
+    class LegacySignal
+    {
+    public:
+      using signature_type = unsigned int(const cell_iterator &,
+                                          const CellStatus);
+      using combiner_type  = CellWeightSum<unsigned int>;
+
+      using slot_function_type = boost::function<signature_type>;
+      using slot_type =
+        boost::signals2::slot<signature_type, slot_function_type>;
+
+      LegacySignal(
+        boost::signals2::signal<signature_type, combiner_type> &new_signal)
+        : new_signal(new_signal)
+      {}
+
+      ~LegacySignal()
+      {
+        base_weight.disconnect();
+      }
+
+      DEAL_II_DEPRECATED_EARLY
+      boost::signals2::connection
+      connect(
+        const slot_type &                 slot,
+        boost::signals2::connect_position position = boost::signals2::at_back)
+      {
+        if (base_weight.connected() == false)
+          {
+            base_weight = new_signal.connect(
+              [](const cell_iterator &, const CellStatus) -> unsigned int {
+                return 1000;
+              });
+            Assert(base_weight.connected() && new_signal.num_slots() == 1,
+                   ExcInternalError());
+          }
+
+        return new_signal.connect(slot, position);
+      }
+
+      DEAL_II_DEPRECATED_EARLY
+      std::size_t
+      num_slots() const
+      {
+        return new_signal.num_slots() -
+               static_cast<std::size_t>(base_weight.connected());
+      }
+
+      DEAL_II_DEPRECATED_EARLY
+      bool
+      empty() const
+      {
+        if (num_slots() == 0)
+          {
+            Assert(new_signal.num_slots() == 0, ExcInternalError());
+            return true;
+          }
+        return false;
+      }
+
+      template <typename S>
+      DEAL_II_DEPRECATED_EARLY void
+      disconnect(const S &connection)
+      {
+        new_signal.disconnect(connection);
+
+        if (num_slots() == 0)
+          {
+            Assert(base_weight.connected() && new_signal.num_slots() == 1,
+                   ExcInternalError());
+            new_signal.disconnect(base_weight);
+          }
+      }
+
+      DEAL_II_DEPRECATED_EARLY
+      unsigned int
+      operator()(const cell_iterator &iterator, const CellStatus status)
+      {
+        return new_signal(iterator, status);
+      }
+
+    private:
+      boost::signals2::connection base_weight;
+
+      boost::signals2::signal<signature_type, combiner_type> &new_signal;
+    };
+#endif
+
+    /**
+     * @copydoc weight
+     *
+     * As a reference a value of 1000 is added for every cell to the total
+     * weight. This means a signal return value of 1000 (resulting in a weight
+     * of 2000) means that it is twice as expensive for a process to handle this
+     * particular cell.
+     *
+     * @deprecated Use the weight signal instead which omits the base weight.
+     * You can invoke the old behavior by connecting a function to the signal
+     * that returns the base weight like this:
+     * @code{.cc}
+     * triangulation.signals.weight.connect(
+     *   [](const typename Triangulation<dim>::cell_iterator &,
+     *      const typename Triangulation<dim>::CellStatus)
+     *     -> unsigned int { return 1000; });
+     * @endcode
+     */
+    LegacySignal cell_weight;
 
     /**
      * This signal is triggered at the beginning of execution of the
index 2218f2ee42301aa43926bfe65e45c066e7fc8bb8..2271efcf674e48c0039ddd1683e50d81ef147cb7 100644 (file)
@@ -50,7 +50,7 @@ namespace parallel
       &weighting_function)
   {
     connection.disconnect();
-    connection = dof_handler.get_triangulation().signals.cell_weight.connect(
+    connection = dof_handler.get_triangulation().signals.weight.connect(
       make_weighting_callback(dof_handler, weighting_function));
   }
 
index 617c9ec3b62f79834fccef5f00ad598fc495a200..773d56fe1aea5258c26955ef3826db32967708bf 100644 (file)
@@ -1453,7 +1453,7 @@ namespace
   {
   public:
     /**
-     * This constructor assumes the cell_weights are already sorted in the
+     * This constructor assumes the @p cell_weights are already sorted in the
      * order that p4est will encounter the cells, and they do not contain
      * ghost cells or artificial cells.
      */
@@ -3351,7 +3351,7 @@ namespace parallel
         {
           // partition the new mesh between all processors. If cell weights
           // have not been given balance the number of cells.
-          if (this->signals.cell_weight.num_slots() == 0)
+          if (this->signals.weight.empty())
             dealii::internal::p4est::functions<dim>::partition(
               parallel_forest,
               /* prepare coarsening */ 1,
@@ -3529,7 +3529,7 @@ namespace parallel
                         (parallel_forest->mpisize + 1));
         }
 
-      if (this->signals.cell_weight.num_slots() == 0)
+      if (this->signals.weight.empty())
         {
           // no cell weights given -- call p4est's 'partition' without a
           // callback for cell weights
@@ -4056,16 +4056,16 @@ namespace parallel
 
       // Iterate over p4est and Triangulation relations
       // to find refined/coarsened/kept
-      // cells. Then append cell_weight.
+      // cells. Then append weight.
       // Note that we need to follow the p4est ordering
-      // instead of the deal.II ordering to get the cell_weights
+      // instead of the deal.II ordering to get the weights
       // in the same order p4est will encounter them during repartitioning.
       for (const auto &cell_rel : this->local_cell_relations)
         {
           const auto &cell_it     = cell_rel.first;
           const auto &cell_status = cell_rel.second;
 
-          weights.push_back(this->signals.cell_weight(cell_it, cell_status));
+          weights.push_back(this->signals.weight(cell_it, cell_status));
         }
 
       return weights;
index a9566e6d7aac9923c085dfe3b1c55fb7235a2a8d..4e77809375f3cdb131644dc44c95cb2fc023acc0 100644 (file)
@@ -3816,7 +3816,7 @@ namespace GridTools
     std::vector<unsigned int> cell_weights;
 
     // Get cell weighting if a signal has been attached to the triangulation
-    if (!triangulation.signals.cell_weight.empty())
+    if (!triangulation.signals.weight.empty())
       {
         cell_weights.resize(triangulation.n_active_cells(), 0U);
 
@@ -3826,7 +3826,7 @@ namespace GridTools
         for (const auto &cell : triangulation.active_cell_iterators())
           if (cell->is_locally_owned())
             cell_weights[cell->active_cell_index()] =
-              triangulation.signals.cell_weight(
+              triangulation.signals.weight(
                 cell, Triangulation<dim, spacedim>::CellStatus::CELL_PERSIST);
 
         // If this is a parallel triangulation, we then need to also
@@ -3920,7 +3920,7 @@ namespace GridTools
     std::vector<unsigned int> cell_weights;
 
     // Get cell weighting if a signal has been attached to the triangulation
-    if (!triangulation.signals.cell_weight.empty())
+    if (!triangulation.signals.weight.empty())
       {
         cell_weights.resize(triangulation.n_active_cells(), 0U);
 
@@ -3930,7 +3930,7 @@ namespace GridTools
         for (const auto &cell : triangulation.active_cell_iterators() |
                                   IteratorFilters::LocallyOwnedCell())
           cell_weights[cell->active_cell_index()] =
-            triangulation.signals.cell_weight(
+            triangulation.signals.weight(
               cell, Triangulation<dim, spacedim>::CellStatus::CELL_PERSIST);
 
         // If this is a parallel triangulation, we then need to also
@@ -4060,7 +4060,7 @@ namespace GridTools
                       "are already partitioned implicitly and can not be "
                       "partitioned again explicitly."));
     Assert(n_partitions > 0, ExcInvalidNumberOfPartitions(n_partitions));
-    Assert(triangulation.signals.cell_weight.empty(), ExcNotImplemented());
+    Assert(triangulation.signals.weight.empty(), ExcNotImplemented());
 
     // signal that partitioning is going to happen
     triangulation.signals.pre_partition();
index bdd6110dc5b26441a2338eb583e593d80fed47e0..b965a3a19a5ca58e6c9e6c0a101195efcd2b47ac 100644 (file)
@@ -75,16 +75,12 @@ test()
 
   // repartition the mesh as described above, first in some arbitrary
   // way, and then with all equal weights
-  tr.signals.cell_weight.connect(std::bind(&cell_weight_1<dim>,
-                                           std::placeholders::_1,
-                                           std::placeholders::_2));
+  tr.signals.weight.connect(&cell_weight_1<dim>);
   tr.repartition();
 
-  tr.signals.cell_weight.disconnect_all_slots();
+  tr.signals.weight.disconnect_all_slots();
 
-  tr.signals.cell_weight.connect(std::bind(&cell_weight_2<dim>,
-                                           std::placeholders::_1,
-                                           std::placeholders::_2));
+  tr.signals.weight.connect(&cell_weight_2<dim>);
   tr.repartition();
 
   const auto n_locally_owned_active_cells_per_processor =
index c51f11a9c98aea3e912960dd0cfc379cbb5e5440..3e3ed3678eea84ef3aee3eb072c9dfcb87b0309d 100644 (file)
@@ -64,12 +64,10 @@ test()
 
   // repartition the mesh as described above, first in some arbitrary
   // way, and then with no weights
-  tr.signals.cell_weight.connect(std::bind(&cell_weight_1<dim>,
-                                           std::placeholders::_1,
-                                           std::placeholders::_2));
+  tr.signals.weight.connect(&cell_weight_1<dim>);
   tr.repartition();
 
-  tr.signals.cell_weight.disconnect_all_slots();
+  tr.signals.weight.disconnect_all_slots();
   tr.repartition();
 
   const auto n_locally_owned_active_cells_per_processor =
index e2e0700a662a43a0a9a54c1ad6e25ded97ba2915..b55f488d0f7a4a458df502ae23bac79671c2e3aa 100644 (file)
@@ -55,7 +55,7 @@ test()
 
   GridGenerator::subdivided_hyper_cube(tr, 16);
 
-  tr.signals.cell_weight.connect(&cell_weight<dim>);
+  tr.signals.weight.connect(&cell_weight<dim>);
   tr.refine_global(1);
 
   const auto n_locally_owned_active_cells_per_processor =
index 54082fd8886762a55aaa3e47a43424d4577be44b..1313a8589e568d8b4f175b74adbaf28371cc4466 100644 (file)
@@ -63,7 +63,7 @@ test()
 
   tr.refine_global(1);
 
-  tr.signals.cell_weight.connect(&cell_weight<dim>);
+  tr.signals.weight.connect(&cell_weight<dim>);
 
   tr.repartition();
 
index 65c505e4f9715a2f7b488a4bccda670ae295bf7c..9305ba654dfcadbc8863db561a9a301bad417dec 100644 (file)
@@ -57,7 +57,7 @@ test()
 
   GridGenerator::subdivided_hyper_cube(tr, 16);
 
-  tr.signals.cell_weight.connect(&cell_weight<dim>);
+  tr.signals.weight.connect(&cell_weight<dim>);
 
   tr.refine_global(1);
 
index 823eaec3e3d9e4a8ee2388b39e530f52130c00a5..c0d027b0e5551fa1d325a68e700b1cb63f87de74 100644 (file)
@@ -77,7 +77,7 @@ test()
 
 
   // repartition the mesh; attach different weights to all cells
-  tr.signals.cell_weight.connect(&cell_weight<dim>);
+  tr.signals.weight.connect(&cell_weight<dim>);
 
   tr.repartition();
 
index ded8471d68c74d76a3af0de47fa06cc9f864dc5d..f4865b0fa54946d266364ccda2f09ebee3f6b6bc 100644 (file)
@@ -77,7 +77,7 @@ test()
 
   // repartition the mesh; attach different weights to all cells
   n_global_active_cells = tr.n_global_active_cells();
-  tr.signals.cell_weight.connect(&cell_weight<dim>);
+  tr.signals.weight.connect(&cell_weight<dim>);
   tr.repartition();
 
   const auto n_locally_owned_active_cells_per_processor =
index 544f59feed2b740d962103ecef55db680eb8bd68..2767eef934166de61ebe44a2d0dc3d22a70d4983 100644 (file)
@@ -249,7 +249,7 @@ namespace Step68
   // return value of this function (representing "work for this cell") is
   // calculated based on the number of particles in the current cell.
   // The function is
-  // connected to the cell_weight() signal inside the triangulation, and will be
+  // connected to the `weight` signal inside the triangulation, and will be
   // called once per cell, whenever the triangulation repartitions the domain
   // between ranks (the connection is created inside the
   // generate_particles() function of this class).
@@ -327,7 +327,7 @@ namespace Step68
     // be created once, so we might as well have set it up in the constructor
     // of this class, but for the purpose of this example we want to group the
     // particle related instructions.
-    background_triangulation.signals.cell_weight.connect(
+    background_triangulation.signals.weight.connect(
       [&](
         const typename parallel::distributed::Triangulation<dim>::cell_iterator
           &cell,

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.