Flow123d  DF_patch_fe_data_tables-b678bc1
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<class FV, unsigned int spacedim = 3>
30 class MapScalar {
31 public:
32  /// Empty method.
33  inline void fill_values_vec(FMT_UNUSED FV &fe_values, FMT_UNUSED const ElementValues<spacedim> &elm_values,
34  FMT_UNUSED const FEInternalData &fe_data) {}
35 
36  /// Update shape_values of given FEValues object.
37  inline void update_values(FV &fe_values, FMT_UNUSED const ElementValues<spacedim> &elm_values,
38  const 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.set_shape_value(i, j, fe_data.ref_shape_values[i][j][0] );
44  }
45  }
46 
47  /// Update shape_gradients of given FEValues object.
48  inline void update_gradients(FV &fe_values, const ElementValues<spacedim> &elm_values,
49  const FEInternalData &fe_data) {
50  ASSERT(fe_values.fe_type_ == FEScalar);
51 
52  for (unsigned int i = 0; i < fe_data.n_points; i++)
53  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
54  fe_values.set_shape_gradient(i, j, trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j] );
55  }
56 };
57 
58 
59 /**
60  * @brief Helper class allows update values and gradients of FEValues of FEVectorPiola type
61  */
62 template<class FV, unsigned int spacedim = 3>
63 class MapPiola {
64 public:
65  /// Empty method.
66  inline void fill_values_vec(FMT_UNUSED FV &fe_values, FMT_UNUSED const ElementValues<spacedim> &elm_values,
67  FMT_UNUSED const FEInternalData &fe_data) {}
68 
69  /// Update shape_values of given FEValues object.
70  inline void update_values(FV &fe_values, const ElementValues<spacedim> &elm_values,
71  const FEInternalData &fe_data) {
72  ASSERT(fe_values.fe_type_ == FEVectorPiola);
73 
74  for (unsigned int i = 0; i < fe_data.n_points; i++)
75  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
76  {
77  arma::vec fv_vec = elm_values.jacobian(i)*fe_data.ref_shape_values[i][j]/elm_values.determinant(i);
78  for (unsigned int c=0; c<spacedim; c++)
79  fe_values.set_shape_value(i, j*spacedim+c, fv_vec(c));
80  }
81  }
82 
83  /// Update shape_gradients of given FEValues object.
84  inline void update_gradients(FV &fe_values, const ElementValues<spacedim> &elm_values,
85  const FEInternalData &fe_data) {
86  ASSERT(fe_values.fe_type_ == FEVectorPiola);
87 
88  for (unsigned int i = 0; i < fe_data.n_points; i++)
89  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
90  {
91  arma::mat grads = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j] * trans(elm_values.jacobian(i))
92  / elm_values.determinant(i);
93  for (unsigned int c=0; c<spacedim; c++)
94  fe_values.set_shape_gradient(i, j*spacedim+c, grads.col(c));
95  }
96  }
97 };
98 
99 
100 /**
101  * @brief Helper class allows update values and gradients of FEValues of FEVectorContravariant type
102  */
103 template<class FV, unsigned int spacedim = 3>
105 public:
106  /// Empty method.
107  inline void fill_values_vec(FMT_UNUSED FV &fe_values, FMT_UNUSED const ElementValues<spacedim> &elm_values,
108  FMT_UNUSED const FEInternalData &fe_data) {}
109 
110  /// Update shape_values of given FEValues object.
111  inline void update_values(FV &fe_values, const ElementValues<spacedim> &elm_values,
112  const FEInternalData &fe_data) {
113  ASSERT(fe_values.fe_type_ == FEVectorContravariant);
114 
115  for (unsigned int i = 0; i < fe_data.n_points; i++)
116  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
117  {
118  arma::vec fv_vec = elm_values.jacobian(i) * fe_data.ref_shape_values[i][j];
119  for (unsigned int c=0; c<spacedim; c++)
120  fe_values.set_shape_value(i, j*spacedim+c, fv_vec[c]);
121  }
122  }
123 
124  /// Update shape_gradients of given FEValues object.
125  inline void update_gradients(FV &fe_values, const ElementValues<spacedim> &elm_values,
126  const FEInternalData &fe_data) {
127  ASSERT(fe_values.fe_type_ == FEVectorContravariant);
128 
129  for (unsigned int i = 0; i < fe_data.n_points; i++)
130  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
131  {
132  arma::mat grads = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j] * trans(elm_values.jacobian(i));
133  for (unsigned int c=0; c<spacedim; c++)
134  fe_values.set_shape_gradient(i, j*spacedim+c, grads.col(c));
135  }
136  }
137 };
138 
139 
140 /**
141  * @brief Helper class allows update values and gradients of FEValues of FEVector type
142  */
143 template<class FV, unsigned int spacedim = 3>
144 class MapVector {
145 public:
146  /// Empty method.
147  inline void fill_values_vec(FMT_UNUSED FV &fe_values, FMT_UNUSED const ElementValues<spacedim> &elm_values,
148  FMT_UNUSED const FEInternalData &fe_data) {}
149 
150  /// Update shape_values of given FEValues object.
151  inline void update_values(FV &fe_values, FMT_UNUSED const ElementValues<spacedim> &elm_values,
152  const FEInternalData &fe_data) {
153  ASSERT(fe_values.fe_type_ == FEVector);
154 
155  for (unsigned int i = 0; i < fe_data.n_points; i++)
156  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
157  {
158  arma::vec fv_vec = fe_data.ref_shape_values[i][j];
159  for (unsigned int c=0; c<spacedim; c++)
160  fe_values.set_shape_value(i, j*spacedim+c, fv_vec[c]);
161  }
162  }
163 
164  /// Update shape_gradients of given FEValues object.
165  inline void update_gradients(FV &fe_values, const ElementValues<spacedim> &elm_values,
166  const FEInternalData &fe_data) {
167  ASSERT(fe_values.fe_type_ == FEVector);
168 
169  for (unsigned int i = 0; i < fe_data.n_points; i++)
170  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
171  {
172  arma::mat grads = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j];
173  for (unsigned int c=0; c<spacedim; c++)
174  fe_values.set_shape_gradient(i, j*spacedim+c, grads.col(c));
175  }
176  }
177 };
178 
179 
180 /**
181  * @brief Helper class allows update values and gradients of FEValues of FETensor type
182  */
183 template<class FV, unsigned int spacedim = 3>
184 class MapTensor {
185 public:
186  /// Empty method.
187  inline void fill_values_vec(FMT_UNUSED FV &fe_values, FMT_UNUSED const ElementValues<spacedim> &elm_values,
188  FMT_UNUSED const FEInternalData &fe_data) {}
189 
190  /// Update shape_values of given FEValues object.
191  inline void update_values(FV &fe_values, FMT_UNUSED const ElementValues<spacedim> &elm_values,
192  const FEInternalData &fe_data) {
193  ASSERT(fe_values.fe_type_ == FETensor);
194 
195  for (unsigned int i = 0; i < fe_data.n_points; i++)
196  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
197  {
198  arma::vec fv_vec = fe_data.ref_shape_values[i][j];
199  for (unsigned int c=0; c<spacedim*spacedim; c++)
200  fe_values.set_shape_value(i, j*spacedim*spacedim+c, fv_vec[c]);
201  }
202  }
203 
204  /// Update shape_gradients of given FEValues object.
205  inline void update_gradients(FV &fe_values, const ElementValues<spacedim> &elm_values,
206  const FEInternalData &fe_data) {
207  ASSERT(fe_values.fe_type_ == FETensor);
208 
209  for (unsigned int i = 0; i < fe_data.n_points; i++)
210  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
211  {
212  arma::mat grads = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j];
213  for (unsigned int c=0; c<spacedim*spacedim; c++)
214  fe_values.set_shape_gradient(i, j*spacedim*spacedim+c, grads.col(c));
215  }
216  }
217 };
218 
219 
220 /**
221  * @brief Helper class allows update values and gradients of FEValues of FEMixedSystem type
222  */
223 template<class FV, unsigned int spacedim = 3>
224 class MapSystem {
225 public:
226  /// Fill fe_values_vec of components of mixed system FEValues object.
227  inline void fill_values_vec(FV &fe_values, const ElementValues<spacedim> &elm_values,
228  const FEInternalData &fe_data) {
229  ASSERT(fe_values.fe_type_ == FEMixedSystem);
230 
231  unsigned int comp_offset = 0;
232  for (unsigned int f=0; f<fe_values.fe_sys_dofs_.size(); f++)
233  {
234  // fill fe_values for base FE
235  FEInternalData vec_fe_data(fe_data, fe_values.fe_sys_dofs_[f], comp_offset, fe_values.fe_sys_n_components_[f]);
236  fe_values.fe_values_vec[f].fill_data(elm_values, vec_fe_data);
237 
238  comp_offset += fe_values.fe_sys_n_components_[f];
239  }
240  }
241 
242  /// Update shape_values of given FEValues object.
243  inline void update_values(FV &fe_values, FMT_UNUSED const ElementValues<spacedim> &elm_values,
244  const FEInternalData &fe_data) {
245  ASSERT(fe_values.fe_type_ == FEMixedSystem);
246 
247  arma::vec fv_vec;
248  unsigned int comp_offset = 0;
249  unsigned int shape_offset = 0;
250  for (unsigned int f=0; f<fe_values.fe_sys_dofs_.size(); f++)
251  {
252  // gather fe_values in vectors for FESystem
253  for (unsigned int i=0; i<fe_data.n_points; i++)
254  for (unsigned int n=0; n<fe_values.fe_sys_dofs_[f].size(); n++)
255  for (unsigned int c=0; c<fe_values.fe_sys_n_space_components_[f]; c++)
256  fe_values.set_shape_value( i, shape_offset+fe_values.n_components_*n+comp_offset+c,
257  fe_values.fe_values_vec[f].shape_value(i, n*fe_values.fe_sys_n_space_components_[f]+c) );
258 
259  comp_offset += fe_values.fe_sys_n_space_components_[f];
260  shape_offset += fe_values.fe_sys_dofs_[f].size()*fe_values.n_components_;
261  }
262  }
263 
264  /// Update shape_gradients of given FEValues object.
265  inline void update_gradients(FV &fe_values, FMT_UNUSED const ElementValues<spacedim> &elm_values,
266  const FEInternalData &fe_data) {
267  ASSERT(fe_values.fe_type_ == FEMixedSystem);
268 
269  arma::mat grads;
270  unsigned int comp_offset = 0;
271  unsigned int shape_offset = 0;
272  for (unsigned int f=0; f<fe_values.fe_sys_dofs_.size(); f++)
273  {
274  // gather fe_values in vectors for FESystem
275  for (unsigned int i=0; i<fe_data.n_points; i++)
276  for (unsigned int n=0; n<fe_values.fe_sys_dofs_[f].size(); n++)
277  for (unsigned int c=0; c<fe_values.fe_sys_n_space_components_[f]; c++)
278  fe_values.set_shape_gradient(
279  i,
280  shape_offset+fe_values.n_components_*n+comp_offset+c,
281  fe_values.fe_values_vec[f].shape_grad(i, n*fe_values.fe_sys_n_space_components_[f]+c) );
282 
283  comp_offset += fe_values.fe_sys_n_space_components_[f];
284  shape_offset += fe_values.fe_sys_dofs_[f].size()*fe_values.n_components_;
285  }
286  }
287 };
288 
289 
290 #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:62
std::vector< std::vector< arma::mat > > ref_shape_grads
Precomputed gradients of basis functions at the quadrature points.
Definition: fe_values.hh:89
unsigned int n_points
Number of quadrature points.
Definition: fe_values.hh:92
unsigned int n_dofs
Number of dofs (shape functions).
Definition: fe_values.hh:95
std::vector< std::vector< arma::vec > > ref_shape_values
Precomputed values of basis functions at the quadrature points.
Definition: fe_values.hh:80
Helper class allows update values and gradients of FEValues of FEVectorContravariant type.
void update_gradients(FV &fe_values, const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Update shape_gradients of given FEValues object.
void fill_values_vec(FMT_UNUSED FV &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, FMT_UNUSED const FEInternalData &fe_data)
Empty method.
void update_values(FV &fe_values, const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Update shape_values of given FEValues object.
Helper class allows update values and gradients of FEValues of FEVectorPiola type.
void update_gradients(FV &fe_values, const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Update shape_gradients of given FEValues object.
void update_values(FV &fe_values, const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Update shape_values of given FEValues object.
void fill_values_vec(FMT_UNUSED FV &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, FMT_UNUSED const FEInternalData &fe_data)
Empty method.
Helper class allows update values and gradients of FEValues of FEScalar type.
void fill_values_vec(FMT_UNUSED FV &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, FMT_UNUSED const FEInternalData &fe_data)
Empty method.
void update_values(FV &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Update shape_values of given FEValues object.
void update_gradients(FV &fe_values, const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Update shape_gradients of given FEValues object.
Helper class allows update values and gradients of FEValues of FEMixedSystem type.
void update_gradients(FV &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Update shape_gradients of given FEValues object.
void fill_values_vec(FV &fe_values, const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Fill fe_values_vec of components of mixed system FEValues object.
void update_values(FV &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Update shape_values of given FEValues object.
Helper class allows update values and gradients of FEValues of FETensor type.
void update_gradients(FV &fe_values, const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Update shape_gradients of given FEValues object.
void update_values(FV &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Update shape_values of given FEValues object.
void fill_values_vec(FMT_UNUSED FV &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, FMT_UNUSED const FEInternalData &fe_data)
Empty method.
Helper class allows update values and gradients of FEValues of FEVector type.
void update_gradients(FV &fe_values, const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Update shape_gradients of given FEValues object.
void fill_values_vec(FMT_UNUSED FV &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, FMT_UNUSED const FEInternalData &fe_data)
Empty method.
void update_values(FV &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, const 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