]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Move constrained_dofs to Kokkos::View
authorBruno Turcksin <bruno.turcksin@gmail.com>
Tue, 21 Mar 2023 15:09:59 +0000 (11:09 -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 108274852b53fc1f172172fbcb2944b98b25c646..b3d1238624483b8ca0a6058e0583e57d9c6086ae 100644 (file)
@@ -536,9 +536,10 @@ namespace CUDAWrappers
       JxW;
 
     /**
-     * Pointer to the constrained degrees of freedom.
+     * Kokkos::View to the constrained degrees of freedom.
      */
-    types::global_dof_index *constrained_dofs;
+    Kokkos::View<types::global_dof_index *, MemorySpace::Default::kokkos_space>
+      constrained_dofs;
 
     /**
      * Mask deciding where constraints are set on a given cell.
index 7ebf79bd9135e36ce3ca8ead58f891ab66236148..5dd210b3f147598da23649aad85de268a64ebf1b 100644 (file)
@@ -607,7 +607,6 @@ namespace CUDAWrappers
   MatrixFree<dim, Number>::MatrixFree()
     : my_id(-1)
     , n_dofs(0)
-    , constrained_dofs(nullptr)
     , padding_length(0)
     , dof_handler(nullptr)
   {}
@@ -723,8 +722,6 @@ namespace CUDAWrappers
       Utilities::CUDA::free(constraint_mask_color_ptr);
     constraint_mask.clear();
 
-    Utilities::CUDA::free(constrained_dofs);
-
     internal::used_objects[my_id].store(false);
     my_id = -1;
   }
@@ -744,7 +741,7 @@ namespace CUDAWrappers
            ExcMessage("src and dst vectors have different size."));
     // FIXME When using C++17, we can use KOKKOS_CLASS_LAMBDA and this
     // work-around can be removed.
-    types::global_dof_index *constr_dofs = constrained_dofs;
+    auto               constr_dofs = constrained_dofs;
     const unsigned int size = internal::VectorLocalSize<VectorType>::get(dst);
     const Number *     src_ptr = src.get_values();
     Number *           dst_ptr = dst.get_values();
@@ -775,7 +772,7 @@ namespace CUDAWrappers
     Number *dst_ptr = dst.get_values();
     // FIXME When using C++17, we can use KOKKOS_CLASS_LAMBDA and this
     // work-around can be removed.
-    types::global_dof_index *constr_dofs = constrained_dofs;
+    auto constr_dofs = constrained_dofs;
     // When working with distributed vectors, the constrained dofs are
     // computed for ghosted vectors but we want to set the values of the
     // constrained dofs of non-ghosted vectors.
@@ -1144,17 +1141,17 @@ namespace CUDAWrappers
               }
           }
 
-        cuda_error = cudaMalloc(&constrained_dofs,
-                                n_constrained_dofs *
-                                  sizeof(dealii::types::global_dof_index));
-        AssertCuda(cuda_error);
+        constrained_dofs = Kokkos::View<types::global_dof_index *,
+                                        MemorySpace::Default::kokkos_space>(
+          Kokkos::view_alloc("constrained_dofs", Kokkos::WithoutInitializing),
+          n_constrained_dofs);
 
-        cuda_error = cudaMemcpy(constrained_dofs,
-                                constrained_dofs_host.data(),
-                                n_constrained_dofs *
-                                  sizeof(dealii::types::global_dof_index),
-                                cudaMemcpyHostToDevice);
-        AssertCuda(cuda_error);
+        Kokkos::View<types::global_dof_index *,
+                     MemorySpace::Default::kokkos_space,
+                     Kokkos::MemoryTraits<Kokkos::Unmanaged>>
+          constrained_dofs_host_view(constrained_dofs_host.data(),
+                                     constrained_dofs_host.size());
+        Kokkos::deep_copy(constrained_dofs, constrained_dofs_host_view);
       }
   }
 

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.