Flow123d  DF_mechanic_bench-4968b1b
fe_values_map.hh
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file fe_values_map.hh
15  * @author Jan Stebel, David Flanderka
16  */
17 
18 #ifndef FE_VALUES_MAP_HH_
19 #define FE_VALUES_MAP_HH_
20 
21 #include <vector> // for vector
22 #include "fem/fe_values.hh" // for FEValues
23 #include "fem/element_values.hh" // for ElementValues
24 
25 
26 /**
27  * @brief Helper class allows update values and gradients of FEValues of FEScalar type.
28  */
29 template<unsigned int spacedim = 3>
30 class MapScalar {
31 public:
32  /// Empty method.
34  FMT_UNUSED const typename FEValues<spacedim>::FEInternalData &fe_data) {}
35 
36  /// Update shape_values of given FEValues object.
37  inline void update_values(FEValues<spacedim> &fe_values, FMT_UNUSED const ElementValues<spacedim> &elm_values,
38  const typename FEValues<spacedim>::FEInternalData &fe_data) {
39  ASSERT(fe_values.fe_type_ == FEScalar);
40 
41  for (unsigned int i = 0; i < fe_data.n_points; i++)
42  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
43  fe_values.shape_values[i][j] = fe_data.ref_shape_values[i][j][0];
44  }
45 
46  /// Update shape_gradients of given FEValues object.
47  inline void update_gradients(FEValues<spacedim> &fe_values, const ElementValues<spacedim> &elm_values,
48  const typename FEValues<spacedim>::FEInternalData &fe_data) {
49  ASSERT(fe_values.fe_type_ == FEScalar);
50 
51  for (unsigned int i = 0; i < fe_data.n_points; i++)
52  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
53  fe_values.shape_gradients[i][j] = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j];
54  }
55 };
56 
57 
58 /**
59  * @brief Helper class allows update values and gradients of FEValues of FEVectorPiola type
60  */
61 template<unsigned int spacedim = 3>
62 class MapPiola {
63 public:
64  /// Empty method.
66  FMT_UNUSED const typename FEValues<spacedim>::FEInternalData &fe_data) {}
67 
68  /// Update shape_values of given FEValues object.
69  inline void update_values(FEValues<spacedim> &fe_values, const ElementValues<spacedim> &elm_values,
70  const typename FEValues<spacedim>::FEInternalData &fe_data) {
71  ASSERT(fe_values.fe_type_ == FEVectorPiola);
72 
73  for (unsigned int i = 0; i < fe_data.n_points; i++)
74  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
75  {
76  arma::vec fv_vec = elm_values.jacobian(i)*fe_data.ref_shape_values[i][j]/elm_values.determinant(i);
77  for (unsigned int c=0; c<spacedim; c++)
78  fe_values.shape_values[i][j*spacedim+c] = fv_vec(c);
79  }
80  }
81 
82  /// Update shape_gradients of given FEValues object.
83  inline void update_gradients(FEValues<spacedim> &fe_values, const ElementValues<spacedim> &elm_values,
84  const typename FEValues<spacedim>::FEInternalData &fe_data) {
85  ASSERT(fe_values.fe_type_ == FEVectorPiola);
86 
87  for (unsigned int i = 0; i < fe_data.n_points; i++)
88  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
89  {
90  arma::mat grads = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j] * trans(elm_values.jacobian(i))
91  / elm_values.determinant(i);
92  for (unsigned int c=0; c<spacedim; c++)
93  fe_values.shape_gradients[i][j*spacedim+c] = grads.col(c);
94  }
95  }
96 };
97 
98 
99 /**
100  * @brief Helper class allows update values and gradients of FEValues of FEVectorContravariant type
101  */
102 template<unsigned int spacedim = 3>
104 public:
105  /// Empty method.
107  FMT_UNUSED const typename FEValues<spacedim>::FEInternalData &fe_data) {}
108 
109  /// Update shape_values of given FEValues object.
110  inline void update_values(FEValues<spacedim> &fe_values, const ElementValues<spacedim> &elm_values,
111  const typename FEValues<spacedim>::FEInternalData &fe_data) {
112  ASSERT(fe_values.fe_type_ == FEVectorContravariant);
113 
114  for (unsigned int i = 0; i < fe_data.n_points; i++)
115  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
116  {
117  arma::vec fv_vec = elm_values.jacobian(i) * fe_data.ref_shape_values[i][j];
118  for (unsigned int c=0; c<spacedim; c++)
119  fe_values.shape_values[i][j*spacedim+c] = fv_vec[c];
120  }
121  }
122 
123  /// Update shape_gradients of given FEValues object.
124  inline void update_gradients(FEValues<spacedim> &fe_values, const ElementValues<spacedim> &elm_values,
125  const typename FEValues<spacedim>::FEInternalData &fe_data) {
126  ASSERT(fe_values.fe_type_ == FEVectorContravariant);
127 
128  for (unsigned int i = 0; i < fe_data.n_points; i++)
129  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
130  {
131  arma::mat grads = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j] * trans(elm_values.jacobian(i));
132  for (unsigned int c=0; c<spacedim; c++)
133  fe_values.shape_gradients[i][j*spacedim+c] = grads.col(c);
134  }
135  }
136 };
137 
138 
139 /**
140  * @brief Helper class allows update values and gradients of FEValues of FEVector type
141  */
142 template<unsigned int spacedim = 3>
143 class MapVector {
144 public:
145  /// Empty method.
147  FMT_UNUSED const typename FEValues<spacedim>::FEInternalData &fe_data) {}
148 
149  /// Update shape_values of given FEValues object.
150  inline void update_values(FEValues<spacedim> &fe_values, FMT_UNUSED const ElementValues<spacedim> &elm_values,
151  const typename FEValues<spacedim>::FEInternalData &fe_data) {
152  ASSERT(fe_values.fe_type_ == FEVector);
153 
154  for (unsigned int i = 0; i < fe_data.n_points; i++)
155  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
156  {
157  arma::vec fv_vec = fe_data.ref_shape_values[i][j];
158  for (unsigned int c=0; c<spacedim; c++)
159  fe_values.shape_values[i][j*spacedim+c] = fv_vec[c];
160  }
161  }
162 
163  /// Update shape_gradients of given FEValues object.
164  inline void update_gradients(FEValues<spacedim> &fe_values, const ElementValues<spacedim> &elm_values,
165  const typename FEValues<spacedim>::FEInternalData &fe_data) {
166  ASSERT(fe_values.fe_type_ == FEVector);
167 
168  for (unsigned int i = 0; i < fe_data.n_points; i++)
169  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
170  {
171  arma::mat grads = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j];
172  for (unsigned int c=0; c<spacedim; c++)
173  fe_values.shape_gradients[i][j*spacedim+c] = grads.col(c);
174  }
175  }
176 };
177 
178 
179 /**
180  * @brief Helper class allows update values and gradients of FEValues of FETensor type
181  */
182 template<unsigned int spacedim = 3>
183 class MapTensor {
184 public:
185  /// Empty method.
187  FMT_UNUSED const typename FEValues<spacedim>::FEInternalData &fe_data) {}
188 
189  /// Update shape_values of given FEValues object.
190  inline void update_values(FEValues<spacedim> &fe_values, FMT_UNUSED const ElementValues<spacedim> &elm_values,
191  const typename FEValues<spacedim>::FEInternalData &fe_data) {
192  ASSERT(fe_values.fe_type_ == FETensor);
193 
194  for (unsigned int i = 0; i < fe_data.n_points; i++)
195  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
196  {
197  arma::vec fv_vec = fe_data.ref_shape_values[i][j];
198  for (unsigned int c=0; c<spacedim*spacedim; c++)
199  fe_values.shape_values[i][j*spacedim*spacedim+c] = fv_vec[c];
200  }
201  }
202 
203  /// Update shape_gradients of given FEValues object.
204  inline void update_gradients(FEValues<spacedim> &fe_values, const ElementValues<spacedim> &elm_values,
205  const typename FEValues<spacedim>::FEInternalData &fe_data) {
206  ASSERT(fe_values.fe_type_ == FETensor);
207 
208  for (unsigned int i = 0; i < fe_data.n_points; i++)
209  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
210  {
211  arma::mat grads = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j];
212  for (unsigned int c=0; c<spacedim*spacedim; c++)
213  fe_values.shape_gradients[i][j*spacedim*spacedim+c] = grads.col(c);
214  }
215  }
216 };
217 
218 
219 /**
220  * @brief Helper class allows update values and gradients of FEValues of FEMixedSystem type
221  */
222 template<unsigned int spacedim = 3>
223 class MapSystem {
224 public:
225  /// Fill fe_values_vec of components of mixed system FEValues object.
226  inline void fill_values_vec(FEValues<spacedim> &fe_values, const ElementValues<spacedim> &elm_values,
227  const typename FEValues<spacedim>::FEInternalData &fe_data) {
228  ASSERT(fe_values.fe_type_ == FEMixedSystem);
229 
230  unsigned int comp_offset = 0;
231  for (unsigned int f=0; f<fe_values.fe_sys_dofs_.size(); f++)
232  {
233  // fill fe_values for base FE
234  typename FEValues<spacedim>::FEInternalData vec_fe_data(fe_data, fe_values.fe_sys_dofs_[f], comp_offset, fe_values.fe_sys_n_components_[f]);
235  fe_values.fe_values_vec[f].fill_data(elm_values, vec_fe_data);
236 
237  comp_offset += fe_values.fe_sys_n_components_[f];
238  }
239  }
240 
241  /// Update shape_values of given FEValues object.
242  inline void update_values(FEValues<spacedim> &fe_values, FMT_UNUSED const ElementValues<spacedim> &elm_values,
243  const typename FEValues<spacedim>::FEInternalData &fe_data) {
244  ASSERT(fe_values.fe_type_ == FEMixedSystem);
245 
246  arma::vec fv_vec;
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++)
250  {
251  // gather fe_values in vectors for FESystem
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++)
254  for (unsigned int c=0; c<fe_values.fe_sys_n_space_components_[f]; c++)
255  fe_values.shape_values[i][shape_offset+fe_values.n_components_*n+comp_offset+c] =
256  fe_values.fe_values_vec[f].shape_values[i][n*fe_values.fe_sys_n_space_components_[f]+c];
257 
258  comp_offset += fe_values.fe_sys_n_space_components_[f];
259  shape_offset += fe_values.fe_sys_dofs_[f].size()*fe_values.n_components_;
260  }
261  }
262 
263  /// Update shape_gradients of given FEValues object.
264  inline void update_gradients(FEValues<spacedim> &fe_values, FMT_UNUSED const ElementValues<spacedim> &elm_values,
265  const typename FEValues<spacedim>::FEInternalData &fe_data) {
266  ASSERT(fe_values.fe_type_ == FEMixedSystem);
267 
268  arma::mat grads;
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++)
272  {
273  // gather fe_values in vectors for FESystem
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++)
276  for (unsigned int c=0; c<fe_values.fe_sys_n_space_components_[f]; c++)
277  fe_values.shape_gradients[i][shape_offset+fe_values.n_components_*n+comp_offset+c] =
278  fe_values.fe_values_vec[f].shape_gradients[i][n*fe_values.fe_sys_n_space_components_[f]+c];
279 
280  comp_offset += fe_values.fe_sys_n_space_components_[f];
281  shape_offset += fe_values.fe_sys_dofs_[f].size()*fe_values.n_components_;
282  }
283  }
284 };
285 
286 
287 #endif // FE_VALUES_MAP_HH_
#define ASSERT(expr)
Definition: asserts.hh:351
Class for computation of data on cell and side.
arma::mat jacobian(const unsigned int point_no) const
Return Jacobian matrix at point point_no.
double determinant(const unsigned int point_no) const
Return the relative volume change of the cell (Jacobian determinant).
arma::mat inverse_jacobian(const unsigned int point_no) const
Return inverse Jacobian matrix at point point_no.
Structure for storing the precomputed finite element data.
Definition: fe_values.hh:316
std::vector< std::vector< arma::vec > > ref_shape_values
Precomputed values of basis functions at the quadrature points.
Definition: fe_values.hh:334
unsigned int n_dofs
Number of dofs (shape functions).
Definition: fe_values.hh:349
std::vector< std::vector< arma::mat > > ref_shape_grads
Precomputed gradients of basis functions at the quadrature points.
Definition: fe_values.hh:343
unsigned int n_points
Number of quadrature points.
Definition: fe_values.hh:346
Calculates finite element data on the actual cell.
Definition: fe_values.hh:67
std::vector< unsigned int > fe_sys_n_space_components_
Numbers of components of FESystem sub-elements in real space.
Definition: fe_values.hh:396
std::vector< std::vector< double > > shape_values
Shape functions evaluated at the quadrature points.
Definition: fe_values.hh:399
FEType fe_type_
Type of finite element (scalar, vector, tensor).
Definition: fe_values.hh:387
std::vector< unsigned int > fe_sys_n_components_
Numbers of components of FESystem sub-elements in reference space.
Definition: fe_values.hh:393
std::vector< std::vector< unsigned int > > fe_sys_dofs_
Dof indices of FESystem sub-elements.
Definition: fe_values.hh:390
std::vector< std::vector< arma::vec::fixed< spacedim > > > shape_gradients
Definition: fe_values.hh:403
unsigned int n_components_
Number of components of the FE.
Definition: fe_values.hh:415
std::vector< FEValues< spacedim > > fe_values_vec
Vector of FEValues for sub-elements of FESystem.
Definition: fe_values.hh:412
Helper class allows update values and gradients of FEValues of FEVectorContravariant type.
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 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, const ElementValues< spacedim > &elm_values, const typename FEValues< spacedim >::FEInternalData &fe_data)
Update shape_gradients of given FEValues object.
Helper class allows update values and gradients of FEValues of FEVectorPiola 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.
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 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 FEScalar 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_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_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.
Helper class allows update values and gradients of FEValues of FEMixedSystem type.
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.
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.
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.
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.
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 FEVector 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_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_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.
Class ElementValues calculates data related to transformation of reference cell to actual cell (Jacob...
Class FEValues calculates finite element data on the actual cells such as shape function values,...
@ FEScalar
@ FEMixedSystem
@ FEVectorPiola
@ FETensor
@ FEVectorContravariant
@ FEVector
ArmaMat< double, N, M > mat
Definition: armor.hh:936
ArmaVec< double, N > vec
Definition: armor.hh:933
#define FMT_UNUSED
Definition: posix.h:75