Flow123d  release_2.2.0-914-gf1a3a4f
fe_p.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 fe_p.cc
15  * @brief
16  */
17 
18 // !! implementation of specializations has to be i *.cc file to avoid multiple definition error during linking
19 #include "fe_p.hh"
20 #include "mesh/ref_element.hh"
21 
22 
23 
24 
25 PolynomialSpace::PolynomialSpace(unsigned int degree, unsigned int dim)
26  : FunctionSpace(dim, 1),
27  degree_(degree)
28 {
29 // computes powers of all monomials up to given @p degree
30 // the order is: 1,x,x^2, y, yx,y^2
31 //
32 // TODO: - check and possibly rewrite to be more clear (use sum_degree temporary
33 // - change order of monomials: 1, x, y, xy, x^2 , y^2 (increasing order)
34 // - allow Q polynomials: 1,x, y, xy
35 // - can use tensor products
36 
37  arma::uvec pows(dim);
38  pows.zeros();
39 
40  unsigned int degree_sum=0;
41  unsigned int i_dim;
42 
43 
44  while (true) {
45  powers.push_back(pows);
46 
47  // increment pows
48  for(i_dim=0; i_dim < dim; i_dim++) {
49  if (degree_sum < degree) {
50  pows[i_dim]++;
51  degree_sum++;
52  break;
53  } else { // if degree_sum == degree, we find first non empty power, free it, and raise the next one
54  degree_sum-=pows[i_dim];
55  pows[i_dim]=0;
56  }
57  }
58  if (i_dim == dim) break; // just after pow == (0, 0, .., degree)
59  }
60 }
61 
62 
63 const double PolynomialSpace::basis_value(unsigned int i,
64  const arma::vec &point,
65  unsigned int comp_index
66  ) const
67 {
68  ASSERT(comp_index == 0);
69  OLD_ASSERT(i<=powers.size(), "Index of basis function is out of range.");
70  ASSERT(point.size()==space_dim_);
71 
72  double v = 1;
73  for (unsigned int j=0; j<this->space_dim_; j++)
74  v *= pow(point[j], (int) powers[i][j]);
75 
76  return v;
77 }
78 
79 
80 const arma::vec PolynomialSpace::basis_grad(unsigned int i,
81  const arma::vec &p,
82  unsigned int comp_index
83  ) const
84 {
85  ASSERT(comp_index == 0);
86  OLD_ASSERT(i<=powers.size(), "Index of basis function is out of range.");
87 
88  arma::vec grad(this->space_dim_);
89 
90  for (unsigned int j=0; j<this->space_dim_; j++)
91  {
92  grad[j] = powers[i][j];
93  if (powers[i][j] == 0) continue;
94 
95  for (unsigned int k=0; k<this->space_dim_; k++)
96  {
97  grad[j] *= pow(p[k], (int) (k==j?powers[i][k]-1:powers[i][k]));
98  }
99  }
100  return grad;
101 }
102 
103 
104 
105 
106 
107 
108 
109 
110 
111 
112 
113 template<unsigned int dim, unsigned int spacedim>
115 {
116  if (degree_ == 0)
117  {
118  // we define nodal dof:
119  // coords = barycentric coordinates of the support point,
120  // coefs = 1 (dof value = function value at the point)
121  arma::vec coords = arma::ones<arma::vec>(dim+1)/(dim+1);
122  this->dofs_.push_back(Dof(dim, 0, coords, { 1 }, Value));
123  }
124  else
125  {
126  // Create vector of barycentric coordinates.
127  // First we make vector uvbc which contains barycentric coordinates
128  // multiplied by degree_, so that its values are unsigned ints.
129  // Then by counting the nonzero barycentric coordinates we can decide
130  // whether the dof lies on node, line, triangle or tetrahedron.
132  arma::uvec ubc = arma::zeros<arma::uvec>(dim+1);
133  ubc[dim] = degree_;
134  bool finish = false;
135  do {
136  uvbc.push_back(ubc);
137  if (ubc[dim] > 0)
138  {
139  // by default increment the first coordinate
140  ubc[0] += 1;
141  ubc[dim] -= 1;
142  }
143  else
144  {
145  // if sum of coordinates is maximal (last coordinate is zero)
146  // then find first nonzero coordinate,
147  // set it to zero, and increment the following coordinate.
148  unsigned int c = 0;
149  while (ubc[c] == 0) c++;
150  // if the first nonzero coordinate is the last but one, we reach the end
151  if (c == dim-1) finish = true;
152  else {
153  ubc[dim] = ubc[c]-1;
154  ubc[c] = 0;
155  ubc[c+1] += 1;
156  }
157  }
158  } while (!finish);
159 
160  // define dofs
161  for (auto ubc : uvbc)
162  {
163  // we define nodal dof:
164  // coords = barycentric coordinates of the support point,
165  // coefs = 1 (dof value = function value at the point)
166 
167  // count nonzero coordinates in ubc: 1=>nodal dof, 2=>dof on line etc.
168  arma::uvec nonzeros = find(ubc);
169  // convert "unsigned barycentric coordinates" to real ones
170  arma::vec coords = arma::conv_to<arma::vec>::from(ubc);
171  coords /= degree_;
172 
173  // find index of n-face
174  std::pair<unsigned int, unsigned int> zeros = RefElement<dim>::zeros_positions(coords);
175  unsigned int n_face_idx;
176  switch (dim-zeros.first) {
177  case 0:
178  n_face_idx = RefElement<dim>::template topology_idx<0>(zeros.second);
179  break;
180  case 1:
181  n_face_idx = RefElement<dim>::template topology_idx<1>(zeros.second);
182  break;
183  case 2:
184  n_face_idx = RefElement<dim>::template topology_idx<2>(zeros.second);
185  break;
186  case 3:
187  n_face_idx = RefElement<dim>::template topology_idx<3>(zeros.second);
188  break;
189  }
190  this->dofs_.push_back(Dof(nonzeros.size()-1, n_face_idx, coords, { 1 }, Value));
191  }
192  }
193 }
194 
195 
196 
197 
198 template<unsigned int dim, unsigned int spacedim>
199 FE_P<dim,spacedim>::FE_P(unsigned int degree)
200  : FiniteElement<dim,spacedim>(),
201  degree_(degree)
202 {
203  this->function_space_ = new PolynomialSpace(degree,dim);
204 
205  init_dofs();
206 
207  this->setup_components();
208 
209  this->compute_node_matrix();
210 }
211 
212 
213 
214 
215 
216 
217 
218 
219 
220 
221 
222 
223 
224 
225 
226 template<unsigned int dim, unsigned int spacedim>
228  : FE_P<dim,spacedim>(degree)
229 {
230  // all dofs are "inside" the cell (not shared with neighbours)
231  for (unsigned int i=0; i<this->dofs_.size(); i++)
232  this->dofs_[i].dim = dim;
233 
234  this->setup_components();
235 
236  this->compute_node_matrix();
237 }
238 
239 
240 
241 
242 
243 template class FE_P_disc<1, 3>;
244 template class FE_P_disc<2, 3>;
245 template class FE_P_disc<3, 3>;
const arma::vec basis_grad(unsigned int basis_index, const arma::vec &point, unsigned int comp_index=0) const override
Gradient of the i th basis function at point point.
Definition: fe_p.cc:80
PolynomialSpace(unsigned int degree, unsigned int dim)
Constructor.
Definition: fe_p.cc:25
Space of polynomial functions.
Definition: fe_p.hh:35
std::vector< Dof > dofs_
Set of degrees of freedom (functionals) defining the FE.
const double basis_value(unsigned int basis_index, const arma::vec &point, unsigned int comp_index=0) const override
Value of the i th basis function at point point.
Definition: fe_p.cc:63
static std::pair< unsigned int, unsigned int > zeros_positions(const BaryPoint &barycentric, double tolerance=std::numeric_limits< double >::epsilon()*2)
Definition: ref_element.cc:393
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:346
std::vector< arma::uvec > powers
Coefficients of basis functions.
Definition: fe_p.hh:69
#define OLD_ASSERT(...)
Definition: global_defs.h:131
FE_P_disc(unsigned int degree)
Constructor.
Definition: fe_p.cc:227
virtual void compute_node_matrix()
Initializes the node_matrix for computing the coefficients of the shape functions in the raw basis of...
Definitions of basic Lagrangean finite elements with polynomial shape functions.
void setup_components()
Initialize vectors with information about components of basis functions.
void init_dofs()
Definition: fe_p.cc:114
FE_P(unsigned int degree)
Constructor.
Definition: fe_p.cc:199
Conforming Lagrangean finite element on dim dimensional simplex.
Definition: fe_p.hh:84
const unsigned int dim() const override
Dimension of function space (number of basis functions).
Definition: fe_p.hh:56
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc...
unsigned int space_dim_
Space dimension of function arguments (i.e. 1, 2 or 3).
FunctionSpace * function_space_
Function space defining the FE.
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Definition: dofhandler.hh:30
const unsigned int degree_
Max. degree of polynomials.
Definition: fe_p.hh:61