Reference documentation for deal.II version 9.5.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\}}\)
Loading...
Searching...
No Matches
Classes | Functions
MatrixFreeTools Namespace Reference

Classes

class  ElementActivationAndDeactivationMatrixFree
 

Functions

template<int dim, typename AdditionalData >
void categorize_by_boundary_ids (const Triangulation< dim > &tria, AdditionalData &additional_data)
 
template<int dim, int fe_degree, int n_q_points_1d, int n_components, typename Number , typename VectorizedArrayType , typename VectorType >
void compute_diagonal (const MatrixFree< dim, Number, VectorizedArrayType > &matrix_free, VectorType &diagonal_global, const std::function< void(FEEvaluation< dim, fe_degree, n_q_points_1d, n_components, Number, VectorizedArrayType > &)> &local_vmult, const unsigned int dof_no=0, const unsigned int quad_no=0, const unsigned int first_selected_component=0)
 
template<typename CLASS , int dim, int fe_degree, int n_q_points_1d, int n_components, typename Number , typename VectorizedArrayType , typename VectorType >
void compute_diagonal (const MatrixFree< dim, Number, VectorizedArrayType > &matrix_free, VectorType &diagonal_global, void(CLASS::*cell_operation)(FEEvaluation< dim, fe_degree, n_q_points_1d, n_components, Number, VectorizedArrayType > &) const, const CLASS *owning_class, const unsigned int dof_no=0, const unsigned int quad_no=0, const unsigned int first_selected_component=0)
 
template<int dim, int fe_degree, int n_q_points_1d, int n_components, typename Number , typename VectorizedArrayType , typename MatrixType >
void compute_matrix (const MatrixFree< dim, Number, VectorizedArrayType > &matrix_free, const AffineConstraints< Number > &constraints, MatrixType &matrix, const std::function< void(FEEvaluation< dim, fe_degree, n_q_points_1d, n_components, Number, VectorizedArrayType > &)> &local_vmult, const unsigned int dof_no=0, const unsigned int quad_no=0, const unsigned int first_selected_component=0)
 
template<typename CLASS , int dim, int fe_degree, int n_q_points_1d, int n_components, typename Number , typename VectorizedArrayType , typename MatrixType >
void compute_matrix (const MatrixFree< dim, Number, VectorizedArrayType > &matrix_free, const AffineConstraints< Number > &constraints, MatrixType &matrix, void(CLASS::*cell_operation)(FEEvaluation< dim, fe_degree, n_q_points_1d, n_components, Number, VectorizedArrayType > &) const, const CLASS *owning_class, const unsigned int dof_no=0, const unsigned int quad_no=0, const unsigned int first_selected_component=0)
 

Detailed Description

A namespace for utility functions in the context of matrix-free operator evaluation.

Function Documentation

◆ categorize_by_boundary_ids()

template<int dim, typename AdditionalData >
void MatrixFreeTools::categorize_by_boundary_ids ( const Triangulation< dim > &  tria,
AdditionalData &  additional_data 
)

Modify additional_data so that cells are categorized according to their boundary IDs, making face integrals in the case of cell-centric loop simpler.

◆ compute_diagonal() [1/2]

template<int dim, int fe_degree, int n_q_points_1d, int n_components, typename Number , typename VectorizedArrayType , typename VectorType >
void MatrixFreeTools::compute_diagonal ( const MatrixFree< dim, Number, VectorizedArrayType > &  matrix_free,
VectorType &  diagonal_global,
const std::function< void(FEEvaluation< dim, fe_degree, n_q_points_1d, n_components, Number, VectorizedArrayType > &)> &  local_vmult,
const unsigned int  dof_no = 0,
const unsigned int  quad_no = 0,
const unsigned int  first_selected_component = 0 
)

Compute the diagonal of a linear operator (diagonal_global), given matrix_free and the local cell integral operation local_vmult. The vector is initialized to the right size in the function.

The parameters dof_no, quad_no, and first_selected_component are passed to the constructor of the FEEvaluation that is internally set up.

◆ compute_diagonal() [2/2]

template<typename CLASS , int dim, int fe_degree, int n_q_points_1d, int n_components, typename Number , typename VectorizedArrayType , typename VectorType >
void MatrixFreeTools::compute_diagonal ( const MatrixFree< dim, Number, VectorizedArrayType > &  matrix_free,
VectorType &  diagonal_global,
void(CLASS::*)(FEEvaluation< dim, fe_degree, n_q_points_1d, n_components, Number, VectorizedArrayType > &) const  cell_operation,
const CLASS *  owning_class,
const unsigned int  dof_no = 0,
const unsigned int  quad_no = 0,
const unsigned int  first_selected_component = 0 
)

Same as above but with a class and a function pointer.

◆ compute_matrix() [1/2]

template<int dim, int fe_degree, int n_q_points_1d, int n_components, typename Number , typename VectorizedArrayType , typename MatrixType >
void MatrixFreeTools::compute_matrix ( const MatrixFree< dim, Number, VectorizedArrayType > &  matrix_free,
const AffineConstraints< Number > &  constraints,
MatrixType &  matrix,
const std::function< void(FEEvaluation< dim, fe_degree, n_q_points_1d, n_components, Number, VectorizedArrayType > &)> &  local_vmult,
const unsigned int  dof_no = 0,
const unsigned int  quad_no = 0,
const unsigned int  first_selected_component = 0 
)

Compute the matrix representation of a linear operator (matrix), given matrix_free and the local cell integral operation local_vmult. Constrained entries on the diagonal are set to one.

The parameters dof_no, quad_no, and first_selected_component are passed to the constructor of the FEEvaluation that is internally set up.

◆ compute_matrix() [2/2]

template<typename CLASS , int dim, int fe_degree, int n_q_points_1d, int n_components, typename Number , typename VectorizedArrayType , typename MatrixType >
void MatrixFreeTools::compute_matrix ( const MatrixFree< dim, Number, VectorizedArrayType > &  matrix_free,
const AffineConstraints< Number > &  constraints,
MatrixType &  matrix,
void(CLASS::*)(FEEvaluation< dim, fe_degree, n_q_points_1d, n_components, Number, VectorizedArrayType > &) const  cell_operation,
const CLASS *  owning_class,
const unsigned int  dof_no = 0,
const unsigned int  quad_no = 0,
const unsigned int  first_selected_component = 0 
)

Same as above but with a class and a function pointer.