Description

Laminated thick shell with geometrically exact kinematics, with 4 nodes.

It generalizes the Reissner thick shell theory (in fact each layer requires a ChMaterialShellReissner) by using the Chroscielewski 6-dof field shell theory as discussed in:

Wojciech Witkowski, "4-Node combined shell element with semi-EAS-ANS strain interpolations in 6-parameter shell theories with drilling degrees of freedom", Comp.Mech 2009.

This specific implementation is based on the paper:

Marco Morandini, Pierangelo Masarati, "Implementation and validation of a 4-node shell finite element", IDETC/CIE 2014.

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

        v
        ^
        |
D o-----+-----o C
  |     |     |
--+-----+-----+----> u
  |     |     |
A o-----+-----o B

#include <ChElementShellReissner4.h>

Inheritance diagram for chrono::fea::ChElementShellReissner4:
Collaboration diagram for chrono::fea::ChElementShellReissner4:

Classes

class  Layer
 Definition of a layer. More...
 

Public Types

enum  IntegrationPoint {
  IP_1_1 = 0, IP_1_2 = 1, IP_1_3 = 2, IP_2_1 = 3,
  IP_2_2 = 4, IP_2_3 = 5, IP_3_1 = 6, IP_3_2 = 7,
  IP_3_3 = 8, NUMIP = 4
}
 
enum  ShearStrainEvaluationPoint {
  SSEP_1 = 0, SSEP_2 = 1, SSEP_3 = 2, SSEP_4 = 3,
  NUMSSEP = 4
}
 
enum  NodeName {
  NODE1 = 0, NODE2 = 1, NODE3 = 2, NODE4 = 3,
  NUMNODES = 4
}
 
enum  Deformations { STRAIN = 0, CURVAT = 1, NUMDEFORM = 2 }
 
enum  InnerEASdofs { IDOFS = 7 }
 
using ShapeVector = ChMatrixNM< double, 1, 4 >
 

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 GetNodeNdofs (int n) override
 Get the number of coordinates from the n-th node used by this element.
 
void SetNodes (std::shared_ptr< ChNodeFEAxyzrot > nodeA, std::shared_ptr< ChNodeFEAxyzrot > nodeB, std::shared_ptr< ChNodeFEAxyzrot > nodeC, std::shared_ptr< ChNodeFEAxyzrot > nodeD)
 Specify the nodes of this element. More...
 
virtual std::shared_ptr< ChNodeFEAbaseGetNodeN (int n) override
 Access the n-th node of this element.
 
std::shared_ptr< ChNodeFEAxyzrotGetNodeA () const
 Get a handle to the first node of this element.
 
std::shared_ptr< ChNodeFEAxyzrotGetNodeB () const
 Get a handle to the second node of this element.
 
std::shared_ptr< ChNodeFEAxyzrotGetNodeC () const
 Get a handle to the third node of this element.
 
std::shared_ptr< ChNodeFEAxyzrotGetNodeD () const
 Get a handle to the fourth node of this element.
 
void SetAsNeutral ()
 Sets the neutral rotations of nodes A,B,C,D, at once, assuming the current element position is for zero strain.
 
void AddLayer (double thickness, double theta, std::shared_ptr< ChMaterialShellReissner > 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 GetLengthX () const
 Set the structural damping: this is the Rayleigh "alpha" OBSOLETE create a ChDampingReissnerRayleigh object and add to layer material to have the same effect void SetAlphaDamp(double a) { m_Alpha = a; }
More...
 
double GetLengthY () const
 Get the element length in the Y direction.
 
double GetThickness ()
 Get the total thickness of the shell element (might be sum of multiple layer thicknesses)
 
ChQuaternion GetAvgRot ()
 
void ShapeFunctions (ShapeVector &N, double x, double y)
 Fills the N shape function matrix.
 
void ShapeFunctionsDerivativeX (ShapeVector &Nx, double x, double y)
 Fills the Nx shape function derivative matrix with respect to X.
 
void ShapeFunctionsDerivativeY (ShapeVector &Ny, double x, double y)
 Fills the Ny shape function derivative matrix with respect to Y.
 
ChVector EvaluateGP (int igp)
 
ChVector EvaluatePT (int ipt)
 
virtual unsigned int iGetNumDof (void) const
 Inner EAS dofs.
 
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 ComputeMmatrixGlobal (ChMatrixRef M) override
 Set M as the global mass matrix.
 
virtual void ComputeInternalForces (ChVectorDynamic<> &Fi) override
 Computes the internal forces. More...
 
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...
 
void ComputeInternalJacobians (double Kfactor, double Rfactor)
 Compute Jacobians of the internal forces. More...
 
void ComputeMassMatrix ()
 Compute the mass matrix of the element. 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 EvaluateSectionVelNorm (double U, double V, ChVector<> &Result) override
 Virtual method to plot velocity field distribution. More...
 
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...
 
- Public Member Functions inherited from chrono::fea::ChElementShell
virtual bool IsTriangleShell ()
 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.
 
- 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 ComputeGravityForces (ChVectorDynamic<> &Fg, const ChVector<> &G_acc) override
 Compute the gravitational forces. 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 int GetNdofs_active ()
 Get the actual number of active degrees of freedom. More...
 
virtual int GetNodeNdofs_active (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::ChLoadableUV
virtual bool IsTriangleIntegrationNeeded ()
 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.
 
- 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.
 

Public Attributes

std::vector< std::shared_ptr< ChNodeFEAxyzrot > > 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
 
double m_lenX
 element length in X direction
 
double m_lenY
 element length in Y direction
 
ChMatrixNM< double, 24, 24 > m_MassMatrix
 mass matrix
 
ChMatrixNM< double, 24, 24 > m_JacobianMatrix
 Jacobian matrix (Kfactor*[K] + Rfactor*[R])
 
ChVariablesGenericDiagonalMassmvariables
 
ChVector xa_0 [NUMNODES]
 
ChVector xa [NUMNODES]
 
ChMatrix33 iTa [NUMNODES]
 
ChMatrix33 iTa_i [NUMIP]
 
ChMatrix33 iTa_A [NUMSSEP]
 
ChVector phi_tilde_n [NUMNODES]
 
ChVector phi_tilde_i [NUMIP]
 
ChVector phi_tilde_A [NUMSSEP]
 
ChVector phi_tilde_0
 
ChMatrix33 T0_overline
 
ChMatrix33 T_overline
 
ChMatrix33 T_0_0
 
ChMatrix33 T_0_i [NUMIP]
 
ChMatrix33 T_0_A [NUMSSEP]
 
ChMatrix33 T_0
 
ChMatrix33 T_i [NUMIP]
 
ChMatrix33 T_A [NUMSSEP]
 
ChMatrix33 Phi_Delta_i [NUMIP][NUMNODES]
 
ChMatrix33 Phi_Delta_A [NUMIP][NUMNODES]
 
ChMatrix33 Kappa_delta_i_1 [NUMIP][NUMNODES]
 
ChMatrix33 Kappa_delta_i_2 [NUMIP][NUMNODES]
 
ChMatrix33 Q_i [NUMIP]
 
ChMatrix33 Q_A [NUMSSEP]
 
ChVector k_1_i [NUMIP]
 
ChVector k_2_i [NUMIP]
 
ChVector eps_tilde_1_0_i [NUMIP]
 
ChVector eps_tilde_2_0_i [NUMIP]
 
ChVector eps_tilde_1_0_A [NUMSSEP]
 
ChVector eps_tilde_2_0_A [NUMSSEP]
 
ChVector eps_tilde_1_i [NUMIP]
 
ChVector eps_tilde_2_i [NUMIP]
 
ChVector eps_tilde_1_A [NUMSSEP]
 
ChVector eps_tilde_2_A [NUMSSEP]
 
ChVector k_tilde_1_0_i [NUMIP]
 
ChVector k_tilde_2_0_i [NUMIP]
 
ChVector k_tilde_1_i [NUMIP]
 
ChVector k_tilde_2_i [NUMIP]
 
ChMatrixNM< double, 2, 2 > S_alpha_beta_0
 
ChMatrixNM< double, 2, 2 > S_alpha_beta_i [NUMIP]
 
ChMatrixNM< double, 2, 2 > S_alpha_beta_A [NUMSSEP]
 
double alpha_0
 
double alpha_i [NUMIP]
 
ChMatrixNM< double, 4, 2 > L_alpha_beta_i [NUMIP]
 
ChMatrixNM< double, 4, 2 > L_alpha_beta_A [NUMSSEP]
 
ChMatrixNM< double, 12, 24 > B_overline_i [NUMIP]
 
ChMatrixNM< double, 15, 24 > D_overline_i [NUMIP]
 
ChMatrixNM< double, 15, 15 > G_i [NUMIP]
 
ChMatrixNM< double, 12, IDOFS > P_i [NUMIP]
 
ChMatrixNM< double, IDOFS, IDOFS > K_beta_beta_i [NUMIP]
 
ChVector y_i_1 [NUMIP]
 
ChVector y_i_2 [NUMIP]
 
ChVectorN< double, IDOFS > beta
 
ChVectorN< double, 12 > epsilon_hat
 
ChVectorN< double, 12 > epsilon
 
ChVectorN< double, 12 > stress_i [NUMIP]
 
bool bFirstRes
 

Static Public Attributes

static double xi_i [NUMIP][2]
 
static double w_i [NUMIP] = {1., 1., 1., 1.}
 
static double xi_A [NUMSSEP][2] = {{0., 1.}, {-1., 0.}, {0., -1.}, {1., 0.}}
 
static double xi_n [NUMNODES][2]
 
static double xi_0 [2] = {0., 0.}
 

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::ChElementShellReissner4::AddLayer ( double  thickness,
double  theta,
std::shared_ptr< ChMaterialShellReissner 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::ChElementShellReissner4::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.

◆ ComputeInternalJacobians()

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

◆ ComputeMassMatrix()

void chrono::fea::ChElementShellReissner4::ComputeMassMatrix ( )

Compute the mass matrix of the element.

Note: in this 'basic' implementation, constant section and constant material are assumed

◆ ComputeNF() [1/2]

void chrono::fea::ChElementShellReissner4::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::ChElementShellReissner4::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::ChElementShellReissner4::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::ChElementShellReissner4::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::ChElementShellReissner4::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::ChElementShellReissner4::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::ChElementShellReissner4::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::ChElementShellReissner4::GetDensity ( )
overridevirtual

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

Density is mass per unit surface.

Implements chrono::ChLoadableUVW.

◆ GetLengthX()

double chrono::fea::ChElementShellReissner4::GetLengthX ( ) const
inline

Set the structural damping: this is the Rayleigh "alpha" OBSOLETE create a ChDampingReissnerRayleigh object and add to layer material to have the same effect void SetAlphaDamp(double a) { m_Alpha = a; }

Get the element length in the X direction.

◆ GetStateBlock()

void chrono::fea::ChElementShellReissner4::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 Rx_a Rx_a Rx_a x_b y_b z_b Rx_b Ry_b Rz_b}

Implements chrono::fea::ChElementBase.

◆ SetLayerZreferenceCentered()

void chrono::fea::ChElementShellReissner4::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::ChElementShellReissner4::SetNodes ( std::shared_ptr< ChNodeFEAxyzrot nodeA,
std::shared_ptr< ChNodeFEAxyzrot nodeB,
std::shared_ptr< ChNodeFEAxyzrot nodeC,
std::shared_ptr< ChNodeFEAxyzrot nodeD 
)

Specify the nodes of this element.

The node numbering is in ccw fashion as in the following scheme: v ^ D o--—+--—o C | | | –+--—+--—+-> u | | | A o--—+--—o B

Member Data Documentation

◆ xi_i

double chrono::fea::ChElementShellReissner4::xi_i
static
Initial value:
= {{-1. / std::sqrt(3.), -1. / std::sqrt(3.)},
{1. / std::sqrt(3.), -1. / std::sqrt(3.)},
{1. / std::sqrt(3.), 1. / std::sqrt(3.)},
{-1. / std::sqrt(3.), 1. / std::sqrt(3.)}}

◆ xi_n

double chrono::fea::ChElementShellReissner4::xi_n
static
Initial value:
= {{1., 1.},
{-1., 1.},
{-1., -1.},
{1., -1.}}

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