Description
Elasticity of 6-field Reissner-Mindlin shells (kinematically-exact shell theory as in Witkowski et al.) to be used in a ChMaterialShellReissner.
This class implements material properties for a layer from the Reissner theory, for the case of orthotropic linear elastic material. This is useful for laminated shells. One direction can be made softer than the other. Note that the angle and the thickness are defined when adding a material with this elasticity to a shell finite element (ex. ChElementShellReissner4) as a layer. Previously: ChMaterialShellReissnerOrthotropic
#include <ChMaterialShellReissner.h>
Public Member Functions | |
ChElasticityReissnerOrthotropic (double m_E_x, double m_E_y, double m_nu_xy, double m_G_xy, double m_G_xz, double m_G_yz, double m_alpha=1.0, double m_beta=0.1) | |
Construct an orthotropic material. More... | |
ChElasticityReissnerOrthotropic (double m_E, double m_nu, double m_alpha=1.0, double m_beta=0.1) | |
Construct an orthotropic material as sub case isotropic. More... | |
double | GetYoungModulusX () const |
Return the elasticity moduli, on x. | |
double | GetYoungModulusY () const |
Return the elasticity moduli, on y. | |
double | GetPoissonRatioXY () const |
Return the Poisson ratio, for xy. | |
double | GetPoissonRatioYX () const |
Return the Poisson ratio, for yx (follows xy as it must be nu_yx*E_x = nu_xy*E_y) | |
double | GetShearModulusXY () const |
Return the shear mod, in plane. | |
double | GetShearModulusXZ () const |
Return the shear mod, transverse. | |
double | GetShearModulusYZ () const |
Return the shear mod, transverse. | |
double | GetShearFactor () const |
Return the shear factor. | |
double | GetTorqueFactor () const |
Return the torque factor. | |
virtual void | ComputeStress (ChVector3d &n_u, ChVector3d &n_v, ChVector3d &m_u, ChVector3d &m_v, const ChVector3d &eps_u, const ChVector3d &eps_v, const ChVector3d &kur_u, const ChVector3d &kur_v, const double z_inf, const double z_sup, const double angle) |
The FE code will evaluate this function to compute per-unit-length forces/torques given the u,v strains/curvatures. More... | |
virtual void | ComputeStiffnessMatrix (ChMatrixRef mC, const ChVector3d &eps_u, const ChVector3d &eps_v, const ChVector3d &kur_u, const ChVector3d &kur_v, const double z_inf, const double z_sup, const double angle) |
/// Compute the 12x12 stiffness matrix [Km] , that is [ds/de], the tangent of the constitutive relation stresses/strains. More... | |
Additional Inherited Members | |
Public Attributes inherited from chrono::fea::ChElasticityReissner | |
ChMaterialShellReissner * | section |
Constructor & Destructor Documentation
◆ ChElasticityReissnerOrthotropic() [1/2]
chrono::fea::ChElasticityReissnerOrthotropic::ChElasticityReissnerOrthotropic | ( | double | m_E_x, |
double | m_E_y, | ||
double | m_nu_xy, | ||
double | m_G_xy, | ||
double | m_G_xz, | ||
double | m_G_yz, | ||
double | m_alpha = 1.0 , |
||
double | m_beta = 0.1 |
||
) |
Construct an orthotropic material.
- Parameters
-
m_E_x Young's modulus on x m_E_y Young's modulus on y m_nu_xy Poisson ratio xy (for yx it holds: nu_yx*E_x = nu_xy*E_y) m_G_xy Shear modulus, in plane m_G_xz Shear modulus, transverse m_G_yz Shear modulus, transverse m_alpha shear factor m_beta torque factor
◆ ChElasticityReissnerOrthotropic() [2/2]
chrono::fea::ChElasticityReissnerOrthotropic::ChElasticityReissnerOrthotropic | ( | double | m_E, |
double | m_nu, | ||
double | m_alpha = 1.0 , |
||
double | m_beta = 0.1 |
||
) |
Construct an orthotropic material as sub case isotropic.
- Parameters
-
m_E Young's modulus on x m_nu Poisson ratio m_alpha shear factor m_beta torque factor
Member Function Documentation
◆ ComputeStiffnessMatrix()
|
virtual |
/// Compute the 12x12 stiffness matrix [Km] , that is [ds/de], the tangent of the constitutive relation stresses/strains.
- Parameters
-
mC tangent matrix eps_u strains along u direction eps_v strains along v direction kur_u curvature along u direction kur_v curvature along v direction z_inf layer lower z value (along thickness coord) z_sup layer upper z value (along thickness coord) angle layer angle respect to x (if needed) -not used in this, isotropic
Reimplemented from chrono::fea::ChElasticityReissner.
◆ ComputeStress()
|
virtual |
The FE code will evaluate this function to compute per-unit-length forces/torques given the u,v strains/curvatures.
- Parameters
-
n_u forces along u direction (per unit length) n_v forces along v direction (per unit length) m_u torques along u direction (per unit length) m_v torques along v direction (per unit length) eps_u strains along u direction eps_v strains along v direction kur_u curvature along u direction kur_v curvature along v direction z_inf layer lower z value (along thickness coord) z_sup layer upper z value (along thickness coord) angle layer angle respect to x (if needed) -not used in this, isotropic
Implements chrono::fea::ChElasticityReissner.
The documentation for this class was generated from the following files:
- /builds/uwsbel/chrono/src/chrono/fea/ChMaterialShellReissner.h
- /builds/uwsbel/chrono/src/chrono/fea/ChMaterialShellReissner.cpp