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
Functions
LocalIntegrators::L2 Namespace Reference

Local integrators related to L2-inner products. More...

Functions

template<int dim>
void mass_matrix (FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double factor=1.)
 
template<int dim>
void weighted_mass_matrix (FullMatrix< double > &M, const FEValuesBase< dim > &fe, const std::vector< double > &weights)
 
template<int dim, typename number >
void L2 (Vector< number > &result, const FEValuesBase< dim > &fe, const std::vector< double > &input, const double factor=1.)
 
template<int dim, typename number >
void L2 (Vector< number > &result, const FEValuesBase< dim > &fe, const ArrayView< const std::vector< double > > &input, const double factor=1.)
 
template<int dim>
void jump_matrix (FullMatrix< double > &M11, FullMatrix< double > &M12, FullMatrix< double > &M21, FullMatrix< double > &M22, const FEValuesBase< dim > &fe1, const FEValuesBase< dim > &fe2, const double factor1=1., const double factor2=1.)
 

Detailed Description

Local integrators related to L2-inner products.

Function Documentation

◆ mass_matrix()

template<int dim>
void LocalIntegrators::L2::mass_matrix ( FullMatrix< double > &  M,
const FEValuesBase< dim > &  fe,
const double  factor = 1. 
)

The mass matrix for scalar or vector values finite elements.

\[ \int_Z uv\,dx \quad \text{or} \quad \int_Z \mathbf u\cdot \mathbf v\,dx \]

Likewise, this term can be used on faces, where it computes the integrals

\[ \int_F uv\,ds \quad \text{or} \quad \int_F \mathbf u\cdot \mathbf v\,ds \]

Parameters
MThe mass matrix obtained as result.
feThe FEValues object describing the local trial function space. update_values and update_JxW_values must be set.
factorA constant that multiplies the mass matrix.

Definition at line 58 of file l2.h.

◆ weighted_mass_matrix()

template<int dim>
void LocalIntegrators::L2::weighted_mass_matrix ( FullMatrix< double > &  M,
const FEValuesBase< dim > &  fe,
const std::vector< double > &  weights 
)

The weighted mass matrix for scalar or vector values finite elements.

\[ \int_Z \omega(x) uv\,dx \quad \text{or} \quad \int_Z \omega(x) \mathbf u\cdot \mathbf v\,dx \]

Likewise, this term can be used on faces, where it computes the integrals

\[ \int_F \omega(x) uv\,ds \quad \text{or} \quad \int_F \omega(x) \mathbf u\cdot \mathbf v\,ds \]

Parameters
MThe weighted mass matrix obtained as result.
feThe FEValues object describing the local trial function space. update_values and update_JxW_values must be set.
weightsThe weights, \(\omega(x)\), evaluated at the quadrature points in the finite element (size must be equal to the number of quadrature points in the element).

Definition at line 109 of file l2.h.

◆ L2() [1/2]

template<int dim, typename number >
void LocalIntegrators::L2::L2 ( Vector< number > &  result,
const FEValuesBase< dim > &  fe,
const std::vector< double > &  input,
const double  factor = 1. 
)

L2-inner product for scalar functions.

\[ \int_Z fv\,dx \quad \text{or} \quad \int_F fv\,ds \]

Parameters
resultThe vector obtained as result.
feThe FEValues object describing the local trial function space. update_values and update_JxW_values must be set.
inputThe representation of \(f\) evaluated at the quadrature points in the finite element (size must be equal to the number of quadrature points in the element).
factorA constant that multiplies the result.

Definition at line 160 of file l2.h.

◆ L2() [2/2]

template<int dim, typename number >
void LocalIntegrators::L2::L2 ( Vector< number > &  result,
const FEValuesBase< dim > &  fe,
const ArrayView< const std::vector< double > > &  input,
const double  factor = 1. 
)

L2-inner product for a slice of a vector valued right hand side.

\[ \int_Z \mathbf f\cdot \mathbf v\,dx \quad \text{or} \quad \int_F \mathbf f\cdot \mathbf v\,ds \]

Parameters
resultThe vector obtained as result.
feThe FEValues object describing the local trial function space. update_values and update_JxW_values must be set.
inputThe vector valued representation of \(\mathbf f\) evaluated at the quadrature points in the finite element (size of each component must be equal to the number of quadrature points in the element).
factorA constant that multiplies the result.

Definition at line 190 of file l2.h.

◆ jump_matrix()

template<int dim>
void LocalIntegrators::L2::jump_matrix ( FullMatrix< double > &  M11,
FullMatrix< double > &  M12,
FullMatrix< double > &  M21,
FullMatrix< double > &  M22,
const FEValuesBase< dim > &  fe1,
const FEValuesBase< dim > &  fe2,
const double  factor1 = 1.,
const double  factor2 = 1. 
)

The jump matrix between two cells for scalar or vector values finite elements. Note that the factor \(\gamma\) can be used to implement weighted jumps.

\[ \int_F [\gamma u][\gamma v]\,ds \quad \text{or} \int_F [\gamma \mathbf u]\cdot [\gamma \mathbf v]\,ds \]

Using appropriate weights, this term can be used to penalize violation of conformity in H1.

Note that for the parameters that follow, the external matrix refers to the flux between cells, while the internal matrix refers to entries coupling inside the cell.

Parameters
M11The internal matrix for the first cell obtained as result.
M12The external matrix for the first cell obtained as result.
M21The external matrix for the second cell obtained as result.
M22The internal matrix for the second cell obtained as result.
fe1The FEValues object describing the local trial function space for the first cell. update_values and update_JxW_values must be set.
fe2The FEValues object describing the local trial function space for the second cell. update_values and update_JxW_values must be set.
factor1A constant that multiplies the shape functions for the first cell.
factor2A constant that multiplies the shape functions for the second cell.

Definition at line 238 of file l2.h.