Flow123d  release_3.0.0-973-g92f55e826
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 
27 #include <new> // for operator new[]
28 #include <ostream> // for operator<<
29 #include <string> // for operator<<
30 #include <vector> // for vector
31 #include "fem/update_flags.hh" // for operator&, operator|=
32 #include "system/exceptions.hh" // for ExcAssertMsg::~ExcAsse...
33 
34 template<unsigned int dim> class FESystem;
35 template<unsigned int dim, unsigned int spacedim> class FESideValues;
36 template<unsigned int dim, unsigned int spacedim> class FEValues;
37 template<unsigned int dim, unsigned int spacedim> class FEValuesBase;
38 template<unsigned int dim, unsigned int spacedim> class FEValuesData;
39 template<unsigned int dim> class FE_P_disc;
40 template<unsigned int dim> class Quadrature;
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.
102  arma::vec coords;
103 
104  /// Coefficients of linear combination of function value components.
105  arma::vec coefs;
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 
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  /**
304  * @brief Destructor.
305  */
306  virtual ~FiniteElement() {};
307 
308 protected:
309 
310  /**
311  * @brief Clears all internal structures.
312  */
313  void init(bool primitive = true,
314  FEType type = FEScalar);
315 
316  /**
317  * @brief Initialize vectors with information about components of basis functions.
318  */
319  void setup_components();
320 
321  /**
322  * @brief Decides which additional quantities have to be computed
323  * for each cell.
324  *
325  * @param flags Computed update flags.
326  */
327  virtual UpdateFlags update_each(UpdateFlags flags);
328 
329  /**
330  * @brief Initializes the @p node_matrix for computing the coefficients
331  * of the shape functions in the raw basis of @p functions_space_.
332  * This is done by evaluating the @p dofs_ for the basis function
333  * and by inverting the resulting matrix.
334  */
335  virtual void compute_node_matrix();
336 
337  /**
338  * @brief Indicates whether the basis functions have one or more
339  * nonzero components (scalar FE spaces are always primitive).
340  */
341  inline const bool is_primitive() const
342  { return is_primitive_; }
343 
344  /**
345  * @brief Returns the component index for vector valued finite elements.
346  * @param sys_idx Index of shape function.
347  */
348  inline unsigned int system_to_component_index(unsigned sys_idx) const
349  { return component_indices_[sys_idx]; }
350 
351  /**
352  * @brief Returns the mask of nonzero components for given basis function.
353  * @param sys_idx Index of basis function.
354  */
355  inline const std::vector<bool> &get_nonzero_components(unsigned int sys_idx) const
356  { return nonzero_components_[sys_idx]; }
357 
358 
359 
360  /// Type of FiniteElement.
362 
363  /**
364  * @brief Primitive FE is using componentwise shape functions,
365  * i.e. only one component is nonzero for each shape function.
366  */
368 
369  /// Indices of nonzero components of shape functions (for primitive FE).
371 
372  /// Footprints of nonzero components of shape functions.
374 
375  /**
376  * @brief Matrix that determines the coefficients of the raw basis
377  * functions from the values at the support points.
378  */
379  arma::mat node_matrix;
380 
381  /// Function space defining the FE.
382  std::shared_ptr<FunctionSpace> function_space_;
383 
384  /// Set of degrees of freedom (functionals) defining the FE.
386 
387 
388  friend class FESystem<dim>;
389  friend class FEValuesBase<dim,3>;
390  friend class FEValues<dim,3>;
391  friend class FESideValues<dim,3>;
392  friend class FE_P_disc<dim>;
393 };
394 
395 
396 
397 
398 #endif /* FINITE_ELEMENT_HH_ */
FunctionSpace::~FunctionSpace
virtual ~FunctionSpace()
Definition: finite_element.hh:167
Dof::n_face_idx
unsigned int n_face_idx
Index of n-face to which the dof is associated.
Definition: finite_element.hh:99
FiniteElement::FiniteElement
FiniteElement()
Constructor.
Definition: finite_element.cc:73
FEScalar
@ FEScalar
Definition: finite_element.hh:204
FiniteElement::init
void init(bool primitive=true, FEType type=FEScalar)
Clears all internal structures.
Definition: finite_element.cc:80
FiniteElement::shape_grad
arma::vec::fixed< dim > shape_grad(const unsigned int i, const arma::vec::fixed< dim > &p, const unsigned int comp=0) const
Calculates the comp-th component of the gradient of the i-th shape function at the point p on the ref...
Definition: finite_element.cc:126
string.h
FiniteElement::nonzero_components_
std::vector< std::vector< bool > > nonzero_components_
Footprints of nonzero components of shape functions.
Definition: finite_element.hh:373
Dof::coefs
arma::vec coefs
Coefficients of linear combination of function value components.
Definition: finite_element.hh:105
DofType
DofType
Definition: finite_element.hh:47
std::vector< bool >
FiniteElement::node_matrix
arma::mat node_matrix
Matrix that determines the coefficients of the raw basis functions from the values at the support poi...
Definition: finite_element.hh:379
FiniteElement::dofs_
std::vector< Dof > dofs_
Set of degrees of freedom (functionals) defining the FE.
Definition: finite_element.hh:385
FiniteElement::dof
const Dof & dof(unsigned int i) const
Returns i -th degree of freedom.
Definition: finite_element.hh:300
FEMixedSystem
@ FEMixedSystem
Definition: finite_element.hh:209
FEType
FEType
Definition: finite_element.hh:203
FiniteElement::component_indices_
std::vector< unsigned int > component_indices_
Indices of nonzero components of shape functions (for primitive FE).
Definition: finite_element.hh:370
FiniteElement::~FiniteElement
virtual ~FiniteElement()
Destructor.
Definition: finite_element.hh:306
FunctionSpace::basis_grad
virtual const arma::vec basis_grad(unsigned int basis_index, const arma::vec &point, unsigned int comp_index=0) const =0
Gradient of the i th basis function at point point.
exceptions.hh
FiniteElement::system_to_component_index
unsigned int system_to_component_index(unsigned sys_idx) const
Returns the component index for vector valued finite elements.
Definition: finite_element.hh:348
Dof::dim
unsigned int dim
Association to n-face of given dimension (point, line, triangle, tetrahedron.
Definition: finite_element.hh:96
FEValues
Calculates finite element data on the actual cell.
Definition: fe_values.hh:532
FunctionSpace
Definition: finite_element.hh:125
FiniteElement::n_components
unsigned int n_components() const
Returns numer of components of the basis function.
Definition: finite_element.hh:296
FiniteElement::update_each
virtual UpdateFlags update_each(UpdateFlags flags)
Decides which additional quantities have to be computed for each cell.
Definition: finite_element.cc:145
FunctionSpace::dim
virtual const unsigned int dim() const =0
Dimension of function space (number of basis functions).
FiniteElement::is_primitive_
bool is_primitive_
Primitive FE is using componentwise shape functions, i.e. only one component is nonzero for each shap...
Definition: finite_element.hh:367
FiniteElement::compute_node_matrix
virtual void compute_node_matrix()
Initializes the node_matrix for computing the coefficients of the shape functions in the raw basis of...
Definition: finite_element.cc:97
FiniteElement::type_
FEType type_
Type of FiniteElement.
Definition: finite_element.hh:361
FiniteElement::n_dofs
const unsigned int n_dofs() const
Returns the number of degrees of freedom needed by the finite element.
Definition: finite_element.hh:266
FunctionSpace::space_dim
const unsigned int space_dim() const
Getter for space dimension.
Definition: finite_element.hh:162
FEVectorContravariant
@ FEVectorContravariant
Definition: finite_element.hh:206
FunctionSpace::n_components
const unsigned int n_components() const
Getter for number of components.
Definition: finite_element.hh:165
FETensor
@ FETensor
Definition: finite_element.hh:208
Dof
Definition: finite_element.hh:74
FEVector
@ FEVector
Definition: finite_element.hh:205
FiniteElement
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Definition: discrete_space.hh:25
FEValuesBase
Base class for FEValues and FESideValues.
Definition: fe_values.hh:37
Value
@ Value
Definition: finite_element.hh:47
FunctionSpace::basis_value
virtual const double basis_value(unsigned int basis_index, const arma::vec &point, unsigned int comp_index=0) const =0
Value of the i th basis function at point point.
FunctionSpace::FunctionSpace
FunctionSpace(unsigned int space_dim, unsigned int n_components)
Definition: finite_element.hh:128
update_flags.hh
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
FiniteElement::is_primitive
const bool is_primitive() const
Indicates whether the basis functions have one or more nonzero components (scalar FE spaces are alway...
Definition: finite_element.hh:341
Dof::Dof
Dof(unsigned int dim_, unsigned int n_face_idx_, arma::vec coords_, arma::vec coefs_, DofType type_)
Definition: finite_element.hh:77
FEValuesData
Class FEValuesData holds the arrays of data computed by Mapping and FiniteElement.
Definition: fe_values.hh:86
Dof::coords
arma::vec coords
Barycentric coordinates.
Definition: finite_element.hh:102
Dof::type
DofType type
Definition: finite_element.hh:111
FiniteElement::function_space_
std::shared_ptr< FunctionSpace > function_space_
Function space defining the FE.
Definition: finite_element.hh:382
UpdateFlags
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
Definition: update_flags.hh:67
FunctionSpace::space_dim_
unsigned int space_dim_
Space dimension of function arguments (i.e. 1, 2 or 3).
Definition: finite_element.hh:172
FESideValues
Calculates finite element data on a side.
Definition: fe_values.hh:577
FiniteElement::setup_components
void setup_components()
Initialize vectors with information about components of basis functions.
Definition: finite_element.cc:89
Quadrature
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: fe_values.hh:35
FEVectorPiola
@ FEVectorPiola
Definition: finite_element.hh:207
FunctionSpace::n_components_
unsigned int n_components_
Number of components of function values.
Definition: finite_element.hh:175
FESystem
Compound finite element on dim dimensional simplex.
Definition: fe_system.hh:99
Dof::evaluate
const double evaluate(const FS &function_space, unsigned int basis_idx) const
Evaulate dof for basis function of given function space.
Definition: finite_element.cc:33
FE_P_disc
Discontinuous Lagrangean finite element on dim dimensional simplex.
Definition: fe_p.hh:108
FiniteElement::get_nonzero_components
const std::vector< bool > & get_nonzero_components(unsigned int sys_idx) const
Returns the mask of nonzero components for given basis function.
Definition: finite_element.hh:355
FiniteElement::shape_value
double shape_value(const unsigned int i, const arma::vec::fixed< dim > &p, const unsigned int comp=0) const
Calculates the value of the comp-th component of the i-th shape function at the point p on the refere...
Definition: finite_element.cc:111