Flow123d  release_2.2.0-914-gf1a3a4f
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 <armadillo>
25 #include <vector>
26 #include "fem/update_flags.hh"
27 #include "fem/fe_values_views.hh"
28 #include "mesh/ref_element.hh"
29 #include "mesh/mesh_types.hh"
30 
31 class DOFHandlerBase;
32 template<unsigned int dim> class Quadrature;
33 template<unsigned int dim, unsigned int spacedim> class FiniteElement;
34 template<unsigned int dim, unsigned int spacedim> class FEValuesBase;
35 template<unsigned int dim, unsigned int spacedim> class Mapping;
36 
37 struct FEInternalData;
38 struct MappingInternalData;
39 class SideIter;
40 
41 
42 
43 /**
44  * @brief Class FEValuesData holds the arrays of data computed by
45  * Mapping and FiniteElement.
46  */
47 template<unsigned int dim, unsigned int spacedim>
49 {
50 public:
51 
52  /**
53  * @brief Resize the data arrays.
54  * @param size Number of quadrature points.
55  * @param flags Update flags to be stores.
56  * @param n_comp Number of components of shape values.
57  */
58  void allocate(unsigned int size, UpdateFlags flags, unsigned n_comp);
59 
60 
61 
62  /**
63  * @brief Transformed quadrature weights.
64  *
65  * Values at quadrature points of the product of the Jacobian
66  * determinant of the mapping and the weight at the particular
67  * quadrature point.
68  */
70 
71  /**
72  * @brief Jacobians of the mapping at the quadrature points.
73  */
75 
76  /**
77  * @brief Determinants of Jacobians at quadrature points.
78  */
80 
81  /**
82  * @brief Inverse Jacobians at the quadrature points.
83  */
85 
86  /**
87  * @brief Coordinates of quadrature points in the actual cell coordinate system.
88  */
90 
91  /**
92  * @brief Shape functions evaluated at the quadrature points.
93  */
95 
96  /**
97  * @brief Gradients of shape functions evaluated at the quadrature points.
98  *
99  * Each row of the matrix contains the gradient of one shape function.
100  */
102 
103 // /**
104 // * @brief Shape functions (for vectorial finite elements) evaluated at
105 // * quadrature points.
106 // */
107 // std::vector<std::vector<arma::vec::fixed<spacedim> > > shape_vectors;
108 //
109 // /**
110 // * @brief Gradients of shape functions (for vectorial finite elements).
111 // */
112 // std::vector<std::vector<arma::mat::fixed<spacedim,spacedim> > > shape_grad_vectors;
113 
114  /**
115  * @brief Normal vectors to the element at the quadrature points lying
116  * on a side.
117  */
119 
120  /**
121  * @brief Flags that indicate which finite element quantities are to be computed.
122  */
124 
125  /**
126  * @brief Iterator to the last reinit-ed cell.
127  */
129 
130 };
131 
132 
133 /**
134  * @brief Abstract base class with certain methods independent of the template parameter @p dim.
135  */
136 template<unsigned int spacedim>
138 {
139 public:
140  virtual ~FEValuesSpaceBase() {}
141  /**
142  * @brief Return the value of the @p function_no-th shape function at
143  * the @p point_no-th quadrature point.
144  *
145  * @param function_no Number of the shape function.
146  * @param point_no Number of the quadrature point.
147  */
148  virtual double shape_value(const unsigned int function_no, const unsigned int point_no) = 0;
149 
150  /**
151  * @brief Return the gradient of the @p function_no-th shape function at
152  * the @p point_no-th quadrature point.
153  *
154  * @param function_no Number of the shape function.
155  * @param point_no Number of the quadrature point.
156  */
157  virtual arma::vec::fixed<spacedim> shape_grad(const unsigned int function_no, const unsigned int point_no) = 0;
158 
159  /**
160  * @brief Return the product of Jacobian determinant and the quadrature
161  * weight at given quadrature point.
162  *
163  * @param point_no Number of the quadrature point.
164  */
165  virtual double JxW(const unsigned int point_no) = 0;
166 
167  /**
168  * @brief Returns the normal vector to a side at given quadrature point.
169  *
170  * @param point_no Number of the quadrature point.
171  */
172  virtual arma::vec::fixed<spacedim> normal_vector(unsigned int point_no) = 0;
173 
174  /**
175  * @brief Returns the number of shape functions.
176  */
177  virtual unsigned int n_dofs() = 0;
178 
179 };
180 
181 /**
182  * @brief Base class for FEValues and FESideValues
183  */
184 template<unsigned int dim, unsigned int spacedim>
185 class FEValuesBase : public FEValuesSpaceBase<spacedim>
186 {
187 private:
188 
189  // internal structure that stores all possible views
190  // for scalar and vector-valued components of the FE
191  struct ViewsCache {
194 
195  void resize(FEValuesBase &fv, unsigned int size);
196  };
197 
198 public:
199 
200  /**
201  * @brief Default constructor
202  */
203  FEValuesBase();
204 
205 
206  /**
207  * Correct deallocation of objects created by 'initialize' methods.
208  */
209  virtual ~FEValuesBase();
210 
211 
212  /**
213  * @brief Allocates space for computed data.
214  *
215  * @param _mapping The mapping between reference and actual cell.
216  * @param _quadrature The quadrature rule.
217  * @param _fe The finite element.
218  * @param flags The update flags.
219  */
220  void allocate(Mapping<dim,spacedim> &_mapping,
221  Quadrature<dim> &_quadrature,
223  UpdateFlags flags);
224 
225  /**
226  * @brief Determine quantities to be recomputed on each cell.
227  *
228  * @param flags Flags that indicate what has to be recomputed.
229  */
230  UpdateFlags update_each(UpdateFlags flags);
231 
232  /**
233  * @brief Return the value of the @p function_no-th shape function at
234  * the @p point_no-th quadrature point.
235  *
236  * @param function_no Number of the shape function.
237  * @param point_no Number of the quadrature point.
238  */
239  double shape_value(const unsigned int function_no, const unsigned int point_no);
240 
241 
242  /**
243  * @brief Return the gradient of the @p function_no-th shape function at
244  * the @p point_no-th quadrature point.
245  *
246  * @param function_no Number of the shape function.
247  * @param point_no Number of the quadrature point.
248  */
249  arma::vec::fixed<spacedim> shape_grad(const unsigned int function_no, const unsigned int point_no);
250 
251  /**
252  * @brief Return the value of the @p function_no-th shape function at
253  * the @p point_no-th quadrature point.
254  *
255  * For vectorial finite elements.
256  *
257  * @param function_no Number of the shape function.
258  * @param point_no Number of the quadrature point.
259  */
260  double shape_value_component(const unsigned int function_no,
261  const unsigned int point_no,
262  const unsigned int comp) const;
263 
264  /**
265  * @brief Return the gradient of the @p function_no-th shape function at
266  * the @p point_no-th quadrature point.
267  *
268  * For vectorial finite elements.
269  *
270  * @param function_no Number of the shape function.
271  * @param point_no Number of the quadrature point.
272  */
273  arma::vec::fixed<spacedim> shape_grad_component(const unsigned int function_no,
274  const unsigned int point_no,
275  const unsigned int comp) const;
276 
277  /**
278  * @brief Return the relative volume change of the cell (Jacobian determinant).
279  *
280  * If dim==spacedim then the sign may be negative, otherwise the
281  * result is a positive number.
282  *
283  * @param point_no Number of the quadrature point.
284  */
285  inline double determinant(const unsigned int point_no)
286  {
287  ASSERT_LT_DBG(point_no, quadrature->size());
288  return data.determinants[point_no];
289  }
290 
291  /**
292  * @brief Return the product of Jacobian determinant and the quadrature
293  * weight at given quadrature point.
294  *
295  * @param point_no Number of the quadrature point.
296  */
297  inline double JxW(const unsigned int point_no)
298  {
299  ASSERT_LT_DBG(point_no, quadrature->size());
300  return data.JxW_values[point_no];
301  }
302 
303  /**
304  * @brief Return coordinates of the quadrature point in the actual cell system.
305  *
306  * @param point_no Number of the quadrature point.
307  */
308  inline arma::vec::fixed<spacedim> point(const unsigned int point_no)
309  {
310  ASSERT_LT_DBG(point_no, quadrature->size());
311  return data.points[point_no];
312  }
313 
314  /**
315  * @brief Return coordinates of all quadrature points in the actual cell system.
316  *
317  */
319  {
320  return data.points;
321  }
322 
323 
324  /**
325  * @brief Returns the normal vector to a side at given quadrature point.
326  *
327  * @param point_no Number of the quadrature point.
328  */
329  inline arma::vec::fixed<spacedim> normal_vector(unsigned int point_no)
330  {
331  ASSERT_LT_DBG(point_no, quadrature->size());
332  return data.normal_vectors[point_no];
333  }
334 
335  /**
336  * @brief Accessor to scalar values of multicomponent FE.
337  * @param i Index of scalar component.
338  */
340  {
341  ASSERT_LT_DBG(i, views_cache_.scalars.size());
342  return views_cache_.scalars[i];
343  }
344 
345  /**
346  * @brief Accessor to vector values of multicomponent FE.
347  * @param i Index of first vector component.
348  */
350  {
351  ASSERT_LT_DBG(i, views_cache_.vectors.size());
352  return views_cache_.vectors[i];
353  }
354 
355  /**
356  * @brief Returns the number of quadrature points.
357  */
358  inline unsigned int n_points()
359  {
360  return quadrature->size();
361  }
362 
363  /**
364  * @brief Returns the number of shape functions.
365  */
366  inline unsigned int n_dofs()
367  {
368  return fe->n_dofs();
369  }
370 
371 
372  /**
373  * @brief Returns the quadrature in use.
374  */
376  {
377  return quadrature;
378  }
379 
380  /**
381  * @brief Returns the finite element in use.
382  */
384  {
385  return fe;
386  }
387 
388  /**
389  * @brief Returns the mapping in use.
390  */
392  {
393  return mapping;
394  }
395 
396 protected:
397 
398  /**
399  * @brief The mapping from the reference cell to the actual cell.
400  */
402 
403  /**
404  * @brief The quadrature rule used to calculate integrals.
405  */
407 
408  /**
409  * @brief The used finite element.
410  */
412 
413  /**
414  * @brief Precomputed mapping data.
415  */
417 
418  /**
419  * @brief Precomputed finite element data.
420  */
422 
423  /**
424  * @brief Data computed by the mapping and finite element.
425  */
427 
428  /// Number of components of the FE.
429  unsigned int n_components_;
430 
431  /// Auxiliary storage of FEValuesViews accessors.
433 };
434 
435 
436 
437 
438 /**
439  * @brief Calculates finite element data on the actual cell.
440  *
441  * FEValues takes care of the calculation of finite element data on
442  * the actual cell such as values of shape functions at quadrature
443  * points, gradients of shape functions, Jacobians of the mapping
444  * from the reference cell etc.
445  * @param dim Dimension of the reference cell.
446  * @param spacedim Dimension of the Euclidean space where the actual
447  * cell lives.
448  */
449 template<unsigned int dim, unsigned int spacedim>
450 class FEValues : public FEValuesBase<dim,spacedim>
451 {
452 public:
453 
454  /**
455  * @brief Constructor.
456  *
457  * Initializes structures and calculates
458  * cell-independent data.
459  *
460  * @param _mapping The mapping between the reference and actual cell.
461  * @param _quadrature The quadrature rule.
462  * @param _fe The finite element.
463  * @param _flags The update flags.
464  */
466  Quadrature<dim> &_quadrature,
468  UpdateFlags _flags);
469 
470  /**
471  * @brief Update cell-dependent data (gradients, Jacobians etc.)
472  *
473  * @param cell The actual cell.
474  */
475  void reinit(ElementFullIter &cell);
476 
477 
478 };
479 
480 
481 
482 
483 /**
484  * @brief Calculates finite element data on a side.
485  *
486  * FESideValues takes care of the calculation of finite element data on
487  * a side of the actual cell such as values of shape functions at quadrature
488  * points, gradients of shape functions, Jacobians of the mapping
489  * from the reference cell etc.
490  * @param dim Dimension of the reference cell.
491  * @param spacedim Dimension of the Euclidean space where the actual
492  * cell lives.
493  */
494 template<unsigned int dim, unsigned int spacedim>
495 class FESideValues : public FEValuesBase<dim,spacedim>
496 {
497 
498 public:
499 
500  /**
501  * @brief Constructor.
502  *
503  * Initializes structures and calculates
504  * cell-independent data.
505  *
506  * @param _mapping The mapping between the reference and actual cell.
507  * @param _sub_quadrature The quadrature rule on the side.
508  * @param _fe The finite element.
509  * @param flags The update flags.
510  */
512  Quadrature<dim-1> &_sub_quadrature,
514  UpdateFlags flags);
515 
516  /// Destructor.
517  virtual ~FESideValues();
518 
519  /**
520  * @brief Update cell-dependent data (gradients, Jacobians etc.)
521  *
522  * @param cell The actual cell.
523  * @param sid Number of the side of the cell.
524  */
525  void reinit(ElementFullIter &cell,
526  unsigned int sid);
527 
528 
529 
530 private:
531 
532  /**
533  * @brief Quadrature for the integration on the element sides.
534  */
535  const Quadrature<dim-1> *sub_quadrature;
536 
538 
540 
542 
543 };
544 
545 
546 
547 
548 
549 #endif /* FE_VALUES_HH_ */
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:297
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
Definition: update_flags.hh:67
unsigned int n_dofs()
Returns the number of shape functions.
Definition: fe_values.hh:366
std::vector< arma::mat::fixed< dim, spacedim > > inverse_jacobians
Inverse Jacobians at the quadrature points.
Definition: fe_values.hh:84
MappingInternalData * mapping_data
Precomputed mapping data.
Definition: fe_values.hh:416
UpdateFlags update_flags
Flags that indicate which finite element quantities are to be computed.
Definition: fe_values.hh:123
FiniteElement< dim, spacedim > * get_fe() const
Returns the finite element in use.
Definition: fe_values.hh:383
FEValuesData< dim, spacedim > data
Data computed by the mapping and finite element.
Definition: fe_values.hh:426
std::vector< double > JxW_values
Transformed quadrature weights.
Definition: fe_values.hh:69
FEInternalData * fe_data
Precomputed finite element data.
Definition: fe_values.hh:421
unsigned int n_points()
Returns the number of quadrature points.
Definition: fe_values.hh:358
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
void allocate(unsigned int size, UpdateFlags flags, unsigned n_comp)
Resize the data arrays.
Definition: fe_values.cc:42
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: fe_values.hh:32
vector< FEValuesViews::Vector< dim, spacedim > > vectors
Definition: fe_values.hh:193
unsigned int n_components_
Number of components of the FE.
Definition: fe_values.hh:429
std::vector< arma::vec::fixed< spacedim > > normal_vectors
Shape functions (for vectorial finite elements) evaluated at quadrature points.
Definition: fe_values.hh:118
vector< FEValuesViews::Scalar< dim, spacedim > > scalars
Definition: fe_values.hh:192
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:329
std::vector< arma::mat::fixed< spacedim, dim > > jacobians
Jacobians of the mapping at the quadrature points.
Definition: fe_values.hh:74
std::vector< arma::vec::fixed< spacedim > > points
Coordinates of quadrature points in the actual cell coordinate system.
Definition: fe_values.hh:89
Abstract class for the mapping between reference and actual cell.
Definition: fe_values.hh:35
FiniteElement< dim, spacedim > * fe
The used finite element.
Definition: fe_values.hh:411
virtual ~FEValuesSpaceBase()
Definition: fe_values.hh:140
ElementFullIter * present_cell
Iterator to the last reinit-ed cell.
Definition: fe_values.hh:128
const Quadrature< dim-1 > * sub_quadrature
Quadrature for the integration on the element sides.
Definition: fe_values.hh:535
const FEValuesViews::Vector< dim, spacedim > & vector_view(unsigned int i) const
Accessor to vector values of multicomponent FE.
Definition: fe_values.hh:349
Quadrature< dim > * quadrature
The quadrature rule used to calculate integrals.
Definition: fe_values.hh:406
Mapping< dim, spacedim > * get_mapping() const
Returns the mapping in use.
Definition: fe_values.hh:391
Class FEValuesData holds the arrays of data computed by Mapping and FiniteElement.
Definition: fe_values.hh:48
Quadrature< dim > * get_quadrature() const
Returns the quadrature in use.
Definition: fe_values.hh:375
ViewsCache views_cache_
Auxiliary storage of FEValuesViews accessors.
Definition: fe_values.hh:432
Mapping< dim, spacedim > * mapping
The mapping from the reference cell to the actual cell.
Definition: fe_values.hh:401
std::vector< std::vector< double > > shape_values
Shape functions evaluated at the quadrature points.
Definition: fe_values.hh:94
std::vector< double > determinants
Determinants of Jacobians at quadrature points.
Definition: fe_values.hh:79
Calculates finite element data on the actual cell.
Definition: fe_values.hh:450
vector< arma::vec::fixed< spacedim > > & point_list()
Return coordinates of all quadrature points in the actual cell system.
Definition: fe_values.hh:318
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc...
Structure for storing the precomputed finite element data.
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Definition: dofhandler.hh:30
Base class for FEValues and FESideValues.
Definition: fe_values.hh:34
Abstract base class with certain methods independent of the template parameter dim.
Definition: fe_values.hh:137
std::vector< std::vector< arma::vec::fixed< spacedim > > > shape_gradients
Gradients of shape functions evaluated at the quadrature points.
Definition: fe_values.hh:101
Mapping data that can be precomputed on the actual cell.
Definition: mapping.hh:99
const FEValuesViews::Scalar< dim, spacedim > & scalar_view(unsigned int i) const
Accessor to scalar values of multicomponent FE.
Definition: fe_values.hh:339
double determinant(const unsigned int point_no)
Return the relative volume change of the cell (Jacobian determinant).
Definition: fe_values.hh:285
Calculates finite element data on a side.
Definition: fe_values.hh:495
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:299
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:308