]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Convert AffineConstraintsImplementation::set_zero_parallel
authorDaniel Arndt <arndtd@ornl.gov>
Sat, 31 Dec 2022 17:33:10 +0000 (18:33 +0100)
committerDaniel Arndt <arndtd@ornl.gov>
Sat, 31 Dec 2022 17:40:59 +0000 (18:40 +0100)
include/deal.II/lac/affine_constraints.templates.h
source/lac/affine_constraints.cc
source/lac/affine_constraints_cuda.cc [deleted file]

index 9b2c26cf730bcbee2a5655e3ad84a2a9b84acf2d..c5d27717ff1234b0fcc0cc375193acd7dc5e2699 100644 (file)
@@ -2256,24 +2256,12 @@ namespace internal
       vec.zero_out_ghost_values();
     }
 
-#ifdef DEAL_II_WITH_CUDA
-    template <typename Number>
-    __global__ void
-    set_zero_kernel(const size_type *  constrained_dofs,
-                    const unsigned int n_constrained_dofs,
-                    Number *           dst)
-    {
-      const unsigned int index = threadIdx.x + blockDim.x * blockIdx.x;
-      if (index < n_constrained_dofs)
-        dst[constrained_dofs[index]] = 0;
-    }
-
     template <typename number>
     void
     set_zero_parallel(
-      const std::vector<size_type> &                                 cm,
-      LinearAlgebra::distributed::Vector<number, MemorySpace::CUDA> &vec,
-      size_type                                                      shift = 0)
+      const std::vector<size_type> &                                    cm,
+      LinearAlgebra::distributed::Vector<number, MemorySpace::Default> &vec,
+      size_type shift = 0)
     {
       Assert(shift == 0, ExcNotImplemented());
       (void)shift;
@@ -2285,22 +2273,30 @@ namespace internal
           constrained_local_dofs_host.push_back(
             vec.get_partitioner()->global_to_local(global_index));
 
-      const int  n_constraints = constrained_local_dofs_host.size();
-      size_type *constrained_local_dofs_device;
-      Utilities::CUDA::malloc(constrained_local_dofs_device, n_constraints);
-      Utilities::CUDA::copy_to_dev(constrained_local_dofs_host,
-                                   constrained_local_dofs_device);
-
-      const int n_blocks = 1 + (n_constraints - 1) / CUDAWrappers::block_size;
-      set_zero_kernel<<<n_blocks, CUDAWrappers::block_size>>>(
-        constrained_local_dofs_device, n_constraints, vec.get_values());
-      AssertCudaKernel();
-
-      Utilities::CUDA::free(constrained_local_dofs_device);
+      const int n_constraints = constrained_local_dofs_host.size();
+      Kokkos::View<size_type *, MemorySpace::Default::kokkos_space>
+        constrained_local_dofs_device(
+          Kokkos::view_alloc(Kokkos::WithoutInitializing,
+                             "constrained_local_dofs_device"),
+          n_constraints);
+      Kokkos::deep_copy(constrained_local_dofs_device,
+                        Kokkos::View<size_type *, Kokkos::HostSpace>(
+                          constrained_local_dofs_host.data(),
+                          constrained_local_dofs_host.size()));
+
+      using ExecutionSpace =
+        MemorySpace::Default::kokkos_space::execution_space;
+      ExecutionSpace exec;
+      auto           local_values = vec.get_values();
+      Kokkos::parallel_for(
+        Kokkos::RangePolicy<ExecutionSpace>(exec, 0, n_constraints),
+        KOKKOS_LAMBDA(int i) {
+          local_values[constrained_local_dofs_device[i]] = 0;
+        });
+      exec.fence();
 
       vec.zero_out_ghost_values();
     }
-#endif
 
     template <class VectorType>
     void
index 86ee7dd8566ac5c2d2c0c6d418b04f3a4df1edc1..440bf4a2e9c21151de3fab369c2b867ecd091a5d 100644 (file)
@@ -150,4 +150,30 @@ dealii::AffineConstraints<double>::distribute<
 #  endif
 #endif
 
+namespace internal
+{
+  namespace AffineConstraintsImplementation
+  {
+    template void
+    set_zero_all(
+      const std::vector<types::global_dof_index> &                     cm,
+      LinearAlgebra::distributed::Vector<float, MemorySpace::Default> &vec);
+
+    template void
+    set_zero_all(
+      const std::vector<types::global_dof_index> &                      cm,
+      LinearAlgebra::distributed::Vector<double, MemorySpace::Default> &vec);
+  } // namespace AffineConstraintsImplementation
+} // namespace internal
+
+template void
+AffineConstraints<float>::set_zero<
+  LinearAlgebra::distributed::Vector<float, MemorySpace::Default>>(
+  LinearAlgebra::distributed::Vector<float, MemorySpace::Default> &) const;
+
+template void
+AffineConstraints<double>::set_zero<
+  LinearAlgebra::distributed::Vector<double, MemorySpace::Default>>(
+  LinearAlgebra::distributed::Vector<double, MemorySpace::Default> &) const;
+
 DEAL_II_NAMESPACE_CLOSE
diff --git a/source/lac/affine_constraints_cuda.cc b/source/lac/affine_constraints_cuda.cc
deleted file mode 100644 (file)
index a64f7bb..0000000
+++ /dev/null
@@ -1,50 +0,0 @@
-// ---------------------------------------------------------------------
-//
-// Copyright (C) 2019 by the deal.II authors
-//
-// This file is part of the deal.II library.
-//
-// The deal.II library is free software; you can use it, redistribute
-// it, and/or modify it under the terms of the GNU Lesser General
-// Public License as published by the Free Software Foundation; either
-// version 2.1 of the License, or (at your option) any later version.
-// The full text of the license can be found in the file LICENSE.md at
-// the top level directory of deal.II.
-//
-// ---------------------------------------------------------------------
-
-#include <deal.II/lac/affine_constraints.templates.h>
-
-DEAL_II_NAMESPACE_OPEN
-
-#ifndef DOXYGEN
-
-namespace internal
-{
-  namespace AffineConstraintsImplementation
-  {
-    template void
-    set_zero_all(
-      const std::vector<types::global_dof_index> &                  cm,
-      LinearAlgebra::distributed::Vector<float, MemorySpace::CUDA> &vec);
-
-    template void
-    set_zero_all(
-      const std::vector<types::global_dof_index> &                   cm,
-      LinearAlgebra::distributed::Vector<double, MemorySpace::CUDA> &vec);
-  } // namespace AffineConstraintsImplementation
-} // namespace internal
-
-template void
-AffineConstraints<float>::set_zero<
-  LinearAlgebra::distributed::Vector<float, MemorySpace::CUDA>>(
-  LinearAlgebra::distributed::Vector<float, MemorySpace::CUDA> &) const;
-
-template void
-AffineConstraints<double>::set_zero<
-  LinearAlgebra::distributed::Vector<double, MemorySpace::CUDA>>(
-  LinearAlgebra::distributed::Vector<double, MemorySpace::CUDA> &) const;
-
-#endif // DOXYGEN
-
-DEAL_II_NAMESPACE_CLOSE

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.