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