Flow123d
DF_patch_fe_data_tables-32b3de9
|
Class for computation of data on cell and side. More...
#include <element_values.hh>
Public Member Functions | |
ElementValues (Quadrature &_quadrature, UpdateFlags _flags, unsigned int dim) | |
Constructor. More... | |
virtual | ~ElementValues () |
Correct deallocation of objects created by 'initialize' methods. More... | |
void | reinit (const ElementAccessor< spacedim > &cell) |
Update cell-dependent data (gradients, Jacobians etc.) More... | |
void | reinit (const Side &cell_side) |
Update side-dependent data (Jacobians etc.) More... | |
UpdateFlags | update_each (UpdateFlags flags) |
Determine quantities to be recomputed on each cell. More... | |
double | determinant (const unsigned int point_no) const |
Return the relative volume change of the cell (Jacobian determinant). More... | |
arma::mat | jacobian (const unsigned int point_no) const |
Return Jacobian matrix at point point_no . More... | |
arma::mat | inverse_jacobian (const unsigned int point_no) const |
Return inverse Jacobian matrix at point point_no . More... | |
double | JxW (const unsigned int point_no) const |
Return the product of Jacobian determinant and the quadrature weight at given quadrature point. More... | |
double | side_JxW (const unsigned int point_no) const |
Return the product of side Jacobian determinant and the quadrature weight at given quadrature point. More... | |
arma::vec::fixed< spacedim > | point (const unsigned int point_no) const |
Return coordinates of the quadrature point in the actual cell system. More... | |
const Armor::array & | point_list () const |
Return coordinates of all quadrature points in the actual cell system. More... | |
arma::vec::fixed< spacedim > | normal_vector (unsigned int point_no) |
Returns the normal vector to a side at given quadrature point. More... | |
const ElementAccessor< spacedim > & | cell () const |
Return cell at which the values were reinited. More... | |
const Side & | side () const |
Return cell side where the values were reinited. More... | |
Public Member Functions inherited from RefElementValues< spacedim > | |
RefElementValues (Quadrature &_quadrature, unsigned int dim) | |
~RefElementValues () | |
Correct deallocation of objects created by 'initialize' methods. More... | |
unsigned int | n_points () |
Returns the number of quadrature points. More... | |
void | ref_initialize (Quadrature &_quadrature, unsigned int dim) |
Initialize ref_data or side_ref_data. More... | |
Protected Member Functions | |
template<unsigned int dim> | |
void | fill_data () |
Compute data from reference cell and using MappingP1. More... | |
template<unsigned int dim> | |
void | fill_side_data () |
Calculates the mapping data on a side of a cell. More... | |
Protected Member Functions inherited from RefElementValues< spacedim > | |
RefElementData * | init_ref_data (const Quadrature &q) |
Precompute data on reference element. More... | |
Protected Attributes | |
ElementData< spacedim > | data |
Data computed by the mapping. More... | |
Additional Inherited Members | |
Public Attributes inherited from RefElementValues< spacedim > | |
const unsigned int | dim_ |
Dimension of space of reference cell. More... | |
const unsigned int | n_points_ |
Number of integration points. More... | |
const unsigned int | n_sides_ |
Number of sides in reference cell. More... | |
RefElementData * | ref_data |
Data on reference element. More... | |
std::vector< RefElementData * > | side_ref_data |
Data on reference element (for each side ). More... | |
Class for computation of data on cell and side.
Definition at line 170 of file element_values.hh.
ElementValues< spacedim >::ElementValues | ( | Quadrature & | _quadrature, |
UpdateFlags | _flags, | ||
unsigned int | dim | ||
) |
Constructor.
Initializes structures and calculates cell-independent data. The quadrature can have dimension either dim or dim-1 (for values on side). In the later case also normal vectors and side jacobians are computed.
_quadrature | The quadrature rule. |
_flags | The update flags. |
dim | Dimension of space of reference cell. |
Definition at line 174 of file element_values.cc.
|
inlinevirtual |
Correct deallocation of objects created by 'initialize' methods.
Definition at line 191 of file element_values.hh.
|
inline |
Return cell at which the values were reinited.
Definition at line 299 of file element_values.hh.
|
inline |
Return the relative volume change of the cell (Jacobian determinant).
If dim==spacedim then the sign may be negative, otherwise the result is a positive number.
point_no | Number of the quadrature point. |
Definition at line 225 of file element_values.hh.
|
protected |
Compute data from reference cell and using MappingP1.
Definition at line 290 of file element_values.cc.
|
protected |
Calculates the mapping data on a side of a cell.
Definition at line 376 of file element_values.cc.
|
inline |
Return inverse Jacobian matrix at point point_no
.
Definition at line 239 of file element_values.hh.
|
inline |
Return Jacobian matrix at point point_no
.
Definition at line 232 of file element_values.hh.
|
inline |
Return the product of Jacobian determinant and the quadrature weight at given quadrature point.
point_no | Number of the quadrature point. |
Definition at line 251 of file element_values.hh.
|
inline |
Returns the normal vector to a side at given quadrature point.
point_no | Number of the quadrature point. |
Definition at line 292 of file element_values.hh.
|
inline |
Return coordinates of the quadrature point in the actual cell system.
point_no | Number of the quadrature point. |
Definition at line 274 of file element_values.hh.
|
inline |
Return coordinates of all quadrature points in the actual cell system.
Definition at line 281 of file element_values.hh.
void ElementValues< spacedim >::reinit | ( | const ElementAccessor< spacedim > & | cell | ) |
Update cell-dependent data (gradients, Jacobians etc.)
cell | The actual cell. |
Definition at line 232 of file element_values.cc.
void ElementValues< spacedim >::reinit | ( | const Side & | cell_side | ) |
Update side-dependent data (Jacobians etc.)
cell_side | The actual cell and side. |
Definition at line 261 of file element_values.cc.
|
inline |
Return cell side where the values were reinited.
Definition at line 303 of file element_values.hh.
|
inline |
Return the product of side Jacobian determinant and the quadrature weight at given quadrature point.
point_no | Number of the quadrature point. |
Definition at line 263 of file element_values.hh.
UpdateFlags ElementValues< spacedim >::update_each | ( | UpdateFlags | flags | ) |
Determine quantities to be recomputed on each cell.
flags | Flags that indicate what has to be recomputed. |
Definition at line 148 of file element_values.cc.
|
protected |
Data computed by the mapping.
Definition at line 319 of file element_values.hh.