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

index b3d1238624483b8ca0a6058e0583e57d9c6086ae..80c5218aecfa1f2033d4de14a8611be2ccc8f6df 100644 (file)
@@ -221,7 +221,9 @@ namespace CUDAWrappers
       /**
        * Mask deciding where constraints are set on a given cell.
        */
-      dealii::internal::MatrixFreeFunctions::ConstraintKinds *constraint_mask;
+      Kokkos::View<dealii::internal::MatrixFreeFunctions::ConstraintKinds *,
+                   MemorySpace::Default::kokkos_space>
+        constraint_mask;
 
       /**
        * If true, use graph coloring has been used and we can simply add into
@@ -544,7 +546,9 @@ namespace CUDAWrappers
     /**
      * Mask deciding where constraints are set on a given cell.
      */
-    std::vector<dealii::internal::MatrixFreeFunctions::ConstraintKinds *>
+    std::vector<
+      Kokkos::View<dealii::internal::MatrixFreeFunctions::ConstraintKinds *,
+                   MemorySpace::Default::kokkos_space>>
       constraint_mask;
 
     /**
@@ -782,8 +786,9 @@ namespace CUDAWrappers
     /**
      * Mask deciding where constraints are set on a given cell.
      */
-    std::vector<dealii::internal::MatrixFreeFunctions::ConstraintKinds>
-      constraint_mask;
+    typename Kokkos::View<
+      dealii::internal::MatrixFreeFunctions::ConstraintKinds *,
+      MemorySpace::Default::kokkos_space>::HostMirror constraint_mask;
 
     /**
      * If true, use graph coloring has been used and we can simply add into
@@ -837,9 +842,8 @@ namespace CUDAWrappers
         Kokkos::deep_copy(data_host.JxW, data.JxW);
       }
 
-    data_host.constraint_mask.resize(data_host.n_cells);
-    Utilities::CUDA::copy_to_host(data.constraint_mask,
-                                  data_host.constraint_mask);
+    data_host.constraint_mask = Kokkos::create_mirror(data.constraint_mask);
+    Kokkos::deep_copy(data_host.constraint_mask, data.constraint_mask);
 
     return data_host;
   }
index 5dd210b3f147598da23649aad85de268a64ebf1b..f39594ed3963de751b86455bbd940871a854f417 100644 (file)
@@ -240,8 +240,9 @@ namespace CUDAWrappers
       Kokkos::View<Number **, MemorySpace::Default::kokkos_space> JxW;
       Kokkos::View<Number **[dim][dim], MemorySpace::Default::kokkos_space>
         inv_jacobian;
-      std::vector<dealii::internal::MatrixFreeFunctions::ConstraintKinds>
-        constraint_mask_host;
+      Kokkos::View<dealii::internal::MatrixFreeFunctions::ConstraintKinds *,
+                   MemorySpace::Default::kokkos_space>
+        constraint_mask;
       // Local buffer
       std::vector<types::global_dof_index> local_dof_indices;
       FEValues<dim>                        fe_values;
@@ -380,7 +381,11 @@ namespace CUDAWrappers
             n_cells,
             dofs_per_cell);
 
-      constraint_mask_host.resize(n_cells);
+      // Initialize to zero, i.e., unconstrained cell
+      constraint_mask =
+        Kokkos::View<dealii::internal::MatrixFreeFunctions::ConstraintKinds *,
+                     MemorySpace::Default::kokkos_space>(
+          "constraint_mask_" + std::to_string(color), n_cells);
     }
 
 
@@ -404,6 +409,10 @@ namespace CUDAWrappers
       for (unsigned int i = 0; i < dofs_per_cell; ++i)
         lexicographic_dof_indices[i] = local_dof_indices[lexicographic_inv[i]];
 
+      // FIXME too many deep_copy
+      auto constraint_mask_host =
+        Kokkos::create_mirror_view_and_copy(MemorySpace::Host::kokkos_space{},
+                                            constraint_mask);
       const ArrayView<dealii::internal::MatrixFreeFunctions::ConstraintKinds>
         cell_id_view(constraint_mask_host[cell_id]);
 
@@ -412,6 +421,7 @@ namespace CUDAWrappers
                                       {lexicographic_inv},
                                       lexicographic_dof_indices,
                                       cell_id_view);
+      Kokkos::deep_copy(constraint_mask, constraint_mask_host);
 
       // FIXME too many deep_copy
       auto local_to_global_host =
@@ -491,11 +501,7 @@ namespace CUDAWrappers
           data->inv_jacobian[color] = inv_jacobian;
         }
 
-      alloc_and_copy(
-        &data->constraint_mask[color],
-        ArrayView<const dealii::internal::MatrixFreeFunctions::ConstraintKinds>(
-          constraint_mask_host.data(), constraint_mask_host.size()),
-        n_cells);
+      data->constraint_mask[color] = constraint_mask;
     }
 
 
@@ -718,10 +724,6 @@ namespace CUDAWrappers
   void
   MatrixFree<dim, Number>::free()
   {
-    for (auto &constraint_mask_color_ptr : constraint_mask)
-      Utilities::CUDA::free(constraint_mask_color_ptr);
-    constraint_mask.clear();
-
     internal::used_objects[my_id].store(false);
     my_id = -1;
   }

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.