Go to the documentation of this file.
21 #ifndef ELEMENT_VALUES_HH_
22 #define ELEMENT_VALUES_HH_
66 template<
unsigned int spacedim = 3>
127 template<
unsigned int spacedim>
185 return data.determinants[point_no];
192 return data.jacobians.arma_mat(point_no);
199 return data.inverse_jacobians.arma_mat(point_no);
208 inline double JxW(
const unsigned int point_no)
const
211 return data.JxW_values[point_no];
220 inline double side_JxW(
const unsigned int point_no)
const
223 return data.side_JxW_values[point_no];
231 inline arma::vec::fixed<spacedim>
point(
const unsigned int point_no)
const
234 return data.points.template vec<spacedim>(point_no);
252 return data.normal_vectors.template vec<spacedim>(point_no);
261 {
return data.cell; }
265 {
return data.side; }
275 template<
unsigned int dim>
279 template<
unsigned int dim>
UpdateFlags update_flags
Flags that indicate which finite element quantities are to be computed.
const unsigned int dim_
Dimension of space of reference cell.
UpdateFlags update_each(UpdateFlags flags)
Determine quantities to be recomputed on each cell.
std::vector< double > weights
Quadrature weights.
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc.
Armor::array points
Coordinates (spacedim) of quadrature points in the actual cell coordinate system.
const ElementAccessor< spacedim > & cell() const
Return cell at which the values were reinited.
RefElementData(unsigned int np)
Resize vectors to size np.
RefElementData * init_ref_data(const Quadrature &q)
Precompute data on reference element.
Structure for storing the precomputed element data.
unsigned int n_points()
Returns the number of quadrature points.
std::vector< double > side_JxW_values
JxW values for sides.
std::vector< RefElementData * > side_ref_data
Data on reference element (for each side ).
RefElementData * ref_data
Data on reference element.
const unsigned int dim_
Dimension of space of reference cell.
unsigned int n_points
Number of quadrature points.
void fill_data()
Compute data from reference cell and using MappingP1.
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()
Correct deallocation of objects created by 'initialize' methods.
Class ElementData holds the arrays of data computed by Mapping.
ElementValues(Quadrature &_quadrature, UpdateFlags _flags, unsigned int dim)
Constructor.
#define ASSERT_LT(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
std::vector< double > JxW_values
Transformed quadrature weights.
std::vector< double > determinants
Determinants of Jacobians at quadrature points.
void reinit(const ElementAccessor< spacedim > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
const unsigned int n_points_
Number of integration points.
ArmaMat< double, N, M > mat
const Side & side() const
Return cell side where the values were reinited.
Side side
Iterator to last updated cell side.
void fill_side_data()
Calculates the mapping data on a side of a cell.
Armor::array jacobians
Jacobians (spacedim x dim) of the mapping at the quadrature points.
ElementAccessor< spacedim > cell
Iterator to last updated cell.
arma::mat inverse_jacobian(const unsigned int point_no) const
Return inverse Jacobian matrix at point point_no.
arma::mat jacobian(const unsigned int point_no) const
Return Jacobian matrix at point point_no.
Class for computation of data on cell and side.
const unsigned int n_sides_
Number of sides in reference cell.
double JxW(const unsigned int point_no) const
Return the product of Jacobian determinant and the quadrature weight at given quadrature point.
Armor::array inverse_jacobians
Inverse Jacobians (dim x spacedim) at the quadrature points.
void print()
Print calculated data.
Armor::array normal_vectors
Normal vectors (spacedim) to the element at the quadrature points lying on a side.
ElementData(unsigned int size, UpdateFlags flags, unsigned int dim)
Resize the data arrays.
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
arma::vec::fixed< spacedim > point(const unsigned int point_no) const
Return coordinates of the quadrature point in the actual cell system.
double determinant(const unsigned int point_no) const
Return the relative volume change of the cell (Jacobian determinant).
arma::vec::fixed< spacedim > normal_vector(unsigned int point_no)
Returns the normal vector to a side at given quadrature point.
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
ElementData< spacedim > data
Data computed by the mapping.
Base class for quadrature rules on simplices in arbitrary dimensions.
std::vector< arma::vec > bar_coords
Barycentric coordinates of quadrature points.
const Armor::array & point_list() const
Return coordinates of all quadrature points in the actual cell system.