Flow123d  JS_before_hm-1822-gbb48b12e9
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 "mesh/accessors.hh"
20 #include "fem/fe_values_views.hh"
21 #include "fem/fe_values.hh"
22 #include "fem/finite_element.hh"
23 #include "quadrature/quadrature.hh"
24 
25 using namespace FEValuesViews;
26 
27 template<unsigned int spacedim>
28 double FEValuesViews::Scalar<spacedim>::value(unsigned int function_no, unsigned int point_no) const
29 {
30  ASSERT_LT_DBG( function_no, fe_values_.n_dofs() );
31  ASSERT_LT_DBG( point_no, fe_values_.n_points() );
32  return fe_values_.shape_value_component(function_no, point_no, component_);
33 }
34 
35 template<unsigned int spacedim>
36 arma::vec::fixed<spacedim> FEValuesViews::Scalar<spacedim>::grad(unsigned int function_no, unsigned int point_no) const
37 {
38  ASSERT_LT_DBG( function_no, fe_values_.n_dofs() );
39  ASSERT_LT_DBG( point_no, fe_values_.n_points() );
40  return fe_values_.shape_grad_component(function_no, point_no, component_);
41 }
42 
43 template<unsigned int spacedim>
45 { return fe_values_; }
46 
47 
48 
49 
50 template<unsigned int spacedim>
51 arma::vec::fixed<spacedim> FEValuesViews::Vector<spacedim>::value(unsigned int function_no, unsigned int point_no) const
52 {
53  ASSERT_LT_DBG( function_no, fe_values_.n_dofs() );
54  ASSERT_LT_DBG( point_no, fe_values_.n_points() );
55  arma::vec::fixed<spacedim> v;
56  for (unsigned int c=0; c<spacedim; ++c)
57  v(c) = fe_values_.shape_value_component(function_no, point_no, first_vector_component_+c);
58  return v;
59 }
60 
61 template<unsigned int spacedim>
62 arma::mat::fixed<spacedim,spacedim> FEValuesViews::Vector<spacedim>::grad(unsigned int function_no, unsigned int point_no) const
63 {
64  ASSERT_LT_DBG( function_no, fe_values_.n_dofs() );
65  ASSERT_LT_DBG( point_no, fe_values_.n_points() );
66  arma::mat::fixed<spacedim,spacedim> g;
67  for (unsigned int c=0; c<spacedim; ++c)
68  g.col(c) = fe_values_.shape_grad_component(function_no, point_no, first_vector_component_+c);
69  return g.t();
70 }
71 
72 template<unsigned int spacedim>
73 arma::mat::fixed<spacedim,spacedim> FEValuesViews::Vector<spacedim>::sym_grad(unsigned int function_no, unsigned int point_no) const
74 {
75  ASSERT_LT_DBG( function_no, fe_values_.n_dofs() );
76  ASSERT_LT_DBG( point_no, fe_values_.n_points() );
77  arma::mat::fixed<spacedim,spacedim> g = grad(function_no, point_no);
78  return 0.5*(g+trans(g));
79 }
80 
81 template<unsigned int spacedim>
82 double FEValuesViews::Vector<spacedim>::divergence(unsigned int function_no, unsigned int point_no) const
83 {
84  ASSERT_LT_DBG( function_no, fe_values_.n_dofs() );
85  ASSERT_LT_DBG( point_no, fe_values_.n_points() );
86  double div = 0;
87  for (unsigned int c=0; c<spacedim; ++c)
88  div += fe_values_.shape_grad_component(function_no, point_no, first_vector_component_+c)(first_vector_component_+c);
89  return div;
90 }
91 
92 template<unsigned int spacedim>
94 { return fe_values_; }
95 
96 
97 
98 template<unsigned int spacedim>
99 arma::mat::fixed<spacedim,spacedim> FEValuesViews::Tensor<spacedim>::value(unsigned int function_no, unsigned int point_no) const
100 {
101  ASSERT_LT_DBG( function_no, fe_values_.n_dofs() );
102  ASSERT_LT_DBG( point_no, fe_values_.n_points() );
103  arma::mat::fixed<spacedim,spacedim> v;
104  for (unsigned int c=0; c<spacedim*spacedim; ++c)
105  v(c/spacedim,c%spacedim) = fe_values_.shape_value_component(function_no, point_no, first_tensor_component_+c);
106  return v;
107 }
108 
109 template<unsigned int spacedim>
110 arma::mat::fixed<spacedim,spacedim> FEValuesViews::Tensor<spacedim>::derivative(
111  unsigned int variable_no,
112  unsigned int function_no,
113  unsigned int point_no) const
114 {
115  ASSERT_LT_DBG( variable_no, spacedim );
116  ASSERT_LT_DBG( function_no, fe_values_.n_dofs() );
117  ASSERT_LT_DBG( point_no, fe_values_.n_points() );
118  arma::mat::fixed<spacedim,spacedim> g;
119  for (unsigned int c=0; c<spacedim*spacedim; ++c)
120  g(c/spacedim,c%spacedim) = fe_values_.shape_grad_component(function_no, point_no, first_tensor_component_+c)[first_tensor_component_+variable_no];
121  return g;
122 }
123 
124 template<unsigned int spacedim>
125 arma::vec::fixed<spacedim> FEValuesViews::Tensor<spacedim>::divergence(unsigned int function_no, unsigned int point_no) const
126 {
127  ASSERT_LT_DBG( function_no, fe_values_.n_dofs() );
128  ASSERT_LT_DBG( point_no, fe_values_.n_points() );
129  arma::vec::fixed<spacedim> div;
130  div.zeros();
131  for (unsigned int c=0; c<spacedim*spacedim; ++c)
132  div(c%spacedim) += fe_values_.shape_grad_component(function_no, point_no, first_tensor_component_+c)[first_tensor_component_+c/spacedim];
133  return div;
134 }
135 
136 template<unsigned int spacedim>
138 { return fe_values_; }
139 
140 
141 
142 
143 
144 template class FEValuesViews::Scalar<3>;
145 template class FEValuesViews::Vector<3>;
146 template class FEValuesViews::Tensor<3>;
147 
148 
FEValuesViews
Definition: fe_values_views.hh:35
FEValuesViews::Vector::divergence
double divergence(unsigned int function_no, unsigned int point_no) const
Return divergence of vector-valued shape function.
Definition: fe_values_views.cc:82
fe_values_views.hh
FEValuesViews::Scalar
Definition: fe_values_views.hh:38
FEValuesViews::Vector< 3 >
FEValuesViews::Scalar::grad
arma::vec::fixed< spacedim > grad(unsigned int function_no, unsigned int point_no) const
Return gradient of scalar shape function.
Definition: fe_values_views.cc:36
FEValuesViews::Tensor
Definition: fe_values_views.hh:126
FEValuesViews::Scalar::value
double value(unsigned int function_no, unsigned int point_no) const
Return value of scalar shape function.
Definition: fe_values_views.cc:28
fe_values.hh
Class FEValues calculates finite element data on the actual cells such as shape function values,...
FEValuesViews::Tensor::derivative
arma::mat::fixed< spacedim, spacedim > derivative(unsigned int variable_no, unsigned int function_no, unsigned int point_no) const
Return partial derivative of tensor-valued shape function.
Definition: fe_values_views.cc:110
FEValuesViews::Scalar::base
const FEValues< spacedim > & base() const
Returns the FEValues class.
Definition: fe_values_views.cc:44
FEValues
Calculates finite element data on the actual cell.
Definition: fe_values.hh:66
FEValuesViews::Vector::sym_grad
arma::mat::fixed< spacedim, spacedim > sym_grad(unsigned int function_no, unsigned int point_no) const
Return symmetric gradient of vector-valued shape function.
Definition: fe_values_views.cc:73
accessors.hh
FEValuesViews::Tensor::divergence
arma::vec::fixed< spacedim > divergence(unsigned int function_no, unsigned int point_no) const
Return divergence of tensor-valued shape function.
Definition: fe_values_views.cc:125
FEValuesViews::Vector::base
const FEValues< spacedim > & base() const
Returns the FEValues class.
Definition: fe_values_views.cc:93
finite_element.hh
Abstract class for description of finite elements.
FEValuesViews::Tensor::value
arma::mat::fixed< spacedim, spacedim > value(unsigned int function_no, unsigned int point_no) const
Return value of tensor-valued shape function.
Definition: fe_values_views.cc:99
FEValuesViews::Vector::grad
arma::mat::fixed< spacedim, spacedim > grad(unsigned int function_no, unsigned int point_no) const
Return gradient of vector-valued shape function.
Definition: fe_values_views.cc:62
FEValuesViews::Tensor::base
const FEValues< spacedim > & base() const
Returns the FEValues class.
Definition: fe_values_views.cc:137
quadrature.hh
Basic definitions of numerical quadrature rules.
ASSERT_LT_DBG
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:300
FEValuesViews::Vector::value
arma::vec::fixed< spacedim > value(unsigned int function_no, unsigned int point_no) const
Return value of vector-valued shape function.
Definition: fe_values_views.cc:51