Flow123d  release_3.0.0-1191-g7740876
finite_element.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 finite_element.hh
15  * @brief Abstract class for description of finite elements.
16  * @author Jan Stebel
17  */
18 
19 #ifndef FINITE_ELEMENT_HH_
20 #define FINITE_ELEMENT_HH_
21 
22 #include <armadillo>
23 #include <stdio.h> // for sprintf
24 #include <string.h> // for memcpy
25 #include <algorithm> // for max, min
26 #include <boost/assign/list_of.hpp> // for generic_list, list_of
27 #include <boost/exception/info.hpp> // for error_info::error_info...
28 #include <new> // for operator new[]
29 #include <ostream> // for operator<<
30 #include <string> // for operator<<
31 #include <vector> // for vector
32 #include "fem/update_flags.hh" // for operator&, operator|=
33 #include "system/exceptions.hh" // for ExcAssertMsg::~ExcAsse...
34 
35 template<unsigned int dim> class FESystem;
36 template<unsigned int dim, unsigned int spacedim> class FESideValues;
37 template<unsigned int dim, unsigned int spacedim> class FEValues;
38 template<unsigned int dim, unsigned int spacedim> class FEValuesBase;
39 template<unsigned int dim, unsigned int spacedim> class FEValuesData;
40 template<unsigned int dim> class FE_P_disc;
41 
42 
43 
44 
45 
46 // Possible types are: value, gradient, cell integral, ...
47 enum DofType { Value = 1 };
48 
49 /**
50  * Class Dof is a general description of functionals (degrees of freedom)
51  * determining the finite element. We assume that the Dof is defined as
52  * a linear combination of components of function value at a given point:
53  *
54  * Dof_value = a_1.f_1(x) + ... + a_n.f_n(x),
55  *
56  * where (a_1, ... , a_n) are given by @p coefs, x is the support point
57  * given by barycentric @p coords and (f_1, ..., f_n) is a generally
58  * vector-valued function. For the simplest Dof, i.e. value of a scalar
59  * function at point p, we set
60  *
61  * @p coords = p, @p coefs = { 1 }.
62  * The member @p dim denotes the affiliation of Dof to n-face:
63  * Nodal dofs have dim = 0,
64  * Dofs on lines: dim = 1,
65  * Dofs on triangles: dim = 2,
66  * Dofs in tetrahedron: dim = 3.
67  * It means that when a node, line or triangle is shared by adjacent cells,
68  * also the Dofs on this n-face are shared by the cells. Therefore
69  * for DG finite elements we set for all dofs the highest possible dimension.
70  *
71  * The class implements the method evaluate() which computes the Dof value
72  * for a basis function from given FunctionSpace.
73  */
74 class Dof {
75 public:
76 
77  Dof(unsigned int dim_,
78  unsigned int n_face_idx_,
79  arma::vec coords_,
80  arma::vec coefs_,
81  DofType type_)
82 
83  : dim(dim_),
84  n_face_idx(n_face_idx_),
85  coords(coords_),
86  coefs(coefs_),
87  type(type_)
88  {}
89 
90  /// Evaulate dof for basis function of given function space.
91  template<class FS>
92  const double evaluate(const FS &function_space,
93  unsigned int basis_idx) const;
94 
95  /// Association to n-face of given dimension (point, line, triangle, tetrahedron.
96  unsigned int dim;
97 
98  /// Index of n-face to which the dof is associated.
99  unsigned int n_face_idx;
100 
101  /// Barycentric coordinates.
103 
104  /// Coefficients of linear combination of function value components.
106 
107  /**
108  * Currently we consider only type=Value, possibly we could have Gradient,
109  * CellIntegral, FaceIntegral or other types.
110  */
112 };
113 
114 
115 /**
116  * FunctionSpace is an abstract class that is used to describe finite elements.
117  * It is determined by the dimension of the field of definition (@p space_dim_),
118  * by the dimension of the range (@p n_components_) and by the values and
119  * gradients of a basis functions.
120  *
121  * Combining FunctionSpace with Dof(s), the FiniteElement class constructs the
122  * shape functions, i.e. basis of FunctionSpace for which the Dof(s) attain
123  * the value 0 or 1.
124  */
126 public:
127 
128  FunctionSpace(unsigned int space_dim,
129  unsigned int n_components)
130 
131  : space_dim_(space_dim),
132  n_components_(n_components)
133  {}
134 
135 
136  /**
137  * @brief Value of the @p i th basis function at point @p point.
138  * @param basis_index Index of the basis function.
139  * @param point Point coordinates.
140  * @param comp_index Index of component (>0 for vector-valued functions).
141  */
142  virtual const double basis_value(unsigned int basis_index,
143  const arma::vec &point,
144  unsigned int comp_index = 0
145  ) const = 0;
146 
147  /**
148  * @brief Gradient of the @p i th basis function at point @p point.
149  * @param basis_index Index of the basis function.
150  * @param point Point coordinates.
151  * @param comp_index Index of component (>0 for vector-valued functions).
152  */
153  virtual const arma::vec basis_grad(unsigned int basis_index,
154  const arma::vec &point,
155  unsigned int comp_index = 0
156  ) const = 0;
157 
158  /// Dimension of function space (number of basis functions).
159  virtual const unsigned int dim() const = 0;
160 
161  /// Getter for space dimension.
162  const unsigned int space_dim() const { return space_dim_; }
163 
164  /// Getter for number of components.
165  const unsigned int n_components() const { return n_components_; }
166 
167  virtual ~FunctionSpace() {}
168 
169 protected:
170 
171  /// Space dimension of function arguments (i.e. 1, 2 or 3).
172  unsigned int space_dim_;
173 
174  /// Number of components of function values.
175  unsigned int n_components_;
176 };
177 
178 
179 /**
180  * Types of FiniteElement: scalar, vector-valued, tensor-valued or mixed system.
181  *
182  * The type also indicates how the shape functions and their gradients are transformed
183  * from reference element to arbitrary element. In particular:
184  *
185  * TYPE OBJECT EXPRESSION
186  * -----------------------------------------------------------
187  * FEScalar, FEVector, value ref_value
188  * FETensor grad J^{-T} * ref_grad
189  * -----------------------------------------------------------
190  * FEVectorContravariant value J * ref_value
191  * grad J^{-T} * ref_grad * J^T
192  * -----------------------------------------------------------
193  * FEVectorPiola value J * ref_value / |J|
194  * grad J^{-T} * ref_grad * J^T / |J|
195  * -----------------------------------------------------------
196  * FEMixedSystem value depends on sub-elements
197  * grad depends on sub-elements
198  *
199  * The transformation itself is done in FEValuesBase::fill_..._data() methods.
200  *
201  * Note that we use columnwise gradients, i.e. gradient of each component is a column vector.
202  */
203 enum FEType {
204  FEScalar = 0,
205  FEVector = 1,
208  FETensor = 4,
210 };
211 
212 
213 
214 /**
215  * @brief Abstract class for the description of a general finite element on
216  * a reference simplex in @p dim dimensions.
217  *
218  * The finite element is determined by a @p function_space_ and a set
219  * of @p dofs_. Further it must be set whether the finite element
220  * @p is_primitive_, which means that even if the functions in
221  * @p function_space_ are vector-valued, the basis functions have
222  * only one nonzero component (this is typical for tensor-product FE,
223  * e.g. vector-valued polynomial FE, but not for Raviart-Thomas FE).
224  *
225  * Description of dofs:
226  *
227  * The reference cell consists of lower dimensional entities (points,
228  * lines, triangles). Each dof is associated to one of these
229  * entities. If the entity is shared by 2 or more neighbouring cells
230  * in the mesh then this dof is shared by the finite elements on all
231  * of these cells. If a dof is associated to the cell itself then it
232  * is not shared with neighbouring cells.
233  *
234  *
235  * Shape functions:
236  *
237  * Sometimes it is convenient to describe the function space using
238  * a basis (called the raw basis) that is different from the set of
239  * shape functions for the finite element (the actual basis). For
240  * this reason we define the support points which play the role of
241  * nodal functionals associated to the particular dofs. To convert
242  * between the two bases one can use the @p node_matrix, which is
243  * constructed by the method compute_node_matrix(). In the case of
244  * non-Lagrangean finite elements the dofs are not associated to
245  * nodal functionals but e.g. to derivatives or averages. For that
246  * reason we distinguish between the unit support points which are
247  * uniquely associated to the dofs and the generalized support
248  * points that are auxiliary for the calculation of the dof
249  * functionals.
250  *
251  *
252  */
253 template<unsigned int dim>
254 class FiniteElement {
255 public:
256 
257  /**
258  * @brief Constructor.
259  */
260  FiniteElement();
261 
262  /**
263  * @brief Returns the number of degrees of freedom needed by the finite
264  * element.
265  */
266  inline const unsigned int n_dofs() const
267  { return dofs_.size(); }
268 
269  /**
270  * @brief Calculates the value of the @p comp-th component of
271  * the @p i-th shape function at the
272  * point @p p on the reference element.
273  *
274  * @param i Number of the shape function.
275  * @param p Point of evaluation.
276  * @param comp Number of vector component.
277  */
278  double shape_value(const unsigned int i,
279  const arma::vec::fixed<dim> &p,
280  const unsigned int comp = 0) const;
281 
282  /**
283  * @brief Calculates the @p comp-th component of the gradient
284  * of the @p i-th shape function at the point @p p on the
285  * reference element.
286  *
287  * @param i Number of the shape function.
288  * @param p Point of evaluation.
289  * @param comp Number of vector component.
290  */
291  arma::vec::fixed<dim> shape_grad(const unsigned int i,
292  const arma::vec::fixed<dim> &p,
293  const unsigned int comp = 0) const;
294 
295  /// Returns numer of components of the basis function.
296  inline unsigned int n_components() const
297  { return function_space_->n_components(); }
298 
299  /// Returns @p i -th degree of freedom.
300  inline const Dof &dof(unsigned int i) const
301  { return dofs_[i]; }
302 
303  /// Number of components of FE in a mapped space with dimension @p spacedim.
304  unsigned int n_space_components(unsigned int spacedim);
305 
306  /**
307  * @brief Destructor.
308  */
309  virtual ~FiniteElement() {};
310 
311 protected:
312 
313  /**
314  * @brief Clears all internal structures.
315  */
316  void init(bool primitive = true,
317  FEType type = FEScalar);
318 
319  /**
320  * @brief Initialize vectors with information about components of basis functions.
321  */
322  void setup_components();
323 
324  /**
325  * @brief Decides which additional quantities have to be computed
326  * for each cell.
327  *
328  * @param flags Computed update flags.
329  */
330  virtual UpdateFlags update_each(UpdateFlags flags);
331 
332  /**
333  * @brief Initializes the @p node_matrix for computing the coefficients
334  * of the shape functions in the raw basis of @p functions_space_.
335  * This is done by evaluating the @p dofs_ for the basis function
336  * and by inverting the resulting matrix.
337  */
338  virtual void compute_node_matrix();
339 
340  /**
341  * @brief Indicates whether the basis functions have one or more
342  * nonzero components (scalar FE spaces are always primitive).
343  */
344  inline const bool is_primitive() const
345  { return is_primitive_; }
346 
347  /**
348  * @brief Returns the component index for vector valued finite elements.
349  * @param sys_idx Index of shape function.
350  */
351  inline unsigned int system_to_component_index(unsigned sys_idx) const
352  { return component_indices_[sys_idx]; }
353 
354  /**
355  * @brief Returns the mask of nonzero components for given basis function.
356  * @param sys_idx Index of basis function.
357  */
358  inline const std::vector<bool> &get_nonzero_components(unsigned int sys_idx) const
359  { return nonzero_components_[sys_idx]; }
360 
361 
362 
363  /// Type of FiniteElement.
365 
366  /**
367  * @brief Primitive FE is using componentwise shape functions,
368  * i.e. only one component is nonzero for each shape function.
369  */
371 
372  /// Indices of nonzero components of shape functions (for primitive FE).
374 
375  /// Footprints of nonzero components of shape functions.
377 
378  /**
379  * @brief Matrix that determines the coefficients of the raw basis
380  * functions from the values at the support points.
381  */
383 
384  /// Function space defining the FE.
385  std::shared_ptr<FunctionSpace> function_space_;
386 
387  /// Set of degrees of freedom (functionals) defining the FE.
389 
390 
391  friend class FESystem<dim>;
392  friend class FEValuesBase<dim,3>;
393  friend class FEValues<dim,3>;
394  friend class FESideValues<dim,3>;
395  friend class FE_P_disc<dim>;
396  friend class SubDOFHandlerMultiDim;
397 };
398 
399 
400 
401 
402 #endif /* FINITE_ELEMENT_HH_ */
const unsigned int space_dim() const
Getter for space dimension.
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
Definition: update_flags.hh:67
arma::vec coords
Barycentric coordinates.
Mat< double, N, M > mat
Definition: armor.hh:214
DofType type
FunctionSpace(unsigned int space_dim, unsigned int n_components)
arma::mat node_matrix
Matrix that determines the coefficients of the raw basis functions from the values at the support poi...
unsigned int n_components_
Number of components of function values.
bool is_primitive_
Primitive FE is using componentwise shape functions, i.e. only one component is nonzero for each shap...
unsigned int dim
Association to n-face of given dimension (point, line, triangle, tetrahedron.
std::vector< std::vector< bool > > nonzero_components_
Footprints of nonzero components of shape functions.
Mat< double, N, 1 > vec
Definition: armor.hh:211
FEType
const bool is_primitive() const
Indicates whether the basis functions have one or more nonzero components (scalar FE spaces are alway...
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
virtual ~FiniteElement()
Destructor.
unsigned int n_face_idx
Index of n-face to which the dof is associated.
Discontinuous Lagrangean finite element on dim dimensional simplex.
Definition: fe_p.hh:110
Compound finite element on dim dimensional simplex.
Definition: fe_system.hh:101
std::vector< Dof > dofs_
Set of degrees of freedom (functionals) defining the FE.
unsigned int n_components() const
Returns numer of components of the basis function.
unsigned int system_to_component_index(unsigned sys_idx) const
Returns the component index for vector valued finite elements.
const unsigned int n_components() const
Getter for number of components.
Class FEValuesData holds the arrays of data computed by Mapping and FiniteElement.
Definition: fe_values.hh:86
const double evaluate(const FS &function_space, unsigned int basis_idx) const
Evaulate dof for basis function of given function space.
std::shared_ptr< FunctionSpace > function_space_
Function space defining the FE.
DofType
const Dof & dof(unsigned int i) const
Returns i -th degree of freedom.
const unsigned int n_dofs() const
Returns the number of degrees of freedom needed by the finite element.
Calculates finite element data on the actual cell.
Definition: fe_values.hh:525
Dof(unsigned int dim_, unsigned int n_face_idx_, arma::vec coords_, arma::vec coefs_, DofType type_)
unsigned int space_dim_
Space dimension of function arguments (i.e. 1, 2 or 3).
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Base class for FEValues and FESideValues.
Definition: fe_values.hh:37
arma::vec coefs
Coefficients of linear combination of function value components.
std::vector< unsigned int > component_indices_
Indices of nonzero components of shape functions (for primitive FE).
const std::vector< bool > & get_nonzero_components(unsigned int sys_idx) const
Returns the mask of nonzero components for given basis function.
FEType type_
Type of FiniteElement.
Calculates finite element data on a side.
Definition: fe_values.hh:582
virtual ~FunctionSpace()