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