Flow123d  release_1.8.2-1603-g0109a2b
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  /**
139  * @brief Return the value of the @p function_no-th shape function at
140  * the @p point_no-th quadrature point.
141  *
142  * @param function_no Number of the shape function.
143  * @param point_no Number of the quadrature point.
144  */
145  virtual double shape_value(const unsigned int function_no, const unsigned int point_no) = 0;
146 
147  /**
148  * @brief Return the gradient of the @p function_no-th shape function at
149  * the @p point_no-th quadrature point.
150  *
151  * @param function_no Number of the shape function.
152  * @param point_no Number of the quadrature point.
153  */
154  virtual arma::vec::fixed<spacedim> shape_grad(const unsigned int function_no, const unsigned int point_no) = 0;
155 
156  /**
157  * @brief Return the product of Jacobian determinant and the quadrature
158  * weight at given quadrature point.
159  *
160  * @param point_no Number of the quadrature point.
161  */
162  virtual double JxW(const unsigned int point_no) = 0;
163 
164  /**
165  * @brief Returns the normal vector to a side at given quadrature point.
166  *
167  * @param point_no Number of the quadrature point.
168  */
169  virtual arma::vec::fixed<spacedim> normal_vector(unsigned int point_no) = 0;
170 
171  /**
172  * @brief Returns the number of shape functions.
173  */
174  virtual unsigned int n_dofs() = 0;
175 
176 };
177 
178 /**
179  * @brief Base class for FEValues and FESideValues
180  */
181 template<unsigned int dim, unsigned int spacedim>
182 class FEValuesBase : public FEValuesSpaceBase<spacedim>
183 {
184 public:
185 
186  /**
187  * @brief Default constructor
188  */
189  FEValuesBase();
190 
191 
192  /**
193  * Correct deallocation of objects created by 'initialize' methods.
194  */
195  virtual ~FEValuesBase();
196 
197 
198  /**
199  * @brief Allocates space for computed data.
200  *
201  * @param _mapping The mapping between reference and actual cell.
202  * @param _quadrature The quadrature rule.
203  * @param _fe The finite element.
204  * @param flags The update flags.
205  */
206  void allocate(Mapping<dim,spacedim> &_mapping,
207  Quadrature<dim> &_quadrature,
209  UpdateFlags flags);
210 
211  /**
212  * @brief Determine quantities to be recomputed on each cell.
213  *
214  * @param flags Flags that indicate what has to be recomputed.
215  */
216  UpdateFlags update_each(UpdateFlags flags);
217 
218  /**
219  * @brief Return the value of the @p function_no-th shape function at
220  * the @p point_no-th quadrature point.
221  *
222  * @param function_no Number of the shape function.
223  * @param point_no Number of the quadrature point.
224  */
225  inline double shape_value(const unsigned int function_no, const unsigned int point_no)
226  {
227  return data.shape_values[point_no][function_no];
228  }
229 
230 
231  /**
232  * @brief Return the gradient of the @p function_no-th shape function at
233  * the @p point_no-th quadrature point.
234  *
235  * @param function_no Number of the shape function.
236  * @param point_no Number of the quadrature point.
237  */
238  inline arma::vec::fixed<spacedim> shape_grad(const unsigned int function_no, const unsigned int point_no)
239  {
240  return trans(data.shape_gradients[point_no].row(function_no));
241 
242  }
243 
244  /**
245  * @brief Return the value of the @p function_no-th shape function at
246  * the @p point_no-th quadrature point.
247  *
248  * For vectorial finite elements.
249  *
250  * @param function_no Number of the shape function.
251  * @param point_no Number of the quadrature point.
252  */
253  inline arma::vec::fixed<spacedim> shape_vector(const unsigned int function_no, const unsigned int point_no)
254  {
255  return data.shape_vectors[point_no][function_no];
256  }
257 
258  /**
259  * @brief Return the gradient of the @p function_no-th shape function at
260  * the @p point_no-th quadrature point.
261  *
262  * For vectorial finite elements.
263  *
264  * @param function_no Number of the shape function.
265  * @param point_no Number of the quadrature point.
266  */
267  inline arma::mat::fixed<spacedim,spacedim> shape_grad_vector(const unsigned int function_no, const unsigned int point_no)
268  {
269  return data.shape_grad_vectors[point_no][function_no];
270  }
271 
272  /**
273  * @brief Return the relative volume change of the cell (Jacobian determinant).
274  *
275  * If dim==spacedim then the sign may be negative, otherwise the
276  * result is a positive number.
277  *
278  * @param point_no Number of the quadrature point.
279  */
280  inline double determinant(const unsigned int point_no)
281  {
282  return data.determinants[point_no];
283  }
284 
285  /**
286  * @brief Return the product of Jacobian determinant and the quadrature
287  * weight at given quadrature point.
288  *
289  * @param point_no Number of the quadrature point.
290  */
291  inline double JxW(const unsigned int point_no)
292  {
293  return data.JxW_values[point_no];
294  }
295 
296  /**
297  * @brief Return coordinates of the quadrature point in the actual cell system.
298  *
299  * @param point_no Number of the quadrature point.
300  */
301  inline arma::vec::fixed<spacedim> point(const unsigned int point_no)
302  {
303  return data.points[point_no];
304  }
305 
306  /**
307  * @brief Return coordinates of all quadrature points in the actual cell system.
308  *
309  */
311  {
312  return data.points;
313  }
314 
315 
316  /**
317  * @brief Returns the normal vector to a side at given quadrature point.
318  *
319  * @param point_no Number of the quadrature point.
320  */
321  inline arma::vec::fixed<spacedim> normal_vector(unsigned int point_no)
322  {
323  return data.normal_vectors[point_no];
324  }
325 
326  /**
327  * @brief Returns the number of quadrature points.
328  */
329  inline unsigned int n_points()
330  {
331  return quadrature->size();
332  }
333 
334  /**
335  * @brief Returns the number of shape functions.
336  */
337  inline unsigned int n_dofs()
338  {
339  return fe->n_dofs();
340  }
341 
342 
343  /**
344  * @brief Returns the quadrature in use.
345  */
347  {
348  return quadrature;
349  }
350 
351  /**
352  * @brief Returns the finite element in use.
353  */
355  {
356  return fe;
357  }
358 
359  /**
360  * @brief Returns the mapping in use.
361  */
363  {
364  return mapping;
365  }
366 
367 protected:
368 
369  /**
370  * @brief The mapping from the reference cell to the actual cell.
371  */
373 
374  /**
375  * @brief The quadrature rule used to calculate integrals.
376  */
378 
379  /**
380  * @brief The used finite element.
381  */
383 
384  /**
385  * @brief Precomputed mapping data.
386  */
388 
389  /**
390  * @brief Precomputed finite element data.
391  */
393 
394  /**
395  * @brief Data computed by the mapping and finite element.
396  */
398 };
399 
400 
401 
402 
403 /**
404  * @brief Calculates finite element data on the actual cell.
405  *
406  * FEValues takes care of the calculation of finite element data on
407  * the actual cell such as values of shape functions at quadrature
408  * points, gradients of shape functions, Jacobians of the mapping
409  * from the reference cell etc.
410  * @param dim Dimension of the reference cell.
411  * @param spacedim Dimension of the Euclidean space where the actual
412  * cell lives.
413  */
414 template<unsigned int dim, unsigned int spacedim>
415 class FEValues : public FEValuesBase<dim,spacedim>
416 {
417 public:
418 
419  /**
420  * @brief Constructor.
421  *
422  * Initializes structures and calculates
423  * cell-independent data.
424  *
425  * @param _mapping The mapping between the reference and actual cell.
426  * @param _quadrature The quadrature rule.
427  * @param _fe The finite element.
428  * @param _flags The update flags.
429  */
431  Quadrature<dim> &_quadrature,
433  UpdateFlags _flags);
434 
435  /**
436  * @brief Update cell-dependent data (gradients, Jacobians etc.)
437  *
438  * @param cell The actual cell.
439  */
440  void reinit(ElementFullIter &cell);
441 
442 
443 };
444 
445 
446 
447 
448 /**
449  * @brief Calculates finite element data on a side.
450  *
451  * FESideValues takes care of the calculation of finite element data on
452  * a side of the actual cell such as values of shape functions at quadrature
453  * points, gradients of shape functions, Jacobians of the mapping
454  * from the reference cell etc.
455  * @param dim Dimension of the reference cell.
456  * @param spacedim Dimension of the Euclidean space where the actual
457  * cell lives.
458  */
459 template<unsigned int dim, unsigned int spacedim>
460 class FESideValues : public FEValuesBase<dim,spacedim>
461 {
462 
463 public:
464 
465  /**
466  * @brief Constructor.
467  *
468  * Initializes structures and calculates
469  * cell-independent data.
470  *
471  * @param _mapping The mapping between the reference and actual cell.
472  * @param _sub_quadrature The quadrature rule on the side.
473  * @param _fe The finite element.
474  * @param flags The update flags.
475  */
477  Quadrature<dim-1> &_sub_quadrature,
479  UpdateFlags flags);
480 
481  /// Destructor.
482  virtual ~FESideValues();
483 
484  /**
485  * @brief Update cell-dependent data (gradients, Jacobians etc.)
486  *
487  * @param cell The actual cell.
488  * @param sid Number of the side of the cell.
489  */
490  void reinit(ElementFullIter &cell,
491  unsigned int sid);
492 
493 
494 
495 private:
496 
497  /**
498  * @brief Quadrature for the integration on the element sides.
499  */
500  const Quadrature<dim-1> *sub_quadrature;
501 
503 
505 
507 
508 };
509 
510 
511 
512 
513 
514 #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:291
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:337
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:387
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:354
FEValuesData< dim, spacedim > data
Data computed by the mapping and finite element.
Definition: fe_values.hh:397
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:392
unsigned int n_points()
Returns the number of quadrature points.
Definition: fe_values.hh:329
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:321
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:382
std::vector< arma::mat > shape_gradients
Gradients of shape functions evaluated at the quadrature points.
Definition: fe_values.hh:99
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:500
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:253
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:267
Quadrature< dim > * quadrature
The quadrature rule used to calculate integrals.
Definition: fe_values.hh:377
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:225
Mapping< dim, spacedim > * get_mapping() const
Returns the mapping in use.
Definition: fe_values.hh:362
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:346
Mapping< dim, spacedim > * mapping
The mapping from the reference cell to the actual cell.
Definition: fe_values.hh:372
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:415
vector< arma::vec::fixed< spacedim > > & point_list()
Return coordinates of all quadrature points in the actual cell system.
Definition: fe_values.hh:310
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:238
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:280
Calculates finite element data on a side.
Definition: fe_values.hh:460
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:301