]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Delete code paths for Kokkos versions before 3.4 18289/head
authorDaniel Arndt <arndtd@ornl.gov>
Wed, 26 Mar 2025 13:32:57 +0000 (09:32 -0400)
committerDaniel Arndt <arndtd@ornl.gov>
Wed, 26 Mar 2025 13:32:57 +0000 (09:32 -0400)
include/deal.II/lac/la_parallel_vector.templates.h
include/deal.II/lac/vector_operations_internal.h
include/deal.II/matrix_free/portable_matrix_free.templates.h
include/deal.II/matrix_free/portable_tensor_product_kernels.h
source/base/exceptions.cc

index d6f3aa34a7d67e2c2bd4adc8b54dba2698bc95df..ef8fcdab85ca7fa6b77f19c6254dbadbed4f446d 100644 (file)
@@ -492,9 +492,7 @@ namespace LinearAlgebra
               ::dealii::MemorySpace::Default::kokkos_space::execution_space>(
               exec, 0, size),
             KOKKOS_LAMBDA(size_type i, RealType & update) {
-#if KOKKOS_VERSION < 30400
-              update = fmax(update, fabs(data.values(i)));
-#elif KOKKOS_VERSION < 30700
+#if KOKKOS_VERSION < 30700
               update = Kokkos::Experimental::fmax(
                 update, Kokkos::Experimental::fabs(data.values(i)));
 #else
index 65bb9a15f09583ccbabdece9f14d79df8f32d80f..7760ac83e99e1f1c76d0f6b8322d26ef5a941157 100644 (file)
@@ -2524,9 +2524,7 @@ namespace internal
             ::dealii::MemorySpace::Default::kokkos_space::execution_space>(
             exec, optional_offset, optional_offset + size),
           KOKKOS_LAMBDA(size_type i, Number & update) {
-#if KOKKOS_VERSION < 30400
-            update += std::abs(data.values(i));
-#elif KOKKOS_VERSION < 30700
+#if KOKKOS_VERSION < 30700
             update += Kokkos::Experimental::fabs(data.values(i));
 #else
             update += Kokkos::abs(data.values(i));
@@ -2554,9 +2552,7 @@ namespace internal
             ::dealii::MemorySpace::Default::kokkos_space::execution_space>(
             exec, 0, size),
           KOKKOS_LAMBDA(size_type i, Number & update) {
-#if KOKKOS_VERSION < 30400
-            update += std::pow(fabs(data.values(i)), exp);
-#elif KOKKOS_VERSION < 30700
+#if KOKKOS_VERSION < 30700
             update += Kokkos::Experimental::pow(
               Kokkos::Experimental::fabs(data.values(i)), exp);
 #else
index 2cb7acb9db9958844586de450833ca9e49234948..0070ff58761366f807aad1ad998f373344a1c195 100644 (file)
@@ -664,12 +664,7 @@ namespace Portable
           MemorySpace::Default::kokkos_space::execution_space exec;
           Kokkos::TeamPolicy<
             MemorySpace::Default::kokkos_space::execution_space>
-            team_policy(
-#if KOKKOS_VERSION >= 20900
-              exec,
-#endif
-              n_cells[i],
-              Kokkos::AUTO);
+            team_policy(exec, n_cells[i], Kokkos::AUTO);
 
           Kokkos::parallel_for(
             "dealii::MatrixFree::evaluate_coeff_cell_loop",
@@ -677,22 +672,18 @@ namespace Portable
             KOKKOS_LAMBDA(const Kokkos::TeamPolicy<
                           MemorySpace::Default::kokkos_space::execution_space>::
                             member_type &team_member) {
-              Kokkos::parallel_for(
-#if KOKKOS_VERSION >= 20900
-                Kokkos::TeamVectorRange(team_member, n_q_points),
-#else
-                Kokkos::TeamThreadRange(team_member, n_q_points),
-#endif
-                [&](const int q_point) {
-                  const int cell = team_member.league_rank();
-
-                  Data data{team_member,
-                            cell,
-                            &color_data,
-                            /*shared_data*/ nullptr};
-
-                  func(&data, cell, q_point);
-                });
+              Kokkos::parallel_for(Kokkos::TeamVectorRange(team_member,
+                                                           n_q_points),
+                                   [&](const int q_point) {
+                                     const int cell = team_member.league_rank();
+
+                                     Data data{team_member,
+                                               cell,
+                                               &color_data,
+                                               /*shared_data*/ nullptr};
+
+                                     func(&data, cell, q_point);
+                                   });
             });
         }
   }
@@ -1037,12 +1028,7 @@ namespace Portable
           MemorySpace::Default::kokkos_space::execution_space exec;
           Kokkos::TeamPolicy<
             MemorySpace::Default::kokkos_space::execution_space>
-            team_policy(
-#if KOKKOS_VERSION >= 20900
-              exec,
-#endif
-              n_cells[color],
-              Kokkos::AUTO);
+            team_policy(exec, n_cells[color], Kokkos::AUTO);
 
 
           internal::
@@ -1099,12 +1085,7 @@ namespace Portable
               {
                 Kokkos::TeamPolicy<
                   MemorySpace::Default::kokkos_space::execution_space>
-                  team_policy(
-#if KOKKOS_VERSION >= 20900
-                    exec,
-#endif
-                    n_cells[0],
-                    Kokkos::AUTO);
+                  team_policy(exec, n_cells[0], Kokkos::AUTO);
 
                 internal::ApplyKernel<dim, Number, Functor, false> apply_kernel(
                   func, get_data(0), src, dst);
@@ -1122,12 +1103,7 @@ namespace Portable
               {
                 Kokkos::TeamPolicy<
                   MemorySpace::Default::kokkos_space::execution_space>
-                  team_policy(
-#if KOKKOS_VERSION >= 20900
-                    exec,
-#endif
-                    n_cells[1],
-                    Kokkos::AUTO);
+                  team_policy(exec, n_cells[1], Kokkos::AUTO);
 
                 internal::ApplyKernel<dim, Number, Functor, false> apply_kernel(
                   func, get_data(1), src, dst);
@@ -1150,12 +1126,7 @@ namespace Portable
               {
                 Kokkos::TeamPolicy<
                   MemorySpace::Default::kokkos_space::execution_space>
-                  team_policy(
-#if KOKKOS_VERSION >= 20900
-                    exec,
-#endif
-                    n_cells[2],
-                    Kokkos::AUTO);
+                  team_policy(exec, n_cells[2], Kokkos::AUTO);
 
                 internal::ApplyKernel<dim, Number, Functor, false> apply_kernel(
                   func, get_data(2), src, dst);
@@ -1183,12 +1154,7 @@ namespace Portable
                 {
                   Kokkos::TeamPolicy<
                     MemorySpace::Default::kokkos_space::execution_space>
-                    team_policy(
-#if KOKKOS_VERSION >= 20900
-                      exec,
-#endif
-                      n_cells[i],
-                      Kokkos::AUTO);
+                    team_policy(exec, n_cells[i], Kokkos::AUTO);
 
                   internal::ApplyKernel<dim, Number, Functor, false>
                     apply_kernel(func, get_data(i), src, dst);
@@ -1220,12 +1186,7 @@ namespace Portable
             {
               Kokkos::TeamPolicy<
                 MemorySpace::Default::kokkos_space::execution_space>
-                team_policy(
-#if KOKKOS_VERSION >= 20900
-                  exec,
-#endif
-                  n_cells[i],
-                  Kokkos::AUTO);
+                team_policy(exec, n_cells[i], Kokkos::AUTO);
 
               internal::ApplyKernel<dim, Number, Functor, false> apply_kernel(
                 func, get_data(i), ghosted_src, ghosted_dst);
index dbf673b78d6a9c65983243ea2a2edab2562b8db5..27a7c168fe5d8b1f9c433c833accadb901f05cea 100644 (file)
@@ -58,18 +58,13 @@ namespace Portable
     {
       Assert(dst.size() >= N, ExcInternalError());
       Assert(src.size() >= N, ExcInternalError());
-      Kokkos::parallel_for(
-#if KOKKOS_VERSION >= 20900
-        Kokkos::TeamVectorRange(team_member, N),
-#else
-        Kokkos::TeamThreadRange(team_member, N),
-#endif
-        [&](const int i) {
-          if constexpr (add)
-            Kokkos::atomic_add(&dst(i), src(i));
-          else
-            dst(i) = src(i);
-        });
+      Kokkos::parallel_for(Kokkos::TeamVectorRange(team_member, N),
+                           [&](const int i) {
+                             if constexpr (add)
+                               Kokkos::atomic_add(&dst(i), src(i));
+                             else
+                               dst(i) = src(i);
+                           });
 
       team_member.team_barrier();
     }
index 3518861656e257d362a8c0f5c67040277420317a..aa4399a07a5060bd69e53f133b26f018ef56866b 100644 (file)
@@ -520,25 +520,11 @@ namespace deal_II_exceptions
         }
 #endif
 
-        // Let's abort the program here. On the host, we need to call
-        // std::abort, on devices we need to do something different.
-        // Kokkos::abort() does the right thing in all circumstances.
-
-#if KOKKOS_VERSION < 30200
-      if constexpr (std::is_same_v<Kokkos::DefaultExecutionSpace,
-                                   Kokkos::DefaultHostExecutionSpace>)
-        {
-          // FIXME_KOKKOS Older Kokkos versions don't declare Kokkos::abort as
-          // [[noreturn]]. In case Kokkos is only configured with host backends,
-          // we can just use std::abort instead.
-          std::abort();
-        }
-      else
-#endif
-        {
-          Kokkos::abort(
-            "Abort() was called during dealing with an assertion or exception.");
-        }
+      // Let's abort the program here. On the host, we need to call
+      // std::abort, on devices we need to do something different.
+      // Kokkos::abort() does the right thing in all circumstances.
+      Kokkos::abort(
+        "Abort() was called during dealing with an assertion or exception.");
     }
 
 

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.