]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Remove custom repartitioning from load. 12596/head
authorMarc Fehling <mafehling.git@gmail.com>
Sat, 24 Jul 2021 00:10:09 +0000 (18:10 -0600)
committerMarc Fehling <mafehling.git@gmail.com>
Wed, 28 Jul 2021 19:52:17 +0000 (13:52 -0600)
doc/news/changes/incompatibilities/20210723Fehling [new file with mode: 0644]
include/deal.II/distributed/tria.h
source/distributed/tria.cc

diff --git a/doc/news/changes/incompatibilities/20210723Fehling b/doc/news/changes/incompatibilities/20210723Fehling
new file mode 100644 (file)
index 0000000..276b0a3
--- /dev/null
@@ -0,0 +1,6 @@
+Changed: Custom repartitioning via Triangulation::Signals::cell_weight
+is no longer part of parallel::distributed::Triangulation::load(). If
+desired, call parallel::distributed::Triangulation::repartition()
+manually after load().
+<br>
+(Marc Fehling, 2021/07/23)
index 65860b9fcc976754332afa71bf21516cdacca0ba..f64343ccc1cf9c29aaf7a9c9ac6675fb0fb31970 100644 (file)
@@ -594,17 +594,22 @@ namespace parallel
        * You do not need to load with the same number of MPI processes that
        * you saved with. Rather, if a mesh is loaded with a different number
        * of MPI processes than used at the time of saving, the mesh is
-       * repartitioned appropriately. Cell-based data that was saved with
+       * repartitioned that the number of cells is balanced among all processes.
+       * Individual repartitioning, e.g., based on the number of dofs or
+       * particles per cell, needs to be invoked manually by calling
+       * repartition() afterwards.
+       *
+       * Cell-based data that was saved with
        * DistributedTriangulationBase::DataTransfer::register_data_attach() can
        * be read in with
        * DistributedTriangulationBase::DataTransfer::notify_ready_to_unpack()
        * after calling load().
        *
-       * If you use p4est version > 0.3.4.2 the @p autopartition flag tells
-       * p4est to ignore the partitioning that the triangulation had when it
-       * was saved and make it uniform upon loading. If @p autopartition is
-       * set to false, the triangulation is only repartitioned if needed (i.e.
-       * if a different number of MPI processes is encountered).
+       * The @p autopartition flag tells p4est to ignore the partitioning that
+       * the triangulation had when it was saved and make it uniform upon
+       * loading. If @p autopartition is set to true, the triangulation will
+       * always be repartitioned. If set to false, it is only repartitioned if
+       * needed (i.e., if a different number of MPI processes is encountered).
        */
       virtual void
       load(const std::string &filename,
index e7f0271692ef504e763a61330345cfd335120c34..384da4c46531ff45867082f4c1f9d7a6410da2f0 100644 (file)
@@ -1534,8 +1534,8 @@ namespace parallel
         filename.c_str(),
         this->mpi_communicator,
         0,
-        false,
-        autopartition,
+        0,
+        static_cast<int>(autopartition),
         0,
         this,
         &connectivity);
@@ -1547,49 +1547,10 @@ namespace parallel
           // number of CPUs and so everything works without this call, but
           // this command changes the distribution for some reason, so we
           // will leave it in here.
-          if (this->signals.cell_weight.num_slots() == 0)
-            {
-              // no cell weights given -- call p4est's 'partition' without a
-              // callback for cell weights
-              dealii::internal::p4est::functions<dim>::partition(
-                parallel_forest,
-                /* prepare coarsening */ 1,
-                /* weight_callback */ nullptr);
-            }
-          else
-            {
-              // first, align both dealii and p4est meshes, which is a
-              // requirement of our internal cell weights algorithm
-              try
-                {
-                  copy_local_forest_to_triangulation();
-                }
-              catch (const typename Triangulation<dim>::DistortedCellList &)
-                {
-                  // the underlying triangulation should not be checking for
-                  // distorted cells
-                  Assert(false, ExcInternalError());
-                }
-
-              // get cell weights for a weighted repartitioning.
-              const std::vector<unsigned int> cell_weights = get_cell_weights();
-
-              PartitionWeights<dim, spacedim> partition_weights(cell_weights);
-
-              // attach (temporarily) a pointer to the cell weights through
-              // p4est's user_pointer object
-              Assert(parallel_forest->user_pointer == this, ExcInternalError());
-              parallel_forest->user_pointer = &partition_weights;
-
-              dealii::internal::p4est::functions<dim>::partition(
-                parallel_forest,
-                /* prepare coarsening */ 1,
-                /* weight_callback */
-                &PartitionWeights<dim, spacedim>::cell_weight);
-
-              // reset the user pointer to its previous state
-              parallel_forest->user_pointer = this;
-            }
+          dealii::internal::p4est::functions<dim>::partition(
+            parallel_forest,
+            /* prepare coarsening */ 1,
+            /* weight_callback */ nullptr);
         }
 
       try

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.