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 ComputeGravityForces (ChVectorDynamic<> &Fi, const ChVector<> &G_acc)=0
 Computes the gravitational forces and set values in the Fi vector, with n.rows = n.of dof of element.
 
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 * w * c.
 
virtual void EleIntLoadResidual_F_gravity (ChVectorDynamic<> &R, const ChVector<> &G_acc, const double c)=0
 Adds the contribution of gravity loads, multiplied by a scaling factor c, as: R += M * g * c Note that it is up to the element implementation to build a proper g vector that contains G_acc values in the proper stride (ex. More...
 
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

◆ ComputeMmatrixGlobal()

◆ EleIntLoadResidual_F_gravity()

virtual void chrono::fea::ChElementBase::EleIntLoadResidual_F_gravity ( ChVectorDynamic<> &  R,
const ChVector<> &  G_acc,
const double  c 
)
pure virtual

Adds the contribution of gravity loads, multiplied by a scaling factor c, as: R += M * g * c Note that it is up to the element implementation to build a proper g vector that contains G_acc values in the proper stride (ex.

tetahedrons have 4x copies of G_acc in g). Note that elements can provide fast implementations that do not need to build any internal M matrix, and not even the g vector, for instance if using lumped masses.

Implemented in chrono::fea::ChElementGeneric.

◆ GetNdofs()

◆ GetNodeNdofs()

◆ GetStateBlock()

◆ 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()

◆ 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