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 unsigned int GetNumNodes () override
 Get the number of nodes used by this element, not considering those marked with "nullptr" ex. if at boundary.
 
virtual unsigned int GetNumCoordsPosLevel () override
 Get the number of 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 std::shared_ptr< ChNodeFEAbaseGetNode (unsigned int n) override
 Access the n-th node of this element, not considering those marked with "nullptr" ex. if at boundary n=0..5.
 
std::shared_ptr< ChNodeFEAxyzGetNodeMainTriangle (unsigned 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< ChNodeFEAxyzGetNodeNeighbour (unsigned 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
 Set 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.
 
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 EvaluateSectionDisplacement (const double u, const double v, ChVector3d &u_displ, ChVector3d &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, ChVector3d &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, ChVector3d &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, ChVector3d &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 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
 Get the number of DOFs sub-blocks.
 
virtual unsigned int GetSubBlockOffset (unsigned int nblock) override
 Get the offset of the specified sub-block of DOFs in global vector.
 
virtual unsigned int GetSubBlockSize (unsigned int nblock) override
 Get the size of the specified 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 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 ChVector3d 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
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 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 ComputeGravityForces (ChVectorDynamic<> &Fg, const ChVector3d &G_acc) override
 Compute the gravitational forces. More...
 
virtual void ComputeMmatrixGlobal (ChMatrixRef M) override
 Calculate the mass matrix, expressed in global reference. 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 unsigned int GetNumCoordsPosLevelActive ()
 Get the actual number of active degrees of freedom. More...
 
virtual unsigned int GetNodeNumCoordsPosLevelActive (unsigned int n)
 Get the actual number of active coordinates from the specified node that are used by this element. More...
 
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
 initial element triangle area
 
ChVector3d l0
 initial lengths
 
ChVector3d cM [3]
 cM coefficients
 
ChVector3d cI [3]
 cI coefficients
 
ChVector3d rI
 rI coefficients = RIi/(RMi/RIi) = (1/hIi)/((1/hMi)(1/hIi))
 
ChVector3d phi0
 initial edge bendings
 
ChVector3d phi
 actual edge bendings (last computed)
 
ChVector3d k0
 initial curvature (not needed?)
 
ChVector3d e0
 initial strain
 
ChVector3d k
 actual curvature (last computed)
 
ChVector3d e
 actual strain (last computed)
 
ChVector3d n
 actual stress, membrane (last computed)
 
ChVector3d m
 actual stress, bending (last computed)
 

Additional Inherited Members

- Protected Attributes inherited from chrono::fea::ChElementShell
double mass
 
- Protected Attributes inherited from chrono::fea::ChElementGeneric
ChKRMBlock 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

Set 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()

ChVector3d 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.

◆ EleIntLoadLumpedMass_Md()

void chrono::fea::ChElementShellBST::EleIntLoadLumpedMass_Md ( ChVectorDynamic<> &  Md,
double &  error,
const double  c 
)
overridevirtual

Adds the lumped mass to a Md vector, representing a mass diagonal matrix.

Used by lumped explicit integrators. If mass lumping is impossible or approximate, adds scalar error to "error" parameter. Md += c*diag(M) or Md += c*HRZ(M) or other lumping heuristics This default implementation is VERY INEFFICIENT.

Reimplemented from chrono::fea::ChElementGeneric.

◆ EvaluateSectionDisplacement()

void chrono::fea::ChElementShellBST::EvaluateSectionDisplacement ( const double  u,
const double  v,
ChVector3d u_displ,
ChVector3d 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,
ChVector3d 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,
ChVector3d 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,
ChVector3d 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->GetNumCoordsPosLevel(), 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