Description

Classical Timoshenko beam element with two nodes, and tapered sections.

For this tapered beam element, the averaged section properties are used to formulate the mass, stiffness and damping matrices. Note that there are also ChElementCableANCF if no torsional effects are needed, as in cables.

#include <ChElementBeamTaperedTimoshenko.h>

Inheritance diagram for chrono::fea::ChElementBeamTaperedTimoshenko:
Collaboration diagram for chrono::fea::ChElementBeamTaperedTimoshenko:

Public Member Functions

virtual int GetNnodes () override
 Gets the number of nodes used by this element.
 
virtual int GetNdofs () override
 Gets the number of coordinates in the field used by the referenced nodes. More...
 
virtual int GetNodeNdofs (int n) override
 Get the number of coordinates from the n-th node that are used by this element. More...
 
virtual std::shared_ptr< ChNodeFEAbaseGetNodeN (int n) override
 Access the nth node.
 
virtual void SetNodes (std::shared_ptr< ChNodeFEAxyzrot > nodeA, std::shared_ptr< ChNodeFEAxyzrot > nodeB)
 
void SetTaperedSection (std::shared_ptr< ChBeamSectionTaperedTimoshenkoAdvancedGeneric > my_material)
 Set the tapered section & material of beam element .
 
std::shared_ptr< ChBeamSectionTaperedTimoshenkoAdvancedGenericGetTaperedSection ()
 Get the tapered section & material of the element.
 
std::shared_ptr< ChNodeFEAxyzrotGetNodeA ()
 Get the first node (beginning)
 
std::shared_ptr< ChNodeFEAxyzrotGetNodeB ()
 Get the second node (ending)
 
void SetNodeAreferenceRot (ChQuaternion<> mrot)
 Set the reference rotation of nodeA respect to the element rotation.
 
ChQuaternion GetNodeAreferenceRot ()
 Get the reference rotation of nodeA respect to the element rotation.
 
void SetNodeBreferenceRot (ChQuaternion<> mrot)
 Set the reference rotation of nodeB respect to the element rotation.
 
ChQuaternion GetNodeBreferenceRot ()
 Get the reference rotation of nodeB respect to the element rotation.
 
ChQuaternion GetAbsoluteRotation ()
 Get the absolute rotation of element in space This is not the same of Rotation() , that expresses the accumulated rotation from starting point.
 
ChQuaternion GetRefRotation ()
 Get the original reference rotation of element in space.
 
void SetDisableCorotate (bool md)
 Set this as true to have the beam behave like a non-corotated beam hence do not update the corotated reference. More...
 
void SetForceSymmetricStiffness (bool md)
 Set this as true to force the tangent stiffness matrix to be inexact, but symmetric. More...
 
void SetUseGeometricStiffness (bool md)
 Set this as false to disable the contribution of geometric stiffness to the total tangent stiffness. More...
 
void SetUseRc (bool md)
 Set this as true to include the transformation matrix due to the different elastic center offsets at two ends of beam element with respect to the centerline reference, in which case the connection line of two elastic centers is not parallel to the one of two centerline references at two ends. More...
 
void SetUseRs (bool md)
 Set this as true to include the transformation matrix due to the different shear center offsets at two ends of beam element with respect to the centerline reference, in which case the connection line of two shear centers is not parallel to the one of two centerline references at two ends. More...
 
void ShapeFunctionsTimoshenko (ShapeFunctionGroup &NN, double eta)
 Shape functions for Timoshenko beam. More...
 
virtual void Update () override
 Update: this is called at least at each time step. More...
 
virtual void UpdateRotation () override
 Compute large rotation of element for corotational approach The reference frame of this Euler-Bernoulli beam has X aligned to two nodes and Y parallel to Y of 1st node.
 
virtual void GetStateBlock (ChVectorDynamic<> &mD) override
 Fills the D vector with the current field values at the nodes of the element, with proper ordering. More...
 
void GetField_dt (ChVectorDynamic<> &mD_dt)
 Fills the Ddt vector with the current time derivatives of field values at the nodes of the element, with proper ordering. More...
 
void GetField_dtdt (ChVectorDynamic<> &mD_dtdt)
 Fills the Ddtdt vector with the current time derivatives of field values at the nodes of the element, with proper ordering. More...
 
void ComputeStiffnessMatrix ()
 Computes the local (material) stiffness matrix of the element: K = integral( [B]' * [D] * [B] ), Note: in this 'basic' implementation, constant section and constant material are assumed, so the explicit result of quadrature is used. More...
 
void ComputeDampingMatrix ()
 Computes the local (material) damping matrix of the element: R = beta * integral( [B]' * [D] * [B] ), Note: in this 'basic' implementation, constant section and constant material are assumed, so the explicit result of quadrature is used. More...
 
void ComputeGeometricStiffnessMatrix ()
 Computes the local geometric stiffness Kg of the element. More...
 
void ComputeKiRimatricesLocal (bool inertial_damping, bool inertial_stiffness)
 Compute the inertia stiffness matrix [Ki^] and inertial damping matrix [Ri^] which are due to the gyroscopic effect.
 
virtual void ComputeKRMmatricesGlobal (ChMatrixRef H, double Kfactor, double Rfactor=0, double Mfactor=0) override
 Sets H as the global stiffness matrix K, scaled by Kfactor. More...
 
virtual void ComputeInternalForces (ChVectorDynamic<> &Fi) override
 Computes the internal forces (e.g. More...
 
virtual void ComputeInternalForces (ChVectorDynamic<> &Fi, bool Mfactor, bool Kfactor, bool Rfactor, bool Gfactor)
 Computes the inertial forces, damping forces, centrifugal forces and gyroscopic moments, then you could consider them as applied external forces, if you want to do the static solve when including nodal velocites and accelarations. More...
 
virtual void ComputeGravityForces (ChVectorDynamic<> &Fg, const ChVector<> &G_acc) override
 Compute gravity forces, grouped in the Fg vector, one node after the other. More...
 
virtual void EvaluateSectionDisplacement (const double eta, ChVector<> &u_displ, ChVector<> &u_rotaz) override
 Gets the xyz displacement of a point on the beam line, and the rotation RxRyRz of section plane, at abscyssa 'eta'. More...
 
virtual void EvaluateSectionFrame (const double eta, ChVector<> &point, ChQuaternion<> &rot) override
 Gets the absolute xyz position of a point on the beam line, and the absolute rotation of section plane, at abscissa 'eta'. More...
 
virtual void EvaluateSectionForceTorque (const double eta, ChVector<> &Fforce, ChVector<> &Mtorque) override
 Gets the force (traction x, shear y, shear z) and the torque (torsion on x, bending on y, on bending on z) at a section along the beam line, at abscissa 'eta'. More...
 
virtual void EvaluateSectionStrain (const double eta, ChVector<> &StrainV) override
 Gets the axial and bending strain of the ANCF "cable" element.
 
virtual void EvaluateSectionStrain (const double eta, ChVector<> &StrainV_trans, ChVector<> &StrainV_rot)
 Gets the strains(traction along x, shear along y, along shear z, torsion about x, bending about y, on bending about z) at a section along the beam line, at abscissa 'eta'. More...
 
virtual void EvaluateElementStrainEnergy (ChVector<> &StrainEnergyV_trans, ChVector<> &StrainEnergyV_rot)
 Gets the elastic strain energy(traction along x, shear along y, along shear z, torsion about x, bending about y, on bending about z) in the element.
 
virtual void EvaluateElementDampingEnergy (ChVector<> &DampingEnergyV_trans, ChVector<> &DampingEnergyV_rot)
 Gets the damping dissipated energy(traction along x, shear along y, along shear z, torsion about x, bending about y, on bending about z) in the element.
 
virtual int LoadableGet_ndof_x () override
 Gets the number of DOFs affected by this element (position part)
 
virtual int LoadableGet_ndof_w () override
 Gets the number of DOFs affected by this element (speed part)
 
virtual void LoadableGetStateBlock_x (int block_offset, ChState &mD) override
 Gets all the DOFs packed in a single vector (position part)
 
virtual void LoadableGetStateBlock_w (int block_offset, ChStateDelta &mD) override
 Gets all the DOFs packed in a single vector (speed part)
 
virtual void LoadableStateIncrement (const unsigned int off_x, ChState &x_new, const ChState &x, const unsigned int off_v, const ChStateDelta &Dv) override
 Increment all DOFs using a delta.
 
virtual int Get_field_ncoords () override
 Number of coordinates in the interpolated field, ex=3 for a tetrahedron finite element or a cable, = 1 for a thermal problem, etc.
 
virtual int GetSubBlocks () override
 Get the number of DOFs sub-blocks.
 
virtual unsigned int GetSubBlockOffset (int nblock) override
 Get the offset of the specified sub-block of DOFs in global vector.
 
virtual unsigned int GetSubBlockSize (int nblock) override
 Get the size of the specified sub-block of DOFs in global vector.
 
virtual bool IsSubBlockActive (int nblock) const override
 Check if the specified sub-block of DOFs is active.
 
virtual void LoadableGetVariables (std::vector< ChVariables * > &mvars) override
 Get the pointers to the contained ChVariables, appending to the mvars vector.
 
virtual void ComputeNF (const double U, ChVectorDynamic<> &Qi, double &detJ, const ChVectorDynamic<> &F, ChVectorDynamic<> *state_x, ChVectorDynamic<> *state_w) override
 Evaluate N'*F , where N is some type of shape function evaluated at U coordinates of the line, each ranging in -1..+1 F is a load, N'*F is the resulting generalized load Returns also det[J] with J=[dx/du,..], that might be useful in gauss quadrature. More...
 
virtual void ComputeNF (const double U, const double V, const double W, ChVectorDynamic<> &Qi, double &detJ, const ChVectorDynamic<> &F, ChVectorDynamic<> *state_x, ChVectorDynamic<> *state_w) override
 Evaluate N'*F , where N is some type of shape function evaluated at U,V,W coordinates of the volume, each ranging in -1..+1 F is a load, N'*F is the resulting generalized load Returns also det[J] with J=[dx/du,..], that might be useful in gauss quadrature. More...
 
virtual double GetDensity () override
 This is needed so that it can be accessed by ChLoaderVolumeGravity.
 
- Public Member Functions inherited from chrono::fea::ChElementBeam
double GetMass ()
 The full mass of the beam, (with const. section, density, etc.)
 
double GetRestLength ()
 The rest length of the bar.
 
void SetRestLength (double ml)
 Set the rest length of the bar (usually this should be automatically done when SetupInitial is called on beams element, given the current state, but one might need to override this, ex for precompressed beams etc).
 
- Public Member Functions inherited from chrono::fea::ChElementGeneric
ChKblockGenericKstiffness ()
 Access the proxy to stiffness, for sparse solver.
 
virtual void EleIntLoadResidual_F (ChVectorDynamic<> &R, const double c) override
 (This is a default (a bit unoptimal) book keeping so that in children classes you can avoid implementing this EleIntLoadResidual_F function, unless you need faster code)
 
virtual void EleIntLoadResidual_Mv (ChVectorDynamic<> &R, const ChVectorDynamic<> &w, const double c) override
 (This is a default (VERY UNOPTIMAL) book keeping so that in children classes you can avoid implementing this EleIntLoadResidual_Mv function, unless you need faster code.)
 
virtual void EleIntLoadResidual_F_gravity (ChVectorDynamic<> &R, const ChVector<> &G_acc, const double c) override
 (This is a default (VERY UNOPTIMAL) book keeping so that in children classes you can avoid implementing this EleIntLoadResidual_F_gravity function, unless you need faster code. More...
 
virtual void ComputeMmatrixGlobal (ChMatrixRef M) override
 Returns the global mass matrix. More...
 
virtual void InjectKRMmatrices (ChSystemDescriptor &mdescriptor) override
 Tell to a system descriptor that there are item(s) of type ChKblock in this object (for further passing it to a solver)
 
virtual void KRMmatricesLoad (double Kfactor, double Rfactor, double Mfactor) override
 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.) override
 Adds the internal forces, expressed as nodal forces, into the encapsulated ChVariables, in the 'fb' part: qf+=forces*factor (This is a default (a bit unoptimal) book keeping so that in children classes you can avoid implementing this VariablesFbLoadInternalForces function, unless you need faster code)
 
virtual void VariablesFbIncrementMq () override
 Adds M*q (internal masses multiplied current 'qb') to Fb, ex. More...
 
- Public Member Functions inherited from chrono::fea::ChElementBase
virtual void ComputeNodalMass ()
 Compute element's nodal masses.
 
virtual void EleDoIntegration ()
 This is optionally implemented if there is some internal state that requires integration.
 
- Public Member Functions inherited from chrono::ChLoadableUVW
virtual bool IsTetrahedronIntegrationNeeded ()
 If true, use quadrature over u,v,w in [0..1] range as tetrahedron volumetric coords (with z=1-u-v-w) otherwise use default quadrature over u,v,w in [-1..+1] as box isoparametric coords.
 
virtual bool IsTrianglePrismIntegrationNeeded ()
 If true, use quadrature over u,v in [0..1] range as triangle natural coords (with z=1-u-v), and use linear quadrature over w in [-1..+1], otherwise use default quadrature over u,v,w in [-1..+1] as box isoparametric coords.
 
- Public Member Functions inherited from chrono::fea::ChElementCorotational
ChMatrix33Rotation ()
 Access the cumulative rotation matrix of the element. More...
 

Protected Member Functions

virtual void SetupInitial (ChSystem *system) override
 Initial setup. More...
 
void ComputeTransformMatrix ()
 compute the transformation matrix due to offset and rotation of axes.
 
void ComputeTransformMatrixAtPoint (ChMatrixDynamic<> &mT, const double eta)
 compute the transformation matrix due to offset and rotation of axes, at dimensionless abscissa eta.
 

Protected Attributes

std::vector< std::shared_ptr< ChNodeFEAxyzrot > > nodes
 
std::shared_ptr< ChBeamSectionTaperedTimoshenkoAdvancedGenerictapered_section
 Tapered section & material of beam element.
 
ChMatrixDynamic Km
 local material stiffness matrix
 
ChMatrixDynamic Kg
 local geometric stiffness matrix NORMALIZED by P
 
ChMatrixDynamic M
 local material mass matrix. More...
 
ChMatrixDynamic Rm
 local material damping matrix
 
ChMatrixDynamic Ri
 local inertial-damping (gyroscopic damping) matrix
 
ChMatrixDynamic Ki
 local inertial-stiffness matrix
 
ChQuaternion q_refrotA
 
ChQuaternion q_refrotB
 
ChQuaternion q_element_abs_rot
 
ChQuaternion q_element_ref_rot
 
bool disable_corotate
 
bool force_symmetric_stiffness
 
bool use_geometric_stiffness
 whether include geometric stiffness matrix
 
bool use_Rc
 whether use the transformation matrix for elastic axis orientation
 
bool use_Rs
 whether use the transformation matrix for shear axis orientation
 
ChMatrixDynamic T
 transformation matrix for entire beam element
 
ChMatrixDynamic Rc
 transformation matrix for elastic axis orientation
 
ChMatrixDynamic Rs
 transformation matrix for shear axis orientation
 
- Protected Attributes inherited from chrono::fea::ChElementBeam
double mass
 
double length
 
- Protected Attributes inherited from chrono::fea::ChElementGeneric
ChKblockGeneric Kmatr
 
- Protected Attributes inherited from chrono::fea::ChElementCorotational
ChMatrix33 A
 

Member Function Documentation

◆ ComputeDampingMatrix()

void chrono::fea::ChElementBeamTaperedTimoshenko::ComputeDampingMatrix ( )

Computes the local (material) damping matrix of the element: R = beta * integral( [B]' * [D] * [B] ), Note: in this 'basic' implementation, constant section and constant material are assumed, so the explicit result of quadrature is used.

Also, this local material damping matrix is constant, computed only at the beginning for performance reasons; if you later change some material property, call this or InitialSetup(). Only the stiffness term(beta) is used for this implemented Rayleigh damping model.

◆ ComputeGeometricStiffnessMatrix()

void chrono::fea::ChElementBeamTaperedTimoshenko::ComputeGeometricStiffnessMatrix ( )

Computes the local geometric stiffness Kg of the element.

Note: this->Kg will be set as the geometric stiffness EXCLUDING the multiplication by the P pull force, in fact P multiplication happens in all terms, thus this allows making the Kg as a constant matrix that is computed only at the beginning, and later it is multiplied by P all times the real Kg is needed. If you later change some material property, call this or InitialSetup().

◆ ComputeGravityForces()

void chrono::fea::ChElementBeamTaperedTimoshenko::ComputeGravityForces ( ChVectorDynamic<> &  Fg,
const ChVector<> &  G_acc 
)
overridevirtual

Compute gravity forces, grouped in the Fg vector, one node after the other.

TODO for the lumped mass matrix case, the mM * mG product can be unrolled into few multiplications as mM mostly zero, and same for mG

Reimplemented from chrono::fea::ChElementGeneric.

◆ ComputeInternalForces() [1/2]

void chrono::fea::ChElementBeamTaperedTimoshenko::ComputeInternalForces ( ChVectorDynamic<> &  Fi)
overridevirtual

Computes the internal forces (e.g.

the actual position of nodes is not in relaxed reference position) and set values in the Fi vector. This functionality can be used to output the forces and torques at two nodes directly.

Implements chrono::fea::ChElementBase.

◆ ComputeInternalForces() [2/2]

void chrono::fea::ChElementBeamTaperedTimoshenko::ComputeInternalForces ( ChVectorDynamic<> &  Fi,
bool  Mfactor,
bool  Kfactor,
bool  Rfactor,
bool  Gfactor 
)
virtual

Computes the inertial forces, damping forces, centrifugal forces and gyroscopic moments, then you could consider them as applied external forces, if you want to do the static solve when including nodal velocites and accelarations.

Strictly speaking, it is not static solve any more. We can name it as quasi-static equilibrium solving, just the same as Simpack. It is recommended to use ChStaticNonLinearRheonomicAnalysis to do this kind of quasi-static equilibrium solving in case of rotating beams.

◆ ComputeKRMmatricesGlobal()

void chrono::fea::ChElementBeamTaperedTimoshenko::ComputeKRMmatricesGlobal ( ChMatrixRef  H,
double  Kfactor,
double  Rfactor = 0,
double  Mfactor = 0 
)
overridevirtual

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

Optionally, also superimposes global damping matrix R, scaled by Rfactor, and global mass matrix M multiplied by Mfactor.

< current angular velocity of section of node A, in material frame

< current angular acceleration of section of node A, in material frame

< current acceleration of section of node A, in material frame)

< current angular velocity of section of node B, in material frame

< current angular acceleration of section of node B, in material frame

< current acceleration of section of node B, in material frame

Implements chrono::fea::ChElementBase.

◆ ComputeNF() [1/2]

void chrono::fea::ChElementBeamTaperedTimoshenko::ComputeNF ( const double  U,
ChVectorDynamic<> &  Qi,
double &  detJ,
const ChVectorDynamic<> &  F,
ChVectorDynamic<> *  state_x,
ChVectorDynamic<> *  state_w 
)
overridevirtual

Evaluate N'*F , where N is some type of shape function evaluated at U coordinates of the line, each ranging in -1..+1 F is a load, N'*F is the resulting generalized load Returns also det[J] with J=[dx/du,..], that might be useful in gauss quadrature.

Parameters
Uparametric coordinate in line
QiReturn result of Q = N'*F here
detJReturn det[J] here
FInput F vector, size is =n. field coords.
state_xif != 0, update state (pos. part) to this, then evaluate Q
state_wif != 0, update state (speed part) to this, then evaluate Q

Implements chrono::ChLoadableU.

Reimplemented in chrono::fea::ChElementBeamTaperedTimoshenkoFPM.

◆ ComputeNF() [2/2]

void chrono::fea::ChElementBeamTaperedTimoshenko::ComputeNF ( const double  U,
const double  V,
const double  W,
ChVectorDynamic<> &  Qi,
double &  detJ,
const ChVectorDynamic<> &  F,
ChVectorDynamic<> *  state_x,
ChVectorDynamic<> *  state_w 
)
overridevirtual

Evaluate N'*F , where N is some type of shape function evaluated at U,V,W coordinates of the volume, each ranging in -1..+1 F is a load, N'*F is the resulting generalized load Returns also det[J] with J=[dx/du,..], that might be useful in gauss quadrature.

Parameters
Uparametric coordinate in volume
Vparametric coordinate in volume
Wparametric coordinate in volume
QiReturn result of N'*F here, maybe with offset block_offset
detJReturn det[J] here
FInput F vector, size is = n.field coords.
state_xif != 0, update state (pos. part) to this, then evaluate Q
state_wif != 0, update state (speed part) to this, then evaluate Q

Implements chrono::ChLoadableUVW.

Reimplemented in chrono::fea::ChElementBeamTaperedTimoshenkoFPM.

◆ ComputeStiffnessMatrix()

void chrono::fea::ChElementBeamTaperedTimoshenko::ComputeStiffnessMatrix ( )

Computes the local (material) stiffness matrix of the element: K = integral( [B]' * [D] * [B] ), Note: in this 'basic' implementation, constant section and constant material are assumed, so the explicit result of quadrature is used.

Also, this local material stiffness matrix is constant, computed only at the beginning for performance reasons; if you later change some material property, call this or InitialSetup().

◆ EvaluateSectionDisplacement()

void chrono::fea::ChElementBeamTaperedTimoshenko::EvaluateSectionDisplacement ( const double  eta,
ChVector<> &  u_displ,
ChVector<> &  u_rotaz 
)
overridevirtual

Gets the xyz displacement of a point on the beam line, and the rotation RxRyRz of section plane, at abscyssa 'eta'.

Note, eta=-1 at node1, eta=+1 at node2. Results are not corotated.

Implements chrono::fea::ChElementBeam.

Reimplemented in chrono::fea::ChElementBeamTaperedTimoshenkoFPM.

◆ EvaluateSectionForceTorque()

void chrono::fea::ChElementBeamTaperedTimoshenko::EvaluateSectionForceTorque ( const double  eta,
ChVector<> &  Fforce,
ChVector<> &  Mtorque 
)
overridevirtual

Gets the force (traction x, shear y, shear z) and the torque (torsion on x, bending on y, on bending on z) at a section along the beam line, at abscissa 'eta'.

Note, eta=-1 at node1, eta=+1 at node2. Results are not corotated, and are expressed in the reference system of beam.

Implements chrono::fea::ChElementBeam.

Reimplemented in chrono::fea::ChElementBeamTaperedTimoshenkoFPM.

◆ EvaluateSectionFrame()

void chrono::fea::ChElementBeamTaperedTimoshenko::EvaluateSectionFrame ( const double  eta,
ChVector<> &  point,
ChQuaternion<> &  rot 
)
overridevirtual

Gets the absolute xyz position of a point on the beam line, and the absolute rotation of section plane, at abscissa 'eta'.

Note, eta=-1 at node1, eta=+1 at node2. Results are corotated (expressed in world reference)

Implements chrono::fea::ChElementBeam.

◆ EvaluateSectionStrain()

void chrono::fea::ChElementBeamTaperedTimoshenko::EvaluateSectionStrain ( const double  eta,
ChVector<> &  StrainV_trans,
ChVector<> &  StrainV_rot 
)
virtual

Gets the strains(traction along x, shear along y, along shear z, torsion about x, bending about y, on bending about z) at a section along the beam line, at abscissa 'eta'.

It's evaluated at the elastic center. Note, eta=-1 at node1, eta=+1 at node2. Results are not corotated, and are expressed in the reference system of beam.

Reimplemented in chrono::fea::ChElementBeamTaperedTimoshenkoFPM.

◆ GetField_dt()

void chrono::fea::ChElementBeamTaperedTimoshenko::GetField_dt ( ChVectorDynamic<> &  mD_dt)

Fills the Ddt vector with the current time derivatives of field values at the nodes of the element, with proper ordering.

If the D vector has not the size of this->GetNdofs(), it will be resized. For corotational elements, field is assumed in local reference! Give that this element includes rotations at nodes, this gives: {v_a v_a v_a wx_a wy_a wz_a v_b v_b v_b wx_b wy_b wz_b}

◆ GetField_dtdt()

void chrono::fea::ChElementBeamTaperedTimoshenko::GetField_dtdt ( ChVectorDynamic<> &  mD_dtdt)

Fills the Ddtdt vector with the current time derivatives of field values at the nodes of the element, with proper ordering.

If the D vector has not the size of this->GetNdofs(), it will be resized. For corotational elements, field is assumed in local reference! Give that this element includes rotations at nodes, this gives: {acc_a acc_a acc_a accx_a accy_a accz_a acc_b acc_b acc_b accx_b accy_b accz_b}

◆ GetNdofs()

virtual int chrono::fea::ChElementBeamTaperedTimoshenko::GetNdofs ( )
inlineoverridevirtual

Gets the number of coordinates in the field used by the referenced nodes.

This is for example the size (n.of rows/columns) of the local stiffness matrix.

Implements chrono::fea::ChElementBase.

◆ GetNodeNdofs()

virtual int chrono::fea::ChElementBeamTaperedTimoshenko::GetNodeNdofs ( int  n)
inlineoverridevirtual

Get the number of coordinates from the n-th node that are used by this element.

Note that this may be different from the value returned by GetNodeN(n)->Get_ndof_w();

Implements chrono::fea::ChElementBase.

◆ GetStateBlock()

void chrono::fea::ChElementBeamTaperedTimoshenko::GetStateBlock ( ChVectorDynamic<> &  mD)
overridevirtual

Fills the D vector with the current field values at the nodes of the element, with proper ordering.

If the D vector has not the size of this->GetNdofs(), it will be resized. For corotational elements, field is assumed in local reference! Given that this element includes rotations at nodes, this gives: {x_a y_a z_a Rx_a Ry_a Rz_a x_b y_b z_b Rx_b Ry_b Rz_b}

Implements chrono::fea::ChElementBase.

◆ SetDisableCorotate()

void chrono::fea::ChElementBeamTaperedTimoshenko::SetDisableCorotate ( bool  md)
inline

Set this as true to have the beam behave like a non-corotated beam hence do not update the corotated reference.

Just for benchmarks!

◆ SetForceSymmetricStiffness()

void chrono::fea::ChElementBeamTaperedTimoshenko::SetForceSymmetricStiffness ( bool  md)
inline

Set this as true to force the tangent stiffness matrix to be inexact, but symmetric.

This allows the use of faster solvers. For systems close to the equilibrium, the tangent stiffness would be symmetric anyway.

◆ SetupInitial()

void chrono::fea::ChElementBeamTaperedTimoshenko::SetupInitial ( ChSystem system)
overrideprotectedvirtual

Initial setup.

Precompute mass and matrices that do not change during the simulation, such as the local tangent stiffness Kl of each element, if needed, etc.

Reimplemented from chrono::fea::ChElementBase.

◆ SetUseGeometricStiffness()

void chrono::fea::ChElementBeamTaperedTimoshenko::SetUseGeometricStiffness ( bool  md)
inline

Set this as false to disable the contribution of geometric stiffness to the total tangent stiffness.

By default it is on.

◆ SetUseRc()

void chrono::fea::ChElementBeamTaperedTimoshenko::SetUseRc ( bool  md)
inline

Set this as true to include the transformation matrix due to the different elastic center offsets at two ends of beam element with respect to the centerline reference, in which case the connection line of two elastic centers is not parallel to the one of two centerline references at two ends.

By default it is on. Please refer to ANSYS theory document for more information.

◆ SetUseRs()

void chrono::fea::ChElementBeamTaperedTimoshenko::SetUseRs ( bool  md)
inline

Set this as true to include the transformation matrix due to the different shear center offsets at two ends of beam element with respect to the centerline reference, in which case the connection line of two shear centers is not parallel to the one of two centerline references at two ends.

By default it is on. This transformation matrix could decrease the equivalent torsional stiffness of beam element, and hence generate larger torsional rotation. Please refer to ANSYS theory document for more information.

◆ ShapeFunctionsTimoshenko()

void chrono::fea::ChElementBeamTaperedTimoshenko::ShapeFunctionsTimoshenko ( ShapeFunctionGroup &  NN,
double  eta 
)

Shape functions for Timoshenko beam.

Please refer to the textbook: J. S. Przemieniecki, Theory of Matrix Structural Analysis-Dover Publications (1985).

◆ Update()

void chrono::fea::ChElementBeamTaperedTimoshenko::Update ( )
overridevirtual

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 from chrono::fea::ChElementBase.

Member Data Documentation

◆ M

ChMatrixDynamic chrono::fea::ChElementBeamTaperedTimoshenko::M
protected

local material mass matrix.

It could be lumped or consistent mass matrix, depending on SetLumpedMassMatrix(true/false)


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