Description

Simple beam element with two nodes and ANCF gradient-deficient formulation.

For this 'basic' implementation, constant section and constant material are assumed along the beam coordinate. Torsional stiffness is impossible because of the formulation.
Based on the formulation in
"Analysis of Thin Beams and Cables Using the Absolute Nodal Co-ordinate Formulation", J. Gerstmayr, A. Shabana, Nonlinear Dynamics (2006) 45: 109-130, DOI: 10.1007/s11071-006-1856-1

#include <ChElementCableANCF.h>

Inheritance diagram for chrono::fea::ChElementCableANCF:
Collaboration diagram for chrono::fea::ChElementCableANCF:

Public Types

using ShapeVector = ChMatrixNM< double, 1, 4 >
 

Public Member Functions

virtual unsigned int GetNumNodes () override
 Get the number of nodes used by this element.
 
virtual unsigned int GetNumCoordsPosLevel () override
 Get the number of coordinates in the field used by the referenced nodes.
 
virtual unsigned int GetNumCoordsPosLevelActive () override
 Get the number of active coordinates in the field used by the referenced nodes.
 
virtual unsigned int GetNodeNumCoordsPosLevel (unsigned int n) override
 Get the number of coordinates from the n-th node used by this element.
 
virtual unsigned int GetNodeNumCoordsPosLevelActive (unsigned int n) override
 Get the number of coordinates from the n-th node used by this element.
 
virtual std::shared_ptr< ChNodeFEAbaseGetNode (unsigned int n) override
 Access the nth node.
 
virtual void SetNodes (std::shared_ptr< ChNodeFEAxyzD > nodeA, std::shared_ptr< ChNodeFEAxyzD > nodeB)
 
void SetSection (std::shared_ptr< ChBeamSectionCable > section)
 Set the section & material of beam element. More...
 
std::shared_ptr< ChBeamSectionCableGetSection ()
 Get the section & material of the element.
 
std::shared_ptr< ChNodeFEAxyzDGetNodeA ()
 Get the first node (beginning).
 
std::shared_ptr< ChNodeFEAxyzDGetNodeB ()
 Get the second node (ending).
 
double GetCurrLength ()
 Get element length.
 
virtual void ShapeFunctions (ShapeVector &N, double xi)
 Fills the N shape function matrix with the values of shape functions at abscissa 'xi'. More...
 
virtual void ShapeFunctionsDerivatives (ShapeVector &Nd, double xi)
 Fills the N shape function derivative matrix with the values of shape function derivatives at abscissa 'xi'. More...
 
virtual void ShapeFunctionsDerivatives2 (ShapeVector &Ndd, double xi)
 Fills the N shape function derivative matrix with the values of shape function 2nd derivatives at abscissa 'xi'. More...
 
virtual void Update () override
 Update, called at least at each time step. More...
 
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...
 
virtual void ComputeInternalJacobians (double Kfactor, double Rfactor)
 Computes the stiffness matrix of the element: K = integral( .... More...
 
virtual void ComputeMassMatrix ()
 Computes the mass matrix of the element. More...
 
virtual void ComputeMmatrixGlobal (ChMatrixRef M) override
 Sets M as the global mass matrix.
 
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 and set values in the Fi vector. More...
 
virtual void ComputeGravityForces (ChVectorDynamic<> &Fg, const ChVector3d &G_acc) override
 Compute the generalized force vector due to gravity using the efficient ANCF specific method.
 
virtual void EvaluateSectionDisplacement (const double eta, ChVector3d &u_displ, ChVector3d &u_rotaz) override
 Gets the xyz displacement of a point on the beam line and the rotation RxRyRz of section plane, at abscissa 'eta'. More...
 
virtual void EvaluateSectionFrame (const double eta, ChVector3d &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, ChVector3d &Fforce, ChVector3d &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, ChVector3d &StrainV) override
 Gets the axial and bending strain of the ANCF element torque (torsion on x, bending on y, on bending on z) at a section along the beam line, at abscissa 'eta'. More...
 
void SetAlphaDamp (double a)
 Set structural damping.
 
virtual unsigned int GetLoadableNumCoordsPosLevel () override
 Gets the number of DOFs affected by this element (position part).
 
virtual unsigned int GetLoadableNumCoordsVelLevel () override
 Gets the number of DOFs affected by this element (speed part).
 
virtual void LoadableGetStateBlockPosLevel (int block_offset, ChState &mD) override
 Gets all the DOFs packed in a single vector (position part).
 
virtual void LoadableGetStateBlockVelLevel (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 unsigned int GetNumFieldCoords () override
 Number of coordinates in the interpolated field.
 
virtual unsigned int GetNumSubBlocks () override
 Get the number of DOFs sub-blocks.
 
virtual unsigned int GetSubBlockOffset (unsigned int nblock) override
 Get the offset of the specified sub-block of DOFs in global vector.
 
virtual unsigned int GetSubBlockSize (unsigned int nblock) override
 Get the size of the specified sub-block of DOFs in global vector.
 
virtual bool IsSubBlockActive (unsigned 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,V coordinates of the surface, each ranging in -1..+1. 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. More...
 
virtual double GetDensity () override
 Return the material density for this element.
 
- 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
ChKRMBlockKstiffness ()
 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 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 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 Attributes

bool m_use_damping
 Boolean indicating whether internal damping is added.
 
double m_alpha
 Scaling factor for internal damping.
 

Additional Inherited Members

- Protected Attributes inherited from chrono::fea::ChElementANCF
int m_element_dof
 actual number of degrees of freedom for the element
 
bool m_full_dof
 true if all node variables are active (not fixed)
 
ChArray< int > m_mapping_dof
 indices of active DOFs (set only is some are fixed)
 
- Protected Attributes inherited from chrono::fea::ChElementBeam
double mass
 
double length
 
- Protected Attributes inherited from chrono::fea::ChElementGeneric
ChKRMBlock Kmatr
 

Member Function Documentation

◆ ComputeInternalForces()

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

Computes the internal forces and set values in the Fi vector.

(e.g. the actual position of nodes is not in relaxed reference position).

Implements chrono::fea::ChElementBase.

◆ ComputeInternalJacobians()

void chrono::fea::ChElementCableANCF::ComputeInternalJacobians ( double  Kfactor,
double  Rfactor 
)
virtual

Computes the stiffness matrix of the element: K = integral( ....

), Note: in this 'basic' implementation, constant section and constant material are assumed.

◆ ComputeKRMmatricesGlobal()

void chrono::fea::ChElementCableANCF::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.

Implements chrono::fea::ChElementBase.

◆ ComputeMassMatrix()

void chrono::fea::ChElementCableANCF::ComputeMassMatrix ( )
virtual

Computes the mass matrix of the element.

Note: in this 'basic' implementation, constant section and constant material are assumed.

◆ ComputeNF() [1/2]

void chrono::fea::ChElementCableANCF::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,V coordinates of the surface, 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.

◆ ComputeNF() [2/2]

void chrono::fea::ChElementCableANCF::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.

◆ EvaluateSectionDisplacement()

void chrono::fea::ChElementCableANCF::EvaluateSectionDisplacement ( const double  eta,
ChVector3d u_displ,
ChVector3d u_rotaz 
)
overridevirtual

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

Note, eta=-1 at node1, eta=+1 at node2. Note that 'displ' is the displ.state of 2 nodes, e.g. get it as GetStateBlock() Results are not corotated.

Implements chrono::fea::ChElementBeam.

◆ EvaluateSectionForceTorque()

void chrono::fea::ChElementCableANCF::EvaluateSectionForceTorque ( const double  eta,
ChVector3d Fforce,
ChVector3d 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. Note that 'displ' is the displ.state of 2 nodes, ex. get it as GetStateBlock(). Results are not corotated, and are expressed in the reference system of beam. This is not mandatory for the element to work, but it can be useful for plotting, showing results, etc.

Implements chrono::fea::ChElementBeam.

◆ EvaluateSectionFrame()

void chrono::fea::ChElementCableANCF::EvaluateSectionFrame ( const double  eta,
ChVector3d 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. Note that 'displ' is the displ.state of 2 nodes, e.g. get it as GetStateBlock() Results are corotated (expressed in world reference)

Implements chrono::fea::ChElementBeam.

◆ EvaluateSectionStrain()

void chrono::fea::ChElementCableANCF::EvaluateSectionStrain ( const double  eta,
ChVector3d StrainV 
)
overridevirtual

Gets the axial and bending strain of the ANCF element 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. Note that 'displ' is the displ.state of 2 nodes, ex. get it as GetStateBlock(). Results are not corotated, and are expressed in the reference system of beam. This is not mandatory for the element to work, but it can be useful for plotting, showing results, etc.

Implements chrono::fea::ChElementBeam.

◆ GetStateBlock()

void chrono::fea::ChElementCableANCF::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->GetNumCoordsPosLevel(), it will be resized. {x_a y_a z_a Dx_a Dx_a Dx_a x_b y_b z_b Dx_b Dy_b Dz_b}

Implements chrono::fea::ChElementBase.

◆ SetSection()

void chrono::fea::ChElementCableANCF::SetSection ( std::shared_ptr< ChBeamSectionCable section)
inline

Set the section & material of beam element.

It is a shared property, so it can be shared between other beams.

◆ ShapeFunctions()

void chrono::fea::ChElementCableANCF::ShapeFunctions ( ShapeVector &  N,
double  xi 
)
virtual

Fills the N shape function matrix with the values of shape functions at abscissa 'xi'.

Note, xi=0 at node1, xi=+1 at node2. N should be a 3x12 parse matrix, N = [s1*eye(3) s2*eye(3) s3*eye(3) s4*eye(3)], but is stored here in a compressed form: only the s1 s2 s3 s4 values in a 4x1 column vector.

◆ ShapeFunctionsDerivatives()

void chrono::fea::ChElementCableANCF::ShapeFunctionsDerivatives ( ShapeVector &  Nd,
double  xi 
)
virtual

Fills the N shape function derivative matrix with the values of shape function derivatives at abscissa 'xi'.

Note, xi=0 at node1, xi=+1 at node2. In a compressed form, only four values are stored in a 4x1 column vector.

◆ ShapeFunctionsDerivatives2()

void chrono::fea::ChElementCableANCF::ShapeFunctionsDerivatives2 ( ShapeVector &  Ndd,
double  xi 
)
virtual

Fills the N shape function derivative matrix with the values of shape function 2nd derivatives at abscissa 'xi'.

Note, xi=0 at node1, xi=+1 at node2. In a compressed form, only four values are stored in a 4x1 column vector.

◆ Update()

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

Update, 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 should be implemented in this function.

Reimplemented from chrono::fea::ChElementBase.


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