From 99709360c0b9cae519ad1970c96f2c2d7b1a4627 Mon Sep 17 00:00:00 2001 From: Timo Heister Date: Fri, 21 Mar 2025 14:27:57 -0400 Subject: [PATCH] Portable::MatrixFree support cell_loop() with BlockVector This is a first step for systems. Currently only sequential loops are supported. --- .../matrix_free/portable_matrix_free.h | 105 +++++++++ .../portable_matrix_free.templates.h | 73 +++++-- .../matrix_free_kokkos/cell_loop_block_01.cc | 204 ++++++++++++++++++ .../cell_loop_block_01.output | 4 + 4 files changed, 371 insertions(+), 15 deletions(-) create mode 100644 tests/matrix_free_kokkos/cell_loop_block_01.cc create mode 100644 tests/matrix_free_kokkos/cell_loop_block_01.output diff --git a/include/deal.II/matrix_free/portable_matrix_free.h b/include/deal.II/matrix_free/portable_matrix_free.h index ec1a4f0fda..713e6c707b 100644 --- a/include/deal.II/matrix_free/portable_matrix_free.h +++ b/include/deal.II/matrix_free/portable_matrix_free.h @@ -18,6 +18,7 @@ #include +#include #include #include #include @@ -37,6 +38,7 @@ #include +#include #include @@ -65,6 +67,14 @@ namespace Portable struct SharedData; #endif + /** + * Maximum number of DofHandler supported at the same time in + * Portable::MatrixFree computations. This limit also applies for the number + * of blocks in a BlockVector and the number of FEEvaluation objects that can + * be active in a single cell_loop(). + */ + static constexpr const unsigned int n_max_dof_handlers = 5; + /** * Type for source and destination vectors in device functions like * MatrixFree::cell_loop(). @@ -76,7 +86,89 @@ namespace Portable using DeviceVector = Kokkos::View; + /** + * A block vector used for source and destination vectors in device functions + * like MatrixFree::cell_loop(). + * + * The maximum number of block is limited by the constant @p n_max_dof_handlers. + */ + template + class DeviceBlockVector + { + public: + /** + * Constructor. + */ + DeviceBlockVector() + : n_blocks(0) + {} + /** + * Constructor. + */ + DeviceBlockVector(const DeviceBlockVector &other) = default; + + /** + * Constructor from a DeviceVector. Creates a DeviceBlockVector + * with a single block. + */ + explicit DeviceBlockVector(const DeviceVector &src) + : n_blocks(1) + , blocks{src} + {} + + /** + * Constructor from a LinearAlgebra::distributed::BlockVector. Creates + * a DeviceVector from each block and stores it. + */ + template + DeviceBlockVector( + const LinearAlgebra::distributed::BlockVector &src) + : n_blocks(src.n_blocks()) + { + Assert(src.n_blocks() <= n_max_dof_handlers, + ExcMessage("Portable::MatrixFree is configured with " + + Utilities::to_string(n_max_dof_handlers) + + " but you are passing a BlockVector with " + + Utilities::to_string(src.n_blocks()) + " blocks.")); + + for (int b = 0; b < src.n_blocks(); ++b) + blocks[b] = DeviceVector(src.block(b).get_values(), + src.block(b).locally_owned_size()); + } + + /** + * Access block @p index. + */ + DEAL_II_HOST_DEVICE + DeviceVector & + block(unsigned int index) + { + AssertIndexRange(index, n_blocks); + return blocks[index]; + } + + /** + * Access block @p index. + */ + DEAL_II_HOST_DEVICE const DeviceVector & + block(unsigned int index) const + { + AssertIndexRange(index, n_blocks); + return blocks[index]; + } + + private: + /** + * The number of components / blocks + */ + unsigned int n_blocks; + + /** + * Storage for the blocks + */ + Kokkos::Array, n_max_dof_handlers> blocks; + }; /** * This class collects all the data that is stored for the matrix free @@ -504,6 +596,19 @@ namespace Portable LinearAlgebra::distributed::Vector &dst) const; + /** + * Same as above but for BlockVector. + */ + template + void + distributed_cell_loop( + const Functor &func, + const LinearAlgebra::distributed::BlockVector &src, + LinearAlgebra::distributed::BlockVector + &dst) const; + + /** * Unique ID associated with the object. */ diff --git a/include/deal.II/matrix_free/portable_matrix_free.templates.h b/include/deal.II/matrix_free/portable_matrix_free.templates.h index 1d4135d748..2cb7acb9db 100644 --- a/include/deal.II/matrix_free/portable_matrix_free.templates.h +++ b/include/deal.II/matrix_free/portable_matrix_free.templates.h @@ -18,6 +18,8 @@ #include +#include +#include #include #include @@ -27,6 +29,8 @@ #include #include +#include + #include #include #include @@ -330,7 +334,7 @@ namespace Portable - template + template struct ApplyKernel { using TeamHandle = Kokkos::TeamPolicy< @@ -359,14 +363,28 @@ namespace Portable LinearAlgebra::distributed::Vector &dst) : func(func) , gpu_data(gpu_data) - , src(src.get_values(), src.locally_owned_size()) - , dst(dst.get_values(), dst.locally_owned_size()) + , src(DeviceVector(src.get_values(), src.locally_owned_size())) + , dst(DeviceVector(dst.get_values(), dst.locally_owned_size())) + {} + + ApplyKernel( + Functor func, + const typename MatrixFree::PrecomputedData gpu_data, + const LinearAlgebra::distributed::BlockVector + &src, + LinearAlgebra::distributed::BlockVector + &dst) + : func(func) + , gpu_data(gpu_data) + , src(src) + , dst(dst) {} Functor func; const typename MatrixFree::PrecomputedData gpu_data; - const DeviceVector src; - DeviceVector dst; + const DeviceBlockVector src; + DeviceBlockVector dst; // Provide the shared memory capacity. This function takes the team_size @@ -405,8 +423,16 @@ namespace Portable &gpu_data, &shared_data}; - DeviceVector nonconstdst = dst; - func(&data, src, nonconstdst); + if constexpr (IsBlock) + { + DeviceBlockVector nonconstdst = dst; + func(&data, src, nonconstdst); + } + else + { + DeviceVector nonconstdst = dst.block(0); + func(&data, src.block(0), nonconstdst); + } } }; } // namespace internal @@ -1018,8 +1044,10 @@ namespace Portable n_cells[color], Kokkos::AUTO); - internal::ApplyKernel apply_kernel( - func, get_data(color), src, dst); + + internal:: + ApplyKernel::value> + apply_kernel(func, get_data(color), src, dst); Kokkos::parallel_for("dealii::MatrixFree::serial_cell_loop", team_policy, @@ -1030,6 +1058,21 @@ namespace Portable + template + template + void + MatrixFree::distributed_cell_loop( + const Functor &, + const LinearAlgebra::distributed::BlockVector + &, + LinearAlgebra::distributed::BlockVector &) + const + { + Assert(false, ExcNotImplemented()); + } + + + template template void @@ -1063,7 +1106,7 @@ namespace Portable n_cells[0], Kokkos::AUTO); - internal::ApplyKernel apply_kernel( + internal::ApplyKernel apply_kernel( func, get_data(0), src, dst); Kokkos::parallel_for( @@ -1086,7 +1129,7 @@ namespace Portable n_cells[1], Kokkos::AUTO); - internal::ApplyKernel apply_kernel( + internal::ApplyKernel apply_kernel( func, get_data(1), src, dst); Kokkos::parallel_for( @@ -1114,7 +1157,7 @@ namespace Portable n_cells[2], Kokkos::AUTO); - internal::ApplyKernel apply_kernel( + internal::ApplyKernel apply_kernel( func, get_data(2), src, dst); Kokkos::parallel_for( @@ -1147,8 +1190,8 @@ namespace Portable n_cells[i], Kokkos::AUTO); - internal::ApplyKernel apply_kernel( - func, get_data(i), src, dst); + internal::ApplyKernel + apply_kernel(func, get_data(i), src, dst); Kokkos::parallel_for( "dealii::MatrixFree::distributed_cell_loop_" + @@ -1184,7 +1227,7 @@ namespace Portable n_cells[i], Kokkos::AUTO); - internal::ApplyKernel apply_kernel( + internal::ApplyKernel apply_kernel( func, get_data(i), ghosted_src, ghosted_dst); Kokkos::parallel_for( diff --git a/tests/matrix_free_kokkos/cell_loop_block_01.cc b/tests/matrix_free_kokkos/cell_loop_block_01.cc new file mode 100644 index 0000000000..ba1e413723 --- /dev/null +++ b/tests/matrix_free_kokkos/cell_loop_block_01.cc @@ -0,0 +1,204 @@ +// ------------------------------------------------------------------------ +// +// SPDX-License-Identifier: LGPL-2.1-or-later +// Copyright (C) 2019 - 2024 by the deal.II authors +// +// This file is part of the deal.II library. +// +// Part of the source code is dual licensed under Apache-2.0 WITH +// LLVM-exception OR LGPL-2.1-or-later. Detailed license information +// governing the source code and code contributions can be found in +// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II. +// +// ------------------------------------------------------------------------ + + +// check that Portable::FEEvaluation::submit_dof_value/get_dof_value +// works correctly. + +#include "deal.II/base/memory_space.h" + +#include +#include + +#include "deal.II/lac/la_parallel_block_vector.h" + +#include +#include + +#include "../tests.h" + +#include "Kokkos_Core.hpp" + +template +class MatrixFreeTest +{ +public: + static const unsigned int n_local_dofs = Utilities::pow(fe_degree + 1, dim); + static const unsigned int n_q_points = Utilities::pow(n_q_points_1d, dim); + + MatrixFreeTest(const Portable::MatrixFree &data_in) + : data(data_in){}; + + DEAL_II_HOST_DEVICE void + operator()(const typename Portable::MatrixFree::Data *data, + const Portable::DeviceVector &src, + Portable::DeviceVector &dst) const + { + Portable::FEEvaluation fe_eval( + data); + + fe_eval.read_dof_values(src); + fe_eval.distribute_local_to_global(dst); + } + + void + test() const + { + LinearAlgebra::distributed::Vector dst; + LinearAlgebra::distributed::Vector src; + + data.initialize_dof_vector(dst); + data.initialize_dof_vector(src); + src.add(1.0); + + data.cell_loop(*this, src, dst); + + Kokkos::fence(); + + deallog << "OK:" << dst.linfty_norm() << std::endl; + } + +protected: + const Portable::MatrixFree &data; +}; + + +template +class MatrixFreeTestBlock +{ +public: + static const unsigned int n_local_dofs = Utilities::pow(fe_degree + 1, dim); + static const unsigned int n_q_points = Utilities::pow(n_q_points_1d, dim); + + MatrixFreeTestBlock(const Portable::MatrixFree &data_in) + : data(data_in){}; + + DEAL_II_HOST_DEVICE void + operator()(const typename Portable::MatrixFree::Data *data, + const Portable::DeviceBlockVector &src, + Portable::DeviceBlockVector &dst) const + { + Portable::FEEvaluation fe_eval( + data); + + fe_eval.read_dof_values(src.block(0)); + fe_eval.distribute_local_to_global(dst.block(0)); + } + + void + test() const + { + LinearAlgebra::distributed::BlockVector dst( + 1); + LinearAlgebra::distributed::BlockVector src( + 1); + + data.initialize_dof_vector(dst.block(0)); + data.initialize_dof_vector(src.block(0)); + src.add(1.0); + + data.cell_loop(*this, src, dst); + + Kokkos::fence(); + + deallog << "OK:" << dst.block(0).linfty_norm() << std::endl; + } + +protected: + const Portable::MatrixFree &data; +}; + +template +const unsigned int + MatrixFreeTest::n_local_dofs; + +template +const unsigned int + MatrixFreeTest::n_q_points; + +template +const unsigned int + MatrixFreeTestBlock::n_local_dofs; + +template +const unsigned int + MatrixFreeTestBlock::n_q_points; + + +template +void +do_test(const DoFHandler &dof, + const AffineConstraints &constraints) +{ + Portable::MatrixFree mf_data; + { + const QGauss<1> quad(fe_degree + 1); + typename Portable::MatrixFree::AdditionalData data; + data.mapping_update_flags = update_values | update_gradients | + update_JxW_values | update_quadrature_points; + mf_data.reinit(dof, constraints, quad, data); + } + + deallog << "Testing " << dof.get_fe().get_name() << std::endl; + MatrixFreeTest mf(mf_data); + mf.test(); + MatrixFreeTestBlock mfb(mf_data); + mfb.test(); +} + + +template +void +test() +{ + Triangulation tria; + GridGenerator::hyper_cube(tria); + + // refine first and last cell + tria.begin(tria.n_levels() - 1)->set_refine_flag(); + tria.last()->set_refine_flag(); + tria.execute_coarsening_and_refinement(); + + FE_Q fe(fe_degree); + DoFHandler dof(tria); + dof.distribute_dofs(fe); + + AffineConstraints constraints; + DoFTools::make_hanging_node_constraints(dof, constraints); + constraints.close(); + + do_test(dof, constraints); +} + + + +int +main() +{ + initlog(); + + Kokkos::initialize(); + + test<2, 1>(); + + Kokkos::finalize(); + + return 0; +} diff --git a/tests/matrix_free_kokkos/cell_loop_block_01.output b/tests/matrix_free_kokkos/cell_loop_block_01.output new file mode 100644 index 0000000000..215a02e6a5 --- /dev/null +++ b/tests/matrix_free_kokkos/cell_loop_block_01.output @@ -0,0 +1,4 @@ + +DEAL::Testing FE_Q<2>(1) +DEAL::OK:4.00000 +DEAL::OK:4.00000 -- 2.39.5