Flow123d  release_3.0.0-893-gf7bf019
finite_element.cc
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.cc
15  * @brief Abstract class for description of finite elements.
16  * @author Jan Stebel
17  */
18 
19 #include "system/system.hh"
20 #include "quadrature/quadrature.hh"
21 #include "fem/dofhandler.hh"
22 #include "fem/finite_element.hh"
23 #include "fem/fe_values.hh"
24 
25 
26 
27 using namespace std;
28 
29 
30 
31 
32 
33 template<class FS> const double Dof::evaluate(const FS &function_space,
34  unsigned int basis_idx) const
35 {
36  // Check that FS is derived from FunctionSpace.
37  static_assert(std::is_base_of<FunctionSpace, FS>::value, "FS must be derived from FunctionSpace.");
38 
39  // We cannot evaluate dof on dim-dimensional n-face if the function space lies on lower-dimensional n-face.
40  ASSERT(function_space.space_dim()+1 == coords.size());
41 
42  switch (type)
43  {
44  case Value:
45  {
46  // evaluate basis function and return the linear combination of components
47  arma::vec vec_value(function_space.n_components());
48 
49  // drop off the 0-th barycentric coordinate
50  // in case of space_dim=0 subvec does not work
51  arma::vec f_coords(function_space.space_dim());
52  if (function_space.space_dim() > 0)
53  f_coords = coords.subvec(1,coords.size()-1);
54  for (unsigned int c=0; c<function_space.n_components(); c++)
55  vec_value[c] = function_space.basis_value(basis_idx, f_coords, c);
56  return dot(coefs, vec_value);
57  break;
58  }
59 
60  default:
61  OLD_ASSERT(false, "Dof evaluation not implemented for this type.");
62  }
63  return 0;
64 }
65 
66 
67 
68 
69 
70 
71 
72 template<unsigned int dim>
74  : function_space_(nullptr)
75 {
76  init();
77 }
78 
79 template<unsigned int dim>
80 void FiniteElement<dim>::init(bool primitive, FEType type)
81 {
82  dofs_.clear();
83  is_primitive_ = primitive;
84  type_ = type;
85 }
86 
87 
88 template<unsigned int dim>
90 {
91  component_indices_.resize(dofs_.size(), 0);
92  nonzero_components_.resize(dofs_.size(), { true });
93 }
94 
95 
96 template<unsigned int dim> inline
98 {
99  arma::mat M(dofs_.size(), dofs_.size());
100 
101  for (unsigned int i = 0; i < dofs_.size(); i++)
102  for (unsigned int j = 0; j < dofs_.size(); j++) {
103  M(j, i) = dofs_[i].evaluate(*function_space_, j);
104 
105  }
106  node_matrix = arma::inv(M);
107 }
108 
109 
110 template<unsigned int dim>
111 double FiniteElement<dim>::shape_value(const unsigned int i,
112  const arma::vec::fixed<dim> &p,
113  const unsigned int comp) const
114 {
115  ASSERT_DBG( comp < n_components() );
116  ASSERT_DBG( i < dofs_.size()).error("Index of basis function is out of range.");
117 
118  double value = 0;
119  for (unsigned int j=0; j<function_space_->dim(); j++)
120  value += function_space_->basis_value(j, p, comp) * node_matrix(i,j);
121 
122  return value;
123 }
124 
125 template<unsigned int dim>
126 arma::vec::fixed<dim> FiniteElement<dim>::shape_grad(const unsigned int i,
127  const arma::vec::fixed<dim> &p,
128  const unsigned int comp) const
129 {
130  ASSERT_DBG( comp < n_components() );
131  ASSERT_DBG( i < dofs_.size()).error("Index of basis function is out of range.");
132 
133  // uword is a typedef for an unsigned integer type; it is used for matrix indices as well as all internal counters and loops
134  // sword is a typedef for a signed integer type
135  arma::vec grad((arma::uword)dim);
136  grad.zeros();
137  for (unsigned int j=0; j<function_space_->dim(); j++)
138  grad += function_space_->basis_grad(j, p, comp) * node_matrix(i,j);
139 
140  return grad;
141 }
142 
143 
144 template<unsigned int dim> inline
146 {
147  UpdateFlags f = flags;
148 
149  switch (type_)
150  {
151  case FEScalar:
152  case FEVector:
153  case FETensor:
154  if (flags & update_gradients)
156  break;
158  if (flags & update_values)
159  f |= update_jacobians;
160  if (flags & update_gradients)
162  break;
163  case FEVectorPiola:
164  if (flags & update_values)
166  if (flags & update_gradients)
168  break;
169  default:;
170  }
171 
172  return f;
173 }
174 
175 
176 
177 
178 
179 
180 
181 template class FiniteElement<0>;
182 template class FiniteElement<1>;
183 template class FiniteElement<2>;
184 template class FiniteElement<3>;
185 
186 
FiniteElement()
Constructor.
Shape function values.
Definition: update_flags.hh:87
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
Definition: update_flags.hh:67
Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between loca...
arma::mat node_matrix
Matrix that determines the coefficients of the raw basis functions from the values at the support poi...
Determinant of the Jacobian.
bool is_primitive_
Primitive FE is using componentwise shape functions, i.e. only one component is nonzero for each shap...
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...
std::vector< std::vector< bool > > nonzero_components_
Footprints of nonzero components of shape functions.
Class FEValues calculates finite element data on the actual cells such as shape function values...
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:346
FEType
static constexpr bool value
Definition: json.hpp:87
Volume element.
#define OLD_ASSERT(...)
Definition: global_defs.h:131
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< Dof > dofs_
Set of degrees of freedom (functionals) defining the FE.
virtual UpdateFlags update_each(UpdateFlags flags)
Decides which additional quantities have to be computed for each cell.
unsigned int n_components() const
Returns numer of components of the basis function.
virtual void compute_node_matrix()
Initializes the node_matrix for computing the coefficients of the shape functions in the raw basis of...
Basic definitions of numerical quadrature rules.
Shape function gradients.
Definition: update_flags.hh:95
void init(bool primitive=true, FEType type=FEScalar)
Clears all internal structures.
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.
void setup_components()
Initialize vectors with information about components of basis functions.
#define ASSERT_DBG(expr)
Definition: asserts.hh:349
Abstract class for description of finite elements.
std::vector< unsigned int > component_indices_
Indices of nonzero components of shape functions (for primitive FE).
FEType type_
Type of FiniteElement.