Reference documentation for deal.II version 9.4.1
|
#include <deal.II/matrix_free/cuda_matrix_free.h>
Classes | |
struct | AdditionalData |
struct | Data |
Public Types | |
enum | ParallelizationScheme { parallel_in_elem , parallel_over_elem } |
using | jacobian_type = Tensor< 2, dim, Tensor< 1, dim, Number > > |
using | point_type = Point< dim, Number > |
using | CellFilter = FilteredIterator< typename DoFHandler< dim >::active_cell_iterator > |
Public Member Functions | |
MatrixFree () | |
~MatrixFree () | |
unsigned int | get_padding_length () const |
template<typename IteratorFiltersType > | |
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 | reinit (const Mapping< dim > &mapping, const DoFHandler< dim > &dof_handler, const AffineConstraints< Number > &constraints, const Quadrature< 1 > &quad, const AdditionalData &additional_data=AdditionalData()) |
void | reinit (const DoFHandler< dim > &dof_handler, const AffineConstraints< Number > &constraints, const Quadrature< 1 > &quad, const AdditionalData &additional_data=AdditionalData()) |
Data | get_data (unsigned int color) const |
template<typename Functor , typename VectorType > | |
void | cell_loop (const Functor &func, const VectorType &src, VectorType &dst) const |
template<typename Functor > | |
void | evaluate_coefficients (Functor func) const |
template<typename VectorType > | |
void | copy_constrained_values (const VectorType &src, VectorType &dst) const |
template<typename VectorType > | |
void | set_constrained_values (const Number val, VectorType &dst) const |
void | initialize_dof_vector (LinearAlgebra::CUDAWrappers::Vector< Number > &vec) const |
void | initialize_dof_vector (LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &vec) const |
const std::vector< std::vector< CellFilter > > & | get_colored_graph () const |
const std::shared_ptr< const Utilities::MPI::Partitioner > & | get_vector_partitioner () const |
void | free () |
const DoFHandler< dim > & | get_dof_handler () const |
std::size_t | memory_consumption () const |
Private Member Functions | |
template<typename IteratorFiltersType > | |
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) |
template<typename Functor , typename VectorType > | |
void | serial_cell_loop (const Functor &func, const VectorType &src, VectorType &dst) const |
template<typename Functor > | |
void | distributed_cell_loop (const Functor &func, const LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &src, LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &dst) const |
template<typename Functor > | |
void | distributed_cell_loop (const Functor &func, const LinearAlgebra::CUDAWrappers::Vector< Number > &src, LinearAlgebra::CUDAWrappers::Vector< Number > &dst) const |
template<typename VectorType > | |
void | serial_copy_constrained_values (const VectorType &src, VectorType &dst) const |
void | distributed_copy_constrained_values (const LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &src, LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &dst) const |
void | distributed_copy_constrained_values (const LinearAlgebra::CUDAWrappers::Vector< Number > &src, LinearAlgebra::CUDAWrappers::Vector< Number > &dst) const |
template<typename VectorType > | |
void | serial_set_constrained_values (const Number val, VectorType &dst) const |
void | distributed_set_constrained_values (const Number val, LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > &dst) const |
void | distributed_set_constrained_values (const Number val, LinearAlgebra::CUDAWrappers::Vector< Number > &dst) const |
Private Attributes | |
int | my_id |
ParallelizationScheme | parallelization_scheme |
bool | use_coloring |
bool | overlap_communication_computation |
types::global_dof_index | n_dofs |
unsigned int | fe_degree |
unsigned int | dofs_per_cell |
unsigned int | n_constrained_dofs |
unsigned int | q_points_per_cell |
unsigned int | n_colors |
std::vector< unsigned int > | n_cells |
std::vector< point_type * > | q_points |
std::vector< types::global_dof_index * > | local_to_global |
std::vector< Number * > | inv_jacobian |
std::vector< Number * > | JxW |
types::global_dof_index * | constrained_dofs |
std::vector<::internal::MatrixFreeFunctions::ConstraintKinds * > | constraint_mask |
std::vector< dim3 > | grid_dim |
std::vector< dim3 > | block_dim |
std::shared_ptr< const Utilities::MPI::Partitioner > | partitioner |
unsigned int | cells_per_block |
dim3 | constraint_grid_dim |
dim3 | constraint_block_dim |
unsigned int | padding_length |
std::vector< unsigned int > | row_start |
const DoFHandler< dim > * | dof_handler |
std::vector< std::vector< CellFilter > > | graph |
Friends | |
class | internal::ReinitHelper< dim, Number > |
Related Functions | |
(Note that these are not member functions.) | |
template<int dim> | |
unsigned int | q_point_id_in_cell (const unsigned int n_q_points_1d) |
template<int dim, typename Number > | |
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) |
template<int dim, typename Number > | |
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) |
template<int dim, typename Number > | |
DataHost< dim, Number > | copy_mf_data_to_host (const typename::CUDAWrappers::MatrixFree< dim, Number >::Data &data, const UpdateFlags &update_flags) |
template<int dim, typename Number > | |
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) |
template<int dim, typename Number > | |
Point< dim, Number > | get_quadrature_point_host (const unsigned int cell, const DataHost< dim, Number > &data, const unsigned int i) |
Subscriptor functionality | |
Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class. | |
std::atomic< unsigned int > | counter |
std::map< std::string, unsigned int > | counter_map |
std::vector< std::atomic< bool > * > | validity_pointers |
const std::type_info * | object_info |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
unsigned int | n_subscriptions () const |
template<typename StreamType > | |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
void | check_no_subscribers () const noexcept |
using | map_value_type = decltype(counter_map)::value_type |
using | map_iterator = decltype(counter_map)::iterator |
static std::mutex | mutex |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
This class collects all the data that is stored for the matrix free implementation. The storage scheme is tailored towards several loops performed with the same data, i.e., typically doing many matrix-vector products or residual computations on the same mesh.
This class does not implement any operations involving finite element basis functions, i.e., regarding the operation performed on the cells. For these operations, the class FEEvaluation is designed to use the data collected in this class.
This class implements a loop over all cells (cell_loop()). This loop is scheduled in such a way that cells that share degrees of freedom are not worked on simultaneously, which implies that it is possible to write to vectors in parallel without having to explicitly synchronize access to these vectors and matrices. This class does not implement any shape values, all it does is to cache the respective data. To implement finite element operations, use the class CUDAWrappers::FEEvalutation.
This class traverse the cells in a different order than the usual Triangulation class in deal.II.
Definition at line 92 of file cuda_matrix_free.h.
using CUDAWrappers::MatrixFree< dim, Number >::jacobian_type = Tensor<2, dim, Tensor<1, dim, Number> > |
Definition at line 95 of file cuda_matrix_free.h.
using CUDAWrappers::MatrixFree< dim, Number >::point_type = Point<dim, Number> |
Definition at line 96 of file cuda_matrix_free.h.
using CUDAWrappers::MatrixFree< dim, Number >::CellFilter = FilteredIterator<typename DoFHandler<dim>::active_cell_iterator> |
Definition at line 97 of file cuda_matrix_free.h.
enum CUDAWrappers::MatrixFree::ParallelizationScheme |
Parallelization scheme used: parallel_in_elem (parallelism at the level of degrees of freedom) or parallel_over_elem (parallelism at the level of cells)
Enumerator | |
---|---|
parallel_in_elem | |
parallel_over_elem |
Definition at line 105 of file cuda_matrix_free.h.
CUDAWrappers::MatrixFree< dim, Number >::MatrixFree | ( | ) |
Default constructor.
CUDAWrappers::MatrixFree< dim, Number >::~MatrixFree | ( | ) |
Destructor.
unsigned int CUDAWrappers::MatrixFree< dim, Number >::get_padding_length | ( | ) | const |
Return the length of the padding.
void CUDAWrappers::MatrixFree< dim, Number >::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() |
||
) |
Extracts the information needed to perform loops over cells. The DoFHandler and AffineConstraints objects describe the layout of degrees of freedom, the DoFHandler and the mapping describe the transformation from unit to real cell, and the finite element underlying the DoFHandler together with the quadrature formula describe the local operations. This function takes an IteratorFilters object (predicate) to loop over a subset of the active cells. When using MPI, the predicate should filter out non locally owned cells.
void CUDAWrappers::MatrixFree< dim, Number >::reinit | ( | const Mapping< dim > & | mapping, |
const DoFHandler< dim > & | dof_handler, | ||
const AffineConstraints< Number > & | constraints, | ||
const Quadrature< 1 > & | quad, | ||
const AdditionalData & | additional_data = AdditionalData() |
||
) |
Same as above using Iterators::LocallyOwnedCell() as predicate.
void CUDAWrappers::MatrixFree< dim, Number >::reinit | ( | const DoFHandler< dim > & | dof_handler, |
const AffineConstraints< Number > & | constraints, | ||
const Quadrature< 1 > & | quad, | ||
const AdditionalData & | additional_data = AdditionalData() |
||
) |
Initializes the data structures. Same as above but using a Q1 mapping.
Data CUDAWrappers::MatrixFree< dim, Number >::get_data | ( | unsigned int | color | ) | const |
Return the Data structure associated with color
.
void CUDAWrappers::MatrixFree< dim, Number >::cell_loop | ( | const Functor & | func, |
const VectorType & | src, | ||
VectorType & | dst | ||
) | const |
This method runs the loop over all cells and apply the local operation on each element in parallel. func
is a functor which is applied on each color.
func
needs to define
void CUDAWrappers::MatrixFree< dim, Number >::evaluate_coefficients | ( | Functor | func | ) | const |
This method runs the loop over all cells and apply the local operation on each element in parallel. This function is very similar to cell_loop() but it uses a simpler functor.
func
needs to define
void CUDAWrappers::MatrixFree< dim, Number >::copy_constrained_values | ( | const VectorType & | src, |
VectorType & | dst | ||
) | const |
Copy the values of the constrained entries from src
to dst
. This is used to impose zero Dirichlet boundary condition.
void CUDAWrappers::MatrixFree< dim, Number >::set_constrained_values | ( | const Number | val, |
VectorType & | dst | ||
) | const |
Set the entries in dst
corresponding to constrained values to val
. The main purpose of this function is to set the constrained entries of the source vector used in cell_loop() to zero.
void CUDAWrappers::MatrixFree< dim, Number >::initialize_dof_vector | ( | LinearAlgebra::CUDAWrappers::Vector< Number > & | vec | ) | const |
Initialize a serial vector. The size corresponds to the number of degrees of freedom in the DoFHandler object.
void CUDAWrappers::MatrixFree< dim, Number >::initialize_dof_vector | ( | LinearAlgebra::distributed::Vector< Number, MemorySpace::CUDA > & | vec | ) | const |
Initialize a distributed vector. The local elements correspond to the locally owned degrees of freedom and the ghost elements correspond to the (additional) locally relevant dofs.
const std::vector< std::vector< CellFilter > > & CUDAWrappers::MatrixFree< dim, Number >::get_colored_graph | ( | ) | const |
Return the colored graph of locally owned active cells.
const std::shared_ptr< const Utilities::MPI::Partitioner > & CUDAWrappers::MatrixFree< dim, Number >::get_vector_partitioner | ( | ) | const |
Return the partitioner that represents the locally owned data and the ghost indices where access is needed to for the cell loop. The partitioner is constructed from the locally owned dofs and ghost dofs given by the respective fields. If you want to have specific information about these objects, you can query them with the respective access functions. If you just want to initialize a (parallel) vector, you should usually prefer this data structure as the data exchange information can be reused from one vector to another.
void CUDAWrappers::MatrixFree< dim, Number >::free | ( | ) |
Free all the memory allocated.
const DoFHandler< dim > & CUDAWrappers::MatrixFree< dim, Number >::get_dof_handler | ( | ) | const |
Return the DoFHandler.
std::size_t CUDAWrappers::MatrixFree< dim, Number >::memory_consumption | ( | ) | const |
Return an approximation of the memory consumption of this class in bytes.
|
private |
Initializes the data structures.
|
private |
Helper function. Loop over all the cells and apply the functor on each element in parallel. This function is used when MPI is not used.
|
private |
Helper function. Loop over all the cells and apply the functor on each element in parallel. This function is used when MPI is used.
|
private |
This function should never be called. Calling it results in an internal error. This function exists only because cell_loop needs distributed_cell_loop() to exist for LinearAlgebra::CUDAWrappers::Vector.
|
private |
Helper function. Copy the values of the constrained entries of src
to dst
. This function is used when MPI is not used.
|
private |
Helper function. Copy the values of the constrained entries of src
to dst
. This function is used when MPI is used.
|
private |
This function should never be called. Calling it results in an internal error. This function exists only because copy_constrained_values needs distributed_copy_constrained_values() to exist for LinearAlgebra::CUDAWrappers::Vector.
|
private |
Helper function. Set the constrained entries of dst
to val
. This function is used when MPI is not used.
|
private |
Helper function. Set the constrained entries of dst
to val
. This function is used when MPI is used.
|
private |
This function should never be called. Calling it results in an internal error. This function exists only because set_constrained_values needs distributed_set_constrained_values() to exist for LinearAlgebra::CUDAWrappers::Vector.
|
friend |
Definition at line 659 of file cuda_matrix_free.h.
Compute the quadrature point index in the local cell of a given thread.
Definition at line 731 of file cuda_matrix_free.h.
|
related |
Return the quadrature point index local of a given thread. The index is only unique for a given MPI process.
Definition at line 750 of file cuda_matrix_free.h.
|
related |
Return the quadrature point associated with a given thread.
Definition at line 769 of file cuda_matrix_free.h.
|
related |
Copy data
from the device to the device. update_flags
should be identical to the one used in MatrixFree::AdditionalData.
Definition at line 849 of file cuda_matrix_free.h.
|
related |
This function is the host version of local_q_point_id().
Definition at line 902 of file cuda_matrix_free.h.
|
related |
This function is the host version of get_quadrature_point(). It assumes that the data in MatrixFree<dim, Number>::Data has been copied to the host using copy_mf_data_to_host().
Definition at line 921 of file cuda_matrix_free.h.
|
private |
Unique ID associated with the object.
Definition at line 514 of file cuda_matrix_free.h.
|
private |
Parallelization scheme used, parallelization over degrees of freedom or over cells.
Definition at line 520 of file cuda_matrix_free.h.
|
private |
If true, use graph coloring. Otherwise, use atomic operations. Graph coloring ensures bitwise reproducibility but is slower on Pascal and newer architectures.
Definition at line 527 of file cuda_matrix_free.h.
|
private |
Overlap MPI communications with computation. This requires CUDA-aware MPI and use_coloring must be false.
Definition at line 533 of file cuda_matrix_free.h.
|
private |
Total number of degrees of freedom.
Definition at line 538 of file cuda_matrix_free.h.
|
private |
Degree of the finite element used.
Definition at line 543 of file cuda_matrix_free.h.
|
private |
Number of degrees of freedom per cell.
Definition at line 548 of file cuda_matrix_free.h.
|
private |
Number of constrained degrees of freedom.
Definition at line 553 of file cuda_matrix_free.h.
|
private |
Number of quadrature points per cells.
Definition at line 558 of file cuda_matrix_free.h.
|
private |
Number of colors produced by the graph coloring algorithm.
Definition at line 563 of file cuda_matrix_free.h.
|
private |
Number of cells in each color.
Definition at line 568 of file cuda_matrix_free.h.
|
private |
Vector of pointers to the quadrature points associated to the cells of each color.
Definition at line 574 of file cuda_matrix_free.h.
|
private |
Map the position in the local vector to the position in the global vector.
Definition at line 580 of file cuda_matrix_free.h.
|
private |
Vector of pointer to the inverse Jacobian associated to the cells of each color.
Definition at line 586 of file cuda_matrix_free.h.
|
private |
Vector of pointer to the Jacobian times the weights associated to the cells of each color.
Definition at line 592 of file cuda_matrix_free.h.
|
private |
Pointer to the constrained degrees of freedom.
Definition at line 597 of file cuda_matrix_free.h.
|
private |
Mask deciding where constraints are set on a given cell.
Definition at line 603 of file cuda_matrix_free.h.
|
private |
Grid dimensions associated to the different colors. The grid dimensions are used to launch the CUDA kernels.
Definition at line 609 of file cuda_matrix_free.h.
|
private |
Block dimensions associated to the different colors. The block dimensions are used to launch the CUDA kernels.
Definition at line 615 of file cuda_matrix_free.h.
|
private |
Shared pointer to a Partitioner for distributed Vectors used in cell_loop. When MPI is not used the pointer is null.
Definition at line 621 of file cuda_matrix_free.h.
|
private |
Cells per block (determined by the function cells_per_block_shmem() ).
Definition at line 626 of file cuda_matrix_free.h.
|
private |
Grid dimensions used to launch the CUDA kernels in *_constrained_values-operations.
Definition at line 632 of file cuda_matrix_free.h.
|
private |
Block dimensions used to launch the CUDA kernels in *_constrained_values-operations.
Definition at line 638 of file cuda_matrix_free.h.
|
private |
Length of the padding (closest power of two larger than or equal to the number of thread).
Definition at line 644 of file cuda_matrix_free.h.
|
private |
Row start of each color.
Definition at line 649 of file cuda_matrix_free.h.
|
private |
Pointer to the DoFHandler associated with the object.
Definition at line 654 of file cuda_matrix_free.h.
|
private |
Colored graphed of locally owned active cells.
Definition at line 659 of file cuda_matrix_free.h.