Reference documentation for deal.II version GIT relicensing-249-g48dc7357c7 2024-03-29 12:30:02+00:00
\(\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::Advection Namespace Reference

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.)
 

Detailed Description

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.

Function Documentation

◆ cell_matrix()

template<int dim>
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.

Parameters
MThe advection matrix obtained as result
feThe FEValues object describing the local trial function space. update_values and update_gradients, and update_JxW_values must be set.
fetestThe FEValues object describing the local test function space. update_values and update_gradients must be set.
velocityThe 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.
factoris an optional multiplication factor for the result.

Definition at line 74 of file advection.h.

◆ cell_residual() [1/4]

template<int dim>
void LocalIntegrators::Advection::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. 
)
inline

Scalar advection residual operator in strong form

\[ r_i = \int_Z (\mathbf w \cdot \nabla)u\, v_i \, dx. \]

Warning
This is not the residual consistent with cell_matrix(), but with its transpose.

Definition at line 130 of file advection.h.

◆ cell_residual() [2/4]

template<int dim>
void LocalIntegrators::Advection::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. 
)
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. \]

Warning
This is not the residual consistent with cell_matrix(), but with its transpose.

Definition at line 173 of file advection.h.

◆ cell_residual() [3/4]

template<int dim>
void LocalIntegrators::Advection::cell_residual ( Vector< double > &  result,
const FEValuesBase< dim > &  fe,
const std::vector< double > &  input,
const ArrayView< const std::vector< double > > &  velocity,
double  factor = 1. 
)
inline

Scalar advection residual operator in weak form

\[ r_i = \int_Z (\mathbf w \cdot \nabla)v\, u_i \, dx. \]

Definition at line 215 of file advection.h.

◆ cell_residual() [4/4]

template<int dim>
void LocalIntegrators::Advection::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. 
)
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 255 of file advection.h.

◆ upwind_value_matrix() [1/2]

template<int dim>
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 309 of file advection.h.

◆ upwind_value_residual() [1/2]

template<int dim>
void LocalIntegrators::Advection::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. 
)
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 382 of file advection.h.

◆ upwind_value_residual() [2/2]

template<int dim>
void LocalIntegrators::Advection::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. 
)
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 449 of file advection.h.

◆ upwind_value_matrix() [2/2]

template<int dim>
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 516 of file advection.h.

◆ upwind_face_residual() [1/2]

template<int dim>
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 601 of file advection.h.

◆ upwind_face_residual() [2/2]

template<int dim>
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 678 of file advection.h.