]> https://gitweb.dealii.org/ - dealii.git/commitdiff
use global_cell_index in many places
authorTimo Heister <timo.heister@gmail.com>
Sun, 1 Mar 2020 19:58:34 +0000 (14:58 -0500)
committerTimo Heister <timo.heister@gmail.com>
Wed, 4 Mar 2020 15:54:01 +0000 (10:54 -0500)
include/deal.II/base/parsed_convergence_table.h
include/deal.II/distributed/grid_refinement.h
include/deal.II/distributed/tria_base.h
include/deal.II/grid/grid_refinement.h
include/deal.II/grid/tria.h
source/distributed/grid_refinement.cc
source/distributed/grid_refinement.inst.in
source/distributed/tria_base.cc
source/grid/grid_refinement.cc
source/grid/grid_refinement.inst.in
source/grid/tria.cc

index 76659f449ad9fbc2f7168718cf6952a541450e9e..41033011369f81068021bb3e133fbb95add18a9d 100644 (file)
@@ -566,7 +566,7 @@ ParsedConvergenceTable::error_from_exact(
       AssertDimension(exact.n_components, n_components);
       AssertDimension(dh.get_fe().n_components(), n_components);
 
-      const unsigned int n_active_cells =
+      const types::global_cell_index n_active_cells =
         dh.get_triangulation().n_global_active_cells();
       const unsigned int n_dofs = dh.n_dofs();
 
index 6b4358ec89c18f7878405e6261b0eeb4bc343b49..22fcd228652d05764dec21b6a17064a6b521aa97 100644 (file)
@@ -57,7 +57,7 @@ namespace internal
           number
           compute_threshold(const dealii::Vector<number> &   criteria,
                             const std::pair<double, double> &global_min_and_max,
-                            const types::global_dof_index    n_target_cells,
+                            const types::global_cell_index   n_target_cells,
                             MPI_Comm                         mpi_communicator);
         } // namespace RefineAndCoarsenFixedNumber
 
@@ -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 types::global_dof_index max_n_cells =
-          std::numeric_limits<types::global_dof_index>::max());
+        const double                   top_fraction_of_cells,
+        const double                   bottom_fraction_of_cells,
+        const types::global_cell_index max_n_cells =
+          std::numeric_limits<types::global_cell_index>::max());
 
       /**
        * Like dealii::GridRefinement::refine_and_coarsen_fixed_fraction, but
index 0e6507defc314aecad343ab9cf9f65be0ef248d7..458c72699c046fd95207806e8b9f6f2fdcc5c3a4 100644 (file)
@@ -145,7 +145,7 @@ namespace parallel
      * by each processor. This equals the overall number of active cells in
      * the triangulation.
      */
-    virtual types::global_dof_index
+    virtual types::global_cell_index
     n_global_active_cells() const override;
 
     /**
@@ -244,7 +244,7 @@ namespace parallel
        * The total number of active cells (sum of @p
        * n_locally_owned_active_cells).
        */
-      types::global_dof_index n_global_active_cells;
+      types::global_cell_index n_global_active_cells;
       /**
        * The global number of levels computed as the maximum number of levels
        * taken over all MPI ranks, so <tt>n_levels()<=n_global_levels =
index 3ec8e038627b2adcc756dfc48ecdb1cc175d5d97..af09fb1c5095ef141d23a87c654b640e2eb7db3f 100644 (file)
@@ -87,10 +87,10 @@ namespace GridRefinement
   template <int dim>
   std::pair<double, double>
   adjust_refine_and_coarsen_number_fraction(
-    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);
+    const types::global_cell_index current_n_cells,
+    const types::global_cell_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 621dab4caa5e81def84f726aa62996499be8ecb3..b8902f7225838c14e18cebc40a97184f5cb0106f 100644 (file)
@@ -3006,7 +3006,7 @@ public:
    * may return a value greater than the number of active cells reported by
    * the triangulation object on the current processor.
    */
-  virtual types::global_dof_index
+  virtual types::global_cell_index
   n_global_active_cells() const;
 
 
index fc45b9f954552bebfa4f0d1f45566e7326151e24..4eb262fab715002c599e842cec8b263703d7eaf7 100644 (file)
@@ -228,7 +228,7 @@ namespace internal
           number
           compute_threshold(const dealii::Vector<number> &   criteria,
                             const std::pair<double, double> &global_min_and_max,
-                            const types::global_dof_index    n_target_cells,
+                            const types::global_cell_index   n_target_cells,
                             MPI_Comm                         mpi_communicator)
           {
             double interesting_range[2] = {global_min_and_max.first,
@@ -257,7 +257,7 @@ namespace internal
 
                 // Count how many of our own elements would be above this
                 // threshold:
-                types::global_dof_index my_count =
+                types::global_cell_index my_count =
                   std::count_if(criteria.begin(),
                                 criteria.end(),
                                 [test_threshold](const double c) {
@@ -265,16 +265,8 @@ namespace internal
                                 });
 
                 // Potentially accumulate in a 64bit int to avoid overflow:
-                types::global_dof_index total_count = 0;
-
-                ierr = MPI_Reduce(&my_count,
-                                  &total_count,
-                                  1,
-                                  DEAL_II_DOF_INDEX_MPI_TYPE,
-                                  MPI_SUM,
-                                  master_mpi_rank,
-                                  mpi_communicator);
-                AssertThrowMPI(ierr);
+                types::global_cell_index total_count =
+                  Utilities::MPI::sum(my_count, mpi_communicator);
 
                 // now adjust the range. if we have too many cells, we take the
                 // upper half of the previous range, otherwise the lower half.
@@ -434,9 +426,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 types::global_dof_index max_n_cells)
+        const double                   top_fraction_of_cells,
+        const double                   bottom_fraction_of_cells,
+        const types::global_cell_index max_n_cells)
       {
         Assert(criteria.size() == tria.n_active_cells(),
                ExcDimensionMismatch(criteria.size(), tria.n_active_cells()));
@@ -478,8 +470,8 @@ namespace parallel
           GridRefinement::RefineAndCoarsenFixedNumber::compute_threshold(
             locally_owned_indicators,
             global_min_and_max,
-            static_cast<types::global_dof_index>(adjusted_fractions.first *
-                                                 tria.n_global_active_cells()),
+            static_cast<types::global_cell_index>(adjusted_fractions.first *
+                                                  tria.n_global_active_cells()),
             mpi_communicator);
 
         // compute bottom threshold only if necessary. otherwise use a threshold
@@ -489,8 +481,8 @@ namespace parallel
             GridRefinement::RefineAndCoarsenFixedNumber::compute_threshold(
               locally_owned_indicators,
               global_min_and_max,
-              static_cast<types::global_dof_index>(
-                std::ceil((1 - adjusted_fractions.second) *
+              static_cast<types::global_cell_index>(
+                std::ceil((1. - adjusted_fractions.second) *
                           tria.n_global_active_cells())),
               mpi_communicator);
         else
index 688294715d0850e01ae91faaade780e98031927f..c0fc46f4893decf9a97f45b7e768b210dc92480f 100644 (file)
@@ -34,7 +34,7 @@ for (S : REAL_SCALARS)
               template S
               compute_threshold<S>(const dealii::Vector<S> &,
                                    const std::pair<double, double> &,
-                                   const types::global_dof_index,
+                                   const types::global_cell_index,
                                    MPI_Comm);
             \}
             namespace RefineAndCoarsenFixedFraction
@@ -69,7 +69,7 @@ for (S : REAL_SCALARS; deal_II_dimension : DIMENSIONS)
             const dealii::Vector<S> &,
             const double,
             const double,
-            const types::global_dof_index);
+            const types::global_cell_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 types::global_dof_index);
+            const types::global_cell_index);
 
           template void refine_and_coarsen_fixed_fraction<deal_II_dimension - 1,
                                                           S,
index 189be2a997ab78799a02d9018ccf9d7a5b38a3a3..95ec90f6df3076d277728dcec64bc152f35cdd47 100644 (file)
@@ -131,7 +131,7 @@ namespace parallel
   }
 
   template <int dim, int spacedim>
-  types::global_dof_index
+  types::global_cell_index
   TriangulationBase<dim, spacedim>::n_global_active_cells() const
   {
     return number_cache.n_global_active_cells;
@@ -219,7 +219,7 @@ namespace parallel
 
     // Potentially cast to a 64 bit type before accumulating to avoid overflow:
     number_cache.n_global_active_cells =
-      Utilities::MPI::sum(static_cast<types::global_dof_index>(
+      Utilities::MPI::sum(static_cast<types::global_cell_index>(
                             number_cache.n_locally_owned_active_cells),
                           this->mpi_communicator);
 
index 29dc3bdf706fefc9d534d1d98daad8bcf71e1d4b..93dc2350eb41a43a571960f8c172f471a7f21672 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 types::global_dof_index current_n_cells,
-  const types::global_dof_index max_n_cells,
-  const double                  top_fraction,
-  const double                  bottom_fraction)
+  const types::global_cell_index current_n_cells,
+  const types::global_cell_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<types::global_dof_index>(
+  else if (static_cast<types::global_cell_index>(
              current_n_cells + refine_cells * cell_increase_on_refine -
              coarsen_cells * cell_decrease_on_coarsen) > max_n_cells)
     {
index 8bda2e72e1a70e2d402f252c3fcbe158928e8b3e..14452345cd096070f572c44da7ff4b7d8e9ad1e5 100644 (file)
@@ -99,8 +99,8 @@ for (deal_II_dimension : DIMENSIONS)
   {
     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,
+      deal_II_dimension>(const types::global_cell_index,
+                         const types::global_cell_index,
                          const double,
                          const double);
   }
index 354fb5098559e0c0d4cdcca1af65ef68bc3104fc..87dfd8653cf175d476bf71586247fb045f12c76e 100644 (file)
@@ -12695,7 +12695,7 @@ Triangulation<dim, spacedim>::n_active_cells() const
 }
 
 template <int dim, int spacedim>
-types::global_dof_index
+types::global_cell_index
 Triangulation<dim, spacedim>::n_global_active_cells() const
 {
   return n_active_cells();

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.