28 template<
class FV,
unsigned int spacedim>
31 ASSERT_LT( function_no, fe_values_.n_dofs() );
32 ASSERT_LT( point_no, fe_values_.n_points() );
33 return fe_values_.shape_value_component(function_no, point_no, component_);
36 template<
class FV,
unsigned int spacedim>
39 ASSERT_LT( function_no, fe_values_.n_dofs() );
40 ASSERT_LT( point_no, fe_values_.n_points() );
41 return fe_values_.shape_grad_component(function_no, point_no, component_);
44 template<
class FV,
unsigned int spacedim>
46 {
return fe_values_; }
51 template<
class FV,
unsigned int spacedim>
54 ASSERT_LT( function_no, fe_values_.n_dofs() );
55 ASSERT_LT( point_no, fe_values_.n_points() );
56 arma::vec::fixed<spacedim> v;
57 for (
unsigned int c=0; c<spacedim; ++c)
58 v(c) = fe_values_.shape_value_component(function_no, point_no, first_vector_component_+c);
62 template<
class FV,
unsigned int spacedim>
65 ASSERT_LT( function_no, fe_values_.n_dofs() );
66 ASSERT_LT( point_no, fe_values_.n_points() );
67 arma::mat::fixed<spacedim,spacedim> g;
68 for (
unsigned int c=0; c<spacedim; ++c)
69 g.col(c) = fe_values_.shape_grad_component(function_no, point_no, first_vector_component_+c);
73 template<
class FV,
unsigned int spacedim>
76 ASSERT_LT( function_no, fe_values_.n_dofs() );
77 ASSERT_LT( point_no, fe_values_.n_points() );
78 arma::mat::fixed<spacedim,spacedim> g = grad(function_no, point_no);
79 return 0.5*(g+trans(g));
82 template<
class FV,
unsigned int spacedim>
85 ASSERT_LT( function_no, fe_values_.n_dofs() );
86 ASSERT_LT( point_no, fe_values_.n_points() );
88 for (
unsigned int c=0; c<spacedim; ++c)
89 div += fe_values_.shape_grad_component(function_no, point_no, first_vector_component_+c)(first_vector_component_+c);
93 template<
class FV,
unsigned int spacedim>
95 {
return fe_values_; }
99 template<
class FV,
unsigned int spacedim>
102 ASSERT_LT( function_no, fe_values_.n_dofs() );
104 arma::mat::fixed<spacedim,spacedim> v;
105 for (
unsigned int c=0; c<spacedim*spacedim; ++c)
106 v(c/spacedim,c%spacedim) = fe_values_.shape_value_component(function_no, point_no, first_tensor_component_+c);
110 template<
class FV,
unsigned int spacedim>
112 unsigned int variable_no,
113 unsigned int function_no,
114 unsigned int point_no)
const
117 ASSERT_LT( function_no, fe_values_.n_dofs() );
118 ASSERT_LT( point_no, fe_values_.n_points() );
119 arma::mat::fixed<spacedim,spacedim> g;
120 for (
unsigned int c=0; c<spacedim*spacedim; ++c)
121 g(c/spacedim,c%spacedim) = fe_values_.shape_grad_component(function_no, point_no, first_tensor_component_+c)[first_tensor_component_+variable_no];
125 template<
class FV,
unsigned int spacedim>
128 ASSERT_LT( function_no, fe_values_.n_dofs() );
129 ASSERT_LT( point_no, fe_values_.n_points() );
130 arma::vec::fixed<spacedim> div;
132 for (
unsigned int c=0; c<spacedim*spacedim; ++c)
133 div(c%spacedim) += fe_values_.shape_grad_component(function_no, point_no, first_tensor_component_+c)[first_tensor_component_+c/spacedim];
137 template<
class FV,
unsigned int spacedim>
139 {
return fe_values_; }
#define ASSERT_LT(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
const FV & base() const
Returns the FEValues class.
double value(unsigned int function_no, unsigned int point_no) const
Return value of scalar shape function.
arma::vec::fixed< spacedim > grad(unsigned int function_no, unsigned int point_no) const
Return gradient 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.
arma::mat::fixed< spacedim, spacedim > value(unsigned int function_no, unsigned int point_no) const
Return value of tensor-valued shape function.
arma::vec::fixed< spacedim > divergence(unsigned int function_no, unsigned int point_no) const
Return divergence of tensor-valued shape function.
const FV & base() const
Returns the FEValues class.
arma::mat::fixed< spacedim, spacedim > sym_grad(unsigned int function_no, unsigned int point_no) const
Return symmetric gradient of vector-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 FV & base() const
Returns the FEValues class.
double divergence(unsigned int function_no, unsigned int point_no) const
Return divergence of vector-valued shape function.
arma::vec::fixed< spacedim > value(unsigned int function_no, unsigned int point_no) const
Return value of vector-valued shape function.
Class FEValues calculates finite element data on the actual cells such as shape function values,...
Abstract class for description of finite elements.
Class FEValues calculates finite element data on the actual cells such as shape function values,...
Basic definitions of numerical quadrature rules.