Flow123d  jenkins-Flow123d-linux-release-multijob-282
fe_values.hh
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2007 Technical University of Liberec. All rights reserved.
4  *
5  * Please make a following refer to Flow123d on your project site if you use the program for any purpose,
6  * especially for academic research:
7  * Flow123d, Research Centre: Advanced Remedial Technologies, Technical University of Liberec, Czech Republic
8  *
9  * This program is free software; you can redistribute it and/or modify it under the terms
10  * of the GNU General Public License version 3 as published by the Free Software Foundation.
11  *
12  * This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
13  * without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
14  * See the GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License along with this program; if not,
17  * write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 021110-1307, USA.
18  *
19  *
20  * $Id$
21  * $Revision$
22  * $LastChangedBy$
23  * $LastChangedDate$
24  *
25  * @file
26  * @brief Class FEValues calculates finite element data on the actual
27  * cells such as shape function values, gradients, Jacobian of
28  * the mapping from the reference cell etc.
29  * @author Jan Stebel
30  */
31 
32 #ifndef FE_VALUES_HH_
33 #define FE_VALUES_HH_
34 
35 #include <armadillo>
36 #include <vector>
37 #include "fem/update_flags.hh"
38 #include "mesh/ref_element.hh"
39 #include "mesh/mesh_types.hh"
40 
41 class DOFHandlerBase;
42 template<unsigned int dim> class Quadrature;
43 template<unsigned int dim, unsigned int spacedim> class FEValuesBase;
44 template<unsigned int dim, unsigned int spacedim> class Mapping;
45 
46 struct FEInternalData;
47 struct MappingInternalData;
48 class SideIter;
49 
50 
51 
52 /**
53  * @brief Class FEValuesData holds the arrays of data computed by
54  * Mapping and FiniteElement.
55  */
56 template<unsigned int dim, unsigned int spacedim>
58 {
59 public:
60 
61  /**
62  * @brief Resize the data arrays.
63  * @param size Number of quadrature points.
64  * @param flags Update flags to be stores.
65  * @param is_scalar If true, the structures for scalar values are allocated. Otherwise the vectorial structures are used.
66  */
67  void allocate(unsigned int size, UpdateFlags flags, bool is_scalar = true);
68 
69 
70 
71  /**
72  * @brief Transformed quadrature weights.
73  *
74  * Values at quadrature points of the product of the Jacobian
75  * determinant of the mapping and the weight at the particular
76  * quadrature point.
77  */
79 
80  /**
81  * @brief Jacobians of the mapping at the quadrature points.
82  */
84 
85  /**
86  * @brief Determinants of Jacobians at quadrature points.
87  */
89 
90  /**
91  * @brief Inverse Jacobians at the quadrature points.
92  */
94 
95  /**
96  * @brief Coordinates of quadrature points in the actual cell coordinate system.
97  */
99 
100  /**
101  * @brief Shape functions evaluated at the quadrature points.
102  */
104 
105  /**
106  * @brief Gradients of shape functions evaluated at the quadrature points.
107  *
108  * Each row of the matrix contains the gradient of one shape function.
109  */
111 
112  /**
113  * @brief Shape functions (for vectorial finite elements) evaluated at
114  * quadrature points.
115  */
117 
118  /**
119  * @brief Gradients of shape functions (for vectorial finite elements).
120  */
122 
123  /**
124  * @brief Normal vectors to the element at the quadrature points lying
125  * on a side.
126  */
128 
129  /**
130  * @brief Flags that indicate which finite element quantities are to be computed.
131  */
133 
134  /**
135  * @brief Iterator to the last reinit-ed cell.
136  */
138 
139 };
140 
141 
142 /**
143  * @brief Abstract base class with certain methods independent of the template parameter @p dim.
144  */
145 template<unsigned int spacedim>
147 {
148 public:
149  /**
150  * @brief Return the value of the @p function_no-th shape function at
151  * the @p point_no-th quadrature point.
152  *
153  * @param function_no Number of the shape function.
154  * @param point_no Number of the quadrature point.
155  */
156  virtual const double shape_value(const unsigned int function_no, const unsigned int point_no) = 0;
157 
158  /**
159  * @brief Return the gradient of the @p function_no-th shape function at
160  * the @p point_no-th quadrature point.
161  *
162  * @param function_no Number of the shape function.
163  * @param point_no Number of the quadrature point.
164  */
165  virtual const arma::vec::fixed<spacedim> shape_grad(const unsigned int function_no, const unsigned int point_no) = 0;
166 
167  /**
168  * @brief Return the product of Jacobian determinant and the quadrature
169  * weight at given quadrature point.
170  *
171  * @param point_no Number of the quadrature point.
172  */
173  virtual const double JxW(const unsigned int point_no) = 0;
174 
175  /**
176  * @brief Returns the normal vector to a side at given quadrature point.
177  *
178  * @param point_no Number of the quadrature point.
179  */
180  virtual const arma::vec::fixed<spacedim> normal_vector(unsigned int point_no) = 0;
181 
182  /**
183  * @brief Returns the number of shape functions.
184  */
185  virtual const unsigned int n_dofs() = 0;
186 
187 };
188 
189 /**
190  * @brief Base class for FEValues and FESideValues
191  */
192 template<unsigned int dim, unsigned int spacedim>
193 class FEValuesBase : public FEValuesSpaceBase<spacedim>
194 {
195 public:
196 
197  /**
198  * @brief Default constructor
199  */
200  FEValuesBase();
201 
202 
203  /**
204  * Correct deallocation of objects created by 'initialize' methods.
205  */
206  virtual ~FEValuesBase();
207 
208 
209  /**
210  * @brief Allocates space for computed data.
211  *
212  * @param _mapping The mapping between reference and actual cell.
213  * @param _quadrature The quadrature rule.
214  * @param _fe The finite element.
215  * @param flags The update flags.
216  */
217  void allocate(Mapping<dim,spacedim> &_mapping,
218  Quadrature<dim> &_quadrature,
220  UpdateFlags flags);
221 
222  /**
223  * @brief Determine quantities to be recomputed on each cell.
224  *
225  * @param flags Flags that indicate what has to be recomputed.
226  */
228 
229  /**
230  * @brief Return the value of the @p function_no-th shape function at
231  * the @p point_no-th quadrature point.
232  *
233  * @param function_no Number of the shape function.
234  * @param point_no Number of the quadrature point.
235  */
236  const double shape_value(const unsigned int function_no, const unsigned int point_no);
237 
238  /**
239  * @brief Return the gradient of the @p function_no-th shape function at
240  * the @p point_no-th quadrature point.
241  *
242  * @param function_no Number of the shape function.
243  * @param point_no Number of the quadrature point.
244  */
245  const arma::vec::fixed<spacedim> shape_grad(const unsigned int function_no, const unsigned int point_no);
246 
247  /**
248  * @brief Return the value of the @p function_no-th shape function at
249  * the @p point_no-th quadrature point.
250  *
251  * For vectorial finite elements.
252  *
253  * @param function_no Number of the shape function.
254  * @param point_no Number of the quadrature point.
255  */
256  const arma::vec::fixed<spacedim> shape_vector(const unsigned int function_no, const unsigned int point_no);
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  const arma::mat::fixed<spacedim,spacedim> shape_grad_vector(const unsigned int function_no, const unsigned int point_no);
268 
269  /**
270  * @brief Return the relative volume change of the cell (Jacobian determinant).
271  *
272  * If dim==spacedim then the sign may be negative, otherwise the
273  * result is a positive number.
274  *
275  * @param point_no Number of the quadrature point.
276  */
277  const double determinant(const unsigned int point_no);
278 
279  /**
280  * @brief Return the product of Jacobian determinant and the quadrature
281  * weight at given quadrature point.
282  *
283  * @param point_no Number of the quadrature point.
284  */
285  const double JxW(const unsigned int point_no);
286 
287  /**
288  * @brief Return coordinates of the quadrature point in the actual cell system.
289  *
290  * @param point_no Number of the quadrature point.
291  */
292  const arma::vec::fixed<spacedim> point(const unsigned int point_no);
293 
294  /**
295  * @brief Return coordinates of all quadrature points in the actual cell system.
296  *
297  */
299 
300  /**
301  * @brief Returns the normal vector to a side at given quadrature point.
302  *
303  * @param point_no Number of the quadrature point.
304  */
305  const arma::vec::fixed<spacedim> normal_vector(unsigned int point_no);
306 
307  /**
308  * @brief Returns the number of quadrature points.
309  */
310  const unsigned int n_points();
311 
312  /**
313  * @brief Returns the number of shape functions.
314  */
315  const unsigned int n_dofs();
316 
317  /**
318  * @brief Returns the quadrature in use.
319  */
320  const Quadrature<dim> * get_quadrature() const;
321 
322 protected:
323 
324  /**
325  * @brief The mapping from the reference cell to the actual cell.
326  */
328 
329  /**
330  * @brief The quadrature rule used to calculate integrals.
331  */
333 
334  /**
335  * @brief The used finite element.
336  */
338 
339  /**
340  * @brief Precomputed mapping data.
341  */
343 
344  /**
345  * @brief Precomputed finite element data.
346  */
348 
349  /**
350  * @brief Data computed by the mapping and finite element.
351  */
353 };
354 
355 
356 
357 
358 /**
359  * @brief Calculates finite element data on the actual cell.
360  *
361  * FEValues takes care of the calculation of finite element data on
362  * the actual cell such as values of shape functions at quadrature
363  * points, gradients of shape functions, Jacobians of the mapping
364  * from the reference cell etc.
365  * @param dim Dimension of the reference cell.
366  * @param spacedim Dimension of the Euclidean space where the actual
367  * cell lives.
368  */
369 template<unsigned int dim, unsigned int spacedim>
370 class FEValues : public FEValuesBase<dim,spacedim>
371 {
372 public:
373 
374  /**
375  * @brief Constructor.
376  *
377  * Initializes structures and calculates
378  * cell-independent data.
379  *
380  * @param _mapping The mapping between the reference and actual cell.
381  * @param _quadrature The quadrature rule.
382  * @param _fe The finite element.
383  * @param _flags The update flags.
384  */
386  Quadrature<dim> &_quadrature,
388  UpdateFlags _flags);
389 
390  /**
391  * @brief Update cell-dependent data (gradients, Jacobians etc.)
392  *
393  * @param cell The actual cell.
394  */
395  void reinit(ElementFullIter &cell);
396 
397 
398 };
399 
400 
401 
402 
403 /**
404  * @brief Calculates finite element data on a side.
405  *
406  * FESideValues takes care of the calculation of finite element data on
407  * a side of 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 FESideValues : public FEValuesBase<dim,spacedim>
416 {
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 _sub_quadrature The quadrature rule on the side.
428  * @param _fe The finite element.
429  * @param flags The update flags.
430  */
432  Quadrature<dim-1> &_sub_quadrature,
434  UpdateFlags flags);
435 
436  /// Destructor.
437  virtual ~FESideValues();
438 
439  /**
440  * @brief Update cell-dependent data (gradients, Jacobians etc.)
441  *
442  * @param cell The actual cell.
443  * @param sid Number of the side of the cell.
444  */
445  void reinit(ElementFullIter &cell,
446  unsigned int sid);
447 
448 
449 
450 private:
451 
452  /**
453  * @brief Quadrature for the integration on the element sides.
454  */
455  const Quadrature<dim-1> *sub_quadrature;
456 
458 
460 
462 
463 };
464 
465 
466 
467 
468 
469 #endif /* FE_VALUES_HH_ */
virtual const double JxW(const unsigned int point_no)=0
Return the product of Jacobian determinant and the quadrature weight at given quadrature point...
std::vector< arma::vec > shape_values
Shape functions evaluated at the quadrature points.
Definition: fe_values.hh:103
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
Definition: update_flags.hh:78
const Quadrature< dim > * get_quadrature() const
Returns the quadrature in use.
Definition: fe_values.cc:214
const arma::vec::fixed< spacedim > point(const unsigned int point_no)
Return coordinates of the quadrature point in the actual cell system.
Definition: fe_values.cc:184
std::vector< arma::mat::fixed< dim, spacedim > > inverse_jacobians
Inverse Jacobians at the quadrature points.
Definition: fe_values.hh:93
MappingInternalData * mapping_data
Precomputed mapping data.
Definition: fe_values.hh:342
UpdateFlags update_flags
Flags that indicate which finite element quantities are to be computed.
Definition: fe_values.hh:132
UpdateFlags update_each(UpdateFlags flags)
Determine quantities to be recomputed on each cell.
Definition: fe_values.cc:140
const double JxW(const unsigned int point_no)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point...
Definition: fe_values.cc:178
const arma::vec::fixed< spacedim > normal_vector(unsigned int point_no)
Returns the normal vector to a side at given quadrature point.
Definition: fe_values.cc:196
virtual ~FEValuesBase()
Definition: fe_values.cc:115
const double determinant(const unsigned int point_no)
Return the relative volume change of the cell (Jacobian determinant).
Definition: fe_values.cc:172
FEValuesBase()
Default constructor.
Definition: fe_values.cc:107
FEValuesData< dim, spacedim > data
Data computed by the mapping and finite element.
Definition: fe_values.hh:352
std::vector< std::vector< arma::mat::fixed< spacedim, spacedim > > > shape_grad_vectors
Gradients of shape functions (for vectorial finite elements).
Definition: fe_values.hh:121
std::vector< double > JxW_values
Transformed quadrature weights.
Definition: fe_values.hh:78
FEInternalData * fe_data
Precomputed finite element data.
Definition: fe_values.hh:347
const 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:154
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:42
virtual ~FESideValues()
Destructor.
Definition: fe_values.cc:294
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:127
const unsigned int n_points()
Returns the number of quadrature points.
Definition: fe_values.cc:202
std::vector< arma::mat::fixed< spacedim, dim > > jacobians
Jacobians of the mapping at the quadrature points.
Definition: fe_values.hh:83
std::vector< arma::vec::fixed< spacedim > > points
Coordinates of quadrature points in the actual cell coordinate system.
Definition: fe_values.hh:98
Abstract class for the mapping between reference and actual cell.
Definition: fe_values.hh:44
virtual const arma::vec::fixed< spacedim > normal_vector(unsigned int point_no)=0
Returns the normal vector to a side at given quadrature point.
FiniteElement< dim, spacedim > * fe
The used finite element.
Definition: fe_values.hh:337
std::vector< arma::mat > shape_gradients
Gradients of shape functions evaluated at the quadrature points.
Definition: fe_values.hh:110
ElementFullIter * present_cell
Iterator to the last reinit-ed cell.
Definition: fe_values.hh:137
const Quadrature< dim-1 > * sub_quadrature
Quadrature for the integration on the element sides.
Definition: fe_values.hh:455
MappingInternalData * side_mapping_data[RefElement< dim >::n_sides][RefElement< dim >::n_side_permutations]
Definition: fe_values.hh:459
virtual const arma::vec::fixed< spacedim > shape_grad(const unsigned int function_no, const unsigned int point_no)=0
Return the gradient of the function_no-th shape function at the point_no-th quadrature point...
Quadrature< dim > side_quadrature[RefElement< dim >::n_sides][RefElement< dim >::n_side_permutations]
Definition: fe_values.hh:457
Quadrature< dim > * quadrature
The quadrature rule used to calculate integrals.
Definition: fe_values.hh:332
virtual const unsigned int n_dofs()=0
Returns the number of shape functions.
const vector< arma::vec::fixed< spacedim > > & point_list()
Return coordinates of all quadrature points in the actual cell system.
Definition: fe_values.cc:190
void allocate(unsigned int size, UpdateFlags flags, bool is_scalar=true)
Resize the data arrays.
Definition: fe_values.cc:53
std::vector< std::vector< arma::vec::fixed< spacedim > > > shape_vectors
Shape functions (for vectorial finite elements) evaluated at quadrature points.
Definition: fe_values.hh:116
const 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.cc:160
Class FEValuesData holds the arrays of data computed by Mapping and FiniteElement.
Definition: fe_values.hh:57
FEInternalData * side_fe_data[RefElement< dim >::n_sides][RefElement< dim >::n_side_permutations]
Definition: fe_values.hh:461
virtual const double shape_value(const unsigned int function_no, const unsigned int point_no)=0
Return the value of the function_no-th shape function at the point_no-th quadrature point...
FEValues(Mapping< dim, spacedim > &_mapping, Quadrature< dim > &_quadrature, FiniteElement< dim, spacedim > &_fe, UpdateFlags _flags)
Constructor.
Definition: fe_values.cc:228
void reinit(ElementFullIter &cell, unsigned int sid)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:312
const 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:148
Mapping< dim, spacedim > * mapping
The mapping from the reference cell to the actual cell.
Definition: fe_values.hh:327
void reinit(ElementFullIter &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:244
std::vector< double > determinants
Determinants of Jacobians at quadrature points.
Definition: fe_values.hh:88
Calculates finite element data on the actual cell.
Definition: fe_values.hh:370
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc...
Structure for storing the precomputed finite element data.
FESideValues(Mapping< dim, spacedim > &_mapping, Quadrature< dim-1 > &_sub_quadrature, FiniteElement< dim, spacedim > &_fe, UpdateFlags flags)
Constructor.
Definition: fe_values.cc:269
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Definition: dofhandler.hh:40
Base class for FEValues and FESideValues.
Definition: fe_values.hh:43
Abstract base class with certain methods independent of the template parameter dim.
Definition: fe_values.hh:146
void allocate(Mapping< dim, spacedim > &_mapping, Quadrature< dim > &_quadrature, FiniteElement< dim, spacedim > &_fe, UpdateFlags flags)
Allocates space for computed data.
Definition: fe_values.cc:123
Mapping data that can be precomputed on the actual cell.
Definition: mapping.hh:110
const unsigned int n_dofs()
Returns the number of shape functions.
Definition: fe_values.cc:208
const 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.cc:166
Calculates finite element data on a side.
Definition: fe_values.hh:415