From: Daniel Arndt Date: Sat, 31 Dec 2022 17:33:10 +0000 (+0100) Subject: Convert AffineConstraintsImplementation::set_zero_parallel X-Git-Tag: v9.5.0-rc1~620^2~6 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=4a6dd29ab5fe27942927424e5cc8b27e0cd7d144;p=dealii.git Convert AffineConstraintsImplementation::set_zero_parallel --- diff --git a/include/deal.II/lac/affine_constraints.templates.h b/include/deal.II/lac/affine_constraints.templates.h index 9b2c26cf73..c5d27717ff 100644 --- a/include/deal.II/lac/affine_constraints.templates.h +++ b/include/deal.II/lac/affine_constraints.templates.h @@ -2256,24 +2256,12 @@ namespace internal vec.zero_out_ghost_values(); } -#ifdef DEAL_II_WITH_CUDA - template - __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 void set_zero_parallel( - const std::vector & cm, - LinearAlgebra::distributed::Vector &vec, - size_type shift = 0) + const std::vector & cm, + LinearAlgebra::distributed::Vector &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<<>>( - 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 + constrained_local_dofs_device( + Kokkos::view_alloc(Kokkos::WithoutInitializing, + "constrained_local_dofs_device"), + n_constraints); + Kokkos::deep_copy(constrained_local_dofs_device, + Kokkos::View( + 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(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 void diff --git a/source/lac/affine_constraints.cc b/source/lac/affine_constraints.cc index 86ee7dd856..440bf4a2e9 100644 --- a/source/lac/affine_constraints.cc +++ b/source/lac/affine_constraints.cc @@ -150,4 +150,30 @@ dealii::AffineConstraints::distribute< # endif #endif +namespace internal +{ + namespace AffineConstraintsImplementation + { + template void + set_zero_all( + const std::vector & cm, + LinearAlgebra::distributed::Vector &vec); + + template void + set_zero_all( + const std::vector & cm, + LinearAlgebra::distributed::Vector &vec); + } // namespace AffineConstraintsImplementation +} // namespace internal + +template void +AffineConstraints::set_zero< + LinearAlgebra::distributed::Vector>( + LinearAlgebra::distributed::Vector &) const; + +template void +AffineConstraints::set_zero< + LinearAlgebra::distributed::Vector>( + LinearAlgebra::distributed::Vector &) 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 index a64f7bb119..0000000000 --- a/source/lac/affine_constraints_cuda.cc +++ /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_NAMESPACE_OPEN - -#ifndef DOXYGEN - -namespace internal -{ - namespace AffineConstraintsImplementation - { - template void - set_zero_all( - const std::vector & cm, - LinearAlgebra::distributed::Vector &vec); - - template void - set_zero_all( - const std::vector & cm, - LinearAlgebra::distributed::Vector &vec); - } // namespace AffineConstraintsImplementation -} // namespace internal - -template void -AffineConstraints::set_zero< - LinearAlgebra::distributed::Vector>( - LinearAlgebra::distributed::Vector &) const; - -template void -AffineConstraints::set_zero< - LinearAlgebra::distributed::Vector>( - LinearAlgebra::distributed::Vector &) const; - -#endif // DOXYGEN - -DEAL_II_NAMESPACE_CLOSE