]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Cleanup commits 17617/head
authorRene Gassmoeller <rene.gassmoeller@mailbox.org>
Tue, 27 Aug 2024 14:30:57 +0000 (16:30 +0200)
committerRene Gassmoeller <rene.gassmoeller@mailbox.org>
Tue, 27 Aug 2024 15:00:55 +0000 (17:00 +0200)
source/base/function_lib.cc
source/base/mpi.cc
source/grid/tria.cc
tests/matrix_free/matrix_vector_hessians_cells.cc

index 2e2ec9e1ee47edf8b6171254bdb6e20cd9b955b1..d64a9271839c3db4fd9390ddc09790e6e1b310c2 100644 (file)
@@ -2402,10 +2402,9 @@ namespace Functions
     Tensor<1, 1>
     gradient_interpolate(const Table<1, double> &data_values,
                          const TableIndices<1>  &ix,
-                         const Point<1>         &p_unit,
-                         const Point<1>         &dx)
+                         const Point<1> & /*p_unit*/,
+                         const Point<1> &dx)
     {
-      [[mabe_unused]]p_unit;
       Tensor<1, 1> grad;
       grad[0] = (data_values[ix[0] + 1] - data_values[ix[0]]) / dx[0];
       return grad;
index 68e3dae9bfbc0b17963cefef8a89b5e9dc290156..7857ec2bfb8299e71f03a6afab4e7c8cf150ce15 100644 (file)
@@ -794,8 +794,6 @@ namespace Utilities
     void
     CollectiveMutex::lock(const MPI_Comm comm)
     {
-      [[mabe_unused]]comm;
-
       Assert(
         !locked,
         ExcMessage(
@@ -821,6 +819,8 @@ namespace Utilities
           // nothing to do as blocking barrier already completed
 #  endif
         }
+#else
+      (void)comm;
 #endif
 
       locked = true;
@@ -831,8 +831,6 @@ namespace Utilities
     void
     CollectiveMutex::unlock(const MPI_Comm comm)
     {
-      [[mabe_unused]]comm;
-
       // First check if this function is called during exception handling
       // if so, abort. This can happen if a ScopedLock is destroyed.
       internal::CollectiveMutexImplementation::check_exception();
@@ -857,6 +855,8 @@ namespace Utilities
           AssertThrowMPI(ierr);
 #  endif
         }
+#else
+      (void)comm;
 #endif
 
       locked = false;
index 210ddf0453155eb7cc815c2a84c2c2eb0d978888..0ceb2bd635695a533ebe143eb8bac6e3a04e84fc 100644 (file)
@@ -669,10 +669,10 @@ namespace internal
   template <int dim, int spacedim>
   DEAL_II_CXX20_REQUIRES((concepts::is_valid_dim_spacedim<dim, spacedim>))
   void CellAttachedDataSerializer<dim, spacedim>::save(
-    [[maybe_unused]] const unsigned int global_first_cell,
-    [[maybe_unused]] const unsigned int global_num_cells,
-    const std::string                  &file_basename,
-    [[maybe_unused]] const MPI_Comm    &mpi_communicator) const
+    const unsigned int global_first_cell,
+    const unsigned int global_num_cells,
+    const std::string &file_basename,
+    const MPI_Comm    &mpi_communicator) const
   {
     Assert(sizes_fixed_cumulative.size() > 0,
            ExcMessage("No data has been packed!"));
@@ -848,6 +848,10 @@ namespace internal
     else
 #endif
       {
+        (void)global_first_cell;
+        (void)global_num_cells;
+        (void)mpi_communicator;
+
         //
         // ---------- Fixed size data ----------
         //
@@ -896,13 +900,13 @@ namespace internal
   template <int dim, int spacedim>
   DEAL_II_CXX20_REQUIRES((concepts::is_valid_dim_spacedim<dim, spacedim>))
   void CellAttachedDataSerializer<dim, spacedim>::load(
-    [[maybe_unused]] const unsigned int global_first_cell,
-    [[maybe_unused]] const unsigned int global_num_cells,
-    const unsigned int                  local_num_cells,
-    const std::string                  &file_basename,
-    const unsigned int                  n_attached_deserialize_fixed,
-    const unsigned int                  n_attached_deserialize_variable,
-    const MPI_Comm                     &mpi_communicator)
+    const unsigned int global_first_cell,
+    const unsigned int global_num_cells,
+    const unsigned int local_num_cells,
+    const std::string &file_basename,
+    const unsigned int n_attached_deserialize_fixed,
+    const unsigned int n_attached_deserialize_variable,
+    const MPI_Comm    &mpi_communicator)
   {
     Assert(dest_data_fixed.empty(),
            ExcMessage("Previously loaded data has not been released yet!"));
@@ -1060,6 +1064,8 @@ namespace internal
 #endif
       {
         (void)mpi_communicator;
+        (void)global_first_cell;
+        (void)global_num_cells;
 
         //
         // ---------- Fixed size data ----------
@@ -3473,7 +3479,7 @@ namespace internal
                     return a.vertices < b.vertices;
                   });
 
-        unsigned int counter = 0;
+        [[maybe_unused]] unsigned int counter = 0;
 
         std::vector<unsigned int> key;
         key.reserve(GeometryInfo<structdim>::vertices_per_cell);
@@ -3519,7 +3525,6 @@ namespace internal
             if (subcell_object->boundary_id !=
                 numbers::internal_face_boundary_id)
               {
-                [[mabe_unused]]vertex_locations;
                 AssertThrow(
                   boundary_id != numbers::internal_face_boundary_id,
                   ExcMessage(
@@ -3549,7 +3554,6 @@ namespace internal
         // make sure that all subcelldata entries have been processed
         // TODO: this is not guaranteed, why?
         // AssertDimension(counter, boundary_objects_in.size());
-        [[mabe_unused]]counter;
       }
 
 
@@ -5038,7 +5042,6 @@ namespace internal
                 triangulation.vertices[next_unused_vertex] = line->center(true);
 
                 bool pair_found = false;
-                [[mabe_unused]]pair_found;
                 for (; next_unused_line != endl; ++next_unused_line)
                   if (!next_unused_line->used() &&
                       !(++next_unused_line)->used())
@@ -5808,7 +5811,6 @@ namespace internal
                 // two child lines.  To this end, find a pair of
                 // unused lines
                 bool pair_found = false;
-                [[mabe_unused]]pair_found;
                 for (; next_unused_line != endl; ++next_unused_line)
                   if (!next_unused_line->used() &&
                       !(++next_unused_line)->used())
@@ -6230,10 +6232,7 @@ namespace internal
                 }
 
               for (const unsigned int line : quad->line_indices())
-                {
-                  AssertIsNotUsed(new_lines[line]);
-                  [[mabe_unused]]line;
-                }
+                AssertIsNotUsed(new_lines[line]);
 
               // 2) create new quads (properties are set below). Both triangles
               // and quads are divided in four.
@@ -6252,10 +6251,7 @@ namespace internal
               quad->set_refinement_case(RefinementCase<2>::cut_xy);
 
               for (const auto &quad : new_quads)
-                {
-                  AssertIsNotUsed(quad);
-                  [[mabe_unused]]quad;
-                }
+                AssertIsNotUsed(quad);
 
               // 3) set vertex indices and set new vertex
 
@@ -11928,16 +11924,12 @@ namespace internal
       template <int dim, int spacedim>
       static void
       delete_children(
-        Triangulation<dim, spacedim>                         &triangulation,
-        typename Triangulation<dim, spacedim>::cell_iterator &cell,
-        std::vector<unsigned int>                            &line_cell_count,
-        std::vector<unsigned int>                            &quad_cell_count)
+        Triangulation<dim, spacedim> & /*triangulation*/,
+        typename Triangulation<dim, spacedim>::cell_iterator & /*cell*/,
+        std::vector<unsigned int> & /*line_cell_count*/,
+        std::vector<unsigned int> & /*quad_cell_count*/)
       {
         AssertThrow(false, ExcNotImplemented());
-        [[mabe_unused]]triangulation;
-        [[mabe_unused]]cell;
-        [[mabe_unused]]line_cell_count;
-        [[mabe_unused]]quad_cell_count;
       }
 
       template <int dim, int spacedim>
@@ -11952,10 +11944,9 @@ namespace internal
       template <int dim, int spacedim>
       static void
       prevent_distorted_boundary_cells(
-        Triangulation<dim, spacedim> &triangulation)
+        Triangulation<dim, spacedim> & /*triangulation*/)
       {
         // nothing to do since anisotropy is not supported
-        [[mabe_unused]]triangulation;
       }
 
       template <int dim, int spacedim>
@@ -12293,7 +12284,6 @@ void Triangulation<dim, spacedim>::set_all_manifold_ids_on_boundary(
             }
     }
 
-  [[mabe_unused]]boundary_found;
   Assert(boundary_found, ExcBoundaryIdNotFound(b_id));
 }
 
@@ -16304,7 +16294,6 @@ void Triangulation<dim, spacedim>::load_attached_data(
       // CellStatus::cell_will_persist.
       for (const auto &cell_rel : this->local_cell_relations)
         {
-          [[mabe_unused]]cell_rel;
           Assert((cell_rel.second == // cell_status
                   ::dealii::CellStatus::cell_will_persist),
                  ExcInternalError());
index d8f4340a8a056bb6a2cc2fec9d1e6191f8238e79..c40b13567339671680634910e0567c3a51649c56 100644 (file)
@@ -68,15 +68,14 @@ test_hessians(const dealii::FE_Poly<dim>                    &fe,
   additional_data.mapping_update_flags_inner_faces =
     update_values | update_gradients | update_hessians;
 
-  [[maybe_unused]]MatrixFree<dim, double, VectorizedArrayType> matrix_free;
+  MatrixFree<dim, double, VectorizedArrayType> matrix_free;
   matrix_free.reinit(mapping, dof_handler, constraints, quad, additional_data);
 
   if (Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0)
     deallog << "Working with " << fe.get_name() << " and "
             << dof_handler.n_dofs() << " dofs" << std::endl;
 
-  LinearAlgebra::distributed::Vector<double> dst2;
-  [[maybe_unused]]LinearAlgebra::distributed::Vector<double> src, dst;
+  LinearAlgebra::distributed::Vector<double> src, dst, dst2;
   matrix_free.initialize_dof_vector(src);
   for (auto &v : src)
     v = random_value<double>();
@@ -102,7 +101,7 @@ test_hessians(const dealii::FE_Poly<dim>                    &fe,
   matrix_free.template loop<LinearAlgebra::distributed::Vector<double>,
                             LinearAlgebra::distributed::Vector<double>>(
     [&](
-      const auto &matrix_free, auto &dst, const auto &src, [[maybe_unused]]const auto &range) {
+      const auto &matrix_free, auto &dst, const auto &src, const auto &range) {
       FEEvaluation<dim, -1, 0, 1, double, VectorizedArrayType> fe_eval(
         matrix_free);
       for (unsigned int cell = range.first; cell < range.second; ++cell)
@@ -184,12 +183,14 @@ test_hessians(const dealii::FE_Poly<dim>                    &fe,
             }
         }
     },
-    [&](
-      const auto &matrix_free, auto &dst, const auto &src, const auto &range) {
-    },
-    [&](
-      const auto &matrix_free, auto &dst, const auto &src, const auto &range) {
-    },
+    [&](const auto & /*matrix_free*/,
+        auto & /*dst*/,
+        const auto & /*src*/,
+        const auto & /*range*/) {},
+    [&](const auto & /*matrix_free*/,
+        auto & /*dst*/,
+        const auto & /*src*/,
+        const auto & /*range*/) {},
     dst,
     src,
     true);

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.