Reference documentation for deal.II version 9.0.0
|
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 VectorSlice< const std::vector< 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.) |
Local integrators related to L2-inner products.
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 \]
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 \]
The size of the vector weights
must be equal to the number of quadrature points in the finite element.
void LocalIntegrators::L2::L2 | ( | Vector< number > & | result, |
const FEValuesBase< dim > & | fe, | ||
const std::vector< double > & | input, | ||
const double | factor = 1. |
||
) |
void LocalIntegrators::L2::L2 | ( | Vector< number > & | result, |
const FEValuesBase< dim > & | fe, | ||
const VectorSlice< const std::vector< std::vector< double > > > & | input, | ||
const double | factor = 1. |
||
) |
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.