]> https://gitweb.dealii.org/ - dealii.git/commitdiff
cell_weight: check that global sum of weights is larger than 0.
authorMarc Fehling <mafehling.git@gmail.com>
Sat, 12 Mar 2022 03:14:09 +0000 (20:14 -0700)
committerMarc Fehling <mafehling.git@gmail.com>
Sat, 12 Mar 2022 03:28:05 +0000 (20:28 -0700)
source/distributed/tria.cc
source/grid/grid_tools.cc

index c879f53d9aea0b1151f9a79b466ecc8c05fe4146..617c9ec3b62f79834fccef5f00ad598fc495a200 100644 (file)
@@ -3361,6 +3361,20 @@ namespace parallel
               // get cell weights for a weighted repartitioning.
               const std::vector<unsigned int> cell_weights = get_cell_weights();
 
+#  ifdef DEBUG
+              // verify that the global sum of weights is larger than 0
+              const auto local_weights_sum =
+                std::accumulate(cell_weights.begin(),
+                                cell_weights.end(),
+                                std::uint64_t(0));
+              const auto global_weights_sum =
+                Utilities::MPI::sum(local_weights_sum, this->mpi_communicator);
+              Assert(global_weights_sum > 0,
+                     ExcMessage(
+                       "The global sum of weights over all active cells "
+                       "is zero. Please verify how you generate weights."));
+#  endif
+
               PartitionWeights<dim, spacedim> partition_weights(cell_weights);
 
               // attach (temporarily) a pointer to the cell weights through
@@ -3529,6 +3543,19 @@ namespace parallel
           // get cell weights for a weighted repartitioning.
           const std::vector<unsigned int> cell_weights = get_cell_weights();
 
+#  ifdef DEBUG
+          // verify that the global sum of weights is larger than 0
+          const auto local_weights_sum = std::accumulate(cell_weights.begin(),
+                                                         cell_weights.end(),
+                                                         std::uint64_t(0));
+          const auto global_weights_sum =
+            Utilities::MPI::sum(local_weights_sum, this->mpi_communicator);
+          Assert(global_weights_sum > 0,
+                 ExcMessage(
+                   "The global sum of weights over all active cells "
+                   "is zero. Please verify how you generate weights."));
+#  endif
+
           PartitionWeights<dim, spacedim> partition_weights(cell_weights);
 
           // attach (temporarily) a pointer to the cell weights through
index a88f8047e0de23106efa00663c3d4c0f1522ef05..a9566e6d7aac9923c085dfe3b1c55fb7235a2a8d 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2001 - 2021 by the deal.II authors
+// Copyright (C) 2001 - 2022 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -3832,7 +3832,7 @@ namespace GridTools
         // If this is a parallel triangulation, we then need to also
         // get the weights for all other cells. We have asserted above
         // that this function can't be used for
-        // parallel::distribute::Triangulation objects, so the only
+        // parallel::distributed::Triangulation objects, so the only
         // ones we have to worry about here are
         // parallel::shared::Triangulation
         if (const auto shared_tria =
@@ -3841,6 +3841,16 @@ namespace GridTools
           Utilities::MPI::sum(cell_weights,
                               shared_tria->get_communicator(),
                               cell_weights);
+
+#ifdef DEBUG
+        // verify that the global sum of weights is larger than 0
+        const auto global_cell_sum = std::accumulate(cell_weights.begin(),
+                                                     cell_weights.end(),
+                                                     std::uint64_t(0));
+        Assert(global_cell_sum > 0,
+               ExcMessage("The global sum of weights over all active cells "
+                          "is zero. Please verify how you generate weights."));
+#endif
       }
 
     // Call the other more general function
@@ -3935,6 +3945,16 @@ namespace GridTools
           Utilities::MPI::sum(cell_weights,
                               shared_tria->get_communicator(),
                               cell_weights);
+
+#ifdef DEBUG
+        // verify that the global sum of weights is larger than 0
+        const auto global_cell_sum = std::accumulate(cell_weights.begin(),
+                                                     cell_weights.end(),
+                                                     std::uint64_t(0));
+        Assert(global_cell_sum > 0,
+               ExcMessage("The global sum of weights over all active cells "
+                          "is zero. Please verify how you generate weights."));
+#endif
       }
 
     // Call the other more general function

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.