From 428a74f83e242dee6b6d8926f0e3aa81a5c318e5 Mon Sep 17 00:00:00 2001 From: Daniel Arndt Date: Thu, 4 Jul 2024 11:35:52 -0400 Subject: [PATCH] Fix CppCheck complaints --- include/deal.II/algorithms/any_data.h | 4 ++-- include/deal.II/algorithms/newton.templates.h | 2 +- include/deal.II/algorithms/timestep_control.h | 2 +- include/deal.II/base/parallel.h | 3 ++- include/deal.II/base/table_handler.h | 4 ++-- include/deal.II/base/thread_management.h | 4 ++-- include/deal.II/distributed/cell_data_transfer.h | 4 ++-- .../distributed/cell_data_transfer.templates.h | 4 ++-- include/deal.II/fe/fe_system.h | 9 ++++----- include/deal.II/grid/tria.h | 2 +- include/deal.II/lac/sparsity_pattern.h | 15 ++++++++------- source/algorithms/timestep_control.cc | 2 +- source/base/convergence_table.cc | 10 ++++++---- source/base/exceptions.cc | 2 +- source/base/mpi_compute_index_owner_internal.cc | 6 ++---- source/base/mu_parser_internal.cc | 7 +++---- source/base/parallel.cc | 2 +- source/base/parameter_acceptor.cc | 9 ++------- source/base/parameter_handler.cc | 2 +- source/base/polynomials_barycentric.cc | 3 +-- source/base/quadrature_lib.cc | 4 ++-- source/base/table_handler.cc | 9 +++++---- source/base/utilities.cc | 8 ++++---- source/dofs/dof_handler.cc | 8 -------- source/dofs/dof_tools.cc | 2 +- source/fe/fe_dgq.cc | 4 ++-- 26 files changed, 59 insertions(+), 72 deletions(-) diff --git a/include/deal.II/algorithms/any_data.h b/include/deal.II/algorithms/any_data.h index 6ed4250633..c23fa925e8 100644 --- a/include/deal.II/algorithms/any_data.h +++ b/include/deal.II/algorithms/any_data.h @@ -230,7 +230,7 @@ inline type AnyData::entry(const unsigned int i) { AssertIndexRange(i, size()); - type *p = std::any_cast(&data[i]); + const type *p = std::any_cast(&data[i]); Assert(p != nullptr, ExcTypeMismatch(typeid(type).name(), data[i].type().name())); return *p; @@ -349,7 +349,7 @@ inline type AnyData::entry(const std::string &n) { const unsigned int i = find(n); - type *p = std::any_cast(&data[i]); + const type *p = std::any_cast(&data[i]); Assert(p != 0, ExcTypeMismatch(typeid(type).name(), data[i].type().name())); return *p; } diff --git a/include/deal.II/algorithms/newton.templates.h b/include/deal.II/algorithms/newton.templates.h index 55c81f0a69..6feed39cc7 100644 --- a/include/deal.II/algorithms/newton.templates.h +++ b/include/deal.II/algorithms/newton.templates.h @@ -156,7 +156,7 @@ namespace Algorithms { (*inverse_derivative)(out2, src2); } - catch (SolverControl::NoConvergence &e) + catch (const SolverControl::NoConvergence &e) { deallog << "Inner iteration failed after " << e.last_step << " steps with residual " << e.last_residual << std::endl; diff --git a/include/deal.II/algorithms/timestep_control.h b/include/deal.II/algorithms/timestep_control.h index a9f71f1fca..72aed99d98 100644 --- a/include/deal.II/algorithms/timestep_control.h +++ b/include/deal.II/algorithms/timestep_control.h @@ -84,7 +84,7 @@ namespace Algorithms * the parameters just read. */ void - parse_parameters(ParameterHandler ¶m); + parse_parameters(const ParameterHandler ¶m); /** * Return the left end of the time interval. diff --git a/include/deal.II/base/parallel.h b/include/deal.II/base/parallel.h index 1495f62105..8c20c92583 100644 --- a/include/deal.II/base/parallel.h +++ b/include/deal.II/base/parallel.h @@ -674,7 +674,8 @@ namespace parallel * again. */ void - release_one_partitioner(std::shared_ptr &p); + release_one_partitioner( + const std::shared_ptr &p); private: /** diff --git a/include/deal.II/base/table_handler.h b/include/deal.II/base/table_handler.h index 5a7d502b3b..506e20afa6 100644 --- a/include/deal.II/base/table_handler.h +++ b/include/deal.II/base/table_handler.h @@ -961,7 +961,7 @@ TableHandler::add_value(const std::string &key, const T value) while (columns[key].entries.size() + 1 < max_col_length) { columns[key].entries.push_back(internal::TableEntry(T())); - internal::TableEntry &entry = columns[key].entries.back(); + const internal::TableEntry &entry = columns[key].entries.back(); entry.cache_string(columns[key].scientific, columns[key].precision); columns[key].max_length = std::max(columns[key].max_length, @@ -972,7 +972,7 @@ TableHandler::add_value(const std::string &key, const T value) // now push the value given to this function columns[key].entries.push_back(internal::TableEntry(value)); - internal::TableEntry &entry = columns[key].entries.back(); + const internal::TableEntry &entry = columns[key].entries.back(); entry.cache_string(columns[key].scientific, columns[key].precision); columns[key].max_length = std::max(columns[key].max_length, diff --git a/include/deal.II/base/thread_management.h b/include/deal.II/base/thread_management.h index b69d4611a3..f130ef0c95 100644 --- a/include/deal.II/base/thread_management.h +++ b/include/deal.II/base/thread_management.h @@ -393,7 +393,7 @@ namespace Threads struct maybe_make_ref { static T - act(T &t) + act(const T &t) { return t; } @@ -1432,7 +1432,7 @@ namespace Threads void join_all() const { - for (auto &t : tasks) + for (const auto &t : tasks) t.join(); } diff --git a/include/deal.II/distributed/cell_data_transfer.h b/include/deal.II/distributed/cell_data_transfer.h index 8614979328..8b2e1bef2f 100644 --- a/include/deal.II/distributed/cell_data_transfer.h +++ b/include/deal.II/distributed/cell_data_transfer.h @@ -302,8 +302,8 @@ namespace parallel cell_iterator &cell, const CellStatus status, const boost::iterator_range::const_iterator> - &data_range, - std::vector &all_out); + &data_range, + const std::vector &all_out); }; } // namespace distributed } // namespace parallel diff --git a/include/deal.II/distributed/cell_data_transfer.templates.h b/include/deal.II/distributed/cell_data_transfer.templates.h index ebbb1a856b..7ceffefe97 100644 --- a/include/deal.II/distributed/cell_data_transfer.templates.h +++ b/include/deal.II/distributed/cell_data_transfer.templates.h @@ -308,8 +308,8 @@ namespace parallel cell_iterator &cell, const CellStatus status, const boost::iterator_range::const_iterator> - &data_range, - std::vector &all_out) + &data_range, + const std::vector &all_out) { std::vector cell_data; diff --git a/include/deal.II/fe/fe_system.h b/include/deal.II/fe/fe_system.h index 68c1c724db..cc8ffe3716 100644 --- a/include/deal.II/fe/fe_system.h +++ b/include/deal.II/fe/fe_system.h @@ -572,11 +572,10 @@ public: */ FESystem(FESystem &&other_fe_system) noexcept : FiniteElement(std::move(other_fe_system)) - { - base_elements = std::move(other_fe_system.base_elements); - generalized_support_points_index_table = - std::move(other_fe_system.generalized_support_points_index_table); - } + , base_elements(std::move(other_fe_system.base_elements)) + , generalized_support_points_index_table( + std::move(other_fe_system.generalized_support_points_index_table)) + {} /** * Destructor. diff --git a/include/deal.II/grid/tria.h b/include/deal.II/grid/tria.h index d92c8ff735..7973e4667f 100644 --- a/include/deal.II/grid/tria.h +++ b/include/deal.II/grid/tria.h @@ -4724,7 +4724,7 @@ void Triangulation::load(Archive &ar, const unsigned int) // they are easy enough to rebuild upon re-loading data. do // this here. don't forget to first resize the fields appropriately { - for (auto &level : levels) + for (const auto &level : levels) { level->active_cell_indices.resize(level->refine_flags.size()); level->global_active_cell_indices.resize(level->refine_flags.size()); diff --git a/include/deal.II/lac/sparsity_pattern.h b/include/deal.II/lac/sparsity_pattern.h index 842cce0b5d..7b0e93d358 100644 --- a/include/deal.II/lac/sparsity_pattern.h +++ b/include/deal.II/lac/sparsity_pattern.h @@ -1408,14 +1408,15 @@ SparsityPattern::operator==(const SparsityPattern &sp2) const return false; if (rows > 0) - for (size_type i = 0; i < rows + 1; ++i) - if (rowstart[i] != sp2.rowstart[i]) - return false; + { + for (size_type i = 0; i < rows + 1; ++i) + if (rowstart[i] != sp2.rowstart[i]) + return false; - if (rows > 0) - for (size_type i = 0; i < rowstart[rows]; ++i) - if (colnums[i] != sp2.colnums[i]) - return false; + for (size_type i = 0; i < rowstart[rows]; ++i) + if (colnums[i] != sp2.colnums[i]) + return false; + } return true; } diff --git a/source/algorithms/timestep_control.cc b/source/algorithms/timestep_control.cc index 8f7c66cb1f..9425cda8bc 100644 --- a/source/algorithms/timestep_control.cc +++ b/source/algorithms/timestep_control.cc @@ -60,7 +60,7 @@ namespace Algorithms void - TimestepControl::parse_parameters(ParameterHandler ¶m) + TimestepControl::parse_parameters(const ParameterHandler ¶m) { start(param.get_double("Start")); start_step(param.get_double("First step")); diff --git a/source/base/convergence_table.cc b/source/base/convergence_table.cc index 24387162c3..5647cc68eb 100644 --- a/source/base/convergence_table.cc +++ b/source/base/convergence_table.cc @@ -36,8 +36,9 @@ ConvergenceTable::evaluate_convergence_rates( // the top that don't yet exist set_auto_fill_mode(false); - std::vector &entries = columns[data_column_key].entries; - std::vector &ref_entries = + const std::vector &entries = + columns[data_column_key].entries; + const std::vector &ref_entries = columns[reference_column_key].entries; std::string rate_key = data_column_key + "..."; @@ -132,8 +133,9 @@ ConvergenceTable::evaluate_convergence_rates(const std::string &data_column_key, // the top that don't yet exist set_auto_fill_mode(false); - std::vector &entries = columns[data_column_key].entries; - std::string rate_key = data_column_key + "..."; + const std::vector &entries = + columns[data_column_key].entries; + std::string rate_key = data_column_key + "..."; const unsigned int n = entries.size(); diff --git a/source/base/exceptions.cc b/source/base/exceptions.cc index 50cd791f1d..9b32898f14 100644 --- a/source/base/exceptions.cc +++ b/source/base/exceptions.cc @@ -285,7 +285,7 @@ ExceptionBase::print_stack_trace(std::ostream &out) const std::string functionname = stacktrace_entry.substr(pos_start + 1, pos_end - pos_start - 1); - stacktrace_entry = stacktrace_entry.substr(0, pos_start); + stacktrace_entry.resize(pos_start); stacktrace_entry += ": "; // demangle, and if successful replace old mangled string by diff --git a/source/base/mpi_compute_index_owner_internal.cc b/source/base/mpi_compute_index_owner_internal.cc index 651aebf7f9..c7825de990 100644 --- a/source/base/mpi_compute_index_owner_internal.cc +++ b/source/base/mpi_compute_index_owner_internal.cc @@ -114,10 +114,8 @@ namespace Utilities } else { - if (data_map.find(index) == data_map.end()) - data_map[index] = invalid_index_value; - - return data_map[index]; + return data_map.try_emplace(index, invalid_index_value) + .first->second; } } diff --git a/source/base/mu_parser_internal.cc b/source/base/mu_parser_internal.cc index f8d32dae9b..739a56d01e 100644 --- a/source/base/mu_parser_internal.cc +++ b/source/base/mu_parser_internal.cc @@ -167,10 +167,9 @@ namespace internal // which is initialized with the seed itself static std::map rng_map; - if (rng_map.find(seed) == rng_map.end()) - rng_map[seed] = std::mt19937(static_cast(seed)); - - return uniform_distribution(rng_map[seed]); + return uniform_distribution( + rng_map.try_emplace(seed, std::mt19937(static_cast(seed))) + .first->second); } diff --git a/source/base/parallel.cc b/source/base/parallel.cc index d426446b53..7ccafffb2c 100644 --- a/source/base/parallel.cc +++ b/source/base/parallel.cc @@ -85,7 +85,7 @@ namespace parallel void TBBPartitioner::release_one_partitioner( - std::shared_ptr &p) + const std::shared_ptr &p) { if (p.get() == my_partitioner.get()) { diff --git a/source/base/parameter_acceptor.cc b/source/base/parameter_acceptor.cc index ebaf925239..7610532d8e 100644 --- a/source/base/parameter_acceptor.cc +++ b/source/base/parameter_acceptor.cc @@ -59,12 +59,7 @@ ParameterAcceptor::ParameterAcceptor(const std::string &name) ParameterAcceptor::~ParameterAcceptor() { std::lock_guard l(class_list_mutex); - // Notice that it is possible that the class is no longer in the static list. - // This happens when the clear() method has been called. We must guard that - // this is not the case, and therefore we only remove ourselves from the list - // if we are actually there. - if (class_list.find(this) != class_list.end()) - class_list.erase(this); + class_list.erase(this); } @@ -198,7 +193,7 @@ ParameterAcceptor::get_section_path() const acceptor_it != class_list.rend(); ++acceptor_it) { - auto *const acceptor = *acceptor_it; + const auto *const acceptor = *acceptor_it; if (acceptor->get_acceptor_id() >= get_acceptor_id()) continue; bool has_trailing = acceptor->get_section_name().back() == sep; diff --git a/source/base/parameter_handler.cc b/source/base/parameter_handler.cc index fefa8fc96a..15ccac94b2 100644 --- a/source/base/parameter_handler.cc +++ b/source/base/parameter_handler.cc @@ -315,7 +315,7 @@ namespace current_section.sort(compare); // Now transverse subsections tree recursively. - for (auto &p : current_section) + for (const auto &p : current_section) { if ((is_parameter_node(p.second) == false) && (is_alias_node(p.second) == false)) diff --git a/source/base/polynomials_barycentric.cc b/source/base/polynomials_barycentric.cc index a132c5898a..4fa894813e 100644 --- a/source/base/polynomials_barycentric.cc +++ b/source/base/polynomials_barycentric.cc @@ -156,9 +156,8 @@ BarycentricPolynomials::BarycentricPolynomials( const std::vector &polynomials) : ScalarPolynomialsBase(internal::get_degree(polynomials), polynomials.size()) + , polys(polynomials) { - polys = polynomials; - poly_grads.resize(polynomials.size()); poly_hessians.resize(polynomials.size()); poly_third_derivatives.resize(polynomials.size()); diff --git a/source/base/quadrature_lib.cc b/source/base/quadrature_lib.cc index 577ef8e48b..281b35719e 100644 --- a/source/base/quadrature_lib.cc +++ b/source/base/quadrature_lib.cc @@ -1351,13 +1351,13 @@ QGaussRadauChebyshev<1>::QGaussRadauChebyshev(const unsigned int n, Assert(n > 0, ExcMessage("Need at least one point for quadrature rules.")); std::vector points = internal::QGaussRadauChebyshev::get_quadrature_points(n, end_point); - std::vector weights = + std::vector new_weights = internal::QGaussRadauChebyshev::get_quadrature_weights(n, end_point); for (unsigned int i = 0; i < this->size(); ++i) { this->quadrature_points[i] = Point<1>(points[i]); - this->weights[i] = weights[i]; + this->weights[i] = new_weights[i]; } } diff --git a/source/base/table_handler.cc b/source/base/table_handler.cc index 6207e6426d..bb32bbc9ff 100644 --- a/source/base/table_handler.cc +++ b/source/base/table_handler.cc @@ -90,7 +90,7 @@ namespace internal else ss.setf(std::ios::fixed, std::ios::floatfield); - std::visit([&ss](auto &v) { ss << v; }, value); + std::visit([&ss](const auto &v) { ss << v; }, value); cached_value = ss.str(); if (cached_value.empty()) @@ -155,7 +155,7 @@ TableHandler::Column::pad_column_below(const unsigned int size) while (entries.size() < size) { entries.push_back(entries.back().get_default_constructed_copy()); - internal::TableEntry &entry = entries.back(); + const internal::TableEntry &entry = entries.back(); entry.cache_string(scientific, precision); max_length = std::max(max_length, @@ -218,7 +218,7 @@ TableHandler::start_new_row() while (column.second.entries.size() < max_col_length) { column.second.entries.emplace_back(""); - internal::TableEntry &entry = column.second.entries.back(); + const internal::TableEntry &entry = column.second.entries.back(); entry.cache_string(column.second.scientific, column.second.precision); column.second.max_length = std::max(column.second.max_length, @@ -701,7 +701,8 @@ TableHandler::write_tex(std::ostream &out, const bool with_header) const else out.setf(std::ios::fixed, std::ios::floatfield); - std::visit([&out](auto &v) { out << v; }, column.entries[i].value); + std::visit([&out](const auto &v) { out << v; }, + column.entries[i].value); if (j < n_cols - 1) out << " & "; diff --git a/source/base/utilities.cc b/source/base/utilities.cc index 003a0a0ce8..4485574b33 100644 --- a/source/base/utilities.cc +++ b/source/base/utilities.cc @@ -1012,8 +1012,8 @@ namespace Utilities std::string get_time() { - std::time_t time1 = std::time(nullptr); - std::tm *time = std::localtime(&time1); + std::time_t time1 = std::time(nullptr); + const std::tm *time = std::localtime(&time1); std::ostringstream o; o << time->tm_hour << ":" << (time->tm_min < 10 ? "0" : "") @@ -1028,8 +1028,8 @@ namespace Utilities std::string get_date() { - std::time_t time1 = std::time(nullptr); - std::tm *time = std::localtime(&time1); + std::time_t time1 = std::time(nullptr); + const std::tm *time = std::localtime(&time1); std::ostringstream o; o << time->tm_year + 1900 << "/" << time->tm_mon + 1 << "/" diff --git a/source/dofs/dof_handler.cc b/source/dofs/dof_handler.cc index af18a2aaf9..8c2baebda2 100644 --- a/source/dofs/dof_handler.cc +++ b/source/dofs/dof_handler.cc @@ -56,21 +56,13 @@ namespace internal { std::string policy_name; if (dynamic_cast *>(&policy) || - dynamic_cast *>(&policy)) policy_name = "Policy::Sequential<"; else if (dynamic_cast< - const typename dealii::internal::DoFHandlerImplementation:: - Policy::ParallelDistributed *>(&policy) || - dynamic_cast< const typename dealii::internal::DoFHandlerImplementation:: Policy::ParallelDistributed *>(&policy)) policy_name = "Policy::ParallelDistributed<"; else if (dynamic_cast< - const typename dealii::internal::DoFHandlerImplementation:: - Policy::ParallelShared *>(&policy) || - dynamic_cast< const typename dealii::internal::DoFHandlerImplementation:: Policy::ParallelShared *>(&policy)) policy_name = "Policy::ParallelShared<"; diff --git a/source/dofs/dof_tools.cc b/source/dofs/dof_tools.cc index 3fd1a6641c..f1809ebd2c 100644 --- a/source/dofs/dof_tools.cc +++ b/source/dofs/dof_tools.cc @@ -489,7 +489,7 @@ namespace DoFTools index_per_comp[comp_i].add_index( locally_owned_dofs.nth_index_in_set(i)); } - for (auto &c : index_per_comp) + for (const auto &c : index_per_comp) c.compress(); return index_per_comp; } diff --git a/source/fe/fe_dgq.cc b/source/fe/fe_dgq.cc index 83ea9c889a..ffd2dc6a02 100644 --- a/source/fe/fe_dgq.cc +++ b/source/fe/fe_dgq.cc @@ -870,7 +870,7 @@ FE_DGQArbitraryNodes::get_name() const bool equidistant = true; std::vector points(this->degree + 1); - auto *const polynomial_space = + const auto *const polynomial_space = dynamic_cast *>(this->poly_space.get()); Assert(polynomial_space != nullptr, ExcInternalError()); std::vector lexicographic = @@ -987,7 +987,7 @@ FE_DGQArbitraryNodes::clone() const { // Construct a dummy quadrature formula containing the FE's nodes: std::vector> qpoints(this->degree + 1); - auto *const polynomial_space = + const auto *const polynomial_space = dynamic_cast *>(this->poly_space.get()); Assert(polynomial_space != nullptr, ExcInternalError()); std::vector lexicographic = -- 2.39.5