Flow123d
JS_before_hm-2198-g122e1f2e2
|
Go to the documentation of this file.
19 #ifndef FE_VALUES_VIEWS_HH
20 #define FE_VALUES_VIEWS_HH
25 template<
unsigned int spacedim>
class FEValues;
37 template<
unsigned int spacedim = 3>
52 double value(
unsigned int function_no,
unsigned int point_no)
const;
59 arma::vec::fixed<spacedim>
grad(
unsigned int function_no,
unsigned int point_no)
const;
74 template<
unsigned int spacedim = 3>
89 arma::vec::fixed<spacedim>
value(
unsigned int function_no,
unsigned int point_no)
const;
96 arma::mat::fixed<spacedim,spacedim>
grad(
unsigned int function_no,
unsigned int point_no)
const;
103 arma::mat::fixed<spacedim,spacedim>
sym_grad(
unsigned int function_no,
unsigned int point_no)
const;
110 double divergence(
unsigned int function_no,
unsigned int point_no)
const;
125 template<
unsigned int spacedim = 3>
140 arma::mat::fixed<spacedim,spacedim>
value(
unsigned int function_no,
unsigned int point_no)
const;
148 arma::mat::fixed<spacedim,spacedim>
derivative(
149 unsigned int variable_no,
150 unsigned int function_no,
151 unsigned int point_no)
const;
162 arma::vec::fixed<spacedim>
divergence(
unsigned int function_no,
unsigned int point_no)
const;
unsigned int component_
Index of the scalar component.
const FEValues< spacedim > & fe_values_
Base FEValues class for access to the FE.
double divergence(unsigned int function_no, unsigned int point_no) const
Return divergence of vector-valued shape function.
arma::vec::fixed< spacedim > grad(unsigned int function_no, unsigned int point_no) const
Return gradient of scalar shape function.
double value(unsigned int function_no, unsigned int point_no) const
Return value of scalar shape function.
arma::mat::fixed< spacedim, spacedim > derivative(unsigned int variable_no, unsigned int function_no, unsigned int point_no) const
Return partial derivative of tensor-valued shape function.
const FEValues< spacedim > & base() const
Returns the FEValues class.
unsigned int first_tensor_component_
Index of the first component of the vector.
Calculates finite element data on the actual cell.
Vector(const FEValues< spacedim > &fe_values, unsigned int component)
arma::mat::fixed< spacedim, spacedim > sym_grad(unsigned int function_no, unsigned int point_no) const
Return symmetric gradient of vector-valued shape function.
Scalar(const FEValues< spacedim > &fe_values, unsigned int component)
Tensor(const FEValues< spacedim > &fe_values, unsigned int component)
const FEValues< spacedim > & fe_values_
Base FEValues class for access to the FE.
arma::vec::fixed< spacedim > divergence(unsigned int function_no, unsigned int point_no) const
Return divergence of tensor-valued shape function.
const FEValues< spacedim > & base() const
Returns the FEValues class.
arma::mat::fixed< spacedim, spacedim > value(unsigned int function_no, unsigned int point_no) const
Return value of tensor-valued shape function.
arma::mat::fixed< spacedim, spacedim > grad(unsigned int function_no, unsigned int point_no) const
Return gradient of vector-valued shape function.
const FEValues< spacedim > & base() const
Returns the FEValues class.
arma::vec::fixed< spacedim > value(unsigned int function_no, unsigned int point_no) const
Return value of vector-valued shape function.
const FEValues< spacedim > & fe_values_
Base FEValues class for access to the FE.
unsigned int first_vector_component_
Index of the first component of the vector.