chrono::fea::ChElementBase Class Referenceabstract

Description

Base class for all finite elements, that can be used in the ChMesh physics item.

#include <ChElementBase.h>

Inheritance diagram for chrono::fea::ChElementBase:

Public Member Functions

virtual int GetNnodes ()=0
 Gets the number of nodes used by this element.
 
virtual int GetNdofs ()=0
 Gets the number of coordinates in the field used by the referenced nodes. More...
 
virtual int GetNodeNdofs (int n)=0
 Get the number of coordinates from the n-th node that are used by this element. More...
 
virtual std::shared_ptr< ChNodeFEAbaseGetNodeN (int n)=0
 Access the nth node.
 
virtual void GetStateBlock (ChVectorDynamic<> &mD)=0
 Fills the D vector with the current field values at the nodes of the element, with proper ordering. More...
 
virtual void ComputeMmatrixGlobal (ChMatrixRef M)=0
 Sets M as the mass matrix. More...
 
virtual void ComputeNodalMass ()
 Compute element's nodal masses.
 
virtual void ComputeKRMmatricesGlobal (ChMatrixRef H, double Kfactor, double Rfactor=0, double Mfactor=0)=0
 Sets H as the stiffness matrix K, scaled by Kfactor. More...
 
virtual void ComputeInternalForces (ChVectorDynamic<> &Fi)=0
 Computes the internal forces (ex. More...
 
virtual void Update ()
 Update: this is called at least at each time step. More...
 
virtual void EleDoIntegration ()
 This is optionally implemented if there is some internal state that requires integration.
 
virtual void EleIntLoadResidual_F (ChVectorDynamic<> &R, const double c)
 Adds the internal forces (pasted at global nodes offsets) into a global vector R, multiplied by a scaling factor c, as R += forces * c.
 
virtual void EleIntLoadResidual_Mv (ChVectorDynamic<> &R, const ChVectorDynamic<> &w, const double c)
 Adds the product of element mass M by a vector w (pasted at global nodes offsets) into a global vector R, multiplied by a scaling factor c, as R += M * v * c.
 
virtual void InjectKRMmatrices (ChSystemDescriptor &mdescriptor)=0
 Tell to a system descriptor that there are item(s) of type ChKblock in this object (for further passing it to a solver) Basically does nothing, but inherited classes must specialize this.
 
virtual void KRMmatricesLoad (double Kfactor, double Rfactor, double Mfactor)=0
 Adds the current stiffness K and damping R and mass M matrices in encapsulated ChKblock item(s), if any. More...
 
virtual void VariablesFbLoadInternalForces (double factor=1.)
 Adds the internal forces, expressed as nodal forces, into the encapsulated ChVariables, in the 'fb' part: qf+=forces*factor WILL BE DEPRECATED - see EleIntLoadResidual_F.
 
virtual void VariablesFbIncrementMq ()
 Adds M*q (internal masses multiplied current 'qb') to Fb, ex. More...
 

Friends

class ChMesh
 

Member Function Documentation

◆ ComputeInternalForces()

◆ ComputeKRMmatricesGlobal()

virtual void chrono::fea::ChElementBase::ComputeKRMmatricesGlobal ( ChMatrixRef  H,
double  Kfactor,
double  Rfactor = 0,
double  Mfactor = 0 
)
pure virtual

Sets H as the stiffness matrix K, scaled by Kfactor.

Optionally, also superimposes global damping matrix R, scaled by Rfactor, and mass matrix M, scaled by Mfactor. Matrices are expressed in global reference. Corotational elements can take the local Kl & Rl matrices and rotate them.

Implemented in chrono::fea::ChElementShellReissner4, chrono::fea::ChElementShellBST, chrono::fea::ChElementTetra_4_P, chrono::fea::ChElementShellANCF_8, chrono::fea::ChElementBeamIGA, chrono::fea::ChElementShellANCF, chrono::fea::ChElementBeamANCF, chrono::fea::ChElementBeamEuler, chrono::fea::ChElementHexa_20, chrono::fea::ChElementCableANCF, chrono::fea::ChElementTetra_10, chrono::fea::ChElementHexa_8, chrono::fea::ChElementTetra_4, chrono::fea::ChElementBar, and chrono::fea::ChElementSpring.

◆ ComputeMmatrixGlobal()

virtual void chrono::fea::ChElementBase::ComputeMmatrixGlobal ( ChMatrixRef  M)
pure virtual

◆ GetNdofs()

◆ GetNodeNdofs()

◆ GetStateBlock()

virtual void chrono::fea::ChElementBase::GetStateBlock ( ChVectorDynamic<> &  mD)
pure virtual

◆ KRMmatricesLoad()

virtual void chrono::fea::ChElementBase::KRMmatricesLoad ( double  Kfactor,
double  Rfactor,
double  Mfactor 
)
pure virtual

Adds the current stiffness K and damping R and mass M matrices in encapsulated ChKblock item(s), if any.

The K, R, M matrices are added with scaling values Kfactor, Rfactor, Mfactor.

Implemented in chrono::fea::ChElementGeneric.

◆ Update()

virtual void chrono::fea::ChElementBase::Update ( )
inlinevirtual

Update: this is called at least at each time step.

If the element has to keep updated some auxiliary data, such as the rotation matrices for corotational approach, this is the proper place.

Reimplemented in chrono::fea::ChElementShellReissner4, chrono::fea::ChElementShellBST, chrono::fea::ChElementShellANCF_8, chrono::fea::ChElementShellANCF, chrono::fea::ChElementBeamANCF, chrono::fea::ChElementBeamIGA, chrono::fea::ChElementBeamEuler, chrono::fea::ChElementCableANCF, chrono::fea::ChElementHexahedron, and chrono::fea::ChElementTetrahedron.

◆ VariablesFbIncrementMq()

virtual void chrono::fea::ChElementBase::VariablesFbIncrementMq ( )
inlinevirtual

Adds M*q (internal masses multiplied current 'qb') to Fb, ex.

if qb is initialized with v_old using VariablesQbLoadSpeed, this method can be used in timestepping schemes that do: M*v_new = M*v_old + forces*dt WILL BE DEPRECATED

Reimplemented in chrono::fea::ChElementGeneric.


The documentation for this class was generated from the following file:
  • /builds/uwsbel/chrono/src/chrono/fea/ChElementBase.h