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