]> https://gitweb.dealii.org/ - dealii.git/commitdiff
cell_weight: remove initial weight from distributed applications.
authorMarc Fehling <mafehling.git@gmail.com>
Sat, 12 Mar 2022 03:13:01 +0000 (20:13 -0700)
committerMarc Fehling <mafehling.git@gmail.com>
Sat, 12 Mar 2022 03:28:05 +0000 (20:28 -0700)
doc/news/changes/incompatibilities/20220225Fehling [new file with mode: 0644]
examples/step-68/step-68.cc
examples/step-75/step-75.cc
include/deal.II/distributed/cell_weights.h
include/deal.II/grid/tria.h
source/distributed/tria.cc
tests/particles/step-68.cc

diff --git a/doc/news/changes/incompatibilities/20220225Fehling b/doc/news/changes/incompatibilities/20220225Fehling
new file mode 100644 (file)
index 0000000..4ef1b4f
--- /dev/null
@@ -0,0 +1,12 @@
+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:
+@code{.cc}
+triangulation.signals.cell_weight.connect(
+  [](const typename parallel::distributed::Triangulation<dim>::cell_iterator &,
+     const typename parallel::distributed::Triangulation<dim>::CellStatus)
+    -> unsigned int { return 1000; });
+@endcode
+(Marc Fehling, 2022/02/25)
index 2d7bacd367087b17224a1bbd75ad2436d1110202..a79f2f106f9e1e141020368c9122dc572bed2476 100644 (file)
@@ -1,6 +1,6 @@
 /* ---------------------------------------------------------------------
  *
- * Copyright (C) 2020 - 2021 by the deal.II authors
+ * Copyright (C) 2020 - 2022 by the deal.II authors
  *
  * This file is part of the deal.II library.
  *
@@ -26,7 +26,6 @@
 #include <deal.II/base/parameter_acceptor.h>
 #include <deal.II/base/timer.h>
 
-#include <deal.II/distributed/cell_weights.h>
 #include <deal.II/distributed/solution_transfer.h>
 #include <deal.II/distributed/tria.h>
 
@@ -334,21 +333,18 @@ namespace Step68
     const typename parallel::distributed::Triangulation<dim>::CellStatus status)
     const
   {
-    // We do not assign any weight to cells we do not own (i.e., artificial
-    // or ghost cells)
-    if (!cell->is_locally_owned())
-      return 0;
-
-    // This determines how important particle work is compared to cell
-    // work (by default every cell has a weight of 1000).
-    // We set the weight per particle much higher to indicate that
-    // the particle load is the only one that is important to distribute the
-    // cells in this example. The optimal value of this number depends on the
-    // application and can range from 0 (cheap particle operations,
-    // expensive cell operations) to much larger than 1000 (expensive
-    // particle operations, cheap cell operations, like presumed in this
-    // example).
-    const unsigned int particle_weight = 10000;
+    // First, we introduce a base weight that will be assigned to every cell.
+    const unsigned int base_weight = 1;
+
+    // The following variable then determines how important particle work is
+    // compared to cell work. We set the weight per particle much higher to
+    // indicate that the particle load is the only one that is important to
+    // distribute the cells in this example. The optimal value of this number
+    // depends on the application and can range from 0 (cheap particle
+    // operations, expensive cell operations) to much larger than the base
+    // weight of 1 (expensive particle operations, cheap cell operations, like
+    // presumed in this example).
+    const unsigned int particle_weight = 10;
 
     // This example does not use adaptive refinement, therefore every cell
     // should have the status `CELL_PERSIST`. However this function can also
index 54874562214cdf516ab7bad391cbbae9f484de3e..0847e66ebe4ebc788e8591871c70952592afcea5 100644 (file)
@@ -1,6 +1,6 @@
 /* ---------------------------------------------------------------------
  *
- * Copyright (C) 2021 by the deal.II authors
+ * Copyright (C) 2021 - 2022 by the deal.II authors
  *
  * This file is part of the deal.II library.
  *
@@ -191,7 +191,7 @@ namespace Step75
     double p_refine_fraction  = 0.9;
     double p_coarsen_fraction = 0.9;
 
-    double weighting_factor   = 1e6;
+    double weighting_factor   = 1.;
     double weighting_exponent = 1.;
   };
 
@@ -980,16 +980,13 @@ namespace Step75
     // for hp-adaptation and right before repartitioning for load balancing is
     // about to happen. Functions can be registered that will attach weights in
     // the form that $a (n_\text{dofs})^b$ with a provided pair of parameters
-    // $(a,b)$. We register such a function in the following. Every cell will be
-    // charged with a constant weight at creation, which is a value of 1000 (see
-    // Triangulation::Signals::cell_weight).
+    // $(a,b)$. We register such a function in the following.
     //
     // For load balancing, efficient solvers like the one we use should scale
-    // linearly with the number of degrees of freedom owned. Further, to
-    // increase the impact of the weights we would like to attach, make sure
-    // that the individual weight will exceed this base weight by orders of
-    // magnitude. We set the parameters for cell weighting correspondingly: A
-    // large weighting factor of $10^6$ and an exponent of $1$.
+    // linearly with the number of degrees of freedom owned. We set the
+    // parameters for cell weighting correspondingly: A weighting factor of $1$
+    // and an exponent of $1$ (see the definitions of the `weighting_factor` and
+    // `weighting_exponent` above).
     cell_weights = std::make_unique<parallel::CellWeights<dim>>(
       dof_handler,
       parallel::CellWeights<dim>::ndofs_weighting(
index 7804dc0a6e6ab749738b5151d78fbd3cbdcac0cc..8fce0349f3cf18ce75cb4862035d97caa40b7494 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2018 - 2021 by the deal.II authors
+// Copyright (C) 2018 - 2022 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -53,12 +53,11 @@ namespace parallel
    * The connected weighting function may be changed anytime using the
    * CellWeights::reinit() function. The following code snippet demonstrates how
    * to achieve each cell being weighted by its current number of degrees of
-   * freedom. We chose a factor of `1000` that corresponds to the initial weight
-   * each cell is assigned to upon creation.
+   * freedom.
    * @code
    * parallel::CellWeights<dim, spacedim> cell_weights(
    *   hp_dof_handler,
-   *   parallel::CellWeights<dim, spacedim>::ndofs_weighting({1000, 1}));
+   *   parallel::CellWeights<dim, spacedim>::ndofs_weighting({1, 1}));
    * @endcode
    *
    * On the other hand, you are also able to take care of handling the signal
@@ -69,8 +68,7 @@ namespace parallel
    *   hp_dof_handler.get_triangulation().signals.cell_weight.connect(
    *     parallel::CellWeights<dim, spacedim>::make_weighting_callback(
    *       hp_dof_handler,
-   *       parallel::CellWeights<dim, spacedim>::ndofs_weighting(
-   *         {1000, 1}));
+   *       parallel::CellWeights<dim, spacedim>::ndofs_weighting({1, 1}));
    * @endcode
    *
    * The use of this class is demonstrated in step-75.
@@ -85,6 +83,17 @@ namespace parallel
    * the weight_callback() function. Use CellWeights::reinit() to deregister the
    * weighting function on the old Triangulation and connect it to the new one.
    *
+   * @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,
+   * 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
+   * parallel::distributed::Triangulation::execute_coarsening_and_refinement(),
+   * parallel::distributed::Triangulation::refine_global(),
+   * parallel::distributed::Triangulation::repartition(), or
+   * GridTools::partition_triangulation() for the very first time.
+   *
    * @ingroup distributed
    */
   template <int dim, int spacedim = dim>
@@ -156,7 +165,7 @@ namespace parallel
      * Choose a constant weight @p factor on each cell.
      */
     static WeightingFunction
-    constant_weighting(const unsigned int factor = 1000);
+    constant_weighting(const unsigned int factor = 1);
 
     /**
      * The pair of floating point numbers $(a,b)$ provided via
index 0331b96848b0961f0118fdb2a9a948219d942b7a..7badd095d0723f6d67932feb1b6deb828e322e48 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 1998 - 2021 by the deal.II authors
+// Copyright (C) 1998 - 2022 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -2173,8 +2173,7 @@ public:
 
     /**
      * This signal is triggered for each cell during every automatic or manual
-     * repartitioning. This signal is somewhat special in that it is only
-     * triggered for distributed parallel calculations and only if functions
+     * repartitioning. This signal will only be triggered if functions
      * are connected to it. It is intended to allow a weighted repartitioning
      * of the domain to balance the computational load across processes in a
      * different way than balancing the number of cells. Any connected
@@ -2183,18 +2182,23 @@ public:
      * coarsened or left untouched (see the documentation of the CellStatus
      * enum for more information). The function is expected to return an
      * unsigned integer, which is interpreted as the additional computational
-     * load of this cell. If this cell is going to be coarsened, the signal is
-     * called for the parent cell and you need to provide the weight of the
-     * future parent cell. If this cell is going to be refined the function
-     * should return a weight, which will be equally assigned to every future
-     * child cell of the current cell. 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. If several
-     * functions are connected to this signal, their return values will be
-     * summed to calculate the final weight.
+     * load of this cell.
      *
-     * This function is used in step-68.
+     * In serial and parallel shared applications, partitioning happens after
+     * refinement. So all cells will have the `CELL_PERSIST` status.
+     *
+     * In parallel distributed applications, partitioning happens during
+     * refinement. If this cell is going to be coarsened, the signal is called
+     * for the parent cell and you need to provide the weight of the future
+     * parent cell. If this cell is going to be refined, the function is called
+     * on all children and you should ideally return the same weight for all
+     * children.
+     *
+     * If several functions are connected to this signal, their return values
+     * will be summed to calculate the final weight via `CellWeightSum`.
+     *
+     * This function is used in step-68 and implicitely in step-75 using the
+     * parallel::CellWeights class.
      */
     boost::signals2::signal<unsigned int(const cell_iterator &,
                                          const CellStatus),
index 670741955620bc7da7a38661d653e29365182c88..c879f53d9aea0b1151f9a79b466ecc8c05fe4146 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2008 - 2021 by the deal.II authors
+// Copyright (C) 2008 - 2022 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -4038,47 +4038,7 @@ namespace parallel
           const auto &cell_it     = cell_rel.first;
           const auto &cell_status = cell_rel.second;
 
-          switch (cell_status)
-            {
-              case parallel::distributed::Triangulation<dim,
-                                                        spacedim>::CELL_PERSIST:
-                weights.push_back(1000);
-                weights.back() += this->signals.cell_weight(
-                  cell_it,
-                  parallel::distributed::Triangulation<dim,
-                                                       spacedim>::CELL_PERSIST);
-                break;
-
-              case parallel::distributed::Triangulation<dim,
-                                                        spacedim>::CELL_REFINE:
-              case parallel::distributed::Triangulation<dim,
-                                                        spacedim>::CELL_INVALID:
-                {
-                  // calculate weight of parent cell
-                  unsigned int parent_weight = 1000;
-                  parent_weight += this->signals.cell_weight(
-                    cell_it,
-                    parallel::distributed::Triangulation<dim, spacedim>::
-                      CELL_REFINE);
-                  // assign the weight of the parent cell equally to all
-                  // children
-                  weights.push_back(parent_weight);
-                  break;
-                }
-
-              case parallel::distributed::Triangulation<dim,
-                                                        spacedim>::CELL_COARSEN:
-                weights.push_back(1000);
-                weights.back() += this->signals.cell_weight(
-                  cell_it,
-                  parallel::distributed::Triangulation<dim,
-                                                       spacedim>::CELL_COARSEN);
-                break;
-
-              default:
-                Assert(false, ExcInternalError());
-                break;
-            }
+          weights.push_back(this->signals.cell_weight(cell_it, cell_status));
         }
 
       return weights;
index e4d6bbd44062410217a996004698f173b1641fa7..72713fbc65540d58a9603b5fab1a2316bb54800c 100644 (file)
@@ -1,6 +1,6 @@
 /* ---------------------------------------------------------------------
  *
- * Copyright (C) 2020 - 2021 by the deal.II authors
+ * Copyright (C) 2020 - 2022 by the deal.II authors
  *
  * This file is part of the deal.II library.
  *
@@ -26,7 +26,6 @@
 #include <deal.II/base/parameter_acceptor.h>
 #include <deal.II/base/timer.h>
 
-#include <deal.II/distributed/cell_weights.h>
 #include <deal.II/distributed/solution_transfer.h>
 #include <deal.II/distributed/tria.h>
 
@@ -262,21 +261,18 @@ namespace Step68
     const typename parallel::distributed::Triangulation<dim>::CellStatus status)
     const
   {
-    // We do not assign any weight to cells we do not own (i.e., artificial
-    // or ghost cells)
-    if (!cell->is_locally_owned())
-      return 0;
-
-    // This determines how important particle work is compared to cell
-    // work (by default every cell has a weight of 1000).
-    // We set the weight per particle much higher to indicate that
-    // the particle load is the only one that is important to distribute the
-    // cells in this example. The optimal value of this number depends on the
-    // application and can range from 0 (cheap particle operations,
-    // expensive cell operations) to much larger than 1000 (expensive
-    // particle operations, cheap cell operations, like presumed in this
-    // example).
-    const unsigned int particle_weight = 10000;
+    // First, we introduce a base weight that will be assigned to every cell.
+    const unsigned int base_weight = 1;
+
+    // The following variable then determines how important particle work is
+    // compared to cell work. We set the weight per particle much higher to
+    // indicate that the particle load is the only one that is important to
+    // distribute the cells in this example. The optimal value of this number
+    // depends on the application and can range from 0 (cheap particle
+    // operations, expensive cell operations) to much larger than the base
+    // weight of 1 (expensive particle operations, cheap cell operations, like
+    // presumed in this example).
+    const unsigned int particle_weight = 10;
 
     // This example does not use adaptive refinement, therefore every cell
     // should have the status `CELL_PERSIST`. However this function can also

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.