LifeV
CurrentFE Class Reference

CurrentFE - A primordial class for the assembly of the local matrices/vectors retaining the values on the real cells. More...

#include <CurrentFE.hpp>

+ Inheritance diagram for CurrentFE:
+ Collaboration diagram for CurrentFE:

Public Member Functions

template<typename MeshElementMarked >
void update (const MeshElementMarked &geoele, flag_Type upFlag)
 
template<typename MeshElementMarked >
void computeCellNodes (const MeshElementMarked &geoele)
 

Protected Member Functions

 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

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 ReferenceFEM_refFE
 
const GeometricMapM_geoMap
 
QuadratureRuleM_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
 

Constructor & Destructor

 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...
 

Methods

template<typename MeshElementMarkedType >
void update (const MeshElementMarkedType &geoele, flag_Type upFlag)
 Update method using a given element. More...
 
virtual void update (const std::vector< std::vector< Real > > &pts, flag_Type upFlag)
 Update method using only point coordinates. It used the flags, as defined in this page. 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...
 
virtual Real measure () const
 Return the measure of the current element. 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...
 

Set Methods

void setQuadRule (const QuadratureRule &newQuadRule)
 Setter for the quadrature rule. More...
 

Get Methods

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 ReferenceFErefFE () const
 Getter for the reference FE. More...
 
const GeometricMapgeoMap () const
 Getter for the GeometricMap reference. More...
 
const QuadratureRulequadRule () 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...
 

Get Methods for the FE values

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...
 

Old methods (for backward compatibility, avoid using them)

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)
 

Detailed Description

CurrentFE - A primordial class for the assembly of the local matrices/vectors retaining the values on the real cells.

This class is an essential piece of any finite element solver implemented in LifeV, as it is used each time one wants to assemble the matrix or the right hand side associated to a given discrete problem (in variational form). The role of the class is, given a cell (or similar data), to compute the values of the basis function, their derivatives, the location of the quadrature nodes,...

During the assembly procedure, the code of the solver will typically loops on the cells of a given mesh. Each time a cell is selected, a local matrix is built (using only the degrees of freedom of that cell). This local matrix represents a local contribution that is added in the global matrix (concerning all the degrees of freedom of the problem).

Example: suppose that we are assembling the local matrix of a Laplacian problem. The local contributions are given by $ \int_{K} \nabla \phi_i \cdot \nabla \phi_j $ where $K$ is the considered cell. As we use numerical quadrature, to build the local matrix, we have to provide:

  1. The values of the derivatives of the basis functions in the quadrature nodes
  2. The weights associated to the quadrature and adapted to the considered cell

The CurrentFE has then to compute these values. To precise that we want to be able to access these values, we use flags. Here we need the two flags (see this page for more informations on the flags)

  1. For the derivative of the basis functions : UPDATE_DPHI
  2. For the modified weights : UPDATE_WDET

The code that we would typically use is:

my_currentFE.update(my_mesh.volumeList(i), UPDATE_DPHI | UPDATE_WDET);

This line of code asks the CurrentFE to compute the values specified with the flags and corresponding to the cell given by the i-th volume of the desired mesh. Now, these values are accessible using simple calls:

Real my_value = my_currentFE.dphi(0,0,1);
Real my_weight = my_currentFE.wDetJacobian(1);

These two lines just get some values already stored in the CurrentFE (so it is quite cheap!). Usually, one here uses the elementary operation already defined in the file lifefem/elemOper.hpp.

Note: One can change the quadrature of the CurrentFE at any moment, but this is a quite expensive operation and all updates previously done are lost. Use this feature only when it is really needed!

Author
Samuel Quinodoz samue.nosp@m.l.qu.nosp@m.inodo.nosp@m.z@ep.nosp@m.fl.ch
Date
13 Apr 2010
Version
2.0
Todo:

Put all the remaining public member in private

Put ifdef for the checks (definition of the booleans)

Change the old update for wrappers to the new update

CXXFLAG to disable the boost asserts: -DBOOST_DISABLE_ASSERTS?

Definition at line 243 of file CurrentFE.hpp.

Constructor & Destructor Documentation

◆ CurrentFE() [1/4]

CurrentFE ( const ReferenceFE refFE,
const GeometricMap geoMap,
const QuadratureRule qr 
)

Full constructor with default quadrature.

Parameters
refFEReference finite element used
geoMapGeometric mapping used
qrQuadrature rule for the computations

Definition at line 43 of file CurrentFE.cpp.

+ Here is the caller graph for this function:

◆ CurrentFE() [2/4]

CurrentFE ( 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!

Parameters
refFEReference finite element used
geoMapGeometric mapping used

Definition at line 157 of file CurrentFE.cpp.

+ Here is the caller graph for this function:

◆ ~CurrentFE()

virtual ~CurrentFE ( )
inlinevirtual

Destructor.

Definition at line 270 of file CurrentFE.hpp.

◆ CurrentFE() [3/4]

CurrentFE ( )
protected

◆ CurrentFE() [4/4]

CurrentFE ( const CurrentFE fe)
protected

Definition at line 193 of file CurrentFE.cpp.

+ Here is the caller graph for this function:

Member Function Documentation

◆ update() [1/5]

void update ( const MeshElementMarkedType &  geoele,
flag_Type  upFlag 
)

Update method using a given element.

This method is the most important one in this class. It allows the user to update different values on the selected cell, like the values of the derivatives of the basis functions,...

Parameters
geoeleThe cell that we are looking at
upFlagThe flag to explain the quantities that we want to update

◆ update() [2/5]

void update ( const std::vector< std::vector< Real > > &  pts,
flag_Type  upFlag 
)
virtual

Update method using only point coordinates. It used the flags, as defined in this page.

Reimplemented in CurrentFEManifold.

Definition at line 489 of file CurrentFE.cpp.

◆ update() [3/5]

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.

Definition at line 548 of file CurrentFE.cpp.

◆ measure()

Real measure ( ) const
virtual

Return the measure of the current element.

Reimplemented in CurrentFEManifold.

Definition at line 563 of file CurrentFE.cpp.

+ Here is the caller graph for this function:

◆ barycenter()

void barycenter ( Real x,
Real y,
Real z 
) const

return the barycenter of the element

Definition at line 575 of file CurrentFE.cpp.

◆ diameter()

Real diameter ( ) const

return the diameter of the element in the 1-norm

Definition at line 595 of file CurrentFE.cpp.

+ Here is the caller graph for this function:

◆ diameter2()

Real diameter2 ( ) const

return the diameter of the element in the 2-norm

Definition at line 620 of file CurrentFE.cpp.

◆ coorMap() [1/2]

void coorMap ( Real x,
Real y,
Real z,
Real  xi,
Real  eta,
Real  zeta 
) const

compute the coordinate (x,y,z)= F(xi,eta,zeta), (F: geo mappping) where (xi,eta,zeta) are the coor in the ref element and (x,y,z) the coor in the current element (if the code is compiled in 2D mode then z=0 and zeta is disgarded)

Definition at line 648 of file CurrentFE.cpp.

+ Here is the caller graph for this function:

◆ coorMap() [2/2]

GeoVector coorMap ( const GeoVector P) const

return the coordinates in the current element of the point P given in the reference frame.

Definition at line 671 of file CurrentFE.cpp.

+ Here is the caller graph for this function:

◆ quadRuleVTKexport()

void quadRuleVTKexport ( const std::string &  filename) const

Export the quadrature rule on the current FE.

This method can be used to position of the quadrature nodes in the considered element using VTK format. The quadrature point must be updated. It differs from the quadRule::VTKexport in the sens that here, the quadrature is mapped on the current element, while it is still in the reference element for quadRule::VTKexport.

Definition at line 689 of file CurrentFE.cpp.

◆ setQuadRule()

void setQuadRule ( const QuadratureRule newQuadRule)

Setter for the quadrature rule.

This method can be used to change the quadrature rule used. Notice that it is an expensive method. Use it only when it is really needed.

Parameters
newQuadRuleThe new quadrature rule

Definition at line 719 of file CurrentFE.cpp.

◆ currentId()

UInt currentId ( ) const
inline

Getter for the ID of the current cell.

Definition at line 353 of file CurrentFE.hpp.

◆ currentLocalId()

UInt currentLocalId ( ) const
inline

Getter for the local ID of the current cell.

Definition at line 359 of file CurrentFE.hpp.

+ Here is the caller graph for this function:

◆ nbQuadPt()

UInt nbQuadPt ( ) const
inline

Getter for the number of quadrature nodes.

Definition at line 365 of file CurrentFE.hpp.

◆ nbGeoNode()

UInt nbGeoNode ( ) const
inline

Getter for the number of geometrical nodes.

Definition at line 371 of file CurrentFE.hpp.

◆ nbFEDof()

UInt nbFEDof ( ) const
inline

Getter for the number of nodes.

Definition at line 377 of file CurrentFE.hpp.

◆ refFE()

const ReferenceFE& refFE ( ) const
inline

Getter for the reference FE.

Definition at line 383 of file CurrentFE.hpp.

◆ geoMap()

const GeometricMap& geoMap ( ) const
inline

Getter for the GeometricMap reference.

Definition at line 389 of file CurrentFE.hpp.

◆ quadRule()

const QuadratureRule& quadRule ( ) const
inline

Getter for the quadrature rule.

Definition at line 395 of file CurrentFE.hpp.

+ Here is the caller graph for this function:

◆ nbPattern()

UInt nbPattern ( ) const
inline

Getter for the number of entries in the pattern.

Definition at line 401 of file CurrentFE.hpp.

◆ nbDiag()

UInt nbDiag ( ) const
inline

Getter for the diagonal entries in the pattern.

Definition at line 407 of file CurrentFE.hpp.

+ Here is the caller graph for this function:

◆ nbUpper()

UInt nbUpper ( ) const
inline

Getter for the number of entries in the pattern.

Definition at line 413 of file CurrentFE.hpp.

+ Here is the caller graph for this function:

◆ LIFEV_DEPRECATED()

LIFEV_DEPRECATED ( UInt  ) const
inline

Old getter for the number of local coordinates.

Definition at line 419 of file CurrentFE.hpp.

◆ nbLocalCoor()

UInt nbLocalCoor ( ) const
inline

Getter for the number of local coordinates.

Definition at line 425 of file CurrentFE.hpp.

◆ cellNode()

Real cellNode ( UInt  node,
UInt  coordinate 
) const
inline

Getter for the nodes of the cell.

Definition at line 437 of file CurrentFE.hpp.

◆ quadNode()

Real quadNode ( UInt  node,
UInt  coordinate 
) const
inline

Getter for the quadrature nodes.

Definition at line 444 of file CurrentFE.hpp.

+ Here is the caller graph for this function:

◆ detJacobian()

Real detJacobian ( UInt  quadNode) const
inline

Getter for the determinant of the jacobian in a given quadrature node.

Definition at line 451 of file CurrentFE.hpp.

+ Here is the caller graph for this function:

◆ wDetJacobian()

Real wDetJacobian ( UInt  quadNode) const
inline

Getter for the weighted jacobian determinant.

Definition at line 458 of file CurrentFE.hpp.

+ Here is the caller graph for this function:

◆ tInverseJacobian()

Real tInverseJacobian ( UInt  element1,
UInt  element2,
UInt  quadNode 
) const
inline

Getter for the inverse of the jacobian.

Definition at line 465 of file CurrentFE.hpp.

+ Here is the caller graph for this function:

◆ phi() [1/2]

Real phi ( UInt  node,
UInt  quadNode 
) const
inline

Getter for basis function values (scalar FE case)

Definition at line 472 of file CurrentFE.hpp.

+ Here is the caller graph for this function:

◆ phi() [2/2]

Real phi ( UInt  node,
UInt  component,
UInt  quadNode 
) const
inline

Getter for basis function values (vectorial FE case)

Definition at line 479 of file CurrentFE.hpp.

+ Here is the caller graph for this function:

◆ dphi()

Real dphi ( UInt  node,
UInt  derivative,
UInt  quadNode 
) const
inline

Getter for the derivatives of the basis functions.

Definition at line 486 of file CurrentFE.hpp.

+ Here is the caller graph for this function:

◆ d2phi()

Real d2phi ( UInt  node,
UInt  derivative1,
UInt  derivative2,
UInt  quadNode 
) const
inline

Getter for the second derivatives of the basis functions.

Definition at line 493 of file CurrentFE.hpp.

◆ divPhiRef()

Real divPhiRef ( UInt  node,
UInt  quadNode 
) const
inline

Getter for the divergence of a vectorial FE in the reference frame.

Definition at line 500 of file CurrentFE.hpp.

+ Here is the caller graph for this function:

◆ point()

Real point ( UInt  node,
UInt  coordinate 
) const
inline

Old accessor, use cellNode instead.

Definition at line 513 of file CurrentFE.hpp.

◆ quadPt()

Real quadPt ( UInt  node,
UInt  coordinate 
) const
inline

Old accessor, use quadNode instead.

Definition at line 520 of file CurrentFE.hpp.

+ Here is the caller graph for this function:

◆ weightDet()

Real weightDet ( UInt  quadNode) const
inline

Old accessor, use wDetJacobian instead.

Definition at line 527 of file CurrentFE.hpp.

+ Here is the caller graph for this function:

◆ detJac()

Real detJac ( UInt  quadNode) const
inline

Getter for the determinant of the jacobian in a given quadrature node.

Definition at line 534 of file CurrentFE.hpp.

◆ tInvJac()

Real tInvJac ( UInt  element1,
UInt  element2,
UInt  quadNode 
) const
inline

Old accessor, use iInverseJacobian instead.

Definition at line 541 of file CurrentFE.hpp.

◆ phiDer()

Real phiDer ( UInt  node,
UInt  derivative,
UInt  quadNode 
) const
inline

Old accessor, use dphi instead.

Definition at line 548 of file CurrentFE.hpp.

+ Here is the caller graph for this function:

◆ phiDer2()

Real phiDer2 ( UInt  node,
UInt  derivative1,
UInt  derivative2,
UInt  quadNode 
) const
inline

Old accessor, use d2phi instead.

Definition at line 555 of file CurrentFE.hpp.

◆ computeCellNodes() [1/3]

void computeCellNodes ( const MeshElementMarkedType &  geoele)
protected

Update the nodes of the cell to the current one.

◆ computeCellNodes() [2/3]

void computeCellNodes ( const std::vector< std::vector< Real > > &  pts)
protected

Update only the nodes of the cells to the current one.

Definition at line 829 of file CurrentFE.cpp.

◆ computeQuadNodes()

void computeQuadNodes ( )
protected

Update the location of the quadrature in the current cell.

Definition at line 840 of file CurrentFE.cpp.

+ Here is the caller graph for this function:

◆ computeDphiGeometricMap()

void computeDphiGeometricMap ( )
protected

Compute the values of the derivatives of the mapping in the quadrature nodes.

Definition at line 855 of file CurrentFE.cpp.

+ Here is the caller graph for this function:

◆ computeJacobian()

void computeJacobian ( )
protected

Compute the jacobian in the quadrature nodes.

Definition at line 871 of file CurrentFE.cpp.

+ Here is the caller graph for this function:

◆ computeTInverseJacobian()

void computeTInverseJacobian ( )
protected

Compute the transposed inverse of the jacobian in the quadrature nodes.

Definition at line 896 of file CurrentFE.cpp.

+ Here is the caller graph for this function:

◆ computeDetJacobian()

void computeDetJacobian ( )
protected

Compute the determinant of the jacobian.

Definition at line 973 of file CurrentFE.cpp.

+ Here is the caller graph for this function:

◆ computeWDetJacobian()

void computeWDetJacobian ( )
protected

Compute the determinant of the jacobian in the quadrature nodes.

Definition at line 1031 of file CurrentFE.cpp.

+ Here is the caller graph for this function:

◆ computeDphiRef()

void computeDphiRef ( )
protected

Compute the value of the derivatives in the reference frame.

Definition at line 1043 of file CurrentFE.cpp.

+ Here is the caller graph for this function:

◆ computeDphi()

void computeDphi ( )
protected

Compute the value of the derivatives in the current element.

Definition at line 1059 of file CurrentFE.cpp.

+ Here is the caller graph for this function:

◆ computeD2phiRef()

void computeD2phiRef ( )
protected

Compute the value of the second derivatives in the reference frame.

Definition at line 1085 of file CurrentFE.cpp.

+ Here is the caller graph for this function:

◆ computeD2phi()

void computeD2phi ( )
protected

Compute the value of the derivatives in the current element.

Definition at line 1105 of file CurrentFE.cpp.

+ Here is the caller graph for this function:

◆ computePhiVect()

void computePhiVect ( )
protected

Compute the value of the vectorial FE using Piola transform.

Definition at line 1138 of file CurrentFE.cpp.

+ Here is the caller graph for this function:

◆ coorBackMap()

void coorBackMap ( Real  x,
Real  y,
Real  z,
Real xi,
Real eta,
Real zeta 
) const

compute the coordinate (xi,eta,zeta)=inv(F)(x,y,z)

Definition at line 244 of file CurrentFE.cpp.

◆ pointJacobian()

Real pointJacobian ( Real  hat_x,
Real  hat_y,
Real  hat_z,
int  compx,
int  compzeta 
) const

compute the jacobian at a given point : d x_compx / d zeta_compzeta

Definition at line 350 of file CurrentFE.cpp.

+ Here is the caller graph for this function:

◆ pointInverseJacobian()

Real pointInverseJacobian ( Real  hat_x,
Real  hat_y,
Real  hat_z,
int  compx,
int  compzeta 
) const

compute the inverse jacobian

Definition at line 424 of file CurrentFE.cpp.

◆ pointDetJacobian()

Real pointDetJacobian ( Real  hat_x,
Real  hat_y,
Real  hat_z 
) const

compute the determinant of the Jacobian at a given point

Definition at line 368 of file CurrentFE.cpp.

◆ coorQuadPt()

void coorQuadPt ( Real x,
Real y,
Real z,
int  ig 
) const
inline

return (x,y,z) = the global coordinates of the quadrature point ig in the current element.

Warning
this function is almost obsolete since if you call the function updateFirstDerivQuadPt rather than updateFirstDeriv (for example), the coordinates of the quadrature points have already been computed and may be obtained via quadPt(ig,icoor). This is usually much less expensive since it avoids many calls to coorQuadPt

Definition at line 710 of file CurrentFE.hpp.

+ Here is the caller graph for this function:

◆ patternFirst()

int patternFirst ( int  i) const
inline

patternFirst(i): row index in the element matrix of the i-th term of the pattern

Definition at line 716 of file CurrentFE.hpp.

+ Here is the caller graph for this function:

◆ patternSecond()

int patternSecond ( int  i) const
inline

patternSecond(i): column index in the element matrix of the i-th term of the pattern

Definition at line 721 of file CurrentFE.hpp.

+ Here is the caller graph for this function:

◆ update() [4/5]

void update ( const MeshElementMarkedType &  geoele)

DECLARATION of the update methods.

IMPLEMENTATION of the CurrentFE::update methods.

minimal update: we just identify the id of the current element

Definition at line 832 of file CurrentFE.hpp.

◆ updateJac()

void updateJac ( const MeshElementMarkedType &  geoele)

compute the arrays detJac, weightDet, jacobian on the current element

compute the jacobian and its determinant...

Definition at line 843 of file CurrentFE.hpp.

◆ updateJacQuadPt()

void updateJacQuadPt ( const MeshElementMarkedType &  geoele)

compute the arrays detJac, weightDet, jacobian and quadPt on the current element

Definition at line 861 of file CurrentFE.hpp.

◆ updateFirstDeriv()

void updateFirstDeriv ( const MeshElementMarkedType &  geoele)

compute the arrays detJac, weightDet, jacobian, tInvJac, phiDer on the current element

compute the inverse jacobian...

product InvJac by dPhiRef to compute phiDer

Definition at line 881 of file CurrentFE.hpp.

◆ updateFirstDerivQuadPt()

void updateFirstDerivQuadPt ( const MeshElementMarkedType &  geoele)

compute the arrays detJac, weightDet, jacobian, tInvJac, phiDer and quadPt on the current element

compute the inverse jacobian...

product InvJac by dPhiRef to compute phiDer

and the coordinates of the quadrature points

Definition at line 903 of file CurrentFE.hpp.

◆ updateSecondDeriv()

void updateSecondDeriv ( const MeshElementMarkedType &  geoele)

compute the arrays detJac, weightDet, jacobian, tInvJac, phiDer2 on the current element

compute the inverse jacobian...

compute the second derivative

Definition at line 928 of file CurrentFE.hpp.

◆ updateSecondDerivQuadPt()

void updateSecondDerivQuadPt ( const MeshElementMarkedType &  geoele)

compute the arrays detJac, weightDet, jacobian, tInvJac, phiDer2 on the current element

compute the inverse jacobian...

compute the second derivative

and the coordinates of the quadrature points

Definition at line 950 of file CurrentFE.hpp.

◆ updateFirstSecondDeriv()

void updateFirstSecondDeriv ( const MeshElementMarkedType &  geoele)

compute the arrays detJac, weightDet, jacobian, tInvJac, phiDer, phiDer2 on the current element

compute the inverse jacobian...

compute phiDer and phiDer2

Definition at line 973 of file CurrentFE.hpp.

◆ updateFirstSecondDerivQuadPt()

void updateFirstSecondDerivQuadPt ( const MeshElementMarkedType &  geoele)

compute the arrays detJac, weightDet, jacobian, tInvJac, phiDer, phiDer2 on the current element

compute the inverse jacobian...

compute phiDer and phiDer2

and the coordinates of the quadrature points

Definition at line 996 of file CurrentFE.hpp.

◆ update() [5/5]

void update ( const MeshElementMarked geoele,
flag_Type  upFlag 
)

Definition at line 791 of file CurrentFE.hpp.

◆ computeCellNodes() [3/3]

void computeCellNodes ( const MeshElementMarked geoele)

Definition at line 809 of file CurrentFE.hpp.

Field Documentation

◆ M_nbNode

const UInt M_nbNode
protected

Definition at line 611 of file CurrentFE.hpp.

◆ M_nbLocalCoor

const UInt M_nbLocalCoor
protected

Definition at line 612 of file CurrentFE.hpp.

◆ M_nbDiag

const UInt M_nbDiag
protected

Definition at line 613 of file CurrentFE.hpp.

◆ M_nbUpper

const UInt M_nbUpper
protected

Definition at line 614 of file CurrentFE.hpp.

◆ M_nbPattern

const UInt M_nbPattern
protected

Definition at line 615 of file CurrentFE.hpp.

◆ M_currentId

UInt M_currentId
protected

Definition at line 617 of file CurrentFE.hpp.

◆ M_currentLocalId

UInt M_currentLocalId
protected

Definition at line 618 of file CurrentFE.hpp.

◆ M_nbGeoNode

UInt M_nbGeoNode
protected

Definition at line 620 of file CurrentFE.hpp.

◆ M_nbQuadPt

UInt M_nbQuadPt
protected

Definition at line 621 of file CurrentFE.hpp.

◆ M_refFE

const ReferenceFE* M_refFE
protected

Definition at line 625 of file CurrentFE.hpp.

◆ M_geoMap

const GeometricMap* M_geoMap
protected

Definition at line 626 of file CurrentFE.hpp.

◆ M_quadRule

QuadratureRule* M_quadRule
protected

Definition at line 627 of file CurrentFE.hpp.

◆ M_cellNodes

boost::multi_array<Real, 2> M_cellNodes
protected

Definition at line 633 of file CurrentFE.hpp.

◆ M_quadNodes

boost::multi_array<Real, 2> M_quadNodes
protected

Definition at line 634 of file CurrentFE.hpp.

◆ M_dphiGeometricMap

boost::multi_array<Real, 3> M_dphiGeometricMap
protected

Definition at line 636 of file CurrentFE.hpp.

◆ M_jacobian

boost::multi_array<Real, 3> M_jacobian
protected

Definition at line 637 of file CurrentFE.hpp.

◆ M_detJacobian

boost::multi_array<Real, 1> M_detJacobian
protected

Definition at line 638 of file CurrentFE.hpp.

◆ M_wDetJacobian

boost::multi_array<Real, 1> M_wDetJacobian
protected

Definition at line 639 of file CurrentFE.hpp.

◆ M_tInverseJacobian

boost::multi_array<Real, 3> M_tInverseJacobian
protected

Definition at line 640 of file CurrentFE.hpp.

◆ M_phi

boost::multi_array<Real, 3> M_phi
protected

Definition at line 642 of file CurrentFE.hpp.

◆ M_dphi

boost::multi_array<Real, 3> M_dphi
protected

Definition at line 643 of file CurrentFE.hpp.

◆ M_d2phi

boost::multi_array<Real, 4> M_d2phi
protected

Definition at line 644 of file CurrentFE.hpp.

◆ M_phiVect

boost::multi_array<Real, 3> M_phiVect
protected

Definition at line 645 of file CurrentFE.hpp.

◆ M_dphiRef

boost::multi_array<Real, 3> M_dphiRef
protected

Definition at line 648 of file CurrentFE.hpp.

◆ M_d2phiRef

boost::multi_array<Real, 4> M_d2phiRef
protected

Definition at line 649 of file CurrentFE.hpp.

◆ M_divPhiRef

boost::multi_array<Real, 2> M_divPhiRef
protected

Definition at line 650 of file CurrentFE.hpp.

◆ M_cellNodesUpdated

bool M_cellNodesUpdated
protected

Definition at line 654 of file CurrentFE.hpp.

◆ M_quadNodesUpdated

bool M_quadNodesUpdated
protected

Definition at line 655 of file CurrentFE.hpp.

◆ M_dphiGeometricMapUpdated

bool M_dphiGeometricMapUpdated
protected

Definition at line 657 of file CurrentFE.hpp.

◆ M_jacobianUpdated

bool M_jacobianUpdated
protected

Definition at line 658 of file CurrentFE.hpp.

◆ M_detJacobianUpdated

bool M_detJacobianUpdated
protected

Definition at line 659 of file CurrentFE.hpp.

◆ M_wDetJacobianUpdated

bool M_wDetJacobianUpdated
protected

Definition at line 660 of file CurrentFE.hpp.

◆ M_tInverseJacobianUpdated

bool M_tInverseJacobianUpdated
protected

Definition at line 661 of file CurrentFE.hpp.

◆ M_phiUpdated

bool M_phiUpdated
protected

Definition at line 663 of file CurrentFE.hpp.

◆ M_dphiUpdated

bool M_dphiUpdated
protected

Definition at line 664 of file CurrentFE.hpp.

◆ M_d2phiUpdated

bool M_d2phiUpdated
protected

Definition at line 665 of file CurrentFE.hpp.

◆ M_phiVectUpdated

bool M_phiVectUpdated
protected

Definition at line 666 of file CurrentFE.hpp.

◆ M_dphiRefUpdated

bool M_dphiRefUpdated
protected

Definition at line 668 of file CurrentFE.hpp.

◆ M_d2phiRefUpdated

bool M_d2phiRefUpdated
protected

Definition at line 669 of file CurrentFE.hpp.

◆ M_divPhiRefUpdated

bool M_divPhiRefUpdated
protected

Definition at line 670 of file CurrentFE.hpp.


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