Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
cuda_matrix_free.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2016 - 2020 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #ifndef dealii_cuda_matrix_free_h
18 #define dealii_cuda_matrix_free_h
19 
20 #include <deal.II/base/config.h>
21 
22 #ifdef DEAL_II_COMPILER_CUDA_AWARE
23 
24 # include <deal.II/base/cuda_size.h>
25 # include <deal.II/base/mpi.h>
26 # include <deal.II/base/quadrature.h>
27 # include <deal.II/base/tensor.h>
28 
30 
32 # include <deal.II/fe/mapping.h>
33 # include <deal.II/fe/mapping_q1.h>
34 
36 # include <deal.II/lac/cuda_vector.h>
38 
39 
41 
42 namespace CUDAWrappers
43 {
44  // forward declaration
45 # ifndef DOXYGEN
46  namespace internal
47  {
48  template <int dim, typename Number>
49  class ReinitHelper;
50  }
51 # endif
52 
79  template <int dim, typename Number = double>
80  class MatrixFree : public Subscriptor
81  {
82  public:
85 
92  {
95  };
96 
101  {
110  const bool use_coloring = false,
111  const bool overlap_communication_computation = false)
116  {
117 # ifndef DEAL_II_MPI_WITH_CUDA_SUPPORT
118  AssertThrow(
120  ExcMessage(
121  "Overlapping communication and computation requires CUDA-aware MPI."));
122 # endif
124  AssertThrow(
125  use_coloring == false || overlap_communication_computation == false,
126  ExcMessage(
127  "Overlapping communication and coloring are incompatible options. Only one of them can be enabled."));
128  }
144 
151 
157  };
158 
163  struct Data
164  {
169 
175 
179  Number *inv_jacobian;
180 
184  Number *JxW;
185 
189  unsigned int n_cells;
190 
194  unsigned int padding_length;
195 
199  unsigned int row_start;
200 
204  unsigned int *constraint_mask;
205 
211  };
212 
216  MatrixFree();
217 
221  unsigned int
222  get_padding_length() const;
223 
232  void
233  reinit(const Mapping<dim> & mapping,
234  const DoFHandler<dim> & dof_handler,
235  const AffineConstraints<Number> &constraints,
236  const Quadrature<1> & quad,
237  const AdditionalData & additional_data = AdditionalData());
238 
242  void
243  reinit(const DoFHandler<dim> & dof_handler,
244  const AffineConstraints<Number> &constraints,
245  const Quadrature<1> & quad,
247 
251  Data
252  get_data(unsigned int color) const;
253 
254  // clang-format off
272  // clang-format on
273  template <typename Functor, typename VectorType>
274  void
275  cell_loop(const Functor & func,
276  const VectorType &src,
277  VectorType & dst) const;
278 
294  template <typename Functor>
295  void
296  evaluate_coefficients(Functor func) const;
297 
302  template <typename VectorType>
303  void
304  copy_constrained_values(const VectorType &src, VectorType &dst) const;
305 
311  template <typename VectorType>
312  void
313  set_constrained_values(const Number val, VectorType &dst) const;
314 
319  void
322 
328  void
331 
335  void
336  free();
337 
341  std::size_t
342  memory_consumption() const;
343 
344  private:
348  void
349  internal_reinit(const Mapping<dim> & mapping,
350  const DoFHandler<dim> & dof_handler,
351  const AffineConstraints<Number> &constraints,
352  const Quadrature<1> & quad,
353  std::shared_ptr<const MPI_Comm> comm,
354  const AdditionalData additional_data);
355 
360  template <typename Functor, typename VectorType>
361  void
362  serial_cell_loop(const Functor & func,
363  const VectorType &src,
364  VectorType & dst) const;
365 
370  template <typename Functor>
371  void
373  const Functor & func,
376 
382  template <typename Functor>
383  void
385  const Functor & func,
388 
393  template <typename VectorType>
394  void
396  VectorType & dst) const;
397 
402  void
406 
413  void
417 
422  template <typename VectorType>
423  void
424  serial_set_constrained_values(const Number val, VectorType &dst) const;
425 
430  void
432  const Number val,
434 
441  void
443  const Number val,
445 
451 
458 
464 
469 
473  unsigned int fe_degree;
474 
478  unsigned int dofs_per_cell;
479 
483  unsigned int n_constrained_dofs;
484 
488  unsigned int q_points_per_cell;
489 
493  unsigned int n_colors;
494 
498  std::vector<unsigned int> n_cells;
499 
504  std::vector<point_type *> q_points;
505 
510  std::vector<types::global_dof_index *> local_to_global;
511 
516  std::vector<Number *> inv_jacobian;
517 
522  std::vector<Number *> JxW;
523 
528 
532  std::vector<unsigned int *> constraint_mask;
533 
538  std::vector<dim3> grid_dim;
539 
544  std::vector<dim3> block_dim;
545 
550  std::shared_ptr<const Utilities::MPI::Partitioner> partitioner;
551 
555  unsigned int cells_per_block;
556 
562 
568 
573  unsigned int padding_length;
574 
578  std::vector<unsigned int> row_start;
579 
580  friend class internal::ReinitHelper<dim, Number>;
581  };
582 
583 
584 
585  // TODO find a better place to put these things
586 
590  template <int dim, typename Number>
591  struct SharedData
592  {
596  __device__
597  SharedData(Number *vd, Number *gq[dim])
598  : values(vd)
599  {
600  for (int d = 0; d < dim; ++d)
601  gradients[d] = gq[d];
602  }
603 
607  Number *values;
608 
614  Number *gradients[dim];
615  };
616 
617 
618 
619  // This function determines the number of cells per block, possibly at compile
620  // time (by virtue of being 'constexpr')
621  // TODO this function should be rewritten using meta-programming
622  __host__ __device__ constexpr unsigned int
623  cells_per_block_shmem(int dim, int fe_degree)
624  {
625  /* clang-format off */
626  // We are limiting the number of threads according to the
627  // following formulas:
628  // - in 2D: `threads = cells * (k+1)^d <= 4*CUDAWrappers::warp_size`
629  // - in 3D: `threads = cells * (k+1)^d <= 2*CUDAWrappers::warp_size`
630  return dim==2 ? (fe_degree==1 ? CUDAWrappers::warp_size : // 128
631  fe_degree==2 ? CUDAWrappers::warp_size/4 : // 72
632  fe_degree==3 ? CUDAWrappers::warp_size/8 : // 64
633  fe_degree==4 ? CUDAWrappers::warp_size/8 : // 100
634  1) :
635  dim==3 ? (fe_degree==1 ? CUDAWrappers::warp_size/4 : // 64
636  fe_degree==2 ? CUDAWrappers::warp_size/16 : // 54
637  1) : 1;
638  /* clang-format on */
639  }
640 
641 
642 
648  template <int dim>
649  __device__ inline unsigned int
650  q_point_id_in_cell(const unsigned int n_q_points_1d)
651  {
652  return (dim == 1 ?
653  threadIdx.x % n_q_points_1d :
654  dim == 2 ?
655  threadIdx.x % n_q_points_1d + n_q_points_1d * threadIdx.y :
656  threadIdx.x % n_q_points_1d +
657  n_q_points_1d * (threadIdx.y + n_q_points_1d * threadIdx.z));
658  }
659 
660 
661 
668  template <int dim, typename Number>
669  __device__ inline unsigned int
671  const unsigned int cell,
672  const typename CUDAWrappers::MatrixFree<dim, Number>::Data *data,
673  const unsigned int n_q_points_1d,
674  const unsigned int n_q_points)
675  {
676  return (data->row_start / data->padding_length + cell) * n_q_points +
677  q_point_id_in_cell<dim>(n_q_points_1d);
678  }
679 
680 
681 
687  template <int dim, typename Number>
688  __device__ inline typename CUDAWrappers::MatrixFree<dim, Number>::point_type &
690  const unsigned int cell,
691  const typename CUDAWrappers::MatrixFree<dim, Number>::Data *data,
692  const unsigned int n_q_points_1d)
693  {
694  return *(data->q_points + data->padding_length * cell +
695  q_point_id_in_cell<dim>(n_q_points_1d));
696  }
697 } // namespace CUDAWrappers
698 
700 
701 #endif
702 
703 #endif
CUDAWrappers::MatrixFree::Data::row_start
unsigned int row_start
Definition: cuda_matrix_free.h:199
CUDAWrappers::MatrixFree
Definition: cuda_matrix_free.h:80
CUDAWrappers::MatrixFree::AdditionalData::mapping_update_flags
UpdateFlags mapping_update_flags
Definition: cuda_matrix_free.h:143
CUDAWrappers::MatrixFree::Data::local_to_global
types::global_dof_index * local_to_global
Definition: cuda_matrix_free.h:174
CUDAWrappers::MatrixFree::Data::padding_length
unsigned int padding_length
Definition: cuda_matrix_free.h:194
CUDAWrappers::MatrixFree::q_points
std::vector< point_type * > q_points
Definition: cuda_matrix_free.h:504
update_quadrature_points
@ update_quadrature_points
Transformed quadrature points.
Definition: fe_update_flags.h:122
fe_update_flags.h
CUDAWrappers::MatrixFree::AdditionalData::use_coloring
bool use_coloring
Definition: cuda_matrix_free.h:150
CUDAWrappers::MatrixFree::grid_dim
std::vector< dim3 > grid_dim
Definition: cuda_matrix_free.h:538
CUDAWrappers::MatrixFree::serial_copy_constrained_values
void serial_copy_constrained_values(const VectorType &src, VectorType &dst) const
LinearAlgebra::distributed::Vector
Definition: la_parallel_vector.h:226
CUDAWrappers::MatrixFree::parallel_over_elem
@ parallel_over_elem
Definition: cuda_matrix_free.h:94
cuda_vector.h
CUDAWrappers::MatrixFree::internal_reinit
void internal_reinit(const Mapping< dim > &mapping, const DoFHandler< dim > &dof_handler, const AffineConstraints< Number > &constraints, const Quadrature< 1 > &quad, std::shared_ptr< const MPI_Comm > comm, const AdditionalData additional_data)
CUDAWrappers::MatrixFree::constrained_dofs
types::global_dof_index * constrained_dofs
Definition: cuda_matrix_free.h:527
CUDAWrappers::MatrixFree::cell_loop
void cell_loop(const Functor &func, const VectorType &src, VectorType &dst) const
CUDAWrappers::SharedData::SharedData
SharedData(Number *vd, Number *gq[dim])
Definition: cuda_matrix_free.h:597
cuda_size.h
CUDAWrappers::SharedData::gradients
Number * gradients[dim]
Definition: cuda_matrix_free.h:614
mapping_q1.h
CUDAWrappers::MatrixFree::AdditionalData::parallelization_scheme
ParallelizationScheme parallelization_scheme
Definition: cuda_matrix_free.h:133
CUDAWrappers::MatrixFree::use_coloring
bool use_coloring
Definition: cuda_matrix_free.h:457
VectorType
CUDAWrappers::SharedData
Definition: cuda_matrix_free.h:591
CUDAWrappers::MatrixFree::copy_constrained_values
void copy_constrained_values(const VectorType &src, VectorType &dst) const
CUDAWrappers::MatrixFree::serial_cell_loop
void serial_cell_loop(const Functor &func, const VectorType &src, VectorType &dst) const
CUDAWrappers::MatrixFree::q_point_id_in_cell
unsigned int q_point_id_in_cell(const unsigned int n_q_points_1d)
Definition: cuda_matrix_free.h:650
CUDAWrappers::warp_size
constexpr int warp_size
Definition: cuda_size.h:40
CUDAWrappers::MatrixFree::Data::q_points
point_type * q_points
Definition: cuda_matrix_free.h:168
CUDAWrappers::SharedData::values
Number * values
Definition: cuda_matrix_free.h:607
Physics::Elasticity::Kinematics::d
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
DoFHandler< dim >
la_parallel_vector.h
CUDAWrappers::MatrixFree::n_cells
std::vector< unsigned int > n_cells
Definition: cuda_matrix_free.h:498
CUDAWrappers::MatrixFree::cells_per_block
unsigned int cells_per_block
Definition: cuda_matrix_free.h:555
CUDAWrappers::MatrixFree::parallel_in_elem
@ parallel_in_elem
Definition: cuda_matrix_free.h:93
CUDAWrappers::MatrixFree::n_dofs
types::global_dof_index n_dofs
Definition: cuda_matrix_free.h:468
Subscriptor
Definition: subscriptor.h:62
LinearAlgebra::CUDAWrappers::Vector
Definition: cuda_vector.h:56
CUDAWrappers::MatrixFree::get_padding_length
unsigned int get_padding_length() const
CUDAWrappers::MatrixFree::n_constrained_dofs
unsigned int n_constrained_dofs
Definition: cuda_matrix_free.h:483
CUDAWrappers::MatrixFree::Data::JxW
Number * JxW
Definition: cuda_matrix_free.h:184
CUDAWrappers::MatrixFree::reinit
void reinit(const Mapping< dim > &mapping, const DoFHandler< dim > &dof_handler, const AffineConstraints< Number > &constraints, const Quadrature< 1 > &quad, const AdditionalData &additional_data=AdditionalData())
Mapping< dim >
tensor.h
CUDAWrappers::MatrixFree::Data
Definition: cuda_matrix_free.h:163
CUDAWrappers::MatrixFree::get_quadrature_point
CUDAWrappers::MatrixFree< dim, Number >::point_type & get_quadrature_point(const unsigned int cell, const typename CUDAWrappers::MatrixFree< dim, Number >::Data *data, const unsigned int n_q_points_1d)
Definition: cuda_matrix_free.h:689
CUDAWrappers::MatrixFree::set_constrained_values
void set_constrained_values(const Number val, VectorType &dst) const
CUDAWrappers::MatrixFree::AdditionalData::overlap_communication_computation
bool overlap_communication_computation
Definition: cuda_matrix_free.h:156
CUDAWrappers::MatrixFree::row_start
std::vector< unsigned int > row_start
Definition: cuda_matrix_free.h:578
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
CUDAWrappers
Definition: cuda_size.h:23
CUDAWrappers::MatrixFree::distributed_copy_constrained_values
void distributed_copy_constrained_values(const LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &src, LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &dst) const
mpi.h
CUDAWrappers::MatrixFree::Data::inv_jacobian
Number * inv_jacobian
Definition: cuda_matrix_free.h:179
Tensor
Definition: tensor.h:450
CUDAWrappers::MatrixFree::serial_set_constrained_values
void serial_set_constrained_values(const Number val, VectorType &dst) const
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
CUDAWrappers::MatrixFree::local_to_global
std::vector< types::global_dof_index * > local_to_global
Definition: cuda_matrix_free.h:510
CUDAWrappers::MatrixFree::memory_consumption
std::size_t memory_consumption() const
update_gradients
@ update_gradients
Shape function gradients.
Definition: fe_update_flags.h:84
CUDAWrappers::MatrixFree::q_points_per_cell
unsigned int q_points_per_cell
Definition: cuda_matrix_free.h:488
CUDAWrappers::MatrixFree::ParallelizationScheme
ParallelizationScheme
Definition: cuda_matrix_free.h:91
CUDAWrappers::MatrixFree::Data::use_coloring
bool use_coloring
Definition: cuda_matrix_free.h:210
CUDAWrappers::MatrixFree::evaluate_coefficients
void evaluate_coefficients(Functor func) const
CUDAWrappers::MatrixFree::free
void free()
CUDAWrappers::MatrixFree::Data::constraint_mask
unsigned int * constraint_mask
Definition: cuda_matrix_free.h:204
CUDAWrappers::MatrixFree::partitioner
std::shared_ptr< const Utilities::MPI::Partitioner > partitioner
Definition: cuda_matrix_free.h:550
CUDAWrappers::MatrixFree::local_q_point_id
unsigned int local_q_point_id(const unsigned int cell, const typename CUDAWrappers::MatrixFree< dim, Number >::Data *data, const unsigned int n_q_points_1d, const unsigned int n_q_points)
Definition: cuda_matrix_free.h:670
CUDAWrappers::MatrixFree::dofs_per_cell
unsigned int dofs_per_cell
Definition: cuda_matrix_free.h:478
CUDAWrappers::MatrixFree::distributed_set_constrained_values
void distributed_set_constrained_values(const Number val, LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &dst) const
CUDAWrappers::MatrixFree::distributed_cell_loop
void distributed_cell_loop(const Functor &func, const LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &src, LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &dst) const
CUDAWrappers::MatrixFree::block_dim
std::vector< dim3 > block_dim
Definition: cuda_matrix_free.h:544
CUDAWrappers::MatrixFree::constraint_mask
std::vector< unsigned int * > constraint_mask
Definition: cuda_matrix_free.h:532
UpdateFlags
UpdateFlags
Definition: fe_update_flags.h:66
CUDAWrappers::MatrixFree::parallelization_scheme
ParallelizationScheme parallelization_scheme
Definition: cuda_matrix_free.h:450
CUDAWrappers::MatrixFree::overlap_communication_computation
bool overlap_communication_computation
Definition: cuda_matrix_free.h:463
unsigned int
CUDAWrappers::MatrixFree::Data::n_cells
unsigned int n_cells
Definition: cuda_matrix_free.h:189
AffineConstraints
Definition: affine_constraints.h:180
update_JxW_values
@ update_JxW_values
Transformed quadrature weights.
Definition: fe_update_flags.h:129
CUDAWrappers::MatrixFree::fe_degree
unsigned int fe_degree
Definition: cuda_matrix_free.h:473
mapping.h
CUDAWrappers::MatrixFree::constraint_block_dim
dim3 constraint_block_dim
Definition: cuda_matrix_free.h:567
dof_handler.h
affine_constraints.h
CUDAWrappers::MatrixFree::MatrixFree
MatrixFree()
CUDAWrappers::MatrixFree::get_data
Data get_data(unsigned int color) const
Point
Definition: point.h:111
quadrature.h
CUDAWrappers::MatrixFree::constraint_grid_dim
dim3 constraint_grid_dim
Definition: cuda_matrix_free.h:561
config.h
CUDAWrappers::MatrixFree::AdditionalData
Definition: cuda_matrix_free.h:100
internal
Definition: aligned_vector.h:369
CUDAWrappers::MatrixFree::inv_jacobian
std::vector< Number * > inv_jacobian
Definition: cuda_matrix_free.h:516
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
CUDAWrappers::MatrixFree::initialize_dof_vector
void initialize_dof_vector(LinearAlgebra::CUDAWrappers::Vector< Number > &vec) const
Quadrature< 1 >
CUDAWrappers::MatrixFree::n_colors
unsigned int n_colors
Definition: cuda_matrix_free.h:493
AssertThrow
#define AssertThrow(cond, exc)
Definition: exceptions.h:1531
CUDAWrappers::cells_per_block_shmem
constexpr __host__ unsigned int cells_per_block_shmem(int dim, int fe_degree)
Definition: cuda_matrix_free.h:623
CUDAWrappers::MatrixFree::AdditionalData::AdditionalData
AdditionalData(const ParallelizationScheme parallelization_scheme=parallel_in_elem, const UpdateFlags mapping_update_flags=update_gradients|update_JxW_values|update_quadrature_points, const bool use_coloring=false, const bool overlap_communication_computation=false)
Definition: cuda_matrix_free.h:105
CUDAWrappers::MatrixFree::padding_length
unsigned int padding_length
Definition: cuda_matrix_free.h:573
CUDAWrappers::MatrixFree::JxW
std::vector< Number * > JxW
Definition: cuda_matrix_free.h:522