Reference documentation for deal.II version 9.4.1
|
#include <deal.II/matrix_free/operators.h>
Public Member Functions | |
CellwiseInverseMassMatrix (const FEEvaluationBase< dim, n_components, Number, false, VectorizedArrayType > &fe_eval) | |
void | apply (const AlignedVector< VectorizedArrayType > &inverse_coefficient, const unsigned int n_actual_components, const VectorizedArrayType *in_array, VectorizedArrayType *out_array) const |
void | apply (const VectorizedArrayType *in_array, VectorizedArrayType *out_array) const |
void | transform_from_q_points_to_basis (const unsigned int n_actual_components, const VectorizedArrayType *in_array, VectorizedArrayType *out_array) const |
void | fill_inverse_JxW_values (AlignedVector< VectorizedArrayType > &inverse_jxw) const |
Private Attributes | |
const FEEvaluationBase< dim, n_components, Number, false, VectorizedArrayType > & | fe_eval |
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).
Definition at line 623 of file operators.h.
|
inline |
Constructor. Initializes the shape information from the ShapeInfo field in the class FEEval.
Definition at line 1007 of file operators.h.
|
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 FEEvaluation::dofs_per_cell 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 1085 of file operators.h.
|
inline |
Applies the inverse mass matrix operation on an input array, using the inverse of the JxW values provided by the fe_eval
argument passed to the constructor of this class. Note that the user code must call FEEvaluation::reinit() on the underlying evaluator to make the FEEvaluationBase::JxW() method return the information of the correct cell. It is assumed that the pointers of the input and output arrays are valid over the length FEEvaluation::dofs_per_cell, which is the number of entries processed by this function. The in_array
and out_array
arguments may point to the same memory position.
Definition at line 1062 of file operators.h.
|
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
is equivalent to
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 1124 of file operators.h.
|
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 1029 of file operators.h.
|
private |
A reference to the FEEvaluation object for getting the JxW_values.
Definition at line 725 of file operators.h.