]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Avoid a couple more uses of the DEBUG macro. 18181/head
authorWolfgang Bangerth <bangerth@colostate.edu>
Fri, 28 Feb 2025 16:09:37 +0000 (09:09 -0700)
committerWolfgang Bangerth <bangerth@colostate.edu>
Fri, 28 Feb 2025 17:02:20 +0000 (10:02 -0700)
include/deal.II/base/partitioner.templates.h
include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h

index e582618649a093bf0321a21fcf7f24b7e6a13d91..d0f711f887f77233703acc7a7b8d59aee1ee6ee2 100644 (file)
@@ -311,10 +311,9 @@ namespace Utilities
       // inserted data, therefore the communication is still initialized.
       // Having different code in debug and optimized mode is somewhat
       // dangerous, but it really saves communication so do it anyway
-#    ifndef DEBUG
-      if (vector_operation == VectorOperation::insert)
-        return;
-#    endif
+      if constexpr (running_in_debug_mode() == false)
+        if (vector_operation == VectorOperation::insert)
+          return;
 
       // nothing to do when we neither have import
       // nor ghost indices.
@@ -545,24 +544,23 @@ namespace Utilities
 
       // in optimized mode, no communication was started, so leave the
       // function directly (and only clear ghosts)
-#    ifndef DEBUG
-      if (vector_operation == VectorOperation::insert)
-        {
-          Assert(requests.empty(),
-                 ExcInternalError(
-                   "Did not expect a non-empty communication "
-                   "request when inserting. Check that the same "
-                   "vector_operation argument was passed to "
-                   "import_from_ghosted_array_start as is passed "
-                   "to import_from_ghosted_array_finish."));
-
-          Kokkos::deep_copy(
-            Kokkos::View<Number *, typename MemorySpaceType::kokkos_space>(
-              ghost_array.data(), ghost_array.size()),
-            0);
-          return;
-        }
-#    endif
+      if constexpr (running_in_debug_mode() == false)
+        if (vector_operation == VectorOperation::insert)
+          {
+            Assert(requests.empty(),
+                   ExcInternalError(
+                     "Did not expect a non-empty communication "
+                     "request when inserting. Check that the same "
+                     "vector_operation argument was passed to "
+                     "import_from_ghosted_array_start as is passed "
+                     "to import_from_ghosted_array_finish."));
+
+            Kokkos::deep_copy(
+              Kokkos::View<Number *, typename MemorySpaceType::kokkos_space>(
+                ghost_array.data(), ghost_array.size()),
+              0);
+            return;
+          }
 
       // nothing to do when we neither have import nor ghost indices.
       if (n_ghost_indices() == 0 && n_import_indices() == 0)
index 786148df417b7e46fda029038855a2ccad3bf171..546ce431a25ba1b1629ab86ccfc497bf26bc49e4 100644 (file)
@@ -4628,86 +4628,87 @@ void
 MGTransferMF<dim, Number, MemorySpace>::assert_dof_handler(
   const DoFHandler<dim> &dof_handler_out) const
 {
-#ifndef DEBUG
   (void)dof_handler_out;
-#else
-
-  const auto dof_handler_and_level_in = get_dof_handler_fine();
-  const auto dof_handler_in           = dof_handler_and_level_in.first;
-  const auto level_in                 = dof_handler_and_level_in.second;
-
-  if ((dof_handler_out.n_dofs() == 0) ||  // dummy DoFHandler
-      (dof_handler_in == nullptr) ||      // single level
-      (dof_handler_in == &dof_handler_out // same DoFHandler
-       ))
-    return; // nothing to do
-
-  if (this->perform_plain_copy)
+  if constexpr (running_in_debug_mode())
     {
-      // global-coarsening path: compare indices of cells
+      const auto dof_handler_and_level_in = get_dof_handler_fine();
+      const auto dof_handler_in           = dof_handler_and_level_in.first;
+      const auto level_in                 = dof_handler_and_level_in.second;
+
+      if ((dof_handler_out.n_dofs() == 0) ||  // dummy DoFHandler
+          (dof_handler_in == nullptr) ||      // single level
+          (dof_handler_in == &dof_handler_out // same DoFHandler
+           ))
+        return; // nothing to do
 
-      std::vector<types::global_dof_index> dof_indices_in;
-      std::vector<types::global_dof_index> dof_indices_out;
+      if (this->perform_plain_copy)
+        {
+          // global-coarsening path: compare indices of cells
 
-      internal::loop_over_active_or_level_cells(
-        dof_handler_in->get_triangulation(), level_in, [&](const auto &cell) {
-          const auto cell_id = cell->id();
+          std::vector<types::global_dof_index> dof_indices_in;
+          std::vector<types::global_dof_index> dof_indices_out;
 
-          Assert(
-            dof_handler_out.get_triangulation().contains_cell(cell_id),
-            ExcMessage(
-              "DoFHandler instances used for set up of MGTransferMF and copy_to_mg(), "
-              "copy_from_mg(), or interpolate_to_mg() are not compatible."));
+          internal::loop_over_active_or_level_cells(
+            dof_handler_in->get_triangulation(),
+            level_in,
+            [&](const auto &cell) {
+              const auto cell_id = cell->id();
 
-          if (level_in == numbers::invalid_unsigned_int)
-            {
-              const auto cell_in =
-                cell->as_dof_handler_iterator(*dof_handler_in);
-              const auto cell_out =
-                dof_handler_out.get_triangulation()
-                  .create_cell_iterator(cell_id)
-                  ->as_dof_handler_iterator(dof_handler_out);
+              Assert(
+                dof_handler_out.get_triangulation().contains_cell(cell_id),
+                ExcMessage(
+                  "DoFHandler instances used for set up of MGTransferMF and copy_to_mg(), "
+                  "copy_from_mg(), or interpolate_to_mg() are not compatible."));
 
-              AssertDimension(cell_in->get_fe().n_dofs_per_cell(),
-                              cell_out->get_fe().n_dofs_per_cell());
+              if (level_in == numbers::invalid_unsigned_int)
+                {
+                  const auto cell_in =
+                    cell->as_dof_handler_iterator(*dof_handler_in);
+                  const auto cell_out =
+                    dof_handler_out.get_triangulation()
+                      .create_cell_iterator(cell_id)
+                      ->as_dof_handler_iterator(dof_handler_out);
 
-              dof_indices_in.resize(cell_in->get_fe().n_dofs_per_cell());
-              dof_indices_out.resize(cell_out->get_fe().n_dofs_per_cell());
+                  AssertDimension(cell_in->get_fe().n_dofs_per_cell(),
+                                  cell_out->get_fe().n_dofs_per_cell());
 
-              cell_in->get_dof_indices(dof_indices_in);
-              cell_out->get_dof_indices(dof_indices_out);
-            }
-          else
-            {
-              const auto cell_in =
-                cell->as_dof_handler_level_iterator(*dof_handler_in);
-              const auto cell_out =
-                dof_handler_out.get_triangulation()
-                  .create_cell_iterator(cell_id)
-                  ->as_dof_handler_level_iterator(dof_handler_out);
+                  dof_indices_in.resize(cell_in->get_fe().n_dofs_per_cell());
+                  dof_indices_out.resize(cell_out->get_fe().n_dofs_per_cell());
 
-              AssertDimension(cell_in->get_fe().n_dofs_per_cell(),
-                              cell_out->get_fe().n_dofs_per_cell());
+                  cell_in->get_dof_indices(dof_indices_in);
+                  cell_out->get_dof_indices(dof_indices_out);
+                }
+              else
+                {
+                  const auto cell_in =
+                    cell->as_dof_handler_level_iterator(*dof_handler_in);
+                  const auto cell_out =
+                    dof_handler_out.get_triangulation()
+                      .create_cell_iterator(cell_id)
+                      ->as_dof_handler_level_iterator(dof_handler_out);
 
-              dof_indices_in.resize(cell_in->get_fe().n_dofs_per_cell());
-              dof_indices_out.resize(cell_out->get_fe().n_dofs_per_cell());
+                  AssertDimension(cell_in->get_fe().n_dofs_per_cell(),
+                                  cell_out->get_fe().n_dofs_per_cell());
 
-              cell_in->get_mg_dof_indices(dof_indices_in);
-              cell_out->get_mg_dof_indices(dof_indices_out);
-            }
+                  dof_indices_in.resize(cell_in->get_fe().n_dofs_per_cell());
+                  dof_indices_out.resize(cell_out->get_fe().n_dofs_per_cell());
 
-          Assert(
-            dof_indices_in == dof_indices_out,
-            ExcMessage(
-              "DoFHandler instances used for set up of MGTransferMF and copy_to_mg(), "
-              "copy_from_mg(), or interpolate_to_mg() are not compatible."));
-        });
-    }
-  else if (this->perform_renumbered_plain_copy)
-    {
-      // nothing to do
+                  cell_in->get_mg_dof_indices(dof_indices_in);
+                  cell_out->get_mg_dof_indices(dof_indices_out);
+                }
+
+              Assert(
+                dof_indices_in == dof_indices_out,
+                ExcMessage(
+                  "DoFHandler instances used for set up of MGTransferMF and copy_to_mg(), "
+                  "copy_from_mg(), or interpolate_to_mg() are not compatible."));
+            });
+        }
+      else if (this->perform_renumbered_plain_copy)
+        {
+          // nothing to do
+        }
     }
-#endif
 }
 
 

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.