Reference documentation for deal.II version 9.1.1
|
Local integrators related to advection along a vector field and its DG formulations. More...
Functions | |
template<int dim> | |
void | cell_matrix (FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, const ArrayView< const std::vector< double >> &velocity, const double factor=1.) |
template<int dim> | |
void | cell_residual (Vector< double > &result, const FEValuesBase< dim > &fe, const std::vector< Tensor< 1, dim >> &input, const ArrayView< const std::vector< double >> &velocity, double factor=1.) |
template<int dim> | |
void | cell_residual (Vector< double > &result, const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim >>> &input, const ArrayView< const std::vector< double >> &velocity, double factor=1.) |
template<int dim> | |
void | cell_residual (Vector< double > &result, const FEValuesBase< dim > &fe, const std::vector< double > &input, const ArrayView< const std::vector< double >> &velocity, double factor=1.) |
template<int dim> | |
void | cell_residual (Vector< double > &result, const FEValuesBase< dim > &fe, const ArrayView< const std::vector< double >> &input, const ArrayView< const std::vector< double >> &velocity, double factor=1.) |
template<int dim> | |
void | upwind_value_matrix (FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, const ArrayView< const std::vector< double >> &velocity, double factor=1.) |
template<int dim> | |
void | upwind_value_residual (Vector< double > &result, const FEValuesBase< dim > &fe, const std::vector< double > &input, const std::vector< double > &data, const ArrayView< const std::vector< double >> &velocity, double factor=1.) |
template<int dim> | |
void | upwind_value_residual (Vector< double > &result, const FEValuesBase< dim > &fe, const ArrayView< const std::vector< double >> &input, const ArrayView< const std::vector< double >> &data, const ArrayView< const std::vector< double >> &velocity, double factor=1.) |
template<int dim> | |
void | upwind_value_matrix (FullMatrix< double > &M11, FullMatrix< double > &M12, FullMatrix< double > &M21, FullMatrix< double > &M22, const FEValuesBase< dim > &fe1, const FEValuesBase< dim > &fe2, const FEValuesBase< dim > &fetest1, const FEValuesBase< dim > &fetest2, const ArrayView< const std::vector< double >> &velocity, const double factor=1.) |
template<int dim> | |
void | upwind_face_residual (Vector< double > &result1, Vector< double > &result2, const FEValuesBase< dim > &fe1, const FEValuesBase< dim > &fe2, const std::vector< double > &input1, const std::vector< double > &input2, const ArrayView< const std::vector< double >> &velocity, const double factor=1.) |
template<int dim> | |
void | upwind_face_residual (Vector< double > &result1, Vector< double > &result2, const FEValuesBase< dim > &fe1, const FEValuesBase< dim > &fe2, const ArrayView< const std::vector< double >> &input1, const ArrayView< const std::vector< double >> &input2, const ArrayView< const std::vector< double >> &velocity, const double factor=1.) |
Local integrators related to advection along a vector field and its DG formulations.
All advection operators depend on an advection velocity denoted by w in the formulas below. It is denoted as velocity
in the parameter lists.
The functions cell_matrix() and both upwind_value_matrix() are taking the equation in weak form, that is, the directional derivative is on the test function.
void LocalIntegrators::Advection::cell_matrix | ( | FullMatrix< double > & | M, |
const FEValuesBase< dim > & | fe, | ||
const FEValuesBase< dim > & | fetest, | ||
const ArrayView< const std::vector< double >> & | velocity, | ||
const double | factor = 1. |
||
) |
Advection along the direction w in weak form with derivative on the test function
\[ m_{ij} = \int_Z u_j\,(\mathbf w \cdot \nabla) v_i \, dx. \]
The FiniteElement in fe
may be scalar or vector valued. In the latter case, the advection operator is applied to each component separately.
M | The advection matrix obtained as result |
fe | The FEValues object describing the local trial function space. update_values and update_gradients, and update_JxW_values must be set. |
fetest | The FEValues object describing the local test function space. update_values and update_gradients must be set. |
velocity | The advection velocity, a vector of dimension dim . Each component may either contain a vector of length one, in which case a constant velocity is assumed, or a vector with as many entries as quadrature points if the velocity is not constant. |
factor | is an optional multiplication factor for the result. |
Definition at line 80 of file advection.h.
|
inline |
Scalar advection residual operator in strong form
\[ r_i = \int_Z (\mathbf w \cdot \nabla)u\, v_i \, dx. \]
Definition at line 136 of file advection.h.
|
inline |
Vector-valued advection residual operator in strong form
\[ r_i = \int_Z \bigl((\mathbf w \cdot \nabla) \mathbf u\bigr) \cdot\mathbf v_i \, dx. \]
Definition at line 179 of file advection.h.
|
inline |
Scalar advection residual operator in weak form
\[ r_i = \int_Z (\mathbf w \cdot \nabla)v\, u_i \, dx. \]
Definition at line 221 of file advection.h.
|
inline |
Vector-valued advection residual operator in weak form
\[ r_i = \int_Z \bigl((\mathbf w \cdot \nabla) \mathbf v\bigr) \cdot\mathbf u_i \, dx. \]
Definition at line 261 of file advection.h.
void LocalIntegrators::Advection::upwind_value_matrix | ( | FullMatrix< double > & | M, |
const FEValuesBase< dim > & | fe, | ||
const FEValuesBase< dim > & | fetest, | ||
const ArrayView< const std::vector< double >> & | velocity, | ||
double | factor = 1. |
||
) |
Upwind flux at the boundary for weak advection operator. This is the value of the trial function at the outflow boundary and zero else:
\[ a_{ij} = \int_{\partial\Omega} [\mathbf w\cdot\mathbf n]_+ u_i v_j \, ds \]
The velocity
is provided as an ArrayView, having dim
vectors, one for each velocity component. Each of the vectors must either have only a single entry, if the advection velocity is constant, or have an entry for each quadrature point.
The finite element can have several components, in which case each component is advected by the same velocity.
Definition at line 315 of file advection.h.
|
inline |
Scalar case: Residual for upwind flux at the boundary for weak advection operator. This is the value of the trial function at the outflow boundary and the value of the incoming boundary condition on the inflow boundary:
\[ a_{ij} = \int_{\partial\Omega} (\mathbf w\cdot\mathbf n) \widehat u v_j \, ds \]
Here, the numerical flux \(\widehat u\) is the upwind value at the face, namely the finite element function whose values are given in the argument input
on the outflow boundary. On the inflow boundary, it is the inhomogeneous boundary value in the argument data
.
The velocity
is provided as an ArrayView, having dim
vectors, one for each velocity component. Each of the vectors must either have only a single entry, if the advection velocity is constant, or have an entry for each quadrature point.
The finite element can have several components, in which case each component is advected by the same velocity.
Definition at line 388 of file advection.h.
|
inline |
Vector-valued case: Residual for upwind flux at the boundary for weak advection operator. This is the value of the trial function at the outflow boundary and the value of the incoming boundary condition on the inflow boundary:
\[ a_{ij} = \int_{\partial\Omega} (\mathbf w\cdot\mathbf n) \widehat u v_j \, ds \]
Here, the numerical flux \(\widehat u\) is the upwind value at the face, namely the finite element function whose values are given in the argument input
on the outflow boundary. On the inflow boundary, it is the inhomogeneous boundary value in the argument data
.
The velocity
is provided as an ArrayView, having dim
vectors, one for each velocity component. Each of the vectors must either have only a single entry, if the advection velocity is constant, or have an entry for each quadrature point.
The finite element can have several components, in which case each component is advected by the same velocity.
Definition at line 455 of file advection.h.
void LocalIntegrators::Advection::upwind_value_matrix | ( | FullMatrix< double > & | M11, |
FullMatrix< double > & | M12, | ||
FullMatrix< double > & | M21, | ||
FullMatrix< double > & | M22, | ||
const FEValuesBase< dim > & | fe1, | ||
const FEValuesBase< dim > & | fe2, | ||
const FEValuesBase< dim > & | fetest1, | ||
const FEValuesBase< dim > & | fetest2, | ||
const ArrayView< const std::vector< double >> & | velocity, | ||
const double | factor = 1. |
||
) |
Upwind flux in the interior for weak advection operator. Matrix entries correspond to the upwind value of the trial function, multiplied by the jump of the test functions
\[ a_{ij} = \int_F \left|\mathbf w \cdot \mathbf n\right| u^\uparrow (v^\uparrow-v^\downarrow) \,ds \]
The velocity
is provided as an ArrayView, having dim
vectors, one for each velocity component. Each of the vectors must either have only a single entry, if the advection velocity is constant, or have an entry for each quadrature point.
The finite element can have several components, in which case each component is advected the same way.
Definition at line 522 of file advection.h.
void LocalIntegrators::Advection::upwind_face_residual | ( | Vector< double > & | result1, |
Vector< double > & | result2, | ||
const FEValuesBase< dim > & | fe1, | ||
const FEValuesBase< dim > & | fe2, | ||
const std::vector< double > & | input1, | ||
const std::vector< double > & | input2, | ||
const ArrayView< const std::vector< double >> & | velocity, | ||
const double | factor = 1. |
||
) |
Scalar case: Upwind flux in the interior for weak advection operator. Matrix entries correspond to the upwind value of the trial function, multiplied by the jump of the test functions
\[ a_{ij} = \int_F \left|\mathbf w \cdot \mathbf n\right| u^\uparrow (v^\uparrow-v^\downarrow) \,ds \]
The velocity
is provided as an ArrayView, having dim
vectors, one for each velocity component. Each of the vectors must either have only a single entry, if the advection velocity is constant, or have an entry for each quadrature point.
The finite element can have several components, in which case each component is advected the same way.
Definition at line 607 of file advection.h.
void LocalIntegrators::Advection::upwind_face_residual | ( | Vector< double > & | result1, |
Vector< double > & | result2, | ||
const FEValuesBase< dim > & | fe1, | ||
const FEValuesBase< dim > & | fe2, | ||
const ArrayView< const std::vector< double >> & | input1, | ||
const ArrayView< const std::vector< double >> & | input2, | ||
const ArrayView< const std::vector< double >> & | velocity, | ||
const double | factor = 1. |
||
) |
Vector-valued case: Upwind flux in the interior for weak advection operator. Matrix entries correspond to the upwind value of the trial function, multiplied by the jump of the test functions
\[ a_{ij} = \int_F \left|\mathbf w \cdot \mathbf n\right| u^\uparrow (v^\uparrow-v^\downarrow) \,ds \]
The velocity
is provided as an ArrayView, having dim
vectors, one for each velocity component. Each of the vectors must either have only a single entry, if the advection velocity is constant, or have an entry for each quadrature point.
The finite element can have several components, in which case each component is advected the same way.
Definition at line 684 of file advection.h.