]> https://gitweb.dealii.org/ - dealii.git/commitdiff
minor MGRepartiniong edits 12928/head
authorTimo Heister <timo.heister@gmail.com>
Mon, 8 Nov 2021 23:41:37 +0000 (18:41 -0500)
committerTimo Heister <timo.heister@gmail.com>
Mon, 8 Nov 2021 23:41:37 +0000 (18:41 -0500)
include/deal.II/distributed/repartitioning_policy_tools.h
include/deal.II/multigrid/mg_transfer_global_coarsening.h
source/distributed/repartitioning_policy_tools.cc

index ff283c56fdbf978f5db5dc18350dbcec6de78970..81d412c45cccd7d3c90eff57909cfbb268f17c0a 100644 (file)
@@ -39,7 +39,11 @@ DEAL_II_NAMESPACE_OPEN
 namespace RepartitioningPolicyTools
 {
   /**
-   * A base class of a repartitioning policy.
+   * The base class for repartitioning policies.
+   *
+   * Used in
+   * MGTransferGlobalCoarseningTools::create_geometric_coarsening_sequence().
+   * See the description of RepartitioningPolicyTools for more information.
    */
   template <int dim, int spacedim = dim>
   class Base
@@ -121,8 +125,10 @@ namespace RepartitioningPolicyTools
   };
 
   /**
-   * A policy that allows to specify a minimal number of cells per process. If
-   * a threshold is reached, processes might be left without cells.
+   * A policy that allows to specify a minimal number of cells per
+   * process. If a threshold is reached, processes might be left
+   * without cells. The cells will be distributed evenly among the
+   * remaining processes.
    */
   template <int dim, int spacedim = dim>
   class MinimalGranularityPolicy : public Base<dim, spacedim>
index a02f83a821737753b1bdeca31c6ed97c61f72000..e1a42451cb4c808ac47beaa2a6cae30711546866 100644 (file)
@@ -135,8 +135,8 @@ namespace MGTransferGlobalCoarseningTools
 
   /**
    * Similar to the above function but taking in a constant version of
-   * @p tria and as a consequence not allowing to directly use it for
-   * coarsening, requiring that internally a temporal copy is created.
+   * @p tria. As a consequence, it can not be used for coarsening directly,
+   * so a temporary copy will be created internally.
    */
   template <int dim, int spacedim>
   std::vector<std::shared_ptr<const Triangulation<dim, spacedim>>>
index 9b203f93c2abeb9d067b12dc76aac4051d64573f..be1947389837c3fca519c6a2f8fde656e9dba118 100644 (file)
@@ -219,12 +219,14 @@ namespace RepartitioningPolicyTools
     // a repartitioning kicks in with the aim that all processes that own
     // cells have at least the specified number of cells
 
-    const unsigned int n_global_active_cells = tria_in.n_global_active_cells();
+    const types::global_cell_index n_global_active_cells =
+      tria_in.n_global_active_cells();
 
     const unsigned int n_partitions =
       std::max<unsigned int>(1,
-                             std::min(n_global_active_cells / n_min_cells,
-                                      Utilities::MPI::n_mpi_processes(comm)));
+                             std::min<types::global_cell_index>(
+                               n_global_active_cells / n_min_cells,
+                               Utilities::MPI::n_mpi_processes(comm)));
 
     const unsigned int min_cells = n_global_active_cells / n_partitions;
 
@@ -322,6 +324,7 @@ namespace RepartitioningPolicyTools
                      Utilities::MPI::internal::mpi_type_id(&total_weight),
                      n_subdomains - 1,
                      mpi_communicator);
+    AssertThrowMPI(ierr);
 
     // setup partition
     LinearAlgebra::distributed::Vector<double> partition(partitioner);

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.