Flow123d  release_2.2.0-914-gf1a3a4f
fe_values_views.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_system.cc
15  * @brief Class FESystem for compound finite elements.
16  * @author Jan Stebel
17  */
18 
19 #include "fem/fe_values_views.hh"
20 #include "fem/fe_values.hh"
21 #include "fem/finite_element.hh"
22 #include "quadrature/quadrature.hh"
23 
24 using namespace FEValuesViews;
25 
26 template<unsigned int dim, unsigned int spacedim>
27 double FEValuesViews::Scalar<dim,spacedim>::value(unsigned int function_no, unsigned int point_no) const
28 {
29  ASSERT_LT_DBG( function_no, fe_values_.n_dofs() );
30  ASSERT_LT_DBG( point_no, fe_values_.n_points() );
31  return fe_values_.shape_value_component(function_no, point_no, component_);
32 }
33 
34 template<unsigned int dim, unsigned int spacedim>
35 arma::vec::fixed<spacedim> FEValuesViews::Scalar<dim,spacedim>::grad(unsigned int function_no, unsigned int point_no) const
36 {
37  ASSERT_LT_DBG( function_no, fe_values_.n_dofs() );
38  ASSERT_LT_DBG( point_no, fe_values_.n_points() );
39  return fe_values_.shape_grad_component(function_no, point_no, component_);
40 }
41 
42 template<unsigned int dim, unsigned int spacedim>
44 { return fe_values_; }
45 
46 
47 
48 
49 template<unsigned int dim, unsigned int spacedim>
50 arma::vec::fixed<spacedim> FEValuesViews::Vector<dim,spacedim>::value(unsigned int function_no, unsigned int point_no) const
51 {
52  ASSERT_LT_DBG( function_no, fe_values_.n_dofs() );
53  ASSERT_LT_DBG( point_no, fe_values_.n_points() );
54  arma::vec::fixed<spacedim> v;
55  for (unsigned int c=0; c<spacedim; ++c)
56  v(c) = fe_values_.shape_value_component(function_no, point_no, first_vector_component_+c);
57  return v;
58 }
59 
60 template<unsigned int dim, unsigned int spacedim>
61 arma::mat::fixed<spacedim,spacedim> FEValuesViews::Vector<dim,spacedim>::grad(unsigned int function_no, unsigned int point_no) const
62 {
63  ASSERT_LT_DBG( function_no, fe_values_.n_dofs() );
64  ASSERT_LT_DBG( point_no, fe_values_.n_points() );
65  arma::mat::fixed<spacedim,spacedim> g;
66  for (unsigned int c=0; c<spacedim; ++c)
67  g.col(c) = fe_values_.shape_grad_component(function_no, point_no, first_vector_component_+c);
68  return g.t();
69 }
70 
71 template<unsigned int dim, unsigned int spacedim>
72 arma::mat::fixed<spacedim,spacedim> FEValuesViews::Vector<dim,spacedim>::sym_grad(unsigned int function_no, unsigned int point_no) const
73 {
74  ASSERT_LT_DBG( function_no, fe_values_.n_dofs() );
75  ASSERT_LT_DBG( point_no, fe_values_.n_points() );
76  arma::mat::fixed<spacedim,spacedim> g = grad(function_no, point_no);
77  return 0.5*(g+trans(g));
78 }
79 
80 template<unsigned int dim, unsigned int spacedim>
81 double FEValuesViews::Vector<dim,spacedim>::divergence(unsigned int function_no, unsigned int point_no) const
82 {
83  ASSERT_LT_DBG( function_no, fe_values_.n_dofs() );
84  ASSERT_LT_DBG( point_no, fe_values_.n_points() );
85  double div = 0;
86  for (unsigned int c=0; c<spacedim; ++c)
87  div += fe_values_.shape_grad_component(function_no, point_no, first_vector_component_+c)(first_vector_component_+c);
88  return div;
89 }
90 
91 template<unsigned int dim, unsigned int spacedim>
93 { return fe_values_; }
94 
95 
96 
97 
98 
99 template class FEValuesViews::Scalar<1,3>;
100 template class FEValuesViews::Scalar<2,3>;
101 template class FEValuesViews::Scalar<3,3>;
102 
103 template class FEValuesViews::Vector<1,3>;
104 template class FEValuesViews::Vector<2,3>;
105 template class FEValuesViews::Vector<3,3>;
106 
107 
FEValuesBase< dim, spacedim > & base() const
Returns the FEValuesBase class.
Class FEValues calculates finite element data on the actual cells such as shape function values...
double value(unsigned int function_no, unsigned int point_no) const
Return value of scalar shape function.
arma::mat::fixed< spacedim, spacedim > grad(unsigned int function_no, unsigned int point_no) const
Return gradient of vector-valued shape function.
Basic definitions of numerical quadrature rules.
FEValuesBase< dim, spacedim > & base() const
Returns the FEValuesBase class.
double divergence(unsigned int function_no, unsigned int point_no) const
Return divergence of vector-valued shape function.
arma::vec::fixed< spacedim > value(unsigned int function_no, unsigned int point_no) const
Return value of vector-valued shape function.
arma::vec::fixed< spacedim > grad(unsigned int function_no, unsigned int point_no) const
Return gradient of scalar shape function.
arma::mat::fixed< spacedim, spacedim > sym_grad(unsigned int function_no, unsigned int point_no) const
Return symmetric gradient of vector-valued shape function.
Abstract class for description of finite elements.
Base class for FEValues and FESideValues.
Definition: fe_values.hh:34
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:299