Description

ANCF beam element with three nodes.

The coordinates at each node are the position vector and the 2 position vector gradients in the cross-section.

The node numbering, as follows:

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

where C is the third and central node.

#include <ChElementBeamANCF_3333.h>

Inheritance diagram for chrono::fea::ChElementBeamANCF_3333:
Collaboration diagram for chrono::fea::ChElementBeamANCF_3333:

Public Types

enum  IntFrcMethod { IntFrcMethod::ContInt, IntFrcMethod::PreInt }
 Internal force calculation method. More...
 
using VectorN = ChVectorN< double, NSF >
 
using Vector3N = ChVectorN< double, 3 *NSF >
 
using VectorNIP_D0 = ChVectorN< double, NIP_D0 >
 
using VectorNIP_Dv = ChVectorN< double, NIP_Dv >
 
using Matrix3xN = ChMatrixNM< double, 3, NSF >
 
using MatrixNx3 = ChMatrixNM< double, NSF, 3 >
 
using MatrixNx3c = ChMatrixNM_col< double, NSF, 3 >
 
using MatrixNx6 = ChMatrixNM< double, NSF, 6 >
 
using MatrixNxN = ChMatrixNM< double, NSF, NSF >
 

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 active 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 thicknessY, double thicknessZ)
 Specify the element dimensions.
 
void SetMaterial (std::shared_ptr< ChMaterialBeamANCF > beam_mat)
 Specify the element material.
 
std::shared_ptr< ChMaterialBeamANCFGetMaterial () const
 Return the material.
 
virtual std::shared_ptr< ChNodeFEAbaseGetNode (unsigned 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 (middle) node of this element.
 
void SetAlphaDamp (double a)
 Set the structural damping.
 
double GetLengthX () const
 Get the element length in the xi direction (when there is no deformation of the element)
 
double GetThicknessY ()
 Get the total thickness of the beam element in the eta direction (when there is no deformation of the element)
 
double GetThicknessZ ()
 Get the total thickness of the beam element in the zeta direction (when there is no deformation of the element)
 
void SetIntFrcCalcMethod (IntFrcMethod method)
 Set the calculation method to use for the generalized internal force and its Jacobian calculations. More...
 
IntFrcMethod GetIntFrcCalcMethod ()
 Return the type of calculation method currently set for the generalized internal force and its Jacobian calculations.
 
ChMatrix33 GetGreenLagrangeStrain (const double xi, const double eta, const double zeta)
 Get the Green-Lagrange strain tensor at the normalized element coordinates (xi, eta, zeta) at the current state of the element. More...
 
ChMatrix33 GetPK2Stress (const double xi, const double eta, const double zeta)
 Get the 2nd Piola-Kirchoff stress tensor at the normalized element coordinates (xi, eta, zeta) at the current state of the element. More...
 
double GetVonMissesStress (const double xi, const double eta, const double zeta)
 Get the von Mises stress value at the normalized element coordinates (xi, eta, zeta) at the current state of the element. More...
 
virtual void GetStateBlock (ChVectorDynamic<> &mD) override
 Fill the D vector (column matrix) with the current field values at the nodes of the element, with proper ordering. More...
 
virtual void Update () override
 Update the state of this element.
 
virtual void ComputeMmatrixGlobal (ChMatrixRef M) override
 Set M equal to the global mass matrix.
 
virtual void ComputeNodalMass () override
 Add contribution of element inertia to total nodal masses.
 
virtual void ComputeInternalForces (ChVectorDynamic<> &Fi) override
 Compute the generalized internal force vector for the current nodal coordinates as set the value in the Fi vector.
 
virtual void ComputeKRMmatricesGlobal (ChMatrixRef H, double Kfactor, double Rfactor=0, double Mfactor=0) override
 Set H as a linear combination of M, K, and R. 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 EvaluateSectionStrain (const double, chrono::ChVector3d &) override
 Gets the axial and bending strain of the ANCF "cable" element.
 
virtual void EvaluateSectionForceTorque (const double, chrono::ChVector3d &, chrono::ChVector3d &) 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 EvaluateSectionDisplacement (const double xi, 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 '(xi,0,0)'. More...
 
virtual void EvaluateSectionFrame (const double xi, 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 '(xi,0,0)'. More...
 
void EvaluateSectionPoint (const double xi, ChVector3d &point)
 Gets the absolute xyz position of a point on the beam line specified in normalized coordinates xi = -1 at node A and xi = 1 at node B.
 
void EvaluateSectionVel (const double xi, ChVector3d &Result)
 Gets the absolute xyz velocity of a point on the beam line specified in normalized coordinates xi = -1 at node A and xi = 1 at node B.
 
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 (velocity 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 (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 unsigned int GetNumFieldCoords () 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 unsigned int GetNumSubBlocks () override
 Tell the number of DOFs blocks (ex. =1 for a body, =4 for a tetrahedron, etc.)
 
virtual unsigned int GetSubBlockOffset (unsigned int nblock) override
 Get the offset of the i-th sub-block of DOFs in global vector.
 
virtual unsigned int GetSubBlockSize (unsigned int nblock) override
 Get the size of the i-th 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 xi, 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 xi coordinate of the beam 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 xi, const double eta, const double zeta, 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 xi,eta,zeta 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...
 
ChVector3d ComputeTangent (const double xi)
 Gets the tangent to the beam axis at the parametric coordinate xi. More...
 
- 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 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.
 

Static Public Attributes

static const int NP = 3
 number of Gauss quadrature along beam axis
 
static const int NT = 2
 number of quadrature points through cross section
 
static const int NIP_D0 = NP * NT * NT
 number of Gauss quadrature points excluding the Poisson effect for the Enhanced Continuum Mechanics method
 
static const int NIP_Dv
 number of Gauss quadrature points including the Poisson effect for reduced integration along the beam axis for the Enhanced Continuum Mechanics method (Full Gauss quadrature points along the beam axis (xi direction) and 1 point Gauss quadrature in the eta and zeta directions) More...
 
static const int NIP
 total number of integration points for the Enhanced Continuum Mechanics method More...
 
static const int NSF = 9
 number of shape functions
 

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 Enumeration Documentation

◆ IntFrcMethod

Internal force calculation method.

Enumerator
ContInt 

"Continuous Integration" style method - Typically fastest for this element

PreInt 

"Pre-Integration" style method

Member Function Documentation

◆ ComputeKRMmatricesGlobal()

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

Set H as a linear combination of M, K, and R.

H = Mfactor * [M] + Kfactor * [K] + Rfactor * [R], where [M] is the mass matrix, [K] is the stiffness matrix, and [R] is the damping matrix.

Implements chrono::fea::ChElementBase.

◆ ComputeNF() [1/2]

void chrono::fea::ChElementBeamANCF_3333::ComputeNF ( const double  xi,
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 xi coordinate of the beam 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.

For this ANCF element, only the first 6 entries in F are used in the calculation. The first three entries is the applied force in global coordinates and the second 3 entries is the applied moment in global space.

Parameters
xiparametric coordinate along the beam axis
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_3333::ComputeNF ( const double  xi,
const double  eta,
const double  zeta,
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 xi,eta,zeta 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.

For this ANCF element, only the first 6 entries in F are used in the calculation. The first three entries is the applied force in global coordinates and the second 3 entries is the applied moment in global space.

Parameters
xiparametric coordinate in volume
etaparametric coordinate in volume
zetaparametric 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.

◆ ComputeTangent()

ChVector3d chrono::fea::ChElementBeamANCF_3333::ComputeTangent ( const double  xi)

Gets the tangent to the beam axis at the parametric coordinate xi.

xi = -1 at node A and xi = 1 at node B

◆ EvaluateSectionDisplacement()

virtual void chrono::fea::ChElementBeamANCF_3333::EvaluateSectionDisplacement ( const double  xi,
ChVector3d u_displ,
ChVector3d u_rotaz 
)
inlineoverridevirtual

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

xi = -1 at node A and xi = 1 at node B

Implements chrono::fea::ChElementBeam.

◆ EvaluateSectionForceTorque()

virtual void chrono::fea::ChElementBeamANCF_3333::EvaluateSectionForceTorque ( const double  eta,
chrono::ChVector3d Fforce,
chrono::ChVector3d Mtorque 
)
inlineoverridevirtual

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.

◆ EvaluateSectionFrame()

void chrono::fea::ChElementBeamANCF_3333::EvaluateSectionFrame ( const double  xi,
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 '(xi,0,0)'.

xi = -1 at node A and xi = 1 at node B

Implements chrono::fea::ChElementBeam.

◆ GetDensity()

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

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

Density is the average mass per unit volume in the reference state of the element.

Implements chrono::ChLoadableUVW.

◆ GetGreenLagrangeStrain()

ChMatrix33 chrono::fea::ChElementBeamANCF_3333::GetGreenLagrangeStrain ( const double  xi,
const double  eta,
const double  zeta 
)

Get the Green-Lagrange strain tensor at the normalized element coordinates (xi, eta, zeta) at the current state of the element.

Normalized element coordinates span from -1 to 1.

◆ GetPK2Stress()

ChMatrix33 chrono::fea::ChElementBeamANCF_3333::GetPK2Stress ( const double  xi,
const double  eta,
const double  zeta 
)

Get the 2nd Piola-Kirchoff stress tensor at the normalized element coordinates (xi, eta, zeta) at the current state of the element.

Normalized element coordinates span from -1 to 1.

◆ GetStateBlock()

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

Fill the D vector (column matrix) 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. {Pos_a Dv_a Dw_a Pos_b Dv_b Dw_b Pos_c Dv_c Dw_c}

Implements chrono::fea::ChElementBase.

◆ GetVonMissesStress()

double chrono::fea::ChElementBeamANCF_3333::GetVonMissesStress ( const double  xi,
const double  eta,
const double  zeta 
)

Get the von Mises stress value at the normalized element coordinates (xi, eta, zeta) at the current state of the element.

Normalized element coordinates span from -1 to 1.

◆ SetIntFrcCalcMethod()

void chrono::fea::ChElementBeamANCF_3333::SetIntFrcCalcMethod ( IntFrcMethod  method)

Set the calculation method to use for the generalized internal force and its Jacobian calculations.

This should be set prior to the start of the simulation since can be a significant amount of pre-calculation overhead required.

Member Data Documentation

◆ NIP

const int chrono::fea::ChElementBeamANCF_3333::NIP
static
Initial value:

total number of integration points for the Enhanced Continuum Mechanics method

◆ NIP_Dv

const int chrono::fea::ChElementBeamANCF_3333::NIP_Dv
static
Initial value:
=

number of Gauss quadrature points including the Poisson effect for reduced integration along the beam axis for the Enhanced Continuum Mechanics method (Full Gauss quadrature points along the beam axis (xi direction) and 1 point Gauss quadrature in the eta and zeta directions)


The documentation for this class was generated from the following files:
  • /builds/uwsbel/chrono/src/chrono/fea/ChElementBeamANCF_3333.h
  • /builds/uwsbel/chrono/src/chrono/fea/ChElementBeamANCF_3333.cpp
static const int NIP_D0
number of Gauss quadrature points excluding the Poisson effect for the Enhanced Continuum Mechanics m...
Definition: ChElementBeamANCF_3333.h:74
static const int NP
number of Gauss quadrature along beam axis
Definition: ChElementBeamANCF_3333.h:72
static const int NIP_Dv
number of Gauss quadrature points including the Poisson effect for reduced integration along the beam...
Definition: ChElementBeamANCF_3333.h:76