Reference documentation for deal.II version 9.1.0-pre
Public Types | Public Member Functions | Private Member Functions | Private Attributes | List of all members
MatrixFreeOperators::LaplaceOperator< dim, fe_degree, n_q_points_1d, n_components, VectorType > Class Template Reference

#include <deal.II/matrix_free/operators.h>

Inheritance diagram for MatrixFreeOperators::LaplaceOperator< dim, fe_degree, n_q_points_1d, n_components, VectorType >:
[legend]

Public Types

using value_type = typename Base< dim, VectorType >::value_type
 
using size_type = typename Base< dim, VectorType >::size_type
 
- Public Types inherited from MatrixFreeOperators::Base< dim, VectorType >
using value_type = typename VectorType::value_type
 
using size_type = typename VectorType::size_type
 

Public Member Functions

 LaplaceOperator ()
 
virtual void compute_diagonal ()
 
void set_coefficient (const std::shared_ptr< Table< 2, VectorizedArray< value_type >>> &scalar_coefficient)
 
virtual void clear ()
 
std::shared_ptr< Table< 2, VectorizedArray< value_type > > > get_coefficient ()
 
- Public Member Functions inherited from MatrixFreeOperators::Base< dim, VectorType >
 Base ()
 
virtual ~Base () override=default
 
void initialize (std::shared_ptr< const MatrixFree< dim, 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 >> 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 >> 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 (VectorType &dst, const VectorType &src) const
 
void vmult_interface_up (VectorType &dst, const VectorType &src) const
 
void vmult (VectorType &dst, const VectorType &src) const
 
void Tvmult (VectorType &dst, const VectorType &src) const
 
void vmult_add (VectorType &dst, const VectorType &src) const
 
void Tvmult_add (VectorType &dst, const VectorType &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 (VectorType &vec) const
 
std::shared_ptr< const MatrixFree< dim, value_type > > get_matrix_free () const
 
const std::shared_ptr< DiagonalMatrix< VectorType > > & get_matrix_diagonal_inverse () const
 
const std::shared_ptr< DiagonalMatrix< VectorType > > & get_matrix_diagonal () const
 
void precondition_Jacobi (VectorType &dst, const VectorType &src, const value_type omega) const
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&) noexcept
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&) noexcept
 
void subscribe (const char *identifier=nullptr) const
 
void unsubscribe (const char *identifier=nullptr) 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
 
void local_apply_cell (const MatrixFree< dim, value_type > &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 > &data, VectorType &dst, const unsigned int &, 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, VectorizedArray< value_type > > > scalar_coefficient
 

Additional Inherited Members

- Static Public Member Functions inherited from Subscriptor
static::ExceptionBase & ExcInUse (int arg1, std::string arg2, std::string arg3)
 
static::ExceptionBase & ExcNoSubscriber (std::string arg1, std::string arg2)
 
- Protected Member Functions inherited from MatrixFreeOperators::Base< dim, VectorType >
void preprocess_constraints (VectorType &dst, const VectorType &src) const
 
void postprocess_constraints (VectorType &dst, const VectorType &src) const
 
void set_constrained_entries_to_one (VectorType &dst) const
 
virtual void Tapply_add (VectorType &dst, const VectorType &src) const
 
- Protected Attributes inherited from MatrixFreeOperators::Base< dim, VectorType >
std::shared_ptr< const MatrixFree< dim, value_type > > data
 
std::shared_ptr< DiagonalMatrix< VectorType > > diagonal_entries
 
std::shared_ptr< DiagonalMatrix< VectorType > > inverse_diagonal_entries
 
std::vector< unsigned int > selected_rows
 
std::vector< unsigned int > selected_columns
 

Detailed Description

template<int dim, int fe_degree, int n_q_points_1d = fe_degree + 1, int n_components = 1, typename VectorType = LinearAlgebra::distributed::Vector<double>>
class MatrixFreeOperators::LaplaceOperator< dim, fe_degree, n_q_points_1d, n_components, VectorType >

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.

Author
Denis Davydov, 2016

Definition at line 740 of file operators.h.

Member Typedef Documentation

template<int dim, int fe_degree, int n_q_points_1d = fe_degree + 1, int n_components = 1, typename VectorType = LinearAlgebra::distributed::Vector<double>>
using MatrixFreeOperators::LaplaceOperator< dim, fe_degree, n_q_points_1d, n_components, VectorType >::value_type = typename Base<dim, VectorType>::value_type

Number alias.

Definition at line 746 of file operators.h.

template<int dim, int fe_degree, int n_q_points_1d = fe_degree + 1, int n_components = 1, typename VectorType = LinearAlgebra::distributed::Vector<double>>
using MatrixFreeOperators::LaplaceOperator< dim, fe_degree, n_q_points_1d, n_components, VectorType >::size_type = typename Base<dim, VectorType>::size_type

size_type needed for preconditioner classes.

Definition at line 751 of file operators.h.

Constructor & Destructor Documentation

template<int dim, int fe_degree, int n_q_points_1d, int n_components, typename VectorType >
MatrixFreeOperators::LaplaceOperator< dim, fe_degree, n_q_points_1d, n_components, VectorType >::LaplaceOperator ( )

Constructor.

Definition at line 1746 of file operators.h.

Member Function Documentation

template<int dim, int fe_degree, int n_q_points_1d, int n_components, typename VectorType >
void MatrixFreeOperators::LaplaceOperator< dim, fe_degree, n_q_points_1d, n_components, VectorType >::compute_diagonal ( )
virtual

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.

Implements MatrixFreeOperators::Base< dim, VectorType >.

Definition at line 1812 of file operators.h.

template<int dim, int fe_degree, int n_q_points_1d = fe_degree + 1, int n_components = 1, typename VectorType = LinearAlgebra::distributed::Vector<double>>
void MatrixFreeOperators::LaplaceOperator< dim, fe_degree, n_q_points_1d, n_components, VectorType >::set_coefficient ( const std::shared_ptr< Table< 2, VectorizedArray< value_type >>> &  scalar_coefficient)

Set the heterogeneous scalar coefficient scalar_coefficient to be used at the quadrature points. The Table should be of correct size, consistent with the total number of quadrature points in dim-dimensions, controlled by the n_q_points_1d template parameter. Here, (*scalar_coefficient)(cell,q) corresponds to the value of the coefficient, where cell is an index into a set of cell batches as administered by the MatrixFree framework (which does not work on individual cells, but instead of batches of cells at once), and q is the number of the quadrature point within this batch.

Such tables can be initialized by

std::shared_ptr<Table<2, VectorizedArray<double> > > coefficient;
coefficient = std::make_shared<Table<2, VectorizedArray<double> > >();
{
const unsigned int n_cells = mf_data.n_macro_cells();
const unsigned int n_q_points = fe_eval.n_q_points;
coefficient->reinit(n_cells, n_q_points);
for (unsigned int cell=0; cell<n_cells; ++cell)
{
fe_eval.reinit(cell);
for (unsigned int q=0; q<n_q_points; ++q)
(*coefficient)(cell,q) =
function.value(fe_eval.quadrature_point(q));
}
}

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 1774 of file operators.h.

template<int dim, int fe_degree, int n_q_points_1d, int n_components, typename VectorType >
void MatrixFreeOperators::LaplaceOperator< dim, fe_degree, n_q_points_1d, n_components, VectorType >::clear ( )
virtual

Release all memory and return to a state just like after having called the default constructor.

Reimplemented from MatrixFreeOperators::Base< dim, VectorType >.

Definition at line 1759 of file operators.h.

template<int dim, int fe_degree, int n_q_points_1d, int n_components, typename VectorType >
std::shared_ptr< Table< 2, VectorizedArray< typename LaplaceOperator< dim, fe_degree, n_q_points_1d, n_components, VectorType >::value_type > > > MatrixFreeOperators::LaplaceOperator< dim, fe_degree, n_q_points_1d, n_components, VectorType >::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 1797 of file operators.h.

template<int dim, int fe_degree, int n_q_points_1d, int n_components, typename VectorType >
void MatrixFreeOperators::LaplaceOperator< dim, fe_degree, n_q_points_1d, n_components, VectorType >::apply_add ( VectorType &  dst,
const VectorType &  src 
) const
privatevirtual

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().

Implements MatrixFreeOperators::Base< dim, VectorType >.

Definition at line 1855 of file operators.h.

template<int dim, int fe_degree, int n_q_points_1d = fe_degree + 1, int n_components = 1, typename VectorType = LinearAlgebra::distributed::Vector<double>>
void MatrixFreeOperators::LaplaceOperator< dim, fe_degree, n_q_points_1d, n_components, VectorType >::local_apply_cell ( const MatrixFree< dim, value_type > &  data,
VectorType &  dst,
const VectorType &  src,
const std::pair< unsigned int, unsigned int > &  cell_range 
) const
private

Applies the Laplace operator on a cell.

Definition at line 1926 of file operators.h.

template<int dim, int fe_degree, int n_q_points_1d = fe_degree + 1, int n_components = 1, typename VectorType = LinearAlgebra::distributed::Vector<double>>
void MatrixFreeOperators::LaplaceOperator< dim, fe_degree, n_q_points_1d, n_components, VectorType >::local_diagonal_cell ( const MatrixFree< dim, value_type > &  data,
VectorType &  dst,
const unsigned int &  ,
const std::pair< unsigned int, unsigned int > &  cell_range 
) const
private

Apply diagonal part of the Laplace operator on a cell.

Definition at line 1952 of file operators.h.

template<int dim, int fe_degree, int n_q_points_1d = fe_degree + 1, int n_components = 1, typename VectorType = LinearAlgebra::distributed::Vector<double>>
void MatrixFreeOperators::LaplaceOperator< dim, fe_degree, n_q_points_1d, n_components, VectorType >::do_operation_on_cell ( FEEvaluation< dim, fe_degree, n_q_points_1d, n_components, value_type > &  phi,
const unsigned int  cell 
) const
private

Apply Laplace operator on a cell cell.

Definition at line 1887 of file operators.h.

Member Data Documentation

template<int dim, int fe_degree, int n_q_points_1d = fe_degree + 1, int n_components = 1, typename VectorType = LinearAlgebra::distributed::Vector<double>>
std::shared_ptr<Table<2, VectorizedArray<value_type> > > MatrixFreeOperators::LaplaceOperator< dim, fe_degree, n_q_points_1d, n_components, VectorType >::scalar_coefficient
private

User-provided heterogeneity coefficient.

Definition at line 870 of file operators.h.


The documentation for this class was generated from the following file: