]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Portable::MatrixFree support cell_loop() with BlockVector 18276/head
authorTimo Heister <timo.heister@gmail.com>
Fri, 21 Mar 2025 18:27:57 +0000 (14:27 -0400)
committerTimo Heister <timo.heister@gmail.com>
Tue, 25 Mar 2025 20:04:43 +0000 (16:04 -0400)
This is a first step for systems. Currently only sequential loops are supported.

include/deal.II/matrix_free/portable_matrix_free.h
include/deal.II/matrix_free/portable_matrix_free.templates.h
tests/matrix_free_kokkos/cell_loop_block_01.cc [new file with mode: 0644]
tests/matrix_free_kokkos/cell_loop_block_01.output [new file with mode: 0644]

index ec1a4f0fda79cc0d5b4fc2847c2a66ef7e56178e..713e6c707bf6624b27f16c23e5baa50dd70c1193 100644 (file)
@@ -18,6 +18,7 @@
 
 #include <deal.II/base/config.h>
 
+#include <deal.II/base/exceptions.h>
 #include <deal.II/base/memory_space.h>
 #include <deal.II/base/mpi_stub.h>
 #include <deal.II/base/partitioner.h>
@@ -37,6 +38,7 @@
 
 #include <deal.II/matrix_free/shape_info.h>
 
+#include <Kokkos_Array.hpp>
 #include <Kokkos_Core.hpp>
 
 
@@ -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<Number *, MemorySpace::Default::kokkos_space>;
 
+  /**
+   * 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 <typename Number>
+  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<Number> &src)
+      : n_blocks(1)
+      , blocks{src}
+    {}
+
+    /**
+     * Constructor from a LinearAlgebra::distributed::BlockVector. Creates
+     * a DeviceVector from each block and stores it.
+     */
+    template <typename MemorySpace>
+    DeviceBlockVector(
+      const LinearAlgebra::distributed::BlockVector<Number, MemorySpace> &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<Number>(src.block(b).get_values(),
+                                         src.block(b).locally_owned_size());
+    }
+
+    /**
+     * Access block @p index.
+     */
+    DEAL_II_HOST_DEVICE
+    DeviceVector<Number> &
+    block(unsigned int index)
+    {
+      AssertIndexRange(index, n_blocks);
+      return blocks[index];
+    }
+
+    /**
+     * Access block @p index.
+     */
+    DEAL_II_HOST_DEVICE const DeviceVector<Number>                           &
+    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<DeviceVector<Number>, 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<Number, MemorySpace::Default> &dst)
       const;
 
+    /**
+     * Same as above but for BlockVector.
+     */
+    template <typename Functor>
+    void
+    distributed_cell_loop(
+      const Functor                                                       &func,
+      const LinearAlgebra::distributed::BlockVector<Number,
+                                                    MemorySpace::Default> &src,
+      LinearAlgebra::distributed::BlockVector<Number, MemorySpace::Default>
+        &dst) const;
+
+
     /**
      * Unique ID associated with the object.
      */
index 1d4135d748ed44a34f10dc6ce4ecc2ca1f090468..2cb7acb9db9958844586de450833ca9e49234948 100644 (file)
@@ -18,6 +18,8 @@
 
 #include <deal.II/base/config.h>
 
+#include <deal.II/base/exception_macros.h>
+#include <deal.II/base/exceptions.h>
 #include <deal.II/base/graph_coloring.h>
 #include <deal.II/base/memory_space.h>
 
@@ -27,6 +29,8 @@
 #include <deal.II/fe/fe_values.h>
 #include <deal.II/fe/mapping_q1.h>
 
+#include <deal.II/lac/block_vector_base.h>
+
 #include <deal.II/matrix_free/portable_hanging_nodes_internal.h>
 #include <deal.II/matrix_free/portable_matrix_free.h>
 #include <deal.II/matrix_free/shape_info.h>
@@ -330,7 +334,7 @@ namespace Portable
 
 
 
-    template <int dim, typename Number, typename Functor>
+    template <int dim, typename Number, typename Functor, bool IsBlock>
     struct ApplyKernel
     {
       using TeamHandle = Kokkos::TeamPolicy<
@@ -359,14 +363,28 @@ namespace Portable
         LinearAlgebra::distributed::Vector<Number, MemorySpace::Default> &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<Number>(src.get_values(), src.locally_owned_size()))
+        , dst(DeviceVector<Number>(dst.get_values(), dst.locally_owned_size()))
+      {}
+
+      ApplyKernel(
+        Functor                                                 func,
+        const typename MatrixFree<dim, Number>::PrecomputedData gpu_data,
+        const LinearAlgebra::distributed::BlockVector<Number,
+                                                      MemorySpace::Default>
+          &src,
+        LinearAlgebra::distributed::BlockVector<Number, MemorySpace::Default>
+          &dst)
+        : func(func)
+        , gpu_data(gpu_data)
+        , src(src)
+        , dst(dst)
       {}
 
       Functor                                                 func;
       const typename MatrixFree<dim, Number>::PrecomputedData gpu_data;
-      const DeviceVector<Number>                              src;
-      DeviceVector<Number>                                    dst;
+      const DeviceBlockVector<Number>                         src;
+      DeviceBlockVector<Number>                               dst;
 
 
       // Provide the shared memory capacity. This function takes the team_size
@@ -405,8 +423,16 @@ namespace Portable
                                                     &gpu_data,
                                                     &shared_data};
 
-        DeviceVector<Number> nonconstdst = dst;
-        func(&data, src, nonconstdst);
+        if constexpr (IsBlock)
+          {
+            DeviceBlockVector<Number> nonconstdst = dst;
+            func(&data, src, nonconstdst);
+          }
+        else
+          {
+            DeviceVector<Number> 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<dim, Number, Functor> apply_kernel(
-            func, get_data(color), src, dst);
+
+          internal::
+            ApplyKernel<dim, Number, Functor, IsBlockVector<VectorType>::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 <int dim, typename Number>
+  template <typename Functor>
+  void
+  MatrixFree<dim, Number>::distributed_cell_loop(
+    const Functor &,
+    const LinearAlgebra::distributed::BlockVector<Number, MemorySpace::Default>
+      &,
+    LinearAlgebra::distributed::BlockVector<Number, MemorySpace::Default> &)
+    const
+  {
+    Assert(false, ExcNotImplemented());
+  }
+
+
+
   template <int dim, typename Number>
   template <typename Functor>
   void
@@ -1063,7 +1106,7 @@ namespace Portable
                     n_cells[0],
                     Kokkos::AUTO);
 
-                internal::ApplyKernel<dim, Number, Functor> apply_kernel(
+                internal::ApplyKernel<dim, Number, Functor, false> apply_kernel(
                   func, get_data(0), src, dst);
 
                 Kokkos::parallel_for(
@@ -1086,7 +1129,7 @@ namespace Portable
                     n_cells[1],
                     Kokkos::AUTO);
 
-                internal::ApplyKernel<dim, Number, Functor> apply_kernel(
+                internal::ApplyKernel<dim, Number, Functor, false> apply_kernel(
                   func, get_data(1), src, dst);
 
                 Kokkos::parallel_for(
@@ -1114,7 +1157,7 @@ namespace Portable
                     n_cells[2],
                     Kokkos::AUTO);
 
-                internal::ApplyKernel<dim, Number, Functor> apply_kernel(
+                internal::ApplyKernel<dim, Number, Functor, false> apply_kernel(
                   func, get_data(2), src, dst);
 
                 Kokkos::parallel_for(
@@ -1147,8 +1190,8 @@ namespace Portable
                       n_cells[i],
                       Kokkos::AUTO);
 
-                  internal::ApplyKernel<dim, Number, Functor> apply_kernel(
-                    func, get_data(i), src, dst);
+                  internal::ApplyKernel<dim, Number, Functor, false>
+                    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<dim, Number, Functor> apply_kernel(
+              internal::ApplyKernel<dim, Number, Functor, false> 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 (file)
index 0000000..ba1e413
--- /dev/null
@@ -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 <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/manifold_lib.h>
+
+#include "deal.II/lac/la_parallel_block_vector.h"
+
+#include <deal.II/matrix_free/portable_fe_evaluation.h>
+#include <deal.II/matrix_free/portable_matrix_free.h>
+
+#include "../tests.h"
+
+#include "Kokkos_Core.hpp"
+
+template <int dim,
+          int fe_degree,
+          int n_q_points_1d = fe_degree + 1,
+          typename Number   = double>
+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<dim, Number> &data_in)
+    : data(data_in){};
+
+  DEAL_II_HOST_DEVICE void
+  operator()(const typename Portable::MatrixFree<dim, Number>::Data *data,
+             const Portable::DeviceVector<Number>                   &src,
+             Portable::DeviceVector<Number>                         &dst) const
+  {
+    Portable::FEEvaluation<dim, fe_degree, n_q_points_1d, 1, Number> fe_eval(
+      data);
+
+    fe_eval.read_dof_values(src);
+    fe_eval.distribute_local_to_global(dst);
+  }
+
+  void
+  test() const
+  {
+    LinearAlgebra::distributed::Vector<Number, MemorySpace::Default> dst;
+    LinearAlgebra::distributed::Vector<Number, MemorySpace::Default> 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<dim, Number> &data;
+};
+
+
+template <int dim,
+          int fe_degree,
+          int n_q_points_1d = fe_degree + 1,
+          typename Number   = double>
+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<dim, Number> &data_in)
+    : data(data_in){};
+
+  DEAL_II_HOST_DEVICE void
+  operator()(const typename Portable::MatrixFree<dim, Number>::Data *data,
+             const Portable::DeviceBlockVector<Number>              &src,
+             Portable::DeviceBlockVector<Number>                    &dst) const
+  {
+    Portable::FEEvaluation<dim, fe_degree, n_q_points_1d, 1, Number> 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<Number, MemorySpace::Default> dst(
+      1);
+    LinearAlgebra::distributed::BlockVector<Number, MemorySpace::Default> 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<dim, Number> &data;
+};
+
+template <int dim, int fe_degree, int n_q_points_1d, typename Number>
+const unsigned int
+  MatrixFreeTest<dim, fe_degree, n_q_points_1d, Number>::n_local_dofs;
+
+template <int dim, int fe_degree, int n_q_points_1d, typename Number>
+const unsigned int
+  MatrixFreeTest<dim, fe_degree, n_q_points_1d, Number>::n_q_points;
+
+template <int dim, int fe_degree, int n_q_points_1d, typename Number>
+const unsigned int
+  MatrixFreeTestBlock<dim, fe_degree, n_q_points_1d, Number>::n_local_dofs;
+
+template <int dim, int fe_degree, int n_q_points_1d, typename Number>
+const unsigned int
+  MatrixFreeTestBlock<dim, fe_degree, n_q_points_1d, Number>::n_q_points;
+
+
+template <int dim, int fe_degree, typename number>
+void
+do_test(const DoFHandler<dim>           &dof,
+        const AffineConstraints<double> &constraints)
+{
+  Portable::MatrixFree<dim, number> mf_data;
+  {
+    const QGauss<1> quad(fe_degree + 1);
+    typename Portable::MatrixFree<dim, number>::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<dim, fe_degree, fe_degree + 1, number> mf(mf_data);
+  mf.test();
+  MatrixFreeTestBlock<dim, fe_degree, fe_degree + 1, number> mfb(mf_data);
+  mfb.test();
+}
+
+
+template <int dim, int fe_degree>
+void
+test()
+{
+  Triangulation<dim> 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<dim>       fe(fe_degree);
+  DoFHandler<dim> dof(tria);
+  dof.distribute_dofs(fe);
+
+  AffineConstraints<double> constraints;
+  DoFTools::make_hanging_node_constraints(dof, constraints);
+  constraints.close();
+
+  do_test<dim, fe_degree, double>(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 (file)
index 0000000..215a02e
--- /dev/null
@@ -0,0 +1,4 @@
+
+DEAL::Testing FE_Q<2>(1)
+DEAL::OK:4.00000
+DEAL::OK:4.00000

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.