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