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