]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fix issues reported by cppcheck
authorDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Thu, 26 Apr 2018 14:40:47 +0000 (16:40 +0200)
committerDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Thu, 26 Apr 2018 14:40:47 +0000 (16:40 +0200)
source/base/function_lib.cc
source/base/mpi.cc
source/base/polynomial.cc
source/distributed/grid_refinement.cc
source/distributed/tria.cc
source/dofs/dof_handler_policy.cc
source/grid/grid_in.cc
source/grid/grid_tools.cc
source/grid/manifold_lib.cc
source/lac/trilinos_epetra_vector.cc
source/matrix_free/task_info.cc

index 41d0eedb02382d56df2f90d1096d95154998327d..ae9838e1495ca0729259b39c80869027878fa9d1 100644 (file)
@@ -2212,15 +2212,17 @@ namespace Functions
   {
     Assert(dim==2, ExcNotImplemented());
     AssertDimension(points.size(), values.size());
+#ifndef DEAL_II_HAVE_JN
+    (void) points;
+    (void) values;
+    Assert(false, ExcMessage("The Bessel function jn was not found by CMake."));
+#else
     for (unsigned int k=0; k<points.size(); ++k)
       {
-#ifdef DEAL_II_HAVE_JN
         const double r = points[k].distance(center);
         values[k] = jn(order, r*wave_number);
-#else
-        Assert(false, ExcMessage("The Bessel function jn was not found by CMake."));
-#endif
       }
+#endif
   }
 
 
@@ -2230,11 +2232,11 @@ namespace Functions
                           const unsigned int) const
   {
     Assert(dim==2, ExcNotImplemented());
+#ifdef DEAL_II_HAVE_JN
     const double r = p.distance(center);
     const double co = (r==0.) ? 0. : (p(0)-center(0))/r;
     const double si = (r==0.) ? 0. : (p(1)-center(1))/r;
 
-#ifdef DEAL_II_HAVE_JN
     const double dJn = (order==0)
                        ? (-jn(1, r*wave_number))
                        : (.5*(jn(order-1, wave_number*r) -jn(order+1, wave_number*r)));
@@ -2243,6 +2245,7 @@ namespace Functions
     result[1] = wave_number * si * dJn;
     return result;
 #else
+    (void) p;
     Assert(false, ExcMessage("The Bessel function jn was not found by CMake."));
     return Tensor<1,dim>();
 #endif
@@ -2259,6 +2262,9 @@ namespace Functions
   {
     Assert(dim==2, ExcNotImplemented());
     AssertDimension(points.size(), gradients.size());
+#ifndef DEAL_II_HAVE_JN
+    Assert(false, ExcMessage("The Bessel function jn was not found by CMake."));
+#else
     for (unsigned int k=0; k<points.size(); ++k)
       {
         const Point<dim> &p = points[k];
@@ -2266,18 +2272,14 @@ namespace Functions
         const double co = (r==0.) ? 0. : (p(0)-center(0))/r;
         const double si = (r==0.) ? 0. : (p(1)-center(1))/r;
 
-#ifdef DEAL_II_HAVE_JN
         const double dJn = (order==0)
                            ? (-jn(1, r*wave_number))
                            : (.5*(jn(order-1, wave_number*r) -jn(order+1, wave_number*r)));
-#else
-        const double dJn = 0.;
-        Assert(false, ExcMessage("The Bessel function jn was not found by CMake."));
-#endif
         Tensor<1,dim> &result = gradients[k];
         result[0] = wave_number * co * dJn;
         result[1] = wave_number * si * dJn;
       }
+#endif
   }
 
 
index 6975286e46bdc70903007a7e33dceeae774da9ef..5d0a06e9a271ab07c3370671f0aa5eba1e528484 100644 (file)
@@ -307,7 +307,7 @@ namespace Utilities
                           "in a program since it initializes the MPI system."));
 
 
-      int ierr;
+      int ierr = 0;
 #ifdef DEAL_II_WITH_MPI
       // if we have PETSc, we will initialize it and let it handle MPI.
       // Otherwise, we will do it.
index 87aae95d6e67511d80f37caa20d064fd40903fc4..72952d9762ba5bb52d169e8ed66230b48c113025 100644 (file)
@@ -1005,8 +1005,7 @@ namespace Polynomials
     // check: does the information
     // already exist?
     if (  (recursive_coefficients.size() < k+1) ||
-          ((recursive_coefficients.size() >= k+1) &&
-           (recursive_coefficients[k].get() == nullptr)) )
+          (recursive_coefficients[k].get() == nullptr))
       // no, then generate the
       // respective coefficients
       {
@@ -1040,20 +1039,20 @@ namespace Polynomials
             // later assign it to the
             // coefficients array to
             // make it const
-            std::vector<double> *c0 = new std::vector<double>(2);
-            (*c0)[0] =  1.;
-            (*c0)[1] = -1.;
+            std::vector<double> c0(2);
+            c0[0] =  1.;
+            c0[1] = -1.;
 
-            std::vector<double> *c1 = new std::vector<double>(2);
-            (*c1)[0] = 0.;
-            (*c1)[1] = 1.;
+            std::vector<double> c1(2);
+            c1[0] = 0.;
+            c1[1] = 1.;
 
             // now make these arrays
             // const
             recursive_coefficients[0] =
-              std::unique_ptr<const std::vector<double> >(c0);
+              std_cxx14::make_unique<const std::vector<double> >(std::move(c0));
             recursive_coefficients[1] =
-              std::unique_ptr<const std::vector<double> >(c1);
+              std_cxx14::make_unique<const std::vector<double> >(std::move(c1));
           }
         else if (k==2)
           {
@@ -1061,16 +1060,16 @@ namespace Polynomials
             compute_coefficients(1);
             coefficients_lock.acquire ();
 
-            std::vector<double> *c2 = new std::vector<double>(3);
+            std::vector<double> c2(3);
 
             const double a = 1.; //1./8.;
 
-            (*c2)[0] =   0.*a;
-            (*c2)[1] =  -4.*a;
-            (*c2)[2] =   4.*a;
+            c2[0] =   0.*a;
+            c2[1] =  -4.*a;
+            c2[2] =   4.*a;
 
             recursive_coefficients[2] =
-              std::unique_ptr<const std::vector<double> >(c2);
+              std_cxx14::make_unique<const std::vector<double> >(std::move(c2));
           }
         else
           {
@@ -1086,17 +1085,17 @@ namespace Polynomials
             compute_coefficients(k-1);
             coefficients_lock.acquire ();
 
-            std::vector<double> *ck = new std::vector<double>(k+1);
+            std::vector<double> ck (k+1);
 
             const double a = 1.; //1./(2.*k);
 
-            (*ck)[0] = - a*(*recursive_coefficients[k-1])[0];
+            ck[0] = - a*(*recursive_coefficients[k-1])[0];
 
             for (unsigned int i=1; i<=k-1; ++i)
-              (*ck)[i] = a*( 2.*(*recursive_coefficients[k-1])[i-1]
-                             - (*recursive_coefficients[k-1])[i] );
+              ck[i] = a*( 2.*(*recursive_coefficients[k-1])[i-1]
+                          - (*recursive_coefficients[k-1])[i] );
 
-            (*ck)[k] = a*2.*(*recursive_coefficients[k-1])[k-1];
+            ck[k] = a*2.*(*recursive_coefficients[k-1])[k-1];
             // for even degrees, we need
             // to add a multiple of
             // basis fcn phi_2
@@ -1106,15 +1105,15 @@ namespace Polynomials
                 //for (unsigned int i=1; i<=k; i++)
                 //  b /= 2.*i;
 
-                (*ck)[1] += b*(*recursive_coefficients[2])[1];
-                (*ck)[2] += b*(*recursive_coefficients[2])[2];
+                ck[1] += b*(*recursive_coefficients[2])[1];
+                ck[2] += b*(*recursive_coefficients[2])[2];
               }
             // finally assign the newly
             // created vector to the
             // const pointer in the
             // coefficients array
             recursive_coefficients[k] =
-              std::unique_ptr<const std::vector<double> >(ck);
+              std_cxx14::make_unique<const std::vector<double> >(std::move(ck));
           };
       };
   }
index 978448b9bbbdda8f65c3fdfc48e3e7525e03b2e3..319d39ee7e189956b1f68b9278bad38f967e6e30 100644 (file)
@@ -226,7 +226,7 @@ namespace
     template <typename number>
     number
     compute_threshold (const dealii::Vector<number>   &criteria,
-                       const std::pair<double,double>  global_min_and_max,
+                       const std::pair<double,double> &global_min_and_max,
                        const unsigned int              n_target_cells,
                        MPI_Comm                        mpi_communicator)
     {
@@ -316,7 +316,7 @@ namespace
     template <typename number>
     number
     compute_threshold (const dealii::Vector<number>   &criteria,
-                       const std::pair<double,double>  global_min_and_max,
+                       const std::pair<double,double> &global_min_and_max,
                        const double                    target_error,
                        MPI_Comm                        mpi_communicator)
     {
index 6f1e0e244a03053e1c62e13520ea490dba0ebbec..d653e0421b77c493404a87c4e6875b6dd52e2fb3 100644 (file)
@@ -519,10 +519,9 @@ namespace
                   typename internal::p4est::types<dim>::quadrant &ghost_quadrant,
                   types::subdomain_id                             ghost_owner)
   {
-    int i, child_id;
     int l = ghost_quadrant.level;
 
-    for (i = 0; i < l; i++)
+    for (int i = 0; i < l; i++)
       {
         typename Triangulation<dim,spacedim>::cell_iterator cell (tria, i, dealii_index);
         if (cell->has_children () == false)
@@ -532,7 +531,7 @@ namespace
             return;
           }
 
-        child_id = internal::p4est::functions<dim>::quadrant_ancestor_id (&ghost_quadrant, i + 1);
+        const int child_id = internal::p4est::functions<dim>::quadrant_ancestor_id (&ghost_quadrant, i + 1);
         dealii_index = cell->child_index(child_id);
       }
 
@@ -1223,7 +1222,7 @@ namespace
      * order that p4est will encounter the cells, and they do not contain
      * ghost cells or artificial cells.
      */
-    PartitionWeights (const std::vector<unsigned int> &cell_weights);
+    explicit PartitionWeights (const std::vector<unsigned int> &cell_weights);
 
     /**
      * A callback function that we pass to the p4est data structures when a
index c260432503928a86aed47a7ede91f460af41c4be..93f5cac25e8da85b0616b26d504f8f5d7b5e470a 100644 (file)
@@ -122,8 +122,7 @@ namespace internal
                   void *)
           {
             if (cell->has_children()
-                ||
-                (!cell->has_children() && !cell->is_artificial()))
+                || !cell->is_artificial())
               cell->update_cell_dof_indices_cache ();
           };
 
index 24954a41fca1b98ebd54b649ae62c9d7edc09c49..8db5169fbb0e1ece9005d4c305fe0406862e8cc4 100644 (file)
@@ -1640,7 +1640,6 @@ void GridIn<2>::read_netcdf (const std::string &filename)
 #else
   const unsigned int dim=2;
   const unsigned int spacedim=2;
-  const bool output=false;
   Assert (tria != 0, ExcNoTriangulationSelected());
   // this function assumes the TAU
   // grid format.
@@ -1913,7 +1912,6 @@ void GridIn<3>::read_netcdf (const std::string &filename)
 #else
   const unsigned int dim=3;
   const unsigned int spacedim=3;
-  const bool output=false;
   Assert (tria != 0, ExcNoTriangulationSelected());
   // this function assumes the TAU
   // grid format.
@@ -1926,8 +1924,6 @@ void GridIn<3>::read_netcdf (const std::string &filename)
   NcDim *elements_dim=nc.get_dim("no_of_elements");
   AssertThrow(elements_dim->is_valid(), ExcIO());
   const unsigned int n_cells=elements_dim->size();
-  if (output)
-    std::cout << "n_cell=" << n_cells << std::endl;
   // and n_hexes
   NcDim *hexes_dim=nc.get_dim("no_of_hexaeders");
   AssertThrow(hexes_dim->is_valid(), ExcIO());
@@ -1960,17 +1956,6 @@ void GridIn<3>::read_netcdf (const std::string &filename)
   for (unsigned int i=0; i<vertex_indices.size(); ++i)
     AssertThrow(vertex_indices[i]>=0, ExcIO());
 
-  if (output)
-    {
-      std::cout << "vertex_indices:" << std::endl;
-      for (unsigned int cell=0, v=0; cell<n_cells; ++cell)
-        {
-          for (unsigned int i=0; i<vertices_per_hex; ++i)
-            std::cout << vertex_indices[v++] << " ";
-          std::cout << std::endl;
-        }
-    }
-
   // next we read
   //   double points_xc(no_of_points)
   //   double points_yc(no_of_points)
@@ -1978,8 +1963,6 @@ void GridIn<3>::read_netcdf (const std::string &filename)
   NcDim *vertices_dim=nc.get_dim("no_of_points");
   AssertThrow(vertices_dim->is_valid(), ExcIO());
   const unsigned int n_vertices=vertices_dim->size();
-  if (output)
-    std::cout << "n_vertices=" << n_vertices << std::endl;
 
   NcVar *points_xc=nc.get_var("points_xc");
   NcVar *points_yc=nc.get_var("points_yc");
@@ -1998,19 +1981,9 @@ void GridIn<3>::read_netcdf (const std::string &filename)
               static_cast<int>(n_vertices), ExcIO());
   std::vector<std::vector<double> > point_values(
     3, std::vector<double> (n_vertices));
-  // we switch y and z
-  const bool switch_y_z=false;
   points_xc->get(&*point_values[0].begin(), n_vertices);
-  if (switch_y_z)
-    {
-      points_yc->get(&*point_values[2].begin(), n_vertices);
-      points_zc->get(&*point_values[1].begin(), n_vertices);
-    }
-  else
-    {
-      points_yc->get(&*point_values[1].begin(), n_vertices);
-      points_zc->get(&*point_values[2].begin(), n_vertices);
-    }
+  points_yc->get(&*point_values[1].begin(), n_vertices);
+  points_zc->get(&*point_values[2].begin(), n_vertices);
 
   // and fill the vertices
   std::vector<Point<spacedim> > vertices (n_vertices);
@@ -2052,18 +2025,6 @@ void GridIn<3>::read_netcdf (const std::string &filename)
   std::vector<int> bvertex_indices(n_bquads*vertices_per_quad);
   bvertex_indices_var->get(&*bvertex_indices.begin(), n_bquads, vertices_per_quad);
 
-  if (output)
-    {
-      std::cout << "bquads: ";
-      for (unsigned int i=0; i<n_bquads; ++i)
-        {
-          for (unsigned int v=0; v<GeometryInfo<dim>::vertices_per_face; ++v)
-            std::cout << bvertex_indices[
-                        i*GeometryInfo<dim>::vertices_per_face+v] << " ";
-          std::cout << std::endl;
-        }
-    }
-
   // next we read
   // int boundarymarker_of_surfaces(
   //   no_of_surfaceelements)
index e244e2b21d8ce4ae43a58b194a3e1a8f8067448c..58083093c7659d1d066906a7ca12720f3512ccc3 100644 (file)
@@ -566,7 +566,7 @@ namespace GridTools
     class Shift
     {
     public:
-      Shift (const Tensor<1,spacedim> &shift)
+      explicit Shift (const Tensor<1,spacedim> &shift)
         :
         shift(shift)
       {}
@@ -585,7 +585,7 @@ namespace GridTools
     class Rotate2d
     {
     public:
-      Rotate2d (const double angle)
+      explicit Rotate2d (const double angle)
         :
         angle(angle)
       {}
@@ -633,7 +633,7 @@ namespace GridTools
     class Scale
     {
     public:
-      Scale (const double factor)
+      explicit Scale (const double factor)
         :
         factor(factor)
       {}
index ab2178a4736ff6ddc9ec1aec034ed9e65574d3fa..5d775c71c12b0bd0ff1a09018894823e21dc9d7c 100644 (file)
@@ -615,7 +615,7 @@ get_new_points (const ArrayView<const Point<spacedim>> &surrounding_points,
   // be good enough. This tolerance was chosen such that the first iteration
   // for a at least three time refined HyperShell mesh with radii .5 and 1.
   // doesn't already succeed.
-  if (max_distance < 2e-2 || spacedim<3)
+  if (max_distance < 2e-2)
     {
       for (unsigned int row=0; row<weight_rows; ++row)
         new_points[row] = center + new_candidates[row].first * new_candidates[row].second;
@@ -1469,7 +1469,7 @@ namespace
 
         for (unsigned int line=0; line<GeometryInfo<2>::lines_per_cell; ++line)
           {
-            const double my_weight = line%2 ? chart_point[line/2] : 1-chart_point[line/2];
+            const double my_weight = (line%2) ? chart_point[line/2] : 1-chart_point[line/2];
             const double line_point = chart_point[1-line/2];
 
             // Same manifold or invalid id which will go back to the same
index 47446f04e90dfbe5a410dba12e980764f12f46dd..a73927d6b5e5d6f2215dd133ae65fb7ffd607a26 100644 (file)
@@ -145,8 +145,7 @@ namespace LinearAlgebra
           // Check if the communication pattern already exists and if it can be
           // reused.
           if ((source_stored_elements.size() != V.get_stored_elements().size()) ||
-              ((source_stored_elements.size() == V.get_stored_elements().size()) &&
-               (source_stored_elements != V.get_stored_elements())))
+              (source_stored_elements != V.get_stored_elements()))
             {
               create_epetra_comm_pattern(V.get_stored_elements(),
                                          dynamic_cast<const Epetra_MpiComm &>(vector->Comm()).Comm());
index 97dd3fe79b2d1ba1cc6059a0cb5505b7c49f4c9c..8ef9e0b814935977e26a9b0db477ecf0fe6c1436 100644 (file)
@@ -1413,8 +1413,8 @@ namespace internal
       Assert(!hp_bool || cell_active_fe_index.size() == n_active_cells,
              ExcInternalError());
 
-      unsigned int n_macro_cells_before = 0;
       {
+        unsigned int n_macro_cells_before = 0;
         // Create partitioning within partitions.
 
         // For each block of cells, this variable saves to which partitions
@@ -1426,7 +1426,6 @@ namespace internal
         partition_row_index.resize(partition+1,0);
         cell_partition_data.resize(1,0);
 
-        unsigned int start_up = 0;
         unsigned int counter = 0;
         unsigned int missing_macros;
         for (unsigned int part=0; part<partition; ++part)
@@ -1435,7 +1434,7 @@ namespace internal
             neighbor_list.resize(0);
             bool work = true;
             unsigned int partition_l2 = 0;
-            start_up = partition_size[part];
+            unsigned int start_up = partition_size[part];
             unsigned int partition_counter = 0;
             while (work)
               {
@@ -1680,10 +1679,7 @@ namespace internal
      const std::vector<unsigned int>       &partition_size,
      std::vector<unsigned int>             &partition_color_list)
     {
-
       const unsigned int n_macro_cells = *(cell_partition_data.end()-2);
-      std::vector<unsigned int> neighbor_list;
-      std::vector<unsigned int> neighbor_neighbor_list;
       std::vector<unsigned int> cell_color(n_blocks, n_macro_cells);
       std::vector<bool> color_finder;
 

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.