Flow123d  master-f44eb46
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
ElementValues< spacedim > Class Template Reference

Class for computation of data on cell and side. More...

#include <element_values.hh>

Collaboration diagram for ElementValues< spacedim >:
Collaboration graph
[legend]

Public Member Functions

 ElementValues (Quadrature &_quadrature, UpdateFlags _flags, unsigned int dim)
 Constructor. More...
 
 ~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::arraypoint_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...
 
unsigned int n_points ()
 Returns the number of quadrature points. More...
 
const ElementAccessor< spacedim > & cell () const
 Return cell at which the values were reinited. More...
 
const Sideside () const
 Return cell side where the values were reinited. More...
 

Protected Member Functions

RefElementDatainit_ref_data (const Quadrature &q)
 Precompute data on reference element. More...
 
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 Attributes

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...
 
RefElementDataref_data
 Data on reference element. More...
 
std::vector< RefElementData * > side_ref_data
 Data on reference element (for each side ). More...
 
ElementData< spacedim > data
 Data computed by the mapping. More...
 

Detailed Description

template<unsigned int spacedim>
class ElementValues< spacedim >

Class for computation of data on cell and side.

Definition at line 128 of file element_values.hh.

Constructor & Destructor Documentation

◆ ElementValues()

template<unsigned int spacedim>
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.

Parameters
_quadratureThe quadrature rule.
_flagsThe update flags.
dimDimension of space of reference cell.

Definition at line 159 of file element_values.cc.

◆ ~ElementValues()

template<unsigned int spacedim>
ElementValues< spacedim >::~ElementValues

Correct deallocation of objects created by 'initialize' methods.

Definition at line 205 of file element_values.cc.

Member Function Documentation

◆ cell()

template<unsigned int spacedim>
const ElementAccessor<spacedim>& ElementValues< spacedim >::cell ( ) const
inline

Return cell at which the values were reinited.

Definition at line 260 of file element_values.hh.

◆ determinant()

template<unsigned int spacedim>
double ElementValues< spacedim >::determinant ( const unsigned int  point_no) const
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.

Parameters
point_noNumber of the quadrature point.

Definition at line 182 of file element_values.hh.

Here is the caller graph for this function:

◆ fill_data()

template<unsigned int spacedim>
template<unsigned int dim>
void ElementValues< spacedim >::fill_data
protected

Compute data from reference cell and using MappingP1.

Definition at line 273 of file element_values.cc.

◆ fill_side_data()

template<unsigned int spacedim>
template<unsigned int dim>
void ElementValues< spacedim >::fill_side_data
protected

Calculates the mapping data on a side of a cell.

Definition at line 349 of file element_values.cc.

◆ init_ref_data()

template<unsigned int spacedim>
RefElementData * ElementValues< spacedim >::init_ref_data ( const Quadrature q)
protected

Precompute data on reference element.

Definition at line 97 of file element_values.cc.

Here is the caller graph for this function:

◆ inverse_jacobian()

template<unsigned int spacedim>
arma::mat ElementValues< spacedim >::inverse_jacobian ( const unsigned int  point_no) const
inline

Return inverse Jacobian matrix at point point_no.

Definition at line 196 of file element_values.hh.

Here is the caller graph for this function:

◆ jacobian()

template<unsigned int spacedim>
arma::mat ElementValues< spacedim >::jacobian ( const unsigned int  point_no) const
inline

Return Jacobian matrix at point point_no.

Definition at line 189 of file element_values.hh.

Here is the caller graph for this function:

◆ JxW()

template<unsigned int spacedim>
double ElementValues< spacedim >::JxW ( const unsigned int  point_no) const
inline

Return the product of Jacobian determinant and the quadrature weight at given quadrature point.

Parameters
point_noNumber of the quadrature point.

Definition at line 208 of file element_values.hh.

◆ n_points()

template<unsigned int spacedim>
unsigned int ElementValues< spacedim >::n_points ( )
inline

Returns the number of quadrature points.

Definition at line 256 of file element_values.hh.

◆ normal_vector()

template<unsigned int spacedim>
arma::vec::fixed<spacedim> ElementValues< spacedim >::normal_vector ( unsigned int  point_no)
inline

Returns the normal vector to a side at given quadrature point.

Parameters
point_noNumber of the quadrature point.

Definition at line 249 of file element_values.hh.

◆ point()

template<unsigned int spacedim>
arma::vec::fixed<spacedim> ElementValues< spacedim >::point ( const unsigned int  point_no) const
inline

Return coordinates of the quadrature point in the actual cell system.

Parameters
point_noNumber of the quadrature point.

Definition at line 231 of file element_values.hh.

◆ point_list()

template<unsigned int spacedim>
const Armor::array& ElementValues< spacedim >::point_list ( ) const
inline

Return coordinates of all quadrature points in the actual cell system.

Definition at line 238 of file element_values.hh.

◆ reinit() [1/2]

template<unsigned int spacedim>
void ElementValues< spacedim >::reinit ( const ElementAccessor< spacedim > &  cell)

Update cell-dependent data (gradients, Jacobians etc.)

Parameters
cellThe actual cell.

Definition at line 215 of file element_values.cc.

◆ reinit() [2/2]

template<unsigned int spacedim>
void ElementValues< spacedim >::reinit ( const Side cell_side)

Update side-dependent data (Jacobians etc.)

Parameters
cell_sideThe actual cell and side.

Definition at line 244 of file element_values.cc.

◆ side()

template<unsigned int spacedim>
const Side& ElementValues< spacedim >::side ( ) const
inline

Return cell side where the values were reinited.

Definition at line 264 of file element_values.hh.

◆ side_JxW()

template<unsigned int spacedim>
double ElementValues< spacedim >::side_JxW ( const unsigned int  point_no) const
inline

Return the product of side Jacobian determinant and the quadrature weight at given quadrature point.

Parameters
point_noNumber of the quadrature point.

Definition at line 220 of file element_values.hh.

◆ update_each()

template<unsigned int spacedim>
UpdateFlags ElementValues< spacedim >::update_each ( UpdateFlags  flags)

Determine quantities to be recomputed on each cell.

Parameters
flagsFlags that indicate what has to be recomputed.

Definition at line 133 of file element_values.cc.

Member Data Documentation

◆ data

template<unsigned int spacedim>
ElementData<spacedim> ElementValues< spacedim >::data
protected

Data computed by the mapping.

Definition at line 300 of file element_values.hh.

◆ dim_

template<unsigned int spacedim>
const unsigned int ElementValues< spacedim >::dim_
protected

Dimension of space of reference cell.

Definition at line 285 of file element_values.hh.

◆ n_points_

template<unsigned int spacedim>
const unsigned int ElementValues< spacedim >::n_points_
protected

Number of integration points.

Definition at line 288 of file element_values.hh.

◆ n_sides_

template<unsigned int spacedim>
const unsigned int ElementValues< spacedim >::n_sides_
protected

Number of sides in reference cell.

Definition at line 291 of file element_values.hh.

◆ ref_data

template<unsigned int spacedim>
RefElementData* ElementValues< spacedim >::ref_data
protected

Data on reference element.

Definition at line 294 of file element_values.hh.

◆ side_ref_data

template<unsigned int spacedim>
std::vector<RefElementData*> ElementValues< spacedim >::side_ref_data
protected

Data on reference element (for each side ).

Definition at line 297 of file element_values.hh.


The documentation for this class was generated from the following files: