Description
Base class for most structural elements of 'shell' type.
#include <ChElementShell.h>
Public Member Functions | |
virtual void | EvaluateSectionDisplacement (const double u, const double v, ChVector3d &u_displ, ChVector3d &u_rotaz)=0 |
Gets the xyz displacement of a point on the shell, and the rotation RxRyRz of section reference, at parametric coordinates 'u' and 'v'. More... | |
virtual void | EvaluateSectionFrame (const double u, const double v, ChVector3d &point, ChQuaternion<> &rot)=0 |
Gets the absolute xyz position of a point on the shell, and the absolute rotation of section reference, at parametric coordinates 'u' and 'v'. More... | |
virtual void | EvaluateSectionPoint (const double u, const double v, ChVector3d &point)=0 |
Gets the absolute xyz position of a point on the shell, at parametric coordinates 'u' and 'v'. More... | |
virtual void | EvaluateSectionVelNorm (double U, double V, ChVector3d &Result)=0 |
Virtual method to plot velocity field distribution. More... | |
virtual bool | IsTriangleShell () |
Return false if quadrilateral shell - hence u,v parametric coordinates assumed in -1..+1, return true if triangular shell - hence u,v are triangle natural coordinates assumed in 0..+1. | |
Public Member Functions inherited from chrono::fea::ChElementGeneric | |
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 | |
double | mass |
Protected Attributes inherited from chrono::fea::ChElementGeneric | |
ChKRMBlock | Kmatr |
Member Function Documentation
◆ EvaluateSectionDisplacement()
|
pure virtual |
Gets the xyz displacement of a point on the shell, and the rotation RxRyRz of section reference, at parametric coordinates 'u' and 'v'.
Note, u=-1..+1 , v= -1..+1 parametric coordinates, except if triangular shell, where u=0..+1, v=0..+1, natural triangle coords. Results are not corotated.
Implemented in chrono::fea::ChElementShellANCF_3833, chrono::fea::ChElementShellANCF_3443, chrono::fea::ChElementShellReissner4, chrono::fea::ChElementShellBST, and chrono::fea::ChElementShellANCF_3423.
◆ EvaluateSectionFrame()
|
pure virtual |
Gets the absolute xyz position of a point on the shell, and the absolute rotation of section reference, at parametric coordinates 'u' and 'v'.
Note, u=-1..+1 , v= -1..+1 parametric coordinates, except if triangular shell, where u=0..+1, v=0..+1, natural triangle coords. Results are corotated.
Implemented in chrono::fea::ChElementShellANCF_3833, chrono::fea::ChElementShellANCF_3443, chrono::fea::ChElementShellReissner4, chrono::fea::ChElementShellBST, and chrono::fea::ChElementShellANCF_3423.
◆ EvaluateSectionPoint()
|
pure virtual |
Gets the absolute xyz position of a point on the shell, at parametric coordinates 'u' and 'v'.
Note, u=-1..+1 , v= -1..+1 parametric coordinates, except if triangular shell, where u=0..+1, v=0..+1, natural triangle coords. Results are corotated.
Implemented in chrono::fea::ChElementShellANCF_3833, chrono::fea::ChElementShellANCF_3443, chrono::fea::ChElementShellReissner4, chrono::fea::ChElementShellBST, and chrono::fea::ChElementShellANCF_3423.
◆ EvaluateSectionVelNorm()
|
pure virtual |
Virtual method to plot velocity field distribution.
Note, u=-1..+1 , v= -1..+1 parametric coordinates, except if triangular shell, where u=0..+1, v=0..+1, natural triangle coords.
Implemented in chrono::fea::ChElementShellReissner4, chrono::fea::ChElementShellANCF_3423, chrono::fea::ChElementShellBST, chrono::fea::ChElementShellANCF_3833, and chrono::fea::ChElementShellANCF_3443.
The documentation for this class was generated from the following file:
- /builds/uwsbel/chrono/src/chrono/fea/ChElementShell.h