17#ifndef dealii_cuda_matrix_free_h
18#define dealii_cuda_matrix_free_h
41#include <Kokkos_Core.hpp>
50 namespace MatrixFreeFunctions
62 template <
int dim,
typename Number>
93 template <
int dim,
typename Number =
double>
119#ifndef DEAL_II_MPI_WITH_DEVICE_SUPPORT
123 "Overlapping communication and computation requires CUDA-aware MPI."));
129 "Overlapping communication and coloring are incompatible options. Only one of them can be enabled."));
165 Kokkos::View<point_type **, MemorySpace::Default::kokkos_space>
q_points;
178 Kokkos::View<Number **[dim][dim], MemorySpace::Default::kokkos_space>
184 Kokkos::View<Number **, MemorySpace::Default::kokkos_space>
JxW;
196 Kokkos::View<Number *, MemorySpace::Default::kokkos_space>
shape_values;
201 Kokkos::View<Number *, MemorySpace::Default::kokkos_space>
207 Kokkos::View<Number *, MemorySpace::Default::kokkos_space>
213 Kokkos::View<Number *, MemorySpace::Default::kokkos_space>
243 const unsigned int n_q_points,
244 const unsigned int q_point)
const
256 const unsigned int q_point)
const
283 template <
typename IteratorFiltersType>
289 const IteratorFiltersType & iterator_filter,
336 template <
typename Functor,
typename VectorType>
339 const VectorType &src,
340 VectorType & dst)
const;
357 template <
typename Functor>
365 template <
typename VectorType>
374 template <
typename VectorType>
378#ifdef DEAL_II_WITH_CUDA
401 const std::vector<std::vector<CellFilter>> &
414 const std::shared_ptr<const Utilities::MPI::Partitioner> &
433 template <
typename IteratorFiltersType>
439 const IteratorFiltersType & iterator_filter,
440 const std::shared_ptr<const MPI_Comm> &
comm,
447 template <
typename Functor,
typename VectorType>
450 const VectorType &src,
451 VectorType & dst)
const;
457 template <
typename Functor>
466#ifdef DEAL_II_WITH_CUDA
472 template <
typename Functor>
475 const Functor & func,
537 std::vector<Kokkos::View<point_type **, MemorySpace::Default::kokkos_space>>
553 Kokkos::View<Number **[dim][dim], MemorySpace::Default::kokkos_space>>
560 std::vector<Kokkos::View<Number **, MemorySpace::Default::kokkos_space>>
566 Kokkos::View<types::global_dof_index *, MemorySpace::Default::kokkos_space>
580 Kokkos::View<Number *, MemorySpace::Default::kokkos_space>
shape_values;
590 Kokkos::View<Number *, MemorySpace::Default::kokkos_space>
596 Kokkos::View<Number *, MemorySpace::Default::kokkos_space>
625 std::vector<std::vector<CellFilter>>
graph;
627 friend class internal::ReinitHelper<dim, Number>;
632 template <
int dim,
typename Number>
636 MemorySpace::Default::kokkos_space::execution_space>::member_type;
640 MemorySpace::Default::kokkos_space::execution_space::scratch_memory_space,
641 Kokkos::MemoryTraits<Kokkos::Unmanaged>>;
644 MemorySpace::Default::kokkos_space::execution_space::scratch_memory_space,
645 Kokkos::MemoryTraits<Kokkos::Unmanaged>>;
678 template <
int dim,
typename Number>
684 typename Kokkos::View<Point<dim, Number> **,
699 typename Kokkos::View<Number **[dim][dim],
706 typename Kokkos::View<Number **,
727 typename Kokkos::View<
744 const unsigned int n_q_points,
745 const unsigned int q_point)
const
757 const unsigned int q_point)
const
771 template <
int dim,
typename Number>
772 DataHost<dim, Number>
774 const typename ::CUDAWrappers::MatrixFree<dim, Number>::Data &data,
779 data_host.
n_cells = data.n_cells;
786 data_host.
q_points = Kokkos::create_mirror(data.q_points);
787 Kokkos::deep_copy(data_host.
q_points, data.q_points);
790 data_host.
local_to_global = Kokkos::create_mirror(data.local_to_global);
795 data_host.
inv_jacobian = Kokkos::create_mirror(data.inv_jacobian);
796 Kokkos::deep_copy(data_host.
inv_jacobian, data.inv_jacobian);
801 data_host.
JxW = Kokkos::create_mirror(data.JxW);
802 Kokkos::deep_copy(data_host.
JxW, data.JxW);
805 data_host.
constraint_mask = Kokkos::create_mirror(data.constraint_mask);
816 template <
int dim,
typename Number>
817 inline const std::vector<std::vector<
826 template <
int dim,
typename Number>
827 inline const std::shared_ptr<const Utilities::MPI::Partitioner> &
835 template <
int dim,
typename Number>
std::vector< Kokkos::View< point_type **, MemorySpace::Default::kokkos_space > > q_points
void evaluate_coefficients(Functor func) const
Kokkos::View< Number *, MemorySpace::Default::kokkos_space > constraint_weights
const DoFHandler< dim > * dof_handler
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 reinit(const DoFHandler< dim > &dof_handler, const AffineConstraints< Number > &constraints, const Quadrature< 1 > &quad, const AdditionalData &additional_data=AdditionalData())
std::vector< unsigned int > n_cells
std::vector< Kokkos::View< Number **[dim][dim], MemorySpace::Default::kokkos_space > > inv_jacobian
unsigned int n_constrained_dofs
std::vector< Kokkos::View< Number **, MemorySpace::Default::kokkos_space > > JxW
const DoFHandler< dim > & get_dof_handler() const
const std::shared_ptr< const Utilities::MPI::Partitioner > & get_vector_partitioner() const
Kokkos::View< Number *, MemorySpace::Default::kokkos_space > co_shape_gradients
Kokkos::View< Number *, MemorySpace::Default::kokkos_space > shape_gradients
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
void distributed_cell_loop(const Functor &func, const LinearAlgebra::distributed::Vector< Number, MemorySpace::Default > &src, LinearAlgebra::distributed::Vector< Number, MemorySpace::Default > &dst) const
Kokkos::View< Number *, MemorySpace::Default::kokkos_space > shape_values
unsigned int dofs_per_cell
void initialize_dof_vector(LinearAlgebra::distributed::Vector< Number, MemorySpace::Default > &vec) const
const std::vector< std::vector< CellFilter > > & get_colored_graph() const
Kokkos::View< types::global_dof_index *, MemorySpace::Default::kokkos_space > constrained_dofs
std::vector< std::vector< CellFilter > > graph
std::vector< Kokkos::View< types::global_dof_index **, MemorySpace::Default::kokkos_space > > local_to_global
void copy_constrained_values(const VectorType &src, VectorType &dst) const
std::vector< Kokkos::View<::internal::MatrixFreeFunctions::ConstraintKinds *, MemorySpace::Default::kokkos_space > > constraint_mask
std::shared_ptr< const Utilities::MPI::Partitioner > partitioner
unsigned int padding_length
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
void internal_reinit(const Mapping< dim > &mapping, const DoFHandler< dim > &dof_handler, const AffineConstraints< Number > &constraints, const Quadrature< 1 > &quad, const IteratorFiltersType &iterator_filter, const std::shared_ptr< const MPI_Comm > &comm, const AdditionalData additional_data)
DataHost< dim, Number > copy_mf_data_to_host(const typename::CUDAWrappers::MatrixFree< dim, Number >::Data &data, const UpdateFlags &update_flags)
void cell_loop(const Functor &func, const VectorType &src, VectorType &dst) const
Abstract base class for mapping classes.
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
static ::ExceptionBase & ExcNotInitialized()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
@ update_JxW_values
Transformed quadrature weights.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
#define DEAL_II_HOST_DEVICE
Kokkos::View< Number **[dim][dim], MemorySpace::Default::kokkos_space >::HostMirror inv_jacobian
unsigned int padding_length
Kokkos::View<::internal::MatrixFreeFunctions::ConstraintKinds *, MemorySpace::Default::kokkos_space >::HostMirror constraint_mask
Kokkos::View< types::global_dof_index **, MemorySpace::Default::kokkos_space >::HostMirror local_to_global
Point< dim, Number > get_quadrature_point(const unsigned int cell, const unsigned int q_point) const
unsigned int local_q_point_id(const unsigned int cell, const unsigned int n_q_points, const unsigned int q_point) const
Kokkos::View< Number **, MemorySpace::Default::kokkos_space >::HostMirror JxW
Kokkos::View< Point< dim, Number > **, MemorySpace::Default::kokkos_space >::HostMirror q_points
UpdateFlags mapping_update_flags
bool overlap_communication_computation
AdditionalData(const UpdateFlags mapping_update_flags=update_gradients|update_JxW_values|update_quadrature_points, const bool use_coloring=false, const bool overlap_communication_computation=false)
Kokkos::View< types::global_dof_index **, MemorySpace::Default::kokkos_space > local_to_global
Kokkos::View< Number **[dim][dim], MemorySpace::Default::kokkos_space > inv_jacobian
Kokkos::View<::internal::MatrixFreeFunctions::ConstraintKinds *, MemorySpace::Default::kokkos_space > constraint_mask
CUDAWrappers::MatrixFree< dim, Number >::point_type & get_quadrature_point(const unsigned int cell, const unsigned int q_point) const
Kokkos::View< Number **, MemorySpace::Default::kokkos_space > JxW
Kokkos::View< Number *, MemorySpace::Default::kokkos_space > co_shape_gradients
Kokkos::View< Number *, MemorySpace::Default::kokkos_space > constraint_weights
unsigned int padding_length
unsigned int local_q_point_id(const unsigned int cell, const unsigned int n_q_points, const unsigned int q_point) const
Kokkos::View< Number *, MemorySpace::Default::kokkos_space > shape_gradients
Kokkos::View< Number *, MemorySpace::Default::kokkos_space > shape_values
Kokkos::View< point_type **, MemorySpace::Default::kokkos_space > q_points
Kokkos::View< Number *, MemorySpace::Default::kokkos_space::execution_space::scratch_memory_space, Kokkos::MemoryTraits< Kokkos::Unmanaged > > SharedView1D
SharedData(const TeamHandle &team_member, const SharedView1D &values, const SharedView2D &gradients)
Kokkos::TeamPolicy< MemorySpace::Default::kokkos_space::execution_space >::member_type TeamHandle
Kokkos::View< Number *[dim], MemorySpace::Default::kokkos_space::execution_space::scratch_memory_space, Kokkos::MemoryTraits< Kokkos::Unmanaged > > SharedView2D
::Kokkos::DefaultExecutionSpace::memory_space kokkos_space