OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Functions
AppShop::Matrix::Cell Namespace Reference

Integration of cell matrices. More...

Functions

template<int dim>
void mass (FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double a=1.)
 The mass matrix.
 
template<int dim>
void mass (FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, const unsigned int component, const double a=1.)
 The mass matrix between different elements.
 
template<int dim>
void mass_scaled (FullMatrix< double > &M, const FEValuesBase< dim > &fe, const std::vector< double > &a, bool divide)
 Mass matrix with varying coefficient.
 
template<int dim>
void mass_scaled (FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, const unsigned int select, const std::vector< double > &factor, bool divide=false)
 Mass matrix with varying coefficient for different elements.
 
template<int dim>
void advection (FullMatrix< double > &M, const FEValuesBase< dim > &fe, const std::vector< Point< dim > > &velocity)
 Linear advection in divergence form.
 
template<int dim>
void advection (FullMatrix< double > &M, const FEValuesBase< dim > &fe, const VectorSlice< const std::vector< std::vector< double > > > &velocity, double factor=1.)
 Linear advection in divergence form.
 
template<int dim>
void div_primal (std::vector< MatrixBlock< FullMatrix< double > > > &M, const std::vector< unsigned int > &indices, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, double factor=1.)
 Cell matrix for divergence.
 
template<int dim>
void div (std::vector< MatrixBlock< FullMatrix< double > > > &M, const std::vector< unsigned int > &indices, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, double factor=1.)
 Cell matrix for divergence.
 
template<int dim>
void div_scaled (FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, const std::vector< double > &factor)
 Same as div_cell(), but with a possibly varying coefficient.
 
template<int dim>
void di (FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, unsigned int i)
 Directional derivative in direction i.
 
template<int dim>
void di_scaled (FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, unsigned int d, std::vector< double > &factor)
 Same as di_cell(), but with varying coefficient.
 
template<int dim>
void eps_div (FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest)
 Discretization for the transpose of the symmetric gradient.
 
template<int dim>
void laplacian (FullMatrix< double > &M, const FEValuesBase< dim > &fe, double nu=1.)
 Laplacian in weak form.
 
template<int dim>
void laplacian (FullMatrix< double > &M, const FEValuesBase< dim > &fe, Tensor< 2, dim > nu)
 Laplacian in weak form with a Tensor coefficient to allow anisotropic media

\[ (K \nabla u, \nabla v) \]

.

 
template<int dim>
void laplacian_scaled (FullMatrix< double > &M, const FEValuesBase< dim > &fe, const std::vector< double > &nu, bool divide=false)
 Laplacian with varying coefficient.
 
template<int dim>
void laplacian_scaled (FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double nu, bool divide=false)
 Laplacian with constant coefficient.
 
template<int dim>
void div_div (FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double factor=1.)
 Second order divergence form for vector valued elements.
 

Detailed Description

Integration of cell matrices.

Function Documentation

template<int dim>
void AppShop::Matrix::Cell::advection ( FullMatrix< double > &  M,
const FEValuesBase< dim > &  fe,
const std::vector< Point< dim > > &  velocity 
)

Linear advection in divergence form.

The derivative is on the test function.

\[ -(u \mathbf w, \nabla v) \]

template<int dim>
void AppShop::Matrix::Cell::advection ( FullMatrix< double > &  M,
const FEValuesBase< dim > &  fe,
const VectorSlice< const std::vector< std::vector< double > > > &  velocity,
double  factor = 1. 
)

Linear advection in divergence form.

The derivative is on the test function.

\[ -(u \mathbf w, \nabla v) \]

template<int dim>
void AppShop::Matrix::Cell::di ( FullMatrix< double > &  M,
const FEValuesBase< dim > &  fe,
const FEValuesBase< dim > &  fetest,
unsigned int  i 
)

Directional derivative in direction i.

\[ -(u, \partial_i v) \]

template<int dim>
void AppShop::Matrix::Cell::di_scaled ( FullMatrix< double > &  M,
const FEValuesBase< dim > &  fe,
const FEValuesBase< dim > &  fetest,
unsigned int  d,
std::vector< double > &  factor 
)

Same as di_cell(), but with varying coefficient.

\[ -(\alpha u, \partial_i v) \]

template<int dim>
void AppShop::Matrix::Cell::div ( std::vector< MatrixBlock< FullMatrix< double > > > &  M,
const std::vector< unsigned int > &  indices,
const FEValuesBase< dim > &  fe,
const FEValuesBase< dim > &  fetest,
double  factor = 1. 
)

Cell matrix for divergence.

The derivative is on the test function.

\[ -(u , \nabla v) \]

  • M Result matrix vector
  • indices: the indices in the vector where we find the matrices.
  • fe: finite element space for trial functions.
  • fetest: finite element space for trial functions.
template<int dim>
void AppShop::Matrix::Cell::div_div ( FullMatrix< double > &  M,
const FEValuesBase< dim > &  fe,
const double  factor = 1. 
)

Second order divergence form for vector valued elements.

\[ (\nabla \cdot u, \nabla \cdot v) \]

template<int dim>
void AppShop::Matrix::Cell::div_primal ( std::vector< MatrixBlock< FullMatrix< double > > > &  M,
const std::vector< unsigned int > &  indices,
const FEValuesBase< dim > &  fe,
const FEValuesBase< dim > &  fetest,
double  factor = 1. 
)

Cell matrix for divergence.

The derivative is on the trial function.

\[ (\nabla \cdot u , v) \]

  • M Result matrix vector
  • indices: the indices in the vector where we find the matrices.
  • fe: finite element space for trial functions.
  • fetest: finite element space for trial functions.
template<int dim>
void AppShop::Matrix::Cell::div_scaled ( FullMatrix< double > &  M,
const FEValuesBase< dim > &  fe,
const FEValuesBase< dim > &  fetest,
const std::vector< double > &  factor 
)

Same as div_cell(), but with a possibly varying coefficient.

\[ -( \alpha u, \nabla v) \]

template<int dim>
void AppShop::Matrix::Cell::eps_div ( FullMatrix< double > &  M,
const FEValuesBase< dim > &  fe,
const FEValuesBase< dim > &  fetest 
)

Discretization for the transpose of the symmetric gradient.

\[ -(u, \varepsilon(v)) \]

template<int dim>
void AppShop::Matrix::Cell::laplacian ( FullMatrix< double > &  M,
const FEValuesBase< dim > &  fe,
double  nu = 1. 
)

Laplacian in weak form.

\[ (\nu \nabla u, \nabla v) \]

template<int dim>
void AppShop::Matrix::Cell::laplacian ( FullMatrix< double > &  M,
const FEValuesBase< dim > &  fe,
Tensor< 2, dim nu 
)

Laplacian in weak form with a Tensor coefficient to allow anisotropic media

\[ (K \nabla u, \nabla v) \]

.

template<int dim>
void AppShop::Matrix::Cell::laplacian_scaled ( FullMatrix< double > &  M,
const FEValuesBase< dim > &  fe,
const std::vector< double > &  nu,
bool  divide = false 
)

Laplacian with varying coefficient.

If divide is true, divide gradient u times gradient v by factor.

\[ (\nu \nabla u, \nabla v) \]

template<int dim>
void AppShop::Matrix::Cell::laplacian_scaled ( FullMatrix< double > &  M,
const FEValuesBase< dim > &  fe,
const double  nu,
bool  divide = false 
)

Laplacian with constant coefficient.

If divide is true, divide gradient u times gradient v by factor.

\[ (\nu \nabla u, \nabla v) \]

template<int dim>
void AppShop::Matrix::Cell::mass ( FullMatrix< double > &  M,
const FEValuesBase< dim > &  fe,
const double  a = 1. 
)

The mass matrix.

\[ (a u,v) \]

template<int dim>
void AppShop::Matrix::Cell::mass ( FullMatrix< double > &  M,
const FEValuesBase< dim > &  fe,
const FEValuesBase< dim > &  fetest,
const unsigned int  component,
const double  a = 1. 
)

The mass matrix between different elements.

\[ (a u,v) \]

If component has a nonnegative value and the element is vector-valued, then only the component corresponding to this number is assembled.

template<int dim>
void AppShop::Matrix::Cell::mass_scaled ( FullMatrix< double > &  M,
const FEValuesBase< dim > &  fe,
const std::vector< double > &  a,
bool  divide 
)

Mass matrix with varying coefficient.

\[ (a u,v) \]

If divide is true, then divide by a instead of multiplying.

template<int dim>
void AppShop::Matrix::Cell::mass_scaled ( FullMatrix< double > &  M,
const FEValuesBase< dim > &  fe,
const FEValuesBase< dim > &  fetest,
const unsigned int  select,
const std::vector< double > &  factor,
bool  divide = false 
)

Mass matrix with varying coefficient for different elements.

\[ (a u,v) \]

If divide is true, then divide by a instead of multiplying.