]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Move local_to_global to Kokkos::View
authorBruno Turcksin <bruno.turcksin@gmail.com>
Tue, 21 Mar 2023 13:49:59 +0000 (09:49 -0400)
committerBruno Turcksin <bruno.turcksin@gmail.com>
Thu, 6 Apr 2023 13:11:26 +0000 (13:11 +0000)
include/deal.II/matrix_free/cuda_fe_evaluation.h
include/deal.II/matrix_free/cuda_matrix_free.h
include/deal.II/matrix_free/cuda_matrix_free.templates.h

index db1b3e06380d870a4152353ffab8005b2bb68065..aa1bd4e1fc33eca6462d5235242522bab137ce9d 100644 (file)
@@ -240,10 +240,19 @@ namespace CUDAWrappers
     apply_for_each_quad_point(const Functor &func);
 
   private:
-    types::global_dof_index *local_to_global;
-    unsigned int             n_cells;
-    unsigned int             padding_length;
-    const unsigned int       mf_object_id;
+    // FIXME We would like to use
+    // Kokkos::Subview<Kokkos::View<types::global_dof_index **,
+    // MemorySpace::Default::kokkos_space>, int, decltype(Kokkos::ALL)>
+    // but we get error: incomplete type is not allowed. I cannot reproduce
+    // outside of deal.II. Need to investigate more.
+    Kokkos::Subview<Kokkos::View<types::global_dof_index **,
+                                 MemorySpace::Default::kokkos_space>,
+                    int,
+                    Kokkos::pair<int, int>>
+                       local_to_global;
+    unsigned int       n_cells;
+    unsigned int       padding_length;
+    const unsigned int mf_object_id;
 
     const dealii::internal::MatrixFreeFunctions::ConstraintKinds
       constraint_mask;
@@ -290,7 +299,11 @@ namespace CUDAWrappers
     FEEvaluation(const unsigned int       cell_id,
                  const data_type *        data,
                  SharedData<dim, Number> *shdata)
-    : n_cells(data->n_cells)
+    : local_to_global(Kokkos::subview(
+        data->local_to_global,
+        cell_id,
+        Kokkos::pair<int, int>(0, Utilities::pow(n_q_points_1d, dim))))
+    , n_cells(data->n_cells)
     , padding_length(data->padding_length)
     , mf_object_id(data->id)
     , constraint_mask(data->constraint_mask[cell_id])
@@ -307,8 +320,6 @@ namespace CUDAWrappers
         Kokkos::pair<int, int>(0, Utilities::pow(n_q_points_1d, dim))))
     , values(shdata->values)
   {
-    local_to_global = data->local_to_global + padding_length * cell_id;
-
     for (unsigned int i = 0; i < dim; ++i)
       gradients[i] = shdata->gradients[i];
   }
index d6588813f6966176ca43e674c970b04c17cef16d..108274852b53fc1f172172fbcb2944b98b25c646 100644 (file)
@@ -183,7 +183,9 @@ namespace CUDAWrappers
        * Map the position in the local vector to the position in the global
        * vector.
        */
-      types::global_dof_index *local_to_global;
+      Kokkos::View<types::global_dof_index **,
+                   MemorySpace::Default::kokkos_space>
+        local_to_global;
 
       /**
        * Kokkos::View of the inverse Jacobian.
@@ -514,7 +516,9 @@ namespace CUDAWrappers
      * Map the position in the local vector to the position in the global
      * vector.
      */
-    std::vector<types::global_dof_index *> local_to_global;
+    std::vector<Kokkos::View<types::global_dof_index **,
+                             MemorySpace::Default::kokkos_space>>
+      local_to_global;
 
     /**
      * Vector of Kokkos::View of the inverse Jacobian associated to the cells of
@@ -737,7 +741,9 @@ namespace CUDAWrappers
      * Map the position in the local vector to the position in the global
      * vector.
      */
-    std::vector<types::global_dof_index> local_to_global;
+    typename Kokkos::View<types::global_dof_index **,
+                          MemorySpace::Default::kokkos_space>::HostMirror
+      local_to_global;
 
     /**
      * Kokkos::View of inverse Jacobians on the host.
@@ -815,9 +821,8 @@ namespace CUDAWrappers
         Kokkos::deep_copy(data_host.q_points, data.q_points);
       }
 
-    data_host.local_to_global.resize(n_elements);
-    Utilities::CUDA::copy_to_host(data.local_to_global,
-                                  data_host.local_to_global);
+    data_host.local_to_global = Kokkos::create_mirror(data.local_to_global);
+    Kokkos::deep_copy(data_host.local_to_global, data.local_to_global);
 
     if (update_flags & update_gradients)
       {
index d23c1b82bc10e42da251a9eb2283d9a51512e438..7ebf79bd9135e36ce3ca8ead58f891ab66236148 100644 (file)
@@ -232,8 +232,9 @@ namespace CUDAWrappers
 
     private:
       MatrixFree<dim, Number> *data;
-      // Host data
-      std::vector<types::global_dof_index> local_to_global_host;
+      Kokkos::View<types::global_dof_index **,
+                   MemorySpace::Default::kokkos_space>
+        local_to_global;
       Kokkos::View<Point<dim, Number> **, MemorySpace::Default::kokkos_space>
                                                                   q_points;
       Kokkos::View<Number **, MemorySpace::Default::kokkos_space> JxW;
@@ -348,7 +349,13 @@ namespace CUDAWrappers
         data->block_dim[color] =
           dim3(n_dofs_1d * cells_per_block, n_dofs_1d, n_dofs_1d);
 
-      local_to_global_host.resize(n_cells * padding_length);
+
+      local_to_global = Kokkos::View<types::global_dof_index **,
+                                     MemorySpace::Default::kokkos_space>(
+        Kokkos::view_alloc("local_to_global_" + std::to_string(color),
+                           Kokkos::WithoutInitializing),
+        n_cells,
+        dofs_per_cell);
 
       if (update_flags & update_quadrature_points)
         q_points = Kokkos::View<Point<dim, Number> **,
@@ -406,9 +413,13 @@ namespace CUDAWrappers
                                       lexicographic_dof_indices,
                                       cell_id_view);
 
-      memcpy(&local_to_global_host[cell_id * padding_length],
-             lexicographic_dof_indices.data(),
-             dofs_per_cell * sizeof(types::global_dof_index));
+      // FIXME too many deep_copy
+      auto local_to_global_host =
+        Kokkos::create_mirror_view_and_copy(MemorySpace::Host::kokkos_space{},
+                                            local_to_global);
+      for (unsigned int i = 0; i < dofs_per_cell; ++i)
+        local_to_global_host(cell_id, i) = lexicographic_dof_indices[i];
+      Kokkos::deep_copy(local_to_global, local_to_global_host);
 
       fe_values.reinit(cell);
 
@@ -460,11 +471,7 @@ namespace CUDAWrappers
       const unsigned int n_cells = data->n_cells[color];
 
       // Local-to-global mapping
-      alloc_and_copy(
-        &data->local_to_global[color],
-        ArrayView<const types::global_dof_index>(local_to_global_host.data(),
-                                                 local_to_global_host.size()),
-        n_cells * padding_length);
+      data->local_to_global[color] = local_to_global;
 
       // Quadrature points
       if (update_flags & update_quadrature_points)
@@ -712,10 +719,6 @@ namespace CUDAWrappers
   void
   MatrixFree<dim, Number>::free()
   {
-    for (auto &local_to_global_color_ptr : local_to_global)
-      Utilities::CUDA::free(local_to_global_color_ptr);
-    local_to_global.clear();
-
     for (auto &constraint_mask_color_ptr : constraint_mask)
       Utilities::CUDA::free(constraint_mask_color_ptr);
     constraint_mask.clear();

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.