Flow123d  build_with_4.0.3-c0baa07
fe_values_views.hh
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.hh
15  * @brief Class FESystem for compound finite elements.
16  * @author Jan Stebel
17  */
18 
19 #ifndef FE_VALUES_VIEWS_HH
20 #define FE_VALUES_VIEWS_HH
21 
22 #include <vector>
23 #include <armadillo>
24 
25 template<unsigned int spacedim> class FEValues;
26 
27 
28 
29 /**
30  * FEValuesViews classes realize access to shape functions and their derivatives
31  * for scalar, vector, tensor etc. part of the finite element.
32  *
33  * The classes are motivated by deal.II.
34  */
35 namespace FEValuesViews {
36 
37  template<unsigned int spacedim = 3>
38  class Scalar {
39 
40  public:
41 
42  Scalar(const FEValues<spacedim> &fe_values, unsigned int component)
43  : fe_values_(fe_values),
44  component_(component)
45  {};
46 
47  /**
48  * @brief Return value of scalar shape function.
49  * @param function_no Index of shape function within the FE.
50  * @param point_no Index of quadrature point.
51  */
52  double value(unsigned int function_no, unsigned int point_no) const;
53 
54  /**
55  * @brief Return gradient of scalar shape function.
56  * @param function_no Index of shape function within the FE.
57  * @param point_no Index of quadrature point.
58  */
59  arma::vec::fixed<spacedim> grad(unsigned int function_no, unsigned int point_no) const;
60 
61  /// Returns the FEValues class.
62  const FEValues<spacedim> &base() const;
63 
64  private:
65 
66  /// Base FEValues class for access to the FE.
68 
69  /// Index of the scalar component.
70  unsigned int component_;
71  };
72 
73 
74  template<unsigned int spacedim = 3>
75  class Vector {
76 
77  public:
78 
79  Vector(const FEValues<spacedim> &fe_values, unsigned int component)
80  : fe_values_(fe_values),
81  first_vector_component_(component)
82  {};
83 
84  /**
85  * @brief Return value of vector-valued shape function.
86  * @param function_no Index of shape function within the FE.
87  * @param point_no Index of quadrature point.
88  */
89  arma::vec::fixed<spacedim> value(unsigned int function_no, unsigned int point_no) const;
90 
91  /**
92  * @brief Return gradient of vector-valued shape function.
93  * @param function_no Index of shape function within the FE.
94  * @param point_no Index of quadrature point.
95  */
96  arma::mat::fixed<spacedim,spacedim> grad(unsigned int function_no, unsigned int point_no) const;
97 
98  /**
99  * @brief Return symmetric gradient of vector-valued shape function.
100  * @param function_no Index of shape function within the FE.
101  * @param point_no Index of quadrature point.
102  */
103  arma::mat::fixed<spacedim,spacedim> sym_grad(unsigned int function_no, unsigned int point_no) const;
104 
105  /**
106  * @brief Return divergence of vector-valued shape function.
107  * @param function_no Index of shape function within the FE.
108  * @param point_no Index of quadrature point.
109  */
110  double divergence(unsigned int function_no, unsigned int point_no) const;
111 
112  /// Returns the FEValues class.
113  const FEValues<spacedim> &base() const;
114 
115  private:
116 
117  /// Base FEValues class for access to the FE.
119 
120  /// Index of the first component of the vector.
122  };
123 
124 
125  template<unsigned int spacedim = 3>
126  class Tensor {
127 
128  public:
129 
130  Tensor(const FEValues<spacedim> &fe_values, unsigned int component)
131  : fe_values_(fe_values),
132  first_tensor_component_(component)
133  {};
134 
135  /**
136  * @brief Return value of tensor-valued shape function.
137  * @param function_no Index of shape function within the FE.
138  * @param point_no Index of quadrature point.
139  */
140  arma::mat::fixed<spacedim,spacedim> value(unsigned int function_no, unsigned int point_no) const;
141 
142  /**
143  * @brief Return partial derivative of tensor-valued shape function.
144  * @param variable_no Index of spacial variable w.r. to which we differentiate.
145  * @param function_no Index of shape function within the FE.
146  * @param point_no Index of quadrature point.
147  */
148  arma::mat::fixed<spacedim,spacedim> derivative(
149  unsigned int variable_no,
150  unsigned int function_no,
151  unsigned int point_no) const;
152 
153  /**
154  * @brief Return divergence of tensor-valued shape function.
155  *
156  * The result is a vector whose components are divergences of tensor columns, i.e.
157  * (div T)_i = dT_ji / dx_j.
158  *
159  * @param function_no Index of shape function within the FE.
160  * @param point_no Index of quadrature point.
161  */
162  arma::vec::fixed<spacedim> divergence(unsigned int function_no, unsigned int point_no) const;
163 
164  /// Returns the FEValues class.
165  const FEValues<spacedim> &base() const;
166 
167  private:
168 
169  /// Base FEValues class for access to the FE.
171 
172  /// Index of the first component of the vector.
174 
175  };
176 
177 }
178 
179 
180 
181 
182 
183 
184 #endif
unsigned int component_
Index of the scalar component.
const FEValues< spacedim > & base() const
Returns the FEValues class.
const FEValues< spacedim > & fe_values_
Base FEValues class for access to the FE.
arma::vec::fixed< spacedim > grad(unsigned int function_no, unsigned int point_no) const
Return gradient of scalar shape function.
Scalar(const FEValues< spacedim > &fe_values, unsigned int component)
double value(unsigned int function_no, unsigned int point_no) const
Return value of scalar shape function.
const FEValues< spacedim > & fe_values_
Base FEValues class for access to the FE.
unsigned int first_tensor_component_
Index of the first component of the vector.
Tensor(const FEValues< spacedim > &fe_values, unsigned int component)
arma::mat::fixed< spacedim, spacedim > value(unsigned int function_no, unsigned int point_no) const
Return value of tensor-valued shape function.
const FEValues< spacedim > & base() const
Returns the FEValues class.
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::vec::fixed< spacedim > divergence(unsigned int function_no, unsigned int point_no) const
Return divergence of tensor-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.
arma::mat::fixed< spacedim, spacedim > sym_grad(unsigned int function_no, unsigned int point_no) const
Return symmetric gradient of vector-valued shape function.
Vector(const FEValues< spacedim > &fe_values, unsigned int component)
const FEValues< spacedim > & base() const
Returns the FEValues class.
arma::vec::fixed< spacedim > value(unsigned int function_no, unsigned int point_no) const
Return value of vector-valued shape function.
unsigned int first_vector_component_
Index of the first component of the vector.
double divergence(unsigned int function_no, unsigned int point_no) const
Return divergence of vector-valued shape function.
const FEValues< spacedim > & fe_values_
Base FEValues class for access to the FE.
Calculates finite element data on the actual cell.
Definition: fe_values.hh:67