27 template<
class FV,
unsigned int spacedim>
30 ASSERT_LT( function_no, fe_values_.n_dofs() );
31 ASSERT_LT( point_no, fe_values_.n_points() );
32 return fe_values_.shape_value_component(function_no, point_no, component_);
35 template<
class FV,
unsigned int spacedim>
38 ASSERT_LT( function_no, fe_values_.n_dofs() );
39 ASSERT_LT( point_no, fe_values_.n_points() );
40 return fe_values_.shape_grad_component(function_no, point_no, component_);
43 template<
class FV,
unsigned int spacedim>
45 {
return fe_values_; }
50 template<
class FV,
unsigned int spacedim>
53 ASSERT_LT( function_no, fe_values_.n_dofs() );
54 ASSERT_LT( point_no, fe_values_.n_points() );
55 arma::vec::fixed<spacedim> v;
56 for (
unsigned int c=0; c<spacedim; ++c)
57 v(c) = fe_values_.shape_value_component(function_no, point_no, first_vector_component_+c);
61 template<
class FV,
unsigned int spacedim>
64 ASSERT_LT( function_no, fe_values_.n_dofs() );
65 ASSERT_LT( point_no, fe_values_.n_points() );
66 arma::mat::fixed<spacedim,spacedim> g;
67 for (
unsigned int c=0; c<spacedim; ++c)
68 g.col(c) = fe_values_.shape_grad_component(function_no, point_no, first_vector_component_+c);
72 template<
class FV,
unsigned int spacedim>
75 ASSERT_LT( function_no, fe_values_.n_dofs() );
76 ASSERT_LT( point_no, fe_values_.n_points() );
77 arma::mat::fixed<spacedim,spacedim> g = grad(function_no, point_no);
78 return 0.5*(g+trans(g));
81 template<
class FV,
unsigned int spacedim>
84 ASSERT_LT( function_no, fe_values_.n_dofs() );
85 ASSERT_LT( point_no, fe_values_.n_points() );
87 for (
unsigned int c=0; c<spacedim; ++c)
88 div += fe_values_.shape_grad_component(function_no, point_no, first_vector_component_+c)(first_vector_component_+c);
92 template<
class FV,
unsigned int spacedim>
94 {
return fe_values_; }
98 template<
class FV,
unsigned int spacedim>
101 ASSERT_LT( function_no, fe_values_.n_dofs() );
102 ASSERT_LT( point_no, fe_values_.n_points() );
103 arma::mat::fixed<spacedim,spacedim> v;
104 for (
unsigned int c=0; c<spacedim*spacedim; ++c)
105 v(c/spacedim,c%spacedim) = fe_values_.shape_value_component(function_no, point_no, first_tensor_component_+c);
109 template<
class FV,
unsigned int spacedim>
111 unsigned int variable_no,
112 unsigned int function_no,
113 unsigned int point_no)
const
116 ASSERT_LT( function_no, fe_values_.n_dofs() );
117 ASSERT_LT( point_no, fe_values_.n_points() );
118 arma::mat::fixed<spacedim,spacedim> g;
119 for (
unsigned int c=0; c<spacedim*spacedim; ++c)
120 g(c/spacedim,c%spacedim) = fe_values_.shape_grad_component(function_no, point_no, first_tensor_component_+c)[first_tensor_component_+variable_no];
124 template<
class FV,
unsigned int spacedim>
127 ASSERT_LT( function_no, fe_values_.n_dofs() );
128 ASSERT_LT( point_no, fe_values_.n_points() );
129 arma::vec::fixed<spacedim> div;
131 for (
unsigned int c=0; c<spacedim*spacedim; ++c)
132 div(c%spacedim) += fe_values_.shape_grad_component(function_no, point_no, first_tensor_component_+c)[first_tensor_component_+c/spacedim];
136 template<
class FV,
unsigned int spacedim>
138 {
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.
Basic definitions of numerical quadrature rules.