chrono::fea::ChElementHexaCorot_8 Class Reference

Description

Class for FEA elements of hexahedron type (isoparametric 3D bricks) with 8 nodes.

This element has a linear displacement field.

#include <ChElementHexaCorot_8.h>

Inheritance diagram for chrono::fea::ChElementHexaCorot_8:
Collaboration diagram for chrono::fea::ChElementHexaCorot_8:

Public Types

using ShapeVector = ChMatrixNM< double, 1, 8 >
 

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...
 
double GetVolume ()
 
virtual std::shared_ptr< ChNodeFEAbaseGetNodeN (int n) override
 Access the nth node.
 
virtual std::shared_ptr< ChNodeFEAxyzGetHexahedronNode (int n) override
 Return the specified hexahedron node (0 <= n <= 7).
 
virtual void SetNodes (std::shared_ptr< ChNodeFEAxyz > nodeA, std::shared_ptr< ChNodeFEAxyz > nodeB, std::shared_ptr< ChNodeFEAxyz > nodeC, std::shared_ptr< ChNodeFEAxyz > nodeD, std::shared_ptr< ChNodeFEAxyz > nodeE, std::shared_ptr< ChNodeFEAxyz > nodeF, std::shared_ptr< ChNodeFEAxyz > nodeG, std::shared_ptr< ChNodeFEAxyz > nodeH)
 
virtual void SetDefaultIntegrationRule ()
 
virtual void SetReducedIntegrationRule ()
 
virtual void SetIntegrationRule (int nPoints)
 
void ShapeFunctions (ShapeVector &N, double z0, double z1, double z2)
 Fills the N shape function matrix with the values of shape functions at z0,z1,z2 parametric coordinates, where each zi is in [-1...+1] range. More...
 
virtual void GetStateBlock (ChVectorDynamic<> &mD) override
 Fills the D vector (displacement) with the current field values at the nodes of the element, with proper ordering. More...
 
virtual void ComputeJacobian (ChMatrixDynamic<> &Jacobian, ChMatrixDynamic<> &J1, ChVector<> coord)
 Puts inside 'Jacobian' and 'J1' the Jacobian matrix and the shape functions derivatives matrix of the element. More...
 
virtual void ComputeMatrB (ChMatrixDynamic<> &MatrB, double zeta1, double zeta2, double zeta3, double &JacobianDet)
 Computes the matrix of partial derivatives and puts data in "MatrB" evaluated at natural coordinates zeta1,...,zeta4 . More...
 
virtual void ComputeMatrB (ChGaussPoint *GaussPt, double &JacobianDet)
 Computes the matrix of partial derivatives and puts data in "GaussPt". More...
 
virtual void ComputeStiffnessMatrix ()
 Computes the global STIFFNESS MATRIX of the element: K = sum (w_i * [B]' * [D] * [B]) The number of Gauss Point is defined by SetIntegrationRule function (default: 8 Gp).
 
virtual void Update () override
 Update element at each time step.
 
virtual void UpdateRotation () override
 Given the actual position of the nodes, recompute the cumulative rotation matrix A.
 
ChStrainTensor GetStrain (double z1, double z2, double z3)
 Returns the strain tensor at given parameters. More...
 
ChStressTensor GetStress (double z1, double z2, double z3)
 Returns the stress tensor at given parameters. More...
 
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 (ex. More...
 
void SetMaterial (std::shared_ptr< ChContinuumElastic > my_material)
 Set the material of the element.
 
std::shared_ptr< ChContinuumElasticGetMaterial ()
 
const ChMatrixDynamicGetStiffnessMatrix () const
 Get the StiffnessMatrix.
 
ChGaussPointGetGaussPoint (int N)
 Get the Nth gauss point.
 
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: here the {x,y,z} displacement.
 
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, 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::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 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::fea::ChElementCorotational
ChMatrix33Rotation ()
 Access the cumulative rotation matrix of the element. More...
 
- 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.
 

Additional Inherited Members

- Protected Attributes inherited from chrono::fea::ChElementGeneric
ChKblockGeneric Kmatr
 
- Protected Attributes inherited from chrono::fea::ChElementCorotational
ChMatrix33 A
 

Member Function Documentation

◆ ComputeInternalForces()

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

Computes the internal forces (ex.

the actual position of nodes is not in relaxed reference position) and set values in the Fi vector.

Implements chrono::fea::ChElementBase.

◆ ComputeJacobian()

void chrono::fea::ChElementHexaCorot_8::ComputeJacobian ( ChMatrixDynamic<> &  Jacobian,
ChMatrixDynamic<> &  J1,
ChVector<>  coord 
)
virtual

Puts inside 'Jacobian' and 'J1' the Jacobian matrix and the shape functions derivatives matrix of the element.

The vector "coord" contains the natural coordinates of the integration point. in case of hexahedral elements natural coords vary in the classical range -1 ... +1.

◆ ComputeKRMmatricesGlobal()

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

◆ ComputeMatrB() [1/2]

void chrono::fea::ChElementHexaCorot_8::ComputeMatrB ( ChGaussPoint GaussPt,
double &  JacobianDet 
)
virtual

Computes the matrix of partial derivatives and puts data in "GaussPt".

Stores the determinant of the jacobian in "JacobianDet".

◆ ComputeMatrB() [2/2]

void chrono::fea::ChElementHexaCorot_8::ComputeMatrB ( ChMatrixDynamic<> &  MatrB,
double  zeta1,
double  zeta2,
double  zeta3,
double &  JacobianDet 
)
virtual

Computes the matrix of partial derivatives and puts data in "MatrB" evaluated at natural coordinates zeta1,...,zeta4 .

Also computes determinant of jacobian. note: in case of hexahedral elements natural coord. vary in the range -1 ... +1

◆ ComputeNF()

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

◆ GetNdofs()

virtual int chrono::fea::ChElementHexaCorot_8::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::ChElementHexaCorot_8::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::ChElementHexaCorot_8::GetStateBlock ( ChVectorDynamic<> &  mD)
overridevirtual

Fills the D vector (displacement) 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!

Implements chrono::fea::ChElementBase.

◆ GetStrain()

ChStrainTensor chrono::fea::ChElementHexaCorot_8::GetStrain ( double  z1,
double  z2,
double  z3 
)

Returns the strain tensor at given parameters.

The tensor is in the original undeformed unrotated reference.

◆ GetStress()

ChStressTensor chrono::fea::ChElementHexaCorot_8::GetStress ( double  z1,
double  z2,
double  z3 
)

Returns the stress tensor at given parameters.

The tensor is in the original undeformed unrotated reference.

◆ ShapeFunctions()

void chrono::fea::ChElementHexaCorot_8::ShapeFunctions ( ShapeVector &  N,
double  z0,
double  z1,
double  z2 
)

Fills the N shape function matrix with the values of shape functions at z0,z1,z2 parametric coordinates, where each zi is in [-1...+1] range.

It stores the Ni(z0,z1,z2) values in a 1 row, 8 columns matrix.


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