Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
Public Member Functions | Private Attributes | List of all members
MatrixFreeOperators::CellwiseInverseMassMatrix< dim, fe_degree, n_components, Number > Class Template Reference

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

Public Member Functions

 CellwiseInverseMassMatrix (const FEEvaluationBase< dim, n_components, Number > &fe_eval)
 
void apply (const AlignedVector< VectorizedArray< Number >> &inverse_coefficient, const unsigned int n_actual_components, const VectorizedArray< Number > *in_array, VectorizedArray< Number > *out_array) const
 
void transform_from_q_points_to_basis (const unsigned int n_actual_components, const VectorizedArray< Number > *in_array, VectorizedArray< Number > *out_array) const
 
void fill_inverse_JxW_values (AlignedVector< VectorizedArray< Number >> &inverse_jxw) const
 

Private Attributes

const FEEvaluationBase< dim, n_components, Number > & fe_eval
 
AlignedVector< VectorizedArray< Number > > inverse_shape
 

Detailed Description

template<int dim, int fe_degree, int n_components = 1, typename Number = double>
class MatrixFreeOperators::CellwiseInverseMassMatrix< dim, fe_degree, n_components, Number >

This class implements the operation of the action of the inverse of a mass matrix on an element for the special case of an evaluation object with as many quadrature points as there are cell degrees of freedom. It uses algorithms from FEEvaluation and produces the exact mass matrix for DGQ elements. This algorithm uses tensor products of inverse 1D shape matrices over quadrature points, so the inverse operation is exactly as expensive as applying the forward operator on each cell. Of course, for continuous finite elements this operation does not produce the inverse of a mass operation as the coupling between the elements cannot be considered by this operation.

The equation may contain variable coefficients, so the user is required to provide an array for the inverse of the local coefficient (this class provide a helper method 'fill_inverse_JxW_values' to get the inverse of a constant-coefficient operator).

Author
Martin Kronbichler, 2014

Definition at line 615 of file operators.h.

Constructor & Destructor Documentation

◆ CellwiseInverseMassMatrix()

template<int dim, int fe_degree, int n_components, typename Number >
MatrixFreeOperators::CellwiseInverseMassMatrix< dim, fe_degree, n_components, Number >::CellwiseInverseMassMatrix ( const FEEvaluationBase< dim, n_components, Number > &  fe_eval)
inline

Constructor. Initializes the shape information from the ShapeInfo field in the class FEEval.

Definition at line 917 of file operators.h.

Member Function Documentation

◆ apply()

template<int dim, int fe_degree, int n_components, typename Number >
void MatrixFreeOperators::CellwiseInverseMassMatrix< dim, fe_degree, n_components, Number >::apply ( const AlignedVector< VectorizedArray< Number >> &  inverse_coefficient,
const unsigned int  n_actual_components,
const VectorizedArray< Number > *  in_array,
VectorizedArray< Number > *  out_array 
) const
inline

Applies the inverse mass matrix operation on an input array. It is assumed that the passed input and output arrays are of correct size, namely FEEval::dofs_per_cell * n_components long. The inverse of the local coefficient (also containing the inverse JxW values) must be passed as first argument. Passing more than one component in the coefficient is allowed.

Definition at line 970 of file operators.h.

◆ transform_from_q_points_to_basis()

template<int dim, int fe_degree, int n_components, typename Number >
void MatrixFreeOperators::CellwiseInverseMassMatrix< dim, fe_degree, n_components, Number >::transform_from_q_points_to_basis ( const unsigned int  n_actual_components,
const VectorizedArray< Number > *  in_array,
VectorizedArray< Number > *  out_array 
) const
inline

This operation performs a projection from the data given in quadrature points to the actual basis underlying this object. This projection can also be interpreted as a change of the basis from the Lagrange interpolation polynomials in the quadrature points to the basis underlying the current fe_eval object.

Calling this function on an array as

inverse_mass.transform_from_q_points_to_basis(1, array,
phi.begin_dof_values());

is equivalent to

for (unsigned int q=0; q<phi.n_q_points; ++q)
phi.submit_value(array[q], q);
phi.integrate(true, false);
inverse_mass.apply(coefficients, 1, phi.begin_dof_values(),
phi.begin_dof_values());

provided that coefficients holds the inverse of the quadrature weights and no additional coefficients. This setup highlights the underlying projection, testing a right hand side and applying an inverse mass matrix. This function works both for the scalar case as described in the example or for multiple components that are laid out component by component.

Compared to the more verbose alternative, the given procedure is considerably faster because it can bypass the integrate() step and the first half of the transformation to the quadrature points, reducing the number of tensor product calls from 3*dim*n_components to dim*n_components.

Definition at line 1042 of file operators.h.

◆ fill_inverse_JxW_values()

template<int dim, int fe_degree, int n_components, typename Number >
void MatrixFreeOperators::CellwiseInverseMassMatrix< dim, fe_degree, n_components, Number >::fill_inverse_JxW_values ( AlignedVector< VectorizedArray< Number >> &  inverse_jxw) const
inline

Fills the given array with the inverse of the JxW values, i.e., a mass matrix with coefficient 1. Non-unit coefficients must be multiplied (in inverse form) to this array.

Definition at line 946 of file operators.h.

Member Data Documentation

◆ fe_eval

template<int dim, int fe_degree, int n_components = 1, typename Number = double>
const FEEvaluationBase<dim, n_components, Number>& MatrixFreeOperators::CellwiseInverseMassMatrix< dim, fe_degree, n_components, Number >::fe_eval
private

A reference to the FEEvaluation object for getting the JxW_values.

Definition at line 690 of file operators.h.

◆ inverse_shape

template<int dim, int fe_degree, int n_components = 1, typename Number = double>
AlignedVector<VectorizedArray<Number> > MatrixFreeOperators::CellwiseInverseMassMatrix< dim, fe_degree, n_components, Number >::inverse_shape
private

A structure to hold inverse shape functions

Definition at line 695 of file operators.h.


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