Flow123d
JS_before_hm-2206-gca2a8366f
|
Go to the documentation of this file.
18 #ifndef FE_VALUES_MAP_HH_
19 #define FE_VALUES_MAP_HH_
29 template<
unsigned int spacedim = 3>
41 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
42 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
51 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
52 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
61 template<
unsigned int spacedim = 3>
73 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
74 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
77 for (
unsigned int c=0; c<spacedim; c++)
87 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
88 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
92 for (
unsigned int c=0; c<spacedim; c++)
102 template<
unsigned int spacedim = 3>
114 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
115 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
118 for (
unsigned int c=0; c<spacedim; c++)
128 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
129 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
132 for (
unsigned int c=0; c<spacedim; c++)
142 template<
unsigned int spacedim = 3>
154 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
155 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
158 for (
unsigned int c=0; c<spacedim; c++)
168 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
169 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
172 for (
unsigned int c=0; c<spacedim; c++)
182 template<
unsigned int spacedim = 3>
194 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
195 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
198 for (
unsigned int c=0; c<spacedim*spacedim; c++)
199 fe_values.
shape_values[i][j*spacedim*spacedim+c] = fv_vec[c];
208 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
209 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
212 for (
unsigned int c=0; c<spacedim*spacedim; c++)
222 template<
unsigned int spacedim = 3>
230 unsigned int comp_offset = 0;
231 for (
unsigned int f=0; f<fe_values.
fe_sys_dofs_.size(); f++)
235 fe_values.
fe_values_vec[f].fill_data(elm_values, vec_fe_data);
247 unsigned int comp_offset = 0;
248 unsigned int shape_offset = 0;
249 for (
unsigned int f=0; f<fe_values.
fe_sys_dofs_.size(); f++)
252 for (
unsigned int i=0; i<fe_data.
n_points; i++)
253 for (
unsigned int n=0; n<fe_values.
fe_sys_dofs_[f].size(); n++)
269 unsigned int comp_offset = 0;
270 unsigned int shape_offset = 0;
271 for (
unsigned int f=0; f<fe_values.
fe_sys_dofs_.size(); f++)
274 for (
unsigned int i=0; i<fe_data.
n_points; i++)
275 for (
unsigned int n=0; n<fe_values.
fe_sys_dofs_[f].size(); n++)
287 #endif // FE_VALUES_MAP_HH_
std::vector< std::vector< double > > shape_values
Shape functions evaluated at the quadrature points.
Helper class allows update values and gradients of FEValues of FEVectorPiola type.
void fill_values_vec(FMT_UNUSED FEValues< spacedim > &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, FMT_UNUSED const typename FEValues< spacedim >::FEInternalData &fe_data)
Empty method.
void update_values(FEValues< spacedim > &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, const typename FEValues< spacedim >::FEInternalData &fe_data)
Update shape_values of given FEValues object.
std::vector< FEValues< spacedim > > fe_values_vec
Vector of FEValues for sub-elements of FESystem.
std::vector< std::vector< unsigned int > > fe_sys_dofs_
Dof indices of FESystem sub-elements.
Class FEValues calculates finite element data on the actual cells such as shape function values,...
void update_values(FEValues< spacedim > &fe_values, const ElementValues< spacedim > &elm_values, const typename FEValues< spacedim >::FEInternalData &fe_data)
Update shape_values of given FEValues object.
void update_gradients(FEValues< spacedim > &fe_values, const ElementValues< spacedim > &elm_values, const typename FEValues< spacedim >::FEInternalData &fe_data)
Update shape_gradients of given FEValues object.
void update_gradients(FEValues< spacedim > &fe_values, const ElementValues< spacedim > &elm_values, const typename FEValues< spacedim >::FEInternalData &fe_data)
Update shape_gradients of given FEValues object.
void fill_values_vec(FMT_UNUSED FEValues< spacedim > &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, FMT_UNUSED const typename FEValues< spacedim >::FEInternalData &fe_data)
Empty method.
void fill_values_vec(FEValues< spacedim > &fe_values, const ElementValues< spacedim > &elm_values, const typename FEValues< spacedim >::FEInternalData &fe_data)
Fill fe_values_vec of components of mixed system FEValues object.
Helper class allows update values and gradients of FEValues of FEVector type.
Calculates finite element data on the actual cell.
void update_values(FEValues< spacedim > &fe_values, const ElementValues< spacedim > &elm_values, const typename FEValues< spacedim >::FEInternalData &fe_data)
Update shape_values of given FEValues object.
void update_values(FEValues< spacedim > &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, const typename FEValues< spacedim >::FEInternalData &fe_data)
Update shape_values of given FEValues object.
void fill_values_vec(FMT_UNUSED FEValues< spacedim > &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, FMT_UNUSED const typename FEValues< spacedim >::FEInternalData &fe_data)
Empty method.
void update_gradients(FEValues< spacedim > &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, const typename FEValues< spacedim >::FEInternalData &fe_data)
Update shape_gradients of given FEValues object.
Structure for storing the precomputed finite element data.
ArmaMat< double, N, M > mat
std::vector< std::vector< arma::mat > > ref_shape_grads
Precomputed gradients of basis functions at the quadrature points.
Helper class allows update values and gradients of FEValues of FETensor type.
void update_values(FEValues< spacedim > &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, const typename FEValues< spacedim >::FEInternalData &fe_data)
Update shape_values of given FEValues object.
void update_gradients(FEValues< spacedim > &fe_values, const ElementValues< spacedim > &elm_values, const typename FEValues< spacedim >::FEInternalData &fe_data)
Update shape_gradients of given FEValues object.
arma::mat inverse_jacobian(const unsigned int point_no) const
Return inverse Jacobian matrix at point point_no.
arma::mat jacobian(const unsigned int point_no) const
Return Jacobian matrix at point point_no.
Helper class allows update values and gradients of FEValues of FEMixedSystem type.
std::vector< unsigned int > fe_sys_n_space_components_
Numbers of components of FESystem sub-elements in real space.
Class for computation of data on cell and side.
std::vector< std::vector< arma::vec > > ref_shape_values
Precomputed values of basis functions at the quadrature points.
std::vector< unsigned int > fe_sys_n_components_
Numbers of components of FESystem sub-elements in reference space.
std::vector< std::vector< arma::vec::fixed< spacedim > > > shape_gradients
Helper class allows update values and gradients of FEValues of FEScalar type.
unsigned int n_components_
Number of components of the FE.
void update_values(FEValues< spacedim > &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, const typename FEValues< spacedim >::FEInternalData &fe_data)
Update shape_values of given FEValues object.
double determinant(const unsigned int point_no) const
Return the relative volume change of the cell (Jacobian determinant).
void update_gradients(FEValues< spacedim > &fe_values, const ElementValues< spacedim > &elm_values, const typename FEValues< spacedim >::FEInternalData &fe_data)
Update shape_gradients of given FEValues object.
Class ElementValues calculates data related to transformation of reference cell to actual cell (Jacob...
unsigned int n_dofs
Number of dofs (shape functions).
void fill_values_vec(FMT_UNUSED FEValues< spacedim > &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, FMT_UNUSED const typename FEValues< spacedim >::FEInternalData &fe_data)
Empty method.
void fill_values_vec(FMT_UNUSED FEValues< spacedim > &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, FMT_UNUSED const typename FEValues< spacedim >::FEInternalData &fe_data)
Empty method.
Helper class allows update values and gradients of FEValues of FEVectorContravariant type.
void update_gradients(FEValues< spacedim > &fe_values, const ElementValues< spacedim > &elm_values, const typename FEValues< spacedim >::FEInternalData &fe_data)
Update shape_gradients of given FEValues object.
unsigned int n_points
Number of quadrature points.
FEType fe_type_
Type of finite element (scalar, vector, tensor).