17 #ifndef dealii_cuda_matrix_free_h 18 #define dealii_cuda_matrix_free_h 20 #include <deal.II/base/config.h> 22 #ifdef DEAL_II_COMPILER_CUDA_AWARE 24 # include <deal.II/base/mpi.h> 25 # include <deal.II/base/quadrature.h> 26 # include <deal.II/base/tensor.h> 28 # include <deal.II/dofs/dof_handler.h> 30 # include <deal.II/fe/fe_update_flags.h> 31 # include <deal.II/fe/mapping.h> 32 # include <deal.II/fe/mapping_q1.h> 34 # include <deal.II/lac/affine_constraints.h> 35 # include <deal.II/lac/cuda_vector.h> 36 # include <deal.II/lac/la_parallel_vector.h> 39 DEAL_II_NAMESPACE_OPEN
46 template <
int dim,
typename Number>
76 template <
int dim,
typename Number =
double>
101 , mapping_update_flags(mapping_update_flags)
133 Number * inv_jacobian;
135 unsigned int n_cells;
136 unsigned int padding_length;
137 unsigned int row_start;
138 unsigned int * constraint_mask;
165 const AdditionalData additional_data = AdditionalData());
174 const AdditionalData AdditionalData = AdditionalData());
201 template <
typename Functor,
typename VectorType>
204 const VectorType &src,
205 VectorType & dst)
const;
222 template <
typename Functor>
230 template <
typename VectorType>
239 template <
typename VectorType>
264 std::shared_ptr<const MPI_Comm> comm,
265 const AdditionalData additional_data);
271 template <
typename Functor,
typename VectorType>
274 const VectorType &src,
275 VectorType & dst)
const;
281 template <
typename Functor>
284 const Functor & func,
293 template <
typename Functor>
296 const Functor & func,
304 template <
typename VectorType>
307 VectorType & dst)
const;
333 template <
typename VectorType>
415 std::vector<Number *>
JxW;
421 std::vector<unsigned int *> constraint_mask;
442 unsigned int cells_per_block;
443 dim3 constraint_grid_dim;
444 dim3 constraint_block_dim;
446 unsigned int padding_length;
447 std::vector<unsigned int> row_start;
449 friend class internal::ReinitHelper<dim, Number>;
456 template <
int dim,
typename Number>
460 SharedData(Number *vd, Number *gq[dim])
463 for (
int d = 0; d < dim; ++d)
464 gradients[d] = gq[d];
468 Number *gradients[dim];
476 __host__ __device__ constexpr
unsigned int 477 cells_per_block_shmem(
int dim,
int fe_degree)
480 return dim==2 ? (fe_degree==1 ? 32 :
485 dim==3 ? (fe_degree==1 ? 8 :
499 __device__
inline unsigned int 503 threadIdx.x % n_q_points_1d :
505 threadIdx.x % n_q_points_1d + n_q_points_1d * threadIdx.y :
506 threadIdx.x % n_q_points_1d +
507 n_q_points_1d * (threadIdx.y + n_q_points_1d * threadIdx.z));
518 template <
int dim,
typename Number>
519 __device__
inline unsigned int 521 const unsigned int cell,
523 const unsigned int n_q_points_1d,
524 const unsigned int n_q_points)
526 return (data->row_start / data->padding_length + cell) * n_q_points +
527 q_point_id_in_cell<dim>(n_q_points_1d);
537 template <
int dim,
typename Number>
540 const unsigned int cell,
542 const unsigned int n_q_points_1d)
544 return *(data->q_points + data->padding_length * cell +
545 q_point_id_in_cell<dim>(n_q_points_1d));
549 DEAL_II_NAMESPACE_CLOSE
Transformed quadrature weights.
ParallelizationScheme parallelization_scheme
std::shared_ptr< const Utilities::MPI::Partitioner > partitioner
void distributed_cell_loop(const Functor &func, const LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &src, LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &dst) const
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)
void serial_set_constrained_values(const Number val, VectorType &dst) const
void serial_copy_constrained_values(const VectorType &src, VectorType &dst) const
void distributed_set_constrained_values(const Number val, LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &dst) const
Data get_data(unsigned int color) 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)
unsigned int q_points_per_cell
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 distributed_copy_constrained_values(const LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &src, LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &dst) const
std::vector< point_type * > q_points
void serial_cell_loop(const Functor &func, const VectorType &src, VectorType &dst) const
unsigned int dofs_per_cell
std::vector< Number * > inv_jacobian
std::vector< dim3 > block_dim
void reinit(const Mapping< dim > &mapping, const DoFHandler< dim > &dof_handler, const AffineConstraints< Number > &constraints, const Quadrature< 1 > &quad, const AdditionalData additional_data=AdditionalData())
unsigned int global_dof_index
std::vector< dim3 > grid_dim
std::vector< types::global_dof_index * > local_to_global
void copy_constrained_values(const VectorType &src, VectorType &dst) const
std::size_t memory_consumption() const
Shape function gradients.
unsigned int get_padding_length() const
void set_constrained_values(const Number val, VectorType &dst) const
std::vector< Number * > JxW
std::vector< unsigned int > n_cells
unsigned int q_point_id_in_cell(const unsigned int n_q_points_1d)
void evaluate_coefficients(Functor func) const
void cell_loop(const Functor &func, const VectorType &src, VectorType &dst) const
unsigned int n_constrained_dofs