Description

A Kirchhoff-Love thin shell element of triangular shape.

This is based on the BST Basic Shell Triangle idea of Onate et al., where 3 purely translational nodes are used for the triangle itself, and 3 other purely translational nodes are used from the neighbouring triangles to infer curvature.

The node numbering is as in the following scheme, where 0-1-2 are the nodes of the triangle, and 3-4-5 are the nodes from neighbouring elements in the mesh:

         2     
  4 o----o---o 3
     \   |\  |
      \  | \ |
       \ |  |
       0 o---o 1
          \  |
           \ |
            |
             o 5

#include <ChElementShellBST.h>

Inheritance diagram for chrono::fea::ChElementShellBST:
Collaboration diagram for chrono::fea::ChElementShellBST:

Classes

class  Layer
 Definition of a layer. More...
 

Public Types

using ShapeVector = ChMatrixNM< double, 1, 3 >
 

Public Member Functions

void SetNodes (std::shared_ptr< ChNodeFEAxyz > node0, std::shared_ptr< ChNodeFEAxyz > node1, std::shared_ptr< ChNodeFEAxyz > node2, std::shared_ptr< ChNodeFEAxyz > node3, std::shared_ptr< ChNodeFEAxyz > node4, std::shared_ptr< ChNodeFEAxyz > node5)
 Specify the nodes of this element. More...
 
virtual int GetNnodes () override
 Get the number of nodes used by this element, not considering those marked with "nullptr" ex. if at boundary.
 
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.
 
virtual std::shared_ptr< ChNodeFEAbaseGetNodeN (int n) override
 Access the n-th node of this element, not considering those marked with "nullptr" ex. if at boundary.
 
std::shared_ptr< ChNodeFEAxyzGetNodeTriangleN (int n) const
 Get a handle to the n node of this element, among the three of the triangle part, n=0..2.
 
std::shared_ptr< ChNodeFEAxyzGetNodeNeighbourN (int n) const
 Get a handle to the n node from the three of the neighbouring triangles, n=0..2. Even if nullptr.
 
void SetAsNeutral ()
 Sets the neutral rotations of nodes at once, assuming the current element position is for zero strain.
 
void AddLayer (double thickness, double theta, std::shared_ptr< ChMaterialShellKirchhoff > material)
 Add a layer. More...
 
void SetLayerZreferenceCentered ()
 Impose the reference z level of shell element as centered along the total thickness. More...
 
void SetLayerZreference (double z_from_bottom)
 Impose the reference z level of shell element respect to the lower face of bottom layer Note! Use after you added all layers.
 
size_t GetNumLayers () const
 Get the number of layers.
 
const LayerGetLayer (size_t i) const
 Get a handle to the specified layer.
 
double GetThickness ()
 Get the total thickness of the shell element (might be sum of multiple layer thicknesses)
 
void ShapeFunctions (ShapeVector &N, const double u, const double v)
 Fills the N shape function matrix, evaluated in "natural" triangle coordinates U and V.
 
void ShapeFunctionsDerivativeU (ShapeVector &Nu, const double u, const double v)
 Fills the Nu shape function derivative matrix with respect to U. Note, may remove u,v parms cause constant.
 
void ShapeFunctionsDerivativeV (ShapeVector &Nv, const double u, const double v)
 Fills the Nv shape function derivative matrix with respect to V. Note, may remove u,v parms cause constant.
 
void ShapeFunctionsDerivativeX (ShapeVector &Nx, const ChMatrixNM< double, 2, 2 > &Jux, const double u, const double v)
 Fills the Nx shape function derivative matrix with respect to X local element direction. More...
 
void ShapeFunctionsDerivativeY (ShapeVector &Ny, const ChMatrixNM< double, 2, 2 > &Jux, const double u, const double v)
 Fills the Ny shape function derivative matrix with respect to Y local element direction. More...
 
virtual void GetStateBlock (ChVectorDynamic<> &mD) override
 Fill the D vector with the current field values at thenodes 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 ComputeInternalForces (ChVectorDynamic<> &Fi) override
 Computes the internal forces. More...
 
void ComputeInternalForces_impl (ChVectorDynamic<> &Fi, ChState &state_x, ChStateDelta &state_w, bool used_for_differentiation=false)
 
virtual void Update () override
 Update the state of this element.
 
virtual void EvaluateSectionDisplacement (const double u, const double v, ChVector<> &u_displ, ChVector<> &u_rotaz) override
 Gets the xyz displacement of a point on the shell, and the rotation RxRyRz of section reference, at parametric coordinates 'u' and 'v'. More...
 
virtual void EvaluateSectionFrame (const double u, const double v, ChVector<> &point, ChQuaternion<> &rot) override
 Gets the absolute xyz position of a point on the shell, and the absolute rotation of section reference, at parametric coordinates 'u' and 'v'. More...
 
virtual void EvaluateSectionPoint (const double u, const double v, ChVector<> &point) override
 Gets the absolute xyz position of a point on the shell, at parametric coordinates 'u' and 'v'. More...
 
virtual void EvaluateSectionVelNorm (double U, double V, ChVector<> &Result) override
 Virtual method to plot velocity field distribution. More...
 
virtual bool IsTriangleShell () override
 Return false if quadrilateral shell - hence u,v parametric coordinates assumed in -1..+1, return true if triangular shell - hence u,v are triangle natural coordinates assumed in 0..+1.
 
void ComputeInternalJacobians (double Kfactor, double Rfactor)
 Compute Jacobians of the internal forces. More...
 
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
 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, 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...
 
virtual ChVector ComputeNormal (const double U, const double V) override
 Gets the normal to the surface at the parametric coordinate U,V. More...
 
virtual bool IsTriangleIntegrationNeeded () override
 If true, use quadrature over u,v in [0..1] range as triangle area coords (with z=1-u-v) otherwise use default quadrature over u,v in [-1..+1] as rectangular isoparametric coords.
 
virtual bool IsTrianglePrismIntegrationNeeded () override
 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 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::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.
 

Public Attributes

std::vector< std::shared_ptr< ChNodeFEAxyz > > m_nodes
 element nodes
 
std::vector< Layerm_layers
 element layers
 
std::vector< double > m_layers_z
 layer separation z values (not scaled, default centered tot thickness)
 
double tot_thickness
 total element thickness
 
ChMatrixNM< double, 2, 2 > Jux
 jacobian [d{u,v}/d{x,y}], as inverse of [d{x,y}/d{u,v}]
 
double area
 
ChVector l0
 
ChVector cM [3]
 
ChVector cI [3]
 
ChVector rI
 
ChVector phi0
 
ChVector phi
 
ChVector k0
 
ChVector e0
 
ChVector k
 
ChVector e
 
ChVector n
 
ChVector m
 

Additional Inherited Members

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

Member Function Documentation

◆ AddLayer()

void chrono::fea::ChElementShellBST::AddLayer ( double  thickness,
double  theta,
std::shared_ptr< ChMaterialShellKirchhoff material 
)

Add a layer.

By default, when adding more than one layer, the reference z level of the shell element is centered along the total thickness.

Parameters
thicknesslayer thickness
thetafiber angle (radians)
materiallayer material

◆ ComputeInternalForces()

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

◆ ComputeInternalForces_impl()

void chrono::fea::ChElementShellBST::ComputeInternalForces_impl ( ChVectorDynamic<> &  Fi,
ChState state_x,
ChStateDelta state_w,
bool  used_for_differentiation = false 
)
Parameters
Fivector of internal forces
state_xstate position to evaluate Fi
state_wstate speed to evaluate Fi
used_for_differentiationtrue if called during finite-difference Jacobian approximation

◆ ComputeInternalJacobians()

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

◆ ComputeNF() [1/2]

void chrono::fea::ChElementShellBST::ComputeNF ( const double  U,
const double  V,
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
Vparametric 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 Q
state_wif != 0, update state (speed part) to this, then evaluate Q

Implements chrono::ChLoadableUV.

◆ ComputeNF() [2/2]

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

◆ ComputeNormal()

ChVector chrono::fea::ChElementShellBST::ComputeNormal ( const double  U,
const double  V 
)
overridevirtual

Gets the normal to the surface at the parametric coordinate U,V.

Each coordinate ranging in -1..+1.

Implements chrono::ChLoadableUV.

◆ EvaluateSectionDisplacement()

void chrono::fea::ChElementShellBST::EvaluateSectionDisplacement ( const double  u,
const double  v,
ChVector<> &  u_displ,
ChVector<> &  u_rotaz 
)
overridevirtual

Gets the xyz displacement of a point on the shell, and the rotation RxRyRz of section reference, at parametric coordinates 'u' and 'v'.

Note, u=-1..+1 , v= -1..+1 parametric coordinates, except if triangular shell, where u=0..+1, v=0..+1, natural triangle coords. Results are not corotated.

Implements chrono::fea::ChElementShell.

◆ EvaluateSectionFrame()

void chrono::fea::ChElementShellBST::EvaluateSectionFrame ( const double  u,
const double  v,
ChVector<> &  point,
ChQuaternion<> &  rot 
)
overridevirtual

Gets the absolute xyz position of a point on the shell, and the absolute rotation of section reference, at parametric coordinates 'u' and 'v'.

Note, u=-1..+1 , v= -1..+1 parametric coordinates, except if triangular shell, where u=0..+1, v=0..+1, natural triangle coords. Results are corotated.

Implements chrono::fea::ChElementShell.

◆ EvaluateSectionPoint()

void chrono::fea::ChElementShellBST::EvaluateSectionPoint ( const double  u,
const double  v,
ChVector<> &  point 
)
overridevirtual

Gets the absolute xyz position of a point on the shell, at parametric coordinates 'u' and 'v'.

Note, u=-1..+1 , v= -1..+1 parametric coordinates, except if triangular shell, where u=0..+1, v=0..+1, natural triangle coords. Results are corotated.

Implements chrono::fea::ChElementShell.

◆ EvaluateSectionVelNorm()

void chrono::fea::ChElementShellBST::EvaluateSectionVelNorm ( double  U,
double  V,
ChVector<> &  Result 
)
overridevirtual

Virtual method to plot velocity field distribution.

Note, u=-1..+1 , v= -1..+1 parametric coordinates, except if triangular shell, where u=0..+1, v=0..+1, natural triangle coords.

Implements chrono::fea::ChElementShell.

◆ GetDensity()

double chrono::fea::ChElementShellBST::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::ChElementShellBST::GetStateBlock ( ChVectorDynamic<> &  mD)
overridevirtual

Fill the D vector with the current field values at thenodes of the element, with proper ordering.

If the D vector has not the size of this->GetNdofs_x(), it will be resized. {x_a y_a z_a x_b y_b z_b ....}

Implements chrono::fea::ChElementBase.

◆ SetLayerZreferenceCentered()

void chrono::fea::ChElementShellBST::SetLayerZreferenceCentered ( )

Impose the reference z level of shell element as centered along the total thickness.

This is the default behavior each time you call AddLayer(); Note! Use after you added all layers.

◆ SetNodes()

void chrono::fea::ChElementShellBST::SetNodes ( std::shared_ptr< ChNodeFEAxyz node0,
std::shared_ptr< ChNodeFEAxyz node1,
std::shared_ptr< ChNodeFEAxyz node2,
std::shared_ptr< ChNodeFEAxyz node3,
std::shared_ptr< ChNodeFEAxyz node4,
std::shared_ptr< ChNodeFEAxyz node5 
)

Specify the nodes of this element.

The node numbering is as in the following scheme, where 0-1-2 are the nodes of the triangle, and 3-4-5 are the nodes from neighbouring elements in the mesh: On the boundary, one or more of nodes 3,4,5 might be "nullptr".

     
         2     
  4 o----o---o 3
     \   |\  |
      \  | \ |
       \ |  |
       0 o---o 1
          \  |
           \ |
            |
             o 5


◆ ShapeFunctionsDerivativeX()

void chrono::fea::ChElementShellBST::ShapeFunctionsDerivativeX ( ShapeVector &  Nx,
const ChMatrixNM< double, 2, 2 > &  Jux,
const double  u,
const double  v 
)

Fills the Nx shape function derivative matrix with respect to X local element direction.

The jacobian Jux = [d{u,v}/d{x,y}] must be provided, often from inverse of Jxu = [d{x,y}/d{u,v}], so that [Nx;Ny]=[Jux]*[Nu;Nv]. Note, may remove u,v parms cause constant

◆ ShapeFunctionsDerivativeY()

void chrono::fea::ChElementShellBST::ShapeFunctionsDerivativeY ( ShapeVector &  Ny,
const ChMatrixNM< double, 2, 2 > &  Jux,
const double  u,
const double  v 
)

Fills the Ny shape function derivative matrix with respect to Y local element direction.

The jacobian Jux = [d{u,v}/d{x,y}] must be provided, often from inverse of Jxu = [d{x,y}/d{u,v}], so that [Nx;Ny]=[Jux]*[Nu;Nv]. Note, may remove u,v parms cause constant


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