Description

ANCF beam element with 3 nodes.

This class implements a continuum-based elastic force formulation.

The node numbering, as follows:

              v
              ^
              |
A o-----+-----o-----+-----o B -> u
             /C
            w

where C is the third and central node.

#include <ChElementBeamANCF.h>

Inheritance diagram for chrono::fea::ChElementBeamANCF:
Collaboration diagram for chrono::fea::ChElementBeamANCF:

Public Types

enum  StrainFormulation { CMPoisson, CMNoPoisson }
 Poisson effect selection. More...
 
using ShapeVector = ChMatrixNM< double, 1, 9 >
 

Public Member Functions

virtual int GetNnodes () override
 Get the number of nodes used by this element.
 
virtual int GetNdofs () override
 Get the number of coordinates in the field used by the referenced nodes.
 
virtual int GetNodeNdofs (int n) override
 Get the number of coordinates from the n-th node used by this element.
 
void SetNodes (std::shared_ptr< ChNodeFEAxyzDD > nodeA, std::shared_ptr< ChNodeFEAxyzDD > nodeB, std::shared_ptr< ChNodeFEAxyzDD > nodeC)
 Specify the nodes of this element.
 
void SetDimensions (double lenX, double beam_h, double beam_w)
 Specify the element dimensions.
 
void SetMaterial (std::shared_ptr< ChMaterialBeamANCF > beam_mat)
 Specify the element material.
 
virtual std::shared_ptr< ChNodeFEAbaseGetNodeN (int n) override
 Access the n-th node of this element.
 
std::shared_ptr< ChNodeFEAxyzDDGetNodeA () const
 Get a handle to the first node of this element.
 
std::shared_ptr< ChNodeFEAxyzDDGetNodeB () const
 Get a handle to the second node of this element.
 
std::shared_ptr< ChNodeFEAxyzDDGetNodeC () const
 Get a handle to the third node of this element.
 
std::shared_ptr< ChMaterialBeamANCFGetMaterial () const
 Return the material.
 
void SetGravityOn (bool val)
 Turn gravity on/off.
 
void SetAlphaDamp (double a)
 Set the structural damping.
 
double GetLengthX () const
 Get the element length in the X direction.
 
double GetThicknessY ()
 Get the total thickness of the beam element.
 
double GetThicknessZ ()
 Get the total thickness of the beam element.
 
void ShapeFunctions (ShapeVector &N, double x, double y, double z)
 Fills the N shape function matrix. More...
 
void ShapeFunctionsDerivativeX (ShapeVector &Nx, double x, double y, double z)
 Fills the Nx shape function derivative matrix with respect to X.
 
void ShapeFunctionsDerivativeY (ShapeVector &Ny, double x, double y, double z)
 Fills the Ny shape function derivative matrix with respect to Y.
 
void ShapeFunctionsDerivativeZ (ShapeVector &Nz, double x, double y, double z)
 Fills the Nz shape function derivative matrix with respect to Z.
 
ChVector EvaluateBeamSectionStrains ()
 Return a vector with three strain components. 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 ComputeKRMmatricesGlobal (ChMatrixRef H, double Kfactor, double Rfactor=0, double Mfactor=0) override
 Sets H as the stiffness matrix K, scaled by Kfactor. More...
 
virtual void ComputeMmatrixGlobal (ChMatrixRef M) override
 Returns the global mass matrix. More...
 
virtual void ComputeNodalMass () override
 Add contribution of element inertia to total nodal masses. More...
 
virtual void ComputeInternalForces (ChVectorDynamic<> &Fi) override
 Computes the internal forces. More...
 
virtual void Update () override
 Update the state of this element.
 
virtual void EvaluateSectionStrain (const double, chrono::ChVector< double > &) override
 
virtual void EvaluateSectionForceTorque (const double, chrono::ChVector< double > &, chrono::ChVector< double > &) override
 
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 abscissa '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...
 
void ComputeInternalJacobians (double Kfactor, double Rfactor)
 Compute Jacobians of the internal forces. More...
 
void ComputeMassMatrix ()
 Compute the mass matrix of the element. More...
 
void ComputeGravityForce (const ChVector<> &g_acc)
 Compute the gravitational forces.
 
double Calc_detJ0 (double x, double y, double z)
 
double Calc_detJ0 (double x, double y, double z, ShapeVector &Nx, ShapeVector &Ny, ShapeVector &Nz, ChMatrixNM< double, 1, 3 > &Nx_d0, ChMatrixNM< double, 1, 3 > &Ny_d0, ChMatrixNM< double, 1, 3 > &Nz_d0)
 
void CalcCoordMatrix (ChMatrixNM< double, 9, 3 > &d)
 
void CalcCoordDerivMatrix (ChVectorN< double, 27 > &dt)
 
void SetStrainFormulation (StrainFormulation model)
 Set the strain formulation.
 
StrainFormulation GetStrainFormulation () const
 Get the strain formulation.
 
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 (velocity 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 (velocity 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
 Tell the number of DOFs blocks (ex. =1 for a body, =4 for a tetrahedron, etc.)
 
virtual unsigned int GetSubBlockOffset (int nblock) override
 Get the offset of the i-th sub-block of DOFs in global vector.
 
virtual unsigned int GetSubBlockSize (int nblock) override
 Get the size of the i-th sub-block of DOFs in global vector.
 
void EvaluateSectionVelNorm (double U, ChVector<> &Result)
 
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 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. More...
 
ChVector ComputeTangent (const double U)
 Gets the tangent to the centerline at the parametric coordinate U. More...
 
- Public Member Functions inherited from chrono::fea::ChElementBeam
virtual void EvaluateSectionForceTorque (const double eta, ChVector<> &Fforce, ChVector<> &Mtorque)=0
 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)=0
 Gets the axial and bending strain of the ANCF "cable" element.
 
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 ComputeGravityForces (ChVectorDynamic<> &Fg, const ChVector<> &G_acc) override
 (This is the default implementation (POTENTIALLY INEFFICIENT) so that in children classes you can avoid implementing this ComputeGravityForces() function, unless you need faster code.) 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.
 
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 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.
 

Friends

class BeamANCF_Mass
 
class BeamANCF_Gravity
 
class BeamANCF_Force
 
class BeamANCF_ForceNu
 
class BeamANCF_Jacobian
 
class BeamANCF_JacobianNu
 

Additional Inherited Members

- Protected Attributes inherited from chrono::fea::ChElementBeam
double mass
 
double length
 
- Protected Attributes inherited from chrono::fea::ChElementGeneric
ChKblockGeneric Kmatr
 

Member Enumeration Documentation

◆ StrainFormulation

Poisson effect selection.

Enumerator
CMPoisson 

Continuum-Mechanics formulation, including Poisson effects.

CMNoPoisson 

Continuum-Mechanics formulation, disregarding Poisson effects.

Member Function Documentation

◆ ComputeInternalForces()

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

Implements chrono::fea::ChElementBase.

◆ ComputeInternalJacobians()

void chrono::fea::ChElementBeamANCF::ComputeInternalJacobians ( double  Kfactor,
double  Rfactor 
)

Compute Jacobians of the internal forces.

This function calculates a linear combination of the stiffness (K) and damping (R) matrices, J = Kfactor * K + Rfactor * R for given coefficients Kfactor and Rfactor. This Jacobian will be further combined with the global mass matrix M and included in the global stiffness matrix H in the function ComputeKRMmatricesGlobal().

◆ ComputeKRMmatricesGlobal()

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

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

Optionally, also superimposes global damping matrix R, scaled by Rfactor, and mass matrix M, scaled by Mfactor. Matrices are expressed in global reference. Corotational elements can take the local Kl & Rl matrices and rotate them.

Implements chrono::fea::ChElementBase.

◆ ComputeMassMatrix()

void chrono::fea::ChElementBeamANCF::ComputeMassMatrix ( )

Compute the mass matrix of the element.

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

◆ ComputeMmatrixGlobal()

void chrono::fea::ChElementBeamANCF::ComputeMmatrixGlobal ( ChMatrixRef  M)
overridevirtual

Returns the global mass matrix.

This is the default implementation, POTENTIALLY VERY INEFFICIENT. Children classes may need to override this with a more efficient version.

Reimplemented from chrono::fea::ChElementGeneric.

◆ ComputeNF() [1/2]

void chrono::fea::ChElementBeamANCF::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 surface
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
state_wif != 0, update state (speed part) to this, then evaluate

Implements chrono::ChLoadableU.

◆ ComputeNF() [2/2]

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

◆ ComputeNodalMass()

void chrono::fea::ChElementBeamANCF::ComputeNodalMass ( )
overridevirtual

Add contribution of element inertia to total nodal masses.

This class computes and adds corresponding masses to ElementGeneric member m_TotalMass.

Reimplemented from chrono::fea::ChElementBase.

◆ ComputeTangent()

ChVector chrono::fea::ChElementBeamANCF::ComputeTangent ( const double  U)

Gets the tangent to the centerline at the parametric coordinate U.

Each coordinate ranging in -1..+1.

◆ EvaluateBeamSectionStrains()

ChVector chrono::fea::ChElementBeamANCF::EvaluateBeamSectionStrains ( )

Return a vector with three strain components.

Beta

◆ EvaluateSectionDisplacement()

void chrono::fea::ChElementBeamANCF::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 abscissa 'eta'.

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

Implements chrono::fea::ChElementBeam.

◆ EvaluateSectionFrame()

void chrono::fea::ChElementBeamANCF::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. Note, 'displ' is the displ.state of 2 nodes, ex. get it as GetStateBlock() Results are corotated (expressed in world reference)

Implements chrono::fea::ChElementBeam.

◆ GetDensity()

double chrono::fea::ChElementBeamANCF::GetDensity ( )
overridevirtual

This is needed so that it can be accessed by ChLoaderVolumeGravity.

Density is mass per unit surface.

Implements chrono::ChLoadableUVW.

◆ GetStateBlock()

void chrono::fea::ChElementBeamANCF::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 size is not this->GetNdofs(), it will be resized. For corotational elements, field is assumed in local reference!

Implements chrono::fea::ChElementBase.

◆ ShapeFunctions()

void chrono::fea::ChElementBeamANCF::ShapeFunctions ( ShapeVector &  N,
double  x,
double  y,
double  z 
)

Fills the N shape function matrix.

NOTE! actually N should be a 3row, 27 column sparse matrix, as N = [s1*eye(3) s2*eye(3) s3*eye(3) s4*eye(3)...]; , but here we only store the s1...s9 values in a vector.


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