|
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\}}\)
Go to the documentation of this file.
17 #ifndef dealii_cuda_matrix_free_h
18 #define dealii_cuda_matrix_free_h
22 #ifdef DEAL_II_COMPILER_CUDA_AWARE
48 template <
int dim,
typename Number>
79 template <
int dim,
typename Number =
double>
117 # ifndef DEAL_II_MPI_WITH_CUDA_SUPPORT
121 "Overlapping communication and computation requires CUDA-aware MPI."));
127 "Overlapping communication and coloring are incompatible options. Only one of them can be enabled."));
273 template <
typename Functor,
typename VectorType>
294 template <
typename Functor>
302 template <
typename VectorType>
311 template <
typename VectorType>
353 std::shared_ptr<const MPI_Comm> comm,
360 template <
typename Functor,
typename VectorType>
370 template <
typename Functor>
373 const Functor & func,
382 template <
typename Functor>
385 const Functor & func,
393 template <
typename VectorType>
422 template <
typename VectorType>
522 std::vector<Number *>
JxW;
580 friend class internal::ReinitHelper<dim, Number>;
590 template <
int dim,
typename Number>
600 for (
int d = 0;
d < dim; ++
d)
622 __host__ __device__ constexpr
unsigned int
649 __device__
inline unsigned int
653 threadIdx.x % n_q_points_1d :
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));
668 template <
int dim,
typename Number>
669 __device__
inline unsigned int
671 const unsigned int cell,
673 const unsigned int n_q_points_1d,
674 const unsigned int n_q_points)
677 q_point_id_in_cell<dim>(n_q_points_1d);
687 template <
int dim,
typename Number>
690 const unsigned int cell,
692 const unsigned int n_q_points_1d)
695 q_point_id_in_cell<dim>(n_q_points_1d));
UpdateFlags mapping_update_flags
types::global_dof_index * local_to_global
unsigned int padding_length
std::vector< point_type * > q_points
@ update_quadrature_points
Transformed quadrature points.
std::vector< dim3 > grid_dim
void serial_copy_constrained_values(const VectorType &src, VectorType &dst) const
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)
types::global_dof_index * constrained_dofs
void cell_loop(const Functor &func, const VectorType &src, VectorType &dst) const
SharedData(Number *vd, Number *gq[dim])
ParallelizationScheme parallelization_scheme
void copy_constrained_values(const VectorType &src, VectorType &dst) const
void serial_cell_loop(const Functor &func, const VectorType &src, VectorType &dst) const
unsigned int q_point_id_in_cell(const unsigned int n_q_points_1d)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
std::vector< unsigned int > n_cells
unsigned int cells_per_block
types::global_dof_index n_dofs
unsigned int get_padding_length() const
unsigned int n_constrained_dofs
void reinit(const Mapping< dim > &mapping, const DoFHandler< dim > &dof_handler, const AffineConstraints< Number > &constraints, const Quadrature< 1 > &quad, const AdditionalData &additional_data=AdditionalData())
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)
void set_constrained_values(const Number val, VectorType &dst) const
bool overlap_communication_computation
std::vector< unsigned int > row_start
static ::ExceptionBase & ExcMessage(std::string arg1)
void distributed_copy_constrained_values(const LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &src, LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &dst) const
void serial_set_constrained_values(const Number val, VectorType &dst) const
#define DEAL_II_NAMESPACE_OPEN
std::vector< types::global_dof_index * > local_to_global
std::size_t memory_consumption() const
@ update_gradients
Shape function gradients.
unsigned int q_points_per_cell
void evaluate_coefficients(Functor func) const
unsigned int * constraint_mask
std::shared_ptr< const Utilities::MPI::Partitioner > partitioner
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)
unsigned int dofs_per_cell
void distributed_set_constrained_values(const Number val, LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &dst) const
void distributed_cell_loop(const Functor &func, const LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &src, LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &dst) const
std::vector< dim3 > block_dim
std::vector< unsigned int * > constraint_mask
ParallelizationScheme parallelization_scheme
bool overlap_communication_computation
@ update_JxW_values
Transformed quadrature weights.
dim3 constraint_block_dim
Data get_data(unsigned int color) const
std::vector< Number * > inv_jacobian
#define DEAL_II_NAMESPACE_CLOSE
void initialize_dof_vector(LinearAlgebra::CUDAWrappers::Vector< Number > &vec) const
#define AssertThrow(cond, exc)
constexpr __host__ unsigned int cells_per_block_shmem(int dim, int fe_degree)
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)
unsigned int padding_length
std::vector< Number * > JxW