LifeV
|
A class for a finite element on a manifold. More...
#include <CurrentFEManifold.hpp>
Constructor(s) & Destructor | |
CurrentFEManifold (const ReferenceFE &refFE, const GeometricMap &geoMap) | |
Constructor without quadrature specification. More... | |
CurrentFEManifold (const ReferenceFE &refFE, const GeometricMap &geoMap, const QuadratureRule &qr) | |
Full constructor with default quadrature. More... | |
CurrentFEManifold (const ReferenceFE &refFE, const GeometricMap &geoMap, const QuadratureRule &qr, const Real *refCoor, UInt currentId, Real invArea=1.) | |
Constructor used to create a static reference FE for hybrid elements. More... | |
CurrentFEManifold (const CurrentFEManifold &feManifold) | |
Copy constructor. More... | |
virtual | ~CurrentFEManifold () |
Methods | |
virtual void | update (const std::vector< std::vector< Real > > &pts, flag_Type upFlag) |
Update method using only point coordinates. More... | |
void | coorMap (Real &x, Real &y, Real &z, Real xi, Real eta) const |
Maps a point from the reference element to the current element. More... | |
Real | measure () const |
Compute the measure of the current element. More... | |
template<typename FunctorType > | |
Real | integral (const FunctorType &f) const |
Compute the integral of a function f over the current element. More... | |
template<typename FunctorType > | |
Real | normalIntegral (const FunctorType &f) const |
Compute the integral of the normal component of a vector function f over the current element. More... | |
Get Methods | |
boost::multi_array< Real, 2 > | M_phiGeo |
Geometric map in the geometric nodes (different from M_phi only for Hybrid elements) More... | |
boost::multi_array< Real, 3 > | M_tangents |
Normal and tangent vectors, metric tensor and weighted measure in the quadrature nodes. More... | |
boost::multi_array< Real, 2 > | M_normal |
boost::multi_array< Real, 3 > | M_metric |
boost::multi_array< Real, 1 > | M_detMetric |
boost::multi_array< Real, 3 > | M_inverseMetric |
boost::multi_array< Real, 1 > | M_wRootDetMetric |
bool | M_tangentsUpdated |
Check variables. More... | |
bool | M_normalUpdated |
bool | M_metricUpdated |
bool | M_detMetricUpdated |
bool | M_inverseMetricUpdated |
bool | M_wRootDetMetricUpdated |
Real | phiGeo (UInt iGeoNode, UInt quadNode) const |
Values of the geometric basis functions on quadrature points (different from M_phi(iGeoNode,quadNode) only for Hybrid elements) More... | |
Real | tangent (UInt tangent, UInt coordinate, UInt quadNode) const |
Values of the tangents on the quadrature points. More... | |
Real | normal (UInt coordinate, UInt quadNode) const |
Values of the normal on the quadrature points. More... | |
Real | metric (UInt iCoor, UInt jCoor, UInt quadNode) const |
Metric tensor on the quadrature points. More... | |
Real | detMetric (UInt quadNode) |
Determinant of the metric tensor on the quadrature points. More... | |
Real | inverseMetric (UInt iCoor, UInt jCoor, UInt quadNode) const |
Inverse of the metric tensor on the quadrature points. More... | |
Real | wRootDetMetric (UInt quadNode) const |
Square root of the determinant of the metric times the weight on the quadrature points. More... | |
void | computeTangents () |
Computes the tangent vectors in the quadrature nodes. More... | |
void | computeNormal () |
Computes the normal vectors in the quadrature nodes. More... | |
void | computeMetric () |
Computes the metric in the quadrature nodes. More... | |
void | computeDetMetric () |
Computes the determinant of the metric tensor in the quadrature nodes. More... | |
void | computeInverseMetric () |
Computes the inverse of the metric tensor in the quadrature nodes. More... | |
void | computeWRootDetMetric () |
Computes the square root of the determinant of the metric times the weight in the quadrature nodes. More... | |
Additional Inherited Members | |
Public Member Functions inherited from CurrentFE | |
template<typename MeshElementMarked > | |
void | update (const MeshElementMarked &geoele, flag_Type upFlag) |
template<typename MeshElementMarked > | |
void | computeCellNodes (const MeshElementMarked &geoele) |
CurrentFE (const ReferenceFE &refFE, const GeometricMap &geoMap, const QuadratureRule &qr) | |
Full constructor with default quadrature. More... | |
CurrentFE (const ReferenceFE &refFE, const GeometricMap &geoMap) | |
Constructor without quadrature specification. More... | |
virtual | ~CurrentFE () |
Destructor. More... | |
template<typename MeshElementMarkedType > | |
void | update (const MeshElementMarkedType &geoele, flag_Type upFlag) |
Update method using a given element. More... | |
void | update (const std::vector< GeoVector > &pts, flag_Type upFlag) |
Update method using only point coordinates. It used the flags, as defined in this page. More... | |
void | barycenter (Real &x, Real &y, Real &z) const |
return the barycenter of the element More... | |
Real | diameter () const |
return the diameter of the element in the 1-norm More... | |
Real | diameter2 () const |
return the diameter of the element in the 2-norm More... | |
void | coorMap (Real &x, Real &y, Real &z, Real xi, Real eta, Real zeta) const |
GeoVector | coorMap (const GeoVector &P) const |
void | quadRuleVTKexport (const std::string &filename) const |
Export the quadrature rule on the current FE. More... | |
void | setQuadRule (const QuadratureRule &newQuadRule) |
Setter for the quadrature rule. More... | |
UInt | currentId () const |
Getter for the ID of the current cell. More... | |
UInt | currentLocalId () const |
Getter for the local ID of the current cell. More... | |
UInt | nbQuadPt () const |
Getter for the number of quadrature nodes. More... | |
UInt | nbGeoNode () const |
Getter for the number of geometrical nodes. More... | |
UInt | nbFEDof () const |
Getter for the number of nodes. More... | |
const ReferenceFE & | refFE () const |
Getter for the reference FE. More... | |
const GeometricMap & | geoMap () const |
Getter for the GeometricMap reference. More... | |
const QuadratureRule & | quadRule () const |
Getter for the quadrature rule. More... | |
UInt | nbPattern () const |
Getter for the number of entries in the pattern. More... | |
UInt | nbDiag () const |
Getter for the diagonal entries in the pattern. More... | |
UInt | nbUpper () const |
Getter for the number of entries in the pattern. More... | |
LIFEV_DEPRECATED (UInt) nbCoor() const | |
Old getter for the number of local coordinates. More... | |
UInt | nbLocalCoor () const |
Getter for the number of local coordinates. More... | |
Real | cellNode (UInt node, UInt coordinate) const |
Getter for the nodes of the cell. More... | |
Real | quadNode (UInt node, UInt coordinate) const |
Getter for the quadrature nodes. More... | |
Real | detJacobian (UInt quadNode) const |
Getter for the determinant of the jacobian in a given quadrature node. More... | |
Real | wDetJacobian (UInt quadNode) const |
Getter for the weighted jacobian determinant. More... | |
Real | tInverseJacobian (UInt element1, UInt element2, UInt quadNode) const |
Getter for the inverse of the jacobian. More... | |
Real | phi (UInt node, UInt quadNode) const |
Getter for basis function values (scalar FE case) More... | |
Real | phi (UInt node, UInt component, UInt quadNode) const |
Getter for basis function values (vectorial FE case) More... | |
Real | dphi (UInt node, UInt derivative, UInt quadNode) const |
Getter for the derivatives of the basis functions. More... | |
Real | d2phi (UInt node, UInt derivative1, UInt derivative2, UInt quadNode) const |
Getter for the second derivatives of the basis functions. More... | |
Real | divPhiRef (UInt node, UInt quadNode) const |
Getter for the divergence of a vectorial FE in the reference frame. More... | |
Real | point (UInt node, UInt coordinate) const |
Old accessor, use cellNode instead. More... | |
Real | quadPt (UInt node, UInt coordinate) const |
Old accessor, use quadNode instead. More... | |
Real | weightDet (UInt quadNode) const |
Old accessor, use wDetJacobian instead. More... | |
Real | detJac (UInt quadNode) const |
Getter for the determinant of the jacobian in a given quadrature node. More... | |
Real | tInvJac (UInt element1, UInt element2, UInt quadNode) const |
Old accessor, use iInverseJacobian instead. More... | |
Real | phiDer (UInt node, UInt derivative, UInt quadNode) const |
Old accessor, use dphi instead. More... | |
Real | phiDer2 (UInt node, UInt derivative1, UInt derivative2, UInt quadNode) const |
Old accessor, use d2phi instead. More... | |
void | coorBackMap (Real x, Real y, Real z, Real &xi, Real &eta, Real &zeta) const |
Real | pointJacobian (Real hat_x, Real hat_y, Real hat_z, int compx, int compzeta) const |
Real | pointInverseJacobian (Real hat_x, Real hat_y, Real hat_z, int compx, int compzeta) const |
Real | pointDetJacobian (Real hat_x, Real hat_y, Real hat_z) const |
void | coorQuadPt (Real &x, Real &y, Real &z, int ig) const |
int | patternFirst (int i) const |
patternFirst(i): row index in the element matrix of the i-th term of the pattern More... | |
int | patternSecond (int i) const |
patternSecond(i): column index in the element matrix of the i-th term of the pattern More... | |
template<class MeshElementMarkedType > | |
void | update (const MeshElementMarkedType &geoele) |
DECLARATION of the update methods. More... | |
template<class MeshElementMarkedType > | |
void | updateJac (const MeshElementMarkedType &geoele) |
template<class MeshElementMarkedType > | |
void | updateJacQuadPt (const MeshElementMarkedType &geoele) |
template<class MeshElementMarkedType > | |
void | updateFirstDeriv (const MeshElementMarkedType &geoele) |
template<class MeshElementMarkedType > | |
void | updateFirstDerivQuadPt (const MeshElementMarkedType &geoele) |
template<class MeshElementMarkedType > | |
void | updateSecondDeriv (const MeshElementMarkedType &geoele) |
template<class MeshElementMarkedType > | |
void | updateSecondDerivQuadPt (const MeshElementMarkedType &geoele) |
template<class MeshElementMarkedType > | |
void | updateFirstSecondDeriv (const MeshElementMarkedType &geoele) |
template<class MeshElementMarkedType > | |
void | updateFirstSecondDerivQuadPt (const MeshElementMarkedType &geoele) |
Protected Member Functions inherited from CurrentFE | |
CurrentFE () | |
CurrentFE (const CurrentFE &) | |
template<typename MeshElementMarkedType > | |
void | computeCellNodes (const MeshElementMarkedType &geoele) |
Update the nodes of the cell to the current one. More... | |
void | computeCellNodes (const std::vector< std::vector< Real > > &pts) |
Update only the nodes of the cells to the current one. More... | |
void | computeQuadNodes () |
Update the location of the quadrature in the current cell. More... | |
void | computeDphiGeometricMap () |
Compute the values of the derivatives of the mapping in the quadrature nodes. More... | |
void | computeJacobian () |
Compute the jacobian in the quadrature nodes. More... | |
void | computeTInverseJacobian () |
Compute the transposed inverse of the jacobian in the quadrature nodes. More... | |
void | computeDetJacobian () |
Compute the determinant of the jacobian. More... | |
void | computeWDetJacobian () |
Compute the determinant of the jacobian in the quadrature nodes. More... | |
void | computeDphiRef () |
Compute the value of the derivatives in the reference frame. More... | |
void | computeDphi () |
Compute the value of the derivatives in the current element. More... | |
void | computeD2phiRef () |
Compute the value of the second derivatives in the reference frame. More... | |
void | computeD2phi () |
Compute the value of the derivatives in the current element. More... | |
void | computePhiVect () |
Compute the value of the vectorial FE using Piola transform. More... | |
Protected Attributes inherited from CurrentFE | |
const UInt | M_nbNode |
const UInt | M_nbLocalCoor |
const UInt | M_nbDiag |
const UInt | M_nbUpper |
const UInt | M_nbPattern |
UInt | M_currentId |
UInt | M_currentLocalId |
UInt | M_nbGeoNode |
UInt | M_nbQuadPt |
const ReferenceFE * | M_refFE |
const GeometricMap * | M_geoMap |
QuadratureRule * | M_quadRule |
boost::multi_array< Real, 2 > | M_cellNodes |
boost::multi_array< Real, 2 > | M_quadNodes |
boost::multi_array< Real, 3 > | M_dphiGeometricMap |
boost::multi_array< Real, 3 > | M_jacobian |
boost::multi_array< Real, 1 > | M_detJacobian |
boost::multi_array< Real, 1 > | M_wDetJacobian |
boost::multi_array< Real, 3 > | M_tInverseJacobian |
boost::multi_array< Real, 3 > | M_phi |
boost::multi_array< Real, 3 > | M_dphi |
boost::multi_array< Real, 4 > | M_d2phi |
boost::multi_array< Real, 3 > | M_phiVect |
boost::multi_array< Real, 3 > | M_dphiRef |
boost::multi_array< Real, 4 > | M_d2phiRef |
boost::multi_array< Real, 2 > | M_divPhiRef |
bool | M_cellNodesUpdated |
bool | M_quadNodesUpdated |
bool | M_dphiGeometricMapUpdated |
bool | M_jacobianUpdated |
bool | M_detJacobianUpdated |
bool | M_wDetJacobianUpdated |
bool | M_tInverseJacobianUpdated |
bool | M_phiUpdated |
bool | M_dphiUpdated |
bool | M_d2phiUpdated |
bool | M_phiVectUpdated |
bool | M_dphiRefUpdated |
bool | M_d2phiRefUpdated |
bool | M_divPhiRefUpdated |
A class for a finite element on a manifold.
This class inherits from CurrentFE and adds some functionality related to the manifold, like normals, tangents, etc.
Definition at line 86 of file CurrentFEManifold.hpp.
CurrentFEManifold | ( | const ReferenceFE & | refFE, |
const GeometricMap & | geoMap | ||
) |
Constructor without quadrature specification.
If you use this constructor, you have to use the method CurrentFE::setQuadRule before using any update method!
refFE | Reference finite element used |
geoMap | Geometric mapping used |
Definition at line 55 of file CurrentFEManifold.cpp.
CurrentFEManifold | ( | const ReferenceFE & | refFE, |
const GeometricMap & | geoMap, | ||
const QuadratureRule & | qr | ||
) |
Full constructor with default quadrature.
refFE | Reference finite element used |
geoMap | Geometric mapping used |
qr | Quadrature rule for the computations |
Definition at line 36 of file CurrentFEManifold.cpp.
CurrentFEManifold | ( | const ReferenceFE & | refFE, |
const GeometricMap & | geoMap, | ||
const QuadratureRule & | qr, | ||
const Real * | refCoor, | ||
UInt | currentId, | ||
Real | invArea = 1. |
||
) |
Constructor used to create a static reference FE for hybrid elements.
Definition at line 74 of file CurrentFEManifold.cpp.
CurrentFEManifold | ( | const CurrentFEManifold & | feManifold | ) |
Copy constructor.
Definition at line 127 of file CurrentFEManifold.cpp.
|
inlinevirtual |
Definition at line 119 of file CurrentFEManifold.hpp.
Update method using only point coordinates.
Overrides the method of base class CurrentFE. Actually, it calls the method of the base class and then performs some extra updates if upFlag contains also some boundary updates
pts | The coordinates of the points defining the current element |
upFlag | The update flag which determines what kind of update has to be done (see top of the file) |
Reimplemented from CurrentFE.
Definition at line 172 of file CurrentFEManifold.cpp.
Maps a point from the reference element to the current element.
return the coordinate (x,y,z)= F(xi,eta), (F: geo mappping) where (xi,eta) are the coordinates in the reference element and (x,y,z) are the coor in the current element (if the code is compiled in 2D mode then z=0 and eta is disgarded)
Definition at line 155 of file CurrentFEManifold.cpp.
|
virtual |
Compute the measure of the current element.
Overrides the corresponding method in the mother class
Reimplemented from CurrentFE.
Definition at line 160 of file CurrentFEManifold.cpp.
Real integral | ( | const FunctorType & | f | ) | const |
Compute the integral of a function f over the current element.
f | The function to be integrated. It must have the signature Real f(Real x, Real y, Real z) |
Definition at line 266 of file CurrentFEManifold.hpp.
Real normalIntegral | ( | const FunctorType & | f | ) | const |
Compute the integral of the normal component of a vector function f over the current element.
f | The function to be integrated. It must have the signature void f(Real x, Real y, Real z, Real* result) where result is an array of size equal to the dimension of the ambient space. The function should NOT create the array, but instead should assume that the array is already created. The dimension of the ambient spase should be equal to M_nbLocalCoor+1 |
Definition at line 279 of file CurrentFEManifold.hpp.
Values of the geometric basis functions on quadrature points (different from M_phi(iGeoNode,quadNode) only for Hybrid elements)
Definition at line 178 of file CurrentFEManifold.hpp.
Values of the tangents on the quadrature points.
Definition at line 184 of file CurrentFEManifold.hpp.
Values of the normal on the quadrature points.
Definition at line 191 of file CurrentFEManifold.hpp.
Metric tensor on the quadrature points.
Definition at line 198 of file CurrentFEManifold.hpp.
Determinant of the metric tensor on the quadrature points.
Definition at line 205 of file CurrentFEManifold.hpp.
Inverse of the metric tensor on the quadrature points.
Definition at line 212 of file CurrentFEManifold.hpp.
Square root of the determinant of the metric times the weight on the quadrature points.
Definition at line 219 of file CurrentFEManifold.hpp.
|
protected |
Computes the tangent vectors in the quadrature nodes.
Definition at line 217 of file CurrentFEManifold.cpp.
|
protected |
Computes the normal vectors in the quadrature nodes.
Definition at line 236 of file CurrentFEManifold.cpp.
|
protected |
Computes the metric in the quadrature nodes.
Definition at line 281 of file CurrentFEManifold.cpp.
|
protected |
Computes the determinant of the metric tensor in the quadrature nodes.
Definition at line 310 of file CurrentFEManifold.cpp.
|
protected |
Computes the inverse of the metric tensor in the quadrature nodes.
Definition at line 331 of file CurrentFEManifold.cpp.
|
protected |
Computes the square root of the determinant of the metric times the weight in the quadrature nodes.
Definition at line 359 of file CurrentFEManifold.cpp.
|
protected |
Geometric map in the geometric nodes (different from M_phi only for Hybrid elements)
Definition at line 246 of file CurrentFEManifold.hpp.
|
protected |
Normal and tangent vectors, metric tensor and weighted measure in the quadrature nodes.
Definition at line 249 of file CurrentFEManifold.hpp.
|
protected |
Definition at line 250 of file CurrentFEManifold.hpp.
|
protected |
Definition at line 251 of file CurrentFEManifold.hpp.
|
protected |
Definition at line 252 of file CurrentFEManifold.hpp.
|
protected |
Definition at line 253 of file CurrentFEManifold.hpp.
|
protected |
Definition at line 254 of file CurrentFEManifold.hpp.
|
protected |
Check variables.
Definition at line 257 of file CurrentFEManifold.hpp.
|
protected |
Definition at line 258 of file CurrentFEManifold.hpp.
|
protected |
Definition at line 259 of file CurrentFEManifold.hpp.
|
protected |
Definition at line 260 of file CurrentFEManifold.hpp.
|
protected |
Definition at line 261 of file CurrentFEManifold.hpp.
|
protected |
Definition at line 262 of file CurrentFEManifold.hpp.