-
Notifications
You must be signed in to change notification settings - Fork 11
/
abdcal.py
executable file
·375 lines (309 loc) · 11.4 KB
/
abdcal.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
"""
ABD matrix calculator for a given stacking sequence.
The methods in this file will call create a ABD matrix of a composit for given
ply properties and stacking sequence.
A.J.J. Lagerweij
COHMAS Mechanical Engineering KAUST
2020
"""
# Import external package.
import numpy as np
# Import local package.
import deformation
# Determine the stiffness matrix for an in plane stress ply.
def QPlaneStress(El, Et, nult, G):
r"""
Generate the plane stress local stiffness matrix.
Parameters
----------
El : float
Elastic modulus in the longitudional direction.
Et : float
Elastic modulus in the transverse direction.
nult : float
Poisson ratio in longitudional-transverse direction.
G : float
Shear modulus in longitudional-transverse directions.
Returns
-------
Q : matrix
Stiffness matrix in longitudional-transverse directions.
"""
nutl = nult*Et/El
Q = np.matrix([[El/(1-nult*nutl), (nult*Et)/(1-nult*nutl), 0],
[(nutl*El)/(1-nult*nutl), Et/(1-nult*nutl), 0],
[0, 0, G]])
return Q
# Determine the stiffness matrix for an in plane strain ply.
def QPlaneStrain(El, Et, nult, G):
r"""
Generate the plane strain local stiffness matrix.
.. warning::
Not yet implemented.
It raises a `NotImplementedError`.
Parameters
----------
El : float
Elastic modulus in the longitudional direction.
Et : float
Elastic modulus in the transverse direction.
nult : float
Poisson ratio in longitudional-transverse direction.
G : float
Shear modulus in longitudional-transverse directions.
Returns
-------
Q : matrix
Stiffness matrix in longitudional-transverse directions.
"""
raise NotImplementedError(
'QPlaneStrain is undefined. Use plane stress version QPlaneStress.')
return 0
# Rotate the stiffness matrix over a given angle.
def stiffness_rotation(stiffness, angle):
r"""
Rotate the stiffness matrix against given angle.
This rotates the stiffness matrix from local to the global axis sytem.
Use a negative angle to rotate from global to local system.
Parameters
----------
stiffness : matrix
The matrix that must be rotated.
angle : float
The rotation angle in degrees.
Returns
-------
stiffness_rot : matrix
A rotated version of the matrix.
"""
angle = angle * np.pi/180 # convert to radians
m = np.cos(angle)
n = np.sin(angle)
T1 = np.matrix([[m**2, n**2, 2*m*n],
[n**2, m**2, -2*m*n],
[-m*n, m*n, m**2-n**2]])
T2 = np.matrix([[m**2, n**2, m*n],
[n**2, m**2, -m*n],
[-2*m*n, 2*m*n, m**2-n**2]])
stiffness_rot = np.linalg.inv(T1) * stiffness * T2
return stiffness_rot
# Rotate the compliance matrix over a given angle.
def compliance_rotation(compliance, angle):
r"""
Rotate the compliance matrix over a given angle.
This rotates the complianc matrix from local to the global axis sytem.
Use a negative angle to rotate from global to local system.
Parameters
----------
compliance : matrix
The matrix that must be rotated.
angle : float
The rotation angle in degrees.
Returns
-------
stiffness_rot : matrix
A rotated version of the matrix.
"""
angle = angle * np.pi/180 # convert to radians
m = np.cos(angle)
n = np.sin(angle)
T1 = np.matrix([[m**2, n**2, 2*m*n],
[n**2, m**2, -2*m*n],
[-m*n, m*n, m**2-n**2]])
T2 = np.matrix([[m**2, n**2, m*n],
[n**2, m**2, -m*n],
[-2*m*n, 2*m*n, m**2-n**2]])
compliance_rot = np.linalg.inv(T2) * compliance * T1
return compliance_rot
# Calculate the stiffness matrix of a thin (no bending) laminate.
def abdthin(Q, angles, thickness, truncate=False):
r"""
ABD matrix calculator for a thin laminate.
In the thin laminate theroy it is assumed that the out of plane stiffness
is negligible and that the layup is symmetric. Hence only the membrane
(A part) of the ABD matrix remains. Top plies should be listed first in the
lists of Q, angles and thickness.
Parameters
----------
Q : list
The stiffness matrix of each ply in its l-t axis system.
angles : list
The rotation of each ply in degrees.
thickness : list
The thickness of each ply.
truncate : bool
Truncates very small numbers when true.
Returns
-------
C : matrix
The stiffness matrix of the thin laminate.
"""
# Create an empty stiffness matrix.
C = np.zeros((3, 3))
# Loop over all plies.
for i in range(len(angles)):
Q_bar = stiffness_rotation(Q[i], angles[i])
C += Q_bar * thickness[i]
# Go from the stresses * length^2 to running loads.
C = np.matrix(C) / np.sum(thickness)
# Truncate very small values.
if truncate is True:
C = np.matrix(np.where(np.abs(C) < np.max(C)*1e-6, 0, C))
return C
# Calculate the ABD matrix of a given laminate.
def abd(Q, angles, thickness, truncate=False):
r"""
Calculate the full ABD matrix of a laminate.
Top plies should be listed first in the lists of Q, angles and thickness.
Parameters
----------
Q : list
The stiffness matrix of each ply in its l-t axis system.
angles : list
The rotation of each ply in degrees.
thickness : list
The thickness of each ply.
truncate : bool
Truncates very small numbers when true.
Returns
-------
ABD : matrix
The stiffness matrix of the thin laminate.
"""
# Calculate the total thickness.
h = sum(thickness) / 2
# Create empty matricces for A B en D.
A = np.zeros((3, 3))
B = np.zeros((3, 3))
D = np.zeros((3, 3))
# Loop over all plies
for i in range(len(angles)):
# Calculate the z coordinates of the top and bottom of the ply.
z_top = np.sum(thickness[:i]) - h
z_bot = np.sum(thickness[:i+1]) - h
# Rotate the local stiffenss matrix.
Q_bar = stiffness_rotation(Q[i], angles[i])
# Calculate the contribution to the A, B and D matrix of this layer.
Ai = Q_bar * (z_bot - z_top)
Bi = 1/2 * Q_bar * (z_bot**2 - z_top**2)
Di = 1/3 * Q_bar * (z_bot**3 - z_top**3)
# Summ this layer to the previous ones.
A = A + Ai
B = B + Bi
D = D + Di
# Compile the entirety of the ABD matrix.
ABD = np.matrix([[A[0, 0], A[0, 1], A[0, 2], B[0, 0], B[0, 1], B[0, 2]],
[A[1, 0], A[1, 1], A[1, 2], B[1, 0], B[1, 1], B[1, 2]],
[A[2, 0], A[2, 1], A[2, 2], B[2, 0], B[2, 1], B[2, 2]],
[B[0, 0], B[0, 1], B[0, 2], D[0, 0], D[0, 1], D[0, 2]],
[B[1, 0], B[1, 1], B[1, 2], D[1, 0], D[1, 1], D[1, 2]],
[B[2, 0], B[2, 1], B[2, 2], D[2, 0], D[2, 1], D[2, 2]]])
# Truncate very small values.
if truncate is True:
ABD = np.matrix(np.where(np.abs(ABD) < np.max(ABD)*1e-6, 0, ABD))
return ABD
# Calculate the global coefficent of thermal expansion of a thin laminate.
def ctethin(Q, angles, thickness, alpha):
r"""
Thin laminate Coefficient of Thermal Expansion calculator.
This funcion calculates the CTE in the x-y axis sytem. It summs up the
rotated CTE of each layer and weights them by layer stiffness and
thickness. Here it is assumed that there is no bending behaviour this is
only true for thin and symmetric layups.
Top plies should be listed first in the lists of Q, angles, thickness and
alpha.
Parameters
----------
Q : list
The stiffness matrix of each ply in its l-t axis system.
angles : list
The rotation of each ply in degrees.
thickness : list
The thickness of each ply.
alpha : list
The coeficcient of termal expansion of each ply in l-t axis system.
Returns
-------
cte : vector
The coefficient of termal expansion of the laminate, in x-y axis sytem.
"""
# Calculate the laminate compliance matrix.
C = abdthin(Q, angles, thickness)
S = matrix_inverse(C)
# Create an empty CTE vector.
cte = np.zeros((3, 1))
# Loop ver all plies.
for i in range(len(angles)):
alpha_bar = deformation.strain_rotation(alpha[i], angles[i])
Q_bar = stiffness_rotation(Q[i], angles[i])
cte += Q_bar * alpha_bar * thickness[i]
# Normalize the cte to be in strain per temperature unit.
cte = cte / sum(thickness)
cte = S * cte
return cte
# Calculate the global coefficent of thermal expansion.
def cte(Q, angles, thickness, alpha):
r"""
Coefficient of Thermal Expansion calculator.
This funcion calculates the CTE in the x-y axis sytem. It summs up the
rotated CTE of each layer and multiplies them by layer stiffness and
thickness. The resulting CTE relates thermal change (:math:`\Delta T`) to
the deformation vector (strains and curvatures).
Top plies should be listed first in the lists of Q, angles,
thickness and alpha.
Parameters
----------
Q : list
The stiffness matrix of each ply in its l-t axis system.
angles : list
The rotation of each ply in degrees.
thickness : list
The thickness of each ply.
alpha : list
The coeficcient of termal expansion of each ply in l-t axis system.
Returns
-------
cte : vector
The coefficient of termal expansion of the laminate, in x-y axis sytem.
"""
# Calculate the total thickness.
h = sum(thickness) / 2
# Create an empty cte * stiffness vector vector.
alphaEt = np.zeros((6, 1))
# Loop ver all plies.
for i in range(len(angles)):
# Calculate the z coordinates of the top and bottom of the ply.
z_top = np.sum(thickness[:i]) - h
z_bot = np.sum(thickness[:i+1]) - h
# rotate the ply cte and Q f
alpha_bar = deformation.strain_rotation(alpha[i], angles[i])
Q_bar = stiffness_rotation(Q[i], angles[i])
# Calculate the contribution if this ply to the total thermal load.
alphaEt[:3] += Q_bar * alpha_bar * (z_bot - z_top)
alphaEt[3:] += 1/2 * Q_bar * alpha_bar * (z_bot**2 - z_top**2)
# Convert alphaEt which links Delta T and running loads to strain and
# curvature.
ABD = abd(Q, angles, thickness)
ABD_inv = matrix_inverse(ABD)
cte = deformation.load_applied(ABD_inv, alphaEt) # stren per unit change T
return cte
# Invert the matrix.
def matrix_inverse(matrix):
r"""
Inverts a matrix.
Parameters
----------
matrix : matrix
The matrix to be inverted.
Returns
-------
inverse : matrix
The inverse of the matrix.
"""
# Check if this ABD matrix is invertible.
if np.linalg.det(matrix) == 0:
error_message = 'Determinant is equal to zero inverting not posible'
raise ValueError(error_message)
# Invert the matrix.
inverse = np.linalg.inv(matrix)
return inverse