Flow123d
fe_p.hh
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2007 Technical University of Liberec. All rights reserved.
4  *
5  * Please make a following refer to Flow123d on your project site if you use the program for any purpose,
6  * especially for academic research:
7  * Flow123d, Research Centre: Advanced Remedial Technologies, Technical University of Liberec, Czech Republic
8  *
9  * This program is free software; you can redistribute it and/or modify it under the terms
10  * of the GNU General Public License version 3 as published by the Free Software Foundation.
11  *
12  * This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
13  * without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
14  * See the GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License along with this program; if not,
17  * write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 021110-1307, USA.
18  *
19  *
20  * $Id$
21  * $Revision$
22  * $LastChangedBy$
23  * $LastChangedDate$
24  *
25  * @file
26  * @brief Definitions of basic Lagrangean finite elements with polynomial shape functions.
27  * @author Jan Stebel
28  */
29 
30 #ifndef FE_P_HH_
31 #define FE_P_HH_
32 
33 #include <vector>
34 
35 #include "system/global_defs.h"
36 #include "system/system.hh"
37 #include "fem/finite_element.hh"
38 
39 
40 /**
41  * @brief Space of polynomial functions.
42  *
43  * This class serves for evaluation of the value and gradient
44  * of a polynomial of order @p degree in @p dim variables.
45  */
46 template<unsigned int degree, unsigned int dim>
48 {
49 public:
50 
51  /**
52  * @brief Constructor.
53  *
54  * Creates the coefficients of the basis.
55  */
57 
58  /**
59  * @brief Value of the @p i th basis function at point @p p.
60  * @param i Number of the basis function.
61  * @param p Point at which the function is evaluated.
62  */
63  const double basis_value(unsigned int i, const arma::vec::fixed<dim> &p) const;
64 
65  /**
66  * @brief Gradient of the @p i th basis function at point @p p.
67  * @param i Number of the basis function.
68  * @param p Point at which the function is evaluated.
69  */
70  const arma::vec::fixed<dim> basis_grad(unsigned int i, const arma::vec::fixed<dim> &p) const;
71 
72 private:
73 
74  /**
75  * @brief Coefficients of basis functions.
76  *
77  * Powers of x, y, z, ... in the i-th basis function are stored
78  * in powers[i].
79  */
81 
82 };
83 
84 
85 
86 /**
87  * @brief Distribution of dofs for polynomial finite elements.
88  *
89  * The class holds the information on the total number of dofs
90  * as well as the number of dofs associated to geometrical entities
91  * such as points, lines, triangles and tetrahedra.
92  * Moreover, some dofs are grouped to pairs, triples or sextuples
93  * which are invariant to rotation/reflection of the element.
94  *
95  * The coordinates of unit support points are provided.
96  * The values at support points uniquely determine the finite
97  * element function.
98  *
99  */
100 template<unsigned int degree, unsigned int dim>
102 {
103 public:
104 
105  /**
106  * @brief Constructor.
107  *
108  * Initializes all variables.
109  */
110  DofDistribution();
111 
112  /// Total number of degrees of freedom at one finite element.
113  unsigned int number_of_dofs;
114 
115  /**
116  * @brief Number of single dofs at one geometrical entity of the given
117  * dimension (point, line, triangle, tetrahedron).
118  */
119  unsigned int number_of_single_dofs[dim + 1];
120 
121  /**
122  * @brief Number of pairs of dofs at one geometrical entity of the given
123  * dimension (applicable to lines and triangles).
124  */
125  unsigned int number_of_pairs[dim + 1];
126 
127  /**
128  * @brief Number of triples of dofs associated to one triangle.
129  */
130  unsigned int number_of_triples[dim + 1];
131 
132  /**
133  * @brief Number of sextuples of dofs associated to one triangle.
134  */
135  unsigned int number_of_sextuples[dim + 1];
136 
137  /**
138  * @brief Support points.
139  *
140  * Support points are points in the reference element where
141  * function values determine the dofs. In case of Lagrangean
142  * finite elements the dof values are precisely the function
143  * values at @p unit_support_points.
144  */
146 
147 };
148 
149 
150 /**
151  * @brief Conforming Lagrangean finite element on @p dim dimensional simplex.
152  *
153  * The finite element functions are continuous across the interfaces.
154  */
155 template <unsigned int degree, unsigned int dim, unsigned int spacedim>
156 class FE_P : public FiniteElement<dim,spacedim>
157 {
165 
166 public:
167  /// Constructor.
168  FE_P();
169 
170  /**
171  * @brief Returns the @p ith basis function evaluated at the point @p p.
172  * @param i Number of the basis function.
173  * @param p Point of evaluation.
174  */
175  double basis_value(const unsigned int i, const arma::vec::fixed<dim> &p) const;
176 
177  /**
178  * @brief Returns the gradient of the @p ith basis function at the point @p p.
179  * @param i Number of the basis function.
180  * @param p Point of evaluation.
181  */
182  arma::vec::fixed<dim> basis_grad(const unsigned int i, const arma::vec::fixed<dim> &p) const;
183 
184  /**
185  * @brief The vector variant of basis_value must be implemented but may not be used.
186  */
187  arma::vec::fixed<dim> basis_vector(const unsigned int i, const arma::vec::fixed<dim> &p) const;
188 
189  /**
190  * @brief The vector variant of basis_grad must be implemented but may not be used.
191  */
192  arma::mat::fixed<dim,dim> basis_grad_vector(const unsigned int i, const arma::vec::fixed<dim> &p) const;
193 
194 private:
195 
196  /// The auxiliary polynomial space.
198 
199  /// The auxiliary dof distribution.
201 };
202 
203 
204 /**
205  * @brief Discontinuous Lagrangean finite element on @p dim dimensional simplex.
206  *
207  * No continuity of the finite element functions across the interfaces is
208  * imposed.
209  */
210 template <unsigned int degree, unsigned int dim, unsigned int spacedim>
211 class FE_P_disc : public FiniteElement<dim,spacedim>
212 {
220 
221 public:
222 
223  /// Constructor.
224  FE_P_disc();
225 
226  /**
227  * @brief Returns the @p ith basis function evaluated at the point @p p.
228  * @param i Number of the basis function.
229  * @param p Point of evaluation.
230  */
231  double basis_value(const unsigned int i, const arma::vec::fixed<dim> &p) const;
232 
233  /**
234  * @brief Returns the gradient of the @p ith basis function at the point @p p.
235  * @param i Number of the basis function.
236  * @param p Point of evaluation.
237  */
238  arma::vec::fixed<dim> basis_grad(const unsigned int i, const arma::vec::fixed<dim> &p) const;
239 
240  /**
241  * @brief The vector variant of basis_value must be implemented but may not be used.
242  */
243  arma::vec::fixed<dim> basis_vector(const unsigned int i, const arma::vec::fixed<dim> &p) const;
244 
245  /**
246  * @brief The vector variant of basis_grad must be implemented but may not be used.
247  */
248  arma::mat::fixed<dim,dim> basis_grad_vector(const unsigned int i, const arma::vec::fixed<dim> &p) const;
249 
250 private:
251 
252  /// The auxiliary polynomial space.
254 
255  /// The auxiliary dof distribution.
257 };
258 
259 
260 
261 template<unsigned int degree, unsigned int dim>
263 {
264 // computes powers of all monomials up to given @p degree
265 // the order is: 1,x,x^2, y, yx,y^2
266 //
267 // TODO: - check and possibly rewrite to be more clear (use sum_degree temporary
268 // - change order of monomials: 1, x, y, xy, x^2 , y^2 (increasing order)
269 // - allow Q polynomials: 1,x, y, xy
270 // - can use tensor products
271 
272  arma::uvec::fixed<dim> pows;
273  pows.zeros();
274 
275  unsigned int degree_sum=0;
276  unsigned int i_dim;
277 
278 
279  while (true) {
280  powers.push_back(pows);
281 
282  // increment pows
283  for(i_dim=0; i_dim < dim; i_dim++) {
284  if (degree_sum < degree) {
285  pows[i_dim]++;
286  degree_sum++;
287  break;
288  } else { // if degree_sum == degree, we find first non empty power, free it, and raise the next one
289  degree_sum-=pows[i_dim];
290  pows[i_dim]=0;
291  }
292  }
293  if (i_dim == dim) break; // just after pow == (0, 0, .., degree)
294  }
295 }
296 
297 template<unsigned int degree, unsigned int dim>
298 const double PolynomialSpace<degree,dim>::basis_value(unsigned int i, const arma::vec::fixed<dim> &p) const
299 {
300  ASSERT(i<=powers.size(), "Index of basis function is out of range.");
301 
302  double v = 1;
303 
304  for (unsigned int j=0; j<dim; j++)
305  v *= pow(p[j], (int) powers[i][j]);
306 
307  return v;
308 }
309 
310 
311 template<unsigned int degree, unsigned int dim>
312 const arma::vec::fixed<dim> PolynomialSpace<degree,dim>::basis_grad(unsigned int i, const arma::vec::fixed<dim> &p) const
313 {
314  ASSERT(i<=powers.size(), "Index of basis function is out of range.");
315 
316  arma::vec::fixed<dim> grad;
317 
318  for (unsigned int j=0; j<dim; j++)
319  {
320  grad[j] = powers[i][j];
321  if (powers[i][j] == 0) continue;
322 
323  for (unsigned int k=0; k<dim; k++)
324  {
325  grad[j] *= pow(p[k], (int) (k==j?powers[i][k]-1:powers[i][k]));
326  }
327  }
328  return grad;
329 }
330 
331 
332 
333 
334 
335 
336 
337 
338 
339 
340 
341 template<unsigned int degree, unsigned int dim, unsigned int spacedim>
343 {
344  this->init();
345 
346  for (int i=0; i<=dim; i++)
347  {
348  number_of_dofs += dof_distribution.number_of_single_dofs[i]
349  +2*dof_distribution.number_of_pairs[i]
350  +3*dof_distribution.number_of_triples[i]
351  +6*dof_distribution.number_of_sextuples[i];
352 
353  number_of_single_dofs[i] = dof_distribution.number_of_single_dofs[i];
354  number_of_pairs[i] = dof_distribution.number_of_pairs[i];
355  number_of_triples[i] = dof_distribution.number_of_triples[i];
356  number_of_sextuples[i] = dof_distribution.number_of_sextuples[i];
357  }
358 
359  for (int i=0; i<dof_distribution.unit_support_points.size(); i++)
360  unit_support_points.push_back(dof_distribution.unit_support_points[i]);
361 
362  order = degree;
363 
364  this->compute_node_matrix();
365 }
366 
367 template<unsigned int degree, unsigned int dim, unsigned int spacedim>
368 double FE_P<degree,dim,spacedim>::basis_value(const unsigned int i, const arma::vec::fixed<dim> &p) const
369 {
370  ASSERT(i <= number_of_dofs, "Index of basis function is out of range.");
371  return poly_space.basis_value(i, p);
372 }
373 
374 template<unsigned int degree, unsigned int dim, unsigned int spacedim>
375 arma::vec::fixed<dim> FE_P<degree,dim,spacedim>::basis_grad(const unsigned int i, const arma::vec::fixed<dim> &p) const
376 {
377  ASSERT(i <= number_of_dofs, "Index of basis function is out of range.");
378  return poly_space.basis_grad(i, p);
379 }
380 
381 template<unsigned int degree, unsigned int dim, unsigned int spacedim>
382 arma::vec::fixed<dim> FE_P<degree,dim,spacedim>::basis_vector(const unsigned int i, const arma::vec::fixed<dim> &p) const
383 {
384  ASSERT(false, "basis_vector() may not be called for scalar finite element.");
385 }
386 
387 template<unsigned int degree, unsigned int dim, unsigned int spacedim>
388 arma::mat::fixed<dim,dim> FE_P<degree,dim,spacedim>::basis_grad_vector(const unsigned int i, const arma::vec::fixed<dim> &p) const
389 {
390  ASSERT(false, "basis_grad_vector() may not be called for scalar finite element.");
391 }
392 
393 
394 
395 
396 
397 
398 
399 
400 
401 
402 
403 
404 
405 
406 
407 template<unsigned int degree, unsigned int dim, unsigned int spacedim>
409 {
410  this->init();
411 
412  number_of_dofs += dof_distribution.number_of_dofs;
413 
414  number_of_single_dofs[dim] = number_of_dofs;
415 
416  for (unsigned int i=0; i<dof_distribution.unit_support_points.size(); i++)
417  unit_support_points.push_back(dof_distribution.unit_support_points[i]);
418 
419  order = degree;
420 
421  this->compute_node_matrix();
422 }
423 
424 template<unsigned int degree, unsigned int dim, unsigned int spacedim>
425 double FE_P_disc<degree,dim,spacedim>::basis_value(const unsigned int i, const arma::vec::fixed<dim> &p) const
426 {
427  ASSERT(i <= number_of_dofs, "Index of basis function is out of range.");
428  return poly_space.basis_value(i, p);
429 }
430 
431 template<unsigned int degree, unsigned int dim, unsigned int spacedim>
432 arma::vec::fixed<dim> FE_P_disc<degree,dim,spacedim>::basis_grad(const unsigned int i, const arma::vec::fixed<dim> &p) const
433 {
434  ASSERT(i <= number_of_dofs, "Index of basis function is out of range.");
435  return poly_space.basis_grad(i, p);
436 }
437 
438 template<unsigned int degree, unsigned int dim, unsigned int spacedim>
439 arma::vec::fixed<dim> FE_P_disc<degree,dim,spacedim>::basis_vector(const unsigned int i, const arma::vec::fixed<dim> &p) const
440 {
441  ASSERT(false, "basis_vector() may not be called for scalar finite element.");
442 }
443 
444 template<unsigned int degree, unsigned int dim, unsigned int spacedim>
445 arma::mat::fixed<dim,dim> FE_P_disc<degree,dim,spacedim>::basis_grad_vector(const unsigned int i, const arma::vec::fixed<dim> &p) const
446 {
447  ASSERT(false, "basis_grad_vector() may not be called for scalar finite element.");
448 }
449 
450 
451 
452 
453 
454 
455 
456 
457 
458 
459 
460 
461 #ifndef DOXYGEN_SHOULD_SKIP_THIS
462 
463 /****** Template specializations ******/
464 
465 /*** 1D finite elements ***/
466 
467 // P0 constant element
468 template<>
470 
471 // P1 linear element
472 template<>
474 
475 // P2 quadratic element
476 template<>
478 
479 // P3 cubic element
480 template<>
482 
483 
484 /*** 2D finite elements ***/
485 
486 // P0 constant element
487 template<>
489 
490 
491 // P1 linear element
492 template<>
494 
495 // P2 quadratic element
496 template<>
498 
499 // P3 cubic element
500 template<>
502 
503 
504 
505 /*** 3D finite elements ***/
506 
507 // P0 constant element
508 template<>
510 
511 
512 // P1 linear element
513 template<>
515 
516 
517 // P2 quadratic element
518 template<>
520 
521 // P3 cubic element
522 template<>
524 
525 #endif /* DOXYGEN_SHOULD_SKIP_THIS */
526 
527 
528 
529 
530 #endif /* FE_P_HH_ */