Property definitions#
Beam property module (pyfe3d.beamprop
)#
- class pyfe3d.beamprop.BeamProp#
General class to represent beam properties.
About shear correction factors, usually one should only apply the shear correction factor to the shear modulus
G
. It is, however, possible to have more complex scenarios where the different cross section parametersAy
,Az
,J
have different corrections applied due to geometric properties. It is the user’s responsability to apply these correction factors properly to achieve a general representation of the beam element behaviour by means of this property class. No shear correction factors or geometric correction factors are applied internally during the calculation of the structural matrices.- Attributes:
- A,float
Area of the cross section
- E,float
Young modulus
- G,float
Shear modulus
- Iyy,float
Second moment of area about the y axis \(\int_y \int_z z^2 dy dz\)
- Izz,float
Second moment of area about the z axis \(\int_y \int_z y^2 dy dz\)
- Iyz,float
Product moment of area \(\int_y \int_z y z dy dz\)
- J,float
Torsion stiffness or torsion constant, see https://en.wikipedia.org/wiki/Torsion_constant.
- Ay,float
Integral \(\int_y \int_z y dy dz\)
- Az,float
Integral \(\int_y \int_z z dy dz\)
- intrho,float
Integral \(\int_{y_e} \int_{z_e} \rho(y, z) dy dz\), where \(\rho\) Is the density
- intrhoy,float
Integral \(\int_y \int_z y \rho(y, z) dy dz\)
- intrhoz,float
Integral \(\int_y \int_z z \rho(y, z) dy dz\)
- intrhoy2,float
Integral \(\int_y \int_z y^2 \rho(y, z) dy dz\)
- intrhoz2,float
Integral \(\int_y \int_z z^2 \rho(y, z) dy dz\)
- intrhoyz,float
Integral \(\int_y \int_z y z \rho(y, z) dy dz\)
Notes
For beams with homogeneous material along the cross section some of the quantities above can be calculated as follows, assuming that
rho
is the material density:intrho = A*rho
intrhoy2 = Izz*rho
intrhoz2 = Iyy*rho
intrhoyz = Iyz*rho
- A#
A: ‘double’
- Ay#
Ay: ‘double’
- Az#
Az: ‘double’
- E#
E: ‘double’
- G#
G: ‘double’
- Iyy#
Iyy: ‘double’
- Iyz#
Iyz: ‘double’
- Izz#
Izz: ‘double’
- J#
J: ‘double’
- intrho#
intrho: ‘double’
- intrhoy#
intrhoy: ‘double’
- intrhoy2#
intrhoy2: ‘double’
- intrhoyz#
intrhoyz: ‘double’
- intrhoz#
intrhoz: ‘double’
- intrhoz2#
intrhoz2: ‘double’
Shell property module (pyfe3d.shellprop
)#
Highly based on the composites. module.
- class pyfe3d.shellprop.GradABDE#
Container to store the gradients of the ABDE matrices with respect to the lamination parameters
- Attributes:
- gradAij, gradBij, gradDij, gradEijtuple of 2D np.array objects
The shapes of these gradient matrices are:
gradAij: (6, 5) gradBij: (6, 5) gradDij: (6, 5) gradEij: (3, 3)
They contain the gradients of each laminate stiffness with respect to the thickness and respective lamination parameters. The rows and columns correspond to:
Methods
calc_LP_grad
(self, double thickness, ...)Gradients of the shell stiffnesses with respect to the thickness and lamination parameters
- calc_LP_grad(self, double thickness, MatLamina mat, LaminationParameters lp) void #
Gradients of the shell stiffnesses with respect to the thickness and lamination parameters
- Parameters:
- thicknessfloat
The total thickness of the laminate
- mat
MatLamina
object Material object
- lp
LaminationParameters
object The container class with all lamination parameters already defined
- Returns:
- None
The attributes of the object are updated.
- gradAij#
gradAij: ‘double[:, ::1]’
- gradBij#
gradBij: ‘double[:, ::1]’
- gradDij#
gradDij: ‘double[:, ::1]’
- gradEij#
gradEij: ‘double[:, ::1]’
- class pyfe3d.shellprop.Lamina#
- Attributes:
Methods
get_constitutive_matrix
(self)Return the constitutive matrix
Return displacement transformation matrix from lamina to laminate
Return stress transformation matrix from laminate to lamina
Return stress transformation matrix from lamina to laminate
rebuild
(self)Update constitutive matrices
- cos2t#
cos2t: ‘double’
- cos4t#
cos4t: ‘double’
- cost#
cost: ‘double’
- get_constitutive_matrix(self) double[:, ::1] #
Return the constitutive matrix
- get_transf_matrix_displ_to_laminate(self) double[:, ::1] #
Return displacement transformation matrix from lamina to laminate
- get_transf_matrix_stress_to_lamina(self) double[:, ::1] #
Return stress transformation matrix from laminate to lamina
- get_transf_matrix_stress_to_laminate(self) double[:, ::1] #
Return stress transformation matrix from lamina to laminate
- h#
h: ‘double’
- matlamina#
matlamina: pyfe3d.shellprop.MatLamina
- plyid#
plyid: ‘int’
- q11L#
q11L: ‘double’
- q12L#
q12L: ‘double’
- q16L#
q16L: ‘double’
- q22L#
q22L: ‘double’
- q26L#
q26L: ‘double’
- q44L#
q44L: ‘double’
- q45L#
q45L: ‘double’
- q55L#
q55L: ‘double’
- q66L#
q66L: ‘double’
- rebuild(self) void #
Update constitutive matrices
Reference:
Reddy, J. N., Mechanics of Laminated Composite Plates and Shells - Theory and Analysys. Second Edition. CRC PRESS, 2004.
- sin2t#
sin2t: ‘double’
- sin4t#
sin4t: ‘double’
- sint#
sint: ‘double’
- thetadeg#
thetadeg: ‘double’
- class pyfe3d.shellprop.LaminationParameters#
Lamination parameters
- Attributes:
- xiA1, xiA2, xiA3, xiA4float
Lamination parameters \(\xi_{Ai}\) (in-plane)
- xiB1, xiB2, xiB3, xiB4float
Lamination parameters \(\xi_{Bi}\) (in-plane coupling with bending)
- xiD1, xiD2, xiD3, xiD4float
Lamination parameters \(\xi_{Di}\) (bending)
- xiE1, xiE2float
Lamination parameters \(\xi_{Ei}\) (transverse shear)
- xiA1#
xiA1: ‘double’
- xiA2#
xiA2: ‘double’
- xiA3#
xiA3: ‘double’
- xiA4#
xiA4: ‘double’
- xiB1#
xiB1: ‘double’
- xiB2#
xiB2: ‘double’
- xiB3#
xiB3: ‘double’
- xiB4#
xiB4: ‘double’
- xiD1#
xiD1: ‘double’
- xiD2#
xiD2: ‘double’
- xiD3#
xiD3: ‘double’
- xiD4#
xiD4: ‘double’
- xiE1#
xiE1: ‘double’
- xiE2#
xiE2: ‘double’
- xiE3#
xiE3: ‘double’
- xiE4#
xiE4: ‘double’
- class pyfe3d.shellprop.MatLamina#
Orthotropic material lamina
- Attributes:
e1
floate1: ‘double’
e2
floate2: ‘double’
g12
floatg12: ‘double’
g13
floatg13: ‘double’
g23
floatg23: ‘double’
nu12
nu12: ‘double’
nu13
nu13: ‘double’
nu23
nu23: ‘double’
nu21
nu21: ‘double’
nu31
nu31: ‘double’
nu32
nu32: ‘double’
rho
rho: ‘double’
a1
a1: ‘double’
a2
a2: ‘double’
a3
a3: ‘double’
tref
tref: ‘double’
- st1,st2
allowable tensile stresses for directions 1 and 2
- sc1,sc2
allowable compressive stresses for directions 1 and 2
ss12
ss12: ‘double’
q11
q11: ‘double’
q12
q12: ‘double’
q13
q13: ‘double’
q21
q21: ‘double’
q22
q22: ‘double’
q23
q23: ‘double’
q31
q31: ‘double’
q32
q32: ‘double’
q33
q33: ‘double’
q44
q44: ‘double’
q55
q55: ‘double’
q66
q66: ‘double’
- ci
lamina stiffness constants
- ui
lamina material invariants
Methods
get_constitutive_matrix
(self)Return the constitutive matrix
get_invariant_matrix
(self)Return the invariant matrix
rebuild
(self)Update constitutive and invariant terms
Trace-normalize the lamina properties for plane stress
Notes
For isotropic materials when the user defines \(\nu\) and \(E\), \(G\) will be recaculated based on equation: \(G = E/(2 \times (1+\nu))\); in a lower priority if the user defines \(\nu\) and \(G\), \(E\) will be recaculated based on equation: \(E = 2 \times (1+\nu) \times G\).
- a1#
a1: ‘double’
- a2#
a2: ‘double’
- a3#
a3: ‘double’
- c11#
c11: ‘double’
- c12#
c12: ‘double’
- c13#
c13: ‘double’
- c22#
c22: ‘double’
- c23#
c23: ‘double’
- c33#
c33: ‘double’
- c44#
c44: ‘double’
- c55#
c55: ‘double’
- c66#
c66: ‘double’
- e1#
e1: ‘double’
- e2#
e2: ‘double’
- e3#
e3: ‘double’
- g12#
g12: ‘double’
- g13#
g13: ‘double’
- g23#
g23: ‘double’
- get_constitutive_matrix(self) double[:, ::1] #
Return the constitutive matrix
- get_invariant_matrix(self) double[:, ::1] #
Return the invariant matrix
- nu12#
nu12: ‘double’
- nu13#
nu13: ‘double’
- nu21#
nu21: ‘double’
- nu23#
nu23: ‘double’
- nu31#
nu31: ‘double’
- nu32#
nu32: ‘double’
- q11#
q11: ‘double’
- q12#
q12: ‘double’
- q13#
q13: ‘double’
- q21#
q21: ‘double’
- q22#
q22: ‘double’
- q23#
q23: ‘double’
- q31#
q31: ‘double’
- q32#
q32: ‘double’
- q33#
q33: ‘double’
- q44#
q44: ‘double’
- q55#
q55: ‘double’
- q66#
q66: ‘double’
- rebuild(self) void #
Update constitutive and invariant terms
Reference:
Reddy, J. N., Mechanics of laminated composite plates and shells. Theory and analysis. Second Edition. CRC Press, 2004.
- rho#
rho: ‘double’
- sc1#
sc1: ‘double’
- sc2#
sc2: ‘double’
- ss12#
ss12: ‘double’
- st1#
st1: ‘double’
- st2#
st2: ‘double’
- trace_normalize_plane_stress(self) void #
Trace-normalize the lamina properties for plane stress
Modify the original
MatLamina
object with a trace-normalization performed after calculating the trace according to Eq. 1 of reference:Melo, J. D. D., Bi, J., and Tsai, S. W., 2017, “A Novel Invariant-Based Design Approach to Carbon Fiber Reinforced Laminates,” Compos. Struct., 159, pp. 44–52.
The trace calculated as \(tr = Q_{11} + Q_{22} + 2Q_{66}\). The universal in-plane stress stiffness components \(Q_{11},Q_{12},Q_{22},Q_{44},Q_{55},Q_{66}\) are divided by \(tr\), and the invariants \(U_1,U_2,U_3,U_4,U_5,U_6,U_7\) are calculated with the normalized stiffnesses, such they also become trace-normalized invariants. These can be accessed using the
u1,u2,u3,u4,u5,u6,u7
attributes.
- tref#
tref: ‘double’
- u1#
u1: ‘double’
- u2#
u2: ‘double’
- u3#
u3: ‘double’
- u4#
u4: ‘double’
- u5#
u5: ‘double’
- u6#
u6: ‘double’
- u7#
u7: ‘double’
- class pyfe3d.shellprop.ShellProp#
- Attributes:
plies
listplies: list
stack
liststack: list
h
floath: ‘double’
offset
floatoffset: ‘double’
- e1, e2float
Equivalent laminate moduli in directions 1 and 2
g12
floatg12: ‘double’
- nu12, nu21float
Equivalent laminate Poisson ratios in the 12 and 21 directions
- scf_k13, scf_k23float
Shear correction factor in the 13 and 23 directions
intrho
floatintrho: ‘double’
intrhoz
floatintrhoz: ‘double’
intrhoz2
floatintrhoz2: ‘double’
Methods
calc_constitutive_matrix
(self)Calculate the laminate constitutive terms
Calculate the equivalent laminate properties
Calculate the lamination parameters.
calc_scf
(self)Update shear correction factors of the
ShellProp
objectforce_balanced
(self)Force a balanced laminate
force_orthotropic
(self)Force an orthotropic laminate
force_symmetric
(self)Force a symmetric laminate
- A11#
A11: ‘double’
- A12#
A12: ‘double’
- A16#
A16: ‘double’
- A22#
A22: ‘double’
- A26#
A26: ‘double’
- A66#
A66: ‘double’
- B11#
B11: ‘double’
- B12#
B12: ‘double’
- B16#
B16: ‘double’
- B22#
B22: ‘double’
- B26#
B26: ‘double’
- B66#
B66: ‘double’
- D11#
D11: ‘double’
- D12#
D12: ‘double’
- D16#
D16: ‘double’
- D22#
D22: ‘double’
- D26#
D26: ‘double’
- D66#
D66: ‘double’
- E44#
E44: ‘double’
- E45#
E45: ‘double’
- E55#
E55: ‘double’
- calc_constitutive_matrix(self) void #
Calculate the laminate constitutive terms
This is the commonly called
ABD
matrix withshape=(6, 6)
when the classical laminated plate theory is used, or theABDE
matrix when the first-order shear deformation theory is used, containing the transverse shear terms.
- calc_equivalent_properties(self) void #
Calculate the equivalent laminate properties
The following attributes are updated:
e1
,e2
,g12
,`u12
,nu21
- calc_lamination_parameters(self) LaminationParameters #
Calculate the lamination parameters.
The following attributes are calculated:
xiA
,xiB
,xiD
,xiE
- calc_scf(self) void #
Update shear correction factors of the
ShellProp
objectReference:
Vlachoutsis, S. “Shear correction factors for plates and shells”, Int. Journal for Numerical Methods in Engineering, Vol. 33, 1537-1552, 1992.
http://onlinelibrary.wiley.com/doi/10.1002/nme.1620330712/full
Using “one shear correction factor” (see reference), assuming:
constant G13, G23, E1, E2, nu12, nu21 within each ply
g1 calculated using z at the middle of each ply
zn1 =
ShellProp
offset
attribute
- Returns:
- k13, k23tuple
Shear correction factors. Also updates attributes:
scf_k13
andscf_k23
.
- e1#
e1: ‘double’
- e2#
e2: ‘double’
- force_balanced(self) void #
Force a balanced laminate
The attributes \(A_{16}\), \(A_{26}\), \(B_{16}\), \(B_{26}\) are set to zero to force a balanced laminate.
- force_orthotropic(self) void #
Force an orthotropic laminate
The attributes \(A_{16}\), \(A_{26}\), \(B_{16}\), \(B_{26}\), \(D_{16}\), \(D_{26}\) are set to zero to force an orthotropic laminate.
- force_symmetric(self) void #
Force a symmetric laminate
The \(B_{ij}\) terms of the constitutive matrix are set to zero.
- g12#
g12: ‘double’
- h#
h: ‘double’
- intrho#
intrho: ‘double’
- intrhoz#
intrhoz: ‘double’
- intrhoz2#
intrhoz2: ‘double’
- nu12#
nu12: ‘double’
- nu21#
nu21: ‘double’
- offset#
offset: ‘double’
- plies#
plies: list
- scf_k13#
scf_k13: ‘double’
- scf_k23#
scf_k23: ‘double’
- stack#
stack: list
- pyfe3d.shellprop.force_balanced_LP(LaminationParameters lp) LaminationParameters #
Force balanced lamination parameters
The lamination parameters \(\xi_{A2}\) and \(\xi_{A4}\) are set to null to force a balanced laminate.
- pyfe3d.shellprop.force_orthotropic_LP(LaminationParameters lp) LaminationParameters #
Force orthotropic lamination parameters
The lamination parameters \(\xi_{A2}\), \(\xi_{A4}\), \(\xi_{B2}\), \(\xi_{B4}\), \(\xi_{D2}\) and \(\xi_{D4}\) are set to null to force an orthotropic laminate. The \(\xi_{D2}\) and \(\xi_{D4}\) are related to the bend-twist coupling and become often very small for balanced laminates with a large amount of plies.
- pyfe3d.shellprop.force_symmetric_LP(LaminationParameters lp) LaminationParameters #
Force symmetric lamination parameters
The lamination parameters \(\xi_{Bi}\) are set to null to force a symmetric laminate.
- pyfe3d.shellprop.shellprop_from_LaminationParameters(double thickness, MatLamina mat, LaminationParameters lp) ShellProp #
Return a
ShellProp
object based in the thickness, material and lamination parameters- Parameters:
- thicknessfloat
The total thickness of the laminate
- mat
MatLamina
object Material object
- lp
LaminationParameters
object The container class with all lamination parameters already defined
- Returns:
- lam
ShellProp
laminate with the ABD and ABDE matrices already calculated
- lam
- pyfe3d.shellprop.shellprop_from_lamination_parameters(double thickness, MatLamina matlamina, double xiA1, double xiA2, double xiA3, double xiA4, double xiB1, double xiB2, double xiB3, double xiB4, double xiD1, double xiD2, double xiD3, double xiD4, double xiE1=0, double xiE2=0) ShellProp #
Return a
ShellProp
object based in the thickness, material and lamination parametersNote that \(\xi_{E1}\) and \(\xi_{E2}\) are optional and usually equal to zero, becoming important only when the transverse shear modulus is different in the two directions, i.e. when \(G_{13} \ne G{23}\).
- Parameters:
- thicknessfloat
The total thickness of the plate
- matlamina
MatLamina
object Material object
- xiAj, xiBj, xiDj, xiEjfloat
The 14 lamination parameters according to the first-order shear deformation theory: \(\xi_{A1} \cdots \xi_{A4}\), \(\xi_{B1} \cdots \xi_{B4}\), \(\xi_{D1} \cdots \xi_{D4}\), \(\xi_{E1}\) and \(\xi_{E2}\)
- Returns:
- lam
ShellProp
Shell property with the ABD and ABDE matrices already calculated
- lam
Shell property utils module (pyfe3d.shellprop_utils
)#
Highly based on the composites. module.
- pyfe3d.shellprop_utils.isotropic_plate(thickness, E, nu, offset=0.0, calc_scf=True, rho=0.0)#
Read data for an isotropic plate
ShellProp
object is returned based on the inputs given.- Parameters:
- thicknessfloat
Plate thickness.
- Efloat
Young modulus.
- nufloat, optional
Poisson’s ratio.
- rhofloat, optional
Material density
- offsetfloat, optional
Offset along the normal axis about the mid-surface, which influences the extension-bending coupling (B matrix).
- calc_scfbool, optional
If True, use
ShellProp.calc_scf()
to compute shear correction factors, otherwise the default value of 5/6 is used.
- pyfe3d.shellprop_utils.laminated_plate(stack, plyt=None, laminaprop=None, rho=0.0, plyts=None, laminaprops=None, rhos=None, offset=0.0, calc_scf=True)#
Read a laminate stacking sequence data.
ShellProp
object is returned based on the inputs given.- Parameters:
- stacklist
Angles of the stacking sequence in degrees.
- plytfloat, optional
When all plies have the same thickness,
plyt
can be supplied.- laminaproptuple, optional
When all plies have the same material properties,
laminaprop
can be supplied.- rhofloat, optional
Uniform material density to be used for all plies.
- plytslist, optional
A list of floats with the thickness of each ply.
- laminapropslist, optional
A list of tuples with a laminaprop for each ply.
- rhoslist, optional
A list of floats with the material density of each ply.
- offsetfloat, optional
Offset along the normal axis about the mid-surface, which influences the laminate properties.
- calc_scfbool, optional
If True, use
ShellProp.calc_scf()
to compute shear correction factors, otherwise the default value of 5/6 is used
Notes
plyt
orplyts
must be suppliedlaminaprop
orlaminaprops
must be suppliedFor orthotropic plies, the
laminaprop
should be:laminaprop = (E11, E22, nu12, G12, G13, G23)
For isotropic plies, the
laminaprop
should be:laminaprop = (E, nu)
- pyfe3d.shellprop_utils.read_laminaprop(laminaprop, rho=0)#
Returns a
MatLamina
object based on an inputlaminaprop
tuple- Parameters:
- laminaproplist or tuple
For the most general case of tri-axial stress, use a tuple containing the folliwing entries:
laminaprop = (e1, e2, nu12, g12, g13, g23, e3, nu13, nu23)
For isotropic materials aiming calculations with tri-axial stresses, use:
g = e/(2*(1+nu)) laminaprop = (e, e, nu, g, g, g, e, nu, nu)
For othotropic materials with in-plane stresses the user can only supply:
laminaprop = (e1, e2, nu12, g12, g13, g23)
For isotropic materials with in-plane stresses the user can only supply:
laminaprop = (e, nu) # new
symbol
value
e1
Young Module in direction 1
e2
Young Module in direction 2
nu12
12 Poisson’s ratio
g12
12 Shear Modulus
g13
13 Shear Modulus
g23
13 Shear Modulus
e3
Young Module in direction 3
nu13
13 Poisson’s ratio
nu23
23 Poisson’s ratio
- rhofloat, optional
Material density
- Returns:
- matlamMatLamina
A
MatLamina
object.