Flow123d
release_2.2.0-914-gf1a3a4f
|
Affine mapping between reference and actual cell. More...
#include <mapping_p1.hh>
Public Types | |
typedef arma::vec::fixed< dim+1 > | BaryPoint |
typedef arma::vec::fixed< spacedim > | RealPoint |
typedef arma::mat::fixed< spacedim, dim+1 > | ElementMap |
Public Member Functions | |
MappingP1 () | |
Constructor. More... | |
MappingInternalData * | initialize (const Quadrature< dim > &q, UpdateFlags flags) |
Initializes the structures and computes static data. More... | |
UpdateFlags | update_each (UpdateFlags flags) |
Determines which additional quantities have to be computed. More... | |
void | fill_fe_values (const typename DOFHandlerBase::CellIterator &cell, const Quadrature< dim > &q, MappingInternalData &data, FEValuesData< dim, spacedim > &fv_data) |
Calculates the mapping data on the actual cell. More... | |
void | fill_fe_side_values (const typename DOFHandlerBase::CellIterator &cell, unsigned int sid, const Quadrature< dim > &q, MappingInternalData &data, FEValuesData< dim, spacedim > &fv_data) |
Calculates the mapping data on a side of a cell. More... | |
ElementMap | element_map (const Element &elm) const |
BaryPoint | project_real_to_unit (const RealPoint &point, const ElementMap &map) const |
RealPoint | project_unit_to_real (const BaryPoint &point, const ElementMap &map) const |
BaryPoint | clip_to_element (BaryPoint &barycentric) |
Public Member Functions inherited from Mapping< dim, spacedim > | |
virtual | ~Mapping () |
Destructor. More... | |
Private Attributes | |
arma::mat::fixed< dim+1, dim > | grad |
Auxiliary matrix of gradients of shape functions (used for computation of the Jacobian). More... | |
Affine mapping between reference and actual cell.
Class MappingP1 implements the affine transformation of the reference cell onto the actual cell.
dim | Dimension of the cells. |
spacedim | Dimension of the Euclidean space. |
Definition at line 54 of file mapping_p1.hh.
typedef arma::vec::fixed<dim+1> MappingP1< dim, spacedim >::BaryPoint |
Definition at line 58 of file mapping_p1.hh.
typedef arma::mat::fixed<spacedim, dim+1> MappingP1< dim, spacedim >::ElementMap |
Definition at line 60 of file mapping_p1.hh.
typedef arma::vec::fixed<spacedim> MappingP1< dim, spacedim >::RealPoint |
Definition at line 59 of file mapping_p1.hh.
Constructor.
Definition at line 31 of file mapping_p1.cc.
auto MappingP1< dim, spacedim >::clip_to_element | ( | BaryPoint & | barycentric | ) |
Clip a point given by barycentric cocordinates to the element. If the point is out of the element the closest point projection to the element surface is used.
Definition at line 297 of file mapping_p1.cc.
auto MappingP1< dim, spacedim >::element_map | ( | const Element & | elm | ) | const |
Map from reference element (barycentric coords) to global coord system. Matrix(3, dim+1) M: x_real = M * x_bary; M columns are real coordinates of nodes.
Definition at line 265 of file mapping_p1.cc.
|
virtual |
Calculates the mapping data on a side of a cell.
cell | The actual cell. |
sid | Number of the side. |
q | The quadrature rule with points on the side. |
data | Precomputed mapping data. |
fv_data | Data to be computed. |
Implements Mapping< dim, spacedim >.
Definition at line 162 of file mapping_p1.cc.
|
virtual |
Calculates the mapping data on the actual cell.
cell | The actual cell. |
q | Quadrature rule. |
data | Precomputed mapping data. |
fv_data | Data to be computed. |
Implements Mapping< dim, spacedim >.
Definition at line 89 of file mapping_p1.cc.
|
virtual |
Initializes the structures and computes static data.
q | Quadrature rule. |
flags | Update flags. |
Implements Mapping< dim, spacedim >.
Definition at line 36 of file mapping_p1.cc.
auto MappingP1< dim, spacedim >::project_real_to_unit | ( | const RealPoint & | point, |
const ElementMap & | map | ||
) | const |
Project given point in real coordinates to reference element (barycentic coordinates). Result vector have dimension dim()+1. Use RefElement<dim>::bary_to_local() to get local coordinates.
Definition at line 275 of file mapping_p1.cc.
auto MappingP1< dim, spacedim >::project_unit_to_real | ( | const BaryPoint & | point, |
const ElementMap & | map | ||
) | const |
Project given point from reference element (barycentic coordinates) to real coordinates. Use RefElement<dim>::local_to_bary() to get barycentric coordinates in input.
Definition at line 291 of file mapping_p1.cc.
|
virtual |
Determines which additional quantities have to be computed.
flags | Update flags for required quantities. |
Implements Mapping< dim, spacedim >.
Definition at line 72 of file mapping_p1.cc.
|
private |
Auxiliary matrix of gradients of shape functions (used for computation of the Jacobian).
Definition at line 146 of file mapping_p1.hh.