]> https://gitweb.dealii.org/ - dealii.git/commitdiff
fix 64bit overflow in refine_and_coarsen_fixed_number 9568/head
authorTimo Heister <timo.heister@gmail.com>
Wed, 26 Feb 2020 13:13:52 +0000 (08:13 -0500)
committerTimo Heister <timo.heister@gmail.com>
Wed, 26 Feb 2020 13:13:52 +0000 (08:13 -0500)
The argument max_n_cells (that probably nobody uses) was declared as an
unsigned int, meaning that we can not refine over 4 billion cells with
this function. This is now fixed.

include/deal.II/distributed/grid_refinement.h
include/deal.II/grid/grid_refinement.h
source/distributed/grid_refinement.cc
source/distributed/grid_refinement.inst.in
source/grid/grid_refinement.cc
source/grid/grid_refinement.inst.in

index 4f96a88bc26626c2327aa369e723d5e4032267f1..6b4358ec89c18f7878405e6261b0eeb4bc343b49 100644 (file)
@@ -131,10 +131,10 @@ namespace parallel
       refine_and_coarsen_fixed_number(
         parallel::distributed::Triangulation<dim, spacedim> &tria,
         const dealii::Vector<Number> &                       criteria,
-        const double       top_fraction_of_cells,
-        const double       bottom_fraction_of_cells,
-        const unsigned int max_n_cells =
-          std::numeric_limits<unsigned int>::max());
+        const double                  top_fraction_of_cells,
+        const double                  bottom_fraction_of_cells,
+        const types::global_dof_index max_n_cells =
+          std::numeric_limits<types::global_dof_index>::max());
 
       /**
        * Like dealii::GridRefinement::refine_and_coarsen_fixed_fraction, but
index 52a6fdaf7e8c8b6c67456cfb34c7e22b0067c171..3ec8e038627b2adcc756dfc48ecdb1cc175d5d97 100644 (file)
@@ -87,10 +87,10 @@ namespace GridRefinement
   template <int dim>
   std::pair<double, double>
   adjust_refine_and_coarsen_number_fraction(
-    const unsigned int current_n_cells,
-    const unsigned int max_n_cells,
-    const double       top_fraction_of_cells,
-    const double       bottom_fraction_of_cells);
+    const types::global_dof_index current_n_cells,
+    const types::global_dof_index max_n_cells,
+    const double                  top_fraction_of_cells,
+    const double                  bottom_fraction_of_cells);
 
   /**
    * This function provides a strategy to mark cells for refinement and
index a1625db018f19dc755d2990c1f06b116f38be887..fc45b9f954552bebfa4f0d1f45566e7326151e24 100644 (file)
@@ -434,9 +434,9 @@ namespace parallel
       refine_and_coarsen_fixed_number(
         parallel::distributed::Triangulation<dim, spacedim> &tria,
         const dealii::Vector<Number> &                       criteria,
-        const double       top_fraction_of_cells,
-        const double       bottom_fraction_of_cells,
-        const unsigned int max_n_cells)
+        const double                  top_fraction_of_cells,
+        const double                  bottom_fraction_of_cells,
+        const types::global_dof_index max_n_cells)
       {
         Assert(criteria.size() == tria.n_active_cells(),
                ExcDimensionMismatch(criteria.size(), tria.n_active_cells()));
index 3e4477074f2e0d7f06b615a6cd8777bbb87a5d58..688294715d0850e01ae91faaade780e98031927f 100644 (file)
@@ -69,7 +69,7 @@ for (S : REAL_SCALARS; deal_II_dimension : DIMENSIONS)
             const dealii::Vector<S> &,
             const double,
             const double,
-            const unsigned int);
+            const types::global_dof_index);
 
           template void
           refine_and_coarsen_fixed_fraction<deal_II_dimension,
@@ -100,7 +100,7 @@ for (S : REAL_SCALARS; deal_II_dimension : DIMENSIONS)
             const dealii::Vector<S> &,
             const double,
             const double,
-            const unsigned int);
+            const types::global_dof_index);
 
           template void refine_and_coarsen_fixed_fraction<deal_II_dimension - 1,
                                                           S,
index f5581c3af2240cce7eca8d73ebdb3c7097672a4a..29dc3bdf706fefc9d534d1d98daad8bcf71e1d4b 100644 (file)
@@ -110,10 +110,10 @@ GridRefinement::coarsen(Triangulation<dim, spacedim> &tria,
 template <int dim>
 std::pair<double, double>
 GridRefinement::adjust_refine_and_coarsen_number_fraction(
-  const unsigned int current_n_cells,
-  const unsigned int max_n_cells,
-  const double       top_fraction,
-  const double       bottom_fraction)
+  const types::global_dof_index current_n_cells,
+  const types::global_dof_index max_n_cells,
+  const double                  top_fraction,
+  const double                  bottom_fraction)
 {
   Assert(top_fraction >= 0, ExcInvalidParameterValue());
   Assert(top_fraction <= 1, ExcInvalidParameterValue());
@@ -166,7 +166,7 @@ GridRefinement::adjust_refine_and_coarsen_number_fraction(
   // again, this is true for isotropically
   // refined cells. we take this as an
   // approximation of a mixed refinement.
-  else if (static_cast<unsigned int>(
+  else if (static_cast<types::global_dof_index>(
              current_n_cells + refine_cells * cell_increase_on_refine -
              coarsen_cells * cell_decrease_on_coarsen) > max_n_cells)
     {
index e173c988cfe6f8591248357dee3484aef1e2f66f..8bda2e72e1a70e2d402f252c3fcbe158928e8b3e 100644 (file)
@@ -97,7 +97,10 @@ for (S : REAL_SCALARS; deal_II_dimension : DIMENSIONS)
 
 for (deal_II_dimension : DIMENSIONS)
   {
-    template std::pair<double, double> GridRefinement::
-      adjust_refine_and_coarsen_number_fraction<deal_II_dimension>(
-        const unsigned int, const unsigned int, const double, const double);
+    template std::pair<double, double>
+    GridRefinement::adjust_refine_and_coarsen_number_fraction<
+      deal_II_dimension>(const types::global_dof_index,
+                         const types::global_dof_index,
+                         const double,
+                         const double);
   }

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.