Reference documentation for deal.II version 9.2.0
|
#include <deal.II/matrix_free/operators.h>
Public Types | |
using | value_type = typename Base< dim, VectorType, VectorizedArrayType >::value_type |
using | size_type = typename Base< dim, VectorType, VectorizedArrayType >::size_type |
Public Types inherited from MatrixFreeOperators::Base< dim, LinearAlgebra::distributed::Vector< double >, VectorizedArray< typename LinearAlgebra::distributed::Vector< double > ::value_type > > | |
using | value_type = typename LinearAlgebra::distributed::Vector< double > ::value_type |
using | size_type = typename LinearAlgebra::distributed::Vector< double > ::size_type |
Public Member Functions | |
LaplaceOperator () | |
virtual void | compute_diagonal () override |
void | set_coefficient (const std::shared_ptr< Table< 2, VectorizedArrayType >> &scalar_coefficient) |
virtual void | clear () override |
std::shared_ptr< Table< 2, VectorizedArrayType > > | get_coefficient () |
Public Member Functions inherited from MatrixFreeOperators::Base< dim, LinearAlgebra::distributed::Vector< double >, VectorizedArray< typename LinearAlgebra::distributed::Vector< double > ::value_type > > | |
Base () | |
virtual | ~Base () override=default |
void | initialize (std::shared_ptr< const MatrixFree< dim, value_type, VectorizedArray< typename LinearAlgebra::distributed::Vector< double > ::value_type > >> data, const std::vector< unsigned int > &selected_row_blocks=std::vector< unsigned int >(), const std::vector< unsigned int > &selected_column_blocks=std::vector< unsigned int >()) |
void | initialize (std::shared_ptr< const MatrixFree< dim, value_type, VectorizedArray< typename LinearAlgebra::distributed::Vector< double > ::value_type > >> data, const MGConstrainedDoFs &mg_constrained_dofs, const unsigned int level, const std::vector< unsigned int > &selected_row_blocks=std::vector< unsigned int >()) |
void | initialize (std::shared_ptr< const MatrixFree< dim, value_type, VectorizedArray< typename LinearAlgebra::distributed::Vector< double > ::value_type > >> data_, const std::vector< MGConstrainedDoFs > &mg_constrained_dofs, const unsigned int level, const std::vector< unsigned int > &selected_row_blocks=std::vector< unsigned int >()) |
size_type | m () const |
size_type | n () const |
void | vmult_interface_down (LinearAlgebra::distributed::Vector< double > &dst, const LinearAlgebra::distributed::Vector< double > &src) const |
void | vmult_interface_up (LinearAlgebra::distributed::Vector< double > &dst, const LinearAlgebra::distributed::Vector< double > &src) const |
void | vmult (LinearAlgebra::distributed::Vector< double > &dst, const LinearAlgebra::distributed::Vector< double > &src) const |
void | Tvmult (LinearAlgebra::distributed::Vector< double > &dst, const LinearAlgebra::distributed::Vector< double > &src) const |
void | vmult_add (LinearAlgebra::distributed::Vector< double > &dst, const LinearAlgebra::distributed::Vector< double > &src) const |
void | Tvmult_add (LinearAlgebra::distributed::Vector< double > &dst, const LinearAlgebra::distributed::Vector< double > &src) const |
value_type | el (const unsigned int row, const unsigned int col) const |
virtual std::size_t | memory_consumption () const |
void | initialize_dof_vector (LinearAlgebra::distributed::Vector< double > &vec) const |
std::shared_ptr< const MatrixFree< dim, value_type, VectorizedArray< typename LinearAlgebra::distributed::Vector< double > ::value_type > > > | get_matrix_free () const |
const std::shared_ptr< DiagonalMatrix< LinearAlgebra::distributed::Vector< double > > > & | get_matrix_diagonal_inverse () const |
const std::shared_ptr< DiagonalMatrix< LinearAlgebra::distributed::Vector< double > > > & | get_matrix_diagonal () const |
void | precondition_Jacobi (LinearAlgebra::distributed::Vector< double > &dst, const LinearAlgebra::distributed::Vector< double > &src, const value_type omega) const |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
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) |
Private Member Functions | |
virtual void | apply_add (VectorType &dst, const VectorType &src) const override |
void | local_apply_cell (const MatrixFree< dim, value_type, VectorizedArrayType > &data, VectorType &dst, const VectorType &src, const std::pair< unsigned int, unsigned int > &cell_range) const |
void | local_diagonal_cell (const MatrixFree< dim, value_type, VectorizedArrayType > &data, VectorType &dst, const VectorType &, const std::pair< unsigned int, unsigned int > &cell_range) const |
void | do_operation_on_cell (FEEvaluation< dim, fe_degree, n_q_points_1d, n_components, value_type > &phi, const unsigned int cell) const |
Private Attributes | |
std::shared_ptr< Table< 2, VectorizedArrayType > > | scalar_coefficient |
This class implements the operation of the action of a Laplace matrix, namely \( L_{ij} = \int_\Omega c(\mathbf x) \mathbf \nabla N_i(\mathbf x) \cdot \mathbf \nabla N_j(\mathbf x)\,d \mathbf x\), where \(c(\mathbf x)\) is the scalar heterogeneity coefficient.
Note that this class only supports the non-blocked vector variant of the Base operator because only a single FEEvaluation object is used in the apply function.
Definition at line 815 of file operators.h.
using MatrixFreeOperators::LaplaceOperator< dim, fe_degree, n_q_points_1d, n_components, VectorType, VectorizedArrayType >::value_type = typename Base<dim, VectorType, VectorizedArrayType>::value_type |
Number alias.
Definition at line 822 of file operators.h.
using MatrixFreeOperators::LaplaceOperator< dim, fe_degree, n_q_points_1d, n_components, VectorType, VectorizedArrayType >::size_type = typename Base<dim, VectorType, VectorizedArrayType>::size_type |
size_type needed for preconditioner classes.
Definition at line 828 of file operators.h.
MatrixFreeOperators::LaplaceOperator< dim, fe_degree, n_q_points_1d, n_components, VectorType, VectorizedArrayType >::LaplaceOperator |
Constructor.
Definition at line 1903 of file operators.h.
|
overridevirtual |
The diagonal is approximated by computing a local diagonal matrix per element and distributing it to the global diagonal. This will lead to wrong results on element with hanging nodes but is still an acceptable approximation to be used in preconditioners.
Definition at line 1982 of file operators.h.
void MatrixFreeOperators::LaplaceOperator< dim, fe_degree, n_q_points_1d, n_components, VectorType, VectorizedArrayType >::set_coefficient | ( | const std::shared_ptr< Table< 2, VectorizedArrayType >> & | scalar_coefficient | ) |
Set the heterogeneous scalar coefficient scalar_coefficient
to be used at the quadrature points. The Table needs to have as many rows as there are cell batches in the underlying MatrixFree object, MatrixFree::n_cell_batches(). The number of batches is related to the fact that the matrix-free operators do not work on individual cells, but instead of batches of cells at once due to vectorization. The Table can take two different numbers of columns. One case is to select it equal to the total number of quadrature points in dim
dimensions, which is the dim
th power of the n_q_points_1d
template parameter. Here, (*scalar_coefficient)(cell,q)
corresponds to the value of the coefficient on cell batch cell
and quadrature point index q
. The second supported variant is a Table with a single column, in which case the same variable coefficient value is used at all quadrature points of a cell.
Such tables can be initialized by
where mf_data
is a MatrixFree object and function
is a function which provides the following method VectorizedArray<double> value(const Point<dim, VectorizedArray<double> > &p_vec)
.
If this function is not called, the coefficient is assumed to be unity.
The argument to this function is a shared pointer to such a table. The class stores the shared pointer to this table, not a deep copy and uses it to form the Laplace matrix. Consequently, you can update the table and re-use the current object to obtain the action of a Laplace matrix with this updated coefficient. Alternatively, if the table values are only to be filled once, the original shared pointer can also go out of scope in user code and the clear() command or destructor of this class will delete the table.
Definition at line 1942 of file operators.h.
|
overridevirtual |
Resets all data structures back to the same state as for a newly constructed object.
Reimplemented from MatrixFreeOperators::Base< dim, LinearAlgebra::distributed::Vector< double >, VectorizedArray< typename LinearAlgebra::distributed::Vector< double > ::value_type > >.
Definition at line 1921 of file operators.h.
std::shared_ptr< Table< 2, VectorizedArrayType > > MatrixFreeOperators::LaplaceOperator< dim, fe_degree, n_q_points_1d, n_components, VectorType, VectorizedArrayType >::get_coefficient |
Read/Write access to coefficients to be used in Laplace operator.
The function will throw an error if coefficients are not previously set by set_coefficient() function.
Definition at line 1962 of file operators.h.
|
overrideprivatevirtual |
Applies the laplace matrix operation on an input vector. It is assumed that the passed input and output vector are correctly initialized using initialize_dof_vector().
Definition at line 2032 of file operators.h.
|
private |
Applies the Laplace operator on a cell.
Definition at line 2132 of file operators.h.
|
private |
Apply diagonal part of the Laplace operator on a cell.
Definition at line 2168 of file operators.h.
|
private |
Apply Laplace operator on a cell cell
.
Definition at line 2068 of file operators.h.
|
private |
User-provided heterogeneity coefficient.
Definition at line 955 of file operators.h.