Description

ANCF brick element with eight nodes.

The coordinates at each node are the position vector and 3 position vector gradients.

The node numbering is in ccw fashion as in the following scheme:

Bottom Layer:
        v
        ^
        |
D o-----+-----o C
  |     |     |
--+-----+-----+----> u
  |     |     |
A o-----+-----o B
Top Layer:
        v
        ^
        |
H o-----+-----o G
  |     |     |
--+-----+-----+----> u
  |     |     |
E o-----+-----o F

#include <ChElementHexaANCF_3843.h>

Inheritance diagram for chrono::fea::ChElementHexaANCF_3843:
Collaboration diagram for chrono::fea::ChElementHexaANCF_3843:

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 = ChVectorN< double, NIP >
 
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 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 GetNdofs_active () override
 Get the number of active 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.
 
virtual int GetNodeNdofs_active (int n) override
 Get the number of active coordinates from the n-th node used by this element.
 
void SetNodes (std::shared_ptr< ChNodeFEAxyzDDD > nodeA, std::shared_ptr< ChNodeFEAxyzDDD > nodeB, std::shared_ptr< ChNodeFEAxyzDDD > nodeC, std::shared_ptr< ChNodeFEAxyzDDD > nodeD, std::shared_ptr< ChNodeFEAxyzDDD > nodeE, std::shared_ptr< ChNodeFEAxyzDDD > nodeF, std::shared_ptr< ChNodeFEAxyzDDD > nodeG, std::shared_ptr< ChNodeFEAxyzDDD > nodeH)
 Specify the nodes of this element.
 
void SetDimensions (double lenX, double lenY, double lenZ)
 Specify the element dimensions.
 
void SetMaterial (std::shared_ptr< ChMaterialHexaANCF > brick_mat)
 Specify the element material.
 
std::shared_ptr< ChMaterialHexaANCFGetMaterial () const
 Return the material.
 
virtual std::shared_ptr< ChNodeFEAbaseGetNodeN (int n) override
 Access the n-th node of this element.
 
virtual std::shared_ptr< ChNodeFEAxyzGetHexahedronNode (int n) override
 Return the specified hexahedron node (0 <= n <= 7).
 
std::shared_ptr< ChNodeFEAxyzDDDGetNodeA () const
 Get a handle to the first node of this element.
 
std::shared_ptr< ChNodeFEAxyzDDDGetNodeB () const
 Get a handle to the second node of this element.
 
std::shared_ptr< ChNodeFEAxyzDDDGetNodeC () const
 Get a handle to the third node of this element.
 
std::shared_ptr< ChNodeFEAxyzDDDGetNodeD () const
 Get a handle to the fourth node of this element.
 
std::shared_ptr< ChNodeFEAxyzDDDGetNodeE () const
 Get a handle to the 5th node of this element.
 
std::shared_ptr< ChNodeFEAxyzDDDGetNodeF () const
 Get a handle to the 6th node of this element.
 
std::shared_ptr< ChNodeFEAxyzDDDGetNodeG () const
 Get a handle to the 7th node of this element.
 
std::shared_ptr< ChNodeFEAxyzDDDGetNodeH () const
 Get a handle to the 8th node of this element.
 
void SetAlphaDamp (double a)
 Set the structural damping.
 
double GetLengthX () const
 Get the element length in the X direction.
 
double GetLengthY ()
 Get the element length in the Y direction.
 
double GetLengthZ ()
 Get the element length in the Z direction.
 
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 ChVector<> &G_acc) override
 Compute the generalized force vector due to gravity using the efficient ANCF specific method.
 
virtual void EvaluateElementDisplacement (const double xi, const double eta, const double zeta, ChVector<> &u_displ, ChVector<> &u_rotaz)
 Gets the xyz displacement of a point in the element, and the approximate rotation RxRyRz at that point '(xi,eta,zeta)'.
 
virtual void EvaluateElementFrame (const double xi, const double eta, const double zeta, ChVector<> &point, ChQuaternion<> &rot)
 Gets the absolute xyz position of a point in the element, and the approximate rotation RxRyRz at that point '(xi,eta,zeta)'. More...
 
virtual void EvaluateElementPoint (const double xi, const double eta, const double zeta, ChVector<> &point)
 Gets the absolute xyz position of a point in the element specified in normalized coordinates.
 
virtual void EvaluateElementVel (const double xi, const double eta, const double zeta, ChVector<> &Result)
 Gets the absolute xyz velocity of a point in the element specified in normalized coordinates.
 
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.
 
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 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...
 
- 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
 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 EleIntLoadResidual_F_gravity (ChVectorDynamic<> &R, const ChVector<> &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
 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
 Add 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
 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 = 4
 number of Gauss quadrature along beam axis
 
static const int NIP = NP * NP * NP
 number of Gauss quadrature points
 
static const int NSF = 32
 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::ChElementGeneric
ChKblockGeneric Kmatr
 

Member Enumeration Documentation

◆ IntFrcMethod

Internal force calculation method.

Enumerator
ContInt 

"Continuous Integration" style method - Typically the Fastest for this element

PreInt 

"Pre-Integration" style method

Member Function Documentation

◆ ComputeKRMmatricesGlobal()

void chrono::fea::ChElementHexaANCF_3843::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()

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

◆ EvaluateElementFrame()

void chrono::fea::ChElementHexaANCF_3843::EvaluateElementFrame ( const double  xi,
const double  eta,
const double  zeta,
ChVector<> &  point,
ChQuaternion<> &  rot 
)
virtual

Gets the absolute xyz position of a point in the element, and the approximate rotation RxRyRz at that point '(xi,eta,zeta)'.

Note, nodeA = (xi=-1, eta=-1, zeta=-1) Note, nodeB = (xi=1, eta=-1, zeta=-1) Note, nodeC = (xi=1, eta=1, zeta=-1) Note, nodeD = (xi=-1, eta=1, zeta=-1) Note, nodeE = (xi=-1, eta=-1, zeta=1) Note, nodeF = (xi=1, eta=-1, zeta=1) Note, nodeG = (xi=1, eta=1, zeta=1) Note, nodeH = (xi=-1, eta=1, zeta=1)

◆ GetDensity()

double chrono::fea::ChElementHexaANCF_3843::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::ChElementHexaANCF_3843::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::ChElementHexaANCF_3843::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::ChElementHexaANCF_3843::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->GetNdofs(), it will be resized. {Pos_a Du_a Dv_a Dw_a Pos_b Du_b Dv_b Dw_b Pos_c Du_c Dv_c Dw_c Pos_d Du_d Dv_d Dw_d}

Implements chrono::fea::ChElementBase.

◆ GetVonMissesStress()

double chrono::fea::ChElementHexaANCF_3843::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::ChElementHexaANCF_3843::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.


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