Description
Class for all elements whose stiffness matrix can be seen as an NxN block-matrix split among N nodes.
Most FEA elements inherited from ChElementGeneric need to implement at most the two fundamental methods ComputeKRMmatricesGlobal(), ComputeInternalForces(), and optionally ComputeGravityForces().
#include <ChElementGeneric.h>
Public Member Functions | |
ChKRMBlock & | Kstiffness () |
Access the proxy to stiffness, for sparse solver. | |
virtual void | EleIntLoadResidual_F (ChVectorDynamic<> &R, const double c) override |
Add the internal forces (pasted at global nodes offsets) into a global vector R, multiplied by a scaling factor c, as R += forces * c This default implementation is SLIGHTLY INEFFICIENT. | |
virtual void | EleIntLoadResidual_Mv (ChVectorDynamic<> &R, const ChVectorDynamic<> &w, const double c) override |
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 This default implementation is VERY INEFFICIENT. | |
virtual void | EleIntLoadLumpedMass_Md (ChVectorDynamic<> &Md, double &error, const double c) override |
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) override |
Add the contribution of gravity loads, multiplied by a scaling factor c, as: R += M * g * c This default implementation is VERY INEFFICIENT. More... | |
virtual void | ComputeGravityForces (ChVectorDynamic<> &Fg, const ChVector3d &G_acc) override |
Compute the gravitational forces. More... | |
virtual void | ComputeMmatrixGlobal (ChMatrixRef M) override |
Calculate the mass matrix, expressed in global reference. More... | |
virtual void | InjectKRMMatrices (ChSystemDescriptor &descriptor) override |
Register with the given system descriptor any ChKRMBlock objects associated with this item. | |
virtual void | LoadKRMMatrices (double Kfactor, double Rfactor, double Mfactor) override |
Compute and load current stiffnes (K), damping (R), and mass (M) matrices in encapsulated ChKRMBlock objects. More... | |
virtual void | VariablesFbLoadInternalForces (double factor=1.) override |
Add the internal forces, expressed as nodal forces, into the encapsulated ChVariables. | |
virtual void | VariablesFbIncrementMq () override |
Add M*q (internal masses multiplied current 'qb'). | |
Public Member Functions inherited from chrono::fea::ChElementBase | |
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< ChNodeFEAbase > | GetNode (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 | 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 | 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. | |
Protected Attributes | |
ChKRMBlock | Kmatr |
Member Function Documentation
◆ ComputeGravityForces()
|
overridevirtual |
Compute the gravitational forces.
This default implementation (POTENTIALLY INEFFICIENT) uses a temporary ChLoaderGravity that applies the load to elements only if they are inherited by ChLoadableUVW so it can use GetDensity() and Gauss quadrature.
Implements chrono::fea::ChElementBase.
Reimplemented in chrono::fea::ChElementShellANCF_3833, chrono::fea::ChElementShellANCF_3443, chrono::fea::ChElementHexaANCF_3843, chrono::fea::ChElementBeamTaperedTimoshenko, chrono::fea::ChElementShellANCF_3423, chrono::fea::ChElementBeamANCF_3333, chrono::fea::ChElementBeamANCF_3243, chrono::fea::ChElementBeamEuler, chrono::fea::ChElementBeamIGA, and chrono::fea::ChElementCableANCF.
◆ ComputeMmatrixGlobal()
|
overridevirtual |
Calculate the mass matrix, expressed in global reference.
This default implementation (POTENTIALLY VERY INEFFICIENT) should be overriden by derived classes a more efficient version.
Implements chrono::fea::ChElementBase.
Reimplemented in chrono::fea::ChElementShellReissner4, chrono::fea::ChElementShellANCF_3833, chrono::fea::ChElementShellANCF_3443, chrono::fea::ChElementHexaANCF_3843, chrono::fea::ChElementShellANCF_3423, chrono::fea::ChElementBeamANCF_3333, chrono::fea::ChElementBeamANCF_3243, and chrono::fea::ChElementCableANCF.
◆ EleIntLoadLumpedMass_Md()
|
overridevirtual |
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 This default implementation is VERY INEFFICIENT.
Reimplemented from chrono::fea::ChElementBase.
Reimplemented in chrono::fea::ChElementShellBST.
◆ EleIntLoadResidual_F_gravity()
|
overridevirtual |
Add the contribution of gravity loads, multiplied by a scaling factor c, as: R += M * g * c This default implementation is VERY INEFFICIENT.
This fallback implementation uses a temp ChLoaderGravity that applies the load to elements only if they are inherited by ChLoadableUVW so it can use GetDensity() and Gauss quadrature.
Implements chrono::fea::ChElementBase.
◆ LoadKRMMatrices()
|
overridevirtual |
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.
Implements chrono::fea::ChElementBase.
The documentation for this class was generated from the following files:
- /builds/uwsbel/chrono/src/chrono/fea/ChElementGeneric.h
- /builds/uwsbel/chrono/src/chrono/fea/ChElementGeneric.cpp