Flow123d  JS_before_hm-1001-gfa0c761
fe_values.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.hh
15  * @brief Class FEValues calculates finite element data on the actual
16  * cells such as shape function values, gradients, Jacobian of
17  * the mapping from the reference cell etc.
18  * @author Jan Stebel
19  */
20 
21 #ifndef FE_VALUES_HH_
22 #define FE_VALUES_HH_
23 
24 #include <string.h> // for memcpy
25 #include <algorithm> // for swap
26 #include <new> // for operator new[]
27 #include <string> // for operator<<
28 #include <vector> // for vector
29 #include "fem/element_values.hh" // for ElementValues
30 #include "fem/fe_values_views.hh" // for FEValuesViews
31 #include "mesh/ref_element.hh" // for RefElement
32 #include "mesh/accessors.hh"
33 #include "fem/update_flags.hh" // for UpdateFlags
34 #include "tools/mixed.hh"
36 
37 class Quadrature;
38 template<unsigned int dim> class FiniteElement;
39 
40 
41 
42 
43 
44 
45 
46 
47 
48 /**
49  * @brief Calculates finite element data on the actual cell.
50  *
51  * FEValues takes care of the calculation of finite element data on
52  * the actual cell or on its side, (values of shape functions at quadrature
53  * points, gradients of shape functions, Jacobians of the mapping
54  * from the reference cell etc.).
55  * @param spacedim Dimension of the Euclidean space where the actual
56  * cell lives.
57  */
58 template<unsigned int spacedim = 3>
59 class FEValues
60 {
61 private:
62 
63  // internal structure that stores all possible views
64  // for scalar and vector-valued components of the FE
65  struct ViewsCache {
69 
70  template<unsigned int DIM>
71  void initialize(const FEValues &fv, const FiniteElement<DIM> &fe);
72  };
73 
74 public:
75 
76  /// Default constructor with postponed initialization.
77  FEValues();
78 
79 
80  /// Constructor with initialization of data structures
81  /// (see initialize() for description of parameters).
82  template<unsigned int DIM>
83  FEValues(Quadrature &_quadrature,
84  FiniteElement<DIM> &_fe,
85  UpdateFlags _flags)
86  : FEValues()
87  {
88  initialize(_quadrature, _fe, _flags);
89  }
90 
91 
92  /// Correct deallocation of objects created by 'initialize' methods.
93  ~FEValues();
94 
95 
96  /**
97  * @brief Allocates space for computed data.
98  *
99  * @param n_points Number of quadrature points.
100  * @param _fe The finite element.
101  * @param flags The update flags.
102  */
103  template<unsigned int DIM>
104  void allocate(unsigned int n_points,
105  FiniteElement<DIM> &_fe,
106  UpdateFlags flags);
107 
108  /**
109  * @brief Initialize structures and calculates cell-independent data.
110  *
111  * @param _quadrature The quadrature rule for the cell associated
112  * to given finite element or for the cell side.
113  * @param _fe The finite element.
114  * @param _flags The update flags.
115  */
116  template<unsigned int DIM>
117  void initialize(Quadrature &_quadrature,
118  FiniteElement<DIM> &_fe,
119  UpdateFlags _flags);
120 
121  /**
122  * @brief Update cell-dependent data (gradients, Jacobians etc.)
123  *
124  * @param cell The actual cell.
125  */
126  void reinit(const ElementAccessor<spacedim> &cell);
127 
128  /**
129  * @brief Update cell side-dependent FE data (values, gradients).
130  *
131  * @param cell_side Accessor to cell side.
132  */
133  void reinit(const Side &cell_side);
134 
135  /**
136  * @brief Return the value of the @p function_no-th shape function at
137  * the @p point_no-th quadrature point.
138  *
139  * @param function_no Number of the shape function.
140  * @param point_no Number of the quadrature point.
141  */
142  double shape_value(const unsigned int function_no, const unsigned int point_no);
143 
144 
145  /**
146  * @brief Return the gradient of the @p function_no-th shape function at
147  * the @p point_no-th quadrature point.
148  *
149  * @param function_no Number of the shape function.
150  * @param point_no Number of the quadrature point.
151  */
152  arma::vec::fixed<spacedim> shape_grad(const unsigned int function_no, const unsigned int point_no);
153 
154  /**
155  * @brief Return the value of the @p function_no-th shape function at
156  * the @p point_no-th quadrature point.
157  *
158  * For vectorial finite elements.
159  *
160  * @param function_no Number of the shape function.
161  * @param point_no Number of the quadrature point.
162  */
163  double shape_value_component(const unsigned int function_no,
164  const unsigned int point_no,
165  const unsigned int comp) const;
166 
167  /**
168  * @brief Return the gradient of the @p function_no-th shape function at
169  * the @p point_no-th quadrature point.
170  *
171  * For vectorial finite elements.
172  *
173  * @param function_no Number of the shape function.
174  * @param point_no Number of the quadrature point.
175  */
176  arma::vec::fixed<spacedim> shape_grad_component(const unsigned int function_no,
177  const unsigned int point_no,
178  const unsigned int comp) const;
179 
180  /**
181  * @brief Return the relative volume change of the cell (Jacobian determinant).
182  *
183  * If dim_==spacedim then the sign may be negative, otherwise the
184  * result is a positive number.
185  *
186  * @param point_no Number of the quadrature point.
187  */
188  inline double determinant(const unsigned int point_no)
189  {
190  ASSERT_LT_DBG(point_no, n_points_);
191  return elm_values->determinant(point_no);
192  }
193 
194  /**
195  * @brief Return the product of Jacobian determinant and the quadrature
196  * weight at given quadrature point.
197  *
198  * @param point_no Number of the quadrature point.
199  */
200  inline double JxW(const unsigned int point_no)
201  {
202  ASSERT_LT_DBG(point_no, n_points_);
203  // TODO: This is temporary solution to distinguish JxW on element and side_JxW on side.
204  // In future we should call the appropriate method in elm_values.
205  return (elm_values->cell().is_valid()) ? elm_values->JxW(point_no) : elm_values->side_JxW(point_no);
206  }
207 
208  /**
209  * @brief Return coordinates of the quadrature point in the actual cell system.
210  *
211  * @param point_no Number of the quadrature point.
212  */
213  inline arma::vec::fixed<spacedim> point(const unsigned int point_no)
214  {
215  ASSERT_LT_DBG(point_no, n_points_);
216  return elm_values->point(point_no);
217  }
218 
219  /**
220  * @brief Return coordinates of all quadrature points in the actual cell system.
221  *
222  */
223  inline const Armor::array &point_list() const
224  {
225  return elm_values->point_list();
226  }
227 
228 
229  /**
230  * @brief Returns the normal vector to a side at given quadrature point.
231  *
232  * @param point_no Number of the quadrature point.
233  */
234  inline arma::vec::fixed<spacedim> normal_vector(unsigned int point_no)
235  {
236  ASSERT_LT_DBG(point_no, n_points_);
237  return elm_values->normal_vector(point_no);
238  }
239 
240  /**
241  * @brief Accessor to scalar values of multicomponent FE.
242  * @param i Index of scalar component.
243  */
244  const FEValuesViews::Scalar<spacedim> &scalar_view(unsigned int i) const
245  {
247  return views_cache_.scalars[i];
248  }
249 
250  /**
251  * @brief Accessor to vector values of multicomponent FE.
252  * @param i Index of first vector component.
253  */
254  const FEValuesViews::Vector<spacedim> &vector_view(unsigned int i) const
255  {
257  return views_cache_.vectors[i];
258  }
259 
260  /**
261  * @brief Accessor to tensor values of multicomponent FE.
262  * @param i Index of first tensor component.
263  */
264  const FEValuesViews::Tensor<spacedim> &tensor_view(unsigned int i) const
265  {
267  return views_cache_.tensors[i];
268  }
269 
270  /**
271  * @brief Returns the number of quadrature points.
272  */
273  inline unsigned int n_points() const
274  { return n_points_; }
275 
276  /**
277  * @brief Returns the number of shape functions.
278  */
279  inline unsigned int n_dofs() const
280  {
281  return n_dofs_;
282  }
283 
284  /// Return dimension of reference space.
285  inline unsigned int dim() const
286  { return dim_; }
287 
288 
289 protected:
290 
291  /// Structure for storing the precomputed finite element data.
293  {
294  public:
295 
296  FEInternalData(unsigned int np, unsigned int nd);
297 
298  /// Create a new instance of FEInternalData for a FESystem component or subvector.
299  FEInternalData(const FEInternalData &fe_system_data,
300  const std::vector<unsigned int> &dof_indices,
301  unsigned int first_component_idx,
302  unsigned int ncomponents = 1);
303 
304  /**
305  * @brief Precomputed values of basis functions at the quadrature points.
306  *
307  * Dimensions: (no. of quadrature points)
308  * x (no. of dofs)
309  * x (no. of components in ref. cell)
310  */
312 
313  /**
314  * @brief Precomputed gradients of basis functions at the quadrature points.
315  *
316  * Dimensions: (no. of quadrature points)
317  * x (no. of dofs)
318  * x ((dim_ of. ref. cell)x(no. of components in ref. cell))
319  */
321 
322  /// Number of quadrature points.
323  unsigned int n_points;
324 
325  /// Number of dofs (shape functions).
326  unsigned int n_dofs;
327  };
328 
329 
330 
331 
332 
333  /// Precompute finite element data on reference element.
334  template<unsigned int DIM>
335  std::shared_ptr<FEInternalData> init_fe_data(const FiniteElement<DIM> &fe, const Quadrature &q);
336 
337  /**
338  * @brief Computes the shape function values and gradients on the actual cell
339  * and fills the FEValues structure.
340  *
341  * @param fe_data Precomputed finite element data.
342  */
344 
345  /// Compute shape functions and gradients on the actual cell for scalar FE.
347 
348  /// Compute shape functions and gradients on the actual cell for vectorial FE.
350 
351  /// Compute shape functions and gradients on the actual cell for vectorial FE.
353 
354  /// Compute shape functions and gradients on the actual cell for Raviart-Thomas FE.
356 
357  /// Compute shape functions and gradients on the actual cell for tensorial FE.
359 
360  /// Compute shape functions and gradients on the actual cell for mixed system of FE.
362 
363 
364  /// Dimension of reference space.
365  unsigned int dim_;
366 
367  /// Number of integration points.
368  unsigned int n_points_;
369 
370  /// Number of finite element dofs.
371  unsigned int n_dofs_;
372 
373  /// Type of finite element (scalar, vector, tensor).
375 
376  /// Dof indices of FESystem sub-elements.
378 
379  /// Numbers of components of FESystem sub-elements in reference space.
381 
382  /// Numbers of components of FESystem sub-elements in real space.
384 
385  /// Shape functions evaluated at the quadrature points.
387 
388  /// Gradients of shape functions evaluated at the quadrature points.
389  /// Each row of the matrix contains the gradient of one shape function.
391 
392  /// Flags that indicate which finite element quantities are to be computed.
394 
395  /// Auxiliary object for calculation of element-dependent data.
396  std::shared_ptr<ElementValues<spacedim> > elm_values;
397 
398  /// Vector of FEValues for sub-elements of FESystem.
400 
401  /// Number of components of the FE.
402  unsigned int n_components_;
403 
404  /// Auxiliary storage of FEValuesViews accessors.
405  ViewsCache views_cache_;
406 
407  /// Precomputed finite element data.
408  std::shared_ptr<FEInternalData> fe_data;
409 
410  /// Precomputed FE data (shape functions on reference element) for all sides and permuted quadrature points.
412 };
413 
414 
415 
416 
417 
418 
420  QGauss::array &quadrature,
422  UpdateFlags flags);
423 
424 
425 
426 
427 
428 
429 
430 #endif /* FE_VALUES_HH_ */
std::shared_ptr< ElementValues< spacedim > > elm_values
Auxiliary object for calculation of element-dependent data.
Definition: fe_values.hh:396
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
Definition: update_flags.hh:67
std::vector< std::vector< unsigned int > > fe_sys_dofs_
Dof indices of FESystem sub-elements.
Definition: fe_values.hh:377
const FEValuesViews::Vector< spacedim > & vector_view(unsigned int i) const
Accessor to vector values of multicomponent FE.
Definition: fe_values.hh:254
vector< FEValuesViews::Tensor< spacedim > > tensors
Definition: fe_values.hh:68
arma::vec::fixed< spacedim > shape_grad(const unsigned int function_no, const unsigned int point_no)
Return the gradient of the function_no-th shape function at the point_no-th quadrature point...
Definition: fe_values.cc:277
~FEValues()
Correct deallocation of objects created by &#39;initialize&#39; methods.
Definition: fe_values.cc:130
double determinant(const unsigned int point_no)
Return the relative volume change of the cell (Jacobian determinant).
Definition: fe_values.hh:188
FEValues(Quadrature &_quadrature, FiniteElement< DIM > &_fe, UpdateFlags _flags)
Definition: fe_values.hh:83
FEType fe_type_
Type of finite element (scalar, vector, tensor).
Definition: fe_values.hh:374
std::shared_ptr< FEInternalData > fe_data
Precomputed finite element data.
Definition: fe_values.hh:408
std::vector< std::vector< arma::vec::fixed< spacedim > > > shape_gradients
Definition: fe_values.hh:390
std::shared_ptr< FEInternalData > init_fe_data(const FiniteElement< DIM > &fe, const Quadrature &q)
Precompute finite element data on reference element.
void fill_data(const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Computes the shape function values and gradients on the actual cell and fills the FEValues structure...
Definition: fe_values.cc:515
double shape_value_component(const unsigned int function_no, const unsigned int point_no, const unsigned int comp) const
Return the value of the function_no-th shape function at the point_no-th quadrature point...
Definition: fe_values.cc:286
unsigned int n_dofs() const
Returns the number of shape functions.
Definition: fe_values.hh:279
void fill_system_data(const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for mixed system of FE.
Definition: fe_values.cc:458
Class ElementValues calculates data related to transformation of reference cell to actual cell (Jacob...
Class for computation of data on cell and side.
void fill_vec_data(const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for vectorial FE.
Definition: fe_values.cc:329
FEType
const Armor::array & point_list() const
Return coordinates of all quadrature points in the actual cell system.
Definition: fe_values.hh:223
std::vector< FEValues< 3 > > mixed_fe_values(QGauss::array &quadrature, MixedPtr< FiniteElement > fe, UpdateFlags flags)
Definition: fe_values.cc:582
void initialize(const FEValues &fv, const FiniteElement< DIM > &fe)
Definition: fe_values.cc:68
void reinit(const ElementAccessor< spacedim > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:548
void allocate(unsigned int n_points, FiniteElement< DIM > &_fe, UpdateFlags flags)
Allocates space for computed data.
Definition: fe_values.cc:183
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
vector< FEValuesViews::Scalar< spacedim > > scalars
Definition: fe_values.hh:66
double shape_value(const unsigned int function_no, const unsigned int point_no)
Return the value of the function_no-th shape function at the point_no-th quadrature point...
Definition: fe_values.cc:268
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
std::vector< unsigned int > fe_sys_n_space_components_
Numbers of components of FESystem sub-elements in real space.
Definition: fe_values.hh:383
UpdateFlags update_flags
Flags that indicate which finite element quantities are to be computed.
Definition: fe_values.hh:393
unsigned int n_points
Number of quadrature points.
Definition: fe_values.hh:323
arma::vec::fixed< spacedim > shape_grad_component(const unsigned int function_no, const unsigned int point_no, const unsigned int comp) const
Return the gradient of the function_no-th shape function at the point_no-th quadrature point...
Definition: fe_values.cc:298
unsigned int dim() const
Return dimension of reference space.
Definition: fe_values.hh:285
unsigned int n_points_
Number of integration points.
Definition: fe_values.hh:368
vector< FEValuesViews::Vector< spacedim > > vectors
Definition: fe_values.hh:67
std::vector< std::vector< arma::mat > > ref_shape_grads
Precomputed gradients of basis functions at the quadrature points.
Definition: fe_values.hh:320
const FEValuesViews::Tensor< spacedim > & tensor_view(unsigned int i) const
Accessor to tensor values of multicomponent FE.
Definition: fe_values.hh:264
unsigned int dim_
Dimension of reference space.
Definition: fe_values.hh:365
unsigned int n_dofs
Number of dofs (shape functions).
Definition: fe_values.hh:326
Structure for storing the precomputed finite element data.
Definition: fe_values.hh:292
void fill_vec_piola_data(const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for Raviart-Thomas FE.
Definition: fe_values.cc:393
std::array< QGauss, 4 > array
void fill_vec_contravariant_data(const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for vectorial FE.
Definition: fe_values.cc:361
std::vector< std::vector< shared_ptr< FEInternalData > > > side_fe_data
Precomputed FE data (shape functions on reference element) for all sides and permuted quadrature poin...
Definition: fe_values.hh:411
ViewsCache views_cache_
Auxiliary storage of FEValuesViews accessors.
Definition: fe_values.hh:405
Definitions of particular quadrature rules on simplices.
Calculates finite element data on the actual cell.
Definition: fe_values.hh:59
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc...
unsigned int n_dofs_
Number of finite element dofs.
Definition: fe_values.hh:371
std::vector< unsigned int > fe_sys_n_components_
Numbers of components of FESystem sub-elements in reference space.
Definition: fe_values.hh:380
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
double JxW(const unsigned int point_no)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point...
Definition: fe_values.hh:200
unsigned int n_components_
Number of components of the FE.
Definition: fe_values.hh:402
std::vector< FEValues< spacedim > > fe_values_vec
Vector of FEValues for sub-elements of FESystem.
Definition: fe_values.hh:399
void fill_tensor_data(const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for tensorial FE.
Definition: fe_values.cc:426
arma::vec::fixed< spacedim > point(const unsigned int point_no)
Return coordinates of the quadrature point in the actual cell system.
Definition: fe_values.hh:213
FEValues()
Default constructor with postponed initialization.
Definition: fe_values.cc:122
const FEValuesViews::Scalar< spacedim > & scalar_view(unsigned int i) const
Accessor to scalar values of multicomponent FE.
Definition: fe_values.hh:244
void fill_scalar_data(const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for scalar FE.
Definition: fe_values.cc:310
unsigned int n_points() const
Returns the number of quadrature points.
Definition: fe_values.hh:273
std::vector< std::vector< arma::vec > > ref_shape_values
Precomputed values of basis functions at the quadrature points.
Definition: fe_values.hh:311
std::vector< std::vector< double > > shape_values
Shape functions evaluated at the quadrature points.
Definition: fe_values.hh:386
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:300
arma::vec::fixed< spacedim > normal_vector(unsigned int point_no)
Returns the normal vector to a side at given quadrature point.
Definition: fe_values.hh:234