17#ifndef dealii_cuda_matrix_free_h
18#define dealii_cuda_matrix_free_h
22#ifdef DEAL_II_COMPILER_CUDA_AWARE
50 template <
int dim,
typename Number>
81 template <
int dim,
typename Number =
double>
121# ifndef DEAL_II_MPI_WITH_CUDA_SUPPORT
125 "Overlapping communication and computation requires CUDA-aware MPI."));
131 "Overlapping communication and coloring are incompatible options. Only one of them can be enabled."));
248 template <
typename IteratorFiltersType>
254 const IteratorFiltersType & iterator_filter,
301 template <
typename Functor,
typename VectorType>
304 const VectorType &src,
305 VectorType & dst)
const;
322 template <
typename Functor>
330 template <
typename VectorType>
339 template <
typename VectorType>
363 const std::vector<std::vector<CellFilter>> &
376 const std::shared_ptr<const Utilities::MPI::Partitioner> &
401 template <
typename IteratorFiltersType>
407 const IteratorFiltersType & iterator_filter,
408 std::shared_ptr<const MPI_Comm>
comm,
415 template <
typename Functor,
typename VectorType>
418 const VectorType &src,
419 VectorType & dst)
const;
425 template <
typename Functor>
428 const Functor & func,
437 template <
typename Functor>
440 const Functor & func,
448 template <
typename VectorType>
451 VectorType & dst)
const;
477 template <
typename VectorType>
582 std::vector<Number *>
JxW;
648 std::vector<std::vector<CellFilter>>
graph;
650 friend class internal::ReinitHelper<dim, Number>;
660 template <
int dim,
typename Number>
670 for (
int d = 0;
d < dim; ++
d)
692 __host__ __device__
constexpr unsigned int
719 __device__
inline unsigned int
723 threadIdx.x % n_q_points_1d :
725 threadIdx.x % n_q_points_1d + n_q_points_1d * threadIdx.y :
726 threadIdx.x % n_q_points_1d +
727 n_q_points_1d * (threadIdx.y + n_q_points_1d * threadIdx.z));
738 template <
int dim,
typename Number>
739 __device__
inline unsigned int
741 const unsigned int cell,
743 const unsigned int n_q_points_1d,
744 const unsigned int n_q_points)
747 q_point_id_in_cell<dim>(n_q_points_1d);
757 template <
int dim,
typename Number>
760 const unsigned int cell,
762 const unsigned int n_q_points_1d)
765 q_point_id_in_cell<dim>(n_q_points_1d));
772 template <
int dim,
typename Number>
836 template <
int dim,
typename Number>
839 const typename ::CUDAWrappers::MatrixFree<dim, Number>::Data &data,
844 data_host.
id = data.id;
845 data_host.
n_cells = data.n_cells;
850 const unsigned int n_elements =
854 data_host.
q_points.resize(n_elements);
871 data_host.
JxW.resize(n_elements);
889 template <
int dim,
typename Number>
893 const unsigned int n_q_points,
894 const unsigned int i)
908 template <
int dim,
typename Number>
912 const unsigned int i)
922 template <
int dim,
typename Number>
923 inline const std::vector<std::vector<
932 template <
int dim,
typename Number>
933 inline const std::shared_ptr<const Utilities::MPI::Partitioner> &
941 template <
int dim,
typename Number>
void evaluate_coefficients(Functor func) const
void serial_copy_constrained_values(const VectorType &src, VectorType &dst) const
std::vector< Number * > JxW
std::vector< unsigned int * > constraint_mask
const DoFHandler< dim > * dof_handler
void initialize_dof_vector(LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &vec) const
unsigned int q_points_per_cell
Data get_data(unsigned int color) const
unsigned int get_padding_length() const
void initialize_dof_vector(LinearAlgebra::CUDAWrappers::Vector< Number > &vec) const
bool overlap_communication_computation
void serial_set_constrained_values(const Number val, VectorType &dst) const
std::vector< dim3 > block_dim
unsigned int cells_per_block
void distributed_copy_constrained_values(const LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &src, LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &dst) const
std::vector< unsigned int > n_cells
ParallelizationScheme parallelization_scheme
unsigned int n_constrained_dofs
void internal_reinit(const Mapping< dim > &mapping, const DoFHandler< dim > &dof_handler, const AffineConstraints< Number > &constraints, const Quadrature< 1 > &quad, const IteratorFiltersType &iterator_filter, std::shared_ptr< const MPI_Comm > comm, const AdditionalData additional_data)
const DoFHandler< dim > & get_dof_handler() const
const std::shared_ptr< const Utilities::MPI::Partitioner > & get_vector_partitioner() const
void reinit(const Mapping< dim > &mapping, const DoFHandler< dim > &dof_handler, const AffineConstraints< Number > &constraints, const Quadrature< 1 > &quad, const IteratorFiltersType &iterator_filter, const AdditionalData &additional_data=AdditionalData())
void serial_cell_loop(const Functor &func, const VectorType &src, VectorType &dst) const
std::size_t memory_consumption() const
std::vector< point_type * > q_points
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_set_constrained_values(const Number val, LinearAlgebra::CUDAWrappers::Vector< Number > &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)
unsigned int dofs_per_cell
void distributed_set_constrained_values(const Number val, LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &dst) const
Point< dim, Number > get_quadrature_point_host(const unsigned int cell, const DataHost< dim, Number > &data, const unsigned int i)
void distributed_cell_loop(const Functor &func, const LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &src, LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &dst) const
const std::vector< std::vector< CellFilter > > & get_colored_graph() const
std::vector< dim3 > grid_dim
std::vector< types::global_dof_index * > local_to_global
dim3 constraint_block_dim
std::vector< std::vector< CellFilter > > graph
void copy_constrained_values(const VectorType &src, VectorType &dst) const
std::shared_ptr< const Utilities::MPI::Partitioner > partitioner
unsigned int padding_length
unsigned int q_point_id_in_cell(const unsigned int n_q_points_1d)
void distributed_copy_constrained_values(const LinearAlgebra::CUDAWrappers::Vector< Number > &src, LinearAlgebra::CUDAWrappers::Vector< Number > &dst) const
void reinit(const Mapping< dim > &mapping, const DoFHandler< dim > &dof_handler, const AffineConstraints< Number > &constraints, const Quadrature< 1 > &quad, const AdditionalData &additional_data=AdditionalData())
types::global_dof_index n_dofs
void distributed_cell_loop(const Functor &func, const LinearAlgebra::CUDAWrappers::Vector< Number > &src, LinearAlgebra::CUDAWrappers::Vector< Number > &dst) const
std::vector< unsigned int > row_start
void set_constrained_values(const Number val, VectorType &dst) const
DataHost< dim, Number > copy_mf_data_to_host(const typename::CUDAWrappers::MatrixFree< dim, Number >::Data &data, const UpdateFlags &update_flags)
types::global_dof_index * constrained_dofs
std::vector< Number * > inv_jacobian
unsigned int local_q_point_id_host(const unsigned int cell, const DataHost< dim, Number > &data, const unsigned int n_q_points, const unsigned int i)
void cell_loop(const Functor &func, const VectorType &src, VectorType &dst) const
void reinit(const DoFHandler< dim > &dof_handler, const AffineConstraints< Number > &constraints, const Quadrature< 1 > &quad, const AdditionalData &AdditionalData=AdditionalData())
Abstract base class for mapping classes.
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
@ update_JxW_values
Transformed quadrature weights.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
#define Assert(cond, exc)
static ::ExceptionBase & ExcNotInitialized()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
__host__ constexpr unsigned int cells_per_block_shmem(int dim, int fe_degree)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
void copy_to_host(const ArrayView< const T, MemorySpace::CUDA > &in, ArrayView< T, MemorySpace::Host > &out)
*braid_SplitCommworld & comm
unsigned int padding_length
std::vector< unsigned int > constraint_mask
std::vector< types::global_dof_index > local_to_global
std::vector< Point< dim, Number > > q_points
std::vector< Number > inv_jacobian
std::vector< Number > JxW
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)
ParallelizationScheme parallelization_scheme
UpdateFlags mapping_update_flags
bool overlap_communication_computation
types::global_dof_index * local_to_global
unsigned int * constraint_mask
unsigned int padding_length
SharedData(Number *vd, Number *gq[dim])