Flow123d  DF_patch_fe_data_tables-32b3de9
element_values.hh
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2019 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 element_values.hh
15  * @brief Class ElementValues calculates data related to transformation
16  * of reference cell to actual cell (Jacobian, inverse Jacobian,
17  * determinant, point coordinates etc.).
18  * @author Jan Stebel
19  */
20 
21 #ifndef ELEMENT_VALUES_HH_
22 #define ELEMENT_VALUES_HH_
23 
24 // #include <string.h> // for memcpy
25 // #include <algorithm> // for swap
26 // #include <new> // for operator new[]
27 // #include <string> // for operator<<
28 #include <vector> // for vector
29 #include "mesh/ref_element.hh" // for RefElement
30 #include "mesh/accessors.hh" // for ElementAccessor
31 #include "fem/update_flags.hh" // for UpdateFlags
32 #include "fem/dh_cell_accessor.hh" // for DHCellAccessor, DHCellSide
33 
34 class Quadrature;
35 template<unsigned int spacedim> class PatchPointValues;
36 
37 
38 
39 
40 /**
41  * @brief Structure for storing the precomputed element data.
42  */
44 {
45 public:
46 
47  /// Resize vectors to size @p np.
48  RefElementData(unsigned int np);
49 
50  /// Barycentric coordinates of quadrature points.
52 
53  /// Quadrature weights.
55 
56  /// Number of quadrature points.
57  unsigned int n_points;
58 };
59 
60 
61 
62 
63 /**
64  * @brief Class ElementData holds the arrays of data computed by
65  * Mapping.
66  */
67 template<unsigned int spacedim = 3>
69 {
70 public:
71 
72  /**
73  * @brief Resize the data arrays.
74  * @param size Number of quadrature points.
75  * @param flags Update flags to be stores.
76  */
77  ElementData(unsigned int size, UpdateFlags flags, unsigned int dim);
78 
79  /// Print calculated data.
80  void print();
81 
82  /// Dimension of space of reference cell.
83  const unsigned int dim_;
84 
85  /**
86  * @brief Transformed quadrature weights.
87  *
88  * Values at quadrature points of the product of the Jacobian
89  * determinant of the mapping and the weight at the particular
90  * quadrature point.
91  */
93 
94  /// JxW values for sides.
96 
97  /// Jacobians (spacedim x dim) of the mapping at the quadrature points.
99 
100  /// Determinants of Jacobians at quadrature points.
102 
103  /// Inverse Jacobians (dim x spacedim) at the quadrature points.
105 
106  /// Coordinates (spacedim) of quadrature points in the actual cell coordinate system.
108 
109  /// Normal vectors (spacedim) to the element at the quadrature points lying on a side.
111 
112  /// Flags that indicate which finite element quantities are to be computed.
114 
115  /// Iterator to last updated cell.
117 
118  /// Iterator to last updated cell side.
120 
121 };
122 
123 
124 template<unsigned int spacedim>
126 {
127 public:
128  RefElementValues(Quadrature &_quadrature,
129  unsigned int dim);
130 
131  /// Correct deallocation of objects created by 'initialize' methods.
133 
134  /// Returns the number of quadrature points.
135  inline unsigned int n_points()
136  { return n_points_; }
137 
138  /// Initialize ref_data or side_ref_data
139  void ref_initialize(Quadrature &_quadrature, unsigned int dim);
140 
141 
142  /// Dimension of space of reference cell.
143  const unsigned int dim_;
144 
145  /// Number of integration points.
146  const unsigned int n_points_;
147 
148  /// Number of sides in reference cell.
149  const unsigned int n_sides_;
150 
151  /// Data on reference element.
153 
154  /// Data on reference element (for each side ).
156 
157 protected:
158  /// Precompute data on reference element.
160 
161 
162  friend class PatchPointValues<spacedim>;
163 };
164 
165 
166 /**
167  * @brief Class for computation of data on cell and side.
168  */
169 template<unsigned int spacedim>
170 class ElementValues : public RefElementValues<spacedim>
171 {
172 public:
173 
174  /**
175  * @brief Constructor.
176  *
177  * Initializes structures and calculates
178  * cell-independent data. The quadrature can have dimension
179  * either dim or dim-1 (for values on side). In the later
180  * case also normal vectors and side jacobians are computed.
181  *
182  * @param _quadrature The quadrature rule.
183  * @param _flags The update flags.
184  * @param dim Dimension of space of reference cell.
185  */
186  ElementValues(Quadrature &_quadrature,
187  UpdateFlags _flags,
188  unsigned int dim);
189 
190  /// Correct deallocation of objects created by 'initialize' methods.
191  virtual ~ElementValues()
192  {}
193 
194 
195 
196  /**
197  * @brief Update cell-dependent data (gradients, Jacobians etc.)
198  *
199  * @param cell The actual cell.
200  */
202 
203  /**
204  * @brief Update side-dependent data (Jacobians etc.)
205  *
206  * @param cell_side The actual cell and side.
207  */
208  void reinit(const Side &cell_side);
209 
210  /**
211  * @brief Determine quantities to be recomputed on each cell.
212  *
213  * @param flags Flags that indicate what has to be recomputed.
214  */
216 
217  /**
218  * @brief Return the relative volume change of the cell (Jacobian determinant).
219  *
220  * If dim==spacedim then the sign may be negative, otherwise the
221  * result is a positive number.
222  *
223  * @param point_no Number of the quadrature point.
224  */
225  inline double determinant(const unsigned int point_no) const
226  {
227  ASSERT_LT(point_no, this->n_points_);
228  return data.determinants[point_no];
229  }
230 
231  /// Return Jacobian matrix at point @p point_no.
232  inline arma::mat jacobian(const unsigned int point_no) const
233  {
234  ASSERT_LT(point_no, this->n_points_);
235  return data.jacobians.arma_mat(point_no);
236  }
237 
238  /// Return inverse Jacobian matrix at point @p point_no.
239  inline arma::mat inverse_jacobian(const unsigned int point_no) const
240  {
241  ASSERT_LT(point_no, this->n_points_);
242  return data.inverse_jacobians.arma_mat(point_no);
243  }
244 
245  /**
246  * @brief Return the product of Jacobian determinant and the quadrature
247  * weight at given quadrature point.
248  *
249  * @param point_no Number of the quadrature point.
250  */
251  inline double JxW(const unsigned int point_no) const
252  {
253  ASSERT_LT(point_no, this->n_points_);
254  return data.JxW_values[point_no];
255  }
256 
257  /**
258  * @brief Return the product of side Jacobian determinant and the quadrature
259  * weight at given quadrature point.
260  *
261  * @param point_no Number of the quadrature point.
262  */
263  inline double side_JxW(const unsigned int point_no) const
264  {
265  ASSERT_LT(point_no, this->n_points_);
266  return data.side_JxW_values[point_no];
267  }
268 
269  /**
270  * @brief Return coordinates of the quadrature point in the actual cell system.
271  *
272  * @param point_no Number of the quadrature point.
273  */
274  inline arma::vec::fixed<spacedim> point(const unsigned int point_no) const
275  {
276  ASSERT_LT(point_no, this->n_points_);
277  return data.points.template vec<spacedim>(point_no);
278  }
279 
280  /// Return coordinates of all quadrature points in the actual cell system.
281  inline const Armor::array &point_list() const
282  {
283  return data.points;
284  }
285 
286 
287  /**
288  * @brief Returns the normal vector to a side at given quadrature point.
289  *
290  * @param point_no Number of the quadrature point.
291  */
292  inline arma::vec::fixed<spacedim> normal_vector(unsigned int point_no)
293  {
294  ASSERT_LT(point_no, this->n_points_);
295  return data.normal_vectors.template vec<spacedim>(point_no);
296  }
297 
298  /// Return cell at which the values were reinited.
300  { return data.cell; }
301 
302  /// Return cell side where the values were reinited.
303  const Side &side() const
304  { return data.side; }
305 
306 
307 
308 protected:
309  /// Compute data from reference cell and using MappingP1.
310  template<unsigned int dim>
311  void fill_data();
312 
313  /// Calculates the mapping data on a side of a cell.
314  template<unsigned int dim>
315  void fill_side_data();
316 
317 
318  /// Data computed by the mapping.
320 
321 };
322 
323 #endif /* ELEMENT_VALUES_HH_ */
#define ASSERT_LT(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:301
arma::mat arma_mat(uint i) const
Definition: armor.hh:896
Class ElementData holds the arrays of data computed by Mapping.
UpdateFlags update_flags
Flags that indicate which finite element quantities are to be computed.
Armor::array inverse_jacobians
Inverse Jacobians (dim x spacedim) at the quadrature points.
std::vector< double > determinants
Determinants of Jacobians at quadrature points.
Armor::array jacobians
Jacobians (spacedim x dim) of the mapping at the quadrature points.
ElementAccessor< spacedim > cell
Iterator to last updated cell.
ElementData(unsigned int size, UpdateFlags flags, unsigned int dim)
Resize the data arrays.
void print()
Print calculated data.
std::vector< double > JxW_values
Transformed quadrature weights.
Armor::array normal_vectors
Normal vectors (spacedim) to the element at the quadrature points lying on a side.
Armor::array points
Coordinates (spacedim) of quadrature points in the actual cell coordinate system.
std::vector< double > side_JxW_values
JxW values for sides.
const unsigned int dim_
Dimension of space of reference cell.
Side side
Iterator to last updated cell side.
Class for computation of data on cell and side.
const ElementAccessor< spacedim > & cell() const
Return cell at which the values were reinited.
arma::vec::fixed< spacedim > normal_vector(unsigned int point_no)
Returns the normal vector to a side at given quadrature point.
arma::vec::fixed< spacedim > point(const unsigned int point_no) const
Return coordinates of the quadrature point in the actual cell system.
virtual ~ElementValues()
Correct deallocation of objects created by 'initialize' methods.
arma::mat jacobian(const unsigned int point_no) const
Return Jacobian matrix at point point_no.
double determinant(const unsigned int point_no) const
Return the relative volume change of the cell (Jacobian determinant).
double side_JxW(const unsigned int point_no) const
Return the product of side Jacobian determinant and the quadrature weight at given quadrature point.
ElementValues(Quadrature &_quadrature, UpdateFlags _flags, unsigned int dim)
Constructor.
const Armor::array & point_list() const
Return coordinates of all quadrature points in the actual cell system.
arma::mat inverse_jacobian(const unsigned int point_no) const
Return inverse Jacobian matrix at point point_no.
const Side & side() const
Return cell side where the values were reinited.
ElementData< spacedim > data
Data computed by the mapping.
void fill_side_data()
Calculates the mapping data on a side of a cell.
void reinit(const ElementAccessor< spacedim > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
UpdateFlags update_each(UpdateFlags flags)
Determine quantities to be recomputed on each cell.
double JxW(const unsigned int point_no) const
Return the product of Jacobian determinant and the quadrature weight at given quadrature point.
void fill_data()
Compute data from reference cell and using MappingP1.
Class for storing FE data of quadrature points.
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
Structure for storing the precomputed element data.
RefElementData(unsigned int np)
Resize vectors to size np.
unsigned int n_points
Number of quadrature points.
std::vector< arma::vec > bar_coords
Barycentric coordinates of quadrature points.
std::vector< double > weights
Quadrature weights.
void ref_initialize(Quadrature &_quadrature, unsigned int dim)
Initialize ref_data or side_ref_data.
RefElementData * ref_data
Data on reference element.
const unsigned int n_points_
Number of integration points.
std::vector< RefElementData * > side_ref_data
Data on reference element (for each side ).
const unsigned int n_sides_
Number of sides in reference cell.
RefElementData * init_ref_data(const Quadrature &q)
Precompute data on reference element.
unsigned int n_points()
Returns the number of quadrature points.
const unsigned int dim_
Dimension of space of reference cell.
~RefElementValues()
Correct deallocation of objects created by 'initialize' methods.
RefElementValues(Quadrature &_quadrature, unsigned int dim)
ArmaMat< double, N, M > mat
Definition: armor.hh:936
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc.
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
Definition: update_flags.hh:68