Newer
Older
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
376
377
378
379
380
381
382
import surface_geom_SEM as sgs
import element_stiff_matrix as esm
import numpy as np
from geomdl import exchange
ex = np.array((1, 0, 0), dtype=np.int32)
ey = np.array((0, 1, 0), dtype=np.int32)
ez = np.array((0, 0, 1), dtype=np.int32)
class SurfaceGeneratedSEM(object):
'''
'''
def __init__(self, nurbs_parent_surface, lobatto_pw, node_1_ub,\
node_1_vb, node_3_ub, node_3_vb):
self.parent_nurbs_surf = nurbs_parent_surface
self.lobatto_pw = lobatto_pw
self.node_1_u = node_1_ub
self.node_1_v = node_1_vb
self.node_3_u = node_3_ub
self.node_3_v = node_3_vb
self.thk = nurbs_parent_surface.thk
# The more Pythonic way to deal with attributes is using decorators.
# However, as I won't change the attributes outside of the class I
# do not use it to make the class simpler.
# self._parent_nurbs_surf = nurbs_surface
# self._lobatto_pw = lobatto_pw
# self._node_1_u = node_1_ub
# self._node_1_v = node_1_vb
# self._node_3_u = node_3_ub
# self._node_3_v = node_3_vb
# @property
# def parent_nurbs_surf(self):
# return self._parent_nurbs_surf
# @property
# def lobatto_pw(self):
# return self._lobatto_pw
# @property
# def node_1_u(self):
# return self._node_1_u
# @property
# def node_1_v(self):
# return self._node_1_v
# @property
# def node_3_u(self):
# return self._node_3_u
# @property
# def node_3_v(self):
# return self._node_3_v
# @parent_nurbs_surf.setter
# def parent_nurbs_surf(self, value):
# self._parent_nurbs_surf = value
# @node_1_u.setter
# def node_1_u(self, value):
# self._node_1_u = value
# @node_1_v.setter
# def node_1_v(self, value):
# self._node_1_v = value
# @node_3_u.setter
# def node_3_u(self, value):
# self._node_3_u = value
# @node_3_v.setter
# def node_3_v(self, value):
# self._node_3_v = value
def nodes_physical_coordinate(self):
'''This function generate the physical coordinate of the
lobatto points.
-Output:
A nxnx3 numpy array, in that n is the number of lobatto points
'''
dim = self.lobatto_pw.shape[0]
nodes_coor_mtx = np.zeros((dim, dim, 3))
for i in range(dim):
xi2 = self.lobatto_pw[i, 0]
for j in range(dim):
xi1 = self.lobatto_pw[j, 0]
physical_coor_local = esm.physical_crd_xi(self.parent_nurbs_surf,\
xi1, xi2, self.node_1_u, self.node_1_v,\
self.node_3_u, self.node_3_v)
nodes_coor_mtx[i, j, 0] = physical_coor_local[0]
nodes_coor_mtx[i, j, 1] = physical_coor_local[1]
nodes_coor_mtx[i, j, 2] = physical_coor_local[2]
return nodes_coor_mtx
# def coorsys_director_tanvec(self, nodes_physical_coor, \
# lag_xi1, lag_xi2, der_lag_dxi1, der_lag_dxi2):
# '''
# In this method, the nodal coordinate system at each node
# specified by the Lagrangina functions at that point is calculated.
# -Output:
# A 5x3 np.array cotains the coordinate system and tangent vectors.
# '''
# dim = self.lobatto_pw.shape[0]
# g1 = np.zeros(3)
# g2 = np.zeros(3)
# # lag_xi1 = esm.lagfunc(self.lobatto_pw, xi1)
# # lag_xi2 = esm.lagfunc(self.lobatto_pw, xi2)
# # der_lag_dxi1 = esm.der_lagfunc_dxi(self.lobatto_pw, xi1)
# # der_lag_dxi2 = esm.der_lagfunc_dxi(self.lobatto_pw, xi2)
# for i in range(dim):
# for j in range(dim):
# g1 = g1 + np.array([der_lag_dxi1[j] * lag_xi2[i] * nodes_physical_coor[i, j, 0],\
# der_lag_dxi1[j] * lag_xi2[i] * nodes_physical_coor[i, j, 1],\
# der_lag_dxi1[j] * lag_xi2[i] * nodes_physical_coor[i, j, 2]])
# g2 = g2 + np.array([lag_xi1[j] * der_lag_dxi2[i] * nodes_physical_coor[i, j, 0],\
# lag_xi1[j] * der_lag_dxi2[i] * nodes_physical_coor[i, j, 1],\
# lag_xi1[j] * der_lag_dxi2[i] * nodes_physical_coor[i, j, 2]])
# vdir_not_unit = np.cross(g1, g2)
# vdir_unit = vdir_not_unit / np.linalg.norm(vdir_not_unit)
# if abs(abs(vdir_unit[1])-1) > 0.00001 :
# v1 = (np.cross(ey, vdir_unit))
# v1_unit= v1 / np.linalg.norm(v1)
# v2_unit = np.cross(vdir_unit, v1_unit)
# else:
# if vdir_unit[1] > 0:
# v1_unit = np.array((1, 0, 0))
# v2_unit = np.array((0, 0, -1))
# vdir_unit = np.array((0, 1, 0))
# else:
# v1_unit = np.array((1, 0, 0))
# v2_unit = np.array((0, 0, 1))
# vdir_unit = np.array((0, 1, 0))
# return np.array([v1_unit, v2_unit, vdir_unit, g1, g2])
def coorsys_director_tanvec_allnodes(self):
'''
In this method, the nodal coordinate system at each node
specified by the Lagrangina functions at that point is calculated.
-Output:
A nxnx5x3 np.array cotains the coordinate system and tangent vectors.
'''
nodes_physical_coor = self.nodes_physical_coordinate()
dim = self.lobatto_pw.shape[0]
coor_tan_mtx = np.zeros((dim, dim, 5, 3))
# lag_xi1 = esm.lagfunc(self.lobatto_pw, xi1)
# lag_xi2 = esm.lagfunc(self.lobatto_pw, xi2)
# der_lag_dxi1 = esm.der_lagfunc_dxi(self.lobatto_pw, xi1)
# der_lag_dxi2 = esm.der_lagfunc_dxi(self.lobatto_pw, xi2)
for i_main in range(dim):
xi2 = self.lobatto_pw[i_main, 0]
lag_xi2 = esm.lagfunc(self.lobatto_pw, xi2)
der_lag_dxi2 = esm.der_lagfunc_dxi(self.lobatto_pw, xi2)
for j_main in range(dim):
xi1 = self.lobatto_pw[j_main, 0]
lag_xi1 = esm.lagfunc(self.lobatto_pw, xi1)
der_lag_dxi1 = esm.der_lagfunc_dxi(self.lobatto_pw, xi1)
g1 = np.zeros(3)
g2 = np.zeros(3)
for i in range(dim):
for j in range(dim):
g1 = g1 + np.array([der_lag_dxi1[j] * lag_xi2[i] * nodes_physical_coor[i, j, 0],\
der_lag_dxi1[j] * lag_xi2[i] * nodes_physical_coor[i, j, 1],\
der_lag_dxi1[j] * lag_xi2[i] * nodes_physical_coor[i, j, 2]])
g2 = g2 + np.array([lag_xi1[j] * der_lag_dxi2[i] * nodes_physical_coor[i, j, 0],\
lag_xi1[j] * der_lag_dxi2[i] * nodes_physical_coor[i, j, 1],\
lag_xi1[j] * der_lag_dxi2[i] * nodes_physical_coor[i, j, 2]])
vdir_not_unit = np.cross(g1, g2)
vdir_unit = vdir_not_unit / np.linalg.norm(vdir_not_unit)
if abs(abs(vdir_unit[1])-1) > 0.00001 :
v1 = (np.cross(ey, vdir_unit))
v1_unit= v1 / np.linalg.norm(v1)
v2_unit = np.cross(vdir_unit, v1_unit)
else:
if vdir_unit[1] > 0:
v1_unit = np.array((1, 0, 0))
v2_unit = np.array((0, 0, -1))
vdir_unit = np.array((0, 1, 0))
else:
v1_unit = np.array((1, 0, 0))
v2_unit = np.array((0, 0, 1))
vdir_unit = np.array((0, 1, 0))
coor_tan_mtx[i_main, j_main, 0] = v1_unit
coor_tan_mtx[i_main, j_main, 1] = v2_unit
coor_tan_mtx[i_main, j_main, 2] = vdir_unit
coor_tan_mtx[i_main, j_main, 3] = g1
coor_tan_mtx[i_main, j_main, 4] = g2
return coor_tan_mtx
def jacobian_mtx(self, coorsys_tanvec_mtx, row_num, col_num,\
xi3, lag_xi1, lag_xi2, der_lag_dxi1, der_lag_dxi2):
'''
In this method, the Jacobian matrix is calculated. Unlike our developed
method, in isoparametric elements, we need the director vectors at
all nodes to iterpolate it, which is necessary for calculation of the
derivatives of the director at a specific point. This is the reason that
the coorsys_tanvec_matx which is calculated from
"coorsys_director_tanvec_allnodes" method is imported.
-Output:
A 3x3 np.array which is Jacobian matrix.
'''
vdir_unit = coorsys_tanvec_mtx[row_num, col_num, 2]
g1 = coorsys_tanvec_mtx[row_num, col_num, 3]
g2 = coorsys_tanvec_mtx[row_num, col_num, 4]
dim = self.lobatto_pw.shape[0]
dvdir_unit_dxi1 = np.zeros(3)
dvdir_unit_dxi2 = np.zeros(3)
# with open('director_unit_isoprm.dat', 'a') as du:
# np.savetxt(du, coorsys_tanvec_mtx[row_num, col_num, 2])
# # du.write('\n')
for i in range(dim):
for j in range(dim):
dvdir_unit_dxi1 = dvdir_unit_dxi1 + \
np.array([der_lag_dxi1[j] * lag_xi2[i] * coorsys_tanvec_mtx[i, j, 2, 0],\
der_lag_dxi1[j] * lag_xi2[i] * coorsys_tanvec_mtx[i, j, 2, 1],\
der_lag_dxi1[j] * lag_xi2[i] * coorsys_tanvec_mtx[i, j, 2, 2]])
dvdir_unit_dxi2 = dvdir_unit_dxi2 +\
np.array([lag_xi1[j] * der_lag_dxi2[i] * coorsys_tanvec_mtx[i, j, 2, 0],\
lag_xi1[j] * der_lag_dxi2[i] * coorsys_tanvec_mtx[i, j, 2, 1],\
lag_xi1[j] * der_lag_dxi2[i] * coorsys_tanvec_mtx[i, j, 2, 2]])
dx_dxi1 = g1[0] + xi3/2 * self.thk * dvdir_unit_dxi1[0]
dy_dxi1 = g1[1] + xi3/2 * self.thk * dvdir_unit_dxi1[1]
dz_dxi1 = g1[2] + xi3/2 * self.thk * dvdir_unit_dxi1[2]
dx_dxi2 = g2[0] + xi3/2 * self.thk * dvdir_unit_dxi2[0]
dy_dxi2 = g2[1] + xi3/2 * self.thk * dvdir_unit_dxi2[1]
dz_dxi2 = g2[2] + xi3/2 * self.thk * dvdir_unit_dxi2[2]
dx_dxi3 = 1/2 * self.thk * vdir_unit[0]
dy_dxi3 = 1/2 * self.thk * vdir_unit[1]
dz_dxi3 = 1/2 * self.thk * vdir_unit[2]
# dx_dxi1 = dvdir_unit_dxi1[0]
# dy_dxi1 = dvdir_unit_dxi1[1]
# dz_dxi1 = dvdir_unit_dxi1[2]
# dx_dxi2 = dvdir_unit_dxi2[0]
# dy_dxi2 = dvdir_unit_dxi2[1]
# dz_dxi2 = dvdir_unit_dxi2[2]
# dx_dxi3 = vdir_unit[0]
# dy_dxi3 = vdir_unit[1]
# dz_dxi3 = vdir_unit[2]
jac = np.array(((dx_dxi1, dy_dxi1, dz_dxi1),(dx_dxi2, dy_dxi2, dz_dxi2),\
(dx_dxi3, dy_dxi3, dz_dxi3)))
return jac
def area(self):
'''This function is only for the test of the
validity of derivatives of mapping and jacobian'''
node_coor = self.nodes_physical_coordinate()
coorsys_tanvec_mtx = self.coorsys_director_tanvec_allnodes()
dim = self.lobatto_pw.shape[0]
local_area = 0
for i in range(dim):
xi2 = self.lobatto_pw[i, 0]
w2 = self.lobatto_pw[i, 1]
lag_xi2 = esm.lagfunc(self.lobatto_pw, xi2)
der_lag_dxi2 = esm.der_lagfunc_dxi(self.lobatto_pw, xi2)
for j in range(dim):
xi1 = self.lobatto_pw[j, 0]
w1 = self.lobatto_pw[j, 1]
lag_xi1 = esm.lagfunc(self.lobatto_pw, xi1)
der_lag_dxi1 = esm.der_lagfunc_dxi(self.lobatto_pw, xi1)
jac_mtx = self.jacobian_mtx( coorsys_tanvec_mtx, i, j, 0, lag_xi1, lag_xi2, der_lag_dxi1, der_lag_dxi2)
vector_1 = jac_mtx[0,:]
vector_2 = jac_mtx[1,:]
local_area = local_area +\
np.linalg.norm((np.cross(vector_1,vector_2)))* w1 * w2
return local_area
def volume(self, number_gauss_point):
'''This function is only for the test of the
validity of derivatives of mapping and jacobian'''
node_coor = self.nodes_physical_coordinate()
coorsys_tanvec_mtx = self.coorsys_director_tanvec_allnodes()
dim = self.lobatto_pw.shape[0]
gauss_pw_nparray = np.array(esm.gauss_pw[number_gauss_point-2])
dim_gauss = gauss_pw_nparray.shape[0]
local_volume = 0
for i in range(dim):
xi2 = self.lobatto_pw[i, 0]
w2 = lobatto_pw[i, 1]
lag_xi2 = esm.lagfunc(self.lobatto_pw, xi2)
der_lag_dxi2 = esm.der_lagfunc_dxi(self.lobatto_pw, xi2)
for j in range(dim):
xi1 = self.lobatto_pw[j, 0]
w1 = lobatto_pw[j, 1]
lag_xi1 = esm.lagfunc(self.lobatto_pw, xi1)
for k in range(dim_gauss):
xi3 = gauss_pw_nparray[k, 0]
w3 =gauss_pw_nparray[k, 1]
der_lag_dxi1 = esm.der_lagfunc_dxi(self.lobatto_pw, xi1)
jac_mtx = self.jacobian_mtx( coorsys_tanvec_mtx, i, j, xi3,\
lag_xi1, lag_xi2, der_lag_dxi1, der_lag_dxi2)
local_volume = local_volume +\
np.linalg.det(jac_mtx)* w1 * w2 * w3
return local_volume
if __name__=="__main__":
thk = 0.5
data = exchange.import_json("double-curved-free-form_5B.json") #curved_beam_lineload_2.json double-curved-free-form.json
nurbs_surface = sgs.SurfaceGeo(data, 0, thk)
lobatto_pw_all = esm.lbto_pw("node_weight_all.dat")
i_main = 9
if i_main == 1:
lobatto_pw = lobatto_pw_all[1:3,:]
else:
index = np.argwhere(lobatto_pw_all==i_main)
lobatto_pw = lobatto_pw_all[index[0, 0] + 1:\
index[0, 0] + (i_main+1) +1, :]
node_1_ub = 0
node_1_vb = 0
node_3_ub = 1
node_3_vb = 1
surf_lag = SurfaceGeneratedSEM(nurbs_surface, lobatto_pw, node_1_ub,\
node_1_vb, node_3_ub, node_3_vb)
coorsys_tanvec_mtx = surf_lag.coorsys_director_tanvec_allnodes()
dim = lobatto_pw.shape[0]
xi3 = 0
with open("jacobian_isoprm.dat", "w") as f:
pass
for i in range(dim):
xi2 = surf_lag.lobatto_pw[i, 0]
w2 = surf_lag.lobatto_pw[i, 1]
lag_xi2 = esm.lagfunc(surf_lag.lobatto_pw, xi2)
der_lag_dxi2 = esm.der_lagfunc_dxi(surf_lag.lobatto_pw, xi2)
for j in range(dim):
xi1 = surf_lag.lobatto_pw[j, 0]
w1 = surf_lag.lobatto_pw[j, 1]
lag_xi1 = esm.lagfunc(surf_lag.lobatto_pw, xi1)
der_lag_dxi1 = esm.der_lagfunc_dxi(surf_lag.lobatto_pw, xi1)
jac_mtx = surf_lag.jacobian_mtx( coorsys_tanvec_mtx, i, j, xi3,\
lag_xi1, lag_xi2, der_lag_dxi1, der_lag_dxi2)
with open ('jacobian_isoprm.dat', 'a') as jac_file:
np.savetxt(jac_file, jac_mtx)
jac_file.write("\n")
print("End")
# print(surf_lag.thk)
# print(surf_lag.nodes_physical_coordinate())
# node_coor = surf_lag.nodes_physical_coordinate()
# print(surf_lag.coorsys_director_tanvec_allnodes())
# print(surf_lag.area())
# print(surf_lag.volume(2))
# print('End')