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 unsigned int GetNumNodes ()=0
 Get the number of nodes used by this element.
 
virtual unsigned int GetNumCoordsPosLevel ()=0
 Get the number of coordinates in the field used by the referenced nodes. More...
 
virtual unsigned int GetNumCoordsPosLevelActive ()
 Get the actual number of active degrees of freedom. More...
 
virtual unsigned int GetNodeNumCoordsPosLevel (unsigned int n)=0
 Get the number of coordinates from the specified node that are used by this element. More...
 
virtual unsigned int GetNodeNumCoordsPosLevelActive (unsigned int n)
 Get the actual number of active coordinates from the specified node that are used by this element. More...
 
virtual std::shared_ptr< ChNodeFEAbaseGetNode (unsigned int n)=0
 Access the nth node.
 
virtual void GetStateBlock (ChVectorDynamic<> &mD)=0
 Fill the D vector with the current field values at the nodes of the element, with proper ordering. More...
 
virtual void ComputeMmatrixGlobal (ChMatrixRef M)=0
 Calculate the mass matrix, expressed in global reference.
 
virtual void ComputeNodalMass ()
 Compute element's nodal masses.
 
virtual void ComputeKRMmatricesGlobal (ChMatrixRef H, double Kfactor, double Rfactor=0, double Mfactor=0)=0
 Set H as the stiffness matrix K, scaled by Kfactor. More...
 
virtual void ComputeInternalForces (ChVectorDynamic<> &Fi)=0
 Compute the internal forces. More...
 
virtual void ComputeGravityForces (ChVectorDynamic<> &Fi, const ChVector3d &G_acc)=0
 Compute the gravitational forces. More...
 
virtual void Update ()
 Update, 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)
 Add 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)
 Add 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 EleIntLoadLumpedMass_Md (ChVectorDynamic<> &Md, double &error, const double c)
 Adds the lumped mass to a Md vector, representing a mass diagonal matrix. More...
 
virtual void EleIntLoadResidual_F_gravity (ChVectorDynamic<> &R, const ChVector3d &G_acc, const double c)=0
 Add 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 &descriptor)=0
 Register with the given system descriptor any ChKRMBlock objects associated with this item.
 
virtual void LoadKRMMatrices (double Kfactor, double Rfactor, double Mfactor)=0
 Compute and load current stiffnes (K), damping (R), and mass (M) matrices in encapsulated ChKRMBlock objects. More...
 
virtual void VariablesFbLoadInternalForces (double factor=1.0)
 Add the internal forces, expressed as nodal forces, into the encapsulated ChVariables. More...
 
virtual void VariablesFbIncrementMq ()
 Add M*q (internal masses multiplied current 'qb'). More...
 

Friends

class ChMesh
 

Member Function Documentation

◆ ComputeGravityForces()

virtual void chrono::fea::ChElementBase::ComputeGravityForces ( ChVectorDynamic<> &  Fi,
const ChVector3d G_acc 
)
pure virtual

◆ ComputeInternalForces()

◆ ComputeKRMmatricesGlobal()

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

◆ EleIntLoadLumpedMass_Md()

virtual void chrono::fea::ChElementBase::EleIntLoadLumpedMass_Md ( ChVectorDynamic<> &  Md,
double &  error,
const double  c 
)
inlinevirtual

Adds the lumped mass to a Md vector, representing a mass diagonal matrix.

Used by lumped explicit integrators. If mass lumping is impossible or approximate, adds scalar error to "error" parameter. Md += c*diag(M) or Md += c*HRZ(M) or other lumping heuristics

Reimplemented in chrono::fea::ChElementShellBST, and chrono::fea::ChElementGeneric.

◆ EleIntLoadResidual_F_gravity()

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

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

◆ GetNodeNumCoordsPosLevel()

◆ GetNodeNumCoordsPosLevelActive()

virtual unsigned int chrono::fea::ChElementBase::GetNodeNumCoordsPosLevelActive ( unsigned int  n)
inlinevirtual

Get the actual number of active coordinates from the specified node that are used by this element.

The default implementation returns the full number of coordinates for this element, but some elements may have nodes with fixed variables.

Reimplemented in chrono::fea::ChElementShellANCF_3443, chrono::fea::ChElementShellANCF_3833, chrono::fea::ChElementBeamANCF_3243, chrono::fea::ChElementHexaANCF_3843, chrono::fea::ChElementBeamANCF_3333, chrono::fea::ChElementShellANCF_3423, chrono::fea::ChElementCableANCF, chrono::fea::ChElementHexaANCF_3813_9, and chrono::fea::ChElementHexaANCF_3813.

◆ GetNumCoordsPosLevel()

◆ GetNumCoordsPosLevelActive()

virtual unsigned int chrono::fea::ChElementBase::GetNumCoordsPosLevelActive ( )
inlinevirtual

Get the actual number of active degrees of freedom.

The default implementation returns the full number of coordinates for this element, but some elements may have nodes with fixed variables.

Reimplemented in chrono::fea::ChElementShellANCF_3443, chrono::fea::ChElementShellANCF_3833, chrono::fea::ChElementBeamANCF_3243, chrono::fea::ChElementHexaANCF_3843, chrono::fea::ChElementBeamANCF_3333, chrono::fea::ChElementShellANCF_3423, chrono::fea::ChElementCableANCF, chrono::fea::ChElementHexaANCF_3813_9, and chrono::fea::ChElementHexaANCF_3813.

◆ GetStateBlock()

◆ LoadKRMMatrices()

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

Compute and load current stiffnes (K), damping (R), and mass (M) matrices in encapsulated ChKRMBlock objects.

The resulting KRM blocks represent linear combinations of the K, R, and M matrices, with the specified coefficients Kfactor, Rfactor,and Mfactor, respectively.

Implemented in chrono::fea::ChElementGeneric.

◆ Update()

◆ VariablesFbIncrementMq()

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

Add M*q (internal masses multiplied current 'qb').

Update fb. For example, 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.

◆ VariablesFbLoadInternalForces()

virtual void chrono::fea::ChElementBase::VariablesFbLoadInternalForces ( double  factor = 1.0)
inlinevirtual

Add the internal forces, expressed as nodal forces, into the encapsulated ChVariables.

Update the 'fb' part: qf+=forces*factor WILL BE DEPRECATED - see EleIntLoadResidual_F

Reimplemented in chrono::fea::ChElementGeneric.


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