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