]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Avoid C-style cast in grid
authorDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Wed, 28 Nov 2018 21:56:46 +0000 (22:56 +0100)
committerDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Wed, 28 Nov 2018 21:56:46 +0000 (22:56 +0100)
include/deal.II/grid/filtered_iterator.h
source/grid/grid_generator.cc
source/grid/grid_in.cc
source/grid/grid_out.cc
source/grid/grid_refinement.cc
source/grid/grid_tools.cc
source/grid/tria.cc
source/grid/tria_levels.cc

index d4f7189d1c0b2fdd0950f3f1c74584dd2443cf73..c7b6b8e771e48f9898f76d524598de01f055c989 100644 (file)
@@ -960,7 +960,7 @@ inline FilteredIterator<BaseIterator>::FilteredIterator(
     // address of fi, GCC would not cast fi to the base class of type
     // BaseIterator but tries to go through constructing a new
     // BaseIterator with an Accessor.
-  BaseIterator(*(BaseIterator *)(&fi))
+  BaseIterator(*static_cast<const BaseIterator *>(&fi))
   , predicate(fi.predicate->clone())
 {}
 
index 8ed2339b288e66c87fe5e45d44961edcc9bdbf07..392680e204485ae3ab0eff78fe3a56e7a1e722a5 100644 (file)
@@ -1176,22 +1176,21 @@ namespace GridGenerator
       {
         case 1:
           for (unsigned int x = 0; x <= repetitions[0]; ++x)
-            points.push_back(p1 + (double)x * delta[0]);
+            points.push_back(p1 + x * delta[0]);
           break;
 
         case 2:
           for (unsigned int y = 0; y <= repetitions[1]; ++y)
             for (unsigned int x = 0; x <= repetitions[0]; ++x)
-              points.push_back(p1 + (double)x * delta[0] +
-                               (double)y * delta[1]);
+              points.push_back(p1 + x * delta[0] + y * delta[1]);
           break;
 
         case 3:
           for (unsigned int z = 0; z <= repetitions[2]; ++z)
             for (unsigned int y = 0; y <= repetitions[1]; ++y)
               for (unsigned int x = 0; x <= repetitions[0]; ++x)
-                points.push_back(p1 + (double)x * delta[0] +
-                                 (double)y * delta[1] + (double)z * delta[2]);
+                points.push_back(p1 + x * delta[0] + y * delta[1] +
+                                 z * delta[2]);
           break;
 
         default:
@@ -1806,22 +1805,21 @@ namespace GridGenerator
       {
         case 1:
           for (unsigned int x = 0; x <= repetitions[0]; ++x)
-            points.push_back(p1 + (double)x * delta[0]);
+            points.push_back(p1 + x * delta[0]);
           break;
 
         case 2:
           for (unsigned int y = 0; y <= repetitions[1]; ++y)
             for (unsigned int x = 0; x <= repetitions[0]; ++x)
-              points.push_back(p1 + (double)x * delta[0] +
-                               (double)y * delta[1]);
+              points.push_back(p1 + x * delta[0] + y * delta[1]);
           break;
 
         case 3:
           for (unsigned int z = 0; z <= repetitions[2]; ++z)
             for (unsigned int y = 0; y <= repetitions[1]; ++y)
               for (unsigned int x = 0; x <= repetitions[0]; ++x)
-                points.push_back(p1 + (double)x * delta[0] +
-                                 (double)y * delta[1] + (double)z * delta[2]);
+                points.push_back(p1 + x * delta[0] + y * delta[1] +
+                                 z * delta[2]);
           break;
 
         default:
@@ -5040,7 +5038,7 @@ namespace GridGenerator
                       treated_vertices[vv] = true;
                       for (unsigned int i = 0; i <= Nz; ++i)
                         {
-                          double d = ((double)i) * L / ((double)Nz);
+                          double d = i * L / Nz;
                           switch (vv - i * 16)
                             {
                               case 1:
index 323547464403926228ca12e5525b11e6b9c8f0ef..8e091d0a48e35b5cb8063119bd183b998fcdae72 100644 (file)
@@ -2961,7 +2961,7 @@ GridIn<dim, spacedim>::read_assimp(const std::string &filename,
                   cells[valid_cell].vertices[f] =
                     mFaces[i].mIndices[f] + v_offset;
                 }
-              cells[valid_cell].material_id = (types::material_id)m;
+              cells[valid_cell].material_id = m;
               ++valid_cell;
             }
           else
index ed724e5607a34c79ebe765476eb573a03c934780..1b13dc9b7a4cb6d459790f3c2da7f71abc1d7b25 100644 (file)
@@ -875,7 +875,7 @@ GridOut::write_dx(const Triangulation<dim, spacedim> &tria,
       out << "object \"material\" class array type int rank 0 items " << n_cells
           << " data follows" << '\n';
       for (cell = tria.begin_active(); cell != endc; ++cell)
-        out << ' ' << (unsigned int)cell->material_id();
+        out << ' ' << cell->material_id();
       out << '\n' << "attribute \"dep\" string \"connections\"" << '\n' << '\n';
 
       out << "object \"level\" class array type int rank 0 items " << n_cells
@@ -944,10 +944,17 @@ GridOut::write_dx(const Triangulation<dim, spacedim> &tria,
           << " data follows" << '\n';
       for (cell = tria.begin_active(); cell != endc; ++cell)
         {
-          // Little trick to get -1
-          // for the interior
+          // Little trick to get -1 for the interior
           for (unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; ++f)
-            out << ' ' << (int)(signed char)cell->face(f)->boundary_id();
+            {
+              const types::boundary_id boundary_id =
+                cell->face(f)->boundary_id();
+              out << ' '
+                  << (boundary_id == static_cast<types::boundary_id>(-1) ?
+                        -1 :
+                        static_cast<std::make_signed<types::boundary_id>::type>(
+                          boundary_id));
+            }
           out << '\n';
         }
       out << "attribute \"dep\" string \"connections\"" << '\n' << '\n';
@@ -1451,8 +1458,8 @@ GridOut::write_xfig(const Triangulation<2> &tria,
             cell->vertex(GeometryInfo<dim>::ucd_to_deal[k % nv]);
           for (unsigned int d = 0; d < static_cast<unsigned int>(dim); ++d)
             {
-              int val = (int)(1200 * xfig_flags.scaling(d) *
-                              (p(d) - xfig_flags.offset(d)));
+              int val =
+                1200 * xfig_flags.scaling(d) * (p(d) - xfig_flags.offset(d));
               out << '\t' << ((d == 0) ? val : -val);
             }
           out << std::endl;
@@ -1471,8 +1478,7 @@ GridOut::write_xfig(const Triangulation<2> &tria,
                 out << "2 1 "
                     // with line style and thickness
                     << xfig_flags.boundary_style << ' '
-                    << xfig_flags.boundary_thickness << ' '
-                    << (1 + (unsigned int)bi);
+                    << xfig_flags.boundary_thickness << ' ' << 1 + bi;
                 // Fill color
                 out << " -1 ";
                 // Depth 100 less than cells
@@ -1496,8 +1502,8 @@ GridOut::write_xfig(const Triangulation<2> &tria,
                     for (unsigned int d = 0; d < static_cast<unsigned int>(dim);
                          ++d)
                       {
-                        int val = (int)(1200 * xfig_flags.scaling(d) *
-                                        (p(d) - xfig_flags.offset(d)));
+                        int val = 1200 * xfig_flags.scaling(d) *
+                                  (p(d) - xfig_flags.offset(d));
                         out << '\t' << ((d == 0) ? val : -val);
                       }
                     out << std::endl;
@@ -1609,9 +1615,7 @@ GridOut::write_svg(const Triangulation<2, 2> &tria, std::ostream &out) const
   // bounding box of the given triangulation and check
   // the cells for material id, level number, subdomain id
   // (, and level subdomain id).
-  for (Triangulation<2, 2>::cell_iterator cell = tria.begin();
-       cell != tria.end();
-       ++cell)
+  for (const auto &cell : tria.cell_iterators())
     {
       for (unsigned int vertex_index = 0; vertex_index < 4; vertex_index++)
         {
@@ -1626,13 +1630,13 @@ GridOut::write_svg(const Triangulation<2, 2> &tria, std::ostream &out) const
             y_max = cell->vertex(vertex_index)[1];
         }
 
-      if ((unsigned int)cell->level() < min_level)
+      if (static_cast<unsigned int>(cell->level()) < min_level)
         min_level = cell->level();
-      if ((unsigned int)cell->level() > max_level)
+      if (static_cast<unsigned int>(cell->level()) > max_level)
         max_level = cell->level();
 
-      materials[(unsigned int)cell->material_id()] = 1;
-      levels[(unsigned int)cell->level()]          = 1;
+      materials[cell->material_id()] = 1;
+      levels[cell->level()]          = 1;
       if (cell->active())
         subdomains[cell->subdomain_id() + 2] = 1;
       level_subdomains[cell->level_subdomain_id() + 2] = 1;
@@ -1795,7 +1799,8 @@ GridOut::write_svg(const Triangulation<2, 2> &tria, std::ostream &out) const
   if (svg_flags.convert_level_number_to_height)
     {
       point[2] = svg_flags.level_height_factor *
-                 ((float)tria.begin()->level() / (float)n_levels) *
+                 (static_cast<float>(tria.begin()->level()) /
+                  static_cast<float>(n_levels)) *
                  std::max(x_dimension, y_dimension);
     }
 
@@ -1818,9 +1823,10 @@ GridOut::write_svg(const Triangulation<2, 2> &tria, std::ostream &out) const
 
       if (svg_flags.convert_level_number_to_height)
         {
-          point[2] = svg_flags.level_height_factor *
-                     ((float)cell->level() / (float)n_levels) *
-                     std::max(x_dimension, y_dimension);
+          point[2] =
+            svg_flags.level_height_factor *
+            (static_cast<float>(cell->level()) / static_cast<float>(n_levels)) *
+            std::max(x_dimension, y_dimension);
         }
 
       projection_decomposition = GridOut::svg_project_point(point,
@@ -1896,7 +1902,7 @@ GridOut::write_svg(const Triangulation<2, 2> &tria, std::ostream &out) const
       if (y_min_perspective > projection_decomposition[1])
         y_min_perspective = projection_decomposition[1];
 
-      if ((unsigned int)cell->level() == min_level)
+      if (static_cast<unsigned int>(cell->level()) == min_level)
         min_level_min_vertex_distance = cell->minimum_vertex_distance();
     }
 
@@ -2139,10 +2145,10 @@ GridOut::write_svg(const Triangulation<2, 2> &tria, std::ostream &out) const
               switch (svg_flags.coloring)
                 {
                   case GridOutFlags::Svg::material_id:
-                    out << (unsigned int)cell->material_id();
+                    out << cell->material_id();
                     break;
                   case GridOutFlags::Svg::level_number:
-                    out << (unsigned int)cell->level();
+                    out << static_cast<unsigned int>(cell->level());
                     break;
                   case GridOutFlags::Svg::subdomain_id:
                     if (cell->active())
@@ -2169,7 +2175,8 @@ GridOut::write_svg(const Triangulation<2, 2> &tria, std::ostream &out) const
           if (svg_flags.convert_level_number_to_height)
             {
               point[2] = svg_flags.level_height_factor *
-                         ((float)cell->level() / (float)n_levels) *
+                         (static_cast<float>(cell->level()) /
+                          static_cast<float>(n_levels)) *
                          std::max(x_dimension, y_dimension);
             }
 
@@ -2307,7 +2314,8 @@ GridOut::write_svg(const Triangulation<2, 2> &tria, std::ostream &out) const
               if (svg_flags.convert_level_number_to_height)
                 {
                   point[2] = svg_flags.level_height_factor *
-                             ((float)cell->level() / (float)n_levels) *
+                             (static_cast<float>(cell->level()) /
+                              static_cast<float>(n_levels)) *
                              std::max(x_dimension, y_dimension);
                 }
 
@@ -2326,11 +2334,8 @@ GridOut::write_svg(const Triangulation<2, 2> &tria, std::ostream &out) const
                                            camera_focus);
 
               const unsigned int font_size_this_cell =
-                static_cast<unsigned int>(.5 +
-                                          cell_label_font_size *
-                                            std::pow(.5,
-                                                     (float)cell->level() - 4. +
-                                                       3.5 * distance_factor));
+                .5 + cell_label_font_size *
+                       std::pow(.5, cell->level() - 4. + 3.5 * distance_factor);
 
               out << "  <text"
                   << " x=\""
@@ -2366,7 +2371,7 @@ GridOut::write_svg(const Triangulation<2, 2> &tria, std::ostream &out) const
                   if (svg_flags.label_level_number ||
                       svg_flags.label_cell_index)
                     out << ',';
-                  out << (int)cell->material_id();
+                  out << static_cast<int>(cell->material_id());
                 }
 
               if (svg_flags.label_subdomain_id)
@@ -2408,7 +2413,8 @@ GridOut::write_svg(const Triangulation<2, 2> &tria, std::ostream &out) const
                       if (svg_flags.convert_level_number_to_height)
                         {
                           point[2] = svg_flags.level_height_factor *
-                                     ((float)cell->level() / (float)n_levels) *
+                                     (static_cast<float>(cell->level()) /
+                                      static_cast<float>(n_levels)) *
                                      std::max(x_dimension, y_dimension);
                         }
 
@@ -2445,7 +2451,8 @@ GridOut::write_svg(const Triangulation<2, 2> &tria, std::ostream &out) const
                       if (svg_flags.convert_level_number_to_height)
                         {
                           point[2] = svg_flags.level_height_factor *
-                                     ((float)cell->level() / (float)n_levels) *
+                                     (static_cast<float>(cell->level()) /
+                                      static_cast<float>(n_levels)) *
                                      std::max(x_dimension, y_dimension);
                         }
 
@@ -3235,12 +3242,10 @@ GridOut::write_mesh_per_processor_as_vtu(
           patch.vertices[vertex] = cell->vertex(vertex);
           patch.data(0, vertex)  = cell->level();
           if (!cell->has_children())
-            patch.data(1, vertex) =
-              (double)static_cast<int>(cell->subdomain_id());
+            patch.data(1, vertex) = static_cast<int>(cell->subdomain_id());
           else
             patch.data(1, vertex) = -1.0;
-          patch.data(2, vertex) =
-            (double)static_cast<int>(cell->level_subdomain_id());
+          patch.data(2, vertex) = static_cast<int>(cell->level_subdomain_id());
           patch.data(3, vertex) = tria.locally_owned_subdomain();
         }
 
index 785679955304fbe87ba73a078db8d5621f9d07ac..0c6af2b6cc22c760425d3087eb2eda3215fa51ca 100644 (file)
@@ -419,9 +419,9 @@ GridRefinement::refine_and_coarsen_optimize(Triangulation<dim, spacedim> &tria,
       expected_error_reduction +=
         (1 - std::pow(2., -1. * order)) * criteria(cell_indices[M]);
 
-      const double cost =
-        std::pow(((std::pow(2., dim) - 1) * (1 + M) + N), (double)order / dim) *
-        (original_error - expected_error_reduction);
+      const double cost = std::pow(((std::pow(2., dim) - 1) * (1 + M) + N),
+                                   static_cast<double>(order) / dim) *
+                          (original_error - expected_error_reduction);
       if (cost <= min_cost)
         {
           min_cost = cost;
index 24def8abd2c574f2075b210c4b526668c22f2cc4..2a85cf49a75d6f0b99d64517cdbf8b8f9c010ebb 100644 (file)
@@ -2812,8 +2812,8 @@ namespace GridTools
       if (cell->active())
         {
           while (current_cell_idx >=
-                 std::floor((long)n_active_cells * (current_proc_idx + 1) /
-                            n_partitions))
+                 std::floor(static_cast<long>(n_active_cells) *
+                            (current_proc_idx + 1) / n_partitions))
             ++current_proc_idx;
           cell->set_subdomain_id(current_proc_idx);
           ++current_cell_idx;
@@ -3674,12 +3674,12 @@ namespace GridTools
           if (cell->face(f)->at_boundary())
             for (unsigned int e = 0; e < GeometryInfo<dim>::lines_per_face; ++e)
               {
-                const auto bid = cell->face(f)->line(e)->boundary_id();
-                const auto ind = std::find(src_boundary_ids.begin(),
-                                           src_boundary_ids.end(),
-                                           bid) -
-                                 src_boundary_ids.begin();
-                if ((unsigned int)ind < src_boundary_ids.size())
+                const auto         bid = cell->face(f)->line(e)->boundary_id();
+                const unsigned int ind = std::find(src_boundary_ids.begin(),
+                                                   src_boundary_ids.end(),
+                                                   bid) -
+                                         src_boundary_ids.begin();
+                if (ind < src_boundary_ids.size())
                   cell->face(f)->line(e)->set_manifold_id(
                     dst_manifold_ids[ind]);
               }
@@ -3692,12 +3692,12 @@ namespace GridTools
       for (unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; ++f)
         if (cell->face(f)->at_boundary())
           {
-            const auto bid = cell->face(f)->boundary_id();
-            const auto ind =
+            const auto         bid = cell->face(f)->boundary_id();
+            const unsigned int ind =
               std::find(src_boundary_ids.begin(), src_boundary_ids.end(), bid) -
               src_boundary_ids.begin();
 
-            if ((unsigned int)ind < src_boundary_ids.size())
+            if (ind < src_boundary_ids.size())
               {
                 // assign the manifold id
                 cell->face(f)->set_manifold_id(dst_manifold_ids[ind]);
@@ -3710,11 +3710,11 @@ namespace GridTools
                    ++e)
                 {
                   const auto bid = cell->face(f)->line(e)->boundary_id();
-                  const auto ind = std::find(src_boundary_ids.begin(),
-                                             src_boundary_ids.end(),
-                                             bid) -
-                                   src_boundary_ids.begin();
-                  if ((unsigned int)ind < src_boundary_ids.size())
+                  const unsigned int ind = std::find(src_boundary_ids.begin(),
+                                                     src_boundary_ids.end(),
+                                                     bid) -
+                                           src_boundary_ids.begin();
+                  if (ind < src_boundary_ids.size())
                     cell->face(f)->line(e)->set_boundary_id(
                       reset_boundary_ids[ind]);
                 }
index 36fef9772d2110ee35a5d4526a025ccbe9b335f0..34ece8f21c07a63b4da11ac227902da97181a8cb 100644 (file)
@@ -1195,7 +1195,8 @@ namespace internal
       types::boundary_id,
       << "The input data for creating a triangulation contained "
       << "information about a line with indices " << arg1 << " and " << arg2
-      << " that is described to have boundary indicator " << (int)arg3
+      << " that is described to have boundary indicator "
+      << static_cast<int>(arg3)
       << ". However, this is an internal line not located on the "
       << "boundary. You cannot assign a boundary indicator to it." << std::endl
       << std::endl
@@ -1228,7 +1229,8 @@ namespace internal
       << "The input data for creating a triangulation contained "
       << "information about a quad with indices " << arg1 << ", " << arg2
       << ", " << arg3 << ", and " << arg4
-      << " that is described to have boundary indicator " << (int)arg5
+      << " that is described to have boundary indicator "
+      << static_cast<int>(arg5)
       << ". However, this is an internal quad not located on the "
       << "boundary. You cannot assign a boundary indicator to it." << std::endl
       << std::endl
@@ -1493,9 +1495,10 @@ namespace internal
         // access any of these numbers, we can do this in the
         // background
         Threads::Task<void> update_lines = Threads::new_task(
-          (void (*)(const Triangulation<dim, spacedim> &,
-                    const unsigned int,
-                    internal::TriangulationImplementation::NumberCache<1> &))(
+          static_cast<
+            void (*)(const Triangulation<dim, spacedim> &,
+                     const unsigned int,
+                     internal::TriangulationImplementation::NumberCache<1> &)>(
             &compute_number_cache<dim, spacedim>),
           triangulation,
           level_objects,
@@ -1599,9 +1602,10 @@ namespace internal
         // don't access any of these numbers, we can do this in the
         // background
         Threads::Task<void> update_quads_and_lines = Threads::new_task(
-          (void (*)(const Triangulation<dim, spacedim> &,
-                    const unsigned int,
-                    internal::TriangulationImplementation::NumberCache<2> &))(
+          static_cast<
+            void (*)(const Triangulation<dim, spacedim> &,
+                     const unsigned int,
+                     internal::TriangulationImplementation::NumberCache<2> &)>(
             &compute_number_cache<dim, spacedim>),
           triangulation,
           level_objects,
index 54a1fbe410c6bc14e0cc8120caca15cd51dea11f..ba07d372e273aadcf7878481e5f1c34351e893d8 100644 (file)
@@ -73,7 +73,7 @@ namespace internal
           else
             direction_flags.clear();
 
-          parents.reserve((int)(total_cells + 1) / 2);
+          parents.reserve((total_cells + 1) / 2);
           parents.insert(parents.end(),
                          (total_cells + 1) / 2 - parents.size(),
                          -1);
@@ -165,7 +165,7 @@ namespace internal
           else
             direction_flags.clear();
 
-          parents.reserve((int)(total_cells + 1) / 2);
+          parents.reserve((total_cells + 1) / 2);
           parents.insert(parents.end(),
                          (total_cells + 1) / 2 - parents.size(),
                          -1);

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.