]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Improve documentation of functors and replace functor by Functor 7519/head
authorBruno Turcksin <bruno.turcksin@gmail.com>
Fri, 14 Dec 2018 14:48:26 +0000 (14:48 +0000)
committerBruno Turcksin <bruno.turcksin@gmail.com>
Fri, 14 Dec 2018 14:50:00 +0000 (14:50 +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 d63a8624a28cbd11d00bcd449e59a3d1b2c9d643..8bb8b6c4efa7c027c0d1ee9e20f84d21c5af9222 100644 (file)
@@ -160,12 +160,21 @@ namespace CUDAWrappers
     __device__ void
     submit_gradient(const gradient_type &grad_in, const unsigned int q_point);
 
+    // clang-format off
     /**
-     * Apply the function @p func on every quadrature point.
+     * Apply the functor @p func on every quadrature point.
+     *
+     * @p func needs to define
+     * \code
+     * __device__ void operator()(
+     *   CUDAWrappers::FEEvaluation<dim, fe_degree, n_q_points_1d, n_components, Number> *fe_eval,
+     *   const unsigned int                                                               q_point) const;
+     * \endcode
      */
-    template <typename functor>
+    // clang-format on
+    template <typename Functor>
     __device__ void
-    apply_quad_point_operations(const functor &func);
+    apply_quad_point_operations(const Functor &func);
 
   private:
     types::global_dof_index *local_to_global;
@@ -427,10 +436,10 @@ namespace CUDAWrappers
             int n_q_points_1d,
             int n_components_,
             typename Number>
-  template <typename functor>
+  template <typename Functor>
   __device__ void
   FEEvaluation<dim, fe_degree, n_q_points_1d, n_components_, Number>::
-    apply_quad_point_operations(const functor &func)
+    apply_quad_point_operations(const Functor &func)
   {
     const unsigned int q_point =
       (dim == 1 ?
index 7bdec6767650e80cd4363bce1b6555b274892b64..67af1c1d7a3efa66d31378faa22e27340e838055 100644 (file)
@@ -208,13 +208,28 @@ namespace CUDAWrappers
     Data
     get_data(unsigned int color) const;
 
+    // clang-format off
     /**
      * This method runs the loop over all cells and apply the local operation on
      * each element in parallel. @p func is a functor which is applied on each color.
+     *
+     * @p func needs to define
+     * \code
+     * __device__ void operator()(
+     *   const unsigned int                                          cell,
+     *   const typename CUDAWrappers::MatrixFree<dim, Number>::Data *gpu_data,
+     *   CUDAWrappers::SharedData<dim, Number> *                     shared_data,
+     *   const Number *                                              src,
+     *   Number *                                                    dst) const;
+     *   static const unsigned int n_dofs_1d;
+     *   static const unsigned int n_local_dofs;
+     *   static const unsigned int n_q_points;
+     * \endcode
      */
-    template <typename functor, typename VectorType>
+    // clang-format on
+    template <typename Functor, typename VectorType>
     void
-    cell_loop(const functor &   func,
+    cell_loop(const Functor &   func,
               const VectorType &src,
               VectorType &      dst) const;
 
@@ -282,9 +297,9 @@ namespace CUDAWrappers
      * Helper function. Loop over all the cells and apply the functor on each
      * element in parallel. This function is used when MPI is not used.
      */
-    template <typename functor, typename VectorType>
+    template <typename Functor, typename VectorType>
     void
-    serial_cell_loop(const functor &   func,
+    serial_cell_loop(const Functor &   func,
                      const VectorType &src,
                      VectorType &      dst) const;
 
@@ -292,10 +307,10 @@ namespace CUDAWrappers
      * Helper function. Loop over all the cells and apply the functor on each
      * element in parallel. This function is used when MPI is used.
      */
-    template <typename functor>
+    template <typename Functor>
     void
     distributed_cell_loop(
-      const functor &                                                      func,
+      const Functor &                                                      func,
       const LinearAlgebra::distributed::Vector<Number, MemorySpace::CUDA> &src,
       LinearAlgebra::distributed::Vector<Number, MemorySpace::CUDA> &dst) const;
 
@@ -304,10 +319,10 @@ namespace CUDAWrappers
      * error. This function exists only because cell_loop needs
      * distributed_cell_loop() to exist for LinearAlgebra::CUDAWrappers::Vector.
      */
-    template <typename functor>
+    template <typename Functor>
     void
     distributed_cell_loop(
-      const functor &                                    func,
+      const Functor &                                    func,
       const LinearAlgebra::CUDAWrappers::Vector<Number> &src,
       LinearAlgebra::CUDAWrappers::Vector<Number> &      dst) const;
 
index ecfd74a09a9266b57a313c4a232564637fc74c89..001d5fa0ce9d38a6de53b403e2f10b98d91a88f7 100644 (file)
@@ -472,34 +472,34 @@ namespace CUDAWrappers
 
 
 
-    template <int dim, typename Number, typename functor>
+    template <int dim, typename Number, typename Functor>
     __global__ void
-    apply_kernel_shmem(functor                                      func,
+    apply_kernel_shmem(Functor                                      func,
                        const typename MatrixFree<dim, Number>::Data gpu_data,
                        const Number *                               src,
                        Number *                                     dst)
     {
       constexpr unsigned int cells_per_block =
-        cells_per_block_shmem(dim, functor::n_dofs_1d - 1);
+        cells_per_block_shmem(dim, Functor::n_dofs_1d - 1);
 
       constexpr unsigned int n_dofs_per_block =
-        cells_per_block * functor::n_local_dofs;
+        cells_per_block * Functor::n_local_dofs;
       constexpr unsigned int n_q_points_per_block =
-        cells_per_block * functor::n_q_points;
+        cells_per_block * Functor::n_q_points;
       // TODO make use of dynamically allocated shared memory
       __shared__ Number values[n_dofs_per_block];
       __shared__ Number gradients[dim][n_q_points_per_block];
 
-      const unsigned int local_cell = threadIdx.x / functor::n_dofs_1d;
+      const unsigned int local_cell = threadIdx.x / Functor::n_dofs_1d;
       const unsigned int cell =
         local_cell + cells_per_block * (blockIdx.x + gridDim.x * blockIdx.y);
 
       Number *gq[dim];
       for (int d = 0; d < dim; ++d)
-        gq[d] = &gradients[d][local_cell * functor::n_q_points];
+        gq[d] = &gradients[d][local_cell * Functor::n_q_points];
 
       SharedData<dim, Number> shared_data(
-        &values[local_cell * functor::n_local_dofs], gq);
+        &values[local_cell * Functor::n_local_dofs], gq);
 
       if (cell < gpu_data.n_cells)
         func(cell, &gpu_data, &shared_data, src, dst);
@@ -741,9 +741,9 @@ namespace CUDAWrappers
 
 
   template <int dim, typename Number>
-  template <typename functor, typename VectorType>
+  template <typename Functor, typename VectorType>
   void
-  MatrixFree<dim, Number>::cell_loop(const functor &   func,
+  MatrixFree<dim, Number>::cell_loop(const Functor &   func,
                                      const VectorType &src,
                                      VectorType &      dst) const
   {
@@ -971,15 +971,15 @@ namespace CUDAWrappers
 
 
   template <int dim, typename Number>
-  template <typename functor, typename VectorType>
+  template <typename Functor, typename VectorType>
   void
-  MatrixFree<dim, Number>::serial_cell_loop(const functor &   func,
+  MatrixFree<dim, Number>::serial_cell_loop(const Functor &   func,
                                             const VectorType &src,
                                             VectorType &      dst) const
   {
     // Execute the loop on the cells
     for (unsigned int i = 0; i < n_colors; ++i)
-      internal::apply_kernel_shmem<dim, Number, functor>
+      internal::apply_kernel_shmem<dim, Number, Functor>
         <<<grid_dim[i], block_dim[i]>>>(func,
                                         get_data(i),
                                         src.get_values(),
@@ -989,10 +989,10 @@ namespace CUDAWrappers
 
 
   template <int dim, typename Number>
-  template <typename functor>
+  template <typename Functor>
   void
   MatrixFree<dim, Number>::distributed_cell_loop(
-    const functor &                                                      func,
+    const Functor &                                                      func,
     const LinearAlgebra::distributed::Vector<Number, MemorySpace::CUDA> &src,
     LinearAlgebra::distributed::Vector<Number, MemorySpace::CUDA> &dst) const
   {
@@ -1005,7 +1005,7 @@ namespace CUDAWrappers
 
     // Execute the loop on the cells
     for (unsigned int i = 0; i < n_colors; ++i)
-      internal::apply_kernel_shmem<dim, Number, functor>
+      internal::apply_kernel_shmem<dim, Number, Functor>
         <<<grid_dim[i], block_dim[i]>>>(func,
                                         get_data(i),
                                         ghosted_src.get_values(),
@@ -1019,10 +1019,10 @@ namespace CUDAWrappers
 
 
   template <int dim, typename Number>
-  template <typename functor>
+  template <typename Functor>
   void
   MatrixFree<dim, Number>::distributed_cell_loop(
-    const functor &,
+    const Functor &,
     const LinearAlgebra::CUDAWrappers::Vector<Number> &,
     LinearAlgebra::CUDAWrappers::Vector<Number> &) const
   {

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.