Reference documentation for deal.II version 9.0.0
|
Local integrators related to the divergence operator and its trace. More...
Functions | |
template<int dim> | |
void | cell_matrix (FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, double factor=1.) |
template<int dim, typename number > | |
void | cell_residual (Vector< number > &result, const FEValuesBase< dim > &fetest, const VectorSlice< const std::vector< std::vector< Tensor< 1, dim > > > > &input, const double factor=1.) |
template<int dim, typename number > | |
void | cell_residual (Vector< number > &result, const FEValuesBase< dim > &fetest, const VectorSlice< const std::vector< std::vector< double > > > &input, const double factor=1.) |
template<int dim> | |
void | gradient_matrix (FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, double factor=1.) |
template<int dim, typename number > | |
void | gradient_residual (Vector< number > &result, const FEValuesBase< dim > &fetest, const std::vector< Tensor< 1, dim > > &input, const double factor=1.) |
template<int dim, typename number > | |
void | gradient_residual (Vector< number > &result, const FEValuesBase< dim > &fetest, const std::vector< double > &input, const double factor=1.) |
template<int dim> | |
void | u_dot_n_matrix (FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, double factor=1.) |
template<int dim, typename number > | |
void | u_dot_n_residual (Vector< number > &result, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, const VectorSlice< const std::vector< std::vector< double > > > &data, double factor=1.) |
template<int dim, typename number > | |
void | u_times_n_residual (Vector< number > &result, const FEValuesBase< dim > &fetest, const std::vector< double > &data, double factor=1.) |
template<int dim> | |
void | u_dot_n_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, double factor=1.) |
template<int dim> | |
void | grad_div_matrix (FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double factor=1.) |
template<int dim, typename number > | |
void | grad_div_residual (Vector< number > &result, const FEValuesBase< dim > &fetest, const VectorSlice< const std::vector< std::vector< Tensor< 1, dim > > > > &input, const double factor=1.) |
template<int dim> | |
void | u_dot_n_jump_matrix (FullMatrix< double > &M11, FullMatrix< double > &M12, FullMatrix< double > &M21, FullMatrix< double > &M22, const FEValuesBase< dim > &fe1, const FEValuesBase< dim > &fe2, double factor=1.) |
template<int dim> | |
double | norm (const FEValuesBase< dim > &fe, const VectorSlice< const std::vector< std::vector< Tensor< 1, dim > > > > &Du) |
Local integrators related to the divergence operator and its trace.
void LocalIntegrators::Divergence::cell_matrix | ( | FullMatrix< double > & | M, |
const FEValuesBase< dim > & | fe, | ||
const FEValuesBase< dim > & | fetest, | ||
double | factor = 1. |
||
) |
Cell matrix for divergence. The derivative is on the trial function.
\[ \int_Z v\nabla \cdot \mathbf u \,dx \]
This is the strong divergence operator and the trial space should be at least Hdiv. The test functions may be discontinuous.
Definition at line 53 of file divergence.h.
void LocalIntegrators::Divergence::cell_residual | ( | Vector< number > & | result, |
const FEValuesBase< dim > & | fetest, | ||
const VectorSlice< const std::vector< std::vector< Tensor< 1, dim > > > > & | input, | ||
const double | factor = 1. |
||
) |
The residual of the divergence operator in strong form.
\[ \int_Z v\nabla \cdot \mathbf u \,dx \]
This is the strong divergence operator and the trial space should be at least Hdiv. The test functions may be discontinuous.
The function cell_matrix() is the Frechet derivative of this function with respect to the test functions.
Definition at line 95 of file divergence.h.
void LocalIntegrators::Divergence::cell_residual | ( | Vector< number > & | result, |
const FEValuesBase< dim > & | fetest, | ||
const VectorSlice< const std::vector< std::vector< double > > > & | input, | ||
const double | factor = 1. |
||
) |
The residual of the divergence operator in weak form.
\[ - \int_Z \nabla v \cdot \mathbf u \,dx \]
This is the weak divergence operator and the test space should be at least H1. The trial functions may be discontinuous.
Definition at line 130 of file divergence.h.
void LocalIntegrators::Divergence::gradient_matrix | ( | FullMatrix< double > & | M, |
const FEValuesBase< dim > & | fe, | ||
const FEValuesBase< dim > & | fetest, | ||
double | factor = 1. |
||
) |
Cell matrix for gradient. The derivative is on the trial function.
\[ \int_Z \nabla u \cdot \mathbf v\,dx \]
This is the strong gradient and the trial space should be at least in H1. The test functions can be discontinuous.
Definition at line 163 of file divergence.h.
void LocalIntegrators::Divergence::gradient_residual | ( | Vector< number > & | result, |
const FEValuesBase< dim > & | fetest, | ||
const std::vector< Tensor< 1, dim > > & | input, | ||
const double | factor = 1. |
||
) |
The residual of the gradient operator in strong form.
\[ \int_Z \mathbf v\cdot\nabla u \,dx \]
This is the strong gradient operator and the trial space should be at least H1. The test functions may be discontinuous.
The function gradient_matrix() is the Frechet derivative of this function with respect to the test functions.
Definition at line 206 of file divergence.h.
void LocalIntegrators::Divergence::gradient_residual | ( | Vector< number > & | result, |
const FEValuesBase< dim > & | fetest, | ||
const std::vector< double > & | input, | ||
const double | factor = 1. |
||
) |
The residual of the gradient operator in weak form.
\[ -\int_Z \nabla\cdot \mathbf v u \,dx \]
This is the weak gradient operator and the test space should be at least Hdiv. The trial functions may be discontinuous.
Definition at line 240 of file divergence.h.
void LocalIntegrators::Divergence::u_dot_n_matrix | ( | FullMatrix< double > & | M, |
const FEValuesBase< dim > & | fe, | ||
const FEValuesBase< dim > & | fetest, | ||
double | factor = 1. |
||
) |
The trace of the divergence operator, namely the product of the normal component of the vector valued trial space and the test space.
\[ \int_F (\mathbf u\cdot \mathbf n) v \,ds \]
Definition at line 271 of file divergence.h.
void LocalIntegrators::Divergence::u_dot_n_residual | ( | Vector< number > & | result, |
const FEValuesBase< dim > & | fe, | ||
const FEValuesBase< dim > & | fetest, | ||
const VectorSlice< const std::vector< std::vector< double > > > & | data, | ||
double | factor = 1. |
||
) |
The trace of the divergence operator, namely the product of the normal component of the vector valued trial space and the test space.
\[ \int_F (\mathbf u\cdot \mathbf n) v \,ds \]
Definition at line 308 of file divergence.h.
void LocalIntegrators::Divergence::u_times_n_residual | ( | Vector< number > & | result, |
const FEValuesBase< dim > & | fetest, | ||
const std::vector< double > & | data, | ||
double | factor = 1. |
||
) |
The trace of the gradient operator, namely the product of the normal component of the vector valued test space and the trial space.
\[ \int_F u (\mathbf v\cdot \mathbf n) \,ds \]
Definition at line 344 of file divergence.h.
void LocalIntegrators::Divergence::u_dot_n_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, | ||
double | factor = 1. |
||
) |
The trace of the divergence operator, namely the product of the jump of the normal component of the vector valued trial function and the mean value of the test function.
\[ \int_F (\mathbf u_1\cdot \mathbf n_1 + \mathbf u_2 \cdot \mathbf n_2) \frac{v_1+v_2}{2} \,ds \]
Definition at line 379 of file divergence.h.
void LocalIntegrators::Divergence::grad_div_matrix | ( | FullMatrix< double > & | M, |
const FEValuesBase< dim > & | fe, | ||
const double | factor = 1. |
||
) |
Definition at line 437 of file divergence.h.
void LocalIntegrators::Divergence::grad_div_residual | ( | Vector< number > & | result, |
const FEValuesBase< dim > & | fetest, | ||
const VectorSlice< const std::vector< std::vector< Tensor< 1, dim > > > > & | input, | ||
const double | factor = 1. |
||
) |
Definition at line 457 of file divergence.h.
void LocalIntegrators::Divergence::u_dot_n_jump_matrix | ( | FullMatrix< double > & | M11, |
FullMatrix< double > & | M12, | ||
FullMatrix< double > & | M21, | ||
FullMatrix< double > & | M22, | ||
const FEValuesBase< dim > & | fe1, | ||
const FEValuesBase< dim > & | fe2, | ||
double | factor = 1. |
||
) |
The jump of the normal component
\[ \int_F (\mathbf u_1\cdot \mathbf n_1 + \mathbf u_2 \cdot \mathbf n_2) (\mathbf v_1\cdot \mathbf n_1 + \mathbf v_2 \cdot \mathbf n_2) \,ds \]
Definition at line 480 of file divergence.h.
double LocalIntegrators::Divergence::norm | ( | const FEValuesBase< dim > & | fe, |
const VectorSlice< const std::vector< std::vector< Tensor< 1, dim > > > > & | Du | ||
) |
The L2-norm of the divergence over the quadrature set determined by the FEValuesBase object.
The vector is expected to consist of dim vectors of length equal to the number of quadrature points. The number of components of the finite element has to be equal to the space dimension.
Definition at line 534 of file divergence.h.